Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef FastSimulation_Event_FBaseSimEvent_H
0002 #define FastSimulation_Event_FBaseSimEvent_H
0003 
0004 // Data Formats
0005 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
0006 #include "DataFormats/Math/interface/Point3D.h"
0007 
0008 // HepPDT Headers
0009 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
0010 
0011 // Famos Headers
0012 #include "CommonTools/BaseParticlePropagator/interface/RawParticle.h"
0013 #include "FastSimDataFormats/NuclearInteractions/interface/FSimVertexType.h"
0014 #include "FastSimDataFormats/NuclearInteractions/interface/FSimVertexTypeFwd.h"
0015 
0016 #include <vector>
0017 
0018 /** FSimEvent special features for FAMOS
0019  *
0020  * \author Patrick Janot, CERN
0021  * \date: 9-Dec-2003
0022  */
0023 
0024 //class FSimEvent;
0025 class FSimTrack;
0026 class FSimVertex;
0027 class KineParticleFilter;
0028 
0029 class SimTrack;
0030 class SimVertex;
0031 
0032 namespace edm {
0033   class ParameterSet;
0034 }
0035 
0036 namespace HepMC {
0037   class GenEvent;
0038   class GenParticle;
0039   class GenVertex;
0040 }  // namespace HepMC
0041 
0042 class FBaseSimEvent {
0043 public:
0044   /// Default constructor
0045   FBaseSimEvent(const edm::ParameterSet& kine);
0046 
0047   ///  usual virtual destructor
0048   ~FBaseSimEvent();
0049 
0050   /// Initialize the particle data table
0051   void initializePdt(const HepPDT::ParticleDataTable* aPdt);
0052 
0053   /// Get the pointer to the particle data table
0054   inline const HepPDT::ParticleDataTable* theTable() const { return pdt; }
0055 
0056   /// fill the FBaseSimEvent from the current HepMC::GenEvent
0057   void fill(const HepMC::GenEvent& hev);
0058 
0059   /// fill the FBaseSimEvent from SimTrack's and SimVert'ices
0060   void fill(const std::vector<SimTrack>&, const std::vector<SimVertex>&);
0061 
0062   /// print the original MCTruth event
0063   void printMCTruth(const HepMC::GenEvent& hev);
0064 
0065   /// Add the particles and their vertices to the list
0066   void addParticles(const HepMC::GenEvent& hev);
0067 
0068   /// print the FBaseSimEvent in an intelligible way
0069   void print() const;
0070 
0071   /// clear the FBaseSimEvent content before the next event
0072   void clear();
0073 
0074   /// Add an id in the vector of charged tracks id's
0075   void addChargedTrack(int id);
0076 
0077   /// Number of tracks
0078   inline unsigned int nTracks() const { return nSimTracks; }
0079 
0080   /// Number of vertices
0081   inline unsigned int nVertices() const { return nSimVertices; }
0082 
0083   /// Number of generator particles
0084   inline unsigned int nGenParts() const { return nGenParticles; }
0085 
0086   /// Number of "reconstructed" charged tracks
0087   inline unsigned int nChargedTracks() const { return nChargedParticleTracks; }
0088 
0089   /// Return track with given Id
0090   inline FSimTrack& track(int id) const;
0091 
0092   /// Return vertex with given Id
0093   inline FSimVertex& vertex(int id) const;
0094 
0095   /// Return vertex with given Id
0096   inline FSimVertexType& vertexType(int id) const;
0097 
0098   /// return "reconstructed" charged tracks index.
0099   int chargedTrack(int id) const;
0100 
0101   /// return embedded track with given id
0102   inline const SimTrack& embdTrack(int i) const;
0103 
0104   /// return embedded vertex with given id
0105   inline const SimVertex& embdVertex(int i) const;
0106 
0107   /// return embedded vertex type with given id
0108   inline const FSimVertexType& embdVertexType(int i) const;
0109 
0110   /// return MC track with a given id
0111   const HepMC::GenParticle* embdGenpart(int i) const;
0112 
0113   /// Add a new track to the Event and to the various lists
0114   int addSimTrack(const RawParticle* p, int iv, int ig = -1, const HepMC::GenVertex* ev = nullptr);
0115 
0116   /// Add a new vertex to the Event and to the various lists
0117   int addSimVertex(const XYZTLorentzVector& decayVertex,
0118                    int im = -1,
0119                    FSimVertexType::VertexType type = FSimVertexType::ANY);
0120 
0121   const KineParticleFilter& filter() const { return *myFilter; }
0122 
0123 protected:
0124   /// The pointer to the vector of FSimTrack's
0125   inline std::vector<FSimTrack>* tracks() const { return theSimTracks; }
0126 
0127   /// The pointer to the vector of FSimVertex's
0128   inline std::vector<FSimVertex>* vertices() const { return theSimVertices; }
0129 
0130   /// The pointer to the vector of GenParticle's
0131   inline std::vector<HepMC::GenParticle*>* genparts() const { return theGenParticles; }
0132 
0133 private:
0134   std::vector<FSimTrack>* theSimTracks;
0135   std::vector<FSimVertex>* theSimVertices;
0136   FSimVertexTypeCollection* theFSimVerticesType;
0137   std::vector<HepMC::GenParticle*>* theGenParticles;
0138 
0139   std::vector<unsigned>* theChargedTracks;
0140 
0141   unsigned int nSimTracks;
0142   unsigned int nSimVertices;
0143   unsigned int nGenParticles;
0144   unsigned int nChargedParticleTracks;
0145 
0146   unsigned int theTrackSize;
0147   unsigned int theVertexSize;
0148   unsigned int theGenSize;
0149   unsigned int theChargedSize;
0150   unsigned int initialSize;
0151 
0152   /// The particle filter
0153   KineParticleFilter* myFilter;
0154 
0155   double sigmaVerteX;
0156   double sigmaVerteY;
0157   double sigmaVerteZ;
0158 
0159   const ParticleDataTable* pdt;
0160 
0161   double lateVertexPosition;
0162 
0163   //  Histos* myHistos;
0164 };
0165 
0166 #include "FastSimulation/Event/interface/FBaseSimEvent.icc"
0167 
0168 #endif  // FBaseSimEvent_H