]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliAnalysisTaskFlowEventforRP.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliAnalysisTaskFlowEventforRP.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 ////////////////////////////////////////////////////
17 // AliAnalysisTaskFlowEventforRP:
18 //
19 // analysis task for filling the flow event
20 // from MCEvent, ESD
21 // and put it in an output stream so the calculated
22 // Reaction Plane can be stored in the AODHeader
23 // when the AOD is made from the ESD 
24 // for cuts the correction framework is used
25 // which also outputs QA histograms to view
26 // the effects of the cuts
27 ////////////////////////////////////////////////////
28
29 #include "Riostream.h" //needed as include
30 #include "TChain.h"
31 #include "TTree.h"
32 #include "TFile.h" //needed as include
33 #include "TList.h"
34 #include "TRandom3.h"
35 #include "TTimeStamp.h"
36
37 // ALICE Analysis Framework
38 #include "AliAnalysisTaskSE.h"
39 #include "AliAnalysisManager.h"
40
41 // ESD interface
42 #include "AliESDEvent.h"
43 #include "AliESDInputHandler.h"
44
45 // AOD interface
46 #include "AliAODEvent.h"
47 #include "AliAODInputHandler.h"
48
49 // Monte Carlo Eventp
50 #include "AliAODHandler.h"
51 #include "AliMCEventHandler.h"
52 #include "AliMCEvent.h"
53
54 // ALICE Correction Framework
55 #include "AliCFManager.h"
56
57 // Interface to Event generators to get Reaction Plane Angle
58 #include "AliGenCocktailEventHeader.h"
59 #include "AliGenHijingEventHeader.h"
60 #include "AliGenGeVSimEventHeader.h"
61 #include "AliGenEposEventHeader.h"
62
63 // Interface to make the Flow Event Simple used in the flow analysis methods
64 #include "AliFlowEvent.h"
65 #include "AliFlowVector.h"
66 #include "AliAnalysisTaskFlowEventforRP.h"
67
68 using std::cout;
69 using std::endl;
70 ClassImp(AliAnalysisTaskFlowEventforRP)
71   
72 //________________________________________________________________________
73 AliAnalysisTaskFlowEventforRP::AliAnalysisTaskFlowEventforRP(const char *name) : 
74   AliAnalysisTaskSE(name), 
75   fAnalysisType("ESD"),
76   fCFManager1(NULL),
77   fCFManager2(NULL),
78   fMinMult(0),
79   fMaxMult(10000000),
80   fMCReactionPlaneAngle(0.)
81     
82 {
83   // Constructor
84   cout<<"AliAnalysisTaskFlowEventforRP::AliAnalysisTaskFlowEventforRP(const char *name)"<<endl;
85
86   // Define input and output slots here
87   // Input slot #0 works with a TChain
88   DefineInput(0, TChain::Class());
89   // Define here the flow event output
90   DefineOutput(0, AliFlowEventSimple::Class());  
91   
92 }
93
94 //________________________________________________________________________
95 AliAnalysisTaskFlowEventforRP::AliAnalysisTaskFlowEventforRP() : 
96   fAnalysisType("ESD"),
97   fCFManager1(NULL),
98   fCFManager2(NULL),
99   fMinMult(0),
100   fMaxMult(10000000),
101   fMCReactionPlaneAngle(0.)
102 {
103   // Constructor
104   cout<<"AliAnalysisTaskFlowEventforRP::AliAnalysisTaskFlowEventforRP()"<<endl;
105 }
106
107
108 //________________________________________________________________________
109 AliAnalysisTaskFlowEventforRP::~AliAnalysisTaskFlowEventforRP()
110 {
111   //
112   // Destructor
113   //
114   
115 }
116
117
118 //________________________________________________________________________
119 void AliAnalysisTaskFlowEventforRP::UserCreateOutputObjects() 
120 {
121   // Called at every worker node to initialize
122   cout<<"AliAnalysisTaskFlowEventforRP::UserCreateOutputObjects()"<<endl;
123
124   if (!(fAnalysisType == "ESD")) {
125     cout<<"WRONG ANALYSIS TYPE! only ESD for this method."<<endl;
126     exit(1);
127   }
128     
129 }
130
131 //________________________________________________________________________
132 void AliAnalysisTaskFlowEventforRP::UserExec(Option_t *) 
133 {
134   // Main loop
135
136   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
137   AliFlowEvent* fEvent = NULL;
138   AliMCEvent* mcEvent  = MCEvent();
139   Double_t fRP = 0.; // the monte carlo reaction plane angle
140   
141   // Fill the FlowEventSimple for ESD input  
142   if (fAnalysisType == "ESD") {
143     if (!(fCFManager1&&fCFManager2))
144       {
145         cout << "ERROR: No pointer to correction framework cuts! " << endl; 
146         return; 
147       }
148     if (!esd)
149       {
150         AliError("ERROR: ESD not available");
151         return;
152       }
153     
154     //check the offline trigger (check if the event has the correct trigger)
155     //AliInfo(Form("ESD has %d tracks", fInputEvent->GetNumberOfTracks()));
156     
157     //check multiplicity
158     if (!fCFManager1->CheckEventCuts(AliCFManager::kEvtRecCuts,esd))
159       {
160         cout << "Event does not pass multiplicity cuts" << endl;
161         return;
162       }
163     
164     // make the flowevent
165     fEvent = new AliFlowEvent(esd,fCFManager1,fCFManager2);
166     
167     if (mcEvent && mcEvent->GenEventHeader()) 
168       {
169         fEvent->SetMCReactionPlaneAngle(mcEvent);
170         fRP = fEvent->GetMCReactionPlaneAngle();
171       }
172     
173     //check final event cuts
174     Int_t mult = fEvent->NumberOfTracks();
175     cout << "FlowEvent has "<<mult<<" tracks"<<endl;
176     if (mult<fMinMult || mult>fMaxMult)
177       {
178         cout << "FlowEvent cut on multiplicity" << endl;
179         return;
180       }
181
182   
183     // get the flow vector     
184     AliFlowVector vQ = fEvent->GetQ();                      
185     Double_t dRP[1] = {0.0};   
186     // Phi is a Double_t, but SetQTheta() needs as input Double_t*, 
187     // an array of doubles. 
188     dRP[0] = vQ.Phi()/2; 
189       
190     cout<<"The reaction plane from MC is "<<fRP<<endl;
191     cout<<"The calculated reaction plane is "<<dRP[0]<<endl;
192
193     
194     // Update the header
195     AliAODHeader* header = dynamic_cast<AliAODHeader*>(AODEvent()->GetHeader());
196     if(!header) AliFatal("Not a standard AOD");
197     header->SetRunNumber(esd->GetRunNumber());
198     header->SetQTheta(dRP,1);
199         
200   }
201     
202   PostData(0,fEvent);
203   
204   
205
206
207 //________________________________________________________________________
208 void AliAnalysisTaskFlowEventforRP::Terminate(Option_t *) 
209 {
210   // Called once at the end of the query -- do not call in case of CAF
211   
212 }
213
214