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 AliFlowAnalysisWithMCEventPlane_cxx
18 #include "Riostream.h" //needed as include
19 #include "TFile.h" //needed as include
20 #include "TProfile.h" //needed as include
21 #include "TComplex.h" //needed as include
26 #include "AliFlowCommonConstants.h" //needed as include
27 #include "AliFlowEventSimple.h"
28 #include "AliFlowTrackSimple.h"
29 #include "AliFlowCommonHist.h"
30 #include "AliFlowCommonHistResults.h"
31 #include "AliFlowAnalysisWithMCEventPlane.h"
35 // AliFlowAnalysisWithMCEventPlane:
36 // Description: Maker to analyze Flow from the generated MC reaction plane.
37 // This class is used to get the real value of the flow
38 // to compare the other methods to when analysing simulated events
39 // author: N. van der Kolk (kolk@nikhef.nl)
41 ClassImp(AliFlowAnalysisWithMCEventPlane)
43 //-----------------------------------------------------------------------
45 AliFlowAnalysisWithMCEventPlane::AliFlowAnalysisWithMCEventPlane():
52 fCommonHistsRes(NULL),
54 fHistProIntFlow(NULL),
55 fHistProDiffFlowPtRP(NULL),
56 fHistProDiffFlowEtaRP(NULL),
57 fHistProDiffFlowPtPOI(NULL),
58 fHistProDiffFlowEtaPOI(NULL)
62 fHistList = new TList();
64 fQsum = new TVector2; // flow vector sum
68 //-----------------------------------------------------------------------
71 AliFlowAnalysisWithMCEventPlane::~AliFlowAnalysisWithMCEventPlane()
78 //-----------------------------------------------------------------------
80 void AliFlowAnalysisWithMCEventPlane::WriteHistograms(TString* outputFileName)
82 //store the final results in output .root file
83 TFile *output = new TFile(outputFileName->Data(),"RECREATE");
84 output->WriteObject(fHistList, "cobjMCEP","SingleKey");
88 //-----------------------------------------------------------------------
90 void AliFlowAnalysisWithMCEventPlane::WriteHistograms(TString outputFileName)
92 //store the final results in output .root file
93 TFile *output = new TFile(outputFileName.Data(),"RECREATE");
94 output->WriteObject(fHistList, "cobjMCEP","SingleKey");
98 //-----------------------------------------------------------------------
99 void AliFlowAnalysisWithMCEventPlane::Init() {
101 //Define all histograms
102 cout<<"---Analysis with the real MC Event Plane---"<<endl;
104 Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt();
105 Double_t dPtMin = AliFlowCommonConstants::GetPtMin();
106 Double_t dPtMax = AliFlowCommonConstants::GetPtMax();
108 Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta();
109 Double_t dEtaMin = AliFlowCommonConstants::GetEtaMin();
110 Double_t dEtaMax = AliFlowCommonConstants::GetEtaMax();
112 fCommonHists = new AliFlowCommonHist("AliFlowCommonHistMCEP");
113 fHistList->Add(fCommonHists);
114 fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsMCEP");
115 fHistList->Add(fCommonHistsRes);
117 fHistRP = new TH1F("Flow_RP_MCEP","Flow_RP_MCEP",100,0.,3.14);
118 fHistRP->SetXTitle("Reaction Plane Angle");
119 fHistRP->SetYTitle("Counts");
120 fHistList->Add(fHistRP);
122 fHistProIntFlow = new TProfile("FlowPro_V_MCEP","FlowPro_V_MCEP",1,0.,1.);
123 fHistProIntFlow->SetLabelSize(0.06);
124 (fHistProIntFlow->GetXaxis())->SetBinLabel(1,"v_{n}{2}");
125 fHistProIntFlow->SetYTitle("");
126 fHistList->Add(fHistProIntFlow);
128 fHistProDiffFlowPtRP = new TProfile("FlowPro_VPtRP_MCEP","FlowPro_VPtRP_MCEP",iNbinsPt,dPtMin,dPtMax);
129 fHistProDiffFlowPtRP->SetXTitle("P_{t}");
130 fHistProDiffFlowPtRP->SetYTitle("");
131 fHistList->Add(fHistProDiffFlowPtRP);
133 fHistProDiffFlowEtaRP = new TProfile("FlowPro_VetaRP_MCEP","FlowPro_VetaRP_MCEP",iNbinsEta,dEtaMin,dEtaMax);
134 fHistProDiffFlowEtaRP->SetXTitle("#eta");
135 fHistProDiffFlowEtaRP->SetYTitle("");
136 fHistList->Add(fHistProDiffFlowEtaRP);
138 fHistProDiffFlowPtPOI = new TProfile("FlowPro_VPtPOI_MCEP","FlowPro_VPtPOI_MCEP",iNbinsPt,dPtMin,dPtMax);
139 fHistProDiffFlowPtPOI->SetXTitle("P_{t}");
140 fHistProDiffFlowPtPOI->SetYTitle("");
141 fHistList->Add(fHistProDiffFlowPtPOI);
143 fHistProDiffFlowEtaPOI = new TProfile("FlowPro_VetaPOI_MCEP","FlowPro_VetaPOI_MCEP",iNbinsEta,dEtaMin,dEtaMax);
144 fHistProDiffFlowEtaPOI->SetXTitle("#eta");
145 fHistProDiffFlowEtaPOI->SetYTitle("");
146 fHistList->Add(fHistProDiffFlowEtaPOI);
148 fEventNumber = 0; //set number of events to zero
152 //-----------------------------------------------------------------------
154 void AliFlowAnalysisWithMCEventPlane::Make(AliFlowEventSimple* anEvent) {
156 //Calculate v2 from the MC reaction plane
159 // get the MC reaction plane angle
160 Double_t aRP = anEvent->GetMCReactionPlaneAngle();
161 //fill control histograms
162 fCommonHists->FillControlHistograms(anEvent);
164 //get the Q vector from the FlowEvent
165 AliFlowVector vQ = anEvent->GetQ();
166 //cout<<"vQ.Mod() = " << vQ.Mod() << endl;
167 //for chi calculation:
169 //cout<<"fQsum.Mod() = "<<fQsum.Mod()<<endl;
171 //cout<<"fQ2sum = "<<fQ2sum<<endl;
179 //Double_t dPi = TMath::Pi();
182 //loop over the tracks of the event
183 Int_t iNumberOfTracks = anEvent->NumberOfTracks();
184 for (Int_t i=0;i<iNumberOfTracks;i++)
186 AliFlowTrackSimple* pTrack = anEvent->GetTrack(i) ;
188 if (pTrack->InRPSelection()){
189 dPhi = pTrack->Phi();
190 dv2 = TMath::Cos(2*(dPhi-aRP));
192 dEta = pTrack->Eta();
193 //no-name int. flow (to be improved):
194 fHistProIntFlow->Fill(0.,dv2);
195 //differential flow (Pt, RP):
196 fHistProDiffFlowPtRP->Fill(dPt,dv2,1.);
197 //differential flow (Eta, RP):
198 fHistProDiffFlowEtaRP->Fill(dEta,dv2,1.);
200 if (pTrack->InPOISelection()) {
201 dPhi = pTrack->Phi();
202 //if (dPhi<0.) dPhi+=2*TMath::Pi();
204 dv2 = TMath::Cos(2*(dPhi-aRP));
206 dEta = pTrack->Eta();
207 //differential flow (Pt, POI):
208 fHistProDiffFlowPtPOI->Fill(dPt,dv2,1.);
209 //differential flow (Eta, POI):
210 fHistProDiffFlowEtaPOI->Fill(dEta,dv2,1.);
216 // cout<<"@@@@@ "<<fEventNumber<<" events processed"<<endl;
220 //--------------------------------------------------------------------
221 void AliFlowAnalysisWithMCEventPlane::Finish() {
223 //*************make histograms etc.
224 if (fDebug) cout<<"AliFlowAnalysisWithMCEventPlane::Terminate()"<<endl;
226 Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt();
227 Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta();
229 // no-name int. flow (to be improved):
230 Double_t dV = fHistProIntFlow->GetBinContent(1);
231 Double_t dErrV = fHistProIntFlow->GetBinError(1); // to be improved (treatment of errors for non-Gaussian distribution needed!)
232 // fill no-name int. flow (to be improved):
233 fCommonHistsRes->FillIntegratedFlow(dV,dErrV);
234 cout<<"dV{MC} is "<<dV<<" +- "<<dErrV<<endl;
237 TH1F* fHistPtRP = fCommonHists->GetHistPtRP();
238 Double_t dYieldPtRP = 0.;
240 Double_t dErrVRP = 0.;
241 Double_t dSumRP = 0.;
242 //differential flow (RP, Pt):
243 Double_t dvPtRP = 0.;
244 Double_t dErrvPtRP = 0.;
245 for(Int_t b=1;b<iNbinsPt;b++)
247 dvPtRP = fHistProDiffFlowPtRP->GetBinContent(b);
248 dErrvPtRP = fHistProDiffFlowPtRP->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
249 fCommonHistsRes->FillDifferentialFlowPtRP(b, dvPtRP, dErrvPtRP);
251 //integrated flow (RP)
252 dYieldPtRP = fHistPtRP->GetBinContent(b);
253 dVRP += dvPtRP*dYieldPtRP;
254 dSumRP += dYieldPtRP;
255 //error on integrated flow
256 dErrVRP += dYieldPtRP*dYieldPtRP*dErrvPtRP*dErrvPtRP;
260 dVRP /= dSumRP; //because pt distribution should be normalised
261 dErrVRP /= (dSumRP*dSumRP);
262 dErrVRP = TMath::Sqrt(dErrVRP);
264 // fill integrated flow (RP):
265 fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrVRP);
266 cout<<"dV{MC} (RP) is "<<dVRP<<" +- "<<dErrVRP<<endl;
268 //differential flow (RP, Eta):
269 Double_t dvEtaRP = 0.;
270 Double_t dErrvEtaRP = 0.;
271 for(Int_t b=1;b<iNbinsEta;b++)
273 dvEtaRP = fHistProDiffFlowEtaRP->GetBinContent(b);
274 dErrvEtaRP = fHistProDiffFlowEtaRP->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
275 fCommonHistsRes->FillDifferentialFlowEtaRP(b, dvEtaRP, dErrvEtaRP);
279 TH1F* fHistPtPOI = fCommonHists->GetHistPtPOI();
280 Double_t dYieldPtPOI = 0.;
282 Double_t dErrVPOI = 0.;
283 Double_t dSumPOI = 0.;
284 Double_t dv2proPtPOI = 0.;
285 Double_t dErrdifcombPtPOI = 0.;
286 Double_t dv2proEtaPOI = 0.;
287 Double_t dErrdifcombEtaPOI = 0.;
289 if(fHistProDiffFlowPtPOI) {
290 for(Int_t b=1;b<iNbinsPt;b++){
291 dv2proPtPOI = fHistProDiffFlowPtPOI->GetBinContent(b);
292 dErrdifcombPtPOI = fHistProDiffFlowPtPOI->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
294 fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2proPtPOI, dErrdifcombPtPOI);
296 //integrated flow (POI)
297 dYieldPtPOI = fHistPtPOI->GetBinContent(b);
298 dVPOI += dv2proPtPOI*dYieldPtPOI;
299 dSumPOI += dYieldPtPOI;
300 //error on integrated flow
301 dErrVPOI += dYieldPtPOI*dYieldPtPOI*dErrdifcombPtPOI*dErrdifcombPtPOI;
303 }//end of for(Int_t b=0;b<iNbinsPt;b++)
304 } else { cout<<"fHistProFlow is NULL"<<endl; }
305 if (dSumPOI != 0. ) {
306 dVPOI /= dSumPOI; //because pt distribution should be normalised
307 dErrVPOI /= (dSumPOI*dSumPOI);
308 dErrVPOI = TMath::Sqrt(dErrVPOI);
310 cout<<"dV{MC} (POI) is "<<dVPOI<<" +- "<<dErrVPOI<<endl;
312 fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrVPOI);
315 if(fHistProDiffFlowEtaPOI)
317 for(Int_t b=1;b<iNbinsEta;b++)
319 dv2proEtaPOI = fHistProDiffFlowEtaPOI->GetBinContent(b);
320 dErrdifcombEtaPOI = fHistProDiffFlowEtaPOI->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
321 //fill common hist results:
322 fCommonHistsRes->FillDifferentialFlowEtaPOI(b, dv2proEtaPOI, dErrdifcombEtaPOI);
327 //cout<<".....finished"<<endl;