Package 'spatstat.local'

Title: Extension to 'spatstat' for Local Composite Likelihood
Description: Extension to the 'spatstat' package, enabling the user to fit point process models to point pattern data by local composite likelihood ('geographically weighted regression').
Authors: Adrian Baddeley [aut, cre]
Maintainer: Adrian Baddeley <[email protected]>
License: GPL (>= 2)
Version: 5.1-0
Built: 2024-09-14 05:38:01 UTC
Source: https://github.com/baddstats/spatstat.local

Help Index


Local Composite Likelihood

Description

Extension of the spatstat package, for fitting spatial point process models by local composite likelihood.

Details

The main functions are

locppm Local likelihood fit of Poisson model
Local pseudolikelihood fit of Gibbs model
locmincon Local minimum contrast fit
of Neyman-Scott or Cox model
loccit Local composite likelihood fit
of Neyman-Scott or Cox model

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.


Cross Validated Bandwidth Selection for Locally Fitted Point Process Model

Description

Uses cross-validation to select a smoothing bandwidth for locally fitting a Cox or cluster point process model.

Usage

bw.loccit(..., use.fft=TRUE, 
               srange = NULL, ns = 9, sigma = NULL,
               fftopt=list(), 
               verbose = TRUE)

Arguments

...

Arguments passed to kppm to fit the homogeneous version of the model.

use.fft

Logical value indicating whether to use a quick-and-dirty approximation based on a first order Taylor expansion.

srange

Range of values of the smoothing parameter sigma to be searched. A numeric vector of length 2 giving the minimum and maximum values of sigma.

ns

Number of values of the smoothing parameter sigma in the range srange to be searched. A positive integer.

sigma

Vector of values of the smoothing parameter to be searched.

fftopt

Developer use only.

verbose

Logical value indicating whether to display progress reports.

Details

This function determines the optimal value of the smoothing parameter sigma to be used in a call to loccit.

The function loccit fits a Cox or cluster point process model to point pattern data by local composite likelihood. The degree of local smoothing is controlled by a smoothing parameter sigma which is an argument to loccit.

For each value of sigma in a search interval, the function bw.loccit fits the model locally and evaluates a cross-validation criterion. The optimal value of sigma is returned.

Value

A numerical value giving the selected bandwidth. The result also belongs to the class "bw.optim" which can be plotted.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

loccit

Examples

X <- redwood[owin(c(0,1), c(-1,-1/2))]
  Ns <- if(interactive()) 16 else 2
  b <- bw.loccit(X, ~1, "Thomas", srange=c(0.07, 0.14), ns=Ns)
  b
  plot(b)

Cross Validated Bandwidth Selection for Locally Fitted Point Process Model

Description

Uses cross-validation to select a smoothing bandwidth for locally fitting a Poisson or Gibbs point process model.

Usage

bw.locppm(...,
          method = c("fft", "exact", "taylor"), 
          srange = NULL, ns = 9, sigma = NULL,
          additive = TRUE, 
          verbose = TRUE)

Arguments

...

Arguments passed to ppm to fit the homogeneous version of the model.

method

Method of calculation. The default method="fft" is much faster than the other choices.

srange

Range of values of the smoothing parameter sigma to be searched. A numeric vector of length 2 giving the minimum and maximum values of sigma.

ns

Number of values of the smoothing parameter sigma in the range srange to be searched. A positive integer.

sigma

Vector of values of the smoothing parameter to be searched. Overrides the values of ns and srange.

additive

Logical value indicating whether to calculate the leverage approximation on the scale of the intensity (additive=TRUE) or the log intensity (additive=FALSE). Applies only when method = "taylor".

verbose

Logical value indicating whether to display progress reports.

Details

This function determines the optimal value of the smoothing parameter sigma to be used in a call to locppm.

The function locppm fits a Poisson or Gibbs point process model to point pattern data by local composite likelihood. The degree of local smoothing is controlled by a smoothing parameter sigma which is an argument to locppm.

This function bw.locppm determines the optimal value of sigma by cross-validation. For each value of sigma in a search interval, the function bw.locppm fits the model locally with smoothing bandwidth sigma, and evaluates the composite likelihood cross-validation criterion LCV(sigma) defined in Baddeley (2016), section 3.2. The value of sigma which maximises LCV(sigma) is returned.

Value

A numerical value giving the selected bandwidth. The result also belongs to the class "bw.optim" which can be plotted.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm

Examples

Ns <- if(interactive()) 16 else 2
  b <- bw.locppm(swedishpines, ~1, srange=c(2.5,4.5), ns=Ns)
  b
  plot(b)

Homogeneity Test for Local Poisson or Gibbs Model

Description

Conducts a Monte Carlo test of homogeneity for a locally-fitted Poisson or Gibbs point process model.

Usage

homtest(X, ..., nsim = 19,
           test = c("residuals", "score", "taylor", "likelihood"),
           locations = c("coarse", "fine", "split"),
           ladjust = NULL,
           use.fft = NULL,
           simul = NULL,
           verbose = TRUE, Xname = NULL)

Arguments

X

A point pattern (object of class "ppp").

...

Additional arguments passed to locppm to determine the locally fitted model, and passed to ppm to determine the homogeneous model.

nsim

Number of simulations for the Monte Carlo test.

test

The local test statistic to be used: either "likelihood" for the local likelihood ratio test statistic, "taylor" for the Taylor approximation to the local likelihood ratio test statistic, "score" for the local score test statistic, or "residuals" for the squared local residuals.

locations, use.fft

Arguments passed to locppm to control the calculation of variances (if method="local"). See Details.

ladjust

Optional argument passed to homtestmap specifying a data-dependent adjustment of the test statistic.

simul

Optional information that determines how to simulate the realisations from the null hypothesis. An expression in the R language (that will be evaluated nsim times to obtain the simulated patterns), or a list that contains the simulated point patterns.

verbose

Logical value indicating whether to print progress reports.

Xname

Optional character string name for the dataset X, to be printed in the test output.

Details

This function performs a Monte Carlo test of the null hypothesis of homogeneity (i.e.\ constant parameter values) for the locally-fitted Poisson point process or Gibbs point process specified by the arguments.

The type of test is controlled by the argument test.

  • test="likelihood": the locally fitted model is computed as locppm(X, ...). The local composite likelihood ratio test statistic of this model is computed at each location, and the mean of this statistic over the window is computed.

  • test="taylor": the locally fitted model is computed as locppm(X, ...). The Taylor approximation to the local composite likelihood ratio test statistic of this model is computed at each location, and the mean of this statistic over the window is computed.

  • test="score": the locally fitted model is computed as locppm(X, ...). The local score test statistic of this model is computed at each location, and the mean of this statistic over the window is computed.

  • method="residuals": the homogeneous model is fitted as ppm(X, ...). The smoothed score residuals of this model are computed at each location, and the mean of the squared norm over the window is computed.

The test statistic is computed for the data pattern X and for each of nsim simulated realisations from the homogeneous model. The Monte Carlo pp-value is computed.

Value

An object of class "htest" containing the test outcome.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

To compute the test statistic only, see homtestmap.

Examples

## Not run: 
   homtest(swedishpines)
   
## End(Not run)

Test Statistic for Homogeneity Test

Description

Compute the test statistic for the test of homogeneity of a locally-fitted Poisson or Gibbs point process model.

Usage

homteststat(object, ..., verbose = FALSE)

homtestmap(object, ...,
           what=c("components", "statistic", "pvalue"),
           test = c("score", "taylor", "likelihood"),
           ladjust=c("none", "moment", "PSS"),
           calibrate=c("chisq", "Satterthwaite", "firstmoment"),
           simple = !is.null(theta0),
           theta0 = NULL,
           poolmoments=NULL,
           sigma = NULL, 
           saveall = FALSE, 
           use.fft = TRUE,
           verbose = TRUE)

## S3 method for class 'homtestmap'
update(object, ..., 
           what=NULL, test=NULL, ladjust=NULL,
           calibrate=NULL, saveall=FALSE, poolmoments=NULL)

Arguments

object

Locally-fitted point process (object of class "locppm") or an object previously computed by homtestmap.

...

For homteststat, arguments passed to homtestmap. For homtestmap, additional unmatched arguments are ignored.

what

Character string (partially matched) indicating whether to return the vector components of the local test statistic, or the value of the local test statistic, or the local pp-values.

test

Character string (partially matched) indicating whether to perform the local score test (test="score"), or the local composite likelihood ratio test approximately (test="taylor") or exactly (test="likelihood").

ladjust

Character string (partially matched) specifying an adjustment to the composite likelihood ratio test statistic.

calibrate

Character string (partially matched) specifying how to calculate pp-values from the test statistic.

simple

Logical value indicating whether to treat the fitted model as a simple null hypothesis (simple=TRUE) or as an estimate of the parameters in a composite null hypothesis (simple=FALSE, the default).

theta0

Coefficient vector specifying a simple null hypothesis.

poolmoments

Logical value indicating how to calculate the reference distribution for the likelihood ratio test statistic (and thus how to calculate pp-values). See Details.

sigma

Smoothing bandwidth.

saveall

Logical value indicating whether to compute a complete set of sufficient statistics and save them as an attribute of the result. See Details.

use.fft

For software testing purposes only. Logical value indicating whether to use data computed by the Fast Fourier Transform.

verbose

Logical value indicating whether to print progress reports.

Details

These functions are used by homtest to perform a Monte Carlo test of the null hypothesis of homogeneity (i.e.\ constant parameter values) for the locally-fitted Poisson point process or Gibbs point process object.

The function homtestmap computes either the local likelihood ratio test statistic or the local score test statistic. If what="statistic", then the result is a scalar-valued function giving the local values of the test statistic. If what="pvalue", the result is a scalar-valued function p(v)p(v) giving the local pp-value at each location vv. If what="components", the result is a vector-valued function T(v)T(v) containing the components of the quadratic form; the squared norm of T(v)T(v) is equal to the desired test statistic at each location vv.

If saveall=TRUE, then a complete set of sufficient statistics is calculated and stored as an attribute of the result. This makes it possible to compute all of the statistics and pp values described above.

The function update.homtestmap, a method for the generic function update, converts an object of class "homtestmap" from one of these formats to another, where possible. Except in trivial cases, this requires that the "homtestmap" object was computed with saveall=TRUE.

The function homteststat computes the mean of the local test statistic or the mean of the local pp-values over the observation window.

To compute the pp-values when test="likelihood" or test="taylor", the values of the local likelihood ratio test statistic are referred to a gamma distribution whose first two moments are estimated from the data. If poolmoments=FALSE, the local estimates of the moments are used; if poolmoments=TRUE, the spatial average of these estimates is used. The default is to use pooling whenever it is theoretically justified, namely when the template model is a stationary point process.

Finer control over the computation is possible using the arguments ... passed to locppm.

Value

For homteststat, a numeric value giving the test statistic.

For homtestmap and update.homtestmap, a spatially-sampled function object (class "ssf"; see ssf). This object also belongs to the special class "homtestmap" which has a print method.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

homtest

Examples

example(locppm)
   plot(H <- homtestmap(fit))
   H

Locally Fitted Cluster or Cox Point Process Model

Description

Fits a Neyman-Scott cluster process or Cox point process model using a locally-weighted composite likelihood.

Usage

loccit(X, trend = ~1,
       clusters = c("Thomas", "MatClust", "Cauchy", "VarGamma", "LGCP"),
       covariates = NULL,
       ...,
       diagnostics = FALSE,
       taylor = FALSE,
       sigma = NULL, f = 1/4,
       clustargs = list(), control = list(), 
       rmax,
       covfunargs=NULL, use.gam=FALSE, nd=NULL, eps=NULL,
       niter=3,
       fftopt = list(),
       verbose = TRUE)

Arguments

X

Point pattern.

trend

Formula (without a left hand side) specifying the form of the logarithm of the intensity.

clusters

Character string determining the cluster model. Partially matched.

covariates

The values of any spatial covariates (other than the Cartesian coordinates) required by the model. A named list of pixel images, functions, windows or numeric constants.

diagnostics

Whether to perform auxiliary calculations in addition to the local estimates of the model parameters.

...

Additional arguments passed to as.mask to control the spatial resolution in the Fast Fourier Transform.

taylor

Logical value indicating whether to fit the model exactly at each spatial location (taylor=FALSE, the default) or to compute a first-order Taylor approximation to the fitted parameters (taylor=TRUE). The Taylor approximation is much faster.

sigma

Standard deviation of Gaussian kernel for local likelihood.

f

Argument passed to bw.frac to compute a value for sigma if it is missing or NULL.

clustargs

List of additional parameters for the cluster model, passed to the function RFcov in the RandomFields package.

control

List of control arguments passed to the generic optimisation function optim.

rmax

Maximum distance between pairs of points that will contribute to the composite likelihood.

covfunargs, use.gam, nd, eps

Arguments passed to ppm to control the intensity model and intensity fitting.

niter

Number of iterations in algorithm if taylor=FALSE.

fftopt

Developer use only.

verbose

Logical. If TRUE, print progress reports.

Details

This function fits a Cox or cluster process model to point pattern data locally, using the local Palm likelihood technique (Baddeley, 2016, section 8).

It can be used in the same way as kppm and effectively performs local fitting of the same model.

Value

An object of class "loccit".

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm

Examples

X <- redwood[owin(c(0,1), c(-1,-1/2))]
   fit <- loccit(X, ~1, "Thomas", nd=5, control=list(maxit=20))
   fit

Locally Fitted Cluster or Cox Point Process Model

Description

Fits a Neyman-Scott cluster process or Cox point process model using local minimum contrast.

Usage

locmincon(..., sigma = NULL, f = 1/4, verbose = TRUE,
                 localstatargs = list(), LocalStats = NULL,
                 tau = NULL)

Arguments

...

Arguments passed to kppm to determine the template homogeneous model.

sigma

Standard deviation of Gaussian kernel for local likelihood.

f

Argument passed to bw.frac to compute a value for sigma if it is missing or NULL.

verbose

Logical. If TRUE, print progress reports.

localstatargs

Optional. List of arguments to be passed to the local statistic localK, localKinhom, localpcf or localpcfinhom.

LocalStats

Optional. Values of the local statistics, if they have already been computed.

tau

Optional. Bandwidth for smoothing the fitted cluster parameters.

Details

The template or homogeneous model is first fitted by kppm. The statistic used to fit the template model is determined (as explained in the help for kppm) by the arguments statistic and trend.

The local version of this statistic is then computed. If statistic="K" and trend=~1 for example, the template model is fitted using the KK function Kest, and the local version is the local KK function localK. The possibilities are:

statistic stationary? template local
"K" yes Kest localK
"K" no Kinhom localKinhom
"pcf" yes pcf localpcf
"pcf" no pcfinhom localpcfinhom

These local functions, one for each data point, are then spatially averaged, using a Gaussian kernel with standard deviation sigma. Finally the model is fitted to each of the averaged local functions to obtain a local fit at each data point.

Value

Object of class "locmincon".

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

loccit

Examples

X <- redwood[owin(c(0,1), c(-1,-1/2))]
   fit <- locmincon(X, ~1, "Thomas", sigma=0.07)
   fit

Locally Fitted Poisson or Gibbs Point Process Model

Description

Fits Poisson or Gibbs point process model using local likelihood or pseudolikelihood.

Usage

locppm(..., sigma = NULL, f = 1/4,
       vcalc = c("none", "t", "hessian", "hom", "lik", "full"),
       locations=c("split", "fine", "coarse"),
       ngrid = NULL, grideps = NULL, verbose = TRUE,
       use.fft=FALSE, fft.algorithm="closepairs")

Arguments

...

Arguments passed to ppm to fit the homogeneous model.

sigma

Standard deviation of Gaussian kernel for local likelihood.

f

Argument passed to bw.frac to compute a value for sigma if it is missing or NULL.

vcalc

Type of variance calculation to be performed. See Details.

locations

Spatial locations for local calculations. See Details.

ngrid

Dimensions of coarse grid, if used. See Details. Incompatible with grideps.

grideps

Grid spacing of coarse grid, if used. See Details. Incompatible with ngrid.

verbose

Logical. If TRUE, print progress reports.

use.fft

Logical value indicating whether to perform computations using the Fast Fourier Transform. With use.fft = TRUE the code runs much faster but some quantities are not computed exactly. See Details.

fft.algorithm

Developer use only.

Details

This function fits a Poisson or Gibbs point process model to point pattern data by local likelihood or local pseudolikelihood respectively.

This command should be used in the same way as ppm. The point pattern data and the specification of the model are given in the leading arguments ... which are passed directly to ppm.

In all cases, the local estimates of the coefficients are computed. However, because the variance calculations are time-consuming, the default is not to perform them. This is controlled by the argument vcalc.

vcalc = "none":

no variance calculations are performed.

vcalc = "t":

the tt statistic for each parameter is computed for the local model.

vcalc = "hessian":

the local Hessian matrix is computed, and its negative inverse is used as a surrogate for the local variance.

vcalc = "hom":

No local fitting is performed. Calculations are performed only for the homogeneous (template) model. The variance of the local parameter estimates under the homogeneous model is computed.

vcalc = "lik":

In addition to the calculations for vcalc="hom" described above, if use.fft=FALSE the algorithm also computes the local composite likelihood ratio test statistic for the test of homogeneity. If use.fft=TRUE then vcalc="lik" is equivalent to vcalc="hom".

vcalc = "full":

all variance calculations are performed for the local model.

The spatial locations, where the model fits and variance calculations are performed, are determined by the argument locations.

locations = "fine":

The calculations are performed at every quadrature point of the model. This can take a very long time.

locations = "coarse":

The calculations are performed at the points of a coarse grid with dimensions specified by ngrid or grideps.

locations = "split":

The fitted coefficients are computed at every quadrature point of the model, but the variance calculations (if any) are performed at a coarse grid of locations, specified by ngrid or grideps. If neither ngrid nor grideps is specified, the default is ngrid=10.

If use.fft=FALSE (the default), all desired quantities are computed exactly, by an iterative algorithm that fits a separate model at each spatial location. This can be quite slow.

If use.fft=TRUE, we only compute quantities that can be obtained using the Fast Fourier Transform, resulting in much faster calculations (sometimes 3 orders of magnitude faster) when locations="fine". Properties of the homogeneous model are computed accurately. Properties of the locally-fitted model are approximated by a first order Taylor expansion.

Value

An object of class "locppm" representing the fitted model.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

methods.locppm, plot.locppm

Examples

fit <- locppm(swedishpines, ~1, sigma=9, nd=20)
   fit

Methods for Local Cluster or Cox Models

Description

Methods for various generic functions, for the class "locmincon" of locally fitted cluster or Cox point process models.

Usage

## S3 method for class 'locmincon'
as.ppp(X, ...)

  ## S3 method for class 'locmincon'
print(x, ...)

Arguments

x, X

A locally-fitted Cox or cluster point process model (object of class "locmincon").

...

Additional arguments

Details

Objects of class "locmincon" represent locally fitted cluster or Cox point process models.

The functions documented here provided methods for this class, for the generic functions as.ppp and print.

Value

as.ppp returns an object of class "ppp".

print returns NULL.

Author(s)

Adrian Baddeley

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locmincon

Examples

example(locmincon)
  fit
  as.ppp(fit)

Methods for Local Gibbs Models

Description

Methods for various generic functions, for the class "locppm" of locally fitted Gibbs point process models.

Usage

## S3 method for class 'locppm'
as.interact(object)

  ## S3 method for class 'locppm'
as.ppm(object)

  ## S3 method for class 'locppm'
coef(object, ..., which = c("local", "homogeneous"))

  ## S3 method for class 'locppm'
confint(object, parm, level = 0.95, ..., which = c("local", "homogeneous"))

  ## S3 method for class 'locppm'
is.poisson(x)

  ## S3 method for class 'locppm'
print(x, ...)

Arguments

object, x

A locally-fitted Gibbs point process model (object of class "locppm").

...

Additional arguments passed to the default method (for confint.locppm) or ignored (by coef.locppm).

which

Character string determining whether to perform calculations for the local Gibbs model (which="local", the default) or the corresponding homogeneous Gibbs model (which="homogeneous").

parm

The parameter or parameters for which a confidence interval is desired. A character string or character vector matching the names of coef(object), or an index or index vector that can be applied to coef(object).

level

Confidence level: a number between 0 and 1.

Details

Objects of class "locppm" represent locally fitted Gibbs point process models.

The functions documented here provided methods for this class, for the generic functions as.interact, as.ppm, coef, confint, is.poisson and print.

For the coef and confint methods, the calculations can be performed either on the locally fitted model or on its homogeneous equivalent, by changing the argument which.

Value

as.interact returns an interaction structure (object of class "interact").

as.ppm returns a fitted Gibbs model (object of class "ppm").

coef and confint return a numeric vector if which="homogeneous" and an object of class "ssf" if which="local".

is.poisson returns a logical value.

print returns NULL.

Author(s)

Adrian Baddeley

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm

Examples

fit <- locppm(swedishpines, ~1, sigma=9, nd=20,
                 vcalc="full", locations="coarse")
   fit
   is.poisson(fit)
   coef(fit)
   coef(fit, which="homogeneous")
   confint(fit)
   confint(fit, which="homogeneous")
   as.ppm(fit)
   as.interact(fit)

Plot a Locally Fitted Cluster or Cox Point Process Model

Description

Plot an object of class "loccit" representing a locally-fitted cluster or Cox point process model.

Usage

## S3 method for class 'loccit'
plot(x, ...,
               what = c("modelpar", "coefs", "lambda"),
               how = c("smoothed", "exact"), which = NULL,
               pre=NULL, post=NULL)

Arguments

x

The model to be plotted. A locally-fitted cluster or Cox point process model (object of class "loccit").

...

Arguments passed to plot.ppp or plot.im to control the plot.

what

Character string determining which quantities to display: "modelpar" for the cluster model parameters, "coefs" for the trend coefficients, or "lambda" for the fitted intensity.

how

Character string determining whether to display the fitted parameter values at the data points (how="exact") or the smoothed fitted parameters as pixel images (how="smoothed").

which

Optional. Which component(s) of the vector-valued quantity to display. An index or index vector. Default is to plot all components.

pre, post

Transformations to apply before and after smoothing.

Details

This is a method for the generic command plot for the class "loccit".

The argument which, if present, specifies which fitted parameters are displayed. It may be any kind of index for a numeric vector.

The quantities are computed at irregularly-placed points. If how="exact" the exact computed values will be displayed as circles centred at the locations where they were computed. If how="smoothed" these values will be kernel-smoothed using Smooth.ppp and displayed as a pixel image.

Value

NULL.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

loccit, plot, plot.default

Examples

X <- redwood[owin(c(0,1), c(-1,-1/2))]
   fitc <- loccit(X, ~1, "Thomas", nd=5, control=list(maxit=20))
   plot(fitc, how="exact")  
   plot(fitc, how="smoothed")

Plot a Locally Fitted Cluster or Cox Point Process Model

Description

Plot an object of class "locmincon" representing a locally-fitted cluster or Cox point process model.

Usage

## S3 method for class 'locmincon'
plot(x, ...,
               how = c("exact", "smoothed"),
               which = NULL, sigma = NULL, do.points = TRUE)

Arguments

x

The model to be plotted. A locally-fitted cluster or Cox point process model (object of class "locmincon" or "loccit").

...

Arguments passed to plot.ppp or plot.im to control the plot.

how

Character string determining whether to display the fitted parameter values at the data points (how="exact") or the smoothed fitted parameters as pixel images (how="smoothed").

which

Optional. Which component(s) of the vector-valued quantity to display. An index or index vector. Default is to plot all components.

sigma

Numeric. Smoothing bandwidth to be used if how="smoothed".

do.points

Logical. Whether to display the original point data as well.

Details

This is a method for the generic command plot for the class "locmincon".

The argument which, if present, specifies which fitted parameters are displayed. It may be any kind of index for a numeric vector.

The quantities are computed at irregularly-placed points. If how="exact" the exact computed values will be displayed as circles centred at the locations where they were computed. If how="smoothed" these values will be kernel-smoothed using Smooth.ppp and displayed as a pixel image.

Value

NULL.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locmincon, methods.locmincon, plot, plot.default

Examples

X <- redwood[owin(c(0,1), c(-1,-1/2))]
   fitm <- locmincon(X, ~1, "Thomas", sigma=0.07)
   plot(fitm, how="smoothed")  
   plot(fitm, how="exact")

Plot a Locally Fitted Poisson or Gibbs Model

Description

Plot an object of class "locppm" representing a locally-fitted Poisson or Gibbs point process model.

Usage

## S3 method for class 'locppm'
plot(x, ..., what = "cg", which = NULL)

## S3 method for class 'locppm'
contour(x, ..., what = "cg", which = NULL)

Arguments

x

A locally-fitted Poisson or Gibbs point process model (object of class "locppm").

...

Arguments passed to plot.ssf to control the plot.

what

What quantity to display. A character string. The default is to display the fitted coefficient vectors.

which

Which component(s) of the vector-valued quantity to display. An index or index vector.

Details

These are methods for the generic commands plot and contour, for the class "locppm".

The argument what specifies what quantity will be displayed:

"cg" Fitted coefficients of local model
"vg" Local variance matrix for Gibbs model
"vh" Local variance matrix for homogeneous model
"tg" tt-statistics based on "coefs" and "vg"

Typically these quantities are vector-valued (matrices are converted to vectors). The argument which, if present, specifies which elements of the vector are displayed. It may be any kind of index for a numeric vector.

The plotting is performed by plot.ssf.

Value

NULL.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm, methods.locppm, plot, plot.default

Examples

fit <- locppm(swedishpines, ~1, sigma=9, nd=20,
             vcalc="hessian", locations="coarse")
   plot(fit)  
   plot(fit, what="Vg")

Prediction for Locally-Fitted Cox or Cluster Model

Description

Computes the fitted intensity of a locally-fitted Cox process or cluster process model.

Usage

## S3 method for class 'loccit'
predict(object, ...)

  ## S3 method for class 'loccit'
fitted(object, ..., new.coef=NULL)

Arguments

object

Locally fitted point process model (object of class "loccit" fitted by loccit).

...

Arguments passed to predict.locppm.

new.coef

New values for the fitted coefficients. A matrix in which each row gives the fitted coefficients at one of the quadrature points of the model.

Details

The fitted intensity is computed.

Value

An object of class "ssf" as described in ssf.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

loccit, predict.locppm.

Examples

X <- redwood[owin(c(0,1), c(-1,-1/2))]
  fit <- loccit(X, ~1, "Thomas", nd=5, control=list(maxit=20))
  lam <- predict(fit)

Prediction of a Locally Fitted Poisson or Gibbs Point Process Model

Description

Computes the fitted intensity of a locally-fitted Poisson point process model, or the fitted intensity, trend or conditional intensity of a locally-fitted Gibbs point process model.

Usage

## S3 method for class 'locppm'
fitted(object, ...,
                        type = c("cif", "trend", "intensity"),
                        new.coef=NULL)

## S3 method for class 'locppm'
predict(object, ...,
                         type = c("cif", "trend", "intensity"),
                         locations=NULL, new.coef=NULL)

Arguments

object

A locally-fitted Poisson or Gibbs point process model (object of class "locppm").

...

Currently ignored.

new.coef

New vector or matrix of values for the model coefficients.

locations

Point pattern of locations where prediction should be computed.

type

Character string (partially matched) specifying the type of predicted value: the conditional intensity "cif" (the default), or the first order trend, or the intensity. For Poisson models all three options are equivalent.

Details

These are methods for the generic functions fitted and predict for the class "locppm" of locally-fitted Gibbs point process models.

The fitted method computes, for each quadrature point v (or in general, at each point v where a local model was fitted), the intensity of the locally-fitted model at v. The result is a numeric vector.

The predict computes the fitted intensity at any specified set of locations, and returns the result as an ssf object.

Value

For fitted.locppm, a numeric vector.

For predict.locppm, an object of class "ssf" as described in ssf.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm

Examples

fit <- locppm(cells, sigma=0.1, use.fft=TRUE)
  lam <- predict(fit)

Sibling Probability of Locally Fitted Cluster Point Process

Description

Computes the sibling probability of a locally fitted cluster point process model.

Usage

## S3 method for class 'loccit'
psib(object)

  ## S3 method for class 'locmincon'
psib(object)

Arguments

object

Fitted cluster point process model (object of class "loccit" or "locmincon").

Details

In a Poisson cluster process, two points are called siblings if they belong to the same cluster, that is, if they had the same parent point. If two points of the process are separated by a distance rr, the probability that they are siblings is p(r)=11/g(r)p(r) = 1 - 1/g(r) where gg is the pair correlation function of the process.

The value p(0)=11/g(0)p(0) = 1 - 1/g(0) is the probability that, if two points of the process are situated very close to each other, they came from the same cluster. This probability is an index of the strength of clustering, with high values suggesting strong clustering.

This concept was proposed in Baddeley, Rubak and Turner (2015, page 479) and Baddeley (2016).

The function psib is generic, with methods for "kppm", "loccit" and "locmincon".

The functions described here are the methods for locally-fitted cluster models of class "loccit" and "locmincon". They compute the spatially-varying sibling probability of the locally-fitted model.

Value

A spatially sampled function (object of class "ssf") giving the spatially-varying sibling probability.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

psib, kppm, loccit

Examples

## Not run: 
   fit <- loccit(redwood, ~1, "Thomas")
   
## End(Not run)
   
   fit
   plot(psib(fit))

Smooth a Locally Fitted Cluster or Cox Point Process Model

Description

Applies kernel smoothing to the fitted cluster parameters of a locally-fitted cluster or Cox point process model.

Usage

## S3 method for class 'locmincon'
Smooth(X, tau = NULL, ...)

Arguments

X

Object of class "locmincon".

tau

Smoothing bandwidth.

...

Additional arguments passed to Smooth.ppp controlling the smoothing and the pixel resolution.

Details

An object of class "locmincon" represents a locally-fitted Cox or cluster point process model. It provides estimates of the cluster parameters at each of the data points of the original point pattern dataset.

The parameter estimates will be smoothed using a Gaussian kernel with standard deviation tau.

Value

A pixel image or a list of pixel images.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locmincon, Smooth.ppp

Examples

fit <- locmincon(redwood)
   Smooth(fit, tau=0.1)

Smooth a locally fitted Gibbs model

Description

Applies kernel smoothing to one of the components of a locally-fitted Gibbs point process model.

Usage

## S3 method for class 'locppm'
Smooth(X, ..., what = "cg")

Arguments

X

A locally-fitted Gibbs point process model (object of class "locppm").

...

Arguments passed to Smooth.ppp to control the smoothing.

what

Component to be smoothed. A character string. The default is to smooth the fitted coefficient vectors.

Details

This function extracts the selected quantity from the fitted object and spatially smooths it using Smooth.ppp. The result is a pixel image or a list of pixel images.

Value

A pixel image or a list of pixel images.

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm, Smooth.ppp

Examples

fit <- locppm(cells, sigma=0.1, use.fft=TRUE)
  plot(Smooth(fit))

Test of Effect in Locally Fitted Point Process Model

Description

Perform a local tt-test for the presence of a covariate effect in a locally fitted Poisson or Gibbs point process model.

Usage

ttestmap(object, term, ...,
         method = c("exact", "hessian", "taylor"),
         grid = FALSE, 
         ngrid = NULL, grideps = NULL,
         verbose = TRUE)

Arguments

object

Locally fitted Poisson or Gibbs point process model (object of class "locppm").

term

Term to be dropped from the model. A character string matching a term in the model formula

...

Ignored.

method

Choice of method to be used to evaluate the tt statistic. See Details.

grid

Logical. If FALSE, compute the test statistic at all quadrature points of the model. If TRUE, compute at a coarse grid of locations.

ngrid

Number of grid points (in each axis direction) for the coarse grid. Incompatible with grideps.

grideps

Spacing (horizontal and vertical) between grid points for the coarse grid. Incompatible with ngrid.

verbose

Logical value indicating whether to print progress reports.

Details

The argument object should be a locally-fitted Poisson or Gibbs point process model (object of class "locppm" created by locppm).

This function computes the local tt test statistic for the test that a particular covariate effect in the model is zero. This is described in Baddeley (2016, sections 3 and 5).

Value

Object of class "ssf".

Author(s)

Adrian Baddeley [email protected].

References

Baddeley, A. (2017) Local composite likelihood for spatial point patterns. Spatial Statistics 22, 261–295. DOI: 10.1016/j.spasta.2017.03.001

Baddeley, A., Rubak, E. and Turner, R. (2015) Spatial Point Patterns: Methodology and Applications with R. Chapman and Hall/CRC Press.

See Also

locppm

Examples

fit <- with(copper,
        locppm(Points, ~D, covariates=list(D=distfun(Lines)), nd=c(7,15)))
 plot(ttestmap(fit, "D"))

Evaluate an Expression for a Locally Fitted Model

Description

Given a locally-fitted Cox or cluster point process model, evaluate an expression involving the fitted cluster parameters.

Usage

## S3 method for class 'locmincon'
with(data, ...)

  ## S3 method for class 'loccit'
with(data, ...)

Arguments

data

An object of class "locmincon" or "loccit" representing a locally-fitted Cox or cluster point process model.

...

Arguments passed to with.default specifying the expression to be evaluated.

Details

These are method for the generic function with for the classes "locmincon" and "loccit".

An object of class "locmincon" or "loccit" represents a locally-fitted Cox or cluster point process model. It contains a data frame which provides estimates of the cluster parameters at each of the data points of the original point pattern dataset.

The expression specified by ... will be evaluated in this dataframe. If the result of evaluation is a data frame with one row for each data point, or a numeric vector with one entry for each data point, then the result will be an object of class "ssf" containing this information. Otherwise, the result will be a numeric vector.

Value

An object of class "ssf" or a numeric vector.

Author(s)

Adrian Baddeley [email protected].

See Also

ssf

Examples

example(locmincon)
  with(fit, kappa * sigma2)
  example(locmincon)
  with(fit, kappa * sigma2)