Example to illustrate high resolution peak searching function (class TSpectrum2).
Author: Miroslav Morhac, Olivier Couet
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, April 17, 2024 at 11:22 AM.
const Int_t nbinsx = 64;
const Int_t nbinsy = 64;
std::vector<Double_t *> source(nbinsx), dest(nbinsx);
for (Int_t i = 0; i < nbinsx; i++) {
source[i] = new Double_t[nbinsy];
dest[i] = new Double_t[nbinsy];
}
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/spectrum/TSpectrum2.root";
TFile *f = TFile::Open(file.Data());
gStyle->SetOptStat(0);
auto search = (TH2F *)f->Get("search1");
TSpectrum2 s;
for (Int_t i = 0; i < nbinsx; i++) {
for (Int_t j = 0; j < nbinsy; j++) {
source[i][j] = search->GetBinContent(i + 1, j + 1);
}
}
Int_t nfound = s.SearchHighRes(source.data(), dest.data(), nbinsx, nbinsy, 2, 2, kFALSE, 3, kFALSE, 1); // 3, 10, 100
printf("Found %d candidate peaks\n", nfound);
Double_t *PositionX = s.GetPositionX();
Double_t *PositionY = s.GetPositionY();
search->Draw("CONT");
TMarker m;
m.SetMarkerStyle(23);
m.SetMarkerColor(kRed);
for (Int_t i = 0; i < nfound; i++) {
printf("posx= %d, posy= %d, value=%d\n", (Int_t)(PositionX[i] + 0.5), (Int_t)(PositionY[i] + 0.5),
(Int_t)source[(Int_t)(PositionX[i] + 0.5)][(Int_t)(PositionY[i] + 0.5)]);
m.DrawMarker(PositionX[i], PositionY[i]);
}
for (Int_t i = 0; i < nbinsx; i++) {
delete[] source[i];
delete[] dest[i];
}
Found 13 candidate peaks posx= 46, posy= 44, value=1962 posx= 39, posy= 4, value=1829 posx= 46, posy= 19, value=1762 posx= 50, posy= 32, value=1745 posx= 23, posy= 56, value=1710 posx= 35, posy= 48, value=1546 posx= 19, posy= 3, value=1524 posx= 39, posy= 58, value=1499 posx= 15, posy= 22, value=1376 posx= 38, posy= 16, value=1170 posx= 4, posy= 21, value=1166 posx= 10, posy= 16, value=1013 posx= 8, posy= 35, value=708
input_line_43:13:1: warning: 'search' shadows a declaration with the same name in the 'std' namespace; use '::search' to reference this declaration auto search = (TH2F *)f->Get("search1"); ^ Info in <TCanvas::MakeDefCanvas>: created default TCanvas with name c1
Draw all canvases
gROOT->GetListOfCanvases()->Draw()