CMS 3D CMS Logo

WZInterestingEventSelector.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 //
4 //
5 //
6 
7 
8 // system include files
9 #include <memory>
10 #include <iostream>
11 #include <fstream>
12 #include <sstream>
13 #include <string>
14 #include <cstdio>
15 #include <cstdlib>
16 // user include files
19 
23 
25 
37 
38 #include "TLorentzVector.h"
39 //
40 // class declaration
41 //
42 
43 using namespace reco;
44 
46 public:
47  struct event
48  {
49  long run;
50  long event;
51  long ls;
52  int nEle;
53  float maxPt;
54  float maxPtEleEta;
55  float maxPtElePhi;
56  float invMass;
57  float met;
58  float metPhi;
59  };
60 
62  ~WZInterestingEventSelector() override;
63 
64 private:
65  bool filter(edm::Event&, const edm::EventSetup&) override;
66  void endJob() override;
67  bool electronSelection( const GsfElectron* eleRef , const math::XYZPoint& bspotPosition);
68  // ----------member data ---------------------------
69 
70  //std::vector<event> interestingEvents_;
71 
72  //Pt cut
73  float ptCut_;
75 
76  //EB ID+ISO cuts
77  float eb_trIsoCut_;
80  float eb_hoeCut_;
81  float eb_seeCut_;
82 
83  //EE ID+ISO cuts
84  float ee_trIsoCut_;
87  float ee_hoeCut_;
88  float ee_seeCut_;
89 
90  //met Cut
91  float metCut_;
92 
93  //invMass Cut
94  float invMassCut_;
95 
99 
100 };
101 
102 //
103 // constructors and destructor
104 //
106 {
107  ptCut_ = iConfig.getParameter<double>("ptCut");
108  missHitCut_ = iConfig.getParameter<int>("missHitsCut");
109 
110  eb_trIsoCut_ = iConfig.getParameter<double>("eb_trIsoCut");
111  eb_ecalIsoCut_ = iConfig.getParameter<double>("eb_ecalIsoCut");
112  eb_hcalIsoCut_ = iConfig.getParameter<double>("eb_hcalIsoCut");
113  eb_hoeCut_ = iConfig.getParameter<double>("eb_hoeCut");
114  eb_seeCut_ = iConfig.getParameter<double>("eb_seeCut");
115 
116  ee_trIsoCut_ = iConfig.getParameter<double>("ee_trIsoCut");
117  ee_ecalIsoCut_ = iConfig.getParameter<double>("ee_ecalIsoCut");
118  ee_hcalIsoCut_ = iConfig.getParameter<double>("ee_hcalIsoCut");
119  ee_hoeCut_ = iConfig.getParameter<double>("ee_hoeCut");
120  ee_seeCut_ = iConfig.getParameter<double>("ee_seeCut");
121 
122  metCut_ = iConfig.getParameter<double>("metCut");
123  invMassCut_ = iConfig.getParameter<double>("invMassCut");
124 
125  electronCollection_ = iConfig.getUntrackedParameter<edm::InputTag>("electronCollection",edm::InputTag("gsfElectrons"));
126  pfMetCollection_ = iConfig.getUntrackedParameter<edm::InputTag>("pfMetCollection",edm::InputTag("pfMet"));
127  offlineBSCollection_ = iConfig.getUntrackedParameter<edm::InputTag>("offlineBSCollection",edm::InputTag("offlineBeamSpot"));
128 
129 }
130 
131 
133 {
134 
135  // do anything here that needs to be done at desctruction time
136  // (e.g. close files, deallocate resources etc.)
137 
138 }
139 
140 
141 //
142 // member functions
143 //
144 
145 // ------------ method called on each new Event ------------
146 
148 
149 // if (interestingEvents_.size()<1)
150 // return;
151 
152 // std::ostringstream oss;
153 // for (unsigned int iEvent=0;iEvent<interestingEvents_.size();++iEvent)
154 // {
155 // oss << "==================================" << std::endl;
156 // oss << "Run: " << interestingEvents_[iEvent].run << " Event: " << interestingEvents_[iEvent].event << " LS: " << interestingEvents_[iEvent].ls << std::endl;
157 // oss << "nGoodEle: " << interestingEvents_[iEvent].nEle << " maxPt " << interestingEvents_[iEvent].maxPt << " maxPtEta " << interestingEvents_[iEvent].maxPtEleEta << " maxPtPhi " << interestingEvents_[iEvent].maxPtElePhi << std::endl;
158 // oss << "invMass " << interestingEvents_[iEvent].invMass << " met " << interestingEvents_[iEvent].met << " metPhi " << interestingEvents_[iEvent].metPhi << std::endl;
159 // }
160 // std::string mailText;
161 // mailText = oss.str();
162 
163 // std::ofstream outputTxt;
164 // outputTxt.open("interestingEvents.txt");
165 // outputTxt << mailText;
166 // outputTxt.close();
167 
168  //Sending email
169 // std::ostringstream subject;
170 // subject << "Interesting events in Run#" << interestingEvents_[0].run;
171 
172 // std::ostringstream command;
173 // command << "cat interestingEvents.txt | mail -s \"" << subject.str() << "\" Paolo.Meridiani@cern.ch";
174 
175 // std::string commandStr = command.str();
176 // char* pch = (char*)malloc( sizeof( char ) *(commandStr.length() +1) );
177 // string::traits_type::copy( pch, commandStr.c_str(), commandStr.length() +1 );
178 // int i=system(pch);
179 
180 }
181 
183 {
184 
185 
186 // if (eleRef->trackerDrivenSeed() && !eleRef->ecalDrivenSeed())
187 // return false;
188 
189 // if (eleRef->ecalDrivenSeed())
190 // {
191 
192  if (eleRef->pt()<ptCut_) return false;
193 
194  if (eleRef->isEB())
195  {
196  if (eleRef->dr03TkSumPt()/eleRef->pt()>eb_trIsoCut_) return false;
197  if (eleRef->dr03EcalRecHitSumEt()/eleRef->pt()>eb_ecalIsoCut_) return false;
198  if (eleRef->dr03HcalTowerSumEt()/eleRef->pt()>eb_hcalIsoCut_) return false;
199  if (eleRef->sigmaIetaIeta()>eb_seeCut_) return false;
200  if (eleRef->hcalOverEcal()>eb_hoeCut_) return false;
201  }
202  else if (eleRef->isEE())
203  {
204  if (eleRef->dr03TkSumPt()/eleRef->pt()>ee_trIsoCut_) return false;
205  if (eleRef->dr03EcalRecHitSumEt()/eleRef->pt()>ee_ecalIsoCut_) return false;
206  if (eleRef->dr03HcalTowerSumEt()/eleRef->pt()>ee_hcalIsoCut_) return false;
207  if (eleRef->sigmaIetaIeta()>ee_seeCut_) return false;
208  if (eleRef->hcalOverEcal()>ee_hoeCut_) return false;
209  }
210 
211  if (eleRef->gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > missHitCut_) return false;
212 
213  return true;
214 }
215 
216 bool
218 {
219  // using namespace edm;
221  iEvent.getByLabel(electronCollection_,gsfElectrons);
222 
223 // edm::Handle<reco::CaloMETCollection> caloMET;
224 // iEvent.getByLabel(edm::InputTag("met"), caloMET);
225 
227  iEvent.getByLabel(pfMetCollection_, pfMET);
228 
229  edm::Handle<reco::BeamSpot> pBeamSpot;
230  iEvent.getByLabel(offlineBSCollection_, pBeamSpot);
231 
232  const reco::BeamSpot *bspot = pBeamSpot.product();
233  const math::XYZPoint& bspotPosition = bspot->position();
234 
235  std::vector<const GsfElectron*> goodElectrons;
236  float ptMax=-999.;
237  const GsfElectron* ptMaxEle=nullptr;
238  for(reco::GsfElectronCollection::const_iterator myEle=gsfElectrons->begin();myEle!=gsfElectrons->end();++myEle)
239  {
240  //Apply a minimal isolated electron selection
241  if (!electronSelection(&(*myEle),bspotPosition)) continue;
242  goodElectrons.push_back(&(*myEle));
243  if (myEle->pt() > ptMax )
244  {
245  ptMax = myEle->pt();
246  ptMaxEle = &(*myEle);
247  }
248  }
249 
250  float maxInv=-999.;
251  TLorentzVector v1;
252  if (ptMaxEle)
253  v1.SetPtEtaPhiM(ptMaxEle->pt(),ptMaxEle->eta(),ptMaxEle->phi(),0);
254  if (goodElectrons.size()>1)
255  {
256  for (unsigned int iEle=0; iEle<goodElectrons.size(); ++iEle)
257  if (goodElectrons[iEle]!=ptMaxEle && (goodElectrons[iEle]->charge() * ptMaxEle->charge() == -1) )
258  {
259  TLorentzVector v2;
260  v2.SetPtEtaPhiM(goodElectrons[iEle]->pt(),goodElectrons[iEle]->eta(),goodElectrons[iEle]->phi(),0.);
261  if ( (v1+v2).M() > maxInv )
262  maxInv = (v1+v2).M();
263  }
264  }
265 
266  //Z filt: Retain event if more then 1 good ele and invMass above threshold (zee)
267  if (goodElectrons.size()>1 && maxInv > invMassCut_)
268  {
269  //interestingEvents_.push_back(thisEvent);
270  return true;
271  }
272 
273  //W filt: Retain event also event with at least 1 good ele and some met
274  if (!goodElectrons.empty() && (pfMET->begin()->et()>metCut_))
275  {
276  //interestingEvents_.push_back(thisEvent);
277  return true;
278  }
279 
280  return false;
281 }
282 
283 //define this as a plug-in
T getParameter(std::string const &) const
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:186
T getUntrackedParameter(std::string const &, T const &) const
double eta() const final
momentum pseudorapidity
double pt() const final
transverse momentum
int charge() const final
electric charge
Definition: LeafCandidate.h:91
bool isEE() const
Definition: GsfElectron.h:357
bool isEB() const
Definition: GsfElectron.h:356
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float sigmaIetaIeta() const
Definition: GsfElectron.h:440
float hcalOverEcal() const
Definition: GsfElectron.h:448
float dr03TkSumPt() const
Definition: GsfElectron.h:552
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
T const * product() const
Definition: Handle.h:74
bool electronSelection(const GsfElectron *eleRef, const math::XYZPoint &bspotPosition)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
float dr03EcalRecHitSumEt() const
Definition: GsfElectron.h:554
fixed size matrix
float dr03HcalTowerSumEt() const
Definition: GsfElectron.h:557
const Point & position() const
position
Definition: BeamSpot.h:62
double phi() const final
momentum azimuthal angle
bool filter(edm::Event &, const edm::EventSetup &) override
Definition: event.py:1
WZInterestingEventSelector(const edm::ParameterSet &)