Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
FitAwmi.C File Reference

Detailed Description

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 35.9048 9
created -9.1 15.9577 4
created -8.5 15.9577 4
created -7.9 27.926 7
created -7.3 23.9365 6
created -6.7 39.8942 10
created -6.1 7.97885 2
created -5.5 7.97885 2
created -4.9 3.98942 1
created -4.3 31.9154 8
created -3.7 27.926 7
created -3.1 23.9365 6
created -2.5 39.8942 10
created -1.9 19.9471 5
created -1.3 35.9048 9
created -0.7 15.9577 4
created -0.1 39.8942 10
created 0.5 23.9365 6
created 1.1 3.98942 1
created 1.7 23.9365 6
created 2.3 11.9683 3
created 2.9 35.9048 9
created 3.5 7.97885 2
created 4.1 3.98942 1
created 4.7 35.9048 9
created 5.3 31.9154 8
created 5.9 19.9471 5
created 6.5 31.9154 8
created 7.1 31.9154 8
created 7.7 15.9577 4
created 8.3 31.9154 8
created 8.9 15.9577 4
created 9.5 31.9154 8
the total number of created peaks = 33 with sigma = 0.1
the total number of found peaks = 33 with sigma = 0.100002 (+-4.26266e-05)
fit chi^2 = 5.55031e-06
found -6.7 (+-0.000325903) 39.8939 (+-0.128547) 10.0001 (+-0.00105492)
found -2.5 (+-0.000326697) 39.894 (+-0.128583) 10.0001 (+-0.00105521)
found -0.0999996 (+-0.000326477) 39.894 (+-0.128572) 10.0001 (+-0.00105513)
found -9.7 (+-0.000344085) 35.9043 (+-0.121962) 9.00005 (+-0.00100088)
found -1.3 (+-0.000344156) 35.9046 (+-0.121975) 9.00012 (+-0.00100099)
found 2.9 (+-0.000343007) 35.9044 (+-0.121928) 9.00006 (+-0.0010006)
found 4.7 (+-0.00034362) 35.9046 (+-0.121957) 9.00011 (+-0.00100084)
found 5.3 (+-0.000366397) 31.9155 (+-0.115053) 8.00018 (+-0.000944178)
found -4.3 (+-0.000364519) 31.9152 (+-0.114984) 8.0001 (+-0.000943611)
found 6.5 (+-0.000366218) 31.9155 (+-0.115045) 8.00017 (+-0.000944118)
found 7.1 (+-0.000365949) 31.9154 (+-0.115035) 8.00015 (+-0.000944035)
found 8.3 (+-0.000365042) 31.9152 (+-0.115) 8.0001 (+-0.000943743)
found 9.5 (+-0.000362392) 31.9154 (+-0.114916) 8.00015 (+-0.000943053)
found -7.9 (+-0.000391147) 27.926 (+-0.107603) 7.00013 (+-0.000883041)
found -3.7 (+-0.000392171) 27.9262 (+-0.107638) 7.00018 (+-0.000883331)
found 0.499997 (+-0.000422315) 23.9367 (+-0.0996211) 6.00014 (+-0.000817539)
found -7.3 (+-0.000424856) 23.937 (+-0.0996923) 6.00022 (+-0.000818124)
found -3.1 (+-0.000424856) 23.937 (+-0.0996923) 6.00022 (+-0.000818124)
found 1.7 (+-0.000420327) 23.9363 (+-0.0995605) 6.00005 (+-0.000817042)
found -1.9 (+-0.000466774) 19.9477 (+-0.0910419) 5.00025 (+-0.000747134)
found 5.9 (+-0.000465984) 19.9476 (+-0.0910211) 5.00021 (+-0.000746963)
found -9.1 (+-0.000520801) 15.9581 (+-0.0814088) 4.00017 (+-0.00066808)
found -0.699999 (+-0.000523105) 15.9584 (+-0.0814569) 4.00025 (+-0.000668475)
found 7.7 (+-0.00052215) 15.9582 (+-0.0814364) 4.00021 (+-0.000668306)
found 8.9 (+-0.00052215) 15.9582 (+-0.0814364) 4.00021 (+-0.000668306)
found -8.5 (+-0.000520127) 15.958 (+-0.0813943) 4.00014 (+-0.000667961)
found 2.3 (+-0.000604286) 11.9688 (+-0.0705482) 3.00019 (+-0.000578953)
found -6.10001 (+-0.000740131) 7.97932 (+-0.0576058) 2.00016 (+-0.000472741)
found 3.49999 (+-0.000737728) 7.97922 (+-0.0575817) 2.00013 (+-0.000472544)
found -5.5 (+-0.000731897) 7.97886 (+-0.0575176) 2.00004 (+-0.000472017)
found 1.1 (+-0.00105879) 3.98997 (+-0.0407997) 1.00016 (+-0.000334822)
found -4.89999 (+-0.0010538) 3.98987 (+-0.040773) 1.00013 (+-0.000334603)
found 4.10002 (+-0.00105481) 3.98993 (+-0.0407793) 1.00014 (+-0.000334655)
#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>
{
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.));
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)
{
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
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);
}
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);
pfit->GetSigma(sigma, sigmaErr);
// current TSpectrumFit needs a sqrt(2) correction factor for sigma
// convert "bin numbers" into "x-axis values"
sigma *= 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;
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;
}
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;
}
#define d(i)
Definition RSha256.hxx:102
#define f(i)
Definition RSha256.hxx:104
#define h(i)
Definition RSha256.hxx:106
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t dest
Option_t Option_t TPoint TPoint const char x1
float xmin
float xmax
#define gROOT
Definition TROOT.h:406
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:233
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:645
A PolyMarker is defined by an array on N points in a 2-D space.
Definition TPolyMarker.h:31
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:615
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:682
Advanced 1-dimensional spectra fitting functions.
Advanced Spectra Processing.
Definition TSpectrum.h:18
Int_t SearchHighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
One-dimensional high-resolution peak search function.
Double_t * GetPositionX() const
Definition TSpectrum.h:58
const Double_t sigma
constexpr Double_t Sqrt2()
Definition TMath.h:86
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
constexpr Double_t TwoPi()
Definition TMath.h:44
Author

Definition in file FitAwmi.C.