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 "Rtypes.h"
00017
00018 class ParticleBase
00019 {
00020 public:
00022 typedef int Charge;
00024 typedef math::XYZTLorentzVectorD LorentzVector;
00026 typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
00028 typedef math::XYZPointD Point;
00030 typedef math::XYZVectorD Vector;
00032 ParticleBase() : cachePolarFixed_( false ), cacheCartesianFixed_( false ) { }
00034 ParticleBase( Charge q, const LorentzVector & p4, const Point & vertex = Point( 0, 0, 0 ),
00035 int pdgId = 0, int status = 0, bool integerCharge = true ) :
00036 qx3_( q ), pt_( p4.pt() ), eta_( p4.eta() ), phi_( p4.phi() ), mass_( p4.mass() ),
00037 vertex_( vertex ), pdgId_( pdgId ), status_( status ),
00038 cachePolarFixed_( false ), cacheCartesianFixed_( false )
00039 {
00040 if ( integerCharge ) qx3_ *= 3;
00041 }
00043 ParticleBase( Charge q, const PolarLorentzVector & p4, const Point & vertex = Point( 0, 0, 0 ),
00044 int pdgId = 0, int status = 0, bool integerCharge = true ) :
00045 qx3_( q ), pt_( p4.pt() ), eta_( p4.eta() ), phi_( p4.phi() ), mass_( p4.mass() ),
00046 vertex_( vertex ), pdgId_( pdgId ), status_( status ),
00047 cachePolarFixed_( false ), cacheCartesianFixed_( false )
00048 {
00049 if ( integerCharge ) qx3_ *= 3;
00050 }
00052 virtual ~ParticleBase() { }
00054 int charge() const
00055 {
00056 return qx3_ / 3;
00057 }
00059 void setCharge( Charge q )
00060 {
00061 qx3_ = q * 3;
00062 }
00064 int threeCharge() const
00065 {
00066 return qx3_;
00067 }
00069 void setThreeCharge( Charge qx3 )
00070 {
00071 qx3_ = qx3;
00072 }
00074 const LorentzVector & p4() const
00075 {
00076 cacheCartesian();
00077 return p4Cartesian_;
00078 }
00080 const PolarLorentzVector & polarP4() const
00081 {
00082 cachePolar();
00083 return p4Polar_;
00084 }
00086 Vector momentum() const
00087 {
00088 cacheCartesian();
00089 return p4Cartesian_.Vect();
00090 }
00093 Vector boostToCM() const
00094 {
00095 cacheCartesian();
00096 return p4Cartesian_.BoostToCM();
00097 }
00099 double p() const
00100 {
00101 cacheCartesian();
00102 return p4Cartesian_.P();
00103 }
00105 double energy() const
00106 {
00107 cacheCartesian();
00108 return p4Cartesian_.E();
00109 }
00111 double et() const
00112 {
00113 cachePolar();
00114 return p4Polar_.Et();
00115 }
00117 double mass() const
00118 {
00119 return mass_;
00120 }
00122 double massSqr() const
00123 {
00124 return mass_ * mass_;
00125 }
00127 double mt() const
00128 {
00129 cachePolar();
00130 return p4Polar_.Mt();
00131 }
00133 double mtSqr() const
00134 {
00135 cachePolar();
00136 return p4Polar_.Mt2();
00137 }
00139 double px() const
00140 {
00141 cacheCartesian();
00142 return p4Cartesian_.Px();
00143 }
00145 double py() const
00146 {
00147 cacheCartesian();
00148 return p4Cartesian_.Py();
00149 }
00151 double pz() const
00152 {
00153 cacheCartesian();
00154 return p4Cartesian_.Pz();
00155 }
00157 double pt() const
00158 {
00159 return pt_;
00160 }
00162 double phi() const
00163 {
00164 return phi_;
00165 }
00167 double theta() const
00168 {
00169 cacheCartesian();
00170 return p4Cartesian_.Theta();
00171 }
00173 double eta() const
00174 {
00175 return eta_;
00176 }
00178 double rapidity() const
00179 {
00180 cachePolar();
00181 return p4Polar_.Rapidity();
00182 }
00184 double y() const
00185 {
00186 return rapidity();
00187 }
00189 void setP4( const LorentzVector & p4 )
00190 {
00191 p4Cartesian_ = p4;
00192 p4Polar_ = p4;
00193 pt_ = p4Polar_.pt();
00194 eta_ = p4Polar_.eta();
00195 phi_ = p4Polar_.phi();
00196 mass_ = p4Polar_.mass();
00197 cachePolarFixed_ = true;
00198 cacheCartesianFixed_ = true;
00199 }
00201 void setP4( const PolarLorentzVector & p4 )
00202 {
00203 p4Polar_ = p4;
00204 pt_ = p4Polar_.pt();
00205 eta_ = p4Polar_.eta();
00206 phi_ = p4Polar_.phi();
00207 mass_ = p4Polar_.mass();
00208 cachePolarFixed_ = true;
00209 cacheCartesianFixed_ = false;
00210 }
00212 void setMass( double m )
00213 {
00214 mass_ = m;
00215 clearCache();
00216 }
00217 void setPz( double pz )
00218 {
00219 cacheCartesian();
00220 p4Cartesian_.SetPz(pz);
00221 p4Polar_ = p4Cartesian_;
00222 pt_ = p4Polar_.pt();
00223 eta_ = p4Polar_.eta();
00224 phi_ = p4Polar_.phi();
00225 mass_ = p4Polar_.mass();
00226 }
00228 const Point & vertex() const
00229 {
00230 return vertex_;
00231 }
00233 double vx() const
00234 {
00235 return vertex_.X();
00236 }
00238 double vy() const
00239 {
00240 return vertex_.Y();
00241 }
00243 double vz() const
00244 {
00245 return vertex_.Z();
00246 }
00248 void setVertex( const Point & vertex )
00249 {
00250 vertex_ = vertex;
00251 }
00253 int pdgId() const
00254 {
00255 return pdgId_;
00256 }
00257
00258 void setPdgId( int pdgId )
00259 {
00260 pdgId_ = pdgId;
00261 }
00263 int status() const
00264 {
00265 return status_;
00266 }
00268 void setStatus( int status )
00269 {
00270 status_ = status;
00271 }
00273 static const unsigned int longLivedTag;
00275 void setLongLived()
00276 {
00277 status_ |= longLivedTag;
00278 }
00280 bool longLived() const
00281 {
00282 return status_ & longLivedTag;
00283 }
00284
00285 protected:
00287 Charge qx3_;
00289 float pt_, eta_, phi_, mass_;
00291 Point vertex_;
00293 int pdgId_;
00295 int status_;
00297 mutable PolarLorentzVector p4Polar_;
00299 mutable LorentzVector p4Cartesian_;
00301 mutable bool cachePolarFixed_, cacheCartesianFixed_;
00303 inline void cachePolar() const
00304 {
00305 if ( cachePolarFixed_ ) return;
00306 p4Polar_ = PolarLorentzVector( pt_, eta_, phi_, mass_ );
00307 cachePolarFixed_ = true;
00308 }
00310 inline void cacheCartesian() const
00311 {
00312 if ( cacheCartesianFixed_ ) return;
00313 cachePolar();
00314 p4Cartesian_ = p4Polar_;
00315 cacheCartesianFixed_ = true;
00316 }
00318 inline void clearCache() const
00319 {
00320 cachePolarFixed_ = false;
00321 cacheCartesianFixed_ = false;
00322 }
00323 };
00324
00325 #endif