Classes | Macros | Typedefs | Functions
verdict.h File Reference

Header file for verdict library that calculates metrics for finite elements. Also see: Main Page. More...

#include "verdict_mangle.h"
Include dependency graph for verdict.h:

Go to the source code of this file.

Classes

struct  HexMetricVals
 A struct to hold values computed by v_hex_quality. More...
 
struct  EdgeMetricVals
 A struct to hold values computed by v_edge_quality. More...
 
struct  KnifeMetricVals
 A struct to hold values computed by v_knife_quality. More...
 
struct  QuadMetricVals
 A struct to hold values computed by v_quad_quality. More...
 
struct  PyramidMetricVals
 A struct to hold values computed by v_pyramid_quality. More...
 
struct  WedgeMetricVals
 A struct to hold values computed by v_wedge_quality. More...
 
struct  TetMetricVals
 A struct to hold values computed by v_tet_quality. More...
 
struct  TriMetricVals
 A struct to hold values computed by v_tri_quality. More...
 

Macros

#define VERDICT_VERSION   120
 
#define BUILD_SHARED_LIBS
 
#define VERDICT_SHARED_LIB
 
#define VERDICT_MANGLE
 
#define VERDICT_DBL_MIN   1.0E-30
 
#define VERDICT_DBL_MAX   1.0E+30
 
#define VERDICT_PI   3.1415926535897932384626
 
#define VERDICT_ABI_IMPORT
 
#define VERDICT_ABI_EXPORT
 
#define C_FUNC_DEF   VERDICT_EXTERN_C VERDICT_ABI_IMPORT
 
Atomic hex quality flags
#define V_HEX_MAX_EDGE_RATIO
 Request that the maximum edge ratio be computed. More...
 
#define V_HEX_SKEW
 Request that skew be computed. More...
 
#define V_HEX_TAPER
 Request that taper be computed. More...
 
#define V_HEX_VOLUME
 Request that volume be computed. More...
 
#define V_HEX_STRETCH
 Request that stretch be computed. More...
 
#define V_HEX_DIAGONAL
 Request that the diagonal be computed. More...
 
#define V_HEX_DIMENSION
 Request that the dimension be computed. More...
 
#define V_HEX_ODDY
 Request that the Oddy quality be computed. More...
 
#define V_HEX_MAX_ASPECT_FROBENIUS
 Request that maximum Frobenius aspect be computed. More...
 
#define V_HEX_CONDITION
 Request that the condition be computed. More...
 
#define V_HEX_JACOBIAN
 Request that the Jacobian be computed. More...
 
#define V_HEX_SCALED_JACOBIAN
 Request that the scaled Jacobian be computed. More...
 
#define V_HEX_SHEAR
 Request that shear be computed. More...
 
#define V_HEX_SHAPE
 Request that shape be computed. More...
 
#define V_HEX_RELATIVE_SIZE_SQUARED
 Request that the square of the relative size be computed. More...
 
#define V_HEX_SHAPE_AND_SIZE
 Request that the product of shape and size be computed. More...
 
#define V_HEX_SHEAR_AND_SIZE
 Request that the product of shear and size be computed. More...
 
#define V_HEX_DISTORTION
 Request that distortion be computed. More...
 
#define V_HEX_EDGE_RATIO
 Request that the edge ratio be computed. More...
 
#define V_HEX_MED_ASPECT_FROBENIUS
 Request that average Frobenius aspect be computed. More...
 
Composite hex quality flags

These flags are bitwise combinations of the atomic flags above.

#define V_HEX_ALL
 Request that all hex metrics be computed. More...
 
#define V_HEX_TRADITIONAL
 Request the traditional hex qualities be computed. More...
 
#define V_HEX_DIAGNOSTIC
 Request hex qualities used to diagnose mesh or solver problems. More...
 
#define V_HEX_ALGEBRAIC
 Request hex qualities that have an algebraic expression be computed. More...
 
#define V_HEX_ROBINSON
 Request hex qualities for Mr. Robinson. More...
 
Atomic tetrahedral quality flags
#define V_TET_RADIUS_RATIO
 Request that the radius ratio be computed. More...
 
#define V_TET_ASPECT_BETA
 Request that the aspect ratio beta be computed. More...
 
#define V_TET_ASPECT_GAMMA
 Request that the aspect ratio gamma be computed. More...
 
#define V_TET_VOLUME
 Request that the volume be computed. More...
 
#define V_TET_CONDITION
 Request that the condition be computed. More...
 
#define V_TET_JACOBIAN
 Request that the Jacobian be computed. More...
 
#define V_TET_SCALED_JACOBIAN
 Request that the scaled Jacobian be computed. More...
 
#define V_TET_SHAPE
 Request that the shape be computed. More...
 
#define V_TET_RELATIVE_SIZE_SQUARED
 Request that the square of the relative size be computed. More...
 
#define V_TET_SHAPE_AND_SIZE
 Request that the product of shape and size be computed. More...
 
#define V_TET_DISTORTION
 Request that the distortion be computed. More...
 
#define V_TET_EDGE_RATIO
 Request that the ratio of edge lengths be computed. More...
 
#define V_TET_ASPECT_RATIO
 Request that the aspect ratio be computed. More...
 
#define V_TET_ASPECT_FROBENIUS
 Request that the Frobenius aspect be computed. More...
 
#define V_TET_MINIMUM_ANGLE
 Request that the minimum dihedral angle be computed. More...
 
#define V_TET_COLLAPSE_RATIO
 Request that the collapse ratio be computed. More...
 
#define V_TET_EQUIVOLUME_SKEW
 Request that the equivolume skew as defined by Fluent be computed. More...
 
#define V_TET_SQUISH_INDEX
 Request that the squish index as defined by Fluent be computed. More...
 
#define V_TET_ALL
 Composite tetrahedral quality flags These flags are bitwise combinations of the atomic flags above. More...
 
#define V_TET_TRADITIONAL
 Request that the tradtional tetrahedral qualities be computed. More...
 
#define V_TET_DIAGNOSTIC
 Request hex qualities used to diagnose mesh or solver problems. More...
 
#define V_TET_ALGEBRAIC
 Request hex qualities that have an algebraic expression be computed. More...
 
Pyramid bit fields
#define V_PYRAMID_VOLUME
 
#define V_PYRAMID_ALL
 
Wedge bit fields
#define V_WEDGE_VOLUME
 
#define V_WEDGE_EDGE_RATIO   2
 
#define V_WEDGE_MAX_ASPECT_FROBENIUS   4
 
#define V_WEDGE_MEAN_ASPECT_FROBENIUS   8
 
#define V_WEDGE_JACOBIAN   16
 
#define V_WEDGE_SCALED_JACOBIAN   32
 
#define V_WEDGE_DISTORTION   64
 
#define V_WEDGE_MAX_STRETCH   128
 
#define V_WEDGE_SHAPE   256
 
#define V_WEDGE_CONDITION   512
 
#define V_WEDGE_ALL
 
Knife bit fields
#define V_KNIFE_VOLUME
 
#define V_KNIFE_ALL
 
Quad bit fields
#define V_QUAD_MAX_EDGE_RATIO
 
#define V_QUAD_SKEW
 
#define V_QUAD_TAPER
 
#define V_QUAD_WARPAGE
 
#define V_QUAD_AREA
 
#define V_QUAD_STRETCH
 
#define V_QUAD_MINIMUM_ANGLE
 
#define V_QUAD_MAXIMUM_ANGLE
 
#define V_QUAD_ODDY
 
#define V_QUAD_CONDITION
 
#define V_QUAD_JACOBIAN
 
#define V_QUAD_SCALED_JACOBIAN
 
#define V_QUAD_SHEAR
 
#define V_QUAD_SHAPE
 
#define V_QUAD_RELATIVE_SIZE_SQUARED
 
#define V_QUAD_SHAPE_AND_SIZE
 
#define V_QUAD_SHEAR_AND_SIZE
 
#define V_QUAD_DISTORTION
 
#define V_QUAD_EDGE_RATIO
 
#define V_QUAD_ASPECT_RATIO
 
#define V_QUAD_RADIUS_RATIO
 
#define V_QUAD_MED_ASPECT_FROBENIUS
 
#define V_QUAD_MAX_ASPECT_FROBENIUS
 
#define V_QUAD_ALL
 
#define V_QUAD_TRADITIONAL
 
#define V_QUAD_DIAGNOSTIC
 
#define V_QUAD_ALGEBRAIC
 
#define V_QUAD_ROBINSON
 
Tri bit fields
#define V_TRI_ASPECT_FROBENIUS
 
#define V_TRI_AREA
 
#define V_TRI_MINIMUM_ANGLE
 
#define V_TRI_MAXIMUM_ANGLE
 
#define V_TRI_CONDITION
 
#define V_TRI_SCALED_JACOBIAN
 
#define V_TRI_SHAPE
 
#define V_TRI_RELATIVE_SIZE_SQUARED
 
#define V_TRI_SHAPE_AND_SIZE
 
#define V_TRI_DISTORTION
 
#define V_TRI_RADIUS_RATIO
 
#define V_TRI_EDGE_RATIO
 
#define V_TRI_ALL
 
#define V_TRI_TRADITIONAL
 
#define V_TRI_DIAGNOSTIC
 
#define V_TRI_ALGEBRAIC
 
#define V_EDGE_LENGTH
 
#define V_EDGE_ALL
 

Typedefs

typedef double(* VerdictFunction) (int, double[][3])
 
typedef int(* ComputeNormal) (double point[3], double normal[3])
 

Functions

C_FUNC_DEF void v_hex_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct HexMetricVals *metric_vals)
 Calculates quality metrics for hexahedral elements. More...
 
C_FUNC_DEF void v_tet_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct TetMetricVals *metric_vals)
 Calculates quality metrics for tetrahedral elements. More...
 
C_FUNC_DEF void v_pyramid_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct PyramidMetricVals *metric_vals)
 Calculates quality metrics for pyramid elements. More...
 
C_FUNC_DEF void v_wedge_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct WedgeMetricVals *metric_vals)
 Calculates quality metrics for wedge elements. More...
 
C_FUNC_DEF void v_knife_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct KnifeMetricVals *metric_vals)
 Calculates quality metrics for knife elements. More...
 
C_FUNC_DEF void v_quad_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct QuadMetricVals *metric_vals)
 Calculates quality metrics for quadralateral elements. More...
 
C_FUNC_DEF void v_tri_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct TriMetricVals *metric_vals)
 Calculates quality metrics for triangle elements. More...
 
C_FUNC_DEF void v_edge_quality (int num_nodes, double coordinates[][3], unsigned int metrics_request_flag, struct EdgeMetricVals *metric_vals)
 Calculates quality metrics for edge elements. More...
 
C_FUNC_DEF void v_set_hex_size (double size)
 Sets average size (volume) of hex, needed for v_hex_relative_size(...) More...
 
C_FUNC_DEF double v_hex_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates hex edge ratio metric. More...
 
C_FUNC_DEF double v_hex_max_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates hex maximum of edge ratio. More...
 
C_FUNC_DEF double v_hex_skew (int num_nodes, double coordinates[][3])
 Calculates hex skew metric. More...
 
C_FUNC_DEF double v_hex_taper (int num_nodes, double coordinates[][3])
 Calculates hex taper metric. More...
 
C_FUNC_DEF double v_hex_volume (int num_nodes, double coordinates[][3])
 Calculates hex volume. More...
 
C_FUNC_DEF double v_hex_stretch (int num_nodes, double coordinates[][3])
 Calculates hex stretch metric. More...
 
C_FUNC_DEF double v_hex_diagonal (int num_nodes, double coordinates[][3])
 Calculates hex diagonal metric. More...
 
C_FUNC_DEF double v_hex_dimension (int num_nodes, double coordinates[][3])
 Calculates hex dimension metric. More...
 
C_FUNC_DEF double v_hex_oddy (int num_nodes, double coordinates[][3])
 Calculates hex oddy metric. More...
 
C_FUNC_DEF double v_hex_med_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates hex condition metric. More...
 
C_FUNC_DEF double v_hex_max_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates hex condition metric. More...
 
C_FUNC_DEF double v_hex_condition (int num_nodes, double coordinates[][3])
 Calculates hex condition metric. This is a synonym for v_hex_max_aspect_frobenius. More...
 
C_FUNC_DEF double v_hex_jacobian (int num_nodes, double coordinates[][3])
 Calculates hex jacobian metric. More...
 
C_FUNC_DEF double v_hex_scaled_jacobian (int num_nodes, double coordinates[][3])
 Calculates hex scaled jacobian metric. More...
 
C_FUNC_DEF double v_hex_shear (int num_nodes, double coordinates[][3])
 Calculates hex shear metric. More...
 
C_FUNC_DEF double v_hex_shape (int num_nodes, double coordinates[][3])
 Calculates hex shape metric. More...
 
C_FUNC_DEF double v_hex_relative_size_squared (int num_nodes, double coordinates[][3])
 Calculates hex relative size metric. More...
 
C_FUNC_DEF double v_hex_shape_and_size (int num_nodes, double coordinates[][3])
 Calculates hex shape-size metric. More...
 
C_FUNC_DEF double v_hex_shear_and_size (int num_nodes, double coordinates[][3])
 Calculates hex shear-size metric. More...
 
C_FUNC_DEF double v_hex_distortion (int num_nodes, double coordinates[][3])
 Calculates hex distortion metric. More...
 
C_FUNC_DEF void v_set_tet_size (double size)
 Sets average size (volume) of tet, needed for v_tet_relative_size(...) More...
 
C_FUNC_DEF double v_tet_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates tet edge ratio metric. More...
 
C_FUNC_DEF double v_tet_radius_ratio (int num_nodes, double coordinates[][3])
 Calculates tet radius ratio metric. More...
 
C_FUNC_DEF double v_tet_aspect_beta (int num_nodes, double coordinates[][3])
 Calculates the radius ratio metric of a positively oriented tet. More...
 
C_FUNC_DEF double v_tet_aspect_ratio (int num_nodes, double coordinates[][3])
 Calculates tet aspect ratio metric. More...
 
C_FUNC_DEF double v_tet_aspect_gamma (int num_nodes, double coordinates[][3])
 Calculates tet aspect gamma metric. More...
 
C_FUNC_DEF double v_tet_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates tet aspect frobenius metric. More...
 
C_FUNC_DEF double v_tet_minimum_angle (int num_nodes, double coordinates[][3])
 Calculates tet minimum dihedral angle. More...
 
C_FUNC_DEF double v_tet_collapse_ratio (int num_nodes, double coordinates[][3])
 Calculates tet collapse ratio metric. More...
 
C_FUNC_DEF double v_tet_volume (int num_nodes, double coordinates[][3])
 Calculates tet volume. More...
 
C_FUNC_DEF double v_tet_condition (int num_nodes, double coordinates[][3])
 Calculates tet condition metric. More...
 
C_FUNC_DEF double v_tet_jacobian (int num_nodes, double coordinates[][3])
 Calculates tet jacobian. More...
 
C_FUNC_DEF double v_tet_scaled_jacobian (int num_nodes, double coordinates[][3])
 Calculates tet scaled jacobian. More...
 
C_FUNC_DEF double v_tet_shape (int num_nodes, double coordinates[][3])
 Calculates tet shape metric. More...
 
C_FUNC_DEF double v_tet_relative_size_squared (int num_nodes, double coordinates[][3])
 Calculates tet relative size metric. More...
 
C_FUNC_DEF double v_tet_shape_and_size (int num_nodes, double coordinates[][3])
 Calculates tet shape-size metric. More...
 
C_FUNC_DEF double v_tet_distortion (int num_nodes, double coordinates[][3])
 Calculates tet distortion metric. More...
 
C_FUNC_DEF double v_pyramid_volume (int num_nodes, double coordinates[][3])
 Calculates pyramid volume. More...
 
C_FUNC_DEF double v_wedge_volume (int num_nodes, double coordinates[][3])
 Calculates wedge volume. More...
 
C_FUNC_DEF double v_wedge_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates wedge edge ratio metric. More...
 
C_FUNC_DEF double v_wedge_max_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates wedge max aspect forbenius. More...
 
C_FUNC_DEF double v_wedge_mean_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates wedge mean aspect forbenius. More...
 
C_FUNC_DEF double v_wedge_jacobian (int num_nodes, double coordinates[][3])
 Calculates wedge jacobian metric. More...
 
C_FUNC_DEF double v_wedge_distortion (int num_nodes, double coordinates[][3])
 Calculates wedge distortion metric. More...
 
C_FUNC_DEF double v_wedge_max_stretch (int num_nodes, double coordinates[][3])
 Calculates the wedge stretch. More...
 
C_FUNC_DEF double v_wedge_scaled_jacobian (int num_nodes, double coordinates[][3])
 Calculates wedge scaled jacobian metric. More...
 
C_FUNC_DEF double v_wedge_shape (int num_nodes, double coordinates[][3])
 Calculates the wedge shape metric. More...
 
C_FUNC_DEF double v_wedge_condition (int num_nodes, double coordinates[][3])
 Calculates wedge max aspect forbenius. More...
 
C_FUNC_DEF double v_knife_volume (int num_nodes, double coordinates[][3])
 Calculates knife volume. More...
 
C_FUNC_DEF double v_edge_length (int num_nodes, double coordinates[][3])
 Calculates edge length. More...
 
C_FUNC_DEF void v_set_quad_size (double size)
 Sets average size (area) of quad, needed for v_quad_relative_size(...) More...
 
C_FUNC_DEF double v_quad_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates quad edge ratio. More...
 
C_FUNC_DEF double v_quad_max_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates quad maximum of edge ratio. More...
 
C_FUNC_DEF double v_quad_aspect_ratio (int num_nodes, double coordinates[][3])
 Calculates quad aspect ratio. More...
 
C_FUNC_DEF double v_quad_radius_ratio (int num_nodes, double coordinates[][3])
 Calculates quad radius ratio. More...
 
C_FUNC_DEF double v_quad_med_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates quad average Frobenius aspect. More...
 
C_FUNC_DEF double v_quad_max_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates quad maximum Frobenius aspect. More...
 
C_FUNC_DEF double v_quad_skew (int num_nodes, double coordinates[][3])
 Calculates quad skew metric. More...
 
C_FUNC_DEF double v_quad_taper (int num_nodes, double coordinates[][3])
 Calculates quad taper metric. More...
 
C_FUNC_DEF double v_quad_warpage (int num_nodes, double coordinates[][3])
 Calculates quad warpage metric. More...
 
C_FUNC_DEF double v_quad_area (int num_nodes, double coordinates[][3])
 Calculates quad area. More...
 
C_FUNC_DEF double v_quad_stretch (int num_nodes, double coordinates[][3])
 Calculates quad strech metric. More...
 
C_FUNC_DEF double v_quad_minimum_angle (int num_nodes, double coordinates[][3])
 Calculates quad's smallest angle. More...
 
C_FUNC_DEF double v_quad_maximum_angle (int num_nodes, double coordinates[][3])
 Calculates quad's largest angle. More...
 
C_FUNC_DEF double v_quad_oddy (int num_nodes, double coordinates[][3])
 Calculates quad oddy metric. More...
 
C_FUNC_DEF double v_quad_condition (int num_nodes, double coordinates[][3])
 Calculates quad condition number metric. More...
 
C_FUNC_DEF double v_quad_jacobian (int num_nodes, double coordinates[][3])
 Calculates quad jacobian. More...
 
C_FUNC_DEF double v_quad_scaled_jacobian (int num_nodes, double coordinates[][3])
 Calculates quad scaled jacobian. More...
 
C_FUNC_DEF double v_quad_shear (int num_nodes, double coordinates[][3])
 Calculates quad shear metric. More...
 
C_FUNC_DEF double v_quad_shape (int num_nodes, double coordinates[][3])
 Calculates quad shape metric. More...
 
C_FUNC_DEF double v_quad_relative_size_squared (int num_nodes, double coordinates[][3])
 Calculates quad relative size metric. More...
 
C_FUNC_DEF double v_quad_shape_and_size (int num_nodes, double coordinates[][3])
 Calculates quad shape-size metric. More...
 
C_FUNC_DEF double v_quad_shear_and_size (int num_nodes, double coordinates[][3])
 Calculates quad shear-size metric. More...
 
C_FUNC_DEF double v_quad_distortion (int num_nodes, double coordinates[][3])
 Calculates quad distortion metric. More...
 
C_FUNC_DEF void v_set_tri_size (double size)
 Sets average size (area) of tri, needed for v_tri_relative_size(...) More...
 
C_FUNC_DEF void v_set_tri_normal_func (ComputeNormal func)
 Sets function pointer to calculate triangle normal wrt surface. More...
 
C_FUNC_DEF double v_tri_edge_ratio (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_aspect_ratio (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_radius_ratio (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_aspect_frobenius (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_area (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_minimum_angle (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_maximum_angle (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_condition (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_scaled_jacobian (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_shear (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_relative_size_squared (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_shape (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_shape_and_size (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 
C_FUNC_DEF double v_tri_distortion (int num_nodes, double coordinates[][3])
 Calculates triangle metric. More...
 

Detailed Description

Header file for verdict library that calculates metrics for finite elements. Also see: Main Page.

verdict.h is the header file for applications/libraries to include to compute quality metrics.

This file is part of VERDICT

Definition in file verdict.h.

Macro Definition Documentation

◆ VERDICT_VERSION

#define VERDICT_VERSION   120

Definition at line 34 of file verdict.h.

◆ BUILD_SHARED_LIBS

#define BUILD_SHARED_LIBS

Definition at line 36 of file verdict.h.

◆ VERDICT_SHARED_LIB

#define VERDICT_SHARED_LIB

Definition at line 38 of file verdict.h.

◆ VERDICT_MANGLE

#define VERDICT_MANGLE

Definition at line 71 of file verdict.h.

◆ VERDICT_DBL_MIN

#define VERDICT_DBL_MIN   1.0E-30

Definition at line 76 of file verdict.h.

◆ VERDICT_DBL_MAX

#define VERDICT_DBL_MAX   1.0E+30

Definition at line 77 of file verdict.h.

◆ VERDICT_PI

#define VERDICT_PI   3.1415926535897932384626

Definition at line 78 of file verdict.h.

◆ VERDICT_ABI_IMPORT

#define VERDICT_ABI_IMPORT

Definition at line 87 of file verdict.h.

◆ VERDICT_ABI_EXPORT

#define VERDICT_ABI_EXPORT

Definition at line 88 of file verdict.h.

◆ C_FUNC_DEF

#define C_FUNC_DEF   VERDICT_EXTERN_C VERDICT_ABI_IMPORT

Definition at line 99 of file verdict.h.

◆ V_HEX_MAX_EDGE_RATIO

#define V_HEX_MAX_EDGE_RATIO

Request that the maximum edge ratio be computed.

See also
v_hex_max_edge_ratio

Definition at line 388 of file verdict.h.

◆ V_HEX_SKEW

#define V_HEX_SKEW

Request that skew be computed.

See also
v_hex_skew

Definition at line 390 of file verdict.h.

◆ V_HEX_TAPER

#define V_HEX_TAPER

Request that taper be computed.

See also
v_hex_taper

Definition at line 392 of file verdict.h.

◆ V_HEX_VOLUME

#define V_HEX_VOLUME

Request that volume be computed.

See also
v_hex_volume

Definition at line 394 of file verdict.h.

◆ V_HEX_STRETCH

#define V_HEX_STRETCH

Request that stretch be computed.

See also
v_hex_stretch

Definition at line 396 of file verdict.h.

◆ V_HEX_DIAGONAL

#define V_HEX_DIAGONAL

Request that the diagonal be computed.

See also
v_hex_diagonal

Definition at line 398 of file verdict.h.

◆ V_HEX_DIMENSION

#define V_HEX_DIMENSION

Request that the dimension be computed.

See also
v_hex_dimension

Definition at line 400 of file verdict.h.

◆ V_HEX_ODDY

#define V_HEX_ODDY

Request that the Oddy quality be computed.

See also
v_hex_oddy

Definition at line 402 of file verdict.h.

◆ V_HEX_MAX_ASPECT_FROBENIUS

#define V_HEX_MAX_ASPECT_FROBENIUS

Request that maximum Frobenius aspect be computed.

See also
v_hex_condition

Definition at line 404 of file verdict.h.

◆ V_HEX_CONDITION

#define V_HEX_CONDITION

Request that the condition be computed.

See also
v_hex_condition

Definition at line 406 of file verdict.h.

◆ V_HEX_JACOBIAN

#define V_HEX_JACOBIAN

Request that the Jacobian be computed.

See also
v_hex_jacobian

Definition at line 408 of file verdict.h.

◆ V_HEX_SCALED_JACOBIAN

#define V_HEX_SCALED_JACOBIAN

Request that the scaled Jacobian be computed.

See also
v_hex_scaled_jacobian

Definition at line 410 of file verdict.h.

◆ V_HEX_SHEAR

#define V_HEX_SHEAR

Request that shear be computed.

See also
v_hex_shear

Definition at line 412 of file verdict.h.

◆ V_HEX_SHAPE

#define V_HEX_SHAPE

Request that shape be computed.

See also
v_hex_shape

Definition at line 414 of file verdict.h.

◆ V_HEX_RELATIVE_SIZE_SQUARED

#define V_HEX_RELATIVE_SIZE_SQUARED

Request that the square of the relative size be computed.

See also
v_hex_relative_size_squared

Definition at line 416 of file verdict.h.

◆ V_HEX_SHAPE_AND_SIZE

#define V_HEX_SHAPE_AND_SIZE

Request that the product of shape and size be computed.

See also
v_hex_shape_and_size

Definition at line 418 of file verdict.h.

◆ V_HEX_SHEAR_AND_SIZE

#define V_HEX_SHEAR_AND_SIZE

Request that the product of shear and size be computed.

See also
v_hex_shear_and_size

Definition at line 420 of file verdict.h.

◆ V_HEX_DISTORTION

#define V_HEX_DISTORTION

Request that distortion be computed.

See also
v_hex_distortion

Definition at line 422 of file verdict.h.

◆ V_HEX_EDGE_RATIO

#define V_HEX_EDGE_RATIO

Request that the edge ratio be computed.

See also
v_hex_edge_ratio

Definition at line 424 of file verdict.h.

◆ V_HEX_MED_ASPECT_FROBENIUS

#define V_HEX_MED_ASPECT_FROBENIUS

Request that average Frobenius aspect be computed.

See also
v_hex_med_aspect_frobenius

Definition at line 426 of file verdict.h.

◆ V_HEX_ALL

#define V_HEX_ALL

Request that all hex metrics be computed.

Definition at line 431 of file verdict.h.

◆ V_HEX_TRADITIONAL

#define V_HEX_TRADITIONAL

Request the traditional hex qualities be computed.

This includes V_HEX_MAX_EDGE_RATIO, V_HEX_SKEW, V_HEX_TAPER V_HEX_STRETCH, V_HEX_DIAGONAL, V_HEX_ODDY, V_HEX_CONDITION, V_HEX_JACOBIAN, V_HEX_SCALED_JACOBIAN, and V_HEX_DIMENSION.

Definition at line 440 of file verdict.h.

◆ V_HEX_DIAGNOSTIC

#define V_HEX_DIAGNOSTIC

Request hex qualities used to diagnose mesh or solver problems.

This is currently a synonym for V_HEX_VOLUME but may include more qualities in the future.

Definition at line 455 of file verdict.h.

◆ V_HEX_ALGEBRAIC

#define V_HEX_ALGEBRAIC

Request hex qualities that have an algebraic expression be computed.

This includes V_HEX_SHAPE, V_HEX_SHEAR, V_HEX_RELATIVE_SIZE_SQUARED, V_HEX_SHAPE_AND_SIZE, and V_HEX_SHEAR_AND_SIZE.

Definition at line 462 of file verdict.h.

◆ V_HEX_ROBINSON

#define V_HEX_ROBINSON

Request hex qualities for Mr. Robinson.

Definition at line 470 of file verdict.h.

◆ V_TET_RADIUS_RATIO

#define V_TET_RADIUS_RATIO

Request that the radius ratio be computed.

See also
v_tet_radius_ratio

Definition at line 480 of file verdict.h.

◆ V_TET_ASPECT_BETA

#define V_TET_ASPECT_BETA

Request that the aspect ratio beta be computed.

See also
v_tet_aspect_beta

Definition at line 482 of file verdict.h.

◆ V_TET_ASPECT_GAMMA

#define V_TET_ASPECT_GAMMA

Request that the aspect ratio gamma be computed.

See also
v_tet_aspect_gamma

Definition at line 484 of file verdict.h.

◆ V_TET_VOLUME

#define V_TET_VOLUME

Request that the volume be computed.

See also
v_tet_volume

Definition at line 486 of file verdict.h.

◆ V_TET_CONDITION

#define V_TET_CONDITION

Request that the condition be computed.

See also
v_tet_condition

Definition at line 488 of file verdict.h.

◆ V_TET_JACOBIAN

#define V_TET_JACOBIAN

Request that the Jacobian be computed.

See also
v_tet_jacobian

Definition at line 490 of file verdict.h.

◆ V_TET_SCALED_JACOBIAN

#define V_TET_SCALED_JACOBIAN

Request that the scaled Jacobian be computed.

See also
v_tet_scaled_jacobian

Definition at line 492 of file verdict.h.

◆ V_TET_SHAPE

#define V_TET_SHAPE

Request that the shape be computed.

See also
v_tet_shape

Definition at line 494 of file verdict.h.

◆ V_TET_RELATIVE_SIZE_SQUARED

#define V_TET_RELATIVE_SIZE_SQUARED

Request that the square of the relative size be computed.

See also
v_tet_relative_size_squared

Definition at line 496 of file verdict.h.

◆ V_TET_SHAPE_AND_SIZE

#define V_TET_SHAPE_AND_SIZE

Request that the product of shape and size be computed.

See also
v_tet_shape_and_size

Definition at line 498 of file verdict.h.

◆ V_TET_DISTORTION

#define V_TET_DISTORTION

Request that the distortion be computed.

See also
v_tet_distortion

Definition at line 500 of file verdict.h.

◆ V_TET_EDGE_RATIO

#define V_TET_EDGE_RATIO

Request that the ratio of edge lengths be computed.

See also
v_tet_edge_ratio

Definition at line 502 of file verdict.h.

◆ V_TET_ASPECT_RATIO

#define V_TET_ASPECT_RATIO

Request that the aspect ratio be computed.

See also
v_tet_aspect_ratio

Definition at line 504 of file verdict.h.

◆ V_TET_ASPECT_FROBENIUS

#define V_TET_ASPECT_FROBENIUS

Request that the Frobenius aspect be computed.

See also
v_tet_aspect_frobenius

Definition at line 506 of file verdict.h.

◆ V_TET_MINIMUM_ANGLE

#define V_TET_MINIMUM_ANGLE

Request that the minimum dihedral angle be computed.

See also
v_tet_minimum_angle

Definition at line 508 of file verdict.h.

◆ V_TET_COLLAPSE_RATIO

#define V_TET_COLLAPSE_RATIO

Request that the collapse ratio be computed.

See also
v_tet_collapse_ratio

Definition at line 510 of file verdict.h.

◆ V_TET_EQUIVOLUME_SKEW

#define V_TET_EQUIVOLUME_SKEW

Request that the equivolume skew as defined by Fluent be computed.

See also
v_tet_equivolume_skew

Definition at line 512 of file verdict.h.

◆ V_TET_SQUISH_INDEX

#define V_TET_SQUISH_INDEX

Request that the squish index as defined by Fluent be computed.

See also
v_tet_squish_index

Definition at line 514 of file verdict.h.

◆ V_TET_ALL

#define V_TET_ALL

Composite tetrahedral quality flags These flags are bitwise combinations of the atomic flags above.

Request that the all tetrahedral qualities be computed.

Definition at line 520 of file verdict.h.

◆ V_TET_TRADITIONAL

#define V_TET_TRADITIONAL

Request that the tradtional tetrahedral qualities be computed.

This includes V_TET_ASPECT_BETA, V_TET_ASPECT_GAMMA, V_TET_CONDITION, V_TET_JACOBIAN, and V_TET_SCALED_JACOBIAN.

Definition at line 527 of file verdict.h.

◆ V_TET_DIAGNOSTIC

#define V_TET_DIAGNOSTIC

Request hex qualities used to diagnose mesh or solver problems.

This is currently a synonym for V_TET_VOLUME but may include more qualities in the future.

Definition at line 536 of file verdict.h.

◆ V_TET_ALGEBRAIC

#define V_TET_ALGEBRAIC

Request hex qualities that have an algebraic expression be computed.

This includes V_TET_SHAPE, V_TET_RELATIVE_SIZE_SQUARED, and V_TET_SHAPE_AND_SIZE.

Definition at line 542 of file verdict.h.

◆ V_PYRAMID_VOLUME

#define V_PYRAMID_VOLUME

Definition at line 553 of file verdict.h.

◆ V_PYRAMID_ALL

#define V_PYRAMID_ALL

Definition at line 555 of file verdict.h.

◆ V_WEDGE_VOLUME

#define V_WEDGE_VOLUME

Definition at line 563 of file verdict.h.

◆ V_WEDGE_EDGE_RATIO

#define V_WEDGE_EDGE_RATIO   2

Definition at line 564 of file verdict.h.

◆ V_WEDGE_MAX_ASPECT_FROBENIUS

#define V_WEDGE_MAX_ASPECT_FROBENIUS   4

Definition at line 565 of file verdict.h.

◆ V_WEDGE_MEAN_ASPECT_FROBENIUS

#define V_WEDGE_MEAN_ASPECT_FROBENIUS   8

Definition at line 566 of file verdict.h.

◆ V_WEDGE_JACOBIAN

#define V_WEDGE_JACOBIAN   16

Definition at line 567 of file verdict.h.

◆ V_WEDGE_SCALED_JACOBIAN

#define V_WEDGE_SCALED_JACOBIAN   32

Definition at line 568 of file verdict.h.

◆ V_WEDGE_DISTORTION

#define V_WEDGE_DISTORTION   64

Definition at line 569 of file verdict.h.

◆ V_WEDGE_MAX_STRETCH

#define V_WEDGE_MAX_STRETCH   128

Definition at line 570 of file verdict.h.

◆ V_WEDGE_SHAPE

#define V_WEDGE_SHAPE   256

Definition at line 571 of file verdict.h.

◆ V_WEDGE_CONDITION

#define V_WEDGE_CONDITION   512

Definition at line 572 of file verdict.h.

◆ V_WEDGE_ALL

#define V_WEDGE_ALL

Definition at line 574 of file verdict.h.

◆ V_KNIFE_VOLUME

#define V_KNIFE_VOLUME

Definition at line 582 of file verdict.h.

◆ V_KNIFE_ALL

#define V_KNIFE_ALL

Definition at line 584 of file verdict.h.

◆ V_QUAD_MAX_EDGE_RATIO

#define V_QUAD_MAX_EDGE_RATIO

Definition at line 592 of file verdict.h.

◆ V_QUAD_SKEW

#define V_QUAD_SKEW

Definition at line 594 of file verdict.h.

◆ V_QUAD_TAPER

#define V_QUAD_TAPER

Definition at line 596 of file verdict.h.

◆ V_QUAD_WARPAGE

#define V_QUAD_WARPAGE

Definition at line 598 of file verdict.h.

◆ V_QUAD_AREA

#define V_QUAD_AREA

Definition at line 600 of file verdict.h.

◆ V_QUAD_STRETCH

#define V_QUAD_STRETCH

Definition at line 602 of file verdict.h.

◆ V_QUAD_MINIMUM_ANGLE

#define V_QUAD_MINIMUM_ANGLE

Definition at line 604 of file verdict.h.

◆ V_QUAD_MAXIMUM_ANGLE

#define V_QUAD_MAXIMUM_ANGLE

Definition at line 606 of file verdict.h.

◆ V_QUAD_ODDY

#define V_QUAD_ODDY

Definition at line 608 of file verdict.h.

◆ V_QUAD_CONDITION

#define V_QUAD_CONDITION

Definition at line 610 of file verdict.h.

◆ V_QUAD_JACOBIAN

#define V_QUAD_JACOBIAN

Definition at line 612 of file verdict.h.

◆ V_QUAD_SCALED_JACOBIAN

#define V_QUAD_SCALED_JACOBIAN

Definition at line 614 of file verdict.h.

◆ V_QUAD_SHEAR

#define V_QUAD_SHEAR

Definition at line 616 of file verdict.h.

◆ V_QUAD_SHAPE

#define V_QUAD_SHAPE

Definition at line 618 of file verdict.h.

◆ V_QUAD_RELATIVE_SIZE_SQUARED

#define V_QUAD_RELATIVE_SIZE_SQUARED

Definition at line 620 of file verdict.h.

◆ V_QUAD_SHAPE_AND_SIZE

#define V_QUAD_SHAPE_AND_SIZE

Definition at line 622 of file verdict.h.

◆ V_QUAD_SHEAR_AND_SIZE

#define V_QUAD_SHEAR_AND_SIZE

Definition at line 624 of file verdict.h.

◆ V_QUAD_DISTORTION

#define V_QUAD_DISTORTION

Definition at line 626 of file verdict.h.

◆ V_QUAD_EDGE_RATIO

#define V_QUAD_EDGE_RATIO

Definition at line 628 of file verdict.h.

◆ V_QUAD_ASPECT_RATIO

#define V_QUAD_ASPECT_RATIO

Definition at line 630 of file verdict.h.

◆ V_QUAD_RADIUS_RATIO

#define V_QUAD_RADIUS_RATIO

Definition at line 632 of file verdict.h.

◆ V_QUAD_MED_ASPECT_FROBENIUS

#define V_QUAD_MED_ASPECT_FROBENIUS

Definition at line 634 of file verdict.h.

◆ V_QUAD_MAX_ASPECT_FROBENIUS

#define V_QUAD_MAX_ASPECT_FROBENIUS

Definition at line 636 of file verdict.h.

◆ V_QUAD_ALL

#define V_QUAD_ALL

Definition at line 638 of file verdict.h.

◆ V_QUAD_TRADITIONAL

#define V_QUAD_TRADITIONAL

Definition at line 640 of file verdict.h.

◆ V_QUAD_DIAGNOSTIC

#define V_QUAD_DIAGNOSTIC

Definition at line 653 of file verdict.h.

◆ V_QUAD_ALGEBRAIC

#define V_QUAD_ALGEBRAIC

Definition at line 656 of file verdict.h.

◆ V_QUAD_ROBINSON

#define V_QUAD_ROBINSON

Definition at line 662 of file verdict.h.

◆ V_TRI_ASPECT_FROBENIUS

#define V_TRI_ASPECT_FROBENIUS

Definition at line 673 of file verdict.h.

◆ V_TRI_AREA

#define V_TRI_AREA

Definition at line 675 of file verdict.h.

◆ V_TRI_MINIMUM_ANGLE

#define V_TRI_MINIMUM_ANGLE

Definition at line 677 of file verdict.h.

◆ V_TRI_MAXIMUM_ANGLE

#define V_TRI_MAXIMUM_ANGLE

Definition at line 679 of file verdict.h.

◆ V_TRI_CONDITION

#define V_TRI_CONDITION

Definition at line 681 of file verdict.h.

◆ V_TRI_SCALED_JACOBIAN

#define V_TRI_SCALED_JACOBIAN

Definition at line 683 of file verdict.h.

◆ V_TRI_SHAPE

#define V_TRI_SHAPE

Definition at line 685 of file verdict.h.

◆ V_TRI_RELATIVE_SIZE_SQUARED

#define V_TRI_RELATIVE_SIZE_SQUARED

Definition at line 687 of file verdict.h.

◆ V_TRI_SHAPE_AND_SIZE

#define V_TRI_SHAPE_AND_SIZE

Definition at line 689 of file verdict.h.

◆ V_TRI_DISTORTION

#define V_TRI_DISTORTION

Definition at line 691 of file verdict.h.

◆ V_TRI_RADIUS_RATIO

#define V_TRI_RADIUS_RATIO

Definition at line 693 of file verdict.h.

◆ V_TRI_EDGE_RATIO

#define V_TRI_EDGE_RATIO

Definition at line 695 of file verdict.h.

◆ V_TRI_ALL

#define V_TRI_ALL

Definition at line 697 of file verdict.h.

◆ V_TRI_TRADITIONAL

#define V_TRI_TRADITIONAL

Definition at line 699 of file verdict.h.

◆ V_TRI_DIAGNOSTIC

#define V_TRI_DIAGNOSTIC

Definition at line 706 of file verdict.h.

◆ V_TRI_ALGEBRAIC

#define V_TRI_ALGEBRAIC

Definition at line 709 of file verdict.h.

◆ V_EDGE_LENGTH

#define V_EDGE_LENGTH

Definition at line 719 of file verdict.h.

◆ V_EDGE_ALL

#define V_EDGE_ALL

Definition at line 721 of file verdict.h.

Typedef Documentation

◆ VerdictFunction

typedef double(* VerdictFunction) (int, double[][3])

Definition at line 111 of file verdict.h.

◆ ComputeNormal

typedef int(* ComputeNormal) (double point[3], double normal[3])

Definition at line 112 of file verdict.h.

Function Documentation

◆ v_hex_quality()

C_FUNC_DEF void v_hex_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct HexMetricVals metric_vals 
)

Calculates quality metrics for hexahedral elements.

◆ v_tet_quality()

C_FUNC_DEF void v_tet_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct TetMetricVals metric_vals 
)

Calculates quality metrics for tetrahedral elements.

◆ v_pyramid_quality()

C_FUNC_DEF void v_pyramid_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct PyramidMetricVals metric_vals 
)

Calculates quality metrics for pyramid elements.

◆ v_wedge_quality()

C_FUNC_DEF void v_wedge_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct WedgeMetricVals metric_vals 
)

Calculates quality metrics for wedge elements.

◆ v_knife_quality()

C_FUNC_DEF void v_knife_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct KnifeMetricVals metric_vals 
)

Calculates quality metrics for knife elements.

◆ v_quad_quality()

C_FUNC_DEF void v_quad_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct QuadMetricVals metric_vals 
)

Calculates quality metrics for quadralateral elements.

◆ v_tri_quality()

C_FUNC_DEF void v_tri_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct TriMetricVals metric_vals 
)

Calculates quality metrics for triangle elements.

◆ v_edge_quality()

C_FUNC_DEF void v_edge_quality ( int  num_nodes,
double  coordinates[][3],
unsigned int  metrics_request_flag,
struct EdgeMetricVals metric_vals 
)

Calculates quality metrics for edge elements.

◆ v_set_hex_size()

C_FUNC_DEF void v_set_hex_size ( double  size)

Sets average size (volume) of hex, needed for v_hex_relative_size(...)

◆ v_hex_edge_ratio()

C_FUNC_DEF double v_hex_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex edge ratio metric.

Hmax / Hmin where Hmax and Hmin are respectively the maximum and the minimum edge lengths

◆ v_hex_max_edge_ratio()

C_FUNC_DEF double v_hex_max_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex maximum of edge ratio.

Maximum edge length ratio at hex center. Reference — L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989.

◆ v_hex_skew()

C_FUNC_DEF double v_hex_skew ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex skew metric.

Maximum |cos A| where A is the angle between edges at hex center. Reference — L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989.

◆ v_hex_taper()

C_FUNC_DEF double v_hex_taper ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex taper metric.

Maximum ratio of lengths derived from opposite edges. Reference — L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989.

◆ v_hex_volume()

C_FUNC_DEF double v_hex_volume ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex volume.

Jacobian at hex center. Reference — L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989.

◆ v_hex_stretch()

C_FUNC_DEF double v_hex_stretch ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex stretch metric.

Sqrt(3) * minimum edge length / maximum diagonal length. Reference — FIMESH code

◆ v_hex_diagonal()

C_FUNC_DEF double v_hex_diagonal ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex diagonal metric.

Minimum diagonal length / maximum diagonal length. Reference — Unknown

◆ v_hex_dimension()

C_FUNC_DEF double v_hex_dimension ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex dimension metric.

Pronto-specific characteristic length for stable time step calculation. Char_length = Volume / 2 grad Volume. Reference — L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989.

◆ v_hex_oddy()

C_FUNC_DEF double v_hex_oddy ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex oddy metric.

◆ v_hex_med_aspect_frobenius()

C_FUNC_DEF double v_hex_med_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex condition metric.

Average Frobenius condition number of the Jacobian matrix at 8 corners.

◆ v_hex_max_aspect_frobenius()

C_FUNC_DEF double v_hex_max_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex condition metric.

Maximum Frobenius condition number of the Jacobian matrix at 8 corners. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_hex_condition()

C_FUNC_DEF double v_hex_condition ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex condition metric. This is a synonym for v_hex_max_aspect_frobenius.

◆ v_hex_jacobian()

C_FUNC_DEF double v_hex_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex jacobian metric.

Minimum pointwise volume of local map at 8 corners & center of hex. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_hex_scaled_jacobian()

C_FUNC_DEF double v_hex_scaled_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex scaled jacobian metric.

Minimum Jacobian divided by the lengths of the 3 edge vectors. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_hex_shear()

C_FUNC_DEF double v_hex_shear ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex shear metric.

3/Mean Ratio of Jacobian Skew matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_hex_shape()

C_FUNC_DEF double v_hex_shape ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex shape metric.

3/Mean Ratio of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_hex_relative_size_squared()

C_FUNC_DEF double v_hex_relative_size_squared ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex relative size metric.

3/Mean Ratio of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_hex_shape_and_size()

C_FUNC_DEF double v_hex_shape_and_size ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex shape-size metric.

Product of Shape and Relative Size. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_hex_shear_and_size()

C_FUNC_DEF double v_hex_shear_and_size ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex shear-size metric.

Product of Shear and Relative Size. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_hex_distortion()

C_FUNC_DEF double v_hex_distortion ( int  num_nodes,
double  coordinates[][3] 
)

Calculates hex distortion metric.

{min(|J|)/actual volume}*parent volume, parent volume = 8 for hex. Reference — SDRC/IDEAS Simulation: Finite Element Modeling–User's Guide

◆ v_set_tet_size()

C_FUNC_DEF void v_set_tet_size ( double  size)

Sets average size (volume) of tet, needed for v_tet_relative_size(...)

◆ v_tet_edge_ratio()

C_FUNC_DEF double v_tet_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet edge ratio metric.

Hmax / Hmin where Hmax and Hmin are respectively the maximum and the minimum edge lengths

◆ v_tet_radius_ratio()

C_FUNC_DEF double v_tet_radius_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet radius ratio metric.

CR / (3.0 * IR) where CR = circumsphere radius, IR = inscribed sphere radius. Reference — V. N. Parthasarathy et al, A comparison of tetrahedron quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261.

◆ v_tet_aspect_beta()

C_FUNC_DEF double v_tet_aspect_beta ( int  num_nodes,
double  coordinates[][3] 
)

Calculates the radius ratio metric of a positively oriented tet.

CR / (3.0 * IR) where CR = circumsphere radius, IR = inscribed sphere radius if the element is positively-oriented. Reference — V. N. Parthasarathy et al, A comparison of tetrahedron quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261.

◆ v_tet_aspect_ratio()

C_FUNC_DEF double v_tet_aspect_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet aspect ratio metric.

Hmax / (2 sqrt(6) r) where Hmax and r respectively denote the greatest edge length and the inradius of the tetrahedron Reference — P. Frey and P.-L. George, Meshing, Hermes (2000).

◆ v_tet_aspect_gamma()

C_FUNC_DEF double v_tet_aspect_gamma ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet aspect gamma metric.

Srms**3 / (8.479670*V) where Srms = sqrt(Sum(Si**2)/6), Si = edge length. Reference — V. N. Parthasarathy et al, A comparison of tetrahedron quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261.

◆ v_tet_aspect_frobenius()

C_FUNC_DEF double v_tet_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet aspect frobenius metric.

Frobenius condition number when the reference element is regular Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_tet_minimum_angle()

C_FUNC_DEF double v_tet_minimum_angle ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet minimum dihedral angle.

Minimum (nonoriented) dihedral angle of a tetrahedron, expressed in degrees.

◆ v_tet_collapse_ratio()

C_FUNC_DEF double v_tet_collapse_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet collapse ratio metric.

Collapse ratio

◆ v_tet_volume()

C_FUNC_DEF double v_tet_volume ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet volume.

(1/6) * Jacobian at corner node. Reference — V. N. Parthasarathy et al, A comparison of tetrahedron quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261.

◆ v_tet_condition()

C_FUNC_DEF double v_tet_condition ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet condition metric.

Condition number of the Jacobian matrix at any corner. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_tet_jacobian()

C_FUNC_DEF double v_tet_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet jacobian.

Minimum pointwise volume at any corner. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_tet_scaled_jacobian()

C_FUNC_DEF double v_tet_scaled_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet scaled jacobian.

Minimum Jacobian divided by the lengths of 3 edge vectors Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_tet_shape()

C_FUNC_DEF double v_tet_shape ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet shape metric.

3/Mean Ratio of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_tet_relative_size_squared()

C_FUNC_DEF double v_tet_relative_size_squared ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet relative size metric.

Min( J, 1/J ), where J is determinant of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_tet_shape_and_size()

C_FUNC_DEF double v_tet_shape_and_size ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet shape-size metric.

Product of Shape and Relative Size. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_tet_distortion()

C_FUNC_DEF double v_tet_distortion ( int  num_nodes,
double  coordinates[][3] 
)

Calculates tet distortion metric.

{min(|J|)/actual volume}*parent volume, parent volume = 1/6 for tet. Reference — SDRC/IDEAS Simulation: Finite Element Modeling–User's Guide

◆ v_pyramid_volume()

C_FUNC_DEF double v_pyramid_volume ( int  num_nodes,
double  coordinates[][3] 
)

Calculates pyramid volume.

◆ v_wedge_volume()

C_FUNC_DEF double v_wedge_volume ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge volume.

◆ v_wedge_edge_ratio()

C_FUNC_DEF double v_wedge_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge edge ratio metric.

Hmax / Hmin where Hmax and Hmin are respectively the maximum and the minimum edge lengths

◆ v_wedge_max_aspect_frobenius()

C_FUNC_DEF double v_wedge_max_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge max aspect forbenius.

max(F_0123, F_1204, F_2015, F_3540, F_4351, F_5432) Reference — Adaptation of hex max aspect frobenius

◆ v_wedge_mean_aspect_frobenius()

C_FUNC_DEF double v_wedge_mean_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge mean aspect forbenius.

1/6 * (F_0123 + F_1204 + F+2015 + F_3540 + F_4351 + F_5432) Reference — Adaptation of hex mean aspect frobenius

◆ v_wedge_jacobian()

C_FUNC_DEF double v_wedge_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge jacobian metric.

min{((L_2 X L_0) * L_3)_k} Reference — Adaptation of Tet jacobian metric.

◆ v_wedge_distortion()

C_FUNC_DEF double v_wedge_distortion ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge distortion metric.

{min(|J|)/actual volume}*parent volume. Reference — Adaptation of Hex distortion metric.

◆ v_wedge_max_stretch()

C_FUNC_DEF double v_wedge_max_stretch ( int  num_nodes,
double  coordinates[][3] 
)

Calculates the wedge stretch.

Minimum of the stretch of each quadrilateral face. Reference – See quadrilateral stretch

◆ v_wedge_scaled_jacobian()

C_FUNC_DEF double v_wedge_scaled_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge scaled jacobian metric.

Reference — Adaptation of Hex and Tet scaled jacobian metric.

◆ v_wedge_shape()

C_FUNC_DEF double v_wedge_shape ( int  num_nodes,
double  coordinates[][3] 
)

Calculates the wedge shape metric.

3 divided by the minimum mean ratio of the Jacobian matrix at each element corner. Reference – Adaptaation of Hex shape metric.

◆ v_wedge_condition()

C_FUNC_DEF double v_wedge_condition ( int  num_nodes,
double  coordinates[][3] 
)

Calculates wedge max aspect forbenius.

max(F_0123, F_1204, F_2015, F_3540, F_4351, F_5432) Reference — Adaptation of hex max aspect frobenius

◆ v_knife_volume()

C_FUNC_DEF double v_knife_volume ( int  num_nodes,
double  coordinates[][3] 
)

Calculates knife volume.

◆ v_edge_length()

C_FUNC_DEF double v_edge_length ( int  num_nodes,
double  coordinates[][3] 
)

Calculates edge length.

◆ v_set_quad_size()

C_FUNC_DEF void v_set_quad_size ( double  size)

Sets average size (area) of quad, needed for v_quad_relative_size(...)

◆ v_quad_edge_ratio()

C_FUNC_DEF double v_quad_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad edge ratio.

edge ratio Reference — P. P. Pebay, Planar Quadrangle Quality Measures, Eng. Comp., 2004, 20(2):157-173

◆ v_quad_max_edge_ratio()

C_FUNC_DEF double v_quad_max_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad maximum of edge ratio.

Maximum edge length ratio at quad center. Reference — J. Robinson, CRE Method of element testing and the Jacobian shape parameters, Eng. Comput., Vol 4, 1987.

◆ v_quad_aspect_ratio()

C_FUNC_DEF double v_quad_aspect_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad aspect ratio.

aspect ratio Reference — P. P. Pebay, Planar Quadrangle Quality Measures, Eng. Comp., 2004, 20(2):157-173

◆ v_quad_radius_ratio()

C_FUNC_DEF double v_quad_radius_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad radius ratio.

radius ratio Reference — P. P. Pebay, Planar Quadrangle Quality Measures, Eng. Comp., 2004, 20(2):157-173

◆ v_quad_med_aspect_frobenius()

C_FUNC_DEF double v_quad_med_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad average Frobenius aspect.

average Frobenius aspect Reference — P. P. Pebay, Planar Quadrangle Quality Measures, Eng. Comp., 2004, 20(2):157-173

◆ v_quad_max_aspect_frobenius()

C_FUNC_DEF double v_quad_max_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad maximum Frobenius aspect.

average Frobenius aspect Reference — P. P. Pebay, Planar Quadrangle Quality Measures, Eng. Comp., 2004, 20(2):157-173

◆ v_quad_skew()

C_FUNC_DEF double v_quad_skew ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad skew metric.

Maximum |cos A| where A is the angle between edges at quad center. Reference — J. Robinson, CRE Method of element testing and the Jacobian shape parameters, Eng. Comput., Vol 4, 1987.

◆ v_quad_taper()

C_FUNC_DEF double v_quad_taper ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad taper metric.

Maximum ratio of lengths derived from opposite edges. Reference — J. Robinson, CRE Method of element testing and the Jacobian shape parameters, Eng. Comput., Vol 4, 1987.

◆ v_quad_warpage()

C_FUNC_DEF double v_quad_warpage ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad warpage metric.

Cosine of Minimum Dihedral Angle formed by Planes Intersecting in Diagonals. Reference — J. Robinson, CRE Method of element testing and the Jacobian shape parameters, Eng. Comput., Vol 4, 1987.

◆ v_quad_area()

C_FUNC_DEF double v_quad_area ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad area.

Jacobian at quad center. Reference — J. Robinson, CRE Method of element testing and the Jacobian shape parameters, Eng. Comput., Vol 4, 1987.

◆ v_quad_stretch()

C_FUNC_DEF double v_quad_stretch ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad strech metric.

Sqrt(2) * minimum edge length / maximum diagonal length. Reference — FIMESH code.

◆ v_quad_minimum_angle()

C_FUNC_DEF double v_quad_minimum_angle ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad's smallest angle.

Smallest included quad angle (degrees). Reference — Unknown.

◆ v_quad_maximum_angle()

C_FUNC_DEF double v_quad_maximum_angle ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad's largest angle.

Largest included quad angle (degrees). Reference — Unknown.

◆ v_quad_oddy()

C_FUNC_DEF double v_quad_oddy ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad oddy metric.

◆ v_quad_condition()

C_FUNC_DEF double v_quad_condition ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad condition number metric.

Maximum condition number of the Jacobian matrix at 4 corners. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_quad_jacobian()

C_FUNC_DEF double v_quad_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad jacobian.

Minimum pointwise volume of local map at 4 corners & center of quad. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_quad_scaled_jacobian()

C_FUNC_DEF double v_quad_scaled_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad scaled jacobian.

Minimum Jacobian divided by the lengths of the 2 edge vectors. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_quad_shear()

C_FUNC_DEF double v_quad_shear ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad shear metric.

2/Condition number of Jacobian Skew matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_quad_shape()

C_FUNC_DEF double v_quad_shape ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad shape metric.

2/Condition number of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_quad_relative_size_squared()

C_FUNC_DEF double v_quad_relative_size_squared ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad relative size metric.

Min( J, 1/J ), where J is determinant of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_quad_shape_and_size()

C_FUNC_DEF double v_quad_shape_and_size ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad shape-size metric.

Product of Shape and Relative Size. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_quad_shear_and_size()

C_FUNC_DEF double v_quad_shear_and_size ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad shear-size metric.

Product of Shear and Relative Size. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_quad_distortion()

C_FUNC_DEF double v_quad_distortion ( int  num_nodes,
double  coordinates[][3] 
)

Calculates quad distortion metric.

{min(|J|)/actual area}*parent area, parent area = 4 for quad. Reference — SDRC/IDEAS Simulation: Finite Element Modeling–User's Guide

◆ v_set_tri_size()

C_FUNC_DEF void v_set_tri_size ( double  size)

Sets average size (area) of tri, needed for v_tri_relative_size(...)

◆ v_set_tri_normal_func()

C_FUNC_DEF void v_set_tri_normal_func ( ComputeNormal  func)

Sets function pointer to calculate triangle normal wrt surface.

◆ v_tri_edge_ratio()

C_FUNC_DEF double v_tri_edge_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

edge ratio Reference — P. P. Pebay & T. J. Baker, Analysis of Triangle Quality Measures, AMS Math. Comp., 2003, 72(244):1817-1839

◆ v_tri_aspect_ratio()

C_FUNC_DEF double v_tri_aspect_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

aspect ratio Reference — P. P. Pebay & T. J. Baker, Analysis of Triangle Quality Measures, AMS Math. Comp., 2003, 72(244):1817-1839

◆ v_tri_radius_ratio()

C_FUNC_DEF double v_tri_radius_ratio ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

radius ratio Reference — P. P. Pebay & T. J. Baker, Analysis of Triangle Quality Measures, AMS Math. Comp., 2003, 72(244):1817-1839

◆ v_tri_aspect_frobenius()

C_FUNC_DEF double v_tri_aspect_frobenius ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Frobenius aspect

◆ v_tri_area()

C_FUNC_DEF double v_tri_area ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Maximum included angle in triangle

◆ v_tri_minimum_angle()

C_FUNC_DEF double v_tri_minimum_angle ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Minimum included angle in triangle

◆ v_tri_maximum_angle()

C_FUNC_DEF double v_tri_maximum_angle ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Maximum included angle in triangle

◆ v_tri_condition()

C_FUNC_DEF double v_tri_condition ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Condition number of the Jacobian matrix. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_tri_scaled_jacobian()

C_FUNC_DEF double v_tri_scaled_jacobian ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Minimum Jacobian divided by the lengths of 2 edge vectors. Reference — P. Knupp, Achieving Finite Element Mesh Quality via Optimization of the Jacobian Matrix Norm and Associated Quantities, Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185.

◆ v_tri_shear()

C_FUNC_DEF double v_tri_shear ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

◆ v_tri_relative_size_squared()

C_FUNC_DEF double v_tri_relative_size_squared ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Min( J, 1/J ), where J is determinant of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_tri_shape()

C_FUNC_DEF double v_tri_shape ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

2/Condition number of weighted Jacobian matrix. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_tri_shape_and_size()

C_FUNC_DEF double v_tri_shape_and_size ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

Product of Shape and Relative Size. Reference — P. Knupp, Algebraic Mesh Quality Metrics for Unstructured Initial Meshes, submitted for publication.

◆ v_tri_distortion()

C_FUNC_DEF double v_tri_distortion ( int  num_nodes,
double  coordinates[][3] 
)

Calculates triangle metric.

{min(|J|)/actual area}*parent area, parent area = 1/2 for triangular element. Reference — SDRC/IDEAS Simulation: Finite Element Modeling–User's Guide