Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
FitAwmi.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// This macro fits the source spectrum using the AWMI algorithm
4/// from the "TSpectrumFit" class ("TSpectrum" class is used to find peaks).
5///
6/// To try this macro, in a ROOT (5 or 6) prompt, do:
7///
8/// ~~~{.cpp}
9/// root > .x FitAwmi.C
10/// ~~~
11///
12/// or:
13///
14/// ~~~{.cpp}
15/// root > .x FitAwmi.C++
16/// root > FitAwmi(); // re-run with another random set of peaks
17/// ~~~
18///
19/// \macro_image
20/// \macro_output
21/// \macro_code
22///
23/// \author
24
25#include "TROOT.h"
26#include "TMath.h"
27#include "TRandom.h"
28#include "TH1.h"
29#include "TF1.h"
30#include "TCanvas.h"
31#include "TSpectrum.h"
32#include "TSpectrumFit.h"
33#include "TPolyMarker.h"
34#include "TList.h"
35
36#include <iostream>
37
39{
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 " << xmin << " " << (area / sigma / TMath::Sqrt(TMath::TwoPi())) << " " << area
59 << std::endl;
60 xmin += 3. * sigma;
61 }
62 std::cout << "the total number of created peaks = " << npeaks << " with sigma = " << sigma << std::endl;
63 return h;
64}
65
66void FitAwmi(void)
67{
68
70
71 TCanvas *cFit = ((TCanvas *)(gROOT->GetListOfCanvases()->FindObject("cFit")));
72 if (!cFit)
73 cFit = new TCanvas("cFit", "cFit", 10, 10, 1000, 700);
74 else
75 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++)
84 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++)
92 FixAmp[i] = FixPos[i] = kFALSE;
93
94 Double_t *Pos, *Amp = new Double_t[nfound]; // ROOT 6
95
96 Pos = s->GetPositionX(); // 0 ... (nbins - 1)
97 for (i = 0; i < nfound; i++) {
98 bin = 1 + Int_t(Pos[i] + 0.5); // the "nearest" bin
99 Amp[i] = h->GetBinContent(bin);
100 }
102 pfit->SetFitParameters(0, (nbins - 1), 1000, 0.1, pfit->kFitOptimChiCounts, pfit->kFitAlphaHalving, pfit->kFitPower2,
103 pfit->kFitTaylorOrderFirst);
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);
115 d->SetNameTitle("d", "");
116 d->Reset("M");
117 for (i = 0; i < nbins; i++)
118 d->SetBinContent(i + 1, source[i]);
119 Double_t x1 = d->GetBinCenter(1), dx = d->GetBinWidth(1);
121 pfit->GetSigma(sigma, sigmaErr);
122
123 // current TSpectrumFit needs a sqrt(2) correction factor for sigma
124 sigma /= TMath::Sqrt2();
126 // convert "bin numbers" into "x-axis values"
127 sigma *= dx;
128 sigmaErr *= dx;
129
130 std::cout << "the total number of found peaks = " << nfound << " with sigma = " << sigma << " (+-" << sigmaErr << ")"
131 << std::endl;
132 std::cout << "fit chi^2 = " << pfit->GetChi() << std::endl;
133 for (i = 0; i < nfound; i++) {
134 bin = 1 + Int_t(Positions[i] + 0.5); // the "nearest" bin
135 Pos[i] = d->GetBinCenter(bin);
136 Amp[i] = d->GetBinContent(bin);
137
138 // convert "bin numbers" into "x-axis values"
139 Positions[i] = x1 + Positions[i] * dx;
140 PositionsErrors[i] *= dx;
141 Areas[i] *= dx;
142 AreasErrors[i] *= dx;
143
144 std::cout << "found " << Positions[i] << " (+-" << PositionsErrors[i] << ") " << Amplitudes[i] << " (+-"
145 << AmplitudesErrors[i] << ") " << Areas[i] << " (+-" << AreasErrors[i] << ")" << std::endl;
146 }
147 d->SetLineColor(kRed);
148 d->SetLineWidth(1);
149 d->Draw("SAME L");
150 TPolyMarker *pm = ((TPolyMarker *)(h->GetListOfFunctions()->FindObject("TPolyMarker")));
151 if (pm) {
152 h->GetListOfFunctions()->Remove(pm);
153 delete pm;
154 }
155 pm = new TPolyMarker(nfound, Pos, Amp);
156 h->GetListOfFunctions()->Add(pm);
157 pm->SetMarkerStyle(23);
158 pm->SetMarkerColor(kRed);
159 pm->SetMarkerSize(1);
160 // cleanup
161 delete pfit;
162 delete[] Amp;
163 delete[] FixAmp;
164 delete[] FixPos;
165 delete s;
166 delete[] dest;
167 delete[] source;
168 return;
169}
#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