Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:53:05

0001 /** \class L1TTkEleFilter
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 "L1TTkEleFilter.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 L1TTkEleFilter::L1TTkEleFilter(const edm::ParameterSet& iConfig)
0031     : HLTFilter(iConfig),
0032       l1TkEleTag1_(iConfig.getParameter<edm::InputTag>("inputTag1")),
0033       l1TkEleTag2_(iConfig.getParameter<edm::InputTag>("inputTag2")),
0034       tkEleToken1_(consumes<TkEleCollection>(l1TkEleTag1_)),
0035       tkEleToken2_(consumes<TkEleCollection>(l1TkEleTag2_)) {
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   trkIsolation_ = iConfig.getParameter<std::vector<double> >("TrkIsolation");
0047   quality1_ = iConfig.getParameter<int>("Quality1");
0048   quality2_ = iConfig.getParameter<int>("Quality2");
0049   qual1IsMask_ = iConfig.getParameter<bool>("Qual1IsMask");
0050   qual2IsMask_ = iConfig.getParameter<bool>("Qual2IsMask");
0051   applyQual1_ = iConfig.getParameter<bool>("ApplyQual1");
0052   applyQual2_ = iConfig.getParameter<bool>("ApplyQual2");
0053 
0054   if (etaBinsForIsolation_.size() != (trkIsolation_.size() + 1))
0055     throw cms::Exception("ConfigurationError")
0056         << "Vector of isolation values should have same size of vector of eta bins plus one.";
0057 }
0058 
0059 L1TTkEleFilter::~L1TTkEleFilter() = default;
0060 
0061 //
0062 // member functions
0063 //
0064 
0065 void L1TTkEleFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0066   edm::ParameterSetDescription desc;
0067 
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("L1TkElectrons1"));
0076   desc.add<edm::InputTag>("inputTag2", edm::InputTag("L1TkElectrons2"));
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("L1TTkEleFilter", desc);
0095 }
0096 
0097 // ------------ method called to produce the data  ------------
0098 bool L1TTkEleFilter::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(l1TkEleTag1_);
0113     filterproduct.addCollectionTag(l1TkEleTag2_);
0114   }
0115 
0116   // Specific filter code
0117 
0118   // get hold of products from Event
0119 
0120   /// Barrel colleciton
0121   Handle<l1t::TkElectronCollection> tkEles1;
0122   iEvent.getByToken(tkEleToken1_, tkEles1);
0123 
0124   /// Endcap collection
0125   Handle<l1t::TkElectronCollection> tkEles2;
0126   iEvent.getByToken(tkEleToken2_, tkEles2);
0127 
0128   int ntrkEle(0);
0129   // Loop over first collection
0130   auto atrkEles(tkEles1->begin());
0131   auto otrkEles(tkEles1->end());
0132   TkEleCollection::const_iterator itkEle;
0133   for (itkEle = atrkEles; itkEle != otrkEles; itkEle++) {
0134     double offlinePt = this->TkEleOfflineEt(itkEle->pt(), itkEle->eta());
0135     bool passQuality(false);
0136     bool passIsolation(false);
0137 
0138     if (applyQual1_) {
0139       if (qual1IsMask_)
0140         passQuality = (itkEle->hwQual() & quality1_);
0141       else
0142         passQuality = (itkEle->hwQual() == quality1_);
0143     } else
0144       passQuality = true;
0145 
0146     // There has to be a better way to do this.
0147     for (unsigned int etabin = 1; etabin != etaBinsForIsolation_.size(); ++etabin) {
0148       if (std::abs(itkEle->eta()) < etaBinsForIsolation_.at(etabin) and
0149           std::abs(itkEle->eta()) > etaBinsForIsolation_.at(etabin - 1) and
0150           itkEle->trkIsol() < trkIsolation_.at(etabin - 1))
0151         passIsolation = true;
0152     }
0153 
0154     if (offlinePt >= min_Pt_ && std::abs(itkEle->eta()) < max_AbsEta1_ && std::abs(itkEle->eta()) >= min_AbsEta1_ &&
0155         passQuality && passIsolation) {
0156       ntrkEle++;
0157       l1t::TkElectronRef ref1(l1t::TkElectronRef(tkEles1, distance(atrkEles, itkEle)));
0158       filterproduct.addObject(trigger::TriggerObjectType::TriggerL1TkEle, ref1);
0159     }
0160   }
0161 
0162   // Loop over second collection. Notice we don't reset ntrkEle
0163   atrkEles = tkEles2->begin();
0164   otrkEles = tkEles2->end();
0165   for (itkEle = atrkEles; itkEle != otrkEles; itkEle++) {
0166     double offlinePt = this->TkEleOfflineEt(itkEle->pt(), itkEle->eta());
0167     bool passQuality(false);
0168     bool passIsolation(false);
0169 
0170     if (applyQual2_) {
0171       if (qual2IsMask_)
0172         passQuality = (itkEle->hwQual() & quality2_);
0173       else
0174         passQuality = (itkEle->hwQual() == quality2_);
0175     } else
0176       passQuality = true;
0177 
0178     for (unsigned int etabin = 1; etabin != etaBinsForIsolation_.size(); ++etabin) {
0179       if (std::abs(itkEle->eta()) < etaBinsForIsolation_.at(etabin) and
0180           std::abs(itkEle->eta()) > etaBinsForIsolation_.at(etabin - 1) and
0181           itkEle->trkIsol() < trkIsolation_.at(etabin - 1))
0182         passIsolation = true;
0183     }
0184 
0185     if (offlinePt >= min_Pt_ && std::abs(itkEle->eta()) <= max_AbsEta2_ && std::abs(itkEle->eta()) >= min_AbsEta2_ &&
0186         passQuality && passIsolation) {
0187       ntrkEle++;
0188       l1t::TkElectronRef ref2(l1t::TkElectronRef(tkEles2, distance(atrkEles, itkEle)));
0189       filterproduct.addObject(trigger::TriggerObjectType::TriggerL1TkEle, ref2);
0190     }
0191   }
0192 
0193   // return with final filter decision
0194   const bool accept(ntrkEle >= min_N_);
0195   return accept;
0196 }
0197 
0198 double L1TTkEleFilter::TkEleOfflineEt(double Et, double Eta) const {
0199   if (std::abs(Eta) < 1.5)
0200     return (barrelScalings_.at(0) + Et * barrelScalings_.at(1) + Et * Et * barrelScalings_.at(2));
0201   else
0202     return (endcapScalings_.at(0) + Et * endcapScalings_.at(1) + Et * Et * endcapScalings_.at(2));
0203 }