VTK  9.6.1
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
53
54#ifndef vtkMeshQuality_h
55#define vtkMeshQuality_h
56
57#include "vtkDataSetAlgorithm.h"
58#include "vtkDeprecation.h" // For deprecation macro
59#include "vtkFiltersVerdictModule.h" // For export macro
60
61#include <functional> // For std::function
62
63VTK_ABI_NAMESPACE_BEGIN
64class vtkCell;
65class vtkDataArray;
66class vtkDoubleArray;
67
68class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
69{
70private:
71 class vtkMeshQualityFunctor;
72
73public:
74 void PrintSelf(ostream& os, vtkIndent indent) override;
77
79
87 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
89
91
98 vtkSetMacro(LinearApproximation, bool);
99 vtkGetMacro(LinearApproximation, bool);
100 vtkBooleanMacro(LinearApproximation, bool);
102
108 static void ComputeAverageCellSize(vtkDataSet* dataset);
109
163 {
164 AREA = 28,
165 ASPECT_FROBENIUS = 3,
166 ASPECT_GAMMA = 27,
167 ASPECT_RATIO = 1,
168 COLLAPSE_RATIO = 7,
169 CONDITION = 9,
170 DIAGONAL = 21,
171 DIMENSION = 22,
172 DISTORTION = 15,
173 EDGE_RATIO = 0,
174 EQUIANGLE_SKEW = 29,
175 EQUIVOLUME_SKEW = 30,
176 JACOBIAN = 25,
177 INRADIUS = 37,
178 MAX_ANGLE = 8,
179 MAX_ASPECT_FROBENIUS = 5,
180 MAX_EDGE_RATIO = 16,
181 MAX_STRETCH = 31,
182 MEAN_ASPECT_FROBENIUS = 32,
183 MEAN_RATIO = 33,
184 MED_ASPECT_FROBENIUS = 4,
185 MIN_ANGLE = 6,
186 NODAL_JACOBIAN_RATIO = 34,
187 NORMALIZED_INRADIUS = 35,
188 ODDY = 23,
189 RADIUS_RATIO = 2,
190 RELATIVE_SIZE_SQUARED = 12,
191 SCALED_JACOBIAN = 10,
192 SHAPE = 13,
193 SHAPE_AND_SIZE = 14,
194 SHEAR = 11,
195 SHEAR_AND_SIZE = 24,
196 SKEW = 17,
197 SQUISH_INDEX = 36,
198 STRETCH = 20,
199 TAPER = 18,
200 VOLUME = 19,
201 WARPAGE = 26,
202 TOTAL_QUALITY_MEASURE_TYPES = 38,
203 NONE = TOTAL_QUALITY_MEASURE_TYPES
204 };
205
209 static const char* QualityMeasureNames[];
210
212
220 virtual void SetTriangleQualityMeasure(int measure)
221 {
222 this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
223 }
285
286
288
301 virtual void SetQuadQualityMeasure(int measure)
302 {
303 this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
304 }
384
385
387
395 virtual void SetTetQualityMeasure(int measure)
396 {
397 this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
398 }
478
479
481
487 virtual void SetPyramidQualityMeasure(int measure)
488 {
489 this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
490 }
512
513
515
522 virtual void SetWedgeQualityMeasure(int measure)
523 {
524 this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
525 }
568
569
571
580 virtual void SetHexQualityMeasure(int measure)
581 {
582 this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
583 }
652
653
659 static double TriangleArea(vtkCell* cell, bool linearApproximation = false);
660
672 static double TriangleAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
673
683 static double TriangleAspectRatio(vtkCell* cell, bool linearApproximation = false);
684
690 static double TriangleCondition(vtkCell* cell, bool linearApproximation = false);
691
698 static double TriangleDistortion(vtkCell* cell, bool linearApproximation = false);
699
709 static double TriangleEdgeRatio(vtkCell* cell, bool linearApproximation = false);
710
716 static double TriangleEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
717
723 static double TriangleMaxAngle(vtkCell* cell, bool linearApproximation = false);
724
730 static double TriangleMinAngle(vtkCell* cell, bool linearApproximation = false);
731
739 static double TriangleNormalizedInradius(vtkCell* cell, bool linearApproximation = false);
740
750 static double TriangleRadiusRatio(vtkCell* cell, bool linearApproximation = false);
751
761 static double TriangleRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
762
768 static double TriangleScaledJacobian(vtkCell* cell, bool linearApproximation = false);
769
775 static double TriangleShape(vtkCell* cell, bool linearApproximation = false);
776
786 static double TriangleShapeAndSize(vtkCell* cell, bool linearApproximation = false);
787
794 static double QuadArea(vtkCell* cell, bool linearApproximation = false);
795
805 static double QuadAspectRatio(vtkCell* cell, bool linearApproximation = false);
806
814 static double QuadCondition(vtkCell* cell, bool linearApproximation = false);
815
823 static double QuadDistortion(vtkCell* cell, bool linearApproximation = false);
824
834 static double QuadEdgeRatio(vtkCell* cell, bool linearApproximation = false);
835
841 static double QuadEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
842
850 static double QuadJacobian(vtkCell* cell, bool linearApproximation = false);
851
857 static double QuadMaxAngle(vtkCell* cell, bool linearApproximation = false);
858
870 static double QuadMaxAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
871
877 static double QuadMaxEdgeRatio(vtkCell* cell, bool linearApproximation = false);
878
890 static double QuadMedAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
891
897 static double QuadMinAngle(vtkCell* cell, bool linearApproximation = false);
898
906 static double QuadOddy(vtkCell* cell, bool linearApproximation = false);
907
920 static double QuadRadiusRatio(vtkCell* cell, bool linearApproximation = false);
921
933 static double QuadRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
934
942 static double QuadScaledJacobian(vtkCell* cell, bool linearApproximation = false);
943
950 static double QuadShape(vtkCell* cell, bool linearApproximation = false);
951
962 static double QuadShapeAndSize(vtkCell* cell, bool linearApproximation = false);
963
970 static double QuadShear(vtkCell* cell, bool linearApproximation = false);
971
982 static double QuadShearAndSize(vtkCell* cell, bool linearApproximation = false);
983
991 static double QuadSkew(vtkCell* cell, bool linearApproximation = false);
992
999 static double QuadStretch(vtkCell* cell, bool linearApproximation = false);
1000
1007 static double QuadTaper(vtkCell* cell, bool linearApproximation = false);
1008
1016 static double QuadWarpage(vtkCell* cell, bool linearApproximation = false);
1017
1030 static double TetAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1031
1039 static double TetAspectGamma(vtkCell* cell, bool linearApproximation = false);
1040
1050 static double TetAspectRatio(vtkCell* cell, bool linearApproximation = false);
1051
1060 static double TetCollapseRatio(vtkCell* cell, bool linearApproximation = false);
1061
1068 static double TetCondition(vtkCell* cell, bool linearApproximation = false);
1069
1078 static double TetDistortion(vtkCell* cell, bool linearApproximation = false);
1079
1089 static double TetEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1090
1096 static double TetEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1097
1103 static double TetEquivolumeSkew(vtkCell* cell, bool linearApproximation = false);
1104
1111 static double TetInradius(vtkCell* cell, bool linearApproximation = false);
1112
1119 static double TetJacobian(vtkCell* cell, bool linearApproximation = false);
1120
1128 static double TetMeanRatio(vtkCell* cell, bool linearApproximation = false);
1129
1135 static double TetMinAngle(vtkCell* cell, bool linearApproximation = false);
1136
1144 static double TetNormalizedInradius(vtkCell* cell, bool linearApproximation = false);
1145
1155 static double TetRadiusRatio(vtkCell* cell, bool linearApproximation = false);
1156
1168 static double TetRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
1169
1177 static double TetScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1178
1185 static double TetShape(vtkCell* cell, bool linearApproximation = false);
1186
1197 static double TetShapeAndSize(vtkCell* cell, bool linearApproximation = false);
1198
1204 static double TetSquishIndex(vtkCell* cell, bool linearApproximation = false);
1205
1212 static double TetVolume(vtkCell* cell, bool linearApproximation = false);
1213
1219 static double PyramidEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1220
1227 static double PyramidJacobian(vtkCell* cell, bool linearApproximation = false);
1228
1235 static double PyramidScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1236
1244 static double PyramidShape(vtkCell* cell, bool linearApproximation = false);
1245
1251 static double PyramidVolume(vtkCell* cell, bool linearApproximation = false);
1252
1259 static double WedgeCondition(vtkCell* cell, bool linearApproximation = false);
1260
1267 static double WedgeDistortion(vtkCell* cell, bool linearApproximation = false);
1268
1276 static double WedgeEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1277
1283 static double WedgeEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1284
1291 static double WedgeJacobian(vtkCell* cell, bool linearApproximation = false);
1292
1299 static double WedgeMaxAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1300
1308 static double WedgeMaxStretch(vtkCell* cell, bool linearApproximation = false);
1309
1317 static double WedgeMeanAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1318
1330 static double WedgeScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1331
1339 static double WedgeShape(vtkCell* cell, bool linearApproximation = false);
1340
1346 static double WedgeVolume(vtkCell* cell, bool linearApproximation = false);
1347
1354 static double HexCondition(vtkCell* cell, bool linearApproximation = false);
1355
1362 static double HexDiagonal(vtkCell* cell, bool linearApproximation = false);
1363
1371 static double HexDimension(vtkCell* cell, bool linearApproximation = false);
1372
1380 static double HexDistortion(vtkCell* cell, bool linearApproximation = false);
1381
1391 static double HexEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1392
1398 static double HexEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1399
1407 static double HexJacobian(vtkCell* cell, bool linearApproximation = false);
1408
1415 static double HexMaxAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1416
1422 static double HexMaxEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1423
1430 static double HexMedAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1431
1438 static double HexNodalJacobianRatio(vtkCell* cell, bool linearApproximation = false);
1439
1447 static double HexOddy(vtkCell* cell, bool linearApproximation = false);
1448
1460 static double HexRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
1461
1469 static double HexScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1470
1477 static double HexShape(vtkCell* cell, bool linearApproximation = false);
1478
1489 static double HexShapeAndSize(vtkCell* cell, bool linearApproximation = false);
1490
1497 static double HexShear(vtkCell* cell, bool linearApproximation = false);
1498
1508 static double HexShearAndSize(vtkCell* cell, bool linearApproximation = false);
1509
1518 static double HexSkew(vtkCell* cell, bool linearApproximation = false);
1519
1526 static double HexStretch(vtkCell* cell, bool linearApproximation = false);
1527
1534 static double HexTaper(vtkCell* cell, bool linearApproximation = false);
1541 static double HexVolume(vtkCell* cell, bool linearApproximation = false);
1542
1553 VTK_DEPRECATED_IN_9_6_0("Please use SetSaveCellQuality instead.")
1554 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1555 VTK_DEPRECATED_IN_9_6_0("Please use GetSaveCellQuality instead.")
1557 VTK_DEPRECATED_IN_9_6_0("Please use SaveCellQualityOn instead.")
1558 void RatioOn() { this->SaveCellQualityOn(); }
1559 VTK_DEPRECATED_IN_9_6_0("Please use SaveCellQualityOff instead.")
1560 void RatioOff() { this->SaveCellQualityOff(); }
1561
1562protected:
1564 ~vtkMeshQuality() override = default;
1565
1567
1576
1577 using CellQualityType = std::function<double(vtkCell*, bool)>;
1591
1592 // Variables used to store the average size (2D: area / 3D: volume)
1593 static double TriangleAverageSize;
1594 static double QuadAverageSize;
1595 static double TetAverageSize;
1596 static double PyramidAverageSize;
1597 static double WedgeAverageSize;
1598 static double HexAverageSize;
1599
1600private:
1601 vtkMeshQuality(const vtkMeshQuality&) = delete;
1602 void operator=(const vtkMeshQuality&) = delete;
1603};
1604
1605VTK_ABI_NAMESPACE_END
1606#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition vtkCell.h:50
abstract class to specify dataset behavior
Definition vtkDataSet.h:57
dynamic, self-adjusting array of double
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
static double QuadScaledJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the scaled Jacobian of a quadrilateral.
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.
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.
CellQualityType GetQuadraticHexQualityMeasureFunction() const
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.
static double QuadDistortion(vtkCell *cell, bool linearApproximation=false)
Calculate the distortion of a quadrilateral.
CellQualityType GetHexQualityMeasureFunction() const
void SetTriangleQualityMeasureToArea()
Set/Get the particular estimator used to function the quality of triangles.
CellQualityType GetPyramidQualityMeasureFunction() const
vtkGetEnumMacro(HexQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetShape(vtkCell *cell, bool linearApproximation=false)
Calculate the shape of a tetrahedron.
static double HexOddy(vtkCell *cell, bool linearApproximation=false)
Calculate the oddy of a hexahedron.
void SetWedgeQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of wedges.
static double TriangleAspectRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the aspect ratio of a triangle.
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 TetScaledJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the scaled Jacobian of a tetrahedron.
void SetPyramidQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of pyramids.
static double QuadAspectRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the aspect ratio of a planar quadrilateral.
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.
CellQualityType GetTetQualityMeasureFunction() const
virtual void SetPyramidQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of pyramids.
static double TriangleRadiusRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the radius ratio of a triangle.
static double TetMeanRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the mean ratio of a tetrahedron.
static double WedgeMaxAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the max aspect Frobenius of a wedge.
static double TetEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the edge ratio of a tetrahedron.
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.
static double TetJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the Jacobian of a tetrahedron.
static double WedgeJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the Jacobian of a wedge.
static double TetRelativeSizeSquared(vtkCell *cell, bool linearApproximation=false)
Calculate the relative size squared of a tetrahedron.
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 HexJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the Jacobian of a hexahedron.
void SetWedgeQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of wedges.
static double QuadAverageSize
CellQualityType GetQuadraticQuadQualityMeasureFunction() const
void SetTriangleQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to function the quality of triangles.
void SetTetQualityMeasureToInradius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
QualityMeasureTypes TriangleQualityMeasure
QualityMeasureTypes
Enum which lists the Quality Measures Types.
CellQualityType GetQuadraticTriangleQualityMeasureFunction() const
static double TetAspectRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the aspect ratio of a tetrahedron.
void SetQuadQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TetAspectGamma(vtkCell *cell, bool linearApproximation=false)
Calculate the aspect gamma of a tetrahedron.
static double HexStretch(vtkCell *cell, bool linearApproximation=false)
Calculate the stretch of a hexahedron.
void SetHexQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TriangleCondition(vtkCell *cell, bool linearApproximation=false)
Calculate the condition number of a triangle.
static double QuadShapeAndSize(vtkCell *cell, bool linearApproximation=false)
Calculate the shape and size of a quadrilateral.
vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to function the quality of triangles.
vtkSetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of wedges.
static double TriangleEquiangleSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equiangle skew of a triangle.
void SetHexQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
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 TriangleDistortion(vtkCell *cell, bool linearApproximation=false)
Calculate the distortion of a triangle.
static double QuadTaper(vtkCell *cell, bool linearApproximation=false)
Calculate the taper of a quadrilateral.
vtkTypeBool GetRatio()
static double TetCollapseRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the collapse ratio of a tetrahedron.
static double TriangleMinAngle(vtkCell *cell, bool linearApproximation=false)
Calculate the minimal (nonoriented) angle of a triangle, expressed in degrees.
static vtkMeshQuality * New()
void SetTriangleQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to function the quality of triangles.
static double QuadSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the skew of a quadrilateral.
virtual void SetQuadQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the edge ratio of a hexahedron.
void SetTetQualityMeasureToSquishIndex()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadRadiusRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the radius ratio of a planar quadrilateral.
static double PyramidVolume(vtkCell *cell, bool linearApproximation=false)
Calculate the volume of a pyramid.
static double HexTaper(vtkCell *cell, bool linearApproximation=false)
Calculate the taper 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 TetSquishIndex(vtkCell *cell, bool linearApproximation=false)
Calculate the squish index of a tetrahedron.
void SetTriangleQualityMeasureToAspectRatio()
Set/Get the particular estimator used to function the quality of triangles.
static double TriangleScaledJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the scaled Jacobian of a triangle.
CellQualityType GetQuadraticTetQualityMeasureFunction() const
void SetHexQualityMeasureToStretch()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetInradius(vtkCell *cell, bool linearApproximation=false)
Calculate the inradius of a tetrahedron.
static double TetAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the Frobenius condition number of the transformation matrix from a regular tetrahedron to a...
void SetTetQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetQuadQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
CellQualityType GetTriQuadraticHexQualityMeasureFunction() const
void SetTetQualityMeasureToNormalizedInradius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadStretch(vtkCell *cell, bool linearApproximation=false)
Calculate the stretch of a quadrilateral.
static double HexMaxAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the maximal Frobenius aspect of the 8 corner tetrahedra of a hexahedron,...
void SetTetQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double HexEquiangleSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equiangle skew of a hexahedron.
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 QuadMaxAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the maximal Frobenius aspect of the 4 corner triangles of a planar quadrilateral,...
QualityMeasureTypes QuadQualityMeasure
CellQualityType GetTriangleQualityMeasureFunction() const
static double HexSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the skew of a hexahedron.
static double PyramidEquiangleSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equiangle skew of a pyramid.
virtual void SaveCellQualityOff()
This variable controls whether or not cell quality is stored as cell data in the resulting mesh or di...
static double QuadEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the edge ratio of a quadrilateral.
static double HexNodalJacobianRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the nodal Jacobian ratio of a hexahedron.
void SetHexQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetAverageSize
static double HexShape(vtkCell *cell, bool linearApproximation=false)
Calculate the shape of a hexahedron.
static double TetCondition(vtkCell *cell, bool linearApproximation=false)
Calculate the condition number of a tetrahedron.
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 WedgeMeanAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the mean aspect Frobenius of a wedge.
void SetTriangleQualityMeasureToNormalizedInradius()
Set/Get the particular estimator used to function the quality of triangles.
static double TriangleRelativeSizeSquared(vtkCell *cell, bool linearApproximation=false)
Calculate the square of the relative size of a triangle.
static double WedgeDistortion(vtkCell *cell, bool linearApproximation=false)
Calculate the distortion of a wedge.
QualityMeasureTypes TetQualityMeasure
void SetTetQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadMinAngle(vtkCell *cell, bool linearApproximation=false)
Calculate the minimal (nonoriented) angle of a quadrilateral, expressed in degrees.
static double QuadRelativeSizeSquared(vtkCell *cell, bool linearApproximation=false)
Calculate the relative size squared of a quadrilateral.
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.
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 QuadShear(vtkCell *cell, bool linearApproximation=false)
Calculate the shear of a quadrilateral.
void SetHexQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double PyramidJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the Jacobian of a pyramid.
void SetQuadQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
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.
static double HexCondition(vtkCell *cell, bool linearApproximation=false)
Calculate the condition of a hexahedron.
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
static double HexShapeAndSize(vtkCell *cell, bool linearApproximation=false)
Calculate the shape and size of a hexahedron.
static double HexRelativeSizeSquared(vtkCell *cell, bool linearApproximation=false)
Calculate the relative size squared of a hexahedron.
static double WedgeVolume(vtkCell *cell, bool linearApproximation=false)
Calculate the volume of a wedge.
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.
CellQualityType GetBiQuadraticQuadQualityMeasureFunction() const
static double HexShearAndSize(vtkCell *cell, bool linearApproximation=false)
Calculate the shear and size 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.
CellQualityType GetBiQuadraticTriangleQualityMeasureFunction() const
static double TetDistortion(vtkCell *cell, bool linearApproximation=false)
Calculate the distortion of a tetrahedron.
void SetWedgeQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of wedges.
void SetTriangleQualityMeasureToCondition()
Set/Get the particular estimator used to function the quality of triangles.
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.
void SetWedgeQualityMeasureToMaxStretch()
Set/Get the particular estimator used to measure the quality of wedges.
~vtkMeshQuality() override=default
void SetTetQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double TriangleShapeAndSize(vtkCell *cell, bool linearApproximation=false)
Calculate the product of shape and relative size of a triangle.
static double WedgeShape(vtkCell *cell, bool linearApproximation=false)
Calculate the shape of a wedge.
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 TriangleArea(vtkCell *cell, bool linearApproximation=false)
Calculate the area of a triangle.
CellQualityType GetWedgeQualityMeasureFunction() const
QualityMeasureTypes HexQualityMeasure
static double HexDistortion(vtkCell *cell, bool linearApproximation=false)
Calculate the distortion of a hexahedron.
void SetHexQualityMeasureToNodalJacobianRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TriangleNormalizedInradius(vtkCell *cell, bool linearApproximation=false)
Calculate the normalized in-radius of a triangle.
static double WedgeEquiangleSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equiangle skew of a wedge.
static double TriangleMaxAngle(vtkCell *cell, bool linearApproximation=false)
Calculate the maximal (nonoriented) angle of a triangle, expressed in degrees.
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.
static void ComputeAverageCellSize(vtkDataSet *dataset)
Compute the average cell size for each cell type in the mesh, and store it in field data arrays (TriA...
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 TetRadiusRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the radius ratio of a tetrahedron.
vtkTypeBool SaveCellQuality
void SetHexQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TetEquivolumeSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equivolume skew of a tetrahedron.
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 QuadWarpage(vtkCell *cell, bool linearApproximation=false)
Calculate the warpage of a quadrilateral.
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.
static double TetEquiangleSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equiangle skew of a tetrahedron.
void SetQuadQualityMeasureToMaxAngle()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double WedgeMaxStretch(vtkCell *cell, bool linearApproximation=false)
Calculate the max stretch of a wedge.
static const char * QualityMeasureNames[]
Array which lists the Quality Measures Names.
static double PyramidShape(vtkCell *cell, bool linearApproximation=false)
Calculate the shape of a pyramid.
static double TetVolume(vtkCell *cell, bool linearApproximation=false)
Calculate the volume of a tetrahedron.
void SetTetQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double TriangleAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the Frobenius condition number of the transformation matrix from an equilateral triangle to...
vtkSetEnumMacro(TetQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadOddy(vtkCell *cell, bool linearApproximation=false)
Calculate the oddy of a quadrilateral.
static double WedgeScaledJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the scaled Jacobian a wedge.
void SetTetQualityMeasureToRadiusRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadMedAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the average Frobenius aspect of the 4 corner triangles of a planar quadrilateral,...
void SetHexQualityMeasureToOddy()
Set/Get the particular estimator used to measure the quality of hexahedra.
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 WedgeEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the edge ratio of a wedge.
static double TriangleShape(vtkCell *cell, bool linearApproximation=false)
Calculate the shape of a triangle.
static double HexShear(vtkCell *cell, bool linearApproximation=false)
Calculate the shear of a hexahedron.
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 PyramidAverageSize
void SetWedgeQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of wedges.
static double HexScaledJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the scaled Jacobian of a hexahedron.
static double QuadMaxAngle(vtkCell *cell, bool linearApproximation=false)
Calculate the maximum (nonoriented) angle of a quadrilateral, expressed in degrees.
void SetPyramidQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of pyramids.
static double TriangleEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the edge ratio of a triangle.
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 HexMaxEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the maximum edge ratio of a hexahedron at its center.
static double QuadShape(vtkCell *cell, bool linearApproximation=false)
Calculate the shear of a quadrilateral.
virtual void SaveCellQualityOn()
This variable controls whether or not cell quality is stored as cell data in the resulting mesh or di...
static double TetNormalizedInradius(vtkCell *cell, bool linearApproximation=false)
Calculate the normalized in-radius of a tetrahedron.
QualityMeasureTypes PyramidQualityMeasure
static double HexDimension(vtkCell *cell, bool linearApproximation=false)
Calculate the dimension of a hexahedron.
static double QuadShearAndSize(vtkCell *cell, bool linearApproximation=false)
Calculate the shear and size of a quadrilateral.
void SetHexQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadEquiangleSkew(vtkCell *cell, bool linearApproximation=false)
Calculate the equiangle skew of a quadrilateral.
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.
vtkSetEnumMacro(QuadQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of quadrilaterals.
vtkSetEnumMacro(HexQualityMeasure, QualityMeasureTypes)
Set/Get the particular estimator used to measure the quality of hexahedra.
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 QuadArea(vtkCell *cell, bool linearApproximation=false)
Calculate the area of a quadrilateral.
static double QuadMaxEdgeRatio(vtkCell *cell, bool linearApproximation=false)
Calculate the maximum edge length ratio of a quadrilateral at quad center.
std::function< double(vtkCell *, bool)> CellQualityType
static double QuadCondition(vtkCell *cell, bool linearApproximation=false)
Calculate the condition number of a quadrilateral.
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.
void SetQuadQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double PyramidScaledJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the Jacobian of a pyramid.
virtual void SetHexQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexMedAspectFrobenius(vtkCell *cell, bool linearApproximation=false)
Calculate the average Frobenius aspect of the 8 corner tetrahedra of a hexahedron,...
CellQualityType GetQuadQualityMeasureFunction() const
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 WedgeCondition(vtkCell *cell, bool linearApproximation=false)
Calculate the condition number of a wedge.
static double WedgeAverageSize
static double HexDiagonal(vtkCell *cell, bool linearApproximation=false)
Calculate the diagonal of a hexahedron.
void SetWedgeQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of wedges.
QualityMeasureTypes WedgeQualityMeasure
static double TetMinAngle(vtkCell *cell, bool linearApproximation=false)
Calculate the minimal (nonoriented) dihedral angle of a tetrahedron, expressed in degrees.
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 HexVolume(vtkCell *cell, bool linearApproximation=false)
Calculate the volume of a hexahedron.
static double TetShapeAndSize(vtkCell *cell, bool linearApproximation=false)
Calculate the shape and size of a tetrahedron.
static double QuadJacobian(vtkCell *cell, bool linearApproximation=false)
Calculate the Jacobian of a quadrilateral.
int vtkTypeBool
Definition vtkABI.h:64
#define vtkDataArray
#define VTK_DEPRECATED_IN_9_6_0(reason)