Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:20

0001 #ifndef FASTSIM_PARTICLEMANAGER_H
0002 #define FASTSIM_PARTICLEMANAGER_H
0003 
0004 #include "DataFormats/Math/interface/LorentzVector.h"
0005 #include "HepMC/GenEvent.h"
0006 #include <vector>
0007 #include <memory>
0008 
0009 #include "SimDataFormats/Track/interface/SimTrack.h"
0010 #include "SimDataFormats/Vertex/interface/SimVertex.h"
0011 
0012 ///////////////////////////////////////////////
0013 // Author: L. Vanelderen, S. Kurz
0014 // Date: 29 May 2017
0015 //////////////////////////////////////////////////////////
0016 
0017 namespace HepPDT {
0018   class ParticleDataTable;
0019 }
0020 
0021 class RandomEngineAndDistribution;
0022 
0023 namespace fastsim {
0024   class Particle;
0025   class ParticleFilter;
0026   class SimplifiedGeometry;
0027 
0028   //! Manages GenParticles and Secondaries from interactions.
0029   /*!
0030         Manages which particle has to be propagated next, this includes GenParticles and secondaries from the interactions.
0031         Furthermore, checks if all necessary information is included with the GenParticles (charge, lifetime), otherwise 
0032         reads information from HepPDT::ParticleDataTable.
0033         Also handles secondaries, including closestChargedDaughter algorithm which is used for FastSim (cheat) tracking: a charged daughter can
0034         continue the track of a charged mother, see addSecondaries(...).
0035     */
0036   class ParticleManager {
0037   public:
0038     //! Constructor.
0039     /*!
0040             \param genEvent Get the GenEvent.
0041             \param particleDataTable Get information about particles, e.g. charge, lifetime.
0042             \param beamPipeRadius Radius of the beampipe.
0043             \param deltaRchargedMother For FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter.
0044             \param particleFilter Selects which particles have to be propagated.
0045             \param simTracks The SimTracks.
0046             \param simVertices The SimVertices.
0047         */
0048     ParticleManager(const HepMC::GenEvent& genEvent,
0049                     const HepPDT::ParticleDataTable& particleDataTable,
0050                     double beamPipeRadius,
0051                     double deltaRchargedMother,
0052                     const ParticleFilter& particleFilter,
0053                     std::vector<SimTrack>& simTracks,
0054                     std::vector<SimVertex>& simVertices,
0055                     bool useFastSimsDecayer);
0056 
0057     //! Default destructor.
0058     ~ParticleManager();
0059 
0060     //! Returns the next particle that has to be propagated (secondary or genParticle).
0061     /*!
0062             Main method of this class. At first consideres particles from the buffer (secondaries) if there are none, then the
0063             next GenParticle is considered. Only returns particles that pass the (kinetic) cuts of the ParticleFilter.
0064             Furthermore, uses the ParticleDataTable to ensure all necessary information about the particle is stored (lifetime, charge).
0065             \return The next particle that has to be propagated.
0066             \sa ParticleFilter
0067         */
0068     std::unique_ptr<Particle> nextParticle(const RandomEngineAndDistribution& random);
0069 
0070     //! Adds secondaries that are produced by any of the interactions (or particle decay) to the buffer.
0071     /*!
0072             Also checks which charged daughter is closest to a charged mother (in deltaR) and assigns the same SimTrack ID.
0073             \param vertexPosition The origin vertex (interaction or particle decay took place here).
0074             \param motherSimTrackId SimTrack ID of the mother particle, necessary for FastSim (cheat) tracking.
0075             \param secondaries All secondaries that where produced in a single particle decay or interaction.
0076         */
0077     void addSecondaries(const math::XYZTLorentzVector& vertexPosition,
0078                         int motherSimTrackId,
0079                         std::vector<std::unique_ptr<Particle> >& secondaries,
0080                         const SimplifiedGeometry* layer = nullptr);
0081 
0082     //! Returns the position of a given SimVertex. Needed for interfacing the code with the old calorimetry.
0083     const SimVertex getSimVertex(unsigned i) { return simVertices_->at(i); }
0084 
0085     //! Returns a given SimTrack. Needed for interfacing the code with the old calorimetry.
0086     const SimTrack getSimTrack(unsigned i) { return simTracks_->at(i); }
0087 
0088     //! Necessary to add an end vertex to a particle.
0089     /*!
0090             Needed if particle is no longer propagated for some reason (e.g. remaining energy below threshold) and no
0091             secondaries where produced at that point.
0092             \return Index of that vertex.
0093         */
0094     unsigned addEndVertex(const Particle* particle);
0095 
0096   private:
0097     //! Add a simVertex (simVertex contains information about the track it was produced).
0098     /*!
0099             Add a origin vertex for any particle.
0100             \param position Position of the vertex.
0101             \param motherIndex Index of the parent's simTrack.
0102             \return Index of that simVertex.
0103         */
0104     unsigned addSimVertex(const math::XYZTLorentzVector& position, int motherIndex);
0105 
0106     //! Add a simTrack (simTrack contains some basic info about the particle, e.g. pdgId).
0107     /*!
0108             Add a simTrack for a given particle and assign an index for that track. This might also be the index of the track
0109             of the mother particle (FastSim cheat tracking).
0110             \param particle Particle that produces that simTrack.
0111             \return Index of that simTrack.
0112         */
0113     unsigned addSimTrack(const Particle* particle);
0114     void exoticRelativesChecker(const HepMC::GenVertex* originVertex, int& hasExoticAssociation, int ngendepth);
0115 
0116     //! Returns next particle from the GenEvent that has to be propagated.
0117     /*!
0118             Tries to get some basic information about the status of the particle from the GenEvent and does some first rejection cuts based on them.
0119         */
0120     std::unique_ptr<Particle> nextGenParticle();
0121 
0122     // data members
0123     const HepMC::GenEvent* const genEvent_;  //!< The GenEvent
0124     HepMC::GenEvent::particle_const_iterator
0125         genParticleIterator_;  //!< Iterator to keep track on which GenParticles where already considered.
0126     const HepMC::GenEvent::particle_const_iterator genParticleEnd_;  //!< The last particle of the GenEvent.
0127     int genParticleIndex_;  //!< Index of particle in the GenEvent (if it is a GenParticle)
0128     const HepPDT::ParticleDataTable* const
0129         particleDataTable_;         //!< Necessary to get information like lifetime and charge of a particle if unknown.
0130     const double beamPipeRadius2_;  //!< (Radius of the beampipe)^2
0131     const double
0132         deltaRchargedMother_;  //!< For FastSim (cheat) tracking: cut on the angle between a charged mother and charged daughter.
0133     const ParticleFilter* const particleFilter_;  //!< (Kinematic) cuts on the particles that have to be propagated.
0134     std::vector<SimTrack>* simTracks_;            //!< The generated SimTrack of this event.
0135     std::vector<SimVertex>* simVertices_;         //!< The generated SimVertices of this event.
0136     bool useFastSimsDecayer_;
0137     double momentumUnitConversionFactor_;  //!< Convert pythia units to GeV (FastSim standard)
0138     double lengthUnitConversionFactor_;    //!< Convert pythia unis to cm (FastSim standard)
0139     double lengthUnitConversionFactor2_;   //!< Convert pythia unis to cm^2 (FastSim standard)
0140     double timeUnitConversionFactor_;      //!< Convert pythia unis to ns (FastSim standard)
0141     std::vector<std::unique_ptr<Particle> >
0142         particleBuffer_;  //!< The vector of all secondaries that are not yet propagated in the event.
0143   };
0144 }  // namespace fastsim
0145 
0146 inline bool isExotic(int pdgid_) {
0147   unsigned int pdgid = std::abs(pdgid_);
0148   return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) ||  // SUSY, R-hadron, and technicolor particles
0149           pdgid == 17 ||                                                // 4th generation lepton
0150           pdgid == 34 ||                                                // W-prime
0151           pdgid == 37 ||                                                // charged Higgs
0152           pdgid == 39);                                                 // bulk graviton
0153 }
0154 
0155 #endif