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
00028 FourVector() : Vector<4>() { }
00029
00030 template<typename V4>
00031 FourVector(const V4& other) {
00032 this->setT(other.t());
00033 this->setX(other.x());
00034 this->setY(other.y());
00035 this->setZ(other.z());
00036 }
00037
00038 FourVector(const Vector<4>& other)
00039 : Vector<4>(other) { }
00040
00041 FourVector(const double t, const double x, const double y, const double z) {
00042 this->setT(t);
00043 this->setX(x);
00044 this->setY(y);
00045 this->setZ(z);
00046 }
00047
00048 virtual ~FourVector() { }
00049
00050 public:
00051
00052 double t() const { return get(0); }
00053 double x() const { return get(1); }
00054 double y() const { return get(2); }
00055 double z() const { return get(3); }
00056 FourVector& setT(const double t) { set(0, t); return *this; }
00057 FourVector& setX(const double x) { set(1, x); return *this; }
00058 FourVector& setY(const double y) { set(2, y); return *this; }
00059 FourVector& setZ(const double z) { set(3, z); return *this; }
00060
00061 double invariant() const {
00062
00063 return (t() + z())*(t() - z()) - x()*x() - y()*y();
00064 }
00065
00066
00067 double angle(const FourVector& v) const {
00068 return vector3().angle( v.vector3() );
00069 }
00070
00071
00072 double angle(const Vector3& v3) const {
00073 return vector3().angle(v3);
00074 }
00075
00076
00077
00078
00079 double polarRadius2() const {
00080 return vector3().polarRadius2();
00081 }
00082
00083
00084 double perp2() const {
00085 return vector3().perp2();
00086 }
00087
00088
00089 double rho2() const {
00090 return vector3().rho2();
00091 }
00092
00093
00094 double polarRadius() const {
00095 return vector3().polarRadius();
00096 }
00097
00098
00099 double perp() const {
00100 return vector3().perp();
00101 }
00102
00103
00104 double rho() const {
00105 return vector3().rho();
00106 }
00107
00108
00109 double azimuthalAngle(const PhiMapping mapping = ZERO_2PI) const {
00110 return vector3().azimuthalAngle(mapping);
00111 }
00112
00113
00114 double phi(const PhiMapping mapping = ZERO_2PI) const {
00115 return vector3().phi(mapping);
00116 }
00117
00118
00119 double polarAngle() const {
00120 return vector3().polarAngle();
00121 }
00122
00123
00124 double theta() const {
00125 return vector3().theta();
00126 }
00127
00128
00129 double pseudorapidity() const {
00130 return vector3().pseudorapidity();
00131 }
00132
00133
00134 double eta() const {
00135 return vector3().eta();
00136 }
00137
00138
00139 Vector3 vector3() const {
00140 return Vector3(get(1), get(2), get(3));
00141 }
00142
00143
00144 public:
00145
00146
00147 double contract(const FourVector& v) const {
00148 const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
00149 return result;
00150 }
00151
00152
00153 double dot(const FourVector& v) const {
00154 return contract(v);
00155 }
00156
00157
00158 double operator*(const FourVector& v) const {
00159 return contract(v);
00160 }
00161
00162
00163 FourVector& operator*=(double a) {
00164 _vec = multiply(a, *this)._vec;
00165 return *this;
00166 }
00167
00168
00169 FourVector& operator/=(double a) {
00170 _vec = multiply(1.0/a, *this)._vec;
00171 return *this;
00172 }
00173
00174 FourVector& operator+=(const FourVector& v) {
00175 _vec = add(*this, v)._vec;
00176 return *this;
00177 }
00178
00179 FourVector& operator-=(const FourVector& v) {
00180 _vec = add(*this, -v)._vec;
00181 return *this;
00182 }
00183
00184 FourVector operator-() const {
00185 FourVector result;
00186 result._vec = -_vec;
00187 return result;
00188 }
00189
00190 };
00191
00192
00193
00194 inline double contract(const FourVector& a, const FourVector& b) {
00195 return a.contract(b);
00196 }
00197
00198
00199 inline double dot(const FourVector& a, const FourVector& b) {
00200 return contract(a, b);
00201 }
00202
00203 inline FourVector multiply(const double a, const FourVector& v) {
00204 FourVector result;
00205 result._vec = a * v._vec;
00206 return result;
00207 }
00208
00209 inline FourVector multiply(const FourVector& v, const double a) {
00210 return multiply(a, v);
00211 }
00212
00213 inline FourVector operator*(const double a, const FourVector& v) {
00214 return multiply(a, v);
00215 }
00216
00217 inline FourVector operator*(const FourVector& v, const double a) {
00218 return multiply(a, v);
00219 }
00220
00221 inline FourVector operator/(const FourVector& v, const double a) {
00222 return multiply(1.0/a, v);
00223 }
00224
00225 inline FourVector add(const FourVector& a, const FourVector& b) {
00226 FourVector result;
00227 result._vec = a._vec + b._vec;
00228 return result;
00229 }
00230
00231 inline FourVector operator+(const FourVector& a, const FourVector& b) {
00232 return add(a, b);
00233 }
00234
00235 inline FourVector operator-(const FourVector& a, const FourVector& b) {
00236 return add(a, -b);
00237 }
00238
00239
00240
00241 inline double invariant(const FourVector& lv) {
00242 return lv.invariant();
00243 }
00244
00245
00246 inline double angle(const FourVector& a, const FourVector& b) {
00247 return a.angle(b);
00248 }
00249
00250
00251 inline double angle(const Vector3& a, const FourVector& b) {
00252 return angle( a, b.vector3() );
00253 }
00254
00255
00256 inline double angle(const FourVector& a, const Vector3& b) {
00257 return a.angle(b);
00258 }
00259
00260
00261 inline double polarRadius2(const FourVector& v) {
00262 return v.polarRadius2();
00263 }
00264
00265 inline double perp2(const FourVector& v) {
00266 return v.perp2();
00267 }
00268
00269 inline double rho2(const FourVector& v) {
00270 return v.rho2();
00271 }
00272
00273
00274 inline double polarRadius(const FourVector& v) {
00275 return v.polarRadius();
00276 }
00277
00278 inline double perp(const FourVector& v) {
00279 return v.perp();
00280 }
00281
00282 inline double rho(const FourVector& v) {
00283 return v.rho();
00284 }
00285
00286
00287 inline double azimuthalAngle(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00288 return v.azimuthalAngle(mapping);
00289 }
00290
00291 inline double phi(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00292 return v.phi(mapping);
00293 }
00294
00295
00296
00297 inline double polarAngle(const FourVector& v) {
00298 return v.polarAngle();
00299 }
00300
00301 inline double theta(const FourVector& v) {
00302 return v.theta();
00303 }
00304
00305
00306 inline double pseudorapidity(const FourVector& v) {
00307 return v.pseudorapidity();
00308 }
00309
00310 inline double eta(const FourVector& v) {
00311 return v.eta();
00312 }
00313
00314
00315
00316
00317
00318
00319
00320
00321 class FourMomentum : public FourVector {
00322 friend FourMomentum multiply(const double a, const FourMomentum& v);
00323 friend FourMomentum multiply(const FourMomentum& v, const double a);
00324 friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
00325 friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
00326
00327 public:
00328 FourMomentum() { }
00329
00330 template<typename V4>
00331 FourMomentum(const V4& other) {
00332 this->setE(other.t());
00333 this->setPx(other.x());
00334 this->setPy(other.y());
00335 this->setPz(other.z());
00336 }
00337
00338 FourMomentum(const Vector<4>& other)
00339 : FourVector(other) { }
00340
00341 FourMomentum(const double E, const double px, const double py, const double pz) {
00342 this->setE(E);
00343 this->setPx(px);
00344 this->setPy(py);
00345 this->setPz(pz);
00346 }
00347
00348 ~FourMomentum() {}
00349
00350 public:
00351
00352 double E() const { return t(); }
00353
00354
00355 Vector3 p() const { return vector3(); }
00356
00357
00358 double px() const { return x(); }
00359
00360
00361 double py() const { return y(); }
00362
00363
00364 double pz() const { return z(); }
00365
00366
00367 FourMomentum& setE(double E) { setT(E); return *this; }
00368
00369
00370 FourMomentum& setPx(double px) { setX(px); return *this; }
00371
00372
00373 FourMomentum& setPy(double py) { setY(py); return *this; }
00374
00375
00376 FourMomentum& setPz(double pz) { setZ(pz); return *this; }
00377
00378
00379 double mass() const {
00380 assert(Rivet::isZero(mass2()) || mass2() > 0);
00381 if (Rivet::isZero(mass2())) {
00382 return 0.0;
00383 } else {
00384 return sqrt(mass2());
00385 }
00386 }
00387
00388
00389 double massT() const {
00390 return mass() * sin(polarAngle());
00391 }
00392
00393
00394 double mass2() const {
00395 return invariant();
00396 }
00397
00398
00399 double massT2() const {
00400 return massT() * massT();
00401 }
00402
00403
00404 double rapidity() const {
00405 return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
00406 }
00407
00408
00409 double pT2() const {
00410 return vector3().polarRadius2();
00411 }
00412
00413
00414 double pT() const {
00415 return sqrt(pT2());
00416 }
00417
00418
00419 double Et2() const {
00420 return Et() * Et();
00421 }
00422
00423
00424 double Et() const {
00425 return E() * sin(polarAngle());
00426 }
00427
00428
00429 Vector3 boostVector() const {
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441 return Vector3(px()/E(), py()/E(), pz()/E());
00442 }
00443
00444
00445 struct byEAscending {
00446 bool operator()(const FourMomentum& left, const FourMomentum& right) const{
00447 double pt2left = left.E();
00448 double pt2right = right.E();
00449 return pt2left < pt2right;
00450 }
00451
00452 bool operator()(const FourMomentum* left, const FourMomentum* right) const{
00453 return (*this)(left, right);
00454 }
00455 };
00456
00457
00458 struct byEDescending {
00459 bool operator()(const FourMomentum& left, const FourMomentum& right) const{
00460 return byEAscending()(right, left);
00461 }
00462
00463 bool operator()(const FourMomentum* left, const FourVector* right) const{
00464 return (*this)(left, right);
00465 }
00466 };
00467
00468
00469
00470 FourMomentum& operator*=(double a) {
00471 _vec = multiply(a, *this)._vec;
00472 return *this;
00473 }
00474
00475
00476 FourMomentum& operator/=(double a) {
00477 _vec = multiply(1.0/a, *this)._vec;
00478 return *this;
00479 }
00480
00481 FourMomentum& operator+=(const FourMomentum& v) {
00482 _vec = add(*this, v)._vec;
00483 return *this;
00484 }
00485
00486 FourMomentum& operator-=(const FourMomentum& v) {
00487 _vec = add(*this, -v)._vec;
00488 return *this;
00489 }
00490
00491 FourMomentum operator-() const {
00492 FourMomentum result;
00493 result._vec = -_vec;
00494 return result;
00495 }
00496
00497
00498 };
00499
00500
00501 inline FourMomentum multiply(const double a, const FourMomentum& v) {
00502 FourMomentum result;
00503 result._vec = a * v._vec;
00504 return result;
00505 }
00506
00507 inline FourMomentum multiply(const FourMomentum& v, const double a) {
00508 return multiply(a, v);
00509 }
00510
00511 inline FourMomentum operator*(const double a, const FourMomentum& v) {
00512 return multiply(a, v);
00513 }
00514
00515 inline FourMomentum operator*(const FourMomentum& v, const double a) {
00516 return multiply(a, v);
00517 }
00518
00519 inline FourMomentum operator/(const FourMomentum& v, const double a) {
00520 return multiply(1.0/a, v);
00521 }
00522
00523 inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
00524 FourMomentum result;
00525 result._vec = a._vec + b._vec;
00526 return result;
00527 }
00528
00529 inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
00530 return add(a, b);
00531 }
00532
00533 inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
00534 return add(a, -b);
00535 }
00536
00537
00538
00539
00540 inline double mass(const FourMomentum& v) {
00541 return v.mass();
00542 }
00543
00544
00545 inline double massT(const FourMomentum& v) {
00546 return v.massT();
00547 }
00548
00549
00550 inline double mass2(const FourMomentum& v) {
00551 return v.mass2();
00552 }
00553
00554
00555 inline double massT2(const FourMomentum& v) {
00556 return v.massT2();
00557 }
00558
00559
00560 inline double rapidity(const FourMomentum& v) {
00561 return v.rapidity();
00562 }
00563
00564
00565 inline double pT2(const FourMomentum& v) {
00566 return v.pT2();
00567 }
00568
00569
00570 inline double pT(const FourMomentum& v) {
00571 return v.pT();
00572 }
00573
00574
00575 inline double Et2(const FourMomentum& v) {
00576 return v.Et2();}
00577
00578
00579 inline double Et(const FourMomentum& v) {
00580 return v.Et();
00581 }
00582
00583
00584 inline Vector3 boostVector(const FourMomentum& v) {
00585 return v.boostVector();
00586 }
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601
00602 inline double deltaR(const FourVector& a, const FourVector& b,
00603 RapScheme scheme = PSEUDORAPIDITY) {
00604 switch (scheme) {
00605 case PSEUDORAPIDITY :
00606 return deltaR(a.vector3(), b.vector3());
00607 case RAPIDITY:
00608 {
00609 const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
00610 const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
00611 if (!ma || !mb) {
00612 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00613 throw std::runtime_error(err);
00614 }
00615 return deltaR(*ma, *mb, scheme);
00616 }
00617 default:
00618 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00619 }
00620 }
00621
00622
00623
00624
00625
00626
00627
00628 inline double deltaR(const FourVector& v,
00629 double eta2, double phi2,
00630 RapScheme scheme = PSEUDORAPIDITY) {
00631 switch (scheme) {
00632 case PSEUDORAPIDITY :
00633 return deltaR(v.vector3(), eta2, phi2);
00634 case RAPIDITY:
00635 {
00636 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00637 if (!mv) {
00638 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00639 throw std::runtime_error(err);
00640 }
00641 return deltaR(*mv, eta2, phi2, scheme);
00642 }
00643 default:
00644 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00645 }
00646 }
00647
00648
00649
00650
00651
00652
00653
00654 inline double deltaR(double eta1, double phi1,
00655 const FourVector& v,
00656 RapScheme scheme = PSEUDORAPIDITY) {
00657 switch (scheme) {
00658 case PSEUDORAPIDITY :
00659 return deltaR(eta1, phi1, v.vector3());
00660 case RAPIDITY:
00661 {
00662 const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00663 if (!mv) {
00664 string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00665 throw std::runtime_error(err);
00666 }
00667 return deltaR(eta1, phi1, *mv, scheme);
00668 }
00669 default:
00670 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00671 }
00672 }
00673
00674
00675
00676
00677
00678
00679
00680 inline double deltaR(const FourMomentum& a, const FourMomentum& b,
00681 RapScheme scheme = PSEUDORAPIDITY) {
00682 switch (scheme) {
00683 case PSEUDORAPIDITY:
00684 return deltaR(a.vector3(), b.vector3());
00685 case RAPIDITY:
00686 return deltaR(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
00687 default:
00688 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00689 }
00690 }
00691
00692
00693
00694
00695
00696
00697 inline double deltaR(const FourMomentum& v,
00698 double eta2, double phi2,
00699 RapScheme scheme = PSEUDORAPIDITY) {
00700 switch (scheme) {
00701 case PSEUDORAPIDITY:
00702 return deltaR(v.vector3(), eta2, phi2);
00703 case RAPIDITY:
00704 return deltaR(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
00705 default:
00706 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00707 }
00708 }
00709
00710
00711
00712
00713
00714
00715
00716 inline double deltaR(double eta1, double phi1,
00717 const FourMomentum& v,
00718 RapScheme scheme = PSEUDORAPIDITY) {
00719 switch (scheme) {
00720 case PSEUDORAPIDITY:
00721 return deltaR(eta1, phi1, v.vector3());
00722 case RAPIDITY:
00723 return deltaR(eta1, phi1, v.rapidity(), v.azimuthalAngle());
00724 default:
00725 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00726 }
00727 }
00728
00729
00730
00731
00732
00733
00734 inline double deltaR(const FourMomentum& a, const FourVector& b,
00735 RapScheme scheme = PSEUDORAPIDITY) {
00736 switch (scheme) {
00737 case PSEUDORAPIDITY:
00738 return deltaR(a.vector3(), b.vector3());
00739 case RAPIDITY:
00740 return deltaR(a.rapidity(), a.azimuthalAngle(), FourMomentum(b).rapidity(), b.azimuthalAngle());
00741 default:
00742 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00743 }
00744 }
00745
00746
00747
00748
00749
00750
00751 inline double deltaR(const FourVector& a, const FourMomentum& b,
00752 RapScheme scheme = PSEUDORAPIDITY) {
00753 switch (scheme) {
00754 case PSEUDORAPIDITY:
00755 return deltaR(a.vector3(), b.vector3());
00756 case RAPIDITY:
00757 return deltaR(FourMomentum(a).rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
00758 default:
00759 throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00760 }
00761 }
00762
00763
00764
00765 inline double deltaR(const FourMomentum& a, const Vector3& b) {
00766 return deltaR(a.vector3(), b);
00767 }
00768
00769
00770
00771 inline double deltaR(const Vector3& a, const FourMomentum& b) {
00772 return deltaR(a, b.vector3());
00773 }
00774
00775
00776
00777 inline double deltaR(const FourVector& a, const Vector3& b) {
00778 return deltaR(a.vector3(), b);
00779 }
00780
00781
00782
00783 inline double deltaR(const Vector3& a, const FourVector& b) {
00784 return deltaR(a, b.vector3());
00785 }
00786
00787
00788
00789
00790
00791
00792
00793 inline double deltaPhi(const FourMomentum& a, const FourMomentum& b) {
00794 return deltaPhi(a.vector3(), b.vector3());
00795 }
00796
00797
00798 inline double deltaPhi(const FourMomentum& v, double phi2) {
00799 return deltaPhi(v.vector3(), phi2);
00800 }
00801
00802
00803 inline double deltaPhi(double phi1, const FourMomentum& v) {
00804 return deltaPhi(phi1, v.vector3());
00805 }
00806
00807
00808 inline double deltaPhi(const FourVector& a, const FourVector& b) {
00809 return deltaPhi(a.vector3(), b.vector3());
00810 }
00811
00812
00813 inline double deltaPhi(const FourVector& v, double phi2) {
00814 return deltaPhi(v.vector3(), phi2);
00815 }
00816
00817
00818 inline double deltaPhi(double phi1, const FourVector& v) {
00819 return deltaPhi(phi1, v.vector3());
00820 }
00821
00822
00823 inline double deltaPhi(const FourVector& a, const FourMomentum& b) {
00824 return deltaPhi(a.vector3(), b.vector3());
00825 }
00826
00827
00828 inline double deltaPhi(const FourMomentum& a, const FourVector& b) {
00829 return deltaPhi(a.vector3(), b.vector3());
00830 }
00831
00832
00833 inline double deltaPhi(const FourVector& a, const Vector3& b) {
00834 return deltaPhi(a.vector3(), b);
00835 }
00836
00837
00838 inline double deltaPhi(const Vector3& a, const FourVector& b) {
00839 return deltaPhi(a, b.vector3());
00840 }
00841
00842
00843 inline double deltaPhi(const FourMomentum& a, const Vector3& b) {
00844 return deltaPhi(a.vector3(), b);
00845 }
00846
00847
00848 inline double deltaPhi(const Vector3& a, const FourMomentum& b) {
00849 return deltaPhi(a, b.vector3());
00850 }
00851
00852
00853
00854
00855
00856
00857
00858 inline double deltaEta(const FourMomentum& a, const FourMomentum& b) {
00859 return deltaEta(a.vector3(), b.vector3());
00860 }
00861
00862
00863 inline double deltaEta(const FourMomentum& v, double eta2) {
00864 return deltaEta(v.vector3(), eta2);
00865 }
00866
00867
00868 inline double deltaEta(double eta1, const FourMomentum& v) {
00869 return deltaEta(eta1, v.vector3());
00870 }
00871
00872
00873 inline double deltaEta(const FourVector& a, const FourVector& b) {
00874 return deltaEta(a.vector3(), b.vector3());
00875 }
00876
00877
00878 inline double deltaEta(const FourVector& v, double eta2) {
00879 return deltaEta(v.vector3(), eta2);
00880 }
00881
00882
00883 inline double deltaEta(double eta1, const FourVector& v) {
00884 return deltaEta(eta1, v.vector3());
00885 }
00886
00887
00888 inline double deltaEta(const FourVector& a, const FourMomentum& b) {
00889 return deltaEta(a.vector3(), b.vector3());
00890 }
00891
00892
00893 inline double deltaEta(const FourMomentum& a, const FourVector& b) {
00894 return deltaEta(a.vector3(), b.vector3());
00895 }
00896
00897
00898 inline double deltaEta(const FourVector& a, const Vector3& b) {
00899 return deltaEta(a.vector3(), b);
00900 }
00901
00902
00903 inline double deltaEta(const Vector3& a, const FourVector& b) {
00904 return deltaEta(a, b.vector3());
00905 }
00906
00907
00908 inline double deltaEta(const FourMomentum& a, const Vector3& b) {
00909 return deltaEta(a.vector3(), b);
00910 }
00911
00912
00913 inline double deltaEta(const Vector3& a, const FourMomentum& b) {
00914 return deltaEta(a, b.vector3());
00915 }
00916
00917
00918
00919
00920
00921
00922
00923
00924
00925
00926 inline const string toString(const FourVector& lv) {
00927 ostringstream out;
00928 out << "(" << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
00929 << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
00930 << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
00931 << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
00932 << ")";
00933 return out.str();
00934 }
00935
00936
00937 inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
00938 out << toString(lv);
00939 return out;
00940 }
00941
00942
00943
00944
00945 }
00946
00947 #endif