Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-05 23:51:51

0001 #ifndef SimTracker_TrackHistory_VertexClassifier_h
0002 #define SimTracker_TrackHistory_VertexClassifier_h
0003 
0004 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
0010 #include "SimTracker/TrackHistory/interface/CMSProcessTypes.h"
0011 #include "SimTracker/TrackHistory/interface/VertexCategories.h"
0012 #include "SimTracker/TrackHistory/interface/VertexHistory.h"
0013 
0014 //! Get track history and classify it in function of their .
0015 class VertexClassifier : public VertexCategories {
0016 public:
0017   //! Type to the associate category
0018   typedef VertexCategories Categories;
0019 
0020   //! Constructor by ParameterSet
0021   VertexClassifier(edm::ParameterSet const &pset, edm::ConsumesCollector);
0022 
0023   virtual ~VertexClassifier() {}
0024 
0025   //! Pre-process event information (for accessing reconstraction information)
0026   virtual void newEvent(edm::Event const &, edm::EventSetup const &);
0027 
0028   //! Classify the RecoVertex in categories.
0029   VertexClassifier const &evaluate(reco::VertexBaseRef const &);
0030 
0031   //! Classify the TrackingVertex in categories.
0032   VertexClassifier const &evaluate(TrackingVertexRef const &);
0033 
0034   //! Classify the RecoVertex in categories.
0035   VertexClassifier const &evaluate(reco::VertexRef const &vertex) { return evaluate(reco::VertexBaseRef(vertex)); }
0036 
0037   //! Returns a reference to the vertex history used in the classification.
0038   VertexHistory const &history() const { return tracer_; }
0039 
0040   static void fillPSetDescription(edm::ParameterSetDescription &desc);
0041 
0042 private:
0043   VertexHistory tracer_;
0044 
0045   const G4toCMSLegacyProcTypeMap g4toCMSProcMap_;
0046 
0047   const edm::InputTag hepMCLabel_;
0048 
0049   double longLivedDecayLength_;
0050   double vertexClusteringDistance_;
0051 
0052   edm::Handle<edm::HepMCProduct> mcInformation_;
0053 
0054   edm::ESHandle<ParticleDataTable> particleDataTable_;
0055   edm::ESGetToken<ParticleDataTable, edm::DefaultRecord> particleDataTableToken_;
0056 
0057   //! Get reconstruction information
0058   void reconstructionInformation(reco::TrackBaseRef const &);
0059 
0060   //! Get all the information related to the simulation details
0061   void simulationInformation();
0062 
0063   //! Get all the information related to decay process
0064   void processesAtGenerator();
0065 
0066   //! Get information about conversion and other interactions
0067   void processesAtSimulation();
0068 
0069   //! Get geometrical information about the vertices
0070   void vertexInformation();
0071 
0072   //! Auxiliary class holding simulated primary vertices
0073   struct GeneratedPrimaryVertex {
0074     GeneratedPrimaryVertex(double x1, double y1, double z1) : x(x1), y(y1), z(z1), ptsq(0), nGenTrk(0) {}
0075 
0076     bool operator<(GeneratedPrimaryVertex const &reference) const { return ptsq < reference.ptsq; }
0077 
0078     double x, y, z;
0079     double ptsq;
0080     int nGenTrk;
0081 
0082     HepMC::FourVector ptot;
0083 
0084     std::vector<int> finalstateParticles;
0085     std::vector<int> simTrackIndex;
0086     std::vector<int> genVertex;
0087   };
0088 
0089   std::vector<GeneratedPrimaryVertex> genpvs_;
0090 
0091   // Auxiliary function to get the generated primary vertex
0092   bool isFinalstateParticle(const HepMC::GenParticle *);
0093   bool isCharged(const HepMC::GenParticle *);
0094   void genPrimaryVertices();
0095 };
0096 
0097 #endif