CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalZmassTask.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Zanalyzer
4 // Class: Zanalyzer
5 //
13 //
14 // Original Author: Vieri Candelise
15 // Created: Wed May 11 14:53:26 CEST 2011
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
24 
25 #include "DQM/Physics/src/EwkDQM.h"
31 #include "TLorentzVector.h"
32 #include "TMath.h"
33 #include <string>
34 #include <cmath>
38 #include <iostream>
41 
44 
46 
47 public:
48  explicit EcalZmassTask (const edm::ParameterSet &);
50 
51 private:
52  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
53  virtual void analyze (const edm::Event &, const edm::EventSetup &) override;
54 
57 
59 
73 
74 };
75 
77  electronCollectionToken_(consumes<reco::GsfElectronCollection>(parameters.getParameter < edm::InputTag > ("electronCollection"))),
78  trackCollectionToken_(consumes<reco::GsfTrackCollection>(parameters.getParameter<edm::InputTag>("trackCollection"))),
79  prefixME_(parameters.getUntrackedParameter < std::string > ("prefixME", ""))
80 {
81 }
82 
83 // ------------ method called for each event ------------
84 void
86  const edm::EventSetup & iSetup)
87 {
88  using namespace edm;
89  Handle < reco::GsfElectronCollection > electronCollection;
90  iEvent.getByToken (electronCollectionToken_, electronCollection);
91  if (!electronCollection.isValid ()) return;
92 
93  //get GSF Tracks
95  iEvent.getByToken (trackCollectionToken_, gsftracks_h);
96 
97  bool isIsolatedBarrel;
98  bool isIDBarrel;
99  bool isConvertedBarrel;
100  bool isIsolatedEndcap;
101  bool isIDEndcap;
102  bool isConvertedEndcap;
103 
104  int elIsAccepted=0;
105  int elIsAcceptedEB=0;
106  int elIsAcceptedEE=0;
107 
108  std::vector<TLorentzVector> LV;
109 
110  for (reco::GsfElectronCollection::const_iterator recoElectron =
111  electronCollection->begin ();
112  recoElectron != electronCollection->end (); recoElectron++)
113  {
114 
115  if (recoElectron->et () <= 25)
116  continue;
117 
118  // Define Isolation variables
119  double IsoTrk = (recoElectron->dr03TkSumPt () / recoElectron->et ());
120  double IsoEcal =
121  (recoElectron->dr03EcalRecHitSumEt () / recoElectron->et ());
122  double IsoHcal =
123  (recoElectron->dr03HcalTowerSumEt () / recoElectron->et ());
124  double HE = (recoElectron->hcalOverEcal ());
125 
126  //Define ID variables
127 
128  float DeltaPhiTkClu = recoElectron->deltaPhiSuperClusterTrackAtVtx ();
129  float DeltaEtaTkClu = recoElectron->deltaEtaSuperClusterTrackAtVtx ();
130  float sigmaIeIe = recoElectron->sigmaIetaIeta ();
131 
132  //Define Conversion Rejection Variables
133 
134  float Dcot = recoElectron->convDcot ();
135  float Dist = recoElectron->convDist ();
136  int NumberOfExpectedInnerHits = recoElectron->gsfTrack()->hitPattern().numberOfHits(reco::HitPattern::MISSING_INNER_HITS);
137 
138  //quality flags
139 
140  isIsolatedBarrel = false;
141  isIDBarrel = false;
142  isConvertedBarrel = false;
143  isIsolatedEndcap = false;
144  isIDEndcap = false;
145  isConvertedEndcap = false;
146 
147 
148  /***** Barrel WP80 Cuts *****/
149 
150  if (fabs (recoElectron->eta ()) <= 1.4442)
151  {
152 
153  /* Isolation */
154  if (IsoTrk < 0.09 && IsoEcal < 0.07 && IsoHcal < 0.10)
155  {
156  isIsolatedBarrel = true;
157  }
158 
159  /* Identification */
160  if (fabs (DeltaEtaTkClu) < 0.004 && fabs (DeltaPhiTkClu) < 0.06
161  && sigmaIeIe < 0.01 && HE < 0.04)
162  {
163  isIDBarrel = true;
164  }
165 
166  /* Conversion Rejection */
167  if ((fabs (Dist) >= 0.02 || fabs (Dcot) >= 0.02)
168  && NumberOfExpectedInnerHits <= 1.0)
169  {
170  isConvertedBarrel = true;
171  }
172  }
173 
174  if (isIsolatedBarrel && isIDBarrel && isConvertedBarrel) {
175  elIsAccepted++;
176  elIsAcceptedEB++;
177  TLorentzVector b_e2(recoElectron->momentum ().x (),recoElectron->momentum ().y (),recoElectron->momentum ().z (), recoElectron->p ());
178  LV.push_back(b_e2);
179  }
180 
181  /***** Endcap WP80 Cuts *****/
182 
183  if (fabs (recoElectron->eta ()) >= 1.5660
184  && fabs (recoElectron->eta ()) <= 2.5000)
185  {
186 
187  /* Isolation */
188  if (IsoTrk < 0.04 && IsoEcal < 0.05 && IsoHcal < 0.025)
189  {
190  isIsolatedEndcap = true;
191  }
192 
193  /* Identification */
194  if (fabs (DeltaEtaTkClu) < 0.007 && fabs (DeltaPhiTkClu) < 0.03
195  && sigmaIeIe < 0.031 && HE < 0.15)
196  {
197  isIDEndcap = true;
198  }
199 
200  /* Conversion Rejection */
201  if ((fabs (Dcot) > 0.02 || fabs (Dist) > 0.02)
202  && NumberOfExpectedInnerHits <= 1.0)
203  {
204  isConvertedEndcap = true;
205  }
206  }
207  if (isIsolatedEndcap && isIDEndcap && isConvertedEndcap) {
208  elIsAccepted++;
209  elIsAcceptedEE++;
210  TLorentzVector e_e2(recoElectron->momentum ().x (),recoElectron->momentum ().y (),recoElectron->momentum ().z (), recoElectron->p ());
211  LV.push_back(e_e2);
212  }
213 
214  }
215 
216  // Calculate the Z invariant masses
217 
218  if (elIsAccepted>1){
219  double e_ee_invMass=0;
220  if (LV.size()==2){
221  TLorentzVector e_pair = LV[0] + LV[1];
222  e_ee_invMass = e_pair.M ();
223  }
224 
225  if (elIsAcceptedEB==2){
226  h_ee_invMass_BB->Fill(e_ee_invMass);
227  }
228  if (elIsAcceptedEE==2){
229  h_ee_invMass_EE->Fill(e_ee_invMass);
230  }
231  if (elIsAcceptedEB==1 && elIsAcceptedEE==1){
232  h_ee_invMass_EB->Fill(e_ee_invMass);
233  }
234 
235  LV.clear();
236 
237  }
238 }
239 
240 void
242 {
243  std::string logTraceName("EcalZmassTask");
244 
245  h_ee_invMass_EB = 0;
246  h_ee_invMass_EE = 0;
247  h_ee_invMass_BB = 0;
248 
249  LogTrace (logTraceName) << "Parameters initialization";
250 
251  iBooker.setCurrentFolder (prefixME_ + "/Zmass"); // Use folder with name of PAG
252 
254  iBooker.book1D ("Z peak - WP80 EB-EE",
255  "Z peak - WP80 EB-EE;InvMass (GeV)", 60, 60.0, 120.0);
257  iBooker.book1D ("Z peak - WP80 EE-EE",
258  "Z peak - WP80 EE-EE;InvMass (Gev)", 60, 60.0, 120.0);
260  iBooker.book1D ("Z peak - WP80 EB-EB",
261  "Z peak - WP80 EB-EB;InvMass (Gev)", 60, 60.0, 120.0);
262 }
263 
264 //define this as a plug-in
MonitorElement * h_e2_eta
MonitorElement * h_ee_invMass_EB
const std::string prefixME_
MonitorElement * h_e1_phi
MonitorElement * h_95_ee_invMass_EE
MonitorElement * h_95_ee_invMass_EB
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * h_e2_et
const edm::EDGetTokenT< reco::GsfElectronCollection > electronCollectionToken_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * h_e1_et
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
void Fill(long long x)
EcalZmassTask(const edm::ParameterSet &)
math::XYZTLorentzVectorD LV
int iEvent
Definition: GenABIO.cc:230
std::vector< GsfTrack > GsfTrackCollection
collection of GsfTracks
Definition: GsfTrackFwd.h:9
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * h_ee_invMass_EE
#define LogTrace(id)
MonitorElement * h_ee_invMass
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:273
MonitorElement * h_95_ee_invMass_BB
MonitorElement * h_ee_invMass_BB
const edm::EDGetTokenT< reco::GsfTrackCollection > trackCollectionToken_
MonitorElement * h_e2_phi
MonitorElement * h_e1_eta
Definition: Run.h:43