]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliHFEcollection.cxx
Changes for #82873: Module debugging broken (Christian)
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliHFEcollection.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, 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 // Collection class for histograms
17 // Stores either histograms or vectors of histograms
18 // 
19 // Author:
20 // Matus Kalisky  <matus.kalisky@cern.ch>
21 //
22
23 //#include <iostream>
24
25 #include <TH1F.h>
26 #include <TH2F.h>
27 #include <THnSparse.h>
28 #include <TProfile.h>
29 #include <TString.h>
30 #include <TBrowser.h>
31 #include <TMath.h>
32
33 #include "AliLog.h"
34 #include "AliHFEcollection.h"
35
36 using namespace std;
37
38
39 ClassImp(AliHFEcollection)
40
41 //___________________________________________________________________
42 AliHFEcollection::AliHFEcollection():
43   TNamed()
44   , fList(NULL)
45 {
46
47   //
48   // default constructor
49   //
50 }
51 //___________________________________________________________________
52 AliHFEcollection::AliHFEcollection(const char* name, const char* title):
53   TNamed(name, title)
54   , fList(NULL)
55 {
56  
57   //
58   // constructor
59   //
60  
61   fList = new THashList();
62   if(fList){
63     fList->SetOwner();
64     fList->SetName(Form("list_%s", name));
65   }
66 }
67 //___________________________________________________________________
68 AliHFEcollection::AliHFEcollection(const AliHFEcollection &c) :
69   TNamed(c)
70   , fList(NULL)
71 {
72
73   //
74   // copy operator
75   // 
76   
77   c.Copy(*this);
78 }
79 //___________________________________________________________________
80 AliHFEcollection &AliHFEcollection::operator=(const AliHFEcollection &ref)
81 {
82   //
83   // Assignment operator
84   //
85   
86   if(this != &ref){
87     ref.Copy(*this);
88   }
89   return *this;
90 }
91 //___________________________________________________________________
92 void AliHFEcollection::Copy(TObject &ref) const {
93
94   //
95   // Performs the copying of the object
96   //
97
98   AliHFEcollection &target = dynamic_cast<AliHFEcollection &>(ref);
99
100   // Clone List Content
101   target.fList = new THashList();          
102   for(Int_t ien = 0; ien < fList->GetEntries(); ien++)
103     target.fList->Add(fList->At(ien)->Clone());
104 }
105 //___________________________________________________________________
106 AliHFEcollection::~AliHFEcollection(){
107
108   //
109   // Destructor
110   //
111   delete fList;
112   AliDebug(1, "DESTRUCTOR");
113 }
114 //___________________________________________________________________
115 Bool_t AliHFEcollection::CreateTH1F(const char* name, const char* title, Int_t nBin, Float_t nMin, Float_t nMax, Int_t logAxis){
116
117   //
118   // Creates a TH1F histogram for the collection
119   //
120
121   if(!fList){
122     AliError("No TList pointer ! ");
123     return kFALSE;
124   }
125   else{
126     fList->Add(new TH1F(name, title, nBin, nMin, nMax));
127     if(logAxis >= 0){
128       BinLogAxis(name, logAxis);
129     }
130     return CheckObject(name);
131   }
132 }
133
134 //___________________________________________________________________
135 Bool_t AliHFEcollection::CreateTH1Farray(const char* name, const char* title, Int_t nBin, const Double_t* xbins){
136
137   //
138   // Creates a TH1F histogram for the collection 2nd version 
139   //
140
141   if(!fList){
142     AliError("No TList pointer ! ");
143     return kFALSE;
144   }
145   else{
146     fList->Add(new TH1F(name, title, nBin, xbins));
147     return CheckObject(name);
148   }
149 }
150
151 //___________________________________________________________________
152 Bool_t AliHFEcollection::CreateTH2F(const char* name, const char* title, Int_t nBinX, Float_t nMinX, Float_t nMaxX, Int_t nBinY, Float_t nMinY, Float_t nMaxY, Int_t logAxis){
153
154   //
155   // Creates a TH2F histogram for the collection
156   //
157
158   if(!fList){
159     AliError("No TList pointer ! ");
160     return kFALSE;
161   }
162   fList->Add(new TH2F(name, title, nBinX, nMinX, nMaxX, nBinY, nMinY, nMaxY));
163   if(logAxis >= 0){
164     BinLogAxis(name, logAxis);
165   }
166   return CheckObject(name); 
167 }
168 //___________________________________________________________________
169 Bool_t AliHFEcollection::CreateTH1Fvector1(Int_t X, const char* name, const char* title, Int_t nBin, Float_t nMin, Float_t nMax, Int_t logAxis){
170
171   //
172   // create a 1 dimensional array of size [X]
173   //
174
175   if(!fList){
176     AliError("No TList pointer ! ");
177     return kFALSE;
178   }
179   if(X <=0){
180     AliError("can not create array with negative or zero size ");
181     return kFALSE;
182   }
183   TString hname;
184   for(Int_t i=0; i<X; ++i){
185     hname = "";
186     hname.Append(Form("%s_[%d]", name, i));
187     //cout<<" -D: name: "<<name.str().c_str()<<endl;
188     //cout<<" -D: nBin: "<<_nBin<<" ,Min: "<<_nMin<<" , Max: "<<_nMax<<endl;
189     CreateTH1F(hname.Data(), title, nBin, nMin, nMax, logAxis);
190     if(!CheckObject(hname.Data())){
191       AliError(Form("Not possible to create object: %s", hname.Data()));
192       return kFALSE;
193     }    
194   }
195   return kTRUE;  
196 }
197 //___________________________________________________________________
198 Bool_t AliHFEcollection::CreateTH2Fvector1(Int_t X, const char* name, const char* title, Int_t nBinX, Float_t nMinX, Float_t nMaxX, Int_t nBinY, Float_t nMinY, Float_t nMaxY, Int_t logAxis){
199
200   //
201   // create a 1 dimensinal array of TH2F histograms with size [X]
202   //
203
204   if(!fList){
205     AliError("No TList pointer !");
206     return kFALSE;
207   }
208   if(X <=0){
209     AliError("can not create array with negative or zero size ");
210     return kFALSE;
211   }
212   TString hname;
213   for(Int_t i=0; i<X; ++i){
214     hname = "";
215     hname.Append(Form("%s_[%d]", name, i));
216     //cout<<" -D: name: "<<name<<endl;
217     //cout<<" -D: nBin: "<<_nBin<<" ,Min: "<<_nMin<<" , Max: "<<_nMax<<endl;
218     CreateTH2F(hname.Data(), title, nBinX, nMinX, nMaxX, nBinY, nMinY, nMaxY, logAxis);
219     if(!CheckObject(hname.Data())){
220       AliError(Form("Not possible to create object: %s", hname.Data()));
221       return kFALSE;
222     }    
223   }
224   return kTRUE;  
225 }
226 //___________________________________________________________________
227 Bool_t AliHFEcollection::CreateTH1Fvector2(Int_t X, Int_t Y, const char* name, const char* title, Int_t nBin, Float_t nMin, Float_t nMax, Int_t logAxis){
228
229   //
230   // create a 2 dimensional array of histograms of size [X, Y]
231   //
232
233   if(!fList){
234     AliError("No TList pointer ! ");
235     return kFALSE;
236   }
237   if(X <=0 || Y <=0){
238     AliError("can not create array with negative or zero size ");
239     return kFALSE;
240   }
241   TString hname;
242   for(Int_t i=0; i<X; ++i){
243     for(Int_t j=0; j<Y; ++j){
244       hname = "";
245       hname.Append(Form("%s_[%d][%d]", name, i, j));
246       //cout<<" -D: name: "<<name.str().c_str()<<endl;
247       //cout<<" -D: nBin: "<<_nBin<<" ,Min: "<<_nMin<<" , Max: "<<_nMax<<endl;
248       CreateTH1F(hname.Data(), title, nBin, nMin, nMax, logAxis);
249       if(!CheckObject(hname.Data())){
250               AliError(Form("Not possible to create object: %s", hname.Data()));
251               return kFALSE;
252       }
253     }
254   }
255   return kTRUE;  
256 }
257 //___________________________________________________________________
258 Bool_t AliHFEcollection::CreateProfile(const char* name, const char* title, Int_t nbins, Double_t xmin, Double_t xmax){
259   
260   //
261   // create a simple TProfile
262   //
263
264   if(!fList){
265     AliError("No TList pointer ! ");
266     return kFALSE;
267   }
268   fList->Add(new TProfile(name, title, nbins, xmin, xmax));
269   return CheckObject(name);
270
271 }
272 //___________________________________________________________________
273 Bool_t AliHFEcollection::CreateTHnSparse(const char* name, const char* title, Int_t dim, Int_t* nbins, Double_t* xmin, Double_t* xmax){
274
275   //
276   // create 'dim' dimensional THnSparse
277   //
278
279   if(!fList){
280     AliError("No TList pointer ! ");
281     return kFALSE;
282   }
283   fList->Add(new THnSparseF(name, title, dim, nbins, xmin, xmax));
284   return CheckObject(name);
285
286 }
287 //___________________________________________________________________
288 TObject* AliHFEcollection::Get(const char* name){ 
289
290   //
291   // Get histogram with the required name
292   // 
293   
294
295   if(!CheckObject(name)){
296     AliWarning(Form("Not possible to return pointer to the object '%s'\n", name));
297     return 0;
298   }
299
300   return fList->FindObject(name);
301   
302 }
303 //___________________________________________________________________
304 Bool_t AliHFEcollection::Fill(const char* name, Double_t v){
305
306   //
307   // fill function for one TH1 histograms
308   //
309
310   if(!CheckObject(name)){
311     AliError(Form("Not possible to fill the object '%s', the object does not exist\n", name));
312     return kFALSE;
313   }
314
315   TH1 *htmp = dynamic_cast<TH1F*>(fList->FindObject(name));
316   // chack the possible object types
317   if(htmp){
318     htmp->Fill(v);
319     return kTRUE;
320   }
321   
322   return kFALSE;
323
324 }
325 //___________________________________________________________________
326 Bool_t AliHFEcollection::Fill(const char* name, Int_t v){
327
328   //
329   // fill function for one TH1 histograms for integer numbers
330   //
331
332   return Fill(name, v*1.0);
333 }
334 //___________________________________________________________________
335 Bool_t AliHFEcollection::Fill(const char* name, Int_t X, Double_t v){
336
337   //
338   // fill function for one dimension arrays of TH1
339   //
340
341   const char* n = Form("%s_[%d]", name, X);
342   TObject *o = Get(n);
343   if(!o){
344     return kFALSE;
345   }
346   Fill(o->GetName(), v);
347   return kTRUE;
348 }
349 //___________________________________________________________________
350 Bool_t AliHFEcollection::Fill(const char* name, Int_t X, Int_t Y, Double_t v){
351
352   //
353   // Fill function fir 2 dimensional TH1 arrays
354   //
355   
356   const char* n = Form("%s_[%d][%d]", name, X, Y);
357   TObject *o = Get(n);
358   if(!o){
359     return kFALSE;
360   }
361   Fill(o->GetName(), v);
362   return kTRUE;
363 }
364 //___________________________________________________________________
365 Bool_t AliHFEcollection::Fill(const char* name, Int_t X, Double_t v1, Double_t v2){
366
367   //
368   // fill function for one dimension array of TH2
369   //
370
371   const char* n = Form("%s_[%d]", name, X);
372   TObject *o = Get(n);
373   if(!o){
374     return kFALSE;
375   }
376   Fill(o->GetName(), v1, v2);
377   
378   return kTRUE;
379 }
380 //___________________________________________________________________
381 Bool_t AliHFEcollection::Fill(const char* name, Double_t v1, Double_t v2){
382
383   //
384   // fill function for TH2 objects
385   //
386
387   if(!CheckObject(name)){
388     AliError(Form("Not possible to fill the object '%s', the object does not exist\n", name));
389     return kFALSE;
390   }
391
392   // chack the possible object types
393   if(fList->FindObject(name)->InheritsFrom("TH2")){
394     TH2 *h2 = dynamic_cast<TH2F*>(fList->FindObject(name));
395     if(h2) h2->Fill(v1, v2);
396     return kTRUE;
397   }  
398   if(fList->FindObject(name)->InheritsFrom("TProfile")){
399     TProfile *pr = dynamic_cast<TProfile*>(fList->FindObject(name));
400     if(pr) pr->Fill(v1, v2);
401     return kTRUE;
402   }  
403   
404   return kFALSE;
405   
406 }
407 //___________________________________________________________________
408 Bool_t AliHFEcollection::Fill(const char* name, Double_t* entry, Double_t weight){
409   //
410   // Fill a THnSparse object
411   //
412
413   if(!CheckObject(name)){
414     AliError(Form("Not possible to fill the object '%s', the object does not exist\n", name));
415     return kFALSE;
416   }
417   
418   THnSparseF *htmp = dynamic_cast<THnSparseF*>(fList->FindObject(name));
419    if(htmp){
420      htmp->Fill(entry, weight);
421      return kTRUE;
422    }
423    return kFALSE;
424
425 }
426 //___________________________________________________________________
427 Bool_t AliHFEcollection::CheckObject(const char* name){
428
429   //
430   // check wheter the creation of the histogram was succesfull
431   //
432   
433   if(!fList){
434     AliError("No TList pointer ! ");
435     return kFALSE;
436   }
437   
438   if(!fList->FindObject(name)){
439     AliWarning(Form("Creating or Finding the object '%s' failed\n", name));
440     return kFALSE;
441   }
442   return kTRUE;
443 }
444 //___________________________________________________________________
445 Bool_t AliHFEcollection::Sumw2(const char* name){
446   //
447   // Set Sumw2 for the given object
448   //
449   if(!CheckObject(name)){
450     return kFALSE;
451   }
452
453   TObject *o = Get(name);
454   THnSparse *htmp = dynamic_cast<THnSparse*>(o);
455   if(htmp){
456     htmp->Sumw2();
457   }
458   return kTRUE;
459 }
460 //___________________________________________________________________
461 Bool_t AliHFEcollection::BinLogAxis(const char* name, Int_t dim){
462
463   // 
464   // converts the axis (defined by the dimension) of THx or THnSparse
465   // object to Log scale. Number of bins and bin min and bin max are preserved
466   //
467
468
469   if(!CheckObject(name)){
470     return kFALSE;
471   }
472
473   TObject *o = Get(name);
474   TAxis *axis = NULL;
475   if(o->InheritsFrom("TH1")){
476     TH1 *h1 = dynamic_cast<TH1F*>(o);
477     if(h1) axis = h1->GetXaxis();
478   }
479   if(o->InheritsFrom("TH2")){
480     TH2 *h2 = dynamic_cast<TH2F*>(o);
481     if(h2){
482       if(0 == dim){
483         axis = h2->GetXaxis();
484       }
485       else if(1 == dim){
486         axis = h2->GetYaxis();
487       }
488       else{
489          AliError("Only dim = 0 or 1 possible for TH2F");
490       }
491     }
492   }
493   if(o->InheritsFrom("THnSparse")){
494     THnSparse *hs = dynamic_cast<THnSparse*>(o);
495     if(hs) axis = hs->GetAxis(dim);
496   }
497   
498   if(!axis){
499     AliError(Form("Axis '%d' could not be identified in the object '%s'\n", dim, name));
500     return kFALSE;
501   }
502   
503   Int_t bins = axis->GetNbins();
504
505   Double_t from = axis->GetXmin();
506   if(from <= 0){
507     AliError(Form(" Log binning not possible for object '%s'because the '%d' axis starts from '%f\n'", name, dim, from));
508     return kFALSE;
509   }
510   Double_t to = axis->GetXmax();
511   Double_t *newBins = new Double_t[bins+1];
512   newBins[0] = from;
513   Double_t factor = TMath::Power(to/from, 1./bins);
514   for(Int_t i=1; i<=bins; ++i){
515     newBins[i] = factor * newBins[i-1];
516   }
517   axis->Set(bins, newBins);
518   delete[] newBins;
519
520   return kTRUE;
521
522
523 }
524 //___________________________________________________________________
525 Long64_t AliHFEcollection::Merge(TCollection *list){
526
527   //
528   // Merge the collections
529   //
530   if(!list)
531     return 0;
532   if(list->IsEmpty())
533     return 1;
534   
535   TIter it(list);
536   TObject *o = NULL;
537   Int_t index = 0;
538   TList templist;       // Create temporary list containing all the lists to merge
539   while((o = it())){
540     AliHFEcollection *coll = dynamic_cast<AliHFEcollection *>(o);
541     if(!coll) continue; 
542     templist.Add(coll->fList);
543     index++;
544   }
545   fList->Merge(&templist);
546   return index + 1;
547 }
548 //____________________________________________________________________
549 void AliHFEcollection::Browse(TBrowser *b)
550 {
551
552   //
553   // Browse the content of the directory.
554   //
555
556    if (b) {
557       TObject *obj = 0;
558       TIter nextin(fList);
559
560       //Add objects that are only in memory
561       while ((obj = nextin())) {
562          b->Add(obj, obj->GetName());
563       }
564    }
565 }