]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnFunction.cxx
Update to the analysis on MC data only.
[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 <TString.h>
23
24 #include "AliLog.h"
25
26 #include "AliRsnDaughter.h"
27 #include "AliRsnEvent.h"
28 #include "AliRsnPairDef.h"
29 #include "AliRsnMother.h"
30 #include "AliRsnValue.h"
31
32 #include "AliRsnFunction.h"
33
34 ClassImp(AliRsnFunction)
35
36 //________________________________________________________________________________________
37 AliRsnFunction::AliRsnFunction(Bool_t useTH1) :
38     TObject(),
39     fPairDef(0x0),
40     fAxisList("AliRsnValue", 0),
41     fPair(0x0),
42     fEvent(0x0),
43     fUseTH1(useTH1),
44     fSize(0),
45     fH1(0x0),
46     fHSparse(0x0)
47 {
48 //
49 // Constructor.
50 //
51 }
52
53 //________________________________________________________________________________________
54 AliRsnFunction::AliRsnFunction(const AliRsnFunction &copy) :
55     TObject(copy),
56     fPairDef(copy.fPairDef),
57     fAxisList(copy.fAxisList),
58     fPair(copy.fPair),
59     fEvent(copy.fEvent),
60     fUseTH1(copy.fUseTH1),
61     fSize(copy.fSize),
62     fH1(0x0),
63     fHSparse(0x0)
64 {
65 //
66 // Copy constructor.
67 //
68 }
69
70 //________________________________________________________________________________________
71 const AliRsnFunction& AliRsnFunction::operator=(const AliRsnFunction& copy)
72 {
73 //
74 // Assignment operator.
75 //
76
77   //SetName(copy.GetName());
78   //SetTitle(copy.GetTitle());
79
80   fPairDef = copy.fPairDef;
81   fPair = copy.fPair;
82   fEvent = copy.fEvent;
83   fUseTH1 = copy.fUseTH1;
84   fSize = copy.fSize;
85
86   if (fH1) delete fH1;
87   fH1 = 0x0;
88   
89   if (fHSparse) delete fHSparse;
90   fHSparse = 0x0;
91
92   return (*this);
93 }
94
95 //________________________________________________________________________________________
96 const char* AliRsnFunction::GetName() const
97 {
98 //
99 // Defines the name of this object according to
100 // the function type and binning
101 //
102
103   TString name("");
104
105   TObjArrayIter next(&fAxisList);
106   AliRsnValue *axis = 0;
107
108   while ((axis = (AliRsnValue*)next())) {
109     if (name.Length() > 1) name += '_';
110     name += axis->GetName();
111   }
112
113   return name.Data();
114 }
115
116 //________________________________________________________________________________________
117 void AliRsnFunction::AddAxis(AliRsnValue *const axis)
118 {
119   AliDebug(AliLog::kDebug+2,"<-");
120   Int_t size = fAxisList.GetEntries();
121   new(fAxisList[size]) AliRsnValue(*axis);
122   AliDebug(AliLog::kDebug+2,"->");
123   
124   if (fAxisList.GetEntries() > 3)
125   {
126     AliWarning("A TH1-type output cannot add more than 3 axes: switching to THnSparse -- THIS COULD CAUSE VERY LARGE FILES!!!");
127     fUseTH1 = kFALSE;
128   }
129 }
130
131 //________________________________________________________________________________________
132 TH1* AliRsnFunction::CreateHistogram(const char *histoName, const char *histoTitle)
133 {
134 //
135 // Creates and returns the histogram defined using
136 // arguments fo name and title, and the first histoDef for binning.
137 // Variable-sized histogram binning is always called, due to use of histoDef,
138 // even if the bins are equal, since they are defined in this class.
139 // Eventually present histoDef's in other slots of array (1, 2) are ignored.
140 //
141 // This version produces a THnSparseF.
142 //
143
144   fSize = fAxisList.GetEntries();
145   if (!fSize) {
146     AliError("No axes defined");
147     return 0x0;
148   }
149   else if (fSize < 1 || fSize > 3)
150   {
151     AliError("Too few or too many axes defined");
152     return 0x0;
153   }
154
155   Int_t    *nbins = new Int_t   [fSize];
156   Double_t *min   = new Double_t[fSize];
157   Double_t *max   = new Double_t[fSize];
158
159   // retrieve binnings for main and secondary axes
160   AliRsnValue *fcnAxis = 0;
161   for (Int_t i = 0; i < fSize; i++) {
162     fcnAxis = (AliRsnValue*)fAxisList.At(i);
163     if (!fcnAxis) {
164       nbins[i] = 0;
165       min[i]   = 0.0;
166       max[i]   = 0.0;
167       AliError("Empty axis");
168       continue;
169     }
170     nbins[i] = fcnAxis->GetNBins();
171     min[i]   = fcnAxis->GetMin();
172     max[i]   = fcnAxis->GetMax();
173   }
174
175   // create histogram depending on the number of axes
176   switch (fSize)
177   {
178     case 1:
179       fH1 = new TH1F(histoName, histoTitle, nbins[0], min[0], max[0]);
180       break;
181     case 2:
182       fH1 = new TH2F(histoName, histoTitle, nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
183       break;
184     case 3:
185       fH1 = new TH3F(histoName, histoTitle, nbins[0], min[0], max[0], nbins[1], min[1], max[1], nbins[2], min[2], max[2]);
186       break;
187   }
188   fH1->Sumw2();
189
190   return fH1;
191 }
192
193 //________________________________________________________________________________________
194 THnSparseF* AliRsnFunction::CreateHistogramSparse(const char *histoName, const char *histoTitle)
195 {
196 //
197 // Creates and returns the histogram defined using
198 // arguments fo name and title, and the first histoDef for binning.
199 // Variable-sized histogram binning is always called, due to use of histoDef,
200 // even if the bins are equal, since they are defined in this class.
201 // Eventually present histoDef's in other slots of array (1, 2) are ignored.
202 //
203 // This version produces a THnSparseF.
204 //
205
206   fSize = fAxisList.GetEntries();
207   if (!fSize) {
208     AliError("No axes defined");
209     return 0x0;
210   }
211
212   Int_t    *nbins = new Int_t   [fSize];
213   Double_t *min   = new Double_t[fSize];
214   Double_t *max   = new Double_t[fSize];
215
216   // retrieve binnings for main and secondary axes
217   AliRsnValue *fcnAxis = 0;
218   for (Int_t i = 0; i < fSize; i++) {
219     fcnAxis = (AliRsnValue*)fAxisList.At(i);
220     if (!fcnAxis) {
221       nbins[i] = 0;
222       min[i]   = 0.0;
223       max[i]   = 0.0;
224       AliError("Empty axis");
225       continue;
226     }
227     nbins[i] = fcnAxis->GetNBins();
228     min[i]   = fcnAxis->GetMin();
229     max[i]   = fcnAxis->GetMax();
230   }
231
232   Int_t size = fAxisList.GetEntries();
233   if (!size) {
234     AliError("No axes defined");
235     return 0x0;
236   }
237
238   // create histogram
239   fHSparse = new THnSparseF(histoName, histoTitle, size, nbins, min, max);
240   fHSparse->Sumw2();
241   
242   // clean heap
243   delete [] nbins;
244   delete [] min;
245   delete [] max;
246
247   return fHSparse;
248 }
249
250
251 //________________________________________________________________________________________
252 Bool_t AliRsnFunction::Fill()
253 {
254 //
255 // Fill function histogram with values computed from given input object.
256 //
257
258   AliDebug(AliLog::kDebug +2,"->");
259
260   Int_t  i;
261   Double_t *values = new Double_t[fSize];
262
263   AliRsnValue *fcnAxis = 0;
264   for (i = 0; i < fSize; i++) {
265     fcnAxis = (AliRsnValue*)fAxisList.At(i);
266     if (!fcnAxis) {
267       values[i] = 0.0;
268       continue;
269     }
270     if (fcnAxis->Eval(fPair, fPairDef, fEvent)) values[i] = fcnAxis->GetValue();
271   }
272   
273   // fill histogram
274   if (fUseTH1)
275   {
276     // check presence of output histogram
277     if (!fH1) {
278       AliError("Required a TH1 whish is not initialized");
279       return kFALSE;
280     }
281     
282     // fill according to dimensions
283     switch (fSize)
284     {
285       case 1:
286         {
287           TH1F *h1 = (TH1F*)fH1;
288           h1->Fill(values[0]);
289         }
290         break;
291       case 2:
292         {
293           TH2F *h2 = (TH2F*)fH1;
294           h2->Fill(values[0], values[1]);
295         }
296         break;
297       case 3:
298         {
299           TH3F *h3 = (TH3F*)fH1;
300           h3->Fill(values[0], values[1], values[2]);
301         }
302         break;
303       default:
304         AliError(Form("Wrong size : %d", fSize));
305         return kFALSE;
306     }
307   }
308   else
309   {
310     // check presence of output histogram
311     if (!fHSparse) {
312       AliError("Required a THnSparseF which is not initialized");
313       return kFALSE;
314     }
315     
316     fHSparse->Fill(values);
317   }
318   
319   delete [] values;
320
321   AliDebug(AliLog::kDebug +2,"->");
322   return kTRUE;
323 }