SMILX
1.01
|
Members for atomic (one image at a time) operations. More...
Functions | |
static vtkSmartPointer< vtkImageData > | milx::Image< TImage >::ConvertITKImageToVTKImage (itk::SmartPointer< TImage > img) |
Converts a ITK image object to an VTK image object. You MUST DeepCopy the result as the ITK smartpointer is not aware of the VTK smartpointer. | |
static vtkSmartPointer< vtkImageData > | milx::Image< TImage >::ConvertITKVectorImageToVTKImage (itk::SmartPointer< TImage > img) |
Converts a ITK Vector image object to an VTK image object. You MUST DeepCopy the result as the ITK smartpointer is not aware of the VTK smartpointer. | |
template<typename TRefImage , class TPrecision > | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ApplyOrientationToITKImage (itk::SmartPointer< TImage > img, itk::SmartPointer< TRefImage > refImage, const bool labelledImage, const bool flipY=true, const bool ignoreDirection=false) |
static vtkSmartPointer< vtkImageData > | milx::Image< TImage >::ApplyOrientationToVTKImage (vtkSmartPointer< vtkImageData > img, itk::SmartPointer< TImage > refImage, vtkSmartPointer< vtkMatrix4x4 > &transformMatrix, const bool labelledImage, const bool flipY=true) |
Applies orientation/direction and origin to a VTK image from a reference image. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ConvertVTKImageToITKImage (vtkSmartPointer< vtkImageData > img) |
Converts a VTK image object to an ITK image object. | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::CastImage (itk::SmartPointer< TImage > img) |
Casts an image from one type to another. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::DuplicateImage (itk::SmartPointer< TImage > img) |
Duplicates the image into a new image. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::DuplicateImage (const TImage *img) |
Duplicates the image into a new image. | |
template<typename TVector > | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ImportVectorToImage (vnl_vector< TVector > &vec, typename TImage::SizeType size, itk::SmartPointer< TImage > image=NULL) |
Imports a VNL vector to an ITK image object. More... | |
template<typename TMatrix > | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ImportMatrixToImage (vnl_matrix< TMatrix > &matrix, itk::SmartPointer< TImage > image=NULL) |
Imports a VNL matrix to an ITK image object. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::BlankImage (typename TImage::PixelType value, typename TImage::SizeType imgSize) |
Creates an empty image filled with value given. | |
template<typename TSubImage > | |
static itk::SmartPointer< TSubImage > | milx::Image< TImage >::ExtractSubImage (itk::SmartPointer< TImage > img, typename TImage::RegionType imgRegion) |
Extract a sub image (such as a slice) from the current image given the region. | |
template<typename TImageSlice > | |
static itk::SmartPointer< TImageSlice > | milx::Image< TImage >::ExtractSlice (itk::SmartPointer< TImage > img, int *extent) |
Extract a sub image (such as a slice) from the current image given the extent (typically obtained from a VTK image actor etc.). More... | |
template<typename TImageComponent > | |
static itk::SmartPointer< TImageComponent > | milx::Image< TImage >::ExtractComponent (itk::SmartPointer< TImage > img, int component) |
Extract a component from the current vector image given the component index. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ResizeImage (itk::SmartPointer< TImage > img, typename TImage::SizeType imgSize) |
Resizes current image using current spacing. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ResizeImage (itk::SmartPointer< TImage > img, typename TImage::SizeType imgSize, typename TImage::SpacingType outputSpacing) |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ResizeImage (itk::SmartPointer< TImage > img, typename TImage::SizeType imgSize, typename TImage::SpacingType outputSpacing, typename TImage::PointType outputOrigin, typename TImage::DirectionType outputDirection) |
Resizes current image using given spacing, origin and direction. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::SubsampleImage (itk::SmartPointer< TImage > img, typename TImage::SizeType factors) |
Subsamples or shrinks the current image using the factors provided. More... | |
template<typename TOutImage , typename TTransform , typename TPrecision > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::TransformImage (itk::SmartPointer< TImage > img, itk::SmartPointer< TOutImage > refImg, itk::SmartPointer< TTransform > transf, const bool inverse, const int interp=1) |
Transform the image into a new reference image space given the transform. More... | |
template<typename TOutImage , typename TTransform , typename TPrecision > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::TransformImage (itk::SmartPointer< TImage > img, itk::SmartPointer< TTransform > transf, const bool inverse, const int interp=1) |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::ResampleImage (itk::SmartPointer< TImage > img, itk::SmartPointer< TOutImage > refImg, const bool linearInterp=false) |
Resample the image into a new reference image space given. | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::ResampleLabel (itk::SmartPointer< TImage > img, itk::SmartPointer< TOutImage > refImg) |
Resample the a labelled image into a new reference image space given ensuring no interpolator artefacts. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::AddImages (itk::SmartPointer< TImage > img1, itk::SmartPointer< TImage > img2) |
Adds image 1 to image 2, returning the result. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::SubtractImages (itk::SmartPointer< TImage > img1, itk::SmartPointer< TImage > img2) |
Subtracts image 2 from image 1, returning the result. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::DifferenceImages (itk::SmartPointer< TImage > img1, itk::SmartPointer< TImage > img2) |
Same as SubtractImages(). | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::MultiplyImages (itk::SmartPointer< TImage > img1, itk::SmartPointer< TImage > img2) |
Multiplies (element-wise) image 2 from image 1, returning the result. | |
template<class TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::ScaleImage (itk::SmartPointer< TImage > img, float scaling) |
Scales the image intensities by scaling factor and returns the result. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ScaleVectorImage (itk::SmartPointer< TImage > img, float scaling, int numberOfComponents) |
Scales each component of the vector image intensities by scaling factor and returns the result. | |
static void | milx::Image< TImage >::Information (itk::SmartPointer< TImage > img) |
Prints the information of the image to standard output. | |
static double | milx::Image< TImage >::ImageMaximum (itk::SmartPointer< TImage > img) |
Returns the maximum pixel value of an image. | |
static double | milx::Image< TImage >::ImageMinimum (itk::SmartPointer< TImage > img) |
Returns the minimum pixel value of an image. | |
static std::string | milx::Image< TImage >::ImageOrientation (itk::SmartPointer< TImage > img) |
Returns the orientation flag of an image. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::CheckerBoard (itk::SmartPointer< TImage > img, itk::SmartPointer< TImage > imgToCheckerBoard, const int numberOfSquares=0) |
Generates a checker board pattern where each box alternates between two images. Ideal for comparison among images. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::RescaleIntensity (itk::SmartPointer< TImage > img, float minValue, float maxValue) |
static itk::SmartPointer< TImage > | milx::Image< TImage >::InvertIntensity (itk::SmartPointer< TImage > img, float maxValue) |
Generates an image with the intensities reversed based on the max pixel value. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::HistogramEqualisation (itk::SmartPointer< TImage > img, float alpha=0.3, float beta=0.3, float radius=5) |
Generates an image with the intensities after histogram equalisation. Defaults to classic histogram equalisation. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::GradientMagnitude (itk::SmartPointer< TImage > img) |
Generates an image with the gradient magnitude of the given image. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::SobelEdges (itk::SmartPointer< TImage > img) |
Generates an image with Sobel edges, i.e. uses the Sobel edge detection on the input image. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::Laplacian (itk::SmartPointer< TImage > img) |
Generates an Laplacian image from the input image. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::CannyEdges (itk::SmartPointer< TImage > img, float variance, float lowerThreshold, float upperThreshold) |
Generates an Canny edge image from the input image. | |
static itk::SmartPointer< itk::Image< float, TImage::ImageDimension > > | milx::Image< TImage >::Normalization (itk::SmartPointer< TImage > img) |
Generates an image which is statistically normalised so having pixel values between -1 and 1. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::MatchInformation (itk::SmartPointer< TImage > img, itk::SmartPointer< TImage > imgToMatch, bool originOnly=false) |
Changes the image info to match that of provided image. By default, only the spacing, region and origin are changed. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::CopyInformation (itk::SmartPointer< TImage > img, itk::SmartPointer< TImage > imgToMatch, bool originOnly=false) |
Changes the image info to match that of provided image. By default, only the spacing, region and origin are changed. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::MatchHistogram (itk::SmartPointer< TImage > img, itk::SmartPointer< TImage > imgToMatch, const int bins=128) |
Changes the image gray levels to match histogram of image provided. | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::DistanceMap (itk::SmartPointer< TImage > img, const bool binaryImage=false, const bool signedDistances=false, const bool computeInsideObject=false, const bool squaredDistance=false) |
Generates a distance map image, where the distances are from the object boundary to the outside or vice versa. More... | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::FlipImage (itk::SmartPointer< TImage > img, bool xAxis=false, bool yAxis=true, bool zAxis=false, bool aboutOrigin=true) |
Generates an image with axes flipped. Currently flips the y-axis as default, but can do other axes also in combination. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::PadImageByConstant (itk::SmartPointer< TImage > img, size_t xAxis, size_t yAxis, size_t zAxis, typename TImage::PixelType value) |
Generates an image padded by extending each axes by size given in both directions. The created areas have the value given by value. | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::AnisotropicDiffusion (itk::SmartPointer< TImage > img, const int iterations, const float timestep) |
Generates the gradient anisotropic diffusion (smoothing) of an image using the number of iterations. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::Bilateral (itk::SmartPointer< TImage > img, const float sigmaRange, const float sigmaSpatial) |
Generates the (non-linear) Bilateral smoothing of an image using the sigmas provided. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::GaussianSmooth (itk::SmartPointer< TImage > img, const float variance) |
Generates the Gaussian smoothing (by convolution) of an image using the variance provided. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::Median (itk::SmartPointer< TImage > img, const int radius) |
Generates the (non-linear) median filtering (smoothing) of an image of given neighbourhood radius. | |
template<typename TMaskImage > | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::MaskImage (itk::SmartPointer< TImage > img, itk::SmartPointer< TMaskImage > maskImg) |
Mask an image by given binary image. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::RelabelImage (itk::SmartPointer< TImage > labelledImg) |
Returns a labelled image with labelled values relabelled consecutively based on connectivity. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ThresholdAboveImage (itk::SmartPointer< TImage > img, float outsideValue, float aboveValue) |
Generates an image with the intensities above a certain level thresholded (capped) to the outsideValue. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ThresholdBelowImage (itk::SmartPointer< TImage > img, float outsideValue, float aboveValue) |
Generates an image with the intensities below a certain level thresholded (capped) to the outsideValue. | |
static itk::SmartPointer< TImage > | milx::Image< TImage >::ThresholdImage (itk::SmartPointer< TImage > img, float outsideValue, float belowValue, float aboveValue) |
Generates an image with the intensities below and above a certain level thresholded (capped) to the outsideValue. | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::BinaryThresholdImage (itk::SmartPointer< TImage > img, float outsideValue, float insideValue, float belowValue, float aboveValue) |
Generates a binary image thresholded at the intensities below and above a certain level to the outsideValue. | |
static double | milx::Image< TImage >::OtsuThreshold (itk::SmartPointer< TImage > img, const int bins) |
Returns the Otsu threshold of an image of given the number of histogram bins. | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::OtsuThresholdImage (itk::SmartPointer< TImage > img, const int bins) |
Generates the Otsu threshold of an image of given the number of histogram bins. More... | |
template<typename TOutImage > | |
static itk::SmartPointer< TOutImage > | milx::Image< TImage >::OtsuMultipleThresholdImage (itk::SmartPointer< TImage > img, const int bins, const int noOfLabels=1) |
Generates the multiple Otsu threshold of an image of given the number of histogram bins. More... | |
template<typename TScalarImage > | |
static itk::SmartPointer< TScalarImage > | milx::Image< TImage >::VectorMagnitude (itk::SmartPointer< TImage > img) |
Generates an image comprised of the magnitude of the vectors in the image. | |
void | milx::Model::Clean () |
Removes duplicate points etc. from model. More... | |
void | milx::Model::Triangulate () |
Ensures the polygons in the model are triangulated. Useful when some filters required triangulation before use. More... | |
void | milx::Model::MatchInformation (vtkSmartPointer< vtkPolyData > matchMesh, const bool rescale) |
Matches the centroid and centroid scale of the model to the one given. | |
void | milx::Model::IterativeClosestPointsRegistration (vtkSmartPointer< vtkPolyData > fixedMesh, const bool similarity=false, const int maxPoints=1024) |
Computes the ICP registration of the current mesh to the one provided. | |
void | milx::Model::LandmarkBasedRegistration (vtkSmartPointer< vtkPolyData > fixedMesh, const bool similarity=false) |
Computes the landmark transform registration of the current mesh to the one provided. Requires point correspondence and same number of points. | |
void | milx::Model::Decimate (const float reduceFactor) |
Decimate the mesh (in terms of points) using the Decimate PRO. Generally use QuadricDecimate() instead. More... | |
void | milx::Model::QuadricDecimate (const float reduceFactor) |
Decimate the mesh (in terms of points) using the Quadric algorithm. More... | |
void | milx::Model::ClusterDecimate () |
Decimate the mesh (in terms of points) using the Quadric Clustering algorithm. More... | |
void | milx::Model::LaplacianSmoothing (const int iterations) |
Smooth the mesh (in terms of points) using the Laplacian algorithm. It is unstable, use WindowedSincSmoothing() instead unless massive non-edge preserving smoothing is required. | |
void | milx::Model::WindowedSincSmoothing (const int iterations, const float passBand=0.1) |
Smooth the mesh (in terms of points) using the Taubin (1995) optimal filter. | |
void | milx::Model::Curvature (const bool meanCurvature=true) |
Compute the curvature (using Laplace-Beltrami operator) of the mesh. | |
void | milx::Model::Gradient () |
Compute the gradient (using Laplace-Beltrami operator) of the scalars on the mesh. | |
void | milx::Model::DistanceMap (vtkSmartPointer< vtkPoints > loopPoints) |
Compute the distance map of the mesh points given a loop (as points) assumed to be on the mesh surface. | |
void | milx::Model::GaussianWeightScalars (float stddev, float clampValue) |
Evaluate the Gaussian function with the given standard deviation for all values of the scalars on the mesh. More... | |
void | milx::Model::Rotate (const bool xAxis, const bool yAxis, const bool zAxis, const float angle, const coordinate centre) |
Rotate the mesh along the selected axes from the centre by angle. More... | |
void | milx::Model::Flip (const bool xAxis, const bool yAxis, const bool zAxis) |
Flip the mesh along the selected axes. More... | |
void | milx::Model::Clip (const coordinateType belowValue, const coordinateType aboveValue) |
Clip the mesh based on a scalar threshold of the mesh scalars. More... | |
void | milx::Model::Threshold (const coordinateType belowVal, const coordinateType aboveVal, const coordinateType outsideVal) |
Threshold the scalars on the mesh based on a levels provided. | |
void | milx::Model::Threshold (const coordinateType belowVal, const coordinateType aboveVal, const coordinateType insideVal, const coordinateType outsideVal) |
Threshold the scalars on the mesh based on a levels provided. Inside value means a binary threshold is effectively done. | |
void | milx::Model::Mask (vtkSmartPointer< vtkPolyData > maskMesh) |
Mask the scalars on the mesh based on values of maskMesh being >= 1. | |
void | milx::Model::IsoSurface (vtkSmartPointer< vtkImageData > img, double minValue) |
vtkSmartPointer< vtkImageData > | milx::Model::Voxelise (const unsigned char insideValue, double *spacing, double *bounds=NULL, const size_t padVoxels=1) |
void | milx::Model::DelaunayGraph3D () |
Compute the Delaunay 3D graph of the points in the mesh. More... | |
void | milx::Model::Delaunay2DTriangulation () |
Compute the Delaunay 2D triangluation of the points in the mesh. More... | |
void | milx::Model::DelaunayTriangulation () |
Compute the Delaunay 3D triangluation of the points in the mesh. | |
void | milx::Model::Append (vtkSmartPointer< vtkPolyData > model) |
Appends the model provided to the member into the current model. | |
static void | milx::Model::Scale (vtkSmartPointer< vtkPolyData > model, const coordinateType scale, const bool useCentroid=true) |
Scales the coordinates of each point in the model provided by scale. More... | |
void | milx::Model::Scale (const coordinateType scale, const bool useCentroid=true) |
Members for atomic (one image at a time) operations.
Members for atomic (one model at a time) operations.
|
static |
Applies orientation/direction and origin to a VTK image from a reference image.
Flip y-axis to comply with VTK coordinate system? The transformMatrix, where the extraneous transforms (from the image orientation, such as flipping) will be stored, need not be pre-allocated.
Definition at line 1034 of file milxImage.h.
void milx::Model::Clean | ( | ) |
Removes duplicate points etc. from model.
Triangulate first
Definition at line 286 of file milxModel.cxx.
void milx::Model::Clip | ( | const coordinateType | belowValue, |
const coordinateType | aboveValue | ||
) |
Clip the mesh based on a scalar threshold of the mesh scalars.
connect raw threshold output
Definition at line 675 of file milxModel.cxx.
void milx::Model::ClusterDecimate | ( | ) |
Decimate the mesh (in terms of points) using the Quadric Clustering algorithm.
This is a one-step mesh simiplification operation. See QuadricDecimate() for similar approach with more control.
Decimate using Quadric Decimation
Definition at line 454 of file milxModel.cxx.
|
inlinestatic |
Changes the image info to match that of provided image. By default, only the spacing, region and origin are changed.
Set the originOnly argument true to only match the origin. Function is alias for MatchInformation().
Definition at line 513 of file milxImage.h.
void milx::Model::Decimate | ( | const float | reduceFactor | ) |
Decimate the mesh (in terms of points) using the Decimate PRO. Generally use QuadricDecimate() instead.
Decimate using Quadric Decimation
Definition at line 416 of file milxModel.cxx.
void milx::Model::Delaunay2DTriangulation | ( | ) |
Compute the Delaunay 2D triangluation of the points in the mesh.
Generate Delaunay tessellation
Definition at line 885 of file milxModel.cxx.
void milx::Model::DelaunayGraph3D | ( | ) |
Compute the Delaunay 3D graph of the points in the mesh.
Generate Delaunay Graph
Now we just pretty the mesh up with tubed edges and balls at the vertices.
Tubes for the edges
Definition at line 853 of file milxModel.cxx.
|
static |
Generates a distance map image, where the distances are from the object boundary to the outside or vice versa.
The inside is considered as having negative distances. Outside is treated as having positive distances. Use the signed argument to get signed distances using the Maurer method (TPAMI, 2003). Using the computeInside argument to compute distance within the object rather than outside. If binary image is set to true, integer only computation is done if possible.
Definition at line 2446 of file milxImage.h.
|
static |
Extract a component from the current vector image given the component index.
No checks are done.
Definition at line 1358 of file milxImage.h.
|
static |
Extract a sub image (such as a slice) from the current image given the extent (typically obtained from a VTK image actor etc.).
Note that extent must be of length 2*dimension. No checks are done.
Pick the hyper plane
Definition at line 1324 of file milxImage.h.
void milx::Model::Flip | ( | const bool | xAxis, |
const bool | yAxis, | ||
const bool | zAxis | ||
) |
Flip the mesh along the selected axes.
connect raw threshold output
Definition at line 643 of file milxModel.cxx.
void milx::Model::GaussianWeightScalars | ( | float | stddev, |
float | clampValue | ||
) |
Evaluate the Gaussian function with the given standard deviation for all values of the scalars on the mesh.
Clamp value is the lowest allowed value, so that the scalars will never have a value below this value.
Definition at line 584 of file milxModel.cxx.
|
static |
Generates an image with the intensities after histogram equalisation. Defaults to classic histogram equalisation.
Use parameters to change equalisation type/style. See reference "Adaptive Image Contrast Enhancement using Generalizations of Histogram Equalization." J.Alex Stark. IEEE Transactions on Image Processing, May 2000.
Definition at line 2258 of file milxImage.h.
|
static |
Imports a VNL matrix to an ITK image object.
Assumes the vector is not empty and that if the image is provided, it is setup correctly (spacing, origin etc.).
Definition at line 1234 of file milxImage.h.
|
static |
Imports a VNL vector to an ITK image object.
Assumes the vector is not empty and that if the image is provided, it is setup correctly (spacing, origin etc.)
Definition at line 1195 of file milxImage.h.
void milx::Model::IsoSurface | ( | vtkSmartPointer< vtkImageData > | img, |
double | minValue | ||
) |
Iso surface uses marching cubes to generate surface. The minValue is essentially a threshold below exclusive of its value.
Definition at line 744 of file milxModel.cxx.
|
static |
Changes the image info to match that of provided image. By default, only the spacing, region and origin are changed.
Set the originOnly argument true to only match the origin.
Definition at line 2386 of file milxImage.h.
|
static |
Generates the multiple Otsu threshold of an image of given the number of histogram bins.
Uses binary threshold filter internally to threshold image above the multiple Otsu thresholds found.
Definition at line 3033 of file milxImage.h.
|
static |
Generates the Otsu threshold of an image of given the number of histogram bins.
Uses binary threshold filter internally to threshold image above the Otsu threshold found.
Definition at line 3009 of file milxImage.h.
void milx::Model::QuadricDecimate | ( | const float | reduceFactor | ) |
Decimate the mesh (in terms of points) using the Quadric algorithm.
Decimate using Quadric Decimation
Definition at line 435 of file milxModel.cxx.
|
static |
Resizes current image using current spacing.
Resizes current image using given spacing.
Definition at line 1379 of file milxImage.h.
void milx::Model::Rotate | ( | const bool | xAxis, |
const bool | yAxis, | ||
const bool | zAxis, | ||
const float | angle, | ||
const coordinate | centre | ||
) |
Rotate the mesh along the selected axes from the centre by angle.
connect raw threshold output
Definition at line 608 of file milxModel.cxx.
|
static |
Scales the coordinates of each point in the model provided by scale.
Scales the coordinates of each point in the internal model provided by scale.
useCentroid enables if the scaling is to be done wrt the centroid of the model. Note: The scale is just applied to the model provided and is not copied internally like the filters! Note: Does not keep track of Previous Result.
useCentroid enables if the scaling is to be done wrt the centroid of the model. Note: Does not keep track of Previous Result.
Definition at line 1476 of file milxModel.cxx.
|
static |
Subsamples or shrinks the current image using the factors provided.
Using a factor of 2 reduces the size of the dimension in the image by half. factors must match image dimension, this is why SizeType is used.
Definition at line 1459 of file milxImage.h.
|
static |
Transform the image into a new reference image space given the transform.
Supports Affine, Versor Rigid 3D, Rigid 3D, Centered Euler 3D and Euler 3D transforms Interpolation variable: 0 - NN, 1 - Linear, 2 - BSpline
Up cast transform to correct type
Definition at line 1484 of file milxImage.h.
|
static |
Up cast transform to correct type
Definition at line 1652 of file milxImage.h.
void milx::Model::Triangulate | ( | ) |
Ensures the polygons in the model are triangulated. Useful when some filters required triangulation before use.
Triangulate first
Definition at line 305 of file milxModel.cxx.
vtkSmartPointer< vtkImageData > milx::Model::Voxelise | ( | const unsigned char | insideValue, |
double * | spacing, | ||
double * | bounds = NULL , |
||
const size_t | padVoxels = 1 |
||
) |
Voxelises mode. Returns NULL is failed.
Definition at line 766 of file milxModel.cxx.