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

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

Public Member Functions

real(dp) function kernel_density_h_1d_dp (ix, silverman, mask)
 
real(sp) function kernel_density_h_1d_sp (ix, silverman, mask)
 

Detailed Description

Approximates the bandwith h of the kernel.

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}{3}^{0.2} n^{-0.2} \sigma_x \]

where n is the number of given data points and \( \sigma \) is the standard deviation of the data.
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_density and kernel_cumdensity.

Example

! given data, e.g. temperature
x = (/ 26.1_dp, 24.5_dp, 24.8_dp, 24.5_dp, 24.1_dp /)
! estimate bandwidth via cross-validation (time-consuming)
h = kernel_density_h(x,silverman=.false.)
! estimate bandwidth with Silverman''s rule of thumb (default)
h = kernel_density_h(x,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_i \) 1D-array with data points
[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_density_h(:) — approximated bandwith h for kernel smoother
Note
Data need to be one-dimensional. Multi-dimensional data handling not implemented yet.
Authors
Juliane Mai, Matthias Cuntz
Date
Mar 2013

Definition at line 262 of file mo_kernel.f90.

Member Function/Subroutine Documentation

◆ kernel_density_h_1d_dp()

real(dp) function mo_kernel::kernel_density_h::kernel_density_h_1d_dp ( real(dp), dimension(:), intent(in)  ix,
logical, intent(in), optional  silverman,
logical, dimension(:), intent(in), optional  mask 
)

Definition at line 1056 of file mo_kernel.f90.

◆ kernel_density_h_1d_sp()

real(sp) function mo_kernel::kernel_density_h::kernel_density_h_1d_sp ( real(sp), dimension(:), intent(in)  ix,
logical, intent(in), optional  silverman,
logical, dimension(:), intent(in), optional  mask 
)

Definition at line 1107 of file mo_kernel.f90.


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