]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithMCEventPlane.cxx
3f1df2c14c97409d9362afd4af0aa40865d1427c
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowAnalysisWithMCEventPlane.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
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 **************************************************************************/
15
16 #define AliFlowAnalysisWithMCEventPlane_cxx
17  
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
22 #include "TList.h"
23
24 class TH1F;
25
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"
32
33 class AliFlowVector;
34
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)
40
41 ClassImp(AliFlowAnalysisWithMCEventPlane)
42
43   //-----------------------------------------------------------------------
44  
45  AliFlowAnalysisWithMCEventPlane::AliFlowAnalysisWithMCEventPlane():
46    fQsum(NULL),
47    fQ2sum(0),
48    fEventNumber(0),
49    fDebug(kFALSE),
50    fHistList(NULL),
51    fCommonHists(NULL),
52    fCommonHistsRes(NULL),
53    fHistRP(NULL),
54    fHistProIntFlow(NULL),
55    fHistProDiffFlowPtRP(NULL),
56    fHistProDiffFlowEtaRP(NULL),
57    fHistProDiffFlowPtPOI(NULL),
58    fHistProDiffFlowEtaPOI(NULL)
59 {
60
61   // Constructor.
62   fHistList = new TList();
63
64   fQsum = new TVector2;        // flow vector sum
65 }
66
67  
68  //-----------------------------------------------------------------------
69
70
71  AliFlowAnalysisWithMCEventPlane::~AliFlowAnalysisWithMCEventPlane() 
72  {
73    //destructor
74    delete fHistList;
75    delete fQsum;
76  }
77  
78 //-----------------------------------------------------------------------
79
80 void AliFlowAnalysisWithMCEventPlane::WriteHistograms(TString* outputFileName)
81 {
82  //store the final results in output .root file
83  TFile *output = new TFile(outputFileName->Data(),"RECREATE");
84  output->WriteObject(fHistList, "cobjMCEP","SingleKey");
85  delete output;
86 }
87
88 //-----------------------------------------------------------------------
89
90 void AliFlowAnalysisWithMCEventPlane::WriteHistograms(TString outputFileName)
91 {
92  //store the final results in output .root file
93  TFile *output = new TFile(outputFileName.Data(),"RECREATE");
94  output->WriteObject(fHistList, "cobjMCEP","SingleKey");
95  delete output;
96 }
97
98 //-----------------------------------------------------------------------
99 void AliFlowAnalysisWithMCEventPlane::Init() {
100
101   //Define all histograms
102   cout<<"---Analysis with the real MC Event Plane---"<<endl;
103
104   Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt();
105   Double_t dPtMin = AliFlowCommonConstants::GetPtMin();      
106   Double_t dPtMax = AliFlowCommonConstants::GetPtMax();
107   
108   Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta();
109   Double_t dEtaMin = AliFlowCommonConstants::GetEtaMin();            
110   Double_t dEtaMax = AliFlowCommonConstants::GetEtaMax();  
111
112   fCommonHists = new AliFlowCommonHist("AliFlowCommonHistMCEP");
113   fHistList->Add(fCommonHists);
114   fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsMCEP");
115   fHistList->Add(fCommonHistsRes);
116   
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);
121   
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);
127  
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);  
132   
133   fHistProDiffFlowEtaRP = new TProfile("FlowPro_VetaRP_MCEP","FlowPro_VetaRP_MCEP",iNbinsEta,dEtaMin,dEtaMax);
134   fHistProDiffFlowEtaRP->SetXTitle("#eta");
135   fHistProDiffFlowEtaRP->SetYTitle("");
136   fHistList->Add(fHistProDiffFlowEtaRP);
137   
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);  
142   
143   fHistProDiffFlowEtaPOI = new TProfile("FlowPro_VetaPOI_MCEP","FlowPro_VetaPOI_MCEP",iNbinsEta,dEtaMin,dEtaMax);
144   fHistProDiffFlowEtaPOI->SetXTitle("#eta");
145   fHistProDiffFlowEtaPOI->SetYTitle("");
146   fHistList->Add(fHistProDiffFlowEtaPOI);          
147  
148   fEventNumber = 0;  //set number of events to zero
149       
150
151  
152 //-----------------------------------------------------------------------
153  
154 void AliFlowAnalysisWithMCEventPlane::Make(AliFlowEventSimple* anEvent) {
155
156   //Calculate v2 from the MC reaction plane
157   if (anEvent) {
158   
159     // get the MC reaction plane angle
160     Double_t aRP = anEvent->GetMCReactionPlaneAngle();  
161     //fill control histograms     
162     fCommonHists->FillControlHistograms(anEvent);
163
164     //get the Q vector from the FlowEvent
165     AliFlowVector vQ = anEvent->GetQ(); 
166     //cout<<"vQ.Mod() = " << vQ.Mod() << endl;
167     //for chi calculation:
168     *fQsum += vQ;
169     //cout<<"fQsum.Mod() = "<<fQsum.Mod()<<endl;
170     fQ2sum += vQ.Mod2();
171     //cout<<"fQ2sum = "<<fQ2sum<<endl;
172         
173     fHistRP->Fill(aRP);   
174     
175     Double_t dPhi = 0.;
176     Double_t dv2  = 0.;
177     Double_t dPt  = 0.;
178     Double_t dEta = 0.;
179     //Double_t dPi = TMath::Pi();                   
180                                                              
181     //calculate flow
182     //loop over the tracks of the event
183     Int_t iNumberOfTracks = anEvent->NumberOfTracks(); 
184     for (Int_t i=0;i<iNumberOfTracks;i++) 
185       {
186         AliFlowTrackSimple* pTrack = anEvent->GetTrack(i) ; 
187         if (pTrack){
188           if (pTrack->InRPSelection()){
189             dPhi = pTrack->Phi();
190             dv2  = TMath::Cos(2*(dPhi-aRP));
191             dPt  = pTrack->Pt();
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.);
199           }
200           if (pTrack->InPOISelection()) {
201             dPhi = pTrack->Phi();
202             //if (dPhi<0.) dPhi+=2*TMath::Pi();
203             //calculate flow v2:
204             dv2  = TMath::Cos(2*(dPhi-aRP));
205             dPt  = pTrack->Pt();
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.); 
211           }           
212         }//track selected
213       }//loop over tracks
214           
215     fEventNumber++;
216     //    cout<<"@@@@@ "<<fEventNumber<<" events processed"<<endl;
217   }
218 }
219
220   //--------------------------------------------------------------------    
221 void AliFlowAnalysisWithMCEventPlane::Finish() {
222    
223   //*************make histograms etc. 
224   if (fDebug) cout<<"AliFlowAnalysisWithMCEventPlane::Terminate()"<<endl;
225    
226   Int_t iNbinsPt  = AliFlowCommonConstants::GetNbinsPt();  
227   Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta(); 
228          
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;
235   
236   //RP:
237   TH1F* fHistPtRP = fCommonHists->GetHistPtRP(); 
238   Double_t dYieldPtRP = 0.;
239   Double_t dVRP = 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++)
246   {
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);
250    if(fHistPtRP){
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;
257       }
258   }
259   if (dSumRP != 0. ) {
260     dVRP /= dSumRP;  //because pt distribution should be normalised
261     dErrVRP /= (dSumRP*dSumRP);
262     dErrVRP = TMath::Sqrt(dErrVRP); 
263   }
264   // fill integrated flow (RP):
265   fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrVRP);
266   cout<<"dV{MC} (RP) is  "<<dVRP<<" +- "<<dErrVRP<<endl;
267   
268   //differential flow (RP, Eta): 
269   Double_t dvEtaRP = 0.;           
270   Double_t dErrvEtaRP = 0.;
271   for(Int_t b=1;b<iNbinsEta;b++)
272   {
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);
276   }
277                                                                                                                                    
278   //POI:
279   TH1F* fHistPtPOI = fCommonHists->GetHistPtPOI(); 
280   Double_t dYieldPtPOI = 0.;
281   Double_t dVPOI = 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.;   
288   //Pt:
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!)
293       //fill TH1D
294       fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2proPtPOI, dErrdifcombPtPOI); 
295       if (fHistPtPOI){
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;
302       }
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); 
309   }
310   cout<<"dV{MC} (POI) is "<<dVPOI<<" +- "<<dErrVPOI<<endl;
311
312   fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrVPOI);
313   
314   //Eta:
315   if(fHistProDiffFlowEtaPOI)
316   {
317    for(Int_t b=1;b<iNbinsEta;b++)
318    {
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); 
323    }
324   }   
325   
326   cout<<endl;                     
327   //cout<<".....finished"<<endl;
328 }
329
330  
331