CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTEgammaAllCombMassFilter.cc
Go to the documentation of this file.
1 
10 
12 
16 
19 
20 //
21 // constructors and destructor
22 //
24 {
25  firstLegLastFilterTag_ = iConfig.getParameter<edm::InputTag>("firstLegLastFilter");
26  secondLegLastFilterTag_= iConfig.getParameter<edm::InputTag>("secondLegLastFilter");
27  minMass_ = iConfig.getParameter<double> ("minMass");
28  firstLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(firstLegLastFilterTag_);
29  secondLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(secondLegLastFilterTag_);
30 }
31 
33 
34 void
38  desc.add<edm::InputTag>("firstLegLastFilter",edm::InputTag("firstFilter"));
39  desc.add<edm::InputTag>("secondLegLastFilter",edm::InputTag("secondFilter"));
40  desc.add<double>("minMass",-1.0);
41  descriptions.add("hltEgammaAllCombMassFilter",desc);
42 }
43 
44 
45 // ------------ method called to produce the data ------------
46 
48 {
49  //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)
50 
51  //trigger::TriggerObjectType firstLegTrigType;
52  std::vector<math::XYZTLorentzVector> firstLegP4s;
53 
54  //trigger::TriggerObjectType secondLegTrigType;
55  std::vector<math::XYZTLorentzVector> secondLegP4s;
56 
58 
59  getP4OfLegCands(iEvent,firstLegLastFilterToken_,firstLegP4s);
60  getP4OfLegCands(iEvent,secondLegLastFilterToken_,secondLegP4s);
61 
62  bool accept=false;
63  for(size_t firstLegNr=0;firstLegNr<firstLegP4s.size();firstLegNr++){
64  for(size_t secondLegNr=0;secondLegNr<secondLegP4s.size();secondLegNr++){
65  math::XYZTLorentzVector pairP4 = firstLegP4s[firstLegNr] + secondLegP4s[secondLegNr];
66  double mass = pairP4.M();
67  if(mass>=minMass_) accept=true;
68  }
69  }
70  for(size_t firstLegNr=0;firstLegNr<firstLegP4s.size();firstLegNr++){
71  for(size_t secondLegNr=0;secondLegNr<firstLegP4s.size();secondLegNr++){
72  math::XYZTLorentzVector pairP4 = firstLegP4s[firstLegNr] + firstLegP4s[secondLegNr];
73  double mass = pairP4.M();
74  if(mass>=minMass_) accept=true;
75  }
76  }
77  for(size_t firstLegNr=0;firstLegNr<secondLegP4s.size();firstLegNr++){
78  for(size_t secondLegNr=0;secondLegNr<secondLegP4s.size();secondLegNr++){
79  math::XYZTLorentzVector pairP4 = secondLegP4s[firstLegNr] + secondLegP4s[secondLegNr];
80  double mass = pairP4.M();
81  if(mass>=minMass_) accept=true;
82  }
83  }
84 
85  return accept;
86 }
87 
88 void HLTEgammaAllCombMassFilter::getP4OfLegCands(const edm::Event& iEvent, const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs>& filterToken, std::vector<math::XYZTLorentzVector>& p4s)
89 {
91  iEvent.getByToken(filterToken,filterOutput);
92 
93  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
94  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > phoCands;
95  filterOutput->getObjects(trigger::TriggerPhoton,phoCands);
96  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
97  filterOutput->getObjects(trigger::TriggerCluster,clusCands);
98  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
99  filterOutput->getObjects(trigger::TriggerElectron,eleCands);
100 
101  if(!phoCands.empty()){ //its photons
102  for(size_t candNr=0;candNr<phoCands.size();candNr++){
103  p4s.push_back(phoCands[candNr]->p4());
104  }
105  }else if(!clusCands.empty()){ //try trigger cluster (should never be this, at the time of writing (17/1/11) this would indicate an error)
106  for(size_t candNr=0;candNr<clusCands.size();candNr++){
107  p4s.push_back(clusCands[candNr]->p4());
108  }
109  }else if(!eleCands.empty()){
110  for(size_t candNr=0;candNr<eleCands.size();candNr++){
111  p4s.push_back(eleCands[candNr]->p4());
112  }
113  }
114 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > firstLegLastFilterToken_
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > secondLegLastFilterToken_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
HLTEgammaAllCombMassFilter(const edm::ParameterSet &)
double p4[4]
Definition: TauolaWrapper.h:92
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void getP4OfLegCands(const edm::Event &iEvent, const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > &filterToken, std::vector< math::XYZTLorentzVector > &p4s)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)