]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDEventPlaneTask.cxx
Renamed some member functions for more logical names
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDEventPlaneTask.cxx
1 //
2 // Calculate the FMD eventplane 
3 //
4 // Inputs:
5 //  - AliAODEvent
6 //
7 // Outputs:
8 //  - AnalysisResults.root
9 //
10 #include <TList.h>
11 #include <TMath.h>
12 #include "TH2D.h"
13 #include "AliLog.h"
14 #include "TAxis.h"
15 #include "AliAnalysisManager.h"
16 #include "AliFMDEventPlaneTask.h"
17 #include "AliAODHandler.h"
18 #include "AliAODInputHandler.h"
19 #include "AliAODForwardMult.h"
20 #include "AliAODEvent.h"
21 #include "AliAODForwardEP.h"
22
23 ClassImp(AliFMDEventPlaneTask)
24 #if 0
25 ; // For emacs 
26 #endif
27
28 AliFMDEventPlaneTask::AliFMDEventPlaneTask()
29   : AliAnalysisTaskSE(),
30     fSumList(0),         // Sum list
31     fOutputList(0),      // Output list
32     fAOD(0),             // AOD input event
33     fMC(0),              // MC flag
34     fEventPlaneFinder(), // EP finder
35     fZvertex(1111),      // Z vertex 
36     fCent(-1),           // Centrality
37     fHistCent(),         // Diagnostics histogram
38     fHistVertexSel(),    // Diagnostics histogram
39     fHistVertexAll()     // Diagnostics histogram
40 {
41   // 
42   // Default constructor
43   //
44   DGUARD(fDebug, 3,"Default CTOR of AliFMDEventPlaneTask");
45 }
46 //_____________________________________________________________________
47 AliFMDEventPlaneTask::AliFMDEventPlaneTask(const char* name) 
48   : AliAnalysisTaskSE(name),
49     fSumList(0),                     // Sum list
50     fOutputList(0),                  // Output list
51     fAOD(0),                         // AOD input event
52     fMC(0),                          // MC flag
53     fEventPlaneFinder("eventPlane"), // EP finder
54     fZvertex(1111),                  // Z vertex 
55     fCent(-1),                       // Centrality
56     fHistCent(0),                    // Diagnostics histogram
57     fHistVertexSel(0),               // Diagnostics histogram
58     fHistVertexAll(0)                // Diagnostics histogram
59
60 {
61   // 
62   // Constructor
63   //
64   // Parameters:
65   //  name: Name of task
66   //
67   DGUARD(fDebug, 3,"Named CTOR of AliFMDEventPlaneTask: %s", name);
68
69   DefineOutput(1, TList::Class());
70   DefineOutput(2, TList::Class());
71
72 }
73 //_____________________________________________________________________
74 AliFMDEventPlaneTask::AliFMDEventPlaneTask(const AliFMDEventPlaneTask& o)
75   : AliAnalysisTaskSE(o),
76     fSumList(o.fSumList),                    // Sumlist
77     fOutputList(o.fOutputList),              // Output list
78     fAOD(o.fAOD),                            // AOD input event
79     fMC(o.fMC),                              // MC flag
80     fEventPlaneFinder(o.fEventPlaneFinder),  // EP finder
81     fZvertex(o.fZvertex),                    // Z vertex 
82     fCent(o.fCent),                          // Centrality
83     fHistCent(o.fHistCent),                  // Diagnostics histogram
84     fHistVertexSel(o.fHistVertexSel),        // Diagnostics histogram
85     fHistVertexAll(o.fHistVertexAll)         // Diagnostics histogram
86 {
87   // 
88   // Copy constructor 
89   // 
90   // Parameters:
91   //    o Object to copy from 
92   //
93   DGUARD(fDebug, 3,"Copy CTOR of AliFMDEventPlaneTask");
94 }
95 //_____________________________________________________________________
96 AliFMDEventPlaneTask&
97 AliFMDEventPlaneTask::operator=(const AliFMDEventPlaneTask& o)
98 {
99   // 
100   // Assignment operator 
101   //
102   DGUARD(fDebug,3,"Assignment of AliFMDEventPlaneTask");
103   if (&o == this) return *this;
104   fSumList           = o.fSumList;
105   fOutputList        = o.fOutputList;
106   fAOD               = o.fAOD;
107   fMC                = o.fMC;
108   fEventPlaneFinder  = o.fEventPlaneFinder;
109   fZvertex           = o.fZvertex;
110   fCent              = o.fCent;
111   fHistCent          = o.fHistCent;
112   fHistVertexSel     = o.fHistVertexSel;
113   fHistVertexAll     = o.fHistVertexAll;
114
115   return *this;
116 }
117 //_____________________________________________________________________
118 void AliFMDEventPlaneTask::UserCreateOutputObjects()
119 {
120   //
121   // Create output objects
122   //
123   DGUARD(fDebug,1,"Create user objects of AliFMDEventPlaneTask");
124   if (!fSumList)
125     fSumList = new TList();
126   fSumList->SetName("Sums");
127   fSumList->SetOwner();
128
129   // Diagnostics histograms
130   fHistCent          = new TH1D("hCent", "Centralities", 100, 0, 100);
131   fHistVertexSel     = new TH1D("hVertexSel", "Selectec vertices", 40, -20, 20);
132   fHistVertexAll     = new TH1D("hVertexAll", "All vertices", 40, -20, 20);
133
134   fSumList->Add(fHistCent);
135   fSumList->Add(fHistVertexSel);
136   fSumList->Add(fHistVertexAll);
137
138   // Init of EventPlaneFinder
139   TAxis* pe = new TAxis(200, -4., 6.);
140   fEventPlaneFinder.CreateOutputObjects(fSumList);
141   fEventPlaneFinder.SetupForData(*pe);
142
143   PostData(1, fSumList);
144
145 }
146 //_____________________________________________________________________
147 void AliFMDEventPlaneTask::UserExec(Option_t */*option*/)
148 {
149   // 
150   // Called each event
151   //
152   // Parameters:
153   //  option: Not used
154   //
155   DGUARD(fDebug,1,"Process an event in AliFMDEventPlaneTask");
156
157   // Reset data members
158   fCent = -1;
159   fZvertex = 1111;
160
161   // Get input event
162   fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
163   if (!fAOD) return;
164
165   AliAODForwardMult* aodfmult = 
166     static_cast<AliAODForwardMult*>(fAOD->FindListObject("Forward"));
167
168   if (!aodfmult) return;
169   if (!AODCheck(aodfmult)) return;
170
171   if (fAOD->GetRunNumber() != fEventPlaneFinder.GetRunNumber())
172     fEventPlaneFinder.SetRunNumber(fAOD->GetRunNumber());
173
174   AliAODForwardEP aodep;
175   TH2D fmdHist = aodfmult->GetHistogram();
176
177   fEventPlaneFinder.FindEventplane(fAOD, aodep, &fmdHist, 0);
178
179   PostData(1, fSumList);
180
181 }
182 //_____________________________________________________________________
183 void AliFMDEventPlaneTask::Terminate(Option_t */*option*/)
184 {
185   //
186   // Terminate - Called after all events
187   //
188   // Parameters:
189   //  option: Not used
190   //
191   DGUARD(fDebug,1,"Process merged output of AliFMDEventPlaneTask");
192
193   // Reinitiate lists if Terminate is called separately!
194   fSumList = dynamic_cast<TList*> (GetOutputData(1));
195   if(!fSumList) {
196     AliError("Could not retrieve TList fSumList"); 
197     return; 
198   }
199  
200   if (!fOutputList)
201     fOutputList = new TList();
202   fOutputList->SetName("Results");
203   fOutputList->SetOwner();
204
205  // Calculations can be done here: Currently there are none
206
207   PostData(2, fOutputList);
208
209 }
210 // _____________________________________________________________________
211 Bool_t AliFMDEventPlaneTask::AODCheck(const AliAODForwardMult* aodfm) 
212 {
213   // 
214   // Function to check that and AOD event meets the cuts
215   //
216   // Parameters: 
217   //  AliAODForwardMult: forward mult object with trigger and vertex info
218   //
219   DGUARD(fDebug,2,"Check AOD in AliFMDEventPlaneTask");
220
221   if (!aodfm->IsTriggerBits(AliAODForwardMult::kOffline)) return kFALSE;
222
223   fCent = (Double_t)aodfm->GetCentrality();
224   if (0. >= fCent || fCent >= 80.) return kFALSE;
225   fHistCent->Fill(fCent);
226
227   fZvertex = aodfm->GetIpZ();
228   fHistVertexAll->Fill(fZvertex);
229   if (TMath::Abs(fZvertex) >= 10.) return kFALSE;
230   fHistVertexSel->Fill(fZvertex);
231
232   return kTRUE;
233
234 }
235 //_____________________________________________________________________
236 //
237 //
238 // EOF