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 SimTrack getSimTrack (unsigned i)
 Returns a given SimTrack. Needed for interfacing the code with the old calorimetry. 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, bool useFastSimsDecayer)
 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...
 
void exoticRelativesChecker (const HepMC::GenVertex *originVertex, int &hasExoticAssociation, int ngendepth)
 
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...
 
bool useFastSimsDecayer_
 

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 36 of file ParticleManager.h.

Constructor & Destructor Documentation

◆ ParticleManager()

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,
bool  useFastSimsDecayer 
)

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(), genWeightsTable_cfi::genEvent, genEvent_, lengthUnitConversionFactor_, position, and timeUnitConversionFactor_.

25  : genEvent_(&genEvent),
26  genParticleIterator_(genEvent_->particles_begin()),
27  genParticleEnd_(genEvent_->particles_end()),
29  particleDataTable_(&particleDataTable),
36  // prepare unit convsersions
37  // --------------------------------------------
38  // | | hepmc | cms |
39  // --------------------------------------------
40  // | length | genEvent_->length_unit | cm |
41  // | momentum | genEvent_->momentum_unit | GeV |
42  // | time | length unit (t*c) | ns |
43  // --------------------------------------------
44  ,
45  momentumUnitConversionFactor_(conversion_factor(genEvent_->momentum_unit(), HepMC::Units::GEV)),
46  lengthUnitConversionFactor_(conversion_factor(genEvent_->length_unit(), HepMC::Units::LengthUnit::CM)),
49 
50 {
51  // add the main vertex from the signal event to the simvertex collection
52  if (genEvent.vertices_begin() != genEvent_->vertices_end()) {
53  const HepMC::FourVector& position = (*genEvent.vertices_begin())->position();
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:12
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:289
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). ...

◆ ~ParticleManager()

fastsim::ParticleManager::~ParticleManager ( )

Default destructor.

Definition at line 62 of file ParticleManager.cc.

62 {}

Member Function Documentation

◆ addEndVertex()

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 180 of file ParticleManager.cc.

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

180  {
181  return this->addSimVertex(particle->position(), particle->simTrackIndex());
182 }
unsigned addSimVertex(const math::XYZTLorentzVector &position, int motherIndex)
Add a simVertex (simVertex contains information about the track it was produced). ...

◆ addSecondaries()

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 121 of file ParticleManager.cc.

References muonTagProbeFilters_cff::distMin, heavyIonCSV_trainingSettings::idx, and eostools::move().

Referenced by FastSimProducer::createFSimTrack().

124  {
125  // vertex must be within the accepted volume
126  if (!particleFilter_->acceptsVtx(vertexPosition)) {
127  return;
128  }
129 
130  // no need to create vertex in case no particles are produced
131  if (secondaries.empty()) {
132  return;
133  }
134 
135  // add simVertex
136  unsigned simVertexIndex = addSimVertex(vertexPosition, parentSimTrackIndex);
137 
138  // closest charged daughter continues the track of the mother particle
139  // simplified tracking algorithm for fastSim
140  double distMin = 99999.;
141  int idx = -1;
142  int idxMin = -1;
143  for (auto& secondary : secondaries) {
144  idx++;
145  if (secondary->getMotherDeltaR() != -1) {
146  if (secondary->getMotherDeltaR() > deltaRchargedMother_) {
147  // larger than max requirement on deltaR
148  secondary->resetMother();
149  } else {
150  if (secondary->getMotherDeltaR() < distMin) {
151  distMin = secondary->getMotherDeltaR();
152  idxMin = idx;
153  }
154  }
155  }
156  }
157 
158  // add secondaries to buffer
159  idx = -1;
160  for (auto& secondary : secondaries) {
161  idx++;
162  if (idxMin != -1) {
163  // reset all but the particle with the lowest deltaR (which is at idxMin)
164  if (secondary->getMotherDeltaR() != -1 && idx != idxMin) {
165  secondary->resetMother();
166  }
167  }
168 
169  // set origin vertex
170  secondary->setSimVertexIndex(simVertexIndex);
171  //
172  if (layer) {
173  secondary->setOnLayer(layer->isForward(), layer->index());
174  }
175  // ...and add particle to buffer
176  particleBuffer_.push_back(std::move(secondary));
177  }
178 }
bool acceptsVtx(const math::XYZTLorentzVector &originVertexPosition) const
Vertex within tracker volume.
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.
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
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:511

◆ addSimTrack()

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 190 of file ParticleManager.cc.

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

190  {
191  int simTrackIndex = simTracks_->size();
192  simTracks_->emplace_back(
193  particle->pdgId(), particle->momentum(), particle->simVertexIndex(), particle->genParticleIndex());
194  simTracks_->back().setTrackId(simTrackIndex);
195  return simTrackIndex;
196 }
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.

◆ addSimVertex()

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 184 of file ParticleManager.cc.

References position.

Referenced by ParticleManager().

184  {
185  int simVertexIndex = simVertices_->size();
186  simVertices_->emplace_back(position.Vect(), position.T(), parentSimTrackIndex, simVertexIndex);
187  return simVertexIndex;
188 }
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.

◆ exoticRelativesChecker()

void fastsim::ParticleManager::exoticRelativesChecker ( const HepMC::GenVertex *  originVertex,
int &  hasExoticAssociation,
int  ngendepth = 0 
)
private

Definition at line 291 of file ParticleManager.cc.

References funct::abs(), GenParticle::GenParticle, and isExotic().

293  {
294  if (ngendepth > 99 || exoticRelativeId_ == -1 || isExotic(std::abs(exoticRelativeId_)))
295  return;
296  ngendepth += 1;
297  std::vector<HepMC::GenParticle*>::const_iterator relativesIterator_ = originVertex->particles_in_const_begin();
298  std::vector<HepMC::GenParticle*>::const_iterator relativesIteratorEnd_ = originVertex->particles_in_const_end();
299  for (; relativesIterator_ != relativesIteratorEnd_; ++relativesIterator_) {
300  const HepMC::GenParticle& genRelative = **relativesIterator_;
301  if (isExotic(std::abs(genRelative.pdg_id()))) {
302  exoticRelativeId_ = genRelative.pdg_id();
303  if (ngendepth == 100)
304  exoticRelativeId_ = -1;
305  return;
306  }
307  const HepMC::GenVertex* vertex_ = genRelative.production_vertex();
308  if (!vertex_)
309  return;
310  exoticRelativesChecker(vertex_, exoticRelativeId_, ngendepth);
311  }
312  return;
313 }
bool isExotic(int pdgid_)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void exoticRelativesChecker(const HepMC::GenVertex *originVertex, int &hasExoticAssociation, int ngendepth)

◆ getSimTrack()

const SimTrack fastsim::ParticleManager::getSimTrack ( unsigned  i)
inline

Returns a given SimTrack. Needed for interfacing the code with the old calorimetry.

Definition at line 86 of file ParticleManager.h.

References mps_fire::i, and simTracks_.

Referenced by FastSimProducer::createFSimTrack().

86 { return simTracks_->at(i); }
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.

◆ getSimVertex()

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 83 of file ParticleManager.h.

References mps_fire::i, and simVertices_.

Referenced by FastSimProducer::createFSimTrack().

83 { return simVertices_->at(i); }
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.

◆ nextGenParticle()

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 198 of file ParticleManager.cc.

References funct::abs(), GenParticle::GenParticle, isExotic(), and fastsim::Constants::speedOfLight.

198  {
199  // only consider particles that start in the beam pipe and end outside the beam pipe
200  // try to get the decay time from pythia
201  // use hepmc units
202  // make the link simtrack to simvertex
203  // try not to change the simvertex structure
204 
205  // loop over gen particles
207  // some handy pointers and references
208  const HepMC::GenParticle& particle = **genParticleIterator_;
209  const HepMC::GenVertex* productionVertex = particle.production_vertex();
210  const HepMC::GenVertex* endVertex = particle.end_vertex();
211 
212  // skip incoming particles
213  if (!productionVertex) {
214  continue;
215  }
216  if (std::abs(particle.pdg_id()) < 10 || std::abs(particle.pdg_id()) == 21) {
217  continue;
218  }
219  // particles which do not descend from exotics must be produced within the beampipe
220  int exoticRelativeId = 0;
221  const bool producedWithinBeamPipe =
222  productionVertex->position().perp2() * lengthUnitConversionFactor2_ < beamPipeRadius2_;
223  if (!producedWithinBeamPipe && useFastSimsDecayer_) {
224  exoticRelativesChecker(productionVertex, exoticRelativeId, 0);
225  if (!isExotic(exoticRelativeId)) {
226  continue;
227  }
228  }
229 
230  // FastSim will not make hits out of particles that decay before reaching the beam pipe
231  const bool decayedWithinBeamPipe =
232  endVertex && endVertex->position().perp2() * lengthUnitConversionFactor2_ < beamPipeRadius2_;
233  if (decayedWithinBeamPipe) {
234  continue;
235  }
236 
237  // SM particles that descend from exotics and cross the beam pipe radius should make hits but not be decayed
238  if (producedWithinBeamPipe && !decayedWithinBeamPipe && useFastSimsDecayer_) {
239  exoticRelativesChecker(productionVertex, exoticRelativeId, 0);
240  }
241 
242  // make the particle
243  std::unique_ptr<Particle> newParticle(
244  new Particle(particle.pdg_id(),
245  math::XYZTLorentzVector(productionVertex->position().x() * lengthUnitConversionFactor_,
246  productionVertex->position().y() * lengthUnitConversionFactor_,
247  productionVertex->position().z() * lengthUnitConversionFactor_,
248  productionVertex->position().t() * timeUnitConversionFactor_),
249  math::XYZTLorentzVector(particle.momentum().x() * momentumUnitConversionFactor_,
250  particle.momentum().y() * momentumUnitConversionFactor_,
251  particle.momentum().z() * momentumUnitConversionFactor_,
252  particle.momentum().e() * momentumUnitConversionFactor_)));
253  newParticle->setGenParticleIndex(genParticleIndex_);
254  if (isExotic(exoticRelativeId)) {
255  newParticle->setMotherPdgId(exoticRelativeId);
256  }
257 
258  // try to get the life time of the particle from the genEvent
259  if (endVertex) {
260  double labFrameLifeTime =
261  (endVertex->position().t() - productionVertex->position().t()) * timeUnitConversionFactor_;
262  newParticle->setRemainingProperLifeTimeC(labFrameLifeTime / newParticle->gamma() *
264  }
265 
266  // Find production vertex if it already exists. Otherwise create new vertex
267  // Possible to recreate the whole GenEvent using SimTracks/SimVertices (see FBaseSimEvent::fill(..))
268  bool foundVtx = false;
269  for (const auto& simVtx : *simVertices_) {
270  if (std::abs(simVtx.position().x() - newParticle->position().x()) < 1E-3 &&
271  std::abs(simVtx.position().y() - newParticle->position().y()) < 1E-3 &&
272  std::abs(simVtx.position().z() - newParticle->position().z()) < 1E-3) {
273  newParticle->setSimVertexIndex(simVtx.vertexId());
274  foundVtx = true;
275  break;
276  }
277  }
278  if (!foundVtx)
279  newParticle->setSimVertexIndex(addSimVertex(newParticle->position(), -1));
280 
281  // iterator/index has to be increased in case of return (is not done by the loop then)
284  // and return
285  return newParticle;
286  }
287 
288  return std::unique_ptr<Particle>();
289 }
double momentumUnitConversionFactor_
Convert pythia units to GeV (FastSim standard)
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:12
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
bool isExotic(int pdgid_)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void exoticRelativesChecker(const HepMC::GenVertex *originVertex, int &hasExoticAssociation, int ngendepth)
double lengthUnitConversionFactor_
Convert pythia unis to cm (FastSim standard)
double timeUnitConversionFactor_
Convert pythia unis to ns (FastSim standard)
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). ...

◆ nextParticle()

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 MillePedeFileConverter_cfg::e, RandomEngineAndDistribution::flatShoot(), dqm-mbProfile::log, eostools::move(), LHEGenericFilter_cfi::ParticleID, and ApeEstimator_cff::width.

64  {
65  std::unique_ptr<fastsim::Particle> particle;
66 
67  // retrieve particle from buffer
68  if (!particleBuffer_.empty()) {
69  particle = std::move(particleBuffer_.back());
70  particleBuffer_.pop_back();
71  }
72  // or from genParticle list
73  else {
74  particle = nextGenParticle();
75  if (!particle)
76  return nullptr;
77  }
78 
79  // if filter does not accept, skip particle
80  if (!particleFilter_->accepts(*particle)) {
81  return nextParticle(random);
82  }
83 
84  // lifetime or charge of particle are not yet set
85  if (!particle->remainingProperLifeTimeIsSet() || !particle->chargeIsSet()) {
86  // retrieve the particle data
87  const HepPDT::ParticleData* particleData = particleDataTable_->particle(HepPDT::ParticleID(particle->pdgId()));
88  if (!particleData) {
89  // in very few events the Decayer (pythia) produces high mass resonances that are for some reason not present in the table (even though they should technically be)
90  // they have short lifetimes, so decay them right away (charge and lifetime cannot be taken from table)
91  particle->setRemainingProperLifeTimeC(0.);
92  particle->setCharge(0.);
93  }
94 
95  // set lifetime
96  if (!particle->remainingProperLifeTimeIsSet()) {
97  // The lifetime is 0. in the Pythia Particle Data Table! Calculate from width instead (ct=hbar/width).
98  // ct=particleData->lifetime().value();
99  double width = particleData->totalWidth().value();
100  if (width > 1.0e-35) {
101  particle->setRemainingProperLifeTimeC(-log(random.flatShoot()) * 6.582119e-25 / width / 10.); // ct in cm
102  } else {
103  particle->setStable();
104  }
105  }
106 
107  // set charge
108  if (!particle->chargeIsSet()) {
109  particle->setCharge(particleData->charge());
110  }
111  }
112 
113  // add corresponding simTrack to simTrack collection
114  unsigned simTrackIndex = addSimTrack(particle.get());
115  particle->setSimTrackIndex(simTrackIndex);
116 
117  // and return
118  return particle;
119 }
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.
double flatShoot(double xmin=0.0, double xmax=1.0) const
def move(src, dest)
Definition: eostools.py:511
bool accepts(const Particle &particle) const
Check all if all criteria are fullfilled.

Member Data Documentation

◆ beamPipeRadius2_

const double fastsim::ParticleManager::beamPipeRadius2_
private

(Radius of the beampipe)^2

Definition at line 130 of file ParticleManager.h.

◆ deltaRchargedMother_

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.

◆ genEvent_

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

The GenEvent.

Definition at line 123 of file ParticleManager.h.

Referenced by ParticleManager().

◆ genParticleEnd_

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

The last particle of the GenEvent.

Definition at line 126 of file ParticleManager.h.

◆ genParticleIndex_

int fastsim::ParticleManager::genParticleIndex_
private

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

Definition at line 127 of file ParticleManager.h.

◆ genParticleIterator_

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

Iterator to keep track on which GenParticles where already considered.

Definition at line 125 of file ParticleManager.h.

◆ lengthUnitConversionFactor2_

double fastsim::ParticleManager::lengthUnitConversionFactor2_
private

Convert pythia unis to cm^2 (FastSim standard)

Definition at line 139 of file ParticleManager.h.

◆ lengthUnitConversionFactor_

double fastsim::ParticleManager::lengthUnitConversionFactor_
private

Convert pythia unis to cm (FastSim standard)

Definition at line 138 of file ParticleManager.h.

Referenced by ParticleManager().

◆ momentumUnitConversionFactor_

double fastsim::ParticleManager::momentumUnitConversionFactor_
private

Convert pythia units to GeV (FastSim standard)

Definition at line 137 of file ParticleManager.h.

◆ particleBuffer_

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 142 of file ParticleManager.h.

◆ particleDataTable_

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

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

Definition at line 129 of file ParticleManager.h.

◆ particleFilter_

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.

◆ simTracks_

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

The generated SimTrack of this event.

Definition at line 134 of file ParticleManager.h.

Referenced by getSimTrack().

◆ simVertices_

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

The generated SimVertices of this event.

Definition at line 135 of file ParticleManager.h.

Referenced by getSimVertex().

◆ timeUnitConversionFactor_

double fastsim::ParticleManager::timeUnitConversionFactor_
private

Convert pythia unis to ns (FastSim standard)

Definition at line 140 of file ParticleManager.h.

Referenced by ParticleManager().

◆ useFastSimsDecayer_

bool fastsim::ParticleManager::useFastSimsDecayer_
private

Definition at line 136 of file ParticleManager.h.