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