Mesh Oriented datABase
(version 5.5.1)
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
*/
79
struct
HexMetricVals
80
{
81
/** \sa v_hex_edge_ratio */
82
double
edge_ratio
;
83
/** \sa v_hex_max_edge_ratio */
84
double
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 */
100
double
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 */
106
double
scaled_jacobian
;
107
/** \sa v_hex_shear */
108
double
shear
;
109
/** \sa v_hex_shape */
110
double
shape
;
111
/** \sa v_hex_relative_size */
112
double
relative_size_squared
;
113
/** \sa v_hex_shape_and_size */
114
double
shape_and_size
;
115
/** \sa v_hex_shear_and_size */
116
double
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
*/
127
struct
EdgeMetricVals
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
*/
138
struct
KnifeMetricVals
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
164
struct
QuadMetricVals
165
{
166
/** \sa v_quad_edge_ratio function */
167
double
edge_ratio
;
168
/** \sa v_quad_max_edge_ratio function */
169
double
max_edge_ratio
;
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 */
175
double
med_aspect_frobenius
;
176
/** \sa v_quad_max_aspect_frobenius function */
177
double
max_aspect_frobenius
;
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 */
189
double
minimum_angle
;
190
/** \sa v_quad_largest_angle function */
191
double
maximum_angle
;
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 */
199
double
scaled_jacobian
;
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 */
205
double
relative_size_squared
;
206
/** \sa v_quad_shape_and_size function */
207
double
shape_and_size
;
208
/** \sa v_quad_shear_and_size function */
209
double
shear_and_size
;
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
*/
220
struct
PyramidMetricVals
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
*/
231
struct
WedgeMetricVals
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
*/
242
struct
TetMetricVals
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 */
255
double
aspect_frobenius
;
256
/** \sa v_tet_minimum_angle */
257
double
minimum_angle
;
258
/** \sa v_tet_collapse_ratio*/
259
double
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 */
267
double
scaled_jacobian
;
268
/** \sa v_tet_shape */
269
double
shape
;
270
/** \sa v_tet_relative_size */
271
double
relative_size_squared
;
272
/** \sa v_tet_shape_and_size*/
273
double
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
*/
284
struct
TriMetricVals
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 */
293
double
aspect_frobenius
;
294
/** \sa v_tri_area*/
295
double
area
;
296
/** \sa v_tri_minimum_angle*/
297
double
minimum_angle
;
298
/** \sa v_tri_maximum_angle */
299
double
maximum_angle
;
300
/** \sa v_tri_condition */
301
double
condition
;
302
/** \sa v_tri_scaled_jacobian */
303
double
scaled_jacobian
;
304
/** \sa v_tri_shape */
305
double
shape
;
306
/** \sa v_tri_relative_size_squared */
307
double
relative_size_squared
;
308
/** \sa v_tri_shape_and_size */
309
double
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
/*! \mainpage
469
Verdict is a library used to calculate metrics on the following type of elements:
470
471
\li Hexahedra
472
\li Tetrahedra
473
\li Pryamid
474
\li Wedge
475
\li Knife
476
\li Quadrilateral
477
\li Triangle
478
\li Edge
479
480
Verdict calculates individual or multiple metrics on a single element.
481
The v_*_quality(...) functions allow for efficient calculations of
482
multiple metrics on a single element. Individual metrics may be
483
calculated on a single element as well.
484
485
\section jack Using Verdict
486
487
The v_*_quality functions take the following parameters:
488
489
\param num_nodes Number of nodes in the element.
490
\param coordinates 2D array containing x,y,z coordinate data of the nodes.
491
\param metrics_request_flag Bitfield to define which metrics to calculate.
492
\param metric_vals Struct to hold the metric calculated values.
493
494
All other functions take these parameters below and return the calculated
495
metric value.
496
497
\param num_nodes Number of nodes in the element.
498
\param coordinates 2D array containing x,y,z coordinate data of the nodes.
499
500
501
\par Setting the metrics_request_flag:
502
In order to use v_*_quality functions you must know how to set the bitfield argument
503
correctly. To calculate aspect ratio, condition number, shape and shear metrics on a triangle, set
504
the "metrics_request_flag" like so:
505
506
\code
507
unsigned int metrics_request_flag = 0;
508
metrics_request_flag += V_TRI_ASPECT_FROBENIUS;
509
metrics_request_flag += V_CONDITION;
510
metrics_request_flag += V_SHAPE;
511
metrics_request_flag += V_SHEAR;
512
\endcode
513
514
The bitwise field can also be set for many metrics at once using #deinfed numbers. V_HEX_ALL,
515
V_HEX_DIAGNOSTIC, V_TRI_ALGEBRAIC are examples.
516
517
Below is an example of how use Verdict's functions:
518
519
520
Example: \code
521
QuadMetricVals quad_metrics = {0};
522
unsigned long metrics_flag = 0;
523
metrics_flag += V_QUAD_SHAPE;
524
metrics_flag += V_QUAD_DISTORTION;
525
metrics_flag += V_QUAD_AREA;
526
double quad_nodes[4][3];
527
528
//node 1
529
quad_node[0][0] = 0; //x
530
quad_node[0][1] = 0; //y
531
quad_node[0][2] = 0; //z
532
533
//node 2
534
quad_node[1][0] = 1;
535
quad_node[1][1] = 0.1;
536
quad_node[1][2] = 0.1;
537
538
//node 3
539
quad_node[2][0] = 0.9;
540
quad_node[2][1] = 0.9;
541
quad_node[2][2] = -0.1;
542
543
//node 4
544
quad_node[3][0] = -0.05;
545
quad_node[3][1] = 1;
546
quad_node[3][2] = 0;
547
548
//calculate multiple metrics with one call
549
v_quad_quality( 4, quad_nodes, metrics_flag, quad_metrics );
550
double my_shape = quad_metrics.shape;
551
double my_distortion = quad_metrics.distortion;
552
double my_area = quad_metrics.area;
553
554
//calculate an individual metric
555
double my_relative_size = v_quad_relative_size( 4, quad_nodes );
556
\endcode
557
558
*/
559
560
//! Calculates quality metrics for hexahedral elements.
561
C_FUNC_DEF
void
v_hex_quality
(
int
num_nodes,
562
double
coordinates[][3],
563
unsigned
int
metrics_request_flag,
564
struct
HexMetricVals
* metric_vals );
565
566
//! Calculates quality metrics for tetrahedral elements.
567
C_FUNC_DEF
void
v_tet_quality
(
int
num_nodes,
568
double
coordinates[][3],
569
unsigned
int
metrics_request_flag,
570
struct
TetMetricVals
* metric_vals );
571
572
//! Calculates quality metrics for pyramid elements.
573
C_FUNC_DEF
void
v_pyramid_quality
(
int
num_nodes,
574
double
coordinates[][3],
575
unsigned
int
metrics_request_flag,
576
struct
PyramidMetricVals
* metric_vals );
577
578
//! Calculates quality metrics for wedge elements.
579
C_FUNC_DEF
void
v_wedge_quality
(
int
num_nodes,
580
double
coordinates[][3],
581
unsigned
int
metrics_request_flag,
582
struct
WedgeMetricVals
* metric_vals );
583
584
//! Calculates quality metrics for knife elements.
585
C_FUNC_DEF
void
v_knife_quality
(
int
num_nodes,
586
double
coordinates[][3],
587
unsigned
int
metrics_request_flag,
588
struct
KnifeMetricVals
* metric_vals );
589
590
//! Calculates quality metrics for quadrilateral elements.
591
C_FUNC_DEF
void
v_quad_quality
(
int
num_nodes,
592
double
coordinates[][3],
593
unsigned
int
metrics_request_flag,
594
struct
QuadMetricVals
* metric_vals );
595
596
//! Calculates quality metrics for triangle elements.
597
C_FUNC_DEF
void
v_tri_quality
(
int
num_nodes,
598
double
coordinates[][3],
599
unsigned
int
metrics_request_flag,
600
struct
TriMetricVals
* metric_vals );
601
602
//! Calculates quality metrics for edge elements.
603
C_FUNC_DEF
void
v_edge_quality
(
int
num_nodes,
604
double
coordinates[][3],
605
unsigned
int
metrics_request_flag,
606
struct
EdgeMetricVals
* metric_vals );
607
608
/* individual quality functions for hex elements */
609
610
//! Sets average size (volume) of hex, needed for v_hex_relative_size(...)
611
C_FUNC_DEF
void
v_set_hex_size
(
double
size
);
612
613
//! Calculates hex edge ratio metric.
614
/** Hmax / Hmin where Hmax and Hmin are respectively the maximum and the
615
minimum edge lengths */
616
C_FUNC_DEF
double
v_hex_edge_ratio
(
int
num_nodes,
double
coordinates[][3] );
617
618
//! Calculates hex maximum of edge ratio
619
/**Maximum edge length ratio at hex center.
620
Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
621
Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
622
C_FUNC_DEF
double
v_hex_max_edge_ratio
(
int
num_nodes,
double
coordinates[][3] );
623
624
//! Calculates hex skew metric.
625
/** Maximum |cos A| where A is the angle between edges at hex center.
626
Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
627
Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
628
C_FUNC_DEF
double
v_hex_skew
(
int
num_nodes,
double
coordinates[][3] );
629
630
//! Calculates hex taper metric
631
/** Maximum ratio of lengths derived from opposite edges.
632
Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
633
Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
634
C_FUNC_DEF
double
v_hex_taper
(
int
num_nodes,
double
coordinates[][3] );
635
636
//! Calculates hex volume
637
/** Jacobian at hex center.
638
Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
639
Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
640
C_FUNC_DEF
double
v_hex_volume
(
int
num_nodes,
double
coordinates[][3] );
641
642
//! Calculates hex stretch metric
643
/** Sqrt(3) * minimum edge length / maximum diagonal length.
644
Reference --- FIMESH code */
645
C_FUNC_DEF
double
v_hex_stretch
(
int
num_nodes,
double
coordinates[][3] );
646
647
//! Calculates hex diagonal metric
648
/** Minimum diagonal length / maximum diagonal length.
649
Reference --- Unknown */
650
C_FUNC_DEF
double
v_hex_diagonal
(
int
num_nodes,
double
coordinates[][3] );
651
652
//! Calculates hex dimension metric
653
/** Pronto-specific characteristic length for stable time step calculation.
654
Char_length = Volume / 2 grad Volume.
655
Reference --- L.M. Taylor, and D.P. Flanagan, Pronto3D - A Three Dimensional Transient
656
Solid Dynamics Program, SAND87-1912, Sandia National Laboratories, 1989. */
657
C_FUNC_DEF
double
v_hex_dimension
(
int
num_nodes,
double
coordinates[][3] );
658
659
//! Calculates hex oddy metric
660
C_FUNC_DEF
double
v_hex_oddy
(
int
num_nodes,
double
coordinates[][3] );
661
662
//! Calculates hex condition metric
663
/** Average Frobenius condition number of the Jacobian matrix at 8 corners. */
664
C_FUNC_DEF
double
v_hex_med_aspect_frobenius
(
int
num_nodes,
double
coordinates[][3] );
665
666
//! Calculates hex condition metric
667
/** Maximum Frobenius condition number of the Jacobian matrix at 8 corners.
668
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
669
Optimization of the Jacobian Matrix Norm and Associated Quantities,
670
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
671
C_FUNC_DEF
double
v_hex_max_aspect_frobenius
(
int
num_nodes,
double
coordinates[][3] );
672
C_FUNC_DEF
double
v_hex_condition
(
int
num_nodes,
double
coordinates[][3] );
673
674
//! Calculates hex jacobian metric
675
/** Minimum pointwise volume of local map at 8 corners & center of hex.
676
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
677
Optimization of the Jacobian Matrix Norm and Associated Quantities,
678
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
679
C_FUNC_DEF
double
v_hex_jacobian
(
int
num_nodes,
double
coordinates[][3] );
680
681
//! Calculates hex scaled jacobian metric
682
/** Minimum Jacobian divided by the lengths of the 3 edge vectors.
683
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
684
Optimization of the Jacobian Matrix Norm and Associated Quantities,
685
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
686
C_FUNC_DEF
double
v_hex_scaled_jacobian
(
int
num_nodes,
double
coordinates[][3] );
687
688
//! Calculates hex shear metric
689
/** 3/Mean Ratio of Jacobian Skew matrix.
690
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
691
Unstructured Initial Meshes, submitted for publication. */
692
C_FUNC_DEF
double
v_hex_shear
(
int
num_nodes,
double
coordinates[][3] );
693
694
//! Calculates hex shape metric.
695
/** 3/Mean Ratio of weighted Jacobian matrix.
696
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
697
Unstructured Initial Meshes, submitted for publication. */
698
C_FUNC_DEF
double
v_hex_shape
(
int
num_nodes,
double
coordinates[][3] );
699
700
//! Calculates hex relative size metric.
701
/** 3/Mean Ratio of weighted Jacobian matrix.
702
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
703
Unstructured Initial Meshes, submitted for publication. */
704
C_FUNC_DEF
double
v_hex_relative_size_squared
(
int
num_nodes,
double
coordinates[][3] );
705
706
//! Calculates hex shape-size metric.
707
/** Product of Shape and Relative Size.
708
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
709
Unstructured Initial Meshes, submitted for publication. */
710
C_FUNC_DEF
double
v_hex_shape_and_size
(
int
num_nodes,
double
coordinates[][3] );
711
712
//! Calculates hex shear-size metric
713
/** Product of Shear and Relative Size.
714
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
715
Unstructured Initial Meshes, submitted for publication. */
716
C_FUNC_DEF
double
v_hex_shear_and_size
(
int
num_nodes,
double
coordinates[][3] );
717
718
//! Calculates hex distortion metric
719
/** {min(|J|)/actual volume}*parent volume, parent volume = 8 for hex.
720
Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
721
C_FUNC_DEF
double
v_hex_distortion
(
int
num_nodes,
double
coordinates[][3] );
722
723
/* individual quality functions for tet elements */
724
725
//! Sets average size (volume) of tet, needed for v_tet_relative_size(...)
726
C_FUNC_DEF
void
v_set_tet_size
(
double
size
);
727
728
//! Calculates tet edge ratio metric.
729
/** Hmax / Hmin where Hmax and Hmin are respectively the maximum and the
730
minimum edge lengths */
731
C_FUNC_DEF
double
v_tet_edge_ratio
(
int
num_nodes,
double
coordinates[][3] );
732
733
//! Calculates tet radius ratio metric.
734
/** CR / (3.0 * IR) where CR = circumsphere radius, IR = inscribed sphere radius.
735
Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
736
quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
737
C_FUNC_DEF
double
v_tet_radius_ratio
(
int
num_nodes,
double
coordinates[][3] );
738
739
//! Calculates the radius ratio metric of a positively oriented tet.
740
/** CR / (3.0 * IR) where CR = circumsphere radius, IR = inscribed sphere radius
741
if the element is positively-oriented.
742
Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
743
quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
744
C_FUNC_DEF
double
v_tet_aspect_beta
(
int
num_nodes,
double
coordinates[][3] );
745
746
//! Calculates tet aspect ratio metric.
747
/** Hmax / (2 sqrt(6) r) where Hmax and r respectively denote the greatest edge
748
length and the inradius of the tetrahedron
749
Reference --- P. Frey and P.-L. George, Meshing, Hermes (2000). */
750
C_FUNC_DEF
double
v_tet_aspect_ratio
(
int
num_nodes,
double
coordinates[][3] );
751
752
//! Calculates tet aspect gamma metric.
753
/** Srms**3 / (8.479670*V) where Srms = sqrt(Sum(Si**2)/6), Si = edge length.
754
Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
755
quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
756
C_FUNC_DEF
double
v_tet_aspect_gamma
(
int
num_nodes,
double
coordinates[][3] );
757
758
//! Calculates tet aspect frobenius metric.
759
/** Frobenius condition number when the reference element is regular
760
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
761
Optimization of the Jacobian Matrix Norm and Associated Quantities,
762
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
763
C_FUNC_DEF
double
v_tet_aspect_frobenius
(
int
num_nodes,
double
coordinates[][3] );
764
765
//! Calculates tet minimum dihedral angle.
766
/** Minimum (nonoriented) dihedral angle of a tetrahedron, expressed in degrees. */
767
C_FUNC_DEF
double
v_tet_minimum_angle
(
int
num_nodes,
double
coordinates[][3] );
768
769
//! Calculates tet collapse ratio metric.
770
/** Collapse ratio */
771
C_FUNC_DEF
double
v_tet_collapse_ratio
(
int
num_nodes,
double
coordinates[][3] );
772
773
//! Calculates tet volume.
774
/** (1/6) * Jacobian at corner node.
775
Reference --- V. N. Parthasarathy et al, A comparison of tetrahedron
776
quality measures, Finite Elem. Anal. Des., Vol 15(1993), 255-261. */
777
C_FUNC_DEF
double
v_tet_volume
(
int
num_nodes,
double
coordinates[][3] );
778
779
//! Calculates tet condition metric.
780
/** Condition number of the Jacobian matrix at any corner.
781
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
782
Optimization of the Jacobian Matrix Norm and Associated Quantities,
783
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
784
C_FUNC_DEF
double
v_tet_condition
(
int
num_nodes,
double
coordinates[][3] );
785
786
//! Calculates tet jacobian.
787
/** Minimum pointwise volume at any corner.
788
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
789
Optimization of the Jacobian Matrix Norm and Associated Quantities,
790
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
791
C_FUNC_DEF
double
v_tet_jacobian
(
int
num_nodes,
double
coordinates[][3] );
792
793
//! Calculates tet scaled jacobian.
794
/** Minimum Jacobian divided by the lengths of 3 edge vectors
795
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
796
Optimization of the Jacobian Matrix Norm and Associated Quantities,
797
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
798
C_FUNC_DEF
double
v_tet_scaled_jacobian
(
int
num_nodes,
double
coordinates[][3] );
799
800
//! Calculates tet shape metric.
801
/** 3/Mean Ratio of weighted Jacobian matrix.
802
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
803
Unstructured Initial Meshes, submitted for publication. */
804
C_FUNC_DEF
double
v_tet_shape
(
int
num_nodes,
double
coordinates[][3] );
805
806
//! Calculates tet relative size metric.
807
/** Min( J, 1/J ), where J is determinant of weighted Jacobian matrix.
808
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
809
Unstructured Initial Meshes, submitted for publication. */
810
C_FUNC_DEF
double
v_tet_relative_size_squared
(
int
num_nodes,
double
coordinates[][3] );
811
812
//! Calculates tet shape-size metric.
813
/** Product of Shape and Relative Size.
814
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
815
Unstructured Initial Meshes, submitted for publication. */
816
C_FUNC_DEF
double
v_tet_shape_and_size
(
int
num_nodes,
double
coordinates[][3] );
817
818
//! Calculates tet distortion metric.
819
/** {min(|J|)/actual volume}*parent volume, parent volume = 1/6 for tet.
820
Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
821
C_FUNC_DEF
double
v_tet_distortion
(
int
num_nodes,
double
coordinates[][3] );
822
823
/* individual quality functions for pyramid elements */
824
825
//! Calculates pyramid volume.
826
C_FUNC_DEF
double
v_pyramid_volume
(
int
num_nodes,
double
coordinates[][3] );
827
828
/* individual quality functions for wedge elements */
829
830
//! Calculates wedge volume.
831
C_FUNC_DEF
double
v_wedge_volume
(
int
num_nodes,
double
coordinates[][3] );
832
833
/* individual quality functions for knife elements */
834
835
//! Calculates knife volume.
836
C_FUNC_DEF
double
v_knife_volume
(
int
num_nodes,
double
coordinates[][3] );
837
838
/* individual quality functions for edge elements */
839
840
//! Calculates edge length.
841
C_FUNC_DEF
double
v_edge_length
(
int
num_nodes,
double
coordinates[][3] );
842
843
/* individual quality functions for quad elements */
844
//! Sets average size (area) of quad, needed for v_quad_relative_size(...)
845
C_FUNC_DEF
void
v_set_quad_size
(
double
size
);
846
847
//! Calculates quad edge ratio
848
/** edge ratio
849
Reference --- P. P. Pebay, Planar Quadrangle Quality
850
Measures, Eng. Comp., 2004, 20(2):157-173 */
851
C_FUNC_DEF
double
v_quad_edge_ratio
(
int
num_nodes,
double
coordinates[][3] );
852
853
//! Calculates quad maximum of edge ratio.
854
/** Maximum edge length ratio at quad center.
855
Reference --- J. Robinson, CRE Method of element testing and the
856
Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
857
C_FUNC_DEF
double
v_quad_max_edge_ratio
(
int
num_nodes,
double
coordinates[][3] );
858
859
//! Calculates quad aspect ratio
860
/** aspect ratio
861
Reference --- P. P. Pebay, Planar Quadrangle Quality
862
Measures, Eng. Comp., 2004, 20(2):157-173 */
863
C_FUNC_DEF
double
v_quad_aspect_ratio
(
int
num_nodes,
double
coordinates[][3] );
864
865
//! Calculates quad radius ratio
866
/** radius ratio
867
Reference --- P. P. Pebay, Planar Quadrangle Quality
868
Measures, Eng. Comp., 2004, 20(2):157-173 */
869
C_FUNC_DEF
double
v_quad_radius_ratio
(
int
num_nodes,
double
coordinates[][3] );
870
871
//! Calculates quad average Frobenius aspect
872
/** average Frobenius aspect
873
Reference --- P. P. Pebay, Planar Quadrangle Quality
874
Measures, Eng. Comp., 2004, 20(2):157-173 */
875
C_FUNC_DEF
double
v_quad_med_aspect_frobenius
(
int
num_nodes,
double
coordinates[][3] );
876
877
//! Calculates quad maximum Frobenius aspect
878
/** average Frobenius aspect
879
Reference --- P. P. Pebay, Planar Quadrangle Quality
880
Measures, Eng. Comp., 2004, 20(2):157-173 */
881
C_FUNC_DEF
double
v_quad_max_aspect_frobenius
(
int
num_nodes,
double
coordinates[][3] );
882
883
//! Calculates quad skew metric.
884
/** Maximum |cos A| where A is the angle between edges at quad center.
885
Reference --- J. Robinson, CRE Method of element testing and the
886
Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
887
C_FUNC_DEF
double
v_quad_skew
(
int
num_nodes,
double
coordinates[][3] );
888
889
//! Calculates quad taper metric.
890
/** Maximum ratio of lengths derived from opposite edges.
891
Reference --- J. Robinson, CRE Method of element testing and the
892
Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
893
C_FUNC_DEF
double
v_quad_taper
(
int
num_nodes,
double
coordinates[][3] );
894
895
//! Calculates quad warpage metric.
896
/** Cosine of Minimum Dihedral Angle formed by Planes Intersecting in Diagonals.
897
Reference --- J. Robinson, CRE Method of element testing and the
898
Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
899
C_FUNC_DEF
double
v_quad_warpage
(
int
num_nodes,
double
coordinates[][3] );
900
901
//! Calculates quad area.
902
/** Jacobian at quad center.
903
Reference --- J. Robinson, CRE Method of element testing and the
904
Jacobian shape parameters, Eng. Comput., Vol 4, 1987. */
905
C_FUNC_DEF
double
v_quad_area
(
int
num_nodes,
double
coordinates[][3] );
906
907
//! Calculates quad strech metric.
908
/** Sqrt(2) * minimum edge length / maximum diagonal length.
909
Reference --- FIMESH code. */
910
C_FUNC_DEF
double
v_quad_stretch
(
int
num_nodes,
double
coordinates[][3] );
911
912
//! Calculates quad's smallest angle.
913
/** Smallest included quad angle (degrees).
914
Reference --- Unknown. */
915
C_FUNC_DEF
double
v_quad_minimum_angle
(
int
num_nodes,
double
coordinates[][3] );
916
917
//! Calculates quad's largest angle.
918
/** Largest included quad angle (degrees).
919
Reference --- Unknown. */
920
C_FUNC_DEF
double
v_quad_maximum_angle
(
int
num_nodes,
double
coordinates[][3] );
921
922
//! Calculates quad oddy metric.
923
C_FUNC_DEF
double
v_quad_oddy
(
int
num_nodes,
double
coordinates[][3] );
924
925
//! Calculates quad condition number metric.
926
/** Maximum condition number of the Jacobian matrix at 4 corners.
927
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
928
Optimization of the Jacobian Matrix Norm and Associated Quantities,
929
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
930
C_FUNC_DEF
double
v_quad_condition
(
int
num_nodes,
double
coordinates[][3] );
931
932
//! Calculates quad jacobian.
933
/** Minimum pointwise volume of local map at 4 corners & center of quad.
934
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
935
Optimization of the Jacobian Matrix Norm and Associated Quantities,
936
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
937
C_FUNC_DEF
double
v_quad_jacobian
(
int
num_nodes,
double
coordinates[][3] );
938
939
//! Calculates quad scaled jacobian.
940
/** Minimum Jacobian divided by the lengths of the 2 edge vectors.
941
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
942
Optimization of the Jacobian Matrix Norm and Associated Quantities,
943
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
944
C_FUNC_DEF
double
v_quad_scaled_jacobian
(
int
num_nodes,
double
coordinates[][3] );
945
946
//! Calculates quad shear metric.
947
/** 2/Condition number of Jacobian Skew matrix.
948
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
949
Unstructured Initial Meshes, submitted for publication. */
950
C_FUNC_DEF
double
v_quad_shear
(
int
num_nodes,
double
coordinates[][3] );
951
952
//! Calculates quad shape metric.
953
/** 2/Condition number of weighted Jacobian matrix.
954
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
955
Unstructured Initial Meshes, submitted for publication. */
956
C_FUNC_DEF
double
v_quad_shape
(
int
num_nodes,
double
coordinates[][3] );
957
958
//! Calculates quad relative size metric.
959
/** Min( J, 1/J ), where J is determinant of weighted Jacobian matrix.
960
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
961
Unstructured Initial Meshes, submitted for publication. */
962
C_FUNC_DEF
double
v_quad_relative_size_squared
(
int
num_nodes,
double
coordinates[][3] );
963
964
//! Calculates quad shape-size metric.
965
/** Product of Shape and Relative Size.
966
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
967
Unstructured Initial Meshes, submitted for publication. */
968
C_FUNC_DEF
double
v_quad_shape_and_size
(
int
num_nodes,
double
coordinates[][3] );
969
970
//! Calculates quad shear-size metric.
971
/** Product of Shear and Relative Size.
972
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
973
Unstructured Initial Meshes, submitted for publication. */
974
C_FUNC_DEF
double
v_quad_shear_and_size
(
int
num_nodes,
double
coordinates[][3] );
975
976
//! Calculates quad distortion metric.
977
/** {min(|J|)/actual area}*parent area, parent area = 4 for quad.
978
Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
979
C_FUNC_DEF
double
v_quad_distortion
(
int
num_nodes,
double
coordinates[][3] );
980
981
/* individual quality functions for tri elements */
982
983
//! Sets average size (area) of tri, needed for v_tri_relative_size(...)
984
C_FUNC_DEF
void
v_set_tri_size
(
double
size
);
985
986
//! Sets fuction pointer to calculate tri normal wrt surface
987
C_FUNC_DEF
void
v_set_tri_normal_func
(
ComputeNormal
func );
988
989
//! Calculates tri metric.
990
/** edge ratio
991
Reference --- P. P. Pebay & T. J. Baker, Analysis of Triangle Quality
992
Measures, AMS Math. Comp., 2003, 72(244):1817-1839 */
993
C_FUNC_DEF
double
v_tri_edge_ratio
(
int
num_nodes,
double
coordinates[][3] );
994
995
//! Calculates tri metric.
996
/** aspect ratio
997
Reference --- P. P. Pebay & T. J. Baker, Analysis of Triangle Quality
998
Measures, AMS Math. Comp., 2003, 72(244):1817-1839 */
999
C_FUNC_DEF
double
v_tri_aspect_ratio
(
int
num_nodes,
double
coordinates[][3] );
1000
1001
//! Calculates tri metric.
1002
/** radius ratio
1003
Reference --- P. P. Pebay & T. J. Baker, Analysis of Triangle Quality
1004
Measures, AMS Math. Comp., 2003, 72(244):1817-1839 */
1005
C_FUNC_DEF
double
v_tri_radius_ratio
(
int
num_nodes,
double
coordinates[][3] );
1006
1007
//! Calculates tri metric.
1008
/** Frobenius aspect */
1009
C_FUNC_DEF
double
v_tri_aspect_frobenius
(
int
num_nodes,
double
coordinates[][3] );
1010
1011
//! Calculates tri metric.
1012
/** Maximum included angle in triangle */
1013
C_FUNC_DEF
double
v_tri_area
(
int
num_nodes,
double
coordinates[][3] );
1014
1015
//! Calculates tri metric.
1016
/** Minimum included angle in triangle */
1017
C_FUNC_DEF
double
v_tri_minimum_angle
(
int
num_nodes,
double
coordinates[][3] );
1018
1019
//! Calculates tri metric.
1020
/** Maximum included angle in triangle */
1021
C_FUNC_DEF
double
v_tri_maximum_angle
(
int
num_nodes,
double
coordinates[][3] );
1022
1023
//! Calculates tri metric.
1024
/** Condition number of the Jacobian matrix.
1025
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
1026
Optimization of the Jacobian Matrix Norm and Associated Quantities,
1027
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
1028
C_FUNC_DEF
double
v_tri_condition
(
int
num_nodes,
double
coordinates[][3] );
1029
1030
//! Calculates tri metric.
1031
/** Minimum Jacobian divided by the lengths of 2 edge vectors.
1032
Reference --- P. Knupp, Achieving Finite Element Mesh Quality via
1033
Optimization of the Jacobian Matrix Norm and Associated Quantities,
1034
Intl. J. Numer. Meth. Engng. 2000, 48:1165-1185. */
1035
C_FUNC_DEF
double
v_tri_scaled_jacobian
(
int
num_nodes,
double
coordinates[][3] );
1036
1037
//! Calculates tri metric.
1038
/** Min( J, 1/J ), where J is determinant of weighted Jacobian matrix.
1039
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
1040
Unstructured Initial Meshes, submitted for publication. */
1041
C_FUNC_DEF
double
v_tri_relative_size_squared
(
int
num_nodes,
double
coordinates[][3] );
1042
1043
//! Calculates tri metric.
1044
/** 2/Condition number of weighted Jacobian matrix.
1045
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
1046
Unstructured Initial Meshes, submitted for publication. */
1047
C_FUNC_DEF
double
v_tri_shape
(
int
num_nodes,
double
coordinates[][3] );
1048
1049
//! Calculates tri metric.
1050
/** Product of Shape and Relative Size.
1051
Reference --- P. Knupp, Algebraic Mesh Quality Metrics for
1052
Unstructured Initial Meshes, submitted for publication. */
1053
C_FUNC_DEF
double
v_tri_shape_and_size
(
int
num_nodes,
double
coordinates[][3] );
1054
1055
//! Calculates tri metric.
1056
/** {min(|J|)/actual area}*parent area, parent area = 1/2 for triangular element.
1057
Reference --- SDRC/IDEAS Simulation: Finite Element Modeling--User's Guide */
1058
C_FUNC_DEF
double
v_tri_distortion
(
int
num_nodes,
double
coordinates[][3] );
1059
1060
#endif
/* VERDICT_INC_LIB */
src
moab
verdict.h
Generated on Tue Oct 29 2024 02:05:52 for Mesh Oriented datABase by
1.9.1.