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"
36 ClassImp(AliRsnListOutput)
38 //________________________________________________________________________________________
39 AliRsnListOutput::AliRsnListOutput(const char *name, AliRsnListOutput::EOut type) :
48 fCheckHistRange(kTRUE),
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),
68 fCheckHistRange(copy.fCheckHistRange),
73 // Since the pointer objects must be initialized in a second step,
74 // they are never copied, and then they are initialized to zero.
78 //________________________________________________________________________________________
79 AliRsnListOutput &AliRsnListOutput::operator=(const AliRsnListOutput ©)
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.
88 TNamed::operator=(copy);
91 fSkipFailed = copy.fSkipFailed;
94 fValues = copy.fValues;
95 fNValues = copy.fNValues;
98 fCheckHistRange = copy.fCheckHistRange;
106 //__________________________________________________________________________________________________
107 AliRsnListOutput::~AliRsnListOutput()
111 // Deletes the output objects.
117 //__________________________________________________________________________________________________
118 void AliRsnListOutput::Reset()
121 // Clear all output objects. In general, only one will be initialized at a time.
127 //_____________________________________________________________________________
128 void AliRsnListOutput::AddValue(AliRsnValue *value)
131 // Adds a value computation object to the list.
134 fValues.AddLast(value);
138 //________________________________________________________________________________________
139 Bool_t AliRsnListOutput::Init(const char *prefix, TList *list)
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>
151 // all output objects are cleared
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();
158 AliError("Need at least 1 value");
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;
166 // resize the output array
167 fArray.Set(fNValues);
170 Bool_t isPair=kFALSE;
173 TString name(GetName());
174 if (!name.CompareTo("pair")) isPair = kTRUE;
175 if (isPair) name = "";
176 else name.Prepend(".");
177 name.Prepend(prefix);
179 // TString name(Form("%s.%s", prefix, GetName()));
180 AliRsnValue *val = 0x0;
181 for (i = 0; i < fNValues; i++) {
184 AliError(Form("Slot %d in value list is NULL", i));
189 name += val->GetName();
194 TObject *object = 0x0;
196 // initialize appropriate output object
197 // and, if successful, insert into the list
200 // name.Append("_hist");
201 object = CreateHistogram(name.Data());
204 // name.Append("_hsparse");
205 object = CreateHistogramSparse(name.Data());
208 // name.Append("_cf");
209 object = CreateCFContainer(name.Data());
212 AliWarning("Wrong type output or initialization failure");
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()));
219 fIndex = fList->IndexOf(object);
226 //________________________________________________________________________________________
227 TH1 *AliRsnListOutput::CreateHistogram(const char *name)
230 // Initialize the 'default' TH1 output object.
231 // In case one of the expected axes is NULL, the initialization fails.
234 // we expect to have maximum 3 axes in this case
235 Int_t i, nbins[3] = {0, 0, 0};
237 for (i = 0; i < fNValues; i++) {
238 AliRsnValue *val = GetValue(i);
240 AliError(Form("Expected axis %d is NULL", i));
243 nbins[i] = GetValue(i)->GetArray().GetSize() - 1;
244 array[i] = GetValue(i)->GetArray();
249 // create histogram depending on the number of axes
252 hist = new TH1F(name, "", nbins[0], array[0].GetArray());
255 hist = new TH2F(name, "", nbins[0], array[0].GetArray(), nbins[1], array[1].GetArray());
258 hist = new TH3F(name, "", nbins[0], array[0].GetArray(), nbins[1], array[1].GetArray(), nbins[2], array[2].GetArray());
261 //AliError(Form("Wrong number of dimensions: %d", fNValues))
265 if (hist) hist->Sumw2();
269 //________________________________________________________________________________________
270 THnSparseF *AliRsnListOutput::CreateHistogramSparse(const char *name)
273 // Initialize the THnSparse output object.
274 // In case one of the expected axes is NULL, the initialization fails.
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();
288 THnSparseF *hist = new THnSparseF(name, "", fNValues, nbins);
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++) {
295 if (val) hist->GetAxis(i)->SetName(val->GetName());
296 hist->GetAxis(i)->Set(nbins[i], array[i].GetArray());
306 //________________________________________________________________________________________
307 AliCFContainer *AliRsnListOutput::CreateCFContainer(const char *name)
310 // Initialize the AliCFContainer output object.
311 // In case one of the expected axes is NULL, the initialization fails.
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();
325 AliCFContainer *cont = new AliCFContainer(name, "", fSteps, fNValues, nbins);
327 // set the bin limits for each axis
328 for (i = 0; i < fNValues; i++) {
329 cont->SetBinLimits(i, array[i].GetArray());
339 //________________________________________________________________________________________
340 Bool_t AliRsnListOutput::Fill(TObject *target, Int_t step)
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.
352 Bool_t globalOK = kTRUE;
353 AliRsnValue *val = 0x0;
354 for (i = 0; i < fNValues; i++) {
357 AliError("NULL value found");
360 globalOK = globalOK && val->Eval(target);
361 fArray[i] = (Double_t)val->GetComputedValue();
363 if (!globalOK && fSkipFailed) return kFALSE;
366 if (!fList || fIndex < 0) {
367 AliError("List not initialized");
370 TObject *obj = fList->At(fIndex);
372 AliError("Null pointer");
377 //AliInfo(Form("[%s] Object index, name, type = %d, %s (%s)", GetName(), fIndex, obj->GetName(), obj->ClassName()));
380 if (obj->IsA() == TH1F::Class()) {
381 TH1F *h = (TH1F *)obj;
384 } else if (obj->IsA() == TH2F::Class()) {
385 TH2F *h = (TH2F *)obj;
386 h->Fill(fArray[0], fArray[1]);
388 } else if (obj->IsA() == TH3F::Class()) {
389 TH3F *h = (TH3F *)obj;
390 h->Fill(fArray[0], fArray[1], fArray[2]);
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;
399 h->Fill(fArray.GetArray());
401 } else if (obj->InheritsFrom(AliCFContainer::Class())) {
402 AliCFContainer *c = (AliCFContainer *)obj;
403 c->Fill(fArray.GetArray(), step);
406 AliError(Form("Not handled class '%s'", obj->ClassName()));
411 //________________________________________________________________________________________
412 Bool_t AliRsnListOutput::Fill(AliRsnEvent *ev, AliRsnDaughter *d)
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.
421 if (!fList || fIndex < 0) {
422 AliError("List not initialized");
425 TObject *obj = fList->At(fIndex);
427 AliError("Null pointer");
432 TIter next(&fValues);
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);
445 values[i] = (Double_t)val->GetComputedValue();
446 if (!globalOK && fSkipFailed) return kFALSE;
451 if (obj->IsA() == TH1F::Class()) {
452 TH1F *h = (TH1F *)obj;
455 } else if (obj->IsA() == TH2F::Class()) {
456 TH2F *h = (TH2F *)obj;
457 h->Fill(values[0], values[1]);
459 } else if (obj->IsA() == TH3F::Class()) {
460 TH3F *h = (TH3F *)obj;
461 h->Fill(values[0], values[1], values[2]);
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;
473 AliError(Form("Not handled class '%s'", obj->ClassName()));