Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SearchHR1.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// Example to illustrate high resolution peak searching function (class TSpectrum).
4///
5/// \macro_output
6/// \macro_image
7/// \macro_code
8///
9/// \authors Miroslav Morhac, Olivier Couet
10
11void SearchHR1()
12{
13 Double_t fPositionX[100];
14 Double_t fPositionY[100];
15 Int_t fNPeaks = 0;
16 Int_t i, nfound, bin;
17 const Int_t nbins = 1024;
18 Double_t xmin = 0;
19 Double_t xmax = nbins;
20 Double_t a;
21 Double_t source[nbins], dest[nbins];
22 gROOT->ForceStyle();
23
24 TString dir = gROOT->GetTutorialDir();
25 TString file = dir + "/legacy/spectrum/TSpectrum.root";
26 TFile *f = new TFile(file.Data());
27 TH1F *h = (TH1F *)f->Get("back2");
28 h->SetTitle("High resolution peak searching, number of iterations = 3");
29 h->GetXaxis()->SetRange(1, nbins);
30 TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
31 h->Draw("L");
32
33 for (i = 0; i < nbins; i++)
34 source[i] = h->GetBinContent(i + 1);
35
36 h->Draw("L");
37
38 TSpectrum *s = new TSpectrum();
39
40 nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
42 for (i = 0; i < nfound; i++) {
43 a = xpeaks[i];
44 bin = 1 + Int_t(a + 0.5);
45 fPositionX[i] = h->GetBinCenter(bin);
46 fPositionY[i] = h->GetBinContent(bin);
47 }
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);
57 pm->SetMarkerColor(kRed);
58 pm->SetMarkerSize(1.3);
59
60 for (i = 0; i < nbins; i++)
61 d->SetBinContent(i + 1, dest[i]);
62 d->SetLineColor(kRed);
63 d->Draw("SAME");
64
65 printf("Found %d candidate peaks\n", nfound);
66 for (i = 0; i < nfound; i++)
67 printf("posx= %f, posy= %f\n", fPositionX[i], fPositionY[i]);
68}
#define d(i)
Definition RSha256.hxx:102
#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
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:53
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