1 // Object holding the Energy loss fit 'correction'
3 // These are generated from Monte-Carlo or real ESDs.
4 #include "AliFMDCorrELossFit.h"
5 #include "AliForwardUtil.h"
8 #include <TVirtualPad.h>
16 Double_t AliFMDCorrELossFit::ELossFit::fgMaxRelError = .2;
17 Double_t AliFMDCorrELossFit::ELossFit::fgLeastWeight = 1e-5;
18 Double_t AliFMDCorrELossFit::ELossFit::fgMaxChi2nu = 20;
20 //____________________________________________________________________
21 AliFMDCorrELossFit::ELossFit::ELossFit()
43 // Default constructor
47 //____________________________________________________________________
48 AliFMDCorrELossFit::ELossFit::ELossFit(Int_t quality, const TF1& f)
49 : fN(f.GetNpar() > AliForwardUtil::ELossFitter::kN ?
50 f.GetParameter(AliForwardUtil::ELossFitter::kN) :
53 fChi2(f.GetChisquare()),
54 fC(f.GetParameter(AliForwardUtil::ELossFitter::kC)),
55 fDelta(f.GetParameter(AliForwardUtil::ELossFitter::kDelta)),
56 fXi(f.GetParameter(AliForwardUtil::ELossFitter::kXi)),
57 fSigma(f.GetParameter(AliForwardUtil::ELossFitter::kSigma)),
58 fSigmaN(f.GetParameter(AliForwardUtil::ELossFitter::kSigmaN)),
60 fEC(f.GetParError(AliForwardUtil::ELossFitter::kC)),
61 fEDelta(f.GetParError(AliForwardUtil::ELossFitter::kDelta)),
62 fEXi(f.GetParError(AliForwardUtil::ELossFitter::kXi)),
63 fESigma(f.GetParError(AliForwardUtil::ELossFitter::kSigma)),
64 fESigmaN(f.GetParError(AliForwardUtil::ELossFitter::kSigmaN)),
72 // Construct from a function
75 // quality Quality flag
79 fA = new Double_t[fN];
80 fEA = new Double_t[fN];
81 for (Int_t i = 0; i < fN-1; i++) {
82 fA[i] = f.GetParameter(AliForwardUtil::ELossFitter::kA+i);
83 fEA[i] = f.GetParError(AliForwardUtil::ELossFitter::kA+i);
89 //____________________________________________________________________
90 AliFMDCorrELossFit::ELossFit::ELossFit(Int_t quality,UShort_t n,
91 Double_t chi2, UShort_t nu,
92 Double_t c, Double_t ec,
93 Double_t delta, Double_t edelta,
94 Double_t xi, Double_t exi,
95 Double_t sigma, Double_t esigma,
96 Double_t sigman, Double_t esigman,
97 Double_t* a, Double_t* ea)
119 // Constructor with full parameter set
122 // quality Quality flag
123 // n @f$ N@f$ - Number of fitted peaks
124 // chi2 @f$ \chi^2 @f$
125 // nu @f$ \nu @f$ - number degrees of freedom
126 // c @f$ C@f$ - scale constant
127 // ec @f$ \delta C@f$ - error on @f$ C@f$
128 // delta @f$ \Delta@f$ - Most probable value
129 // edelta @f$ \delta\Delta@f$ - error on @f$\Delta@f$
130 // xi @f$ \xi@f$ - width
131 // exi @f$ \delta\xi@f$ - error on @f$\xi@f$
132 // sigma @f$ \sigma@f$ - Width of Gaussian
133 // esigma @f$ \delta\sigma@f$ - error on @f$\sigma@f$
134 // sigman @f$ \sigma_n@f$ - Noise width
135 // esigman @f$ \delta\sigma_n@f$ - error on @f$\sigma_n@f$
136 // a Array of @f$ N-1@f$ weights @f$ a_i@f$ for
138 // ea Array of @f$ N-1@f$ error on the weights @f$ a_i@f$ for
142 fA = new Double_t[fN];
143 fEA = new Double_t[fN];
144 for (Int_t i = 0; i < fN-1; i++) {
151 //____________________________________________________________________
152 AliFMDCorrELossFit::ELossFit::ELossFit(const ELossFit& o)
167 fESigmaN(o.fESigmaN),
169 fQuality(o.fQuality),
178 // o Object to copy from
181 fA = new Double_t[fN];
182 fEA = new Double_t[fN];
183 for (Int_t i = 0; i < fN-1; i++) {
191 //____________________________________________________________________
192 AliFMDCorrELossFit::ELossFit&
193 AliFMDCorrELossFit::ELossFit::operator=(const ELossFit& o)
196 // Assignment operator
199 // o Object to assign from
202 // Reference to this object
216 fESigmaN = o.fESigmaN;
217 fQuality = o.fQuality;
221 if (fA) delete [] fA;
222 if (fEA) delete [] fEA;
226 if (fN <= 0) return *this;
227 fA = new Double_t[fN];
228 fEA = new Double_t[fN];
229 for (Int_t i = 0; i < fN; i++) {
237 //____________________________________________________________________
238 AliFMDCorrELossFit::ELossFit::~ELossFit()
241 if (fEA) delete[] fEA;
245 //____________________________________________________________________
247 AliFMDCorrELossFit::ELossFit::FindMaxWeight(Double_t maxRelError,
248 Double_t leastWeight,
252 // Find the maximum weight to use. The maximum weight is the
253 // largest i for which
255 // - @f$ i \leq \max{N}@f$
256 // - @f$ a_i > \min{a}@f$
257 // - @f$ \delta a_i/a_i > \delta_{max}@f$
260 // maxRelError @f$ \min{a}@f$
261 // leastWeight @f$ \delta_{max}@f$
262 // maxN @f$ \max{N}@f$
265 // The largest index @f$ i@f$ for which the above
266 // conditions hold. Will never return less than 1.
268 Int_t n = TMath::Min(maxN, UShort_t(fN-1));
270 // fN is one larger than we have data
271 for (Int_t i = 0; i < n-1; i++, m++) {
272 if (fA[i] < leastWeight) break;
273 if (fEA[i] / fA[i] > maxRelError) break;
278 //____________________________________________________________________
280 AliFMDCorrELossFit::ELossFit::Evaluate(Double_t x,
286 // f_N(x;\Delta,\xi,\sigma') =
287 // \sum_{i=1}^{n} a_i f(x;\Delta_i,\xi_i,\sigma_i')
290 // (see AliForwardUtil::NLandauGaus) for the maximum @f$ N @f$
291 // that fulfills the requirements
294 // x Where to evaluate
295 // maxN @f$ \max{N}@f$
298 // @f$ f_N(x;\Delta,\xi,\sigma')@f$
300 return AliForwardUtil::NLandauGaus(x, fDelta, fXi, fSigma, fSigmaN,
301 TMath::Min(maxN, UShort_t(fN)), fA);
304 //____________________________________________________________________
306 AliFMDCorrELossFit::ELossFit::EvaluateWeighted(Double_t x,
312 // f_W(x;\Delta,\xi,\sigma') =
313 // \frac{\sum_{i=1}^{n} i a_i f_i(x;\Delta,\xi,\sigma')}{
314 // f_N(x;\Delta,\xi,\sigma')} =
315 // \frac{\sum_{i=1}^{n} i a_i f(x;\Delta_i,\xi_i,\sigma_i')}{
316 // \sum_{i=1}^{n} a_i f(x;\Delta_i,\xi_i,\sigma_i')}
318 // where @f$ n@f$ fulfills the requirements (see FindMaxWeight).
320 // If the denominator is zero, then 1 is returned.
322 // See also AliForwardUtil::ILandauGaus and AliForwardUtil::NLandauGaus
323 // for more information on the evaluated functions.
326 // x Where to evaluate
327 // maxN @f$ \max{N}@f$
330 // @f$ f_W(x;\Delta,\xi,\sigma')@f$.
332 UShort_t n = TMath::Min(maxN, UShort_t(fN-1));
335 for (Int_t i = 1; i <= n; i++) {
336 Double_t a = (i == 1 ? 1 : fA[i-1]);
337 if (fA[i-1] < 0) break;
338 Double_t f = AliForwardUtil::ILandauGaus(x,fDelta,fXi,fSigma,fSigmaN,i);
342 if (den <= 0) return 1;
347 #define OUTPAR(N,V,E) \
348 std::setprecision(9) << \
349 std::setw(12) << N << ": " << \
350 std::setw(14) << V << " +/- " << \
351 std::setw(14) << E << " (" << \
352 std::setprecision(-1) << \
353 std::setw(5) << 100*(V>0?E/V:1) << "%)\n"
356 //____________________________________________________________________
358 AliFMDCorrELossFit::ELossFit::Compare(const TObject* o) const
361 // Compare to another ELossFit object.
363 // - +1, if this quality is better (larger) than other objects quality
364 // - -1, if this quality is worse (smaller) than other objects quality
365 // - +1, if this @f$|\chi^2/\nu-1|@f$ is smaller than the same for other
366 // - -1, if this @f$|\chi^2/\nu-1|@f$ is larger than the same for other
370 // o Other object to compare to
372 const ELossFit* other = static_cast<const ELossFit*>(o);
373 if (this->fQuality < other->fQuality) return -1;
374 if (this->fQuality > other->fQuality) return +1;
375 Double_t chi2nu = (fNu == 0 ? 1000 : fChi2 / fNu);
376 Double_t oChi2nu = (other->fNu == 0 ? 1000 : other->fChi2 / other->fNu);
377 if (TMath::Abs(chi2nu-1) < TMath::Abs(oChi2nu-1)) return +1;
378 if (TMath::Abs(chi2nu-1) > TMath::Abs(oChi2nu-1)) return -1;
382 //____________________________________________________________________
384 AliFMDCorrELossFit::ELossFit::Print(Option_t*) const
387 // Information to standard output
392 std::cout << GetName() << ":\n"
393 << " chi^2/nu = " << fChi2 << "/" << fNu << " = "
394 << (fNu == 0 ? 999 : fChi2 / fNu) << "\n"
395 << " Quality: " << fQuality << "\n"
396 << " NParticles: " << fN << " (" << FindMaxWeight() << ")\n"
397 << OUTPAR("Delta", fDelta, fEDelta)
398 << OUTPAR("xi", fXi, fEXi)
399 << OUTPAR("sigma", fSigma, fESigma)
400 << OUTPAR("sigma_n", fSigmaN, fESigmaN);
401 for (Int_t i = 0; i < fN-1; i++)
402 std::cout << OUTPAR(Form("a%d", i+2), fA[i], fEA[i]);
403 std::cout << std::flush;
406 //____________________________________________________________________
408 AliFMDCorrELossFit::ELossFit::GetName() const
411 // Get the name of this object
417 return Form("FMD%d%c_etabin%03d", fDet, fRing, fBin);
420 //____________________________________________________________________
422 AliFMDCorrELossFit::ELossFit::Browse(TBrowser* b)
425 // Browse this object
430 Draw(b ? b->GetDrawOption() : "comp");
435 //____________________________________________________________________
437 AliFMDCorrELossFit::ELossFit::Draw(Option_t* option)
444 // - COMP Draw components too
449 if (opt.Contains("COMP")) {
450 opt.ReplaceAll("COMP","");
453 if (!opt.Contains("SAME")) {
458 TF1* tot = AliForwardUtil::MakeNLandauGaus(1,
463 tot->SetLineColor(kBlack);
464 tot->SetLineWidth(2);
465 tot->SetLineStyle(1);
466 tot->SetTitle(GetName());
467 Double_t max = tot->GetMaximum();
469 if (!opt.Contains("SAME")) {
470 TH1* frame = new TH1F(GetName(),
471 Form("FMD%d%c, eta bin %d",fDet,fRing,fBin),
473 frame->SetMinimum(max/10000);
474 frame->SetMaximum(max*1.4);
475 frame->SetXTitle("#Delta / #Delta_{mip}");
479 tot->DrawCopy(opt.Data());
489 Int_t maxW = FindMaxWeight();
490 for (Int_t i=1; i <= fN; i++) {
491 TF1* f = AliForwardUtil::MakeILandauGaus((i == 1 ? 1 : fA[i-2]),
496 f->SetLineStyle(i > maxW ? 2 : 1);
497 min = TMath::Min(f->GetMaximum(), min);
498 f->DrawCopy(opt.Data());
502 tot->GetHistogram()->SetMaximum(max);
503 tot->GetHistogram()->SetMinimum(min);
504 tot->GetHistogram()->GetYaxis()->SetRangeUser(min, max);
510 //____________________________________________________________________
511 #define CHECKPAR(V,E,T) ((V > 0) && (E / V < T))
513 //____________________________________________________________________
515 AliFMDCorrELossFit::ELossFit::CalculateQuality(Double_t maxChi2nu,
516 Double_t maxRelError,
517 Double_t leastWeight)
520 // Calculate the quality
523 if (fNu > 0 && fChi2 / fNu < maxChi2nu) qual += 4;;
524 if (CHECKPAR(fDelta, fEDelta, maxRelError)) qual++;
525 if (CHECKPAR(fXi, fEXi, maxRelError)) qual++;
526 if (CHECKPAR(fSigma, fESigma, maxRelError)) qual++;
527 if (CHECKPAR(fSigmaN, fESigmaN, maxRelError)) qual++;
528 qual += FindMaxWeight(1.5*maxRelError, leastWeight, fN);
532 //____________________________________________________________________
533 AliFMDCorrELossFit::AliFMDCorrELossFit()
540 // Default constructor
542 fRings.SetOwner(kTRUE);
543 fEtaAxis.SetTitle("#eta");
544 fEtaAxis.SetName("etaAxis");
545 fRings.SetName("rings");
548 //____________________________________________________________________
549 AliFMDCorrELossFit::AliFMDCorrELossFit(const AliFMDCorrELossFit& o)
552 fEtaAxis(o.fEtaAxis.GetNbins(),o.fEtaAxis.GetXmin(),o.fEtaAxis.GetXmax()),
559 // o Object to copy from
561 fEtaAxis.SetTitle("#eta");
562 fEtaAxis.SetName("etaAxis");
565 //____________________________________________________________________
566 AliFMDCorrELossFit::~AliFMDCorrELossFit()
574 //____________________________________________________________________
576 AliFMDCorrELossFit::operator=(const AliFMDCorrELossFit& o)
579 // Assignment operator
582 // o Object to assign from
585 // Reference to this object
589 SetEtaAxis(o.fEtaAxis.GetNbins(), o.fEtaAxis.GetXmin(), o.fEtaAxis.GetXmax());
593 //____________________________________________________________________
595 AliFMDCorrELossFit::FindEtaBin(Double_t eta) const
598 // Find the eta bin corresponding to the given eta
604 // Bin (in @f$[1,N_{bins}]@f$) corresponding to the given
605 // eta, or 0 if out of range.
607 if (fEtaAxis.GetXmin() == fEtaAxis.GetXmax() || fEtaAxis.GetNbins() == 0) {
608 AliWarning("No eta axis defined");
611 Int_t bin = const_cast<TAxis&>(fEtaAxis).FindBin(eta);
612 if (bin <= 0 || bin > fEtaAxis.GetNbins()) return 0;
616 //____________________________________________________________________
618 AliFMDCorrELossFit::SetFit(UShort_t d, Char_t r, Int_t etaBin, ELossFit* fit)
621 // Set the fit parameters from a function
626 // etaBin Eta (bin number, 1->nBins)
627 // f ELoss fit result - note, the object will take ownership
629 TObjArray* ringArray = GetOrMakeRingArray(d, r);
631 AliError(Form("Failed to make ring array for FMD%d%c", d, r));
634 if (etaBin <= 0 || etaBin >= fEtaAxis.GetNbins()+1) {
635 AliError(Form("bin=%d is out of range [%d,%d]",
636 etaBin, 1, fEtaAxis.GetNbins()));
639 // AliInfo(Form("Adding fit %p at %3d", fit, etaBin));
640 ringArray->AddAtAndExpand(fit, etaBin);
644 //____________________________________________________________________
646 AliFMDCorrELossFit::SetFit(UShort_t d, Char_t r, Double_t eta, ELossFit* fit)
649 // Set the fit parameters from a function
655 // f ELoss fit result - note, the object will take ownership
657 Int_t bin = FindEtaBin(eta);
659 AliError(Form("eta=%f is out of range [%f,%f]",
660 eta, fEtaAxis.GetXmin(), fEtaAxis.GetXmax()));
664 return SetFit(d, r, bin, fit);
666 //____________________________________________________________________
668 AliFMDCorrELossFit::SetFit(UShort_t d, Char_t r,
670 Int_t quality,UShort_t n,
671 Double_t chi2, UShort_t nu,
672 Double_t c, Double_t ec,
673 Double_t delta, Double_t edelta,
674 Double_t xi, Double_t exi,
675 Double_t sigma, Double_t esigma,
676 Double_t sigman, Double_t esigman,
677 Double_t* a, Double_t* ea)
680 // Set the fit parameters from a function
686 // quality Quality flag
687 // n @f$ N@f$ - Number of fitted peaks
688 // chi2 @f$ \chi^2 @f$
689 // nu @f$ \nu @f$ - number degrees of freedom
690 // c @f$ C@f$ - scale constant
691 // ec @f$ \delta C@f$ - error on @f$ C@f$
692 // delta @f$ \Delta@f$ - most probable value
693 // edelta @f$ \delta\Delta@f$ - error on @f$\Delta@f$
694 // xi @f$ \xi@f$ - Landau width
695 // exi @f$ \delta\xi@f$ - error on @f$\xi@f$
696 // sigma @f$ \sigma@f$ - Gaussian width
697 // esigma @f$ \delta\sigma@f$ - error on @f$\sigma@f$
698 // sigman @f$ \sigma_n@f$ - Noise width
699 // esigman @f$ \delta\sigma_n@f$ - error on @f$\sigma_n@f$
700 // a Array of @f$ N-1@f$ weights @f$ a_i@f$ for
702 // ea Array of @f$ N-1@f$ errors on weights @f$ a_i@f$ for
705 ELossFit* e = new ELossFit(quality, n,
713 if (!SetFit(d, r, eta, e)) {
719 //____________________________________________________________________
721 AliFMDCorrELossFit::SetFit(UShort_t d, Char_t r, Double_t eta,
722 Int_t quality, const TF1& f)
725 // Set the fit parameters from a function
731 // quality Quality flag
732 // f Function from fit
734 ELossFit* e = new ELossFit(quality, f);
735 if (!SetFit(d, r, eta, e)) {
741 //____________________________________________________________________
742 AliFMDCorrELossFit::ELossFit*
743 AliFMDCorrELossFit::FindFit(UShort_t d, Char_t r, Int_t etabin) const
746 // Find the fit corresponding to the specified parameters
751 // etabin Eta bin (1 based)
754 // Fit parameters or null in case of problems
756 TObjArray* ringArray = GetRingArray(d, r);
758 AliError(Form("Failed to make ring array for FMD%d%c", d, r));
761 if (etabin <= 0 || etabin >= fEtaAxis.GetNbins()) {
762 AliError(Form("Eta bin=%3d out of bounds [%d,%d] for FMD%d%c",
763 etabin, 1, fEtaAxis.GetNbins(), d, r));
766 if (etabin > ringArray->GetEntriesFast()) {
767 AliError(Form("Eta bin=%3d out of bounds [%d,%d] for FMD%d%c",
768 etabin, 1, ringArray->GetEntriesFast(), d, r));
771 else if (etabin >= ringArray->GetEntriesFast()) {
772 AliWarning(Form("Eta bin=%3d out of bounds by +1 [%d,%d] for FMD%d%c, "
773 "trying %3d", etabin, 1, ringArray->GetEntriesFast(), d, r,
777 else if (!ringArray->At(etabin)) {
778 AliWarning(Form("Eta bin=%d has no fit for FMD%d%c, trying %03d",
779 etabin, d, r, etabin+1));
782 return static_cast<ELossFit*>(ringArray->At(etabin));
784 //____________________________________________________________________
785 AliFMDCorrELossFit::ELossFit*
786 AliFMDCorrELossFit::FindFit(UShort_t d, Char_t r, Double_t eta) const
789 // Find the fit corresponding to the specified parameters
797 // Fit parameters or null in case of problems
799 Int_t etabin = FindEtaBin(eta);
800 return FindFit(d, r, etabin);
802 //____________________________________________________________________
804 AliFMDCorrELossFit::GetRingArray(UShort_t d, Char_t r) const
807 // Get the ring array corresponding to the specified ring
814 // Pointer to ring array, or null in case of problems
818 case 1: idx = 0; break;
819 case 2: idx = (r == 'i' || r == 'I') ? 1 : 2; break;
820 case 3: idx = (r == 'o' || r == 'I') ? 3 : 4; break;
822 if (idx < 0 || idx >= fRings.GetEntriesFast()) return 0;
823 return static_cast<TObjArray*>(fRings.At(idx));
825 //____________________________________________________________________
827 AliFMDCorrELossFit::GetOrMakeRingArray(UShort_t d, Char_t r)
830 // Get the ring array corresponding to the specified ring
837 // Pointer to ring array, or newly created container
841 case 1: idx = 0; break;
842 case 2: idx = (r == 'i' || r == 'I') ? 1 : 2; break;
843 case 3: idx = (r == 'o' || r == 'I') ? 3 : 4; break;
845 if (idx < 0) return 0;
846 if (idx >= fRings.GetEntriesFast() || !fRings.At(idx)) {
847 TObjArray* a = new TObjArray(0);
848 // TOrdCollection* a = new TOrdCollection(fEtaAxis.GetNbins());
849 a->SetName(Form("FMD%d%c", d, r));
851 fRings.AddAtAndExpand(a, idx);
853 return static_cast<TObjArray*>(fRings.At(idx));
857 TH1D* MakeHist(const TAxis& axis, const char* name, const char* title,
860 TH1D* h = new TH1D(Form("%s_%s", name, title),
861 Form("%s %s", name, title),
862 axis.GetNbins(), axis.GetXmin(), axis.GetXmax());
864 h->SetMarkerStyle(20);
865 h->SetMarkerColor(color);
866 h->SetMarkerSize(0.5);
867 h->SetFillColor(color);
868 h->SetFillStyle(3001);
869 h->SetLineColor(color);
874 //____________________________________________________________________
876 AliFMDCorrELossFit::Draw(Option_t* option)
882 // option Options. Possible values are
883 // - err Plot error bars
885 TString opt(Form("nostack %s", option));
887 Bool_t rel = (opt.Contains("rel"));
888 Bool_t err = (opt.Contains("err"));
889 if (rel) opt.ReplaceAll("rel","");
890 if (err) opt.ReplaceAll("err","");
891 Int_t nRings = fRings.GetEntriesFast();
893 for (Int_t i = 0; i < nRings; i++) {
894 if (!fRings.At(i)) continue;
895 TObjArray* a = static_cast<TObjArray*>(fRings.At(i));
896 Int_t nFits = a->GetEntriesFast();
898 for (Int_t j = 0; j < nFits; j++) {
899 ELossFit* fit = static_cast<ELossFit*>(a->At(j));
901 maxN = TMath::Max(maxN, UShort_t(fit->fN));
904 // AliInfo(Form("Maximum N is %d", maxN));
905 Int_t nPad = 7+maxN-1; // 7 regular params, and maxN-1 weights
906 TVirtualPad* pad = gPad;
907 pad->Divide(2, (nPad+1)/2, 0.1, 0, 0);
917 stacks.AddAt(chi2nu= new THStack("chi2", "#chi^{2}/#nu"), 0);
918 stacks.AddAt(c = new THStack("c", "C"), 1);
919 stacks.AddAt(delta = new THStack("delta", "#Delta_{mp}"), 2);
920 stacks.AddAt(xi = new THStack("xi", "#xi"), 3);
921 stacks.AddAt(sigma = new THStack("sigma", "#sigma"), 4);
922 stacks.AddAt(sigman= new THStack("sigman", "#sigma_{n}"), 5);
923 stacks.AddAt(n = new THStack("n", "N"), 6);
924 for (Int_t i = 1; i <= maxN; i++) {
925 stacks.AddAt(new THStack(Form("a_%02d", i+1), Form("a_{%d}", i+1)), 6+i);
928 for (Int_t i = 0; i < nRings; i++) {
929 if (!fRings.At(i)) continue;
930 TObjArray* a = static_cast<TObjArray*>(fRings.At(i));
931 Int_t nFits = a->GetEntriesFast();
934 case 0: color = kRed+2; break;
935 case 1: color = kGreen+2; break;
936 case 2: color = kGreen-2; break;
937 case 3: color = kBlue+2; break;
938 case 4: color = kBlue-2; break;
941 TH1D* hChi = MakeHist(fEtaAxis,a->GetName(), "chi2", color);
942 TH1D* hC = MakeHist(fEtaAxis,a->GetName(), "c", color);
943 TH1D* hDelta = MakeHist(fEtaAxis,a->GetName(), "delta", color);
944 TH1D* hXi = MakeHist(fEtaAxis,a->GetName(), "xi", color);
945 TH1D* hSigma = MakeHist(fEtaAxis,a->GetName(), "sigma", color);
946 TH1D* hSigmaN = MakeHist(fEtaAxis,a->GetName(), "sigman", color);
947 TH1D* hN = MakeHist(fEtaAxis,a->GetName(), "n", color);
949 const char* ho = (rel || !err ? "hist" : "e");
950 chi2nu->Add(hChi, "hist"); // 0
951 c ->Add(hC, ho); // 1
952 delta ->Add(hDelta, ho); // 2
953 xi ->Add(hXi, ho); // 3
954 sigma ->Add(hSigma, ho); // 4
955 sigman->Add(hSigmaN,ho); // 5
956 n ->Add(hN, "hist"); // 6
957 hChi->SetFillColor(color);
958 hChi->SetFillStyle(3001);
959 hN->SetFillColor(color);
960 hN->SetFillStyle(3001);
962 for (Int_t k = 1; k <= maxN; k++) {
963 hA[k-1] = MakeHist(fEtaAxis,a->GetName(), Form("a%02d",k+1), color);
964 static_cast<THStack*>(stacks.At(6+k))->Add(hA[k-1], ho);
967 for (Int_t j = 0; j < nFits; j++) {
968 ELossFit* f = static_cast<ELossFit*>(a->At(j));
972 Int_t nW = f->FindMaxWeight();
973 hChi ->SetBinContent(b, (f->fNu <= 0 ? 0 : f->fChi2 / f->fNu));
974 hN ->SetBinContent(b, nW);
977 hC ->SetBinContent(b, f->fC >0 ?f->fEC /f->fC :0);
978 hDelta ->SetBinContent(b, f->fDelta >0 ?f->fEDelta /f->fDelta :0);
979 hXi ->SetBinContent(b, f->fXi >0 ?f->fEXi /f->fXi :0);
980 hSigma ->SetBinContent(b, f->fSigma >0 ?f->fESigma /f->fSigma :0);
981 hSigmaN->SetBinContent(b, f->fSigmaN>0 ?f->fESigmaN/f->fSigmaN :0);
984 hC ->SetBinContent(b, f->fC);
985 hDelta ->SetBinContent(b, f->fDelta);
986 hXi ->SetBinContent(b, f->fXi);
987 hSigma ->SetBinContent(b, f->fSigma);
988 hSigmaN->SetBinContent(b, f->fSigmaN);
989 hC ->SetBinError(b, f->fEC);
990 hDelta ->SetBinError(b, f->fEDelta);
991 hXi ->SetBinError(b, f->fEXi);
992 hSigma ->SetBinError(b, f->fESigma);
993 hSigmaN->SetBinError(b, f->fESigmaN);
995 for (Int_t k = 0; k < f->fN-1 && k < maxN; k++) {
997 hA[k]->SetBinContent(b, f->fA[k] > 0 ? f->fEA[k] / f->fA[k] : 0);
999 hA[k]->SetBinContent(b, f->fA[k]);
1000 hA[k]->SetBinError(b, f->fEA[k]);
1005 Int_t nPad2 = (nPad+1) / 2;
1006 for (Int_t i = 0; i < nPad; i++) {
1007 Int_t iPad = 1 + i/nPad2 + 2 * (i % nPad2);
1008 TVirtualPad* p = pad->cd(iPad);
1009 p->SetLeftMargin(.15);
1014 if (rel && i != 0 && i != 6 && i != 5 && i != 4) p->SetLogy();
1016 THStack* stack = static_cast<THStack*>(stacks.At(i));
1017 // AliInfo(Form("Drawing %s (%d) in pad %d", stack->GetName(), i, iPad));
1018 stack->Draw(opt.Data());
1020 TString tit(stack->GetTitle());
1021 if (rel && i != 0 && i != 6)
1022 tit = Form("#delta %s/%s", tit.Data(), tit.Data());
1023 TH1* hist = stack->GetHistogram();
1024 TAxis* yaxis = hist->GetYaxis();
1025 yaxis->SetTitle(tit.Data());
1026 yaxis->SetTitleSize(0.15);
1027 yaxis->SetLabelSize(0.08);
1028 yaxis->SetTitleOffset(0.35);
1029 yaxis->SetNdivisions(5);
1031 TAxis* xaxis = stack->GetHistogram()->GetXaxis();
1032 xaxis->SetTitle("#eta");
1033 xaxis->SetTitleSize(0.15);
1034 xaxis->SetLabelSize(0.08);
1035 xaxis->SetTitleOffset(0.35);
1036 xaxis->SetNdivisions(10);
1041 case 0: break; // chi^2/nu
1043 case 2: stack->SetMinimum(0.4); break; // Delta
1044 case 3: stack->SetMinimum(0.02);break; // xi
1045 case 4: stack->SetMinimum(0.05);break; // sigma
1046 case 5: break; // sigmaN
1048 stack->SetMinimum(-.1);
1049 stack->SetMaximum(maxN+.5);
1051 default: break; // Weights
1054 stack->DrawClone(opt.Data());
1059 //____________________________________________________________________
1061 AliFMDCorrELossFit::Print(Option_t* option) const
1063 TString opt(option);
1065 Int_t nRings = fRings.GetEntriesFast();
1066 bool recurse = opt.Contains("R");
1068 std::cout << "Low cut in fit range: " << fLowCut << "\n"
1069 << "Eta axis: " << fEtaAxis.GetNbins()
1070 << " bins, range [" << fEtaAxis.GetXmin() << ","
1071 << fEtaAxis.GetXmax() << "]" << std::endl;
1073 for (Int_t i = 0; i < nRings; i++) {
1074 if (!fRings.At(i)) continue;
1075 TObjArray* a = static_cast<TObjArray*>(fRings.At(i));
1076 Int_t nFits = a->GetEntriesFast();
1078 std::cout << a->GetName() << " [" << nFits << " entries]"
1079 << (recurse ? ":\n" : "\t");
1080 Int_t min = fEtaAxis.GetNbins()+1;
1082 for (Int_t j = 0; j < nFits; j++) {
1083 if (!a->At(j)) continue;
1085 min = TMath::Min(j, min);
1086 max = TMath::Max(j, max);
1089 std::cout << "Bin # " << j << "\t";
1090 ELossFit* fit = static_cast<ELossFit*>(a->At(j));
1095 std::cout << " bin range: " << std::setw(3) << min
1096 << "-" << std::setw(3) << max << " " << std::setw(3)
1097 << (max-min+1) << " bins" << std::endl;
1101 //____________________________________________________________________
1103 AliFMDCorrELossFit::Browse(TBrowser* b)
1106 // Browse this object
1117 //____________________________________________________________________