2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
8 #include "PixelCalibAlgs/MultiHisto.h"
9 #include "CxxUtils/checker_macros.h"
16 ////////////////////////////////////////////////////////////////////////////////
18 // Building up & filling & destroing the thing
20 ////////////////////////////////////////////////////////////////////////////////
22 // creating the MultiHisto
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){
28 // this should be added!!! if( ! histo_model->InheritsFrom(TH1F::Class()) ) return -1;
30 m_model = new ht(histo_model);
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);
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() );
49 MultiHisto<ht>::MultiHisto(const ht &histo_model,const std::vector<std::string> &nameDiv,
50 const std::vector<std::vector <float> > &binsDiv):
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);
57 // this should be added!!! if( ! histo_model->InheritsFrom(TH1F::Class()) ) return -1;
59 m_model = new ht(histo_model);
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);
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() );
80 MultiHisto<ht>::~MultiHisto(){
83 // filling the histograms
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);
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);
98 // reading the histograms from a file
100 int MultiHisto<ht>::FillFromFile(TDirectory *histodir){
102 int nhistos = m_myhistos.size();
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);
118 ////////////////////////////////////////////////////////////////////////////////
120 // Writing the histos
122 ////////////////////////////////////////////////////////////////////////////////
124 // without counting bins
126 int MultiHisto<ht>::Write(const char* name, Int_t option, Int_t bufsize) const {
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);
136 for(int i=0; i < nhistos; i++)
137 if(m_myhistos[i]->Write(0,option,bufsize)) writtenhistos++;
140 return writtenhistos;
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);
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();
158 ////////////////////////////////////////////////////////////////////////////////
160 // Drawing the histos
162 ////////////////////////////////////////////////////////////////////////////////
165 // draw without canvas and without color
167 void MultiHisto<ht>::Draw(const Option_t* option){
169 std::vector<TCanvas *> c1;
173 // draw without canvas
175 void MultiHisto<ht>::Draw(const int color,const Option_t* option){
177 std::vector<TCanvas *> c1;
178 Draw(c1,color,option);
182 // draw with canvases
184 void MultiHisto<ht>::Draw(std::vector <TCanvas*> &c1, const int color,const Option_t* option){
187 while(c1.size() < GetNhistos()){
188 TCanvas *c = new TCanvas();
192 for(unsigned int i = 0 ; i < GetNhistos() ; i++){
193 ht* swap = GetHisto(i);
194 c1[i]->SetName(swap->GetName());
196 swap->SetLineColor(color);
204 ////////////////////////////////////////////////////////////////////////////////
206 // Getting global index
208 ////////////////////////////////////////////////////////////////////////////////
210 // from a detailed index set
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];
218 if(globalindex >= int(m_myhistos.size())) globalindex = -1;
222 // from a set of parametrs
224 int MultiHisto<ht>::GetGlobalIndex(const std::vector<double> &pars) const{
225 std::vector<int> indexes = GetDivisionsIndexes(pars);
226 return GetGlobalIndex(indexes);
230 int MultiHisto<ht>::GetGlobalIndex(const std::vector<float> &pars) const{
231 std::vector<int> indexes = GetDivisionsIndexes(pars);
232 return GetGlobalIndex(indexes);
235 ////////////////////////////////////////////////////////////////////////////////
237 // Getting detailed indexes
239 ////////////////////////////////////////////////////////////////////////////////
241 // from a global index
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];
251 indxfinal.push_back(curr_index);
252 for(int i = indx.size()-1; i >= 0; i--) indxfinal.push_back(indx[i]);
256 // from a set of parameters
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));
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);
276 std::vector<int> MultiHisto<ht>::GetDivisionsIndexes(const std::vector<float> &pars) const{
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));
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);
295 ////////////////////////////////////////////////////////////////////////////////
297 // Composing titles and names
299 ////////////////////////////////////////////////////////////////////////////////
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();
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]);
325 ///////////////////////////////////////////////////////////////////////////////
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 );
340 basename << low << " #leq " << div_name << " < " << high;
341 return basename.str();
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]);
355 ////////////////////////////////////////////////////////////////////////////////
359 ////////////////////////////////////////////////////////////////////////////////
362 const char * MultiHisto<ht>::GetName() const{
363 return m_model->GetName();
367 const char * MultiHisto<ht>::GetTitle() const{
368 return m_model->GetTitle();
373 unsigned int MultiHisto<ht>::GetNhistos() const{
374 return m_myhistos.size();
379 ht* MultiHisto<ht>::GetHisto(const std::vector<int> &indexes){
380 int i = GetGlobalIndex(indexes);
385 ht* MultiHisto<ht>::GetHisto(int i){
387 if(i >= 0 && i < int(m_myhistos.size())) current = m_myhistos[i];