ATLAS Offline Software
MultiHisto.icc
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef MultiHisto_cxx
6 #define MultiHisto_cxx
7 
8 #include "PixelCalibAlgs/MultiHisto.h"
9 #include "CxxUtils/checker_macros.h"
10 #include "sstream"
11 #include <vector>
12 #include <string>
13 #include <TROOT.h>
14 #include <TCanvas.h>
15 
16 ////////////////////////////////////////////////////////////////////////////////
17 //
18 // Building up & filling & destroing the thing
19 //
20 ////////////////////////////////////////////////////////////////////////////////
21 
22 // creating the MultiHisto
23 template <class ht>
24 MultiHisto<ht>::MultiHisto(const ht &histo_model,const std::vector<std::string> &nameDiv,
25  const std::vector<int> &nDiv, const std::vector<double *> &binsDiv):
26  m_div_Nbins(nDiv), m_div_names(nameDiv){
27 
28 // this should be added!!! if( ! histo_model->InheritsFrom(TH1F::Class()) ) return -1;
29 
30  m_model = new ht(histo_model);
31  int nhistos = 1;
32  for(unsigned int i=0; i < m_div_Nbins.size(); i++){
33  nhistos*=m_div_Nbins[i];
34  double *tmp = new double[m_div_Nbins[i]+1];
35  for(int j=0; j < m_div_Nbins[i]+1; j++) tmp[j] = (binsDiv[i])[j];
36  m_div_bins.push_back(tmp);
37  }
38  std::vector<int> indexes;
39  for(int i = 0; i < nhistos; i++){
40  ht *clone = (ht *)m_model->Clone();
41  m_myhistos.push_back(clone);
42  indexes = GetDivisionsIndexes(i);
43  m_myhistos[i]->SetTitle( ComposeTitle(indexes).c_str() );
44  m_myhistos[i]->SetName( ComposeName(indexes).c_str() );
45  }
46 }
47 
48 template <class ht>
49 MultiHisto<ht>::MultiHisto(const ht &histo_model,const std::vector<std::string> &nameDiv,
50  const std::vector<std::vector <float> > &binsDiv):
51  m_div_names(nameDiv){
52 
53  int nBinnages = m_div_names.size();
54  for(int i = 0; i < nBinnages ; i++)
55  m_div_Nbins.push_back((binsDiv[i]).size()-1);
56 
57 // this should be added!!! if( ! histo_model->InheritsFrom(TH1F::Class()) ) return -1;
58 
59  m_model = new ht(histo_model);
60  int nhistos = 1;
61  for(unsigned int i=0; i < m_div_Nbins.size(); i++){
62  nhistos*=m_div_Nbins[i];
63  double *tmp = new double[m_div_Nbins[i]+1];
64  for(int j=0; j < m_div_Nbins[i]+1; j++) tmp[j] = (binsDiv[i])[j];
65  m_div_bins.push_back(tmp);
66  }
67  std::vector<int> indexes;
68  for(int i = 0; i < nhistos; i++){
69  ht *clone = (ht *)m_model->Clone();
70  m_myhistos.push_back(clone);
71  indexes = GetDivisionsIndexes(i);
72  m_myhistos[i]->SetTitle( ComposeTitle(indexes).c_str() );
73  m_myhistos[i]->SetName( ComposeName(indexes).c_str() );
74  }
75 }
76 
77 
78 // deleting the thing
79 template <class ht>
80 MultiHisto<ht>::~MultiHisto(){
81 }
82 
83 // filling the histograms
84 template <class ht>
85 int MultiHisto<ht>::Fill(double xvar, double yvar,const std::vector<double> &pars){
86  int index = GetGlobalIndex(pars);
87  if( index >= 0 ) m_myhistos[index]->Fill(xvar,yvar);
88  return index;
89 }
90 
91 template <class ht>
92 int MultiHisto<ht>::Fill(double xvar, double yvar,const std::vector<float> &pars){
93  int index = GetGlobalIndex(pars);
94  if( index >= 0 ) m_myhistos[index]->Fill(xvar,yvar);
95  return index;
96 }
97 
98 // reading the histograms from a file
99 template <class ht>
100 int MultiHisto<ht>::FillFromFile(TDirectory *histodir){
101 
102  int nhistos = m_myhistos.size();
103  int readhistos = 0;
104  TDirectory *olddir = gDirectory;
105  if(histodir) histodir->cd();
106  for(int i=0; i < nhistos; i++){
107  std::vector <int> indexes = GetDivisionsIndexes(i);
108  std::string name = ComposeName(indexes);
109  if( ( m_myhistos[i] = (ht *) gDirectory->Get(name.c_str()) ) ){
110  m_myhistos[i]->SetDirectory(gROOT);
111  readhistos++;
112  }
113  }
114  olddir->cd();
115  return readhistos;
116 }
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 //
120 // Writing the histos
121 //
122 ////////////////////////////////////////////////////////////////////////////////
123 
124 // without counting bins
125 template <class ht>
126 int MultiHisto<ht>::Write(const char* name, Int_t option, Int_t bufsize) const {
127 
128  int nhistos = m_myhistos.size();
129  int writtenhistos = -1;
130  TDirectory *olddir = gDirectory;
131  if(name == 0) name = m_model->GetName();
132  TDirectory *histodir = gDirectory->mkdir(name);
133  if(histodir){
134  histodir->cd();
135  writtenhistos++;
136  for(int i=0; i < nhistos; i++)
137  if(m_myhistos[i]->Write(0,option,bufsize)) writtenhistos++;
138  olddir->cd();
139  }
140  return writtenhistos;
141 }
142 
143 template <class ht>
144 int MultiHisto<ht>::Write (const char* name, Int_t option, Int_t bufsize) {
145  return const_cast<const MultiHisto<ht>*>(this)->Write(name,option,bufsize);
146 }
147 
148 
149 // counting bins
150 template <class ht>
151 int MultiHisto<ht>::WriteAndCount(int *nhistos, int *nbins) const{
152  int nwritten = Write();
153  *nhistos += nwritten;
154  *nbins += nwritten*m_model->GetNbinsY()*m_model->GetNbinsX();
155  return nwritten;
156 }
157 
158 ////////////////////////////////////////////////////////////////////////////////
159 //
160 // Drawing the histos
161 //
162 ////////////////////////////////////////////////////////////////////////////////
163 
164 
165 // draw without canvas and without color
166 template <class ht>
167 void MultiHisto<ht>::Draw(const Option_t* option){
168 
169  std::vector<TCanvas *> c1;
170  Draw(1,option);
171 }
172 
173 // draw without canvas
174 template <class ht>
175 void MultiHisto<ht>::Draw(const int color,const Option_t* option){
176 
177  std::vector<TCanvas *> c1;
178  Draw(c1,color,option);
179 }
180 
181 
182 // draw with canvases
183 template <class ht>
184 void MultiHisto<ht>::Draw(std::vector <TCanvas*> &c1, const int color,const Option_t* option){
185 
186 
187  while(c1.size() < GetNhistos()){
188  TCanvas *c = new TCanvas();
189  c1.push_back(c);
190  }
191 
192  for(unsigned int i = 0 ; i < GetNhistos() ; i++){
193  ht* swap = GetHisto(i);
194  c1[i]->SetName(swap->GetName());
195  c1[i]->Print();
196  swap->SetLineColor(color);
197  swap->Draw(option);
198  }
199 }
200 
201 
202 
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 //
206 // Getting global index
207 //
208 ////////////////////////////////////////////////////////////////////////////////
209 
210 // from a detailed index set
211 template <class ht>
212 int MultiHisto<ht>::GetGlobalIndex(const std::vector<int> &indexes) const{
213  int globalindex = 0, factor = 1;
214  for(int i = indexes.size()-1; i >= 0; i--){
215  globalindex += (indexes[i])*factor;
216  factor *= m_div_Nbins[i];
217  }
218  if(globalindex >= int(m_myhistos.size())) globalindex = -1;
219  return globalindex;
220 }
221 
222 // from a set of parametrs
223 template <class ht>
224 int MultiHisto<ht>::GetGlobalIndex(const std::vector<double> &pars) const{
225  std::vector<int> indexes = GetDivisionsIndexes(pars);
226  return GetGlobalIndex(indexes);
227 }
228 
229 template <class ht>
230 int MultiHisto<ht>::GetGlobalIndex(const std::vector<float> &pars) const{
231  std::vector<int> indexes = GetDivisionsIndexes(pars);
232  return GetGlobalIndex(indexes);
233 }
234 
235 ////////////////////////////////////////////////////////////////////////////////
236 //
237 // Getting detailed indexes
238 //
239 ////////////////////////////////////////////////////////////////////////////////
240 
241 // from a global index
242 template <class ht>
243 std::vector<int> MultiHisto<ht>::GetDivisionsIndexes(int globalindex) const{
244  int curr_index = globalindex, offset =0;
245  std::vector<int> indx, indxfinal;
246  for(int iNdivs = m_div_Nbins.size()-1; iNdivs != 0; iNdivs--){
247  offset = curr_index%m_div_Nbins[iNdivs];
248  indx.push_back(offset);
249  curr_index = curr_index/m_div_Nbins[iNdivs];
250  }
251  indxfinal.push_back(curr_index);
252  for(int i = indx.size()-1; i >= 0; i--) indxfinal.push_back(indx[i]);
253  return indxfinal;
254 }
255 
256 // from a set of parameters
257 template <class ht>
258 std::vector<int> MultiHisto<ht>::GetDivisionsIndexes(const std::vector<double> &pars) const{
259  std::vector<int> indexes;
260  for(unsigned int i = 0; i < m_div_Nbins.size(); i++)
261  for(int j = 0; j < m_div_Nbins[i]; j++)
262  if(pars[i] >= (m_div_bins[i])[j]
263  && pars[i] < (m_div_bins[i])[j+1]){
264  indexes.push_back(int(j));
265  break;
266  }
267  unsigned int size = indexes.size();
268  if(size < m_div_Nbins.size()){
269  for(unsigned int i = 0; i < size; i++) indexes.pop_back();
270  indexes.push_back(-1);
271  }
272  return indexes;
273 }
274 
275 template <class ht>
276 std::vector<int> MultiHisto<ht>::GetDivisionsIndexes(const std::vector<float> &pars) const{
277 
278  std::vector<int> indexes;
279  for(unsigned int i = 0; i < m_div_Nbins.size(); i++)
280  for(int j = 0; j < m_div_Nbins[i]; j++)
281  if(pars[i] >= (m_div_bins[i])[j]
282  && pars[i] < (m_div_bins[i])[j+1]){
283  indexes.push_back(int(j));
284  break;
285  }
286  unsigned int size = indexes.size();
287  if(size < m_div_Nbins.size()){
288  for(unsigned int i = 0; i < size; i++) indexes.pop_back();
289  indexes.push_back(-1);
290  }
291  return indexes;
292 
293 }
294 
295 ////////////////////////////////////////////////////////////////////////////////
296 //
297 // Composing titles and names
298 //
299 ////////////////////////////////////////////////////////////////////////////////
300 
301 template <class ht>
302 std::string MultiHisto<ht>::NameString(int idiv, int index){
303  std::ostringstream basename;
304  std::string name = m_div_names[idiv];
305  size_t found = name.find('#');
306  if(found != std::string::npos) name = name.substr(found+1);
307  found = name.find(' ');
308  if(found != std::string::npos) name.resize(0,found);
309  basename << name << (m_div_bins[idiv])[index] << "-" <<
310  (m_div_bins[idiv])[index+1];
311  return basename.str();
312 }
313 
314 
315 
316 template <class ht>
317 std::string MultiHisto<ht>::ComposeName(const std::vector<int> &indexes){
318  std::string name = std::string(m_model->GetName());
319  for(unsigned int i = 0; i < m_div_Nbins.size(); i++)
320  if(indexes[i] >= 0 && indexes[i] < m_div_Nbins[i])
321  name += "_" + NameString(i,indexes[i]);
322  return name;
323 }
324 
325 ///////////////////////////////////////////////////////////////////////////////
326 
327 template <class ht>
328 std::string MultiHisto<ht>::TitleString(int idiv, int index){
329  std::ostringstream basename;
330  std::string div_name = m_div_names[idiv];
331  double low = (m_div_bins[idiv])[index];
332  double high = (m_div_bins[idiv])[index+1];
333  if( (high - low) == 1 && div_name.find("ayer") != std::string::npos)
334  basename << div_name << " " <<ceil( low );
335  else if( (high - low) == 1 &&
336  div_name.find("uster") != std::string::npos &&
337  div_name.find("ize") != std::string::npos)
338  basename << div_name << " = " << ceil( low );
339  else
340  basename << low << " #leq " << div_name << " < " << high;
341  return basename.str();
342 }
343 
344 
345 template <class ht>
346 std::string MultiHisto<ht>::ComposeTitle(const std::vector<int> &indexes){
347  std::string title = std::string(m_model->GetTitle());
348  for(unsigned int i = 0; i < m_div_Nbins.size(); i++)
349  if(indexes[i] >= 0 && indexes[i] < m_div_Nbins[i])
350  title += " - " + TitleString(i,indexes[i]);
351  return title;
352 }
353 
354 
355 ////////////////////////////////////////////////////////////////////////////////
356 //
357 // Accessing Data
358 //
359 ////////////////////////////////////////////////////////////////////////////////
360 
361 template <class ht>
362 const char * MultiHisto<ht>::GetName() const{
363  return m_model->GetName();
364 }
365 
366 template <class ht>
367 const char * MultiHisto<ht>::GetTitle() const{
368  return m_model->GetTitle();
369 }
370 
371 
372 template <class ht>
373 unsigned int MultiHisto<ht>::GetNhistos() const{
374  return m_myhistos.size();
375 }
376 
377 
378 template <class ht>
379 ht* MultiHisto<ht>::GetHisto(const std::vector<int> &indexes){
380  int i = GetGlobalIndex(indexes);
381  return GetHisto(i);
382 }
383 
384 template <class ht>
385 ht* MultiHisto<ht>::GetHisto(int i){
386  ht *current = NULL;
387  if(i >= 0 && i < int(m_myhistos.size())) current = m_myhistos[i];
388  return current;
389 }
390 
391 
392 
393 #endif