/************************************************************************* * Copyright(c) 1998-2008, 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. * **************************************************************************/ /* $Log$ */ //#define AliFlowAnalysisWithLYZEventPlane_cxx #include "Riostream.h" //needed as include #include "TComplex.h" //needed as include #include "TProfile.h" //needed as include class TH1F; class TFile; class TList; class TVector2; #include "AliFlowLYZConstants.h" //needed as include #include "AliFlowCommonConstants.h" //needed as include #include "AliFlowEventSimple.h" #include "AliFlowTrackSimple.h" #include "AliFlowCommonHist.h" #include "AliFlowCommonHistResults.h" #include "AliFlowLYZEventPlane.h" #include "AliFlowAnalysisWithLYZEventPlane.h" class AliFlowVector; // AliFlowAnalysisWithLYZEventPlane: // // Class to do flow analysis with the event plane from the LYZ method // // author: N. van der Kolk (kolk@nikhef.nl) ClassImp(AliFlowAnalysisWithLYZEventPlane) //----------------------------------------------------------------------- AliFlowAnalysisWithLYZEventPlane::AliFlowAnalysisWithLYZEventPlane(): fHistList(NULL), fSecondRunList(NULL), fSecondReDtheta(NULL), fSecondImDtheta(NULL), fFirstr0theta(NULL), fHistProVetaRP(NULL), fHistProVetaPOI(NULL), fHistProVPtRP(NULL), fHistProVPtPOI(NULL), fHistProWr(NULL), fHistProWrCorr(NULL), fHistQsumforChi(NULL), fHistDeltaPhi(NULL), fHistDeltaPhi2(NULL), fHistDeltaPhihere(NULL), fHistPhiEP(NULL), fHistPhiEPhere(NULL), fHistPhiLYZ(NULL), fHistPhiLYZ2(NULL), fCommonHists(NULL), fCommonHistsRes(NULL), fEventNumber(0), fQsum(NULL), fQ2sum(0) { // Constructor. fQsum = new TVector2(); // flow vector sum fHistList = new TList(); fSecondRunList = new TList(); } //----------------------------------------------------------------------- AliFlowAnalysisWithLYZEventPlane::~AliFlowAnalysisWithLYZEventPlane() { //destructor delete fQsum; delete fHistList; delete fSecondRunList; } //----------------------------------------------------------------------- void AliFlowAnalysisWithLYZEventPlane::WriteHistograms(TString* outputFileName) { //store the final results in output .root file TFile *output = new TFile(outputFileName->Data(),"RECREATE"); //output->WriteObject(fHistList, "cobjLYZEP","SingleKey"); fHistList->SetName("cobjLYZEP"); fHistList->Write(fHistList->GetName(), TObject::kSingleKey); delete output; } //----------------------------------------------------------------------- void AliFlowAnalysisWithLYZEventPlane::WriteHistograms(TString outputFileName) { //store the final results in output .root file TFile *output = new TFile(outputFileName.Data(),"RECREATE"); //output->WriteObject(fHistList, "cobjLYZEP","SingleKey"); fHistList->SetName("cobjLYZEP"); fHistList->Write(fHistList->GetName(), TObject::kSingleKey); delete output; } //----------------------------------------------------------------------- void AliFlowAnalysisWithLYZEventPlane::Init() { //Initialise all histograms cout<<"---Analysis with Lee Yang Zeros Event Plane Method---"<FindObject("Second_FlowPro_ReDtheta_LYZ"); fHistList->Add(fSecondReDtheta); fSecondImDtheta = (TProfile*)fSecondRunList->FindObject("Second_FlowPro_ImDtheta_LYZ"); fHistList->Add(fSecondImDtheta); fFirstr0theta = (TProfile*)fSecondRunList->FindObject("First_FlowPro_r0theta_LYZ"); fHistList->Add(fFirstr0theta); //warnings if (!fSecondReDtheta) {cout<<"fSecondReDtheta is NULL!"<Add(fCommonHists); fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsLYZEP"); fHistList->Add(fCommonHistsRes); Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt(); Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta(); Double_t dPtMin = AliFlowCommonConstants::GetPtMin(); Double_t dPtMax = AliFlowCommonConstants::GetPtMax(); Double_t dEtaMin = AliFlowCommonConstants::GetEtaMin(); Double_t dEtaMax = AliFlowCommonConstants::GetEtaMax(); fHistProVetaRP = new TProfile("FlowPro_VetaRP_LYZEP","FlowPro_VetaRP_LYZEP",iNbinsEta,dEtaMin,dEtaMax); fHistProVetaRP->SetXTitle("rapidity"); fHistProVetaRP->SetYTitle("v_{2}(#eta) for RP selection"); fHistList->Add(fHistProVetaRP); fHistProVetaPOI = new TProfile("FlowPro_VetaPOI_LYZEP","FlowPro_VetaPOI_LYZEP",iNbinsEta,dEtaMin,dEtaMax); fHistProVetaPOI->SetXTitle("rapidity"); fHistProVetaPOI->SetYTitle("v_{2}(#eta) for POI selection"); fHistList->Add(fHistProVetaPOI); fHistProVPtRP = new TProfile("FlowPro_VPtRP_LYZEP","FlowPro_VPtRP_LYZEP",iNbinsPt,dPtMin,dPtMax); fHistProVPtRP->SetXTitle("Pt"); fHistProVPtRP->SetYTitle("v_{2}(p_{T}) for RP selection"); fHistList->Add(fHistProVPtRP); fHistProVPtPOI = new TProfile("FlowPro_VPtPOI_LYZEP","FlowPro_VPtPOI_LYZEP",iNbinsPt,dPtMin,dPtMax); fHistProVPtPOI->SetXTitle("p_{T}"); fHistProVPtPOI->SetYTitle("v_{2}(p_{T}) for POI selection"); fHistList->Add(fHistProVPtPOI); fHistProWr = new TProfile("FlowPro_Wr_LYZEP","FlowPro_Wr_LYZEP",100,0.,0.25); fHistProWr->SetXTitle("Q"); fHistProWr->SetYTitle("Wr"); fHistList->Add(fHistProWr); fHistQsumforChi = new TH1F("Flow_QsumforChi_LYZEP","Flow_QsumforChi_LYZEP",3,-1.,2.); fHistQsumforChi->SetXTitle("Qsum.X , Qsum.Y, Q2sum"); fHistQsumforChi->SetYTitle("value"); fHistList->Add(fHistQsumforChi); fHistDeltaPhi = new TH1F("Flow_DeltaPhi_LYZEP","Flow_DeltaPhi_LYZEP",100,0.,3.14); fHistDeltaPhi->SetXTitle("DeltaPhi"); fHistDeltaPhi->SetYTitle("Counts"); fHistList->Add(fHistDeltaPhi); fHistPhiLYZ = new TH1F("Flow_PhiLYZ_LYZEP","Flow_PhiLYZ_LYZEP",100,0.,3.14); fHistPhiLYZ->SetXTitle("Phi from LYZ"); fHistPhiLYZ->SetYTitle("Counts"); fHistList->Add(fHistPhiLYZ); fHistPhiEP = new TH1F("Flow_PhiEP_LYZEP","Flow_PhiEP_LYZEP",100,0.,3.14); fHistPhiEP->SetXTitle("Phi from EP"); fHistPhiEP->SetYTitle("Counts"); fHistList->Add(fHistPhiEP); fEventNumber = 0; //set number of events to zero } //----------------------------------------------------------------------- void AliFlowAnalysisWithLYZEventPlane::Make(AliFlowEventSimple* anEvent, AliFlowLYZEventPlane* aLYZEP) { //Get the event plane and weight for each event if (anEvent) { //fill control histograms fCommonHists->FillControlHistograms(anEvent); //get the Q vector from the FlowEvent AliFlowVector vQ = anEvent->GetQ(); if (vQ.X()== 0. && vQ.Y()== 0. ) { cout<<"Q vector is NULL!"<SetBinContent(1,fQsum->X()); fHistQsumforChi->SetBinContent(2,fQsum->Y()); fQ2sum += vQ.Mod2(); fHistQsumforChi->SetBinContent(3,fQ2sum); //cout<<"fQ2sum = "<CalculateRPandW(vQ); Double_t dWR = aLYZEP->GetWR(); Double_t dRP = aLYZEP->GetPsi(); //fHistProWr->Fill(vQ.Mod(),dWR); //this weight is always positive fHistPhiLYZ->Fill(dRP); //plot difference between event plane from EP-method and LYZ-method Double_t dRPEP = vQ.Phi()/2; //gives distribution from (0 to pi) //Double_t dRPEP = 0.5*TMath::ATan2(vQ.Y(),vQ.X()); //gives distribution from (-pi/2 to pi/2) //cout<<"dRPEP = "<< dRPEP <Fill(dRPEP); Double_t dDeltaPhi = dRPEP - dRP; if (dDeltaPhi < 0.) { dDeltaPhi += TMath::Pi(); } //to shift distribution from (-pi/2 to pi/2) to (0 to pi) //cout<<"dDeltaPhi = "<Fill(dDeltaPhi); //Flip sign of WR Double_t dLow = TMath::Pi()/4.; Double_t dHigh = 3.*(TMath::Pi()/4.); if ((dDeltaPhi > dLow) && (dDeltaPhi < dHigh)){ dRP -= (TMath::Pi()/2); dWR = -dWR; cerr<<"*** dRP modified ***"<Fill(vQ.Mod(),dWR); //corrected weight //calculate flow for RP and POI selections //loop over the tracks of the event Int_t iNumberOfTracks = anEvent->NumberOfTracks(); for (Int_t i=0;iGetTrack(i) ; if (pTrack){ Double_t dPhi = pTrack->Phi(); //if (dPhi<0.) fPhi+=2*TMath::Pi(); Double_t dPt = pTrack->Pt(); Double_t dEta = pTrack->Eta(); //calculate flow v2: Double_t dv2 = dWR * TMath::Cos(2*(dPhi-dRP)); if (pTrack->InRPSelection()) { //fill histograms for RP selection fHistProVetaRP -> Fill(dEta,dv2); fHistProVPtRP -> Fill(dPt,dv2); } if (pTrack->InPOISelection()) { //fill histograms for POI selection fHistProVetaPOI -> Fill(dEta,dv2); fHistProVPtPOI -> Fill(dPt,dv2); } }//track }//loop over tracks fEventNumber++; // cout<<"@@@@@ "< (outputListHistos->FindObject("AliFlowCommonHistLYZEP")); AliFlowCommonHistResults *pCommonHistResults = dynamic_cast (outputListHistos->FindObject("AliFlowCommonHistResultsLYZEP")); TProfile* pHistProR0theta = dynamic_cast (outputListHistos->FindObject("First_FlowPro_r0theta_LYZ")); TProfile* pHistProVetaRP = dynamic_cast (outputListHistos->FindObject("FlowPro_VetaRP_LYZEP")); TProfile* pHistProVetaPOI = dynamic_cast (outputListHistos->FindObject("FlowPro_VetaPOI_LYZEP")); TProfile* pHistProVPtRP = dynamic_cast (outputListHistos->FindObject("FlowPro_VPtRP_LYZEP")); TProfile* pHistProVPtPOI = dynamic_cast (outputListHistos->FindObject("FlowPro_VPtPOI_LYZEP")); TH1F* pHistQsumforChi = dynamic_cast (outputListHistos->FindObject("Flow_QsumforChi_LYZEP")); if (pCommonHist && pCommonHistResults && pHistProR0theta && pHistProVetaRP && pHistProVetaPOI && pHistProVPtRP && pHistProVPtPOI && pHistQsumforChi ) { this -> SetCommonHists(pCommonHist); this -> SetCommonHistsRes(pCommonHistResults); this -> SetFirstr0theta(pHistProR0theta); this -> SetHistProVetaRP(pHistProVetaRP); this -> SetHistProVetaPOI(pHistProVetaPOI); this -> SetHistProVPtRP(pHistProVPtRP); this -> SetHistProVPtPOI(pHistProVPtPOI); this -> SetHistQsumforChi(pHistQsumforChi); } } else { cout<<"WARNING: Histograms needed to run Finish() are not accessable!"<GetHistMultOrig()->GetEntries()); //cout<<"number of events processed is "<Set(fHistQsumforChi->GetBinContent(1),fHistQsumforChi->GetBinContent(2)); SetQ2sum(fHistQsumforChi->GetBinContent(3)); //calculate dV the mean of dVtheta Double_t dVtheta = 0; Double_t dV = 0; for (Int_t theta=0;thetaGetBinContent(theta+1); if (dR0!=0.) { dVtheta = dJ01/dR0 ;} dV += dVtheta; } dV /= iNtheta; //calculate dChi Double_t dSigma2 = 0; Double_t dChi= 0; if (fEventNumber!=0) { *fQsum /= fEventNumber; //cerr<<"fQsum->X() = "<X()<Y() = "<Y()<0) dChi = dV/TMath::Sqrt(dSigma2); else dChi = -1.; fCommonHistsRes->FillChiRP(dChi); // recalculate statistical errors on integrated flow //combining 5 theta angles to 1 relative error BP eq. 89 Double_t dRelErr2comb = 0.; Int_t iEvts = fEventNumber; if (iEvts!=0) { for (Int_t theta=0;thetaFillIntegratedFlow(dV, dVErr); cout<<"*************************************"<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* fHistPtRP = fCommonHists->GetHistPtRP(); //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); cout<<"dV(RP) = "<GetHistPtPOI(); //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); //cerr<<"dNprime = "<GetBinEntries(b)); //cerr<<"iNprime = "<FillIntegratedFlowPOI(dVPOI,dErrV); cout<<"dV(POI) = "<