CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/SimDataFormats/TrackingAnalysis/src/TrackingParticle.cc

Go to the documentation of this file.
00001 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
00002 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00003 
00004 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00005 
00006 #include <FWCore/MessageLogger/interface/MessageLogger.h>
00007 
00008 const unsigned int TrackingParticle::longLivedTag = 65536;
00009 
00010 TrackingParticle::TrackingParticle()
00011 {
00012         // No operation
00013 }
00014 
00015 TrackingParticle::TrackingParticle( const SimTrack& simtrk, const TrackingVertexRef& parentVertex )
00016 {
00017         addG4Track( simtrk );
00018         setParentVertex( parentVertex );
00019 }
00020 
00021 TrackingParticle::~TrackingParticle()
00022 {
00023 }
00024 
00025 int TrackingParticle::pdgId() const
00026 {
00027         if( genParticles_.empty() ) return g4Tracks_.at( 0 ).type();
00028         else return (*genParticles_.begin())->pdgId();
00029 }
00030 
00031 EncodedEventId TrackingParticle::eventId() const
00032 {
00033         return g4Tracks_.at( 0 ).eventId();
00034 }
00035 
00036 void TrackingParticle::addGenParticle( const reco::GenParticleRef& ref )
00037 {
00038         genParticles_.push_back( ref );
00039 }
00040 
00041 void TrackingParticle::addG4Track( const SimTrack& t )
00042 {
00043         g4Tracks_.push_back( t );
00044 }
00045 
00046 TrackingParticle::genp_iterator TrackingParticle::genParticle_begin() const
00047 {
00048         return genParticles_.begin();
00049 }
00050 
00051 TrackingParticle::genp_iterator TrackingParticle::genParticle_end() const
00052 {
00053         return genParticles_.end();
00054 }
00055 
00056 TrackingParticle::g4t_iterator TrackingParticle::g4Track_begin() const
00057 {
00058         return g4Tracks_.begin();
00059 }
00060 
00061 TrackingParticle::g4t_iterator TrackingParticle::g4Track_end() const
00062 {
00063         return g4Tracks_.end();
00064 }
00065 
00066 void TrackingParticle::setParentVertex( const TrackingVertexRef& ref )
00067 {
00068         parentVertex_=ref;
00069 }
00070 
00071 void TrackingParticle::addDecayVertex( const TrackingVertexRef& ref )
00072 {
00073         decayVertices_.push_back( ref );
00074 }
00075 
00076 void TrackingParticle::clearParentVertex()
00077 {
00078         parentVertex_=TrackingVertexRef();
00079 }
00080 
00081 void TrackingParticle::clearDecayVertices()
00082 {
00083         decayVertices_.clear();
00084 }
00085 
00086 const reco::GenParticleRefVector& TrackingParticle::genParticles() const
00087 {
00088         return genParticles_;
00089 }
00090 
00091 const std::vector<SimTrack>& TrackingParticle::g4Tracks() const
00092 {
00093         return g4Tracks_;
00094 }
00095 
00096 const TrackingVertexRef& TrackingParticle::parentVertex() const
00097 {
00098         return parentVertex_;
00099 }
00100 
00101 const TrackingVertexRefVector& TrackingParticle::decayVertices() const
00102 {
00103         return decayVertices_;
00104 }
00105 
00106 tv_iterator TrackingParticle::decayVertices_begin() const
00107 {
00108         return decayVertices_.begin();
00109 }
00110 
00111 tv_iterator TrackingParticle::decayVertices_end() const
00112 {
00113         return decayVertices_.end();
00114 }
00115 
00116 int TrackingParticle::charge() const
00117 {
00118         return g4Tracks_.at( 0 ).charge();
00119 }
00120 
00121 int TrackingParticle::threeCharge() const
00122 {
00123         return g4Tracks_.at( 0 ).charge()*3;
00124 }
00125 
00126 const TrackingParticle::LorentzVector& TrackingParticle::p4() const
00127 {
00128         return g4Tracks_.at( 0 ).momentum();
00129 }
00130 
00131 TrackingParticle::Vector TrackingParticle::momentum() const
00132 {
00133         return p4().Vect();
00134 }
00135 
00136 TrackingParticle::Vector TrackingParticle::boostToCM() const
00137 {
00138         return p4().BoostToCM();
00139 }
00140 
00141 double TrackingParticle::p() const
00142 {
00143         return p4().P();
00144 }
00145 
00146 double TrackingParticle::energy() const
00147 {
00148         return p4().E();
00149 }
00150 
00151 double TrackingParticle::et() const
00152 {
00153         return p4().Et();
00154 }
00155 
00156 double TrackingParticle::mass() const
00157 {
00158         return p4().M();
00159 }
00160 
00161 double TrackingParticle::massSqr() const
00162 {
00163         return pow( mass(), 2 );
00164 }
00165 
00166 double TrackingParticle::mt() const
00167 {
00168         return p4().Mt();
00169 }
00170 
00171 double TrackingParticle::mtSqr() const
00172 {
00173         return p4().Mt2();
00174 }
00175 
00176 double TrackingParticle::px() const
00177 {
00178         return p4().Px();
00179 }
00180 
00181 double TrackingParticle::py() const
00182 {
00183         return p4().Py();
00184 }
00185 
00186 double TrackingParticle::pz() const
00187 {
00188         return p4().Pz();
00189 }
00190 
00191 double TrackingParticle::pt() const
00192 {
00193         return p4().Pt();
00194 }
00195 
00196 double TrackingParticle::phi() const
00197 {
00198         return p4().Phi();
00199 }
00200 
00201 double TrackingParticle::theta() const
00202 {
00203         return p4().Theta();
00204 }
00205 
00206 double TrackingParticle::eta() const
00207 {
00208         return p4().Eta();
00209 }
00210 
00211 double TrackingParticle::rapidity() const
00212 {
00213         return p4().Rapidity();
00214 }
00215 
00216 double TrackingParticle::y() const
00217 {
00218         return rapidity();
00219 }
00220 
00221 TrackingParticle::Point TrackingParticle::vertex() const
00222 {
00223         return Point( vx(), vy(), vz() );
00224 }
00225 
00226 double TrackingParticle::vx() const
00227 {
00228         const TrackingVertex& r=( *parentVertex_);
00229         return r.position().X();
00230 }
00231 
00232 double TrackingParticle::vy() const
00233 {
00234         const TrackingVertex& r=( *parentVertex_);
00235         return r.position().Y();
00236 }
00237 
00238 double TrackingParticle::vz() const
00239 {
00240         const TrackingVertex& r=( *parentVertex_);
00241         return r.position().Z();
00242 }
00243 
00244 int TrackingParticle::status() const
00245 {
00246         if( genParticles_.empty() ) return -99; // Use the old invalid status flag that used to be set by TrackingTruthProducer.
00247         else return (*genParticles_.begin())->status();
00248 }
00249 
00250 bool TrackingParticle::longLived() const
00251 {
00252         return status()&longLivedTag;
00253 }
00254 
00255 int TrackingParticle::numberOfHits() const
00256 {
00257     return numberOfHits_;
00258 }
00259 
00260 int TrackingParticle::numberOfTrackerHits() const
00261 {
00262     return numberOfTrackerHits_;
00263 }
00264 
00265 int TrackingParticle::matchedHit() const
00266 {
00267         edm::LogWarning("TrackingParticle") << "The method matchedHit() has been deprecated. Use numberOfTrackerLayers() instead.";
00268         return numberOfTrackerLayers_;
00269 }
00270 
00271 int TrackingParticle::numberOfTrackerLayers() const
00272 {
00273         return numberOfTrackerLayers_;
00274 }
00275 
00276 void TrackingParticle::setNumberOfHits( int numberOfHits )
00277 {
00278     numberOfHits_=numberOfHits;
00279 }
00280 
00281 void TrackingParticle::setNumberOfTrackerHits( int numberOfTrackerHits )
00282 {
00283     numberOfTrackerHits_=numberOfTrackerHits;
00284 }
00285 
00286 void TrackingParticle::setNumberOfTrackerLayers( const int numberOfTrackerLayers )
00287 {
00288         numberOfTrackerLayers_=numberOfTrackerLayers;
00289 }
00290 
00291 std::ostream& operator<< (std::ostream& s, TrackingParticle const & tp)
00292 {
00293     s << "TP momentum, q, ID, & Event #: "
00294     << tp.p4()                      << " " << tp.charge() << " "   << tp.pdgId() << " "
00295     << tp.eventId().bunchCrossing() << "." << tp.eventId().event() << std::endl;
00296 
00297     for (TrackingParticle::genp_iterator hepT = tp.genParticle_begin(); hepT !=  tp.genParticle_end(); ++hepT)
00298     {
00299         s << " HepMC Track Momentum " << (*hepT)->momentum().rho() << std::endl;
00300     }
00301 
00302     for (TrackingParticle::g4t_iterator g4T = tp.g4Track_begin(); g4T !=  tp.g4Track_end(); ++g4T)
00303     {
00304         s << " Geant Track Momentum  " << g4T->momentum() << std::endl;
00305         s << " Geant Track ID & type " << g4T->trackId() << " " << g4T->type() << std::endl;
00306         if (g4T->type() !=  tp.pdgId())
00307         {
00308             s << " Mismatch b/t TrackingParticle and Geant types" << std::endl;
00309         }
00310     }
00311     // Loop over decay vertices
00312     s << " TP Vertex " << tp.vertex() << std::endl;
00313     s << " Source vertex: " << tp.parentVertex()->position() << std::endl;
00314     s << " " << tp.decayVertices().size() << " Decay vertices" << std::endl;
00315     for (tv_iterator iTV = tp.decayVertices_begin(); iTV != tp.decayVertices_end(); ++iTV)
00316     {
00317         s << " Decay vertices:      " << (**iTV).position() << std::endl;
00318     }
00319 
00320     return s;
00321 }