Logo ROOT  
Reference Guide
FitAwmi.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// \notebook
4///
5/// This macro fits the source spectrum using the AWMI algorithm
6/// from the "TSpectrumFit" class ("TSpectrum" class is used to find peaks).
7///
8/// To try this macro, in a ROOT (5 or 6) prompt, do:
9///
10/// ~~~{.cpp}
11/// root > .x FitAwmi.C
12/// ~~~
13///
14/// or:
15///
16/// ~~~{.cpp}
17/// root > .x FitAwmi.C++
18/// root > FitAwmi(); // re-run with another random set of peaks
19/// ~~~
20///
21/// \macro_image
22/// \macro_output
23/// \macro_code
24///
25/// \author
26
27#include "TROOT.h"
28#include "TMath.h"
29#include "TRandom.h"
30#include "TH1.h"
31#include "TF1.h"
32#include "TCanvas.h"
33#include "TSpectrum.h"
34#include "TSpectrumFit.h"
35#include "TPolyMarker.h"
36#include "TList.h"
37
38#include <iostream>
39
40TH1F *FitAwmi_Create_Spectrum(void) {
41 Int_t nbins = 1000;
42 Double_t xmin = -10., xmax = 10.;
43 delete gROOT->FindObject("h"); // prevent "memory leak"
44 TH1F *h = new TH1F("h", "simulated spectrum", nbins, xmin, xmax);
45 h->SetStats(kFALSE);
46 TF1 f("f", "TMath::Gaus(x, [0], [1], 1)", xmin, xmax);
47 // f.SetParNames("mean", "sigma");
48 gRandom->SetSeed(0); // make it really random
49 // create well separated peaks with exactly known means and areas
50 // note: TSpectrumFit assumes that all peaks have the same sigma
51 Double_t sigma = (xmax - xmin) / Double_t(nbins) * Int_t(gRandom->Uniform(2., 6.));
52 Int_t npeaks = 0;
53 while (xmax > (xmin + 6. * sigma)) {
54 npeaks++;
55 xmin += 3. * sigma; // "mean"
56 f.SetParameters(xmin, sigma);
57 Double_t area = 1. * Int_t(gRandom->Uniform(1., 11.));
58 h->Add(&f, area, ""); // "" ... or ... "I"
59 std::cout << "created "
60 << xmin << " "
61 << (area / sigma / TMath::Sqrt(TMath::TwoPi())) << " "
62 << area << std::endl;
63 xmin += 3. * sigma;
64 }
65 std::cout << "the total number of created peaks = " << npeaks
66 << " with sigma = " << sigma << std::endl;
67 return h;
68}
69
70void FitAwmi(void) {
71
72 TH1F *h = FitAwmi_Create_Spectrum();
73
74 TCanvas *cFit = ((TCanvas *)(gROOT->GetListOfCanvases()->FindObject("cFit")));
75 if (!cFit) cFit = new TCanvas("cFit", "cFit", 10, 10, 1000, 700);
76 else cFit->Clear();
77 h->Draw("L");
78 Int_t i, nfound, bin;
79 Int_t nbins = h->GetNbinsX();
80
81 Double_t *source = new Double_t[nbins];
82 Double_t *dest = new Double_t[nbins];
83
84 for (i = 0; i < nbins; i++) source[i] = h->GetBinContent(i + 1);
85 TSpectrum *s = new TSpectrum(); // note: default maxpositions = 100
86 // searching for candidate peaks positions
87 nfound = s->SearchHighRes(source, dest, nbins, 2., 2., kFALSE, 10000, kFALSE, 0);
88 // filling in the initial estimates of the input parameters
89 Bool_t *FixPos = new Bool_t[nfound];
90 Bool_t *FixAmp = new Bool_t[nfound];
91 for(i = 0; i < nfound; i++) FixAmp[i] = FixPos[i] = kFALSE;
92
93 Double_t *Pos, *Amp = new Double_t[nfound]; // ROOT 6
94
95 Pos = s->GetPositionX(); // 0 ... (nbins - 1)
96 for (i = 0; i < nfound; i++) {
97 bin = 1 + Int_t(Pos[i] + 0.5); // the "nearest" bin
98 Amp[i] = h->GetBinContent(bin);
99 }
100 TSpectrumFit *pfit = new TSpectrumFit(nfound);
101 pfit->SetFitParameters(0, (nbins - 1), 1000, 0.1, pfit->kFitOptimChiCounts,
102 pfit->kFitAlphaHalving, pfit->kFitPower2,
104 pfit->SetPeakParameters(2., kFALSE, Pos, FixPos, Amp, FixAmp);
105 // pfit->SetBackgroundParameters(source[0], kFALSE, 0., kFALSE, 0., kFALSE);
106 pfit->FitAwmi(source);
107 Double_t *Positions = pfit->GetPositions();
108 Double_t *PositionsErrors = pfit->GetPositionsErrors();
109 Double_t *Amplitudes = pfit->GetAmplitudes();
110 Double_t *AmplitudesErrors = pfit->GetAmplitudesErrors();
111 Double_t *Areas = pfit->GetAreas();
112 Double_t *AreasErrors = pfit->GetAreasErrors();
113 delete gROOT->FindObject("d"); // prevent "memory leak"
114 TH1F *d = new TH1F(*h); d->SetNameTitle("d", ""); d->Reset("M");
115 for (i = 0; i < nbins; i++) d->SetBinContent(i + 1, source[i]);
116 Double_t x1 = d->GetBinCenter(1), dx = d->GetBinWidth(1);
117 Double_t sigma, sigmaErr;
118 pfit->GetSigma(sigma, sigmaErr);
119
120 // current TSpectrumFit needs a sqrt(2) correction factor for sigma
121 sigma /= TMath::Sqrt2(); sigmaErr /= TMath::Sqrt2();
122 // convert "bin numbers" into "x-axis values"
123 sigma *= dx; sigmaErr *= dx;
124
125 std::cout << "the total number of found peaks = " << nfound
126 << " with sigma = " << sigma << " (+-" << sigmaErr << ")"
127 << std::endl;
128 std::cout << "fit chi^2 = " << pfit->GetChi() << std::endl;
129 for (i = 0; i < nfound; i++) {
130 bin = 1 + Int_t(Positions[i] + 0.5); // the "nearest" bin
131 Pos[i] = d->GetBinCenter(bin);
132 Amp[i] = d->GetBinContent(bin);
133
134 // convert "bin numbers" into "x-axis values"
135 Positions[i] = x1 + Positions[i] * dx;
136 PositionsErrors[i] *= dx;
137 Areas[i] *= dx;
138 AreasErrors[i] *= dx;
139
140 std::cout << "found "
141 << Positions[i] << " (+-" << PositionsErrors[i] << ") "
142 << Amplitudes[i] << " (+-" << AmplitudesErrors[i] << ") "
143 << Areas[i] << " (+-" << AreasErrors[i] << ")"
144 << std::endl;
145 }
146 d->SetLineColor(kRed); d->SetLineWidth(1);
147 d->Draw("SAME L");
148 TPolyMarker *pm = ((TPolyMarker*)(h->GetListOfFunctions()->FindObject("TPolyMarker")));
149 if (pm) {
150 h->GetListOfFunctions()->Remove(pm);
151 delete pm;
152 }
153 pm = new TPolyMarker(nfound, Pos, Amp);
154 h->GetListOfFunctions()->Add(pm);
155 pm->SetMarkerStyle(23);
156 pm->SetMarkerColor(kRed);
157 pm->SetMarkerSize(1);
158 // cleanup
159 delete pfit;
160 delete [] Amp;
161 delete [] FixAmp;
162 delete [] FixPos;
163 delete s;
164 delete [] dest;
165 delete [] source;
166 return;
167}
#define d(i)
Definition: RSha256.hxx:102
#define f(i)
Definition: RSha256.hxx:104
#define h(i)
Definition: RSha256.hxx:106
static const double x1[5]
int Int_t
Definition: RtypesCore.h:43
const Bool_t kFALSE
Definition: RtypesCore.h:90
bool Bool_t
Definition: RtypesCore.h:61
double Double_t
Definition: RtypesCore.h:57
@ kRed
Definition: Rtypes.h:64
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
#define gROOT
Definition: TROOT.h:406
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:41
The Canvas class.
Definition: TCanvas.h:27
void Clear(Option_t *option="")
Remove all primitives from the canvas.
Definition: TCanvas.cxx:720
1-Dim function class
Definition: TF1.h:210
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
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:597
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:635
Advanced 1-dimensional spectra fitting functions.
Definition: TSpectrumFit.h:18
Double_t GetChi() const
Definition: TSpectrumFit.h:121
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:
Double_t * GetAmplitudesErrors() const
Definition: TSpectrumFit.h:117
void FitAwmi(Double_t *source)
This function fits the source spectrum.
Double_t * GetAreasErrors() const
Definition: TSpectrumFit.h:119
void GetSigma(Double_t &sigma, Double_t &sigmaErr)
This function gets the sigma parameter and its error.
@ kFitTaylorOrderFirst
Definition: TSpectrumFit.h:82
Double_t * GetAreas() const
Definition: TSpectrumFit.h:118
Double_t * GetAmplitudes() const
Definition: TSpectrumFit.h:116
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:
Double_t * GetPositionsErrors() const
Definition: TSpectrumFit.h:123
Double_t * GetPositions() const
Definition: TSpectrumFit.h:122
Advanced Spectra Processing.
Definition: TSpectrum.h:18
const Double_t sigma
static constexpr double s
constexpr Double_t Sqrt2()
Definition: TMath.h:89
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
constexpr Double_t TwoPi()
Definition: TMath.h:45
#define dest(otri, vertexptr)
Definition: triangle.c:1040