#include <SimTracker/TrackHistory/interface/TrackHistory.h>
Public Types | |
typedef std::vector< const HepMC::GenParticle * > | GenParticleTrail |
GenParticle trail type. | |
typedef std::vector< const HepMC::GenVertex * > | GenVertexTrail |
GenVertex trail type. | |
typedef std::vector < TrackingParticleRef > | SimParticleTrail |
SimParticle trail type. | |
typedef std::vector < TrackingVertexRef > | SimVertexTrail |
SimVertex trail type. | |
Public Member Functions | |
void | depth (int d) |
Set the depth of the history. | |
bool | evaluate (reco::TrackBaseRef) |
Evaluate reco::Track history using a given association. | |
bool | evaluate (TrackingParticleRef tpr) |
Evaluate track history using a TrackingParticleRef. | |
const HepMC::GenParticle * | genParticle () const |
Returns a pointer to most primitive status 1 or 2 particle. | |
const GenParticleTrail & | genParticleTrail () const |
Return all generated particle in the history. | |
const GenVertexTrail & | genVertexTrail () const |
Return all generated vertex in the history. | |
void | newEvent (const edm::Event &, const edm::EventSetup &) |
Pre-process event information (for accessing reconstruction information). | |
const TrackingParticleRef & | simParticle () const |
Return the initial tracking particle from the history. | |
const SimParticleTrail & | simParticleTrail () const |
Return all the simulated particle in the history. | |
const SimVertexTrail & | simVertexTrail () const |
Return all the simulated vertices in the history. | |
TrackHistory (const edm::ParameterSet &) | |
Constructor by pset. | |
Protected Attributes | |
int | depth_ |
GenParticleTrail | genParticleTrail_ |
GenVertexTrail | genVertexTrail_ |
SimParticleTrail | simParticleTrail_ |
SimVertexTrail | simVertexTrail_ |
Private Member Functions | |
TrackingParticleRef | match (reco::TrackBaseRef) |
void | resetTrails (TrackingParticleRef tpr) |
void | traceGenHistory (const HepMC::GenParticle *) |
bool | traceSimHistory (TrackingParticleRef, int) |
Private Attributes | |
reco::RecoToSimCollection | association_ |
bool | bestMatchByMaxValue_ |
bool | newEvent_ |
std::string | trackAssociator_ |
edm::InputTag | trackingTruth_ |
edm::InputTag | trackProducer_ |
Definition at line 27 of file TrackHistory.h.
typedef std::vector<const HepMC::GenParticle *> TrackHistory::GenParticleTrail |
typedef std::vector<const HepMC::GenVertex *> TrackHistory::GenVertexTrail |
typedef std::vector<TrackingParticleRef> TrackHistory::SimParticleTrail |
typedef std::vector<TrackingVertexRef> TrackHistory::SimVertexTrail |
TrackHistory::TrackHistory | ( | const edm::ParameterSet & | config | ) |
Constructor by pset.
Definition at line 15 of file TrackHistory.cc.
References bestMatchByMaxValue_, depth_, edm::ParameterSet::getUntrackedParameter(), trackAssociator_, trackingTruth_, and trackProducer_.
00018 { 00019 // Default depth 00020 depth_ = -1; 00021 00022 // Name of the track collection 00023 trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" ); 00024 00025 // Name of the traking pariticle collection 00026 trackingTruth_ = config.getUntrackedParameter<edm::InputTag> ( "trackingTruth" ); 00027 00028 // Track association record 00029 trackAssociator_ = config.getUntrackedParameter<std::string> ( "trackAssociator" ); 00030 00031 // Association by max. value 00032 bestMatchByMaxValue_ = config.getUntrackedParameter<bool> ( "bestMatchByMaxValue" ); 00033 }
Set the depth of the history.
Definition at line 62 of file TrackHistory.h.
References depth_.
Referenced by RecoTOA::analyze(), TruthTOA::analyze(), and TrackClassifier::TrackClassifier().
bool TrackHistory::evaluate | ( | reco::TrackBaseRef | tr | ) |
Evaluate reco::Track history using a given association.
Definition at line 150 of file TrackHistory.cc.
References evaluate(), edm::Ref< C, T, F >::isNull(), match(), and tpr.
00151 { 00152 TrackingParticleRef tpr( match(tr) ); 00153 00154 if ( !tpr.isNull() ) 00155 { 00156 evaluate(tpr); 00157 return true; 00158 } 00159 00160 return false; 00161 }
bool TrackHistory::evaluate | ( | TrackingParticleRef | tpr | ) | [inline] |
Evaluate track history using a TrackingParticleRef.
Definition at line 75 of file TrackHistory.h.
References depth_, resetTrails(), and traceSimHistory().
Referenced by RecoTOA::analyze(), VertexTHA::analyze(), TruthTOA::analyze(), evaluate(), TrackClassifier::evaluate(), and GenTrackMatcher::produce().
00076 { 00077 resetTrails(tpr); 00078 return traceSimHistory(tpr, depth_); 00079 }
const HepMC::GenParticle* TrackHistory::genParticle | ( | ) | const [inline] |
Returns a pointer to most primitive status 1 or 2 particle.
Definition at line 96 of file TrackHistory.h.
References genParticleTrail_.
Referenced by RecoTOA::analyze(), TruthTOA::analyze(), TrackClassifier::hadronFlavor(), and GenTrackMatcher::produce().
00097 { 00098 if ( genParticleTrail_.empty() ) return 0; 00099 return genParticleTrail_[genParticleTrail_.size()-1]; 00100 }
const GenParticleTrail& TrackHistory::genParticleTrail | ( | ) | const [inline] |
Return all generated particle in the history.
Definition at line 121 of file TrackHistory.h.
References genParticleTrail_.
Referenced by SimpleTHA::analyze(), and TrackClassifier::vertexInformation().
00122 { 00123 return genParticleTrail_; 00124 }
const GenVertexTrail& TrackHistory::genVertexTrail | ( | ) | const [inline] |
Return all generated vertex in the history.
Definition at line 115 of file TrackHistory.h.
References genVertexTrail_.
Referenced by SimpleTHA::analyze(), and TrackClassifier::decayProcesses().
00116 { 00117 return genVertexTrail_; 00118 }
TrackingParticleRef TrackHistory::match | ( | reco::TrackBaseRef | tr | ) | [private] |
Definition at line 164 of file TrackHistory.cc.
References association_, bestMatchByMaxValue_, edm::AssociationMap< Tag >::end(), edm::AssociationMap< Tag >::find(), i, m, edm::second(), tp, and tpr.
Referenced by evaluate().
00165 { 00166 TrackingParticleRef tpr; 00167 00168 reco::RecoToSimCollection::const_iterator pos = association_.find(tr); 00169 if (pos == association_.end()) 00170 return tpr; 00171 00172 const std::vector<std::pair<TrackingParticleRef, double> > &tp = pos->val; 00173 00174 double m = bestMatchByMaxValue_ ? -1e30 : 1e30; 00175 00176 for (std::size_t i = 0; i < tp.size(); i++) 00177 { 00178 if (bestMatchByMaxValue_ ? (tp[i].second > m) : (tp[i].second < m)) 00179 { 00180 tpr = tp[i].first; 00181 m = tp[i].second; 00182 } 00183 } 00184 00185 return tpr; 00186 }
void TrackHistory::newEvent | ( | const edm::Event & | event, | |
const edm::EventSetup & | setup | |||
) |
Pre-process event information (for accessing reconstruction information).
Definition at line 36 of file TrackHistory.cc.
References association_, edm::EventSetup::get(), trackAssociator_, trackingTruth_, and trackProducer_.
Referenced by RecoTOA::analyze(), VertexTHA::analyze(), TrackClassifier::newEvent(), and GenTrackMatcher::produce().
00039 { 00040 // Track collection 00041 edm::Handle<edm::View<reco::Track> > trackCollection; 00042 event.getByLabel(trackProducer_, trackCollection); 00043 00044 // Tracking particle information 00045 edm::Handle<TrackingParticleCollection> TPCollection; 00046 event.getByLabel(trackingTruth_, TPCollection); 00047 00048 // Get the track associator 00049 edm::ESHandle<TrackAssociatorBase> associator; 00050 setup.get<TrackAssociatorRecord>().get(trackAssociator_, associator); 00051 00052 // Get the map between recotracks and tp 00053 association_ = associator->associateRecoToSim (trackCollection, TPCollection, &event); 00054 }
void TrackHistory::resetTrails | ( | TrackingParticleRef | tpr | ) | [inline, private] |
Definition at line 139 of file TrackHistory.h.
References genParticleTrail_, genVertexTrail_, simParticleTrail_, and simVertexTrail_.
Referenced by evaluate().
00140 { 00141 // save the initial particle in the trail 00142 simParticleTrail_.clear(); 00143 simParticleTrail_.push_back(tpr); 00144 00145 // clear the remaining trails 00146 simVertexTrail_.clear(); 00147 genVertexTrail_.clear(); 00148 genParticleTrail_.clear(); 00149 }
const TrackingParticleRef& TrackHistory::simParticle | ( | ) | const [inline] |
Return the initial tracking particle from the history.
Definition at line 90 of file TrackHistory.h.
References simParticleTrail_.
Referenced by TrackClassifier::reconstructionInformation(), and TrackClassifier::simulationInformation().
00091 { 00092 return simParticleTrail_[0]; 00093 }
const SimParticleTrail& TrackHistory::simParticleTrail | ( | ) | const [inline] |
Return all the simulated particle in the history.
Definition at line 109 of file TrackHistory.h.
References simParticleTrail_.
Referenced by SimpleTHA::analyze(), TrackClassifier::conversionInteraction(), TrackClassifier::qualityInformation(), and TrackClassifier::vertexInformation().
00110 { 00111 return simParticleTrail_; 00112 }
const SimVertexTrail& TrackHistory::simVertexTrail | ( | ) | const [inline] |
Return all the simulated vertices in the history.
Definition at line 103 of file TrackHistory.h.
References simVertexTrail_.
Referenced by VertexTHA::analyze(), and SimpleTHA::analyze().
00104 { 00105 return simVertexTrail_; 00106 }
void TrackHistory::traceGenHistory | ( | const HepMC::GenParticle * | gpp | ) | [private] |
Definition at line 57 of file TrackHistory.cc.
References funct::abs(), depth_, genParticleTrail_, and genVertexTrail_.
Referenced by traceSimHistory().
00058 { 00059 // Third stop criteria: status abs(depth_) particles after the hadronization. 00060 // The after hadronization is done by detecting the pdg_id pythia code from 88 to 99 00061 if ( gpp->status() <= abs(depth_) && (gpp->pdg_id() < 88 || gpp->pdg_id() > 99) ) 00062 { 00063 genParticleTrail_.push_back(gpp); 00064 // Get the producer vertex. 00065 HepMC::GenVertex * vertex = gpp->production_vertex(); 00066 // Verify if has a vertex associated 00067 if ( vertex ) 00068 { 00069 genVertexTrail_.push_back(vertex); 00070 if ( vertex->particles_in_size() ) // Verify if the vertex has incoming particles 00071 traceGenHistory( *(vertex->particles_in_const_begin()) ); 00072 } 00073 } 00074 }
bool TrackHistory::traceSimHistory | ( | TrackingParticleRef | tpr, | |
int | depth | |||
) | [private] |
Definition at line 77 of file TrackHistory.cc.
References depth_, lat::endl(), find(), edm::Ref< C, T, F >::isNonnull(), itd, its, LogDebug, simParticleTrail_, simVertexTrail_, and traceGenHistory().
Referenced by evaluate().
00078 { 00079 // first stop condition: if the required depth is reached 00080 if ( depth == depth_ && depth_ >= 0 ) return true; 00081 00082 // sencond stop condition: if a gen particle is associated to the TP 00083 if ( !tpr->genParticle().empty() ) 00084 { 00085 LogDebug("TrackHistory") << "Particle " << tpr->pdgId() << " has a GenParicle image." << std::endl; 00086 traceGenHistory(&(**(tpr->genParticle_begin()))); 00087 return true; 00088 } 00089 00090 // get a reference to the TP's parent vertex 00091 TrackingVertexRef parentVertex = tpr->parentVertex(); 00092 00093 // verify if the parent vertex exists 00094 if ( parentVertex.isNonnull() ) 00095 { 00096 // save the vertex in the trail 00097 simVertexTrail_.push_back(parentVertex); 00098 00099 if ( !parentVertex->sourceTracks().empty() ) 00100 { 00101 LogDebug("TrackHistory") << "No GenParticle image for " << tpr->pdgId() << " moving on to the parent particle." << std::endl; 00102 00103 // select the original source in case of combined vertices 00104 bool flag = false; 00105 TrackingVertex::tp_iterator itd, its; 00106 00107 for (its = parentVertex->sourceTracks_begin(); its != parentVertex->sourceTracks_end(); its++) 00108 { 00109 for (itd = parentVertex->daughterTracks_begin(); itd != parentVertex->daughterTracks_end(); itd++) 00110 if (itd != its) 00111 { 00112 flag = true; 00113 break; 00114 } 00115 if (flag) 00116 break; 00117 } 00118 00119 // verify if the new particle is not in the trail (looping partiles) 00120 if ( 00121 std::find( 00122 simParticleTrail_.begin(), 00123 simParticleTrail_.end(), 00124 *its 00125 ) != simParticleTrail_.end() 00126 ) 00127 { 00128 LogDebug("TrackHistory") << "WARNING: Looping track found." << std::endl; 00129 return false; 00130 } 00131 00132 // save particle in the trail 00133 simParticleTrail_.push_back(*its); 00134 return traceSimHistory (*its, --depth); 00135 } 00136 else 00137 { 00138 LogDebug("TrackHistory") << "WARNING: Source track for tracking vertex cannot be found." << std::endl; 00139 } 00140 } 00141 else 00142 { 00143 LogDebug("TrackHistory") << " WARNING: Parent vertex for tracking particle cannot be found."; 00144 } 00145 00146 return false; 00147 }
bool TrackHistory::bestMatchByMaxValue_ [private] |
int TrackHistory::depth_ [protected] |
Definition at line 128 of file TrackHistory.h.
Referenced by depth(), evaluate(), traceGenHistory(), traceSimHistory(), and TrackHistory().
GenParticleTrail TrackHistory::genParticleTrail_ [protected] |
Definition at line 131 of file TrackHistory.h.
Referenced by genParticle(), genParticleTrail(), resetTrails(), and traceGenHistory().
GenVertexTrail TrackHistory::genVertexTrail_ [protected] |
Definition at line 130 of file TrackHistory.h.
Referenced by genVertexTrail(), resetTrails(), and traceGenHistory().
bool TrackHistory::newEvent_ [private] |
Definition at line 137 of file TrackHistory.h.
SimParticleTrail TrackHistory::simParticleTrail_ [protected] |
Definition at line 133 of file TrackHistory.h.
Referenced by resetTrails(), simParticle(), simParticleTrail(), and traceSimHistory().
SimVertexTrail TrackHistory::simVertexTrail_ [protected] |
Definition at line 132 of file TrackHistory.h.
Referenced by resetTrails(), simVertexTrail(), and traceSimHistory().
std::string TrackHistory::trackAssociator_ [private] |
edm::InputTag TrackHistory::trackingTruth_ [private] |
edm::InputTag TrackHistory::trackProducer_ [private] |