verdict.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Module: verdict.h.in
4 
5  Copyright (c) 2006 Sandia Corporation.
6  All rights reserved.
7  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
8 
9  This software is distributed WITHOUT ANY WARRANTY; without even
10  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11  PURPOSE. See the above copyright notice for more information.
12 
13 =========================================================================*/
14 
15 
27 // .SECTION Thanks
28 // Prior to its inclusion within VTK, this code was developed by the CUBIT
29 // project at Sandia National Laboratories.
30 
31 #ifndef __verdict_h
32 #define __verdict_h
33 
34 #define VERDICT_VERSION 120
35 
36 #define BUILD_SHARED_LIBS
37 #ifdef BUILD_SHARED_LIBS
38 # define VERDICT_SHARED_LIB
39 #endif
40 
71 #define VERDICT_MANGLE
72 #ifdef VERDICT_MANGLE
73 # include "verdict_mangle.h"
74 #endif
75 
76 #define VERDICT_DBL_MIN 1.0E-30
77 #define VERDICT_DBL_MAX 1.0E+30
78 #define VERDICT_PI 3.1415926535897932384626
79 
80 #if defined(_WIN32) || defined (__CYGWIN__)
81 # define VERDICT_ABI_IMPORT __declspec(dllimport)
82 # define VERDICT_ABI_EXPORT __declspec(dllexport)
83 #elif __GNUC__ >= 4
84 # define VERDICT_ABI_IMPORT __attribute__ ((visibility("default")))
85 # define VERDICT_ABI_EXPORT __attribute__ ((visibility("default")))
86 #else
87 # define VERDICT_ABI_IMPORT
88 # define VERDICT_ABI_EXPORT
89 #endif
90 
91 #ifdef __cplusplus
92 # define VERDICT_EXTERN_C extern "C"
93 #endif
94 
95 #if defined(VERDICT_SHARED_LIB)
96 # ifdef verdict_EXPORTS
97 # define C_FUNC_DEF VERDICT_EXTERN_C VERDICT_ABI_EXPORT
98 # else
99 # define C_FUNC_DEF VERDICT_EXTERN_C VERDICT_ABI_IMPORT
100 # endif
101 #else
102 # define C_FUNC_DEF VERDICT_EXTERN_C
103 #endif
104 
105 
106 /* typedef for the user if they want to use
107  * function pointers */
108 #ifdef __cplusplus
109 extern "C" {
110 #endif
111  typedef double(*VerdictFunction)(int, double[][3]);
112  typedef int(*ComputeNormal)(double point[3], double normal[3]);
113 #ifdef __cplusplus
114 }
115 #endif
116 
117 
124 {
126  double edge_ratio ;
128  double max_edge_ratio ;
130  double skew ;
132  double taper ;
134  double volume ;
136  double stretch ;
138  double diagonal ;
140  double dimension ;
142  double oddy ;
148  double condition ;
150  double jacobian ;
154  double shear ;
156  double shape ;
160  double shape_and_size ;
162  double shear_and_size ;
164  double distortion;
165 };
166 
173 {
174  double length ;
175 };
176 
177 
184 {
185  double volume ;
186 };
187 
210 {
212  double edge_ratio ;
214  double max_edge_ratio ;
216  double aspect_ratio ;
218  double radius_ratio ;
224  double skew ;
226  double taper ;
228  double warpage ;
230  double area ;
232  double stretch ;
234  double minimum_angle ;
236  double maximum_angle ;
238  double oddy ;
240  double condition ;
242  double jacobian ;
246  double shear ;
248  double shape ;
252  double shape_and_size ;
254  double shear_and_size ;
256  double distortion;
257 };
258 
265 {
266  double volume ;
267 };
268 
275 {
277  double volume ;
279  double edge_ratio ;
285  double jacobian;
287  double distortion;
289  double max_stretch;
293  double shape;
295  double condition;
296 };
297 
304 {
306  double edge_ratio;
308  double radius_ratio;
310  double aspect_beta;
312  double aspect_ratio ;
314  double aspect_gamma ;
318  double minimum_angle ;
322  double volume ;
324  double condition ;
326  double jacobian ;
330  double shape ;
334  double shape_and_size ;
336  double distortion;
340  double squish_index;
341 
342 };
343 
350 {
352  double edge_ratio ;
354  double aspect_ratio ;
356  double radius_ratio ;
360  double area ;
362  double minimum_angle ;
364  double maximum_angle ;
366  double condition ;
370  double shape ;
374  double shape_and_size ;
376  double distortion;
377 };
378 
379 
380 
381 /* definition of bit fields to determine which metrics to calculate */
382 
385 
386 
388 #define V_HEX_MAX_EDGE_RATIO 1
389 
390 #define V_HEX_SKEW 2
391 
392 #define V_HEX_TAPER 4
393 
394 #define V_HEX_VOLUME 8
395 
396 #define V_HEX_STRETCH 16
397 
398 #define V_HEX_DIAGONAL 32
399 
400 #define V_HEX_DIMENSION 64
401 
402 #define V_HEX_ODDY 128
403 
404 #define V_HEX_MAX_ASPECT_FROBENIUS 256
405 
406 #define V_HEX_CONDITION 256
407 
408 #define V_HEX_JACOBIAN 512
409 
410 #define V_HEX_SCALED_JACOBIAN 1024
411 
412 #define V_HEX_SHEAR 2048
413 
414 #define V_HEX_SHAPE 4096
415 
416 #define V_HEX_RELATIVE_SIZE_SQUARED 8192
417 
418 #define V_HEX_SHAPE_AND_SIZE 16384
419 
420 #define V_HEX_SHEAR_AND_SIZE 32768
421 
422 #define V_HEX_DISTORTION 65536
423 
424 #define V_HEX_EDGE_RATIO 131072
425 
426 #define V_HEX_MED_ASPECT_FROBENIUS 262144
427 
430 
431 #define V_HEX_ALL 524287
432 
440 #define V_HEX_TRADITIONAL \
441  V_HEX_MAX_EDGE_RATIO + \
442  V_HEX_SKEW + \
443  V_HEX_TAPER + \
444  V_HEX_STRETCH + \
445  V_HEX_DIAGONAL + \
446  V_HEX_ODDY + \
447  V_HEX_CONDITION + \
448  V_HEX_JACOBIAN + \
449  V_HEX_SCALED_JACOBIAN + \
450  V_HEX_DIMENSION
451 
455 #define V_HEX_DIAGNOSTIC V_HEX_VOLUME
456 
462 #define V_HEX_ALGEBRAIC \
463  V_HEX_SHAPE + \
464  V_HEX_SHEAR + \
465  V_HEX_RELATIVE_SIZE_SQUARED + \
466  V_HEX_SHAPE_AND_SIZE + \
467  V_HEX_SHEAR_AND_SIZE
468 
470 #define V_HEX_ROBINSON \
471  V_HEX_SKEW + \
472  V_HEX_TAPER
473 
474 
475 
478 
479 
480 #define V_TET_RADIUS_RATIO 1
481 
482 #define V_TET_ASPECT_BETA 2
483 
484 #define V_TET_ASPECT_GAMMA 4
485 
486 #define V_TET_VOLUME 8
487 
488 #define V_TET_CONDITION 16
489 
490 #define V_TET_JACOBIAN 32
491 
492 #define V_TET_SCALED_JACOBIAN 64
493 
494 #define V_TET_SHAPE 128
495 
496 #define V_TET_RELATIVE_SIZE_SQUARED 256
497 
498 #define V_TET_SHAPE_AND_SIZE 512
499 
500 #define V_TET_DISTORTION 1024
501 
502 #define V_TET_EDGE_RATIO 2048
503 
504 #define V_TET_ASPECT_RATIO 4096
505 
506 #define V_TET_ASPECT_FROBENIUS 8192
507 
508 #define V_TET_MINIMUM_ANGLE 16384
509 
510 #define V_TET_COLLAPSE_RATIO 32768
511 
512 #define V_TET_EQUIVOLUME_SKEW 65536
513 
514 #define V_TET_SQUISH_INDEX 131072
515 
518 
519 
520 #define V_TET_ALL 262143
521 
527 #define V_TET_TRADITIONAL V_TET_ASPECT_BETA + \
528  V_TET_ASPECT_GAMMA + \
529  V_TET_CONDITION + \
530  V_TET_JACOBIAN + \
531  V_TET_SCALED_JACOBIAN
532 
536 #define V_TET_DIAGNOSTIC V_TET_VOLUME
537 
542 #define V_TET_ALGEBRAIC V_TET_SHAPE + \
543  V_TET_RELATIVE_SIZE_SQUARED + \
544  V_TET_SHAPE_AND_SIZE
545 
546 
547 
548 
551 
552 
553 #define V_PYRAMID_VOLUME 1
554 
555 #define V_PYRAMID_ALL \
556  V_PYRAMID_VOLUME
557 
558 
561 
562 
563 #define V_WEDGE_VOLUME 1
564 #define V_WEDGE_EDGE_RATIO 2
565 #define V_WEDGE_MAX_ASPECT_FROBENIUS 4
566 #define V_WEDGE_MEAN_ASPECT_FROBENIUS 8
567 #define V_WEDGE_JACOBIAN 16
568 #define V_WEDGE_SCALED_JACOBIAN 32
569 #define V_WEDGE_DISTORTION 64
570 #define V_WEDGE_MAX_STRETCH 128
571 #define V_WEDGE_SHAPE 256
572 #define V_WEDGE_CONDITION 512
573 
574 #define V_WEDGE_ALL \
575  V_WEDGE_VOLUME
576 
577 
580 
581 
582 #define V_KNIFE_VOLUME 1
583 
584 #define V_KNIFE_ALL \
585  V_KNIFE_VOLUME
586 
587 
590 
591 
592 #define V_QUAD_MAX_EDGE_RATIO 1
593 
594 #define V_QUAD_SKEW 2
595 
596 #define V_QUAD_TAPER 4
597 
598 #define V_QUAD_WARPAGE 8
599 
600 #define V_QUAD_AREA 16
601 
602 #define V_QUAD_STRETCH 32
603 
604 #define V_QUAD_MINIMUM_ANGLE 64
605 
606 #define V_QUAD_MAXIMUM_ANGLE 128
607 
608 #define V_QUAD_ODDY 256
609 
610 #define V_QUAD_CONDITION 512
611 
612 #define V_QUAD_JACOBIAN 1024
613 
614 #define V_QUAD_SCALED_JACOBIAN 2048
615 
616 #define V_QUAD_SHEAR 4096
617 
618 #define V_QUAD_SHAPE 8192
619 
620 #define V_QUAD_RELATIVE_SIZE_SQUARED 16384
621 
622 #define V_QUAD_SHAPE_AND_SIZE 32768
623 
624 #define V_QUAD_SHEAR_AND_SIZE 65536
625 
626 #define V_QUAD_DISTORTION 131072
627 
628 #define V_QUAD_EDGE_RATIO 262144
629 
630 #define V_QUAD_ASPECT_RATIO 524288
631 
632 #define V_QUAD_RADIUS_RATIO 1048576
633 
634 #define V_QUAD_MED_ASPECT_FROBENIUS 2097152
635 
636 #define V_QUAD_MAX_ASPECT_FROBENIUS 4194304
637 
638 #define V_QUAD_ALL 8388607
639 
640 #define V_QUAD_TRADITIONAL \
641  V_QUAD_MAX_EDGE_RATIO + \
642  V_QUAD_SKEW + \
643  V_QUAD_TAPER + \
644  V_QUAD_WARPAGE + \
645  V_QUAD_STRETCH + \
646  V_QUAD_MINIMUM_ANGLE + \
647  V_QUAD_MAXIMUM_ANGLE + \
648  V_QUAD_ODDY + \
649  V_QUAD_CONDITION + \
650  V_QUAD_JACOBIAN + \
651  V_QUAD_SCALED_JACOBIAN
652 
653 #define V_QUAD_DIAGNOSTIC \
654  V_QUAD_AREA
655 
656 #define V_QUAD_ALGEBRAIC \
657  V_QUAD_SHEAR + \
658  V_QUAD_SHAPE + \
659  V_QUAD_RELATIVE_SIZE_SQUARED + \
660  V_QUAD_SHAPE_AND_SIZE
661 
662 #define V_QUAD_ROBINSON \
663  V_QUAD_MAX_EDGE_RATIO + \
664  V_QUAD_SKEW + \
665  V_QUAD_TAPER
666 
667 
668 
671 
672 
673 #define V_TRI_ASPECT_FROBENIUS 1
674 
675 #define V_TRI_AREA 2
676 
677 #define V_TRI_MINIMUM_ANGLE 4
678 
679 #define V_TRI_MAXIMUM_ANGLE 8
680 
681 #define V_TRI_CONDITION 16
682 
683 #define V_TRI_SCALED_JACOBIAN 32
684 
685 #define V_TRI_SHAPE 64
686 
687 #define V_TRI_RELATIVE_SIZE_SQUARED 128
688 
689 #define V_TRI_SHAPE_AND_SIZE 256
690 
691 #define V_TRI_DISTORTION 512
692 
693 #define V_TRI_RADIUS_RATIO 1024
694 
695 #define V_TRI_EDGE_RATIO 2048
696 
697 #define V_TRI_ALL 4095
698 
699 #define V_TRI_TRADITIONAL \
700  V_TRI_ASPECT_FROBENIUS + \
701  V_TRI_MINIMUM_ANGLE + \
702  V_TRI_MAXIMUM_ANGLE + \
703  V_TRI_CONDITION + \
704  V_TRI_SCALED_JACOBIAN
705 
706 #define V_TRI_DIAGNOSTIC \
707  V_TRI_AREA
708 
709 #define V_TRI_ALGEBRAIC \
710  V_TRI_SHAPE + \
711  V_TRI_SHAPE_AND_SIZE + \
712  V_TRI_RELATIVE_SIZE_SQUARED
713 
714 
717 
718 
719 #define V_EDGE_LENGTH 1
720 
721 #define V_EDGE_ALL \
722  V_EDGE_LENGTH
723 
724 
725 
866  C_FUNC_DEF void v_hex_quality( int num_nodes, double coordinates[][3],
868  unsigned int metrics_request_flag, struct HexMetricVals *metric_vals );
869 
871  C_FUNC_DEF void v_tet_quality( int num_nodes, double coordinates[][3],
872  unsigned int metrics_request_flag, struct TetMetricVals *metric_vals );
873 
875  C_FUNC_DEF void v_pyramid_quality( int num_nodes, double coordinates[][3],
876  unsigned int metrics_request_flag, struct PyramidMetricVals *metric_vals );
877 
879  C_FUNC_DEF void v_wedge_quality( int num_nodes, double coordinates[][3],
880  unsigned int metrics_request_flag, struct WedgeMetricVals *metric_vals );
881 
883  C_FUNC_DEF void v_knife_quality( int num_nodes, double coordinates[][3],
884  unsigned int metrics_request_flag, struct KnifeMetricVals *metric_vals );
885 
887  C_FUNC_DEF void v_quad_quality( int num_nodes, double coordinates[][3],
888  unsigned int metrics_request_flag, struct QuadMetricVals *metric_vals );
889 
891  C_FUNC_DEF void v_tri_quality( int num_nodes, double coordinates[][3],
892  unsigned int metrics_request_flag, struct TriMetricVals *metric_vals );
893 
895  C_FUNC_DEF void v_edge_quality( int num_nodes, double coordinates[][3],
896  unsigned int metrics_request_flag, struct EdgeMetricVals *metric_vals );
897 
898 
899 
900 /* individual quality functions for hex elements */
901 
903  C_FUNC_DEF void v_set_hex_size( double size );
904 
906 
908  C_FUNC_DEF double v_hex_edge_ratio( int num_nodes, double coordinates[][3] );
909 
911 
914  C_FUNC_DEF double v_hex_max_edge_ratio( int num_nodes, double coordinates[][3] );
915 
917 
920  C_FUNC_DEF double v_hex_skew( int num_nodes, double coordinates[][3] );
921 
923 
926  C_FUNC_DEF double v_hex_taper( int num_nodes, double coordinates[][3] );
927 
929 
932  C_FUNC_DEF double v_hex_volume( int num_nodes, double coordinates[][3] );
933 
935 
937  C_FUNC_DEF double v_hex_stretch( int num_nodes, double coordinates[][3] );
938 
940 
942  C_FUNC_DEF double v_hex_diagonal( int num_nodes, double coordinates[][3] );
943 
945 
949  C_FUNC_DEF double v_hex_dimension( int num_nodes, double coordinates[][3] );
950 
952  C_FUNC_DEF double v_hex_oddy( int num_nodes, double coordinates[][3] );
953 
955 
956  C_FUNC_DEF double v_hex_med_aspect_frobenius( int num_nodes, double coordinates[][3] );
957 
959 
963  C_FUNC_DEF double v_hex_max_aspect_frobenius( int num_nodes, double coordinates[][3] );
965  C_FUNC_DEF double v_hex_condition( int num_nodes, double coordinates[][3] );
966 
968 
972  C_FUNC_DEF double v_hex_jacobian( int num_nodes, double coordinates[][3] );
973 
975 
979  C_FUNC_DEF double v_hex_scaled_jacobian( int num_nodes, double coordinates[][3] );
980 
982 
985  C_FUNC_DEF double v_hex_shear( int num_nodes, double coordinates[][3] );
986 
988 
991  C_FUNC_DEF double v_hex_shape( int num_nodes, double coordinates[][3] );
992 
994 
997  C_FUNC_DEF double v_hex_relative_size_squared( int num_nodes, double coordinates[][3] );
998 
1000 
1003  C_FUNC_DEF double v_hex_shape_and_size( int num_nodes, double coordinates[][3] );
1004 
1006 
1009  C_FUNC_DEF double v_hex_shear_and_size( int num_nodes, double coordinates[][3] );
1010 
1012 
1014  C_FUNC_DEF double v_hex_distortion( int num_nodes, double coordinates[][3] );
1015 
1016 /* individual quality functions for tet elements */
1017 
1019  C_FUNC_DEF void v_set_tet_size( double size );
1020 
1022 
1024  C_FUNC_DEF double v_tet_edge_ratio( int num_nodes, double coordinates[][3] );
1025 
1027 
1030  C_FUNC_DEF double v_tet_radius_ratio( int num_nodes, double coordinates[][3] );
1031 
1033 
1037  C_FUNC_DEF double v_tet_aspect_beta( int num_nodes, double coordinates[][3] );
1038 
1040 
1043  C_FUNC_DEF double v_tet_aspect_ratio( int num_nodes, double coordinates[][3] );
1044 
1046 
1049  C_FUNC_DEF double v_tet_aspect_gamma( int num_nodes, double coordinates[][3] );
1050 
1052 
1056  C_FUNC_DEF double v_tet_aspect_frobenius( int num_nodes, double coordinates[][3] );
1057 
1059 
1060  C_FUNC_DEF double v_tet_minimum_angle( int num_nodes, double coordinates[][3] );
1061 
1063 
1064  C_FUNC_DEF double v_tet_collapse_ratio( int num_nodes, double coordinates[][3] );
1065 
1067 
1070  C_FUNC_DEF double v_tet_volume( int num_nodes, double coordinates[][3] );
1071 
1073 
1077  C_FUNC_DEF double v_tet_condition( int num_nodes, double coordinates[][3] );
1078 
1080 
1084  C_FUNC_DEF double v_tet_jacobian( int num_nodes, double coordinates[][3] );
1085 
1087 
1091  C_FUNC_DEF double v_tet_scaled_jacobian( int num_nodes, double coordinates[][3] );
1092 
1094 
1097  C_FUNC_DEF double v_tet_shape( int num_nodes, double coordinates[][3] );
1098 
1100 
1103  C_FUNC_DEF double v_tet_relative_size_squared( int num_nodes, double coordinates[][3] );
1104 
1106 
1109  C_FUNC_DEF double v_tet_shape_and_size( int num_nodes, double coordinates[][3] );
1110 
1112 
1114  C_FUNC_DEF double v_tet_distortion( int num_nodes, double coordinates[][3] );
1115 
1116 /* individual quality functions for pyramid elements */
1117 
1119  C_FUNC_DEF double v_pyramid_volume( int num_nodes, double coordinates[][3] );
1120 
1121 
1122 /* individual quality functions for wedge elements */
1123 
1125  C_FUNC_DEF double v_wedge_volume( int num_nodes, double coordinates[][3] );
1126 
1128 
1130  C_FUNC_DEF double v_wedge_edge_ratio( int num_nodes, double coordinates[][3] );
1131 
1133 
1135  C_FUNC_DEF double v_wedge_max_aspect_frobenius( int num_nodes, double coordinates[][3] );
1136 
1138 
1140  C_FUNC_DEF double v_wedge_mean_aspect_frobenius( int num_nodes, double coordinates[][3] );
1141 
1143 
1145  C_FUNC_DEF double v_wedge_jacobian( int num_nodes, double coordinates[][3]);
1146 
1148 
1150  C_FUNC_DEF double v_wedge_distortion( int num_nodes, double coordinates[][3] );
1151 
1153 
1155  C_FUNC_DEF double v_wedge_max_stretch( int num_nodes, double coordinates[][3] );
1156 
1158 
1159  C_FUNC_DEF double v_wedge_scaled_jacobian( int num_nodes, double coordinates[][3] );
1160 
1162 
1165  C_FUNC_DEF double v_wedge_shape( int num_nodes, double coordinates[][3] );
1166 
1168 
1170  C_FUNC_DEF double v_wedge_condition( int num_nodes, double coordinates[][3] );
1171 
1172 
1173 /* individual quality functions for knife elements */
1174 
1176  C_FUNC_DEF double v_knife_volume( int num_nodes, double coordinates[][3] );
1177 
1178 
1179 /* individual quality functions for edge elements */
1180 
1182  C_FUNC_DEF double v_edge_length( int num_nodes, double coordinates[][3] );
1183 
1184 
1185 /* individual quality functions for quad elements */
1187  C_FUNC_DEF void v_set_quad_size( double size );
1188 
1190 
1193  C_FUNC_DEF double v_quad_edge_ratio( int num_nodes, double coordinates[][3] );
1194 
1196 
1199  C_FUNC_DEF double v_quad_max_edge_ratio( int num_nodes, double coordinates[][3] );
1200 
1202 
1205  C_FUNC_DEF double v_quad_aspect_ratio( int num_nodes, double coordinates[][3] );
1206 
1208 
1211  C_FUNC_DEF double v_quad_radius_ratio( int num_nodes, double coordinates[][3] );
1212 
1214 
1217  C_FUNC_DEF double v_quad_med_aspect_frobenius( int num_nodes, double coordinates[][3] );
1218 
1220 
1223  C_FUNC_DEF double v_quad_max_aspect_frobenius( int num_nodes, double coordinates[][3] );
1224 
1226 
1229  C_FUNC_DEF double v_quad_skew( int num_nodes, double coordinates[][3] );
1230 
1232 
1235  C_FUNC_DEF double v_quad_taper( int num_nodes, double coordinates[][3] );
1236 
1238 
1241  C_FUNC_DEF double v_quad_warpage( int num_nodes, double coordinates[][3] );
1242 
1244 
1247  C_FUNC_DEF double v_quad_area( int num_nodes, double coordinates[][3] );
1248 
1250 
1252  C_FUNC_DEF double v_quad_stretch( int num_nodes, double coordinates[][3] );
1253 
1255 
1257  C_FUNC_DEF double v_quad_minimum_angle( int num_nodes, double coordinates[][3] );
1258 
1260 
1262  C_FUNC_DEF double v_quad_maximum_angle( int num_nodes, double coordinates[][3] );
1263 
1265  C_FUNC_DEF double v_quad_oddy( int num_nodes, double coordinates[][3] );
1266 
1268 
1272  C_FUNC_DEF double v_quad_condition( int num_nodes, double coordinates[][3] );
1273 
1275 
1279  C_FUNC_DEF double v_quad_jacobian( int num_nodes, double coordinates[][3] );
1280 
1282 
1286  C_FUNC_DEF double v_quad_scaled_jacobian( int num_nodes, double coordinates[][3] );
1287 
1289 
1292  C_FUNC_DEF double v_quad_shear( int num_nodes, double coordinates[][3] );
1293 
1295 
1298  C_FUNC_DEF double v_quad_shape( int num_nodes, double coordinates[][3] );
1299 
1301 
1304  C_FUNC_DEF double v_quad_relative_size_squared( int num_nodes, double coordinates[][3] );
1305 
1307 
1310  C_FUNC_DEF double v_quad_shape_and_size( int num_nodes, double coordinates[][3] );
1311 
1313 
1316  C_FUNC_DEF double v_quad_shear_and_size( int num_nodes, double coordinates[][3] );
1317 
1319 
1321  C_FUNC_DEF double v_quad_distortion( int num_nodes, double coordinates[][3] );
1322 
1323 
1324 /* individual quality functions for triangle elements */
1325 
1327  C_FUNC_DEF void v_set_tri_size( double size );
1328 
1331 
1333 
1336  C_FUNC_DEF double v_tri_edge_ratio( int num_nodes, double coordinates[][3] );
1337 
1339 
1342  C_FUNC_DEF double v_tri_aspect_ratio( int num_nodes, double coordinates[][3] );
1343 
1345 
1348  C_FUNC_DEF double v_tri_radius_ratio( int num_nodes, double coordinates[][3] );
1349 
1351 
1352  C_FUNC_DEF double v_tri_aspect_frobenius( int num_nodes, double coordinates[][3] );
1353 
1355 
1356  C_FUNC_DEF double v_tri_area( int num_nodes, double coordinates[][3] );
1357 
1359 
1360  C_FUNC_DEF double v_tri_minimum_angle( int num_nodes, double coordinates[][3] );
1361 
1363 
1364  C_FUNC_DEF double v_tri_maximum_angle( int num_nodes, double coordinates[][3] );
1365 
1367 
1371  C_FUNC_DEF double v_tri_condition( int num_nodes, double coordinates[][3] );
1372 
1374 
1378  C_FUNC_DEF double v_tri_scaled_jacobian( int num_nodes, double coordinates[][3] );
1379 
1381 
1382  C_FUNC_DEF double v_tri_shear( int num_nodes, double coordinates[][3] );
1383 
1385 
1388  C_FUNC_DEF double v_tri_relative_size_squared( int num_nodes, double coordinates[][3] );
1389 
1391 
1394  C_FUNC_DEF double v_tri_shape( int num_nodes, double coordinates[][3] );
1395 
1397 
1400  C_FUNC_DEF double v_tri_shape_and_size( int num_nodes, double coordinates[][3] );
1401 
1403 
1405  C_FUNC_DEF double v_tri_distortion( int num_nodes, double coordinates[][3] );
1406 
1407 
1408 
1409 #endif /* __verdict_h */
1410 
1411 
1412 
C_FUNC_DEF double v_knife_volume(int num_nodes, double coordinates[][3])
Calculates knife volume.
C_FUNC_DEF double v_tet_volume(int num_nodes, double coordinates[][3])
Calculates tet volume.
C_FUNC_DEF double v_tri_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates triangle metric.
A struct to hold values computed by v_tri_quality.
Definition: verdict.h:349
double scaled_jacobian
Definition: verdict.h:328
C_FUNC_DEF double v_quad_warpage(int num_nodes, double coordinates[][3])
Calculates quad warpage metric.
C_FUNC_DEF void v_set_tri_normal_func(ComputeNormal func)
Sets function pointer to calculate triangle normal wrt surface.
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.
A struct to hold values computed by v_pyramid_quality.
Definition: verdict.h:264
double length
Definition: verdict.h:174
double jacobian
Definition: verdict.h:242
C_FUNC_DEF double v_tri_distortion(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double shape
Definition: verdict.h:293
C_FUNC_DEF double v_tri_aspect_ratio(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double shear_and_size
Definition: verdict.h:162
C_FUNC_DEF double v_quad_edge_ratio(int num_nodes, double coordinates[][3])
Calculates quad edge ratio.
int
double minimum_angle
Definition: verdict.h:362
double taper
Definition: verdict.h:226
C_FUNC_DEF double v_quad_distortion(int num_nodes, double coordinates[][3])
Calculates quad distortion metric.
A struct to hold values computed by v_quad_quality.
Definition: verdict.h:209
C_FUNC_DEF double v_tri_maximum_angle(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double aspect_ratio
Definition: verdict.h:312
C_FUNC_DEF double v_tri_shear(int num_nodes, double coordinates[][3])
Calculates triangle metric.
C_FUNC_DEF double v_quad_oddy(int num_nodes, double coordinates[][3])
Calculates quad oddy metric.
double area
Definition: verdict.h:230
double area
Definition: verdict.h:360
double aspect_frobenius
Definition: verdict.h:316
double condition
Definition: verdict.h:148
double edge_ratio
Definition: verdict.h:306
double scaled_jacobian
Definition: verdict.h:244
double shape_and_size
Definition: verdict.h:334
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.
C_FUNC_DEF double v_hex_diagonal(int num_nodes, double coordinates[][3])
Calculates hex diagonal metric.
double max_stretch
Definition: verdict.h:289
C_FUNC_DEF double v_tet_shape_and_size(int num_nodes, double coordinates[][3])
Calculates tet shape-size metric.
double max_aspect_frobenius
Definition: verdict.h:222
double volume
Definition: verdict.h:322
C_FUNC_DEF double v_hex_oddy(int num_nodes, double coordinates[][3])
Calculates hex oddy metric.
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.
double minimum_angle
Definition: verdict.h:234
C_FUNC_DEF double v_tri_shape(int num_nodes, double coordinates[][3])
Calculates triangle metric.
C_FUNC_DEF double v_hex_dimension(int num_nodes, double coordinates[][3])
Calculates hex dimension metric.
C_FUNC_DEF double v_quad_skew(int num_nodes, double coordinates[][3])
Calculates quad skew metric.
C_FUNC_DEF double v_quad_condition(int num_nodes, double coordinates[][3])
Calculates quad condition number metric.
C_FUNC_DEF double v_wedge_condition(int num_nodes, double coordinates[][3])
Calculates wedge max aspect forbenius.
C_FUNC_DEF double v_quad_shape_and_size(int num_nodes, double coordinates[][3])
Calculates quad shape-size metric.
C_FUNC_DEF double v_wedge_jacobian(int num_nodes, double coordinates[][3])
Calculates wedge jacobian metric.
C_FUNC_DEF double v_wedge_scaled_jacobian(int num_nodes, double coordinates[][3])
Calculates wedge scaled jacobian metric.
double max_aspect_frobenius
Definition: verdict.h:146
A struct to hold values computed by v_edge_quality.
Definition: verdict.h:172
double stretch
Definition: verdict.h:136
C_FUNC_DEF double v_wedge_mean_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates wedge mean aspect forbenius.
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.
double edge_ratio
Definition: verdict.h:126
C_FUNC_DEF double v_quad_jacobian(int num_nodes, double coordinates[][3])
Calculates quad jacobian.
double radius_ratio
Definition: verdict.h:308
C_FUNC_DEF double v_wedge_edge_ratio(int num_nodes, double coordinates[][3])
Calculates wedge edge ratio metric.
C_FUNC_DEF void v_set_tet_size(double size)
Sets average size (volume) of tet, needed for v_tet_relative_size(...)
C_FUNC_DEF double v_tet_minimum_angle(int num_nodes, double coordinates[][3])
Calculates tet minimum dihedral angle.
double dimension
Definition: verdict.h:140
C_FUNC_DEF double v_pyramid_volume(int num_nodes, double coordinates[][3])
Calculates pyramid volume.
A struct to hold values computed by v_knife_quality.
Definition: verdict.h:183
double scaled_jacobian
Definition: verdict.h:152
double squish_index
Definition: verdict.h:340
double maximum_angle
Definition: verdict.h:364
C_FUNC_DEF double v_hex_scaled_jacobian(int num_nodes, double coordinates[][3])
Calculates hex scaled jacobian metric.
double max_aspect_frobenius
Definition: verdict.h:281
C_FUNC_DEF double v_quad_stretch(int num_nodes, double coordinates[][3])
Calculates quad strech metric.
C_FUNC_DEF double v_quad_shape(int num_nodes, double coordinates[][3])
Calculates quad shape metric.
C_FUNC_DEF double v_hex_max_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates hex condition metric.
C_FUNC_DEF double v_tet_condition(int num_nodes, double coordinates[][3])
Calculates tet condition metric.
double jacobian
Definition: verdict.h:150
C_FUNC_DEF double v_quad_scaled_jacobian(int num_nodes, double coordinates[][3])
Calculates quad scaled jacobian.
double shear
Definition: verdict.h:246
C_FUNC_DEF double v_tri_edge_ratio(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double aspect_frobenius
Definition: verdict.h:358
double volume
Definition: verdict.h:134
C_FUNC_DEF double v_wedge_shape(int num_nodes, double coordinates[][3])
Calculates the wedge shape metric.
C_FUNC_DEF double v_hex_jacobian(int num_nodes, double coordinates[][3])
Calculates hex jacobian metric.
double relative_size_squared
Definition: verdict.h:332
double aspect_beta
Definition: verdict.h:310
C_FUNC_DEF double v_wedge_volume(int num_nodes, double coordinates[][3])
Calculates wedge volume.
C_FUNC_DEF double v_edge_length(int num_nodes, double coordinates[][3])
Calculates edge length.
double skew
Definition: verdict.h:130
double equivolume_skew
Definition: verdict.h:338
double shape_and_size
Definition: verdict.h:252
C_FUNC_DEF double v_quad_shear_and_size(int num_nodes, double coordinates[][3])
Calculates quad shear-size metric.
double oddy
Definition: verdict.h:238
double scaled_jacobian
Definition: verdict.h:291
double volume
Definition: verdict.h:277
double med_aspect_frobenius
Definition: verdict.h:144
C_FUNC_DEF void v_set_tri_size(double size)
Sets average size (area) of tri, needed for v_tri_relative_size(...)
double
double distortion
Definition: verdict.h:287
C_FUNC_DEF double v_tet_aspect_beta(int num_nodes, double coordinates[][3])
Calculates the radius ratio metric of a positively oriented tet.
#define C_FUNC_DEF
Definition: verdict.h:99
C_FUNC_DEF double v_tet_scaled_jacobian(int num_nodes, double coordinates[][3])
Calculates tet scaled jacobian.
C_FUNC_DEF double v_wedge_max_stretch(int num_nodes, double coordinates[][3])
Calculates the wedge stretch.
double shape_and_size
Definition: verdict.h:160
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.
double jacobian
Definition: verdict.h:326
C_FUNC_DEF double v_hex_stretch(int num_nodes, double coordinates[][3])
Calculates hex stretch metric.
C_FUNC_DEF void v_set_hex_size(double size)
Sets average size (volume) of hex, needed for v_hex_relative_size(...)
double stretch
Definition: verdict.h:232
C_FUNC_DEF double v_tet_radius_ratio(int num_nodes, double coordinates[][3])
Calculates tet radius ratio metric.
C_FUNC_DEF double v_hex_edge_ratio(int num_nodes, double coordinates[][3])
Calculates hex edge ratio metric.
C_FUNC_DEF double v_hex_shear(int num_nodes, double coordinates[][3])
Calculates hex shear metric.
C_FUNC_DEF double v_hex_skew(int num_nodes, double coordinates[][3])
Calculates hex skew metric.
double condition
Definition: verdict.h:240
double relative_size_squared
Definition: verdict.h:372
C_FUNC_DEF double v_quad_area(int num_nodes, double coordinates[][3])
Calculates quad area.
double edge_ratio
Definition: verdict.h:279
C_FUNC_DEF double v_hex_distortion(int num_nodes, double coordinates[][3])
Calculates hex distortion metric.
double med_aspect_frobenius
Definition: verdict.h:220
C_FUNC_DEF double v_hex_volume(int num_nodes, double coordinates[][3])
Calculates hex volume.
C_FUNC_DEF double v_tet_aspect_gamma(int num_nodes, double coordinates[][3])
Calculates tet aspect gamma metric.
double shape
Definition: verdict.h:248
double condition
Definition: verdict.h:366
double volume
Definition: verdict.h:185
double taper
Definition: verdict.h:132
C_FUNC_DEF double v_quad_maximum_angle(int num_nodes, double coordinates[][3])
Calculates quad's largest angle.
double shear
Definition: verdict.h:154
double maximum_angle
Definition: verdict.h:236
C_FUNC_DEF double v_quad_minimum_angle(int num_nodes, double coordinates[][3])
Calculates quad's smallest angle.
A struct to hold values computed by v_tet_quality.
Definition: verdict.h:303
C_FUNC_DEF double v_quad_radius_ratio(int num_nodes, double coordinates[][3])
Calculates quad radius ratio.
A struct to hold values computed by v_wedge_quality.
Definition: verdict.h:274
C_FUNC_DEF double v_tri_minimum_angle(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double max_edge_ratio
Definition: verdict.h:214
double distortion
Definition: verdict.h:256
double relative_size_squared
Definition: verdict.h:158
C_FUNC_DEF double v_quad_aspect_ratio(int num_nodes, double coordinates[][3])
Calculates quad aspect ratio.
A struct to hold values computed by v_hex_quality.
Definition: verdict.h:123
C_FUNC_DEF double v_quad_taper(int num_nodes, double coordinates[][3])
Calculates quad taper metric.
double aspect_ratio
Definition: verdict.h:216
C_FUNC_DEF double v_quad_relative_size_squared(int num_nodes, double coordinates[][3])
Calculates quad relative size metric.
C_FUNC_DEF double v_quad_shear(int num_nodes, double coordinates[][3])
Calculates quad shear metric.
double edge_ratio
Definition: verdict.h:352
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.
double distortion
Definition: verdict.h:164
C_FUNC_DEF double v_tet_shape(int num_nodes, double coordinates[][3])
Calculates tet shape metric.
C_FUNC_DEF double v_tet_edge_ratio(int num_nodes, double coordinates[][3])
Calculates tet edge ratio metric.
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.
double aspect_ratio
Definition: verdict.h:354
C_FUNC_DEF void v_set_quad_size(double size)
Sets average size (area) of quad, needed for v_quad_relative_size(...)
double scaled_jacobian
Definition: verdict.h:368
C_FUNC_DEF double v_hex_shape_and_size(int num_nodes, double coordinates[][3])
Calculates hex shape-size metric.
C_FUNC_DEF double v_hex_relative_size_squared(int num_nodes, double coordinates[][3])
Calculates hex relative size metric.
double minimum_angle
Definition: verdict.h:318
C_FUNC_DEF double v_tri_shape_and_size(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double(* VerdictFunction)(int, double[][3])
Definition: verdict.h:111
double condition
Definition: verdict.h:295
double edge_ratio
Definition: verdict.h:212
double shape
Definition: verdict.h:156
double mean_aspect_frobenius
Definition: verdict.h:283
C_FUNC_DEF double v_quad_max_edge_ratio(int num_nodes, double coordinates[][3])
Calculates quad maximum of edge ratio.
C_FUNC_DEF double v_tet_relative_size_squared(int num_nodes, double coordinates[][3])
Calculates tet relative size metric.
point
C_FUNC_DEF double v_tet_distortion(int num_nodes, double coordinates[][3])
Calculates tet distortion metric.
double relative_size_squared
Definition: verdict.h:250
double max_edge_ratio
Definition: verdict.h:128
C_FUNC_DEF double v_hex_med_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates hex condition metric.
double shear_and_size
Definition: verdict.h:254
C_FUNC_DEF double v_tet_aspect_ratio(int num_nodes, double coordinates[][3])
Calculates tet aspect ratio metric.
C_FUNC_DEF double v_tri_area(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double radius_ratio
Definition: verdict.h:356
C_FUNC_DEF double v_wedge_distortion(int num_nodes, double coordinates[][3])
Calculates wedge distortion metric.
C_FUNC_DEF double v_tri_radius_ratio(int num_nodes, double coordinates[][3])
Calculates triangle metric.
C_FUNC_DEF double v_hex_shear_and_size(int num_nodes, double coordinates[][3])
Calculates hex shear-size metric.
double jacobian
Definition: verdict.h:285
double warpage
Definition: verdict.h:228
C_FUNC_DEF double v_hex_max_edge_ratio(int num_nodes, double coordinates[][3])
Calculates hex maximum of edge ratio.
C_FUNC_DEF double v_tet_jacobian(int num_nodes, double coordinates[][3])
Calculates tet jacobian.
double oddy
Definition: verdict.h:142
C_FUNC_DEF double v_quad_med_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates quad average Frobenius aspect.
double collapse_ratio
Definition: verdict.h:320
C_FUNC_DEF double v_tri_scaled_jacobian(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double aspect_gamma
Definition: verdict.h:314
double diagonal
Definition: verdict.h:138
double shape_and_size
Definition: verdict.h:374
double radius_ratio
Definition: verdict.h:218
double distortion
Definition: verdict.h:376
C_FUNC_DEF double v_hex_shape(int num_nodes, double coordinates[][3])
Calculates hex shape metric.
C_FUNC_DEF double v_hex_taper(int num_nodes, double coordinates[][3])
Calculates hex taper metric.
double distortion
Definition: verdict.h:336
C_FUNC_DEF double v_tri_condition(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double skew
Definition: verdict.h:224
C_FUNC_DEF double v_tri_relative_size_squared(int num_nodes, double coordinates[][3])
Calculates triangle metric.
double shape
Definition: verdict.h:370
C_FUNC_DEF double v_quad_max_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates quad maximum Frobenius aspect.
C_FUNC_DEF double v_tet_collapse_ratio(int num_nodes, double coordinates[][3])
Calculates tet collapse ratio metric.
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.
double condition
Definition: verdict.h:324
double shape
Definition: verdict.h:330
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.
C_FUNC_DEF double v_wedge_max_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates wedge max aspect forbenius.
int(* ComputeNormal)(double point[3], double normal[3])
Definition: verdict.h:112
C_FUNC_DEF double v_tet_aspect_frobenius(int num_nodes, double coordinates[][3])
Calculates tet aspect frobenius metric.