91#ifndef vtkMeshQuality_h
92#define vtkMeshQuality_h
95#include "vtkFiltersVerdictModule.h"
100#define VTK_QUALITY_EDGE_RATIO 0
101#define VTK_QUALITY_ASPECT_RATIO 1
102#define VTK_QUALITY_RADIUS_RATIO 2
103#define VTK_QUALITY_ASPECT_FROBENIUS 3
104#define VTK_QUALITY_MED_ASPECT_FROBENIUS 4
105#define VTK_QUALITY_MAX_ASPECT_FROBENIUS 5
106#define VTK_QUALITY_MIN_ANGLE 6
107#define VTK_QUALITY_COLLAPSE_RATIO 7
108#define VTK_QUALITY_MAX_ANGLE 8
109#define VTK_QUALITY_CONDITION 9
110#define VTK_QUALITY_SCALED_JACOBIAN 10
111#define VTK_QUALITY_SHEAR 11
112#define VTK_QUALITY_RELATIVE_SIZE_SQUARED 12
113#define VTK_QUALITY_SHAPE 13
114#define VTK_QUALITY_SHAPE_AND_SIZE 14
115#define VTK_QUALITY_DISTORTION 15
116#define VTK_QUALITY_MAX_EDGE_RATIO 16
117#define VTK_QUALITY_SKEW 17
118#define VTK_QUALITY_TAPER 18
119#define VTK_QUALITY_VOLUME 19
120#define VTK_QUALITY_STRETCH 20
121#define VTK_QUALITY_DIAGONAL 21
122#define VTK_QUALITY_DIMENSION 22
123#define VTK_QUALITY_ODDY 23
124#define VTK_QUALITY_SHEAR_AND_SIZE 24
125#define VTK_QUALITY_JACOBIAN 25
126#define VTK_QUALITY_WARPAGE 26
127#define VTK_QUALITY_ASPECT_GAMMA 27
128#define VTK_QUALITY_AREA 28
129#define VTK_QUALITY_ASPECT_BETA 29
159 vtkSetMacro(TriangleQualityMeasure,
int);
160 vtkGetMacro(TriangleQualityMeasure,
int);
225 vtkSetMacro(QuadQualityMeasure,
int);
226 vtkGetMacro(QuadQualityMeasure,
int);
290 vtkSetMacro(TetQualityMeasure,
int);
291 vtkGetMacro(TetQualityMeasure,
int);
337 vtkSetMacro(HexQualityMeasure,
int);
338 vtkGetMacro(HexQualityMeasure,
int);
772 if (!((cv != 0) ^ (this->Volume != 0)))
780 this->CompatibilityModeOn();
817 if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
821 this->CompatibilityMode = cm;
823 if (this->CompatibilityMode)
854 static double CurrentTriNormal[3];
abstract class to specify cell behavior
abstract superclass for arrays of numeric data
Superclass for algorithms that produce output of the same type as input.
a simple class to control print indentation
Calculate functions of quality of the elements of a mesh.
void SetQuadQualityMeasureToSkew()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double QuadWarpage(vtkCell *cell)
void SetHexQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadTaper(vtkCell *cell)
void SetTriangleQualityMeasureToArea()
Set/Get the particular estimator used to function the quality of triangles.
static double HexOddy(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double TetAspectGamma(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
void SetTetQualityMeasureToMinAngle()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetHexQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexShear(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double HexScaledJacobian(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetTriangleQualityMeasureToRadiusRatio()
Set/Get the particular estimator used to function the quality of triangles.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
void SetQuadQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToDiagonal()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetHexQualityMeasureToMedAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTetQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of tetrahedra.
virtual void SetRatio(vtkTypeBool r)
These methods are deprecated.
static double QuadJacobian(vtkCell *cell)
void SetTriangleQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a triangle.
void SetHexQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadOddy(vtkCell *cell)
static double QuadAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a planar quadrilateral.
static double TriangleAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the Frobenius condition number of the transformation matr...
static double QuadScaledJacobian(vtkCell *cell)
static double HexMaxEdgeRatio(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexShapeAndSize(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double QuadShear(vtkCell *cell)
void SetHexQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToAspectRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleShapeAndSize(vtkCell *cell)
This is a static function used to calculate the product of shape and relative size of a triangle.
static int GetCurrentTriangleNormal(double point[3], double normal[3])
A function called by some VERDICT triangle quality functions to test for inverted triangles.
static double HexTaper(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static vtkMeshQuality * New()
vtkDataArray * CellNormals
static double TetAspectBeta(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TetJacobian(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double HexDistortion(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double TetCollapseRatio(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TriangleDistortion(vtkCell *cell)
This is a static function used to calculate the distortion of a triangle.
static double HexVolume(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadArea(vtkCell *cell)
void SetTriangleQualityMeasureToAspectRatio()
Set/Get the particular estimator used to function the quality of triangles.
void SetHexQualityMeasureToStretch()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a tetrahedron.
void SetQuadQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexRelativeSizeSquared(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetTetQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetQuadQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexShape(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TriangleScaledJacobian(vtkCell *cell)
This is a static function used to calculate the scaled Jacobian of a triangle.
void SetHexQualityMeasureToDimension()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetRelativeSizeSquared(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TriangleShape(vtkCell *cell)
This is a static function used to calculate the shape of a triangle.
void SetTetQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadMedAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the average Frobenius aspect of the 4 corner triangles of...
static double TriangleCondition(vtkCell *cell)
This is a static function used to calculate the condition number of a triangle.
void SetQuadQualityMeasureToTaper()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TetMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) dihedral angle of a tetrahedron...
static double TetShapeandSize(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
void SetTetQualityMeasureToAspectRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTriangleQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to function the quality of triangles.
void SetHexQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double QuadCondition(vtkCell *cell)
void SetQuadQualityMeasureToRadiusRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetQuadQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToSkew()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToShearAndSize()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
~vtkMeshQuality() override
void SetQuadQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a triangle.
static double HexJacobian(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetQuadQualityMeasureToMedAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TetDistortion(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double QuadShapeAndSize(vtkCell *cell)
static double TriangleMaxAngle(vtkCell *cell)
This is a static function used to calculate the maximal (nonoriented) angle of a triangle,...
void SetTriangleQualityMeasureToCondition()
Set/Get the particular estimator used to function the quality of triangles.
static double TetCondition(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TriangleRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a triangle.
static double TetScaledJacobian(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double QuadEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a quadrilateral.
void SetTetQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadShearAndSize(vtkCell *cell)
static double HexDimension(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToTaper()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTriangleQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexStretch(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double QuadMaxAngle(vtkCell *cell)
static double QuadStretch(vtkCell *cell)
void SetTriangleQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to function the quality of triangles.
virtual void SetCompatibilityMode(vtkTypeBool cm)
CompatibilityMode governs whether, when both a quality function and cell volume are to be stored as c...
static double HexEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a hexahedron.
void SetTetQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetHexQualityMeasureToShearAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 4 corner triangles of...
vtkTypeBool GetVolume()
These methods are deprecated.
vtkTypeBool SaveCellQuality
void SetHexQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadMaxEdgeRatios(vtkCell *cell)
static double TetAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the Frobenius condition number of the transformation matr...
static double HexMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetQuadQualityMeasureToArea()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetQuadQualityMeasureToWarpage()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexShearAndSize(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double HexSkew(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetTetQualityMeasureToCollapseRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTriangleQualityMeasureToMaxAngle()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToMaxAngle()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
vtkTypeBool CompatibilityMode
void SetTetQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadShape(vtkCell *cell)
void SetQuadQualityMeasureToMaxEdgeRatios()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTetQualityMeasureToRadiusRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetHexQualityMeasureToOddy()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexMedAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the average Frobenius aspect of the 8 corner tetrahedra o...
void SetTetQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetHexQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TriangleRelativeSizeSquared(vtkCell *cell)
This is a static function used to calculate the square of the relative size of a triangle.
static double QuadRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a planar quadrilateral.
void SetHexQualityMeasureToMaxEdgeRatios()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadRelativeSizeSquared(vtkCell *cell)
static double TetVolume(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double QuadSkew(vtkCell *cell)
static double QuadDistortion(vtkCell *cell)
void SetTetQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetQuadQualityMeasureToOddy()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTetQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
int TriangleQualityMeasure
void SetHexQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTetQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetQuadQualityMeasureToMinAngle()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a triangle,...
static double HexDiagonal(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetTriangleQualityMeasureToShape()
Set/Get the particular estimator used to function the quality of triangles.
void SetTriangleQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to function the quality of triangles.
static double TriangleArea(vtkCell *cell)
This is a static function used to calculate the area of a triangle.
void SetTetQualityMeasureToAspectGamma()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a quadrilateral,...
void SetQuadQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTriangleQualityMeasureToDistortion()
Set/Get the particular estimator used to function the quality of triangles.
void SetTriangleQualityMeasureToMinAngle()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToStretch()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TetRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a tetrahedron.
void SetTetQualityMeasureToAspectBeta()
Set/Get the particular estimator used to measure the quality of tetrahedra.
virtual void SetVolume(vtkTypeBool cv)
These methods are deprecated.
static double TetAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a tetrahedron.
static double TetShape(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
void SetQuadQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexCondition(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
virtual void Modified()
Update the modification time for this object.
#define VTK_QUALITY_SHAPE_AND_SIZE
#define VTK_QUALITY_STRETCH
#define VTK_QUALITY_MAX_ANGLE
#define VTK_QUALITY_MIN_ANGLE
#define VTK_QUALITY_MAX_ASPECT_FROBENIUS
#define VTK_QUALITY_SHEAR_AND_SIZE
#define VTK_QUALITY_RELATIVE_SIZE_SQUARED
#define VTK_QUALITY_JACOBIAN
#define VTK_QUALITY_SHEAR
#define VTK_QUALITY_ASPECT_GAMMA
#define VTK_QUALITY_VOLUME
#define VTK_QUALITY_EDGE_RATIO
#define VTK_QUALITY_DISTORTION
#define VTK_QUALITY_SHAPE
#define VTK_QUALITY_WARPAGE
#define VTK_QUALITY_RADIUS_RATIO
#define VTK_QUALITY_MED_ASPECT_FROBENIUS
#define VTK_QUALITY_CONDITION
#define VTK_QUALITY_DIAGONAL
#define VTK_QUALITY_SCALED_JACOBIAN
#define VTK_QUALITY_COLLAPSE_RATIO
#define VTK_QUALITY_TAPER
#define VTK_QUALITY_DIMENSION
#define VTK_QUALITY_ASPECT_BETA
#define VTK_QUALITY_MAX_EDGE_RATIO
#define VTK_QUALITY_ASPECT_RATIO
#define VTK_QUALITY_ASPECT_FROBENIUS