Title: | Spatial Bayesian Methods for Task Functional MRI Studies |
---|---|
Description: | Performs a spatial Bayesian general linear model (GLM) for task functional magnetic resonance imaging (fMRI) data on the cortical surface. Additional models include group analysis and inference to detect thresholded areas of activation. Includes direct support for the 'CIFTI' neuroimaging file format. For more information see A. F. Mejia, Y. R. Yue, D. Bolin, F. Lindgren, M. A. Lindquist (2020) <doi:10.1080/01621459.2019.1611582> and D. Spencer, Y. R. Yue, D. Bolin, S. Ryan, A. F. Mejia (2022) <doi:10.1016/j.neuroimage.2022.118908>. |
Authors: | Amanda Mejia [aut, cre], Damon Pham [ctb] , David Bolin [ctb], Yu (Ryan) Yue [ctb], Daniel Spencer [aut] , Sarah Ryan [ctb] |
Maintainer: | Amanda Mejia <[email protected]> |
License: | GPL-3 |
Version: | 0.8.0 |
Built: | 2024-11-23 05:07:19 UTC |
Source: | https://github.com/mandymejia/bayesfmri |
Performs a spatial Bayesian general linear model (GLM) for task functional magnetic resonance imaging (fMRI) data on the cortical surface. Additional models include group analysis and inference to detect thresholded areas of activation. Includes direct support for the 'CIFTI' neuroimaging file format. For more information see A. F. Mejia, Y. R. Yue, D. Bolin, F. Lindgren, M. A. Lindquist (2020) doi:10.1080/01621459.2019.1611582 and D. Spencer, Y. R. Yue, D. Bolin, S. Ryan, A. F. Mejia (2022) doi:10.1016/j.neuroimage.2022.118908.
Maintainer: Amanda Mejia [email protected]
Authors:
Daniel Spencer [email protected] (ORCID)
Other contributors:
Damon Pham [email protected] (ORCID) [contributor]
David Bolin [email protected] [contributor]
Yu (Ryan) Yue [email protected] [contributor]
Sarah Ryan [email protected] [contributor]
Useful links:
Report bugs at https://github.com/mandymejia/BayesfMRI/issues
Perform the EM algorithm of the Bayesian GLM fitting
.findTheta(theta, spde, y, X, QK, Psi, A, Ns, tol, verbose = FALSE)
.findTheta(theta, spde, y, X, QK, Psi, A, Ns, tol, verbose = FALSE)
theta |
the vector of initial values for theta |
spde |
a list containing the sparse matrix elements Cmat, Gmat, and GtCinvG |
y |
the vector of response values |
X |
the sparse matrix of the data values |
QK |
a sparse matrix of the prior precision found using the initial values of the hyperparameters |
Psi |
a sparse matrix representation of the basis function mapping the data locations to the mesh vertices |
A |
a precomputed matrix crossprod(X%*%Psi) |
Ns |
the number of columns for the random matrix used in the Hutchinson estimator |
tol |
a value for the tolerance used for a stopping rule (compared to
the squared norm of the differences between |
verbose |
(logical) Should intermediate output be displayed? |
Get the prewhitening matrix for a single data location
.getSqrtInvCpp(AR_coefs, nTime, avg_var)
.getSqrtInvCpp(AR_coefs, nTime, avg_var)
AR_coefs |
a length-p vector where p is the AR order |
nTime |
(integer) the length of the time series that is being prewhitened |
avg_var |
a scalar value of the residual variances of the AR model |
Find the initial values of kappa2 and phi
.initialKP(theta, spde, w, n_sess, tol, verbose)
.initialKP(theta, spde, w, n_sess, tol, verbose)
theta |
a vector of length two containing the range and scale parameters kappa2 and phi, in that order |
spde |
a list containing the sparse matrix elements Cmat, Gmat, and GtCinvG |
w |
the beta_hat estimates for a single task |
n_sess |
the number of sessions |
tol |
the stopping rule tolerance |
verbose |
(logical) Should intermediate output be displayed? |
Find the log of the determinant of Q_tilde
.logDetQt(kappa2, in_list, n_sess)
.logDetQt(kappa2, in_list, n_sess)
kappa2 |
a scalar |
in_list |
a list with elements Cmat, Gmat, and GtCinvG |
n_sess |
the integer number of sessions |
Identify areas of activation for each field from the result of BayesGLM
or fit_bayesglm
.
activations( x, Bayes = TRUE, gamma = NULL, alpha = 0.05, correction = c("FWER", "FDR", "none"), fields = NULL, sessions = NULL, verbose = 1 ) id_activations( x, Bayes = TRUE, gamma = NULL, alpha = 0.05, correction = c("FWER", "FDR", "none"), fields = NULL, sessions = NULL, verbose = 1 )
activations( x, Bayes = TRUE, gamma = NULL, alpha = 0.05, correction = c("FWER", "FDR", "none"), fields = NULL, sessions = NULL, verbose = 1 ) id_activations( x, Bayes = TRUE, gamma = NULL, alpha = 0.05, correction = c("FWER", "FDR", "none"), fields = NULL, sessions = NULL, verbose = 1 )
x |
Result of |
Bayes |
Use spatial Bayesian modeling to identify activations based on
the joint posterior distribution? Default: |
gamma |
Activation threshold, for example |
alpha |
Significance level. Default: |
correction |
For the classical method only: Type of multiple comparisons
correction: |
fields |
The field(s) to identify activations for. Give either the name(s)
as a character vector, or the numerical indices. If |
sessions |
The session(s) to identify activations for. Give either the
name(s) as a character vector, or the numerical indices. If |
verbose |
|
An "act_BGLM"
or "act_fit_bglm"
object, a
list which indicates the activated locations along with related information.
aic
aic |
(For prewhitening) Use the Akaike information criterion (AIC) to
select AR model orders between |
ar_order
ar_order |
(For prewhitening) The order of the autoregressive (AR) model
to use for prewhitening. If For multi-session modeling, note that a single AR model is used; its coefficients will be the average estimate from each session. |
ar_smooth
ar_smooth |
(For prewhitening) The FWHM parameter for spatially
smoothing the coefficient estimates for the AR model to use for
prewhitening. Recall that
|
Bayes
Bayes |
Perform spatial Bayesian modeling? Default: |
Performs spatial Bayesian GLM for task fMRI activation with CIFTI-format data. The cortex is modeled as a surface mesh, and subcortical structures are modeled as distinct volumetric regions. Includes the pre-processing steps of nuisance regression, prewhitening, scaling, and variance normalization. Supports both single- and multi-session analysis. Can also compute just the classical (spatially-independent)
BayesGLM( BOLD, brainstructures = c("left", "right"), subROI = c("Amygdala-L", "Amygdala-R", "Caudate-L", "Caudate-R", "Hippocampus-L", "Hippocampus-R", "Thalamus-L", "Thalamus-R"), design, nuisance = NULL, hpf = NULL, TR = NULL, surfL = NULL, surfR = NULL, resamp_res = 10000, nbhd_order = 1, buffer = c(1, 1, 3, 4, 4), session_names = NULL, scale_BOLD = c("mean", "sd", "none"), Bayes = TRUE, hyperpriors = c("informative", "default"), ar_order = 6, ar_smooth = 5, aic = FALSE, n_threads = 4, return_INLA = c("trimmed", "full", "minimal"), verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
BayesGLM( BOLD, brainstructures = c("left", "right"), subROI = c("Amygdala-L", "Amygdala-R", "Caudate-L", "Caudate-R", "Hippocampus-L", "Hippocampus-R", "Thalamus-L", "Thalamus-R"), design, nuisance = NULL, hpf = NULL, TR = NULL, surfL = NULL, surfR = NULL, resamp_res = 10000, nbhd_order = 1, buffer = c(1, 1, 3, 4, 4), session_names = NULL, scale_BOLD = c("mean", "sd", "none"), Bayes = TRUE, hyperpriors = c("informative", "default"), ar_order = 6, ar_smooth = 5, aic = FALSE, n_threads = 4, return_INLA = c("trimmed", "full", "minimal"), verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
BOLD |
fMRI timeseries data in CIFTI format ("*.dtseries.nii").
For single-session analysis this can be a file path to a CIFTI file or a
If |
brainstructures |
Character vector indicating which brain structure(s)
of |
subROI |
Which subcortical ROIs should be analyzed? |
design |
A numeric matrix or |
nuisance |
(Optional) A Detrending/high-pass filtering is accomplished by adding DCT bases to the
nuisance matrix; see the parameters |
hpf |
Add DCT bases to Using at least two DCT bases is as sufficient for detrending as using linear
and quadratic drift terms in the nuisance matrix. So if DCT detrending is
being used here, there is no need to add linear and quadratic drift terms to
|
TR |
Temporal resolution of the data, in seconds. |
surfL , surfR
|
For cortex spatial model. Left and right cortex surface
geometry in GIFTI format ("*.surf.gii"). These can be a file path to
a GIFTI file or a Surfaces can alternatively be provided through the |
resamp_res |
For cortex spatial model. The number of vertices to which
each cortical surface should be resampled, or For computational feasibility, a value of |
nbhd_order |
For volumetric model. What order neighborhood around data
locations to keep? |
buffer |
For volumetric model. The number of extra voxel layers around
the bounding box. Set to |
session_names |
The names of the task-fMRI |
scale_BOLD |
Controls scaling the BOLD response at each location.
|
Bayes |
Perform spatial Bayesian modeling? Default: |
hyperpriors |
Should informative or default non-informative hyperpriors be assumed on SPDE hyperparameters? |
ar_order |
(For prewhitening) The order of the autoregressive (AR) model
to use for prewhitening. If For multi-session modeling, note that a single AR model is used; its coefficients will be the average estimate from each session. |
ar_smooth |
(For prewhitening) The FWHM parameter for spatially
smoothing the coefficient estimates for the AR model to use for
prewhitening. Recall that
|
aic |
(For prewhitening) Use the Akaike information criterion (AIC) to
select AR model orders between |
n_threads |
The maximum number of threads to use for parallel
computations: prewhitening parameter estimation, and the inla-program model
estimation. Default: |
return_INLA |
Return the INLA model object? (It can be large.) Use
|
verbose |
|
meanTol , varTol
|
Tolerance for mean and variance of each data location.
Locations which do not meet these thresholds are masked out of the analysis.
Default: |
To use BayesGLM
, the design matrix must first be constructed
with make_design
.
An object of class "BayesGLM"
: a list with elements
The field coefficients for the Bayesian model.
The field coefficients for the classical model.
The entire list of GLM results, except for parameters estimated for the classical model.
Parameters estimated for the classical model from the GLM.
data.frame
summarizing the spatial features of each brain structure modeled.
data.frame
with the name
and nTime
of each BOLD session.
data.frame
with the name
, related task
, and HRF_order
of each field.
This function uses a system wrapper for the 'wb_command' executable. The user must first download and install the Connectome Workbench, available from https://www.humanconnectome.org/software/get-connectome-workbench .
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
INLA computation times increase greatly when the number of columns in the
design matrix exceeds five: when there are more than five tasks, or more
than three tasks each with a temporal derivative modeled as a field. In
cases like the latter, we recommend modeling the temporal derivatives as
nuisance signals using the option dHRF_as="nuisance"
, rather than
modeling the temporal derivatives as fields.
Performs group-level Bayesian GLM estimation and inference using the joint approach described in Mejia et al. (2020).
BayesGLM2( results, contrasts = NULL, quantiles = NULL, excursion_type = NULL, contrast_names = NULL, gamma = 0, alpha = 0.05, nsamp_theta = 50, nsamp_beta = 100, num_cores = NULL, verbose = 1 ) BayesGLM_group( results, contrasts = NULL, quantiles = NULL, excursion_type = NULL, gamma = 0, alpha = 0.05, nsamp_theta = 50, nsamp_beta = 100, num_cores = NULL, verbose = 1 )
BayesGLM2( results, contrasts = NULL, quantiles = NULL, excursion_type = NULL, contrast_names = NULL, gamma = 0, alpha = 0.05, nsamp_theta = 50, nsamp_beta = 100, num_cores = NULL, verbose = 1 ) BayesGLM_group( results, contrasts = NULL, quantiles = NULL, excursion_type = NULL, gamma = 0, alpha = 0.05, nsamp_theta = 50, nsamp_beta = 100, num_cores = NULL, verbose = 1 )
results |
Either (1) a length |
contrasts |
(Optional) A list of contrast vectors that specify the
group-level summaries of interest. If Each contrast vector is length For a single session/subject, the contrast vector for the first field would be:
so the full contrast vector for the group average over all sessions/subjects for the first field would be:
To obtain the group average for the first field, for just the first session, input zeros for the remaining sessions:
To obtain the group mean difference between two sessions (
To obtain the mean over sessions of the first field, just for the first subject:
|
quantiles |
(Optional) Vector of posterior quantiles to return in addition to the posterior mean. |
excursion_type |
(For inference only) The type of excursion function for
the contrast (">", "<", "!="), or a vector thereof (each element
corresponding to one contrast). If |
contrast_names |
(Optional) Names of contrasts. |
gamma |
(For inference only) Activation threshold for the excursion set,
or a vector thereof (each element corresponding to one contrast). Default:
|
alpha |
(For inference only) Significance level for activation for the
excursion set, or a vector thereof (each element corresponding to one
contrast). Default: |
nsamp_theta |
Number of theta values to sample from posterior. Default:
|
nsamp_beta |
Number of beta vectors to sample conditional on each theta
value sampled. Default: |
num_cores |
The number of cores to use for sampling betas in parallel. If
|
verbose |
|
A list containing the estimates, PPMs and areas of activation for each contrast.
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
BOLD
BOLD |
fMRI timeseries data in CIFTI format ("*.dtseries.nii").
For single-session analysis this can be a file path to a CIFTI file or a
If |
brainstructures
brainstructures |
Character vector indicating which brain structure(s)
of |
buffer
buffer |
For volumetric model. The number of extra voxel layers around
the bounding box. Set to |
Take the central derivative of numeric vectors by averaging the forward and backward differences.
cderiv(x)
cderiv(x)
x |
A numeric matrix, or a vector which will be converted to a single-column matrix. |
A matrix or vector the same dimensions as x
, with the
derivative taken for each column of x
. The first and last rows may
need to be deleted, depending on the application.
x <- cderiv(seq(5)) stopifnot(all(x == c(.5, 1, 1, 1, .5)))
x <- cderiv(seq(5)) stopifnot(all(x == c(.5, 1, 1, 1, .5)))
Connectome Workbench
This function uses a system wrapper for the 'wb_command' executable. The user must first download and install the Connectome Workbench, available from https://www.humanconnectome.org/software/get-connectome-workbench .
contrasts
contrasts |
List of contrast vectors to be passed to |
design
design |
A numeric matrix or |
EM
EM |
(logical) Should the EM implementation of the Bayesian GLM be used?
Default: |
faces
faces |
An |
field_names
field_names |
(Optional) Names of fields represented in design matrix. |
Performs spatial Bayesian GLM for task fMRI activation
fit_bayesglm( BOLD, design, nuisance = NULL, spatial, scale_BOLD = c("mean", "sd", "none"), Bayes = TRUE, hyperpriors = c("informative", "default"), ar_order = 6, ar_smooth = 5, aic = FALSE, n_threads = 4, return_INLA = c("trimmed", "full", "minimal"), verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
fit_bayesglm( BOLD, design, nuisance = NULL, spatial, scale_BOLD = c("mean", "sd", "none"), Bayes = TRUE, hyperpriors = c("informative", "default"), ar_order = 6, ar_smooth = 5, aic = FALSE, n_threads = 4, return_INLA = c("trimmed", "full", "minimal"), verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
BOLD , design , nuisance
|
Session-length list of numeric matrices/arrays, each with volumes along the first dimension. |
spatial |
Gives the spatial information:
For voxel data, a list of six: |
scale_BOLD |
Controls scaling the BOLD response at each location.
|
Bayes |
Perform spatial Bayesian modeling? Default: |
hyperpriors |
Should informative or default non-informative hyperpriors be assumed on SPDE hyperparameters? |
ar_order |
(For prewhitening) The order of the autoregressive (AR) model
to use for prewhitening. If For multi-session modeling, note that a single AR model is used; its coefficients will be the average estimate from each session. |
ar_smooth |
(For prewhitening) The FWHM parameter for spatially
smoothing the coefficient estimates for the AR model to use for
prewhitening. Recall that
|
aic |
(For prewhitening) Use the Akaike information criterion (AIC) to
select AR model orders between |
n_threads |
The maximum number of threads to use for parallel
computations: prewhitening parameter estimation, and the inla-program model
estimation. Default: |
return_INLA |
Return the INLA model object? (It can be large.) Use
|
verbose |
|
meanTol , varTol
|
Tolerance for mean, variance and SNR of each data location.
Locations which do not meet these thresholds are masked out of the analysis.
Default: |
A "BayesGLM"
object: a list with elements
The full result of the call to INLA::inla
.
The estimated coefficients for the Bayesian model.
Results from the classical model: field estimates, field standard error estimates, residuals, degrees of freedom, and the mask.
The model mesh including only the locations analyzed, i.e. within mask
, without missing values, and meeting meanTol
and varTol
.
The original mesh provided.
A mask of mesh_orig
indicating the locations inside mesh
.
The design matrix, after centering and scaling, but before any nuisance regression or prewhitening.
The names of the fields.
The names of the sessions.
Hyperparameter posterior densities.
Theta estimates from the Bayesian model.
For joint group modeling.
For joint group modeling.
For joint group modeling.
For joint group modeling: The BOLD data after any centering, scaling, nuisance regression, or prewhitening.
For joint group modeling: The design matrix after any centering, scaling, nuisance regression, or prewhitening.
Vectors of values across locations: phi
(AR coefficients averaged across sessions), sigma_sq
(residual variance averaged across sessions), and AIC (the maximum across sessions).
match.call() for this function call.
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
hpf
hpf |
Add DCT bases to Using at least two DCT bases is as sufficient for detrending as using linear
and quadratic drift terms in the nuisance matrix. So if DCT detrending is
being used here, there is no need to add linear and quadratic drift terms to
|
Calculate the HRF from a time vector and parameters, or its derivative with respect to delay or dispersion.
HRF_calc( t, deriv = 0, a1 = 6, b1 = 1, a2 = 16/6 * a1 * sqrt(b1), b2 = b1, c = 1/6, o = 0 )
HRF_calc( t, deriv = 0, a1 = 6, b1 = 1, a2 = 16/6 * a1 * sqrt(b1), b2 = b1, c = 1/6, o = 0 )
t |
time vector (in units of seconds) |
deriv |
|
a1 |
delay of response. Default: |
b1 |
response dispersion. Default: |
a2 |
delay of undershoot. Default: |
b2 |
dispersion of undershoot. Default: |
c |
scale of undershoot. Default: |
o |
onset of response. Default: |
HRF vector (or dHRF, or d2HRF) corresponding to time vector t
Calculate the HRF from a time vector and parameters. Optionally compute the first or second derivative of the HRF instead. Form of HRF is similar to SPM but here the response and undershoot are scaled so the difference of the HRFs peaks at 1 and -c
HRF_main(t, a1 = 6, b1 = 1, a2 = NULL, b2 = NULL, c = 1/6, o = 0)
HRF_main(t, a1 = 6, b1 = 1, a2 = NULL, b2 = NULL, c = 1/6, o = 0)
t |
time vector (in seconds). Must be equally spaced. |
a1 |
delay of response. Default: |
b1 |
response dispersion. Default: |
a2 |
delay of undershoot. Default: |
b2 |
dispersion of undershoot. Default: |
c |
scale of undershoot. Default: |
o |
onset of response (in seconds). Default: |
HRF vector corresponding to time vector t
Calculate the HRF from a time vector and parameters. Optionally compute the first or second derivative of the HRF instead.
HRF96(t, deriv = 0, a1 = 6, b1 = 0.9, a2 = 12, b2 = 0.9, c = 0.35)
HRF96(t, deriv = 0, a1 = 6, b1 = 0.9, a2 = 12, b2 = 0.9, c = 0.35)
t |
time vector |
deriv |
|
a1 |
delay of response. Default: |
b1 |
response dispersion. Default: |
a2 |
delay of undershoot. Default: |
b2 |
dispersion of undershoot. Default: |
c |
scale of undershoot. Default: |
HRF vector (or dHRF, or d2HRF) corresponding to time
upsample <- 100 HRF96(seq(0, 30, by=1/upsample))
upsample <- 100 HRF96(seq(0, 30, by=1/upsample))
INLA
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
INLA Latent Fields
INLA computation times increase greatly when the number of columns in the
design matrix exceeds five: when there are more than five tasks, or more
than three tasks each with a temporal derivative modeled as a field. In
cases like the latter, we recommend modeling the temporal derivatives as
nuisance signals using the option dHRF_as="nuisance"
, rather than
modeling the temporal derivatives as fields.
Make the design matrix for the GLM, from the task information.
make_design( EVs, nTime, TR, dHRF = 0, upsample = 100, onset = NULL, offset = NULL, scale_design = TRUE, onsets_sep = FALSE, offsets_sep = FALSE, verbose = TRUE, ... )
make_design( EVs, nTime, TR, dHRF = 0, upsample = 100, onset = NULL, offset = NULL, scale_design = TRUE, onsets_sep = FALSE, offsets_sep = FALSE, verbose = TRUE, ... )
EVs |
The explanatory variables i.e. the task stimulus information,
from which a design matrix will be constructed. This is a list where each
entry represents a task as a matrix of onsets (first column) and durations
(second column) for each stimuli (each row) of the task, in seconds. List
names should be the task names. An example of a properly-formatted |
nTime |
the number of timepoints (volumes) in the task fMRI data. |
TR |
the temporal resolution of the data, in seconds. |
dHRF |
Controls the extent of HRF derivatives modeling. Set to If there are several tasks and |
upsample |
Upsample factor for convolving stimulus boxcar or stick
function with canonical HRF. Default: |
onset , offset
|
Add task regressors indicating the onset and/or offset of
each event block? Provide the names of the tasks as a character vector. All
onsets (or offsets) across the specified tasks will be represented by one
additional column in the design matrix. The task names must match the names
of Onsets/offset modeling is only compatible with a block design experiment.
An error will be raised if the events in |
scale_design |
Scale the columns of the design matrix? Default:
|
onsets_sep , offsets_sep
|
Model the onsets ( |
verbose |
Print diagnostic messages? Default: |
... |
Additional arguments to |
A "BfMRI_design"
object: a list with elements
The volumes by fields design matrix. Column names are field names.
The name of each task from the provided onsets.
The input dHRF
parameter.
Additional HRF modeling results.
Mask out data locations that are invalid (missing data, low mean, or low variance) for any session.
make_mask(BOLD, meanTol = 1e-06, varTol = 1e-06, verbose = TRUE)
make_mask(BOLD, meanTol = 1e-06, varTol = 1e-06, verbose = TRUE)
BOLD |
A session-length list of |
meanTol , varTol
|
Tolerance for mean and variance of each data location.
Locations which do not meet these thresholds are masked out of the analysis.
Defaults: |
verbose |
Print messages counting how many locations are removed?
Default: |
A logical vector indicating locations that are valid across all sessions.
nT <- 30 nV <- 400 BOLD1 <- matrix(rnorm(nT*nV), nrow=nT) BOLD1[,seq(30,50)] <- NA BOLD2 <- matrix(rnorm(nT*nV), nrow=nT) BOLD2[,65] <- BOLD2[,65] / 1e10 BOLD <- list(sess1=BOLD1, sess2=BOLD2) make_mask(BOLD)
nT <- 30 nV <- 400 BOLD1 <- matrix(rnorm(nT*nV), nrow=nT) BOLD1[,seq(30,50)] <- NA BOLD2 <- matrix(rnorm(nT*nV), nrow=nT) BOLD2[,65] <- BOLD2[,65] / 1e10 BOLD <- list(sess1=BOLD1, sess2=BOLD2) make_mask(BOLD)
Make INLA triangular mesh from faces and vertices
make_mesh(vertices, faces)
make_mesh(vertices, faces)
vertices |
A |
faces |
An |
INLA triangular mesh
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
mask: vertices
mask |
A length |
max_threads
max_threads |
The maximum number of threads to use in the inla-program
for model estimation. |
mean and variance tolerance
meanTol , varTol
|
Tolerance for mean and variance of each data location.
Locations which do not meet these thresholds are masked out of the analysis.
Default: |
mesh: either
mesh |
An |
mesh: INLA only
mesh |
An |
Performs classical Bayesian GLM for task fMRI activation with CIFTI-format data, evaluating multiple design matrices. Includes the pre-processing steps of nuisance regression. Supports single-session analysis only.
multiGLM( BOLD, design, brainstructures = c("left", "right"), TR = NULL, resamp_res = 10000, hpf = NULL, nuisance = NULL, design_canonical = NULL, verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
multiGLM( BOLD, design, brainstructures = c("left", "right"), TR = NULL, resamp_res = 10000, hpf = NULL, nuisance = NULL, design_canonical = NULL, verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
BOLD |
fMRI timeseries data in CIFTI format ("*.dtseries.nii").
For single-session analysis this can be a file path to a CIFTI file or a
If |
design |
A 3D numeric array that is locations by fields by designs. |
brainstructures |
Character vector indicating which brain structure(s)
of |
TR |
Temporal resolution of the data, in seconds. |
resamp_res |
For cortex spatial model. The number of vertices to which
each cortical surface should be resampled, or For computational feasibility, a value of |
hpf |
Add DCT bases to Using at least two DCT bases is as sufficient for detrending as using linear
and quadratic drift terms in the nuisance matrix. So if DCT detrending is
being used here, there is no need to add linear and quadratic drift terms to
|
nuisance |
(Optional) A Detrending/high-pass filtering is accomplished by adding DCT bases to the
nuisance matrix; see the parameters |
design_canonical |
TO DO |
verbose |
|
meanTol , varTol
|
Tolerance for mean and variance of each data location.
Locations which do not meet these thresholds are masked out of the analysis.
Default: |
An object of class "mGLM"
: a list with elements
data.frame
summarizing the spatial features of each brain structure modeled.
data.frame
with the name
, related task
, and HRF_order
of each field.
This function uses a system wrapper for the 'wb_command' executable. The user must first download and install the Connectome Workbench, available from https://www.humanconnectome.org/software/get-connectome-workbench .
Performs classical GLM for task fMRI activation, comparing multiple designs
multiGLM_fun( BOLD, design, nuisance = NULL, design_canonical = NULL, verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
multiGLM_fun( BOLD, design, nuisance = NULL, design_canonical = NULL, verbose = 1, meanTol = 1e-06, varTol = 1e-06 )
BOLD , design , nuisance
|
Session-length list of numeric matrices/arrays, each with volumes along the first dimension. |
design_canonical |
TO DO |
verbose |
|
meanTol , varTol
|
Tolerance for mean, variance and SNR of each data location.
Locations which do not meet these thresholds are masked out of the analysis.
Default: |
A "CompareGLM"
object: a list with elements
The estimated coefficients for the Bayesian model.
A mask of mesh_orig
indicating the locations inside mesh
.
The design matrix, after centering and scaling, but before any nuisance regression or prewhitening.
The names of the fields.
The names of the sessions.
Hyperparameter posterior densities.
Theta estimates from the Bayesian model.
For joint group modeling.
For joint group modeling.
For joint group modeling.
For joint group modeling: The BOLD data after any centering, scaling, nuisance regression, or prewhitening.
For joint group modeling: The design matrix after any centering, scaling, nuisance regression, or prewhitening.
Vectors of values across locations: phi
(AR coefficients averaged across sessions), sigma_sq
(residual variance averaged across sessions), and AIC (the maximum across sessions).
match.call() for this function call.
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
n_threads
n_threads |
The maximum number of threads to use for parallel
computations: prewhitening parameter estimation, and the inla-program model
estimation. Default: |
nbhd_order
nbhd_order |
For volumetric model. What order neighborhood around data
locations to keep? |
nuisance
nuisance |
(Optional) A Detrending/high-pass filtering is accomplished by adding DCT bases to the
nuisance matrix; see the parameters |
Plot design matrix
Plot design with lineplot
Plot design with imageplot
plot_design(design, method = c("lineplot", "imageplot"), ...) plot_design_line( design, colors = "Set1", linetype = "solid", linewidth = 0.7, alpha = 0.8 ) plot_design_image(design)
plot_design(design, method = c("lineplot", "imageplot"), ...) plot_design_line( design, colors = "Set1", linetype = "solid", linewidth = 0.7, alpha = 0.8 ) plot_design_image(design)
design |
The timepoints by fields design matrix or data.frame. |
method |
|
... |
Additional arguments to |
colors |
The name of a ColorBrewer palette (see
RColorBrewer::brewer.pal.info and colorbrewer2.org), the name of a
viridisLite palette, or a character vector of colors. Default:
|
linetype , linewidth , alpha
|
Parameters for |
A ggplot
A ggplot
A ggplot
view_xifti
to plot a "act_BGLM"
objectS3 method: use view_xifti
to plot a "act_BGLM"
object
## S3 method for class 'act_BGLM' plot(x, idx = NULL, title = NULL, session = NULL, ...)
## S3 method for class 'act_BGLM' plot(x, idx = NULL, title = NULL, session = NULL, ...)
x |
An object of class "act_BGLM" |
idx |
Which field should be plotted? Give the numeric indices or the
names. |
title |
If NULL, the field names associated with idx will be used. |
session |
Which session should be plotted? |
... |
Additional arguments to |
Result of the call to ciftiTools::view_cifti_surface
.
view_xifti
to plot a "BGLM"
objectS3 method: use view_xifti
to plot a "BGLM"
object
## S3 method for class 'BfMRI_design' plot(x, ...)
## S3 method for class 'BfMRI_design' plot(x, ...)
x |
An object of class "BfMRI_design". |
... |
Additional arguments to |
Result of the call to plot_design
view_xifti
to plot a "BGLM"
objectS3 method: use view_xifti
to plot a "BGLM"
object
## S3 method for class 'BGLM' plot( x, Bayes = NULL, idx = NULL, title = NULL, session = NULL, zlim = c(-1, 1), ... )
## S3 method for class 'BGLM' plot( x, Bayes = NULL, idx = NULL, title = NULL, session = NULL, zlim = c(-1, 1), ... )
x |
An object of class "BGLM" |
Bayes |
|
idx |
Which field should be plotted? Give the numeric indices or the
names. |
title |
If NULL, the field names associated with idx will be used. |
session |
Which session should be plotted? |
zlim |
Overrides the |
... |
Additional arguments to |
Result of the call to ciftiTools::view_cifti
.
view_xifti
to plot a "BGLM2"
objectS3 method: use view_xifti
to plot a "BGLM2"
object
## S3 method for class 'BGLM2' plot(x, idx = NULL, what = c("contrasts", "activations"), zlim = c(-1, 1), ...)
## S3 method for class 'BGLM2' plot(x, idx = NULL, what = c("contrasts", "activations"), zlim = c(-1, 1), ...)
x |
An object of class "BGLM2" |
idx |
Which contrast should be plotted? Give the numeric indices or the
names. |
what |
Estimates of the |
zlim |
Overrides the |
... |
Additional arguments to |
Result of the call to ciftiTools::view_cifti
.
view_xifti
to plot a "prev_BGLM"
objectS3 method: use view_xifti
to plot a "prev_BGLM"
object
## S3 method for class 'prev_BGLM' plot( x, idx = NULL, session = NULL, drop_zeros = NULL, colors = "plasma", zlim = c(0, 1), ... )
## S3 method for class 'prev_BGLM' plot( x, idx = NULL, session = NULL, drop_zeros = NULL, colors = "plasma", zlim = c(0, 1), ... )
x |
An object of class "prev_BGLM" |
idx |
Which task should be plotted? Give the numeric indices or the
names. |
session |
Which session should be plotted? |
drop_zeros |
Color locations without any activation across all results
(zero prevalence) the same color as the medial wall? Default: |
colors , zlim
|
See |
... |
Additional arguments to |
Result of the call to ciftiTools::view_cifti_surface
.
Activations prevalence.
prevalence(act_list, gamma_idx = 1)
prevalence(act_list, gamma_idx = 1)
act_list |
List of activations from |
gamma_idx |
If activations at multiple thresholds were computed, which threshold should be used for prevalence? Default: the first (lowest). |
A list containing the prevalences of activation, as a proportion of
the results from act_list
.
resamp_res
resamp_res |
For cortex spatial model. The number of vertices to which
each cortical surface should be resampled, or For computational feasibility, a value of |
return_INLA
return_INLA |
Return the INLA model object? (It can be large.) Use
|
Scale the BOLD timeseries
scale_BOLD(BOLD, scale = c("mean", "sd", "none"), v_means = NULL)
scale_BOLD(BOLD, scale = c("mean", "sd", "none"), v_means = NULL)
BOLD |
fMRI data as a locations by time ( |
scale |
Option for scaling the BOLD response. |
v_means |
Original means of the BOLD data. ONLY provide if data has already been centered. \code{"mean"} scaling will scale the data to percent local signal change. \code{"sd"} scaling will scale the data by local standard deviation. \code{"none"} will only center the data, not scale it. |
Scale to units of percent local signal change and centers
scale_BOLD
scale_BOLD |
Controls scaling the BOLD response at each location.
|
session_names
session_names |
The names of the task-fMRI |
"act_BGLM"
objectSummary method for class "act_BGLM"
## S3 method for class 'act_BGLM' summary(object, ...) ## S3 method for class 'summary.act_BGLM' print(x, ...) ## S3 method for class 'act_BGLM' print(x, ...)
## S3 method for class 'act_BGLM' summary(object, ...) ## S3 method for class 'summary.act_BGLM' print(x, ...) ## S3 method for class 'act_BGLM' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.act_BGLM"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"act_fit_bglm"
objectSummary method for class "act_fit_bglm"
## S3 method for class 'act_fit_bglm' summary(object, ...) ## S3 method for class 'summary.act_fit_bglm' print(x, ...) ## S3 method for class 'act_fit_bglm' print(x, ...)
## S3 method for class 'act_fit_bglm' summary(object, ...) ## S3 method for class 'summary.act_fit_bglm' print(x, ...) ## S3 method for class 'act_fit_bglm' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.act_fit_bglm"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"BfMRI_design"
objectSummary method for class "BfMRI_design"
## S3 method for class 'BfMRI_design' summary(object, ...) ## S3 method for class 'summary.BfMRI_design' print(x, ...) ## S3 method for class 'BfMRI_design' print(x, ...)
## S3 method for class 'BfMRI_design' summary(object, ...) ## S3 method for class 'summary.BfMRI_design' print(x, ...) ## S3 method for class 'BfMRI_design' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.BfMRI_design"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"BGLM"
objectSummary method for class "BGLM"
## S3 method for class 'BGLM' summary(object, ...) ## S3 method for class 'summary.BGLM' print(x, ...) ## S3 method for class 'BGLM' print(x, ...)
## S3 method for class 'BGLM' summary(object, ...) ## S3 method for class 'summary.BGLM' print(x, ...) ## S3 method for class 'BGLM' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.BGLM"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"BGLM2"
objectSummary method for class "BGLM2"
## S3 method for class 'BGLM2' summary(object, ...) ## S3 method for class 'summary.BGLM2' print(x, ...) ## S3 method for class 'BGLM2' print(x, ...)
## S3 method for class 'BGLM2' summary(object, ...) ## S3 method for class 'summary.BGLM2' print(x, ...) ## S3 method for class 'BGLM2' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.BGLM2"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"fit_bglm"
objectSummary method for class "fit_bglm"
## S3 method for class 'fit_bglm' summary(object, ...) ## S3 method for class 'summary.fit_bglm' print(x, ...) ## S3 method for class 'fit_bglm' print(x, ...)
## S3 method for class 'fit_bglm' summary(object, ...) ## S3 method for class 'summary.fit_bglm' print(x, ...) ## S3 method for class 'fit_bglm' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.fit_bglm"
object, a list summarizing the properties
of object
.
NULL
, invisibly.
NULL
, invisibly.
"fit_bglm2"
objectSummary method for class "fit_bglm2"
## S3 method for class 'fit_bglm2' summary(object, ...) ## S3 method for class 'summary.fit_bglm2' print(x, ...) ## S3 method for class 'fit_bglm2' print(x, ...)
## S3 method for class 'fit_bglm2' summary(object, ...) ## S3 method for class 'summary.fit_bglm2' print(x, ...) ## S3 method for class 'fit_bglm2' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.fit_bglm2"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"prev_BGLM"
objectSummary method for class "prev_BGLM"
Summary method for class "prev_BGLM"
## S3 method for class 'prev_BGLM' summary(object, ...) ## S3 method for class 'summary.prev_BGLM' print(x, ...) ## S3 method for class 'prev_BGLM' print(x, ...) ## S3 method for class 'prev_BGLM' summary(object, ...) ## S3 method for class 'summary.prev_BGLM' print(x, ...) ## S3 method for class 'prev_BGLM' print(x, ...)
## S3 method for class 'prev_BGLM' summary(object, ...) ## S3 method for class 'summary.prev_BGLM' print(x, ...) ## S3 method for class 'prev_BGLM' print(x, ...) ## S3 method for class 'prev_BGLM' summary(object, ...) ## S3 method for class 'summary.prev_BGLM' print(x, ...) ## S3 method for class 'prev_BGLM' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.prev_BGLM"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
A "summary.prev_BGLM"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
"prev_fit_bglm"
objectSummary method for class "prev_fit_bglm"
Summary method for class "prev_fit_bglm"
## S3 method for class 'prev_fit_bglm' summary(object, ...) ## S3 method for class 'summary.prev_fit_bglm' print(x, ...) ## S3 method for class 'prev_fit_bglm' print(x, ...) ## S3 method for class 'prev_fit_bglm' summary(object, ...) ## S3 method for class 'summary.prev_fit_bglm' print(x, ...) ## S3 method for class 'prev_fit_bglm' print(x, ...)
## S3 method for class 'prev_fit_bglm' summary(object, ...) ## S3 method for class 'summary.prev_fit_bglm' print(x, ...) ## S3 method for class 'prev_fit_bglm' print(x, ...) ## S3 method for class 'prev_fit_bglm' summary(object, ...) ## S3 method for class 'summary.prev_fit_bglm' print(x, ...) ## S3 method for class 'prev_fit_bglm' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
Object of class |
A "summary.prev_fit_bglm"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
A "summary.prev_fit_bglm"
object, a list summarizing the
properties of object
.
NULL
, invisibly.
NULL
, invisibly.
surfaces
surfL , surfR
|
For cortex spatial model. Left and right cortex surface
geometry in GIFTI format ("*.surf.gii"). These can be a file path to
a GIFTI file or a Surfaces can alternatively be provided through the |
trim_INLA
trim_INLA |
(logical) should the |
verbose
verbose |
|
Compute surface areas of each vertex in a triangular mesh.
vertex_areas(mesh)
vertex_areas(mesh)
mesh |
An |
Vector of areas
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.
vertices
vertices |
A |
Construct a triangular mesh from a 3D volumetric mask
vol2spde(mask, res, nbhd_order = 1, buffer = c(1, 1, 3, 4, 4))
vol2spde(mask, res, nbhd_order = 1, buffer = c(1, 1, 3, 4, 4))
mask |
An array of 0s and 1s representing a volumetric mask |
res |
The spatial resolution in each direction, in mm. For example, c(2,2,2) indicates 2mm isotropic voxels. |
nbhd_order |
For volumetric data, what order neighborhood around data locations to keep? (0 = no neighbors, 1 = 1st-order neighbors, 2 = 1st- and 2nd-order neighbors, etc.). Smaller values will provide greater computational efficiency at the cost of higher variance around the edge of the data. |
buffer |
For volumetric data, size of extra voxels layers around the bounding box, in terms of voxels. Set to NULL for no buffer. |
An inla.spde2 object.
This function requires the INLA
package, which is not a CRAN package.
See https://www.r-inla.org/download-install for easy installation instructions.