LKrigSetup: Create or update the LatticeKrig model object (LKinfo) for...

Description Usage Arguments Details Value Author(s) Examples

View source: R/LKrigSetup.R

Description

This function combines some input arguments with defaults for other to create a list describing the LatticeKrig spatial model. A key to specifying the LatticeKrig spatial model is specifying the geometry, e.g. LKRectangle for a 2-d rectangular domain. Each geometry has some parameters that control the basic model setup and these are included through the ... arguments of this function. See the help for this argument below for some examples.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
LKrigSetup(x = NULL, nlevel = NULL, alpha = NA, alphaObject =
                 NULL, nu = NULL, a.wght = NA, a.wghtObject = NULL, NC
                 = NULL, NC.buffer = NULL, normalize = TRUE, lambda =
                 NA, sigma = NA, rho = NA, rho.object = NULL,
                 latticeInfo = NULL, basisInfo = NULL, LKGeometry =
                 "LKRectangle", distance.type = "Euclidean",
                 BasisFunction = "WendlandFunction", overlap = 2.5, V =
                 NULL, BasisType = "Radial", fixedFunction =
                 "LKrigDefaultFixedFunction", fixedFunctionArgs =
                 list(m = 2), collapseFixedEffect = FALSE, max.points =
                 NULL, mean.neighbor = 50, choleskyMemory = NULL,
                 verbose = FALSE, noCheck = FALSE, returnCall = FALSE,
                 dense = FALSE, ...)
                 
LatticeKrigEasyDefaults(argList, nlevel, x)

LKinfoUpdate(LKinfo, ... )

Arguments

argList

Argument supplied to the top level LatticeKrig function.

alpha

A vector of length nlevel with the relative variances for the different multi-resolution levels.

alphaObject

For non-stationary models an object to be used with the predict function to give the alpha values at the process locations. Typically this is a list of "predict" objects. See nonstationaryModels for examples how to use this option.

a.wght

This parameter controls the correlation range in the SAR model. In most cases, a scalar value and for the 2-d model by default greater than 4. If a vector this can specify an anisotropic set of weights. To specify a.wght parameters that are different for each level they should be in the form of a list of length nlevel. E.g. a.wght = list( 4.5, 5, 10) will specify different a.wghts for 3 different levels. The setup function will check that the values are in a valid range for a geometry and the length of the list agrees with the number of levels.

The details of how this is connected to the covariance function varies based on the geometry. However, qualitatively this is related to a range parameter. For the LKRectangle geometry and a stationary model, at level k the center point has weight 1 with the 4 nearest neighbors given weight -1/a.wght[k]. In this case a.wght must be greater than 4 for the fields to be stationary and following Lindgren and Rue the range parameter is approximately 1/sqrt(a.wght-4).

a.wghtObject

For non-stationary models an object to be used with the predict function to give the a.wght values at the lattice locations. See nonstationaryModels for examples how to used option.

basisInfo

A list with extra components the object used to describe the multi-resolution basis. Usually this will not be needed for standard models.

BasisType

A character string indicating the type of basis function. Currently this is either "Radial" or "Tensor".

choleskyMemory

A list that will be used in the spam call to do the Cholesky decomposition. See the memory argument in chol.spam.

collapseFixedEffect

If FALSE the fixed part of the model is found separately for each replicated data set. If TRUE the estimate is polled across replicates.This is largely a modeling decision whether variation among the replicate fields is due to the spatial component or also include variation in the fixed effects across replicates – guess they are not really fixed then!

dense

If FALSE sparse linear algebra is used for the computations . If TRUE the matrices are made "dense" (zeroes are filled in) and the ordinary Lapack functions are used for the linear algebra. This option is primarily for testing and timing sparse verses standard algorithms.

distance.type

A text string indicate type distance to use between spatial locations when evaluating the basis functions. Default is "Euclidean". Other choices when locations are in degrees longitude and latitude are "Chordal" and "GreatCircle". The default radius is in miles (3963.34) but can be reset using the "Radius" attribute to the text string.

e.g. dtype<- "Chordal"; attr( dtype, which="Radius") <- 6371

will set radius to kilometers.

fixedFunction

A text string that is the name of the function used to find the fixed part of the spatial model based on the locations. The default is a linear (m=2) polynomial in the spatial coordinates. See LKrigDefaultFixedFunction for more details. Set this to NULL if you do not want to include a fixed effect in the spatial model.

fixedFunctionArgs

A list containing arguments to supply when evaluating the fixed function.

lambda

The "noise to signal ratio" or also known as the smoothing parameter it is the parameter lambda = sigma^2/rho. If specified then sigma and rho typically are estimated in LKrig by maximum likelihood. If lambda is not specified then it is set as lambda = sigma^2/ rho. Note that to evaluate the spatial process model, e.g. using the function LKrig.cov, a value of lambda is not needed and this argument can default to NA.

latticeInfo

Part or all of the object used to describe the Markov random field lattice. In the standard cases this list is created in the setup function and need not be specified. See LKrigSetupLattice for details. Note that the contents of this list is concatenated to any additional components supplied by LKrigSetupLattice.

LKGeometry

A text string that gives the names of the model geometry. The default is "LKrectangle" assuming the spatial domain is a rectangle. Other common choices are "LKInterval" (1 d problem ) and "LKBox" (3 d problem)

LKinfo

A list that has class "LKinfo".

mean.neighbor

The average number of nonzero points when each basis function is evaluated at a set of points points in the spatial domain.

max.points

This is a parameter for the nearest neighbor distance function that sets the maximum array size for the nonzero distances among points. e.g. with 100 points each with 20 nonzero neighbors max.points needs to be 2000 = 100*20. Specifically if the total number of nonzero values when each basis function is evaluated at all the spatial locations. The easier way to specify space is by using the mean.neighbor argument.

NC

For regular grids of lattice points the maximum number of lattice grid points for a spatial coordinate and at the coarsest level of resolution. For a example, for a square region, (and LKGeometry = "LKRectangle") NC=5 results in a 5X5 = 25 lattice points at the first level. Note that the default is that lattice points are also taken to have the same spacing in every dimension.

NC.buffer

Number of extra lattice points added outside the spatial domain for regular grids of lattice points. This helps to reduce boundary effects from the SAR model. NC.buffer=5 and NC=5 for a square region will result in (5+ 2*5) X (5+ 2*5) = 225 lattice locations at the coarsest level of resolution. Note that this number by default is fixed for finer resolutions and so does not contribute as much to the total number of lattice points.

nlevel

Number of levels in multi-resolution. Note that each subsequent level increases the number of basis functions within the spatial domain size by a factor of roughly 4.

noCheck

If FALSE do not make any checks on the consistency of the different parts of the final LKinfo object. e.q. values of a.wght within the range of a stationary Markov random field.

normalize

If TRUE the basis functions will be normalized to give a marginal variance of one for each level of multi-resolution. (Normalizing by levels makes it easier to interpret the alpha weights.)

nu

A smoothness parameter that controls relative sizes of alpha. Currently this parameter only makes sense for the 2D rectangular domain (LKRectangle)

overlap

Controls the overlap among the radial basis functions and should be in units of the lattice spacing. For the rectangular case the default of 2.5 means that the support of the Wendland basis function will overlap 2.5 lattice units in each direction. See LKrig.basis for how overlap adjusts scaling in the basis function formula.

rho

A scalar, the sill or marginal variance of the process.

BasisFunction

Text string giving the 1-d form for basis function.

rho.object

A prediction object to specify part of the marginal variance for the process. Specifically the form is VAR(g(x1))= rho*h(x1) Calling predict(rho.object,x1) should return a vector with the values of h at the (arbitrary) spatial locations in x1. If omitted assumed to be the constant one – the usual case.

returnCall

If TRUE the call to LKrigSetup is also included as part of the LKinfo object.

sigma

The measurement error standard deviation. Also called the nugget variance in geostatistics.

V

See entry in LKrig.

verbose

If TRUE print out intermediate information.

x

Spatial locations that define the spatial domain for prediction. This is only used to determine ranges of the grid for the basis functions so, for example, for a rectangular domain only two points are required that bound the rest of the data locations. E.g. x= rbind( c( 0,0), c(1,1)) will set the domain to be the unit square.

...

Specific arguments that will be included in the setupArgs list and also passed to LKrigSetupLattice. For LKinfoUpdate these specify the components of LKinfo to update. For example for a rectangular domain (the default) the argument NC is needed to specify the initial lattice size NC.buffer sets the number of extra points included in the lattice and defaults to 5. NC for the LKSphere geometry is the initial level of the icosahedral grid.

Details

Many of the functions within LKrigSetup are overloaded to adapt to the LKGeometry class. This makes it easy to add new geometries or other models to the LatticeKrig framework. The required components of this object (see below) outline how the latticeKrig model is structured and what should be common features independent of the geometry. The key components are latticeInfo that contains the information used to generate the spatial autoregressive matrix on the lattice (see LKrigSAR) and the lattice centers (see LKrigLatticeCenters). The component basisInfo used to generate the radial basis function (see LKrig.basis)

The function LKrigEasyDefaults is used in the top level function LattticeKrig to make the logic of different default choices more readable and reduces the clutter in this function. Its main purpose is to find a reasonable choice for NC when this is not specified.

The function LKinfoUpdate is more of a utility used for clarity that allows one to update the LKinfo object with particular components without having to recreate the entire object. This function is used in the MLE search when just values of alpha or a.wght are being varied.

Value

An object with class "LKinfo" and also the additional class given by LKGeometry. The required components are:

nlevel

Number of levels

alpha

alpha parameters as a list that has nlevel components and possibly some attributes.

a.wght

a.wght parameters as a list that has nlevel components and possibly some attributes.

nu

nu parameter

normalize

A logical indicating whether to normalize.

lambda

Value of lambda.

sigma

Value of sigma.

rho

Value for rho.

rho.object

Value for rho.object.

latticeInfo

A list with specific multi-resolution lattice information

setupArgs

All arguments passed in the call and any in in ...

basisInfo

A list with basis information.

call

The actual call used to create this object.

Author(s)

Doug Nychka

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
  data(ozone2)
  # find the ranges of the  data, this is the same as passing
  # the entire set of observation locations and is more compact 
  sDomain<-apply( ozone2$lon.lat, 2,"range")
  LKinfo0<- LKrigSetup( sDomain, NC=10, nlevel=2, alpha=c(1,.5),
                       a.wght = 5)
  print( LKinfo0)
  
  #Gigantic buffer added note extra basis functions. 
  LKinfo<- LKrigSetup( sDomain, NC=10, NC.buffer= 15, nlevel=2, 
  alpha=c(1,.5),a.wght = 5)
  print( LKinfo)
  
  LKinfo2<- LKinfoUpdate( LKinfo,  a.wght=4.1, NC=12)
  LKinfo3<- LKrigSetup( sDomain, NC=12, nlevel=2, alpha=c(1,.5),
                        a.wght=4.1)
# LKinfo2 and LKinfo3 should be the same.                         
 

Example output

Loading required package: spam
Loading required package: dotCall64
Loading required package: grid
Spam version 2.2-2 (2019-03-07) is loaded.
Type 'help( Spam)' or 'demo( spam)' for a short introduction 
and overview of this package.
Help for individual functions is also obtained by adding the
suffix '.spam' to the function name, e.g. 'help( chol.spam)'.

Attaching package: 'spam'

The following objects are masked from 'package:base':

    backsolve, forwardsolve

Loading required package: fields
Loading required package: maps
See https://github.com/NCAR/Fields for
 an extensive vignette, other supplements and source code 
Classes for this object are:  LKinfo LKRectangle
The second class usually will indicate the geometry
     e.g.  2-d rectangle is  LKRectangle
 
Some details on spatial autoregression flags:
stationary:  TRUE TRUE
first order (by level):  TRUE TRUE
isotropic:  TRUE TRUE
 
Ranges of locations in raw scale:
        [,1]   [,2]
[1,] -93.572 36.791
[2,] -82.960 44.453
 
Logical (collapseFixedEffect) if fixed effects will be pooled: FALSE
 
Number of levels: 2
delta scalings: 1.179111 0.5895556
with an overlap parameter of  2.5
alpha:  1 0.5
 
a.wght:  5 5
 
Basis  type: Radial using  WendlandFunction  and Euclidean  distance.
Basis functions will be normalized
 
Total number of basis functions  1007
 Level Basis size      
     1        340 20 17
     2        667 29 23
 
Lambda value:  NA
Classes for this object are:  LKinfo LKRectangle
The second class usually will indicate the geometry
     e.g.  2-d rectangle is  LKRectangle
 
Some details on spatial autoregression flags:
stationary:  TRUE TRUE
first order (by level):  TRUE TRUE
isotropic:  TRUE TRUE
 
Ranges of locations in raw scale:
        [,1]   [,2]
[1,] -93.572 36.791
[2,] -82.960 44.453
 
Logical (collapseFixedEffect) if fixed effects will be pooled: FALSE
 
Number of levels: 2
delta scalings: 1.179111 0.5895556
with an overlap parameter of  2.5
alpha:  1 0.5
 
a.wght:  5 5
 
Basis  type: Radial using  WendlandFunction  and Euclidean  distance.
Basis functions will be normalized
 
Total number of basis functions  3587
 Level Basis size      
     1       1480 40 37
     2       2107 49 43
 
Lambda value:  NA

LatticeKrig documentation built on Nov. 9, 2019, 5:07 p.m.