library: libGraf #include "TGraphSmooth.h" |
TGraphSmooth
class description - source file - inheritance tree (.pdf)
public:
TGraphSmooth()
TGraphSmooth(const char* name)
TGraphSmooth(const TGraphSmooth&)
virtual ~TGraphSmooth()
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")
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)
void Approxin(TGraph* grin, Int_t iKind, Double_t& Ylow, Double_t& Yhigh, Int_t rule, Int_t iTies)
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)
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)
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)
static TClass* Class()
virtual TClass* IsA() const
void Lowess(Double_t* x, Double_t* y, Int_t n, Double_t* ys, Double_t span, Int_t iter, Double_t delta)
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)
TGraphSmooth& operator=(const TGraphSmooth&)
static void Psort(Double_t* x, Int_t n, Int_t k)
static void Rank(Int_t n, Double_t* a, Int_t* index, Int_t* rank, Bool_t down = kTRUE)
static Int_t Rcmp(Double_t x, Double_t y)
virtual void ShowMembers(TMemberInspector& insp, char* parent)
void Smoothin(TGraph* grin)
TGraph* SmoothKern(TGraph* grin, Option_t* option = "normal", Double_t bandwidth = 0.5, Int_t nout = 100, Double_t* xout = 0)
TGraph* SmoothLowess(TGraph* grin, Option_t* option = "", Double_t span = 0.67, Int_t iter = 3, Double_t delta = 0)
TGraph* SmoothSuper(TGraph* grin, Option_t* option = "", Double_t bass = 0, Double_t span = 0, Bool_t isPeriodic = kFALSE, Double_t* w = 0)
virtual void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
protected:
Int_t fNin Number of input points
Int_t fNout Number of output points
TGraph* fGin Input graph
TGraph* fGout Output graph
Double_t fMinX Minimum value of array X
Double_t fMaxX Maximum value of array X
TGraphSmooth
A helper class to smooth TGraph
see examples in $ROOTSYS/tutorials/motorcycle.C and approx.C
______________________________________________________________________
TGraphSmooth(): TNamed()
*-*-*-*-*-*-*-*-*Default GraphSmooth constructor *-*-*-*-*-*-*-*-*-*-*
===============================
TGraphSmooth(const char *name): TNamed(name,"")
*-*-*-*-*-*-*-*-*GraphSmooth constructor *-*-*-*-*-*-*-*-*-*-*-*-*-*-*
=======================
~TGraphSmooth()
*-*-*-*-*-*-*-*-*GraphSmooth destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
======================
void Smoothin(TGraph *grin)
*-*-*-*-*-*-*-*-*Sort input data points*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
======================
TGraph* SmoothKern(TGraph *grin, Option_t *option,
Double_t bandwidth, Int_t nout, Double_t *xout)
*-*-*-*-*-*-*-*-*Smooth data with Kernel smoother*-*-*-*-*-*-*-*-*-*-*
================================
Smooth grin with the Nadaraya-Watson kernel regression estimate.
Arguments:
grin: input graph
option: the kernel to be used: "box", "normal"
bandwidth: the bandwidth. The kernels are scaled so that their quartiles
(viewed as probability densities) are at +/- 0.25*bandwidth.
nout: If xout is not specified, interpolation takes place at equally
spaced points spanning the interval [min(x), max(x)], where
nout = max(nout, number of input data).
xout: an optional set of values at which to evaluate the fit
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*-*-*-*-*-*-*-*-*-*-*
*-* =================================
Based on R function ksmooth: Translated to C++ by C. Stratowa
(R source file: ksmooth.c by B.D.Ripley Copyright (C) 1998)
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
TGraph* SmoothLowess(TGraph *grin, Option_t *option ,
Double_t span, Int_t iter, Double_t delta)
*-*-*-*-*-*-*-*-*Smooth data with Lowess smoother*-*-*-*-*-*-*-*-*-*-*
================================
This function performs the computations for the LOWESS smoother
(see the reference below). Lowess returns the output points
x and y which give the coordinates of the smooth.
Arguments:
grin: Input graph
span: the smoother span. This gives the proportion of points in the plot
which influence the smooth at each value.
Larger values give more smoothness.
iter: the number of robustifying iterations which should be performed.
Using smaller values of iter will make lowess run faster.
delta: values of x which lie within delta of each other replaced by a
single value in the output from lowess.
For delta = 0, delta will be calculated.
References:
Cleveland, W. S. (1979) Robust locally weighted regression and smoothing
scatterplots. J. Amer. Statist. Assoc. 74, 829-836.
Cleveland, W. S. (1981) LOWESS: A program for smoothing scatterplots
by robust locally weighted regression.
The American Statistician, 35, 54.
==================
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*-*-*-*-*-*-*-*-*-*-*-*-*-*
==========================
Based on R function clowess: Translated to C++ by C. Stratowa
(R source file: lowess.c by R Development Core Team (C) 1999-2001)
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
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. Stratowa
(R source file: lowess.c by R Development Core Team (C) 1999-2001)
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
TGraph* SmoothSuper(TGraph *grin, Option_t *,
Double_t bass, Double_t span, Bool_t isPeriodic, Double_t *w)
*-*-*-*-*-*-*-*-*Smooth data with Super smoother*-*-*-*-*-*-*-*-*-*-*-*
===============================
Smooth the (x, y) values by Friedman's ``super smoother''.
Arguments:
grin: graph for smoothing
span: the fraction of the observations in the span of the running lines
smoother, or 0 to choose this by leave-one-out cross-validation.
bass: controls the smoothness of the fitted curve.
Values of up to 10 indicate increasing smoothness.
isPeriodic: if TRUE, the x values are assumed to be in [0, 1]
and of period 1.
w: case weights
Details:
supsmu is a running lines smoother which chooses between three spans for
the lines. The running lines smoothers are symmetric, with k/2 data points
each side of the predicted point, and values of k as 0.5 * n, 0.2 * n and
0.05 * n, where n is the number of data points. If span is specified,
a single smoother with span span * n is used.
The best of the three smoothers is chosen by cross-validation for each
prediction. The best spans are then smoothed by a running lines smoother
and the final prediction chosen by linear interpolation.
The FORTRAN code says: ``For small samples (n < 40) or if there are
substantial serial correlations between observations close in x - value,
then a prespecified fixed span smoother (span > 0) should be used.
Reasonable span values are 0.2 to 0.4.''
References:
Friedman, J. H. (1984) SMART User's Guide.
Laboratory for Computational Statistics,
Stanford University Technical Report No. 1.
Friedman, J. H. (1984) A variable span scatterplot smoother.
Laboratory for Computational Statistics,
Stanford University Technical Report No. 5.
==================
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)
*-*-*-*-*-*-*-*-*Friedmannīs super smoother *-*-*-*-*-*-*-*-*-*-*-*-*-*
==========================
super smoother (Friedman, 1984).
version 10/10/84
coded and copywrite (c) 1984 by:
Jerome H. Friedman
department of statistics
and
stanford linear accelerator center
stanford university
all rights reserved.
input:
n : number of observations (x,y - pairs).
x(n) : ordered abscissa values.
y(n) : corresponding ordinate (response) values.
w(n) : weight for each (x,y) observation.
iper : periodic variable flag.
iper=1 => x is ordered interval variable.
iper=2 => x is a periodic variable with values
in the range (0.0,1.0) and period 1.0.
span : smoother span (fraction of observations in window).
span=0.0 => automatic (variable) span selection.
alpha : controls high frequency (small span) penality
used with automatic span selection (bass tone control).
(alpha.le.0.0 or alpha.gt.10.0 => no effect.)
output:
smo(n) : smoothed ordinate (response) values.
scratch:
sc(n,7) : internal working storage.
note:
for small samples (n < 40) or if there are substantial serial
correlations between observations close in x - value, then
a prespecified fixed span smoother (span > 0) should be
used. reasonable span values are 0.2 to 0.4.
current implementation:
Based on R function supsmu: Translated to C++ by C. Stratowa
(R source file: ppr.f by B.D.Ripley Copyright (C) 1994-97)
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
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. Stratowa
(R source file: ppr.f by B.D.Ripley Copyright (C) 1994-97)
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
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*-*-*-*
==============================================
TGraph* Approx(TGraph *grin, Option_t *option, Int_t nout, Double_t *xout,
Double_t yleft, Double_t yright, Int_t rule, Double_t f, Option_t *ties)
*-*-*-*-*-*-*-*-*Approximate data points*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
=======================
Arguments:
grin: graph giving the coordinates of the points to be interpolated.
Alternatively a single plotting structure can be specified:
option: specifies the interpolation method to be used.
Choices are "linear" (iKind = 1) or "constant" (iKind = 2).
nout: If xout is not specified, interpolation takes place at n equally
spaced points spanning the interval [min(x), max(x)], where
nout = max(nout, number of input data).
xout: an optional set of values specifying where interpolation is to
take place.
yleft: the value to be returned when input x values less than min(x).
The default is defined by the value of rule given below.
yright: the value to be returned when input x values greater than max(x).
The default is defined by the value of rule given below.
rule: an integer describing how interpolation is to take place outside
the interval [min(x), max(x)]. If rule is 0 then the given yleft
and yright values are returned, if it is 1 then 0 is returned
for such points and if it is 2, the value at the closest data
extreme is used.
f: For method="constant" a number between 0 and 1 inclusive,
indicating a compromise between left- and right-continuous step
functions. If y0 and y1 are the values to the left and right of
the point then the value is y0*f+y1*(1-f) so that f=0 is
right-continuous and f=1 is left-continuous
ties: Handling of tied x values. An integer describing a function with
a single vector argument returning a single number result:
ties = "ordered" (iTies = 0): input x are "ordered"
ties = "mean" (iTies = 1): function "mean"
ties = "min" (iTies = 2): function "min"
ties = "max" (iTies = 3): function "max"
Details:
At least two complete (x, y) pairs are required.
If there are duplicated (tied) x values and ties is a function it is
applied to the y values for each distinct x value. Useful functions in
this context include mean, min, and max.
If ties="ordered" the x values are assumed to be already ordered. The
first y value will be used for interpolation to the left and the last
one for interpolation to the right.
Value:
approx returns a graph with components x and y, containing n coordinates
which interpolate the given data points according to the method (and rule)
desired.
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*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-* ==========================
Approximate y(v), given (x,y)[i], i = 0,..,n-1
Based on R function approx1: Translated to C++ by Christian Stratowa
(R source file: approx.c by R Development Core Team (C) 1999-2001)
Int_t Rcmp(Double_t x, Double_t y)
static function
if (ISNAN(x)) return 1;
if (ISNAN(y)) return -1;
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_sort.c by R Development Core Team (C) 1999-2001)
void Rank(Int_t n, Double_t *a, Int_t *index, Int_t *rank, Bool_t down)
static function
Inline Functions
TClass* Class()
TClass* IsA() const
void ShowMembers(TMemberInspector& insp, char* parent)
void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
TGraphSmooth TGraphSmooth(const TGraphSmooth&)
TGraphSmooth& operator=(const TGraphSmooth&)
Author: Christian Stratowa 30/09/2001
Last update: root/graf:$Name: $:$Id: TGraphSmooth.cxx,v 1.9 2005/09/05 07:25:22 brun Exp $
ROOT page - Class index - Class Hierarchy - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.