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=========================================================================*/
62#ifndef vtkMeshQuality_h
63#define vtkMeshQuality_h
64
65#include "vtkDataSetAlgorithm.h"
66#include "vtkFiltersVerdictModule.h" // For export macro
67
68class vtkCell;
69class 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
102class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
103{
104public:
105 void PrintSelf(ostream& os, vtkIndent indent) override;
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
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
804protected:
806 ~vtkMeshQuality() override;
807
809
813 static int GetCurrentTriangleNormal(double point[3], double normal[3]);
814
820
823
825 static double CurrentTriNormal[3];
826
827private:
828 vtkMeshQuality(const vtkMeshQuality&) = delete;
829 void operator=(const vtkMeshQuality&) = delete;
830};
831
832#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:58
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
Superclass for algorithms that produce output of the same type as input.
a simple class to control print indentation
Definition: vtkIndent.h:34
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