Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:57

0001 #ifndef RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
0002 #define RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
0003 
0004 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0005 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0006 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
0007 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "DataFormats/VertexReco/interface/Vertex.h"
0013 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0014 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "FWCore/Framework/interface/ConsumesCollector.h"
0017 
0018 class GlobalTrackingRegionWithVerticesProducer : public TrackingRegionProducer {
0019 public:
0020   GlobalTrackingRegionWithVerticesProducer(const edm::ParameterSet& cfg, edm::ConsumesCollector&& iC) {
0021     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
0022 
0023     thePtMin = regionPSet.getParameter<double>("ptMin");
0024     theOriginRadius = regionPSet.getParameter<double>("originRadius");
0025     theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
0026     token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
0027     thePrecise = regionPSet.getParameter<bool>("precise");
0028     theUseMS = regionPSet.getParameter<bool>("useMultipleScattering");
0029 
0030     theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
0031     theFixedError = regionPSet.getParameter<double>("fixedError");
0032 
0033     theMaxNVertices = regionPSet.getParameter<int>("maxNVertices");
0034 
0035     theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
0036     theUseFakeVertices = regionPSet.getParameter<bool>("useFakeVertices");
0037     theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
0038     token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("VertexCollection"));
0039 
0040     //information for Heavy ion region scaling
0041     theOriginRScaling = regionPSet.getParameter<bool>("originRScaling4BigEvts");
0042     thePtMinScaling = regionPSet.getParameter<bool>("ptMinScaling4BigEvts");
0043     theHalfLengthScaling = regionPSet.getParameter<bool>("halfLengthScaling4BigEvts");
0044     theAllowEmpty = regionPSet.getParameter<bool>("allowEmpty");
0045     theMinOriginR = regionPSet.getParameter<double>("minOriginR");
0046     theMaxPtMin = regionPSet.getParameter<double>("maxPtMin");
0047     theMinHalfLength = regionPSet.getParameter<double>("minHalfLength");
0048     theScalingStart = regionPSet.getParameter<double>("scalingStartNPix");
0049     theScalingEnd = regionPSet.getParameter<double>("scalingEndNPix");
0050     edm::InputTag pixelClustersForScaling = regionPSet.getParameter<edm::InputTag>("pixelClustersForScaling");
0051     if (theOriginRScaling || thePtMinScaling || theHalfLengthScaling)
0052       token_pc = iC.consumes<edmNew::DetSetVector<SiPixelCluster> >(pixelClustersForScaling);
0053 
0054     if (theUseMS) {
0055       token_msmaker = iC.esConsumes();
0056     }
0057   }
0058 
0059   ~GlobalTrackingRegionWithVerticesProducer() override {}
0060 
0061   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0062     edm::ParameterSetDescription desc;
0063 
0064     desc.add<bool>("precise", true);
0065     desc.add<bool>("useMultipleScattering", false);
0066     desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
0067     desc.add<bool>("useFixedError", true);
0068     desc.add<double>("originRadius", 0.2);
0069     desc.add<double>("sigmaZVertex", 3.0);
0070     desc.add<double>("fixedError", 0.2);
0071     desc.add<edm::InputTag>("VertexCollection", edm::InputTag("firstStepPrimaryVertices"));
0072     desc.add<double>("ptMin", 0.9);
0073     desc.add<bool>("useFoundVertices", true);
0074     desc.add<bool>("useFakeVertices", false);
0075     desc.add<int>("maxNVertices", -1)->setComment("-1 for all vertices");
0076     desc.add<double>("nSigmaZ", 4.0);
0077     desc.add<edm::InputTag>("pixelClustersForScaling", edm::InputTag("siPixelClusters"));
0078     desc.add<bool>("originRScaling4BigEvts", false);
0079     desc.add<bool>("ptMinScaling4BigEvts", false);
0080     desc.add<bool>("halfLengthScaling4BigEvts", false);
0081     desc.add<bool>("allowEmpty", false);
0082     desc.add<double>("minOriginR", 0);
0083     desc.add<double>("maxPtMin", 1000);
0084     desc.add<double>("minHalfLength", 0);
0085     desc.add<double>("scalingStartNPix", 0.0);
0086     desc.add<double>("scalingEndNPix", 1.0);
0087 
0088     // Only for backwards-compatibility
0089     edm::ParameterSetDescription descRegion;
0090     descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
0091 
0092     descriptions.add("globalTrackingRegionWithVertices", descRegion);
0093   }
0094 
0095   std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
0096                                                         const edm::EventSetup& es) const override {
0097     std::vector<std::unique_ptr<TrackingRegion> > result;
0098 
0099     GlobalPoint theOrigin;
0100     edm::Handle<reco::BeamSpot> bsHandle;
0101     ev.getByToken(token_beamSpot, bsHandle);
0102     double bsSigmaZ;
0103     if (bsHandle.isValid()) {
0104       const reco::BeamSpot& bs = *bsHandle;
0105       bsSigmaZ = theNSigmaZ * bs.sigmaZ();
0106       theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
0107     } else {
0108       throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
0109     }
0110 
0111     const MultipleScatteringParametrisationMaker* msmaker = nullptr;
0112     if (theUseMS) {
0113       msmaker = &es.getData(token_msmaker);
0114     }
0115 
0116     if (theUseFoundVertices) {
0117       edm::Handle<reco::VertexCollection> vertexCollection;
0118       ev.getByToken(token_vertex, vertexCollection);
0119 
0120       edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusterDSV;
0121       bool doScaling = theOriginRScaling || thePtMinScaling || theHalfLengthScaling;
0122       if (doScaling)
0123         ev.getByToken(token_pc, pixelClusterDSV);
0124 
0125       for (reco::VertexCollection::const_iterator iV = vertexCollection->begin(); iV != vertexCollection->end(); iV++) {
0126         if (!iV->isValid())
0127           continue;
0128         if (iV->isFake() && !(theUseFakeVertices && theUseFixedError))
0129           continue;
0130         GlobalPoint theOrigin_ = GlobalPoint(iV->x(), iV->y(), iV->z());
0131 
0132         //scaling origin radius, half length, min pt for high-occupancy HI events to keep timing reasonable
0133         if (doScaling) {
0134           //Use the unscaled radius unless one of the two conditions below is met
0135           double scaledOriginRadius = theOriginRadius;
0136           double scaledHalfLength = theFixedError;
0137           double scaledPtMin = thePtMin;
0138 
0139           //calculate nPixels (adapted from TkSeedGenerator/src/ClusterChecker.cc)
0140           double nPix = 0;
0141 
0142           const edmNew::DetSetVector<SiPixelCluster>& input = *pixelClusterDSV;
0143           nPix = input.dataSize();
0144 
0145           //first condition is for high occupancy, second makes sure we won't divide by zero or a negative number
0146           if ((nPix > theScalingEnd) || ((theScalingEnd - theScalingStart) <= 0)) {
0147             if (theOriginRScaling)
0148               scaledOriginRadius = theMinOriginR;  // sets parameters to minimum value from PSet
0149             if (theHalfLengthScaling)
0150               scaledHalfLength = theMinHalfLength;
0151             if (thePtMinScaling)
0152               scaledPtMin = theMaxPtMin;
0153           }
0154           //second condition - scale radius linearly by Npix in the region from ScalingStart to ScalingEnd
0155           else if ((nPix <= theScalingEnd) && (nPix > theScalingStart)) {
0156             float slopeFactor = (nPix - theScalingStart) / (theScalingEnd - theScalingStart);
0157             if (theOriginRScaling)
0158               scaledOriginRadius = theOriginRadius - (theOriginRadius - theMinOriginR) * slopeFactor;
0159             if (theHalfLengthScaling)
0160               scaledHalfLength = theFixedError - (theFixedError - theMinHalfLength) * slopeFactor;
0161             if (thePtMinScaling)
0162               scaledPtMin = thePtMin - (thePtMin - theMaxPtMin) * slopeFactor;
0163           }
0164           //if region has 0 size, return 'result' empty, otherwise make a tracking region
0165           if (scaledOriginRadius != 0 && scaledHalfLength != 0) {
0166             result.push_back(std::make_unique<GlobalTrackingRegion>(
0167                 scaledPtMin, theOrigin_, scaledOriginRadius, scaledHalfLength, thePrecise, theUseMS, msmaker));
0168           }
0169         }  //end of region scaling code, pp behavior below
0170 
0171         else {
0172           double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError()) * theSigmaZVertex);
0173           result.push_back(std::make_unique<GlobalTrackingRegion>(
0174               thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise, theUseMS, msmaker));
0175           if (theMaxNVertices >= 0 && result.size() >= static_cast<unsigned>(theMaxNVertices))
0176             break;
0177         }
0178       }
0179 
0180       if (result.empty() && !(theOriginRScaling || thePtMinScaling || theHalfLengthScaling || theAllowEmpty)) {
0181         result.push_back(std::make_unique<GlobalTrackingRegion>(
0182             thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS, msmaker));
0183       }
0184     } else {
0185       result.push_back(std::make_unique<GlobalTrackingRegion>(
0186           thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS, msmaker));
0187     }
0188 
0189     return result;
0190   }
0191 
0192 private:
0193   double thePtMin;
0194   double theOriginRadius;
0195   double theNSigmaZ;
0196   edm::InputTag theBeamSpotTag;
0197 
0198   double theSigmaZVertex;
0199   double theFixedError;
0200   int theMaxNVertices;
0201   bool thePrecise;
0202   bool theUseMS;
0203 
0204   bool theUseFoundVertices;
0205   bool theUseFakeVertices;
0206   bool theUseFixedError;
0207   edm::EDGetTokenT<reco::VertexCollection> token_vertex;
0208   edm::EDGetTokenT<reco::BeamSpot> token_beamSpot;
0209   edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> token_msmaker;
0210 
0211   //HI-related variables
0212   bool theOriginRScaling;
0213   bool thePtMinScaling;
0214   bool theHalfLengthScaling;
0215   bool theAllowEmpty;
0216   double theMinOriginR;
0217   double theMaxPtMin;
0218   double theMinHalfLength;
0219   double theScalingStart;
0220   double theScalingEnd;
0221   edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > token_pc;
0222 };
0223 
0224 #endif