Working with 4D NeuroVectors (NeuroVec)
2025-09-01
Source:vignettes/NeuroVector.Rmd
NeuroVector.Rmd
Working with neuroimaging time-series data
The neuroim2
package contains data structures and
functions for reading, accessing, and processing 4-dimensional
neuroimaging data.
In this vignette, we’ll introduce NeuroVec
(4D images)
and related helpers you’ll use most often:
- Read/write on-disk images (
read_vec
,write_vec
) - Spatial metadata via
NeuroSpace
(dimensions, spacing, origin) - Voxel- and ROI-based access (
series
,series_roi
,coords
) - Common transforms (z-scoring with
scale_series
, concatenation withconcat
) - Dense vs. sparse representations (
DenseNeuroVec
,SparseNeuroVec
,as.sparse
)
Reading a four-dimensional NifTI image with read_vec
Here we read in an example image. This file is 4D in the package data
(64×64×25×4), so read_vec()
returns a NeuroVec
with four volumes (timepoints).
library(purrr)
library(ggplot2)
file_name <- system.file("extdata", "global_mask_v4.nii", package="neuroim2")
vec <- read_vec(file_name)
dim(vec)
#> [1] 64 64 25 4
vec
#>
#> DenseNeuroVec (3.13 bytes MB)
#>
#> - Spatial Info ---------------------------
#> | Dimensions : 64 x 64 x 25 (4 timepoints)
#> | Total Voxels : 102,400
#> | Spacing : 3.5 x 3.5 x 3.7
#>
#> - Properties ---------------------------
#> | Origin : 112 x -108 x -46.2
#> | Orientation : Right-to-Left Posterior-to-Anterior Inferior-to-Superior
#>
#> - Statistics ---------------------------
#> Mean +/- SD : 0.288 +/- 0.453
#>
#> Label: /home/runner/work/_temp/Library/neuroim2/extdata/global_mask_v4.nii
Now imagine we have a set of images. We can read multiple files with
read_vec
. Passing multiple paths returns a
NeuroVecSeq
(a sequence of vectors) rather than a single
concatenated 4D vector.
file_name <- system.file("extdata", "global_mask_v4.nii", package="neuroim2")
vec_multi <- read_vec(c(file_name, file_name, file_name))
dim(vec_multi)
#> [1] 64 64 25 12
vec2 <- read_vec(rep(file_name, 10))
vec2
#>
#> NeuroVecSeq (10 vectors)
#>
#> += Sequence Info ---------------------------
#> | Length : 10
#> | Total Time : 40 points
#>
#> += Spatial Info ---------------------------
#> | Dimensions : 64 x 64 x 25
#> | Spacing : 3.5 x 3.5 x 3.7
#> | Origin : 112 x -108 x -46.2
#> | Orientation : Right-to-Left Posterior-to-Anterior Inferior-to-Superior
#>
#> += Vector Details --------------------------
#> 1. DenseNeuroVec (4 timepoints)
#> 2. DenseNeuroVec (4 timepoints)
#> 3. DenseNeuroVec (4 timepoints)
#> 4. DenseNeuroVec (4 timepoints)
#> 5. DenseNeuroVec (4 timepoints)
#> 6. DenseNeuroVec (4 timepoints)
#> 7. DenseNeuroVec (4 timepoints)
#> 8. DenseNeuroVec (4 timepoints)
#> 9. DenseNeuroVec (4 timepoints)
#> 10. DenseNeuroVec (4 timepoints)
To extract a subset of volumes from a 4D vector, use
sub_vector
directly:
# Extract a subset of volumes (first 3 timepoints)
vec_1_3 <- sub_vector(vec, 1:3)
dim(vec_1_3)
#> [1] 64 64 25 3
vec_1_3
#>
#> DenseNeuroVec (2.35 bytes MB)
#>
#> - Spatial Info ---------------------------
#> | Dimensions : 64 x 64 x 25 (3 timepoints)
#> | Total Voxels : 102,400
#> | Spacing : 3.5 x 3.5 x 3.7
#>
#> - Properties ---------------------------
#> | Origin : 112 x -108 x -46.2
#> | Orientation : Right-to-Left Posterior-to-Anterior Inferior-to-Superior
#>
#> - Statistics ---------------------------
#> Mean +/- SD : 0.288 +/- 0.453
#>
#> Label: none
Extracting time-series data using the series
and
series_roi
functions
To get the time-series at voxel (1,1,1) we can use the
series
function:
series(vec_1_3, 1,1,1)
#> [1] 0 0 0
We can extract a 4d region of interest with the
series_roi
as follows:
file_name <- system.file("extdata", "global_mask_v4.nii", package="neuroim2")
vol <- read_vol(file_name)
roi <- spherical_roi(vol, c(12,12,12), radius=8)
rvec1 <- series_roi(vec, roi)
## or alternatively as a pipeline
rvec2 <- read_vol(file_name) %>% spherical_roi(c(12,12,12), radius=8) %>% series_roi(vec,.)
rvec2
#>
#> === ROIVec Object ===
#>
#> - Structure
#> Points: 49
#> Features: 3 (147 total)
#> Memory: 9.6 Kb
#>
#> - Spatial Properties
#> Parent Space: 64 x 64 x 25 x 4
#> Centroid: [13.0, 13.0, 13.0 mm]
#>
#> - Value Properties
#> Range: [0.00, 0.00]
#>
#> ======================================
#>
#> Access Methods:
#> . Get Points: coords(object)
#> . Get Values: as.matrix(object)
#> . Subset: object[1:10, ]
## we can extract the ROI values with the `values` method.
assertthat::assert_that(all(values(rvec1) == values(rvec2)))
#> [1] TRUE
assertthat::assert_that(all(coords(rvec1) == coords(rvec2)))
#> [1] TRUE
We can also extract an ROI using 1d indices:
r1 <- series_roi(vec, 1:100)
r1
#>
#> === ROIVec Object ===
#>
#> - Structure
#> Points: 100
#> Features: 3 (300 total)
#> Memory: 11.2 Kb
#>
#> - Spatial Properties
#> Parent Space: 64 x 64 x 25 x 4
#> Centroid: [27.5, 1.4, 1.0 mm]
#>
#> - Value Properties
#> Range: [0.00, 0.00]
#>
#> ======================================
#>
#> Access Methods:
#> . Get Points: coords(object)
#> . Get Values: as.matrix(object)
#> . Subset: object[1:10, ]
Or we can extract a plain matrix using the series
function:
We can also use coordinate indexing using voxel coordinates. First we load a binary mask with the same spatial dimensions as our NeuroVec:
mask <- read_vol(system.file("extdata", "global_mask_v4.nii", package="neuroim2"))
Now we convert indices to voxels and extract a matrix of values at the specified locations:
vox <- index_to_grid(mask, 1:100)
r3 <- series(vec, vox)
dim(r3)
#> [1] 4 100
And the same using series_roi
:
r4 <- series_roi(vec,vox)
r4
#>
#> === ROIVec Object ===
#>
#> - Structure
#> Points: 100
#> Features: 3 (300 total)
#> Memory: 12.4 Kb
#>
#> - Spatial Properties
#> Parent Space: 64 x 64 x 25 x 4
#> Centroid: [27.5, 1.4, 1.0 mm]
#>
#> - Value Properties
#> Range: [0.00, 0.00]
#>
#> ======================================
#>
#> Access Methods:
#> . Get Points: coords(object)
#> . Get Values: as.matrix(object)
#> . Subset: object[1:10, ]
Inspecting spatial metadata with NeuroSpace
Every NeuroVec
carries a NeuroSpace
describing its geometry.
sp <- space(vec)
sp # dimensions, spacing, origin, axes
#>
#> NeuroSpace Object
#>
#> >> Dimensions
#> Grid Size: 64 x 64 x 25 x 4
#> Memory: 5.9 KB
#>
#> >> Spatial Properties
#> Spacing: 3.50 x 3.50 x 3.70 mm
#> Origin: 112.00 x -108.00 x -46.20 mm
#>
#> >> Anatomical Orientation
#> X: Right-to-Left | Y: Posterior-to-Anterior | Z: Inferior-to-Superior
#>
#> >> World Transformation
#> Forward (Voxel to World):
#> -3.500 0.000 -0.000 112.000
#> 0.000 3.500 -0.000 -108.000
#> 0.000 0.000 3.700 -46.200
#> 0.000 0.000 0.000 1.000
#> Inverse (World to Voxel):
#> -0.286 -0.000 -0.000 32.000
#> 0.000 0.286 0.000 30.857
#> 0.000 0.000 0.270 12.486
#> 0.000 0.000 0.000 1.000
#>
#> >> Bounding Box
#> Min Corner: -108.5, -108.0, -46.2 mm
#> Max Corner: 112.0, 112.5, 42.6 mm
#>
#> ==================================================
dim(vec) # 4D dims: X × Y × Z × T
#> [1] 64 64 25 4
spacing(vec) # voxel spacing (mm)
#> [1] 3.5 3.5 3.7
origin(vec) # image origin
#> [1] 112.0 -108.0 -46.2
ndim(vec) # == 4 for time series
#> [1] 4
The default mask for dense vectors is “all voxels are valid”:
m <- mask(vec)
m
#>
#> === NeuroVol Object ===
#>
#> * Basic Information
#> Type: LogicalNeuroVol
#> Dimensions: 64 x 64 x 25 (406.6 Kb)
#> Total Voxels: 102,400
#>
#> * Data Properties
#> Value Range: [1.00, 1.00]
#>
#> * Spatial Properties
#> Spacing: 3.50 x 3.50 x 3.70 mm
#> Origin: 112.0, -108.0, -46.2 mm
#> Axes: Left-to-Right x Posterior-to-Anterior x Inferior-to-Superior
#>
#> ======================================
#>
#> Access Methods:
#> . Get Slice: slice(object, zlevel=10)
#> . Get Value: object[i, j, k]
#> . Plot: plot(object) # shows multiple slices
Creating a NeuroVec from an in-memory array
You can build a NeuroVec
directly from arrays or
matrices:
set.seed(1)
dims <- c(24, 24, 24, 5)
arr <- array(rnorm(prod(dims)), dims)
sp4 <- NeuroSpace(dims, spacing = c(2,2,2))
dvec <- NeuroVec(arr, sp4)
dim(dvec)
#> [1] 24 24 24 5
You can also start from a matrix (voxels × time or time × voxels)
using DenseNeuroVec
:
Time-series transforms: z-scoring and summary volumes
Z-score each voxel’s time-series (center and scale across time):
vec_z <- scale_series(dvec, center = TRUE, scale = TRUE)
dim(vec_z)
#> [1] 24 24 24 5
Compute a mean volume across time and return a 3D
NeuroVol
:
M <- as.matrix(dvec) # voxels × time
vmean <- rowMeans(M) # per-voxel mean
mean3d <- NeuroVol(vmean, drop_dim(space(dvec)))
mean3d
#>
#> === NeuroVol Object ===
#>
#> * Basic Information
#> Type: DenseNeuroVol
#> Dimensions: 24 x 24 x 24 (114.6 Kb)
#> Total Voxels: 13,824
#>
#> * Data Properties
#> Value Range: [-1.64, 1.76]
#>
#> * Spatial Properties
#> Spacing: 2.00 x 2.00 x 2.00 mm
#> Origin: 0.0, 0.0, 0.0 mm
#> Axes: Left-to-Right x Posterior-to-Anterior x Inferior-to-Superior
#>
#> ======================================
#>
#> Access Methods:
#> . Get Slice: slice(object, zlevel=10)
#> . Get Value: object[i, j, k]
#> . Plot: plot(object) # shows multiple slices
Dense ↔︎ sparse workflows
Sparse representations store only voxels within a mask. This is handy for large ROIs or brain masks.
# Build a random mask and convert a dense vec to sparse
mask_arr <- array(runif(prod(dims[1:3])) > 0.7, dims[1:3])
mask_vol <- LogicalNeuroVol(mask_arr, drop_dim(sp4))
svec <- as.sparse(dvec, mask_vol) # SparseNeuroVec with explicit mask
svec # note the stored mask and cardinality
#>
#> SparseNeuroVec
#>
#> += Spatial Info ---------------------------
#> | Dimensions : 24 x 24 x 24
#> | Time Points : 5
#> | Spacing : 2 x 2 x 2
#> | Origin : 0 x 0 x 0
#>
#> +- Sparse Info ----------------------------
#> | Cardinality : 4128
#>
#> += Memory Usage --------------------------
#> Size : 306.02 KB
# Convert back to dense if needed
dense_again <- as.dense(svec)
all.equal(dim(dense_again), dim(dvec))
#> [1] TRUE
Tip: For file-backed or memory-mapped vectors, convert to
DenseNeuroVec
via a matrix if you need dense-only
operations:
dv_dense <- DenseNeuroVec(as.matrix(vec), space(vec))
Writing vectors to disk
You can write NeuroVec
and NeuroVol
objects
as NIfTI files:
tmp_vec <- tempfile(fileext = ".nii.gz")
write_vec(vec_1_3, tmp_vec)
file.exists(tmp_vec)
#> [1] TRUE
unlink(tmp_vec)
Putting it together with an ROI
Combine ROI extraction with time-series transforms:
vol3d <- read_vol(system.file("extdata", "global_mask_v4.nii", package="neuroim2"))
roi <- spherical_roi(vol3d, c(12,12,12), radius = 6)
rts <- series_roi(vec, roi) # ROIVec (T × N with coords)
# z-score each column (voxel) across time, then average within ROI
mat_roi <- values(rts) # T × N
mat_z <- base::scale(mat_roi, center=TRUE, scale=TRUE)
roi_mean <- rowMeans(mat_z)
length(roi_mean) # matches time dimension
#> [1] 4
That’s the core workflow for 4D data in neuroim2
: load
or create a NeuroVec
, inspect its NeuroSpace
,
access time-series via voxels or ROIs, apply simple transforms, and
optionally move between dense and sparse representations.