Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:56

0001 #include "RecoMuon/GlobalTrackingTools/interface/MuonTrackingRegionByPtBuilder.h"
0002 
0003 //-------------------------------

0004 // Collaborating Class Headers --

0005 //-------------------------------

0006 
0007 #include "DataFormats/Common/interface/Handle.h"
0008 
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "MagneticField/Engine/interface/MagneticField.h"
0015 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0016 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
0017 #include "RecoTracker/MeasurementDet/interface/MeasurementTrackerEvent.h"
0018 #include "RecoTracker/Record/interface/TrackerMultipleScatteringRecord.h"
0019 #include "RecoTracker/TkMSParametrization/interface/MultipleScatteringParametrisationMaker.h"
0020 
0021 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0022 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
0023 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0024 
0025 //

0026 // constructor

0027 //

0028 void MuonTrackingRegionByPtBuilder::build(const edm::ParameterSet& par, edm::ConsumesCollector& iC) {
0029   // Adjust errors on dz

0030   theNsigmaDz = par.getParameter<double>("Rescale_Dz");
0031 
0032   // Flag to switch to use Vertices instead of BeamSpot

0033   useVertex = par.getParameter<bool>("UseVertex");
0034 
0035   // Flag to use fixed limits for Eta, Phi, Z, pT

0036   useFixedZ = par.getParameter<bool>("Z_fixed");
0037   useFixedPt = par.getParameter<bool>("Pt_fixed");
0038 
0039   // Minimum value for pT

0040   thePtMin = par.getParameter<double>("Pt_min");
0041 
0042   // The static region size along the Z direction

0043   theHalfZ = par.getParameter<double>("DeltaZ");
0044 
0045   // The transverse distance of the region from the BS/PV

0046   theDeltaR = par.getParameter<double>("DeltaR");
0047 
0048   // Parameterized ROIs depending on the pt

0049   ptRanges_ = par.getParameter<std::vector<double>>("ptRanges");
0050   if (ptRanges_.size() < 2) {
0051     edm::LogError("MuonTrackingRegionByPtBuilder") << "Size of ptRanges does not be less than 2." << std::endl;
0052   }
0053 
0054   deltaEtas_ = par.getParameter<std::vector<double>>("deltaEtas");
0055   if (deltaEtas_.size() != ptRanges_.size() - 1) {
0056     edm::LogError("MuonTrackingRegionByPtBuilder")
0057         << "Size of deltaEtas does not match number of pt bins." << std::endl;
0058   }
0059 
0060   deltaPhis_ = par.getParameter<std::vector<double>>("deltaPhis");
0061   if (deltaPhis_.size() != ptRanges_.size() - 1) {
0062     edm::LogError("MuonTrackingRegionByPtBuilder")
0063         << "Size of deltaPhis does not match number of pt bins." << std::endl;
0064   }
0065 
0066   // Maximum number of regions to build when looping over Muons

0067   theMaxRegions = par.getParameter<int>("maxRegions");
0068 
0069   // Flag to use precise??

0070   thePrecise = par.getParameter<bool>("precise");
0071 
0072   // perigee reference point

0073   theOnDemand = RectangularEtaPhiTrackingRegion::intToUseMeasurementTracker(par.getParameter<int>("OnDemand"));
0074   if (theOnDemand != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever) {
0075     theMeasurementTrackerToken =
0076         iC.consumes<MeasurementTrackerEvent>(par.getParameter<edm::InputTag>("MeasurementTrackerName"));
0077   }
0078 
0079   // Vertex collection and Beam Spot

0080   beamSpotToken = iC.consumes<reco::BeamSpot>(par.getParameter<edm::InputTag>("beamSpot"));
0081   vertexCollectionToken = iC.consumes<reco::VertexCollection>(par.getParameter<edm::InputTag>("vertexCollection"));
0082 
0083   // Input muon collection

0084   inputCollectionToken = iC.consumes<reco::TrackCollection>(par.getParameter<edm::InputTag>("input"));
0085 
0086   bfieldToken = iC.esConsumes();
0087   if (thePrecise) {
0088     msmakerToken = iC.esConsumes();
0089   }
0090 }
0091 
0092 //

0093 // Member function to be compatible with TrackingRegionProducerFactory: create many ROI for many tracks

0094 //

0095 std::vector<std::unique_ptr<TrackingRegion>> MuonTrackingRegionByPtBuilder::regions(const edm::Event& ev,
0096                                                                                     const edm::EventSetup& es) const {
0097   std::vector<std::unique_ptr<TrackingRegion>> result;
0098 
0099   edm::Handle<reco::TrackCollection> tracks;
0100   ev.getByToken(inputCollectionToken, tracks);
0101 
0102   int nRegions = 0;
0103   for (auto it = tracks->cbegin(), ed = tracks->cend(); it != ed && nRegions < theMaxRegions; ++it) {
0104     result.push_back(region(*it, ev, es));
0105     nRegions++;
0106   }
0107 
0108   return result;
0109 }
0110 
0111 //

0112 // Call region on Track from TrackRef

0113 //

0114 std::unique_ptr<RectangularEtaPhiTrackingRegion> MuonTrackingRegionByPtBuilder::region(
0115     const reco::TrackRef& track) const {
0116   return region(*track);
0117 }
0118 
0119 //

0120 // ToDo: Not sure if this is needed?

0121 //

0122 void MuonTrackingRegionByPtBuilder::setEvent(const edm::Event& event, const edm::EventSetup& es) {
0123   theEvent = &event;
0124   theEventSetup = &es;
0125 }
0126 
0127 //

0128 //  Main member function called to create the ROI

0129 //

0130 std::unique_ptr<RectangularEtaPhiTrackingRegion> MuonTrackingRegionByPtBuilder::region(
0131     const reco::Track& staTrack, const edm::Event& ev, const edm::EventSetup& es) const {
0132   // get track momentum/direction at vertex

0133   const math::XYZVector& mom = staTrack.momentum();
0134   GlobalVector dirVector(mom.x(), mom.y(), mom.z());
0135   double pt = staTrack.pt();
0136 
0137   // Fix for StandAlone tracks with low momentum

0138   const math::XYZVector& innerMomentum = staTrack.innerMomentum();
0139   GlobalVector forSmallMomentum(innerMomentum.x(), innerMomentum.y(), innerMomentum.z());
0140   if (staTrack.p() <= 1.5) {
0141     pt = std::abs(forSmallMomentum.perp());
0142   }
0143 
0144   // initial vertex position - in the following it is replaced with beamspot/vertexing

0145   GlobalPoint vertexPos(0.0, 0.0, 0.0);
0146   // standard 15.9, if useVertex than use error from  vertex

0147   double deltaZ = theHalfZ;
0148 
0149   // retrieve beam spot information

0150   edm::Handle<reco::BeamSpot> bs;
0151   bool bsHandleFlag = ev.getByToken(beamSpotToken, bs);
0152 
0153   // check the validity, otherwise vertexing

0154   if (bsHandleFlag && bs.isValid() && !useVertex) {
0155     vertexPos = GlobalPoint(bs->x0(), bs->y0(), bs->z0());
0156     deltaZ = useFixedZ ? theHalfZ : bs->sigmaZ() * theNsigmaDz;
0157   } else {
0158     // get originZPos from list of reconstructed vertices (first or all)

0159     edm::Handle<reco::VertexCollection> vertexCollection;
0160     bool vtxHandleFlag = ev.getByToken(vertexCollectionToken, vertexCollection);
0161     // check if there exists at least one reconstructed vertex

0162     if (vtxHandleFlag && !vertexCollection->empty()) {
0163       // use the first vertex in the collection and assume it is the primary event vertex

0164       reco::VertexCollection::const_iterator vtx = vertexCollection->begin();
0165       if (!vtx->isFake() && vtx->isValid()) {
0166         vertexPos = GlobalPoint(vtx->x(), vtx->y(), vtx->z());
0167         deltaZ = useFixedZ ? theHalfZ : vtx->zError() * theNsigmaDz;
0168       }
0169     }
0170   }
0171 
0172   // set delta eta and delta phi depending on the pt

0173   auto region_dEta = deltaEtas_.at(0);
0174   auto region_dPhi = deltaPhis_.at(0);
0175   if (pt < ptRanges_.back()) {
0176     auto lowEdge = std::upper_bound(ptRanges_.begin(), ptRanges_.end(), pt);
0177     region_dEta = deltaEtas_.at(lowEdge - ptRanges_.begin() - 1);
0178     region_dPhi = deltaPhis_.at(lowEdge - ptRanges_.begin() - 1);
0179   }
0180 
0181   float deltaR = theDeltaR;
0182   double minPt = useFixedPt ? thePtMin : std::max(thePtMin, pt * 0.6);
0183 
0184   const MeasurementTrackerEvent* measurementTracker = nullptr;
0185   if (!theMeasurementTrackerToken.isUninitialized()) {
0186     edm::Handle<MeasurementTrackerEvent> hmte;
0187     ev.getByToken(theMeasurementTrackerToken, hmte);
0188     measurementTracker = hmte.product();
0189   }
0190 
0191   const auto& bfield = es.getData(bfieldToken);
0192   const MultipleScatteringParametrisationMaker* msmaker = nullptr;
0193   if (thePrecise) {
0194     msmaker = &es.getData(msmakerToken);
0195   }
0196 
0197   auto region = std::make_unique<RectangularEtaPhiTrackingRegion>(dirVector,
0198                                                                   vertexPos,
0199                                                                   minPt,
0200                                                                   deltaR,
0201                                                                   deltaZ,
0202                                                                   region_dEta,
0203                                                                   region_dPhi,
0204                                                                   bfield,
0205                                                                   msmaker,
0206                                                                   thePrecise,
0207                                                                   theOnDemand,
0208                                                                   measurementTracker);
0209 
0210   LogDebug("MuonTrackingRegionByPtBuilder")
0211       << "the region parameters are:\n"
0212       << "\n dirVector: " << dirVector << "\n vertexPos: " << vertexPos << "\n minPt: " << minPt
0213       << "\n deltaR:" << deltaR << "\n deltaZ:" << deltaZ << "\n region_dEta:" << region_dEta
0214       << "\n region_dPhi:" << region_dPhi << "\n on demand parameter: " << static_cast<int>(theOnDemand);
0215 
0216   return region;
0217 }
0218 
0219 void MuonTrackingRegionByPtBuilder::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0220   edm::ParameterSetDescription desc;
0221   desc.add<double>("DeltaR", 0.2);
0222   desc.add<edm::InputTag>("beamSpot", edm::InputTag(""));
0223   desc.add<int>("OnDemand", -1);
0224   desc.add<edm::InputTag>("vertexCollection", edm::InputTag(""));
0225   desc.add<edm::InputTag>("MeasurementTrackerName", edm::InputTag(""));
0226   desc.add<bool>("UseVertex", false);
0227   desc.add<double>("Rescale_Dz", 3.0);
0228   desc.add<bool>("Pt_fixed", false);
0229   desc.add<bool>("Z_fixed", true);
0230   desc.add<double>("Pt_min", 1.5);
0231   desc.add<double>("DeltaZ", 15.9);
0232   desc.add<std::vector<double>>("ptRanges", {0., 1.e9});
0233   desc.add<std::vector<double>>("deltaEtas", {0.2});
0234   desc.add<std::vector<double>>("deltaPhis", {0.15});
0235   desc.add<int>("maxRegions", 1);
0236   desc.add<bool>("precise", true);
0237   desc.add<edm::InputTag>("input", edm::InputTag(""));
0238   descriptions.add("MuonTrackingRegionByPtBuilder", desc);
0239 
0240   descriptions.setComment(
0241       "Build a TrackingRegion around a standalone muon. Options to define region around beamspot or primary vertex and "
0242       "dynamic regions are included.");
0243 }