]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliForwardQATask.cxx
Major refactoring of the code.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliForwardQATask.cxx
1 // 
2 // Calculate the multiplicity in the forward regions event-by-event 
3 // 
4 // Inputs: 
5 //   - AliESDEvent 
6 //
7 // Outputs: 
8 //   - AliAODForwardMult 
9 // 
10 // Histograms 
11 //   
12 // Corrections used 
13 //
14 #include "AliForwardQATask.h"
15 #include "AliForwardUtil.h"
16 #include "AliTriggerAnalysis.h"
17 #include "AliPhysicsSelection.h"
18 #include "AliLog.h"
19 #include "AliESDEvent.h"
20 #include "AliAODHandler.h"
21 #include "AliMultiplicity.h"
22 #include "AliInputEventHandler.h"
23 #include "AliForwardCorrectionManager.h"
24 #include "AliAnalysisManager.h"
25 #include "AliAODForwardMult.h"
26 #include <TH1.h>
27 #include <TDirectory.h>
28 #include <TTree.h>
29 #include <TROOT.h>
30 #include <TStopwatch.h>
31
32 //====================================================================
33 AliForwardQATask::AliForwardQATask()
34   : AliBaseESDTask(),
35     fEnableLowFlux(false), 
36     fESDFMD(),
37     fHistos(),
38     fEventInspector(),
39     fEnergyFitter(),
40     fSharingFilter(),
41     fDensityCalculator()
42 {
43   // 
44   // Constructor
45   //
46   fCloneList = true;
47 }
48
49 //____________________________________________________________________
50 AliForwardQATask::AliForwardQATask(const char* name)
51   : AliBaseESDTask(name, "", &(AliForwardCorrectionManager::Instance())),
52     fEnableLowFlux(false), 
53     fESDFMD(),
54     fHistos(),
55     fEventInspector("event"),
56     fEnergyFitter("energy"),
57     fSharingFilter("sharing"), 
58     fDensityCalculator("density")
59 {
60   // 
61   // Constructor 
62   // 
63   // Parameters:
64   //    name Name of task 
65   //
66   fEnergyFitter.SetNParticles(1); // Just find the 1st peak 
67   fEnergyFitter.SetDoMakeObject(false); 
68   fEnergyFitter.SetUseIncreasingBins(true);
69   fEnergyFitter.SetDoFits(kTRUE);
70   fEnergyFitter.SetLowCut(0.4);
71   fEnergyFitter.SetFitRangeBinWidth(4);
72   fEnergyFitter.SetMinEntries(1000);
73   fCloneList = true;
74 }
75 //____________________________________________________________________
76 void
77 AliForwardQATask::SetDebug(Int_t dbg)
78 {
79   // 
80   // Set debug level 
81   // 
82   // Parameters:
83   //    dbg Debug level
84   //
85   AliBaseESDTask::        SetDebug(dbg);
86   GetEnergyFitter()     .SetDebug(dbg);
87   GetSharingFilter()    .SetDebug(dbg);
88   GetDensityCalculator().SetDebug(dbg);
89 }
90
91 //____________________________________________________________________
92 TAxis*
93 AliForwardQATask::DefaultEtaAxis() const
94 {
95   static TAxis* a = new TAxis(240, -6, 6);
96   return a;
97 }
98 //____________________________________________________________________
99 TAxis*
100 AliForwardQATask::DefaultVertexAxis() const
101 {
102   static TAxis* a = new TAxis(10, -10, 10);
103   return a;
104 }
105
106 //____________________________________________________________________
107 Bool_t
108 AliForwardQATask::Book()
109 {
110   // 
111   // Create output objects 
112   // 
113   //
114   UInt_t what = AliForwardCorrectionManager::kAll;
115   what ^= AliForwardCorrectionManager::kDoubleHit;
116   what ^= AliForwardCorrectionManager::kVertexBias;
117   what ^= AliForwardCorrectionManager::kAcceptance;
118   what ^= AliForwardCorrectionManager::kMergingEfficiency;
119   what ^= AliForwardCorrectionManager::kELossFits;
120   fNeededCorrections = what;
121   fExtraCorrections  = AliForwardCorrectionManager::kELossFits;
122   
123   fEnergyFitter.CreateOutputObjects(fList);
124   fSharingFilter.CreateOutputObjects(fList);
125   fDensityCalculator.CreateOutputObjects(fList);
126   
127   return true;
128 }
129 //____________________________________________________________________
130 Bool_t
131 AliForwardQATask::PreData(const TAxis& /*vertex*/, const TAxis& eta)
132 {
133   // 
134   // Initialise the sub objects and stuff.  Called on first event 
135   // 
136   //
137   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
138   // We allow fall-back queries so that we may proceed in case we have no 
139   // valid corrections 
140   if (!fcm.GetELossFit()) { 
141     AliWarning("No energy loss fits");
142     
143     // Fall-back values if we do not have the energy loss fits 
144     AliFMDMultCuts& sfLCuts = GetSharingFilter().GetLCuts();
145     if (sfLCuts.GetMethod() != AliFMDMultCuts::kFixed) { 
146       Double_t cut = 0.3;
147       AliWarningF("Using fixed cut @ %f for the lower bound "
148                   "of the sharing filter", cut);
149       sfLCuts.SetMultCuts(cut);
150     }
151     AliFMDMultCuts& sfHCuts = GetSharingFilter().GetHCuts();
152     if (sfHCuts.GetMethod() != AliFMDMultCuts::kFixed) { 
153       Double_t cut = 100;
154       AliWarningF("Using fixed cut @ %f for the upper bound "
155                   "of the sharing filter", cut);
156       sfHCuts.SetMultCuts(cut);
157     }
158     AliFMDMultCuts& dcCuts  = GetDensityCalculator().GetCuts();
159     if (dcCuts.GetMethod() != AliFMDMultCuts::kFixed) { 
160       Double_t cut = 0.3;
161       AliWarningF("Using fixed cut @ %f for the lower bound "
162                   "of the density calculator", cut);
163       dcCuts.SetMultCuts(cut);
164     }
165   }
166   else 
167     fcm.GetELossFit()->CacheBins(GetDensityCalculator().GetMinQuality());
168
169   fHistos.Init(eta);
170
171   // GetEventInspector().SetupForData(vertex);
172   GetEnergyFitter()     .SetupForData(eta);
173   GetSharingFilter()    .SetupForData(eta);
174   GetDensityCalculator().SetupForData(eta);
175
176   return true;
177 }
178
179 //____________________________________________________________________
180 Bool_t
181 AliForwardQATask::PreEvent()
182 {
183   // Clear stuff 
184   fHistos.Clear();
185   fESDFMD.Clear();
186   return true;
187 }
188 //____________________________________________________________________
189 Bool_t
190 AliForwardQATask::Event(AliESDEvent& esd)
191 {
192   // 
193   // Process each event 
194   // 
195   // Parameters:
196   //    option Not used
197   //  
198   DGUARD(fDebug,1,"Process the input event");
199
200   if (fFirstEvent) { 
201     // If the first event flag wasn't cleared in the above call to
202     // GetESDEvent, we should not do anything, since nothing has been
203     // initialised yet, so we opt out here (with a warning) 
204     AliWarning("Nothing has been initialized yet, opt'ing out");
205     return false;
206   }
207
208   Bool_t   lowFlux   = kFALSE;
209   UInt_t   triggers  = 0;
210   UShort_t ivz       = 0;
211   TVector3 ip;
212   Double_t cent      = -1;
213   UShort_t nClusters = 0;
214   UInt_t   found     = fEventInspector.Process(&esd, triggers, lowFlux, 
215                                                ivz, ip, cent, nClusters);
216   
217   Bool_t ok = true;
218   if (found & AliFMDEventInspector::kNoEvent)    ok = false;
219   if (found & AliFMDEventInspector::kNoTriggers) ok = false;
220   if (found & AliFMDEventInspector::kNoSPD)      ok = false;
221   if (found & AliFMDEventInspector::kNoFMD)      ok = false;
222   if (found & AliFMDEventInspector::kNoVertex)   ok = false;
223   if (triggers & AliAODForwardMult::kPileUp)     ok = false;
224   if (triggers & AliAODForwardMult::kA)          ok = false;
225   if (triggers & AliAODForwardMult::kC)          ok = false;
226   if (triggers & AliAODForwardMult::kE)          ok = false;
227   if (!(triggers & AliAODForwardMult::kOffline)) ok = false;
228   if (found & AliFMDEventInspector::kBadVertex)  ok = false;
229   if (!ok) { 
230     DMSG(fDebug,2,"Event failed selection: %s", 
231          AliFMDEventInspector::CodeString(found));
232     return false;
233   }
234   DMSG(fDebug,2,"Event triggers: %s", 
235        AliAODForwardMult::GetTriggerString(triggers));
236
237   // We we do not want to use low flux specific code, we disable it here. 
238   if (!fEnableLowFlux) lowFlux = false;
239
240   // Get FMD data 
241   AliESDFMD* esdFMD = esd.GetFMDData();
242   
243   // Run the energy loss fitter 
244   if (!fEnergyFitter.Accumulate(*esdFMD, cent, 
245                                 triggers & AliAODForwardMult::kEmpty)) {
246     AliWarning("Energy fitter failed");
247     return false;
248   }
249   
250   //  // Apply the sharing filter (or hit merging or clustering if you like)
251   if (!fSharingFilter.Filter(*esdFMD, lowFlux, fESDFMD, ip.Z())) { 
252     AliWarning("Sharing filter failed!");
253     return false;
254   }
255  
256   // Calculate the inclusive charged particle density 
257   if (!fDensityCalculator.Calculate(fESDFMD, fHistos, lowFlux, cent, ip)) { 
258     // if (!fDensityCalculator.Calculate(*esdFMD, fHistos, ivz, lowFlux)) { 
259     AliWarning("Density calculator failed!");
260     return false;
261   }
262   
263   return true;
264 }
265
266 //____________________________________________________________________
267 Bool_t
268 AliForwardQATask::Finalize()
269 {
270   // 
271   // End of job
272   // 
273   // Parameters:
274   //    option Not used 
275   //
276   if (fDebug) AliInfo("In Forwards terminate");
277   TStopwatch swt;
278   swt.Start();
279
280   // Get our histograms from the container 
281   TH1I* hEventsTr    = 0;
282   TH1I* hEventsTrVtx = 0;
283   TH1I* hEventsAcc   = 0;
284   TH1I* hTriggers    = 0;
285   if (!fEventInspector.FetchHistograms(fList, 
286                                        hEventsTr, 
287                                        hEventsTrVtx, 
288                                        hEventsAcc,
289                                        hTriggers)) { 
290     AliError(Form("Didn't get histograms from event selector "
291                   "(hEventsTr=%p,hEventsTrVtx=%p,hEventsAcc=%p)", 
292                   hEventsTr, hEventsTrVtx,hEventsAcc));
293     return false;
294   }
295
296   TStopwatch swf;
297   swf.Start();
298   fEnergyFitter.Fit(fResults);
299   swf.Stop();
300   AliInfoF("Fitting took %d real-time seconds, and %f CPU seconds", 
301            Int_t(swf.RealTime()), swf.CpuTime());
302
303   fSharingFilter.Terminate(fList,fResults,Int_t(hEventsTr->Integral()));
304   fDensityCalculator.Terminate(fList,fResults,Int_t(hEventsTrVtx->Integral()));
305
306   if (fDebug) AliInfoF("Posting post processing results to %s", 
307                        fResults->GetName());
308   swt.Stop();
309   AliInfoF("Finalize took %d real-time seconds, and %f CPU seconds", 
310            Int_t(swt.RealTime()), swt.CpuTime());
311
312   return true;
313 }
314
315 //____________________________________________________________________
316 void
317 AliForwardQATask::Print(Option_t* option) const
318 {
319   // 
320   // Print information 
321   // 
322   // Parameters:
323   //    option Not used
324   //
325   AliBaseESDTask::Print(option);
326   gROOT->IncreaseDirLevel();
327   GetEnergyFitter()     .Print(option);
328   GetSharingFilter()    .Print(option);
329   GetDensityCalculator().Print(option);
330   gROOT->DecreaseDirLevel();
331 }
332
333 //
334 // EOF
335 //