File indexing completed on 2025-12-15 10:30:26
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #ifndef ROOT_TVector3
0012 #define ROOT_TVector3
0013
0014 #include "TError.h"
0015 #include "TVector2.h"
0016 #include "TMatrix.h"
0017 #include "TMath.h"
0018
0019 class TRotation;
0020
0021
0022 class TVector3 : public TObject {
0023
0024 public:
0025
0026 typedef Double_t Scalar;
0027
0028 TVector3();
0029
0030 TVector3(Double_t x, Double_t y, Double_t z);
0031
0032
0033 TVector3(const Double_t *);
0034 TVector3(const Float_t *);
0035
0036
0037 TVector3(const TVector3 &) noexcept;
0038
0039
0040 ~TVector3() override = default;
0041
0042
0043
0044 constexpr std::size_t size() const { return 3; }
0045
0046 Double_t operator () (int) const;
0047 inline Double_t operator [] (int) const;
0048
0049
0050 Double_t & operator () (int);
0051 inline Double_t & operator [] (int);
0052
0053
0054 inline Double_t x() const;
0055 inline Double_t y() const;
0056 inline Double_t z() const;
0057 inline Double_t X() const;
0058 inline Double_t Y() const;
0059 inline Double_t Z() const;
0060 inline Double_t Px() const;
0061 inline Double_t Py() const;
0062 inline Double_t Pz() const;
0063
0064
0065 inline void SetX(Double_t);
0066 inline void SetY(Double_t);
0067 inline void SetZ(Double_t);
0068 inline void SetXYZ(Double_t x, Double_t y, Double_t z);
0069 void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi);
0070 void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi);
0071
0072 inline void GetXYZ(Double_t *carray) const;
0073 inline void GetXYZ(Float_t *carray) const;
0074
0075
0076
0077 Double_t Phi() const;
0078
0079
0080 Double_t Theta() const;
0081
0082
0083 inline Double_t CosTheta() const;
0084
0085
0086 inline Double_t Mag2() const;
0087
0088
0089 Double_t Mag() const { return TMath::Sqrt(Mag2()); }
0090
0091
0092 void SetPhi(Double_t);
0093
0094
0095 void SetTheta(Double_t);
0096
0097
0098 inline void SetMag(Double_t);
0099
0100
0101 inline Double_t Perp2() const;
0102
0103
0104 inline Double_t Pt() const;
0105 Double_t Perp() const;
0106
0107
0108 inline void SetPerp(Double_t);
0109
0110
0111 inline Double_t Perp2(const TVector3 &) const;
0112
0113
0114 inline Double_t Pt(const TVector3 &) const;
0115 Double_t Perp(const TVector3 &) const;
0116
0117
0118 inline Double_t DeltaPhi(const TVector3 &) const;
0119 Double_t DeltaR(const TVector3 &) const;
0120 inline Double_t DrEtaPhi(const TVector3 &) const;
0121 inline TVector2 EtaPhiVector() const;
0122 void SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi);
0123
0124 inline TVector3 & operator = (const TVector3 &);
0125
0126
0127 inline Bool_t operator == (const TVector3 &) const;
0128 inline Bool_t operator != (const TVector3 &) const;
0129
0130
0131 inline TVector3 & operator += (const TVector3 &);
0132
0133
0134 inline TVector3 & operator -= (const TVector3 &);
0135
0136
0137 inline TVector3 operator - () const;
0138
0139
0140 inline TVector3 & operator *= (Double_t);
0141
0142
0143 TVector3 Unit() const;
0144
0145
0146 inline TVector3 Orthogonal() const;
0147
0148
0149 inline Double_t Dot(const TVector3 &) const;
0150
0151
0152 inline TVector3 Cross(const TVector3 &) const;
0153
0154
0155 Double_t Angle(const TVector3 &) const;
0156
0157
0158 Double_t PseudoRapidity() const;
0159
0160
0161 inline Double_t Eta() const;
0162
0163 void RotateX(Double_t);
0164
0165
0166 void RotateY(Double_t);
0167
0168
0169 void RotateZ(Double_t);
0170
0171
0172 void RotateUz(const TVector3&);
0173
0174
0175 void Rotate(Double_t, const TVector3 &);
0176
0177
0178 TVector3 & operator *= (const TRotation &);
0179 TVector3 & Transform(const TRotation &);
0180
0181
0182 inline TVector2 XYvector() const;
0183
0184 void Print(Option_t* option="") const override;
0185
0186 private:
0187
0188 Double_t fX, fY, fZ;
0189
0190
0191 ClassDefOverride(TVector3,3)
0192
0193
0194 friend class TLorentzVector;
0195 };
0196
0197 TVector3 operator + (const TVector3 &, const TVector3 &);
0198
0199
0200 TVector3 operator - (const TVector3 &, const TVector3 &);
0201
0202
0203 Double_t operator * (const TVector3 &, const TVector3 &);
0204
0205
0206 TVector3 operator * (const TVector3 &, Double_t a);
0207 TVector3 operator * (Double_t a, const TVector3 &);
0208
0209
0210 TVector3 operator * (const TMatrix &, const TVector3 &);
0211
0212
0213 inline Double_t & TVector3::operator[] (int i) { return operator()(i); }
0214 inline Double_t TVector3::operator[] (int i) const { return operator()(i); }
0215
0216 inline Double_t TVector3::x() const { return fX; }
0217 inline Double_t TVector3::y() const { return fY; }
0218 inline Double_t TVector3::z() const { return fZ; }
0219 inline Double_t TVector3::X() const { return fX; }
0220 inline Double_t TVector3::Y() const { return fY; }
0221 inline Double_t TVector3::Z() const { return fZ; }
0222 inline Double_t TVector3::Px() const { return fX; }
0223 inline Double_t TVector3::Py() const { return fY; }
0224 inline Double_t TVector3::Pz() const { return fZ; }
0225
0226 inline void TVector3::SetX(Double_t xx) { fX = xx; }
0227 inline void TVector3::SetY(Double_t yy) { fY = yy; }
0228 inline void TVector3::SetZ(Double_t zz) { fZ = zz; }
0229
0230 inline void TVector3::SetXYZ(Double_t xx, Double_t yy, Double_t zz) {
0231 fX = xx;
0232 fY = yy;
0233 fZ = zz;
0234 }
0235
0236 inline void TVector3::GetXYZ(Double_t *carray) const {
0237 carray[0] = fX;
0238 carray[1] = fY;
0239 carray[2] = fZ;
0240 }
0241
0242 inline void TVector3::GetXYZ(Float_t *carray) const {
0243 carray[0] = fX;
0244 carray[1] = fY;
0245 carray[2] = fZ;
0246 }
0247
0248
0249
0250 inline TVector3::TVector3()
0251 : fX(0.0), fY(0.0), fZ(0.0) {}
0252
0253 inline TVector3::TVector3(const TVector3 & p) noexcept : TObject(p),
0254 fX(p.fX), fY(p.fY), fZ(p.fZ) {}
0255
0256 inline TVector3::TVector3(Double_t xx, Double_t yy, Double_t zz)
0257 : fX(xx), fY(yy), fZ(zz) {}
0258
0259 inline TVector3::TVector3(const Double_t * x0)
0260 : fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
0261
0262 inline TVector3::TVector3(const Float_t * x0)
0263 : fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
0264
0265
0266 inline Double_t TVector3::operator () (int i) const {
0267 switch(i) {
0268 case 0:
0269 return fX;
0270 case 1:
0271 return fY;
0272 case 2:
0273 return fZ;
0274 default:
0275 Error("operator()(i)", "bad index (%d) returning 0",i);
0276 }
0277 return 0.;
0278 }
0279
0280 inline Double_t & TVector3::operator () (int i) {
0281 switch(i) {
0282 case 0:
0283 return fX;
0284 case 1:
0285 return fY;
0286 case 2:
0287 return fZ;
0288 default:
0289 Error("operator()(i)", "bad index (%d) returning &fX",i);
0290 }
0291 return fX;
0292 }
0293
0294 inline TVector3 & TVector3::operator = (const TVector3 & p) {
0295 fX = p.fX;
0296 fY = p.fY;
0297 fZ = p.fZ;
0298 return *this;
0299 }
0300
0301 inline Bool_t TVector3::operator == (const TVector3& v) const {
0302 return (v.fX==fX && v.fY==fY && v.fZ==fZ) ? kTRUE : kFALSE;
0303 }
0304
0305 inline Bool_t TVector3::operator != (const TVector3& v) const {
0306 return (v.fX!=fX || v.fY!=fY || v.fZ!=fZ) ? kTRUE : kFALSE;
0307 }
0308
0309 inline TVector3& TVector3::operator += (const TVector3 & p) {
0310 fX += p.fX;
0311 fY += p.fY;
0312 fZ += p.fZ;
0313 return *this;
0314 }
0315
0316 inline TVector3& TVector3::operator -= (const TVector3 & p) {
0317 fX -= p.fX;
0318 fY -= p.fY;
0319 fZ -= p.fZ;
0320 return *this;
0321 }
0322
0323 inline TVector3 TVector3::operator - () const {
0324 return TVector3(-fX, -fY, -fZ);
0325 }
0326
0327 inline TVector3& TVector3::operator *= (Double_t a) {
0328 fX *= a;
0329 fY *= a;
0330 fZ *= a;
0331 return *this;
0332 }
0333
0334 inline Double_t TVector3::Dot(const TVector3 & p) const {
0335 return fX*p.fX + fY*p.fY + fZ*p.fZ;
0336 }
0337
0338 inline TVector3 TVector3::Cross(const TVector3 & p) const {
0339 return TVector3(fY*p.fZ-p.fY*fZ, fZ*p.fX-p.fZ*fX, fX*p.fY-p.fX*fY);
0340 }
0341
0342 inline Double_t TVector3::Mag2() const { return fX*fX + fY*fY + fZ*fZ; }
0343
0344
0345 inline TVector3 TVector3::Orthogonal() const {
0346 Double_t xx = fX < 0.0 ? -fX : fX;
0347 Double_t yy = fY < 0.0 ? -fY : fY;
0348 Double_t zz = fZ < 0.0 ? -fZ : fZ;
0349 if (xx < yy) {
0350 return xx < zz ? TVector3(0,fZ,-fY) : TVector3(fY,-fX,0);
0351 } else {
0352 return yy < zz ? TVector3(-fZ,0,fX) : TVector3(fY,-fX,0);
0353 }
0354 }
0355
0356 inline Double_t TVector3::Perp2() const { return fX*fX + fY*fY; }
0357
0358
0359 inline Double_t TVector3::Pt() const { return Perp(); }
0360
0361 inline Double_t TVector3::Perp2(const TVector3 & p) const {
0362 Double_t tot = p.Mag2();
0363 Double_t ss = Dot(p);
0364 Double_t per = Mag2();
0365 if (tot > 0.0) per -= ss*ss/tot;
0366 if (per < 0) per = 0;
0367 return per;
0368 }
0369
0370 inline Double_t TVector3::Pt(const TVector3 & p) const {
0371 return Perp(p);
0372 }
0373
0374 inline Double_t TVector3::CosTheta() const {
0375 Double_t ptot = Mag();
0376 return ptot == 0.0 ? 1.0 : fZ/ptot;
0377 }
0378
0379 inline void TVector3::SetMag(Double_t ma) {
0380 Double_t factor = Mag();
0381 if (factor == 0) {
0382 Warning("SetMag","zero vector can't be stretched");
0383 } else {
0384 factor = ma/factor;
0385 SetX(fX*factor);
0386 SetY(fY*factor);
0387 SetZ(fZ*factor);
0388 }
0389 }
0390
0391 inline void TVector3::SetPerp(Double_t r) {
0392 Double_t p = Perp();
0393 if (p != 0.0) {
0394 fX *= r/p;
0395 fY *= r/p;
0396 }
0397 }
0398
0399 inline Double_t TVector3::DeltaPhi(const TVector3 & v) const {
0400 return TVector2::Phi_mpi_pi(Phi()-v.Phi());
0401 }
0402
0403 inline Double_t TVector3::Eta() const {
0404 return PseudoRapidity();
0405 }
0406
0407 inline Double_t TVector3::DrEtaPhi(const TVector3 & v) const{
0408 return DeltaR(v);
0409 }
0410
0411
0412 inline TVector2 TVector3::EtaPhiVector() const {
0413 return TVector2 (Eta(),Phi());
0414 }
0415
0416 inline TVector2 TVector3::XYvector() const {
0417 return TVector2(fX,fY);
0418 }
0419
0420
0421 #endif