1 /*************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /********************************
17 * integrated flow estimate by *
18 * fitting q-distribution *
20 * author: Ante Bilandzic *
23 * based on the macro written *
24 * by Sergei Voloshin *
25 *******************************/
27 #define AliFittingFunctionsForQDistribution_cxx
29 #include "Riostream.h"
33 #include "TParticle.h"
42 #include "AliFlowEventSimple.h"
43 #include "AliFlowTrackSimple.h"
44 #include "AliFittingQDistribution.h"
45 #include "AliFlowCommonConstants.h"
46 #include "AliFlowCommonHistResults.h"
47 #include "AliFittingFunctionsForQDistribution.h"
49 ClassImp(AliFittingFunctionsForQDistribution)
51 //================================================================================================================_
53 AliFittingFunctionsForQDistribution::AliFittingFunctionsForQDistribution():
55 fQDistributionFQD(NULL),
63 AliFittingFunctionsForQDistribution::~AliFittingFunctionsForQDistribution()
68 AliFittingFunctionsForQDistribution::AliFittingFunctionsForQDistribution(TProfile *AvMult, TH1D *QDistribution, TH1D *intFlowRes, TH1D *sigma2, AliFlowCommonHistResults *chr):
70 fQDistributionFQD(QDistribution),
71 fIntFlowResFQD(intFlowRes),
78 //================================================================================================================
80 void AliFittingFunctionsForQDistribution::Calculate()
82 //fitting q-distribution
84 Int_t n=2;//harmonic (to be improved)
86 Double_t AvM = fAvMultFQD->GetBinContent(1);
88 Int_t nEvts = (Int_t)(fAvMultFQD->GetBinEntries(1));
90 Double_t qmin=(fQDistributionFQD->GetXaxis())->GetXmin();
91 Double_t qmax=(fQDistributionFQD->GetXaxis())->GetXmax();
92 Double_t bin=fQDistributionFQD->GetBinWidth(4);//assuming that all bins have the same width
93 Double_t ent=fQDistributionFQD->GetEntries();
94 Double_t norm=bin*ent;//assuming that all bins have the same width
96 TF1 *fittingFun = new TF1("fittingFun","[2]*(x/[1])*exp(-(x*x+[0]*[0])/(2.*[1]))*TMath::BesselI0(x*[0]/[1])",qmin,qmax);
98 fittingFun->SetParNames("V","sigma","norm");
99 fittingFun->SetParameters(0.1*pow(AvM,0.5),0.5,norm);
101 fittingFun->SetParLimits(0,0.,10.);//to be improved (limits)
102 fittingFun->SetParLimits(1,0.,5.5);//to be improved (limits)
103 fittingFun->FixParameter(2,norm);
105 Double_t v=0.,errorv=0.,sigma2=0.,errorsigma2=0.;
107 if(fQDistributionFQD->GetEntries()>50)//to be improved (only a pragmatic fix)
109 fQDistributionFQD->Fit("fittingFun","NQ","",qmin,qmax);
112 v = fittingFun->GetParameter(0)/pow(AvM,0.5);
113 errorv = fittingFun->GetParError(0)/pow(AvM,0.5);
115 sigma2 = fittingFun->GetParameter(1);
116 errorsigma2 = fittingFun->GetParError(1);
120 cout<<"***************************************"<<endl;
121 cout<<"***************************************"<<endl;
122 cout<<" integrated flow by fitting "<<endl;
123 cout<<" q-distribution: "<<endl;
125 cout<<" v_"<<n<<"{FQD} = "<<v<<" +/- "<<errorv<<endl;
126 cout<<" sigma^2 = "<<sigma2<<" +/- "<<errorsigma2<<endl;
127 //cout<<"vm = "<<v*pow(AvM,0.5)<<endl;
129 cout<<" nEvts = "<<nEvts<<", AvM = "<<AvM<<endl;
130 cout<<"***************************************"<<endl;
131 cout<<"***************************************"<<endl;
132 fIntFlowResFQD->SetBinContent(1,v);
133 fIntFlowResFQD->SetBinError(1,errorv);
134 fSigma2->SetBinContent(1,sigma2);
135 fSigma2->SetBinError(1,errorsigma2);
138 fchrFQD->FillIntegratedFlow(v,errorv);
139 fchrFQD->FillChi(v*pow(AvM,0.5));
141 fchrFQD->FillIntegratedFlowRP(v,errorv);
142 fchrFQD->FillChiRP(v*pow(AvM,0.5));
146 //fQDistributionFQD->Draw("");
147 //fittingFun->Draw("SAME");
149 }//end of Calculate()