]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowAnalysisWithMCEventPlane.cxx
spring cleanup
[u/mrichter/AliRoot.git] / PWG2 / FLOW / 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
23 class TH1F;
24
25 #include "AliFlowCommonConstants.h"    //needed as include
26 #include "AliFlowEventSimple.h"
27 #include "AliFlowTrackSimple.h"
28 #include "AliFlowCommonHist.h"
29 #include "AliFlowCommonHistResults.h"
30 #include "AliFlowAnalysisWithMCEventPlane.h"
31
32 class AliFlowVector;
33
34 // AliFlowAnalysisWithMCEventPlane:
35 // Description: Maker to analyze Flow from the generated MC reaction plane.
36 //              This class is used to get the real value of the flow 
37 //              to compare the other methods to when analysing simulated events
38 // author: N. van der Kolk (kolk@nikhef.nl)
39
40 ClassImp(AliFlowAnalysisWithMCEventPlane)
41
42   //-----------------------------------------------------------------------
43  
44  AliFlowAnalysisWithMCEventPlane::AliFlowAnalysisWithMCEventPlane():
45    fQsum(NULL),
46    fQ2sum(0),
47    fEventNumber(0),
48    fDebug(kFALSE),
49    fHistFileName(0),
50    fHistFile(0),
51    fCommonHists(NULL),
52    fCommonHistsRes(NULL),
53    fHistProFlow(NULL),
54    fHistRP(NULL)
55
56 {
57
58   // Constructor.
59   fQsum = new TVector2;        // flow vector sum
60 }
61
62  
63  //-----------------------------------------------------------------------
64
65
66  AliFlowAnalysisWithMCEventPlane::~AliFlowAnalysisWithMCEventPlane() 
67  {
68    //destructor
69    delete fQsum;
70  }
71  
72
73 //-----------------------------------------------------------------------
74 void AliFlowAnalysisWithMCEventPlane::Init() {
75
76   //Define all histograms
77   cout<<"---Analysis with the real MC Event Plane---"<<endl;
78
79   Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt();
80   Double_t  dPtMin = AliFlowCommonConstants::GetPtMin();             
81   Double_t  dPtMax = AliFlowCommonConstants::GetPtMax();
82
83   // analysis file (output)
84   fHistFile = new TFile(fHistFileName.Data(),"RECREATE") ;
85
86   fCommonHists = new AliFlowCommonHist("MC");
87   fCommonHistsRes = new AliFlowCommonHistResults("MC");
88
89   fHistProFlow = new TProfile("FlowPro_VPt_MC","FlowPro_VPt_MC",iNbinsPt,dPtMin,dPtMax);
90   fHistProFlow->SetXTitle("Pt");
91   fHistProFlow->SetYTitle("v2 (%)");
92
93   fHistRP = new TH1F("Flow_RP_MC","Flow_RP_MC",100,0.,3.14);
94   fHistRP->SetXTitle("Reaction Plane Angle");
95   fHistRP->SetYTitle("Counts");
96
97  
98   fEventNumber = 0;  //set number of events to zero
99       
100
101  
102 //-----------------------------------------------------------------------
103  
104 void AliFlowAnalysisWithMCEventPlane::Make(AliFlowEventSimple* anEvent, Double_t aRP) {
105
106   //Calculate v2 from the MC reaction plane
107   if (anEvent) {
108          
109     //fill control histograms     
110     fCommonHists->FillControlHistograms(anEvent);
111
112     //get the Q vector from the FlowEvent
113     AliFlowVector vQ = anEvent->GetQ(); 
114     //cout<<"vQ.Mod() = " << vQ.Mod() << endl;
115     //for chi calculation:
116     *fQsum += vQ;
117     //cout<<"fQsum.Mod() = "<<fQsum.Mod()<<endl;
118     fQ2sum += vQ.Mod2();
119     cout<<"fQ2sum = "<<fQ2sum<<endl;
120         
121     fHistRP->Fill(aRP);   
122               
123     //calculate flow
124     //loop over the tracks of the event
125     Int_t iNumberOfTracks = anEvent->NumberOfTracks(); 
126     for (Int_t i=0;i<iNumberOfTracks;i++) 
127       {
128         AliFlowTrackSimple* pTrack = anEvent->GetTrack(i) ; 
129         if (pTrack){
130           if (pTrack->UseForDifferentialFlow()) {
131             Double_t dPhi = pTrack->Phi();
132             //if (dPhi<0.) dPhi+=2*TMath::Pi();
133             //calculate flow v2:
134             Double_t dv2 = TMath::Cos(2*(dPhi-aRP));
135             Double_t dPt = pTrack->Pt();
136             //fill histogram
137             fHistProFlow->Fill(dPt,100*dv2);  
138           }  
139         }//track selected
140       }//loop over tracks
141           
142     fEventNumber++;
143     cout<<"@@@@@ "<<fEventNumber<<" events processed"<<endl;
144   }
145 }
146
147   //--------------------------------------------------------------------    
148 void AliFlowAnalysisWithMCEventPlane::Finish() {
149    
150   //*************make histograms etc. 
151   if (fDebug) cout<<"AliFlowAnalysisWithMCEventPlane::Terminate()"<<endl;
152      
153   Int_t iNbinsPt = AliFlowCommonConstants::GetNbinsPt();
154     
155   TH1F* fHistPtDiff = fCommonHists->GetfHistPtDiff();
156   Double_t dV = 0.;
157   Double_t dErrV = 0.;
158   Double_t dSum = 0.;
159   for(Int_t b=0;b<iNbinsPt;b++){
160     Double_t dv2pro = 0.;
161     Double_t dErrdifcomb = 0.; //in case error from profile is correct
162     if(fHistProFlow) {
163       dv2pro = fHistProFlow->GetBinContent(b);
164       dErrdifcomb = fHistProFlow->GetBinError(b);
165       //fill TH1D
166       fCommonHistsRes->FillDifferentialFlow(b, dv2pro, dErrdifcomb); 
167       if (fHistPtDiff){
168         //integrated flow
169         Double_t dYield = fHistPtDiff->GetBinContent(b);
170         dV += dv2pro/100*dYield ;
171         dSum += dYield;
172         //error on integrated flow
173         dErrV += dYield*dYield*(dErrdifcomb/100)*(dErrdifcomb/100);
174       }
175     }
176   }
177   dV /= dSum;  //because pt distribution should be normalised
178   dErrV /= dSum*dSum;
179   dErrV = TMath::Sqrt(dErrV);
180   cout<<"dV is "<<dV<<" +- "<<dErrV<<endl;
181   fCommonHistsRes->FillIntegratedFlow(dV,dErrV); 
182
183   // write to file
184   fHistFile->Write();
185           
186   cout<<".....finished"<<endl;
187  }
188
189