Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-03-05 03:16:38

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 #include "L1GTOptionalParam.h"
0022 
0023 #include <cinttypes>
0024 #include <memory>
0025 #include <vector>
0026 #include <set>
0027 
0028 #include <ap_int.h>
0029 
0030 using namespace l1t;
0031 
0032 class L1GTDoubleObjectCond : public edm::global::EDFilter<> {
0033 public:
0034   explicit L1GTDoubleObjectCond(const edm::ParameterSet&);
0035   ~L1GTDoubleObjectCond() override = default;
0036 
0037   static void fillDescriptions(edm::ConfigurationDescriptions&);
0038 
0039 private:
0040   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0041 
0042   const L1GTScales scales_;
0043 
0044   const L1GTSingleCollectionCut collection1Cuts_;
0045   const L1GTSingleCollectionCut collection2Cuts_;
0046 
0047   const bool enable_sanity_checks_;
0048   const bool inv_mass_checks_;
0049 
0050   const L1GTCorrelationalCut deltaCuts_;
0051 
0052   const std::optional<unsigned int> minQualityScoreSum_;
0053   const std::optional<unsigned int> maxQualityScoreSum_;
0054 
0055   const edm::EDGetTokenT<P2GTCandidateCollection> token1_;
0056   const edm::EDGetTokenT<P2GTCandidateCollection> token2_;
0057   const edm::EDGetTokenT<P2GTCandidateCollection> primVertToken_;
0058 };
0059 
0060 L1GTDoubleObjectCond::L1GTDoubleObjectCond(const edm::ParameterSet& config)
0061     : scales_(config.getParameter<edm::ParameterSet>("scales")),
0062       collection1Cuts_(config.getParameterSet("collection1"), config, scales_),
0063       collection2Cuts_(config.getParameterSet("collection2"), config, scales_),
0064       enable_sanity_checks_(config.getUntrackedParameter<bool>("sanity_checks")),
0065       inv_mass_checks_(config.getUntrackedParameter<bool>("inv_mass_checks")),
0066       deltaCuts_(config, config, scales_, enable_sanity_checks_, inv_mass_checks_),
0067       minQualityScoreSum_(getOptionalParam<unsigned int>("minQualityScoreSum", config)),
0068       maxQualityScoreSum_(getOptionalParam<unsigned int>("maxQualityScoreSum", config)),
0069       token1_(consumes<P2GTCandidateCollection>(collection1Cuts_.tag())),
0070       token2_(collection1Cuts_.tag() == collection2Cuts_.tag()
0071                   ? token1_
0072                   : consumes<P2GTCandidateCollection>(collection2Cuts_.tag())),
0073       primVertToken_(consumes<P2GTCandidateCollection>(config.getParameter<edm::InputTag>("primVertTag"))) {
0074   produces<P2GTCandidateVectorRef>(collection1Cuts_.tag().instance());
0075 
0076   if (!(collection1Cuts_.tag() == collection2Cuts_.tag())) {
0077     produces<P2GTCandidateVectorRef>(collection2Cuts_.tag().instance());
0078   }
0079 
0080   if (inv_mass_checks_) {
0081     produces<InvariantMassErrorCollection>();
0082   }
0083 
0084   if ((minQualityScoreSum_ || maxQualityScoreSum_) && !(collection1Cuts_.tag() == collection2Cuts_.tag())) {
0085     throw cms::Exception("Configuration") << "A qualityScore sum can only be calculated within one collection.";
0086   }
0087 }
0088 
0089 void L1GTDoubleObjectCond::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0090   edm::ParameterSetDescription desc;
0091 
0092   edm::ParameterSetDescription collection1Desc;
0093   L1GTSingleCollectionCut::fillPSetDescription(collection1Desc);
0094   desc.add<edm::ParameterSetDescription>("collection1", collection1Desc);
0095 
0096   edm::ParameterSetDescription collection2Desc;
0097   L1GTSingleCollectionCut::fillPSetDescription(collection2Desc);
0098   desc.add<edm::ParameterSetDescription>("collection2", collection2Desc);
0099 
0100   desc.add<edm::InputTag>("primVertTag");
0101 
0102   desc.addUntracked<bool>("sanity_checks", false);
0103   desc.addUntracked<bool>("inv_mass_checks", false);
0104 
0105   L1GTCorrelationalCut::fillPSetDescription(desc);
0106   L1GTCorrelationalCut::fillLUTDescriptions(desc);
0107 
0108   desc.addOptional<unsigned int>("minQualityScoreSum");
0109   desc.addOptional<unsigned int>("maxQualityScoreSum");
0110 
0111   edm::ParameterSetDescription scalesDesc;
0112   L1GTScales::fillPSetDescription(scalesDesc);
0113   desc.add<edm::ParameterSetDescription>("scales", scalesDesc);
0114 
0115   descriptions.addWithDefaultLabel(desc);
0116 }
0117 
0118 bool L1GTDoubleObjectCond::filter(edm::StreamID, edm::Event& event, const edm::EventSetup& setup) const {
0119   edm::Handle<P2GTCandidateCollection> col1 = event.getHandle(token1_);
0120   edm::Handle<P2GTCandidateCollection> col2 = event.getHandle(token2_);
0121   edm::Handle<P2GTCandidateCollection> primVertCol = event.getHandle(primVertToken_);
0122 
0123   bool condition_result = false;
0124 
0125   std::set<std::size_t> triggeredIdcs1;
0126   std::set<std::size_t> triggeredIdcs2;
0127 
0128   InvariantMassErrorCollection massErrors;
0129 
0130   for (std::size_t idx1 = 0; idx1 < col1->size(); ++idx1) {
0131     bool single1Pass = collection1Cuts_.checkObject(col1->at(idx1));
0132     single1Pass &= collection1Cuts_.checkPrimaryVertices(col1->at(idx1), *primVertCol);
0133 
0134     for (std::size_t idx2 = 0; idx2 < col2->size(); ++idx2) {
0135       // If we're looking at the same collection then we shouldn't use the same object in one comparison.
0136       if (col1.product() == col2.product() && idx1 == idx2) {
0137         continue;
0138       }
0139 
0140       bool pass = single1Pass;
0141       pass &= collection2Cuts_.checkObject(col2->at(idx2));
0142       pass &= collection2Cuts_.checkPrimaryVertices(col2->at(idx2), *primVertCol);
0143       pass &= deltaCuts_.checkObjects(col1->at(idx1), col2->at(idx2), massErrors);
0144 
0145       if (minQualityScoreSum_ || maxQualityScoreSum_) {
0146         unsigned int qualityScoreSum =
0147             col1->at(idx1).hwQualityScore().to_uint() + col2->at(idx2).hwQualityScore().to_uint();
0148 
0149         pass &= minQualityScoreSum_ ? qualityScoreSum > minQualityScoreSum_ : true;
0150         pass &= maxQualityScoreSum_ ? qualityScoreSum < maxQualityScoreSum_ : true;
0151       }
0152 
0153       condition_result |= pass;
0154 
0155       if (pass) {
0156         triggeredIdcs1.emplace(idx1);
0157         if (col1.product() != col2.product()) {
0158           triggeredIdcs2.emplace(idx2);
0159         } else {
0160           triggeredIdcs1.emplace(idx2);
0161         }
0162       }
0163     }
0164   }
0165 
0166   condition_result &= collection1Cuts_.checkCollection(*col1);
0167   condition_result &= collection2Cuts_.checkCollection(*col2);
0168 
0169   if (condition_result) {
0170     std::unique_ptr<P2GTCandidateVectorRef> triggerCol1 = std::make_unique<P2GTCandidateVectorRef>();
0171 
0172     for (std::size_t idx : triggeredIdcs1) {
0173       triggerCol1->push_back(P2GTCandidateRef(col1, idx));
0174     }
0175     event.put(std::move(triggerCol1), collection1Cuts_.tag().instance());
0176 
0177     if (col1.product() != col2.product()) {
0178       std::unique_ptr<P2GTCandidateVectorRef> triggerCol2 = std::make_unique<P2GTCandidateVectorRef>();
0179 
0180       for (std::size_t idx : triggeredIdcs2) {
0181         triggerCol2->push_back(P2GTCandidateRef(col2, idx));
0182       }
0183       event.put(std::move(triggerCol2), collection2Cuts_.tag().instance());
0184     }
0185   }
0186 
0187   if (inv_mass_checks_) {
0188     event.put(std::make_unique<InvariantMassErrorCollection>(std::move(massErrors)), "");
0189   }
0190 
0191   return condition_result;
0192 }
0193 
0194 DEFINE_FWK_MODULE(L1GTDoubleObjectCond);