CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTEgammaDoubleEtPhiFilter.cc
Go to the documentation of this file.
1 
9 
11 
15 
17 
19 public:
21  return lhs->et() > rhs->et();
22  }
23 };
24 
25 //
26 // constructors and destructor
27 //
29 {
30  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
31  etcut1_ = iConfig.getParameter<double> ("etcut1");
32  etcut2_ = iConfig.getParameter<double> ("etcut2");
33  min_Acop_ = iConfig.getParameter<double> ("MinAcop");
34  max_Acop_ = iConfig.getParameter<double> ("MaxAcop");
35  min_EtBalance_ = iConfig.getParameter<double> ("MinEtBalance");
36  max_EtBalance_ = iConfig.getParameter<double> ("MaxEtBalance");
37  npaircut_ = iConfig.getParameter<int> ("npaircut");
38  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
39 }
40 
42 
43 void
47  desc.add<edm::InputTag>("candTag",edm::InputTag("hltDoubleL1MatchFilter"));
48  desc.add<double>("etcut1", 6.0);
49  desc.add<double>("etcut2", 6.0);
50  desc.add<double>("MinAcop", -0.1);
51  desc.add<double>("MaxAcop", 0.6);
52  desc.add<double>("MinEtBalance", -1.0);
53  desc.add<double>("MaxEtBalance", 10.0);
54  desc.add<int>("npaircut", 1);
55  descriptions.add("hltEgammaDoubleEtPhiFilter",desc);
56 }
57 
58 // ------------ method called to produce the data ------------
59 bool
61 {
62  using namespace trigger;
63 
64  // Ref to Candidate object to be recorded in filter object
66 
68  iEvent.getByToken (candToken_,PrevFilterOutput);
69 
70  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > mysortedrecoecalcands;
71  PrevFilterOutput->getObjects(TriggerCluster, mysortedrecoecalcands);
72  if(mysortedrecoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton,mysortedrecoecalcands); //we dont know if its type trigger cluster or trigger photon
73  // Sort the list
74  std::sort(mysortedrecoecalcands.begin(), mysortedrecoecalcands.end(), EgammaHLTEtSortCriterium());
76 
77  int n(0);
78  for (unsigned int i=0; i<mysortedrecoecalcands.size(); i++) {
79  ref1 = mysortedrecoecalcands[i];
80  if( ref1->et() >= etcut1_){
81 
82  for (unsigned int j=i+1; j<mysortedrecoecalcands.size(); j++) {
83  ref2 = mysortedrecoecalcands[j];
84  if( ref2->et() >= etcut2_ ){
85 
86  // Acoplanarity
87  double acop = fabs(ref1->phi()-ref2->phi());
88  if (acop>M_PI) acop = 2*M_PI - acop;
89  acop = M_PI - acop;
90  LogDebug("HLTEgammaDoubleEtPhiFilter") << " ... 1-2 acop= " << acop;
91 
92  if ((acop>=min_Acop_) && (acop<=max_Acop_))
93  {
94  // Et balance
95  double etbalance = fabs(ref1->et()-ref2->et());
96  if ((etbalance>=min_EtBalance_) && (etbalance<=max_EtBalance_))
97  {
98  filterproduct.addObject(TriggerCluster, ref1);
99  filterproduct.addObject(TriggerCluster, ref2);
100  n++;
101  }
102  }
103  }
104  }
105  }
106  }
107 
108 
109  // filter decision
110  bool accept(n>=npaircut_);
111 
112  return accept;
113 }
114 
115 
116 
#define LogDebug(id)
T getParameter(std::string const &) const
bool operator()(edm::Ref< reco::RecoEcalCandidateCollection > lhs, edm::Ref< reco::RecoEcalCandidateCollection > rhs)
int i
Definition: DBlmapReader.cc:9
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;)
int iEvent
Definition: GenABIO.cc:230
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define M_PI
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
HLTEgammaDoubleEtPhiFilter(const edm::ParameterSet &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)