Bayesian Scalar-on-Function Regression with `refundBayes::sofr_bayes`
2026-03-03
Source:vignettes/sofr_bayes_vignette.Rmd
sofr_bayes_vignette.RmdIntroduction
This vignette provides a detailed guide to the
sofr_bayes() function in the refundBayes
package, which fits Bayesian Scalar-on-Function Regression (SoFR) models
using Stan. The function is designed with a syntax similar to
mgcv::gam, making it accessible to users familiar with the
frequentist approach while providing full Bayesian posterior
inference.
The methodology follows the tutorial by Jiang, Crainiceanu, and Cui (2025), Tutorial on Bayesian Functional Regression Using Stan, published in Statistics in Medicine.
Install the refundBayes Package
The refundBayes package can be installed from
GitHub:
library(remotes)
remotes::install_github("https://github.com/ZirenJiang/refundBayes")Statistical Model
The SoFR Model
Scalar-on-Function Regression (SoFR) models the relationship between a scalar outcome and one or more functional predictors (curves or trajectories observed over a continuum), along with optional scalar covariates.
For subject , let be the scalar outcome, be a vector of scalar predictors, and for be a functional predictor observed at time points over a domain . The SoFR model assumes that the distribution of belongs to an exponential family with mean , and the linear predictor has the following structure:
where:
- is the overall intercept,
- is the unknown functional coefficient that characterizes the effect of the functional predictor on the outcome,
- is a vector of scalar regression coefficients,
- is the domain of the functional predictor. This domain is not restricted to ; it is determined by the actual time points in the data (e.g., for minute-level 24-hour activity data, or for a generic index set).
The integral is approximated using a Riemann sum over the observed time points.
Basis Expansion and Penalized Splines
The functional coefficient is represented nonparametrically using a set of pre-specified spline basis functions :
With this expansion, the linear predictor becomes:
where is a vector with entries:
Here are the Riemann sum integration weights and are the time points at which the functional predictor is observed.
The Role of tmat, lmat, and
wmat
The Riemann sum approximation
to the integral
is constructed directly from three user-supplied matrices in the
data argument:
tmat(an matrix): contains the time points at which the functional predictor is observed. The -th entry equals . The range of values intmatdetermines the domain of integration . For example, if the functional predictor is observed at minutes within a day, thentmathas entries ranging from to and . There is no requirement that the domain be rescaled to .lmat(an matrix): contains the integration weights for the Riemann sum approximation. The -th entry equals . For equally spaced time points with spacing , every entry oflmatequals . For unevenly spaced time points,lmatreflects the varying widths of the integration intervals. These weights, together withtmat, fully specify how the numerical integration is performed and over what domain.wmat(an matrix): contains the functional predictor values. The -th row contains the observed values for subject .
In the formula
s(tmat, by = lmat * wmat, bs = "cc", k = 10), the
mgcv infrastructure uses tmat to construct the
spline basis
at the observed time points, and the by = lmat * wmat
argument provides the element-wise product
that enters the Riemann sum. This means the basis functions are
evaluated on the scale of tmat, and the integration weights
in lmat ensure that the discrete sum correctly approximates
the integral over the actual domain
— regardless of whether it is
,
,
or any other interval.
Note that for all subjects, the time points are assumed to be
identical so that
for all
.
Thus every row of tmat is the same, and every row of
lmat is the same. The matrices are replicated across rows
to match the mgcv syntax, which expects all terms in the
formula to have the same dimensions.
Smoothness Penalty
To induce smoothness on , a quadratic penalty on the spline coefficients is applied. The penalty is based on the integrated squared second derivative of :
where is the penalty matrix. In the Bayesian framework, this penalty is equivalent to placing a multivariate normal prior on the spline coefficients:
where is the smoothing parameter that controls the smoothness of , and is estimated from the data.
The sofr_bayes() Function
Usage
sofr_bayes(
formula,
data,
family = gaussian(),
joint_FPCA = NULL,
intercept = TRUE,
runStan = TRUE,
niter = 3000,
nwarmup = 1000,
nchain = 3,
ncores = 1
)Arguments
| Argument | Description |
|---|---|
formula |
Functional regression formula, using the same syntax as
mgcv::gam. Functional predictors are specified using the
s() term with by = lmat * wmat to encode the
Riemann sum integration (see Example below). |
data |
A data frame containing all scalar and functional variables used in the model. |
family |
Distribution of the outcome variable. Currently
supports gaussian() and binomial(). Default is
gaussian(). |
joint_FPCA |
A logical (TRUE/FALSE) vector
of the same length as the number of functional predictors, indicating
whether to jointly model FPCA for each functional predictor. Default is
NULL, which sets all entries to FALSE (no
joint FPCA). |
intercept |
Logical. Whether to include an intercept term in the
linear predictor. Default is TRUE. |
runStan |
Logical. Whether to run the Stan program. If
FALSE, the function only generates the Stan code and data
without sampling. This is useful for inspecting or modifying the
generated Stan code. Default is TRUE. |
niter |
Total number of Bayesian posterior sampling iterations
(including warmup). Default is 3000. |
nwarmup |
Number of warmup (burn-in) iterations. These samples
are discarded and not used for inference. Default is
1000. |
nchain |
Number of Markov chains for posterior sampling.
Multiple chains help assess convergence. Default is 3. |
ncores |
Number of CPU cores to use when executing the chains in
parallel. Default is 1. |
Return Value
The function returns a list of class "refundBayes"
containing the following elements:
| Element | Description |
|---|---|
stanfit |
The Stan fit object (class stanfit). Can
be used for convergence diagnostics, traceplots, and additional
summaries via the rstan package. |
spline_basis |
Basis functions used to reconstruct the functional coefficients from the posterior samples. |
stancode |
A character string containing the generated Stan model code. |
standata |
A list containing the data passed to the Stan model. |
int |
A vector of posterior samples for the intercept term
.
NULL if intercept = FALSE. |
scalar_coef |
A matrix of posterior samples for scalar coefficients
,
where each row is one posterior sample and each column corresponds to
one scalar predictor. NULL if no scalar predictors are
included. |
func_coef |
A list of posterior samples for functional coefficients. Each element is a matrix where each row is one posterior sample and each column corresponds to one location on the functional domain. |
family |
The distribution family used for the outcome. |
Formula Syntax
The formula follows the mgcv::gam syntax. The key
component for specifying functional predictors is:
s(tmat, by = lmat * wmat, bs = "cc", k = 10)where:
-
tmat: an matrix of time points. Each row contains the same observation times (replicated across subjects). The values intmatdefine the domain over which the functional predictor is observed and the integral is computed. For example,tmatmay contain values from to , from to , or from to — the integration domain adapts accordingly. -
lmat: an matrix of Riemann sum weights. The -th entry equals . These weights control the numerical integration and should be consistent with the spacing of the time points intmat. For equally spaced time points,lmatis a constant matrix; for irregular spacing, it reflects the actual gaps between consecutive time points. -
wmat: an matrix of functional predictor values. The -th row contains the observed values for subject . -
bs: the type of spline basis. Common choices include"cr"(cubic regression splines) and"cc"(cyclic cubic regression splines, suitable for periodic data). -
k: the number of basis functions (maximum degrees of freedom for the spline).
Scalar predictors are included as standard formula terms (e.g.,
y ~ X1 + s(tmat, by = lmat * wmat, bs = "cr", k = 10)).
Example: Bayesian SoFR with Binary Outcome
We demonstrate the sofr_bayes() function using a
simulated example dataset with a binary outcome and one functional
predictor.
Load Data
## Alternative sample data
# data.SoFR <- readRDS("data/example_data_sofr.rds")
## Load the example data
set.seed(123)
n <- 100
M <- 50
tgrid <- seq(0, 1, length.out = M)
dt <- tgrid[2] - tgrid[1]
tmat <- matrix(rep(tgrid, each = n), nrow = n)
lmat <- matrix(dt, nrow = n, ncol = M)
wmat <- t(apply(matrix(rnorm(n * M), n, M), 1, cumsum)) / sqrt(M)
beta_true <- sin(2 * pi * tgrid)
X1 <- rnorm(n)
eta <- 0.5 * X1 + wmat %*% (beta_true * dt)
prob <- plogis(eta)
y <- rbinom(n, 1, prob)
data.SoFR <- data.frame(y = y, X1 = X1)
data.SoFR$tmat <- tmat
data.SoFR$lmat <- lmat
data.SoFR$wmat <- wmatThe example dataset data.SoFR contains:
-
y: a binary outcome variable, -
X1: a scalar predictor, -
tmat: the time point matrix (defines the domain ), -
lmat: the Riemann sum weight matrix (defines the integration weights over ), -
wmat: the functional predictor matrix.
Fit the Bayesian SoFR Model
library(refundBayes)
refundBayes_SoFR <- refundBayes::sofr_bayes(
y ~ X1 + s(tmat, by = lmat * wmat, bs = "cc", k = 10),
data = data.SoFR,
family = binomial(),
runStan = TRUE,
niter = 1500,
nwarmup = 500,
nchain = 3,
ncores = 3
)In this call:
- The formula specifies a binary outcome
ywith one scalar predictorX1and one functional predictor encoded vias(tmat, by = lmat * wmat, bs = "cc", k = 10). - The spline basis is evaluated at the time points stored in
tmat, and the integration over the domain (determined bytmat) is approximated using the weights inlmat. -
bs = "cc"uses cyclic cubic regression splines, which are appropriate when the functional predictor is periodic (e.g., physical activity measured over a 24-hour cycle). -
k = 10specifies 10 basis functions. In practice, 30–40 basis functions are often sufficient for moderately smooth functional data on dense grids. -
family = binomial()specifies a logistic regression for the binary outcome. - The sampler runs 3 chains in parallel, each with 1500 total iterations (500 warmup + 1000 posterior samples).
Visualization
The plot() method for refundBayes objects
displays the estimated functional coefficient
along with pointwise 95% credible intervals:
Extracting Posterior Summaries
Posterior summaries of the functional coefficient can be computed
directly from the func_coef element:
## Posterior mean of the functional coefficient
mean_curve <- apply(refundBayes_SoFR$func_coef[[1]], 2, mean)
## Pointwise 95% credible interval
upper_curve <- apply(refundBayes_SoFR$func_coef[[1]], 2,
function(x) quantile(x, prob = 0.975))
lower_curve <- apply(refundBayes_SoFR$func_coef[[1]], 2,
function(x) quantile(x, prob = 0.025))The posterior samples in func_coef[[1]] are stored as a
matrix, where
is the number of posterior samples and
is the number of time points on the functional domain.
Comparison with Frequentist Results
The Bayesian results can be compared with frequentist estimates
obtained via mgcv::gam:
library(mgcv)
## Fit frequentist SoFR using mgcv
fit_freq <- gam(
y ~ s(tmat, by = lmat * wmat, bs = "cc", k = 10) + X1,
data = data.SoFR,
family = "binomial"
)
## Extract frequentist estimates
freq_result <- plot(fit_freq)The functional coefficient estimates from the Bayesian and frequentist approaches are generally comparable in shape and magnitude, though the Bayesian credible intervals tend to be slightly wider due to accounting for uncertainty in the smoothing parameter.
Inspecting the Generated Stan Code
Setting runStan = FALSE allows you to inspect or modify
the Stan code before running the model:
## Generate Stan code without running the sampler
sofr_code <- refundBayes::sofr_bayes(
y ~ X1 + s(tmat, by = lmat * wmat, bs = "cc", k = 10),
data = data.SoFR,
family = binomial(),
runStan = FALSE
)
## Print the generated Stan code
cat(sofr_code$stancode)Practical Recommendations
-
Number of basis functions (
k): For illustrative purposes,k = 10is often used. In practice, 30–40 basis functions are recommended for moderately smooth functional data observed on dense grids. -
Spline type (
bs): Use"cr"(cubic regression splines) for general functional predictors. Use"cc"(cyclic cubic regression splines) when the functional predictor is periodic (e.g., 24-hour activity patterns). -
Number of iterations: Ensure sufficient posterior
samples for reliable inference. A common setup is
niter = 3000withnwarmup = 1000, but more complex models may require additional iterations. -
Convergence diagnostics: After fitting, examine
traceplots and
statistics using the
rstanpackage (e.g.,rstan::traceplot(refundBayes_SoFR$stanfit)) to ensure that the Markov chains have converged. -
Joint FPCA: When functional predictors are measured
with substantial noise, consider setting
joint_FPCA = TRUEfor the relevant predictor to jointly estimate FPCA scores and regression coefficients.
References
- Jiang, Z., Crainiceanu, C., and Cui, E. (2025). Tutorial on Bayesian Functional Regression Using Stan. Statistics in Medicine, 44(20–22), e70265.
- Crainiceanu, C. M., Goldsmith, J., Leroux, A., and Cui, E. (2024). Functional Data Analysis with R. CRC Press.
- Crainiceanu, C. M. and Goldsmith, A. J. (2010). Bayesian Functional Data Analysis Using WinBUGS. Journal of Statistical Software, 32(11), 1–33.
- Wood, S. (2001). mgcv: GAMs and Generalized Ridge Regression for R. R News, 1(2), 20–25.
- Carpenter, B., Gelman, A., Hoffman, M. D., et al. (2017). Stan: A Probabilistic Programming Language. Journal of Statistical Software, 76(1), 1–32.