]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FORWARD/analysis2/AliForwardQATask.cxx
Added ignores
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliForwardQATask.cxx
CommitLineData
96624385 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>
e2213ed5 30#include <TStopwatch.h>
96624385 31
32//====================================================================
33AliForwardQATask::AliForwardQATask()
34 : AliAnalysisTaskSE(),
35 fEnableLowFlux(false),
36 fFirstEvent(true),
37 fCorrManager(0),
38 fESDFMD(),
39 fHistos(),
40 fEventInspector(),
41 fEnergyFitter(),
42 fSharingFilter(),
43 fDensityCalculator(),
44 fList(0)
45{
46 //
47 // Constructor
48 //
49}
50
51//____________________________________________________________________
52AliForwardQATask::AliForwardQATask(const char* name)
53 : AliAnalysisTaskSE(name),
54 fEnableLowFlux(false),
55 fFirstEvent(true),
56 fCorrManager(0),
57 fESDFMD(),
58 fHistos(),
59 fEventInspector("event"),
60 fEnergyFitter("energy"),
61 fSharingFilter("sharing"),
62 fDensityCalculator("density"),
63 fList(0)
64{
65 //
66 // Constructor
67 //
68 // Parameters:
69 // name Name of task
70 //
71 DefineOutput(1, TList::Class());
72 DefineOutput(2, TList::Class());
73 fCorrManager = &AliForwardCorrectionManager::Instance();
74 fEnergyFitter.SetNParticles(1); // Just find the 1st peak
75 fEnergyFitter.SetDoMakeObject(false);
76 fEnergyFitter.SetUseIncreasingBins(true);
77 fEnergyFitter.SetDoFits(kTRUE);
78 fEnergyFitter.SetLowCut(0.4);
79 fEnergyFitter.SetFitRangeBinWidth(4);
80 fEnergyFitter.SetMinEntries(1000);
81}
82
83//____________________________________________________________________
84AliForwardQATask::AliForwardQATask(const AliForwardQATask& o)
85 : AliAnalysisTaskSE(o),
86 fEnableLowFlux(o.fEnableLowFlux),
87 fFirstEvent(o.fFirstEvent),
88 fCorrManager(o.fCorrManager),
89 fESDFMD(o.fESDFMD),
90 fHistos(o.fHistos),
91 fEventInspector(o.fEventInspector),
92 fEnergyFitter(o.fEnergyFitter),
93 fSharingFilter(o.fSharingFilter),
94 fDensityCalculator(o.fDensityCalculator),
95 fList(o.fList)
96{
97 //
98 // Copy constructor
99 //
100 // Parameters:
101 // o Object to copy from
102 //
103 DefineOutput(1, TList::Class());
104 DefineOutput(2, TList::Class());
105}
106
107//____________________________________________________________________
108AliForwardQATask&
109AliForwardQATask::operator=(const AliForwardQATask& o)
110{
111 //
112 // Assignment operator
113 //
114 // Parameters:
115 // o Object to assign from
116 //
117 // Return:
118 // Reference to this object
119 //
120 AliAnalysisTaskSE::operator=(o);
121
122 fEnableLowFlux = o.fEnableLowFlux;
123 fFirstEvent = o.fFirstEvent;
124 fCorrManager = o.fCorrManager;
125 fEventInspector = o.fEventInspector;
126 fEnergyFitter = o.fEnergyFitter;
127 fSharingFilter = o.fSharingFilter;
128 fDensityCalculator = o.fDensityCalculator;
129 fHistos = o.fHistos;
130 fList = o.fList;
131
132 return *this;
133}
134
135//____________________________________________________________________
136void
137AliForwardQATask::SetDebug(Int_t dbg)
138{
139 //
140 // Set debug level
141 //
142 // Parameters:
143 // dbg Debug level
144 //
145 fEventInspector.SetDebug(dbg);
146 fEnergyFitter.SetDebug(dbg);
147 fSharingFilter.SetDebug(dbg);
148 fDensityCalculator.SetDebug(dbg);
149}
150
151//____________________________________________________________________
152Bool_t
153AliForwardQATask::CheckCorrections(UInt_t what) const
154{
155 //
156 // Check if all needed corrections are there and accounted for. If not,
157 // do a Fatal exit
158 //
159 // Parameters:
160 // what Which corrections is needed
161 //
162 // Return:
163 // true if all present, false otherwise
164 //
165
166 AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
167 // Check that we have the energy loss fits, needed by
168 // AliFMDSharingFilter
169 // AliFMDDensityCalculator
170 if (what & AliForwardCorrectionManager::kELossFits && !fcm.GetELossFit()) {
171 AliFatal(Form("No energy loss fits"));
172 return false;
173 }
174 return true;
175}
176
177//____________________________________________________________________
178Bool_t
179AliForwardQATask::ReadCorrections(const TAxis*& pe,
180 const TAxis*& pv,
181 Bool_t mc)
182{
183 //
184 // Read corrections
185 //
186 //
187 UInt_t what = AliForwardCorrectionManager::kAll;
188 what ^= AliForwardCorrectionManager::kDoubleHit;
189 what ^= AliForwardCorrectionManager::kVertexBias;
190 what ^= AliForwardCorrectionManager::kAcceptance;
191 what ^= AliForwardCorrectionManager::kMergingEfficiency;
192
193 AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
194 if (!fcm.Init(GetEventInspector().GetCollisionSystem(),
195 GetEventInspector().GetEnergy(),
196 GetEventInspector().GetField(),
197 mc,
198 what)) return false;
199 if (!CheckCorrections(what)) return false;
200
201 // Sett our persistency pointer
202 // fCorrManager = &fcm;
203
204 // Get the eta axis from the secondary maps - if read in
205 if (!pe) {
206 pe = fcm.GetEtaAxis();
207 if (!pe) AliFatal("No eta axis defined");
208 }
209 // Get the vertex axis from the secondary maps - if read in
210 if (!pv) {
211 pv = fcm.GetVertexAxis();
212 if (!pv) AliFatal("No vertex axis defined");
213 }
214
215 return true;
216}
217
218//____________________________________________________________________
219AliESDEvent*
220AliForwardQATask::GetESDEvent()
221{
222 //
223 // Get the ESD event. IF this is the first event, initialise
224 //
225 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
226 if (!esd) {
227 AliWarning("No ESD event found for input event");
228 return 0;
229 }
230
231 // On the first event, initialize the parameters
232 if (fFirstEvent && esd->GetESDRun()) {
233 GetEventInspector().ReadRunDetails(esd);
234
422a78c8 235 AliInfoF("Initializing with parameters from the ESD:\n"
236 " AliESDEvent::GetBeamEnergy() ->%f\n"
237 " AliESDEvent::GetBeamType() ->%s\n"
238 " AliESDEvent::GetCurrentL3() ->%f\n"
239 " AliESDEvent::GetMagneticField()->%f\n"
240 " AliESDEvent::GetRunNumber() ->%d\n",
241 esd->GetBeamEnergy(),
242 esd->GetBeamType(),
243 esd->GetCurrentL3(),
244 esd->GetMagneticField(),
245 esd->GetRunNumber());
96624385 246
247 fFirstEvent = false;
248
422a78c8 249 if (!InitializeSubs()) {
250 AliWarning("Initialisation of sub algorithms failed!");
251 return 0;
252 }
96624385 253 }
254 return esd;
255}
256//____________________________________________________________________
422a78c8 257Bool_t
96624385 258AliForwardQATask::InitializeSubs()
259{
260 //
261 // Initialise the sub objects and stuff. Called on first event
262 //
263 //
264 const TAxis* pe = 0;
265 const TAxis* pv = 0;
266
422a78c8 267 if (!ReadCorrections(pe,pv)) return false;
96624385 268
269 fHistos.Init(*pe);
270
271 fEventInspector.Init(*pv);
272 fEnergyFitter.Init(*pe);
273 fSharingFilter.Init();
274 fDensityCalculator.Init(*pe);
275
276 this->Print();
422a78c8 277
278 return true;
96624385 279}
280
281//____________________________________________________________________
282void
283AliForwardQATask::UserCreateOutputObjects()
284{
285 //
286 // Create output objects
287 //
288 //
289 fList = new TList;
290 fList->SetOwner();
291
292 fEventInspector.DefineOutput(fList);
293 fEnergyFitter.DefineOutput(fList);
294 fSharingFilter.DefineOutput(fList);
295 fDensityCalculator.DefineOutput(fList);
296
297 PostData(1, fList);
298}
299//____________________________________________________________________
300void
301AliForwardQATask::UserExec(Option_t*)
302{
303 //
304 // Process each event
305 //
306 // Parameters:
307 // option Not used
308 //
309
310 // static Int_t cnt = 0;
311 // cnt++;
312 // Get the input data
313 AliESDEvent* esd = GetESDEvent();
422a78c8 314 if (!esd) {
315 AliWarning("Got no ESD event");
316 return;
317 }
96624385 318
319 // Clear stuff
320 fHistos.Clear();
321 fESDFMD.Clear();
322
323 Bool_t lowFlux = kFALSE;
324 UInt_t triggers = 0;
325 UShort_t ivz = 0;
326 Double_t vz = 0;
327 Double_t cent = -1;
328 UShort_t nClusters = 0;
329 UInt_t found = fEventInspector.Process(esd, triggers, lowFlux,
330 ivz, vz, cent, nClusters);
331
332 if (found & AliFMDEventInspector::kNoEvent) return;
333 if (found & AliFMDEventInspector::kNoTriggers) return;
334 if (found & AliFMDEventInspector::kNoSPD) return;
335 if (found & AliFMDEventInspector::kNoFMD) return;
336 if (found & AliFMDEventInspector::kNoVertex) return;
337 if (triggers & AliAODForwardMult::kPileUp) return;
338 if (found & AliFMDEventInspector::kBadVertex) return;
339
340 // We we do not want to use low flux specific code, we disable it here.
341 if (!fEnableLowFlux) lowFlux = false;
342
343 // Get FMD data
344 AliESDFMD* esdFMD = esd->GetFMDData();
345
346 // Run the energy loss fitter
347 if (!fEnergyFitter.Accumulate(*esdFMD, cent,
348 triggers & AliAODForwardMult::kEmpty)) {
349 AliWarning("Energy fitter failed");
350 return;
351 }
352
353 // // Apply the sharing filter (or hit merging or clustering if you like)
354 if (!fSharingFilter.Filter(*esdFMD, lowFlux, fESDFMD)) {
355 AliWarning("Sharing filter failed!");
356 return;
357 }
358
359 // Calculate the inclusive charged particle density
360 if (!fDensityCalculator.Calculate(fESDFMD, fHistos, ivz, lowFlux)) {
361 // if (!fDensityCalculator.Calculate(*esdFMD, fHistos, ivz, lowFlux)) {
362 AliWarning("Density calculator failed!");
363 return;
364 }
365 PostData(1, fList);
366}
367
368//____________________________________________________________________
369void
370AliForwardQATask::Terminate(Option_t*)
371{
372 //
373 // End of job
374 //
375 // Parameters:
376 // option Not used
377 //
e2213ed5 378 TStopwatch swt;
379 swt.Start();
96624385 380
381 TList* list = dynamic_cast<TList*>(GetOutputData(1));
382 if (!list) {
383 AliError(Form("No output list defined (%p)", GetOutputData(1)));
384 if (GetOutputData(1)) GetOutputData(1)->Print();
385 return;
386 }
387
388 // Get our histograms from the container
389 TH1I* hEventsTr = 0;//static_cast<TH1I*>(list->FindObject("nEventsTr"));
390 TH1I* hEventsTrVtx = 0;//static_cast<TH1I*>(list->FindObject("nEventsTrVtx"));
391 TH1I* hTriggers = 0;
392 if (!fEventInspector.FetchHistograms(list, hEventsTr,
393 hEventsTrVtx, hTriggers)) {
394 AliError(Form("Didn't get histograms from event selector "
395 "(hEventsTr=%p,hEventsTrVtx=%p)",
396 hEventsTr, hEventsTrVtx));
397 return;
398 }
399
e2213ed5 400 TStopwatch swf;
401 swf.Start();
96624385 402 fEnergyFitter.Fit(list);
e2213ed5 403 swf.Stop();
404 AliInfoF("Fitting took %d real-time seconds, and %f CPU seconds",
405 Int_t(swf.RealTime()), swf.CpuTime());
96624385 406
407 fSharingFilter.ScaleHistograms(list,Int_t(hEventsTr->Integral()));
408 fDensityCalculator.ScaleHistograms(list,Int_t(hEventsTrVtx->Integral()));
409
410 // Make a deep copy and post that as output 2
411 TList* list2 = static_cast<TList*>(list->Clone(Form("%sResults",
412 list->GetName())));
413 list2->SetOwner();
414 PostData(2, list2);
415
e2213ed5 416 swt.Stop();
417 AliInfoF("Terminate took %d real-time seconds, and %f CPU seconds",
418 Int_t(swt.RealTime()), swt.CpuTime());
419
96624385 420}
421
422//____________________________________________________________________
423void
424AliForwardQATask::Print(Option_t* option) const
425{
426 //
427 // Print information
428 //
429 // Parameters:
430 // option Not used
431 //
432
433 std::cout << ClassName() << ": " << GetName() << "\n"
434 << " Enable low flux code: " << (fEnableLowFlux ? "yes" : "no")
435 << "\n"
436 << " Off-line trigger mask: 0x"
437 << std::hex << std::setfill('0')
438 << std::setw (8) << fOfflineTriggerMask
439 << std::dec << std::setfill (' ') << std::endl;
440 gROOT->IncreaseDirLevel();
441 if (fCorrManager) fCorrManager->Print();
442 else
443 std::cout << " Correction manager not set yet" << std::endl;
444 GetEventInspector() .Print(option);
445 GetEnergyFitter() .Print(option);
446 GetSharingFilter() .Print(option);
447 gROOT->DecreaseDirLevel();
448}
449
450//
451// EOF
452//