ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TSpectrum.cxx
Go to the documentation of this file.
1 // @(#)root/spectrum:$Id$
2 // Author: Miroslav Morhac 27/05/99
3 
4 #include "TSpectrum.h"
5 #include "TPolyMarker.h"
6 #include "TVirtualPad.h"
7 #include "TList.h"
8 #include "TH1.h"
9 #include "TMath.h"
10 
11 /** \class TSpectrum
12  \ingroup Spectrum
13  \brief Advanced Spectra Processing
14  \author Miroslav Morhac
15 
16  This class contains advanced spectra processing functions for:
17 
18  - One-dimensional background estimation
19  - One-dimensional smoothing
20  - One-dimensional deconvolution
21  - One-dimensional peak search
22 
23  The algorithms in this class have been published in the following references:
24 
25  1. M.Morhac et al.: Background elimination methods for multidimensional coincidence gamma-ray spectra. Nuclear Instruments and Methods in Physics Research A 401 (1997) 113-132.
26  2. M.Morhac et al.: Efficient one- and two-dimensional Gold deconvolution and its application to gamma-ray spectra decomposition. Nuclear Instruments and Methods in Physics Research A 401 (1997) 385-408.
27  3. M.Morhac et al.: Identification of peaks in multidimensional coincidence gamma-ray spectra. Nuclear Instruments and Methods in Research Physics A 443(2000), 108-125.
28 
29  These NIM papers are also available as doc or ps files from:
30 
31  - [Spectrum.doc](ftp://root.cern.ch/root/Spectrum.doc)
32  - [SpectrumDec.ps.gz](ftp://root.cern.ch/root/SpectrumDec.ps.gz)
33  - [SpectrumSrc.ps.gz](ftp://root.cern.ch/root/SpectrumSrc.ps.gz)
34  - [SpectrumBck.ps.gz](ftp://root.cern.ch/root/SpectrumBck.ps.gz)
35 */
36 
39 
40 #define PEAK_WINDOW 1024
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 /// Constructor.
45 
46 TSpectrum::TSpectrum() :TNamed("Spectrum", "Miroslav Morhac peak finder")
47 {
48  Int_t n = 100;
49  fMaxPeaks = n;
50  fPosition = new Double_t[n];
51  fPositionX = new Double_t[n];
52  fPositionY = new Double_t[n];
53  fResolution = 1;
54  fHistogram = 0;
55  fNPeaks = 0;
56 }
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 /// - maxpositions: maximum number of peaks
60 /// - resolution: determines resolution of the neighbouring peaks
61 /// default value is 1 correspond to 3 sigma distance
62 /// between peaks. Higher values allow higher resolution
63 /// (smaller distance between peaks.
64 /// May be set later through SetResolution.
65 
66 TSpectrum::TSpectrum(Int_t maxpositions, Double_t resolution)
67  :TNamed("Spectrum", "Miroslav Morhac peak finder")
68 {
69  Int_t n = maxpositions;
70  if (n <= 0) n = 1;
71  fMaxPeaks = n;
72  fPosition = new Double_t[n];
73  fPositionX = new Double_t[n];
74  fPositionY = new Double_t[n];
75  fHistogram = 0;
76  fNPeaks = 0;
77  SetResolution(resolution);
78 }
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// Destructor.
82 
84 {
85  delete [] fPosition;
86  delete [] fPositionX;
87  delete [] fPositionY;
88  delete fHistogram;
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Static function: Set average window of searched peaks
93 /// (see TSpectrum::SearchHighRes).
94 
96 {
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 /// Static function: Set max number of decon iterations in deconvolution
102 /// operation (see TSpectrum::SearchHighRes).
103 
105 {
106  fgIterations = n;
107 }
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// One-dimensional background estimation function.
111 ///
112 /// This function calculates the background spectrum in the input histogram h.
113 /// The background is returned as a histogram.
114 ///
115 /// #### Parameters:
116 ///
117 /// - h: input 1-d histogram
118 /// - numberIterations, (default value = 20)
119 /// Increasing numberIterations make the result smoother and lower.
120 /// - option: may contain one of the following options:
121 ///
122 /// - to set the direction parameter
123 /// "BackIncreasingWindow". By default the direction is BackDecreasingWindow
124 /// - filterOrder-order of clipping filter, (default "BackOrder2")
125 /// -possible values= "BackOrder4"
126 /// "BackOrder6"
127 /// "BackOrder8"
128 /// - "nosmoothing"- if selected, the background is not smoothed
129 /// By default the background is smoothed.
130 /// - smoothWindow-width of smoothing window, (default is "BackSmoothing3")
131 /// -possible values= "BackSmoothing5"
132 /// "BackSmoothing7"
133 /// "BackSmoothing9"
134 /// "BackSmoothing11"
135 /// "BackSmoothing13"
136 /// "BackSmoothing15"
137 /// - "Compton" if selected the estimation of Compton edge
138 /// will be included.
139 /// - "same" : if this option is specified, the resulting background
140 /// histogram is superimposed on the picture in the current pad.
141 ///
142 /// NOTE that the background is only evaluated in the current range of h.
143 /// ie, if h has a bin range (set via `h->GetXaxis()->SetRange(binmin,binmax)`,
144 /// the returned histogram will be created with the same number of bins
145 /// as the input histogram h, but only bins from `binmin` to `binmax` will be filled
146 /// with the estimated background.
147 
148 TH1 *TSpectrum::Background(const TH1 * h, int numberIterations,
149  Option_t * option)
150 {
151  if (h == 0) return 0;
152  Int_t dimension = h->GetDimension();
153  if (dimension > 1) {
154  Error("Search", "Only implemented for 1-d histograms");
155  return 0;
156  }
157  TString opt = option;
158  opt.ToLower();
159 
160  //set options
161  Int_t direction = kBackDecreasingWindow;
162  if (opt.Contains("backincreasingwindow")) direction = kBackIncreasingWindow;
163  Int_t filterOrder = kBackOrder2;
164  if (opt.Contains("backorder4")) filterOrder = kBackOrder4;
165  if (opt.Contains("backorder6")) filterOrder = kBackOrder6;
166  if (opt.Contains("backorder8")) filterOrder = kBackOrder8;
167  Bool_t smoothing = kTRUE;
168  if (opt.Contains("nosmoothing")) smoothing = kFALSE;
169  Int_t smoothWindow = kBackSmoothing3;
170  if (opt.Contains("backsmoothing5")) smoothWindow = kBackSmoothing5;
171  if (opt.Contains("backsmoothing7")) smoothWindow = kBackSmoothing7;
172  if (opt.Contains("backsmoothing9")) smoothWindow = kBackSmoothing9;
173  if (opt.Contains("backsmoothing11")) smoothWindow = kBackSmoothing11;
174  if (opt.Contains("backsmoothing13")) smoothWindow = kBackSmoothing13;
175  if (opt.Contains("backsmoothing15")) smoothWindow = kBackSmoothing15;
176  Bool_t compton = kFALSE;
177  if (opt.Contains("compton")) compton = kTRUE;
178 
179  Int_t first = h->GetXaxis()->GetFirst();
180  Int_t last = h->GetXaxis()->GetLast();
181  Int_t size = last-first+1;
182  Int_t i;
183  Double_t * source = new Double_t[size];
184  for (i = 0; i < size; i++) source[i] = h->GetBinContent(i + first);
185 
186  //find background (source is input and in output contains the background
187  Background(source,size,numberIterations, direction, filterOrder,smoothing,
188  smoothWindow,compton);
189 
190  //create output histogram containing background
191  //only bins in the range of the input histogram are filled
192  Int_t nch = strlen(h->GetName());
193  char *hbname = new char[nch+20];
194  snprintf(hbname,nch+20,"%s_background",h->GetName());
195  TH1 *hb = (TH1*)h->Clone(hbname);
196  hb->Reset();
197  hb->GetListOfFunctions()->Delete();
198  hb->SetLineColor(2);
199  for (i=0; i< size; i++) hb->SetBinContent(i+first,source[i]);
200  hb->SetEntries(size);
201 
202  //if option "same is specified, draw the result in the pad
203  if (opt.Contains("same")) {
204  if (gPad) delete gPad->GetPrimitive(hbname);
205  hb->Draw("same");
206  }
207  delete [] source;
208  delete [] hbname;
209  return hb;
210 }
211 
212 ////////////////////////////////////////////////////////////////////////////////
213 /// Print the array of positions.
214 
216 {
217  printf("\nNumber of positions = %d\n",fNPeaks);
218  for (Int_t i=0;i<fNPeaks;i++) {
219  printf(" x[%d] = %g, y[%d] = %g\n",i,fPositionX[i],i,fPositionY[i]);
220  }
221 }
222 
223 ////////////////////////////////////////////////////////////////////////////////
224 /// One-dimensional peak search function
225 ///
226 /// This function searches for peaks in source spectrum in hin
227 /// The number of found peaks and their positions are written into
228 /// the members fNpeaks and fPositionX.
229 /// The search is performed in the current histogram range.
230 ///
231 /// #### Parameters:
232 ///
233 /// - hin: pointer to the histogram of source spectrum
234 /// - sigma: sigma of searched peaks, for details we refer to manual
235 /// - threshold: (default=0.05) peaks with amplitude less than
236 /// threshold*highest_peak are discarded. 0<threshold<1
237 ///
238 /// By default, the background is removed before deconvolution.
239 /// Specify the option "nobackground" to not remove the background.
240 ///
241 /// By default the "Markov" chain algorithm is used.
242 /// Specify the option "noMarkov" to disable this algorithm
243 /// Note that by default the source spectrum is replaced by a new spectrum
244 ///
245 /// By default a polymarker object is created and added to the list of
246 /// functions of the histogram. The histogram is drawn with the specified
247 /// option and the polymarker object drawn on top of the histogram.
248 /// The polymarker coordinates correspond to the npeaks peaks found in
249 /// the histogram.
250 ///
251 /// A pointer to the polymarker object can be retrieved later via:
252 /// ~~~ {.cpp}
253 /// TList *functions = hin->GetListOfFunctions();
254 /// TPolyMarker *pm = (TPolyMarker*)functions->FindObject("TPolyMarker");
255 /// ~~~
256 /// Specify the option "goff" to disable the storage and drawing of the
257 /// polymarker.
258 ///
259 /// To disable the final drawing of the histogram with the search results (in case
260 /// you want to draw it yourself) specify "nodraw" in the options parameter.
261 
263  Double_t threshold)
264 {
265  if (hin == 0) return 0;
266  Int_t dimension = hin->GetDimension();
267  if (dimension > 2) {
268  Error("Search", "Only implemented for 1-d and 2-d histograms");
269  return 0;
270  }
271  if (threshold <=0 || threshold >= 1) {
272  Warning("Search","threshold must 0<threshold<1, threshold=0.05 assumed");
273  threshold = 0.05;
274  }
275  TString opt = option;
276  opt.ToLower();
278  if (opt.Contains("nobackground")) {
279  background = kFALSE;
280  opt.ReplaceAll("nobackground","");
281  }
282  Bool_t markov = kTRUE;
283  if (opt.Contains("nomarkov")) {
284  markov = kFALSE;
285  opt.ReplaceAll("nomarkov","");
286  }
287  Bool_t draw = kTRUE;
288  if (opt.Contains("nodraw")) {
289  draw = kFALSE;
290  opt.ReplaceAll("nodraw","");
291  }
292  if (dimension == 1) {
293  Int_t first = hin->GetXaxis()->GetFirst();
294  Int_t last = hin->GetXaxis()->GetLast();
295  Int_t size = last-first+1;
296  Int_t i, bin, npeaks;
297  Double_t * source = new Double_t[size];
298  Double_t * dest = new Double_t[size];
299  for (i = 0; i < size; i++) source[i] = hin->GetBinContent(i + first);
300  if (sigma < 1) {
301  sigma = size/fMaxPeaks;
302  if (sigma < 1) sigma = 1;
303  if (sigma > 8) sigma = 8;
304  }
305  npeaks = SearchHighRes(source, dest, size, sigma, 100*threshold,
306  background, fgIterations, markov, fgAverageWindow);
307 
308  for (i = 0; i < npeaks; i++) {
309  bin = first + Int_t(fPositionX[i] + 0.5);
310  fPositionX[i] = hin->GetBinCenter(bin);
311  fPositionY[i] = hin->GetBinContent(bin);
312  }
313  delete [] source;
314  delete [] dest;
315 
316  if (opt.Contains("goff"))
317  return npeaks;
318  if (!npeaks) return 0;
319  TPolyMarker * pm =
320  (TPolyMarker*)hin->GetListOfFunctions()->FindObject("TPolyMarker");
321  if (pm) {
322  hin->GetListOfFunctions()->Remove(pm);
323  delete pm;
324  }
325  pm = new TPolyMarker(npeaks, fPositionX, fPositionY);
326  hin->GetListOfFunctions()->Add(pm);
327  pm->SetMarkerStyle(23);
328  pm->SetMarkerColor(kRed);
329  pm->SetMarkerSize(1.3);
330  opt.ReplaceAll(" ","");
331  opt.ReplaceAll(",","");
332  if (draw)
333  ((TH1*)hin)->Draw(opt.Data());
334  return npeaks;
335  }
336  return 0;
337 }
338 
339 ////////////////////////////////////////////////////////////////////////////////
340 /// resolution: determines resolution of the neighbouring peaks
341 /// default value is 1 correspond to 3 sigma distance
342 /// between peaks. Higher values allow higher resolution
343 /// (smaller distance between peaks.
344 /// May be set later through SetResolution.
345 
347 {
348  if (resolution > 1)
349  fResolution = resolution;
350  else
351  fResolution = 1;
352 }
353 
354 ////////////////////////////////////////////////////////////////////////////////
355 /// This function calculates background spectrum from source spectrum.
356 /// The result is placed in the vector pointed by spe1945ctrum pointer.
357 /// The goal is to separate the useful information (peaks) from useless
358 /// information (background).
359 ///
360 /// - method is based on Sensitive Nonlinear Iterative Peak (SNIP) clipping
361 /// algorithm.
362 /// - new value in the channel "i" is calculated
363 ///
364 /// \f[
365 /// v_p(i) = min \left\{ v_{p-1}(i)^{\frac{\left[v_{p-1}(i+p)+v_{p-1}(i-p)\right]}{2}} \right\}
366 /// \f]
367 ///
368 /// where p = 1, 2, ..., numberIterations. In fact it represents second order
369 /// difference filter (-1,2,-1).
370 ///
371 /// One can also change the
372 /// direction of the change of the clipping window, the order of the clipping
373 /// filter, to include smoothing, to set width of smoothing window and to include
374 /// the estimation of Compton edges. On successful completion it returns 0. On
375 /// error it returns pointer to the string describing error.
376 ///
377 /// #### Parameters:
378 ///
379 /// - spectrum: pointer to the vector of source spectrum
380 /// - ssize: length of the spectrum vector
381 /// - numberIterations: maximal width of clipping window,
382 /// - direction: direction of change of clipping window.
383 /// Possible values: kBackIncreasingWindow, kBackDecreasingWindow
384 /// - filterOrder: order of clipping filter.
385 /// Possible values: kBackOrder2, kBackOrder4, kBackOrder6, kBackOrder8
386 /// - smoothing: logical variable whether the smoothing operation in the
387 /// estimation of background will be included.
388 /// Possible values: kFALSE, kTRUE
389 /// - smoothWindow: width of smoothing window.
390 /// Possible values: kBackSmoothing3, kBackSmoothing5, kBackSmoothing7,
391 /// kBackSmoothing9, kBackSmoothing11, kBackSmoothing13, kBackSmoothing15.
392 /// - compton: logical variable whether the estimation of Compton edge will be
393 /// included. Possible values: kFALSE, kTRUE.
394 ///
395 /// #### References:
396 ///
397 /// 1. C. G Ryan et al.: SNIP, a statistics-sensitive background treatment for the
398 /// quantitative analysis of PIXE spectra in geoscience applications. NIM, B34
399 /// (1988), 396-402.
400 ///
401 /// 2. M. Morhac;, J. Kliman, V. Matouoek, M. Veselsky, I. Turzo:
402 /// Background elimination methods for multidimensional gamma-ray spectra. NIM,
403 /// A401 (1997) 113-132.
404 ///
405 /// 3. D. D. Burgess, R. J. Tervo: Background estimation for gamma-ray
406 /// spectroscopy. NIM 214 (1983), 431-434.
407 ///
408 /// ### Example 1 script Background_incr.c:
409 ///
410 /// \image html TSpectrum_Background_incr.jpg Fig. 1 Example of the estimation of background for number of iterations=6. Original spectrum is shown in black color, estimated background in red color.
411 ///
412 /// #### Script:
413 ///
414 /// ~~~ {.cpp}
415 /// // Example to illustrate the background estimator (class TSpectrum).
416 /// // To execute this example, do
417 /// // root > .x Background_incr.C
418 ///
419 /// #include <TSpectrum>
420 ///
421 /// void Background_incr() {
422 /// Int_t i;
423 /// Double_t nbins = 256;
424 /// Double_t xmin = 0;
425 /// Double_t xmax = nbins;
426 /// Double_t * source = new Double_t[nbins];
427 /// TH1F *back = new TH1F("back","",nbins,xmin,xmax);
428 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
429 /// TFile *f = new TFile("spectra/TSpectrum.root");
430 /// back=(TH1F*) f->Get("back1;1");
431 /// TCanvas *Background = gROOT->GetListOfCanvases()->FindObject("Background");
432 /// if (!Background) Background =
433 /// new TCanvas("Background",
434 /// "Estimation of background with increasing window",
435 /// 10,10,1000,700);
436 /// back->Draw("L");
437 /// TSpectrum *s = new TSpectrum();
438 /// for (i = 0; i < nbins; i++) source[i]=back->GetBinContent(i + 1);
439 /// s->Background(source,nbins,6,kBackIncreasingWindow,kBackOrder2,kFALSE,
440 /// kBackSmoothing3,kFALSE);
441 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
442 /// d->SetLineColor(kRed);
443 /// d->Draw("SAME L");
444 /// }
445 /// ~~~
446 ///
447 /// ### Example 2 script Background_decr.c:
448 ///
449 /// In Fig. 1. one can notice that at the edges of the peaks the estimated
450 /// background goes under the peaks. An alternative approach is to decrease the
451 /// clipping window from a given value numberIterations to the value of one, which
452 /// is presented in this example.
453 ///
454 /// \image html TSpectrum_Background_decr.jpg Fig. 2 Example of the estimation of background for numberIterations=6 using decreasing clipping window algorithm. Original spectrum is shown in black color, estimated background in red color.
455 ///
456 /// #### Script:
457 ///
458 /// ~~~ {.cpp}
459 /// // Example to illustrate the background estimator (class TSpectrum).
460 /// // To execute this example, do
461 /// // root > .x Background_decr.C
462 ///
463 /// #include <TSpectrum>
464 ///
465 /// void Background_decr() {
466 /// Int_t i;
467 /// Double_t nbins = 256;
468 /// Double_t xmin = 0;
469 /// Double_t xmax = nbins;
470 /// Double_t * source = new Double_t[nbins];
471 /// TH1F *back = new TH1F("back","",nbins,xmin,xmax);
472 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
473 /// TFile *f = new TFile("spectra/TSpectrum.root");
474 /// back=(TH1F*) f->Get("back1;1");
475 /// TCanvas *Background = gROOT->GetListOfCanvases()->FindObject("Background");
476 /// if (!Background) Background =
477 /// new TCanvas("Background","Estimation of background with decreasing window",
478 /// 10,10,1000,700);
479 /// back->Draw("L");
480 /// TSpectrum *s = new TSpectrum();
481 /// for (i = 0; i < nbins; i++) source[i]=back->GetBinContent(i + 1);
482 /// s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kFALSE,
483 /// kBackSmoothing3,kFALSE);
484 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
485 /// d->SetLineColor(kRed);
486 /// d->Draw("SAME L");
487 /// }
488 /// ~~~
489 ///
490 /// ### Example 3 script Background_width.c:
491 ///
492 /// The question is how to choose the width of the clipping window, i.e.,
493 /// numberIterations parameter. The influence of this parameter on the estimated
494 /// background is illustrated in Fig. 3.
495 ///
496 /// \image html TSpectrum_Background_width.jpg Fig. 3 Example of the influence of clipping window width on the estimated background for numberIterations=4 (red line), 6 (blue line) 8 (green line) using decreasing clipping window algorithm.
497 ///
498 /// in general one should set this parameter so that the value
499 /// 2*numberIterations+1 was greater than the widths of preserved objects (peaks).
500 ///
501 /// #### Script:
502 ///
503 /// ~~~ {.cpp}
504 /// // Example to illustrate the influence of the clipping window width on the
505 /// // estimated background. To execute this example, do:
506 /// // root > .x Background_width.C
507 ///
508 /// #include <TSpectrum>
509 ///
510 /// void Background_width() {
511 /// Int_t i;
512 /// Double_t nbins = 256;
513 /// Double_t xmin = 0;
514 /// Double_t xmax = nbins;
515 /// Double_t * source = new Double_t[nbins];
516 /// TH1F *h = new TH1F("h","",nbins,xmin,xmax);
517 /// TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
518 /// TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
519 /// TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
520 /// TFile *f = new TFile("spectra/TSpectrum.root");
521 /// h=(TH1F*) f->Get("back1;1");
522 /// TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
523 /// if (!background) background = new TCanvas("background",
524 /// "Influence of clipping window width on the estimated background",
525 /// 10,10,1000,700);
526 /// h->Draw("L");
527 /// TSpectrum *s = new TSpectrum();
528 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
529 /// s->Background(source,nbins,4,kBackDecreasingWindow,kBackOrder2,kFALSE,
530 /// kBackSmoothing3,kFALSE);
531 /// for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
532 /// d1->SetLineColor(kRed);
533 /// d1->Draw("SAME L");
534 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
535 /// s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kFALSE,
536 /// kBackSmoothing3,kFALSE);
537 /// for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
538 /// d2->SetLineColor(kBlue);
539 /// d2->Draw("SAME L");
540 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
541 /// s->Background(source,nbins,8,kBackDecreasingWindow,kBackOrder2,kFALSE,
542 /// kBackSmoothing3,kFALSE);
543 /// for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,source[i]);
544 /// d3->SetLineColor(kGreen);
545 /// d3->Draw("SAME L");
546 /// }
547 /// ~~~
548 ///
549 /// ### Example 4 script Background_width2.c:
550 ///
551 /// another example for very complex spectrum is given in Fig. 4.
552 ///
553 /// \image html TSpectrum_Background_width2.jpg Fig. 4 Example of the influence of clipping window width on the estimated background for numberIterations=10 (red line), 20 (blue line), 30 (green line) and 40 (magenta line) using decreasing clipping window algorithm.
554 ///
555 /// #### Script:
556 ///
557 /// ~~~ {.cpp}
558 /// // Example to illustrate the influence of the clipping window width on the
559 /// // estimated background. To execute this example, do:
560 /// // root > .x Background_width2.C
561 ///
562 /// #include <TSpectrum>
563 ///
564 /// void Background_width2() {
565 /// Int_t i;
566 /// Double_t nbins = 4096;
567 /// Double_t xmin = 0;
568 /// Double_t xmax = 4096;
569 /// Double_t * source = new Double_t[nbins];
570 /// TH1F *h = new TH1F("h","",nbins,xmin,xmax);
571 /// TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
572 /// TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
573 /// TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
574 /// TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
575 /// TFile *f = new TFile("spectra/TSpectrum.root");
576 /// h=(TH1F*) f->Get("back2;1");
577 /// TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
578 /// if (!background) background = new TCanvas("background",
579 /// "Influence of clipping window width on the estimated background",
580 /// 10,10,1000,700);
581 /// h->SetAxisRange(0,1000);
582 /// h->SetMaximum(20000);
583 /// h->Draw("L");
584 /// TSpectrum *s = new TSpectrum();
585 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
586 /// s->Background(source,nbins,10,kBackDecreasingWindow,kBackOrder2,kFALSE,
587 /// kBackSmoothing3,kFALSE);
588 /// for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
589 /// d1->SetLineColor(kRed);
590 /// d1->Draw("SAME L");
591 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
592 /// s->Background(source,nbins,20,kBackDecreasingWindow,kBackOrder2,kFALSE,
593 /// kBackSmoothing3,kFALSE);
594 /// for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
595 /// d2->SetLineColor(kBlue);
596 /// d2->Draw("SAME L");
597 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
598 /// s->Background(source,nbins,30,kBackDecreasingWindow,kBackOrder2,kFALSE,
599 /// kBackSmoothing3,kFALSE);
600 /// for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,source[i]);
601 /// d3->SetLineColor(kGreen);
602 /// d3->Draw("SAME L");
603 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
604 /// s->Background(source,nbins,10,kBackDecreasingWindow,kBackOrder2,kFALSE,
605 /// kBackSmoothing3,kFALSE);
606 /// for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,source[i]);
607 /// d4->SetLineColor(kMagenta);
608 /// d4->Draw("SAME L");
609 /// }
610 /// ~~~
611 ///
612 /// ### Example 5 script Background_order.c:
613 ///
614 /// Second order difference filter removes linear (quasi-linear) background and
615 /// preserves symmetrical peaks. However if the shape of the background is more
616 /// complex one can employ higher-order clipping filters (see example in Fig. 5)
617 ///
618 /// \image html TSpectrum_Background_order.jpg Fig. 5 Example of the influence of clipping filter difference order on the estimated background for fNnumberIterations=40, 2-nd order red line, 4-th order blue line, 6-th order green line and 8-th order magenta line, and using decreasing clipping window algorithm.
619 ///
620 /// #### Script:
621 ///
622 /// ~~~ {.cpp}
623 /// // Example to illustrate the influence of the clipping filter difference order
624 /// // on the estimated background. To execute this example, do
625 /// // root > .x Background_order.C
626 ///
627 /// #include <TSpectrum>
628 ///
629 /// void Background_order() {
630 /// Int_t i;
631 /// Double_t nbins = 4096;
632 /// Double_t xmin = 0;
633 /// Double_t xmax = 4096;
634 /// Double_t * source = new Double_t[nbins];
635 /// TH1F *h = new TH1F("h","",nbins,xmin,xmax);
636 /// TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
637 /// TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
638 /// TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
639 /// TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
640 /// TFile *f = new TFile("spectra/TSpectrum.root");
641 /// h=(TH1F*) f->Get("back2;1");
642 /// TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
643 /// if (!background) background = new TCanvas("background",
644 /// "Influence of clipping filter difference order on the estimated background",
645 /// 10,10,1000,700);
646 /// h->SetAxisRange(1220,1460);
647 /// h->SetMaximum(11000);
648 /// h->Draw("L");
649 /// TSpectrum *s = new TSpectrum();
650 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
651 /// s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder2,kFALSE,
652 /// kBackSmoothing3,kFALSE);
653 /// for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
654 /// d1->SetLineColor(kRed);
655 /// d1->Draw("SAME L");
656 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
657 /// s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder4,kFALSE,
658 /// kBackSmoothing3,kFALSE);
659 /// for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
660 /// d2->SetLineColor(kBlue);
661 /// d2->Draw("SAME L");
662 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
663 /// s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder6,kFALSE,
664 /// kBackSmoothing3,kFALSE);
665 /// for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,source[i]);
666 /// d3->SetLineColor(kGreen);
667 /// d3->Draw("SAME L");
668 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
669 /// s->Background(source,nbins,40,kBackDecreasingWindow,kBackOrder8,kFALSE,
670 /// kBackSmoothing3,kFALSE);
671 /// for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,source[i]);
672 /// d4->SetLineColor(kMagenta);
673 /// d4->Draw("SAME L");
674 /// }
675 /// ~~~
676 ///
677 /// ### Example 6 script Background_smooth.c:
678 ///
679 /// The estimate of the background can be influenced by noise present in the
680 /// spectrum. We proposed the algorithm of the background estimate with
681 /// simultaneous smoothing. In the original algorithm without smoothing, the
682 /// estimated background snatches the lower spikes in the noise. Consequently,
683 /// the areas of peaks are biased by this error.
684 ///
685 /// \image html TSpectrum_Background_smooth1.jpg Fig. 7 Principle of background estimation algorithm with simultaneous smoothing.
686 /// \image html TSpectrum_Background_smooth2.jpg Fig. 8 Illustration of non-smoothing (red line) and smoothing algorithm of background estimation (blue line).
687 ///
688 /// #### Script:
689 ///
690 /// ~~~ {.cpp}
691 /// // Example to illustrate the background estimator (class TSpectrum) including
692 /// // Compton edges. To execute this example, do:
693 /// // root > .x Background_smooth.C
694 ///
695 /// #include <TSpectrum>
696 ///
697 /// void Background_smooth() {
698 /// Int_t i;
699 /// Double_t nbins = 4096;
700 /// Double_t xmin = 0;
701 /// Double_t xmax = nbins;
702 /// Double_t * source = new Double_t[nbins];
703 /// TH1F *h = new TH1F("h","",nbins,xmin,xmax);
704 /// TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
705 /// TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
706 /// TFile *f = new TFile("spectra/TSpectrum.root");
707 /// h=(TH1F*) f->Get("back4;1");
708 /// TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
709 /// if (!background) background = new TCanvas("background",
710 /// "Estimation of background with noise",10,10,1000,700);
711 /// h->SetAxisRange(3460,3830);
712 /// h->Draw("L");
713 /// TSpectrum *s = new TSpectrum();
714 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
715 /// s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kFALSE,
716 /// kBackSmoothing3,kFALSE);
717 /// for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
718 /// d1->SetLineColor(kRed);
719 /// d1->Draw("SAME L");
720 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
721 /// s->Background(source,nbins,6,kBackDecreasingWindow,kBackOrder2,kTRUE,
722 /// kBackSmoothing3,kFALSE);
723 /// for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,source[i]);
724 /// d2->SetLineColor(kBlue);
725 /// d2->Draw("SAME L");
726 /// }
727 /// ~~~
728 ///
729 /// ### Example 8 script Background_compton.c:
730 ///
731 /// Sometimes it is necessary to include also the Compton edges into the estimate of
732 /// the background. In Fig. 8 we present the example of the synthetic spectrum
733 /// with Compton edges. The background was estimated using the 8-th order filter
734 /// with the estimation of the Compton edges using decreasing
735 /// clipping window algorithm (numberIterations=10) with smoothing
736 /// (smoothingWindow=5).
737 ///
738 /// \image html TSpectrum_Background_compton.jpg Fig. 8 Example of the estimate of the background with Compton edges (red line) for numberIterations=10, 8-th order difference filter, using decreasing clipping window algorithm and smoothing (smoothingWindow=5).
739 ///
740 /// #### Script:
741 ///
742 /// ~~~ {.cpp}
743 /// // Example to illustrate the background estimator (class TSpectrum) including
744 /// // Compton edges. To execute this example, do:
745 /// // root > .x Background_compton.C
746 ///
747 /// #include <TSpectrum>
748 ///
749 /// void Background_compton() {
750 /// Int_t i;
751 /// Double_t nbins = 512;
752 /// Double_t xmin = 0;
753 /// Double_t xmax = nbins;
754 /// Double_t * source = new Double_t[nbins];
755 /// TH1F *h = new TH1F("h","",nbins,xmin,xmax);
756 /// TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
757 /// TFile *f = new TFile("spectra/TSpectrum.root");
758 /// h=(TH1F*) f->Get("back3;1");
759 /// TCanvas *background = gROOT->GetListOfCanvases()->FindObject("background");
760 /// if (!background) background = new TCanvas("background",
761 /// "Estimation of background with Compton edges under peaks",10,10,1000,700);
762 /// h->Draw("L");
763 /// TSpectrum *s = new TSpectrum();
764 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
765 /// s->Background(source,nbins,10,kBackDecreasingWindow,kBackOrder8,kTRUE,
766 /// kBackSmoothing5,,kTRUE);
767 /// for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,source[i]);
768 /// d1->SetLineColor(kRed);
769 /// d1->Draw("SAME L");
770 /// }
771 /// ~~~
772 
773 const char *TSpectrum::Background(Double_t *spectrum, int ssize,
774  int numberIterations,
775  int direction, int filterOrder,
776  bool smoothing,int smoothWindow,
777  bool compton)
778 {
779  int i, j, w, bw, b1, b2, priz;
780  Double_t 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;
781  if (ssize <= 0)
782  return "Wrong Parameters";
783  if (numberIterations < 1)
784  return "Width of Clipping Window Must Be Positive";
785  if (ssize < 2 * numberIterations + 1)
786  return "Too Large Clipping Window";
787  if (smoothing == kTRUE && smoothWindow != kBackSmoothing3 && smoothWindow != kBackSmoothing5 && smoothWindow != kBackSmoothing7 && smoothWindow != kBackSmoothing9 && smoothWindow != kBackSmoothing11 && smoothWindow != kBackSmoothing13 && smoothWindow != kBackSmoothing15)
788  return "Incorrect width of smoothing window";
789  Double_t *working_space = new Double_t[2 * ssize];
790  for (i = 0; i < ssize; i++){
791  working_space[i] = spectrum[i];
792  working_space[i + ssize] = spectrum[i];
793  }
794  bw=(smoothWindow-1)/2;
795  if (direction == kBackIncreasingWindow)
796  i = 1;
797  else if(direction == kBackDecreasingWindow)
798  i = numberIterations;
799  if (filterOrder == kBackOrder2) {
800  do{
801  for (j = i; j < ssize - i; j++) {
802  if (smoothing == kFALSE){
803  a = working_space[ssize + j];
804  b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
805  if (b < a)
806  a = b;
807  working_space[j] = a;
808  }
809 
810  else if (smoothing == kTRUE){
811  a = working_space[ssize + j];
812  av = 0;
813  men = 0;
814  for (w = j - bw; w <= j + bw; w++){
815  if ( w >= 0 && w < ssize){
816  av += working_space[ssize + w];
817  men +=1;
818  }
819  }
820  av = av / men;
821  b = 0;
822  men = 0;
823  for (w = j - i - bw; w <= j - i + bw; w++){
824  if ( w >= 0 && w < ssize){
825  b += working_space[ssize + w];
826  men +=1;
827  }
828  }
829  b = b / men;
830  c = 0;
831  men = 0;
832  for (w = j + i - bw; w <= j + i + bw; w++){
833  if ( w >= 0 && w < ssize){
834  c += working_space[ssize + w];
835  men +=1;
836  }
837  }
838  c = c / men;
839  b = (b + c) / 2;
840  if (b < a)
841  av = b;
842  working_space[j]=av;
843  }
844  }
845  for (j = i; j < ssize - i; j++)
846  working_space[ssize + j] = working_space[j];
847  if (direction == kBackIncreasingWindow)
848  i+=1;
849  else if(direction == kBackDecreasingWindow)
850  i-=1;
851  }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
852  }
853 
854  else if (filterOrder == kBackOrder4) {
855  do{
856  for (j = i; j < ssize - i; j++) {
857  if (smoothing == kFALSE){
858  a = working_space[ssize + j];
859  b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
860  c = 0;
861  ai = i / 2;
862  c -= working_space[ssize + j - (Int_t) (2 * ai)] / 6;
863  c += 4 * working_space[ssize + j - (Int_t) ai] / 6;
864  c += 4 * working_space[ssize + j + (Int_t) ai] / 6;
865  c -= working_space[ssize + j + (Int_t) (2 * ai)] / 6;
866  if (b < c)
867  b = c;
868  if (b < a)
869  a = b;
870  working_space[j] = a;
871  }
872 
873  else if (smoothing == kTRUE){
874  a = working_space[ssize + j];
875  av = 0;
876  men = 0;
877  for (w = j - bw; w <= j + bw; w++){
878  if ( w >= 0 && w < ssize){
879  av += working_space[ssize + w];
880  men +=1;
881  }
882  }
883  av = av / men;
884  b = 0;
885  men = 0;
886  for (w = j - i - bw; w <= j - i + bw; w++){
887  if ( w >= 0 && w < ssize){
888  b += working_space[ssize + w];
889  men +=1;
890  }
891  }
892  b = b / men;
893  c = 0;
894  men = 0;
895  for (w = j + i - bw; w <= j + i + bw; w++){
896  if ( w >= 0 && w < ssize){
897  c += working_space[ssize + w];
898  men +=1;
899  }
900  }
901  c = c / men;
902  b = (b + c) / 2;
903  ai = i / 2;
904  b4 = 0, men = 0;
905  for (w = j - (Int_t)(2 * ai) - bw; w <= j - (Int_t)(2 * ai) + bw; w++){
906  if (w >= 0 && w < ssize){
907  b4 += working_space[ssize + w];
908  men +=1;
909  }
910  }
911  b4 = b4 / men;
912  c4 = 0, men = 0;
913  for (w = j - (Int_t)ai - bw; w <= j - (Int_t)ai + bw; w++){
914  if (w >= 0 && w < ssize){
915  c4 += working_space[ssize + w];
916  men +=1;
917  }
918  }
919  c4 = c4 / men;
920  d4 = 0, men = 0;
921  for (w = j + (Int_t)ai - bw; w <= j + (Int_t)ai + bw; w++){
922  if (w >= 0 && w < ssize){
923  d4 += working_space[ssize + w];
924  men +=1;
925  }
926  }
927  d4 = d4 / men;
928  e4 = 0, men = 0;
929  for (w = j + (Int_t)(2 * ai) - bw; w <= j + (Int_t)(2 * ai) + bw; w++){
930  if (w >= 0 && w < ssize){
931  e4 += working_space[ssize + w];
932  men +=1;
933  }
934  }
935  e4 = e4 / men;
936  b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
937  if (b < b4)
938  b = b4;
939  if (b < a)
940  av = b;
941  working_space[j]=av;
942  }
943  }
944  for (j = i; j < ssize - i; j++)
945  working_space[ssize + j] = working_space[j];
946  if (direction == kBackIncreasingWindow)
947  i+=1;
948  else if(direction == kBackDecreasingWindow)
949  i-=1;
950  }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
951  }
952 
953  else if (filterOrder == kBackOrder6) {
954  do{
955  for (j = i; j < ssize - i; j++) {
956  if (smoothing == kFALSE){
957  a = working_space[ssize + j];
958  b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
959  c = 0;
960  ai = i / 2;
961  c -= working_space[ssize + j - (Int_t) (2 * ai)] / 6;
962  c += 4 * working_space[ssize + j - (Int_t) ai] / 6;
963  c += 4 * working_space[ssize + j + (Int_t) ai] / 6;
964  c -= working_space[ssize + j + (Int_t) (2 * ai)] / 6;
965  d = 0;
966  ai = i / 3;
967  d += working_space[ssize + j - (Int_t) (3 * ai)] / 20;
968  d -= 6 * working_space[ssize + j - (Int_t) (2 * ai)] / 20;
969  d += 15 * working_space[ssize + j - (Int_t) ai] / 20;
970  d += 15 * working_space[ssize + j + (Int_t) ai] / 20;
971  d -= 6 * working_space[ssize + j + (Int_t) (2 * ai)] / 20;
972  d += working_space[ssize + j + (Int_t) (3 * ai)] / 20;
973  if (b < d)
974  b = d;
975  if (b < c)
976  b = c;
977  if (b < a)
978  a = b;
979  working_space[j] = a;
980  }
981 
982  else if (smoothing == kTRUE){
983  a = working_space[ssize + j];
984  av = 0;
985  men = 0;
986  for (w = j - bw; w <= j + bw; w++){
987  if ( w >= 0 && w < ssize){
988  av += working_space[ssize + w];
989  men +=1;
990  }
991  }
992  av = av / men;
993  b = 0;
994  men = 0;
995  for (w = j - i - bw; w <= j - i + bw; w++){
996  if ( w >= 0 && w < ssize){
997  b += working_space[ssize + w];
998  men +=1;
999  }
1000  }
1001  b = b / men;
1002  c = 0;
1003  men = 0;
1004  for (w = j + i - bw; w <= j + i + bw; w++){
1005  if ( w >= 0 && w < ssize){
1006  c += working_space[ssize + w];
1007  men +=1;
1008  }
1009  }
1010  c = c / men;
1011  b = (b + c) / 2;
1012  ai = i / 2;
1013  b4 = 0, men = 0;
1014  for (w = j - (Int_t)(2 * ai) - bw; w <= j - (Int_t)(2 * ai) + bw; w++){
1015  if (w >= 0 && w < ssize){
1016  b4 += working_space[ssize + w];
1017  men +=1;
1018  }
1019  }
1020  b4 = b4 / men;
1021  c4 = 0, men = 0;
1022  for (w = j - (Int_t)ai - bw; w <= j - (Int_t)ai + bw; w++){
1023  if (w >= 0 && w < ssize){
1024  c4 += working_space[ssize + w];
1025  men +=1;
1026  }
1027  }
1028  c4 = c4 / men;
1029  d4 = 0, men = 0;
1030  for (w = j + (Int_t)ai - bw; w <= j + (Int_t)ai + bw; w++){
1031  if (w >= 0 && w < ssize){
1032  d4 += working_space[ssize + w];
1033  men +=1;
1034  }
1035  }
1036  d4 = d4 / men;
1037  e4 = 0, men = 0;
1038  for (w = j + (Int_t)(2 * ai) - bw; w <= j + (Int_t)(2 * ai) + bw; w++){
1039  if (w >= 0 && w < ssize){
1040  e4 += working_space[ssize + w];
1041  men +=1;
1042  }
1043  }
1044  e4 = e4 / men;
1045  b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
1046  ai = i / 3;
1047  b6 = 0, men = 0;
1048  for (w = j - (Int_t)(3 * ai) - bw; w <= j - (Int_t)(3 * ai) + bw; w++){
1049  if (w >= 0 && w < ssize){
1050  b6 += working_space[ssize + w];
1051  men +=1;
1052  }
1053  }
1054  b6 = b6 / men;
1055  c6 = 0, men = 0;
1056  for (w = j - (Int_t)(2 * ai) - bw; w <= j - (Int_t)(2 * ai) + bw; w++){
1057  if (w >= 0 && w < ssize){
1058  c6 += working_space[ssize + w];
1059  men +=1;
1060  }
1061  }
1062  c6 = c6 / men;
1063  d6 = 0, men = 0;
1064  for (w = j - (Int_t)ai - bw; w <= j - (Int_t)ai + bw; w++){
1065  if (w >= 0 && w < ssize){
1066  d6 += working_space[ssize + w];
1067  men +=1;
1068  }
1069  }
1070  d6 = d6 / men;
1071  e6 = 0, men = 0;
1072  for (w = j + (Int_t)ai - bw; w <= j + (Int_t)ai + bw; w++){
1073  if (w >= 0 && w < ssize){
1074  e6 += working_space[ssize + w];
1075  men +=1;
1076  }
1077  }
1078  e6 = e6 / men;
1079  f6 = 0, men = 0;
1080  for (w = j + (Int_t)(2 * ai) - bw; w <= j + (Int_t)(2 * ai) + bw; w++){
1081  if (w >= 0 && w < ssize){
1082  f6 += working_space[ssize + w];
1083  men +=1;
1084  }
1085  }
1086  f6 = f6 / men;
1087  g6 = 0, men = 0;
1088  for (w = j + (Int_t)(3 * ai) - bw; w <= j + (Int_t)(3 * ai) + bw; w++){
1089  if (w >= 0 && w < ssize){
1090  g6 += working_space[ssize + w];
1091  men +=1;
1092  }
1093  }
1094  g6 = g6 / men;
1095  b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
1096  if (b < b6)
1097  b = b6;
1098  if (b < b4)
1099  b = b4;
1100  if (b < a)
1101  av = b;
1102  working_space[j]=av;
1103  }
1104  }
1105  for (j = i; j < ssize - i; j++)
1106  working_space[ssize + j] = working_space[j];
1107  if (direction == kBackIncreasingWindow)
1108  i+=1;
1109  else if(direction == kBackDecreasingWindow)
1110  i-=1;
1111  }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
1112  }
1113 
1114  else if (filterOrder == kBackOrder8) {
1115  do{
1116  for (j = i; j < ssize - i; j++) {
1117  if (smoothing == kFALSE){
1118  a = working_space[ssize + j];
1119  b = (working_space[ssize + j - i] + working_space[ssize + j + i]) / 2.0;
1120  c = 0;
1121  ai = i / 2;
1122  c -= working_space[ssize + j - (Int_t) (2 * ai)] / 6;
1123  c += 4 * working_space[ssize + j - (Int_t) ai] / 6;
1124  c += 4 * working_space[ssize + j + (Int_t) ai] / 6;
1125  c -= working_space[ssize + j + (Int_t) (2 * ai)] / 6;
1126  d = 0;
1127  ai = i / 3;
1128  d += working_space[ssize + j - (Int_t) (3 * ai)] / 20;
1129  d -= 6 * working_space[ssize + j - (Int_t) (2 * ai)] / 20;
1130  d += 15 * working_space[ssize + j - (Int_t) ai] / 20;
1131  d += 15 * working_space[ssize + j + (Int_t) ai] / 20;
1132  d -= 6 * working_space[ssize + j + (Int_t) (2 * ai)] / 20;
1133  d += working_space[ssize + j + (Int_t) (3 * ai)] / 20;
1134  e = 0;
1135  ai = i / 4;
1136  e -= working_space[ssize + j - (Int_t) (4 * ai)] / 70;
1137  e += 8 * working_space[ssize + j - (Int_t) (3 * ai)] / 70;
1138  e -= 28 * working_space[ssize + j - (Int_t) (2 * ai)] / 70;
1139  e += 56 * working_space[ssize + j - (Int_t) ai] / 70;
1140  e += 56 * working_space[ssize + j + (Int_t) ai] / 70;
1141  e -= 28 * working_space[ssize + j + (Int_t) (2 * ai)] / 70;
1142  e += 8 * working_space[ssize + j + (Int_t) (3 * ai)] / 70;
1143  e -= working_space[ssize + j + (Int_t) (4 * ai)] / 70;
1144  if (b < e)
1145  b = e;
1146  if (b < d)
1147  b = d;
1148  if (b < c)
1149  b = c;
1150  if (b < a)
1151  a = b;
1152  working_space[j] = a;
1153  }
1154 
1155  else if (smoothing == kTRUE){
1156  a = working_space[ssize + j];
1157  av = 0;
1158  men = 0;
1159  for (w = j - bw; w <= j + bw; w++){
1160  if ( w >= 0 && w < ssize){
1161  av += working_space[ssize + w];
1162  men +=1;
1163  }
1164  }
1165  av = av / men;
1166  b = 0;
1167  men = 0;
1168  for (w = j - i - bw; w <= j - i + bw; w++){
1169  if ( w >= 0 && w < ssize){
1170  b += working_space[ssize + w];
1171  men +=1;
1172  }
1173  }
1174  b = b / men;
1175  c = 0;
1176  men = 0;
1177  for (w = j + i - bw; w <= j + i + bw; w++){
1178  if ( w >= 0 && w < ssize){
1179  c += working_space[ssize + w];
1180  men +=1;
1181  }
1182  }
1183  c = c / men;
1184  b = (b + c) / 2;
1185  ai = i / 2;
1186  b4 = 0, men = 0;
1187  for (w = j - (Int_t)(2 * ai) - bw; w <= j - (Int_t)(2 * ai) + bw; w++){
1188  if (w >= 0 && w < ssize){
1189  b4 += working_space[ssize + w];
1190  men +=1;
1191  }
1192  }
1193  b4 = b4 / men;
1194  c4 = 0, men = 0;
1195  for (w = j - (Int_t)ai - bw; w <= j - (Int_t)ai + bw; w++){
1196  if (w >= 0 && w < ssize){
1197  c4 += working_space[ssize + w];
1198  men +=1;
1199  }
1200  }
1201  c4 = c4 / men;
1202  d4 = 0, men = 0;
1203  for (w = j + (Int_t)ai - bw; w <= j + (Int_t)ai + bw; w++){
1204  if (w >= 0 && w < ssize){
1205  d4 += working_space[ssize + w];
1206  men +=1;
1207  }
1208  }
1209  d4 = d4 / men;
1210  e4 = 0, men = 0;
1211  for (w = j + (Int_t)(2 * ai) - bw; w <= j + (Int_t)(2 * ai) + bw; w++){
1212  if (w >= 0 && w < ssize){
1213  e4 += working_space[ssize + w];
1214  men +=1;
1215  }
1216  }
1217  e4 = e4 / men;
1218  b4 = (-b4 + 4 * c4 + 4 * d4 - e4) / 6;
1219  ai = i / 3;
1220  b6 = 0, men = 0;
1221  for (w = j - (Int_t)(3 * ai) - bw; w <= j - (Int_t)(3 * ai) + bw; w++){
1222  if (w >= 0 && w < ssize){
1223  b6 += working_space[ssize + w];
1224  men +=1;
1225  }
1226  }
1227  b6 = b6 / men;
1228  c6 = 0, men = 0;
1229  for (w = j - (Int_t)(2 * ai) - bw; w <= j - (Int_t)(2 * ai) + bw; w++){
1230  if (w >= 0 && w < ssize){
1231  c6 += working_space[ssize + w];
1232  men +=1;
1233  }
1234  }
1235  c6 = c6 / men;
1236  d6 = 0, men = 0;
1237  for (w = j - (Int_t)ai - bw; w <= j - (Int_t)ai + bw; w++){
1238  if (w >= 0 && w < ssize){
1239  d6 += working_space[ssize + w];
1240  men +=1;
1241  }
1242  }
1243  d6 = d6 / men;
1244  e6 = 0, men = 0;
1245  for (w = j + (Int_t)ai - bw; w <= j + (Int_t)ai + bw; w++){
1246  if (w >= 0 && w < ssize){
1247  e6 += working_space[ssize + w];
1248  men +=1;
1249  }
1250  }
1251  e6 = e6 / men;
1252  f6 = 0, men = 0;
1253  for (w = j + (Int_t)(2 * ai) - bw; w <= j + (Int_t)(2 * ai) + bw; w++){
1254  if (w >= 0 && w < ssize){
1255  f6 += working_space[ssize + w];
1256  men +=1;
1257  }
1258  }
1259  f6 = f6 / men;
1260  g6 = 0, men = 0;
1261  for (w = j + (Int_t)(3 * ai) - bw; w <= j + (Int_t)(3 * ai) + bw; w++){
1262  if (w >= 0 && w < ssize){
1263  g6 += working_space[ssize + w];
1264  men +=1;
1265  }
1266  }
1267  g6 = g6 / men;
1268  b6 = (b6 - 6 * c6 + 15 * d6 + 15 * e6 - 6 * f6 + g6) / 20;
1269  ai = i / 4;
1270  b8 = 0, men = 0;
1271  for (w = j - (Int_t)(4 * ai) - bw; w <= j - (Int_t)(4 * ai) + bw; w++){
1272  if (w >= 0 && w < ssize){
1273  b8 += working_space[ssize + w];
1274  men +=1;
1275  }
1276  }
1277  b8 = b8 / men;
1278  c8 = 0, men = 0;
1279  for (w = j - (Int_t)(3 * ai) - bw; w <= j - (Int_t)(3 * ai) + bw; w++){
1280  if (w >= 0 && w < ssize){
1281  c8 += working_space[ssize + w];
1282  men +=1;
1283  }
1284  }
1285  c8 = c8 / men;
1286  d8 = 0, men = 0;
1287  for (w = j - (Int_t)(2 * ai) - bw; w <= j - (Int_t)(2 * ai) + bw; w++){
1288  if (w >= 0 && w < ssize){
1289  d8 += working_space[ssize + w];
1290  men +=1;
1291  }
1292  }
1293  d8 = d8 / men;
1294  e8 = 0, men = 0;
1295  for (w = j - (Int_t)ai - bw; w <= j - (Int_t)ai + bw; w++){
1296  if (w >= 0 && w < ssize){
1297  e8 += working_space[ssize + w];
1298  men +=1;
1299  }
1300  }
1301  e8 = e8 / men;
1302  f8 = 0, men = 0;
1303  for (w = j + (Int_t)ai - bw; w <= j + (Int_t)ai + bw; w++){
1304  if (w >= 0 && w < ssize){
1305  f8 += working_space[ssize + w];
1306  men +=1;
1307  }
1308  }
1309  f8 = f8 / men;
1310  g8 = 0, men = 0;
1311  for (w = j + (Int_t)(2 * ai) - bw; w <= j + (Int_t)(2 * ai) + bw; w++){
1312  if (w >= 0 && w < ssize){
1313  g8 += working_space[ssize + w];
1314  men +=1;
1315  }
1316  }
1317  g8 = g8 / men;
1318  h8 = 0, men = 0;
1319  for (w = j + (Int_t)(3 * ai) - bw; w <= j + (Int_t)(3 * ai) + bw; w++){
1320  if (w >= 0 && w < ssize){
1321  h8 += working_space[ssize + w];
1322  men +=1;
1323  }
1324  }
1325  h8 = h8 / men;
1326  i8 = 0, men = 0;
1327  for (w = j + (Int_t)(4 * ai) - bw; w <= j + (Int_t)(4 * ai) + bw; w++){
1328  if (w >= 0 && w < ssize){
1329  i8 += working_space[ssize + w];
1330  men +=1;
1331  }
1332  }
1333  i8 = i8 / men;
1334  b8 = ( -b8 + 8 * c8 - 28 * d8 + 56 * e8 - 56 * f8 - 28 * g8 + 8 * h8 - i8)/70;
1335  if (b < b8)
1336  b = b8;
1337  if (b < b6)
1338  b = b6;
1339  if (b < b4)
1340  b = b4;
1341  if (b < a)
1342  av = b;
1343  working_space[j]=av;
1344  }
1345  }
1346  for (j = i; j < ssize - i; j++)
1347  working_space[ssize + j] = working_space[j];
1348  if (direction == kBackIncreasingWindow)
1349  i += 1;
1350  else if(direction == kBackDecreasingWindow)
1351  i -= 1;
1352  }while((direction == kBackIncreasingWindow && i <= numberIterations) || (direction == kBackDecreasingWindow && i >= 1));
1353  }
1354 
1355  if (compton == kTRUE) {
1356  for (i = 0, b2 = 0; i < ssize; i++){
1357  b1 = b2;
1358  a = working_space[i], b = spectrum[i];
1359  j = i;
1360  if (TMath::Abs(a - b) >= 1) {
1361  b1 = i - 1;
1362  if (b1 < 0)
1363  b1 = 0;
1364  yb1 = working_space[b1];
1365  for (b2 = b1 + 1, c = 0, priz = 0; priz == 0 && b2 < ssize; b2++){
1366  a = working_space[b2], b = spectrum[b2];
1367  c = c + b - yb1;
1368  if (TMath::Abs(a - b) < 1) {
1369  priz = 1;
1370  yb2 = b;
1371  }
1372  }
1373  if (b2 == ssize)
1374  b2 -= 1;
1375  yb2 = working_space[b2];
1376  if (yb1 <= yb2){
1377  for (j = b1, c = 0; j <= b2; j++){
1378  b = spectrum[j];
1379  c = c + b - yb1;
1380  }
1381  if (c > 1){
1382  c = (yb2 - yb1) / c;
1383  for (j = b1, d = 0; j <= b2 && j < ssize; j++){
1384  b = spectrum[j];
1385  d = d + b - yb1;
1386  a = c * d + yb1;
1387  working_space[ssize + j] = a;
1388  }
1389  }
1390  }
1391 
1392  else{
1393  for (j = b2, c = 0; j >= b1; j--){
1394  b = spectrum[j];
1395  c = c + b - yb2;
1396  }
1397  if (c > 1){
1398  c = (yb1 - yb2) / c;
1399  for (j = b2, d = 0;j >= b1 && j >= 0; j--){
1400  b = spectrum[j];
1401  d = d + b - yb2;
1402  a = c * d + yb2;
1403  working_space[ssize + j] = a;
1404  }
1405  }
1406  }
1407  i=b2;
1408  }
1409  }
1410  }
1411 
1412  for (j = 0; j < ssize; j++)
1413  spectrum[j] = working_space[ssize + j];
1414  delete[]working_space;
1415  return 0;
1416 }
1417 
1418 ////////////////////////////////////////////////////////////////////////////////
1419 /// One-dimensional markov spectrum smoothing function
1420 ///
1421 /// This function calculates smoothed spectrum from source spectrum based on
1422 /// Markov chain method. The result is placed in the array pointed by source
1423 /// pointer. On successful completion it returns 0. On error it returns pointer
1424 /// to the string describing error.
1425 ///
1426 /// #### Parameters:
1427 ///
1428 /// - source: pointer to the array of source spectrum
1429 /// - ssize: length of source array
1430 /// - averWindow: width of averaging smoothing window
1431 ///
1432 /// The goal of this function is the suppression of the statistical fluctuations.
1433 /// The algorithm is based on discrete Markov chain, which has very simple
1434 /// invariant distribution:
1435 ///
1436 /// \f[
1437 /// U_2 = \frac{p_{1,2}}{p_{2,1}}U_1, U_3 = \frac{p_{2,3}}{p_{3,2}}U_2U_1, ... , U_n = \frac{p_{n-1,n}}{p_{n,n-1}}U_{n-1}...U_2U_1
1438 /// \f]
1439 /// \f$ U_1\f$ being defined from the normalization condition
1440 /// \f$ \sum_{i=1}^{n} U_i=1\f$. \f$ n \f$ is the length of the smoothed spectrum and
1441 /// \f[
1442 /// p_{i,i\pm 1} = A_i\sum_{k=1}^{m} exp\left[ \frac{y(i\pm k)-y(i)}{y(i\pm k)+y(i)}\right]
1443 /// \f]
1444 ///
1445 /// #### Reference:
1446 ///
1447 /// 1. Z.K. Silagadze, A new algorithm for automatic photopeak searches.
1448 /// NIM A 376 (1996), 451.
1449 ///
1450 /// ### Example 14 - script Smoothing.c
1451 ///
1452 /// \image html TSpectrum_Smoothing1.jpg Fig. 23 Original noisy spectrum
1453 /// \image html TSpectrum_Smoothing2.jpg Fig. 24 Smoothed spectrum m=3
1454 /// \image html TSpectrum_Smoothing3.jpg Fig. 25 Smoothed spectrum
1455 /// \image html TSpectrum_Smoothing4.jpg Fig. 26 Smoothed spectrum m=10
1456 ///
1457 /// #### Script:
1458 ///
1459 /// ~~~ {.cpp}
1460 /// // Example to illustrate smoothing using Markov algorithm (class TSpectrum).
1461 /// // To execute this example, do
1462 /// // root > .x Smoothing.C
1463 ///
1464 /// void Smoothing() {
1465 /// Int_t i;
1466 /// Double_t nbins = 1024;
1467 /// Double_t xmin = 0;
1468 /// Double_t xmax = nbins;
1469 /// Double_t * source = new Double_t[nbins];
1470 /// TH1F *h = new TH1F("h","Smoothed spectrum for m=3",nbins,xmin,xmax);
1471 /// TFile *f = new TFile("spectra/TSpectrum.root");
1472 /// h=(TH1F*) f->Get("smooth1;1");
1473 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
1474 /// TCanvas *Smooth1 = gROOT->GetListOfCanvases()->FindObject("Smooth1");
1475 /// if (!Smooth1) Smooth1 = new TCanvas("Smooth1","Smooth1",10,10,1000,700);
1476 /// TSpectrum *s = new TSpectrum();
1477 /// s->SmoothMarkov(source,1024,3); //3, 7, 10
1478 /// for (i = 0; i < nbins; i++) h->SetBinContent(i + 1,source[i]);
1479 /// h->SetAxisRange(330,880);
1480 /// h->Draw("L");
1481 /// }
1482 /// ~~~
1483 
1484 const char* TSpectrum::SmoothMarkov(Double_t *source, int ssize, int averWindow)
1485 {
1486  int xmin, xmax, i, l;
1487  Double_t a, b, maxch;
1488  Double_t nom, nip, nim, sp, sm, area = 0;
1489  if(averWindow <= 0)
1490  return "Averaging Window must be positive";
1491  Double_t *working_space = new Double_t[ssize];
1492  xmin = 0,xmax = ssize - 1;
1493  for(i = 0, maxch = 0; i < ssize; i++){
1494  working_space[i]=0;
1495  if(maxch < source[i])
1496  maxch = source[i];
1497 
1498  area += source[i];
1499  }
1500  if(maxch == 0) {
1501  delete [] working_space;
1502  return 0 ;
1503  }
1504 
1505  nom = 1;
1506  working_space[xmin] = 1;
1507  for(i = xmin; i < xmax; i++){
1508  nip = source[i] / maxch;
1509  nim = source[i + 1] / maxch;
1510  sp = 0,sm = 0;
1511  for(l = 1; l <= averWindow; l++){
1512  if((i + l) > xmax)
1513  a = source[xmax] / maxch;
1514 
1515  else
1516  a = source[i + l] / maxch;
1517  b = a - nip;
1518  if(a + nip <= 0)
1519  a = 1;
1520 
1521  else
1522  a = TMath::Sqrt(a + nip);
1523  b = b / a;
1524  b = TMath::Exp(b);
1525  sp = sp + b;
1526  if((i - l + 1) < xmin)
1527  a = source[xmin] / maxch;
1528 
1529  else
1530  a = source[i - l + 1] / maxch;
1531  b = a - nim;
1532  if(a + nim <= 0)
1533  a = 1;
1534  else
1535  a = TMath::Sqrt(a + nim);
1536  b = b / a;
1537  b = TMath::Exp(b);
1538  sm = sm + b;
1539  }
1540  a = sp / sm;
1541  a = working_space[i + 1] = working_space[i] * a;
1542  nom = nom + a;
1543  }
1544  for(i = xmin; i <= xmax; i++){
1545  working_space[i] = working_space[i] / nom;
1546  }
1547  for(i = 0; i < ssize; i++)
1548  source[i] = working_space[i] * area;
1549  delete[]working_space;
1550  return 0;
1551 }
1552 
1553 ////////////////////////////////////////////////////////////////////////////////
1554 /// One-dimensional deconvolution function
1555 ///
1556 /// This function calculates deconvolution from source spectrum according to
1557 /// response spectrum using Gold deconvolution algorithm. The result is placed
1558 /// in the vector pointed by source pointer. On successful completion it
1559 /// returns 0. On error it returns pointer to the string describing error. If
1560 /// desired after every numberIterations one can apply boosting operation
1561 /// (exponential function with exponent given by boost coefficient) and repeat
1562 /// it numberRepetitions times.
1563 ///
1564 /// #### Parameters:
1565 ///
1566 /// - source: pointer to the vector of source spectrum
1567 /// - response: pointer to the vector of response spectrum
1568 /// - ssize: length of source and response spectra
1569 /// - numberIterations, for details we refer to the reference given below
1570 /// - numberRepetitions, for repeated boosted deconvolution
1571 /// - boost, boosting coefficient
1572 ///
1573 /// The goal of this function is the improvement of the resolution in spectra,
1574 /// decomposition of multiplets. The mathematical formulation of
1575 /// the convolution system is:
1576 ///
1577 /// \f[
1578 /// y(i) = \sum_{k=0}^{N-1} h(i-k)x(k), i=0,1,2,...,N-1
1579 /// \f]
1580 ///
1581 /// where h(i) is the impulse response function, x, y are input and output
1582 /// vectors, respectively, N is the length of x and h vectors. In matrix form
1583 /// we have:
1584 /**
1585  \f[
1586  \begin{bmatrix}
1587  y(0) \\
1588  y(1) \\
1589  \dots \\
1590  y(2N-2) \\
1591  y(2N-1)
1592  \end{bmatrix}
1593  =
1594  \begin{bmatrix}
1595  h(0) & 0 & 0 & \dots & 0 \\
1596  h(1) & h(0) & 0 & \dots & \dots \\
1597  \dots & h(1) & h(0) & \dots & \dots \\
1598  \dots & \dots & h(1) & \dots & \dots \\
1599  \dots & \dots & \dots & \dots & \dots \\
1600  h(N-1) & \dots & \dots &\dots & 0 \\
1601  0 & h(N-1) & \dots & \dots & h(0) \\
1602  0 & 0 & h(N-1) & \dots & h(1) \\
1603  \dots & \dots & \dots & \dots & \dots \\
1604  0 & 0 & 0 & \dots & h(N-1)
1605  \end{bmatrix}
1606  \begin{bmatrix}
1607  x(0) \\
1608  x(1) \\
1609  \dots \\
1610  x(N-1)
1611  \end{bmatrix}
1612  \f]
1613 */
1614 /// Let us assume that we know the response and the output vector (spectrum) of
1615 /// the above given system. The deconvolution represents solution of the
1616 /// overdetermined system of linear equations, i.e., the calculation of the
1617 /// vector x. From numerical stability point of view the operation of
1618 /// deconvolution is extremely critical (ill-posed problem) as well as time
1619 /// consuming operation. The Gold deconvolution algorithm proves to work very
1620 /// well, other methods (Fourier, VanCittert etc) oscillate. It is suitable to
1621 /// process positive definite data (e.g. histograms).
1622 ///
1623 /// #### Gold deconvolution algorithm:
1624 /**
1625  \f[
1626  y = Hx \\
1627  H^{T}y = H^{T}Hx \\
1628  y^{'} = H^{'}x \\
1629  x_{i}^{(k+1)} = \frac{y_{i}^{'}}{\sum_{m=0}^{N-1}H_{im}^{'}x_{m}{(k)}}x_{i}{(k)}, i=0,1,2,...,N-1 \\
1630  k = 1,2,3,...,L
1631  x^{0} = [1,1, ..., 1]^T
1632  \f]
1633 */
1634 /// Where L is given number of iterations (numberIterations parameter).
1635 ///
1636 /// #### Boosted deconvolution:
1637 ///
1638 /// 1. Set the initial solution:
1639 /// \f$ x^{(0)} = [1,1,...,1]^{T} \f$
1640 /// 2. Set required number of repetitions R and iterations L.
1641 /// 3. Set r = 1.
1642 /// 4. Using Gold deconvolution algorithm for k=1,2,...,L find
1643 /// \f$ x^{(L)} \f$
1644 /// 5. If r = R stop calculation, else
1645 ///
1646 /// 1. Apply boosting operation, i.e., set
1647 /// \f$ x^{(0)}(i) = [x^{(L)}(i)]^{p} \f$
1648 /// i=0,1,...N-1 and p is boosting coefficient >0.
1649 /// 2. r = r + 1
1650 /// 3. continue in 4.
1651 ///
1652 /// #### References:
1653 ///
1654 /// 1. Gold R., ANL-6984, Argonne National Laboratories, Argonne Ill, 1964.
1655 /// 2. Coote G.E., Iterative smoothing and deconvolution of one- and two-dimensional
1656 /// elemental distribution data, NIM B 130 (1997) 118.
1657 /// 3. M. Morhac;, J. Kliman, V. Matouoek, M. Veselsky,
1658 /// I. Turzo: Efficient one- and two-dimensional Gold deconvolution and
1659 /// its application to gamma-ray spectra decomposition. NIM, A401 (1997) 385-408.
1660 /// 4. Morhac; M., Matouoek V., Kliman J., Efficient algorithm of multidimensional
1661 /// deconvolution and its application to nuclear data processing, Digital Signal
1662 /// Processing 13 (2003) 144.
1663 ///
1664 /// ### Example 8 - script Deconvolution.c :
1665 ///
1666 /// response function (usually peak) should be shifted left to the first
1667 /// non-zero channel (bin) (see Fig. 9)
1668 ///
1669 /// \image html TSpectrum_Deconvolution1.jpg Fig. 9 Response spectrum.
1670 /// \image html TSpectrum_Deconvolution2.jpg Fig. 10 Principle how the response matrix is composed inside of the Deconvolution function.
1671 /// \image html TSpectrum_Deconvolution3.jpg Fig. 11 Example of Gold deconvolution. The original source spectrum is drawn with black color, the spectrum after the deconvolution (10000 iterations) with red color.
1672 ///
1673 /// #### Script:
1674 ///
1675 /// ~~~ {.cpp}
1676 /// // Example to illustrate deconvolution function (class TSpectrum).
1677 /// // To execute this example, do
1678 /// // root > .x Deconvolution.C
1679 ///
1680 /// #include <TSpectrum>
1681 ///
1682 /// void Deconvolution() {
1683 /// Int_t i;
1684 /// Double_t nbins = 256;
1685 /// Double_t xmin = 0;
1686 /// Double_t xmax = nbins;
1687 /// Double_t * source = new Double_t[nbins];
1688 /// Double_t * response = new Double_t[nbins];
1689 /// TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
1690 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
1691 /// TFile *f = new TFile("spectra/TSpectrum.root");
1692 /// h=(TH1F*) f->Get("decon1;1");
1693 /// TFile *fr = new TFile("spectra/TSpectrum.root");
1694 /// d=(TH1F*) fr->Get("decon_response;1");
1695 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
1696 /// for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
1697 /// TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
1698 /// if (!Decon1) Decon1 = new TCanvas("Decon1","Decon1",10,10,1000,700);
1699 /// h->Draw("L");
1700 /// TSpectrum *s = new TSpectrum();
1701 /// s->Deconvolution(source,response,256,1000,1,1);
1702 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
1703 /// d->SetLineColor(kRed);
1704 /// d->Draw("SAME L");
1705 /// }
1706 /// ~~~
1707 ///
1708 /// ### Examples of Gold deconvolution method:
1709 ///
1710 /// First let us study the influence of the number of iterations on the
1711 /// deconvolved spectrum (Fig. 12).
1712 ///
1713 /// \image html TSpectrum_Deconvolution_wide1.jpg Fig. 12 Study of Gold deconvolution algorithm.The original source spectrum is drawn with black color, spectrum after 100 iterations with red color, spectrum after 1000 iterations with blue color, spectrum after 10000 iterations with green color and spectrum after 100000 iterations with magenta color.
1714 ///
1715 /// For relatively narrow peaks in the above given example the Gold
1716 /// deconvolution method is able to decompose overlapping peaks practically to
1717 /// delta - functions. In the next example we have chosen a synthetic data
1718 /// (spectrum, 256 channels) consisting of 5 very closely positioned, relatively
1719 /// wide peaks (sigma =5), with added noise (Fig. 13). Thin lines represent
1720 /// pure Gaussians (see Table 1); thick line is a resulting spectrum with
1721 /// additive noise (10% of the amplitude of small peaks).
1722 ///
1723 /// \image html TSpectrum_Deconvolution_wide2.jpg Fig. 13 Testing example of synthetic spectrum composed of 5 Gaussians with added noise.
1724 ///
1725 /// | Peak # | Position | Height | Area |
1726 /// |----------|----------|--------|--------|
1727 /// | 1 | 50 | 500 | 10159 |
1728 /// | 2 | 70 | 3000 | 60957 |
1729 /// | 3 | 80 | 1000 | 20319 |
1730 /// | 4 | 100 | 5000 | 101596 |
1731 /// | 5 | 110 | 500 | 10159 |
1732 ///
1733 /// Table 1 Positions, heights and areas of peaks in the spectrum shown in Fig. 13.
1734 ///
1735 /// In ideal case, we should obtain the result given in Fig. 14. The areas of
1736 /// the Gaussian components of the spectrum are concentrated completely to
1737 /// delta-functions. When solving the overdetermined system of linear equations
1738 /// with data from Fig. 13 in the sense of minimum least squares criterion
1739 /// without any regularisation we obtain the result with large oscillations
1740 /// (Fig. 15). From mathematical point of view, it is the optimal solution in
1741 /// the unconstrained space of independent variables. From physical point of
1742 /// view we are interested only in a meaningful solution. Therefore, we have to
1743 /// employ regularisation techniques (e.g. Gold deconvolution) and/or to
1744 /// confine the space of allowed solutions to subspace of positive solutions.
1745 ///
1746 /// \image html TSpectrum_Deconvolution_wide3.jpg Fig. 14 The same spectrum like in Fig. 13, outlined bars show the contents of present components (peaks).
1747 /// \image html TSpectrum_Deconvolution_wide4.jpg Fig. 15 Least squares solution of the system of linear equations without regularisation.
1748 ///
1749 /// ### Example 9 - script Deconvolution_wide.c
1750 ///
1751 /// When we employ Gold deconvolution algorithm we obtain the result given in
1752 /// Fig. 16. One can observe that the resulting spectrum is smooth. On the
1753 /// other hand the method is not able to decompose completely the peaks in the
1754 /// spectrum.
1755 ///
1756 /// \image html TSpectrum_Deconvolution_wide5.jpg Fig 16 Example of Gold deconvolution for closely positioned wide peaks. The original source spectrum is drawn with black color, the spectrum after the deconvolution (10000 iterations) with red color.
1757 ///
1758 /// #### Script:
1759 ///
1760 /// ~~~ {.cpp}
1761 /// // Example to illustrate deconvolution function (class TSpectrum).
1762 /// // To execute this example, do
1763 /// // root > .x Deconvolution_wide.C
1764 ///
1765 /// #include <TSpectrum>
1766 ///
1767 /// void Deconvolution_wide() {
1768 /// Int_t i;
1769 /// Double_t nbins = 256;
1770 /// Double_t xmin = 0;
1771 /// Double_t xmax = nbins;
1772 /// Double_t * source = new Double_t[nbins];
1773 /// Double_t * response = new Double_t[nbins];
1774 /// TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
1775 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
1776 /// TFile *f = new TFile("spectra/TSpectrum.root");
1777 /// h=(TH1F*) f->Get("decon3;1");
1778 /// TFile *fr = new TFile("spectra/TSpectrum.root");
1779 /// d=(TH1F*) fr->Get("decon_response_wide;1");
1780 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
1781 /// for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
1782 /// TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
1783 /// if (!Decon1) Decon1 = new TCanvas("Decon1",
1784 /// "Deconvolution of closely positioned overlapping peaks using Gold deconvolution method",10,10,1000,700);
1785 /// h->SetMaximum(30000);
1786 /// h->Draw("L");
1787 /// TSpectrum *s = new TSpectrum();
1788 /// s->Deconvolution(source,response,256,10000,1,1);
1789 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
1790 /// d->SetLineColor(kRed);
1791 /// d->Draw("SAME L");
1792 /// }
1793 /// ~~~
1794 ///
1795 /// ### Example 10 - script Deconvolution_wide_boost.c :
1796 ///
1797 /// Further let us employ boosting operation into deconvolution (Fig. 17).
1798 ///
1799 /// \image html TSpectrum_Deconvolution_wide6.jpg Fig. 17 The original source spectrum is drawn with black color, the spectrum after the deconvolution with red color. Number of iterations = 200, number of repetitions = 50 and boosting coefficient = 1.2.
1800 ///
1801 /// | Peak # | Original/Estimated (max) position | Original/Estimated area |
1802 /// |---------|------------------------------------|-------------------------|
1803 /// | 1 | 50/49 | 10159/10419 |
1804 /// | 2 | 70/70 | 60957/58933 |
1805 /// | 3 | 80/79 | 20319/19935 |
1806 /// | 4 | 100/100 | 101596/105413 |
1807 /// | 5 | 110/117 | 10159/6676 |
1808 ///
1809 /// Table 2 Results of the estimation of peaks in spectrum shown in Fig. 17.
1810 ///
1811 /// One can observe that peaks are decomposed practically to delta functions.
1812 /// Number of peaks is correct, positions of big peaks as well as their areas
1813 /// are relatively well estimated. However there is a considerable error in
1814 /// the estimation of the position of small right hand peak.
1815 ///
1816 /// #### Script:
1817 ///
1818 /// ~~~ {.cpp}
1819 /// // Example to illustrate deconvolution function (class TSpectrum).
1820 /// // To execute this example, do
1821 /// // root > .x Deconvolution_wide_boost.C
1822 ///
1823 /// #include <TSpectrum>
1824 ///
1825 /// void Deconvolution_wide_boost() {
1826 /// Int_t i;
1827 /// Double_t nbins = 256;
1828 /// Double_t xmin = 0;
1829 /// Double_t xmax = nbins;
1830 /// Double_t * source = new Double_t[nbins];
1831 /// Double_t * response = new Double_t[nbins];
1832 /// TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
1833 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
1834 /// TFile *f = new TFile("spectra/TSpectrum.root");
1835 /// h=(TH1F*) f->Get("decon3;1");
1836 /// TFile *fr = new TFile("spectra/TSpectrum.root");
1837 /// d=(TH1F*) fr->Get("decon_response_wide;1");
1838 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
1839 /// for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
1840 /// TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
1841 /// if (!Decon1) Decon1 = new TCanvas("Decon1",
1842 /// "Deconvolution of closely positioned overlapping peaks using boosted Gold deconvolution method",10,10,1000,700);
1843 /// h->SetMaximum(110000);
1844 /// h->Draw("L");
1845 /// TSpectrum *s = new TSpectrum();
1846 /// s->Deconvolution(source,response,256,200,50,1.2);
1847 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
1848 /// d->SetLineColor(kRed);
1849 /// d->Draw("SAME L");
1850 /// }
1851 /// ~~~
1852 
1853 const char *TSpectrum::Deconvolution(Double_t *source, const Double_t *response,
1854  int ssize, int numberIterations,
1855  int numberRepetitions, Double_t boost )
1856 {
1857  if (ssize <= 0)
1858  return "Wrong Parameters";
1859 
1860  if (numberRepetitions <= 0)
1861  return "Wrong Parameters";
1862 
1863  // working_space-pointer to the working vector
1864  // (its size must be 4*ssize of source spectrum)
1865  Double_t *working_space = new Double_t[4 * ssize];
1866  int i, j, k, lindex, posit, lh_gold, l, repet;
1867  Double_t lda, ldb, ldc, area, maximum;
1868  area = 0;
1869  lh_gold = -1;
1870  posit = 0;
1871  maximum = 0;
1872 
1873 //read response vector
1874  for (i = 0; i < ssize; i++) {
1875  lda = response[i];
1876  if (lda != 0)
1877  lh_gold = i + 1;
1878  working_space[i] = lda;
1879  area += lda;
1880  if (lda > maximum) {
1881  maximum = lda;
1882  posit = i;
1883  }
1884  }
1885  if (lh_gold == -1) {
1886  delete [] working_space;
1887  return "ZERO RESPONSE VECTOR";
1888  }
1889 
1890 //read source vector
1891  for (i = 0; i < ssize; i++)
1892  working_space[2 * ssize + i] = source[i];
1893 
1894 // create matrix at*a and vector at*y
1895  for (i = 0; i < ssize; i++){
1896  lda = 0;
1897  for (j = 0; j < ssize; j++){
1898  ldb = working_space[j];
1899  k = i + j;
1900  if (k < ssize){
1901  ldc = working_space[k];
1902  lda = lda + ldb * ldc;
1903  }
1904  }
1905  working_space[ssize + i] = lda;
1906  lda = 0;
1907  for (k = 0; k < ssize; k++){
1908  l = k - i;
1909  if (l >= 0){
1910  ldb = working_space[l];
1911  ldc = working_space[2 * ssize + k];
1912  lda = lda + ldb * ldc;
1913  }
1914  }
1915  working_space[3 * ssize + i]=lda;
1916  }
1917 
1918 // move vector at*y
1919  for (i = 0; i < ssize; i++){
1920  working_space[2 * ssize + i] = working_space[3 * ssize + i];
1921  }
1922 
1923 //initialization of resulting vector
1924  for (i = 0; i < ssize; i++)
1925  working_space[i] = 1;
1926 
1927  //**START OF ITERATIONS**
1928  for (repet = 0; repet < numberRepetitions; repet++) {
1929  if (repet != 0) {
1930  for (i = 0; i < ssize; i++)
1931  working_space[i] = TMath::Power(working_space[i], boost);
1932  }
1933  for (lindex = 0; lindex < numberIterations; lindex++) {
1934  for (i = 0; i < ssize; i++) {
1935  if (working_space[2 * ssize + i] > 0.000001
1936  && working_space[i] > 0.000001) {
1937  lda = 0;
1938  for (j = 0; j < lh_gold; j++) {
1939  ldb = working_space[j + ssize];
1940  if (j != 0){
1941  k = i + j;
1942  ldc = 0;
1943  if (k < ssize)
1944  ldc = working_space[k];
1945  k = i - j;
1946  if (k >= 0)
1947  ldc += working_space[k];
1948  }
1949 
1950  else
1951  ldc = working_space[i];
1952  lda = lda + ldb * ldc;
1953  }
1954  ldb = working_space[2 * ssize + i];
1955  if (lda != 0)
1956  lda = ldb / lda;
1957 
1958  else
1959  lda = 0;
1960  ldb = working_space[i];
1961  lda = lda * ldb;
1962  working_space[3 * ssize + i] = lda;
1963  }
1964  }
1965  for (i = 0; i < ssize; i++)
1966  working_space[i] = working_space[3 * ssize + i];
1967  }
1968  }
1969 
1970 //shift resulting spectrum
1971  for (i = 0; i < ssize; i++) {
1972  lda = working_space[i];
1973  j = i + posit;
1974  j = j % ssize;
1975  working_space[ssize + j] = lda;
1976  }
1977 
1978 //write back resulting spectrum
1979  for (i = 0; i < ssize; i++)
1980  source[i] = area * working_space[ssize + i];
1981  delete[]working_space;
1982  return 0;
1983 }
1984 
1985 
1986 ////////////////////////////////////////////////////////////////////////////////
1987 /// One-dimensional deconvolution function.
1988 ///
1989 /// This function calculates deconvolution from source spectrum according to
1990 /// response spectrum using Richardson-Lucy deconvolution algorithm. The result
1991 /// is placed in the vector pointed by source pointer. On successful completion
1992 /// it returns 0. On error it returns pointer to the string describing error.
1993 /// If desired after every numberIterations one can apply boosting operation
1994 /// (exponential function with exponent given by boost coefficient) and repeat
1995 /// it numberRepetitions times (see Gold deconvolution).
1996 ///
1997 /// #### Parameters:
1998 ///
1999 /// - source: pointer to the vector of source spectrum
2000 /// - response: pointer to the vector of response spectrum
2001 /// - ssize: length of source and response spectra
2002 /// - numberIterations, for details we refer to the reference given above
2003 /// - numberRepetitions, for repeated boosted deconvolution
2004 /// - boost, boosting coefficient
2005 ///
2006 /// ### Richardson-Lucy deconvolution algorithm:
2007 ///
2008 /// For discrete systems it has the form:
2009 /**
2010  \f[
2011  x^{(n)}(i) = x^{(n-1)}(i) \sum_{j=0}^{N-1}h(i,j)\frac{y(j)}{\sum_{k=0}^{M-1}h(j,k)x^{(n-1)}(k)} \\
2012  i \in \left<0,M-1\right>
2013  \f]
2014 */
2015 /// for positive input data and response matrix this iterative method forces
2016 /// the deconvoluted spectra to be non-negative. The Richardson-Lucy
2017 /// iteration converges to the maximum likelihood solution for Poisson statistics
2018 /// in the data.
2019 ///
2020 /// #### References:
2021 ///
2022 /// 1. Abreu M.C. et al., A four-dimensional deconvolution method to correct NA38
2023 /// experimental data, NIM A 405 (1998) 139.
2024 /// 2. Lucy L.B., A.J. 79 (1974) 745.
2025 /// 3. Richardson W.H., J. Opt. Soc. Am. 62 (1972) 55.
2026 ///
2027 /// ### Examples of Richardson-Lucy deconvolution method:
2028 ///
2029 /// ### Example 11 - script DeconvolutionRL_wide.c :
2030 ///
2031 /// When we employ Richardson-Lucy deconvolution algorithm to our data from
2032 /// Fig. 13 we obtain the result given in Fig. 18. One can observe improvements
2033 /// as compared to the result achieved by Gold deconvolution. Nevertheless it is
2034 /// unable to decompose the multiplet.
2035 ///
2036 /// \image html TSpectrum_DeconvolutionRL_wide1.jpg Fig. 18 Example of Richardson-Lucy deconvolution for closely positioned wide peaks. The original source spectrum is drawn with black color, the spectrum after the deconvolution (10000 iterations) with red color.
2037 ///
2038 /// #### Script:
2039 ///
2040 /// ~~~ {.cpp}
2041 /// // Example to illustrate deconvolution function (class TSpectrum).
2042 /// // To execute this example, do
2043 /// // root > .x DeconvolutionRL_wide.C
2044 ///
2045 /// #include <TSpectrum>
2046 ///
2047 /// void DeconvolutionRL_wide() {
2048 /// Int_t i;
2049 /// Double_t nbins = 256;
2050 /// Double_t xmin = 0;
2051 /// Double_t xmax = nbins;
2052 /// Double_t * source = new Double_t[nbins];
2053 /// Double_t * response = new Double_t[nbins];
2054 /// TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
2055 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
2056 /// TFile *f = new TFile("spectra/TSpectrum.root");
2057 /// h=(TH1F*) f->Get("decon3;1");
2058 /// TFile *fr = new TFile("spectra/TSpectrum.root");
2059 /// d=(TH1F*) fr->Get("decon_response_wide;1");
2060 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2061 /// for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
2062 /// TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
2063 /// if (!Decon1) Decon1 = new TCanvas("Decon1",
2064 /// "Deconvolution of closely positioned overlapping peaks using Richardson-Lucy deconvolution method",
2065 /// 10,10,1000,700);
2066 /// h->SetMaximum(30000);
2067 /// h->Draw("L");
2068 /// TSpectrum *s = new TSpectrum();
2069 /// s->DeconvolutionRL(source,response,256,10000,1,1);
2070 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
2071 /// d->SetLineColor(kRed);
2072 /// d->Draw("SAME L");
2073 /// }
2074 /// ~~~
2075 ///
2076 /// ### Example 12 - script DeconvolutionRL_wide_boost.c :
2077 ///
2078 /// Further let us employ boosting operation into deconvolution (Fig. 19).
2079 ///
2080 /// \image html TSpectrum_DeconvolutionRL_wide2.jpg Fig. 19 The original source spectrum is drawn with black color, the spectrum after the deconvolution with red color. Number of iterations = 200, number of repetitions = 50 and boosting coefficient = 1.2.
2081 ///
2082 /// | Peak # | Original/Estimated (max) position | Original/Estimated area |
2083 /// |--------|-----------------------------------|--------------------------|
2084 /// | 1 | 50/51 | 10159/11426 |
2085 /// | 2 | 70/71 | 60957/65003 |
2086 /// | 3 | 80/81 | 20319/12813 |
2087 /// | 4 | 100/100 | 101596/101851 |
2088 /// | 5 | 110/111 | 10159/8920 |
2089 ///
2090 /// Table 3 Results of the estimation of peaks in spectrum shown in Fig. 19.
2091 ///
2092 /// One can observe improvements in the estimation of peak positions as compared
2093 /// to the results achieved by Gold deconvolution.
2094 ///
2095 /// #### Script:
2096 ///
2097 /// ~~~ {.cpp}
2098 /// // Example to illustrate deconvolution function (class TSpectrum).
2099 /// // To execute this example, do
2100 /// // root > .x DeconvolutionRL_wide_boost.C
2101 ///
2102 /// #include <TSpectrum>
2103 ///
2104 /// void DeconvolutionRL_wide_boost() {
2105 /// Int_t i;
2106 /// Double_t nbins = 256;
2107 /// Double_t xmin = 0;
2108 /// Double_t xmax = nbins;
2109 /// Double_t * source = new Double_t[nbins];
2110 /// Double_t * response = new Double_t[nbins];
2111 /// TH1F *h = new TH1F("h","Deconvolution",nbins,xmin,xmax);
2112 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
2113 /// TFile *f = new TFile("spectra/TSpectrum.root");
2114 /// h=(TH1F*) f->Get("decon3;1");
2115 /// TFile *fr = new TFile("spectra/TSpectrum.root");
2116 /// d=(TH1F*) fr->Get("decon_response_wide;1");
2117 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2118 /// for (i = 0; i < nbins; i++) response[i]=d->GetBinContent(i + 1);
2119 /// TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
2120 /// if (!Decon1) Decon1 = new TCanvas("Decon1",
2121 /// "Deconvolution of closely positioned overlapping peaks using boosted Richardson-Lucy deconvolution method",
2122 /// 10,10,1000,700);
2123 /// h->SetMaximum(110000);
2124 /// h->Draw("L");
2125 /// TSpectrum *s = new TSpectrum();
2126 /// s->DeconvolutionRL(source,response,256,200,50,1.2);
2127 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,source[i]);
2128 /// d->SetLineColor(kRed);
2129 /// d->Draw("SAME L");
2130 /// }
2131 /// ~~~
2132 
2133 const char *TSpectrum::DeconvolutionRL(Double_t *source, const Double_t *response,
2134  int ssize, int numberIterations,
2135  int numberRepetitions, Double_t boost )
2136 {
2137  if (ssize <= 0)
2138  return "Wrong Parameters";
2139 
2140  if (numberRepetitions <= 0)
2141  return "Wrong Parameters";
2142 
2143  // working_space-pointer to the working vector
2144  // (its size must be 4*ssize of source spectrum)
2145  Double_t *working_space = new Double_t[4 * ssize];
2146  int i, j, k, lindex, posit, lh_gold, repet, kmin, kmax;
2147  Double_t lda, ldb, ldc, maximum;
2148  lh_gold = -1;
2149  posit = 0;
2150  maximum = 0;
2151 
2152 //read response vector
2153  for (i = 0; i < ssize; i++) {
2154  lda = response[i];
2155  if (lda != 0)
2156  lh_gold = i + 1;
2157  working_space[ssize + i] = lda;
2158  if (lda > maximum) {
2159  maximum = lda;
2160  posit = i;
2161  }
2162  }
2163  if (lh_gold == -1) {
2164  delete [] working_space;
2165  return "ZERO RESPONSE VECTOR";
2166  }
2167 
2168 //read source vector
2169  for (i = 0; i < ssize; i++)
2170  working_space[2 * ssize + i] = source[i];
2171 
2172 //initialization of resulting vector
2173  for (i = 0; i < ssize; i++){
2174  if (i <= ssize - lh_gold)
2175  working_space[i] = 1;
2176 
2177  else
2178  working_space[i] = 0;
2179 
2180  }
2181  //**START OF ITERATIONS**
2182  for (repet = 0; repet < numberRepetitions; repet++) {
2183  if (repet != 0) {
2184  for (i = 0; i < ssize; i++)
2185  working_space[i] = TMath::Power(working_space[i], boost);
2186  }
2187  for (lindex = 0; lindex < numberIterations; lindex++) {
2188  for (i = 0; i <= ssize - lh_gold; i++){
2189  lda = 0;
2190  if (working_space[i] > 0){//x[i]
2191  for (j = i; j < i + lh_gold; j++){
2192  ldb = working_space[2 * ssize + j];//y[j]
2193  if (j < ssize){
2194  if (ldb > 0){//y[j]
2195  kmax = j;
2196  if (kmax > lh_gold - 1)
2197  kmax = lh_gold - 1;
2198  kmin = j + lh_gold - ssize;
2199  if (kmin < 0)
2200  kmin = 0;
2201  ldc = 0;
2202  for (k = kmax; k >= kmin; k--){
2203  ldc += working_space[ssize + k] * working_space[j - k];//h[k]*x[j-k]
2204  }
2205  if (ldc > 0)
2206  ldb = ldb / ldc;
2207 
2208  else
2209  ldb = 0;
2210  }
2211  ldb = ldb * working_space[ssize + j - i];//y[j]*h[j-i]/suma(h[j][k]x[k])
2212  }
2213  lda += ldb;
2214  }
2215  lda = lda * working_space[i];
2216  }
2217  working_space[3 * ssize + i] = lda;
2218  }
2219  for (i = 0; i < ssize; i++)
2220  working_space[i] = working_space[3 * ssize + i];
2221  }
2222  }
2223 
2224 //shift resulting spectrum
2225  for (i = 0; i < ssize; i++) {
2226  lda = working_space[i];
2227  j = i + posit;
2228  j = j % ssize;
2229  working_space[ssize + j] = lda;
2230  }
2231 
2232 //write back resulting spectrum
2233  for (i = 0; i < ssize; i++)
2234  source[i] = working_space[ssize + i];
2235  delete[]working_space;
2236  return 0;
2237 }
2238 
2239 
2240 ////////////////////////////////////////////////////////////////////////////////
2241 /// One-dimensional unfolding function
2242 ///
2243 /// This function unfolds source spectrum according to response matrix columns.
2244 /// The result is placed in the vector pointed by source pointer.
2245 /// The coefficients of the resulting vector represent contents of the columns
2246 /// (weights) in the input vector. On successful completion it returns 0. On
2247 /// error it returns pointer to the string describing error. If desired after
2248 /// every numberIterations one can apply boosting operation (exponential
2249 /// function with exponent given by boost coefficient) and repeat it
2250 /// numberRepetitions times. For details we refer to [1].
2251 ///
2252 /// #### Parameters:
2253 ///
2254 /// - source: pointer to the vector of source spectrum
2255 /// - respMatrix: pointer to the matrix of response spectra
2256 /// - ssizex: length of source spectrum and # of columns of the response
2257 /// matrix. ssizex must be >= ssizey.
2258 /// - ssizey: length of destination spectrum and # of rows of the response
2259 /// matrix.
2260 /// - numberIterations: number of iterations
2261 /// - numberRepetitions: number of repetitions for boosted deconvolution.
2262 /// It must be greater or equal to one.
2263 /// - boost: boosting coefficient, applies only if numberRepetitions is
2264 /// greater than one.
2265 ///
2266 /// ### Unfolding:
2267 ///
2268 /// The goal is the decomposition of spectrum to a given set of component spectra.
2269 ///
2270 /// The mathematical formulation of the discrete linear system is:
2271 ///
2272 /// \f[
2273 /// y(i) = \sum_{k=0}^{N_y-1} h(i,k)x(k), i = 0,1,2,...,N_x-1
2274 /// \f]
2275 /**
2276  \f[
2277  \begin{bmatrix}
2278  y(0) \\
2279  y(1) \\
2280  \dots \\
2281  y(N_x-1)
2282  \end{bmatrix}
2283  =
2284  \begin{bmatrix}
2285  h(0,0) & h(0,1) & \dots & h(0,N_y-1) \\
2286  h(1,0) & h(1,1) & \dots & h(1,N_y-1) \\
2287  \dots \\
2288  h(N_x-1,0) & h(N_x-1,1) & \dots & h(N_x-1,N_y-1)
2289  \end{bmatrix}
2290  \begin{bmatrix}
2291  x(0) \\
2292  x(1) \\
2293  \dots \\
2294  x(N_x-1)
2295  \end{bmatrix}
2296  \f]
2297 */
2298 /// #### References:
2299 ///
2300 /// 1. Jandel M., Morhac; M., Kliman J., Krupa L., Matouoek
2301 /// V., Hamilton J. H., Ramaya A. V.:
2302 /// Decomposition of continuum gamma-ray spectra using synthesised response matrix.
2303 /// NIM A 516 (2004), 172-183.
2304 ///
2305 /// ### Example of unfolding:
2306 ///
2307 /// ### Example 13 - script Unfolding.c:
2308 ///
2309 /// \image html TSpectrum_Unfolding3.gif Fig. 20 Response matrix composed of neutron spectra of pure chemical elements.
2310 /// \image html TSpectrum_Unfolding2.jpg Fig. 21 Source neutron spectrum to be decomposed
2311 /// \image html TSpectrum_Unfolding3.jpg Fig. 22 Spectrum after decomposition, contains 10 coefficients, which correspond to contents of chemical components (dominant 8-th and 10-th components, i.e. O, Si)
2312 ///
2313 /// #### Script:
2314 ///
2315 /// ~~~ {.cpp}
2316 /// // Example to illustrate unfolding function (class TSpectrum).
2317 /// // To execute this example, do
2318 /// // root > .x Unfolding.C
2319 ///
2320 /// #include <TSpectrum>
2321 ///
2322 /// void Unfolding() {
2323 /// Int_t i, j;
2324 /// Int_t nbinsx = 2048;
2325 /// Int_t nbinsy = 10;
2326 /// Double_t xmin = 0;
2327 /// Double_t xmax = nbinsx;
2328 /// Double_t ymin = 0;
2329 /// Double_t ymax = nbinsy;
2330 /// Double_t * source = new Double_t[nbinsx];
2331 /// Double_t ** response = new Double_t *[nbinsy];
2332 /// for (i=0;i<nbinsy;i++) response[i]=new Double_t[nbinsx];
2333 /// TH1F *h = new TH1F("h","",nbinsx,xmin,xmax);
2334 /// TH1F *d = new TH1F("d","Decomposition - unfolding",nbinsx,xmin,xmax);
2335 /// TH2F *decon_unf_resp = new TH2F("decon_unf_resp","Root File",nbinsy,ymin,ymax,nbinsx,xmin,xmax);
2336 /// TFile *f = new TFile("spectra/TSpectrum.root");
2337 /// h=(TH1F*) f->Get("decon_unf_in;1");
2338 /// TFile *fr = new TFile("spectra/TSpectrum.root");
2339 /// decon_unf_resp = (TH2F*) fr->Get("decon_unf_resp;1");
2340 /// for (i = 0; i < nbinsx; i++) source[i] = h->GetBinContent(i + 1);
2341 /// for (i = 0; i < nbinsy; i++){
2342 /// for (j = 0; j< nbinsx; j++){
2343 /// response[i][j] = decon_unf_resp->GetBinContent(i + 1, j + 1);
2344 /// }
2345 /// }
2346 /// TCanvas *Decon1 = gROOT->GetListOfCanvases()->FindObject("Decon1");
2347 /// if (!Decon1) Decon1 = new TCanvas("Decon1","Decon1",10,10,1000,700);
2348 /// h->Draw("L");
2349 /// TSpectrum *s = new TSpectrum();
2350 /// s->Unfolding(source,response,nbinsx,nbinsy,1000,1,1);
2351 /// for (i = 0; i < nbinsy; i++) d->SetBinContent(i + 1,source[i]);
2352 /// d->SetLineColor(kRed);
2353 /// d->SetAxisRange(0,nbinsy);
2354 /// d->Draw("");
2355 /// }
2356 /// ~~~
2357 
2358 const char *TSpectrum::Unfolding(Double_t *source,
2359  const Double_t **respMatrix,
2360  int ssizex, int ssizey,
2361  int numberIterations,
2362  int numberRepetitions, Double_t boost)
2363 {
2364  int i, j, k, lindex, lhx = 0, repet;
2365  Double_t lda, ldb, ldc, area;
2366  if (ssizex <= 0 || ssizey <= 0)
2367  return "Wrong Parameters";
2368  if (ssizex < ssizey)
2369  return "Sizex must be greater than sizey)";
2370  if (numberIterations <= 0)
2371  return "Number of iterations must be positive";
2372  Double_t *working_space =
2373  new Double_t[ssizex * ssizey + 2 * ssizey * ssizey + 4 * ssizex];
2374 
2375 /*read response matrix*/
2376  for (j = 0; j < ssizey && lhx != -1; j++) {
2377  area = 0;
2378  lhx = -1;
2379  for (i = 0; i < ssizex; i++) {
2380  lda = respMatrix[j][i];
2381  if (lda != 0) {
2382  lhx = i + 1;
2383  }
2384  working_space[j * ssizex + i] = lda;
2385  area = area + lda;
2386  }
2387  if (lhx != -1) {
2388  for (i = 0; i < ssizex; i++)
2389  working_space[j * ssizex + i] /= area;
2390  }
2391  }
2392  if (lhx == -1) {
2393  delete [] working_space;
2394  return ("ZERO COLUMN IN RESPONSE MATRIX");
2395  }
2396 
2397 /*read source vector*/
2398  for (i = 0; i < ssizex; i++)
2399  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
2400  source[i];
2401 
2402 /*create matrix at*a + at*y */
2403  for (i = 0; i < ssizey; i++) {
2404  for (j = 0; j < ssizey; j++) {
2405  lda = 0;
2406  for (k = 0; k < ssizex; k++) {
2407  ldb = working_space[ssizex * i + k];
2408  ldc = working_space[ssizex * j + k];
2409  lda = lda + ldb * ldc;
2410  }
2411  working_space[ssizex * ssizey + ssizey * i + j] = lda;
2412  }
2413  lda = 0;
2414  for (k = 0; k < ssizex; k++) {
2415  ldb = working_space[ssizex * i + k];
2416  ldc =
2417  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
2418  k];
2419  lda = lda + ldb * ldc;
2420  }
2421  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
2422  lda;
2423  }
2424 
2425 /*move vector at*y*/
2426  for (i = 0; i < ssizey; i++)
2427  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
2428  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
2429 
2430 /*create matrix at*a*at*a + vector at*a*at*y */
2431  for (i = 0; i < ssizey; i++) {
2432  for (j = 0; j < ssizey; j++) {
2433  lda = 0;
2434  for (k = 0; k < ssizey; k++) {
2435  ldb = working_space[ssizex * ssizey + ssizey * i + k];
2436  ldc = working_space[ssizex * ssizey + ssizey * j + k];
2437  lda = lda + ldb * ldc;
2438  }
2439  working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j] =
2440  lda;
2441  }
2442  lda = 0;
2443  for (k = 0; k < ssizey; k++) {
2444  ldb = working_space[ssizex * ssizey + ssizey * i + k];
2445  ldc =
2446  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex +
2447  k];
2448  lda = lda + ldb * ldc;
2449  }
2450  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] =
2451  lda;
2452  }
2453 
2454 /*move at*a*at*y*/
2455  for (i = 0; i < ssizey; i++)
2456  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i] =
2457  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
2458 
2459 /*initialization in resulting vector */
2460  for (i = 0; i < ssizey; i++)
2461  working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = 1;
2462 
2463  /***START OF ITERATIONS***/
2464  for (repet = 0; repet < numberRepetitions; repet++) {
2465  if (repet != 0) {
2466  for (i = 0; i < ssizey; i++)
2467  working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] = TMath::Power(working_space[ssizex * ssizey + 2 * ssizey * ssizey + i], boost);
2468  }
2469  for (lindex = 0; lindex < numberIterations; lindex++) {
2470  for (i = 0; i < ssizey; i++) {
2471  lda = 0;
2472  for (j = 0; j < ssizey; j++) {
2473  ldb =
2474  working_space[ssizex * ssizey + ssizey * ssizey + ssizey * i + j];
2475  ldc = working_space[ssizex * ssizey + 2 * ssizey * ssizey + j];
2476  lda = lda + ldb * ldc;
2477  }
2478  ldb =
2479  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 2 * ssizex + i];
2480  if (lda != 0) {
2481  lda = ldb / lda;
2482  }
2483 
2484  else
2485  lda = 0;
2486  ldb = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
2487  lda = lda * ldb;
2488  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i] = lda;
2489  }
2490  for (i = 0; i < ssizey; i++)
2491  working_space[ssizex * ssizey + 2 * ssizey * ssizey + i] =
2492  working_space[ssizex * ssizey + 2 * ssizey * ssizey + 3 * ssizex + i];
2493  }
2494  }
2495 
2496 /*write back resulting spectrum*/
2497  for (i = 0; i < ssizex; i++) {
2498  if (i < ssizey)
2499  source[i] = working_space[ssizex * ssizey + 2 * ssizey * ssizey + i];
2500 
2501  else
2502  source[i] = 0;
2503  }
2504  delete[]working_space;
2505  return 0;
2506 }
2507 
2508 
2509 ////////////////////////////////////////////////////////////////////////////////
2510 /// One-dimensional high-resolution peak search function
2511 ///
2512 /// This function searches for peaks in source spectrum. It is based on
2513 /// deconvolution method. First the background is removed (if desired), then
2514 /// Markov smoothed spectrum is calculated (if desired), then the response
2515 /// function is generated according to given sigma and deconvolution is
2516 /// carried out. The order of peaks is arranged according to their heights in
2517 /// the spectrum after background elimination. The highest peak is the first in
2518 /// the list. On success it returns number of found peaks.
2519 ///
2520 /// #### Parameters:
2521 ///
2522 /// - source: pointer to the vector of source spectrum.
2523 /// - destVector: pointer to the vector of resulting deconvolved spectrum.
2524 /// - ssize: length of source spectrum.
2525 /// - sigma: sigma of searched peaks, for details we refer to manual.
2526 /// - threshold: threshold value in % for selected peaks, peaks with
2527 /// amplitude less than threshold*highest_peak/100
2528 /// are ignored, see manual.
2529 /// - backgroundRemove: logical variable, set if the removal of
2530 /// background before deconvolution is desired.
2531 /// - deconIterations-number of iterations in deconvolution operation.
2532 /// - markov: logical variable, if it is true, first the source spectrum
2533 /// is replaced by new spectrum calculated using Markov
2534 /// chains method.
2535 /// - averWindow: averaging window of searched peaks, for details
2536 /// we refer to manual (applies only for Markov method).
2537 ///
2538 /// ### Peaks searching:
2539 ///
2540 /// The goal of this function is to identify automatically the peaks in spectrum
2541 /// with the presence of the continuous background and statistical
2542 /// fluctuations - noise.
2543 ///
2544 /// The common problems connected with correct peak identification are:
2545 ///
2546 /// - non-sensitivity to noise, i.e., only statistically
2547 /// relevant peaks should be identified.
2548 /// - non-sensitivity of the algorithm to continuous
2549 /// background.
2550 /// - ability to identify peaks close to the edges of the
2551 /// spectrum region. Usually peak finders fail to detect them.
2552 /// - resolution, decomposition of Double_tts and multiplets.
2553 /// The algorithm should be able to recognise close positioned peaks.
2554 /// - ability to identify peaks with different sigma.
2555 ///
2556 /// \image html TSpectrum_Searching1.jpg Fig. 27 An example of one-dimensional synthetic spectrum with found peaks denoted by markers.
2557 ///
2558 /// #### References:
2559 ///
2560 /// 1. M.A. Mariscotti: A method for identification of peaks in the presence of
2561 /// background and its application to spectrum analysis. NIM 50 (1967),
2562 /// 309-320.
2563 /// 2. M. Morhac;, J. Kliman, V. Matouoek, M. Veselsky,
2564 /// I. Turzo.:Identification of peaks in
2565 /// multidimensional coincidence gamma-ray spectra. NIM, A443 (2000) 108-125.
2566 /// 3. Z.K. Silagadze, A new algorithm for automatic photopeak searches. NIM
2567 /// A 376 (1996), 451.
2568 ///
2569 /// Examples of peak searching method:
2570 ///
2571 /// The SearchHighRes function provides users with the possibility to vary the
2572 /// input parameters and with the access to the output deconvolved data in the
2573 /// destination spectrum. Based on the output data one can tune the parameters.
2574 ///
2575 /// ### Example 15 - script SearchHR1.c:
2576 ///
2577 /// \image html TSpectrum_Searching1.jpg Fig. 28 One-dimensional spectrum with found peaks denoted by markers, 3 iterations steps in the deconvolution.
2578 /// \image html TSpectrum_Searching2.jpg Fig. 29 One-dimensional spectrum with found peaks denoted by markers, 8 iterations steps in the deconvolution.
2579 ///
2580 /// #### Script:
2581 ///
2582 /// ~~~ {.cpp}
2583 /// // Example to illustrate high resolution peak searching function (class TSpectrum).
2584 /// // To execute this example, do
2585 /// // root > .x SearchHR1.C
2586 ///
2587 /// #include <TSpectrum>
2588 ///
2589 /// void SearchHR1() {
2590 /// Double_t fPositionX[100];
2591 /// Double_t fPositionY[100];
2592 /// Int_t fNPeaks = 0;
2593 /// Int_t i,nfound,bin;
2594 /// Double_t nbins = 1024,a;
2595 /// Double_t xmin = 0;
2596 /// Double_t xmax = nbins;
2597 /// Double_t * source = new Double_t[nbins];
2598 /// Double_t * dest = new Double_t[nbins];
2599 /// TH1F *h = new TH1F("h","High resolution peak searching, number of iterations = 3",nbins,xmin,xmax);
2600 /// TH1F *d = new TH1F("d","",nbins,xmin,xmax);
2601 /// TFile *f = new TFile("spectra/TSpectrum.root");
2602 /// h=(TH1F*) f->Get("search2;1");
2603 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2604 /// TCanvas *Search = gROOT->GetListOfCanvases()->FindObject("Search");
2605 /// if (!Search) Search = new TCanvas("Search","Search",10,10,1000,700);
2606 /// h->SetMaximum(4000);
2607 /// h->Draw("L");
2608 /// TSpectrum *s = new TSpectrum();
2609 /// nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
2610 /// Double_t *xpeaks = s->GetPositionX();
2611 /// for (i = 0; i < nfound; i++) {
2612 /// a=xpeaks[i];
2613 /// bin = 1 + Int_t(a + 0.5);
2614 /// fPositionX[i] = h->GetBinCenter(bin);
2615 /// fPositionY[i] = h->GetBinContent(bin);
2616 /// }
2617 /// TPolyMarker * pm = (TPolyMarker*)h->GetListOfFunctions()->FindObject("TPolyMarker");
2618 /// if (pm) {
2619 /// h->GetListOfFunctions()->Remove(pm);
2620 /// delete pm;
2621 /// }
2622 /// pm = new TPolyMarker(nfound, fPositionX, fPositionY);
2623 /// h->GetListOfFunctions()->Add(pm);
2624 /// pm->SetMarkerStyle(23);
2625 /// pm->SetMarkerColor(kRed);
2626 /// pm->SetMarkerSize(1.3);
2627 /// for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,dest[i]);
2628 /// d->SetLineColor(kRed);
2629 /// d->Draw("SAME");
2630 /// printf("Found %d candidate peaks\n",nfound);
2631 /// for(i=0;i<nfound;i++)
2632 /// printf("posx= %d, posy= %d\n",fPositionX[i], fPositionY[i]);
2633 /// }
2634 /// ~~~
2635 ///
2636 /// ### Example 16 - script SearchHR3.c:
2637 ///
2638 /// | Peak # | Position | Sigma |
2639 /// |--------|----------|-------|
2640 /// | 1 | 118 | 26 |
2641 /// | 2 | 162 | 41 |
2642 /// | 3 | 310 | 4 |
2643 /// | 4 | 330 | 8 |
2644 /// | 5 | 482 | 22 |
2645 /// | 6 | 491 | 26 |
2646 /// | 7 | 740 | 21 |
2647 /// | 8 | 852 | 15 |
2648 /// | 9 | 954 | 12 |
2649 /// | 10 | 989 | 13 |
2650 ///
2651 /// Table 4 Positions and sigma of peaks in the following examples.
2652 ///
2653 /// \image html TSpectrum_Searching3.jpg Fig. 30 Influence of number of iterations (3-red, 10-blue, 100- green, 1000-magenta), sigma=8, smoothing width=3.
2654 /// \image html TSpectrum_Searching4.jpg Fig. 31 Influence of sigma (3-red, 8-blue, 20- green, 43-magenta), num. iter.=10, sm. width=3.
2655 /// \image html TSpectrum_Searching5.jpg Fig. 32 Influence smoothing width (0-red, 3-blue, 7- green, 20-magenta), num. iter.=10, sigma=8.
2656 ///
2657 /// #### Script:
2658 ///
2659 /// ~~~ {.cpp}
2660 /// // Example to illustrate the influence of number of iterations in deconvolution in high resolution peak searching function (class TSpectrum).
2661 /// // To execute this example, do
2662 /// // root > .x SearchHR3.C
2663 ///
2664 /// #include <TSpectrum>
2665 ///
2666 /// void SearchHR3() {
2667 /// Double_t fPositionX[100];
2668 /// Double_t fPositionY[100];
2669 /// Int_t fNPeaks = 0;
2670 /// Int_t i,nfound,bin;
2671 /// Double_t nbins = 1024,a;
2672 /// Double_t xmin = 0;
2673 /// Double_t xmax = nbins;
2674 /// Double_t * source = new Double_t[nbins];
2675 /// Double_t * dest = new Double_t[nbins];
2676 /// TH1F *h = new TH1F("h","Influence of # of iterations in deconvolution in peak searching",nbins,xmin,xmax);
2677 /// TH1F *d1 = new TH1F("d1","",nbins,xmin,xmax);
2678 /// TH1F *d2 = new TH1F("d2","",nbins,xmin,xmax);
2679 /// TH1F *d3 = new TH1F("d3","",nbins,xmin,xmax);
2680 /// TH1F *d4 = new TH1F("d4","",nbins,xmin,xmax);
2681 /// TFile *f = new TFile("spectra/TSpectrum.root");
2682 /// h=(TH1F*) f->Get("search3;1");
2683 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2684 /// TCanvas *Search = gROOT->GetListOfCanvases()->FindObject("Search");
2685 /// if (!Search) Search = new TCanvas("Search","Search",10,10,1000,700);
2686 /// h->SetMaximum(1300);
2687 /// h->Draw("L");
2688 /// TSpectrum *s = new TSpectrum();
2689 /// nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
2690 /// Double_t *xpeaks = s->GetPositionX();
2691 /// for (i = 0; i < nfound; i++) {
2692 /// a=xpeaks[i];
2693 /// bin = 1 + Int_t(a + 0.5);
2694 /// fPositionX[i] = h->GetBinCenter(bin);
2695 /// fPositionY[i] = h->GetBinContent(bin);
2696 /// }
2697 /// TPolyMarker * pm = (TPolyMarker*)h->GetListOfFunctions()->FindObject("TPolyMarker");
2698 /// if (pm) {
2699 /// h->GetListOfFunctions()->Remove(pm);
2700 /// delete pm;
2701 /// }
2702 /// pm = new TPolyMarker(nfound, fPositionX, fPositionY);
2703 /// h->GetListOfFunctions()->Add(pm);
2704 /// pm->SetMarkerStyle(23);
2705 /// pm->SetMarkerColor(kRed);
2706 /// pm->SetMarkerSize(1.3);
2707 /// for (i = 0; i < nbins; i++) d1->SetBinContent(i + 1,dest[i]);
2708 /// h->Draw("");
2709 /// d1->SetLineColor(kRed);
2710 /// d1->Draw("SAME");
2711 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2712 /// s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 10, kTRUE, 3);
2713 /// for (i = 0; i < nbins; i++) d2->SetBinContent(i + 1,dest[i]);
2714 /// d2->SetLineColor(kBlue);
2715 /// d2->Draw("SAME");
2716 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2717 /// s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 100, kTRUE, 3);
2718 /// for (i = 0; i < nbins; i++) d3->SetBinContent(i + 1,dest[i]);
2719 /// d3->SetLineColor(kGreen);
2720 /// d3->Draw("SAME");
2721 /// for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);
2722 /// s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 1000, kTRUE, 3);
2723 /// for (i = 0; i < nbins; i++) d4->SetBinContent(i + 1,dest[i]);
2724 /// d4->SetLineColor(kMagenta);
2725 /// d4->Draw("SAME");
2726 /// printf("Found %d candidate peaks\n",nfound);
2727 /// }
2728 /// ~~~
2729 
2730 Int_t TSpectrum::SearchHighRes(Double_t *source,Double_t *destVector, int ssize,
2731  Double_t sigma, Double_t threshold,
2732  bool backgroundRemove,int deconIterations,
2733  bool markov, int averWindow)
2734 {
2735  int i, j, numberIterations = (Int_t)(7 * sigma + 0.5);
2736  Double_t a, b, c;
2737  int k, lindex, posit, imin, imax, jmin, jmax, lh_gold, priz;
2738  Double_t lda, ldb, ldc, area, maximum, maximum_decon;
2739  int xmin, xmax, l, peak_index = 0, size_ext = ssize + 2 * numberIterations, shift = numberIterations, bw = 2, w;
2740  Double_t maxch;
2741  Double_t nom, nip, nim, sp, sm, plocha = 0;
2742  Double_t m0low=0,m1low=0,m2low=0,l0low=0,l1low=0,detlow,av,men;
2743  if (sigma < 1) {
2744  Error("SearchHighRes", "Invalid sigma, must be greater than or equal to 1");
2745  return 0;
2746  }
2747 
2748  if(threshold<=0 || threshold>=100){
2749  Error("SearchHighRes", "Invalid threshold, must be positive and less than 100");
2750  return 0;
2751  }
2752 
2753  j = (Int_t) (5.0 * sigma + 0.5);
2754  if (j >= PEAK_WINDOW / 2) {
2755  Error("SearchHighRes", "Too large sigma");
2756  return 0;
2757  }
2758 
2759  if (markov == true) {
2760  if (averWindow <= 0) {
2761  Error("SearchHighRes", "Averaging window must be positive");
2762  return 0;
2763  }
2764  }
2765 
2766  if(backgroundRemove == true){
2767  if(ssize < 2 * numberIterations + 1){
2768  Error("SearchHighRes", "Too large clipping window");
2769  return 0;
2770  }
2771  }
2772 
2773  k = int(2 * sigma+0.5);
2774  if(k >= 2){
2775  for(i = 0;i < k;i++){
2776  a = i,b = source[i];
2777  m0low += 1,m1low += a,m2low += a * a,l0low += b,l1low += a * b;
2778  }
2779  detlow = m0low * m2low - m1low * m1low;
2780  if(detlow != 0)
2781  l1low = (-l0low * m1low + l1low * m0low) / detlow;
2782 
2783  else
2784  l1low = 0;
2785  if(l1low > 0)
2786  l1low=0;
2787  }
2788 
2789  else{
2790  l1low = 0;
2791  }
2792 
2793  i = (Int_t)(7 * sigma + 0.5);
2794  i = 2 * i;
2795  Double_t *working_space = new Double_t [7 * (ssize + i)];
2796  for (j=0;j<7 * (ssize + i);j++) working_space[j] = 0;
2797  for(i = 0; i < size_ext; i++){
2798  if(i < shift){
2799  a = i - shift;
2800  working_space[i + size_ext] = source[0] + l1low * a;
2801  if(working_space[i + size_ext] < 0)
2802  working_space[i + size_ext]=0;
2803  }
2804 
2805  else if(i >= ssize + shift){
2806  a = i - (ssize - 1 + shift);
2807  working_space[i + size_ext] = source[ssize - 1];
2808  if(working_space[i + size_ext] < 0)
2809  working_space[i + size_ext]=0;
2810  }
2811 
2812  else
2813  working_space[i + size_ext] = source[i - shift];
2814  }
2815 
2816  if(backgroundRemove == true){
2817  for(i = 1; i <= numberIterations; i++){
2818  for(j = i; j < size_ext - i; j++){
2819  if(markov == false){
2820  a = working_space[size_ext + j];
2821  b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
2822  if(b < a)
2823  a = b;
2824 
2825  working_space[j]=a;
2826  }
2827 
2828  else{
2829  a = working_space[size_ext + j];
2830  av = 0;
2831  men = 0;
2832  for (w = j - bw; w <= j + bw; w++){
2833  if ( w >= 0 && w < size_ext){
2834  av += working_space[size_ext + w];
2835  men +=1;
2836  }
2837  }
2838  av = av / men;
2839  b = 0;
2840  men = 0;
2841  for (w = j - i - bw; w <= j - i + bw; w++){
2842  if ( w >= 0 && w < size_ext){
2843  b += working_space[size_ext + w];
2844  men +=1;
2845  }
2846  }
2847  b = b / men;
2848  c = 0;
2849  men = 0;
2850  for (w = j + i - bw; w <= j + i + bw; w++){
2851  if ( w >= 0 && w < size_ext){
2852  c += working_space[size_ext + w];
2853  men +=1;
2854  }
2855  }
2856  c = c / men;
2857  b = (b + c) / 2;
2858  if (b < a)
2859  av = b;
2860  working_space[j]=av;
2861  }
2862  }
2863  for(j = i; j < size_ext - i; j++)
2864  working_space[size_ext + j] = working_space[j];
2865  }
2866  for(j = 0;j < size_ext; j++){
2867  if(j < shift){
2868  a = j - shift;
2869  b = source[0] + l1low * a;
2870  if (b < 0) b = 0;
2871  working_space[size_ext + j] = b - working_space[size_ext + j];
2872  }
2873 
2874  else if(j >= ssize + shift){
2875  a = j - (ssize - 1 + shift);
2876  b = source[ssize - 1];
2877  if (b < 0) b = 0;
2878  working_space[size_ext + j] = b - working_space[size_ext + j];
2879  }
2880 
2881  else{
2882  working_space[size_ext + j] = source[j - shift] - working_space[size_ext + j];
2883  }
2884  }
2885  for(j = 0;j < size_ext; j++){
2886  if(working_space[size_ext + j] < 0) working_space[size_ext + j] = 0;
2887  }
2888  }
2889 
2890  for(i = 0; i < size_ext; i++){
2891  working_space[i + 6*size_ext] = working_space[i + size_ext];
2892  }
2893 
2894  if(markov == true){
2895  for(j = 0; j < size_ext; j++)
2896  working_space[2 * size_ext + j] = working_space[size_ext + j];
2897  xmin = 0,xmax = size_ext - 1;
2898  for(i = 0, maxch = 0; i < size_ext; i++){
2899  working_space[i] = 0;
2900  if(maxch < working_space[2 * size_ext + i])
2901  maxch = working_space[2 * size_ext + i];
2902  plocha += working_space[2 * size_ext + i];
2903  }
2904  if(maxch == 0) {
2905  delete [] working_space;
2906  return 0;
2907  }
2908 
2909  nom = 1;
2910  working_space[xmin] = 1;
2911  for(i = xmin; i < xmax; i++){
2912  nip = working_space[2 * size_ext + i] / maxch;
2913  nim = working_space[2 * size_ext + i + 1] / maxch;
2914  sp = 0,sm = 0;
2915  for(l = 1; l <= averWindow; l++){
2916  if((i + l) > xmax)
2917  a = working_space[2 * size_ext + xmax] / maxch;
2918 
2919  else
2920  a = working_space[2 * size_ext + i + l] / maxch;
2921 
2922  b = a - nip;
2923  if(a + nip <= 0)
2924  a=1;
2925 
2926  else
2927  a = TMath::Sqrt(a + nip);
2928 
2929  b = b / a;
2930  b = TMath::Exp(b);
2931  sp = sp + b;
2932  if((i - l + 1) < xmin)
2933  a = working_space[2 * size_ext + xmin] / maxch;
2934 
2935  else
2936  a = working_space[2 * size_ext + i - l + 1] / maxch;
2937 
2938  b = a - nim;
2939  if(a + nim <= 0)
2940  a = 1;
2941 
2942  else
2943  a = TMath::Sqrt(a + nim);
2944 
2945  b = b / a;
2946  b = TMath::Exp(b);
2947  sm = sm + b;
2948  }
2949  a = sp / sm;
2950  a = working_space[i + 1] = working_space[i] * a;
2951  nom = nom + a;
2952  }
2953  for(i = xmin; i <= xmax; i++){
2954  working_space[i] = working_space[i] / nom;
2955  }
2956  for(j = 0; j < size_ext; j++)
2957  working_space[size_ext + j] = working_space[j] * plocha;
2958  for(j = 0; j < size_ext; j++){
2959  working_space[2 * size_ext + j] = working_space[size_ext + j];
2960  }
2961  if(backgroundRemove == true){
2962  for(i = 1; i <= numberIterations; i++){
2963  for(j = i; j < size_ext - i; j++){
2964  a = working_space[size_ext + j];
2965  b = (working_space[size_ext + j - i] + working_space[size_ext + j + i]) / 2.0;
2966  if(b < a)
2967  a = b;
2968  working_space[j] = a;
2969  }
2970  for(j = i; j < size_ext - i; j++)
2971  working_space[size_ext + j] = working_space[j];
2972  }
2973  for(j = 0; j < size_ext; j++){
2974  working_space[size_ext + j] = working_space[2 * size_ext + j] - working_space[size_ext + j];
2975  }
2976  }
2977  }
2978 //deconvolution starts
2979  area = 0;
2980  lh_gold = -1;
2981  posit = 0;
2982  maximum = 0;
2983 //generate response vector
2984  for(i = 0; i < size_ext; i++){
2985  lda = (Double_t)i - 3 * sigma;
2986  lda = lda * lda / (2 * sigma * sigma);
2987  j = (Int_t)(1000 * TMath::Exp(-lda));
2988  lda = j;
2989  if(lda != 0)
2990  lh_gold = i + 1;
2991 
2992  working_space[i] = lda;
2993  area = area + lda;
2994  if(lda > maximum){
2995  maximum = lda;
2996  posit = i;
2997  }
2998  }
2999 //read source vector
3000  for(i = 0; i < size_ext; i++)
3001  working_space[2 * size_ext + i] = TMath::Abs(working_space[size_ext + i]);
3002 //create matrix at*a(vector b)
3003  i = lh_gold - 1;
3004  if(i > size_ext)
3005  i = size_ext;
3006 
3007  imin = -i,imax = i;
3008  for(i = imin; i <= imax; i++){
3009  lda = 0;
3010  jmin = 0;
3011  if(i < 0)
3012  jmin = -i;
3013  jmax = lh_gold - 1 - i;
3014  if(jmax > (lh_gold - 1))
3015  jmax = lh_gold - 1;
3016 
3017  for(j = jmin;j <= jmax; j++){
3018  ldb = working_space[j];
3019  ldc = working_space[i + j];
3020  lda = lda + ldb * ldc;
3021  }
3022  working_space[size_ext + i - imin] = lda;
3023  }
3024 //create vector p
3025  i = lh_gold - 1;
3026  imin = -i,imax = size_ext + i - 1;
3027  for(i = imin; i <= imax; i++){
3028  lda = 0;
3029  for(j = 0; j <= (lh_gold - 1); j++){
3030  ldb = working_space[j];
3031  k = i + j;
3032  if(k >= 0 && k < size_ext){
3033  ldc = working_space[2 * size_ext + k];
3034  lda = lda + ldb * ldc;
3035  }
3036 
3037  }
3038  working_space[4 * size_ext + i - imin] = lda;
3039  }
3040 //move vector p
3041  for(i = imin; i <= imax; i++)
3042  working_space[2 * size_ext + i - imin] = working_space[4 * size_ext + i - imin];
3043 //initialization of resulting vector
3044  for(i = 0; i < size_ext; i++)
3045  working_space[i] = 1;
3046 //START OF ITERATIONS
3047  for(lindex = 0; lindex < deconIterations; lindex++){
3048  for(i = 0; i < size_ext; i++){
3049  if(TMath::Abs(working_space[2 * size_ext + i]) > 0.00001 && TMath::Abs(working_space[i]) > 0.00001){
3050  lda=0;
3051  jmin = lh_gold - 1;
3052  if(jmin > i)
3053  jmin = i;
3054 
3055  jmin = -jmin;
3056  jmax = lh_gold - 1;
3057  if(jmax > (size_ext - 1 - i))
3058  jmax=size_ext-1-i;
3059 
3060  for(j = jmin; j <= jmax; j++){
3061  ldb = working_space[j + lh_gold - 1 + size_ext];
3062  ldc = working_space[i + j];
3063  lda = lda + ldb * ldc;
3064  }
3065  ldb = working_space[2 * size_ext + i];
3066  if(lda != 0)
3067  lda = ldb / lda;
3068 
3069  else
3070  lda = 0;
3071 
3072  ldb = working_space[i];
3073  lda = lda * ldb;
3074  working_space[3 * size_ext + i] = lda;
3075  }
3076  }
3077  for(i = 0; i < size_ext; i++){
3078  working_space[i] = working_space[3 * size_ext + i];
3079  }
3080  }
3081 //shift resulting spectrum
3082  for(i=0;i<size_ext;i++){
3083  lda = working_space[i];
3084  j = i + posit;
3085  j = j % size_ext;
3086  working_space[size_ext + j] = lda;
3087  }
3088 //write back resulting spectrum
3089  maximum = 0, maximum_decon = 0;
3090  j = lh_gold - 1;
3091  for(i = 0; i < size_ext - j; i++){
3092  if(i >= shift && i < ssize + shift){
3093  working_space[i] = area * working_space[size_ext + i + j];
3094  if(maximum_decon < working_space[i])
3095  maximum_decon = working_space[i];
3096  if(maximum < working_space[6 * size_ext + i])
3097  maximum = working_space[6 * size_ext + i];
3098  }
3099 
3100  else
3101  working_space[i] = 0;
3102  }
3103  lda=1;
3104  if(lda>threshold)
3105  lda=threshold;
3106  lda=lda/100;
3107 
3108 //searching for peaks in deconvolved spectrum
3109  for(i = 1; i < size_ext - 1; i++){
3110  if(working_space[i] > working_space[i - 1] && working_space[i] > working_space[i + 1]){
3111  if(i >= shift && i < ssize + shift){
3112  if(working_space[i] > lda*maximum_decon && working_space[6 * size_ext + i] > threshold * maximum / 100.0){
3113  for(j = i - 1, a = 0, b = 0; j <= i + 1; j++){
3114  a += (Double_t)(j - shift) * working_space[j];
3115  b += working_space[j];
3116  }
3117  a = a / b;
3118  if(a < 0)
3119  a = 0;
3120 
3121  if(a >= ssize)
3122  a = ssize - 1;
3123  if(peak_index == 0){
3124  fPositionX[0] = a;
3125  peak_index = 1;
3126  }
3127 
3128  else{
3129  for(j = 0, priz = 0; j < peak_index && priz == 0; j++){
3130  if(working_space[6 * size_ext + shift + (Int_t)a] > working_space[6 * size_ext + shift + (Int_t)fPositionX[j]])
3131  priz = 1;
3132  }
3133  if(priz == 0){
3134  if(j < fMaxPeaks){
3135  fPositionX[j] = a;
3136  }
3137  }
3138 
3139  else{
3140  for(k = peak_index; k >= j; k--){
3141  if(k < fMaxPeaks){
3142  fPositionX[k] = fPositionX[k - 1];
3143  }
3144  }
3145  fPositionX[j - 1] = a;
3146  }
3147  if(peak_index < fMaxPeaks)
3148  peak_index += 1;
3149  }
3150  }
3151  }
3152  }
3153  }
3154 
3155  for(i = 0; i < ssize; i++) destVector[i] = working_space[i + shift];
3156  delete [] working_space;
3157  fNPeaks = peak_index;
3158  if(peak_index == fMaxPeaks)
3159  Warning("SearchHighRes", "Peak buffer full");
3160  return fNPeaks;
3161 }
3162 
3163 ////////////////////////////////////////////////////////////////////////////////
3164 /// Old name of SearcHighRes introduced for back compatibility.
3165 /// This function will be removed after the June 2006 release
3166 
3167 Int_t TSpectrum::Search1HighRes(Double_t *source,Double_t *destVector, int ssize,
3168  Double_t sigma, Double_t threshold,
3169  bool backgroundRemove,int deconIterations,
3170  bool markov, int averWindow)
3171 {
3172 
3173 
3174  return SearchHighRes(source,destVector,ssize,sigma,threshold,backgroundRemove,
3175  deconIterations,markov,averWindow);
3176 }
3177 
3178 ////////////////////////////////////////////////////////////////////////////////
3179 /// Static function, interface to TSpectrum::Search.
3180 
3182 {
3183 
3184  TSpectrum s;
3185  return s.Search(hist,sigma,option,threshold);
3186 }
3187 
3188 ////////////////////////////////////////////////////////////////////////////////
3189 /// Static function, interface to TSpectrum::Background.
3190 
3191 TH1 *TSpectrum::StaticBackground(const TH1 *hist,Int_t niter, Option_t *option)
3192 {
3193  TSpectrum s;
3194  return s.Background(hist,niter,option);
3195 }
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:429
TList * GetListOfFunctions() const
Definition: TH1.h:244
Int_t Search1HighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
Old name of SearcHighRes introduced for back compatibility.
Definition: TSpectrum.cxx:3167
float xmin
Definition: THbookFile.cxx:93
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4629
void f6()
Definition: na49.C:122
return c
const char Option_t
Definition: RtypesCore.h:62
Definition: Rtypes.h:61
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
virtual Int_t GetDimension() const
Definition: TH1.h:283
TH1 * h
Definition: legend2.C:5
static TH1 * StaticBackground(const TH1 *hist, Int_t niter=20, Option_t *option="")
Static function, interface to TSpectrum::Background.
Definition: TSpectrum.cxx:3191
const char * DeconvolutionRL(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
Definition: TSpectrum.cxx:2133
Basic string class.
Definition: TString.h:137
ClassImp(TSpectrum) TSpectrum
Constructor.
Definition: TSpectrum.cxx:41
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1075
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:497
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TH1.cxx:2565
const char * Data() const
Definition: TString.h:349
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
Double_t * fPositionX
[fNPeaks] X position of peaks
Definition: TSpectrum.h:30
int d
Definition: tornado.py:11
Int_t SearchHighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
One-dimensional high-resolution peak search function.
Definition: TSpectrum.cxx:2730
virtual void SetMarkerColor(Color_t mcolor=1)
Definition: TAttMarker.h:51
const Double_t sigma
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
virtual void Delete(Option_t *option="")
Delete this object.
Definition: TObject.cxx:228
virtual Double_t GetBinCenter(Int_t bin) const
return bin center for 1D historam Better to use h1.GetXaxis().GetBinCenter(bin)
Definition: TH1.cxx:8470
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
One-dimensional peak search function.
Definition: TSpectrum.cxx:262
virtual ~TSpectrum()
Destructor.
Definition: TSpectrum.cxx:83
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:675
LVector boost(const LVector &v, const BoostVector &b)
Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost The only re...
Definition: VectorUtil.h:329
TH1 * fHistogram
resulting histogram
Definition: TSpectrum.h:33
Double_t * fPositionY
[fNPeaks] Y position of peaks
Definition: TSpectrum.h:31
const char * Deconvolution(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
Definition: TSpectrum.cxx:1853
bool first
Definition: line3Dfit.C:48
tuple w
Definition: qtexample.py:51
TLine * l
Definition: textangle.C:4
static void SetAverageWindow(Int_t w=3)
Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes).
Definition: TSpectrum.cxx:95
static Int_t fgAverageWindow
Average window of searched peaks.
Definition: TSpectrum.h:34
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
float xmax
Definition: THbookFile.cxx:93
void f8()
Definition: na49.C:163
#define PEAK_WINDOW
Definition: TSpectrum.cxx:40
static void SetDeconIterations(Int_t n=3)
Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::Search...
Definition: TSpectrum.cxx:104
virtual void SetMarkerSize(Size_t msize=1)
Definition: TAttMarker.h:54
Double_t Exp(Double_t x)
Definition: TMath.h:495
A PolyMarker is defined by an array on N points in a 2-D space.
Definition: TPolyMarker.h:37
double Double_t
Definition: RtypesCore.h:55
virtual TH1 * Background(const TH1 *hist, Int_t niter=20, Option_t *option="")
One-dimensional background estimation function.
Definition: TSpectrum.cxx:148
static Int_t fgIterations
Maximum number of decon iterations (default=3)
Definition: TSpectrum.h:35
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Int_t fMaxPeaks
Maximum number of peaks to be found.
Definition: TSpectrum.h:27
The TH1 histogram class.
Definition: TH1.h:80
void SetResolution(Double_t resolution=1)
resolution: determines resolution of the neighbouring peaks default value is 1 correspond to 3 sigma ...
Definition: TSpectrum.cxx:346
Advanced Spectra Processing.
Definition: TSpectrum.h:20
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:440
#define hbname
Definition: THbookFile.cxx:123
Double_t * fPosition
[fNPeaks] array of current peak positions
Definition: TSpectrum.h:29
virtual void Print(Option_t *option="") const
Print the array of positions.
Definition: TSpectrum.cxx:215
virtual void Add(TObject *obj)
Definition: TList.h:81
#define dest(otri, vertexptr)
Definition: triangle.c:1040
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
static Int_t StaticSearch(const TH1 *hist, Double_t sigma=2, Option_t *option="goff", Double_t threshold=0.05)
Static function, interface to TSpectrum::Search.
Definition: TSpectrum.cxx:3181
#define gPad
Definition: TVirtualPad.h:288
Int_t fNPeaks
number of peaks found
Definition: TSpectrum.h:28
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
const Bool_t kTRUE
Definition: Rtypes.h:91
const char * SmoothMarkov(Double_t *source, Int_t ssize, Int_t averWindow)
One-dimensional markov spectrum smoothing function.
Definition: TSpectrum.cxx:1484
const Int_t n
Definition: legend1.C:16
const char * Unfolding(Double_t *source, const Double_t **respMatrix, Int_t ssizex, Int_t ssizey, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional unfolding function.
Definition: TSpectrum.cxx:2358
TH1F * background
Definition: fithist.C:4
Double_t fResolution
resolution of the neighboring peaks
Definition: TSpectrum.h:32
TAxis * GetXaxis()
Definition: TH1.h:319
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904