CMS 3D CMS Logo

ECALMultifitAnalyzer_HI.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: DQMOffline/ECALMultifitAnalyzer_HI
4 // Class: ECALMultifitAnalyzer_HI
5 //
12 //
13 // Original Author: R. Alex Barbieri
14 // Created: Tue, 17 Nov 2015 17:18:50 GMT
15 // Modified for DQM: Raghav Kunnawalkam Elayavalli
16 // Wednesday 18 Nov 2015 9:45:AM GVA time
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 #include <vector>
23 
24 // user include files
28 
31 
33 
37 
43 
48 
50 
51 #include "TH2F.h"
52 
53 //
54 // class declaration
55 //
56 
58 public:
61 
62 private:
63  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
64  void analyze(const edm::Event &, const edm::EventSetup &) override;
65 
70 
76 
89 
102 };
103 
104 //
105 // constructors and destructor
106 //
108  :
109 
110  recoPhotonsCollection_(consumes<std::vector<reco::Photon>>(iConfig.getParameter<edm::InputTag>("recoPhotonSrc"))),
111  caloJetToken_(consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("recoJetSrc"))),
112  RecHitCollection_EB_(consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("RecHitCollection_EB"))),
113  RecHitCollection_EE_(consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("RecHitCollection_EE"))),
114  mRechitEnergyThreshold(iConfig.getParameter<double>("rechitEnergyThreshold")),
115  mRecoPhotonPtThreshold(iConfig.getParameter<double>("recoPhotonPtThreshold")),
116  mRecoJetPtThreshold(iConfig.getParameter<double>("recoJetPtThreshold")),
117  mDeltaRPhotonThreshold(iConfig.getParameter<double>("deltaRPhotonThreshold")),
118  mDeltaRJetThreshold(iConfig.getParameter<double>("deltaRJetThreshold")) {}
119 
120 //
121 // member functions
122 //
123 
124 // ------------ method called for each event ------------
126  using namespace edm;
127 
129  iSetup.get<CaloGeometryRecord>().get(geomH);
130  const CaloGeometry *geom = geomH.product();
131 
132  Handle<std::vector<reco::Photon>> recoPhotonsHandle;
133  iEvent.getByToken(recoPhotonsCollection_, recoPhotonsHandle);
134 
135  Handle<reco::CaloJetCollection> recoJetHandle;
136  iEvent.getByToken(caloJetToken_, recoJetHandle);
137 
139  iEvent.getByToken(RecHitCollection_EB_, ebHandle);
140 
142  iEvent.getByToken(RecHitCollection_EE_, eeHandle);
143 
144  for (EcalRecHitCollection::const_iterator hit = ebHandle->begin(); hit != ebHandle->end(); ++hit) {
145  eb_chi2->Fill(hit->chi2());
146  eb_errors->Fill(hit->energyError());
147  double eta = geom->getGeometry(hit->detid())->getPosition().eta();
148  double phi = geom->getGeometry(hit->detid())->getPosition().phi();
149  eb_chi2_eta->Fill(eta, hit->chi2());
150  eb_errors_eta->Fill(eta, hit->energyError());
151  if (hit->energy() > mRechitEnergyThreshold) {
152  eb_chi2_e5->Fill(hit->chi2());
153  eb_errors_e5->Fill(hit->energyError());
154  eb_chi2_e5_eta->Fill(eta, hit->chi2());
155  eb_errors_e5_eta->Fill(eta, hit->energyError());
156  }
157 
158  for (std::vector<reco::Photon>::const_iterator pho = recoPhotonsHandle->begin(); pho != recoPhotonsHandle->end();
159  ++pho) {
160  if (pho->et() < mRecoPhotonPtThreshold)
161  continue;
162  double dr = reco::deltaR(eta, phi, pho->eta(), pho->phi());
163  if (dr < mDeltaRPhotonThreshold) {
164  eb_chi2_photon15->Fill(hit->chi2());
165  eb_errors_photon15->Fill(hit->energyError());
166  }
167  }
168  for (std::vector<reco::CaloJet>::const_iterator jet = recoJetHandle->begin(); jet != recoJetHandle->end(); ++jet) {
169  if (jet->pt() < mRecoJetPtThreshold)
170  continue;
171  double dr = reco::deltaR(eta, phi, jet->eta(), jet->phi());
172  if (dr < mDeltaRJetThreshold) {
173  eb_chi2_jet30->Fill(hit->chi2());
174  eb_errors_jet30->Fill(hit->energyError());
175  }
176  }
177  }
178 
179  for (EcalRecHitCollection::const_iterator hit = eeHandle->begin(); hit != eeHandle->end(); ++hit) {
180  ee_chi2->Fill(hit->chi2());
181  ee_errors->Fill(hit->energyError());
182  double eta = geom->getGeometry(hit->detid())->getPosition().eta();
183  double phi = geom->getGeometry(hit->detid())->getPosition().phi();
184  ee_chi2_eta->Fill(eta, hit->chi2());
185  ee_errors_eta->Fill(eta, hit->energyError());
186  if (hit->energy() > mRechitEnergyThreshold) {
187  ee_chi2_e5->Fill(hit->chi2());
188  ee_errors_e5->Fill(hit->energyError());
189  ee_chi2_e5_eta->Fill(eta, hit->chi2());
190  ee_errors_e5_eta->Fill(eta, hit->energyError());
191  }
192 
193  for (std::vector<reco::Photon>::const_iterator pho = recoPhotonsHandle->begin(); pho != recoPhotonsHandle->end();
194  ++pho) {
195  if (pho->et() < mRecoPhotonPtThreshold)
196  continue;
197  double dr = reco::deltaR(eta, phi, pho->eta(), pho->phi());
198  if (dr < mDeltaRPhotonThreshold) {
199  ee_chi2_photon15->Fill(hit->chi2());
200  ee_errors_photon15->Fill(hit->energyError());
201  }
202  }
203  for (std::vector<reco::CaloJet>::const_iterator jet = recoJetHandle->begin(); jet != recoJetHandle->end(); ++jet) {
204  if (jet->pt() < mRecoJetPtThreshold)
205  continue;
206  double dr = reco::deltaR(eta, phi, jet->eta(), jet->phi());
207  if (dr < mDeltaRJetThreshold) {
208  ee_chi2_jet30->Fill(hit->chi2());
209  ee_errors_jet30->Fill(hit->energyError());
210  }
211  }
212  }
213 }
214 
215 // ------------ method called once each job just before starting event loop
216 // ------------
218  iBooker.setCurrentFolder("EcalCalibration/HIN_MultiFitAnalyzer");
219 
220  eb_chi2 = nullptr;
221  eb_chi2_eta = nullptr;
222  eb_chi2_e5 = nullptr;
223  eb_chi2_e5_eta = nullptr;
224  eb_errors = nullptr;
225  eb_errors_eta = nullptr;
226  eb_errors_e5 = nullptr;
227  eb_errors_e5_eta = nullptr;
228  eb_chi2_photon15 = nullptr;
229  eb_errors_photon15 = nullptr;
230  eb_chi2_jet30 = nullptr;
231  eb_errors_jet30 = nullptr;
232 
233  ee_chi2 = nullptr;
234  ee_chi2_eta = nullptr;
235  ee_chi2_e5 = nullptr;
236  ee_chi2_e5_eta = nullptr;
237  ee_errors = nullptr;
238  ee_errors_eta = nullptr;
239  ee_errors_e5 = nullptr;
240  ee_errors_e5_eta = nullptr;
241  ee_chi2_photon15 = nullptr;
242  ee_errors_photon15 = nullptr;
243  ee_chi2_jet30 = nullptr;
244  ee_errors_jet30 = nullptr;
245 
246  const int nBins = 500;
247  const float maxChi2 = 70;
248  const float maxError = 0.5;
249 
250  TH2F *hProfile_Chi2 = new TH2F("hProfile_Chi2", "", nBins, -5, 5, nBins, 0, maxChi2);
251  TH2F *hProfile_Err = new TH2F("hProfile_Err", "", nBins, -5, 5, nBins, 0, maxError);
252 
253  eb_chi2 = iBooker.book1D("rechit_eb_chi2", "Rechit eb_chi2;chi2 fit value;", nBins, 0, maxChi2);
254  eb_chi2_eta = iBooker.book2D("rechit_eb_eta_Vs_mean_Chi2", hProfile_Chi2);
255  eb_chi2_e5 = iBooker.book1D(Form("rechit_eb_chi2_e%2.0f", mRechitEnergyThreshold),
256  Form("Rechit eb_chi2, e>%2.0fGeV;chi2 fit value;", mRechitEnergyThreshold),
257  nBins,
258  0,
259  maxChi2);
260  eb_chi2_e5_eta = iBooker.book2D(Form("rechit_eb_chi2_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Chi2);
261  eb_errors = iBooker.book1D("rechit_eb_errors", "Rechit eb_errors;error on the energy;", nBins, 0, maxError);
262  eb_errors_eta = iBooker.book2D("rechit_eb_errors_eta", hProfile_Err);
263  eb_errors_e5 = iBooker.book1D(Form("rechit_eb_errors_e%2.0f", mRechitEnergyThreshold),
264  "Rechit eb_errors, e>5GeV;error on the energy;",
265  nBins,
266  0,
267  maxError);
268  eb_errors_e5_eta = iBooker.book2D(Form("rechit_eb_errors_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Err);
269  eb_chi2_photon15 = iBooker.book1D(Form("rechit_eb_chi2_photon%2.0f", mRecoPhotonPtThreshold),
270  "Rechit eb_chi2 near photons;chi2 fit value;",
271  nBins,
272  0,
273  maxChi2);
274  eb_errors_photon15 = iBooker.book1D(Form("rechit_eb_errors_photon%2.0f", mRecoPhotonPtThreshold),
275  "Rechit eb_errors near photons;error on the energy;",
276  nBins,
277  0,
278  maxError);
279  eb_chi2_jet30 = iBooker.book1D(Form("rechit_eb_chi2_jet%2.0f", mRecoJetPtThreshold),
280  "Rechit eb_chi2 near jets;chi2 fit value;",
281  nBins,
282  0,
283  maxChi2);
284  eb_errors_jet30 = iBooker.book1D(Form("rechit_eb_errors_jet%2.0f", mRecoJetPtThreshold),
285  "Rechit eb_errors near jets;error on the energy;",
286  nBins,
287  0,
288  maxError);
289 
290  ee_chi2 = iBooker.book1D("rechit_ee_chi2", "Rechit ee_chi2;chi2 fit value;", nBins, 0, maxChi2);
291  ee_chi2_eta = iBooker.book2D("rechit_ee_chi2_eta", hProfile_Chi2);
292  ee_chi2_e5 = iBooker.book1D(Form("rechit_ee_chi2_e%2.0f", mRechitEnergyThreshold),
293  "Rechit ee_chi2, e>5GeV;chi2 fit value;",
294  nBins,
295  0,
296  maxChi2);
297  ee_chi2_e5_eta = iBooker.book2D(Form("rechit_ee_chi2_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Chi2);
298  ee_errors = iBooker.book1D("rechit_ee_errors", "Rechit ee_errors;error on the energy;", nBins, 0, maxError);
299  ee_errors_eta = iBooker.book2D("rechit_ee_errors_eta", hProfile_Err);
300  ee_errors_e5 = iBooker.book1D(Form("rechit_ee_errors_e%2.0f", mRechitEnergyThreshold),
301  "Rechit ee_errors, e>5GeV;error on the energy;",
302  nBins,
303  0,
304  maxError);
305  ee_errors_e5_eta = iBooker.book2D(Form("rechit_ee_errors_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Err);
306  ee_chi2_photon15 = iBooker.book1D(Form("rechit_ee_chi2_photon%2.0f", mRecoPhotonPtThreshold),
307  "Rechit ee_chi2 near photons;chi2 fit value;",
308  nBins,
309  0,
310  maxChi2);
311  ee_errors_photon15 = iBooker.book1D(Form("rechit_ee_errors_photon%2.0f", mRecoPhotonPtThreshold),
312  "Rechit ee_errors near photons;error on the energy;",
313  nBins,
314  0,
315  maxError);
316  ee_chi2_jet30 = iBooker.book1D(Form("rechit_ee_chi2_jet%2.0f", mRecoJetPtThreshold),
317  "Rechit ee_chi2 near jets;chi2 fit value;",
318  nBins,
319  0,
320  maxChi2);
321  ee_errors_jet30 = iBooker.book1D(Form("rechit_ee_errors_jet%2.0f", mRecoJetPtThreshold),
322  "Rechit ee_errors near jets;error on the energy;",
323  nBins,
324  0,
325  maxError);
326 
327  delete hProfile_Chi2;
328  delete hProfile_Err;
329 }
330 
331 // define this as a plug-in
edm::EDGetTokenT< EcalRecHitCollection > RecHitCollection_EB_
Definition: Photon.py:1
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< EcalRecHit >::const_iterator const_iterator
ECALMultifitAnalyzer_HI(const edm::ParameterSet &)
edm::EDGetTokenT< reco::CaloJetCollection > caloJetToken_
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
const_iterator end() const
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:71
edm::EDGetTokenT< EcalRecHitCollection > RecHitCollection_EE_
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:85
T const * product() const
Definition: ESHandle.h:86
edm::EDGetTokenT< std::vector< reco::Photon > > recoPhotonsCollection_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
const_iterator begin() const
Definition: Run.h:45
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects