GSI Object Oriented Online Offline (Go4)  GO4-5.3.2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
TGo4FitPeakFinder.cxx
Go to the documentation of this file.
1 // $Id: TGo4FitPeakFinder.cxx 478 2009-10-29 12:26:09Z linev $
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 für 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 #include "TGo4FitModelPolynom.h"
27 
29 
30 }
31 
32 TGo4FitPeakFinder::TGo4FitPeakFinder(const char* Name, const char* DataName, Bool_t ClearModels, Int_t PolOrder) :
33  TGo4FitterAction(Name, "Peak finder action"),
34  fiPeakFinderType(0), fxDataName(DataName), fbClearModels(ClearModels),
35  fbUsePolynom(PolOrder>=0), fiPolynomOrder(PolOrder>=0 ? PolOrder : 0),
36  fd0MinWidth(5.), fd0MaxWidth(30.), fd0MaxAmplFactor(0.2),
37  fd1LineWidth(5.),
38  fd2NoiseFactor(2.), fd2NoiseMinimum(5.), fi2ChannelSum(2) {
39 }
40 
42 }
43 
45  if (PolynomOrder<0) SetUsePolynom(kFALSE);
46  else {
47  SetUsePolynom(kTRUE);
48  SetPolynomOrder(PolynomOrder);
49  }
50 }
51 
52 void TGo4FitPeakFinder::SetupForFirst(Double_t MaxAmplFactor, Double_t MinWidth, Double_t MaxWidth) {
54  Set0MaxAmplFactor(MaxAmplFactor);
55  Set0MinWidth(MinWidth);
56  Set0MaxWidth(MaxWidth);
57 }
58 
59 void TGo4FitPeakFinder::SetupForSecond(Double_t LineWidth) {
61  Set1LineWidth(LineWidth);
62 }
63 
64 void TGo4FitPeakFinder::SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum) {
66  Set2NoiseFactor(NoiseFactor);
67  Set2NoiseMinimum(NoiseMinimum);
68  Set2ChannelSum(ChannelSum);
69 }
70 
72  TGo4Fitter* fitter = dynamic_cast<TGo4Fitter*> (Fitter);
73  if (fitter==0) return;
74 
75  TGo4FitData* data = fitter->FindData(GetDataName());
76  if (data==0) return;
77 
78  if (GetClearModels())
79  fitter->DeleteModelsAssosiatedTo(data->GetName());
80 
81  if (GetPeakFinderType()==0)
82  SergeyLinevPeakFinder(fitter,
83  data,
84  GetUsePolynom() ? GetPolynomOrder() : -1,
86  Get0MinWidth(),
87  Get0MaxWidth());
88 
89  if (GetPeakFinderType()==1)
90  ROOTPeakFinder(fitter,
91  data,
92  GetUsePolynom() ? GetPolynomOrder() : -1,
93  Get1LineWidth());
94 
95  if (GetPeakFinderType()==2)
96  HansEsselPeakFinder(fitter,
97  data,
98  500,
100  Get2NoiseFactor(),
102  GetUsePolynom() ? GetPolynomOrder() : -1);
103 }
104 
105 void TGo4FitPeakFinder::Print(Option_t* option) const {
106  TGo4FitterAction::Print(option);
107 }
108 /*
109 void TGo4FitPeakFinder::ROOTPeakFinder(TGo4Fitter* fitter, TGo4FitData* data, Int_t PolynomOrder, Double_t Sigma) {
110  if ((fitter==0) || (data==0)) return;
111 
112  TGo4FitDataIter* iter = data->MakeIter();
113  if (iter==0) return;
114 
115  Int_t size = iter->CountPoints();
116 
117  if ((size<10) || (iter->ScalesSize()!=1)) { delete iter; return; }
118 
119  TArrayF Bins(size), Scales(size);
120 
121  Int_t nbin = 0;
122  if (iter->Reset()) do {
123  Bins[nbin] = data->GetAmplValue() * iter->Value();
124  Scales[nbin] = iter->x();
125  nbin++;
126  } while (iter->Next());
127 
128  delete iter;
129 
130  if (PolynomOrder>=0)
131  fitter->AddPolynomX(data->GetName(), "Pol", PolynomOrder);
132 
133  TSpectrum sp(100);
134  sp.Search1(Bins.GetArray(), size, Sigma);
135 
136  for(Int_t n=0;n<sp.GetNPeaks();n++) {
137  Double_t dindx = sp.GetPositionX()[n];
138  Int_t left = TMath::Nint(dindx-Sigma/2.);
139  Int_t right = TMath::Nint(dindx+Sigma/2.);
140 
141  Double_t pos = Scales[TMath::Nint(dindx)];
142  Double_t width = TMath::Abs(Scales[right]-Scales[left]);
143 
144  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",n), pos, width);
145  }
146 }
147 */
148 
149 void TGo4FitPeakFinder::ROOTPeakFinder(TGo4Fitter* fitter, TGo4FitData* data, Int_t PolynomOrder, Double_t Sigma) {
150  if ((fitter==0) || (data==0)) return;
151 
152  TGo4FitDataIter* iter = data->MakeIter();
153  if (iter==0) return;
154 
155  Int_t size = iter->CountPoints();
156 
157  if ((size<10) || (iter->ScalesSize()!=1)) { delete iter; return; }
158 
159  TArrayD Bins(size), Scales(size), HScales(size+1);
160 
161  Int_t nbin = 0;
162  if (iter->Reset()) do {
163  Bins[nbin] = data->GetAmplValue() * iter->Value();
164  Double_t pos = iter->x();
165  Double_t width = iter->HasWidths() ? iter->Widths()[0] : 1.;
166  Scales[nbin] = pos;
167  HScales[nbin] = pos - width/2.;
168  HScales[nbin+1] = pos + width/2;
169  nbin++;
170  } while (iter->Next());
171 
172  delete iter;
173 
174  if (PolynomOrder>=0)
175  fitter->AddPolynomX(data->GetName(), "Pol", PolynomOrder);
176 
177  TH1D histo("ROOTPeakFinder", "ROOTPeakFinder spectrum", size, HScales.GetArray());
178  histo.SetDirectory(0);
179  for (Int_t i=0;i<size;i++)
180  histo.SetBinContent(i+1, Bins[i]);
181 
182  TSpectrum sp(100);
183  sp.Search(&histo, Sigma);
184 // sp.Search1(Bins.GetArray(), size, Sigma);
185 
186  for(Int_t n=0;n<sp.GetNPeaks();n++) {
187  Double_t pos = sp.GetPositionX()[n];
188  Double_t width = Sigma;
189 
190 
191 // this is old code for ROOT 3.*, may be usefull for somebody
192 /*
193  Double_t dindx = sp.GetPositionX()[n];
194  Int_t left = TMath::Nint(dindx-Sigma/2.);
195  Int_t right = TMath::Nint(dindx+Sigma/2.);
196  if (left<0) left = 0;
197  if (right>=size) right = size-1;
198 
199  Double_t pos = Scales[TMath::Nint(dindx)];
200  Double_t width = TMath::Abs(Scales[right]-Scales[left]);
201 */
202  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",n), pos, width);
203  }
204 }
205 
206 
207 Double_t TGo4FitPeakFinder::CalcPolynom(const TArrayD& Coef, Double_t x) {
208  Double_t power = 1., sum = 0.;
209  for(Int_t n=0;n<Coef.GetSize();n++)
210  { sum+=Coef[n]*power; power*=x; }
211  return sum;
212 }
213 
214 void TGo4FitPeakFinder::DefinePolynom(Int_t size, Double_t* bins, Double_t* scales, TArrayD& Coef,
215  Double_t* weight, Double_t* backgr, Char_t* use) {
216  Int_t maxorder = Coef.GetSize()-1;
217 
218  TMatrixD matr(maxorder+1,maxorder+1);
219  TVectorD vector(maxorder+1);
220 
221  for(Int_t order1=0;order1<=maxorder;order1++)
222  for(Int_t order2=order1;order2<=maxorder;order2++) {
223  Double_t sum = 0;
224  for(Int_t i=0;i<size;i++)
225  if((use==0) || use[i]) {
226  Double_t zn = TMath::Power(scales[i],order1+order2);
227  if (weight) zn*=weight[i];
228  sum+=zn;
229  }
230  matr(order1,order2) = sum;
231  matr(order2,order1) = sum;
232  }
233 
234  for(Int_t order1=0;order1<=maxorder;order1++) {
235  Double_t sum = 0;
236  for(Int_t i=0;i<size;i++)
237  if((use==0) || use[i]) {
238  Double_t zn = bins[i];
239  if (backgr) zn-=backgr[i];
240  zn*=TMath::Power(scales[i],order1);
241  if (weight) zn*=weight[i];
242  sum+=zn;
243  }
244  vector(order1) = sum;
245  }
246 
247  matr.Invert();
248  vector*=matr;
249 
250  for(Int_t order=0;order<=maxorder;order++)
251  Coef[order]=vector(order);
252 }
253 
254 void TGo4FitPeakFinder::DefinePolynomEx(Int_t size, Double_t* bins, Double_t* scales, Double_t* weight, Double_t* backgr,
255  Int_t lbound, Int_t rbound, TArrayD& Coef) {
256  TArrayC use(size); use.Reset(0);
257  if (lbound<0) lbound = 0;
258  if (rbound>=size) rbound=size-1;
259  for(Int_t i=lbound;i<=rbound;i++) use[i]=1;
260  DefinePolynom(size, bins, scales, Coef, weight, backgr, use.GetArray());
261 }
262 
263 Bool_t FindValue(Double_t* bins, Double_t* backgr, Char_t* use, Int_t size, Int_t& p, Int_t dir, Double_t value) {
264  p += dir;
265  while ((p>=0) && (p<size) && use[p]) {
266  Double_t zn = bins[p]-backgr[p];
267  if (zn<value) return kTRUE;
268  p+=dir;
269  }
270  return kFALSE;
271 }
272 
274  TGo4FitData* data,
275  Int_t PolOrder,
276  Double_t AmplThreshold,
277  Double_t MinWidth,
278  Double_t MaxWidth) {
279  if ((fitter==0) || (data==0)) return;
280 
281  TGo4FitDataIter* iter = data->MakeIter();
282  if (iter==0) return;
283 
284  Int_t size = iter->CountPoints();
285 
286  if ((size<10) || (iter->ScalesSize()!=1)) { delete iter; return; }
287 
288  TArrayD Bins(size), Scales(size), Weights(size), Background(size);
289 
290  Weights.Reset(1.); Background.Reset(0.);
291 
292  Double_t* bins = Bins.GetArray();
293  Double_t* scales = Scales.GetArray();
294  Double_t* weights = Weights.GetArray();
295  Double_t* backgr = Background.GetArray();
296 
297  Int_t nbin = 0;
298  if (iter->Reset()) do {
299  bins[nbin] = data->GetAmplValue() * iter->Value();
300  if (iter->StandardDeviation()>0) weights[nbin] = 1./iter->StandardDeviation();
301  scales[nbin] = iter->x();
302  nbin++;
303  } while (iter->Next());
304 
305  delete iter;
306 
307  if ((PolOrder>=0) && (PolOrder<(size+1)*5)) {
308  Int_t nsegments = (PolOrder+1)*2;
309  TArrayC usage(size); usage.Reset(0);
310  for (Int_t nseg=0; nseg<nsegments; nseg++) {
311  Int_t lbound = Int_t((nseg+0.)/nsegments*size);
312  Int_t rbound = Int_t((nseg+1.)/nsegments*size-1.);
313  Int_t nselect = (rbound-lbound)/10;
314  if (nselect<1) nselect = 1;
315  Int_t nsel = 0;
316  while (nsel++<nselect) {
317  Int_t pmin=-1;
318  for(Int_t i=lbound;i<=rbound;i++)
319  if(!usage[i])
320  if ((pmin<0) || (bins[i]<bins[pmin])) pmin = i;
321  usage[pmin] = 1;
322  }
323  }
324 
325  TArrayD Coef(PolOrder+1);
326  DefinePolynom(size, bins, scales, Coef, weights, 0, usage.GetArray());
327 
328  for(Int_t i=0;i<size;i++)
329  backgr[i] = CalcPolynom(Coef, scales[i]);
330 
331  fitter->AddPolynomX(data->GetName(), "Pol", Coef);
332  }
333 
334  Double_t max=0.;
335  for(Int_t i=0;i<size;i++) {
336  Double_t zn = bins[i]-backgr[i];
337  if(zn>max) max = zn;
338  }
339  AmplThreshold = AmplThreshold*max;
340 
341  TArrayC usage(size); usage.Reset(1);
342  Char_t* use = usage.GetArray();
343 
344  Bool_t validline = kFALSE;
345 
346  do {
347 
348  validline = kFALSE;
349 
350  Int_t pmax=-1; Double_t max=0.;
351  for(Int_t i=0;i<size;i++)
352  if (use[i]) {
353  Double_t zn = bins[i]-backgr[i];
354  if((pmax<0) || (zn>max)) { pmax = i; max = zn; }
355  }
356 
357  Double_t sigma = TMath::Sqrt(1./weights[pmax]);
358 
359  if ((max<2.*sigma) || (max<AmplThreshold)) break;
360 
361  Int_t lbound = pmax, rbound = pmax;
362 
363  Bool_t lok = FindValue(bins, backgr, use, size, lbound, -1, max*0.5);
364  Bool_t rok = FindValue(bins, backgr, use, size, rbound, +1, max*0.5);
365 
366  if(lok || rok) {
367 
368  validline = kTRUE;
369 
370  Int_t w = (lok && rok) ? (rbound-lbound)/4 : ( lok ? (pmax-lbound)/2 : (rbound-pmax)/2 ) ;
371  if (w<2) w=2;
372 
373  TArrayD Coef1(3), Coef2(2);
374  DefinePolynomEx(size, bins, scales, weights, backgr, pmax-w, pmax+w, Coef1);
375  Double_t middle = -0.5*Coef1[1]/Coef1[2];
376  Double_t amplitude = CalcPolynom(Coef1, middle);
377 
378  Double_t left = 0, right = 0;
379  if (lok) {
380  DefinePolynomEx(size, bins, scales, weights, backgr, lbound-w, lbound+w, Coef2);
381  left = (amplitude*.5-Coef2[0])/Coef2[1];
382  }
383 
384  if (rok) {
385  DefinePolynomEx(size, bins, scales, weights, backgr, rbound-w, rbound+w, Coef2);
386  right = (amplitude*.5-Coef2[0])/Coef2[1];
387  }
388 
389  if (!lok) left = 2*middle-right; else
390  if (!rok) right = 2*middle-left;
391 
392  Double_t width = (right-left)*0.5;
393 
394  lbound = pmax-4*w; if (lbound<0) lbound=0;
395  rbound = pmax+4*w; if (rbound>=size) rbound=size-1;
396  for(Int_t i=lbound;i<=rbound;i++) use[i]=0;
397 
398  if(lok && rok && (width>=MinWidth) && (width<=MaxWidth) && (amplitude>AmplThreshold) && (amplitude>2*sigma))
399  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",0), middle, width, amplitude);
400  }
401 
402  } while(validline);
403 }
404 
405 #include "f_find_peaks.c"
406 
408  TGo4FitData* data,
409  Int_t MaxNumPeaks,
410  Int_t ChannelSum,
411  Double_t NoiseFactor,
412  Double_t NoiseMinimum,
413  Int_t MinimasOrder) {
414  if ((fitter==0) || (data==0) || (MaxNumPeaks<1)) return;
415 
416  TGo4FitDataIter* iter = data->MakeIter();
417  if (iter==0) return;
418 
419  Int_t size = iter->CountPoints();
420 
421  if ((size<10) || (iter->ScalesSize()!=1)) { delete iter; return; }
422 
423  TArrayD Bins(size), Scales(size);
424  Double_t* bins = Bins.GetArray();
425  Double_t* scales = Scales.GetArray();
426 
427  Int_t nbin = 0;
428  if (iter->Reset()) do {
429  bins[nbin] = data->GetAmplValue() * iter->Value();
430  scales[nbin] = iter->x();
431  nbin++;
432  } while (iter->Next());
433 
434  delete iter;
435 
436 
437  int NumPeaks = 0;
438  TArrayI Minimas(MaxNumPeaks), Maximas(MaxNumPeaks), MaximaWidths(MaxNumPeaks);
439  TArrayI PeaksLeft(MaxNumPeaks), PeaksRight(MaxNumPeaks), PeaksMaximum(MaxNumPeaks);
440  TArrayD MinimaValues(MaxNumPeaks), MinimaScales(MaxNumPeaks), MaximaIntegrals(MaxNumPeaks);
441 
442  go4fit_find_peaks(bins,
443  TYPE__DOUBLE,
444  0,
445  size,
446  ChannelSum,
447  NoiseFactor,
448  NoiseMinimum,
449  0,
450  MaxNumPeaks,
451  &NumPeaks,
452  Minimas.GetArray(),
453  MinimaValues.GetArray(),
454  Maximas.GetArray(),
455  MaximaWidths.GetArray(),
456  MaximaIntegrals.GetArray(),
457  PeaksLeft.GetArray(),
458  PeaksMaximum.GetArray(),
459  PeaksRight.GetArray(),
460  0,
461  0);
462 
463  if (NumPeaks<=0) return;
464 
465  Int_t NumMinimas = (NumPeaks<MaxNumPeaks) ? NumPeaks+1 : MaxNumPeaks;
466  for(int n=0; n<NumMinimas; n++) {
467  MinimaScales[n] = scales[Minimas[n]];
468  }
469 
470  if ((MinimasOrder>=0) && (NumMinimas>1)) {
471  Int_t order = MinimasOrder;
472  if (order>=NumMinimas) order = NumMinimas-1;
473  TArrayD Coef(order+1);
474  DefinePolynom(NumMinimas, MinimaValues.GetArray(), MinimaScales.GetArray(), Coef);
475  fitter->AddPolynomX(data->GetName(), "Pol", Coef);
476  }
477 
478  for(int n=0;n<NumPeaks;n++) {
479  Double_t pos = scales[PeaksMaximum[n]];
480  Double_t width = 0.;
481  if (PeaksRight[n]-PeaksLeft[n]<1)
482  width = scales[PeaksMaximum[n]] - scales[PeaksMaximum[n]-1];
483  else
484  width = TMath::Abs(scales[PeaksRight[n]] - scales[PeaksLeft[n]])/2.;
485  Double_t ampl = bins[PeaksMaximum[n]];
486 
487  fitter->AddGauss1(data->GetName(), fitter->FindNextName("Gauss",0), pos, width, ampl);
488  }
489 }
490 
491 
492 
Double_t Get1LineWidth()
static void DefinePolynom(Int_t size, Double_t *bins, Double_t *scales, TArrayD &Coef, Double_t *weight=0, Double_t *backgr=0, Char_t *use=0)
void SetupForSecond(Double_t LineWidth)
Bool_t HasWidths() const
Definition: TGo4FitData.h:499
Double_t x() const
Definition: TGo4FitData.h:484
void Set2NoiseFactor(Double_t factor)
virtual TGo4FitDataIter * MakeIter()
Definition: TGo4FitData.h:165
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:225
void SergeyLinevPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolOrder, Double_t AmplThreshold, Double_t MinWidth, Double_t MaxWidth)
void Set2NoiseMinimum(Double_t min)
Double_t Get2NoiseFactor()
void SetPeakFinderType(Int_t typ)
virtual Bool_t Reset(Bool_t UseRanges=kTRUE)
void DeleteModelsAssosiatedTo(const char *DataName)
Definition: TGo4Fitter.cxx:398
const Double_t * Widths() const
Definition: TGo4FitData.h:504
const char * GetDataName()
Double_t Get0MaxAmplFactor()
#define TYPE__DOUBLE
Definition: f_find_peaks.c:19
virtual void Print(Option_t *option) 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:347
void Set2ChannelSum(Int_t sum)
void Print(Option_t *option) 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)
Int_t CountPoints(Bool_t UseRanges=kTRUE)
Double_t Get2NoiseMinimum()
virtual Bool_t Next(Bool_t UseRanges=kTRUE)
TString FindNextName(const char *Head, Int_t start, Bool_t isModel=kTRUE)
virtual void DoAction(TGo4FitterAbstract *Fitter)
void SetupPolynomialBackground(Int_t PolynomOrder)
void SetUsePolynom(Bool_t use)
void SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum)
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 usage(const char *subtopic=0)
void Set1LineWidth(Double_t width)
Double_t Value() const
Definition: TGo4FitData.h:514
static Double_t CalcPolynom(const TArrayD &Coef, Double_t x)
void ROOTPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolynomOrder, Double_t Sigma)
TGo4FitData * FindData(const char *DataName)
Definition: TGo4Fitter.cxx:113
void SetPolynomOrder(Int_t order)
Double_t StandardDeviation() const
Definition: TGo4FitData.h:519
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)
Int_t ScalesSize() const
Definition: TGo4FitData.h:474