File indexing completed on 2024-09-12 04:16:39
0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/global/EDFilter.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "DataFormats/L1Trigger/interface/P2GTCandidate.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006
0007 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "FWCore/Utilities/interface/EDGetToken.h"
0010 #include "DataFormats/Common/interface/Ref.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012
0013 #include "L1Trigger/Phase2L1GT/interface/L1GTInvariantMassError.h"
0014
0015 #include "L1Trigger/Phase2L1GT/interface/L1GTScales.h"
0016
0017 #include "L1GTOptionalParam.h"
0018 #include "L1GTSingleCollectionCut.h"
0019 #include "L1GTCorrelationalCut.h"
0020 #include "L1GTSingleInOutLUT.h"
0021
0022 #include <cinttypes>
0023 #include <memory>
0024 #include <vector>
0025 #include <set>
0026
0027 #include <ap_int.h>
0028
0029 using namespace l1t;
0030
0031 class L1GTDoubleObjectCond : public edm::global::EDFilter<> {
0032 public:
0033 explicit L1GTDoubleObjectCond(const edm::ParameterSet&);
0034 ~L1GTDoubleObjectCond() override = default;
0035
0036 static void fillDescriptions(edm::ConfigurationDescriptions&);
0037
0038 private:
0039 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0040
0041 const L1GTScales scales_;
0042
0043 const L1GTSingleCollectionCut collection1Cuts_;
0044 const L1GTSingleCollectionCut collection2Cuts_;
0045
0046 const bool enable_sanity_checks_;
0047 const bool inv_mass_checks_;
0048
0049 const L1GTCorrelationalCut deltaCuts_;
0050
0051 const edm::EDGetTokenT<P2GTCandidateCollection> token1_;
0052 const edm::EDGetTokenT<P2GTCandidateCollection> token2_;
0053 const edm::EDGetTokenT<P2GTCandidateCollection> primVertToken_;
0054 };
0055
0056 L1GTDoubleObjectCond::L1GTDoubleObjectCond(const edm::ParameterSet& config)
0057 : scales_(config.getParameter<edm::ParameterSet>("scales")),
0058 collection1Cuts_(config.getParameterSet("collection1"), config, scales_),
0059 collection2Cuts_(config.getParameterSet("collection2"), config, scales_),
0060 enable_sanity_checks_(config.getUntrackedParameter<bool>("sanity_checks")),
0061 inv_mass_checks_(config.getUntrackedParameter<bool>("inv_mass_checks")),
0062 deltaCuts_(config, config, scales_, enable_sanity_checks_, inv_mass_checks_),
0063 token1_(consumes<P2GTCandidateCollection>(collection1Cuts_.tag())),
0064 token2_(collection1Cuts_.tag() == collection2Cuts_.tag()
0065 ? token1_
0066 : consumes<P2GTCandidateCollection>(collection2Cuts_.tag())),
0067 primVertToken_(consumes<P2GTCandidateCollection>(config.getParameter<edm::InputTag>("primVertTag"))) {
0068 produces<P2GTCandidateVectorRef>(collection1Cuts_.tag().instance());
0069
0070 if (!(collection1Cuts_.tag() == collection2Cuts_.tag())) {
0071 produces<P2GTCandidateVectorRef>(collection2Cuts_.tag().instance());
0072 }
0073
0074 if (inv_mass_checks_) {
0075 produces<InvariantMassErrorCollection>();
0076 }
0077 }
0078
0079 void L1GTDoubleObjectCond::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0080 edm::ParameterSetDescription desc;
0081
0082 edm::ParameterSetDescription collection1Desc;
0083 L1GTSingleCollectionCut::fillPSetDescription(collection1Desc);
0084 desc.add<edm::ParameterSetDescription>("collection1", collection1Desc);
0085
0086 edm::ParameterSetDescription collection2Desc;
0087 L1GTSingleCollectionCut::fillPSetDescription(collection2Desc);
0088 desc.add<edm::ParameterSetDescription>("collection2", collection2Desc);
0089
0090 desc.add<edm::InputTag>("primVertTag");
0091
0092 desc.addUntracked<bool>("sanity_checks", false);
0093 desc.addUntracked<bool>("inv_mass_checks", false);
0094
0095 L1GTCorrelationalCut::fillPSetDescription(desc);
0096 L1GTCorrelationalCut::fillLUTDescriptions(desc);
0097
0098 edm::ParameterSetDescription scalesDesc;
0099 L1GTScales::fillPSetDescription(scalesDesc);
0100 desc.add<edm::ParameterSetDescription>("scales", scalesDesc);
0101
0102 descriptions.addWithDefaultLabel(desc);
0103 }
0104
0105 bool L1GTDoubleObjectCond::filter(edm::StreamID, edm::Event& event, const edm::EventSetup& setup) const {
0106 edm::Handle<P2GTCandidateCollection> col1 = event.getHandle(token1_);
0107 edm::Handle<P2GTCandidateCollection> col2 = event.getHandle(token2_);
0108 edm::Handle<P2GTCandidateCollection> primVertCol = event.getHandle(primVertToken_);
0109
0110 bool condition_result = false;
0111
0112 std::set<std::size_t> triggeredIdcs1;
0113 std::set<std::size_t> triggeredIdcs2;
0114
0115 InvariantMassErrorCollection massErrors;
0116
0117 for (std::size_t idx1 = 0; idx1 < col1->size(); ++idx1) {
0118 bool single1Pass = collection1Cuts_.checkObject(col1->at(idx1));
0119 single1Pass &= collection1Cuts_.checkPrimaryVertices(col1->at(idx1), *primVertCol);
0120
0121 for (std::size_t idx2 = 0; idx2 < col2->size(); ++idx2) {
0122
0123 if (col1.product() == col2.product() && idx1 == idx2) {
0124 continue;
0125 }
0126
0127 bool pass = single1Pass;
0128 pass &= collection2Cuts_.checkObject(col2->at(idx2));
0129 pass &= collection2Cuts_.checkPrimaryVertices(col2->at(idx2), *primVertCol);
0130 pass &= deltaCuts_.checkObjects(col1->at(idx1), col2->at(idx2), massErrors);
0131
0132 condition_result |= pass;
0133
0134 if (pass) {
0135 triggeredIdcs1.emplace(idx1);
0136 if (col1.product() != col2.product()) {
0137 triggeredIdcs2.emplace(idx2);
0138 } else {
0139 triggeredIdcs1.emplace(idx2);
0140 }
0141 }
0142 }
0143 }
0144
0145 condition_result &= collection1Cuts_.checkCollection(*col1);
0146 condition_result &= collection2Cuts_.checkCollection(*col2);
0147
0148 if (condition_result) {
0149 std::unique_ptr<P2GTCandidateVectorRef> triggerCol1 = std::make_unique<P2GTCandidateVectorRef>();
0150
0151 for (std::size_t idx : triggeredIdcs1) {
0152 triggerCol1->push_back(P2GTCandidateRef(col1, idx));
0153 }
0154 event.put(std::move(triggerCol1), collection1Cuts_.tag().instance());
0155
0156 if (col1.product() != col2.product()) {
0157 std::unique_ptr<P2GTCandidateVectorRef> triggerCol2 = std::make_unique<P2GTCandidateVectorRef>();
0158
0159 for (std::size_t idx : triggeredIdcs2) {
0160 triggerCol2->push_back(P2GTCandidateRef(col2, idx));
0161 }
0162 event.put(std::move(triggerCol2), collection2Cuts_.tag().instance());
0163 }
0164 }
0165
0166 if (inv_mass_checks_) {
0167 event.put(std::make_unique<InvariantMassErrorCollection>(std::move(massErrors)), "");
0168 }
0169
0170 return condition_result;
0171 }
0172
0173 DEFINE_FWK_MODULE(L1GTDoubleObjectCond);