#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 |
Description: Histogram handling for Hcal Test Beam 2002 studies
Usage: Sets up histograms and stores in a file
Definition at line 31 of file HcalTB02Histo.h.
HcalTB02Histo::HcalTB02Histo | ( | const edm::ParameterSet & | ps | ) |
Definition at line 26 of file HcalTB02Histo.cc.
References fileName, edm::ParameterSet::getUntrackedParameter(), h, edm::Service< T >::isAvailable(), rt_histoProf, rt_tbTimes, rt_TransProf, indexGen::title, and verbose.
: rt_tbTimes(0), rt_TransProf(0) { verbose = ps.getUntrackedParameter<bool>("Verbose", false); edm::LogInfo("HcalTBSim") << "HcalTB02Histo:: initialised with o/p file " << fileName << " verbosity " << verbose; // Book histograms edm::Service<TFileService> tfile; if ( !tfile.isAvailable() ) throw cms::Exception("BadConfig") << "TFileService unavailable: " << "please add it to config file"; char title[80]; for (int ilayer=0; ilayer<19; ilayer++) { sprintf(title, "Scint. Energy in Layer %d ", ilayer); TH1D *h = tfile->make<TH1D>(title, title, 500, 0., 1.5); rt_histoProf.push_back(h); edm::LogInfo("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " <<title; } sprintf(title, "All Hit Time slices"); rt_tbTimes = tfile->make<TH1D>(title, title, 100, 0., 100.); edm::LogInfo("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title; sprintf(title, "Transv. Shower Profile"); rt_TransProf = tfile->make<TH2D>(title, title, 100, 0., 1., 1000, 0., 0.5); edm::LogInfo("HcalTBSim") << "HcalTB02Histo:: Initialise Histo " << title; }
HcalTB02Histo::~HcalTB02Histo | ( | ) | [virtual] |
Definition at line 55 of file HcalTB02Histo.cc.
{ edm::LogInfo("HcalTBSim") << " Deleting HcalTB02Histo"; }
void HcalTB02Histo::fillAllTime | ( | float | v | ) |
Definition at line 63 of file HcalTB02Histo.cc.
References LogDebug, rt_tbTimes, and v.
Referenced by HcalTB02Analysis::update().
{ LogDebug("HcalTBSim") << "HcalTB02Histo::Fill Time histo with " << v; if (rt_tbTimes) { rt_tbTimes->Fill(v); } }
void HcalTB02Histo::fillProfile | ( | int | ilayer, |
float | value | ||
) |
Definition at line 80 of file HcalTB02Histo.cc.
References LogDebug, rt_histoProf, and relativeConstraints::value.
Referenced by HcalTB02Analysis::update().
{ if (ilayer < (int)(rt_histoProf.size())) { rt_histoProf[ilayer]->Fill(value); LogDebug("HcalTBSim") << "HcalTB02Histo::Fill profile " << ilayer << " with " << value; } }
void HcalTB02Histo::fillTransProf | ( | float | u, |
float | v | ||
) |
Definition at line 71 of file HcalTB02Histo.cc.
References LogDebug, rt_TransProf, and v.
Referenced by HcalTB02Analysis::update().
{ LogDebug("HcalTBSim") << "HcalTB02Histo:::Fill Shower Transv. Profile histo" << " with " << u << " and " << v; if (rt_TransProf) { rt_TransProf->Fill(u,v); } }
float HcalTB02Histo::getMean | ( | int | ilayer | ) |
Definition at line 89 of file HcalTB02Histo.cc.
References rt_histoProf.
{ if (ilayer < (int)(rt_histoProf.size())) { return rt_histoProf[ilayer]->GetMean(); } else { return 0; } }
float HcalTB02Histo::getRMS | ( | int | ilayer | ) |
Definition at line 98 of file HcalTB02Histo.cc.
References rt_histoProf.
{ if (ilayer < (int)(rt_histoProf.size())) { return rt_histoProf[ilayer]->GetRMS(); } else { return 0; } }
std::string HcalTB02Histo::fileName [private] |
Definition at line 49 of file HcalTB02Histo.h.
Referenced by HcalTB02Histo().
std::vector<TH1D *> HcalTB02Histo::rt_histoProf [private] |
Definition at line 54 of file HcalTB02Histo.h.
Referenced by fillProfile(), getMean(), getRMS(), and HcalTB02Histo().
TH1D* HcalTB02Histo::rt_tbTimes [private] |
Definition at line 52 of file HcalTB02Histo.h.
Referenced by fillAllTime(), and HcalTB02Histo().
TH2D* HcalTB02Histo::rt_TransProf [private] |
Definition at line 53 of file HcalTB02Histo.h.
Referenced by fillTransProf(), and HcalTB02Histo().
bool HcalTB02Histo::verbose [private] |
Definition at line 50 of file HcalTB02Histo.h.
Referenced by HcalTB02Histo().