CMS 3D CMS Logo

ZDCDigiStudy.cc
Go to the documentation of this file.
1 // Package: ZDCDigiStudy
3 // Class: ZDCDigiStudy
4 //
5 /*
6  Description:
7  This code has been developed to be a check for the ZDC sim. In 2009, it was found that the ZDC Simulation was unrealistic and needed repair. The aim of this code is to show the user the input and output of a ZDC MinBias simulation.
8 
9  Implementation:
10  First a MinBias simulation should be run, it could be pythia,hijin,or hydjet. This will output a .root file which should have information about recoGenParticles, hcalunsuppresseddigis. Use this .root file as the input into the cfg.py which is found in the main directory of this package. This output will be another .root file which is meant to be viewed in a TBrowser
11 
12 */
13 //
14 // Original Author: Jaime Gomez (U. of Maryland) with SIGNIFICANT assistance of Dr. Jefferey Temple (U. of Maryland)
15 //
16 //
17 // Created: Summer 2012
19 
22 
24 #include "CLHEP/Units/GlobalSystemOfUnits.h"
25 
26 //#define EDM_ML_DEBUG
27 
29  : zdcHits(ps.getUntrackedParameter<std::string>("HitCollection", "ZdcHits")),
30  outFile_(ps.getUntrackedParameter<std::string>("outputFile", "zdcHitStudy.root")),
31  verbose_(ps.getUntrackedParameter<bool>("Verbose", false)),
32  checkHit_(true),
33  tok_zdc_(consumes<ZDCDigiCollection>(edm::InputTag("simHcalUnsuppressedDigis"))) {
34  edm::LogVerbatim("ZDCDigiStudy") << " Hits: " << zdcHits << " / " << checkHit_ << " Output: " << outFile_;
35 }
36 
38  ib.setCurrentFolder("ZDCDigiValidation");
39  // run histos only since there is dqmEndRun processing.
40  ib.setScope(MonitorElementData::Scope::RUN);
41 
42  //Histograms for Hits
44  //# Below we are filling the histograms made in the .h file. The syntax is as follows: #
45  //# plot_code_name = dbe_->TypeofPlot[(1,2,3)-D,(F,I,D)]("Name as it will appear","Title",axis options); #
46  //# They will be stored in the TFile subdirectory set by : dbe_->setCurrentFolder("FolderIwant") #
47  //# axis options are like (#ofbins,min,max) #
49 
50  if (checkHit_) {
52 
54  ib.setCurrentFolder("ZDCDigiValidation/ZDC_Digis/1D_fC");
55  meZdcfCPHAD = ib.book1D("PHAD_TotalfC", "PZDC_HAD_TotalfC", 1000, -50, 10000);
56  meZdcfCPHAD->setAxisTitle("Counts", 2);
57  meZdcfCPHAD->setAxisTitle("fC", 1);
59  meZdcfCPTOT = ib.book1D("PZDC_TotalfC", "PZDC_TotalfC", 1000, -50, 20000);
60  meZdcfCPTOT->setAxisTitle("Counts", 2);
61  meZdcfCPTOT->setAxisTitle("fC", 1);
63  meZdcfCNHAD = ib.book1D("NHAD_TotalfC", "NZDC_HAD_TotalfC", 1000, -50, 10000);
64  meZdcfCNHAD->setAxisTitle("Counts", 2);
65  meZdcfCNHAD->setAxisTitle("fC", 1);
67  meZdcfCNTOT = ib.book1D("NZDC_TotalfC", "NZDC_TotalfC", 1000, -50, 20000);
68  meZdcfCNTOT->setAxisTitle("Counts", 2);
69  meZdcfCNTOT->setAxisTitle("fC", 1);
71 
73  ib.setCurrentFolder("ZDCDigiValidation/ZDC_Digis/fCvsTS/PZDC");
74 
76  meZdcPEM1fCvsTS = ib.book1D("PEM1_fCvsTS", "P-EM1_AveragefC_vsTS", 10, 0, 9);
77  meZdcPEM1fCvsTS->setAxisTitle("fC", 2);
78  meZdcPEM1fCvsTS->setAxisTitle("TS", 1);
81  meZdcPEM2fCvsTS = ib.book1D("PEM2_fCvsTS", "P-EM2_AveragefC_vsTS", 10, 0, 9);
82  meZdcPEM2fCvsTS->setAxisTitle("fC", 2);
83  meZdcPEM2fCvsTS->setAxisTitle("TS", 1);
86  meZdcPEM3fCvsTS = ib.book1D("PEM3_fCvsTS", "P-EM3_AveragefC_vsTS", 10, 0, 9);
87  meZdcPEM3fCvsTS->setAxisTitle("fC", 2);
88  meZdcPEM3fCvsTS->setAxisTitle("TS", 1);
91  meZdcPEM4fCvsTS = ib.book1D("PEM4_fCvsTS", "P-EM4_AveragefC_vsTS", 10, 0, 9);
92  meZdcPEM4fCvsTS->setAxisTitle("fC", 2);
93  meZdcPEM4fCvsTS->setAxisTitle("TS", 1);
96  meZdcPEM5fCvsTS = ib.book1D("PEM5_fCvsTS", "P-EM5_AveragefC_vsTS", 10, 0, 9);
97  meZdcPEM5fCvsTS->setAxisTitle("fC", 2);
98  meZdcPEM5fCvsTS->setAxisTitle("TS", 1);
101  meZdcPHAD1fCvsTS = ib.book1D("PHAD1_fCvsTS", "P-HAD1_AveragefC_vsTS", 10, 0, 9);
102  meZdcPHAD1fCvsTS->setAxisTitle("fC", 2);
103  meZdcPHAD1fCvsTS->setAxisTitle("TS", 1);
106  meZdcPHAD2fCvsTS = ib.book1D("PHAD2_fCvsTS", "P-HAD2_AveragefC_vsTS", 10, 0, 9);
107  meZdcPHAD2fCvsTS->setAxisTitle("fC", 2);
108  meZdcPHAD2fCvsTS->setAxisTitle("TS", 1);
111  meZdcPHAD3fCvsTS = ib.book1D("PHAD3_fCvsTS", "P-HAD3_AveragefC_vsTS", 10, 0, 9);
112  meZdcPHAD3fCvsTS->setAxisTitle("fC", 2);
113  meZdcPHAD3fCvsTS->setAxisTitle("TS", 1);
116  meZdcPHAD4fCvsTS = ib.book1D("PHAD4_fCvsTS", "P-HAD4_AveragefC_vsTS", 10, 0, 9);
117  meZdcPHAD4fCvsTS->setAxisTitle("fC", 2);
118  meZdcPHAD4fCvsTS->setAxisTitle("TS", 1);
120  ib.setCurrentFolder("ZDCDigiValidation/ZDC_Digis/fCvsTS/NZDC");
121 
123  meZdcNEM1fCvsTS = ib.book1D("NEM1_fCvsTS", "N-EM1_AveragefC_vsTS", 10, 0, 9);
124  meZdcNEM1fCvsTS->setAxisTitle("fC", 2);
125  meZdcNEM1fCvsTS->setAxisTitle("TS", 1);
128  meZdcNEM2fCvsTS = ib.book1D("NEM2_fCvsTS", "N-EM2_AveragefC_vsTS", 10, 0, 9);
129  meZdcNEM2fCvsTS->setAxisTitle("fC", 2);
130  meZdcNEM2fCvsTS->setAxisTitle("TS", 1);
133  meZdcNEM3fCvsTS = ib.book1D("NEM3_fCvsTS", "N-EM3_AveragefC_vsTS", 10, 0, 9);
134  meZdcNEM3fCvsTS->setAxisTitle("fC", 2);
135  meZdcNEM3fCvsTS->setAxisTitle("TS", 1);
138  meZdcNEM4fCvsTS = ib.book1D("NEM4_fCvsTS", "N-EM4_AveragefC_vsTS", 10, 0, 9);
139  meZdcNEM4fCvsTS->setAxisTitle("fC", 2);
140  meZdcNEM4fCvsTS->setAxisTitle("TS", 1);
143  meZdcNEM5fCvsTS = ib.book1D("NEM5_fCvsTS", "N-EM5_AveragefC_vsTS", 10, 0, 9);
144  meZdcNEM5fCvsTS->setAxisTitle("fC", 2);
145  meZdcNEM5fCvsTS->setAxisTitle("TS", 1);
148  meZdcNHAD1fCvsTS = ib.book1D("NHAD1_fCvsTS", "N-HAD1_AveragefC_vsTS", 10, 0, 9);
149  meZdcNHAD1fCvsTS->setAxisTitle("fC", 2);
150  meZdcNHAD1fCvsTS->setAxisTitle("TS", 1);
153  meZdcNHAD2fCvsTS = ib.book1D("NHAD2_fCvsTS", "N-HAD2_AveragefC_vsTS", 10, 0, 9);
154  meZdcNHAD2fCvsTS->setAxisTitle("fC", 2);
155  meZdcNHAD2fCvsTS->setAxisTitle("TS", 1);
158  meZdcNHAD3fCvsTS = ib.book1D("NHAD3_fCvsTS", "N-HAD3_AveragefC_vsTS", 10, 0, 9);
159  meZdcNHAD3fCvsTS->setAxisTitle("fC", 2);
160  meZdcNHAD3fCvsTS->setAxisTitle("TS", 1);
163  meZdcNHAD4fCvsTS = ib.book1D("NHAD4_fCvsTS", "N-HAD4_AveragefC_vsTS", 10, 0, 9);
164  meZdcNHAD4fCvsTS->setAxisTitle("fC", 2);
165  meZdcNHAD4fCvsTS->setAxisTitle("TS", 1);
167 
169  ib.setCurrentFolder("ZDCDigiValidation/ZDC_Digis/2D_EMvHAD");
172  meZdcfCPEMvHAD = ib.book2D("PEMvPHAD", "PZDC_EMvHAD", 200, -25, 12000, 200, -25, 15000);
173  meZdcfCPEMvHAD->setAxisTitle("SumEM_fC", 2);
174  meZdcfCPEMvHAD->setAxisTitle("SumHAD_fC", 1);
175  meZdcfCPEMvHAD->setOption("colz");
177  meZdcfCNEMvHAD = ib.book2D("NEMvNHAD", "NZDC_EMvHAD", 1000, -25, 12000, 1000, -25, 15000);
178  meZdcfCNEMvHAD->setAxisTitle("SumEM_fC", 2);
179  meZdcfCNEMvHAD->setAxisTitle("SumHAD_fC", 1);
180  meZdcfCNEMvHAD->setOption("colz");
182  }
183 }
184 
185 /*void ZDCDigiStudy::endJob() {
186  if (dbe_ && outFile_.size() > 0) dbe_->save(outFile_);
187  }*/
188 
189 //void ZDCDigiStudy::analyze(const edm::Event& e, const edm::EventSetup& ) {
192 
193  using namespace edm;
194  bool gotZDCDigis = true;
195 
196  Handle<ZDCDigiCollection> zdchandle;
197  if (!(iEvent.getByToken(tok_zdc_, zdchandle))) {
198  gotZDCDigis = false; //this is a boolean set up to check if there are ZDCDigis in the input root file
199  }
200  if (!(zdchandle.isValid())) {
201  gotZDCDigis = false; //if it is not there, leave it false
202  }
203 
204  double totalPHADCharge = 0;
205  double totalNHADCharge = 0;
206  double totalPEMCharge = 0;
207  double totalNEMCharge = 0;
208  double totalPCharge = 0;
209  double totalNCharge = 0;
210 
212  if (gotZDCDigis == true) {
213  for (ZDCDigiCollection::const_iterator zdc = zdchandle->begin(); zdc != zdchandle->end(); ++zdc) {
214  const ZDCDataFrame digi = (const ZDCDataFrame)(*zdc);
215 #ifdef EDM_ML_DEBUG
216  edm::LogVerbatim("ZDCDigiStudy") << "CHANNEL = " << zdc->id().channel();
217 #endif
218 
220  if (digi.id().section() == 2) { // require HAD
221  if (digi.id().zside() == 1) { // require POS
222  for (int i = 0; i < digi.size(); ++i) // loop over all 10 TS because each digi has 10 entries
223  {
224  if (digi.id().channel() == 1) { //here i specify PHAD1
226  i, digi.sample(i).nominal_fC()); //filling the plot name with the nominal fC value for each TS
227  if (i == 0)
228  meZdcPHAD1fCvsTS->Fill(-1, 1); // on first iteration of loop, increment underflow bin
229  } //NEW AVERAGE Thingy
230  if (digi.id().channel() == 2) {
232  if (i == 0)
233  meZdcPHAD2fCvsTS->Fill(-1, 1);
234  }
235  if (digi.id().channel() == 3) {
237  if (i == 0)
238  meZdcPHAD3fCvsTS->Fill(-1, 1);
239  }
240  if (digi.id().channel() == 4) {
242  if (i == 0)
243  meZdcPHAD4fCvsTS->Fill(-1, 1);
244  }
245  if (i == 4 || i == 5 || i == 6)
246  totalPHADCharge += digi.sample(i).nominal_fC();
247  } // loop over all (10) TS for the given digi
248  } else {
249  for (int i = 0; i < digi.size(); ++i) {
250  if (digi.id().channel() == 1) {
252  if (i == 0)
253  meZdcNHAD1fCvsTS->Fill(-1, 1);
254  }
255  if (digi.id().channel() == 2) {
257  if (i == 0)
258  meZdcNHAD2fCvsTS->Fill(-1, 1);
259  }
260  if (digi.id().channel() == 3) {
262  if (i == 0)
263  meZdcNHAD3fCvsTS->Fill(-1, 1);
264  }
265  if (digi.id().channel() == 4) {
267  if (i == 0)
268  meZdcNHAD4fCvsTS->Fill(-1, 1);
269  }
270  if (i == 4 || i == 5 || i == 6)
271  totalNHADCharge += digi.sample(i).nominal_fC();
272  } //loop over all 10 TS
273  } //Requires NHAd
274  } //Requires HAD sections
276  if (digi.id().section() ==
277  1) { //require EM....here i do the smae thing that i did above but now for P/N EM sections
278  if (digi.id().zside() == 1) { //require pos
279  for (int i = 0; i < digi.size(); ++i) {
280  if (digi.id().channel() == 1) {
282  if (i == 0)
283  meZdcPEM1fCvsTS->Fill(-1, 1);
284  }
285  if (digi.id().channel() == 2) {
287  if (i == 0)
288  meZdcPEM2fCvsTS->Fill(-1, 1);
289  }
290  if (digi.id().channel() == 3) {
292  if (i == 0)
293  meZdcPEM3fCvsTS->Fill(-1, 1);
294  }
295  if (digi.id().channel() == 4) {
297  if (i == 0)
298  meZdcPEM4fCvsTS->Fill(-1, 1);
299  }
300  if (digi.id().channel() == 5) {
302  if (i == 0)
303  meZdcPEM5fCvsTS->Fill(-1, 1);
304  }
305  if (i == 4 || i == 5 || i == 6)
306  totalPEMCharge += digi.sample(i).nominal_fC();
307  }
308  } else {
309  for (int i = 0; i < digi.size(); ++i) {
310  if (digi.id().channel() == 1) {
312  if (i == 0)
313  meZdcNEM1fCvsTS->Fill(-1, 1);
314  }
315  if (digi.id().channel() == 2) {
317  if (i == 0)
318  meZdcNEM2fCvsTS->Fill(-1, 1);
319  }
320  if (digi.id().channel() == 3) {
322  if (i == 0)
323  meZdcNEM3fCvsTS->Fill(-1, 1);
324  }
325  if (digi.id().channel() == 4) {
327  if (i == 0)
328  meZdcNEM4fCvsTS->Fill(-1, 1);
329  }
330  if (digi.id().channel() == 5) {
332  if (i == 0)
333  meZdcNEM5fCvsTS->Fill(-1, 1);
334  }
335  if (i == 4 || i == 5 || i == 6)
336  totalNEMCharge += digi.sample(i).nominal_fC();
337  }
338  }
339  }
340 
341  totalPCharge = totalPHADCharge + (0.1) * totalPEMCharge;
342  totalNCharge = totalNHADCharge + (0.1) * totalNEMCharge;
343 
344 #ifdef EDM_ML_DEBUG
345  edm::LogVerbatim("ZDCDigiStudy") << "CHANNEL = " << digi.id().channel();
346  for (int i = 0; i < digi.size(); ++i)
347  edm::LogVerbatim("ZDCDigiStudy") << "SAMPLE = " << i << " ADC = " << digi.sample(i).adc()
348  << " fC = " << digi.sample(i).nominal_fC();
349 #endif
350  // digi[i] should be the sample as digi.sample(i), I think
351  } // loop on all (22) ZDC digis
352  }
354 
355  // Now fill total charge histogram
356  meZdcfCPEMvHAD->Fill(totalPCharge, totalPEMCharge);
357  meZdcfCNEMvHAD->Fill(totalNCharge, totalNEMCharge);
358  meZdcfCPHAD->Fill(totalPHADCharge);
359  meZdcfCNHAD->Fill(totalNHADCharge);
360  meZdcfCNTOT->Fill(totalNCharge);
361  meZdcfCPTOT->Fill(totalPCharge);
362 }
363 
365 
367  int nevents =
369  ->GetBinContent(
370  0); //grab the number of digis that were read in and stored in the underflow bin, and call them Nevents
372  ->Scale(
373  1. /
374  nevents); // divide histogram by nevents thereby creating an average..it was done this way so that in DQM when everything is done in parallel and added at the end then the average will add appropriately
375 
376  int nevents1 = (meZdcPHAD2fCvsTS->getTH1F())->GetBinContent(0);
377  (meZdcPHAD2fCvsTS->getTH1F())->Scale(1. / nevents1);
378 
379  int nevents2 = (meZdcPHAD3fCvsTS->getTH1F())->GetBinContent(0);
380  (meZdcPHAD3fCvsTS->getTH1F())->Scale(1. / nevents2);
381 
382  int nevents3 = (meZdcPHAD4fCvsTS->getTH1F())->GetBinContent(0);
383  (meZdcPHAD4fCvsTS->getTH1F())->Scale(1. / nevents3);
384 
385  int nevents4 = (meZdcNHAD1fCvsTS->getTH1F())->GetBinContent(0);
386  (meZdcNHAD1fCvsTS->getTH1F())->Scale(1. / nevents4);
387 
388  int nevents5 = (meZdcNHAD2fCvsTS->getTH1F())->GetBinContent(0);
389  (meZdcNHAD2fCvsTS->getTH1F())->Scale(1. / nevents5);
390 
391  int nevents6 = (meZdcNHAD3fCvsTS->getTH1F())->GetBinContent(0);
392  (meZdcNHAD3fCvsTS->getTH1F())->Scale(1. / nevents6);
393 
394  int nevents7 = (meZdcNHAD4fCvsTS->getTH1F())->GetBinContent(0);
395  (meZdcNHAD4fCvsTS->getTH1F())->Scale(1. / nevents7);
396 
397  int nevents8 = (meZdcPEM1fCvsTS->getTH1F())->GetBinContent(0);
398  (meZdcPEM1fCvsTS->getTH1F())->Scale(1. / nevents8);
399 
400  int nevents9 = (meZdcPEM2fCvsTS->getTH1F())->GetBinContent(0);
401  (meZdcPEM2fCvsTS->getTH1F())->Scale(1. / nevents9);
402 
403  int nevents10 = (meZdcPEM3fCvsTS->getTH1F())->GetBinContent(0);
404  (meZdcPEM3fCvsTS->getTH1F())->Scale(1. / nevents10);
405 
406  int nevents11 = (meZdcPEM4fCvsTS->getTH1F())->GetBinContent(0);
407  (meZdcPEM4fCvsTS->getTH1F())->Scale(1. / nevents11);
408 
409  int nevents12 = (meZdcPEM5fCvsTS->getTH1F())->GetBinContent(0);
410  (meZdcPEM5fCvsTS->getTH1F())->Scale(1. / nevents12);
411 
412  int nevents13 = (meZdcNEM1fCvsTS->getTH1F())->GetBinContent(0);
413  (meZdcNEM1fCvsTS->getTH1F())->Scale(1. / nevents13);
414 
415  int nevents14 = (meZdcNEM2fCvsTS->getTH1F())->GetBinContent(0);
416  (meZdcNEM2fCvsTS->getTH1F())->Scale(1. / nevents14);
417 
418  int nevents15 = (meZdcNEM3fCvsTS->getTH1F())->GetBinContent(0);
419  (meZdcNEM3fCvsTS->getTH1F())->Scale(1. / nevents15);
420 
421  int nevents16 = (meZdcNEM4fCvsTS->getTH1F())->GetBinContent(0);
422  (meZdcNEM4fCvsTS->getTH1F())->Scale(1. / nevents16);
423 
424  int nevents17 = (meZdcNEM5fCvsTS->getTH1F())->GetBinContent(0);
425  (meZdcNEM5fCvsTS->getTH1F())->Scale(1. / nevents17);
426 }
427 
428 //define this as a plug-in
Log< level::Info, true > LogVerbatim
const HcalQIESample & sample(int i) const
access a sample
Definition: ZDCDataFrame.h:39
virtual void setOption(const char *option)
constexpr double nominal_fC() const
get the nominal FC (no calibrations applied)
Definition: HcalQIESample.h:45
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ZDCDigiStudy.cc:37
std::vector< T >::const_iterator const_iterator
MonitorElement * meZdcfCPTOT
Definition: ZDCDigiStudy.h:73
MonitorElement * meZdcPEM1fCvsTS
Definition: ZDCDigiStudy.h:78
MonitorElement * meZdcfCPEMvHAD
Definition: ZDCDigiStudy.h:76
void dqmEndRun(const edm::Run &run, const edm::EventSetup &c) override
MonitorElement * meZdcPHAD2fCvsTS
Definition: ZDCDigiStudy.h:84
MonitorElement * meZdcfCNTOT
Definition: ZDCDigiStudy.h:75
MonitorElement * meZdcPHAD3fCvsTS
Definition: ZDCDigiStudy.h:85
MonitorElement * meZdcfCPHAD
Definition: ZDCDigiStudy.h:72
void Fill(long long x)
const std::string zdcHits
Definition: ZDCDigiStudy.h:64
MonitorElement * meZdcNHAD2fCvsTS
Definition: ZDCDigiStudy.h:93
void analyze(const edm::Event &e, const edm::EventSetup &c) override
int iEvent
Definition: GenABIO.cc:224
MonitorElement * meZdcPEM3fCvsTS
Definition: ZDCDigiStudy.h:80
const bool checkHit_
Definition: ZDCDigiStudy.h:65
MonitorElement * meZdcNHAD3fCvsTS
Definition: ZDCDigiStudy.h:94
MonitorElement * meZdcPEM2fCvsTS
Definition: ZDCDigiStudy.h:79
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * meZdcNEM2fCvsTS
Definition: ZDCDigiStudy.h:88
MonitorElement * meZdcNEM4fCvsTS
Definition: ZDCDigiStudy.h:90
MonitorElement * meZdcNEM3fCvsTS
Definition: ZDCDigiStudy.h:89
const_iterator begin() const
const HcalZDCDetId & id() const
Definition: ZDCDataFrame.h:22
MonitorElement * meZdcPHAD1fCvsTS
Definition: ZDCDigiStudy.h:83
int size() const
total number of samples in the digi
Definition: ZDCDataFrame.h:26
const_iterator end() const
MonitorElement * meZdcNHAD4fCvsTS
Definition: ZDCDigiStudy.h:95
MonitorElement * meZdcNHAD1fCvsTS
Definition: ZDCDigiStudy.h:92
const edm::EDGetTokenT< ZDCDigiCollection > tok_zdc_
Definition: ZDCDigiStudy.h:67
bool isValid() const
Definition: HandleBase.h:70
MonitorElement * meZdcfCNHAD
Definition: ZDCDigiStudy.h:74
constexpr Section section() const
get the section
Definition: HcalZDCDetId.h:92
HLT enums.
constexpr int adc() const
get the ADC sample
Definition: HcalQIESample.h:43
MonitorElement * meZdcfCNEMvHAD
Definition: ZDCDigiStudy.h:77
constexpr int32_t channel() const
get the channel
Definition: HcalZDCDetId.h:112
MonitorElement * meZdcPHAD4fCvsTS
Definition: ZDCDigiStudy.h:86
constexpr int32_t zside() const
get the z-side of the cell (1/-1)
Definition: HcalZDCDetId.h:90
ZDCDigiStudy(const edm::ParameterSet &ps)
Definition: ZDCDigiStudy.cc:28
MonitorElement * meZdcPEM4fCvsTS
Definition: ZDCDigiStudy.h:81
MonitorElement * meZdcPEM5fCvsTS
Definition: ZDCDigiStudy.h:82
MonitorElement * meZdcNEM1fCvsTS
Definition: ZDCDigiStudy.h:87
const std::string outFile_
Definition: ZDCDigiStudy.h:64
Definition: Run.h:45
MonitorElement * meZdcNEM5fCvsTS
Definition: ZDCDigiStudy.h:91
ib
Definition: cuy.py:661
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)