Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:28

0001 /** \class HLTSinglet
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Martin Grunewald
0007  *
0008  */
0009 
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 
0013 #include "DataFormats/Common/interface/Handle.h"
0014 #include "DataFormats/Common/interface/Ref.h"
0015 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"
0016 #include "DataFormats/L1Trigger/interface/L1EtMissParticle.h"
0017 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
0018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0019 
0020 #include "HLTSinglet.h"
0021 
0022 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0023 
0024 // extract the candidate type
0025 template <typename T>
0026 int getObjectType(const T&) {
0027   return 0;
0028 }
0029 
0030 // specialize for type l1extra::L1EmParticle
0031 template <typename T>
0032 int getObjectType(const l1extra::L1EmParticle& candidate) {
0033   switch (candidate.type()) {
0034     case l1extra::L1EmParticle::kIsolated:
0035       return trigger::TriggerL1IsoEG;
0036     case l1extra::L1EmParticle::kNonIsolated:
0037       return trigger::TriggerL1NoIsoEG;
0038     default:
0039       return 0;
0040   }
0041 }
0042 
0043 // specialize for type l1extra::L1EtMissParticle
0044 template <typename T>
0045 int getObjectType(const l1extra::L1EtMissParticle& candidate) {
0046   switch (candidate.type()) {
0047     case l1extra::L1EtMissParticle::kMET:
0048       return trigger::TriggerL1ETM;
0049     case l1extra::L1EtMissParticle::kMHT:
0050       return trigger::TriggerL1HTM;
0051     default:
0052       return 0;
0053   }
0054 }
0055 
0056 // specialize for type l1extra::L1JetParticle
0057 template <typename T>
0058 int getObjectType(const l1extra::L1JetParticle& candidate) {
0059   switch (candidate.type()) {
0060     case l1extra::L1JetParticle::kCentral:
0061       return trigger::TriggerL1CenJet;
0062     case l1extra::L1JetParticle::kForward:
0063       return trigger::TriggerL1ForJet;
0064     case l1extra::L1JetParticle::kTau:
0065       return trigger::TriggerL1TauJet;
0066     default:
0067       return 0;
0068   }
0069 }
0070 
0071 //
0072 // constructors and destructor
0073 //
0074 template <typename T>
0075 HLTSinglet<T>::HLTSinglet(const edm::ParameterSet& iConfig)
0076     : HLTFilter(iConfig),
0077       inputTag_(iConfig.template getParameter<edm::InputTag>("inputTag")),
0078       inputToken_(consumes<std::vector<T>>(inputTag_)),
0079       triggerType_(iConfig.template getParameter<int>("triggerType")),
0080       min_N_(iConfig.template getParameter<int>("MinN")),
0081       min_E_(iConfig.template getParameter<double>("MinE")),
0082       min_Pt_(iConfig.template getParameter<double>("MinPt")),
0083       min_Mass_(iConfig.template getParameter<double>("MinMass")),
0084       max_Mass_(iConfig.template getParameter<double>("MaxMass")),
0085       min_Eta_(iConfig.template getParameter<double>("MinEta")),
0086       max_Eta_(iConfig.template getParameter<double>("MaxEta")) {
0087   LogDebug("") << "Input/ptcut/etacut/ncut : " << inputTag_.encode() << " " << min_E_ << " " << min_Pt_ << " "
0088                << min_Mass_ << " " << max_Mass_ << " " << min_Eta_ << " " << max_Eta_ << " " << min_N_;
0089 }
0090 
0091 template <typename T>
0092 HLTSinglet<T>::~HLTSinglet() = default;
0093 
0094 template <typename T>
0095 void HLTSinglet<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0096   edm::ParameterSetDescription desc;
0097   makeHLTFilterDescription(desc);
0098   desc.add<edm::InputTag>("inputTag", edm::InputTag("hltCollection"));
0099   desc.add<int>("triggerType", 0);
0100   desc.add<double>("MinE", -1.0);
0101   desc.add<double>("MinPt", -1.0);
0102   desc.add<double>("MinMass", -1.0);
0103   desc.add<double>("MaxMass", -1.0);
0104   desc.add<double>("MinEta", -1.0);
0105   desc.add<double>("MaxEta", -1.0);
0106   desc.add<int>("MinN", 1);
0107   descriptions.add(defaultModuleLabel<HLTSinglet<T>>(), desc);
0108 }
0109 
0110 //
0111 // member functions
0112 //
0113 
0114 // ------------ method called to produce the data  ------------
0115 template <typename T>
0116 bool HLTSinglet<T>::hltFilter(edm::Event& iEvent,
0117                               const edm::EventSetup& iSetup,
0118                               trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0119   using namespace std;
0120   using namespace edm;
0121   using namespace reco;
0122   using namespace trigger;
0123 
0124   typedef vector<T> TCollection;
0125   typedef Ref<TCollection> TRef;
0126 
0127   // All HLT filters must create and fill an HLT filter object,
0128   // recording any reconstructed physics objects satisfying (or not)
0129   // this HLT filter, and place it in the Event.
0130 
0131   // The filter object
0132   if (saveTags())
0133     filterproduct.addCollectionTag(inputTag_);
0134 
0135   // Ref to Candidate object to be recorded in filter object
0136   TRef ref;
0137 
0138   // get hold of collection of objects
0139   Handle<TCollection> objects;
0140   iEvent.getByToken(inputToken_, objects);
0141 
0142   // look at all objects, check cuts and add to filter object
0143   int n(0);
0144   typename TCollection::const_iterator i(objects->begin());
0145   for (; i != objects->end(); i++) {
0146     if ((i->energy() >= min_E_) && (i->pt() >= min_Pt_) && (i->mass() >= min_Mass_) &&
0147         ((max_Mass_ < 0.0) || (i->mass() <= max_Mass_)) && ((min_Eta_ < 0.0) || (std::abs(i->eta()) >= min_Eta_)) &&
0148         ((max_Eta_ < 0.0) || (std::abs(i->eta()) <= max_Eta_))) {
0149       n++;
0150       ref = TRef(objects, distance(objects->begin(), i));
0151       int tid = getObjectType<T>(*i);
0152       if (tid == 0)
0153         tid = triggerType_;
0154       filterproduct.addObject(tid, ref);
0155     }
0156   }
0157 
0158   // filter decision
0159   bool accept(n >= min_N_);
0160 
0161   return accept;
0162 }