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 =========================================================================*/
74 #ifndef vtkMeshQuality_h
75 #define vtkMeshQuality_h
76 
77 #include "vtkDataSetAlgorithm.h"
78 #include "vtkDeprecation.h" // For deprecation
79 #include "vtkFiltersVerdictModule.h" // For export macro
80 
81 class vtkCell;
82 class vtkDataArray;
83 class vtkDoubleArray;
84 class vtkMeshQualityFunctor;
85 
86 class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
87 {
88 private:
89  friend class vtkMeshQualityFunctor;
90 
91 public:
92  void PrintSelf(ostream& os, vtkIndent indent) override;
94  static vtkMeshQuality* New();
95 
97 
103  vtkSetMacro(SaveCellQuality, vtkTypeBool);
104  vtkGetMacro(SaveCellQuality, vtkTypeBool);
105  vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
107 
109 
116  vtkSetMacro(LinearApproximation, bool);
117  vtkGetMacro(LinearApproximation, bool);
118  vtkBooleanMacro(LinearApproximation, bool);
120 
125  {
126  AREA = 28,
127  ASPECT_FROBENIUS = 3,
128  ASPECT_GAMMA = 27,
129  ASPECT_RATIO = 1,
130  COLLAPSE_RATIO = 7,
131  CONDITION = 9,
132  DIAGONAL = 21,
133  DIMENSION = 22,
134  DISTORTION = 15,
135  EDGE_RATIO = 0,
136  EQUIANGLE_SKEW = 29,
137  EQUIVOLUME_SKEW = 30,
138  JACOBIAN = 25,
139  MAX_ANGLE = 8,
140  MAX_ASPECT_FROBENIUS = 5,
141  MAX_EDGE_RATIO = 16,
142  MAX_STRETCH = 31,
143  MEAN_ASPECT_FROBENIUS = 32,
144  MEAN_RATIO = 33,
145  MED_ASPECT_FROBENIUS = 4,
146  MIN_ANGLE = 6,
147  NODAL_JACOBIAN_RATIO = 34,
148  NORMALIZED_INRADIUS = 35,
149  ODDY = 23,
150  RADIUS_RATIO = 2,
151  RELATIVE_SIZE_SQUARED = 12,
152  SCALED_JACOBIAN = 10,
153  SHAPE = 13,
154  SHAPE_AND_SIZE = 14,
155  SHEAR = 11,
156  SHEAR_AND_SIZE = 24,
157  SKEW = 17,
158  SQUISH_INDEX = 36,
159  STRETCH = 20,
160  TAPER = 18,
161  VOLUME = 19,
162  WARPAGE = 26,
163  TOTAL_QUALITY_MEASURE_TYPES = 37,
164  NONE = TOTAL_QUALITY_MEASURE_TYPES
165  };
166 
170  static const char* QualityMeasureNames[];
171 
173 
180  vtkSetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
181  virtual void SetTriangleQualityMeasure(int measure)
182  {
183  this->SetTriangleQualityMeasure(static_cast<QualityMeasureTypes>(measure));
184  }
185  vtkGetEnumMacro(TriangleQualityMeasure, QualityMeasureTypes);
187  {
188  this->SetTriangleQualityMeasure(QualityMeasureTypes::AREA);
189  }
191  {
192  this->SetTriangleQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
193  }
195  {
196  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
197  }
199  {
200  this->SetTriangleQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
201  }
203  {
204  this->SetTriangleQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
205  }
207  {
208  this->SetTriangleQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
209  }
211  {
212  this->SetTriangleQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
213  }
215  {
216  this->SetTriangleQualityMeasure(QualityMeasureTypes::CONDITION);
217  }
219  {
220  this->SetTriangleQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
221  }
223  {
224  this->SetTriangleQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
225  }
227  {
228  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE);
229  }
231  {
232  this->SetTriangleQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
233  }
235  {
236  this->SetTriangleQualityMeasure(QualityMeasureTypes::DISTORTION);
237  }
239  {
240  this->SetTriangleQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
241  }
243  {
244  this->SetTriangleQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
245  }
247 
249 
261  vtkSetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
262  virtual void SetQuadQualityMeasure(int measure)
263  {
264  this->SetQuadQualityMeasure(static_cast<QualityMeasureTypes>(measure));
265  }
266  vtkGetEnumMacro(QuadQualityMeasure, QualityMeasureTypes);
268  {
269  this->SetQuadQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
270  }
272  {
273  this->SetQuadQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
274  }
276  {
277  this->SetQuadQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
278  }
280  {
281  this->SetQuadQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
282  }
284  {
285  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
286  }
288  {
289  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
290  }
291  void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(QualityMeasureTypes::SKEW); }
292  void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(QualityMeasureTypes::TAPER); }
294  {
295  this->SetQuadQualityMeasure(QualityMeasureTypes::WARPAGE);
296  }
297  void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(QualityMeasureTypes::AREA); }
299  {
300  this->SetQuadQualityMeasure(QualityMeasureTypes::STRETCH);
301  }
303  {
304  this->SetQuadQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
305  }
307  {
308  this->SetQuadQualityMeasure(QualityMeasureTypes::MAX_ANGLE);
309  }
310  void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(QualityMeasureTypes::ODDY); }
312  {
313  this->SetQuadQualityMeasure(QualityMeasureTypes::CONDITION);
314  }
316  {
317  this->SetQuadQualityMeasure(QualityMeasureTypes::JACOBIAN);
318  }
320  {
321  this->SetQuadQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
322  }
323  void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR); }
324  void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE); }
326  {
327  this->SetQuadQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
328  }
330  {
331  this->SetQuadQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
332  }
334  {
335  this->SetQuadQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
336  }
338  {
339  this->SetQuadQualityMeasure(QualityMeasureTypes::DISTORTION);
340  }
342  {
343  this->SetQuadQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
344  }
346 
348 
355  vtkSetEnumMacro(TetQualityMeasure, QualityMeasureTypes);
356  virtual void SetTetQualityMeasure(int measure)
357  {
358  this->SetTetQualityMeasure(static_cast<QualityMeasureTypes>(measure));
359  }
360  vtkGetEnumMacro(TetQualityMeasure, QualityMeasureTypes);
362  {
363  this->SetTetQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
364  }
366  {
367  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_RATIO);
368  }
370  {
371  this->SetTetQualityMeasure(QualityMeasureTypes::RADIUS_RATIO);
372  }
374  {
375  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_FROBENIUS);
376  }
378  {
379  this->SetTetQualityMeasure(QualityMeasureTypes::MIN_ANGLE);
380  }
382  {
383  this->SetTetQualityMeasure(QualityMeasureTypes::COLLAPSE_RATIO);
384  }
386  {
387  this->SetTetQualityMeasure(QualityMeasureTypes::ASPECT_GAMMA);
388  }
389  void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(QualityMeasureTypes::VOLUME); }
391  {
392  this->SetTetQualityMeasure(QualityMeasureTypes::CONDITION);
393  }
395  {
396  this->SetTetQualityMeasure(QualityMeasureTypes::JACOBIAN);
397  }
399  {
400  this->SetTetQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
401  }
402  void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE); }
404  {
405  this->SetTetQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
406  }
408  {
409  this->SetTetQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
410  }
412  {
413  this->SetTetQualityMeasure(QualityMeasureTypes::DISTORTION);
414  }
416  {
417  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
418  }
420  {
421  this->SetTetQualityMeasure(QualityMeasureTypes::EQUIVOLUME_SKEW);
422  }
424  {
425  this->SetTetQualityMeasure(QualityMeasureTypes::MEAN_RATIO);
426  }
428  {
429  this->SetTetQualityMeasure(QualityMeasureTypes::NORMALIZED_INRADIUS);
430  }
432  {
433  this->SetTetQualityMeasure(QualityMeasureTypes::SQUISH_INDEX);
434  }
436 
438 
443  vtkSetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
444  virtual void SetPyramidQualityMeasure(int measure)
445  {
446  this->SetPyramidQualityMeasure(static_cast<QualityMeasureTypes>(measure));
447  }
448  vtkGetEnumMacro(PyramidQualityMeasure, QualityMeasureTypes);
450  {
451  this->SetPyramidQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
452  }
454  {
455  this->SetPyramidQualityMeasure(QualityMeasureTypes::JACOBIAN);
456  }
458  {
459  this->SetPyramidQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
460  }
462  {
463  this->SetPyramidQualityMeasure(QualityMeasureTypes::SHAPE);
464  }
466  {
467  this->SetPyramidQualityMeasure(QualityMeasureTypes::VOLUME);
468  }
470 
472 
478  vtkSetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
479  virtual void SetWedgeQualityMeasure(int measure)
480  {
481  this->SetWedgeQualityMeasure(static_cast<QualityMeasureTypes>(measure));
482  }
483  vtkGetEnumMacro(WedgeQualityMeasure, QualityMeasureTypes);
485  {
486  this->SetWedgeQualityMeasure(QualityMeasureTypes::CONDITION);
487  }
489  {
490  this->SetWedgeQualityMeasure(QualityMeasureTypes::DISTORTION);
491  }
493  {
494  this->SetWedgeQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
495  }
497  {
498  this->SetWedgeQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
499  }
501  {
502  this->SetWedgeQualityMeasure(QualityMeasureTypes::JACOBIAN);
503  }
505  {
506  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
507  }
509  {
510  this->SetWedgeQualityMeasure(QualityMeasureTypes::MAX_STRETCH);
511  }
513  {
514  this->SetWedgeQualityMeasure(QualityMeasureTypes::MEAN_ASPECT_FROBENIUS);
515  }
517  {
518  this->SetWedgeQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
519  }
520  void SetWedgeQualityMeasureToShape() { this->SetWedgeQualityMeasure(QualityMeasureTypes::SHAPE); }
522  {
523  this->SetWedgeQualityMeasure(QualityMeasureTypes::VOLUME);
524  }
526 
528 
536  vtkSetEnumMacro(HexQualityMeasure, QualityMeasureTypes);
537  virtual void SetHexQualityMeasure(int measure)
538  {
539  this->SetHexQualityMeasure(static_cast<QualityMeasureTypes>(measure));
540  }
541  vtkGetEnumMacro(HexQualityMeasure, QualityMeasureTypes);
543  {
544  this->SetHexQualityMeasure(QualityMeasureTypes::EDGE_RATIO);
545  }
547  {
548  this->SetHexQualityMeasure(QualityMeasureTypes::MED_ASPECT_FROBENIUS);
549  }
551  {
552  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_ASPECT_FROBENIUS);
553  }
555  {
556  this->SetHexQualityMeasure(QualityMeasureTypes::MAX_EDGE_RATIO);
557  }
558  void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(QualityMeasureTypes::SKEW); }
559  void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(QualityMeasureTypes::TAPER); }
560  void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(QualityMeasureTypes::VOLUME); }
561  void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(QualityMeasureTypes::STRETCH); }
563  {
564  this->SetHexQualityMeasure(QualityMeasureTypes::DIAGONAL);
565  }
567  {
568  this->SetHexQualityMeasure(QualityMeasureTypes::DIMENSION);
569  }
570  void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(QualityMeasureTypes::ODDY); }
572  {
573  this->SetHexQualityMeasure(QualityMeasureTypes::CONDITION);
574  }
576  {
577  this->SetHexQualityMeasure(QualityMeasureTypes::JACOBIAN);
578  }
580  {
581  this->SetHexQualityMeasure(QualityMeasureTypes::SCALED_JACOBIAN);
582  }
583  void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR); }
584  void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE); }
586  {
587  this->SetHexQualityMeasure(QualityMeasureTypes::RELATIVE_SIZE_SQUARED);
588  }
590  {
591  this->SetHexQualityMeasure(QualityMeasureTypes::SHAPE_AND_SIZE);
592  }
594  {
595  this->SetHexQualityMeasure(QualityMeasureTypes::SHEAR_AND_SIZE);
596  }
598  {
599  this->SetHexQualityMeasure(QualityMeasureTypes::DISTORTION);
600  }
602  {
603  this->SetHexQualityMeasure(QualityMeasureTypes::EQUIANGLE_SKEW);
604  }
606  {
607  this->SetHexQualityMeasure(QualityMeasureTypes::NODAL_JACOBIAN_RATIO);
608  }
610 
614  static double TriangleArea(vtkCell* cell);
615 
623  static double TriangleEdgeRatio(vtkCell* cell);
624 
632  static double TriangleAspectRatio(vtkCell* cell);
633 
641  static double TriangleRadiusRatio(vtkCell* cell);
642 
652  static double TriangleAspectFrobenius(vtkCell* cell);
653 
657  static double TriangleMinAngle(vtkCell* cell);
658 
662  static double TriangleMaxAngle(vtkCell* cell);
663 
667  static double TriangleCondition(vtkCell* cell);
668 
672  static double TriangleScaledJacobian(vtkCell* cell);
673 
680  static double TriangleRelativeSizeSquared(vtkCell* cell);
681 
685  static double TriangleShape(vtkCell* cell);
686 
693  static double TriangleShapeAndSize(vtkCell* cell);
694 
698  static double TriangleDistortion(vtkCell* cell);
699 
703  static double TriangleEquiangleSkew(vtkCell* cell);
704 
710  static double TriangleNormalizedInradius(vtkCell* cell);
711 
719  static double QuadEdgeRatio(vtkCell* cell);
720 
728  static double QuadAspectRatio(vtkCell* cell);
729 
740  static double QuadRadiusRatio(vtkCell* cell);
741 
751  static double QuadMedAspectFrobenius(vtkCell* cell);
752 
762  static double QuadMaxAspectFrobenius(vtkCell* cell);
763 
767  static double QuadMinAngle(vtkCell* cell);
768 
772  static double QuadMaxEdgeRatio(vtkCell* cell);
773 
779  static double QuadSkew(vtkCell* cell);
780 
785  static double QuadTaper(vtkCell* cell);
786 
792  static double QuadWarpage(vtkCell* cell);
793 
798  static double QuadArea(vtkCell* cell);
799 
804  static double QuadStretch(vtkCell* cell);
805 
809  static double QuadMaxAngle(vtkCell* cell);
810 
816  static double QuadOddy(vtkCell* cell);
817 
823  static double QuadCondition(vtkCell* cell);
824 
830  static double QuadJacobian(vtkCell* cell);
831 
837  static double QuadScaledJacobian(vtkCell* cell);
838 
843  static double QuadShear(vtkCell* cell);
844 
849  static double QuadShape(vtkCell* cell);
850 
859  static double QuadRelativeSizeSquared(vtkCell* cell);
860 
868  static double QuadShapeAndSize(vtkCell* cell);
869 
877  static double QuadShearAndSize(vtkCell* cell);
878 
884  static double QuadDistortion(vtkCell* cell);
885 
889  static double QuadEquiangleSkew(vtkCell* cell);
890 
898  static double TetEdgeRatio(vtkCell* cell);
899 
907  static double TetAspectRatio(vtkCell* cell);
908 
916  static double TetRadiusRatio(vtkCell* cell);
917 
928  static double TetAspectFrobenius(vtkCell* cell);
929 
933  static double TetMinAngle(vtkCell* cell);
934 
941  static double TetCollapseRatio(vtkCell* cell);
942 
948  static double TetAspectGamma(vtkCell* cell);
949 
954  static double TetVolume(vtkCell* cell);
955 
960  static double TetCondition(vtkCell* cell);
961 
966  static double TetJacobian(vtkCell* cell);
967 
973  static double TetScaledJacobian(vtkCell* cell);
974 
979  static double TetShape(vtkCell* cell);
980 
989  static double TetRelativeSizeSquared(vtkCell* cell);
990 
998  static double TetShapeAndSize(vtkCell* cell);
999 
1005  static double TetDistortion(vtkCell* cell);
1006 
1010  static double TetEquiangleSkew(vtkCell* cell);
1011 
1015  static double TetEquivolumeSkew(vtkCell* cell);
1016 
1022  static double TetMeanRatio(vtkCell* cell);
1023 
1029  static double TetNormalizedInradius(vtkCell* cell);
1030 
1034  static double TetSquishIndex(vtkCell* cell);
1035 
1039  static double PyramidEquiangleSkew(vtkCell* cell);
1040 
1045  static double PyramidJacobian(vtkCell* cell);
1046 
1051  static double PyramidScaledJacobian(vtkCell* cell);
1052 
1058  static double PyramidShape(vtkCell* cell);
1059 
1063  static double PyramidVolume(vtkCell* cell);
1064 
1069  static double WedgeCondition(vtkCell* cell);
1070 
1075  static double WedgeDistortion(vtkCell* cell);
1076 
1082  static double WedgeEdgeRatio(vtkCell* cell);
1083 
1087  static double WedgeEquiangleSkew(vtkCell* cell);
1088 
1093  static double WedgeJacobian(vtkCell* cell);
1094 
1099  static double WedgeMaxAspectFrobenius(vtkCell* cell);
1100 
1106  static double WedgeMaxStretch(vtkCell* cell);
1107 
1113  static double WedgeMeanAspectFrobenius(vtkCell* cell);
1114 
1124  static double WedgeScaledJacobian(vtkCell* cell);
1125 
1131  static double WedgeShape(vtkCell* cell);
1132 
1136  static double WedgeVolume(vtkCell* cell);
1137 
1145  static double HexEdgeRatio(vtkCell* cell);
1146 
1151  static double HexMedAspectFrobenius(vtkCell* cell);
1152 
1157  static double HexMaxAspectFrobenius(vtkCell* cell);
1158 
1162  static double HexMaxEdgeRatio(vtkCell* cell);
1163 
1169  static double HexSkew(vtkCell* cell);
1170 
1175  static double HexTaper(vtkCell* cell);
1176 
1181  static double HexVolume(vtkCell* cell);
1182 
1187  static double HexStretch(vtkCell* cell);
1188 
1193  static double HexDiagonal(vtkCell* cell);
1194 
1200  static double HexDimension(vtkCell* cell);
1201 
1207  static double HexOddy(vtkCell* cell);
1208 
1213  static double HexCondition(vtkCell* cell);
1214 
1220  static double HexJacobian(vtkCell* cell);
1221 
1227  static double HexScaledJacobian(vtkCell* cell);
1228 
1233  static double HexShear(vtkCell* cell);
1234 
1239  static double HexShape(vtkCell* cell);
1240 
1249  static double HexRelativeSizeSquared(vtkCell* cell);
1250 
1258  static double HexShapeAndSize(vtkCell* cell);
1259 
1267  static double HexShearAndSize(vtkCell* cell);
1268 
1274  static double HexDistortion(vtkCell* cell);
1275 
1279  static double HexEquiangleSkew(vtkCell* cell);
1280 
1285  static double HexNodalJacobianRatio(vtkCell* cell);
1286 
1297  virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
1298  vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
1299  vtkBooleanMacro(Ratio, vtkTypeBool);
1300 
1302 
1319  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1320  virtual void SetVolume(vtkTypeBool cv)
1321  {
1322  if (!((cv != 0) ^ (this->Volume != 0)))
1323  {
1324  return;
1325  }
1326  this->Modified();
1327  this->Volume = cv;
1328  if (this->Volume)
1329  {
1330  this->CompatibilityMode = 1;
1331  }
1332  }
1333  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1334  vtkTypeBool GetVolume() { return this->Volume; }
1335  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1336  void VolumeOn()
1337  {
1338  if (!this->Volume)
1339  {
1340  this->Volume = 1;
1341  this->Modified();
1342  }
1343  }
1344  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1345  void VolumeOff()
1346  {
1347  if (this->Volume)
1348  {
1349  this->Volume = 0;
1350  this->Modified();
1351  }
1352  }
1354 
1356 
1383  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1384  virtual void SetCompatibilityMode(vtkTypeBool cm)
1385  {
1386  if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
1387  {
1388  return;
1389  }
1390  this->CompatibilityMode = cm;
1391  this->Modified();
1392  if (this->CompatibilityMode)
1393  {
1394  this->Volume = 1;
1395  this->TetQualityMeasure = QualityMeasureTypes::RADIUS_RATIO;
1396  }
1397  }
1398  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1399  vtkGetMacro(CompatibilityMode, vtkTypeBool);
1400  VTK_DEPRECATED_IN_9_2_0("Deprecating compatibility mode for this filter")
1401  void CompatibilityModeOn()
1402  {
1403  if (!this->CompatibilityMode)
1404  {
1405  this->CompatibilityMode = 1;
1406  this->Modified();
1407  }
1408  }
1409  VTK_DEPRECATED_IN_9_2_0("Part of deprecating compatibility mode for this filter")
1410  void CompatibilityModeOff()
1411  {
1412  if (this->CompatibilityMode)
1413  {
1414  this->CompatibilityMode = 0;
1415  this->Modified();
1416  }
1417  }
1419 
1420 protected:
1421  vtkMeshQuality();
1422  ~vtkMeshQuality() override = default;
1423 
1425 
1434 
1435  using CellQualityType = double (*)(vtkCell*);
1436  CellQualityType GetTriangleQualityMeasureFunction();
1437  CellQualityType GetQuadQualityMeasureFunction();
1438  CellQualityType GetTetQualityMeasureFunction();
1439  CellQualityType GetPyramidQualityMeasureFunction();
1440  CellQualityType GetWedgeQualityMeasureFunction();
1441  CellQualityType GetHexQualityMeasureFunction();
1442 
1443  // VTK_DEPRECATED_IN_9_2_0 Those 2 attributes need to be removed, and instance in the code as
1444  // well.
1447 
1448  // Variables used to store the average size (2D: area / 3D: volume)
1449  static double TriangleAverageSize;
1450  static double QuadAverageSize;
1451  static double TetAverageSize;
1452  static double PyramidAverageSize;
1453  static double WedgeAverageSize;
1454  static double HexAverageSize;
1455 
1456 private:
1457  vtkMeshQuality(const vtkMeshQuality&) = delete;
1458  void operator=(const vtkMeshQuality&) = delete;
1459 };
1460 
1461 #endif // vtkMeshQuality_h
QualityMeasureTypes
Enum which lists the Quality Measures Types.
void SetQuadQualityMeasureToEdgeRatio()
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
void SetQuadQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToMedAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTetQualityMeasureToMinAngle()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetHexQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToRadiusRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
QualityMeasureTypes WedgeQualityMeasure
void SetTriangleQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to function the quality of triangles.
void SetWedgeQualityMeasureToJacobian()
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.
void SetTriangleQualityMeasureToAspectRatio()
Set/Get the particular estimator used to function the quality of triangles.
void SetHexQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetPyramidQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of pyramids.
void SetTetQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTetQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTetQualityMeasureToMeanRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTetQualityMeasureToEquiangleSkew()
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.
void SetQuadQualityMeasureToOddy()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
Store vtkAlgorithm input/output information.
void SetHexQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetWedgeQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToMaxAngle()
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.
void SetWedgeQualityMeasureToMeanAspectFrobenius()
Set/Get the particular estimator used to measure the quality of wedges.
void SetHexQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToMaxEdgeRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetPyramidQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of pyramids.
void SetTetQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of tetrahedra.
#define VTK_DEPRECATED_IN_9_2_0(reason)
void SetTriangleQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to function the quality of triangles.
void SetPyramidQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of pyramids.
void SetTetQualityMeasureToNormalizedInradius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
double(*)(vtkCell *) CellQualityType
void SetTetQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of tetrahedra.
QualityMeasureTypes PyramidQualityMeasure
void SetTetQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTriangleQualityMeasureToDistortion()
Set/Get the particular estimator used to function the quality of triangles.
static double WedgeAverageSize
void SetTriangleQualityMeasureToCondition()
Set/Get the particular estimator used to function the quality of triangles.
virtual void SetPyramidQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of pyramids.
void SetQuadQualityMeasureToMinAngle()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetQuadQualityMeasureToMedAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetQuadQualityMeasureToShearAndSize()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToNodalJacobianRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetHexQualityMeasureToSkew()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
vtkTypeBool Volume
void SetTriangleQualityMeasureToShape()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToAspectRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTriangleQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to function the quality of triangles.
void SetTetQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetWedgeQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of wedges.
void SetTetQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of tetrahedra.
dynamic, self-adjusting array of double
void SetQuadQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetQuadQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
int vtkTypeBool
Definition: vtkABI.h:69
vtkTypeBool SaveCellQuality
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetWedgeQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of wedges.
void SetHexQualityMeasureToDiagonal()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetWedgeQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of wedges.
abstract class to specify cell behavior
Definition: vtkCell.h:60
QualityMeasureTypes TriangleQualityMeasure
void SetQuadQualityMeasureToStretch()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTriangleQualityMeasureToRadiusRatio()
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.
void SetTetQualityMeasureToCollapseRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetPyramidQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of pyramids.
void SetTetQualityMeasureToAspectFrobenius()
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.
a simple class to control print indentation
Definition: vtkIndent.h:39
void SetHexQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTetQualityMeasureToVolume()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetWedgeQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToEquiangleSkew()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTriangleQualityMeasureToMinAngle()
Set/Get the particular estimator used to function the quality of triangles.
QualityMeasureTypes HexQualityMeasure
void SetTetQualityMeasureToRadiusRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
virtual void SetHexQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
void SetQuadQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
virtual void SetQuadQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of quadrilaterals.
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called within ProcessRequest when a request asks the algorithm to do its work.
void SetTriangleQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to function the quality of triangles.
virtual void Modified()
Update the modification time for this object.
void SetWedgeQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToArea()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetHexQualityMeasureToStretch()
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 SetHexQualityMeasureToMaxEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetHexQualityMeasureToShearAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
virtual void SetRatio(vtkTypeBool r)
These methods are deprecated.
void SetHexQualityMeasureToTaper()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetHexQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetTetQualityMeasureToSquishIndex()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetHexQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexAverageSize
vtkTypeBool GetRatio()
static double TetAverageSize
void SetHexQualityMeasureToOddy()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetWedgeQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToSkew()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
virtual void SetTetQualityMeasure(int measure)
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetTetQualityMeasureToAspectRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetWedgeQualityMeasureToMaxStretch()
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToShape()
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.
QualityMeasureTypes TetQualityMeasure
static double QuadAverageSize
void SetTriangleQualityMeasureToMaxAngle()
Set/Get the particular estimator used to function the quality of triangles.
Store zero or more vtkInformation instances.
void SetTetQualityMeasureToAspectGamma()
Set/Get the particular estimator used to measure the quality of tetrahedra.
QualityMeasureTypes QuadQualityMeasure
Calculate functions of quality of the elements of a mesh.
Superclass for algorithms that produce output of the same type as input.
void SetHexQualityMeasureToDimension()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetPyramidQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of pyramids.
vtkTypeBool CompatibilityMode
virtual void SetTriangleQualityMeasure(int measure)
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToTaper()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTetQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static vtkDataSetAlgorithm * New()
void SetTriangleQualityMeasureToNormalizedInradius()
Set/Get the particular estimator used to function the quality of triangles.
static double PyramidAverageSize
void SetWedgeQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of wedges.
void SetQuadQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
void SetTetQualityMeasureToEquivolumeSkew()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetHexQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of hexahedra.