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