]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/AliFMDEnergyFitterTask.cxx
corrected obvious typo. which value to chose between 0.045 and 0.050 is a bit matter...
[u/mrichter/AliRoot.git] / PWG2 / 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 }
48
49 //____________________________________________________________________
50 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask(const char* name)
51   : AliAnalysisTaskSE(name), 
52     fFirstEvent(true),
53     fEventInspector("event"),
54     fEnergyFitter("energy"),
55     fList(0),
56     fbLow(0),
57     fbHigh(100)
58 {
59   // 
60   // Constructor 
61   // 
62   // Parameters:
63   //    name Name of task 
64   //
65   DefineOutput(1, TList::Class());
66   DefineOutput(2, TList::Class());
67 }
68
69 //____________________________________________________________________
70 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask(const AliFMDEnergyFitterTask& o)
71   : AliAnalysisTaskSE(o),
72     fFirstEvent(o.fFirstEvent),
73     fEventInspector(o.fEventInspector),
74     fEnergyFitter(o.fEnergyFitter),
75     fList(o.fList),
76     fbLow(o.fbLow),
77     fbHigh(o.fbHigh)
78 {
79   // 
80   // Copy constructor 
81   // 
82   // Parameters:
83   //    o Object to copy from 
84   //
85   DefineOutput(1, TList::Class());
86   DefineOutput(2, TList::Class());
87 }
88
89 //____________________________________________________________________
90 AliFMDEnergyFitterTask&
91 AliFMDEnergyFitterTask::operator=(const AliFMDEnergyFitterTask& o)
92 {
93   // 
94   // Assignment operator 
95   // 
96   // Parameters:
97   //    o Object to assign from 
98   // 
99   // Return:
100   //    Reference to this object 
101   //
102   AliAnalysisTaskSE::operator=(o);
103
104   fFirstEvent        = o.fFirstEvent;
105   fEventInspector    = o.fEventInspector;
106   fEnergyFitter      = o.fEnergyFitter;
107   fList              = o.fList;
108   fbLow              = o.fbLow;
109   fbHigh             = o.fbHigh;
110
111   return *this;
112 }
113
114 //____________________________________________________________________
115 void
116 AliFMDEnergyFitterTask::SetDebug(Int_t dbg)
117 {
118   // 
119   // Set the debug level 
120   // 
121   // Parameters:
122   //    dbg Debug level
123   //
124   fEventInspector.SetDebug(dbg);
125   fEnergyFitter.SetDebug(dbg);
126 }
127
128 //____________________________________________________________________
129 void
130 AliFMDEnergyFitterTask::Init()
131 {
132   // 
133   // Initialize the task 
134   // 
135   //
136   fFirstEvent = true;
137 }
138
139 //____________________________________________________________________
140 void
141 AliFMDEnergyFitterTask::InitializeSubs()
142 {
143   // 
144   // Initialise the sub objects and stuff.  Called on first event 
145   // 
146   //
147   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
148   UShort_t sys = fEventInspector.GetCollisionSystem();
149   UShort_t sNN = fEventInspector.GetEnergy();
150   Short_t  fld = fEventInspector.GetField();
151   fcm.Init(sys, sNN, fld, 0);
152   TAxis eAxis(0,0,0);
153   TAxis vAxis(10,-10,10);
154   fEnergyFitter.Init(eAxis);
155   fEventInspector.Init(vAxis);
156
157 }
158
159 //____________________________________________________________________
160 void
161 AliFMDEnergyFitterTask::UserCreateOutputObjects()
162 {
163   // 
164   // Create output objects 
165   // 
166   //
167   fList = new TList;
168   fList->SetOwner();
169
170   fEventInspector.DefineOutput(fList);
171   fEnergyFitter.DefineOutput(fList);
172
173   PostData(1, fList);
174 }
175 //____________________________________________________________________
176 void
177 AliFMDEnergyFitterTask::UserExec(Option_t*)
178 {
179   // 
180   // Process each event 
181   // 
182   // Parameters:
183   //    option Not used
184   //  
185
186   // static Int_t cnt = 0;
187   // cnt++;
188   // Get the input data 
189   
190   AliMCEvent* mcevent = MCEvent();
191   if(mcevent) {
192     AliHeader* header            = mcevent->Header();
193     AliGenHijingEventHeader* hijingHeader = 
194       dynamic_cast<AliGenHijingEventHeader*>(header->GenEventHeader());
195     if(hijingHeader) {
196       Float_t b = hijingHeader->ImpactParameter();
197       if(b<fbLow || b>fbHigh) return;
198       else
199         std::cout<<"Selecting event with impact parameter "<<b<<std::endl;
200     }
201     
202   }
203   
204   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
205   // AliInfo(Form("Event # %6d (esd=%p)", cnt, esd));
206   if (!esd) { 
207     AliWarning("No ESD event found for input event");
208     return;
209   }
210
211   // On the first event, initialize the parameters 
212   if (fFirstEvent && esd->GetESDRun()) { 
213     fEventInspector.ReadRunDetails(esd);
214     
215     AliInfo(Form("Initializing with parameters from the ESD:\n"
216                  "         AliESDEvent::GetBeamEnergy()   ->%f\n"
217                  "         AliESDEvent::GetBeamType()     ->%s\n"
218                  "         AliESDEvent::GetCurrentL3()    ->%f\n"
219                  "         AliESDEvent::GetMagneticField()->%f\n"
220                  "         AliESDEvent::GetRunNumber()    ->%d\n",
221                  esd->GetBeamEnergy(), 
222                  esd->GetBeamType(),
223                  esd->GetCurrentL3(), 
224                  esd->GetMagneticField(),
225                  esd->GetRunNumber()));
226
227               
228
229     // AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
230     // pars->SetParametersFromESD(esd);
231     // pars->PrintStatus();
232     fFirstEvent = false;
233
234     InitializeSubs();
235   }
236   Bool_t   lowFlux   = kFALSE;
237   UInt_t   triggers  = 0;
238   UShort_t ivz       = 0;
239   Double_t vz        = 0;
240   Double_t cent      = 0;
241   UShort_t nClusters = 0;
242   UInt_t   found     = fEventInspector.Process(esd, triggers, lowFlux, 
243                                                ivz, vz, cent, nClusters);
244   if (found & AliFMDEventInspector::kNoEvent)    return;
245   if (found & AliFMDEventInspector::kNoTriggers) return;
246   if (found & AliFMDEventInspector::kNoSPD)     return;
247   if (found & AliFMDEventInspector::kNoFMD)     return;
248   if (found & AliFMDEventInspector::kNoVertex)  return;
249   if (found & AliFMDEventInspector::kBadVertex) return;
250   
251   //  if(cent > 0) {
252   //  if( cent < 40 || cent >50 ) return;
253   //  else std::cout<<"selecting event with cent "<<cent<<std::endl;
254   // }
255   
256   // Get FMD data 
257   AliESDFMD* esdFMD = esd->GetFMDData();
258   // Do the energy stuff 
259   if (!fEnergyFitter.Accumulate(*esdFMD, cent, 
260                                 triggers & AliAODForwardMult::kEmpty)){
261     AliWarning("Energy fitter failed");
262     return;
263   }
264   PostData(1, fList);
265 }
266
267 //____________________________________________________________________
268 void
269 AliFMDEnergyFitterTask::Terminate(Option_t*)
270 {
271   // 
272   // End of job
273   // 
274   // Parameters:
275   //    option Not used 
276   //
277   AliInfo(Form("Running terminate of %s", GetName()));
278   TList* list = dynamic_cast<TList*>(GetOutputData(1));
279   if (!list) {
280     AliError(Form("No output list defined (%p)", GetOutputData(1)));
281     if (GetOutputData(1)) GetOutputData(1)->Print();
282     return;
283   }
284   
285   AliInfo("Fitting energy loss spectra");
286   fEnergyFitter.Fit(list);
287
288   // Make a deep copy and post that as output 2 
289   TList* list2 = static_cast<TList*>(list->Clone(Form("%sResults", 
290                                                       list->GetName())));
291   list2->SetOwner();
292   PostData(2, list2);
293 }
294
295 //____________________________________________________________________
296 void
297 AliFMDEnergyFitterTask::Print(Option_t*) const
298 {
299   // 
300   // Print information 
301   // 
302   // Parameters:
303   //    option Not used
304   //
305 }
306
307 //
308 // EOF
309 //