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 Thursday, February 29, 2024 at 12:27 PM.
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("search4");
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, 5, kTRUE, 3, kFALSE, 3);
printf("Found %d candidate peaks\n", nfound);
Double_t *PositionX = s.GetPositionX();
Double_t *PositionY = s.GetPositionY();
search->Draw("COL");
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 5 candidate peaks posx= 32, posy= 18, value=2164 posx= 29, posy= 42, value=1959 posx= 24, posy= 26, value=1667 posx= 16, posy= 29, value=1631 posx= 20, posy= 46, value=1334
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("search4"); ^ Info in <TCanvas::MakeDefCanvas>: created default TCanvas with name c1
Draw all canvases
gROOT->GetListOfCanvases()->Draw()