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