VTK
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 =========================================================================*/
65 #ifndef vtkMeshQuality_h
66 #define vtkMeshQuality_h
67 
68 #include "vtkFiltersVerdictModule.h" // For export macro
69 #include "vtkDataSetAlgorithm.h"
70 
71 class vtkCell;
72 class vtkDataArray;
73 
74 #define VTK_QUALITY_EDGE_RATIO 0
75 #define VTK_QUALITY_ASPECT_RATIO 1
76 #define VTK_QUALITY_RADIUS_RATIO 2
77 #define VTK_QUALITY_ASPECT_FROBENIUS 3
78 #define VTK_QUALITY_MED_ASPECT_FROBENIUS 4
79 #define VTK_QUALITY_MAX_ASPECT_FROBENIUS 5
80 #define VTK_QUALITY_MIN_ANGLE 6
81 #define VTK_QUALITY_COLLAPSE_RATIO 7
82 #define VTK_QUALITY_MAX_ANGLE 8
83 #define VTK_QUALITY_CONDITION 9
84 #define VTK_QUALITY_SCALED_JACOBIAN 10
85 #define VTK_QUALITY_SHEAR 11
86 #define VTK_QUALITY_RELATIVE_SIZE_SQUARED 12
87 #define VTK_QUALITY_SHAPE 13
88 #define VTK_QUALITY_SHAPE_AND_SIZE 14
89 #define VTK_QUALITY_DISTORTION 15
90 #define VTK_QUALITY_MAX_EDGE_RATIO 16
91 #define VTK_QUALITY_SKEW 17
92 #define VTK_QUALITY_TAPER 18
93 #define VTK_QUALITY_VOLUME 19
94 #define VTK_QUALITY_STRETCH 20
95 #define VTK_QUALITY_DIAGONAL 21
96 #define VTK_QUALITY_DIMENSION 22
97 #define VTK_QUALITY_ODDY 23
98 #define VTK_QUALITY_SHEAR_AND_SIZE 24
99 #define VTK_QUALITY_JACOBIAN 25
100 #define VTK_QUALITY_WARPAGE 26
101 #define VTK_QUALITY_ASPECT_GAMMA 27
102 #define VTK_QUALITY_AREA 28
103 #define VTK_QUALITY_ASPECT_BETA 29
104 
105 class VTKFILTERSVERDICT_EXPORT vtkMeshQuality : public vtkDataSetAlgorithm
106 {
107 public:
108  void PrintSelf(ostream& os, vtkIndent indent);
110  static vtkMeshQuality* New();
111 
113 
119  vtkSetMacro(SaveCellQuality,int);
120  vtkGetMacro(SaveCellQuality,int);
121  vtkBooleanMacro(SaveCellQuality,int);
123 
125 
133  vtkSetMacro(TriangleQualityMeasure,int);
134  vtkGetMacro(TriangleQualityMeasure,int);
136  {
137  this->SetTriangleQualityMeasure( VTK_QUALITY_AREA );
138  }
140  {
141  this->SetTriangleQualityMeasure( VTK_QUALITY_EDGE_RATIO );
142  }
144  {
145  this->SetTriangleQualityMeasure( VTK_QUALITY_ASPECT_RATIO );
146  }
148  {
149  this->SetTriangleQualityMeasure( VTK_QUALITY_RADIUS_RATIO );
150  }
152  {
153  this->SetTriangleQualityMeasure( VTK_QUALITY_ASPECT_FROBENIUS );
154  }
156  {
157  this->SetTriangleQualityMeasure( VTK_QUALITY_MIN_ANGLE );
158  }
160  {
161  this->SetTriangleQualityMeasure( VTK_QUALITY_MAX_ANGLE );
162  }
164  {
165  this->SetTriangleQualityMeasure( VTK_QUALITY_CONDITION );
166  }
168  {
169  this->SetTriangleQualityMeasure( VTK_QUALITY_SCALED_JACOBIAN );
170  }
172  {
173  this->SetTriangleQualityMeasure( VTK_QUALITY_RELATIVE_SIZE_SQUARED );
174  }
176  {
177  this->SetTriangleQualityMeasure( VTK_QUALITY_SHAPE );
178  }
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);
208  {
209  this->SetQuadQualityMeasure( VTK_QUALITY_EDGE_RATIO );
210  }
212  {
213  this->SetQuadQualityMeasure( VTK_QUALITY_ASPECT_RATIO );
214  }
216  {
217  this->SetQuadQualityMeasure( VTK_QUALITY_RADIUS_RATIO );
218  }
220  {
221  this->SetQuadQualityMeasure( VTK_QUALITY_MED_ASPECT_FROBENIUS );
222  }
224  {
225  this->SetQuadQualityMeasure( VTK_QUALITY_MAX_ASPECT_FROBENIUS );
226  }
228  {
229  this->SetQuadQualityMeasure( VTK_QUALITY_MAX_EDGE_RATIO );
230  }
232  {
233  this->SetQuadQualityMeasure( VTK_QUALITY_SKEW );
234  }
236  {
237  this->SetQuadQualityMeasure( VTK_QUALITY_TAPER );
238  }
240  {
241  this->SetQuadQualityMeasure( VTK_QUALITY_WARPAGE );
242  }
244  {
245  this->SetQuadQualityMeasure( VTK_QUALITY_AREA );
246  }
248  {
249  this->SetQuadQualityMeasure( VTK_QUALITY_STRETCH );
250  }
252  {
253  this->SetQuadQualityMeasure( VTK_QUALITY_MIN_ANGLE );
254  }
256  {
257  this->SetQuadQualityMeasure( VTK_QUALITY_MAX_ANGLE );
258  }
260  {
261  this->SetQuadQualityMeasure( VTK_QUALITY_ODDY );
262  }
264  {
265  this->SetQuadQualityMeasure( VTK_QUALITY_CONDITION );
266  }
268  {
269  this->SetQuadQualityMeasure( VTK_QUALITY_JACOBIAN );
270  }
272  {
273  this->SetQuadQualityMeasure( VTK_QUALITY_SCALED_JACOBIAN );
274  }
276  {
277  this->SetQuadQualityMeasure( VTK_QUALITY_SHEAR );
278  }
280  {
281  this->SetQuadQualityMeasure( VTK_QUALITY_SHAPE );
282  }
284  {
285  this->SetQuadQualityMeasure( VTK_QUALITY_RELATIVE_SIZE_SQUARED );
286  }
288  {
289  this->SetQuadQualityMeasure( VTK_QUALITY_SHAPE_AND_SIZE );
290  }
292  {
293  this->SetQuadQualityMeasure( VTK_QUALITY_SHEAR_AND_SIZE );
294  }
296  {
297  this->SetQuadQualityMeasure( VTK_QUALITY_DISTORTION );
298  }
300 
302 
312  vtkSetMacro(TetQualityMeasure,int);
313  vtkGetMacro(TetQualityMeasure,int);
315  {
316  this->SetTetQualityMeasure( VTK_QUALITY_EDGE_RATIO );
317  }
319  {
320  this->SetTetQualityMeasure( VTK_QUALITY_ASPECT_RATIO );
321  }
323  {
324  this->SetTetQualityMeasure( VTK_QUALITY_RADIUS_RATIO );
325  }
327  {
328  this->SetTetQualityMeasure( VTK_QUALITY_ASPECT_FROBENIUS );
329  }
331  {
332  this->SetTetQualityMeasure( VTK_QUALITY_MIN_ANGLE );
333  }
335  {
336  this->SetTetQualityMeasure( VTK_QUALITY_COLLAPSE_RATIO );
337  }
339  {
340  this->SetTetQualityMeasure( VTK_QUALITY_ASPECT_BETA );
341  }
343  {
344  this->SetTetQualityMeasure( VTK_QUALITY_ASPECT_GAMMA );
345  }
347  {
348  this->SetTetQualityMeasure( VTK_QUALITY_VOLUME );
349  }
351  {
352  this->SetTetQualityMeasure( VTK_QUALITY_CONDITION );
353  }
355  {
356  this->SetTetQualityMeasure( VTK_QUALITY_JACOBIAN );
357  }
359  {
360  this->SetTetQualityMeasure( VTK_QUALITY_SCALED_JACOBIAN );
361  }
363  {
364  this->SetTetQualityMeasure( VTK_QUALITY_SHAPE );
365  }
367  {
368  this->SetTetQualityMeasure( VTK_QUALITY_RELATIVE_SIZE_SQUARED );
369  }
371  {
372  this->SetTetQualityMeasure( VTK_QUALITY_SHAPE_AND_SIZE );
373  }
375  {
376  this->SetTetQualityMeasure( VTK_QUALITY_DISTORTION );
377  }
379 
381 
392  vtkSetMacro(HexQualityMeasure,int);
393  vtkGetMacro(HexQualityMeasure,int);
395  {
396  this->SetHexQualityMeasure( VTK_QUALITY_EDGE_RATIO );
397  }
399  {
400  this->SetHexQualityMeasure( VTK_QUALITY_MED_ASPECT_FROBENIUS );
401  }
403  {
404  this->SetHexQualityMeasure( VTK_QUALITY_MAX_ASPECT_FROBENIUS );
405  }
407  {
408  this->SetHexQualityMeasure( VTK_QUALITY_MAX_EDGE_RATIO );
409  }
411  {
412  this->SetHexQualityMeasure( VTK_QUALITY_SKEW );
413  }
415  {
416  this->SetHexQualityMeasure( VTK_QUALITY_TAPER );
417  }
419  {
420  this->SetHexQualityMeasure( VTK_QUALITY_VOLUME );
421  }
423  {
424  this->SetHexQualityMeasure( VTK_QUALITY_STRETCH );
425  }
427  {
428  this->SetHexQualityMeasure( VTK_QUALITY_DIAGONAL );
429  }
431  {
432  this->SetHexQualityMeasure( VTK_QUALITY_DIMENSION );
433  }
435  {
436  this->SetHexQualityMeasure( VTK_QUALITY_ODDY );
437  }
439  {
440  this->SetHexQualityMeasure( VTK_QUALITY_CONDITION );
441  }
443  {
444  this->SetHexQualityMeasure( VTK_QUALITY_JACOBIAN );
445  }
447  {
448  this->SetHexQualityMeasure( VTK_QUALITY_SCALED_JACOBIAN );
449  }
451  {
452  this->SetHexQualityMeasure( VTK_QUALITY_SHEAR );
453  }
455  {
456  this->SetHexQualityMeasure( VTK_QUALITY_SHAPE );
457  }
459  {
460  this->SetHexQualityMeasure( VTK_QUALITY_RELATIVE_SIZE_SQUARED );
461  }
463  {
464  this->SetHexQualityMeasure( VTK_QUALITY_SHAPE_AND_SIZE );
465  }
467  {
468  this->SetHexQualityMeasure( VTK_QUALITY_SHEAR_AND_SIZE );
469  }
471  {
472  this->SetHexQualityMeasure( VTK_QUALITY_DISTORTION );
473  }
475 
482  static double TriangleArea( vtkCell* cell );
483 
494  static double TriangleEdgeRatio( vtkCell* cell );
495 
506  static double TriangleAspectRatio( vtkCell* cell );
507 
518  static double TriangleRadiusRatio( vtkCell* cell );
519 
532  static double TriangleAspectFrobenius( vtkCell* cell );
533 
541  static double TriangleMinAngle( vtkCell* cell );
542 
550  static double TriangleMaxAngle( vtkCell* cell );
551 
559  static double TriangleCondition( vtkCell* cell );
560 
567  static double TriangleScaledJacobian( vtkCell* cell );
568 
575  static double TriangleRelativeSizeSquared( vtkCell* cell );
576 
583  static double TriangleShape( vtkCell* cell );
584 
591  static double TriangleShapeAndSize( vtkCell* cell );
592 
599  static double TriangleDistortion( vtkCell* cell );
600 
611  static double QuadEdgeRatio( vtkCell* cell );
612 
624  static double QuadAspectRatio( vtkCell* cell );
625 
640  static double QuadRadiusRatio( vtkCell* cell );
641 
655  static double QuadMedAspectFrobenius( vtkCell* cell );
656 
670  static double QuadMaxAspectFrobenius( vtkCell* cell );
671 
679  static double QuadMinAngle( vtkCell* cell );
680 
681  static double QuadMaxEdgeRatios( vtkCell* cell );
682  static double QuadSkew( vtkCell* cell );
683  static double QuadTaper( vtkCell* cell );
684  static double QuadWarpage( vtkCell* cell );
685  static double QuadArea( vtkCell* cell );
686  static double QuadStretch( vtkCell* cell );
687  static double QuadMaxAngle( vtkCell* cell );
688  static double QuadOddy( vtkCell* cell );
689  static double QuadCondition( vtkCell* cell );
690  static double QuadJacobian( vtkCell* cell );
691  static double QuadScaledJacobian( vtkCell* cell );
692  static double QuadShear( vtkCell* cell );
693  static double QuadShape( vtkCell* cell );
694  static double QuadRelativeSizeSquared( vtkCell* cell );
695  static double QuadShapeAndSize( vtkCell* cell );
696  static double QuadShearAndSize( vtkCell* cell );
697  static double QuadDistortion( vtkCell* cell );
698 
709  static double TetEdgeRatio( vtkCell* cell );
710 
721  static double TetAspectRatio( vtkCell* cell );
722 
733  static double TetRadiusRatio( vtkCell* cell );
734 
748  static double TetAspectFrobenius( vtkCell* cell );
749 
757  static double TetMinAngle( vtkCell* cell );
758 
760 
769  static double TetCollapseRatio( vtkCell* cell );
770  static double TetAspectBeta( vtkCell* cell );
771  static double TetAspectGamma( vtkCell* cell );
772  static double TetVolume( vtkCell* cell );
773  static double TetCondition( vtkCell* cell );
774  static double TetJacobian( vtkCell* cell );
775  static double TetScaledJacobian( vtkCell* cell );
776  static double TetShape( vtkCell* cell );
777  static double TetRelativeSizeSquared( vtkCell* cell );
778  static double TetShapeandSize( vtkCell* cell );
779  static double TetDistortion( vtkCell* cell );
781 
792  static double HexEdgeRatio( vtkCell* cell );
793 
802  static double HexMedAspectFrobenius( vtkCell* cell );
803 
805 
813  static double HexMaxAspectFrobenius( vtkCell* cell );
814  static double HexMaxEdgeRatio( vtkCell* cell );
815  static double HexSkew( vtkCell* cell );
816  static double HexTaper( vtkCell* cell );
817  static double HexVolume( vtkCell* cell );
818  static double HexStretch( vtkCell* cell );
819  static double HexDiagonal( vtkCell* cell );
820  static double HexDimension( vtkCell* cell );
821  static double HexOddy( vtkCell* cell );
822  static double HexCondition( vtkCell* cell );
823  static double HexJacobian( vtkCell* cell );
824  static double HexScaledJacobian( vtkCell* cell );
825  static double HexShear( vtkCell* cell );
826  static double HexShape( vtkCell* cell );
827  static double HexRelativeSizeSquared( vtkCell* cell );
828  static double HexShapeAndSize( vtkCell* cell );
829  static double HexShearAndSize( vtkCell* cell );
830  static double HexDistortion( vtkCell* cell );
832 
843  virtual void SetRatio( int r ) { this->SetSaveCellQuality( r ); }
844  int GetRatio() { return this->GetSaveCellQuality(); }
845  vtkBooleanMacro(Ratio,int);
846 
848 
865  virtual void SetVolume( int cv )
866  {
867  if ( ! ((cv != 0) ^ (this->Volume != 0)) )
868  {
869  return;
870  }
871  this->Modified();
872  this->Volume = cv;
873  if ( this->Volume )
874  {
875  this->CompatibilityModeOn();
876  }
877  }
878  int GetVolume()
879  {
880  return this->Volume;
881  }
882  vtkBooleanMacro(Volume,int);
884 
886 
913  virtual void SetCompatibilityMode( int cm )
914  {
915  if ( !((cm != 0) ^ (this->CompatibilityMode != 0)) )
916  {
917  return;
918  }
919  this->CompatibilityMode = cm;
920  this->Modified();
921  if ( this->CompatibilityMode )
922  {
923  this->Volume = 1;
924  this->TetQualityMeasure = VTK_QUALITY_RADIUS_RATIO;
925  }
926  }
927  vtkGetMacro(CompatibilityMode,int);
928  vtkBooleanMacro(CompatibilityMode,int);
930 
931 protected:
934 
936 
940  static int GetCurrentTriangleNormal( double point[3], double normal[3] );
941 
947 
949  int Volume;
950 
952  static double CurrentTriNormal[3];
953 
954 private:
955  vtkMeshQuality( const vtkMeshQuality& ) VTK_DELETE_FUNCTION;
956  void operator = ( const vtkMeshQuality& ) VTK_DELETE_FUNCTION;
957 };
958 
959 #endif // vtkMeshQuality_h
abstract class to specify cell behavior
Definition: vtkCell.h:60
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
Superclass for algorithms that produce output of the same type as input.
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Calculate functions of quality of the elements of a mesh.
void SetQuadQualityMeasureToSkew()
static double QuadWarpage(vtkCell *cell)
void SetHexQualityMeasureToScaledJacobian()
static double QuadTaper(vtkCell *cell)
void SetTriangleQualityMeasureToArea()
static double HexOddy(vtkCell *cell)
static double TetAspectGamma(vtkCell *cell)
void SetTetQualityMeasureToMinAngle()
void SetHexQualityMeasureToDistortion()
static double HexShear(vtkCell *cell)
static double HexScaledJacobian(vtkCell *cell)
void SetTriangleQualityMeasureToRadiusRatio()
void SetQuadQualityMeasureToRelativeSizeSquared()
void SetHexQualityMeasureToDiagonal()
void SetHexQualityMeasureToMedAspectFrobenius()
void SetTetQualityMeasureToDistortion()
virtual void SetRatio(int r)
These methods are deprecated.
static double QuadJacobian(vtkCell *cell)
void SetTriangleQualityMeasureToScaledJacobian()
void SetQuadQualityMeasureToMaxAspectFrobenius()
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
static double TriangleEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a triangle.
void SetHexQualityMeasureToMaxAspectFrobenius()
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...
virtual void SetCompatibilityMode(int cm)
CompatibilityMode governs whether, when both a quality function and cell volume are to be stored as c...
static double QuadScaledJacobian(vtkCell *cell)
static double HexMaxEdgeRatio(vtkCell *cell)
void SetHexQualityMeasureToJacobian()
static double HexShapeAndSize(vtkCell *cell)
static double QuadShear(vtkCell *cell)
void SetHexQualityMeasureToShape()
void SetQuadQualityMeasureToAspectRatio()
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.
static double HexTaper(vtkCell *cell)
vtkDataArray * CellNormals
static double TetAspectBeta(vtkCell *cell)
static double TetJacobian(vtkCell *cell)
static double HexDistortion(vtkCell *cell)
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)
void SetHexQualityMeasureToShear()
static double QuadArea(vtkCell *cell)
void SetTriangleQualityMeasureToAspectRatio()
void SetHexQualityMeasureToStretch()
static double TetEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a tetrahedron.
void SetQuadQualityMeasureToCondition()
static double HexRelativeSizeSquared(vtkCell *cell)
void SetTetQualityMeasureToShape()
void SetQuadQualityMeasureToShear()
static double HexShape(vtkCell *cell)
void SetHexQualityMeasureToCondition()
static double TriangleScaledJacobian(vtkCell *cell)
This is a static function used to calculate the scaled Jacobian of a triangle.
void SetHexQualityMeasureToDimension()
static double TetRelativeSizeSquared(vtkCell *cell)
static double TriangleShape(vtkCell *cell)
This is a static function used to calculate the shape of a triangle.
void SetTetQualityMeasureToVolume()
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.
virtual void SetVolume(int cv)
These methods are deprecated.
void SetQuadQualityMeasureToTaper()
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)
void SetTetQualityMeasureToAspectRatio()
void SetTriangleQualityMeasureToRelativeSizeSquared()
void SetHexQualityMeasureToRelativeSizeSquared()
void SetQuadQualityMeasureToEdgeRatio()
static double QuadCondition(vtkCell *cell)
void SetQuadQualityMeasureToRadiusRatio()
void SetQuadQualityMeasureToShapeAndSize()
void SetHexQualityMeasureToSkew()
void SetQuadQualityMeasureToShearAndSize()
void SetQuadQualityMeasureToJacobian()
static double TriangleAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a triangle.
static double HexJacobian(vtkCell *cell)
void SetQuadQualityMeasureToMedAspectFrobenius()
static double TetDistortion(vtkCell *cell)
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()
static double TetCondition(vtkCell *cell)
static double TriangleRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a triangle.
static double TetScaledJacobian(vtkCell *cell)
static double QuadEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a quadrilateral.
void SetTetQualityMeasureToAspectFrobenius()
static double QuadShearAndSize(vtkCell *cell)
static double HexDimension(vtkCell *cell)
void SetHexQualityMeasureToTaper()
void SetTriangleQualityMeasureToShapeAndSize()
void SetQuadQualityMeasureToShape()
static double HexStretch(vtkCell *cell)
static double QuadMaxAngle(vtkCell *cell)
static double QuadStretch(vtkCell *cell)
void SetTriangleQualityMeasureToAspectFrobenius()
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called within ProcessRequest when a request asks the algorithm to do its work.
static double HexEdgeRatio(vtkCell *cell)
This is a static function used to calculate the edge ratio of a hexahedron.
void SetTetQualityMeasureToRelativeSizeSquared()
void SetHexQualityMeasureToShearAndSize()
static double QuadMaxAspectFrobenius(vtkCell *cell)
This is a static function used to calculate the maximal Frobenius aspect of the 4 corner triangles of...
void SetHexQualityMeasureToShapeAndSize()
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()
void SetQuadQualityMeasureToWarpage()
static double HexShearAndSize(vtkCell *cell)
static double HexSkew(vtkCell *cell)
void SetTetQualityMeasureToCollapseRatio()
void SetTriangleQualityMeasureToMaxAngle()
void SetQuadQualityMeasureToMaxAngle()
void SetTetQualityMeasureToJacobian()
static double QuadShape(vtkCell *cell)
void SetQuadQualityMeasureToMaxEdgeRatios()
static vtkMeshQuality * New()
void SetTetQualityMeasureToRadiusRatio()
void SetHexQualityMeasureToOddy()
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()
void SetHexQualityMeasureToVolume()
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()
static double QuadRelativeSizeSquared(vtkCell *cell)
static double TetVolume(vtkCell *cell)
static double QuadSkew(vtkCell *cell)
static double QuadDistortion(vtkCell *cell)
void SetTetQualityMeasureToCondition()
void SetQuadQualityMeasureToOddy()
void SetTetQualityMeasureToEdgeRatio()
void SetHexQualityMeasureToEdgeRatio()
void SetTetQualityMeasureToScaledJacobian()
void SetQuadQualityMeasureToMinAngle()
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)
void SetTriangleQualityMeasureToShape()
void SetTriangleQualityMeasureToEdgeRatio()
static double TriangleArea(vtkCell *cell)
This is a static function used to calculate the area of a triangle.
void SetTetQualityMeasureToAspectGamma()
static double QuadMinAngle(vtkCell *cell)
This is a static function used to calculate the minimal (nonoriented) angle of a quadrilateral,...
void SetQuadQualityMeasureToScaledJacobian()
void SetTriangleQualityMeasureToDistortion()
void SetTriangleQualityMeasureToMinAngle()
void SetQuadQualityMeasureToStretch()
static double TetRadiusRatio(vtkCell *cell)
This is a static function used to calculate the radius ratio of a tetrahedron.
void SetTetQualityMeasureToAspectBeta()
static double TetAspectRatio(vtkCell *cell)
This is a static function used to calculate the aspect ratio of a tetrahedron.
static double TetShape(vtkCell *cell)
void SetQuadQualityMeasureToDistortion()
static double HexCondition(vtkCell *cell)
virtual void Modified()
Update the modification time for this object.
@ point
Definition: vtkX3D.h:236
#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
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.