Title: | Estimate Brain Networks and Connectivity with ICA and Empirical Priors |
---|---|
Description: | Implements the template ICA (independent components analysis) model proposed in Mejia et al. (2020) <doi:10.1080/01621459.2019.1679638> and the spatial template ICA model proposed in proposed in Mejia et al. (2022) <doi:10.1080/10618600.2022.2104289>. Both models estimate subject-level brain as deviations from known population-level networks, which are estimated using standard ICA algorithms. Both models employ an expectation-maximization algorithm for estimation of the latent brain networks and unknown model parameters. Includes direct support for 'CIFTI', 'GIFTI', and 'NIFTI' neuroimaging file formats. |
Authors: | Amanda Mejia [aut, cre], Damon Pham [aut] , Daniel Spencer [ctb] , Mary Beth Nebel [ctb] |
Maintainer: | Amanda Mejia <[email protected]> |
License: | GPL-3 |
Version: | 0.9.1 |
Built: | 2024-11-25 16:44:49 UTC |
Source: | CRAN |
Identify areas of activation in each independent component map from the result of (spatial) template ICA.
activations( tICA, u = NULL, z = NULL, alpha = 0.01, type = c(">", "abs >", "<", "!="), method_p = "BH", verbose = FALSE, which.ICs = NULL, deviation = FALSE )
activations( tICA, u = NULL, z = NULL, alpha = 0.01, type = c(">", "abs >", "<", "!="), method_p = "BH", verbose = FALSE, which.ICs = NULL, deviation = FALSE )
tICA |
Fitted (spatial) template ICA object from |
u , z
|
Set a threshold value for activation? A threshold value can be
specified directly with |
alpha |
Significance level for hypothesis testing. Default: |
type |
Type of region: |
method_p |
If the input is a |
verbose |
If |
which.ICs |
Indices of ICs for which to identify activations. If
|
deviation |
If |
A list containing activation maps for each IC, the joint and marginal PPMs for each IC, and the parameters used for computing activation. If the input represented CIFTI- or NIFTI-format data, then the activations maps will be formatted accordingly.
Use summary
to obtain information about the activations results.
For CIFTI-format activations, use plot
to visualize the activation
maps.
## Not run: activations(tICA_result, alpha=.05, deviation=TRUE) ## End(Not run)
## Not run: activations(tICA_result, alpha=.05, deviation=TRUE) ## End(Not run)
Bdiag m2
bdiag_m2(mat, N)
bdiag_m2(mat, N)
mat |
a k x k 'matrix' |
N |
how many times to repeat |
a sparse (Nk x Nk) matrix of class "dgCMatrix"
.
Cholesky-based FC sampling
Chol_samp_fun(Chol_vals, p, M, chol_diag, chol_offdiag, Chol_mat_blank)
Chol_samp_fun(Chol_vals, p, M, chol_diag, chol_offdiag, Chol_mat_blank)
Chol_vals |
Matrix of Cholesky factorizations (upper triangular values) for all template sessions (nN*nM x nChol) |
p |
Pivot/reordering applied to FC matrices prior to Cholesky factorization |
M |
Number of samples to draw |
chol_diag |
Indices of diagonal upper triangular elements |
chol_offdiag |
Indices of off-diagonal upper triangular elements |
Chol_mat_blank |
A nLxnL matrix indexing the upper triangular elements |
Performs dimension reduction and prewhitening based on probabilistic PCA using SVD. If dimensionality is not specified, it is estimated using the method described in Minka (2008).
dim_reduce(X, Q = NULL, Q_max = 100)
dim_reduce(X, Q = NULL, Q_max = 100)
X |
A numeric matrix, with each column being a centered timeseries.
For fMRI data, |
Q |
Number of latent dimensions to estimate. If |
Q_max |
Maximal number of principal components for automatic
dimensionality selection with PESEL. Default: |
A list containing the dimension-reduced data (data_reduced
, a
matrix), prewhitening/dimension reduction matrix (
H
,
a matrix) and its (pseudo-)inverse (
Hinv
, a
matrix), the noise variance (
sigma_sq
), the correlation matrix of the
dimension-reduced data (C_diag
, a matrix), and the
dimensionality (
Q
).
nT <- 30 nV <- 400 nQ <- 7 X <- matrix(rnorm(nV*nQ), nrow=nV) %*% diag(seq(nQ, 1)) %*% matrix(rnorm(nQ*nT), nrow=nQ) dim_reduce(X, Q=nQ)
nT <- 30 nV <- 400 nQ <- 7 X <- matrix(rnorm(nV*nQ), nrow=nV) %*% diag(seq(nQ, 1)) %*% matrix(rnorm(nQ*nT), nrow=nQ) dim_reduce(X, Q=nQ)
EM Algorithms for Template ICA Models
EM_templateICA.spatial( template_mean, template_var, meshes, BOLD, theta0, C_diag, H, Hinv, maxiter = 100, usePar = FALSE, epsilon = 0.001, reduce_dim = TRUE, verbose = FALSE ) EM_templateICA.independent( template_mean, template_var, BOLD, theta0, C_diag, H, Hinv, maxiter = 100, epsilon = 0.001, reduce_dim = FALSE, usePar = FALSE, verbose )
EM_templateICA.spatial( template_mean, template_var, meshes, BOLD, theta0, C_diag, H, Hinv, maxiter = 100, usePar = FALSE, epsilon = 0.001, reduce_dim = TRUE, verbose = FALSE ) EM_templateICA.independent( template_mean, template_var, BOLD, theta0, C_diag, H, Hinv, maxiter = 100, epsilon = 0.001, reduce_dim = FALSE, usePar = FALSE, verbose )
template_mean |
( |
template_var |
( |
meshes |
|
BOLD |
( |
theta0 |
(list) initial guess at parameter values: A ( |
C_diag |
( |
H , Hinv
|
For dimension reduction
of the spatial template ICA model, which assumes that all IC's have the
same smoothness parameter, |
maxiter |
Maximum number of EM iterations. Default: 100. |
usePar |
Parallelize the computation? Default: |
epsilon |
Smallest proportion change between iterations. Default: 0.001. |
reduce_dim |
Reduce the temporal dimension of the data using PCA?
Default: |
verbose |
If |
EM_templateICA.spatial
implements the expectation-maximization
(EM) algorithm described in Mejia et al. (2019+) for estimating the
subject-level ICs and unknown parameters in the template ICA model with
spatial priors on subject effects.
In both models, if original fMRI timeseries has covariance
, the prewhitened timeseries achieved by premultiplying
by (
) matrix
from PCA has diagonal covariance
, so C_diag is
.
A list: theta (list of final parameter estimates), subICmean
(estimates of subject-level ICs), subICvar (variance of subject-level ICs,
for non-spatial model) or subjICcov (covariance matrix of subject-level ICs,
for spatial model – note that only diagonal and values for neighbors are
computed), and success (flag indicating convergence (TRUE
) or not
(FALSE
))
Universally estimate IW dof parameter nu based on method of moments, so that no empirical variance is under-estimated
estimate_nu(var_FC, mean_FC)
estimate_nu(var_FC, mean_FC)
var_FC |
Empirical between-subject variance of covariance matrices (QxQ) |
mean_FC |
Empirical mean of covariance matrices (QxQ) |
estimate for nu
Estimate IW dof parameter nu based on method of moments
estimate_nu_matrix(var_FC, mean_FC)
estimate_nu_matrix(var_FC, mean_FC)
var_FC |
Empirical between-subject variance of covariance matrices (QxQ) |
mean_FC |
Empirical mean of covariance matrices (QxQ) |
QxQ matrix of estimates for nu
Estimate template for Template ICA based on fMRI data
estimate_template( BOLD, BOLD2 = NULL, GICA, mask = NULL, inds = NULL, scale = c("local", "global", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, brainstructures = "all", resamp_res = NULL, keep_S = FALSE, keep_FC = FALSE, FC = TRUE, FC_nPivots = 100, FC_nSamp = 50000, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE ) estimate_template.cifti( BOLD, BOLD2 = NULL, GICA, inds = NULL, scale = c("local", "global", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, brainstructures = "all", resamp_res = resamp_res, keep_S = FALSE, keep_FC = FALSE, FC = FALSE, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE ) estimate_template.gifti( BOLD, BOLD2 = NULL, GICA, inds = NULL, scale = c("local", "global", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, brainstructures = "all", keep_S = FALSE, keep_FC = FALSE, FC = FALSE, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE ) estimate_template.nifti( BOLD, BOLD2 = NULL, GICA, inds = NULL, scale = c("local", "global", "none"), TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, mask = NULL, keep_S = FALSE, keep_FC = FALSE, FC = FALSE, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE )
estimate_template( BOLD, BOLD2 = NULL, GICA, mask = NULL, inds = NULL, scale = c("local", "global", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, brainstructures = "all", resamp_res = NULL, keep_S = FALSE, keep_FC = FALSE, FC = TRUE, FC_nPivots = 100, FC_nSamp = 50000, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE ) estimate_template.cifti( BOLD, BOLD2 = NULL, GICA, inds = NULL, scale = c("local", "global", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, brainstructures = "all", resamp_res = resamp_res, keep_S = FALSE, keep_FC = FALSE, FC = FALSE, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE ) estimate_template.gifti( BOLD, BOLD2 = NULL, GICA, inds = NULL, scale = c("local", "global", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, brainstructures = "all", keep_S = FALSE, keep_FC = FALSE, FC = FALSE, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE ) estimate_template.nifti( BOLD, BOLD2 = NULL, GICA, inds = NULL, scale = c("local", "global", "none"), TR = NULL, hpf = 0.01, GSR = FALSE, Q2 = 0, Q2_max = NULL, mask = NULL, keep_S = FALSE, keep_FC = FALSE, FC = FALSE, varTol = 1e-06, maskTol = 0.1, missingTol = 0.1, usePar = FALSE, wb_path = NULL, verbose = TRUE )
BOLD , BOLD2
|
Vector of subject-level fMRI data in one of the following
formats: CIFTI file paths, If |
GICA |
Group ICA maps in a format compatible with New: can also be a parcellation in CIFTI format (other formats to be implemented in the future). The parcellation should have the same locations and one column, with integer values indicating the parcel to which each location belongs to. Each parcel is modeled as a brain map; instead of the first step of dual regression, the medial timecourse of each parcel is used. |
mask |
Required if |
inds |
Numeric indices of the group ICs to include in the template. If
If |
scale |
|
scale_sm_surfL , scale_sm_surfR , scale_sm_FWHM
|
Only applies if
If If To create a |
TR |
The temporal resolution of the data, i.e. the time between volumes,
in seconds. |
hpf |
The frequency at which to apply a highpass filter to the data
during pre-processing, in Hertz. Default: The highpass filter serves to detrend the data, since low-frequency variance is associated with noise. Highpass filtering is accomplished by nuisance regression of discrete cosine transform (DCT) bases. Note the |
GSR |
Center BOLD across columns (each image)? This
is equivalent to performing global signal regression. Default:
|
Q2 , Q2_max
|
Obtain dual regression estimates after denoising? Denoising is based on modeling and removing nuisance ICs. It may result in a cleaner estimate for smaller datasets, but it may be unnecessary (and time-consuming) for larger datasets. Set If |
brainstructures |
Only applies if the entries of |
resamp_res |
Only applies if the entries of |
keep_S |
Keep the DR estimates of S? If |
keep_FC |
Keep the DR estimates of the FC cor(A)? If |
FC |
Include the functional connectivity template? Default: |
FC_nPivots |
Number of pivots to use in Cholesky-based FC template estimation. Set to zero to skip Cholesky-based FC template estimation. Default: 100. |
FC_nSamp |
Number of FC matrix samples to generate across all pivots. This should be a multiple of FC_nPivots. |
varTol |
Tolerance for variance of each data location. For each scan,
locations which do not meet this threshold are masked out of the analysis.
Default: |
maskTol |
For computing the dual regression results for each subject:
tolerance for number of locations masked out due to low
variance or missing values. If more than this many locations are masked out,
a subject is skipped without calculating dual regression. If |
missingTol |
For computing the variance decomposition across all subjects:
tolerance for number of subjects masked out due to low variance or missing
values at a given location. If more than this many subjects are masked out,
the location's value will be |
usePar , wb_path
|
Parallelize the DR computations over subjects? Default:
|
verbose |
Display progress updates? Default: |
All fMRI data (entries in BOLD
and BOLD2
, and GICA
) must
be in the same spatial resolution.
A list: the template
and var_decomp
with entries in
matrix format; the mask
of locations without template values due to
too many low variance or missing values; the function params
such as
the type of scaling and detrending performed; the dat_struct
which can be
used to convert template
and var_decomp
to "xifti"
or
"nifti"
objects if the BOLD
format was CIFTI or NIFTI data;
and DR results if isTRUE(keep_S)
and/or isTRUE(keep_FC)
.
Use summary
to print a description of the template results, and
for CIFTI-format data use plot
to plot the template mean and variance
estimates. Use export_template
to save the templates to
individual RDS, CIFTI, or NIFTI files (depending on the BOLD
format).
nT <- 21 nV <- 140 nQ <- 6 mU <- matrix(rnorm(nV*nQ), nrow=nV) mS <- mU %*% diag(seq(nQ, 1)) %*% matrix(rnorm(nQ*nT), nrow=nQ) BOLD <- list(B1=mS, B2=mS, B3=mS) BOLD <- lapply(BOLD, function(x){x + rnorm(nV*nT, sd=.05)}) GICA <- mU estimate_template(BOLD=BOLD, GICA=mU, FC_nSamp=2000) ## Not run: estimate_template( run1_cifti_fnames, run2_cifti_fnames, gICA_cifti_fname, brainstructures="all", scale="global", TR=0.71, Q2=NULL, varTol=10 ) ## End(Not run)
nT <- 21 nV <- 140 nQ <- 6 mU <- matrix(rnorm(nV*nQ), nrow=nV) mS <- mU %*% diag(seq(nQ, 1)) %*% matrix(rnorm(nQ*nT), nrow=nQ) BOLD <- list(B1=mS, B2=mS, B3=mS) BOLD <- lapply(BOLD, function(x){x + rnorm(nV*nT, sd=.05)}) GICA <- mU estimate_template(BOLD=BOLD, GICA=mU, FC_nSamp=2000) ## Not run: estimate_template( run1_cifti_fnames, run2_cifti_fnames, gICA_cifti_fname, brainstructures="all", scale="global", TR=0.71, Q2=NULL, varTol=10 ) ## End(Not run)
Estimate FC template
estimate_template_FC(FC0, nu_adjust = 1)
estimate_template_FC(FC0, nu_adjust = 1)
FC0 |
The FC estimates from |
nu_adjust |
Factor by which to adjust estimate of nu. Values < 1 will inflate the template variance to avoid an over-informative prior on FC. |
Estimate variance decomposition and templates from DR estimates.
estimate_template_from_DR(DR, LV = NULL)
estimate_template_from_DR(DR, LV = NULL)
DR |
the test/retest dual regression estimates, as an array with
dimensions ( |
LV |
A length-two integer vector giving the dimensions |
List of two elements: the templates and the variance decomposition.
There are two version of the variance template: varUB
gives the
unbiased variance estimate, and varNN
gives the upwardly-biased
non-negative variance estimate. Values in varUB
will need to be
clamped above zero before using in templateICA
.
Estimation of effective sample size
estimate.ESS(mesh, Y, ind = NULL, trace = FALSE)
estimate.ESS(mesh, Y, ind = NULL, trace = FALSE)
mesh |
INLA mesh |
Y |
data |
ind |
index of the data locations in the mesh |
trace |
If |
The functions computes the effective sample size as
as in Bretherton et al. (1999), Journal of Climate.
Estimate of the effective sample size
Export the templates (mean and variance) as separate files for
visualization or processing outside of templateICAr
.
export_template( x, out_fname = NULL, var_method = c("non-negative", "unbiased") )
export_template( x, out_fname = NULL, var_method = c("non-negative", "unbiased") )
x |
The result of |
out_fname |
Use |
var_method |
|
If is.null(out_fname)
, the templates in data matrix,
"xifti"
, or "nifti"
format, to match the format of the
original BOLD data. Otherwise, the paths to the new files specified by
out_fname
. If template includes functional connectivity components,
the FC template and its mean and variance will be included.
## Not run: tm <- estimate_template(cii1_fnames, cii2_fnames, gICA_fname) export_template(tm, out_fname="my_template", var_method="unbiased") ## End(Not run)
## Not run: tm <- estimate_template(cii1_fnames, cii2_fnames, gICA_fname) export_template(tm, out_fname="my_template", var_method="unbiased") ## End(Not run)
Compute inverse covariance matrix for AR process (up to a constant scaling factor)
getInvCovAR(ar, ntime)
getInvCovAR(ar, ntime)
ar |
vector of p AR parameters |
ntime |
number of time points in timeseries |
inverse covariance matrix for AR process (up to a constant scaling factor)
Perform group ICA based on CIFTI data
groupICA.cifti( cifti_fnames, subjects = NULL, brainstructures = "all", num_PCs = 100, num_ICs, max_rows_GPCA = 10000, GSR = FALSE, scale = c("local", "global", "none"), scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, verbose = TRUE, out_fname = NULL, surfL = NULL, surfR = NULL, smooth = TRUE, smooth_surf_FWHM = 5, smooth_vol_FWHM = 5, smooth_zeroes_as_NA = FALSE, smooth_subcortical_merged = FALSE )
groupICA.cifti( cifti_fnames, subjects = NULL, brainstructures = "all", num_PCs = 100, num_ICs, max_rows_GPCA = 10000, GSR = FALSE, scale = c("local", "global", "none"), scale_sm_FWHM = 2, TR = NULL, hpf = 0.01, verbose = TRUE, out_fname = NULL, surfL = NULL, surfR = NULL, smooth = TRUE, smooth_surf_FWHM = 5, smooth_vol_FWHM = 5, smooth_zeroes_as_NA = FALSE, smooth_subcortical_merged = FALSE )
cifti_fnames |
Vector of file paths of CIFTI-format fMRI timeseries (*.dtseries.nii) for which to calculate group ICA |
subjects |
Use this argument if some subjects have more than one scan.
This is a numeric or factor vector the same length as |
brainstructures |
Character vector indicating which brain structure(s)
to obtain: |
num_PCs |
Maximum number of PCs to retain in initial subject-level PCA |
num_ICs |
Number of independent components to identify. |
max_rows_GPCA |
The maximum number of rows for the matrix on which group
level PCA will be performed. This matrix is the result of temporal concatenation
and contains |
GSR , scale , scale_sm_FWHM , TR , hpf
|
Center BOLD columns, scale by the
standard deviation, and detrend voxel timecourses? See
|
verbose |
If |
out_fname |
(Optional) If not specified, a xifti object will be returned but the GICA maps will not be written to a CIFTI file. |
surfL |
(Optional) File path to a surface GIFTI for the left cortex. If provided, will be part of xifti result object for visualization in R. Will also be used to perform any smoothing. |
surfR |
(Optional) File path to a surface GIFTI for the right cortex. If provided, will be part of xifti result object for visualization in R. Will also be used to perform any smoothing. |
smooth |
Smooth the CIFTI data prior to reading in each file? Default:
|
smooth_surf_FWHM , smooth_vol_FWHM , smooth_zeroes_as_NA , smooth_subcortical_merged
|
See |
out_fname
if a file was written, or the GICA as a "xifti"
object
if not.
Compute theoretical Inverse-Wishart variance of covariance matrix elements
IW_var(nu, p, xbar_ij, xbar_ii, xbar_jj)
IW_var(nu, p, xbar_ij, xbar_ii, xbar_jj)
nu |
Inverse Wishart degrees of freedom parameter |
p |
Matrix dimension for IW distribution |
xbar_ij |
Empirical mean of covariance matrices at element (i,j) |
xbar_ii |
Empirical mean of covariance matrices at the ith diagonal element |
xbar_jj |
Empirical mean of covariance matrices at the jth diagonal element |
Theoretical IW variance for covariance element (i,j)
Compute theoretical Inverse-Wishart variance of correlation matrix elements
IW_var_cor(nu, p, xbar_ij)
IW_var_cor(nu, p, xbar_ij)
nu |
Inverse Wishart degrees of freedom parameter |
p |
Matrix dimension for IW distribution |
xbar_ij |
Empirical mean of covariance matrices at element (i,j) |
Theoretical IW variance for correlation element (i,j)
Compute likelihood in SPDE model for ESS estimation
lik(theta, Y, G, C, ind = NULL)
lik(theta, Y, G, C, ind = NULL)
theta |
Value of hyperparameters |
Y |
Data vector |
G |
SPDE G matrix |
C |
SPDE C matrix |
ind |
Indices of data locations in the mesh |
Log likelihood value
"surf"
objectCreate INLA mesh and observation weight matrix based on a "surf"
object
make_mesh(surf = NULL, inds_data = NULL, inds_mesh = NULL)
make_mesh(surf = NULL, inds_data = NULL, inds_mesh = NULL)
surf |
|
inds_data |
Subset of vertices to include in analysis, e.g. non-medial wall locations. |
inds_mesh |
Subset of vertices to retain in mesh, e.g. non-medial wall
locations. Must be a superset of |
List containing INLA mesh, observation weight matrix A for translating between mesh locations and original data locations, the brain mask used to create the mesh, and the number of original and mesh data locations.
Create INLA mesh and observation weight matrix based on a binary brain mask
make_mesh_2D(mask)
make_mesh_2D(mask)
mask |
Brain mask (matrix of 0 and 1 or |
This function requires the INLA
package, which is not a CRAN
package. See https://www.r-inla.org/download-install for easy
installation instructions.
List containing INLA mesh, observation weight matrix A for translating between mesh locations and original data locations, the brain mask used to create the mesh, and the number of original and mesh data locations.
Center the data across space and/or time, detrend, and scale, in that order. For dual regression, row centering is required and column centering is not recommended. Scaling and detrending depend on the user preference.
norm_BOLD( BOLD, center_rows = TRUE, center_cols = FALSE, scale = c("local", "global", "none"), scale_sm_xifti = NULL, scale_sm_FWHM = 2, scale_sm_xifti_mask = NULL, TR = NULL, hpf = 0.01 )
norm_BOLD( BOLD, center_rows = TRUE, center_cols = FALSE, scale = c("local", "global", "none"), scale_sm_xifti = NULL, scale_sm_FWHM = 2, scale_sm_xifti_mask = NULL, TR = NULL, hpf = 0.01 )
BOLD |
fMRI numeric data matrix ( |
center_rows , center_cols
|
Center BOLD data across rows (each data
location's time series) or columns (each time point's image)? Default:
|
scale |
|
scale_sm_xifti , scale_sm_FWHM
|
Only applies if |
scale_sm_xifti_mask |
For local scaling with smoothing, the data must be unmasked to be mapped back to the surface. So if the data are masked, provide the mask here. |
TR |
The temporal resolution of the data, i.e. the time between volumes,
in seconds. |
hpf |
The frequency at which to apply a highpass filter to the data
during pre-processing, in Hertz. Default: The highpass filter serves to detrend the data, since low-frequency variance is associated with noise. Highpass filtering is accomplished by nuisance regression of discrete cosine transform (DCT) bases. Note the |
Normalized BOLD data matrix ()
Orthonormalizes a square, invertible matrix
orthonorm(X)
orthonorm(X)
X |
A square matrix to be orthonormalized. |
Y is orthonormal if $YY'=Y'Y=I$. Orthonormalization of X is given by $X (X'X)^(-.5)$.
X after orthonormalization
Plot template
## S3 method for class 'template.cifti' plot( x, stat = c("both", "mean", "sd", "var"), var_method = c("non-negative", "unbiased"), ... )
## S3 method for class 'template.cifti' plot( x, stat = c("both", "mean", "sd", "var"), var_method = c("non-negative", "unbiased"), ... )
x |
The template from |
stat |
|
var_method |
|
... |
Additional arguments to |
The plot
Plot template
## S3 method for class 'template.gifti' plot( x, stat = c("both", "mean", "sd", "var"), var_method = c("non-negative", "unbiased"), ... )
## S3 method for class 'template.gifti' plot( x, stat = c("both", "mean", "sd", "var"), var_method = c("non-negative", "unbiased"), ... )
x |
The template from |
stat |
|
var_method |
|
... |
Additional arguments to |
The plot
Plot template
## S3 method for class 'template.matrix' plot(x, ...)
## S3 method for class 'template.matrix' plot(x, ...)
x |
The template from |
... |
Additional arguments |
The plot
Based on oro.nifti::image
.
## S3 method for class 'template.nifti' plot( x, stat = c("mean", "sd", "var"), plane = c("axial", "sagittal", "coronal"), n_slices = 9, slices = NULL, var_method = c("non-negative", "unbiased"), ... )
## S3 method for class 'template.nifti' plot( x, stat = c("mean", "sd", "var"), plane = c("axial", "sagittal", "coronal"), n_slices = 9, slices = NULL, var_method = c("non-negative", "unbiased"), ... )
x |
The template from |
stat |
|
plane , n_slices , slices
|
Anatomical plane and which slice indices to show. Default: 9 axial slices. |
var_method |
|
... |
Additional arguments to |
Consider using struct_template
to obtain the 3D volumes to plot with a different
viewer function (e.g. from oro.nifti
) if desired.
The plot
Plot activations
## S3 method for class 'tICA_act.cifti' plot(x, stat = c("active", "pvals", "pvals_adj", "tstats", "se"), ...)
## S3 method for class 'tICA_act.cifti' plot(x, stat = c("active", "pvals", "pvals_adj", "tstats", "se"), ...)
x |
The activations from |
stat |
|
... |
Additional arguments to |
The activations plot
Plot template
## S3 method for class 'tICA.cifti' plot(x, stat = c("mean", "se", "both"), ...)
## S3 method for class 'tICA.cifti' plot(x, stat = c("mean", "se", "both"), ...)
x |
The result of |
stat |
|
... |
Additional arguments to |
The plot
This feature is not supported yet.
## S3 method for class 'tICA.matrix' plot(x, ...)
## S3 method for class 'tICA.matrix' plot(x, ...)
x |
The result of |
... |
Additional arguments |
Nothing, because an error is raised.
Plot template
## S3 method for class 'tICA.nifti' plot( x, stat = c("mean", "se"), plane = c("axial", "sagittal", "coronal"), n_slices = 9, slices = NULL, ... )
## S3 method for class 'tICA.nifti' plot( x, stat = c("mean", "se"), plane = c("axial", "sagittal", "coronal"), n_slices = 9, slices = NULL, ... )
x |
The result of |
stat |
|
plane , n_slices , slices
|
Anatomical plane and which slice indices to show. Default: 9 axial slices. |
... |
Additional arguments |
The plot
Estimate residual autocorrelation for prewhitening
pw_estimate(A, ar_order, aic = FALSE)
pw_estimate(A, ar_order, aic = FALSE)
A |
Estimated A matrix (T x Q) |
ar_order , aic
|
Order of the AR model used to prewhiten the data at each location.
If |
Estimated AR coefficients and residual variance at every vertex
Resample a CIFTI template to a new spatial resolution.
resample_template(x, resamp_res, verbose = FALSE)
resample_template(x, resamp_res, verbose = FALSE)
x |
The |
resamp_res |
The new resampling resolution. |
verbose |
Give occasional updates? Default: |
The resampled "template.cifti"
object.
Compute matrix square root of X'X
sqrt_XtX(X, inverse = FALSE)
sqrt_XtX(X, inverse = FALSE)
X |
A numerical matrix |
inverse |
if inverse=TRUE, compute inverse of square root |
A matrix equalling the (inverse) matrix square root of X'X
"template.cifti"
objectSummary method for class "template.cifti"
## S3 method for class 'template.cifti' summary(object, ...) ## S3 method for class 'summary.template.cifti' print(x, ...) ## S3 method for class 'template.cifti' print(x, ...)
## S3 method for class 'template.cifti' summary(object, ...) ## S3 method for class 'summary.template.cifti' print(x, ...) ## S3 method for class 'template.cifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The template from |
A list summarizing the template: data dimensions, options used for template estimation, etc.
Nothing, invisibly.
Nothing, invisibly.
"template.gifti"
objectSummary method for class "template.gifti"
## S3 method for class 'template.gifti' summary(object, ...) ## S3 method for class 'summary.template.gifti' print(x, ...) ## S3 method for class 'template.gifti' print(x, ...)
## S3 method for class 'template.gifti' summary(object, ...) ## S3 method for class 'summary.template.gifti' print(x, ...) ## S3 method for class 'template.gifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The template from |
A list summarizing the template: data dimensions, options used for template estimation, etc.
Nothing, invisibly.
Nothing, invisibly.
"template.matrix"
objectSummary method for class "template.matrix"
## S3 method for class 'template.matrix' summary(object, ...) ## S3 method for class 'summary.template.matrix' print(x, ...) ## S3 method for class 'template.matrix' print(x, ...)
## S3 method for class 'template.matrix' summary(object, ...) ## S3 method for class 'summary.template.matrix' print(x, ...) ## S3 method for class 'template.matrix' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The template from |
A list summarizing the template: data dimensions, options used for template estimation, etc.
Nothing, invisibly.
Nothing, invisibly.
"template.nifti"
objectSummary method for class "template.nifti"
## S3 method for class 'template.nifti' summary(object, ...) ## S3 method for class 'summary.template.nifti' print(x, ...) ## S3 method for class 'template.nifti' print(x, ...)
## S3 method for class 'template.nifti' summary(object, ...) ## S3 method for class 'summary.template.nifti' print(x, ...) ## S3 method for class 'template.nifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The template from |
A list summarizing the template: data dimensions, options used for template estimation, etc.
Nothing, invisibly.
Nothing, invisibly.
"tICA_act.cifti"
objectSummary method for class "tICA_act.cifti"
## S3 method for class 'tICA_act.cifti' summary(object, ...) ## S3 method for class 'summary.tICA_act.cifti' print(x, ...) ## S3 method for class 'tICA_act.cifti' print(x, ...)
## S3 method for class 'tICA_act.cifti' summary(object, ...) ## S3 method for class 'summary.tICA_act.cifti' print(x, ...) ## S3 method for class 'tICA_act.cifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The activations from |
A list summarizing the data and results for the activations analysis.
Nothing, invisibly.
Nothing, invisibly.
"tICA_act.matrix"
objectSummary method for class "tICA_act.matrix"
## S3 method for class 'tICA_act.matrix' summary(object, ...) ## S3 method for class 'summary.tICA_act.matrix' print(x, ...) ## S3 method for class 'tICA_act.matrix' print(x, ...)
## S3 method for class 'tICA_act.matrix' summary(object, ...) ## S3 method for class 'summary.tICA_act.matrix' print(x, ...) ## S3 method for class 'tICA_act.matrix' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The activations from |
A list summarizing the data and results for the activations analysis.
Nothing, invisibly.
Nothing, invisibly.
"tICA_act.nifti"
objectSummary method for class "tICA_act.nifti"
## S3 method for class 'tICA_act.nifti' summary(object, ...) ## S3 method for class 'summary.tICA_act.nifti' print(x, ...) ## S3 method for class 'tICA_act.nifti' print(x, ...)
## S3 method for class 'tICA_act.nifti' summary(object, ...) ## S3 method for class 'summary.tICA_act.nifti' print(x, ...) ## S3 method for class 'tICA_act.nifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The activations from |
A list summarizing the data and results for the activations analysis.
Nothing, invisibly.
Nothing, invisibly.
"tICA.cifti"
objectSummary method for class "tICA.cifti"
## S3 method for class 'tICA.cifti' summary(object, ...) ## S3 method for class 'summary.tICA.cifti' print(x, ...) ## S3 method for class 'tICA.cifti' print(x, ...)
## S3 method for class 'tICA.cifti' summary(object, ...) ## S3 method for class 'summary.tICA.cifti' print(x, ...) ## S3 method for class 'tICA.cifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The result of |
A list summarizing of the results of the templateICA analysis.
Nothing, invisibly.
Nothing, invisibly.
"tICA.matrix"
objectSummary method for class "tICA.matrix"
## S3 method for class 'tICA.matrix' summary(object, ...) ## S3 method for class 'summary.tICA.matrix' print(x, ...) ## S3 method for class 'tICA.matrix' print(x, ...)
## S3 method for class 'tICA.matrix' summary(object, ...) ## S3 method for class 'summary.tICA.matrix' print(x, ...) ## S3 method for class 'tICA.matrix' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The template from |
A list summarizing of the results of the templateICA analysis.
Nothing, invisibly.
Nothing, invisibly.
"tICA.nifti"
objectSummary method for class "tICA.nifti"
## S3 method for class 'tICA.nifti' summary(object, ...) ## S3 method for class 'summary.tICA.nifti' print(x, ...) ## S3 method for class 'tICA.nifti' print(x, ...)
## S3 method for class 'tICA.nifti' summary(object, ...) ## S3 method for class 'summary.tICA.nifti' print(x, ...) ## S3 method for class 'tICA.nifti' print(x, ...)
object |
Object of class |
... |
further arguments passed to or from other methods. |
x |
The template from |
A list summarizing of the results of the templateICA analysis.
Nothing, invisibly.
Nothing, invisibly.
Perform template independent component analysis (ICA) using variational Bayes (VB) or expectation-maximization (EM).
templateICA( BOLD, template, tvar_method = c("non-negative", "unbiased"), scale = c("template", "global", "local", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = "template", nuisance = NULL, scrub = NULL, TR = NULL, hpf = "template", GSR = "template", Q2 = "template", Q2_max = "template", brainstructures = "template", mask = NULL, time_inds = NULL, varTol = "template", spatial_model = NULL, resamp_res = NULL, rm_mwall = TRUE, reduce_dim = FALSE, method_FC = c("VB1", "VB2", "none"), maxiter = 100, miniter = 3, epsilon = 0.001, kappa_init = 0.2, usePar = TRUE, PW = FALSE, verbose = TRUE )
templateICA( BOLD, template, tvar_method = c("non-negative", "unbiased"), scale = c("template", "global", "local", "none"), scale_sm_surfL = NULL, scale_sm_surfR = NULL, scale_sm_FWHM = "template", nuisance = NULL, scrub = NULL, TR = NULL, hpf = "template", GSR = "template", Q2 = "template", Q2_max = "template", brainstructures = "template", mask = NULL, time_inds = NULL, varTol = "template", spatial_model = NULL, resamp_res = NULL, rm_mwall = TRUE, reduce_dim = FALSE, method_FC = c("VB1", "VB2", "none"), maxiter = 100, miniter = 3, epsilon = 0.001, kappa_init = 0.2, usePar = TRUE, PW = FALSE, verbose = TRUE )
BOLD |
Vector of subject-level fMRI data in one of the following
formats: CIFTI file paths, If multiple BOLD data are provided, they will be independently centered, scaled, detrended (if applicable), and denoised (if applicable). Then they will be concatenated together followed by computing the initial dual regression estimate. |
template |
Template estimates in a format compatible with |
tvar_method |
Which calculation of the template variance to use:
|
scale |
|
scale_sm_surfL , scale_sm_surfR , scale_sm_FWHM
|
Only applies if
If If To create a |
nuisance |
(Optional) Signals to regress from the data, given as a
numeric matrix with the same number of rows as there are volumes in the
|
scrub |
(Optional) A numeric vector of integers indicating the indices
of volumes to scrub from the BOLD data. (List the volumes to remove, not the
ones to keep.) If multiple |
TR , hpf
|
These arguments control detrending. Note that if multiple |
GSR |
Center BOLD across columns (each image)? This
is equivalent to performing global signal regression. Default:
|
Q2 , Q2_max
|
Denoise the BOLD data? Denoising is based on modeling and removing nuisance ICs. It may result in a cleaner estimate for smaller datasets, but it may be unnecessary (and time-consuming) for larger datasets. Set If The defaults for both arguments is |
brainstructures |
Only applies if the entries of |
mask |
Required if and only if the entries of |
time_inds |
Subset of fMRI BOLD volumes to include in analysis, as a
vector of numeric integers. If |
varTol |
Tolerance for variance of each data location. For each scan,
locations which do not meet this threshold are masked out of the analysis.
Default: |
spatial_model |
Should spatial modeling be performed? If If If If |
resamp_res |
Only applies if |
rm_mwall |
Only applies if |
reduce_dim |
Reduce the temporal dimension of the data using PCA?
Default: |
method_FC |
Variational Bayes (VB) method for FC template ICA model:
|
maxiter |
Maximum number of EM or VB iterations. Default: |
miniter |
Minimum number of EM or VB iterations. Default: |
epsilon |
Smallest proportion change between iterations. Default:
|
kappa_init |
Starting value for kappa. Default: |
usePar |
Parallelize the computation? Default: |
PW |
Prewhiten to account for residual autocorrelation? Default: |
verbose |
If |
A (spatial) template ICA object, which is a list containing:
subjICmean
, the estimated independent components
S;
subjICse
, the standard errors of S; the
mask
of locations without template values due to too many low
variance or missing values; the nuisance
design matrix or matrices if
applicable; and the function params
such as the type of scaling and
detrending performed.
If BOLD
represented CIFTI or NIFTI data, subjICmean
and
subjICse
will be formatted as "xifti"
or "nifti"
objects, respectively.
## Not run: tm <- estimate_template(cii1_fnames, cii2_fnames, gICA_fname) templateICA(newcii_fname, tm, spatial_model=TRUE, resamp_res=2000) ## End(Not run)
## Not run: tm <- estimate_template(cii1_fnames, cii2_fnames, gICA_fname) templateICA(newcii_fname, tm, spatial_model=TRUE, resamp_res=2000) ## End(Not run)
Parameter Estimates in EM Algorithm for Template ICA Model
UpdateTheta_templateICA.spatial( template_mean, template_var, meshes, BOLD, theta, C_diag, H, Hinv, s0_vec, D, Dinv_s0, verbose = FALSE, return_MAP = FALSE, update = c("all", "kappa", "A") ) UpdateTheta_templateICA.independent( template_mean, template_var, BOLD, theta, C_diag, H, Hinv, update_nu0sq = TRUE, return_MAP = FALSE, verbose = TRUE )
UpdateTheta_templateICA.spatial( template_mean, template_var, meshes, BOLD, theta, C_diag, H, Hinv, s0_vec, D, Dinv_s0, verbose = FALSE, return_MAP = FALSE, update = c("all", "kappa", "A") ) UpdateTheta_templateICA.independent( template_mean, template_var, BOLD, theta, C_diag, H, Hinv, update_nu0sq = TRUE, return_MAP = FALSE, verbose = TRUE )
template_mean |
( |
template_var |
( |
meshes |
|
BOLD |
( |
theta |
(list) current parameter estimates |
C_diag |
|
H , Hinv
|
For dimension reduction |
s0_vec |
Vectorized template means |
D |
Sparse diagonal matrix of template standard deviations |
Dinv_s0 |
The inverse of D times s0_vec |
verbose |
If |
return_MAP |
If |
update |
Which parameters to update. Either |
update_nu0sq |
For non-spatial model: updating |
An updated list of parameter estimates, theta, OR if
return_MAP=TRUE
, the posterior mean and precision of the latent fields
Transform upper-triangular elements to matrix form
UT2mat(x, diag = TRUE)
UT2mat(x, diag = TRUE)
x |
Vector of upper triangular values |
diag |
Are diagonal values included in x? Default: |
Compute the error between empirical and theoretical variance of covariance matrix elements
var_sq_err(nu, p, var_ij, xbar_ij, xbar_ii, xbar_jj)
var_sq_err(nu, p, var_ij, xbar_ij, xbar_ii, xbar_jj)
nu |
Inverse Wishart degrees of freedom parameter |
p |
Matrix dimension for IW distribution |
var_ij |
Empirical between-subject variance of covariance matrices at element (i,j) |
xbar_ij |
Empirical mean of covariance matrices at element (i,j) |
xbar_ii |
Empirical mean of covariance matrices at the ith diagonal element |
xbar_jj |
Empirical mean of covariance matrices at the jth diagonal element |
Squared difference between the empirical and theoretical IW variance of covariance matrices at element (i,j)
Compute the overall error between empirical and theoretical variance of CORRELATION matrix elements
var_sq_err_constrained(nu, p, var, xbar, M = 10000)
var_sq_err_constrained(nu, p, var, xbar, M = 10000)
nu |
Inverse Wishart degrees of freedom parameter |
p |
Matrix dimension for IW distribution |
var |
Empirical between-subject variances of CORRELATION matrix (upper triangle) |
xbar |
Empirical mean of CORRELATION matrix (upper triangle) |
M |
Penalty to assign if theoretical variance is smaller than empirical variance |
Sum of squared difference between the empirical and theoretical IW variance of CORRELATION matrix, but with constraint that theoretical variances must not be smaller than empirical variances