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

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  // 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  , momentumUnitConversionFactor_(conversion_factor( genEvent_->momentum_unit(), HepMC::Units::GEV ))
45  , lengthUnitConversionFactor_(conversion_factor(genEvent_->length_unit(),HepMC::Units::LengthUnit::CM))
48  , fixLongLivedBug_(fixLongLivedBug)
49 {
50 
51  // add the main vertex from the signal event to the simvertex collection
52  if(genEvent.vertices_begin() != genEvent_->vertices_end())
53  {
54  const HepMC::FourVector & position = (*genEvent.vertices_begin())->position();
56  position.y()*lengthUnitConversionFactor_,
57  position.z()*lengthUnitConversionFactor_,
58  position.t()*timeUnitConversionFactor_)
59  ,-1);
60  }
61 }
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 63 of file ParticleManager.cc.

63 {}

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

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

201 {
202  return this->addSimVertex(particle->position(), particle->simTrackIndex());
203 }
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 133 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().

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

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

Referenced by nextParticle().

219 {
220  int simTrackIndex = simTracks_->size();
221  simTracks_->emplace_back(particle->pdgId(),
222  particle->momentum(),
223  particle->simVertexIndex(),
224  particle->genParticleIndex());
225  simTracks_->back().setTrackId(simTrackIndex);
226  return simTrackIndex;
227 }
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 205 of file ParticleManager.cc.

References simVertices_.

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

208 {
209  int simVertexIndex = simVertices_->size();
210  simVertices_->emplace_back(position.Vect(),
211  position.T(),
212  parentSimTrackIndex,
213  simVertexIndex);
214  return simVertexIndex;
215 }
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 318 of file ParticleManager.cc.

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

Referenced by nextGenParticle().

320  {
321  if (ngendepth > 99 || exoticRelativeId_ == -1 || isExotic(fixLongLivedBug_, std::abs(exoticRelativeId_)))
322  return;
323  ngendepth += 1;
324  std::vector<HepMC::GenParticle*>::const_iterator relativesIterator_ = originVertex->particles_in_const_begin();
325  std::vector<HepMC::GenParticle*>::const_iterator relativesIteratorEnd_ = originVertex->particles_in_const_end();
326  for (; relativesIterator_ != relativesIteratorEnd_; ++relativesIterator_) {
327  const HepMC::GenParticle& genRelative = **relativesIterator_;
328  if (isExotic(fixLongLivedBug_, std::abs(genRelative.pdg_id()))) {
329  exoticRelativeId_ = genRelative.pdg_id();
330  if (ngendepth == 100)
331  exoticRelativeId_ = -1;
332  return;
333  }
334  const HepMC::GenVertex* vertex_ = genRelative.production_vertex();
335  if (!vertex_)
336  return;
337  exoticRelativesChecker(vertex_, exoticRelativeId_, ngendepth);
338  }
339  return;
340 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isExotic(bool fixLongLivedBug, int pdgid_)
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 92 of file ParticleManager.h.

References objects.autophobj::motherIndex, and position.

Referenced by FastSimProducer::createFSimTrack().

92 { 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 89 of file ParticleManager.h.

Referenced by FastSimProducer::createFSimTrack().

89 { 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 229 of file ParticleManager.cc.

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

Referenced by nextParticle().

230 {
231  // only consider particles that start in the beam pipe and end outside the beam pipe
232  // try to get the decay time from pythia
233  // use hepmc units
234  // make the link simtrack to simvertex
235  // try not to change the simvertex structure
236 
237  // loop over gen particles
239  {
240  // some handy pointers and references
241  const HepMC::GenParticle & particle = **genParticleIterator_;
242  const HepMC::GenVertex * productionVertex = particle.production_vertex();
243  const HepMC::GenVertex * endVertex = particle.end_vertex();
244  // skip incoming particles
245  if(!productionVertex){
246  continue;
247  }
248 
249  if (std::abs(particle.pdg_id()) < 10 || std::abs(particle.pdg_id()) == 21) {
250  continue;
251  }
252 
253  // particles which do not descend from exotics must be produced within the beampipe
254  int exoticRelativeId = 0;
255  const bool producedWithinBeamPipe = productionVertex->position().perp2() * lengthUnitConversionFactor2_ < beamPipeRadius2_;
256  if (!producedWithinBeamPipe) //
257  {
258  exoticRelativesChecker(productionVertex, exoticRelativeId, 0);
259  if (!isExotic(fixLongLivedBug_, exoticRelativeId)) {
260  continue;
261  }
262  }
263  const bool decayedWithinBeamPipe = endVertex && endVertex->position().perp2() * lengthUnitConversionFactor2_ < beamPipeRadius2_;
264  // FastSim will not make hits out of particles that decay before reaching the beam pipe
265  if(decayedWithinBeamPipe)
266  {
267  continue;
268  }
269 
270  // SM particles that descend from exotics and cross the beam pipe radius should make hits but not be decayed, by default it will duplicate FastSim hits for long lived particles and so anything produced without activating fixLongLivedBug_ is physically wrong
271  if (fixLongLivedBug_ && producedWithinBeamPipe && !decayedWithinBeamPipe){
272  exoticRelativesChecker(productionVertex, exoticRelativeId, 0);
273  }
274 
275  // make the particle
276  std::unique_ptr<Particle> newParticle(
277  new Particle(particle.pdg_id(),
278  math::XYZTLorentzVector(productionVertex->position().x()*lengthUnitConversionFactor_,
279  productionVertex->position().y()*lengthUnitConversionFactor_,
280  productionVertex->position().z()*lengthUnitConversionFactor_,
281  productionVertex->position().t()*timeUnitConversionFactor_),
283  particle.momentum().y()*momentumUnitConversionFactor_,
284  particle.momentum().z()*momentumUnitConversionFactor_,
285  particle.momentum().e()*momentumUnitConversionFactor_)));
286  newParticle->setGenParticleIndex(genParticleIndex_);
287  if (isExotic(fixLongLivedBug_, exoticRelativeId)) {
288  newParticle->setMotherPdgId(exoticRelativeId);
289  }
290  // try to get the life time of the particle from the genEvent
291  if(endVertex)
292  {
293  double labFrameLifeTime = (endVertex->position().t() - productionVertex->position().t())*timeUnitConversionFactor_;
294  newParticle->setRemainingProperLifeTimeC(labFrameLifeTime / newParticle->gamma() * fastsim::Constants::speedOfLight);
295  }
296 
297  // Find production vertex if it already exists. Otherwise create new vertex
298  // Possible to recreate the whole GenEvent using SimTracks/SimVertices (see FBaseSimEvent::fill(..))
299  bool foundVtx = false;
300  for(const auto& simVtx : *simVertices_){
301  if(std::abs(simVtx.position().x() - newParticle->position().x()) < 1E-3 && std::abs(simVtx.position().y() - newParticle->position().y()) < 1E-3 && std::abs(simVtx.position().z() - newParticle->position().z()) < 1E-3){
302  newParticle->setSimVertexIndex(simVtx.vertexId());
303  foundVtx = true;
304  break;
305  }
306  }
307  if(!foundVtx) newParticle->setSimVertexIndex(addSimVertex(newParticle->position(), -1));
308 
309  // iterator/index has to be increased in case of return (is not done by the loop then)
311  // and return
312  return newParticle;
313  }
314 
315  return std::unique_ptr<Particle>();
316 }
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isExotic(bool fixLongLivedBug, int pdgid_)
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). ...
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 65 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  {
72  particle = std::move(particleBuffer_.back());
73  particleBuffer_.pop_back();
74  }
75  // or from genParticle list
76  else
77  {
78  particle = nextGenParticle();
79  if(!particle) return nullptr;
80  }
81 
82  // if filter does not accept, skip particle
83  if(!particleFilter_->accepts(*particle))
84  {
85  return nextParticle(random);
86  }
87 
88  // lifetime or charge of particle are not yet set
89  if(!particle->remainingProperLifeTimeIsSet() || !particle->chargeIsSet())
90  {
91  // retrieve the particle data
92  const HepPDT::ParticleData * particleData = particleDataTable_->particle(HepPDT::ParticleID(particle->pdgId()));
93  if(!particleData)
94  {
95  // 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)
96  // they have short lifetimes, so decay them right away (charge and lifetime cannot be taken from table)
97  particle->setRemainingProperLifeTimeC(0.);
98  particle->setCharge(0.);
99  }
100 
101  // set lifetime
102  if(!particle->remainingProperLifeTimeIsSet())
103  {
104  // The lifetime is 0. in the Pythia Particle Data Table! Calculate from width instead (ct=hbar/width).
105  // ct=particleData->lifetime().value();
106  double width = particleData->totalWidth().value();
107  if(width > 1.0e-35)
108  {
109  particle->setRemainingProperLifeTimeC(-log(random.flatShoot())*6.582119e-25/width/10.); // ct in cm
110  }
111  else
112  {
113  particle->setStable();
114  }
115  }
116 
117  // set charge
118  if(!particle->chargeIsSet())
119  {
120  particle->setCharge(particleData->charge());
121  }
122  }
123 
124  // add corresponding simTrack to simTrack collection
125  unsigned simTrackIndex = addSimTrack(particle.get());
126  particle->setSimTrackIndex(simTrackIndex);
127 
128  // and return
129  return particle;
130 }
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 137 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 138 of file ParticleManager.h.

Referenced by addSecondaries().

bool fastsim::ParticleManager::fixLongLivedBug_
private

Definition at line 147 of file ParticleManager.h.

Referenced by exoticRelativesChecker(), and nextGenParticle().

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

The GenEvent.

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

Referenced by nextGenParticle().

double fastsim::ParticleManager::lengthUnitConversionFactor2_
private

Convert pythia unis to cm^2 (FastSim standard)

Definition at line 144 of file ParticleManager.h.

Referenced by nextGenParticle().

double fastsim::ParticleManager::lengthUnitConversionFactor_
private

Convert pythia unis to cm (FastSim standard)

Definition at line 143 of file ParticleManager.h.

Referenced by nextGenParticle(), and ParticleManager().

double fastsim::ParticleManager::momentumUnitConversionFactor_
private

Convert pythia units to GeV (FastSim standard)

Definition at line 142 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 146 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 136 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 139 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 140 of file ParticleManager.h.

Referenced by addSimTrack().

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

The generated SimVertices of this event.

Definition at line 141 of file ParticleManager.h.

Referenced by addSimVertex(), and nextGenParticle().

double fastsim::ParticleManager::timeUnitConversionFactor_
private

Convert pythia unis to ns (FastSim standard)

Definition at line 145 of file ParticleManager.h.

Referenced by nextGenParticle(), and ParticleManager().