CMS 3D CMS Logo

ParticleManager.cc
Go to the documentation of this file.
2 
3 #include "HepMC/GenEvent.h"
4 #include "HepMC/Units.h"
5 #include "HepPDT/ParticleDataTable.hh"
6 
11 
16 
18  const HepPDT::ParticleDataTable& particleDataTable,
19  double beamPipeRadius,
20  double deltaRchargedMother,
21  const fastsim::ParticleFilter& particleFilter,
22  std::vector<SimTrack>& simTracks,
23  std::vector<SimVertex>& simVertices,
24  bool fixLongLivedBug,
25  bool useFastSimsDecayer)
26  : genEvent_(&genEvent),
27  genParticleIterator_(genEvent_->particles_begin()),
28  genParticleEnd_(genEvent_->particles_end()),
29  genParticleIndex_(1),
30  particleDataTable_(&particleDataTable),
31  beamPipeRadius2_(beamPipeRadius * beamPipeRadius),
32  deltaRchargedMother_(deltaRchargedMother),
33  particleFilter_(&particleFilter),
34  simTracks_(&simTracks),
35  simVertices_(&simVertices),
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)),
49  lengthUnitConversionFactor2_(lengthUnitConversionFactor_ * lengthUnitConversionFactor_),
50  timeUnitConversionFactor_(lengthUnitConversionFactor_ / fastsim::Constants::speedOfLight)
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 }
63 
65 
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 }
122 
124  int parentSimTrackIndex,
125  std::vector<std::unique_ptr<Particle> >& secondaries,
126  const SimplifiedGeometry* layer) {
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 }
181 
183  return this->addSimVertex(particle->position(), particle->simTrackIndex());
184 }
185 
187  int simVertexIndex = simVertices_->size();
188  simVertices_->emplace_back(position.Vect(), position.T(), parentSimTrackIndex, simVertexIndex);
189  return simVertexIndex;
190 }
191 
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 }
199 
200 std::unique_ptr<fastsim::Particle> fastsim::ParticleManager::nextGenParticle() {
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 }
292 
293 void fastsim::ParticleManager::exoticRelativesChecker(const HepMC::GenVertex* originVertex,
294  int& exoticRelativeId_,
295  int ngendepth = 0) {
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 }
Implementation of a generic detector layer (base class for forward/barrel layers).
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:142
double flatShoot(double xmin=0.0, double xmax=1.0) const
HepPDT::ParticleDataTable ParticleDataTable
const double deltaRchargedMother_
For FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter...
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.
int genParticleIndex() const
Return index of the particle in the genParticle vector.
Definition: Particle.h:167
double momentumUnitConversionFactor_
Convert Pythia units to GeV (FastSim standard).
std::vector< std::unique_ptr< Particle > > particleBuffer_
The vector of all secondaries that are not yet propagated in the event.
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
virtual bool isForward() const =0
Returns false/true depending if the object is a (non-abstract) barrel/forward layer.
HepMC::GenEvent::particle_const_iterator genParticleIterator_
Iterator to keep track on which GenParticles where already considered.
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:136
TRandom random
Definition: MVATrainer.cc:138
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).
unsigned addSimTrack(const Particle *particle)
Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId).
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int simVertexIndex() const
Return index of the origin vertex.
Definition: Particle.h:161
const double beamPipeRadius2_
(Radius of the beampipe)^2
bool isExotic(int pdgid_)
~ParticleManager()
Default destructor.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void addSecondaries(const math::XYZTLorentzVector &vertexPosition, int motherSimTrackId, std::vector< std::unique_ptr< Particle > > &secondaries, const SimplifiedGeometry *layer=nullptr)
Adds secondaries that are produced by any of the interactions (or particle decay) to the buffer...
std::unique_ptr< Particle > nextParticle(const RandomEngineAndDistribution &random)
Returns the next particle that has to be propagated (secondary or genParticle).
void exoticRelativesChecker(const HepMC::GenVertex *originVertex, int &hasExoticAssociation, int ngendepth)
HepPDT::ParticleData ParticleData
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:155
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.
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.
std::unique_ptr< Particle > nextGenParticle()
Returns next particle from the GenEvent that has to be propagated.
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).
int index() const
Return index of this layer (layers are numbered according to their position in the detector)...
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< SimVertex > * simVertices_
The generated SimVertices of this event.
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:145
(Kinematic) cuts on the particles that are propagated.
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:19
int genParticleIndex_
Index of particle in the GenEvent (if it is a GenParticle)
unsigned addSimVertex(const math::XYZTLorentzVector &position, int motherIndex)
Add a simVertex (simVertex contains information about the track it was produced). ...
unsigned addEndVertex(const Particle *particle)
Necessary to add an end vertex to a particle.
def move(src, dest)
Definition: eostools.py:511
bool acceptsVtx(const math::XYZTLorentzVector &originVertexPosition) const
Vertex within tracker volume.