VTK  9.4.20251007
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 "vtkFiltersVerdictModule.h" // For export macro
59
60VTK_ABI_NAMESPACE_BEGIN
61class vtkCell;
62class vtkDataArray;
63class vtkDoubleArray;
64class vtkMeshQualityFunctor;
65
66class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
67{
68private:
70
71public:
72 void PrintSelf(ostream& os, vtkIndent indent) override;
75
77
85 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
87
89
96 vtkSetMacro(LinearApproximation, bool);
97 vtkGetMacro(LinearApproximation, bool);
98 vtkBooleanMacro(LinearApproximation, bool);
100
146
150 static const char* QualityMeasureNames[];
151
153
161 virtual void SetTriangleQualityMeasure(int measure)
162 {
163 this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
164 }
226
227
229
242 virtual void SetQuadQualityMeasure(int measure)
243 {
244 this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
245 }
325
326
328
336 virtual void SetTetQualityMeasure(int measure)
337 {
338 this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
339 }
415
416
418
424 virtual void SetPyramidQualityMeasure(int measure)
425 {
426 this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
427 }
449
450
452
459 virtual void SetWedgeQualityMeasure(int measure)
460 {
461 this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
462 }
505
506
508
517 virtual void SetHexQualityMeasure(int measure)
518 {
519 this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
520 }
589
590
594 static double TriangleArea(vtkCell* cell);
595
603 static double TriangleEdgeRatio(vtkCell* cell);
604
612 static double TriangleAspectRatio(vtkCell* cell);
613
621 static double TriangleRadiusRatio(vtkCell* cell);
622
632 static double TriangleAspectFrobenius(vtkCell* cell);
633
637 static double TriangleMinAngle(vtkCell* cell);
638
642 static double TriangleMaxAngle(vtkCell* cell);
643
647 static double TriangleCondition(vtkCell* cell);
648
652 static double TriangleScaledJacobian(vtkCell* cell);
653
661
665 static double TriangleShape(vtkCell* cell);
666
673 static double TriangleShapeAndSize(vtkCell* cell);
674
678 static double TriangleDistortion(vtkCell* cell);
679
683 static double TriangleEquiangleSkew(vtkCell* cell);
684
691
699 static double QuadEdgeRatio(vtkCell* cell);
700
708 static double QuadAspectRatio(vtkCell* cell);
709
720 static double QuadRadiusRatio(vtkCell* cell);
721
731 static double QuadMedAspectFrobenius(vtkCell* cell);
732
742 static double QuadMaxAspectFrobenius(vtkCell* cell);
743
747 static double QuadMinAngle(vtkCell* cell);
748
752 static double QuadMaxEdgeRatio(vtkCell* cell);
753
759 static double QuadSkew(vtkCell* cell);
760
765 static double QuadTaper(vtkCell* cell);
766
772 static double QuadWarpage(vtkCell* cell);
773
778 static double QuadArea(vtkCell* cell);
779
784 static double QuadStretch(vtkCell* cell);
785
789 static double QuadMaxAngle(vtkCell* cell);
790
796 static double QuadOddy(vtkCell* cell);
797
803 static double QuadCondition(vtkCell* cell);
804
810 static double QuadJacobian(vtkCell* cell);
811
817 static double QuadScaledJacobian(vtkCell* cell);
818
823 static double QuadShear(vtkCell* cell);
824
829 static double QuadShape(vtkCell* cell);
830
839 static double QuadRelativeSizeSquared(vtkCell* cell);
840
848 static double QuadShapeAndSize(vtkCell* cell);
849
857 static double QuadShearAndSize(vtkCell* cell);
858
864 static double QuadDistortion(vtkCell* cell);
865
869 static double QuadEquiangleSkew(vtkCell* cell);
870
878 static double TetEdgeRatio(vtkCell* cell);
879
887 static double TetAspectRatio(vtkCell* cell);
888
896 static double TetRadiusRatio(vtkCell* cell);
897
908 static double TetAspectFrobenius(vtkCell* cell);
909
913 static double TetMinAngle(vtkCell* cell);
914
921 static double TetCollapseRatio(vtkCell* cell);
922
928 static double TetAspectGamma(vtkCell* cell);
929
934 static double TetVolume(vtkCell* cell);
935
940 static double TetCondition(vtkCell* cell);
941
946 static double TetJacobian(vtkCell* cell);
947
953 static double TetScaledJacobian(vtkCell* cell);
954
959 static double TetShape(vtkCell* cell);
960
969 static double TetRelativeSizeSquared(vtkCell* cell);
970
978 static double TetShapeAndSize(vtkCell* cell);
979
985 static double TetDistortion(vtkCell* cell);
986
990 static double TetEquiangleSkew(vtkCell* cell);
991
995 static double TetEquivolumeSkew(vtkCell* cell);
996
1002 static double TetMeanRatio(vtkCell* cell);
1003
1009 static double TetNormalizedInradius(vtkCell* cell);
1010
1014 static double TetSquishIndex(vtkCell* cell);
1015
1019 static double PyramidEquiangleSkew(vtkCell* cell);
1020
1025 static double PyramidJacobian(vtkCell* cell);
1026
1031 static double PyramidScaledJacobian(vtkCell* cell);
1032
1038 static double PyramidShape(vtkCell* cell);
1039
1043 static double PyramidVolume(vtkCell* cell);
1044
1049 static double WedgeCondition(vtkCell* cell);
1050
1055 static double WedgeDistortion(vtkCell* cell);
1056
1062 static double WedgeEdgeRatio(vtkCell* cell);
1063
1067 static double WedgeEquiangleSkew(vtkCell* cell);
1068
1073 static double WedgeJacobian(vtkCell* cell);
1074
1079 static double WedgeMaxAspectFrobenius(vtkCell* cell);
1080
1086 static double WedgeMaxStretch(vtkCell* cell);
1087
1093 static double WedgeMeanAspectFrobenius(vtkCell* cell);
1094
1104 static double WedgeScaledJacobian(vtkCell* cell);
1105
1111 static double WedgeShape(vtkCell* cell);
1112
1116 static double WedgeVolume(vtkCell* cell);
1117
1125 static double HexEdgeRatio(vtkCell* cell);
1126
1131 static double HexMedAspectFrobenius(vtkCell* cell);
1132
1137 static double HexMaxAspectFrobenius(vtkCell* cell);
1138
1142 static double HexMaxEdgeRatio(vtkCell* cell);
1143
1149 static double HexSkew(vtkCell* cell);
1150
1155 static double HexTaper(vtkCell* cell);
1156
1161 static double HexVolume(vtkCell* cell);
1162
1167 static double HexStretch(vtkCell* cell);
1168
1173 static double HexDiagonal(vtkCell* cell);
1174
1180 static double HexDimension(vtkCell* cell);
1181
1187 static double HexOddy(vtkCell* cell);
1188
1193 static double HexCondition(vtkCell* cell);
1194
1200 static double HexJacobian(vtkCell* cell);
1201
1207 static double HexScaledJacobian(vtkCell* cell);
1208
1213 static double HexShear(vtkCell* cell);
1214
1219 static double HexShape(vtkCell* cell);
1220
1229 static double HexRelativeSizeSquared(vtkCell* cell);
1230
1238 static double HexShapeAndSize(vtkCell* cell);
1239
1247 static double HexShearAndSize(vtkCell* cell);
1248
1254 static double HexDistortion(vtkCell* cell);
1255
1259 static double HexEquiangleSkew(vtkCell* cell);
1260
1265 static double HexNodalJacobianRatio(vtkCell* cell);
1266
1277 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1279 vtkBooleanMacro(Ratio, vtkTypeBool);
1280
1281protected:
1283 ~vtkMeshQuality() override = default;
1284
1286
1295
1296 using CellQualityType = double (*)(vtkCell*);
1303
1304 // Variables used to store the average size (2D: area / 3D: volume)
1305 static double TriangleAverageSize;
1306 static double QuadAverageSize;
1307 static double TetAverageSize;
1308 static double PyramidAverageSize;
1309 static double WedgeAverageSize;
1310 static double HexAverageSize;
1311
1312private:
1313 vtkMeshQuality(const vtkMeshQuality&) = delete;
1314 void operator=(const vtkMeshQuality&) = delete;
1315};
1316
1317VTK_ABI_NAMESPACE_END
1318#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition vtkCell.h:51
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.
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