VTK  9.3.0
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
63 #ifndef vtkMeshQuality_h
64 #define vtkMeshQuality_h
65 
66 #include "vtkDataSetAlgorithm.h"
67 #include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_2_0
68 #include "vtkFiltersVerdictModule.h" // For export macro
69 
70 VTK_ABI_NAMESPACE_BEGIN
71 class vtkCell;
72 class vtkDataArray;
73 class vtkDoubleArray;
74 class vtkMeshQualityFunctor;
75 
76 class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
77 {
78 private:
79  friend class vtkMeshQualityFunctor;
80 
81 public:
82  void PrintSelf(ostream& os, vtkIndent indent) override;
84  static vtkMeshQuality* New();
85 
87 
93  vtkSetMacro(SaveCellQuality, vtkTypeBool);
94  vtkGetMacro(SaveCellQuality, vtkTypeBool);
95  vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
97 
99 
106  vtkSetMacro(LinearApproximation, bool);
107  vtkGetMacro(LinearApproximation, bool);
108  vtkBooleanMacro(LinearApproximation, bool);
110 
115  {
116  AREA = 28,
117  ASPECT_FROBENIUS = 3,
118  ASPECT_GAMMA = 27,
119  ASPECT_RATIO = 1,
120  COLLAPSE_RATIO = 7,
121  CONDITION = 9,
122  DIAGONAL = 21,
123  DIMENSION = 22,
124  DISTORTION = 15,
125  EDGE_RATIO = 0,
126  EQUIANGLE_SKEW = 29,
127  EQUIVOLUME_SKEW = 30,
128  JACOBIAN = 25,
129  MAX_ANGLE = 8,
130  MAX_ASPECT_FROBENIUS = 5,
131  MAX_EDGE_RATIO = 16,
132  MAX_STRETCH = 31,
133  MEAN_ASPECT_FROBENIUS = 32,
134  MEAN_RATIO = 33,
135  MED_ASPECT_FROBENIUS = 4,
136  MIN_ANGLE = 6,
137  NODAL_JACOBIAN_RATIO = 34,
138  NORMALIZED_INRADIUS = 35,
139  ODDY = 23,
140  RADIUS_RATIO = 2,
141  RELATIVE_SIZE_SQUARED = 12,
142  SCALED_JACOBIAN = 10,
143  SHAPE = 13,
144  SHAPE_AND_SIZE = 14,
145  SHEAR = 11,
146  SHEAR_AND_SIZE = 24,
147  SKEW = 17,
148  SQUISH_INDEX = 36,
149  STRETCH = 20,
150  TAPER = 18,
151  VOLUME = 19,
152  WARPAGE = 26,
153  TOTAL_QUALITY_MEASURE_TYPES = 37,
154  NONE = TOTAL_QUALITY_MEASURE_TYPES
155  };
156 
160  static const char* QualityMeasureNames[];
161 
163 
170  vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
171  virtual void SetTriangleQualityMeasure(int measure)
172  {
173  this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
174  }
175  vtkGetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
177  {
178  this->SetTriangleQualityMeasure(QualityMeasureTypes::AREA);
179  }
181  {
182  this->SetTriangleQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
183  }
185  {
186  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
187  }
189  {
190  this->SetTriangleQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
191  }
193  {
194  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
195  }
197  {
198  this->SetTriangleQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
199  }
201  {
202  this->SetTriangleQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
203  }
205  {
206  this->SetTriangleQualityMeasure(QualityMeasureTypes::CONDITION);
207  }
209  {
210  this->SetTriangleQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
211  }
213  {
214  this->SetTriangleQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
215  }
217  {
218  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE);
219  }
221  {
222  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
223  }
225  {
226  this->SetTriangleQualityMeasure(QualityMeasureTypes::DISTORTION);
227  }
229  {
230  this->SetTriangleQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
231  }
233  {
234  this->SetTriangleQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
235  }
237 
239 
251  vtkSetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
252  virtual void SetQuadQualityMeasure(int measure)
253  {
254  this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
255  }
256  vtkGetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
258  {
259  this->SetQuadQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
260  }
262  {
263  this->SetQuadQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
264  }
266  {
267  this->SetQuadQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
268  }
270  {
271  this->SetQuadQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
272  }
274  {
275  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
276  }
278  {
279  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
280  }
281  void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(QualityMeasureTypes::SKEW); }
282  void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(QualityMeasureTypes::TAPER); }
284  {
285  this->SetQuadQualityMeasure(QualityMeasureTypes::WARPAGE);
286  }
287  void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(QualityMeasureTypes::AREA); }
289  {
290  this->SetQuadQualityMeasure(QualityMeasureTypes::STRETCH);
291  }
293  {
294  this->SetQuadQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
295  }
297  {
298  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
299  }
300  void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(QualityMeasureTypes::ODDY); }
302  {
303  this->SetQuadQualityMeasure(QualityMeasureTypes::CONDITION);
304  }
306  {
307  this->SetQuadQualityMeasure(QualityMeasureTypes::JACOBIAN);
308  }
310  {
311  this->SetQuadQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
312  }
313  void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR); }
314  void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE); }
316  {
317  this->SetQuadQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
318  }
320  {
321  this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
322  }
324  {
325  this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
326  }
328  {
329  this->SetQuadQualityMeasure(QualityMeasureTypes::DISTORTION);
330  }
332  {
333  this->SetQuadQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
334  }
336 
338 
346  virtual void SetTetQualityMeasure(int measure)
347  {
348  this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
349  }
352  {
353  this->SetTetQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
354  }
356  {
357  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
358  }
360  {
361  this->SetTetQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
362  }
364  {
365  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
366  }
368  {
369  this->SetTetQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
370  }
372  {
373  this->SetTetQualityMeasure(QualityMeasureTypes::COLLAPSE_RATIO);
374  }
376  {
377  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_GAMMA);
378  }
379  void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(QualityMeasureTypes::VOLUME); }
381  {
382  this->SetTetQualityMeasure(QualityMeasureTypes::CONDITION);
383  }
385  {
386  this->SetTetQualityMeasure(QualityMeasureTypes::JACOBIAN);
387  }
389  {
390  this->SetTetQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
391  }
392  void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE); }
394  {
395  this->SetTetQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
396  }
398  {
399  this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
400  }
402  {
403  this->SetTetQualityMeasure(QualityMeasureTypes::DISTORTION);
404  }
406  {
407  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
408  }
410  {
411  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIVOLUME_SKEW);
412  }
414  {
415  this->SetTetQualityMeasure(QualityMeasureTypes::MEAN_RATIO);
416  }
418  {
419  this->SetTetQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
420  }
422  {
423  this->SetTetQualityMeasure(QualityMeasureTypes::SQUISH_INDEX);
424  }
426 
428 
433  vtkSetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
434  virtual void SetPyramidQualityMeasure(int measure)
435  {
436  this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
437  }
438  vtkGetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
440  {
441  this->SetPyramidQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
442  }
444  {
445  this->SetPyramidQualityMeasure(QualityMeasureTypes::JACOBIAN);
446  }
448  {
449  this->SetPyramidQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
450  }
452  {
453  this->SetPyramidQualityMeasure(QualityMeasureTypes::SHAPE);
454  }
456  {
457  this->SetPyramidQualityMeasure(QualityMeasureTypes::VOLUME);
458  }
460 
462 
468  vtkSetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
469  virtual void SetWedgeQualityMeasure(int measure)
470  {
471  this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
472  }
473  vtkGetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
475  {
476  this->SetWedgeQualityMeasure(QualityMeasureTypes::CONDITION);
477  }
479  {
480  this->SetWedgeQualityMeasure(QualityMeasureTypes::DISTORTION);
481  }
483  {
484  this->SetWedgeQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
485  }
487  {
488  this->SetWedgeQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
489  }
491  {
492  this->SetWedgeQualityMeasure(QualityMeasureTypes::JACOBIAN);
493  }
495  {
496  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
497  }
499  {
500  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_STRETCH);
501  }
503  {
504  this->SetWedgeQualityMeasure(QualityMeasureTypes::MEAN_ASPECT_FROBENIUS);
505  }
507  {
508  this->SetWedgeQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
509  }
510  void SetWedgeQualityMeasureToShape() { this->SetWedgeQualityMeasure(QualityMeasureTypes::SHAPE); }
512  {
513  this->SetWedgeQualityMeasure(QualityMeasureTypes::VOLUME);
514  }
516 
518 
527  virtual void SetHexQualityMeasure(int measure)
528  {
529  this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
530  }
533  {
534  this->SetHexQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
535  }
537  {
538  this->SetHexQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
539  }
541  {
542  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
543  }
545  {
546  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
547  }
548  void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(QualityMeasureTypes::SKEW); }
549  void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(QualityMeasureTypes::TAPER); }
550  void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(QualityMeasureTypes::VOLUME); }
551  void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(QualityMeasureTypes::STRETCH); }
553  {
554  this->SetHexQualityMeasure(QualityMeasureTypes::DIAGONAL);
555  }
557  {
558  this->SetHexQualityMeasure(QualityMeasureTypes::DIMENSION);
559  }
560  void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(QualityMeasureTypes::ODDY); }
562  {
563  this->SetHexQualityMeasure(QualityMeasureTypes::CONDITION);
564  }
566  {
567  this->SetHexQualityMeasure(QualityMeasureTypes::JACOBIAN);
568  }
570  {
571  this->SetHexQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
572  }
573  void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR); }
574  void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE); }
576  {
577  this->SetHexQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
578  }
580  {
581  this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
582  }
584  {
585  this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
586  }
588  {
589  this->SetHexQualityMeasure(QualityMeasureTypes::DISTORTION);
590  }
592  {
593  this->SetHexQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
594  }
596  {
597  this->SetHexQualityMeasure(QualityMeasureTypes::NODAL_JACOBIAN_RATIO);
598  }
600 
604  static double TriangleArea(vtkCell* cell);
605 
613  static double TriangleEdgeRatio(vtkCell* cell);
614 
622  static double TriangleAspectRatio(vtkCell* cell);
623 
631  static double TriangleRadiusRatio(vtkCell* cell);
632 
642  static double TriangleAspectFrobenius(vtkCell* cell);
643 
647  static double TriangleMinAngle(vtkCell* cell);
648 
652  static double TriangleMaxAngle(vtkCell* cell);
653 
657  static double TriangleCondition(vtkCell* cell);
658 
662  static double TriangleScaledJacobian(vtkCell* cell);
663 
670  static double TriangleRelativeSizeSquared(vtkCell* cell);
671 
675  static double TriangleShape(vtkCell* cell);
676 
683  static double TriangleShapeAndSize(vtkCell* cell);
684 
688  static double TriangleDistortion(vtkCell* cell);
689 
693  static double TriangleEquiangleSkew(vtkCell* cell);
694 
700  static double TriangleNormalizedInradius(vtkCell* cell);
701 
709  static double QuadEdgeRatio(vtkCell* cell);
710 
718  static double QuadAspectRatio(vtkCell* cell);
719 
730  static double QuadRadiusRatio(vtkCell* cell);
731 
741  static double QuadMedAspectFrobenius(vtkCell* cell);
742 
752  static double QuadMaxAspectFrobenius(vtkCell* cell);
753 
757  static double QuadMinAngle(vtkCell* cell);
758 
762  static double QuadMaxEdgeRatio(vtkCell* cell);
763 
769  static double QuadSkew(vtkCell* cell);
770 
775  static double QuadTaper(vtkCell* cell);
776 
782  static double QuadWarpage(vtkCell* cell);
783 
788  static double QuadArea(vtkCell* cell);
789 
794  static double QuadStretch(vtkCell* cell);
795 
799  static double QuadMaxAngle(vtkCell* cell);
800 
806  static double QuadOddy(vtkCell* cell);
807 
813  static double QuadCondition(vtkCell* cell);
814 
820  static double QuadJacobian(vtkCell* cell);
821 
827  static double QuadScaledJacobian(vtkCell* cell);
828 
833  static double QuadShear(vtkCell* cell);
834 
839  static double QuadShape(vtkCell* cell);
840 
849  static double QuadRelativeSizeSquared(vtkCell* cell);
850 
858  static double QuadShapeAndSize(vtkCell* cell);
859 
867  static double QuadShearAndSize(vtkCell* cell);
868 
874  static double QuadDistortion(vtkCell* cell);
875 
879  static double QuadEquiangleSkew(vtkCell* cell);
880 
888  static double TetEdgeRatio(vtkCell* cell);
889 
897  static double TetAspectRatio(vtkCell* cell);
898 
906  static double TetRadiusRatio(vtkCell* cell);
907 
918  static double TetAspectFrobenius(vtkCell* cell);
919 
923  static double TetMinAngle(vtkCell* cell);
924 
931  static double TetCollapseRatio(vtkCell* cell);
932 
938  static double TetAspectGamma(vtkCell* cell);
939 
944  static double TetVolume(vtkCell* cell);
945 
950  static double TetCondition(vtkCell* cell);
951 
956  static double TetJacobian(vtkCell* cell);
957 
963  static double TetScaledJacobian(vtkCell* cell);
964 
969  static double TetShape(vtkCell* cell);
970 
979  static double TetRelativeSizeSquared(vtkCell* cell);
980 
988  static double TetShapeAndSize(vtkCell* cell);
989 
995  static double TetDistortion(vtkCell* cell);
996 
1000  static double TetEquiangleSkew(vtkCell* cell);
1001 
1005  static double TetEquivolumeSkew(vtkCell* cell);
1006 
1012  static double TetMeanRatio(vtkCell* cell);
1013 
1019  static double TetNormalizedInradius(vtkCell* cell);
1020 
1024  static double TetSquishIndex(vtkCell* cell);
1025 
1029  static double PyramidEquiangleSkew(vtkCell* cell);
1030 
1035  static double PyramidJacobian(vtkCell* cell);
1036 
1041  static double PyramidScaledJacobian(vtkCell* cell);
1042 
1048  static double PyramidShape(vtkCell* cell);
1049 
1053  static double PyramidVolume(vtkCell* cell);
1054 
1059  static double WedgeCondition(vtkCell* cell);
1060 
1065  static double WedgeDistortion(vtkCell* cell);
1066 
1072  static double WedgeEdgeRatio(vtkCell* cell);
1073 
1077  static double WedgeEquiangleSkew(vtkCell* cell);
1078 
1083  static double WedgeJacobian(vtkCell* cell);
1084 
1089  static double WedgeMaxAspectFrobenius(vtkCell* cell);
1090 
1096  static double WedgeMaxStretch(vtkCell* cell);
1097 
1103  static double WedgeMeanAspectFrobenius(vtkCell* cell);
1104 
1114  static double WedgeScaledJacobian(vtkCell* cell);
1115 
1121  static double WedgeShape(vtkCell* cell);
1122 
1126  static double WedgeVolume(vtkCell* cell);
1127 
1135  static double HexEdgeRatio(vtkCell* cell);
1136 
1141  static double HexMedAspectFrobenius(vtkCell* cell);
1142 
1147  static double HexMaxAspectFrobenius(vtkCell* cell);
1148 
1152  static double HexMaxEdgeRatio(vtkCell* cell);
1153 
1159  static double HexSkew(vtkCell* cell);
1160 
1165  static double HexTaper(vtkCell* cell);
1166 
1171  static double HexVolume(vtkCell* cell);
1172 
1177  static double HexStretch(vtkCell* cell);
1178 
1183  static double HexDiagonal(vtkCell* cell);
1184 
1190  static double HexDimension(vtkCell* cell);
1191 
1197  static double HexOddy(vtkCell* cell);
1198 
1203  static double HexCondition(vtkCell* cell);
1204 
1210  static double HexJacobian(vtkCell* cell);
1211 
1217  static double HexScaledJacobian(vtkCell* cell);
1218 
1223  static double HexShear(vtkCell* cell);
1224 
1229  static double HexShape(vtkCell* cell);
1230 
1239  static double HexRelativeSizeSquared(vtkCell* cell);
1240 
1248  static double HexShapeAndSize(vtkCell* cell);
1249 
1257  static double HexShearAndSize(vtkCell* cell);
1258 
1264  static double HexDistortion(vtkCell* cell);
1265 
1269  static double HexEquiangleSkew(vtkCell* cell);
1270 
1275  static double HexNodalJacobianRatio(vtkCell* cell);
1276 
1287  virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1288  vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
1289  vtkBooleanMacro(Ratio, vtkTypeBool);
1290 
1292 
1309  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1310  virtual void SetVolume(vtkTypeBool cv)
1311  {
1312  if (!((cv != 0) ^ (this->Volume != 0)))
1313  {
1314  return;
1315  }
1316  this->Modified();
1317  this->Volume = cv;
1318  if (this->Volume)
1319  {
1320  this->CompatibilityMode = 1;
1321  }
1322  }
1323  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1324  vtkTypeBool GetVolume() { return this->Volume; }
1325  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1326  void VolumeOn()
1327  {
1328  if (!this->Volume)
1329  {
1330  this->Volume = 1;
1331  this->Modified();
1332  }
1333  }
1334  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1335  void VolumeOff()
1336  {
1337  if (this->Volume)
1338  {
1339  this->Volume = 0;
1340  this->Modified();
1341  }
1342  }
1344 
1346 
1373  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1374  virtual void SetCompatibilityMode(vtkTypeBool cm)
1375  {
1376  if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
1377  {
1378  return;
1379  }
1380  this->CompatibilityMode = cm;
1381  this->Modified();
1382  if (this->CompatibilityMode)
1383  {
1384  this->Volume = 1;
1385  this->TetQualityMeasure = QualityMeasureTypes::RADIUS_RATIO;
1386  }
1387  }
1388  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1389  vtkGetMacro(CompatibilityMode, vtkTypeBool);
1390  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1391  void CompatibilityModeOn()
1392  {
1393  if (!this->CompatibilityMode)
1394  {
1395  this->CompatibilityMode = 1;
1396  this->Modified();
1397  }
1398  }
1399  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1400  void CompatibilityModeOff()
1401  {
1402  if (this->CompatibilityMode)
1403  {
1404  this->CompatibilityMode = 0;
1405  this->Modified();
1406  }
1407  }
1409 
1410 protected:
1412  ~vtkMeshQuality() override = default;
1413 
1415 
1424 
1425  using CellQualityType = double (*)(vtkCell*);
1432 
1433  // VTK_DEPRECATED_IN_9_2_0 Those 2 attributes need to be removed, and instance in the code as
1434  // well.
1437 
1438  // Variables used to store the average size (2D: area / 3D: volume)
1439  static double TriangleAverageSize;
1440  static double QuadAverageSize;
1441  static double TetAverageSize;
1442  static double PyramidAverageSize;
1443  static double WedgeAverageSize;
1444  static double HexAverageSize;
1445 
1446 private:
1447  vtkMeshQuality(const vtkMeshQuality&) = delete;
1448  void operator=(const vtkMeshQuality&) = delete;
1449 };
1450 
1451 VTK_ABI_NAMESPACE_END
1452 #endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:59
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
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:38
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.
vtkTypeBool Volume
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.
vtkTypeBool CompatibilityMode
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.
virtual void Modified()
Update the modification time for this object.
@ mode
Definition: vtkX3D.h:247
int vtkTypeBool
Definition: vtkABI.h:64
#define VTK_DEPRECATED_IN_9_2_0(reason)