This macro fits the source spectrum using the AWMI algorithm from the "TSpectrumFit" class ("TSpectrum" class is used to find peaks).
To try this macro, in a ROOT (5 or 6) prompt, do:
root > .x FitAwmi.C
or:
root > .x FitAwmi.C++
root > FitAwmi(); // re-run with another random set of peaks
Author:
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, March 19, 2024 at 07:19 PM.
Definition of a helper function:
%%cpp -d
#include "TROOT.h"
#include "TMath.h"
#include "TRandom.h"
#include "TH1.h"
#include "TF1.h"
#include "TCanvas.h"
#include "TSpectrum.h"
#include "TSpectrumFit.h"
#include "TPolyMarker.h"
#include "TList.h"
#include <iostream>
TH1F *FitAwmi_Create_Spectrum(void) {
Int_t nbins = 1000;
Double_t xmin = -10., xmax = 10.;
delete gROOT->FindObject("h"); // prevent "memory leak"
TH1F *h = new TH1F("h", "simulated spectrum", nbins, xmin, xmax);
h->SetStats(kFALSE);
TF1 f("f", "TMath::Gaus(x, [0], [1], 1)", xmin, xmax);
// f.SetParNames("mean", "sigma");
gRandom->SetSeed(0); // make it really random
// create well separated peaks with exactly known means and areas
// note: TSpectrumFit assumes that all peaks have the same sigma
Double_t sigma = (xmax - xmin) / Double_t(nbins) * Int_t(gRandom->Uniform(2., 6.));
Int_t npeaks = 0;
while (xmax > (xmin + 6. * sigma)) {
npeaks++;
xmin += 3. * sigma; // "mean"
f.SetParameters(xmin, sigma);
Double_t area = 1. * Int_t(gRandom->Uniform(1., 11.));
h->Add(&f, area, ""); // "" ... or ... "I"
std::cout << "created "
<< xmin << " "
<< (area / sigma / TMath::Sqrt(TMath::TwoPi())) << " "
<< area << std::endl;
xmin += 3. * sigma;
}
std::cout << "the total number of created peaks = " << npeaks
<< " with sigma = " << sigma << std::endl;
return h;
}
Arguments are defined.
TH1F *h = FitAwmi_Create_Spectrum();
TCanvas *cFit = ((TCanvas *)(gROOT->GetListOfCanvases()->FindObject("cFit")));
if (!cFit) cFit = new TCanvas("cFit", "cFit", 10, 10, 1000, 700);
else cFit->Clear();
h->Draw("L");
Int_t i, nfound, bin;
Int_t nbins = h->GetNbinsX();
Double_t *source = new Double_t[nbins];
Double_t *dest = new Double_t[nbins];
for (i = 0; i < nbins; i++) source[i] = h->GetBinContent(i + 1);
TSpectrum *s = new TSpectrum(); // note: default maxpositions = 100
created -9.7 7.97885 2 created -9.1 23.9365 6 created -8.5 11.9683 3 created -7.9 11.9683 3 created -7.3 39.8942 10 created -6.7 31.9154 8 created -6.1 15.9577 4 created -5.5 23.9365 6 created -4.9 15.9577 4 created -4.3 31.9154 8 created -3.7 27.926 7 created -3.1 3.98942 1 created -2.5 11.9683 3 created -1.9 19.9471 5 created -1.3 7.97885 2 created -0.7 39.8942 10 created -0.1 3.98942 1 created 0.5 39.8942 10 created 1.1 15.9577 4 created 1.7 23.9365 6 created 2.3 31.9154 8 created 2.9 7.97885 2 created 3.5 7.97885 2 created 4.1 3.98942 1 created 4.7 31.9154 8 created 5.3 19.9471 5 created 5.9 23.9365 6 created 6.5 31.9154 8 created 7.1 11.9683 3 created 7.7 23.9365 6 created 8.3 23.9365 6 created 8.9 11.9683 3 created 9.5 35.9048 9 the total number of created peaks = 33 with sigma = 0.1
searching for candidate peaks positions
nfound = s->SearchHighRes(source, dest, nbins, 2., 2., kFALSE, 10000, kFALSE, 0);
filling in the initial estimates of the input parameters
Bool_t *FixPos = new Bool_t[nfound];
Bool_t *FixAmp = new Bool_t[nfound];
for(i = 0; i < nfound; i++) FixAmp[i] = FixPos[i] = kFALSE;
Double_t *Pos, *Amp = new Double_t[nfound]; // ROOT 6
Pos = s->GetPositionX(); // 0 ... (nbins - 1)
for (i = 0; i < nfound; i++) {
bin = 1 + Int_t(Pos[i] + 0.5); // the "nearest" bin
Amp[i] = h->GetBinContent(bin);
}
TSpectrumFit *pfit = new TSpectrumFit(nfound);
pfit->SetFitParameters(0, (nbins - 1), 1000, 0.1, pfit->kFitOptimChiCounts,
pfit->kFitAlphaHalving, pfit->kFitPower2,
pfit->kFitTaylorOrderFirst);
pfit->SetPeakParameters(2., kFALSE, Pos, FixPos, Amp, FixAmp);
pfit->SetBackgroundParameters(source[0], kFALSE, 0., kFALSE, 0., kFALSE);
pfit->FitAwmi(source);
Double_t *Positions = pfit->GetPositions();
Double_t *PositionsErrors = pfit->GetPositionsErrors();
Double_t *Amplitudes = pfit->GetAmplitudes();
Double_t *AmplitudesErrors = pfit->GetAmplitudesErrors();
Double_t *Areas = pfit->GetAreas();
Double_t *AreasErrors = pfit->GetAreasErrors();
delete gROOT->FindObject("d"); // prevent "memory leak"
TH1F *d = new TH1F(*h); d->SetNameTitle("d", ""); d->Reset("M");
for (i = 0; i < nbins; i++) d->SetBinContent(i + 1, source[i]);
Double_t x1 = d->GetBinCenter(1), dx = d->GetBinWidth(1);
Double_t sigma, sigmaErr;
pfit->GetSigma(sigma, sigmaErr);
current TSpectrumFit needs a sqrt(2) correction factor for sigma
sigma /= TMath::Sqrt2(); sigmaErr /= TMath::Sqrt2();
convert "bin numbers" into "x-axis values"
sigma *= dx; sigmaErr *= dx;
std::cout << "the total number of found peaks = " << nfound
<< " with sigma = " << sigma << " (+-" << sigmaErr << ")"
<< std::endl;
std::cout << "fit chi^2 = " << pfit->GetChi() << std::endl;
for (i = 0; i < nfound; i++) {
bin = 1 + Int_t(Positions[i] + 0.5); // the "nearest" bin
Pos[i] = d->GetBinCenter(bin);
Amp[i] = d->GetBinContent(bin);
// convert "bin numbers" into "x-axis values"
Positions[i] = x1 + Positions[i] * dx;
PositionsErrors[i] *= dx;
Areas[i] *= dx;
AreasErrors[i] *= dx;
std::cout << "found "
<< Positions[i] << " (+-" << PositionsErrors[i] << ") "
<< Amplitudes[i] << " (+-" << AmplitudesErrors[i] << ") "
<< Areas[i] << " (+-" << AreasErrors[i] << ")"
<< std::endl;
}
d->SetLineColor(kRed); d->SetLineWidth(1);
d->Draw("SAME L");
TPolyMarker *pm = ((TPolyMarker*)(h->GetListOfFunctions()->FindObject("TPolyMarker")));
if (pm) {
h->GetListOfFunctions()->Remove(pm);
delete pm;
}
pm = new TPolyMarker(nfound, Pos, Amp);
h->GetListOfFunctions()->Add(pm);
pm->SetMarkerStyle(23);
pm->SetMarkerColor(kRed);
pm->SetMarkerSize(1);
the total number of found peaks = 33 with sigma = 0.100002 (+-4.8236e-05) fit chi^2 = 6.23939e-06 found -7.3 (+-0.00034623) 39.894 (+-0.136325) 10.0001 (+-0.00111875) found -0.7 (+-0.000344047) 39.8936 (+-0.136226) 10 (+-0.00111793) found 0.500001 (+-0.000344648) 39.8937 (+-0.136253) 10.0001 (+-0.00111816) found 9.5 (+-0.000361833) 35.9047 (+-0.129213) 9.00015 (+-0.00106039) found -6.7 (+-0.000388365) 31.9155 (+-0.121982) 8.00018 (+-0.00100104) found -4.3 (+-0.000387796) 31.9153 (+-0.121959) 8.00014 (+-0.00100085) found 2.3 (+-0.000386828) 31.9152 (+-0.121923) 8.0001 (+-0.00100056) found 4.7 (+-0.000386015) 31.9151 (+-0.121894) 8.00008 (+-0.00100032) found 6.5 (+-0.000387239) 31.9152 (+-0.121938) 8.00012 (+-0.00100068) found -3.7 (+-0.000413717) 27.9259 (+-0.114057) 7.00012 (+-0.000936011) found -9.1 (+-0.000446386) 23.9363 (+-0.105579) 6.00006 (+-0.000866435) found -5.5 (+-0.000447769) 23.9365 (+-0.105619) 6.0001 (+-0.000866758) found 1.7 (+-0.000449005) 23.9367 (+-0.105656) 6.00016 (+-0.000867063) found 5.9 (+-0.000449373) 23.9368 (+-0.105667) 6.00017 (+-0.000867152) found 7.7 (+-0.000448025) 23.9366 (+-0.105627) 6.00012 (+-0.000866824) found 8.3 (+-0.000448025) 23.9366 (+-0.105627) 6.00012 (+-0.000866824) found -1.9 (+-0.00048953) 19.947 (+-0.0963929) 5.00006 (+-0.000791047) found 5.3 (+-0.000493414) 19.9475 (+-0.0964892) 5.00018 (+-0.000791837) found -6.1 (+-0.000552827) 15.9581 (+-0.0863271) 4.00018 (+-0.000708442) found -4.9 (+-0.000552827) 15.9581 (+-0.0863271) 4.00018 (+-0.000708442) found 1.1 (+-0.00055349) 15.9582 (+-0.0863415) 4.00021 (+-0.00070856) found -8.5 (+-0.000637212) 11.9685 (+-0.0747441) 3.00012 (+-0.000613387) found 7.1 (+-0.000640261) 11.9688 (+-0.0747921) 3.00018 (+-0.000613781) found 8.9 (+-0.0006407) 11.9688 (+-0.0747994) 3.0002 (+-0.00061384) found -7.9 (+-0.000639044) 11.9687 (+-0.0747745) 3.00017 (+-0.000613636) found -2.5 (+-0.00063428) 11.9684 (+-0.0747013) 3.00008 (+-0.000613035) found -1.29999 (+-0.000788236) 7.97947 (+-0.0611132) 2.0002 (+-0.000501525) found 2.89999 (+-0.000783561) 7.97922 (+-0.0610634) 2.00013 (+-0.000501116) found -9.69999 (+-0.000780307) 7.97901 (+-0.0610246) 2.00008 (+-0.000500798) found 3.5 (+-0.000776001) 7.97886 (+-0.0609836) 2.00004 (+-0.000500461) found -0.1 (+-0.00113189) 3.99039 (+-0.0433143) 1.00026 (+-0.000355458) found -3.10001 (+-0.00111868) 3.98987 (+-0.0432365) 1.00013 (+-0.00035482) found 4.10001 (+-0.0011173) 3.98987 (+-0.04323) 1.00013 (+-0.000354767)
cleanup
delete pfit;
delete [] Amp;
delete [] FixAmp;
delete [] FixPos;
delete s;
delete [] dest;
delete [] source;
return;
Draw all canvases
gROOT->GetListOfCanvases()->Draw()