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

Approximates the bandwith h of the kernel for regression. More...

Public Member Functions

real(dp) function, dimension(size(ix, 2)) kernel_regression_h_2d_dp (ix, iy, silverman, mask)
 
real(sp) function, dimension(size(ix, 2)) kernel_regression_h_2d_sp (ix, iy, silverman, mask)
 
real(dp) function kernel_regression_h_1d_dp (ix, iy, silverman, mask)
 
real(sp) function kernel_regression_h_1d_sp (ix, iy, silverman, mask)
 

Detailed Description

Approximates the bandwith h of the kernel for regression.

Approximates the bandwith h of the kernel for a given dataset x either using Silverman''s rule-of-thumb or a cross-validation method.
By default the bandwidth h is approximated by 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_density_h.
The dataset x can be single or double precision. The result will have the same numerical precision.
The result of this function can be used as the optional input for kernel_regression.
The code is adapted from the kernel_regression.py of the UFZ CHS Python library.

Example

! 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.)

-> 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]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]logical, optional :: mask(:)mask x values at calculation.
Returns
real(sp/dp), allocatable :: kernel_regression_h(:) — approximated bandwith h for kernel regression
number of bandwidths equals number of independent variables
Authors
Matthias Cuntz, Juliane Mai
Date
Mar 2013

Definition at line 401 of file mo_kernel.f90.

Member Function/Subroutine Documentation

◆ kernel_regression_h_1d_dp()

real(dp) function mo_kernel::kernel_regression_h::kernel_regression_h_1d_dp ( real(dp), dimension(:), intent(in)  ix,
real(dp), dimension(:), intent(in)  iy,
logical, intent(in), optional  silverman,
logical, dimension(:), intent(in), optional  mask 
)

Definition at line 1570 of file mo_kernel.f90.

◆ kernel_regression_h_1d_sp()

real(sp) function mo_kernel::kernel_regression_h::kernel_regression_h_1d_sp ( real(sp), dimension(:), intent(in)  ix,
real(sp), dimension(:), intent(in)  iy,
logical, intent(in), optional  silverman,
logical, dimension(:), intent(in), optional  mask 
)

Definition at line 1627 of file mo_kernel.f90.

◆ kernel_regression_h_2d_dp()

real(dp) function, dimension(size(ix,2)) mo_kernel::kernel_regression_h::kernel_regression_h_2d_dp ( real(dp), dimension(:,:), intent(in)  ix,
real(dp), dimension(:), intent(in)  iy,
logical, intent(in), optional  silverman,
logical, dimension(:), intent(in), optional  mask 
)

Definition at line 1684 of file mo_kernel.f90.

◆ kernel_regression_h_2d_sp()

real(sp) function, dimension(size(ix,2)) mo_kernel::kernel_regression_h::kernel_regression_h_2d_sp ( real(sp), dimension(:,:), intent(in)  ix,
real(sp), dimension(:), intent(in)  iy,
logical, intent(in), optional  silverman,
logical, dimension(:), intent(in), optional  mask 
)

Definition at line 1742 of file mo_kernel.f90.


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