Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
SearchHR1.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_spectrum
3/// \notebook
4/// Example to illustrate 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 SearchHR1() {
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+"/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++) source[i]=h->GetBinContent(i + 1);
34
35 h->Draw("L");
36
37 TSpectrum *s = new TSpectrum();
38
39 nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
40 Double_t *xpeaks = s->GetPositionX();
41 for (i = 0; i < nfound; i++) {
42 a=xpeaks[i];
43 bin = 1 + Int_t(a + 0.5);
44 fPositionX[i] = h->GetBinCenter(bin);
45 fPositionY[i] = h->GetBinContent(bin);
46 }
47
48 TPolyMarker * pm = (TPolyMarker*)h->GetListOfFunctions()->FindObject("TPolyMarker");
49 if (pm) {
50 h->GetListOfFunctions()->Remove(pm);
51 delete pm;
52 }
53 pm = new TPolyMarker(nfound, fPositionX, fPositionY);
54 h->GetListOfFunctions()->Add(pm);
55 pm->SetMarkerStyle(23);
57 pm->SetMarkerSize(1.3);
58
59 for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,dest[i]);
60 d->SetLineColor(kRed);
61 d->Draw("SAME");
62
63 printf("Found %d candidate peaks\n",nfound);
64 for( i=0;i<nfound;i++) printf("posx= %f, posy= %f\n",fPositionX[i], fPositionY[i]);
65}
#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
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 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:623
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6747
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:420
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