Logo ROOT   6.08/07
Reference Guide
peaks.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_spectrum
3 /// \notebook
4 /// Getting Contours From TH2D.
5 /// Illustrates how to find peaks in histograms.
6 ///
7 /// This script generates a random number of gaussian peaks
8 /// on top of a linear background.
9 /// The position of the peaks is found via TSpectrum and injected
10 /// as initial values of parameters to make a global fit.
11 /// The background is computed and drawn on top of the original histogram.
12 ///
13 /// To execute this example, do:
14 ///
15 /// ~~~{.cpp}
16 /// root > .x peaks.C (generate 10 peaks by default)
17 /// root > .x peaks.C++ (use the compiler)
18 /// root > .x peaks.C++(30) (generates 30 peaks)
19 /// ~~~
20 ///
21 /// To execute only the first part of the script (without fitting)
22 /// specify a negative value for the number of peaks, eg
23 ///
24 /// ~~~{.cpp}
25 /// root > .x peaks.C(-20)
26 /// ~~~
27 ///
28 /// \macro_image
29 /// \macro_code
30 ///
31 /// \author Rene Brun
32 
33 #include "TCanvas.h"
34 #include "TMath.h"
35 #include "TH1.h"
36 #include "TF1.h"
37 #include "TRandom.h"
38 #include "TSpectrum.h"
39 #include "TVirtualFitter.h"
40 
41 Int_t npeaks = 30;
42 Double_t fpeaks(Double_t *x, Double_t *par) {
43  Double_t result = par[0] + par[1]*x[0];
44  for (Int_t p=0;p<npeaks;p++) {
45  Double_t norm = par[3*p+2];
46  Double_t mean = par[3*p+3];
47  Double_t sigma = par[3*p+4];
48  result += norm*TMath::Gaus(x[0],mean,sigma);
49  }
50  return result;
51 }
52 void peaks(Int_t np=10) {
53  npeaks = TMath::Abs(np);
54  TH1F *h = new TH1F("h","test",500,0,1000);
55  //generate n peaks at random
56  Double_t par[3000];
57  par[0] = 0.8;
58  par[1] = -0.6/1000;
59  Int_t p;
60  for (p=0;p<npeaks;p++) {
61  par[3*p+2] = 1;
62  par[3*p+3] = 10+gRandom->Rndm()*980;
63  par[3*p+4] = 3+2*gRandom->Rndm();
64  }
65  TF1 *f = new TF1("f",fpeaks,0,1000,2+3*npeaks);
66  f->SetNpx(1000);
67  f->SetParameters(par);
68  TCanvas *c1 = new TCanvas("c1","c1",10,10,1000,900);
69  c1->Divide(1,2);
70  c1->cd(1);
71  h->FillRandom("f",200000);
72  h->Draw();
73  TH1F *h2 = (TH1F*)h->Clone("h2");
74  //Use TSpectrum to find the peak candidates
75  TSpectrum *s = new TSpectrum(2*npeaks);
76  Int_t nfound = s->Search(h,2,"",0.10);
77  printf("Found %d candidate peaks to fit\n",nfound);
78  //Estimate background using TSpectrum::Background
79  TH1 *hb = s->Background(h,20,"same");
80  if (hb) c1->Update();
81  if (np <0) return;
82 
83  //estimate linear background using a fitting method
84  c1->cd(2);
85  TF1 *fline = new TF1("fline","pol1",0,1000);
86  h->Fit("fline","qn");
87  //Loop on all found peaks. Eliminate peaks at the background level
88  par[0] = fline->GetParameter(0);
89  par[1] = fline->GetParameter(1);
90  npeaks = 0;
91  Double_t *xpeaks = s->GetPositionX();
92  for (p=0;p<nfound;p++) {
93  Double_t xp = xpeaks[p];
94  Int_t bin = h->GetXaxis()->FindBin(xp);
95  Double_t yp = h->GetBinContent(bin);
96  if (yp-TMath::Sqrt(yp) < fline->Eval(xp)) continue;
97  par[3*npeaks+2] = yp;
98  par[3*npeaks+3] = xp;
99  par[3*npeaks+4] = 3;
100  npeaks++;
101  }
102  printf("Found %d useful peaks to fit\n",npeaks);
103  printf("Now fitting: Be patient\n");
104  TF1 *fit = new TF1("fit",fpeaks,0,1000,2+3*npeaks);
105  //we may have more than the default 25 parameters
106  TVirtualFitter::Fitter(h2,10+3*npeaks);
107  fit->SetParameters(par);
108  fit->SetNpx(1000);
109  h2->Fit("fit");
110 }
double par[1]
Definition: unuranDistr.cxx:38
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:439
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3156
return c1
Definition: legend1.C:41
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:302
TH1 * h
Definition: legend2.C:5
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4638
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
int Int_t
Definition: RtypesCore.h:41
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t x[n]
Definition: legend1.C:17
const Double_t sigma
Double_t * GetPositionX() const
Definition: TSpectrum.h:60
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:512
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
One-dimensional peak search function.
Definition: TSpectrum.cxx:266
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH1.cxx:3293
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2851
R__EXTERN TRandom * gRandom
Definition: TRandom.h:66
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculate a gaussian function with mean and sigma.
Definition: TMath.cxx:452
virtual Int_t FindBin(Double_t x)
Find bin number corresponding to abscissa x.
Definition: TAxis.cxx:279
The Canvas class.
Definition: TCanvas.h:41
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition: TF1.cxx:1196
double f(double x)
double Double_t
Definition: RtypesCore.h:55
virtual TH1 * Background(const TH1 *hist, Int_t niter=20, Option_t *option="")
One-dimensional background estimation function.
Definition: TSpectrum.cxx:152
The TH1 histogram class.
Definition: TH1.h:80
Advanced Spectra Processing.
Definition: TSpectrum.h:20
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1089
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:359
1-Dim function class
Definition: TF1.h:149
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TH1.cxx:2544
double result[121]
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2183
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3563
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
TAxis * GetXaxis()
Definition: TH1.h:324