SMILX  1.01
Functions
Collection operations

Members for operating of a collection of models. More...

Functions

static vtkSmartPointer< vtkMatrix4x4 > milx::Model::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...
 
void milx::Model::ConcatenateCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Concatenates all the models in the collection together into one model. More...
 
void milx::Model::ScaleCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType scale)
 Scales the coordinates of each point in all the models in the collection. More...
 
void milx::Model::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 milx::Model::LaplacianCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const size_t iterations)
 Smoothes all points in all the models in the collection using the Laplacian algorithm. More...
 
void milx::Model::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 milx::Model::DecimateCollection (vtkSmartPointer< vtkPolyDataCollection > collection, const coordinateType factor)
 Decimates all points in all the models in the collection using the Quadric algorithm algorithm. More...
 
void milx::Model::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 milx::Model::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 milx::Model::ScalarDifferenceCollection (vtkSmartPointer< vtkPolyDataCollection > collection1, vtkSmartPointer< vtkPolyDataCollection > collection2, vtkSmartPointer< vtkPolyDataCollection > &results)
 
void milx::Model::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 milx::Model::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 milx::Model::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 milx::Model::ScalarRemoveCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Removes the scalars over the collection assuming the meshes have the same number of points. More...
 
void milx::Model::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 milx::Model::MeanSquaredErrorCollection (vtkSmartPointer< vtkPolyDataCollection > collection)
 Computes the MSE of all models in collection accummulatively.
 
vtkSmartPointer< vtkPolyDataCollection > milx::Model::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 > milx::Model::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...
 

Detailed Description

Members for operating of a collection of models.

Function Documentation

◆ ClipCollection()

void milx::Model::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.

The collection is modified in place.

< clip the model inplace

Definition at line 1784 of file milxModel.cxx.

◆ ConcatenateCollection()

void milx::Model::ConcatenateCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection)

Concatenates all the models in the collection together into one model.

This member calls Append() over the members of the collection to give the concatenated result. The order of the models in the collection will effect the order of the point IDs of the resultant model.

< concatenate the model

Definition at line 1620 of file milxModel.cxx.

◆ DecimateCollection()

void milx::Model::DecimateCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection,
const coordinateType  factor 
)

Decimates all points in all the models in the collection using the Quadric algorithm algorithm.

The collection is modified in place.

< scale the model inplace

Definition at line 1689 of file milxModel.cxx.

◆ FlipCollection()

void milx::Model::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.

The collection is modified in place.

< scale the model inplace

Definition at line 1672 of file milxModel.cxx.

◆ IterativeClosestPointsAlignCollection()

vtkSmartPointer< vtkPolyDataCollection > milx::Model::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.

The entire collection is then registered to the first mesh in a least squares sense. The member returns a collection of aligned meshes and the mesh internally stored is the mean mesh. The rigid argument makes the alignment rigid body, by default similarity is used instead.

Definition at line 2064 of file milxModel.cxx.

◆ LaplacianCollection()

void milx::Model::LaplacianCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection,
const size_t  iterations 
)

Smoothes all points in all the models in the collection using the Laplacian algorithm.

The collection is modified in place.

< scale the model inplace

Definition at line 1655 of file milxModel.cxx.

◆ ProcrustesAlign()

vtkSmartPointer< vtkMatrix4x4 > milx::Model::ProcrustesAlign ( vtkSmartPointer< vtkPolyData >  source,
vtkSmartPointer< vtkPolyData >  target,
bool  rigid = false 
)
static

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.

The transformation matrix is returned, it is NULL if failure was encountered. The meshes are aligned in a least squares sense. The member returns the aligned mesh and is internally allocated. The scalars are preserved from the source. The rigid argument makes the alignment rigid body, by default similarity is used instead.

Definition at line 1594 of file milxModel.cxx.

◆ ProcrustesAlignCollection()

vtkSmartPointer< vtkPolyDataCollection > milx::Model::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.

The entire collection is then registered to one another in a least squares sense. The member returns a collection of aligned meshes and the mesh internally stored is the mean mesh. The mean mesh scalars are the mean of the scalars of all the meshes. The rigid argument makes the alignment rigid body, by default similarity is used instead.

Definition at line 1970 of file milxModel.cxx.

◆ ScalarCopyCollection()

void milx::Model::ScalarCopyCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection)

Copies the scalars over the collection (from the first mesh) assuming the meshes have the same number of points.

The collection is modified in place.

< inplace

Definition at line 1930 of file milxModel.cxx.

◆ ScalarDifferenceCollection() [1/2]

void milx::Model::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.

The result is effectively a sum of differences across the collection with respect to the first model.

< threshold the model inplace

Definition at line 1753 of file milxModel.cxx.

◆ ScalarDifferenceCollection() [2/2]

void milx::Model::ScalarDifferenceCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection1,
vtkSmartPointer< vtkPolyDataCollection >  collection2,
vtkSmartPointer< vtkPolyDataCollection > &  results 
)

< diff the model scalars

Definition at line 1800 of file milxModel.cxx.

◆ ScalarRemoveCollection()

void milx::Model::ScalarRemoveCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection)

Removes the scalars over the collection assuming the meshes have the same number of points.

The collection is modified in place.

< inplace

Definition at line 1917 of file milxModel.cxx.

◆ ScalarStatisticsCollection()

void milx::Model::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.

Use Result() to get result, which is populated by replicating the first mesh with arrays for each statistic.

Definition at line 1831 of file milxModel.cxx.

◆ ScalarThresholdCollection()

void milx::Model::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.

The collection is modified in place.

< threshold the model inplace

Definition at line 1771 of file milxModel.cxx.

◆ ScaleCollection()

void milx::Model::ScaleCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection,
const coordinateType  scale 
)

Scales the coordinates of each point in all the models in the collection.

The collection is modified in place.

< scale the model inplace

Definition at line 1629 of file milxModel.cxx.

◆ SmoothCollection()

void milx::Model::SmoothCollection ( vtkSmartPointer< vtkPolyDataCollection >  collection,
const size_t  iterations 
)

Smoothes all points in all the models in the collection using the Windowed Sinc algorithm.

The collection is modified in place.

< scale the model inplace

Definition at line 1638 of file milxModel.cxx.

◆ SplitCollection()

void milx::Model::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.

The order of the models in the components needs to match the order of the models in the collection.

Definition at line 1706 of file milxModel.cxx.