VTK  9.1.0
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=========================================================================*/
71#ifndef vtkMeshQuality_h
72#define vtkMeshQuality_h
73
74#include "vtkDataSetAlgorithm.h"
75#include "vtkFiltersVerdictModule.h" // For export macro
76
77class vtkCell;
78class vtkDataArray;
79
80#define VTK_QUALITY_EDGE_RATIO 0
81#define VTK_QUALITY_ASPECT_RATIO 1
82#define VTK_QUALITY_RADIUS_RATIO 2
83#define VTK_QUALITY_ASPECT_FROBENIUS 3
84#define VTK_QUALITY_MED_ASPECT_FROBENIUS 4
85#define VTK_QUALITY_MAX_ASPECT_FROBENIUS 5
86#define VTK_QUALITY_MIN_ANGLE 6
87#define VTK_QUALITY_COLLAPSE_RATIO 7
88#define VTK_QUALITY_MAX_ANGLE 8
89#define VTK_QUALITY_CONDITION 9
90#define VTK_QUALITY_SCALED_JACOBIAN 10
91#define VTK_QUALITY_SHEAR 11
92#define VTK_QUALITY_RELATIVE_SIZE_SQUARED 12
93#define VTK_QUALITY_SHAPE 13
94#define VTK_QUALITY_SHAPE_AND_SIZE 14
95#define VTK_QUALITY_DISTORTION 15
96#define VTK_QUALITY_MAX_EDGE_RATIO 16
97#define VTK_QUALITY_SKEW 17
98#define VTK_QUALITY_TAPER 18
99#define VTK_QUALITY_VOLUME 19
100#define VTK_QUALITY_STRETCH 20
101#define VTK_QUALITY_DIAGONAL 21
102#define VTK_QUALITY_DIMENSION 22
103#define VTK_QUALITY_ODDY 23
104#define VTK_QUALITY_SHEAR_AND_SIZE 24
105#define VTK_QUALITY_JACOBIAN 25
106#define VTK_QUALITY_WARPAGE 26
107#define VTK_QUALITY_ASPECT_GAMMA 27
108#define VTK_QUALITY_AREA 28
109#define VTK_QUALITY_ASPECT_BETA 29
110
111class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
112{
113public:
114 void PrintSelf(ostream& os, vtkIndent indent) override;
117
119
125 vtkSetMacro(SaveCellQuality, vtkTypeBool);
126 vtkGetMacro(SaveCellQuality, vtkTypeBool);
127 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
129
131
139 vtkSetMacro(TriangleQualityMeasure, int);
140 vtkGetMacro(TriangleQualityMeasure, int);
141 void SetTriangleQualityMeasureToArea() { this->SetTriangleQualityMeasure(VTK_QUALITY_AREA); }
143 {
144 this->SetTriangleQualityMeasure(VTK_QUALITY_EDGE_RATIO);
145 }
147 {
148 this->SetTriangleQualityMeasure(VTK_QUALITY_ASPECT_RATIO);
149 }
151 {
152 this->SetTriangleQualityMeasure(VTK_QUALITY_RADIUS_RATIO);
153 }
155 {
156 this->SetTriangleQualityMeasure(VTK_QUALITY_ASPECT_FROBENIUS);
157 }
159 {
160 this->SetTriangleQualityMeasure(VTK_QUALITY_MIN_ANGLE);
161 }
163 {
164 this->SetTriangleQualityMeasure(VTK_QUALITY_MAX_ANGLE);
165 }
167 {
168 this->SetTriangleQualityMeasure(VTK_QUALITY_CONDITION);
169 }
171 {
172 this->SetTriangleQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
173 }
175 {
176 this->SetTriangleQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
177 }
178 void SetTriangleQualityMeasureToShape() { this->SetTriangleQualityMeasure(VTK_QUALITY_SHAPE); }
180 {
181 this->SetTriangleQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
182 }
184 {
185 this->SetTriangleQualityMeasure(VTK_QUALITY_DISTORTION);
186 }
188
190
205 vtkSetMacro(QuadQualityMeasure, int);
206 vtkGetMacro(QuadQualityMeasure, int);
207 void SetQuadQualityMeasureToEdgeRatio() { this->SetQuadQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
209 {
210 this->SetQuadQualityMeasure(VTK_QUALITY_ASPECT_RATIO);
211 }
213 {
214 this->SetQuadQualityMeasure(VTK_QUALITY_RADIUS_RATIO);
215 }
217 {
218 this->SetQuadQualityMeasure(VTK_QUALITY_MED_ASPECT_FROBENIUS);
219 }
221 {
222 this->SetQuadQualityMeasure(VTK_QUALITY_MAX_ASPECT_FROBENIUS);
223 }
225 {
226 this->SetQuadQualityMeasure(VTK_QUALITY_MAX_EDGE_RATIO);
227 }
228 void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(VTK_QUALITY_SKEW); }
229 void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(VTK_QUALITY_TAPER); }
230 void SetQuadQualityMeasureToWarpage() { this->SetQuadQualityMeasure(VTK_QUALITY_WARPAGE); }
231 void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(VTK_QUALITY_AREA); }
232 void SetQuadQualityMeasureToStretch() { this->SetQuadQualityMeasure(VTK_QUALITY_STRETCH); }
233 void SetQuadQualityMeasureToMinAngle() { this->SetQuadQualityMeasure(VTK_QUALITY_MIN_ANGLE); }
234 void SetQuadQualityMeasureToMaxAngle() { this->SetQuadQualityMeasure(VTK_QUALITY_MAX_ANGLE); }
235 void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(VTK_QUALITY_ODDY); }
236 void SetQuadQualityMeasureToCondition() { this->SetQuadQualityMeasure(VTK_QUALITY_CONDITION); }
237 void SetQuadQualityMeasureToJacobian() { this->SetQuadQualityMeasure(VTK_QUALITY_JACOBIAN); }
239 {
240 this->SetQuadQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
241 }
242 void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(VTK_QUALITY_SHEAR); }
243 void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(VTK_QUALITY_SHAPE); }
245 {
246 this->SetQuadQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
247 }
249 {
250 this->SetQuadQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
251 }
253 {
254 this->SetQuadQualityMeasure(VTK_QUALITY_SHEAR_AND_SIZE);
255 }
256 void SetQuadQualityMeasureToDistortion() { this->SetQuadQualityMeasure(VTK_QUALITY_DISTORTION); }
258
260
270 vtkSetMacro(TetQualityMeasure, int);
271 vtkGetMacro(TetQualityMeasure, int);
272 void SetTetQualityMeasureToEdgeRatio() { this->SetTetQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
273 void SetTetQualityMeasureToAspectRatio() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_RATIO); }
274 void SetTetQualityMeasureToRadiusRatio() { this->SetTetQualityMeasure(VTK_QUALITY_RADIUS_RATIO); }
276 {
277 this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_FROBENIUS);
278 }
279 void SetTetQualityMeasureToMinAngle() { this->SetTetQualityMeasure(VTK_QUALITY_MIN_ANGLE); }
281 {
282 this->SetTetQualityMeasure(VTK_QUALITY_COLLAPSE_RATIO);
283 }
284 void SetTetQualityMeasureToAspectBeta() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_BETA); }
285 void SetTetQualityMeasureToAspectGamma() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_GAMMA); }
286 void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(VTK_QUALITY_VOLUME); }
287 void SetTetQualityMeasureToCondition() { this->SetTetQualityMeasure(VTK_QUALITY_CONDITION); }
288 void SetTetQualityMeasureToJacobian() { this->SetTetQualityMeasure(VTK_QUALITY_JACOBIAN); }
290 {
291 this->SetTetQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
292 }
293 void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(VTK_QUALITY_SHAPE); }
295 {
296 this->SetTetQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
297 }
299 {
300 this->SetTetQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
301 }
302 void SetTetQualityMeasureToDistortion() { this->SetTetQualityMeasure(VTK_QUALITY_DISTORTION); }
304
306
317 vtkSetMacro(HexQualityMeasure, int);
318 vtkGetMacro(HexQualityMeasure, int);
319 void SetHexQualityMeasureToEdgeRatio() { this->SetHexQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
321 {
322 this->SetHexQualityMeasure(VTK_QUALITY_MED_ASPECT_FROBENIUS);
323 }
325 {
326 this->SetHexQualityMeasure(VTK_QUALITY_MAX_ASPECT_FROBENIUS);
327 }
329 {
330 this->SetHexQualityMeasure(VTK_QUALITY_MAX_EDGE_RATIO);
331 }
332 void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(VTK_QUALITY_SKEW); }
333 void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(VTK_QUALITY_TAPER); }
334 void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(VTK_QUALITY_VOLUME); }
335 void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(VTK_QUALITY_STRETCH); }
336 void SetHexQualityMeasureToDiagonal() { this->SetHexQualityMeasure(VTK_QUALITY_DIAGONAL); }
337 void SetHexQualityMeasureToDimension() { this->SetHexQualityMeasure(VTK_QUALITY_DIMENSION); }
338 void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(VTK_QUALITY_ODDY); }
339 void SetHexQualityMeasureToCondition() { this->SetHexQualityMeasure(VTK_QUALITY_CONDITION); }
340 void SetHexQualityMeasureToJacobian() { this->SetHexQualityMeasure(VTK_QUALITY_JACOBIAN); }
342 {
343 this->SetHexQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
344 }
345 void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(VTK_QUALITY_SHEAR); }
346 void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(VTK_QUALITY_SHAPE); }
348 {
349 this->SetHexQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
350 }
352 {
353 this->SetHexQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
354 }
356 {
357 this->SetHexQualityMeasure(VTK_QUALITY_SHEAR_AND_SIZE);
358 }
359 void SetHexQualityMeasureToDistortion() { this->SetHexQualityMeasure(VTK_QUALITY_DISTORTION); }
361
368 static double TriangleArea(vtkCell* cell);
369
380 static double TriangleEdgeRatio(vtkCell* cell);
381
392 static double TriangleAspectRatio(vtkCell* cell);
393
404 static double TriangleRadiusRatio(vtkCell* cell);
405
418 static double TriangleAspectFrobenius(vtkCell* cell);
419
427 static double TriangleMinAngle(vtkCell* cell);
428
436 static double TriangleMaxAngle(vtkCell* cell);
437
445 static double TriangleCondition(vtkCell* cell);
446
453 static double TriangleScaledJacobian(vtkCell* cell);
454
462
469 static double TriangleShape(vtkCell* cell);
470
476 static double TriangleShapeAndSize(vtkCell* cell);
477
484 static double TriangleDistortion(vtkCell* cell);
485
496 static double QuadEdgeRatio(vtkCell* cell);
497
509 static double QuadAspectRatio(vtkCell* cell);
510
525 static double QuadRadiusRatio(vtkCell* cell);
526
540 static double QuadMedAspectFrobenius(vtkCell* cell);
541
555 static double QuadMaxAspectFrobenius(vtkCell* cell);
556
564 static double QuadMinAngle(vtkCell* cell);
565
566 static double QuadMaxEdgeRatios(vtkCell* cell);
567 static double QuadSkew(vtkCell* cell);
568 static double QuadTaper(vtkCell* cell);
569 static double QuadWarpage(vtkCell* cell);
570 static double QuadArea(vtkCell* cell);
571 static double QuadStretch(vtkCell* cell);
572 static double QuadMaxAngle(vtkCell* cell);
573 static double QuadOddy(vtkCell* cell);
574 static double QuadCondition(vtkCell* cell);
575 static double QuadJacobian(vtkCell* cell);
576 static double QuadScaledJacobian(vtkCell* cell);
577 static double QuadShear(vtkCell* cell);
578 static double QuadShape(vtkCell* cell);
579 static double QuadRelativeSizeSquared(vtkCell* cell);
580 static double QuadShapeAndSize(vtkCell* cell);
581 static double QuadShearAndSize(vtkCell* cell);
582 static double QuadDistortion(vtkCell* cell);
583
594 static double TetEdgeRatio(vtkCell* cell);
595
606 static double TetAspectRatio(vtkCell* cell);
607
618 static double TetRadiusRatio(vtkCell* cell);
619
633 static double TetAspectFrobenius(vtkCell* cell);
634
642 static double TetMinAngle(vtkCell* cell);
643
645
654 static double TetCollapseRatio(vtkCell* cell);
655 static double TetAspectBeta(vtkCell* cell);
656 static double TetAspectGamma(vtkCell* cell);
657 static double TetVolume(vtkCell* cell);
658 static double TetCondition(vtkCell* cell);
659 static double TetJacobian(vtkCell* cell);
660 static double TetScaledJacobian(vtkCell* cell);
661 static double TetShape(vtkCell* cell);
662 static double TetRelativeSizeSquared(vtkCell* cell);
663 static double TetShapeandSize(vtkCell* cell);
664 static double TetDistortion(vtkCell* cell);
666
677 static double HexEdgeRatio(vtkCell* cell);
678
687 static double HexMedAspectFrobenius(vtkCell* cell);
688
690
698 static double HexMaxAspectFrobenius(vtkCell* cell);
699 static double HexMaxEdgeRatio(vtkCell* cell);
700 static double HexSkew(vtkCell* cell);
701 static double HexTaper(vtkCell* cell);
702 static double HexVolume(vtkCell* cell);
703 static double HexStretch(vtkCell* cell);
704 static double HexDiagonal(vtkCell* cell);
705 static double HexDimension(vtkCell* cell);
706 static double HexOddy(vtkCell* cell);
707 static double HexCondition(vtkCell* cell);
708 static double HexJacobian(vtkCell* cell);
709 static double HexScaledJacobian(vtkCell* cell);
710 static double HexShear(vtkCell* cell);
711 static double HexShape(vtkCell* cell);
712 static double HexRelativeSizeSquared(vtkCell* cell);
713 static double HexShapeAndSize(vtkCell* cell);
714 static double HexShearAndSize(vtkCell* cell);
715 static double HexDistortion(vtkCell* cell);
717
728 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
729 vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
730 vtkBooleanMacro(Ratio, vtkTypeBool);
731
733
750 virtual void SetVolume(vtkTypeBool cv)
751 {
752 if (!((cv != 0) ^ (this->Volume != 0)))
753 {
754 return;
755 }
756 this->Modified();
757 this->Volume = cv;
758 if (this->Volume)
759 {
760 this->CompatibilityModeOn();
761 }
762 }
763 vtkTypeBool GetVolume() { return this->Volume; }
764 vtkBooleanMacro(Volume, vtkTypeBool);
766
768
796 {
797 if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
798 {
799 return;
800 }
801 this->CompatibilityMode = cm;
802 this->Modified();
803 if (this->CompatibilityMode)
804 {
805 this->Volume = 1;
806 this->TetQualityMeasure = VTK_QUALITY_RADIUS_RATIO;
807 }
808 }
809 vtkGetMacro(CompatibilityMode, vtkTypeBool);
810 vtkBooleanMacro(CompatibilityMode, vtkTypeBool);
812
813protected:
815 ~vtkMeshQuality() override;
816
818
822 static int GetCurrentTriangleNormal(double point[3], double normal[3]);
823
829
832
834 static double CurrentTriNormal[3];
835
836private:
837 vtkMeshQuality(const vtkMeshQuality&) = delete;
838 void operator=(const vtkMeshQuality&) = delete;
839};
840
841#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:67
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
Superclass for algorithms that produce output of the same type as input.
a simple class to control print indentation
Definition: vtkIndent.h:43
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.
static double QuadWarpage(vtkCell *cell)
void SetHexQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadTaper(vtkCell *cell)
void SetTriangleQualityMeasureToArea()
Set/Get the particular estimator used to function the quality of triangles.
static double HexOddy(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double TetAspectGamma(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
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)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double HexScaledJacobian(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
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.
static double QuadJacobian(vtkCell *cell)
void SetTriangleQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to function the quality of triangles.
void SetQuadQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a triangle.
void SetHexQualityMeasureToMaxAspectFrobenius()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadOddy(vtkCell *cell)
static double QuadAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a planar quadrilateral.
static double TriangleAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the Frobenius condition number of the transformation matr...
static double QuadScaledJacobian(vtkCell *cell)
static double HexMaxEdgeRatio(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double HexShapeAndSize(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double QuadShear(vtkCell *cell)
void SetHexQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of hexahedra.
void SetQuadQualityMeasureToAspectRatio()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleShapeAndSize(vtkCell *cell)
This is a static function used to calculate the product of shape and relative size of a triangle.
static int GetCurrentTriangleNormal(double point[3], double normal[3])
A function called by some VERDICT triangle quality functions to test for inverted triangles.
vtkTypeBool GetRatio()
static double HexTaper(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static vtkMeshQuality * New()
vtkDataArray * CellNormals
static double TetAspectBeta(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TetJacobian(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double HexDistortion(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double TetCollapseRatio(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TriangleDistortion(vtkCell *cell)
This is a static function used to calculate the distortion of a triangle.
static double HexVolume(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadArea(vtkCell *cell)
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.
static double TetEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a tetrahedron.
void SetQuadQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexRelativeSizeSquared(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetTetQualityMeasureToShape()
Set/Get the particular estimator used to measure the quality of tetrahedra.
void SetQuadQualityMeasureToShear()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexShape(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetHexQualityMeasureToCondition()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double TriangleScaledJacobian(vtkCell *cell)
This is a static function used to 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)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TriangleShape(vtkCell *cell)
This is a static function used to 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)
This is a static function used to calculate the average Frobenius aspect of the 4 corner triangles of...
static double TriangleCondition(vtkCell *cell)
This is a static function used to calculate the condition number of a triangle.
void SetQuadQualityMeasureToTaper()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TetMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) dihedral angle of a tetrahedron...
static double TetShapeandSize(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
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.
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)
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.
~vtkMeshQuality() override
void SetQuadQualityMeasureToJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TriangleAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a triangle.
static double HexJacobian(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetQuadQualityMeasureToMedAspectFrobenius()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double TetDistortion(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double QuadShapeAndSize(vtkCell *cell)
static double TriangleMaxAngle(vtkCell *cell)
This is a static function used to calculate the maximal (nonoriented) angle of a triangle,...
void SetTriangleQualityMeasureToCondition()
Set/Get the particular estimator used to function the quality of triangles.
static double TetCondition(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double TriangleRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a triangle.
static double TetScaledJacobian(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double QuadEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a quadrilateral.
void SetTetQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadShearAndSize(vtkCell *cell)
static double HexDimension(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
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)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double QuadMaxAngle(vtkCell *cell)
static double QuadStretch(vtkCell *cell)
void SetTriangleQualityMeasureToAspectFrobenius()
Set/Get the particular estimator used to function the quality of triangles.
virtual void SetCompatibilityMode(vtkTypeBool cm)
CompatibilityMode governs whether, when both a quality function and cell volume are to be stored as c...
static double HexEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a hexahedron.
void SetTetQualityMeasureToRelativeSizeSquared()
Set/Get the particular estimator used to measure the quality of tetrahedra.
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)
This is a static function used to calculate the maximal Frobenius aspect of the 4 corner triangles of...
vtkTypeBool GetVolume()
These methods are deprecated.
vtkTypeBool SaveCellQuality
void SetHexQualityMeasureToShapeAndSize()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadMaxEdgeRatios(vtkCell *cell)
static double TetAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the Frobenius condition number of the transformation matr...
static double HexMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
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 HexShearAndSize(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
static double HexSkew(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
void SetTetQualityMeasureToCollapseRatio()
Set/Get the particular estimator used to measure the quality of tetrahedra.
vtkTypeBool Volume
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)
void SetQuadQualityMeasureToMaxEdgeRatios()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
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)
This is a static function used to calculate the average Frobenius aspect of the 8 corner tetrahedra o...
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)
This is a static function used to calculate the square of the relative size of a triangle.
static double QuadRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a planar quadrilateral.
void SetHexQualityMeasureToMaxEdgeRatios()
Set/Get the particular estimator used to measure the quality of hexahedra.
static double QuadRelativeSizeSquared(vtkCell *cell)
static double TetVolume(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
static double QuadSkew(vtkCell *cell)
static double QuadDistortion(vtkCell *cell)
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.
void SetHexQualityMeasureToEdgeRatio()
Set/Get the particular estimator used to measure the quality of hexahedra.
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 TriangleMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a triangle,...
static double HexDiagonal(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
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)
This is a static function used to calculate the area of a triangle.
void SetTetQualityMeasureToAspectGamma()
Set/Get the particular estimator used to measure the quality of tetrahedra.
static double QuadMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a quadrilateral,...
void SetQuadQualityMeasureToScaledJacobian()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
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)
This is a static function used to calculate the radius ratio of a tetrahedron.
void SetTetQualityMeasureToAspectBeta()
Set/Get the particular estimator used to measure the quality of tetrahedra.
virtual void SetVolume(vtkTypeBool cv)
These methods are deprecated.
static double TetAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a tetrahedron.
static double TetShape(vtkCell *cell)
This is a static function used to calculate the collapse ratio of a tetrahedron.
void SetQuadQualityMeasureToDistortion()
Set/Get the particular estimator used to measure the quality of quadrilaterals.
static double HexCondition(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 8 corner tetrahedra o...
virtual void Modified()
Update the modification time for this object.
@ point
Definition: vtkX3D.h:242
int vtkTypeBool
Definition: vtkABI.h:69
#define VTK_QUALITY_AREA
#define VTK_QUALITY_SHAPE_AND_SIZE
#define VTK_QUALITY_STRETCH
#define VTK_QUALITY_MAX_ANGLE
#define VTK_QUALITY_MIN_ANGLE
#define VTK_QUALITY_MAX_ASPECT_FROBENIUS
#define VTK_QUALITY_SHEAR_AND_SIZE
#define VTK_QUALITY_RELATIVE_SIZE_SQUARED
#define VTK_QUALITY_JACOBIAN
#define VTK_QUALITY_SHEAR
#define VTK_QUALITY_ASPECT_GAMMA
#define VTK_QUALITY_VOLUME
#define VTK_QUALITY_EDGE_RATIO
#define VTK_QUALITY_DISTORTION
#define VTK_QUALITY_SHAPE
#define VTK_QUALITY_WARPAGE
#define VTK_QUALITY_RADIUS_RATIO
#define VTK_QUALITY_MED_ASPECT_FROBENIUS
#define VTK_QUALITY_CONDITION
#define VTK_QUALITY_DIAGONAL
#define VTK_QUALITY_ODDY
#define VTK_QUALITY_SCALED_JACOBIAN
#define VTK_QUALITY_COLLAPSE_RATIO
#define VTK_QUALITY_TAPER
#define VTK_QUALITY_DIMENSION
#define VTK_QUALITY_ASPECT_BETA
#define VTK_QUALITY_MAX_EDGE_RATIO
#define VTK_QUALITY_ASPECT_RATIO
#define VTK_QUALITY_SKEW
#define VTK_QUALITY_ASPECT_FROBENIUS