]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDEnergyFitterTask.cxx
Many changes in one.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDEnergyFitterTask.cxx
1 // 
2 // Histogram and fit the energy loss distributions for the FMD
3 // 
4 // Inputs: 
5 //   - AliESDEvent 
6 //
7 // Outputs: 
8 //   - None
9 // 
10 // Histograms:
11 //   
12 // Corrections used:
13 //   - None
14 // 
15 // 
16 //
17 #include "AliFMDEnergyFitterTask.h"
18 #include "AliLog.h"
19 #include "AliESDEvent.h"
20 #include "AliAODForwardMult.h"
21 #include "AliForwardCorrectionManager.h"
22 #include "AliAnalysisManager.h"
23 #include "AliAnalysisDataSlot.h"
24 #include "AliAnalysisDataContainer.h"
25 #include <TH1.h>
26 #include <TDirectory.h>
27 #include <TTree.h>
28 #include <TFile.h>
29 #include "AliMCEvent.h"
30 #include "AliGenHijingEventHeader.h"
31 #include "AliHeader.h"
32 #include <iostream>
33
34 //====================================================================
35 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask()
36   : AliAnalysisTaskSE(),
37     fFirstEvent(true),
38     fEventInspector(),
39     fEnergyFitter(),
40     fList(0),
41     fbLow(0),
42     fbHigh(100)
43 {
44   // 
45   // Constructor
46   //
47   DGUARD(fDebug, 3,"Default CTOR of AliFMDEnergyFitterTask");
48 }
49
50 //____________________________________________________________________
51 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask(const char* name)
52   : AliAnalysisTaskSE(name), 
53     fFirstEvent(true),
54     fEventInspector("event"),
55     fEnergyFitter("energy"),
56     fList(0),
57     fbLow(0),
58     fbHigh(100)
59 {
60   // 
61   // Constructor 
62   // 
63   // Parameters:
64   //    name Name of task 
65   //
66   DGUARD(fDebug, 3,"Named CTOR of AliFMDEnergyFitterTask: %s", name);
67   DefineOutput(1, TList::Class());
68   DefineOutput(2, TList::Class());
69   fBranchNames = 
70     "ESD:AliESDRun.,AliESDHeader.,AliMultiplicity.,"
71     "AliESDFMD.,SPDVertex.,PrimaryVertex.";
72 }
73
74 //____________________________________________________________________
75 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask(const AliFMDEnergyFitterTask& o)
76   : AliAnalysisTaskSE(o),
77     fFirstEvent(o.fFirstEvent),
78     fEventInspector(o.fEventInspector),
79     fEnergyFitter(o.fEnergyFitter),
80     fList(o.fList),
81     fbLow(o.fbLow),
82     fbHigh(o.fbHigh)
83 {
84   // 
85   // Copy constructor 
86   // 
87   // Parameters:
88   //    o Object to copy from 
89   //
90   DGUARD(fDebug, 3,"COPY CTOR of AliFMDEnergyFitterTask");
91   DefineOutput(1, TList::Class());
92   DefineOutput(2, TList::Class());
93 }
94
95 //____________________________________________________________________
96 AliFMDEnergyFitterTask&
97 AliFMDEnergyFitterTask::operator=(const AliFMDEnergyFitterTask& o)
98 {
99   // 
100   // Assignment operator 
101   // 
102   // Parameters:
103   //    o Object to assign from 
104   // 
105   // Return:
106   //    Reference to this object 
107   //
108   fDebug = o.fDebug;
109   DGUARD(fDebug,3,"Assignment of AliFMDEnergyFitterTask");
110   if (&o == this) return *this; 
111   AliAnalysisTaskSE::operator=(o);
112
113   fFirstEvent        = o.fFirstEvent;
114   fEventInspector    = o.fEventInspector;
115   fEnergyFitter      = o.fEnergyFitter;
116   fList              = o.fList;
117   fbLow              = o.fbLow;
118   fbHigh             = o.fbHigh;
119
120   return *this;
121 }
122
123 //____________________________________________________________________
124 void
125 AliFMDEnergyFitterTask::SetDebug(Int_t dbg)
126 {
127   // 
128   // Set the debug level 
129   // 
130   // Parameters:
131   //    dbg Debug level
132   //
133   fDebug = dbg;
134   fEventInspector.SetDebug(dbg);
135   fEnergyFitter.SetDebug(dbg);
136 }
137
138 //____________________________________________________________________
139 void
140 AliFMDEnergyFitterTask::Init()
141 {
142   // 
143   // Initialize the task 
144   // 
145   //
146   DGUARD(fDebug,1,"Initialize of AliFMDEnergyFitterTask");
147   fFirstEvent = true;
148 }
149
150 //____________________________________________________________________
151 void
152 AliFMDEnergyFitterTask::InitializeSubs()
153 {
154   // 
155   // Initialise the sub objects and stuff.  Called on first event 
156   // 
157   //
158   DGUARD(fDebug,1,"Initialize subs of AliFMDEnergyFitterTask");
159   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
160   UShort_t sys = fEventInspector.GetCollisionSystem();
161   UShort_t sNN = fEventInspector.GetEnergy();
162   Short_t  fld = fEventInspector.GetField();
163   fcm.Init(sys, sNN, fld, 0);
164   TAxis eAxis(0,0,0);
165   TAxis vAxis(10,-10,10);
166   fEnergyFitter.Init(eAxis);
167   fEventInspector.Init(vAxis);
168
169 }
170
171 //____________________________________________________________________
172 void
173 AliFMDEnergyFitterTask::UserCreateOutputObjects()
174 {
175   // 
176   // Create output objects 
177   // 
178   //
179   DGUARD(fDebug,1,"Create output objects of AliFMDEnergyFitterTask");
180   fList = new TList;
181   fList->SetOwner();
182
183   fEventInspector.DefineOutput(fList);
184   fEnergyFitter.DefineOutput(fList);
185
186   PostData(1, fList);
187 }
188 //____________________________________________________________________
189 void
190 AliFMDEnergyFitterTask::UserExec(Option_t*)
191 {
192   // 
193   // Process each event 
194   // 
195   // Parameters:
196   //    option Not used
197   //  
198
199   // static Int_t cnt = 0;
200   // cnt++;
201   // Get the input data 
202   DGUARD(fDebug,3,"Analyse event of AliFMDEnergyFitterTask");
203   
204   AliMCEvent* mcevent = MCEvent();
205   if(mcevent) {
206     AliHeader* header            = mcevent->Header();
207     AliGenHijingEventHeader* hijingHeader = 
208       dynamic_cast<AliGenHijingEventHeader*>(header->GenEventHeader());
209     if(hijingHeader) {
210       Float_t b = hijingHeader->ImpactParameter();
211       if(b<fbLow || b>fbHigh) return;
212       else
213         std::cout<<"Selecting event with impact parameter "<<b<<std::endl;
214     }
215     
216   }
217   
218   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
219   // AliInfo(Form("Event # %6d (esd=%p)", cnt, esd));
220   if (!esd) { 
221     AliWarning("No ESD event found for input event");
222     return;
223   }
224
225   // On the first event, initialize the parameters 
226   if (fFirstEvent && esd->GetESDRun()) { 
227     fEventInspector.ReadRunDetails(esd);
228     
229     AliInfo(Form("Initializing with parameters from the ESD:\n"
230                  "         AliESDEvent::GetBeamEnergy()   ->%f\n"
231                  "         AliESDEvent::GetBeamType()     ->%s\n"
232                  "         AliESDEvent::GetCurrentL3()    ->%f\n"
233                  "         AliESDEvent::GetMagneticField()->%f\n"
234                  "         AliESDEvent::GetRunNumber()    ->%d\n",
235                  esd->GetBeamEnergy(), 
236                  esd->GetBeamType(),
237                  esd->GetCurrentL3(), 
238                  esd->GetMagneticField(),
239                  esd->GetRunNumber()));
240
241               
242
243     // AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
244     // pars->SetParametersFromESD(esd);
245     // pars->PrintStatus();
246     fFirstEvent = false;
247
248     InitializeSubs();
249   }
250   Bool_t   lowFlux   = kFALSE;
251   UInt_t   triggers  = 0;
252   UShort_t ivz       = 0;
253   TVector3 ip;
254   Double_t cent      = 0;
255   UShort_t nClusters = 0;
256   UInt_t   found     = fEventInspector.Process(esd, triggers, lowFlux, 
257                                                ivz, ip, cent, nClusters);
258   if (found & AliFMDEventInspector::kNoEvent)    return;
259   if (found & AliFMDEventInspector::kNoTriggers) return;
260   if (found & AliFMDEventInspector::kNoSPD)     return;
261   if (found & AliFMDEventInspector::kNoFMD)     return;
262   if (found & AliFMDEventInspector::kNoVertex)  return;
263   if (found & AliFMDEventInspector::kBadVertex) return;
264   
265   //  if(cent > 0) {
266   //  if( cent < 40 || cent >50 ) return;
267   //  else std::cout<<"selecting event with cent "<<cent<<std::endl;
268   // }
269   
270   // Get FMD data 
271   AliESDFMD* esdFMD = esd->GetFMDData();
272   // Do the energy stuff 
273   if (!fEnergyFitter.Accumulate(*esdFMD, cent, 
274                                 triggers & AliAODForwardMult::kEmpty)){
275     AliWarning("Energy fitter failed");
276     return;
277   }
278   PostData(1, fList);
279 }
280
281 //____________________________________________________________________
282 void
283 AliFMDEnergyFitterTask::FinishTaskOutput()
284 {
285   if (!fList) { 
286     Warning("FinishTaskOutput", "No list defined");
287   }
288   // else {
289   //   fList->ls();
290   // }
291   // if (fDebug)
292   // gDebug = 1;
293   AliAnalysisTaskSE::FinishTaskOutput();
294 }
295
296 //____________________________________________________________________
297 void
298 AliFMDEnergyFitterTask::Terminate(Option_t*)
299 {
300   // 
301   // End of job
302   // 
303   // Parameters:
304   //    option Not used 
305   //
306   DGUARD(fDebug,1,"Processing merged output of AliFMDEnergyFitterTask");
307
308   TList* list = dynamic_cast<TList*>(GetOutputData(1));
309   if (!list) {
310     AliError(Form("No output list defined (%p)", GetOutputData(1)));
311     if (GetOutputData(1)) GetOutputData(1)->Print();
312     return;
313   }
314   
315   AliInfo("Fitting energy loss spectra");
316   fEnergyFitter.Fit(list);
317
318   // Make a deep copy and post that as output 2 
319   TList* list2 = static_cast<TList*>(list->Clone(Form("%sResults", 
320                                                       list->GetName())));
321   list2->SetOwner();
322   PostData(2, list2);
323 }
324
325 //____________________________________________________________________
326 void
327 AliFMDEnergyFitterTask::Print(Option_t*) const
328 {
329   // 
330   // Print information 
331   // 
332   // Parameters:
333   //    option Not used
334   //
335 }
336
337 //
338 // EOF
339 //