]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDCorrELossFit.cxx
remove buggy histos
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDCorrELossFit.cxx
1 // Object holding the Energy loss fit 'correction' 
2 // 
3 // These are generated from Monte-Carlo or real ESDs. 
4 //
5 #include "AliFMDCorrELossFit.h"
6 #include "AliForwardUtil.h"
7 #include <TF1.h>
8 #include <TBrowser.h>
9 #include <TVirtualPad.h>
10 #include <THStack.h>
11 #include <TH1D.h>
12 #include <AliLog.h>
13 #include <TMath.h>
14 #include <TList.h>
15 #include <iostream>
16 #include <iomanip>
17
18 Double_t AliFMDCorrELossFit::ELossFit::fgMaxRelError = .12;
19 Double_t AliFMDCorrELossFit::ELossFit::fgLeastWeight = 1e-5;
20 Double_t AliFMDCorrELossFit::ELossFit::fgMaxChi2nu   = 20;
21
22 //____________________________________________________________________
23 AliFMDCorrELossFit::ELossFit::ELossFit()
24   : fN(0),
25     fNu(0),
26     fChi2(0),
27     fC(0),
28     fDelta(0),
29     fXi(0),
30     fSigma(0),
31     fSigmaN(0),
32     fA(0),
33     fEC(0),
34     fEDelta(0),
35     fEXi(0),
36     fESigma(0),
37     fESigmaN(0),
38     fEA(0),
39     fQuality(0), 
40     fDet(0), 
41     fRing('\0'),
42     fBin(0)
43 {
44   //
45   // Default constructor 
46   // 
47   //
48 }
49 //____________________________________________________________________
50 AliFMDCorrELossFit::ELossFit::ELossFit(Int_t quality, const TF1& f)
51   : fN(f.GetNpar() > AliForwardUtil::ELossFitter::kN ? 
52        Int_t(f.GetParameter(AliForwardUtil::ELossFitter::kN)) : 
53        1),
54     fNu(f.GetNDF()),
55     fChi2(f.GetChisquare()),
56     fC(f.GetParameter(AliForwardUtil::ELossFitter::kC)),
57     fDelta(f.GetParameter(AliForwardUtil::ELossFitter::kDelta)),
58     fXi(f.GetParameter(AliForwardUtil::ELossFitter::kXi)),
59     fSigma(f.GetParameter(AliForwardUtil::ELossFitter::kSigma)),
60     fSigmaN(f.GetParameter(AliForwardUtil::ELossFitter::kSigmaN)),
61     fA(0),
62     fEC(f.GetParError(AliForwardUtil::ELossFitter::kC)),
63     fEDelta(f.GetParError(AliForwardUtil::ELossFitter::kDelta)),
64     fEXi(f.GetParError(AliForwardUtil::ELossFitter::kXi)),
65     fESigma(f.GetParError(AliForwardUtil::ELossFitter::kSigma)),
66     fESigmaN(f.GetParError(AliForwardUtil::ELossFitter::kSigmaN)),
67     fEA(0),
68     fQuality(quality),
69     fDet(0), 
70     fRing('\0'),
71     fBin(0)
72 {
73   // 
74   // Construct from a function
75   // 
76   // Parameters:
77   //    quality Quality flag
78   //    f       Function
79   //
80   if (fN <= 0) return;
81   fA  = new Double_t[fN];
82   fEA = new Double_t[fN];
83   for (Int_t i = 0; i < fN-1; i++) { 
84     fA[i]  = f.GetParameter(AliForwardUtil::ELossFitter::kA+i);
85     fEA[i] = f.GetParError(AliForwardUtil::ELossFitter::kA+i);
86   }
87   fA[fN-1]  = -9999;
88   fEA[fN-1] = -9999;
89 }
90
91 //____________________________________________________________________
92 AliFMDCorrELossFit::ELossFit::ELossFit(Int_t     quality,UShort_t  n, 
93                                        Double_t  chi2,   UShort_t  nu, 
94                                        Double_t  c,      Double_t  ec, 
95                                        Double_t  delta,  Double_t  edelta, 
96                                        Double_t  xi,     Double_t  exi,
97                                        Double_t  sigma,  Double_t  esigma, 
98                                        Double_t  sigman, Double_t  esigman, 
99                                        const Double_t* a,const Double_t* ea)
100   : fN(n),
101     fNu(nu),
102     fChi2(chi2),
103     fC(c),
104     fDelta(delta),
105     fXi(xi),
106     fSigma(sigma),
107     fSigmaN(sigman),
108     fA(0),
109     fEC(ec),
110     fEDelta(edelta),
111     fEXi(exi),
112     fESigma(esigma),
113     fESigmaN(esigman),
114     fEA(0),
115     fQuality(quality),
116     fDet(0), 
117     fRing('\0'),
118     fBin(0)
119 {
120   // 
121   // Constructor with full parameter set
122   // 
123   // Parameters:
124   //    quality   Quality flag
125   //    n         @f$ N@f$ - Number of fitted peaks
126   //    chi2      @f$ \chi^2 @f$
127   //    nu        @f$ \nu @f$ - number degrees of freedom
128   //    c         @f$ C@f$ - scale constant
129   //    ec        @f$ \delta C@f$ - error on @f$ C@f$ 
130   //    delta     @f$ \Delta@f$ - Most probable value             
131   //    edelta    @f$ \delta\Delta@f$ - error on @f$\Delta@f$ 
132   //    xi        @f$ \xi@f$ - width  
133   //    exi       @f$ \delta\xi@f$ - error on @f$\xi@f$ 
134   //    sigma     @f$ \sigma@f$ - Width of Gaussian                
135   //    esigma    @f$ \delta\sigma@f$ - error on @f$\sigma@f$ 
136   //    sigman    @f$ \sigma_n@f$ - Noise width           
137   //    esigman   @f$ \delta\sigma_n@f$ - error on @f$\sigma_n@f$ 
138   //    a         Array of @f$ N-1@f$ weights @f$ a_i@f$ for 
139   //                  @f$ i=2,\ldots@f$ 
140   //    ea        Array of @f$ N-1@f$ error on the weights @f$ a_i@f$ for 
141   //                  @f$ i=2,\ldots@f$ 
142   //
143   if (fN <= 0) return;
144   fA  = new Double_t[fN];
145   fEA = new Double_t[fN];
146   for (Int_t i = 0; i < fN-1; i++) { 
147     fA[i]  = a[i];
148     fEA[i] = ea[i];
149   }
150   fA[fN-1]  = -9999;
151   fEA[fN-1] = -9999;
152 }
153 //____________________________________________________________________
154 AliFMDCorrELossFit::ELossFit::ELossFit(const ELossFit& o)
155   : TObject(o), 
156     fN(o.fN),
157     fNu(o.fNu),
158     fChi2(o.fChi2),
159     fC(o.fC),
160     fDelta(o.fDelta),
161     fXi(o.fXi),
162     fSigma(o.fSigma),
163     fSigmaN(o.fSigmaN),
164     fA(0),
165     fEC(o.fEC),
166     fEDelta(o.fEDelta),
167     fEXi(o.fEXi),
168     fESigma(o.fESigma),
169     fESigmaN(o.fESigmaN),
170     fEA(0),
171     fQuality(o.fQuality),
172     fDet(o.fDet), 
173     fRing(o.fRing),
174     fBin(o.fBin)
175 {
176   // 
177   // Copy constructor 
178   // 
179   // Parameters:
180   //    o Object to copy from 
181   //
182   if (fN <= 0) return;
183   fA  = new Double_t[fN];
184   fEA = new Double_t[fN];
185   for (Int_t i = 0; i < fN-1; i++) { 
186     fA[i]  = o.fA[i];
187     fEA[i] = o.fEA[i];
188   }
189   fA[fN-1]  = -9999;
190   fEA[fN-1] = -9999;
191 }
192
193 //____________________________________________________________________
194 AliFMDCorrELossFit::ELossFit&
195 AliFMDCorrELossFit::ELossFit::operator=(const ELossFit& o)
196 {
197   // 
198   // Assignment operator 
199   // 
200   // Parameters:
201   //    o Object to assign from 
202   // 
203   // Return:
204   //    Reference to this object 
205   //
206   if (&o == this) return *this; 
207   fN       = o.fN;
208   fNu      = o.fNu;
209   fChi2    = o.fChi2;
210   fC       = o.fC;
211   fDelta   = o.fDelta;
212   fXi      = o.fXi;
213   fSigma   = o.fSigma;
214   fSigmaN  = o.fSigmaN;
215   fEC      = o.fEC;
216   fEDelta  = o.fEDelta;
217   fEXi     = o.fEXi;
218   fESigma  = o.fESigma;
219   fESigmaN = o.fESigmaN;
220   fQuality = o.fQuality;
221   fDet     = o.fDet; 
222   fRing    = o.fRing;
223   fBin     = o.fBin;
224   if (fA)  delete [] fA;
225   if (fEA) delete [] fEA; 
226   fA  = 0;
227   fEA = 0;
228
229   if (fN <= 0) return *this;
230   fA  = new Double_t[fN];
231   fEA = new Double_t[fN];
232   for (Int_t i = 0; i < fN; i++) { 
233     fA[i]  = o.fA[i];
234     fEA[i] = o.fEA[i];
235   }
236
237   return *this;
238 }
239
240 //____________________________________________________________________
241 AliFMDCorrELossFit::ELossFit::~ELossFit()
242 {
243   if (fA)  delete[] fA;
244   if (fEA) delete[] fEA;
245 }
246
247
248 //____________________________________________________________________
249 Int_t
250 AliFMDCorrELossFit::ELossFit::FindMaxWeight(Double_t maxRelError, 
251                                             Double_t leastWeight, 
252                                             UShort_t  maxN) const
253 {
254   // 
255   // Find the maximum weight to use.  The maximum weight is the
256   // largest i for which 
257   // 
258   // - @f$ i \leq \max{N}@f$ 
259   // - @f$ a_i > \min{a}@f$ 
260   // - @f$ \delta a_i/a_i > \delta_{max}@f$ 
261   // 
262   // Parameters:
263   //    maxRelError @f$ \min{a}@f$ 
264   //    leastWeight @f$ \delta_{max}@f$ 
265   //    maxN        @f$ \max{N}@f$      
266   // 
267   // Return:
268   //    The largest index @f$ i@f$ for which the above
269   // conditions hold.  Will never return less than 1. 
270   //
271   Int_t n = TMath::Min(maxN, UShort_t(fN-1));
272   Int_t m = 1;
273   // fN is one larger than we have data 
274   for (Int_t i = 0; i < n-1; i++, m++) { 
275     if (fA[i] < leastWeight)  break;
276     if (fEA[i] / fA[i] > maxRelError) break;
277   }
278   return m;
279 }
280
281 //____________________________________________________________________
282 Double_t 
283 AliFMDCorrELossFit::ELossFit::Evaluate(Double_t x, 
284                                        UShort_t maxN) const
285 {
286   // 
287   // Evaluate 
288   // @f[ 
289   //  f_N(x;\Delta,\xi,\sigma') = 
290   //     \sum_{i=1}^{n} a_i f(x;\Delta_i,\xi_i,\sigma_i')
291   // @f] 
292   //
293   // (see AliForwardUtil::NLandauGaus) for the maximum @f$ N @f$
294   // that fulfills the requirements 
295   // 
296   // Parameters:
297   //    x           Where to evaluate 
298   //    maxN      @f$ \max{N}@f$    
299   // 
300   // Return:
301   //    @f$ f_N(x;\Delta,\xi,\sigma')@f$ 
302   //
303   return AliForwardUtil::NLandauGaus(x, fDelta, fXi, fSigma, fSigmaN, 
304                                      TMath::Min(maxN, UShort_t(fN)), fA);
305 }
306
307 //____________________________________________________________________
308 Double_t 
309 AliFMDCorrELossFit::ELossFit::EvaluateWeighted(Double_t x, 
310                                                UShort_t maxN) const
311 {                                                                       
312   // 
313   // Evaluate 
314   // @f[ 
315   //   f_W(x;\Delta,\xi,\sigma') = 
316   //   \frac{\sum_{i=1}^{n} i a_i f_i(x;\Delta,\xi,\sigma')}{
317   //     f_N(x;\Delta,\xi,\sigma')} = 
318   //   \frac{\sum_{i=1}^{n} i a_i f(x;\Delta_i,\xi_i,\sigma_i')}{
319   //     \sum_{i=1}^{n} a_i f(x;\Delta_i,\xi_i,\sigma_i')}
320   // @f] 
321   // where @f$ n@f$ fulfills the requirements (see FindMaxWeight). 
322   //
323   // If the denominator is zero, then 1 is returned. 
324   //
325   // See also AliForwardUtil::ILandauGaus and AliForwardUtil::NLandauGaus
326   // for more information on the evaluated functions. 
327   // 
328   // Parameters:
329   //    x           Where to evaluate 
330   //    maxN      @f$ \max{N}@f$      
331   // 
332   // Return:
333   //    @f$ f_W(x;\Delta,\xi,\sigma')@f$.  
334   //
335   UShort_t n   = TMath::Min(maxN, UShort_t(fN-1));
336   Double_t num = 0;
337   Double_t den = 0;
338   for (Int_t i = 1; i <= n; i++) {
339     Double_t a = (i == 1 ? 1 : fA[i-1]);
340     if (fA[i-1] < 0) break;
341     Double_t f = AliForwardUtil::ILandauGaus(x,fDelta,fXi,fSigma,fSigmaN,i);
342     num += i * a * f;
343     den += a * f;
344   }
345   if (den <= 0) return 1;
346   return num / den;
347 }
348
349
350 #define OUTPAR(N,V,E)                   \
351   std::setprecision(9) <<               \
352   std::setw(12) << N << ": " <<         \
353   std::setw(14) << V << " +/- " <<      \
354   std::setw(14) << E << " (" <<         \
355   std::setprecision(-1) <<              \
356   std::setw(5)  << 100*(V>0?E/V:1) << "%)\n"
357   
358
359 //____________________________________________________________________
360 Int_t
361 AliFMDCorrELossFit::ELossFit::Compare(const TObject* o) const
362 {
363   // 
364   // Compare to another ELossFit object. 
365   // 
366   // - +1, if this quality is better (larger) than other objects quality
367   // - -1, if this quality is worse (smaller) than other objects quality
368   // - +1, if this @f$|\chi^2/\nu-1|@f$ is smaller than the same for other
369   // - -1, if this @f$|\chi^2/\nu-1|@f$ is larger than the same for other
370   // - 0 otherwise 
371   // 
372   // Parameters:
373   //    o Other object to compare to 
374   //
375   const ELossFit* other = static_cast<const ELossFit*>(o);
376   if (this->fQuality < other->fQuality) return -1;
377   if (this->fQuality > other->fQuality) return +1;
378   Double_t chi2nu  = (fNu == 0 ? 1000 : fChi2 / fNu);
379   Double_t oChi2nu = (other->fNu == 0 ? 1000 : other->fChi2 / other->fNu);
380   if (TMath::Abs(chi2nu-1) < TMath::Abs(oChi2nu-1)) return +1;
381   if (TMath::Abs(chi2nu-1) > TMath::Abs(oChi2nu-1)) return -1;
382   return 0;
383 }
384
385 //____________________________________________________________________
386 void
387 AliFMDCorrELossFit::ELossFit::Print(Option_t*) const
388 {
389   // 
390   // Information to standard output 
391   // 
392   // Parameters:
393   //    option Not used 
394   //
395   std::cout << GetName() << ":\n"
396             << " chi^2/nu = " << fChi2 << "/" << fNu << " = " 
397             << (fNu == 0 ? 999 : fChi2 / fNu) << "\n"
398             << "     Quality:   " << fQuality << "\n" 
399             << "  NParticles:   " << fN << "  (" << FindMaxWeight() << ")\n"
400             << OUTPAR("Delta", fDelta, fEDelta) 
401             << OUTPAR("xi", fXi, fEXi)
402             << OUTPAR("sigma", fSigma, fESigma)
403             << OUTPAR("sigma_n", fSigmaN, fESigmaN);
404   for (Int_t i = 0; i < fN-1; i++) 
405     std::cout << OUTPAR(Form("a%d", i+2), fA[i], fEA[i]);
406   std::cout << std::flush;
407 }
408
409 //____________________________________________________________________
410 const Char_t*
411 AliFMDCorrELossFit::ELossFit::GetName() const 
412 {
413   // 
414   // Get the name of this object 
415   // 
416   // 
417   // Return:
418   //    
419   //
420   return Form("FMD%d%c_etabin%03d", fDet, fRing, fBin);
421 }
422
423 //____________________________________________________________________
424 void
425 AliFMDCorrELossFit::ELossFit::Browse(TBrowser* b)
426 {
427   // 
428   // Browse this object 
429   // 
430   // Parameters:
431   //    b Browser
432   //
433   Draw(b ? b->GetDrawOption() : "comp");
434   gPad->SetLogy();
435   gPad->Update();
436 }
437      
438 //____________________________________________________________________
439 void
440 AliFMDCorrELossFit::ELossFit::Draw(Option_t* option)
441 {
442   // 
443   // Draw this fit 
444   // 
445   // Parameters:
446   //    option Options 
447   //  - COMP  Draw components too 
448   //
449   TString opt(option);
450   opt.ToUpper();
451   bool comp = false;
452   if (opt.Contains("COMP")) { 
453     opt.ReplaceAll("COMP","");
454     comp = true;
455   }
456   if (!opt.Contains("SAME")) { 
457     gPad->Clear();
458   }
459
460   TObjArray cleanup;
461   TF1* tot = AliForwardUtil::MakeNLandauGaus(1, 
462                                              fDelta, fXi, 
463                                              fSigma, fSigmaN, 
464                                              fN,     fA, 
465                                              0.01,   10);
466   tot->SetLineColor(kBlack);
467   tot->SetLineWidth(2);
468   tot->SetLineStyle(1);
469   tot->SetTitle(GetName());
470   Double_t max = tot->GetMaximum();
471
472   if (!opt.Contains("SAME")) {
473     TH1* frame = new TH1F(GetName(), 
474                           Form("FMD%d%c, eta bin %d",fDet,fRing,fBin),
475                           100, 0, 10);
476     frame->SetMinimum(max/10000);
477     frame->SetMaximum(max*1.4);
478     frame->SetXTitle("#Delta / #Delta_{mip}");
479     frame->Draw();
480     opt.Append(" SAME");
481   }
482   tot->DrawCopy(opt.Data());
483   cleanup.Add(tot);
484
485   if (!comp) { 
486     gPad->cd();
487     return;
488   }
489
490   Double_t min = max;
491   opt.Append(" same");
492   Int_t maxW = FindMaxWeight();
493   for (Int_t i=1; i <= fN; i++) { 
494     TF1* f = AliForwardUtil::MakeILandauGaus((i == 1 ? 1 : fA[i-2]), 
495                                              fDelta, fXi, 
496                                              fSigma, fSigmaN, 
497                                              i,      0.01, 10);
498     f->SetLineWidth(2);
499     f->SetLineStyle(i > maxW ? 2 : 1);
500     min = TMath::Min(f->GetMaximum(), min);
501     f->DrawCopy(opt.Data());
502     cleanup.Add(f);
503   }
504   min /= 100;
505   tot->GetHistogram()->SetMaximum(max);
506   tot->GetHistogram()->SetMinimum(min);
507   tot->GetHistogram()->GetYaxis()->SetRangeUser(min, max);
508
509   gPad->cd();
510 }
511
512                                                  
513 //____________________________________________________________________
514 #define CHECKPAR(V,E,T) ((V > 0) && (E / V < T))
515
516 //____________________________________________________________________
517 Double_t
518 AliFMDCorrELossFit::ELossFit::GetLowerBound(Double_t f) const
519 {
520   // 
521   // Return 
522   //    Delta * f
523   return f * fDelta;
524 }
525 //____________________________________________________________________
526 Double_t
527 AliFMDCorrELossFit::ELossFit::GetLowerBound(Double_t f, 
528                                             Bool_t includeSigma) const
529 {
530   // 
531   // Return 
532   //    Delta - f * (xi + sigma)
533   return fDelta - f * (fXi + (includeSigma ? fSigma : 0));
534 }
535
536 //____________________________________________________________________
537 void 
538 AliFMDCorrELossFit::ELossFit::CalculateQuality(Double_t maxChi2nu, 
539                                                Double_t maxRelError, 
540                                                Double_t leastWeight)
541 {
542   // 
543   // Calculate the quality 
544   //
545   Int_t qual = 0;
546   if (fNu > 0 && fChi2 / fNu < maxChi2nu) qual += 4;;
547   if (CHECKPAR(fDelta,  fEDelta,  maxRelError)) qual++;
548   if (CHECKPAR(fXi,     fEXi,     maxRelError)) qual++;
549   if (CHECKPAR(fSigma,  fESigma,  maxRelError)) qual++;
550   if (CHECKPAR(fSigmaN, fESigmaN, maxRelError)) qual++;
551   qual += FindMaxWeight(1.5*maxRelError, leastWeight, fN);
552   fQuality = qual;
553 }
554
555 //____________________________________________________________________
556 AliFMDCorrELossFit::AliFMDCorrELossFit()
557   : TObject(), 
558     fRings(), 
559     fEtaAxis(0,0,0), 
560     fLowCut(0)
561 {
562   // 
563   // Default constructor 
564   //
565   fRings.SetOwner(kTRUE);
566   fEtaAxis.SetTitle("#eta");
567   fEtaAxis.SetName("etaAxis");
568   fRings.SetName("rings");
569 }
570
571 //____________________________________________________________________
572 AliFMDCorrELossFit::AliFMDCorrELossFit(const AliFMDCorrELossFit& o)
573   : TObject(o), 
574     fRings(o.fRings),
575     fEtaAxis(o.fEtaAxis.GetNbins(),o.fEtaAxis.GetXmin(),o.fEtaAxis.GetXmax()), 
576     fLowCut(0)
577 {
578   // 
579   // Copy constructor 
580   // 
581   // Parameters:
582   //    o Object to copy from 
583   //
584   fEtaAxis.SetTitle("#eta");
585   fEtaAxis.SetName("etaAxis");
586 }
587
588 //____________________________________________________________________
589 AliFMDCorrELossFit::~AliFMDCorrELossFit()
590 {
591   // 
592   // Destructor 
593   //
594   fRings.Clear();
595 }
596
597 //____________________________________________________________________
598 AliFMDCorrELossFit&
599 AliFMDCorrELossFit::operator=(const AliFMDCorrELossFit& o)
600 {
601   // 
602   // Assignment operator 
603   // 
604   // Parameters:
605   //    o Object to assign from 
606   // 
607   // Return:
608   //    Reference to this object 
609   //
610   if (&o == this) return *this; 
611   fRings = o.fRings;
612   fLowCut = o.fLowCut;
613   SetEtaAxis(o.fEtaAxis.GetNbins(), o.fEtaAxis.GetXmin(), o.fEtaAxis.GetXmax());
614
615   return *this;
616 }
617 //____________________________________________________________________
618 Int_t
619 AliFMDCorrELossFit::FindEtaBin(Double_t eta) const
620 {
621   // 
622   // Find the eta bin corresponding to the given eta 
623   // 
624   // Parameters:
625   //    eta  Eta value 
626   // 
627   // Return:
628   //    Bin (in @f$[1,N_{bins}]@f$) corresponding to the given
629   // eta, or 0 if out of range.
630   //
631   if (TMath::Abs(fEtaAxis.GetXmin() - fEtaAxis.GetXmax()) < 1e-6 
632       || fEtaAxis.GetNbins() == 0) {
633     AliWarning("No eta axis defined");
634     return -1;
635   }
636   Int_t bin = const_cast<TAxis&>(fEtaAxis).FindBin(eta);
637   if (bin <= 0 || bin > fEtaAxis.GetNbins()) return 0;
638   return bin;
639 }
640
641 //____________________________________________________________________
642 Bool_t
643 AliFMDCorrELossFit::SetFit(UShort_t d, Char_t r, Int_t etaBin, ELossFit* fit)
644 {
645   // 
646   // Set the fit parameters from a function 
647   // 
648   // Parameters:
649   //    d       Detector
650   //    r       Ring 
651   //    etaBin  Eta (bin number, 1->nBins)
652   //    f       ELoss fit result - note, the object will take ownership
653   //  
654   TObjArray* ringArray = GetOrMakeRingArray(d, r);
655   if (!ringArray) { 
656     AliError(Form("Failed to make ring array for FMD%d%c", d, r));
657     return kFALSE;
658   }
659   if (etaBin <= 0 || etaBin >= fEtaAxis.GetNbins()+1) { 
660     AliError(Form("bin=%d is out of range [%d,%d]", 
661                   etaBin, 1, fEtaAxis.GetNbins()));
662     return kFALSE;
663   }
664   // AliInfo(Form("Adding fit %p at %3d", fit, etaBin));
665   ringArray->AddAtAndExpand(fit, etaBin);
666   return kTRUE;
667 }
668
669 //____________________________________________________________________
670 Bool_t
671 AliFMDCorrELossFit::SetFit(UShort_t d, Char_t r, Double_t eta, ELossFit* fit)
672 {
673   // 
674   // Set the fit parameters from a function 
675   // 
676   // Parameters:
677   //    d    Detector
678   //    r    Ring 
679   //    eta  Eta 
680   //    f    ELoss fit result - note, the object will take ownership
681   //  
682   Int_t bin = FindEtaBin(eta);
683   if (bin <= 0) { 
684     AliError(Form("eta=%f is out of range [%f,%f]", 
685                   eta, fEtaAxis.GetXmin(), fEtaAxis.GetXmax()));
686     return kFALSE;
687   }
688
689   return SetFit(d, r, bin, fit);
690 }
691 //____________________________________________________________________
692 Bool_t
693 AliFMDCorrELossFit::SetFit(UShort_t  d,      Char_t    r, 
694                            Double_t  eta, 
695                            Int_t     quality,UShort_t  n, 
696                            Double_t  chi2,   UShort_t  nu, 
697                            Double_t  c,      Double_t  ec, 
698                            Double_t  delta,  Double_t  edelta, 
699                            Double_t  xi,     Double_t  exi,
700                            Double_t  sigma,  Double_t  esigma, 
701                            Double_t  sigman, Double_t  esigman, 
702                            Double_t* a,      Double_t* ea)
703 {
704   // 
705   // Set the fit parameters from a function 
706   // 
707   // Parameters:
708   //    d         Detector number
709   //    r         Ring identifier 
710   //    eta       Eta value
711   //    quality   Quality flag
712   //    n         @f$ N@f$ - Number of fitted peaks
713   //    chi2      @f$ \chi^2 @f$
714   //    nu        @f$ \nu @f$ - number degrees of freedom
715   //    c         @f$ C@f$ - scale constant
716   //    ec        @f$ \delta C@f$ - error on @f$ C@f$ 
717   //    delta     @f$ \Delta@f$ - most probable value
718   //    edelta    @f$ \delta\Delta@f$ - error on @f$\Delta@f$ 
719   //    xi        @f$ \xi@f$ - Landau width               
720   //    exi       @f$ \delta\xi@f$ - error on @f$\xi@f$ 
721   //    sigma     @f$ \sigma@f$ - Gaussian width
722   //    esigma    @f$ \delta\sigma@f$ - error on @f$\sigma@f$ 
723   //    sigman    @f$ \sigma_n@f$ - Noise width           
724   //    esigman   @f$ \delta\sigma_n@f$ - error on @f$\sigma_n@f$ 
725   //    a         Array of @f$ N-1@f$ weights @f$ a_i@f$ for 
726   //                  @f$ i=2,\ldots@f$ 
727   //    ea        Array of @f$ N-1@f$ errors on weights @f$ a_i@f$ for 
728   //                  @f$ i=2,\ldots@f$ 
729   //
730   ELossFit* e = new ELossFit(quality, n, 
731                              chi2,    nu,
732                              c,       ec,
733                              delta,   edelta,
734                              xi,      exi,
735                              sigma,   esigma,
736                              sigman,  esigman,
737                              a,       ea);
738   if (!SetFit(d, r, eta, e)) { 
739     delete e;
740     return kFALSE;
741   }
742   return kTRUE;
743 }
744 //____________________________________________________________________
745 Bool_t
746 AliFMDCorrELossFit::SetFit(UShort_t  d, Char_t r, Double_t eta, 
747                            Int_t quality, const TF1& f)
748 {
749   // 
750   // Set the fit parameters from a function 
751   // 
752   // Parameters:
753   //    d        Detector
754   //    r        Ring 
755   //    eta      Eta 
756   //    quality  Quality flag
757   //    f        Function from fit 
758   //  
759   ELossFit* e = new ELossFit(quality, f);
760   if (!SetFit(d, r, eta, e)) { 
761     delete e;
762     return kFALSE;
763   }
764   return kTRUE;
765 }
766 //____________________________________________________________________
767 AliFMDCorrELossFit::ELossFit*
768 AliFMDCorrELossFit::GetFit(UShort_t  d, Char_t r, Int_t etabin) const
769 {
770   // 
771   // Get the fit corresponding to the specified parameters 
772   // 
773   // Parameters:
774   //    d      Detector 
775   //    r      Ring 
776   //    etabin Eta bin (1 based)
777   // 
778   // Return:
779   //    Fit parameters or null in case of problems 
780   //
781   TObjArray* ringArray = GetRingArray(d, r);
782   if (!ringArray)                                   return 0;
783   if (etabin <= 0 || etabin >= fEtaAxis.GetNbins()) return 0;
784   if      (etabin > ringArray->GetEntriesFast())    return 0;
785   else if (etabin >= ringArray->GetEntriesFast())   etabin--;
786   else if (!ringArray->At(etabin))                  etabin++;
787   return static_cast<ELossFit*>(ringArray->At(etabin));
788 }
789 //____________________________________________________________________
790 AliFMDCorrELossFit::ELossFit*
791 AliFMDCorrELossFit::GetFit(UShort_t  d, Char_t r, Double_t eta) const
792 {
793   // 
794   // Find the fit corresponding to the specified parameters 
795   // 
796   // Parameters:
797   //    d   Detector 
798   //    r   Ring 
799   //    eta Eta value 
800   // 
801   // Return:
802   //    Fit parameters or null in case of problems 
803   //
804   Int_t etabin = FindEtaBin(eta);
805   return GetFit(d, r, etabin);
806 }
807
808 //____________________________________________________________________
809 AliFMDCorrELossFit::ELossFit*
810 AliFMDCorrELossFit::FindFit(UShort_t  d, Char_t r, Int_t etabin) const
811 {
812   // 
813   // Find the fit corresponding to the specified parameters 
814   // 
815   // Parameters:
816   //    d      Detector 
817   //    r      Ring 
818   //    etabin Eta bin (1 based)
819   // 
820   // Return:
821   //    Fit parameters or null in case of problems 
822   //
823   TObjArray* ringArray = GetRingArray(d, r);
824   if (!ringArray) { 
825     AliError(Form("Failed to make ring array for FMD%d%c", d, r));
826     return 0;
827   }
828   if (etabin <= 0 || etabin >= fEtaAxis.GetNbins()) { 
829     // AliError(Form("Eta bin=%3d out of bounds [%d,%d] for FMD%d%c", 
830     //               etabin, 1, fEtaAxis.GetNbins(), d, r));
831     return 0;
832   }
833   if (etabin > ringArray->GetEntriesFast()) { 
834     // AliError(Form("Eta bin=%3d out of bounds [%d,%d] for FMD%d%c", 
835     //               etabin, 1, ringArray->GetEntriesFast(), d, r));
836     return 0;
837   }
838   else if (etabin >= ringArray->GetEntriesFast()) { 
839     // AliWarning(Form("Eta bin=%3d out of bounds by +1 [%d,%d] for FMD%d%c, " 
840     //              "trying %3d", etabin, 1, ringArray->GetEntriesFast(), d, r,
841     //              etabin-1));
842     etabin--;
843   }
844   else if (!ringArray->At(etabin)) { 
845     // AliWarning(Form("Eta bin=%d has no fit for FMD%d%c, trying %03d", 
846     //              etabin, d, r, etabin+1));
847     etabin++;
848   }
849   return static_cast<ELossFit*>(ringArray->At(etabin));
850 }
851 //____________________________________________________________________
852 AliFMDCorrELossFit::ELossFit*
853 AliFMDCorrELossFit::FindFit(UShort_t  d, Char_t r, Double_t eta) const
854 {
855   // 
856   // Find the fit corresponding to the specified parameters 
857   // 
858   // Parameters:
859   //    d   Detector 
860   //    r   Ring 
861   //    eta Eta value 
862   // 
863   // Return:
864   //    Fit parameters or null in case of problems 
865   //
866   Int_t etabin = FindEtaBin(eta);
867   return FindFit(d, r, etabin);
868 }
869 //____________________________________________________________________
870 TObjArray*
871 AliFMDCorrELossFit::GetRingArray(UShort_t d, Char_t r) const
872 {
873   // 
874   // Get the ring array corresponding to the specified ring
875   // 
876   // Parameters:
877   //    d Detector 
878   //    r Ring 
879   // 
880   // Return:
881   //    Pointer to ring array, or null in case of problems
882   //
883   Int_t idx = -1;
884   switch (d) { 
885   case 1:   idx = 0; break;
886   case 2:   idx = (r == 'i' || r == 'I') ? 1 : 2; break;
887   case 3:   idx = (r == 'o' || r == 'I') ? 3 : 4; break;
888   }
889   if (idx < 0 || idx >= fRings.GetEntriesFast()) return 0;
890   return static_cast<TObjArray*>(fRings.At(idx));
891 }
892 //____________________________________________________________________
893 TObjArray*
894 AliFMDCorrELossFit::GetOrMakeRingArray(UShort_t d, Char_t r)
895 {
896   // 
897   // Get the ring array corresponding to the specified ring
898   // 
899   // Parameters:
900   //    d Detector 
901   //    r Ring 
902   // 
903   // Return:
904   //    Pointer to ring array, or newly created container 
905   //
906   Int_t idx = -1;
907   switch (d) { 
908   case 1:   idx = 0; break;
909   case 2:   idx = (r == 'i' || r == 'I') ? 1 : 2; break;
910   case 3:   idx = (r == 'o' || r == 'I') ? 3 : 4; break;
911   }
912   if (idx < 0) return 0;
913   if (idx >= fRings.GetEntriesFast() || !fRings.At(idx)) {
914     TObjArray* a = new TObjArray(0);
915     // TOrdCollection* a = new TOrdCollection(fEtaAxis.GetNbins());
916     a->SetName(Form("FMD%d%c", d, r));
917     a->SetOwner();
918     fRings.AddAtAndExpand(a, idx);
919   }
920   return static_cast<TObjArray*>(fRings.At(idx));
921 }
922
923 //____________________________________________________________________
924 Double_t
925 AliFMDCorrELossFit::GetLowerBound(UShort_t  d, Char_t r, Int_t etabin,
926                                   Double_t f) const
927 {
928   ELossFit* fit = GetFit(d, r, etabin);
929   if (!fit) return -1024;
930   return fit->GetLowerBound(f);
931 }
932 //____________________________________________________________________
933 Double_t
934 AliFMDCorrELossFit::GetLowerBound(UShort_t  d, Char_t r, Double_t eta,
935                                   Double_t f) const
936 {
937   Int_t bin = FindEtaBin(eta);
938   if (bin <= 0) return -1024;
939   return GetLowerBound(d, r, Int_t(bin), f);
940 }
941 //____________________________________________________________________
942 Double_t
943 AliFMDCorrELossFit::GetLowerBound(UShort_t  d, Char_t r, Int_t etabin,
944                                   Double_t f, Bool_t showErrors, 
945                                   Bool_t includeSigma) const
946 {
947   ELossFit* fit = GetFit(d, r, etabin);
948   if (!fit) { 
949     if (showErrors) {
950       AliWarning(Form("No fit for FMD%d%c @ etabin=%d", d, r, etabin));
951     }
952     return -1024;
953   }
954   return fit->GetLowerBound(f, includeSigma);
955 }
956
957 //____________________________________________________________________
958 Double_t
959 AliFMDCorrELossFit::GetLowerBound(UShort_t  d, Char_t r, Double_t eta,
960                                   Double_t f, Bool_t showErrors, 
961                                   Bool_t includeSigma) const
962 {
963   Int_t bin = FindEtaBin(eta);
964   if (bin <= 0) { 
965     if (showErrors)
966       AliError(Form("eta=%f out of bounds for FMD%d%c", eta, d, r));
967     return -1024;
968   }
969   return GetLowerBound(d, r, bin, f, showErrors, includeSigma);
970 }
971
972 //____________________________________________________________________
973 namespace { 
974   TH1D* MakeHist(const TAxis& axis, const char* name, const char* title, 
975                  Int_t color)
976   {
977     TH1D* h = new TH1D(Form("%s_%s", name, title), 
978                        Form("%s %s", name, title), 
979                        axis.GetNbins(), axis.GetXmin(), axis.GetXmax());
980     h->SetDirectory(0);
981     h->SetMarkerStyle(20);
982     h->SetMarkerColor(color);
983     h->SetMarkerSize(0.5);
984     h->SetFillColor(color);
985     h->SetFillStyle(3001);
986     h->SetLineColor(color);
987     return h;
988   }
989 }
990
991   
992
993 #define IDX2RING(I) (i == 0 || i == 1 || i == 3 ? 'I' : 'O')
994 #define IDX2DET(I)  (i == 0 ? 1 : (i == 1 || i == 2 ? 2 : 3))
995 //____________________________________________________________________
996 TList*
997 AliFMDCorrELossFit::GetStacks(Bool_t err, Bool_t rel, UShort_t maxN) const
998 {
999   // Get a list of THStacks 
1000   Int_t nRings = fRings.GetEntriesFast();
1001   Int_t nPad   = 6+maxN-1; // 7 regular params, and maxN-1 weights
1002
1003   enum { 
1004     kChi2nu = 0, 
1005     kC      = 1, 
1006     kDelta  = 2, 
1007     kXi     = 3, 
1008     kSigma  = 4, 
1009     kN      = 5 
1010   };
1011   
1012   THStack* sChi2nu;
1013   THStack* sC;
1014   THStack* sDelta;
1015   THStack* sXi;
1016   THStack* sSigma;
1017   // THStack* sigman;
1018   THStack* n;
1019   TList* stacks = new TList;
1020   stacks->AddAt(sChi2nu= new THStack("chi2",   "#chi^{2}/#nu"), kChi2nu);
1021   stacks->AddAt(sC     = new THStack("c",       "C"),           kC);
1022   stacks->AddAt(sDelta = new THStack("delta",  "#Delta_{mp}"),  kDelta);
1023   stacks->AddAt(sXi    = new THStack("xi",     "#xi"),          kXi);
1024   stacks->AddAt(sSigma = new THStack("sigma",  "#sigma"),       kSigma);
1025   //stacks->AddAt(sigman= new THStack("sigman", "#sigma_{n}"),   5);
1026   stacks->AddAt(n     = new THStack("n",      "N"),            kN);
1027   for (Int_t i = 1; i <= maxN; i++) {
1028     stacks->AddAt(new THStack(Form("a_%02d", i+1), Form("a_{%d}", i+1)), kN+i);
1029   }
1030   
1031   TArrayD min(nPad);
1032   TArrayD max(nPad);
1033   min.Reset(100000);
1034   max.Reset(-100000);
1035
1036   for (Int_t i = 0; i < nRings; i++) { 
1037     if (!fRings.At(i)) continue;
1038     TObjArray* a     = static_cast<TObjArray*>(fRings.At(i));
1039     Int_t      nFits = a->GetEntriesFast();
1040     Int_t      color = AliForwardUtil::RingColor(IDX2DET(i), IDX2RING(i));
1041
1042     TH1D* hChi    = MakeHist(fEtaAxis,a->GetName(), "chi2",   color);
1043     TH1D* hC      = MakeHist(fEtaAxis,a->GetName(), "c",      color);
1044     TH1D* hDelta  = MakeHist(fEtaAxis,a->GetName(), "delta",  color);
1045     TH1D* hXi     = MakeHist(fEtaAxis,a->GetName(), "xi",     color);
1046     TH1D* hSigma  = MakeHist(fEtaAxis,a->GetName(), "sigma",  color);
1047     // TH1D* hSigmaN = MakeHist(fEtaAxis,a->GetName(), "sigman", color);
1048     TH1D* hN      = MakeHist(fEtaAxis,a->GetName(), "n",      color);
1049     TH1D* hA[maxN];
1050     const char* ho = (rel || !err ? "hist" : "e");
1051     sChi2nu->Add(hChi,   "hist"); // 0
1052     sC     ->Add(hC,     ho);     // 1
1053     sDelta ->Add(hDelta, ho);     // 2
1054     sXi    ->Add(hXi,    ho);     // 3
1055     sSigma ->Add(hSigma, ho);     // 4
1056     // sigman->Add(hSigmaN,ho);     // 5
1057     n     ->Add(hN,     "hist"); // 5
1058     hChi->SetFillColor(color);
1059     hChi->SetFillStyle(3001);
1060     hN->SetFillColor(color);
1061     hN->SetFillStyle(3001);
1062
1063     for (Int_t k = 1; k <= maxN; k++) { 
1064       hA[k-1] = MakeHist(fEtaAxis,a->GetName(), Form("a%02d",k+1), color);
1065       static_cast<THStack*>(stacks->At(kN+k))->Add(hA[k-1], ho);
1066     }
1067                           
1068     for (Int_t j = 0; j < nFits; j++) {
1069       ELossFit* f = static_cast<ELossFit*>(a->At(j));
1070       if (!f) continue;
1071
1072       Int_t     b      = f->fBin;
1073       Int_t     nW     = f->FindMaxWeight();
1074       Double_t  vChi2nu = (f->fNu <= 0 ? 0 : f->fChi2 / f->fNu);
1075       Double_t  vC      = (rel ? (f->fC     >0 ?f->fEC     /f->fC      :0) 
1076                            : f->fC);
1077       Double_t  vDelta  = (rel ? (f->fDelta >0 ?f->fEDelta /f->fDelta  :0) 
1078                            : f->fDelta);
1079       Double_t  vXi     = (rel ? (f->fXi    >0 ?f->fEXi    /f->fXi     :0) 
1080                            : f->fXi);
1081       Double_t  vSigma  = (rel ? (f->fSigma >0 ?f->fESigma /f->fSigma  :0) 
1082                            : f->fSigma);
1083       // Double_t  sigman = (rel ? (f->fSigmaN>0 ?f->fESigmaN/f->fSigmaN :0) 
1084       //                     : f->SigmaN); 
1085       hChi   ->SetBinContent(b, vChi2nu);
1086       hN     ->SetBinContent(b, nW);
1087       hC     ->SetBinContent(b, vC);
1088       hDelta ->SetBinContent(b, vDelta);
1089       hXi    ->SetBinContent(b, vXi);
1090       hSigma ->SetBinContent(b, vSigma);
1091       if (vChi2nu > 1e-12) {
1092         min[kChi2nu] = TMath::Min(min[kChi2nu], vChi2nu);
1093         max[kChi2nu] = TMath::Max(max[kChi2nu], vChi2nu);
1094       }
1095       if (vC > 1e-12) {
1096         min[kC] = TMath::Min(min[kC], vC);
1097         max[kC] = TMath::Max(max[kC], vC);
1098       }
1099       if (vDelta > 1e-12) {
1100         min[kDelta] = TMath::Min(min[kDelta], vDelta);
1101         max[kDelta] = TMath::Max(max[kDelta], vDelta);
1102       }
1103       if (vXi > 1e-12) {
1104         min[kXi] = TMath::Min(min[kXi], vXi);
1105         max[kXi] = TMath::Max(max[kXi], vXi);
1106       }
1107       if (vSigma > 1e-12) {
1108         min[kSigma] = TMath::Min(min[kSigma], vSigma);
1109         max[kSigma] = TMath::Max(max[kSigma], vSigma);
1110       }
1111       if (nW > 1e-12) { 
1112         min[kN] = TMath::Min(min[kN], Double_t(nW));
1113         max[kN] = TMath::Max(max[kN], Double_t(nW));
1114       }
1115       // hSigmaN->SetBinContent(b,  sigman);
1116       if (!rel) {
1117         hC     ->SetBinError(b, f->fEC);
1118         hDelta ->SetBinError(b, f->fEDelta);
1119         hXi    ->SetBinError(b, f->fEXi);
1120         hSigma ->SetBinError(b, f->fESigma);
1121         // hSigmaN->SetBinError(b, f->fESigmaN);
1122       }
1123       for (Int_t k = 0; k < f->fN-1 && k < maxN; k++) { 
1124         Double_t vA = (rel ? (f->fA[k] > 0 ? f->fEA[k] / f->fA[k] : 0) 
1125                        : f->fA[k]);
1126         hA[k]->SetBinContent(b, vA);
1127         if (!rel) hA[k]->SetBinError(b, f->fEA[k]);
1128         if (vA > 1e-12) {
1129           min[kN+1+k] = TMath::Min(min[kN+1+k], vA);
1130           max[kN+1+k] = TMath::Max(max[kN+1+k], vA);
1131         }
1132       }
1133     }
1134   }
1135   return stacks;
1136 }
1137
1138 //____________________________________________________________________
1139 void
1140 AliFMDCorrELossFit::Draw(Option_t* option)
1141 {
1142   // 
1143   // Draw this object 
1144   // 
1145   // Parameters:
1146   //    option Options.  Possible values are 
1147   //  - err Plot error bars 
1148   //
1149   TString opt(Form("nostack %s", option));
1150   opt.ToLower();
1151   Bool_t  rel = (opt.Contains("relative"));
1152   Bool_t  err = (opt.Contains("error"));
1153   if (rel) opt.ReplaceAll("relative","");
1154   if (err) opt.ReplaceAll("error","");
1155
1156   UShort_t maxN   = 0;
1157   Int_t nRings = fRings.GetEntriesFast();
1158   for (Int_t i = 0; i < nRings; i++) { 
1159     if (!fRings.At(i)) continue;
1160     TObjArray* a     = static_cast<TObjArray*>(fRings.At(i));
1161     Int_t      nFits = a->GetEntriesFast();
1162
1163     for (Int_t j = 0; j < nFits; j++) {
1164       ELossFit* fit = static_cast<ELossFit*>(a->At(j));
1165       if (!fit) continue;
1166       maxN          = TMath::Max(maxN, UShort_t(fit->fN));
1167     }
1168   }
1169   // AliInfo(Form("Maximum N is %d", maxN));
1170   Int_t nPad = 6+maxN-1; // 7 regular params, and maxN-1 weights
1171   TVirtualPad* pad = gPad;
1172   pad->Divide(2, (nPad+1)/2, 0.1, 0, 0);
1173
1174   TList* stacks = GetStacks(err, rel, maxN);
1175
1176   Int_t nPad2 = (nPad+1) / 2;
1177   for (Int_t i = 0; i < nPad; i++) {
1178     Int_t iPad = 1 + i/nPad2 + 2 * (i % nPad2);
1179     TVirtualPad* p = pad->cd(iPad);
1180     p->SetLeftMargin(.15);
1181     p->SetFillColor(0);
1182     p->SetFillStyle(0);
1183     p->SetGridx();
1184     p->SetGridy();
1185     if (rel && i != 0 && i != 6 && i != 5 && i != 4) p->SetLogy();
1186
1187
1188     THStack* stack = static_cast<THStack*>(stacks->At(i));
1189
1190     // Double_t powMax = TMath::Log10(max[i]);
1191     // Double_t powMin = min[i] <= 0 ? powMax : TMath::Log10(min[i]);
1192     // if (powMax-powMin > 2. && min[i] != 0) p->SetLogy();
1193
1194     // stack->SetMinimum(min[i]);
1195     // stack->SetMaximum(max[i]);
1196     stack->Draw(opt.Data());
1197
1198     TString tit(stack->GetTitle());
1199     if (rel && i != 0 && i != 6)
1200       tit = Form("#delta %s/%s", tit.Data(), tit.Data());
1201     TH1*   hist  = stack->GetHistogram();
1202     TAxis* yaxis = hist->GetYaxis();
1203     yaxis->SetTitle(tit.Data());
1204     yaxis->SetTitleSize(0.15);
1205     yaxis->SetLabelSize(0.08);
1206     yaxis->SetTitleOffset(0.35);
1207     yaxis->SetTitleFont(132);
1208     yaxis->SetLabelFont(132);
1209     yaxis->SetNdivisions(5);
1210
1211
1212     TAxis* xaxis = stack->GetHistogram()->GetXaxis();
1213     xaxis->SetTitle("#eta");
1214     xaxis->SetTitleSize(0.15);
1215     xaxis->SetLabelSize(0.08);
1216     xaxis->SetTitleOffset(0.35);
1217     xaxis->SetTitleFont(132);
1218     xaxis->SetLabelFont(132);
1219     xaxis->SetNdivisions(10);
1220
1221     stack->Draw(opt.Data());
1222   }
1223   pad->cd();      
1224 }
1225
1226 //____________________________________________________________________
1227 void
1228 AliFMDCorrELossFit::Print(Option_t* option) const
1229 {
1230   // 
1231   // Print this object.  
1232   // 
1233   // Parameters:
1234   //    option Options 
1235   //   - R   Print recursive  
1236   //
1237   //
1238   TString opt(option);
1239   opt.ToUpper();
1240   Int_t nRings = fRings.GetEntriesFast();
1241   bool recurse = opt.Contains("R");
1242
1243   std::cout << "Low cut in fit range: " << fLowCut << "\n"
1244             << "Eta axis:             " << fEtaAxis.GetNbins() 
1245             << " bins, range [" << fEtaAxis.GetXmin() << "," 
1246             << fEtaAxis.GetXmax() << "]" << std::endl;
1247   
1248   for (Int_t i = 0; i < nRings; i++) { 
1249     if (!fRings.At(i)) continue;
1250     TObjArray* a     = static_cast<TObjArray*>(fRings.At(i));
1251     Int_t      nFits = a->GetEntriesFast();
1252
1253     std::cout << a->GetName() << " [" << nFits << " entries]" 
1254               << (recurse ? ":\n" : "\t");
1255     Int_t min = fEtaAxis.GetNbins()+1;
1256     Int_t max = 0;
1257     for (Int_t j = 0; j < nFits; j++) {
1258       if (!a->At(j)) continue;
1259       
1260       min = TMath::Min(j, min);
1261       max = TMath::Max(j, max);
1262
1263       if (recurse) {
1264         std::cout << "Bin # " << j << "\t";
1265         ELossFit* fit = static_cast<ELossFit*>(a->At(j));
1266         fit->Print(option);
1267       }
1268     }
1269     if (!recurse) 
1270       std::cout << " bin range: " << std::setw(3) << min 
1271                 << "-" << std::setw(3) << max << " " << std::setw(3) 
1272                 << (max-min+1) << " bins" << std::endl;
1273   }
1274 }
1275
1276 //____________________________________________________________________
1277 void
1278 AliFMDCorrELossFit::Browse(TBrowser* b)
1279 {
1280   // 
1281   // Browse this object 
1282   // 
1283   // Parameters:
1284   //    b 
1285   //
1286   b->Add(&fRings);
1287   b->Add(&fEtaAxis);
1288 }
1289
1290
1291
1292 //____________________________________________________________________
1293 //
1294 // EOF
1295 //