CMS 3D CMS Logo

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 // system include files
20 #include <memory>
21 
22 // user include files
23 
24 #include "DQM/Physics/src/EwkDQM.h"
35 #include "TLorentzVector.h"
36 #include "TMath.h"
37 #include <cmath>
38 #include <iostream>
39 #include <string>
40 
43 
44 class EcalZmassTask : public DQMEDAnalyzer {
45 public:
46  explicit EcalZmassTask(const edm::ParameterSet &);
47  ~EcalZmassTask() override {}
48 
49 private:
50  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
51  void analyze(const edm::Event &, const edm::EventSetup &) override;
52 
55 
57 
71 };
72 
75  consumes<reco::GsfElectronCollection>(parameters.getParameter<edm::InputTag>("electronCollection"))),
77  consumes<reco::GsfTrackCollection>(parameters.getParameter<edm::InputTag>("trackCollection"))),
78  prefixME_(parameters.getUntrackedParameter<std::string>("prefixME", "")) {}
79 
80 // ------------ method called for each event ------------
82  using namespace edm;
84  iEvent.getByToken(electronCollectionToken_, electronCollection);
85  if (!electronCollection.isValid())
86  return;
87 
88  // get GSF Tracks
90  iEvent.getByToken(trackCollectionToken_, gsftracks_h);
91 
92  bool isIsolatedBarrel;
93  bool isIDBarrel;
94  bool isConvertedBarrel;
95  bool isIsolatedEndcap;
96  bool isIDEndcap;
97  bool isConvertedEndcap;
98 
99  int elIsAccepted = 0;
100  int elIsAcceptedEB = 0;
101  int elIsAcceptedEE = 0;
102 
103  std::vector<TLorentzVector> LV;
104 
105  for (reco::GsfElectronCollection::const_iterator recoElectron = electronCollection->begin();
106  recoElectron != electronCollection->end();
107  recoElectron++) {
108  if (recoElectron->et() <= 25)
109  continue;
110 
111  // Define Isolation variables
112  double IsoTrk = (recoElectron->dr03TkSumPt() / recoElectron->et());
113  double IsoEcal = (recoElectron->dr03EcalRecHitSumEt() / recoElectron->et());
114  double IsoHcal = (recoElectron->dr03HcalTowerSumEt() / recoElectron->et());
115  double HE = (recoElectron->hcalOverEcal());
116 
117  // Define ID variables
118 
119  float DeltaPhiTkClu = recoElectron->deltaPhiSuperClusterTrackAtVtx();
120  float DeltaEtaTkClu = recoElectron->deltaEtaSuperClusterTrackAtVtx();
121  float sigmaIeIe = recoElectron->sigmaIetaIeta();
122 
123  // Define Conversion Rejection Variables
124 
125  float Dcot = recoElectron->convDcot();
126  float Dist = recoElectron->convDist();
127  int NumberOfExpectedInnerHits =
128  recoElectron->gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
129 
130  // quality flags
131 
132  isIsolatedBarrel = false;
133  isIDBarrel = false;
134  isConvertedBarrel = false;
135  isIsolatedEndcap = false;
136  isIDEndcap = false;
137  isConvertedEndcap = false;
138 
139  /***** Barrel WP80 Cuts *****/
140 
141  if (fabs(recoElectron->eta()) <= 1.4442) {
142  /* Isolation */
143  if (IsoTrk < 0.09 && IsoEcal < 0.07 && IsoHcal < 0.10) {
144  isIsolatedBarrel = true;
145  }
146 
147  /* Identification */
148  if (fabs(DeltaEtaTkClu) < 0.004 && fabs(DeltaPhiTkClu) < 0.06 && sigmaIeIe < 0.01 && HE < 0.04) {
149  isIDBarrel = true;
150  }
151 
152  /* Conversion Rejection */
153  if ((fabs(Dist) >= 0.02 || fabs(Dcot) >= 0.02) && NumberOfExpectedInnerHits <= 1.0) {
154  isConvertedBarrel = true;
155  }
156  }
157 
158  if (isIsolatedBarrel && isIDBarrel && isConvertedBarrel) {
159  elIsAccepted++;
160  elIsAcceptedEB++;
161  TLorentzVector b_e2(
162  recoElectron->momentum().x(), recoElectron->momentum().y(), recoElectron->momentum().z(), recoElectron->p());
163  LV.push_back(b_e2);
164  }
165 
166  /***** Endcap WP80 Cuts *****/
167 
168  if (fabs(recoElectron->eta()) >= 1.5660 && fabs(recoElectron->eta()) <= 2.5000) {
169  /* Isolation */
170  if (IsoTrk < 0.04 && IsoEcal < 0.05 && IsoHcal < 0.025) {
171  isIsolatedEndcap = true;
172  }
173 
174  /* Identification */
175  if (fabs(DeltaEtaTkClu) < 0.007 && fabs(DeltaPhiTkClu) < 0.03 && sigmaIeIe < 0.031 && HE < 0.15) {
176  isIDEndcap = true;
177  }
178 
179  /* Conversion Rejection */
180  if ((fabs(Dcot) > 0.02 || fabs(Dist) > 0.02) && NumberOfExpectedInnerHits <= 1.0) {
181  isConvertedEndcap = true;
182  }
183  }
184  if (isIsolatedEndcap && isIDEndcap && isConvertedEndcap) {
185  elIsAccepted++;
186  elIsAcceptedEE++;
187  TLorentzVector e_e2(
188  recoElectron->momentum().x(), recoElectron->momentum().y(), recoElectron->momentum().z(), recoElectron->p());
189  LV.push_back(e_e2);
190  }
191  }
192 
193  // Calculate the Z invariant masses
194 
195  if (elIsAccepted > 1) {
196  double e_ee_invMass = 0;
197  if (LV.size() == 2) {
198  TLorentzVector e_pair = LV[0] + LV[1];
199  e_ee_invMass = e_pair.M();
200  }
201 
202  if (elIsAcceptedEB == 2) {
203  h_ee_invMass_BB->Fill(e_ee_invMass);
204  }
205  if (elIsAcceptedEE == 2) {
206  h_ee_invMass_EE->Fill(e_ee_invMass);
207  }
208  if (elIsAcceptedEB == 1 && elIsAcceptedEE == 1) {
209  h_ee_invMass_EB->Fill(e_ee_invMass);
210  }
211 
212  LV.clear();
213  }
214 }
215 
217  std::string logTraceName("EcalZmassTask");
218 
219  h_ee_invMass_EB = nullptr;
220  h_ee_invMass_EE = nullptr;
221  h_ee_invMass_BB = nullptr;
222 
223  LogTrace(logTraceName) << "Parameters initialization";
224 
225  iBooker.setCurrentFolder(prefixME_ + "/Zmass"); // Use folder with name of PAG
226 
227  h_ee_invMass_EB = iBooker.book1D("Z peak - WP80 EB-EE", "Z peak - WP80 EB-EE;InvMass (GeV)", 60, 60.0, 120.0);
228  h_ee_invMass_EE = iBooker.book1D("Z peak - WP80 EE-EE", "Z peak - WP80 EE-EE;InvMass (Gev)", 60, 60.0, 120.0);
229  h_ee_invMass_BB = iBooker.book1D("Z peak - WP80 EB-EB", "Z peak - WP80 EB-EB;InvMass (Gev)", 60, 60.0, 120.0);
230 }
231 
232 // define this as a plug-in
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
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:525
math::XYZTLorentzVectorD LV
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
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
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void Fill(long long x)
EcalZmassTask(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< GsfTrack > GsfTrackCollection
collection of GsfTracks
Definition: GsfTrackFwd.h:9
bool isValid() const
Definition: HandleBase.h:70
MonitorElement * h_ee_invMass_EE
#define LogTrace(id)
MonitorElement * h_ee_invMass
void analyze(const edm::Event &, const edm::EventSetup &) override
~EcalZmassTask() override
MonitorElement * h_95_ee_invMass_BB
MonitorElement * h_ee_invMass_BB
fixed size matrix
HLT enums.
const edm::EDGetTokenT< reco::GsfTrackCollection > trackCollectionToken_
MonitorElement * h_e2_phi
MonitorElement * h_e1_eta
Definition: Run.h:45