VTK  9.3.20240328
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 "vtkFiltersVerdictModule.h" // For export macro
88 
89 VTK_ABI_NAMESPACE_BEGIN
90 class vtkCell;
91 class vtkDataArray;
92 class vtkDoubleArray;
93 class vtkMeshQualityFunctor;
94 
95 class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
96 {
97 private:
98  friend class vtkMeshQualityFunctor;
99 
100 public:
101  void PrintSelf(ostream& os, vtkIndent indent) override;
103  static vtkMeshQuality* New();
104 
106 
112  vtkSetMacro(SaveCellQuality, vtkTypeBool);
113  vtkGetMacro(SaveCellQuality, vtkTypeBool);
114  vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
116 
118 
125  vtkSetMacro(LinearApproximation, bool);
126  vtkGetMacro(LinearApproximation, bool);
127  vtkBooleanMacro(LinearApproximation, bool);
129 
134  {
135  AREA = 28,
136  ASPECT_FROBENIUS = 3,
137  ASPECT_GAMMA = 27,
138  ASPECT_RATIO = 1,
139  COLLAPSE_RATIO = 7,
140  CONDITION = 9,
141  DIAGONAL = 21,
142  DIMENSION = 22,
143  DISTORTION = 15,
144  EDGE_RATIO = 0,
145  EQUIANGLE_SKEW = 29,
146  EQUIVOLUME_SKEW = 30,
147  JACOBIAN = 25,
148  MAX_ANGLE = 8,
149  MAX_ASPECT_FROBENIUS = 5,
150  MAX_EDGE_RATIO = 16,
151  MAX_STRETCH = 31,
152  MEAN_ASPECT_FROBENIUS = 32,
153  MEAN_RATIO = 33,
154  MED_ASPECT_FROBENIUS = 4,
155  MIN_ANGLE = 6,
156  NODAL_JACOBIAN_RATIO = 34,
157  NORMALIZED_INRADIUS = 35,
158  ODDY = 23,
159  RADIUS_RATIO = 2,
160  RELATIVE_SIZE_SQUARED = 12,
161  SCALED_JACOBIAN = 10,
162  SHAPE = 13,
163  SHAPE_AND_SIZE = 14,
164  SHEAR = 11,
165  SHEAR_AND_SIZE = 24,
166  SKEW = 17,
167  SQUISH_INDEX = 36,
168  STRETCH = 20,
169  TAPER = 18,
170  VOLUME = 19,
171  WARPAGE = 26,
172  TOTAL_QUALITY_MEASURE_TYPES = 37,
173  NONE = TOTAL_QUALITY_MEASURE_TYPES
174  };
175 
179  static const char* QualityMeasureNames[];
180 
182 
189  vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
190  virtual void SetTriangleQualityMeasure(int measure)
191  {
192  this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
193  }
194  vtkGetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
196  {
197  this->SetTriangleQualityMeasure(QualityMeasureTypes::AREA);
198  }
200  {
201  this->SetTriangleQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
202  }
204  {
205  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
206  }
208  {
209  this->SetTriangleQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
210  }
212  {
213  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
214  }
216  {
217  this->SetTriangleQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
218  }
220  {
221  this->SetTriangleQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
222  }
224  {
225  this->SetTriangleQualityMeasure(QualityMeasureTypes::CONDITION);
226  }
228  {
229  this->SetTriangleQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
230  }
232  {
233  this->SetTriangleQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
234  }
236  {
237  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE);
238  }
240  {
241  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
242  }
244  {
245  this->SetTriangleQualityMeasure(QualityMeasureTypes::DISTORTION);
246  }
248  {
249  this->SetTriangleQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
250  }
252  {
253  this->SetTriangleQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
254  }
256 
258 
270  vtkSetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
271  virtual void SetQuadQualityMeasure(int measure)
272  {
273  this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
274  }
275  vtkGetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
277  {
278  this->SetQuadQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
279  }
281  {
282  this->SetQuadQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
283  }
285  {
286  this->SetQuadQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
287  }
289  {
290  this->SetQuadQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
291  }
293  {
294  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
295  }
297  {
298  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
299  }
300  void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(QualityMeasureTypes::SKEW); }
301  void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(QualityMeasureTypes::TAPER); }
303  {
304  this->SetQuadQualityMeasure(QualityMeasureTypes::WARPAGE);
305  }
306  void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(QualityMeasureTypes::AREA); }
308  {
309  this->SetQuadQualityMeasure(QualityMeasureTypes::STRETCH);
310  }
312  {
313  this->SetQuadQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
314  }
316  {
317  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
318  }
319  void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(QualityMeasureTypes::ODDY); }
321  {
322  this->SetQuadQualityMeasure(QualityMeasureTypes::CONDITION);
323  }
325  {
326  this->SetQuadQualityMeasure(QualityMeasureTypes::JACOBIAN);
327  }
329  {
330  this->SetQuadQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
331  }
332  void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR); }
333  void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE); }
335  {
336  this->SetQuadQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
337  }
339  {
340  this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
341  }
343  {
344  this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
345  }
347  {
348  this->SetQuadQualityMeasure(QualityMeasureTypes::DISTORTION);
349  }
351  {
352  this->SetQuadQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
353  }
355 
357 
365  virtual void SetTetQualityMeasure(int measure)
366  {
367  this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
368  }
371  {
372  this->SetTetQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
373  }
375  {
376  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
377  }
379  {
380  this->SetTetQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
381  }
383  {
384  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
385  }
387  {
388  this->SetTetQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
389  }
391  {
392  this->SetTetQualityMeasure(QualityMeasureTypes::COLLAPSE_RATIO);
393  }
395  {
396  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_GAMMA);
397  }
398  void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(QualityMeasureTypes::VOLUME); }
400  {
401  this->SetTetQualityMeasure(QualityMeasureTypes::CONDITION);
402  }
404  {
405  this->SetTetQualityMeasure(QualityMeasureTypes::JACOBIAN);
406  }
408  {
409  this->SetTetQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
410  }
411  void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE); }
413  {
414  this->SetTetQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
415  }
417  {
418  this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
419  }
421  {
422  this->SetTetQualityMeasure(QualityMeasureTypes::DISTORTION);
423  }
425  {
426  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
427  }
429  {
430  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIVOLUME_SKEW);
431  }
433  {
434  this->SetTetQualityMeasure(QualityMeasureTypes::MEAN_RATIO);
435  }
437  {
438  this->SetTetQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
439  }
441  {
442  this->SetTetQualityMeasure(QualityMeasureTypes::SQUISH_INDEX);
443  }
445 
447 
452  vtkSetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
453  virtual void SetPyramidQualityMeasure(int measure)
454  {
455  this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
456  }
457  vtkGetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
459  {
460  this->SetPyramidQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
461  }
463  {
464  this->SetPyramidQualityMeasure(QualityMeasureTypes::JACOBIAN);
465  }
467  {
468  this->SetPyramidQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
469  }
471  {
472  this->SetPyramidQualityMeasure(QualityMeasureTypes::SHAPE);
473  }
475  {
476  this->SetPyramidQualityMeasure(QualityMeasureTypes::VOLUME);
477  }
479 
481 
487  vtkSetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
488  virtual void SetWedgeQualityMeasure(int measure)
489  {
490  this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
491  }
492  vtkGetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
494  {
495  this->SetWedgeQualityMeasure(QualityMeasureTypes::CONDITION);
496  }
498  {
499  this->SetWedgeQualityMeasure(QualityMeasureTypes::DISTORTION);
500  }
502  {
503  this->SetWedgeQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
504  }
506  {
507  this->SetWedgeQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
508  }
510  {
511  this->SetWedgeQualityMeasure(QualityMeasureTypes::JACOBIAN);
512  }
514  {
515  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
516  }
518  {
519  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_STRETCH);
520  }
522  {
523  this->SetWedgeQualityMeasure(QualityMeasureTypes::MEAN_ASPECT_FROBENIUS);
524  }
526  {
527  this->SetWedgeQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
528  }
529  void SetWedgeQualityMeasureToShape() { this->SetWedgeQualityMeasure(QualityMeasureTypes::SHAPE); }
531  {
532  this->SetWedgeQualityMeasure(QualityMeasureTypes::VOLUME);
533  }
535 
537 
546  virtual void SetHexQualityMeasure(int measure)
547  {
548  this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
549  }
552  {
553  this->SetHexQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
554  }
556  {
557  this->SetHexQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
558  }
560  {
561  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
562  }
564  {
565  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
566  }
567  void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(QualityMeasureTypes::SKEW); }
568  void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(QualityMeasureTypes::TAPER); }
569  void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(QualityMeasureTypes::VOLUME); }
570  void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(QualityMeasureTypes::STRETCH); }
572  {
573  this->SetHexQualityMeasure(QualityMeasureTypes::DIAGONAL);
574  }
576  {
577  this->SetHexQualityMeasure(QualityMeasureTypes::DIMENSION);
578  }
579  void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(QualityMeasureTypes::ODDY); }
581  {
582  this->SetHexQualityMeasure(QualityMeasureTypes::CONDITION);
583  }
585  {
586  this->SetHexQualityMeasure(QualityMeasureTypes::JACOBIAN);
587  }
589  {
590  this->SetHexQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
591  }
592  void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR); }
593  void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE); }
595  {
596  this->SetHexQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
597  }
599  {
600  this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
601  }
603  {
604  this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
605  }
607  {
608  this->SetHexQualityMeasure(QualityMeasureTypes::DISTORTION);
609  }
611  {
612  this->SetHexQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
613  }
615  {
616  this->SetHexQualityMeasure(QualityMeasureTypes::NODAL_JACOBIAN_RATIO);
617  }
619 
623  static double TriangleArea(vtkCell* cell);
624 
632  static double TriangleEdgeRatio(vtkCell* cell);
633 
641  static double TriangleAspectRatio(vtkCell* cell);
642 
650  static double TriangleRadiusRatio(vtkCell* cell);
651 
661  static double TriangleAspectFrobenius(vtkCell* cell);
662 
666  static double TriangleMinAngle(vtkCell* cell);
667 
671  static double TriangleMaxAngle(vtkCell* cell);
672 
676  static double TriangleCondition(vtkCell* cell);
677 
681  static double TriangleScaledJacobian(vtkCell* cell);
682 
689  static double TriangleRelativeSizeSquared(vtkCell* cell);
690 
694  static double TriangleShape(vtkCell* cell);
695 
702  static double TriangleShapeAndSize(vtkCell* cell);
703 
707  static double TriangleDistortion(vtkCell* cell);
708 
712  static double TriangleEquiangleSkew(vtkCell* cell);
713 
719  static double TriangleNormalizedInradius(vtkCell* cell);
720 
728  static double QuadEdgeRatio(vtkCell* cell);
729 
737  static double QuadAspectRatio(vtkCell* cell);
738 
749  static double QuadRadiusRatio(vtkCell* cell);
750 
760  static double QuadMedAspectFrobenius(vtkCell* cell);
761 
771  static double QuadMaxAspectFrobenius(vtkCell* cell);
772 
776  static double QuadMinAngle(vtkCell* cell);
777 
781  static double QuadMaxEdgeRatio(vtkCell* cell);
782 
788  static double QuadSkew(vtkCell* cell);
789 
794  static double QuadTaper(vtkCell* cell);
795 
801  static double QuadWarpage(vtkCell* cell);
802 
807  static double QuadArea(vtkCell* cell);
808 
813  static double QuadStretch(vtkCell* cell);
814 
818  static double QuadMaxAngle(vtkCell* cell);
819 
825  static double QuadOddy(vtkCell* cell);
826 
832  static double QuadCondition(vtkCell* cell);
833 
839  static double QuadJacobian(vtkCell* cell);
840 
846  static double QuadScaledJacobian(vtkCell* cell);
847 
852  static double QuadShear(vtkCell* cell);
853 
858  static double QuadShape(vtkCell* cell);
859 
868  static double QuadRelativeSizeSquared(vtkCell* cell);
869 
877  static double QuadShapeAndSize(vtkCell* cell);
878 
886  static double QuadShearAndSize(vtkCell* cell);
887 
893  static double QuadDistortion(vtkCell* cell);
894 
898  static double QuadEquiangleSkew(vtkCell* cell);
899 
907  static double TetEdgeRatio(vtkCell* cell);
908 
916  static double TetAspectRatio(vtkCell* cell);
917 
925  static double TetRadiusRatio(vtkCell* cell);
926 
937  static double TetAspectFrobenius(vtkCell* cell);
938 
942  static double TetMinAngle(vtkCell* cell);
943 
950  static double TetCollapseRatio(vtkCell* cell);
951 
957  static double TetAspectGamma(vtkCell* cell);
958 
963  static double TetVolume(vtkCell* cell);
964 
969  static double TetCondition(vtkCell* cell);
970 
975  static double TetJacobian(vtkCell* cell);
976 
982  static double TetScaledJacobian(vtkCell* cell);
983 
988  static double TetShape(vtkCell* cell);
989 
998  static double TetRelativeSizeSquared(vtkCell* cell);
999 
1007  static double TetShapeAndSize(vtkCell* cell);
1008 
1014  static double TetDistortion(vtkCell* cell);
1015 
1019  static double TetEquiangleSkew(vtkCell* cell);
1020 
1024  static double TetEquivolumeSkew(vtkCell* cell);
1025 
1031  static double TetMeanRatio(vtkCell* cell);
1032 
1038  static double TetNormalizedInradius(vtkCell* cell);
1039 
1043  static double TetSquishIndex(vtkCell* cell);
1044 
1048  static double PyramidEquiangleSkew(vtkCell* cell);
1049 
1054  static double PyramidJacobian(vtkCell* cell);
1055 
1060  static double PyramidScaledJacobian(vtkCell* cell);
1061 
1067  static double PyramidShape(vtkCell* cell);
1068 
1072  static double PyramidVolume(vtkCell* cell);
1073 
1078  static double WedgeCondition(vtkCell* cell);
1079 
1084  static double WedgeDistortion(vtkCell* cell);
1085 
1091  static double WedgeEdgeRatio(vtkCell* cell);
1092 
1096  static double WedgeEquiangleSkew(vtkCell* cell);
1097 
1102  static double WedgeJacobian(vtkCell* cell);
1103 
1108  static double WedgeMaxAspectFrobenius(vtkCell* cell);
1109 
1115  static double WedgeMaxStretch(vtkCell* cell);
1116 
1122  static double WedgeMeanAspectFrobenius(vtkCell* cell);
1123 
1133  static double WedgeScaledJacobian(vtkCell* cell);
1134 
1140  static double WedgeShape(vtkCell* cell);
1141 
1145  static double WedgeVolume(vtkCell* cell);
1146 
1154  static double HexEdgeRatio(vtkCell* cell);
1155 
1160  static double HexMedAspectFrobenius(vtkCell* cell);
1161 
1166  static double HexMaxAspectFrobenius(vtkCell* cell);
1167 
1171  static double HexMaxEdgeRatio(vtkCell* cell);
1172 
1178  static double HexSkew(vtkCell* cell);
1179 
1184  static double HexTaper(vtkCell* cell);
1185 
1190  static double HexVolume(vtkCell* cell);
1191 
1196  static double HexStretch(vtkCell* cell);
1197 
1202  static double HexDiagonal(vtkCell* cell);
1203 
1209  static double HexDimension(vtkCell* cell);
1210 
1216  static double HexOddy(vtkCell* cell);
1217 
1222  static double HexCondition(vtkCell* cell);
1223 
1229  static double HexJacobian(vtkCell* cell);
1230 
1236  static double HexScaledJacobian(vtkCell* cell);
1237 
1242  static double HexShear(vtkCell* cell);
1243 
1248  static double HexShape(vtkCell* cell);
1249 
1258  static double HexRelativeSizeSquared(vtkCell* cell);
1259 
1267  static double HexShapeAndSize(vtkCell* cell);
1268 
1276  static double HexShearAndSize(vtkCell* cell);
1277 
1283  static double HexDistortion(vtkCell* cell);
1284 
1288  static double HexEquiangleSkew(vtkCell* cell);
1289 
1294  static double HexNodalJacobianRatio(vtkCell* cell);
1295 
1306  virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1307  vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
1308  vtkBooleanMacro(Ratio, vtkTypeBool);
1309 
1310 protected:
1312  ~vtkMeshQuality() override = default;
1313 
1315 
1324 
1325  using CellQualityType = double (*)(vtkCell*);
1332 
1333  // Variables used to store the average size (2D: area / 3D: volume)
1334  static double TriangleAverageSize;
1335  static double QuadAverageSize;
1336  static double TetAverageSize;
1337  static double PyramidAverageSize;
1338  static double WedgeAverageSize;
1339  static double HexAverageSize;
1340 
1341 private:
1342  vtkMeshQuality(const vtkMeshQuality&) = delete;
1343  void operator=(const vtkMeshQuality&) = delete;
1344 };
1345 
1346 VTK_ABI_NAMESPACE_END
1347 #endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:130
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:154
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.
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.
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.
static vtkMeshQuality * New()
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