CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/SimDataFormats/TrackingAnalysis/interface/ParticleBase.h

Go to the documentation of this file.
00001 #ifndef SimDataFormats_ParticleBase_h
00002 #define SimDataFormats_ParticleBase_h
00003 
00013 #include "DataFormats/Math/interface/Point3D.h"
00014 #include "DataFormats/Math/interface/Vector3D.h"
00015 #include "DataFormats/Math/interface/LorentzVector.h"
00016 #include "DataFormats/Common/interface/BoolCache.h"
00017 #include "Rtypes.h"
00018 
00019 class ParticleBase
00020 {
00021 public:
00023     typedef int Charge;
00025     typedef math::XYZTLorentzVectorD LorentzVector;
00027     typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
00029     typedef math::XYZPointD Point;
00031     typedef math::XYZVectorD Vector;
00033     ParticleBase() : cachePolarFixed_( false ) { }
00035     ParticleBase( Charge q, const LorentzVector & p4, const Point & vertex = Point( 0, 0, 0 ),
00036                   int pdgId = 0, int status = 0, bool integerCharge = true ) :
00037             qx3_( q ), pt_( p4.pt() ), eta_( p4.eta() ), phi_( p4.phi() ), mass_( p4.mass() ),
00038             vertex_( vertex ), pdgId_( pdgId ), status_( status ),
00039             cachePolarFixed_( false ), cacheCartesianFixed_( false )
00040     {
00041         if ( integerCharge ) qx3_ *= 3;
00042     }
00044     ParticleBase( Charge q, const PolarLorentzVector & p4, const Point & vertex = Point( 0, 0, 0 ),
00045                   int pdgId = 0, int status = 0, bool integerCharge = true ) :
00046             qx3_( q ), pt_( p4.pt() ), eta_( p4.eta() ), phi_( p4.phi() ), mass_( p4.mass() ),
00047             vertex_( vertex ), pdgId_( pdgId ), status_( status ),
00048             cachePolarFixed_( false ), cacheCartesianFixed_( false )
00049     {
00050         if ( integerCharge ) qx3_ *= 3;
00051     }
00053     virtual ~ParticleBase() { }
00055     int charge() const
00056     {
00057         return qx3_ / 3;
00058     }
00060     void setCharge( Charge q )
00061     {
00062         qx3_ = q * 3;
00063     }
00065     int threeCharge() const
00066     {
00067         return qx3_;
00068     }
00070     void setThreeCharge( Charge qx3 )
00071     {
00072         qx3_ = qx3;
00073     }
00075     const LorentzVector & p4() const
00076     {
00077         cacheCartesian();
00078         return p4Cartesian_;
00079     }
00081     const PolarLorentzVector & polarP4() const
00082     {
00083         cachePolar();
00084         return p4Polar_;
00085     }
00087     Vector momentum() const
00088     {
00089         cacheCartesian();
00090         return p4Cartesian_.Vect();
00091     }
00094     Vector boostToCM() const
00095     {
00096         cacheCartesian();
00097         return p4Cartesian_.BoostToCM();
00098     }
00100     double p() const
00101     {
00102         cacheCartesian();
00103         return p4Cartesian_.P();
00104     }
00106     double energy() const
00107     {
00108         cacheCartesian();
00109         return p4Cartesian_.E();
00110     }
00112     double et() const
00113     {
00114         cachePolar();
00115         return p4Polar_.Et();
00116     }
00118     double mass() const
00119     {
00120         return mass_;
00121     }
00123     double massSqr() const
00124     {
00125         return mass_ * mass_;
00126     }
00128     double mt() const
00129     {
00130         cachePolar();
00131         return p4Polar_.Mt();
00132     }
00134     double mtSqr() const
00135     {
00136         cachePolar();
00137         return p4Polar_.Mt2();
00138     }
00140     double px() const
00141     {
00142         cacheCartesian();
00143         return p4Cartesian_.Px();
00144     }
00146     double py() const
00147     {
00148         cacheCartesian();
00149         return p4Cartesian_.Py();
00150     }
00152     double pz() const
00153     {
00154         cacheCartesian();
00155         return p4Cartesian_.Pz();
00156     }
00158     double pt() const
00159     {
00160         return pt_;
00161     }
00163     double phi() const
00164     {
00165         return phi_;
00166     }
00168     double theta() const
00169     {
00170         cacheCartesian();
00171         return p4Cartesian_.Theta();
00172     }
00174     double eta() const
00175     {
00176         return eta_;
00177     }
00179     double rapidity() const
00180     {
00181         cachePolar();
00182         return p4Polar_.Rapidity();
00183     }
00185     double y() const
00186     {
00187         return rapidity();
00188     }
00190     void setP4( const LorentzVector & p4 )
00191     {
00192         p4Cartesian_ = p4;
00193         p4Polar_ = p4;
00194         pt_ = p4Polar_.pt();
00195         eta_ = p4Polar_.eta();
00196         phi_ = p4Polar_.phi();
00197         mass_ = p4Polar_.mass();
00198         cachePolarFixed_ = true;
00199         cacheCartesianFixed_ = true;
00200     }
00202     void setP4( const PolarLorentzVector & p4 )
00203     {
00204         p4Polar_ = p4;
00205         pt_ = p4Polar_.pt();
00206         eta_ = p4Polar_.eta();
00207         phi_ = p4Polar_.phi();
00208         mass_ = p4Polar_.mass();
00209         cachePolarFixed_ = true;
00210         cacheCartesianFixed_ = false;
00211     }
00213     void setMass( double m )
00214     {
00215         mass_ = m;
00216         clearCache();
00217     }
00218     void setPz( double pz )
00219     {
00220         cacheCartesian();
00221         p4Cartesian_.SetPz(pz);
00222         p4Polar_ = p4Cartesian_;
00223         pt_ = p4Polar_.pt();
00224         eta_ = p4Polar_.eta();
00225         phi_ = p4Polar_.phi();
00226         mass_ = p4Polar_.mass();
00227     }
00229     const Point & vertex() const
00230     {
00231         return vertex_;
00232     }
00234     double vx() const
00235     {
00236         return vertex_.X();
00237     }
00239     double vy() const
00240     {
00241         return vertex_.Y();
00242     }
00244     double vz() const
00245     {
00246         return vertex_.Z();
00247     }
00249     void setVertex( const Point & vertex )
00250     {
00251         vertex_ = vertex;
00252     }
00254     int pdgId() const
00255     {
00256         return pdgId_;
00257     }
00258     // set PDG identifier
00259     void setPdgId( int pdgId )
00260     {
00261         pdgId_ = pdgId;
00262     }
00264     int status() const
00265     {
00266         return status_;
00267     }
00269     void setStatus( int status )
00270     {
00271         status_ = status;
00272     }
00274     static const unsigned int longLivedTag;
00276     void setLongLived()
00277     {
00278         status_ |= longLivedTag;
00279     }
00281     bool longLived() const
00282     {
00283         return status_ & longLivedTag;
00284     }
00285 
00286 protected:
00288     Charge qx3_;
00290     float pt_, eta_, phi_, mass_;
00292     Point vertex_;
00294     int pdgId_;
00296     int status_;
00298     mutable PolarLorentzVector p4Polar_;
00300     mutable LorentzVector p4Cartesian_;
00302     mutable  edm::BoolCache cachePolarFixed_, cacheCartesianFixed_;
00304     inline void cachePolar() const
00305     {
00306         if ( cachePolarFixed_ ) return;
00307         p4Polar_ = PolarLorentzVector( pt_, eta_, phi_, mass_ );
00308         cachePolarFixed_ = true;
00309     }
00311     inline void cacheCartesian() const
00312     {
00313         if ( cacheCartesianFixed_ ) return;
00314         cachePolar();
00315         p4Cartesian_ = p4Polar_;
00316         cacheCartesianFixed_ = true;
00317     }
00319     inline void clearCache() const
00320     {
00321         cachePolarFixed_ = false;
00322         cacheCartesianFixed_ = false;
00323     }
00324 };
00325 
00326 #endif