CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  relaxed_ = iConfig.getUntrackedParameter<bool> ("relaxed",true) ;
37  L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand");
38  L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand");
39  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
40 }
41 
43 
44 void
48  desc.add<edm::InputTag>("candTag",edm::InputTag("hltTrackIsolFilter"));
49  desc.add<edm::InputTag>("L1IsoCand",edm::InputTag("hltL1IsoRecoEcalCandidate"));
50  desc.add<edm::InputTag>("L1NonIsoCand",edm::InputTag("hltL1NonIsoRecoEcalCandidate"));
51  desc.addUntracked<bool>("relaxed",true);
52  desc.add<double>("etcut1", 30.0);
53  desc.add<double>("etcut2", 20.0);
54  desc.add<int>("npaircut", 1);
55  descriptions.add("hltEgammaDoubleEtFilter",desc);
56 }
57 
58 // ------------ method called to produce the data ------------
59 bool
61 {
62  using namespace trigger;
63 
64  // The filter object
65  if (saveTags()) {
66  filterproduct.addCollectionTag(L1IsoCollTag_);
67  if (relaxed_) filterproduct.addCollectionTag(L1NonIsoCollTag_);
68  }
69  // Ref to Candidate object to be recorded in filter object
71  iEvent.getByToken (candToken_,PrevFilterOutput);
72 
73  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > mysortedrecoecalcands;
74  PrevFilterOutput->getObjects(TriggerPhoton, mysortedrecoecalcands);
75  if(mysortedrecoecalcands.empty()) PrevFilterOutput->getObjects(TriggerCluster,mysortedrecoecalcands); //we dont know if its type trigger cluster or trigger photon
76 
77  // look at all candidates, check cuts and add to filter object
78  int n(0);
79 
80  // Sort the list
81  std::sort(mysortedrecoecalcands.begin(), mysortedrecoecalcands.end(), EgammaHLTEtSortCriterium());
83  for (unsigned int i=0; i<mysortedrecoecalcands.size(); i++) {
84  ref1 = mysortedrecoecalcands[i];
85  if( ref1->et() >= etcut1_){
86  for (unsigned int j=i+1; j<mysortedrecoecalcands.size(); j++) {
87  ref2 = mysortedrecoecalcands[j];
88  if( ref2->et() >= etcut2_ ){
89  filterproduct.addObject(TriggerPhoton, ref1);
90  filterproduct.addObject(TriggerPhoton, ref2);
91  n++;
92  }
93  }
94  }
95  }
96 
97 
98  // filter decision
99  bool accept(n>=npaircut_);
100 
101  return accept;
102 }
103 
104 
105 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool operator()(edm::Ref< reco::RecoEcalCandidateCollection > lhs, edm::Ref< reco::RecoEcalCandidateCollection > rhs)
int i
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const 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 &)