Skip to contents

rSPDE inlabru mapper

Usage

bru_get_mapper.inla_rspde(model, ...)

ibm_n.bru_mapper_inla_rspde(mapper, ...)

ibm_values.bru_mapper_inla_rspde(mapper, ...)

ibm_jacobian.bru_mapper_inla_rspde(mapper, input, ...)

Arguments

model

An inla_rspde for which to construct or extract a mapper

...

Arguments passed on to other methods

mapper

A bru_mapper_inla_rspde object

input

The values for which to produce a mapping matrix

Examples

# \donttest{
# devel version
if (requireNamespace("INLA", quietly = TRUE) &&
  requireNamespace("inlabru", quietly = TRUE)) {
  library(INLA)
  library(inlabru)

  set.seed(123)
  m <- 100
  loc_2d_mesh <- matrix(runif(m * 2), m, 2)
  mesh_2d <- inla.mesh.2d(
    loc = loc_2d_mesh,
    cutoff = 0.05,
    max.edge = c(0.1, 0.5)
  )
  sigma <- 1
  range <- 0.2
  nu <- 0.8
  kappa <- sqrt(8 * nu) / range
  op <- matern.operators(
    mesh = mesh_2d, nu = nu,
    range = range, sigma = sigma, m = 2,
    parameterization = "matern"
  )
  u <- simulate(op)
  A <- inla.spde.make.A(
    mesh = mesh_2d,
    loc = loc_2d_mesh
  )
  sigma.e <- 0.1
  y <- A %*% u + rnorm(m) * sigma.e
  y <- as.vector(y)

  data_df <- data.frame(
    y = y, x1 = loc_2d_mesh[, 1],
    x2 = loc_2d_mesh[, 2]
  )
  coordinates(data_df) <- c("x1", "x2")
  rspde_model <- rspde.matern(
    mesh = mesh_2d,
    nu_upper_bound = 2
  )

  cmp <- y ~ Intercept(1) +
    field(coordinates, model = rspde_model)


  rspde_fit <- bru(cmp, data = data_df)
  summary(rspde_fit)
}
#> Loading required package: sp
#> This is INLA_24.04.25-1 built 2024-04-25 17:05:50 UTC.
#>  - See www.r-inla.org/contact-us for how to get help.
#>  - List available models/likelihoods/etc with inla.list.models()
#>  - Use inla.doc(<NAME>) to access documentation
#> Loading required package: fmesher
#> inlabru version: 2.10.1
#> INLA version: 24.04.25-1
#> Components:
#> Intercept: main = linear(1), group = exchangeable(1L), replicate = iid(1L)
#> field: main = cgeneric(coordinates), group = exchangeable(1L), replicate = iid(1L)
#> Likelihoods:
#>   Family: 'gaussian'
#>     Data class: 'SpatialPointsDataFrame'
#>     Predictor: y ~ .
#> Time used:
#>     Pre = 0.437, Running = 7.08, Post = 0.669, Total = 8.18 
#> Fixed effects:
#>            mean   sd 0.025quant 0.5quant 0.975quant  mode kld
#> Intercept 0.286 0.18     -0.072    0.286      0.648 0.286   0
#> 
#> Random effects:
#>   Name	  Model
#>     field CGeneric
#> 
#> Model hyperparameters:
#>                                           mean     sd 0.025quant 0.5quant
#> Precision for the Gaussian observations 168.95 87.633      50.30   151.65
#> Theta1 for field                         -3.29  0.465      -4.20    -3.29
#> Theta2 for field                          2.81  0.246       2.34     2.81
#> Theta3 for field                         -1.38  0.194      -1.76    -1.38
#>                                         0.975quant   mode
#> Precision for the Gaussian observations     385.43 118.47
#> Theta1 for field                             -2.37  -3.31
#> Theta2 for field                              3.31   2.79
#> Theta3 for field                             -1.00  -1.37
#> 
#> Deviance Information Criterion (DIC) ...............: -123.05
#> Deviance Information Criterion (DIC, saturated) ....: 199.02
#> Effective number of parameters .....................: 96.25
#> 
#> Watanabe-Akaike information criterion (WAIC) ...: -144.90
#> Effective number of parameters .................: 55.61
#> 
#> Marginal log-Likelihood:  -109.91 
#>  is computed 
#> Posterior summaries for the linear predictor and the fitted values are computed
#> (Posterior marginals needs also 'control.compute=list(return.marginals.predictor=TRUE)')
#> 
# devel.tag
# }