%extend itkPyBuffer@PyBufferTypes@{ %pythoncode %{ def GetArrayViewFromImage(image, keep_axes=False, update=True): """Get a NumPy array view of a ITK Image. When *keep_axes* is *False*, the NumPy array will have C-order indexing. This is the reverse of how indices are specified in ITK, i.e. k,j,i versus i,j,k. However C-order indexing is expected by most algorithms in NumPy / SciPy. """ if image.GetBufferPointer() is None: return None if update: # Ensure the image regions and image pixel buffer have been updated # correctly source = image.GetSource() if source: source.UpdateLargestPossibleRegion() itksize = image.GetLargestPossibleRegion().GetSize() dim = len(itksize) shape = [int(itksize[idx]) for idx in range(dim)] if(image.GetNumberOfComponentsPerPixel() > 1): shape = [image.GetNumberOfComponentsPerPixel(), ] + shape if keep_axes == False: shape.reverse() pixelType = "@PixelType@" numpy_dtype = _get_numpy_pixelid(pixelType) memview = itkPyBuffer@PyBufferTypes@._GetArrayViewFromImage(image) ndarr_view = np.asarray(memview).view(dtype = numpy_dtype).reshape(shape).view(np.ndarray) itk_view = NDArrayITKBase(ndarr_view, image) return itk_view GetArrayViewFromImage = staticmethod(GetArrayViewFromImage) def GetArrayFromImage(image, keep_axes=False, update=True): """Get a NumPy ndarray from an ITK Image. When *keep_axes* is *False*, the NumPy array will have C-order indexing. This is the reverse of how indices are specified in ITK, i.e. k,j,i versus i,j,k. However C-order indexing is expected by most algorithms in NumPy / SciPy. This is a deep copy of the image buffer and is completely safe and without potential side effects. """ if image.GetBufferPointer() is None: return None arrayView = itkPyBuffer@PyBufferTypes@.GetArrayViewFromImage(image, keep_axes, update) # perform deep copy of the image buffer arr = np.array(arrayView, copy=True) return arr GetArrayFromImage = staticmethod(GetArrayFromImage) def GetImageViewFromArray(ndarr, is_vector=False, need_contiguous=True): """Get an ITK Image view of a NumPy array. If is_vector is True, then a 3D array will be treated as a 2D vector image, otherwise it will be treated as a 3D image. If the array uses Fortran-order indexing, i.e. i,j,k, the Image Size will have the same dimensions as the array shape. If the array uses C-order indexing, i.e. k,j,i, the image Size will have the dimensions reversed from the array shape. Therefore, since the *np.transpose* operator on a 2D array simply inverts the indexing scheme, the Image representation will be the same for an array and its transpose. If flipping is desired, see *np.reshape*. By default, a warning is issued if this function is called on a non-contiguous array, since a copy is performed and care must be taken to keep a reference to the copied array. This warning can be suppressed with need_contiguous=False """ assert ndarr.ndim in (1, 2, 3, 4, 5), \ "Only arrays of 1, 2, 3, 4 or 5 dimensions are supported." if not ndarr.flags['C_CONTIGUOUS'] and not ndarr.flags['F_CONTIGUOUS']: ndarr = np.ascontiguousarray(ndarr) if need_contiguous: import warnings msg = """Because the input array was not contiguous, the returned image is not a view of the array passed to this function. Instead, it is a view of the member "base" of the returned image! If that member is ever garbage collected, this view becomes invalid.""" warnings.warn(msg) if is_vector: if ndarr.flags['C_CONTIGUOUS']: imgview = itkPyBuffer@PyBufferTypes@._GetImageViewFromArray(ndarr, ndarr.shape[-2::-1], ndarr.shape[-1]) else: imgview = itkPyBuffer@PyBufferTypes@._GetImageViewFromArray(ndarr, ndarr.shape[-1:0:-1], ndarr.shape[0]) else: imgview = itkPyBuffer@PyBufferTypes@._GetImageViewFromArray(ndarr, ndarr.shape[::-1], 1) # Keep a reference imgview._SetBase(ndarr) return imgview GetImageViewFromArray = staticmethod(GetImageViewFromArray) def GetImageFromArray(ndarr, is_vector=False): """Get an ITK Image of a NumPy array. This is a deep copy of the NumPy array buffer and is completely safe without potential side effects. If is_vector is True, then a 3D array will be treated as a 2D vector image, otherwise it will be treated as a 3D image. If the array uses Fortran-order indexing, i.e. i,j,k, the Image Size will have the same dimensions as the array shape. If the array uses C-order indexing, i.e. k,j,i, the image Size will have the dimensions reversed from the array shape. Therefore, since the *np.transpose* operator on a 2D array simply inverts the indexing scheme, the Image representation will be the same for an array and its transpose. If flipping is desired, see *np.reshape*. """ # Create a temporary image view of the array imageView = itkPyBuffer@PyBufferTypes@.GetImageViewFromArray(ndarr, is_vector, need_contiguous=False) # Duplicate the image to let it manage its own memory buffer import itk duplicator = itk.ImageDuplicator.New(imageView) duplicator.Update() return duplicator.GetOutput() GetImageFromArray = staticmethod(GetImageFromArray) %} };