]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithFittingQDistribution.cxx
fixing coding violations from FC
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowAnalysisWithFittingQDistribution.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2008, 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  * integrated flow estimate by  *
18  *   fitting q-distribution     * 
19  *                              *
20  * author: Ante Bilandzic       * 
21  *          (anteb@nikhef.nl)   *
22  *                              *  
23  *  based on the macro written  *
24  *     by Sergei Voloshin       *
25  *******************************/  
26
27 #define AliFlowAnalysisWithFittingQDistribution_cxx
28
29 #include "Riostream.h"
30 #include "AliFlowCommonConstants.h"
31 #include "AliFlowCommonHist.h"
32 #include "AliFlowCommonHistResults.h"
33 #include "TChain.h"
34 #include "TFile.h"
35 #include "TList.h" 
36 #include "TF1.h"
37 #include "TLegend.h"
38 #include "TParticle.h"
39 #include "TProfile.h"
40 #include "AliFlowEventSimple.h"
41 #include "AliFlowTrackSimple.h"
42 #include "AliFlowAnalysisWithFittingQDistribution.h"
43
44 class TH1;
45 class TGraph;
46 class TPave;
47 class TLatex;
48 class TMarker;
49 class TObjArray;
50 class TList;
51 class TSystem;
52 class TROOT;
53 class AliFlowVector;
54 class TVector;
55
56 //================================================================================================================
57
58 ClassImp(AliFlowAnalysisWithFittingQDistribution)
59
60 AliFlowAnalysisWithFittingQDistribution::AliFlowAnalysisWithFittingQDistribution():  
61  fHistList(NULL),
62  fCommonHists(NULL),
63  fCommonHistsResults(NULL),
64  fnBinsPhi(0),
65  fPhiMin(0),
66  fPhiMax(0),
67  fPhiBinWidth(0),
68  fnBinsPt(0),
69  fPtMin(0),
70  fPtMax(0),
71  fPtBinWidth(0),
72  fnBinsEta(0),
73  fEtaMin(0),
74  fEtaMax(0),
75  fEtaBinWidth(0),
76  fHarmonic(2),
77  fAnalysisLabel(NULL),
78  fWeightsList(NULL),
79  fUsePhiWeights(kFALSE),
80  fUsePtWeights(kFALSE),
81  fUseEtaWeights(kFALSE),
82  fUseParticleWeights(NULL),
83  fPhiWeights(NULL),
84  fPtWeights(NULL),
85  fEtaWeights(NULL),
86  // fitting parameters with default values harwired here (use dedicated macro fqd.C to change them):
87  fFittingParameters(NULL), 
88  fTreshold(5),
89  fvStart(0.05),
90  fvMin(0.0),
91  fvMax(0.25),
92  fSigma2Start(0.75),
93  fSigma2Min(0.5), 
94  fSigma2Max(2.5),
95  fPlotResults(kFALSE),
96  // rest:
97  fLegend(NULL)
98  {
99   // constructor 
100   
101   // base list to hold all output objects:
102   fHistList = new TList();
103   fHistList->SetName("cobjFQD");
104   fHistList->SetOwner(kTRUE);
105   
106   // analysis label;
107   fAnalysisLabel = new TString();
108  
109   // list to hold histograms with phi, pt and eta weights:      
110   fWeightsList = new TList();
111
112   // initialize all arrays:  
113   this->InitializeArrays();
114
115  } // end of constructor
116  
117  
118 //================================================================================================================
119  
120
121 AliFlowAnalysisWithFittingQDistribution::~AliFlowAnalysisWithFittingQDistribution()
122 {
123  // desctructor
124  delete fHistList; 
125 }
126
127
128 //================================================================================================================
129
130
131 void AliFlowAnalysisWithFittingQDistribution::Init()
132 {
133  // access constants and book everything
134  
135  //save old value and prevent histograms from being added to directory
136  //to avoid name clashes in case multiple analaysis objects are used
137  //in an analysis
138  Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
139  TH1::AddDirectory(kFALSE);
140  
141  // access constants:
142  this->AccessConstants();
143  
144  // booking:
145  this->BookCommonHistograms();
146  this->BookAndFillWeightsHistograms();
147  this->BookEverythingForDistributions();
148  
149  // store fitting parameters:  
150  this->StoreFittingParameters();
151  
152  // nest lists:
153  fWeightsList->SetName("Weights");
154  fWeightsList->SetOwner(kTRUE);   
155  fHistList->Add(fWeightsList);
156  
157  // set harmonic in common control histograms (to be improved (should I do this somewhere else?)):
158  (fCommonHists->GetHarmonic())->Fill(0.5,fHarmonic); 
159  
160  //restore old status
161  TH1::AddDirectory(oldHistAddStatus);
162 } // end of void AliFlowAnalysisWithFittingQDistribution::Init()
163
164
165 //================================================================================================================
166
167
168 void AliFlowAnalysisWithFittingQDistribution::Make(AliFlowEventSimple* anEvent)
169 {
170  // loop over data
171  
172  // a) fill the common control histograms;
173  // b) loop over data and calculate non-weighted and weighted Q-vector and sum of particle weights;
174  // c) fill histograms for q-distribution;
175  // d) reset e-b-e quantities.
176  
177  // a) fill the common control histograms:
178  fCommonHists->FillControlHistograms(anEvent); 
179  
180  Double_t dPhi = 0.; // azimuthal angle in the laboratory frame
181  Double_t dPt  = 0.; // transverse momentum
182  Double_t dEta = 0.; // pseudorapidity
183
184  Double_t wPhi = 1.; // phi weight
185  Double_t wPt  = 1.; // pt weight
186  Double_t wEta = 1.; // eta weight
187  
188  Int_t nPrim = anEvent->NumberOfTracks();  // nPrim = total number of primary tracks, i.e. nPrim = nRP + nPOI + rest, where:
189                                            // nRP   = # of particles used to determine the reaction plane;
190                                            // nPOI  = # of particles of interest for a detailed flow analysis;
191                                            // rest  = # of particles which are not niether RPs nor POIs.  
192  
193  Int_t n = fHarmonic; // shortcut for the harmonic
194  
195  Double_t dReQ[2] = {0.}; // real part of Q-vector [0=particle weights not used, 1=particle weights used]
196  Double_t dImQ[2] = {0.}; // imaginary part of Q-vector [0=particle weights not used, 1=particle weights used]
197  Double_t dSumOfParticleWeights[2] = {0.}; // [0=particle weights not used, 1=particle weights used] 
198                                                                                                                                
199  AliFlowTrackSimple *aftsTrack = NULL;                                          
200  
201  // b) loop over data and calculate non-weighted and weighted Q-vector and sum of particle weights:                                          
202  for(Int_t i=0;i<nPrim;i++) 
203  { 
204   aftsTrack=anEvent->GetTrack(i);
205   if(aftsTrack)
206   {
207    if(!(aftsTrack->InRPSelection())) continue; // consider only tracks which are RPs
208      
209    dPhi = aftsTrack->Phi();
210    dPt  = aftsTrack->Pt();
211    dEta = aftsTrack->Eta();
212   
213    if(fUsePhiWeights && fPhiWeights && fnBinsPhi) // determine phi weight for this particle:
214    {
215     wPhi = fPhiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*fnBinsPhi/TMath::TwoPi())));
216    }
217    if(fUsePtWeights && fPtWeights && fnBinsPt) // determine pt weight for this particle:
218    {
219     wPt = fPtWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-fPtMin)/fPtBinWidth))); 
220    }            
221    if(fUseEtaWeights && fEtaWeights && fEtaBinWidth) // determine eta weight for this particle: 
222    {
223     wEta = fEtaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-fEtaMin)/fEtaBinWidth))); 
224    } 
225     
226    // calculate real and imaginary part of non-weighted and weighted Q-vector and sum of particle weights for this event:
227    for(Int_t pW=0;pW<1+(Int_t)(fUsePhiWeights||fUsePtWeights||fUseEtaWeights);pW++) // particle weights not used (0) or used (1)
228    {
229     // Q-vector:
230     dReQ[pW]+=pow(wPhi*wPt*wEta,pW)*TMath::Cos(n*dPhi); 
231     dImQ[pW]+=pow(wPhi*wPt*wEta,pW)*TMath::Sin(n*dPhi);
232     // sum of particle weights:
233     dSumOfParticleWeights[pW] += pow(wPhi*wPt*wEta,pW); // if pW = 0, this sum gives # of RPs, i.e. multiplicity
234    } 
235    
236   } // end of if(aftsTrack)
237  } // end of for(Int_t i=0;i<nPrim;i++)                                      
238                                            
239  // c) fill histograms for q-distribution:
240  // calculating first q = Q\sqrt{sum of particle weights} (Remark: if particle weights are unit than sum of particle weights = multiplicity)
241  Double_t q=0;                                          
242  for(Int_t pW=0;pW<1+(Int_t)(fUsePhiWeights||fUsePtWeights||fUseEtaWeights);pW++) // particle weights not used (0) or used (1)
243  {
244   if(dSumOfParticleWeights[pW])
245   {
246    q = pow(dReQ[pW]*dReQ[pW]+dImQ[pW]*dImQ[pW],0.5)/pow(dSumOfParticleWeights[pW],0.5);
247    // fill histograms:
248    fqDistribution[pW]->Fill(q,1.);
249    fSumOfParticleWeights[pW]->Fill(dSumOfParticleWeights[pW],1.);
250   }
251  } 
252  
253  // d) reset e-b-e quantities:
254  for(Int_t pW=0;pW<1+(Int_t)(fUsePhiWeights||fUsePtWeights||fUseEtaWeights);pW++) // particle weights not used (0) or used (1)
255  {
256   dReQ[pW] = 0.;
257   dImQ[pW] = 0.;
258   dSumOfParticleWeights[pW] = 0.;
259  }
260
261 } // end of Make()
262
263
264 //================================================================================================================
265
266
267 void AliFlowAnalysisWithFittingQDistribution::Finish(Bool_t doFit)
268 {
269  // calculate the final results
270  
271  // a) acces the constants and all fitting paremeters;
272  // b) access the flags for particle weights;
273  // c) do final fit;
274  // d) fill common hist results;
275  // e) print on the screen the final results.
276  
277  // a) access the constants and all fitting paremeters:
278  this->AccessConstants();
279  this->AccessFittingParameters();
280  
281  // b) access the flags for particle weights: 
282  fUsePhiWeights = (Bool_t)fUseParticleWeights->GetBinContent(1); 
283  fUsePtWeights = (Bool_t)fUseParticleWeights->GetBinContent(2); 
284  fUseEtaWeights = (Bool_t)fUseParticleWeights->GetBinContent(3);
285
286  // to be improved (moved somewhere else):
287  if(fPlotResults)
288  {
289   fLegend = new TLegend(0.6,0.55,0.85,0.7); 
290  }
291
292  // c) do final fit:             
293  if(doFit) 
294  {
295   // particle weights not used:
296   // 1) sigma^2 not fitted (fixed to 0.5):
297   this->DoFit(kFALSE,kFALSE);
298   // 2) sigma^2 fitted:
299   this->DoFit(kFALSE,kTRUE);
300   // particle weights used:
301   if(fUsePhiWeights||fUsePtWeights||fUseEtaWeights)      
302   {
303    // 1) sigma^2 not fitted (fixed to 0.5):
304    this->DoFit(kTRUE,kFALSE);  
305    // 2) sigma^2 fitted:
306    this->DoFit(kTRUE,kTRUE);  
307   }
308   
309   // d) fill common hist results (by default fill results obtained with sigma^2 fitted):
310   if(fUsePhiWeights||fUsePtWeights||fUseEtaWeights)
311   {
312    this->FillCommonHistResultsIntFlow(kTRUE,kTRUE); 
313   } else
314     {
315      this->FillCommonHistResultsIntFlow(kFALSE,kTRUE);    
316     } 
317   
318   // e) print on the screen the final results:
319   this->PrintFinalResultsForIntegratedFlow();  
320   
321  } // end of if(doFit)
322    
323 } // end of void AliFlowAnalysisWithFittingQDistribution::Finish(Bool_t doFit)
324
325
326 //================================================================================================================
327
328
329 void AliFlowAnalysisWithFittingQDistribution::GetOutputHistograms(TList *outputListHistos) 
330 {
331  // get pointers to all output histograms (called before Finish()) 
332  
333  if(outputListHistos)   
334  {   
335   // 1.) common control histograms and common histograms for final results:
336   TString commonHistName = "AliFlowCommonHistFQD";
337   commonHistName += fAnalysisLabel->Data();
338   AliFlowCommonHist *commonHist = dynamic_cast<AliFlowCommonHist*>(outputListHistos->FindObject(commonHistName.Data()));
339   if(commonHist) this->SetCommonHists(commonHist); 
340   
341   TString commonHistResName = "AliFlowCommonHistResultsFQD";
342   commonHistResName += fAnalysisLabel->Data();
343   AliFlowCommonHistResults *commonHistRes = dynamic_cast<AliFlowCommonHistResults*>
344                                             (outputListHistos->FindObject(commonHistResName.Data()));
345   if(commonHistRes) this->SetCommonHistsResults(commonHistRes); 
346   
347   // 2.) weights: 
348   TList *weightsList = dynamic_cast<TList*>(outputListHistos->FindObject("Weights"));
349   if(weightsList) this->SetWeightsList(weightsList);
350   Bool_t bUsePhiWeights = kFALSE;
351   Bool_t bUsePtWeights = kFALSE;
352   Bool_t bUseEtaWeights = kFALSE;
353   TString fUseParticleWeightsName = "fUseParticleWeightsFQD";
354   fUseParticleWeightsName += fAnalysisLabel->Data();
355   TProfile *useParticleWeights = dynamic_cast<TProfile*>(weightsList->FindObject(fUseParticleWeightsName.Data()));
356   if(useParticleWeights)
357   {
358    this->SetUseParticleWeights(useParticleWeights);  
359    bUsePhiWeights = (Int_t)useParticleWeights->GetBinContent(1);
360    bUsePtWeights = (Int_t)useParticleWeights->GetBinContent(2);
361    bUseEtaWeights = (Int_t)useParticleWeights->GetBinContent(3);
362   }
363   
364   // 3.) distributions and 4.) final results of fitting:
365   TString pWeightsFlag[2] = {"pWeights not used","pWeights used"};
366   TString sigmaFlag[2] = {"#sigma^{2} not fitted","#sigma^{2} fitted"};
367   
368   // q-distribution:
369   TString qDistributionName = "fqDistribution";
370   qDistributionName += fAnalysisLabel->Data();
371   // sum of particle weights:
372   TString sumOfParticleWeightsName = "fSumOfParticleWeightsName"; 
373   sumOfParticleWeightsName += fAnalysisLabel->Data();
374   // final results for integrated flow:
375   TString intFlowName = "fIntFlowFQD";
376   intFlowName += fAnalysisLabel->Data();
377   // sigma^2:
378   TString sigma2Name = "fSigma2";
379   sigma2Name += fAnalysisLabel->Data();
380   // chi^2:
381   TString chi2Name = "fChi2";
382   chi2Name += fAnalysisLabel->Data();
383   
384   TH1D *qDistribution[2] = {NULL};
385   TH1D *sumOfParticleWeights[2] = {NULL};
386   TH1D *intFlow[2][2] = {{NULL}};
387   TH1D *sigma2[2][2] = {{NULL}};
388   TH1D *chi2[2][2] = {{NULL}};
389    
390   for(Int_t pW=0;pW<1+(Int_t)(bUsePhiWeights||bUsePtWeights||bUseEtaWeights);pW++)
391   {
392    // q-distribution:
393    qDistribution[pW] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("%s, %s",qDistributionName.Data(),pWeightsFlag[pW].Data())));
394    if(qDistribution[pW])
395    {
396     this->SetqDistribution(qDistribution[pW],pW);
397    } else
398      {
399       cout<<"WARNING: qDistribution[pW] is NULL in AFAWFQD::GOH() !!!!"<<endl;
400       cout<<"pW = "<<pW<<endl;
401      }
402    // sum of particle weights:
403    sumOfParticleWeights[pW] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("%s, %s",sumOfParticleWeightsName.Data(),pWeightsFlag[pW].Data())));
404    if(sumOfParticleWeights[pW])
405    {
406     this->SetSumOfParticleWeights(sumOfParticleWeights[pW],pW);
407    } else
408      {
409       cout<<"WARNING: sumOfParticleWeights[pW] is NULL in AFAWFQD::GOH() !!!!"<<endl;
410       cout<<"pW = "<<pW<<endl;
411      }
412    // final results:
413    for(Int_t f=0;f<2;f++)
414    {
415     // final results for integrated flow:
416     intFlow[pW][f] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("%s, %s, %s",intFlowName.Data(),pWeightsFlag[pW].Data(),sigmaFlag[f].Data())));
417     if(intFlow[pW][f])
418     {
419      this->SetIntFlow(intFlow[pW][f],pW,f);
420     } else 
421       {
422        cout<<"WARNING: intFlow[pW][f] is NULL in AFAWFQD::GOH() !!!!"<<endl;
423        cout<<"pW = "<<pW<<endl;
424        cout<<"f  = "<<f<<endl;
425       }
426     // sigma^2:
427     sigma2[pW][f] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("%s, %s, %s",sigma2Name.Data(),pWeightsFlag[pW].Data(),sigmaFlag[f].Data())));
428     if(sigma2[pW][f])
429     {
430      this->SetSigma2(sigma2[pW][f],pW,f);
431     } else 
432       {
433        cout<<"WARNING: sigma2[pW][f] is NULL in AFAWFQD::GOH() !!!!"<<endl;
434        cout<<"pW = "<<pW<<endl;
435        cout<<"f  = "<<f<<endl;
436       } 
437     // chi^2:
438     chi2[pW][f] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("%s, %s, %s",chi2Name.Data(),pWeightsFlag[pW].Data(),sigmaFlag[f].Data())));
439     if(chi2[pW][f])
440     {
441      this->SetChi2(chi2[pW][f],pW,f);
442     } else 
443       {
444        cout<<"WARNING: chi2[pW][f] is NULL in AFAWFQD::GOH() !!!!"<<endl;
445        cout<<"pW = "<<pW<<endl;
446        cout<<"f  = "<<f<<endl;
447       }  
448       
449    } // end of for(Int_t f=0;f<2;f++)
450   } // end of for(Int_t pW=0;pW<1+(Int_t)(bUsePhiWeights||bUsePtWeights||bUseEtaWeights);pW++)
451   
452   // 5.) fitting parameters:
453   // q-distribution:
454   TString fittingParametersName = "fFittingParameters";
455   fittingParametersName += fAnalysisLabel->Data();
456   TProfile *fittingParameters = NULL;
457   fittingParameters = dynamic_cast<TProfile*>(outputListHistos->FindObject(fittingParametersName.Data()));
458   if(fittingParameters)
459   {
460    this->SetFittingParameters(fittingParameters);
461   } else
462     {
463      cout<<"WARNING:fittingParameters is NULL in AFAWFQD::GOH() !!!!"<<endl;
464     }
465   
466  } else // to if(outputListHistos)
467    {
468     cout<<"WARNING: outputListHistos is NULL in AFAWFQD::GOH() !!!!"<<endl;
469     exit(0);
470    } 
471     
472    
473 } // end of void AliFlowAnalysisWithFittingQDistribution::GetOutputHistograms(TList *outputListHistos) 
474
475
476 //================================================================================================================
477
478
479 void AliFlowAnalysisWithFittingQDistribution::WriteHistograms(TString* outputFileName)
480 {
481  //store the final results in output .root file
482  TFile *output = new TFile(outputFileName->Data(),"RECREATE");
483  //output->WriteObject(fHistList, "cobjFQD","SingleKey");
484  fHistList->SetName("cobjFQD");
485  fHistList->SetOwner(kTRUE);
486  fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
487  delete output;
488 }
489
490 //================================================================================================================
491
492 void AliFlowAnalysisWithFittingQDistribution::WriteHistograms(TString outputFileName)
493 {
494  //store the final results in output .root file
495  TFile *output = new TFile(outputFileName.Data(),"RECREATE");
496  //output->WriteObject(fHistList, "cobjFQD","SingleKey");
497  fHistList->SetName("cobjFQD");
498  fHistList->SetOwner(kTRUE);
499  fHistList->Write(fHistList->GetName(), TObject::kSingleKey); 
500  delete output;
501 }
502
503
504 //================================================================================================================
505
506
507 void AliFlowAnalysisWithFittingQDistribution::WriteHistograms(TDirectoryFile *outputFileName)
508 {
509  //store the final results in output .root file
510  fHistList->SetName("cobjFQD");
511  fHistList->SetOwner(kTRUE);
512  outputFileName->Add(fHistList);
513  outputFileName->Write(outputFileName->GetName(), TObject::kSingleKey);
514 }
515
516
517 //================================================================================================================
518
519
520 void AliFlowAnalysisWithFittingQDistribution::InitializeArrays()
521 {
522  // initialize all arrays
523  
524  for(Int_t pW=0;pW<2;pW++) // particle weights not used (0) or used (1)
525  {
526   fSumOfParticleWeights[pW] = NULL;
527   fqDistribution[pW] = NULL; 
528   for(Int_t f=0;f<2;f++) // sigma^2 not fitted (0) or fitted (1)
529   {
530    fIntFlow[pW][f] = NULL;
531    fSigma2[pW][f] = NULL;
532    fChi2[pW][f] = NULL;
533   }
534  } 
535
536 } // end of void AliFlowAnalysisWithFittingQDistribution::InitializeArrays()
537
538
539 //================================================================================================================
540
541
542 void AliFlowAnalysisWithFittingQDistribution::BookCommonHistograms()
543 {
544  // book common histograms
545  
546  // common control histogram: 
547  TString commonHistName = "AliFlowCommonHistFQD";
548  commonHistName += fAnalysisLabel->Data();
549  fCommonHists = new AliFlowCommonHist(commonHistName.Data());
550  fHistList->Add(fCommonHists);  
551
552  // common histograms for final results:
553  TString commonHistResName = "AliFlowCommonHistResultsFQD";
554  commonHistResName += fAnalysisLabel->Data();
555  fCommonHistsResults = new AliFlowCommonHistResults(commonHistResName.Data());
556  fHistList->Add(fCommonHistsResults); 
557
558 } // end of void AliFlowAnalysisWithFittingQDistribution::BookCommonHistograms()
559
560
561 //================================================================================================================
562
563
564 void AliFlowAnalysisWithFittingQDistribution::BookAndFillWeightsHistograms()
565 {
566  // book and fill histograms which hold phi, pt and eta weights
567
568  if(!fWeightsList)
569  {
570   cout<<"WARNING: fWeightsList is NULL in AFAWFQD::BAFWH() !!!!"<<endl;
571   exit(0);  
572  }
573     
574  TString fUseParticleWeightsName = "fUseParticleWeightsFQD";
575  fUseParticleWeightsName += fAnalysisLabel->Data();
576  fUseParticleWeights = new TProfile(fUseParticleWeightsName.Data(),"0 = particle weight not used, 1 = particle weight used ",3,0,3);
577  fUseParticleWeights->SetLabelSize(0.08);
578  (fUseParticleWeights->GetXaxis())->SetBinLabel(1,"w_{#phi}");
579  (fUseParticleWeights->GetXaxis())->SetBinLabel(2,"w_{p_{T}}");
580  (fUseParticleWeights->GetXaxis())->SetBinLabel(3,"w_{#eta}");
581  fUseParticleWeights->Fill(0.5,(Int_t)fUsePhiWeights);
582  fUseParticleWeights->Fill(1.5,(Int_t)fUsePtWeights);
583  fUseParticleWeights->Fill(2.5,(Int_t)fUseEtaWeights);
584  fWeightsList->Add(fUseParticleWeights); 
585   
586  if(fUsePhiWeights)
587  {
588   if(fWeightsList->FindObject("phi_weights"))
589   {
590    fPhiWeights = dynamic_cast<TH1F*>(fWeightsList->FindObject("phi_weights"));
591    if(fPhiWeights->GetBinWidth(1) != fPhiBinWidth)
592    {
593     cout<<"WARNING: fPhiWeights->GetBinWidth(1) != fPhiBinWidth in AFAWFQD::BAFWH() !!!!        "<<endl;
594     cout<<"         This indicates inconsistent binning in phi histograms throughout the code."<<endl;
595     exit(0);
596    }
597   } else 
598     {
599      cout<<"WARNING: fWeightsList->FindObject(\"phi_weights\") is NULL in AFAWFQD::BAFWH() !!!!"<<endl;
600      exit(0);
601     }
602  } // end of if(fUsePhiWeights)
603  
604  if(fUsePtWeights) 
605  {
606   if(fWeightsList->FindObject("pt_weights"))
607   {
608    fPtWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("pt_weights"));
609    if(fPtWeights->GetBinWidth(1) != fPtBinWidth)
610    {
611     cout<<"WARNING: fPtWeights->GetBinWidth(1) != fPtBinWidth in AFAWFQD::BAFWH() !!!!         "<<endl;
612     cout<<"         This indicates insconsistent binning in pt histograms throughout the code."<<endl;
613     exit(0);
614    }
615   } else 
616     {
617      cout<<"WARNING: fWeightsList->FindObject(\"pt_weights\") is NULL in AFAWFQD::BAFWH() !!!!"<<endl;
618      exit(0);
619     }
620  } // end of if(fUsePtWeights)    
621
622  if(fUseEtaWeights) 
623  {
624   if(fWeightsList->FindObject("eta_weights"))
625   {
626    fEtaWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("eta_weights"));
627    if(fEtaWeights->GetBinWidth(1) != fEtaBinWidth)
628    {
629     cout<<"WARNING: fEtaWeights->GetBinWidth(1) != fEtaBinWidth in AFAWFQD::BAFWH() !!!!        "<<endl;
630     cout<<"         This indicates insconsistent binning in eta histograms throughout the code."<<endl;
631     exit(0);
632    }
633   } else 
634     {
635      cout<<"WARNING: fUseEtaWeights && fWeightsList->FindObject(\"eta_weights\") is NULL in AFAWFQD::BAFWH() !!!!"<<endl;
636      exit(0);
637     }
638  } // end of if(fUseEtaWeights)
639  
640 } // end of AliFlowAnalysisWithFittingQDistribution::BookAndFillWeightsHistograms()
641
642
643 //================================================================================================================================
644
645
646 void AliFlowAnalysisWithFittingQDistribution::AccessConstants()
647 {
648  // access needed common constants from AliFlowCommonConstants
649  
650  fnBinsPhi = AliFlowCommonConstants::GetMaster()->GetNbinsPhi();
651  fPhiMin = AliFlowCommonConstants::GetMaster()->GetPhiMin();         
652  fPhiMax = AliFlowCommonConstants::GetMaster()->GetPhiMax();
653  if(fnBinsPhi) fPhiBinWidth = (fPhiMax-fPhiMin)/fnBinsPhi;  
654  fnBinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
655  fPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();           
656  fPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
657  if(fnBinsPt) fPtBinWidth = (fPtMax-fPtMin)/fnBinsPt;  
658  fnBinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
659  fEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();         
660  fEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
661  if(fnBinsEta) fEtaBinWidth = (fEtaMax-fEtaMin)/fnBinsEta;  
662  
663 } // end of void AliFlowAnalysisWithFittingQDistribution::AccessConstants()
664
665
666 //================================================================================================================================
667
668
669 void AliFlowAnalysisWithFittingQDistribution::BookEverythingForDistributions()
670 {
671  // book histograms for distributions
672  
673  TString pWeightsFlag[2] = {"pWeights not used","pWeights used"};
674  TString sigmaFlag[2] = {"#sigma^{2} not fitted","#sigma^{2} fitted"};
675  // q-distribution:
676  TString fqDistributionName = "fqDistribution";
677  fqDistributionName += fAnalysisLabel->Data();
678  // sum of particle weights:
679  TString fSumOfParticleWeightsName = "fSumOfParticleWeightsName";
680  fSumOfParticleWeightsName += fAnalysisLabel->Data();
681  // final results for integrated flow:
682  TString fIntFlowName = "fIntFlowFQD";
683  fIntFlowName += fAnalysisLabel->Data();
684  // sigma^2:
685  TString fSigma2Name = "fSigma2";
686  fSigma2Name += fAnalysisLabel->Data();
687  // chi^2:
688  TString fChi2Name = "fChi2";
689  fChi2Name += fAnalysisLabel->Data();
690  
691  for(Int_t pW=0;pW<1+(Int_t)(fUsePhiWeights||fUsePtWeights||fUseEtaWeights);pW++) // particle weights not used (0) or used (1)
692  {
693   // q-distribution:
694   fqDistribution[pW] = new TH1D(Form("%s, %s",fqDistributionName.Data(),pWeightsFlag[pW].Data()),"q-distribution",10000,0,1000);
695   fqDistribution[pW]->SetXTitle("q_{n}=Q_{n}/#sqrt{M}");
696   fqDistribution[pW]->SetYTitle("Counts");
697   fHistList->Add(fqDistribution[pW]);
698   // sum of particle weights: 
699   fSumOfParticleWeights[pW] = new TH1D(Form("%s, %s",fSumOfParticleWeightsName.Data(),pWeightsFlag[pW].Data()),"Sum of particle weights",10000,0,10000);
700   fSumOfParticleWeights[pW]->SetXTitle("#sum_{i=1}^{N} w_{i}");
701   fSumOfParticleWeights[pW]->SetYTitle("Counts");
702   fHistList->Add(fSumOfParticleWeights[pW]);
703   
704   for(Int_t f=0;f<2;f++) // sigma^2 not fitted (0) or fitted (1)
705   {
706    // final results for integrated flow:
707    fIntFlow[pW][f] = new TH1D(Form("%s, %s, %s",fIntFlowName.Data(),pWeightsFlag[pW].Data(),sigmaFlag[f].Data()),"Integrated Flow",1,0,1);
708    fIntFlow[pW][f]->SetLabelSize(0.08);
709    (fIntFlow[pW][f]->GetXaxis())->SetBinLabel(1,"v_{n}");
710    fHistList->Add(fIntFlow[pW][f]);
711    // sigma^2:
712    fSigma2[pW][f] = new TH1D(Form("%s, %s, %s",fSigma2Name.Data(),pWeightsFlag[pW].Data(),sigmaFlag[f].Data()),"#sigma^{2}",1,0,1);
713    fSigma2[pW][f]->SetLabelSize(0.08);
714    (fSigma2[pW][f]->GetXaxis())->SetBinLabel(1,"#sigma^{2}");
715    fHistList->Add(fSigma2[pW][f]);
716    // chi^2:
717    fChi2[pW][f] = new TH1D(Form("%s, %s, %s",fChi2Name.Data(),pWeightsFlag[pW].Data(),sigmaFlag[f].Data()),"#chi^{2} (Minuit)",1,0,1);
718    fChi2[pW][f]->SetLabelSize(0.08);
719    (fChi2[pW][f]->GetXaxis())->SetLabelOffset(0.01);
720    (fChi2[pW][f]->GetXaxis())->SetBinLabel(1,"#chi^{2}");
721    fHistList->Add(fChi2[pW][f]);
722   } // end of for(Int_t f=0;f<2;f++) // sigma^2 not fitted or fitted
723   
724  } // end of for(Int_t pW=0;pW<1+(Int_t)(fUsePhiWeights||fUsePtWeights||fUseEtaWeights);pW++) // particle weights not used (0) or used (1)
725  
726  // book profile fFittingParameters which will hold all fitting parameters:
727  TString fFittingParametersName = "fFittingParameters";
728  fFittingParametersName += fAnalysisLabel->Data(); 
729  fFittingParameters = new TProfile(fFittingParametersName.Data(),"Parameters for fitting q-distribution",8,0,8);
730  fFittingParameters->SetLabelSize(0.05);
731  (fFittingParameters->GetXaxis())->SetBinLabel(1,"treshold");
732  (fFittingParameters->GetXaxis())->SetBinLabel(2,"starting v_{n}");
733  (fFittingParameters->GetXaxis())->SetBinLabel(3,"min. v_{n}");
734  (fFittingParameters->GetXaxis())->SetBinLabel(4,"max. v_{n}");
735  (fFittingParameters->GetXaxis())->SetBinLabel(5,"starting #sigma^{2}");
736  (fFittingParameters->GetXaxis())->SetBinLabel(6,"min. #sigma^{2}");
737  (fFittingParameters->GetXaxis())->SetBinLabel(7,"max. #sigma^{2}");
738  (fFittingParameters->GetXaxis())->SetBinLabel(8,"plot or not?");
739  fHistList->Add(fFittingParameters);
740  
741 } // end of void AliFlowAnalysisWithFittingQDistribution::BookEverythingForDistributions()
742
743
744 //================================================================================================================================
745
746
747 void AliFlowAnalysisWithFittingQDistribution::DoFit(Bool_t useParticleWeights, Bool_t sigma2Fitted)
748 {
749  // do final fit for q-distribution
750  
751  // shortcuts for flags:
752  Int_t pW = (Int_t)(useParticleWeights);
753  Int_t s2F = (Int_t)(sigma2Fitted);
754  
755  if(!(fqDistribution[pW] && fSumOfParticleWeights[pW] && fIntFlow[pW][s2F] && fSigma2[pW][s2F] && fChi2[pW][s2F])) 
756  { 
757   cout<<"WARNING: fqDistribution[pW] && fSumOfParticleWeights[pW] && fIntFlow[pW][s2F] && fSigma2[pW][s2F] && fChi2[pW][s2F] is NULL in AFAWFQD::DoFit() !!!!"<<endl;
758   cout<<"pW  = "<<pW<<endl;
759   cout<<"s2F = "<<s2F<<endl;
760   exit(0);
761  }
762  
763  // average multiplicity and number of events:
764  Double_t AvM = fSumOfParticleWeights[pW]->GetMean(1);
765  //Int_t nEvts = (Int_t)fSumOfParticleWeights[pW]->GetEntries();
766  
767  // start fitting from the bin with at least fTreshold entries, 
768  // finish fitting at the bin with at least fTreshold entries:
769  Int_t binMin = fqDistribution[pW]->FindFirstBinAbove(fTreshold);  
770  Int_t binMax = fqDistribution[pW]->FindLastBinAbove(fTreshold);
771  Double_t binWidth = fqDistribution[pW]->GetBinWidth(4); // assuming that all bins have the same width 
772  if(binWidth == 0) 
773  {
774   cout<<"WARNING: binWidth == 0 in AFAWFQD::DoFit()"<<endl;
775   exit(0);
776  }
777  Double_t qmin = (binMin-1)*binWidth; 
778  Double_t qmax = (binMax)*binWidth;
779  Double_t ent = 0.; // number of entries between binMin and binMax:
780  for(Int_t b=binMin;b<=binMax;b++)
781  {
782   ent += fqDistribution[pW]->GetBinContent(b);
783  }
784  Double_t norm = binWidth*ent; // norm (assuming that all bins have the same width)
785
786  // fitting function:
787  TF1 *fittingFun = new TF1("fittingFun","[2]*(x/[1])*exp(-(x*x+[0]*[0])/(2.*[1]))*TMath::BesselI0(x*[0]/[1])",qmin,qmax); 
788  
789  fittingFun->SetParNames("v*sqrt{sum of particle weights}","sigma^2","norm");
790  fittingFun->SetParameters(fvStart*pow(AvM,0.5),fSigma2Start,norm);         
791  fittingFun->SetParLimits(0,fvMin*pow(AvM,0.5),fvMax*pow(AvM,0.5)); 
792  
793  if(s2F == 0)
794  {
795   fittingFun->FixParameter(1,0.5);
796  } else
797    {
798     fittingFun->SetParLimits(1,fSigma2Min,fSigma2Max);          
799    }
800  fittingFun->FixParameter(2,norm);  
801
802  // fitting (do it only if # of entries >50): // to be improved (this is only a pragmatics fix to avoid TMinuit crash)
803  if(ent > 50)
804  {
805   fqDistribution[pW]->Fit("fittingFun","NQ","",qmin,qmax);
806  }
807  // results:
808  Double_t v = 0.; // integrated flow
809  Double_t vError = 0.; // error of integrated flow 
810  Double_t sigma2 = 0.; // sigma^2
811  Double_t sigma2Error = 0.; // error of sigma^2
812  Double_t chi2 = 0; // chi^2 from Minuit
813  
814  if(AvM)
815  { 
816   // integrated flow:
817   v = fittingFun->GetParameter(0)/pow(AvM,0.5);
818   vError = fittingFun->GetParError(0)/pow(AvM,0.5);
819   fIntFlow[pW][s2F]->SetBinContent(1,v); // s2F is shortcut for "sigma^2 fitted"
820   fIntFlow[pW][s2F]->SetBinError(1,vError); // s2F is shortcut for "sigma^2 fitted"
821  }
822  
823  if(s2F == 0) // sigma^2 not fitted, but fixed to 0.5
824  {
825   // sigma^2:
826   sigma2 = 0.5;
827   fSigma2[pW][0]->SetBinContent(1,sigma2);  
828   fSigma2[pW][0]->SetBinError(1,0.);
829   // chi^2:
830   chi2 = fittingFun->GetChisquare();
831   fChi2[pW][0]->SetBinContent(1,chi2);  
832   //fChi2[pW][0]->SetBinError(1,0.);  
833  } else // sigma^2 fitted
834    {
835     // sigma^2:
836     sigma2 = fittingFun->GetParameter(1);
837     sigma2Error = fittingFun->GetParError(1);
838     fSigma2[pW][1]->SetBinContent(1,sigma2);  
839     fSigma2[pW][1]->SetBinError(1,sigma2Error);    
840     // chi^2:
841     chi2 = fittingFun->GetChisquare();
842     fChi2[pW][1]->SetBinContent(1,chi2);  
843     //fChi2[pW][1]->SetBinError(1,0.);  
844    }
845  
846  if(fPlotResults && !(fUsePhiWeights||fUsePtWeights||fUseEtaWeights)) // to be improved (plot also the plot when particle weights are used)
847  {
848   // set ranges: // to be improved (there is certainly a better way to implement this)
849   Int_t firstNonEmptyBin = fqDistribution[pW]->FindFirstBinAbove(0);
850   Double_t lowRange = fqDistribution[pW]->GetBinLowEdge(firstNonEmptyBin);
851   Int_t lastNonEmptyBin = fqDistribution[pW]->FindLastBinAbove(0);
852   Double_t upperRange = fqDistribution[pW]->GetBinLowEdge(lastNonEmptyBin+10);
853   (fqDistribution[pW]->GetXaxis())->SetRangeUser(lowRange,upperRange);
854   
855   if(s2F == 0)
856   { 
857    // to be improved (there is certainly a better way to implement this)
858    fqDistribution[pW]->SetFillColor(16);  
859    fqDistribution[pW]->SetTitle("Fitted q-distribution");
860    fqDistribution[pW]->Draw("");
861    fLegend->AddEntry(fqDistribution[pW],"q-distribution","f");
862    TF1 *fittingFunTemp = (TF1*)fittingFun->Clone("fittingFunTemp");
863    fittingFunTemp->SetLineColor(4); // 4 = blue color
864    fittingFunTemp->Draw("SAME"); 
865    fLegend->AddEntry("fittingFunTemp","#sigma^{2} fixed","l");
866    fLegend->Draw("SAME");      
867   } else
868     {    
869      fittingFun->SetLineColor(2); // 2 = red color   
870      fittingFun->Draw("SAME");
871      fLegend->AddEntry("fittingFun","#sigma^{2} fitted","l");    
872     } 
873  } // end of if(fPlotResults)
874
875 } // end of void AliFlowAnalysisWithFittingQDistribution::DoFit(Bool_t useParticleWeights)
876
877
878 //================================================================================================================================ 
879
880
881 void AliFlowAnalysisWithFittingQDistribution::FillCommonHistResultsIntFlow(Bool_t useParticleWeights, Bool_t sigma2Fitted)
882 {
883  // fill in AliFlowCommonHistResults histograms relevant for 'NONAME' integrated flow (to be improved (name))
884  
885  // shortcuts for the flags:
886  Int_t pW = (Int_t)(useParticleWeights); // particle weights not used (0) or used (1)
887  Int_t s2F = (Int_t)(sigma2Fitted); // 0 = sigma^2 not fitted (but fixed to 0.5), 1 = sigma^2 fitted
888  
889  if(!fIntFlow[pW][s2F])
890  {
891   cout<<"WARNING: fIntFlow[pW][s2F] is NULL in AFAWFQD::FCHRIF() !!!!"<<endl;
892   cout<<"pW  = "<<pW<<endl;
893   cout<<"s2F = "<<s2F<<endl;
894   exit(0); 
895  }  
896  
897  if(!fSumOfParticleWeights[pW])
898  {
899   cout<<"WARNING: fSumOfParticleWeights[pW] is NULL in AFAWFQD::FCHRIF() !!!!"<<endl;
900   cout<<"pW = "<<pW<<endl;
901   exit(0);
902  }
903  
904  if(!(fCommonHistsResults))
905  {
906   cout<<"WARNING: fCommonHistsResults is NULL in AFAWFQD::FCHRIF() !!!!"<<endl; 
907   exit(0);
908  }
909   
910  // fill integrated flow:
911  Double_t v = fIntFlow[pW][s2F]->GetBinContent(1); 
912  Double_t vError = fIntFlow[pW][s2F]->GetBinError(1);
913  
914  fCommonHistsResults->FillIntegratedFlow(v,vError);   
915  
916  // fill chi (this chi stands for resolution, not to be confused with chi2 used before):
917  Double_t AvM = fSumOfParticleWeights[pW]->GetMean(1);
918  Double_t chi = AvM*pow(v,2); 
919  if(chi>=0)
920  {
921   fCommonHistsResults->FillChi(pow(chi,0.5));   
922   fCommonHistsResults->FillChiRP(pow(chi,0.5));   
923  }
924    
925 } // end of void AliFlowAnalysisWithFittingQDistribution::FillCommonHistResultsIntFlow(Bool_t useParticleWeights, Bool_t sigma2NotFixed) 
926
927
928 //================================================================================================================================ 
929
930
931 void AliFlowAnalysisWithFittingQDistribution::PrintFinalResultsForIntegratedFlow()
932 {
933  // print the final results for integrated flow on the screen
934  
935  // shortcuts: pW  = particle weights 
936  //            s2F = sigma^2 fitted 
937  
938  for(Int_t pW=0;pW<1+(Int_t)(fUsePhiWeights||fUsePtWeights||fUseEtaWeights);pW++)
939  {
940   if(!fSumOfParticleWeights[pW])
941   {
942    cout<<"WARNING: fSumOfParticleWeights[pW] is NULL in AFAWFQD::FCHRIF() !!!!"<<endl;
943    cout<<"pW = "<<pW<<endl;
944    exit(0);
945   }
946   for(Int_t s2F=0;s2F<2;s2F++)
947   {
948    if(!fIntFlow[pW][s2F])
949    {
950     cout<<"WARNING: fIntFlow[pW][s2F] is NULL in AFAWFQD::FCHRIF() !!!!"<<endl;
951     cout<<"pW  = "<<pW<<endl;
952     cout<<"s2F = "<<s2F<<endl;
953     exit(0); 
954    }
955   }  
956  }  
957  
958  if(!(fCommonHistsResults))
959  {
960   cout<<"WARNING: fCommonHistsResults is NULL in AFAWFQD::FCHRIF() !!!!"<<endl; 
961   exit(0);
962  }
963  
964  // shortcut for the harmonic:
965  Int_t n = 0;
966  if(fCommonHists->GetHarmonic())
967  {
968   n = (Int_t)(fCommonHists->GetHarmonic())->GetBinContent(1); 
969  }
970  
971  // printing:
972  cout<<" "<<endl;
973  cout<<"***************************************"<<endl;
974  cout<<"***************************************"<<endl;
975  cout<<"      integrated flow by fitting "<<endl;
976  cout<<"           q-distribution:      "<<endl;
977  if(fUsePhiWeights||fUsePtWeights||fUseEtaWeights)
978  {
979   cout<<"           (with weights)       "<<endl;
980  } else
981    {
982     cout<<"          (without weights)       "<<endl;
983    }   
984  cout<<endl;
985
986  if(fUsePhiWeights||fUsePtWeights||fUseEtaWeights)
987  {
988   cout<<"1.) sigma^2 not fitted: "<<endl;
989   cout<<"  v_"<<n<<"{FQD} = "<<fIntFlow[1][0]->GetBinContent(1)<<" +/- "<<fIntFlow[1][0]->GetBinError(1)<<endl;
990   cout<<"  sigma^2 = 0.5 +/- 0 "<<endl; 
991   cout<<"  chi^2 = "<<fChi2[1][0]->GetBinContent(1)<<" (Minuit)"<<endl; 
992   cout<<" "<<endl;   
993   cout<<"2.) sigma^2 fitted: "<<endl;
994   cout<<"  v_"<<n<<"{FQD} = "<<fIntFlow[1][1]->GetBinContent(1)<<" +/- "<<fIntFlow[1][1]->GetBinError(1)<<endl;
995   cout<<"  sigma^2 = "<<fSigma2[1][1]->GetBinContent(1)<<" +/- "<<fSigma2[1][1]->GetBinError(1)<<endl; 
996   cout<<"  chi^2 = "<<fChi2[1][1]->GetBinContent(1)<<" (Minuit)"<<endl; 
997   cout<<" "<<endl; 
998   cout<<"      nEvts = "<<fSumOfParticleWeights[1]->GetEntries()<<", AvM = "<<fSumOfParticleWeights[1]->GetMean()<<endl;
999   cout<<" "<<endl;
1000  } else
1001    { 
1002     cout<<"1.) sigma^2 not fitted: "<<endl;
1003     cout<<endl;
1004     cout<<"  v_"<<n<<"{FQD} = "<<fIntFlow[0][0]->GetBinContent(1)<<" +/- "<<fIntFlow[0][0]->GetBinError(1)<<endl;
1005     cout<<"  sigma^2 = 0.5 +/- 0 "<<endl; 
1006     cout<<"  chi^2 = "<<fChi2[0][0]->GetBinContent(1)<<" (Minuit)"<<endl;  
1007     cout<<" "<<endl;   
1008     cout<<"2.) sigma^2 fitted: "<<endl;
1009     cout<<endl;
1010     cout<<"  v_"<<n<<"{FQD} = "<<fIntFlow[0][1]->GetBinContent(1)<<" +/- "<<fIntFlow[0][1]->GetBinError(1)<<endl;
1011     cout<<"  sigma^2 = "<<fSigma2[0][1]->GetBinContent(1)<<" +/- "<<fSigma2[0][1]->GetBinError(1)<<endl; 
1012     cout<<"  chi^2 = "<<fChi2[0][1]->GetBinContent(1)<<" (Minuit)"<<endl; 
1013     cout<<" "<<endl;  
1014     cout<<"      nEvts = "<<fSumOfParticleWeights[0]->GetEntries()<<", AvM = "<<fSumOfParticleWeights[0]->GetMean()<<endl;
1015     cout<<" "<<endl;
1016    }
1017     
1018  cout<<"***************************************"<<endl;
1019  cout<<"***************************************"<<endl; 
1020  cout<<endl;
1021   
1022 } // end of void AliFlowAnalysisWithFittingQDistribution::PrintFinalResultsForIntegratedFlow()
1023
1024
1025 //================================================================================================================================ 
1026
1027
1028 void AliFlowAnalysisWithFittingQDistribution::StoreFittingParameters()
1029 {
1030  // store fitting parameters in profile fFittingParameters
1031  
1032  // Binning of fFittingParameters is organized as follows:
1033  // 1st bin: fTreshold
1034  // 2nd bin: fvStart
1035  // 3rd bin: fvMin
1036  // 4th bin: fvMax
1037  // 5th bin: fSigma2Start
1038  // 6th bin: fSigma2Min
1039  // 7th bin: fSigma2Max
1040  // 8th bin: fPlotResults
1041  
1042  if(!fFittingParameters)
1043  {
1044   cout<<"WARNING: fFittingParameters is NULL in AFAWFQD::SFP() !!!!"<<endl;
1045   exit(0);
1046  }
1047  
1048  fFittingParameters->Reset();
1049  fFittingParameters->Fill(0.5,fTreshold);
1050  fFittingParameters->Fill(1.5,fvStart);
1051  fFittingParameters->Fill(2.5,fvMin);
1052  fFittingParameters->Fill(3.5,fvMax);
1053  fFittingParameters->Fill(4.5,fSigma2Start);
1054  fFittingParameters->Fill(5.5,fSigma2Min);
1055  fFittingParameters->Fill(6.5,fSigma2Max);
1056  fFittingParameters->Fill(7.5,(Int_t)fPlotResults);
1057  
1058 } // end of void AliFlowAnalysisWithFittingQDistribution::StoreFittingParameters()
1059
1060
1061 //================================================================================================================================ 
1062
1063
1064 void AliFlowAnalysisWithFittingQDistribution::AccessFittingParameters()
1065 {
1066  // access fitting parameters:
1067  
1068  if(!fFittingParameters)
1069  {
1070   cout<<"WARNING: fFittingParameters is NULL in AFAWFQD::AFP() !!!!"<<endl;
1071   exit(0);
1072  }
1073  
1074  fTreshold = fFittingParameters->GetBinContent(1);
1075  fvStart = fFittingParameters->GetBinContent(2);
1076  fvMin = fFittingParameters->GetBinContent(3);
1077  fvMax = fFittingParameters->GetBinContent(4);
1078  fSigma2Start = fFittingParameters->GetBinContent(5);
1079  fSigma2Min = fFittingParameters->GetBinContent(6);
1080  fSigma2Max = fFittingParameters->GetBinContent(7);
1081  fPlotResults = (Bool_t) fFittingParameters->GetBinContent(8);
1082  
1083 } // end of void AliFlowAnalysisWithFittingQDistribution::AccessFittingParameters()