https://github.com/cran/spatstat
Raw File
Tip revision: bd9c4530fb27ff8cdce3e3ab11283555cd489a10 authored by Adrian Baddeley on 27 September 2006, 21:21:44 UTC
version 1.9-6
Tip revision: bd9c453
Kest.Rd
\name{Kest}
\alias{Kest}
\title{K-function}
\description{
Estimates the reduced second moment function \eqn{K(r)} 
from a point pattern in a window of arbitrary shape.
}
\synopsis{
  Kest(X, \dots, r=NULL, breaks=NULL, 
     correction=c("border", "isotropic", "Ripley", "translate"),
    nlarge=3000)
}
\usage{
  Kest(X, \dots, r, correction=c("border", "isotropic", "Ripley", "translate"),
     nlarge=3000)
}
\arguments{
  \item{X}{The observed point pattern, 
    from which an estimate of \eqn{K(r)} will be computed.
    An object of class \code{"ppp"}, or data
    in any format acceptable to \code{\link{as.ppp}()}.
  }
  \item{\dots}{Ignored.}
  \item{r}{
    Optional. Vector of values for the argument \eqn{r} at which \eqn{K(r)} 
    should be evaluated. There is a sensible default.
  }
  \item{correction}{
    Optional. A character vector containing any selection of the
    options \code{"border"}, \code{"bord.modif"},
    \code{"isotropic"}, \code{"Ripley"} or \code{"translate"}.
    It specifies the edge correction(s) to be applied.
  }
  \item{nlarge}{
    Optional. Efficiency threshold.
    If the number of points exceeds \code{nlarge}, then only the
    border correction will be computed, using a fast algorithm.
  }
}
\value{
  An object of class \code{"fv"}, see \code{\link{fv.object}},
  which can be plotted directly using \code{\link{plot.fv}}.

  Essentially a data frame containing columns
  \item{r}{the vector of values of the argument \eqn{r} 
    at which the function \eqn{K} has been  estimated
  }
  \item{theo}{the theoretical value \eqn{K(r) = \pi r^2}{K(r) = pi * r^2}
    for a stationary Poisson process
  }
  together with columns named 
  \code{"border"}, \code{"bord.modif"},
  \code{"iso"} and/or \code{"trans"},
  according to the selected edge corrections. These columns contain
  estimates of the function \eqn{K(r)} obtained by the edge corrections
  named.
}
\details{
  The \eqn{K} function (variously called ``Ripley's K-function''
  and the ``reduced second moment function'')
  of a stationary point process \eqn{X} is defined so that
  \eqn{\lambda K(r)}{lambda K(r)} equals the expected number of
  additional random points within a distance \eqn{r} of a
  typical random point of \eqn{X}. Here \eqn{\lambda}{lambda}
  is the intensity of the process,
  i.e. the expected number of points of \eqn{X} per unit area.
  The \eqn{K} function is determined by the 
  second order moment properties of \eqn{X}.
 
  An estimate of \eqn{K} derived from a spatial point pattern dataset
  can be used in exploratory data analysis and formal inference
  about the pattern (Cressie, 1991; Diggle, 1983; Ripley, 1988).
  In exploratory analyses, the estimate of \eqn{K} is a useful statistic 
  summarising aspects of inter-point ``dependence'' and ``clustering''.
  For inferential purposes, the estimate of \eqn{K} is usually compared to the 
  true value of \eqn{K} for a completely random (Poisson) point process,
  which is \eqn{K(r) = \pi r^2}{K(r) = pi * r^2}.
  Deviations between the empirical and theoretical \eqn{K} curves
  may suggest spatial clustering or spatial regularity.
 
  This routine \code{Kest} estimates the \eqn{K} function
  of a stationary point process, given observation of the process
  inside a known, bounded window. 
  The argument \code{X} is interpreted as a point pattern object 
  (of class \code{"ppp"}, see \code{\link{ppp.object}}) and can
  be supplied in any of the formats recognised by
  \code{\link{as.ppp}()}.

  The estimation of \eqn{K} is hampered by edge effects arising from 
  the unobservability of points of the random pattern outside the window. 
  An edge correction is needed to reduce bias (Baddeley, 1998; Ripley, 1988). 
  The corrections implemented here are
  \describe{
    \item{border}{the border method or
      ``reduced sample'' estimator (see Ripley, 1988). This is
      the least efficient (statistically) and the fastest to compute.
      It can be computed for a window of arbitrary shape.
    }
    \item{isotropic/Ripley}{Ripley's isotropic correction
      (see Ripley, 1988; Ohser, 1983).
      This is implemented for rectangular and polygonal windows
      (not for binary masks).
    }
    \item{translate}{Translation correction (Ohser, 1983).
      Implemented for all window geometries, but slow for
      complex windows. 
    }
  }
  Note that the estimator assumes the process is stationary (spatially
  homogeneous). For inhomogeneous point patterns, see
  \code{\link{Kinhom}}.

  If the point pattern \code{X} contains more than about 3000 points,
  the isotropic and translation edge corrections can be computationally
  prohibitive. The computations for the border method are much faster,
  and are statistically efficient when there are large numbers of
  points. Accordingly, if the number of points in \code{X} exceeds
  the threshold \code{nlarge}, then only the border correction will be
  computed. Setting \code{nlarge=Inf} will prevent this from happening.
  Setting \code{nlarge=0} is equivalent to selecting only the border
  correction with \code{correction="border"}.

  The estimator \code{Kest} ignores marks.
  Its counterparts for multitype point patterns
  are \code{\link{Kcross}}, \code{\link{Kdot}},
  and for general marked point patterns
  see \code{\link{Kmulti}}. 

  Some writers, particularly Stoyan (1994, 1995) advocate the use of
  the ``pair correlation function''
  \deqn{
    g(r) = \frac{K'(r)}{2\pi r}
  }{
    g(r) = K'(r)/ ( 2 * pi * r) 
  }
  where \eqn{K'(r)} is the derivative of \eqn{K(r)}.
  See \code{\link{pcf}} on how to estimate this function.
}
\references{
Baddeley, A.J. Spatial sampling and censoring.
     In O.E. Barndorff-Nielsen, W.S. Kendall and
     M.N.M. van Lieshout (eds) 
     \emph{Stochastic Geometry: Likelihood and Computation}.
     Chapman and Hall, 1998.
     Chapter 2, pages 37--78.
  
  Cressie, N.A.C. \emph{Statistics for spatial data}.
    John Wiley and Sons, 1991.

  Diggle, P.J. \emph{Statistical analysis of spatial point patterns}.
  Academic Press, 1983.

  Ohser, J. (1983)
  On estimators for the reduced second moment measure of
  point processes. \emph{Mathematische Operationsforschung und
  Statistik, series Statistics}, \bold{14}, 63 -- 71.
    
  Ripley, B.D. \emph{Statistical inference for spatial processes}.
  Cambridge University Press, 1988.

  Stoyan, D, Kendall, W.S. and Mecke, J. (1995)
  \emph{Stochastic geometry and its applications}.
  2nd edition. Springer Verlag.

  Stoyan, D. and Stoyan, H. (1994)
  Fractals, random shapes and point fields:
  methods of geometrical statistics.
  John Wiley and Sons.
} 
\section{Warnings}{
  The estimator of \eqn{K(r)} is approximately unbiased for each fixed \eqn{r}.
  Bias increases with \eqn{r} and depends on the window geometry.
  For a rectangular window it is prudent to restrict the \eqn{r} values to
  a maximum of \eqn{1/4} of the smaller side length of the rectangle.
  Bias may become appreciable for point patterns consisting of 
  fewer than 15 points.
 
  While \eqn{K(r)} is always a non-decreasing function, the estimator 
  of \eqn{K} is not guaranteed to be non-decreasing. This is rarely 
  a problem in practice.
}
\seealso{
  \code{\link{Fest}},
  \code{\link{Gest}},
  \code{\link{Jest}},
  \code{\link{pcf}},
  \code{\link{reduced.sample}},
  \code{\link{Kcross}},
  \code{\link{Kdot}},
  \code{\link{Kinhom}},
  \code{\link{Kmulti}}
}
\examples{
 pp <- runifpoint(50)
 K <- Kest(pp)
 data(cells)
 K <- Kest(cells, correction="isotropic")
 plot(K)
 plot(K, main="K function for cells")
 # plot the L function
 plot(K, sqrt(iso/pi) ~ r)
 plot(K, sqrt(./pi) ~ r, ylab="L(r)", main="L function for cells")
}
\author{Adrian Baddeley
  \email{adrian@maths.uwa.edu.au}
  \url{http://www.maths.uwa.edu.au/~adrian/}
  and Rolf Turner
  \email{rolf@math.unb.ca}
  \url{http://www.math.unb.ca/~rolf}
}
\keyword{spatial}
\keyword{nonparametric}

 
 
back to top