| Title: | Spatial Null Models and Transforms for Brain Map Comparison |
|---|---|
| Description: | Implements spatial null models and coordinate-space transformations for statistical comparison of brain maps, following the framework described in Markello et al. (2022) <doi:10.1038/s41592-022-01625-w>. Provides variogram-matching surrogates (Burt et al. 2020), Moran spectral randomization (Wagner & Dray 2015), and spin-based permutation tests (Alexander-Bloch et al. 2018). Includes an R interface to the 'neuromaps' annotation registry for browsing, downloading, and comparing brain map annotations from the Open Science Framework ('OSF'). Integrates with 'ciftiTools' for coordinate-space transforms. |
| Authors: | Athanasia Mo Mowinckel [aut, cre, cph] (ORCID: <https://orcid.org/0000-0002-5756-0223>) |
| Maintainer: | Athanasia Mo Mowinckel <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 0.2.2 |
| Built: | 2026-05-30 08:21:04 UTC |
| Source: | https://github.com/LCBC-UiO/neuromapr |
Reads a FreeSurfer .annot file and writes a GIFTI label file.
annot_to_gifti(annot_path, output_path = NULL)annot_to_gifti(annot_path, output_path = NULL)
annot_path |
Path to FreeSurfer |
output_path |
Output GIFTI path. If |
The output file path (invisibly).
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
## Not run: annot_to_gifti("lh.aparc.annot") annot_to_gifti("lh.aparc.annot", "lh.aparc.label.gii") ## End(Not run)## Not run: annot_to_gifti("lh.aparc.annot") annot_to_gifti("lh.aparc.annot", "lh.aparc.label.gii") ## End(Not run)
Verifies that wb_command is available. If wb_path is NULL,
checks ciftiTools default and system PATH.
check_wb_command(wb_path = NULL)check_wb_command(wb_path = NULL)
wb_path |
Optional explicit path to |
Path to wb_command executable.
## Not run: check_wb_command() ## End(Not run)## Not run: check_wb_command() ## End(Not run)
Removes the session-level cache of the neuromaps annotation registry,
forcing a fresh download on the next call to neuromaps_available() or
fetch_neuromaps_annotation().
clear_neuromaps_cache()clear_neuromaps_cache()
NULL, invisibly.
clear_neuromaps_cache()clear_neuromaps_cache()
Computes the correlation between two brain maps and optionally tests significance using a spatial null model to account for spatial autocorrelation.
compare_maps( x, y, method = c("pearson", "spearman"), null_method = NULL, n_perm = 1000L, nulls = NULL, distmat = NULL, coords = NULL, seed = NULL, na.rm = TRUE, verbose = TRUE, ... ) ## S3 method for class 'neuromaps_enhanced_comparison' print(x, ...)compare_maps( x, y, method = c("pearson", "spearman"), null_method = NULL, n_perm = 1000L, nulls = NULL, distmat = NULL, coords = NULL, seed = NULL, na.rm = TRUE, verbose = TRUE, ... ) ## S3 method for class 'neuromaps_enhanced_comparison' print(x, ...)
x, y
|
Numeric vectors of brain map values, or file paths to GIFTI/NIfTI files. |
method |
Correlation method: |
null_method |
Optional null model method for empirical p-values.
One of |
n_perm |
Integer number of null permutations. |
nulls |
Pre-computed null_distribution object for |
distmat |
Distance matrix (passed to null model if needed). |
coords |
Coordinate list (passed to spin null models if needed). |
seed |
Optional integer seed for reproducibility. |
na.rm |
Logical, remove NA values before computing correlation. |
verbose |
Logical, print progress messages. |
... |
Additional arguments passed to |
A neuromaps_enhanced_comparison object (inherits
neuromaps_comparison) with additional fields p_null, null_method,
null_r, and n_perm.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
x <- rnorm(50) y <- x + rnorm(50) compare_maps(x, y, verbose = FALSE)x <- rnorm(50) y <- x + rnorm(50) compare_maps(x, y, verbose = FALSE)
Download brain map annotation files from the neuromaps OSF repository. Surface annotations return two files (left and right hemispheres), volume annotations return one.
fetch_neuromaps_annotation( source, desc, space, density = NULL, resolution = NULL, hemisphere = NULL, data_dir = neuromaps_cache_dir(), overwrite = FALSE, verbose = TRUE )fetch_neuromaps_annotation( source, desc, space, density = NULL, resolution = NULL, hemisphere = NULL, data_dir = neuromaps_cache_dir(), overwrite = FALSE, verbose = TRUE )
source |
Data source identifier (e.g. |
desc |
Map descriptor key (e.g. |
space |
Coordinate space (e.g. |
density |
Surface vertex density (e.g. |
resolution |
Volume voxel resolution (e.g. |
hemisphere |
Hemisphere ( |
data_dir |
Directory for cached downloads. Defaults to
|
overwrite |
Re-download even if cached file exists. |
verbose |
Print progress messages. |
Character vector of downloaded file path(s).
fetch_neuromaps_annotation("abagen", "genepc1", "fsaverage", density = "10k")fetch_neuromaps_annotation("abagen", "genepc1", "fsaverage", density = "10k")
Reads a FreeSurfer morphometry file (.curv, .thickness, .sulc) and
writes a GIFTI func file.
fsmorph_to_gifti(morph_path, output_path = NULL)fsmorph_to_gifti(morph_path, output_path = NULL)
morph_path |
Path to FreeSurfer morphometry file. |
output_path |
Output GIFTI path. If |
The output file path (invisibly).
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
## Not run: fsmorph_to_gifti("lh.thickness") fsmorph_to_gifti("lh.thickness", "lh.thickness.func.gii") ## End(Not run)## Not run: fsmorph_to_gifti("lh.thickness") fsmorph_to_gifti("lh.thickness", "lh.thickness.func.gii") ## End(Not run)
Dispatches to the appropriate null model method for generating spatially-constrained surrogate brain maps.
generate_nulls( data, method = c("burt2020", "moran", "spin_vasa", "spin_hungarian", "alexander_bloch", "baum", "cornblath", "burt2018"), n_perm = 1000L, distmat = NULL, coords = NULL, parcellation = NULL, seed = NULL, ... )generate_nulls( data, method = c("burt2020", "moran", "spin_vasa", "spin_hungarian", "alexander_bloch", "baum", "cornblath", "burt2018"), n_perm = 1000L, distmat = NULL, coords = NULL, parcellation = NULL, seed = NULL, ... )
data |
Numeric vector of brain map values. |
method |
Character string specifying the null model method. |
n_perm |
Integer number of null permutations to generate. |
distmat |
Distance matrix (required for |
coords |
List with |
parcellation |
Integer vector of parcel labels (required for |
seed |
Optional integer seed for reproducibility. |
... |
Additional arguments passed to the specific null method
(e.g. |
A null_distribution object.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
data <- rnorm(100) distmat <- as.matrix(dist(seq_len(100))) nd <- generate_nulls(data, method = "moran", distmat = distmat, n_perm = 10L)data <- rnorm(100) distmat <- as.matrix(dist(seq_len(100))) nd <- generate_nulls(data, method = "moran", distmat = distmat, n_perm = 10L)
Finds the centroid of each parcel using one of three methods.
get_parcel_centroids( vertices, labels, method = c("average", "surface", "geodesic"), faces = NULL )get_parcel_centroids( vertices, labels, method = c("average", "surface", "geodesic"), faces = NULL )
vertices |
Numeric matrix (n x 3) of vertex coordinates. |
labels |
Integer vector of parcel labels. |
method |
Centroid method: |
faces |
Integer matrix (m x 3) of face indices. Required for
|
Numeric matrix (n_parcels x 3) with rownames set to parcel labels.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
vertices <- matrix(rnorm(30), ncol = 3) labels <- c(1L, 1L, 1L, 2L, 2L, 2L, 3L, 3L, 3L, 3L) get_parcel_centroids(vertices, labels, method = "average")vertices <- matrix(rnorm(30), ncol = 3) labels <- c(1L, 1L, 1L, 2L, 2L, 2L, 3L, 3L, 3L, 3L) get_parcel_centroids(vertices, labels, method = "average")
Builds a graph from a triangular mesh and computes shortest-path (Dijkstra) distances between vertices.
get_surface_distance(vertices, faces, source_vertices = NULL)get_surface_distance(vertices, faces, source_vertices = NULL)
vertices |
Numeric matrix (n x 3) of vertex coordinates. |
faces |
Integer matrix (m x 3) of face indices (1-indexed). |
source_vertices |
Optional integer vector of source vertex indices.
If |
Numeric distance matrix. If source_vertices is provided, returns
a length(source_vertices) x n matrix; otherwise an n x n matrix.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
vertices <- matrix( c(0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1), nrow = 4, byrow = TRUE ) faces <- matrix(c(1L, 2L, 3L, 2L, 3L, 4L), nrow = 2, byrow = TRUE) get_surface_distance(vertices, faces)vertices <- matrix( c(0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1), nrow = 4, byrow = TRUE ) faces <- matrix(c(1L, 2L, 3L, 2L, 3L, 4L), nrow = 2, byrow = TRUE) get_surface_distance(vertices, faces)
Extracts unique edges from triangular faces, computes Euclidean edge weights, and returns an igraph graph object suitable for geodesic distance computation.
make_surf_graph(vertices, faces)make_surf_graph(vertices, faces)
vertices |
Numeric matrix (n x 3) of vertex coordinates. |
faces |
Integer matrix (m x 3) of face indices (1-indexed). |
An igraph graph object with weighted edges.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
vertices <- matrix( c(0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1), nrow = 4, byrow = TRUE ) faces <- matrix(c(1L, 2L, 3L, 2L, 3L, 4L), nrow = 2, byrow = TRUE) g <- make_surf_graph(vertices, faces)vertices <- matrix( c(0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1), nrow = 4, byrow = TRUE ) faces <- matrix(c(1L, 2L, 3L, 2L, 3L, 4L), nrow = 2, byrow = TRUE) g <- make_surf_graph(vertices, faces)
Query the neuromaps registry to see which brain map annotations are available for download. Data is fetched from the neuromaps project's GitHub repository and cached for the session.
neuromaps_available( source = NULL, desc = NULL, space = NULL, density = NULL, resolution = NULL, hemisphere = NULL, tags = NULL, format = NULL, refresh = FALSE, fixed = FALSE )neuromaps_available( source = NULL, desc = NULL, space = NULL, density = NULL, resolution = NULL, hemisphere = NULL, tags = NULL, format = NULL, refresh = FALSE, fixed = FALSE )
source |
Data source identifier (e.g. |
desc |
Map descriptor key (e.g. |
space |
Coordinate space (e.g. |
density |
Surface vertex density (e.g. |
resolution |
Volume voxel resolution (e.g. |
hemisphere |
Hemisphere ( |
tags |
Character vector of tags. All must match (AND logic). |
format |
Filter by format ( |
refresh |
Logical. If |
fixed |
Logical. If |
All string filter parameters (source, desc, space, density,
resolution, hemisphere, format) are treated as R regular
expressions and matched with grepl(). For example,
source = "^beliveau$" matches exactly, while source = "bel" matches
any source containing "bel". Set fixed = TRUE for literal string
matching. The tags parameter always uses exact matching (AND logic).
When used in fetch_neuromaps_annotation(), source, desc, and space
are exact matches.
A tibble of available annotations with columns: source, desc, space, den, res, hemi, format, fname, full_desc, tags, N, age.
neuromaps_available() neuromaps_available(source = "beliveau") neuromaps_available(tags = "pet")neuromaps_available() neuromaps_available(source = "beliveau") neuromaps_available(tags = "pet")
Original vertex-level spin test: rotates coordinates and assigns each rotated vertex the value of its nearest original vertex (no optimal matching).
null_alexander_bloch( data, coords, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )null_alexander_bloch( data, coords, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )
data |
Numeric vector of brain map values. |
coords |
List with |
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
rotation |
Rotation generation method: |
A null_distribution object.
Alexander-Bloch AF et al. (2018) NeuroImage 175:111-120. doi:10.1016/j.neuroimage.2018.04.023
coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) data <- rnorm(20) nd <- null_alexander_bloch(data, coords, n_perm = 10L, seed = 1L)coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) data <- rnorm(20) nd <- null_alexander_bloch(data, coords, n_perm = 10L, seed = 1L)
Spin-based null model with maximum-overlap parcel reassignment. After rotating vertex coordinates, each original parcel is assigned the value of the rotated parcel with the most vertex overlap.
null_baum( data, coords, parcellation, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )null_baum( data, coords, parcellation, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )
data |
Numeric vector of brain map values. |
coords |
List with |
parcellation |
Integer vector of parcel labels for all vertices.
|
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
rotation |
Rotation generation method: |
A null_distribution object.
Baum GL et al. (2020) PNAS 117:21854-21861. doi:10.1073/pnas.2005518117
coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) parcellation <- c(rep(1L, 5), rep(2L, 5), rep(3L, 5), rep(4L, 5)) data <- c(1.0, 2.0, 3.0, 4.0) nd <- null_baum(data, coords, parcellation, n_perm = 10L, seed = 1L)coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) parcellation <- c(rep(1L, 5), rep(2L, 5), rep(3L, 5), rep(4L, 5)) data <- c(1.0, 2.0, 3.0, 4.0) nd <- null_baum(data, coords, parcellation, n_perm = 10L, seed = 1L)
Generates surrogate brain maps using a spatial autoregressive (SAR) model. Estimates spatial autocorrelation and distance decay parameters, then generates surrogates by solving the SAR equation and rank-matching to the original data.
null_burt2018(data, distmat, n_perm = 1000L, seed = NULL)null_burt2018(data, distmat, n_perm = 1000L, seed = NULL)
data |
Numeric vector of brain map values. |
distmat |
Distance matrix between parcels/vertices. |
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
A null_distribution object.
Burt JB et al. (2018) Nature Neuroscience 21:1251-1259. doi:10.1038/s41593-018-0195-0
data <- rnorm(50) distmat <- as.matrix(dist(matrix(rnorm(100), 50, 2))) nd <- null_burt2018(data, distmat, n_perm = 10L, seed = 1L)data <- rnorm(50) distmat <- as.matrix(dist(matrix(rnorm(100), 50, 2))) nd <- null_burt2018(data, distmat, n_perm = 10L, seed = 1L)
Generates spatially-constrained surrogate brain maps by matching the empirical variogram of the original data through smoothed random permutations.
null_burt2020( data, distmat, n_perm = 1000L, seed = NULL, ns = 500L, nh = 25L, pv = 25, knn = 1000L, deltas = seq(0.1, 0.9, by = 0.1), kernel = c("exponential", "gaussian", "uniform"), resample = FALSE )null_burt2020( data, distmat, n_perm = 1000L, seed = NULL, ns = 500L, nh = 25L, pv = 25, knn = 1000L, deltas = seq(0.1, 0.9, by = 0.1), kernel = c("exponential", "gaussian", "uniform"), resample = FALSE )
data |
Numeric vector of brain map values. |
distmat |
Distance matrix between parcels/vertices. |
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
ns |
Integer, subsample size for variogram computation. |
nh |
Integer, number of distance bins for variogram. |
pv |
Numeric, percentile cutoff for maximum distance in variogram. |
knn |
Integer, number of nearest neighbors for smoothing. |
deltas |
Numeric vector of smoothing levels (fractions of |
kernel |
Smoothing kernel function. |
resample |
Logical. If |
A null_distribution object.
Burt JB et al. (2020) NeuroImage 220:117038. doi:10.1016/j.neuroimage.2020.117038
data <- rnorm(50) distmat <- as.matrix(dist(matrix(rnorm(100), 50, 2))) nd <- null_burt2020(data, distmat, n_perm = 10L, seed = 1L)data <- rnorm(50) distmat <- as.matrix(dist(matrix(rnorm(100), 50, 2))) nd <- null_burt2020(data, distmat, n_perm = 10L, seed = 1L)
Spin-based null model where each rotated vertex receives the label of its nearest non-medial-wall original vertex, then parcels are reassigned by majority vote among the resulting vertex labels.
null_cornblath( data, coords, parcellation, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )null_cornblath( data, coords, parcellation, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )
data |
Numeric vector of brain map values. |
coords |
List with |
parcellation |
Integer vector of parcel labels for all vertices.
|
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
rotation |
Rotation generation method: |
A null_distribution object.
Cornblath EJ et al. (2020) Communications Biology 3:590. doi:10.1038/s42003-020-01296-5
coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) parcellation <- c(rep(1L, 5), rep(2L, 5), rep(3L, 5), rep(4L, 5)) data <- c(1.0, 2.0, 3.0, 4.0) nd <- null_cornblath(data, coords, parcellation, n_perm = 10L, seed = 1L)coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) parcellation <- c(rep(1L, 5), rep(2L, 5), rep(3L, 5), rep(4L, 5)) data <- c(1.0, 2.0, 3.0, 4.0) nd <- null_cornblath(data, coords, parcellation, n_perm = 10L, seed = 1L)
Create a null distribution object
new_null_distribution(nulls, method, observed, params = list()) ## S3 method for class 'null_distribution' print(x, ...) ## S3 method for class 'null_distribution' summary(object, ...) ## S3 method for class 'null_distribution' as.matrix(x, ...) ## S3 method for class 'null_distribution' plot(x, parcel = 1L, ...)new_null_distribution(nulls, method, observed, params = list()) ## S3 method for class 'null_distribution' print(x, ...) ## S3 method for class 'null_distribution' summary(object, ...) ## S3 method for class 'null_distribution' as.matrix(x, ...) ## S3 method for class 'null_distribution' plot(x, parcel = 1L, ...)
nulls |
Numeric matrix (n x n_perm) of surrogate values. |
method |
Character string identifying the null model method. |
observed |
Numeric vector of original data values. |
params |
Named list of algorithm parameters. |
x |
A |
... |
Ignored. |
object |
A |
parcel |
Integer index of the parcel to plot. |
A null_distribution object.
nulls <- matrix(rnorm(30), nrow = 3, ncol = 10) nd <- new_null_distribution(nulls, "test", observed = c(1, 2, 3)) print(nd) summary(nd)nulls <- matrix(rnorm(30), nrow = 3, ncol = 10) nd <- new_null_distribution(nulls, "test", observed = c(1, 2, 3)) print(nd) summary(nd)
Generates spatially-constrained surrogate brain maps using Moran's eigenvector maps (MEMs) for spectral randomization.
null_moran( data, distmat, n_perm = 1000L, seed = NULL, procedure = c("pair", "singleton"), kernel = c("inverse_distance", "exponential", "gaussian", "bisquare"), tol = 1e-06 )null_moran( data, distmat, n_perm = 1000L, seed = NULL, procedure = c("pair", "singleton"), kernel = c("inverse_distance", "exponential", "gaussian", "bisquare"), tol = 1e-06 )
data |
Numeric vector of brain map values. |
distmat |
Distance matrix between parcels/vertices. |
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
procedure |
Character, either |
kernel |
Weight matrix kernel: |
tol |
Numeric tolerance for eigenvalue comparison. |
A null_distribution object.
Wagner HH, Dray S (2015) Methods in Ecology and Evolution 6:1169-1178. doi:10.1111/2041-210X.12407
data <- rnorm(50) distmat <- as.matrix(dist(matrix(rnorm(100), 50, 2))) nd <- null_moran(data, distmat, n_perm = 10L, seed = 1L)data <- rnorm(50) distmat <- as.matrix(dist(matrix(rnorm(100), 50, 2))) nd <- null_moran(data, distmat, n_perm = 10L, seed = 1L)
Generate spatially-constrained null distributions using spin-based permutation of spherical coordinates.
null_spin_vasa( data, coords, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") ) null_spin_hungarian( data, coords, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )null_spin_vasa( data, coords, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") ) null_spin_hungarian( data, coords, n_perm = 1000L, seed = NULL, rotation = c("euler", "rodrigues") )
data |
Numeric vector of brain map values. |
coords |
List with |
n_perm |
Integer number of null permutations to generate. |
seed |
Optional integer seed for reproducibility. |
rotation |
Rotation generation method: |
A null_distribution object.
Alexander-Bloch AF et al. (2018) NeuroImage 175:111-120. doi:10.1016/j.neuroimage.2018.04.023
Vasa F et al. (2018) Cerebral Cortex 28:3293-3303. doi:10.1093/cercor/bhx195
Markello RD, Misic B (2021) NeuroImage 236:118052. doi:10.1016/j.neuroimage.2021.118052
coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) data <- rnorm(20) nd <- null_spin_vasa(data, coords, n_perm = 10L, seed = 1L)coords <- list(lh = matrix(rnorm(30), 10, 3), rh = matrix(rnorm(30), 10, 3)) data <- rnorm(20) nd <- null_spin_vasa(data, coords, n_perm = 10L, seed = 1L)
High-level wrapper that reads data and parcellation from file paths or vectors, then aggregates vertices into parcels.
parcellate(data, parcellation, summary_func = mean)parcellate(data, parcellation, summary_func = mean)
data |
Numeric vector or file path to a GIFTI/NIfTI brain map. |
parcellation |
Integer vector of labels or file path to a GIFTI label file. |
summary_func |
Function to summarise each parcel (default: mean). |
Named numeric vector of parcel-level values.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
data <- c(1.0, 2.0, 3.0, 4.0) labels <- c(1L, 1L, 2L, 2L) parcellate(data, labels)data <- c(1.0, 2.0, 3.0, 4.0) labels <- c(1L, 1L, 2L, 2L) parcellate(data, labels)
Expands parcel-level values to a vertex-level vector using parcellation labels.
parcels_to_vertices(parcel_data, labels, fill = NA_real_)parcels_to_vertices(parcel_data, labels, fill = NA_real_)
parcel_data |
Named numeric vector of parcel values (names match labels). |
labels |
Integer vector of parcel labels. |
fill |
Value for medial wall vertices (default: |
Numeric vector of length(labels).
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
parcel_data <- c("1" = 10, "2" = 20) labels <- c(1L, 1L, 2L, 2L, 0L) parcels_to_vertices(parcel_data, labels)parcel_data <- c("1" = 10, "2" = 20) labels <- c(1L, 1L, 2L, 2L, 0L) parcels_to_vertices(parcel_data, labels)
Computes a user-specified metric between two vectors and tests significance using either spatially-constrained null surrogates or simple random permutation.
permtest_metric( x, y, metric_func = stats::cor, n_perm = 1000L, seed = NULL, null_method = NULL, distmat = NULL, coords = NULL, parcellation = NULL, ... )permtest_metric( x, y, metric_func = stats::cor, n_perm = 1000L, seed = NULL, null_method = NULL, distmat = NULL, coords = NULL, parcellation = NULL, ... )
x, y
|
Numeric vectors. |
metric_func |
Function taking |
n_perm |
Integer number of permutations. |
seed |
Optional integer seed for reproducibility. |
null_method |
Optional null model method passed to |
distmat |
Distance matrix (passed to |
coords |
Coordinate list (passed to |
parcellation |
Integer vector (passed to |
... |
Additional arguments passed to |
List with $observed, $null_values, $p_value, and $n_perm.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
x <- rnorm(100) y <- x + rnorm(100) result <- permtest_metric(x, y, n_perm = 99L, seed = 1L) result$observed result$p_valuex <- rnorm(100) y <- x + rnorm(100) result <- permtest_metric(x, y, n_perm = 99L, seed = 1L) result$observed result$p_value
Reads GIFTI (.func.gii) or NIfTI (.nii.gz) files and returns the
values as a numeric vector. Used internally by compare_maps() when
file paths are passed instead of numeric vectors.
read_brain_map_values(path)read_brain_map_values(path)
path |
Path to a |
A numeric vector of map values.
## Not run: read_brain_map_values("cortical_thickness.func.gii") read_brain_map_values("brain_volume.nii.gz") ## End(Not run)## Not run: read_brain_map_values("cortical_thickness.func.gii") read_brain_map_values("brain_volume.nii.gz") ## End(Not run)
Aligns two brain maps into the same coordinate space and density before comparison. Supports multiple strategies for choosing the target space.
resample_images( src, trg, src_space = c("fsaverage", "fsLR"), trg_space = c("fsaverage", "fsLR"), strategy = c("downsample_only", "transform_to_src", "transform_to_trg", "transform_to_alt"), alt_space = NULL, alt_density = NULL, hemisphere = c("left", "right"), area_surf_current = NULL, area_surf_new = NULL, wb_path = NULL, verbose = TRUE )resample_images( src, trg, src_space = c("fsaverage", "fsLR"), trg_space = c("fsaverage", "fsLR"), strategy = c("downsample_only", "transform_to_src", "transform_to_trg", "transform_to_alt"), alt_space = NULL, alt_density = NULL, hemisphere = c("left", "right"), area_surf_current = NULL, area_surf_new = NULL, wb_path = NULL, verbose = TRUE )
src |
Character, file path to the source GIFTI. |
trg |
Character, file path to the target GIFTI. |
src_space |
Source coordinate space ( |
trg_space |
Target coordinate space ( |
strategy |
Resampling strategy. One of |
alt_space |
Alternative space for |
alt_density |
Alternative density for |
hemisphere |
Which hemispheres: |
area_surf_current |
Path to a current-resolution area-correction
surface (e.g. midthickness). Passed to |
area_surf_new |
Path to a target-resolution area-correction surface.
Passed to |
wb_path |
Path to |
verbose |
Logical, print progress messages. |
List with $src and $trg file paths in the aligned space.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
## Not run: resample_images("src.func.gii", "trg.func.gii", src_space = "fsaverage", trg_space = "fsaverage" ) ## End(Not run)## Not run: resample_images("src.func.gii", "trg.func.gii", src_space = "fsaverage", trg_space = "fsaverage" ) ## End(Not run)
Resamples GIFTI surface files between fsaverage and fsLR coordinate
spaces using Connectome Workbench via ciftiTools.
transform_to_space( paths, target_space = c("fsLR", "fsaverage"), target_density = "32k", hemisphere = c("left", "right"), method = c("adaptive", "barycentric"), area_surf_current = NULL, area_surf_new = NULL, wb_path = NULL, verbose = TRUE )transform_to_space( paths, target_space = c("fsLR", "fsaverage"), target_density = "32k", hemisphere = c("left", "right"), method = c("adaptive", "barycentric"), area_surf_current = NULL, area_surf_new = NULL, wb_path = NULL, verbose = TRUE )
paths |
Character vector of GIFTI file paths to transform. |
target_space |
Target coordinate space: |
target_density |
Target mesh density (e.g., |
hemisphere |
Which hemispheres to transform: |
method |
Resampling method: |
area_surf_current |
Path to the current-resolution area-correction
surface (e.g. midthickness). Recommended when |
area_surf_new |
Path to the target-resolution area-correction surface.
If |
wb_path |
Path to |
verbose |
Logical, print progress messages. |
Character vector of transformed file paths.
Robinson EC et al. (2014) NeuroImage 100:414-426. doi:10.1016/j.neuroimage.2014.07.033
Robinson EC et al. (2018) NeuroImage 167:150-165. doi:10.1016/j.neuroimage.2017.10.037
## Not run: transform_to_space("map.func.gii", target_space = "fsLR") ## End(Not run)## Not run: transform_to_space("map.func.gii", target_space = "fsLR") ## End(Not run)
Inverse of parcellate(): maps parcel-level values back to vertices.
unparcellate(parcel_data, parcellation, fill = NA_real_)unparcellate(parcel_data, parcellation, fill = NA_real_)
parcel_data |
Named numeric vector of parcel values. |
parcellation |
Integer vector of labels or file path to a GIFTI label file. |
fill |
Value for medial wall vertices (default: |
Numeric vector of vertex-level values.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
parcel_data <- c("1" = 10, "2" = 20) labels <- c(1L, 1L, 2L, 2L, 0L) unparcellate(parcel_data, labels)parcel_data <- c("1" = 10, "2" = 20) labels <- c(1L, 1L, 2L, 2L, 0L) unparcellate(parcel_data, labels)
Each triangle distributes one-third of its area to each of its three vertices. Triangle area is computed via the cross-product formula.
vertex_areas(vertices, faces)vertex_areas(vertices, faces)
vertices |
Numeric matrix (n x 3) of vertex coordinates. |
faces |
Integer matrix (m x 3) of face indices (1-indexed). |
Numeric vector of length nrow(vertices).
vertices <- matrix(c(0, 1, 0, 0, 0, 1, 0, 0, 0), nrow = 3, byrow = TRUE) faces <- matrix(c(1L, 2L, 3L), nrow = 1) vertex_areas(vertices, faces)vertices <- matrix(c(0, 1, 0, 0, 0, 1, 0, 0, 0), nrow = 3, byrow = TRUE) faces <- matrix(c(1L, 2L, 3L), nrow = 1) vertex_areas(vertices, faces)
Summarises vertex-level data by parcellation labels.
vertices_to_parcels(data, labels, summary_func = mean)vertices_to_parcels(data, labels, summary_func = mean)
data |
Numeric vector of vertex-level values. |
labels |
Integer vector of parcel labels. |
summary_func |
Function to summarise each parcel (default: mean). |
Named numeric vector of parcel-level values.
Markello RD et al. (2022) Nature Methods 19:1472-1480. doi:10.1038/s41592-022-01625-w
data <- c(1.0, 2.0, 3.0, 4.0) labels <- c(1L, 1L, 2L, 2L) vertices_to_parcels(data, labels)data <- c(1.0, 2.0, 3.0, 4.0) labels <- c(1L, 1L, 2L, 2L) vertices_to_parcels(data, labels)