actor

apply_affine(aff, pts) Apply affine matrix aff to points pts.
axes([scale, colorx, colory, colorz, opacity]) Create an actor with the coordinate’s system axes where red = x, green = y, blue = z.
colormap_lookup_table([scale_range, …]) Lookup table for the colormap.
contour_from_roi(data[, affine, color, opacity]) Generates surface actor from a binary ROI.
create_colormap(v[, name, auto]) Create colors from a specific colormap and return it as an array of shape (N,3) where every row gives the corresponding r,g,b value.
dots(points[, color, opacity, dot_size]) Create one or more 3d points
label([text, pos, scale, color]) Create a label actor.
line(lines[, colors, opacity, linewidth, …]) Create an actor for one or more lines.
lines_to_vtk_polydata(lines[, colors]) Create a vtkPolyData with lines and colors
numpy_to_vtk_colors(colors) Numpy color array to a vtk color array
numpy_to_vtk_points(points) Numpy points array to a vtk points array
odf_slicer(odfs[, affine, mask, sphere, …]) Slice spherical fields in native or world coordinates
peak_slicer(peaks_dirs[, peaks_values, …]) Visualize peak directions as given from peaks_from_model
point(points, colors[, opacity, …]) Visualize points as sphere glyphs
scalar_bar([lookup_table, title]) Default scalar bar actor for a given colormap (colorbar)
set_input(vtk_object, inp) Generic input function which takes into account VTK 5 or 6
set_polydata_triangles(polydata, triangles) set polydata triangles with a numpy array (ndarrays Nx3 int)
set_polydata_vertices(polydata, vertices) set polydata vertices with a numpy array (ndarrays Nx3 int)
slicer(data[, affine, value_range, opacity, …]) Cuts 3D scalar or rgb volumes into 2D images
sphere(centers, colors[, radii, theta, phi, …]) Visualize one or many spheres with different colors and radii
streamtube(lines[, colors, opacity, …]) Uses streamtubes to visualize polylines
tensor_slicer(evals, evecs[, affine, mask, …]) Slice many tensors as ellipsoids in native or world coordinates

apply_affine

fury.actor.apply_affine(aff, pts)[source]

Apply affine matrix aff to points pts.

Returns result of application of aff to the right of pts. The coordinate dimension of pts should be the last. For the 3D case, aff will be shape (4,4) and pts will have final axis length 3 - maybe it will just be N by 3. The return value is the transformed points, in this case:: res = np.dot(aff[:3,:3], pts.T) + aff[:3,3:4] transformed_pts = res.T This routine is more general than 3D, in that aff can have any shape (N,N), and pts can have any shape, as long as the last dimension is for the coordinates, and is therefore length N-1.

Parameters:
aff : (N, N) array-like

Homogenous affine, for 3D points, will be 4 by 4. Contrary to first appearance, the affine will be applied on the left of pts.

pts : (…, N-1) array-like

Points, where the last dimension contains the coordinates of each point. For 3D, the last dimension will be length 3.

Returns:
transformed_pts : (…, N-1) array

transformed points

Notes

Copied from nibabel to remove dependency.

Examples

>>> aff = np.array([[0,2,0,10],[3,0,0,11],[0,0,4,12],[0,0,0,1]])
>>> pts = np.array([[1,2,3],[2,3,4],[4,5,6],[6,7,8]])
>>> apply_affine(aff, pts) #doctest: +ELLIPSIS
array([[14, 14, 24],
       [16, 17, 28],
       [20, 23, 36],
       [24, 29, 44]]...)
Just to show that in the simple 3D case, it is equivalent to:
>>> (np.dot(aff[:3,:3], pts.T) + aff[:3,3:4]).T #doctest: +ELLIPSIS
array([[14, 14, 24],
       [16, 17, 28],
       [20, 23, 36],
       [24, 29, 44]]...)
But `pts` can be a more complicated shape:
>>> pts = pts.reshape((2,2,3))
>>> apply_affine(aff, pts) #doctest: +ELLIPSIS
array([[[14, 14, 24],
        [16, 17, 28]],
<BLANKLINE>
       [[20, 23, 36],
        [24, 29, 44]]]...)

axes

fury.actor.axes(scale=(1, 1, 1), colorx=(1, 0, 0), colory=(0, 1, 0), colorz=(0, 0, 1), opacity=1)[source]

Create an actor with the coordinate’s system axes where red = x, green = y, blue = z.

Parameters:
scale : tuple (3,)

Axes size e.g. (100, 100, 100). Default is (1, 1, 1).

colorx : tuple (3,)

x-axis color. Default red (1, 0, 0).

colory : tuple (3,)

y-axis color. Default green (0, 1, 0).

colorz : tuple (3,)

z-axis color. Default blue (0, 0, 1).

opacity : float, optional

Takes values from 0 (fully transparent) to 1 (opaque). Default is 1.

Returns:
vtkAssembly

colormap_lookup_table

fury.actor.colormap_lookup_table(scale_range=(0, 1), hue_range=(0.8, 0), saturation_range=(1, 1), value_range=(0.8, 0.8))[source]

Lookup table for the colormap.

Parameters:
scale_range : tuple

It can be anything e.g. (0, 1) or (0, 255). Usually it is the mininum and maximum value of your data. Default is (0, 1).

hue_range : tuple of floats

HSV values (min 0 and max 1). Default is (0.8, 0).

saturation_range : tuple of floats

HSV values (min 0 and max 1). Default is (1, 1).

value_range : tuple of floats

HSV value (min 0 and max 1). Default is (0.8, 0.8).

Returns:
lookup_table : vtkLookupTable

contour_from_roi

fury.actor.contour_from_roi(data, affine=None, color=array([1, 0, 0]), opacity=1)[source]

Generates surface actor from a binary ROI.

The color and opacity of the surface can be customized.

Parameters:
data : array, shape (X, Y, Z)

An ROI file that will be binarized and displayed.

affine : array, shape (4, 4)

Grid to space (usually RAS 1mm) transformation matrix. Default is None. If None then the identity matrix is used.

color : (1, 3) ndarray

RGB values in [0,1].

opacity : float

Opacity of surface between 0 and 1.

Returns:
contour_assembly : vtkAssembly

ROI surface object displayed in space coordinates as calculated by the affine parameter.

create_colormap

fury.actor.create_colormap(v, name='plasma', auto=True)[source]

Create colors from a specific colormap and return it as an array of shape (N,3) where every row gives the corresponding r,g,b value. The colormaps we use are similar with those of matplotlib.

Parameters:
v : (N,) array

vector of values to be mapped in RGB colors according to colormap

name : str.

Name of the colormap. Currently implemented: ‘jet’, ‘blues’, ‘accent’, ‘bone’ and matplotlib colormaps if you have matplotlib installed. For example, we suggest using ‘plasma’, ‘viridis’ or ‘inferno’. ‘jet’ is popular but can be often misleading and we will deprecate it the future.

auto : bool,

if auto is True then v is interpolated to [0, 10] from v.min() to v.max()

Notes

FURY supports a few colormaps for those who do not use Matplotlib, for more colormaps consider downloading Matplotlib (see matplotlib.org).

dots

fury.actor.dots(points, color=(1, 0, 0), opacity=1, dot_size=5)[source]

Create one or more 3d points

Parameters:
points : ndarray, (N, 3)
color : tuple (3,)
opacity : float, optional

Takes values from 0 (fully transparent) to 1 (opaque)

dot_size : int
Returns:
vtkActor

See also

fury.actor.point

label

fury.actor.label(text='Origin', pos=(0, 0, 0), scale=(0.2, 0.2, 0.2), color=(1, 1, 1))[source]

Create a label actor.

This actor will always face the camera

Parameters:
text : str

Text for the label.

pos : (3,) array_like, optional

Left down position of the label.

scale : (3,) array_like

Changes the size of the label.

color : (3,) array_like

Label color as (r,g,b) tuple.

Returns:
l : vtkActor object

Label.

Examples

>>> from fury import window, actor
>>> ren = window.Renderer()
>>> l = actor.label(text='Hello')
>>> ren.add(l)
>>> #window.show(ren)

line

fury.actor.line(lines, colors=None, opacity=1, linewidth=1, spline_subdiv=None, lod=True, lod_points=10000, lod_points_size=3, lookup_colormap=None)[source]

Create an actor for one or more lines.

Parameters:
lines : list of arrays
colors : array (N, 3), list of arrays, tuple (3,), array (K,), None

If None then a standard orientation colormap is used for every line. If one tuple of color is used. Then all streamlines will have the same colour. If an array (N, 3) is given, where N is equal to the number of lines. Then every line is coloured with a different RGB color. If a list of RGB arrays is given then every point of every line takes a different color. If an array (K, ) is given, where K is the number of points of all lines then these are considered as the values to be used by the colormap. If an array (L, ) is given, where L is the number of streamlines then these are considered as the values to be used by the colormap per streamline. If an array (X, Y, Z) or (X, Y, Z, 3) is given then the values for the colormap are interpolated automatically using trilinear interpolation.

opacity : float, optional

Takes values from 0 (fully transparent) to 1 (opaque). Default is 1.

linewidth : float, optional

Line thickness. Default is 1.

spline_subdiv : int, optional

Number of splines subdivision to smooth streamtubes. Default is None which means no subdivision.

lod : bool

Use vtkLODActor(level of detail) rather than vtkActor. Default is True. Level of detail actors do not render the full geometry when the frame rate is low.

lod_points : int

Number of points to be used when LOD is in effect. Default is 10000.

lod_points_size : int

Size of points when lod is in effect. Default is 3.

lookup_colormap : bool, optional

Add a default lookup table to the colormap. Default is None which calls fury.actor.colormap_lookup_table().

Returns:
v : vtkActor or vtkLODActor object

Line.

Examples

>>> from fury import actor, window
>>> ren = window.Renderer()
>>> lines = [np.random.rand(10, 3), np.random.rand(20, 3)]
>>> colors = np.random.rand(2, 3)
>>> c = actor.line(lines, colors)
>>> ren.add(c)
>>> #window.show(ren)

lines_to_vtk_polydata

fury.actor.lines_to_vtk_polydata(lines, colors=None)[source]

Create a vtkPolyData with lines and colors

Parameters:
lines : list

list of N curves represented as 2D ndarrays

colors : array (N, 3), list of arrays, tuple (3,), array (K,), None

If None then a standard orientation colormap is used for every line. If one tuple of color is used. Then all streamlines will have the same colour. If an array (N, 3) is given, where N is equal to the number of lines. Then every line is coloured with a different RGB color. If a list of RGB arrays is given then every point of every line takes a different color. If an array (K, 3) is given, where K is the number of points of all lines then every point is colored with a different RGB color. If an array (K,) is given, where K is the number of points of all lines then these are considered as the values to be used by the colormap. If an array (L,) is given, where L is the number of streamlines then these are considered as the values to be used by the colormap per streamline. If an array (X, Y, Z) or (X, Y, Z, 3) is given then the values for the colormap are interpolated automatically using trilinear interpolation.

Returns:
poly_data : vtkPolyData
is_colormap : bool, true if the input color array was a colormap

numpy_to_vtk_colors

fury.actor.numpy_to_vtk_colors(colors)[source]

Numpy color array to a vtk color array

Parameters:
colors: ndarray
Returns:
vtk_colors : vtkDataArray

Notes

If colors are not already in UNSIGNED_CHAR you may need to multiply by 255.

Examples

>>> import numpy as np
>>> from fury.utils import numpy_to_vtk_colors
>>> rgb_array = np.random.rand(100, 3)
>>> vtk_colors = numpy_to_vtk_colors(255 * rgb_array)

numpy_to_vtk_points

fury.actor.numpy_to_vtk_points(points)[source]

Numpy points array to a vtk points array

Parameters:
points : ndarray
Returns:
vtk_points : vtkPoints()

odf_slicer

fury.actor.odf_slicer(odfs, affine=None, mask=None, sphere=None, scale=2.2, norm=True, radial_scale=True, opacity=1.0, colormap='plasma', global_cm=False)[source]

Slice spherical fields in native or world coordinates

Parameters:
odfs : ndarray

4D array of spherical functions

affine : array

4x4 transformation array from native coordinates to world coordinates

mask : ndarray

3D mask

sphere : Sphere

a sphere

scale : float

Distance between spheres.

norm : bool

Normalize sphere_values.

radial_scale : bool

Scale sphere points according to odf values.

opacity : float

Takes values from 0 (fully transparent) to 1 (opaque). Default is 1.

colormap : None or str

If None then white color is used. Otherwise the name of colormap is given. Matplotlib colormaps are supported (e.g., ‘inferno’).

global_cm : bool

If True the colormap will be applied in all ODFs. If False it will be applied individually at each voxel (default False).

Returns:
actor : vtkActor

Spheres

peak_slicer

fury.actor.peak_slicer(peaks_dirs, peaks_values=None, mask=None, affine=None, colors=(1, 0, 0), opacity=1.0, linewidth=1, lod=False, lod_points=10000, lod_points_size=3)[source]

Visualize peak directions as given from peaks_from_model

Parameters:
peaks_dirs : ndarray

Peak directions. The shape of the array can be (M, 3) or (X, M, 3) or (X, Y, M, 3) or (X, Y, Z, M, 3)

peaks_values : ndarray

Peak values. The shape of the array can be (M, ) or (X, M) or (X, Y, M) or (X, Y, Z, M)

affine : array

4x4 transformation array from native coordinates to world coordinates

mask : ndarray

3D mask

colors : tuple or None

Default red color. If None then every peak gets an orientation color in similarity to a DEC map.

opacity : float, optional

Takes values from 0 (fully transparent) to 1 (opaque)

linewidth : float, optional

Line thickness. Default is 1.

lod : bool

Use vtkLODActor(level of detail) rather than vtkActor. Default is False. Level of detail actors do not render the full geometry when the frame rate is low.

lod_points : int

Number of points to be used when LOD is in effect. Default is 10000.

lod_points_size : int

Size of points when lod is in effect. Default is 3.

Returns:
vtkActor

point

fury.actor.point(points, colors, opacity=1.0, point_radius=0.1, theta=8, phi=8)[source]

Visualize points as sphere glyphs

Parameters:
points : ndarray, shape (N, 3)
colors : ndarray (N,3) or tuple (3,)
point_radius : float
theta : int
phi : int
opacity : float, optional

Takes values from 0 (fully transparent) to 1 (opaque)

Returns:
vtkActor

Examples

>>> from fury import window, actor
>>> ren = window.Renderer()
>>> pts = np.random.rand(5, 3)
>>> point_actor = actor.point(pts, window.colors.coral)
>>> ren.add(point_actor)
>>> # window.show(ren)

scalar_bar

fury.actor.scalar_bar(lookup_table=None, title=' ')[source]

Default scalar bar actor for a given colormap (colorbar)

Parameters:
lookup_table : vtkLookupTable or None

If None then colormap_lookup_table is called with default options.

title : str
Returns:
scalar_bar : vtkScalarBarActor

set_input

fury.actor.set_input(vtk_object, inp)[source]

Generic input function which takes into account VTK 5 or 6

Parameters:
vtk_object: vtk object
inp: vtkPolyData or vtkImageData or vtkAlgorithmOutput
Returns:
vtk_object

Notes

This can be used in the following way::
from fury.utils import set_input poly_mapper = set_input(vtk.vtkPolyDataMapper(), poly_data)

set_polydata_triangles

fury.actor.set_polydata_triangles(polydata, triangles)[source]

set polydata triangles with a numpy array (ndarrays Nx3 int)

Parameters:
polydata : vtkPolyData
triangles : array (N, 3)

triangles, represented as 2D ndarrays (Nx3)

set_polydata_vertices

fury.actor.set_polydata_vertices(polydata, vertices)[source]

set polydata vertices with a numpy array (ndarrays Nx3 int)

Parameters:
polydata : vtkPolyData
vertices : vertices, represented as 2D ndarrays (Nx3)

slicer

fury.actor.slicer(data, affine=None, value_range=None, opacity=1.0, lookup_colormap=None, interpolation='linear', picking_tol=0.025)[source]

Cuts 3D scalar or rgb volumes into 2D images

Parameters:
data : array, shape (X, Y, Z) or (X, Y, Z, 3)

A grayscale or rgb 4D volume as a numpy array.

affine : array, shape (4, 4)

Grid to space (usually RAS 1mm) transformation matrix. Default is None. If None then the identity matrix is used.

value_range : None or tuple (2,)

If None then the values will be interpolated from (data.min(), data.max()) to (0, 255). Otherwise from (value_range[0], value_range[1]) to (0, 255).

opacity : float, optional

Opacity of 0 means completely transparent and 1 completely visible.

lookup_colormap : vtkLookupTable

If None (default) then a grayscale map is created.

interpolation : string

If ‘linear’ (default) then linear interpolation is used on the final texture mapping. If ‘nearest’ then nearest neighbor interpolation is used on the final texture mapping.

picking_tol : float

The tolerance for the vtkCellPicker, specified as a fraction of rendering window size.

Returns:
image_actor : ImageActor

An object that is capable of displaying different parts of the volume as slices. The key method of this object is display_extent where one can input grid coordinates and display the slice in space (or grid) coordinates as calculated by the affine parameter.

sphere

fury.actor.sphere(centers, colors, radii=1.0, theta=16, phi=16, vertices=None, faces=None)[source]

Visualize one or many spheres with different colors and radii

Parameters:
centers : ndarray, shape (N, 3)
colors : ndarray (N,3) or (N, 4) or tuple (3,) or tuple (4,)

RGB or RGBA (for opacity) R, G, B and A should be at the range [0, 1]

radii : float or ndarray, shape (N,)
theta : int
phi : int
vertices : ndarray, shape (N, 3)
faces : ndarray, shape (M, 3)

If faces is None then a sphere is created based on theta and phi angles If not then a sphere is created with the provided vertices and faces.

Returns:
vtkActor

Examples

>>> from fury import window, actor
>>> ren = window.Renderer()
>>> centers = np.random.rand(5, 3)
>>> sphere_actor = actor.sphere(centers, window.colors.coral)
>>> ren.add(sphere_actor)
>>> # window.show(ren)

streamtube

fury.actor.streamtube(lines, colors=None, opacity=1, linewidth=0.1, tube_sides=9, lod=True, lod_points=10000, lod_points_size=3, spline_subdiv=None, lookup_colormap=None)[source]

Uses streamtubes to visualize polylines

Parameters:
lines : list

list of N curves represented as 2D ndarrays

colors : array (N, 3), list of arrays, tuple (3,), array (K,), None

If None then a standard orientation colormap is used for every line. If one tuple of color is used. Then all streamlines will have the same colour. If an array (N, 3) is given, where N is equal to the number of lines. Then every line is coloured with a different RGB color. If a list of RGB arrays is given then every point of every line takes a different color. If an array (K, ) is given, where K is the number of points of all lines then these are considered as the values to be used by the colormap. If an array (L, ) is given, where L is the number of streamlines then these are considered as the values to be used by the colormap per streamline. If an array (X, Y, Z) or (X, Y, Z, 3) is given then the values for the colormap are interpolated automatically using trilinear interpolation.

opacity : float

Takes values from 0 (fully transparent) to 1 (opaque). Default is 1.

linewidth : float

Default is 0.01.

tube_sides : int

Default is 9.

lod : bool

Use vtkLODActor(level of detail) rather than vtkActor. Default is True. Level of detail actors do not render the full geometry when the frame rate is low.

lod_points : int

Number of points to be used when LOD is in effect. Default is 10000.

lod_points_size : int

Size of points when lod is in effect. Default is 3.

spline_subdiv : int

Number of splines subdivision to smooth streamtubes. Default is None.

lookup_colormap : vtkLookupTable

Add a default lookup table to the colormap. Default is None which calls fury.actor.colormap_lookup_table().

Notes

Streamtubes can be heavy on GPU when loading many streamlines and therefore, you may experience slow rendering time depending on system GPU. A solution to this problem is to reduce the number of points in each streamline. In Dipy we provide an algorithm that will reduce the number of points on the straighter parts of the streamline but keep more points on the curvier parts. This can be used in the following way:

from dipy.tracking.distances import approx_polygon_track
lines = [approx_polygon_track(line, 0.2) for line in lines]

Alternatively we suggest using the line actor which is much more efficient.

Examples

>>> import numpy as np
>>> from fury import actor, window
>>> ren = window.Renderer()
>>> lines = [np.random.rand(10, 3), np.random.rand(20, 3)]
>>> colors = np.random.rand(2, 3)
>>> c = actor.streamtube(lines, colors)
>>> ren.add(c)
>>> #window.show(ren)

tensor_slicer

fury.actor.tensor_slicer(evals, evecs, affine=None, mask=None, sphere=None, scale=2.2, norm=True, opacity=1.0, scalar_colors=None)[source]

Slice many tensors as ellipsoids in native or world coordinates

Parameters:
evals : (3,) or (X, 3) or (X, Y, 3) or (X, Y, Z, 3) ndarray

eigenvalues

evecs : (3, 3) or (X, 3, 3) or (X, Y, 3, 3) or (X, Y, Z, 3, 3) ndarray

eigenvectors

affine : array

4x4 transformation array from native coordinates to world coordinates*

mask : ndarray

3D mask

sphere : Sphere

a sphere

scale : float

Distance between spheres.

norm : bool

Normalize sphere_values.

opacity : float

Takes values from 0 (fully transparent) to 1 (opaque). Default is 1.

scalar_colors : (3,) or (X, 3) or (X, Y, 3) or (X, Y, Z, 3) ndarray

RGB colors used to show the tensors Default None, color the ellipsoids using color_fa

Returns:
actor : vtkActor

Ellipsoid