CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  : genEvent_(&genEvent),
25  genParticleIterator_(genEvent_->particles_begin()),
26  genParticleEnd_(genEvent_->particles_end()),
27  genParticleIndex_(1),
28  particleDataTable_(&particleDataTable),
29  beamPipeRadius2_(beamPipeRadius * beamPipeRadius),
30  deltaRchargedMother_(deltaRchargedMother),
31  particleFilter_(&particleFilter),
32  simTracks_(&simTracks),
33  simVertices_(&simVertices)
34  // prepare unit convsersions
35  // --------------------------------------------
36  // | | hepmc | cms |
37  // --------------------------------------------
38  // | length | genEvent_->length_unit | cm |
39  // | momentum | genEvent_->momentum_unit | GeV |
40  // | time | length unit (t*c) | ns |
41  // --------------------------------------------
42  ,
43  momentumUnitConversionFactor_(conversion_factor(genEvent_->momentum_unit(), HepMC::Units::GEV)),
44  lengthUnitConversionFactor_(conversion_factor(genEvent_->length_unit(), HepMC::Units::LengthUnit::CM)),
45  lengthUnitConversionFactor2_(lengthUnitConversionFactor_ * lengthUnitConversionFactor_),
46  timeUnitConversionFactor_(lengthUnitConversionFactor_ / fastsim::Constants::speedOfLight)
47 
48 {
49  // add the main vertex from the signal event to the simvertex collection
50  if (genEvent.vertices_begin() != genEvent_->vertices_end()) {
51  const HepMC::FourVector& position = (*genEvent.vertices_begin())->position();
53  position.y() * lengthUnitConversionFactor_,
54  position.z() * lengthUnitConversionFactor_,
55  position.t() * timeUnitConversionFactor_),
56  -1);
57  }
58 }
59 
61 
62 std::unique_ptr<fastsim::Particle> fastsim::ParticleManager::nextParticle(const RandomEngineAndDistribution& random) {
63  std::unique_ptr<fastsim::Particle> particle;
64 
65  // retrieve particle from buffer
66  if (!particleBuffer_.empty()) {
67  particle = std::move(particleBuffer_.back());
68  particleBuffer_.pop_back();
69  }
70  // or from genParticle list
71  else {
72  particle = nextGenParticle();
73  if (!particle)
74  return nullptr;
75  }
76 
77  // if filter does not accept, skip particle
78  if (!particleFilter_->accepts(*particle)) {
79  return nextParticle(random);
80  }
81 
82  // lifetime or charge of particle are not yet set
83  if (!particle->remainingProperLifeTimeIsSet() || !particle->chargeIsSet()) {
84  // retrieve the particle data
85  const HepPDT::ParticleData* particleData = particleDataTable_->particle(HepPDT::ParticleID(particle->pdgId()));
86  if (!particleData) {
87  // 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)
88  // they have short lifetimes, so decay them right away (charge and lifetime cannot be taken from table)
89  particle->setRemainingProperLifeTimeC(0.);
90  particle->setCharge(0.);
91  }
92 
93  // set lifetime
94  if (!particle->remainingProperLifeTimeIsSet()) {
95  // The lifetime is 0. in the Pythia Particle Data Table! Calculate from width instead (ct=hbar/width).
96  // ct=particleData->lifetime().value();
97  double width = particleData->totalWidth().value();
98  if (width > 1.0e-35) {
99  particle->setRemainingProperLifeTimeC(-log(random.flatShoot()) * 6.582119e-25 / width / 10.); // ct in cm
100  } else {
101  particle->setStable();
102  }
103  }
104 
105  // set charge
106  if (!particle->chargeIsSet()) {
107  particle->setCharge(particleData->charge());
108  }
109  }
110 
111  // add corresponding simTrack to simTrack collection
112  unsigned simTrackIndex = addSimTrack(particle.get());
113  particle->setSimTrackIndex(simTrackIndex);
114 
115  // and return
116  return particle;
117 }
118 
120  int parentSimTrackIndex,
121  std::vector<std::unique_ptr<Particle> >& secondaries,
122  const SimplifiedGeometry* layer) {
123  // vertex must be within the accepted volume
124  if (!particleFilter_->acceptsVtx(vertexPosition)) {
125  return;
126  }
127 
128  // no need to create vertex in case no particles are produced
129  if (secondaries.empty()) {
130  return;
131  }
132 
133  // add simVertex
134  unsigned simVertexIndex = addSimVertex(vertexPosition, parentSimTrackIndex);
135 
136  // closest charged daughter continues the track of the mother particle
137  // simplified tracking algorithm for fastSim
138  double distMin = 99999.;
139  int idx = -1;
140  int idxMin = -1;
141  for (auto& secondary : secondaries) {
142  idx++;
143  if (secondary->getMotherDeltaR() != -1) {
144  if (secondary->getMotherDeltaR() > deltaRchargedMother_) {
145  // larger than max requirement on deltaR
146  secondary->resetMother();
147  } else {
148  if (secondary->getMotherDeltaR() < distMin) {
149  distMin = secondary->getMotherDeltaR();
150  idxMin = idx;
151  }
152  }
153  }
154  }
155 
156  // add secondaries to buffer
157  idx = -1;
158  for (auto& secondary : secondaries) {
159  idx++;
160  if (idxMin != -1) {
161  // reset all but the particle with the lowest deltaR (which is at idxMin)
162  if (secondary->getMotherDeltaR() != -1 && idx != idxMin) {
163  secondary->resetMother();
164  }
165  }
166 
167  // set origin vertex
168  secondary->setSimVertexIndex(simVertexIndex);
169  //
170  if (layer) {
171  secondary->setOnLayer(layer->isForward(), layer->index());
172  }
173  // ...and add particle to buffer
174  particleBuffer_.push_back(std::move(secondary));
175  }
176 }
177 
179  return this->addSimVertex(particle->position(), particle->simTrackIndex());
180 }
181 
183  int simVertexIndex = simVertices_->size();
184  simVertices_->emplace_back(position.Vect(), position.T(), parentSimTrackIndex, simVertexIndex);
185  return simVertexIndex;
186 }
187 
189  int simTrackIndex = simTracks_->size();
190  simTracks_->emplace_back(
191  particle->pdgId(), particle->momentum(), particle->simVertexIndex(), particle->genParticleIndex());
192  simTracks_->back().setTrackId(simTrackIndex);
193  return simTrackIndex;
194 }
195 
196 std::unique_ptr<fastsim::Particle> fastsim::ParticleManager::nextGenParticle() {
197  // only consider particles that start in the beam pipe and end outside the beam pipe
198  // try to get the decay time from pythia
199  // use hepmc units
200  // make the link simtrack to simvertex
201  // try not to change the simvertex structure
202 
203  // loop over gen particles
204  for (; genParticleIterator_ != genParticleEnd_; ++genParticleIterator_, ++genParticleIndex_) {
205  // some handy pointers and references
206  const HepMC::GenParticle& particle = **genParticleIterator_;
207  const HepMC::GenVertex* productionVertex = particle.production_vertex();
208  const HepMC::GenVertex* endVertex = particle.end_vertex();
209 
210  // skip incoming particles
211  if (!productionVertex) {
212  continue;
213  }
214  if (std::abs(particle.pdg_id()) < 10 || std::abs(particle.pdg_id()) == 21) {
215  continue;
216  }
217  // particles which do not descend from exotics must be produced within the beampipe
218  int exoticRelativeId = 0;
219  const bool producedWithinBeamPipe =
220  productionVertex->position().perp2() * lengthUnitConversionFactor2_ < beamPipeRadius2_;
221  if (!producedWithinBeamPipe) {
222  exoticRelativesChecker(productionVertex, exoticRelativeId, 0);
223  if (!isExotic(exoticRelativeId)) {
224  continue;
225  }
226  }
227 
228  // FastSim will not make hits out of particles that decay before reaching the beam pipe
229  const bool decayedWithinBeamPipe =
230  endVertex && endVertex->position().perp2() * lengthUnitConversionFactor2_ < beamPipeRadius2_;
231  if (decayedWithinBeamPipe) {
232  continue;
233  }
234 
235  // SM particles that descend from exotics and cross the beam pipe radius should make hits but not be decayed
236  if (producedWithinBeamPipe && !decayedWithinBeamPipe) {
237  exoticRelativesChecker(productionVertex, exoticRelativeId, 0);
238  }
239 
240  // make the particle
241  std::unique_ptr<Particle> newParticle(
242  new Particle(particle.pdg_id(),
243  math::XYZTLorentzVector(productionVertex->position().x() * lengthUnitConversionFactor_,
244  productionVertex->position().y() * lengthUnitConversionFactor_,
245  productionVertex->position().z() * lengthUnitConversionFactor_,
246  productionVertex->position().t() * timeUnitConversionFactor_),
247  math::XYZTLorentzVector(particle.momentum().x() * momentumUnitConversionFactor_,
248  particle.momentum().y() * momentumUnitConversionFactor_,
249  particle.momentum().z() * momentumUnitConversionFactor_,
250  particle.momentum().e() * momentumUnitConversionFactor_)));
251  newParticle->setGenParticleIndex(genParticleIndex_);
252  if (isExotic(exoticRelativeId)) {
253  newParticle->setMotherPdgId(exoticRelativeId);
254  }
255 
256  // try to get the life time of the particle from the genEvent
257  if (endVertex) {
258  double labFrameLifeTime =
259  (endVertex->position().t() - productionVertex->position().t()) * timeUnitConversionFactor_;
260  newParticle->setRemainingProperLifeTimeC(labFrameLifeTime / newParticle->gamma() *
262  }
263 
264  // Find production vertex if it already exists. Otherwise create new vertex
265  // Possible to recreate the whole GenEvent using SimTracks/SimVertices (see FBaseSimEvent::fill(..))
266  bool foundVtx = false;
267  for (const auto& simVtx : *simVertices_) {
268  if (std::abs(simVtx.position().x() - newParticle->position().x()) < 1E-3 &&
269  std::abs(simVtx.position().y() - newParticle->position().y()) < 1E-3 &&
270  std::abs(simVtx.position().z() - newParticle->position().z()) < 1E-3) {
271  newParticle->setSimVertexIndex(simVtx.vertexId());
272  foundVtx = true;
273  break;
274  }
275  }
276  if (!foundVtx)
277  newParticle->setSimVertexIndex(addSimVertex(newParticle->position(), -1));
278 
279  // iterator/index has to be increased in case of return (is not done by the loop then)
280  ++genParticleIterator_;
281  ++genParticleIndex_;
282  // and return
283  return newParticle;
284  }
285 
286  return std::unique_ptr<Particle>();
287 }
288 
289 void fastsim::ParticleManager::exoticRelativesChecker(const HepMC::GenVertex* originVertex,
290  int& exoticRelativeId_,
291  int ngendepth = 0) {
292  if (ngendepth > 99 || exoticRelativeId_ == -1 || isExotic(std::abs(exoticRelativeId_)))
293  return;
294  ngendepth += 1;
295  std::vector<HepMC::GenParticle*>::const_iterator relativesIterator_ = originVertex->particles_in_const_begin();
296  std::vector<HepMC::GenParticle*>::const_iterator relativesIteratorEnd_ = originVertex->particles_in_const_end();
297  for (; relativesIterator_ != relativesIteratorEnd_; ++relativesIterator_) {
298  const HepMC::GenParticle& genRelative = **relativesIterator_;
299  if (isExotic(std::abs(genRelative.pdg_id()))) {
300  exoticRelativeId_ = genRelative.pdg_id();
301  if (ngendepth == 100)
302  exoticRelativeId_ = -1;
303  return;
304  }
305  const HepMC::GenVertex* vertex_ = genRelative.production_vertex();
306  if (!vertex_)
307  return;
308  exoticRelativesChecker(vertex_, exoticRelativeId_, ngendepth);
309  }
310  return;
311 }
static std::vector< std::string > checklist log
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:140
double flatShoot(double xmin=0.0, double xmax=1.0) const
HepPDT::ParticleDataTable ParticleDataTable
int genParticleIndex() const
Return index of the particle in the genParticle vector.
Definition: Particle.h:165
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:12
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
const HepMC::GenEvent *const genEvent_
The GenEvent.
unsigned addSimTrack(const Particle *particle)
Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId).
constexpr std::array< uint8_t, layerIndexSize > layer
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:159
bool isExotic(int pdgid_)
ParticleManager(const HepMC::GenEvent &genEvent, const HepPDT::ParticleDataTable &particleDataTable, double beamPipeRadius, double deltaRchargedMother, const ParticleFilter &particleFilter, std::vector< SimTrack > &simTracks, std::vector< SimVertex > &simVertices)
Constructor.
~ParticleManager()
Default destructor.
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
tuple genEvent
Definition: MCTruth.py:34
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:153
virtual bool isForward() const =0
Returns false/true depending if the object is a (non-abstract) barrel/forward layer.
std::unique_ptr< Particle > nextGenParticle()
Returns next particle from the GenEvent that has to be propagated.
double lengthUnitConversionFactor_
Convert pythia unis to cm (FastSim standard)
double timeUnitConversionFactor_
Convert pythia unis 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:289
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
(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:16
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.