io#
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 [1].
References#
| 
 | 
Module: io.gradients#
| 
 | Read b-values and b-vectors from disk. | 
Module: io.image#
| 
 | Load only the data array from a nifti file. | 
| 
 | Load data and other information from a nifti file. | 
| 
 | Save a data array into a nifti file. | 
| 
 | Save Quality Assurance metrics. | 
Module: io.peaks#
| 
 | Load a PeaksAndMetrics HDF5 file (PAM5). | 
| 
 | Load a PeaksAndMetrics HDF5 file (PAM5). | 
| 
 | Save PeaksAndMetrics object attributes in a PAM5 file (HDF5). | 
| 
 | Save all important attributes of object PeaksAndMetrics in a PAM5 file (HDF5). | 
| 
 | Save SH, directions, indices and values of peaks to Nifti. | 
| 
 | Save SH, directions, indices and values of peaks to Nifti. | 
| 
 | Return SH, directions, indices and values of peaks to pam5. | 
| 
 | Convert diffusion tensor to pam5. | 
Module: io.pickles#
Load and save pickles
| 
 | Save dix to fname as pickle. | 
| 
 | Load object from pickle file fname. | 
Module: io.stateful_tractogram#
| 
 | Enum to simplify future change to convention | 
| 
 | Enum to simplify future change to convention | 
| 
 | Class for stateful representation of collections of streamlines Object designed to be identical no matter the file format (trk, tck, vtk, fib, dpy). | 
Module: io.streamline#
| 
 | Save the stateful tractogram in any format (trk/tck/vtk/vtp/fib/dpy) | 
| 
 | Load the stateful tractogram from any format (trk/tck/vtk/vtp/fib/dpy) | 
| 
 | Generate a loading function that performs a file extension check to restrict the user to a single file format. | 
| 
 | Generate a saving function that performs a file extension check to restrict the user to a single file format. | 
| 
 | Load the stateful tractogram of the .trk format | 
| 
 | Load the stateful tractogram of the .tck format | 
| 
 | Load the stateful tractogram of the .trx format | 
| 
 | Load the stateful tractogram of the .vtk format | 
| 
 | Load the stateful tractogram of the .vtp format | 
| 
 | Load the stateful tractogram of the .fib format | 
| 
 | Load the stateful tractogram of the .dpy format | 
| 
 | Save the stateful tractogram of the .trk format | 
| 
 | Save the stateful tractogram of the .tck format | 
| 
 | Save the stateful tractogram of the .trx format | 
| 
 | Save the stateful tractogram of the .vtk format | 
| 
 | Save the stateful tractogram of the .vtp format | 
| 
 | Save the stateful tractogram of the .fib format | 
| 
 | Save the stateful tractogram of the .dpy format | 
Module: io.surface#
| 
 | Load pial file. | 
| 
 | Load gifti file. | 
Module: io.utils#
Utility functions for file formats
| 
 | Returns a Nifti1Image with a symmetric matrix intent | 
| 
 | reshapes the input to have 5 dimensions, adds extra dimensions just before the last dimension | 
| 
 | Create a nifti-compliant directional-encoded color FA image. | 
| 
 | Convert a nifti-compliant directional-encoded color FA image into a nifti image with RGB encoded in floating point resolution. | 
| 
 | Validate basic data type and value of spatial attribute. | 
| 
 | Returns the clean basename and extension of a file. | 
| 
 | Will compare the spatial attribute of 2 references. | 
| 
 | Will compare the spatial attribute of 2 references | 
| 
 | Write a standard trk/tck header from spatial attribute | 
| 
 | Write a standard nifti header from spatial attribute | 
| 
 | Saves the given input dataframe to .h5 file | 
| 
 | Helper function that handles inputs that can be paths, nifti img or arrays | 
Module: io.vtk#
| 
 | Load a vtk polydata to a supported format file. | 
| 
 | Save a vtk polydata to a supported format file. | 
| 
 | Save streamlines as vtk polydata to a supported format file. | 
| 
 | Load streamlines from vtk polydata. | 
Dpy#
- class dipy.io.dpy.Dpy(fname, *, mode='r', compression=0)[source]#
- Bases: - object- Methods - read one track each time - read the entire tractography - read_tracksi(indices)- read tracks with specific indices - write_track(track)- write on track each time - write_tracks(tracks)- write many tracks together - close - version 
read_bvals_bvecs#
- dipy.io.gradients.read_bvals_bvecs(fbvals, fbvecs)[source]#
- 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:
- bvalsarray, (N,) or None
- bvecsarray, (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)[source]#
- 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. 
 
- Returns:
- data: np.ndarray or nib.ArrayProxy
 
 - See also 
load_nifti#
- dipy.io.image.load_nifti(fname, *, return_img=False, return_voxsize=False, return_coords=False, as_ndarray=True)[source]#
- 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. 
- return_voxsize: bool, optional
- Whether to return the nifti header zooms. 
- return_coordsbool, optional
- Whether to return the nifti header aff2axcodes. 
- 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. 
 
- 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))
 
 - See also 
save_nifti#
- dipy.io.image.save_nifti(fname, data, affine, *, hdr=None, dtype=None)[source]#
- 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. 
 
- Returns:
- None
 
 
save_qa_metric#
load_peaks#
- dipy.io.peaks.load_peaks(fname, *, verbose=False)[source]#
- Load a PeaksAndMetrics HDF5 file (PAM5). - dipy.io.peaks.load_peaks is deprecated, Please usedipy.io.peaks.load_pam instead - deprecated from version: 1.10 
- Will raise <class ‘dipy.utils.deprecator.ExpiredDeprecationError’> as of version: 1.12 
 - Parameters:
- fnamestring
- Filename of PAM5 file. 
- verbosebool
- Print summary information about the loaded file. 
 
- Returns:
- pamPeaksAndMetrics object
 
 
load_pam#
save_peaks#
- dipy.io.peaks.save_peaks(fname, pam, *, affine=None, verbose=False)[source]#
- Save PeaksAndMetrics object attributes in a PAM5 file (HDF5). - dipy.io.peaks.save_peaks is deprecated, Please usedipy.io.peaks.save_pam instead - deprecated from version: 1.10.0 
- Will raise <class ‘dipy.utils.deprecator.ExpiredDeprecationError’> as of version: 1.12.0 
 - 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. 
 
 
save_pam#
- dipy.io.peaks.save_pam(fname, pam, *, affine=None, verbose=False)[source]#
- Save all important attributes of object PeaksAndMetrics in a PAM5 file (HDF5). - Parameters:
- fnamestr
- Filename of PAM5 file. 
- pamPeaksAndMetrics
- Object holding peaks information and metrics. 
- affinendarray, optional
- The 4x4 matrix transforming the date from native to world coordinates. PeaksAndMetrics should have that attribute but if not it can be provided here. 
- verbosebool, optional
- 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=None, reshape_dirs=False)[source]#
- Save SH, directions, indices and values of peaks to Nifti. - dipy.io.peaks.peaks_to_niftis is deprecated, Please use dipy.io.peaks.pam_to_niftis instead - deprecated from version: 1.10.0 
- Will raise <class ‘dipy.utils.deprecator.ExpiredDeprecationError’> as of version: 1.12.0 
 - Parameters:
- pamPeaksAndMetrics
- Object holding peaks information and metrics. 
- fname_shmstr
- Spherical Harmonics coefficients filename. 
- fname_dirsstr
- Peaks direction filename. 
- fname_valuesstr
- Peaks values filename. 
- fname_indicesstr
- Peaks indices filename. 
- fname_gfastr, optional
- Generalized FA filename. 
- reshape_dirsbool, optional
- If True, reshape peaks for visualization. 
 
 
pam_to_niftis#
- dipy.io.peaks.pam_to_niftis(pam, *, fname_peaks_dir='peaks_dirs.nii.gz', fname_peaks_values='peaks_values.nii.gz', fname_peaks_indices='peaks_indices.nii.gz', fname_shm='shm.nii.gz', fname_gfa='gfa.nii.gz', fname_sphere='sphere.txt', fname_b='B.nii.gz', fname_qa='qa.nii.gz', reshape_dirs=False)[source]#
- Save SH, directions, indices and values of peaks to Nifti. - Parameters:
- pamPeaksAndMetrics
- Object holding peaks information and metrics. 
- fname_peaks_dirstr, optional
- Peaks direction filename. 
- fname_peaks_valuesstr, optional
- Peaks values filename. 
- fname_peaks_indicesstr, optional
- Peaks indices filename. 
- fname_shmstr, optional
- Spherical Harmonics coefficients filename. It will be saved if available. 
- fname_gfastr, optional
- Generalized FA filename. It will be saved if available. 
- fname_spherestr, optional
- Sphere vertices filename. It will be saved if available. 
- fname_bstr, optional
- B Matrix filename. Matrix that transforms spherical harmonics to spherical function. It will be saved if available. 
- fname_qastr, optional
- Quantitative Anisotropy filename. It will be saved if available. 
- reshape_dirsbool, optional
- If True, Reshape and convert to float32 a set of peaks for visualisation with mrtrix or the fibernavigator. 
 
 
niftis_to_pam#
- dipy.io.peaks.niftis_to_pam(affine, peak_dirs, peak_values, peak_indices, *, shm_coeff=None, sphere=None, gfa=None, B=None, qa=None, odf=None, total_weight=None, ang_thr=None, pam_file=None)[source]#
- Return SH, directions, indices and values of peaks to pam5. - Parameters:
- affinearray, (4, 4)
- The matrix defining the affine transform. 
- peak_dirsndarray
- The direction of each peak. 
- peak_valuesndarray
- The value of the peaks. 
- peak_indicesndarray
- Indices (in sphere vertices) of the peaks in each voxel. 
- shm_coeffarray, optional
- Spherical harmonics coefficients. 
- sphereSphere class instance, optional
- The sphere providing discrete directions for evaluation. 
- gfandarray, optional
- Generalized FA volume. 
- Bndarray, optional
- Matrix that transforms spherical harmonics to spherical function. 
- qaarray, optional
- Quantitative Anisotropy in each voxel. 
- odfndarray, optional
- SH coefficients for the ODF spherical function. 
- total_weightfloat, optional
- Total weight of the peaks. 
- ang_thrfloat, optional
- Angular threshold of the peaks. 
- pam_filestr, optional
- Filename of the desired pam file. 
 
- Returns:
- pamPeaksAndMetrics
- Object holding peak_dirs, shm_coeffs and other attributes. 
 
 
tensor_to_pam#
- dipy.io.peaks.tensor_to_pam(evals, evecs, affine, *, shm_coeff=None, sphere=None, gfa=None, B=None, qa=None, odf=None, total_weight=None, ang_thr=None, pam_file=None, npeaks=5, generate_peaks_indices=True)[source]#
- Convert diffusion tensor to pam5. - Parameters:
- evalsndarray
- Eigenvalues of a diffusion tensor. shape should be (…,3). 
- evecsndarray
- Eigen vectors from the tensor model. 
- affinearray, (4, 4)
- The matrix defining the affine transform. 
- shm_coeffarray, optional
- Spherical harmonics coefficients. 
- sphereSphere class instance, optional
- The sphere providing discrete directions for evaluation. 
- gfandarray, optional
- Generalized FA volume. 
- Bndarray, optional
- Matrix that transforms spherical harmonics to spherical function. 
- qaarray, optional
- Quantitative Anisotropy in each voxel. 
- odfndarray, optional
- SH coefficients for the ODF spherical function. 
- pam_filestr, optional
- Filename of the desired pam file. 
- npeaksint, optional
- Maximum number of peaks found. 
- generate_peaks_indicesbool, optional
- total_weightfloat, optional
- Total weight of the peaks. 
- ang_thrfloat, optional
- Angular threshold of the peaks. 
 
- Returns:
- pamPeaksAndMetrics
- Object holding peaks information and metrics. 
 
 
save_pickle#
- dipy.io.pickles.save_pickle(fname, dix)[source]#
- Save dix to fname as pickle. - Parameters:
- fnamestr
- filename to save object e.g. a dictionary 
- dixstr
- dictionary or other object 
 
 - See also - 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) 
load_pickle#
Space#
Origin#
StatefulTractogram#
- class dipy.io.stateful_tractogram.StatefulTractogram(streamlines, reference, space, *, origin=Origin.NIFTI, data_per_point=None, data_per_streamline=None)[source]#
- 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. - Attributes:
- affine
- Getter for the reference affine 
- data_per_point
- Getter for data_per_point 
- data_per_streamline
- Getter for data_per_streamline 
- dimensions
- Getter for the reference dimensions 
- dtype_dict
- Getter for dtype_dict 
- origin
- Getter for origin standard 
- space
- Getter for the current space 
- space_attributes
- Getter for spatial attribute 
- streamlines
- Partially safe getter for streamlines 
- voxel_order
- Getter for the reference voxel order 
- voxel_sizes
- Getter for the reference voxel sizes 
 
 - Methods - are_compatible(sft_1, sft_2)- Compatibility verification of two StatefulTractogram to ensure space, origin, data_per_point and data_per_streamline consistency - Compute the bounding box of the streamlines in their current state - from_sft(streamlines, sft, *[, ...])- Create an instance of StatefulTractogram from another instance of StatefulTractogram. - Return a list of the data_per_point attribute names - Return a list of the data_per_streamline attribute names - Safe getter for streamlines (for slicing) - Verify that the bounding box is valid in voxel space. - remove_invalid_streamlines(*[, epsilon])- Remove streamlines with invalid coordinates from the object. - Safe function to shift streamlines so the center of voxel is the origin - 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 affine#
- Getter for the reference affine 
 - static are_compatible(sft_1, sft_2)[source]#
- Compatibility verification of two StatefulTractogram to ensure space, origin, data_per_point and data_per_streamline consistency 
 - compute_bounding_box()[source]#
- 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)[source]#
- 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 
- —–
 
 
 - is_bbox_in_vox_valid()[source]#
- 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)[source]#
- 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_origin(target_origin)[source]#
- Safe function to change streamlines to a particular origin standard False means NIFTI (center) and True means TrackVis (corner) 
 - to_space(target_space)[source]#
- Safe function to transform streamlines to a particular space using an enum and update state 
 - property voxel_order#
- Getter for the reference voxel order 
 - property voxel_sizes#
- Getter for the reference voxel sizes 
 
save_tractogram#
- dipy.io.streamline.save_tractogram(sft, filename, *, bbox_valid_check=True)[source]#
- 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)[source]#
- 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)[source]#
- 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)[source]#
- 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_trx#
- dipy.io.streamline.load_trx(filename, reference, *, to_space=Space.RASMM, to_origin=Origin.NIFTI, bbox_valid_check=True, trk_header_check=True)#
- Load the stateful tractogram of the .trx 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_trx#
- dipy.io.streamline.save_trx(sft, filename, bbox_valid_check=True)#
- Save the stateful tractogram of the .trx 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 
 
 
load_pial#
- dipy.io.surface.load_pial(fname, *, return_meta=False)[source]#
- Load pial file. - Parameters:
- fnamestr
- Absolute path of the file. 
- return_metabool, optional
- Whether to read the metadata of the file or not, by default False. 
 
- Returns:
- tuple
- (vertices, faces) if return_meta=False. Otherwise, (vertices, faces, metadata). 
 
 
load_gifti#
nifti1_symmat#
- dipy.io.utils.nifti1_symmat(image_data, *args, **kwargs)[source]#
- 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 
- *args
- Passed to Nifti1Image 
- **kwargs
- Passed to Nifti1Image 
 
- Returns:
- imageNifti1Image
- 5d, extra dimensions added before the last. Has symmetric matrix intent code 
 
 
make5d#
decfa#
- dipy.io.utils.decfa(img_orig, *, scale=False)[source]#
- 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. 
 
 - Notes - For a description of this format, see: - https://nifti.nimh.nih.gov/nifti-1/documentation/nifti1fields/nifti1fields_pages/datatype.html 
decfa_to_float#
- dipy.io.utils.decfa_to_float(img_orig)[source]#
- 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:
- imgNifti1Image class instance with float dtype.
 
 - Notes - For a description of this format, see: - https://nifti.nimh.nih.gov/nifti-1/documentation/nifti1fields/nifti1fields_pages/datatype.html 
is_reference_info_valid#
- dipy.io.utils.is_reference_info_valid(affine, dimensions, voxel_sizes, voxel_order)[source]#
- 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 verifies 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 
 
 
split_name_with_gz#
get_reference_info#
- dipy.io.utils.get_reference_info(reference)[source]#
- Will compare the spatial attribute of 2 references. - Parameters:
- referenceNifti or Trk filename, Nifti1Image or TrkFile, Nifti1Header or
- trk.header (dict), TrxFile or trx.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’ 
 
 
 
is_header_compatible#
- dipy.io.utils.is_header_compatible(reference_1, reference_2)[source]#
- Will compare the spatial attribute of 2 references - Parameters:
- reference_1Nifti or Trk filename, Nifti1Image or TrkFile,
- Nifti1Header or trk.header (dict) Reference that provides the spatial attribute. 
- reference_2Nifti or Trk filename, Nifti1Image or TrkFile,
- Nifti1Header or trk.header (dict) Reference that provides the spatial attribute. 
 
- Returns:
- outputbool
- Does all the spatial attribute match 
 
 
create_tractogram_header#
create_nifti_header#
save_buan_profiles_hdf5#
- dipy.io.utils.save_buan_profiles_hdf5(fname, dt, *, key=None)[source]#
- 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 
- keystr, optional
- Key to retrieve the contents in the HDF5 file. The file rootname will be used if not provided. 
 
 
read_img_arr_or_path#
- dipy.io.utils.read_img_arr_or_path(data, *, affine=None)[source]#
- 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, affinendarray and 4x4 array
 
 
load_polydata#
save_polydata#
save_vtk_streamlines#
- dipy.io.vtk.save_vtk_streamlines(streamlines, filename, *, to_lps=True, binary=False)[source]#
- 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)[source]#
- 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