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

Approximates the cumulative density function (CDF). More...

Public Member Functions

real(dp) function, dimension(:), allocatable kernel_cumdensity_1d_dp (ix, h, silverman, xout, romberg, nintegrate, epsint, mask, nodata)
 
real(sp) function, dimension(:), allocatable kernel_cumdensity_1d_sp (ix, h, silverman, xout, romberg, nintegrate, epsint, mask, nodata)
 

Detailed Description

Approximates the cumulative density function (CDF).

Approximates the cumulative density function (CDF) to a given 1D data set using a Gaussian kernel.

The bandwith of the kernel can be pre-determined using the function kernel_density_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}{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.
If the CDF 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

! 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.)
! calc cumulative density with the estimated bandwidth h at given output points xout
cdf = kernel_cumdensity(x, h=h, xout=xout)
! gives cumulative density at 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_i \) 1D-array with data points
[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(:)If present, the CDF will be approximated at this arguments, otherwise the CDF is approximated at x.
[in]logical, optional :: rombergIf .true., use Romberg converging integration If .true., use 5-point Newton-Cotes with fixed nintegrate points
[in]integer(i4), optional :: nintegrate2**nintegrate if Romberg or nintegrate number of sampling points for integration between output points. Default: 10 if Romberg, 101 otherwise.
[in]real(sp/dp), optional :: epsintmaximum relative error for Romberg integration. Default: 1e-6.
[in]logical, optional :: mask(:)mask x 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_cumdensity(:) — smoothed CDF at either x or xout
Note
Data need to be one-dimensional. Multi-dimensional data handling not implemented yet.
Author
Juliane Mai
Date
Mar 2013
Author
Matthias Cuntz
Date
Mar 2013
May 2013
  • sort -> qsort
Author
Matthias Cuntz
Date
Mar 2016
  • Romberg integration in cumdensity

Definition at line 128 of file mo_kernel.f90.

Member Function/Subroutine Documentation

◆ kernel_cumdensity_1d_dp()

real(dp) function, dimension(:), allocatable mo_kernel::kernel_cumdensity::kernel_cumdensity_1d_dp ( real(dp), dimension(:), intent(in)  ix,
real(dp), intent(in), optional  h,
logical, intent(in), optional  silverman,
real(dp), dimension(:), intent(in), optional  xout,
logical, intent(in), optional  romberg,
integer(i4), intent(in), optional  nintegrate,
real(dp), intent(in), optional  epsint,
logical, dimension(:), intent(in), optional  mask,
real(dp), intent(in), optional  nodata 
)

Definition at line 460 of file mo_kernel.f90.

◆ kernel_cumdensity_1d_sp()

real(sp) function, dimension(:), allocatable mo_kernel::kernel_cumdensity::kernel_cumdensity_1d_sp ( real(sp), dimension(:), intent(in)  ix,
real(sp), intent(in), optional  h,
logical, intent(in), optional  silverman,
real(sp), dimension(:), intent(in), optional  xout,
logical, intent(in), optional  romberg,
integer(i4), intent(in), optional  nintegrate,
real(sp), intent(in), optional  epsint,
logical, dimension(:), intent(in), optional  mask,
real(sp), intent(in), optional  nodata 
)

Definition at line 659 of file mo_kernel.f90.


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