CMS 3D CMS Logo

HLTElectronPFMTFilter.cc
Go to the documentation of this file.
3 
4 template <typename T>
6 {
7  // MHT parameters
8  inputMetTag_ = iConfig.getParameter< edm::InputTag > ("inputMetTag");
9  minMht_ = iConfig.getParameter<double> ("minMht");
10 
11  // Electron parameters
12  inputEleTag_ = iConfig.getParameter< edm::InputTag > ("inputEleTag");
13  lowerMTCut_ = iConfig.getParameter<double> ("lowerMTCut");
14  upperMTCut_ = iConfig.getParameter<double> ("upperMTCut");
15  minN_ = iConfig.getParameter<int>("minN");
16  l1EGTag_ = iConfig.getParameter< edm::InputTag > ("l1EGCand");
17 
18  inputMetToken_ = consumes<reco::METCollection> (inputMetTag_);
19  inputEleToken_ = consumes<trigger::TriggerFilterObjectWithRefs> (inputEleTag_);
20 }
21 
22 template <typename T>
24 
25 template <typename T>
29  desc.add<edm::InputTag>("inputMetTag", edm::InputTag("hltPFMHT"));
30  desc.add<edm::InputTag>("inputEleTag", edm::InputTag("hltEle25CaloIdVTTrkIdTCaloIsoTTrkIsoTTrackIsolFilter"));
31  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
32  desc.add<int>("minN", 0);
33  desc.add<double>("minMht", 0.0);
34  desc.add<double>("lowerMTCut", 0.0);
35  desc.add<double>("upperMTCut", 9999.0);
36  descriptions.add(defaultModuleLabel<HLTElectronPFMTFilter<T>>(), desc);
37 }
38 
39 template <typename T>
41 {
42  using namespace std;
43  using namespace edm;
44  using namespace reco;
45  using namespace trigger;
46  // The filter object
47  if (saveTags()) {
48  filterproduct.addCollectionTag(inputMetTag_);
49  filterproduct.addCollectionTag(l1EGTag_);
50  }
51 
52  // Get the Met collection
54  iEvent.getByToken(inputMetToken_,pfMHT);
55 
56  // Sanity check:
57  if(!pfMHT.isValid()) {
58  edm::LogError("HLTElectronPFMTFilter") << "missing input Met collection!";
59  }
60 
61  const METCollection *metcol = pfMHT.product();
62  const MET *met;
63  met = &(metcol->front());
64 
66  iEvent.getByToken (inputEleToken_,PrevFilterOutput);
67 
68  int nW = 0;
69 
70  vector< Ref< vector<T> > > refEleCollection ;
71 
72  PrevFilterOutput->getObjects(TriggerElectron, refEleCollection);
73  int trigger_type = trigger::TriggerElectron;
74  if(refEleCollection.empty()){
75  PrevFilterOutput->getObjects(TriggerCluster, refEleCollection);
76  trigger_type = trigger::TriggerCluster;
77  if(refEleCollection.empty()){
78  PrevFilterOutput->getObjects(TriggerPhoton, refEleCollection);
79  trigger_type = trigger::TriggerPhoton;
80  }
81  }
82 
83 
84  TLorentzVector pMET(met->px(), met->py(),0.0,sqrt(met->px()*met->px() + met->py()*met->py()));
85 
86  for (unsigned int i=0; i<refEleCollection.size(); i++) {
87  TLorentzVector pThisEle(refEleCollection.at(i)->px(), refEleCollection.at(i)->py(),
88  0.0, refEleCollection.at(i)->et() );
89  TLorentzVector pTot = pMET + pThisEle;
90  double mass = pTot.M();
91 
92  if(mass>=lowerMTCut_ && mass<=upperMTCut_ && pMET.E()>= minMht_){
93  nW++;
94  filterproduct.addObject(trigger_type, refEleCollection.at(i));
95  }
96  }
97 
98  // filter decision
99  const bool accept(nW>=minN_);
100  return accept;
101 
102 }
103 
104 
105 // declare the specialisations as framework plugins
107 
110 
T getParameter(std::string const &) const
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputEleToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
double px() const final
x coordinate of momentum vector
std::string defaultModuleLabel()
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: MET.h:42
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< reco::METCollection > inputMetToken_
Collection of MET.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
HLTElectronPFMTFilter(const edm::ParameterSet &)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
~HLTElectronPFMTFilter() override
T const * product() const
Definition: Handle.h:74
met
===> hadronic RAZOR
double py() const final
y coordinate of momentum vector
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.