86 for (i=0;i<
fNin;i++) {
94 for (i=0;i<
fNin;i++) {
123 if (opt.
Contains(
"normal")) kernel = 2;
134 index =
new Int_t[nout];
147 if (index) {
delete [] index; index = 0;}
173 while ((imin < n) && (x[imin] < xp[0] - cutoff))
176 for (
Int_t j=0;j<np;j++) {
181 for (
Int_t i=imin;i<
n;i++) {
182 if (x[i] < x0 - cutoff) imin = i;
183 if (x[i] > x0 + cutoff)
break;
185 if (kernel == 1) w = 1;
254 Int_t i, iiter, j, last, m1,
m2, nleft, nright,
ns;
255 Double_t alpha,
c1, c9, cmad, cut, d1, d2, denom,
r;
276 while (iiter <= iter+1) {
285 d1 = x[i] - x[nleft];
286 d2 = x[nright+1] - x[i];
299 Lowest(&x[1], &y[1], n, x[i], ys[i], nleft, nright,
300 res, iterg1, rw, ok);
301 if (!ok) ys[i] = y[i];
305 denom = x[i]-x[last];
308 for(j = last+1; j < i; j++) {
309 alpha = (x[j]-x[last])/denom;
310 ys[j] = alpha*ys[i] + (1.-alpha)*ys[last];
318 cut = x[last] + delta;
319 for (i = last+1; i <=
n; i++) {
322 if (x[i] == x[last]) {
334 res[i] = y[i+1] - ys[i+1];
349 cmad = 3.*(rw[m1]+rw[
m2]);
356 for(i=0 ; i<
n ; i++) {
361 rw[i] = (1.-(r/cmad)*(r/cmad))*(1.-(r/cmad)*(r/cmad));
402 d = (r/
h)*(r/h)*(r/
h);
403 w[j] = (1.-
d)*(1.- d)*(1.-
d);
408 }
else if (x[j] > xs)
420 for(j=nleft ; j<=nrt ; j++)
425 for(j=nleft ; j<=nrt ; j++)
429 for(j=nleft ; j<=nrt ; j++)
430 c += w[j]*(x[j]-a)*(x[j]-
a);
434 for(j=nleft; j <= nrt; j++)
435 w[j] *= (b*(x[j]-a) + 1.);
439 for(j=nleft; j <= nrt; j++)
485 if (span < 0 || span > 1) {
486 std::cout <<
"Error: Span must be between 0 and 1" << std::endl;
495 if (fMinX < 0 || fMaxX > 1) {
496 std::cout <<
"Error: x must be between 0 and 1 for periodic smooth" << std::endl;
505 for (i=0; i<
fNout; i++) {
511 for (i=0; i<
fNin; i++) {
512 if (w == 0) weight[i] = 1;
513 else weight[i] = w[i];
517 Int_t nTmp = (fNin+1)*8;
519 for (i=0; i<nTmp; i++) {
584 Double_t spans[3] = { 0.05, 0.2, 0.5 };
607 if (sw > 0.0) a = sy / sw;
608 for (j=1;j<=
n;++j) smo[j] = a;
615 while (scale <= 0.0) {
625 if (iper == 2 && (x[1] < 0.0 || x[n] > 1.0)) {
628 if (jper < 1 || jper > 2) {
632 BDRsmooth(n, &x[1], &y[1], &w[1], span, jper, vsmlsq,
633 &smo[1], &sc[sc_offset]);
638 for (i = 1; i <= 3; ++i) {
639 BDRsmooth(n, &x[1], &y[1], &w[1], spans[i - 1], jper, vsmlsq,
640 &sc[((i<<1)-1)*n + 1], &sc[n*7 + 1]);
641 BDRsmooth(n, &x[1], &sc[n*7 + 1], &w[1], spans[1], -jper, vsmlsq,
642 &sc[(i<<1)*n + 1], &h[1]);
645 for (j=1; j<=
n; ++j) {
647 for (i=1; i<=3; ++i) {
648 if (sc[j + (i<<1)*n] < resmin) {
649 resmin = sc[j + (i<<1)*n];
650 sc[j + n*7] = spans[i-1];
654 if (alpha>0.0 && alpha<=10.0 && resmin<sc[j + n*6] && resmin>0.0) {
658 sc[j + n*7] += (spans[2] - sc[j + n*7]) *
TMath::Power(d1, d2);
662 BDRsmooth(n, &x[1], &sc[n*7 + 1], &w[1], spans[1], -jper, vsmlsq,
663 &sc[(n<<1) + 1], &h[1]);
665 for (j=1; j<=
n; ++j) {
666 if (sc[j + (n<<1)] <= spans[0]) {
667 sc[j + (n<<1)] = spans[0];
669 if (sc[j + (n<<1)] >= spans[2]) {
670 sc[j + (n<<1)] = spans[2];
672 f = sc[j + (n<<1)] - spans[1];
674 f = -f / (spans[1] - spans[0]);
675 sc[j + (n<<2)] = (1.0 - f) * sc[j + n*3] + f * sc[j +
n];
677 f /= spans[2] - spans[1];
678 sc[j + (n<<2)] = (1.0 - f) * sc[j + n*3] + f * sc[j + n*5];
682 BDRsmooth(n, &x[1], &sc[(n<<2) + 1], &w[1], spans[0], -jper, vsmlsq,
698 Int_t i, j, j0, in, out, it, jper, ibw;
718 ibw = (
Int_t)(span * 0.5 * n + 0.5);
724 for (i=1; i<=it; ++i) {
738 xm = (fbo * xm + wt * xti) / fbw;
739 ym = (fbo * ym + wt * y[j]) / fbw;
743 tmp = fbw * wt * (xti - xm) / fbo;
745 var += tmp * (xti - xm);
746 cvar += tmp * (y[j] - ym);
749 for (j=1; j<=
n; ++j) {
752 if (!(jper != 2 && (out < 1 || in > n))) {
771 tmp = fbo * wt * (xto - xm) / fbw;
773 var -= tmp * (xto - xm);
774 cvar -= tmp * (y[out] - ym);
776 xm = (fbo * xm - wt * xto) / fbw;
777 ym = (fbo * ym - wt * y[out]) / fbw;
783 xm = (fbo * xm + wt * xti) / fbw;
784 ym = (fbo * ym + wt * y[in]) / fbw;
788 tmp = fbw * wt * (xti - xm) / fbo;
790 var += tmp * (xti - xm);
791 cvar += tmp * (y[in] - ym);
798 smo[j] = a * (x[j] - xm) + ym;
832 if (x[j + 1] > x[j]) {
846 for (i=j0; i<=j; ++i) {
869 for (i=0;i<
fNin;i++) {
886 for (i=1;i<fNin+1;i++) {
888 vMin = vMean = vMax = yin[index[i-1]];
889 while ((i < fNin) && (rank[index[i]] == rank[index[i-1]])) {
890 vMean += yin[index[i]];
891 vMax = (vMax < yin[index[i]]) ? yin[index[i]] : vMax;
892 vMin = (vMin > yin[index[i]]) ? yin[index[i]] : vMin;
898 x[k] = xin[index[i-1]];
899 if (ndup == 1) {y[k++] = yin[index[i-1]];}
919 for (i=0;i<
fNin;i++) {
1004 if (opt.
Contains(
"linear")) iKind = 1;
1005 else if (opt.
Contains(
"constant")) iKind = 2;
1007 if (f < 0 || f > 1) {
1008 std::cout <<
"Error: Invalid f value" << std::endl;
1024 std::cout <<
"Error: Method not known: " << ties << std::endl;
1031 Approxin(grin, iKind, ylow, yhigh, rule, iTies);
1045 if (xout == 0) x =
fMinX + i*delta;
1067 if(v < x[i])
return ylow;
1068 if(v > x[j])
return yhigh;
1072 Int_t ij = (i + j)/2;
1073 if(v < x[ij]) j = ij;
1078 if(v == x[j])
return y[j];
1079 if(v == x[i])
return y[i];
1082 return y[i] + (y[j] - y[i]) * ((v - x[i])/(x[j] - x[i]));
1084 return y[i] * (1-
f) + y[j] * f;
1096 if (x < y)
return -1;
1097 if (x > y)
return 1;
1111 for (pL = 0, pR = n - 1; pL < pR; ) {
1113 for(i = pL, j = pR; i <= j;) {
1116 if (i <= j) { w = x[i]; x[i++] = x[j]; x[j--] = w; }
1138 for (
Int_t i=0;i<
n;i++) {
1139 if ((i > 0) && (a[index[i]] == a[index[i-1]])) {
1140 rank[index[i]] = i-1;
1143 rank[index[i]] = i-k;
static void BDRsupsmu(Int_t n, Double_t *x, Double_t *y, Double_t *w, Int_t iper, Double_t span, Double_t alpha, Double_t *smo, Double_t *sc)
Friedmanns super smoother (Friedman, 1984).
TGraph * SmoothKern(TGraph *grin, Option_t *option="normal", Double_t bandwidth=0.5, Int_t nout=100, Double_t *xout=0)
Smooth data with Kernel smoother.
virtual ~TGraphSmooth()
GraphSmooth destructor.
void Approxin(TGraph *grin, Int_t iKind, Double_t &Ylow, Double_t &Yhigh, Int_t rule, Int_t iTies)
Sort data points and eliminate double x values.
Short_t Min(Short_t a, Short_t b)
void ToLower()
Change string to lower-case.
A helper class to smooth TGraph.
TGraph * SmoothSuper(TGraph *grin, Option_t *option="", Double_t bass=0, Double_t span=0, Bool_t isPeriodic=kFALSE, Double_t *w=0)
Smooth data with Super smoother.
static void Lowest(Double_t *x, Double_t *y, Int_t n, Double_t &xs, Double_t &ys, Int_t nleft, Int_t nright, Double_t *w, Bool_t userw, Double_t *rw, Bool_t &ok)
Fit value at x[i] Based on R function lowest: Translated to C++ by C.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
TGraph * SmoothLowess(TGraph *grin, Option_t *option="", Double_t span=0.67, Int_t iter=3, Double_t delta=0)
Smooth data with Lowess smoother.
The TNamed class is the base class for all named ROOT classes.
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
void Lowess(Double_t *x, Double_t *y, Int_t n, Double_t *ys, Double_t span, Int_t iter, Double_t delta)
Lowess regression smoother.
static Int_t Rcmp(Double_t x, Double_t y)
Static function if (ISNAN(x)) return 1; if (ISNAN(y)) return -1;.
static void BDRsmooth(Int_t n, Double_t *x, Double_t *y, Double_t *w, Double_t span, Int_t iper, Double_t vsmlsq, Double_t *smo, Double_t *acvr)
Function for super smoother Based on R function supsmu: Translated to C++ by C.
static constexpr double m2
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
static void Rank(Int_t n, Double_t *a, Int_t *index, Int_t *rank, Bool_t down=kTRUE)
static function
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
TGraph * Approx(TGraph *grin, Option_t *option="linear", Int_t nout=50, Double_t *xout=0, Double_t yleft=0, Double_t yright=0, Int_t rule=0, Double_t f=0, Option_t *ties="mean")
Approximate data points.
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
void Smoothin(TGraph *grin)
Sort input data points.
static void Psort(Double_t *x, Int_t n, Int_t k)
Static function based on R function rPsort: adapted to C++ by Christian Stratowa (R source file: R_so...
Short_t Max(Short_t a, Short_t b)
A Graph is a graphics object made of two arrays X and Y with npoints each.
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
A TGraphErrors is a TGraph with error bars.
Double_t Sqrt(Double_t x)
virtual void Set(Int_t n)
Set number of points in the graph Existing coordinates are preserved New coordinates above fNpoints a...
static Double_t Approx1(Double_t v, Double_t f, Double_t *x, Double_t *y, Int_t n, Int_t iKind, Double_t Ylow, Double_t Yhigh)
Approximate one data point.
static constexpr double ns
static void BDRksmooth(Double_t *x, Double_t *y, Int_t n, Double_t *xp, Double_t *yp, Int_t np, Int_t kernel, Double_t bw)
Smooth data with specified kernel.