CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DataFormats/Candidate/interface/Particle.h

Go to the documentation of this file.
00001 #ifndef Candidate_Particle_h
00002 #define Candidate_Particle_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 namespace reco {
00020   
00021   class Particle {
00022   public:
00024     typedef int Charge;
00026     typedef math::XYZTLorentzVector LorentzVector;
00028     typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
00030     typedef math::XYZPoint Point;
00032     typedef math::XYZVector Vector;
00034     Particle() :
00035       qx3_(0), pt_(0), eta_(0), phi_(0), mass_(0), 
00036       vertex_(0, 0, 0), pdgId_(0), status_(0),
00037       cachePolarFixed_( false ) { }
00039     Particle( Charge q, const LorentzVector & p4, const Point & vertex = Point( 0, 0, 0 ),
00040               int pdgId = 0, int status = 0, bool integerCharge = true ) : 
00041       qx3_( q ), pt_( p4.pt() ), eta_( p4.eta() ), phi_( p4.phi() ), mass_( p4.mass() ),
00042       vertex_( vertex ), pdgId_( pdgId ), status_( status ),
00043       cachePolarFixed_( false ), cacheCartesianFixed_( false ) { 
00044       if ( integerCharge ) qx3_ *= 3;
00045     }
00047     Particle( Charge q, const PolarLorentzVector & p4, const Point & vertex = Point( 0, 0, 0 ),
00048               int pdgId = 0, int status = 0, bool integerCharge = true ) : 
00049       qx3_( q ), pt_( p4.pt() ), eta_( p4.eta() ), phi_( p4.phi() ), mass_( p4.mass() ),
00050       vertex_( vertex ), pdgId_( pdgId ), status_( status ),
00051       cachePolarFixed_( false ), cacheCartesianFixed_( false ) { 
00052       if ( integerCharge ) qx3_ *= 3;
00053     }
00055     virtual ~Particle() { }
00057     int charge() const { return qx3_ / 3; }
00059     void setCharge( Charge q ) { qx3_ = q * 3; }
00061     int threeCharge() const { return qx3_; }
00063     void setThreeCharge( Charge qx3 ) { qx3_ = qx3; }
00065     const LorentzVector & p4() const { cacheCartesian(); return p4Cartesian_; }
00067     const PolarLorentzVector & polarP4() const { cachePolar(); return p4Polar_; }
00069     Vector momentum() const { cacheCartesian(); return p4Cartesian_.Vect(); }
00072     Vector boostToCM() const { cacheCartesian(); return p4Cartesian_.BoostToCM(); }
00074     double p() const { cacheCartesian(); return p4Cartesian_.P(); }
00076     double energy() const { cacheCartesian(); return p4Cartesian_.E(); }  
00078     double et() const { cachePolar(); return p4Polar_.Et(); }  
00080     double mass() const { return mass_; }
00082     double massSqr() const { return mass_ * mass_; }
00084     double mt() const { cachePolar(); return p4Polar_.Mt(); }
00086     double mtSqr() const { cachePolar(); return p4Polar_.Mt2(); }
00088     double px() const { cacheCartesian(); return p4Cartesian_.Px(); }
00090     double py() const { cacheCartesian(); return p4Cartesian_.Py(); }
00092     double pz() const { cacheCartesian(); return p4Cartesian_.Pz(); }
00094     double pt() const { return pt_; }
00096     double phi() const { return phi_; }
00098     double theta() const { cacheCartesian(); return p4Cartesian_.Theta(); }
00100     double eta() const { return eta_; }
00102     double rapidity() const { cachePolar(); return p4Polar_.Rapidity(); }
00104     double y() const { return rapidity(); }
00106     void setP4( const LorentzVector & p4 ) { 
00107       p4Cartesian_ = p4;
00108       p4Polar_ = p4;
00109       pt_ = p4Polar_.pt();
00110       eta_ = p4Polar_.eta();
00111       phi_ = p4Polar_.phi();
00112       mass_ = p4Polar_.mass();
00113       cachePolarFixed_ = true;
00114       cacheCartesianFixed_ = true;      
00115     }
00117     void setP4( const PolarLorentzVector & p4 ) { 
00118       p4Polar_ = p4;
00119       pt_ = p4Polar_.pt();
00120       eta_ = p4Polar_.eta();
00121       phi_ = p4Polar_.phi();
00122       mass_ = p4Polar_.mass();
00123       cachePolarFixed_ = true;
00124       cacheCartesianFixed_ = false;            
00125     }
00127     void setMass( double m ) { 
00128       mass_ = m; 
00129       clearCache(); 
00130     }
00131     void setPz( double pz ) {
00132       cacheCartesian();
00133       p4Cartesian_.SetPz(pz);
00134       p4Polar_ = p4Cartesian_;
00135       pt_ = p4Polar_.pt();
00136       eta_ = p4Polar_.eta();
00137       phi_ = p4Polar_.phi();
00138       mass_ = p4Polar_.mass();
00139     }
00141     const Point & vertex() const { return vertex_; }
00143     double vx() const { return vertex_.X(); }
00145     double vy() const { return vertex_.Y(); }
00147     double vz() const { return vertex_.Z(); }
00149     void setVertex( const Point & vertex ) { vertex_ = vertex; }
00151     int pdgId() const { return pdgId_; }
00152     // set PDG identifier
00153     void setPdgId( int pdgId ) { pdgId_ = pdgId; }
00155     int status() const { return status_; }
00157     void setStatus( int status ) { status_ = status; }
00159     static const unsigned int longLivedTag;
00161     void setLongLived() { status_ |= longLivedTag; }
00163     bool longLived() const { return status_ & longLivedTag; }
00164 
00165   protected:
00167     Charge qx3_;   
00169     float pt_, eta_, phi_, mass_;
00171     Point vertex_;
00173     int pdgId_;
00175     int status_;
00177     mutable PolarLorentzVector p4Polar_;
00179     mutable LorentzVector p4Cartesian_;
00181     mutable  edm::BoolCache cachePolarFixed_, cacheCartesianFixed_;
00183     inline void cachePolar() const { 
00184       if ( cachePolarFixed_ ) return;
00185       p4Polar_ = PolarLorentzVector( pt_, eta_, phi_, mass_ );
00186       cachePolarFixed_ = true;
00187     }
00189     inline void cacheCartesian() const { 
00190       if ( cacheCartesianFixed_ ) return;
00191       cachePolar();
00192       p4Cartesian_ = p4Polar_;
00193       cacheCartesianFixed_ = true;
00194     }
00196     inline void clearCache() const { 
00197       cachePolarFixed_ = false;
00198       cacheCartesianFixed_ = false;
00199     }
00200   };
00201 
00202 }
00203 
00204 #endif