CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ParticleManager.h
Go to the documentation of this file.
1 #ifndef FASTSIM_PARTICLEMANAGER_H
2 #define FASTSIM_PARTICLEMANAGER_H
3 
5 #include "HepMC/GenEvent.h"
6 #include <vector>
7 #include <memory>
8 
11 
13 // Author: L. Vanelderen, S. Kurz
14 // Date: 29 May 2017
16 
17 namespace HepPDT {
18  class ParticleDataTable;
19 }
20 
22 
23 namespace fastsim {
24  class Particle;
25  class ParticleFilter;
26  class SimplifiedGeometry;
27 
29 
37  public:
39 
49  const HepPDT::ParticleDataTable& particleDataTable,
50  double beamPipeRadius,
51  double deltaRchargedMother,
52  const ParticleFilter& particleFilter,
53  std::vector<SimTrack>& simTracks,
54  std::vector<SimVertex>& simVertices);
55 
58 
60 
67  std::unique_ptr<Particle> nextParticle(const RandomEngineAndDistribution& random);
68 
70 
76  void addSecondaries(const math::XYZTLorentzVector& vertexPosition,
77  int motherSimTrackId,
78  std::vector<std::unique_ptr<Particle> >& secondaries,
79  const SimplifiedGeometry* layer = nullptr);
80 
82  const SimVertex getSimVertex(unsigned i) { return simVertices_->at(i); }
83 
85  const SimTrack getSimTrack(unsigned i) { return simTracks_->at(i); }
86 
88 
93  unsigned addEndVertex(const Particle* particle);
94 
95  private:
97 
103  unsigned addSimVertex(const math::XYZTLorentzVector& position, int motherIndex);
104 
106 
112  unsigned addSimTrack(const Particle* particle);
113  void exoticRelativesChecker(const HepMC::GenVertex* originVertex, int& hasExoticAssociation, int ngendepth);
114 
116 
119  std::unique_ptr<Particle> nextGenParticle();
120 
121  // data members
122  const HepMC::GenEvent* const genEvent_;
123  HepMC::GenEvent::particle_const_iterator
125  const HepMC::GenEvent::particle_const_iterator genParticleEnd_;
127  const HepPDT::ParticleDataTable* const
129  const double beamPipeRadius2_;
130  const double
133  std::vector<SimTrack>* simTracks_;
134  std::vector<SimVertex>* simVertices_;
139  std::vector<std::unique_ptr<Particle> >
141  };
142 } // namespace fastsim
143 
144 inline bool isExotic(int pdgid_) {
145  unsigned int pdgid = std::abs(pdgid_);
146  return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) || // SUSY, R-hadron, and technicolor particles
147  pdgid == 17 || // 4th generation lepton
148  pdgid == 34 || // W-prime
149  pdgid == 37 || // charged Higgs
150  pdgid == 39); // bulk graviton
151 }
152 
153 #endif
Implementation of a generic detector layer (base class for forward/barrel layers).
Manages GenParticles and Secondaries from interactions.
HepPDT::ParticleDataTable ParticleDataTable
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)
std::vector< std::unique_ptr< Particle > > particleBuffer_
The vector of all secondaries that are not yet propagated in the event.
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)
const SimVertex getSimVertex(unsigned i)
Returns the position of a given SimVertex. Needed for interfacing the code with the old calorimetry...
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
const SimTrack getSimTrack(unsigned i)
Returns a given SimTrack. Needed for interfacing the code with the old calorimetry.
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_)
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.
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)
std::vector< SimTrack > * simTracks_
The generated SimTrack of this event.
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)
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.
(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
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.