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 //
11 //
12 // Original Author: R. Alex Barbieri
13 // Created: Tue, 17 Nov 2015 17:18:50 GMT
14 // Modified for DQM: Raghav Kunnawalkam Elayavalli
15 // Wednesday 18 Nov 2015 9:45:AM GVA time
16 //
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  recoPhotonsCollection_(consumes<std::vector<reco::Photon> > (iConfig.getParameter<edm::InputTag>("recoPhotonSrc"))),
110  caloJetToken_(consumes<reco::CaloJetCollection> (iConfig.getParameter<edm::InputTag> ("recoJetSrc"))),
111  RecHitCollection_EB_(consumes<EcalRecHitCollection> (iConfig.getParameter<edm::InputTag>("RecHitCollection_EB"))),
112  RecHitCollection_EE_(consumes<EcalRecHitCollection> (iConfig.getParameter<edm::InputTag>("RecHitCollection_EE"))),
113  mRechitEnergyThreshold (iConfig.getParameter<double> ("rechitEnergyThreshold")),
114  mRecoPhotonPtThreshold (iConfig.getParameter<double> ("recoPhotonPtThreshold")),
115  mRecoJetPtThreshold (iConfig.getParameter<double> ("recoJetPtThreshold")),
116  mDeltaRPhotonThreshold (iConfig.getParameter<double> ("deltaRPhotonThreshold")),
117  mDeltaRJetThreshold (iConfig.getParameter<double> ("deltaRJetThreshold"))
118 {
119 
120 }
121 
122 //
123 // member functions
124 //
125 
126 // ------------ method called for each event ------------
128 {
129 
130  using namespace edm;
131 
133  iSetup.get<CaloGeometryRecord>().get(geomH);
134  const CaloGeometry* geom = geomH.product();
135 
136  Handle<std::vector<reco::Photon> > recoPhotonsHandle;
137  iEvent.getByToken(recoPhotonsCollection_, recoPhotonsHandle);
138 
139  Handle<reco::CaloJetCollection> recoJetHandle;
140  iEvent.getByToken(caloJetToken_, recoJetHandle);
141 
143  iEvent.getByToken(RecHitCollection_EB_, ebHandle);
144 
146  iEvent.getByToken(RecHitCollection_EE_, eeHandle);
147 
148  for(EcalRecHitCollection::const_iterator hit = ebHandle->begin(); hit != ebHandle->end(); ++hit) {
149  eb_chi2->Fill(hit->chi2() );
150  eb_errors->Fill(hit->energyError() );
151  double eta = geom->getGeometry(hit->detid())->getPosition().eta();
152  double phi = geom->getGeometry(hit->detid())->getPosition().phi();
153  eb_chi2_eta->Fill(eta, hit->chi2() );
154  eb_errors_eta->Fill(eta, hit->energyError() );
155  if(hit->energy() > mRechitEnergyThreshold)
156  {
157  eb_chi2_e5->Fill(hit->chi2() );
158  eb_errors_e5->Fill(hit->energyError() );
159  eb_chi2_e5_eta->Fill(eta, hit->chi2() );
160  eb_errors_e5_eta->Fill(eta, hit->energyError() );
161  }
162 
163  for (std::vector<reco::Photon>::const_iterator pho = recoPhotonsHandle->begin(); pho != recoPhotonsHandle->end(); ++pho) {
164  if(pho->et() < mRecoPhotonPtThreshold ) continue;
165  double dr = reco::deltaR(eta, phi, pho->eta(), pho->phi());
166  if(dr < mDeltaRPhotonThreshold)
167  {
168  eb_chi2_photon15->Fill(hit->chi2() );
169  eb_errors_photon15->Fill(hit->energyError() );
170  }
171  }
172  for (std::vector<reco::CaloJet>::const_iterator jet = recoJetHandle->begin(); jet != recoJetHandle->end(); ++jet) {
173  if(jet->pt() < mRecoJetPtThreshold ) continue;
174  double dr = reco::deltaR(eta, phi, jet->eta(), jet->phi());
175  if(dr < mDeltaRJetThreshold)
176  {
177  eb_chi2_jet30->Fill(hit->chi2() );
178  eb_errors_jet30->Fill(hit->energyError() );
179  }
180  }
181  }
182 
183  for(EcalRecHitCollection::const_iterator hit = eeHandle->begin(); hit != eeHandle->end(); ++hit) {
184  ee_chi2->Fill(hit->chi2() );
185  ee_errors->Fill(hit->energyError() );
186  double eta = geom->getGeometry(hit->detid())->getPosition().eta();
187  double phi = geom->getGeometry(hit->detid())->getPosition().phi();
188  ee_chi2_eta->Fill(eta, hit->chi2() );
189  ee_errors_eta->Fill(eta, hit->energyError() );
190  if(hit->energy() > mRechitEnergyThreshold)
191  {
192  ee_chi2_e5->Fill(hit->chi2() );
193  ee_errors_e5->Fill(hit->energyError() );
194  ee_chi2_e5_eta->Fill(eta, hit->chi2() );
195  ee_errors_e5_eta->Fill(eta, hit->energyError() );
196  }
197 
198  for (std::vector<reco::Photon>::const_iterator pho = recoPhotonsHandle->begin(); pho != recoPhotonsHandle->end(); ++pho) {
199  if(pho->et() < mRecoPhotonPtThreshold ) continue;
200  double dr = reco::deltaR(eta, phi, pho->eta(), pho->phi());
201  if(dr < mDeltaRPhotonThreshold)
202  {
203  ee_chi2_photon15->Fill(hit->chi2() );
204  ee_errors_photon15->Fill(hit->energyError() );
205  }
206  }
207  for (std::vector<reco::CaloJet>::const_iterator jet = recoJetHandle->begin(); jet != recoJetHandle->end(); ++jet) {
208  if(jet->pt() < mRecoJetPtThreshold ) continue;
209  double dr = reco::deltaR(eta, phi, jet->eta(), jet->phi());
210  if(dr < mDeltaRJetThreshold)
211  {
212  ee_chi2_jet30->Fill(hit->chi2() );
213  ee_errors_jet30->Fill(hit->energyError() );
214  }
215  }
216  }
217 }
218 
219 
220 // ------------ method called once each job just before starting event loop ------------
221 void
223 {
224 
225  iBooker.setCurrentFolder ("EcalCalibration/HIN_MultiFitAnalyzer");
226 
227  eb_chi2 = nullptr;
228  eb_chi2_eta = nullptr;
229  eb_chi2_e5 = nullptr;
230  eb_chi2_e5_eta = nullptr;
231  eb_errors = nullptr;
232  eb_errors_eta = nullptr;
233  eb_errors_e5 = nullptr;
234  eb_errors_e5_eta = nullptr;
235  eb_chi2_photon15 = nullptr;
236  eb_errors_photon15 = nullptr;
237  eb_chi2_jet30 = nullptr;
238  eb_errors_jet30 = nullptr;
239 
240  ee_chi2 = nullptr;
241  ee_chi2_eta = nullptr;
242  ee_chi2_e5 = nullptr;
243  ee_chi2_e5_eta = nullptr;
244  ee_errors = nullptr;
245  ee_errors_eta = nullptr;
246  ee_errors_e5 = nullptr;
247  ee_errors_e5_eta = nullptr;
248  ee_chi2_photon15 = nullptr;
249  ee_errors_photon15 = nullptr;
250  ee_chi2_jet30 = nullptr;
251  ee_errors_jet30 = nullptr;
252 
253  const int nBins = 500;
254  const float maxChi2 = 70;
255  const float maxError = 0.5;
256 
257  TH2F * hProfile_Chi2 = new TH2F("hProfile_Chi2","",nBins, -5, 5, nBins, 0, maxChi2);
258  TH2F * hProfile_Err = new TH2F("hProfile_Err","",nBins, -5, 5, nBins, 0, maxError);
259 
260  eb_chi2 = iBooker.book1D("rechit_eb_chi2","Rechit eb_chi2;chi2 fit value;",nBins,0,maxChi2);
261  eb_chi2_eta = iBooker.book2D("rechit_eb_eta_Vs_mean_Chi2", hProfile_Chi2);
262  eb_chi2_e5 = iBooker.book1D(Form("rechit_eb_chi2_e%2.0f", mRechitEnergyThreshold),Form("Rechit eb_chi2, e>%2.0fGeV;chi2 fit value;", mRechitEnergyThreshold),nBins,0,maxChi2);
263  eb_chi2_e5_eta = iBooker.book2D(Form("rechit_eb_chi2_e%2.0f_eta", mRechitEnergyThreshold),hProfile_Chi2);
264  eb_errors = iBooker.book1D("rechit_eb_errors","Rechit eb_errors;error on the energy;",nBins,0,maxError);
265  eb_errors_eta = iBooker.book2D("rechit_eb_errors_eta",hProfile_Err);
266  eb_errors_e5 = iBooker.book1D(Form("rechit_eb_errors_e%2.0f", mRechitEnergyThreshold),"Rechit eb_errors, e>5GeV;error on the energy;",nBins,0,maxError);
267  eb_errors_e5_eta = iBooker.book2D(Form("rechit_eb_errors_e%2.0f_eta", mRechitEnergyThreshold),hProfile_Err);
268  eb_chi2_photon15 = iBooker.book1D(Form("rechit_eb_chi2_photon%2.0f", mRecoPhotonPtThreshold),"Rechit eb_chi2 near photons;chi2 fit value;",nBins,0,maxChi2);
269  eb_errors_photon15 = iBooker.book1D(Form("rechit_eb_errors_photon%2.0f", mRecoPhotonPtThreshold),"Rechit eb_errors near photons;error on the energy;",nBins,0,maxError);
270  eb_chi2_jet30 = iBooker.book1D(Form("rechit_eb_chi2_jet%2.0f", mRecoJetPtThreshold),"Rechit eb_chi2 near jets;chi2 fit value;",nBins,0,maxChi2);
271  eb_errors_jet30 = iBooker.book1D(Form("rechit_eb_errors_jet%2.0f", mRecoJetPtThreshold),"Rechit eb_errors near jets;error on the energy;",nBins,0,maxError);
272 
273  ee_chi2 = iBooker.book1D("rechit_ee_chi2","Rechit ee_chi2;chi2 fit value;",nBins,0,maxChi2);
274  ee_chi2_eta = iBooker.book2D("rechit_ee_chi2_eta",hProfile_Chi2);
275  ee_chi2_e5 = iBooker.book1D(Form("rechit_ee_chi2_e%2.0f", mRechitEnergyThreshold),"Rechit ee_chi2, e>5GeV;chi2 fit value;",nBins,0,maxChi2);
276  ee_chi2_e5_eta = iBooker.book2D(Form("rechit_ee_chi2_e%2.0f_eta", mRechitEnergyThreshold),hProfile_Chi2);
277  ee_errors = iBooker.book1D("rechit_ee_errors","Rechit ee_errors;error on the energy;",nBins,0,maxError);
278  ee_errors_eta = iBooker.book2D("rechit_ee_errors_eta",hProfile_Err);
279  ee_errors_e5 = iBooker.book1D(Form("rechit_ee_errors_e%2.0f", mRechitEnergyThreshold),"Rechit ee_errors, e>5GeV;error on the energy;",nBins,0,maxError);
280  ee_errors_e5_eta = iBooker.book2D(Form("rechit_ee_errors_e%2.0f_eta", mRechitEnergyThreshold),hProfile_Err);
281  ee_chi2_photon15 = iBooker.book1D(Form("rechit_ee_chi2_photon%2.0f", mRecoPhotonPtThreshold),"Rechit ee_chi2 near photons;chi2 fit value;",nBins,0,maxChi2);
282  ee_errors_photon15 = iBooker.book1D(Form("rechit_ee_errors_photon%2.0f", mRecoPhotonPtThreshold),"Rechit ee_errors near photons;error on the energy;",nBins,0,maxError);
283  ee_chi2_jet30 = iBooker.book1D(Form("rechit_ee_chi2_jet%2.0f", mRecoJetPtThreshold),"Rechit ee_chi2 near jets;chi2 fit value;",nBins,0,maxChi2);
284  ee_errors_jet30 = iBooker.book1D(Form("rechit_ee_errors_jet%2.0f", mRecoJetPtThreshold),"Rechit ee_errors near jets;error on the energy;",nBins,0,maxError);
285 
286  delete hProfile_Chi2;
287  delete hProfile_Err;
288 
289 }
290 
291 
292 //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:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< EcalRecHit >::const_iterator const_iterator
ECALMultifitAnalyzer_HI(const edm::ParameterSet &)
edm::EDGetTokenT< reco::CaloJetCollection > caloJetToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
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
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:63
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:44
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects