ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
peaks2.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_spectrum
3 /// Example to illustrate the 2-d peak finder (class TSpectrum2).
4 ///
5 /// This script generates a random number of 2-d gaussian peaks
6 /// The position of the peaks is found via TSpectrum2
7 /// To execute this example, do:
8 ///
9 /// ~~~ {.cpp}
10 /// root > .x peaks2.C (generate up to 50 peaks by default)
11 /// root > .x peaks2.C(10) (generate up to 10 peaks)
12 /// root > .x peaks2.C+(200) (generate up to 200 peaks via ACLIC)
13 /// ~~~
14 ///
15 /// The script will iterate generating a new histogram having
16 /// between 5 and the maximun number of peaks specified.
17 /// Double Click on the bottom right corner of the pad to go to a new spectrum
18 /// To Quit, select the "quit" item in the canvas "File" menu
19 ///
20 /// \macro_image
21 /// \macro_code
22 ///
23 /// \author Rene Brun
24 
25 #include "TSpectrum2.h"
26 #include "TCanvas.h"
27 #include "TRandom.h"
28 #include "TH2.h"
29 #include "TF2.h"
30 #include "TMath.h"
31 #include "TROOT.h"
32 
33 TSpectrum2 *s;
34 TH2F *h2 = 0;
35 Int_t npeaks = 30;
36 Double_t fpeaks2(Double_t *x, Double_t *par) {
37  Double_t result = 0.1;
38  for (Int_t p=0;p<npeaks;p++) {
39  Double_t norm = par[5*p+0];
40  Double_t mean1 = par[5*p+1];
41  Double_t sigma1 = par[5*p+2];
42  Double_t mean2 = par[5*p+3];
43  Double_t sigma2 = par[5*p+4];
44  result += norm*TMath::Gaus(x[0],mean1,sigma1)*TMath::Gaus(x[1],mean2,sigma2);
45  }
46  return result;
47 }
48 void findPeak2() {
49  printf("Generating histogram with %d peaks\n",npeaks);
50  Int_t nbinsx = 200;
51  Int_t nbinsy = 200;
52  Double_t xmin = 0;
53  Double_t xmax = (Double_t)nbinsx;
54  Double_t ymin = 0;
55  Double_t ymax = (Double_t)nbinsy;
56  Double_t dx = (xmax-xmin)/nbinsx;
57  Double_t dy = (ymax-ymin)/nbinsy;
58  delete h2;
59  h2 = new TH2F("h2","test",nbinsx,xmin,xmax,nbinsy,ymin,ymax);
60  h2->SetStats(0);
61  //generate n peaks at random
62  Double_t par[3000];
63  Int_t p;
64  for (p=0;p<npeaks;p++) {
65  par[5*p+0] = gRandom->Uniform(0.2,1);
66  par[5*p+1] = gRandom->Uniform(xmin,xmax);
67  par[5*p+2] = gRandom->Uniform(dx,5*dx);
68  par[5*p+3] = gRandom->Uniform(ymin,ymax);
69  par[5*p+4] = gRandom->Uniform(dy,5*dy);
70  }
71  TF2 *f2 = new TF2("f2",fpeaks2,xmin,xmax,ymin,ymax,5*npeaks);
72  f2->SetNpx(100);
73  f2->SetNpy(100);
74  f2->SetParameters(par);
75  TCanvas *c1 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("c1");
76  if (!c1) c1 = new TCanvas("c1","c1",10,10,1000,700);
77  h2->FillRandom("f2",500000);
78 
79  //now the real stuff: Finding the peaks
80  Int_t nfound = s->Search(h2,2,"col");
81 
82  //searching good and ghost peaks (approximation)
83  Int_t pf,ngood = 0;
84  Double_t *xpeaks = s->GetPositionX();
85  Double_t *ypeaks = s->GetPositionY();
86  for (p=0;p<npeaks;p++) {
87  for (pf=0;pf<nfound;pf++) {
88  Double_t diffx = TMath::Abs(xpeaks[pf] - par[5*p+1]);
89  Double_t diffy = TMath::Abs(ypeaks[pf] - par[5*p+3]);
90  if (diffx < 2*dx && diffy < 2*dy) ngood++;
91  }
92  }
93  if (ngood > nfound) ngood = nfound;
94  //Search ghost peaks (approximation)
95  Int_t nghost = 0;
96  for (pf=0;pf<nfound;pf++) {
97  Int_t nf=0;
98  for (p=0;p<npeaks;p++) {
99  Double_t diffx = TMath::Abs(xpeaks[pf] - par[5*p+1]);
100  Double_t diffy = TMath::Abs(ypeaks[pf] - par[5*p+3]);
101  if (diffx < 2*dx && diffy < 2*dy) nf++;
102  }
103  if (nf == 0) nghost++;
104  }
105  c1->Update();
106 
107  s->Print();
108  printf("Gener=%d, Found=%d, Good=%d, Ghost=%d\n",npeaks,nfound,ngood,nghost);
109  printf("\nDouble click in the bottom right corner of the pad to continue\n");
110  c1->WaitPrimitive();
111 }
112 void peaks2(Int_t maxpeaks=50) {
113  s = new TSpectrum2(2*maxpeaks);
114  for (int i=0; i<10; ++i) {
115  npeaks = (Int_t)gRandom->Uniform(5,maxpeaks);
116  findPeak2();
117  }
118 }
119 
120 
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
float xmin
Definition: THbookFile.cxx:93
virtual void FillRandom(const char *fname, Int_t ntimes=5000)
Fill histogram following distribution in function fname.
Definition: TH2.cxx:592
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3116
float ymin
Definition: THbookFile.cxx:93
Advanced 2-dimensional spectra processing.
Definition: TSpectrum2.h:20
#define gROOT
Definition: TROOT.h:344
int Int_t
Definition: RtypesCore.h:41
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t * GetPositionX() const
Definition: TSpectrum2.h:46
Double_t * GetPositionY() const
Definition: TSpectrum2.h:47
TH2D * h2
Definition: fit2dHist.C:45
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
virtual void SetNpy(Int_t npy=100)
Set the number of points used to draw the function.
Definition: TF2.cxx:907
A 2-Dim function with parameters.
Definition: TF2.h:33
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
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:453
virtual void Print(Option_t *option="") const
Print the array of positions.
Definition: TSpectrum2.cxx:162
The Canvas class.
Definition: TCanvas.h:48
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:606
virtual TObject * WaitPrimitive(const char *pname="", const char *emode="")
Loop and sleep until a primitive with name=pname is found in the pad.
Definition: TPad.cxx:6231
double result[121]
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2179
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
This function searches for peaks in source spectrum in hin The number of found peaks and their positi...
Definition: TSpectrum2.cxx:202