]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnListOutput.cxx
Fix for K0 pdg code (F.Blanco) + option for filling histos (M.Vala)
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnListOutput.cxx
1 //
2 // Class AliRsnListOutput
3 //
4 // This class defines a base classe to implement a Output
5 // which uses the internal RSN package event format (AliRsnEvent).
6 // It contains some default flags which turn out to be useful:
7 //  - a flag to select only the "true" pairs (tracks from same resonance)
8 //  - a flag to know if the computation is done over two events (mixing)
9 //
10 // Any kind of analysis object should be implemented as inheriting from this
11 // because the AliRsnAnalyzer which executes the analysis will accept a collection
12 // of such objects, in order to have a unique format of processing method
13 //
14 // The user who implements a kind of computation type should inherit from
15 // this class and override the virtual Outputs defined in it, which
16 // initialize the final output histogram and define how to process data.
17 //
18 //
19 // author: A. Pulvirenti             (email: alberto.pulvirenti@ct.infn.it)
20 //
21
22 #include <Riostream.h>
23 #include <TList.h>
24 #include <TCollection.h>
25
26 #include "AliLog.h"
27 #include "AliCFContainer.h"
28
29 #include "AliRsnValue.h"
30 #include "AliRsnValueDaughter.h"
31 #include "AliRsnValueEvent.h"
32 #include "AliRsnLoop.h"
33
34 #include "AliRsnListOutput.h"
35
36 ClassImp(AliRsnListOutput)
37
38 //________________________________________________________________________________________
39 AliRsnListOutput::AliRsnListOutput(const char *name, AliRsnListOutput::EOut type) :
40    TNamed(name, ""),
41    fSkipFailed(kTRUE),
42    fType(type),
43    fSteps(0),
44    fValues(0),
45    fNValues(0),
46    fList(0x0),
47    fIndex(-1),
48    fArray(0)
49 {
50 //
51 // Constructor.
52 // Requires a name for this object (which will be used to name the output object)
53 // and the definition of the output type from the built-in enumeration.
54 //
55 }
56
57 //________________________________________________________________________________________
58 AliRsnListOutput::AliRsnListOutput(const AliRsnListOutput &copy) :
59    TNamed(copy),
60    fSkipFailed(copy.fSkipFailed),
61    fType(copy.fType),
62    fSteps(copy.fSteps),
63    fValues(copy.fValues),
64    fNValues(copy.fNValues),
65    fList(copy.fList),
66    fIndex(copy.fIndex),
67    fArray(0)
68 {
69 //
70 // Copy constructor.
71 // Since the pointer objects must be initialized in a second step,
72 // they are never copied, and then they are initialized to zero.
73 //
74 }
75
76 //________________________________________________________________________________________
77 AliRsnListOutput &AliRsnListOutput::operator=(const AliRsnListOutput &copy)
78 {
79 //
80 // Assignment operator.
81 // Same consideration as the copiy constructor. In this case, there is
82 // the possibility to have the output objects alreasy initialized, but
83 // they are anyway reset.
84 //
85
86    TNamed::operator=(copy);
87    if (this == &copy)
88       return *this;
89    fSkipFailed = copy.fSkipFailed;
90    fType = copy.fType;
91    fSteps = copy.fSteps;
92    fValues = copy.fValues;
93    fNValues = copy.fNValues;
94    fList = copy.fList;
95    fIndex = copy.fIndex;
96    fArray = copy.fArray;
97
98    Reset();
99
100    return (*this);
101 }
102
103 //__________________________________________________________________________________________________
104 AliRsnListOutput::~AliRsnListOutput()
105 {
106 //
107 // Destructor.
108 // Deletes the output objects.
109 //
110
111    Reset();
112 }
113
114 //__________________________________________________________________________________________________
115 void AliRsnListOutput::Reset()
116 {
117 //
118 // Clear all output objects. In general, only one will be initialized at a time.
119 //
120
121    fList = 0x0;
122 }
123
124 //_____________________________________________________________________________
125 void AliRsnListOutput::AddValue(AliRsnValue *value)
126 {
127 //
128 // Adds a value computation object to the list.
129 //
130
131    fValues.AddLast(value);
132 }
133
134
135 //________________________________________________________________________________________
136 Bool_t AliRsnListOutput::Init(const char *prefix, TList *list)
137 {
138 //
139 // Initializes the output for this object.
140 // What kind of output depends on the 'fType' data member,
141 // and in case it is a CF container, also on the 'fSteps'.
142 // The object is named with the following criterion:
143 // <prefix>_<name>_<type>_<varList>
144 //
145
146    Int_t i;
147
148    // all output objects are cleared
149    Reset();
150
151    // count values and set dimension of arrays
152    // do also some checks for a good match between type and output
153    fNValues = fValues.GetEntries();
154    if (fNValues < 1) {
155       AliError("Need at least 1 value");
156       return kFALSE;
157    }
158    if (fType == kHistoDefault && fNValues > 3) {
159       AliInfo(Form("NValues = %d > 3 --> cannot use a normal histogram, need to use a sparse", fNValues));
160       fType = kHistoSparse;
161    }
162
163    // resize the output array
164    fArray.Set(fNValues);
165
166
167    Bool_t isPair=kFALSE;
168
169    // create the name
170    TString name(GetName());
171    if (!name.CompareTo("pair")) isPair = kTRUE;
172    if (isPair) name = "";
173    else name.Prepend(".");
174    name.Prepend(prefix);
175
176 //    TString name(Form("%s.%s", prefix, GetName()));
177    AliRsnValue *val = 0x0;
178    for (i = 0; i < fNValues; i++) {
179       val = GetValue(i);
180       if (!val) {
181          AliError(Form("Slot %d in value list is NULL", i));
182          return kFALSE;
183       }
184       if (!isPair) {
185          name += '_';
186          name += val->GetName();
187       }
188    }
189
190    // allowed objects
191    TObject *object = 0x0;
192
193    // initialize appropriate output object
194    // and, if successful, insert into the list
195    switch (fType) {
196       case kHistoDefault:
197 //          name.Append("_hist");
198          object = CreateHistogram(name.Data());
199          break;
200       case kHistoSparse:
201 //          name.Append("_hsparse");
202          object = CreateHistogramSparse(name.Data());
203          break;
204       case kCFContainer:
205 //          name.Append("_cf");
206          object = CreateCFContainer(name.Data());
207          break;
208       default:
209          AliWarning("Wrong type output or initialization failure");
210    }
211
212    if (object) {
213       //AliInfo(Form("[%s]: initializing output '%s' (obj name = '%s') with %d values and format %d [%s]", GetName(), name.Data(), object->GetName(), fNValues, fType, object->ClassName()));
214       fList = list;
215       fList->Add(object);
216       fIndex = fList->IndexOf(object);
217       return kTRUE;
218    }
219
220    return kFALSE;
221 }
222
223 //________________________________________________________________________________________
224 TH1 *AliRsnListOutput::CreateHistogram(const char *name)
225 {
226 //
227 // Initialize the 'default' TH1 output object.
228 // In case one of the expected axes is NULL, the initialization fails.
229 //
230
231    // we expect to have maximum 3 axes in this case
232    Int_t i, nbins[3] = {0, 0, 0};
233    TArrayD  array[3];
234    for (i = 0; i < fNValues; i++) {
235       AliRsnValue *val = GetValue(i);
236       if (!val) {
237          AliError(Form("Expected axis %d is NULL", i));
238          return 0x0;
239       }
240       nbins[i] = GetValue(i)->GetArray().GetSize() - 1;
241       array[i] = GetValue(i)->GetArray();
242    }
243
244    TH1 *hist = 0x0;
245
246    // create histogram depending on the number of axes
247    switch (fNValues) {
248       case 1:
249          hist = new TH1F(name, "", nbins[0], array[0].GetArray());
250          break;
251       case 2:
252          hist = new TH2F(name, "", nbins[0], array[0].GetArray(), nbins[1], array[1].GetArray());
253          break;
254       case 3:
255          hist = new TH3F(name, "", nbins[0], array[0].GetArray(), nbins[1], array[1].GetArray(), nbins[2], array[2].GetArray());
256          break;
257       default:
258          //AliError(Form("Wrong number of dimensions: %d", fNValues))
259          return 0x0;
260    }
261
262    if (hist) hist->Sumw2();
263    return hist;
264 }
265
266 //________________________________________________________________________________________
267 THnSparseF *AliRsnListOutput::CreateHistogramSparse(const char *name)
268 {
269 //
270 // Initialize the THnSparse output object.
271 // In case one of the expected axes is NULL, the initialization fails.
272 //
273
274    // retrieve binnings and sizes of all axes
275    // since the check for null values is done in Init(),
276    // we assume that here they must all be well defined
277    Int_t i, *nbins = new Int_t[fNValues];
278    TArrayD  *array = new TArrayD[fNValues];
279    for (i = 0; i < fNValues; i++) {
280       nbins[i] = GetValue(i)->GetArray().GetSize() - 1;
281       array[i] = GetValue(i)->GetArray();
282    }
283
284    // create histogram
285    THnSparseF *hist = new THnSparseF(name, "", fNValues, nbins);
286    hist->Sumw2();
287
288    // update the various axes using the definitions given in the array of axes here
289    AliRsnValue *val = 0x0;
290    for (i = 0; i < fNValues; i++) {
291       val = GetValue(i);
292       if (val) hist->GetAxis(i)->SetName(val->GetName());
293       hist->GetAxis(i)->Set(nbins[i], array[i].GetArray());
294    }
295
296    // clear heap
297    delete [] nbins;
298    delete [] array;
299
300    return hist;
301 }
302
303 //________________________________________________________________________________________
304 AliCFContainer *AliRsnListOutput::CreateCFContainer(const char *name)
305 {
306 //
307 // Initialize the AliCFContainer output object.
308 // In case one of the expected axes is NULL, the initialization fails.
309 //
310
311    // retrieve binnings and sizes of all axes
312    // since the check for null values is done in Init(),
313    // we assume that here they must all be well defined
314    Int_t i, *nbins = new Int_t[fNValues];
315    TArrayD  *array = new TArrayD[fNValues];
316    for (i = 0; i < fNValues; i++) {
317       nbins[i] = GetValue(i)->GetArray().GetSize() - 1;
318       array[i] = GetValue(i)->GetArray();
319    }
320
321    // create object
322    AliCFContainer *cont = new AliCFContainer(name, "", fSteps, fNValues, nbins);
323
324    // set the bin limits for each axis
325    for (i = 0; i < fNValues; i++) {
326       cont->SetBinLimits(i, array[i].GetArray());
327    }
328
329    // clear heap
330    delete [] nbins;
331    delete [] array;
332
333    return cont;
334 }
335
336 //________________________________________________________________________________________
337 Bool_t AliRsnListOutput::Fill(TObject *target, Int_t step)
338 {
339 //
340 // Uses the passed argument to compute all values.
341 // If all computations were successful, fill the output
342 // Second argument (step) is needed only in case of CF containers.
343 // Return value is the AND of all computation successes.
344 //
345
346    Int_t  i;
347
348    // do computations
349    Bool_t globalOK = kTRUE;
350    AliRsnValue *val = 0x0;
351    for (i = 0; i < fNValues; i++) {
352       val = GetValue(i);
353       if (!val) {
354          AliError("NULL value found");
355          return kFALSE;
356       }
357       globalOK = globalOK && val->Eval(target);
358       fArray[i] = (Double_t)val->GetComputedValue();
359    }
360    if (!globalOK && fSkipFailed) return kFALSE;
361
362    // retrieve object
363    if (!fList || fIndex < 0) {
364       AliError("List not initialized");
365       return kFALSE;
366    }
367    TObject *obj = fList->At(fIndex);
368    if (!obj) {
369       AliError("Null pointer");
370       return kFALSE;
371    }
372
373    // check
374    //AliInfo(Form("[%s] Object index, name, type = %d, %s (%s)", GetName(), fIndex, obj->GetName(), obj->ClassName()));
375
376    // fill output
377    if (obj->IsA() == TH1F::Class()) {
378       TH1F *h = (TH1F *)obj;
379       h->Fill(fArray[0]);
380       return kTRUE;
381    } else if (obj->IsA() == TH2F::Class()) {
382       TH2F *h = (TH2F *)obj;
383       h->Fill(fArray[0], fArray[1]);
384       return kTRUE;
385    } else if (obj->IsA() == TH3F::Class()) {
386       TH3F *h = (TH3F *)obj;
387       h->Fill(fArray[0], fArray[1], fArray[2]);
388       return kTRUE;
389    } else if (obj->InheritsFrom(THnSparse::Class())) {
390       THnSparseF *h = (THnSparseF *)obj;
391       for (Int_t iAxis = 0; iAxis<h->GetNdimensions(); iAxis++) {
392          if (fArray.At(iAxis)>h->GetAxis(iAxis)->GetXmax() || fArray.At(iAxis)<h->GetAxis(iAxis)->GetXmin()) return kFALSE;
393       }
394       h->Fill(fArray.GetArray());
395       return kTRUE;
396    } else if (obj->InheritsFrom(AliCFContainer::Class())) {
397       AliCFContainer *c = (AliCFContainer *)obj;
398       c->Fill(fArray.GetArray(), step);
399       return kTRUE;
400    } else {
401       AliError(Form("Not handled class '%s'", obj->ClassName()));
402       return kFALSE;
403    }
404 }
405
406 //________________________________________________________________________________________
407 Bool_t AliRsnListOutput::Fill(AliRsnEvent *ev, AliRsnDaughter *d)
408 {
409 //
410 // Uses the passed argument to compute all values.
411 // If all computations were successful, fill the output
412 // Return value is the AND of all computation successes.
413 //
414
415    // retrieve object
416    if (!fList || fIndex < 0) {
417       AliError("List not initialized");
418       return kFALSE;
419    }
420    TObject *obj = fList->At(fIndex);
421    if (!obj) {
422       AliError("Null pointer");
423       return kFALSE;
424    }
425
426
427    TIter next(&fValues);
428    AliRsnValue *val;
429    Int_t i=0;
430    Bool_t globalOK = kTRUE;
431    Double_t values[fValues.GetEntries()];
432    while ((val = (AliRsnValue *)next())) {
433       if (val->InheritsFrom(AliRsnValueDaughter::Class())) {
434          if (!d && fSkipFailed) return kFALSE;
435          globalOK = globalOK && val->Eval(d);
436       } else if (val->InheritsFrom(AliRsnValueEvent::Class())) {
437          if (!ev && fSkipFailed) return kFALSE;
438          globalOK = globalOK && val->Eval(ev);
439       }
440       values[i] = (Double_t)val->GetComputedValue();
441       if (!globalOK && fSkipFailed) return kFALSE;
442       i++;
443    }
444
445    // fill output
446    if (obj->IsA() == TH1F::Class()) {
447       TH1F *h = (TH1F *)obj;
448       h->Fill(values[0]);
449       return kTRUE;
450    } else if (obj->IsA() == TH2F::Class()) {
451       TH2F *h = (TH2F *)obj;
452       h->Fill(values[0], values[1]);
453       return kTRUE;
454    } else if (obj->IsA() == TH3F::Class()) {
455       TH3F *h = (TH3F *)obj;
456       h->Fill(values[0], values[1], values[2]);
457       return kTRUE;
458    } else if (obj->InheritsFrom(THnSparse::Class())) {
459       THnSparseF *h = (THnSparseF *)obj;
460       for (Int_t iAxis = 0; iAxis<h->GetNdimensions(); iAxis++) {
461          if (values[iAxis]>h->GetAxis(iAxis)->GetXmax() || values[iAxis]<h->GetAxis(iAxis)->GetXmin()) return kFALSE;
462       }
463       h->Fill(values);
464       return kTRUE;
465    } else {
466       AliError(Form("Not handled class '%s'", obj->ClassName()));
467       return kFALSE;
468    }
469 }