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 #define AliFlowAnalysisWithScalarProduct_cxx
18 #include "Riostream.h"
27 #include "AliFlowCommonConstants.h"
28 #include "AliFlowEventSimple.h"
29 #include "AliFlowTrackSimple.h"
30 #include "AliFlowCommonHist.h"
31 #include "AliFlowCommonHistResults.h"
32 #include "AliFlowAnalysisWithScalarProduct.h"
33 #include "AliFlowVector.h"
35 //////////////////////////////////////////////////////////////////////////////
36 // AliFlowAnalysisWithScalarProduct:
38 // Maker to analyze Flow with the Scalar product method.
40 // authors: N. van der Kolk (kolk@nikhef.nl), A. Bilandzic (anteb@nikhef.nl)
41 //////////////////////////////////////////////////////////////////////////////
43 ClassImp(AliFlowAnalysisWithScalarProduct)
45 //-----------------------------------------------------------------------
46 AliFlowAnalysisWithScalarProduct::AliFlowAnalysisWithScalarProduct():
51 fUsePhiWeights(kFALSE),
52 fPhiWeightsSub0(NULL),
53 fPhiWeightsSub1(NULL),
56 fHistProUQetaRP(NULL),
57 fHistProUQetaPOI(NULL),
59 fHistProUQPtPOI(NULL),
61 fHistProQaQbNorm(NULL),
62 fApplyCorrectionForNUA(kFALSE),
63 fHistProQaQbReImNorm(NULL),
64 fHistProNonIsotropicTermsQ(NULL),
65 fHistSumOfLinearWeights(NULL),
66 fHistSumOfQuadraticWeights(NULL),
67 fHistProUQQaQbPtRP(NULL),
68 fHistProUQQaQbEtaRP(NULL),
69 fHistProUQQaQbPtPOI(NULL),
70 fHistProUQQaQbEtaPOI(NULL),
72 fCommonHistsRes(NULL),
76 fHistResolution(NULL),
78 fHistQaNormvsMa(NULL),
80 fHistQbNormvsMb(NULL),
85 fWeightsList = new TList();
86 fHistList = new TList();
89 for(Int_t i=0;i<3;i++)
91 fHistSumOfWeightsPtRP[i] = NULL;
92 fHistSumOfWeightsEtaRP[i] = NULL;
93 fHistSumOfWeightsPtPOI[i] = NULL;
94 fHistSumOfWeightsEtaPOI[i] = NULL;
96 for(Int_t rp=0;rp<2;rp++)
98 for(Int_t pe=0;pe<2;pe++)
100 for(Int_t sc=0;sc<2;sc++)
102 fHistProNonIsotropicTermsU[rp][pe][sc] = NULL;
107 //-----------------------------------------------------------------------
108 AliFlowAnalysisWithScalarProduct::~AliFlowAnalysisWithScalarProduct()
115 //-----------------------------------------------------------------------
116 void AliFlowAnalysisWithScalarProduct::WriteHistograms(TString* outputFileName)
118 //store the final results in output .root file
120 TFile *output = new TFile(outputFileName->Data(),"RECREATE");
121 //output->WriteObject(fHistList, "cobjSP","SingleKey");
122 fHistList->SetName("cobjSP");
123 fHistList->SetOwner(kTRUE);
124 fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
128 //-----------------------------------------------------------------------
129 void AliFlowAnalysisWithScalarProduct::WriteHistograms(TString outputFileName)
131 //store the final results in output .root file
133 TFile *output = new TFile(outputFileName.Data(),"RECREATE");
134 //output->WriteObject(fHistList, "cobjSP","SingleKey");
135 fHistList->SetName("cobjSP");
136 fHistList->SetOwner(kTRUE);
137 fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
141 //-----------------------------------------------------------------------
142 void AliFlowAnalysisWithScalarProduct::WriteHistograms(TDirectoryFile *outputFileName)
144 //store the final results in output .root file
145 fHistList->SetName("cobjSP");
146 fHistList->SetOwner(kTRUE);
147 outputFileName->Add(fHistList);
148 outputFileName->Write(outputFileName->GetName(), TObject::kSingleKey);
151 //-----------------------------------------------------------------------
152 void AliFlowAnalysisWithScalarProduct::Init() {
154 //Define all histograms
155 cout<<"---Analysis with the Scalar Product Method--- Init"<<endl;
157 //save old value and prevent histograms from being added to directory
158 //to avoid name clashes in case multiple analaysis objects are used
160 Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
161 TH1::AddDirectory(kFALSE);
163 Int_t iNbinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
164 Double_t dPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();
165 Double_t dPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
166 Int_t iNbinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
167 Double_t dEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();
168 Double_t dEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
170 fHistProFlags = new TProfile("Flow_Flags_SP","Flow_Flags_SP",1,0,1,"s");
171 fHistProFlags->GetXaxis()->SetBinLabel(1,"fApplyCorrectionForNUA");
172 fHistList->Add(fHistProFlags);
174 fHistProUQetaRP = new TProfile("Flow_UQetaRP_SP","Flow_UQetaRP_SP",iNbinsEta,dEtaMin,dEtaMax,"s");
175 fHistProUQetaRP->SetXTitle("{eta}");
176 fHistProUQetaRP->SetYTitle("<uQ>");
177 fHistList->Add(fHistProUQetaRP);
179 fHistProUQetaPOI = new TProfile("Flow_UQetaPOI_SP","Flow_UQetaPOI_SP",iNbinsEta,dEtaMin,dEtaMax,"s");
180 fHistProUQetaPOI->SetXTitle("{eta}");
181 fHistProUQetaPOI->SetYTitle("<uQ>");
182 fHistList->Add(fHistProUQetaPOI);
184 fHistProUQPtRP = new TProfile("Flow_UQPtRP_SP","Flow_UQPtRP_SP",iNbinsPt,dPtMin,dPtMax,"s");
185 fHistProUQPtRP->SetXTitle("p_t (GeV)");
186 fHistProUQPtRP->SetYTitle("<uQ>");
187 fHistList->Add(fHistProUQPtRP);
189 fHistProUQPtPOI = new TProfile("Flow_UQPtPOI_SP","Flow_UQPtPOI_SP",iNbinsPt,dPtMin,dPtMax,"s");
190 fHistProUQPtPOI->SetXTitle("p_t (GeV)");
191 fHistProUQPtPOI->SetYTitle("<uQ>");
192 fHistList->Add(fHistProUQPtPOI);
194 fHistProQaQb = new TProfile("Flow_QaQb_SP","Flow_QaQb_SP", 1, 0.5, 1.5,"s");
195 fHistProQaQb->SetYTitle("<QaQb>");
196 fHistList->Add(fHistProQaQb);
198 fHistProQaQbNorm = new TProfile("FlowPro_QaQbNorm_SP","FlowPro_QaQbNorm_SP", 1, 0.5, 1.5,"s");
199 fHistProQaQbNorm->SetYTitle("<QaQb/MaMb>");
200 fHistList->Add(fHistProQaQbNorm);
202 fHistProQaQbReImNorm = new TProfile("FlowPro_QaQbReImNorm_SP","FlowPro_QaQbReImNorm_SP", 4, 0.5, 4.5,"s");
203 fHistProQaQbReImNorm->GetXaxis()->SetBinLabel(1,"#LT#LTsin(#phi_{a})#GT#GT");
204 fHistProQaQbReImNorm->GetXaxis()->SetBinLabel(2,"#LT#LTcos(#phi_{a})#GT#GT");
205 fHistProQaQbReImNorm->GetXaxis()->SetBinLabel(3,"#LT#LTsin(#phi_{b})#GT#GT");
206 fHistProQaQbReImNorm->GetXaxis()->SetBinLabel(4,"#LT#LTcos(#phi_{b})#GT#GT");
207 fHistList->Add(fHistProQaQbReImNorm);
209 fHistProNonIsotropicTermsQ = new TProfile("FlowPro_NonIsotropicTermsQ_SP","FlowPro_NonIsotropicTermsQ_SP", 2, 0.5, 2.5,"s");
210 fHistProNonIsotropicTermsQ->GetXaxis()->SetBinLabel(1,"#LT#LTsin(#phi_{a+b})#GT#GT");
211 fHistProNonIsotropicTermsQ->GetXaxis()->SetBinLabel(2,"#LT#LTcos(#phi_{a+b})#GT#GT");
212 fHistList->Add(fHistProNonIsotropicTermsQ);
214 TString rpPoi[2] = {"RP","POI"};
215 TString ptEta[2] = {"Pt","Eta"};
216 TString sinCos[2] = {"sin","cos"};
217 Int_t nBinsPtEta[2] = {iNbinsPt,iNbinsEta};
218 Double_t minPtEta[2] = {dPtMin,dEtaMin};
219 Double_t maxPtEta[2] = {dPtMax,dEtaMax};
220 for(Int_t rp=0;rp<2;rp++)
222 for(Int_t pe=0;pe<2;pe++)
224 for(Int_t sc=0;sc<2;sc++)
226 fHistProNonIsotropicTermsU[rp][pe][sc] = new TProfile(Form("Flow_NonIsotropicTerms_%s_%s_%s_SP",rpPoi[rp].Data(),ptEta[pe].Data(),sinCos[sc].Data()),Form("Flow_NonIsotropicTerms_%s_%s_%s_SP",rpPoi[rp].Data(),ptEta[pe].Data(),sinCos[sc].Data()),nBinsPtEta[pe],minPtEta[pe],maxPtEta[pe]);
227 fHistList->Add(fHistProNonIsotropicTermsU[rp][pe][sc]);
232 fHistSumOfLinearWeights = new TH1D("Flow_SumOfLinearWeights_SP","Flow_SumOfLinearWeights_SP",1,-0.5, 0.5);
233 fHistSumOfLinearWeights -> SetYTitle("sum (*)");
234 fHistSumOfLinearWeights -> SetXTitle("sum (Ma*Mb)");
235 fHistList->Add(fHistSumOfLinearWeights);
237 fHistSumOfQuadraticWeights = new TH1D("Flow_SumOfQuadraticWeights_SP","Flow_SumOfQuadraticWeights_SP",1,-0.5, 0.5);
238 fHistSumOfQuadraticWeights -> SetYTitle("sum (*)");
239 fHistSumOfQuadraticWeights -> SetXTitle("sum (Ma*Mb)^2");
240 fHistList->Add(fHistSumOfQuadraticWeights);
242 fHistProUQQaQbPtRP = new TProfile("Flow_UQQaQbPtRP_SP","Flow_UQQaQbPtRP_SP",iNbinsPt,dPtMin,dPtMax);
243 fHistProUQQaQbPtRP -> SetYTitle("<*>");
244 fHistProUQQaQbPtRP -> SetXTitle("<Qu QaQb>");
245 fHistList->Add(fHistProUQQaQbPtRP);
247 fHistProUQQaQbEtaRP = new TProfile("Flow_UQQaQbEtaRP_SP","Flow_UQQaQbEtaRP_SP",iNbinsEta,dEtaMin,dEtaMax);
248 fHistProUQQaQbEtaRP -> SetYTitle("<*>");
249 fHistProUQQaQbEtaRP -> SetXTitle("<Qu QaQb>");
250 fHistList->Add(fHistProUQQaQbEtaRP);
252 fHistProUQQaQbPtPOI = new TProfile("Flow_UQQaQbPtPOI_SP","Flow_UQQaQbPtPOI_SP",iNbinsPt,dPtMin,dPtMax);
253 fHistProUQQaQbPtPOI -> SetYTitle("<*>");
254 fHistProUQQaQbPtPOI -> SetXTitle("<Qu QaQb>");
255 fHistList->Add(fHistProUQQaQbPtPOI);
257 fHistProUQQaQbEtaPOI = new TProfile("Flow_UQQaQbEtaPOI_SP","Flow_UQQaQbEtaPOI_SP",iNbinsEta,dEtaMin,dEtaMax);
258 fHistProUQQaQbEtaPOI -> SetYTitle("<*>");
259 fHistProUQQaQbEtaPOI -> SetXTitle("<Qu QaQb>");
260 fHistList->Add(fHistProUQQaQbEtaPOI);
262 TString weightFlag[3] = {"w_Qu_","w_Qu^2_","w_QuQaQb_"};
263 for(Int_t i=0;i<3;i++)
265 fHistSumOfWeightsPtRP[i] = new TH1D(Form("Flow_SumOfWeights%sPtRP_SP",weightFlag[i].Data()),
266 Form("Flow_SumOfWeights%sPtRP_SP",weightFlag[i].Data()),iNbinsPt,dPtMin,dPtMax);
267 fHistSumOfWeightsPtRP[i] -> SetYTitle("sum (*)");
268 fHistSumOfWeightsPtRP[i] -> SetXTitle("p_{T}");
269 fHistList->Add(fHistSumOfWeightsPtRP[i]);
271 fHistSumOfWeightsEtaRP[i] = new TH1D(Form("Flow_SumOfWeights%sEtaRP_SP",weightFlag[i].Data()),
272 Form("Flow_SumOfWeights%sEtaRP_SP",weightFlag[i].Data()),iNbinsEta,dEtaMin,dEtaMax);
273 fHistSumOfWeightsEtaRP[i] -> SetYTitle("sum (*)");
274 fHistSumOfWeightsEtaRP[i] -> SetXTitle("#eta");
275 fHistList->Add(fHistSumOfWeightsEtaRP[i]);
277 fHistSumOfWeightsPtPOI[i] = new TH1D(Form("Flow_SumOfWeights%sPtPOI_SP",weightFlag[i].Data()),
278 Form("Flow_SumOfWeights%sPtPOI_SP",weightFlag[i].Data()),iNbinsPt,dPtMin,dPtMax);
279 fHistSumOfWeightsPtPOI[i] -> SetYTitle("sum (*)");
280 fHistSumOfWeightsPtPOI[i] -> SetXTitle("p_{T}");
281 fHistList->Add(fHistSumOfWeightsPtPOI[i]);
283 fHistSumOfWeightsEtaPOI[i] = new TH1D(Form("Flow_SumOfWeights%sEtaPOI_SP",weightFlag[i].Data()),
284 Form("Flow_SumOfWeights%sEtaPOI_SP",weightFlag[i].Data()),iNbinsEta,dEtaMin,dEtaMax);
285 fHistSumOfWeightsEtaPOI[i] -> SetYTitle("sum (*)");
286 fHistSumOfWeightsEtaPOI[i] -> SetXTitle("#eta");
287 fHistList->Add(fHistSumOfWeightsEtaPOI[i]);
290 fCommonHists = new AliFlowCommonHist("AliFlowCommonHistSP");
291 fHistList->Add(fCommonHists);
292 fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsSP");
293 fHistList->Add(fCommonHistsRes);
295 fHistQaQb = new TH1D("Flow_QaQb_SP","Flow_QaQb_SP",200,-100.,100.);
296 fHistQaQb -> SetYTitle("dN/dQaQb");
297 fHistQaQb -> SetXTitle("QaQb");
298 fHistList->Add(fHistQaQb);
300 fHistQaQbNorm = new TH1D("Flow_QaQbNorm_SP","Flow_QaQbNorm_SP",44,-1.1,1.1);
301 fHistQaQbNorm -> SetYTitle("dN/d(QaQb/MaMb)");
302 fHistQaQbNorm -> SetXTitle("QaQb/MaMb");
303 fHistList->Add(fHistQaQbNorm);
305 fHistQaQbCos = new TH1D("Flow_QaQbCos_SP","Flow_QaQbCos_SP",63,0.,TMath::Pi());
306 fHistQaQbCos -> SetYTitle("dN/d(#phi)");
307 fHistQaQbCos -> SetXTitle("#phi");
308 fHistList->Add(fHistQaQbCos);
310 fHistResolution = new TH1D("Flow_resolution_SP","Flow_resolution_SP",100,-1.0,1.0);
311 fHistResolution -> SetYTitle("dN/d(cos(2(#phi_a - #phi_b))");
312 fHistResolution -> SetXTitle("cos(2*(#phi_a - #phi_b))");
313 fHistList->Add(fHistResolution);
315 fHistQaNorm = new TH1D("Flow_QaNorm_SP","Flow_QaNorm_SP",22,0.,1.1);
316 fHistQaNorm -> SetYTitle("dN/d(|Qa/Ma|)");
317 fHistQaNorm -> SetXTitle("|Qa/Ma|");
318 fHistList->Add(fHistQaNorm);
320 fHistQaNormvsMa = new TH2D("Flow_QaNormvsMa_SP","Flow_QaNormvsMa_SP",100,0.,100.,22,0.,1.1);
321 fHistQaNormvsMa -> SetYTitle("|Qa/Ma|");
322 fHistQaNormvsMa -> SetXTitle("Ma");
323 fHistList->Add(fHistQaNormvsMa);
325 fHistQbNorm = new TH1D("Flow_QbNorm_SP","Flow_QbNorm_SP",22,0.,1.1);
326 fHistQbNorm -> SetYTitle("dN/d(|Qb/Mb|)");
327 fHistQbNorm -> SetXTitle("|Qb/Mb|");
328 fHistList->Add(fHistQbNorm);
330 fHistQbNormvsMb = new TH2D("Flow_QbNormvsMb_SP","Flow_QbNormvsMb_SP",100,0.,100.,22,0.,1.1);
331 fHistQbNormvsMb -> SetYTitle("|Qb/Mb|");
332 fHistQbNormvsMb -> SetXTitle("|Mb|");
333 fHistList->Add(fHistQbNormvsMb);
335 fHistMavsMb = new TH2D("Flow_MavsMb_SP","Flow_MavsMb_SP",100,0.,100.,100,0.,100.);
336 fHistMavsMb -> SetYTitle("Ma");
337 fHistMavsMb -> SetXTitle("Mb");
338 fHistList->Add(fHistMavsMb);
344 cout<<"WARNING: fWeightsList is NULL in the Scalar Product method."<<endl;
347 if(fWeightsList->FindObject("phi_weights_sub0")) {
348 fPhiWeightsSub0 = dynamic_cast<TH1F*>
349 (fWeightsList->FindObject("phi_weights_sub0"));
350 fHistList->Add(fPhiWeightsSub0);
352 cout<<"WARNING: histogram with phi weights is not accessible in Scalar Product"<<endl;
355 if(fWeightsList->FindObject("phi_weights_sub1")) {
356 fPhiWeightsSub1 = dynamic_cast<TH1F*>
357 (fWeightsList->FindObject("phi_weights_sub1"));
358 fHistList->Add(fPhiWeightsSub1);
360 cout<<"WARNING: histogram with phi weights is not accessible in Scalar Product"<<endl;
364 } // end of if(fUsePhiWeights)
366 fEventNumber = 0; //set number of events to zero
368 //store all boolean flags needed in Finish():
371 TH1::AddDirectory(oldHistAddStatus);
374 //-----------------------------------------------------------------------
375 void AliFlowAnalysisWithScalarProduct::Make(AliFlowEventSimple* anEvent) {
377 //Calculate flow based on <QaQb/MaMb> = <v^2>
382 //get Q vectors for the eta-subevents
383 AliFlowVector* vQarray = new AliFlowVector[2];
384 if (fUsePhiWeights) {
385 anEvent->Get2Qsub(vQarray,2,fWeightsList,kTRUE);
387 anEvent->Get2Qsub(vQarray);
390 AliFlowVector vQa = vQarray[0];
392 AliFlowVector vQb = vQarray[1];
394 //check that the subevents are not empty:
395 Double_t dMa = vQa.GetMult();
396 Double_t dMb = vQb.GetMult();
397 if (dMa != 0 && dMb != 0) {
400 //request that the subevent multiplicities are not too different
401 Double_t dRelDiff = TMath::Abs((dMa - dMb)/(dMa + dMb));
402 if (dRelDiff < fRelDiffMsub) {
404 //fill control histograms
405 fCommonHists->FillControlHistograms(anEvent);
407 //fill some SP control histograms
408 fHistProQaQb -> Fill(1.,vQa*vQb,1.); //Fill with weight 1 -> Weight with MaMb????
409 fHistQaQbCos ->Fill(TMath::ACos((vQa/vQa.Mod())*(vQb/vQb.Mod()))); //Acos(Qa*Qb) = angle
410 fHistResolution -> Fill(TMath::Cos( vQa.Phi()- vQb.Phi() )); //vQa.Phi() returns 2*phi
411 fHistQaQb -> Fill(vQa*vQb);
412 fHistMavsMb -> Fill(dMb,dMa);
414 //get total Q vector = the sum of subevent a and subevent b
415 AliFlowVector vQ = vQa + vQb;
416 //needed to correct for non-uniform acceptance:
419 fHistProNonIsotropicTermsQ->Fill(1.,vQ.Y()/(dMa+dMb),dMa+dMb);
420 fHistProNonIsotropicTermsQ->Fill(2.,vQ.X()/(dMa+dMb),dMa+dMb);
422 //weight the Q vectors for the subevents by the multiplicity
423 //Note: Weight Q only in the particle loop when it is clear if it should be (m-1) or M
424 Double_t dQXa = vQa.X()/dMa;
425 Double_t dQYa = vQa.Y()/dMa;
428 Double_t dQXb = vQb.X()/dMb;
429 Double_t dQYb = vQb.Y()/dMb;
432 //scalar product of the two subevents
433 Double_t dQaQb = (vQa*vQb);
434 fHistProQaQbNorm -> Fill(1.,dQaQb,dMa*dMb); //Fill (QaQb/MaMb) with weight (MaMb).
435 //needed for the error calculation:
436 fHistSumOfLinearWeights -> Fill(0.,dMa*dMb);
437 fHistSumOfQuadraticWeights -> Fill(0.,pow(dMa*dMb,2.));
438 //needed for correcting non-uniform acceptance:
439 fHistProQaQbReImNorm->Fill(1.,dQYa,dMa); // to get <<sin(phi_a)>>
440 fHistProQaQbReImNorm->Fill(2.,dQXa,dMa); // to get <<cos(phi_a)>>
441 fHistProQaQbReImNorm->Fill(3.,dQYb,dMb); // to get <<sin(phi_b)>>
442 fHistProQaQbReImNorm->Fill(4.,dQXb,dMb); // to get <<cos(phi_b)>>
444 //fill some SP control histograms
445 fHistQaQbNorm ->Fill(vQa*vQb);
446 fHistQaNorm ->Fill(vQa.Mod());
447 fHistQaNormvsMa->Fill(dMa,vQa.Mod());
448 fHistQbNorm ->Fill(vQb.Mod());
449 fHistQbNormvsMb->Fill(dMb,vQb.Mod());
451 //loop over the tracks of the event
452 AliFlowTrackSimple* pTrack = NULL;
453 Int_t iNumberOfTracks = anEvent->NumberOfTracks();
454 Double_t dMq = vQ.GetMult();
456 for (Int_t i=0;i<iNumberOfTracks;i++)
458 pTrack = anEvent->GetTrack(i) ;
460 Double_t dPhi = pTrack->Phi();
461 //set default phi weight to 1
463 //phi weight of pTrack
464 if(fUsePhiWeights && fPhiWeightsSub0 && fPhiWeightsSub1) {
465 if (pTrack->InSubevent(0) ) {
466 Int_t iNBinsPhiSub0 = fPhiWeightsSub0->GetNbinsX();
467 dW = fPhiWeightsSub0->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*iNBinsPhiSub0/TMath::TwoPi())));
469 else if ( pTrack->InSubevent(1)) {
470 Int_t iNBinsPhiSub1 = fPhiWeightsSub1->GetNbinsX();
471 dW = fPhiWeightsSub1->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*iNBinsPhiSub1/TMath::TwoPi())));
473 //bin = 1 + value*nbins/range
474 //TMath::Floor rounds to the lower integer
479 Double_t dUX = TMath::Cos(2*dPhi);
480 Double_t dUY = TMath::Sin(2*dPhi);
482 Double_t dModulus = vU.Mod();
483 if (dModulus!=0.) vU.Set(dUX/dModulus,dUY/dModulus); // make length 1
484 else cerr<<"dModulus is zero!"<<endl;
486 //redefine the Q vector and devide by its multiplicity
490 //subtract particle from the flowvector if used to define it
491 if (pTrack->InSubevent(0) || pTrack->InSubevent(1)) {
492 dQmX = (vQ.X() - dW*(pTrack->Weight())*dUX)/(dMq-1);
493 dQmY = (vQ.Y() - dW*(pTrack->Weight())*dUY)/(dMq-1);
496 //dUQ = scalar product of vU and vQm
497 Double_t dUQ = (vU * vQm);
498 Double_t dPt = pTrack->Pt();
499 Double_t dEta = pTrack->Eta();
500 //fill the profile histograms
501 if (pTrack->InRPSelection()) {
502 fHistProUQetaRP -> Fill(dEta,dUQ,(dMq-1)); //Fill (Qu/(Mq-1)) with weight (Mq-1)
503 //needed for the error calculation:
504 fHistProUQQaQbEtaRP -> Fill(dEta,dUQ*dQaQb,(dMq-1)*dMa*dMb); //Fill [Qu/(Mq-1)]*[QaQb/MaMb] with weight (Mq-1)MaMb
505 fHistProUQPtRP -> Fill(dPt,dUQ,(dMq-1)); //Fill (Qu/(Mq-1)) with weight (Mq-1)
506 fHistProUQQaQbPtRP -> Fill(dPt,dUQ*dQaQb,(dMq-1)*dMa*dMb); //Fill [Qu/(Mq-1)]*[QaQb/MaMb] with weight (Mq-1)MaMb
508 fHistSumOfWeightsEtaRP[0]->Fill(dEta,(dMq-1)); // sum of Mq-1
509 fHistSumOfWeightsEtaRP[1]->Fill(dEta,pow((dMq-1),2.));// sum of (Mq-1)^2
510 fHistSumOfWeightsEtaRP[2]->Fill(dEta,(dMq-1)*dMa*dMb);// sum of (Mq-1)*MaMb
511 fHistSumOfWeightsPtRP[0]->Fill(dPt,(dMq-1)); // sum of Mq-1
512 fHistSumOfWeightsPtRP[1]->Fill(dPt,pow((dMq-1),2.)); // sum of (Mq-1)^2
513 fHistSumOfWeightsPtRP[2]->Fill(dPt,(dMq-1)*dMa*dMb); // sum of (Mq-1)*MaMb
514 //nonisotropic terms:
515 fHistProNonIsotropicTermsU[0][0][0]->Fill(dPt,dUY,1.);
516 fHistProNonIsotropicTermsU[0][0][1]->Fill(dPt,dUX,1.);
517 fHistProNonIsotropicTermsU[0][1][0]->Fill(dEta,dUY,1.);
518 fHistProNonIsotropicTermsU[0][1][1]->Fill(dEta,dUX,1.);
520 if (pTrack->InPOISelection()) {
521 fHistProUQetaPOI -> Fill(dEta,dUQ,(dMq-1));//Fill (Qu/(Mq-1)) with weight (Mq-1)
522 //needed for the error calculation:
523 fHistProUQQaQbEtaPOI -> Fill(dEta,dUQ*dQaQb,(dMq-1)*dMa*dMb); //Fill [Qu/(Mq-1)]*[QaQb/MaMb] with weight (Mq-1)MaMb
524 fHistProUQPtPOI -> Fill(dPt,dUQ,(dMq-1)); //Fill (Qu/(Mq-1)) with weight (Mq-1)
525 fHistProUQQaQbPtPOI -> Fill(dPt,dUQ*dQaQb,(dMq-1)*dMa*dMb); //Fill [Qu/(Mq-1)]*[QaQb/MaMb] with weight (Mq-1)MaMb
527 fHistSumOfWeightsEtaPOI[0]->Fill(dEta,(dMq-1)); // sum of Mq-1
528 fHistSumOfWeightsEtaPOI[1]->Fill(dEta,pow((dMq-1),2.));// sum of (Mq-1)^2
529 fHistSumOfWeightsEtaPOI[2]->Fill(dEta,(dMq-1)*dMa*dMb);// sum of (Mq-1)*MaMb
530 fHistSumOfWeightsPtPOI[0]->Fill(dPt,(dMq-1)); // sum of Mq-1
531 fHistSumOfWeightsPtPOI[1]->Fill(dPt,pow((dMq-1),2.)); // sum of (Mq-1)^2
532 fHistSumOfWeightsPtPOI[2]->Fill(dPt,(dMq-1)*dMa*dMb); // sum of (Mq-1)*MaMb
533 //nonisotropic terms:
534 fHistProNonIsotropicTermsU[1][0][0]->Fill(dPt,dUY,1.);
535 fHistProNonIsotropicTermsU[1][0][1]->Fill(dPt,dUX,1.);
536 fHistProNonIsotropicTermsU[1][1][0]->Fill(dEta,dUY,1.);
537 fHistProNonIsotropicTermsU[1][1][1]->Fill(dEta,dUX,1.);
540 } else { //do not subtract the particle from the flowvector
545 //dUQ = scalar product of vU and vQm
546 Double_t dUQ = (vU * vQm);
547 Double_t dPt = pTrack->Pt();
548 Double_t dEta = pTrack->Eta();
549 //fill the profile histograms
550 if (pTrack->InRPSelection()) {
551 fHistProUQetaRP -> Fill(dEta,dUQ,dMq); //Fill (Qu/Mq) with weight Mq
552 //needed for the error calculation:
553 fHistProUQQaQbEtaRP -> Fill(dEta,dUQ*dQaQb,dMq*dMa*dMb); //Fill [Qu/Mq]*[QaQb/MaMb] with weight Mq*MaMb
554 fHistProUQPtRP -> Fill(dPt,dUQ,dMq); //Fill (Qu/Mq) with weight Mq
555 fHistProUQQaQbPtRP -> Fill(dPt,dUQ*dQaQb,dMq*dMa*dMb); //Fill [Qu/Mq]*[QaQb/MaMb] with weight Mq*MaMb
557 fHistSumOfWeightsEtaRP[0]->Fill(dEta,dMq); // sum of Mq
558 fHistSumOfWeightsEtaRP[1]->Fill(dEta,pow(dMq,2.));// sum of Mq^2
559 fHistSumOfWeightsEtaRP[2]->Fill(dEta,dMq*dMa*dMb);// sum of Mq*MaMb
560 fHistSumOfWeightsPtRP[0]->Fill(dPt,dMq); // sum of Mq
561 fHistSumOfWeightsPtRP[1]->Fill(dPt,pow(dMq,2.)); // sum of Mq^2
562 fHistSumOfWeightsPtRP[2]->Fill(dPt,dMq*dMa*dMb); // sum of Mq*MaMb
563 //nonisotropic terms:
564 fHistProNonIsotropicTermsU[0][0][0]->Fill(dPt,dUY,1.);
565 fHistProNonIsotropicTermsU[0][0][1]->Fill(dPt,dUX,1.);
566 fHistProNonIsotropicTermsU[0][1][0]->Fill(dEta,dUY,1.);
567 fHistProNonIsotropicTermsU[0][1][1]->Fill(dEta,dUX,1.);
569 if (pTrack->InPOISelection()) {
570 fHistProUQetaPOI -> Fill(dEta,dUQ,dMq); //Fill (Qu/Mq) with weight Mq
571 //needed for the error calculation:
572 fHistProUQQaQbEtaPOI -> Fill(dEta,dUQ*dQaQb,dMq*dMa*dMb); //Fill [Qu/Mq]*[QaQb/MaMb] with weight Mq*MaMb
573 fHistProUQPtPOI -> Fill(dPt,dUQ,dMq); //Fill (Qu/Mq) with weight Mq
574 fHistProUQQaQbPtPOI -> Fill(dPt,dUQ*dQaQb,dMq*dMa*dMb); //Fill [Qu/Mq]*[QaQb/MaMb] with weight Mq*MaMb
576 fHistSumOfWeightsEtaPOI[0]->Fill(dEta,dMq); // sum of Mq
577 fHistSumOfWeightsEtaPOI[1]->Fill(dEta,pow(dMq,2.));// sum of Mq^2
578 fHistSumOfWeightsEtaPOI[2]->Fill(dEta,dMq*dMa*dMb);// sum of Mq*MaMb
579 fHistSumOfWeightsPtPOI[0]->Fill(dPt,dMq); // sum of Mq
580 fHistSumOfWeightsPtPOI[1]->Fill(dPt,pow(dMq,2.)); // sum of Mq^2
581 fHistSumOfWeightsPtPOI[2]->Fill(dPt,dMq*dMa*dMb); // sum of Mq*MaMb
582 //nonisotropic terms:
583 fHistProNonIsotropicTermsU[1][0][0]->Fill(dPt,dUY,1.);
584 fHistProNonIsotropicTermsU[1][0][1]->Fill(dPt,dUX,1.);
585 fHistProNonIsotropicTermsU[1][1][0]->Fill(dEta,dUY,1.);
586 fHistProNonIsotropicTermsU[1][1][1]->Fill(dEta,dUX,1.);
588 }//track not in subevents
596 } //difference Ma and Mb
598 }// subevents not empty
605 //--------------------------------------------------------------------
606 void AliFlowAnalysisWithScalarProduct::GetOutputHistograms(TList *outputListHistos){
608 //get pointers to all output histograms (called before Finish())
610 if (outputListHistos) {
611 //Get the common histograms from the output list
612 AliFlowCommonHist *pCommonHist = dynamic_cast<AliFlowCommonHist*>
613 (outputListHistos->FindObject("AliFlowCommonHistSP"));
614 AliFlowCommonHistResults *pCommonHistResults = dynamic_cast<AliFlowCommonHistResults*>
615 (outputListHistos->FindObject("AliFlowCommonHistResultsSP"));
616 TProfile* pHistProQaQb = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_QaQb_SP"));
617 TProfile* pHistProQaQbNorm = dynamic_cast<TProfile*>(outputListHistos->FindObject("FlowPro_QaQbNorm_SP"));
618 TProfile* pHistProQaQbReImNorm = dynamic_cast<TProfile*>(outputListHistos->FindObject("FlowPro_QaQbReImNorm_SP"));
619 TProfile* pHistProNonIsotropicTermsQ = dynamic_cast<TProfile*>(outputListHistos->FindObject("FlowPro_NonIsotropicTermsQ_SP"));
620 TH1D* pHistSumOfLinearWeights = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_SumOfLinearWeights_SP"));
621 TH1D* pHistSumOfQuadraticWeights = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_SumOfQuadraticWeights_SP"));
623 TProfile* pHistProFlags = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_Flags_SP"));
624 TProfile* pHistProUQetaRP = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQetaRP_SP"));
625 TProfile* pHistProUQetaPOI = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQetaPOI_SP"));
626 TProfile* pHistProUQPtRP = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQPtRP_SP"));
627 TProfile* pHistProUQPtPOI = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQPtPOI_SP"));
628 TProfile* pHistProUQQaQbPtRP = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQQaQbPtRP_SP"));
629 TProfile* pHistProUQQaQbEtaRP = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQQaQbEtaRP_SP"));
630 TProfile* pHistProUQQaQbPtPOI = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQQaQbPtPOI_SP"));
631 TProfile* pHistProUQQaQbEtaPOI = dynamic_cast<TProfile*>(outputListHistos->FindObject("Flow_UQQaQbEtaPOI_SP"));
632 TString weightFlag[3] = {"w_Qu_","w_Qu^2_","w_QuQaQb_"};
633 TH1D* pHistSumOfWeightsPtRP[3] = {NULL};
634 TH1D* pHistSumOfWeightsEtaRP[3] = {NULL};
635 TH1D* pHistSumOfWeightsPtPOI[3] = {NULL};
636 TH1D* pHistSumOfWeightsEtaPOI[3] = {NULL};
637 for(Int_t i=0;i<3;i++)
639 pHistSumOfWeightsPtRP[i] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("Flow_SumOfWeights%sPtRP_SP",weightFlag[i].Data())));
640 pHistSumOfWeightsEtaRP[i] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("Flow_SumOfWeights%sEtaRP_SP",weightFlag[i].Data())));
641 pHistSumOfWeightsPtPOI[i] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("Flow_SumOfWeights%sPtPOI_SP",weightFlag[i].Data())));
642 pHistSumOfWeightsEtaPOI[i] = dynamic_cast<TH1D*>(outputListHistos->FindObject(Form("Flow_SumOfWeights%sEtaPOI_SP",weightFlag[i].Data())));
644 TString rpPoi[2] = {"RP","POI"};
645 TString ptEta[2] = {"Pt","Eta"};
646 TString sinCos[2] = {"sin","cos"};
647 TProfile *pHistProNonIsotropicTermsU[2][2][2] = {{{NULL}}};
648 for(Int_t rp=0;rp<2;rp++)
650 for(Int_t pe=0;pe<2;pe++)
652 for(Int_t sc=0;sc<2;sc++)
654 pHistProNonIsotropicTermsU[rp][pe][sc] = dynamic_cast<TProfile*>(outputListHistos->FindObject(Form("Flow_NonIsotropicTerms_%s_%s_%s_SP",rpPoi[rp].Data(),ptEta[pe].Data(),sinCos[sc].Data())));
659 TH1D* pHistQaQb = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_QaQb_SP"));
660 TH1D* pHistQaQbNorm = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_QaQbNorm_SP"));
661 TH1D* pHistQaQbCos = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_QaQbCos_SP"));
662 TH1D* pHistResolution = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_resolution_SP"));
663 TH1D* pHistQaNorm = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_QaNorm_SP"));
664 TH2D* pHistQaNormvsMa = dynamic_cast<TH2D*>(outputListHistos->FindObject("Flow_QaNormvsMa_SP"));
665 TH1D* pHistQbNorm = dynamic_cast<TH1D*>(outputListHistos->FindObject("Flow_QbNorm_SP"));
666 TH2D* pHistQbNormvsMb = dynamic_cast<TH2D*>(outputListHistos->FindObject("Flow_QbNormvsMb_SP"));
667 TH2D* pHistMavsMb = dynamic_cast<TH2D*>(outputListHistos->FindObject("Flow_MavsMb_SP"));
669 //pass the pointers to the task
670 if (pCommonHist && pCommonHistResults && pHistProQaQb && pHistProQaQbNorm &&
671 pHistProQaQbReImNorm && pHistProNonIsotropicTermsQ &&
672 pHistSumOfLinearWeights && pHistSumOfQuadraticWeights &&
673 pHistQaQb && pHistQaQbNorm && pHistQaQbCos && pHistResolution &&
674 pHistQaNorm && pHistQaNormvsMa && pHistQbNorm && pHistQbNormvsMb &&
675 pHistMavsMb && pHistProFlags &&
676 pHistProUQetaRP && pHistProUQetaPOI &&
677 pHistProUQPtRP && pHistProUQPtPOI &&
678 pHistProUQQaQbPtRP && pHistProUQQaQbEtaRP &&
679 pHistProUQQaQbPtPOI && pHistProUQQaQbEtaPOI) {
680 this -> SetCommonHists(pCommonHist);
681 this -> SetCommonHistsRes(pCommonHistResults);
682 this -> SetHistProQaQb(pHistProQaQb);
683 this -> SetHistProQaQbNorm(pHistProQaQbNorm);
684 this -> SetHistProQaQbReImNorm(pHistProQaQbReImNorm);
685 this -> SetHistProNonIsotropicTermsQ(pHistProNonIsotropicTermsQ);
686 this -> SetHistSumOfLinearWeights(pHistSumOfLinearWeights);
687 this -> SetHistSumOfQuadraticWeights(pHistSumOfQuadraticWeights);
688 this -> SetHistQaQb(pHistQaQb);
689 this -> SetHistQaQbNorm(pHistQaQbNorm);
690 this -> SetHistQaQbCos(pHistQaQbCos);
691 this -> SetHistResolution(pHistResolution);
692 this -> SetHistQaNorm(pHistQaNorm);
693 this -> SetHistQaNormvsMa(pHistQaNormvsMa);
694 this -> SetHistQbNorm(pHistQbNorm);
695 this -> SetHistQbNormvsMb(pHistQbNormvsMb);
696 this -> SetHistMavsMb(pHistMavsMb);
697 this -> SetHistProFlags(pHistProFlags);
698 this -> SetHistProUQetaRP(pHistProUQetaRP);
699 this -> SetHistProUQetaPOI(pHistProUQetaPOI);
700 this -> SetHistProUQPtRP(pHistProUQPtRP);
701 this -> SetHistProUQPtPOI(pHistProUQPtPOI);
702 this -> SetHistProUQQaQbPtRP(pHistProUQQaQbPtRP);
703 this -> SetHistProUQQaQbEtaRP(pHistProUQQaQbEtaRP);
704 this -> SetHistProUQQaQbPtPOI(pHistProUQQaQbPtPOI);
705 this -> SetHistProUQQaQbEtaPOI(pHistProUQQaQbEtaPOI);
707 for(Int_t i=0;i<3;i++) {
708 if(pHistSumOfWeightsPtRP[i]) this -> SetHistSumOfWeightsPtRP(pHistSumOfWeightsPtRP[i],i);
709 if(pHistSumOfWeightsEtaRP[i]) this -> SetHistSumOfWeightsEtaRP(pHistSumOfWeightsEtaRP[i],i);
710 if(pHistSumOfWeightsPtPOI[i]) this -> SetHistSumOfWeightsPtPOI(pHistSumOfWeightsPtPOI[i],i);
711 if(pHistSumOfWeightsEtaPOI[i]) this -> SetHistSumOfWeightsEtaPOI(pHistSumOfWeightsEtaPOI[i],i);
713 for(Int_t rp=0;rp<2;rp++)
715 for(Int_t pe=0;pe<2;pe++)
717 for(Int_t sc=0;sc<2;sc++)
719 if(pHistProNonIsotropicTermsU[rp][pe][sc]) this->SetHistProNonIsotropicTermsU(pHistProNonIsotropicTermsU[rp][pe][sc],rp,pe,sc);
725 cout<<"WARNING: Histograms needed to run Finish() in SP are not accessable!"<<endl; }
727 } // end of if(outputListHistos)
730 //--------------------------------------------------------------------
731 void AliFlowAnalysisWithScalarProduct::Finish() {
733 //calculate flow and fill the AliFlowCommonHistResults
734 if (fDebug) cout<<"AliFlowAnalysisWithScalarProduct::Finish()"<<endl;
736 //access all boolean flags needed in Finish():
741 for(Int_t b=1;b<=10;b++)
743 cout<<"b = "<<b<<": "<<fHistProUQPtRP->GetBinContent(b)<<endl;
747 cout<<"*************************************"<<endl;
748 cout<<"*************************************"<<endl;
749 cout<<" Integrated flow from "<<endl;
750 cout<<" Scalar product "<<endl;
753 Int_t iNbinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
754 Int_t iNbinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
756 //Calculate the event plane resolution
757 //----------------------------------
758 Double_t dCos2phi = fHistResolution->GetMean();
760 Double_t dResolution = TMath::Sqrt(2*dCos2phi);
761 cout<<"An estimate of the event plane resolution is: "<<dResolution<<endl;
765 //Calculate reference flow (noname)
766 //----------------------------------
767 //weighted average over (QaQb/MaMb) with weight (MaMb)
768 Double_t dQaQb = fHistProQaQbNorm->GetBinContent(1);
769 //non-isotropic terms:
770 Double_t dImQa = fHistProQaQbReImNorm->GetBinContent(1);
771 Double_t dReQa = fHistProQaQbReImNorm->GetBinContent(2);
772 Double_t dImQb = fHistProQaQbReImNorm->GetBinContent(3);
773 Double_t dReQb = fHistProQaQbReImNorm->GetBinContent(4);
775 if(fApplyCorrectionForNUA)
777 dQaQb = dQaQb - dImQa*dImQb - dReQa*dReQb;
783 dV = TMath::Sqrt(dQaQb);
785 //statistical error of dQaQb:
786 // statistical error = term1 * spread * term2:
787 // term1 = sqrt{sum_{i=1}^{N} w^2}/(sum_{i=1}^{N} w)
788 // term2 = 1/sqrt(1-term1^2)
789 Double_t dSpreadQaQb = fHistProQaQbNorm->GetBinError(1);
790 Double_t dSumOfLinearWeights = fHistSumOfLinearWeights->GetBinContent(1);
791 Double_t dSumOfQuadraticWeights = fHistSumOfQuadraticWeights->GetBinContent(1);
792 Double_t dTerm1 = 0.;
793 Double_t dTerm2 = 0.;
794 if(dSumOfLinearWeights) {
795 dTerm1 = pow(dSumOfQuadraticWeights,0.5)/dSumOfLinearWeights;
797 if(1.-pow(dTerm1,2.) > 0.) {
798 dTerm2 = 1./pow(1-pow(dTerm1,2.),0.5);
800 Double_t dStatErrorQaQb = dTerm1 * dSpreadQaQb * dTerm2;
801 //calculate the statistical error
804 dVerr = (1./(2.*pow(dQaQb,0.5)))*dStatErrorQaQb;
806 fCommonHistsRes->FillIntegratedFlow(dV,dVerr);
807 cout<<"dV = "<<dV<<" +- "<<dVerr<<endl;
809 //Calculate differential flow and integrated flow (RP, POI)
810 //---------------------------------------------------------
811 //v as a function of eta for RP selection
812 for(Int_t b=1;b<iNbinsEta+1;b++) {
813 Double_t duQpro = fHistProUQetaRP->GetBinContent(b);
814 if(fApplyCorrectionForNUA)
817 - fHistProNonIsotropicTermsU[0][1][1]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(2)
818 - fHistProNonIsotropicTermsU[0][1][0]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(1);
820 Double_t dv2pro = -999.;
821 if (dV!=0.) { dv2pro = duQpro/dV; }
822 //calculate the statistical error
823 Double_t dv2ProErr = CalculateStatisticalError(b, dStatErrorQaQb, fHistProUQetaRP, fHistProUQQaQbEtaRP, fHistSumOfWeightsEtaRP);
825 fCommonHistsRes->FillDifferentialFlowEtaRP(b, dv2pro, dv2ProErr);
829 //v as a function of eta for POI selection
830 for(Int_t b=1;b<iNbinsEta+1;b++) {
831 Double_t duQpro = fHistProUQetaPOI->GetBinContent(b);
832 if(fApplyCorrectionForNUA)
835 - fHistProNonIsotropicTermsU[1][1][1]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(2)
836 - fHistProNonIsotropicTermsU[1][1][0]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(1);
838 Double_t dv2pro = -999.;
839 if (dV!=0.) { dv2pro = duQpro/dV; }
840 //calculate the statistical error
841 Double_t dv2ProErr = CalculateStatisticalError(b, dStatErrorQaQb, fHistProUQetaPOI, fHistProUQQaQbEtaPOI, fHistSumOfWeightsEtaPOI);
844 fCommonHistsRes->FillDifferentialFlowEtaPOI(b, dv2pro, dv2ProErr);
848 //v as a function of Pt for RP selection
849 TH1F* fHistPtRP = fCommonHists->GetHistPtRP(); //for calculating integrated flow
851 Double_t dSumRP = 0.;
852 Double_t dErrVRP =0.;
854 for(Int_t b=1;b<iNbinsPt+1;b++) {
855 Double_t duQpro = fHistProUQPtRP->GetBinContent(b);
856 if(fApplyCorrectionForNUA)
859 - fHistProNonIsotropicTermsU[0][0][1]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(2)
860 - fHistProNonIsotropicTermsU[0][0][0]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(1);
862 Double_t dv2pro = -999.;
863 if (dV!=0.) { dv2pro = duQpro/dV; }
864 //calculate the statistical error
865 Double_t dv2ProErr = CalculateStatisticalError(b, dStatErrorQaQb, fHistProUQPtRP, fHistProUQQaQbPtRP, fHistSumOfWeightsPtRP);
868 fCommonHistsRes->FillDifferentialFlowPtRP(b, dv2pro, dv2ProErr);
870 //calculate integrated flow for RP selection
872 Double_t dYieldPt = fHistPtRP->GetBinContent(b);
873 dVRP += dv2pro*dYieldPt;
875 dErrVRP += dYieldPt*dYieldPt*dv2ProErr*dv2ProErr;
876 } else { cout<<"fHistPtRP is NULL"<<endl; }
880 dVRP /= dSumRP; //the pt distribution should be normalised
881 dErrVRP /= (dSumRP*dSumRP);
882 dErrVRP = TMath::Sqrt(dErrVRP);
884 fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrVRP);
885 cout<<"dV(RP) = "<<dVRP<<" +- "<<dErrVRP<<endl;
888 //v as a function of Pt for POI selection
889 TH1F* fHistPtPOI = fCommonHists->GetHistPtPOI(); //for calculating integrated flow
891 Double_t dSumPOI = 0.;
892 Double_t dErrVPOI =0.;
894 for(Int_t b=1;b<iNbinsPt+1;b++) {
895 Double_t duQpro = fHistProUQPtPOI->GetBinContent(b);
896 if(fApplyCorrectionForNUA)
899 - fHistProNonIsotropicTermsU[1][0][1]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(2)
900 - fHistProNonIsotropicTermsU[1][0][0]->GetBinContent(b)*fHistProNonIsotropicTermsQ->GetBinContent(1);
902 Double_t dv2pro = -999.;
903 if (dV!=0.) { dv2pro = duQpro/dV; }
904 //calculate the statistical error
905 Double_t dv2ProErr = CalculateStatisticalError(b, dStatErrorQaQb, fHistProUQPtPOI, fHistProUQQaQbPtPOI, fHistSumOfWeightsPtPOI);
908 fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2pro, dv2ProErr);
910 //calculate integrated flow for POI selection
912 Double_t dYieldPt = fHistPtPOI->GetBinContent(b);
913 dVPOI += dv2pro*dYieldPt;
915 dErrVPOI += dYieldPt*dYieldPt*dv2ProErr*dv2ProErr;
916 } else { cout<<"fHistPtPOI is NULL"<<endl; }
920 dVPOI /= dSumPOI; //the pt distribution should be normalised
921 dErrVPOI /= (dSumPOI*dSumPOI);
922 dErrVPOI = TMath::Sqrt(dErrVPOI);
924 fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrVPOI);
925 cout<<"dV(POI) = "<<dVPOI<<" +- "<<dErrVPOI<<endl;
928 cout<<"*************************************"<<endl;
929 cout<<"*************************************"<<endl;
931 //cout<<".....finished"<<endl;
935 //--------------------------------------------------------------------
936 Double_t AliFlowAnalysisWithScalarProduct::CalculateStatisticalError(Int_t b, Double_t aStatErrorQaQb, TProfile* pHistProUQ, TProfile* pHistProUQQaQb, TH1D** pHistSumOfWeights) {
937 //calculate the statistical error for differential flow for bin b
938 Double_t duQproSpread = pHistProUQ->GetBinError(b);
939 Double_t sumOfMq = pHistSumOfWeights[0]->GetBinContent(b);
940 Double_t sumOfMqSquared = pHistSumOfWeights[1]->GetBinContent(b);
941 Double_t dQaQb = fHistProQaQbNorm->GetBinContent(1);
942 Double_t dTerm1 = 0.;
943 Double_t dTerm2 = 0.;
945 dTerm1 = (pow(sumOfMqSquared,0.5)/sumOfMq);
947 if(1.-pow(dTerm1,2.)>0.) {
948 dTerm2 = 1./pow(1.-pow(dTerm1,2.),0.5);
950 Double_t duQproErr = dTerm1*duQproSpread*dTerm2;
952 Double_t dTerm1Cov = pHistSumOfWeights[2]->GetBinContent(b);
953 Double_t dTerm2Cov = fHistSumOfLinearWeights->GetBinContent(1);
954 Double_t dTerm3Cov = sumOfMq;
955 Double_t dWeightedCovariance = 0.;
956 if(dTerm2Cov*dTerm3Cov>0.) {
957 Double_t dDenominator = 1.-dTerm1Cov/(dTerm2Cov*dTerm3Cov);
958 Double_t dPrefactor = dTerm1Cov/(dTerm2Cov*dTerm3Cov);
959 if(dDenominator!=0) {
960 Double_t dCovariance = (pHistProUQQaQb->GetBinContent(b)-dQaQb*pHistProUQ->GetBinContent(b))/dDenominator;
961 dWeightedCovariance = dCovariance*dPrefactor;
965 Double_t dv2ProErr = 0.; // final statitical error:
967 Double_t dv2ProErrorSquared = (1./4.)*pow(dQaQb,-3.)*
968 (pow(pHistProUQ->GetBinContent(b),2.)*pow(aStatErrorQaQb,2.)
969 + 4.*pow(dQaQb,2.)*pow(duQproErr,2.)
970 - 4.*dQaQb*pHistProUQ->GetBinContent(b)*dWeightedCovariance);
971 if(dv2ProErrorSquared>0.) dv2ProErr = pow(dv2ProErrorSquared,0.5);
978 //--------------------------------------------------------------------
980 void AliFlowAnalysisWithScalarProduct::StoreFlags()
982 // Store all boolean flags needed in Finish() in profile fHistProFlags.
984 // Apply correction for non-uniform acceptance or not:
985 fHistProFlags->Fill(0.5,fApplyCorrectionForNUA);
989 //--------------------------------------------------------------------
991 void AliFlowAnalysisWithScalarProduct::AccessFlags()
993 // Access all boolean flags needed in Finish() from profile fHistProFlags.
995 // Apply correction for non-uniform acceptance or not:
996 fApplyCorrectionForNUA = (Bool_t) fHistProFlags->GetBinContent(1);
1000 //--------------------------------------------------------------------