Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:21

0001 #ifndef RecoHI_HiTracking_HITrackingRegionForPrimaryVtxProducer_H
0002 #define RecoHI_HiTracking_HITrackingRegionForPrimaryVtxProducer _H
0003 
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006 
0007 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
0008 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
0009 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/ConsumesCollector.h"
0013 
0014 #include "MagneticField/Engine/interface/MagneticField.h"
0015 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0016 #include "RecoTracker/Record/interface/TrackerMultipleScatteringRecord.h"
0017 #include "RecoTracker/TkMSParametrization/interface/MultipleScatteringParametrisationMaker.h"
0018 
0019 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0020 #include "FWCore/Utilities/interface/InputTag.h"
0021 
0022 #include "DataFormats/VertexReco/interface/Vertex.h"
0023 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0024 
0025 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0026 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0027 #include "DataFormats/Common/interface/DetSetAlgorithm.h"
0028 
0029 #include "DataFormats/Common/interface/DetSetVector.h"
0030 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0031 
0032 #include "TMath.h"
0033 
0034 class HITrackingRegionForPrimaryVtxProducer : public TrackingRegionProducer {
0035 public:
0036   HITrackingRegionForPrimaryVtxProducer(const edm::ParameterSet& cfg, edm::ConsumesCollector&& iC)
0037       : theTtopoToken(iC.esConsumes()), theFieldToken(iC.esConsumes()) {
0038     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
0039     thePtMin = regionPSet.getParameter<double>("ptMin");
0040     theOriginRadius = regionPSet.getParameter<double>("originRadius");
0041     theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
0042     theBeamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot");
0043     theBeamSpotToken = iC.consumes<reco::BeamSpot>(theBeamSpotTag);
0044     thePrecise = regionPSet.getParameter<bool>("precise");
0045     theSiPixelRecHits = regionPSet.getParameter<edm::InputTag>("siPixelRecHits");
0046     theSiPixelRecHitsToken = iC.consumes<SiPixelRecHitCollection>(theSiPixelRecHits);
0047     doVariablePtMin = regionPSet.getParameter<bool>("doVariablePtMin");
0048     double xDir = regionPSet.getParameter<double>("directionXCoord");
0049     double yDir = regionPSet.getParameter<double>("directionYCoord");
0050     double zDir = regionPSet.getParameter<double>("directionZCoord");
0051     theDirection = GlobalVector(xDir, yDir, zDir);
0052 
0053     // for using vertex instead of beamspot
0054     theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
0055     theFixedError = regionPSet.getParameter<double>("fixedError");
0056     theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
0057     theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
0058     vertexCollName = regionPSet.getParameter<edm::InputTag>("VertexCollection");
0059     vertexCollToken = iC.consumes<reco::VertexCollection>(vertexCollName);
0060 
0061     if (thePrecise) {
0062       theMSMakerToken = iC.esConsumes();
0063     }
0064   }
0065 
0066   ~HITrackingRegionForPrimaryVtxProducer() override = default;
0067 
0068   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0069     edm::ParameterSetDescription desc;
0070 
0071     desc.add<double>("ptMin", 0.7);
0072     desc.add<bool>("doVariablePtMin", true);
0073     desc.add<double>("originRadius", 0.2);
0074     desc.add<double>("nSigmaZ", 3.0);
0075     desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
0076     desc.add<bool>("precise", true);
0077     desc.add<bool>("useMultipleScattering", false);
0078     desc.add<bool>("useFakeVertices", false);
0079     desc.add<edm::InputTag>("siPixelRecHits", edm::InputTag("siPixelRecHits"));
0080     desc.add<double>("directionXCoord", 1.0);
0081     desc.add<double>("directionYCoord", 1.0);
0082     desc.add<double>("directionZCoord", 0.0);
0083     desc.add<bool>("useFoundVertices", true);
0084     desc.add<edm::InputTag>("VertexCollection", edm::InputTag("hiPixelClusterVertex"));
0085     desc.add<bool>("useFixedError", true);
0086     desc.add<double>("fixedError", 3.0);
0087     desc.add<double>("sigmaZVertex", 3.0);
0088 
0089     // Only for backwards-compatibility
0090     edm::ParameterSetDescription descRegion;
0091     descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
0092 
0093     descriptions.add("hiTrackingRegionFromClusterVtx", descRegion);
0094   }
0095 
0096   int estimateMultiplicity(const edm::Event& ev, const edm::EventSetup& es) const {
0097     //rechits
0098     edm::Handle<SiPixelRecHitCollection> recHitColl;
0099     ev.getByToken(theSiPixelRecHitsToken, recHitColl);
0100 
0101     auto const& ttopo = es.getData(theTtopoToken);
0102 
0103     std::vector<const TrackingRecHit*> theChosenHits;
0104     edmNew::copyDetSetRange(*recHitColl, theChosenHits, ttopo.pxbDetIdLayerComparator(1));
0105     return theChosenHits.size();
0106   }
0107 
0108   std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
0109                                                         const edm::EventSetup& es) const override {
0110     int estMult = estimateMultiplicity(ev, es);
0111 
0112     // from MC relating first layer pixel hits to "findable" sim tracks with pt>1 GeV
0113     float cc = -38.6447;
0114     float bb = 0.0581765;
0115     float aa = 1.34306e-06;
0116 
0117     float estTracks = aa * estMult * estMult + bb * estMult + cc;
0118 
0119     LogTrace("heavyIonHLTVertexing") << "[HIVertexing]";
0120     LogTrace("heavyIonHLTVertexing") << " [HIVertexing: hits in the 1. layer:" << estMult << "]";
0121     LogTrace("heavyIonHLTVertexing") << " [HIVertexing: estimated number of tracks:" << estTracks << "]";
0122 
0123     float regTracking = 60.;  //if we have more tracks -> regional tracking
0124     float etaB = 10.;
0125     float phiB = TMath::Pi() / 2.;
0126 
0127     float decEta = estTracks / 90.;
0128     etaB = 2.5 / decEta;
0129 
0130     if (estTracks > regTracking) {
0131       LogTrace("heavyIonHLTVertexing") << " [HIVertexing: Regional Tracking]";
0132       LogTrace("heavyIonHLTVertexing") << "  [Regional Tracking: eta range: -" << etaB << ", " << etaB << "]";
0133       LogTrace("heavyIonHLTVertexing") << "  [Regional Tracking: phi range: -" << phiB << ", " << phiB << "]";
0134       LogTrace("heavyIonHLTVertexing") << "  [Regional Tracking: factor of decrease: " << decEta * 2.
0135                                        << "]";  // 2:from phi
0136     }
0137 
0138     float minpt = thePtMin;
0139     float varPtCutoff = 1500;  //cutoff
0140     if (doVariablePtMin && estMult < varPtCutoff) {
0141       minpt = 0.075;
0142       if (estMult > 0)
0143         minpt += estMult * (thePtMin - 0.075) / varPtCutoff;  // lower ptMin linearly with pixel hit multiplicity
0144     }
0145 
0146     // tracking region selection
0147     std::vector<std::unique_ptr<TrackingRegion> > result;
0148     double halflength;
0149     GlobalPoint origin;
0150     edm::Handle<reco::BeamSpot> bsHandle;
0151     ev.getByToken(theBeamSpotToken, bsHandle);
0152     if (bsHandle.isValid()) {
0153       const reco::BeamSpot& bs = *bsHandle;
0154       origin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
0155       halflength = theNSigmaZ * bs.sigmaZ();
0156 
0157       if (theUseFoundVertices) {
0158         edm::Handle<reco::VertexCollection> vertexCollection;
0159         ev.getByToken(vertexCollToken, vertexCollection);
0160 
0161         for (reco::VertexCollection::const_iterator iV = vertexCollection->begin(); iV != vertexCollection->end();
0162              iV++) {
0163           if (iV->isFake() || !iV->isValid())
0164             continue;
0165           origin = GlobalPoint(bs.x0(), bs.y0(), iV->z());
0166           halflength = (theUseFixedError ? theFixedError : (iV->zError()) * theSigmaZVertex);
0167         }
0168       }
0169 
0170       if (estTracks > regTracking) {  // regional tracking
0171         const auto& field = es.getData(theFieldToken);
0172         const MultipleScatteringParametrisationMaker* msmaker = nullptr;
0173         if (thePrecise) {
0174           msmaker = &es.getData(theMSMakerToken);
0175         }
0176 
0177         result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
0178             theDirection, origin, thePtMin, theOriginRadius, halflength, etaB, phiB, field, msmaker, thePrecise));
0179       } else {  // global tracking
0180         LogTrace("heavyIonHLTVertexing") << " [HIVertexing: Global Tracking]";
0181         result.push_back(
0182             std::make_unique<GlobalTrackingRegion>(minpt, origin, theOriginRadius, halflength, thePrecise));
0183       }
0184     }
0185     return result;
0186   }
0187 
0188 private:
0189   double thePtMin;
0190   double theOriginRadius;
0191   double theNSigmaZ;
0192   edm::InputTag theBeamSpotTag;
0193   edm::EDGetTokenT<reco::BeamSpot> theBeamSpotToken;
0194   bool thePrecise;
0195   GlobalVector theDirection;
0196   edm::InputTag theSiPixelRecHits;
0197   edm::EDGetTokenT<SiPixelRecHitCollection> theSiPixelRecHitsToken;
0198   bool doVariablePtMin;
0199 
0200   double theSigmaZVertex;
0201   double theFixedError;
0202   bool theUseFoundVertices;
0203   bool theUseFixedError;
0204   edm::InputTag vertexCollName;
0205   edm::EDGetTokenT<reco::VertexCollection> vertexCollToken;
0206   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> theTtopoToken;
0207   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> theFieldToken;
0208   edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> theMSMakerToken;
0209 };
0210 
0211 #endif