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