0.6.2-dev0
FORCES
FORtran lib for Comp. Env. Sys.
Loading...
Searching...
No Matches
mo_kernel::kernel_regression Interface Reference

Multi-dimensional non-parametric kernel regression. More...

Public Member Functions

real(dp) function, dimension(:), allocatable kernel_regression_2d_dp (ix, iy, h, silverman, xout, mask, nodata)
 
real(sp) function, dimension(:), allocatable kernel_regression_2d_sp (ix, iy, h, silverman, xout, mask, nodata)
 
real(dp) function, dimension(:), allocatable kernel_regression_1d_dp (ix, iy, h, silverman, xout, mask, nodata)
 
real(sp) function, dimension(:), allocatable kernel_regression_1d_sp (ix, iy, h, silverman, xout, mask, nodata)
 

Detailed Description

Multi-dimensional non-parametric kernel regression.

Multi-dimensional non-parametric kernel regression using a Gaussian kernel.

The bandwith of the kernel can be pre-determined using the function kernel_regression_h and specified by the optional argument h. If h is not given the default method to approximate the bandwith h is Silverman''s rule-of-thumb

\[ h = \frac{4}{d+2}^{1/(d+4)} n^{-1/(d+4)} \sigma_{x_i} \]

where \( n \) is the number of given data points, \( d \) is the number of dimensions, and \( \sigma_{x_i} \) is the standard deviation of the data of dimension \( i \).
If the optional argument silverman is set to false, the cross-validation method described by Scott et al. (2005) is applied. Therefore, the likelihood of a given h is maximized using the Nelder-Mead algorithm nelminrange. For large data sets this might be time consuming and should be performed aforehand using the function kernel_regression_h.
The dataset (x,y) can be single or double precision. The result will have the same numerical precision.
If the regression for other datapoints than x is needed the optional argument xout has to be specified. The result will than be of the same size and precision as xout.
Example

the code is adapted from the kernel_regression.py of the jams_python library.
! given data, e.g. temperature
x(:,1) = (/ 26.1_dp, 24.5_dp, 24.8_dp, 14.5_dp, 2.1_dp /)
x(:,2) = (/ 2.1_dp, 4.5_dp, 6.8_dp, 4.8_dp, 0.1_dp /)
y = (/ 28.2_dp, 29.0_dp, 31.6_dp, 19.3_dp, 2.2_dp /)
! estimate bandwidth via cross-validation (time-consuming)
h = kernel_regression_h(x,y,silverman=.false.)
! estimate bandwidth with Silverman''s rule of thumb (default)
h = kernel_regression_h(x,y,silverman=.true.)
fit_y = kernel_regression(x, y, h=h, silverman=.false., xout=xout)
! gives fitted values at either xout values, if specified, or at x values, if xout is not present
! if bandwith h is given : silverman (true/false) is ignored
! if silverman=.true. and h not present : bandwith will be estimated using Silerman''s rule of thumb
! if silverman=.false. and h not present : bandwith will be estimated using cross-validation approach

-> see also example in test directory

Literature

  1. Scott, D. W., & Sain, S. R. (2005). Multi-dimensional Density Estimation. Handbook of Statistics, 24, 229-261. doi:10.1016/S0169-7161(04)24009-3
  2. Haerdle, W., & Mueller, M. (2000). Multivariate and semiparametric kernel regression. In M. G. Schimek (Ed.), Smoothing and regression: Approaches, computation, and application (pp. 357-392). Hoboken, NJ, USA: John Wiley & Sons, Inc. doi:10.1002/9781118150658.ch12s, Cambridge University Press, UK, 1996
Parameters
[in]real(sp/dp) :: x(:)/x(:,:)1D or ND array with independent variables
[in]real(sp/dp) :: y(:)1D array of dependent variables y(i) = f(x(i:)) with unknown f
[in]real(sp/dp), optional :: hBandwith of kernel.
If present, argument silverman is ignored. If not present, the bandwith will be approximated first.
[in]logical, optional :: silvermanBy default Silverman''s Rule-of-thumb will be used to approximate the bandwith of the kernel (silverman=true). If silverman=false the Cross-Validation approach is used to estimate the bandwidth.
[in]real(sp/dp), optional :: xout(:)/xout(:,:)If present, the fitted values will be returned for this independent variables, otherwise the fitted values at x are returned.
[in]logical, optional :: mask(:)mask y values at calculation.
if not xout given, then kernel estimates will have nodata value.
[in]real(sp/dp), optional :: nodataif mask and not xout given, then masked data will have nodata kernel estimate.
Returns
real(sp/dp), allocatable :: kernel_regression(:) — fitted values at eighter x or xout
Authors
Matthias Cuntz, Juliane Mai
Date
Mar 2013

Definition at line 339 of file mo_kernel.f90.

Member Function/Subroutine Documentation

◆ kernel_regression_1d_dp()

real(dp) function, dimension(:), allocatable mo_kernel::kernel_regression::kernel_regression_1d_dp ( real(dp), dimension(:), intent(in)  ix,
real(dp), dimension(:), intent(in)  iy,
real(dp), intent(in), optional  h,
logical, intent(in), optional  silverman,
real(dp), dimension(:), intent(in), optional  xout,
logical, dimension(:), intent(in), optional  mask,
real(dp), intent(in), optional  nodata 
)

Definition at line 1161 of file mo_kernel.f90.

◆ kernel_regression_1d_sp()

real(sp) function, dimension(:), allocatable mo_kernel::kernel_regression::kernel_regression_1d_sp ( real(sp), dimension(:), intent(in)  ix,
real(sp), dimension(:), intent(in)  iy,
real(sp), intent(in), optional  h,
logical, intent(in), optional  silverman,
real(sp), dimension(:), intent(in), optional  xout,
logical, dimension(:), intent(in), optional  mask,
real(sp), intent(in), optional  nodata 
)

Definition at line 1259 of file mo_kernel.f90.

◆ kernel_regression_2d_dp()

real(dp) function, dimension(:), allocatable mo_kernel::kernel_regression::kernel_regression_2d_dp ( real(dp), dimension(:,:), intent(in)  ix,
real(dp), dimension(:), intent(in)  iy,
real(dp), dimension(:), intent(in), optional  h,
logical, intent(in), optional  silverman,
real(dp), dimension(:,:), intent(in), optional  xout,
logical, dimension(:), intent(in), optional  mask,
real(dp), intent(in), optional  nodata 
)

Definition at line 1357 of file mo_kernel.f90.

◆ kernel_regression_2d_sp()

real(sp) function, dimension(:), allocatable mo_kernel::kernel_regression::kernel_regression_2d_sp ( real(sp), dimension(:,:), intent(in)  ix,
real(sp), dimension(:), intent(in)  iy,
real(sp), dimension(:), intent(in), optional  h,
logical, intent(in), optional  silverman,
real(sp), dimension(:,:), intent(in), optional  xout,
logical, dimension(:), intent(in), optional  mask,
real(sp), intent(in), optional  nodata 
)

Definition at line 1462 of file mo_kernel.f90.


The documentation for this interface was generated from the following file: