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 fixLongLivedBug, 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...
 
bool fixLongLivedBug_
 
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 units to cm^2 (FastSim standard). More...
 
double lengthUnitConversionFactor_
 Convert Pythia units 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 units 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

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

26  : genEvent_(&genEvent),
27  genParticleIterator_(genEvent_->particles_begin()),
28  genParticleEnd_(genEvent_->particles_end()),
30  particleDataTable_(&particleDataTable),
31  beamPipeRadius2_(beamPipeRadius * beamPipeRadius),
32  deltaRchargedMother_(deltaRchargedMother),
33  particleFilter_(&particleFilter),
36  fixLongLivedBug_(fixLongLivedBug),
37  useFastSimsDecayer_(useFastSimsDecayer)
38  // prepare unit convsersions
39  // --------------------------------------------
40  // | | hepmc | cms |
41  // --------------------------------------------
42  // | length | genEvent_->length_unit | cm |
43  // | momentum | genEvent_->momentum_unit | GeV |
44  // | time | length unit (t*c) | ns |
45  // --------------------------------------------
46  ,
47  momentumUnitConversionFactor_(conversion_factor(genEvent_->momentum_unit(), HepMC::Units::GEV)),
48  lengthUnitConversionFactor_(conversion_factor(genEvent_->length_unit(), HepMC::Units::LengthUnit::CM)),
51 
52 {
53  // add the main vertex from the signal event to the simvertex collection
54  if (genEvent.vertices_begin() != genEvent_->vertices_end()) {
55  const HepMC::FourVector& position = (*genEvent.vertices_begin())->position();
57  position.y() * lengthUnitConversionFactor_,
58  position.z() * lengthUnitConversionFactor_,
59  position.t() * timeUnitConversionFactor_),
60  -1);
61  }
62 }
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 units 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 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.
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 64 of file ParticleManager.cc.

64 {}

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

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

182  {
183  return this->addSimVertex(particle->position(), particle->simTrackIndex());
184 }
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 123 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().

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

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

Referenced by nextParticle().

192  {
193  int simTrackIndex = simTracks_->size();
194  simTracks_->emplace_back(
195  particle->pdgId(), particle->momentum(), particle->simVertexIndex(), particle->genParticleIndex());
196  simTracks_->back().setTrackId(simTrackIndex);
197  return simTrackIndex;
198 }
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 186 of file ParticleManager.cc.

References simVertices_.

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

186  {
187  int simVertexIndex = simVertices_->size();
188  simVertices_->emplace_back(position.Vect(), position.T(), parentSimTrackIndex, simVertexIndex);
189  return simVertexIndex;
190 }
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
void fastsim::ParticleManager::exoticRelativesChecker ( const HepMC::GenVertex *  originVertex,
int &  hasExoticAssociation,
int  ngendepth = 0 
)
private

Definition at line 293 of file ParticleManager.cc.

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

Referenced by nextGenParticle().

295  {
296  if (ngendepth > 99 || exoticRelativeId_ == -1 || isExotic(std::abs(exoticRelativeId_)))
297  return;
298  ngendepth += 1;
299  std::vector<HepMC::GenParticle*>::const_iterator relativesIterator_ = originVertex->particles_in_const_begin();
300  std::vector<HepMC::GenParticle*>::const_iterator relativesIteratorEnd_ = originVertex->particles_in_const_end();
301  for (; relativesIterator_ != relativesIteratorEnd_; ++relativesIterator_) {
302  const HepMC::GenParticle& genRelative = **relativesIterator_;
303  if (isExotic(std::abs(genRelative.pdg_id()))) {
304  exoticRelativeId_ = genRelative.pdg_id();
305  if (ngendepth == 100)
306  exoticRelativeId_ = -1;
307  return;
308  }
309  const HepMC::GenVertex* vertex_ = genRelative.production_vertex();
310  if (!vertex_)
311  return;
312  exoticRelativesChecker(vertex_, exoticRelativeId_, ngendepth);
313  }
314  return;
315 }
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)
const SimTrack fastsim::ParticleManager::getSimTrack ( unsigned  i)
inline

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

Definition at line 87 of file ParticleManager.h.

References objects.autophobj::motherIndex, and position.

Referenced by FastSimProducer::createFSimTrack().

87 { return simTracks_->at(i); }
std::vector< SimTrack > * simTracks_
The generated SimTrack 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 84 of file ParticleManager.h.

Referenced by FastSimProducer::createFSimTrack().

84 { 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 200 of file ParticleManager.cc.

References funct::abs(), addSimVertex(), beamPipeRadius2_, exoticRelativesChecker(), fixLongLivedBug_, GenParticle::GenParticle, genParticleEnd_, genParticleIndex_, genParticleIterator_, isExotic(), lengthUnitConversionFactor2_, lengthUnitConversionFactor_, momentumUnitConversionFactor_, simVertices_, fastsim::Constants::speedOfLight, timeUnitConversionFactor_, and useFastSimsDecayer_.

Referenced by nextParticle().

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

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

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

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

Referenced by addSecondaries().

bool fastsim::ParticleManager::fixLongLivedBug_
private

Definition at line 137 of file ParticleManager.h.

Referenced by nextGenParticle().

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

The GenEvent.

Definition at line 124 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 127 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 128 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 126 of file ParticleManager.h.

Referenced by nextGenParticle().

double fastsim::ParticleManager::lengthUnitConversionFactor2_
private

Convert Pythia units to cm^2 (FastSim standard).

Definition at line 141 of file ParticleManager.h.

Referenced by nextGenParticle().

double fastsim::ParticleManager::lengthUnitConversionFactor_
private

Convert Pythia units to cm (FastSim standard).

Definition at line 140 of file ParticleManager.h.

Referenced by nextGenParticle(), and ParticleManager().

double fastsim::ParticleManager::momentumUnitConversionFactor_
private

Convert Pythia units to GeV (FastSim standard).

Definition at line 139 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 144 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 134 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 135 of file ParticleManager.h.

Referenced by addSimTrack().

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

The generated SimVertices of this event.

Definition at line 136 of file ParticleManager.h.

Referenced by addSimVertex(), and nextGenParticle().

double fastsim::ParticleManager::timeUnitConversionFactor_
private

Convert Pythia units to ns (FastSim standard).

Definition at line 142 of file ParticleManager.h.

Referenced by nextGenParticle(), and ParticleManager().

bool fastsim::ParticleManager::useFastSimsDecayer_
private

Definition at line 138 of file ParticleManager.h.

Referenced by nextGenParticle().