CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
fastsim::ParticleManager Class Reference

Manages GenParticles and Secondaries from interactions. More...

#include <ParticleManager.h>

Public Member Functions

unsigned addEndVertex (const Particle *particle)
 Necessary to add an end vertex to a particle. More...
 
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. More...
 
const SimVertex getSimVertex (unsigned i)
 Returns the position of a given SimVertex. Needed for interfacing the code with the old calorimetry. More...
 
std::unique_ptr< ParticlenextParticle (const RandomEngineAndDistribution &random)
 Returns the next particle that has to be propagated (secondary or genParticle). More...
 
 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. More...
 
 ~ParticleManager ()
 Default destructor. More...
 

Private Member Functions

unsigned addSimTrack (const Particle *particle)
 Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId). More...
 
unsigned addSimVertex (const math::XYZTLorentzVector &position, int motherIndex)
 Add a simVertex (simVertex contains information about the track it was produced). More...
 
std::unique_ptr< ParticlenextGenParticle ()
 Returns next particle from the GenEvent that has to be propagated. More...
 

Private Attributes

const double beamPipeRadius2_
 (Radius of the beampipe)^2 More...
 
const double deltaRchargedMother_
 For FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter. More...
 
const HepMC::GenEvent *const genEvent_
 The GenEvent. More...
 
const HepMC::GenEvent::particle_const_iterator genParticleEnd_
 The last particle of the GenEvent. More...
 
int genParticleIndex_
 Index of particle in the GenEvent (if it is a GenParticle) More...
 
HepMC::GenEvent::particle_const_iterator genParticleIterator_
 Iterator to keep track on which GenParticles where already considered. More...
 
double lengthUnitConversionFactor2_
 Convert pythia unis to cm^2 (FastSim standard) More...
 
double lengthUnitConversionFactor_
 Convert pythia unis to cm (FastSim standard) More...
 
double momentumUnitConversionFactor_
 Convert pythia units to GeV (FastSim standard) More...
 
std::vector< std::unique_ptr< Particle > > particleBuffer_
 The vector of all secondaries that are not yet propagated in the event. More...
 
const HepPDT::ParticleDataTable *const particleDataTable_
 Necessary to get information like lifetime and charge of a particle if unknown. More...
 
const ParticleFilter *const particleFilter_
 (Kinematic) cuts on the particles that have to be propagated. More...
 
std::vector< SimTrack > * simTracks_
 The generated SimTrack of this event. More...
 
std::vector< SimVertex > * simVertices_
 The generated SimVertices of this event. More...
 
double timeUnitConversionFactor_
 Convert pythia unis to ns (FastSim standard) More...
 

Detailed Description

Manages GenParticles and Secondaries from interactions.

Manages which particle has to be propagated next, this includes GenParticles and secondaries from the interactions. Furthermore, checks if all necessary information is included with the GenParticles (charge, lifetime), otherwise reads information from HepPDT::ParticleDataTable. Also handles secondaries, including closestChargedDaughter algorithm which is used for FastSim (cheat) tracking: a charged daughter can continue the track of a charged mother, see addSecondaries(...).

Definition at line 39 of file ParticleManager.h.

Constructor & Destructor Documentation

fastsim::ParticleManager::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.

Parameters
genEventGet the GenEvent.
particleDataTableGet information about particles, e.g. charge, lifetime.
beamPipeRadiusRadius of the beampipe.
deltaRchargedMotherFor FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter.
particleFilterSelects which particles have to be propagated.
simTracksThe SimTracks.
simVerticesThe SimVertices.

Definition at line 17 of file ParticleManager.cc.

References addSimVertex(), genEvent_, lengthUnitConversionFactor_, position, and timeUnitConversionFactor_.

26  , genParticleIterator_(genEvent_->particles_begin())
27  , genParticleEnd_(genEvent_->particles_end())
29  , particleDataTable_(&particleDataTable)
30  , beamPipeRadius2_(beamPipeRadius*beamPipeRadius)
31  , deltaRchargedMother_(deltaRchargedMother)
32  , particleFilter_(&particleFilter)
34  , simVertices_(&simVertices)
35  // prepare unit convsersions
36  // --------------------------------------------
37  // | | hepmc | cms |
38  // --------------------------------------------
39  // | length | genEvent_->length_unit | cm |
40  // | momentum | genEvent_->momentum_unit | GeV |
41  // | time | length unit (t*c) | ns |
42  // --------------------------------------------
43  , momentumUnitConversionFactor_(conversion_factor( genEvent_->momentum_unit(), HepMC::Units::GEV ))
44  , lengthUnitConversionFactor_(conversion_factor(genEvent_->length_unit(),HepMC::Units::LengthUnit::CM))
47 
48 {
49 
50  // add the main vertex from the signal event to the simvertex collection
51  if(genEvent.vertices_begin() != genEvent_->vertices_end())
52  {
53  const HepMC::FourVector & position = (*genEvent.vertices_begin())->position();
55  position.y()*lengthUnitConversionFactor_,
56  position.z()*lengthUnitConversionFactor_,
57  position.t()*timeUnitConversionFactor_)
58  ,-1);
59  }
60 }
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)
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
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)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const double beamPipeRadius2_
(Radius of the beampipe)^2
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.
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). ...
fastsim::ParticleManager::~ParticleManager ( )

Default destructor.

Definition at line 62 of file ParticleManager.cc.

62 {}

Member Function Documentation

unsigned fastsim::ParticleManager::addEndVertex ( const Particle particle)

Necessary to add an end vertex to a particle.

Needed if particle is no longer propagated for some reason (e.g. remaining energy below threshold) and no secondaries where produced at that point.

Returns
Index of that vertex.

Definition at line 196 of file ParticleManager.cc.

References addSimVertex(), fastsim::Particle::position(), and fastsim::Particle::simTrackIndex().

197 {
198  return this->addSimVertex(particle->position(), particle->simTrackIndex());
199 }
unsigned addSimVertex(const math::XYZTLorentzVector &position, int motherIndex)
Add a simVertex (simVertex contains information about the track it was produced). ...
void fastsim::ParticleManager::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.

Also checks which charged daughter is closest to a charged mother (in deltaR) and assigns the same SimTrack ID.

Parameters
vertexPositionThe origin vertex (interaction or particle decay took place here).
motherSimTrackIdSimTrack ID of the mother particle, necessary for FastSim (cheat) tracking.
secondariesAll secondaries that where produced in a single particle decay or interaction.

Definition at line 129 of file ParticleManager.cc.

References fastsim::ParticleFilter::acceptsVtx(), addSimVertex(), deltaRchargedMother_, training_settings::idx, fastsim::SimplifiedGeometry::index(), fastsim::SimplifiedGeometry::isForward(), eostools::move(), particleBuffer_, and particleFilter_.

Referenced by FastSimProducer::createFSimTrack().

134 {
135 
136  // vertex must be within the accepted volume
137  if(!particleFilter_->acceptsVtx(vertexPosition))
138  {
139  return;
140  }
141 
142  // no need to create vertex in case no particles are produced
143  if(secondaries.empty()){
144  return;
145  }
146 
147  // add simVertex
148  unsigned simVertexIndex = addSimVertex(vertexPosition,parentSimTrackIndex);
149 
150  // closest charged daughter continues the track of the mother particle
151  // simplified tracking algorithm for fastSim
152  double distMin = 99999.;
153  int idx = -1;
154  int idxMin = -1;
155  for(auto & secondary : secondaries)
156  {
157  idx++;
158  if(secondary->getMotherDeltaR() != -1){
159  if(secondary->getMotherDeltaR() > deltaRchargedMother_){
160  // larger than max requirement on deltaR
161  secondary->resetMother();
162  }else{
163  if(secondary->getMotherDeltaR() < distMin){
164  distMin = secondary->getMotherDeltaR();
165  idxMin = idx;
166  }
167  }
168  }
169  }
170 
171  // add secondaries to buffer
172  idx = -1;
173  for(auto & secondary : secondaries)
174  {
175  idx++;
176  if(idxMin != -1){
177  // reset all but the particle with the lowest deltaR (which is at idxMin)
178  if(secondary->getMotherDeltaR() != -1 && idx != idxMin){
179  secondary->resetMother();
180  }
181  }
182 
183  // set origin vertex
184  secondary->setSimVertexIndex(simVertexIndex);
185  //
186  if(layer)
187  {
188  secondary->setOnLayer(layer->isForward(), layer->index());
189  }
190  // ...and add particle to buffer
191  particleBuffer_.push_back(std::move(secondary));
192  }
193 
194 }
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.
std::vector< std::unique_ptr< Particle > > particleBuffer_
The vector of all secondaries that are not yet propagated in the event.
unsigned addSimVertex(const math::XYZTLorentzVector &position, int motherIndex)
Add a simVertex (simVertex contains information about the track it was produced). ...
def move(src, dest)
Definition: eostools.py:510
bool acceptsVtx(const math::XYZTLorentzVector &originVertexPosition) const
Vertex within tracker volume.
unsigned fastsim::ParticleManager::addSimTrack ( const Particle particle)
private

Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId).

Add a simTrack for a given particle and assign an index for that track. This might also be the index of the track of the mother particle (FastSim cheat tracking).

Parameters
particleParticle that produces that simTrack.
Returns
Index of that simTrack.

Definition at line 213 of file ParticleManager.cc.

References fastsim::Particle::genParticleIndex(), fastsim::Particle::getMotherDeltaR(), fastsim::Particle::getMotherSimTrackIndex(), fastsim::Particle::momentum(), fastsim::Particle::pdgId(), simTracks_, and fastsim::Particle::simVertexIndex().

Referenced by nextParticle().

214 {
215  int simTrackIndex;
216  // Again: FastSim cheat tracking -> continue track of mother
217  if(particle->getMotherDeltaR() != -1){
218  simTrackIndex = particle->getMotherSimTrackIndex();
219  }
220  // or create new SimTrack
221  else
222  {
223  simTrackIndex = simTracks_->size();
224  simTracks_->emplace_back(particle->pdgId(),
225  particle->momentum(),
226  particle->simVertexIndex(),
227  particle->genParticleIndex());
228  simTracks_->back().setTrackId(simTrackIndex);
229  }
230  return simTrackIndex;
231 }
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.
unsigned fastsim::ParticleManager::addSimVertex ( const math::XYZTLorentzVector position,
int  motherIndex 
)
private

Add a simVertex (simVertex contains information about the track it was produced).

Add a origin vertex for any particle.

Parameters
positionPosition of the vertex.
motherIndexIndex of the parent's simTrack.
Returns
Index of that simVertex.

Definition at line 201 of file ParticleManager.cc.

References simVertices_.

Referenced by addEndVertex(), addSecondaries(), and ParticleManager().

204 {
205  int simVertexIndex = simVertices_->size();
206  simVertices_->emplace_back(position.Vect(),
207  position.T(),
208  parentSimTrackIndex,
209  simVertexIndex);
210  return simVertexIndex;
211 }
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
const SimVertex fastsim::ParticleManager::getSimVertex ( unsigned  i)
inline

Returns the position of a given SimVertex. Needed for interfacing the code with the old calorimetry.

Definition at line 88 of file ParticleManager.h.

References objects.autophobj::motherIndex, and position.

Referenced by FastSimProducer::createFSimTrack().

88 { return simVertices_->at(i); }
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
std::unique_ptr< fastsim::Particle > fastsim::ParticleManager::nextGenParticle ( )
private

Returns next particle from the GenEvent that has to be propagated.

Tries to get some basic information about the status of the particle from the GenEvent and does some first rejection cuts based on them.

Definition at line 233 of file ParticleManager.cc.

References beamPipeRadius2_, GenParticle::GenParticle, genParticleEnd_, genParticleIndex_, genParticleIterator_, lengthUnitConversionFactor2_, lengthUnitConversionFactor_, momentumUnitConversionFactor_, fastsim::Constants::speedOfLight, and timeUnitConversionFactor_.

Referenced by nextParticle().

234 {
235  // only consider particles that start in the beam pipe and end outside the beam pipe
236  // try to get the decay time from pythia
237  // use hepmc units
238  // make the link simtrack to simvertex
239  // try not to change the simvertex structure
240 
241  // loop over gen particles
243  {
244  // some handy pointers and references
245  const HepMC::GenParticle & particle = **genParticleIterator_;
246  const HepMC::GenVertex * productionVertex = particle.production_vertex();
247  const HepMC::GenVertex * endVertex = particle.end_vertex();
248 
249  // skip incoming particles
250  if(!productionVertex){
251  continue;
252  }
253 
254  // particle must be produced within the beampipe
255  if(productionVertex->position().perp2()*lengthUnitConversionFactor2_ > beamPipeRadius2_)
256  {
257  continue;
258  }
259 
260  // particle must not decay before it reaches the beam pipe
261  if(endVertex && endVertex->position().perp2()*lengthUnitConversionFactor2_ < beamPipeRadius2_)
262  {
263  continue;
264  }
265 
266  // make the particle
267  std::unique_ptr<Particle> newParticle(
268  new Particle(particle.pdg_id(),
269  math::XYZTLorentzVector(productionVertex->position().x()*lengthUnitConversionFactor_,
270  productionVertex->position().y()*lengthUnitConversionFactor_,
271  productionVertex->position().z()*lengthUnitConversionFactor_,
272  productionVertex->position().t()*timeUnitConversionFactor_),
274  particle.momentum().y()*momentumUnitConversionFactor_,
275  particle.momentum().z()*momentumUnitConversionFactor_,
276  particle.momentum().e()*momentumUnitConversionFactor_)));
277  newParticle->setGenParticleIndex(genParticleIndex_);
278 
279  // try to get the life time of the particle from the genEvent
280  if(endVertex)
281  {
282  double labFrameLifeTime = (endVertex->position().t() - productionVertex->position().t())*timeUnitConversionFactor_;
283  newParticle->setRemainingProperLifeTimeC(labFrameLifeTime / newParticle->gamma() * fastsim::Constants::speedOfLight);
284  }
285 
286  // TODO: The products of a b-decay should point to that vertex and not to the primary vertex!
287  // Seems like this information has to be taken from the genEvent. How to do this? Is this really neccessary?
288  // See FBaseSimEvent::fill(..)
289  newParticle->setSimVertexIndex(0);
290 
291  // iterator/index has to be increased in case of return (is not done by the loop then)
293  // and return
294  return newParticle;
295  }
296 
297  return std::unique_ptr<Particle>();
298 }
double momentumUnitConversionFactor_
Convert pythia units to GeV (FastSim standard)
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
HepMC::GenEvent::particle_const_iterator genParticleIterator_
Iterator to keep track on which GenParticles where already considered.
const HepMC::GenEvent::particle_const_iterator genParticleEnd_
The last particle of the GenEvent.
double lengthUnitConversionFactor2_
Convert pythia unis to cm^2 (FastSim standard)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const double beamPipeRadius2_
(Radius of the beampipe)^2
double lengthUnitConversionFactor_
Convert pythia unis to cm (FastSim standard)
double timeUnitConversionFactor_
Convert pythia unis to ns (FastSim standard)
int genParticleIndex_
Index of particle in the GenEvent (if it is a GenParticle)
std::unique_ptr< fastsim::Particle > fastsim::ParticleManager::nextParticle ( const RandomEngineAndDistribution random)

Returns the next particle that has to be propagated (secondary or genParticle).

Main method of this class. At first consideres particles from the buffer (secondaries) if there are none, then the next GenParticle is considered. Only returns particles that pass the (kinetic) cuts of the ParticleFilter. Furthermore, uses the ParticleDataTable to ensure all necessary information about the particle is stored (lifetime, charge).

Returns
The next particle that has to be propagated.
See also
ParticleFilter

Definition at line 64 of file ParticleManager.cc.

References fastsim::ParticleFilter::accepts(), addSimTrack(), MillePedeFileConverter_cfg::e, Exception, RandomEngineAndDistribution::flatShoot(), cmsBatch::log, eostools::move(), nextGenParticle(), particleBuffer_, particleDataTable_, particleFilter_, source_particleGun_cfi::ParticleID, and ApeEstimator_cff::width.

65 {
66  std::unique_ptr<fastsim::Particle> particle;
67 
68  // retrieve particle from buffer
69  if(!particleBuffer_.empty())
70  {
71  particle = std::move(particleBuffer_.back());
72  particleBuffer_.pop_back();
73  }
74  // or from genParticle list
75  else
76  {
77  particle = nextGenParticle();
78  if(!particle) return nullptr;
79  }
80 
81  // if filter does not accept, skip particle
82  if(!particleFilter_->accepts(*particle))
83  {
84  return nextParticle(random);
85  }
86 
87  // lifetime or charge of particle are not yet set
88  if(!particle->remainingProperLifeTimeIsSet() || !particle->chargeIsSet())
89  {
90  // retrieve the particle data
91  const HepPDT::ParticleData * particleData = particleDataTable_->particle(HepPDT::ParticleID(particle->pdgId()));
92  if(!particleData)
93  {
94  throw cms::Exception("fastsim::ParticleManager") << "unknown pdg id: " << particle->pdgId() << std::endl;
95  }
96 
97  // set lifetime
98  if(!particle->remainingProperLifeTimeIsSet())
99  {
100  // The lifetime is 0. in the Pythia Particle Data Table! Calculate from width instead (ct=hbar/width).
101  // ct=particleData->lifetime().value();
102  double width = particleData->totalWidth().value();
103  if(width > 1.0e-35)
104  {
105  particle->setRemainingProperLifeTimeC(-log(random.flatShoot())*6.582119e-25/width/10.); // ct in cm
106  }
107  else
108  {
109  particle->setStable();
110  }
111  }
112 
113  // set charge
114  if(!particle->chargeIsSet())
115  {
116  particle->setCharge(particleData->charge());
117  }
118  }
119 
120  // add corresponding simTrack to simTrack collection
121  unsigned simTrackIndex = addSimTrack(particle.get());
122  particle->setSimTrackIndex(simTrackIndex);
123 
124  // and return
125  return particle;
126 }
double flatShoot(double xmin=0.0, double xmax=1.0) const
bool accepts(const Particle &particle) const
Check all if all criteria are fullfilled.
const ParticleFilter *const particleFilter_
(Kinematic) cuts on the particles that have to be propagated.
std::vector< std::unique_ptr< Particle > > particleBuffer_
The vector of all secondaries that are not yet propagated in the event.
unsigned addSimTrack(const Particle *particle)
Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId).
std::unique_ptr< Particle > nextParticle(const RandomEngineAndDistribution &random)
Returns the next particle that has to be propagated (secondary or genParticle).
HepPDT::ParticleData ParticleData
std::unique_ptr< Particle > nextGenParticle()
Returns next particle from the GenEvent that has to be propagated.
const HepPDT::ParticleDataTable *const particleDataTable_
Necessary to get information like lifetime and charge of a particle if unknown.
def move(src, dest)
Definition: eostools.py:510

Member Data Documentation

const double fastsim::ParticleManager::beamPipeRadius2_
private

(Radius of the beampipe)^2

Definition at line 131 of file ParticleManager.h.

Referenced by nextGenParticle().

const double fastsim::ParticleManager::deltaRchargedMother_
private

For FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter.

Definition at line 132 of file ParticleManager.h.

Referenced by addSecondaries().

const HepMC::GenEvent* const fastsim::ParticleManager::genEvent_
private

The GenEvent.

Definition at line 126 of file ParticleManager.h.

Referenced by ParticleManager().

const HepMC::GenEvent::particle_const_iterator fastsim::ParticleManager::genParticleEnd_
private

The last particle of the GenEvent.

Definition at line 128 of file ParticleManager.h.

Referenced by nextGenParticle().

int fastsim::ParticleManager::genParticleIndex_
private

Index of particle in the GenEvent (if it is a GenParticle)

Definition at line 129 of file ParticleManager.h.

Referenced by nextGenParticle().

HepMC::GenEvent::particle_const_iterator fastsim::ParticleManager::genParticleIterator_
private

Iterator to keep track on which GenParticles where already considered.

Definition at line 127 of file ParticleManager.h.

Referenced by nextGenParticle().

double fastsim::ParticleManager::lengthUnitConversionFactor2_
private

Convert pythia unis to cm^2 (FastSim standard)

Definition at line 138 of file ParticleManager.h.

Referenced by nextGenParticle().

double fastsim::ParticleManager::lengthUnitConversionFactor_
private

Convert pythia unis to cm (FastSim standard)

Definition at line 137 of file ParticleManager.h.

Referenced by nextGenParticle(), and ParticleManager().

double fastsim::ParticleManager::momentumUnitConversionFactor_
private

Convert pythia units to GeV (FastSim standard)

Definition at line 136 of file ParticleManager.h.

Referenced by nextGenParticle().

std::vector<std::unique_ptr<Particle> > fastsim::ParticleManager::particleBuffer_
private

The vector of all secondaries that are not yet propagated in the event.

Definition at line 140 of file ParticleManager.h.

Referenced by addSecondaries(), and nextParticle().

const HepPDT::ParticleDataTable* const fastsim::ParticleManager::particleDataTable_
private

Necessary to get information like lifetime and charge of a particle if unknown.

Definition at line 130 of file ParticleManager.h.

Referenced by nextParticle().

const ParticleFilter* const fastsim::ParticleManager::particleFilter_
private

(Kinematic) cuts on the particles that have to be propagated.

Definition at line 133 of file ParticleManager.h.

Referenced by addSecondaries(), and nextParticle().

std::vector<SimTrack>* fastsim::ParticleManager::simTracks_
private

The generated SimTrack of this event.

Definition at line 134 of file ParticleManager.h.

Referenced by addSimTrack().

std::vector<SimVertex>* fastsim::ParticleManager::simVertices_
private

The generated SimVertices of this event.

Definition at line 135 of file ParticleManager.h.

Referenced by addSimVertex().

double fastsim::ParticleManager::timeUnitConversionFactor_
private

Convert pythia unis to ns (FastSim standard)

Definition at line 139 of file ParticleManager.h.

Referenced by nextGenParticle(), and ParticleManager().