--- title: "Intrinsic models in the rSPDE package" author: "David Bolin" date: "`r Sys.Date()`" output: rmarkdown::html_vignette vignette: > %\VignetteIndexEntry{Intrinsic models in the rSPDE package} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} bibliography: references.bib --- ```{r setup, include = FALSE} knitr::opts_chunk$set( collapse = TRUE, comment = "#>" ) set.seed(1) ``` ```{r inla_link, include = FALSE} inla_link <- function() { sprintf("[%s](%s)", "`R-INLA`", "https://www.r-inla.org") } ``` ## Introduction In this vignette we provide a brief introduction to the intrinsic models implemented in the `rSPDE` package. ## A fractional intrinsic model A basic intrinsic model which is implemented in `rSPDE` is defined as $$ (-\Delta)^{\beta/2}(\tau u) = \mathcal{W}, $$ where $\beta > d/2$ and $d$ is the dimension of the spatial domain. To illustrate these models, we begin by defining a mesh over $[0,2]\times [0, 2]$: ```{r, message=FALSE, warning=FALSE, fig.align='center'} library(fmesher) bnd <- fm_segm(rbind(c(0, 0), c(2, 0), c(2, 2), c(0, 2)), is.bnd = TRUE) mesh_2d <- fm_mesh_2d( boundary = bnd, cutoff = 0.02, max.edge = c(0.1) ) plot(mesh_2d, main = "") ``` We now use the `intrinsic.operators()` function to construct the `rSPDE` representation of the general model. ```{r, message=FALSE} library(rSPDE) tau <- 0.2 beta <- 1.8 fem <- fm_fem(mesh_2d) op <- intrinsic.operators(tau = tau, beta = beta, mesh = mesh_2d, m = 2) ``` To see that the `rSPDE` model is approximating the true model, we can compare the variogram of the approximation (implemented in the function `variogram` in the model object) with the true variogram (implemented in `variogram.intrinsic.spde()`) as follows. ```{r, message=FALSE} point <- matrix(c(1,1),1,2) Gamma <- op$variogram(point) vario <- variogram.intrinsic.spde(point, mesh_2d$loc[,1:2], tau = tau, beta = beta, L = 2, d = 2) d = sqrt((mesh_2d$loc[,1]-point[1])^2 + (mesh_2d$loc[,2]-point[2])^2) plot(d, Gamma, xlim = c(0,0.7), ylim = c(0,3), ylab = "variogram(h)", xlab = "h") points(d,vario,col=2) ``` If we want to increase the accuracy, we can either use a finer mesh or increase the order of the rational approximation through the argument `m` in `intrinsic.operators`. The default value of `m` is 1. We can now use the `simulate` function to simulate a realization of the field $u$: ```{r} u <- simulate(op,nsim = 1) proj <- fm_evaluator(mesh_2d, dims = c(100, 100)) field <- fm_evaluate(proj, field = as.vector(u)) field.df <- data.frame(x1 = proj$lattice$loc[,1], x2 = proj$lattice$loc[,2], y = as.vector(field)) library(ggplot2) library(viridis) ggplot(field.df, aes(x = x1, y = x2, fill = y)) + geom_raster() + scale_fill_viridis() ``` By default, the field is simulated with a zero-integral constraint. ## Fitting the model with `R-INLA` Let us now consider a simple Gaussian linear model where the spatial field $u(\mathbf{s})$ is observed at $m$ locations, $\{\mathbf{s}_1 , \ldots , \mathbf{s}_m \}$ under Gaussian measurement noise. For each $i = 1,\ldots,m,$ we have $$ \begin{align} y_i &= u(\mathbf{s}_i)+\varepsilon_i\\ \end{align}, $$ where $\varepsilon_1,\ldots,\varepsilon_{m}$ are iid normally distributed with mean 0 and standard deviation 0.1. To generate a data set `y` from this model, we first draw some observation locations at random in the domain and then use the `spde.make.A()` functions (that wraps the functions `fm_basis()`, `fm_block()` and `fm_row_kron()` of the `fmesher` package) to construct the observation matrix which can be used to evaluate the simulated field $u$ at the observation locations. After this we simply add the measurment noise. ```{r} n_loc <- 1000 loc_2d_mesh <- matrix(2*runif(n_loc * 2), n_loc, 2) A <- spde.make.A( mesh = mesh_2d, loc = loc_2d_mesh ) sigma.e <- 0.1 y <- A %*% u + rnorm(n_loc) * sigma.e ``` The generated data can be seen in the following image. ```{r,fig.align = "center", echo=TRUE} df <- data.frame(x1 = as.double(loc_2d_mesh[, 1]), x2 = as.double(loc_2d_mesh[, 2]), y = as.double(y)) ggplot(df, aes(x = x1, y = x2, col = y)) + geom_point() + scale_color_viridis() ``` We will now fit the model using our `r inla_link()` implementation of the rational SPDE approach. Further details on this implementation can be found in [R-INLA implementation of the rational SPDE approach](rspde_inla.html). ```{r} library(INLA) rspde.order <- 2 mesh.index <- rspde.make.index(name = "field", mesh = mesh_2d, rspde.order = rspde.order) Abar <- rspde.make.A(mesh = mesh_2d, loc = loc_2d_mesh, rspde.order = rspde.order) st.dat <- inla.stack(data = list(y = as.vector(y)), A = Abar, effects = mesh.index) ``` We now create the model object. ```{r} rspde_model <- rspde.intrinsic(mesh = mesh_2d, rspde.order = rspde.order) ``` Finally, we create the formula and fit the model to the data: ```{r message=FALSE, warning=FALSE} f <- y ~ -1 + f(field, model = rspde_model) rspde_fit <- inla(f, data = inla.stack.data(st.dat), family = "gaussian", control.predictor = list(A = inla.stack.A(st.dat))) ``` To compare the estimated parameters to the true parameters, we can do the following: ```{r} result_fit <- rspde.result(rspde_fit, "field", rspde_model) summary(result_fit) tau <- op$tau nu <- op$beta - 1 #beta = nu + d/2 result_df <- data.frame( parameter = c("tau", "nu", "sigma.e"), true = c(tau, nu, sigma.e), mean = c(result_fit$summary.tau$mean,result_fit$summary.nu$mean, sqrt(1/rspde_fit$summary.hyperpar[1,1])), mode = c(result_fit$summary.tau$mode, result_fit$summary.nu$mode, sqrt(1/rspde_fit$summary.hyperpar[1,6])) ) print(result_df) ``` ## Extreme value models When used for extreme value statistics, one might want to use a particular form of the mean value of the latent field $u$, which is zero at one location $k$ and is given by the diagonal of $Q_{-k,-k}^{-1}$ for the remaining locations. This option can be specified via the `mean.correction` argument of `rspde.intrinsic`: ```{r} rspde_model2 <- rspde.intrinsic(mesh = mesh_2d, rspde.order = rspde.order, mean.correction = TRUE) ``` We can then fit this model as before: ```{r message=FALSE, warning=FALSE} f <- y ~ -1 + f(field, model = rspde_model2) rspde_fit <- inla(f, data = inla.stack.data(st.dat), family = "gaussian", control.predictor = list(A = inla.stack.A(st.dat))) ``` To see the posterior distributions of the parameters we can do: ```{r} result_fit <- rspde.result(rspde_fit, "field", rspde_model2) posterior_df_fit <- gg_df(result_fit) ggplot(posterior_df_fit) + geom_line(aes(x = x, y = y)) + facet_wrap(~parameter, scales = "free") + labs(y = "Density") ``` ## An example with replicates Let us redo the previous example with replicated data to illustrate that replicates are handled in the same way as any other `rSPDE` model. We start by generating some data with 200 observations per replicate ```{r} set.seed(1) tau <- 0.2 beta <- 1.9 op <- intrinsic.operators(tau = tau, beta = beta, mesh = mesh_2d) n.rep <- 5 m <- 1000 loc_2d_mesh <- matrix(2*runif(m * 2), m, 2) A <- spde.make.A( mesh = mesh_2d, loc = loc_2d_mesh, index = rep(1:m, times = n.rep), repl = rep(1:n.rep, each = m) ) u <- simulate(op, nsim = n.rep) y <- as.vector(A %*% as.vector(u)) + rnorm(m * n.rep) * 0.1 ``` We now create the stack, A matrix and index and fit the model: ```{r} Abar.rep <- rspde.make.A( mesh = mesh_2d, loc = loc_2d_mesh, index = rep(1:m, times = n.rep), repl = rep(1:n.rep, each = m) ) mesh.index.rep <- rspde.make.index( name = "field", mesh = mesh_2d, n.repl = n.rep ) st.dat.rep <- inla.stack( data = list(y = y), A = Abar.rep, effects = mesh.index.rep ) rspde_model.rep <- rspde.intrinsic(mesh = mesh_2d, prior.nu.dist = "beta") f.rep <- y ~ -1 + f(field, model = rspde_model.rep, replicate = field.repl ) rspde_fit.rep <- inla(f.rep, data = inla.stack.data(st.dat.rep), family = "gaussian", control.predictor = list(A = inla.stack.A(st.dat.rep)) ) ``` We then compare with the true parameter estimates as before ```{r} result_fit <- rspde.result(rspde_fit.rep, "field", rspde_model.rep) summary(result_fit) tau <- op$tau nu <- op$beta - 1 #beta = nu + d/2 result_df <- data.frame( parameter = c("tau", "nu", "sigma.e"), true = c(tau, nu, sigma.e), mean = c(result_fit$summary.tau$mean,result_fit$summary.nu$mean, sqrt(1/rspde_fit.rep$summary.hyperpar[1,1])), mode = c(result_fit$summary.tau$mode, result_fit$summary.nu$mode, sqrt(1/rspde_fit.rep$summary.hyperpar[1,6])) ) print(result_df) ``` To see the posterior distributions of the parameters we can do: ```{r} result_fit <- rspde.result(rspde_fit.rep, "field", rspde_model.rep) posterior_df_fit <- gg_df(result_fit) ggplot(posterior_df_fit) + geom_line(aes(x = x, y = y)) + facet_wrap(~parameter, scales = "free") + labs(y = "Density") ``` # A more general model The `rSPDE` package also contains a partial implementation of a more general intrinsic model, which we refer to as an intrinsic Matérn model. The model is defined as $$ (-\Delta)^{\beta/2}(\kappa^2-\Delta)^{\alpha/2}(\tau u) = \mathcal{W}, $$ where $\alpha + \beta > d/2$ and $d$ is the dimension of the spatial domain. These models are handled by performing two rational approximations, one for each fractional operator. To illustrate this model, we consider the same mesh as before and use the `intrinsic.matern.operators()` function to construct the `rSPDE` representation of the general model. ```{r, message=FALSE} bnd <- fm_segm(rbind(c(0, 0), c(2, 0), c(2, 2), c(0, 2)), is.bnd = TRUE) mesh_2d <- fm_mesh_2d( boundary = bnd, cutoff = 0.01, max.edge = c(0.05) ) kappa <- 10 tau <- 0.0025 alpha <- 2 beta <- 1 op <- intrinsic.matern.operators(kappa = kappa, tau = tau, alpha = alpha, beta = beta, mesh = mesh_2d) ``` To see that the `rSPDE` model is approximating the true model, we can compare the variogram of the approximation with the true variogram (implemented in `variogram.intrinsic.spde()`) as follows. ```{r, message=FALSE} point <- matrix(c(1,1),1,2) Gamma <- op$variogram(point) vario <- variogram.intrinsic.spde(point, mesh_2d$loc[,1:2], kappa = kappa, alpha = alpha, tau = tau, beta = beta, L = 2, d = 2) d = sqrt((mesh_2d$loc[,1]-point[1])^2 + (mesh_2d$loc[,2]-point[2])^2) plot(d, Gamma, xlim = c(0,0.5), ylim = c(0,4), ylab = "variogram(h)", xlab = "h") lines(sort(d),sort(vario),col=2, lwd = 2) ``` We can now use the `simulate` function to simulate a realization of the field $u$: ```{r} u <- simulate(op,nsim = 1, use_kl = FALSE) proj <- fm_evaluator(mesh_2d, dims = c(100, 100)) field <- fm_evaluate(proj, field = as.vector(u)) field.df <- data.frame(x1 = proj$lattice$loc[,1], x2 = proj$lattice$loc[,2], y = as.vector(field)) library(ggplot2) library(viridis) ggplot(field.df, aes(x = x1, y = x2, fill = y)) + geom_raster() + scale_fill_viridis() ``` By default, the field is simulated with a zero-integral constraint. ## Fitting the model with `R-INLA` We will now fit the model using our `r inla_link()` implementation of the rational SPDE approach. Further details on this implementation can be found in [R-INLA implementation of the rational SPDE approach](rspde_inla.html). We begin by simulating some data as before. ```{r} n_loc <- 2000 loc_2d_mesh <- matrix(2*runif(n_loc * 2), n_loc, 2) A <- spde.make.A( mesh = mesh_2d, loc = loc_2d_mesh ) sigma.e <- 0.1 y <- A %*% u + rnorm(n_loc) * sigma.e ``` The generated data can be seen in the following image. ```{r,fig.align = "center", echo=TRUE} df <- data.frame(x1 = as.double(loc_2d_mesh[, 1]), x2 = as.double(loc_2d_mesh[, 2]), y = as.double(y)) ggplot(df, aes(x = x1, y = x2, col = y)) + geom_point() + scale_color_viridis() ``` To fit the model, we create the $A$ matrix, the index, and the `inla.stack` object. For now, these more general models can only be estimated with $\beta = 1$ and $\alpha = 1$ or $\alpha = 2$. For these non-fractional models, we can use the standard INLA functions to make the required elements. ```{r} mesh.index <- inla.spde.make.index(name = "field", n.spde = mesh_2d$n) st.dat <- inla.stack(data = list(y = as.vector(y)), A = A, effects = mesh.index) ``` We now create the model object. ```{r} rspde_model <- rspde.intrinsic.matern(mesh = mesh_2d, alpha = alpha) ``` Finally, we create the formula and fit the model to the data: ```{r message=FALSE, warning=FALSE} f <- y ~ -1 + f(field, model = rspde_model) rspde_fit <- inla(f, data = inla.stack.data(st.dat), family = "gaussian", control.predictor = list(A = inla.stack.A(st.dat))) ``` We can get a summary of the fit: ```{r} summary(rspde_fit) ``` To get a summary of the fit of the random field only, we can do the following: ```{r} result_fit <- rspde.result(rspde_fit, "field", rspde_model) summary(result_fit) tau <- op$tau result_df <- data.frame( parameter = c("tau", "kappa"), true = c(tau, kappa), mean = c(result_fit$summary.tau$mean, result_fit$summary.kappa$mean), mode = c(result_fit$summary.tau$mode, result_fit$summary.kappa$mode) ) print(result_df) ``` ## Kriging with `R-INLA` implementation Let us now obtain predictions (i.e., do kriging) of the latent field on a dense grid in the region. We begin by creating the grid of locations where we want to compute the predictions. To this end, we can use the `rspde.mesh.projector()` function. This function has the same arguments as the function `inla.mesh.projector()` the only difference being that the rSPDE version also has an argument `nu` and an argument `rspde.order`. Thus, we proceed in the same fashion as we would in `r inla_link()`'s standard SPDE implementation: ```{r create_proj_grid} projgrid <- inla.mesh.projector(mesh_2d, xlim = c(0, 2), ylim = c(0, 2) ) ``` This lattice contains 100 × 100 locations (the default). Let us now calculate the predictions jointly with the estimation. To this end, first, we begin by linking the prediction coordinates to the mesh nodes through an $A$ matrix ```{r A_prd} A.prd <- projgrid$proj$A ``` We now make a stack for the prediction locations. We have no data at the prediction locations, so we set `y= NA`. We then join this stack with the estimation stack. ```{r stk.prd} ef.prd <- list(c(mesh.index)) st.prd <- inla.stack( data = list(y = NA), A = list(A.prd), tag = "prd", effects = ef.prd ) st.all <- inla.stack(st.dat, st.prd) ``` Doing the joint estimation takes a while, and we therefore turn off the computation of certain things that we are not interested in, such as the marginals for the random effect. We will also use a simplified integration strategy (actually only using the posterior mode of the hyper-parameters) through the command `control.inla = list(int.strategy = "eb")`, i.e. empirical Bayes: ```{r fit_prd, message=FALSE, warning=FALSE} rspde_fitprd <- inla(f, family = "Gaussian", data = inla.stack.data(st.all), control.predictor = list( A = inla.stack.A(st.all), compute = TRUE, link = 1 ), control.compute = list( return.marginals = FALSE, return.marginals.predictor = FALSE ), control.inla = list(int.strategy = "eb") ) ``` We then extract the indices to the prediction nodes and then extract the mean and the standard deviation of the response: ```{r stk.mean.sd} id.prd <- inla.stack.index(st.all, "prd")$data m.prd <- matrix(rspde_fitprd$summary.fitted.values$mean[id.prd], 100, 100) sd.prd <- matrix(rspde_fitprd$summary.fitted.values$sd[id.prd], 100, 100) ``` Finally, we plot the results. First the mean: ```{r plot_pred, echo=TRUE, fig.align='center'} field.pred.df <- data.frame(x1 = projgrid$lattice$loc[,1], x2 = projgrid$lattice$loc[,2], y = as.vector(m.prd)) ggplot(field.pred.df, aes(x = x1, y = x2, fill = y)) + geom_raster() + scale_fill_viridis() ``` Then, the marginal standard deviations: ```{r plot_pred_sd, fig.align='center', echo=TRUE} field.pred.sd.df <- data.frame(x1 = proj$lattice$loc[,1], x2 = proj$lattice$loc[,2], sd = as.vector(sd.prd)) ggplot(field.pred.sd.df, aes(x = x1, y = x2, fill = sd)) + geom_raster() + scale_fill_viridis() ``` # Using intrinsic models without `R-INLA` Currently, the more general model is only implemented in `R-INLA` using fixed integer values of the smoothness parameters. However, all intrinsic models are implemented in `rSPDE` in full generality. In this section, we illustrate the `rSPDE` interface. Let us test a model in one dimension. Let us start with generating the model ```{r} L = 20 x <- seq(from = 0, to = L, length.out = 101) mesh <- fm_mesh_1d(x) beta <- 1.1 alpha <- 0 kappa <- 10 tau <- 10 op <- intrinsic.matern.operators(kappa = kappa, tau = tau, alpha = alpha, beta = beta, mesh = mesh, d = 1) vario <- variogram.intrinsic.spde(c(L/2), mesh$loc, tau = tau, beta = beta, alpha = alpha, kappa = kappa, L = L, d = 1) plot(x, vario, type = "l", col = 2, lwd = 2) points(x,op$variogram(L/2),col=1) ``` We now generate some data. The option to use a mean value correction for extremes models is also implemented, so we generate some data using this. ```{r} n.rep <- 100 u <- simulate(op,nsim = n.rep, integral.constraint = FALSE, use_kl = TRUE) drift <- op$mean_correction() u <- u + matrix(rep(drift, times = n.rep), nrow = op$n, ncol= n.rep) sigma.e <- 0.01 n.obs <- 300 obs.loc <- runif(n = n.obs, min = 0, max = L) A <- rSPDE.A1d(x, obs.loc) Y <- as.matrix(A %*% u + sigma.e * matrix(rnorm(n.obs*n.rep),n.obs,n.rep)) ``` Let us now show how to do kriging prediction for this model. ```{r} A <- make_A(op, loc = obs.loc) A.krig <- make_A(op, loc = x) u.krig <- predict(op, A = A, Aprd = A.krig, Y = Y[,1], sigma.e = sigma.e, compute.variances = TRUE ) plot(obs.loc, Y[,1], ylab = "u(x)", xlab = "x", main = "Data and prediction", ylim = c( min(c(min(u.krig$mean - 2 * sqrt(u.krig$variance)),min(u[,1]))), max(c(max(u.krig$mean + 2 * sqrt(u.krig$variance)), max(u[,1]))) ) ) lines(x,u[,1],col=3) lines(x, u.krig$mean) lines(x, u.krig$mean + 2 * sqrt(u.krig$variance), col = 2) lines(x, u.krig$mean - 2 * sqrt(u.krig$variance), col = 2) ``` We now use `rspde_lme` to fit the parameters based on this data. Since we generated data with `alpha=0`, we specify this in the function to indicate that this parameter should not be fitted but kept fixed at `alpha=0` by setting `fix_alpha=0` in `model_options`. We also specify `mean_correction=TRUE` to indicate that we should use the mean value correction when fitting. ```{r} data = data.frame(y = c(Y), loc = rep(obs.loc, n.rep), rep = rep(1:n.rep, each = n.obs)) fit <- rspde_lme(y ~ -1, loc = "loc", repl = "rep", data = data, model = op, mean_correction = TRUE, parallel = TRUE, model_options = list(fix_alpha = 0)) rbind(c(fit$coeff$random_effects[c("beta", "tau")], fit$coeff$measurement_error), c(beta, tau, sigma.e)) ``` ## An example with estimated alpha and beta parameters In the previous example, we fixed the alpha parameter and only estimated beta. Now, let us demonstrate how to estimate both alpha and beta simultaneously. We will set up a new model with different parameter values: ```{r} L = 20 x <- seq(from = 0, to = L, length.out = 101) mesh <- fm_mesh_1d(x) beta <- 1.2 alpha <- 0.3 kappa <- 15 tau <- 7 op <- intrinsic.matern.operators(kappa = kappa, tau = tau, alpha = alpha, beta = beta, mesh = mesh, d = 1) vario <- variogram.intrinsic.spde(c(L/2), mesh$loc, tau = tau, beta = beta, alpha = alpha, kappa = kappa, L = L, d = 1) plot(x, vario, type = "l", col = 2, lwd = 2) points(x, op$variogram(L/2), col = 1) ``` We can note here that the variogram of the approximate model is not particularly close to the variogram of the true continuous model. The reason for this is that the value of `alpha` is very small, and we therefore need a larger order of the rational approximation than the default value of 2. We can adjust the orders of the rational approximations through the `m_alpha` and `m_beta` values in `intrinsic.matern.operators`. Let us increase the value of `m_alpha` and decrease the value of `m_beta`. ```{r} op <- intrinsic.matern.operators(kappa = kappa, tau = tau, alpha = alpha, beta = beta, mesh = mesh, d = 1, m_alpha = 6, m_beta = 1) vario <- variogram.intrinsic.spde(c(L/2), mesh$loc, tau = tau, beta = beta, alpha = alpha, kappa = kappa, L = L, d = 1) plot(x, vario, type = "l", col = 2, lwd = 2) points(x, op$variogram(L/2), col = 1) ``` We now have a better approximation. Similar to the previous example, we will generate data with the mean value correction for extremes models: ```{r} n.rep <- 100 u <- simulate(op, nsim = n.rep, integral.constraint = FALSE, use_kl = TRUE) drift <- op$mean_correction() u <- u + matrix(rep(drift, times = n.rep), nrow = op$n, ncol = n.rep) sigma.e <- 0.015 n.obs <- 300 obs.loc <- runif(n = n.obs, min = 0, max = L) A <- rSPDE.A1d(x, obs.loc) Y <- as.matrix(A %*% u + sigma.e * matrix(rnorm(n.obs*n.rep), n.obs, n.rep)) ``` Let's visualize the data and predictions for this model: ```{r} A <- make_A(op, loc = obs.loc) A.krig <- make_A(op, loc = x) u.krig <- predict(op, A = A, Aprd = A.krig, Y = Y[,1], sigma.e = sigma.e, compute.variances = TRUE ) plot(obs.loc, Y[,1], ylab = "u(x)", xlab = "x", main = "Data and prediction with fractional alpha and beta", ylim = c( min(c(min(u.krig$mean - 2 * sqrt(u.krig$variance)), min(u[,1]))), max(c(max(u.krig$mean + 2 * sqrt(u.krig$variance)), max(u[,1]))) ) ) lines(x, u[,1], col = 3) lines(x, u.krig$mean) lines(x, u.krig$mean + 2 * sqrt(u.krig$variance), col = 2) lines(x, u.krig$mean - 2 * sqrt(u.krig$variance), col = 2) ``` Now, we will use `rspde_lme` to fit the parameters but this time we will not fix alpha, allowing both alpha and beta to be estimated. Unlike the previous example where we set `fix_alpha=0`, we do not include this constraint: ```{r} data = data.frame(y = c(Y), loc = rep(obs.loc, n.rep), rep = rep(1:n.rep, each = n.obs)) op <- intrinsic.matern.operators(kappa = kappa, tau = tau, alpha = 1.3, beta = 1.05, mesh = mesh, d = 1, m_alpha = 3, m_beta = 1) fit <- rspde_lme(y ~ -1, loc = "loc", repl = "rep", data = data, model = op, mean_correction = TRUE, parallel = FALSE) # Compare estimated and true parameter values rbind(c(fit$coeff$random_effects[c("alpha", "beta", "tau", "kappa")], fit$coeff$measurement_error), c(alpha, beta, tau, kappa, sigma.e)) ```