CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTEgammaCombMassFilter.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("hltEgammaCombMassFilter",desc);
42 }
43 
44 // ------------ method called to produce the data ------------
46 {
47  //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)
48 
49  //trigger::TriggerObjectType firstLegTrigType;
50  std::vector<math::XYZTLorentzVector> firstLegP4s;
51 
52  //trigger::TriggerObjectType secondLegTrigType;
53  std::vector<math::XYZTLorentzVector> secondLegP4s;
54 
56 
57  getP4OfLegCands(iEvent,firstLegLastFilterToken_,firstLegP4s);
58  getP4OfLegCands(iEvent,secondLegLastFilterToken_,secondLegP4s);
59 
60  bool accept=false;
61  for(size_t firstLegNr=0;firstLegNr<firstLegP4s.size();firstLegNr++){
62  for(size_t secondLegNr=0;secondLegNr<secondLegP4s.size();secondLegNr++){
63  math::XYZTLorentzVector pairP4 = firstLegP4s[firstLegNr] + secondLegP4s[secondLegNr];
64  double mass = pairP4.M();
65  if(mass>=minMass_) accept=true;
66  }
67  }
68 
69  return accept;
70 }
71 
72 void HLTEgammaCombMassFilter::getP4OfLegCands(const edm::Event& iEvent, const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs>& filterToken, std::vector<math::XYZTLorentzVector>& p4s)
73 {
75  iEvent.getByToken(filterToken,filterOutput);
76 
77  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
78  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > phoCands;
79  filterOutput->getObjects(trigger::TriggerPhoton,phoCands);
80  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
81  filterOutput->getObjects(trigger::TriggerCluster,clusCands);
82  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
83  filterOutput->getObjects(trigger::TriggerElectron,eleCands);
84 
85  if(!phoCands.empty()){ //its photons
86  for(size_t candNr=0;candNr<phoCands.size();candNr++){
87  p4s.push_back(phoCands[candNr]->p4());
88  }
89  }else if(!clusCands.empty()){ //try trigger cluster (should never be this, at the time of writing (17/1/11) this would indicate an error)
90  for(size_t candNr=0;candNr<clusCands.size();candNr++){
91  p4s.push_back(clusCands[candNr]->p4());
92  }
93  }else if(!eleCands.empty()){
94  for(size_t candNr=0;candNr<eleCands.size();candNr++){
95  p4s.push_back(eleCands[candNr]->p4());
96  }
97  }
98 }
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
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
HLTEgammaCombMassFilter(const edm::ParameterSet &)
double p4[4]
Definition: TauolaWrapper.h:92
static void getP4OfLegCands(const edm::Event &iEvent, const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > &filterToken, std::vector< math::XYZTLorentzVector > &p4s)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > firstLegLastFilterToken_
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)