Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef HLTrigger_btau_L3MumuTrackingRegion_H
0002 #define HLTrigger_btau_L3MumuTrackingRegion_H
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 
0007 #include "MagneticField/Engine/interface/MagneticField.h"
0008 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0009 #include "RecoTracker/Record/interface/TrackerMultipleScatteringRecord.h"
0010 #include "RecoTracker/TkMSParametrization/interface/MultipleScatteringParametrisationMaker.h"
0011 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
0012 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
0013 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
0014 #include "RecoTracker/MeasurementDet/interface/MeasurementTrackerEvent.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "DataFormats/VertexReco/interface/Vertex.h"
0018 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0019 #include "DataFormats/TrackReco/interface/Track.h"
0020 
0021 class L3MumuTrackingRegion : public TrackingRegionProducer {
0022 public:
0023   L3MumuTrackingRegion(const edm::ParameterSet& cfg, edm::ConsumesCollector&& iC)
0024       : theFieldToken(iC.esConsumes()), theMSMakerToken(iC.esConsumes()) {
0025     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
0026 
0027     theVertexTag = regionPSet.getParameter<edm::InputTag>("vertexSrc");
0028     theVertex = (theVertexTag.label().length() > 1);
0029     theInputTrkTag = regionPSet.getParameter<edm::InputTag>("TrkSrc");
0030     useVtxTks = regionPSet.getParameter<bool>("UseVtxTks");
0031 
0032     if (theVertex)
0033       theVertexToken = iC.consumes<reco::VertexCollection>(theVertexTag);
0034     if (!(theVertex && useVtxTks))
0035       theInputTrkToken = iC.consumes<reco::TrackCollection>(theInputTrkTag);
0036 
0037     thePtMin = regionPSet.getParameter<double>("ptMin");
0038     theOriginRadius = regionPSet.getParameter<double>("originRadius");
0039     theOriginHalfLength = regionPSet.getParameter<double>("originHalfLength");
0040     theOriginZPos = regionPSet.getParameter<double>("vertexZDefault");
0041 
0042     theDeltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
0043     theDeltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
0044     if (regionPSet.exists("searchOpt")) {
0045       m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
0046     } else {
0047       m_searchOpt = false;
0048     }
0049     m_howToUseMeasurementTracker = RectangularEtaPhiTrackingRegion::stringToUseMeasurementTracker(
0050         regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
0051     if (m_howToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever) {
0052       theMeasurementTrackerToken =
0053           iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<edm::InputTag>("measurementTrackerName"));
0054     }
0055   }
0056 
0057   ~L3MumuTrackingRegion() override = default;
0058 
0059   std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
0060                                                         const edm::EventSetup& es) const override {
0061     std::vector<std::unique_ptr<TrackingRegion> > result;
0062 
0063     const MeasurementTrackerEvent* measurementTracker = nullptr;
0064     if (!theMeasurementTrackerToken.isUninitialized()) {
0065       edm::Handle<MeasurementTrackerEvent> hmte;
0066       ev.getByToken(theMeasurementTrackerToken, hmte);
0067       measurementTracker = hmte.product();
0068     }
0069     const auto& field = es.getData(theFieldToken);
0070     const auto& msmaker = es.getData(theMSMakerToken);
0071 
0072     // optional constraint for vertex
0073     // get highest Pt pixel vertex (if existing)
0074     double deltaZVertex = theOriginHalfLength;
0075     double originz = theOriginZPos;
0076     if (theVertex) {
0077       edm::Handle<reco::VertexCollection> vertices;
0078       ev.getByToken(theVertexToken, vertices);
0079       const reco::VertexCollection vertCollection = *(vertices.product());
0080       reco::VertexCollection::const_iterator ci = vertCollection.begin();
0081       if (!vertCollection.empty()) {
0082         originz = ci->z();
0083       } else {
0084         originz = theOriginZPos;
0085         deltaZVertex = 15.;
0086       }
0087       if (useVtxTks) {
0088         for (ci = vertCollection.begin(); ci != vertCollection.end(); ci++)
0089           for (reco::Vertex::trackRef_iterator trackIt = ci->tracks_begin(); trackIt != ci->tracks_end(); trackIt++) {
0090             reco::TrackRef iTrk = (*trackIt).castTo<reco::TrackRef>();
0091             GlobalVector dirVector((iTrk)->px(), (iTrk)->py(), (iTrk)->pz());
0092             result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(dirVector,
0093                                                                                GlobalPoint(0, 0, float(ci->z())),
0094                                                                                thePtMin,
0095                                                                                theOriginRadius,
0096                                                                                deltaZVertex,
0097                                                                                theDeltaEta,
0098                                                                                theDeltaPhi,
0099                                                                                field,
0100                                                                                &msmaker,
0101                                                                                true,
0102                                                                                m_howToUseMeasurementTracker,
0103                                                                                measurementTracker,
0104                                                                                m_searchOpt));
0105           }
0106         return result;
0107       }
0108     }
0109 
0110     edm::Handle<reco::TrackCollection> trks;
0111     if (!theInputTrkToken.isUninitialized())
0112       ev.getByToken(theInputTrkToken, trks);
0113     for (reco::TrackCollection::const_iterator iTrk = trks->begin(); iTrk != trks->end(); iTrk++) {
0114       GlobalVector dirVector((iTrk)->px(), (iTrk)->py(), (iTrk)->pz());
0115       result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(dirVector,
0116                                                                          GlobalPoint(0, 0, float(originz)),
0117                                                                          thePtMin,
0118                                                                          theOriginRadius,
0119                                                                          deltaZVertex,
0120                                                                          theDeltaEta,
0121                                                                          theDeltaPhi,
0122                                                                          field,
0123                                                                          &msmaker,
0124                                                                          true,
0125                                                                          m_howToUseMeasurementTracker,
0126                                                                          measurementTracker,
0127                                                                          m_searchOpt));
0128     }
0129     return result;
0130   }
0131 
0132 private:
0133   edm::InputTag theVertexTag;
0134   bool theVertex;
0135   edm::EDGetTokenT<reco::VertexCollection> theVertexToken;
0136   edm::InputTag theInputTrkTag;
0137   edm::EDGetTokenT<reco::TrackCollection> theInputTrkToken;
0138 
0139   bool useVtxTks;
0140 
0141   double thePtMin;
0142   double theOriginRadius;
0143   double theOriginHalfLength;
0144   double theOriginZPos;
0145 
0146   double theDeltaEta;
0147   double theDeltaPhi;
0148   edm::EDGetTokenT<MeasurementTrackerEvent> theMeasurementTrackerToken;
0149   RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker;
0150   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> theFieldToken;
0151   edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> theMSMakerToken;
0152   bool m_searchOpt;
0153 };
0154 
0155 #endif