Mesh Oriented datABase  (version 5.6.0)
An array-based unstructured mesh library
verdict.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Module: $RCSfile: verdict.h,v $
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 /*! \file verdict.h
16  \brief Header file for verdict library that calculates metrics for finite elements.
17  Also see: \ref index "Main Page"
18  *
19  * verdict.h is the header file for applications/libraries to include
20  * to compute quality metrics.
21  *
22  * This file is part of VERDICT
23  *
24 */
25 
26 /* note: changes to this file must be propagated to verdict.h.in so that
27  Verdict can be compiled using CMake and autoconf. */
28 
29 #ifndef VERDICT_INC_LIB
30 #define VERDICT_INC_LIB
31 
32 #ifndef VERDICT_VERSION
33 #define VERDICT_VERSION 120
34 #endif
35 
36 #define VERDICT_DBL_MIN 1.0E-30
37 #define VERDICT_DBL_MAX 1.0E+30
38 #define VERDICT_PI 3.1415926535897932384626
39 
40 #ifdef __cplusplus
41 #if defined( WIN32 ) && defined( VERDICT_SHARED_LIB )
42 #ifdef verdict_EXPORTS
43 #define C_FUNC_DEF extern "C" __declspec( dllexport )
44 #else
45 #define C_FUNC_DEF extern "C" __declspec( dllimport )
46 #endif
47 #else
48 #define C_FUNC_DEF extern "C"
49 #endif
50 #else
51 #if defined( WIN32 ) && defined( VERDICT_SHARED_LIB )
52 #ifdef verdict_EXPORTS
53 #define C_FUNC_DEF __declspec( dllexport )
54 #else
55 #define C_FUNC_DEF __declspec( dllimport )
56 #endif
57 #else
58 #define C_FUNC_DEF
59 #endif
60 #endif
61 
62 /* typedef for the user if they want to use
63  * function pointers */
64 #ifdef __cplusplus
65 extern "C" {
66 #endif
67 typedef double ( *VerdictFunction )( int, double[][3] );
68 typedef int ( *ComputeNormal )( double point[3], double normal[3] );
69 #ifdef __cplusplus
70 }
71 #endif
72 
73 /** HexMetricVals is a <em>struct</em> used to return calculated metrics
74  when calling the function <em>v_hex_quality(...)</em>
75 
76  HexMetricVals is used just as QuadMetricVals.
77  For an example, see Detailed Description in QuadMetricVals
78 */
80 {
81  /** \sa v_hex_edge_ratio */
82  double edge_ratio;
83  /** \sa v_hex_max_edge_ratio */
85  /** \sa v_hex_skew */
86  double skew;
87  /** \sa v_hex_taper */
88  double taper;
89  /** \sa v_hex_volume */
90  double volume;
91  /** \sa v_hex_stretch */
92  double stretch;
93  /** \sa v_hex_diagonal */
94  double diagonal;
95  /** \sa v_hex_dimension */
96  double dimension;
97  /** \sa v_hex_oddy */
98  double oddy;
99  /** \sa v_hex_med_aspect_frobenius */
101  /** \sa v_hex_condition */
102  double condition;
103  /** \sa v_hex_jacobian */
104  double jacobian;
105  /** \sa v_hex_scaled_jacobian */
107  /** \sa v_hex_shear */
108  double shear;
109  /** \sa v_hex_shape */
110  double shape;
111  /** \sa v_hex_relative_size */
113  /** \sa v_hex_shape_and_size */
115  /** \sa v_hex_shear_and_size */
117  /** \sa v_hex_distortion */
118  double distortion;
119 };
120 
121 /** EdgeMetricVals is a <em>struct</em> used to return calculated metrics
122  when calling the function <em>v_edge_quality(...)</em>
123 
124  EdgeMetricVals is used just as QuadMetricVals.
125  For an example, see Detailed Description in QuadMetricVals
126 */
128 {
129  double length;
130 };
131 
132 /** KnifeMetricVals is a <em>struct</em> used to return calculated metrics
133  when calling the function <em>v_knife_quality(...)</em>
134 
135  KnifeMetricVals is used just as QuadMetricVals.
136  For an example, see Detailed Description in QuadMetricVals
137 */
139 {
140  double volume;
141 };
142 
143 /** QuadMetricVals is a <em>struct</em> used to return calculated metrics
144  when calling the function <em>v_quad_quality(...)</em>
145 
146  The following is an example of how this struct is used with Verdict.
147 
148  Example: \code
149  QuadMetricVals quad_metrics = {0};
150  unsigned long metrics_flag = 0;
151  metrics_flag += V_QUAD_SHAPE;
152  metrics_flag += V_QUAD_DISTORTION;
153  metrics_flag += V_QUAD_AREA;
154  double quad_nodes[4][3];
155  get_quad_nodes( quad_nodes ); //some user-defined function to load
156  //xyz coordinate info. into array
157  v_quad_quality( 4, quad_nodes, metrics_flag, quad_metrics );
158  double my_shape = quad_metrics.shape;
159  double my_distortion = quad_metrics.distortion;
160  double my_area = quad_metrics.area; \endcode
161 
162  */
163 
165 {
166  /** \sa v_quad_edge_ratio function */
167  double edge_ratio;
168  /** \sa v_quad_max_edge_ratio function */
170  /** \sa v_quad_aspect_ratio function */
171  double aspect_ratio;
172  /** \sa v_quad_radius_ratio function */
173  double radius_ratio;
174  /** \sa v_quad_med_aspect_frobenius function */
176  /** \sa v_quad_max_aspect_frobenius function */
178  /** \sa v_quad_skew function */
179  double skew;
180  /** \sa v_quad_taper function */
181  double taper;
182  /** \sa v_quad_warpage function */
183  double warpage;
184  /** \sa v_quad_area function */
185  double area;
186  /** \sa v_quad_stretch function */
187  double stretch;
188  /** \sa v_quad_smallest_angle function */
190  /** \sa v_quad_largest_angle function */
192  /** \sa v_quad_oddy function */
193  double oddy;
194  /** \sa v_quad_condition function */
195  double condition;
196  /** \sa v_quad_jacobian function */
197  double jacobian;
198  /** \sa v_quad_scaled_jacobian function */
200  /** \sa v_quad_shear function */
201  double shear;
202  /** \sa v_quad_shape function */
203  double shape;
204  /** \sa v_quad_relative_size_squared function */
206  /** \sa v_quad_shape_and_size function */
208  /** \sa v_quad_shear_and_size function */
210  /** \sa v_quad_distortion function */
211  double distortion;
212 };
213 
214 /** PyramidMetricVals is a <em>struct</em> used to return calculated metrics
215  when calling the function <em>v_pyramid_quality(...)</em>
216 
217  PyramidMetricVals is used just as QuadMetricVals.
218  For an example, see Detailed Description in QuadMetricVals
219 */
221 {
222  double volume;
223 };
224 
225 /** WedgeMetricVals is a <em>struct</em> used to return calculated metrics
226  when calling the function <em>v_wedge_quality(...)</em>
227 
228  WedgeMetricVals is used just as QuadMetricVals.
229  For an example, see Detailed Description in QuadMetricVals
230 */
232 {
233  double volume;
234 };
235 
236 /** TetMetricVals is a <em>struct</em> used to return calculated metrics
237  when calling the function <em>v_tet_quality(...)</em>
238 
239  TetMetricVals is used just as QuadMetricVals.
240  For an example, see Detailed Description in QuadMetricVals
241 */
243 {
244  /** \sa v_tet_edge_ratio*/
245  double edge_ratio;
246  /** \sa v_tet_radius_ratio*/
247  double radius_ratio;
248  /** \sa v_tet_aspect_beta*/
249  double aspect_beta;
250  /** \sa v_tet_aspect_ratio */
251  double aspect_ratio;
252  /** \sa v_tet_aspect_gamma */
253  double aspect_gamma;
254  /** \sa v_tet_aspect_frobenius */
256  /** \sa v_tet_minimum_angle */
258  /** \sa v_tet_collapse_ratio*/
260  /** \sa v_tet_volume */
261  double volume;
262  /** \sa v_tet_condition */
263  double condition;
264  /** \sa v_tet_jacobian */
265  double jacobian;
266  /** \sa v_tet_scaled_jacobian */
268  /** \sa v_tet_shape */
269  double shape;
270  /** \sa v_tet_relative_size */
272  /** \sa v_tet_shape_and_size*/
274  /** \sa v_tet_distortion */
275  double distortion;
276 };
277 
278 /** TriMetricVals is a <em>struct</em> used to return calculated metrics
279  when calling the function <em>v_tri_quality(...)</em>
280 
281  TriMetricVals is used just as QuadMetricVals.
282  For an example, see Detailed Description in QuadMetricVals
283 */
285 {
286  /** \sa v_tri_edge_ratio */
287  double edge_ratio;
288  /** \sa v_tri_aspect_ratio */
289  double aspect_ratio;
290  /** \sa v_tri_radius_ratio */
291  double radius_ratio;
292  /** \sa v_tri_aspect_frobenius */
294  /** \sa v_tri_area*/
295  double area;
296  /** \sa v_tri_minimum_angle*/
298  /** \sa v_tri_maximum_angle */
300  /** \sa v_tri_condition */
301  double condition;
302  /** \sa v_tri_scaled_jacobian */
304  /** \sa v_tri_shape */
305  double shape;
306  /** \sa v_tri_relative_size_squared */
308  /** \sa v_tri_shape_and_size */
310  /** \sa v_tri_distortion */
311  double distortion;
312 };
313 
314 /* definition of bit fields to determine which metrics to calculate */
315 
316 //! \name Hex bit fields
317 //!
318 //@{
319 
320 #define V_HEX_MAX_EDGE_RATIO 1 /*!< \hideinitializer */
321 #define V_HEX_SKEW 2 /*!< \hideinitializer */
322 #define V_HEX_TAPER 4 /*!< \hideinitializer */
323 #define V_HEX_VOLUME 8 /*!< \hideinitializer */
324 #define V_HEX_STRETCH 16 /*!< \hideinitializer */
325 #define V_HEX_DIAGONAL 32 /*!< \hideinitializer */
326 #define V_HEX_DIMENSION 64 /*!< \hideinitializer */
327 #define V_HEX_ODDY 128 /*!< \hideinitializer */
328 #define V_HEX_MAX_ASPECT_FROBENIUS 256 /*!< \hideinitializer */
329 #define V_HEX_CONDITION 256 /*!< \hideinitializer */
330 #define V_HEX_JACOBIAN 512 /*!< \hideinitializer */
331 #define V_HEX_SCALED_JACOBIAN 1024 /*!< \hideinitializer */
332 #define V_HEX_SHEAR 2048 /*!< \hideinitializer */
333 #define V_HEX_SHAPE 4096 /*!< \hideinitializer */
334 #define V_HEX_RELATIVE_SIZE_SQUARED 8192 /*!< \hideinitializer */
335 #define V_HEX_SHAPE_AND_SIZE 16384 /*!< \hideinitializer */
336 #define V_HEX_SHEAR_AND_SIZE 32768 /*!< \hideinitializer */
337 #define V_HEX_DISTORTION 65536 /*!< \hideinitializer */
338 #define V_HEX_EDGE_RATIO 131072 /*!< \hideinitializer */
339 #define V_HEX_MED_ASPECT_FROBENIUS 262144 /*!< \hideinitializer */
340 #define V_HEX_ALL 524287 /*!< \hideinitializer */
341 /*!< \hideinitializer */
342 #define V_HEX_TRADITIONAL \
343  ( V_HEX_MAX_EDGE_RATIO + V_HEX_SKEW + V_HEX_TAPER + V_HEX_STRETCH + V_HEX_DIAGONAL + V_HEX_ODDY + \
344  V_HEX_CONDITION + V_HEX_JACOBIAN + V_HEX_SCALED_JACOBIAN + V_HEX_DIMENSION )
345 
346 /*!< \hideinitializer */
347 #define V_HEX_DIAGNOSTIC V_HEX_VOLUME
348 /*!< \hideinitializer */
349 #define V_HEX_ALGEBRAIC \
350  ( V_HEX_SHAPE + V_HEX_SHEAR + V_HEX_RELATIVE_SIZE_SQUARED + V_HEX_SHAPE_AND_SIZE + V_HEX_SHEAR_AND_SIZE )
351 /*!< \hideinitializer */
352 #define V_HEX_ROBINSON ( V_HEX_SKEW + V_HEX_TAPER )
353 //@}
354 
355 //! \name Tet bit fields
356 //!
357 //@{
358 #define V_TET_RADIUS_RATIO 1 /*!< \hideinitializer */
359 #define V_TET_ASPECT_BETA 2 /*!< \hideinitializer */
360 #define V_TET_ASPECT_GAMMA 4 /*!< \hideinitializer */
361 #define V_TET_VOLUME 8 /*!< \hideinitializer */
362 #define V_TET_CONDITION 16 /*!< \hideinitializer */
363 #define V_TET_JACOBIAN 32 /*!< \hideinitializer */
364 #define V_TET_SCALED_JACOBIAN 64 /*!< \hideinitializer */
365 #define V_TET_SHAPE 128 /*!< \hideinitializer */
366 #define V_TET_RELATIVE_SIZE_SQUARED 256 /*!< \hideinitializer */
367 #define V_TET_SHAPE_AND_SIZE 512 /*!< \hideinitializer */
368 #define V_TET_DISTORTION 1024 /*!< \hideinitializer */
369 #define V_TET_EDGE_RATIO 2048 /*!< \hideinitializer */
370 #define V_TET_ASPECT_RATIO 4096 /*!< \hideinitializer */
371 #define V_TET_ASPECT_FROBENIUS 8192 /*!< \hideinitializer */
372 #define V_TET_MINIMUM_ANGLE 16384 /*!< \hideinitializer */
373 #define V_TET_COLLAPSE_RATIO 32768 /*!< \hideinitializer */
374 #define V_TET_ALL 65535 /*!< \hideinitializer */
375 /*!< \hideinitializer */
376 #define V_TET_TRADITIONAL \
377  ( V_TET_ASPECT_BETA + V_TET_ASPECT_GAMMA + V_TET_CONDITION + V_TET_JACOBIAN + V_TET_SCALED_JACOBIAN )
378 /*!< \hideinitializer */
379 #define V_TET_DIAGNOSTIC V_TET_VOLUME
380 /*!< \hideinitializer */
381 #define V_TET_ALGEBRAIC ( V_TET_SHAPE + V_TET_RELATIVE_SIZE_SQUARED + V_TET_SHAPE_AND_SIZE )
382 //@}
383 
384 //! \name Pyramid bit fields
385 //!
386 //@{
387 #define V_PYRAMID_VOLUME 1 /*!< \hideinitializer */
388 //@}
389 
390 //! \name Wedge bit fields
391 //!
392 //@{
393 #define V_WEDGE_VOLUME 1 /*!< \hideinitializer */
394 //@}
395 
396 //! \name Knife bit fields
397 //!
398 //@{
399 #define V_KNIFE_VOLUME 1 /*!< \hideinitializer */
400 //@}
401 
402 //! \name Quad bit fields
403 //!
404 //@{
405 #define V_QUAD_MAX_EDGE_RATIO 1 /*!< \hideinitializer */
406 #define V_QUAD_SKEW 2 /*!< \hideinitializer */
407 #define V_QUAD_TAPER 4 /*!< \hideinitializer */
408 #define V_QUAD_WARPAGE 8 /*!< \hideinitializer */
409 #define V_QUAD_AREA 16 /*!< \hideinitializer */
410 #define V_QUAD_STRETCH 32 /*!< \hideinitializer */
411 #define V_QUAD_MINIMUM_ANGLE 64 /*!< \hideinitializer */
412 #define V_QUAD_MAXIMUM_ANGLE 128 /*!< \hideinitializer */
413 #define V_QUAD_ODDY 256 /*!< \hideinitializer */
414 #define V_QUAD_CONDITION 512 /*!< \hideinitializer */
415 #define V_QUAD_JACOBIAN 1024 /*!< \hideinitializer */
416 #define V_QUAD_SCALED_JACOBIAN 2048 /*!< \hideinitializer */
417 #define V_QUAD_SHEAR 4096 /*!< \hideinitializer */
418 #define V_QUAD_SHAPE 8192 /*!< \hideinitializer */
419 #define V_QUAD_RELATIVE_SIZE_SQUARED 16384 /*!< \hideinitializer */
420 #define V_QUAD_SHAPE_AND_SIZE 32768 /*!< \hideinitializer */
421 #define V_QUAD_SHEAR_AND_SIZE 65536 /*!< \hideinitializer */
422 #define V_QUAD_DISTORTION 131072 /*!< \hideinitializer */
423 #define V_QUAD_EDGE_RATIO 262144 /*!< \hideinitializer */
424 #define V_QUAD_ASPECT_RATIO 524288 /*!< \hideinitializer */
425 #define V_QUAD_RADIUS_RATIO 1048576 /*!< \hideinitializer */
426 #define V_QUAD_MED_ASPECT_FROBENIUS 2097152 /*!< \hideinitializer */
427 #define V_QUAD_MAX_ASPECT_FROBENIUS 4194304 /*!< \hideinitializer */
428 #define V_QUAD_ALL 8388607 /*!< \hideinitializer */
429 /*!< \hideinitializer */
430 #define V_QUAD_TRADITIONAL \
431  ( V_QUAD_MAX_EDGE_RATIO + V_QUAD_SKEW + V_QUAD_TAPER + V_QUAD_WARPAGE + V_QUAD_STRETCH + V_QUAD_MINIMUM_ANGLE + \
432  V_QUAD_MAXIMUM_ANGLE + V_QUAD_ODDY + V_QUAD_CONDITION + V_QUAD_JACOBIAN + V_QUAD_SCALED_JACOBIAN )
433 /*!< \hideinitializer */
434 #define V_QUAD_DIAGNOSTIC V_QUAD_AREA
435 /*!< \hideinitializer */
436 #define V_QUAD_ALGEBRAIC ( V_QUAD_SHEAR + V_QUAD_SHAPE + V_QUAD_RELATIVE_SIZE_SQUARED + V_QUAD_SHAPE_AND_SIZE )
437 /*!< \hideinitializer */
438 #define V_QUAD_ROBINSON ( V_QUAD_MAX_EDGE_RATIO + V_QUAD_SKEW + V_QUAD_TAPER )
439 //@}
440 
441 //! \name Tri bit fields
442 //!
443 //@{
444 #define V_TRI_ASPECT_FROBENIUS 1 /*!< \hideinitializer */
445 #define V_TRI_AREA 2 /*!< \hideinitializer */
446 #define V_TRI_MINIMUM_ANGLE 4 /*!< \hideinitializer */
447 #define V_TRI_MAXIMUM_ANGLE 8 /*!< \hideinitializer */
448 #define V_TRI_CONDITION 16 /*!< \hideinitializer */
449 #define V_TRI_SCALED_JACOBIAN 32 /*!< \hideinitializer */
450 #define V_TRI_SHAPE 64 /*!< \hideinitializer */
451 #define V_TRI_RELATIVE_SIZE_SQUARED 128 /*!< \hideinitializer */
452 #define V_TRI_SHAPE_AND_SIZE 256 /*!< \hideinitializer */
453 #define V_TRI_DISTORTION 512 /*!< \hideinitializer */
454 #define V_TRI_RADIUS_RATIO 1024 /*!< \hideinitializer */
455 #define V_TRI_EDGE_RATIO 2048 /*!< \hideinitializer */
456 #define V_TRI_ALL 4095 /*!< \hideinitializer */
457 /*!< \hideinitializer */
458 #define V_TRI_TRADITIONAL \
459  ( V_TRI_ASPECT_FROBENIUS + V_TRI_MINIMUM_ANGLE + V_TRI_MAXIMUM_ANGLE + V_TRI_CONDITION + V_TRI_SCALED_JACOBIAN )
460 /*!< \hideinitializer */
461 #define V_TRI_DIAGNOSTIC V_TRI_AREA
462 /*!< \hideinitializer */
463 #define V_TRI_ALGEBRAIC ( V_TRI_SHAPE + V_TRI_SHAPE_AND_SIZE + V_TRI_RELATIVE_SIZE_SQUARED )
464 
465 #define V_EDGE_LENGTH 1 /*!< \hideinitializer */
466  //@}
467 
468 /*! \ref verdict "Verdict Mesh Quality Metrics" are available in the MOAB User's Guide.
469  This library provides comprehensive mesh quality assessment capabilities for various element types.
470 */
471 
472 //! Calculates quality metrics for hexahedral elements.
473 C_FUNC_DEF void v_hex_quality( int num_nodes,
474  double coordinates[][3],
475  unsigned int metrics_request_flag,
476  struct HexMetricVals* metric_vals );
477 
478 //! Calculates quality metrics for tetrahedral elements.
479 C_FUNC_DEF void v_tet_quality( int num_nodes,
480  double coordinates[][3],
481  unsigned int metrics_request_flag,
482  struct TetMetricVals* metric_vals );
483 
484 //! Calculates quality metrics for pyramid elements.
485 C_FUNC_DEF void v_pyramid_quality( int num_nodes,
486  double coordinates[][3],
487  unsigned int metrics_request_flag,
488  struct PyramidMetricVals* metric_vals );
489 
490 //! Calculates quality metrics for wedge elements.
491 C_FUNC_DEF void v_wedge_quality( int num_nodes,
492  double coordinates[][3],
493  unsigned int metrics_request_flag,
494  struct WedgeMetricVals* metric_vals );
495 
496 //! Calculates quality metrics for knife elements.
497 C_FUNC_DEF void v_knife_quality( int num_nodes,
498  double coordinates[][3],
499  unsigned int metrics_request_flag,
500  struct KnifeMetricVals* metric_vals );
501 
502 //! Calculates quality metrics for quadrilateral elements.
503 C_FUNC_DEF void v_quad_quality( int num_nodes,
504  double coordinates[][3],
505  unsigned int metrics_request_flag,
506  struct QuadMetricVals* metric_vals );
507 
508 //! Calculates quality metrics for triangle elements.
509 C_FUNC_DEF void v_tri_quality( int num_nodes,
510  double coordinates[][3],
511  unsigned int metrics_request_flag,
512  struct TriMetricVals* metric_vals );
513 
514 //! Calculates quality metrics for edge elements.
515 C_FUNC_DEF void v_edge_quality( int num_nodes,
516  double coordinates[][3],
517  unsigned int metrics_request_flag,
518  struct EdgeMetricVals* metric_vals );
519 
520 /* individual quality functions for hex elements */
521 
522 //! Sets average size (volume) of hex, needed for v_hex_relative_size(...)
523 C_FUNC_DEF void v_set_hex_size( double size );
524 
525 //! Calculates hex edge ratio metric.
526 /** Hmax / Hmin where Hmax and Hmin are respectively the maximum and the
527  minimum edge lengths */
528 C_FUNC_DEF double v_hex_edge_ratio( int num_nodes, double coordinates[][3] );
529 
530 //! Calculates hex maximum of edge ratio
531 /**Maximum edge length ratio at hex center.
532  Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
533  Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
534 C_FUNC_DEF double v_hex_max_edge_ratio( int num_nodes, double coordinates[][3] );
535 
536 //! Calculates hex skew metric.
537 /** Maximum |cos A| where A is the angle between edges at hex center.
538  Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
539  Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
540 C_FUNC_DEF double v_hex_skew( int num_nodes, double coordinates[][3] );
541 
542 //! Calculates hex taper metric
543 /** Maximum ratio of lengths derived from opposite edges.
544  Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
545  Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
546 C_FUNC_DEF double v_hex_taper( int num_nodes, double coordinates[][3] );
547 
548 //! Calculates hex volume
549 /** Jacobian at hex center.
550  Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
551  Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
552 C_FUNC_DEF double v_hex_volume( int num_nodes, double coordinates[][3] );
553 
554 //! Calculates hex stretch metric
555 /** Sqrt(3) * minimum edge length / maximum diagonal length.
556  Reference --- FIMESH code */
557 C_FUNC_DEF double v_hex_stretch( int num_nodes, double coordinates[][3] );
558 
559 //! Calculates hex diagonal metric
560 /** Minimum diagonal length / maximum diagonal length.
561  Reference --- Unknown */
562 C_FUNC_DEF double v_hex_diagonal( int num_nodes, double coordinates[][3] );
563 
564 //! Calculates hex dimension metric
565 /** Pronto-specific characteristic length for stable time step calculation.
566  Char_length = Volume / 2 grad Volume.
567  Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
568  Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
569 C_FUNC_DEF double v_hex_dimension( int num_nodes, double coordinates[][3] );
570 
571 //! Calculates hex oddy metric
572 C_FUNC_DEF double v_hex_oddy( int num_nodes, double coordinates[][3] );
573 
574 //! Calculates hex condition metric
575 /** Average Frobenius condition number of the Jacobian matrix at 8 corners. */
576 C_FUNC_DEF double v_hex_med_aspect_frobenius( int num_nodes, double coordinates[][3] );
577 
578 //! Calculates hex condition metric
579 /** Maximum Frobenius condition number of the Jacobian matrix at 8 corners.
580  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
581  Optimization of the Jacobian Matrix Norm and Associated Quantities,
582  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
583 C_FUNC_DEF double v_hex_max_aspect_frobenius( int num_nodes, double coordinates[][3] );
584 C_FUNC_DEF double v_hex_condition( int num_nodes, double coordinates[][3] );
585 
586 //! Calculates hex jacobian metric
587 /** Minimum pointwise volume of local map at 8 corners & center of hex.
588  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
589  Optimization of the Jacobian Matrix Norm and Associated Quantities,
590  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
591 C_FUNC_DEF double v_hex_jacobian( int num_nodes, double coordinates[][3] );
592 
593 //! Calculates hex scaled jacobian metric
594 /** Minimum Jacobian divided by the lengths of the 3 edge vectors.
595  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
596  Optimization of the Jacobian Matrix Norm and Associated Quantities,
597  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
598 C_FUNC_DEF double v_hex_scaled_jacobian( int num_nodes, double coordinates[][3] );
599 
600 //! Calculates hex shear metric
601 /** 3/Mean Ratio of Jacobian Skew matrix.
602  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
603  Unstructured Initial Meshes, submitted for publication. */
604 C_FUNC_DEF double v_hex_shear( int num_nodes, double coordinates[][3] );
605 
606 //! Calculates hex shape metric.
607 /** 3/Mean Ratio of weighted Jacobian matrix.
608  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
609  Unstructured Initial Meshes, submitted for publication. */
610 C_FUNC_DEF double v_hex_shape( int num_nodes, double coordinates[][3] );
611 
612 //! Calculates hex relative size metric.
613 /** 3/Mean Ratio of weighted Jacobian matrix.
614  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
615  Unstructured Initial Meshes, submitted for publication. */
616 C_FUNC_DEF double v_hex_relative_size_squared( int num_nodes, double coordinates[][3] );
617 
618 //! Calculates hex shape-size metric.
619 /** Product of Shape and Relative Size.
620  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
621  Unstructured Initial Meshes, submitted for publication. */
622 C_FUNC_DEF double v_hex_shape_and_size( int num_nodes, double coordinates[][3] );
623 
624 //! Calculates hex shear-size metric
625 /** Product of Shear and Relative Size.
626  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
627  Unstructured Initial Meshes, submitted for publication. */
628 C_FUNC_DEF double v_hex_shear_and_size( int num_nodes, double coordinates[][3] );
629 
630 //! Calculates hex distortion metric
631 /** {min(|J|)/actual volume}*parent volume, parent volume = 8 for hex.
632  Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
633 C_FUNC_DEF double v_hex_distortion( int num_nodes, double coordinates[][3] );
634 
635 /* individual quality functions for tet elements */
636 
637 //! Sets average size (volume) of tet, needed for v_tet_relative_size(...)
638 C_FUNC_DEF void v_set_tet_size( double size );
639 
640 //! Calculates tet edge ratio metric.
641 /** Hmax / Hmin where Hmax and Hmin are respectively the maximum and the
642  minimum edge lengths */
643 C_FUNC_DEF double v_tet_edge_ratio( int num_nodes, double coordinates[][3] );
644 
645 //! Calculates tet radius ratio metric.
646 /** CR / (3.0 * IR) where CR = circumsphere radius, IR = inscribed sphere radius.
647  Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
648  quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
649 C_FUNC_DEF double v_tet_radius_ratio( int num_nodes, double coordinates[][3] );
650 
651 //! Calculates the radius ratio metric of a positively oriented tet.
652 /** CR / (3.0 * IR) where CR = circumsphere radius, IR = inscribed sphere radius
653  if the element is positively-oriented.
654  Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
655  quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
656 C_FUNC_DEF double v_tet_aspect_beta( int num_nodes, double coordinates[][3] );
657 
658 //! Calculates tet aspect ratio metric.
659 /** Hmax / (2 sqrt(6) r) where Hmax and r respectively denote the greatest edge
660  length and the inradius of the tetrahedron
661  Reference --- P. Frey and P.-L. George, Meshing, Hermes (2000). */
662 C_FUNC_DEF double v_tet_aspect_ratio( int num_nodes, double coordinates[][3] );
663 
664 //! Calculates tet aspect gamma metric.
665 /** Srms**3 / (8.479670*V) where Srms = sqrt(Sum(Si**2)/6), Si = edge length.
666  Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
667  quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
668 C_FUNC_DEF double v_tet_aspect_gamma( int num_nodes, double coordinates[][3] );
669 
670 //! Calculates tet aspect frobenius metric.
671 /** Frobenius condition number when the reference element is regular
672  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
673  Optimization of the Jacobian Matrix Norm and Associated Quantities,
674  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
675 C_FUNC_DEF double v_tet_aspect_frobenius( int num_nodes, double coordinates[][3] );
676 
677 //! Calculates tet minimum dihedral angle.
678 /** Minimum (nonoriented) dihedral angle of a tetrahedron, expressed in degrees. */
679 C_FUNC_DEF double v_tet_minimum_angle( int num_nodes, double coordinates[][3] );
680 
681 //! Calculates tet collapse ratio metric.
682 /** Collapse ratio */
683 C_FUNC_DEF double v_tet_collapse_ratio( int num_nodes, double coordinates[][3] );
684 
685 //! Calculates tet volume.
686 /** (1/6) * Jacobian at corner node.
687  Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
688  quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
689 C_FUNC_DEF double v_tet_volume( int num_nodes, double coordinates[][3] );
690 
691 //! Calculates tet condition metric.
692 /** Condition number of the Jacobian matrix at any corner.
693  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
694  Optimization of the Jacobian Matrix Norm and Associated Quantities,
695  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
696 C_FUNC_DEF double v_tet_condition( int num_nodes, double coordinates[][3] );
697 
698 //! Calculates tet jacobian.
699 /** Minimum pointwise volume at any corner.
700  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
701  Optimization of the Jacobian Matrix Norm and Associated Quantities,
702  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
703 C_FUNC_DEF double v_tet_jacobian( int num_nodes, double coordinates[][3] );
704 
705 //! Calculates tet scaled jacobian.
706 /** Minimum Jacobian divided by the lengths of 3 edge vectors
707  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
708  Optimization of the Jacobian Matrix Norm and Associated Quantities,
709  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
710 C_FUNC_DEF double v_tet_scaled_jacobian( int num_nodes, double coordinates[][3] );
711 
712 //! Calculates tet shape metric.
713 /** 3/Mean Ratio of weighted Jacobian matrix.
714  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
715  Unstructured Initial Meshes, submitted for publication. */
716 C_FUNC_DEF double v_tet_shape( int num_nodes, double coordinates[][3] );
717 
718 //! Calculates tet relative size metric.
719 /** Min( J, 1/J ), where J is determinant of weighted Jacobian matrix.
720  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
721  Unstructured Initial Meshes, submitted for publication. */
722 C_FUNC_DEF double v_tet_relative_size_squared( int num_nodes, double coordinates[][3] );
723 
724 //! Calculates tet shape-size metric.
725 /** Product of Shape and Relative Size.
726  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
727  Unstructured Initial Meshes, submitted for publication. */
728 C_FUNC_DEF double v_tet_shape_and_size( int num_nodes, double coordinates[][3] );
729 
730 //! Calculates tet distortion metric.
731 /** {min(|J|)/actual volume}*parent volume, parent volume = 1/6 for tet.
732  Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
733 C_FUNC_DEF double v_tet_distortion( int num_nodes, double coordinates[][3] );
734 
735 /* individual quality functions for pyramid elements */
736 
737 //! Calculates pyramid volume.
738 C_FUNC_DEF double v_pyramid_volume( int num_nodes, double coordinates[][3] );
739 
740 /* individual quality functions for wedge elements */
741 
742 //! Calculates wedge volume.
743 C_FUNC_DEF double v_wedge_volume( int num_nodes, double coordinates[][3] );
744 
745 /* individual quality functions for knife elements */
746 
747 //! Calculates knife volume.
748 C_FUNC_DEF double v_knife_volume( int num_nodes, double coordinates[][3] );
749 
750 /* individual quality functions for edge elements */
751 
752 //! Calculates edge length.
753 C_FUNC_DEF double v_edge_length( int num_nodes, double coordinates[][3] );
754 
755 /* individual quality functions for quad elements */
756 //! Sets average size (area) of quad, needed for v_quad_relative_size(...)
757 C_FUNC_DEF void v_set_quad_size( double size );
758 
759 //! Calculates quad edge ratio
760 /** edge ratio
761  Reference --- P. P. Pebay, Planar Quadrangle Quality
762  Measures, Eng. Comp., 2004, 20(2):157-173 */
763 C_FUNC_DEF double v_quad_edge_ratio( int num_nodes, double coordinates[][3] );
764 
765 //! Calculates quad maximum of edge ratio.
766 /** Maximum edge length ratio at quad center.
767  Reference --- J. Robinson, CRE Method of element testing and the
768  Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
769 C_FUNC_DEF double v_quad_max_edge_ratio( int num_nodes, double coordinates[][3] );
770 
771 //! Calculates quad aspect ratio
772 /** aspect ratio
773  Reference --- P. P. Pebay, Planar Quadrangle Quality
774  Measures, Eng. Comp., 2004, 20(2):157-173 */
775 C_FUNC_DEF double v_quad_aspect_ratio( int num_nodes, double coordinates[][3] );
776 
777 //! Calculates quad radius ratio
778 /** radius ratio
779  Reference --- P. P. Pebay, Planar Quadrangle Quality
780  Measures, Eng. Comp., 2004, 20(2):157-173 */
781 C_FUNC_DEF double v_quad_radius_ratio( int num_nodes, double coordinates[][3] );
782 
783 //! Calculates quad average Frobenius aspect
784 /** average Frobenius aspect
785  Reference --- P. P. Pebay, Planar Quadrangle Quality
786  Measures, Eng. Comp., 2004, 20(2):157-173 */
787 C_FUNC_DEF double v_quad_med_aspect_frobenius( int num_nodes, double coordinates[][3] );
788 
789 //! Calculates quad maximum Frobenius aspect
790 /** average Frobenius aspect
791  Reference --- P. P. Pebay, Planar Quadrangle Quality
792  Measures, Eng. Comp., 2004, 20(2):157-173 */
793 C_FUNC_DEF double v_quad_max_aspect_frobenius( int num_nodes, double coordinates[][3] );
794 
795 //! Calculates quad skew metric.
796 /** Maximum |cos A| where A is the angle between edges at quad center.
797  Reference --- J. Robinson, CRE Method of element testing and the
798  Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
799 C_FUNC_DEF double v_quad_skew( int num_nodes, double coordinates[][3] );
800 
801 //! Calculates quad taper metric.
802 /** Maximum ratio of lengths derived from opposite edges.
803  Reference --- J. Robinson, CRE Method of element testing and the
804  Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
805 C_FUNC_DEF double v_quad_taper( int num_nodes, double coordinates[][3] );
806 
807 //! Calculates quad warpage metric.
808 /** Cosine of Minimum Dihedral Angle formed by Planes Intersecting in Diagonals.
809  Reference --- J. Robinson, CRE Method of element testing and the
810  Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
811 C_FUNC_DEF double v_quad_warpage( int num_nodes, double coordinates[][3] );
812 
813 //! Calculates quad area.
814 /** Jacobian at quad center.
815  Reference --- J. Robinson, CRE Method of element testing and the
816  Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
817 C_FUNC_DEF double v_quad_area( int num_nodes, double coordinates[][3] );
818 
819 //! Calculates quad strech metric.
820 /** Sqrt(2) * minimum edge length / maximum diagonal length.
821  Reference --- FIMESH code. */
822 C_FUNC_DEF double v_quad_stretch( int num_nodes, double coordinates[][3] );
823 
824 //! Calculates quad's smallest angle.
825 /** Smallest included quad angle (degrees).
826  Reference --- Unknown. */
827 C_FUNC_DEF double v_quad_minimum_angle( int num_nodes, double coordinates[][3] );
828 
829 //! Calculates quad's largest angle.
830 /** Largest included quad angle (degrees).
831  Reference --- Unknown. */
832 C_FUNC_DEF double v_quad_maximum_angle( int num_nodes, double coordinates[][3] );
833 
834 //! Calculates quad oddy metric.
835 C_FUNC_DEF double v_quad_oddy( int num_nodes, double coordinates[][3] );
836 
837 //! Calculates quad condition number metric.
838 /** Maximum condition number of the Jacobian matrix at 4 corners.
839  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
840  Optimization of the Jacobian Matrix Norm and Associated Quantities,
841  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
842 C_FUNC_DEF double v_quad_condition( int num_nodes, double coordinates[][3] );
843 
844 //! Calculates quad jacobian.
845 /** Minimum pointwise volume of local map at 4 corners & center of quad.
846  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
847  Optimization of the Jacobian Matrix Norm and Associated Quantities,
848  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
849 C_FUNC_DEF double v_quad_jacobian( int num_nodes, double coordinates[][3] );
850 
851 //! Calculates quad scaled jacobian.
852 /** Minimum Jacobian divided by the lengths of the 2 edge vectors.
853  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
854  Optimization of the Jacobian Matrix Norm and Associated Quantities,
855  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
856 C_FUNC_DEF double v_quad_scaled_jacobian( int num_nodes, double coordinates[][3] );
857 
858 //! Calculates quad shear metric.
859 /** 2/Condition number of Jacobian Skew matrix.
860  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
861  Unstructured Initial Meshes, submitted for publication. */
862 C_FUNC_DEF double v_quad_shear( int num_nodes, double coordinates[][3] );
863 
864 //! Calculates quad shape metric.
865 /** 2/Condition number of weighted Jacobian matrix.
866  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
867  Unstructured Initial Meshes, submitted for publication. */
868 C_FUNC_DEF double v_quad_shape( int num_nodes, double coordinates[][3] );
869 
870 //! Calculates quad relative size metric.
871 /** Min( J, 1/J ), where J is determinant of weighted Jacobian matrix.
872  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
873  Unstructured Initial Meshes, submitted for publication. */
874 C_FUNC_DEF double v_quad_relative_size_squared( int num_nodes, double coordinates[][3] );
875 
876 //! Calculates quad shape-size metric.
877 /** Product of Shape and Relative Size.
878  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
879  Unstructured Initial Meshes, submitted for publication. */
880 C_FUNC_DEF double v_quad_shape_and_size( int num_nodes, double coordinates[][3] );
881 
882 //! Calculates quad shear-size metric.
883 /** Product of Shear and Relative Size.
884  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
885  Unstructured Initial Meshes, submitted for publication. */
886 C_FUNC_DEF double v_quad_shear_and_size( int num_nodes, double coordinates[][3] );
887 
888 //! Calculates quad distortion metric.
889 /** {min(|J|)/actual area}*parent area, parent area = 4 for quad.
890  Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
891 C_FUNC_DEF double v_quad_distortion( int num_nodes, double coordinates[][3] );
892 
893 /* individual quality functions for tri elements */
894 
895 //! Sets average size (area) of tri, needed for v_tri_relative_size(...)
896 C_FUNC_DEF void v_set_tri_size( double size );
897 
898 //! Sets fuction pointer to calculate tri normal wrt surface
900 
901 //! Calculates tri metric.
902 /** edge ratio
903  Reference --- P. P. Pebay & T. J. Baker, Analysis of Triangle Quality
904  Measures, AMS Math. Comp., 2003, 72(244):1817-1839 */
905 C_FUNC_DEF double v_tri_edge_ratio( int num_nodes, double coordinates[][3] );
906 
907 //! Calculates tri metric.
908 /** aspect ratio
909  Reference --- P. P. Pebay & T. J. Baker, Analysis of Triangle Quality
910  Measures, AMS Math. Comp., 2003, 72(244):1817-1839 */
911 C_FUNC_DEF double v_tri_aspect_ratio( int num_nodes, double coordinates[][3] );
912 
913 //! Calculates tri metric.
914 /** radius ratio
915  Reference --- P. P. Pebay & T. J. Baker, Analysis of Triangle Quality
916  Measures, AMS Math. Comp., 2003, 72(244):1817-1839 */
917 C_FUNC_DEF double v_tri_radius_ratio( int num_nodes, double coordinates[][3] );
918 
919 //! Calculates tri metric.
920 /** Frobenius aspect */
921 C_FUNC_DEF double v_tri_aspect_frobenius( int num_nodes, double coordinates[][3] );
922 
923 //! Calculates tri metric.
924 /** Maximum included angle in triangle */
925 C_FUNC_DEF double v_tri_area( int num_nodes, double coordinates[][3] );
926 
927 //! Calculates tri metric.
928 /** Minimum included angle in triangle */
929 C_FUNC_DEF double v_tri_minimum_angle( int num_nodes, double coordinates[][3] );
930 
931 //! Calculates tri metric.
932 /** Maximum included angle in triangle */
933 C_FUNC_DEF double v_tri_maximum_angle( int num_nodes, double coordinates[][3] );
934 
935 //! Calculates tri metric.
936 /** Condition number of the Jacobian matrix.
937  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
938  Optimization of the Jacobian Matrix Norm and Associated Quantities,
939  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
940 C_FUNC_DEF double v_tri_condition( int num_nodes, double coordinates[][3] );
941 
942 //! Calculates tri metric.
943 /** Minimum Jacobian divided by the lengths of 2 edge vectors.
944  Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
945  Optimization of the Jacobian Matrix Norm and Associated Quantities,
946  Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
947 C_FUNC_DEF double v_tri_scaled_jacobian( int num_nodes, double coordinates[][3] );
948 
949 //! Calculates tri metric.
950 /** Min( J, 1/J ), where J is determinant of weighted Jacobian matrix.
951  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
952  Unstructured Initial Meshes, submitted for publication. */
953 C_FUNC_DEF double v_tri_relative_size_squared( int num_nodes, double coordinates[][3] );
954 
955 //! Calculates tri metric.
956 /** 2/Condition number of weighted Jacobian matrix.
957  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
958  Unstructured Initial Meshes, submitted for publication. */
959 C_FUNC_DEF double v_tri_shape( int num_nodes, double coordinates[][3] );
960 
961 //! Calculates tri metric.
962 /** Product of Shape and Relative Size.
963  Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
964  Unstructured Initial Meshes, submitted for publication. */
965 C_FUNC_DEF double v_tri_shape_and_size( int num_nodes, double coordinates[][3] );
966 
967 //! Calculates tri metric.
968 /** {min(|J|)/actual area}*parent area, parent area = 1/2 for triangular element.
969  Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
970 C_FUNC_DEF double v_tri_distortion( int num_nodes, double coordinates[][3] );
971 
972 #endif /* VERDICT_INC_LIB */