Skip to content

Commit 9895eb0

Browse files
author
Jakub Kaczmarzyk
committed
add documentation and fix style
1 parent 76e9aed commit 9895eb0

File tree

1 file changed

+65
-3
lines changed

1 file changed

+65
-3
lines changed

nibabel/processing.py

Lines changed: 65 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -314,13 +314,74 @@ def smooth_image(img,
314314

315315

316316
def _transform_range(x, new_min, new_max):
317+
"""Transform data to a new range, while maintaining ratios.
318+
319+
Parameters
320+
----------
321+
x : array-like
322+
The data to transform.
323+
new_min, new_max : scalar
324+
The minimum and maximum of the output array.
325+
326+
Returns
327+
-------
328+
transformed : array-like
329+
A copy of the transformed data.
330+
331+
Examples
332+
--------
333+
>>> _transform_range([2, 4, 6], -1, 1)
334+
array([-1., 0., 1.])
335+
"""
317336
x = np.asarray(x)
318337
x_min, x_max = x.min(), x.max()
319338
return (((x - x_min) * (new_max - new_min)) / (x_max - x_min)) + new_min
320339

321340

322-
def conform(from_img, out_shape=(256, 256, 256),
323-
voxel_size=(1.0, 1.0, 1.0), order=3, cval=0.0, out_class=Nifti1Image):
341+
def conform(from_img,
342+
out_shape=(256, 256, 256),
343+
voxel_size=(1.0, 1.0, 1.0),
344+
order=3,
345+
cval=0.0,
346+
out_class=Nifti1Image):
347+
""" Resample image to ``out_shape`` with voxels of size ``voxel_size``.
348+
349+
Using the default arguments, this function is meant to replicate most parts
350+
of FreeSurfer's ``mri_convert --conform`` command. The output image is also
351+
reoriented to RAS.
352+
353+
Parameters
354+
----------
355+
from_img : object
356+
Object having attributes ``dataobj``, ``affine``, ``header`` and
357+
``shape``. If `out_class` is not None, ``img.__class__`` should be able
358+
to construct an image from data, affine and header.
359+
out_shape : sequence, optional
360+
The shape of the output volume. Default is (256, 256, 256).
361+
voxel_size : sequence, optional
362+
The size in millimeters of the voxels in the resampled output. Default
363+
is 1mm isotropic.
364+
order : int, optional
365+
The order of the spline interpolation, default is 3. The order has to
366+
be in the range 0-5 (see ``scipy.ndimage.affine_transform``)
367+
mode : str, optional
368+
Points outside the boundaries of the input are filled according
369+
to the given mode ('constant', 'nearest', 'reflect' or 'wrap').
370+
Default is 'constant' (see ``scipy.ndimage.affine_transform``)
371+
cval : scalar, optional
372+
Value used for points outside the boundaries of the input if
373+
``mode='constant'``. Default is 0.0 (see
374+
``scipy.ndimage.affine_transform``)
375+
out_class : None or SpatialImage class, optional
376+
Class of output image. If None, use ``from_img.__class__``.
377+
378+
Returns
379+
-------
380+
out_img : object
381+
Image of instance specified by `out_class`, containing data output from
382+
resampling `from_img` into axes aligned to the output space of
383+
``from_img.affine``
384+
"""
324385
# Create fake image of the image we want to resample to.
325386
hdr = Nifti1Header()
326387
hdr.set_data_shape(out_shape)
@@ -329,7 +390,8 @@ def conform(from_img, out_shape=(256, 256, 256),
329390
to_img = Nifti1Image(np.empty(out_shape), affine=dst_aff, header=hdr)
330391
# Resample input image.
331392
out_img = resample_from_to(
332-
from_img=from_img, to_vox_map=to_img, order=order, mode="constant", cval=cval, out_class=out_class)
393+
from_img=from_img, to_vox_map=to_img, order=order, mode="constant",
394+
cval=cval, out_class=out_class)
333395
# Cast to uint8.
334396
data = out_img.get_fdata()
335397
data = _transform_range(data, new_min=0.0, new_max=255.0)

0 commit comments

Comments
 (0)