Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-09-19 00:12:00

0001 /** \class L1TTkEmFilter
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Martin Grunewald
0007  *  \author Simone Gennai
0008  *  \author Thiago Tomei
0009  *
0010  */
0011 
0012 #include "L1TTkEmFilter.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 
0015 #include "DataFormats/Common/interface/Handle.h"
0016 
0017 #include "DataFormats/Common/interface/Ref.h"
0018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0019 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0020 
0021 #include "FWCore/Framework/interface/EventSetup.h"
0022 #include "FWCore/Framework/interface/EventSetupRecord.h"
0023 #include "FWCore/Utilities/interface/Exception.h"
0024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0025 
0026 //
0027 // constructors and destructor
0028 //
0029 
0030 L1TTkEmFilter::L1TTkEmFilter(const edm::ParameterSet& iConfig)
0031     : HLTFilter(iConfig),
0032       l1TkEmTag1_(iConfig.getParameter<edm::InputTag>("inputTag1")),
0033       l1TkEmTag2_(iConfig.getParameter<edm::InputTag>("inputTag2")),
0034       tkEmToken1_(consumes<TkEmCollection>(l1TkEmTag1_)),
0035       tkEmToken2_(consumes<TkEmCollection>(l1TkEmTag2_)) {
0036   min_Pt_ = iConfig.getParameter<double>("MinPt");
0037   min_N_ = iConfig.getParameter<int>("MinN");
0038   min_AbsEta1_ = iConfig.getParameter<double>("MinAbsEta1");
0039   max_AbsEta1_ = iConfig.getParameter<double>("MaxAbsEta1");
0040   min_AbsEta2_ = iConfig.getParameter<double>("MinAbsEta2");
0041   max_AbsEta2_ = iConfig.getParameter<double>("MaxAbsEta2");
0042   scalings_ = iConfig.getParameter<edm::ParameterSet>("Scalings");
0043   barrelScalings_ = scalings_.getParameter<std::vector<double> >("barrel");
0044   endcapScalings_ = scalings_.getParameter<std::vector<double> >("endcap");
0045   etaBinsForIsolation_ = iConfig.getParameter<std::vector<double> >("EtaBinsForIsolation");
0046   etaBinsForIsolation_ = iConfig.getParameter<std::vector<double> >("EtaBinsForIsolation");
0047   trkIsolation_ = iConfig.getParameter<std::vector<double> >("TrkIsolation");
0048   quality1_ = iConfig.getParameter<int>("Quality1");
0049   quality2_ = iConfig.getParameter<int>("Quality2");
0050   qual1IsMask_ = iConfig.getParameter<bool>("Qual1IsMask");
0051   qual2IsMask_ = iConfig.getParameter<bool>("Qual2IsMask");
0052   applyQual1_ = iConfig.getParameter<bool>("ApplyQual1");
0053   applyQual2_ = iConfig.getParameter<bool>("ApplyQual2");
0054 
0055   if (etaBinsForIsolation_.size() != (trkIsolation_.size() + 1))
0056     throw cms::Exception("ConfigurationError")
0057         << "Vector of isolation values should have same size of vector of eta bins plus one.";
0058 }
0059 
0060 L1TTkEmFilter::~L1TTkEmFilter() = default;
0061 
0062 //
0063 // member functions
0064 //
0065 
0066 void L1TTkEmFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0067   edm::ParameterSetDescription desc;
0068   makeHLTFilterDescription(desc);
0069   desc.add<double>("MinPt", -1.0);
0070   desc.add<double>("MinAbsEta1", 0.0);
0071   desc.add<double>("MaxAbsEta1", 1.479);
0072   desc.add<double>("MinAbsEta2", 1.479);
0073   desc.add<double>("MaxAbsEta2", 5.0);
0074   desc.add<int>("MinN", 1);
0075   desc.add<edm::InputTag>("inputTag1", edm::InputTag("L1TkEms1"));
0076   desc.add<edm::InputTag>("inputTag2", edm::InputTag("L1TkEms2"));
0077   desc.add<std::vector<double> >("EtaBinsForIsolation", {0.0, 1.479});
0078   desc.add<std::vector<double> >("TrkIsolation",
0079                                  {
0080                                      99999.9,
0081                                  });
0082   desc.add<int>("Quality1", 0);
0083   desc.add<int>("Quality2", 0);
0084   desc.add<bool>("Qual1IsMask", false);
0085   desc.add<bool>("Qual2IsMask", false);
0086   desc.add<bool>("ApplyQual1", false);
0087   desc.add<bool>("ApplyQual2", false);
0088 
0089   edm::ParameterSetDescription descScalings;
0090   descScalings.add<std::vector<double> >("barrel", {0.0, 1.0, 0.0});
0091   descScalings.add<std::vector<double> >("endcap", {0.0, 1.0, 0.0});
0092   desc.add<edm::ParameterSetDescription>("Scalings", descScalings);
0093 
0094   descriptions.add("L1TTkEmFilter", desc);
0095 }
0096 
0097 // ------------ method called to produce the data  ------------
0098 bool L1TTkEmFilter::hltFilter(edm::Event& iEvent,
0099                               const edm::EventSetup& iSetup,
0100                               trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0101   using namespace std;
0102   using namespace edm;
0103   using namespace reco;
0104   using namespace trigger;
0105 
0106   // All HLT filters must create and fill an HLT filter object,
0107   // recording any reconstructed physics objects satisfying (or not)
0108   // this HLT filter, and place it in the Event.
0109 
0110   // The filter object
0111   if (saveTags()) {
0112     filterproduct.addCollectionTag(l1TkEmTag1_);
0113     filterproduct.addCollectionTag(l1TkEmTag2_);
0114   }
0115 
0116   // Specific filter code
0117 
0118   // get hold of products from Event
0119 
0120   /// Barrel collection
0121   Handle<l1t::TkEmCollection> tkEms1;
0122   iEvent.getByToken(tkEmToken1_, tkEms1);
0123 
0124   /// Endcap collection
0125   Handle<l1t::TkEmCollection> tkEms2;
0126   iEvent.getByToken(tkEmToken2_, tkEms2);
0127 
0128   int ntrkEm(0);
0129   // Loop over first collection
0130   auto atrkEms(tkEms1->begin());
0131   auto otrkEms(tkEms1->end());
0132   TkEmCollection::const_iterator itkEm;
0133   for (itkEm = atrkEms; itkEm != otrkEms; itkEm++) {
0134     double offlinePt = this->TkEmOfflineEt(itkEm->pt(), itkEm->eta());
0135     bool passQuality(false);
0136     bool passIsolation(false);
0137 
0138     if (applyQual1_) {
0139       if (qual1IsMask_) {
0140         passQuality = (itkEm->hwQual() & quality1_);
0141       } else {
0142         passQuality = (itkEm->hwQual() == quality1_);
0143       }
0144     } else
0145       passQuality = true;
0146 
0147     // There has to be a better way to do this.
0148     for (unsigned int etabin = 1; etabin != etaBinsForIsolation_.size(); ++etabin) {
0149       if (std::abs(itkEm->eta()) < etaBinsForIsolation_.at(etabin) and
0150           std::abs(itkEm->eta()) > etaBinsForIsolation_.at(etabin - 1) and
0151           itkEm->trkIsol() < trkIsolation_.at(etabin - 1))
0152         passIsolation = true;
0153     }
0154 
0155     if (offlinePt >= min_Pt_ && std::abs(itkEm->eta()) < max_AbsEta1_ && std::abs(itkEm->eta()) >= min_AbsEta1_ &&
0156         passQuality && passIsolation) {
0157       ntrkEm++;
0158       l1t::TkEmRef ref1(l1t::TkEmRef(tkEms1, distance(atrkEms, itkEm)));
0159       filterproduct.addObject(trigger::TriggerObjectType::TriggerL1TkEm, ref1);
0160     }
0161   }
0162 
0163   // Loop over second collection. Notice we don't reset ntrkEm
0164   atrkEms = tkEms2->begin();
0165   otrkEms = tkEms2->end();
0166   for (itkEm = atrkEms; itkEm != otrkEms; itkEm++) {
0167     double offlinePt = this->TkEmOfflineEt(itkEm->pt(), itkEm->eta());
0168     bool passQuality(false);
0169     bool passIsolation(false);
0170 
0171     if (applyQual2_) {
0172       if (qual2IsMask_) {
0173         passQuality = (itkEm->hwQual() & quality2_);
0174       } else {
0175         passQuality = (itkEm->hwQual() == quality2_);
0176       }
0177     } else
0178       passQuality = true;
0179 
0180     for (unsigned int etabin = 1; etabin != etaBinsForIsolation_.size(); ++etabin) {
0181       if (std::abs(itkEm->eta()) < etaBinsForIsolation_.at(etabin) and
0182           std::abs(itkEm->eta()) > etaBinsForIsolation_.at(etabin - 1) and
0183           itkEm->trkIsol() < trkIsolation_.at(etabin - 1))
0184         passIsolation = true;
0185     }
0186 
0187     if (offlinePt >= min_Pt_ && std::abs(itkEm->eta()) <= max_AbsEta2_ && std::abs(itkEm->eta()) >= min_AbsEta2_ &&
0188         passQuality && passIsolation) {
0189       ntrkEm++;
0190       l1t::TkEmRef ref2(l1t::TkEmRef(tkEms2, distance(atrkEms, itkEm)));
0191       filterproduct.addObject(trigger::TriggerObjectType::TriggerL1TkEm, ref2);
0192     }
0193   }
0194 
0195   // return with final filter decision
0196   const bool accept(ntrkEm >= min_N_);
0197   return accept;
0198 }
0199 
0200 double L1TTkEmFilter::TkEmOfflineEt(double Et, double Eta) const {
0201   if (std::abs(Eta) < 1.5)
0202     return (barrelScalings_.at(0) + Et * barrelScalings_.at(1) + Et * Et * barrelScalings_.at(2));
0203   else
0204     return (endcapScalings_.at(0) + Et * endcapScalings_.at(1) + Et * Et * endcapScalings_.at(2));
0205 }