File indexing completed on 2024-04-06 12:18:42
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "FWCore/Framework/interface/global/EDFilter.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017
0018 #include "DataFormats/Common/interface/Ref.h"
0019 #include "DataFormats/Common/interface/Handle.h"
0020 #include "DataFormats/Common/interface/DetSet.h"
0021 #include "DataFormats/Common/interface/DetSetVector.h"
0022
0023 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0024 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0025
0026 #include "DataFormats/CTPPSDetId/interface/CTPPSPixelDetId.h"
0027
0028 #include "DataFormats/CTPPSReco/interface/CTPPSPixelLocalTrack.h" // pixel
0029
0030 #include <unordered_map>
0031
0032 class HLTPPSCalFilter : public edm::global::EDFilter<> {
0033 public:
0034 explicit HLTPPSCalFilter(const edm::ParameterSet&);
0035 ~HLTPPSCalFilter() 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 const edm::InputTag pixelLocalTrackInputTag_;
0042 const edm::EDGetTokenT<edm::DetSetVector<CTPPSPixelLocalTrack>> pixelLocalTrackToken_;
0043
0044 const int minTracks_;
0045 const int maxTracks_;
0046
0047 const bool do_express_;
0048 };
0049
0050 void HLTPPSCalFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0051 edm::ParameterSetDescription desc;
0052
0053 desc.add<edm::InputTag>("pixelLocalTrackInputTag", edm::InputTag("ctppsPixelLocalTracks"))
0054 ->setComment("input tag of the pixel local track collection");
0055
0056 desc.add<int>("minTracks", 1)->setComment("minimum number of tracks in pot");
0057 desc.add<int>("maxTracks", -1)->setComment("maximum number of tracks in pot");
0058
0059 desc.add<bool>("do_express", true)->setComment("toggle on filter type; true for Express, false for Prompt");
0060
0061 desc.add<int>("triggerType", trigger::TriggerTrack);
0062
0063 descriptions.add("hltPPSCalFilter", desc);
0064 }
0065
0066 HLTPPSCalFilter::HLTPPSCalFilter(const edm::ParameterSet& iConfig)
0067 : pixelLocalTrackInputTag_(iConfig.getParameter<edm::InputTag>("pixelLocalTrackInputTag")),
0068 pixelLocalTrackToken_(consumes<edm::DetSetVector<CTPPSPixelLocalTrack>>(pixelLocalTrackInputTag_)),
0069 minTracks_(iConfig.getParameter<int>("minTracks")),
0070 maxTracks_(iConfig.getParameter<int>("maxTracks")),
0071 do_express_(iConfig.getParameter<bool>("do_express")) {}
0072
0073 bool HLTPPSCalFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0074
0075 const auto valid_trks = [](const auto& trk) { return trk.isValid(); };
0076
0077
0078 std::unordered_map<uint32_t, bool> pixel_filter_;
0079
0080
0081 pixel_filter_[CTPPSPixelDetId(0, 0, 3)] = false;
0082 pixel_filter_[CTPPSPixelDetId(0, 2, 3)] = false;
0083 pixel_filter_[CTPPSPixelDetId(1, 0, 3)] = false;
0084 pixel_filter_[CTPPSPixelDetId(1, 2, 3)] = false;
0085
0086
0087 edm::Handle<edm::DetSetVector<CTPPSPixelLocalTrack>> pixelTracks;
0088 iEvent.getByToken(pixelLocalTrackToken_, pixelTracks);
0089
0090 for (const auto& rpv : *pixelTracks) {
0091 if (pixel_filter_.count(rpv.id) == 0) {
0092 continue;
0093 }
0094
0095 pixel_filter_.at(rpv.id) = true;
0096
0097
0098 const auto ntrks = std::count_if(rpv.begin(), rpv.end(), valid_trks);
0099
0100 if (minTracks_ > 0 && ntrks < minTracks_)
0101 pixel_filter_.at(rpv.id) = false;
0102 if (maxTracks_ > 0 && ntrks > maxTracks_)
0103 pixel_filter_.at(rpv.id) = false;
0104 }
0105
0106
0107 if (do_express_) {
0108 return (pixel_filter_.at(CTPPSPixelDetId(0, 0, 3)) && pixel_filter_.at(CTPPSPixelDetId(0, 2, 3))) ||
0109 (pixel_filter_.at(CTPPSPixelDetId(1, 0, 3)) && pixel_filter_.at(CTPPSPixelDetId(1, 2, 3)));
0110 } else {
0111 return (pixel_filter_.at(CTPPSPixelDetId(0, 0, 3)) || pixel_filter_.at(CTPPSPixelDetId(0, 2, 3))) ||
0112 (pixel_filter_.at(CTPPSPixelDetId(1, 0, 3)) || pixel_filter_.at(CTPPSPixelDetId(1, 2, 3)));
0113 }
0114
0115 return false;
0116 }
0117
0118
0119 #include "FWCore/Framework/interface/MakerMacros.h"
0120 DEFINE_FWK_MODULE(HLTPPSCalFilter);