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