Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SearchHR3.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// \notebook
4/// Example to illustrate the influence of number of iterations in deconvolution
5/// in high resolution peak searching function (class TSpectrum).
6///
7/// \macro_output
8/// \macro_image
9/// \macro_code
10///
11/// \authors Miroslav Morhac, Olivier Couet
12
13void SearchHR3() {
14 Double_t fPositionX[100];
15 Double_t fPositionY[100];
16 Int_t fNPeaks = 0;
17 Int_t i,nfound,bin;
18 const Int_t nbins = 1024;
19 Double_t xmin = 0;
20 Double_t xmax = nbins;
21 Double_t a;
22 Double_t source[nbins], dest[nbins];
23 gROOT->ForceStyle();
24
25 TString dir = gROOT->GetTutorialDir();
26 TString file = dir+"/spectrum/TSpectrum.root";
27 TFile *f = new TFile(file.Data());
28 TH1F *h = (TH1F*) f->Get("back2");
29 h->SetTitle("Influence of # of iterations in deconvolution in peak searching");
30 h->GetXaxis()->SetRange(1,nbins);
31
32 TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
33 TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
34 TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
35 TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
36
37 TSpectrum *s = new TSpectrum();
38
39 for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
40 nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
41
42 Double_t *xpeaks = s->GetPositionX();
43 for (i = 0; i < nfound; i++) {
44 a=xpeaks[i];
45 bin = 1 + Int_t(a + 0.5);
46 fPositionX[i] = h->GetBinCenter(bin);
47 fPositionY[i] = h->GetBinContent(bin);
48 }
49 TPolyMarker * pm = (TPolyMarker*)h->GetListOfFunctions()->FindObject("TPolyMarker");
50 if (pm) {
51 h->GetListOfFunctions()->Remove(pm);
52 delete pm;
53 }
54 pm = new TPolyMarker(nfound, fPositionX, fPositionY);
55 h->GetListOfFunctions()->Add(pm);
56 pm->SetMarkerStyle(23);
58 pm->SetMarkerSize(1.3);
59
60 h->Draw("L");
61
62 for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,dest[i]);
63 d1->SetLineColor(kRed);
64 d1->Draw("SAME");
65
66 for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
67 s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 10, kTRUE, 3);
68 for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,dest[i]);
70 d2->Draw("SAME");
71
72 for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
73 s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 100, kTRUE, 3);
74 for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,dest[i]);
76 d3->Draw("SAME");
77
78 for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
79 s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 1000, kTRUE, 3);
80 for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,dest[i]);
82 d4->Draw("SAME");
83
84 printf("Found %d candidate peaks\n",nfound);
85}
#define f(i)
Definition RSha256.hxx:104
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kBlue
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
float xmin
float xmax
#define gROOT
Definition TROOT.h:406
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
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
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
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9190
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:403
A PolyMarker is defined by an array on N points in a 2-D space.
Definition TPolyMarker.h:31
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
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376