]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithMCEventPlane.cxx
fix coding violations
[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    fHistProFlow(NULL),
54    fHistRP(NULL),
55    fHistProIntFlow(NULL),
56    fHistProDiffFlowPtRP(NULL),
57    fHistProDiffFlowEtaRP(NULL),
58    fHistProDiffFlowPtPOI(NULL),
59    fHistProDiffFlowEtaPOI(NULL)
60 {
61
62   // Constructor.
63   fHistList = new TList();
64
65   fQsum = new TVector2;        // flow vector sum
66 }
67
68  
69  //-----------------------------------------------------------------------
70
71
72  AliFlowAnalysisWithMCEventPlane::~AliFlowAnalysisWithMCEventPlane() 
73  {
74    //destructor
75    delete fHistList;
76    delete fQsum;
77  }
78  
79 //-----------------------------------------------------------------------
80
81 void AliFlowAnalysisWithMCEventPlane::WriteHistograms(TString* outputFileName)
82 {
83  //store the final results in output .root file
84  TFile *output = new TFile(outputFileName->Data(),"RECREATE");
85  output->WriteObject(fHistList, "cobjMCEP","SingleKey");
86  delete output;
87 }
88
89 //-----------------------------------------------------------------------
90
91 void AliFlowAnalysisWithMCEventPlane::WriteHistograms(TString outputFileName)
92 {
93  //store the final results in output .root file
94  TFile *output = new TFile(outputFileName.Data(),"RECREATE");
95  output->WriteObject(fHistList, "cobjMCEP","SingleKey");
96  delete output;
97 }
98
99 //-----------------------------------------------------------------------
100 void AliFlowAnalysisWithMCEventPlane::Init() {
101
102   //Define all histograms
103   cout<<"---Analysis with the real MC Event Plane---"<<endl;
104
105   Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt();
106   Double_t dPtMin = AliFlowCommonConstants::GetPtMin();      
107   Double_t dPtMax = AliFlowCommonConstants::GetPtMax();
108   
109   Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta();
110   Double_t dEtaMin = AliFlowCommonConstants::GetEtaMin();            
111   Double_t dEtaMax = AliFlowCommonConstants::GetEtaMax();  
112
113   fCommonHists = new AliFlowCommonHist("AliFlowCommonHistMCEP");
114   fHistList->Add(fCommonHists);
115   fCommonHistsRes = new AliFlowCommonHistResults("AliFlowCommonHistResultsMCEP");
116   fHistList->Add(fCommonHistsRes);
117   
118   fHistProFlow = new TProfile("FlowPro_VPt_MCEP","FlowPro_VPt_MCEP",iNbinsPt,dPtMin,dPtMax);
119   fHistProFlow->SetXTitle("P_{t}");
120   fHistProFlow->SetYTitle("v_{2}");
121   fHistList->Add(fHistProFlow);
122
123   fHistRP = new TH1F("Flow_RP_MCEP","Flow_RP_MCEP",100,0.,3.14);
124   fHistRP->SetXTitle("Reaction Plane Angle");
125   fHistRP->SetYTitle("Counts");
126   fHistList->Add(fHistRP);
127   
128   fHistProIntFlow = new TProfile("fHistProIntFlow","fHistProIntFlow",1,0.,1.);
129   fHistProIntFlow->SetLabelSize(0.06);
130   (fHistProIntFlow->GetXaxis())->SetBinLabel(1,"v_{n}{2}");
131   fHistProIntFlow->SetYTitle("");
132   fHistList->Add(fHistProIntFlow);
133  
134   fHistProDiffFlowPtRP = new TProfile("fHistProDiffFlowPtRP","fHistProDiffFlowPtRP",iNbinsPt,dPtMin,dPtMax);
135   fHistProDiffFlowPtRP->SetXTitle("P_{t}");
136   fHistProDiffFlowPtRP->SetYTitle("");
137   fHistList->Add(fHistProDiffFlowPtRP);  
138   
139   fHistProDiffFlowEtaRP = new TProfile("fHistProDiffFlowEtaRP","fHistProDiffFlowEtaRP",iNbinsEta,dEtaMin,dEtaMax);
140   fHistProDiffFlowEtaRP->SetXTitle("#eta");
141   fHistProDiffFlowEtaRP->SetYTitle("");
142   fHistList->Add(fHistProDiffFlowEtaRP);
143   
144   fHistProDiffFlowPtPOI = new TProfile("fHistProDiffFlowPtPOI","fHistProDiffFlowPtPOI",iNbinsPt,dPtMin,dPtMax);
145   fHistProDiffFlowPtPOI->SetXTitle("P_{t}");
146   fHistProDiffFlowPtPOI->SetYTitle("");
147   fHistList->Add(fHistProDiffFlowPtPOI);  
148   
149   fHistProDiffFlowEtaPOI = new TProfile("fHistProDiffFlowEtaPOI","fHistProDiffFlowEtaPOI",iNbinsEta,dEtaMin,dEtaMax);
150   fHistProDiffFlowEtaPOI->SetXTitle("#eta");
151   fHistProDiffFlowEtaPOI->SetYTitle("");
152   fHistList->Add(fHistProDiffFlowEtaPOI);          
153  
154   fEventNumber = 0;  //set number of events to zero
155       
156
157  
158 //-----------------------------------------------------------------------
159  
160 void AliFlowAnalysisWithMCEventPlane::Make(AliFlowEventSimple* anEvent) {
161
162   //Calculate v2 from the MC reaction plane
163   if (anEvent) {
164   
165     // get the MC reaction plane angle
166     Double_t aRP = anEvent->GetMCReactionPlaneAngle();  
167     //fill control histograms     
168     fCommonHists->FillControlHistograms(anEvent);
169
170     //get the Q vector from the FlowEvent
171     AliFlowVector vQ = anEvent->GetQ(); 
172     //cout<<"vQ.Mod() = " << vQ.Mod() << endl;
173     //for chi calculation:
174     *fQsum += vQ;
175     //cout<<"fQsum.Mod() = "<<fQsum.Mod()<<endl;
176     fQ2sum += vQ.Mod2();
177     //cout<<"fQ2sum = "<<fQ2sum<<endl;
178         
179     fHistRP->Fill(aRP);   
180     
181     Double_t dPhi = 0.;
182     Double_t dv2  = 0.;
183     Double_t dPt  = 0.;
184     Double_t dEta = 0.;
185     //Double_t dPi = TMath::Pi();                   
186                                                              
187     //calculate flow
188     //loop over the tracks of the event
189     Int_t iNumberOfTracks = anEvent->NumberOfTracks(); 
190     for (Int_t i=0;i<iNumberOfTracks;i++) 
191       {
192         AliFlowTrackSimple* pTrack = anEvent->GetTrack(i) ; 
193         if (pTrack){
194           if (pTrack->UseForIntegratedFlow()){
195             dPhi = pTrack->Phi();
196             dv2 = TMath::Cos(2*(dPhi-aRP));
197             dPt = pTrack->Pt();
198             dEta = pTrack->Eta();
199             //no-name int. flow (to be improved):
200             fHistProIntFlow->Fill(0.,dv2);
201             //differential flow (Pt, RP):
202             fHistProDiffFlowPtRP->Fill(dPt,dv2,1.);
203             //differential flow (Eta, RP):
204             fHistProDiffFlowEtaRP->Fill(dEta,dv2,1.);
205           }
206           if (pTrack->UseForDifferentialFlow()) {
207             dPhi = pTrack->Phi();
208             //if (dPhi<0.) dPhi+=2*TMath::Pi();
209             //calculate flow v2:
210             dv2 = TMath::Cos(2*(dPhi-aRP));
211             dPt = pTrack->Pt();
212             dEta = pTrack->Eta();
213             //differential flow (Pt, POI):
214             fHistProDiffFlowPtPOI->Fill(dPt,dv2,1.);
215             //differential flow (Eta, POI):
216             fHistProDiffFlowEtaPOI->Fill(dEta,dv2,1.); 
217           }           
218         }//track selected
219       }//loop over tracks
220           
221     fEventNumber++;
222     //    cout<<"@@@@@ "<<fEventNumber<<" events processed"<<endl;
223   }
224 }
225
226   //--------------------------------------------------------------------    
227 void AliFlowAnalysisWithMCEventPlane::Finish() {
228    
229   //*************make histograms etc. 
230   if (fDebug) cout<<"AliFlowAnalysisWithMCEventPlane::Terminate()"<<endl;
231    
232   Int_t iNbinsPt  = AliFlowCommonConstants::GetNbinsPt();  
233   Int_t iNbinsEta = AliFlowCommonConstants::GetNbinsEta(); 
234          
235   // no-name int. flow (to be improved):
236   Double_t dV = fHistProIntFlow->GetBinContent(1);  
237   Double_t dErrV = fHistProIntFlow->GetBinError(1); // to be improved (treatment of errors for non-Gaussian distribution needed!)  
238   // fill no-name int. flow (to be improved):
239   fCommonHistsRes->FillIntegratedFlow(dV,dErrV);
240   cout<<"dV{MC} is       "<<dV<<" +- "<<dErrV<<endl;
241   
242   //RP:
243   TH1F* fHistPtRP = fCommonHists->GetHistPtInt(); // to be improved (change "int" and "diff" to RP and POI in common control histos)
244   Double_t dYieldPtRP = 0.;
245   Double_t dVRP = 0.;
246   Double_t dErrVRP = 0.;
247   Double_t dSumRP = 0.;
248   //differential flow (RP, Pt): 
249   Double_t dvPtRP = 0.;           
250   Double_t dErrvPtRP = 0.;
251   for(Int_t b=1;b<iNbinsPt;b++)
252   {
253    dvPtRP    = fHistProDiffFlowPtRP->GetBinContent(b);
254    dErrvPtRP = fHistProDiffFlowPtRP->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
255    fCommonHistsRes->FillDifferentialFlowPtRP(b, dvPtRP, dErrvPtRP);
256    if(fHistPtRP){
257         //integrated flow (RP)
258         dYieldPtRP = fHistPtRP->GetBinContent(b);
259         dVRP += dvPtRP*dYieldPtRP;
260         dSumRP += dYieldPtRP;
261         //error on integrated flow
262         dErrVRP += dYieldPtRP*dYieldPtRP*dErrvPtRP*dErrvPtRP;
263       }
264   }
265   if (dSumRP != 0. ) {
266     dVRP /= dSumRP;  //because pt distribution should be normalised
267     dErrVRP /= (dSumRP*dSumRP);
268     dErrVRP = TMath::Sqrt(dErrVRP); 
269   }
270   // fill integrated flow (RP):
271   fCommonHistsRes->FillIntegratedFlowRP(dVRP,dErrVRP);
272   cout<<"dV{MC} (RP) is  "<<dVRP<<" +- "<<dErrVRP<<endl;
273   
274   //differential flow (RP, Eta): 
275   Double_t dvEtaRP = 0.;           
276   Double_t dErrvEtaRP = 0.;
277   for(Int_t b=1;b<iNbinsEta;b++)
278   {
279    dvEtaRP    = fHistProDiffFlowEtaRP->GetBinContent(b);
280    dErrvEtaRP = fHistProDiffFlowEtaRP->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
281    fCommonHistsRes->FillDifferentialFlowEtaRP(b, dvEtaRP, dErrvEtaRP);
282   }
283                                                                                                                                    
284   //POI:
285   TH1F* fHistPtPOI = fCommonHists->GetHistPtDiff(); // to be improved (change "int" and "diff" to RP and POI in common control histos)
286   Double_t dYieldPtPOI = 0.;
287   Double_t dVPOI = 0.;
288   Double_t dErrVPOI = 0.;
289   Double_t dSumPOI = 0.;
290   Double_t dv2proPtPOI = 0.;
291   Double_t dErrdifcombPtPOI = 0.; 
292   Double_t dv2proEtaPOI = 0.;
293   Double_t dErrdifcombEtaPOI = 0.;   
294   //Pt:
295   if(fHistProFlow && fHistProDiffFlowPtPOI) {//to be removed (fHistProFlow)
296     for(Int_t b=1;b<iNbinsPt;b++){
297       //dv2pro = fHistProFlow->GetBinContent(b);//to be removed
298       //dErrdifcomb = fHistProFlow->GetBinError(b);//to be removed
299       dv2proPtPOI = fHistProDiffFlowPtPOI->GetBinContent(b);
300       dErrdifcombPtPOI = fHistProDiffFlowPtPOI->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
301       //fill TH1D
302       fCommonHistsRes->FillDifferentialFlowPtPOI(b, dv2proPtPOI, dErrdifcombPtPOI); 
303       if (fHistPtPOI){
304         //integrated flow (POI)
305         dYieldPtPOI = fHistPtPOI->GetBinContent(b);
306         dVPOI += dv2proPtPOI*dYieldPtPOI;
307         dSumPOI += dYieldPtPOI;
308         //error on integrated flow
309         dErrVPOI += dYieldPtPOI*dYieldPtPOI*dErrdifcombPtPOI*dErrdifcombPtPOI;
310       }
311     }//end of for(Int_t b=0;b<iNbinsPt;b++)  
312   } else { cout<<"fHistProFlow is NULL"<<endl; }
313   if (dSumPOI != 0. ) {
314     dVPOI /= dSumPOI;  //because pt distribution should be normalised
315     dErrVPOI /= (dSumPOI*dSumPOI);
316     dErrVPOI = TMath::Sqrt(dErrVPOI); 
317   }
318   cout<<"dV{MC} (POI) is "<<dVPOI<<" +- "<<dErrVPOI<<endl;
319
320   fCommonHistsRes->FillIntegratedFlowPOI(dVPOI,dErrVPOI);
321   
322   //Eta:
323   if(fHistProDiffFlowEtaPOI)
324   {
325    for(Int_t b=1;b<iNbinsEta;b++)
326    {
327     dv2proEtaPOI = fHistProDiffFlowEtaPOI->GetBinContent(b);
328     dErrdifcombEtaPOI = fHistProDiffFlowEtaPOI->GetBinError(b);//to be improved (treatment of errors for non-Gaussian distribution needed!)
329     //fill common hist results:
330     fCommonHistsRes->FillDifferentialFlowEtaPOI(b, dv2proEtaPOI, dErrdifcombEtaPOI); 
331    }
332   }   
333   
334   cout<<endl;                     
335   cout<<".....finished"<<endl;
336 }
337
338  
339