VTK  9.5.20251102
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
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
90VTK_ABI_NAMESPACE_BEGIN
91class vtkCell;
92class vtkDataArray;
93class vtkDoubleArray;
94class vtkMeshQualityFunctor;
95
96class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
97{
98private:
99 friend class vtkMeshQualityFunctor;
100
101public:
102 void PrintSelf(ostream& os, vtkIndent indent) override;
105
107
113 vtkSetMacro(SaveCellQuality, vtkTypeBool);
114 vtkGetMacro(SaveCellQuality, vtkTypeBool);
115 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
117
119
126 vtkSetMacro(LinearApproximation, bool);
127 vtkGetMacro(LinearApproximation, bool);
128 vtkBooleanMacro(LinearApproximation, bool);
130
146 {
147 AREA = 28,
148 ASPECT_FROBENIUS = 3,
149 ASPECT_GAMMA = 27,
150 ASPECT_RATIO = 1,
151 COLLAPSE_RATIO = 7,
152 CONDITION = 9,
153 DIAGONAL = 21,
154 DIMENSION = 22,
155 DISTORTION = 15,
156 EDGE_RATIO = 0,
157 EQUIANGLE_SKEW = 29,
158 EQUIVOLUME_SKEW = 30,
159 JACOBIAN = 25,
160 MAX_ANGLE = 8,
161 MAX_ASPECT_FROBENIUS = 5,
162 MAX_EDGE_RATIO = 16,
163 MAX_STRETCH = 31,
164 MEAN_ASPECT_FROBENIUS = 32,
165 MEAN_RATIO = 33,
166 MED_ASPECT_FROBENIUS = 4,
167 MIN_ANGLE = 6,
168 NODAL_JACOBIAN_RATIO = 34,
169 NORMALIZED_INRADIUS = 35,
170 ODDY = 23,
171 RADIUS_RATIO = 2,
172 RELATIVE_SIZE_SQUARED = 12,
173 SCALED_JACOBIAN = 10,
174 SHAPE = 13,
175 SHAPE_AND_SIZE = 14,
176 SHEAR = 11,
177 SHEAR_AND_SIZE = 24,
178 SKEW = 17,
179 SQUISH_INDEX = 36,
180 STRETCH = 20,
181 TAPER = 18,
182 VOLUME = 19,
183 WARPAGE = 26,
184 TOTAL_QUALITY_MEASURE_TYPES = 37,
185 NONE = TOTAL_QUALITY_MEASURE_TYPES
186 };
187
191 static const char* QualityMeasureNames[];
192
194
201 vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
202 virtual void SetTriangleQualityMeasure(int measure)
203 {
204 this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
205 }
206 vtkGetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
208 {
209 this->SetTriangleQualityMeasure(QualityMeasureTypes::AREA);
210 }
212 {
213 this->SetTriangleQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
214 }
216 {
217 this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
218 }
220 {
221 this->SetTriangleQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
222 }
224 {
225 this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
226 }
228 {
229 this->SetTriangleQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
230 }
232 {
233 this->SetTriangleQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
234 }
236 {
237 this->SetTriangleQualityMeasure(QualityMeasureTypes::CONDITION);
238 }
240 {
241 this->SetTriangleQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
242 }
244 {
245 this->SetTriangleQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
246 }
248 {
249 this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE);
250 }
252 {
253 this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
254 }
256 {
257 this->SetTriangleQualityMeasure(QualityMeasureTypes::DISTORTION);
258 }
260 {
261 this->SetTriangleQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
262 }
264 {
265 this->SetTriangleQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
266 }
268
270
283 virtual void SetQuadQualityMeasure(int measure)
284 {
285 this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
286 }
289 {
290 this->SetQuadQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
291 }
293 {
294 this->SetQuadQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
295 }
297 {
298 this->SetQuadQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
299 }
301 {
302 this->SetQuadQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
303 }
305 {
306 this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
307 }
309 {
310 this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
311 }
312 void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(QualityMeasureTypes::SKEW); }
313 void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(QualityMeasureTypes::TAPER); }
315 {
316 this->SetQuadQualityMeasure(QualityMeasureTypes::WARPAGE);
317 }
318 void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(QualityMeasureTypes::AREA); }
320 {
321 this->SetQuadQualityMeasure(QualityMeasureTypes::STRETCH);
322 }
324 {
325 this->SetQuadQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
326 }
328 {
329 this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
330 }
331 void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(QualityMeasureTypes::ODDY); }
333 {
334 this->SetQuadQualityMeasure(QualityMeasureTypes::CONDITION);
335 }
337 {
338 this->SetQuadQualityMeasure(QualityMeasureTypes::JACOBIAN);
339 }
341 {
342 this->SetQuadQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
343 }
344 void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR); }
345 void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE); }
347 {
348 this->SetQuadQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
349 }
351 {
352 this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
353 }
355 {
356 this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
357 }
359 {
360 this->SetQuadQualityMeasure(QualityMeasureTypes::DISTORTION);
361 }
363 {
364 this->SetQuadQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
365 }
367
369
377 virtual void SetTetQualityMeasure(int measure)
378 {
379 this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
380 }
383 {
384 this->SetTetQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
385 }
387 {
388 this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
389 }
391 {
392 this->SetTetQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
393 }
395 {
396 this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
397 }
399 {
400 this->SetTetQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
401 }
403 {
404 this->SetTetQualityMeasure(QualityMeasureTypes::COLLAPSE_RATIO);
405 }
407 {
408 this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_GAMMA);
409 }
410 void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(QualityMeasureTypes::VOLUME); }
412 {
413 this->SetTetQualityMeasure(QualityMeasureTypes::CONDITION);
414 }
416 {
417 this->SetTetQualityMeasure(QualityMeasureTypes::JACOBIAN);
418 }
420 {
421 this->SetTetQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
422 }
423 void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE); }
425 {
426 this->SetTetQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
427 }
429 {
430 this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
431 }
433 {
434 this->SetTetQualityMeasure(QualityMeasureTypes::DISTORTION);
435 }
437 {
438 this->SetTetQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
439 }
441 {
442 this->SetTetQualityMeasure(QualityMeasureTypes::EQUIVOLUME_SKEW);
443 }
445 {
446 this->SetTetQualityMeasure(QualityMeasureTypes::MEAN_RATIO);
447 }
449 {
450 this->SetTetQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
451 }
453 {
454 this->SetTetQualityMeasure(QualityMeasureTypes::SQUISH_INDEX);
455 }
457
459
464 vtkSetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
465 virtual void SetPyramidQualityMeasure(int measure)
466 {
467 this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
468 }
469 vtkGetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
471 {
472 this->SetPyramidQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
473 }
475 {
476 this->SetPyramidQualityMeasure(QualityMeasureTypes::JACOBIAN);
477 }
479 {
480 this->SetPyramidQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
481 }
483 {
484 this->SetPyramidQualityMeasure(QualityMeasureTypes::SHAPE);
485 }
487 {
488 this->SetPyramidQualityMeasure(QualityMeasureTypes::VOLUME);
489 }
491
493
500 virtual void SetWedgeQualityMeasure(int measure)
501 {
502 this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
503 }
506 {
507 this->SetWedgeQualityMeasure(QualityMeasureTypes::CONDITION);
508 }
510 {
511 this->SetWedgeQualityMeasure(QualityMeasureTypes::DISTORTION);
512 }
514 {
515 this->SetWedgeQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
516 }
518 {
519 this->SetWedgeQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
520 }
522 {
523 this->SetWedgeQualityMeasure(QualityMeasureTypes::JACOBIAN);
524 }
526 {
527 this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
528 }
530 {
531 this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_STRETCH);
532 }
534 {
535 this->SetWedgeQualityMeasure(QualityMeasureTypes::MEAN_ASPECT_FROBENIUS);
536 }
538 {
539 this->SetWedgeQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
540 }
541 void SetWedgeQualityMeasureToShape() { this->SetWedgeQualityMeasure(QualityMeasureTypes::SHAPE); }
543 {
544 this->SetWedgeQualityMeasure(QualityMeasureTypes::VOLUME);
545 }
547
549
558 virtual void SetHexQualityMeasure(int measure)
559 {
560 this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
561 }
564 {
565 this->SetHexQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
566 }
568 {
569 this->SetHexQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
570 }
572 {
573 this->SetHexQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
574 }
576 {
577 this->SetHexQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
578 }
579 void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(QualityMeasureTypes::SKEW); }
580 void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(QualityMeasureTypes::TAPER); }
581 void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(QualityMeasureTypes::VOLUME); }
582 void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(QualityMeasureTypes::STRETCH); }
584 {
585 this->SetHexQualityMeasure(QualityMeasureTypes::DIAGONAL);
586 }
588 {
589 this->SetHexQualityMeasure(QualityMeasureTypes::DIMENSION);
590 }
591 void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(QualityMeasureTypes::ODDY); }
593 {
594 this->SetHexQualityMeasure(QualityMeasureTypes::CONDITION);
595 }
597 {
598 this->SetHexQualityMeasure(QualityMeasureTypes::JACOBIAN);
599 }
601 {
602 this->SetHexQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
603 }
604 void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR); }
605 void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE); }
607 {
608 this->SetHexQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
609 }
611 {
612 this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
613 }
615 {
616 this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
617 }
619 {
620 this->SetHexQualityMeasure(QualityMeasureTypes::DISTORTION);
621 }
623 {
624 this->SetHexQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
625 }
627 {
628 this->SetHexQualityMeasure(QualityMeasureTypes::NODAL_JACOBIAN_RATIO);
629 }
631
635 static double TriangleArea(vtkCell* cell);
636
644 static double TriangleEdgeRatio(vtkCell* cell);
645
653 static double TriangleAspectRatio(vtkCell* cell);
654
662 static double TriangleRadiusRatio(vtkCell* cell);
663
673 static double TriangleAspectFrobenius(vtkCell* cell);
674
678 static double TriangleMinAngle(vtkCell* cell);
679
683 static double TriangleMaxAngle(vtkCell* cell);
684
688 static double TriangleCondition(vtkCell* cell);
689
693 static double TriangleScaledJacobian(vtkCell* cell);
694
702
706 static double TriangleShape(vtkCell* cell);
707
714 static double TriangleShapeAndSize(vtkCell* cell);
715
719 static double TriangleDistortion(vtkCell* cell);
720
724 static double TriangleEquiangleSkew(vtkCell* cell);
725
732
740 static double QuadEdgeRatio(vtkCell* cell);
741
749 static double QuadAspectRatio(vtkCell* cell);
750
761 static double QuadRadiusRatio(vtkCell* cell);
762
772 static double QuadMedAspectFrobenius(vtkCell* cell);
773
783 static double QuadMaxAspectFrobenius(vtkCell* cell);
784
788 static double QuadMinAngle(vtkCell* cell);
789
793 static double QuadMaxEdgeRatio(vtkCell* cell);
794
800 static double QuadSkew(vtkCell* cell);
801
806 static double QuadTaper(vtkCell* cell);
807
813 static double QuadWarpage(vtkCell* cell);
814
819 static double QuadArea(vtkCell* cell);
820
825 static double QuadStretch(vtkCell* cell);
826
830 static double QuadMaxAngle(vtkCell* cell);
831
837 static double QuadOddy(vtkCell* cell);
838
844 static double QuadCondition(vtkCell* cell);
845
851 static double QuadJacobian(vtkCell* cell);
852
858 static double QuadScaledJacobian(vtkCell* cell);
859
864 static double QuadShear(vtkCell* cell);
865
870 static double QuadShape(vtkCell* cell);
871
880 static double QuadRelativeSizeSquared(vtkCell* cell);
881
889 static double QuadShapeAndSize(vtkCell* cell);
890
898 static double QuadShearAndSize(vtkCell* cell);
899
905 static double QuadDistortion(vtkCell* cell);
906
910 static double QuadEquiangleSkew(vtkCell* cell);
911
919 static double TetEdgeRatio(vtkCell* cell);
920
928 static double TetAspectRatio(vtkCell* cell);
929
937 static double TetRadiusRatio(vtkCell* cell);
938
949 static double TetAspectFrobenius(vtkCell* cell);
950
954 static double TetMinAngle(vtkCell* cell);
955
962 static double TetCollapseRatio(vtkCell* cell);
963
969 static double TetAspectGamma(vtkCell* cell);
970
975 static double TetVolume(vtkCell* cell);
976
981 static double TetCondition(vtkCell* cell);
982
987 static double TetJacobian(vtkCell* cell);
988
994 static double TetScaledJacobian(vtkCell* cell);
995
1000 static double TetShape(vtkCell* cell);
1001
1010 static double TetRelativeSizeSquared(vtkCell* cell);
1011
1019 static double TetShapeAndSize(vtkCell* cell);
1020
1026 static double TetDistortion(vtkCell* cell);
1027
1031 static double TetEquiangleSkew(vtkCell* cell);
1032
1036 static double TetEquivolumeSkew(vtkCell* cell);
1037
1043 static double TetMeanRatio(vtkCell* cell);
1044
1050 static double TetNormalizedInradius(vtkCell* cell);
1051
1055 static double TetSquishIndex(vtkCell* cell);
1056
1060 static double PyramidEquiangleSkew(vtkCell* cell);
1061
1066 static double PyramidJacobian(vtkCell* cell);
1067
1072 static double PyramidScaledJacobian(vtkCell* cell);
1073
1079 static double PyramidShape(vtkCell* cell);
1080
1084 static double PyramidVolume(vtkCell* cell);
1085
1090 static double WedgeCondition(vtkCell* cell);
1091
1096 static double WedgeDistortion(vtkCell* cell);
1097
1103 static double WedgeEdgeRatio(vtkCell* cell);
1104
1108 static double WedgeEquiangleSkew(vtkCell* cell);
1109
1114 static double WedgeJacobian(vtkCell* cell);
1115
1120 static double WedgeMaxAspectFrobenius(vtkCell* cell);
1121
1127 static double WedgeMaxStretch(vtkCell* cell);
1128
1134 static double WedgeMeanAspectFrobenius(vtkCell* cell);
1135
1145 static double WedgeScaledJacobian(vtkCell* cell);
1146
1152 static double WedgeShape(vtkCell* cell);
1153
1157 static double WedgeVolume(vtkCell* cell);
1158
1166 static double HexEdgeRatio(vtkCell* cell);
1167
1172 static double HexMedAspectFrobenius(vtkCell* cell);
1173
1178 static double HexMaxAspectFrobenius(vtkCell* cell);
1179
1183 static double HexMaxEdgeRatio(vtkCell* cell);
1184
1190 static double HexSkew(vtkCell* cell);
1191
1196 static double HexTaper(vtkCell* cell);
1197
1202 static double HexVolume(vtkCell* cell);
1203
1208 static double HexStretch(vtkCell* cell);
1209
1214 static double HexDiagonal(vtkCell* cell);
1215
1221 static double HexDimension(vtkCell* cell);
1222
1228 static double HexOddy(vtkCell* cell);
1229
1234 static double HexCondition(vtkCell* cell);
1235
1241 static double HexJacobian(vtkCell* cell);
1242
1248 static double HexScaledJacobian(vtkCell* cell);
1249
1254 static double HexShear(vtkCell* cell);
1255
1260 static double HexShape(vtkCell* cell);
1261
1270 static double HexRelativeSizeSquared(vtkCell* cell);
1271
1279 static double HexShapeAndSize(vtkCell* cell);
1280
1288 static double HexShearAndSize(vtkCell* cell);
1289
1295 static double HexDistortion(vtkCell* cell);
1296
1300 static double HexEquiangleSkew(vtkCell* cell);
1301
1306 static double HexNodalJacobianRatio(vtkCell* cell);
1307
1318 VTK_DEPRECATED_IN_9_6_0("Please use SetSaveCellQuality instead.")
1319 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1320 VTK_DEPRECATED_IN_9_6_0("Please use GetSaveCellQuality instead.")
1321 vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
1322 VTK_DEPRECATED_IN_9_6_0("Please use SaveCellQualityOn instead.")
1323 void RatioOn() { this->SaveCellQualityOn(); }
1324 VTK_DEPRECATED_IN_9_6_0("Please use SaveCellQualityOff instead.")
1325 void RatioOff() { this->SaveCellQualityOff(); }
1326
1327protected:
1329 ~vtkMeshQuality() override = default;
1330
1332
1341
1342 using CellQualityType = double (*)(vtkCell*);
1349
1350 // Variables used to store the average size (2D: area / 3D: volume)
1351 static double TriangleAverageSize;
1352 static double QuadAverageSize;
1353 static double TetAverageSize;
1354 static double PyramidAverageSize;
1355 static double WedgeAverageSize;
1356 static double HexAverageSize;
1357
1358private:
1359 vtkMeshQuality(const vtkMeshQuality&) = delete;
1360 void operator=(const vtkMeshQuality&) = delete;
1361};
1362
1363VTK_ABI_NAMESPACE_END
1364#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition vtkCell.h:129
abstract superclass for arrays of numeric data
Superclass for algorithms that produce output of the same type as input.
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.
Calculate functions of quality of the elements of a mesh.
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.
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.
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.
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.
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.
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.
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 VTK_DEPRECATED_IN_9_6_0(reason)