Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:05

0001 
0002 #ifndef TrackClassifier_h
0003 #define TrackClassifier_h
0004 
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/Utilities/interface/InputTag.h"
0008 
0009 #include "MagneticField/Engine/interface/MagneticField.h"
0010 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0011 
0012 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0013 
0014 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
0015 
0016 #include "SimTracker/TrackHistory/interface/CMSProcessTypes.h"
0017 #include "SimTracker/TrackHistory/interface/TrackCategories.h"
0018 #include "SimTracker/TrackHistory/interface/TrackHistory.h"
0019 #include "SimTracker/TrackHistory/interface/TrackQuality.h"
0020 
0021 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
0022 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0023 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
0024 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0025 
0026 class TrackTopology;
0027 
0028 //! Get track history and classify it in function of their .
0029 class TrackClassifier : public TrackCategories {
0030 public:
0031   //! Type to the associate category
0032   typedef TrackCategories Categories;
0033 
0034   //! Constructor by ParameterSet
0035   TrackClassifier(edm::ParameterSet const &, edm::ConsumesCollector &&);
0036 
0037   //! Pre-process event information (for accessing reconstraction information)
0038   void newEvent(edm::Event const &, edm::EventSetup const &);
0039 
0040   //! Classify the RecoTrack in categories.
0041   TrackClassifier const &evaluate(reco::TrackBaseRef const &);
0042 
0043   //! Classify the TrackingParticle in categories.
0044   TrackClassifier const &evaluate(TrackingParticleRef const &);
0045 
0046   //! Classify the RecoTrack in categories.
0047   TrackClassifier const &evaluate(reco::TrackRef const &track) { return evaluate(reco::TrackBaseRef(track)); }
0048 
0049   //! Returns a reference to the track history used in the classification.
0050   TrackHistory const &history() const { return tracer_; }
0051 
0052   //! Returns a reference to the track quality used in the classification.
0053   TrackQuality const &quality() const { return quality_; }
0054 
0055 private:
0056   const edm::InputTag hepMCLabel_;
0057   const edm::InputTag beamSpotLabel_;
0058 
0059   double badPull_;
0060   double longLivedDecayLength_;
0061   double vertexClusteringSqDistance_;
0062   unsigned int numberOfInnerLayers_;
0063   unsigned int minTrackerSimHits_;
0064 
0065   TrackHistory tracer_;
0066 
0067   TrackQuality quality_;
0068 
0069   const G4toCMSLegacyProcTypeMap g4toCMSProcMap_;
0070 
0071   edm::ESHandle<MagneticField> magneticField_;
0072   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magneticFieldToken_;
0073 
0074   edm::Handle<edm::HepMCProduct> mcInformation_;
0075 
0076   edm::ESHandle<ParticleDataTable> particleDataTable_;
0077   edm::ESGetToken<ParticleDataTable, PDTRecord> particleDataTableToken_;
0078 
0079   edm::ESHandle<TransientTrackBuilder> transientTrackBuilder_;
0080   edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> transientTrackBuilderToken_;
0081 
0082   edm::Handle<reco::BeamSpot> beamSpot_;
0083 
0084   const TrackerTopology *tTopo_;
0085   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoHandToken_;
0086 
0087   //! Classify all the tracks by their association and reconstruction
0088   //! information
0089   void reconstructionInformation(reco::TrackBaseRef const &);
0090 
0091   //! Get all the information related to the simulation details
0092   void simulationInformation();
0093 
0094   //! Classify all the tracks by their reconstruction quality
0095   void qualityInformation(reco::TrackBaseRef const &);
0096 
0097   //! Get hadron flavor of the initial hadron
0098   void hadronFlavor();
0099 
0100   //! Get all the information related to decay process
0101   void processesAtGenerator();
0102 
0103   //! Get information about conversion and other interactions
0104   void processesAtSimulation();
0105 
0106   //! Get geometrical information about the vertices
0107   void vertexInformation();
0108 
0109   //! Auxiliary class holding simulated primary vertices
0110   struct GeneratedPrimaryVertex {
0111     GeneratedPrimaryVertex(double x1, double y1, double z1) : x(x1), y(y1), z(z1), ptsq(0), nGenTrk(0) {}
0112 
0113     bool operator<(GeneratedPrimaryVertex const &reference) const { return ptsq < reference.ptsq; }
0114 
0115     double x, y, z;
0116     double ptsq;
0117     int nGenTrk;
0118 
0119     HepMC::FourVector ptot;
0120 
0121     std::vector<int> finalstateParticles;
0122     std::vector<int> simTrackIndex;
0123     std::vector<int> genVertex;
0124   };
0125 
0126   std::vector<GeneratedPrimaryVertex> genpvs_;
0127 
0128   // Auxiliary function to get the generated primary vertex
0129   bool isFinalstateParticle(const HepMC::GenParticle *);
0130   bool isCharged(const HepMC::GenParticle *);
0131   void genPrimaryVertices();
0132 };
0133 
0134 #endif