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