Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
Deconvolution_wide_boost.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// \notebook
4/// Example to illustrate deconvolution function (class TSpectrum).
5///
6/// \macro_image
7/// \macro_code
8///
9/// \authors Miroslav Morhac, Olivier Couet
10
11void Deconvolution_wide_boost()
12{
13 Int_t i;
14 const Int_t nbins = 256;
15 Double_t xmin = 0;
16 Double_t xmax = nbins;
17 Double_t source[nbins];
18 Double_t response[nbins];
19 gROOT->ForceStyle();
20
21 TH1F *h = new TH1F("h", "Deconvolution", nbins, xmin, xmax);
22 TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
23
24 TString dir = gROOT->GetTutorialDir();
25 TString file = dir + "/legacy/spectrum/TSpectrum.root";
26 TFile *f = new TFile(file.Data());
27 h = (TH1F *)f->Get("decon3");
28 h->SetTitle("Deconvolution of closely positioned overlapping peaks using boosted Gold deconvolution method");
29 d = (TH1F *)f->Get("decon_response_wide");
30
31 for (i = 0; i < nbins; i++)
32 source[i] = h->GetBinContent(i + 1);
33 for (i = 0; i < nbins; i++)
34 response[i] = d->GetBinContent(i + 1);
35
36 h->SetMaximum(200000);
37 h->Draw("L");
38 TSpectrum *s = new TSpectrum();
39 s->Deconvolution(source, response, 256, 200, 50, 1.2);
40
41 for (i = 0; i < nbins; i++)
42 d->SetBinContent(i + 1, source[i]);
43 d->SetLineColor(kRed);
44 d->Draw("SAME L");
45}
#define d(i)
Definition RSha256.hxx:102
#define f(i)
Definition RSha256.hxx:104
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
float xmin
float xmax
#define gROOT
Definition TROOT.h:406
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:634
Advanced Spectra Processing.
Definition TSpectrum.h:18
const char * Deconvolution(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376