CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackClassifier.h
Go to the documentation of this file.
1 
2 #ifndef TrackClassifier_h
3 #define TrackClassifier_h
4 
8 
11 
13 
15 
20 
25 
26 class TrackTopology;
27 
30 public:
33 
36 
38  void newEvent(edm::Event const &, edm::EventSetup const &);
39 
42 
45 
48 
50  TrackHistory const &history() const { return tracer_; }
51 
53  TrackQuality const &quality() const { return quality_; }
54 
55 private:
58 
59  double badPull_;
62  unsigned int numberOfInnerLayers_;
63  unsigned int minTrackerSimHits_;
64 
66 
68 
70 
73 
75 
78 
81 
83 
86 
90 
92  void simulationInformation();
93 
96 
98  void hadronFlavor();
99 
101  void processesAtGenerator();
102 
104  void processesAtSimulation();
105 
107  void vertexInformation();
108 
111  GeneratedPrimaryVertex(double x1, double y1, double z1) : x(x1), y(y1), z(z1), ptsq(0), nGenTrk(0) {}
112 
113  bool operator<(GeneratedPrimaryVertex const &reference) const { return ptsq < reference.ptsq; }
114 
115  double x, y, z;
116  double ptsq;
117  int nGenTrk;
118 
119  HepMC::FourVector ptot;
120 
121  std::vector<int> finalstateParticles;
122  std::vector<int> simTrackIndex;
123  std::vector<int> genVertex;
124  };
125 
126  std::vector<GeneratedPrimaryVertex> genpvs_;
127 
128  // Auxiliary function to get the generated primary vertex
130  bool isCharged(const HepMC::GenParticle *);
131  void genPrimaryVertices();
132 };
133 
134 #endif
const edm::InputTag beamSpotLabel_
const edm::InputTag hepMCLabel_
TrackQuality const & quality() const
Returns a reference to the track quality used in the classification.
GeneratedPrimaryVertex(double x1, double y1, double z1)
TrackClassifier const & evaluate(reco::TrackRef const &track)
Classify the RecoTrack in categories.
const G4toCMSLegacyProcTypeMap g4toCMSProcMap_
edm::Handle< edm::HepMCProduct > mcInformation_
const TrackerTopology * tTopo_
unsigned int numberOfInnerLayers_
TrackHistory const & history() const
Returns a reference to the track history used in the classification.
unsigned int minTrackerSimHits_
std::vector< GeneratedPrimaryVertex > genpvs_
bool operator<(GeneratedPrimaryVertex const &reference) const
void processesAtSimulation()
Get information about conversion and other interactions.
void simulationInformation()
Get all the information related to the simulation details.
edm::ESGetToken< ParticleDataTable, PDTRecord > particleDataTableToken_
void qualityInformation(reco::TrackBaseRef const &)
Classify all the tracks by their reconstruction quality.
void reconstructionInformation(reco::TrackBaseRef const &)
TrackClassifier const & evaluate(reco::TrackBaseRef const &)
Classify the RecoTrack in categories.
TrackClassifier(edm::ParameterSet const &, edm::ConsumesCollector &&)
Constructor by ParameterSet.
bool isCharged(const HepMC::GenParticle *)
edm::ESHandle< TransientTrackBuilder > transientTrackBuilder_
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackBuilderToken_
Auxiliary class holding simulated primary vertices.
This class traces the simulated and generated history of a given track.
Definition: TrackHistory.h:17
void hadronFlavor()
Get hadron flavor of the initial hadron.
void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
edm::Handle< reco::BeamSpot > beamSpot_
TrackQuality quality_
Get track history and classify it in function of their .
This class analyses the reconstruction quality for a given track.
Definition: TrackQuality.h:29
bool isFinalstateParticle(const HepMC::GenParticle *)
void vertexInformation()
Get geometrical information about the vertices.
void processesAtGenerator()
Get all the information related to decay process.
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoHandToken_
edm::ESHandle< MagneticField > magneticField_
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
double longLivedDecayLength_
double vertexClusteringSqDistance_
TrackCategories Categories
Type to the associate category.
edm::ESHandle< ParticleDataTable > particleDataTable_
TrackHistory tracer_