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