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