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=========================================================================*/
91#ifndef vtkMeshQuality_h
92#define vtkMeshQuality_h
93
94#include "vtkDataSetAlgorithm.h"
95#include "vtkFiltersVerdictModule.h" // For export macro
96
97class vtkCell;
98class vtkDataArray;
99
100#define VTK_QUALITY_EDGE_RATIO 0
101#define VTK_QUALITY_ASPECT_RATIO 1
102#define VTK_QUALITY_RADIUS_RATIO 2
103#define VTK_QUALITY_ASPECT_FROBENIUS 3
104#define VTK_QUALITY_MED_ASPECT_FROBENIUS 4
105#define VTK_QUALITY_MAX_ASPECT_FROBENIUS 5
106#define VTK_QUALITY_MIN_ANGLE 6
107#define VTK_QUALITY_COLLAPSE_RATIO 7
108#define VTK_QUALITY_MAX_ANGLE 8
109#define VTK_QUALITY_CONDITION 9
110#define VTK_QUALITY_SCALED_JACOBIAN 10
111#define VTK_QUALITY_SHEAR 11
112#define VTK_QUALITY_RELATIVE_SIZE_SQUARED 12
113#define VTK_QUALITY_SHAPE 13
114#define VTK_QUALITY_SHAPE_AND_SIZE 14
115#define VTK_QUALITY_DISTORTION 15
116#define VTK_QUALITY_MAX_EDGE_RATIO 16
117#define VTK_QUALITY_SKEW 17
118#define VTK_QUALITY_TAPER 18
119#define VTK_QUALITY_VOLUME 19
120#define VTK_QUALITY_STRETCH 20
121#define VTK_QUALITY_DIAGONAL 21
122#define VTK_QUALITY_DIMENSION 22
123#define VTK_QUALITY_ODDY 23
124#define VTK_QUALITY_SHEAR_AND_SIZE 24
125#define VTK_QUALITY_JACOBIAN 25
126#define VTK_QUALITY_WARPAGE 26
127#define VTK_QUALITY_ASPECT_GAMMA 27
128#define VTK_QUALITY_AREA 28
129#define VTK_QUALITY_ASPECT_BETA 29
130
131class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
132{
133public:
134 void PrintSelf(ostream& os, vtkIndent indent) override;
137
139
145 vtkSetMacro(SaveCellQuality, vtkTypeBool);
146 vtkGetMacro(SaveCellQuality, vtkTypeBool);
147 vtkBooleanMacro(SaveCellQuality, vtkTypeBool);
149
151
159 vtkSetMacro(TriangleQualityMeasure, int);
160 vtkGetMacro(TriangleQualityMeasure, int);
161 void SetTriangleQualityMeasureToArea() { this->SetTriangleQualityMeasure(VTK_QUALITY_AREA); }
163 {
164 this->SetTriangleQualityMeasure(VTK_QUALITY_EDGE_RATIO);
165 }
167 {
168 this->SetTriangleQualityMeasure(VTK_QUALITY_ASPECT_RATIO);
169 }
171 {
172 this->SetTriangleQualityMeasure(VTK_QUALITY_RADIUS_RATIO);
173 }
175 {
176 this->SetTriangleQualityMeasure(VTK_QUALITY_ASPECT_FROBENIUS);
177 }
179 {
180 this->SetTriangleQualityMeasure(VTK_QUALITY_MIN_ANGLE);
181 }
183 {
184 this->SetTriangleQualityMeasure(VTK_QUALITY_MAX_ANGLE);
185 }
187 {
188 this->SetTriangleQualityMeasure(VTK_QUALITY_CONDITION);
189 }
191 {
192 this->SetTriangleQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
193 }
195 {
196 this->SetTriangleQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
197 }
198 void SetTriangleQualityMeasureToShape() { this->SetTriangleQualityMeasure(VTK_QUALITY_SHAPE); }
200 {
201 this->SetTriangleQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
202 }
204 {
205 this->SetTriangleQualityMeasure(VTK_QUALITY_DISTORTION);
206 }
208
210
225 vtkSetMacro(QuadQualityMeasure, int);
226 vtkGetMacro(QuadQualityMeasure, int);
227 void SetQuadQualityMeasureToEdgeRatio() { this->SetQuadQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
229 {
230 this->SetQuadQualityMeasure(VTK_QUALITY_ASPECT_RATIO);
231 }
233 {
234 this->SetQuadQualityMeasure(VTK_QUALITY_RADIUS_RATIO);
235 }
237 {
238 this->SetQuadQualityMeasure(VTK_QUALITY_MED_ASPECT_FROBENIUS);
239 }
241 {
242 this->SetQuadQualityMeasure(VTK_QUALITY_MAX_ASPECT_FROBENIUS);
243 }
245 {
246 this->SetQuadQualityMeasure(VTK_QUALITY_MAX_EDGE_RATIO);
247 }
248 void SetQuadQualityMeasureToSkew() { this->SetQuadQualityMeasure(VTK_QUALITY_SKEW); }
249 void SetQuadQualityMeasureToTaper() { this->SetQuadQualityMeasure(VTK_QUALITY_TAPER); }
250 void SetQuadQualityMeasureToWarpage() { this->SetQuadQualityMeasure(VTK_QUALITY_WARPAGE); }
251 void SetQuadQualityMeasureToArea() { this->SetQuadQualityMeasure(VTK_QUALITY_AREA); }
252 void SetQuadQualityMeasureToStretch() { this->SetQuadQualityMeasure(VTK_QUALITY_STRETCH); }
253 void SetQuadQualityMeasureToMinAngle() { this->SetQuadQualityMeasure(VTK_QUALITY_MIN_ANGLE); }
254 void SetQuadQualityMeasureToMaxAngle() { this->SetQuadQualityMeasure(VTK_QUALITY_MAX_ANGLE); }
255 void SetQuadQualityMeasureToOddy() { this->SetQuadQualityMeasure(VTK_QUALITY_ODDY); }
256 void SetQuadQualityMeasureToCondition() { this->SetQuadQualityMeasure(VTK_QUALITY_CONDITION); }
257 void SetQuadQualityMeasureToJacobian() { this->SetQuadQualityMeasure(VTK_QUALITY_JACOBIAN); }
259 {
260 this->SetQuadQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
261 }
262 void SetQuadQualityMeasureToShear() { this->SetQuadQualityMeasure(VTK_QUALITY_SHEAR); }
263 void SetQuadQualityMeasureToShape() { this->SetQuadQualityMeasure(VTK_QUALITY_SHAPE); }
265 {
266 this->SetQuadQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
267 }
269 {
270 this->SetQuadQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
271 }
273 {
274 this->SetQuadQualityMeasure(VTK_QUALITY_SHEAR_AND_SIZE);
275 }
276 void SetQuadQualityMeasureToDistortion() { this->SetQuadQualityMeasure(VTK_QUALITY_DISTORTION); }
278
280
290 vtkSetMacro(TetQualityMeasure, int);
291 vtkGetMacro(TetQualityMeasure, int);
292 void SetTetQualityMeasureToEdgeRatio() { this->SetTetQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
293 void SetTetQualityMeasureToAspectRatio() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_RATIO); }
294 void SetTetQualityMeasureToRadiusRatio() { this->SetTetQualityMeasure(VTK_QUALITY_RADIUS_RATIO); }
296 {
297 this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_FROBENIUS);
298 }
299 void SetTetQualityMeasureToMinAngle() { this->SetTetQualityMeasure(VTK_QUALITY_MIN_ANGLE); }
301 {
302 this->SetTetQualityMeasure(VTK_QUALITY_COLLAPSE_RATIO);
303 }
304 void SetTetQualityMeasureToAspectBeta() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_BETA); }
305 void SetTetQualityMeasureToAspectGamma() { this->SetTetQualityMeasure(VTK_QUALITY_ASPECT_GAMMA); }
306 void SetTetQualityMeasureToVolume() { this->SetTetQualityMeasure(VTK_QUALITY_VOLUME); }
307 void SetTetQualityMeasureToCondition() { this->SetTetQualityMeasure(VTK_QUALITY_CONDITION); }
308 void SetTetQualityMeasureToJacobian() { this->SetTetQualityMeasure(VTK_QUALITY_JACOBIAN); }
310 {
311 this->SetTetQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
312 }
313 void SetTetQualityMeasureToShape() { this->SetTetQualityMeasure(VTK_QUALITY_SHAPE); }
315 {
316 this->SetTetQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
317 }
319 {
320 this->SetTetQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
321 }
322 void SetTetQualityMeasureToDistortion() { this->SetTetQualityMeasure(VTK_QUALITY_DISTORTION); }
324
326
337 vtkSetMacro(HexQualityMeasure, int);
338 vtkGetMacro(HexQualityMeasure, int);
339 void SetHexQualityMeasureToEdgeRatio() { this->SetHexQualityMeasure(VTK_QUALITY_EDGE_RATIO); }
341 {
342 this->SetHexQualityMeasure(VTK_QUALITY_MED_ASPECT_FROBENIUS);
343 }
345 {
346 this->SetHexQualityMeasure(VTK_QUALITY_MAX_ASPECT_FROBENIUS);
347 }
349 {
350 this->SetHexQualityMeasure(VTK_QUALITY_MAX_EDGE_RATIO);
351 }
352 void SetHexQualityMeasureToSkew() { this->SetHexQualityMeasure(VTK_QUALITY_SKEW); }
353 void SetHexQualityMeasureToTaper() { this->SetHexQualityMeasure(VTK_QUALITY_TAPER); }
354 void SetHexQualityMeasureToVolume() { this->SetHexQualityMeasure(VTK_QUALITY_VOLUME); }
355 void SetHexQualityMeasureToStretch() { this->SetHexQualityMeasure(VTK_QUALITY_STRETCH); }
356 void SetHexQualityMeasureToDiagonal() { this->SetHexQualityMeasure(VTK_QUALITY_DIAGONAL); }
357 void SetHexQualityMeasureToDimension() { this->SetHexQualityMeasure(VTK_QUALITY_DIMENSION); }
358 void SetHexQualityMeasureToOddy() { this->SetHexQualityMeasure(VTK_QUALITY_ODDY); }
359 void SetHexQualityMeasureToCondition() { this->SetHexQualityMeasure(VTK_QUALITY_CONDITION); }
360 void SetHexQualityMeasureToJacobian() { this->SetHexQualityMeasure(VTK_QUALITY_JACOBIAN); }
362 {
363 this->SetHexQualityMeasure(VTK_QUALITY_SCALED_JACOBIAN);
364 }
365 void SetHexQualityMeasureToShear() { this->SetHexQualityMeasure(VTK_QUALITY_SHEAR); }
366 void SetHexQualityMeasureToShape() { this->SetHexQualityMeasure(VTK_QUALITY_SHAPE); }
368 {
369 this->SetHexQualityMeasure(VTK_QUALITY_RELATIVE_SIZE_SQUARED);
370 }
372 {
373 this->SetHexQualityMeasure(VTK_QUALITY_SHAPE_AND_SIZE);
374 }
376 {
377 this->SetHexQualityMeasure(VTK_QUALITY_SHEAR_AND_SIZE);
378 }
379 void SetHexQualityMeasureToDistortion() { this->SetHexQualityMeasure(VTK_QUALITY_DISTORTION); }
381
388 static double TriangleArea(vtkCell* cell);
389
400 static double TriangleEdgeRatio(vtkCell* cell);
401
412 static double TriangleAspectRatio(vtkCell* cell);
413
424 static double TriangleRadiusRatio(vtkCell* cell);
425
438 static double TriangleAspectFrobenius(vtkCell* cell);
439
447 static double TriangleMinAngle(vtkCell* cell);
448
456 static double TriangleMaxAngle(vtkCell* cell);
457
465 static double TriangleCondition(vtkCell* cell);
466
473 static double TriangleScaledJacobian(vtkCell* cell);
474
482
489 static double TriangleShape(vtkCell* cell);
490
496 static double TriangleShapeAndSize(vtkCell* cell);
497
504 static double TriangleDistortion(vtkCell* cell);
505
516 static double QuadEdgeRatio(vtkCell* cell);
517
529 static double QuadAspectRatio(vtkCell* cell);
530
545 static double QuadRadiusRatio(vtkCell* cell);
546
560 static double QuadMedAspectFrobenius(vtkCell* cell);
561
575 static double QuadMaxAspectFrobenius(vtkCell* cell);
576
584 static double QuadMinAngle(vtkCell* cell);
585
586 static double QuadMaxEdgeRatios(vtkCell* cell);
587 static double QuadSkew(vtkCell* cell);
588 static double QuadTaper(vtkCell* cell);
589 static double QuadWarpage(vtkCell* cell);
590 static double QuadArea(vtkCell* cell);
591 static double QuadStretch(vtkCell* cell);
592 static double QuadMaxAngle(vtkCell* cell);
593 static double QuadOddy(vtkCell* cell);
594 static double QuadCondition(vtkCell* cell);
595 static double QuadJacobian(vtkCell* cell);
596 static double QuadScaledJacobian(vtkCell* cell);
597 static double QuadShear(vtkCell* cell);
598 static double QuadShape(vtkCell* cell);
599 static double QuadRelativeSizeSquared(vtkCell* cell);
600 static double QuadShapeAndSize(vtkCell* cell);
601 static double QuadShearAndSize(vtkCell* cell);
602 static double QuadDistortion(vtkCell* cell);
603
614 static double TetEdgeRatio(vtkCell* cell);
615
626 static double TetAspectRatio(vtkCell* cell);
627
638 static double TetRadiusRatio(vtkCell* cell);
639
653 static double TetAspectFrobenius(vtkCell* cell);
654
662 static double TetMinAngle(vtkCell* cell);
663
665
674 static double TetCollapseRatio(vtkCell* cell);
675 static double TetAspectBeta(vtkCell* cell);
676 static double TetAspectGamma(vtkCell* cell);
677 static double TetVolume(vtkCell* cell);
678 static double TetCondition(vtkCell* cell);
679 static double TetJacobian(vtkCell* cell);
680 static double TetScaledJacobian(vtkCell* cell);
681 static double TetShape(vtkCell* cell);
682 static double TetRelativeSizeSquared(vtkCell* cell);
683 static double TetShapeandSize(vtkCell* cell);
684 static double TetDistortion(vtkCell* cell);
686
697 static double HexEdgeRatio(vtkCell* cell);
698
707 static double HexMedAspectFrobenius(vtkCell* cell);
708
710
718 static double HexMaxAspectFrobenius(vtkCell* cell);
719 static double HexMaxEdgeRatio(vtkCell* cell);
720 static double HexSkew(vtkCell* cell);
721 static double HexTaper(vtkCell* cell);
722 static double HexVolume(vtkCell* cell);
723 static double HexStretch(vtkCell* cell);
724 static double HexDiagonal(vtkCell* cell);
725 static double HexDimension(vtkCell* cell);
726 static double HexOddy(vtkCell* cell);
727 static double HexCondition(vtkCell* cell);
728 static double HexJacobian(vtkCell* cell);
729 static double HexScaledJacobian(vtkCell* cell);
730 static double HexShear(vtkCell* cell);
731 static double HexShape(vtkCell* cell);
732 static double HexRelativeSizeSquared(vtkCell* cell);
733 static double HexShapeAndSize(vtkCell* cell);
734 static double HexShearAndSize(vtkCell* cell);
735 static double HexDistortion(vtkCell* cell);
737
748 virtual void SetRatio(vtkTypeBool r) { this->SetSaveCellQuality(r); }
749 vtkTypeBool GetRatio() { return this->GetSaveCellQuality(); }
750 vtkBooleanMacro(Ratio, vtkTypeBool);
751
753
770 virtual void SetVolume(vtkTypeBool cv)
771 {
772 if (!((cv != 0) ^ (this->Volume != 0)))
773 {
774 return;
775 }
776 this->Modified();
777 this->Volume = cv;
778 if (this->Volume)
779 {
780 this->CompatibilityModeOn();
781 }
782 }
783 vtkTypeBool GetVolume() { return this->Volume; }
784 vtkBooleanMacro(Volume, vtkTypeBool);
786
788
816 {
817 if (!((cm != 0) ^ (this->CompatibilityMode != 0)))
818 {
819 return;
820 }
821 this->CompatibilityMode = cm;
822 this->Modified();
823 if (this->CompatibilityMode)
824 {
825 this->Volume = 1;
826 this->TetQualityMeasure = VTK_QUALITY_RADIUS_RATIO;
827 }
828 }
829 vtkGetMacro(CompatibilityMode, vtkTypeBool);
830 vtkBooleanMacro(CompatibilityMode, vtkTypeBool);
832
833protected:
835 ~vtkMeshQuality() override;
836
838
842 static int GetCurrentTriangleNormal(double point[3], double normal[3]);
843
849
852
854 static double CurrentTriNormal[3];
855
856private:
857 vtkMeshQuality(const vtkMeshQuality&) = delete;
858 void operator=(const vtkMeshQuality&) = delete;
859};
860
861#endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:147
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
Superclass for algorithms that produce output of the same type as input.
a simple class to control print indentation
Definition: vtkIndent.h:113
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