22 #include "TSpectrum.h"
33 fiPeakFinderType(0), fxDataName(DataName), fbClearModels(ClearModels),
34 fbUsePolynom(PolOrder>=0), fiPolynomOrder(PolOrder>=0 ? PolOrder : 0),
35 fd0MinWidth(5.), fd0MaxWidth(30.), fd0MaxAmplFactor(0.2),
37 fd2NoiseFactor(2.), fd2NoiseMinimum(5.), fi2ChannelSum(2) {
72 if (fitter==0)
return;
149 if ((fitter==0) || (data==0))
return;
156 if ((size<10) || (iter->
ScalesSize()!=1)) {
delete iter;
return; }
158 TArrayD Bins(size), Scales(size), HScales(size+1);
161 if (iter->
Reset())
do {
163 Double_t pos = iter->
x();
166 HScales[nbin] = pos - width/2.;
167 HScales[nbin+1] = pos + width/2;
169 }
while (iter->
Next());
174 fitter->
AddPolynomX(data->GetName(),
"Pol", PolynomOrder);
176 TH1D histo(
"ROOTPeakFinder",
"ROOTPeakFinder spectrum", size, HScales.GetArray());
177 histo.SetDirectory(0);
178 for (Int_t i=0;i<size;i++)
179 histo.SetBinContent(i+1, Bins[i]);
182 sp.Search(&histo, Sigma);
185 for(Int_t n=0;n<sp.GetNPeaks();n++) {
186 Double_t pos = sp.GetPositionX()[n];
187 Double_t width = Sigma;
207 Double_t power = 1., sum = 0.;
208 for(Int_t n=0;n<Coef.GetSize();n++)
209 { sum+=Coef[n]*power; power*=x; }
214 Double_t* weight, Double_t* backgr, Char_t* use) {
215 Int_t maxorder = Coef.GetSize()-1;
217 TMatrixD matr(maxorder+1,maxorder+1);
218 TVectorD vector(maxorder+1);
220 for(Int_t order1=0;order1<=maxorder;order1++)
221 for(Int_t order2=order1;order2<=maxorder;order2++) {
223 for(Int_t i=0;i<size;i++)
224 if((use==0) || use[i]) {
225 Double_t zn = TMath::Power(scales[i],order1+order2);
226 if (weight) zn*=weight[i];
229 matr(order1,order2) = sum;
230 matr(order2,order1) = sum;
233 for(Int_t order1=0;order1<=maxorder;order1++) {
235 for(Int_t i=0;i<size;i++)
236 if((use==0) || use[i]) {
237 Double_t zn = bins[i];
238 if (backgr) zn-=backgr[i];
239 zn*=TMath::Power(scales[i],order1);
240 if (weight) zn*=weight[i];
243 vector(order1) = sum;
249 for(Int_t order=0;order<=maxorder;order++)
250 Coef[order]=vector(order);
254 Int_t lbound, Int_t rbound, TArrayD& Coef) {
255 TArrayC use(size); use.Reset(0);
256 if (lbound<0) lbound = 0;
257 if (rbound>=size) rbound=size-1;
258 for(Int_t i=lbound;i<=rbound;i++) use[i]=1;
259 DefinePolynom(size, bins, scales, Coef, weight, backgr, use.GetArray());
262 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 while ((p>=0) && (p<size) && use[p]) {
265 Double_t zn = bins[p]-backgr[p];
266 if (zn<value)
return kTRUE;
275 Double_t AmplThreshold,
278 if ((fitter==0) || (data==0))
return;
285 if ((size<10) || (iter->
ScalesSize()!=1)) {
delete iter;
return; }
287 TArrayD Bins(size), Scales(size), Weights(size), Background(size);
289 Weights.Reset(1.); Background.Reset(0.);
291 Double_t* bins = Bins.GetArray();
292 Double_t* scales = Scales.GetArray();
293 Double_t* weights = Weights.GetArray();
294 Double_t* backgr = Background.GetArray();
297 if (iter->
Reset())
do {
300 scales[nbin] = iter->
x();
302 }
while (iter->
Next());
306 if ((PolOrder>=0) && (PolOrder<(size+1)*5)) {
307 Int_t nsegments = (PolOrder+1)*2;
308 TArrayC
usage(size); usage.Reset(0);
309 for (Int_t nseg=0; nseg<nsegments; nseg++) {
310 Int_t lbound = Int_t((nseg+0.)/nsegments*size);
311 Int_t rbound = Int_t((nseg+1.)/nsegments*size-1.);
312 Int_t nselect = (rbound-lbound)/10;
313 if (nselect<1) nselect = 1;
315 while (nsel++<nselect) {
317 for(Int_t i=lbound;i<=rbound;i++)
319 if ((pmin<0) || (bins[i]<bins[pmin])) pmin = i;
324 TArrayD Coef(PolOrder+1);
325 DefinePolynom(size, bins, scales, Coef, weights, 0, usage.GetArray());
327 for(Int_t i=0;i<size;i++)
334 for(Int_t i = 0; i < size; i++) {
335 Double_t zn = bins[i]-backgr[i];
336 if(zn > dmax) dmax = zn;
338 AmplThreshold = AmplThreshold*dmax;
340 TArrayC
usage(size); usage.Reset(1);
341 Char_t* use = usage.GetArray();
343 Bool_t validline = kFALSE;
351 for(Int_t i=0;i<size;i++)
353 Double_t zn = bins[i]-backgr[i];
354 if((pmax<0) || (zn>max)) { pmax = i; max = zn; }
357 Double_t sigma = TMath::Sqrt(1./weights[pmax]);
359 if ((max<2.*sigma) || (max<AmplThreshold))
break;
361 Int_t lbound = pmax, rbound = pmax;
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);
370 Int_t w = (lok && rok) ? (rbound-lbound)/4 : ( lok ? (pmax-lbound)/2 : (rbound-pmax)/2 ) ;
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];
378 Double_t left = 0, right = 0;
380 DefinePolynomEx(size, bins, scales, weights, backgr, lbound-w, lbound+w, Coef2);
381 left = (amplitude*.5-Coef2[0])/Coef2[1];
385 DefinePolynomEx(size, bins, scales, weights, backgr, rbound-w, rbound+w, Coef2);
386 right = (amplitude*.5-Coef2[0])/Coef2[1];
389 if (!lok) left = 2*middle-right;
else
390 if (!rok) right = 2*middle-left;
392 Double_t width = (right-left)*0.5;
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;
398 if(lok && rok && (width>=MinWidth) && (width<=MaxWidth) && (amplitude>AmplThreshold) && (amplitude>2*sigma))
411 Double_t NoiseFactor,
412 Double_t NoiseMinimum,
413 Int_t MinimasOrder) {
414 if ((fitter==0) || (data==0) || (MaxNumPeaks<1))
return;
421 if ((size<10) || (iter->
ScalesSize()!=1)) {
delete iter;
return; }
423 TArrayD Bins(size), Scales(size);
424 Double_t* bins = Bins.GetArray();
425 Double_t* scales = Scales.GetArray();
428 if (iter->
Reset())
do {
430 scales[nbin] = iter->
x();
432 }
while (iter->
Next());
438 TArrayI Minimas(MaxNumPeaks), Maximas(MaxNumPeaks), MaximaWidths(MaxNumPeaks);
439 TArrayI PeaksLeft(MaxNumPeaks), PeaksRight(MaxNumPeaks), PeaksMaximum(MaxNumPeaks);
440 TArrayD MinimaValues(MaxNumPeaks), MinimaScales(MaxNumPeaks), MaximaIntegrals(MaxNumPeaks);
453 MinimaValues.GetArray(),
455 MaximaWidths.GetArray(),
456 MaximaIntegrals.GetArray(),
457 PeaksLeft.GetArray(),
458 PeaksMaximum.GetArray(),
459 PeaksRight.GetArray(),
463 if (NumPeaks<=0)
return;
465 Int_t NumMinimas = (NumPeaks<MaxNumPeaks) ? NumPeaks+1 : MaxNumPeaks;
466 for(
int n=0; n<NumMinimas; n++) {
467 MinimaScales[n] = scales[Minimas[n]];
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);
478 for(
int n=0;n<NumPeaks;n++) {
479 Double_t pos = scales[PeaksMaximum[n]];
481 if (PeaksRight[n]-PeaksLeft[n]<1)
482 width = scales[PeaksMaximum[n]] - scales[PeaksMaximum[n]-1];
484 width = TMath::Abs(scales[PeaksRight[n]] - scales[PeaksLeft[n]])/2.;
485 Double_t ampl = bins[PeaksMaximum[n]];
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)
Int_t GetPeakFinderType() const
void Set2NoiseFactor(Double_t factor)
virtual TGo4FitDataIter * MakeIter()
Double_t Get0MinWidth() 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.)
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)
void Set2NoiseMinimum(Double_t min)
void SetPeakFinderType(Int_t typ)
Bool_t GetUsePolynom() const
virtual Bool_t Reset(Bool_t UseRanges=kTRUE)
Double_t Get2NoiseMinimum() const
void DeleteModelsAssosiatedTo(const char *DataName)
const Double_t * Widths() const
Int_t Get2ChannelSum() const
const char * GetDataName()
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)
void Set2ChannelSum(Int_t sum)
void Print(Option_t *option) const
Double_t Get2NoiseFactor() const
virtual ~TGo4FitPeakFinder()
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)
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)
Double_t Get1LineWidth() const
void SetUsePolynom(Bool_t use)
void SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum)
void usage(const char *subtopic=0)
Double_t Get0MaxWidth() const
void Set1LineWidth(Double_t width)
static Double_t CalcPolynom(const TArrayD &Coef, Double_t x)
Double_t Get0MaxAmplFactor() const
void ROOTPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolynomOrder, Double_t Sigma)
TGo4FitData * FindData(const char *DataName)
void SetPolynomOrder(Int_t order)
Double_t StandardDeviation() const
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 GetPolynomOrder() const