]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnFunction.cxx
d417ea27c350244152a358dc0745a0afdda19537
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnFunction.cxx
1 //
2 // Class AliRsnFunction
3 //
4 // This class defines a base classe to implement a function
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 functions 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 <TString.h>
24
25 #include "AliLog.h"
26
27 #include "AliRsnDaughter.h"
28 #include "AliRsnEvent.h"
29 #include "AliRsnPairDef.h"
30 #include "AliRsnCut.h"
31
32 #include "AliRsnFunction.h"
33
34 ClassImp(AliRsnFunction)
35
36 //________________________________________________________________________________________
37 AliRsnFunction::AliRsnFunction() :
38   TNamed(),
39   fFcnType(kFcnTypes),
40   fPairDef(0x0),
41   fTrack(0x0),
42   fPair(0x0),
43   fEvent(0x0),
44   fHistogram(0x0)
45 {
46 //
47 // Constructor for 1D functions.
48 // Requires only the binning of the output function,
49 // which is stored as 'main' histoDef in fHistoDef[0]
50 //
51
52   fBinType[0] = kNoBins;
53   fBinType[1] = kNoBins;
54
55   fHistoDef[0] = 0x0;
56   fHistoDef[1] = 0x0;
57   fHistoDef[2] = 0x0;
58 }
59
60 //________________________________________________________________________________________
61 AliRsnFunction::AliRsnFunction
62 (EFcnType fcnType, AliRsnHistoDef *hd) :
63   TNamed(),
64   fFcnType(fcnType),
65   fPairDef(0x0),
66   fTrack(0x0),
67   fPair(0x0),
68   fEvent(0x0),
69   fHistogram(0x0)
70 {
71 //
72 // Constructor for 1D functions.
73 // Requires only the binning of the output function,
74 // which is stored as 'main' histoDef in fHistoDef[0]
75 //
76
77   fBinType[0] = kNoBins;
78   fBinType[1] = kNoBins;
79
80   fHistoDef[0] = hd;
81   fHistoDef[1] = 0x0;
82   fHistoDef[2] = 0x0;
83
84   DefineName();
85 }
86
87 //________________________________________________________________________________________
88 AliRsnFunction::AliRsnFunction
89 (EFcnType fcnType, EBinType binType, AliRsnHistoDef *hdMain, AliRsnHistoDef *hdBin) :
90   TNamed(),
91   fFcnType(fcnType),
92   fPairDef(0x0),
93   fTrack(0x0),
94   fPair(0x0),
95   fEvent(0x0),
96   fHistogram(0x0)
97 {
98 //
99 // Constructor for 2D functions.
100 // Requires the binning of the output function,
101 // which is stored as 'main' histoDef in fHistoDef[0],
102 // and a definition for a secondary binning, stored in fHistoDef[1]
103 //
104
105   fBinType[0] = binType;
106   fBinType[1] = kNoBins;
107
108   fHistoDef[0] = hdMain;
109   fHistoDef[1] = hdBin;
110   fHistoDef[2] = 0x0;
111
112   DefineName();
113 }
114
115 //________________________________________________________________________________________
116 AliRsnFunction::AliRsnFunction
117 (EFcnType fcnType, EBinType binType1, EBinType binType2,
118  AliRsnHistoDef *hdMain, AliRsnHistoDef *hdBin1, AliRsnHistoDef *hdBin2) :
119   fFcnType(fcnType),
120   fPairDef(0x0),
121   fTrack(0x0),
122   fPair(0x0),
123   fEvent(0x0),
124   fHistogram(0x0)
125 {
126 //
127 // Constructor for 3D functions.
128 // Requires the binning of the output function,
129 // which is stored as 'main' histoDef in fHistoDef[0],
130 // and a definition for two secondary binnings, stored in fHistoDef[1,2]
131 //
132
133   fBinType[0] = binType1;
134   fBinType[1] = binType2;
135
136   fHistoDef[0] = hdMain;
137   fHistoDef[1] = hdBin1;
138   fHistoDef[2] = hdBin2;
139
140   DefineName();
141 }    
142
143 //________________________________________________________________________________________
144 AliRsnFunction::AliRsnFunction(const AliRsnFunction &copy) :
145   TNamed(copy),
146   fFcnType(copy.fFcnType),
147   fPairDef(copy.fPairDef),
148   fTrack(copy.fTrack),
149   fPair(copy.fPair),
150   fEvent(copy.fEvent),
151   fHistogram(0x0)
152 {
153 //
154 // Copy constructor.
155 //
156
157   Int_t i;
158   for (i = 0; i < 3; i++) {
159     fHistoDef[i] = copy.fHistoDef[i];
160     if (i < 2) fBinType[i] = copy.fBinType[i];
161   }
162
163   DefineName();
164 }
165
166 //________________________________________________________________________________________
167 const AliRsnFunction& AliRsnFunction::operator=(const AliRsnFunction& copy)
168 {
169 //
170 // Assignment operator.
171 //
172
173   SetName(copy.GetName());
174   SetTitle(copy.GetTitle());
175
176   fFcnType = copy.fFcnType;
177
178   fPairDef = copy.fPairDef;
179
180   Int_t i;
181   for (i = 0; i < 3; i++) {
182     fHistoDef[i] = copy.fHistoDef[i];
183     if (i < 2) fBinType[i] = copy.fBinType[i];
184   }
185
186   fTrack = copy.fTrack;
187   fPair = copy.fPair;
188   fEvent = copy.fEvent;
189
190   if (fHistogram) delete fHistogram;
191   fHistogram = 0x0;
192
193   DefineName();
194
195   return (*this);
196 }
197
198 //________________________________________________________________________________________
199 const char* AliRsnFunction::FcnName()
200 {
201 //
202 // Defines the name of this object according to
203 // the function type and binning
204 //
205
206   switch (fFcnType)
207   {
208     case kTrackPt:
209       return  "TRKPT";
210       break;
211     case kTrackEta:
212       return  "TRKETA";
213       break;
214     case kInvMass:
215       return  "IM";
216       break;
217     case kInvMassMC:
218       return  "IMMC";
219       break;
220     case kResolution:
221       return  "RES";
222       break;
223     case kPairPt:
224       return  "PT";
225       break;
226     case kPairEta:
227       return  "ETA";
228       break;
229     case kEventMult:
230       return  "MULT";
231       break;
232     default:
233       return  "UNDEF";
234   }
235 }
236
237 //________________________________________________________________________________________
238 void AliRsnFunction::DefineName()
239 {
240 //
241 // Defines the name of this object according to
242 // the function type and binning
243 //
244
245   Int_t  dim = CheckDim();
246
247   switch (dim)
248   {
249     case 1:
250       SetName(FcnName());
251       break;
252     case 2:
253       SetName(Form("%s_%s", FcnName(), BinName(fBinType[0])));
254       break;
255     case 3:
256       SetName(Form("%s_%s_%s", FcnName(), BinName(fBinType[0]), BinName(fBinType[1])));
257       break;
258     default:
259       SetName("UNDEF");
260   }
261 }
262
263 //________________________________________________________________________________________
264 Double_t AliRsnFunction::Eval()
265 {
266 //
267 // Compute value for functions with 'event' argument type
268 //
269
270   Double_t value;
271
272   switch (fFcnType)
273   {
274     case kTrackPt:
275       return fTrack->Pt();
276     case kTrackEta:
277       return fTrack->Eta();
278     case kInvMass:
279       return fPair->GetInvMass(fPairDef->GetMass(0), fPairDef->GetMass(1));
280     case kInvMassMC:
281       return fPair->GetInvMassMC(fPairDef->GetMass(0), fPairDef->GetMass(1));
282     case kResolution:
283       value  = fPair->GetInvMass(fPairDef->GetMass(0), fPairDef->GetMass(1));
284       value -= fPair->GetInvMassMC(fPairDef->GetMass(0), fPairDef->GetMass(1));
285       value /= fPair->GetInvMassMC(fPairDef->GetMass(0), fPairDef->GetMass(1));
286       return value;
287     case kPairPt:
288       return fPair->GetPt();
289     case kPairEta:
290       return fPair->GetEta();
291     case kEventMult:
292       return fEvent->GetMultiplicity();
293     default:
294       AliWarning("Function type not supported");
295       return -999.0;
296   }
297 }
298
299 //_________________________________________________________________________________________________
300 Bool_t AliRsnFunction::CheckInput(Option_t *option)
301 {
302 //
303 // Checks if the argument type is coherent with
304 // the function type required
305 //
306
307   TString opt(option);
308   opt.ToUpper();
309
310   if (opt.Contains("TRACK")) {
311     if (!fTrack) {
312       AliError("Input track object is NULL");
313       return kFALSE;
314     }
315   }
316
317   if (opt.Contains("PAIR")) {
318     if (!fPair) {
319       AliError("Input pair object is NULL");
320       return kFALSE;
321     }
322   }
323
324   if (opt.Contains("EVENT")) {
325     if (!fEvent) {
326       AliError("Input event object is NULL");
327       return kFALSE;
328     }
329   }
330
331   return kTRUE;
332 }
333
334 //________________________________________________________________________________________
335 TH1* AliRsnFunction::CreateHistogram(const char *histoName, const char *histoTitle)
336 {
337 //
338 // Creates and returns the histogram defined using
339 // arguments fo name and title, and the first histoDef for binning.
340 // Variable-sized histogram binning is always called, due to use of histoDef,
341 // even if the bins are equal, since they are defined in this class.
342 // Eventually present histoDef's in other slots of array (1, 2) are ignored.
343 //
344
345   // first binning is required
346   if (!fHistoDef[0]) return 0;
347
348   // retrieve binnings for main and secondary axes
349   Int_t    i, nbins[3] = {0, 0, 0};
350   Double_t min[3] = {0., 0., 0.}, max[3] = {0., 0., 0.};
351   for (i = 0; i < 3; i++)
352   {
353     if (fHistoDef[i])
354     {
355       nbins[i] = fHistoDef[i]->GetNBins();
356       min[i] = fHistoDef[i]->GetMin();
357       max[i] = fHistoDef[i]->GetMax();
358     }
359   }
360     
361   // define the kind of output according to the number of histoDefs
362   if (fHistogram) delete fHistogram;
363   if (!nbins[1] && !nbins[2]) {
364     fHistogram = new TH1D(histoName, histoTitle, nbins[0], min[0], max[0]);
365     fHistogram->SetXTitle(FcnName());
366   }
367   else if (nbins[1] > 0 && !nbins[2]) {
368     fHistogram = new TH2D(histoName, histoTitle, nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
369     fHistogram->SetXTitle(FcnName());
370     fHistogram->SetYTitle(BinName(fBinType[0]));
371   }
372   else {
373     fHistogram = new TH3D(histoName, histoTitle, nbins[0], min[0], max[0], nbins[1], min[1], max[1], nbins[2], min[2], max[2]);
374     fHistogram->SetXTitle(FcnName());
375     fHistogram->SetYTitle(BinName(fBinType[0]));
376     fHistogram->SetZTitle(BinName(fBinType[1]));
377   }
378   
379   fHistogram->Sumw2();
380
381   return fHistogram;
382 }
383
384 //________________________________________________________________________________________
385 Bool_t AliRsnFunction::Fill()
386 {
387 //
388 // Fill function histogram with values computed from given input object.
389 //
390   AliDebug(AliLog::kDebug +2,"->");
391
392   // checks coherence between fcn type and passed argument
393   switch (fFcnType) {
394     case kTrackPt:
395     case kTrackEta:
396       if (!CheckInput("TRACK")) return kFALSE;
397       break;
398     case kInvMass:
399     case kInvMassMC:
400     case kResolution:
401     case kPairPt:
402     case kPairEta:
403       if (!CheckInput("PAIR")) return kFALSE;
404       break;
405     case kEventMult:
406       if (!CheckInput("EVENT")) return kFALSE;
407       break;
408     default:
409       AliError(Form("Input type %d not defined", (Int_t)fFcnType));
410       return kFALSE;
411   }
412
413   // check presence of output histogram
414   if (!fHistogram) {
415     AliError("Histogram is not yet initialized");
416     return kFALSE;
417   }
418
419   // compute value and stores into histogram
420   Int_t    dim = CheckDim();
421   Double_t mainValue, binValue[2];
422
423   TH1D *h1 = dynamic_cast<TH1D*>(fHistogram);
424   TH2D *h2 = dynamic_cast<TH2D*>(fHistogram);
425   TH3D *h3 = dynamic_cast<TH3D*>(fHistogram);
426   
427   mainValue = Eval();
428
429   switch (dim)
430   {
431     case 1:
432       if (h1) h1->Fill(mainValue);
433       break;
434     case 2:
435       binValue[0] = BinValue(fBinType[0]);
436       if (h2) h2->Fill(mainValue, binValue[0]);
437       break;
438     case 3:
439       binValue[0] = BinValue(fBinType[0]);
440       binValue[1] = BinValue(fBinType[1]);
441       if (h3) h3->Fill(mainValue, binValue[0], binValue[1]);
442       break;
443     default:
444       AliError("Wrong number of dimensions in the histogram. Check HD initialization");
445       return kFALSE;
446   }
447
448   AliDebug(AliLog::kDebug +2,"->");
449   return kTRUE;
450 }
451
452 //________________________________________________________________________________________
453 Double_t AliRsnFunction::BinValue(EBinType binType)
454 {
455 //
456 // Computes the value for binning from the argument.
457 // For each kind of binning type, the object is expected
458 // to be of a given type, otherwise an error is raised.
459 //
460
461   // checks coherence between bin type and passed argument
462   switch (binType) {
463     case kBinPairPt:
464       if (!CheckInput("PAIR")) return 0.0;
465       return fPair->GetPt();
466     case kBinPairEta:
467       if (!CheckInput("PAIR")) return 0.0;
468       return fPair->GetEta();
469     case kBinEventMult:
470       if (!CheckInput("EVENT")) return 0.0;
471       return fEvent->GetMultiplicity();
472     default:
473       AliError(Form("%s: Binning type not defined", GetName()));
474       return 0.0;
475   }
476 }
477
478 //________________________________________________________________________________________
479 Int_t AliRsnFunction::CheckDim()
480 {
481 //
482 // Checks number of dimensions.
483 // Makes sure that eventual binnings are coherent and well defined
484 //
485
486   if (!fHistoDef[0]) return 0;
487   if (fHistoDef[0] && !fHistoDef[1] && !fHistoDef[2]) return 1;
488   if (fHistoDef[0] && fHistoDef[1] && !fHistoDef[2]) return 2;
489
490   return 3;
491 }
492
493 //________________________________________________________________________________________
494 const char* AliRsnFunction::BinName(EBinType binType)
495 {
496 //
497 // Defines the name of binning
498 //
499
500   switch (binType)
501   {
502     case kBinPairPt:
503       return "PT";
504       break;
505     case kBinPairEta:
506       return "ETA";
507       break;
508     case kBinEventMult:
509       return "MULT";
510       break;
511     default:
512       return "UNDEF";
513   }
514 }