]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliDielectronHistos.cxx
filtering updates
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliDielectronHistos.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 //
17 // Generic Histogram container with support for groups and filling of groups by passing
18 // a vector of data
19 //
20 // Authors: 
21 //   Jens Wiechula <Jens.Wiechula@cern.ch> 
22 //   Julian Book   <Julian.Book@cern.ch> 
23 // 
24
25 #include <TH1.h>
26 #include <TH1F.h>
27 #include <TH2.h>
28 #include <TH3.h>
29 #include <THnBase.h>
30 #include <THn.h>
31 #include <THnSparse.h>
32 #include <TProfile.h>
33 #include <TProfile2D.h>
34 #include <TProfile3D.h>
35 #include <TCollection.h>
36 #include <THashList.h>
37 #include <TString.h>
38 #include <TObjString.h>
39 #include <TObjArray.h>
40 #include <TFile.h>
41 #include <TError.h>
42 #include <TCanvas.h>
43 #include <TMath.h>
44 #include <TROOT.h>
45 #include <TLegend.h>
46 #include <TKey.h>
47 #include <TAxis.h>
48 #include <TVirtualPS.h>
49 #include <TVectorD.h>
50
51 #include "AliDielectronHelper.h"
52 #include "AliDielectronVarManager.h"
53 #include "AliDielectronHistos.h"
54
55 ClassImp(AliDielectronHistos)
56
57
58 AliDielectronHistos::AliDielectronHistos() :
59 //   TCollection(),
60   TNamed("AliDielectronHistos","Dielectron Histogram Container"),
61   fHistoList(),
62   fList(0x0),
63   fUsedVars(new TBits(AliDielectronVarManager::kNMaxValues)),
64   fReservedWords(new TString)
65 {
66   //
67   // Default constructor
68   //
69   fHistoList.SetOwner(kTRUE);
70   fHistoList.SetName("Dielectron_Histos");
71 }
72
73 //_____________________________________________________________________________
74 AliDielectronHistos::AliDielectronHistos(const char* name, const char* title) :
75 //   TCollection(),
76   TNamed(name, title),
77   fHistoList(),
78   fList(0x0),
79   fUsedVars(new TBits(AliDielectronVarManager::kNMaxValues)),
80   fReservedWords(new TString)
81 {
82   //
83   // TNamed constructor
84   //
85   fHistoList.SetOwner(kTRUE);
86   fHistoList.SetName(name);
87 }
88
89 //_____________________________________________________________________________
90 AliDielectronHistos::~AliDielectronHistos()
91 {
92   //
93   // Destructor
94   //
95   fHistoList.Clear();
96   if (fUsedVars) delete fUsedVars;
97   if (fList) fList->Clear();
98   delete fReservedWords;
99 }
100
101 //_____________________________________________________________________________
102 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
103                                       UInt_t valTypeP,
104                                       Int_t nbinsX, Double_t xmin, Double_t xmax,
105                                       UInt_t valTypeX, Bool_t logBinX, TString option,
106                                       UInt_t valTypeW)
107 {
108   //
109   // Default histogram creation 1D case
110   //
111
112   TVectorD *binLimX=0x0;
113   
114   if (logBinX) {
115     binLimX=AliDielectronHelper::MakeLogBinning(nbinsX, xmin, xmax);
116   } else {
117     binLimX=AliDielectronHelper::MakeLinBinning(nbinsX, xmin, xmax);
118   }
119   UserProfile(histClass,name,title,valTypeP,binLimX,valTypeX,option,valTypeW);
120 }
121
122 //_____________________________________________________________________________
123 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
124                                       UInt_t valTypeP,
125                                       Int_t nbinsX, Double_t xmin, Double_t xmax,
126                                       Int_t nbinsY, Double_t ymin, Double_t ymax,
127                                       UInt_t valTypeX, UInt_t valTypeY,
128                                       Bool_t logBinX, Bool_t logBinY, TString option,
129                                       UInt_t valTypeW)
130 {
131   //
132   // Default histogram creation 2D case
133   //
134   if (!IsHistogramOk(histClass,name)) return;
135
136   TVectorD *binLimX=0x0;
137   TVectorD *binLimY=0x0;
138   
139   if (logBinX) {
140     binLimX=AliDielectronHelper::MakeLogBinning(nbinsX, xmin, xmax);
141   } else {
142     binLimX=AliDielectronHelper::MakeLinBinning(nbinsX, xmin, xmax);
143   }
144   if (logBinY) {
145     binLimY=AliDielectronHelper::MakeLogBinning(nbinsY, ymin, ymax);
146   } else {
147     binLimY=AliDielectronHelper::MakeLinBinning(nbinsY, ymin, ymax);
148   }
149   UserProfile(histClass,name,title,valTypeP,binLimX,binLimY,valTypeX,valTypeY,option,valTypeW);
150 }
151
152
153 //_____________________________________________________________________________
154 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
155                                       UInt_t valTypeP,
156                                       Int_t nbinsX, Double_t xmin, Double_t xmax,
157                                       Int_t nbinsY, Double_t ymin, Double_t ymax,
158                                       Int_t nbinsZ, Double_t zmin, Double_t zmax,
159                                       UInt_t valTypeX, UInt_t valTypeY, UInt_t valTypeZ,
160                                       Bool_t logBinX, Bool_t logBinY, Bool_t logBinZ, TString option,
161                                       UInt_t valTypeW)
162 {
163   //
164   // Default histogram creation 3D case
165   //
166   if (!IsHistogramOk(histClass,name)) return;
167
168   TVectorD *binLimX=0x0;
169   TVectorD *binLimY=0x0;
170   TVectorD *binLimZ=0x0;
171   
172   if (logBinX) {
173     binLimX=AliDielectronHelper::MakeLogBinning(nbinsX, xmin, xmax);
174   } else {
175     binLimX=AliDielectronHelper::MakeLinBinning(nbinsX, xmin, xmax);
176   }
177   
178   if (logBinY) {
179     binLimY=AliDielectronHelper::MakeLogBinning(nbinsY, ymin, ymax);
180   } else {
181     binLimY=AliDielectronHelper::MakeLinBinning(nbinsY, ymin, ymax);
182   }
183   
184   if (logBinZ) {
185     binLimZ=AliDielectronHelper::MakeLogBinning(nbinsZ, zmin, zmax);
186   } else {
187     binLimZ=AliDielectronHelper::MakeLinBinning(nbinsZ, zmin, zmax);
188   }
189
190   UserProfile(histClass,name,title,valTypeP,binLimX,binLimY,binLimZ,valTypeX,valTypeY,valTypeZ,option,valTypeW);
191 }
192
193 //_____________________________________________________________________________
194 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
195                                       UInt_t valTypeP,
196                                       const char* binning,
197                                       UInt_t valTypeX, TString option,
198                                       UInt_t valTypeW)
199 {
200   //
201   // Histogram creation 1D case with arbitraty binning
202   //
203
204   TVectorD *binLimX=AliDielectronHelper::MakeArbitraryBinning(binning);
205   UserProfile(histClass,name,title,valTypeP,binLimX,valTypeX,option,valTypeW);
206 }
207
208 //_____________________________________________________________________________
209 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
210                                       UInt_t valTypeP,
211                                       const TVectorD * const binsX,
212                                       UInt_t valTypeX/*=kNoAutoFill*/, TString option,
213                                       UInt_t valTypeW)
214 {
215   //
216   // Histogram creation 1D case with arbitraty binning X
217   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
218   //
219
220   Bool_t isOk=kTRUE;
221   isOk&=IsHistogramOk(histClass,name);
222   isOk&=(binsX!=0x0);
223   TH1 *hist=0x0;
224
225   if (isOk){
226     if(valTypeP==kNoProfile)
227       hist=new TH1F(name,title,binsX->GetNrows()-1,binsX->GetMatrixArray());
228     else {
229       TString opt=""; Double_t pmin=0., pmax=0.;
230       if(!option.IsNull()) {
231         TObjArray *arr=option.Tokenize(";");
232         arr->SetOwner();
233         opt=((TObjString*)arr->At(0))->GetString();
234         if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
235         if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
236         delete arr;
237       }
238       hist=new TProfile(name,title,binsX->GetNrows()-1,binsX->GetMatrixArray());
239       ((TProfile*)hist)->BuildOptions(pmin,pmax,opt.Data());
240       //      printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
241     }
242
243     // store variales in axes
244     UInt_t valType[20] = {0};
245     valType[0]=valTypeX;     valType[1]=valTypeP;
246     StoreVariables(hist, valType);
247     hist->SetUniqueID(valTypeW); // store weighting variable
248
249     // store which variables are used
250     for(Int_t i=0; i<4; i++)   fUsedVars->SetBitNumber(valType[i],kTRUE);
251     fUsedVars->SetBitNumber(valTypeW,kTRUE);
252
253     // adapt the name and title of the histogram in case they are empty
254     AdaptNameTitle(hist, histClass);
255
256     Bool_t isReserved=fReservedWords->Contains(histClass);
257     if(valTypeX==kNoAutoFill) hist->SetUniqueID(valTypeX);
258     if (isReserved)
259       UserHistogramReservedWords(histClass, hist, 999);
260     else
261       UserHistogram(histClass, hist, 999);
262   }
263   
264   delete binsX;
265 }
266
267 //_____________________________________________________________________________
268 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
269                                       UInt_t valTypeP,
270                                       const TVectorD * const binsX, const TVectorD * const binsY,
271                                       UInt_t valTypeX/*=kNoAutoFill*/, UInt_t valTypeY/*=0*/, TString option,
272                                       UInt_t valTypeW)
273 {
274   //
275   // Histogram creation 2D case with arbitraty binning X
276   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
277   //
278
279   Bool_t isOk=kTRUE;
280   isOk&=IsHistogramOk(histClass,name);
281   isOk&=(binsX!=0x0);
282   isOk&=(binsY!=0x0);
283   TH1 *hist=0x0;
284
285   if (isOk){
286     if(valTypeP==kNoProfile) {
287       hist=new TH2F(name,title,
288                     binsX->GetNrows()-1,binsX->GetMatrixArray(),
289                     binsY->GetNrows()-1,binsY->GetMatrixArray()); 
290     }
291     else  {
292       TString opt=""; Double_t pmin=0., pmax=0.;
293       if(!option.IsNull()) {
294         TObjArray *arr=option.Tokenize(";");
295         arr->SetOwner();
296         opt=((TObjString*)arr->At(0))->GetString();
297         if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
298         if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
299         delete arr;
300       }
301       hist=new TProfile2D(name,title,
302                           binsX->GetNrows()-1,binsX->GetMatrixArray(),
303                           binsY->GetNrows()-1,binsY->GetMatrixArray());
304       ((TProfile2D*)hist)->BuildOptions(pmin,pmax,opt.Data());
305       //      printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
306     }
307
308     // store variales in axes
309     UInt_t valType[20] = {0};
310     valType[0]=valTypeX;     valType[1]=valTypeY;     valType[2]=valTypeP;
311     StoreVariables(hist, valType);
312     hist->SetUniqueID(valTypeW); // store weighting variable
313
314     // store which variables are used
315     for(Int_t i=0; i<4; i++)   fUsedVars->SetBitNumber(valType[i],kTRUE);
316     fUsedVars->SetBitNumber(valTypeW,kTRUE);
317
318     // adapt the name and title of the histogram in case they are empty
319     AdaptNameTitle(hist, histClass);
320
321     Bool_t isReserved=fReservedWords->Contains(histClass);
322     if(valTypeX==kNoAutoFill) hist->SetUniqueID(valTypeX);
323     if (isReserved)
324       UserHistogramReservedWords(histClass, hist, 999);
325     else
326       UserHistogram(histClass, hist, 999);
327   }
328   
329   delete binsX;
330   delete binsY;
331   
332 }
333
334 //_____________________________________________________________________________
335 void AliDielectronHistos::UserProfile(const char* histClass,const char *name, const char* title,
336                                       UInt_t valTypeP,
337                                       const TVectorD * const binsX, const TVectorD * const binsY, const TVectorD * const binsZ,
338                                       UInt_t valTypeX/*=kNoAutoFill*/, UInt_t valTypeY/*=0*/, UInt_t valTypeZ/*=0*/, TString option,
339                                       UInt_t valTypeW)
340 {
341   //
342   // Histogram creation 3D case with arbitraty binning X
343   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
344   //
345
346   Bool_t isOk=kTRUE;
347   isOk&=IsHistogramOk(histClass,name);
348   isOk&=(binsX!=0x0);
349   isOk&=(binsY!=0x0);
350   isOk&=(binsZ!=0x0);
351   TH1 *hist=0x0;
352
353   if (isOk) {
354     if(valTypeP==kNoProfile) {
355       hist=new TH3F(name,title,
356                     binsX->GetNrows()-1,binsX->GetMatrixArray(),
357                     binsY->GetNrows()-1,binsY->GetMatrixArray(),
358                     binsZ->GetNrows()-1,binsZ->GetMatrixArray());
359     }
360     else {
361       TString opt=""; Double_t pmin=0., pmax=0.;
362       if(!option.IsNull()) {
363         TObjArray *arr=option.Tokenize(";");
364         arr->SetOwner();
365         opt=((TObjString*)arr->At(0))->GetString();
366         if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
367         if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
368         delete arr;
369       }
370       hist=new TProfile3D(name,title,
371                           binsX->GetNrows()-1,binsX->GetMatrixArray(),
372                           binsY->GetNrows()-1,binsY->GetMatrixArray(),
373                           binsZ->GetNrows()-1,binsZ->GetMatrixArray());
374       ((TProfile3D*)hist)->BuildOptions(pmin,pmax,opt.Data());
375       //      printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
376     }
377
378     // store variales in axes
379     UInt_t valType[20] = {0};
380     valType[0]=valTypeX;     valType[1]=valTypeY;     valType[2]=valTypeZ;     valType[3]=valTypeP;
381     StoreVariables(hist, valType);
382     hist->SetUniqueID(valTypeW); // store weighting variable
383
384     // store which variables are used
385     for(Int_t i=0; i<4; i++)   fUsedVars->SetBitNumber(valType[i],kTRUE);
386     fUsedVars->SetBitNumber(valTypeW,kTRUE);
387
388     // adapt the name and title of the histogram in case they are empty
389     AdaptNameTitle(hist, histClass);
390
391     Bool_t isReserved=fReservedWords->Contains(histClass);
392     if(valTypeX==kNoAutoFill) hist->SetUniqueID(valTypeX);
393     if (isReserved)
394       UserHistogramReservedWords(histClass, hist, 999);
395     else
396       UserHistogram(histClass, hist, 999);
397   }
398   
399   delete binsX;
400   delete binsY;
401   delete binsZ;
402 }
403
404 //_____________________________________________________________________________
405 void AliDielectronHistos::UserHistogram(const char* histClass, Int_t ndim, Int_t *bins, Double_t *mins, Double_t *maxs, UInt_t *vars, UInt_t valTypeW)
406 {
407   //
408   // Histogram creation 4-n dimension only with linear binning
409   //
410
411   Bool_t isOk=kTRUE;
412   isOk&=(ndim<21 && ndim>3);
413   if(!isOk) { Warning("UserHistogram","Array sizes should be between 3 and 20. Not adding Histogram to '%s'.", histClass); return; }
414
415   // set automatic histo name
416   TString name;
417   for(Int_t iv=0; iv < ndim; iv++)
418     name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
419   name.Resize(name.Length()-1);
420
421   isOk&=IsHistogramOk(histClass,name);
422
423   THnD *hist;
424   if (isOk) {
425     hist=new THnD(name.Data(),"", ndim, bins, mins, maxs);
426
427     // store variales in axes
428     StoreVariables(hist, vars);
429     hist->SetUniqueID(valTypeW); // store weighting variable
430
431     // store which variables are used
432     for(Int_t i=0; i<20; i++)   fUsedVars->SetBitNumber(vars[i],kTRUE);
433     fUsedVars->SetBitNumber(valTypeW,kTRUE);
434
435     Bool_t isReserved=fReservedWords->Contains(histClass);
436     if (isReserved)
437       UserHistogramReservedWords(histClass, hist, 999);
438     else
439       UserHistogram(histClass, hist, 999);
440
441   }
442 }
443
444 //_____________________________________________________________________________
445 void AliDielectronHistos::UserHistogram(const char* histClass, Int_t ndim, TObjArray *limits, UInt_t *vars, UInt_t valTypeW)
446 {
447   //
448   // Histogram creation n>3 dimension only with non-linear binning
449   //
450
451   Bool_t isOk=kTRUE;
452   isOk&=(ndim<21 && ndim>3);
453   if(!isOk) { Warning("UserHistogram","Array sizes should be between 3 and 20. Not adding Histogram to '%s'.", histClass); return; }
454   isOk&=(ndim==limits->GetEntriesFast());
455   if(!isOk) return;
456
457   // set automatic histo name
458   TString name;
459   for(Int_t iv=0; iv < ndim; iv++)
460     name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
461   name.Resize(name.Length()-1);
462
463   isOk&=IsHistogramOk(histClass,name);
464
465   THnD *hist;
466   Int_t bins[ndim];
467   if (isOk) {
468     // get number of bins
469     for(Int_t idim=0 ;idim<ndim; idim++) {
470       TVectorD *vec = (TVectorD*) limits->At(idim);
471       bins[idim]=vec->GetNrows()-1;
472     }
473
474     hist=new THnD(name.Data(),"", ndim, bins, 0x0, 0x0);
475
476     // set binning
477     for(Int_t idim=0 ;idim<ndim; idim++) {
478       TVectorD *vec = (TVectorD*) limits->At(idim);
479       hist->SetBinEdges(idim,vec->GetMatrixArray());
480     }
481
482     // store variales in axes
483     StoreVariables(hist, vars);
484     hist->SetUniqueID(valTypeW); // store weighting variable
485
486     // store which variables are used
487     for(Int_t i=0; i<20; i++)   fUsedVars->SetBitNumber(vars[i],kTRUE);
488     fUsedVars->SetBitNumber(valTypeW,kTRUE);
489
490     Bool_t isReserved=fReservedWords->Contains(histClass);
491     if (isReserved)
492       UserHistogramReservedWords(histClass, hist, 999);
493     else
494       UserHistogram(histClass, hist, 999);
495
496   }
497 }
498
499 //_____________________________________________________________________________
500 void AliDielectronHistos::UserSparse(const char* histClass, Int_t ndim, Int_t *bins, Double_t *mins, Double_t *maxs, UInt_t *vars,
501                                      UInt_t valTypeW)
502 {
503   //
504   // THnSparse creation with linear binning
505   //
506
507   Bool_t isOk=kTRUE;
508
509   // set automatic histo name
510   TString name;
511   for(Int_t iv=0; iv < ndim; iv++)
512     name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
513   name.Resize(name.Length()-1);
514
515   isOk&=IsHistogramOk(histClass,name);
516
517   THnSparseD *hist;
518   if (isOk) {
519     hist=new THnSparseD(name.Data(),"", ndim, bins, mins, maxs);
520
521     // store variales in axes
522     StoreVariables(hist, vars);
523     hist->SetUniqueID(valTypeW); // store weighting variable
524
525     // store which variables are used
526     for(Int_t i=0; i<20; i++)   fUsedVars->SetBitNumber(vars[i],kTRUE);
527     fUsedVars->SetBitNumber(valTypeW,kTRUE);
528
529     Bool_t isReserved=fReservedWords->Contains(histClass);
530     if (isReserved)
531       UserHistogramReservedWords(histClass, hist, 999);
532     else
533       UserHistogram(histClass, hist, 999);
534
535   }
536 }
537
538 //_____________________________________________________________________________
539 void AliDielectronHistos::UserSparse(const char* histClass, Int_t ndim, TObjArray *limits, UInt_t *vars, UInt_t valTypeW)
540 {
541   //
542   // THnSparse creation with non-linear binning
543   //
544
545   Bool_t isOk=kTRUE;
546   isOk&=(ndim==limits->GetEntriesFast());
547   if(!isOk) return;
548
549   // set automatic histo name
550   TString name;
551   for(Int_t iv=0; iv < ndim; iv++)
552     name+=Form("%s_",AliDielectronVarManager::GetValueName(vars[iv]));
553   name.Resize(name.Length()-1);
554
555   isOk&=IsHistogramOk(histClass,name);
556
557   THnSparseD *hist;
558   Int_t bins[ndim];
559   if (isOk) {
560     // get number of bins
561     for(Int_t idim=0 ;idim<ndim; idim++) {
562       TVectorD *vec = (TVectorD*) limits->At(idim);
563       bins[idim]=vec->GetNrows()-1;
564     }
565
566     hist=new THnSparseD(name.Data(),"", ndim, bins, 0x0, 0x0);
567
568     // set binning
569     for(Int_t idim=0 ;idim<ndim; idim++) {
570       TVectorD *vec = (TVectorD*) limits->At(idim);
571       hist->SetBinEdges(idim,vec->GetMatrixArray());
572     }
573
574     // store variales in axes
575     StoreVariables(hist, vars);
576     hist->SetUniqueID(valTypeW); // store weighting variable
577
578     // store which variables are used
579     for(Int_t i=0; i<20; i++)   fUsedVars->SetBitNumber(vars[i],kTRUE);
580     fUsedVars->SetBitNumber(valTypeW,kTRUE);
581
582     Bool_t isReserved=fReservedWords->Contains(histClass);
583     if (isReserved)
584       UserHistogramReservedWords(histClass, hist, 999);
585     else
586       UserHistogram(histClass, hist, 999);
587
588   }
589 }
590
591 //_____________________________________________________________________________
592 void AliDielectronHistos::UserHistogram(const char* histClass, TObject* hist, UInt_t valTypes)
593 {
594   //
595   // Add any type of user histogram
596   //
597
598   //special case for the calss Pair. where histograms will be created for all pair classes
599   Bool_t isReserved=fReservedWords->Contains(histClass);
600   if (isReserved) {
601     UserHistogramReservedWords(histClass, hist, valTypes);
602     return;
603   }
604
605   if (!IsHistogramOk(histClass,hist->GetName())) return;
606   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
607   //  hist->SetDirectory(0);
608
609   // store variables axis
610   UInt_t valType[20] = {0};
611   // incase valTypes is given old way of extracting variables
612   if(valTypes!=999) {
613     valType[0]=valTypes%1000;          //last three digits
614     valType[1]=valTypes/1000%1000;     //second last three digits
615     valType[2]=valTypes/1000000%1000;  //third last three digits
616     hist->SetUniqueID(valTypes);
617   }
618   else {
619     // extract variables from axis
620     FillVarArray(hist, valType);
621     StoreVariables(hist, valType);
622     hist->SetUniqueID(valType[19]); // store weighting variable
623   }
624
625   classTable->Add(hist);
626 }
627
628 //_____________________________________________________________________________
629 void AliDielectronHistos::AddClass(const char* histClass)
630 {
631   //
632   // Add a class of histograms
633   // Several classes can be added by separating them by a ';' e.g. 'class1;class2;class3'
634   //
635   TString hists(histClass);
636   TObjArray *arr=hists.Tokenize(";");
637   TIter next(arr);
638   TObject *o=0;
639   while ( (o=next()) ){
640     if (fHistoList.FindObject(o->GetName())){
641       Warning("AddClass","Cannot create class '%s' it already exists.",histClass);
642       continue;
643     }
644     if (fReservedWords->Contains(o->GetName())){
645       Error("AddClass","Pair is a reserved word, please use another name");
646       continue;
647     }
648     THashList *table=new THashList;
649     table->SetOwner(kTRUE);
650     table->SetName(o->GetName());
651     fHistoList.Add(table);
652   }
653   delete arr;
654 }
655
656 //_____________________________________________________________________________
657 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval)
658 {
659   //
660   // Fill function 1D case
661   //
662   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
663   TH1* hist=0;
664   if (!classTable || !(hist=(TH1*)classTable->FindObject(name)) ){
665     Warning("Fill","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
666     return;
667   }
668   hist->Fill(xval);
669 }
670
671 //_____________________________________________________________________________
672 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval, Double_t yval)
673 {
674   //
675   // Fill function 2D case
676   //
677   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
678   TH2* hist=0;
679   if (!classTable || !(hist=(TH2*)classTable->FindObject(name)) ){
680     Warning("UserHistogram","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
681     return;
682   }
683   hist->Fill(xval,yval);
684 }
685
686 //_____________________________________________________________________________
687 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval, Double_t yval, Double_t zval)
688 {
689   //
690   // Fill function 3D case
691   //
692   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
693   TH3* hist=0;
694   if (!classTable || !(hist=(TH3*)classTable->FindObject(name)) ){
695     Warning("UserHistogram","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
696     return;
697   }
698   hist->Fill(xval,yval,zval);
699 }
700
701 //_____________________________________________________________________________
702 void AliDielectronHistos::FillClass(const char* histClass, Int_t nValues, const Double_t *values)
703 {
704   //
705   // Fill class 'histClass' (by name)
706   //
707
708   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
709   if (!classTable){
710     Warning("FillClass","Cannot fill class '%s' its not defined. nValues %d",histClass,nValues);
711     return;
712   }
713
714   TIter nextHist(classTable);
715   TObject *obj=0;
716   while ( (obj=(TObject*)nextHist()) )  FillValues(obj, values);
717
718   return;
719 }
720
721 //_____________________________________________________________________________
722 // void AliDielectronHistos::FillClass(const char* histClass, const TVectorD &vals)
723 // {
724 //   //
725 //   //
726 //   //
727 //   FillClass(histClass, vals.GetNrows(), vals.GetMatrixArray());
728 // }
729
730 //_____________________________________________________________________________
731 void AliDielectronHistos::UserHistogramReservedWords(const char* histClass, const TObject *hist, UInt_t valTypes)
732 {
733   //
734   // Creation of histogram for all pair types
735   //
736   TString title(hist->GetTitle());
737   // Same Event Like Sign
738   TIter nextClass(&fHistoList);
739   THashList *l=0;
740   while ( (l=static_cast<THashList*>(nextClass())) ){
741     TString name(l->GetName());
742     if (name.Contains(histClass)){
743       TObject *h=hist->Clone();
744       // Tobject has no function SetDirectory, didn't we need this???
745       //      h->SetDirectory(0);
746       ((TH1*)h)->SetTitle(Form("%s %s",title.Data(),l->GetName()));
747
748       UserHistogram(l->GetName(),h,valTypes);
749     }
750   }
751   delete hist;
752 }
753
754 //_____________________________________________________________________________
755 void AliDielectronHistos::DumpToFile(const char* file)
756 {
757   //
758   // Dump the histogram list to a newly created root file
759   //
760   TFile f(file,"recreate");
761   fHistoList.Write(fHistoList.GetName(),TObject::kSingleKey);
762   f.Close();
763 }
764
765 //_____________________________________________________________________________
766 TObject* AliDielectronHistos::GetHist(const char* histClass, const char* name) const
767 {
768   //
769   // return object 'name' in 'histClass'
770   //
771   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
772   if (!classTable) return 0x0;
773   return classTable->FindObject(name);
774 }
775
776 //_____________________________________________________________________________
777 TH1* AliDielectronHistos::GetHistogram(const char* histClass, const char* name) const
778 {
779   //
780   // return histogram 'name' in 'histClass'
781   //
782   return ((TH1*) GetHist(histClass, name));
783 }
784
785 //_____________________________________________________________________________
786 TObject* AliDielectronHistos::GetHist(const char* cutClass, const char* histClass, const char* name) const
787 {
788   //
789   // return object from list of list of histograms
790   // this function is thought for retrieving histograms if a list of AliDielectronHistos is set
791   //
792   
793   if (!fList) return 0x0;
794   THashList *h=dynamic_cast<THashList*>(fList->FindObject(cutClass));
795   if (!h)return 0x0;
796   THashList *classTable=dynamic_cast<THashList*>(h->FindObject(histClass));
797   if (!classTable) return 0x0;
798   return classTable->FindObject(name);
799 }
800
801 //_____________________________________________________________________________
802 TH1* AliDielectronHistos::GetHistogram(const char* cutClass, const char* histClass, const char* name) const
803 {
804   //
805   // return histogram from list of list of histograms
806   // this function is thought for retrieving histograms if a list of AliDielectronHistos is set
807   //
808   return ((TH1*) GetHist(cutClass, histClass, name));
809 }
810
811 //_____________________________________________________________________________
812 void AliDielectronHistos::Draw(const Option_t* option)
813 {
814   //
815   // Draw histograms
816   //
817
818   TString drawStr(option);
819   TObjArray *arr=drawStr.Tokenize(";");
820   arr->SetOwner();
821   TIter nextOpt(arr);
822
823   TString drawClasses;
824   TObjString *ostr=0x0;
825
826   TString currentOpt;
827   TString testOpt;
828   while ( (ostr=(TObjString*)nextOpt()) ){
829     currentOpt=ostr->GetString();
830     currentOpt.Remove(TString::kBoth,'\t');
831     currentOpt.Remove(TString::kBoth,' ');
832
833     testOpt="classes=";
834     if ( currentOpt.Contains(testOpt.Data()) ){
835       drawClasses=currentOpt(testOpt.Length(),currentOpt.Length());
836     }
837   }
838
839   delete arr;
840   drawStr.ToLower();
841   //optionsfList
842 //   Bool_t same=drawOpt.Contains("same"); //FIXME not yet implemented
843
844   TCanvas *c=0x0;
845   if (gVirtualPS) {
846     if (!gPad){
847       Error("Draw","When writing to a file you have to create a canvas before opening the file!!!");
848       return;
849     }
850     c=gPad->GetCanvas();
851     c->cd();
852 //     c=new TCanvas;
853   }
854   
855   TIter nextClass(&fHistoList);
856   THashList *classTable=0;
857 //   Bool_t first=kTRUE;
858   while ( (classTable=(THashList*)nextClass()) ){
859     //test classes option
860     if (!drawClasses.IsNull() && !drawClasses.Contains(classTable->GetName())) continue;
861     //optimised division
862     Int_t nPads = classTable->GetEntries();
863     Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
864     Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
865
866     //create canvas
867     if (!gVirtualPS){
868       TString canvasName;
869       canvasName.Form("c%s_%s",GetName(),classTable->GetName());
870       c=(TCanvas*)gROOT->FindObject(canvasName.Data());
871       if (!c) c=new TCanvas(canvasName.Data(),Form("%s: %s",GetName(),classTable->GetName()));
872       c->Clear();
873     } else {
874 //       if (first){
875 //         first=kFALSE;
876 //         if (nPads>1) gVirtualPS->NewPage();
877 //       } else {
878         if (nPads>1) c->Clear();
879 //       }
880     }
881     if (nCols>1||nRows>1) c->Divide(nCols,nRows);
882     
883     //loop over histograms and draw them
884     TIter nextHist(classTable);
885     Int_t iPad=0;
886     TH1 *h=0;
887     while ( (h=(TH1*)nextHist()) ){
888       TString drawOpt;
889       if ( (h->InheritsFrom(TH2::Class())) ) drawOpt="colz";
890       if (nCols>1||nRows>1) c->cd(++iPad);
891       if ( TMath::Abs(h->GetXaxis()->GetBinWidth(1)-h->GetXaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogx();
892       if ( TMath::Abs(h->GetYaxis()->GetBinWidth(1)-h->GetYaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogy();
893       if ( TMath::Abs(h->GetZaxis()->GetBinWidth(1)-h->GetZaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogz();
894       TString histOpt=h->GetOption();
895       histOpt.ToLower();
896       if (histOpt.Contains("logx")) gPad->SetLogx();
897       if (histOpt.Contains("logy")) gPad->SetLogy();
898       if (histOpt.Contains("logz")) gPad->SetLogz();
899       histOpt.ReplaceAll("logx","");
900       histOpt.ReplaceAll("logy","");
901       histOpt.ReplaceAll("logz","");
902       h->Draw(drawOpt.Data());
903     }
904     if (gVirtualPS) {
905       c->Update();
906     }
907     
908   }
909 //   if (gVirtualPS) delete c;
910 }
911
912 //_____________________________________________________________________________
913 void AliDielectronHistos::Print(const Option_t* option) const
914 {
915   //
916   // Print classes and histograms
917   //
918   TString optString(option);
919
920   if (optString.IsNull()) PrintStructure();
921
922
923
924 }
925
926 //_____________________________________________________________________________
927 void AliDielectronHistos::PrintStructure() const
928 {
929   //
930   // Print classes and histograms in the class to stdout
931   //
932   if (!fList){
933     TIter nextClass(&fHistoList);
934     THashList *classTable=0;
935     while ( (classTable=(THashList*)nextClass()) ){
936       TIter nextHist(classTable);
937       TObject *o=0;
938       printf("+ %s\n",classTable->GetName());
939       while ( (o=nextHist()) )
940         printf("| ->%s\n",o->GetName());
941     }
942   } else {
943     TIter nextCutClass(fList);
944     THashList *cutClass=0x0;
945     while ( (cutClass=(THashList*)nextCutClass()) ) {
946       printf("+ %s\n",cutClass->GetName());
947       TIter nextClass(cutClass);
948       THashList *classTable=0;
949       while ( (classTable=(THashList*)nextClass()) ){
950         TIter nextHist(classTable);
951         TObject *o=0;
952         printf("|  + %s\n",classTable->GetName());
953         while ( (o=nextHist()) )
954           printf("|  | ->%s\n",o->GetName());
955       }
956       
957     }
958   }
959 }
960
961 //_____________________________________________________________________________
962 void AliDielectronHistos::SetHistogramList(THashList &list, Bool_t setOwner/*=kTRUE*/)
963 {
964   //
965   // set histogram classes and histograms to this instance. It will take onwnership!
966   //
967   ResetHistogramList();
968   TString name(GetName());
969   if (name == "AliDielectronHistos") SetName(list.GetName());
970   TIter next(&list);
971   TObject *o;
972   while ( (o=next()) ){
973     fHistoList.Add(o);
974   }
975   if (setOwner){
976     list.SetOwner(kFALSE);
977     fHistoList.SetOwner(kTRUE);
978   } else {
979     fHistoList.SetOwner(kFALSE);
980   }
981 }
982
983 //_____________________________________________________________________________
984 Bool_t AliDielectronHistos::SetCutClass(const char* cutClass)
985 {
986   //
987   // Assign histogram list according to cutClass
988   //
989
990   if (!fList) return kFALSE;
991   ResetHistogramList();
992   THashList *h=dynamic_cast<THashList*>(fList->FindObject(cutClass));
993   if (!h) {
994     Warning("SetCutClass","cutClass '%s' not found", cutClass);
995     return kFALSE;
996   }
997   SetHistogramList(*h,kFALSE);
998   return kTRUE;
999 }
1000
1001 //_____________________________________________________________________________
1002 Bool_t AliDielectronHistos::IsHistogramOk(const char* histClass, const char* name)
1003 {
1004   //
1005   // check whether the histogram class exists and the histogram itself does not exist yet
1006   //
1007   Bool_t isReserved=fReservedWords->Contains(histClass);
1008   if (!fHistoList.FindObject(histClass)&&!isReserved){
1009     Warning("IsHistogramOk","Cannot create histogram. Class '%s' not defined. Please create it using AddClass before.",histClass);
1010     return kFALSE;
1011   }
1012   if (GetHist(histClass,name)){
1013     Warning("IsHistogramOk","Cannot create histogram '%s' in class '%s': It already exists!",name,histClass);
1014     return kFALSE;
1015   }
1016   return kTRUE;
1017 }
1018
1019 // //_____________________________________________________________________________
1020 // TIterator* AliDielectronHistos::MakeIterator(Bool_t dir) const
1021 // {
1022 //   //
1023 //   //
1024 //   //
1025 //   return new TListIter(&fHistoList, dir);
1026 // }
1027
1028 //_____________________________________________________________________________
1029 void AliDielectronHistos::ReadFromFile(const char* file)
1030 {
1031   //
1032   // Read histos from file
1033   //
1034   TFile f(file);
1035   TIter nextKey(f.GetListOfKeys());
1036   TKey *key=0;
1037   while ( (key=(TKey*)nextKey()) ){
1038     TObject *o=f.Get(key->GetName());
1039     THashList *list=dynamic_cast<THashList*>(o);
1040     if (!list) continue;
1041     SetHistogramList(*list);
1042     break;
1043   }
1044   f.Close();
1045 }
1046
1047 //_____________________________________________________________________________
1048 void AliDielectronHistos::DrawSame(const char* histName, const Option_t *opt)
1049 {
1050   //
1051   // Draw all histograms with the same name into one canvas
1052   // if option contains 'leg' a legend will be created with the class name as caption
1053   // if option contains 'can' a new canvas is created
1054   //
1055
1056   TString optString(opt);
1057   optString.ToLower();
1058   Bool_t optLeg=optString.Contains("leg");
1059   Bool_t optCan=optString.Contains("can");
1060
1061   TLegend *leg=0;
1062   TCanvas *c=0;
1063   if (optCan){
1064     c=(TCanvas*)gROOT->FindObject(Form("c%s",histName));
1065     if (!c) c=new TCanvas(Form("c%s",histName),Form("All '%s' histograms",histName));
1066     c->Clear();
1067     c->cd();
1068   }
1069
1070   if (optLeg) leg=new TLegend(.8,.3,.99,.9);
1071   
1072   Int_t i=0;
1073   TIter next(&fHistoList);
1074   THashList *classTable=0;
1075   Double_t max=-1e10;
1076   TH1 *hFirst=0x0;
1077   while ( (classTable=(THashList*)next()) ){
1078     if ( TH1 *h=(TH1*)classTable->FindObject(histName) ){
1079       if (i==0) hFirst=h;
1080       h->SetLineColor(i+1);
1081       h->SetMarkerColor(i+1);
1082       h->Draw(i>0?"same":"");
1083       if (leg) leg->AddEntry(h,classTable->GetName(),"lp");
1084       ++i;
1085       max=TMath::Max(max,h->GetMaximum());
1086     }
1087   }
1088   if (leg){
1089     leg->SetFillColor(10);
1090     leg->SetY1(.9-i*.05);
1091     leg->Draw();
1092   }
1093   if (hFirst&&(hFirst->GetYaxis()->GetXmax()<max)){
1094     hFirst->SetMaximum(max);
1095   }
1096 }
1097
1098 //_____________________________________________________________________________
1099 void AliDielectronHistos::SetReservedWords(const char* words)
1100 {
1101   //
1102   // set reserved words
1103   //
1104
1105   (*fReservedWords)=words;
1106 }
1107
1108 //_____________________________________________________________________________
1109 void AliDielectronHistos::StoreVariables(TObject *obj, UInt_t valType[20])
1110 {
1111   //
1112   //
1113   //
1114   if (!obj) return;
1115   if      (obj->InheritsFrom(TH1::Class()))         StoreVariables(static_cast<TH1*>(obj), valType);
1116   else if (obj->InheritsFrom(THnBase::Class()))     StoreVariables(static_cast<THnBase*>(obj), valType);
1117
1118   return;
1119
1120 }
1121
1122
1123 //_____________________________________________________________________________
1124 void AliDielectronHistos::StoreVariables(TH1 *obj, UInt_t valType[20])
1125 {
1126   //
1127   // store variables in the axis (special for TProfile3D)
1128   //
1129
1130   Int_t dim   = obj->GetDimension();
1131
1132   // dimension correction for profiles
1133   if(obj->IsA() == TProfile::Class() || obj->IsA() == TProfile2D::Class() || obj->IsA() == TProfile3D::Class()) {
1134     dim++;
1135   }
1136
1137   switch( dim ) {
1138   case 4:
1139     obj->SetUniqueID(valType[3]); // Tprofile3D variable
1140   case 3:
1141     obj->GetZaxis()->SetUniqueID(valType[2]);
1142   case 2:
1143     obj->GetYaxis()->SetUniqueID(valType[1]);
1144   case 1:
1145     obj->GetXaxis()->SetUniqueID(valType[0]);
1146   }
1147
1148   return;
1149 }
1150
1151 //_____________________________________________________________________________
1152 void AliDielectronHistos::StoreVariables(THnBase *obj, UInt_t valType[20])
1153 {
1154   //
1155   // store variables in the axis
1156   //
1157
1158   Int_t dim = obj->GetNdimensions();
1159
1160   for(Int_t it=0; it<dim; it++) {
1161     obj->GetAxis(it)->SetUniqueID(valType[it]);
1162     obj->GetAxis(it)->SetName(Form("%s", AliDielectronVarManager::GetValueName(valType[it])));
1163     obj->GetAxis(it)->SetTitle(Form("%s %s", AliDielectronVarManager::GetValueLabel(valType[it]), AliDielectronVarManager::GetValueUnit(valType[it])));
1164   }
1165   obj->Sumw2();
1166   return;
1167 }
1168
1169 //_____________________________________________________________________________
1170 void AliDielectronHistos::FillValues(TObject *obj, const Double_t *values)
1171 {
1172   //
1173   //
1174   //
1175   if (!obj) return;
1176   if      (obj->InheritsFrom(TH1::Class()))       FillValues(static_cast<TH1*>(obj), values);
1177   else if (obj->InheritsFrom(THnBase::Class()))   FillValues(static_cast<THnBase*>(obj), values);
1178
1179   return;
1180
1181 }
1182
1183 //_____________________________________________________________________________
1184 void AliDielectronHistos::FillValues(TH1 *obj, const Double_t *values)
1185 {
1186   //
1187   // fill values for TH1 inherted classes
1188   //
1189
1190   Int_t dim   = obj->GetDimension();
1191   Bool_t bprf = kFALSE;
1192   //  UInt_t nValues = (UInt_t) AliDielectronVarManager::kNMaxValues;
1193
1194   UInt_t valueTypes=obj->GetUniqueID();
1195   if (valueTypes==(UInt_t)AliDielectronHistos::kNoAutoFill) return;
1196   Bool_t weight = (valueTypes!=kNoWeights);
1197
1198   if(obj->IsA() == TProfile::Class() || obj->IsA() == TProfile2D::Class() || obj->IsA() == TProfile3D::Class())
1199     bprf=kTRUE;
1200
1201   UInt_t value1=obj->GetXaxis()->GetUniqueID();
1202   UInt_t value2=obj->GetYaxis()->GetUniqueID();
1203   UInt_t value3=obj->GetZaxis()->GetUniqueID();
1204   UInt_t value4=obj->GetUniqueID();            // get profile var stored in the unique ID
1205
1206   // ask for inclusive trigger map variables
1207   if(value1!=AliDielectronVarManager::kTriggerInclONL && value1!=AliDielectronVarManager::kTriggerInclOFF &&
1208      value2!=AliDielectronVarManager::kTriggerInclONL && value2!=AliDielectronVarManager::kTriggerInclOFF &&
1209      value3!=AliDielectronVarManager::kTriggerInclONL && value3!=AliDielectronVarManager::kTriggerInclOFF &&
1210      value4!=AliDielectronVarManager::kTriggerInclONL && value4!=AliDielectronVarManager::kTriggerInclOFF ) {
1211     // no trigger map variable selected
1212     switch ( dim ) {
1213     case 1:
1214       if(!bprf && !weight)     obj->Fill(values[value1]);                 // histograms
1215       else if(!bprf && weight) obj->Fill(values[value1], values[value4]); // weighted histograms
1216       else if(bprf && !weight) ((TProfile*)obj)->Fill(values[value1],values[value2]);   // profiles
1217       else                     ((TProfile*)obj)->Fill(values[value1],values[value2], values[value4]); // weighted profiles
1218       break;
1219     case 2:
1220       if(!bprf && !weight)     obj->Fill(values[value1], values[value2]);                 // histograms
1221       else if(!bprf && weight) ((TH2*)obj)->Fill(values[value1], values[value2], values[value4]); // weighted histograms
1222       else if(bprf && !weight) ((TProfile2D*)obj)->Fill(values[value1], values[value2], values[value3]); // profiles
1223       else                     ((TProfile2D*)obj)->Fill(values[value1], values[value2], values[value3], values[value4]); // weighted profiles
1224       break;
1225     case 3:
1226       if(!bprf && !weight)     ((TH3*)obj)->Fill(values[value1], values[value2], values[value3]);                 // histograms
1227       else if(!bprf && weight) ((TH3*)obj)->Fill(values[value1], values[value2], values[value3], values[value4]); // weighted histograms
1228       else if(bprf && !weight) ((TProfile3D*)obj)->Fill(values[value1], values[value2], values[value3], values[value4]); // profiles
1229       else                     printf(" WARNING: weighting NOT possible yet for TProfile3Ds ! \n");
1230       break;
1231     }
1232   }
1233   else {
1234     // fill inclusive trigger map variables
1235     if(weight) return;
1236     switch ( dim ) {
1237     case 1:
1238       for(Int_t i=0; i<30; i++) { if(TESTBIT((UInt_t)values[value1],i)) obj->Fill(i); }
1239       break;
1240     case 2:
1241       if((value1==AliDielectronVarManager::kTriggerInclOFF && value2==AliDielectronVarManager::kTriggerInclONL) ||
1242          (value1==AliDielectronVarManager::kTriggerInclONL && value2==AliDielectronVarManager::kTriggerInclOFF) ) {
1243         for(Int_t i=0; i<30; i++) {
1244           if((UInt_t)values[value1]==BIT(i)) {
1245             for(Int_t i2=0; i2<30; i2++) {
1246               if((UInt_t)values[value2]==BIT(i2)) {
1247                 obj->Fill(i, i2);
1248               } // bit fired
1249             } //loop 2
1250           }//bit fired
1251         } // loop 1
1252       }
1253       else if(value1==AliDielectronVarManager::kTriggerInclONL || value1==AliDielectronVarManager::kTriggerInclOFF) {
1254         for(Int_t i=0; i<30; i++) { if(TESTBIT((UInt_t)values[value1],i)) obj->Fill(i, values[value2]); }
1255       }
1256       else if(value2==AliDielectronVarManager::kTriggerInclONL || value2==AliDielectronVarManager::kTriggerInclOFF) {
1257         for(Int_t i=0; i<30; i++) { if(TESTBIT((UInt_t)values[value2],i)) obj->Fill(values[value1], i); }
1258       }
1259       else //makes no sense
1260         return;
1261       break;
1262     default: return;
1263     }
1264
1265   } //end: trigger filling
1266
1267
1268   return;
1269 }
1270
1271 //_____________________________________________________________________________
1272 void AliDielectronHistos::FillValues(THnBase *obj, const Double_t *values)
1273 {
1274   //
1275   // fill values for THn inherted classes
1276   //
1277
1278   const Int_t dim   = obj->GetNdimensions();
1279
1280   UInt_t valueTypes=obj->GetUniqueID();
1281   if (valueTypes==(UInt_t)AliDielectronHistos::kNoAutoFill) return;
1282   Bool_t weight = (valueTypes!=kNoWeights);
1283
1284   UInt_t value4=obj->GetUniqueID();            // weight variable
1285
1286   Double_t fill[dim];
1287   for(Int_t it=0; it<dim; it++)   fill[it] = values[obj->GetAxis(it)->GetUniqueID()];
1288   if(!weight) obj->Fill(fill);
1289   else obj->Fill(fill, values[value4]);
1290
1291
1292   return;
1293 }
1294
1295 //_____________________________________________________________________________
1296 void AliDielectronHistos::FillVarArray(TObject *obj, UInt_t *valType)
1297 {
1298   //
1299   // extract variables stored in the axis (special for TProfile3D)
1300   //
1301
1302
1303   if (!obj) return;
1304   //  printf(" fillvararray %s \n",obj->GetName());
1305
1306   if (obj->InheritsFrom(TH1::Class())) {
1307     valType[0]=((TH1*)obj)->GetXaxis()->GetUniqueID();
1308     valType[1]=((TH1*)obj)->GetYaxis()->GetUniqueID();
1309     valType[2]=((TH1*)obj)->GetZaxis()->GetUniqueID();
1310     valType[3]=((TH1*)obj)->GetUniqueID();  // tprofile var stored in unique ID
1311   }
1312   else if (obj->InheritsFrom(THnBase::Class())) {
1313     for(Int_t it=0; it<((THn*)obj)->GetNdimensions(); it++)
1314       valType[it]=((THn*)obj)->GetAxis(it)->GetUniqueID();
1315   }
1316   valType[19]=obj->GetUniqueID(); //weights
1317   return;
1318 }
1319
1320 //_____________________________________________________________________________
1321 void AliDielectronHistos::AdaptNameTitle(TH1 *hist, const char* histClass) {
1322
1323   //
1324   // adapt name and title of the histogram
1325   //
1326
1327   Int_t dim            = hist->GetDimension();
1328   TString currentName  = hist->GetName();
1329   TString currentTitle = hist->GetTitle();
1330
1331
1332   Bool_t bname  = (currentName.IsNull());
1333   Bool_t btitle = (currentTitle.IsNull());
1334   Bool_t bprf   = kFALSE;
1335   if(hist->IsA() == TProfile::Class() || hist->IsA() == TProfile2D::Class() || hist->IsA() == TProfile3D::Class())
1336     bprf=kTRUE;
1337
1338   // tprofile options
1339   Double_t pmin=0., pmax=0.;
1340   TString option = "", calcrange="";
1341   Bool_t bStdOpt=kTRUE;
1342   if(bprf) {
1343     switch( dim ) {
1344     case 3:
1345       option = ((TProfile3D*)hist)->GetErrorOption();
1346       pmin   = ((TProfile3D*)hist)->GetTmin();
1347       pmax   = ((TProfile3D*)hist)->GetTmax();
1348       break;
1349     case 2:
1350       option = ((TProfile2D*)hist)->GetErrorOption();
1351       pmin   = ((TProfile2D*)hist)->GetZmin();
1352       pmax   = ((TProfile2D*)hist)->GetZmax();
1353       break;
1354     case 1:
1355       option = ((TProfile*)hist)->GetErrorOption();
1356       pmin   = ((TProfile*)hist)->GetYmin();
1357       pmax   = ((TProfile*)hist)->GetYmax();
1358       break;
1359     }
1360     if(option.Contains("s",TString::kIgnoreCase)) bStdOpt=kFALSE;
1361     if(pmin!=pmax) calcrange=Form("#cbar_{%+.*f}^{%+.*f}",GetPrecision(pmin),pmin,GetPrecision(pmax),pmax);
1362   }
1363
1364   UInt_t varx = hist->GetXaxis()->GetUniqueID();
1365   UInt_t vary = hist->GetYaxis()->GetUniqueID();
1366   UInt_t varz = hist->GetZaxis()->GetUniqueID();
1367   UInt_t varp = hist->GetUniqueID();
1368   Bool_t weight = (varp!=kNoWeights);
1369
1370   // store titles in the axis
1371   if(btitle) {
1372     switch( dim ) {
1373     case 3:
1374       hist->GetXaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varx),
1375                                      Form("%s %s",
1376                                       AliDielectronVarManager::GetValueLabel(varx),
1377                                       AliDielectronVarManager::GetValueUnit(varx))
1378                                  );
1379       hist->GetYaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(vary),
1380                                      Form("%s %s",
1381                                       AliDielectronVarManager::GetValueLabel(vary),
1382                                       AliDielectronVarManager::GetValueUnit(vary))
1383                                  );
1384       hist->GetZaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varz),
1385                                      Form("%s %s",
1386                                       AliDielectronVarManager::GetValueLabel(varz),
1387                                       AliDielectronVarManager::GetValueUnit(varz))
1388                                  );
1389       if(bprf)
1390         hist->SetNameTitle(AliDielectronVarManager::GetValueName(varp),
1391                            Form("%s  %s%s%s%s %s",
1392                             hist->GetTitle(),
1393                             (bStdOpt ? "#LT" : "RMS("),
1394                             AliDielectronVarManager::GetValueLabel(varp),
1395                             (bStdOpt ? "#GT" : ")"),
1396                             calcrange.Data(),
1397                             AliDielectronVarManager::GetValueUnit(varp))
1398                        );
1399       break;
1400     case 2:
1401       hist->GetXaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varx),
1402                                      Form("%s %s",
1403                                       AliDielectronVarManager::GetValueLabel(varx),
1404                                       AliDielectronVarManager::GetValueUnit(varx))
1405                                  );
1406       hist->GetYaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(vary),
1407                                      Form("%s %s",
1408                                       AliDielectronVarManager::GetValueLabel(vary),
1409                                       AliDielectronVarManager::GetValueUnit(vary))
1410                                  );
1411       hist->GetZaxis()->SetTitle(Form("#%ss",histClass));
1412       if(bprf)
1413         hist->GetZaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varz),
1414                                        Form("%s%s%s%s %s",
1415                                         (bStdOpt ? "#LT" : "RMS("),
1416                                         AliDielectronVarManager::GetValueLabel(varz),
1417                                         (bStdOpt ? "#GT" : ")"),
1418                                         calcrange.Data(),
1419                                         AliDielectronVarManager::GetValueUnit(varz))
1420                                    );
1421       if(weight) hist->GetZaxis()->SetTitle(Form("%s weighted %s",
1422                                                  AliDielectronVarManager::GetValueLabel(varp),
1423                                                  hist->GetZaxis()->GetTitle() )
1424                                             );
1425
1426       break;
1427     case 1:
1428       hist->GetXaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(varx),
1429                                      Form("%s %s",
1430                                       AliDielectronVarManager::GetValueLabel(varx),
1431                                       AliDielectronVarManager::GetValueUnit(varx))
1432                                  );
1433       hist->GetYaxis()->SetTitle(Form("#%ss",histClass));
1434       if(bprf)
1435         hist->GetYaxis()->SetNameTitle(AliDielectronVarManager::GetValueName(vary),
1436                                        Form("%s%s%s%s %s",
1437                                         (bStdOpt ? "#LT" : "RMS("),
1438                                         AliDielectronVarManager::GetValueLabel(vary), 
1439                                         (bStdOpt ? "#GT" : ")"),
1440                                         calcrange.Data(),
1441                                         AliDielectronVarManager::GetValueUnit(vary))
1442                                    );
1443       if(weight) hist->GetYaxis()->SetTitle(Form("%s weighted %s",
1444                                                  AliDielectronVarManager::GetValueLabel(varp),
1445                                                  hist->GetYaxis()->GetTitle() )
1446                                             );
1447
1448       break;
1449     }
1450
1451     // create an unique name
1452     if(bname)
1453       switch(dim) {
1454       case 3:
1455         currentName+=Form("%s_",AliDielectronVarManager::GetValueName(varx));
1456         currentName+=Form("%s_",AliDielectronVarManager::GetValueName(vary));
1457         currentName+=Form("%s",AliDielectronVarManager::GetValueName(varz));
1458         if(bprf)   currentName+=Form("-%s%s",AliDielectronVarManager::GetValueName(varp),(bStdOpt ? "avg" : "rms"));
1459         if(weight) currentName+=Form("-wght%s",AliDielectronVarManager::GetValueName(varp));
1460         break;
1461       case 2:
1462         currentName+=Form("%s_",AliDielectronVarManager::GetValueName(varx));
1463         currentName+=Form("%s",AliDielectronVarManager::GetValueName(vary));
1464         if(bprf)   currentName+=Form("-%s%s",AliDielectronVarManager::GetValueName(varz),(bStdOpt ? "avg" : "rms"));
1465         if(weight) currentName+=Form("-wght%s",AliDielectronVarManager::GetValueName(varp));
1466         break;
1467       case 1:
1468         currentName+=Form("%s",AliDielectronVarManager::GetValueName(varx));
1469         if(bprf)   currentName+=Form("-%s%s",AliDielectronVarManager::GetValueName(vary),(bStdOpt ? "avg" : "rms"));
1470         if(weight) currentName+=Form("-wght%s",AliDielectronVarManager::GetValueName(varp));
1471         break;
1472       }
1473     // to differentiate btw. leg and pair histos
1474     if(!strcmp(histClass,"Pair")) currentName.Prepend("p");
1475     hist->SetName(currentName.Data());
1476   }
1477
1478 }
1479
1480 Int_t AliDielectronHistos::GetPrecision(Double_t value) {
1481
1482   //
1483   // computes the precision of a double
1484   // usefull for axis ranges etc
1485   //
1486
1487   Bool_t bfnd     = kFALSE;
1488   Int_t precision = 0;
1489   value*=1e6;
1490   while(!bfnd) {
1491     //    printf(" value %f precision %d bfnd %d \n",TMath::Abs(value*TMath::Power(10,precision)), precision, bfnd);
1492     bfnd = ((TMath::Abs(value*TMath::Power(10,precision))/1e6  -  TMath::Floor(TMath::Abs(value*TMath::Power(10,precision))/1e6)) != 0.0
1493             ? kFALSE
1494             : kTRUE);
1495     if(!bfnd) precision++;
1496   }
1497
1498   //  printf("precision for %f found to be %d \n", value, precision);
1499   return precision;
1500
1501 }
1502
1503
1504