CMS 3D CMS Logo

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

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

Public Member Functions

void fillEdep (double etots, double eecals, double ehcals, double etotq, double eecalq, double ehcalq)
 
void fillLongProf (const std::vector< double > &es, const std::vector< double > &eq)
 
void fillPrimary (double energy, double eta, double phi)
 
void fillTrnsProf (const std::vector< double > &es1, const std::vector< double > &eq1, const std::vector< double > &es2, const std::vector< double > &eq2)
 
 HcalTB04Histo (const edm::ParameterSet &ps)
 
virtual ~HcalTB04Histo ()
 

Private Attributes

TH1D * edecQ
 
TH1D * edecS
 
TH2D * edehQ
 
TH2D * edehS
 
TH1D * edepQ
 
TH1D * edepS
 
TH1D * edhcQ
 
TH1D * edhcS
 
double eHcalMax
 
double eTotMax
 
TH1D * iEta
 
TH1D * iniE
 
TH1D * iPhi
 
TProfile * latqe
 
TProfile * latqf
 
TProfile * latse
 
TProfile * latsf
 
TProfile * lngq
 
TProfile * lngs
 
bool verbose
 

Detailed Description

Description: Histogram handling for Hcal Test Beam 2004 studies

Usage: Sets up histograms and stores in a file

Definition at line 31 of file HcalTB04Histo.h.

Constructor & Destructor Documentation

◆ HcalTB04Histo()

HcalTB04Histo::HcalTB04Histo ( const edm::ParameterSet ps)

Definition at line 27 of file HcalTB04Histo.cc.

References edecQ, edecS, edehQ, edehS, edepQ, edepS, edhcQ, edhcS, edm::ParameterSet::getUntrackedParameter(), iEta, iniE, iPhi, latqe, latqf, latse, latsf, lngq, lngs, compare::tfile, and verbose.

28  : iniE(nullptr),
29  iEta(nullptr),
30  iPhi(nullptr),
31  edepS(nullptr),
32  edecS(nullptr),
33  edhcS(nullptr),
34  edepQ(nullptr),
35  edecQ(nullptr),
36  edhcQ(nullptr),
37  edehS(nullptr),
38  edehQ(nullptr),
39  latse(nullptr),
40  latqe(nullptr),
41  latsf(nullptr),
42  latqf(nullptr),
43  lngs(nullptr),
44  lngq(nullptr) {
45  verbose = ps.getUntrackedParameter<bool>("Verbose", false);
46  double em1 = ps.getUntrackedParameter<double>("ETtotMax", 400.);
47  double em2 = ps.getUntrackedParameter<double>("EHCalMax", 4.0);
48 
49  // Book histograms
51 
52  if (!tfile.isAvailable())
53  throw cms::Exception("BadConfig") << "TFileService unavailable: "
54  << "please add it to config file";
55  iniE = tfile->make<TH1D>("iniE", "Incident Energy (GeV)", 4000, 0., em1);
56  iEta = tfile->make<TH1D>("iEta", "Eta at incidence ", 300, 0., 3.);
57  iPhi = tfile->make<TH1D>("iPhi", "Phi at incidence ", 300, -1., 1.);
58  edepS = tfile->make<TH1D>("edepS", "Energy deposit == Total (Simhit)", 4000, 0., em1);
59  edecS = tfile->make<TH1D>("edecS", "Energy deposit == ECal (Simhit)", 4000, 0., em1);
60  edhcS = tfile->make<TH1D>("edhcS", "Energy deposit == HCal (Simhit)", 4000, 0., em2);
61  edepQ = tfile->make<TH1D>("edepQ", "Energy deposit == Total (QIE)", 4000, 0., em1);
62  edecQ = tfile->make<TH1D>("edecQ", "Energy deposit == ECal (QIE)", 4000, 0., em1);
63  edhcQ = tfile->make<TH1D>("edhcQ", "Energy deposit == HCal (QIE)", 4000, 0., em2);
64  edehS = tfile->make<TH2D>("edehS", "Hcal vs Ecal (Simhit)", 100, 0., em1, 100, 0., em2);
65  edehQ = tfile->make<TH2D>("edehQ", "Hcal vs Ecal (QIE)", 100, 0., em1, 100, 0., em2);
66  latse = tfile->make<TProfile>("latse", "Lat Prof (Eta Sim)", 10, 0., 10.);
67  latqe = tfile->make<TProfile>("latqe", "Lat Prof (Eta QIE)", 10, 0., 10.);
68  latsf = tfile->make<TProfile>("latsf", "Lat Prof (Phi Sim)", 10, 0., 10.);
69  latqf = tfile->make<TProfile>("latqf", "Lat Prof (Phi QIE)", 10, 0., 10.);
70  lngs = tfile->make<TProfile>("lngs", "Long. Prof (Sim)", 20, 0., 20.);
71  lngq = tfile->make<TProfile>("lngq", "Long. Prof (QIE)", 20, 0., 20.);
72 }
TProfile * lngq
Definition: HcalTB04Histo.h:54
TProfile * lngs
Definition: HcalTB04Histo.h:54
T getUntrackedParameter(std::string const &, T const &) const
TProfile * latse
Definition: HcalTB04Histo.h:54
TProfile * latqf
Definition: HcalTB04Histo.h:54
Definition: tfile.py:1
TProfile * latsf
Definition: HcalTB04Histo.h:54
TProfile * latqe
Definition: HcalTB04Histo.h:54

◆ ~HcalTB04Histo()

HcalTB04Histo::~HcalTB04Histo ( )
virtual

Definition at line 74 of file HcalTB04Histo.cc.

74 {}

Member Function Documentation

◆ fillEdep()

void HcalTB04Histo::fillEdep ( double  etots,
double  eecals,
double  ehcals,
double  etotq,
double  eecalq,
double  ehcalq 
)

Definition at line 89 of file HcalTB04Histo.cc.

References edecQ, edecS, edehQ, edehS, edepQ, edepS, edhcQ, and edhcS.

Referenced by HcalTB04Analysis::finalAnalysis().

89  {
90 #ifdef EDM_ML_DEBUG
91  edm::LogVerbatim("HcalTBSim") << "HcalTB04Histo:::fillEdep: Simulated Total " << etots << " ECal " << eecals
92  << " HCal " << ehcals << " Digitised Total " << etotq << " ECal " << eecalq << " HCal "
93  << ehcalq;
94 #endif
95  edepS->Fill(etots);
96  edecS->Fill(eecals);
97  edhcS->Fill(ehcals);
98  edepQ->Fill(etotq);
99  edecQ->Fill(eecalq);
100  edhcQ->Fill(ehcalq);
101  edehS->Fill(eecals, ehcals);
102  edehQ->Fill(eecalq, ehcalq);
103 }
Log< level::Info, true > LogVerbatim

◆ fillLongProf()

void HcalTB04Histo::fillLongProf ( const std::vector< double > &  es,
const std::vector< double > &  eq 
)

Definition at line 135 of file HcalTB04Histo.cc.

References mps_fire::i, lngq, lngs, SiStripPI::min, and dqmiodumpmetadata::n.

Referenced by HcalTB04Analysis::finalAnalysis().

135  {
136 #ifdef EDM_ML_DEBUG
137  unsigned int n = std::min(es.size(), eq.size());
138  for (unsigned int i = 0; i < n; i++)
139  edm::LogVerbatim("HcalTBSim") << "HcalTB04Histo::fillLongProf [" << i << "] Sim " << es[i] << " Dig " << eq[i];
140 #endif
141  for (unsigned int i = 0; i < (es.size()); i++) {
142  double lay = i + 0.5;
143  lngs->Fill(lay, es[i]);
144  }
145  for (unsigned int i = 0; i < (eq.size()); i++) {
146  double lay = i + 0.5;
147  lngq->Fill(lay, eq[i]);
148  }
149 }
TProfile * lngq
Definition: HcalTB04Histo.h:54
TProfile * lngs
Definition: HcalTB04Histo.h:54

◆ fillPrimary()

void HcalTB04Histo::fillPrimary ( double  energy,
double  eta,
double  phi 
)

Definition at line 80 of file HcalTB04Histo.cc.

References HCALHighEnergyHPDFilter_cfi::energy, PVValHelper::eta, iEta, iniE, iPhi, and phi.

Referenced by HcalTB04Analysis::finalAnalysis().

80  {
81 #ifdef EDM_ML_DEBUG
82  edm::LogVerbatim("HcalTBSim") << "HcalTB04Histo::fillPrimary: Energy " << energy << " Eta " << eta << " Phi " << phi;
83 #endif
84  iniE->Fill(energy);
85  iEta->Fill(eta);
86  iPhi->Fill(phi);
87 }
Log< level::Info, true > LogVerbatim

◆ fillTrnsProf()

void HcalTB04Histo::fillTrnsProf ( const std::vector< double > &  es1,
const std::vector< double > &  eq1,
const std::vector< double > &  es2,
const std::vector< double > &  eq2 
)

Definition at line 105 of file HcalTB04Histo.cc.

References mps_fire::i, latqe, latqf, latse, latsf, SiStripPI::min, and dqmiodumpmetadata::n.

Referenced by HcalTB04Analysis::finalAnalysis().

108  {
109 #ifdef EDM_ML_DEBUG
110  unsigned int n1 = std::min(es1.size(), eq1.size());
111  unsigned int n2 = std::min(es2.size(), eq2.size());
112  unsigned int n = std::min(n1, n2);
113  for (unsigned int i = 0; i < n; i++)
114  edm::LogVerbatim("HcalTBSim") << "HcalTB04Histo::fillTrnsProf [" << i << "] SimEta " << es1[i] << " DigEta "
115  << eq1[i] << " SimPhi " << es2[i] << " DigPhi " << eq2[i];
116 #endif
117  for (unsigned int i = 0; i < (es1.size()); i++) {
118  double tow = i + 0.5;
119  latse->Fill(tow, es1[i]);
120  }
121  for (unsigned int i = 0; i < (eq1.size()); i++) {
122  double tow = i + 0.5;
123  latqe->Fill(tow, eq1[i]);
124  }
125  for (unsigned int i = 0; i < (es2.size()); i++) {
126  double tow = i + 0.5;
127  latsf->Fill(tow, es2[i]);
128  }
129  for (unsigned int i = 0; i < (eq2.size()); i++) {
130  double tow = i + 0.5;
131  latqf->Fill(tow, eq2[i]);
132  }
133 }
TProfile * latse
Definition: HcalTB04Histo.h:54
TProfile * latqf
Definition: HcalTB04Histo.h:54
TProfile * latsf
Definition: HcalTB04Histo.h:54
TProfile * latqe
Definition: HcalTB04Histo.h:54

Member Data Documentation

◆ edecQ

TH1D * HcalTB04Histo::edecQ
private

Definition at line 52 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edecS

TH1D * HcalTB04Histo::edecS
private

Definition at line 52 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edehQ

TH2D * HcalTB04Histo::edehQ
private

Definition at line 53 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edehS

TH2D* HcalTB04Histo::edehS
private

Definition at line 53 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edepQ

TH1D * HcalTB04Histo::edepQ
private

Definition at line 52 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edepS

TH1D* HcalTB04Histo::edepS
private

Definition at line 52 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edhcQ

TH1D * HcalTB04Histo::edhcQ
private

Definition at line 52 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ edhcS

TH1D * HcalTB04Histo::edhcS
private

Definition at line 52 of file HcalTB04Histo.h.

Referenced by fillEdep(), and HcalTB04Histo().

◆ eHcalMax

double HcalTB04Histo::eHcalMax
private

Definition at line 49 of file HcalTB04Histo.h.

◆ eTotMax

double HcalTB04Histo::eTotMax
private

Definition at line 49 of file HcalTB04Histo.h.

◆ iEta

TH1D * HcalTB04Histo::iEta
private

Definition at line 51 of file HcalTB04Histo.h.

Referenced by fillPrimary(), and HcalTB04Histo().

◆ iniE

TH1D* HcalTB04Histo::iniE
private

Definition at line 51 of file HcalTB04Histo.h.

Referenced by fillPrimary(), and HcalTB04Histo().

◆ iPhi

TH1D * HcalTB04Histo::iPhi
private

Definition at line 51 of file HcalTB04Histo.h.

Referenced by fillPrimary(), and HcalTB04Histo().

◆ latqe

TProfile * HcalTB04Histo::latqe
private

Definition at line 54 of file HcalTB04Histo.h.

Referenced by fillTrnsProf(), and HcalTB04Histo().

◆ latqf

TProfile * HcalTB04Histo::latqf
private

Definition at line 54 of file HcalTB04Histo.h.

Referenced by fillTrnsProf(), and HcalTB04Histo().

◆ latse

TProfile* HcalTB04Histo::latse
private

Definition at line 54 of file HcalTB04Histo.h.

Referenced by fillTrnsProf(), and HcalTB04Histo().

◆ latsf

TProfile * HcalTB04Histo::latsf
private

Definition at line 54 of file HcalTB04Histo.h.

Referenced by fillTrnsProf(), and HcalTB04Histo().

◆ lngq

TProfile * HcalTB04Histo::lngq
private

Definition at line 54 of file HcalTB04Histo.h.

Referenced by fillLongProf(), and HcalTB04Histo().

◆ lngs

TProfile * HcalTB04Histo::lngs
private

Definition at line 54 of file HcalTB04Histo.h.

Referenced by fillLongProf(), and HcalTB04Histo().

◆ verbose

bool HcalTB04Histo::verbose
private

Definition at line 48 of file HcalTB04Histo.h.

Referenced by HcalTB04Histo().