]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowTasks/AliAnalysisTaskFlowEvent.cxx
973af3176f65daa5bc1ab6093b042c25d9979c96
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowTasks / AliAnalysisTaskFlowEvent.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 // AliAnalysisTaskFlowEvent:
18 //
19 // analysis task for filling the flow event
20 // from MCEvent, ESD, AOD ....
21 // and put it in an output stream so it can
22 // be used by the various flow analysis methods
23 // for cuts the correction framework is used
24 // which also outputs QA histograms to view
25 // the effects of the cuts
26 ////////////////////////////////////////////////////
27
28 #include "Riostream.h" //needed as include
29 #include "TChain.h"
30 #include "TTree.h"
31 #include "TFile.h" //needed as include
32 #include "TList.h"
33 #include "TRandom3.h"
34 #include "TTimeStamp.h"
35
36 // ALICE Analysis Framework
37 #include "AliAnalysisManager.h"
38 #include "AliAnalysisTaskSE.h"
39
40 // ESD interface
41 #include "AliESDEvent.h"
42 #include "AliESDInputHandler.h"
43
44 // AOD interface
45 #include "AliAODEvent.h"
46 #include "AliAODInputHandler.h"
47
48 // Monte Carlo Event
49 #include "AliMCEventHandler.h"
50 #include "AliMCEvent.h"
51
52 // ALICE Correction Framework
53 #include "AliCFManager.h"
54
55 // Interface to Event generators to get Reaction Plane Angle
56 #include "AliGenCocktailEventHeader.h"
57 #include "AliGenHijingEventHeader.h"
58 #include "AliGenGeVSimEventHeader.h"
59 #include "AliGenEposEventHeader.h"
60
61 // Interface to make the Flow Event Simple used in the flow analysis methods
62 #include "AliFlowEvent.h"
63
64 #include "AliAnalysisTaskFlowEvent.h"
65
66 ClassImp(AliAnalysisTaskFlowEvent)
67
68 //________________________________________________________________________
69 AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent() :
70   AliAnalysisTaskSE(),
71   //  fOutputFile(NULL),
72   fAnalysisType("ESD"),
73   fCFManager1(NULL),
74   fCFManager2(NULL),
75   fQAInt(NULL),
76   fQADiff(NULL),
77   fMinMult(0),
78   fMaxMult(10000000),
79   fMinA(-1.0),
80   fMaxA(-0.01),
81   fMinB(0.01),
82   fMaxB(1.0),
83   fQA(kFALSE),
84   fMCReactionPlaneAngle(0.),
85   fCount(0),
86   fNoOfLoops(1),
87   fEllipticFlowValue(0.),
88   fSigmaEllipticFlowValue(0.),
89   fMultiplicityOfEvent(1000000000),
90   fSigmaMultiplicityOfEvent(0),
91   fMyTRandom3(NULL),
92   fbAfterburnerOn(kFALSE)
93 {
94   // Constructor
95   cout<<"AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent()"<<endl;
96 }
97
98 //________________________________________________________________________
99 AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent(const char *name, Bool_t on, UInt_t iseed) :
100   AliAnalysisTaskSE(name),
101   //  fOutputFile(NULL),
102   fAnalysisType("ESD"),
103   fCFManager1(NULL),
104   fCFManager2(NULL),
105   fQAInt(NULL),
106   fQADiff(NULL),
107   fMinMult(0),
108   fMaxMult(10000000),
109   fMinA(-1.0),
110   fMaxA(-0.01),
111   fMinB(0.01),
112   fMaxB(1.0),
113   fQA(on),
114   fMCReactionPlaneAngle(0.),
115   fCount(0),
116   fNoOfLoops(1),
117   fEllipticFlowValue(0.),
118   fSigmaEllipticFlowValue(0.),
119   fMultiplicityOfEvent(1000000000),
120   fSigmaMultiplicityOfEvent(0),
121   fMyTRandom3(NULL),
122   fbAfterburnerOn(kFALSE)
123 {
124   // Constructor
125   cout<<"AliAnalysisTaskFlowEvent::AliAnalysisTaskFlowEvent(const char *name, Bool_t on, UInt_t iseed)"<<endl;
126   fMyTRandom3 = new TRandom3(iseed);
127   gRandom->SetSeed(fMyTRandom3->Integer(65539));
128
129
130   // Define output slots here
131   // Define here the flow event output
132   DefineOutput(1, AliFlowEventSimple::Class());
133   if(on)
134   {
135     DefineOutput(2, TList::Class());
136     DefineOutput(3, TList::Class());
137   }
138   // and for testing open an output file
139   //  fOutputFile = new TFile("FlowEvents.root","RECREATE");
140
141 }
142
143 //________________________________________________________________________
144 AliAnalysisTaskFlowEvent::~AliAnalysisTaskFlowEvent()
145 {
146   //
147   // Destructor
148   //
149   if (fMyTRandom3) delete fMyTRandom3;
150   // objects in the output list are deleted
151   // by the TSelector dtor (I hope)
152
153 }
154
155 //________________________________________________________________________
156 void AliAnalysisTaskFlowEvent::UserCreateOutputObjects()
157 {
158   // Called at every worker node to initialize
159   cout<<"AliAnalysisTaskFlowEvent::CreateOutputObjects()"<<endl;
160
161   if (!(fAnalysisType == "AOD" || fAnalysisType == "ESD" || fAnalysisType == "ESDMCkineESD"  || fAnalysisType == "ESDMCkineMC" || fAnalysisType == "MC"))
162   {
163     cout<<"WRONG ANALYSIS TYPE! only ESD, ESDMCkineESD, ESDMCkineMC, AOD and MC are allowed."<<endl;
164     exit(1);
165   }
166 }
167
168 //________________________________________________________________________
169 void AliAnalysisTaskFlowEvent::UserExec(Option_t *)
170 {
171   // Main loop
172   // Called for each event
173   AliFlowEvent* flowEvent = NULL;
174   AliMCEvent*  mcEvent = MCEvent();                              // from TaskSE
175   AliESDEvent* myESD = dynamic_cast<AliESDEvent*>(InputEvent()); // from TaskSE
176   AliAODEvent* myAOD = dynamic_cast<AliAODEvent*>(InputEvent()); // from TaskSE
177
178   // Make the FlowEvent for MC input
179   if (fAnalysisType == "MC")
180   {
181     // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
182     // This handler can return the current MC event
183     if (!(fCFManager1&&fCFManager2))
184     {
185       cout << "ERROR: No pointer to correction framework cuts! " << endl;
186       return;
187     }
188     if (!mcEvent)
189     {
190       Printf("ERROR: Could not retrieve MC event");
191       return;
192     }
193
194     fCFManager1->SetMCEventInfo(mcEvent);
195     fCFManager2->SetMCEventInfo(mcEvent);
196     
197     Printf("Number of MC particles: %d", mcEvent->GetNumberOfTracks());
198     // analysis
199     flowEvent = new AliFlowEvent(mcEvent,fCFManager1,fCFManager2);
200   }
201   // Make the FlowEvent for ESD input
202   else if (fAnalysisType == "ESD")
203   {
204     if (!(fCFManager1&&fCFManager2))
205     {
206       cout << "ERROR: No pointer to correction framework cuts! " << endl;
207       return;
208     }
209     if (!myESD)
210     {
211       Printf("ERROR: ESD not available");
212       return;
213     }
214     //check the offline trigger (check if the event has the correct trigger)
215     Printf("There are %d tracks in this event", fInputEvent->GetNumberOfTracks());
216     // analysis
217     flowEvent = new AliFlowEvent(myESD,fCFManager1,fCFManager2);
218   }
219   // Make the FlowEvent for ESD input combined with MC info
220   else if (fAnalysisType == "ESDMCkineESD" || fAnalysisType == "ESDMCkineMC" )
221   {
222     if (!(fCFManager1&&fCFManager2))
223     {
224       cout << "ERROR: No pointer to correction framework cuts! " << endl;
225       return;
226     }
227     if (!myESD)
228     {
229       Printf("ERROR: ESD not available");
230       return;
231     }
232     Printf("There are %d tracks in this event", fInputEvent->GetNumberOfTracks());
233
234     if (!mcEvent)
235     {
236       Printf("ERROR: Could not retrieve MC event");
237       return;
238     }
239
240     fCFManager1->SetMCEventInfo(mcEvent);
241     fCFManager2->SetMCEventInfo(mcEvent);
242
243     if (fAnalysisType == "ESDMCkineESD")
244     {
245       flowEvent = new AliFlowEvent(myESD, mcEvent, AliFlowEvent::kESDkine, fCFManager1, fCFManager2 );
246     }
247     else if (fAnalysisType == "ESDMCkineMC")
248     {
249       flowEvent = new AliFlowEvent(myESD, mcEvent, AliFlowEvent::kMCkine, fCFManager1, fCFManager2 );
250     }
251   }
252   // Make the FlowEventSimple for AOD input
253   else if (fAnalysisType == "AOD")
254   {
255     if (!myAOD)
256     {
257       Printf("ERROR: AOD not available");
258       return;
259     }
260     Printf("There are %d tracks in this event", myAOD->GetNumberOfTracks());
261     flowEvent = new AliFlowEvent(myAOD);
262   }
263
264   //check event cuts
265   Int_t mult = flowEvent->NumberOfTracks();
266   if (mult<fMinMult && mult>fMaxMult) return;
267
268   //tag subevents
269   flowEvent->TagSubeventsInEta(fMinA,fMaxA,fMinB,fMaxB);
270
271   ////TODO afterburner
272   //if (fbAfterburnerOn && fMyTRandom3) {
273   //  // set the new value of the values using a after burner
274   //  cout << "settings for afterburner in TaskFlowEvent.cxx:" << endl;
275   //  cout << "fCount = " << fCount << endl;
276   //  cout << "fNoOfLoops = " << fNoOfLoops << endl;
277   //  cout << "fEllipticFlowValue = " << fEllipticFlowValue << endl;
278   //  cout << "fSigmaEllipticFlowValue = " << fSigmaEllipticFlowValue << endl;
279   //  cout << "fMultiplicityOflowEvent = " << fMultiplicityOflowEvent << endl;
280   //  cout << "fSigmaMultiplicityOflowEvent = " << fSigmaMultiplicityOflowEvent << endl;
281   //  Double_t xRPangle;
282   //  if (!mcEvent) { xRPangle = TMath::TwoPi()*(fMyTRandom3->Rndm()); }
283   //  else { xRPangle = fRP; }
284   //  Double_t xNewFlowValue = fMyTRandom3->Gaus(fEllipticFlowValue,fSigmaEllipticFlowValue);
285   //  Int_t nNewMultOflowEvent =  TMath::Nint(fMyTRandom3->Gaus(fMultiplicityOflowEvent,fSigmaMultiplicityOflowEvent));
286   //  cout << "xRPangle = " << xRPangle << endl;
287   //  cout << "xNewFlowValue = " << xNewFlowValue << endl;
288   //  cout << "nNewMultOflowEvent = " << nNewMultOflowEvent << endl;
289   //  cout << "settings for after burner" << endl;
290   //  flowEventMaker->SetMCReactionPlaneAngle(xRPangle);
291   //  flowEventMaker->SetNoOfLoops(fNoOfLoops);
292   //  flowEventMaker->SetEllipticFlowValue(xNewFlowValue);
293   //  flowEventMaker->SetMultiplicityOflowEvent(nNewMultOflowEvent);
294   //  //end settings afterburner
295   //}
296
297   // if monte carlo event get reaction plane from monte carlo (depends on generator)
298   if (mcEvent && mcEvent->GenEventHeader()) flowEvent->SetMCReactionPlaneAngle(mcEvent);
299
300   //fListHistos->Print();
301   //  fOutputFile->WriteObject(flowEvent,"myFlowEventSimple");
302   PostData(1,flowEvent);
303   if (fQA)
304   {
305     PostData(2,fQAInt);
306     PostData(3,fQADiff);
307   }
308 }
309
310 //________________________________________________________________________
311 void AliAnalysisTaskFlowEvent::Terminate(Option_t *)
312 {
313   // Called once at the end of the query -- do not call in case of CAF
314 }
315
316