Logo ROOT  
Reference Guide
FitAwmi.C File Reference

Detailed Description

View in nbviewer Open in SWAN

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
created -9.7 31.9154 8
created -9.1 31.9154 8
created -8.5 7.97885 2
created -7.9 3.98942 1
created -7.3 23.9365 6
created -6.7 11.9683 3
created -6.1 11.9683 3
created -5.5 15.9577 4
created -4.9 11.9683 3
created -4.3 19.9471 5
created -3.7 27.926 7
created -3.1 3.98942 1
created -2.5 3.98942 1
created -1.9 3.98942 1
created -1.3 11.9683 3
created -0.7 23.9365 6
created -0.1 31.9154 8
created 0.5 27.926 7
created 1.1 31.9154 8
created 1.7 7.97885 2
created 2.3 39.8942 10
created 2.9 7.97885 2
created 3.5 39.8942 10
created 4.1 7.97885 2
created 4.7 27.926 7
created 5.3 3.98942 1
created 5.9 3.98942 1
created 6.5 27.926 7
created 7.1 35.9048 9
created 7.7 35.9048 9
created 8.3 39.8942 10
created 8.9 7.97885 2
created 9.5 3.98942 1
the total number of created peaks = 33 with sigma = 0.1
the total number of found peaks = 33 with sigma = 0.100002 (+-1.79137e-05)
fit chi^2 = 7.98115e-07
found 2.3 (+-0.000123213) 39.8937 (+-0.0487284) 10.0001 (+-0.000399889)
found 3.5 (+-0.000123213) 39.8937 (+-0.0487284) 10 (+-0.000399889)
found 8.3 (+-0.000123765) 39.894 (+-0.0487545) 10.0001 (+-0.000400103)
found 7.1 (+-0.000131032) 35.9049 (+-0.0462766) 9.00021 (+-0.000379768)
found 7.7 (+-0.000131217) 35.9051 (+-0.0462849) 9.00024 (+-0.000379837)
found -9.7 (+-0.000138788) 31.9152 (+-0.043619) 8.0001 (+-0.000357959)
found -9.1 (+-0.000138503) 31.9153 (+-0.0436122) 8.00013 (+-0.000357903)
found -0.0999997 (+-0.000138888) 31.9154 (+-0.0436263) 8.00017 (+-0.000358019)
found 1.1 (+-0.00013843) 31.9152 (+-0.0436092) 8.00012 (+-0.000357879)
found -3.7 (+-0.000147697) 27.9257 (+-0.0407834) 7.00008 (+-0.000334689)
found 0.5 (+-0.000148887) 27.9262 (+-0.040823) 7.00021 (+-0.000335014)
found 4.7 (+-0.000147289) 27.9256 (+-0.0407699) 7.00004 (+-0.000334578)
found 6.5 (+-0.000148042) 27.9259 (+-0.0407957) 7.00013 (+-0.00033479)
found -7.3 (+-0.00015939) 23.9363 (+-0.0377538) 6.00005 (+-0.000309827)
found -0.699998 (+-0.000160434) 23.9366 (+-0.0377837) 6.00014 (+-0.000310072)
found -4.3 (+-0.00017589) 19.9472 (+-0.0344952) 5.00013 (+-0.000283084)
found -5.5 (+-0.000196304) 15.9577 (+-0.0308462) 4.00008 (+-0.000253139)
found -6.7 (+-0.000227901) 11.9685 (+-0.0267324) 3.00012 (+-0.00021938)
found -6.1 (+-0.000227456) 11.9684 (+-0.0267254) 3.00009 (+-0.000219322)
found -4.9 (+-0.000227975) 11.9685 (+-0.0267334) 3.00012 (+-0.000219387)
found -1.3 (+-0.000227058) 11.9684 (+-0.0267204) 3.00009 (+-0.000219281)
found 1.7 (+-0.000282724) 7.97963 (+-0.0218663) 2.00023 (+-0.000179445)
found 2.9 (+-0.000283154) 7.97973 (+-0.0218712) 2.00026 (+-0.000179486)
found -8.50001 (+-0.000279535) 7.97916 (+-0.0218327) 2.00012 (+-0.00017917)
found 4.1 (+-0.000282482) 7.97957 (+-0.0218635) 2.00022 (+-0.000179423)
found 8.89999 (+-0.00027995) 7.97927 (+-0.0218376) 2.00014 (+-0.00017921)
found -3.10001 (+-0.000397888) 3.98976 (+-0.0154522) 1.0001 (+-0.000126808)
found 5.29999 (+-0.000397888) 3.98977 (+-0.0154522) 1.00011 (+-0.000126808)
found -7.89999 (+-0.000398723) 3.98976 (+-0.015456) 1.0001 (+-0.000126839)
found 5.90001 (+-0.000397888) 3.98977 (+-0.0154522) 1.00011 (+-0.000126808)
found 9.5 (+-0.000390856) 3.9895 (+-0.0154179) 1.00004 (+-0.000126527)
found -2.5 (+-0.00039348) 3.98945 (+-0.0154274) 1.00003 (+-0.000126605)
found -1.9 (+-0.00039562) 3.98955 (+-0.0154389) 1.00005 (+-0.000126699)
#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;
}
void FitAwmi(void) {
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
// 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->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);
// cleanup
delete pfit;
delete [] Amp;
delete [] FixAmp;
delete [] FixPos;
delete s;
delete [] dest;
delete [] source;
return;
}
Author

Definition in file FitAwmi.C.

TSpectrumFit::GetPositions
Double_t * GetPositions() const
Definition: TSpectrumFit.h:122
f
#define f(i)
Definition: RSha256.hxx:122
dest
#define dest(otri, vertexptr)
Definition: triangle.c:1040
TSpectrumFit::SetFitParameters
void SetFitParameters(Int_t xmin, Int_t xmax, Int_t numberIterations, Double_t alpha, Int_t statisticType, Int_t alphaOptim, Int_t power, Int_t fitTaylor)
This function sets the following fitting parameters:
Definition: TSpectrumFit.cxx:2608
xmax
float xmax
Definition: THbookFile.cxx:95
TMath::Sqrt
Double_t Sqrt(Double_t x)
Definition: TMath.h:691
TRandom.h
TGeant4Unit::s
static constexpr double s
Definition: TGeant4SystemOfUnits.h:168
Int_t
int Int_t
Definition: RtypesCore.h:45
TRandom::Uniform
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:635
TSpectrumFit::GetAreas
Double_t * GetAreas() const
Definition: TSpectrumFit.h:118
TList.h
TCanvas::Clear
void Clear(Option_t *option="")
Remove all primitives from the canvas.
Definition: TCanvas.cxx:723
TCanvas.h
Bool_t
bool Bool_t
Definition: RtypesCore.h:63
TSpectrum.h
x1
static const double x1[5]
Definition: RooGaussKronrodIntegrator1D.cxx:346
TROOT.h
TSpectrumFit::GetAreasErrors
Double_t * GetAreasErrors() const
Definition: TSpectrumFit.h:119
TSpectrumFit::GetAmplitudesErrors
Double_t * GetAmplitudesErrors() const
Definition: TSpectrumFit.h:117
TSpectrumFit::GetAmplitudes
Double_t * GetAmplitudes() const
Definition: TSpectrumFit.h:116
TSpectrumFit::kFitTaylorOrderFirst
@ kFitTaylorOrderFirst
Definition: TSpectrumFit.h:82
TSpectrumFit::GetSigma
void GetSigma(Double_t &sigma, Double_t &sigmaErr)
This function gets the sigma parameter and its error.
Definition: TSpectrumFit.cxx:2725
xmin
float xmin
Definition: THbookFile.cxx:95
h
#define h(i)
Definition: RSha256.hxx:124
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:92
TSpectrumFit::kFitAlphaHalving
@ kFitAlphaHalving
Definition: TSpectrumFit.h:74
gRandom
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
TRandom::SetSeed
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:597
TSpectrumFit::GetChi
Double_t GetChi() const
Definition: TSpectrumFit.h:121
kRed
@ kRed
Definition: Rtypes.h:66
TSpectrumFit::SetPeakParameters
void SetPeakParameters(Double_t sigma, Bool_t fixSigma, const Double_t *positionInit, const Bool_t *fixPosition, const Double_t *ampInit, const Bool_t *fixAmp)
This function sets the following fitting parameters of peaks:
Definition: TSpectrumFit.cxx:2656
TMath::TwoPi
constexpr Double_t TwoPi()
Definition: TMath.h:50
sigma
const Double_t sigma
Definition: h1analysisProxy.h:11
TSpectrumFit::FitAwmi
void FitAwmi(Double_t *source)
This function fits the source spectrum.
Definition: TSpectrumFit.cxx:820
Double_t
double Double_t
Definition: RtypesCore.h:59
TF1.h
TSpectrumFit::kFitOptimChiCounts
@ kFitOptimChiCounts
Definition: TSpectrumFit.h:71
TCanvas
Definition: TCanvas.h:23
TH1F
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:572
TPolyMarker
Definition: TPolyMarker.h:31
TMath::Sqrt2
constexpr Double_t Sqrt2()
Definition: TMath.h:94
d
#define d(i)
Definition: RSha256.hxx:120
TPolyMarker.h
TSpectrumFit::kFitPower2
@ kFitPower2
Definition: TSpectrumFit.h:76
TSpectrumFit
Advanced 1-dimensional spectra fitting functions.
Definition: TSpectrumFit.h:18
TSpectrumFit::GetPositionsErrors
Double_t * GetPositionsErrors() const
Definition: TSpectrumFit.h:123
TF1
1-Dim function class
Definition: TF1.h:212
TSpectrum
Advanced Spectra Processing.
Definition: TSpectrum.h:18
TH1.h
TMath.h
TSpectrumFit.h
gROOT
#define gROOT
Definition: TROOT.h:406