CMS 3D CMS Logo

DigiInvestigatorHistogramMaker.cc
Go to the documentation of this file.
7 #include "TProfile.h"
8 #include "TH1F.h"
9 
11 
13  : _hitname(),
14  _nbins(500),
15  m_maxLS(100),
16  m_LSfrac(4),
17  _scalefact(),
18  _runHisto(true),
19  _fillHisto(false),
20  _binmax(),
21  _labels(),
22  _rhm(iC),
23  _fhm(iC, true),
24  _nmultvsorbrun(),
25  _nmultvsbxrun(),
26  _nmultvsbxfill(),
27  _nmult() {}
28 
31  : _hitname(iConfig.getUntrackedParameter<std::string>("hitName", "digi")),
32  _nbins(iConfig.getUntrackedParameter<int>("numberOfBins", 500)),
33  m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin", 100)),
34  m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction", 4)),
35  _scalefact(iConfig.getUntrackedParameter<int>("scaleFactor", 5)),
36  _runHisto(iConfig.getUntrackedParameter<bool>("runHisto", true)),
37  _fillHisto(iConfig.getUntrackedParameter<bool>("fillHisto", false)),
38  _labels(),
39  _rhm(iC),
40  _fhm(iC, true),
41  _nmultvsorbrun(),
42  _nmultvsbxrun(),
43  _nmultvsbxfill(),
44  _nmult(),
45  _subdirs() {
46  std::vector<edm::ParameterSet> wantedsubds(iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >(
47  "wantedSubDets", std::vector<edm::ParameterSet>()));
48 
49  for (std::vector<edm::ParameterSet>::iterator ps = wantedsubds.begin(); ps != wantedsubds.end(); ++ps) {
50  _labels[ps->getParameter<unsigned int>("detSelection")] = ps->getParameter<std::string>("detLabel");
51  _binmax[ps->getParameter<unsigned int>("detSelection")] = ps->getParameter<int>("binMax");
52  }
53 }
54 
56  for (std::map<unsigned int, std::string>::const_iterator lab = _labels.begin(); lab != _labels.end(); lab++) {
57  const unsigned int i = lab->first;
58  const std::string slab = lab->second;
59 
60  delete _subdirs[i];
61  }
62 }
63 
65  const std::map<unsigned int, std::string>& labels) {
66  _labels = labels;
67  book(dirname);
68 }
69 
72  TFileDirectory subev = tfserv->mkdir(dirname);
73 
74  SiStripTKNumbers trnumb;
75 
76  edm::LogInfo("NumberOfBins") << "Number of Bins: " << _nbins;
77  edm::LogInfo("NumberOfMaxLS") << "Max number of LS before rebinning: " << m_maxLS;
78  edm::LogInfo("StartingLSFrac") << "Fraction of LS in one bin before rebinning: " << m_LSfrac;
79  edm::LogInfo("ScaleFactors") << "x-axis range scale factor: " << _scalefact;
80  edm::LogInfo("BinMaxValue") << "Setting bin max values";
81 
82  for (std::map<unsigned int, std::string>::const_iterator lab = _labels.begin(); lab != _labels.end(); lab++) {
83  const unsigned int i = lab->first;
84  const std::string slab = lab->second;
85 
86  if (_binmax.find(i) == _binmax.end()) {
87  edm::LogVerbatim("NotConfiguredBinMax")
88  << "Bin max for " << lab->second << " not configured: " << trnumb.nstrips(i) << " used";
89  _binmax[i] = trnumb.nstrips(i);
90  }
91 
92  edm::LogVerbatim("BinMaxValue") << "Bin max for " << lab->second << " is " << _binmax[i];
93  }
94 
95  for (std::map<unsigned int, std::string>::const_iterator lab = _labels.begin(); lab != _labels.end(); ++lab) {
96  const int i = lab->first;
97  const std::string slab = lab->second;
98 
99  char name[200];
100  char title[500];
101 
102  _subdirs[i] = new TFileDirectory(subev.mkdir(slab));
103 
104  if (_subdirs[i]) {
105  sprintf(name, "n%sdigi", slab.c_str());
106  sprintf(title, "%s %s multiplicity", slab.c_str(), _hitname.c_str());
107  _nmult[i] = _subdirs[i]->make<TH1F>(name, title, _nbins, 0., (1 + _binmax[i] / (_scalefact * _nbins)) * _nbins);
108  _nmult[i]->GetXaxis()->SetTitle("Number of Hits");
109  _nmult[i]->GetYaxis()->SetTitle("Events");
110 
111  if (_runHisto) {
112  sprintf(name, "n%sdigivsorbrun", slab.c_str());
113  sprintf(title, "%s %s mean multiplicity vs orbit", slab.c_str(), _hitname.c_str());
115  sprintf(name, "n%sdigivsbxrun", slab.c_str());
116  sprintf(title, "%s %s mean multiplicity vs BX", slab.c_str(), _hitname.c_str());
117  _nmultvsbxrun[i] = _rhm.makeTProfile(name, title, 3564, -0.5, 3563.5);
118  }
119  if (_fillHisto) {
120  sprintf(name, "n%sdigivsbxfill", slab.c_str());
121  sprintf(title, "%s %s mean multiplicity vs BX", slab.c_str(), _hitname.c_str());
122  _nmultvsbxfill[i] = _fhm.makeTProfile(name, title, 3564, -0.5, 3563.5);
123  }
124  }
125  }
126 }
127 
129  // char runname[100];
130  // sprintf(runname,"run_%d",nrun);
131 
133 
134  // currdir = &(*tfserv);
135  // _rhm.beginRun(nrun,*currdir);
136 
137  _rhm.beginRun(iRun, tfserv->tFileDirectory());
138  _fhm.beginRun(iRun, tfserv->tFileDirectory());
139 
140  for (std::map<unsigned int, std::string>::const_iterator lab = _labels.begin(); lab != _labels.end(); ++lab) {
141  const int i = lab->first;
142  const std::string slab = lab->second;
143 
144  // char name[200];
145  // char title[500];
146 
147  // TFileDirectory subd =_subdirs[i]->mkdir(runname);
148 
149  // sprintf(name,"n%sdigivsorbrun",slab.c_str());
150  // sprintf(title,"%s %s mean multiplicity vs orbit",slab.c_str(),_hitname.c_str());
151  // _nmultvsorbrun[i] = subd.make<TProfile>(name,title,_norbbin,0.5,11223*_norbbin+0.5);
152  if (_runHisto) {
153  if (*_nmultvsorbrun[i]) {
154  (*_nmultvsorbrun[i])->GetXaxis()->SetTitle("time [orbit#]");
155  (*_nmultvsorbrun[i])->GetYaxis()->SetTitle("Hits");
156  (*_nmultvsorbrun[i])->SetCanExtend(TH1::kXaxis);
157  }
158  if (*_nmultvsbxrun[i]) {
159  (*_nmultvsbxrun[i])->GetXaxis()->SetTitle("BX#");
160  (*_nmultvsbxrun[i])->GetYaxis()->SetTitle("Mean Number of Hits");
161  }
162  }
163  if (_fillHisto) {
164  if (*_nmultvsbxfill[i]) {
165  (*_nmultvsbxfill[i])->GetXaxis()->SetTitle("BX#");
166  (*_nmultvsbxfill[i])->GetYaxis()->SetTitle("Mean Number of Hits");
167  }
168  }
169  }
170 }
171 
172 void DigiInvestigatorHistogramMaker::fill(const edm::Event& iEvent, const std::map<unsigned int, int>& ndigi) {
173  for (std::map<unsigned int, int>::const_iterator digi = ndigi.begin(); digi != ndigi.end(); digi++) {
174  if (_labels.find(digi->first) != _labels.end()) {
175  const unsigned int i = digi->first;
176 
177  _nmult[i]->Fill(digi->second);
178  if (_runHisto) {
179  if (_nmultvsorbrun[i] && *_nmultvsorbrun[i])
180  (*_nmultvsorbrun[i])->Fill(iEvent.orbitNumber(), digi->second);
181  if (_nmultvsbxrun[i] && *_nmultvsbxrun[i])
182  (*_nmultvsbxrun[i])->Fill(iEvent.bunchCrossing() % 3564, digi->second);
183  }
184  if (_fillHisto) {
185  if (_nmultvsbxfill[i] && *_nmultvsbxfill[i])
186  (*_nmultvsbxfill[i])->Fill(iEvent.bunchCrossing() % 3564, digi->second);
187  }
188  }
189  }
190 }
Log< level::Info, true > LogVerbatim
int nstrips(const DetId &detid) const
void book(const std::string dirname, const std::map< unsigned int, std::string > &labels)
std::map< unsigned int, TProfile ** > _nmultvsbxrun
std::map< unsigned int, std::string > _labels
T getUntrackedParameter(std::string const &, T const &) const
std::map< unsigned int, TProfile ** > _nmultvsbxfill
TFileDirectory & tFileDirectory()
Definition: TFileService.h:42
int iEvent
Definition: GenABIO.cc:224
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
TProfile ** makeTProfile(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
std::map< unsigned int, TProfile ** > _nmultvsorbrun
Log< level::Info, false > LogInfo
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
void beginRun(const edm::Run &iRun)
void fill(const edm::Event &iEvent, const std::map< unsigned int, int > &ndigi)
std::map< unsigned int, TFileDirectory * > _subdirs
std::map< unsigned int, TH1F * > _nmult
DigiInvestigatorHistogramMaker(edm::ConsumesCollector &&iC)
Definition: Run.h:45