CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HcalTB02Histo Class Reference

#include <SimG4CMS/HcalTestBeam/interface/HcalTB02Histo.h>

Public Member Functions

void fillAllTime (float v)
 
void fillProfile (int ilayer, float value)
 
void fillTransProf (float u, float v)
 
float getMean (int ilayer)
 
float getRMS (int ilayer)
 
 HcalTB02Histo (const edm::ParameterSet &ps)
 
virtual ~HcalTB02Histo ()
 

Private Attributes

std::string fileName
 
std::vector< TH1D * > rt_histoProf
 
TH1D * rt_tbTimes
 
TH2D * rt_TransProf
 
bool verbose
 

Detailed Description

Description: Histogram handling for Hcal Test Beam 2002 studies

Usage: Sets up histograms and stores in a file

Definition at line 30 of file HcalTB02Histo.h.

Constructor & Destructor Documentation

◆ HcalTB02Histo()

HcalTB02Histo::HcalTB02Histo ( const edm::ParameterSet ps)

Definition at line 26 of file HcalTB02Histo.cc.

References fileName, edm::ParameterSet::getUntrackedParameter(), h, rt_histoProf, rt_tbTimes, rt_TransProf, compare::tfile, runGCPTkAlMap::title, and verbose.

26  : 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 }
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
std::vector< TH1D * > rt_histoProf
Definition: HcalTB02Histo.h:50
Definition: tfile.py:1
TH2D * rt_TransProf
Definition: HcalTB02Histo.h:49
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
std::string fileName
Definition: HcalTB02Histo.h:45

◆ ~HcalTB02Histo()

HcalTB02Histo::~HcalTB02Histo ( )
virtual

Definition at line 59 of file HcalTB02Histo.cc.

59  {
60 #ifdef EDM_ML_DEBUG
61  edm::LogVerbatim("HcalTBSim") << " Deleting HcalTB02Histo";
62 #endif
63 }
Log< level::Info, true > LogVerbatim

Member Function Documentation

◆ fillAllTime()

void HcalTB02Histo::fillAllTime ( float  v)

Definition at line 69 of file HcalTB02Histo.cc.

References rt_tbTimes, and findQualityFiles::v.

69  {
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 }
Log< level::Info, true > LogVerbatim
TH1D * rt_tbTimes
Definition: HcalTB02Histo.h:48

◆ fillProfile()

void HcalTB02Histo::fillProfile ( int  ilayer,
float  value 
)

Definition at line 88 of file HcalTB02Histo.cc.

References rt_histoProf, and relativeConstraints::value.

88  {
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 }
Log< level::Info, true > LogVerbatim
std::vector< TH1D * > rt_histoProf
Definition: HcalTB02Histo.h:50
Definition: value.py:1

◆ fillTransProf()

void HcalTB02Histo::fillTransProf ( float  u,
float  v 
)

Definition at line 78 of file HcalTB02Histo.cc.

References rt_TransProf, and findQualityFiles::v.

78  {
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 }
Log< level::Info, true > LogVerbatim
TH2D * rt_TransProf
Definition: HcalTB02Histo.h:49

◆ getMean()

float HcalTB02Histo::getMean ( int  ilayer)

Definition at line 97 of file HcalTB02Histo.cc.

References rt_histoProf.

97  {
98  if (ilayer < (int)(rt_histoProf.size())) {
99  return rt_histoProf[ilayer]->GetMean();
100  } else {
101  return 0;
102  }
103 }
std::vector< TH1D * > rt_histoProf
Definition: HcalTB02Histo.h:50

◆ getRMS()

float HcalTB02Histo::getRMS ( int  ilayer)

Definition at line 105 of file HcalTB02Histo.cc.

References rt_histoProf.

105  {
106  if (ilayer < (int)(rt_histoProf.size())) {
107  return rt_histoProf[ilayer]->GetRMS();
108  } else {
109  return 0;
110  }
111 }
std::vector< TH1D * > rt_histoProf
Definition: HcalTB02Histo.h:50

Member Data Documentation

◆ fileName

std::string HcalTB02Histo::fileName
private

Definition at line 45 of file HcalTB02Histo.h.

Referenced by HcalTB02Histo().

◆ rt_histoProf

std::vector<TH1D *> HcalTB02Histo::rt_histoProf
private

Definition at line 50 of file HcalTB02Histo.h.

Referenced by fillProfile(), getMean(), getRMS(), and HcalTB02Histo().

◆ rt_tbTimes

TH1D* HcalTB02Histo::rt_tbTimes
private

Definition at line 48 of file HcalTB02Histo.h.

Referenced by fillAllTime(), and HcalTB02Histo().

◆ rt_TransProf

TH2D* HcalTB02Histo::rt_TransProf
private

Definition at line 49 of file HcalTB02Histo.h.

Referenced by fillTransProf(), and HcalTB02Histo().

◆ verbose

bool HcalTB02Histo::verbose
private

Definition at line 46 of file HcalTB02Histo.h.

Referenced by HcalTB02Histo().