Module: io.bvectxt
read_bvec_file (filename[, atol])
|
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead |
ornt_mapping (ornt1, ornt2)
|
Calculate the mapping needing to get from orn1 to orn2. |
reorient_vectors (bvecs, current_ornt, new_ornt)
|
Change the orientation of gradients or other vectors. |
reorient_on_axis (bvecs, current_ornt, new_ornt)
|
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead |
orientation_from_string (string_ornt)
|
Return an array representation of an ornt string. |
orientation_to_string (ornt)
|
Return a string representation of a 3d ornt. |
Module: io.dpy
A class for handling large tractography datasets.
It is built using the h5py which in turn implement
key features of the HDF5 (hierarchical data format) API .
Dpy (fname[, mode, compression])
|
|
Module: io.image
load_nifti_data (fname[, as_ndarray])
|
Load only the data array from a nifti file. |
load_nifti (fname[, return_img, ...])
|
Load data and other information from a nifti file. |
save_nifti (fname, data, affine[, hdr, dtype])
|
Save a data array into a nifti file. |
save_qa_metric (fname, xopt, fopt)
|
Save Quality Assurance metrics. |
Module: io.peaks
load_peaks (fname[, verbose])
|
Load a PeaksAndMetrics HDF5 file (PAM5) |
save_peaks (fname, pam[, affine, verbose])
|
Save all important attributes of object PeaksAndMetrics in a PAM5 file (HDF5). |
peaks_to_niftis (pam, fname_shm, fname_dirs, ...)
|
Save SH, directions, indices and values of peaks to Nifti. |
Module: io.pickles
Load and save pickles
Module: io.stateful_tractogram
Space (value)
|
Enum to simplify future change to convention |
Origin (value)
|
Enum to simplify future change to convention |
StatefulTractogram (streamlines, reference, space)
|
Class for stateful representation of collections of streamlines Object designed to be identical no matter the file format (trk, tck, vtk, fib, dpy). |
logger
|
Instances of the Logger class represent a single logging channel. |
set_sft_logger_level (log_level)
|
Change the logger of the StatefulTractogram to one on the following: DEBUG, INFO, WARNING, CRITICAL, ERROR |
Module: io.streamline
save_tractogram (sft, filename[, ...])
|
Save the stateful tractogram in any format (trk/tck/vtk/vtp/fib/dpy) |
load_tractogram (filename, reference[, ...])
|
Load the stateful tractogram from any format (trk/tck/vtk/vtp/fib/dpy) |
load_generator (ttype)
|
Generate a loading function that performs a file extension check to restrict the user to a single file format. |
save_generator (ttype)
|
Generate a saving function that performs a file extension check to restrict the user to a single file format. |
load_trk (filename, reference[, to_space, ...])
|
Load the stateful tractogram of the .trk format |
load_tck (filename, reference[, to_space, ...])
|
Load the stateful tractogram of the .tck format |
load_vtk (filename, reference[, to_space, ...])
|
Load the stateful tractogram of the .vtk format |
load_vtp (filename, reference[, to_space, ...])
|
Load the stateful tractogram of the .vtp format |
load_fib (filename, reference[, to_space, ...])
|
Load the stateful tractogram of the .fib format |
load_dpy (filename, reference[, to_space, ...])
|
Load the stateful tractogram of the .dpy format |
save_trk (sft, filename[, bbox_valid_check])
|
Save the stateful tractogram of the .trk format |
save_tck (sft, filename[, bbox_valid_check])
|
Save the stateful tractogram of the .tck format |
save_vtk (sft, filename[, bbox_valid_check])
|
Save the stateful tractogram of the .vtk format |
save_vtp (sft, filename[, bbox_valid_check])
|
Save the stateful tractogram of the .vtp format |
save_fib (sft, filename[, bbox_valid_check])
|
Save the stateful tractogram of the .fib format |
save_dpy (sft, filename[, bbox_valid_check])
|
Save the stateful tractogram of the .dpy format |
Module: io.utils
Utility functions for file formats
nifti1_symmat (image_data, *args, **kwargs)
|
Returns a Nifti1Image with a symmetric matrix intent |
make5d (data)
|
reshapes the input to have 5 dimensions, adds extra dimensions just before the last dimension |
decfa (img_orig[, scale])
|
Create a nifti-compliant directional-encoded color FA image. |
decfa_to_float (img_orig)
|
Convert a nifti-compliant directional-encoded color FA image into a nifti image with RGB encoded in floating point resolution. |
is_reference_info_valid (affine, dimensions, ...)
|
Validate basic data type and value of spatial attribute. |
get_reference_info (reference)
|
Will compare the spatial attribute of 2 references |
is_header_compatible (reference_1, reference_2)
|
Will compare the spatial attribute of 2 references |
create_tractogram_header (tractogram_type, ...)
|
Write a standard trk/tck header from spatial attribute |
create_nifti_header (affine, dimensions, ...)
|
Write a standard nifti header from spatial attribute |
save_buan_profiles_hdf5 (fname, dt)
|
Saves the given input dataframe to .h5 file |
read_img_arr_or_path (data[, affine])
|
Helper function that handles inputs that can be paths, nifti img or arrays |
Module: io.vtk
load_polydata (file_name)
|
Load a vtk polydata to a supported format file. |
save_polydata (polydata, file_name[, binary, ...])
|
Save a vtk polydata to a supported format file. |
save_vtk_streamlines (streamlines, filename)
|
Save streamlines as vtk polydata to a supported format file. |
load_vtk_streamlines (filename[, to_lps])
|
Load streamlines from vtk polydata. |
read_bvec_file
-
dipy.io.bvectxt.read_bvec_file(filename, atol=0.001)
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead
Read gradient table information from a pair of files with extensions
.bvec and .bval. The bval file should have one row of values
representing the bvalues of each volume in the dwi data set. The bvec
file should have three rows, where the rows are the x, y, and z
components of the normalized gradient direction for each of the
volumes.
Parameters
- filename :
The path to the either the bvec or bval file
- atolfloat, optional
The tolerance used to check all the gradient directions are
normalized. Default is .001
ornt_mapping
-
dipy.io.bvectxt.ornt_mapping(ornt1, ornt2)
Calculate the mapping needing to get from orn1 to orn2.
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead
reorient_vectors
-
dipy.io.bvectxt.reorient_vectors(bvecs, current_ornt, new_ornt, axis=0)
Change the orientation of gradients or other vectors.
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead
Moves vectors, storted along axis, from current_ornt to new_ornt. For
example the vector [x, y, z] in “RAS” will be [-x, -y, z] in “LPS”.
R: Right
A: Anterior
S: Superior
L: Left
P: Posterior
I: Inferior
reorient_on_axis
-
dipy.io.bvectxt.reorient_on_axis(bvecs, current_ornt, new_ornt, axis=0)
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead
orientation_from_string
-
dipy.io.bvectxt.orientation_from_string(string_ornt)
Return an array representation of an ornt string.
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead
orientation_to_string
-
dipy.io.bvectxt.orientation_to_string(ornt)
Return a string representation of a 3d ornt.
dipy.io.bvectxt module is deprecated, Please use dipy.core.gradients module instead
-
class dipy.io.dpy.Dpy(fname, mode='r', compression=0)
Bases: object
-
__init__(fname, mode='r', compression=0)
Advanced storage system for tractography based on HDF5
Parameters
fname : str, full filename
mode : ‘r’ read
‘w’ write
‘r+’ read and write only if file already exists
compression : 0 no compression to 9 maximum compression
Examples
>>> import os
>>> from tempfile import mkstemp #temp file
>>> from dipy.io.dpy import Dpy
>>> def dpy_example():
... fd,fname = mkstemp()
... fname += '.dpy'#add correct extension
... dpw = Dpy(fname,'w')
... A=np.ones((5,3))
... B=2*A.copy()
... C=3*A.copy()
... dpw.write_track(A)
... dpw.write_track(B)
... dpw.write_track(C)
... dpw.close()
... dpr = Dpy(fname,'r')
... dpr.read_track()
... dpr.read_track()
... dpr.read_tracksi([0, 1, 2, 0, 0, 2])
... dpr.close()
... os.remove(fname) #delete file from disk
>>> dpy_example()
-
close()
-
read_track()
read one track each time
-
read_tracks()
read the entire tractography
-
read_tracksi(indices)
read tracks with specific indices
-
version()
-
write_track(track)
write on track each time
-
write_tracks(tracks)
write many tracks together
read_bvals_bvecs
-
dipy.io.gradients.read_bvals_bvecs(fbvals, fbvecs)
Read b-values and b-vectors from disk.
Parameters
- fbvalsstr
Full path to file with b-values. None to not read bvals.
- fbvecsstr
Full path of file with b-vectors. None to not read bvecs.
Returns
bvals : array, (N,) or None
bvecs : array, (N, 3) or None
Notes
Files can be either ‘.bvals’/’.bvecs’ or ‘.txt’ or ‘.npy’ (containing
arrays stored with the appropriate values).
load_nifti_data
-
dipy.io.image.load_nifti_data(fname, as_ndarray=True)
Load only the data array from a nifti file.
Parameters
- fnamestr
Full path to the file.
- as_ndarray: bool, optional
convert nibabel ArrayProxy to a numpy.ndarray.
If you want to save memory and delay this casting, just turn this
option to False (default: True)
Returns
data: np.ndarray or nib.ArrayProxy
load_nifti
-
dipy.io.image.load_nifti(fname, return_img=False, return_voxsize=False, return_coords=False, as_ndarray=True)
Load data and other information from a nifti file.
Parameters
- fnamestr
Full path to a nifti file.
- return_imgbool, optional
Whether to return the nibabel nifti img object. Default: False
- return_voxsize: bool, optional
Whether to return the nifti header zooms. Default: False
- return_coordsbool, optional
Whether to return the nifti header aff2axcodes. Default: False
- as_ndarray: bool, optional
convert nibabel ArrayProxy to a numpy.ndarray.
If you want to save memory and delay this casting, just turn this
option to False (default: True)
Returns
A tuple, with (at the most, if all keyword args are set to True):
(data, img.affine, img, vox_size, nib.aff2axcodes(img.affine))
save_nifti
-
dipy.io.image.save_nifti(fname, data, affine, hdr=None, dtype=None)
Save a data array into a nifti file.
Parameters
- fnamestr
The full path to the file to be saved.
- datandarray
The array with the data to save.
- affine4x4 array
The affine transform associated with the file.
- hdrnifti header, optional
May contain additional information to store in the file header.
save_qa_metric
-
dipy.io.image.save_qa_metric(fname, xopt, fopt)
Save Quality Assurance metrics.
Parameters
- fname: string
File name to save the metric values.
- xopt: numpy array
The metric containing the
optimal parameters for
image registration.
- fopt: int
The distance between the registered images.
load_peaks
-
dipy.io.peaks.load_peaks(fname, verbose=False)
Load a PeaksAndMetrics HDF5 file (PAM5)
Parameters
- fnamestring
Filename of PAM5 file.
- verbosebool
Print summary information about the loaded file.
Returns
pam : PeaksAndMetrics object
save_peaks
-
dipy.io.peaks.save_peaks(fname, pam, affine=None, verbose=False)
Save all important attributes of object PeaksAndMetrics in a PAM5 file
(HDF5).
Parameters
- fnamestring
Filename of PAM5 file
- pamPeaksAndMetrics
Object holding peak_dirs, shm_coeffs and other attributes
- affinearray
The 4x4 matrix transforming the date from native to world coordinates.
PeaksAndMetrics should have that attribute but if not it can be
provided here. Default None.
- verbosebool
Print summary information about the saved file.
peaks_to_niftis
-
dipy.io.peaks.peaks_to_niftis(pam, fname_shm, fname_dirs, fname_values, fname_indices, fname_gfa, reshape_dirs=False)
Save SH, directions, indices and values of peaks to Nifti.
save_pickle
-
dipy.io.pickles.save_pickle(fname, dix)
Save dix to fname as pickle.
Parameters
- fnamestr
filename to save object e.g. a dictionary
- dixstr
dictionary or other object
Examples
>>> import os
>>> from tempfile import mkstemp
>>> fd, fname = mkstemp() # make temporary file (opened, attached to fh)
>>> d={0:{'d':1}}
>>> save_pickle(fname, d)
>>> d2=load_pickle(fname)
We remove the temporary file we created for neatness
>>> os.close(fd) # the file is still open, we need to close the fh
>>> os.remove(fname)
See Also
dipy.io.pickles.load_pickle
load_pickle
-
dipy.io.pickles.load_pickle(fname)
Load object from pickle file fname.
Parameters
- fnamestr
filename to load dict or other python object
Returns
- dixobject
dictionary or other object
Examples
dipy.io.pickles.save_pickle
-
class dipy.io.stateful_tractogram.Space(value)
Bases: Enum
Enum to simplify future change to convention
-
__init__()
-
RASMM = 'rasmm'
-
VOX = 'vox'
-
VOXMM = 'voxmm'
-
class dipy.io.stateful_tractogram.Origin(value)
Bases: Enum
Enum to simplify future change to convention
-
__init__()
-
NIFTI = 'center'
-
TRACKVIS = 'corner'
-
class dipy.io.stateful_tractogram.StatefulTractogram(streamlines, reference, space, origin=Origin.NIFTI, data_per_point=None, data_per_streamline=None)
Bases: object
Class for stateful representation of collections of streamlines
Object designed to be identical no matter the file format
(trk, tck, vtk, fib, dpy). Facilitate transformation between space and
data manipulation for each streamline / point.
-
__init__(streamlines, reference, space, origin=Origin.NIFTI, data_per_point=None, data_per_streamline=None)
Create a strict, state-aware, robust tractogram
Parameters
- streamlineslist or ArraySequence
Streamlines of the tractogram
- referenceNifti or Trk filename, Nifti1Image or TrkFile,
Nifti1Header, trk.header (dict) or another Stateful Tractogram
Reference that provides the spatial attributes.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- spaceEnum (dipy.io.stateful_tractogram.Space)
Current space in which the streamlines are (vox, voxmm or rasmm)
After tracking the space is VOX, after loading with nibabel
the space is RASMM
- originEnum (dipy.io.stateful_tractogram.Origin), optional
Current origin in which the streamlines are (center or corner)
After loading with nibabel the origin is CENTER
- data_per_pointdict, optional
Dictionary in which each key has X items, each items has Y_i items
X being the number of streamlines
Y_i being the number of points on streamlines #i
- data_per_streamlinedict, optional
Dictionary in which each key has X items
X being the number of streamlines
Notes
Very important to respect the convention, verify that streamlines
match the reference and are effectively in the right space.
Any change to the number of streamlines, data_per_point or
data_per_streamline requires particular verification.
In a case of manipulation not allowed by this object, use Nibabel
directly and be careful.
-
property affine
Getter for the reference affine
-
static are_compatible(sft_1, sft_2)
Compatibility verification of two StatefulTractogram to ensure space,
origin, data_per_point and data_per_streamline consistency
-
compute_bounding_box()
Compute the bounding box of the streamlines in their current state
Returns
- outputndarray
8 corners of the XYZ aligned box, all zeros if no streamlines
-
property data_per_point
Getter for data_per_point
-
property data_per_streamline
Getter for data_per_streamline
-
property dimensions
Getter for the reference dimensions
-
property dtype_dict
Getter for dtype_dict
-
static from_sft(streamlines, sft, data_per_point=None, data_per_streamline=None)
Create an instance of StatefulTractogram from another instance
of StatefulTractogram.
Parameters
- streamlineslist or ArraySequence
Streamlines of the tractogram
- sftStatefulTractogram,
The other StatefulTractogram to copy the space_attribute AND
state from.
- data_per_pointdict, optional
Dictionary in which each key has X items, each items has Y_i items
X being the number of streamlines
Y_i being the number of points on streamlines #i
- data_per_streamlinedict, optional
Dictionary in which each key has X items
X being the number of streamlines
-
get_data_per_point_keys()
Return a list of the data_per_point attribute names
-
get_data_per_streamline_keys()
Return a list of the data_per_streamline attribute names
-
get_streamlines_copy()
Safe getter for streamlines (for slicing)
-
is_bbox_in_vox_valid()
Verify that the bounding box is valid in voxel space.
Negative coordinates or coordinates above the volume dimensions
are considered invalid in voxel space.
Returns
- outputbool
Are the streamlines within the volume of the associated reference
-
property origin
Getter for origin standard
-
remove_invalid_streamlines(epsilon=0.001)
Remove streamlines with invalid coordinates from the object.
Will also remove the data_per_point and data_per_streamline.
Invalid coordinates are any X,Y,Z values above the reference
dimensions or below zero
Parameters
- epsilonfloat (optional)
Epsilon value for the bounding box verification.
Default is 1e-6.
Returns
- outputtuple
Tuple of two list, indices_to_remove, indices_to_keep
-
property space
Getter for the current space
-
property space_attributes
Getter for spatial attribute
-
property streamlines
Partially safe getter for streamlines
-
to_center()
Safe function to shift streamlines so the center of voxel is
the origin
-
to_corner()
Safe function to shift streamlines so the corner of voxel is
the origin
-
to_origin(target_origin)
Safe function to change streamlines to a particular origin standard
False means NIFTI (center) and True means TrackVis (corner)
-
to_rasmm()
Safe function to transform streamlines and update state
-
to_space(target_space)
Safe function to transform streamlines to a particular space using
an enum and update state
-
to_vox()
Safe function to transform streamlines and update state
-
to_voxmm()
Safe function to transform streamlines and update state
-
property voxel_order
Getter for the reference voxel order
-
property voxel_sizes
Getter for the reference voxel sizes
logger
-
dipy.io.stateful_tractogram.logger()
Instances of the Logger class represent a single logging channel. A
“logging channel” indicates an area of an application. Exactly how an
“area” is defined is up to the application developer. Since an
application can have any number of areas, logging channels are identified
by a unique string. Application areas can be nested (e.g. an area
of “input processing” might include sub-areas “read CSV files”, “read
XLS files” and “read Gnumeric files”). To cater for this natural nesting,
channel names are organized into a namespace hierarchy where levels are
separated by periods, much like the Java or Python package namespace. So
in the instance given above, channel names might be “input” for the upper
level, and “input.csv”, “input.xls” and “input.gnu” for the sub-levels.
There is no arbitrary limit to the depth of nesting.
set_sft_logger_level
-
dipy.io.stateful_tractogram.set_sft_logger_level(log_level)
Change the logger of the StatefulTractogram
to one on the following: DEBUG, INFO, WARNING, CRITICAL, ERROR
Parameters
- log_levelstr
Log level for the StatefulTractogram only
save_tractogram
-
dipy.io.streamline.save_tractogram(sft, filename, bbox_valid_check=True)
Save the stateful tractogram in any format (trk/tck/vtk/vtp/fib/dpy)
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
load_tractogram
-
dipy.io.streamline.load_tractogram(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram from any format (trk/tck/vtk/vtp/fib/dpy)
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
load_generator
-
dipy.io.streamline.load_generator(ttype)
Generate a loading function that performs a file extension
check to restrict the user to a single file format.
Parameters
- ttypestring
Extension of the file format that requires a loader
Returns
- outputfunction
Function (load_tractogram) that handle only one file format
save_generator
-
dipy.io.streamline.save_generator(ttype)
Generate a saving function that performs a file extension
check to restrict the user to a single file format.
Parameters
- ttypestring
Extension of the file format that requires a saver
Returns
- outputfunction
Function (save_tractogram) that handle only one file format
load_trk
-
dipy.io.streamline.load_trk(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram of the .trk format
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
load_tck
-
dipy.io.streamline.load_tck(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram of the .tck format
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
load_vtk
-
dipy.io.streamline.load_vtk(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram of the .vtk format
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
load_vtp
-
dipy.io.streamline.load_vtp(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram of the .vtp format
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
load_fib
-
dipy.io.streamline.load_fib(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram of the .fib format
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
load_dpy
-
dipy.io.streamline.load_dpy(filename, reference, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)
Load the stateful tractogram of the .dpy format
Parameters
- filenamestring
Filename with valid extension
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict), or ‘same’ if the input is a trk file.
Reference that provides the spatial attribute.
Typically a nifti-related object from the native diffusion used for
streamlines generation
- to_spaceEnum (dipy.io.stateful_tractogram.Space)
Space to which the streamlines will be transformed after loading
- to_originEnum (dipy.io.stateful_tractogram.Origin)
- Origin to which the streamlines will be transformed after loading
NIFTI standard, default (center of the voxel)
TRACKVIS standard (corner of the voxel)
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
- trk_header_checkbool
Verification that the reference has the same header as the spatial
attributes as the input tractogram when a Trk is loaded
Returns
- outputStatefulTractogram
The tractogram to load (must have been saved properly)
save_trk
-
dipy.io.streamline.save_trk(sft, filename, bbox_valid_check=True)
Save the stateful tractogram of the .trk format
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
save_tck
-
dipy.io.streamline.save_tck(sft, filename, bbox_valid_check=True)
Save the stateful tractogram of the .tck format
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
save_vtk
-
dipy.io.streamline.save_vtk(sft, filename, bbox_valid_check=True)
Save the stateful tractogram of the .vtk format
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
save_vtp
-
dipy.io.streamline.save_vtp(sft, filename, bbox_valid_check=True)
Save the stateful tractogram of the .vtp format
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
save_fib
-
dipy.io.streamline.save_fib(sft, filename, bbox_valid_check=True)
Save the stateful tractogram of the .fib format
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
save_dpy
-
dipy.io.streamline.save_dpy(sft, filename, bbox_valid_check=True)
Save the stateful tractogram of the .dpy format
Parameters
- sftStatefulTractogram
The stateful tractogram to save
- filenamestring
Filename with valid extension
- bbox_valid_checkbool
Verification for negative voxel coordinates or values above the
volume dimensions. Default is True, to enforce valid file.
Returns
- outputbool
True if the saving operation was successful
nifti1_symmat
-
dipy.io.utils.nifti1_symmat(image_data, *args, **kwargs)
Returns a Nifti1Image with a symmetric matrix intent
Parameters
- image_dataarray-like
should have lower triangular elements of a symmetric matrix along the
last dimension
all other arguments and keywords are passed to Nifti1Image
Returns
- imageNifti1Image
5d, extra dimensions added before the last. Has symmetric matrix intent
code
make5d
-
dipy.io.utils.make5d(data)
reshapes the input to have 5 dimensions, adds extra dimensions just
before the last dimension
decfa
-
dipy.io.utils.decfa(img_orig, scale=False)
Create a nifti-compliant directional-encoded color FA image.
Parameters
- img_origNifti1Image class instance.
Contains encoding of the DEC FA image with a 4D volume of data, where
the elements on the last dimension represent R, G and B components.
- scale: bool.
Whether to scale the incoming data from the 0-1 to the 0-255 range
expected in the output.
Returns
- imgNifti1Image class instance with dtype set to store tuples of
uint8 in (R, G, B) order.
decfa_to_float
-
dipy.io.utils.decfa_to_float(img_orig)
Convert a nifti-compliant directional-encoded color FA image into a
nifti image with RGB encoded in floating point resolution.
Parameters
- img_origNifti1Image class instance.
Contains encoding of the DEC FA image with a 3D volume of data, where
each element is a (R, G, B) tuple in uint8.
Returns
img : Nifti1Image class instance with float dtype.
is_reference_info_valid
-
dipy.io.utils.is_reference_info_valid(affine, dimensions, voxel_sizes, voxel_order)
Validate basic data type and value of spatial attribute.
Does not ensure that voxel_sizes and voxel_order are self-coherent with
the affine.
Only verify the following:
affine is of the right type (float) and dimension (4,4)
affine contain values in the rotation part
dimensions is of right type (int) and length (3)
voxel_sizes is of right type (float) and length (3)
voxel_order is of right type (str) and length (3)
The listed parameters are what is expected, provide something else and this
function should fail (cover common mistakes).
Parameters
- affine: ndarray (4,4)
Transformation of VOX to RASMM
- dimensions: ndarray (3,), int16
Volume shape for each axis
- voxel_sizes: ndarray (3,), float32
Size of voxel for each axis
- voxel_order: string
Typically ‘RAS’ or ‘LPS’
Returns
- outputbool
Does the input represent a valid ‘state’ of spatial attribute
get_reference_info
-
dipy.io.utils.get_reference_info(reference)
Will compare the spatial attribute of 2 references
Parameters
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
trk.header (dict)
Reference that provides the spatial attribute.
Returns
- outputtuple
affine ndarray (4,4), np.float32, transformation of VOX to RASMM
dimensions ndarray (3,), int16, volume shape for each axis
voxel_sizes ndarray (3,), float32, size of voxel for each axis
voxel_order, string, Typically ‘RAS’ or ‘LPS’
save_buan_profiles_hdf5
-
dipy.io.utils.save_buan_profiles_hdf5(fname, dt)
Saves the given input dataframe to .h5 file
Parameters
- fnamestring
file name for saving the hdf5 file
- dtPandas DataFrame
DataFrame to be saved as .h5 file
read_img_arr_or_path
-
dipy.io.utils.read_img_arr_or_path(data, affine=None)
Helper function that handles inputs that can be paths, nifti img or arrays
Parameters
- dataarray or nib.Nifti1Image or str.
Either as a 3D/4D array or as a nifti image object, or as
a string containing the full path to a nifti file.
- affine4x4 array, optional.
Must be provided for data provided as an array. If provided together
with Nifti1Image or str data, this input will over-ride the affine
that is stored in the data input. Default: use the affine stored
in data.
Returns
data, affine : ndarray and 4x4 array
load_polydata
-
dipy.io.vtk.load_polydata(file_name)
Load a vtk polydata to a supported format file.
Supported file formats are OBJ, VTK, VTP, FIB, PLY, STL and XML
Parameters
file_name : string
Returns
output : vtkPolyData
save_polydata
-
dipy.io.vtk.save_polydata(polydata, file_name, binary=False, color_array_name=None)
Save a vtk polydata to a supported format file.
Save formats can be VTK, VTP, FIB, PLY, STL and XML.
Parameters
polydata : vtkPolyData
file_name : string
save_vtk_streamlines
-
dipy.io.vtk.save_vtk_streamlines(streamlines, filename, to_lps=True, binary=False)
Save streamlines as vtk polydata to a supported format file.
File formats can be OBJ, VTK, VTP, FIB, PLY, STL and XML
Parameters
- streamlineslist
list of 2D arrays or ArraySequence
- filenamestring
output filename (.obj, .vtk, .fib, .ply, .stl and .xml)
- to_lpsbool
Default to True, will follow the vtk file convention for streamlines
Will be supported by MITKDiffusion and MI-Brain
- binarybool
save the file as binary
load_vtk_streamlines
-
dipy.io.vtk.load_vtk_streamlines(filename, to_lps=True)
Load streamlines from vtk polydata.
Load formats can be VTK, FIB
Parameters
- filenamestring
input filename (.vtk or .fib)
- to_lpsbool
Default to True, will follow the vtk file convention for streamlines
Will be supported by MITKDiffusion and MI-Brain
Returns
- outputlist
list of 2D arrays