CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
VertexHistory Class Reference

This class traces the simulated and generated history of a given track. More...

#include <VertexHistory.h>

Inheritance diagram for VertexHistory:
HistoryBase

Public Member Functions

bool evaluate (TrackingVertexRef tvr)
 Evaluate track history using a TrackingParticleRef. More...
 
bool evaluate (reco::VertexBaseRef)
 Evaluate reco::Vertex history using a given association. More...
 
void newEvent (const edm::Event &, const edm::EventSetup &)
 Pre-process event information (for accessing reconstruction information) More...
 
double quality () const
 Return the quality of the match. More...
 
const reco::VertexBaseRefrecoVertex () const
 Return a reference to the reconstructed track. More...
 
 VertexHistory (const edm::ParameterSet &)
 Constructor by pset. More...
 
- Public Member Functions inherited from HistoryBase
void depth (int d)
 Set the depth of the history. More...
 
const HepMC::GenParticle * genParticle () const
 Returns a pointer to most primitive status 1 or 2 particle. More...
 
GenParticleTrail const & genParticleTrail () const
 Return all generated particle in the history. More...
 
GenVertexTrail const & genVertexTrail () const
 Return all generated vertex in the history. More...
 
 HistoryBase ()
 
const TrackingParticleRefsimParticle () const
 Return the initial tracking particle from the history. More...
 
SimParticleTrail const & simParticleTrail () const
 Return all the simulated particle in the history. More...
 
const TrackingVertexRefsimVertex () const
 Return the initial tracking vertex from the history. More...
 
SimVertexTrail const & simVertexTrail () const
 Return all the simulated vertices in the history. More...
 

Private Attributes

bool bestMatchByMaxValue_
 
bool enableRecoToSim_
 
bool enableSimToReco_
 
double quality_
 
reco::VertexRecoToSimCollection recoToSim_
 
reco::VertexBaseRef recovertex_
 
reco::VertexSimToRecoCollection simToReco_
 
std::string trackAssociator_
 
edm::InputTag trackingTruth_
 
edm::InputTag trackProducer_
 
std::string vertexAssociator_
 
edm::InputTag vertexProducer_
 

Additional Inherited Members

- Public Types inherited from HistoryBase
typedef std::vector< const
HepMC::GenParticle * > 
GenParticleTrail
 GenParticle trail type. More...
 
typedef std::vector< const
HepMC::GenVertex * > 
GenVertexTrail
 GenVertex trail type. More...
 
typedef std::set< const
HepMC::GenVertex * > 
GenVertexTrailHelper
 GenVertex trail helper type. More...
 
typedef std::vector
< TrackingParticleRef
SimParticleTrail
 SimParticle trail type. More...
 
typedef std::vector
< TrackingVertexRef
SimVertexTrail
 SimVertex trail type. More...
 
- Protected Member Functions inherited from HistoryBase
bool evaluate (TrackingParticleRef tpr)
 Evaluate track history using a TrackingParticleRef. More...
 
bool evaluate (TrackingVertexRef tvr)
 Evaluate track history using a TrackingParticleRef. More...
 
- Protected Attributes inherited from HistoryBase
GenParticleTrail genParticleTrail_
 
GenVertexTrail genVertexTrail_
 
GenVertexTrailHelper genVertexTrailHelper_
 
SimParticleTrail simParticleTrail_
 
SimVertexTrail simVertexTrail_
 

Detailed Description

This class traces the simulated and generated history of a given track.

Definition at line 17 of file VertexHistory.h.

Constructor & Destructor Documentation

VertexHistory::VertexHistory ( const edm::ParameterSet config)

Constructor by pset.

Definition at line 6 of file VertexHistory.cc.

References bestMatchByMaxValue_, enableRecoToSim_, enableSimToReco_, edm::ParameterSet::getUntrackedParameter(), quality_, trackAssociator_, trackingTruth_, trackProducer_, vertexAssociator_, and vertexProducer_.

8  : HistoryBase()
9 {
10  // Name of the track collection
11  trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
12 
13  // Name of the track collection
14  vertexProducer_ = config.getUntrackedParameter<edm::InputTag> ( "vertexProducer" );
15 
16  // Name of the traking pariticle collection
17  trackingTruth_ = config.getUntrackedParameter<edm::InputTag> ( "trackingTruth" );
18 
19  // Track association record
20  trackAssociator_ = config.getUntrackedParameter<std::string> ( "trackAssociator" );
21 
22  // Track association record
23  vertexAssociator_ = config.getUntrackedParameter<std::string> ( "vertexAssociator" );
24 
25  // Association by max. value
26  bestMatchByMaxValue_ = config.getUntrackedParameter<bool> ( "bestMatchByMaxValue" );
27 
28  // Enable RecoToSim association
29  enableRecoToSim_ = config.getUntrackedParameter<bool> ( "enableRecoToSim" );
30 
31  // Enable SimToReco association
32  enableSimToReco_ = config.getUntrackedParameter<bool> ( "enableSimToReco" );
33 
34  quality_ = 0.;
35 }
T getUntrackedParameter(std::string const &, T const &) const
std::string trackAssociator_
Definition: VertexHistory.h:88
bool bestMatchByMaxValue_
Definition: VertexHistory.h:76
edm::InputTag trackProducer_
Definition: VertexHistory.h:82
edm::InputTag trackingTruth_
Definition: VertexHistory.h:86
bool enableRecoToSim_
Definition: VertexHistory.h:78
std::string vertexAssociator_
Definition: VertexHistory.h:90
bool enableSimToReco_
Definition: VertexHistory.h:78
edm::InputTag vertexProducer_
Definition: VertexHistory.h:84

Member Function Documentation

bool VertexHistory::evaluate ( TrackingVertexRef  tvr)
inline

Evaluate track history using a TrackingParticleRef.

Definition at line 40 of file VertexHistory.h.

References bestMatchByMaxValue_, enableSimToReco_, HistoryBase::evaluate(), match(), quality_, recovertex_, query::result, and simToReco_.

Referenced by VertexClassifier::evaluate().

41  {
42  if ( enableSimToReco_ )
43  {
44 
45  std::pair<reco::VertexBaseRef, double> result = match(tvr, simToReco_, bestMatchByMaxValue_);
46  recovertex_ = result.first;
47  quality_ = result.second;
48  }
49  return HistoryBase::evaluate(tvr);
50  }
bool bestMatchByMaxValue_
Definition: VertexHistory.h:76
reco::VertexSimToRecoCollection simToReco_
Definition: VertexHistory.h:96
reco::VertexBaseRef recovertex_
Definition: VertexHistory.h:92
tuple result
Definition: query.py:137
bool evaluate(TrackingParticleRef tpr)
Evaluate track history using a TrackingParticleRef.
Definition: HistoryBase.h:114
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
bool enableSimToReco_
Definition: VertexHistory.h:78
bool VertexHistory::evaluate ( reco::VertexBaseRef  tv)

Evaluate reco::Vertex history using a given association.

Definition at line 93 of file VertexHistory.cc.

References bestMatchByMaxValue_, enableRecoToSim_, HistoryBase::evaluate(), edm::Ref< C, T, F >::isNull(), match(), quality_, recoToSim_, recovertex_, and query::result.

94 {
95 
96  if ( !enableRecoToSim_ ) return false;
97 
98  std::pair<TrackingVertexRef, double> result = match(tv, recoToSim_, bestMatchByMaxValue_);
99 
100  TrackingVertexRef tvr( result.first );
101  quality_ = result.second;
102 
103  if ( !tvr.isNull() )
104  {
106 
107  recovertex_ = tv;
108 
109  return true;
110  }
111 
112  return false;
113 }
reco::VertexRecoToSimCollection recoToSim_
Definition: VertexHistory.h:94
bool bestMatchByMaxValue_
Definition: VertexHistory.h:76
reco::VertexBaseRef recovertex_
Definition: VertexHistory.h:92
tuple result
Definition: query.py:137
bool enableRecoToSim_
Definition: VertexHistory.h:78
bool evaluate(TrackingParticleRef tpr)
Evaluate track history using a TrackingParticleRef.
Definition: HistoryBase.h:114
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
void VertexHistory::newEvent ( const edm::Event event,
const edm::EventSetup setup 
)

Pre-process event information (for accessing reconstruction information)

Definition at line 38 of file VertexHistory.cc.

References enableRecoToSim_, enableSimToReco_, edm::EventSetup::get(), recoToSim_, simToReco_, ecalTPGAnalyzer_cfg::TPCollection, trackAssociator_, trackingTruth_, trackProducer_, vertexAssociator_, GoodVertex_cfg::vertexCollection, and vertexProducer_.

Referenced by VertexClassifier::newEvent().

41 {
43  {
44 
45  // Track collection
46  edm::Handle<edm::View<reco::Track> > trackCollection;
47  event.getByLabel(trackProducer_, trackCollection);
48 
49  // Tracking particle information
51  event.getByLabel(trackingTruth_, TPCollection);
52 
53  // Get the track associator
54  edm::ESHandle<TrackAssociatorBase> trackAssociator;
55  setup.get<TrackAssociatorRecord>().get(trackAssociator_, trackAssociator);
56 
57  // Vertex collection
59  event.getByLabel(vertexProducer_, vertexCollection);
60 
61  // Tracking particle information
63  event.getByLabel(trackingTruth_, TVCollection);
64 
65  // Get the track associator
66  edm::ESHandle<VertexAssociatorBase> vertexAssociator;
67  setup.get<VertexAssociatorRecord>().get(vertexAssociator_, vertexAssociator);
68 
69  if ( enableRecoToSim_ )
70  {
71  // Get the map between recovertex -> simvertex
73  trackRecoToSim = trackAssociator->associateRecoToSim(trackCollection, TPCollection, &event);
74 
75  // Calculate the map between recovertex -> simvertex
76  recoToSim_ = vertexAssociator->associateRecoToSim(vertexCollection, TVCollection, event, trackRecoToSim);
77  }
78 
79  if ( enableSimToReco_ )
80  {
81  // Get the map between recovertex <- simvertex
83  trackSimToReco = trackAssociator->associateSimToReco (trackCollection, TPCollection, &event);
84 
85  // Calculate the map between recovertex <- simvertex
86  simToReco_ = vertexAssociator->associateSimToReco(vertexCollection, TVCollection, event, trackSimToReco);
87  }
88 
89  }
90 }
reco::VertexRecoToSimCollection recoToSim_
Definition: VertexHistory.h:94
std::string trackAssociator_
Definition: VertexHistory.h:88
edm::InputTag trackProducer_
Definition: VertexHistory.h:82
reco::VertexSimToRecoCollection simToReco_
Definition: VertexHistory.h:96
tuple vertexCollection
edm::InputTag trackingTruth_
Definition: VertexHistory.h:86
bool enableRecoToSim_
Definition: VertexHistory.h:78
const T & get() const
Definition: EventSetup.h:55
std::string vertexAssociator_
Definition: VertexHistory.h:90
bool enableSimToReco_
Definition: VertexHistory.h:78
edm::InputTag vertexProducer_
Definition: VertexHistory.h:84
double VertexHistory::quality ( void  ) const
inline

Return the quality of the match.

Definition at line 68 of file VertexHistory.h.

References quality_.

Referenced by SVTagInfoValidationAnalyzer::analyze(), and recoBSVTagInfoValidationAnalyzer::analyze().

69  {
70  return quality_;
71  }
const reco::VertexBaseRef& VertexHistory::recoVertex ( ) const
inline

Return a reference to the reconstructed track.

Definition at line 62 of file VertexHistory.h.

References recovertex_.

Referenced by SVTagInfoValidationAnalyzer::analyze(), recoBSVTagInfoValidationAnalyzer::analyze(), and VertexClassifier::evaluate().

63  {
64  return recovertex_;
65  }
reco::VertexBaseRef recovertex_
Definition: VertexHistory.h:92

Member Data Documentation

bool VertexHistory::bestMatchByMaxValue_
private

Definition at line 76 of file VertexHistory.h.

Referenced by evaluate(), and VertexHistory().

bool VertexHistory::enableRecoToSim_
private

Definition at line 78 of file VertexHistory.h.

Referenced by evaluate(), newEvent(), and VertexHistory().

bool VertexHistory::enableSimToReco_
private

Definition at line 78 of file VertexHistory.h.

Referenced by evaluate(), newEvent(), and VertexHistory().

double VertexHistory::quality_
private

Definition at line 80 of file VertexHistory.h.

Referenced by evaluate(), quality(), and VertexHistory().

reco::VertexRecoToSimCollection VertexHistory::recoToSim_
private

Definition at line 94 of file VertexHistory.h.

Referenced by evaluate(), and newEvent().

reco::VertexBaseRef VertexHistory::recovertex_
private

Definition at line 92 of file VertexHistory.h.

Referenced by evaluate(), and recoVertex().

reco::VertexSimToRecoCollection VertexHistory::simToReco_
private

Definition at line 96 of file VertexHistory.h.

Referenced by evaluate(), and newEvent().

std::string VertexHistory::trackAssociator_
private

Definition at line 88 of file VertexHistory.h.

Referenced by newEvent(), and VertexHistory().

edm::InputTag VertexHistory::trackingTruth_
private

Definition at line 86 of file VertexHistory.h.

Referenced by newEvent(), and VertexHistory().

edm::InputTag VertexHistory::trackProducer_
private

Definition at line 82 of file VertexHistory.h.

Referenced by newEvent(), and VertexHistory().

std::string VertexHistory::vertexAssociator_
private

Definition at line 90 of file VertexHistory.h.

Referenced by newEvent(), and VertexHistory().

edm::InputTag VertexHistory::vertexProducer_
private

Definition at line 84 of file VertexHistory.h.

Referenced by newEvent(), and VertexHistory().