Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
DeconvolutionRL_wide.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 DeconvolutionRL_wide() {
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 TString dir = gROOT->GetTutorialDir();
21 TString file = dir+"/spectrum/TSpectrum.root";
22 TFile *f = new TFile(file.Data());
23 TH1F* h = (TH1F*) f->Get("decon3");
24 h->SetTitle("Deconvolution of closely positioned overlapping peaks using Richardson-Lucy deconvolution method");
25 TH1F* d = (TH1F*) f->Get("decon_response_wide");
26
27 for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
28 for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
29
30 h->SetMaximum(30000);
31 h->Draw("L");
32 TSpectrum *s = new TSpectrum();
33 s->DeconvolutionRL(source,response,256,10000,1,1);
34
35 for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
36 d->SetLineColor(kRed);
37 d->Draw("SAME L");
38}
#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:621
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6686
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5029
Advanced Spectra Processing.
Definition TSpectrum.h:18
const char * DeconvolutionRL(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