]> git.uio.no Git - u/mrichter/AliRoot.git/blame - UNICOR/AliDHN.cxx
Code cleaned up, warnings and violations reduced.
[u/mrichter/AliRoot.git] / UNICOR / AliDHN.cxx
CommitLineData
7148817a 1// Author: Dariusz Miskowiec <mailto:d.miskowiec@gsi.de> 2007
2
3//=============================================================================
4// multidimensional histogram
5// Physically, the data are kept in one single one-dimensional histogram.
6// Projecting on 1, 2, and n-1 dimensions is implemented.
7// The histogram can be saved on root file as the one-dimensional data
8// histogram and the axes, thus eternal forward compatibility is ensured.
9//=============================================================================
10
5a6d201c 11#include <cmath>
7148817a 12#include <TFile.h>
13#include <TDirectory.h>
14#include <TAxis.h>
15#include <TH2.h>
7148817a 16#include "AliDHN.h"
17
18ClassImp(AliDHN)
19
20//=============================================================================
21AliDHN::AliDHN(Char_t *nam, Int_t ndim, TAxis **ax)
22 : TH1D(nam, nam, Albins(ndim,ax), 0.5, Albins(ndim,ax)+0.5),
23 fNdim(ndim)
24{
25 // Constructor for building from scratch.
5a6d201c 26
7148817a 27 // Above, we just have managed to create a single one-dimensional histogram
5a6d201c 28 // with number of bins equal to the product of the numbers of bins in all
7148817a 29 // dimensions. For easy inspection the histogram range was set to -0.5,n-0.5.
30
31 for (int i=0; i<fNdim; i++) ax[i]->Copy(fAxis[i]);
32
33 // for speed reasons, number of bins of each axis is stored on fNbins
34 // and the running product of the latter on fMbins
35 // so fMbins = {...,fNbins[fNdim-2]*fNbins[fNdim-1],fNbins[fNdim-1],1}
36
37 fMbins[fNdim-1] = 1;
38 for (int i=0; i<fNdim; i++) fNbins[i] = fAxis[i].GetNbins();
39 for (int i=fNdim-1; i>0; i--) fMbins[i-1] = fMbins[i]*fNbins[i];
40 printf("%d-dimensional histogram %s with %d bins created\n",fNdim,nam,GetNbinsX());
41}
42//=============================================================================
43AliDHN::AliDHN(Char_t *filnam, Char_t *nam)
44 : TH1D(*((TH1D*) TFile::Open(filnam,"read")->GetDirectory(nam)->Get("histo"))),
45 fNdim(0)
46{
47 // Constructor for reading from file.
48
49 TFile *f = TFile::Open(filnam,"read");
50 f->cd(nam);
51 TAxis *ax[fMaxNdim];
52 char axnam[1000];
53 for (fNdim=0; fNdim<fMaxNdim; fNdim++) {
54 sprintf(axnam,"axis%d",fNdim);
55 ax[fNdim] = (TAxis *) gDirectory->Get(axnam);
56 if (ax[fNdim]) ax[fNdim]->Copy(fAxis[fNdim]);
57 else break;
58 }
59 f->Close();
60
61 fMbins[fNdim-1] = 1;
62 for (int i=0; i<fNdim; i++) fNbins[i] = fAxis[i].GetNbins();
63 for (int i=fNdim-1; i>0; i--) fMbins[i-1] = fMbins[i]*fNbins[i];
64
65 if (GetNbinsX()!=Albins(fNdim,ax)) {
66 printf("number of bins of histo %d differs from product of nbins of axes %d\n",
67 GetNbinsX(),Albins(fNdim,ax));
68 printf("bombing\n");
69 exit(-1);
70 }
71 printf("%d-dimensional histogram %s with %d bins read from file %s\n",
72 fNdim,nam,GetNbinsX(),filnam);
73}
74//=============================================================================
75Int_t AliDHN::Albins(Int_t n, TAxis **ax)
76{
77 // Calculate product of nbins of ax[0]...ax[n-1]
78 // (= total number of bins of the multidimensional histogram to be made).
79
80 Int_t nbins = 1;
81 // while (n--) nbins *= ax[n]->GetNbins();
82 for (int i=0; i<n; i++) nbins *= ax[i]->GetNbins();
83 return nbins;
84}
85//=============================================================================
86Int_t AliDHN::MulToOne(Int_t *k) const
87{
88 // Calculate the 1-dim index n from n-dim indices k[fNdim].
89 // Valid k[i] should be between 0 and fNbins[i]-1.
90 // Valid result will be between 0 and GetNbinsX()-1.
91 // Return -1 if under- or overflow in any dimension.
92
93 Int_t n = 0;
94 for (int i=0; i<fNdim; i++) {
95 if (k[i]<0) return -1;
96 if (k[i]>=fNbins[i]) return -1;
97 n += fMbins[i]*k[i];
98 }
99 return n;
100}
101//=============================================================================
102Int_t AliDHN::MulToOne(Double_t *x)
103{
104 // Calculate the 1-dim index n from n-dim vector x, representing the
105 // abscissa of the n-dim histogram. The result will be between 0 and
106 // GetNbinsX()-1.
107
108 Int_t k[fMaxNdim];
109 for (int i=0; i<fNdim; i++) k[i] = fAxis[i].FindBin(x[i])-1;
110 return MulToOne(k);
111}
112//=============================================================================
5a6d201c 113void AliDHN::OneToMul(Int_t n, Int_t *k) const
7148817a 114{
115 // Calculate the n-dim indices k[fNdim] from 1-dim index n.
116 // Valid n should be between 0 and GetNbinsX()-1.
117 // Valid results will be between 0 and fNbins[i]-1.
118
119 div_t idi; // integer division structure
120 for (int i=0; i<fNdim; i++) {
121 idi = div(n,fMbins[i]);
122 k[i] = idi.quot; // quotient
123 n = idi.rem; // remainder
124 }
125}
126//=============================================================================
5a6d201c 127Int_t AliDHN::Fill(Double_t *xx, Double_t w)
7148817a 128{
129 // Fill the histogram. The array xx holds the abscissa information, w is the
130 // weigth. The 1-dim histogram is filled using the standard Fill method
131 // (direct access to the arrays was tried and was not faster).
132
133 int nbin = MulToOne(xx);
5a6d201c 134 if (nbin == -1) return 0;
135 return TH1D::Fill(nbin+1,w);
7148817a 136}
137//=============================================================================
5a6d201c 138Int_t AliDHN::Fill(Double_t x0, Double_t x1, ...)
139{
140 // Fill the histogram. Arguments are passed as doubles rather than array.
141 va_list ap;
142 Double_t xx[fMaxNdim] = {x0, x1};
143 va_start(ap,x1);
144 for (int i=2; i<fNdim; i++) xx[i] = va_arg(ap,Double_t);
145 Double_t weigth = va_arg(ap,Double_t);
146 va_end(ap);
147 return Fill(xx,weigth);
148}
149//=============================================================================
150Int_t AliDHN::Write() const
7148817a 151{
152 // Save the 1-dim histo and the axes in a subdirectory on file. This might
153 // not be the most elegant way but it is very simple and backward and forward
154 // compatible.
155
156 Int_t nbytes = 0;
157 TH1D histo(*this);
158 histo.SetName("histo"); // hadd bug fix
159 TDirectory *dest = gDirectory->mkdir(GetName());
160 if (dest) {
161 dest->cd();
162 nbytes += histo.Write("histo");
163 char axnam[1000];
164 for (int i=0; i<fNdim; i++) {
165 sprintf(axnam,"axis%d",i);
166 nbytes += fAxis[i].Write(axnam);
167 }
168 printf(" %s stored in %s\n",GetName(),dest->GetPath());
169 dest->cd("..");
170 }
171 return nbytes;
172}
173//=============================================================================
174AliDHN *AliDHN::ProjectAlong(char *nam, Int_t dim, Int_t first, Int_t last)
175{
176 // Reduce dimension dim by summing up its bins between first and last.
177 // Use root convention: bin=1 is the first bin, bin=nbins is the last.
178 // Return the resulting fNdim-1 dimensional histogram.
179
180 if (dim<0 || dim>fNdim-1) return 0;
181 if (first<0) first = 1;
182 if (last<0) last = fNbins[dim];
183
184 // create new (reduced) histogram
185
186 TAxis *ax[fMaxNdim-1];
187 int n=0;
188 for (int i=0; i<fNdim; i++) if (i!=dim) ax[n++] = GetAxis(i);
189 char *name = strlen(nam)? nam : Form("%s_wo%d",GetName(),dim);
190 char *eame = Form("%s_error",name);
191 AliDHN *his = new AliDHN(name,fNdim-1,ax); // result histogram
192 AliDHN *eis = new AliDHN(eame,fNdim-1,ax); // temporary storage for errors
193
194 // sum up the content and errors squared
195
196 Int_t *k = new Int_t[fNdim]; // old hist multiindex
197 Int_t *m = new Int_t[fNdim-1]; // new hist multiindex
198 for (int i=0; i<GetNbinsX(); i++) {
199 OneToMul(i,k);
200 if (k[dim]+1<first) continue;
201 if (k[dim]+1>last) continue;
5a6d201c 202 n = 0;
7148817a 203 for (int j=0; j<fNdim; j++) if (j!=dim) m[n++] = k[j];
204 n = his->MulToOne(m);
205 his->AddBinContent(n+1,GetBinContent(i+1));
206 eis->AddBinContent(n+1,GetBinError(i+1)*GetBinError(i+1));
207 }
208
209 // combine content and errors in one histogram
210
211 for (int i=0; i<his->GetNbinsX(); i++) {
5a6d201c 212 his->SetBinError(i+1,sqrt(eis->GetBinContent(i+1)));
7148817a 213 }
214
215 his->SetLineColor(this->GetLineColor());
216 his->SetFillColor(this->GetFillColor());
217 his->SetMarkerColor(this->GetMarkerColor());
218 his->SetMarkerStyle(this->GetMarkerStyle());
219
220 // some cleanup
221
222 delete eis;
223 return his;
224}
225//=============================================================================
5a6d201c 226TH1D *AliDHN::ProjectOn(char *nam, Int_t dim, Int_t *first, Int_t *last) const
7148817a 227{
228 // Project on dimension dim. Use only bins between first[i] and last[i].
229 // Use root convention: bin=1 is the first bin, bin=nbins is the last.
230 // first[i]=0 means from the first bin
231 // last[i]=0 means till the last bin
232
233 if (dim<0 || dim>fNdim-1) return 0;
234
235 // calculate the projection; lowest index 0
236
237 double *yy = new double[fNbins[dim]]; // value
238 double *ey = new double[fNbins[dim]]; // error
239 for (int i=0; i<fNbins[dim]; i++) yy[i]=0;
240 for (int i=0; i<fNbins[dim]; i++) ey[i]=0;
241 Int_t *k = new Int_t[fNdim];
242 for (int i=0; i<GetNbinsX(); i++) {
243 OneToMul(i,k);
244 int isgood = 1; // bin within the range?
245 for (int j=0; j<fNdim; j++) {
246 if (first) if (first[j]) if (k[j]+1<first[j]) {isgood = 0; break;}
247 if (last) if (last[j]) if (k[j]+1>last[j]) {isgood = 0; break;}
248 }
249 if (isgood) {
250 yy[k[dim]]+=GetBinContent(i+1);
251 ey[k[dim]]+=GetBinError(i+1)*GetBinError(i+1);
252 }
253 }
254
255 // make the projection histogram
256 // use name nam if specified; otherwise generate one
257
258 TH1D *his;
259 char *name = strlen(nam)? nam : Form("%s_proj%d",GetName(),dim);
260 if (fAxis[dim].IsVariableBinSize())
261 his = new TH1D(name,name,fNbins[dim],fAxis[dim].GetXbins()->GetArray());
262 else
263 his = new TH1D(name,name,fNbins[dim],fAxis[dim].GetXmin(),fAxis[dim].GetXmax());
264 his->SetXTitle(fAxis[dim].GetTitle());
265 // or his->GetXaxis()->ImportAttributes(ax);
266 his->Sumw2();
267 his->SetLineColor(this->GetLineColor());
268 his->SetFillColor(this->GetFillColor());
269 his->SetMarkerColor(this->GetMarkerColor());
270 his->SetMarkerStyle(this->GetMarkerStyle());
271 for (int i=0; i<his->GetNbinsX(); i++) {
272 his->SetBinContent(i+1,yy[i]);
5a6d201c 273 his->SetBinError(i+1,sqrt(ey[i]));
7148817a 274 }
275
276 // some cleanup
277
278 delete [] yy;
279 delete [] ey;
280 delete [] k;
281 // if (name!=nam) delete [] name;
282
283 return his;
284}
285//=============================================================================
286TH1D *AliDHN::ProjectOn(char *nam, Int_t dim, Double_t *first, Double_t *last)
287{
288 // Project on dimension dim. Use only bins between first[i] and last[i].
289
290 if (dim<0 || dim>fNdim-1) return 0;
291 Int_t kfirst[fMaxNdim];
292 Int_t klast[fMaxNdim];
293 for (int i=0; i<fNdim; i++) {
294 kfirst[i] = fAxis[i].FindBin(first[i]);
295 klast[i] = fAxis[i].FindBin(last[i]);
296 }
297 return ProjectOn(nam,dim,kfirst,klast);
298}
299//=============================================================================
5a6d201c 300TH2D *AliDHN::ProjectOn(char *nam, Int_t dim0, Int_t dim1, Int_t *first, Int_t *last) const
7148817a 301{
302 // Project on dim1 vs dim0. Use only bins between first[i] and last[i].
303 // Use root convention: bin=1 is the first bin, bin=nbins is the last.
304 // first[i]=0 means from the first bin
305 // last[i]=0 means till the last bin
306
307 if (dim0<0 || dim0>fNdim-1) return 0;
308 if (dim1<0 || dim1>fNdim-1) return 0;
309
310 // calculate the projection
311
312 double **yy = new double*[fNbins[dim0]]; // value
313 double **ey = new double*[fNbins[dim0]]; // error
314 for (int i=0; i<fNbins[dim0]; i++) yy[i] = new double[fNbins[dim1]];
315 for (int i=0; i<fNbins[dim0]; i++) ey[i] = new double[fNbins[dim1]];
316 for (int i=0; i<fNbins[dim0]; i++) for (int j=0; j<fNbins[dim1]; j++) yy[i][j]=0;
317 for (int i=0; i<fNbins[dim0]; i++) for (int j=0; j<fNbins[dim1]; j++) ey[i][j]=0;
318 Int_t *k = new Int_t[fNdim];
319 for (int i=0; i<GetNbinsX(); i++) {
320 OneToMul(i,k);
321 int isgood = 1; // bin within the range?
322 for (int j=0; j<fNdim; j++) {
323 if (first) if (first[j]) if (k[j]+1<first[j]) {isgood = 0; break;}
324 if (last) if (last[j]) if (k[j]+1>last[j]) {isgood = 0; break;}
325 }
326 if (isgood) {
327 yy[k[dim0]][k[dim1]]+=GetBinContent(i+1);
328 ey[k[dim0]][k[dim1]]+=GetBinError(i+1)*GetBinError(i+1);
329 }
330 }
331
332 // make the projection histogram
333 // use name nam if specified; otherwise generate one
334
335 TH2D *his=0;
336 char *name = strlen(nam)? nam : Form("%s_proj%dvs%d",GetName(),dim1,dim0);
337 if (fAxis[dim0].IsVariableBinSize() && fAxis[dim1].IsVariableBinSize())
338 his = new TH2D(name,name,
339 fNbins[dim0],fAxis[dim0].GetXbins()->GetArray(),
340 fNbins[dim1],fAxis[dim1].GetXbins()->GetArray());
341 else if (!fAxis[dim0].IsVariableBinSize() && fAxis[dim1].IsVariableBinSize())
342 his = new TH2D(name,name,
343 fNbins[dim0],fAxis[dim0].GetXmin(),fAxis[dim0].GetXmax(),
344 fNbins[dim1],fAxis[dim1].GetXbins()->GetArray());
345 else if (fAxis[dim0].IsVariableBinSize() && !fAxis[dim1].IsVariableBinSize())
346 his = new TH2D(name,name,
347 fNbins[dim0],fAxis[dim0].GetXbins()->GetArray(),
348 fNbins[dim1],fAxis[dim1].GetXmin(),fAxis[dim1].GetXmax());
349 else if (!fAxis[dim0].IsVariableBinSize() && !fAxis[dim1].IsVariableBinSize())
350 his = new TH2D(name,name,
351 fNbins[dim0],fAxis[dim0].GetXmin(),fAxis[dim0].GetXmax(),
352 fNbins[dim1],fAxis[dim1].GetXmin(),fAxis[dim1].GetXmax());
353 his->SetXTitle(fAxis[dim0].GetTitle());
354 his->SetYTitle(fAxis[dim1].GetTitle());
355 // or his->GetXaxis()->ImportAttributes(ax0);
356 // or his->GetYaxis()->ImportAttributes(ax1);
357 his->Sumw2();
358 his->SetLineColor(this->GetLineColor());
359 his->SetFillColor(this->GetFillColor());
360 his->SetMarkerColor(this->GetMarkerColor());
361 his->SetMarkerStyle(this->GetMarkerStyle());
362 for (int i=0; i<his->GetNbinsX(); i++)
363 for (int j=0; j<his->GetNbinsY(); j++) {
364 his->SetBinContent(i+1,j+1,yy[i][j]);
5a6d201c 365 his->SetBinError(i+1,j+1,sqrt(ey[i][j]));
7148817a 366 }
367
368 // some cleanup
369
370 for (int i=0; i<fNbins[dim0]; i++) delete [] yy[i];
371 for (int i=0; i<fNbins[dim0]; i++) delete [] ey[i];
372 delete [] yy;
373 delete [] ey;
374 delete [] k;
375 // if (name!=nam) delete [] name;
376
377 return his;
378}
379//=============================================================================