/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ #include "Riostream.h" //needed as include #include "TObject.h" //needed as include #include "AliFlowCommonConstants.h" //needed as include #include "AliFlowLYZConstants.h" //needed as include #include "AliFlowAnalysisWithLeeYangZeros.h" #include "AliFlowLYZHist1.h" //needed as include #include "AliFlowLYZHist2.h" //needed as include #include "AliFlowCommonHist.h" //needed as include #include "AliFlowCommonHistResults.h" //needed as include #include "AliFlowEventSimple.h" //needed as include #include "AliFlowTrackSimple.h" //needed as include class AliFlowVector; #include "TMath.h" //needed as include #include "TFile.h" //needed as include #include "TList.h" class TComplex; class TProfile; class TH1F; class TH1D; //Description: Maker to analyze Flow using the LeeYangZeros method // Equation numbers are from Big Paper (BP): Nucl. Phys. A 727, 373 (2003) // Practical Guide (PG): J. Phys. G: Nucl. Part. Phys. 30, S1213 (2004) // Adapted from StFlowLeeYangZerosMaker.cxx // by Markus Oldenberg and Art Poskanzer, LBNL // with advice from Jean-Yves Ollitrault and Nicolas Borghini // //Author: Naomi van der Kolk (kolk@nikhef.nl) ClassImp(AliFlowAnalysisWithLeeYangZeros) //----------------------------------------------------------------------- AliFlowAnalysisWithLeeYangZeros::AliFlowAnalysisWithLeeYangZeros(): fQsum(NULL), fQ2sum(0), fEventNumber(0), fFirstRun(kTRUE), fUseSum(kTRUE), fDoubleLoop(kFALSE), fDebug(kFALSE), fHistList(NULL), fFirstRunList(NULL), fHistProVtheta(NULL), fHistProVetaRP(NULL), fHistProVetaPOI(NULL), fHistProVPtRP(NULL), fHistProVPtPOI(NULL), fHistProR0theta(NULL), fHistProReDenom(NULL), fHistProImDenom(NULL), fHistProReDtheta(NULL), fHistProImDtheta(NULL), fHistQsumforChi(NULL), fCommonHists(NULL), fCommonHistsRes(NULL) { //default constructor if (fDebug) cout<<"****AliFlowAnalysisWithLeeYangZeros::AliFlowAnalysisWithLeeYangZeros default constructor****"<Data(),"RECREATE"); if (GetFirstRun()) { output->WriteObject(fHistList, "cobjLYZ1","SingleKey"); } else { output->WriteObject(fHistList, "cobjLYZ2","SingleKey"); } delete output; } //----------------------------------------------------------------------- void AliFlowAnalysisWithLeeYangZeros::WriteHistograms(TString outputFileName) { //store the final results in output .root file TFile *output = new TFile(outputFileName.Data(),"RECREATE"); if (GetFirstRun()) { output->WriteObject(fHistList, "cobjLYZ1","SingleKey"); } else { output->WriteObject(fHistList, "cobjLYZ2","SingleKey"); } delete output; } //----------------------------------------------------------------------- Bool_t AliFlowAnalysisWithLeeYangZeros::Init() { //init method if (fDebug) cout<<"****AliFlowAnalysisWithLeeYangZeros::Init()****"<Add(fCommonHists); fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsLYZ1"); fHistList->Add(fCommonHistsRes); } else { fCommonHists = new AliFlowCommonHist("AliFlowCommonHistLYZ2"); fHistList->Add(fCommonHists); fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsLYZ2"); fHistList->Add(fCommonHistsRes); } fHistQsumforChi = new TH1F("Flow_QsumforChi_LYZ","Flow_QsumforChi_LYZ",3,-1.,2.); fHistQsumforChi->SetXTitle("Qsum.X , Qsum.Y, Q2sum"); fHistQsumforChi->SetYTitle("value"); fHistList->Add(fHistQsumforChi); //for first loop over events if (fFirstRun){ fHistProR0theta = new TProfile("First_FlowPro_r0theta_LYZ","First_FlowPro_r0theta_LYZ",iNtheta,-0.5,iNtheta-0.5); fHistProR0theta->SetXTitle("#theta"); fHistProR0theta->SetYTitle("r_{0}^{#theta}"); fHistList->Add(fHistProR0theta); fHistProVtheta = new TProfile("First_FlowPro_Vtheta_LYZ","First_FlowPro_Vtheta_LYZ",iNtheta,-0.5,iNtheta-0.5); fHistProVtheta->SetXTitle("#theta"); fHistProVtheta->SetYTitle("V_{n}^{#theta}"); fHistList->Add(fHistProVtheta); //class AliFlowLYZHist1 defines the histograms: fHistProGtheta, fHistProReGtheta, fHistProImGtheta, fHistProR0theta for (Int_t theta=0;thetaAdd(fHist1[theta]); } } //for second loop over events else { fHistProReDenom = new TProfile("Second_FlowPro_ReDenom_LYZ","Second_FlowPro_ReDenom_LYZ" , iNtheta, -0.5, iNtheta-0.5); fHistProReDenom->SetXTitle("#theta"); fHistProReDenom->SetYTitle("Re(Q^{#theta}e^{ir_{0}^{#theta}Q^{#theta}})"); fHistList->Add(fHistProReDenom); fHistProImDenom = new TProfile("Second_FlowPro_ImDenom_LYZ","Second_FlowPro_ImDenom_LYZ" , iNtheta, -0.5, iNtheta-0.5); fHistProImDenom->SetXTitle("#theta"); fHistProImDenom->SetYTitle("Im(Q^{#theta}e^{ir_{0}^{#theta}Q^{#theta}})"); fHistList->Add(fHistProImDenom); fHistProVetaRP = new TProfile("Second_FlowPro_VetaRP_LYZ","Second_FlowPro_VetaRP_LYZ",iNbinsEta,dEtaMin,dEtaMax); fHistProVetaRP->SetXTitle("rapidity"); fHistProVetaRP->SetYTitle("v_{2}(#eta) for RP selection"); fHistList->Add(fHistProVetaRP); fHistProVetaPOI = new TProfile("Second_FlowPro_VetaPOI_LYZ","Second_FlowPro_VetaPOI_LYZ",iNbinsEta,dEtaMin,dEtaMax); fHistProVetaPOI->SetXTitle("rapidity"); fHistProVetaPOI->SetYTitle("v_{2}(#eta) for POI selection"); fHistList->Add(fHistProVetaPOI); fHistProVPtRP = new TProfile("Second_FlowPro_VPtRP_LYZ","Second_FlowPro_VPtRP_LYZ",iNbinsPt,dPtMin,dPtMax); fHistProVPtRP->SetXTitle("Pt"); fHistProVPtRP->SetYTitle("v_{2}(p_{T}) for RP selection"); fHistList->Add(fHistProVPtRP); fHistProVPtPOI = new TProfile("Second_FlowPro_VPtPOI_LYZ","Second_FlowPro_VPtPOI_LYZ",iNbinsPt,dPtMin,dPtMax); fHistProVPtPOI->SetXTitle("p_{T}"); fHistProVPtPOI->SetYTitle("v_{2}(p_{T}) for POI selection"); fHistList->Add(fHistProVPtPOI); fHistProReDtheta = new TProfile("Second_FlowPro_ReDtheta_LYZ","Second_FlowPro_ReDtheta_LYZ",iNtheta, -0.5, iNtheta-0.5); fHistProReDtheta->SetXTitle("#theta"); fHistProReDtheta->SetYTitle("Re(D^{#theta})"); fHistList->Add(fHistProReDtheta); fHistProImDtheta = new TProfile("Second_FlowPro_ImDtheta_LYZ","Second_FlowPro_ImDtheta_LYZ",iNtheta, -0.5, iNtheta-0.5); fHistProImDtheta->SetXTitle("#theta"); fHistProImDtheta->SetYTitle("Im(D^{#theta})"); fHistList->Add(fHistProImDtheta); //class AliFlowLYZHist2 defines the histograms: for (Int_t theta=0;thetaAdd(fHist2RP[theta]); TString namePOI = "AliFlowLYZHist2POI_"; namePOI += theta; fHist2POI[theta]=new AliFlowLYZHist2(theta, "POI", namePOI); fHistList->Add(fHist2POI[theta]); } //read histogram fHistProR0theta from the first run list if (fFirstRunList) { fHistProR0theta = (TProfile*)fFirstRunList->FindObject("First_FlowPro_r0theta_LYZ"); if (!fHistProR0theta) {cout<<"fHistProR0theta has a NULL pointer!"<Add(fHistProR0theta); } else { cout<<"list is NULL pointer!"<FillControlHistograms(anEvent); FillFromFlowEvent(anEvent); } else { fCommonHists->FillControlHistograms(anEvent); SecondFillFromFlowEvent(anEvent); } } else { cout<<"##### FlowLeeYangZero: Stack pointer null"<GetHistMultOrig()->GetEntries()); //cout<<"number of events processed is "<Set(fHistQsumforChi->GetBinContent(1),fHistQsumforChi->GetBinContent(2)); SetQ2sum(fHistQsumforChi->GetBinContent(3)); if (fFirstRun){ Double_t dR0 = 0; Double_t dVtheta = 0; Double_t dv = 0; Double_t dV = 0; for (Int_t theta=0;thetaFillGtheta(); dR0 = fHist1[theta]->GetR0(); //calculate integrated flow if (dR0!=0) dVtheta = dJ01/dR0; //for estimating systematic error resulting from d0 Double_t dBinsize = (AliFlowLYZConstants::fgMax)/(AliFlowLYZConstants::kNbins); Double_t dVplus = dJ01/(dR0+dBinsize); Double_t dVmin = dJ01/(dR0-dBinsize); dv = dVtheta; Double_t dvplus = dVplus; Double_t dvmin = dVmin; if (fDebug) cout<<"dv = "<FillIntegratedFlow(dv2pro, dv2Err); if (fDebug) cout<<"****histograms filled****"<GetBinContent(theta+1); //histogram starts at bin 1 if (fDebug) cerr<<"dR0 = "< Vn^theta = j01/r0^theta if (!fHistProReDenom && !fHistProImDenom) { cout << "Hist pointer fDenom in file does not exist" <GetBinContent(theta+1); dImDenom = fHistProImDenom->GetBinContent(theta+1); cDenom(dReDenom,dImDenom); //for new method and use by others (only with the sum generating function): if (fUseSum) { cDtheta = dR0*cDenom/dJ01; fHistProReDtheta->Fill(theta,cDtheta.Re()); fHistProImDtheta->Fill(theta,cDtheta.Im()); } cDenom *= TComplex::Power(i, m-1); //cerr<<"TComplex::Power(i, m-1) = "<GetBinCenter(be); cNumerRP = fHist2RP[theta]->GetNumerEta(be); if (cNumerRP.Rho()==0) { if (fDebug) cerr<<"WARNING: modulus of cNumerRP is zero in Finish()"<Fill(dEtaRP,dVetaRP); } //loop over bins be //v as a function of eta for POI selection for (Int_t be=1;be<=iNbinsEta;be++) { dEtaPOI = fHist2POI[theta]->GetBinCenter(be); cNumerPOI = fHist2POI[theta]->GetNumerEta(be); if (cNumerPOI.Rho()==0) { if (fDebug) cerr<<"WARNING: modulus of cNumerPOI is zero in Finish()"<Fill(dEtaPOI,dVetaPOI); } //loop over bins be //v as a function of Pt for RP selection for (Int_t bp=1;bp<=iNbinsPt;bp++) { dPtRP = fHist2RP[theta]->GetBinCenterPt(bp); cNumerRP = fHist2RP[theta]->GetNumerPt(bp); if (cNumerRP.Rho()==0) { if (fDebug) cerr<<"WARNING: modulus of cNumerRP is zero"<Fill(dPtRP,dVPtRP); } //loop over bins bp //v as a function of Pt for POI selection for (Int_t bp=1;bp<=iNbinsPt;bp++) { dPtPOI = fHist2POI[theta]->GetBinCenterPt(bp); cNumerPOI = fHist2POI[theta]->GetNumerPt(bp); if (cNumerPOI.Rho()==0) { if (fDebug) cerr<<"WARNING: modulus of cNumerPOI is zero"<Fill(dPtPOI,dVPtPOI); } //loop over bins bp } } }//end of loop over theta //calculate the average of fVtheta = fV dV /= iNtheta; //sigma2 and chi (for statistical error calculations) Double_t dSigma2 = 0; Double_t dChi= 0; if (fEventNumber!=0) { *fQsum /= fEventNumber; //cerr<<"fQsum.X() = "<GetBinContent(b); Double_t dNprime = fCommonHists->GetEntriesInEtaBinRP(b); Double_t dErrdifcomb = 0.; //set error to zero Double_t dErr2difcomb = 0.; //set error to zero //calculate error if (dNprime!=0.) { for (Int_t theta=0;thetaFillDifferentialFlowEtaRP(b, dv2pro, dErrdifcomb); } //loop over bins b //v as a function of eta for POI selection for(Int_t b=0;bGetBinContent(b); Double_t dNprime = fCommonHists->GetEntriesInEtaBinPOI(b); Double_t dErrdifcomb = 0.; //set error to zero Double_t dErr2difcomb = 0.; //set error to zero //calculate error if (dNprime!=0.) { for (Int_t theta=0;thetaFillDifferentialFlowEtaPOI(b, dv2pro, dErrdifcomb); } //loop over bins b //v as a function of Pt for RP selection TH1F* fHistPtInt = fCommonHists->GetHistPtInt(); //for calculating integrated flow Double_t dVRP = 0.; Double_t dSum = 0.; Double_t dErrV =0.; for(Int_t b=0;bGetBinContent(b); Double_t dNprime = fCommonHists->GetEntriesInPtBinRP(b); Double_t dErrdifcomb = 0.; //set error to zero Double_t dErr2difcomb = 0.; //set error to zero //calculate error if (dNprime!=0.) { for (Int_t theta=0;thetaFillIntegratedFlowRP(dVRP,dErrV); //v as a function of Pt for POI selection TH1F* fHistPtDiff = fCommonHists->GetHistPtDiff(); //for calculating integrated flow Double_t dVPOI = 0.; dSum = 0.; dErrV =0.; for(Int_t b=0;bGetBinContent(b); Double_t dNprime = fCommonHists->GetEntriesInPtBinPOI(b); Double_t dErrdifcomb = 0.; //set error to zero Double_t dErr2difcomb = 0.; //set error to zero //calculate error if (dNprime!=0.) { for (Int_t theta=0;thetaFillIntegratedFlowPOI(dVPOI,dErrV); } //secondrun cout<<"----LYZ analysis finished....----"<GetQ(); //weight by the multiplicity Double_t dQX = 0; Double_t dQY = 0; if (vQ.GetMult() != 0) { dQX = vQ.X()/vQ.GetMult(); dQY = vQ.Y()/vQ.GetMult(); } vQ.Set(dQX,dQY); //for chi calculation: *fQsum += vQ; fHistQsumforChi->SetBinContent(1,fQsum->X()); fHistQsumforChi->SetBinContent(2,fQsum->Y()); fQ2sum += vQ.Mod2(); fHistQsumforChi->SetBinContent(3,fQ2sum); //cerr<<"fQ2sum = "<GetBinCenter(bin); //bincentre of bins in histogram //FIXED??? //if (theta == 0) cerr<<"cerr::dR = "< 100.) break; } fHist1[theta]->Fill(dR,cGtheta); //fill real and imaginary part of cGtheta } //loop over bins } //loop over theta return kTRUE; } //----------------------------------------------------------------------- Bool_t AliFlowAnalysisWithLeeYangZeros::SecondFillFromFlowEvent(AliFlowEventSimple* anEvent) { //for differential flow if (fDebug) cout<<"****AliFlowAnalysisWithLeeYangZeros::SecondFillFromFlowEvent()****"<GetQ(); //weight by the multiplicity Double_t dQX = 0.; Double_t dQY = 0.; if (vQ.GetMult() != 0) { dQX = vQ.X()/vQ.GetMult(); dQY = vQ.Y()/vQ.GetMult(); } vQ.Set(dQX,dQY); //for chi calculation: *fQsum += vQ; fHistQsumforChi->SetBinContent(1,fQsum->X()); fHistQsumforChi->SetBinContent(2,fQsum->Y()); fQ2sum += vQ.Mod2(); fHistQsumforChi->SetBinContent(3,fQ2sum); for (Int_t theta=0;thetaFill(dEta,dPt,cNumerRP); } } if (pTrack->UseForDifferentialFlow()) { //POI selection dCosTermPOI = cos(m*dOrder*(dPhi-dTheta)); cNumerPOI = dCosTermPOI*(TComplex::Exp(cExpo)); if (cNumerPOI.Rho()==0) {cerr<<"WARNING: modulus of cNumerPOI is zero in SecondFillFromFlowEvent"<Fill(dEta,dPt,cNumerPOI); } } else { cout << "fHist2 pointer mising" <Fill(theta,cDenom.Re()); //fill the real part of fDenom fHistProImDenom->Fill(theta,cDenom.Im()); //fill the imaginary part of fDenom } else { cout << "Pointers to cDenom mising" << endl; } }//end of loop over theta return kTRUE; } //----------------------------------------------------------------------- Double_t AliFlowAnalysisWithLeeYangZeros::GetQtheta(AliFlowVector aQ, Double_t aTheta) { //calculate Qtheta. Qtheta is the sum over all particles of cos(dOrder*(dPhi-dTheta)) BP eq. 3 if (fDebug) cout<<"****AliFlowAnalysisWithLeeYangZeros::GetQtheta()****"<NumberOfTracks(); for (Int_t i=0;iGetTrack(i) ; if (pTrack){ if (pTrack->UseForIntegratedFlow()) { Double_t dPhi = pTrack->Phi(); Double_t dGIm = aR * dWgt*cos(dOrder*(dPhi - aTheta)); TComplex cGi(1., dGIm); cG *= cGi; //product over all tracks } } else {cerr << "no particle pointer !!!"<NumberOfTracks(); Int_t iNtheta = AliFlowLYZConstants::kTheta; Double_t dTheta = ((double)theta/iNtheta)*TMath::Pi()/dOrder; //for the denominator (use all RP selected particles) for (Int_t i=0;iGetTrack(i) ; if (pTrack){ if (pTrack->UseForIntegratedFlow()) { Double_t dPhi = pTrack->Phi(); Double_t dCosTerm = dWgt*cos(dOrder*(dPhi - dTheta)); //GetGr0theta Double_t dGIm = aR0 * dCosTerm; TComplex cGi(1., dGIm); TComplex cCosTermComplex(1., aR0*dCosTerm); cG *= cGi; //product over all tracks //GetdGr0theta cdGr0 +=(dCosTerm / cCosTermComplex); //sum over all tracks } } //if particle else {cerr << "no particle!!!"<GetTrack(i) ; if (pTrack){ Double_t dEta = pTrack->Eta(); Double_t dPt = pTrack->Pt(); Double_t dPhi = pTrack->Phi(); Double_t dCosTerm = cos(dOrder*(dPhi-dTheta)); TComplex cCosTermComplex(1.,aR0*dCosTerm); //RP selection if (pTrack->UseForIntegratedFlow()) { TComplex cNumerRP = cG*dCosTerm/cCosTermComplex; //PG Eq. 9 fHist2RP[theta]->Fill(dEta,dPt,cNumerRP); } //POI selection if (pTrack->UseForDifferentialFlow()) { TComplex cNumerPOI = cG*dCosTerm/cCosTermComplex; //PG Eq. 9 fHist2POI[theta]->Fill(dEta,dPt,cNumerPOI); } } //if particle else {cerr << "no particle pointer!!!"<