Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SearchHR3.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// Example to illustrate the influence of number of iterations in deconvolution
4/// in high resolution peak searching function (class TSpectrum).
5///
6/// \macro_output
7/// \macro_image
8/// \macro_code
9///
10/// \authors Miroslav Morhac, Olivier Couet
11
12void SearchHR3()
13{
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 + "/legacy/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++)
40 source[i] = h->GetBinContent(i + 1);
41 nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
42
44 for (i = 0; i < nfound; i++) {
45 a = xpeaks[i];
46 bin = 1 + Int_t(a + 0.5);
47 fPositionX[i] = h->GetBinCenter(bin);
48 fPositionY[i] = h->GetBinContent(bin);
49 }
50 TPolyMarker *pm = (TPolyMarker *)h->GetListOfFunctions()->FindObject("TPolyMarker");
51 if (pm) {
52 h->GetListOfFunctions()->Remove(pm);
53 delete pm;
54 }
55 pm = new TPolyMarker(nfound, fPositionX, fPositionY);
56 h->GetListOfFunctions()->Add(pm);
57 pm->SetMarkerStyle(23);
58 pm->SetMarkerColor(kRed);
59 pm->SetMarkerSize(1.3);
60
61 h->Draw("L");
62
63 for (i = 0; i < nbins; i++)
64 d1->SetBinContent(i + 1, dest[i]);
65 d1->SetLineColor(kRed);
66 d1->Draw("SAME");
67
68 for (i = 0; i < nbins; i++)
69 source[i] = h->GetBinContent(i + 1);
70 s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 10, kTRUE, 3);
71 for (i = 0; i < nbins; i++)
72 d2->SetBinContent(i + 1, dest[i]);
73 d2->SetLineColor(kBlue);
74 d2->Draw("SAME");
75
76 for (i = 0; i < nbins; i++)
77 source[i] = h->GetBinContent(i + 1);
78 s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 100, kTRUE, 3);
79 for (i = 0; i < nbins; i++)
80 d3->SetBinContent(i + 1, dest[i]);
81 d3->SetLineColor(kGreen);
82 d3->Draw("SAME");
83
84 for (i = 0; i < nbins; i++)
85 source[i] = h->GetBinContent(i + 1);
86 s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 1000, kTRUE, 3);
87 for (i = 0; i < nbins; i++)
88 d4->SetBinContent(i + 1, dest[i]);
89 d4->SetLineColor(kMagenta);
90 d4->Draw("SAME");
91
92 printf("Found %d candidate peaks\n", nfound);
93}
#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:93
@ kRed
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kBlue
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
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
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6716
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