CMS 3D CMS Logo

HLTEgammaDoubleLegCombFilter.cc
Go to the documentation of this file.
1 
10 
12 
16 
20 
22 
23 //
24 // constructors and destructor
25 //
27 {
28  firstLegLastFilterTag_ = iConfig.getParameter<edm::InputTag>("firstLegLastFilter");
29  secondLegLastFilterTag_= iConfig.getParameter<edm::InputTag>("secondLegLastFilter");
30  nrRequiredFirstLeg_ = iConfig.getParameter<int> ("nrRequiredFirstLeg");
31  nrRequiredSecondLeg_ = iConfig.getParameter<int> ("nrRequiredSecondLeg");
32  nrRequiredUniqueLeg_ = iConfig.getParameter<int> ("nrRequiredUniqueLeg");
33  maxMatchDR_ = iConfig.getParameter<double> ("maxMatchDR");
34  firstLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(firstLegLastFilterTag_);
35  secondLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(secondLegLastFilterTag_);
36 }
37 
38 void
42  desc.add<edm::InputTag>("firstLegLastFilter",edm::InputTag("firstFilter"));
43  desc.add<edm::InputTag>("secondLegLastFilter",edm::InputTag("secondFilter"));
44  desc.add<int>("nrRequiredFirstLeg",0);
45  desc.add<int>("nrRequiredSecondLeg",0);
46  desc.add<int>("nrRequiredUniqueLeg",0);
47  desc.add<double>("maxMatchDR",0.01);
48  descriptions.add("hltEgammaDoubleLegCombFilter",desc);
49 }
50 
52 
53 
54 // ------------ method called to produce the data ------------
56 {
57  //right, issue 1, we dont know if this is a TriggerElectron, TriggerPhoton, TriggerCluster (should never be a TriggerCluster btw as that implies the 4-vectors are not stored in AOD)
58 
59  //trigger::TriggerObjectType firstLegTrigType;
60  std::vector<math::XYZPoint> firstLegP3s;
61 
62  //trigger::TriggerObjectType secondLegTrigType;
63  std::vector<math::XYZPoint> secondLegP3s;
64 
65  getP3OfLegCands(iEvent,firstLegLastFilterToken_,firstLegP3s);
66  getP3OfLegCands(iEvent,secondLegLastFilterToken_,secondLegP3s);
67 
68  std::vector<std::pair<int,int> > matchedCands;
69  matchCands(firstLegP3s,secondLegP3s,matchedCands);
70 
71 
72  int nr1stLegOnly=0;
73  int nr2ndLegOnly=0;
74  int nrBoth=0;;
75 
76  for(auto & matchedCand : matchedCands){
77  if(matchedCand.first>=0){ //we found a first leg cand
78  if(matchedCand.second>=0) nrBoth++;//we also found a second leg cand
79  else nr1stLegOnly++; //we didnt find a second leg cand
80  }else if(matchedCand.second>=0) nr2ndLegOnly++; //we found a second leg cand but we didnt find a first leg
81 
82  }
83 
84  bool accept=false;
85  if(nr1stLegOnly + nr2ndLegOnly + nrBoth >= nrRequiredUniqueLeg_) {
86  if(nr1stLegOnly>=nrRequiredFirstLeg_ && nr2ndLegOnly>=nrRequiredSecondLeg_) accept=true;
87  else{
88  int nrNeededFirstLeg = std::max(0,nrRequiredFirstLeg_ - nr1stLegOnly);
89  int nrNeededSecondLeg = std::max(0,nrRequiredSecondLeg_ - nr2ndLegOnly);
90 
91  if(nrBoth >= nrNeededFirstLeg + nrNeededSecondLeg) accept = true;
92  }
93  }
94 
95  return accept;
96 }
97 
98 //-1 if no match is found
99 void HLTEgammaDoubleLegCombFilter::matchCands(const std::vector<math::XYZPoint>& firstLegP3s,const std::vector<math::XYZPoint>& secondLegP3s,std::vector<std::pair<int,int> >&matchedCands) const
100 {
101  std::vector<size_t> matched2ndLegs;
102  for(size_t firstLegNr=0;firstLegNr<firstLegP3s.size();firstLegNr++){
103  int matchedNr = -1;
104  for(size_t secondLegNr=0;secondLegNr<secondLegP3s.size() && matchedNr==-1;secondLegNr++){
105  if(reco::deltaR2(firstLegP3s[firstLegNr],secondLegP3s[secondLegNr])<maxMatchDR_*maxMatchDR_) matchedNr=secondLegNr;
106  }
107  matchedCands.push_back(std::make_pair(firstLegNr,matchedNr));
108  if(matchedNr>=0) matched2ndLegs.push_back(static_cast<size_t>(matchedNr));
109  }
110  std::sort(matched2ndLegs.begin(),matched2ndLegs.end());
111 
112  for(size_t secondLegNr=0;secondLegNr<secondLegP3s.size();secondLegNr++){
113  if(!std::binary_search(matched2ndLegs.begin(),matched2ndLegs.end(),secondLegNr)){ //wasnt matched already
114  matchedCands.push_back(std::make_pair<int,int>(-1, static_cast<int>(secondLegNr)));
115  }
116  }
117 }
118 
119 //we use position and p3 interchangably here, we only use eta/phi so its alright
121 {
123  iEvent.getByToken(filterToken,filterOutput);
124 
125  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
126  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > phoCands;
127  filterOutput->getObjects(trigger::TriggerPhoton,phoCands);
128  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
129  filterOutput->getObjects(trigger::TriggerCluster,clusCands);
130  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
131  filterOutput->getObjects(trigger::TriggerElectron,eleCands);
132  std::vector<edm::Ref<reco::CaloJetCollection> > jetCands;
133  filterOutput->getObjects(trigger::TriggerJet,jetCands);
134 
135  if(!phoCands.empty()){ //its photons
136  for(auto & phoCand : phoCands){
137  p3s.push_back(phoCand->superCluster()->position());
138  }
139  }else if(!clusCands.empty()){ //try trigger cluster (should never be this, at the time of writing (17/1/11) this would indicate an error)
140  for(auto & clusCand : clusCands){
141  p3s.push_back(clusCand->superCluster()->position());
142  }
143  }else if(!eleCands.empty()){
144  for(auto & eleCand : eleCands){
145  p3s.push_back(eleCand->superCluster()->position());
146  }
147  }else if(!jetCands.empty()){
148  for(auto & jetCand : jetCands){
150  p3.SetX(jetCand->p4().x());
151  p3.SetY(jetCand->p4().y());
152  p3.SetZ(jetCand->p4().z());
153  p3s.push_back(p3);
154  }
155  }
156 }
157 
158 // 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
~HLTEgammaDoubleLegCombFilter() override
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > firstLegLastFilterToken_
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
HLTEgammaDoubleLegCombFilter(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void getP3OfLegCands(const edm::Event &iEvent, const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > &filterToken, std::vector< math::XYZPoint > &p3s)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void matchCands(const std::vector< math::XYZPoint > &firstLegP3s, const std::vector< math::XYZPoint > &secondLegP3s, std::vector< std::pair< int, int > > &matchedCands) const
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > secondLegLastFilterToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double p3[4]
Definition: TauolaWrapper.h:91