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 object 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]
  )
  rspde_model <- rspde.matern(
    mesh = mesh_2d,
    nu_upper_bound = 2
  )

  cmp <- y ~ Intercept(1) +
    field(cbind(x1,x2), model = rspde_model)


  rspde_fit <- bru(cmp, data = data_df)
  summary(rspde_fit)
}
#> This is INLA_24.12.11 built 2024-12-11 19:58:26 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.12.0
#> INLA version: 24.12.11
#> Components:
#> Intercept: main = linear(1), group = exchangeable(1L), replicate = iid(1L), NULL
#> field: main = cgeneric(cbind(x1, x2)), group = exchangeable(1L), replicate = iid(1L), NULL
#> Likelihoods:
#>   Family: 'gaussian'
#>     Tag: ''
#>     Data class: 'data.frame'
#>     Response class: 'numeric'
#>     Predictor: y ~ .
#>     Used components: effects[Intercept, field], latent[]
#> Time used:
#>     Pre = 0.41, Running = 6.84, Post = 0.382, Total = 7.64 
#> Fixed effects:
#>            mean    sd 0.025quant 0.5quant 0.975quant  mode kld
#> Intercept 0.286 0.172     -0.054    0.285      0.629 0.285   0
#> 
#> Random effects:
#>   Name	  Model
#>     field CGeneric
#> 
#> Model hyperparameters:
#>                                            mean     sd 0.025quant 0.5quant
#> Precision for the Gaussian observations 168.733 90.699     49.488  149.940
#> Theta1 for field                         -3.854  0.511     -4.989   -3.807
#> Theta2 for field                          2.881  0.216      2.460    2.879
#> Theta3 for field                         -0.123  0.277     -0.589   -0.145
#>                                         0.975quant    mode
#> Precision for the Gaussian observations    395.673 115.211
#> Theta1 for field                            -3.014  -3.577
#> Theta2 for field                             3.311   2.873
#> Theta3 for field                             0.488  -0.257
#> 
#> Deviance Information Criterion (DIC) ...............: -125.22
#> Deviance Information Criterion (DIC, saturated) ....: 196.31
#> Effective number of parameters .....................: 93.72
#> 
#> Watanabe-Akaike information criterion (WAIC) ...: -145.94
#> Effective number of parameters .................: 54.28
#> 
#> Marginal log-Likelihood:  -109.52 
#>  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
# }