GSI Object Oriented Online Offline (Go4)  GO4-6.3.0
TGo4FitPeakFinder.cxx
Go to the documentation of this file.
1 // $Id$
2 //-----------------------------------------------------------------------
3 // The GSI Online Offline Object Oriented (Go4) Project
4 // Experiment Data Processing at EE department, GSI
5 //-----------------------------------------------------------------------
6 // Copyright (C) 2000- GSI Helmholtzzentrum fuer Schwerionenforschung GmbH
7 // Planckstr. 1, 64291 Darmstadt, Germany
8 // Contact: http://go4.gsi.de
9 //-----------------------------------------------------------------------
10 // This software can be used under the license agreements as stated
11 // in Go4License.txt file which is part of the distribution.
12 //-----------------------------------------------------------------------
13 
14 #include "TGo4FitPeakFinder.h"
15 
16 #include "TMath.h"
17 #include "TH1D.h"
18 #include "TMatrixD.h"
19 #include "TVectorD.h"
20 #include "TArrayD.h"
21 #include "TArrayI.h"
22 #include "TSpectrum.h"
23 
24 #include "TGo4Fitter.h"
25 #include "TGo4FitData.h"
26 
28 
29 TGo4FitPeakFinder::TGo4FitPeakFinder(const char *Name, const char *DataName, Bool_t ClearModels, Int_t PolOrder)
30  : TGo4FitterAction(Name, "Peak finder action"), fiPeakFinderType(0), fxDataName(DataName),
31  fbClearModels(ClearModels), fbUsePolynom(PolOrder >= 0), fiPolynomOrder(PolOrder >= 0 ? PolOrder : 0),
34 {
35 }
36 
38 
40 {
41  if (PolynomOrder < 0)
42  SetUsePolynom(kFALSE);
43  else {
44  SetUsePolynom(kTRUE);
45  SetPolynomOrder(PolynomOrder);
46  }
47 }
48 
49 void TGo4FitPeakFinder::SetupForFirst(Double_t MaxAmplFactor, Double_t MinWidth, Double_t MaxWidth)
50 {
52  Set0MaxAmplFactor(MaxAmplFactor);
53  Set0MinWidth(MinWidth);
54  Set0MaxWidth(MaxWidth);
55 }
56 
57 void TGo4FitPeakFinder::SetupForSecond(Double_t LineWidth)
58 {
60  Set1LineWidth(LineWidth);
61 }
62 
63 void TGo4FitPeakFinder::SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum)
64 {
66  Set2NoiseFactor(NoiseFactor);
67  Set2NoiseMinimum(NoiseMinimum);
68  Set2ChannelSum(ChannelSum);
69 }
70 
72 {
73  TGo4Fitter *fitter = dynamic_cast<TGo4Fitter *>(Fitter);
74  if (!fitter)
75  return;
76 
77  TGo4FitData *data = fitter->FindData(GetDataName());
78  if (!data)
79  return;
80 
81  if (GetClearModels())
82  fitter->DeleteModelsAssosiatedTo(data->GetName());
83 
84  if (GetPeakFinderType() == 0)
86  Get0MaxWidth());
87 
88  if (GetPeakFinderType() == 1)
89  ROOTPeakFinder(fitter, data, GetUsePolynom() ? GetPolynomOrder() : -1, Get1LineWidth());
90 
91  if (GetPeakFinderType() == 2)
93  GetUsePolynom() ? GetPolynomOrder() : -1);
94 }
95 
96 void TGo4FitPeakFinder::Print(Option_t *option) const
97 {
99 }
100 
101 void TGo4FitPeakFinder::ROOTPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolynomOrder, Double_t Sigma)
102 {
103  if (!fitter || !data)
104  return;
105 
106  auto iter = data->MakeIter();
107  if (!iter)
108  return;
109 
110  Int_t size = iter->CountPoints();
111 
112  if ((size < 10) || (iter->ScalesSize() != 1))
113  return;
114 
115  TArrayD Bins(size), Scales(size), HScales(size + 1);
116 
117  Int_t nbin = 0;
118  if (iter->Reset())
119  do {
120  Bins[nbin] = data->GetAmplValue() * iter->Value();
121  Double_t pos = iter->x();
122  Double_t width = iter->HasWidths() ? iter->Widths()[0] : 1.;
123  Scales[nbin] = pos;
124  HScales[nbin] = pos - width / 2.;
125  HScales[nbin + 1] = pos + width / 2;
126  nbin++;
127  } while (iter->Next());
128 
129  if (PolynomOrder >= 0)
130  fitter->AddPolynomX(data->GetName(), "Pol", PolynomOrder);
131 
132  TH1D histo("ROOTPeakFinder", "ROOTPeakFinder spectrum", size, HScales.GetArray());
133  histo.SetDirectory(nullptr);
134  for (Int_t i = 0; i < size; i++)
135  histo.SetBinContent(i + 1, Bins[i]);
136 
137  TSpectrum sp(100);
138  sp.Search(&histo, Sigma);
139  // sp.Search1(Bins.GetArray(), size, Sigma);
140 
141  for (Int_t n = 0; n < sp.GetNPeaks(); n++) {
142  Double_t pos = sp.GetPositionX()[n];
143  Double_t width = Sigma;
144 
145  // this is old code for ROOT 3.*, may be usefull for somebody
146  /*
147  Double_t dindx = sp.GetPositionX()[n];
148  Int_t left = TMath::Nint(dindx-Sigma/2.);
149  Int_t right = TMath::Nint(dindx+Sigma/2.);
150  if (left<0) left = 0;
151  if (right>=size) right = size-1;
152 
153  Double_t pos = Scales[TMath::Nint(dindx)];
154  Double_t width = TMath::Abs(Scales[right]-Scales[left]);
155  */
156  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss", n), pos, width);
157  }
158 }
159 
160 Double_t TGo4FitPeakFinder::CalcPolynom(const TArrayD &Coef, Double_t x)
161 {
162  Double_t power = 1., sum = 0.;
163  for (Int_t n = 0; n < Coef.GetSize(); n++) {
164  sum += Coef[n] * power;
165  power *= x;
166  }
167  return sum;
168 }
169 
170 void TGo4FitPeakFinder::DefinePolynom(Int_t size, Double_t *bins, Double_t *scales, TArrayD &Coef, Double_t *weight,
171  Double_t *backgr, Char_t *use)
172 {
173  Int_t maxorder = Coef.GetSize() - 1;
174 
175  TMatrixD matr(maxorder + 1, maxorder + 1);
176  TVectorD vector(maxorder + 1);
177 
178  for (Int_t order1 = 0; order1 <= maxorder; order1++)
179  for (Int_t order2 = order1; order2 <= maxorder; order2++) {
180  Double_t sum = 0;
181  for (Int_t i = 0; i < size; i++)
182  if (!use || use[i]) {
183  Double_t zn = TMath::Power(scales[i], order1 + order2);
184  if (weight)
185  zn *= weight[i];
186  sum += zn;
187  }
188  matr(order1, order2) = sum;
189  matr(order2, order1) = sum;
190  }
191 
192  for (Int_t order1 = 0; order1 <= maxorder; order1++) {
193  Double_t sum = 0;
194  for (Int_t i = 0; i < size; i++)
195  if (!use || use[i]) {
196  Double_t zn = bins[i];
197  if (backgr)
198  zn -= backgr[i];
199  zn *= TMath::Power(scales[i], order1);
200  if (weight)
201  zn *= weight[i];
202  sum += zn;
203  }
204  vector(order1) = sum;
205  }
206 
207  matr.Invert();
208  vector *= matr;
209 
210  for (Int_t order = 0; order <= maxorder; order++)
211  Coef[order] = vector(order);
212 }
213 
214 void TGo4FitPeakFinder::DefinePolynomEx(Int_t size, Double_t *bins, Double_t *scales, Double_t *weight,
215  Double_t *backgr, Int_t lbound, Int_t rbound, TArrayD &Coef)
216 {
217  TArrayC use(size);
218  use.Reset(0);
219  if (lbound < 0)
220  lbound = 0;
221  if (rbound >= size)
222  rbound = size - 1;
223  for (Int_t i = lbound; i <= rbound; i++)
224  use[i] = 1;
225  DefinePolynom(size, bins, scales, Coef, weight, backgr, use.GetArray());
226 }
227 
228 Bool_t FindValue(Double_t *bins, Double_t *backgr, Char_t *use, Int_t size, Int_t &p, Int_t dir, Double_t value)
229 {
230  p += dir;
231  while ((p >= 0) && (p < size) && use[p]) {
232  Double_t zn = bins[p] - backgr[p];
233  if (zn < value)
234  return kTRUE;
235  p += dir;
236  }
237  return kFALSE;
238 }
239 
241  Double_t AmplThreshold, Double_t MinWidth, Double_t MaxWidth)
242 {
243  if (!fitter || !data)
244  return;
245 
246  auto iter = data->MakeIter();
247  if (!iter)
248  return;
249 
250  Int_t size = iter->CountPoints();
251 
252  if ((size < 10) || (iter->ScalesSize() != 1))
253  return;
254 
255  TArrayD Bins(size), Scales(size), Weights(size), Background(size);
256 
257  Weights.Reset(1.);
258  Background.Reset(0.);
259 
260  Double_t *bins = Bins.GetArray();
261  Double_t *scales = Scales.GetArray();
262  Double_t *weights = Weights.GetArray();
263  Double_t *backgr = Background.GetArray();
264 
265  Int_t nbin = 0;
266  if (iter->Reset())
267  do {
268  bins[nbin] = data->GetAmplValue() * iter->Value();
269  if (iter->StandardDeviation() > 0)
270  weights[nbin] = 1. / iter->StandardDeviation();
271  scales[nbin] = iter->x();
272  nbin++;
273  } while (iter->Next());
274 
275  if ((PolOrder >= 0) && (PolOrder < (size + 1) * 5)) {
276  Int_t nsegments = (PolOrder + 1) * 2;
277  TArrayC usage(size);
278  usage.Reset(0);
279  for (Int_t nseg = 0; nseg < nsegments; nseg++) {
280  Int_t lbound = Int_t((nseg + 0.) / nsegments * size);
281  Int_t rbound = Int_t((nseg + 1.) / nsegments * size - 1.);
282  Int_t nselect = (rbound - lbound) / 10;
283  if (nselect < 1)
284  nselect = 1;
285  Int_t nsel = 0;
286  while (nsel++ < nselect) {
287  Int_t pmin = -1;
288  for (Int_t i = lbound; i <= rbound; i++)
289  if (!usage[i])
290  if ((pmin < 0) || (bins[i] < bins[pmin]))
291  pmin = i;
292  usage[pmin] = 1;
293  }
294  }
295 
296  TArrayD Coef(PolOrder + 1);
297  DefinePolynom(size, bins, scales, Coef, weights, nullptr, usage.GetArray());
298 
299  for (Int_t i = 0; i < size; i++)
300  backgr[i] = CalcPolynom(Coef, scales[i]);
301 
302  fitter->AddPolynomX(data->GetName(), "Pol", Coef);
303  }
304 
305  Double_t dmax = 0.;
306  for (Int_t i = 0; i < size; i++) {
307  Double_t zn = bins[i] - backgr[i];
308  if (zn > dmax)
309  dmax = zn;
310  }
311  AmplThreshold = AmplThreshold * dmax;
312 
313  TArrayC usage(size);
314  usage.Reset(1);
315  Char_t *use = usage.GetArray();
316 
317  Bool_t validline = kFALSE;
318 
319  do {
320 
321  validline = kFALSE;
322 
323  Int_t pmax = -1;
324  Double_t max = 0.;
325  for (Int_t i = 0; i < size; i++)
326  if (use[i]) {
327  Double_t zn = bins[i] - backgr[i];
328  if ((pmax < 0) || (zn > max)) {
329  pmax = i;
330  max = zn;
331  }
332  }
333 
334  Double_t sigma = TMath::Sqrt(1. / weights[pmax]);
335 
336  if ((max < 2. * sigma) || (max < AmplThreshold))
337  break;
338 
339  Int_t lbound = pmax, rbound = pmax;
340 
341  Bool_t lok = FindValue(bins, backgr, use, size, lbound, -1, max * 0.5);
342  Bool_t rok = FindValue(bins, backgr, use, size, rbound, +1, max * 0.5);
343 
344  if (lok || rok) {
345 
346  validline = kTRUE;
347 
348  Int_t w = (lok && rok) ? (rbound - lbound) / 4 : (lok ? (pmax - lbound) / 2 : (rbound - pmax) / 2);
349  if (w < 2)
350  w = 2;
351 
352  TArrayD Coef1(3), Coef2(2);
353  DefinePolynomEx(size, bins, scales, weights, backgr, pmax - w, pmax + w, Coef1);
354  Double_t middle = -0.5 * Coef1[1] / Coef1[2];
355  Double_t amplitude = CalcPolynom(Coef1, middle);
356 
357  Double_t left = 0, right = 0;
358  if (lok) {
359  DefinePolynomEx(size, bins, scales, weights, backgr, lbound - w, lbound + w, Coef2);
360  left = (amplitude * .5 - Coef2[0]) / Coef2[1];
361  }
362 
363  if (rok) {
364  DefinePolynomEx(size, bins, scales, weights, backgr, rbound - w, rbound + w, Coef2);
365  right = (amplitude * .5 - Coef2[0]) / Coef2[1];
366  }
367 
368  if (!lok)
369  left = 2 * middle - right;
370  else if (!rok)
371  right = 2 * middle - left;
372 
373  Double_t width = (right - left) * 0.5;
374 
375  lbound = pmax - 4 * w;
376  if (lbound < 0)
377  lbound = 0;
378  rbound = pmax + 4 * w;
379  if (rbound >= size)
380  rbound = size - 1;
381  for (Int_t i = lbound; i <= rbound; i++)
382  use[i] = 0;
383 
384  if (lok && rok && (width >= MinWidth) && (width <= MaxWidth) && (amplitude > AmplThreshold) &&
385  (amplitude > 2 * sigma))
386  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss", 0), middle, width, amplitude);
387  }
388 
389  } while (validline);
390 }
391 
392 #include "f_find_peaks.c"
393 
394 void TGo4FitPeakFinder::HansEsselPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t MaxNumPeaks, Int_t ChannelSum,
395  Double_t NoiseFactor, Double_t NoiseMinimum, Int_t MinimasOrder)
396 {
397  if (!fitter || !data || (MaxNumPeaks < 1))
398  return;
399 
400  auto iter = data->MakeIter();
401  if (!iter)
402  return;
403 
404  Int_t size = iter->CountPoints();
405 
406  if ((size < 10) || (iter->ScalesSize() != 1))
407  return;
408 
409  TArrayD Bins(size), Scales(size);
410  Double_t *bins = Bins.GetArray();
411  Double_t *scales = Scales.GetArray();
412 
413  Int_t nbin = 0;
414  if (iter->Reset())
415  do {
416  bins[nbin] = data->GetAmplValue() * iter->Value();
417  scales[nbin] = iter->x();
418  nbin++;
419  } while (iter->Next());
420 
421  int NumPeaks = 0;
422  TArrayI Minimas(MaxNumPeaks), Maximas(MaxNumPeaks), MaximaWidths(MaxNumPeaks);
423  TArrayI PeaksLeft(MaxNumPeaks), PeaksRight(MaxNumPeaks), PeaksMaximum(MaxNumPeaks);
424  TArrayD MinimaValues(MaxNumPeaks), MinimaScales(MaxNumPeaks), MaximaIntegrals(MaxNumPeaks);
425 
426  go4fit_find_peaks(bins, TYPE__DOUBLE, 0, size, ChannelSum, NoiseFactor, NoiseMinimum, nullptr, MaxNumPeaks, &NumPeaks,
427  Minimas.GetArray(), MinimaValues.GetArray(), Maximas.GetArray(), MaximaWidths.GetArray(),
428  MaximaIntegrals.GetArray(), PeaksLeft.GetArray(), PeaksMaximum.GetArray(), PeaksRight.GetArray(),
429  nullptr, nullptr);
430 
431  if (NumPeaks <= 0)
432  return;
433 
434  Int_t NumMinimas = (NumPeaks < MaxNumPeaks) ? NumPeaks + 1 : MaxNumPeaks;
435  for (int n = 0; n < NumMinimas; n++) {
436  MinimaScales[n] = scales[Minimas[n]];
437  }
438 
439  if ((MinimasOrder >= 0) && (NumMinimas > 1)) {
440  Int_t order = MinimasOrder;
441  if (order >= NumMinimas)
442  order = NumMinimas - 1;
443  TArrayD Coef(order + 1);
444  DefinePolynom(NumMinimas, MinimaValues.GetArray(), MinimaScales.GetArray(), Coef);
445  fitter->AddPolynomX(data->GetName(), "Pol", Coef);
446  }
447 
448  for (int n = 0; n < NumPeaks; n++) {
449  Double_t pos = scales[PeaksMaximum[n]];
450  Double_t width = 0.;
451  if (PeaksRight[n] - PeaksLeft[n] < 1)
452  width = scales[PeaksMaximum[n]] - scales[PeaksMaximum[n] - 1];
453  else
454  width = TMath::Abs(scales[PeaksRight[n]] - scales[PeaksLeft[n]]) / 2.;
455  Double_t ampl = bins[PeaksMaximum[n]];
456 
457  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss", 0), pos, width, ampl);
458  }
459 }
Bool_t GetUsePolynom() const
void SetupForSecond(Double_t LineWidth)
void Set2NoiseFactor(Double_t factor)
Double_t Get0MinWidth() const
void Print(Option_t *option="") const override
Bool_t GetClearModels() const
Int_t GetPolynomOrder() const
void DoAction(TGo4FitterAbstract *Fitter) override
Int_t GetPeakFinderType() const
static void DefinePolynomEx(Int_t size, Double_t *bins, Double_t *scales, Double_t *weight, Double_t *backgr, Int_t lbound, Int_t rbound, TArrayD &Coef)
void AddPolynomX(const char *DataName, const char *NamePrefix, Int_t MaxOrder=1, Int_t GroupIndex=0, Double_t lrange=0., Double_t rrange=0.)
Definition: TGo4Fitter.cxx:235
static void DefinePolynom(Int_t size, Double_t *bins, Double_t *scales, TArrayD &Coef, Double_t *weight=nullptr, Double_t *backgr=nullptr, Char_t *use=nullptr)
Int_t Get2ChannelSum() const
Double_t Get2NoiseMinimum() const
void SergeyLinevPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolOrder, Double_t AmplThreshold, Double_t MinWidth, Double_t MaxWidth)
void go4fit_find_peaks(void *pfData, int lType, int lFirstChan, int lPoints, int lAver, double dDeltaFact, double dDeltaMin, double *dNoise, int lPeaksMax, int *plPeaks, int *plMinima, double *pdMinima, int *plMean, int *plWidth, double *pdIntegral, int *plLeft, int *plMaximum, int *plRight, float *pfLowNoise, float *pfHighNoise)
Definition: f_find_peaks.c:28
void Set2NoiseMinimum(Double_t min)
void SetPeakFinderType(Int_t typ)
Double_t Get0MaxWidth() const
void DeleteModelsAssosiatedTo(const char *DataName)
Definition: TGo4Fitter.cxx:440
Double_t Get1LineWidth() const
void Set0MaxWidth(Double_t max)
void Set0MaxAmplFactor(Double_t factor)
void SetupForFirst(Double_t MaxAmplFactor, Double_t MinWidth, Double_t MaxWidth)
TGo4FitModelGauss1 * AddGauss1(const char *DataName, const char *ModelName, Double_t iPosition, Double_t iWidth, Double_t iAmpl=1., Int_t Axis=0)
Definition: TGo4Fitter.cxx:382
void Set2ChannelSum(Int_t sum)
const char * GetDataName() const
Bool_t FindValue(Double_t *bins, Double_t *backgr, Char_t *use, Int_t size, Int_t &p, Int_t dir, Double_t value)
#define TYPE__DOUBLE
Definition: f_find_peaks.c:19
TString FindNextName(const char *Head, Int_t start, Bool_t isModel=kTRUE)
void SetupPolynomialBackground(Int_t PolynomOrder)
virtual std::unique_ptr< TGo4FitDataIter > MakeIter()
Definition: TGo4FitData.h:175
void SetUsePolynom(Bool_t use)
void SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum)
void Print(Option_t *option="") const override
void Set1LineWidth(Double_t width)
static Double_t CalcPolynom(const TArrayD &Coef, Double_t x)
Double_t Get2NoiseFactor() const
void ROOTPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolynomOrder, Double_t Sigma)
Double_t Get0MaxAmplFactor() const
TGo4FitData * FindData(const char *DataName)
Definition: TGo4Fitter.cxx:114
void SetPolynomOrder(Int_t order)
void HansEsselPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t MaxNumPeaks=50, Int_t ChannelSum=1, Double_t NoiseFactor=2., Double_t NoiseMinimum=10., Int_t MinimasOrder=-1)
void Set0MinWidth(Double_t min)
void usage(const char *subtopic=nullptr)