// @(#)root/spectrum:$Id$ // Author: Miroslav Morhac 27/05/99 #include "TSpectrum.h" #include "TPolyMarker.h" #include "TVirtualPad.h" #include "TList.h" #include "TH1.h" #include "TMath.h" //////////////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////////////// /** \class TSpectrum \ingroup Hist \brief Advanced Spectra Processing ## Advanced spectra processing This class contains advanced spectra processing functions for: - One-dimensional background estimation - One-dimensional smoothing - One-dimensional deconvolution - One-dimensional peak search Author: Miroslav Morhac Institute of Physics Slovak Academy of Sciences Dubravska cesta 9, 842 28 BRATISLAVA SLOVAKIA email:fyzimiro@savba.sk, fax:+421 7 54772479 The original code in C has been repackaged as a C++ class by R.Brun. The algorithms in this class have been published in the following references: 1. M.Morhac et al.: Background elimination methods for multidimensional coincidence gamma-ray spectra. Nuclear Instruments and Methods in Physics Research A 401 (1997) 113-132. 2. M.Morhac et al.: Efficient one- and two-dimensional Gold deconvolution and its application to gamma-ray spectra decomposition. Nuclear Instruments and Methods in Physics Research A 401 (1997) 385-408. 3. M.Morhac et al.: Identification of peaks in multidimensional coincidence gamma-ray spectra. Nuclear Instruments and Methods in Research Physics A 443(2000), 108-125. These NIM papers are also available as doc or ps files from: - [Spectrum.doc](ftp://root.cern.ch/root/Spectrum.doc) - [SpectrumDec.ps.gz](ftp://root.cern.ch/root/SpectrumDec.ps.gz) - [SpectrumSrc.ps.gz](ftp://root.cern.ch/root/SpectrumSrc.ps.gz) - [SpectrumBck.ps.gz](ftp://root.cern.ch/root/SpectrumBck.ps.gz) */ Int_t TSpectrum::fgIterations = 3; Int_t TSpectrum::fgAverageWindow = 3; #define PEAK_WINDOW 1024 ClassImp(TSpectrum) //////////////////////////////////////////////////////////////////////////////// TSpectrum::TSpectrum() :TNamed("Spectrum", "Miroslav Morhac peak finder") { /* Begin_Html Constructor. End_Html */ Int_t n = 100; fMaxPeaks = n; fPosition = new Double_t[n]; fPositionX = new Double_t[n]; fPositionY = new Double_t[n]; fResolution = 1; fHistogram = 0; fNPeaks = 0; } //////////////////////////////////////////////////////////////////////////////// TSpectrum::TSpectrum(Int_t maxpositions, Double_t resolution) :TNamed("Spectrum", "Miroslav Morhac peak finder") { /* Begin_Html End_Html */ Int_t n = maxpositions; if (n <= 0) n = 1; fMaxPeaks = n; fPosition = new Double_t[n]; fPositionX = new Double_t[n]; fPositionY = new Double_t[n]; fHistogram = 0; fNPeaks = 0; SetResolution(resolution); } //////////////////////////////////////////////////////////////////////////////// TSpectrum::~TSpectrum() { /* Begin_Html Destructor. End_Html */ delete [] fPosition; delete [] fPositionX; delete [] fPositionY; delete fHistogram; } //////////////////////////////////////////////////////////////////////////////// void TSpectrum::SetAverageWindow(Int_t w) { /* Begin_Html Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes). End_Html */ fgAverageWindow = w; } //////////////////////////////////////////////////////////////////////////////// void TSpectrum::SetDeconIterations(Int_t n) { /* Begin_Html Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::SearchHighRes). End_Html */ fgIterations = n; } //////////////////////////////////////////////////////////////////////////////// TH1 *TSpectrum::Background(const TH1 * h, int numberIterations, Option_t * option) { /* Begin_Html One-dimensional background estimation function.

This function calculates the background spectrum in the input histogram h. The background is returned as a histogram.

Function parameters:

NOTE that the background is only evaluated in the current range of h. ie, if h has a bin range (set via h->GetXaxis()->SetRange(binmin,binmax), the returned histogram will be created with the same number of bins as the input histogram h, but only bins from binmin to binmax will be filled with the estimated background. End_Html */ if (h == 0) return 0; Int_t dimension = h->GetDimension(); if (dimension > 1) { Error("Search", "Only implemented for 1-d histograms"); return 0; } TString opt = option; opt.ToLower(); //set options Int_t direction = kBackDecreasingWindow; if (opt.Contains("backincreasingwindow")) direction = kBackIncreasingWindow; Int_t filterOrder = kBackOrder2; if (opt.Contains("backorder4")) filterOrder = kBackOrder4; if (opt.Contains("backorder6")) filterOrder = kBackOrder6; if (opt.Contains("backorder8")) filterOrder = kBackOrder8; Bool_t smoothing = kTRUE; if (opt.Contains("nosmoothing")) smoothing = kFALSE; Int_t smoothWindow = kBackSmoothing3; if (opt.Contains("backsmoothing5")) smoothWindow = kBackSmoothing5; if (opt.Contains("backsmoothing7")) smoothWindow = kBackSmoothing7; if (opt.Contains("backsmoothing9")) smoothWindow = kBackSmoothing9; if (opt.Contains("backsmoothing11")) smoothWindow = kBackSmoothing11; if (opt.Contains("backsmoothing13")) smoothWindow = kBackSmoothing13; if (opt.Contains("backsmoothing15")) smoothWindow = kBackSmoothing15; Bool_t compton = kFALSE; if (opt.Contains("compton")) compton = kTRUE; Int_t first = h->GetXaxis()->GetFirst(); Int_t last = h->GetXaxis()->GetLast(); Int_t size = last-first+1; Int_t i; Double_t * source = new Double_t[size]; for (i = 0; i < size; i++) source[i] = h->GetBinContent(i + first); //find background (source is input and in output contains the background Background(source,size,numberIterations, direction, filterOrder,smoothing, smoothWindow,compton); //create output histogram containing backgound //only bins in the range of the input histogram are filled Int_t nch = strlen(h->GetName()); char *hbname = new char[nch+20]; snprintf(hbname,nch+20,"%s_background",h->GetName()); TH1 *hb = (TH1*)h->Clone(hbname); hb->Reset(); hb->GetListOfFunctions()->Delete(); hb->SetLineColor(2); for (i=0; i< size; i++) hb->SetBinContent(i+first,source[i]); hb->SetEntries(size); //if option "same is specified, draw the result in the pad if (opt.Contains("same")) { if (gPad) delete gPad->GetPrimitive(hbname); hb->Draw("same"); } delete [] source; delete [] hbname; return hb; } //////////////////////////////////////////////////////////////////////////////// void TSpectrum::Print(Option_t *) const { /* Begin_Html Print the array of positions. End_Html */ printf("\nNumber of positions = %d\n",fNPeaks); for (Int_t i=0;iOne-dimensional peak search function

This function searches for peaks in source spectrum in hin The number of found peaks and their positions are written into the members fNpeaks and fPositionX. The search is performed in the current histogram range.

Function parameters: