CMS 3D CMS Logo

HLTEgammaEtFilterPairs.cc
Go to the documentation of this file.
1 
9 
11 
15 
17 
18 //
19 // constructors and destructor
20 //
22 {
23  inputTag_ = iConfig.getParameter< edm::InputTag > ("inputTag");
24  etcutEB1_ = iConfig.getParameter<double> ("etcut1EB");
25  etcutEE1_ = iConfig.getParameter<double> ("etcut1EE");
26  etcutEB2_ = iConfig.getParameter<double> ("etcut2EB");
27  etcutEE2_ = iConfig.getParameter<double> ("etcut2EE");
28  l1EGTag_ = iConfig.getParameter< edm::InputTag > ("l1EGCand");
29  inputToken_ = consumes<trigger::TriggerFilterObjectWithRefs> (inputTag_);
30 }
31 
32 void
36  desc.add<edm::InputTag>("inputTag", edm::InputTag("HLTEgammaL1MatchFilter"));
37  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
38  desc.add<double>("etcut1EB", 1.0);
39  desc.add<double>("etcut1EE", 1.0);
40  desc.add<double>("etcut2EB", 1.0);
41  desc.add<double>("etcut2EE", 1.0);
42  descriptions.add("hltEgammaEtFilterPairs", desc);
43 }
44 
46 
47 
48 // ------------ method called to produce the data ------------
49 bool
51 {
52  using namespace trigger;
53  // The filter object
54  if (saveTags()) {
55  filterproduct.addCollectionTag(l1EGTag_);
56  }
57 
59  iEvent.getByToken (inputToken_, PrevFilterOutput);
60 
61  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands; // vref with your specific C++ collection type
62  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
63  if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton, recoecalcands);
64  // they list should be interpreted as pairs:
65  // <recoecalcands[0],recoecalcands[1]>
66  // <recoecalcands[2],recoecalcands[3]>
67  // <recoecalcands[4],recoecalcands[5]>
68  // .......
69 
70  // Should I check that the size of recoecalcands is even ?
71  int n(0);
72 
73  for (unsigned int i=0; i<recoecalcands.size(); i=i+2) {
74 
77  // std::cout<<"EtFilter: 1) Et Eta phi: "<<r1->et()<<" "<<r1->eta()<<" "<<r1->phi()<<" 2) Et eta phi: "<<r2->et()<<" "<<r2->eta()<<" "<<r2->phi()<<std::endl;
78  bool first = (fabs(r1->eta()) < 1.479 && r1->et() >= etcutEB1_) || (fabs(r1->eta()) >= 1.479 && r1->et() >= etcutEE1_);
79  bool second = (fabs(r2->eta()) < 1.479 && r2->et() >= etcutEB2_) || (fabs(r2->eta()) >= 1.479 && r2->et() >= etcutEE2_);
80 
81  if ( first && second ) {
82  n++;
83  filterproduct.addObject(TriggerCluster,r1 );
84  filterproduct.addObject(TriggerCluster,r2 );
85  }
86  }
87 
88 
89  // filter decision
90  bool accept(n>=1);
91 
92  return accept;
93 }
94 
95 // declare this class as a framework plugin
T getParameter(std::string const &) const
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
~HLTEgammaEtFilterPairs() override
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken_
bool saveTags() const
Definition: HLTFilter.h:45
HLTEgammaEtFilterPairs(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)