Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:08

0001 //
0002 // Class:           HITRegionalPixelSeedGenerator
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "MagneticField/Engine/interface/MagneticField.h"
0010 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0011 #include "RecoTracker/Record/interface/TrackerMultipleScatteringRecord.h"
0012 #include "RecoTracker/TkMSParametrization/interface/MultipleScatteringParametrisationMaker.h"
0013 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
0014 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
0015 #include "DataFormats/VertexReco/interface/Vertex.h"
0016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0017 #include "DataFormats/TrackReco/interface/Track.h"
0018 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0019 #include "DataFormats/Math/interface/Vector3D.h"
0020 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
0021 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
0022 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0023 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidate.h"
0024 
0025 // Math
0026 #include "Math/GenVector/VectorUtil.h"
0027 #include "Math/GenVector/PxPyPzE4D.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 
0031 #include "DataFormats/Common/interface/Ref.h"
0032 #include "DataFormats/JetReco/interface/Jet.h"
0033 
0034 class HITRegionalPixelSeedGenerator : public TrackingRegionProducer {
0035 public:
0036   explicit HITRegionalPixelSeedGenerator(const edm::ParameterSet& conf, edm::ConsumesCollector&& iC)
0037       : m_regionPSet(conf.getParameter<edm::ParameterSet>("RegionPSet")),
0038         ptmin(m_regionPSet.getParameter<double>("ptMin")),
0039         originradius(m_regionPSet.getParameter<double>("originRadius")),
0040         halflength(m_regionPSet.getParameter<double>("originHalfLength")),
0041         etaCenter_(m_regionPSet.getParameter<double>("etaCenter")),
0042         phiCenter_(m_regionPSet.getParameter<double>("phiCenter")),
0043         deltaTrackEta(m_regionPSet.getParameter<double>("deltaEtaTrackRegion")),
0044         deltaTrackPhi(m_regionPSet.getParameter<double>("deltaPhiTrackRegion")),
0045         deltaL1JetEta(m_regionPSet.getParameter<double>("deltaEtaL1JetRegion")),
0046         deltaL1JetPhi(m_regionPSet.getParameter<double>("deltaPhiL1JetRegion")),
0047         usejets_(m_regionPSet.getParameter<bool>("useL1Jets")),
0048         usetracks_(m_regionPSet.getParameter<bool>("useTracks")),
0049         fixedReg_(m_regionPSet.getParameter<bool>("fixedReg")),
0050         useIsoTracks_(m_regionPSet.getParameter<bool>("useIsoTracks")),
0051         token_bfield(iC.esConsumes()),
0052         token_msmaker(iC.esConsumes()) {
0053     edm::LogVerbatim("HITRegionalPixelSeedGenerator") << "Enter the HITRegionalPixelSeedGenerator";
0054 
0055     if (usetracks_)
0056       token_trks = iC.consumes<reco::TrackCollection>(m_regionPSet.getParameter<edm::InputTag>("trackSrc"));
0057     if (usetracks_ || useIsoTracks_ || fixedReg_ || usejets_)
0058       token_vertex = iC.consumes<reco::VertexCollection>(m_regionPSet.getParameter<edm::InputTag>("vertexSrc"));
0059     if (useIsoTracks_)
0060       token_isoTrack =
0061           iC.consumes<trigger::TriggerFilterObjectWithRefs>(m_regionPSet.getParameter<edm::InputTag>("isoTrackSrc"));
0062     if (usejets_)
0063       token_l1jet =
0064           iC.consumes<l1extra::L1JetParticleCollection>(m_regionPSet.getParameter<edm::InputTag>("l1tjetSrc"));
0065   }
0066 
0067   ~HITRegionalPixelSeedGenerator() override = default;
0068 
0069   std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
0070     std::vector<std::unique_ptr<TrackingRegion> > result;
0071     float originz = 0.;
0072 
0073     double deltaZVertex = halflength;
0074 
0075     auto const& bfield = es.getData(token_bfield);
0076     auto const& msmaker = es.getData(token_msmaker);
0077 
0078     if (usetracks_) {
0079       const edm::Handle<reco::TrackCollection>& tracks = e.getHandle(token_trks);
0080 
0081       const reco::VertexCollection& vertCollection = e.get(token_vertex);
0082       reco::VertexCollection::const_iterator ci = vertCollection.begin();
0083 
0084       if (!vertCollection.empty()) {
0085         originz = ci->z();
0086       } else {
0087         deltaZVertex = 15.;
0088       }
0089 
0090       GlobalVector globalVector(0, 0, 1);
0091       if (tracks->empty())
0092         return result;
0093 
0094       reco::TrackCollection::const_iterator itr = tracks->begin();
0095       for (; itr != tracks->end(); itr++) {
0096         GlobalVector ptrVec((itr)->px(), (itr)->py(), (itr)->pz());
0097         globalVector = ptrVec;
0098 
0099         result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
0100                                                                            GlobalPoint(0, 0, originz),
0101                                                                            ptmin,
0102                                                                            originradius,
0103                                                                            deltaZVertex,
0104                                                                            deltaTrackEta,
0105                                                                            deltaTrackPhi,
0106                                                                            bfield,
0107                                                                            &msmaker));
0108       }
0109     }
0110 
0111     if (useIsoTracks_) {
0112       const edm::Handle<trigger::TriggerFilterObjectWithRefs>& isotracks = e.getHandle(token_isoTrack);
0113 
0114       std::vector<edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
0115 
0116       isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
0117 
0118       const reco::VertexCollection& vertCollection = e.get(token_vertex);
0119       reco::VertexCollection::const_iterator ci = vertCollection.begin();
0120 
0121       if (!vertCollection.empty()) {
0122         originz = ci->z();
0123       } else {
0124         deltaZVertex = 15.;
0125       }
0126 
0127       GlobalVector globalVector(0, 0, 1);
0128       if (isoPixTrackRefs.empty())
0129         return result;
0130 
0131       for (uint32_t p = 0; p < isoPixTrackRefs.size(); p++) {
0132         GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),
0133                             (isoPixTrackRefs[p]->track())->py(),
0134                             (isoPixTrackRefs[p]->track())->pz());
0135         globalVector = ptrVec;
0136 
0137         result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
0138                                                                            GlobalPoint(0, 0, originz),
0139                                                                            ptmin,
0140                                                                            originradius,
0141                                                                            deltaZVertex,
0142                                                                            deltaTrackEta,
0143                                                                            deltaTrackPhi,
0144                                                                            bfield,
0145                                                                            &msmaker));
0146       }
0147     }
0148 
0149     if (usejets_) {
0150       const edm::Handle<l1extra::L1JetParticleCollection>& jets = e.getHandle(token_l1jet);
0151       const reco::VertexCollection& vertCollection = e.get(token_vertex);
0152       reco::VertexCollection::const_iterator ci = vertCollection.begin();
0153       if (!vertCollection.empty()) {
0154         originz = ci->z();
0155       } else {
0156         deltaZVertex = 15.;
0157       }
0158 
0159       GlobalVector globalVector(0, 0, 1);
0160       if (jets->empty())
0161         return result;
0162 
0163       for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++) {
0164         GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
0165         GlobalPoint vertex(0, 0, originz);
0166 
0167         result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
0168             jetVector, vertex, ptmin, originradius, deltaZVertex, deltaL1JetEta, deltaL1JetPhi, bfield, &msmaker));
0169       }
0170     }
0171     if (fixedReg_) {
0172       GlobalVector fixedVector(cos(phiCenter_) * sin(2 * atan(exp(-etaCenter_))),
0173                                sin(phiCenter_) * sin(2 * atan(exp(-etaCenter_))),
0174                                cos(2 * atan(exp(-etaCenter_))));
0175       GlobalPoint vertex(0, 0, originz);
0176 
0177       const reco::VertexCollection& vertCollection = e.get(token_vertex);
0178       if (!vertCollection.empty()) {
0179         //      reco::VertexCollection::const_iterator ci = vertCollection.begin();
0180         //      originz = ci->z();
0181       } else {
0182         deltaZVertex = 15.;
0183       }
0184 
0185       result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
0186           fixedVector, vertex, ptmin, originradius, deltaZVertex, deltaL1JetEta, deltaL1JetPhi, bfield, &msmaker));
0187     }
0188 
0189     return result;
0190   }
0191 
0192 private:
0193   const edm::ParameterSet m_regionPSet;
0194   const float ptmin;
0195   const float originradius;
0196   const float halflength;
0197   const double etaCenter_;
0198   const double phiCenter_;
0199   const float deltaTrackEta;
0200   const float deltaTrackPhi;
0201   const float deltaL1JetEta;
0202   const float deltaL1JetPhi;
0203   const bool usejets_;
0204   const bool usetracks_;
0205   const bool fixedReg_;
0206   const bool useIsoTracks_;
0207   edm::EDGetTokenT<reco::TrackCollection> token_trks;
0208   edm::EDGetTokenT<reco::VertexCollection> token_vertex;
0209   edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> token_isoTrack;
0210   edm::EDGetTokenT<l1extra::L1JetParticleCollection> token_l1jet;
0211   const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> token_bfield;
0212   const edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> token_msmaker;
0213 };
0214 
0215 #include "FWCore/PluginManager/interface/ModuleDef.h"
0216 #include "FWCore/Framework/interface/MakerMacros.h"
0217 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducerFactory.h"
0218 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
0219 #
0220 DEFINE_EDM_PLUGIN(TrackingRegionProducerFactory, HITRegionalPixelSeedGenerator, "HITRegionalPixelSeedGenerator");