VTK  9.2.6
vtkMeshQuality.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkMeshQuality.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14  Copyright 2003-2008 Sandia Corporation.
15  Under the terms of Contract DE-AC04-94AL85000, there is a non-exclusive
16  license for use of this work by or on behalf of the
17  U.S. Government. Redistribution and use in source and binary forms, with
18  or without modification, are permitted provided that this Notice and any
19  statement of authorship are reproduced on all copies.
20 
21  Contact: dcthomp@sandia.gov,pppebay@sandia.gov
22 
23 =========================================================================*/
103 #ifndef vtkMeshQuality_h
104 #define vtkMeshQuality_h
105 
106 #include "vtkDataSetAlgorithm.h"
107 #include "vtkDeprecation.h" // For deprecation
108 #include "vtkFiltersVerdictModule.h" // For export macro
109 
110 class vtkCell;
111 class vtkDataArray;
112 class vtkDoubleArray;
113 class vtkMeshQualityFunctor;
114 
115 class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
116 {
117 private:
118  friend class vtkMeshQualityFunctor;
119 
120 public:
121  void PrintSelf(ostream& os, vtkIndent indent) override;
123  static vtkMeshQuality* New();
124 
126 
132  vtkSetMacro(SaveCellQuality, vtkTypeBool);
133  vtkGetMacro(SaveCellQuality, vtkTypeBool);
134  vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
136 
138 
145  vtkSetMacro(LinearApproximation, bool);
146  vtkGetMacro(LinearApproximation, bool);
147  vtkBooleanMacro(LinearApproximation, bool);
149 
154  {
155  AREA = 28,
156  ASPECT_FROBENIUS = 3,
157  ASPECT_GAMMA = 27,
158  ASPECT_RATIO = 1,
159  COLLAPSE_RATIO = 7,
160  CONDITION = 9,
161  DIAGONAL = 21,
162  DIMENSION = 22,
163  DISTORTION = 15,
164  EDGE_RATIO = 0,
165  EQUIANGLE_SKEW = 29,
166  EQUIVOLUME_SKEW = 30,
167  JACOBIAN = 25,
168  MAX_ANGLE = 8,
169  MAX_ASPECT_FROBENIUS = 5,
170  MAX_EDGE_RATIO = 16,
171  MAX_STRETCH = 31,
172  MEAN_ASPECT_FROBENIUS = 32,
173  MEAN_RATIO = 33,
174  MED_ASPECT_FROBENIUS = 4,
175  MIN_ANGLE = 6,
176  NODAL_JACOBIAN_RATIO = 34,
177  NORMALIZED_INRADIUS = 35,
178  ODDY = 23,
179  RADIUS_RATIO = 2,
180  RELATIVE_SIZE_SQUARED = 12,
181  SCALED_JACOBIAN = 10,
182  SHAPE = 13,
183  SHAPE_AND_SIZE = 14,
184  SHEAR = 11,
185  SHEAR_AND_SIZE = 24,
186  SKEW = 17,
187  SQUISH_INDEX = 36,
188  STRETCH = 20,
189  TAPER = 18,
190  VOLUME = 19,
191  WARPAGE = 26,
192  TOTAL_QUALITY_MEASURE_TYPES = 37,
193  NONE = TOTAL_QUALITY_MEASURE_TYPES
194  };
195 
199  static const char* QualityMeasureNames[];
200 
202 
209  vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
210  virtual void SetTriangleQualityMeasure(int measure)
211  {
212  this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
213  }
214  vtkGetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
216  {
217  this->SetTriangleQualityMeasure(QualityMeasureTypes::AREA);
218  }
220  {
221  this->SetTriangleQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
222  }
224  {
225  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
226  }
228  {
229  this->SetTriangleQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
230  }
232  {
233  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
234  }
236  {
237  this->SetTriangleQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
238  }
240  {
241  this->SetTriangleQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
242  }
244  {
245  this->SetTriangleQualityMeasure(QualityMeasureTypes::CONDITION);
246  }
248  {
249  this->SetTriangleQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
250  }
252  {
253  this->SetTriangleQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
254  }
256  {
257  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE);
258  }
260  {
261  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
262  }
264  {
265  this->SetTriangleQualityMeasure(QualityMeasureTypes::DISTORTION);
266  }
268  {
269  this->SetTriangleQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
270  }
272  {
273  this->SetTriangleQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
274  }
276 
278 
290  vtkSetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
291  virtual void SetQuadQualityMeasure(int measure)
292  {
293  this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
294  }
295  vtkGetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
297  {
298  this->SetQuadQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
299  }
301  {
302  this->SetQuadQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
303  }
305  {
306  this->SetQuadQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
307  }
309  {
310  this->SetQuadQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
311  }
313  {
314  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
315  }
317  {
318  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
319  }
320  void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(QualityMeasureTypes::SKEW); }
321  void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(QualityMeasureTypes::TAPER); }
323  {
324  this->SetQuadQualityMeasure(QualityMeasureTypes::WARPAGE);
325  }
326  void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(QualityMeasureTypes::AREA); }
328  {
329  this->SetQuadQualityMeasure(QualityMeasureTypes::STRETCH);
330  }
332  {
333  this->SetQuadQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
334  }
336  {
337  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
338  }
339  void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(QualityMeasureTypes::ODDY); }
341  {
342  this->SetQuadQualityMeasure(QualityMeasureTypes::CONDITION);
343  }
345  {
346  this->SetQuadQualityMeasure(QualityMeasureTypes::JACOBIAN);
347  }
349  {
350  this->SetQuadQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
351  }
352  void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR); }
353  void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE); }
355  {
356  this->SetQuadQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
357  }
359  {
360  this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
361  }
363  {
364  this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
365  }
367  {
368  this->SetQuadQualityMeasure(QualityMeasureTypes::DISTORTION);
369  }
371  {
372  this->SetQuadQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
373  }
375 
377 
385  virtual void SetTetQualityMeasure(int measure)
386  {
387  this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
388  }
391  {
392  this->SetTetQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
393  }
395  {
396  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
397  }
399  {
400  this->SetTetQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
401  }
403  {
404  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
405  }
407  {
408  this->SetTetQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
409  }
411  {
412  this->SetTetQualityMeasure(QualityMeasureTypes::COLLAPSE_RATIO);
413  }
415  {
416  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_GAMMA);
417  }
418  void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(QualityMeasureTypes::VOLUME); }
420  {
421  this->SetTetQualityMeasure(QualityMeasureTypes::CONDITION);
422  }
424  {
425  this->SetTetQualityMeasure(QualityMeasureTypes::JACOBIAN);
426  }
428  {
429  this->SetTetQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
430  }
431  void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE); }
433  {
434  this->SetTetQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
435  }
437  {
438  this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
439  }
441  {
442  this->SetTetQualityMeasure(QualityMeasureTypes::DISTORTION);
443  }
445  {
446  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
447  }
449  {
450  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIVOLUME_SKEW);
451  }
453  {
454  this->SetTetQualityMeasure(QualityMeasureTypes::MEAN_RATIO);
455  }
457  {
458  this->SetTetQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
459  }
461  {
462  this->SetTetQualityMeasure(QualityMeasureTypes::SQUISH_INDEX);
463  }
465 
467 
472  vtkSetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
473  virtual void SetPyramidQualityMeasure(int measure)
474  {
475  this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
476  }
477  vtkGetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
479  {
480  this->SetPyramidQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
481  }
483  {
484  this->SetPyramidQualityMeasure(QualityMeasureTypes::JACOBIAN);
485  }
487  {
488  this->SetPyramidQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
489  }
491  {
492  this->SetPyramidQualityMeasure(QualityMeasureTypes::SHAPE);
493  }
495  {
496  this->SetPyramidQualityMeasure(QualityMeasureTypes::VOLUME);
497  }
499 
501 
507  vtkSetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
508  virtual void SetWedgeQualityMeasure(int measure)
509  {
510  this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
511  }
512  vtkGetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
514  {
515  this->SetWedgeQualityMeasure(QualityMeasureTypes::CONDITION);
516  }
518  {
519  this->SetWedgeQualityMeasure(QualityMeasureTypes::DISTORTION);
520  }
522  {
523  this->SetWedgeQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
524  }
526  {
527  this->SetWedgeQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
528  }
530  {
531  this->SetWedgeQualityMeasure(QualityMeasureTypes::JACOBIAN);
532  }
534  {
535  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
536  }
538  {
539  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_STRETCH);
540  }
542  {
543  this->SetWedgeQualityMeasure(QualityMeasureTypes::MEAN_ASPECT_FROBENIUS);
544  }
546  {
547  this->SetWedgeQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
548  }
549  void SetWedgeQualityMeasureToShape() { this->SetWedgeQualityMeasure(QualityMeasureTypes::SHAPE); }
551  {
552  this->SetWedgeQualityMeasure(QualityMeasureTypes::VOLUME);
553  }
555 
557 
566  virtual void SetHexQualityMeasure(int measure)
567  {
568  this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
569  }
572  {
573  this->SetHexQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
574  }
576  {
577  this->SetHexQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
578  }
580  {
581  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
582  }
584  {
585  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
586  }
587  void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(QualityMeasureTypes::SKEW); }
588  void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(QualityMeasureTypes::TAPER); }
589  void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(QualityMeasureTypes::VOLUME); }
590  void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(QualityMeasureTypes::STRETCH); }
592  {
593  this->SetHexQualityMeasure(QualityMeasureTypes::DIAGONAL);
594  }
596  {
597  this->SetHexQualityMeasure(QualityMeasureTypes::DIMENSION);
598  }
599  void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(QualityMeasureTypes::ODDY); }
601  {
602  this->SetHexQualityMeasure(QualityMeasureTypes::CONDITION);
603  }
605  {
606  this->SetHexQualityMeasure(QualityMeasureTypes::JACOBIAN);
607  }
609  {
610  this->SetHexQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
611  }
612  void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR); }
613  void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE); }
615  {
616  this->SetHexQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
617  }
619  {
620  this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
621  }
623  {
624  this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
625  }
627  {
628  this->SetHexQualityMeasure(QualityMeasureTypes::DISTORTION);
629  }
631  {
632  this->SetHexQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
633  }
635  {
636  this->SetHexQualityMeasure(QualityMeasureTypes::NODAL_JACOBIAN_RATIO);
637  }
639 
643  static double TriangleArea(vtkCell* cell);
644 
652  static double TriangleEdgeRatio(vtkCell* cell);
653 
661  static double TriangleAspectRatio(vtkCell* cell);
662 
670  static double TriangleRadiusRatio(vtkCell* cell);
671 
681  static double TriangleAspectFrobenius(vtkCell* cell);
682 
686  static double TriangleMinAngle(vtkCell* cell);
687 
691  static double TriangleMaxAngle(vtkCell* cell);
692 
696  static double TriangleCondition(vtkCell* cell);
697 
701  static double TriangleScaledJacobian(vtkCell* cell);
702 
709  static double TriangleRelativeSizeSquared(vtkCell* cell);
710 
714  static double TriangleShape(vtkCell* cell);
715 
722  static double TriangleShapeAndSize(vtkCell* cell);
723 
727  static double TriangleDistortion(vtkCell* cell);
728 
732  static double TriangleEquiangleSkew(vtkCell* cell);
733 
739  static double TriangleNormalizedInradius(vtkCell* cell);
740 
748  static double QuadEdgeRatio(vtkCell* cell);
749 
757  static double QuadAspectRatio(vtkCell* cell);
758 
769  static double QuadRadiusRatio(vtkCell* cell);
770 
780  static double QuadMedAspectFrobenius(vtkCell* cell);
781 
791  static double QuadMaxAspectFrobenius(vtkCell* cell);
792 
796  static double QuadMinAngle(vtkCell* cell);
797 
801  static double QuadMaxEdgeRatio(vtkCell* cell);
802 
808  static double QuadSkew(vtkCell* cell);
809 
814  static double QuadTaper(vtkCell* cell);
815 
821  static double QuadWarpage(vtkCell* cell);
822 
827  static double QuadArea(vtkCell* cell);
828 
833  static double QuadStretch(vtkCell* cell);
834 
838  static double QuadMaxAngle(vtkCell* cell);
839 
845  static double QuadOddy(vtkCell* cell);
846 
852  static double QuadCondition(vtkCell* cell);
853 
859  static double QuadJacobian(vtkCell* cell);
860 
866  static double QuadScaledJacobian(vtkCell* cell);
867 
872  static double QuadShear(vtkCell* cell);
873 
878  static double QuadShape(vtkCell* cell);
879 
888  static double QuadRelativeSizeSquared(vtkCell* cell);
889 
897  static double QuadShapeAndSize(vtkCell* cell);
898 
906  static double QuadShearAndSize(vtkCell* cell);
907 
913  static double QuadDistortion(vtkCell* cell);
914 
918  static double QuadEquiangleSkew(vtkCell* cell);
919 
927  static double TetEdgeRatio(vtkCell* cell);
928 
936  static double TetAspectRatio(vtkCell* cell);
937 
945  static double TetRadiusRatio(vtkCell* cell);
946 
957  static double TetAspectFrobenius(vtkCell* cell);
958 
962  static double TetMinAngle(vtkCell* cell);
963 
970  static double TetCollapseRatio(vtkCell* cell);
971 
977  static double TetAspectGamma(vtkCell* cell);
978 
983  static double TetVolume(vtkCell* cell);
984 
989  static double TetCondition(vtkCell* cell);
990 
995  static double TetJacobian(vtkCell* cell);
996 
1002  static double TetScaledJacobian(vtkCell* cell);
1003 
1008  static double TetShape(vtkCell* cell);
1009 
1018  static double TetRelativeSizeSquared(vtkCell* cell);
1019 
1027  static double TetShapeAndSize(vtkCell* cell);
1028 
1034  static double TetDistortion(vtkCell* cell);
1035 
1039  static double TetEquiangleSkew(vtkCell* cell);
1040 
1044  static double TetEquivolumeSkew(vtkCell* cell);
1045 
1051  static double TetMeanRatio(vtkCell* cell);
1052 
1058  static double TetNormalizedInradius(vtkCell* cell);
1059 
1063  static double TetSquishIndex(vtkCell* cell);
1064 
1068  static double PyramidEquiangleSkew(vtkCell* cell);
1069 
1074  static double PyramidJacobian(vtkCell* cell);
1075 
1080  static double PyramidScaledJacobian(vtkCell* cell);
1081 
1087  static double PyramidShape(vtkCell* cell);
1088 
1092  static double PyramidVolume(vtkCell* cell);
1093 
1098  static double WedgeCondition(vtkCell* cell);
1099 
1104  static double WedgeDistortion(vtkCell* cell);
1105 
1111  static double WedgeEdgeRatio(vtkCell* cell);
1112 
1116  static double WedgeEquiangleSkew(vtkCell* cell);
1117 
1122  static double WedgeJacobian(vtkCell* cell);
1123 
1128  static double WedgeMaxAspectFrobenius(vtkCell* cell);
1129 
1135  static double WedgeMaxStretch(vtkCell* cell);
1136 
1142  static double WedgeMeanAspectFrobenius(vtkCell* cell);
1143 
1153  static double WedgeScaledJacobian(vtkCell* cell);
1154 
1160  static double WedgeShape(vtkCell* cell);
1161 
1165  static double WedgeVolume(vtkCell* cell);
1166 
1174  static double HexEdgeRatio(vtkCell* cell);
1175 
1180  static double HexMedAspectFrobenius(vtkCell* cell);
1181 
1186  static double HexMaxAspectFrobenius(vtkCell* cell);
1187 
1191  static double HexMaxEdgeRatio(vtkCell* cell);
1192 
1198  static double HexSkew(vtkCell* cell);
1199 
1204  static double HexTaper(vtkCell* cell);
1205 
1210  static double HexVolume(vtkCell* cell);
1211 
1216  static double HexStretch(vtkCell* cell);
1217 
1222  static double HexDiagonal(vtkCell* cell);
1223 
1229  static double HexDimension(vtkCell* cell);
1230 
1236  static double HexOddy(vtkCell* cell);
1237 
1242  static double HexCondition(vtkCell* cell);
1243 
1249  static double HexJacobian(vtkCell* cell);
1250 
1256  static double HexScaledJacobian(vtkCell* cell);
1257 
1262  static double HexShear(vtkCell* cell);
1263 
1268  static double HexShape(vtkCell* cell);
1269 
1278  static double HexRelativeSizeSquared(vtkCell* cell);
1279 
1287  static double HexShapeAndSize(vtkCell* cell);
1288 
1296  static double HexShearAndSize(vtkCell* cell);
1297 
1303  static double HexDistortion(vtkCell* cell);
1304 
1308  static double HexEquiangleSkew(vtkCell* cell);
1309 
1314  static double HexNodalJacobianRatio(vtkCell* cell);
1315 
1326  virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1327  vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
1328  vtkBooleanMacro(Ratio, vtkTypeBool);
1329 
1331 
1348  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1349  virtual void SetVolume(vtkTypeBool cv)
1350  {
1351  if (!((cv != 0) ^ (this->Volume != 0)))
1352  {
1353  return;
1354  }
1355  this->Modified();
1356  this->Volume = cv;
1357  if (this->Volume)
1358  {
1359  this->CompatibilityMode = 1;
1360  }
1361  }
1362  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1363  vtkTypeBool GetVolume() { return this->Volume; }
1364  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1365  void VolumeOn()
1366  {
1367  if (!this->Volume)
1368  {
1369  this->Volume = 1;
1370  this->Modified();
1371  }
1372  }
1373  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1374  void VolumeOff()
1375  {
1376  if (this->Volume)
1377  {
1378  this->Volume = 0;
1379  this->Modified();
1380  }
1381  }
1383 
1385 
1412  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1413  virtual void SetCompatibilityMode(vtkTypeBool cm)
1414  {
1415  if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
1416  {
1417  return;
1418  }
1419  this->CompatibilityMode = cm;
1420  this->Modified();
1421  if (this->CompatibilityMode)
1422  {
1423  this->Volume = 1;
1424  this->TetQualityMeasure = QualityMeasureTypes::RADIUS_RATIO;
1425  }
1426  }
1427  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1428  vtkGetMacro(CompatibilityMode, vtkTypeBool);
1429  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1430  void CompatibilityModeOn()
1431  {
1432  if (!this->CompatibilityMode)
1433  {
1434  this->CompatibilityMode = 1;
1435  this->Modified();
1436  }
1437  }
1438  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1439  void CompatibilityModeOff()
1440  {
1441  if (this->CompatibilityMode)
1442  {
1443  this->CompatibilityMode = 0;
1444  this->Modified();
1445  }
1446  }
1448 
1449 protected:
1451  ~vtkMeshQuality() override = default;
1452 
1454 
1463 
1464  using CellQualityType = double (*)(vtkCell*);
1471 
1472  // VTK_DEPRECATED_IN_9_2_0 Those 2 attributes need to be removed, and instance in the code as
1473  // well.
1476 
1477  // Variables used to store the average size (2D: area / 3D: volume)
1478  static double TriangleAverageSize;
1479  static double QuadAverageSize;
1480  static double TetAverageSize;
1481  static double PyramidAverageSize;
1482  static double WedgeAverageSize;
1483  static double HexAverageSize;
1484 
1485 private:
1486  vtkMeshQuality(const vtkMeshQuality&) = delete;
1487  void operator=(const vtkMeshQuality&) = delete;
1488 };
1489 
1490 #endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:150
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
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:119
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:253
int vtkTypeBool
Definition: vtkABI.h:69
#define VTK_DEPRECATED_IN_9_2_0(reason)