CMS 3D CMS Logo

HLTEgammaDoubleEtFilter.cc
Go to the documentation of this file.
1 
9 
11 
13 
17 
19 
21 public:
23  return lhs->et() > rhs->et();
24  }
25 };
26 
27 //
28 // constructors and destructor
29 //
31 {
32  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
33  etcut1_ = iConfig.getParameter<double> ("etcut1");
34  etcut2_ = iConfig.getParameter<double> ("etcut2");
35  npaircut_ = iConfig.getParameter<int> ("npaircut");
36  l1EGTag_ = iConfig.getParameter< edm::InputTag > ("l1EGCand");
37  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs> (candTag_);
38 }
39 
41 
42 void
46  desc.add<edm::InputTag>("candTag", edm::InputTag("hltTrackIsolFilter"));
47  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
48  desc.add<double>("etcut1", 30.0);
49  desc.add<double>("etcut2", 20.0);
50  desc.add<int>("npaircut", 1);
51  descriptions.add("hltEgammaDoubleEtFilter", desc);
52 }
53 
54 // ------------ method called to produce the data ------------
55 bool
57 {
58  using namespace trigger;
59 
60  // The filter object
61  if (saveTags()) {
62  filterproduct.addCollectionTag(l1EGTag_);
63  }
64  // Ref to Candidate object to be recorded in filter object
66  iEvent.getByToken (candToken_, PrevFilterOutput);
67 
68  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > mysortedrecoecalcands;
69  PrevFilterOutput->getObjects(TriggerPhoton, mysortedrecoecalcands);
70  if(mysortedrecoecalcands.empty()) PrevFilterOutput->getObjects(TriggerCluster,mysortedrecoecalcands); //we dont know if its type trigger cluster or trigger photon
71 
72  // look at all candidates, check cuts and add to filter object
73  int n(0);
74 
75  // Sort the list
76  std::sort(mysortedrecoecalcands.begin(), mysortedrecoecalcands.end(), EgammaHLTEtSortCriterium());
78  for (unsigned int i=0; i<mysortedrecoecalcands.size(); i++) {
79  ref1 = mysortedrecoecalcands[i];
80  if( ref1->et() >= etcut1_){
81  for (unsigned int j=i+1; j<mysortedrecoecalcands.size(); j++) {
82  ref2 = mysortedrecoecalcands[j];
83  if( ref2->et() >= etcut2_ ){
84  filterproduct.addObject(TriggerPhoton, ref1);
85  filterproduct.addObject(TriggerPhoton, ref2);
86  n++;
87  }
88  }
89  }
90  }
91 
92 
93  // filter decision
94  bool accept(n>=npaircut_);
95 
96  return accept;
97 }
98 
99 
100 
101 
102 // declare this class as a framework plugin
T getParameter(std::string const &) const
bool operator()(edm::Ref< reco::RecoEcalCandidateCollection > lhs, edm::Ref< reco::RecoEcalCandidateCollection > rhs)
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
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>)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
~HLTEgammaDoubleEtFilter() override
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45
HLTEgammaDoubleEtFilter(const edm::ParameterSet &)