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 
13 // Author: L. Vanelderen, S. Kurz
14 // Date: 29 May 2017
16 
17 namespace HepPDT {
18  class ParticleDataTable;
19 }
20 
22 
23 namespace fastsim {
24  class Particle;
25  class ParticleFilter;
26  class SimplifiedGeometry;
27 
29 
37  public:
39 
49  const HepPDT::ParticleDataTable& particleDataTable,
50  double beamPipeRadius,
51  double deltaRchargedMother,
52  const ParticleFilter& particleFilter,
53  std::vector<SimTrack>& simTracks,
54  std::vector<SimVertex>& simVertices,
55  bool fixLongLivedBug,
56  bool useFastSimsDecayer);
57 
59  ~ParticleManager();
60 
62 
69  std::unique_ptr<Particle> nextParticle(const RandomEngineAndDistribution& random);
70 
72 
78  void addSecondaries(const math::XYZTLorentzVector& vertexPosition,
79  int motherSimTrackId,
80  std::vector<std::unique_ptr<Particle> >& secondaries,
81  const SimplifiedGeometry* layer = nullptr);
82 
84  const SimVertex getSimVertex(unsigned i) { return simVertices_->at(i); }
85 
87  const SimTrack getSimTrack(unsigned i) { return simTracks_->at(i); }
88 
90 
95  unsigned addEndVertex(const Particle* particle);
96 
97  private:
99 
105  unsigned addSimVertex(const math::XYZTLorentzVector& position, int motherIndex);
106 
108 
114  unsigned addSimTrack(const Particle* particle);
115  void exoticRelativesChecker(const HepMC::GenVertex* originVertex, int& hasExoticAssociation, int ngendepth);
116 
118 
121  std::unique_ptr<Particle> nextGenParticle();
122 
123  // data members
124  const HepMC::GenEvent* const genEvent_;
125  HepMC::GenEvent::particle_const_iterator
127  const HepMC::GenEvent::particle_const_iterator genParticleEnd_;
129  const HepPDT::ParticleDataTable* const
131  const double beamPipeRadius2_;
132  const double
135  std::vector<SimTrack>* simTracks_;
136  std::vector<SimVertex>* simVertices_;
143  std::vector<std::unique_ptr<Particle>>
145 
146  };
147 } // namespace fastsim
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  pdgid == 39); // bulk graviton
156 }
157 
158 #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 units 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 units 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 units 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)