// root.cern.ch/root/htmldoc/guides/spectrum/Spectrum.html#dimensional-spectra
// Example to illustrate the influence of the clipping window width on the
// estimated background.
//
// Authors: Miroslav Morhac, Olivier Couet
void Background_width2 ()
{
Int_t i;
Int_t nbins = 4096;
Double_t xmin = 0;
Double_t xmax = 4096;
Double_t* source = new Double_t [nbins];
gROOT->ForceStyle ();
TFile* f = new TFile ("TSpectrum.root");
TH1F* hdata = (TH1F*) f->Get ("back2");
hdata->SetTitle ("Influence of clipping window width on the estimated background");
hdata->SetAxisRange (0, 1000);
hdata->SetMaximum (8000);
hdata->Draw ("L");
TSpectrum* Analyser1dim = new TSpectrum ();
/*** Case for width of clipping window = 10 ***/
for (i = 0; i < nbins; i++) source[i] = hdata->GetBinContent (i+1);
Analyser1dim->Background (
source, nbins, 10, TSpectrum::kBackDecreasingWindow,
TSpectrum::kBackOrder2, kFALSE,
TSpectrum::kBackSmoothing3, kFALSE );
/*
const char* TSpectrum::Background (
Double_t* spectrum, pointer to the array of source spectrum
Int_t ssize, length of the spectrum vector
Int_t numberIterations, maximal width of clipping window,
Int_t direction, direction of change of clipping window.
Possible values: kBackIncreasingWindow, kBackDecreasingWindow
Int_t filterOrder, order of clipping filter. Possible values:
kBackOrder2, kBackOrder4, kBackOrder6, kBackOrder8
bool smoothing, logical variable whether the smoothing operation in the
estimation of background will be included. (0/1)
Int_t smoothWindow, width of smoothing window. Possible values:
kBackSmoothing3, kBackSmoothing5, ..., kBackSmoothing15
bool compton whether the estimation of Compton edge will be included. (0/1)
)
See: root.cern.ch/doc/v608/classTSpectrum.html#a7830480612d1ee301964e6c79319f444
*/
TH1F* hwindow10 = new TH1F ("hwindow10", "", nbins, xmin, xmax);
for (i = 0; i < nbins; i++) hwindow10->SetBinContent (i+1, source[i]);
hwindow10->SetLineColor (kRed);
hwindow10->Draw ("SAME L");
/*** Case for width of clipping window = 20 ***/
for (i = 0; i < nbins; i++) source[i] = hdata->GetBinContent (i+1);
Analyser1dim->Background (
source, nbins, 20, TSpectrum::kBackDecreasingWindow,
TSpectrum::kBackOrder2, kFALSE,
TSpectrum::kBackSmoothing3, kFALSE );
TH1F* hwindow20 = new TH1F ("hwindow20", "", nbins, xmin, xmax);
for (i = 0; i < nbins; i++) hwindow20->SetBinContent (i+1, source[i]);
hwindow20->SetLineColor (kBlue);
hwindow20->Draw ("SAME L");
/*** Case for width of clipping window = 30 ***/
for (i = 0; i < nbins; i++) source[i] = hdata->GetBinContent (i+1);
Analyser1dim->Background (
source, nbins, 30, TSpectrum::kBackDecreasingWindow,
TSpectrum::kBackOrder2, kFALSE,
TSpectrum::kBackSmoothing3, kFALSE );
TH1F* hwindow30 = new TH1F ("hwindow30", "", nbins, xmin, xmax);
for (i = 0; i < nbins; i++) hwindow30->SetBinContent (i+1, source[i]);
hwindow30->SetLineColor (kGreen);
hwindow30->Draw ("SAME L");
/*** Case for width of clipping window = 40 ***/
for (i = 0; i < nbins; i++) source[i] = hdata->GetBinContent (i+1);
Analyser1dim->Background (
source, nbins, 40, TSpectrum::kBackDecreasingWindow,
TSpectrum::kBackOrder2, kFALSE,
TSpectrum::kBackSmoothing3, kFALSE );
TH1F* hwindow40 = new TH1F ("hwindow40", "", nbins, xmin, xmax);
for (i = 0; i < nbins; i++) hwindow40->SetBinContent (i+1, source[i]);
hwindow40->SetLineColor (kMagenta);
hwindow40->Draw ("SAME L");
gStyle->SetOptStat (0);
TLegend* legend = new TLegend (0.1, 0.7, 0.48, 0.9);
legend->AddEntry (hdata , "Data" , "l");
legend->AddEntry (hwindow10, "Width = 10", "l");
legend->AddEntry (hwindow20, "Width = 20", "l");
legend->AddEntry (hwindow30, "Width = 30", "l");
legend->AddEntry (hwindow40, "Width = 40", "l");
legend->Draw();
}