Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
Src4.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 TSpectrum2).
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \authors Miroslav Morhac, Olivier Couet
11 
12 #include <TSpectrum2.h>
13 
14 void Src4() {
15  Int_t i, j, nfound;
16  const Int_t nbinsx = 64;
17  const Int_t nbinsy = 64;
18  Double_t xmin = 0;
19  Double_t xmax = (Double_t)nbinsx;
20  Double_t ymin = 0;
21  Double_t ymax = (Double_t)nbinsy;
22  Double_t** source = new Double_t*[nbinsx];
23  for (i=0;i<nbinsx;i++)
24  source[i]=new Double_t[nbinsy];
25  Double_t** dest = new Double_t*[nbinsx];
26  for (i=0;i<nbinsx;i++)
27  dest[i]=new Double_t[nbinsy];
28  TString dir = gROOT->GetTutorialDir();
29  TString file = dir+"/spectrum/TSpectrum2.root";
30  TFile *f = new TFile(file.Data());
31  gStyle->SetOptStat(0);
32  auto search = (TH2F*) f->Get("search2");
33  auto s = new TSpectrum2();
34  for (i = 0; i < nbinsx; i++){
35  for (j = 0; j < nbinsy; j++){
36  source[i][j] = search->GetBinContent(i + 1,j + 1);
37  }
38  }
39  nfound = s->SearchHighRes(source, dest, nbinsx, nbinsy, 3, 5, kFALSE, 10, kTRUE, 3);
40  printf("Found %d candidate peaks\n",nfound);
41  Double_t *PositionX = s->GetPositionX();
42  Double_t *PositionY = s->GetPositionY();
43  search->Draw("CONT");
44  auto m = new TMarker();
45  m->SetMarkerStyle(23);
46  m->SetMarkerColor(kRed);
47  for (i=0;i<nfound;i++) {
48  printf("posx= %d, posy= %d, value=%d\n",(Int_t)(PositionX[i]+0.5), (Int_t)(PositionY[i]+0.5),
49  (Int_t)source[(Int_t)(PositionX[i]+0.5)][(Int_t)(PositionY[i]+0.5)]);
50  m->DrawMarker(PositionX[i],PositionY[i]);
51  }
52 }