| Title: | S7 Data Structures for Diffusion MRI Tractography |
|---|---|
| Description: | Provides three S7 classes — streamline, bundle, and bundle_set — for representing diffusion MRI tractography data in R, together with a concise set of methods for computing shape descriptors (arc-length, curvature, torsion, sinuosity), the Hausdorff distance between streamlines, arc-length reparametrization of streamlines and bundles onto uniform grids, combination of streamlines or bundles into a single bundle, combination of bundles from multiple subjects or sessions into a bundle_set, and coercion to and from the dwiFiber S4 class of the 'dti' package. See Dell'Acqua, F., Descoteaux, M. and Leemans, A. (2024) "Handbook of Diffusion MR Tractography" <doi:10.1016/C2018-0-02520-7> for more about the mathematical and computational underpinnings of diffusion MRI tractography. |
| Authors: | Aymeric Stamm [aut, cre] (ORCID: <https://orcid.org/0000-0002-8725-3654>) |
| Maintainer: | Aymeric Stamm <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 0.1.2 |
| Built: | 2026-05-31 09:59:01 UTC |
| Source: | https://github.com/tractoverse/fiber |
add_shape_descriptors() is an S7 generic that computes a number of shape
descriptors for each streamline object and stores them in the
@streamline_data or @point_data slots as appropriate, with methods
available for the following classes:
This function provides a convenient way to compute shape descriptors and
attach them to streamline or bundle objects. See the documentation for
each individual shape descriptor function (e.g. get_euclidean_length(),
get_curvilinear_length(), get_sinuosity(), get_curvature(),
get_torsion()) for more details on how each descriptor is computed.
add_shape_descriptors( x, descriptors = c("euclidean_length", "curvilinear_length", "sinuosity", "curvature", "torsion") )add_shape_descriptors( x, descriptors = c("euclidean_length", "curvilinear_length", "sinuosity", "curvature", "torsion") )
x |
A streamline or bundle object. |
descriptors |
A character vector of shape descriptors to add. Defaults
to all available descriptors: |
An object of the same class as x with the specified shape descriptors
added to the @streamline_data or @point_data slots of each streamline.
# add multiple shape descriptors to a single streamline pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) sl <- add_shape_descriptors( sl, descriptors = c("euclidean_length", "curvilinear_length", "sinuosity") ) # add multiple shape descriptors to a bundle sl1 <- streamline(points = pts) pts2 <- matrix(runif(60), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl2 <- streamline(points = pts2) b <- bundle(streamlines = list(sl1, sl2)) b <- add_shape_descriptors( b, descriptors = c("euclidean_length", "curvilinear_length", "sinuosity") )# add multiple shape descriptors to a single streamline pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) sl <- add_shape_descriptors( sl, descriptors = c("euclidean_length", "curvilinear_length", "sinuosity") ) # add multiple shape descriptors to a bundle sl1 <- streamline(points = pts) pts2 <- matrix(runif(60), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl2 <- streamline(points = pts2) b <- bundle(streamlines = list(sl1, sl2)) b <- add_shape_descriptors( b, descriptors = c("euclidean_length", "curvilinear_length", "sinuosity") )
add_shape_descriptors() method for bundle objectsAdds multiple shape descriptors to every streamline inside a bundle.
x |
A bundle object. |
descriptors |
A character vector of shape descriptors to add. Defaults
to all available descriptors: |
A bundle with the specified shape descriptors added to the
@streamline_data or @point_data slots of each streamline as appropriate.
add_shape_descriptors() method for streamline objectsAdds multiple shape descriptors to a single streamline object.
x |
A streamline object. |
descriptors |
A character vector of shape descriptors to add. Defaults
to all available descriptors: |
A streamline with the specified shape descriptors added to the
@streamline_data or @point_data slots as appropriate.
as_bundle() converts a supported object into a bundle.
as_bundle(x, ...)as_bundle(x, ...)
x |
An object to coerce. |
... |
Additional arguments (currently unused). |
Currently supported input classes:
streamline: wrapped in a single-element bundle (lossless).
bundle: returned unchanged.
dwiFiber (from dti): each fiber becomes a streamline. The
per-point direction vectors (columns 4–6 of @fibers) are stored as
@point_data$direction_x, @point_data$direction_y, and
@point_data$direction_z. Tracking metadata (method, minfa,
maxangle) are stored in @bundle_data.
A bundle object.
as_streamline(), as_dwifiber()
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- as_bundle(sl) b@n_streamlines # 1pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- as_bundle(sl) b@n_streamlines # 1
as_bundle_set() converts a supported object into a bundle_set.
as_bundle_set(x, ...)as_bundle_set(x, ...)
x |
An object to coerce. |
... |
Additional arguments passed to methods (e.g. |
Currently supported input classes:
bundle_set: returned unchanged.
bundle: wrapped in a single-element bundle_set. An optional name
argument sets the element name (defaults to "bundle_1").
A bundle_set object.
bundle_set(), bind_bundle_sets()
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b <- bundle(streamlines = list(streamline(points = pts))) bs <- as_bundle_set(b, name = "sub-01") bs@n_bundles # 1pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b <- bundle(streamlines = list(streamline(points = pts))) bs <- as_bundle_set(b, name = "sub-01") bs@n_bundles # 1
as_bundle() method for bundle objectsas_bundle() method for bundle objects
x |
A bundle object. |
... |
Additional arguments (currently unused). |
x unchanged.
as_bundle() method for streamline objectsas_bundle() method for streamline objects
x |
A streamline object. |
... |
Additional arguments (currently unused). |
A bundle containing x as its sole streamline.
dwiFiber objectas_dwifiber() converts a streamline or bundle to the S4 class
dwiFiber from the dti package.
as_dwifiber(x, ...)as_dwifiber(x, ...)
x |
A streamline or bundle object. |
... |
Additional arguments (currently unused). |
Per-point direction vectors are taken from @point_data$direction_x,
@point_data$direction_y, and @point_data$direction_z when
present; otherwise they are estimated via finite differences of the
coordinates (forward difference at the first point, backward difference at
the last, central differences in between), then unit-normalised.
Bundle-level metadata stored in @bundle_data under the keys method,
minfa, maxangle, ddim, ddim0, voxelext, orientation, rotation,
level, and source are transferred to the corresponding dwiFiber / dwi
slots when present. MRI-acquisition metadata that cannot be recovered from a
fiber object (gradient directions, b-values, etc.) are filled with neutral
placeholders.
An S4 object of class dwiFiber (from dti).
if (requireNamespace("dti", quietly = TRUE)) { pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- bundle(streamlines = list(sl)) dfi <- as_dwifiber(b) class(dfi) # "dwiFiber" }if (requireNamespace("dti", quietly = TRUE)) { pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- bundle(streamlines = list(sl)) dfi <- as_dwifiber(b) class(dfi) # "dwiFiber" }
as_dwifiber() method for bundle objectsas_dwifiber() method for bundle objects
x |
A bundle object. |
... |
Additional arguments (currently unused). |
An S4 dwiFiber object.
as_dwifiber() method for streamline objectsas_dwifiber() method for streamline objects
x |
A streamline object. |
... |
Additional arguments (currently unused). |
An S4 dwiFiber object.
as_streamline() converts a supported object into a streamline.
as_streamline(x, ...)as_streamline(x, ...)
x |
An object to coerce. |
... |
Additional arguments (currently unused). |
Currently supported input classes:
streamline: returned unchanged.
dwiFiber (from dti): the object must contain exactly one
fiber. For multi-fiber objects use as_bundle() instead.
A streamline object.
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) identical(as_streamline(sl), sl) # TRUE — identity coercionpts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) identical(as_streamline(sl), sl) # TRUE — identity coercion
as_streamline() method for bundle objectsas_streamline() method for bundle objects
x |
A bundle object containing exactly one streamline. |
... |
Additional arguments (currently unused). |
The sole streamline inside x.
as_streamline() method for streamline objectsas_streamline() method for streamline objects
x |
A streamline object. |
... |
Additional arguments (currently unused). |
x unchanged.
Accepts any mix of named bundle objects (passed as name = bundle) or
bundle_set objects. All bundles are collected into a flat named list and
wrapped in a new bundle_set.
bind_bundle_sets(..., set_data = NULL)bind_bundle_sets(..., set_data = NULL)
... |
Named bundle objects or bundle_set objects to combine. Each bare bundle argument must be named so that its label in the resulting set is unambiguous. |
set_data |
A named list of set-level metadata to attach to the
resulting bundle_set. Defaults to the |
A bundle_set containing all input bundles.
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b1 <- bundle(streamlines = list(streamline(points = pts))) b2 <- bundle(streamlines = list(streamline(points = pts))) # two named bare bundles bs <- bind_bundle_sets("sub-01" = b1, "sub-02" = b2) bs@n_bundles # 2 bs@bundle_names # c("sub-01", "sub-02") # combine two bundle_sets bs1 <- bundle_set(list("sub-01" = b1)) bs2 <- bundle_set(list("sub-02" = b2)) bs_all <- bind_bundle_sets(bs1, bs2) bs_all@n_bundles # 2pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b1 <- bundle(streamlines = list(streamline(points = pts))) b2 <- bundle(streamlines = list(streamline(points = pts))) # two named bare bundles bs <- bind_bundle_sets("sub-01" = b1, "sub-02" = b2) bs@n_bundles # 2 bs@bundle_names # c("sub-01", "sub-02") # combine two bundle_sets bs1 <- bundle_set(list("sub-01" = b1)) bs2 <- bundle_set(list("sub-02" = b2)) bs_all <- bind_bundle_sets(bs1, bs2) bs_all@n_bundles # 2
Accepts any mix of streamline and bundle objects. All streamlines are
collected into a flat list and wrapped in a new bundle. bundle_data
from the first bundle argument (if any) is preserved; pass your own via
the bundle_data argument to override.
bind_bundles(..., bundle_data = NULL)bind_bundles(..., bundle_data = NULL)
... |
One or more streamline or bundle objects. |
bundle_data |
A named list of bundle-level metadata to attach to the
resulting bundle. Defaults to an empty list (or the |
A bundle containing all input streamlines.
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl1 <- streamline(points = pts) sl2 <- streamline(points = pts) b1 <- bundle(streamlines = list(sl1)) b2 <- bundle(streamlines = list(sl2)) # combine two bundles b_all <- bind_bundles(b1, b2) b_all@n_streamlines # 2 # mix a bundle and a loose streamline b_mixed <- bind_bundles(b1, sl2) b_mixed@n_streamlines # 2pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl1 <- streamline(points = pts) sl2 <- streamline(points = pts) b1 <- bundle(streamlines = list(sl1)) b2 <- bundle(streamlines = list(sl2)) # combine two bundles b_all <- bind_bundles(b1, b2) b_all@n_streamlines # 2 # mix a bundle and a loose streamline b_mixed <- bind_bundles(b1, sl2) b_mixed@n_streamlines # 2
A bundle is an ordered collection of streamline objects representing a
tractogram or white-matter bundle. It stores two compartments:
@streamlines — a list of streamline objects.
@bundle_data — a named list of bundle-level metadata (arbitrary R
objects, e.g. the affine transform used during tracking).
bundle(streamlines = list(), bundle_data = list())bundle(streamlines = list(), bundle_data = list())
streamlines |
A list of streamline objects. |
bundle_data |
A named list of bundle-level metadata. |
A bundle S7 object.
The following methods are defined for bundle objects:
format(x, ...): Returns a compact character string such as
<bundle [2 streamlines | 10–20 pts/streamline]>.
print(x, ...): Prints the formatted string to the console and
invisibly returns x.
length(x): Returns the number of streamlines (equivalent to
x@n_streamlines).
x[[i]]: Extracts the i-th streamline from the bundle.
x[i]: Returns a new bundle containing only the selected
streamlines, preserving @bundle_data.
@n_streamlinesAn integer scalar giving the number of streamlines in the bundle (read-only).
@bundle_attributesA character vector of the names of the bundle-level attributes (read-only).
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- bundle(streamlines = list(sl)) b@n_streamlines # 1 b@bundle_attributes # NULL (no bundle-level attributes) # bundle_data is stored b2 <- bundle( streamlines = list(sl), bundle_data = list(subject = "sub-01") ) b2@bundle_data$subject # "sub-01" # format(), print(), length() and indexing methods format(b2) print(b2) length(b2) # 1 b2[[1]] # first streamline # subsetting preserves bundle_data pts2 <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl2 <- streamline(points = pts2) b3 <- bundle(streamlines = list(sl, sl2), bundle_data = list(subject = "sub-01")) b3[1]@n_streamlines # 1, bundle_data preservedpts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- bundle(streamlines = list(sl)) b@n_streamlines # 1 b@bundle_attributes # NULL (no bundle-level attributes) # bundle_data is stored b2 <- bundle( streamlines = list(sl), bundle_data = list(subject = "sub-01") ) b2@bundle_data$subject # "sub-01" # format(), print(), length() and indexing methods format(b2) print(b2) length(b2) # 1 b2[[1]] # first streamline # subsetting preserves bundle_data pts2 <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl2 <- streamline(points = pts2) b3 <- bundle(streamlines = list(sl, sl2), bundle_data = list(subject = "sub-01")) b3[1]@n_streamlines # 1, bundle_data preserved
A bundle_set is a named collection of bundle objects, designed for
multi-subject or multi-session studies where each element represents one
subject's (or session's) tractogram. It stores two compartments:
@bundles — a named list of bundle objects. Names typically encode
subject or session identifiers (e.g. "sub-01", "sub-02").
@set_data — a named list of set-level metadata (arbitrary R objects,
e.g. study name, atlas used, acquisition protocol).
bundle_set(bundles = list(), set_data = list())bundle_set(bundles = list(), set_data = list())
bundles |
A named list of bundle objects. |
set_data |
A named list of set-level metadata. |
A bundle_set S7 object.
The following methods are defined for bundle_set objects:
format(x, ...): Returns a compact character string.
print(x, ...): Prints the formatted string to the console and
invisibly returns x.
length(x): Returns the number of bundles.
x[[i]]: Extracts the i-th (or named) bundle from the set.
x[i]: Returns a new bundle_set containing only the selected bundles,
preserving @set_data.
names(x): Returns the names of the bundles.
@n_bundlesAn integer scalar giving the number of bundles in the set (read-only).
@bundle_namesA character vector of the names of the bundles (read-only).
@set_attributesA character vector of the names of the set-level attributes (read-only).
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b1 <- bundle(streamlines = list(streamline(points = pts)), bundle_data = list(subject = "sub-01")) b2 <- bundle(streamlines = list(streamline(points = pts)), bundle_data = list(subject = "sub-02")) bs <- bundle_set(bundles = list("sub-01" = b1, "sub-02" = b2)) bs@n_bundles # 2 bs@bundle_names # c("sub-01", "sub-02") bs[["sub-01"]] # first bundlepts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b1 <- bundle(streamlines = list(streamline(points = pts)), bundle_data = list(subject = "sub-01")) b2 <- bundle(streamlines = list(streamline(points = pts)), bundle_data = list(subject = "sub-02")) bs <- bundle_set(bundles = list("sub-01" = b1, "sub-02" = b2)) bs@n_bundles # 2 bs@bundle_names # c("sub-01", "sub-02") bs[["sub-01"]] # first bundle
compute_hausdorff_distance() is an S7 generic that computes the symmetric
Hausdorff distance between streamline objects based on their 3-D
coordinate matrices, with methods available for the following classes:
The four dispatch cases are:
streamline + streamline: returns a single numeric scalar — the
symmetric Hausdorff distance between the two streamlines.
bundle + missing: returns a symmetric numeric distance matrix of
dimension , where is the number of streamlines in
the bundle, giving all pairwise Hausdorff distances.
bundle + streamline: returns a numeric vector of length
giving the Hausdorff distance from y to each streamline in x.
bundle + bundle: returns a symmetric numeric distance matrix of
dimension , treating the
concatenation of all streamlines from x and y as one collection.
compute_hausdorff_distance(x, y = NULL)compute_hausdorff_distance(x, y = NULL)
x |
A streamline or bundle object. |
y |
A streamline or bundle object, or |
A non-negative numeric scalar when both x and y are streamlines.
A dist object of size when x is a bundle
and y is NULL or a bundle (use as.matrix() to expand to a full
matrix).
A numeric vector of length when x is a bundle and y is
a streamline.
pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") sl1 <- streamline(points = pts1) pts2 <- matrix(runif(30), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl2 <- streamline(points = pts2) # streamline x streamline -> scalar compute_hausdorff_distance(sl1, sl2) # bundle x missing -> pairwise dist object b <- bundle(streamlines = list(sl1, sl2)) compute_hausdorff_distance(b) as.matrix(compute_hausdorff_distance(b)) # bundle x streamline -> vector compute_hausdorff_distance(b, sl1) # bundle x bundle -> combined pairwise matrix b2 <- bundle(streamlines = list(sl2)) compute_hausdorff_distance(b, b2)pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") sl1 <- streamline(points = pts1) pts2 <- matrix(runif(30), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl2 <- streamline(points = pts2) # streamline x streamline -> scalar compute_hausdorff_distance(sl1, sl2) # bundle x missing -> pairwise dist object b <- bundle(streamlines = list(sl1, sl2)) compute_hausdorff_distance(b) as.matrix(compute_hausdorff_distance(b)) # bundle x streamline -> vector compute_hausdorff_distance(b, sl1) # bundle x bundle -> combined pairwise matrix b2 <- bundle(streamlines = list(sl2)) compute_hausdorff_distance(b, b2)
compute_hausdorff_distance() method for bundle objectsDispatches to one of three behaviours depending on y:
x |
A bundle object. |
y |
|
y = NULL: pairwise distances within x as a dist
object of size , computed in C++ via a single linear loop.
y is a streamline: numeric vector of distances from y to each
streamline in x.
y is a bundle: dist object for the concatenation
of all streamlines from x and y.
A dist object of size x@n_streamlines when y is
NULL or a bundle. The lower triangle stores all pairwise symmetric
Hausdorff distances (computed in C++). Use as.matrix() to obtain the
full matrix.
A numeric vector of length x@n_streamlines when y is a streamline.
pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") pts2 <- matrix(runif(30), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl1 <- streamline(points = pts1) sl2 <- streamline(points = pts2) b <- bundle(streamlines = list(sl1, sl2)) # pairwise dist object (size 2) compute_hausdorff_distance(b) as.matrix(compute_hausdorff_distance(b)) # distances from sl1 to each streamline in b compute_hausdorff_distance(b, sl1)pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") pts2 <- matrix(runif(30), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl1 <- streamline(points = pts1) sl2 <- streamline(points = pts2) b <- bundle(streamlines = list(sl1, sl2)) # pairwise dist object (size 2) compute_hausdorff_distance(b) as.matrix(compute_hausdorff_distance(b)) # distances from sl1 to each streamline in b compute_hausdorff_distance(b, sl1)
compute_hausdorff_distance() method for two streamline objectscompute_hausdorff_distance() method for two streamline objects
x |
A streamline object. |
y |
A streamline object. |
A non-negative numeric scalar equal to
, where
is the
directed Hausdorff distance. The core computation is performed in C++
via hausdorff_distance_cpp().
pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") pts2 <- matrix(runif(30), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl1 <- streamline(points = pts1) sl2 <- streamline(points = pts2) compute_hausdorff_distance(sl1, sl2)pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") pts2 <- matrix(runif(30), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl1 <- streamline(points = pts1) sl2 <- streamline(points = pts2) compute_hausdorff_distance(sl1, sl2)
compute_hausdorff_distance() catch-all methodcompute_hausdorff_distance() catch-all method
Does not return a value; always throws an error for unsupported input types.
get_curvature() is function that computes the curvature of a
streamline object. The curvature at each point
along the arc-length abscissa is computed using cubic smoothing
splines (3 degrees of freedom per component).
get_curvature(x)get_curvature(x)
x |
A streamline object. |
A non-negative numeric vector of length x@n_points giving the
curvature at each sampled point along the streamline.
Higher values indicate sharper bending at that location.
pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_curvature(sl)pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_curvature(sl)
get_curvilinear_length() is a function that computes the total arc-length of a
streamline object as the sum of Euclidean segment lengths between
consecutive points.
get_curvilinear_length(x)get_curvilinear_length(x)
x |
A streamline object. |
A non-negative numeric scalar.
pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_curvilinear_length(sl)pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_curvilinear_length(sl)
get_euclidean_length() is a function that computes the Euclidean
(straight-line) distance of a streamline object.
get_euclidean_length(x)get_euclidean_length(x)
x |
A streamline object. |
A non-negative numeric scalar.
pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_euclidean_length(sl)pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_euclidean_length(sl)
get_sinuosity() is a function that computes the ratio of curvilinear
length to Euclidean length for a streamline object, with a value of
1 indicating a perfectly straight streamline and larger values indicating
greater curviness.
get_sinuosity(x)get_sinuosity(x)
x |
A streamline object. |
A numeric scalar .
pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_sinuosity(sl)pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_sinuosity(sl)
get_torsion() is function that computes the torsion of a
streamline object. The torsion at each point
along the arc-length abscissa is computed using cubic smoothing
splines (4 degrees of freedom per component).
get_torsion(x)get_torsion(x)
x |
A streamline object. |
A numeric vector of length x@n_points giving the torsion
at each sampled point along the streamline. Positive values
indicate right-handed twisting; negative values indicate left-handed
twisting; zero indicates a planar curve at that location.
pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_torsion(sl)pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) get_torsion(sl)
Test whether an object is a bundle
is_bundle(x)is_bundle(x)
x |
An object. |
TRUE if x is of class bundle, otherwise FALSE.
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- bundle(streamlines = list(sl)) is_bundle(b) # TRUE is_bundle(sl) # FALSEpts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) b <- bundle(streamlines = list(sl)) is_bundle(b) # TRUE is_bundle(sl) # FALSE
Test whether an object is a bundle_set
is_bundle_set(x)is_bundle_set(x)
x |
An object. |
TRUE if x is of class bundle_set, otherwise FALSE.
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b <- bundle(streamlines = list(streamline(points = pts))) bs <- bundle_set(bundles = list("sub-01" = b)) is_bundle_set(bs) # TRUE is_bundle_set(b) # FALSEpts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) b <- bundle(streamlines = list(streamline(points = pts))) bs <- bundle_set(bundles = list("sub-01" = b)) is_bundle_set(bs) # TRUE is_bundle_set(b) # FALSE
Test whether an object is a streamline
is_streamline(x)is_streamline(x)
x |
An object. |
TRUE if x is of class streamline, otherwise FALSE.
pts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) is_streamline(sl) # TRUE is_streamline(42) # FALSEpts <- matrix(runif(15), ncol = 3, dimnames = list(NULL, c("X", "Y", "Z"))) sl <- streamline(points = pts) is_streamline(sl) # TRUE is_streamline(42) # FALSE
reparametrize() is an S7 generic that resamples the 3-D coordinates (and
any @point_data attributes) of a tractography object onto a uniform
arc-length grid using linear interpolation, with methods available for the
following classes:
reparametrize(x, n_points = NULL)reparametrize(x, n_points = NULL)
x |
A streamline or bundle object. |
n_points |
Number of equally-spaced arc-length points to use.
Pass |
An object of the same class as x reparametrized onto the new
grid.
# reparametrize a single streamline to 10 points pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) sl_reparam <- reparametrize(sl, n_points = 10) # reparametrize a bundle to the mean number of points across its streamlines sl1 <- streamline(points = pts) pts2 <- matrix(runif(60), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl2 <- streamline(points = pts2) b <- bundle(streamlines = list(sl1, sl2)) bundle_reparam <- reparametrize(b)# reparametrize a single streamline to 10 points pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts) sl_reparam <- reparametrize(sl, n_points = 10) # reparametrize a bundle to the mean number of points across its streamlines sl1 <- streamline(points = pts) pts2 <- matrix(runif(60), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") sl2 <- streamline(points = pts2) b <- bundle(streamlines = list(sl1, sl2)) bundle_reparam <- reparametrize(b)
reparametrize() method for bundle objectsResamples every streamline inside a bundle onto a common uniform
arc-length grid. See reparametrize() for the full parameter documentation.
x |
A bundle object. |
n_points |
Number of equally-spaced arc-length points to use.
Pass |
A bundle reparametrized onto the new grid. Every streamline in the
returned bundle has exactly n_points rows in @points (defaulting to
the rounded mean number of points across all streamlines when n_points
is NULL).
pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") pts2 <- matrix(runif(60), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") b <- bundle(streamlines = list(streamline(points = pts1), streamline(points = pts2))) b_reparam <- reparametrize(b, n_points = 15) b_reparam[[1]]@n_points # 15 b_reparam[[2]]@n_points # 15pts1 <- matrix(runif(30), ncol = 3) colnames(pts1) <- c("X", "Y", "Z") pts2 <- matrix(runif(60), ncol = 3) colnames(pts2) <- c("X", "Y", "Z") b <- bundle(streamlines = list(streamline(points = pts1), streamline(points = pts2))) b_reparam <- reparametrize(b, n_points = 15) b_reparam[[1]]@n_points # 15 b_reparam[[2]]@n_points # 15
reparametrize() method for streamline objectsResamples a single streamline onto a uniform arc-length grid. See
reparametrize() for the full parameter documentation.
x |
A streamline object. |
n_points |
Number of equally-spaced arc-length points to use.
Pass |
A streamline reparametrized onto the new grid. The returned
object has the same class as the input but with @points resampled to
exactly n_points rows and all @point_data vectors resampled
correspondingly via linear interpolation.
pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts, point_data = list(FA = runif(10))) sl_reparam <- reparametrize(sl, n_points = 20) sl_reparam@n_points # 20pts <- matrix(runif(30), ncol = 3) colnames(pts) <- c("X", "Y", "Z") sl <- streamline(points = pts, point_data = list(FA = runif(10))) sl_reparam <- reparametrize(sl, n_points = 20) sl_reparam@n_points # 20
A streamline represents a single fibre tract. It stores three data
compartments that mirror the conceptual levels found in tractography file
formats:
@points — an numeric matrix whose columns are named
"X", "Y", and "Z", holding the ordered 3-D coordinates of the
points along the tract.
@point_data — a named list of numeric vectors, each of length ,
holding per-point scalar attributes (e.g. fractional anisotropy sampled
at every point).
@streamline_data — a named list of numeric scalars (length-1 vectors)
holding per-streamline attributes (e.g. a tract-level weight or mean FA).
streamline(points = NULL, point_data = list(), streamline_data = list())streamline(points = NULL, point_data = list(), streamline_data = list())
points |
A numeric matrix with columns |
point_data |
A named list of per-point numeric vectors. |
streamline_data |
A named list of per-streamline numeric scalars. |
A streamline S7 object.
The following methods are defined for streamline objects:
format(x, ...): Returns a compact character string such as
<streamline [10 pts] | point: FA>.
print(x, ...): Prints the formatted string to the console and
invisibly returns x.
@n_pointsAn integer scalar giving the number of points in the streamline (read-only).
@point_attributesA character vector of the names of the per-point attributes (read-only).
@streamline_attributesA character vector of the names of the per-streamline attributes (read-only).
# Create a streamline with 5 points and some attributes sl <- streamline( points = matrix( c(0, 0, 0, 1, 0, 0, 1, 1, 0, 1, 1, 1, 0, 1, 1), ncol = 3, byrow = TRUE, dimnames = list(NULL, c("X", "Y", "Z")) ), point_data = list(FA = c(0.5, 0.6, 0.7, 0.8, 0.9)), streamline_data = list(mean_FA = 0.7) ) sl@n_points # 5 sl@point_attributes # "FA" sl@streamline_attributes # "mean_FA" # format() and print() methods format(sl) # "<streamline [5 pts] | point: FA | streamline: mean_FA>" print(sl)# Create a streamline with 5 points and some attributes sl <- streamline( points = matrix( c(0, 0, 0, 1, 0, 0, 1, 1, 0, 1, 1, 1, 0, 1, 1), ncol = 3, byrow = TRUE, dimnames = list(NULL, c("X", "Y", "Z")) ), point_data = list(FA = c(0.5, 0.6, 0.7, 0.8, 0.9)), streamline_data = list(mean_FA = 0.7) ) sl@n_points # 5 sl@point_attributes # "FA" sl@streamline_attributes # "mean_FA" # format() and print() methods format(sl) # "<streamline [5 pts] | point: FA | streamline: mean_FA>" print(sl)