CMS 3D CMS Logo

ParticleManager.h
Go to the documentation of this file.
1 #ifndef FASTSIM_PARTICLEMANAGER_H
2 #define FASTSIM_PARTICLEMANAGER_H
3 
5 #include "HepMC/GenEvent.h"
6 #include <vector>
7 #include <memory>
8 
11 
12 
14 // Author: L. Vanelderen, S. Kurz
15 // Date: 29 May 2017
17 
18 
19 namespace HepPDT
20 {
21  class ParticleDataTable;
22 }
23 
25 
26 namespace fastsim {
27  class Particle;
28  class ParticleFilter;
29  class SimplifiedGeometry;
30 
32 
40  {
41  public:
43 
53  const HepMC::GenEvent & genEvent,
54  const HepPDT::ParticleDataTable & particleDataTable,
55  double beamPipeRadius,
56  double deltaRchargedMother,
57  const ParticleFilter & particleFilter,
58  std::vector<SimTrack> & simTracks,
59  std::vector<SimVertex> & simVertices);
60 
62  ~ParticleManager();
63 
65 
72  std::unique_ptr<Particle> nextParticle(const RandomEngineAndDistribution & random);
73 
75 
81  void addSecondaries(
82  const math::XYZTLorentzVector & vertexPosition,
83  int motherSimTrackId,
84  std::vector<std::unique_ptr<Particle> > & secondaries,
85  const SimplifiedGeometry * layer = nullptr);
86 
88  const SimVertex getSimVertex(unsigned i) { return simVertices_->at(i); }
89 
91  const SimTrack getSimTrack(unsigned i) { return simTracks_->at(i); }
92 
94 
99  unsigned addEndVertex(const Particle * particle);
100 
101  private:
103 
109  unsigned addSimVertex(
111  int motherIndex);
112 
114 
121  unsigned addSimTrack(const Particle * particle);
122  void exoticRelativesChecker(const HepMC::GenVertex* originVertex, int& hasExoticAssociation, int ngendepth);
123 
125 
128  std::unique_ptr<Particle> nextGenParticle();
129 
130  // data members
131  const HepMC::GenEvent * const genEvent_;
132  HepMC::GenEvent::particle_const_iterator genParticleIterator_;
133  const HepMC::GenEvent::particle_const_iterator genParticleEnd_;
136  const double beamPipeRadius2_;
137  const double deltaRchargedMother_;
139  std::vector<SimTrack> * simTracks_;
140  std::vector<SimVertex> * simVertices_;
145  std::vector<std::unique_ptr<Particle> > particleBuffer_;
146  };
147 }
148 
149 inline bool isExotic(int pdgid_) {
150  unsigned int pdgid = std::abs(pdgid_);
151  return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) || // SUSY, R-hadron, and technicolor particles
152  pdgid == 17 || // 4th generation lepton
153  pdgid == 34 || // W-prime
154  pdgid == 37); // charged Higgs
155 }
156 
157 #endif
Implementation of a generic detector layer (base class for forward/barrel layers).
Manages GenParticles and Secondaries from interactions.
HepPDT::ParticleDataTable ParticleDataTable
const double deltaRchargedMother_
For FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter...
const ParticleFilter *const particleFilter_
(Kinematic) cuts on the particles that have to be propagated.
double momentumUnitConversionFactor_
Convert pythia units to GeV (FastSim standard)
std::vector< std::unique_ptr< Particle > > particleBuffer_
The vector of all secondaries that are not yet propagated in the event.
HepMC::GenEvent::particle_const_iterator genParticleIterator_
Iterator to keep track on which GenParticles where already considered.
TRandom random
Definition: MVATrainer.cc:138
const HepMC::GenEvent *const genEvent_
The GenEvent.
const HepMC::GenEvent::particle_const_iterator genParticleEnd_
The last particle of the GenEvent.
double lengthUnitConversionFactor2_
Convert pythia unis to cm^2 (FastSim standard)
const SimVertex getSimVertex(unsigned i)
Returns the position of a given SimVertex. Needed for interfacing the code with the old calorimetry...
const SimTrack getSimTrack(unsigned i)
Returns a given SimTrack. Needed for interfacing the code with the old calorimetry.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const double beamPipeRadius2_
(Radius of the beampipe)^2
bool isExotic(int pdgid_)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.
double lengthUnitConversionFactor_
Convert pythia unis to cm (FastSim standard)
const HepPDT::ParticleDataTable *const particleDataTable_
Necessary to get information like lifetime and charge of a particle if unknown.
double timeUnitConversionFactor_
Convert pythia unis to ns (FastSim standard)
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
(Kinematic) cuts on the particles that are propagated.
int genParticleIndex_
Index of particle in the GenEvent (if it is a GenParticle)