VTK  9.5.20251126
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 "vtkDeprecation.h" // For deprecation macro
88#include "vtkFiltersVerdictModule.h" // For export macro
89
90#include <functional> // For std::function
91
92VTK_ABI_NAMESPACE_BEGIN
93class vtkCell;
94class vtkDataArray;
95class vtkDoubleArray;
96
97class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
98{
99private:
100 class vtkMeshQualityFunctor;
101
102public:
103 void PrintSelf(ostream& os, vtkIndent indent) override;
106
108
116 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
118
120
127 vtkSetMacro(LinearApproximation, bool);
128 vtkGetMacro(LinearApproximation, bool);
129 vtkBooleanMacro(LinearApproximation, bool);
131
137 static void ComputeAverageCellSize(vtkDataSet* dataset);
138
234
238 static const char* QualityMeasureNames[];
239
241
249 virtual void SetTriangleQualityMeasure(int measure)
250 {
251 this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
252 }
314
315
317
330 virtual void SetQuadQualityMeasure(int measure)
331 {
332 this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
333 }
413
414
416
424 virtual void SetTetQualityMeasure(int measure)
425 {
426 this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
427 }
507
508
510
516 virtual void SetPyramidQualityMeasure(int measure)
517 {
518 this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
519 }
541
542
544
551 virtual void SetWedgeQualityMeasure(int measure)
552 {
553 this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
554 }
597
598
600
609 virtual void SetHexQualityMeasure(int measure)
610 {
611 this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
612 }
681
682
688 static double TriangleArea(vtkCell* cell, bool linearApproximation = false);
689
701 static double TriangleAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
702
712 static double TriangleAspectRatio(vtkCell* cell, bool linearApproximation = false);
713
719 static double TriangleCondition(vtkCell* cell, bool linearApproximation = false);
720
727 static double TriangleDistortion(vtkCell* cell, bool linearApproximation = false);
728
738 static double TriangleEdgeRatio(vtkCell* cell, bool linearApproximation = false);
739
745 static double TriangleEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
746
752 static double TriangleMaxAngle(vtkCell* cell, bool linearApproximation = false);
753
759 static double TriangleMinAngle(vtkCell* cell, bool linearApproximation = false);
760
768 static double TriangleNormalizedInradius(vtkCell* cell, bool linearApproximation = false);
769
779 static double TriangleRadiusRatio(vtkCell* cell, bool linearApproximation = false);
780
790 static double TriangleRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
791
797 static double TriangleScaledJacobian(vtkCell* cell, bool linearApproximation = false);
798
804 static double TriangleShape(vtkCell* cell, bool linearApproximation = false);
805
815 static double TriangleShapeAndSize(vtkCell* cell, bool linearApproximation = false);
816
823 static double QuadArea(vtkCell* cell, bool linearApproximation = false);
824
834 static double QuadAspectRatio(vtkCell* cell, bool linearApproximation = false);
835
843 static double QuadCondition(vtkCell* cell, bool linearApproximation = false);
844
852 static double QuadDistortion(vtkCell* cell, bool linearApproximation = false);
853
863 static double QuadEdgeRatio(vtkCell* cell, bool linearApproximation = false);
864
870 static double QuadEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
871
879 static double QuadJacobian(vtkCell* cell, bool linearApproximation = false);
880
886 static double QuadMaxAngle(vtkCell* cell, bool linearApproximation = false);
887
899 static double QuadMaxAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
900
906 static double QuadMaxEdgeRatio(vtkCell* cell, bool linearApproximation = false);
907
919 static double QuadMedAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
920
926 static double QuadMinAngle(vtkCell* cell, bool linearApproximation = false);
927
935 static double QuadOddy(vtkCell* cell, bool linearApproximation = false);
936
949 static double QuadRadiusRatio(vtkCell* cell, bool linearApproximation = false);
950
962 static double QuadRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
963
971 static double QuadScaledJacobian(vtkCell* cell, bool linearApproximation = false);
972
979 static double QuadShape(vtkCell* cell, bool linearApproximation = false);
980
991 static double QuadShapeAndSize(vtkCell* cell, bool linearApproximation = false);
992
999 static double QuadShear(vtkCell* cell, bool linearApproximation = false);
1000
1011 static double QuadShearAndSize(vtkCell* cell, bool linearApproximation = false);
1012
1020 static double QuadSkew(vtkCell* cell, bool linearApproximation = false);
1021
1028 static double QuadStretch(vtkCell* cell, bool linearApproximation = false);
1029
1036 static double QuadTaper(vtkCell* cell, bool linearApproximation = false);
1037
1045 static double QuadWarpage(vtkCell* cell, bool linearApproximation = false);
1046
1059 static double TetAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1060
1068 static double TetAspectGamma(vtkCell* cell, bool linearApproximation = false);
1069
1079 static double TetAspectRatio(vtkCell* cell, bool linearApproximation = false);
1080
1089 static double TetCollapseRatio(vtkCell* cell, bool linearApproximation = false);
1090
1097 static double TetCondition(vtkCell* cell, bool linearApproximation = false);
1098
1107 static double TetDistortion(vtkCell* cell, bool linearApproximation = false);
1108
1118 static double TetEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1119
1125 static double TetEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1126
1132 static double TetEquivolumeSkew(vtkCell* cell, bool linearApproximation = false);
1133
1140 static double TetInradius(vtkCell* cell, bool linearApproximation = false);
1141
1148 static double TetJacobian(vtkCell* cell, bool linearApproximation = false);
1149
1157 static double TetMeanRatio(vtkCell* cell, bool linearApproximation = false);
1158
1164 static double TetMinAngle(vtkCell* cell, bool linearApproximation = false);
1165
1173 static double TetNormalizedInradius(vtkCell* cell, bool linearApproximation = false);
1174
1184 static double TetRadiusRatio(vtkCell* cell, bool linearApproximation = false);
1185
1197 static double TetRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
1198
1206 static double TetScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1207
1214 static double TetShape(vtkCell* cell, bool linearApproximation = false);
1215
1226 static double TetShapeAndSize(vtkCell* cell, bool linearApproximation = false);
1227
1233 static double TetSquishIndex(vtkCell* cell, bool linearApproximation = false);
1234
1241 static double TetVolume(vtkCell* cell, bool linearApproximation = false);
1242
1248 static double PyramidEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1249
1256 static double PyramidJacobian(vtkCell* cell, bool linearApproximation = false);
1257
1264 static double PyramidScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1265
1273 static double PyramidShape(vtkCell* cell, bool linearApproximation = false);
1274
1280 static double PyramidVolume(vtkCell* cell, bool linearApproximation = false);
1281
1288 static double WedgeCondition(vtkCell* cell, bool linearApproximation = false);
1289
1296 static double WedgeDistortion(vtkCell* cell, bool linearApproximation = false);
1297
1305 static double WedgeEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1306
1312 static double WedgeEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1313
1320 static double WedgeJacobian(vtkCell* cell, bool linearApproximation = false);
1321
1328 static double WedgeMaxAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1329
1337 static double WedgeMaxStretch(vtkCell* cell, bool linearApproximation = false);
1338
1346 static double WedgeMeanAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1347
1359 static double WedgeScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1360
1368 static double WedgeShape(vtkCell* cell, bool linearApproximation = false);
1369
1375 static double WedgeVolume(vtkCell* cell, bool linearApproximation = false);
1376
1383 static double HexCondition(vtkCell* cell, bool linearApproximation = false);
1384
1391 static double HexDiagonal(vtkCell* cell, bool linearApproximation = false);
1392
1400 static double HexDimension(vtkCell* cell, bool linearApproximation = false);
1401
1409 static double HexDistortion(vtkCell* cell, bool linearApproximation = false);
1410
1420 static double HexEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1421
1427 static double HexEquiangleSkew(vtkCell* cell, bool linearApproximation = false);
1428
1436 static double HexJacobian(vtkCell* cell, bool linearApproximation = false);
1437
1444 static double HexMaxAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1445
1451 static double HexMaxEdgeRatio(vtkCell* cell, bool linearApproximation = false);
1452
1459 static double HexMedAspectFrobenius(vtkCell* cell, bool linearApproximation = false);
1460
1467 static double HexNodalJacobianRatio(vtkCell* cell, bool linearApproximation = false);
1468
1476 static double HexOddy(vtkCell* cell, bool linearApproximation = false);
1477
1489 static double HexRelativeSizeSquared(vtkCell* cell, bool linearApproximation = false);
1490
1498 static double HexScaledJacobian(vtkCell* cell, bool linearApproximation = false);
1499
1506 static double HexShape(vtkCell* cell, bool linearApproximation = false);
1507
1518 static double HexShapeAndSize(vtkCell* cell, bool linearApproximation = false);
1519
1526 static double HexShear(vtkCell* cell, bool linearApproximation = false);
1527
1537 static double HexShearAndSize(vtkCell* cell, bool linearApproximation = false);
1538
1547 static double HexSkew(vtkCell* cell, bool linearApproximation = false);
1548
1555 static double HexStretch(vtkCell* cell, bool linearApproximation = false);
1556
1563 static double HexTaper(vtkCell* cell, bool linearApproximation = false);
1570 static double HexVolume(vtkCell* cell, bool linearApproximation = false);
1571
1582 VTK_DEPRECATED_IN_9_6_0("Please use SetSaveCellQuality instead.")
1583 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1584 VTK_DEPRECATED_IN_9_6_0("Please use GetSaveCellQuality instead.")
1586 VTK_DEPRECATED_IN_9_6_0("Please use SaveCellQualityOn instead.")
1587 void RatioOn() { this->SaveCellQualityOn(); }
1588 VTK_DEPRECATED_IN_9_6_0("Please use SaveCellQualityOff instead.")
1589 void RatioOff() { this->SaveCellQualityOff(); }
1590
1591protected:
1593 ~vtkMeshQuality() override = default;
1594
1596
1605
1606 using CellQualityType = std::function<double(vtkCell*, bool)>;
1620
1621 // Variables used to store the average size (2D: area / 3D: volume)
1622 static double TriangleAverageSize;
1623 static double QuadAverageSize;
1624 static double TetAverageSize;
1625 static double PyramidAverageSize;
1626 static double WedgeAverageSize;
1627 static double HexAverageSize;
1628
1629private:
1630 vtkMeshQuality(const vtkMeshQuality&) = delete;
1631 void operator=(const vtkMeshQuality&) = delete;
1632};
1633
1634VTK_ABI_NAMESPACE_END
1635#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition vtkCell.h:129
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
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.
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)