103 if (!fitter || !data)
112 if ((size < 10) || (iter->ScalesSize() != 1))
115 TArrayD Bins(size), Scales(size), HScales(size + 1);
121 Double_t pos = iter->x();
122 Double_t width = iter->HasWidths() ? iter->Widths()[0] : 1.;
124 HScales[nbin] = pos - width / 2.;
125 HScales[nbin + 1] = pos + width / 2;
127 }
while (iter->Next());
129 if (PolynomOrder >= 0)
130 fitter->
AddPolynomX(data->GetName(),
"Pol", PolynomOrder);
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]);
138 sp.Search(&histo, Sigma);
141 for (Int_t n = 0; n < sp.GetNPeaks(); n++) {
142 Double_t pos = sp.GetPositionX()[n];
143 Double_t width = Sigma;
162 Double_t power = 1., sum = 0.;
163 for (Int_t n = 0; n < Coef.GetSize(); n++) {
164 sum += Coef[n] * power;
171 Double_t *backgr, Char_t *use)
173 Int_t maxorder = Coef.GetSize() - 1;
175 TMatrixD matr(maxorder + 1, maxorder + 1);
176 TVectorD vector(maxorder + 1);
178 for (Int_t order1 = 0; order1 <= maxorder; order1++)
179 for (Int_t order2 = order1; order2 <= maxorder; order2++) {
181 for (Int_t i = 0; i < size; i++)
182 if (!use || use[i]) {
183 Double_t zn = TMath::Power(scales[i], order1 + order2);
188 matr(order1, order2) = sum;
189 matr(order2, order1) = sum;
192 for (Int_t order1 = 0; order1 <= maxorder; order1++) {
194 for (Int_t i = 0; i < size; i++)
195 if (!use || use[i]) {
196 Double_t zn = bins[i];
199 zn *= TMath::Power(scales[i], order1);
204 vector(order1) = sum;
210 for (Int_t order = 0; order <= maxorder; order++)
211 Coef[order] = vector(order);
215 Double_t *backgr, Int_t lbound, Int_t rbound, TArrayD &Coef)
223 for (Int_t i = lbound; i <= rbound; i++)
225 DefinePolynom(size, bins, scales, Coef, weight, backgr, use.GetArray());
228Bool_t
FindValue(Double_t *bins, Double_t *backgr, Char_t *use, Int_t size, Int_t &p, Int_t dir, Double_t value)
231 while ((p >= 0) && (p < size) && use[p]) {
232 Double_t zn = bins[p] - backgr[p];
241 Double_t AmplThreshold, Double_t MinWidth, Double_t MaxWidth)
243 if (!fitter || !data)
252 if ((size < 10) || (iter->ScalesSize() != 1))
255 TArrayD Bins(size), Scales(size), Weights(size), Background(size);
258 Background.Reset(0.);
260 Double_t *bins = Bins.GetArray();
261 Double_t *scales = Scales.GetArray();
262 Double_t *weights = Weights.GetArray();
263 Double_t *backgr = Background.GetArray();
269 if (iter->StandardDeviation() > 0)
270 weights[nbin] = 1. / iter->StandardDeviation();
271 scales[nbin] = iter->x();
273 }
while (iter->Next());
275 if ((PolOrder >= 0) && (PolOrder < (size + 1) * 5)) {
276 Int_t nsegments = (PolOrder + 1) * 2;
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;
286 while (nsel++ < nselect) {
288 for (Int_t i = lbound; i <= rbound; i++)
290 if ((pmin < 0) || (bins[i] < bins[pmin]))
296 TArrayD Coef(PolOrder + 1);
299 for (Int_t i = 0; i < size; i++)
306 for (Int_t i = 0; i < size; i++) {
307 Double_t zn = bins[i] - backgr[i];
311 AmplThreshold = AmplThreshold * dmax;
315 Char_t *use =
usage.GetArray();
317 Bool_t validline = kFALSE;
325 for (Int_t i = 0; i < size; i++)
327 Double_t zn = bins[i] - backgr[i];
328 if ((pmax < 0) || (zn > max)) {
334 Double_t sigma = TMath::Sqrt(1. / weights[pmax]);
336 if ((max < 2. * sigma) || (max < AmplThreshold))
339 Int_t lbound = pmax, rbound = pmax;
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);
348 Int_t w = (lok && rok) ? (rbound - lbound) / 4 : (lok ? (pmax - lbound) / 2 : (rbound - pmax) / 2);
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];
357 Double_t left = 0, right = 0;
359 DefinePolynomEx(size, bins, scales, weights, backgr, lbound - w, lbound + w, Coef2);
360 left = (amplitude * .5 - Coef2[0]) / Coef2[1];
364 DefinePolynomEx(size, bins, scales, weights, backgr, rbound - w, rbound + w, Coef2);
365 right = (amplitude * .5 - Coef2[0]) / Coef2[1];
369 left = 2 * middle - right;
371 right = 2 * middle - left;
373 Double_t width = (right - left) * 0.5;
375 lbound = pmax - 4 * w;
378 rbound = pmax + 4 * w;
381 for (Int_t i = lbound; i <= rbound; i++)
384 if (lok && rok && (width >= MinWidth) && (width <= MaxWidth) && (amplitude > AmplThreshold) &&
385 (amplitude > 2 * sigma))
395 Double_t NoiseFactor, Double_t NoiseMinimum, Int_t MinimasOrder)
397 if (!fitter || !data || (MaxNumPeaks < 1))
406 if ((size < 10) || (iter->ScalesSize() != 1))
409 TArrayD Bins(size), Scales(size);
410 Double_t *bins = Bins.GetArray();
411 Double_t *scales = Scales.GetArray();
417 scales[nbin] = iter->x();
419 }
while (iter->Next());
422 TArrayI Minimas(MaxNumPeaks), Maximas(MaxNumPeaks), MaximaWidths(MaxNumPeaks);
423 TArrayI PeaksLeft(MaxNumPeaks), PeaksRight(MaxNumPeaks), PeaksMaximum(MaxNumPeaks);
424 TArrayD MinimaValues(MaxNumPeaks), MinimaScales(MaxNumPeaks), MaximaIntegrals(MaxNumPeaks);
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(),
434 Int_t NumMinimas = (NumPeaks < MaxNumPeaks) ? NumPeaks + 1 : MaxNumPeaks;
435 for (
int n = 0; n < NumMinimas; n++) {
436 MinimaScales[n] = scales[Minimas[n]];
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);
448 for (
int n = 0; n < NumPeaks; n++) {
449 Double_t pos = scales[PeaksMaximum[n]];
451 if (PeaksRight[n] - PeaksLeft[n] < 1)
452 width = scales[PeaksMaximum[n]] - scales[PeaksMaximum[n] - 1];
454 width = TMath::Abs(scales[PeaksRight[n]] - scales[PeaksLeft[n]]) / 2.;
455 Double_t ampl = bins[PeaksMaximum[n]];
void usage(const char *subtopic=nullptr)
Bool_t FindValue(Double_t *bins, Double_t *backgr, Char_t *use, Int_t size, Int_t &p, Int_t dir, Double_t value)
Double_t GetAmplValue()
Return value of amplitude parameter.
Int_t CountPoints(Bool_t UseRanges=kTRUE)
Counts total number of points in data object.
Basic abstract class for representing data, which should be fitted.
virtual std::unique_ptr< TGo4FitDataIter > MakeIter()
Creates iterator for data object.
void Print(Option_t *option="") const override
void Set0MinWidth(Double_t min)
Double_t Get0MaxAmplFactor() const
Double_t Get1LineWidth() const
virtual ~TGo4FitPeakFinder()
Destroys TGo4FitPeakFinder object.
void SetupForThird(Double_t NoiseFactor, Double_t NoiseMinimum, Int_t ChannelSum)
const char * GetDataName() const
void Print(Option_t *option="") const override
Print information on standard output.
void ROOTPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolynomOrder, Double_t Sigma)
void DoAction(TGo4FitterAbstract *Fitter) override
Bool_t GetUsePolynom() 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)
Int_t GetPeakFinderType() const
void SetPolynomOrder(Int_t order)
void SetupForSecond(Double_t LineWidth)
Bool_t GetClearModels() const
Double_t Get2NoiseFactor() const
void SergeyLinevPeakFinder(TGo4Fitter *fitter, TGo4FitData *data, Int_t PolOrder, Double_t AmplThreshold, Double_t MinWidth, Double_t MaxWidth)
Perform simple peak finder algorithm.
void Set0MaxWidth(Double_t max)
Double_t Get0MaxWidth() const
void Set0MaxAmplFactor(Double_t factor)
static Double_t CalcPolynom(const TArrayD &Coef, Double_t x)
void SetPeakFinderType(Int_t typ)
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)
Hans Essel (c) peak finder.
Int_t Get2ChannelSum() const
void Set2NoiseFactor(Double_t factor)
Double_t Get2NoiseMinimum() const
TGo4FitPeakFinder()
Default constructor.
void Set1LineWidth(Double_t width)
void Set2ChannelSum(Int_t sum)
void SetUsePolynom(Bool_t use)
void SetupPolynomialBackground(Int_t PolynomOrder)
Double_t Get0MinWidth() const
void Set2NoiseMinimum(Double_t min)
Double_t fd0MaxAmplFactor
Int_t GetPolynomOrder() const
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)
void SetupForFirst(Double_t MaxAmplFactor, Double_t MinWidth, Double_t MaxWidth)
TGo4FitterAction()
Default constructor.
Central class of Go4Fit package.
void DeleteModelsAssosiatedTo(const char *DataName)
Remove models associated with specific data.
TGo4FitModelGauss1 * AddGauss1(const char *DataName, const char *ModelName, Double_t iPosition, Double_t iWidth, Double_t iAmpl=1., Int_t Axis=0)
Add 1-dim gaussian model to fitter.
void AddPolynomX(const char *DataName, const char *NamePrefix, Int_t MaxOrder=1, Int_t GroupIndex=0, Double_t lrange=0., Double_t rrange=0.)
Construct 1-dim polynom for specified data object for x scale.
TString FindNextName(const char *Head, Int_t start, Bool_t isModel=kTRUE)
TGo4FitData * FindData(const char *DataName)
Return data object with given name.
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)