CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
VertexClassifier Class Reference

Get track history and classify it in function of their . More...

#include <VertexClassifier.h>

Inheritance diagram for VertexClassifier:
VertexCategories VertexClassifierByProxy< Collection > VertexClassifierByProxy< reco::SecondaryVertexTagInfoCollection >

Classes

struct  GeneratedPrimaryVertex
 Auxiliary class holding simulated primary vertices. More...
 

Public Types

typedef VertexCategories Categories
 Type to the associate category. More...
 
- Public Types inherited from VertexCategories
enum  Category {
  Fake = 0, Reconstructed = Fake, SignalEvent, BWeakDecay,
  CWeakDecay, TauDecay, KsDecay, LambdaDecay,
  JpsiDecay, XiDecay, OmegaDecay, SigmaPlusDecay,
  SigmaMinusDecay, LongLivedDecay, KnownProcess, UndefinedProcess,
  UnknownProcess, PrimaryProcess, HadronicProcess, DecayProcess,
  ComptonProcess, AnnihilationProcess, EIoniProcess, HIoniProcess,
  MuIoniProcess, PhotonProcess, MuPairProdProcess, ConversionsProcess,
  EBremProcess, SynchrotronRadiationProcess, MuBremProcess, MuNuclProcess,
  PrimaryVertex, SecondaryVertex, TertiaryVertex, TierciaryVertex = TertiaryVertex,
  Unknown
}
 Categories available to vertexes. More...
 
typedef std::vector< bool > Flags
 Main types associated to the class. More...
 

Public Member Functions

VertexClassifier const & evaluate (reco::VertexBaseRef const &)
 Classify the RecoVertex in categories. More...
 
VertexClassifier const & evaluate (TrackingVertexRef const &)
 Classify the TrackingVertex in categories. More...
 
VertexClassifier const & evaluate (reco::VertexRef const &vertex)
 Classify the RecoVertex in categories. More...
 
VertexHistory const & history () const
 Returns a reference to the vertex history used in the classification. More...
 
virtual void newEvent (edm::Event const &, edm::EventSetup const &)
 Pre-process event information (for accessing reconstraction information) More...
 
 VertexClassifier (edm::ParameterSet const &pset, edm::ConsumesCollector)
 Constructor by ParameterSet. More...
 
virtual ~VertexClassifier ()
 
- Public Member Functions inherited from VertexCategories
const Flagsflags () const
 Returns flags with the category descriptions. More...
 
bool is (Category category) const
 Returns track flag for a given category. More...
 
 VertexCategories ()
 Void constructor. More...
 

Private Member Functions

void genPrimaryVertices ()
 
bool isCharged (const HepMC::GenParticle *)
 
bool isFinalstateParticle (const HepMC::GenParticle *)
 
void processesAtGenerator ()
 Get all the information related to decay process. More...
 
void processesAtSimulation ()
 Get information about conversion and other interactions. More...
 
void reconstructionInformation (reco::TrackBaseRef const &)
 Get reconstruction information. More...
 
void simulationInformation ()
 Get all the information related to the simulation details. More...
 
void vertexInformation ()
 Get geometrical information about the vertices. More...
 

Private Attributes

const G4toCMSLegacyProcTypeMap g4toCMSProcMap_
 
std::vector< GeneratedPrimaryVertexgenpvs_
 
const edm::InputTag hepMCLabel_
 
double longLivedDecayLength_
 
edm::Handle< edm::HepMCProductmcInformation_
 
edm::ESHandle< ParticleDataTableparticleDataTable_
 
edm::ESGetToken< ParticleDataTable, edm::DefaultRecordparticleDataTableToken_
 
VertexHistory tracer_
 
double vertexClusteringDistance_
 

Additional Inherited Members

- Static Public Attributes inherited from VertexCategories
static const char *const Names []
 Name of the different categories. More...
 
- Protected Member Functions inherited from VertexCategories
void reset ()
 Reset the categories flags. More...
 
void unknownVertex ()
 
- Protected Attributes inherited from VertexCategories
Flags flags_
 Flag containers. More...
 

Detailed Description

Get track history and classify it in function of their .

Definition at line 18 of file VertexClassifier.h.

Member Typedef Documentation

◆ Categories

Type to the associate category.

Definition at line 21 of file VertexClassifier.h.

Constructor & Destructor Documentation

◆ VertexClassifier()

VertexClassifier::VertexClassifier ( edm::ParameterSet const &  pset,
edm::ConsumesCollector  collector 
)

Constructor by ParameterSet.

Definition at line 18 of file VertexClassifier.cc.

References edm::ConsumesCollector::consumes(), HistoryBase::depth(), hepMCLabel_, longLivedDecayLength_, tracer_, and vertexClusteringDistance_.

19  : VertexCategories(),
20  tracer_(config, collector),
21  hepMCLabel_(config.getUntrackedParameter<edm::InputTag>("hepMC")),
22  particleDataTableToken_(collector.esConsumes()) {
24  // Set the history depth after hadronization
25  tracer_.depth(-2);
26 
27  // Set the minimum decay length for detecting long decays
28  longLivedDecayLength_ = config.getUntrackedParameter<double>("longLivedDecayLength");
29 
30  // Set the distance for clustering vertices
31  vertexClusteringDistance_ = config.getUntrackedParameter<double>("vertexClusteringDistance");
32 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::ESGetToken< ParticleDataTable, edm::DefaultRecord > particleDataTableToken_
VertexHistory tracer_
const edm::InputTag hepMCLabel_
Definition: config.py:1
double longLivedDecayLength_
void depth(int d)
Set the depth of the history.
Definition: HistoryBase.h:49
double vertexClusteringDistance_
VertexCategories()
Void constructor.

◆ ~VertexClassifier()

virtual VertexClassifier::~VertexClassifier ( )
inlinevirtual

Definition at line 26 of file VertexClassifier.h.

26 {}

Member Function Documentation

◆ evaluate() [1/3]

VertexClassifier const & VertexClassifier::evaluate ( reco::VertexBaseRef const &  vertex)

Classify the RecoVertex in categories.

Definition at line 48 of file VertexClassifier.cc.

References VertexHistory::evaluate(), VertexCategories::Fake, VertexCategories::flags_, processesAtGenerator(), processesAtSimulation(), VertexCategories::reset(), simulationInformation(), tracer_, VertexCategories::unknownVertex(), bphysicsOniaDQM_cfi::vertex, and vertexInformation().

Referenced by VertexHistoryAnalyzer::analyze(), and VertexClassifierByProxy< reco::SecondaryVertexTagInfoCollection >::evaluate().

48  {
49  // Initializing the category vector
50  reset();
51 
52  // Associate and evaluate the vertex history (check for fakes)
53  if (tracer_.evaluate(vertex)) {
54  // Get all the information related to the simulation details
56 
57  // Get all the information related to decay process
59 
60  // Get information about conversion and other interactions
62 
63  // Get geometrical information about the vertices
65 
66  // Check for unkown classification
67  unknownVertex();
68  } else
69  flags_[Fake] = true;
70 
71  return *this;
72 }
VertexHistory tracer_
void vertexInformation()
Get geometrical information about the vertices.
void reset()
Reset the categories flags.
bool evaluate(TrackingVertexRef tvr)
Evaluate track history using a TrackingParticleRef.
Definition: VertexHistory.h:38
Flags flags_
Flag containers.
void processesAtGenerator()
Get all the information related to decay process.
void processesAtSimulation()
Get information about conversion and other interactions.
void simulationInformation()
Get all the information related to the simulation details.

◆ evaluate() [2/3]

VertexClassifier const & VertexClassifier::evaluate ( TrackingVertexRef const &  vertex)

Classify the TrackingVertex in categories.

Definition at line 74 of file VertexClassifier.cc.

References VertexHistory::evaluate(), VertexCategories::flags_, edm::RefToBase< T >::isNonnull(), processesAtGenerator(), processesAtSimulation(), VertexCategories::Reconstructed, VertexHistory::recoVertex(), VertexCategories::reset(), simulationInformation(), tracer_, VertexCategories::unknownVertex(), bphysicsOniaDQM_cfi::vertex, and vertexInformation().

74  {
75  // Initializing the category vector
76  reset();
77 
78  // Trace the history for the given TP
80 
81  // Check for a reconstructed track
83  flags_[Reconstructed] = true;
84  else
85  flags_[Reconstructed] = false;
86 
87  // Get all the information related to the simulation details
89 
90  // Get all the information related to decay process
92 
93  // Get information about conversion and other interactions
95 
96  // Get geometrical information about the vertices
98 
99  // Check for unkown classification
100  unknownVertex();
101 
102  return *this;
103 }
VertexHistory tracer_
const reco::VertexBaseRef & recoVertex() const
Return a reference to the reconstructed track.
Definition: VertexHistory.h:56
void vertexInformation()
Get geometrical information about the vertices.
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:303
void reset()
Reset the categories flags.
bool evaluate(TrackingVertexRef tvr)
Evaluate track history using a TrackingParticleRef.
Definition: VertexHistory.h:38
Flags flags_
Flag containers.
void processesAtGenerator()
Get all the information related to decay process.
void processesAtSimulation()
Get information about conversion and other interactions.
void simulationInformation()
Get all the information related to the simulation details.

◆ evaluate() [3/3]

VertexClassifier const& VertexClassifier::evaluate ( reco::VertexRef const &  vertex)
inline

Classify the RecoVertex in categories.

Definition at line 38 of file VertexClassifier.h.

References evaluate(), and bphysicsOniaDQM_cfi::vertex.

Referenced by evaluate().

VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.

◆ genPrimaryVertices()

void VertexClassifier::genPrimaryVertices ( )
private

Definition at line 382 of file VertexClassifier.cc.

References funct::abs(), HLT_2023v12_cff::distance, spr::find(), genpvs_, edm::HepMCProduct::GetEvent(), isCharged(), isFinalstateParticle(), visualization-live-secondInstance_cfg::m, mcInformation_, parents, funct::pow(), MetAnalyzer::pv(), jetUpdater_cfi::sort, mathSSE::sqrt(), and vertexClusteringDistance_.

Referenced by newEvent().

382  {
383  genpvs_.clear();
384 
385  const HepMC::GenEvent *event = mcInformation_->GetEvent();
386 
387  if (event) {
388  // Loop over the different GenVertex
389  for (HepMC::GenEvent::vertex_const_iterator ivertex = event->vertices_begin(); ivertex != event->vertices_end();
390  ++ivertex) {
391  bool hasParentVertex = false;
392 
393  // Loop over the parents looking to see if they are coming from a
394  // production vertex
395  for (HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(HepMC::parents);
396  iparent != (*ivertex)->particles_end(HepMC::parents);
397  ++iparent)
398  if ((*iparent)->production_vertex()) {
399  hasParentVertex = true;
400  break;
401  }
402 
403  // Reject those vertices with parent vertices
404  if (hasParentVertex)
405  continue;
406 
407  // Get the position of the vertex
408  HepMC::FourVector pos = (*ivertex)->position();
409 
410  double const mm = 0.1;
411 
412  GeneratedPrimaryVertex pv(pos.x() * mm, pos.y() * mm, pos.z() * mm);
413 
414  std::vector<GeneratedPrimaryVertex>::iterator ientry = genpvs_.begin();
415 
416  // Search for a VERY close vertex in the list
417  for (; ientry != genpvs_.end(); ++ientry) {
418  double distance = sqrt(pow(pv.x - ientry->x, 2) + pow(pv.y - ientry->y, 2) + pow(pv.z - ientry->z, 2));
420  break;
421  }
422 
423  // Check if there is not a VERY close vertex and added to the list
424  if (ientry == genpvs_.end())
425  ientry = genpvs_.insert(ientry, pv);
426 
427  // Add the vertex barcodes to the new or existent vertices
428  ientry->genVertex.push_back((*ivertex)->barcode());
429 
430  // Collect final state descendants
431  for (HepMC::GenVertex::particle_iterator idecendants = (*ivertex)->particles_begin(HepMC::descendants);
432  idecendants != (*ivertex)->particles_end(HepMC::descendants);
433  ++idecendants) {
434  if (isFinalstateParticle(*idecendants))
435  if (find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) ==
436  ientry->finalstateParticles.end()) {
437  ientry->finalstateParticles.push_back((*idecendants)->barcode());
438  HepMC::FourVector m = (*idecendants)->momentum();
439 
440  ientry->ptot.setPx(ientry->ptot.px() + m.px());
441  ientry->ptot.setPy(ientry->ptot.py() + m.py());
442  ientry->ptot.setPz(ientry->ptot.pz() + m.pz());
443  ientry->ptot.setE(ientry->ptot.e() + m.e());
444  ientry->ptsq += m.perp() * m.perp();
445 
446  if (m.perp() > 0.8 && std::abs(m.pseudoRapidity()) < 2.5 && isCharged(*idecendants))
447  ientry->nGenTrk++;
448  }
449  }
450  }
451  }
452 
453  std::sort(genpvs_.begin(), genpvs_.end());
454 }
edm::Handle< edm::HepMCProduct > mcInformation_
TPRegexp parents
Definition: eve_filter.cc:21
std::vector< GeneratedPrimaryVertex > genpvs_
bool isCharged(const HepMC::GenParticle *)
bool isFinalstateParticle(const HepMC::GenParticle *)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
T sqrt(T t)
Definition: SSEVec.h:19
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
double vertexClusteringDistance_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Definition: event.py:1

◆ history()

VertexHistory const& VertexClassifier::history ( ) const
inline

Returns a reference to the vertex history used in the classification.

Definition at line 41 of file VertexClassifier.h.

References tracer_.

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

41 { return tracer_; }
VertexHistory tracer_

◆ isCharged()

bool VertexClassifier::isCharged ( const HepMC::GenParticle *  p)
private

Definition at line 372 of file VertexClassifier.cc.

References AlCaHLTBitMon_ParallelJobs::p, and particleDataTable_.

Referenced by genPrimaryVertices().

372  {
373  const ParticleData *part = particleDataTable_->particle(p->pdg_id());
374  if (part)
375  return part->charge() != 0;
376  else {
377  // the new/improved particle table doesn't know anti-particles
378  return particleDataTable_->particle(-p->pdg_id()) != nullptr;
379  }
380 }
edm::ESHandle< ParticleDataTable > particleDataTable_
HepPDT::ParticleData ParticleData
part
Definition: HCALResponse.h:20

◆ isFinalstateParticle()

bool VertexClassifier::isFinalstateParticle ( const HepMC::GenParticle *  p)
private

Definition at line 368 of file VertexClassifier.cc.

References AlCaHLTBitMon_ParallelJobs::p.

Referenced by genPrimaryVertices().

368  {
369  return !p->end_vertex() && p->status() == 1;
370 }

◆ newEvent()

void VertexClassifier::newEvent ( edm::Event const &  event,
edm::EventSetup const &  setup 
)
virtual

Pre-process event information (for accessing reconstraction information)

Reimplemented in VertexClassifierByProxy< Collection >, and VertexClassifierByProxy< reco::SecondaryVertexTagInfoCollection >.

Definition at line 34 of file VertexClassifier.cc.

References genPrimaryVertices(), hepMCLabel_, mcInformation_, VertexHistory::newEvent(), particleDataTable_, particleDataTableToken_, singleTopDQM_cfi::setup, and tracer_.

Referenced by VertexHistoryAnalyzer::analyze(), and VertexClassifierByProxy< reco::SecondaryVertexTagInfoCollection >::newEvent().

34  {
35  // Get the new event information for the tracer
37 
38  // Get hepmc of the event
39  event.getByLabel(hepMCLabel_, mcInformation_);
40 
41  // Get the partivle data table
43 
44  // Create the list of primary vertices associated to the event
46 }
edm::Handle< edm::HepMCProduct > mcInformation_
edm::ESGetToken< ParticleDataTable, edm::DefaultRecord > particleDataTableToken_
VertexHistory tracer_
void newEvent(const edm::Event &, const edm::EventSetup &)
Pre-process event information (for accessing reconstruction information)
const edm::InputTag hepMCLabel_
edm::ESHandle< ParticleDataTable > particleDataTable_
Definition: event.py:1

◆ processesAtGenerator()

void VertexClassifier::processesAtGenerator ( )
private

Get all the information related to decay process.

Definition at line 112 of file VertexClassifier.cc.

References funct::abs(), VertexCategories::BWeakDecay, VertexCategories::CWeakDecay, VertexCategories::flags_, HistoryBase::genVertexTrail(), VertexCategories::JpsiDecay, VertexCategories::KsDecay, VertexCategories::LambdaDecay, VertexCategories::LongLivedDecay, longLivedDecayLength_, VertexCategories::OmegaDecay, parents, particleDataTable_, EgammaObjectsElectrons_cfi::particleID, LHEGenericFilter_cfi::ParticleID, EgammaValidation_cff::pdgid, VertexCategories::SigmaMinusDecay, VertexCategories::SigmaPlusDecay, tracer_, update, bphysicsOniaDQM_cfi::vertex, and VertexCategories::XiDecay.

Referenced by evaluate().

112  {
113  // Get the generated vetices from track history
114  VertexHistory::GenVertexTrail const &genVertexTrail = tracer_.genVertexTrail();
115 
116  // Loop over the generated vertices
117  for (VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin(); ivertex != genVertexTrail.end();
118  ++ivertex) {
119  // Get the pointer to the vertex by removing the const-ness (no const methos
120  // in HepMC::GenVertex)
121  HepMC::GenVertex *vertex = const_cast<HepMC::GenVertex *>(*ivertex);
122 
123  // Loop over the sources looking for specific decays
124  for (HepMC::GenVertex::particle_iterator iparent = vertex->particles_begin(HepMC::parents);
125  iparent != vertex->particles_end(HepMC::parents);
126  ++iparent) {
127  // Collect the pdgid of the parent
128  int pdgid = std::abs((*iparent)->pdg_id());
129  // Get particle type
131 
132  // Check if the particle type is valid one
133  if (particleID.isValid()) {
134  // Get particle data
135  ParticleData const *particleData = particleDataTable_->particle(particleID);
136  // Check if the particle exist in the table
137  if (particleData) {
138  // Check if their life time is bigger than longLivedDecayLength_
139  if (particleData->lifetime() > longLivedDecayLength_) {
140  // Check for B, C weak decays and long lived decays
141  update(flags_[BWeakDecay], particleID.hasBottom());
142  update(flags_[CWeakDecay], particleID.hasCharm());
143  update(flags_[LongLivedDecay], true);
144  }
145  // Check Tau, Ks and Lambda decay
146  update(flags_[TauDecay], pdgid == 15);
147  update(flags_[KsDecay], pdgid == 310);
148  update(flags_[LambdaDecay], pdgid == 3122);
149  update(flags_[JpsiDecay], pdgid == 443);
150  update(flags_[XiDecay], pdgid == 3312);
151  update(flags_[OmegaDecay], pdgid == 3334);
152  update(flags_[SigmaPlusDecay], pdgid == 3222);
153  update(flags_[SigmaMinusDecay], pdgid == 3112);
154  }
155  }
156  }
157  }
158 }
VertexHistory tracer_
TPRegexp parents
Definition: eve_filter.cc:21
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:58
Flags flags_
Flag containers.
double longLivedDecayLength_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
Definition: HistoryBase.h:24
edm::ESHandle< ParticleDataTable > particleDataTable_
HepPDT::ParticleData ParticleData
#define update(a, b)

◆ processesAtSimulation()

void VertexClassifier::processesAtSimulation ( )
private

Get information about conversion and other interactions.

Definition at line 160 of file VertexClassifier.cc.

References funct::abs(), CMS::Annihilation, VertexCategories::AnnihilationProcess, VertexCategories::BWeakDecay, CMS::Compton, VertexCategories::ComptonProcess, CMS::Conversions, VertexCategories::ConversionsProcess, VertexCategories::CWeakDecay, CMS::Decay, VertexCategories::DecayProcess, CMS::EBrem, VertexCategories::EBremProcess, CMS::EIoni, VertexCategories::EIoniProcess, RemoveAddSevLevel::flag, VertexCategories::flags_, g4toCMSProcMap_, CMS::Hadronic, VertexCategories::HadronicProcess, CMS::HIoni, VertexCategories::HIoniProcess, VertexCategories::JpsiDecay, VertexCategories::KnownProcess, VertexCategories::KsDecay, VertexCategories::LambdaDecay, VertexCategories::LongLivedDecay, longLivedDecayLength_, CMS::MuBrem, VertexCategories::MuBremProcess, CMS::MuIoni, VertexCategories::MuIoniProcess, CMS::MuNucl, VertexCategories::MuNuclProcess, CMS::MuPairProd, VertexCategories::MuPairProdProcess, VertexCategories::OmegaDecay, particleDataTable_, EgammaObjectsElectrons_cfi::particleID, LHEGenericFilter_cfi::ParticleID, EgammaValidation_cff::pdgid, CMS::Photon, VertexCategories::PhotonProcess, CMS::Primary, VertexCategories::PrimaryProcess, LaserDQM_cfg::process, G4toCMSLegacyProcTypeMap::processId(), VertexCategories::SigmaMinusDecay, VertexCategories::SigmaPlusDecay, HistoryBase::simVertexTrail(), CMS::SynchrotronRadiation, VertexCategories::SynchrotronRadiationProcess, tracer_, CMS::Undefined, VertexCategories::UndefinedProcess, CMS::Unknown, VertexCategories::UnknownProcess, update, and VertexCategories::XiDecay.

Referenced by evaluate().

160  {
161  VertexHistory::SimVertexTrail const &simVertexTrail = tracer_.simVertexTrail();
162 
163  for (VertexHistory::SimVertexTrail::const_iterator ivertex = simVertexTrail.begin(); ivertex != simVertexTrail.end();
164  ++ivertex) {
165  // pdgid of the real source parent vertex
166  int pdgid = 0;
167 
168  // select the original source in case of combined vertices
169  bool flag = false;
171 
172  for (its = (*ivertex)->sourceTracks_begin(); its != (*ivertex)->sourceTracks_end(); ++its) {
173  for (itd = (*ivertex)->daughterTracks_begin(); itd != (*ivertex)->daughterTracks_end(); ++itd)
174  if (itd != its) {
175  flag = true;
176  break;
177  }
178  if (flag)
179  break;
180  }
181  // Collect the pdgid of the original source track
182  if (its != (*ivertex)->sourceTracks_end())
183  pdgid = std::abs((*its)->pdgId());
184  else
185  pdgid = 0;
186 
187  // Geant4 process type is selected using first Geant4 vertex assigned to
188  // the TrackingVertex
189  unsigned int processG4 = 0;
190 
191  if ((*ivertex)->nG4Vertices() > 0) {
192  processG4 = (*(*ivertex)->g4Vertices_begin()).processType();
193  }
194 
195  unsigned int process = g4toCMSProcMap_.processId(processG4);
196 
197  // Flagging all the different processes
199 
217 
218  // Loop over the simulated particles
219  for (TrackingVertex::tp_iterator iparticle = (*ivertex)->daughterTracks_begin();
220  iparticle != (*ivertex)->daughterTracks_end();
221  ++iparticle) {
222  if ((*iparticle)->numberOfTrackerLayers()) {
223  // Special treatment for decays
224  if (process == CMS::Decay) {
225  // Get particle type
227  // Check if the particle type is valid one
228  if (particleID.isValid()) {
229  // Get particle data
230  ParticleData const *particleData = particleDataTable_->particle(particleID);
231  // Check if the particle exist in the table
232  if (particleData) {
233  // Check if their life time is bigger than 1e-14
234  if (particleDataTable_->particle(particleID)->lifetime() > longLivedDecayLength_) {
235  // Check for B, C weak decays and long lived decays
236  update(flags_[BWeakDecay], particleID.hasBottom());
237  update(flags_[CWeakDecay], particleID.hasCharm());
238  update(flags_[LongLivedDecay], true);
239  }
240  // Check Tau, Ks and Lambda decay
241  update(flags_[TauDecay], pdgid == 15);
242  update(flags_[KsDecay], pdgid == 310);
243  update(flags_[LambdaDecay], pdgid == 3122);
244  update(flags_[JpsiDecay], pdgid == 443);
245  update(flags_[XiDecay], pdgid == 3312);
246  update(flags_[OmegaDecay], pdgid == 3334);
247  update(flags_[SigmaPlusDecay], pdgid == 3222);
248  update(flags_[SigmaMinusDecay], pdgid == 3112);
249  }
250  }
251  }
252  }
253  }
254  }
255 }
VertexHistory tracer_
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:52
const unsigned int processId(unsigned int g4ProcessId) const
Definition: Utils.cc:59
Flags flags_
Flag containers.
double longLivedDecayLength_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::ESHandle< ParticleDataTable > particleDataTable_
const G4toCMSLegacyProcTypeMap g4toCMSProcMap_
HepPDT::ParticleData ParticleData
#define update(a, b)
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:33

◆ reconstructionInformation()

void VertexClassifier::reconstructionInformation ( reco::TrackBaseRef const &  )
private

Get reconstruction information.

◆ simulationInformation()

void VertexClassifier::simulationInformation ( )
private

Get all the information related to the simulation details.

Definition at line 105 of file VertexClassifier.cc.

References EncodedEventId::bunchCrossing(), EncodedEventId::event(), VertexCategories::flags_, VertexCategories::SignalEvent, HistoryBase::simVertex(), and tracer_.

Referenced by evaluate().

105  {
106  // Get the event id for the initial TP.
107  EncodedEventId eventId = tracer_.simVertex()->eventId();
108  // Check for signal events
109  flags_[SignalEvent] = !eventId.bunchCrossing() && !eventId.event();
110 }
VertexHistory tracer_
int event() const
get the contents of the subdetector field (should be protected?)
Flags flags_
Flag containers.
int bunchCrossing() const
get the detector field from this detid
const TrackingVertexRef & simVertex() const
Return the initial tracking vertex from the history.
Definition: HistoryBase.h:70

◆ vertexInformation()

void VertexClassifier::vertexInformation ( )
private

Get geometrical information about the vertices.

Definition at line 257 of file VertexClassifier.cc.

References caloTruthCellsNtuples_cff::Clusters, bsc_activity_cfg::clusters, HLT_2023v12_cff::distance, VertexCategories::flags_, genpvs_, HistoryBase::genVertexTrail(), AlCaHLTBitMon_ParallelJobs::p, funct::pow(), VertexCategories::PrimaryVertex, VertexCategories::SecondaryVertex, HistoryBase::simVertexTrail(), mathSSE::sqrt(), VertexCategories::TertiaryVertex, tracer_, vertexClusteringDistance_, VertexClassifier::GeneratedPrimaryVertex::x, VertexClassifier::GeneratedPrimaryVertex::y, and VertexClassifier::GeneratedPrimaryVertex::z.

Referenced by evaluate().

257  {
258  // Helper class for clusterization
259  typedef std::multimap<double, HepMC::ThreeVector> Clusters;
260  typedef std::pair<double, HepMC::ThreeVector> ClusterPair;
261 
263 
264  // Get the main primary vertex from the list
265  GeneratedPrimaryVertex const &genpv = genpvs_.back();
266 
267  // Get the generated history of the tracks
268  const VertexHistory::GenVertexTrail &genVertexTrail = tracer_.genVertexTrail();
269 
270  // Unit transformation from mm to cm
271  double const mm = 0.1;
272 
273  // Loop over the generated vertexes
274  for (VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin(); ivertex != genVertexTrail.end();
275  ++ivertex) {
276  // Check vertex exist
277  if (*ivertex) {
278  // Measure the distance2 respecto the primary vertex
279  HepMC::ThreeVector p = (*ivertex)->point3d();
280  double distance =
281  sqrt(pow(p.x() * mm - genpv.x, 2) + pow(p.y() * mm - genpv.y, 2) + pow(p.z() * mm - genpv.z, 2));
282 
283  // If there is not any clusters add the first vertex.
284  if (clusters.empty()) {
285  clusters.insert(ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)));
286  continue;
287  }
288 
289  // Check if there is already a cluster in the given distance from primary
290  // vertex
291  Clusters::const_iterator icluster = clusters.lower_bound(distance - vertexClusteringDistance_);
292 
293  if (icluster == clusters.upper_bound(distance + vertexClusteringDistance_)) {
294  clusters.insert(ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)));
295  continue;
296  }
297 
298  bool cluster = false;
299 
300  // Looping over the vertex clusters of a given distance from primary
301  // vertex
302  for (; icluster != clusters.upper_bound(distance + vertexClusteringDistance_); ++icluster) {
303  double difference = sqrt(pow(p.x() * mm - icluster->second.x(), 2) + pow(p.y() * mm - icluster->second.y(), 2) +
304  pow(p.z() * mm - icluster->second.z(), 2));
305 
306  if (difference < vertexClusteringDistance_) {
307  cluster = true;
308  break;
309  }
310  }
311 
312  if (!cluster)
313  clusters.insert(ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)));
314  }
315  }
316 
317  const VertexHistory::SimVertexTrail &simVertexTrail = tracer_.simVertexTrail();
318 
319  // Loop over the generated particles
320  for (VertexHistory::SimVertexTrail::const_reverse_iterator ivertex = simVertexTrail.rbegin();
321  ivertex != simVertexTrail.rend();
322  ++ivertex) {
323  // Look for those with production vertex
324  TrackingVertex::LorentzVector p = (*ivertex)->position();
325 
326  double distance = sqrt(pow(p.x() - genpv.x, 2) + pow(p.y() - genpv.y, 2) + pow(p.z() - genpv.z, 2));
327 
328  // If there is not any clusters add the first vertex.
329  if (clusters.empty()) {
330  clusters.insert(ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())));
331  continue;
332  }
333 
334  // Check if there is already a cluster in the given distance from primary
335  // vertex
336  Clusters::const_iterator icluster = clusters.lower_bound(distance - vertexClusteringDistance_);
337 
338  if (icluster == clusters.upper_bound(distance + vertexClusteringDistance_)) {
339  clusters.insert(ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())));
340  continue;
341  }
342 
343  bool cluster = false;
344 
345  // Looping over the vertex clusters of a given distance from primary vertex
346  for (; icluster != clusters.upper_bound(distance + vertexClusteringDistance_); ++icluster) {
347  double difference = sqrt(pow(p.x() - icluster->second.x(), 2) + pow(p.y() - icluster->second.y(), 2) +
348  pow(p.z() - icluster->second.z(), 2));
349 
350  if (difference < vertexClusteringDistance_) {
351  cluster = true;
352  break;
353  }
354  }
355 
356  if (!cluster)
357  clusters.insert(ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())));
358  }
359 
360  if (clusters.size() == 1)
361  flags_[PrimaryVertex] = true;
362  else if (clusters.size() == 2)
363  flags_[SecondaryVertex] = true;
364  else
365  flags_[TertiaryVertex] = true;
366 }
VertexHistory tracer_
std::vector< GeneratedPrimaryVertex > genpvs_
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
Definition: HistoryBase.h:58
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
Definition: HistoryBase.h:52
math::XYZTLorentzVectorD LorentzVector
Flags flags_
Flag containers.
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
Definition: HistoryBase.h:24
double vertexClusteringDistance_
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Definition: HistoryBase.h:33
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

Member Data Documentation

◆ g4toCMSProcMap_

const G4toCMSLegacyProcTypeMap VertexClassifier::g4toCMSProcMap_
private

Definition at line 46 of file VertexClassifier.h.

Referenced by processesAtSimulation().

◆ genpvs_

std::vector<GeneratedPrimaryVertex> VertexClassifier::genpvs_
private

Definition at line 90 of file VertexClassifier.h.

Referenced by genPrimaryVertices(), and vertexInformation().

◆ hepMCLabel_

const edm::InputTag VertexClassifier::hepMCLabel_
private

Definition at line 48 of file VertexClassifier.h.

Referenced by newEvent(), and VertexClassifier().

◆ longLivedDecayLength_

double VertexClassifier::longLivedDecayLength_
private

◆ mcInformation_

edm::Handle<edm::HepMCProduct> VertexClassifier::mcInformation_
private

Definition at line 53 of file VertexClassifier.h.

Referenced by genPrimaryVertices(), and newEvent().

◆ particleDataTable_

edm::ESHandle<ParticleDataTable> VertexClassifier::particleDataTable_
private

◆ particleDataTableToken_

edm::ESGetToken<ParticleDataTable, edm::DefaultRecord> VertexClassifier::particleDataTableToken_
private

Definition at line 56 of file VertexClassifier.h.

Referenced by newEvent().

◆ tracer_

VertexHistory VertexClassifier::tracer_
private

◆ vertexClusteringDistance_

double VertexClassifier::vertexClusteringDistance_
private

Definition at line 51 of file VertexClassifier.h.

Referenced by genPrimaryVertices(), VertexClassifier(), and vertexInformation().