]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliForwardMultiplicityBase.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliForwardMultiplicityBase.cxx
1 //====================================================================
2 // 
3 // Base class for classes that calculate the multiplicity in the
4 // forward regions event-by-event
5 // 
6 // Inputs: 
7 //   - AliESDEvent 
8 //
9 // Outputs: 
10 //   - AliAODForwardMult 
11 // 
12 // Histograms 
13 //   
14 // Corrections used 
15 #include "AliForwardMultiplicityBase.h"
16 #include "AliForwardCorrectionManager.h"
17 #include "AliForwardUtil.h"
18 #include "AliFMDCorrELossFit.h"
19 #include "AliLog.h"
20 #include "AliAODHandler.h"
21 #include "AliInputEventHandler.h"
22 #include "AliAnalysisManager.h"
23 #include "AliFMDEventInspector.h"
24 #include "AliFMDESDFixer.h"
25 #include "AliFMDSharingFilter.h"
26 #include "AliFMDDensityCalculator.h"
27 #include "AliFMDCorrector.h"
28 #include "AliFMDHistCollector.h"
29 #include "AliFMDEventPlaneFinder.h"
30 #include "AliESDEvent.h"
31 #include <TROOT.h>
32 #include <TSystem.h>
33 #include <TAxis.h>
34 #include <TProfile.h>
35 #include <THStack.h>
36 #include <iostream>
37 #include <iomanip>
38 #define AOD_SLOT 3
39 #ifdef POST_AOD
40 # define DEFINE(N) DefineOutput(N,AliAODForwardMult::Class())
41 # define POST(N)   PostData(N,fAODFMD)
42 #else
43 # define DEFINE(N) do { } while(false)
44 # define POST(N)   do { } while(false)
45 #endif
46
47 //====================================================================
48 AliForwardMultiplicityBase::AliForwardMultiplicityBase(const char* name) 
49   : AliBaseESDTask(name, "AliForwardMultiplicityBase", 
50                    &(AliForwardCorrectionManager::Instance())),
51     fEnableLowFlux(false), 
52     fStorePerRing(false),
53     fHData(0),
54     fHistos(),
55     fAODFMD(false),
56     fAODEP(false),
57     fRingSums(),
58     fDoTiming(false),
59   fHTiming(0),
60   fHStatus(0)
61 {
62   DGUARD(fDebug, 3,"Named CTOR of AliForwardMultiplicityBase %s",name);
63
64   DEFINE(AOD_SLOT);
65 }
66
67
68 //____________________________________________________________________
69 void
70 AliForwardMultiplicityBase::SetDebug(Int_t dbg)
71 {
72   // 
73   // Set debug level 
74   // 
75   // Parameters:
76   //    dbg debug level
77   //
78   AliBaseESDTask::       SetDebug(dbg);
79   GetSharingFilter()    .SetDebug(dbg);
80   GetDensityCalculator().SetDebug(dbg);
81   GetCorrections()      .SetDebug(dbg);
82   GetHistCollector()    .SetDebug(dbg);
83   GetEventPlaneFinder() .SetDebug(dbg);
84 }
85
86 //____________________________________________________________________
87 Bool_t
88 AliForwardMultiplicityBase::Book()
89 {
90   // 
91   // Create output objects 
92   // 
93   //
94   DGUARD(fDebug,1,"Create user ouput");
95   UInt_t what = AliForwardCorrectionManager::kAll;
96   if (!fEnableLowFlux)
97     what ^= AliForwardCorrectionManager::kDoubleHit;
98   if (!GetESDFixer().IsUseNoiseCorrection()) 
99     what ^= AliForwardCorrectionManager::kNoiseGain;
100   if (!GetCorrections().IsUseVertexBias())
101     what ^= AliForwardCorrectionManager::kVertexBias;
102   if (!GetCorrections().IsUseAcceptance())
103     what ^= AliForwardCorrectionManager::kAcceptance;
104   if (!GetCorrections().IsUseMergingEfficiency())
105     what ^= AliForwardCorrectionManager::kMergingEfficiency;
106   fNeededCorrections = what;
107
108   GetESDFixer()         .CreateOutputObjects(fList);
109   GetSharingFilter()    .CreateOutputObjects(fList);
110   GetDensityCalculator().CreateOutputObjects(fList);
111   GetCorrections()      .CreateOutputObjects(fList);
112   GetHistCollector()    .CreateOutputObjects(fList);
113   GetEventPlaneFinder() .CreateOutputObjects(fList);
114
115   TAxis tmp(1, 0, 1);
116   fHistos.Init(tmp);
117
118   if (fDebug > 1) fDoTiming = true;
119   if (fDoTiming) { 
120     fHTiming = new TProfile("timing", "Timing of task", 
121                             kTimingTotal, 0.5, kTimingTotal+.5);
122     fHTiming->SetDirectory(0);
123     fHTiming->SetFillColor(kRed+1);
124     fHTiming->SetFillStyle(3001);
125     fHTiming->SetMarkerStyle(20);
126     fHTiming->SetMarkerColor(kBlack);
127     fHTiming->SetLineColor(kBlack);
128     fHTiming->SetXTitle("Part");
129     fHTiming->SetYTitle("#LTt_{part}#GT [CPU]");
130     fHTiming->SetStats(0);
131     TAxis* xaxis = fHTiming->GetXaxis();
132     xaxis->SetBinLabel(kTimingEventInspector,   
133                        GetEventInspector()   .GetName());
134     xaxis->SetBinLabel(kTimingSharingFilter,    
135                        GetSharingFilter()    .GetName());
136     xaxis->SetBinLabel(kTimingDensityCalculator,        
137                        GetDensityCalculator().GetName());
138     xaxis->SetBinLabel(kTimingCorrections,      
139                        GetCorrections()      .GetName());
140     xaxis->SetBinLabel(kTimingHistCollector,    
141                        GetHistCollector()    .GetName());
142     xaxis->SetBinLabel(kTimingEventPlaneFinder, 
143                        GetEventPlaneFinder() .GetName());
144     xaxis->SetBinLabel(kTimingTotal, "Total");
145     fList->Add(fHTiming);
146   }
147   fHStatus = new TH1I("status", "Status of task",15, .5, 15.5);
148   fHStatus->SetDirectory(0);
149   fHStatus->SetFillColor(kCyan+2);
150   fHStatus->SetFillStyle(3001);
151   fHStatus->SetLineColor(kBlack);
152   fHStatus->SetMarkerStyle(20);
153   fHStatus->SetMarkerColor(kBlack);
154   fHStatus->SetYTitle("Events");
155   TAxis* a = fHStatus->GetXaxis();
156   a->SetBinLabel(1,"No event");
157   a->SetBinLabel(2,"No triggers");
158   a->SetBinLabel(3,"No SPD (not used)");
159   a->SetBinLabel(4,"No FMD");
160   a->SetBinLabel(5,"No Vertex");
161   a->SetBinLabel(6,"Pile-up");
162   a->SetBinLabel(7,"IP_{z} out of range");
163   a->SetBinLabel(8,"Merger failed");
164   a->SetBinLabel(9,"N_{ch} estimator failed");
165   a->SetBinLabel(10,"#Phi_{R} estimator failed");
166   a->SetBinLabel(11,"Too many outliers");
167   a->SetBinLabel(12,"Corrector failed");
168   a->SetBinLabel(13,"Histogram collector failed");
169   a->SetBinLabel(14,"Not added");
170   a->SetBinLabel(15,"All through");
171   fList->Add(fHStatus);
172
173   POST(AOD_SLOT);
174   return true;
175 }
176 //____________________________________________________________________
177 void
178 AliForwardMultiplicityBase::CreateBranches(AliAODHandler* ah)
179 {
180   TObject* obj   = &fAODFMD; ah->AddBranch("AliAODForwardMult", &obj);
181   TObject* epobj = &fAODEP;  ah->AddBranch("AliAODForwardEP", &epobj);
182
183   if (!fStorePerRing) return;
184   
185   AliWarning("Per-ring histograms in AOD\n"
186              "*********************************************************\n"
187              "* For each event 5 additional 2D histogram are stored   *\n"
188              "* in separate branches of the AODs.  This will increase *\n"
189              "* the size of the AODs - proceed with caution           *\n"
190              "*********************************************************");
191   TObject* hists[] = { fHistos.fFMD1i, 
192                        fHistos.fFMD2i, fHistos.fFMD2o, 
193                        fHistos.fFMD3i, fHistos.fFMD3o };
194   for (Int_t i = 0; i < 5; i++) { 
195     ah->AddBranch("TH2D", &(hists[i]));
196   }
197 }
198
199 //____________________________________________________________________
200 Bool_t
201 AliForwardMultiplicityBase::PreData(const TAxis& vertex, const TAxis& eta)
202 {
203   // 
204   // Initialise the sub objects and stuff.  Called on first event 
205   // 
206   //
207   DGUARD(fDebug,1,"Initialize sub-algorithms");
208
209   // Force this here so we select the proper quality 
210   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
211
212   UInt_t what = fNeededCorrections;
213   // Check that we have the energy loss fits, needed by 
214   //   AliFMDSharingFilter 
215   //   AliFMDDensityCalculator 
216   if (what & AliForwardCorrectionManager::kELossFits && !fcm.GetELossFit()) 
217     AliFatal(Form("No energy loss fits"));
218   
219   // Check that we have the double hit correction - (optionally) used by 
220   //  AliFMDDensityCalculator 
221   if (what & AliForwardCorrectionManager::kDoubleHit && !fcm.GetDoubleHit()) 
222     AliFatal("No double hit corrections"); 
223
224   // Check that we have the secondary maps, needed by 
225   //   AliFMDCorrector 
226   //   AliFMDHistCollector
227   if (what & AliForwardCorrectionManager::kSecondaryMap && 
228       !fcm.GetSecondaryMap()) 
229     AliFatal("No secondary corrections");
230
231   // Check that we have the vertex bias correction, needed by 
232   //   AliFMDCorrector 
233   if (what & AliForwardCorrectionManager::kVertexBias && !fcm.GetVertexBias())
234     AliFatal("No event vertex bias corrections");
235
236   // Check that we have the merging efficiencies, optionally used by 
237   //   AliFMDCorrector 
238   if (what & AliForwardCorrectionManager::kMergingEfficiency && 
239       !fcm.GetMergingEfficiency()) 
240     AliFatal("No merging efficiencies");
241
242   // Check that we have the acceptance correction, needed by 
243   //   AliFMDCorrector 
244   if (what & AliForwardCorrectionManager::kAcceptance && !fcm.GetAcceptance())
245     AliFatal("No acceptance corrections");
246
247   // const AliFMDCorrELossFit* fits = fcm.GetELossFit();
248   // fits->CacheBins(GetDensityCalculator().GetMinQuality());
249
250   InitMembers(eta,vertex);
251   
252   GetDensityCalculator().SetupForData(eta);
253   GetSharingFilter()    .SetupForData(eta);
254   GetCorrections()      .SetupForData(eta);
255   GetHistCollector()    .SetupForData(vertex,eta);
256
257   GetEventPlaneFinder() .SetRunNumber(GetEventInspector().GetRunNumber());
258   GetEventPlaneFinder() .SetupForData(eta);
259   
260   fAODFMD.SetBit(AliAODForwardMult::kSecondary, 
261                  GetCorrections().IsUseSecondaryMap());
262   fAODFMD.SetBit(AliAODForwardMult::kVertexBias, 
263                  GetCorrections().IsUseVertexBias());
264   fAODFMD.SetBit(AliAODForwardMult::kAcceptance, 
265                  GetCorrections().IsUseAcceptance());
266   fAODFMD.SetBit(AliAODForwardMult::kMergingEfficiency, 
267                  GetCorrections().IsUseMergingEfficiency());
268   fAODFMD.SetBit(AliAODForwardMult::kSum, 
269                  GetHistCollector().GetMergeMethod() == 
270                  AliFMDHistCollector::kSum);
271   return true;
272 }
273
274 //____________________________________________________________________
275 void
276 AliForwardMultiplicityBase::InitMembers(const TAxis& eta, const TAxis& /*pv*/)
277 {
278   fHistos.ReInit(eta);
279   fAODFMD.Init(eta);
280   fAODEP.Init(eta);
281   fRingSums.Init(eta);
282
283   fHData = static_cast<TH2D*>(fAODFMD.GetHistogram().Clone("d2Ndetadphi"));
284   fHData->SetStats(0);
285   fHData->SetDirectory(0);
286   fList->Add(fHData);
287
288   TList* rings = new TList;
289   rings->SetName("ringSums");
290   rings->SetOwner();
291   fList->Add(rings);
292
293   rings->Add(fRingSums.Get(1, 'I'));
294   rings->Add(fRingSums.Get(2, 'I'));
295   rings->Add(fRingSums.Get(2, 'O'));
296   rings->Add(fRingSums.Get(3, 'I'));
297   rings->Add(fRingSums.Get(3, 'O'));
298   fRingSums.Get(1, 'I')->SetMarkerColor(AliForwardUtil::RingColor(1, 'I'));
299   fRingSums.Get(2, 'I')->SetMarkerColor(AliForwardUtil::RingColor(2, 'I'));
300   fRingSums.Get(2, 'O')->SetMarkerColor(AliForwardUtil::RingColor(2, 'O'));
301   fRingSums.Get(3, 'I')->SetMarkerColor(AliForwardUtil::RingColor(3, 'I'));
302   fRingSums.Get(3, 'O')->SetMarkerColor(AliForwardUtil::RingColor(3, 'O'));
303 }
304 //____________________________________________________________________
305 Bool_t
306 AliForwardMultiplicityBase::PostEvent()
307 {
308   POST(AOD_SLOT);
309   return true;
310 }
311 //____________________________________________________________________
312 Bool_t
313 AliForwardMultiplicityBase::Finalize()
314 {
315   // 
316   // End of job
317   // 
318   // Parameters:
319   //    option Not used 
320   //
321   DGUARD(fDebug,1,"Processing the merged results");
322
323   TList* list   = fList;
324   TList* output = fResults;
325   Double_t nTr = 0, nTrVtx = 0, nAcc = 0;
326   MakeSimpledNdeta(list, output, nTr, nTrVtx, nAcc);
327   AliInfoF("\n"
328            "\t# events w/trigger:                 %f\n"
329            "\t# events w/trigger+vertex:          %f\n"
330            "\t# events accepted by cuts:          %f", 
331            nTr, nTrVtx, nAcc);
332
333   EstimatedNdeta(list, output);
334
335   GetSharingFilter()    .Terminate(list,output,Int_t(nTr));
336   GetDensityCalculator().Terminate(list,output,Int_t(nTrVtx));
337   GetCorrections()      .Terminate(list,output,Int_t(nTrVtx));
338
339   TProfile* timing = static_cast<TProfile*>(list->FindObject("timing"));
340   if (timing) { 
341     TProfile* p = static_cast<TProfile*>(timing->Clone());
342     p->SetDirectory(0);
343     p->Scale(100. / p->GetBinContent(p->GetNbinsX()));
344     p->SetYTitle("#LTt_{part}#GT/#LTt_{total}#GT [%]");
345     p->SetTitle("Relative timing of task");
346     output->Add(p);
347   }
348   return true;
349 }
350
351 //____________________________________________________________________
352 void
353 AliForwardMultiplicityBase::EstimatedNdeta(const TList* input, 
354                                            TList*       output) const
355 {
356   MakeRingdNdeta(input, "ringSums", output, "ringResults");
357 }
358
359 //____________________________________________________________________
360 Bool_t
361 AliForwardMultiplicityBase::MakeSimpledNdeta(const TList* input, 
362                                              TList*       output,
363                                              Double_t&    nTr, 
364                                              Double_t&    nTrVtx, 
365                                              Double_t&    nAcc)
366 {
367   // Get our histograms from the container 
368   TH1I* hEventsTr    = 0;
369   TH1I* hEventsTrVtx = 0;
370   TH1I* hEventsAcc   = 0;
371   TH1I* hTriggers    = 0;
372   TH1*  hStatus      = dynamic_cast<TH1*>(input->FindObject("status"));
373   if (!GetEventInspector().FetchHistograms(input, 
374                                            hEventsTr, 
375                                            hEventsTrVtx, 
376                                            hEventsAcc,
377                                            hTriggers)) { 
378     AliError(Form("Didn't get histograms from event selector "
379                   "(hEventsTr=%p,hEventsTrVtx=%p,hEventsAcc=%p,hTriggers=%p)", 
380                   hEventsTr, hEventsTrVtx, hEventsAcc, hTriggers));
381     input->ls();
382     return false;
383   }
384   nTr              = hEventsTr->Integral();
385   nTrVtx           = hEventsTrVtx->Integral();
386   nAcc             = hEventsAcc->Integral();
387   Double_t vtxEff  = nTrVtx / nTr;
388   Double_t vtxEff2 = 0;
389   if (hStatus) {
390     Double_t nTrg    = hStatus->Integral(3,15);
391     Double_t nTrgVtx = hStatus->Integral(6,15);
392     vtxEff2          = (nTrg > 0 ? nTrgVtx / nTrg : 0);
393   }
394     
395   TH2D*   hData   = static_cast<TH2D*>(input->FindObject("d2Ndetadphi"));
396   if (!hData) { 
397     AliError(Form("Couldn't get our summed histogram from output "
398                   "list %s (d2Ndetadphi=%p)", input->GetName(), hData));
399     input->ls();
400     return false;
401   }
402
403   Int_t nY      = hData->GetNbinsY();
404   TH1D* dNdeta  = hData->ProjectionX("dNdeta",  1,     nY, "e");
405   TH1D* dNdeta_ = hData->ProjectionX("dNdeta_", 1,     nY, "e");
406   TH1D* norm    = hData->ProjectionX("norm",    0,     0,  "");
407   TH1D* phi     = hData->ProjectionX("phi",     nY+1,  nY+1,  "");
408   dNdeta->SetTitle("dN_{ch}/d#eta in the forward regions");
409   dNdeta->SetYTitle("#frac{1}{N}#frac{dN_{ch}}{d#eta}");
410   dNdeta->SetMarkerColor(kRed+1);
411   dNdeta->SetMarkerStyle(20);
412   dNdeta->SetDirectory(0);
413
414   dNdeta_->SetTitle("dN_{ch}/d#eta in the forward regions");
415   dNdeta_->SetYTitle("#frac{1}{N}#frac{dN_{ch}}{d#eta}");
416   dNdeta_->SetMarkerColor(kMagenta+1);
417   dNdeta_->SetMarkerStyle(21);
418   dNdeta_->SetDirectory(0);
419
420   norm->SetTitle("Normalization to  #eta coverage");
421   norm->SetYTitle("#eta coverage");
422   norm->SetLineColor(kBlue+1);
423   norm->SetMarkerColor(kBlue+1);
424   norm->SetMarkerStyle(21);
425   norm->SetFillColor(kBlue+1);
426   norm->SetFillStyle(3005);
427   norm->SetDirectory(0);
428
429   phi->SetTitle("Normalization to  #phi acceptance");
430   phi->SetYTitle("#phi acceptance");
431   phi->SetLineColor(kGreen+1);
432   phi->SetMarkerColor(kGreen+1);
433   phi->SetMarkerStyle(20);
434   phi->SetFillColor(kGreen+1);
435   phi->SetFillStyle(3004);
436   // phi->Scale(1. / nAcc);
437   phi->SetDirectory(0);
438
439   // dNdeta->Divide(norm);
440   dNdeta->Divide(phi);
441   dNdeta->SetStats(0);
442   dNdeta->Scale(vtxEff, "width");
443
444   dNdeta_->Divide(norm);
445   dNdeta_->SetStats(0);
446   dNdeta_->Scale(vtxEff, "width");
447
448   AliInfoF("All\n"
449            "\tNormalization eta:  %d\n"
450            "\tNormalization phi:  %d\n"
451            "\tVertex efficiency:  %f (%f)",
452            Int_t(norm->GetMaximum()), Int_t(phi->GetMaximum()), 
453            vtxEff, vtxEff2);
454   output->Add(dNdeta);
455   output->Add(dNdeta_);
456   output->Add(norm);
457   output->Add(phi);
458
459   return true;
460 }
461
462                                              
463 //____________________________________________________________________
464 void
465 AliForwardMultiplicityBase::MakeRingdNdeta(const TList* input, 
466                                            const char*  inName,
467                                            TList*       output,
468                                            const char*  outName,
469                                            Int_t        style) const
470 {
471   // Make dN/deta for each ring found in the input list.  
472   // 
473   // A stack of all the dN/deta is also made for easy drawing. 
474   // 
475   // Note, that the distributions are normalised to the number of
476   // observed events only - they should be corrected for 
477   DGUARD(fDebug,3,"Make first-shot ring dN/deta");
478
479   if (!input) return;
480   TList* list = static_cast<TList*>(input->FindObject(inName));
481   if (!list) { 
482     AliWarning(Form("No list %s found in %s", inName, input->GetName()));
483     return;
484   }
485   
486   TList* out = new TList;
487   out->SetName(outName);
488   out->SetOwner();
489   output->Add(out);
490
491   THStack*     dndetaRings = new THStack("all", "dN/d#eta per ring");
492   const char*  names[]     = { "FMD1I", "FMD2I", "FMD2O", "FMD3I", "FMD3O", 0 };
493   const char** ptr         = names;
494   
495   while (*ptr) { 
496     TList* thisList = new TList;
497     thisList->SetOwner();
498     thisList->SetName(*ptr);
499     out->Add(thisList);
500
501     TH2D* h = static_cast<TH2D*>(list->FindObject(Form("%s_cache", *ptr)));
502     if (!h) { 
503       AliWarning(Form("Didn't find %s_cache in %s", *ptr, list->GetName()));
504       ptr++;
505       continue;
506     }
507     TH2D* sumPhi = static_cast<TH2D*>(h->Clone("sum_phi"));
508     sumPhi->SetDirectory(0);
509     thisList->Add(sumPhi);
510
511     TH2D* sumEta = static_cast<TH2D*>(h->Clone("sum_eta"));
512     sumEta->SetDirectory(0);
513     thisList->Add(sumEta);
514     
515     Int_t nY   = sumEta->GetNbinsY();
516     TH1D* etaCov =static_cast<TH1D*>(h->ProjectionX("etaCov", 0,    0,    ""));
517     TH1D* phiAcc =static_cast<TH1D*>(h->ProjectionX("phiAcc", nY+1, nY+1, ""));
518
519     etaCov->SetTitle("Normalization to #eta coverage");
520     etaCov->SetYTitle("#eta coverage");
521     etaCov->SetMarkerColor(kBlue+1);
522     etaCov->SetFillColor(kBlue+1);
523     etaCov->SetFillStyle(3005);
524     etaCov->SetDirectory(0);
525     
526     phiAcc->SetTitle("Normalization to #phi acceptance");
527     phiAcc->SetYTitle("#phi acceptance");
528     phiAcc->SetMarkerColor(kGreen+1);
529     phiAcc->SetFillColor(kGreen+1);
530     phiAcc->SetFillStyle(3004);
531     // phiAcc->Scale(1. / nAcc);
532     phiAcc->SetDirectory(0);
533
534     // Double_t s = (etaCov->GetMaximum() > 0 ? 1. / etaCov->GetMaximum() : 1);
535     for (Int_t i = 1; i <= sumEta->GetNbinsX(); i++) { 
536       for (Int_t j = 1; j <= nY; j++) { 
537         Double_t c = sumEta->GetBinContent(i, j);
538         Double_t e = sumEta->GetBinError(i, j);
539         Double_t a = etaCov->GetBinContent(i);
540         Double_t p = phiAcc->GetBinContent(i);
541         // Double_t t = p; // * a
542         sumEta->SetBinContent(i, j, a <= 0 ? 0 : c / a);
543         sumEta->SetBinError(  i, j, a <= 0 ? 0 : e / a);
544         sumPhi->SetBinContent(i, j, p <= 0 ? 0 : c / p);
545         sumPhi->SetBinError(  i, j, p <= 0 ? 0 : e / p);
546       }
547     }
548     // etaCov->Scale(s);
549     // phiAcc->Scale(s);
550
551     TH1D* resPhi  =static_cast<TH1D*>(sumPhi->ProjectionX("dndeta_phi",
552                                                           1,nY,"e"));
553     resPhi->SetMarkerStyle(style);
554     resPhi->SetDirectory(0);
555     resPhi->Scale(1, "width");
556
557     TH1D* resEta  =static_cast<TH1D*>(sumEta->ProjectionX("dndeta_eta",
558                                                           1,nY,"e"));
559     resEta->SetMarkerStyle(style);
560     resEta->SetDirectory(0);
561     resEta->Scale(1, "width");
562
563     thisList->Add(resEta);
564     thisList->Add(etaCov);
565     thisList->Add(resPhi);
566     thisList->Add(phiAcc);
567     dndetaRings->Add(resPhi);
568
569     AliInfoF("%s\n"
570              "\tNormalization eta:  %d\n"
571              "\tNormalization phi:  %d\n",
572              *ptr, 
573              Int_t(etaCov->GetMaximum()), Int_t(phiAcc->GetMaximum()));
574
575     ptr++;
576   }
577   out->Add(dndetaRings);
578 }
579 #define PFB(N,FLAG)                             \
580   do {                                                                  \
581     AliForwardUtil::PrintName(N);                                       \
582     std::cout << std::boolalpha << (FLAG) << std::noboolalpha << std::endl; \
583   } while(false)
584 //____________________________________________________________________
585 void
586 AliForwardMultiplicityBase::Print(Option_t* option) const
587 {
588   // 
589   // Print information 
590   // 
591   // Parameters:
592   //    option Not used
593   //
594   AliBaseESDTask::Print(option);
595   gROOT->IncreaseDirLevel();
596   PFB("Enable low flux code", fEnableLowFlux);
597   PFB("Store per-ring hists", fStorePerRing);
598   PFB("Make timing histogram", fDoTiming);
599   // gROOT->IncreaseDirLevel();
600   GetESDFixer()         .Print(option);        
601   GetSharingFilter()    .Print(option);
602   GetDensityCalculator().Print(option);
603   GetCorrections()      .Print(option);
604   GetHistCollector()    .Print(option);
605   GetEventPlaneFinder() .Print(option);
606   // gROOT->DecreaseDirLevel();
607   gROOT->DecreaseDirLevel();
608 }
609
610
611 //
612 // EOF
613 //