CMS 3D CMS Logo

HcalTB02Histo.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestBeam
4 // Class : HcalTB02Histo
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author:
10 // Created: Sun May 21 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 #include <iostream>
15 #include <cmath>
16 
17 // user include files
22 
23 //#define EDM_ML_DEBUG
24 //
25 // constructors and destructor
26 HcalTB02Histo::HcalTB02Histo(const edm::ParameterSet& ps) : rt_tbTimes(nullptr), rt_TransProf(nullptr) {
27  verbose = ps.getUntrackedParameter<bool>("Verbose", false);
28 #ifdef EDM_ML_DEBUG
29  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: initialised with o/p file " << fileName << " verbosity " << verbose;
30 #endif
31  // Book histograms
33 
34  if (!tfile.isAvailable())
35  throw cms::Exception("BadConfig") << "TFileService unavailable: "
36  << "please add it to config file";
37 
38  char title[80];
39  for (int ilayer = 0; ilayer < 19; ilayer++) {
40  sprintf(title, "Scint. Energy in Layer %d ", ilayer);
41  TH1D* h = tfile->make<TH1D>(title, title, 500, 0., 1.5);
42  rt_histoProf.push_back(h);
43 #ifdef EDM_ML_DEBUG
44  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title;
45 #endif
46  }
47  sprintf(title, "All Hit Time slices");
48  rt_tbTimes = tfile->make<TH1D>(title, title, 100, 0., 100.);
49 #ifdef EDM_ML_DEBUG
50  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title;
51 #endif
52  sprintf(title, "Transv. Shower Profile");
53  rt_TransProf = tfile->make<TH2D>(title, title, 100, 0., 1., 1000, 0., 0.5);
54 #ifdef EDM_ML_DEBUG
55  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title;
56 #endif
57 }
58 
60 #ifdef EDM_ML_DEBUG
61  edm::LogVerbatim("HcalTBSim") << " Deleting HcalTB02Histo";
62 #endif
63 }
64 
65 //
66 // member functions
67 //
68 
70 #ifdef EDM_ML_DEBUG
71  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo::Fill Time histo with " << v;
72 #endif
73  if (rt_tbTimes) {
74  rt_tbTimes->Fill(v);
75  }
76 }
77 
78 void HcalTB02Histo::fillTransProf(float u, float v) {
79 #ifdef EDM_ML_DEBUG
80  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:::Fill Shower Transv. Profile histo"
81  << " with " << u << " and " << v;
82 #endif
83  if (rt_TransProf) {
84  rt_TransProf->Fill(u, v);
85  }
86 }
87 
88 void HcalTB02Histo::fillProfile(int ilayer, float value) {
89  if (ilayer < (int)(rt_histoProf.size())) {
90  rt_histoProf[ilayer]->Fill(value);
91 #ifdef EDM_ML_DEBUG
92  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo::Fill profile " << ilayer << " with " << value;
93 #endif
94  }
95 }
96 
97 float HcalTB02Histo::getMean(int ilayer) {
98  if (ilayer < (int)(rt_histoProf.size())) {
99  return rt_histoProf[ilayer]->GetMean();
100  } else {
101  return 0;
102  }
103 }
104 
105 float HcalTB02Histo::getRMS(int ilayer) {
106  if (ilayer < (int)(rt_histoProf.size())) {
107  return rt_histoProf[ilayer]->GetRMS();
108  } else {
109  return 0;
110  }
111 }
Log< level::Info, true > LogVerbatim
void fillTransProf(float u, float v)
T getUntrackedParameter(std::string const &, T const &) const
std::vector< TH1D * > rt_histoProf
Definition: HcalTB02Histo.h:50
Definition: tfile.py:1
void fillAllTime(float v)
HcalTB02Histo(const edm::ParameterSet &ps)
Definition: value.py:1
float getMean(int ilayer)
TH2D * rt_TransProf
Definition: HcalTB02Histo.h:49
void fillProfile(int ilayer, float value)
TH1D * rt_tbTimes
Definition: HcalTB02Histo.h:48
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
virtual ~HcalTB02Histo()
std::string fileName
Definition: HcalTB02Histo.h:45
float getRMS(int ilayer)