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