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,
53  std::vector<SimTrack>& simTracks,
54  std::vector<SimVertex>& simVertices,
55  bool useFastSimsDecayer);
56 
59 
61 
68  std::unique_ptr<Particle> nextParticle(const RandomEngineAndDistribution& random);
69 
71 
77  void addSecondaries(const math::XYZTLorentzVector& vertexPosition,
78  int motherSimTrackId,
79  std::vector<std::unique_ptr<Particle> >& secondaries,
80  const SimplifiedGeometry* layer = nullptr);
81 
83  const SimVertex getSimVertex(unsigned i) { return simVertices_->at(i); }
84 
86  const SimTrack getSimTrack(unsigned i) { return simTracks_->at(i); }
87 
89 
94  unsigned addEndVertex(const Particle* particle);
95 
96  private:
98 
105 
107 
113  unsigned addSimTrack(const Particle* particle);
114  void exoticRelativesChecker(const HepMC::GenVertex* originVertex, int& hasExoticAssociation, int ngendepth);
115 
117 
120  std::unique_ptr<Particle> nextGenParticle();
121 
122  // data members
123  const HepMC::GenEvent* const genEvent_;
124  HepMC::GenEvent::particle_const_iterator
126  const HepMC::GenEvent::particle_const_iterator genParticleEnd_;
128  const HepPDT::ParticleDataTable* const
130  const double beamPipeRadius2_;
131  const double
134  std::vector<SimTrack>* simTracks_;
135  std::vector<SimVertex>* simVertices_;
141  std::vector<std::unique_ptr<Particle> >
143  };
144 } // namespace fastsim
145 
146 inline bool isExotic(int pdgid_) {
147  unsigned int pdgid = std::abs(pdgid_);
148  return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) || // SUSY, R-hadron, and technicolor particles
149  pdgid == 17 || // 4th generation lepton
150  pdgid == 34 || // W-prime
151  pdgid == 37 || // charged Higgs
152  pdgid == 39); // bulk graviton
153 }
154 
155 #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.
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...
unsigned addSimTrack(const Particle *particle)
Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId).
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_)
~ParticleManager()
Default destructor.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void addSecondaries(const math::XYZTLorentzVector &vertexPosition, int motherSimTrackId, std::vector< std::unique_ptr< Particle > > &secondaries, const SimplifiedGeometry *layer=nullptr)
Adds secondaries that are produced by any of the interactions (or particle decay) to the buffer...
std::unique_ptr< Particle > nextParticle(const RandomEngineAndDistribution &random)
Returns the next particle that has to be propagated (secondary or genParticle).
void exoticRelativesChecker(const HepMC::GenVertex *originVertex, int &hasExoticAssociation, int ngendepth)
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.
std::unique_ptr< Particle > nextGenParticle()
Returns next particle from the GenEvent that has to be propagated.
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:289
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
(Kinematic) cuts on the particles that are propagated.
ParticleManager(const HepMC::GenEvent &genEvent, const HepPDT::ParticleDataTable &particleDataTable, double beamPipeRadius, double deltaRchargedMother, const ParticleFilter &particleFilter, std::vector< SimTrack > &simTracks, std::vector< SimVertex > &simVertices, bool useFastSimsDecayer)
Constructor.
int genParticleIndex_
Index of particle in the GenEvent (if it is a GenParticle)
unsigned addSimVertex(const math::XYZTLorentzVector &position, int motherIndex)
Add a simVertex (simVertex contains information about the track it was produced). ...
unsigned addEndVertex(const Particle *particle)
Necessary to add an end vertex to a particle.