SMILX  1.01
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
milx::Model Class Reference

Represents a model (i.e. a model with cells and scalar values) and their common operations. Also allows batch operations on collection of models. More...

#include <milxModel.h>

Inheritance diagram for milx::Model:
milx::DeformableModel

Public Member Functions

 Model ()
 Standard constructor.
 
 Model (vtkSmartPointer< vtkPolyData > model)
 Constructor that copies the input model.
 
virtual ~Model ()
 Standard Destructor.
 
void AddInput (vtkSmartPointer< vtkPolyData > model)
 Assigns and coninually appends meshes provided to the class.
 
void AddArray (vtkSmartPointer< vtkDataArray > array)
 Appends an array to the model.
 
void SetInput (vtkSmartPointer< vtkPolyData > model)
 Assigns the input model to class. More...
 
void SetInputData (vtkSmartPointer< vtkPolyData > model)
 Assigns the input model to class. Same as SetInput and meant for VTK version > 5.
 
void SetTransform (vtkSmartPointer< vtkAbstractTransform > newTransform)
 Transforms the model by given transform.
 
void SetGraph (vtkSmartPointer< vtkMutableUndirectedGraph > graph)
 Converts a graph to polydata/mesh ready for display.
 
void SetPoints (vtkSmartPointer< vtkPoints > modelPoints)
 Assign points to the model.
 
void SetPolys (vtkSmartPointer< vtkCellArray > modelPolys)
 Assign polygons to the model via the array given.
 
void SetVectors (vtkSmartPointer< vtkDataArray > modelVectors)
 Sets the vector field of the mesh (vectors per vertex of the surface)
 
void SetScalars (vtkSmartPointer< vtkDataArray > modelScalars)
 Sets the scalar field of the mesh (values per vertex of the surface)
 
void SetActiveScalars (std::string nameOfArray)
 Sets the array of name as the scalar field of the mesh (values per vertex of the surface)
 
vtkSmartPointer< vtkPolyData > & Result ()
 Returns the current model, i.e. the result of the latest operation. More...
 
vtkSmartPointer< vtkPolyData > & PreviousResult ()
 Returns the previous model, i.e. the result of the penultimate operation. More...
 
vtkSmartPointer< vtkPolyData > & GetInput ()
 Returns the input model, i.e. the very most initial model. More...
 
vtkSmartPointer< vtkPolyData > & GetOutput ()
 Returns the current model, i.e. the result of the latest operation ITK/VTK style.
 
vtkAlgorithmOutput * GetProducerPort ()
 Returns the current model, i.e. the result of the latest operation ITK/VTK style.
 
vtkIdType GetNumberOfPoints ()
 Returns the total number of points of the model currently held.
 
vtkIdType GetNumberOfArrays ()
 Returns the total number of point-based arrays in the model.
 
vtkSmartPointer< vtkTransform > & GetTransform ()
 Returns the current net transform applied to the model, i.e. concatenations of all (rigid/linear) transforms.
 
vtkSmartPointer< vtkPoints > GetPoints ()
 Returns the points of the model.
 
vtkSmartPointer< vtkDataArray > GetScalars ()
 Returns the active scalars of the model.
 
void GetScalarRange (double *range)
 Returns the scalar range of the model.
 
vtkSmartPointer< vtkDataArray > GetVectors ()
 Returns the active vectors of the model.
 
vtkSmartPointer< vtkDataArray > GetNormals ()
 Returns the normal vectors of the model.
 
vtkSmartPointer< vtkDataArray > GetTensors ()
 Returns the active tensors (3x3 matrix per vertex) of the model.
 
void Update ()
 Updates the model to be the most current.
 
void Clean ()
 Removes duplicate points etc. from model. More...
 
void Triangulate ()
 Ensures the polygons in the model are triangulated. Useful when some filters required triangulation before use. More...
 
void MatchInformation (vtkSmartPointer< vtkPolyData > matchMesh, const bool rescale)
 Matches the centroid and centroid scale of the model to the one given.
 
void 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 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 Decimate (const float reduceFactor)
 Decimate the mesh (in terms of points) using the Decimate PRO. Generally use QuadricDecimate() instead. More...
 
void QuadricDecimate (const float reduceFactor)
 Decimate the mesh (in terms of points) using the Quadric algorithm. More...
 
void ClusterDecimate ()
 Decimate the mesh (in terms of points) using the Quadric Clustering algorithm. More...
 
void 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 WindowedSincSmoothing (const int iterations, const float passBand=0.1)
 Smooth the mesh (in terms of points) using the Taubin (1995) optimal filter.
 
void Curvature (const bool meanCurvature=true)
 Compute the curvature (using Laplace-Beltrami operator) of the mesh.
 
void Gradient ()
 Compute the gradient (using Laplace-Beltrami operator) of the scalars on the mesh.
 
void 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 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 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 Flip (const bool xAxis, const bool yAxis, const bool zAxis)
 Flip the mesh along the selected axes. More...
 
void Clip (const coordinateType belowValue, const coordinateType aboveValue)
 Clip the mesh based on a scalar threshold of the mesh scalars. More...
 
void Threshold (const coordinateType belowVal, const coordinateType aboveVal, const coordinateType outsideVal)
 Threshold the scalars on the mesh based on a levels provided.
 
void 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 Mask (vtkSmartPointer< vtkPolyData > maskMesh)
 Mask the scalars on the mesh based on values of maskMesh being >= 1.
 
void IsoSurface (vtkSmartPointer< vtkImageData > img, double minValue)
 
vtkSmartPointer< vtkImageData > Voxelise (const unsigned char insideValue, double *spacing, double *bounds=NULL, const size_t padVoxels=1)
 
void DelaunayGraph3D ()
 Compute the Delaunay 3D graph of the points in the mesh. More...
 
void Delaunay2DTriangulation ()
 Compute the Delaunay 2D triangluation of the points in the mesh. More...
 
void DelaunayTriangulation ()
 Compute the Delaunay 3D triangluation of the points in the mesh.
 
void Append (vtkSmartPointer< vtkPolyData > model)
 Appends the model provided to the member into the current model.
 
void Scale (const coordinateType scale, const bool useCentroid=true)
 
void GenerateVertices ()
 Generate (only) vertices for display from point data. Good for when only points in mesh.
 
void GenerateVertexScalars ()
 Generate vertex scalars for display from point data. This colours the mesh by point ids.
 
void GenerateVerticesAs2DGlyphs (const GlyphType glyphType)
 Generate (only) vertices (as 2D glyphs like crosses etc.) for display from point data. Good for when only points in mesh.
 
void GenerateTubes ()
 Generate tubes along edges of the mesh.
 
void GeneratePointModel (const double newScale)
 Generate spheres of scaling for each point of the mesh.
 
void GenerateSpheres (const double newScale)
 Generate spheres of scaling for each point of the mesh. Same as GeneratePointModel().
 
void GenerateSampledPoints (const double distance)
 Generate sampled points (only points) for the mesh at given spacing on triangles.
 
void GenerateVectorField (const double newScale=0.0, const bool useNormals=false)
 Generate vector field for vectors present (as an array) in the mesh at given scaling. More...
 
void GenerateTensorField (const double newScale=0.0)
 Generate tensor field for tensors present (as an array) in the mesh at given scaling. More...
 
void GenerateHedgehog (const double newScale=0.0)
 Generate a line field for vectors present (as an array) in the mesh at given scaling. More...
 
void GenerateNormals (const int pointNormals=0)
 Generate normal vectors for each point on the mesh. If pointNormals is (0,1,2), then normals are for (points,cells,both).
 
void GenerateKMeansClustering (int numberOfClusters)
 Generate k-means clustering for the point positions in the mesh.
 
void GenerateQuantisedPoints (float quantiseFactor)
 Generate points of the mesh to have integral cordinates.
 
void GenerateCappedBoundaries ()
 Generate capped boundaries for open meshes. This closes off open ends of a mesh.
 
void GenerateRegions ()
 Generate connected regions for the mesh labelled by scalar values.
 
void GenerateElevation (double x, double y, double z)
 Generate scalar field for the mesh proportional to the elevation (dot product) with a given vector/direction.
 
void GenerateElevation ()
 Generate scalar field for the mesh proportional to the elevation in the z-direction.
 
void ResetScalars ()
 Resets the scalars to zero of the current model. More...
 
void FillScalars (const coordinateType value)
 Sets all the scalars of the current model to value. More...
 
void RemoveScalars ()
 Removes all the scalars from the current model. Same as ClearScalars()
 
void RemoveArrays ()
 Removes all the arrays from the current model. Same as ClearArrayss()
 
void ConcatenateCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Concatenates all the models in the collection together into one model. More...
 
void ScaleCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType scale)
 Scales the coordinates of each point in all the models in the collection. More...
 
void SmoothCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const size_t iterations)
 Smoothes all points in all the models in the collection using the Windowed Sinc algorithm. More...
 
void LaplacianCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const size_t iterations)
 Smoothes all points in all the models in the collection using the Laplacian algorithm. More...
 
void FlipCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const bool xAxis, const bool yAxis, const bool zAxis)
 Flips all points in all the models in the collection using the axis provided. More...
 
void DecimateCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType factor)
 Decimates all points in all the models in the collection using the Quadric algorithm algorithm. More...
 
void SplitCollection (vtkSmartPointer< vtkPolyDataCollection > collection, vtkSmartPointer< vtkPolyDataCollection > components, std::vector< vtkSmartPointer< vtkPolyDataCollection > > &splitCollections)
 Splits each of the models in the collection into a collection of models. More...
 
void ScalarDifferenceCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Computes the difference in scalars cummulatively over the collection wrt the first mesh assuming the meshes have the same number of points. More...
 
void ScalarDifferenceCollection (vtkSmartPointer< vtkPolyDataCollection > collection1, vtkSmartPointer< vtkPolyDataCollection > collection2, vtkSmartPointer< vtkPolyDataCollection > &results)
 
void ScalarThresholdCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType aboveVal, const coordinateType belowVal)
 Computes the Threshold of scalar values of meshes in place over the collection wrt the first mesh assuming the meshes have the same number of points. More...
 
void ClipCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType aboveVal, const coordinateType belowVal)
 Computes the Clipping of models via the Threshold of scalar values in place over the collection. More...
 
void ScalarStatisticsCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Computes the statistics (min, max, variance etc.) of the scalars over the collection assuming the meshes have the same number of points. More...
 
void ScalarRemoveCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Removes the scalars over the collection assuming the meshes have the same number of points. More...
 
void ScalarCopyCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Copies the scalars over the collection (from the first mesh) assuming the meshes have the same number of points. More...
 
double MeanSquaredErrorCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Computes the MSE of all models in collection accummulatively.
 
vtkSmartPointer< vtkPolyDataCollection > ProcrustesAlignCollection (vtkSmartPointer< vtkPolyDataCollection > collection, bool rigid=false)
 Aligns the collection to the mean mesh (computed internally) of the collection assuming the meshes have the same number of points and correspondence. More...
 
vtkSmartPointer< vtkPolyDataCollection > IterativeClosestPointsAlignCollection (vtkSmartPointer< vtkPolyDataCollection > collection, bool rigid=false, vtkSmartPointer< vtkTransformCollection > tCollection=0)
 Aligns the collection to the mean mesh (computed internally) of the collection without needing correspondence. More...
 

Static Public Member Functions

static void Scale (vtkSmartPointer< vtkPolyData > model, const coordinateType scale, const bool useCentroid=true)
 Scales the coordinates of each point in the model provided by scale. More...
 
static void ClearScalars (vtkSmartPointer< vtkPolyData > &mesh)
 Removes all the scalars from the current model. More...
 
static void ClearArrays (vtkSmartPointer< vtkPolyData > &mesh)
 Removes all the arrays from the current model. More...
 
static void ScalarDifference (vtkSmartPointer< vtkPolyData > model1, vtkSmartPointer< vtkPolyData > model2, vtkSmartPointer< vtkPolyData > result)
 Calculates differences within scalars of two meshes assuming the meshes have the same number of points. More...
 
static void ThresholdScalars (vtkSmartPointer< vtkPolyData > model, const coordinateType aboveVal, const coordinateType belowVal, const coordinateType outsideVal)
 Thresholds scalar values of mesh in place. More...
 
static void BinaryThresholdScalars (vtkSmartPointer< vtkPolyData > model, const coordinateType aboveVal, const coordinateType belowVal, const coordinateType insideVal=1.0, const coordinateType outsideVal=0.0)
 Thresholds scalar values of mesh in place. Scalars outside the given range is set to outsideVal and scalars inside the range set to insideVal. More...
 
static void MaskScalars (vtkSmartPointer< vtkPolyData > model1, vtkSmartPointer< vtkPolyData > model2)
 Masks the scalars in model1 with scalars in model2. Only values corresponding to >= 1 in model2 are kept in model1.
 
static vtkSmartPointer< vtkMatrix4x4 > ProcrustesAlign (vtkSmartPointer< vtkPolyData > source, vtkSmartPointer< vtkPolyData > target, bool rigid=false)
 Aligns the source mesh to target mesh assuming the meshes have the same number of points and correspondence. The transform is applied in-place to the source mesh. More...
 

Protected Member Functions

bool IsCurrentModel ()
 
void UpdateModelState (vtkSmartPointer< vtkPolyDataAlgorithm > filter)
 
void UpdateModelStateFromModel (vtkSmartPointer< vtkPolyData > mesh)
 

Protected Attributes

vtkSmartPointer< vtkPolyData > InputModel
 Holds the initial model in the pipeline.
 
vtkSmartPointer< vtkPolyData > CurrentModel
 Holds the current model in the pipeline.
 
vtkSmartPointer< vtkPolyData > PreviousModel
 Holds the previous model in the pipeline.
 
vtkSmartPointer< vtkTransform > CurrentTransform
 transform applied to model
 
bool InternalInPlaceOperation
 Used by the collection members to assign via pointers rather than deep copys.
 

Detailed Description

Represents a model (i.e. a model with cells and scalar values) and their common operations. Also allows batch operations on collection of models.

The class members are divided into two groups, atmoic and collection members. Atomic operations operate on one model at a time, whereas collection members use the atomic operations in a loop to batch that operation over multiple models.

Unlike the milx::Image class, this class tracks the current state of the data using pointers which can be retrieved using the Result() (or GetOutput()) and PreviousResult() members. This class is used extensively throughout the milxQtModel class and in milxSurfaceDiagnosticApp application. See the usage examples given below:

Processing Pipeline

model.SetInput(surface);
model.Triangulate();
model.LaplacianSmoothing(250); //Kill some of the steps caused by MCs
model.WindowedSincSmoothing(20); //Preserve features
model.Clean();
model.Update();

ICP Registration

...
vtkSmartPointer<vtkPolyData> mesh = m_RawSurfaces->GetNextItem();
milx::Model initialisation(m_AtlasSurface);
initialisation.MatchInformation(mesh, true); //Get correct approx scale
initialisation.IterativeClosestPointsRegistration(mesh, false); //Doesn't handle different FoV's

Landmark based registration

milx::Model model(m_AtlasSurface);
model.LandmarkBasedRegistration(mesh, true);
mesh->DeepCopy(model.GetOutput());

Voxelise model

milx::Model model1(surface1);
vtkSmartPointer<vtkImageData> voxelisedModel1 = model1.Voxelise(255, labelSpacing.GetDataPointer(), bigBounds);

Multi-use example

//Read labels and generate iso surface
//Clip mesh to ensure same FoV
vtkSmartPointer<vtkFloatArray> weights = vtkSmartPointer<vtkFloatArray>::New();
milx::Model model(surface);
model.MarkSurfaceInsideImage<FloatImageType>(surface, distanceMap, weights); //sets weights as scalars in here
model.Clip(1.0, 1.0);
//Save distances as scalars on mesh
const bool absValues = true;
vtkSmartPointer<vtkFloatArray> scalars = model.SurfaceScalarsFromImage<FloatImageType>(model.GetOutput(), distanceMap, absValues);
scalars->SetName("Surface Distance Errors");
model.SetScalars(scalars);
//Restore correspondence to copy scalars over to full (unclipped) mesh, since the mesh is clipped to ensure same FoV
cout << "Regions of Surface Outside the Image are marked with -1." << endl;
weights->FillComponent(0, -1.0);
for(int j = 0; j < model.GetOutput()->GetNumberOfPoints(); j ++)
{
vtkIdType index = surface->FindPoint(model.GetOutput()->GetPoint(j));
weights->SetTuple1(index, scalars->GetTuple1(j));
}

Definition at line 113 of file milxModel.h.


The documentation for this class was generated from the following files: