Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Adding variance to scil_visualize_fodf #671

Merged
merged 8 commits into from
Mar 7, 2023
Merged
Show file tree
Hide file tree
Changes from 7 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
55 changes: 47 additions & 8 deletions scilpy/viz/scene_utils.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@
import numpy as np

import vtk
from dipy.reconst.shm import sh_to_sf_matrix
from dipy.reconst.shm import sh_to_sf_matrix, sh_to_sf
from fury import window, actor
from fury.colormap import distinguishable_colormap
from PIL import Image
Expand Down Expand Up @@ -122,7 +122,8 @@ def set_display_extent(slicer_actor, orientation, volume_shape, slice_index):

def create_odf_slicer(sh_fodf, orientation, slice_index, mask, sphere,
nb_subdivide, sh_order, sh_basis, full_basis,
scale, radial_scale, norm, colormap):
scale, radial_scale, norm, colormap, sh_variance=None,
variance_k=1, variance_color=(255, 255, 255)):
"""
Create a ODF slicer actor displaying a fODF slice. The input volume is a
3-dimensional grid containing the SH coefficients of the fODF for each
Expand Down Expand Up @@ -158,6 +159,12 @@ def create_odf_slicer(sh_fodf, orientation, slice_index, mask, sphere,
If True, enables normalization of ODF slicer.
colormap : str
Colormap for the ODF slicer. If None, a RGB colormap is used.
sh_variance : np.ndarray, optional
Spherical harmonics of the variance fODF data.
variance_k : float, optional
Factor that multiplies sqrt(variance).
variance_color : tuple, optional
Color of the variance fODF data, in RGB.

Returns
-------
Expand All @@ -172,14 +179,46 @@ def create_odf_slicer(sh_fodf, orientation, slice_index, mask, sphere,
B_mat = sh_to_sf_matrix(sphere, sh_order, sh_basis,
full_basis, return_inv=False)

odf_actor = actor.odf_slicer(sh_fodf, mask=mask, norm=norm,
radial_scale=radial_scale,
sphere=sphere,
colormap=colormap,
scale=scale, B_matrix=B_mat)
var_actor = None

if sh_variance is not None:
fodf = sh_to_sf(sh_fodf, sphere, sh_order, sh_basis,
full_basis=full_basis)
fodf_var = sh_to_sf(sh_variance, sphere, sh_order, sh_basis,
full_basis=full_basis)
fodf_uncertainty = fodf + variance_k * np.sqrt(np.clip(fodf_var, 0,
None))
# normalise fodf and variance
if norm:
maximums = np.abs(np.append(fodf, fodf_uncertainty, axis=-1))\
.max(axis=-1)
fodf[maximums > 0] /= maximums[maximums > 0][..., None]
fodf_uncertainty[maximums > 0] /= maximums[maximums > 0][..., None]

odf_actor = actor.odf_slicer(fodf, mask=mask, norm=False,
radial_scale=radial_scale,
sphere=sphere, scale=scale,
colormap=colormap)

var_actor = actor.odf_slicer(fodf_uncertainty, mask=mask, norm=False,
radial_scale=radial_scale,
sphere=sphere, scale=scale,
colormap=variance_color)
var_actor.GetProperty().SetDiffuse(0.0)
var_actor.GetProperty().SetAmbient(1.0)
var_actor.GetProperty().SetFrontfaceCulling(True)
else:
odf_actor = actor.odf_slicer(sh_fodf, mask=mask, norm=norm,
radial_scale=radial_scale,
sphere=sphere,
colormap=colormap,
scale=scale, B_matrix=B_mat)
set_display_extent(odf_actor, orientation, sh_fodf.shape[:3], slice_index)
if var_actor is not None:
set_display_extent(var_actor, orientation,
fodf_uncertainty.shape[:3], slice_index)

return odf_actor
return odf_actor, var_actor


def _get_affine_for_texture(orientation, offset):
Expand Down
49 changes: 41 additions & 8 deletions scripts/scil_visualize_fodf.py
Original file line number Diff line number Diff line change
Expand Up @@ -146,6 +146,20 @@ def _build_arg_parser():
help='Length of the peaks segments. '
'[%(default)s]')

# fODF variance options
p.add_argument('--variance',
help='FODF variance file. For the visualization of fodf '
'uncertainty, this variance will be used as follow: '
'mean + k * sqrt(variance), where mean is the input '
'fodf (in_fodf) and k is the scaling factor '
'(variance_k).')
p.add_argument('--variance_k', default=1, type=float,
help='Scaling factor (k) for the computation of the fodf '
'uncertainty. [%(default)s]')
p.add_argument('--var_color', nargs=3, type=float, default=(1, 1, 1),
help='Color of variance outline. Must be RGB values scaled '
'between 0 and 1. [%(default)s]')

return p


Expand Down Expand Up @@ -222,6 +236,16 @@ def _get_data_from_inputs(args):
peak_vals =\
nib.nifti1.load(args.peaks_values).get_fdata(dtype=np.float32)
data['peaks_values'] = peak_vals
if args.variance:
assert_same_resolution([args.variance, args.in_fodf])
variance = nib.nifti1.load(args.variance).get_fdata(dtype=np.float32)
if len(variance.shape) == 3:
variance = np.reshape(variance, variance.shape + (1,))
if variance.shape != fodf.shape:
raise ValueError('Dimensions mismatch between fODF {0} and '
'variance {1}.'
.format(fodf.shape, variance.shape))
data['variance'] = variance

return data

Expand All @@ -246,17 +270,26 @@ def main():
else:
color_rgb = None

variance = data['variance'] if args.variance else None
var_color = np.asarray(args.var_color) * 255
# Instantiate the ODF slicer actor
odf_actor = create_odf_slicer(data['fodf'], args.axis_name,
args.slice_index, mask, sph,
args.sph_subdivide, sh_order,
args.sh_basis, full_basis,
args.scale,
not args.radial_scale_off,
not args.norm_off,
args.colormap or color_rgb)
odf_actor, var_actor = create_odf_slicer(data['fodf'], args.axis_name,
args.slice_index, mask, sph,
args.sph_subdivide, sh_order,
args.sh_basis, full_basis,
args.scale,
not args.radial_scale_off,
not args.norm_off,
args.colormap or color_rgb,
sh_variance=variance,
variance_k=args.variance_k,
variance_color=var_color)
actors.append(odf_actor)

# Instantiate a variance slicer actor if a variance image is supplied
if 'variance' in data:
actors.append(var_actor)

# Instantiate a texture slicer actor if a background image is supplied
if 'bg' in data:
bg_actor = create_texture_slicer(data['bg'],
Expand Down