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