CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: initialised with o/p file " << fileName << " verbosity " << verbose;
29 
30  // Book histograms
32 
33  if (!tfile.isAvailable())
34  throw cms::Exception("BadConfig") << "TFileService unavailable: "
35  << "please add it to config file";
36 
37  char title[80];
38  for (int ilayer = 0; ilayer < 19; ilayer++) {
39  sprintf(title, "Scint. Energy in Layer %d ", ilayer);
40  TH1D* h = tfile->make<TH1D>(title, title, 500, 0., 1.5);
41  rt_histoProf.push_back(h);
42 #ifdef EDM_ML_DEBUG
43  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title;
44 #endif
45  }
46  sprintf(title, "All Hit Time slices");
47  rt_tbTimes = tfile->make<TH1D>(title, title, 100, 0., 100.);
48 #ifdef EDM_ML_DEBUG
49  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title;
50 #endif
51  sprintf(title, "Transv. Shower Profile");
52  rt_TransProf = tfile->make<TH2D>(title, title, 100, 0., 1., 1000, 0., 0.5);
53 #ifdef EDM_ML_DEBUG
54  edm::LogVerbatim("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title;
55 #endif
56 }
57 
59 #ifdef EDM_ML_DEBUG
60  edm::LogVerbatim("HcalTBSim") << " Deleting HcalTB02Histo";
61 #endif
62 }
63 
64 //
65 // member functions
66 //
67 
69  LogDebug("HcalTBSim") << "HcalTB02Histo::Fill Time histo with " << v;
70  if (rt_tbTimes) {
71  rt_tbTimes->Fill(v);
72  }
73 }
74 
75 void HcalTB02Histo::fillTransProf(float u, float v) {
76  LogDebug("HcalTBSim") << "HcalTB02Histo:::Fill Shower Transv. Profile histo"
77  << " with " << u << " and " << v;
78  if (rt_TransProf) {
79  rt_TransProf->Fill(u, v);
80  }
81 }
82 
83 void HcalTB02Histo::fillProfile(int ilayer, float value) {
84  if (ilayer < (int)(rt_histoProf.size())) {
85  rt_histoProf[ilayer]->Fill(value);
86  LogDebug("HcalTBSim") << "HcalTB02Histo::Fill profile " << ilayer << " with " << value;
87  }
88 }
89 
90 float HcalTB02Histo::getMean(int ilayer) {
91  if (ilayer < (int)(rt_histoProf.size())) {
92  return rt_histoProf[ilayer]->GetMean();
93  } else {
94  return 0;
95  }
96 }
97 
98 float HcalTB02Histo::getRMS(int ilayer) {
99  if (ilayer < (int)(rt_histoProf.size())) {
100  return rt_histoProf[ilayer]->GetRMS();
101  } else {
102  return 0;
103  }
104 }
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
void fillTransProf(float u, float v)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
std::vector< TH1D * > rt_histoProf
Definition: HcalTB02Histo.h:50
void fillAllTime(float v)
bool isAvailable() const
Definition: Service.h:40
HcalTB02Histo(const edm::ParameterSet &ps)
float getMean(int ilayer)
TH2D * rt_TransProf
Definition: HcalTB02Histo.h:49
void fillProfile(int ilayer, float value)
TH1D * rt_tbTimes
Definition: HcalTB02Histo.h:48
tuple tfile
Definition: compare.py:324
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)
#define LogDebug(id)