Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/global/EDProducer.h"
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "FWCore/Utilities/interface/InputTag.h"
0011 #include "FWCore/Utilities/interface/EDGetToken.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "RecoTracker/PixelTrackFitting/interface/PixelTrackFilter.h"
0015 #include "RecoHI/HiTracking/interface/HIProtoTrackFilter.h"
0016 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0017 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0018 #include "DataFormats/Common/interface/DetSetAlgorithm.h"
0019 
0020 #include "DataFormats/Common/interface/DetSetVector.h"
0021 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0022 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0023 
0024 class HIProtoTrackFilterProducer : public edm::global::EDProducer<> {
0025 public:
0026   explicit HIProtoTrackFilterProducer(const edm::ParameterSet& iConfig);
0027   ~HIProtoTrackFilterProducer() override;
0028 
0029   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0030 
0031 private:
0032   void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0033 
0034   edm::EDGetTokenT<reco::BeamSpot> theBeamSpotToken;
0035   edm::EDGetTokenT<SiPixelRecHitCollection> theSiPixelRecHitsToken;
0036   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> theTtopoToken;
0037   double theTIPMax;
0038   double theChi2Max, thePtMin;
0039   bool doVariablePtMin;
0040 };
0041 
0042 HIProtoTrackFilterProducer::HIProtoTrackFilterProducer(const edm::ParameterSet& iConfig)
0043     : theBeamSpotToken(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
0044       theSiPixelRecHitsToken(consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("siPixelRecHits"))),
0045       theTtopoToken(esConsumes()),
0046       theTIPMax(iConfig.getParameter<double>("tipMax")),
0047       theChi2Max(iConfig.getParameter<double>("chi2")),
0048       thePtMin(iConfig.getParameter<double>("ptMin")),
0049       doVariablePtMin(iConfig.getParameter<bool>("doVariablePtMin")) {
0050   produces<PixelTrackFilter>();
0051 }
0052 
0053 void HIProtoTrackFilterProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0054   edm::ParameterSetDescription desc;
0055 
0056   desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
0057   desc.add<edm::InputTag>("siPixelRecHits", edm::InputTag("siPixelRecHits"));
0058   desc.add<double>("ptMin", 1.0);
0059   desc.add<double>("tipMax", 1.0);
0060   desc.add<double>("chi2", 1000);
0061   desc.add<bool>("doVariablePtMin", true);
0062 
0063   descriptions.add("hiProtoTrackFilter", desc);
0064 }
0065 
0066 HIProtoTrackFilterProducer::~HIProtoTrackFilterProducer() {}
0067 
0068 void HIProtoTrackFilterProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0069   // Get the beam spot
0070   edm::Handle<reco::BeamSpot> bsHandle;
0071   iEvent.getByToken(theBeamSpotToken, bsHandle);
0072   const reco::BeamSpot* beamSpot = bsHandle.product();
0073 
0074   if (beamSpot) {
0075     edm::LogInfo("HeavyIonVertexing") << "[HIProtoTrackFilterProducer] Proto track selection based on beamspot"
0076                                       << "\n   (x,y,z) = (" << beamSpot->x0() << "," << beamSpot->y0() << ","
0077                                       << beamSpot->z0() << ")";
0078   } else {
0079     edm::EDConsumerBase::Labels labels;
0080     labelsForToken(theBeamSpotToken, labels);
0081     edm::LogError("HeavyIonVertexing")  // this can be made a warning when operator() is fixed
0082         << "No beamspot found with tag '" << labels.module << "'";
0083   }
0084 
0085   // Estimate multiplicity
0086   edm::Handle<SiPixelRecHitCollection> recHitColl;
0087   iEvent.getByToken(theSiPixelRecHitsToken, recHitColl);
0088 
0089   auto const& ttopo = iSetup.getData(theTtopoToken);
0090 
0091   std::vector<const TrackingRecHit*> theChosenHits;
0092   edmNew::copyDetSetRange(*recHitColl, theChosenHits, ttopo.pxbDetIdLayerComparator(1));
0093   float estMult = theChosenHits.size();
0094 
0095   double variablePtMin = thePtMin;
0096   if (doVariablePtMin) {
0097     // parameterize ptMin such that a roughly constant number of selected prototracks passed are to vertexing
0098     float varPtCutoff = 1500;  //cutoff for variable ptMin
0099     if (estMult < varPtCutoff) {
0100       variablePtMin = 0.075;
0101       if (estMult > 0)
0102         variablePtMin = (13. - (varPtCutoff / estMult)) / 12.;
0103       if (variablePtMin < 0.075)
0104         variablePtMin = 0.075;  // don't lower the cut past 75 MeV
0105     }
0106     LogTrace("heavyIonHLTVertexing") << "   [HIProtoTrackFilterProducer: variablePtMin: " << variablePtMin << "]";
0107   }
0108 
0109   auto impl = std::make_unique<HIProtoTrackFilter>(beamSpot, theTIPMax, theChi2Max, variablePtMin);
0110   auto prod = std::make_unique<PixelTrackFilter>(std::move(impl));
0111   iEvent.put(std::move(prod));
0112 }
0113 
0114 DEFINE_FWK_MODULE(HIProtoTrackFilterProducer);