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