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