00001 #ifndef RIVET_MATH_VECTOR4
00002 #define RIVET_MATH_VECTOR4
00003
00004 #include "Rivet/Math/MathHeader.hh"
00005 #include "Rivet/Math/MathUtils.hh"
00006 #include "Rivet/Math/VectorN.hh"
00007 #include "Rivet/Math/Vector3.hh"
00008
00009 namespace Rivet {
00010
00011
00012 class FourVector;
00013 class FourMomentum;
00014 class LorentzTransform;
00015 typedef FourVector Vector4;
00016 FourVector transform(const LorentzTransform& lt, const FourVector& v4);
00017
00018
00019
00020 class FourVector : public Vector<4> {
00021 friend FourVector multiply(const double a, const FourVector& v);
00022 friend FourVector multiply(const FourVector& v, const double a);
00023 friend FourVector add(const FourVector& a, const FourVector& b);
00024 friend FourVector transform(const LorentzTransform& lt, const FourVector& v4);
00025
00026 public:
00027 FourVector() : Vector<4>() { }
00028
00029 template<typename V4>
00030 FourVector(const V4& other) {
00031 this->setT(other.t());
00032 this->setX(other.x());
00033 this->setY(other.y());
00034 this->setZ(other.z());
00035 }
00036
00037 FourVector(const Vector<4>& other)
00038 : Vector<4>(other) { }
00039
00040 FourVector(const double t, const double x, const double y, const double z) {
00041 this->setT(t);
00042 this->setX(x);
00043 this->setY(y);
00044 this->setZ(z);
00045 }
00046
00047 virtual ~FourVector() { }
00048
00049 public:
00050 double t() const { return get(0); }
00051 double x() const { return get(1); }
00052 double y() const { return get(2); }
00053 double z() const { return get(3); }
00054 FourVector& setT(const double t) { set(0, t); return *this; }
00055 FourVector& setX(const double x) { set(1, x); return *this; }
00056 FourVector& setY(const double y) { set(2, y); return *this; }
00057 FourVector& setZ(const double z) { set(3, z); return *this; }
00058
00059 double invariant() const {
00060
00061 return (t() + z())*(t() - z()) - x()*x() - y()*y();
00062 }
00063
00064 double angle(const FourVector& v) const {
00065 return vector3().angle( v.vector3() );
00066 }
00067
00068 double angle(const Vector3& v3) const {
00069 return vector3().angle(v3);
00070 }
00071
00072 double polarRadius2() const {
00073 return vector3().polarRadius2();
00074 }
00075
00076
00077 double perp2() const {
00078 return vector3().perp2();
00079 }
00080
00081
00082 double rho2() const {
00083 return vector3().rho2();
00084 }
00085
00086 double polarRadius() const {
00087 return vector3().polarRadius();
00088 }
00089
00090
00091 double perp() const {
00092 return vector3().perp();
00093 }
00094
00095
00096 double rho() const {
00097 return vector3().rho();
00098 }
00099
00100
00101 double azimuthalAngle(const PhiMapping mapping = ZERO_2PI) const {
00102 return vector3().azimuthalAngle(mapping);
00103 }
00104
00105
00106 double phi(const PhiMapping mapping = ZERO_2PI) const {
00107 return vector3().phi(mapping);
00108 }
00109
00110
00111 double polarAngle() const {
00112 return vector3().polarAngle();
00113 }
00114
00115
00116 double theta() const {
00117 return vector3().theta();
00118 }
00119
00120 double pseudorapidity() const {
00121 return vector3().pseudorapidity();
00122 }
00123
00124
00125 double eta() const {
00126 return vector3().eta();
00127 }
00128
00129
00130 Vector3 vector3() const {
00131 return Vector3(get(1), get(2), get(3));
00132 }
00133
00134 public:
00135
00136 double contract(const FourVector& v) const {
00137 const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
00138 return result;
00139 }
00140
00141
00142 double dot(const FourVector& v) const {
00143 return contract(v);
00144 }
00145
00146
00147 double operator*(const FourVector& v) const {
00148 return contract(v);
00149 }
00150
00151
00152 FourVector& operator*=(double a) {
00153 _vec = multiply(a, *this)._vec;
00154 return *this;
00155 }
00156
00157
00158 FourVector& operator/=(double a) {
00159 _vec = multiply(1.0/a, *this)._vec;
00160 return *this;
00161 }
00162
00163 FourVector& operator+=(const FourVector& v) {
00164 _vec = add(*this, v)._vec;
00165 return *this;
00166 }
00167
00168 FourVector& operator-=(const FourVector& v) {
00169 _vec = add(*this, -v)._vec;
00170 return *this;
00171 }
00172
00173 FourVector operator-() const {
00174 FourVector result;
00175 result._vec = -_vec;
00176 return result;
00177 }
00178
00179 };
00180
00181
00182
00183 inline double contract(const FourVector& a, const FourVector& b) {
00184 return a.contract(b);
00185 }
00186
00187
00188 inline double dot(const FourVector& a, const FourVector& b) {
00189 return contract(a, b);
00190 }
00191
00192 inline FourVector multiply(const double a, const FourVector& v) {
00193 FourVector result;
00194 result._vec = a * v._vec;
00195 return result;
00196 }
00197
00198 inline FourVector multiply(const FourVector& v, const double a) {
00199 return multiply(a, v);
00200 }
00201
00202 inline FourVector operator*(const double a, const FourVector& v) {
00203 return multiply(a, v);
00204 }
00205
00206 inline FourVector operator*(const FourVector& v, const double a) {
00207 return multiply(a, v);
00208 }
00209
00210 inline FourVector operator/(const FourVector& v, const double a) {
00211 return multiply(1.0/a, v);
00212 }
00213
00214 inline FourVector add(const FourVector& a, const FourVector& b) {
00215 FourVector result;
00216 result._vec = a._vec + b._vec;
00217 return result;
00218 }
00219
00220 inline FourVector operator+(const FourVector& a, const FourVector& b) {
00221 return add(a, b);
00222 }
00223
00224 inline FourVector operator-(const FourVector& a, const FourVector& b) {
00225 return add(a, -b);
00226 }
00227
00228
00229
00230 inline double invariant(const FourVector& lv) {
00231 return lv.invariant();
00232 }
00233
00234
00235 inline double angle(const FourVector& a, const FourVector& b) {
00236 return a.angle(b);
00237 }
00238
00239
00240 inline double angle(const Vector3& a, const FourVector& b) {
00241 return angle( a, b.vector3() );
00242 }
00243
00244
00245 inline double angle(const FourVector& a, const Vector3& b) {
00246 return a.angle(b);
00247 }
00248
00249
00250 inline double polarRadius2(const FourVector& v) {
00251 return v.polarRadius2();
00252 }
00253
00254 inline double perp2(const FourVector& v) {
00255 return v.perp2();
00256 }
00257
00258 inline double rho2(const FourVector& v) {
00259 return v.rho2();
00260 }
00261
00262
00263 inline double polarRadius(const FourVector& v) {
00264 return v.polarRadius();
00265 }
00266
00267 inline double perp(const FourVector& v) {
00268 return v.perp();
00269 }
00270
00271 inline double rho(const FourVector& v) {
00272 return v.rho();
00273 }
00274
00275
00276 inline double azimuthalAngle(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00277 return v.azimuthalAngle(mapping);
00278 }
00279
00280 inline double phi(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00281 return v.phi(mapping);
00282 }
00283
00284
00285
00286 inline double polarAngle(const FourVector& v) {
00287 return v.polarAngle();
00288 }
00289
00290 inline double theta(const FourVector& v) {
00291 return v.theta();
00292 }
00293
00294
00295 inline double pseudorapidity(const FourVector& v) {
00296 return v.pseudorapidity();
00297 }
00298
00299 inline double eta(const FourVector& v) {
00300 return v.eta();
00301 }
00302
00303
00304
00305
00306
00307
00308
00309
00310 class FourMomentum : public FourVector {
00311 public:
00312 FourMomentum() { }
00313
00314 template<typename V4>
00315 FourMomentum(const V4& other) {
00316 this->setE(other.t());
00317 this->setPx(other.x());
00318 this->setPy(other.y());
00319 this->setPz(other.z());
00320 }
00321
00322 FourMomentum(const Vector<4>& other)
00323 : FourVector(other) { }
00324
00325 FourMomentum(const double E, const double px, const double py, const double pz) {
00326 this->setE(E);
00327 this->setPx(px);
00328 this->setPy(py);
00329 this->setPz(pz);
00330 }
00331
00332 ~FourMomentum() {}
00333
00334 public:
00335
00336 double E() const { return t(); }
00337
00338
00339 Vector3 p() const { return vector3(); }
00340
00341
00342 double px() const { return x(); }
00343
00344
00345 double py() const { return y(); }
00346
00347
00348 double pz() const { return z(); }
00349
00350
00351 FourMomentum& setE(double E) { setT(E); return *this; }
00352
00353
00354 FourMomentum& setPx(double px) { setX(px); return *this; }
00355
00356
00357 FourMomentum& setPy(double py) { setY(py); return *this; }
00358
00359
00360 FourMomentum& setPz(double pz) { setZ(pz); return *this; }
00361
00362
00363 double mass2() const {
00364 return invariant();
00365 }
00366
00367
00368 double mass() const {
00369 assert(Rivet::isZero(mass2()) || mass2() > 0);
00370 return sqrt(mass2());
00371 }
00372
00373
00374 double rapidity() const {
00375 return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
00376 }
00377
00378
00379 double pT2() const {
00380 return vector3().polarRadius2();
00381 }
00382
00383
00384 double pT() const {
00385 return sqrt(pT2());
00386 }
00387
00388
00389 double Et2() const {
00390 return Et() * Et();
00391 }
00392
00393
00394 double Et() const {
00395 return E() * sin(polarAngle());
00396 }
00397
00398
00399 Vector3 boostVector() const {
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411 return Vector3(px()/E(), py()/E(), pz()/E());
00412 }
00413
00414
00415
00416 struct byEAscending{
00417 bool operator()(const FourMomentum &left, const FourMomentum &right) const{
00418 double pt2left = left.E();
00419 double pt2right = right.E();
00420 return pt2left < pt2right;
00421 }
00422
00423 bool operator()(const FourMomentum *left, const FourMomentum *right) const{
00424 return (*this)(left, right);
00425 }
00426 };
00427
00428
00429
00430 struct byEDescending{
00431 bool operator()(const FourMomentum &left, const FourMomentum &right) const{
00432 return byEAscending()(right, left);
00433 }
00434
00435 bool operator()(const FourMomentum *left, const FourVector *right) const{
00436 return (*this)(left, right);
00437 }
00438 };
00439
00440 };
00441
00442
00443
00444 inline double mass2(const FourMomentum& v) {
00445 return v.mass2();
00446 }
00447
00448
00449 inline double mass(const FourMomentum& v) {
00450 return v.mass();
00451 }
00452
00453
00454 inline double rapidity(const FourMomentum& v) {
00455 return v.rapidity();
00456 }
00457
00458
00459 inline double pT2(const FourMomentum& v) {
00460 return v.pT2();
00461 }
00462
00463
00464 inline double pT(const FourMomentum& v) {
00465 return v.pT();
00466 }
00467
00468
00469 inline double Et2(const FourMomentum& v) {
00470 return v.Et2();}
00471
00472
00473 inline double Et(const FourMomentum& v) {
00474 return v.Et();
00475 }
00476
00477
00478 inline Vector3 boostVector(const FourMomentum& v) {
00479 return v.boostVector();
00480 }
00481
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493 inline double deltaR(const FourVector& a, const FourVector& b,
00494 DeltaRScheme scheme = PSEUDORAPIDITY) {
00495 switch (scheme) {
00496 case PSEUDORAPIDITY :
00497 return deltaR(a.vector3(), b.vector3());
00498 case RAPIDITY:
00499 {
00500 const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
00501 const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
00502 if (!ma || !mb) {
00503 string err = "deltaR with scheme RAPIDITY, can be called with FourMomenta only";
00504 throw std::runtime_error(err);
00505 }
00506 return deltaR(*ma, *mb, scheme);
00507 }
00508 default:
00509 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00510 }
00511 }
00512
00513
00514 inline double deltaR(const FourVector& v,
00515 double eta2, double phi2,
00516 DeltaRScheme scheme = PSEUDORAPIDITY) {
00517 switch (scheme) {
00518 case PSEUDORAPIDITY :
00519 return deltaR(v.vector3(), eta2, phi2);
00520 case RAPIDITY:
00521 {
00522 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00523 if (!mv) {
00524 string err = "deltaR with scheme RAPIDITY, can be called with FourMomenta only";
00525 throw std::runtime_error(err);
00526 }
00527 return deltaR(*mv, eta2, phi2, scheme);
00528 }
00529 default:
00530 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00531 }
00532 }
00533
00534
00535 inline double deltaR(double eta1, double phi1,
00536 const FourVector& v,
00537 DeltaRScheme scheme = PSEUDORAPIDITY) {
00538 switch (scheme) {
00539 case PSEUDORAPIDITY :
00540 return deltaR(eta1, phi1, v.vector3());
00541 case RAPIDITY:
00542 {
00543 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00544 if (!mv) {
00545 string err = "deltaR with scheme RAPIDITY, can be called with FourMomenta only";
00546 throw std::runtime_error(err);
00547 }
00548 return deltaR(eta1, phi1, *mv, scheme);
00549 }
00550 default:
00551 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00552 }
00553 }
00554
00555
00556
00557
00558
00559
00560
00561 inline double deltaR(const FourMomentum& a, const FourMomentum& b,
00562 DeltaRScheme scheme = PSEUDORAPIDITY) {
00563 switch (scheme) {
00564 case PSEUDORAPIDITY:
00565 return deltaR(a.vector3(), b.vector3());
00566 case RAPIDITY:
00567 return deltaR(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
00568 default:
00569 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00570 }
00571 }
00572
00573 inline double deltaR(const FourMomentum& v,
00574 double eta2, double phi2,
00575 DeltaRScheme scheme = PSEUDORAPIDITY) {
00576 switch (scheme) {
00577 case PSEUDORAPIDITY:
00578 return deltaR(v.vector3(), eta2, phi2);
00579 case RAPIDITY:
00580 return deltaR(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
00581 default:
00582 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00583 }
00584 }
00585
00586
00587 inline double deltaR(double eta1, double phi1,
00588 const FourMomentum& v,
00589 DeltaRScheme scheme = PSEUDORAPIDITY) {
00590 switch (scheme) {
00591 case PSEUDORAPIDITY:
00592 return deltaR(eta1, phi1, v.vector3());
00593 case RAPIDITY:
00594 return deltaR(eta1, phi1, v.rapidity(), v.azimuthalAngle());
00595 default:
00596 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00597 }
00598 }
00599
00600
00601
00602
00603
00604
00605 inline const string toString(const FourVector& lv) {
00606 ostringstream out;
00607 out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
00608 << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
00609 << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
00610 << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
00611 << ")";
00612 return out.str();
00613 }
00614
00615
00616 inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
00617 out << toString(lv);
00618 return out;
00619 }
00620
00621
00622 }
00623
00624 #endif