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/// Example to illustrate deconvolution function (class TSpectrum).
4///
5/// \macro_image
6/// \macro_code
7///
8/// \authors Miroslav Morhac, Olivier Couet
9
11{
12 Int_t i;
13 const Int_t nbins = 256;
14 Double_t xmin = 0;
15 Double_t xmax = nbins;
16 Double_t source[nbins];
17 Double_t response[nbins];
18 gROOT->ForceStyle();
19
20 TH1F *h = new TH1F("h", "Deconvolution", nbins, xmin, xmax);
21 TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
22
23 TString dir = gROOT->GetTutorialDir();
24 TString file = dir + "/legacy/spectrum/TSpectrum.root";
25 TFile *f = new TFile(file.Data());
26 h = (TH1F *)f->Get("decon3");
27 h->SetTitle("Deconvolution of closely positioned overlapping peaks using boosted Gold deconvolution method");
28 d = (TH1F *)f->Get("decon_response_wide");
29
30 for (i = 0; i < nbins; i++)
31 source[i] = h->GetBinContent(i + 1);
32 for (i = 0; i < nbins; i++)
33 response[i] = d->GetBinContent(i + 1);
34
35 h->SetMaximum(200000);
36 h->Draw("L");
37 TSpectrum *s = new TSpectrum();
38 s->Deconvolution(source, response, 256, 200, 50, 1.2);
39
40 for (i = 0; i < nbins; i++)
41 d->SetBinContent(i + 1, source[i]);
42 d->SetLineColor(kRed);
43 d->Draw("SAME L");
44}
#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
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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:131
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:645
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