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 Wednesday, April 17, 2024 at 11:21 AM.
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 27.926 7 created -9.1 39.8942 10 created -8.5 35.9048 9 created -7.9 7.97885 2 created -7.3 27.926 7 created -6.7 35.9048 9 created -6.1 27.926 7 created -5.5 35.9048 9 created -4.9 19.9471 5 created -4.3 23.9365 6 created -3.7 3.98942 1 created -3.1 11.9683 3 created -2.5 31.9154 8 created -1.9 35.9048 9 created -1.3 39.8942 10 created -0.7 3.98942 1 created -0.1 35.9048 9 created 0.5 39.8942 10 created 1.1 39.8942 10 created 1.7 23.9365 6 created 2.3 39.8942 10 created 2.9 23.9365 6 created 3.5 31.9154 8 created 4.1 35.9048 9 created 4.7 39.8942 10 created 5.3 15.9577 4 created 5.9 35.9048 9 created 6.5 27.926 7 created 7.1 3.98942 1 created 7.7 7.97885 2 created 8.3 35.9048 9 created 8.9 19.9471 5 created 9.5 23.9365 6 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 (+-3.45197e-05) fit chi^2 = 4.20486e-06 found -9.1 (+-0.000285097) 39.8943 (+-0.111954) 10.0002 (+-0.000918748) found -1.3 (+-0.0002837) 39.894 (+-0.111891) 10.0001 (+-0.000918234) found 0.5 (+-0.000285485) 39.8944 (+-0.111973) 10.0002 (+-0.000918906) found 1.1 (+-0.000285063) 39.8943 (+-0.111952) 10.0002 (+-0.000918737) found 2.3 (+-0.000284525) 39.8941 (+-0.111926) 10.0002 (+-0.000918519) found 4.7 (+-0.000284582) 39.8941 (+-0.111929) 10.0002 (+-0.000918548) found -8.5 (+-0.000299772) 35.9047 (+-0.106179) 9.00015 (+-0.000871355) found -6.7 (+-0.000300467) 35.9048 (+-0.106206) 9.00018 (+-0.000871581) found -5.5 (+-0.000300114) 35.9047 (+-0.106191) 9.00015 (+-0.000871455) found -1.9 (+-0.000301044) 35.905 (+-0.106232) 9.00023 (+-0.000871794) found -0.0999979 (+-0.000299353) 35.9047 (+-0.106163) 9.00014 (+-0.000871225) found 4.1 (+-0.000301044) 35.905 (+-0.106232) 9.00023 (+-0.000871794) found 5.9 (+-0.000299903) 35.9047 (+-0.106182) 9.00014 (+-0.000871382) found 8.3 (+-0.000299004) 35.9045 (+-0.106144) 9.00009 (+-0.000871073) found -2.5 (+-0.000318401) 31.9154 (+-0.100122) 8.00015 (+-0.000821652) found 3.5 (+-0.000319117) 31.9155 (+-0.100149) 8.00019 (+-0.000821875) found -6.1 (+-0.000342094) 27.9264 (+-0.0937144) 7.00023 (+-0.000769066) found 6.5 (+-0.000339804) 27.9259 (+-0.0936391) 7.00013 (+-0.000768449) found -9.7 (+-0.000341082) 27.926 (+-0.0936721) 7.00013 (+-0.000768719) found -7.3 (+-0.000340331) 27.926 (+-0.093655) 7.00014 (+-0.000768579) found 1.7 (+-0.000370399) 23.9371 (+-0.0867906) 6.00026 (+-0.000712246) found 2.9 (+-0.000370011) 23.937 (+-0.0867785) 6.00023 (+-0.000712147) found -4.3 (+-0.000366495) 23.9364 (+-0.0866759) 6.00008 (+-0.000711305) found 9.5 (+-0.000364856) 23.9366 (+-0.0866404) 6.00014 (+-0.000711014) found -4.9 (+-0.000405291) 19.9475 (+-0.0792168) 5.00019 (+-0.000650092) found 8.9 (+-0.000405291) 19.9475 (+-0.0792168) 5.0002 (+-0.000650092) found 5.3 (+-0.000455309) 15.9584 (+-0.0708998) 4.00025 (+-0.000581838) found -3.1 (+-0.000521977) 11.9685 (+-0.0613453) 3.00012 (+-0.000503429) found -7.9 (+-0.000647915) 7.97953 (+-0.0501783) 2.00021 (+-0.000411788) found 7.70001 (+-0.000642116) 7.97922 (+-0.0501189) 2.00013 (+-0.0004113) found -0.700002 (+-0.000928359) 3.99034 (+-0.0355527) 1.00025 (+-0.000291763) found -3.70001 (+-0.000917288) 3.98982 (+-0.0354876) 1.00012 (+-0.000291228) found 7.09999 (+-0.000916259) 3.98982 (+-0.0354828) 1.00012 (+-0.000291189)
cleanup
delete pfit;
delete [] Amp;
delete [] FixAmp;
delete [] FixPos;
delete s;
delete [] dest;
delete [] source;
return;
Draw all canvases
gROOT->GetListOfCanvases()->Draw()