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