2 // Class AliRsnListOutput
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)
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
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.
19 // author: A. Pulvirenti (email: alberto.pulvirenti@ct.infn.it)
22 #include <Riostream.h>
24 #include <TCollection.h>
27 #include "AliCFContainer.h"
29 #include "AliRsnValue.h"
30 #include "AliRsnValueDaughter.h"
31 #include "AliRsnValueEvent.h"
32 #include "AliRsnLoop.h"
34 #include "AliRsnListOutput.h"
37 ClassImp(AliRsnListOutput)
39 //________________________________________________________________________________________
40 AliRsnListOutput::AliRsnListOutput(const char *name, AliRsnListOutput::EOut type) :
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.
58 //________________________________________________________________________________________
59 AliRsnListOutput::AliRsnListOutput(const AliRsnListOutput ©) :
61 fSkipFailed(copy.fSkipFailed),
64 fValues(copy.fValues),
65 fNValues(copy.fNValues),
72 // Since the pointer objects must be initialized in a second step,
73 // they are never copied, and then they are initialized to zero.
77 //________________________________________________________________________________________
78 AliRsnListOutput &AliRsnListOutput::operator=(const AliRsnListOutput ©)
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.
87 TNamed::operator=(copy);
90 fSkipFailed = copy.fSkipFailed;
93 fValues = copy.fValues;
94 fNValues = copy.fNValues;
104 //__________________________________________________________________________________________________
105 AliRsnListOutput::~AliRsnListOutput()
109 // Deletes the output objects.
115 //__________________________________________________________________________________________________
116 void AliRsnListOutput::Reset()
119 // Clear all output objects. In general, only one will be initialized at a time.
125 //_____________________________________________________________________________
126 void AliRsnListOutput::AddValue(AliRsnValue *value)
129 // Adds a value computation object to the list.
132 fValues.AddLast(value);
136 //________________________________________________________________________________________
137 Bool_t AliRsnListOutput::Init(const char *prefix, TList *list)
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>
149 // all output objects are cleared
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();
156 AliError("Need at least 1 value");
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;
164 // resize the output array
165 fArray.Set(fNValues);
168 TString name(Form("%s_%s", prefix, GetName()));
169 AliRsnValue *val = 0x0;
170 for (i = 0; i < fNValues; i++) {
173 AliError(Form("Slot %d in value list is NULL", i));
177 name += val->GetName();
181 TObject *object = 0x0;
183 // initialize appropriate output object
184 // and, if successful, insert into the list
187 name.Append("_hist");
188 object = CreateHistogram(name.Data());
191 name.Append("_hsparse");
192 object = CreateHistogramSparse(name.Data());
196 object = CreateCFContainer(name.Data());
199 AliWarning("Wrong type output or initialization failure");
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()));
206 fIndex = fList->IndexOf(object);
213 //________________________________________________________________________________________
214 TH1 *AliRsnListOutput::CreateHistogram(const char *name)
217 // Initialize the 'default' TH1 output object.
218 // In case one of the expected axes is NULL, the initialization fails.
221 // we expect to have maximum 3 axes in this case
222 Int_t i, nbins[3] = {0, 0, 0};
224 for (i = 0; i < fNValues; i++) {
225 AliRsnValue *val = GetValue(i);
227 AliError(Form("Expected axis %d is NULL", i));
230 nbins[i] = GetValue(i)->GetArray().GetSize() - 1;
231 array[i] = GetValue(i)->GetArray();
236 // create histogram depending on the number of axes
239 hist = new TH1F(name, "", nbins[0], array[0].GetArray());
242 hist = new TH2F(name, "", nbins[0], array[0].GetArray(), nbins[1], array[1].GetArray());
245 hist = new TH3F(name, "", nbins[0], array[0].GetArray(), nbins[1], array[1].GetArray(), nbins[2], array[2].GetArray());
248 //AliError(Form("Wrong number of dimensions: %d", fNValues))
252 if (hist) hist->Sumw2();
256 //________________________________________________________________________________________
257 THnSparseF *AliRsnListOutput::CreateHistogramSparse(const char *name)
260 // Initialize the THnSparse output object.
261 // In case one of the expected axes is NULL, the initialization fails.
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();
275 THnSparseF *hist = new THnSparseF(name, "", fNValues, nbins);
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());
290 //________________________________________________________________________________________
291 AliCFContainer *AliRsnListOutput::CreateCFContainer(const char *name)
294 // Initialize the AliCFContainer output object.
295 // In case one of the expected axes is NULL, the initialization fails.
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();
309 AliCFContainer *cont = new AliCFContainer(name, "", fSteps, fNValues, nbins);
311 // set the bin limits for each axis
312 for (i = 0; i < fNValues; i++) {
313 cont->SetBinLimits(i, array[i].GetArray());
323 //________________________________________________________________________________________
324 Bool_t AliRsnListOutput::Fill(TObject *target, Int_t step)
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.
336 Bool_t globalOK = kTRUE;
337 AliRsnValue *val = 0x0;
338 for (i = 0; i < fNValues; i++) {
341 AliError("NULL value found");
344 globalOK = globalOK && val->Eval(target);
345 fArray[i] = (Double_t)val->GetComputedValue();
347 if (!globalOK && fSkipFailed) return kFALSE;
350 if (!fList || fIndex < 0) {
351 AliError("List not initialized");
354 TObject *obj = fList->At(fIndex);
356 AliError("Null pointer");
361 //AliInfo(Form("[%s] Object index, name, type = %d, %s (%s)", GetName(), fIndex, obj->GetName(), obj->ClassName()));
364 if (obj->IsA() == TH1F::Class()) {
365 TH1F *h = (TH1F *)obj;
368 } else if (obj->IsA() == TH2F::Class()) {
369 TH2F *h = (TH2F *)obj;
370 h->Fill(fArray[0], fArray[1]);
372 } else if (obj->IsA() == TH3F::Class()) {
373 TH3F *h = (TH3F *)obj;
374 h->Fill(fArray[0], fArray[1], fArray[2]);
376 } else if (obj->InheritsFrom(THnSparse::Class())) {
377 THnSparseF *h = (THnSparseF *)obj;
378 h->Fill(fArray.GetArray());
380 } else if (obj->InheritsFrom(AliCFContainer::Class())) {
381 AliCFContainer *c = (AliCFContainer *)obj;
382 c->Fill(fArray.GetArray(), step);
385 AliError(Form("Not handled class '%s'", obj->ClassName()));
390 //________________________________________________________________________________________
391 Bool_t AliRsnListOutput::Fill(AliRsnEvent *ev, AliRsnDaughter *d)
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.
400 if (!fList || fIndex < 0) {
401 AliError("List not initialized");
404 TObject *obj = fList->At(fIndex);
406 AliError("Null pointer");
411 TIter next(&fValues);
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);
424 values[i] = (Double_t)val->GetComputedValue();
425 if (!globalOK && fSkipFailed) return kFALSE;
430 if (obj->IsA() == TH1F::Class()) {
431 TH1F *h = (TH1F *)obj;
434 } else if (obj->IsA() == TH2F::Class()) {
435 TH2F *h = (TH2F *)obj;
436 h->Fill(values[0], values[1]);
438 } else if (obj->IsA() == TH3F::Class()) {
439 TH3F *h = (TH3F *)obj;
440 h->Fill(values[0], values[1], values[2]);
442 } else if (obj->InheritsFrom(THnSparse::Class())) {
443 THnSparseF *h = (THnSparseF *)obj;
447 AliError(Form("Not handled class '%s'", obj->ClassName()));