TSpectrum.cxx

Go to the documentation of this file.
00001 // @(#)root/spectrum:$Id: TSpectrum.cxx 35945 2010-09-30 16:05:07Z brun $
00002 // Author: Miroslav Morhac   27/05/99
00003 
00004 #include "TSpectrum.h"
00005 #include "TPolyMarker.h"
00006 #include "TVirtualPad.h"
00007 #include "TList.h"
00008 #include "TH1.h"
00009 #include "TMath.h"
00010 
00011 
00012 //______________________________________________________________________________
00013 /* Begin_Html
00014 <center><h2>Advanced Spectra Processing</h2></center>
00015 This class contains advanced spectra processing functions for:
00016 <ul>
00017 <li> One-dimensional background estimation
00018 <li> One-dimensional smoothing
00019 <li> One-dimensional deconvolution
00020 <li> One-dimensional peak search
00021 </ul>
00022 <p>
00023 Author:
00024 <br>
00025 <br> Miroslav Morhac
00026 <br> Institute of Physics
00027 <br> Slovak Academy of Sciences
00028 <br> Dubravska cesta 9, 842 28 BRATISLAVA
00029 <br> SLOVAKIA
00030 <br> email:fyzimiro@savba.sk,    fax:+421 7 54772479
00031 <br>
00032 <br>
00033 
00034 The original code in C has been repackaged as a C++ class by R.Brun.
00035 <p>
00036 The algorithms in this class have been published in the following references:
00037 <ol>
00038 <li> M.Morhac et al.: Background elimination methods for
00039 multidimensional coincidence gamma-ray spectra. Nuclear
00040 Instruments and Methods in Physics Research A 401 (1997) 113-132.
00041 <li> M.Morhac et al.: Efficient one- and two-dimensional Gold
00042 deconvolution and its application to gamma-ray spectra
00043 decomposition. Nuclear Instruments and Methods in Physics
00044 Research A 401 (1997) 385-408.
00045 <li> M.Morhac et al.: Identification of peaks in multidimensional
00046 coincidence gamma-ray spectra. Nuclear Instruments and Methods in
00047 Research Physics A  443(2000), 108-125.
00048 </ol>
00049 These NIM papers are also available as doc or ps files from:
00050 <ul>
00051 <li> <A href="ftp://root.cern.ch/root/Spectrum.doc">Spectrum.doc</A><br>
00052 <li> <A href="ftp://root.cern.ch/root/SpectrumDec.ps.gz">SpectrumDec.ps.gz</A><br>
00053 <li> <A href="ftp://root.cern.ch/root/SpectrumSrc.ps.gz">SpectrumSrc.ps.gz</A><br>
00054 <li> <A href="ftp://root.cern.ch/root/SpectrumBck.ps.gz">SpectrumBck.ps.gz</A><br>
00055 </ul>
00056 End_Html */
00057 
00058 Int_t TSpectrum::fgIterations    = 3;
00059 Int_t TSpectrum::fgAverageWindow = 3;
00060 
00061 #define PEAK_WINDOW 1024
00062 ClassImp(TSpectrum)
00063 
00064 
00065 //______________________________________________________________________________
00066 TSpectrum::TSpectrum() :TNamed("Spectrum", "Miroslav Morhac peak finder")
00067 {
00068    /* Begin_Html
00069    Constructor.
00070    End_Html */
00071 
00072    Int_t n = 100;
00073    fMaxPeaks  = n;
00074    fPosition   = new Float_t[n];
00075    fPositionX  = new Float_t[n];
00076    fPositionY  = new Float_t[n];
00077    fResolution = 1;
00078    fHistogram  = 0;
00079    fNPeaks     = 0;
00080 }
00081 
00082 
00083 //______________________________________________________________________________
00084 TSpectrum::TSpectrum(Int_t maxpositions, Float_t resolution)
00085           :TNamed("Spectrum", "Miroslav Morhac peak finder")
00086 {
00087    /* Begin_Html
00088    <ul>
00089    <li> maxpositions: maximum number of peaks
00090    <li> resolution: determines resolution of the neighboring peaks
00091                    default value is 1 correspond to 3 sigma distance
00092                    between peaks. Higher values allow higher resolution
00093                    (smaller distance between peaks.
00094                    May be set later through SetResolution.
00095    </ul>
00096    End_Html */
00097 
00098    Int_t n = maxpositions;
00099    if (n <= 0) n = 1;
00100    fMaxPeaks  = n;
00101    fPosition  = new Float_t[n];
00102    fPositionX = new Float_t[n];
00103    fPositionY = new Float_t[n];
00104    fHistogram = 0;
00105    fNPeaks    = 0;
00106    SetResolution(resolution);
00107 }
00108 
00109 
00110 //______________________________________________________________________________
00111 TSpectrum::~TSpectrum()
00112 {
00113    /* Begin_Html
00114    Destructor.
00115    End_Html */
00116 
00117    delete [] fPosition;
00118    delete [] fPositionX;
00119    delete [] fPositionY;
00120    delete    fHistogram;
00121 }
00122 
00123 
00124 //______________________________________________________________________________
00125 void TSpectrum::SetAverageWindow(Int_t w)
00126 {
00127    /* Begin_Html
00128    Static function: Set average window of searched peaks
00129    (see TSpectrum::SearchHighRes).
00130    End_Html */
00131 
00132    fgAverageWindow = w;
00133 }
00134 
00135 
00136 //______________________________________________________________________________
00137 void TSpectrum::SetDeconIterations(Int_t n)
00138 {
00139    /* Begin_Html
00140    Static function: Set max number of decon iterations in deconvolution
00141    operation (see TSpectrum::SearchHighRes).
00142    End_Html */
00143 
00144    fgIterations = n;
00145 }
00146 
00147 
00148 //______________________________________________________________________________
00149 TH1 *TSpectrum::Background(const TH1 * h, int numberIterations,
00150                            Option_t * option)
00151 {
00152    /* Begin_Html
00153    <b>One-dimensional background estimation function.</b>
00154    <p>
00155    This function calculates the background spectrum in the input histogram h.
00156    The background is returned as a histogram.
00157    <p>
00158    Function parameters:
00159    <ul>
00160    <li> h: input 1-d histogram
00161    <li> numberIterations, (default value = 20)
00162       Increasing numberIterations make the result smoother and lower.
00163    <li> option: may contain one of the following options:
00164       <ul>
00165       <li> to set the direction parameter
00166       "BackIncreasingWindow". By default the direction is BackDecreasingWindow
00167       <li> filterOrder-order of clipping filter,  (default "BackOrder2")
00168                   -possible values= "BackOrder4"
00169                                     "BackOrder6"
00170                                     "BackOrder8"
00171       <li> "nosmoothing"- if selected, the background is not smoothed
00172            By default the background is smoothed.
00173       <li> smoothWindow-width of smoothing window, (default is "BackSmoothing3")
00174                   -possible values= "BackSmoothing5"
00175                                     "BackSmoothing7"
00176                                     "BackSmoothing9"
00177                                     "BackSmoothing11"
00178                                     "BackSmoothing13"
00179                                     "BackSmoothing15"
00180       <li> "Compton" if selected the estimation of Compton edge
00181                   will be included.
00182       <li> "same" : if this option is specified, the resulting background
00183                  histogram is superimposed on the picture in the current pad.
00184       </ul>
00185    </ul>
00186    NOTE that the background is only evaluated in the current range of h.
00187    ie, if h has a bin range (set via h->GetXaxis()->SetRange(binmin,binmax),
00188    the returned histogram will be created with the same number of bins
00189    as the input histogram h, but only bins from binmin to binmax will be filled
00190    with the estimated background.
00191    End_Html */
00192 
00193    if (h == 0) return 0;
00194    Int_t dimension = h->GetDimension();
00195    if (dimension > 1) {
00196       Error("Search", "Only implemented for 1-d histograms");
00197       return 0;
00198    }
00199    TString opt = option;
00200    opt.ToLower();
00201 
00202    //set options
00203    Int_t direction = kBackDecreasingWindow;
00204    if (opt.Contains("backincreasingwindow")) direction = kBackIncreasingWindow;
00205    Int_t filterOrder = kBackOrder2;
00206    if (opt.Contains("backorder4")) filterOrder = kBackOrder4;
00207    if (opt.Contains("backorder6")) filterOrder = kBackOrder6;
00208    if (opt.Contains("backorder8")) filterOrder = kBackOrder8;
00209    Bool_t smoothing = kTRUE;
00210    if (opt.Contains("nosmoothing")) smoothing = kFALSE;
00211    Int_t smoothWindow = kBackSmoothing3;
00212    if (opt.Contains("backsmoothing5"))  smoothWindow = kBackSmoothing5;
00213    if (opt.Contains("backsmoothing7"))  smoothWindow = kBackSmoothing7;
00214    if (opt.Contains("backsmoothing9"))  smoothWindow = kBackSmoothing9;
00215    if (opt.Contains("backsmoothing11")) smoothWindow = kBackSmoothing11;
00216    if (opt.Contains("backsmoothing13")) smoothWindow = kBackSmoothing13;
00217    if (opt.Contains("backsmoothing15")) smoothWindow = kBackSmoothing15;
00218    Bool_t compton = kFALSE;
00219    if (opt.Contains("compton")) compton = kTRUE;
00220 
00221    Int_t first = h->GetXaxis()->GetFirst();
00222    Int_t last  = h->GetXaxis()->GetLast();
00223    Int_t size = last-first+1;
00224    Int_t i;
00225    Float_t * source = new float[size];
00226    for (i = 0; i < size; i++) source[i] = h->GetBinContent(i + first);
00227 
00228    //find background (source is input and in output contains the background
00229    Background(source,size,numberIterations, direction, filterOrder,smoothing,
00230               smoothWindow,compton);
00231 
00232    //create output histogram containing backgound
00233    //only bins in the range of the input histogram are filled
00234    Int_t nch = strlen(h->GetName());
00235    char *hbname = new char[nch+20];
00236    snprintf(hbname,nch+20,"%s_background",h->GetName());
00237    TH1 *hb = (TH1*)h->Clone(hbname);
00238    hb->Reset();
00239    hb->GetListOfFunctions()->Delete();
00240    hb->SetLineColor(2);
00241    for (i=0; i< size; i++) hb->SetBinContent(i+first,source[i]);
00242    hb->SetEntries(size);
00243 
00244    //if option "same is specified, draw the result in the pad
00245    if (opt.Contains("same")) {
00246       if (gPad) delete gPad->GetPrimitive(hbname);
00247       hb->Draw("same");
00248    }
00249    delete [] source;
00250    delete [] hbname;
00251    return hb;
00252 }
00253 
00254 
00255 //______________________________________________________________________________
00256 void TSpectrum::Print(Option_t *) const
00257 {
00258    /* Begin_Html
00259    Print the array of positions.
00260    End_Html */
00261 
00262    printf("\nNumber of positions = %d\n",fNPeaks);
00263    for (Int_t i=0;i<fNPeaks;i++) {
00264       printf(" x[%d] = %g, y[%d] = %g\n",i,fPositionX[i],i,fPositionY[i]);
00265    }
00266 }
00267 
00268 
00269 //______________________________________________________________________________
00270 Int_t TSpectrum::Search(const TH1 * hin, Double_t sigma, Option_t * option,
00271                         Double_t threshold)
00272 {
00273    /* Begin_Html
00274    <b>One-dimensional peak search function</b>
00275    <p>
00276    This function searches for peaks in source spectrum in hin
00277    The number of found peaks and their positions are written into
00278    the members fNpeaks and fPositionX.
00279    The search is performed in the current histogram range.
00280    <p>
00281    Function parameters:
00282    <ul>
00283    <li> hin:       pointer to the histogram of source spectrum
00284    <li> sigma:   sigma of searched peaks, for details we refer to manual
00285    <li> threshold: (default=0.05)  peaks with amplitude less than
00286        threshold*highest_peak are discarded.  0<threshold<1
00287    </ul>
00288    By default, the background is removed before deconvolution.
00289    Specify the option "nobackground" to not remove the background.
00290    <p>
00291    By default the "Markov" chain algorithm is used.
00292    Specify the option "noMarkov" to disable this algorithm
00293    Note that by default the source spectrum is replaced by a new spectrum
00294    <p>
00295    By default a polymarker object is created and added to the list of
00296    functions of the histogram. The histogram is drawn with the specified
00297    option and the polymarker object drawn on top of the histogram.
00298    The polymarker coordinates correspond to the npeaks peaks found in
00299    the histogram.
00300    <p>
00301    A pointer to the polymarker object can be retrieved later via:
00302    <pre>
00303     TList *functions = hin->GetListOfFunctions();
00304     TPolyMarker *pm = (TPolyMarker*)functions->FindObject("TPolyMarker");
00305    </pre>
00306    Specify the option "goff" to disable the storage and drawing of the
00307    polymarker.
00308    End_Html */
00309 
00310    if (hin == 0) return 0;
00311    Int_t dimension = hin->GetDimension();
00312    if (dimension > 2) {
00313       Error("Search", "Only implemented for 1-d and 2-d histograms");
00314       return 0;
00315    }
00316    if (threshold <=0 || threshold >= 1) {
00317       Warning("Search","threshold must 0<threshold<1, threshol=0.05 assumed");
00318       threshold = 0.05;
00319    }
00320    TString opt = option;
00321    opt.ToLower();
00322    Bool_t background = kTRUE;
00323    if (opt.Contains("nobackground")) {
00324       background = kFALSE;
00325       opt.ReplaceAll("nobackground","");
00326    }
00327    Bool_t markov = kTRUE;
00328    if (opt.Contains("nomarkov")) {
00329       markov = kFALSE;
00330       opt.ReplaceAll("nomarkov","");
00331    }
00332    if (dimension == 1) {
00333       Int_t first = hin->GetXaxis()->GetFirst();
00334       Int_t last  = hin->GetXaxis()->GetLast();
00335       Int_t size = last-first+1;
00336       Int_t i, bin, npeaks;
00337       Float_t * source = new float[size];
00338       Float_t * dest   = new float[size];
00339       for (i = 0; i < size; i++) source[i] = hin->GetBinContent(i + first);
00340       if (sigma < 1) {
00341          sigma = size/fMaxPeaks;
00342          if (sigma < 1) sigma = 1;
00343          if (sigma > 8) sigma = 8;
00344       }
00345       npeaks = SearchHighRes(source, dest, size, sigma, 100*threshold,
00346                              background, fgIterations, markov, fgAverageWindow);
00347 
00348       for (i = 0; i < npeaks; i++) {
00349          bin = first + Int_t(fPositionX[i] + 0.5);
00350          fPositionX[i] = hin->GetBinCenter(bin);
00351          fPositionY[i] = hin->GetBinContent(bin);
00352       }
00353       delete [] source;
00354       delete [] dest;
00355 
00356       if (opt.Contains("goff"))
00357          return npeaks;
00358       if (!npeaks) return 0;
00359       TPolyMarker * pm =
00360          (TPolyMarker*)hin->GetListOfFunctions()->FindObject("TPolyMarker");
00361       if (pm) {
00362          hin->GetListOfFunctions()->Remove(pm);
00363          delete pm;
00364       }
00365       pm = new TPolyMarker(npeaks, fPositionX, fPositionY);
00366       hin->GetListOfFunctions()->Add(pm);
00367       pm->SetMarkerStyle(23);
00368       pm->SetMarkerColor(kRed);
00369       pm->SetMarkerSize(1.3);
00370       opt.ReplaceAll(" ","");
00371       opt.ReplaceAll(",","");
00372       ((TH1*)hin)->Draw(opt.Data());
00373       return npeaks;
00374    }
00375    return 0;
00376 }
00377 
00378 
00379 //______________________________________________________________________________
00380 void TSpectrum::SetResolution(Float_t resolution)
00381 {
00382    /* Begin_Html
00383   resolution: determines resolution of the neighboring peaks
00384               default value is 1 correspond to 3 sigma distance
00385               between peaks. Higher values allow higher resolution
00386               (smaller distance between peaks.
00387               May be set later through SetResolution.
00388    End_Html */
00389 
00390    if (resolution > 1)
00391       fResolution = resolution;
00392    else
00393       fResolution = 1;
00394 }
00395 
00396 
00397 //______________________________________________________________________________
00398 const char *TSpectrum::Background(float *spectrum, int ssize,
00399                                           int numberIterations,
00400                                           int direction, int filterOrder,
00401                                           bool smoothing,int smoothWindow,
00402                                           bool compton)
00403 {
00404 /* Begin_Html
00405 This function calculates background spectrum from source spectrum.
00406 The result is placed in the vector pointed by spe1945ctrum pointer.
00407 The goal is to separate the useful information (peaks) from useless
00408 information (background).
00409 
00410 <ul>
00411 <li> method is based on Sensitive Nonlinear Iterative Peak (SNIP) clipping
00412      algorithm.
00413 <li> new value in the channel "i" is calculated
00414 </ul>
00415 
00416 <img width=486 height=72 src="gif/TSpectrum_Background.gif">
00417 
00418 where p = 1, 2, ..., numberIterations. In fact it represents second order
00419 difference filter (-1,2,-1).
00420 
00421 One can also change the
00422 direction of the change of the clipping window, the order of the clipping
00423 filter, to include smoothing, to set width of smoothing window and to include
00424 the estimation of Compton edges. On successful completion it returns 0. On
00425 error it returns pointer to the string describing error.
00426 
00427 <h4>Parameters:</h4>
00428 <ul>
00429 <li> spectrum: pointer to the vector of source spectrum
00430 <li> ssize: length of the spectrum vector
00431 <li> numberIterations: maximal width of clipping window,
00432 <li> direction:  direction of change of clipping window.
00433      Possible values: kBackIncreasingWindow, kBackDecreasingWindow
00434 <li> filterOrder: order of clipping filter.
00435      Possible values: kBackOrder2, kBackOrder4, kBackOrder6, kBackOrder8
00436 <li> smoothing: logical variable whether the smoothing operation in the
00437      estimation of background will be included.
00438      Possible values: kFALSE, kTRUE
00439 <li> smoothWindow: width of smoothing window.
00440      Possible values: kBackSmoothing3, kBackSmoothing5, kBackSmoothing7,
00441      kBackSmoothing9, kBackSmoothing11, kBackSmoothing13, kBackSmoothing15.
00442 <li> compton: logical variable whether the estimation of Compton edge will be
00443      included. Possible values: kFALSE, kTRUE.
00444 </ul>
00445 
00446 
00447 <h4>References:</h4>
00448 <ol>
00449 <li> C. G Ryan et al.: SNIP, a statistics-sensitive background treatment for the
00450 quantitative analysis of PIXE spectra in geoscience applications. NIM, B34
00451 (1988), 396-402.
00452 
00453 <li> M. Morhá&#269;, J. Kliman, V. Matoušek, M. Veselský, I. Turzo:
00454 Background elimination methods for multidimensional gamma-ray spectra. NIM,
00455 A401 (1997) 113-132.
00456 
00457 <li> D. D. Burgess, R. J. Tervo: Background estimation for gamma-ray
00458 spectroscopy. NIM 214 (1983), 431-434.
00459 </ol>
00460 
00461 Example 1 script Background_incr.c:
00462 <p>
00463 <img width=601 height=407 src="gif/TSpectrum_Background_incr.jpg">
00464 <p>
00465 Figure 1 Example of the estimation of background for number of iterations=6.
00466 Original spectrum is shown in black color, estimated background in red color.
00467 <p>
00468 Script:
00469 <pre>
00470 // Example to illustrate the background estimator (class TSpectrum).
00471 // To execute this example, do
00472 // root > .x Background_incr.C
00473 
00474 #include <TSpectrum>
00475 
00476 void Background_incr() {
00477    Int_t i;
00478    Double_t nbins = 256;
00479    Double_t xmin  = 0;
00480    Double_t xmax  = (Double_t)nbins;
00481    Float_t * source = new float[nbins];
00482    TH1F *back = new TH1F("back","",nbins,xmin,xmax);
00483    TH1F *d = new TH1F("d","",nbins,xmin,xmax);
00484    TFile *f = new TFile("spectra\\TSpectrum.root");
00485    back=(TH1F*) f->Get("back1;1");
00486    TCanvas *Background = gROOT->GetListOfCanvases()->FindObject("Background");
00487    if (!Background) Background =
00488      new TCanvas("Background",
00489                  "Estimation of background with increasing window",
00490                  10,10,1000,700);
00491    back->Draw("L");
00492    TSpectrum *s = new TSpectrum();
00493    for (i = 0; i < nbins; i++) source[i]=back->GetBinContent(i + 1);
00494    s->Background(source,nbins,6,kBackIncreasingWindow,kBackOrder2,kFALSE,
00495                  kBackSmoothing3,kFALSE);
00496    for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
00497    d->SetLineColor(kRed);
00498    d->Draw("SAME L");
00499 }
00500 </pre>
00501 
00502 Example 2 script Background_decr.c:
00503 <p>
00504 In Figure 1. one can notice that at the edges of the peaks the estimated
00505 background goes under the peaks. An alternative approach is to decrease the
00506 clipping window from a given value numberIterations to the value of one, which
00507 is presented in this example.
00508 <p>
00509 <img width=601 height=407 src="gif/TSpectrum_Background_decr.jpg">
00510 <p>
00511 Figure 2 Example of the estimation of background for numberIterations=6 using
00512 decreasing clipping window algorithm. Original spectrum is shown in black
00513 color, estimated background in red color.
00514 <p>
00515 Script:
00516 
00517 <pre>
00518 // Example to illustrate the background estimator (class TSpectrum).
00519 // To execute this example, do
00520 // root > .x Background_decr.C
00521 
00522 #include <TSpectrum>
00523 
00524 void Background_decr() {
00525    Int_t i;
00526    Double_t nbins = 256;
00527    Double_t xmin  = 0;
00528    Double_t xmax  = (Double_t)nbins;
00529    Float_t * source = new float[nbins];
00530    TH1F *back = new TH1F("back","",nbins,xmin,xmax);
00531    TH1F *d = new TH1F("d","",nbins,xmin,xmax);
00532    TFile *f = new TFile("spectra\\TSpectrum.root");
00533    back=(TH1F*) f->Get("back1;1");
00534    TCanvas *Background = gROOT->GetListOfCanvases()->FindObject("Background");
00535    if (!Background) Background =
00536      new TCanvas("Background","Estimation of background with decreasing window",
00537                  10,10,1000,700);
00538    back->Draw("L");
00539    TSpectrum *s = new TSpectrum();
00540    for (i = 0; i < nbins; i++) source[i]=back->GetBinContent(i + 1);
00541    s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kFALSE,
00542                  kBackSmoothing3,kFALSE);
00543    for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
00544    d->SetLineColor(kRed);
00545    d->Draw("SAME L");
00546 }
00547 </pre>
00548 
00549 Example 3 script Background_width.c:
00550 <p>
00551 The question is how to choose the width of the clipping window, i.e.,
00552 numberIterations parameter. The influence of this parameter on the estimated
00553 background is illustrated in Figure 3.
00554 <p>
00555 <img width=601 height=407 src="gif/TSpectrum_Background_width.jpg">
00556 <p>
00557 Figure 3 Example of the influence of clipping window width on the estimated
00558 background for numberIterations=4 (red line), 6 (blue line) 8 (green line) using
00559 decreasing clipping window algorithm.
00560 
00561 <p>
00562 in general one should set this parameter so that the value
00563 2*numberIterations+1 was greater than the widths of preserved objects (peaks).
00564 
00565 <p>
00566 Script:
00567 
00568 <pre>
00569 // Example to illustrate the influence of the clipping window width on the
00570 // estimated background. To execute this example, do:
00571 // root > .x Background_width.C
00572 
00573 #include <TSpectrum>
00574 
00575 void Background_width() {
00576    Int_t i;
00577    Double_t nbins = 256;
00578    Double_t xmin  = 0;
00579    Double_t xmax  = (Double_t)nbins;
00580    Float_t * source = new float[nbins];
00581    TH1F *h = new TH1F("h","",nbins,xmin,xmax);
00582    TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
00583    TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
00584    TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
00585    TFile *f = new TFile("spectra\\TSpectrum.root");
00586    h=(TH1F*) f->Get("back1;1");
00587    TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
00588    if (!background) background = new TCanvas("background",
00589    "Influence of clipping window width on the estimated background",
00590    10,10,1000,700);
00591    h->Draw("L");
00592    TSpectrum *s = new TSpectrum();
00593    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00594    s->Background(source,nbins,4,kBackDecreasingWindow,kBackOrder2,kFALSE,
00595    kBackSmoothing3,kFALSE);
00596    for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
00597    d1->SetLineColor(kRed);
00598    d1->Draw("SAME L");
00599    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00600    s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kFALSE,
00601    kBackSmoothing3,kFALSE);
00602    for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
00603    d2->SetLineColor(kBlue);
00604    d2->Draw("SAME L");
00605    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00606    s->Background(source,nbins,8,kBackDecreasingWindow,kBackOrder2,kFALSE,
00607    kBackSmoothing3,kFALSE);
00608    for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,source[i]);
00609    d3->SetLineColor(kGreen);
00610    d3->Draw("SAME L");
00611 }
00612 </pre>
00613 
00614 Example 4 script Background_width2.c:
00615 <p>
00616 another example for very complex spectrum is given in Figure 4.
00617 <p>
00618 <img width=601 height=407 src="gif/TSpectrum_Background_width2.jpg">
00619 <p>
00620 Figure 4 Example of the influence of clipping window width on the estimated
00621 background for numberIterations=10 (red line), 20 (blue line), 30 (green line)
00622 and 40 (magenta line) using decreasing clipping window algorithm.
00623 
00624 <p>
00625 Script:
00626 
00627 <pre>
00628 // Example to illustrate the influence of the clipping window width on the
00629 // estimated background. To execute this example, do:
00630 // root > .x Background_width2.C
00631 
00632 #include <TSpectrum>
00633 
00634 void Background_width2() {
00635    Int_t i;
00636    Double_t nbins = 4096;
00637    Double_t xmin  = 0;
00638    Double_t xmax  = (Double_t)4096;
00639    Float_t * source = new float[nbins];
00640    TH1F *h = new TH1F("h","",nbins,xmin,xmax);
00641    TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
00642    TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
00643    TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
00644    TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
00645    TFile *f = new TFile("spectra\\TSpectrum.root");
00646    h=(TH1F*) f->Get("back2;1");
00647    TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
00648    if (!background) background = new TCanvas("background",
00649    "Influence of clipping window width on the estimated background",
00650    10,10,1000,700);
00651    h->SetAxisRange(0,1000);
00652    h->SetMaximum(20000);
00653    h->Draw("L");
00654    TSpectrum *s = new TSpectrum();
00655    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00656    s->Background(source,nbins,10,kBackDecreasingWindow,kBackOrder2,kFALSE,
00657    kBackSmoothing3,kFALSE);
00658    for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
00659    d1->SetLineColor(kRed);
00660    d1->Draw("SAME L");
00661    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00662    s->Background(source,nbins,20,kBackDecreasingWindow,kBackOrder2,kFALSE,
00663    kBackSmoothing3,kFALSE);
00664    for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
00665    d2->SetLineColor(kBlue);
00666    d2->Draw("SAME L");
00667    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00668    s->Background(source,nbins,30,kBackDecreasingWindow,kBackOrder2,kFALSE,
00669    kBackSmoothing3,kFALSE);
00670    for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,source[i]);
00671    d3->SetLineColor(kGreen);
00672    d3->Draw("SAME L");
00673    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00674    s->Background(source,nbins,10,kBackDecreasingWindow,kBackOrder2,kFALSE,
00675    kBackSmoothing3,kFALSE);
00676    for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,source[i]);
00677    d4->SetLineColor(kMagenta);
00678    d4->Draw("SAME L");
00679 }
00680 </pre>
00681 
00682 Example 5 script Background_order.c:
00683 <p>
00684 Second order difference filter removes linear (quasi-linear) background and
00685 preserves symmetrical peaks. However if the shape of the background is more
00686 complex one can employ higher-order clipping filters (see example in Figure 5)
00687 <p>
00688 <img width=601 height=407 src="gif/TSpectrum_Background_order.jpg">
00689 <p>
00690 Figure 5 Example of the influence of clipping filter difference order on the
00691 estimated background for fNnumberIterations=40, 2-nd order red line, 4-th order
00692 blue line, 6-th order green line and 8-th order magenta line, and using
00693 decreasing clipping window algorithm.
00694 <p>
00695 Script:
00696 <pre>
00697 // Example to illustrate the influence of the clipping filter difference order
00698 // on the estimated background. To execute this example, do
00699 // root > .x Background_order.C
00700 
00701 #include <TSpectrum>
00702 
00703 void Background_order() {
00704    Int_t i;
00705    Double_t nbins = 4096;
00706    Double_t xmin  = 0;
00707    Double_t xmax  = (Double_t)4096;
00708    Float_t * source = new float[nbins];
00709    TH1F *h = new TH1F("h","",nbins,xmin,xmax);
00710    TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
00711    TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
00712    TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
00713    TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
00714    TFile *f = new TFile("spectra\\TSpectrum.root");
00715    h=(TH1F*) f->Get("back2;1");
00716    TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
00717    if (!background) background = new TCanvas("background",
00718    "Influence of clipping filter difference order on the estimated background",
00719    10,10,1000,700);
00720    h->SetAxisRange(1220,1460);
00721    h->SetMaximum(11000);
00722    h->Draw("L");
00723    TSpectrum *s = new TSpectrum();
00724    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00725    s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder2,kFALSE,
00726    kBackSmoothing3,kFALSE);
00727    for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
00728    d1->SetLineColor(kRed);
00729    d1->Draw("SAME L");
00730    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00731    s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder4,kFALSE,
00732    kBackSmoothing3,kFALSE);
00733    for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
00734    d2->SetLineColor(kBlue);
00735    d2->Draw("SAME L");
00736    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00737    s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder6,kFALSE,
00738    kBackSmoothing3,kFALSE);
00739    for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,source[i]);
00740    d3->SetLineColor(kGreen);
00741    d3->Draw("SAME L");
00742    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00743    s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder8,kFALSE,
00744    kBackSmoothing3,kFALSE);
00745    for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,source[i]);
00746    d4->SetLineColor(kMagenta);
00747    d4->Draw("SAME L");
00748 }
00749 </pre>
00750 
00751 Example 6 script Background_smooth.c:
00752 <p>
00753 The estimate of the background can be influenced by noise present in the
00754 spectrum.  We proposed the algorithm of the background estimate with
00755 simultaneous smoothing.  In the original algorithm without smoothing, the
00756 estimated background snatches the lower spikes in the noise. Consequently,
00757 the areas of peaks are biased by this error.
00758 <p>
00759 <img width=554 height=104 src="gif/TSpectrum_Background_smooth1.jpg">
00760 <p>
00761 Figure 7 Principle of background estimation algorithm with simultaneous
00762 smoothing.
00763 <p>
00764 <img width=601 height=407 src="gif/TSpectrum_Background_smooth2.jpg">
00765 <p>
00766 Figure 8 Illustration of non-smoothing (red line) and smoothing algorithm of
00767 background estimation (blue line).
00768 
00769 <p>
00770 
00771 Script:
00772 
00773 <pre>
00774 // Example to illustrate the background estimator (class TSpectrum) including
00775 // Compton edges. To execute this example, do:
00776 // root > .x Background_smooth.C
00777 
00778 #include <TSpectrum>
00779 
00780 void Background_smooth() {
00781    Int_t i;
00782    Double_t nbins = 4096;
00783    Double_t xmin  = 0;
00784    Double_t xmax  = (Double_t)nbins;
00785    Float_t * source = new float[nbins];
00786    TH1F *h = new TH1F("h","",nbins,xmin,xmax);
00787    TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
00788    TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
00789    TFile *f = new TFile("spectra\\TSpectrum.root");
00790    h=(TH1F*) f->Get("back4;1");
00791    TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
00792    if (!background) background = new TCanvas("background",
00793    "Estimation of background with noise",10,10,1000,700);
00794    h->SetAxisRange(3460,3830);
00795    h->Draw("L");
00796    TSpectrum *s = new TSpectrum();
00797    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00798    s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kFALSE,
00799    kBackSmoothing3,kFALSE);
00800    for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
00801    d1->SetLineColor(kRed);
00802    d1->Draw("SAME L");
00803    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00804    s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kTRUE,
00805    kBackSmoothing3,kFALSE);
00806    for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
00807    d2->SetLineColor(kBlue);
00808    d2->Draw("SAME L");
00809 }
00810 </pre>
00811 
00812 Example 8 script Background_compton.c:
00813 <p>
00814 Sometimes it is necessary to include also the Compton edges into the estimate of
00815 the background. In Figure 8 we present the example of the synthetic spectrum
00816 with Compton edges. The background was estimated using the 8-th order filter
00817 with the estimation of the Compton edges using decreasing
00818 clipping window algorithm (numberIterations=10) with smoothing
00819 (smoothingWindow=5).
00820 <p>
00821 <img width=601 height=407 src="gif/TSpectrum_Background_compton.jpg">
00822 <p>
00823 Figure 8 Example of the estimate of the background with Compton edges (red
00824 line) for numberIterations=10, 8-th order difference filter, using decreasing
00825 clipping window algorithm and smoothing (smoothingWindow=5).
00826 <p>
00827 Script:
00828 
00829 <pre>
00830 // Example to illustrate the background estimator (class TSpectrum) including
00831 // Compton edges. To execute this example, do:
00832 // root > .x Background_compton.C
00833 
00834 #include <TSpectrum>
00835 
00836 void Background_compton() {
00837    Int_t i;
00838    Double_t nbins = 512;
00839    Double_t xmin  = 0;
00840    Double_t xmax  = (Double_t)nbins;
00841    Float_t * source = new float[nbins];
00842    TH1F *h = new TH1F("h","",nbins,xmin,xmax);
00843    TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
00844    TFile *f = new TFile("spectra\\TSpectrum.root");
00845    h=(TH1F*) f->Get("back3;1");
00846    TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
00847    if (!background) background = new TCanvas("background",
00848    "Estimation of background with Compton edges under peaks",10,10,1000,700);
00849    h->Draw("L");
00850    TSpectrum *s = new TSpectrum();
00851    for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
00852    s->Background(source,nbins,10,kBackDecreasingWindow,kBackOrder8,kTRUE,
00853    kBackSmoothing5,,kTRUE);
00854    for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
00855    d1->SetLineColor(kRed);
00856    d1->Draw("SAME L");
00857 }
00858 </pre>
00859 
00860 End_Html */
00861 
00862    int i, j, w, bw, b1, b2, priz;
00863    float a, b, c, d, e, yb1, yb2, ai, av, men, b4, c4, d4, e4, b6, c6, d6, e6, f6, g6, b8, c8, d8, e8, f8, g8, h8, i8;
00864    if (ssize <= 0)
00865       return "Wrong Parameters";
00866    if (numberIterations < 1)
00867       return "Width of Clipping Window Must Be Positive";
00868    if (ssize < 2 * numberIterations + 1)
00869       return "Too Large Clipping Window";
00870    if (smoothing == kTRUE && smoothWindow != kBackSmoothing3 && smoothWindow != kBackSmoothing5 && smoothWindow != kBackSmoothing7 && smoothWindow != kBackSmoothing9 && smoothWindow != kBackSmoothing11 && smoothWindow != kBackSmoothing13 && smoothWindow != kBackSmoothing15)
00871       return "Incorrect width of smoothing window";
00872    float *working_space = new float[2 * ssize];
00873    for (i = 0; i < ssize; i++){
00874       working_space[i] = spectrum[i];
00875       working_space[i + ssize] = spectrum[i];
00876    }
00877    bw=(smoothWindow-1)/2;
00878    if (direction == kBackIncreasingWindow)
00879       i = 1;
00880    else if(direction == kBackDecreasingWindow)
00881       i = numberIterations;
00882    if (filterOrder == kBackOrder2) {
00883       do{
00884          for (j = i; j < ssize - i; j++) {
00885             if (smoothing == kFALSE){
00886                a = working_space[ssize + j];
00887                b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
00888                if (b < a)
00889                   a = b;
00890                working_space[j] = a;
00891             }
00892 
00893             else if (smoothing == kTRUE){
00894                a = working_space[ssize + j];
00895                av = 0;
00896                men = 0;
00897                for (w = j - bw; w <= j + bw; w++){
00898                   if ( w >= 0 && w < ssize){
00899                      av += working_space[ssize + w];
00900                      men +=1;
00901                   }
00902                }
00903                av = av / men;
00904                b = 0;
00905                men = 0;
00906                for (w = j - i - bw; w <= j - i + bw; w++){
00907                   if ( w >= 0 && w < ssize){
00908                      b += working_space[ssize + w];
00909                      men +=1;
00910                   }
00911                }
00912                b = b / men;
00913                c = 0;
00914                men = 0;
00915                for (w = j + i - bw; w <= j + i + bw; w++){
00916                   if ( w >= 0 && w < ssize){
00917                      c += working_space[ssize + w];
00918                      men +=1;
00919                   }
00920                }
00921                c = c / men;
00922                b = (b + c) / 2;
00923                if (b < a)
00924                   av = b;
00925                working_space[j]=av;
00926             }
00927          }
00928          for (j = i; j < ssize - i; j++)
00929             working_space[ssize + j] = working_space[j];
00930          if (direction == kBackIncreasingWindow)
00931             i+=1;
00932          else if(direction == kBackDecreasingWindow)
00933             i-=1;
00934       }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
00935    }
00936 
00937    else if (filterOrder == kBackOrder4) {
00938       do{
00939          for (j = i; j < ssize - i; j++) {
00940             if (smoothing == kFALSE){
00941                a = working_space[ssize + j];
00942                b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
00943                c = 0;
00944                ai = i / 2;
00945                c -= working_space[ssize + j - (int) (2 * ai)] / 6;
00946                c += 4 * working_space[ssize + j - (int) ai] / 6;
00947                c += 4 * working_space[ssize + j + (int) ai] / 6;
00948                c -= working_space[ssize + j + (int) (2 * ai)] / 6;
00949                if (b < c)
00950                   b = c;
00951                if (b < a)
00952                   a = b;
00953                working_space[j] = a;
00954             }
00955 
00956             else if (smoothing == kTRUE){
00957                a = working_space[ssize + j];
00958                av = 0;
00959                men = 0;
00960                for (w = j - bw; w <= j + bw; w++){
00961                   if ( w >= 0 && w < ssize){
00962                      av += working_space[ssize + w];
00963                      men +=1;
00964                   }
00965                }
00966                av = av / men;
00967                b = 0;
00968                men = 0;
00969                for (w = j - i - bw; w <= j - i + bw; w++){
00970                   if ( w >= 0 && w < ssize){
00971                      b += working_space[ssize + w];
00972                      men +=1;
00973                   }
00974                }
00975                b = b / men;
00976                c = 0;
00977                men = 0;
00978                for (w = j + i - bw; w <= j + i + bw; w++){
00979                   if ( w >= 0 && w < ssize){
00980                      c += working_space[ssize + w];
00981                      men +=1;
00982                   }
00983                }
00984                c = c / men;
00985                b = (b + c) / 2;
00986                ai = i / 2;
00987                b4 = 0, men = 0;
00988                for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
00989                   if (w >= 0 && w < ssize){
00990                      b4 += working_space[ssize + w];
00991                      men +=1;
00992                   }
00993                }
00994                b4 = b4 / men;
00995                c4 = 0, men = 0;
00996                for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
00997                   if (w >= 0 && w < ssize){
00998                      c4 += working_space[ssize + w];
00999                      men +=1;
01000                   }
01001                }
01002                c4 = c4 / men;
01003                d4 = 0, men = 0;
01004                for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
01005                   if (w >= 0 && w < ssize){
01006                      d4 += working_space[ssize + w];
01007                      men +=1;
01008                   }
01009                }
01010                d4 = d4 / men;
01011                e4 = 0, men = 0;
01012                for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
01013                   if (w >= 0 && w < ssize){
01014                      e4 += working_space[ssize + w];
01015                      men +=1;
01016                   }
01017                }
01018                e4 = e4 / men;
01019                b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
01020                if (b < b4)
01021                   b = b4;
01022                if (b < a)
01023                   av = b;
01024                working_space[j]=av;
01025             }
01026          }
01027          for (j = i; j < ssize - i; j++)
01028             working_space[ssize + j] = working_space[j];
01029          if (direction == kBackIncreasingWindow)
01030             i+=1;
01031          else if(direction == kBackDecreasingWindow)
01032             i-=1;
01033       }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
01034    }
01035 
01036    else if (filterOrder == kBackOrder6) {
01037       do{
01038          for (j = i; j < ssize - i; j++) {
01039             if (smoothing == kFALSE){
01040                a = working_space[ssize + j];
01041                b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
01042                c = 0;
01043                ai = i / 2;
01044                c -= working_space[ssize + j - (int) (2 * ai)] / 6;
01045                c += 4 * working_space[ssize + j - (int) ai] / 6;
01046                c += 4 * working_space[ssize + j + (int) ai] / 6;
01047                c -= working_space[ssize + j + (int) (2 * ai)] / 6;
01048                d = 0;
01049                ai = i / 3;
01050                d += working_space[ssize + j - (int) (3 * ai)] / 20;
01051                d -= 6 * working_space[ssize + j - (int) (2 * ai)] / 20;
01052                d += 15 * working_space[ssize + j - (int) ai] / 20;
01053                d += 15 * working_space[ssize + j + (int) ai] / 20;
01054                d -= 6 * working_space[ssize + j + (int) (2 * ai)] / 20;
01055                d += working_space[ssize + j + (int) (3 * ai)] / 20;
01056                if (b < d)
01057                   b = d;
01058                if (b < c)
01059                   b = c;
01060                if (b < a)
01061                   a = b;
01062                working_space[j] = a;
01063             }
01064 
01065             else if (smoothing == kTRUE){
01066                a = working_space[ssize + j];
01067                av = 0;
01068                men = 0;
01069                for (w = j - bw; w <= j + bw; w++){
01070                   if ( w >= 0 && w < ssize){
01071                      av += working_space[ssize + w];
01072                      men +=1;
01073                   }
01074                }
01075                av = av / men;
01076                b = 0;
01077                men = 0;
01078                for (w = j - i - bw; w <= j - i + bw; w++){
01079                   if ( w >= 0 && w < ssize){
01080                      b += working_space[ssize + w];
01081                      men +=1;
01082                   }
01083                }
01084                b = b / men;
01085                c = 0;
01086                men = 0;
01087                for (w = j + i - bw; w <= j + i + bw; w++){
01088                   if ( w >= 0 && w < ssize){
01089                      c += working_space[ssize + w];
01090                      men +=1;
01091                   }
01092                }
01093                c = c / men;
01094                b = (b + c) / 2;
01095                ai = i / 2;
01096                b4 = 0, men = 0;
01097                for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
01098                   if (w >= 0 && w < ssize){
01099                      b4 += working_space[ssize + w];
01100                      men +=1;
01101                   }
01102                }
01103                b4 = b4 / men;
01104                c4 = 0, men = 0;
01105                for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
01106                   if (w >= 0 && w < ssize){
01107                      c4 += working_space[ssize + w];
01108                      men +=1;
01109                   }
01110                }
01111                c4 = c4 / men;
01112                d4 = 0, men = 0;
01113                for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
01114                   if (w >= 0 && w < ssize){
01115                      d4 += working_space[ssize + w];
01116                      men +=1;
01117                   }
01118                }
01119                d4 = d4 / men;
01120                e4 = 0, men = 0;
01121                for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
01122                   if (w >= 0 && w < ssize){
01123                      e4 += working_space[ssize + w];
01124                      men +=1;
01125                   }
01126                }
01127                e4 = e4 / men;
01128                b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
01129                ai = i / 3;
01130                b6 = 0, men = 0;
01131                for (w = j - (int)(3 * ai) - bw; w <= j - (int)(3 * ai) + bw; w++){
01132                   if (w >= 0 && w < ssize){
01133                      b6 += working_space[ssize + w];
01134                      men +=1;
01135                   }
01136                }
01137                b6 = b6 / men;
01138                c6 = 0, men = 0;
01139                for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
01140                   if (w >= 0 && w < ssize){
01141                      c6 += working_space[ssize + w];
01142                      men +=1;
01143                   }
01144                }
01145                c6 = c6 / men;
01146                d6 = 0, men = 0;
01147                for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
01148                   if (w >= 0 && w < ssize){
01149                      d6 += working_space[ssize + w];
01150                      men +=1;
01151                   }
01152                }
01153                d6 = d6 / men;
01154                e6 = 0, men = 0;
01155                for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
01156                   if (w >= 0 && w < ssize){
01157                      e6 += working_space[ssize + w];
01158                      men +=1;
01159                   }
01160                }
01161                e6 = e6 / men;
01162                f6 = 0, men = 0;
01163                for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
01164                   if (w >= 0 && w < ssize){
01165                      f6 += working_space[ssize + w];
01166                      men +=1;
01167                   }
01168                }
01169                f6 = f6 / men;
01170                g6 = 0, men = 0;
01171                for (w = j + (int)(3 * ai) - bw; w <= j + (int)(3 * ai) + bw; w++){
01172                   if (w >= 0 && w < ssize){
01173                      g6 += working_space[ssize + w];
01174                      men +=1;
01175                   }
01176                }
01177                g6 = g6 / men;
01178                b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
01179                if (b < b6)
01180                   b = b6;
01181                if (b < b4)
01182                   b = b4;
01183                if (b < a)
01184                   av = b;
01185                working_space[j]=av;
01186             }
01187          }
01188          for (j = i; j < ssize - i; j++)
01189             working_space[ssize + j] = working_space[j];
01190          if (direction == kBackIncreasingWindow)
01191             i+=1;
01192          else if(direction == kBackDecreasingWindow)
01193             i-=1;
01194       }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
01195    }
01196 
01197    else if (filterOrder == kBackOrder8) {
01198       do{
01199          for (j = i; j < ssize - i; j++) {
01200             if (smoothing == kFALSE){
01201                a = working_space[ssize + j];
01202                b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
01203                c = 0;
01204                ai = i / 2;
01205                c -= working_space[ssize + j - (int) (2 * ai)] / 6;
01206                c += 4 * working_space[ssize + j - (int) ai] / 6;
01207                c += 4 * working_space[ssize + j + (int) ai] / 6;
01208                c -= working_space[ssize + j + (int) (2 * ai)] / 6;
01209                d = 0;
01210                ai = i / 3;
01211                d += working_space[ssize + j - (int) (3 * ai)] / 20;
01212                d -= 6 * working_space[ssize + j - (int) (2 * ai)] / 20;
01213                d += 15 * working_space[ssize + j - (int) ai] / 20;
01214                d += 15 * working_space[ssize + j + (int) ai] / 20;
01215                d -= 6 * working_space[ssize + j + (int) (2 * ai)] / 20;
01216                d += working_space[ssize + j + (int) (3 * ai)] / 20;
01217                e = 0;
01218                ai = i / 4;
01219                e -= working_space[ssize + j - (int) (4 * ai)] / 70;
01220                e += 8 * working_space[ssize + j - (int) (3 * ai)] / 70;
01221                e -= 28 * working_space[ssize + j - (int) (2 * ai)] / 70;
01222                e += 56 * working_space[ssize + j - (int) ai] / 70;
01223                e += 56 * working_space[ssize + j + (int) ai] / 70;
01224                e -= 28 * working_space[ssize + j + (int) (2 * ai)] / 70;
01225                e += 8 * working_space[ssize + j + (int) (3 * ai)] / 70;
01226                e -= working_space[ssize + j + (int) (4 * ai)] / 70;
01227                if (b < e)
01228                   b = e;
01229                if (b < d)
01230                   b = d;
01231                if (b < c)
01232                   b = c;
01233                if (b < a)
01234                   a = b;
01235                working_space[j] = a;
01236             }
01237 
01238             else if (smoothing == kTRUE){
01239                a = working_space[ssize + j];
01240                av = 0;
01241                men = 0;
01242                for (w = j - bw; w <= j + bw; w++){
01243                   if ( w >= 0 && w < ssize){
01244                      av += working_space[ssize + w];
01245                      men +=1;
01246                   }
01247                }
01248                av = av / men;
01249                b = 0;
01250                men = 0;
01251                for (w = j - i - bw; w <= j - i + bw; w++){
01252                   if ( w >= 0 && w < ssize){
01253                      b += working_space[ssize + w];
01254                      men +=1;
01255                   }
01256                }
01257                b = b / men;
01258                c = 0;
01259                men = 0;
01260                for (w = j + i - bw; w <= j + i + bw; w++){
01261                   if ( w >= 0 && w < ssize){
01262                      c += working_space[ssize + w];
01263                      men +=1;
01264                   }
01265                }
01266                c = c / men;
01267                b = (b + c) / 2;
01268                ai = i / 2;
01269                b4 = 0, men = 0;
01270                for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
01271                   if (w >= 0 && w < ssize){
01272                      b4 += working_space[ssize + w];
01273                      men +=1;
01274                   }
01275                }
01276                b4 = b4 / men;
01277                c4 = 0, men = 0;
01278                for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
01279                   if (w >= 0 && w < ssize){
01280                      c4 += working_space[ssize + w];
01281                      men +=1;
01282                   }
01283                }
01284                c4 = c4 / men;
01285                d4 = 0, men = 0;
01286                for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
01287                   if (w >= 0 && w < ssize){
01288                      d4 += working_space[ssize + w];
01289                      men +=1;
01290                   }
01291                }
01292                d4 = d4 / men;
01293                e4 = 0, men = 0;
01294                for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
01295                   if (w >= 0 && w < ssize){
01296                      e4 += working_space[ssize + w];
01297                      men +=1;
01298                   }
01299                }
01300                e4 = e4 / men;
01301                b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
01302                ai = i / 3;
01303                b6 = 0, men = 0;
01304                for (w = j - (int)(3 * ai) - bw; w <= j - (int)(3 * ai) + bw; w++){
01305                   if (w >= 0 && w < ssize){
01306                      b6 += working_space[ssize + w];
01307                      men +=1;
01308                   }
01309                }
01310                b6 = b6 / men;
01311                c6 = 0, men = 0;
01312                for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
01313                   if (w >= 0 && w < ssize){
01314                      c6 += working_space[ssize + w];
01315                      men +=1;
01316                   }
01317                }
01318                c6 = c6 / men;
01319                d6 = 0, men = 0;
01320                for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
01321                   if (w >= 0 && w < ssize){
01322                      d6 += working_space[ssize + w];
01323                      men +=1;
01324                   }
01325                }
01326                d6 = d6 / men;
01327                e6 = 0, men = 0;
01328                for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
01329                   if (w >= 0 && w < ssize){
01330                      e6 += working_space[ssize + w];
01331                      men +=1;
01332                   }
01333                }
01334                e6 = e6 / men;
01335                f6 = 0, men = 0;
01336                for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
01337                   if (w >= 0 && w < ssize){
01338                      f6 += working_space[ssize + w];
01339                      men +=1;
01340                   }
01341                }
01342                f6 = f6 / men;
01343                g6 = 0, men = 0;
01344                for (w = j + (int)(3 * ai) - bw; w <= j + (int)(3 * ai) + bw; w++){
01345                   if (w >= 0 && w < ssize){
01346                      g6 += working_space[ssize + w];
01347                      men +=1;
01348                   }
01349                }
01350                g6 = g6 / men;
01351                b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
01352                ai = i / 4;
01353                b8 = 0, men = 0;
01354                for (w = j - (int)(4 * ai) - bw; w <= j - (int)(4 * ai) + bw; w++){
01355                   if (w >= 0 && w < ssize){
01356                      b8 += working_space[ssize + w];
01357                      men +=1;
01358                   }
01359                }
01360                b8 = b8 / men;
01361                c8 = 0, men = 0;
01362                for (w = j - (int)(3 * ai) - bw; w <= j - (int)(3 * ai) + bw; w++){
01363                   if (w >= 0 && w < ssize){
01364                      c8 += working_space[ssize + w];
01365                      men +=1;
01366                   }
01367                }
01368                c8 = c8 / men;
01369                d8 = 0, men = 0;
01370                for (w = j - (int)(2 * ai) - bw; w <= j - (int)(2 * ai) + bw; w++){
01371                   if (w >= 0 && w < ssize){
01372                      d8 += working_space[ssize + w];
01373                      men +=1;
01374                   }
01375                }
01376                d8 = d8 / men;
01377                e8 = 0, men = 0;
01378                for (w = j - (int)ai - bw; w <= j - (int)ai + bw; w++){
01379                   if (w >= 0 && w < ssize){
01380                      e8 += working_space[ssize + w];
01381                      men +=1;
01382                   }
01383                }
01384                e8 = e8 / men;
01385                f8 = 0, men = 0;
01386                for (w = j + (int)ai - bw; w <= j + (int)ai + bw; w++){
01387                   if (w >= 0 && w < ssize){
01388                      f8 += working_space[ssize + w];
01389                      men +=1;
01390                   }
01391                }
01392                f8 = f8 / men;
01393                g8 = 0, men = 0;
01394                for (w = j + (int)(2 * ai) - bw; w <= j + (int)(2 * ai) + bw; w++){
01395                   if (w >= 0 && w < ssize){
01396                      g8 += working_space[ssize + w];
01397                      men +=1;
01398                   }
01399                }
01400                g8 = g8 / men;
01401                h8 = 0, men = 0;
01402                for (w = j + (int)(3 * ai) - bw; w <= j + (int)(3 * ai) + bw; w++){
01403                   if (w >= 0 && w < ssize){
01404                      h8 += working_space[ssize + w];
01405                      men +=1;
01406                   }
01407                }
01408                h8 = h8 / men;
01409                i8 = 0, men = 0;
01410                for (w = j + (int)(4 * ai) - bw; w <= j + (int)(4 * ai) + bw; w++){
01411                   if (w >= 0 && w < ssize){
01412                      i8 += working_space[ssize + w];
01413                      men +=1;
01414                   }
01415                }
01416                i8 = i8 / men;
01417                b8 = ( -b8 + 8 * c8 - 28 * d8 + 56 * e8 - 56 * f8 - 28 * g8 + 8 * h8 - i8)/70;
01418                if (b < b8)
01419                   b = b8;
01420                if (b < b6)
01421                   b = b6;
01422                if (b < b4)
01423                   b = b4;
01424                if (b < a)
01425                   av = b;
01426                working_space[j]=av;
01427             }
01428          }
01429          for (j = i; j < ssize - i; j++)
01430             working_space[ssize + j] = working_space[j];
01431          if (direction == kBackIncreasingWindow)
01432             i += 1;
01433          else if(direction == kBackDecreasingWindow)
01434             i -= 1;
01435       }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
01436    }
01437 
01438    if (compton == kTRUE) {
01439       for (i = 0, b2 = 0; i < ssize; i++){
01440          b1 = b2;
01441          a = working_space[i], b = spectrum[i];
01442          j = i;
01443          if (TMath::Abs(a - b) >= 1) {
01444             b1 = i - 1;
01445             if (b1 < 0)
01446                b1 = 0;
01447             yb1 = working_space[b1];
01448             for (b2 = b1 + 1, c = 0, priz = 0; priz == 0 && b2 < ssize; b2++){
01449                a = working_space[b2], b = spectrum[b2];
01450                c = c + b - yb1;
01451                if (TMath::Abs(a - b) < 1) {
01452                   priz = 1;
01453                   yb2 = b;
01454                }
01455             }
01456             if (b2 == ssize)
01457                b2 -= 1;
01458             yb2 = working_space[b2];
01459             if (yb1 <= yb2){
01460                for (j = b1, c = 0; j <= b2; j++){
01461                   b = spectrum[j];
01462                   c = c + b - yb1;
01463                }
01464                if (c > 1){
01465                   c = (yb2 - yb1) / c;
01466                   for (j = b1, d = 0; j <= b2 && j < ssize; j++){
01467                      b = spectrum[j];
01468                      d = d + b - yb1;
01469                      a = c * d + yb1;
01470                      working_space[ssize + j] = a;
01471                   }
01472                }
01473             }
01474 
01475             else{
01476                for (j = b2, c = 0; j >= b1; j--){
01477                   b = spectrum[j];
01478                   c = c + b - yb2;
01479                }
01480                if (c > 1){
01481                   c = (yb1 - yb2) / c;
01482                   for (j = b2, d = 0;j >= b1 && j >= 0; j--){
01483                      b = spectrum[j];
01484                      d = d + b - yb2;
01485                      a = c * d + yb2;
01486                      working_space[ssize + j] = a;
01487                   }
01488                }
01489             }
01490             i=b2;
01491          }
01492       }
01493    }
01494 
01495    for (j = 0; j < ssize; j++)
01496       spectrum[j] = working_space[ssize + j];
01497    delete[]working_space;
01498    return 0;
01499 }
01500 
01501 
01502 //______________________________________________________________________________
01503 const char* TSpectrum::SmoothMarkov(float *source, int ssize, int averWindow)
01504 {
01505    /* Begin_Html
01506    <b>One-dimensional markov spectrum smoothing function</b>
01507    <p>
01508    This function calculates smoothed spectrum from source spectrum based on
01509    Markov chain method. The result is placed in the array pointed by source
01510    pointer. On successful completion it returns 0. On error it returns pointer
01511    to the string describing error.
01512    <p>
01513    Function parameters:
01514    <ul>
01515    <li> source: pointer to the array of source spectrum
01516    <li> ssize: length of source array
01517    <li> averWindow: width of averaging smoothing window
01518    </ul>
01519    The goal of this function is the suppression of the statistical fluctuations.
01520    The algorithm is based on discrete Markov chain, which has very simple
01521    invariant distribution:
01522    <img width=551 height=63 src="gif/TSpectrum_Smoothing1.gif">
01523    <p>
01524    <img width=28 height=36 src="gif/TSpectrum_Smoothing2.gif"> being defined
01525    from the normalization condition
01526    <img width=70 height=52 src="gif/TSpectrum_Smoothing3.gif">.
01527    n is the length of the smoothed spectrum and
01528    <img width=258 height=60 src="gif/TSpectrum_Smoothing4.gif">
01529    <p>
01530    Reference:
01531    <ol>
01532    <li> Z.K. Silagadze, A new algorithm for automatic photopeak searches.
01533    NIM A 376 (1996), 451.
01534    </ol>
01535    <p>
01536    Example 14 - script Smoothing.c
01537    <p>
01538    <img width=296 height=182 src="gif/TSpectrum_Smoothing1.jpg">
01539    Fig. 23 Original noisy spectrum
01540    <p>
01541    <img width=296 height=182 src="gif/TSpectrum_Smoothing2.jpg">
01542    Fig. 24 Smoothed spectrum m=3
01543    <p>
01544    <img width=299 height=184 src="gif/TSpectrum_Smoothing3.jpg">
01545    Fig. 25 Smoothed spectrum
01546    <p>
01547    <img width=299 height=184 src="gif/TSpectrum_Smoothing4.jpg">
01548    Fig.26 Smoothed spectrum m=10
01549    <p>
01550    Script:
01551    <pre>
01552    // Example to illustrate smoothing using Markov algorithm (class TSpectrum).
01553    // To execute this example, do
01554    // root > .x Smoothing.C
01555 
01556    void Smoothing() {
01557       Int_t i;
01558       Double_t nbins = 1024;
01559       Double_t xmin  = 0;
01560       Double_t xmax  = (Double_t)nbins;
01561       Float_t * source = new float[nbins];
01562       TH1F *h = new TH1F("h","Smoothed spectrum for m=3",nbins,xmin,xmax);
01563       TFile *f = new TFile("spectra\\TSpectrum.root");
01564       h=(TH1F*) f->Get("smooth1;1");
01565       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
01566       TCanvas *Smooth1 = gROOT->GetListOfCanvases()->FindObject("Smooth1");
01567       if (!Smooth1) Smooth1 = new TCanvas("Smooth1","Smooth1",10,10,1000,700);
01568       TSpectrum *s = new TSpectrum();
01569       s->SmoothMarkov(source,1024,3);  //3, 7, 10
01570       for (i = 0; i < nbins; i++) h->SetBinContent(i + 1,source[i]);
01571       h->SetAxisRange(330,880);
01572       h->Draw("L");
01573    }
01574    </pre>
01575    End_Html */
01576 
01577    int xmin, xmax, i, l;
01578    float a, b, maxch;
01579    float nom, nip, nim, sp, sm, area = 0;
01580    if(averWindow <= 0)
01581       return "Averaging Window must be positive";
01582    float *working_space = new float[ssize];
01583    xmin = 0,xmax = ssize - 1;
01584    for(i = 0, maxch = 0; i < ssize; i++){
01585       working_space[i]=0;
01586       if(maxch < source[i])
01587          maxch = source[i];
01588 
01589       area += source[i];
01590    }
01591    if(maxch == 0) {
01592       delete [] working_space;
01593       return 0 ;
01594    }
01595 
01596    nom = 1;
01597    working_space[xmin] = 1;
01598    for(i = xmin; i < xmax; i++){
01599       nip = source[i] / maxch;
01600       nim = source[i + 1] / maxch;
01601       sp = 0,sm = 0;
01602       for(l = 1; l <= averWindow; l++){
01603          if((i + l) > xmax)
01604             a = source[xmax] / maxch;
01605 
01606          else
01607             a = source[i + l] / maxch;
01608          b = a - nip;
01609          if(a + nip <= 0)
01610             a = 1;
01611 
01612          else
01613             a = TMath::Sqrt(a + nip);
01614          b = b / a;
01615          b = TMath::Exp(b);
01616          sp = sp + b;
01617          if((i - l + 1) < xmin)
01618             a = source[xmin] / maxch;
01619 
01620          else
01621             a = source[i - l + 1] / maxch;
01622          b = a - nim;
01623          if(a + nim <= 0)
01624             a = 1;
01625          else
01626             a = TMath::Sqrt(a + nim);
01627          b = b / a;
01628          b = TMath::Exp(b);
01629          sm = sm + b;
01630       }
01631       a = sp / sm;
01632       a = working_space[i + 1] = working_space[i] * a;
01633       nom = nom + a;
01634    }
01635    for(i = xmin; i <= xmax; i++){
01636       working_space[i] = working_space[i] / nom;
01637    }
01638    for(i = 0; i < ssize; i++)
01639       source[i] = working_space[i] * area;
01640    delete[]working_space;
01641    return 0;
01642 }
01643 
01644 
01645 //______________________________________________________________________________
01646 const char *TSpectrum::Deconvolution(float *source, const float *response,
01647                                       int ssize, int numberIterations,
01648                                       int numberRepetitions, double boost )
01649 {
01650    /* Begin_Html
01651    <b>One-dimensional deconvolution function</b>
01652    <p>
01653    This function calculates deconvolution from source spectrum according to
01654    response spectrum using Gold deconvolution algorithm. The result is placed
01655    in the vector pointed by source pointer. On successful completion it
01656    returns 0. On error it returns pointer to the string describing error. If
01657    desired after every numberIterations one can apply boosting operation
01658    (exponential function with exponent given by boost coefficient) and repeat
01659    it numberRepetitions times.
01660    <p>
01661    Function parameters:
01662    <ul>
01663    <li>source:  pointer to the vector of source spectrum
01664    <li>response:     pointer to the vector of response spectrum
01665    <li>ssize:    length of source and response spectra
01666    numberIterations, for details we refer to the reference given below
01667    numberRepetitions, for repeated boosted deconvolution
01668    boost, boosting coefficient
01669    </ul>
01670    The goal of this function is the improvement of the resolution in spectra,
01671    decomposition of multiplets. The mathematical formulation of
01672    the convolution system is:
01673    <p>
01674    <img width=585 height=84 src="gif/TSpectrum_Deconvolution1.gif">
01675    <p>
01676    where h(i) is the impulse response function, x, y are input and output
01677    vectors, respectively, N is the length of x and h vectors. In matrix form
01678    we have:
01679    <p>
01680    <img width=597 height=360 src="gif/TSpectrum_Deconvolution2.gif">
01681    <p>
01682    Let us assume that we know the response and the output vector (spectrum) of
01683    the above given system. The deconvolution represents solution of the
01684    overdetermined system of linear equations, i.e., the calculation of the
01685    vector <b>x</b>. From numerical stability point of view the operation of
01686    deconvolution is extremely critical (ill-posed problem) as well as time
01687    consuming operation. The Gold deconvolution algorithm proves to work very
01688    well, other methods (Fourier, VanCittert etc) oscillate. It is suitable to
01689    process positive definite data (e.g. histograms).
01690    <p>
01691    <b>Gold deconvolution algorithm:</b>
01692    <p>
01693    <img width=551 height=233 src="gif/TSpectrum_Deconvolution3.gif">
01694    <p>
01695    Where L is given number of iterations (numberIterations parameter).
01696    <p>
01697    <b>Boosted deconvolution:</b>
01698    <ol>
01699    <li> Set the initial solution:
01700         End_Html Begin_Latex x^{(0)} = [1,1,...,1]^{T} End_Latex Begin_Html
01701    <li> Set required number of repetitions R and iterations L.
01702    <li> Set r = 1.
01703    <li>Using Gold deconvolution algorithm for k=1,2,...,L find
01704        End_Html Begin_Latex x^{(L)} End_Latex Begin_Html
01705    <li> If r = R stop calculation, else
01706       <ol>
01707       <li> Apply boosting operation, i.e., set
01708            End_Html Begin_Latex x^{(0)}(i) = [x^{(L)}(i)]^{p} End_Latex Begin_Html
01709            i=0,1,...N-1 and p is boosting coefficient &gt;0.
01710       <li> r = r + 1
01711       <li> continue in 4.
01712       </ol>
01713    </ol>
01714    <p>
01715    <b>References:</b>
01716    <ol>
01717    <li> Gold R., ANL-6984, Argonne National Laboratories, Argonne Ill, 1964.
01718    <li> Coote G.E., Iterative smoothing and deconvolution of one- and two-dimensional
01719         elemental distribution data, NIM B 130 (1997) 118.
01720    <li> M. Morhá&#269;, J. Kliman, V.  Matoušek, M. Veselský,
01721         I. Turzo: Efficient one- and two-dimensional Gold deconvolution and
01722         its application to gamma-ray spectra decomposition. NIM, A401 (1997) 385-408.
01723    <li> Morhá&#269; M., Matoušek V., Kliman J., Efficient algorithm of multidimensional
01724         deconvolution and its application to nuclear data processing, Digital Signal
01725         Processing 13 (2003) 144.
01726    </ol>
01727    <p>
01728    <i>Example 8 - script Deconvolution.c :</i>
01729    <p>
01730    response function (usually peak) should be shifted left to the first
01731    non-zero channel (bin) (see Figure 9)
01732    <p>
01733    <img width=600 height=340 src="gif/TSpectrum_Deconvolution1.jpg">
01734    <p>
01735    Figure 9 Response spectrum.
01736    <p>
01737    <img width=946 height=407 src="gif/TSpectrum_Deconvolution2.jpg">
01738    <p>
01739    Figure 10 Principle how the response matrix is composed inside of the
01740    Deconvolution function.
01741    <img width=601 height=407 src="gif/TSpectrum_Deconvolution3.jpg">
01742    <p>
01743    Figure 11 Example of Gold deconvolution. The original source spectrum is
01744    drawn with black color, the spectrum after the deconvolution (10000
01745    iterations) with red color.
01746    <p>
01747    Script:
01748    <p>
01749    <pre>
01750    // Example to illustrate deconvolution function (class TSpectrum).
01751    // To execute this example, do
01752    // root > .x Deconvolution.C
01753 
01754    #include <TSpectrum>
01755 
01756    void Deconvolution() {
01757       Int_t i;
01758       Double_t nbins = 256;
01759       Double_t xmin  = 0;
01760       Double_t xmax  = (Double_t)nbins;
01761       Float_t * source = new float[nbins];
01762       Float_t * response = new float[nbins];
01763       TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
01764       TH1F *d = new TH1F("d","",nbins,xmin,xmax);
01765       TFile *f = new TFile("spectra\\TSpectrum.root");
01766       h=(TH1F*) f->Get("decon1;1");
01767       TFile *fr = new TFile("spectra\\TSpectrum.root");
01768       d=(TH1F*) fr->Get("decon_response;1");
01769       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
01770       for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
01771       TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
01772       if (!Decon1) Decon1 = new TCanvas("Decon1","Decon1",10,10,1000,700);
01773       h->Draw("L");
01774       TSpectrum *s = new TSpectrum();
01775       s->Deconvolution(source,response,256,1000,1,1);
01776       for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
01777       d->SetLineColor(kRed);
01778       d->Draw("SAME L");
01779    }
01780    </pre>
01781    <p>
01782    <b>Examples of Gold deconvolution method:</b>
01783    <p>
01784    First let us study the influence of the number of iterations on the
01785    deconvolved spectrum (Figure 12).
01786    <p>
01787    <img width=602 height=409 src="gif/TSpectrum_Deconvolution_wide1.jpg">
01788    <p>
01789    Figure 12 Study of Gold deconvolution algorithm.The original source spectrum
01790    is drawn with black color, spectrum after 100 iterations with red color,
01791    spectrum after 1000 iterations with blue color, spectrum after 10000
01792    iterations with green color and spectrum after 100000 iterations with
01793    magenta color.
01794    <p>
01795    For relatively narrow peaks in the above given example the Gold
01796    deconvolution method is able to decompose overlapping peaks practically to
01797    delta - functions. In the next example we have chosen a synthetic data
01798    (spectrum, 256 channels) consisting of 5 very closely positioned, relatively
01799    wide peaks (sigma =5), with added noise (Figure 13). Thin lines represent
01800    pure Gaussians (see Table 1); thick line is a resulting spectrum with
01801    additive noise (10% of the amplitude of small peaks).
01802    <p>
01803    <img width=600 height=367 src="gif/TSpectrum_Deconvolution_wide2.jpg">
01804    <p>
01805    Figure 13 Testing example of synthetic spectrum composed of 5 Gaussians with
01806    added noise.
01807    <p>
01808    <table border=solid><tr>
01809    <td> Peak # </td><td> Position </td><td> Height </td><td> Area   </td>
01810    </tr><tr>
01811    <td> 1      </td><td> 50       </td><td> 500    </td><td> 10159  </td>
01812    </tr><tr>
01813    <td> 2      </td><td> 70       </td><td> 3000   </td><td> 60957  </td>
01814    </tr><tr>
01815    <td> 3      </td><td> 80       </td><td> 1000   </td><td> 20319  </td>
01816    </tr><tr>
01817    <td> 4      </td><td> 100      </td><td> 5000   </td><td> 101596 </td>
01818    </tr><tr>
01819    <td> 5      </td><td> 110      </td><td> 500    </td><td> 10159  </td>
01820    </tr></table>
01821    <p>
01822    Table 1 Positions, heights and areas of peaks in the spectrum shown in
01823    Figure 13.
01824    <p>
01825    In ideal case, we should obtain the result given in Figure 14. The areas of
01826    the Gaussian components of the spectrum are concentrated completely to
01827    delta-functions. When solving the overdetermined system of linear equations
01828    with data from Figure 13 in the sense of minimum least squares criterion
01829    without any regularization we obtain the result with large oscillations
01830    (Figure 15). From mathematical point of view, it is the optimal solution in
01831    the unconstrained space of independent variables. From physical point of
01832    view we are interested only in a meaningful solution. Therefore, we have to
01833    employ regularization techniques (e.g. Gold deconvolution) and/or to
01834    confine the space of allowed solutions to subspace of positive solutions.
01835    <p>
01836    <img width=589 height=189 src="gif/TSpectrum_Deconvolution_wide3.jpg">
01837    <p>
01838    Figure 14 The same spectrum like in Figure 13, outlined bars show the
01839    contents of present components (peaks).
01840    <img width=585 height=183 src="gif/TSpectrum_Deconvolution_wide4.jpg">
01841    <p>
01842    Figure 15 Least squares solution of the system of linear equations without
01843    regularization.
01844    <p>
01845    <i>Example 9 - script Deconvolution_wide.c</i>
01846    <p>
01847    When we employ Gold deconvolution algorithm we obtain the result given in
01848    Fig. 16. One can observe that the resulting spectrum is smooth. On the
01849    other hand the method is not able to decompose completely the peaks in the
01850    spectrum.
01851    <p>
01852    <img width=601 height=407 src="gif/TSpectrum_Deconvolution_wide5.jpg">
01853    Figure 16 Example of Gold deconvolution for closely positioned wide peaks.
01854    The original source spectrum is drawn with black color, the spectrum after
01855    the deconvolution (10000 iterations) with red color.
01856    <p>
01857    Script:
01858    <p>
01859    <pre>
01860    // Example to illustrate deconvolution function (class TSpectrum).
01861    // To execute this example, do
01862    // root > .x Deconvolution_wide.C
01863 
01864    #include <TSpectrum>
01865 
01866    void Deconvolution_wide() {
01867       Int_t i;
01868       Double_t nbins = 256;
01869       Double_t xmin  = 0;
01870       Double_t xmax  = (Double_t)nbins;
01871       Float_t * source = new float[nbins];
01872       Float_t * response = new float[nbins];
01873       TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
01874       TH1F *d = new TH1F("d","",nbins,xmin,xmax);
01875       TFile *f = new TFile("spectra\\TSpectrum.root");
01876       h=(TH1F*) f->Get("decon3;1");
01877       TFile *fr = new TFile("spectra\\TSpectrum.root");
01878       d=(TH1F*) fr->Get("decon_response_wide;1");
01879       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
01880       for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
01881       TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
01882       if (!Decon1) Decon1 = new TCanvas("Decon1",
01883       "Deconvolution of closely positioned overlapping peaks using Gold deconvolution method",10,10,1000,700);
01884       h->SetMaximum(30000);
01885       h->Draw("L");
01886       TSpectrum *s = new TSpectrum();
01887       s->Deconvolution(source,response,256,10000,1,1);
01888       for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
01889       d->SetLineColor(kRed);
01890       d->Draw("SAME L");
01891    }
01892    </pre>
01893    <p>
01894    <i>Example 10 - script Deconvolution_wide_boost.c :</i>
01895    <p>
01896    Further let us employ boosting operation into deconvolution (Fig. 17).
01897    <p>
01898    <img width=601 height=407 src="gif/TSpectrum_Deconvolution_wide6.jpg">
01899    <p>
01900    Figure 17 The original source spectrum is drawn with black color, the
01901    spectrum after the deconvolution with red color. Number of iterations = 200,
01902    number of repetitions = 50 and boosting coefficient = 1.2.
01903    <p>
01904    <table border=solid><tr>
01905    <td> Peak # </td> <td> Original/Estimated (max) position </td> <td> Original/Estimated area </td>
01906    </tr> <tr>
01907    <td> 1 </td> <td> 50/49 </td> <td> 10159/10419 </td>
01908    </tr> <tr>
01909    <td> 2 </td> <td> 70/70 </td> <td> 60957/58933 </td>
01910    </tr> <tr>
01911    <td> 3 </td> <td> 80/79 </td> <td> 20319/19935 </td>
01912    </tr> <tr>
01913    <td> 4 </td> <td> 100/100 </td> <td> 101596/105413 </td>
01914    </tr> <tr>
01915    <td> 5 </td> <td> 110/117 </td> <td> 10159/6676 </td>
01916    </tr> </table>
01917    <p>
01918    Table 2 Results of the estimation of peaks in spectrum shown in Figure 17.
01919    <p>
01920    One can observe that peaks are decomposed practically to delta functions.
01921    Number of peaks is correct, positions of big peaks as well as their areas
01922    are relatively well estimated. However there is a considerable error in
01923    the estimation of the position of small right hand peak.
01924    <p>
01925    Script:
01926    <p>
01927    <pre>
01928    // Example to illustrate deconvolution function (class TSpectrum).
01929    // To execute this example, do
01930    // root > .x Deconvolution_wide_boost.C
01931 
01932    #include <TSpectrum>
01933 
01934    void Deconvolution_wide_boost() {
01935       Int_t i;
01936       Double_t nbins = 256;
01937       Double_t xmin  = 0;
01938       Double_t xmax  = (Double_t)nbins;
01939       Float_t * source = new float[nbins];
01940       Float_t * response = new float[nbins];
01941       TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
01942       TH1F *d = new TH1F("d","",nbins,xmin,xmax);
01943       TFile *f = new TFile("spectra\\TSpectrum.root");
01944       h=(TH1F*) f->Get("decon3;1");
01945       TFile *fr = new TFile("spectra\\TSpectrum.root");
01946       d=(TH1F*) fr->Get("decon_response_wide;1");
01947       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
01948       for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
01949       TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
01950       if (!Decon1) Decon1 = new TCanvas("Decon1",
01951       "Deconvolution of closely positioned overlapping peaks using boosted Gold deconvolution method",10,10,1000,700);
01952       h->SetMaximum(110000);
01953       h->Draw("L");
01954       TSpectrum *s = new TSpectrum();
01955       s->Deconvolution(source,response,256,200,50,1.2);
01956       for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
01957       d->SetLineColor(kRed);
01958       d->Draw("SAME L");
01959    }
01960    </pre>
01961    End_Html */
01962 
01963    if (ssize <= 0)
01964       return "Wrong Parameters";
01965 
01966    if (numberRepetitions <= 0)
01967       return "Wrong Parameters";
01968 
01969        //   working_space-pointer to the working vector
01970        //   (its size must be 4*ssize of source spectrum)
01971    double *working_space = new double[4 * ssize];
01972    int i, j, k, lindex, posit, lh_gold, l, repet;
01973    double lda, ldb, ldc, area, maximum;
01974    area = 0;
01975    lh_gold = -1;
01976    posit = 0;
01977    maximum = 0;
01978 
01979 //read response vector
01980    for (i = 0; i < ssize; i++) {
01981       lda = response[i];
01982       if (lda != 0)
01983          lh_gold = i + 1;
01984       working_space[i] = lda;
01985       area += lda;
01986       if (lda > maximum) {
01987          maximum = lda;
01988          posit = i;
01989       }
01990    }
01991    if (lh_gold == -1) {
01992       delete [] working_space;
01993       return "ZERO RESPONSE VECTOR";
01994    }
01995 
01996 //read source vector
01997    for (i = 0; i < ssize; i++)
01998       working_space[2 * ssize + i] = source[i];
01999 
02000 // create matrix at*a and vector at*y
02001    for (i = 0; i < ssize; i++){
02002       lda = 0;
02003       for (j = 0; j < ssize; j++){
02004          ldb = working_space[j];
02005          k = i + j;
02006          if (k < ssize){
02007             ldc = working_space[k];
02008             lda = lda + ldb * ldc;
02009          }
02010       }
02011       working_space[ssize + i] = lda;
02012       lda = 0;
02013       for (k = 0; k < ssize; k++){
02014          l = k - i;
02015          if (l >= 0){
02016             ldb = working_space[l];
02017             ldc = working_space[2 * ssize + k];
02018             lda = lda + ldb * ldc;
02019          }
02020       }
02021       working_space[3 * ssize + i]=lda;
02022    }
02023 
02024 // move vector at*y
02025    for (i = 0; i < ssize; i++){
02026       working_space[2 * ssize + i] = working_space[3 * ssize + i];
02027    }
02028 
02029 //initialization of resulting vector
02030    for (i = 0; i < ssize; i++)
02031       working_space[i] = 1;
02032 
02033        //**START OF ITERATIONS**
02034    for (repet = 0; repet < numberRepetitions; repet++) {
02035       if (repet != 0) {
02036          for (i = 0; i < ssize; i++)
02037             working_space[i] = TMath::Power(working_space[i], boost);
02038       }
02039       for (lindex = 0; lindex < numberIterations; lindex++) {
02040          for (i = 0; i < ssize; i++) {
02041             if (working_space[2 * ssize + i] > 0.000001
02042                  && working_space[i] > 0.000001) {
02043                lda = 0;
02044                for (j = 0; j < lh_gold; j++) {
02045                   ldb = working_space[j + ssize];
02046                   if (j != 0){
02047                      k = i + j;
02048                      ldc = 0;
02049                      if (k < ssize)
02050                         ldc = working_space[k];
02051                      k = i - j;
02052                      if (k >= 0)
02053                         ldc += working_space[k];
02054                   }
02055 
02056                   else
02057                      ldc = working_space[i];
02058                   lda = lda + ldb * ldc;
02059                }
02060                ldb = working_space[2 * ssize + i];
02061                if (lda != 0)
02062                   lda = ldb / lda;
02063 
02064                else
02065                   lda = 0;
02066                ldb = working_space[i];
02067                lda = lda * ldb;
02068                working_space[3 * ssize + i] = lda;
02069             }
02070          }
02071          for (i = 0; i < ssize; i++)
02072             working_space[i] = working_space[3 * ssize + i];
02073       }
02074    }
02075 
02076 //shift resulting spectrum
02077    for (i = 0; i < ssize; i++) {
02078       lda = working_space[i];
02079       j = i + posit;
02080       j = j % ssize;
02081       working_space[ssize + j] = lda;
02082    }
02083 
02084 //write back resulting spectrum
02085    for (i = 0; i < ssize; i++)
02086       source[i] = area * working_space[ssize + i];
02087    delete[]working_space;
02088    return 0;
02089 }
02090 
02091 
02092 //______________________________________________________________________________
02093 const char *TSpectrum::DeconvolutionRL(float *source, const float *response,
02094                                       int ssize, int numberIterations,
02095                                       int numberRepetitions, double boost )
02096 {
02097    /* Begin_Html
02098    <b>One-dimensional deconvolution function.</b>
02099    <p>
02100    This function calculates deconvolution from source spectrum according to
02101    response spectrum using Richardson-Lucy deconvolution algorithm. The result
02102    is placed in the vector pointed by source pointer. On successful completion
02103    it returns 0. On error it returns pointer to the string describing error.
02104    If desired after every numberIterations one can apply boosting operation
02105    (exponential function with exponent given by boost coefficient) and repeat
02106    it numberRepetitions times (see Gold deconvolution).
02107    <p>
02108    Function parameters:
02109    <ul>
02110    <li> source:  pointer to the vector of source spectrum
02111    <li> response:     pointer to the vector of response spectrum
02112    <li> ssize:    length of source and response spectra
02113    numberIterations, for details we refer to the reference given above
02114    numberRepetitions, for repeated boosted deconvolution
02115    boost, boosting coefficient
02116    </ul>
02117    <p>
02118    <b>Richardson-Lucy deconvolution algorithm:</b>
02119    <p>
02120    For discrete systems it has the form:
02121    <p>
02122    <img width=438 height=98 src="gif/TSpectrum_DeconvolutionRL1.gif">
02123    <p>
02124    <img width=124 height=39 src="gif/TSpectrum_DeconvolutionRL2.gif">
02125    <p>
02126    for positive input data and response matrix this iterative method forces
02127    the deconvoluted spectra to be non-negative. The Richardson-Lucy
02128    iteration converges to the maximum likelihood solution for Poisson statistics
02129    in the data.
02130    <p>
02131    <b>References:</b>
02132    <ol>
02133    <li> Abreu M.C. et al., A four-dimensional deconvolution method to correct NA38
02134    experimental data, NIM A 405 (1998) 139.
02135    <li> Lucy L.B., A.J. 79 (1974) 745.
02136    <li> Richardson W.H., J. Opt. Soc. Am. 62 (1972) 55.
02137    </ol>
02138    <p>
02139    <b>Examples of Richardson-Lucy deconvolution method:</b>
02140    <p>
02141    <i>Example 11 - script DeconvolutionRL_wide.c :</i>
02142    <p>
02143    When we employ Richardson-Lucy deconvolution algorithm to our data from
02144    Fig. 13 we obtain the result given in Fig. 18. One can observe improvements
02145    as compared to the result achieved by Gold deconvolution. Neverthless it is
02146    unable to decompose the multiplet.
02147    <p>
02148    <img width=601 height=407 src="gif/TSpectrum_DeconvolutionRL_wide1.jpg">
02149    Figure 18 Example of Richardson-Lucy deconvolution for closely positioned
02150    wide peaks. The original source spectrum is drawn with black color, the
02151    spectrum after the deconvolution (10000 iterations) with red color.
02152    <p>
02153    Script:
02154    <p>
02155    <pre>
02156    // Example to illustrate deconvolution function (class TSpectrum).
02157    // To execute this example, do
02158    // root > .x DeconvolutionRL_wide.C
02159 
02160    #include <TSpectrum>
02161 
02162    void DeconvolutionRL_wide() {
02163       Int_t i;
02164       Double_t nbins = 256;
02165       Double_t xmin  = 0;
02166       Double_t xmax  = (Double_t)nbins;
02167       Float_t * source = new float[nbins];
02168       Float_t * response = new float[nbins];
02169       TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
02170       TH1F *d = new TH1F("d","",nbins,xmin,xmax);
02171       TFile *f = new TFile("spectra\\TSpectrum.root");
02172       h=(TH1F*) f->Get("decon3;1");
02173       TFile *fr = new TFile("spectra\\TSpectrum.root");
02174       d=(TH1F*) fr->Get("decon_response_wide;1");
02175       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02176       for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
02177       TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
02178       if (!Decon1) Decon1 = new TCanvas("Decon1",
02179       "Deconvolution of closely positioned overlapping peaks using Richardson-Lucy deconvolution method",
02180       10,10,1000,700);
02181       h->SetMaximum(30000);
02182       h->Draw("L");
02183       TSpectrum *s = new TSpectrum();
02184       s->DeconvolutionRL(source,response,256,10000,1,1);
02185       for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
02186       d->SetLineColor(kRed);
02187       d->Draw("SAME L");
02188    }
02189    </pre>
02190    <p>
02191    <i>Example 12 - script DeconvolutionRL_wide_boost.c :</i>
02192    <p>
02193    Further let us employ boosting operation into deconvolution (Fig. 19).
02194    <img width=601 height=407 src="gif/TSpectrum_DeconvolutionRL_wide2.jpg">
02195    <p>
02196    Figure 19 The original source spectrum is drawn with black color, the
02197    spectrum after the deconvolution with red color. Number of iterations = 200,
02198    number of repetitions = 50 and boosting coefficient = 1.2.
02199    <p>
02200    <table border=solid>
02201    <tr><td> Peak # </td><td> Original/Estimated (max) position </td><td> Original/Estimated area </td></tr>
02202    <tr><td> 1 </td><td> 50/51 </td><td> 10159/11426 </td></tr>
02203    <tr><td> 2 </td><td> 70/71 </td><td> 60957/65003 </td></tr>
02204    <tr><td> 3 </td><td> 80/81 </td><td> 20319/12813 </td></tr>
02205    <tr><td> 4 </td><td> 100/100 </td><td> 101596/101851 </td></tr>
02206    <tr><td> 5 </td><td> 110/111 </td><td> 10159/8920 </td></tr>
02207    </table>
02208    <p>
02209    Table 3 Results of the estimation of peaks in spectrum shown in Figure 19.
02210    <p>
02211    One can observe improvements in the estimation of peak positions as compared
02212    to the results achieved by Gold deconvolution.
02213    <p>
02214    Script:
02215    <pre>
02216    // Example to illustrate deconvolution function (class TSpectrum).
02217    // To execute this example, do
02218    // root > .x DeconvolutionRL_wide_boost.C
02219 
02220    #include <TSpectrum>
02221 
02222    void DeconvolutionRL_wide_boost() {
02223       Int_t i;
02224       Double_t nbins = 256;
02225       Double_t xmin  = 0;
02226       Double_t xmax  = (Double_t)nbins;
02227       Float_t * source = new float[nbins];
02228       Float_t * response = new float[nbins];
02229       TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
02230       TH1F *d = new TH1F("d","",nbins,xmin,xmax);
02231       TFile *f = new TFile("spectra\\TSpectrum.root");
02232       h=(TH1F*) f->Get("decon3;1");
02233       TFile *fr = new TFile("spectra\\TSpectrum.root");
02234       d=(TH1F*) fr->Get("decon_response_wide;1");
02235       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02236       for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
02237       TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
02238       if (!Decon1) Decon1 = new TCanvas("Decon1",
02239       "Deconvolution of closely positioned overlapping peaks using boosted Richardson-Lucy deconvolution method",
02240       10,10,1000,700);
02241       h->SetMaximum(110000);
02242       h->Draw("L");
02243       TSpectrum *s = new TSpectrum();
02244       s->DeconvolutionRL(source,response,256,200,50,1.2);
02245       for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
02246       d->SetLineColor(kRed);
02247       d->Draw("SAME L");
02248    }
02249    </pre>
02250    End_Html */
02251 
02252    if (ssize <= 0)
02253       return "Wrong Parameters";
02254 
02255    if (numberRepetitions <= 0)
02256       return "Wrong Parameters";
02257 
02258        //   working_space-pointer to the working vector
02259        //   (its size must be 4*ssize of source spectrum)
02260    double *working_space = new double[4 * ssize];
02261    int i, j, k, lindex, posit, lh_gold, repet, kmin, kmax;
02262    double lda, ldb, ldc, maximum;
02263    lh_gold = -1;
02264    posit = 0;
02265    maximum = 0;
02266 
02267 //read response vector
02268    for (i = 0; i < ssize; i++) {
02269       lda = response[i];
02270       if (lda != 0)
02271          lh_gold = i + 1;
02272       working_space[ssize + i] = lda;
02273       if (lda > maximum) {
02274          maximum = lda;
02275          posit = i;
02276       }
02277    }
02278    if (lh_gold == -1) {
02279       delete [] working_space;
02280       return "ZERO RESPONSE VECTOR";
02281    }
02282 
02283 //read source vector
02284    for (i = 0; i < ssize; i++)
02285       working_space[2 * ssize + i] = source[i];
02286 
02287 //initialization of resulting vector
02288    for (i = 0; i < ssize; i++){
02289       if (i <= ssize - lh_gold)
02290          working_space[i] = 1;
02291 
02292       else
02293          working_space[i] = 0;
02294 
02295    }
02296        //**START OF ITERATIONS**
02297    for (repet = 0; repet < numberRepetitions; repet++) {
02298       if (repet != 0) {
02299          for (i = 0; i < ssize; i++)
02300             working_space[i] = TMath::Power(working_space[i], boost);
02301       }
02302       for (lindex = 0; lindex < numberIterations; lindex++) {
02303          for (i = 0; i <= ssize - lh_gold; i++){
02304             lda = 0;
02305             if (working_space[i] > 0){//x[i]
02306                for (j = i; j < i + lh_gold; j++){
02307                   ldb = working_space[2 * ssize + j];//y[j]
02308                   if (j < ssize){
02309                      if (ldb > 0){//y[j]
02310                         kmax = j;
02311                         if (kmax > lh_gold - 1)
02312                            kmax = lh_gold - 1;
02313                         kmin = j + lh_gold - ssize;
02314                         if (kmin < 0)
02315                            kmin = 0;
02316                         ldc = 0;
02317                         for (k = kmax; k >= kmin; k--){
02318                            ldc += working_space[ssize + k] * working_space[j - k];//h[k]*x[j-k]
02319                         }
02320                         if (ldc > 0)
02321                            ldb = ldb / ldc;
02322 
02323                         else
02324                            ldb = 0;
02325                      }
02326                      ldb = ldb * working_space[ssize + j - i];//y[j]*h[j-i]/suma(h[j][k]x[k])
02327                   }
02328                   lda += ldb;
02329                }
02330                lda = lda * working_space[i];
02331             }
02332             working_space[3 * ssize + i] = lda;
02333          }
02334          for (i = 0; i < ssize; i++)
02335             working_space[i] = working_space[3 * ssize + i];
02336       }
02337    }
02338 
02339 //shift resulting spectrum
02340    for (i = 0; i < ssize; i++) {
02341       lda = working_space[i];
02342       j = i + posit;
02343       j = j % ssize;
02344       working_space[ssize + j] = lda;
02345    }
02346 
02347 //write back resulting spectrum
02348    for (i = 0; i < ssize; i++)
02349       source[i] = working_space[ssize + i];
02350    delete[]working_space;
02351    return 0;
02352 }
02353 
02354 
02355 //______________________________________________________________________________
02356 const char *TSpectrum::Unfolding(float *source,
02357                                  const float **respMatrix,
02358                                  int ssizex, int ssizey,
02359                                  int numberIterations,
02360                                  int numberRepetitions, double boost)
02361 {
02362    /* Begin_Html
02363    <b>One-dimensional unfolding function</b>
02364    <p>
02365    This function unfolds source spectrum according to response matrix columns.
02366    The result is placed in the vector pointed by source pointer.
02367    The coefficients of the resulting vector represent contents of the columns
02368    (weights) in the input vector. On successful completion it returns 0. On
02369    error it returns pointer to the string describing error. If desired after
02370    every numberIterations one can apply boosting operation (exponential
02371    function with exponent given by boost coefficient) and repeat it
02372    numberRepetitions times. For details we refer to [1].
02373    <p>
02374    Function parameters:
02375    <ul>
02376    <li> source: pointer to the vector of source spectrum
02377    <li> respMatrix: pointer to the matrix of response spectra
02378    <li> ssizex: length of source spectrum and # of columns of the response
02379         matrix. ssizex must be >= ssizey.
02380    <li> ssizey: length of destination spectrum and # of rows of the response
02381         matrix.
02382    <li> numberIterations: number of iterations
02383    <li> numberRepetitions: number of repetitions for boosted deconvolution.
02384         It must be greater or equal to one.
02385    <li> boost: boosting coefficient, applies only if numberRepetitions is
02386         greater than one.
02387    </ul>
02388    <p>
02389    <b>Unfolding:</b>
02390    <p>
02391    The goal is the decomposition of spectrum to a given set of component
02392    spectra.
02393    <p>
02394    The mathematical formulation of the discrete linear system is:
02395    <p>
02396    <img width=588 height=89 src="gif/TSpectrum_Unfolding1.gif">
02397    <p>
02398    <img width=597 height=228 src="gif/TSpectrum_Unfolding2.gif">
02399    <p>
02400    <b>References:</b>
02401    <ol>
02402    <li> Jandel M., Morhá&#269; M., Kliman J., Krupa L., Matoušek
02403    V., Hamilton J. H., Ramaya A. V.:
02404    Decomposition of continuum gamma-ray spectra using synthetized response matrix.
02405    NIM A 516 (2004), 172-183.
02406    </ol>
02407    <p>
02408    <b>Example of unfolding:</b>
02409    <p>
02410    <i>Example 13 - script Unfolding.c:</i>
02411    <p>
02412    <img width=442 height=648 src="gif/TSpectrum_Unfolding3.gif">
02413    <p>
02414    Fig. 20 Response matrix composed of neutron spectra of pure
02415    chemical elements.
02416    <img width=604 height=372 src="gif/TSpectrum_Unfolding2.jpg">
02417    <p>
02418    Fig. 21 Source neutron spectrum to be decomposed
02419    <P>
02420    <img width=600 height=360 src="gif/TSpectrum_Unfolding3.jpg">
02421    <p>
02422    Fig. 22 Spectrum after decomposition, contains 10 coefficients, which
02423    correspond to contents of chemical components (dominant 8-th and 10-th
02424    components, i.e. O, Si)
02425    <p>
02426    Script:
02427    <pre>
02428    // Example to illustrate unfolding function (class TSpectrum).
02429    // To execute this example, do
02430    // root > .x Unfolding.C
02431 
02432    #include <TSpectrum>
02433 
02434    void Unfolding() {
02435       Int_t i, j;
02436       Int_t nbinsx = 2048;
02437       Int_t nbinsy = 10;
02438       Double_t xmin  = 0;
02439       Double_t xmax  = (Double_t)nbinsx;
02440       Double_t ymin  = 0;
02441       Double_t ymax  = (Double_t)nbinsy;
02442       Float_t * source = new float[nbinsx];
02443       Float_t ** response = new float *[nbinsy];
02444       for (i=0;i<nbinsy;i++) response[i]=new float[nbinsx];
02445       TH1F *h = new TH1F("h","",nbinsx,xmin,xmax);
02446       TH1F *d = new TH1F("d","Decomposition - unfolding",nbinsx,xmin,xmax);
02447       TH2F *decon_unf_resp = new TH2F("decon_unf_resp","Root File",nbinsy,ymin,ymax,nbinsx,xmin,xmax);
02448       TFile *f = new TFile("spectra\\TSpectrum.root");
02449       h=(TH1F*) f->Get("decon_unf_in;1");
02450       TFile *fr = new TFile("spectra\\TSpectrum.root");
02451       decon_unf_resp = (TH2F*) fr->Get("decon_unf_resp;1");
02452       for (i = 0; i < nbinsx; i++) source[i] = h->GetBinContent(i + 1);
02453       for (i = 0; i < nbinsy; i++){
02454          for (j = 0; j< nbinsx; j++){
02455             response[i][j] = decon_unf_resp->GetBinContent(i + 1, j + 1);
02456          }
02457       }
02458       TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
02459       if (!Decon1) Decon1 = new TCanvas("Decon1","Decon1",10,10,1000,700);
02460       h->Draw("L");
02461       TSpectrum *s = new TSpectrum();
02462       s->Unfolding(source,response,nbinsx,nbinsy,1000,1,1);
02463       for (i = 0; i < nbinsy; i++) d->SetBinContent(i + 1,source[i]);
02464       d->SetLineColor(kRed);
02465       d->SetAxisRange(0,nbinsy);
02466       d->Draw("");
02467    }
02468    </pre>
02469    End_Html */
02470 
02471    int i, j, k, lindex, lhx = 0, repet;
02472    double lda, ldb, ldc, area;
02473    if (ssizex <= 0 || ssizey <= 0)
02474       return "Wrong Parameters";
02475    if (ssizex < ssizey)
02476       return "Sizex must be greater than sizey)";
02477    if (numberIterations <= 0)
02478       return "Number of iterations must be positive";
02479    double *working_space =
02480        new double[ssizex * ssizey + 2 * ssizey * ssizey + 4 * ssizex];
02481 
02482 /*read response matrix*/
02483    for (j = 0; j < ssizey && lhx != -1; j++) {
02484       area = 0;
02485       lhx = -1;
02486       for (i = 0; i < ssizex; i++) {
02487          lda = respMatrix[j][i];
02488          if (lda != 0) {
02489             lhx = i + 1;
02490          }
02491          working_space[j * ssizex + i] = lda;
02492          area = area + lda;
02493       }
02494       if (lhx != -1) {
02495          for (i = 0; i < ssizex; i++)
02496             working_space[j * ssizex + i] /= area;
02497       }
02498    }
02499    if (lhx == -1) {
02500       delete [] working_space;
02501       return ("ZERO COLUMN IN RESPONSE MATRIX");
02502    }
02503 
02504 /*read source vector*/
02505    for (i = 0; i < ssizex; i++)
02506       working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
02507           source[i];
02508 
02509 /*create matrix at*a + at*y */
02510    for (i = 0; i < ssizey; i++) {
02511       for (j = 0; j < ssizey; j++) {
02512          lda = 0;
02513          for (k = 0; k < ssizex; k++) {
02514             ldb = working_space[ssizex * i + k];
02515             ldc = working_space[ssizex * j + k];
02516             lda = lda + ldb * ldc;
02517          }
02518          working_space[ssizex * ssizey + ssizey * i + j] = lda;
02519       }
02520       lda = 0;
02521       for (k = 0; k < ssizex; k++) {
02522          ldb = working_space[ssizex * i + k];
02523          ldc =
02524              working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
02525                            k];
02526          lda = lda + ldb * ldc;
02527       }
02528       working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
02529           lda;
02530    }
02531 
02532 /*move vector at*y*/
02533    for (i = 0; i < ssizey; i++)
02534       working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
02535           working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
02536 
02537 /*create matrix at*a*at*a + vector at*a*at*y */
02538    for (i = 0; i < ssizey; i++) {
02539       for (j = 0; j < ssizey; j++) {
02540          lda = 0;
02541          for (k = 0; k < ssizey; k++) {
02542             ldb = working_space[ssizex * ssizey + ssizey * i + k];
02543             ldc = working_space[ssizex * ssizey + ssizey * j + k];
02544             lda = lda + ldb * ldc;
02545          }
02546          working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j] =
02547              lda;
02548       }
02549       lda = 0;
02550       for (k = 0; k < ssizey; k++) {
02551          ldb = working_space[ssizex * ssizey + ssizey * i + k];
02552          ldc =
02553              working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
02554                            k];
02555          lda = lda + ldb * ldc;
02556       }
02557       working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
02558           lda;
02559    }
02560 
02561 /*move at*a*at*y*/
02562    for (i = 0; i < ssizey; i++)
02563       working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
02564           working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
02565 
02566 /*initialization in resulting vector */
02567    for (i = 0; i < ssizey; i++)
02568       working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = 1;
02569 
02570         /***START OF ITERATIONS***/
02571    for (repet = 0; repet < numberRepetitions; repet++) {
02572       if (repet != 0) {
02573          for (i = 0; i < ssizey; i++)
02574             working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = TMath::Power(working_space[ssizex * ssizey + 2 * ssizey * ssizey + i], boost);
02575       }
02576       for (lindex = 0; lindex < numberIterations; lindex++) {
02577          for (i = 0; i < ssizey; i++) {
02578             lda = 0;
02579             for (j = 0; j < ssizey; j++) {
02580                ldb =
02581                    working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j];
02582                ldc = working_space[ssizex * ssizey + 2 * ssizey * ssizey + j];
02583                lda = lda + ldb * ldc;
02584             }
02585             ldb =
02586                 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i];
02587             if (lda != 0) {
02588                lda = ldb / lda;
02589             }
02590 
02591             else
02592                lda = 0;
02593             ldb = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
02594             lda = lda * ldb;
02595             working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] = lda;
02596          }
02597          for (i = 0; i < ssizey; i++)
02598             working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] =
02599                 working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
02600       }
02601    }
02602 
02603 /*write back resulting spectrum*/
02604    for (i = 0; i < ssizex; i++) {
02605       if (i < ssizey)
02606          source[i] = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
02607 
02608       else
02609          source[i] = 0;
02610    }
02611    delete[]working_space;
02612    return 0;
02613 }
02614 
02615 
02616 //______________________________________________________________________________
02617 Int_t TSpectrum::SearchHighRes(float *source,float *destVector, int ssize,
02618                                      float sigma, double threshold,
02619                                      bool backgroundRemove,int deconIterations,
02620                                      bool markov, int averWindow)
02621 {
02622    /* Begin_Html
02623    <b>One-dimensional high-resolution peak search function</b>
02624    <p>
02625    This function searches for peaks in source spectrum. It is based on
02626    deconvolution method. First the background is removed (if desired), then
02627    Markov smoothed spectrum is calculated (if desired), then the response
02628    function is generated according to given sigma and deconvolution is
02629    carried out. The order of peaks is arranged according to their heights in
02630    the spectrum after background elimination. The highest peak is the first in
02631    the list. On success it returns number of found peaks.
02632    <p>
02633    <b>Function parameters:</b>
02634    <ul>
02635    <li> source: pointer to the vector of source spectrum.
02636    <li> destVector: pointer to the vector of resulting deconvolved spectrum.
02637    <li> ssize: length of source spectrum.
02638    <li> sigma: sigma of searched peaks, for details we refer to manual.
02639    <li> threshold: threshold value in % for selected peaks, peaks with
02640         amplitude less than threshold*highest_peak/100
02641         are ignored, see manual.
02642    <li> backgroundRemove: logical variable, set if the removal of
02643         background before deconvolution is desired.
02644    <li> deconIterations-number of iterations in deconvolution operation.
02645    <li> markov: logical variable, if it is true, first the source spectrum
02646         is replaced by new spectrum calculated using Markov
02647         chains method.
02648    <li> averWindow: averanging window of searched peaks, for details
02649         we refer to manual (applies only for Markov method).
02650    </ul>
02651    <p>
02652    <b>Peaks searching:</b>
02653    <p>
02654    The goal of this function is to identify automatically the peaks in spectrum
02655    with the presence of the continuous background and statistical
02656    fluctuations - noise.
02657    <p>
02658    The common problems connected with correct peak identification are:
02659    <ul>
02660    <li> non-sensitivity to noise, i.e., only statistically
02661      relevant peaks should be identified.
02662    <li> non-sensitivity of the algorithm to continuous
02663      background.
02664    <li> ability to identify peaks close to the edges of the
02665      spectrum region. Usually peak finders fail to detect them.
02666    <li> resolution, decomposition of doublets and multiplets.
02667      The algorithm should be able to recognize close positioned peaks.
02668    <li> ability to identify peaks with different sigma.
02669    </ul>
02670    <img width=600 height=375 src="gif/TSpectrum_Searching1.jpg">
02671    <p>
02672    Fig. 27 An example of one-dimensional synthetic spectrum with found peaks
02673    denoted by markers.
02674    <p>
02675    <b>References:</b>
02676    <ol>
02677    <li> M.A. Mariscotti: A method for identification of peaks in the presence of
02678    background and its application to spectrum analysis. NIM 50 (1967),
02679    309-320.
02680    <li> M. Morhá&#269;, J. Kliman, V.  Matoušek, M. Veselský,
02681    I. Turzo.:Identification of peaks in
02682    multidimensional coincidence gamma-ray spectra. NIM, A443 (2000) 108-125.
02683    <li> Z.K. Silagadze, A new algorithm for automatic photopeak searches. NIM
02684    A 376 (1996), 451.
02685    </ol>
02686    <p>
02687    <b>Examples of peak searching method:</b>
02688    <p>
02689    The SearchHighRes function provides users with the possibility to vary the
02690    input parameters and with the access to the output deconvolved data in the
02691    destination spectrum. Based on the output data one can tune the parameters.
02692    <p>
02693    Example 15 - script SearchHR1.c:
02694    <img width=600 height=321 src="gif/TSpectrum_Searching1.jpg">
02695    <p>
02696    Fig. 28 One-dimensional spectrum with found peaks denoted by markers, 3
02697    iterations steps in the deconvolution.
02698    <p>
02699    <img width=600 height=323 src="gif/TSpectrum_Searching2.jpg">
02700    Fig. 29 One-dimensional spectrum with found peaks denoted by markers, 8
02701    iterations steps in the deconvolution.
02702    <p>
02703    Script:
02704    <pre>
02705    // Example to illustrate high resolution peak searching function (class TSpectrum).
02706    // To execute this example, do
02707    // root > .x SearchHR1.C
02708 
02709    #include <TSpectrum>
02710 
02711    void SearchHR1() {
02712       Float_t fPositionX[100];
02713       Float_t fPositionY[100];
02714       Int_t fNPeaks = 0;
02715       Int_t i,nfound,bin;
02716       Double_t nbins = 1024,a;
02717       Double_t xmin  = 0;
02718       Double_t xmax  = (Double_t)nbins;
02719       Float_t * source = new float[nbins];
02720       Float_t * dest = new float[nbins];
02721       TH1F *h = new TH1F("h","High resolution peak searching, number of iterations = 3",nbins,xmin,xmax);
02722       TH1F *d = new TH1F("d","",nbins,xmin,xmax);
02723       TFile *f = new TFile("spectra\\TSpectrum.root");
02724       h=(TH1F*) f->Get("search2;1");
02725       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02726       TCanvas *Search = gROOT->GetListOfCanvases()->FindObject("Search");
02727       if (!Search) Search = new TCanvas("Search","Search",10,10,1000,700);
02728       h->SetMaximum(4000);
02729       h->Draw("L");
02730       TSpectrum *s = new TSpectrum();
02731       nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
02732       Float_t *xpeaks = s->GetPositionX();
02733       for (i = 0; i < nfound; i++) {
02734          a=xpeaks[i];
02735          bin = 1 + Int_t(a + 0.5);
02736          fPositionX[i] = h->GetBinCenter(bin);
02737          fPositionY[i] = h->GetBinContent(bin);
02738       }
02739       TPolyMarker * pm = (TPolyMarker*)h->GetListOfFunctions()->FindObject("TPolyMarker");
02740       if (pm) {
02741          h->GetListOfFunctions()->Remove(pm);
02742          delete pm;
02743       }
02744       pm = new TPolyMarker(nfound, fPositionX, fPositionY);
02745       h->GetListOfFunctions()->Add(pm);
02746       pm->SetMarkerStyle(23);
02747       pm->SetMarkerColor(kRed);
02748       pm->SetMarkerSize(1.3);
02749       for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,dest[i]);
02750       d->SetLineColor(kRed);
02751       d->Draw("SAME");
02752       printf("Found %d candidate peaks\n",nfound);
02753       for(i=0;i<nfound;i++)
02754          printf("posx= %d, posy= %d\n",fPositionX[i], fPositionY[i]);
02755       }
02756    </pre>
02757    <p>
02758    Example 16 - script SearchHR3.c:
02759    <p>
02760    <table border=solid>
02761    <tr><td> Peak # </td><td> Position </td><td> Sigma </td></tr>
02762    <tr><td> 1      </td><td> 118      </td><td> 26    </td></tr>
02763    <tr><td> 2      </td><td> 162      </td><td> 41    </td></tr>
02764    <tr><td> 3      </td><td> 310      </td><td> 4     </td></tr>
02765    <tr><td> 4      </td><td> 330      </td><td> 8     </td></tr>
02766    <tr><td> 5      </td><td> 482      </td><td> 22    </td></tr>
02767    <tr><td> 6      </td><td> 491      </td><td> 26    </td></tr>
02768    <tr><td> 7      </td><td> 740      </td><td> 21    </td></tr>
02769    <tr><td> 8      </td><td> 852      </td><td> 15    </td></tr>
02770    <tr><td> 9      </td><td> 954      </td><td> 12    </td></tr>
02771    <tr><td> 10     </td><td> 989      </td><td> 13    </td></tr>
02772    </table>
02773    <p>
02774    Table 4 Positions and sigma of peaks in the following examples.
02775    <p>
02776    <img width=600 height=328 src="gif/TSpectrum_Searching3.jpg">
02777    <p>
02778    Fig. 30 Influence of number of iterations (3-red, 10-blue, 100- green,
02779    1000-magenta), sigma=8, smoothing width=3.
02780    <p>
02781    <img width=600 height=321 src="gif/TSpectrum_Searching4.jpg">
02782    <p>
02783    Fig. 31 Influence of sigma (3-red, 8-blue, 20- green, 43-magenta),
02784    num. iter.=10, sm. width=3.
02785    <p>
02786    <img width=600 height=323 src="gif/TSpectrum_Searching5.jpg"></p>
02787    <p>
02788    Fig. 32 Influence smoothing width (0-red, 3-blue, 7- green, 20-magenta), num.
02789    iter.=10, sigma=8.
02790    <p>
02791    Script:
02792    <pre>
02793    // Example to illustrate the influence of number of iterations in deconvolution in high resolution peak searching function (class TSpectrum).
02794    // To execute this example, do
02795    // root > .x SearchHR3.C
02796 
02797    #include <TSpectrum>
02798 
02799    void SearchHR3() {
02800       Float_t fPositionX[100];
02801       Float_t fPositionY[100];
02802       Int_t fNPeaks = 0;
02803       Int_t i,nfound,bin;
02804       Double_t nbins = 1024,a;
02805       Double_t xmin  = 0;
02806       Double_t xmax  = (Double_t)nbins;
02807       Float_t * source = new float[nbins];
02808       Float_t * dest = new float[nbins];
02809       TH1F *h = new TH1F("h","Influence of # of iterations in deconvolution in peak searching",nbins,xmin,xmax);
02810       TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
02811       TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
02812       TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
02813       TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
02814       TFile *f = new TFile("spectra\\TSpectrum.root");
02815       h=(TH1F*) f->Get("search3;1");
02816       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02817       TCanvas *Search = gROOT->GetListOfCanvases()->FindObject("Search");
02818       if (!Search) Search = new TCanvas("Search","Search",10,10,1000,700);
02819       h->SetMaximum(1300);
02820       h->Draw("L");
02821       TSpectrum *s = new TSpectrum();
02822       nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
02823       Float_t *xpeaks = s->GetPositionX();
02824       for (i = 0; i < nfound; i++) {
02825          a=xpeaks[i];
02826          bin = 1 + Int_t(a + 0.5);
02827          fPositionX[i] = h->GetBinCenter(bin);
02828          fPositionY[i] = h->GetBinContent(bin);
02829       }
02830       TPolyMarker * pm = (TPolyMarker*)h->GetListOfFunctions()->FindObject("TPolyMarker");
02831       if (pm) {
02832          h->GetListOfFunctions()->Remove(pm);
02833          delete pm;
02834       }
02835       pm = new TPolyMarker(nfound, fPositionX, fPositionY);
02836       h->GetListOfFunctions()->Add(pm);
02837       pm->SetMarkerStyle(23);
02838       pm->SetMarkerColor(kRed);
02839       pm->SetMarkerSize(1.3);
02840       for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,dest[i]);
02841       h->Draw("");
02842       d1->SetLineColor(kRed);
02843       d1->Draw("SAME");
02844       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02845       s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 10, kTRUE, 3);
02846       for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,dest[i]);
02847       d2->SetLineColor(kBlue);
02848       d2->Draw("SAME");
02849       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02850       s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 100, kTRUE, 3);
02851       for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,dest[i]);
02852       d3->SetLineColor(kGreen);
02853       d3->Draw("SAME");
02854       for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
02855       s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 1000, kTRUE, 3);
02856       for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,dest[i]);
02857       d4->SetLineColor(kMagenta);
02858       d4->Draw("SAME");
02859       printf("Found %d candidate peaks\n",nfound);
02860    }
02861    </pre>
02862    End_Html */
02863 
02864    int i, j, numberIterations = (int)(7 * sigma + 0.5);
02865    double a, b, c;
02866    int k, lindex, posit, imin, imax, jmin, jmax, lh_gold, priz;
02867    double lda, ldb, ldc, area, maximum, maximum_decon;
02868    int xmin, xmax, l, peak_index = 0, size_ext = ssize + 2 * numberIterations, shift = numberIterations, bw = 2, w;
02869    double maxch;
02870    double nom, nip, nim, sp, sm, plocha = 0;
02871    double m0low=0,m1low=0,m2low=0,l0low=0,l1low=0,detlow,av,men;
02872    if (sigma < 1) {
02873       Error("SearchHighRes", "Invalid sigma, must be greater than or equal to 1");
02874       return 0;
02875    }
02876 
02877    if(threshold<=0 || threshold>=100){
02878       Error("SearchHighRes", "Invalid threshold, must be positive and less than 100");
02879       return 0;
02880    }
02881 
02882    j = (int) (5.0 * sigma + 0.5);
02883    if (j >= PEAK_WINDOW / 2) {
02884       Error("SearchHighRes", "Too large sigma");
02885       return 0;
02886    }
02887 
02888    if (markov == true) {
02889       if (averWindow <= 0) {
02890          Error("SearchHighRes", "Averanging window must be positive");
02891          return 0;
02892       }
02893    }
02894 
02895    if(backgroundRemove == true){
02896       if(ssize < 2 * numberIterations + 1){
02897          Error("SearchHighRes", "Too large clipping window");
02898          return 0;
02899       }
02900    }
02901 
02902    k = int(2 * sigma+0.5);
02903    if(k >= 2){
02904       for(i = 0;i < k;i++){
02905          a = i,b = source[i];
02906          m0low += 1,m1low += a,m2low += a * a,l0low += b,l1low += a * b;
02907       }
02908       detlow = m0low * m2low - m1low * m1low;
02909       if(detlow != 0)
02910          l1low = (-l0low * m1low + l1low * m0low) / detlow;
02911 
02912       else
02913          l1low = 0;
02914       if(l1low > 0)
02915          l1low=0;
02916    }
02917 
02918    else{
02919       l1low = 0;
02920    }
02921 
02922    i = (int)(7 * sigma + 0.5);
02923    i = 2 * i;
02924    double *working_space = new double [7 * (ssize + i)];
02925    for (j=0;j<7 * (ssize + i);j++) working_space[j] = 0;
02926    for(i = 0; i < size_ext; i++){
02927       if(i < shift){
02928          a = i - shift;
02929          working_space[i + size_ext] = source[0] + l1low * a;
02930          if(working_space[i + size_ext] < 0)
02931             working_space[i + size_ext]=0;
02932       }
02933 
02934       else if(i >= ssize + shift){
02935          a = i - (ssize - 1 + shift);
02936          working_space[i + size_ext] = source[ssize - 1];
02937          if(working_space[i + size_ext] < 0)
02938             working_space[i + size_ext]=0;
02939       }
02940 
02941       else
02942          working_space[i + size_ext] = source[i - shift];
02943    }
02944 
02945    if(backgroundRemove == true){
02946       for(i = 1; i <= numberIterations; i++){
02947          for(j = i; j < size_ext - i; j++){
02948             if(markov == false){
02949                a = working_space[size_ext + j];
02950                b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
02951                if(b < a)
02952                   a = b;
02953 
02954                working_space[j]=a;
02955             }
02956 
02957             else{
02958                a = working_space[size_ext + j];
02959                av = 0;
02960                men = 0;
02961                for (w = j - bw; w <= j + bw; w++){
02962                   if ( w >= 0 && w < size_ext){
02963                      av += working_space[size_ext + w];
02964                      men +=1;
02965                   }
02966                }
02967                av = av / men;
02968                b = 0;
02969                men = 0;
02970                for (w = j - i - bw; w <= j - i + bw; w++){
02971                   if ( w >= 0 && w < size_ext){
02972                      b += working_space[size_ext + w];
02973                      men +=1;
02974                   }
02975                }
02976                b = b / men;
02977                c = 0;
02978                men = 0;
02979                for (w = j + i - bw; w <= j + i + bw; w++){
02980                   if ( w >= 0 && w < size_ext){
02981                      c += working_space[size_ext + w];
02982                      men +=1;
02983                   }
02984                }
02985                c = c / men;
02986                b = (b + c) / 2;
02987                if (b < a)
02988                   av = b;
02989                working_space[j]=av;
02990             }
02991          }
02992          for(j = i; j < size_ext - i; j++)
02993             working_space[size_ext + j] = working_space[j];
02994       }
02995       for(j = 0;j < size_ext; j++){
02996          if(j < shift){
02997                   a = j - shift;
02998                   b = source[0] + l1low * a;
02999                   if (b < 0) b = 0;
03000             working_space[size_ext + j] = b - working_space[size_ext + j];
03001          }
03002 
03003          else if(j >= ssize + shift){
03004                   a = j - (ssize - 1 + shift);
03005                   b = source[ssize - 1];
03006                   if (b < 0) b = 0;
03007             working_space[size_ext + j] = b - working_space[size_ext + j];
03008          }
03009 
03010          else{
03011             working_space[size_ext + j] = source[j - shift] - working_space[size_ext + j];
03012          }
03013       }
03014       for(j = 0;j < size_ext; j++){
03015          if(working_space[size_ext + j] < 0) working_space[size_ext + j] = 0;
03016       }
03017    }
03018 
03019    for(i = 0; i < size_ext; i++){
03020       working_space[i + 6*size_ext] = working_space[i + size_ext];
03021    }
03022 
03023    if(markov == true){
03024       for(j = 0; j < size_ext; j++)
03025          working_space[2 * size_ext + j] = working_space[size_ext + j];
03026       xmin = 0,xmax = size_ext - 1;
03027       for(i = 0, maxch = 0; i < size_ext; i++){
03028          working_space[i] = 0;
03029          if(maxch < working_space[2 * size_ext + i])
03030             maxch = working_space[2 * size_ext + i];
03031          plocha += working_space[2 * size_ext + i];
03032       }
03033       if(maxch == 0) {
03034          delete [] working_space;
03035          return 0;
03036       }
03037 
03038       nom = 1;
03039       working_space[xmin] = 1;
03040       for(i = xmin; i < xmax; i++){
03041          nip = working_space[2 * size_ext + i] / maxch;
03042          nim = working_space[2 * size_ext + i + 1] / maxch;
03043          sp = 0,sm = 0;
03044          for(l = 1; l <= averWindow; l++){
03045             if((i + l) > xmax)
03046                a = working_space[2 * size_ext + xmax] / maxch;
03047 
03048             else
03049                a = working_space[2 * size_ext + i + l] / maxch;
03050 
03051             b = a - nip;
03052             if(a + nip <= 0)
03053                a=1;
03054 
03055             else
03056                a = TMath::Sqrt(a + nip);
03057 
03058             b = b / a;
03059             b = TMath::Exp(b);
03060             sp = sp + b;
03061             if((i - l + 1) < xmin)
03062                a = working_space[2 * size_ext + xmin] / maxch;
03063 
03064             else
03065                a = working_space[2 * size_ext + i - l + 1] / maxch;
03066 
03067             b = a - nim;
03068             if(a + nim <= 0)
03069                a = 1;
03070 
03071             else
03072                a = TMath::Sqrt(a + nim);
03073 
03074             b = b / a;
03075             b = TMath::Exp(b);
03076             sm = sm + b;
03077          }
03078          a = sp / sm;
03079          a = working_space[i + 1] = working_space[i] * a;
03080          nom = nom + a;
03081       }
03082       for(i = xmin; i <= xmax; i++){
03083          working_space[i] = working_space[i] / nom;
03084       }
03085       for(j = 0; j < size_ext; j++)
03086          working_space[size_ext + j] = working_space[j] * plocha;
03087       for(j = 0; j < size_ext; j++){
03088          working_space[2 * size_ext + j] = working_space[size_ext + j];
03089       }
03090       if(backgroundRemove == true){
03091          for(i = 1; i <= numberIterations; i++){
03092             for(j = i; j < size_ext - i; j++){
03093                a = working_space[size_ext + j];
03094                b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
03095                if(b < a)
03096                   a = b;
03097                working_space[j] = a;
03098             }
03099             for(j = i; j < size_ext - i; j++)
03100                working_space[size_ext + j] = working_space[j];
03101          }
03102          for(j = 0; j < size_ext; j++){
03103             working_space[size_ext + j] = working_space[2 * size_ext + j] - working_space[size_ext + j];
03104          }
03105       }
03106    }
03107 //deconvolution starts
03108    area = 0;
03109    lh_gold = -1;
03110    posit = 0;
03111    maximum = 0;
03112 //generate response vector
03113    for(i = 0; i < size_ext; i++){
03114       lda = (double)i - 3 * sigma;
03115       lda = lda * lda / (2 * sigma * sigma);
03116       j = (int)(1000 * TMath::Exp(-lda));
03117       lda = j;
03118       if(lda != 0)
03119          lh_gold = i + 1;
03120 
03121       working_space[i] = lda;
03122       area = area + lda;
03123       if(lda > maximum){
03124          maximum = lda;
03125          posit = i;
03126       }
03127    }
03128 //read source vector
03129    for(i = 0; i < size_ext; i++)
03130       working_space[2 * size_ext + i] = TMath::Abs(working_space[size_ext + i]);
03131 //create matrix at*a(vector b)
03132    i = lh_gold - 1;
03133    if(i > size_ext)
03134       i = size_ext;
03135 
03136    imin = -i,imax = i;
03137    for(i = imin; i <= imax; i++){
03138       lda = 0;
03139       jmin = 0;
03140       if(i < 0)
03141          jmin = -i;
03142       jmax = lh_gold - 1 - i;
03143       if(jmax > (lh_gold - 1))
03144          jmax = lh_gold - 1;
03145 
03146       for(j = jmin;j <= jmax; j++){
03147          ldb = working_space[j];
03148          ldc = working_space[i + j];
03149          lda = lda + ldb * ldc;
03150       }
03151       working_space[size_ext + i - imin] = lda;
03152    }
03153 //create vector p
03154    i = lh_gold - 1;
03155    imin = -i,imax = size_ext + i - 1;
03156    for(i = imin; i <= imax; i++){
03157       lda = 0;
03158       for(j = 0; j <= (lh_gold - 1); j++){
03159          ldb = working_space[j];
03160          k = i + j;
03161          if(k >= 0 && k < size_ext){
03162             ldc = working_space[2 * size_ext + k];
03163             lda = lda + ldb * ldc;
03164          }
03165 
03166       }
03167       working_space[4 * size_ext + i - imin] = lda;
03168    }
03169 //move vector p
03170    for(i = imin; i <= imax; i++)
03171       working_space[2 * size_ext + i - imin] = working_space[4 * size_ext + i - imin];
03172 //initialization of resulting vector
03173    for(i = 0; i < size_ext; i++)
03174       working_space[i] = 1;
03175 //START OF ITERATIONS
03176    for(lindex = 0; lindex < deconIterations; lindex++){
03177       for(i = 0; i < size_ext; i++){
03178          if(TMath::Abs(working_space[2 * size_ext + i]) > 0.00001 && TMath::Abs(working_space[i]) > 0.00001){
03179             lda=0;
03180             jmin = lh_gold - 1;
03181             if(jmin > i)
03182                jmin = i;
03183 
03184             jmin = -jmin;
03185             jmax = lh_gold - 1;
03186             if(jmax > (size_ext - 1 - i))
03187                jmax=size_ext-1-i;
03188 
03189             for(j = jmin; j <= jmax; j++){
03190                ldb = working_space[j + lh_gold - 1 + size_ext];
03191                ldc = working_space[i + j];
03192                lda = lda + ldb * ldc;
03193             }
03194             ldb = working_space[2 * size_ext + i];
03195             if(lda != 0)
03196                lda = ldb / lda;
03197 
03198             else
03199                lda = 0;
03200 
03201             ldb = working_space[i];
03202             lda = lda * ldb;
03203             working_space[3 * size_ext + i] = lda;
03204          }
03205       }
03206       for(i = 0; i < size_ext; i++){
03207          working_space[i] = working_space[3 * size_ext + i];
03208       }
03209    }
03210 //shift resulting spectrum
03211    for(i=0;i<size_ext;i++){
03212       lda = working_space[i];
03213       j = i + posit;
03214       j = j % size_ext;
03215       working_space[size_ext + j] = lda;
03216    }
03217 //write back resulting spectrum
03218    maximum = 0, maximum_decon = 0;
03219    j = lh_gold - 1;
03220    for(i = 0; i < size_ext - j; i++){
03221       if(i >= shift && i < ssize + shift){
03222          working_space[i] = area * working_space[size_ext + i + j];
03223          if(maximum_decon < working_space[i])
03224             maximum_decon = working_space[i];
03225          if(maximum < working_space[6 * size_ext + i])
03226             maximum = working_space[6 * size_ext + i];
03227       }
03228 
03229       else
03230          working_space[i] = 0;
03231    }
03232    lda=1;
03233    if(lda>threshold)
03234       lda=threshold;
03235    lda=lda/100;
03236 
03237 //searching for peaks in deconvolved spectrum
03238    for(i = 1; i < size_ext - 1; i++){
03239       if(working_space[i] > working_space[i - 1] && working_space[i] > working_space[i + 1]){
03240          if(i >= shift && i < ssize + shift){
03241             if(working_space[i] > lda*maximum_decon && working_space[6 * size_ext + i] > threshold * maximum / 100.0){
03242                for(j = i - 1, a = 0, b = 0; j <= i + 1; j++){
03243                   a += (double)(j - shift) * working_space[j];
03244                   b += working_space[j];
03245                }
03246                a = a / b;
03247                if(a < 0)
03248                   a = 0;
03249 
03250                if(a >= ssize)
03251                   a = ssize - 1;
03252                if(peak_index == 0){
03253                   fPositionX[0] = a;
03254                   peak_index = 1;
03255                }
03256 
03257                else{
03258                   for(j = 0, priz = 0; j < peak_index && priz == 0; j++){
03259                      if(working_space[6 * size_ext + shift + (int)a] > working_space[6 * size_ext + shift + (int)fPositionX[j]])
03260                         priz = 1;
03261                   }
03262                   if(priz == 0){
03263                      if(j < fMaxPeaks){
03264                         fPositionX[j] = a;
03265                      }
03266                   }
03267 
03268                   else{
03269                      for(k = peak_index; k >= j; k--){
03270                         if(k < fMaxPeaks){
03271                            fPositionX[k] = fPositionX[k - 1];
03272                         }
03273                      }
03274                      fPositionX[j - 1] = a;
03275                   }
03276                   if(peak_index < fMaxPeaks)
03277                      peak_index += 1;
03278                }
03279             }
03280          }
03281       }
03282    }
03283 
03284    for(i = 0; i < ssize; i++) destVector[i] = working_space[i + shift];
03285    delete [] working_space;
03286    fNPeaks = peak_index;
03287    if(peak_index == fMaxPeaks)
03288       Warning("SearchHighRes", "Peak buffer full");
03289    return fNPeaks;
03290 }
03291 
03292 
03293 //______________________________________________________________________________
03294 Int_t TSpectrum::Search1HighRes(float *source,float *destVector, int ssize,
03295                                      float sigma, double threshold,
03296                                      bool backgroundRemove,int deconIterations,
03297                                      bool markov, int averWindow)
03298 {
03299    /* Begin_Html
03300    Old name of SearcHighRes introduced for back compatibility.
03301    This function will be removed after the June 2006 release
03302    End_Html */
03303 
03304    return SearchHighRes(source,destVector,ssize,sigma,threshold,backgroundRemove,
03305                         deconIterations,markov,averWindow);
03306 }
03307 
03308 
03309 //______________________________________________________________________________
03310 Int_t TSpectrum::StaticSearch(const TH1 *hist, Double_t sigma, Option_t *option, Double_t threshold)
03311 {
03312    /* Begin_Html
03313    Static function, interface to TSpectrum::Search.
03314    End_Html */
03315 
03316    TSpectrum s;
03317    return s.Search(hist,sigma,option,threshold);
03318 }
03319 
03320 
03321 //______________________________________________________________________________
03322 TH1 *TSpectrum::StaticBackground(const TH1 *hist,Int_t niter, Option_t *option)
03323 {
03324    /* Begin_Html
03325    Static function, interface to TSpectrum::Background.
03326    End_Html */
03327 
03328    TSpectrum s;
03329    return s.Background(hist,niter,option);
03330 }

Generated on Tue Jul 5 14:24:38 2011 for ROOT_528-00b_version by  doxygen 1.5.1