CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackClassifier.cc
Go to the documentation of this file.
1 
2 #include <cmath>
3 #include <cstdlib>
4 #include <iostream>
5 
6 #include "HepPDT/ParticleID.hh"
7 
9 
10 #define update(a, b) \
11  do { \
12  (a) = (a) | (b); \
13  } while (0)
14 
16  : TrackCategories(),
17  hepMCLabel_(config.getUntrackedParameter<edm::InputTag>("hepMC")),
18  beamSpotLabel_(config.getUntrackedParameter<edm::InputTag>("beamSpot")),
19  tracer_(config, std::move(collector)),
20  quality_(config, collector),
21  magneticFieldToken_(collector.esConsumes<MagneticField, IdealMagneticFieldRecord>()),
22  particleDataTableToken_(collector.esConsumes()),
23  transientTrackBuilderToken_(collector.esConsumes<TransientTrackBuilder, TransientTrackRecord>()),
24  tTopoHandToken_(collector.esConsumes<TrackerTopology, TrackerTopologyRcd>()) {
25  collector.consumes<edm::HepMCProduct>(hepMCLabel_);
26  collector.consumes<reco::BeamSpot>(beamSpotLabel_);
27 
28  // Set the history depth after hadronization
29  tracer_.depth(-2);
30 
31  // Set the maximum d0pull for the bad category
32  badPull_ = config.getUntrackedParameter<double>("badPull");
33 
34  // Set the minimum decay length for detecting long decays
35  longLivedDecayLength_ = config.getUntrackedParameter<double>("longLivedDecayLength");
36 
37  // Set the distance for clustering vertices
38  float vertexClusteringDistance = config.getUntrackedParameter<double>("vertexClusteringDistance");
39  vertexClusteringSqDistance_ = vertexClusteringDistance * vertexClusteringDistance;
40 
41  // Set the number of innermost layers to check for bad hits
42  numberOfInnerLayers_ = config.getUntrackedParameter<unsigned int>("numberOfInnerLayers");
43 
44  // Set the minimum number of simhits in the tracker
45  minTrackerSimHits_ = config.getUntrackedParameter<unsigned int>("minTrackerSimHits");
46 }
47 
49  // Get the new event information for the tracer
50  tracer_.newEvent(event, setup);
51 
52  // Get the new event information for the track quality analyser
53  quality_.newEvent(event, setup);
54 
55  // Get hepmc of the event
56  event.getByLabel(hepMCLabel_, mcInformation_);
57 
58  // Magnetic field
60 
61  // Get the partivle data table
63 
64  // get the beam spot
65  event.getByLabel(beamSpotLabel_, beamSpot_);
66 
67  // Transient track builder
69 
70  // Create the list of primary vertices associated to the event
72 
73  // Retrieve tracker topology from geometry
74  tTopo_ = &setup.getData(tTopoHandToken_);
75 }
76 
78  // Initializing the category vector
79  reset();
80 
81  // Associate and evaluate the track history (check for fakes)
82  if (tracer_.evaluate(track)) {
83  // Classify all the tracks by their association and reconstruction
84  // information
86 
87  // Get all the information related to the simulation details
89 
90  // Analyse the track reconstruction quality
91  qualityInformation(track);
92 
93  // Get hadron flavor of the initial hadron
94  hadronFlavor();
95 
96  // Get all the information related to decay process
98 
99  // Get information about conversion and other interactions
101 
102  // Get geometrical information about the vertices
104 
105  // Check for unkown classification
106  unknownTrack();
107  } else
108  flags_[Fake] = true;
109 
110  return *this;
111 }
112 
114  // Initializing the category vector
115  reset();
116 
117  // Trace the history for the given TP
118  tracer_.evaluate(track);
119 
120  // Collect the associated reco track
122 
123  // If there is a reco truck then evaluate the simulated history
124  if (recotrack.isNonnull()) {
125  flags_[Reconstructed] = true;
126  // Classify all the tracks by their association and reconstruction
127  // information
128  reconstructionInformation(recotrack);
129  // Analyse the track reconstruction quality
130  qualityInformation(recotrack);
131  } else
132  flags_[Reconstructed] = false;
133 
134  // Get all the information related to the simulation details
136 
137  // Get hadron flavor of the initial hadron
138  hadronFlavor();
139 
140  // Get all the information related to decay process
142 
143  // Get information about conversion and other interactions
145 
146  // Get geometrical information about the vertices
148 
149  // Check for unkown classification
150  unknownTrack();
151 
152  return *this;
153 }
154 
157 
158  // Compute tracking particle parameters at point of closest approach to the
159  // beamline
160 
161  const SimTrack *assocTrack = &(*tpr->g4Track_begin());
162 
163  FreeTrajectoryState ftsAtProduction(
164  GlobalPoint(tpr->vertex().x(), tpr->vertex().y(), tpr->vertex().z()),
165  GlobalVector(assocTrack->momentum().x(), assocTrack->momentum().y(), assocTrack->momentum().z()),
166  TrackCharge(track->charge()),
168 
169  try {
170  TSCPBuilderNoMaterial tscpBuilder;
171  TrajectoryStateClosestToPoint tsAtClosestApproach =
172  tscpBuilder(ftsAtProduction, GlobalPoint(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()));
173 
174  GlobalVector v =
175  tsAtClosestApproach.theState().position() - GlobalPoint(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0());
176  GlobalVector p = tsAtClosestApproach.theState().momentum();
177 
178  // Simulated dxy
179  double dxySim = -v.x() * sin(p.phi()) + v.y() * cos(p.phi());
180 
181  // Simulated dz
182  double dzSim = v.z() - (v.x() * p.x() + v.y() * p.y()) * p.z() / p.perp2();
183 
184  // Calculate the dxy pull
185  double dxyPull =
186  std::abs(track->dxy(reco::TrackBase::Point(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0())) - dxySim) /
187  track->dxyError();
188 
189  // Calculate the dx pull
190  double dzPull =
191  std::abs(track->dz(reco::TrackBase::Point(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0())) - dzSim) /
192  track->dzError();
193 
194  // Return true if d0Pull > badD0Pull sigmas
195  flags_[Bad] = (dxyPull > badPull_ || dzPull > badPull_);
196 
197  } catch (cms::Exception const &) {
198  flags_[Bad] = true;
199  }
200 }
201 
203  // Get the event id for the initial TP.
204  EncodedEventId eventId = tracer_.simParticle()->eventId();
205  // Check for signal events
206  flags_[SignalEvent] = !eventId.bunchCrossing() && !eventId.event();
207  // Check for muons
208  flags_[Muon] = (abs(tracer_.simParticle()->pdgId()) == 13);
209  // Check for the number of psimhit in tracker
210  flags_[TrackerSimHits] = tracer_.simParticle()->numberOfTrackerLayers() >= (int)minTrackerSimHits_;
211 }
212 
214  // run the hit-by-hit reconstruction quality analysis
216 
218 
219  // check the innermost layers for bad hits
220  for (unsigned int i = 0; i < maxLayers; i++) {
222 
223  // check all hits in that layer
224  for (unsigned int j = 0; j < layer.hits.size(); j++) {
225  const TrackQuality::Layer::Hit &hit = layer.hits[j];
226 
227  // In those cases the bad hit was used by track reconstruction
229  flags_[BadInnerHits] = true;
230  else if (hit.state == TrackQuality::Layer::Shared)
231  flags_[SharedInnerHits] = true;
232  }
233  }
234 }
235 
237  // Get the initial hadron from the recoGenParticleTrail
238  const reco::GenParticle *particle = tracer_.recoGenParticle();
239 
240  // Check for the initial hadron
241  if (particle) {
242  HepPDT::ParticleID pid(particle->pdgId());
243  flags_[Bottom] = pid.hasBottom();
244  flags_[Charm] = pid.hasCharm();
245  flags_[Light] = !pid.hasCharm() && !pid.hasBottom();
246  }
247 }
248 
250  // pdgid of the "in" particle to the production vertex
251  int pdgid = 0;
252 
253  // Get the generated particles from track history (reco::GenParticle in the
254  // recoGenParticleTrail)
255  TrackHistory::RecoGenParticleTrail const &recoGenParticleTrail = tracer_.recoGenParticleTrail();
256 
257  // Loop over the generated particles (reco::GenParticle in the
258  // recoGenParticleTrail)
259  for (TrackHistory::RecoGenParticleTrail::const_iterator iparticle = recoGenParticleTrail.begin();
260  iparticle != recoGenParticleTrail.end();
261  ++iparticle) {
262  pdgid = std::abs((*iparticle)->pdgId());
263  // Get particle type
264  HepPDT::ParticleID particleID(pdgid);
265 
266  // Check if the particle type is valid one
267  if (particleID.isValid()) {
268  // Get particle data
269  ParticleData const *particleData = particleDataTable_->particle(particleID);
270  // Check if the particle exist in the table
271  if (particleData) {
272  // Check if their life time is bigger than longLivedDecayLength_
273  if (particleData->lifetime() > longLivedDecayLength_)
274  update(flags_[LongLivedDecay], true);
275  // Check for B and C weak decays
276  update(flags_[BWeakDecay], particleID.hasBottom());
277  update(flags_[CWeakDecay], particleID.hasCharm());
278  // Check for B and C pure leptonic decay
279  std::set<int> daughterIds;
280  size_t ndau = (*iparticle)->numberOfDaughters();
281  for (size_t i = 0; i < ndau; ++i) {
282  daughterIds.insert((*iparticle)->daughter(i)->pdgId());
283  }
284  update(flags_[FromBWeakDecayMuon], particleID.hasBottom() && (daughterIds.find(13) != daughterIds.end()));
285  update(flags_[FromCWeakDecayMuon], particleID.hasCharm() && (daughterIds.find(13) != daughterIds.end()));
286  }
287  // Check Tau, Ks and Lambda decay
288  update(flags_[ChargePionDecay], pdgid == 211);
289  update(flags_[ChargeKaonDecay], pdgid == 321);
290  update(flags_[TauDecay], pdgid == 15);
291  update(flags_[KsDecay], pdgid == 310);
292  update(flags_[LambdaDecay], pdgid == 3122);
293  update(flags_[JpsiDecay], pdgid == 443);
294  update(flags_[XiDecay], pdgid == 3312);
295  update(flags_[SigmaPlusDecay], pdgid == 3222);
296  update(flags_[SigmaMinusDecay], pdgid == 3112);
297  }
298  }
299  // Decays in flight
302  update(flags_[DecayOnFlightMuon], (flags_[FromChargePionMuon] || flags_[FromChargeKaonMuon]));
303 }
304 
306  TrackHistory::SimParticleTrail const &simParticleTrail = tracer_.simParticleTrail();
307 
308  // Loop over the simulated particles
309  for (TrackHistory::SimParticleTrail::const_iterator iparticle = simParticleTrail.begin();
310  iparticle != simParticleTrail.end();
311  ++iparticle) {
312  // pdgid of the real source parent vertex
313  int pdgid = 0;
314 
315  // Get a reference to the TP's parent vertex
316  TrackingVertexRef const &parentVertex = (*iparticle)->parentVertex();
317 
318  // Look for the original source track
319  if (parentVertex.isNonnull()) {
320  // select the original source in case of combined vertices
321  bool flag = false;
323 
324  for (its = parentVertex->sourceTracks_begin(); its != parentVertex->sourceTracks_end(); ++its) {
325  for (itd = parentVertex->daughterTracks_begin(); itd != parentVertex->daughterTracks_end(); ++itd)
326  if (itd != its) {
327  flag = true;
328  break;
329  }
330  if (flag)
331  break;
332  }
333 
334  // Collect the pdgid of the original source track
335  if (its != parentVertex->sourceTracks_end())
336  pdgid = std::abs((*its)->pdgId());
337  else
338  pdgid = 0;
339  }
340 
341  unsigned int processG4 = 0;
342 
343  // Check existence of SimVerteces assigned
344  if (parentVertex->nG4Vertices() > 0) {
345  processG4 = (*(parentVertex->g4Vertices_begin())).processType();
346  }
347 
348  unsigned int process = g4toCMSProcMap_.processId(processG4);
349 
350  // Flagging all the different processes
351  update(flags_[KnownProcess], process != CMS::Undefined && process != CMS::Unknown && process != CMS::Primary);
352 
357  update(flags_[DecayProcess], process == CMS::Decay);
360  update(flags_[EIoniProcess], process == CMS::EIoni);
361  update(flags_[HIoniProcess], process == CMS::HIoni);
362  update(flags_[MuIoniProcess], process == CMS::MuIoni);
363  update(flags_[PhotonProcess], process == CMS::Photon);
366  update(flags_[EBremProcess], process == CMS::EBrem);
368  update(flags_[MuBremProcess], process == CMS::MuBrem);
369  update(flags_[MuNuclProcess], process == CMS::MuNucl);
370 
371  // Get particle type
372  HepPDT::ParticleID particleID(pdgid);
373 
374  // Check if the particle type is valid one
375  if (particleID.isValid()) {
376  // Get particle data
377  ParticleData const *particleData = particleDataTable_->particle(particleID);
378  // Special treatment for decays
379  if (process == CMS::Decay) {
380  // Check if the particle exist in the table
381  if (particleData) {
382  // Check if their life time is bigger than 1e-14
383  if (particleDataTable_->particle(particleID)->lifetime() > longLivedDecayLength_)
384  update(flags_[LongLivedDecay], true);
385 
386  // Check for B and C weak decays
387  update(flags_[BWeakDecay], particleID.hasBottom());
388  update(flags_[CWeakDecay], particleID.hasCharm());
389 
390  // Check for B or C pure leptonic decays
391  int daughtId = abs((*iparticle)->pdgId());
392  update(flags_[FromBWeakDecayMuon], particleID.hasBottom() && daughtId == 13);
393  update(flags_[FromCWeakDecayMuon], particleID.hasCharm() && daughtId == 13);
394  }
395  // Check decays
396  update(flags_[ChargePionDecay], pdgid == 211);
397  update(flags_[ChargeKaonDecay], pdgid == 321);
398  update(flags_[TauDecay], pdgid == 15);
399  update(flags_[KsDecay], pdgid == 310);
400  update(flags_[LambdaDecay], pdgid == 3122);
401  update(flags_[JpsiDecay], pdgid == 443);
402  update(flags_[XiDecay], pdgid == 3312);
403  update(flags_[OmegaDecay], pdgid == 3334);
404  update(flags_[SigmaPlusDecay], pdgid == 3222);
405  update(flags_[SigmaMinusDecay], pdgid == 3112);
406  }
407  }
408  }
409  // Decays in flight
412  update(flags_[DecayOnFlightMuon], flags_[FromChargePionMuon] || flags_[FromChargeKaonMuon]);
413 }
414 
416  // Get the main primary vertex from the list
417  GeneratedPrimaryVertex const &genpv = genpvs_.back();
418 
419  // Get the generated history of the tracks
420  const TrackHistory::GenParticleTrail &genParticleTrail = tracer_.genParticleTrail();
421 
422  // Vertex counter
423  int counter = 0;
424 
425  // Unit transformation from mm to cm
426  double const mm = 0.1;
427 
428  double oldX = genpv.x;
429  double oldY = genpv.y;
430  double oldZ = genpv.z;
431 
432  // Loop over the generated particles
433  for (TrackHistory::GenParticleTrail::const_reverse_iterator iparticle = genParticleTrail.rbegin();
434  iparticle != genParticleTrail.rend();
435  ++iparticle) {
436  // Look for those with production vertex
437  HepMC::GenVertex *parent = (*iparticle)->production_vertex();
438  if (parent) {
439  HepMC::ThreeVector p = parent->point3d();
440 
441  double distance2 = pow(p.x() * mm - genpv.x, 2) + pow(p.y() * mm - genpv.y, 2) + pow(p.z() * mm - genpv.z, 2);
442  double difference2 = pow(p.x() * mm - oldX, 2) + pow(p.y() * mm - oldY, 2) + pow(p.z() * mm - oldZ, 2);
443 
444  // std::cout << "Distance2 : " << distance2 << " (" << p.x() * mm << ","
445  // << p.y() * mm << "," << p.z() * mm << ")" << std::endl; std::cout <<
446  // "Difference2 : " << difference2 << std::endl;
447 
448  if (difference2 > vertexClusteringSqDistance_) {
449  if (distance2 > vertexClusteringSqDistance_)
450  counter++;
451  oldX = p.x() * mm;
452  oldY = p.y() * mm;
453  oldZ = p.z() * mm;
454  }
455  }
456  }
457 
458  const TrackHistory::SimParticleTrail &simParticleTrail = tracer_.simParticleTrail();
459 
460  // Loop over the generated particles
461  for (TrackHistory::SimParticleTrail::const_reverse_iterator iparticle = simParticleTrail.rbegin();
462  iparticle != simParticleTrail.rend();
463  ++iparticle) {
464  // Look for those with production vertex
465  TrackingParticle::Point p = (*iparticle)->vertex();
466 
467  double distance2 = pow(p.x() - genpv.x, 2) + pow(p.y() - genpv.y, 2) + pow(p.z() - genpv.z, 2);
468  double difference2 = pow(p.x() - oldX, 2) + pow(p.y() - oldY, 2) + pow(p.z() - oldZ, 2);
469 
470  // std::cout << "Distance2 : " << distance2 << " (" << p.x() << "," << p.y()
471  // << "," << p.z() << ")" << std::endl; std::cout << "Difference2 : " <<
472  // difference2 << std::endl;
473 
474  if (difference2 > vertexClusteringSqDistance_) {
475  if (distance2 > vertexClusteringSqDistance_)
476  counter++;
477  oldX = p.x();
478  oldY = p.y();
479  oldZ = p.z();
480  }
481  }
482 
483  if (!counter)
484  flags_[PrimaryVertex] = true;
485  else if (counter == 1)
486  flags_[SecondaryVertex] = true;
487  else
488  flags_[TertiaryVertex] = true;
489 }
490 
491 bool TrackClassifier::isFinalstateParticle(const HepMC::GenParticle *p) { return !p->end_vertex() && p->status() == 1; }
492 
494  const ParticleData *part = particleDataTable_->particle(p->pdg_id());
495  if (part)
496  return part->charge() != 0;
497  else {
498  // the new/improved particle table doesn't know anti-particles
499  return particleDataTable_->particle(-p->pdg_id()) != nullptr;
500  }
501 }
502 
504  genpvs_.clear();
505 
506  const HepMC::GenEvent *event = mcInformation_->GetEvent();
507 
508  if (event) {
509  int idx = 0;
510 
511  // Loop over the different GenVertex
512  for (HepMC::GenEvent::vertex_const_iterator ivertex = event->vertices_begin(); ivertex != event->vertices_end();
513  ++ivertex) {
514  bool hasParentVertex = false;
515 
516  // Loop over the parents looking to see if they are coming from a
517  // production vertex
518  for (HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(HepMC::parents);
519  iparent != (*ivertex)->particles_end(HepMC::parents);
520  ++iparent)
521  if ((*iparent)->production_vertex()) {
522  hasParentVertex = true;
523  break;
524  }
525 
526  // Reject those vertices with parent vertices
527  if (hasParentVertex)
528  continue;
529 
530  // Get the position of the vertex
531  HepMC::FourVector pos = (*ivertex)->position();
532 
533  double const mm = 0.1;
534 
535  GeneratedPrimaryVertex pv(pos.x() * mm, pos.y() * mm, pos.z() * mm);
536 
537  std::vector<GeneratedPrimaryVertex>::iterator ientry = genpvs_.begin();
538 
539  // Search for a VERY close vertex in the list
540  for (; ientry != genpvs_.end(); ++ientry) {
541  double distance2 = pow(pv.x - ientry->x, 2) + pow(pv.y - ientry->y, 2) + pow(pv.z - ientry->z, 2);
542  if (distance2 < vertexClusteringSqDistance_)
543  break;
544  }
545 
546  // Check if there is not a VERY close vertex and added to the list
547  if (ientry == genpvs_.end())
548  ientry = genpvs_.insert(ientry, pv);
549 
550  // Add the vertex barcodes to the new or existent vertices
551  ientry->genVertex.push_back((*ivertex)->barcode());
552 
553  // Collect final state descendants
554  for (HepMC::GenVertex::particle_iterator idecendants = (*ivertex)->particles_begin(HepMC::descendants);
555  idecendants != (*ivertex)->particles_end(HepMC::descendants);
556  ++idecendants) {
557  if (isFinalstateParticle(*idecendants))
558  if (find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) ==
559  ientry->finalstateParticles.end()) {
560  ientry->finalstateParticles.push_back((*idecendants)->barcode());
561  HepMC::FourVector m = (*idecendants)->momentum();
562 
563  ientry->ptot.setPx(ientry->ptot.px() + m.px());
564  ientry->ptot.setPy(ientry->ptot.py() + m.py());
565  ientry->ptot.setPz(ientry->ptot.pz() + m.pz());
566  ientry->ptot.setE(ientry->ptot.e() + m.e());
567  ientry->ptsq += m.perp() * m.perp();
568 
569  if (m.perp() > 0.8 && std::abs(m.pseudoRapidity()) < 2.5 && isCharged(*idecendants))
570  ientry->nGenTrk++;
571  }
572  }
573  idx++;
574  }
575  }
576 
577  std::sort(genpvs_.begin(), genpvs_.end());
578 }
unsigned int numberOfLayers() const
Return the number of layers with simulated and/or reconstructed hits.
Definition: TrackQuality.h:74
T getUntrackedParameter(std::string const &, T const &) const
const edm::InputTag beamSpotLabel_
const edm::InputTag hepMCLabel_
void newEvent(const edm::Event &, const edm::EventSetup &)
Pre-process event information (for accessing reconstruction information)
Definition: TrackHistory.cc:32
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
int event() const
get the contents of the subdetector field (should be protected?)
TPRegexp parents
Definition: eve_filter.cc:21
const FreeTrajectoryState & theState() const
double dxyError() const
error on dxy
Definition: TrackBase.h:769
const reco::GenParticle * recoGenParticle() const
Definition: HistoryBase.h:82
const G4toCMSLegacyProcTypeMap g4toCMSProcMap_
constexpr uint32_t maxLayers
const TrackingParticleRef & simParticle() const
Return the initial tracking particle from the history.
Definition: HistoryBase.h:67
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::Handle< edm::HepMCProduct > mcInformation_
const TrackerTopology * tTopo_
T y() const
Definition: PV3DBase.h:60
unsigned int numberOfInnerLayers_
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:301
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
unsigned int minTrackerSimHits_
std::vector< GeneratedPrimaryVertex > genpvs_
void processesAtSimulation()
Get information about conversion and other interactions.
constexpr std::array< uint8_t, layerIndexSize > layer
std::vector< const HepMC::GenParticle * > GenParticleTrail
HepMC::GenParticle trail type.
Definition: HistoryBase.h:15
int pdgId() const final
PDG identifier.
math::XYZPointD Point
point in the space
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void simulationInformation()
Get all the information related to the simulation details.
edm::ESGetToken< ParticleDataTable, PDTRecord > particleDataTableToken_
int TrackCharge
Definition: TrackCharge.h:4
std::vector< const reco::GenParticle * > RecoGenParticleTrail
reco::GenParticle trail type.
Definition: HistoryBase.h:18
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_
void evaluate(SimParticleTrail const &, reco::TrackBaseRef const &, const TrackerTopology *tTopo)
Compute information about the track reconstruction quality.
Auxiliary class holding simulated primary vertices.
int bunchCrossing() const
get the detector field from this detid
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
math::XYZPoint Point
point in the space
Definition: TrackBase.h:80
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)
const unsigned int processId(unsigned int g4ProcessId) const
Definition: Utils.cc:59
const Layer & layer(unsigned int index) const
Return information about the given layer by index.
Definition: TrackQuality.h:77
edm::Handle< reco::BeamSpot > beamSpot_
HepPDT::ParticleData ParticleData
GlobalVector momentum() const
TrackQuality quality_
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
double dzError() const
error on dz
Definition: TrackBase.h:778
GlobalPoint position() const
Get track history and classify it in function of their .
part
Definition: HCALResponse.h:20
T const * product() const
Definition: ESHandle.h:86
bool isFinalstateParticle(const HepMC::GenParticle *)
SimParticleTrail const & simParticleTrail() const
Return all the simulated particle in the history.
Definition: HistoryBase.h:55
void vertexInformation()
Get geometrical information about the vertices.
tuple config
parse the configuration file
void processesAtGenerator()
Get all the information related to decay process.
bool evaluate(TrackingParticleRef tpr)
Evaluate track history using a TrackingParticleRef.
Definition: TrackHistory.h:37
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoHandToken_
#define update(a, b)
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:49
edm::ESHandle< MagneticField > magneticField_
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
static std::atomic< unsigned int > counter
std::vector< Hit > hits
Definition: TrackQuality.h:56
Flags flags_
Flag containers.
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
int charge() const
track electric charge
Definition: TrackBase.h:596
const reco::TrackBaseRef & recoTrack() const
Return a reference to the reconstructed track.
Definition: TrackHistory.h:55
void newEvent(const edm::Event &, const edm::EventSetup &)
Pre-process event information (for accessing reconstruction information)
double longLivedDecayLength_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
void reset()
Reset the categories flags.
double vertexClusteringSqDistance_
tuple process
Definition: LaserDQM_cfg.py:3
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
RecoGenParticleTrail const & recoGenParticleTrail() const
Return all reco::GenParticle in the history.
Definition: HistoryBase.h:64
T x() const
Definition: PV3DBase.h:59
edm::ESHandle< ParticleDataTable > particleDataTable_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
GenParticleTrail const & genParticleTrail() const
Return all generated particle (HepMC::GenParticle) in the history.
Definition: HistoryBase.h:61
TrackHistory tracer_
Global3DVector GlobalVector
Definition: GlobalVector.h:10
std::vector< TrackingParticleRef > SimParticleTrail
SimParticle trail type.
Definition: HistoryBase.h:30