]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDEnergyFitterTask.cxx
Major refactoring of the code.
[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 "AliESDFMD.h"
21 #include "AliMCEvent.h"
22 #include "AliAODForwardMult.h"
23 #include "AliAnalysisManager.h"
24 #include "AliForwardCorrectionManager.h"
25 #include <TH1.h>
26 #include <TDirectory.h>
27 #include <TTree.h>
28 #include <TFile.h>
29 #include <TROOT.h>
30 #include <iostream>
31
32 //====================================================================
33 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask()
34   : AliBaseESDTask(),
35     fEventInspector(),
36     fEnergyFitter(),
37     fbLow(0),
38     fbHigh(100)
39 {
40   // 
41   // Constructor
42   //
43   DGUARD(fDebug, 3,"Default CTOR of AliFMDEnergyFitterTask");
44   fCloneList = true;
45 }
46
47 //____________________________________________________________________
48 AliFMDEnergyFitterTask::AliFMDEnergyFitterTask(const char* name)
49   : AliBaseESDTask(name, "", &(AliForwardCorrectionManager::Instance())), 
50     fEventInspector("event"),
51     fEnergyFitter("energy"),
52     fbLow(0),
53     fbHigh(100)
54 {
55   // 
56   // Constructor 
57   // 
58   // Parameters:
59   //    name Name of task 
60   //
61   DGUARD(fDebug, 3,"Named CTOR of AliFMDEnergyFitterTask: %s", name);
62   fCloneList = true;
63 }
64
65
66 //____________________________________________________________________
67 void
68 AliFMDEnergyFitterTask::SetDebug(Int_t dbg)
69 {
70   // 
71   // Set the debug level 
72   // 
73   // Parameters:
74   //    dbg Debug level
75   //
76   AliBaseESDTask::SetDebug(dbg);
77   fEnergyFitter.SetDebug(dbg);
78 }
79 //____________________________________________________________________
80 TAxis*
81 AliFMDEnergyFitterTask::DefaultEtaAxis() const
82 {
83   static TAxis* a = new TAxis(0, 0, 0);
84   return a;
85 }
86 //____________________________________________________________________
87 TAxis*
88 AliFMDEnergyFitterTask::DefaultVertexAxis() const
89 {
90   static TAxis* a = new TAxis(10, -10, 10);
91   return a;
92 }
93
94 //____________________________________________________________________
95 Bool_t
96 AliFMDEnergyFitterTask::Book()
97 {
98   // 
99   // Create output objects 
100   // 
101   //
102   DGUARD(fDebug,1,"Create output objects of AliFMDEnergyFitterTask");
103
104   // We don't need any corrections for this task 
105   fNeededCorrections = 0;
106   fExtraCorrections  = 0;
107
108   fEnergyFitter.CreateOutputObjects(fList);
109
110   return true;
111 }
112 //____________________________________________________________________
113 Bool_t
114 AliFMDEnergyFitterTask::PreData(const TAxis& /*vertex*/, const TAxis& eta)
115 {
116   // 
117   // Initialise the sub objects and stuff.  Called on first event 
118   // 
119   //
120   DGUARD(fDebug,1,"Initialize subs of AliFMDEnergyFitterTask");
121
122   fEnergyFitter.SetupForData(eta);
123
124   return true;
125 }
126
127 //____________________________________________________________________
128 Bool_t
129 AliFMDEnergyFitterTask::Event(AliESDEvent& esd)
130 {
131   // 
132   // Process each event 
133   // 
134   // Parameters:
135   //    option Not used
136   //  
137
138   // static Int_t cnt = 0;
139   // cnt++;
140   // Get the input data 
141   DGUARD(fDebug,3,"Analyse event of AliFMDEnergyFitterTask");
142   // --- Read in the data --------------------------------------------
143   LoadBranches();
144
145   Bool_t   lowFlux   = kFALSE;
146   UInt_t   triggers  = 0;
147   UShort_t ivz       = 0;
148   TVector3 ip;
149   Double_t cent      = 0;
150   UShort_t nClusters = 0;
151   UInt_t   found     = fEventInspector.Process(&esd, triggers, lowFlux, 
152                                                ivz, ip, cent, nClusters);
153   if (found & AliFMDEventInspector::kNoEvent)    return false;
154   if (found & AliFMDEventInspector::kNoTriggers) return false;
155   if (found & AliFMDEventInspector::kNoSPD)      return false;
156   if (found & AliFMDEventInspector::kNoFMD)      return false;
157   if (found & AliFMDEventInspector::kNoVertex)   return false;
158   if (found & AliFMDEventInspector::kBadVertex)  return false;
159
160   // do not process pile-up, A, C, and E events 
161   if (triggers & AliAODForwardMult::kPileUp)     return false;
162   if (triggers & AliAODForwardMult::kA)          return false;
163   if (triggers & AliAODForwardMult::kC)          return false;
164   if (triggers & AliAODForwardMult::kE)          return false;
165   
166   // We want only the events found by off-line 
167   if (!(triggers & AliAODForwardMult::kOffline)) return false;
168
169   //  if(cent > 0) {
170   //  if( cent < 40 || cent >50 ) return;
171   //  else std::cout<<"selecting event with cent "<<cent<<std::endl;
172   // }
173   
174   // Get FMD data 
175   AliESDFMD* esdFMD = esd.GetFMDData();
176   // Do the energy stuff 
177   if (!fEnergyFitter.Accumulate(*esdFMD, cent, 
178                                 triggers & AliAODForwardMult::kEmpty)){
179     AliWarning("Energy fitter failed");
180     return false;
181   }
182
183   return true;
184 }
185
186 //____________________________________________________________________
187 Bool_t
188 AliFMDEnergyFitterTask::Finalize()
189 {
190   // 
191   // End of job
192   // 
193   // Parameters:
194   //    option Not used 
195   //
196   DGUARD(fDebug,1,"Processing merged output of AliFMDEnergyFitterTask");
197
198   AliInfo("Fitting energy loss spectra");
199   fEnergyFitter.Fit(fResults);
200
201   return true;
202 }
203
204 //____________________________________________________________________
205 void
206 AliFMDEnergyFitterTask::Print(Option_t* option) const
207 {
208   // 
209   // Print information 
210   // 
211   // Parameters:
212   //    option Not used
213   //
214   AliBaseESDTask::Print(option);
215   gROOT->IncreaseDirLevel();
216   fEnergyFitter.Print(option);
217   gROOT->DecreaseDirLevel();
218 }
219
220 //
221 // EOF
222 //