Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:42

0001 // <author>Laurent Forthomme</author>
0002 // <email>lforthom@cern.ch</email>
0003 // <created>2020-02-17</created>
0004 // <description>
0005 // HLT filter module to select events with tracks in the PPS detector
0006 // </description>
0007 
0008 #include "FWCore/Framework/interface/global/EDFilter.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/Utilities/interface/InputTag.h"
0015 
0016 #include "DataFormats/Common/interface/Ref.h"
0017 #include "DataFormats/Common/interface/Handle.h"
0018 #include "DataFormats/Common/interface/DetSet.h"
0019 #include "DataFormats/Common/interface/DetSetVector.h"
0020 
0021 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0022 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0023 
0024 #include "DataFormats/CTPPSDetId/interface/CTPPSPixelDetId.h"
0025 
0026 #include "DataFormats/CTPPSReco/interface/CTPPSPixelLocalTrack.h"    // pixel
0027 #include "DataFormats/CTPPSReco/interface/TotemRPLocalTrack.h"       // strip
0028 #include "DataFormats/CTPPSReco/interface/CTPPSDiamondLocalTrack.h"  // diamond
0029 
0030 #include <unordered_map>
0031 
0032 class HLTPPSPerPotTrackFilter : public edm::global::EDFilter<> {
0033 public:
0034   explicit HLTPPSPerPotTrackFilter(const edm::ParameterSet&);
0035   ~HLTPPSPerPotTrackFilter() override = default;
0036 
0037   static void fillDescriptions(edm::ConfigurationDescriptions&);
0038   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0039 
0040 private:
0041   edm::EDGetTokenT<edm::DetSetVector<CTPPSPixelLocalTrack>> pixelLocalTrackToken_;
0042   edm::EDGetTokenT<edm::DetSetVector<TotemRPLocalTrack>> stripLocalTrackToken_;
0043   edm::EDGetTokenT<edm::DetSetVector<CTPPSDiamondLocalTrack>> diamondLocalTrackToken_;
0044 
0045   struct PerPotFilter {
0046     int minTracks;
0047     int maxTracks;
0048   };
0049   std::unordered_map<uint32_t, PerPotFilter> pixel_filter_;
0050   std::unordered_map<uint32_t, PerPotFilter> strip_filter_;
0051   std::unordered_map<uint32_t, PerPotFilter> diam_filter_;
0052 
0053   // Helper tool to count valid tracks
0054   static constexpr auto valid_trks_ = [](const auto& trk) { return trk.isValid(); };
0055 };
0056 
0057 void HLTPPSPerPotTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0058   edm::ParameterSetDescription desc;
0059 
0060   desc.add<edm::InputTag>("pixelLocalTrackInputTag", edm::InputTag("ctppsPixelLocalTracks"))
0061       ->setComment("input tag of the pixel local track collection");
0062   desc.add<edm::InputTag>("stripLocalTrackInputTag", edm::InputTag("totemRPLocalTrackFitter"))
0063       ->setComment("input tag of the strip local track collection");
0064   desc.add<edm::InputTag>("diamondLocalTrackInputTag", edm::InputTag("ctppsDiamondLocalTracks"))
0065       ->setComment("input tag of the diamond local track collection");
0066 
0067   edm::ParameterSetDescription filterValid;
0068   filterValid.add<unsigned int>("detid", 0)->setComment("station/pot raw DetId");
0069   filterValid.add<int>("minTracks", -1)->setComment("minimum number of tracks in pot");
0070   filterValid.add<int>("maxTracks", -1)->setComment("maximum number of tracks in pot");
0071 
0072   std::vector<edm::ParameterSet> vPixelDefault;
0073   auto& near_pix45 = vPixelDefault.emplace_back();
0074   near_pix45.addParameter<unsigned int>("detid", CTPPSPixelDetId(0, 2, 2));  // arm-station-rp
0075   near_pix45.addParameter<int>("minTracks", 2);
0076   near_pix45.addParameter<int>("maxTracks", -1);
0077   auto& far_pix45 = vPixelDefault.emplace_back();
0078   far_pix45.addParameter<unsigned int>("detid", CTPPSPixelDetId(0, 2, 3));  // arm-station-rp
0079   far_pix45.addParameter<int>("minTracks", 2);
0080   far_pix45.addParameter<int>("maxTracks", -1);
0081   auto& near_pix56 = vPixelDefault.emplace_back();
0082   near_pix56.addParameter<unsigned int>("detid", CTPPSPixelDetId(1, 2, 2));  // arm-station-rp
0083   near_pix56.addParameter<int>("minTracks", 2);
0084   near_pix56.addParameter<int>("maxTracks", -1);
0085   auto& far_pix56 = vPixelDefault.emplace_back();
0086   far_pix56.addParameter<unsigned int>("detid", CTPPSPixelDetId(1, 2, 3));  // arm-station-rp
0087   far_pix56.addParameter<int>("minTracks", 2);
0088   far_pix56.addParameter<int>("maxTracks", -1);
0089   desc.addVPSet("pixelFilter", filterValid, vPixelDefault);
0090 
0091   std::vector<edm::ParameterSet> vStripDefault;
0092   desc.addVPSet("stripFilter", filterValid, vStripDefault);
0093 
0094   std::vector<edm::ParameterSet> vDiamDefault;
0095   desc.addVPSet("diamondFilter", filterValid, vDiamDefault);
0096 
0097   desc.add<int>("triggerType", trigger::TriggerTrack);
0098 
0099   descriptions.add("hltPPSPerPotTrackFilter", desc);
0100 }
0101 
0102 HLTPPSPerPotTrackFilter::HLTPPSPerPotTrackFilter(const edm::ParameterSet& iConfig) {
0103   // First define pixels (2017+) selection
0104   const auto& pixelVPset = iConfig.getParameter<std::vector<edm::ParameterSet>>("pixelFilter");
0105   if (!pixelVPset.empty()) {
0106     pixelLocalTrackToken_ = consumes<edm::DetSetVector<CTPPSPixelLocalTrack>>(
0107         iConfig.getParameter<edm::InputTag>("pixelLocalTrackInputTag"));
0108     for (const auto& pset : pixelVPset)
0109       pixel_filter_[pset.getParameter<unsigned int>("detid")] =
0110           PerPotFilter{pset.getParameter<int>("minTracks"), pset.getParameter<int>("maxTracks")};
0111   }
0112   // Then define strips (2016-17) selection
0113   const auto& stripVPset = iConfig.getParameter<std::vector<edm::ParameterSet>>("stripFilter");
0114   if (!stripVPset.empty()) {
0115     stripLocalTrackToken_ =
0116         consumes<edm::DetSetVector<TotemRPLocalTrack>>(iConfig.getParameter<edm::InputTag>("stripLocalTrackInputTag"));
0117     for (const auto& pset : stripVPset)
0118       strip_filter_[pset.getParameter<unsigned int>("detid")] =
0119           PerPotFilter{pset.getParameter<int>("minTracks"), pset.getParameter<int>("maxTracks")};
0120   }
0121   // Finally define diamond (2016+) selection
0122   const auto& diamVPset = iConfig.getParameter<std::vector<edm::ParameterSet>>("diamondFilter");
0123   if (!diamVPset.empty()) {
0124     diamondLocalTrackToken_ = consumes<edm::DetSetVector<CTPPSDiamondLocalTrack>>(
0125         iConfig.getParameter<edm::InputTag>("diamondLocalTrackInputTag"));
0126     for (const auto& pset : diamVPset)
0127       diam_filter_[pset.getParameter<unsigned int>("detid")] =
0128           PerPotFilter{pset.getParameter<int>("minTracks"), pset.getParameter<int>("maxTracks")};
0129   }
0130 }
0131 
0132 bool HLTPPSPerPotTrackFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0133   // First filter on pixels
0134   if (!pixel_filter_.empty()) {
0135     edm::Handle<edm::DetSetVector<CTPPSPixelLocalTrack>> pixelTracks;
0136     iEvent.getByToken(pixelLocalTrackToken_, pixelTracks);
0137 
0138     for (const auto& rpv : *pixelTracks) {
0139       if (pixel_filter_.count(rpv.id) == 0)
0140         continue;
0141       const auto& fltr = pixel_filter_.at(rpv.id);
0142 
0143       const auto ntrks = std::count_if(rpv.begin(), rpv.end(), valid_trks_);
0144       if (fltr.minTracks > 0 && ntrks < fltr.minTracks)
0145         return false;
0146       if (fltr.maxTracks > 0 && ntrks > fltr.maxTracks)
0147         return false;
0148     }
0149   }
0150 
0151   // Then filter on strips
0152   if (!strip_filter_.empty()) {
0153     edm::Handle<edm::DetSetVector<TotemRPLocalTrack>> stripTracks;
0154     iEvent.getByToken(stripLocalTrackToken_, stripTracks);
0155 
0156     for (const auto& rpv : *stripTracks) {
0157       if (strip_filter_.count(rpv.id) == 0)
0158         continue;
0159       const auto& fltr = strip_filter_.at(rpv.id);
0160 
0161       const auto ntrks = std::count_if(rpv.begin(), rpv.end(), valid_trks_);
0162       if (fltr.minTracks > 0 && ntrks < fltr.minTracks)
0163         return false;
0164       if (fltr.maxTracks > 0 && ntrks > fltr.maxTracks)
0165         return false;
0166     }
0167   }
0168 
0169   // Finally filter on diamond
0170   if (!diam_filter_.empty()) {
0171     edm::Handle<edm::DetSetVector<CTPPSDiamondLocalTrack>> diamondTracks;
0172     iEvent.getByToken(diamondLocalTrackToken_, diamondTracks);
0173 
0174     for (const auto& rpv : *diamondTracks) {
0175       if (diam_filter_.count(rpv.id) == 0)
0176         continue;
0177       const auto& fltr = diam_filter_.at(rpv.id);
0178 
0179       const auto ntrks = std::count_if(rpv.begin(), rpv.end(), valid_trks_);
0180       if (fltr.minTracks > 0 && ntrks < fltr.minTracks)
0181         return false;
0182       if (fltr.maxTracks > 0 && ntrks > fltr.maxTracks)
0183         return false;
0184     }
0185   }
0186 
0187   return true;
0188 }
0189 
0190 // define as a framework module
0191 #include "FWCore/Framework/interface/MakerMacros.h"
0192 DEFINE_FWK_MODULE(HLTPPSPerPotTrackFilter);