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