VTK  9.0.1
vtkMeshQuality.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkMeshQuality.h
5  Language: C++
6 
7  Copyright 2003-2006 Sandia Corporation.
8  Under the terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9  license for use of this work by or on behalf of the
10  U.S. Government. Redistribution and use in source and binary forms, with
11  or without modification, are permitted provided that this Notice and any
12  statement of authorship are reproduced on all copies.
13 
14  Contact: dcthomp@sandia.gov,pppebay@sandia.gov
15 
16 =========================================================================*/
62 #ifndef vtkMeshQuality_h
63 #define vtkMeshQuality_h
64 
65 #include "vtkDataSetAlgorithm.h"
66 #include "vtkFiltersVerdictModule.h" // For export macro
67 
68 class vtkCell;
69 class vtkDataArray;
70 
71 #define VTK_QUALITY_EDGE_RATIO 0
72 #define VTK_QUALITY_ASPECT_RATIO 1
73 #define VTK_QUALITY_RADIUS_RATIO 2
74 #define VTK_QUALITY_ASPECT_FROBENIUS 3
75 #define VTK_QUALITY_MED_ASPECT_FROBENIUS 4
76 #define VTK_QUALITY_MAX_ASPECT_FROBENIUS 5
77 #define VTK_QUALITY_MIN_ANGLE 6
78 #define VTK_QUALITY_COLLAPSE_RATIO 7
79 #define VTK_QUALITY_MAX_ANGLE 8
80 #define VTK_QUALITY_CONDITION 9
81 #define VTK_QUALITY_SCALED_JACOBIAN 10
82 #define VTK_QUALITY_SHEAR 11
83 #define VTK_QUALITY_RELATIVE_SIZE_SQUARED 12
84 #define VTK_QUALITY_SHAPE 13
85 #define VTK_QUALITY_SHAPE_AND_SIZE 14
86 #define VTK_QUALITY_DISTORTION 15
87 #define VTK_QUALITY_MAX_EDGE_RATIO 16
88 #define VTK_QUALITY_SKEW 17
89 #define VTK_QUALITY_TAPER 18
90 #define VTK_QUALITY_VOLUME 19
91 #define VTK_QUALITY_STRETCH 20
92 #define VTK_QUALITY_DIAGONAL 21
93 #define VTK_QUALITY_DIMENSION 22
94 #define VTK_QUALITY_ODDY 23
95 #define VTK_QUALITY_SHEAR_AND_SIZE 24
96 #define VTK_QUALITY_JACOBIAN 25
97 #define VTK_QUALITY_WARPAGE 26
98 #define VTK_QUALITY_ASPECT_GAMMA 27
99 #define VTK_QUALITY_AREA 28
100 #define VTK_QUALITY_ASPECT_BETA 29
101 
102 class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
103 {
104 public:
105  void PrintSelf(ostream& os, vtkIndent indent) override;
107  static vtkMeshQuality* New();
108 
110 
116  vtkSetMacro(SaveCellQuality, vtkTypeBool);
117  vtkGetMacro(SaveCellQuality, vtkTypeBool);
118  vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
120 
122 
130  vtkSetMacro(TriangleQualityMeasure, int);
131  vtkGetMacro(TriangleQualityMeasure, int);
132  void SetTriangleQualityMeasureToArea() { this->SetTriangleQualityMeasure(VTK_QUALITY_AREA); }
134  {
135  this->SetTriangleQualityMeasure(VTK_QUALITY_EDGE_RATIO);
136  }
138  {
139  this->SetTriangleQualityMeasure(VTK_QUALITY_ASPECT_RATIO);
140  }
142  {
143  this->SetTriangleQualityMeasure(VTK_QUALITY_RADIUS_RATIO);
144  }
146  {
147  this->SetTriangleQualityMeasure(VTK_QUALITY_ASPECT_FROBENIUS);
148  }
150  {
151  this->SetTriangleQualityMeasure(VTK_QUALITY_MIN_ANGLE);
152  }
154  {
155  this->SetTriangleQualityMeasure(VTK_QUALITY_MAX_ANGLE);
156  }
158  {
159  this->SetTriangleQualityMeasure(VTK_QUALITY_CONDITION);
160  }
162  {
163  this->SetTriangleQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
164  }
166  {
167  this->SetTriangleQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
168  }
169  void SetTriangleQualityMeasureToShape() { this->SetTriangleQualityMeasure(VTK_QUALITY_SHAPE); }
171  {
172  this->SetTriangleQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
173  }
175  {
176  this->SetTriangleQualityMeasure(VTK_QUALITY_DISTORTION);
177  }
179 
181 
196  vtkSetMacro(QuadQualityMeasure, int);
197  vtkGetMacro(QuadQualityMeasure, int);
198  void SetQuadQualityMeasureToEdgeRatio() { this->SetQuadQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
200  {
201  this->SetQuadQualityMeasure(VTK_QUALITY_ASPECT_RATIO);
202  }
204  {
205  this->SetQuadQualityMeasure(VTK_QUALITY_RADIUS_RATIO);
206  }
208  {
209  this->SetQuadQualityMeasure(VTK_QUALITY_MED_ASPECT_FROBENIUS);
210  }
212  {
213  this->SetQuadQualityMeasure(VTK_QUALITY_MAX_ASPECT_FROBENIUS);
214  }
216  {
217  this->SetQuadQualityMeasure(VTK_QUALITY_MAX_EDGE_RATIO);
218  }
219  void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(VTK_QUALITY_SKEW); }
220  void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(VTK_QUALITY_TAPER); }
221  void SetQuadQualityMeasureToWarpage() { this->SetQuadQualityMeasure(VTK_QUALITY_WARPAGE); }
222  void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(VTK_QUALITY_AREA); }
223  void SetQuadQualityMeasureToStretch() { this->SetQuadQualityMeasure(VTK_QUALITY_STRETCH); }
224  void SetQuadQualityMeasureToMinAngle() { this->SetQuadQualityMeasure(VTK_QUALITY_MIN_ANGLE); }
225  void SetQuadQualityMeasureToMaxAngle() { this->SetQuadQualityMeasure(VTK_QUALITY_MAX_ANGLE); }
226  void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(VTK_QUALITY_ODDY); }
227  void SetQuadQualityMeasureToCondition() { this->SetQuadQualityMeasure(VTK_QUALITY_CONDITION); }
228  void SetQuadQualityMeasureToJacobian() { this->SetQuadQualityMeasure(VTK_QUALITY_JACOBIAN); }
230  {
231  this->SetQuadQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
232  }
233  void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(VTK_QUALITY_SHEAR); }
234  void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(VTK_QUALITY_SHAPE); }
236  {
237  this->SetQuadQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
238  }
240  {
241  this->SetQuadQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
242  }
244  {
245  this->SetQuadQualityMeasure(VTK_QUALITY_SHEAR_AND_SIZE);
246  }
247  void SetQuadQualityMeasureToDistortion() { this->SetQuadQualityMeasure(VTK_QUALITY_DISTORTION); }
249 
251 
261  vtkSetMacro(TetQualityMeasure, int);
262  vtkGetMacro(TetQualityMeasure, int);
263  void SetTetQualityMeasureToEdgeRatio() { this->SetTetQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
264  void SetTetQualityMeasureToAspectRatio() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_RATIO); }
265  void SetTetQualityMeasureToRadiusRatio() { this->SetTetQualityMeasure(VTK_QUALITY_RADIUS_RATIO); }
267  {
268  this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_FROBENIUS);
269  }
270  void SetTetQualityMeasureToMinAngle() { this->SetTetQualityMeasure(VTK_QUALITY_MIN_ANGLE); }
272  {
273  this->SetTetQualityMeasure(VTK_QUALITY_COLLAPSE_RATIO);
274  }
275  void SetTetQualityMeasureToAspectBeta() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_BETA); }
276  void SetTetQualityMeasureToAspectGamma() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_GAMMA); }
277  void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(VTK_QUALITY_VOLUME); }
278  void SetTetQualityMeasureToCondition() { this->SetTetQualityMeasure(VTK_QUALITY_CONDITION); }
279  void SetTetQualityMeasureToJacobian() { this->SetTetQualityMeasure(VTK_QUALITY_JACOBIAN); }
281  {
282  this->SetTetQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
283  }
284  void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(VTK_QUALITY_SHAPE); }
286  {
287  this->SetTetQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
288  }
290  {
291  this->SetTetQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
292  }
293  void SetTetQualityMeasureToDistortion() { this->SetTetQualityMeasure(VTK_QUALITY_DISTORTION); }
295 
297 
308  vtkSetMacro(HexQualityMeasure, int);
309  vtkGetMacro(HexQualityMeasure, int);
310  void SetHexQualityMeasureToEdgeRatio() { this->SetHexQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
312  {
313  this->SetHexQualityMeasure(VTK_QUALITY_MED_ASPECT_FROBENIUS);
314  }
316  {
317  this->SetHexQualityMeasure(VTK_QUALITY_MAX_ASPECT_FROBENIUS);
318  }
320  {
321  this->SetHexQualityMeasure(VTK_QUALITY_MAX_EDGE_RATIO);
322  }
323  void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(VTK_QUALITY_SKEW); }
324  void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(VTK_QUALITY_TAPER); }
325  void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(VTK_QUALITY_VOLUME); }
326  void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(VTK_QUALITY_STRETCH); }
327  void SetHexQualityMeasureToDiagonal() { this->SetHexQualityMeasure(VTK_QUALITY_DIAGONAL); }
328  void SetHexQualityMeasureToDimension() { this->SetHexQualityMeasure(VTK_QUALITY_DIMENSION); }
329  void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(VTK_QUALITY_ODDY); }
330  void SetHexQualityMeasureToCondition() { this->SetHexQualityMeasure(VTK_QUALITY_CONDITION); }
331  void SetHexQualityMeasureToJacobian() { this->SetHexQualityMeasure(VTK_QUALITY_JACOBIAN); }
333  {
334  this->SetHexQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
335  }
336  void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(VTK_QUALITY_SHEAR); }
337  void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(VTK_QUALITY_SHAPE); }
339  {
340  this->SetHexQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
341  }
343  {
344  this->SetHexQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
345  }
347  {
348  this->SetHexQualityMeasure(VTK_QUALITY_SHEAR_AND_SIZE);
349  }
350  void SetHexQualityMeasureToDistortion() { this->SetHexQualityMeasure(VTK_QUALITY_DISTORTION); }
352 
359  static double TriangleArea(vtkCell* cell);
360 
371  static double TriangleEdgeRatio(vtkCell* cell);
372 
383  static double TriangleAspectRatio(vtkCell* cell);
384 
395  static double TriangleRadiusRatio(vtkCell* cell);
396 
409  static double TriangleAspectFrobenius(vtkCell* cell);
410 
418  static double TriangleMinAngle(vtkCell* cell);
419 
427  static double TriangleMaxAngle(vtkCell* cell);
428 
436  static double TriangleCondition(vtkCell* cell);
437 
444  static double TriangleScaledJacobian(vtkCell* cell);
445 
452  static double TriangleRelativeSizeSquared(vtkCell* cell);
453 
460  static double TriangleShape(vtkCell* cell);
461 
467  static double TriangleShapeAndSize(vtkCell* cell);
468 
475  static double TriangleDistortion(vtkCell* cell);
476 
487  static double QuadEdgeRatio(vtkCell* cell);
488 
500  static double QuadAspectRatio(vtkCell* cell);
501 
516  static double QuadRadiusRatio(vtkCell* cell);
517 
531  static double QuadMedAspectFrobenius(vtkCell* cell);
532 
546  static double QuadMaxAspectFrobenius(vtkCell* cell);
547 
555  static double QuadMinAngle(vtkCell* cell);
556 
557  static double QuadMaxEdgeRatios(vtkCell* cell);
558  static double QuadSkew(vtkCell* cell);
559  static double QuadTaper(vtkCell* cell);
560  static double QuadWarpage(vtkCell* cell);
561  static double QuadArea(vtkCell* cell);
562  static double QuadStretch(vtkCell* cell);
563  static double QuadMaxAngle(vtkCell* cell);
564  static double QuadOddy(vtkCell* cell);
565  static double QuadCondition(vtkCell* cell);
566  static double QuadJacobian(vtkCell* cell);
567  static double QuadScaledJacobian(vtkCell* cell);
568  static double QuadShear(vtkCell* cell);
569  static double QuadShape(vtkCell* cell);
570  static double QuadRelativeSizeSquared(vtkCell* cell);
571  static double QuadShapeAndSize(vtkCell* cell);
572  static double QuadShearAndSize(vtkCell* cell);
573  static double QuadDistortion(vtkCell* cell);
574 
585  static double TetEdgeRatio(vtkCell* cell);
586 
597  static double TetAspectRatio(vtkCell* cell);
598 
609  static double TetRadiusRatio(vtkCell* cell);
610 
624  static double TetAspectFrobenius(vtkCell* cell);
625 
633  static double TetMinAngle(vtkCell* cell);
634 
636 
645  static double TetCollapseRatio(vtkCell* cell);
646  static double TetAspectBeta(vtkCell* cell);
647  static double TetAspectGamma(vtkCell* cell);
648  static double TetVolume(vtkCell* cell);
649  static double TetCondition(vtkCell* cell);
650  static double TetJacobian(vtkCell* cell);
651  static double TetScaledJacobian(vtkCell* cell);
652  static double TetShape(vtkCell* cell);
653  static double TetRelativeSizeSquared(vtkCell* cell);
654  static double TetShapeandSize(vtkCell* cell);
655  static double TetDistortion(vtkCell* cell);
657 
668  static double HexEdgeRatio(vtkCell* cell);
669 
678  static double HexMedAspectFrobenius(vtkCell* cell);
679 
681 
689  static double HexMaxAspectFrobenius(vtkCell* cell);
690  static double HexMaxEdgeRatio(vtkCell* cell);
691  static double HexSkew(vtkCell* cell);
692  static double HexTaper(vtkCell* cell);
693  static double HexVolume(vtkCell* cell);
694  static double HexStretch(vtkCell* cell);
695  static double HexDiagonal(vtkCell* cell);
696  static double HexDimension(vtkCell* cell);
697  static double HexOddy(vtkCell* cell);
698  static double HexCondition(vtkCell* cell);
699  static double HexJacobian(vtkCell* cell);
700  static double HexScaledJacobian(vtkCell* cell);
701  static double HexShear(vtkCell* cell);
702  static double HexShape(vtkCell* cell);
703  static double HexRelativeSizeSquared(vtkCell* cell);
704  static double HexShapeAndSize(vtkCell* cell);
705  static double HexShearAndSize(vtkCell* cell);
706  static double HexDistortion(vtkCell* cell);
708 
719  virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
720  vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
721  vtkBooleanMacro(Ratio, vtkTypeBool);
722 
724 
741  virtual void SetVolume(vtkTypeBool cv)
742  {
743  if (!((cv != 0) ^ (this->Volume != 0)))
744  {
745  return;
746  }
747  this->Modified();
748  this->Volume = cv;
749  if (this->Volume)
750  {
751  this->CompatibilityModeOn();
752  }
753  }
754  vtkTypeBool GetVolume() { return this->Volume; }
755  vtkBooleanMacro(Volume, vtkTypeBool);
757 
759 
787  {
788  if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
789  {
790  return;
791  }
792  this->CompatibilityMode = cm;
793  this->Modified();
794  if (this->CompatibilityMode)
795  {
796  this->Volume = 1;
797  this->TetQualityMeasure = VTK_QUALITY_RADIUS_RATIO;
798  }
799  }
800  vtkGetMacro(CompatibilityMode, vtkTypeBool);
801  vtkBooleanMacro(CompatibilityMode, vtkTypeBool);
803 
804 protected:
806  ~vtkMeshQuality() override;
807 
809 
813  static int GetCurrentTriangleNormal(double point[3], double normal[3]);
814 
820 
823 
825  static double CurrentTriNormal[3];
826 
827 private:
828  vtkMeshQuality(const vtkMeshQuality&) = delete;
829  void operator=(const vtkMeshQuality&) = delete;
830 };
831 
832 #endif // vtkMeshQuality_h
vtkMeshQuality::QuadWarpage
static double QuadWarpage(vtkCell *cell)
vtkMeshQuality::TetCondition
static double TetCondition(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToScaledJacobian
void SetTetQualityMeasureToScaledJacobian()
Definition: vtkMeshQuality.h:280
VTK_QUALITY_AREA
#define VTK_QUALITY_AREA
Definition: vtkMeshQuality.h:99
vtkMeshQuality::SetTriangleQualityMeasureToAspectRatio
void SetTriangleQualityMeasureToAspectRatio()
Definition: vtkMeshQuality.h:137
vtkMeshQuality::SetTriangleQualityMeasureToCondition
void SetTriangleQualityMeasureToCondition()
Definition: vtkMeshQuality.h:157
vtkMeshQuality::HexEdgeRatio
static double HexEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a hexahedron.
VTK_QUALITY_ODDY
#define VTK_QUALITY_ODDY
Definition: vtkMeshQuality.h:94
VTK_QUALITY_SHAPE
#define VTK_QUALITY_SHAPE
Definition: vtkMeshQuality.h:84
vtkMeshQuality::QuadJacobian
static double QuadJacobian(vtkCell *cell)
vtkMeshQuality::HexShapeAndSize
static double HexShapeAndSize(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToMinAngle
void SetTetQualityMeasureToMinAngle()
Definition: vtkMeshQuality.h:270
vtkMeshQuality::SetTetQualityMeasureToRelativeSizeSquared
void SetTetQualityMeasureToRelativeSizeSquared()
Definition: vtkMeshQuality.h:285
vtkMeshQuality::TetScaledJacobian
static double TetScaledJacobian(vtkCell *cell)
vtkMeshQuality::QuadAspectRatio
static double QuadAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a planar quadrilateral.
VTK_QUALITY_ASPECT_RATIO
#define VTK_QUALITY_ASPECT_RATIO
Definition: vtkMeshQuality.h:72
vtkMeshQuality::SetQuadQualityMeasureToTaper
void SetQuadQualityMeasureToTaper()
Definition: vtkMeshQuality.h:220
VTK_QUALITY_SHEAR
#define VTK_QUALITY_SHEAR
Definition: vtkMeshQuality.h:82
VTK_QUALITY_MAX_ANGLE
#define VTK_QUALITY_MAX_ANGLE
Definition: vtkMeshQuality.h:79
vtkMeshQuality::TriangleQualityMeasure
int TriangleQualityMeasure
Definition: vtkMeshQuality.h:816
vtkMeshQuality::SetVolume
virtual void SetVolume(vtkTypeBool cv)
These methods are deprecated.
Definition: vtkMeshQuality.h:741
vtkMeshQuality::SetQuadQualityMeasureToRelativeSizeSquared
void SetQuadQualityMeasureToRelativeSizeSquared()
Definition: vtkMeshQuality.h:235
VTK_QUALITY_DIAGONAL
#define VTK_QUALITY_DIAGONAL
Definition: vtkMeshQuality.h:92
vtkMeshQuality::SetTriangleQualityMeasureToArea
void SetTriangleQualityMeasureToArea()
Definition: vtkMeshQuality.h:132
vtkMeshQuality::SetTetQualityMeasureToEdgeRatio
void SetTetQualityMeasureToEdgeRatio()
Definition: vtkMeshQuality.h:263
vtkMeshQuality::HexMedAspectFrobenius
static double HexMedAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the average Frobenius aspect of the 8 corner tetrahedra o...
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:36
VTK_QUALITY_CONDITION
#define VTK_QUALITY_CONDITION
Definition: vtkMeshQuality.h:80
vtkMeshQuality::QuadMaxAspectFrobenius
static double QuadMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 4 corner triangles of...
VTK_QUALITY_RADIUS_RATIO
#define VTK_QUALITY_RADIUS_RATIO
Definition: vtkMeshQuality.h:73
vtkMeshQuality::QuadOddy
static double QuadOddy(vtkCell *cell)
vtkMeshQuality::SetTriangleQualityMeasureToRadiusRatio
void SetTriangleQualityMeasureToRadiusRatio()
Definition: vtkMeshQuality.h:141
vtkMeshQuality::SetQuadQualityMeasureToAspectRatio
void SetQuadQualityMeasureToAspectRatio()
Definition: vtkMeshQuality.h:199
vtkMeshQuality::SetHexQualityMeasureToScaledJacobian
void SetHexQualityMeasureToScaledJacobian()
Definition: vtkMeshQuality.h:332
vtkObject::Modified
virtual void Modified()
Update the modification time for this object.
vtkMeshQuality::TriangleRelativeSizeSquared
static double TriangleRelativeSizeSquared(vtkCell *cell)
This is a static function used to calculate the square of the relative size of a triangle.
vtkMeshQuality::SetHexQualityMeasureToEdgeRatio
void SetHexQualityMeasureToEdgeRatio()
Definition: vtkMeshQuality.h:310
vtkMeshQuality::TriangleEdgeRatio
static double TriangleEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a triangle.
vtkMeshQuality::QuadMedAspectFrobenius
static double QuadMedAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the average Frobenius aspect of the 4 corner triangles of...
VTK_QUALITY_VOLUME
#define VTK_QUALITY_VOLUME
Definition: vtkMeshQuality.h:90
vtkMeshQuality::SetQuadQualityMeasureToWarpage
void SetQuadQualityMeasureToWarpage()
Definition: vtkMeshQuality.h:221
vtkMeshQuality::TetAspectGamma
static double TetAspectGamma(vtkCell *cell)
vtkMeshQuality::TriangleScaledJacobian
static double TriangleScaledJacobian(vtkCell *cell)
This is a static function used to calculate the scaled Jacobian of a triangle.
vtkMeshQuality::QuadShearAndSize
static double QuadShearAndSize(vtkCell *cell)
VTK_QUALITY_DISTORTION
#define VTK_QUALITY_DISTORTION
Definition: vtkMeshQuality.h:86
vtkMeshQuality::QuadMaxAngle
static double QuadMaxAngle(vtkCell *cell)
vtkMeshQuality::SetHexQualityMeasureToShearAndSize
void SetHexQualityMeasureToShearAndSize()
Definition: vtkMeshQuality.h:346
vtkMeshQuality::SetQuadQualityMeasureToShapeAndSize
void SetQuadQualityMeasureToShapeAndSize()
Definition: vtkMeshQuality.h:239
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
vtkMeshQuality::HexQualityMeasure
int HexQualityMeasure
Definition: vtkMeshQuality.h:819
VTK_QUALITY_MIN_ANGLE
#define VTK_QUALITY_MIN_ANGLE
Definition: vtkMeshQuality.h:77
vtkMeshQuality
Calculate functions of quality of the elements of a mesh.
Definition: vtkMeshQuality.h:103
vtkMeshQuality::TriangleCondition
static double TriangleCondition(vtkCell *cell)
This is a static function used to calculate the condition number of a triangle.
vtkMeshQuality::SetQuadQualityMeasureToScaledJacobian
void SetQuadQualityMeasureToScaledJacobian()
Definition: vtkMeshQuality.h:229
vtkMeshQuality::TetAspectFrobenius
static double TetAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the Frobenius condition number of the transformation matr...
VTK_QUALITY_MED_ASPECT_FROBENIUS
#define VTK_QUALITY_MED_ASPECT_FROBENIUS
Definition: vtkMeshQuality.h:75
vtkMeshQuality::TetJacobian
static double TetJacobian(vtkCell *cell)
vtkMeshQuality::SetQuadQualityMeasureToMaxAngle
void SetQuadQualityMeasureToMaxAngle()
Definition: vtkMeshQuality.h:225
vtkMeshQuality::SetTetQualityMeasureToCondition
void SetTetQualityMeasureToCondition()
Definition: vtkMeshQuality.h:278
vtkMeshQuality::HexStretch
static double HexStretch(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToAspectBeta
void SetTetQualityMeasureToAspectBeta()
Definition: vtkMeshQuality.h:275
vtkMeshQuality::SetQuadQualityMeasureToMedAspectFrobenius
void SetQuadQualityMeasureToMedAspectFrobenius()
Definition: vtkMeshQuality.h:207
VTK_QUALITY_SKEW
#define VTK_QUALITY_SKEW
Definition: vtkMeshQuality.h:88
vtkMeshQuality::SetHexQualityMeasureToShapeAndSize
void SetHexQualityMeasureToShapeAndSize()
Definition: vtkMeshQuality.h:342
vtkMeshQuality::SetTriangleQualityMeasureToEdgeRatio
void SetTriangleQualityMeasureToEdgeRatio()
Definition: vtkMeshQuality.h:133
vtkMeshQuality::SetTriangleQualityMeasureToDistortion
void SetTriangleQualityMeasureToDistortion()
Definition: vtkMeshQuality.h:174
vtkMeshQuality::SetTriangleQualityMeasureToScaledJacobian
void SetTriangleQualityMeasureToScaledJacobian()
Definition: vtkMeshQuality.h:161
vtkMeshQuality::SetQuadQualityMeasureToMaxAspectFrobenius
void SetQuadQualityMeasureToMaxAspectFrobenius()
Definition: vtkMeshQuality.h:211
VTK_QUALITY_WARPAGE
#define VTK_QUALITY_WARPAGE
Definition: vtkMeshQuality.h:97
vtkMeshQuality::SetHexQualityMeasureToTaper
void SetHexQualityMeasureToTaper()
Definition: vtkMeshQuality.h:324
vtkDataSetAlgorithm
Superclass for algorithms that produce output of the same type as input.
Definition: vtkDataSetAlgorithm.h:49
VTK_QUALITY_ASPECT_FROBENIUS
#define VTK_QUALITY_ASPECT_FROBENIUS
Definition: vtkMeshQuality.h:74
vtkMeshQuality::SetHexQualityMeasureToSkew
void SetHexQualityMeasureToSkew()
Definition: vtkMeshQuality.h:323
vtkMeshQuality::SetHexQualityMeasureToMaxEdgeRatios
void SetHexQualityMeasureToMaxEdgeRatios()
Definition: vtkMeshQuality.h:319
vtkMeshQuality::SetHexQualityMeasureToDistortion
void SetHexQualityMeasureToDistortion()
Definition: vtkMeshQuality.h:350
vtkMeshQuality::QuadShape
static double QuadShape(vtkCell *cell)
vtkMeshQuality::SetQuadQualityMeasureToRadiusRatio
void SetQuadQualityMeasureToRadiusRatio()
Definition: vtkMeshQuality.h:203
vtkMeshQuality::SetHexQualityMeasureToMaxAspectFrobenius
void SetHexQualityMeasureToMaxAspectFrobenius()
Definition: vtkMeshQuality.h:315
vtkMeshQuality::HexTaper
static double HexTaper(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToRadiusRatio
void SetTetQualityMeasureToRadiusRatio()
Definition: vtkMeshQuality.h:265
vtkX3D::point
@ point
Definition: vtkX3D.h:242
vtkMeshQuality::SetQuadQualityMeasureToShearAndSize
void SetQuadQualityMeasureToShearAndSize()
Definition: vtkMeshQuality.h:243
vtkMeshQuality::QuadCondition
static double QuadCondition(vtkCell *cell)
vtkMeshQuality::SetHexQualityMeasureToJacobian
void SetHexQualityMeasureToJacobian()
Definition: vtkMeshQuality.h:331
VTK_QUALITY_STRETCH
#define VTK_QUALITY_STRETCH
Definition: vtkMeshQuality.h:91
vtkMeshQuality::SetHexQualityMeasureToDiagonal
void SetHexQualityMeasureToDiagonal()
Definition: vtkMeshQuality.h:327
vtkMeshQuality::GetVolume
vtkTypeBool GetVolume()
Definition: vtkMeshQuality.h:754
VTK_QUALITY_JACOBIAN
#define VTK_QUALITY_JACOBIAN
Definition: vtkMeshQuality.h:96
vtkMeshQuality::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkMeshQuality::SetQuadQualityMeasureToStretch
void SetQuadQualityMeasureToStretch()
Definition: vtkMeshQuality.h:223
vtkMeshQuality::SetTetQualityMeasureToAspectFrobenius
void SetTetQualityMeasureToAspectFrobenius()
Definition: vtkMeshQuality.h:266
vtkMeshQuality::GetRatio
vtkTypeBool GetRatio()
Definition: vtkMeshQuality.h:720
vtkMeshQuality::SetTetQualityMeasureToShapeAndSize
void SetTetQualityMeasureToShapeAndSize()
Definition: vtkMeshQuality.h:289
VTK_QUALITY_SHEAR_AND_SIZE
#define VTK_QUALITY_SHEAR_AND_SIZE
Definition: vtkMeshQuality.h:95
vtkMeshQuality::SaveCellQuality
vtkTypeBool SaveCellQuality
Definition: vtkMeshQuality.h:815
vtkMeshQuality::SetHexQualityMeasureToVolume
void SetHexQualityMeasureToVolume()
Definition: vtkMeshQuality.h:325
vtkMeshQuality::TetCollapseRatio
static double TetCollapseRatio(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
vtkMeshQuality::SetRatio
virtual void SetRatio(vtkTypeBool r)
These methods are deprecated.
Definition: vtkMeshQuality.h:719
vtkMeshQuality::SetQuadQualityMeasureToShear
void SetQuadQualityMeasureToShear()
Definition: vtkMeshQuality.h:233
vtkMeshQuality::TetMinAngle
static double TetMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) dihedral angle of a tetrahedron...
vtkMeshQuality::SetTetQualityMeasureToJacobian
void SetTetQualityMeasureToJacobian()
Definition: vtkMeshQuality.h:279
vtkMeshQuality::TriangleArea
static double TriangleArea(vtkCell *cell)
This is a static function used to calculate the area of a triangle.
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:57
vtkMeshQuality::SetQuadQualityMeasureToJacobian
void SetQuadQualityMeasureToJacobian()
Definition: vtkMeshQuality.h:228
vtkMeshQuality::SetHexQualityMeasureToShear
void SetHexQualityMeasureToShear()
Definition: vtkMeshQuality.h:336
vtkMeshQuality::SetHexQualityMeasureToOddy
void SetHexQualityMeasureToOddy()
Definition: vtkMeshQuality.h:329
vtkMeshQuality::HexShape
static double HexShape(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToCollapseRatio
void SetTetQualityMeasureToCollapseRatio()
Definition: vtkMeshQuality.h:271
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:34
vtkMeshQuality::TetEdgeRatio
static double TetEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a tetrahedron.
vtkMeshQuality::HexMaxAspectFrobenius
static double HexMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
vtkMeshQuality::QuadArea
static double QuadArea(vtkCell *cell)
vtkMeshQuality::Volume
vtkTypeBool Volume
Definition: vtkMeshQuality.h:822
vtkMeshQuality::SetTetQualityMeasureToVolume
void SetTetQualityMeasureToVolume()
Definition: vtkMeshQuality.h:277
vtkMeshQuality::SetTriangleQualityMeasureToMaxAngle
void SetTriangleQualityMeasureToMaxAngle()
Definition: vtkMeshQuality.h:153
vtkMeshQuality::SetQuadQualityMeasureToMinAngle
void SetQuadQualityMeasureToMinAngle()
Definition: vtkMeshQuality.h:224
vtkMeshQuality::SetQuadQualityMeasureToMaxEdgeRatios
void SetQuadQualityMeasureToMaxEdgeRatios()
Definition: vtkMeshQuality.h:215
vtkMeshQuality::TriangleMaxAngle
static double TriangleMaxAngle(vtkCell *cell)
This is a static function used to calculate the maximal (nonoriented) angle of a triangle,...
vtkMeshQuality::QuadTaper
static double QuadTaper(vtkCell *cell)
vtkMeshQuality::HexShearAndSize
static double HexShearAndSize(vtkCell *cell)
vtkMeshQuality::QuadQualityMeasure
int QuadQualityMeasure
Definition: vtkMeshQuality.h:817
vtkMeshQuality::SetHexQualityMeasureToStretch
void SetHexQualityMeasureToStretch()
Definition: vtkMeshQuality.h:326
vtkMeshQuality::TetQualityMeasure
int TetQualityMeasure
Definition: vtkMeshQuality.h:818
vtkMeshQuality::SetTriangleQualityMeasureToShape
void SetTriangleQualityMeasureToShape()
Definition: vtkMeshQuality.h:169
vtkMeshQuality::HexDimension
static double HexDimension(vtkCell *cell)
vtkMeshQuality::QuadStretch
static double QuadStretch(vtkCell *cell)
vtkMeshQuality::QuadRelativeSizeSquared
static double QuadRelativeSizeSquared(vtkCell *cell)
vtkMeshQuality::SetTriangleQualityMeasureToRelativeSizeSquared
void SetTriangleQualityMeasureToRelativeSizeSquared()
Definition: vtkMeshQuality.h:165
vtkMeshQuality::SetHexQualityMeasureToCondition
void SetHexQualityMeasureToCondition()
Definition: vtkMeshQuality.h:330
VTK_QUALITY_SCALED_JACOBIAN
#define VTK_QUALITY_SCALED_JACOBIAN
Definition: vtkMeshQuality.h:81
vtkMeshQuality::QuadShear
static double QuadShear(vtkCell *cell)
vtkMeshQuality::SetCompatibilityMode
virtual void SetCompatibilityMode(vtkTypeBool cm)
CompatibilityMode governs whether, when both a quality function and cell volume are to be stored as c...
Definition: vtkMeshQuality.h:786
vtkMeshQuality::SetTriangleQualityMeasureToShapeAndSize
void SetTriangleQualityMeasureToShapeAndSize()
Definition: vtkMeshQuality.h:170
vtkMeshQuality::TriangleRadiusRatio
static double TriangleRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a triangle.
vtkMeshQuality::CompatibilityMode
vtkTypeBool CompatibilityMode
Definition: vtkMeshQuality.h:821
vtkMeshQuality::HexScaledJacobian
static double HexScaledJacobian(vtkCell *cell)
vtkMeshQuality::QuadMinAngle
static double QuadMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a quadrilateral,...
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:65
vtkMeshQuality::QuadShapeAndSize
static double QuadShapeAndSize(vtkCell *cell)
vtkMeshQuality::QuadEdgeRatio
static double QuadEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a quadrilateral.
vtkMeshQuality::SetHexQualityMeasureToMedAspectFrobenius
void SetHexQualityMeasureToMedAspectFrobenius()
Definition: vtkMeshQuality.h:311
vtkMeshQuality::HexShear
static double HexShear(vtkCell *cell)
vtkMeshQuality::TriangleDistortion
static double TriangleDistortion(vtkCell *cell)
This is a static function used to calculate the distortion of a triangle.
vtkMeshQuality::HexMaxEdgeRatio
static double HexMaxEdgeRatio(vtkCell *cell)
vtkMeshQuality::HexSkew
static double HexSkew(vtkCell *cell)
vtkMeshQuality::TriangleMinAngle
static double TriangleMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a triangle,...
vtkMeshQuality::TriangleAspectFrobenius
static double TriangleAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the Frobenius condition number of the transformation matr...
vtkDataSetAlgorithm.h
vtkMeshQuality::SetTetQualityMeasureToDistortion
void SetTetQualityMeasureToDistortion()
Definition: vtkMeshQuality.h:293
vtkMeshQuality::SetQuadQualityMeasureToArea
void SetQuadQualityMeasureToArea()
Definition: vtkMeshQuality.h:222
vtkMeshQuality::GetCurrentTriangleNormal
static int GetCurrentTriangleNormal(double point[3], double normal[3])
A function called by some VERDICT triangle quality functions to test for inverted triangles.
VTK_QUALITY_RELATIVE_SIZE_SQUARED
#define VTK_QUALITY_RELATIVE_SIZE_SQUARED
Definition: vtkMeshQuality.h:83
vtkMeshQuality::TetAspectBeta
static double TetAspectBeta(vtkCell *cell)
vtkMeshQuality::TetRadiusRatio
static double TetRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a tetrahedron.
vtkMeshQuality::CellNormals
vtkDataArray * CellNormals
Definition: vtkMeshQuality.h:824
vtkMeshQuality::TetDistortion
static double TetDistortion(vtkCell *cell)
vtkMeshQuality::HexDistortion
static double HexDistortion(vtkCell *cell)
vtkMeshQuality::TriangleShape
static double TriangleShape(vtkCell *cell)
This is a static function used to calculate the shape of a triangle.
vtkMeshQuality::QuadSkew
static double QuadSkew(vtkCell *cell)
vtkMeshQuality::TetVolume
static double TetVolume(vtkCell *cell)
VTK_QUALITY_SHAPE_AND_SIZE
#define VTK_QUALITY_SHAPE_AND_SIZE
Definition: vtkMeshQuality.h:85
vtkMeshQuality::SetQuadQualityMeasureToEdgeRatio
void SetQuadQualityMeasureToEdgeRatio()
Definition: vtkMeshQuality.h:198
vtkMeshQuality::~vtkMeshQuality
~vtkMeshQuality() override
vtkMeshQuality::SetTetQualityMeasureToAspectRatio
void SetTetQualityMeasureToAspectRatio()
Definition: vtkMeshQuality.h:264
vtkMeshQuality::TetRelativeSizeSquared
static double TetRelativeSizeSquared(vtkCell *cell)
vtkMeshQuality::QuadMaxEdgeRatios
static double QuadMaxEdgeRatios(vtkCell *cell)
VTK_QUALITY_TAPER
#define VTK_QUALITY_TAPER
Definition: vtkMeshQuality.h:89
vtkMeshQuality::TriangleAspectRatio
static double TriangleAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a triangle.
VTK_QUALITY_ASPECT_BETA
#define VTK_QUALITY_ASPECT_BETA
Definition: vtkMeshQuality.h:100
vtkMeshQuality::HexJacobian
static double HexJacobian(vtkCell *cell)
vtkMeshQuality::TetShapeandSize
static double TetShapeandSize(vtkCell *cell)
vtkMeshQuality::SetQuadQualityMeasureToSkew
void SetQuadQualityMeasureToSkew()
Definition: vtkMeshQuality.h:219
vtkMeshQuality::SetQuadQualityMeasureToDistortion
void SetQuadQualityMeasureToDistortion()
Definition: vtkMeshQuality.h:247
vtkMeshQuality::HexRelativeSizeSquared
static double HexRelativeSizeSquared(vtkCell *cell)
vtkMeshQuality::QuadRadiusRatio
static double QuadRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a planar quadrilateral.
vtkMeshQuality::TetShape
static double TetShape(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToShape
void SetTetQualityMeasureToShape()
Definition: vtkMeshQuality.h:284
vtkMeshQuality::SetHexQualityMeasureToRelativeSizeSquared
void SetHexQualityMeasureToRelativeSizeSquared()
Definition: vtkMeshQuality.h:338
VTK_QUALITY_DIMENSION
#define VTK_QUALITY_DIMENSION
Definition: vtkMeshQuality.h:93
VTK_QUALITY_EDGE_RATIO
#define VTK_QUALITY_EDGE_RATIO
Definition: vtkMeshQuality.h:71
vtkMeshQuality::QuadScaledJacobian
static double QuadScaledJacobian(vtkCell *cell)
vtkMeshQuality::SetTetQualityMeasureToAspectGamma
void SetTetQualityMeasureToAspectGamma()
Definition: vtkMeshQuality.h:276
vtkMeshQuality::HexDiagonal
static double HexDiagonal(vtkCell *cell)
vtkMeshQuality::HexOddy
static double HexOddy(vtkCell *cell)
vtkMeshQuality::SetHexQualityMeasureToShape
void SetHexQualityMeasureToShape()
Definition: vtkMeshQuality.h:337
vtkMeshQuality::SetQuadQualityMeasureToShape
void SetQuadQualityMeasureToShape()
Definition: vtkMeshQuality.h:234
vtkMeshQuality::QuadDistortion
static double QuadDistortion(vtkCell *cell)
vtkMeshQuality::SetTriangleQualityMeasureToMinAngle
void SetTriangleQualityMeasureToMinAngle()
Definition: vtkMeshQuality.h:149
vtkMeshQuality::HexVolume
static double HexVolume(vtkCell *cell)
vtkMeshQuality::vtkMeshQuality
vtkMeshQuality()
vtkMeshQuality::TriangleShapeAndSize
static double TriangleShapeAndSize(vtkCell *cell)
This is a static function used to calculate the product of shape and relative size of a triangle.
vtkMeshQuality::RequestData
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
VTK_QUALITY_MAX_ASPECT_FROBENIUS
#define VTK_QUALITY_MAX_ASPECT_FROBENIUS
Definition: vtkMeshQuality.h:76
vtkMeshQuality::New
static vtkMeshQuality * New()
VTK_QUALITY_COLLAPSE_RATIO
#define VTK_QUALITY_COLLAPSE_RATIO
Definition: vtkMeshQuality.h:78
VTK_QUALITY_MAX_EDGE_RATIO
#define VTK_QUALITY_MAX_EDGE_RATIO
Definition: vtkMeshQuality.h:87
vtkMeshQuality::SetQuadQualityMeasureToCondition
void SetQuadQualityMeasureToCondition()
Definition: vtkMeshQuality.h:227
VTK_QUALITY_ASPECT_GAMMA
#define VTK_QUALITY_ASPECT_GAMMA
Definition: vtkMeshQuality.h:98
vtkTypeBool
int vtkTypeBool
Definition: vtkABI.h:69
vtkMeshQuality::HexCondition
static double HexCondition(vtkCell *cell)
vtkMeshQuality::SetHexQualityMeasureToDimension
void SetHexQualityMeasureToDimension()
Definition: vtkMeshQuality.h:328
vtkMeshQuality::TetAspectRatio
static double TetAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a tetrahedron.
vtkMeshQuality::SetQuadQualityMeasureToOddy
void SetQuadQualityMeasureToOddy()
Definition: vtkMeshQuality.h:226
vtkMeshQuality::SetTriangleQualityMeasureToAspectFrobenius
void SetTriangleQualityMeasureToAspectFrobenius()
Definition: vtkMeshQuality.h:145