Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:12:19

0001 // @(#)root/eg:$Id$
0002 // Author: Rene Brun , Federico Carminati  26/04/99
0003 
0004 /*************************************************************************
0005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
0006  * All rights reserved.                                                  *
0007  *                                                                       *
0008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
0009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
0010  *************************************************************************/
0011 //////////////////////////////////////////////////////////////////////////
0012 //                                                                      //
0013 // TParticle: defines  equivalent of HEPEVT particle                    //
0014 //////////////////////////////////////////////////////////////////////////
0015 
0016 #ifndef ROOT_TParticle
0017 #define ROOT_TParticle
0018 
0019 #include "TNamed.h"
0020 #include "TAttLine.h"
0021 #include "TAtt3D.h"
0022 #include "TLorentzVector.h"
0023 
0024 class TParticlePDG;
0025 
0026 class TParticle : public TObject, public TAttLine, public TAtt3D {
0027 
0028 
0029 protected:
0030 
0031   Int_t          fPdgCode;              // PDG code of the particle
0032   Int_t          fStatusCode;           // generation status code
0033   Int_t          fMother[2];            // Indices of the mother particles
0034   Int_t          fDaughter[2];          // Indices of the daughter particles
0035   Float_t        fWeight;               // particle weight
0036 
0037   Double_t       fCalcMass;             // Calculated mass
0038 
0039   Double_t       fPx;                   // x component of momentum
0040   Double_t       fPy;                   // y component of momentum
0041   Double_t       fPz;                   // z component of momentum
0042   Double_t       fE;                    // Energy
0043 
0044   Double_t       fVx;                   // x of production vertex
0045   Double_t       fVy;                   // y of production vertex
0046   Double_t       fVz;                   // z of production vertex
0047   Double_t       fVt;                   // t of production vertex
0048 
0049   Double_t       fPolarTheta;           // Polar angle of polarisation
0050   Double_t       fPolarPhi;             // azymutal angle of polarisation
0051 
0052   mutable TParticlePDG* fParticlePDG;   //! reference to the particle record in PDG database
0053   //----------------------------------------------------------------------------
0054   //  functions
0055   //----------------------------------------------------------------------------
0056 public:
0057                                 // ****** constructors and destructor
0058    TParticle();
0059 
0060    TParticle(Int_t pdg, Int_t status,
0061              Int_t mother1, Int_t mother2,
0062              Int_t daughter1, Int_t daughter2,
0063              Double_t px, Double_t py, Double_t pz, Double_t etot,
0064              Double_t vx, Double_t vy, Double_t vz, Double_t time);
0065 
0066    TParticle(Int_t pdg, Int_t status,
0067              Int_t mother1, Int_t mother2,
0068              Int_t daughter1, Int_t daughter2,
0069              const TLorentzVector &p,
0070              const TLorentzVector &v);
0071 
0072    TParticle(const TParticle &part);
0073 
0074    ~TParticle() override;
0075 
0076    TParticle& operator=(const TParticle&);
0077 
0078 //   virtual TString* Name   () const { return fName.Data(); }
0079 //   virtual char*  GetName()   const { return fName.Data(); }
0080 
0081    Double_t       Ek              ()            const { return fE-fCalcMass;                                    }
0082    Int_t          GetStatusCode   ()            const { return fStatusCode;                                     }
0083    Int_t          GetPdgCode      ()            const { return fPdgCode;                                        }
0084    Int_t          GetFirstMother  ()            const { return fMother[0];                                      }
0085    Int_t          GetMother       (Int_t i)     const { return fMother[i];                                      }
0086    Int_t          GetSecondMother ()            const { return fMother[1];                                      }
0087    Bool_t         IsPrimary       ()            const { return fMother[0]>-1 ? kFALSE : kTRUE;                  } //Is this particle primary one?
0088    Int_t          GetFirstDaughter()            const { return fDaughter[0];                                    }
0089    Int_t          GetDaughter     (Int_t i)     const { return fDaughter[i];                                    }
0090    Int_t          GetLastDaughter ()            const { return fDaughter[1];                                    }
0091    Double_t       GetCalcMass     ()            const { return fCalcMass;                                       }
0092    Double_t       GetMass         ()            const;
0093    Int_t          GetNDaughters   ()            const { return fDaughter[1]>0 ? fDaughter[1]-fDaughter[0]+1 : 0;}
0094    Float_t        GetWeight       ()            const { return fWeight;                                         }
0095    void           GetPolarisation(TVector3 &v)  const;
0096    TParticlePDG*  GetPDG          (Int_t mode = 0) const;
0097    Int_t          Beauty          ()            const;
0098    Int_t          Charm           ()            const;
0099    Int_t          Strangeness     ()            const;
0100    Double_t       PhiX            ()            const { return TMath::Pi()+TMath::ATan2(-fPz,-fPy);             } // note that PhiX() returns an angle between 0 and 2pi
0101    Double_t       PhiY            ()            const { return TMath::Pi()+TMath::ATan2(-fPx,-fPz);             } // note that PhiY() returns an angle between 0 and 2pi
0102    Double_t       PhiZ            ()            const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx);             } // note that PhiZ() returns an angle between 0 and 2pi
0103    Double_t       ThetaX          ()            const { return (fPx==0)?TMath::PiOver2():TMath::ACos(fPx/P());  }
0104    Double_t       ThetaY          ()            const { return (fPy==0)?TMath::PiOver2():TMath::ACos(fPy/P());  }
0105    Double_t       ThetaZ          ()            const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P());  }
0106    Double_t       GetPolarTheta   ()            const { return fPolarTheta;                                     }
0107    Double_t       GetPolarPhi     ()            const { return fPolarPhi;                                       }
0108    void           GetPolarisation (Double_t &theta, Double_t &phi) const { theta = fPolarTheta; phi = fPolarPhi; }
0109    void           SetPolarTheta   (Double_t theta) { fPolarTheta = theta; }
0110    void           SetPolarPhi     (Double_t phi) { fPolarPhi = phi; }
0111    void           SetPolarisation (Double_t theta, Double_t phi) { fPolarTheta = theta; fPolarPhi = phi; }
0112    void           Momentum(TLorentzVector &v)   const { v.SetPxPyPzE(fPx,fPy,fPz,fE);                           }
0113    void           ProductionVertex(TLorentzVector &v) const { v.SetXYZT(fVx,fVy,fVz,fVt);                       }
0114 
0115    Double_t       Theta(const TParticle &p) // Returns the angle between momenta of particles
0116    {
0117       Double_t v = P()*p.P();
0118       if (v == 0) v = 1; else v = (fPx*p.Px()+fPy*p.Py()+fPz*p.Pz())/v;
0119       if (v > 1) v = 1; else if (v < -1) v = -1; // just a precaution
0120       return TMath::ACos(v);
0121    }
0122 
0123 // ****** redefine several most oftenly used methods of LORENTZ_VECTOR
0124 
0125    Double_t       Vx              ()            const { return fVx;                                             }
0126    Double_t       Vy              ()            const { return fVy;                                             }
0127    Double_t       Vz              ()            const { return fVz;                                             }
0128    Double_t       T               ()            const { return fVt;                                             }
0129    Double_t       R               ()            const { return TMath::Sqrt(fVx*fVx+fVy*fVy);                    } //Radius of production vertex in cylindrical system
0130    Double_t       Rho             ()            const { return TMath::Sqrt(fVx*fVx+fVy*fVy+fVz*fVz);            } //Radius of production vertex in spherical system
0131    Double_t       Px              ()            const { return fPx;                                             }
0132    Double_t       Py              ()            const { return fPy;                                             }
0133    Double_t       Pz              ()            const { return fPz;                                             }
0134    Double_t       P               ()            const { return TMath::Sqrt(fPx*fPx+fPy*fPy+fPz*fPz);            }
0135    Double_t       Pt              ()            const { return TMath::Sqrt(fPx*fPx+fPy*fPy);                    }
0136    Double_t       Energy          ()            const { return fE;                                              }
0137    Double_t       Eta             ()            const
0138    {
0139       Double_t pmom = P();
0140       if (pmom != TMath::Abs(fPz)) return 0.5*TMath::Log((pmom+fPz)/(pmom-fPz));
0141       else                         return 1.e30;
0142    }
0143    Double_t       Y               ()            const
0144    {
0145       if (fE != TMath::Abs(fPz))   return 0.5*TMath::Log((fE+fPz)/(fE-fPz));
0146       else                         return 1.e30;
0147    }
0148 
0149    Double_t         Phi   () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); }  // note that Phi() returns an angle between 0 and 2pi
0150    Double_t         Theta () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
0151 
0152               // setters
0153 
0154    void           SetFirstMother  (int code)                                               { fMother[0]   = code ; }
0155    void           SetMother  (int i, int code)                                             { fMother[i]   = code ; }
0156    void           SetLastMother  (int code)                                                { fMother[1]   = code ; }
0157    void           SetFirstDaughter(int code)                                               { fDaughter[0] = code ; }
0158    void           SetDaughter(int i, int code)                                             { fDaughter[i] = code ; }
0159    void           SetLastDaughter(int code)                                                { fDaughter[1] = code ; }
0160    void           SetCalcMass(Double_t mass)                                               { fCalcMass=mass;}
0161    void           SetPdgCode(Int_t pdg);
0162    void           SetPolarisation(Double_t polx, Double_t poly, Double_t polz);
0163    void           SetPolarisation(const TVector3& v)                                       {SetPolarisation(v.X(), v.Y(), v.Z());}
0164    void           SetStatusCode(int status)                                                {fStatusCode = status;}
0165    void           SetWeight(Float_t weight = 1)                                            {fWeight = weight; }
0166    void           SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e)           {fPx=px; fPy=py; fPz=pz; fE=e;}
0167    void           SetMomentum(const TLorentzVector& p)                                     {SetMomentum(p.Px(),p.Py(),p.Pz(),p.Energy());}
0168    void           SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t)   {fVx=vx; fVy=vy; fVz=vz; fVt=t;}
0169    void           SetProductionVertex(const TLorentzVector& v)                             {SetProductionVertex(v.X(),v.Y(),v.Z(),v.T());}
0170 
0171              // ****** overloaded functions of TObject
0172 
0173    void      Paint(Option_t *option = "") override;
0174    void      Print(Option_t *option = "") const override;
0175    void      Sizeof3D() const override;
0176    Int_t     DistancetoPrimitive(Int_t px, Int_t py) override;
0177    void      ExecuteEvent(Int_t event, Int_t px, Int_t py) override;
0178    const     char *GetName() const override;
0179    const     char *GetTitle() const override;
0180 
0181 
0182    ClassDefOverride(TParticle,2)  // TParticle vertex particle information
0183 };
0184 
0185 #endif