]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muondep/AliAnalysisMuMuSpectra.cxx
move track selection differences from AOD and ESD to corresponding reader
[u/mrichter/AliRoot.git] / PWG / muondep / AliAnalysisMuMuSpectra.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 // AliAnalysisMuMuSpectra : a class to encapsulate results from MuMu analysis
18 //
19 // Spectra can be merged and converted into histograms
20 //
21 // author: L. Aphecetche (Subatech)
22 //
23
24 #include "AliAnalysisMuMuSpectra.h"
25
26 #include "AliLog.h"
27 #include "AliAnalysisMuMuBinning.h"
28 #include "AliAnalysisMuMuJpsiResult.h"
29 #include "Riostream.h"
30 #include "TH1.h"
31 #include "TList.h"
32 #include "TObjArray.h"
33
34 ClassImp(AliAnalysisMuMuSpectra)
35
36 //______________________________________________________________________________
37 AliAnalysisMuMuSpectra::AliAnalysisMuMuSpectra(const char* name, const char* title) :
38 TNamed(name,title),
39 fBinning(0x0),
40 fBins(0x0),
41 fWeight(1.0)
42 {
43   // default ctor
44 }
45
46 //______________________________________________________________________________
47 AliAnalysisMuMuSpectra::AliAnalysisMuMuSpectra(const AliAnalysisMuMuSpectra& rhs)
48 : TNamed(rhs.GetName(),rhs.GetTitle()),
49 fBinning(0x0),
50 fBins(0x0),
51 fWeight(rhs.Weight())
52 {
53   // copy ctor
54
55   if ( rhs.fBinning )
56   {
57     fBinning = new AliAnalysisMuMuBinning(*rhs.fBinning);
58   }
59
60   TIter next(rhs.fBins);
61   AliAnalysisMuMuResult* bin;
62   
63   while ( ( bin = static_cast<AliAnalysisMuMuResult*>(next()) ) )
64   {
65     if (!fBins)
66     {
67       fBins = new TObjArray;
68       fBins->SetOwner(kTRUE);
69     }
70     fBins->Add(bin);
71   }
72   
73   
74 }
75
76 //______________________________________________________________________________
77 AliAnalysisMuMuSpectra&
78 AliAnalysisMuMuSpectra::operator=(const AliAnalysisMuMuSpectra& rhs)
79 {
80   // assignment operator
81   
82   if (this==&rhs) return *this;
83   
84   delete fBinning;
85   fBinning = 0x0;
86   delete fBins;
87   fBins = 0x0;
88   
89   if ( rhs.fBinning )
90   {
91     fBinning = new AliAnalysisMuMuBinning(*rhs.fBinning);
92   }
93   
94   TIter next(rhs.fBins);
95   AliAnalysisMuMuResult* bin;
96   
97   while ( ( bin = static_cast<AliAnalysisMuMuResult*>(next()) ) )
98   {
99     if (!fBins)
100     {
101       fBins = new TObjArray;
102       fBins->SetOwner(kTRUE);
103     }
104     fBins->Add(bin);
105   }
106   
107   fWeight = rhs.Weight();
108   
109   return *this;
110 }
111
112 //______________________________________________________________________________
113 AliAnalysisMuMuSpectra::~AliAnalysisMuMuSpectra()
114 {
115   // dtor
116   delete fBinning;
117   delete fBins;
118 }
119
120 //______________________________________________________________________________
121 void AliAnalysisMuMuSpectra::AdoptResult(const AliAnalysisMuMuBinning::Range& bin,
122                                          AliAnalysisMuMuResult* result)
123 {
124   // adopt (i.e. we are becoming the owner) a result for a given bin
125   if (!fBinning)
126   {
127     fBinning = new AliAnalysisMuMuBinning;
128     fBins = new TObjArray;
129     fBins->SetOwner(kTRUE);
130   }
131   fBinning->AddBin(bin);
132   fBins->Add(result);
133 }
134
135 //______________________________________________________________________________
136 Bool_t AliAnalysisMuMuSpectra::Correct(const AliAnalysisMuMuSpectra& accEff, const char* particle, const char* subResultName)
137 {
138   /// Correct this spectra by acceff one
139   if (IsEmpty())
140   {
141     AliError("Cannot correct an empty spectra !");
142     return kFALSE;
143   }
144   
145   // check we have the same binning first
146   AliAnalysisMuMuBinning* accEffBins = accEff.Binning();
147   
148   if ( !fBinning->IsEqual(accEffBins) )
149   {
150     AliError("Cannot correct with a spectra which does not have the same binning");
151     return kFALSE;
152   }
153   
154   TObjArray* particles = fBinning->CreateWhatArray();
155   TObjArray* types = fBinning->CreateQuantityArray();
156   
157   if (particles->GetEntries()!=1 || types->GetEntries()!=1 )
158   {
159     delete particles;
160     delete types;
161     return kFALSE;
162   }
163   
164   TObjArray* bins = accEff.BinContentArray();
165
166   
167   for ( Int_t i = 0; i < bins->GetEntries(); ++i )
168   {
169     AliAnalysisMuMuJpsiResult* thisResult = static_cast<AliAnalysisMuMuJpsiResult*>(fBins->At(i));
170     AliAnalysisMuMuJpsiResult* accResult = static_cast<AliAnalysisMuMuJpsiResult*>(bins->At(i));
171 //    AliInfoClass(Form("i=%d",i ));
172 //    StdoutToAliInfoClass(thisResult->Print("full");
173 //                         std::cout << "----" << std::endl;
174 //                         accResult->Print("full"));
175     
176     thisResult->Correct(*accResult,particle,subResultName);
177   }
178   
179   delete particles;
180   delete types;
181   return kTRUE;
182 }
183
184 //______________________________________________________________________________
185 AliAnalysisMuMuResult*
186 AliAnalysisMuMuSpectra::GetResultForBin(const AliAnalysisMuMuBinning::Range& bin) const
187 {
188   /// Get result for a given bin
189   /// Warning: this requires a loop on bins
190   
191   if ( IsEmpty() ) return 0x0;
192   
193   TObjArray* bins = fBinning->CreateBinObjArray();
194   
195   Int_t j(-1);
196   
197   StdoutToAliDebug(1,std::cout << "searching for "; bin.Print());
198   
199   for ( Int_t i = 0; i <= bins->GetLast() && j < 0 ; ++i )
200   {
201     AliAnalysisMuMuBinning::Range* b = static_cast<AliAnalysisMuMuBinning::Range*>(bins->At(i));
202
203     StdoutToAliDebug(1,b->Print(););
204     
205     if ( bin == *b )
206     {
207       j = i;
208     }
209   }
210   
211   delete bins;
212   
213   if (j>=0)
214   {
215     return static_cast<AliAnalysisMuMuResult*>(fBins->At(j));
216   }
217   else
218   {
219     StdoutToAliDebug(1,std::cout << "Could not find result for bin:" << std::endl; bin.Print(););
220   }
221   return 0x0;
222 }
223
224 //______________________________________________________________________________
225 Bool_t AliAnalysisMuMuSpectra::HasValue(const char* what) const
226 {
227     // whether or not our result(s) has a given property
228   if ( IsEmpty() ) return kFALSE;
229   
230   AliAnalysisMuMuResult* r = static_cast<AliAnalysisMuMuResult*>(fBins->First());
231   
232   return r->HasValue(what);
233 }
234
235 //______________________________________________________________________________
236 Bool_t AliAnalysisMuMuSpectra::IsEmpty() const
237 {
238   // whether this spectra is empty or not
239   return ( fBins==0x0 || fBins->GetEntries()<=0 );
240 }
241
242 //______________________________________________________________________________
243 Long64_t AliAnalysisMuMuSpectra::Merge(TCollection* list)
244 {
245   /// Merge method
246   
247   // Merge a list of AliAnalysisMuMuSpectra objects with this
248   // Returns the number of merged objects (including this).
249   
250   if (!list) return 0;
251   
252   if (list->IsEmpty()) return 1;
253   
254   TIter next(list);
255   TObject* currObj;
256   Int_t count(0);
257   
258   TList binningList;
259   
260   // for each bin must do a list of results, and merge that list
261   
262   TObjArray* bins = fBinning->CreateBinObjArray();
263   TIter nextBin(bins);
264   AliAnalysisMuMuBinning::Range* bin;
265
266   Int_t i(0);
267   
268   while ( ( bin = static_cast<AliAnalysisMuMuBinning::Range*>(nextBin()) ) )
269   {
270     next.Reset();
271     
272     TList binList;
273     
274     while ( ( currObj = next() ) )
275     {
276       AliAnalysisMuMuSpectra* spectra = static_cast<AliAnalysisMuMuSpectra*>(currObj);
277
278       if (i==0)
279       {
280         binningList.Add(spectra->Binning());
281   
282         if ( !fBinning->IsEqual(spectra->Binning()) || spectra->BinContentArray()->GetLast() != BinContentArray()->GetLast() )
283         {
284           AliError("Cannot merge spectra with different binning");
285           continue;
286         }
287         ++count;
288       }
289       
290       binList.Add(spectra->GetResultForBin(*bin));
291     }
292     
293     ++i;
294     
295     AliAnalysisMuMuResult* r = static_cast<AliAnalysisMuMuResult*>(GetResultForBin(*bin));
296     r->Merge(&binList);
297     
298   }
299
300   delete bins;
301   
302   return count+1;
303 }
304
305 //_____________________________________________________________________________
306 TH1* AliAnalysisMuMuSpectra::Plot(const char* what, const char* subresult, Bool_t divideByBinWidth) const
307 {
308   // Convert the spectra into an histogram
309   
310   TString swhat(what);
311   swhat.ToUpper();
312   
313   Double_t* bins = fBinning->CreateBinArray();
314   
315   TObjArray* binArray = fBinning->CreateBinObjArray();
316   
317   TH1* h(0x0);
318   
319   AliDebug(1,Form("nbins=%d nresults=%d",binArray->GetEntries(),fBins->GetEntries()));
320   
321   for ( Int_t j = 0; j < TMath::Min(binArray->GetEntries(),fBins->GetEntries()); ++j )
322   {
323     AliAnalysisMuMuJpsiResult* r = static_cast<AliAnalysisMuMuJpsiResult*>(fBins->At(j));
324     
325     if ( strlen(subresult) > 0 && r->SubResults() )
326     {
327       TString sub(subresult);
328       sub.ToUpper();
329       r = static_cast<AliAnalysisMuMuJpsiResult*>(r->SubResult(sub.Data()));
330       if (!r) continue;
331     }
332     
333     const AliAnalysisMuMuBinning::Range& b = r->Bin();
334     
335     if (!h)
336     {
337       h = new TH1F(r->GetName(),r->GetName(),binArray->GetEntries(),bins);
338       h->SetDirectory(0);
339     }
340     
341     Double_t y = r->GetValue(what);
342     Double_t yerr = r->GetErrorStat(what);
343   
344     if ( divideByBinWidth && b.WidthX()>0 )
345     {
346       y /= (b.WidthX());
347       yerr /= (b.WidthX());
348     }
349     
350     if (!TMath::Finite(y)) y = 0.0;
351     if (!TMath::Finite(yerr)) yerr = 0.0;
352
353     std::cout << b.AsString();
354     r->PrintValue(swhat.Data(),"",y,yerr);
355     
356     h->SetBinContent(j+1,y);
357     h->SetBinError(j+1,yerr);
358     
359   }
360   
361   delete binArray;
362   delete[] bins;
363   
364   return h;
365 }
366
367
368 //______________________________________________________________________________
369 void AliAnalysisMuMuSpectra::Print(Option_t* opt) const
370 {
371   // printout
372   
373   if (!IsEmpty())
374   {
375     TString sopt(opt);
376     Int_t nmax = sopt.Atoi();
377     if ( nmax <= 0 ) nmax = fBins->GetEntries();
378     for ( Int_t i = 0; i < nmax; ++i )
379     {
380       AliAnalysisMuMuResult* r = static_cast<AliAnalysisMuMuResult*>(fBins->At(i));
381       if (r) r->Print(opt);
382     }
383   }
384 }
385
386 //______________________________________________________________________________
387 void AliAnalysisMuMuSpectra::Scale(Double_t value)
388 {
389   // scale all bins by value
390
391   TIter next(fBins);
392   AliAnalysisMuMuResult* r;
393
394   while ( ( r = static_cast<AliAnalysisMuMuResult*>(next()) ) )
395   {
396     r->Scale(value);
397   }
398 }
399
400 //______________________________________________________________________________
401 void AliAnalysisMuMuSpectra::SetWeight(Double_t w)
402 {
403   // Set the weight of this spectra
404   fWeight = w;
405   TIter next(fBins);
406   AliAnalysisMuMuResult* r;
407   
408   while ( ( r = static_cast<AliAnalysisMuMuResult*>(next()) ) )
409   {
410     r->SetWeight(Weight());
411   }
412 }