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 **************************************************************************/
20 //#define AliFlowAnalysisWithLYZEventPlane_cxx
22 #include "Riostream.h" //needed as include
23 #include "TComplex.h" //needed as include
24 #include "TProfile.h" //needed as include
31 #include "AliFlowLYZConstants.h" //needed as include
32 #include "AliFlowCommonConstants.h" //needed as include
33 #include "AliFlowEventSimple.h"
34 #include "AliFlowTrackSimple.h"
35 #include "AliFlowCommonHist.h"
36 #include "AliFlowCommonHistResults.h"
37 #include "AliFlowLYZEventPlane.h"
38 #include "AliFlowAnalysisWithLYZEventPlane.h"
42 // AliFlowAnalysisWithLYZEventPlane:
44 // Class to do flow analysis with the event plane from the LYZ method
46 // author: N. van der Kolk (kolk@nikhef.nl)
49 ClassImp(AliFlowAnalysisWithLYZEventPlane)
51 //-----------------------------------------------------------------------
53 AliFlowAnalysisWithLYZEventPlane::AliFlowAnalysisWithLYZEventPlane():
56 fSecondReDtheta(NULL),
57 fSecondImDtheta(NULL),
60 fHistProVetaPOI(NULL),
65 fHistQsumforChi(NULL),
68 fHistDeltaPhihere(NULL),
74 fCommonHistsRes(NULL),
81 fQsum = new TVector2(); // flow vector sum
83 fHistList = new TList();
84 fSecondRunList = new TList();
89 //-----------------------------------------------------------------------
92 AliFlowAnalysisWithLYZEventPlane::~AliFlowAnalysisWithLYZEventPlane()
97 delete fSecondRunList;
101 //-----------------------------------------------------------------------
103 void AliFlowAnalysisWithLYZEventPlane::WriteHistograms(TString* outputFileName)
105 //store the final results in output .root file
107 TFile *output = new TFile(outputFileName->Data(),"RECREATE");
108 //output->WriteObject(fHistList, "cobjLYZEP","SingleKey");
109 fHistList->SetName("cobjLYZEP");
110 fHistList->SetOwner(kTRUE);
111 fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
115 //-----------------------------------------------------------------------
117 void AliFlowAnalysisWithLYZEventPlane::WriteHistograms(TString outputFileName)
119 //store the final results in output .root file
121 TFile *output = new TFile(outputFileName.Data(),"RECREATE");
122 //output->WriteObject(fHistList, "cobjLYZEP","SingleKey");
123 fHistList->SetName("cobjLYZEP");
124 fHistList->SetOwner(kTRUE);
125 fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
129 //-----------------------------------------------------------------------
131 void AliFlowAnalysisWithLYZEventPlane::WriteHistograms(TDirectoryFile *outputFileName)
133 //store the final results in output .root file
134 fHistList->SetName("cobjLYZEP");
135 fHistList->SetOwner(kTRUE);
136 outputFileName->Add(fHistList);
137 outputFileName->Write(outputFileName->GetName(), TObject::kSingleKey);
140 //-----------------------------------------------------------------------
142 void AliFlowAnalysisWithLYZEventPlane::Init() {
144 //Initialise all histograms
145 cout<<"---Analysis with Lee Yang Zeros Event Plane Method---"<<endl;
147 //save old value and prevent histograms from being added to directory
148 //to avoid name clashes in case multiple analaysis objects are used
150 Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
151 TH1::AddDirectory(kFALSE);
154 if (fSecondRunList) {
155 fSecondReDtheta = (TProfile*)fSecondRunList->FindObject("Second_FlowPro_ReDtheta_LYZSUM");
156 fHistList->Add(fSecondReDtheta);
158 fSecondImDtheta = (TProfile*)fSecondRunList->FindObject("Second_FlowPro_ImDtheta_LYZSUM");
159 fHistList->Add(fSecondImDtheta);
161 fFirstr0theta = (TProfile*)fSecondRunList->FindObject("First_FlowPro_r0theta_LYZSUM");
162 fHistList->Add(fFirstr0theta);
165 if (!fSecondReDtheta) {cout<<"fSecondReDtheta is NULL!"<<endl; }
166 if (!fSecondImDtheta) {cout<<"fSecondImDtheta is NULL!"<<endl; }
167 if (!fFirstr0theta) {cout<<"fFirstr0theta is NULL!"<<endl; }
171 fCommonHists = new AliFlowCommonHist("AliFlowCommonHistLYZEP");
172 fHistList->Add(fCommonHists);
174 fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsLYZEP");
175 fHistList->Add(fCommonHistsRes);
177 Int_t iNbinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
178 Int_t iNbinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
179 Double_t dPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();
180 Double_t dPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
181 Double_t dEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();
182 Double_t dEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
184 fHistProVetaRP = new TProfile("FlowPro_VetaRP_LYZEP","FlowPro_VetaRP_LYZEP",iNbinsEta,dEtaMin,dEtaMax);
185 fHistProVetaRP->SetXTitle("rapidity");
186 fHistProVetaRP->SetYTitle("v_{2}(#eta) for RP selection");
187 fHistList->Add(fHistProVetaRP);
189 fHistProVetaPOI = new TProfile("FlowPro_VetaPOI_LYZEP","FlowPro_VetaPOI_LYZEP",iNbinsEta,dEtaMin,dEtaMax);
190 fHistProVetaPOI->SetXTitle("rapidity");
191 fHistProVetaPOI->SetYTitle("v_{2}(#eta) for POI selection");
192 fHistList->Add(fHistProVetaPOI);
194 fHistProVPtRP = new TProfile("FlowPro_VPtRP_LYZEP","FlowPro_VPtRP_LYZEP",iNbinsPt,dPtMin,dPtMax);
195 fHistProVPtRP->SetXTitle("Pt");
196 fHistProVPtRP->SetYTitle("v_{2}(p_{T}) for RP selection");
197 fHistList->Add(fHistProVPtRP);
199 fHistProVPtPOI = new TProfile("FlowPro_VPtPOI_LYZEP","FlowPro_VPtPOI_LYZEP",iNbinsPt,dPtMin,dPtMax);
200 fHistProVPtPOI->SetXTitle("p_{T}");
201 fHistProVPtPOI->SetYTitle("v_{2}(p_{T}) for POI selection");
202 fHistList->Add(fHistProVPtPOI);
204 fHistProWr = new TProfile("FlowPro_Wr_LYZEP","FlowPro_Wr_LYZEP",100,0.,0.25);
205 fHistProWr->SetXTitle("Q");
206 fHistProWr->SetYTitle("Wr");
207 fHistList->Add(fHistProWr);
209 fHistQsumforChi = new TH1F("Flow_QsumforChi_LYZEP","Flow_QsumforChi_LYZEP",3,-1.,2.);
210 fHistQsumforChi->SetXTitle("Qsum.X , Qsum.Y, Q2sum");
211 fHistQsumforChi->SetYTitle("value");
212 fHistList->Add(fHistQsumforChi);
214 fHistDeltaPhi = new TH1F("Flow_DeltaPhi_LYZEP","Flow_DeltaPhi_LYZEP",100,0.,3.14);
215 fHistDeltaPhi->SetXTitle("DeltaPhi");
216 fHistDeltaPhi->SetYTitle("Counts");
217 fHistList->Add(fHistDeltaPhi);
219 fHistPhiLYZ = new TH1F("Flow_PhiLYZ_LYZEP","Flow_PhiLYZ_LYZEP",100,0.,3.14);
220 fHistPhiLYZ->SetXTitle("Phi from LYZ");
221 fHistPhiLYZ->SetYTitle("Counts");
222 fHistList->Add(fHistPhiLYZ);
224 fHistPhiEP = new TH1F("Flow_PhiEP_LYZEP","Flow_PhiEP_LYZEP",100,0.,3.14);
225 fHistPhiEP->SetXTitle("Phi from EP");
226 fHistPhiEP->SetYTitle("Counts");
227 fHistList->Add(fHistPhiEP);
229 fEventNumber = 0; //set number of events to zero
232 TH1::AddDirectory(oldHistAddStatus);
235 //-----------------------------------------------------------------------
237 void AliFlowAnalysisWithLYZEventPlane::Make(AliFlowEventSimple* anEvent, AliFlowLYZEventPlane* aLYZEP) {
239 //Get the event plane and weight for each event
242 //fill control histograms
243 fCommonHists->FillControlHistograms(anEvent);
245 //get the Q vector from the FlowEvent
246 AliFlowVector vQ = anEvent->GetQ();
247 if (vQ.X()== 0. && vQ.Y()== 0. ) { cout<<"Q vector is NULL!"<<endl; }
248 //Weight with the multiplicity
251 if (vQ.GetMult()!=0.) {
252 dQX = vQ.X()/vQ.GetMult();
253 dQY = vQ.Y()/vQ.GetMult();
254 } else {cerr<<"vQ.GetMult() is zero!"<<endl; }
256 //cout<<"vQ("<<dQX<<","<<dQY<<")"<<endl;
258 //for chi calculation:
260 fHistQsumforChi->SetBinContent(1,fQsum->X());
261 fHistQsumforChi->SetBinContent(2,fQsum->Y());
263 fHistQsumforChi->SetBinContent(3,fQ2sum);
264 //cout<<"fQ2sum = "<<fQ2sum<<endl;
266 //call AliFlowLYZEventPlane::CalculateRPandW() here!
267 aLYZEP->CalculateRPandW(vQ);
269 Double_t dWR = aLYZEP->GetWR();
270 Double_t dRP = aLYZEP->GetPsi();
272 //fHistProWr->Fill(vQ.Mod(),dWR); //this weight is always positive
273 fHistPhiLYZ->Fill(dRP);
275 //plot difference between event plane from EP-method and LYZ-method
276 Double_t dRPEP = vQ.Phi()/2; //gives distribution from (0 to pi)
277 //Double_t dRPEP = 0.5*TMath::ATan2(vQ.Y(),vQ.X()); //gives distribution from (-pi/2 to pi/2)
278 //cout<<"dRPEP = "<< dRPEP <<endl;
279 fHistPhiEP->Fill(dRPEP);
281 Double_t dDeltaPhi = dRPEP - dRP;
282 if (dDeltaPhi < 0.) { dDeltaPhi += TMath::Pi(); } //to shift distribution from (-pi/2 to pi/2) to (0 to pi)
283 //cout<<"dDeltaPhi = "<<dDeltaPhi<<endl;
284 fHistDeltaPhi->Fill(dDeltaPhi);
287 Double_t dLow = TMath::Pi()/4.;
288 Double_t dHigh = 3.*(TMath::Pi()/4.);
289 if ((dDeltaPhi > dLow) && (dDeltaPhi < dHigh)){
290 dRP -= (TMath::Pi()/2);
292 cerr<<"*** dRP modified ***"<<endl;
294 fHistProWr->Fill(vQ.Mod(),dWR); //corrected weight
296 //calculate flow for RP and POI selections
297 //loop over the tracks of the event
298 Int_t iNumberOfTracks = anEvent->NumberOfTracks();
299 for (Int_t i=0;i<iNumberOfTracks;i++)
301 AliFlowTrackSimple* pTrack = anEvent->GetTrack(i) ;
303 Double_t dPhi = pTrack->Phi();
304 //if (dPhi<0.) fPhi+=2*TMath::Pi();
305 Double_t dPt = pTrack->Pt();
306 Double_t dEta = pTrack->Eta();
308 Double_t dv2 = dWR * TMath::Cos(2*(dPhi-dRP));
309 if (pTrack->InRPSelection()) {
310 //fill histograms for RP selection
311 fHistProVetaRP -> Fill(dEta,dv2);
312 fHistProVPtRP -> Fill(dPt,dv2);
314 if (pTrack->InPOISelection()) {
315 //fill histograms for POI selection
316 fHistProVetaPOI -> Fill(dEta,dv2);
317 fHistProVPtPOI -> Fill(dPt,dv2);
323 // cout<<"@@@@@ "<<fEventNumber<<" events processed"<<endl;
327 //--------------------------------------------------------------------
328 void AliFlowAnalysisWithLYZEventPlane::GetOutputHistograms(TList *outputListHistos){
329 //get pointers to all output histograms (called before Finish())
330 if (outputListHistos) {
331 //Get the common histograms from the output list
332 AliFlowCommonHist *pCommonHist = dynamic_cast<AliFlowCommonHist*>
333 (outputListHistos->FindObject("AliFlowCommonHistLYZEP"));
334 AliFlowCommonHistResults *pCommonHistResults = dynamic_cast<AliFlowCommonHistResults*>
335 (outputListHistos->FindObject("AliFlowCommonHistResultsLYZEP"));
337 TProfile* pHistProR0theta = dynamic_cast<TProfile*>
338 (outputListHistos->FindObject("First_FlowPro_r0theta_LYZSUM"));
340 TProfile* pHistProVetaRP = dynamic_cast<TProfile*>
341 (outputListHistos->FindObject("FlowPro_VetaRP_LYZEP"));
342 TProfile* pHistProVetaPOI = dynamic_cast<TProfile*>
343 (outputListHistos->FindObject("FlowPro_VetaPOI_LYZEP"));
344 TProfile* pHistProVPtRP = dynamic_cast<TProfile*>
345 (outputListHistos->FindObject("FlowPro_VPtRP_LYZEP"));
346 TProfile* pHistProVPtPOI = dynamic_cast<TProfile*>
347 (outputListHistos->FindObject("FlowPro_VPtPOI_LYZEP"));
349 TH1F* pHistQsumforChi = dynamic_cast<TH1F*>
350 (outputListHistos->FindObject("Flow_QsumforChi_LYZEP"));
352 if (pCommonHist && pCommonHistResults && pHistProR0theta &&
353 pHistProVetaRP && pHistProVetaPOI && pHistProVPtRP &&
354 pHistProVPtPOI && pHistQsumforChi ) {
355 this -> SetCommonHists(pCommonHist);
356 this -> SetCommonHistsRes(pCommonHistResults);
357 this -> SetFirstr0theta(pHistProR0theta);
358 this -> SetHistProVetaRP(pHistProVetaRP);
359 this -> SetHistProVetaPOI(pHistProVetaPOI);
360 this -> SetHistProVPtRP(pHistProVPtRP);
361 this -> SetHistProVPtPOI(pHistProVPtPOI);
362 this -> SetHistQsumforChi(pHistQsumforChi);
365 cout<<"WARNING: Histograms needed to run Finish() are not accessable!"<<endl;
369 //--------------------------------------------------------------------
370 void AliFlowAnalysisWithLYZEventPlane::Finish() {
372 //plot histograms etc.
373 cout<<"AliFlowAnalysisWithLYZEventPlane::Finish()"<<endl;
376 Double_t dJ01 = 2.405;
377 Int_t iNtheta = AliFlowLYZConstants::GetMaster()->GetNtheta();
378 Int_t iNbinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
379 Int_t iNbinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
380 //set the event number
382 SetEventNumber((int)fCommonHists->GetHistMultOrig()->GetEntries());
383 //cout<<"number of events processed is "<<fEventNumber<<endl;
386 //set the sum of Q vectors
387 fQsum->Set(fHistQsumforChi->GetBinContent(1),fHistQsumforChi->GetBinContent(2));
388 SetQ2sum(fHistQsumforChi->GetBinContent(3));
390 //calculate dV the mean of dVtheta
391 Double_t dVtheta = 0;
393 for (Int_t theta=0;theta<iNtheta;theta++) {
394 Double_t dR0 = fFirstr0theta->GetBinContent(theta+1);
395 if (dR0!=0.) { dVtheta = dJ01/dR0 ;}
401 Double_t dSigma2 = 0;
403 if (fEventNumber!=0) {
404 *fQsum /= fEventNumber;
405 //cerr<<"fQsum->X() = "<<fQsum->X()<<endl;
406 //cerr<<"fQsum->Y() = "<<fQsum->Y()<<endl;
407 fQ2sum /= fEventNumber;
408 //cerr<<"fEventNumber = "<<fEventNumber<<endl;
409 //cerr<<"fQ2sum = "<<fQ2sum<<endl;
410 dSigma2 = fQ2sum - TMath::Power(fQsum->X(),2.) - TMath::Power(fQsum->Y(),2.) - TMath::Power(dV,2.); //BP eq. 62
411 //cerr<<"dSigma2"<<dSigma2<<endl;
412 if (dSigma2>0) dChi = dV/TMath::Sqrt(dSigma2);
414 fCommonHistsRes->FillChiRP(dChi);
416 // recalculate statistical errors on integrated flow
417 //combining 5 theta angles to 1 relative error BP eq. 89
418 Double_t dRelErr2comb = 0.;
419 Int_t iEvts = fEventNumber;
421 for (Int_t theta=0;theta<iNtheta;theta++){
422 Double_t dTheta = ((double)theta/iNtheta)*TMath::Pi();
423 Double_t dApluscomb = TMath::Exp((dJ01*dJ01)/(2*dChi*dChi)*
425 Double_t dAmincomb = TMath::Exp(-(dJ01*dJ01)/(2*dChi*dChi)*
427 dRelErr2comb += (1/(2*iEvts*(dJ01*dJ01)*TMath::BesselJ1(dJ01)*
428 TMath::BesselJ1(dJ01)))*
429 (dApluscomb*TMath::BesselJ0(2*dJ01*TMath::Sin(dTheta/2)) +
430 dAmincomb*TMath::BesselJ0(2*dJ01*TMath::Cos(dTheta/2)));
432 dRelErr2comb /= iNtheta;
434 Double_t dRelErrcomb = TMath::Sqrt(dRelErr2comb);
435 Double_t dVErr = dV*dRelErrcomb ;
436 fCommonHistsRes->FillIntegratedFlow(dV, dVErr);
438 cout<<"*************************************"<<endl;
439 cout<<"*************************************"<<endl;
440 cout<<" Integrated flow from "<<endl;
441 cout<<" Lee-Yang Zeroes Event Plane "<<endl;
443 cout<<"dChi = "<<dChi<<endl;
444 cout<<"dV = "<<dV<<" +- "<<dVErr<<endl;
449 //copy content of profile into TH1D, add error and fill the AliFlowCommonHistResults
451 //v as a function of eta for RP selection
452 for(Int_t b=0;b<iNbinsEta;b++) {
453 Double_t dv2pro = fHistProVetaRP->GetBinContent(b);
454 Double_t dNprime = fCommonHists->GetEntriesInEtaBinRP(b);
455 Double_t dErrdifcomb = 0.; //set error to zero
456 Double_t dErr2difcomb = 0.; //set error to zero
459 for (Int_t theta=0;theta<iNtheta;theta++) {
460 Double_t dTheta = ((double)theta/iNtheta)*TMath::Pi();
461 Double_t dApluscomb = TMath::Exp((dJ01*dJ01)/(2*dChi*dChi)*
463 Double_t dAmincomb = TMath::Exp(-(dJ01*dJ01)/(2*dChi*dChi)*
465 dErr2difcomb += (TMath::Cos(dTheta)/(4*dNprime*TMath::BesselJ1(dJ01)*
466 TMath::BesselJ1(dJ01)))*
467 ((dApluscomb*TMath::BesselJ0(2*dJ01*TMath::Sin(dTheta/2))) -
468 (dAmincomb*TMath::BesselJ0(2*dJ01*TMath::Cos(dTheta/2))));
472 if (dErr2difcomb!=0.) {
473 dErr2difcomb /= iNtheta;
474 dErrdifcomb = TMath::Sqrt(dErr2difcomb);
476 else {dErrdifcomb = 0.;}
478 fCommonHistsRes->FillDifferentialFlowEtaRP(b, dv2pro, dErrdifcomb);
482 //v as a function of eta for POI selection
483 for(Int_t b=0;b<iNbinsEta;b++) {
484 Double_t dv2pro = fHistProVetaPOI->GetBinContent(b);
485 Double_t dNprime = fCommonHists->GetEntriesInEtaBinPOI(b);
486 Double_t dErrdifcomb = 0.; //set error to zero
487 Double_t dErr2difcomb = 0.; //set error to zero
490 for (Int_t theta=0;theta<iNtheta;theta++) {
491 Double_t dTheta = ((double)theta/iNtheta)*TMath::Pi();
492 Double_t dApluscomb = TMath::Exp((dJ01*dJ01)/(2*dChi*dChi)*
494 Double_t dAmincomb = TMath::Exp(-(dJ01*dJ01)/(2*dChi*dChi)*
496 dErr2difcomb += (TMath::Cos(dTheta)/(4*dNprime*TMath::BesselJ1(dJ01)*
497 TMath::BesselJ1(dJ01)))*
498 ((dApluscomb*TMath::BesselJ0(2*dJ01*TMath::Sin(dTheta/2))) -
499 (dAmincomb*TMath::BesselJ0(2*dJ01*TMath::Cos(dTheta/2))));
503 if (dErr2difcomb!=0.) {
504 dErr2difcomb /= iNtheta;
505 dErrdifcomb = TMath::Sqrt(dErr2difcomb);
507 else {dErrdifcomb = 0.;}
509 fCommonHistsRes->FillDifferentialFlowEtaPOI(b, dv2pro, dErrdifcomb);
512 //v as a function of Pt for RP selection
513 TH1F* fHistPtRP = fCommonHists->GetHistPtRP(); //for calculating integrated flow
518 for(Int_t b=0;b<iNbinsPt;b++) {
519 Double_t dv2pro = fHistProVPtRP->GetBinContent(b);
520 Double_t dNprime = fCommonHists->GetEntriesInPtBinRP(b);
521 Double_t dErrdifcomb = 0.; //set error to zero
522 Double_t dErr2difcomb = 0.; //set error to zero
525 for (Int_t theta=0;theta<iNtheta;theta++) {
526 Double_t dTheta = ((double)theta/iNtheta)*TMath::Pi();
527 Double_t dApluscomb = TMath::Exp((dJ01*dJ01)/(2*dChi*dChi)*
529 Double_t dAmincomb = TMath::Exp(-(dJ01*dJ01)/(2*dChi*dChi)*
531 dErr2difcomb += (TMath::Cos(dTheta)/(4*dNprime*TMath::BesselJ1(dJ01)*
532 TMath::BesselJ1(dJ01)))*
533 ((dApluscomb*TMath::BesselJ0(2*dJ01*TMath::Sin(dTheta/2))) -
534 (dAmincomb*TMath::BesselJ0(2*dJ01*TMath::Cos(dTheta/2))));
538 if (dErr2difcomb!=0.) {
539 dErr2difcomb /= iNtheta;
540 dErrdifcomb = TMath::Sqrt(dErr2difcomb);
541 //cerr<<"dErrdifcomb = "<<dErrdifcomb<<endl;
543 else {dErrdifcomb = 0.;}
546 fCommonHistsRes->FillDifferentialFlowPtRP(b, dv2pro, dErrdifcomb);
547 //calculate integrated flow for RP selection
549 Double_t dYieldPt = fHistPtRP->GetBinContent(b);
550 dVRP += dv2pro*dYieldPt;
552 dErrV += dYieldPt*dYieldPt*dErrdifcomb*dErrdifcomb;
553 } else { cout<<"fHistPtRP is NULL"<<endl; }
558 dVRP /= dSum; //the pt distribution should be normalised
559 dErrV /= (dSum*dSum);
560 dErrV = TMath::Sqrt(dErrV);
562 fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrV);
564 cout<<"dV(RP) = "<<dVRP<<" +- "<<dErrV<<endl;
568 //v as a function of Pt for POI selection
569 TH1F* fHistPtPOI = fCommonHists->GetHistPtPOI(); //for calculating integrated flow
574 for(Int_t b=0;b<iNbinsPt;b++) {
575 Double_t dv2pro = fHistProVPtPOI->GetBinContent(b);
576 Double_t dNprime = fCommonHists->GetEntriesInPtBinPOI(b);
578 //cerr<<"dNprime = "<<dNprime<<endl;
579 //Int_t iNprime = TMath::Nint(fHistProVPtPOI->GetBinEntries(b));
580 //cerr<<"iNprime = "<<iNprime<<endl;
582 Double_t dErrdifcomb = 0.; //set error to zero
583 Double_t dErr2difcomb = 0.; //set error to zero
586 for (Int_t theta=0;theta<iNtheta;theta++) {
587 Double_t dTheta = ((double)theta/iNtheta)*TMath::Pi();
588 Double_t dApluscomb = TMath::Exp((dJ01*dJ01)/(2*dChi*dChi)*
590 Double_t dAmincomb = TMath::Exp(-(dJ01*dJ01)/(2*dChi*dChi)*
592 dErr2difcomb += (TMath::Cos(dTheta)/(4*dNprime*TMath::BesselJ1(dJ01)*
593 TMath::BesselJ1(dJ01)))*
594 ((dApluscomb*TMath::BesselJ0(2*dJ01*TMath::Sin(dTheta/2))) -
595 (dAmincomb*TMath::BesselJ0(2*dJ01*TMath::Cos(dTheta/2))));
599 if (dErr2difcomb!=0.) {
600 dErr2difcomb /= iNtheta;
601 dErrdifcomb = TMath::Sqrt(dErr2difcomb);
602 //cerr<<"dErrdifcomb = "<<dErrdifcomb<<endl;
604 else {dErrdifcomb = 0.;}
607 fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2pro, dErrdifcomb);
609 //calculate integrated flow for POI selection
611 Double_t dYieldPt = fHistPtPOI->GetBinContent(b);
612 dVPOI += dv2pro*dYieldPt;
614 dErrV += dYieldPt*dYieldPt*dErrdifcomb*dErrdifcomb;
615 } else { cout<<"fHistPtPOI is NULL"<<endl; }
619 dVPOI /= dSum; //the pt distribution should be normalised
620 dErrV /= (dSum*dSum);
621 dErrV = TMath::Sqrt(dErrV);
623 fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrV);
625 cout<<"dV(POI) = "<<dVPOI<<" +- "<<dErrV<<endl;
627 cout<<"*************************************"<<endl;
628 cout<<"*************************************"<<endl;
630 //cout<<".....finished"<<endl;