1 #ifndef FASTSIM_PARTICLEMANAGER_H
2 #define FASTSIM_PARTICLEMANAGER_H
5 #include "HepMC/GenEvent.h"
26 class SimplifiedGeometry;
50 double beamPipeRadius,
51 double deltaRchargedMother,
53 std::vector<SimTrack>& simTracks,
54 std::vector<SimVertex>& simVertices);
78 std::vector<std::unique_ptr<Particle> >& secondaries,
113 void exoticRelativesChecker(
const HepMC::GenVertex* originVertex,
int& hasExoticAssociation,
int ngendepth);
123 HepMC::GenEvent::particle_const_iterator
139 std::vector<std::unique_ptr<Particle> >
145 unsigned int pdgid =
std::abs(pdgid_);
146 return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) ||
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).
constexpr std::array< uint8_t, layerIndexSize > layer
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.
const double beamPipeRadius2_
(Radius of the beampipe)^2
bool isExotic(int pdgid_)
ParticleManager(const HepMC::GenEvent &genEvent, const HepPDT::ParticleDataTable &particleDataTable, double beamPipeRadius, double deltaRchargedMother, const ParticleFilter &particleFilter, std::vector< SimTrack > &simTracks, std::vector< SimVertex > &simVertices)
Constructor.
~ParticleManager()
Default destructor.
Abs< T >::type abs(const T &t)
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]
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
(Kinematic) cuts on the particles that are propagated.
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
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.