]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliDielectronHFhelper.cxx
-iarsene: FMD includes commented out
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliDielectronHFhelper.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 //       Dielectron Histogram framework helper                           //
18 //                                                                       //
19 /*
20
21 A helper class to extract objects(histograms and/or profiles) from
22 a AliDielctronHF array of objects.
23
24
25 How to use it:
26
27   AliDielectronHFhelper *hf = new AliDielectronHFhelper("path/to/the/output/file.root", "ConfigName");
28   // print the structure
29   hf->Print();
30
31   //apply some cuts and print them
32   hf->SetRangeUser("cut1name",cutmin1,cutmax1);
33   hf->SetRangeUser(AliDielectronVarManager::kPt,ptmin,ptmax);
34   hf->PrintCuts();
35
36   // collect 1-,2- or 3-dim histograms or profiles with error option (default:"")
37   TObjArray *arrHists = hf->CollectHistos(AliDielectronVarManager::kM);
38   TObjArray *arrProfs = hf->CollectProfiles("",AliDielectronVarManager::kM,AliDielectronVarManager::kPt);
39
40   // then you are left with an array of histograms for all pair types or MC signals
41
42 */
43 //                                                                       //
44 ///////////////////////////////////////////////////////////////////////////
45
46 #include <TObjArray.h>
47 #include <TKey.h>
48 #include <TList.h>
49 #include <TClass.h>
50 #include <TObject.h>
51 #include <TFile.h>
52 #include <TString.h>
53 #include <TObjString.h>
54 #include <TVectorD.h>
55 #include <TMath.h>
56 #include <TH1.h>
57
58 #include <AliLog.h>
59
60 #include "AliDielectron.h"
61 #include "AliDielectronHFhelper.h"
62 #include "AliDielectronHF.h"
63
64 ClassImp(AliDielectronHFhelper)
65
66 //const char* AliDielectronHFhelper::fCutVars[AliDielectronHFhelper::kMaxCuts] = {"","","","","","","","","",""};
67
68 //________________________________________________________________
69 AliDielectronHFhelper::AliDielectronHFhelper(const char* filename, const char* container) :
70   TNamed(),
71   fMainArr(0x0),
72   fCutVars(0x0),
73   fCutLowLimits(0),
74   fCutUpLimits(0)
75 {
76   //
77   // get HF container(s) from file 'filename'
78   //
79   SetHFArray(filename, container);
80
81 }
82
83 //________________________________________________________________
84 AliDielectronHFhelper::~AliDielectronHFhelper()
85 {
86   //
87   // dtor
88   //
89   if(fMainArr) delete fMainArr;
90   if(fCutVars)     delete fCutVars;
91
92 }
93
94 //________________________________________________________________
95 void AliDielectronHFhelper::SetHFArray(const char* filename, const char* container)
96 {
97   //
98   // get HF container from file
99   //
100
101   TFile *f = TFile::Open(filename);
102
103   TList *l=f->GetListOfKeys();
104   TIter nextKey(l);
105   TKey *k=0x0;
106   while ( (k=static_cast<TKey*>(nextKey())) ){
107
108     TObject *o=k->ReadObj();
109     if (o->IsA()==TList::Class()){
110
111       TList *tlist=(TList*)o;
112
113       TIter next(tlist);
114       TObject *obj=0x0;
115       while ((obj = next())) {
116         TString objname(obj->GetName());
117
118         if( objname.Contains(Form("%s_HF",container)) && obj->IsA()==TObjArray::Class()) {
119           fMainArr = new TObjArray( *(dynamic_cast<TObjArray*>(obj)) );
120           //fMainArr->Print();
121           return;
122         }
123       }
124     }
125   }
126
127 }
128 //________________________________________________________________
129 void AliDielectronHFhelper::SetRangeUser(const char *varname, Double_t min, Double_t max, Bool_t leg)
130 {
131   //
132   // Set range from variable name
133   //
134
135   Int_t size=fCutLowLimits.GetNrows();
136
137   // check if cut is already set
138   for(Int_t icut=0; icut<size; icut++) {
139     TString cutName = fCutVars->At(icut)->GetName();
140     if(!cutName.CompareTo(Form("%s%s",(leg?"Leg":""),varname))) {
141       UnsetRangeUser(varname,leg);
142       SetRangeUser(varname, min, max, leg);
143       return;
144     }
145   }
146
147   if(size>=kMaxCuts) return;
148
149   // arrays
150   if(!fCutVars) {
151     fCutVars = new TObjArray();
152     fCutVars->SetOwner();
153   }
154   fCutLowLimits.ResizeTo(size+1);
155   fCutUpLimits.ResizeTo(size+1);
156
157   // fill
158   TObjString *str = new TObjString(Form("%s%s",(leg?"Leg":""),varname));
159   fCutVars->Add(str);
160   fCutLowLimits(size) = min;
161   fCutUpLimits(size)  = max;
162   AliWarning(Form(" %s [%.2f,%.2f]",fCutVars->At(size)->GetName(),fCutLowLimits(size),fCutUpLimits(size)));
163 }
164
165 //________________________________________________________________
166 void AliDielectronHFhelper::SetRangeUser(AliDielectronVarManager::ValueTypes type, Double_t min, Double_t max, Bool_t leg)
167 {
168   //
169   // Set range from AliDielectronVarManager
170   //
171   SetRangeUser(AliDielectronVarManager::GetValueName(type), min, max, leg);
172 }
173
174 //________________________________________________________________
175 void AliDielectronHFhelper::UnsetRangeUser(const char *varname, Bool_t leg)
176 {
177   //
178   // unset range from variable name
179   //
180
181   Int_t size=fCutLowLimits.GetNrows();
182   TVectorD newlow;
183   TVectorD newup;
184
185   // find cut and build new vectors w/o it
186   Int_t ientries = 0;
187   for(Int_t icut=0; icut<size; icut++) {
188
189     TString cutName = fCutVars->At(icut)->GetName();
190     if(!cutName.CompareTo(Form("%s%s",(leg?"Leg":""),varname))) {
191       fCutVars->AddAt(0x0,icut);
192       continue;
193     }
194
195     // fill new vectors
196     newlow.ResizeTo(ientries+1);
197     newup.ResizeTo(ientries+1);
198     newlow(ientries) = fCutLowLimits(icut);
199     newup(ientries)  = fCutUpLimits(icut);
200
201     ientries++;
202   }
203
204   // adapt new arrays/vectors
205   fCutVars->Compress();
206
207   fCutLowLimits.ResizeTo(ientries);
208   fCutUpLimits.ResizeTo(ientries);
209   for(Int_t icut=0; icut<ientries; icut++) {
210     fCutLowLimits(icut) = newlow(icut);
211     fCutUpLimits(icut)  = newup(icut);
212   }
213   // PrintCuts();
214
215 }
216
217 //________________________________________________________________
218 void AliDielectronHFhelper::UnsetRangeUser(AliDielectronVarManager::ValueTypes type, Bool_t leg)
219 {
220   //
221   // Unset range from AliDielectronVarManager
222   //
223   UnsetRangeUser(AliDielectronVarManager::GetValueName(type), leg);
224 }
225
226 //________________________________________________________________
227 TObjArray* AliDielectronHFhelper::CollectProfiles(TString option,
228                                                   AliDielectronVarManager::ValueTypes varx,
229                                                   AliDielectronVarManager::ValueTypes vary,
230                                                   AliDielectronVarManager::ValueTypes varz,
231                                                   AliDielectronVarManager::ValueTypes vart)
232 {
233   //
234   // collect 1-3 dimensional TProfiles for all kind of pair types or sources
235   //
236
237   // reconstruct the histogram/profile name resp. key
238   Int_t dim    = 0;
239   if(varx < AliDielectronVarManager::kNMaxValues) dim++;
240   if(vary < AliDielectronVarManager::kNMaxValues) dim++;
241   if(varz < AliDielectronVarManager::kNMaxValues) dim++;
242   Bool_t bPairClass=0;
243   if( varx < AliDielectronVarManager::kPairMax ||
244       vary < AliDielectronVarManager::kPairMax ||
245       varz < AliDielectronVarManager::kPairMax ||
246       vart < AliDielectronVarManager::kPairMax  ) bPairClass=kTRUE;
247   Bool_t bprf = !(option.Contains("hist",TString::kIgnoreCase));
248   Bool_t bStdOpt=kTRUE;
249   if(option.Contains("S",TString::kIgnoreCase) || option.Contains("rms",TString::kIgnoreCase)) bStdOpt=kFALSE;
250
251   TString key = "";
252   if(bprf) dim--;
253   switch(dim) {
254   case 3:
255     key+=Form("%s_",AliDielectronVarManager::GetValueName(varx));
256     key+=Form("%s_",AliDielectronVarManager::GetValueName(vary));
257     key+=Form("%s",AliDielectronVarManager::GetValueName(varz));
258     if(bprf) key+=Form("-%s%s",AliDielectronVarManager::GetValueName(vart),(bStdOpt ? "avg" : "rms"));
259     break;
260   case 2:
261     key+=Form("%s_",AliDielectronVarManager::GetValueName(varx));
262     key+=Form("%s",AliDielectronVarManager::GetValueName(vary));
263     if(bprf) key+=Form("-%s%s",AliDielectronVarManager::GetValueName(varz),(bStdOpt ? "avg" : "rms"));
264     break;
265   case 1:
266     key+=Form("%s",AliDielectronVarManager::GetValueName(varx));
267     if(bprf) key+=Form("-%s%s",AliDielectronVarManager::GetValueName(vary),(bStdOpt ? "avg" : "rms"));
268     break;
269   }
270   // to differentiate btw. leg and pair histos
271   if(bPairClass) key.Prepend("p");
272   // prepend HF
273   key.Prepend("HF_");
274
275   //TODO:  printf("--------> KEY: %s \n",key.Data());
276   // get the requested object from the arrays
277   TObjArray *collection = new TObjArray(fMainArr->GetEntriesFast());
278
279   TObjArray *cloneArr = (TObjArray*) fMainArr->Clone("tmpArr");
280   if(!cloneArr) return 0x0;
281   cloneArr->SetOwner(kTRUE);
282
283   // loop over all pair types or sources
284   for(Int_t i=0; i<cloneArr->GetEntriesFast(); i++) {
285     if(!cloneArr->At(i))                             continue;
286     if(!((TObjArray*)cloneArr->At(i))->GetEntries()) continue;
287     TString stepName = cloneArr->At(i)->GetName();
288
289     // find histogram of interest from array of objects
290     TObjArray *arr = (TObjArray*) GetObject(cloneArr->At(i)->GetName(),cloneArr);
291     if(arr) {
292       collection->AddAt(arr->FindObject(key.Data()), i);
293       if(!collection->At(i)) AliError(Form("Object %s does not exist",key.Data()));
294       // modify the key name
295       stepName.ReplaceAll("(","");
296       stepName.ReplaceAll(")","");
297       stepName.ReplaceAll(": ","");
298       stepName.ReplaceAll(" ","_");
299       stepName.ReplaceAll("Signal","");
300       ((TH1*)collection->At(i))->SetName(Form("%s_%s",key.Data(),stepName.Data()));
301     }
302
303   }
304
305   // clean up the clone
306   if(cloneArr) {
307     delete cloneArr;
308     cloneArr=0;
309   }
310
311   return collection;
312 }
313
314 //________________________________________________________________
315 TObject* AliDielectronHFhelper::GetObject(const char *step, TObjArray *histArr)
316 {
317   //
318   // main function to recieve a pair type or MC signal array
319   //
320
321   AliDebug(1,Form(" Step %s selected",step));
322   // TODO: check memory
323
324   TObjArray *stepArr = 0x0; // this is the requested step
325   TObject *hist      = 0x0;
326   if(!histArr) {
327     stepArr = (TObjArray*) fMainArr->FindObject(step)->Clone("tmpArr");
328   }
329   else {
330     stepArr = (TObjArray*) histArr->FindObject(step);
331   }
332
333   if(stepArr) hist   = FindObjects(stepArr);
334   return hist;
335
336 }
337
338 //________________________________________________________________
339 TObject* AliDielectronHFhelper::FindObjects(TObjArray *stepArr)
340 {
341   //
342   // apply cuts and exclude objects from the array
343   //
344
345   // debug
346   // TString title    = stepArr->At(0)->GetTitle();
347   // TObjArray* vars  = title.Tokenize(":");
348   // AliDebug(1,Form(" number of cuts/vars: %d/%d",fCutLowLimits.GetNrows(),vars->GetEntriesFast()));
349
350   // check for missing cuts
351   CheckCuts(stepArr);
352
353   // loop over all cuts
354   for(Int_t icut=0; icut<fCutLowLimits.GetNrows(); icut++) {
355
356     Bool_t bFndBin      = kFALSE; // bin with exact limits found
357     const char *cutvar  = fCutVars->At(icut)->GetName(); // cut variable
358     Double_t min        = fCutLowLimits(icut);           // lower limit
359     Double_t max        = fCutUpLimits(icut);            // upper limit
360     AliDebug(5,Form(" Cut %d: %s [%.2f,%.2f]",icut,cutvar,min,max));
361
362     // loop over the full grid of given step
363     for(Int_t i=0; i<stepArr->GetEntriesFast(); i++) {
364
365       // continue if already empty
366       if(!stepArr->At(i)) continue;
367
368       // collect bins from the name
369       TString title    = stepArr->At(i)->GetName();
370       if(title.IsNull()) continue;
371       AliDebug(5,Form(" %03d object name: %s",i,title.Data()));
372
373       // loop over all variables
374       TObjArray *vars  = title.Tokenize(":");
375       for(Int_t ivar=0; ivar<vars->GetEntriesFast(); ivar++) {
376         TString binvar = vars->At(ivar)->GetName();
377         AliDebug(10,Form(" --> %d check bin var %s",ivar,binvar.Data()));
378
379         // check cuts set by the user, and compare to the bin limits
380         if(binvar.Contains(cutvar)) {
381           // bin limits
382           TObjArray *limits = binvar.Tokenize("#");
383           Double_t binmin = atof(limits->At(1)->GetName()); // lower bin limit
384           Double_t binmax = atof(limits->At(2)->GetName()); // upper bin limit
385           AliDebug(10,Form(" bin %s var %s [%.2f,%.2f]",binvar.Data(),limits->At(0)->GetName(),binmin,binmax));
386           if(limits) delete limits;
387
388           // cut and remove objects from the array
389           if(binmin < min || binmax < min || binmin > max || binmax > max ) {
390             AliDebug(10,Form(" removed, out of range! lower bin,cut: %.2f,%.2f or upper bin,cut: %.2f>%.2f",binmin,min,binmax,max));
391             stepArr->AddAt(0x0,i);
392           }
393           if(bFndBin && !(binmin == min && binmax == max)) {
394             stepArr->AddAt(0x0,i);
395             AliDebug(10,Form(" removed, within range! lower bin,cut: %.2f,%.2f or upper bin,cut: %.2f,%.2f",binmin,min,binmax,max));
396           }
397
398           // did we found a bin with exact the cut ranges
399           // this can happen only once per variable
400           if(binmin==min && binmax==max) bFndBin=kTRUE;
401
402         }
403
404       }
405       // clean up
406       if(vars)   delete vars;
407
408     }
409
410   }
411
412   // compress the array by removing all empty entries
413   stepArr->Compress();
414   AliDebug(1,Form(" Compression: %d objects left",stepArr->GetEntriesFast()));
415
416   // merge left objects
417   TObject* hist = Merge(stepArr);
418   //  if(hist) AliDebug(1,Form(" Merging: %e  entries",hist->GetEntries()));
419   return hist;
420 }
421
422 //________________________________________________________________
423 TObject* AliDielectronHFhelper::Merge(TObjArray *arr)
424 {
425   //
426   // merge left objects to a single one
427   //
428
429   if(arr->GetEntriesFast()<1) { AliError(" No more objects left!"); return 0x0; }
430
431   TObject *final=arr->At(0)->Clone();
432   if(!final) return 0x0;
433
434   TList listH;
435   TString listHargs;
436   listHargs.Form("((TCollection*)0x%lx)", (ULong_t)&listH);
437   Int_t error = 0;
438
439   //  final->Reset("CE");
440   //  final->SetTitle(""); //TODO: change in future
441   for(Int_t i=1; i<arr->GetEntriesFast(); i++) {
442     listH.Add(arr->At(i));
443     //   printf("%d: ent %.0f \n",i,((TH1*)((TObjArray*)arr->At(i))->At(0))->GetEntries());
444     //    final->Add((TH1*)arr->At(i));
445   }
446   //  arr->Clear();
447
448   final->Execute("Merge", listHargs.Data(), &error);
449   return final;
450 }
451
452 //________________________________________________________________
453 void AliDielectronHFhelper::CheckCuts(TObjArray *arr)
454 {
455   //
456   // Compare binning and cut variables. Add necessary cuts (full range, no exclusion)
457   //
458
459   // build array with bin variable, minimum and maximum bin values
460   TString titleFIRST       = arr->First()->GetName();
461   TString titleLAST        = arr->Last()->GetName();
462   TObjArray* binvarsF  = titleFIRST.Tokenize(":#");
463   TObjArray* binvarsL  = titleLAST.Tokenize(":#");
464   Double_t binmin[kMaxCuts]= {0.0};
465   Double_t binmax[kMaxCuts]= {0.0};
466   for(Int_t ivar=0; ivar<binvarsF->GetEntriesFast(); ivar++) {
467
468     TString elementF=binvarsF->At(ivar)->GetName();
469     TString elementL=binvarsL->At(ivar)->GetName();
470     AliDebug(1,Form(" binvar %d: %s,%s",ivar,elementF.Data(),elementL.Data()));
471
472     switch(ivar%3) {
473     case 0: continue; break;
474     case 1: binmin[(int)ivar/3]=atof(elementF.Data()); break;
475     case 2: binmax[(int)ivar/3]=atof(elementL.Data()); break;
476     }
477
478     binvarsF->AddAt(0x0,ivar);
479   }
480   binvarsF->Compress();
481
482   // loop over all vars and cuts, check for missing stuff
483   for(Int_t ivar=0; ivar<binvarsF->GetEntriesFast(); ivar++) {
484
485     TString binvar=binvarsF->At(ivar)->GetName();
486     Bool_t selected=kFALSE;
487
488     AliDebug(1,Form(" check cuts %d %s [%.2f,%.2f]",ivar,binvar.Data(),binmin[ivar],binmax[ivar]));
489     // loop over all cuts and check for missing stuff
490     for(Int_t icut=0; icut<fCutLowLimits.GetNrows(); icut++) {
491       if(binvar.Contains(fCutVars->At(icut)->GetName())) { selected=kTRUE; break; }
492       //      else break;
493     }
494
495     // add missing cut with max limits
496     if(!selected) {
497       AliWarning(Form(" Bin variable %s not covered. Add cut!",binvar.Data()));
498       Bool_t leg = binvar.BeginsWith("Leg");
499       if(leg) binvar.Remove(0,3);
500       SetRangeUser(binvar.Data(),binmin[ivar],binmax[ivar], leg);
501     }
502
503   }
504
505   // clean up
506   if(binvarsF) delete binvarsF;
507   if(binvarsL) delete binvarsL;
508 }
509
510 //________________________________________________________________
511 void AliDielectronHFhelper::Print(const Option_t* /*option*/) const
512 {
513
514   //
515   // Print out object contents
516   //
517   AliInfo(Form(" Container:               %s",fMainArr->GetName()));
518
519   // pairtypes, steps and sources
520   Int_t stepLast=0;
521   AliInfo(Form(" Number of filled steps:  %d",fMainArr->GetEntries()));
522   for(Int_t istep=0; istep<fMainArr->GetEntriesFast(); istep++) {
523     if(!fMainArr->At(istep))                             continue;
524     if(!((TObjArray*)fMainArr->At(istep))->GetEntries()) continue;
525     AliInfo(Form(" step %d: %s",istep,fMainArr->At(istep)->GetName()));
526     stepLast=istep;
527   }
528
529   AliInfo(Form(" Number of objects:    %d",
530                ((TObjArray*) ((TObjArray*)fMainArr->At(stepLast)) ->First())->GetEntriesFast()));
531
532   TString title       = ((TObjArray*)fMainArr->At(stepLast))->First()->GetName();
533   TObjArray* binvars  = title.Tokenize(":");
534   AliInfo(Form(" Number of variables:     %d",binvars->GetEntriesFast()));
535   delete binvars;
536
537   TObjArray* binvars2  = title.Tokenize(":#");
538   for(Int_t ivar=0; ivar<binvars2->GetEntriesFast(); ivar++) {
539     if(ivar%3) continue;
540     AliInfo(Form(" variable %.0f: %s",((Double_t)ivar)/3+1,binvars2->At(ivar)->GetName()));
541   }
542   delete binvars2;
543
544 }
545
546 //________________________________________________________________
547 void AliDielectronHFhelper::PrintCuts()
548 {
549
550   //
551   // Print cuts
552   //
553
554   // loop over all cuts
555   AliInfo(" Selected cuts:");
556   for(Int_t icut=0; icut<fCutLowLimits.GetNrows(); icut++)
557     AliInfo(Form(" %d: %s [%.2f,%.2f]",icut,fCutVars->At(icut)->GetName(),fCutLowLimits(icut),fCutUpLimits(icut)));
558
559 }
560