65 #ifndef vtkMeshQuality_h
66 #define vtkMeshQuality_h
68 #include "vtkFiltersVerdictModule.h"
74 #define VTK_QUALITY_EDGE_RATIO 0
75 #define VTK_QUALITY_ASPECT_RATIO 1
76 #define VTK_QUALITY_RADIUS_RATIO 2
77 #define VTK_QUALITY_ASPECT_FROBENIUS 3
78 #define VTK_QUALITY_MED_ASPECT_FROBENIUS 4
79 #define VTK_QUALITY_MAX_ASPECT_FROBENIUS 5
80 #define VTK_QUALITY_MIN_ANGLE 6
81 #define VTK_QUALITY_COLLAPSE_RATIO 7
82 #define VTK_QUALITY_MAX_ANGLE 8
83 #define VTK_QUALITY_CONDITION 9
84 #define VTK_QUALITY_SCALED_JACOBIAN 10
85 #define VTK_QUALITY_SHEAR 11
86 #define VTK_QUALITY_RELATIVE_SIZE_SQUARED 12
87 #define VTK_QUALITY_SHAPE 13
88 #define VTK_QUALITY_SHAPE_AND_SIZE 14
89 #define VTK_QUALITY_DISTORTION 15
90 #define VTK_QUALITY_MAX_EDGE_RATIO 16
91 #define VTK_QUALITY_SKEW 17
92 #define VTK_QUALITY_TAPER 18
93 #define VTK_QUALITY_VOLUME 19
94 #define VTK_QUALITY_STRETCH 20
95 #define VTK_QUALITY_DIAGONAL 21
96 #define VTK_QUALITY_DIMENSION 22
97 #define VTK_QUALITY_ODDY 23
98 #define VTK_QUALITY_SHEAR_AND_SIZE 24
99 #define VTK_QUALITY_JACOBIAN 25
100 #define VTK_QUALITY_WARPAGE 26
101 #define VTK_QUALITY_ASPECT_GAMMA 27
102 #define VTK_QUALITY_AREA 28
103 #define VTK_QUALITY_ASPECT_BETA 29
120 vtkGetMacro(SaveCellQuality,
int);
134 vtkGetMacro(TriangleQualityMeasure,
int);
206 vtkGetMacro(QuadQualityMeasure,
int);
313 vtkGetMacro(TetQualityMeasure,
int);
393 vtkGetMacro(HexQualityMeasure,
int);
843 virtual void SetRatio(
int r ) { this->SetSaveCellQuality( r ); }
844 int GetRatio() {
return this->GetSaveCellQuality(); }
867 if ( ! ((cv != 0) ^ (this->Volume != 0)) )
875 this->CompatibilityModeOn();
915 if ( !((cm != 0) ^ (this->CompatibilityMode != 0)) )
919 this->CompatibilityMode = cm;
921 if ( this->CompatibilityMode )
927 vtkGetMacro(CompatibilityMode,
int);
952 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()
static double QuadWarpage(vtkCell *cell)
void SetHexQualityMeasureToScaledJacobian()
static double QuadTaper(vtkCell *cell)
void SetTriangleQualityMeasureToArea()
static double HexOddy(vtkCell *cell)
static double TetAspectGamma(vtkCell *cell)
void SetTetQualityMeasureToMinAngle()
void SetHexQualityMeasureToDistortion()
static double HexShear(vtkCell *cell)
static double HexScaledJacobian(vtkCell *cell)
void SetTriangleQualityMeasureToRadiusRatio()
void SetQuadQualityMeasureToRelativeSizeSquared()
void SetHexQualityMeasureToDiagonal()
void SetHexQualityMeasureToMedAspectFrobenius()
void SetTetQualityMeasureToDistortion()
virtual void SetRatio(int r)
These methods are deprecated.
static double QuadJacobian(vtkCell *cell)
void SetTriangleQualityMeasureToScaledJacobian()
void SetQuadQualityMeasureToMaxAspectFrobenius()
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
static double TriangleEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a triangle.
void SetHexQualityMeasureToMaxAspectFrobenius()
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...
virtual void SetCompatibilityMode(int cm)
CompatibilityMode governs whether, when both a quality function and cell volume are to be stored as c...
static double QuadScaledJacobian(vtkCell *cell)
static double HexMaxEdgeRatio(vtkCell *cell)
void SetHexQualityMeasureToJacobian()
static double HexShapeAndSize(vtkCell *cell)
static double QuadShear(vtkCell *cell)
void SetHexQualityMeasureToShape()
void SetQuadQualityMeasureToAspectRatio()
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)
vtkDataArray * CellNormals
static double TetAspectBeta(vtkCell *cell)
static double TetJacobian(vtkCell *cell)
static double HexDistortion(vtkCell *cell)
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)
void SetHexQualityMeasureToShear()
static double QuadArea(vtkCell *cell)
void SetTriangleQualityMeasureToAspectRatio()
void SetHexQualityMeasureToStretch()
static double TetEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a tetrahedron.
void SetQuadQualityMeasureToCondition()
static double HexRelativeSizeSquared(vtkCell *cell)
void SetTetQualityMeasureToShape()
void SetQuadQualityMeasureToShear()
static double HexShape(vtkCell *cell)
void SetHexQualityMeasureToCondition()
static double TriangleScaledJacobian(vtkCell *cell)
This is a static function used to calculate the scaled Jacobian of a triangle.
void SetHexQualityMeasureToDimension()
static double TetRelativeSizeSquared(vtkCell *cell)
static double TriangleShape(vtkCell *cell)
This is a static function used to calculate the shape of a triangle.
void SetTetQualityMeasureToVolume()
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.
virtual void SetVolume(int cv)
These methods are deprecated.
void SetQuadQualityMeasureToTaper()
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)
void SetTetQualityMeasureToAspectRatio()
void SetTriangleQualityMeasureToRelativeSizeSquared()
void SetHexQualityMeasureToRelativeSizeSquared()
void SetQuadQualityMeasureToEdgeRatio()
static double QuadCondition(vtkCell *cell)
void SetQuadQualityMeasureToRadiusRatio()
void SetQuadQualityMeasureToShapeAndSize()
void SetHexQualityMeasureToSkew()
void SetQuadQualityMeasureToShearAndSize()
void SetQuadQualityMeasureToJacobian()
static double TriangleAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a triangle.
static double HexJacobian(vtkCell *cell)
void SetQuadQualityMeasureToMedAspectFrobenius()
static double TetDistortion(vtkCell *cell)
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()
static double TetCondition(vtkCell *cell)
static double TriangleRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a triangle.
static double TetScaledJacobian(vtkCell *cell)
static double QuadEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a quadrilateral.
void SetTetQualityMeasureToAspectFrobenius()
static double QuadShearAndSize(vtkCell *cell)
static double HexDimension(vtkCell *cell)
void SetHexQualityMeasureToTaper()
void SetTriangleQualityMeasureToShapeAndSize()
void SetQuadQualityMeasureToShape()
static double HexStretch(vtkCell *cell)
static double QuadMaxAngle(vtkCell *cell)
static double QuadStretch(vtkCell *cell)
void SetTriangleQualityMeasureToAspectFrobenius()
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called within ProcessRequest when a request asks the algorithm to do its work.
static double HexEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a hexahedron.
void SetTetQualityMeasureToRelativeSizeSquared()
void SetHexQualityMeasureToShearAndSize()
static double QuadMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 4 corner triangles of...
void SetHexQualityMeasureToShapeAndSize()
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()
void SetQuadQualityMeasureToWarpage()
static double HexShearAndSize(vtkCell *cell)
static double HexSkew(vtkCell *cell)
void SetTetQualityMeasureToCollapseRatio()
void SetTriangleQualityMeasureToMaxAngle()
void SetQuadQualityMeasureToMaxAngle()
void SetTetQualityMeasureToJacobian()
static double QuadShape(vtkCell *cell)
void SetQuadQualityMeasureToMaxEdgeRatios()
static vtkMeshQuality * New()
void SetTetQualityMeasureToRadiusRatio()
void SetHexQualityMeasureToOddy()
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()
void SetHexQualityMeasureToVolume()
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()
static double QuadRelativeSizeSquared(vtkCell *cell)
static double TetVolume(vtkCell *cell)
static double QuadSkew(vtkCell *cell)
static double QuadDistortion(vtkCell *cell)
void SetTetQualityMeasureToCondition()
void SetQuadQualityMeasureToOddy()
void SetTetQualityMeasureToEdgeRatio()
int TriangleQualityMeasure
void SetHexQualityMeasureToEdgeRatio()
void SetTetQualityMeasureToScaledJacobian()
void SetQuadQualityMeasureToMinAngle()
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)
void SetTriangleQualityMeasureToShape()
void SetTriangleQualityMeasureToEdgeRatio()
static double TriangleArea(vtkCell *cell)
This is a static function used to calculate the area of a triangle.
void SetTetQualityMeasureToAspectGamma()
static double QuadMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a quadrilateral,...
void SetQuadQualityMeasureToScaledJacobian()
void SetTriangleQualityMeasureToDistortion()
void SetTriangleQualityMeasureToMinAngle()
void SetQuadQualityMeasureToStretch()
static double TetRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a tetrahedron.
void SetTetQualityMeasureToAspectBeta()
static double TetAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a tetrahedron.
static double TetShape(vtkCell *cell)
void SetQuadQualityMeasureToDistortion()
static double HexCondition(vtkCell *cell)
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
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.