VTK  9.5.2
vtkMeshQuality.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright 2003-2008 Sandia Corporation
3// SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-USGov
82
83#ifndef vtkMeshQuality_h
84#define vtkMeshQuality_h
85
86#include "vtkDataSetAlgorithm.h"
87#include "vtkFiltersVerdictModule.h" // For export macro
88
89VTK_ABI_NAMESPACE_BEGIN
90class vtkCell;
91class vtkDataArray;
92class vtkDoubleArray;
93class vtkMeshQualityFunctor;
94
95class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
96{
97private:
99
100public:
101 void PrintSelf(ostream& os, vtkIndent indent) override;
104
106
114 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
116
118
125 vtkSetMacro(LinearApproximation, bool);
126 vtkGetMacro(LinearApproximation, bool);
127 vtkBooleanMacro(LinearApproximation, bool);
129
145 {
146 AREA = 28,
147 ASPECT_FROBENIUS = 3,
148 ASPECT_GAMMA = 27,
149 ASPECT_RATIO = 1,
150 COLLAPSE_RATIO = 7,
151 CONDITION = 9,
152 DIAGONAL = 21,
153 DIMENSION = 22,
154 DISTORTION = 15,
155 EDGE_RATIO = 0,
156 EQUIANGLE_SKEW = 29,
157 EQUIVOLUME_SKEW = 30,
158 JACOBIAN = 25,
159 MAX_ANGLE = 8,
160 MAX_ASPECT_FROBENIUS = 5,
161 MAX_EDGE_RATIO = 16,
162 MAX_STRETCH = 31,
163 MEAN_ASPECT_FROBENIUS = 32,
164 MEAN_RATIO = 33,
165 MED_ASPECT_FROBENIUS = 4,
166 MIN_ANGLE = 6,
167 NODAL_JACOBIAN_RATIO = 34,
168 NORMALIZED_INRADIUS = 35,
169 ODDY = 23,
170 RADIUS_RATIO = 2,
171 RELATIVE_SIZE_SQUARED = 12,
172 SCALED_JACOBIAN = 10,
173 SHAPE = 13,
174 SHAPE_AND_SIZE = 14,
175 SHEAR = 11,
176 SHEAR_AND_SIZE = 24,
177 SKEW = 17,
178 SQUISH_INDEX = 36,
179 STRETCH = 20,
180 TAPER = 18,
181 VOLUME = 19,
182 WARPAGE = 26,
183 TOTAL_QUALITY_MEASURE_TYPES = 37,
184 NONE = TOTAL_QUALITY_MEASURE_TYPES
185 };
186
190 static const char* QualityMeasureNames[];
191
193
201 virtual void SetTriangleQualityMeasure(int measure)
202 {
203 this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
204 }
266
267
269
282 virtual void SetQuadQualityMeasure(int measure)
283 {
284 this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
285 }
365
366
368
376 virtual void SetTetQualityMeasure(int measure)
377 {
378 this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
379 }
455
456
458
464 virtual void SetPyramidQualityMeasure(int measure)
465 {
466 this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
467 }
489
490
492
499 virtual void SetWedgeQualityMeasure(int measure)
500 {
501 this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
502 }
545
546
548
557 virtual void SetHexQualityMeasure(int measure)
558 {
559 this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
560 }
629
630
634 static double TriangleArea(vtkCell* cell);
635
643 static double TriangleEdgeRatio(vtkCell* cell);
644
652 static double TriangleAspectRatio(vtkCell* cell);
653
661 static double TriangleRadiusRatio(vtkCell* cell);
662
672 static double TriangleAspectFrobenius(vtkCell* cell);
673
677 static double TriangleMinAngle(vtkCell* cell);
678
682 static double TriangleMaxAngle(vtkCell* cell);
683
687 static double TriangleCondition(vtkCell* cell);
688
692 static double TriangleScaledJacobian(vtkCell* cell);
693
701
705 static double TriangleShape(vtkCell* cell);
706
713 static double TriangleShapeAndSize(vtkCell* cell);
714
718 static double TriangleDistortion(vtkCell* cell);
719
723 static double TriangleEquiangleSkew(vtkCell* cell);
724
731
739 static double QuadEdgeRatio(vtkCell* cell);
740
748 static double QuadAspectRatio(vtkCell* cell);
749
760 static double QuadRadiusRatio(vtkCell* cell);
761
771 static double QuadMedAspectFrobenius(vtkCell* cell);
772
782 static double QuadMaxAspectFrobenius(vtkCell* cell);
783
787 static double QuadMinAngle(vtkCell* cell);
788
792 static double QuadMaxEdgeRatio(vtkCell* cell);
793
799 static double QuadSkew(vtkCell* cell);
800
805 static double QuadTaper(vtkCell* cell);
806
812 static double QuadWarpage(vtkCell* cell);
813
818 static double QuadArea(vtkCell* cell);
819
824 static double QuadStretch(vtkCell* cell);
825
829 static double QuadMaxAngle(vtkCell* cell);
830
836 static double QuadOddy(vtkCell* cell);
837
843 static double QuadCondition(vtkCell* cell);
844
850 static double QuadJacobian(vtkCell* cell);
851
857 static double QuadScaledJacobian(vtkCell* cell);
858
863 static double QuadShear(vtkCell* cell);
864
869 static double QuadShape(vtkCell* cell);
870
879 static double QuadRelativeSizeSquared(vtkCell* cell);
880
888 static double QuadShapeAndSize(vtkCell* cell);
889
897 static double QuadShearAndSize(vtkCell* cell);
898
904 static double QuadDistortion(vtkCell* cell);
905
909 static double QuadEquiangleSkew(vtkCell* cell);
910
918 static double TetEdgeRatio(vtkCell* cell);
919
927 static double TetAspectRatio(vtkCell* cell);
928
936 static double TetRadiusRatio(vtkCell* cell);
937
948 static double TetAspectFrobenius(vtkCell* cell);
949
953 static double TetMinAngle(vtkCell* cell);
954
961 static double TetCollapseRatio(vtkCell* cell);
962
968 static double TetAspectGamma(vtkCell* cell);
969
974 static double TetVolume(vtkCell* cell);
975
980 static double TetCondition(vtkCell* cell);
981
986 static double TetJacobian(vtkCell* cell);
987
993 static double TetScaledJacobian(vtkCell* cell);
994
999 static double TetShape(vtkCell* cell);
1000
1009 static double TetRelativeSizeSquared(vtkCell* cell);
1010
1018 static double TetShapeAndSize(vtkCell* cell);
1019
1025 static double TetDistortion(vtkCell* cell);
1026
1030 static double TetEquiangleSkew(vtkCell* cell);
1031
1035 static double TetEquivolumeSkew(vtkCell* cell);
1036
1042 static double TetMeanRatio(vtkCell* cell);
1043
1049 static double TetNormalizedInradius(vtkCell* cell);
1050
1054 static double TetSquishIndex(vtkCell* cell);
1055
1059 static double PyramidEquiangleSkew(vtkCell* cell);
1060
1065 static double PyramidJacobian(vtkCell* cell);
1066
1071 static double PyramidScaledJacobian(vtkCell* cell);
1072
1078 static double PyramidShape(vtkCell* cell);
1079
1083 static double PyramidVolume(vtkCell* cell);
1084
1089 static double WedgeCondition(vtkCell* cell);
1090
1095 static double WedgeDistortion(vtkCell* cell);
1096
1102 static double WedgeEdgeRatio(vtkCell* cell);
1103
1107 static double WedgeEquiangleSkew(vtkCell* cell);
1108
1113 static double WedgeJacobian(vtkCell* cell);
1114
1119 static double WedgeMaxAspectFrobenius(vtkCell* cell);
1120
1126 static double WedgeMaxStretch(vtkCell* cell);
1127
1133 static double WedgeMeanAspectFrobenius(vtkCell* cell);
1134
1144 static double WedgeScaledJacobian(vtkCell* cell);
1145
1151 static double WedgeShape(vtkCell* cell);
1152
1156 static double WedgeVolume(vtkCell* cell);
1157
1165 static double HexEdgeRatio(vtkCell* cell);
1166
1171 static double HexMedAspectFrobenius(vtkCell* cell);
1172
1177 static double HexMaxAspectFrobenius(vtkCell* cell);
1178
1182 static double HexMaxEdgeRatio(vtkCell* cell);
1183
1189 static double HexSkew(vtkCell* cell);
1190
1195 static double HexTaper(vtkCell* cell);
1196
1201 static double HexVolume(vtkCell* cell);
1202
1207 static double HexStretch(vtkCell* cell);
1208
1213 static double HexDiagonal(vtkCell* cell);
1214
1220 static double HexDimension(vtkCell* cell);
1221
1227 static double HexOddy(vtkCell* cell);
1228
1233 static double HexCondition(vtkCell* cell);
1234
1240 static double HexJacobian(vtkCell* cell);
1241
1247 static double HexScaledJacobian(vtkCell* cell);
1248
1253 static double HexShear(vtkCell* cell);
1254
1259 static double HexShape(vtkCell* cell);
1260
1269 static double HexRelativeSizeSquared(vtkCell* cell);
1270
1278 static double HexShapeAndSize(vtkCell* cell);
1279
1287 static double HexShearAndSize(vtkCell* cell);
1288
1294 static double HexDistortion(vtkCell* cell);
1295
1299 static double HexEquiangleSkew(vtkCell* cell);
1300
1305 static double HexNodalJacobianRatio(vtkCell* cell);
1306
1317 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1319 vtkBooleanMacro(Ratio, vtkTypeBool);
1320
1321protected:
1323 ~vtkMeshQuality() override = default;
1324
1326
1335
1336 using CellQualityType = double (*)(vtkCell*);
1343
1344 // Variables used to store the average size (2D: area / 3D: volume)
1345 static double TriangleAverageSize;
1346 static double QuadAverageSize;
1347 static double TetAverageSize;
1348 static double PyramidAverageSize;
1349 static double WedgeAverageSize;
1350 static double HexAverageSize;
1351
1352private:
1353 vtkMeshQuality(const vtkMeshQuality&) = delete;
1354 void operator=(const vtkMeshQuality&) = delete;
1355};
1356
1357VTK_ABI_NAMESPACE_END
1358#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition vtkCell.h:130
dynamic, self-adjusting array of double
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
void SetQuadQualityMeasureToSkew()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
virtual void SetWedgeQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of wedges.
static double QuadWarpage(vtkCell *cell)
Calculate the warpage of a quadrilateral.
void SetHexQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTetQualityMeasureToEquivolumeSkew()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadTaper(vtkCell *cell)
Calculate the taper of a quadrilateral.
virtual void SetTetQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetWedgeQualityMeasureToMeanAspectFrobenius()
Set/Get the particular estimator used to measure the quality of wedges.
void SetTriangleQualityMeasureToArea()
Set/Get the particular estimator used to function the quality of triangles.
vtkGetEnumMacro(HexQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetWedgeQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of wedges.
static double TetNormalizedInradius(vtkCell *cell)
Calculate the normalized in-radius of a tetrahedron.
static double HexOddy(vtkCell *cell)
Calculate the oddy of a hexahedron.
void SetWedgeQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of wedges.
virtual void SetTriangleQualityMeasure(int measure)
Set/Get the particular estimator used to function the quality of triangles.
static double WedgeScaledJacobian(vtkCell *cell)
Calculate the scaled Jacobian a wedge.
static double TetAspectGamma(vtkCell *cell)
Calculate the aspect gamma of a tetrahedron.
void SetPyramidQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of pyramids.
static double WedgeMaxStretch(vtkCell *cell)
Calculate the max stretch of a wedge.
void SetPyramidQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of pyramids.
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)
Calculate the shear of a hexahedron.
static double HexScaledJacobian(vtkCell *cell)
Calculate the scaled Jacobian of a hexahedron.
virtual void SetPyramidQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of pyramids.
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.
void SetWedgeQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of wedges.
static double QuadAverageSize
static double QuadJacobian(vtkCell *cell)
Calculate the Jacobian of a quadrilateral.
void SetTriangleQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to function the quality of triangles.
static double PyramidJacobian(vtkCell *cell)
Calculate the Jacobian of a pyramid.
QualityMeasureTypes TriangleQualityMeasure
QualityMeasureTypes
Enum which lists the Quality Measures Types.
void SetQuadQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleNormalizedInradius(vtkCell *cell)
Calculate the normalized in-radius of a triangle.
static double WedgeShape(vtkCell *cell)
Calculate the shape of a wedge.
static double WedgeEdgeRatio(vtkCell *cell)
Calculate the edge ratio of a wedge.
static double TriangleEdgeRatio(vtkCell *cell)
Calculate the edge ratio of a triangle.
void SetHexQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double WedgeMeanAspectFrobenius(vtkCell *cell)
Calculate the mean aspect Frobenius of a wedge.
static double PyramidVolume(vtkCell *cell)
Calculate the volume of a pyramid.
static double QuadOddy(vtkCell *cell)
Calculate the oddy of a quadrilateral.
static double QuadAspectRatio(vtkCell *cell)
Calculate the aspect ratio of a planar quadrilateral.
static double TriangleAspectFrobenius(vtkCell *cell)
Calculate the Frobenius condition number of the transformation matrix from an equilateral triangle to...
vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to function the quality of triangles.
static double TetMeanRatio(vtkCell *cell)
Calculate the mean ratio of a tetrahedron.
static double QuadEquiangleSkew(vtkCell *cell)
Calculate the equiangle skew of a quadrilateral.
static double QuadScaledJacobian(vtkCell *cell)
Calculate the scaled Jacobian of a quadrilateral.
vtkSetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of wedges.
static double HexMaxEdgeRatio(vtkCell *cell)
Calculate the maximum edge ratio of a hexahedron at its center.
void SetHexQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexShapeAndSize(vtkCell *cell)
Calculate the shape and size of a hexahedron.
static double QuadShear(vtkCell *cell)
Calculate the shear of a quadrilateral.
void SetHexQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of hexahedra.
vtkGetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToAspectRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleShapeAndSize(vtkCell *cell)
Calculate the product of shape and relative size of a triangle.
vtkTypeBool GetRatio()
static double HexTaper(vtkCell *cell)
Calculate the taper of a hexahedron.
static vtkMeshQuality * New()
void SetTriangleQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to function the quality of triangles.
CellQualityType GetTriangleQualityMeasureFunction()
virtual void SetQuadQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTetQualityMeasureToSquishIndex()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double TetJacobian(vtkCell *cell)
Calculate the Jacobian of a tetrahedron.
static double HexDistortion(vtkCell *cell)
Calculate the distortion of a hexahedron.
static double TetCollapseRatio(vtkCell *cell)
Calculate the collapse ratio of a tetrahedron.
static double TriangleDistortion(vtkCell *cell)
Calculate the distortion of a triangle.
static double HexVolume(vtkCell *cell)
Calculate the volume of a hexahedron.
void SetQuadQualityMeasureToMaxEdgeRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of hexahedra.
vtkGetEnumMacro(QuadQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double QuadArea(vtkCell *cell)
Calculate the area of a quadrilateral.
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.
void SetTetQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double TetEdgeRatio(vtkCell *cell)
Calculate the edge ratio of a tetrahedron.
void SetQuadQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTetQualityMeasureToNormalizedInradius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double HexRelativeSizeSquared(vtkCell *cell)
Calculate the relative size squared of a hexahedron.
void SetTetQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of tetrahedra.
vtkSetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of pyramids.
void SetQuadQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double PyramidShape(vtkCell *cell)
Calculate the shape of a pyramid.
QualityMeasureTypes QuadQualityMeasure
CellQualityType GetTetQualityMeasureFunction()
static double HexShape(vtkCell *cell)
Calculate the shape of a hexahedron.
void SetHexQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetAverageSize
static double TriangleScaledJacobian(vtkCell *cell)
Calculate the scaled Jacobian of a triangle.
void SetHexQualityMeasureToDimension()
Set/Get the particular estimator used to measure the quality of hexahedra.
virtual vtkTypeBool GetSaveCellQuality()
This variable controls whether or not cell quality is stored as cell data in the resulting mesh or di...
static double TetRelativeSizeSquared(vtkCell *cell)
Calculate the relative size squared of a tetrahedron.
void SetTriangleQualityMeasureToNormalizedInradius()
Set/Get the particular estimator used to function the quality of triangles.
static double TriangleEquiangleSkew(vtkCell *cell)
Calculate the equiangle skew of a triangle.
QualityMeasureTypes TetQualityMeasure
static double TriangleShape(vtkCell *cell)
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)
Calculate the average Frobenius aspect of the 4 corner triangles of a planar quadrilateral,...
static double TriangleCondition(vtkCell *cell)
Calculate the condition number of a triangle.
static double WedgeDistortion(vtkCell *cell)
Calculate the distortion of a wedge.
static double HexEquiangleSkew(vtkCell *cell)
Calculate the equiangle skew of a hexahedron.
void SetQuadQualityMeasureToTaper()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetWedgeQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of wedges.
CellQualityType GetWedgeQualityMeasureFunction()
static double TetMinAngle(vtkCell *cell)
Calculate the minimal (nonoriented) dihedral angle of a tetrahedron, expressed in degrees.
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.
static double WedgeEquiangleSkew(vtkCell *cell)
Calculate the equiangle skew of a wedge.
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)
Calculate the condition number of a quadrilateral.
void SetPyramidQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of pyramids.
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.
static double HexAverageSize
vtkGetEnumMacro(TetQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetQuadQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToMaxEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TriangleAspectRatio(vtkCell *cell)
Calculate the aspect ratio of a triangle.
static double TetSquishIndex(vtkCell *cell)
Calculate the squish index of a tetrahedron.
static double HexJacobian(vtkCell *cell)
Calculate the Jacobian of a hexahedron.
void SetQuadQualityMeasureToMedAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
vtkGetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of pyramids.
static double PyramidEquiangleSkew(vtkCell *cell)
Calculate the equiangle skew of a pyramid.
static double TetDistortion(vtkCell *cell)
Calculate the distortion of a tetrahedron.
static double QuadShapeAndSize(vtkCell *cell)
Calculate the shape and size of a quadrilateral.
void SetWedgeQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of wedges.
static double TriangleMaxAngle(vtkCell *cell)
Calculate the maximal (nonoriented) angle of a triangle, expressed in degrees.
void SetTriangleQualityMeasureToCondition()
Set/Get the particular estimator used to function the quality of triangles.
static double TetCondition(vtkCell *cell)
Calculate the condition number of a tetrahedron.
void SetWedgeQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of wedges.
void SetTetQualityMeasureToMeanRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double TriangleRadiusRatio(vtkCell *cell)
Calculate the radius ratio of a triangle.
void SetWedgeQualityMeasureToMaxStretch()
Set/Get the particular estimator used to measure the quality of wedges.
static double TetScaledJacobian(vtkCell *cell)
Calculate the scaled Jacobian of a tetrahedron.
static double QuadEdgeRatio(vtkCell *cell)
Calculate the edge ratio of a quadrilateral.
~vtkMeshQuality() override=default
void SetTetQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadShearAndSize(vtkCell *cell)
Calculate the shear and size of a quadrilateral.
static double HexDimension(vtkCell *cell)
Calculate the dimension of a hexahedron.
friend class vtkMeshQualityFunctor
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)
Calculate the stretch of a hexahedron.
QualityMeasureTypes HexQualityMeasure
static double QuadMaxAngle(vtkCell *cell)
Calculate the maximum (nonoriented) angle of a quadrilateral, expressed in degrees.
void SetHexQualityMeasureToNodalJacobianRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadStretch(vtkCell *cell)
Calculate the stretch of a quadrilateral.
void SetWedgeQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of wedges.
void SetTriangleQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of hexahedra.
double(*)(vtkCell *) CellQualityType
static double HexEdgeRatio(vtkCell *cell)
Calculate the edge ratio of a hexahedron.
void SetTetQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of tetrahedra.
CellQualityType GetQuadQualityMeasureFunction()
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)
Calculate the maximal Frobenius aspect of the 4 corner triangles of a planar quadrilateral,...
vtkTypeBool SaveCellQuality
void SetHexQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetAspectFrobenius(vtkCell *cell)
Calculate the Frobenius condition number of the transformation matrix from a regular tetrahedron to a...
static double HexMaxAspectFrobenius(vtkCell *cell)
Calculate the maximal Frobenius aspect of the 8 corner tetrahedra of a hexahedron,...
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 TriangleAverageSize
static double HexShearAndSize(vtkCell *cell)
Calculate the shear and size of a hexahedron.
static double HexSkew(vtkCell *cell)
Calculate the skew of a hexahedron.
void SetTetQualityMeasureToCollapseRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double TetEquiangleSkew(vtkCell *cell)
Calculate the equiangle skew of a tetrahedron.
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.
static const char * QualityMeasureNames[]
Array which lists the Quality Measures Names.
void SetTetQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadShape(vtkCell *cell)
Calculate the shear of a quadrilateral.
vtkSetEnumMacro(TetQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of tetrahedra.
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)
Calculate the average Frobenius aspect of the 8 corner tetrahedra of a hexahedron,...
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)
Calculate the square of the relative size of a triangle.
static double QuadRadiusRatio(vtkCell *cell)
Calculate the radius ratio of a planar quadrilateral.
static double QuadRelativeSizeSquared(vtkCell *cell)
Calculate the relative size squared of a quadrilateral.
virtual void SetSaveCellQuality(vtkTypeBool)
This variable controls whether or not cell quality is stored as cell data in the resulting mesh or di...
static double TetEquivolumeSkew(vtkCell *cell)
Calculate the equivolume skew of a tetrahedron.
static double TetVolume(vtkCell *cell)
Calculate the volume of a tetrahedron.
static double QuadSkew(vtkCell *cell)
Calculate the skew of a quadrilateral.
static double PyramidScaledJacobian(vtkCell *cell)
Calculate the Jacobian of a pyramid.
static double PyramidAverageSize
void SetWedgeQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of wedges.
static double QuadDistortion(vtkCell *cell)
Calculate the distortion of a quadrilateral.
void SetPyramidQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of pyramids.
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.
static double WedgeJacobian(vtkCell *cell)
Calculate the Jacobian of a wedge.
QualityMeasureTypes PyramidQualityMeasure
static double QuadMaxEdgeRatio(vtkCell *cell)
Calculate the maximum edge length ratio of a quadrilateral at quad center.
void SetHexQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
CellQualityType GetPyramidQualityMeasureFunction()
static double WedgeCondition(vtkCell *cell)
Calculate the condition number of a wedge.
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 WedgeVolume(vtkCell *cell)
Calculate the volume of a wedge.
vtkSetEnumMacro(QuadQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleMinAngle(vtkCell *cell)
Calculate the minimal (nonoriented) angle of a triangle, expressed in degrees.
vtkSetEnumMacro(HexQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexDiagonal(vtkCell *cell)
Calculate the diagonal of a hexahedron.
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)
Calculate the area of a triangle.
CellQualityType GetHexQualityMeasureFunction()
static double TetShapeAndSize(vtkCell *cell)
Calculate the shape and size of a tetrahedron.
vtkGetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to function the quality of triangles.
void SetTetQualityMeasureToAspectGamma()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadMinAngle(vtkCell *cell)
Calculate the minimal (nonoriented) angle of a quadrilateral, expressed in degrees.
void SetQuadQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
virtual void SetHexQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of hexahedra.
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)
Calculate the radius ratio of a tetrahedron.
static double HexNodalJacobianRatio(vtkCell *cell)
Calculate the nodal Jacobian ratio of a hexahedron.
static double WedgeAverageSize
static double TetAspectRatio(vtkCell *cell)
Calculate the aspect ratio of a tetrahedron.
void SetWedgeQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of wedges.
QualityMeasureTypes WedgeQualityMeasure
static double WedgeMaxAspectFrobenius(vtkCell *cell)
Calculate the max aspect Frobenius of a wedge.
static double TetShape(vtkCell *cell)
Calculate the shape of a tetrahedron.
void SetPyramidQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of pyramids.
void SetQuadQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexCondition(vtkCell *cell)
Calculate the condition of a hexahedron.
int vtkTypeBool
Definition vtkABI.h:64
#define vtkDataArray