1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
2529
2530
2531
2532
2533
2534
2535
2536
2537
2538
2539
2540
2541
2542
2543
2544
2545
2546
2547
2548
2549
2550
2551
2552
2553
2554
2555
2556
2557
2558
2559
2560
2561
2562
2563
2564
2565
2566
2567
2568
2569
2570
2571
2572
2573
2574
2575
2576
2577
2578
2579
2580
2581
2582
2583
2584
2585
2586
2587
2588
2589
2590
2591
2592
2593
2594
2595
2596
2597
2598
2599
2600
2601
2602
2603
2604
2605
2606
2607
2608
2609
2610
2611
2612
2613
2614
2615
2616
2617
2618
2619
2620
2621
2622
2623
2624
2625
2626
2627
2628
2629
2630
2631
2632
2633
2634
2635
2636
2637
2638
2639
2640
2641
2642
2643
2644
2645
2646
2647
2648
2649
2650
2651
2652
2653
2654
2655
2656
2657
2658
2659
2660
2661
2662
2663
2664
2665
2666
2667
2668
2669
2670
2671
2672
2673
2674
2675
2676
2677
2678
2679
2680
2681
2682
2683
2684
2685
2686
2687
2688
2689
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701
2702
2703
2704
2705
2706
2707
2708
2709
2710
2711
2712
2713
2714
2715
2716
2717
2718
2719
2720
2721
2722
2723
2724
2725
2726
2727
2728
2729
2730
2731
2732
2733
2734
2735
2736
2737
2738
2739
2740
2741
2742
2743
2744
2745
2746
2747
2748
2749
2750
2751
2752
2753
2754
2755
2756
2757
2758
2759
2760
2761
2762
2763
2764
2765
2766
2767
2768
2769
2770
2771
2772
2773
2774
2775
2776
2777
2778
2779
2780
2781
2782
2783
2784
2785
2786
2787
2788
2789
2790
2791
2792
2793
2794
2795
2796
2797
2798
2799
2800
2801
2802
2803
2804
2805
2806
2807
2808
2809
2810
2811
2812
2813
2814
2815
2816
2817
2818
2819
2820
2821
2822
2823
2824
2825
2826
2827
2828
2829
2830
2831
2832
2833
2834
2835
2836
2837
2838
2839
2840
2841
2842
2843
2844
2845
2846
2847
2848
2849
2850
2851
2852
2853
2854
2855
2856
2857
2858
2859
2860
2861
2862
2863
2864
2865
2866
2867
2868
2869
2870
2871
2872
2873
2874
2875
2876
2877
2878
2879
2880
2881
2882
2883
2884
2885
2886
2887
2888
2889
2890
2891
2892
2893
2894
2895
2896
2897
2898
2899
2900
2901
2902
2903
2904
2905
2906
2907
2908
2909
2910
2911
2912
2913
2914
2915
2916
2917
2918
2919
2920
2921
2922
2923
2924
2925
2926
2927
2928
2929
2930
2931
2932
2933
2934
2935
2936
2937
2938
2939
2940
2941
2942
2943
2944
2945
2946
2947
2948
2949
2950
2951
2952
2953
2954
2955
2956
2957
2958
2959
2960
2961
2962
2963
2964
2965
2966
2967
2968
2969
2970
2971
2972
2973
2974
2975
2976
2977
2978
2979
2980
2981
2982
2983
2984
2985
2986
2987
2988
2989
2990
2991
2992
2993
2994
2995
2996
2997
2998
2999
3000
3001
3002
3003
3004
3005
3006
3007
3008
3009
3010
3011
3012
3013
3014
3015
3016
3017
3018
3019
3020
3021
3022
3023
3024
3025
3026
3027
3028
3029
3030
3031
3032
3033
3034
3035
3036
3037
3038
3039
3040
3041
3042
3043
3044
3045
3046
3047
3048
3049
3050
3051
3052
3053
3054
3055
3056
3057
3058
3059
3060
3061
3062
3063
3064
3065
3066
3067
3068
3069
3070
3071
3072
3073
3074
3075
3076
3077
3078
3079
3080
3081
3082
3083
3084
3085
3086
3087
3088
3089
3090
3091
3092
3093
3094
3095
3096
3097
3098
3099
3100
3101
3102
3103
3104
3105
3106
3107
3108
3109
3110
3111
3112
3113
3114
3115
3116
3117
3118
3119
3120
3121
3122
3123
3124
3125
3126
3127
3128
3129
3130
3131
3132
3133
3134
3135
3136
3137
3138
3139
3140
3141
3142
3143
3144
3145
3146
3147
3148
3149
3150
3151
3152
3153
3154
3155
3156
3157
3158
3159
3160
3161
3162
3163
3164
3165
3166
3167
3168
3169
3170
3171
3172
3173
3174
3175
3176
3177
3178
3179
3180
3181
3182
3183
3184
3185
3186
3187
3188
3189
3190
3191
3192
3193
3194
3195
3196
3197
3198
3199
3200
3201
3202
3203
3204
3205
3206
3207
3208
3209
3210
3211
3212
3213
3214
3215
3216
3217
3218
3219
3220
3221
3222
3223
3224
3225
3226
3227
3228
3229
3230
3231
3232
3233
3234
3235
3236
3237
3238
3239
3240
3241
3242
3243
3244
3245
3246
3247
3248
3249
3250
3251
3252
3253
3254
3255
3256
3257
3258
3259
3260
3261
3262
3263
3264
3265
3266
3267
3268
3269
3270
3271
3272
3273
3274
3275
3276
3277
3278
3279
3280
3281
3282
3283
3284
3285
3286
3287
3288
3289
3290
3291
3292
3293
3294
3295
3296
3297
3298
3299
3300
3301
3302
3303
3304
3305
3306
3307
3308
3309
3310
3311
3312
3313
3314
3315
3316
3317
3318
3319
3320
3321
3322
3323
3324
3325
3326
3327
3328
3329
3330
3331
3332
3333
3334
3335
3336
3337
3338
3339
3340
3341
3342
3343
3344
3345
3346
3347
3348
3349
3350
3351
3352
3353
3354
3355
3356
3357
3358
3359
3360
3361
3362
3363
3364
3365
3366
3367
3368
3369
3370
3371
3372
3373
3374
3375
3376
3377
3378
3379
3380
3381
3382
3383
3384
3385
3386
3387
3388
3389
3390
3391
3392
3393
3394
3395
3396
3397
3398
3399
3400
3401 | /*=========================================================================
Module: $RCSfile: V_HexMetric.cpp,v $
Copyright (c) 2006 Sandia Corporation.
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
/*
*
* HexMetric.cpp contains quality calculations for hexes
*
* This file is part of VERDICT
*
*/
#define VERDICT_EXPORTS
#include "moab/verdict.h"
#include "VerdictVector.hpp"
#include "V_GaussIntegration.hpp"
#include "verdict_defines.hpp"
#include <memory.h>
#if defined( __BORLANDC__ )
#pragma warn - 8004 /* "assigned a value that is never used" */
#endif
//! the average volume of a hex
double verdict_hex_size = 0;
//! weights based on the average size of a hex
int v_hex_get_weight( VerdictVector& v1, VerdictVector& v2, VerdictVector& v3 )
{
if( verdict_hex_size == 0 ) return 0;
v1.set( 1, 0, 0 );
v2.set( 0, 1, 0 );
v3.set( 0, 0, 1 );
double scale = pow( verdict_hex_size / ( v1 % ( v2 * v3 ) ), 0.33333333333 );
v1 *= scale;
v2 *= scale;
v3 *= scale;
return 1;
}
//! returns the average volume of a hex
C_FUNC_DEF void v_set_hex_size( double size )
{
verdict_hex_size = size;
}
#define make_hex_nodes( coord, pos ) \
for( int mhcii = 0; mhcii < 8; mhcii++ ) \
{ \
( pos )[mhcii].set( ( coord )[mhcii][0], ( coord )[mhcii][1], ( coord )[mhcii][2] ); \
}
#define make_edge_length_squares( edges, lengths ) \
{ \
for( int melii = 0; melii < 12; melii++ ) \
( lengths )[melii] = ( edges )[melii].length_squared(); \
}
//! make VerdictVectors from coordinates
void make_hex_edges( double coordinates[][3], VerdictVector edges[12] )
{
edges[0].set( coordinates[1][0] - coordinates[0][0], coordinates[1][1] - coordinates[0][1],
coordinates[1][2] - coordinates[0][2] );
edges[1].set( coordinates[2][0] - coordinates[1][0], coordinates[2][1] - coordinates[1][1],
coordinates[2][2] - coordinates[1][2] );
edges[2].set( coordinates[3][0] - coordinates[2][0], coordinates[3][1] - coordinates[2][1],
coordinates[3][2] - coordinates[2][2] );
edges[3].set( coordinates[0][0] - coordinates[3][0], coordinates[0][1] - coordinates[3][1],
coordinates[0][2] - coordinates[3][2] );
edges[4].set( coordinates[5][0] - coordinates[4][0], coordinates[5][1] - coordinates[4][1],
coordinates[5][2] - coordinates[4][2] );
edges[5].set( coordinates[6][0] - coordinates[5][0], coordinates[6][1] - coordinates[5][1],
coordinates[6][2] - coordinates[5][2] );
edges[6].set( coordinates[7][0] - coordinates[6][0], coordinates[7][1] - coordinates[6][1],
coordinates[7][2] - coordinates[6][2] );
edges[7].set( coordinates[4][0] - coordinates[7][0], coordinates[4][1] - coordinates[7][1],
coordinates[4][2] - coordinates[7][2] );
edges[8].set( coordinates[4][0] - coordinates[0][0], coordinates[4][1] - coordinates[0][1],
coordinates[4][2] - coordinates[0][2] );
edges[9].set( coordinates[5][0] - coordinates[1][0], coordinates[5][1] - coordinates[1][1],
coordinates[5][2] - coordinates[1][2] );
edges[10].set( coordinates[6][0] - coordinates[2][0], coordinates[6][1] - coordinates[2][1],
coordinates[6][2] - coordinates[2][2] );
edges[11].set( coordinates[7][0] - coordinates[3][0], coordinates[7][1] - coordinates[3][1],
coordinates[7][2] - coordinates[3][2] );
}
/*!
localizes hex coordinates
*/
void localize_hex_coordinates( double coordinates[][3], VerdictVector position[8] )<--- The function 'localize_hex_coordinates' is never used.
{
int ii;
for( ii = 0; ii < 8; ii++ )
{
position[ii].set( coordinates[ii][0], coordinates[ii][1], coordinates[ii][2] );
}
// ... Make centroid of element the center of coordinate system
VerdictVector point_1256 = position[1];
point_1256 += position[2];
point_1256 += position[5];
point_1256 += position[6];
VerdictVector point_0374 = position[0];
point_0374 += position[3];
point_0374 += position[7];
point_0374 += position[4];
VerdictVector centroid = point_1256;
centroid += point_0374;
centroid /= 8.0;
int i;
for( i = 0; i < 8; i++ )
position[i] -= centroid;
// ... Rotate element such that center of side 1-2 and 0-3 define X axis
double DX = point_1256.x() - point_0374.x();
double DY = point_1256.y() - point_0374.y();
double DZ = point_1256.z() - point_0374.z();
double AMAGX = sqrt( DX * DX + DZ * DZ );
double AMAGY = sqrt( DX * DX + DY * DY + DZ * DZ );
double FMAGX = AMAGX == 0.0 ? 1.0 : 0.0;
double FMAGY = AMAGY == 0.0 ? 1.0 : 0.0;
double CZ = DX / ( AMAGX + FMAGX ) + FMAGX;
double SZ = DZ / ( AMAGX + FMAGX );
double CY = AMAGX / ( AMAGY + FMAGY ) + FMAGY;
double SY = DY / ( AMAGY + FMAGY );
double temp;
for( i = 0; i < 8; i++ )
{
temp = CY * CZ * position[i].x() + CY * SZ * position[i].z() + SY * position[i].y();
position[i].y( -SY * CZ * position[i].x() - SY * SZ * position[i].z() + CY * position[i].y() );
position[i].z( -SZ * position[i].x() + CZ * position[i].z() );
position[i].x( temp );
}
// ... Now, rotate about Y
VerdictVector delta = -position[0];
delta -= position[1];
delta += position[2];
delta += position[3];
delta -= position[4];
delta -= position[5];
delta += position[6];
delta += position[7];
DY = delta.y();
DZ = delta.z();
AMAGY = sqrt( DY * DY + DZ * DZ );
FMAGY = AMAGY == 0.0 ? 1.0 : 0.0;
double CX = DY / ( AMAGY + FMAGY ) + FMAGY;
double SX = DZ / ( AMAGY + FMAGY );
for( i = 0; i < 8; i++ )
{
temp = CX * position[i].y() + SX * position[i].z();
position[i].z( -SX * position[i].y() + CX * position[i].z() );
position[i].y( temp );
}
}
double safe_ratio3( const double numerator, const double denominator, const double max_ratio )<--- The function 'safe_ratio3' is never used.
{
// this filter is essential for good running time in practice
double return_value;
const double filter_n = max_ratio * 1.0e-16;
const double filter_d = 1.0e-16;
if( fabs( numerator ) <= filter_n && fabs( denominator ) >= filter_d )
{
return_value = numerator / denominator;
}
else
{
return_value = fabs( numerator ) / max_ratio >= fabs( denominator )
? ( ( numerator >= 0.0 && denominator >= 0.0 ) || ( numerator < 0.0 && denominator < 0.0 )
? max_ratio
: -max_ratio )
: numerator / denominator;
}
return return_value;
}
double safe_ratio( const double numerator, const double denominator )
{
double return_value;
const double filter_n = VERDICT_DBL_MAX;
const double filter_d = VERDICT_DBL_MIN;
if( fabs( numerator ) <= filter_n && fabs( denominator ) >= filter_d )
{
return_value = numerator / denominator;
}
else
{
return_value = VERDICT_DBL_MAX;
}
return return_value;
}
double condition_comp( const VerdictVector& xxi, const VerdictVector& xet, const VerdictVector& xze )
{
double det = xxi % ( xet * xze );
if( det <= VERDICT_DBL_MIN )
{
return VERDICT_DBL_MAX;
}
double term1 = xxi % xxi + xet % xet + xze % xze;<--- Same expression on both sides of '%'. [+]Finding the same expression on both sides of an operator is suspicious and might indicate a cut and paste or logic error. Please examine this code carefully to determine if it is correct.
double term2 =
( ( xxi * xet ) % ( xxi * xet ) ) + ( ( xet * xze ) % ( xet * xze ) ) + ( ( xze * xxi ) % ( xze * xxi ) );<--- Same expression on both sides of '%'. [+]Finding the same expression on both sides of an operator is suspicious and might indicate a cut and paste or logic error. Please examine this code carefully to determine if it is correct.
return sqrt( term1 * term2 ) / det;
}
double oddy_comp( const VerdictVector& xxi, const VerdictVector& xet, const VerdictVector& xze )
{
static const double third = 1.0 / 3.0;
double g11, g12, g13, g22, g23, g33, rt_g;
g11 = xxi % xxi;<--- Same expression on both sides of '%'. [+]Finding the same expression on both sides of an operator is suspicious and might indicate a cut and paste or logic error. Please examine this code carefully to determine if it is correct.
g12 = xxi % xet;
g13 = xxi % xze;
g22 = xet % xet;<--- Same expression on both sides of '%'. [+]Finding the same expression on both sides of an operator is suspicious and might indicate a cut and paste or logic error. Please examine this code carefully to determine if it is correct.
g23 = xet % xze;
g33 = xze % xze;<--- Same expression on both sides of '%'. [+]Finding the same expression on both sides of an operator is suspicious and might indicate a cut and paste or logic error. Please examine this code carefully to determine if it is correct.
rt_g = xxi % ( xet * xze );
double oddy_metric;
if( rt_g > VERDICT_DBL_MIN )
{
double norm_G_squared = g11 * g11 + 2.0 * g12 * g12 + 2.0 * g13 * g13 + g22 * g22 + 2.0 * g23 * g23 + g33 * g33;
double norm_J_squared = g11 + g22 + g33;
oddy_metric = ( norm_G_squared - third * norm_J_squared * norm_J_squared ) / pow( rt_g, 4. * third );
}
else
oddy_metric = VERDICT_DBL_MAX;
return oddy_metric;
}
//! calcualates edge lengths of a hex
double hex_edge_length( int max_min, double coordinates[][3] )
{
double temp[3], edge[12];
int i;
// lengths^2 of edges
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[1][i] - coordinates[0][i];
temp[i] = temp[i] * temp[i];
}
edge[0] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[2][i] - coordinates[1][i];
temp[i] = temp[i] * temp[i];
}
edge[1] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[3][i] - coordinates[2][i];
temp[i] = temp[i] * temp[i];
}
edge[2] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[0][i] - coordinates[3][i];
temp[i] = temp[i] * temp[i];
}
edge[3] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[5][i] - coordinates[4][i];
temp[i] = temp[i] * temp[i];
}
edge[4] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[6][i] - coordinates[5][i];
temp[i] = temp[i] * temp[i];
}
edge[5] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[7][i] - coordinates[6][i];
temp[i] = temp[i] * temp[i];
}
edge[6] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[4][i] - coordinates[7][i];
temp[i] = temp[i] * temp[i];
}
edge[7] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[4][i] - coordinates[0][i];
temp[i] = temp[i] * temp[i];
}
edge[8] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[5][i] - coordinates[1][i];
temp[i] = temp[i] * temp[i];
}
edge[9] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[6][i] - coordinates[2][i];
temp[i] = temp[i] * temp[i];
}
edge[10] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[7][i] - coordinates[3][i];
temp[i] = temp[i] * temp[i];
}
edge[11] = sqrt( temp[0] + temp[1] + temp[2] );
double _edge = edge[0];
if( max_min == 0 )
{
for( i = 1; i < 12; i++ )
_edge = VERDICT_MIN( _edge, edge[i] );
return (double)_edge;
}
else
{
for( i = 1; i < 12; i++ )
_edge = VERDICT_MAX( _edge, edge[i] );
return (double)_edge;
}
}
double diag_length( int max_min, double coordinates[][3] )
{
double temp[3], diag[4];
int i;
// lengths^2 f diag nals
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[6][i] - coordinates[0][i];
temp[i] = temp[i] * temp[i];
}
diag[0] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[4][i] - coordinates[2][i];
temp[i] = temp[i] * temp[i];
}
diag[1] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[7][i] - coordinates[1][i];
temp[i] = temp[i] * temp[i];
}
diag[2] = sqrt( temp[0] + temp[1] + temp[2] );
for( i = 0; i < 3; i++ )
{
temp[i] = coordinates[5][i] - coordinates[3][i];
temp[i] = temp[i] * temp[i];
}
diag[3] = sqrt( temp[0] + temp[1] + temp[2] );
double diagonal = diag[0];
if( max_min == 0 ) // Return min diagonal
{
for( i = 1; i < 4; i++ )
diagonal = VERDICT_MIN( diagonal, diag[i] );
return (double)diagonal;
}
else // Return max diagonal
{
for( i = 1; i < 4; i++ )
diagonal = VERDICT_MAX( diagonal, diag[i] );
return (double)diagonal;
}
}
//! calculates efg values
VerdictVector calc_hex_efg( int efg_index, VerdictVector coordinates[8] )
{
VerdictVector efg;
switch( efg_index )
{
case 1:
efg = coordinates[1];
efg += coordinates[2];
efg += coordinates[5];
efg += coordinates[6];
efg -= coordinates[0];
efg -= coordinates[3];
efg -= coordinates[4];
efg -= coordinates[7];
break;
case 2:
efg = coordinates[2];
efg += coordinates[3];
efg += coordinates[6];
efg += coordinates[7];
efg -= coordinates[0];
efg -= coordinates[1];
efg -= coordinates[4];
efg -= coordinates[5];
break;
case 3:
efg = coordinates[4];
efg += coordinates[5];
efg += coordinates[6];
efg += coordinates[7];
efg -= coordinates[0];
efg -= coordinates[1];
efg -= coordinates[2];
efg -= coordinates[3];
break;
case 12:
efg = coordinates[0];
efg += coordinates[2];
efg += coordinates[4];
efg += coordinates[6];
efg -= coordinates[1];
efg -= coordinates[3];
efg -= coordinates[5];
efg -= coordinates[7];
break;
case 13:
efg = coordinates[0];
efg += coordinates[3];
efg += coordinates[5];
efg += coordinates[6];
efg -= coordinates[1];
efg -= coordinates[2];
efg -= coordinates[4];
efg -= coordinates[7];
break;
case 23:
efg = coordinates[0];
efg += coordinates[1];
efg += coordinates[6];
efg += coordinates[7];
efg -= coordinates[2];
efg -= coordinates[3];
efg -= coordinates[4];
efg -= coordinates[5];
break;
case 123:
efg = coordinates[0];
efg += coordinates[2];
efg += coordinates[5];
efg += coordinates[7];
efg -= coordinates[1];
efg -= coordinates[5];
efg -= coordinates[6];
efg -= coordinates[2];
break;
default:
efg.set( 0, 0, 0 );
}
return efg;
}
/*!
the edge ratio of a hex
NB (P. Pebay 01/23/07):
Hmax / Hmin where Hmax and Hmin are respectively the maximum and the
minimum edge lengths
*/
C_FUNC_DEF double v_hex_edge_ratio( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector edges[12];
make_hex_edges( coordinates, edges );
double a2 = edges[0].length_squared();
double b2 = edges[1].length_squared();
double c2 = edges[2].length_squared();
double d2 = edges[3].length_squared();
double e2 = edges[4].length_squared();
double f2 = edges[5].length_squared();
double g2 = edges[6].length_squared();
double h2 = edges[7].length_squared();
double i2 = edges[8].length_squared();
double j2 = edges[9].length_squared();
double k2 = edges[10].length_squared();
double l2 = edges[11].length_squared();
double mab, mcd, mef, Mab, Mcd, Mef;
double mgh, mij, mkl, Mgh, Mij, Mkl;
if( a2 < b2 )
{
mab = a2;
Mab = b2;
}
else // b2 <= a2
{
mab = b2;
Mab = a2;
}
if( c2 < d2 )
{
mcd = c2;
Mcd = d2;
}
else // d2 <= c2
{
mcd = d2;
Mcd = c2;
}
if( e2 < f2 )
{
mef = e2;
Mef = f2;
}
else // f2 <= e2
{
mef = f2;
Mef = e2;
}
if( g2 < h2 )
{
mgh = g2;
Mgh = h2;
}
else // h2 <= g2
{
mgh = h2;
Mgh = g2;
}
if( i2 < j2 )
{
mij = i2;
Mij = j2;
}
else // j2 <= i2
{
mij = j2;
Mij = i2;
}
if( k2 < l2 )
{
mkl = k2;
Mkl = l2;
}
else // l2 <= k2
{
mkl = l2;
Mkl = k2;
}
double m2;
m2 = mab < mcd ? mab : mcd;
m2 = m2 < mef ? m2 : mef;
m2 = m2 < mgh ? m2 : mgh;
m2 = m2 < mij ? m2 : mij;
m2 = m2 < mkl ? m2 : mkl;
if( m2 < VERDICT_DBL_MIN ) return (double)VERDICT_DBL_MAX;
double M2;
M2 = Mab > Mcd ? Mab : Mcd;
M2 = M2 > Mef ? M2 : Mef;
M2 = M2 > Mgh ? M2 : Mgh;
M2 = M2 > Mij ? M2 : Mij;
M2 = M2 > Mkl ? M2 : Mkl;
m2 = m2 < mef ? m2 : mef;
M2 = Mab > Mcd ? Mab : Mcd;
M2 = M2 > Mef ? M2 : Mef;
double edge_ratio = sqrt( M2 / m2 );
if( edge_ratio > 0 ) return (double)VERDICT_MIN( edge_ratio, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( edge_ratio, -VERDICT_DBL_MAX );
}
/*!
max edge ratio of a hex
Maximum edge length ratio at hex center
*/
C_FUNC_DEF double v_hex_max_edge_ratio( int /*num_nodes*/, double coordinates[][3] )
{
double aspect;
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
double aspect_12, aspect_13, aspect_23;
VerdictVector efg1 = calc_hex_efg( 1, node_pos );
VerdictVector efg2 = calc_hex_efg( 2, node_pos );
VerdictVector efg3 = calc_hex_efg( 3, node_pos );
double mag_efg1 = efg1.length();
double mag_efg2 = efg2.length();
double mag_efg3 = efg3.length();
aspect_12 = safe_ratio( VERDICT_MAX( mag_efg1, mag_efg2 ), VERDICT_MIN( mag_efg1, mag_efg2 ) );
aspect_13 = safe_ratio( VERDICT_MAX( mag_efg1, mag_efg3 ), VERDICT_MIN( mag_efg1, mag_efg3 ) );
aspect_23 = safe_ratio( VERDICT_MAX( mag_efg2, mag_efg3 ), VERDICT_MIN( mag_efg2, mag_efg3 ) );
aspect = VERDICT_MAX( aspect_12, VERDICT_MAX( aspect_13, aspect_23 ) );
if( aspect > 0 ) return (double)VERDICT_MIN( aspect, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( aspect, -VERDICT_DBL_MAX );
}
/*!
skew of a hex
Maximum ||cosA|| where A is the angle between edges at hex center.
*/
C_FUNC_DEF double v_hex_skew( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
double skew_1, skew_2, skew_3;
VerdictVector efg1 = calc_hex_efg( 1, node_pos );
VerdictVector efg2 = calc_hex_efg( 2, node_pos );
VerdictVector efg3 = calc_hex_efg( 3, node_pos );
if( efg1.normalize() <= VERDICT_DBL_MIN ) return VERDICT_DBL_MAX;
if( efg2.normalize() <= VERDICT_DBL_MIN ) return VERDICT_DBL_MAX;
if( efg3.normalize() <= VERDICT_DBL_MIN ) return VERDICT_DBL_MAX;
skew_1 = fabs( efg1 % efg2 );
skew_2 = fabs( efg1 % efg3 );
skew_3 = fabs( efg2 % efg3 );
double skew = ( VERDICT_MAX( skew_1, VERDICT_MAX( skew_2, skew_3 ) ) );
if( skew > 0 ) return (double)VERDICT_MIN( skew, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( skew, -VERDICT_DBL_MAX );
}
/*!
taper of a hex
Maximum ratio of lengths derived from opposite edges.
*/
C_FUNC_DEF double v_hex_taper( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
VerdictVector efg1 = calc_hex_efg( 1, node_pos );
VerdictVector efg2 = calc_hex_efg( 2, node_pos );
VerdictVector efg3 = calc_hex_efg( 3, node_pos );
VerdictVector efg12 = calc_hex_efg( 12, node_pos );
VerdictVector efg13 = calc_hex_efg( 13, node_pos );
VerdictVector efg23 = calc_hex_efg( 23, node_pos );
double taper_1 = fabs( safe_ratio( efg12.length(), VERDICT_MIN( efg1.length(), efg2.length() ) ) );
double taper_2 = fabs( safe_ratio( efg13.length(), VERDICT_MIN( efg1.length(), efg3.length() ) ) );
double taper_3 = fabs( safe_ratio( efg23.length(), VERDICT_MIN( efg2.length(), efg3.length() ) ) );
double taper = (double)VERDICT_MAX( taper_1, VERDICT_MAX( taper_2, taper_3 ) );
if( taper > 0 ) return (double)VERDICT_MIN( taper, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( taper, -VERDICT_DBL_MAX );
}
/*!
volume of a hex
Jacobian at hex center
*/
C_FUNC_DEF double v_hex_volume( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
VerdictVector efg1 = calc_hex_efg( 1, node_pos );
VerdictVector efg2 = calc_hex_efg( 2, node_pos );
VerdictVector efg3 = calc_hex_efg( 3, node_pos );
double volume;
volume = (double)( efg1 % ( efg2 * efg3 ) ) / 64.0;
if( volume > 0 ) return (double)VERDICT_MIN( volume, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( volume, -VERDICT_DBL_MAX );
}
/*!
stretch of a hex
sqrt(3) * minimum edge length / maximum diagonal length
*/
C_FUNC_DEF double v_hex_stretch( int /*num_nodes*/, double coordinates[][3] )
{
static const double HEX_STRETCH_SCALE_FACTOR = sqrt( 3.0 );
double min_edge = hex_edge_length( 0, coordinates );
double max_diag = diag_length( 1, coordinates );
double stretch = HEX_STRETCH_SCALE_FACTOR * safe_ratio( min_edge, max_diag );
if( stretch > 0 ) return (double)VERDICT_MIN( stretch, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( stretch, -VERDICT_DBL_MAX );
}
/*!
diagonal ratio of a hex
Minimum diagonal length / maximum diagonal length
*/
C_FUNC_DEF double v_hex_diagonal( int /*num_nodes*/, double coordinates[][3] )
{
double min_diag = diag_length( 0, coordinates );
double max_diag = diag_length( 1, coordinates );
double diagonal = safe_ratio( min_diag, max_diag );
if( diagonal > 0 ) return (double)VERDICT_MIN( diagonal, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( diagonal, -VERDICT_DBL_MAX );
}
#define SQR( x ) ( ( x ) * ( x ) )
/*!
dimension of a hex
Pronto-specific characteristic length for stable time step calculation.
Char_length = Volume / 2 grad Volume
*/
C_FUNC_DEF double v_hex_dimension( int /*num_nodes*/, double coordinates[][3] )
{
double gradop[9][4];
double x1 = coordinates[0][0];
double x2 = coordinates[1][0];
double x3 = coordinates[2][0];
double x4 = coordinates[3][0];
double x5 = coordinates[4][0];
double x6 = coordinates[5][0];
double x7 = coordinates[6][0];
double x8 = coordinates[7][0];
double y1 = coordinates[0][1];
double y2 = coordinates[1][1];
double y3 = coordinates[2][1];
double y4 = coordinates[3][1];
double y5 = coordinates[4][1];
double y6 = coordinates[5][1];
double y7 = coordinates[6][1];
double y8 = coordinates[7][1];
double z1 = coordinates[0][2];
double z2 = coordinates[1][2];
double z3 = coordinates[2][2];
double z4 = coordinates[3][2];
double z5 = coordinates[4][2];
double z6 = coordinates[5][2];
double z7 = coordinates[6][2];
double z8 = coordinates[7][2];
double z24 = z2 - z4;
double z52 = z5 - z2;
double z45 = z4 - z5;
gradop[1][1] =
( y2 * ( z6 - z3 - z45 ) + y3 * z24 + y4 * ( z3 - z8 - z52 ) + y5 * ( z8 - z6 - z24 ) + y6 * z52 + y8 * z45 ) /
12.0;
double z31 = z3 - z1;
double z63 = z6 - z3;
double z16 = z1 - z6;
gradop[2][1] =
( y3 * ( z7 - z4 - z16 ) + y4 * z31 + y1 * ( z4 - z5 - z63 ) + y6 * ( z5 - z7 - z31 ) + y7 * z63 + y5 * z16 ) /
12.0;
double z42 = z4 - z2;
double z74 = z7 - z4;
double z27 = z2 - z7;
gradop[3][1] =
( y4 * ( z8 - z1 - z27 ) + y1 * z42 + y2 * ( z1 - z6 - z74 ) + y7 * ( z6 - z8 - z42 ) + y8 * z74 + y6 * z27 ) /
12.0;
double z13 = z1 - z3;
double z81 = z8 - z1;
double z38 = z3 - z8;
gradop[4][1] =
( y1 * ( z5 - z2 - z38 ) + y2 * z13 + y3 * ( z2 - z7 - z81 ) + y8 * ( z7 - z5 - z13 ) + y5 * z81 + y7 * z38 ) /
12.0;
double z86 = z8 - z6;
double z18 = z1 - z8;
double z61 = z6 - z1;
gradop[5][1] =
( y8 * ( z4 - z7 - z61 ) + y7 * z86 + y6 * ( z7 - z2 - z18 ) + y1 * ( z2 - z4 - z86 ) + y4 * z18 + y2 * z61 ) /
12.0;
double z57 = z5 - z7;
double z25 = z2 - z5;
double z72 = z7 - z2;
gradop[6][1] =
( y5 * ( z1 - z8 - z72 ) + y8 * z57 + y7 * ( z8 - z3 - z25 ) + y2 * ( z3 - z1 - z57 ) + y1 * z25 + y3 * z72 ) /
12.0;
double z68 = z6 - z8;
double z36 = z3 - z6;
double z83 = z8 - z3;
gradop[7][1] =
( y6 * ( z2 - z5 - z83 ) + y5 * z68 + y8 * ( z5 - z4 - z36 ) + y3 * ( z4 - z2 - z68 ) + y2 * z36 + y4 * z83 ) /
12.0;
double z75 = z7 - z5;
double z47 = z4 - z7;
double z54 = z5 - z4;
gradop[8][1] =
( y7 * ( z3 - z6 - z54 ) + y6 * z75 + y5 * ( z6 - z1 - z47 ) + y4 * ( z1 - z3 - z75 ) + y3 * z47 + y1 * z54 ) /
12.0;
double x24 = x2 - x4;
double x52 = x5 - x2;
double x45 = x4 - x5;
gradop[1][2] =
( z2 * ( x6 - x3 - x45 ) + z3 * x24 + z4 * ( x3 - x8 - x52 ) + z5 * ( x8 - x6 - x24 ) + z6 * x52 + z8 * x45 ) /
12.0;
double x31 = x3 - x1;
double x63 = x6 - x3;
double x16 = x1 - x6;
gradop[2][2] =
( z3 * ( x7 - x4 - x16 ) + z4 * x31 + z1 * ( x4 - x5 - x63 ) + z6 * ( x5 - x7 - x31 ) + z7 * x63 + z5 * x16 ) /
12.0;
double x42 = x4 - x2;
double x74 = x7 - x4;
double x27 = x2 - x7;
gradop[3][2] =
( z4 * ( x8 - x1 - x27 ) + z1 * x42 + z2 * ( x1 - x6 - x74 ) + z7 * ( x6 - x8 - x42 ) + z8 * x74 + z6 * x27 ) /
12.0;
double x13 = x1 - x3;
double x81 = x8 - x1;
double x38 = x3 - x8;
gradop[4][2] =
( z1 * ( x5 - x2 - x38 ) + z2 * x13 + z3 * ( x2 - x7 - x81 ) + z8 * ( x7 - x5 - x13 ) + z5 * x81 + z7 * x38 ) /
12.0;
double x86 = x8 - x6;
double x18 = x1 - x8;
double x61 = x6 - x1;
gradop[5][2] =
( z8 * ( x4 - x7 - x61 ) + z7 * x86 + z6 * ( x7 - x2 - x18 ) + z1 * ( x2 - x4 - x86 ) + z4 * x18 + z2 * x61 ) /
12.0;
double x57 = x5 - x7;
double x25 = x2 - x5;
double x72 = x7 - x2;
gradop[6][2] =
( z5 * ( x1 - x8 - x72 ) + z8 * x57 + z7 * ( x8 - x3 - x25 ) + z2 * ( x3 - x1 - x57 ) + z1 * x25 + z3 * x72 ) /
12.0;
double x68 = x6 - x8;
double x36 = x3 - x6;
double x83 = x8 - x3;
gradop[7][2] =
( z6 * ( x2 - x5 - x83 ) + z5 * x68 + z8 * ( x5 - x4 - x36 ) + z3 * ( x4 - x2 - x68 ) + z2 * x36 + z4 * x83 ) /
12.0;
double x75 = x7 - x5;
double x47 = x4 - x7;
double x54 = x5 - x4;
gradop[8][2] =
( z7 * ( x3 - x6 - x54 ) + z6 * x75 + z5 * ( x6 - x1 - x47 ) + z4 * ( x1 - x3 - x75 ) + z3 * x47 + z1 * x54 ) /
12.0;
double y24 = y2 - y4;
double y52 = y5 - y2;
double y45 = y4 - y5;
gradop[1][3] =
( x2 * ( y6 - y3 - y45 ) + x3 * y24 + x4 * ( y3 - y8 - y52 ) + x5 * ( y8 - y6 - y24 ) + x6 * y52 + x8 * y45 ) /
12.0;
double y31 = y3 - y1;
double y63 = y6 - y3;
double y16 = y1 - y6;
gradop[2][3] =
( x3 * ( y7 - y4 - y16 ) + x4 * y31 + x1 * ( y4 - y5 - y63 ) + x6 * ( y5 - y7 - y31 ) + x7 * y63 + x5 * y16 ) /
12.0;
double y42 = y4 - y2;
double y74 = y7 - y4;
double y27 = y2 - y7;
gradop[3][3] =
( x4 * ( y8 - y1 - y27 ) + x1 * y42 + x2 * ( y1 - y6 - y74 ) + x7 * ( y6 - y8 - y42 ) + x8 * y74 + x6 * y27 ) /
12.0;
double y13 = y1 - y3;
double y81 = y8 - y1;
double y38 = y3 - y8;
gradop[4][3] =
( x1 * ( y5 - y2 - y38 ) + x2 * y13 + x3 * ( y2 - y7 - y81 ) + x8 * ( y7 - y5 - y13 ) + x5 * y81 + x7 * y38 ) /
12.0;
double y86 = y8 - y6;
double y18 = y1 - y8;
double y61 = y6 - y1;
gradop[5][3] =
( x8 * ( y4 - y7 - y61 ) + x7 * y86 + x6 * ( y7 - y2 - y18 ) + x1 * ( y2 - y4 - y86 ) + x4 * y18 + x2 * y61 ) /
12.0;
double y57 = y5 - y7;
double y25 = y2 - y5;
double y72 = y7 - y2;
gradop[6][3] =
( x5 * ( y1 - y8 - y72 ) + x8 * y57 + x7 * ( y8 - y3 - y25 ) + x2 * ( y3 - y1 - y57 ) + x1 * y25 + x3 * y72 ) /
12.0;
double y68 = y6 - y8;
double y36 = y3 - y6;
double y83 = y8 - y3;
gradop[7][3] =
( x6 * ( y2 - y5 - y83 ) + x5 * y68 + x8 * ( y5 - y4 - y36 ) + x3 * ( y4 - y2 - y68 ) + x2 * y36 + x4 * y83 ) /
12.0;
double y75 = y7 - y5;
double y47 = y4 - y7;
double y54 = y5 - y4;
gradop[8][3] =
( x7 * ( y3 - y6 - y54 ) + x6 * y75 + x5 * ( y6 - y1 - y47 ) + x4 * ( y1 - y3 - y75 ) + x3 * y47 + x1 * y54 ) /
12.0;
// calculate element volume and characteristic element aspect ratio
// (used in time step and hourglass control) -
double volume = coordinates[0][0] * gradop[1][1] + coordinates[1][0] * gradop[2][1] +
coordinates[2][0] * gradop[3][1] + coordinates[3][0] * gradop[4][1] +
coordinates[4][0] * gradop[5][1] + coordinates[5][0] * gradop[6][1] +
coordinates[6][0] * gradop[7][1] + coordinates[7][0] * gradop[8][1];
double aspect =
.5 * SQR( volume ) /
( SQR( gradop[1][1] ) + SQR( gradop[2][1] ) + SQR( gradop[3][1] ) + SQR( gradop[4][1] ) + SQR( gradop[5][1] ) +
SQR( gradop[6][1] ) + SQR( gradop[7][1] ) + SQR( gradop[8][1] ) + SQR( gradop[1][2] ) + SQR( gradop[2][2] ) +
SQR( gradop[3][2] ) + SQR( gradop[4][2] ) + SQR( gradop[5][2] ) + SQR( gradop[6][2] ) + SQR( gradop[7][2] ) +
SQR( gradop[8][2] ) + SQR( gradop[1][3] ) + SQR( gradop[2][3] ) + SQR( gradop[3][3] ) + SQR( gradop[4][3] ) +
SQR( gradop[5][3] ) + SQR( gradop[6][3] ) + SQR( gradop[7][3] ) + SQR( gradop[8][3] ) );
return (double)sqrt( aspect );
}
/*!
oddy of a hex
General distortion measure based on left Cauchy-Green Tensor
*/
C_FUNC_DEF double v_hex_oddy( int /*num_nodes*/, double coordinates[][3] )
{
double oddy = 0.0, current_oddy;
VerdictVector xxi, xet, xze;
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
xxi = calc_hex_efg( 1, node_pos );
xet = calc_hex_efg( 2, node_pos );
xze = calc_hex_efg( 3, node_pos );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[1][0] - coordinates[0][0], coordinates[1][1] - coordinates[0][1],
coordinates[1][2] - coordinates[0][2] );
xet.set( coordinates[3][0] - coordinates[0][0], coordinates[3][1] - coordinates[0][1],
coordinates[3][2] - coordinates[0][2] );
xze.set( coordinates[4][0] - coordinates[0][0], coordinates[4][1] - coordinates[0][1],
coordinates[4][2] - coordinates[0][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[2][0] - coordinates[1][0], coordinates[2][1] - coordinates[1][1],
coordinates[2][2] - coordinates[1][2] );
xet.set( coordinates[0][0] - coordinates[1][0], coordinates[0][1] - coordinates[1][1],
coordinates[0][2] - coordinates[1][2] );
xze.set( coordinates[5][0] - coordinates[1][0], coordinates[5][1] - coordinates[1][1],
coordinates[5][2] - coordinates[1][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[3][0] - coordinates[2][0], coordinates[3][1] - coordinates[2][1],
coordinates[3][2] - coordinates[2][2] );
xet.set( coordinates[1][0] - coordinates[2][0], coordinates[1][1] - coordinates[2][1],
coordinates[1][2] - coordinates[2][2] );
xze.set( coordinates[6][0] - coordinates[2][0], coordinates[6][1] - coordinates[2][1],
coordinates[6][2] - coordinates[2][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[0][0] - coordinates[3][0], coordinates[0][1] - coordinates[3][1],
coordinates[0][2] - coordinates[3][2] );
xet.set( coordinates[2][0] - coordinates[3][0], coordinates[2][1] - coordinates[3][1],
coordinates[2][2] - coordinates[3][2] );
xze.set( coordinates[7][0] - coordinates[3][0], coordinates[7][1] - coordinates[3][1],
coordinates[7][2] - coordinates[3][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[7][0] - coordinates[4][0], coordinates[7][1] - coordinates[4][1],
coordinates[7][2] - coordinates[4][2] );
xet.set( coordinates[5][0] - coordinates[4][0], coordinates[5][1] - coordinates[4][1],
coordinates[5][2] - coordinates[4][2] );
xze.set( coordinates[0][0] - coordinates[4][0], coordinates[0][1] - coordinates[4][1],
coordinates[0][2] - coordinates[4][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[4][0] - coordinates[5][0], coordinates[4][1] - coordinates[5][1],
coordinates[4][2] - coordinates[5][2] );
xet.set( coordinates[6][0] - coordinates[5][0], coordinates[6][1] - coordinates[5][1],
coordinates[6][2] - coordinates[5][2] );
xze.set( coordinates[1][0] - coordinates[5][0], coordinates[1][1] - coordinates[5][1],
coordinates[1][2] - coordinates[5][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[5][0] - coordinates[6][0], coordinates[5][1] - coordinates[6][1],
coordinates[5][2] - coordinates[6][2] );
xet.set( coordinates[7][0] - coordinates[6][0], coordinates[7][1] - coordinates[6][1],
coordinates[7][2] - coordinates[6][2] );
xze.set( coordinates[2][0] - coordinates[6][0], coordinates[2][1] - coordinates[6][1],
coordinates[2][2] - coordinates[6][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
xxi.set( coordinates[6][0] - coordinates[7][0], coordinates[6][1] - coordinates[7][1],
coordinates[6][2] - coordinates[7][2] );
xet.set( coordinates[4][0] - coordinates[7][0], coordinates[4][1] - coordinates[7][1],
coordinates[4][2] - coordinates[7][2] );
xze.set( coordinates[3][0] - coordinates[7][0], coordinates[3][1] - coordinates[7][1],
coordinates[3][2] - coordinates[7][2] );
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
if( oddy > 0 ) return (double)VERDICT_MIN( oddy, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( oddy, -VERDICT_DBL_MAX );
}
/*!
the average Frobenius aspect of a hex
NB (P. Pebay 01/20/07):
this metric is calculated by averaging the 8 Frobenius aspects at
each corner of the hex, when the reference corner is right isosceles.
*/
C_FUNC_DEF double v_hex_med_aspect_frobenius( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
double med_aspect_frobenius = 0.;
VerdictVector xxi, xet, xze;
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(1,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(0,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(0,0,1):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
// J(1,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
med_aspect_frobenius += condition_comp( xxi, xet, xze );
med_aspect_frobenius /= 24.;
if( med_aspect_frobenius > 0 ) return (double)VERDICT_MIN( med_aspect_frobenius, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( med_aspect_frobenius, -VERDICT_DBL_MAX );
}
/*!
maximum Frobenius condition number of a hex
Maximum Frobenius condition number of the Jacobian matrix at 8 corners
NB (P. Pebay 01/25/07):
this metric is calculated by taking the maximum of the 8 Frobenius aspects at
each corner of the hex, when the reference corner is right isosceles.
*/
C_FUNC_DEF double v_hex_max_aspect_frobenius( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
double condition = 0.0, current_condition;
VerdictVector xxi, xet, xze;
xxi = calc_hex_efg( 1, node_pos );
xet = calc_hex_efg( 2, node_pos );
xze = calc_hex_efg( 3, node_pos );
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(1,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(0,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(0,0,1):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
// J(1,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
condition /= 3.0;
if( condition > 0 ) return (double)VERDICT_MIN( condition, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( condition, -VERDICT_DBL_MAX );
}
/*!
The maximum Frobenius condition of a hex, a.k.a. condition
NB (P. Pebay 01/25/07):
this method is maintained for backwards compatibility only.
It will become deprecated at some point.
*/
C_FUNC_DEF double v_hex_condition( int /*num_nodes*/, double coordinates[][3] )
{
return v_hex_max_aspect_frobenius( 8, coordinates );
}
/*!
jacobian of a hex
Minimum pointwise volume of local map at 8 corners & center of hex
*/
C_FUNC_DEF double v_hex_jacobian( int /*num_nodes*/, double coordinates[][3] )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
double jacobian = VERDICT_DBL_MAX;
double current_jacobian;
VerdictVector xxi, xet, xze;
xxi = calc_hex_efg( 1, node_pos );
xet = calc_hex_efg( 2, node_pos );
xze = calc_hex_efg( 3, node_pos );
current_jacobian = xxi % ( xet * xze ) / 64.0;
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(1,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(0,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(0,0,1):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
// J(0,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
current_jacobian = xxi % ( xet * xze );
if( current_jacobian < jacobian )
{
jacobian = current_jacobian;
}
if( jacobian > 0 ) return (double)VERDICT_MIN( jacobian, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( jacobian, -VERDICT_DBL_MAX );
}
/*!
scaled jacobian of a hex
Minimum Jacobian divided by the lengths of the 3 edge vectors
*/
C_FUNC_DEF double v_hex_scaled_jacobian( int /*num_nodes*/, double coordinates[][3] )
{
double jacobi, min_norm_jac = VERDICT_DBL_MAX;
// double min_jacobi = VERDICT_DBL_MAX;
double temp_norm_jac, lengths;
double len1_sq, len2_sq, len3_sq;
VerdictVector xxi, xet, xze;
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
xxi = calc_hex_efg( 1, node_pos );
xet = calc_hex_efg( 2, node_pos );
xze = calc_hex_efg( 3, node_pos );
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(1,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(0,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(0,0,1):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac )
min_norm_jac = temp_norm_jac;
else
temp_norm_jac = jacobi;<--- temp_norm_jac is assigned
// J(0,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
jacobi = xxi % ( xet * xze );
// if( jacobi < min_jacobi ) { min_jacobi = jacobi; }
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN )
return (double)VERDICT_DBL_MAX;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
temp_norm_jac = jacobi / lengths;<--- temp_norm_jac is overwritten
if( temp_norm_jac < min_norm_jac ) min_norm_jac = temp_norm_jac;
// else
// temp_norm_jac = jacobi;
if( min_norm_jac > 0 ) return (double)VERDICT_MIN( min_norm_jac, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( min_norm_jac, -VERDICT_DBL_MAX );
}
/*!
shear of a hex
3/Condition number of Jacobian Skew matrix
*/
C_FUNC_DEF double v_hex_shear( int /*num_nodes*/, double coordinates[][3] )
{
double shear;
double min_shear = 1.0;
VerdictVector xxi, xet, xze;
double det, len1_sq, len2_sq, len3_sq, lengths;
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(1,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(0,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(0,0,1):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
// J(0,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
len1_sq = xxi.length_squared();
len2_sq = xet.length_squared();
len3_sq = xze.length_squared();
if( len1_sq <= VERDICT_DBL_MIN || len2_sq <= VERDICT_DBL_MIN || len3_sq <= VERDICT_DBL_MIN ) return 0;
lengths = sqrt( len1_sq * len2_sq * len3_sq );
det = xxi % ( xet * xze );
if( det < VERDICT_DBL_MIN )
{
return 0;
}
shear = det / lengths;
min_shear = VERDICT_MIN( shear, min_shear );
if( min_shear <= VERDICT_DBL_MIN ) min_shear = 0;
if( min_shear > 0 ) return (double)VERDICT_MIN( min_shear, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( min_shear, -VERDICT_DBL_MAX );
}
/*!
shape of a hex
3/Condition number of weighted Jacobian matrix
*/
C_FUNC_DEF double v_hex_shape( int /*num_nodes*/, double coordinates[][3] )
{
double det, shape;
double min_shape = 1.0;
static const double two_thirds = 2.0 / 3.0;
VerdictVector xxi, xet, xze;
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(1,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(0,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(0,0,1):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
// J(1,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
det = xxi % ( xet * xze );
if( det > VERDICT_DBL_MIN )
shape = 3 * pow( det, two_thirds ) / ( xxi % xxi + xet % xet + xze % xze );
else
return 0;
if( shape < min_shape )
{
min_shape = shape;
}
if( min_shape <= VERDICT_DBL_MIN ) min_shape = 0;
if( min_shape > 0 ) return (double)VERDICT_MIN( min_shape, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( min_shape, -VERDICT_DBL_MAX );
}
/*!
relative size of a hex
Min( J, 1/J ), where J is determinant of weighted Jacobian matrix
*/
C_FUNC_DEF double v_hex_relative_size_squared( int /*num_nodes*/, double coordinates[][3] )
{
double size = 0;
double tau;<--- The scope of the variable 'tau' can be reduced. [+]The scope of the variable 'tau' can be reduced. Warning: Be careful when fixing this message, especially when there are inner loops. Here is an example where cppcheck will write that the scope for 'i' can be reduced:
void f(int x)
{
int i = 0;
if (x) {
// it's safe to move 'int i = 0;' here
for (int n = 0; n < 10; ++n) {
// it is possible but not safe to move 'int i = 0;' here
do_something(&i);
}
}
}
When you see this message it is always safe to reduce the variable scope 1 level.
VerdictVector xxi, xet, xze;
double det, det_sum = 0;
v_hex_get_weight( xxi, xet, xze );
// This is the average relative size
double detw = xxi % ( xet * xze );
if( detw < VERDICT_DBL_MIN ) return 0;
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
// J(0,0,0):
xxi = node_pos[1] - node_pos[0];
xet = node_pos[3] - node_pos[0];
xze = node_pos[4] - node_pos[0];
det = xxi % ( xet * xze );
det_sum += det;
// J(1,0,0):
xxi = node_pos[2] - node_pos[1];
xet = node_pos[0] - node_pos[1];
xze = node_pos[5] - node_pos[1];
det = xxi % ( xet * xze );
det_sum += det;
// J(0,1,0):
xxi = node_pos[3] - node_pos[2];
xet = node_pos[1] - node_pos[2];
xze = node_pos[6] - node_pos[2];
det = xxi % ( xet * xze );
det_sum += det;
// J(1,1,0):
xxi = node_pos[0] - node_pos[3];
xet = node_pos[2] - node_pos[3];
xze = node_pos[7] - node_pos[3];
det = xxi % ( xet * xze );
det_sum += det;
// J(0,1,0):
xxi = node_pos[7] - node_pos[4];
xet = node_pos[5] - node_pos[4];
xze = node_pos[0] - node_pos[4];
det = xxi % ( xet * xze );
det_sum += det;
// J(1,0,1):
xxi = node_pos[4] - node_pos[5];
xet = node_pos[6] - node_pos[5];
xze = node_pos[1] - node_pos[5];
det = xxi % ( xet * xze );
det_sum += det;
// J(1,1,1):
xxi = node_pos[5] - node_pos[6];
xet = node_pos[7] - node_pos[6];
xze = node_pos[2] - node_pos[6];
det = xxi % ( xet * xze );
det_sum += det;
// J(1,1,1):
xxi = node_pos[6] - node_pos[7];
xet = node_pos[4] - node_pos[7];
xze = node_pos[3] - node_pos[7];
det = xxi % ( xet * xze );
det_sum += det;
if( det_sum > VERDICT_DBL_MIN )
{
tau = det_sum / ( 8 * detw );
tau = VERDICT_MIN( tau, 1.0 / tau );
size = tau * tau;
}
if( size > 0 ) return (double)VERDICT_MIN( size, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( size, -VERDICT_DBL_MAX );
}
/*!
shape and size of a hex
Product of Shape and Relative Size
*/
C_FUNC_DEF double v_hex_shape_and_size( int num_nodes, double coordinates[][3] )
{
double size = v_hex_relative_size_squared( num_nodes, coordinates );
double shape = v_hex_shape( num_nodes, coordinates );
double shape_size = size * shape;
if( shape_size > 0 ) return (double)VERDICT_MIN( shape_size, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( shape_size, -VERDICT_DBL_MAX );
}
/*!
shear and size of a hex
Product of Shear and Relative Size
*/
C_FUNC_DEF double v_hex_shear_and_size( int num_nodes, double coordinates[][3] )
{
double size = v_hex_relative_size_squared( num_nodes, coordinates );
double shear = v_hex_shear( num_nodes, coordinates );
double shear_size = shear * size;
if( shear_size > 0 ) return (double)VERDICT_MIN( shear_size, VERDICT_DBL_MAX );
return (double)VERDICT_MAX( shear_size, -VERDICT_DBL_MAX );
}
/*!
distortion of a hex
*/
C_FUNC_DEF double v_hex_distortion( int num_nodes, double coordinates[][3] )
{
// use 2x2 gauss points for linear hex and 3x3 for 2nd order hex
int number_of_gauss_points = 0;
if( num_nodes == 8 )
// 2x2 quadrature rule
number_of_gauss_points = 2;
else if( num_nodes == 20 )
// 3x3 quadrature rule
number_of_gauss_points = 3;
int number_dimension = 3;
int total_number_of_gauss_points = number_of_gauss_points * number_of_gauss_points * number_of_gauss_points;
double distortion = VERDICT_DBL_MAX;
// maxTotalNumberGaussPoints =27, maxNumberNodes = 20
// they are defined in GaussIntegration.hpp
// This is used to make these arrays static.
// I tried dynamically allocated arrays but the new and delete
// was very expensive
double shape_function[maxTotalNumberGaussPoints][maxNumberNodes];
double dndy1[maxTotalNumberGaussPoints][maxNumberNodes];
double dndy2[maxTotalNumberGaussPoints][maxNumberNodes];
double dndy3[maxTotalNumberGaussPoints][maxNumberNodes];
double weight[maxTotalNumberGaussPoints];
// create an object of GaussIntegration
GaussIntegration::initialize( number_of_gauss_points, num_nodes, number_dimension );
GaussIntegration::calculate_shape_function_3d_hex();
GaussIntegration::get_shape_func( shape_function[0], dndy1[0], dndy2[0], dndy3[0], weight );
VerdictVector xxi, xet, xze, xin;
double jacobian, minimum_jacobian;
double element_volume = 0.0;
minimum_jacobian = VERDICT_DBL_MAX;
// calculate element volume
int ife, ja;
for( ife = 0; ife < total_number_of_gauss_points; ife++ )
{
xxi.set( 0.0, 0.0, 0.0 );
xet.set( 0.0, 0.0, 0.0 );
xze.set( 0.0, 0.0, 0.0 );
for( ja = 0; ja < num_nodes; ja++ )
{
xin.set( coordinates[ja][0], coordinates[ja][1], coordinates[ja][2] );
xxi += dndy1[ife][ja] * xin;
xet += dndy2[ife][ja] * xin;
xze += dndy3[ife][ja] * xin;
}
jacobian = xxi % ( xet * xze );
if( minimum_jacobian > jacobian ) minimum_jacobian = jacobian;
element_volume += weight[ife] * jacobian;
}
// loop through all nodes
double dndy1_at_node[maxNumberNodes][maxNumberNodes];
double dndy2_at_node[maxNumberNodes][maxNumberNodes];
double dndy3_at_node[maxNumberNodes][maxNumberNodes];
GaussIntegration::calculate_derivative_at_nodes_3d( dndy1_at_node, dndy2_at_node, dndy3_at_node );
int node_id;
for( node_id = 0; node_id < num_nodes; node_id++ )
{
xxi.set( 0.0, 0.0, 0.0 );
xet.set( 0.0, 0.0, 0.0 );
xze.set( 0.0, 0.0, 0.0 );
for( ja = 0; ja < num_nodes; ja++ )
{
xin.set( coordinates[ja][0], coordinates[ja][1], coordinates[ja][2] );
xxi += dndy1_at_node[node_id][ja] * xin;
xet += dndy2_at_node[node_id][ja] * xin;
xze += dndy3_at_node[node_id][ja] * xin;
}
jacobian = xxi % ( xet * xze );
if( minimum_jacobian > jacobian ) minimum_jacobian = jacobian;
}
distortion = minimum_jacobian / element_volume * 8.;
return (double)distortion;
}
/*
C_FUNC_DEF double hex_jac_normjac_oddy_cond( int choices[],
double coordinates[][3],
double answers[4] )
{
//Define variables
int i;
double xxi[3], xet[3], xze[3];
double norm_jacobian = 0.0, current_norm_jac = 0.0;
double jacobian = 0.0, current_jacobian = 0.0;
double oddy = 0.0, current_oddy = 0.0;
double condition = 0.0, current_condition = 0.0;
for( i=0; i<3; i++)
xxi[i] = calc_hex_efg( 2, i, coordinates );
for( i=0; i<3; i++)
xet[i] = calc_hex_efg( 3, i, coordinates );
for( i=0; i<3; i++)
xze[i] = calc_hex_efg( 6, i, coordinates );
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
current_jacobian /= 64.0;
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[1][i] - coordinates[0][i];
xet[i] = coordinates[3][i] - coordinates[0][i];
xze[i] = coordinates[4][i] - coordinates[0][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[1][i] - coordinates[0][i];
xet[i] = coordinates[2][i] - coordinates[1][i];
xze[i] = coordinates[5][i] - coordinates[1][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[2][i] - coordinates[3][i];
xet[i] = coordinates[3][i] - coordinates[0][i];
xze[i] = coordinates[7][i] - coordinates[3][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[5][i] - coordinates[4][i];
xet[i] = coordinates[7][i] - coordinates[4][i];
xze[i] = coordinates[4][i] - coordinates[0][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[2][i] - coordinates[3][i];
xet[i] = coordinates[2][i] - coordinates[1][i];
xze[i] = coordinates[6][i] - coordinates[2][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[5][i] - coordinates[4][i];
xet[i] = coordinates[6][i] - coordinates[5][i];
xze[i] = coordinates[5][i] - coordinates[1][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[6][i] - coordinates[7][i];
xet[i] = coordinates[7][i] - coordinates[4][i];
xze[i] = coordinates[7][i] - coordinates[3][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
}
for( i=0; i<3; i++)
{
xxi[i] = coordinates[6][i] - coordinates[7][i];
xet[i] = coordinates[6][i] - coordinates[5][i];
xze[i] = coordinates[6][i] - coordinates[2][i];
}
// norm jacobian and jacobian
if( choices[0] || choices[1] )
{
current_jacobian = determinant( xxi, xet, xze );
current_norm_jac = normalize_jacobian( current_jacobian,
xxi, xet, xze );
if (current_norm_jac < norm_jacobian) { norm_jacobian = current_norm_jac; }
if (current_jacobian < jacobian) { jacobian = current_jacobian; }
}
// oddy
if( choices[2] )
{
current_oddy = oddy_comp( xxi, xet, xze );
if ( current_oddy > oddy ) { oddy = current_oddy; }
}
// condition
if( choices[3] )
{
current_condition = condition_comp( xxi, xet, xze );
if ( current_condition > condition ) { condition = current_condition; }
condition /= 3.0;
}
answers[0] = jacobian;
answers[1] = norm_jacobian;
answers[2] = oddy;
answers[3] = condition;
return 1.0;
}
*/
/*!
multiple quality metrics of a hex
*/
C_FUNC_DEF void v_hex_quality( int num_nodes,
double coordinates[][3],
unsigned int metrics_request_flag,
HexMetricVals* metric_vals )
{
memset( metric_vals, 0, sizeof( HexMetricVals ) );
// max edge ratio, skew, taper, and volume
if( metrics_request_flag & ( V_HEX_MAX_EDGE_RATIO | V_HEX_SKEW | V_HEX_TAPER ) )
{
VerdictVector node_pos[8];
make_hex_nodes( coordinates, node_pos );
VerdictVector efg1, efg2, efg3;
efg1 = calc_hex_efg( 1, node_pos );
efg2 = calc_hex_efg( 2, node_pos );
efg3 = calc_hex_efg( 3, node_pos );
if( metrics_request_flag & V_HEX_MAX_EDGE_RATIO )
{
double max_edge_ratio_12, max_edge_ratio_13, max_edge_ratio_23;
double mag_efg1 = efg1.length();
double mag_efg2 = efg2.length();
double mag_efg3 = efg3.length();
max_edge_ratio_12 = safe_ratio( VERDICT_MAX( mag_efg1, mag_efg2 ), VERDICT_MIN( mag_efg1, mag_efg2 ) );
max_edge_ratio_13 = safe_ratio( VERDICT_MAX( mag_efg1, mag_efg3 ), VERDICT_MIN( mag_efg1, mag_efg3 ) );
max_edge_ratio_23 = safe_ratio( VERDICT_MAX( mag_efg2, mag_efg3 ), VERDICT_MIN( mag_efg2, mag_efg3 ) );
metric_vals->max_edge_ratio =
(double)VERDICT_MAX( max_edge_ratio_12, VERDICT_MAX( max_edge_ratio_13, max_edge_ratio_23 ) );
}
if( metrics_request_flag & V_HEX_SKEW )
{
VerdictVector vec1 = efg1;
VerdictVector vec2 = efg2;
VerdictVector vec3 = efg3;
if( vec1.normalize() <= VERDICT_DBL_MIN || vec2.normalize() <= VERDICT_DBL_MIN ||
vec3.normalize() <= VERDICT_DBL_MIN )
{
metric_vals->skew = (double)VERDICT_DBL_MAX;
}
else
{
double skewx = fabs( vec1 % vec2 );
double skewy = fabs( vec1 % vec3 );
double skewz = fabs( vec2 % vec3 );
metric_vals->skew = (double)( VERDICT_MAX( skewx, VERDICT_MAX( skewy, skewz ) ) );
}
}
if( metrics_request_flag & V_HEX_TAPER )
{
VerdictVector efg12 = calc_hex_efg( 12, node_pos );
VerdictVector efg13 = calc_hex_efg( 13, node_pos );
VerdictVector efg23 = calc_hex_efg( 23, node_pos );
double taperx = fabs( safe_ratio( efg12.length(), VERDICT_MIN( efg1.length(), efg2.length() ) ) );
double tapery = fabs( safe_ratio( efg13.length(), VERDICT_MIN( efg1.length(), efg3.length() ) ) );
double taperz = fabs( safe_ratio( efg23.length(), VERDICT_MIN( efg2.length(), efg3.length() ) ) );
metric_vals->taper = (double)VERDICT_MAX( taperx, VERDICT_MAX( tapery, taperz ) );
}
}
if( metrics_request_flag & V_HEX_VOLUME )
{
metric_vals->volume = v_hex_volume( 8, coordinates );
}
if( metrics_request_flag &
( V_HEX_JACOBIAN | V_HEX_SCALED_JACOBIAN | V_HEX_CONDITION | V_HEX_SHEAR | V_HEX_SHAPE |
V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE | V_HEX_STRETCH ) )
{
static const double two_thirds = 2.0 / 3.0;
VerdictVector edges[12];
// the length squares
double length_squared[12];<--- Shadow variable
// make vectors from coordinates
make_hex_edges( coordinates, edges );
// calculate the length squares if we need to
if( metrics_request_flag &
( V_HEX_JACOBIAN | V_HEX_SHEAR | V_HEX_SCALED_JACOBIAN | V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE |
V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHEAR_AND_SIZE | V_HEX_STRETCH ) )
make_edge_length_squares( edges, length_squared );
double jacobian = VERDICT_DBL_MAX, scaled_jacobian = VERDICT_DBL_MAX, condition = 0.0, shear = 1.0, shape = 1.0,
oddy = 0.0;
double current_jacobian, current_scaled_jacobian, current_condition, current_shape, detw = 0, det_sum = 0,
current_oddy;
VerdictBoolean rel_size_error = VERDICT_FALSE;
VerdictVector xxi, xet, xze;
// get weights if we need based on average size of a hex
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
v_hex_get_weight( xxi, xet, xze );
detw = xxi % ( xet * xze );
if( detw < VERDICT_DBL_MIN ) rel_size_error = VERDICT_TRUE;
}
xxi = edges[0] - edges[2] + edges[4] - edges[6];
xet = edges[1] - edges[3] + edges[5] - edges[7];
xze = edges[8] + edges[9] + edges[10] + edges[11];
current_jacobian = xxi % ( xet * xze ) / 64.0;
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
current_jacobian *= 64.0;
current_scaled_jacobian =
current_jacobian / sqrt( xxi.length_squared() * xet.length_squared() * xze.length_squared() );
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( xxi, xet, xze );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( xxi, xet, xze );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
// J(0,0,0)
current_jacobian = edges[0] % ( -edges[3] * edges[8] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[0] <= VERDICT_DBL_MIN || length_squared[3] <= VERDICT_DBL_MIN ||
length_squared[8] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[0] * length_squared[3] * length_squared[8] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( edges[0], -edges[3], edges[8] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( edges[0], -edges[3], edges[8] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[0] + length_squared[3] + length_squared[8] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(1,0,0)
current_jacobian = edges[1] % ( -edges[0] * edges[9] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[1] <= VERDICT_DBL_MIN || length_squared[0] <= VERDICT_DBL_MIN ||
length_squared[9] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[1] * length_squared[0] * length_squared[9] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( edges[1], -edges[0], edges[9] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( edges[1], -edges[0], edges[9] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[1] + length_squared[0] + length_squared[9] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(1,1,0)
current_jacobian = edges[2] % ( -edges[1] * edges[10] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[2] <= VERDICT_DBL_MIN || length_squared[1] <= VERDICT_DBL_MIN ||
length_squared[10] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[2] * length_squared[1] * length_squared[10] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( edges[2], -edges[1], edges[10] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( edges[2], -edges[1], edges[10] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[2] + length_squared[1] + length_squared[10] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(0,1,0)
current_jacobian = edges[3] % ( -edges[2] * edges[11] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[3] <= VERDICT_DBL_MIN || length_squared[2] <= VERDICT_DBL_MIN ||
length_squared[11] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[3] * length_squared[2] * length_squared[11] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( edges[3], -edges[2], edges[11] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( edges[3], -edges[2], edges[11] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[3] + length_squared[2] + length_squared[11] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(0,0,1)
current_jacobian = edges[4] % ( -edges[8] * -edges[7] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[4] <= VERDICT_DBL_MIN || length_squared[8] <= VERDICT_DBL_MIN ||
length_squared[7] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[4] * length_squared[8] * length_squared[7] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( edges[4], -edges[8], -edges[7] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( edges[4], -edges[8], -edges[7] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[4] + length_squared[8] + length_squared[7] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(1,0,1)
current_jacobian = -edges[4] % ( edges[5] * -edges[9] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[4] <= VERDICT_DBL_MIN || length_squared[5] <= VERDICT_DBL_MIN ||
length_squared[9] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[4] * length_squared[5] * length_squared[9] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( -edges[4], edges[5], -edges[9] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( -edges[4], edges[5], -edges[9] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[4] + length_squared[5] + length_squared[9] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(1,1,1)
current_jacobian = -edges[5] % ( edges[6] * -edges[10] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[5] <= VERDICT_DBL_MIN || length_squared[6] <= VERDICT_DBL_MIN ||
length_squared[10] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[5] * length_squared[6] * length_squared[10] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( -edges[5], edges[6], -edges[10] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( -edges[5], edges[6], -edges[10] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[5] + length_squared[6] + length_squared[10] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
// J(0,1,1)
current_jacobian = -edges[6] % ( edges[7] * -edges[11] );
if( current_jacobian < jacobian ) jacobian = current_jacobian;
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
det_sum += current_jacobian;
}
if( metrics_request_flag & ( V_HEX_SCALED_JACOBIAN | V_HEX_SHEAR | V_HEX_SHEAR_AND_SIZE ) )
{
if( length_squared[6] <= VERDICT_DBL_MIN || length_squared[7] <= VERDICT_DBL_MIN ||
length_squared[11] <= VERDICT_DBL_MIN )
{
current_scaled_jacobian = VERDICT_DBL_MAX;
}
else
{
current_scaled_jacobian =
current_jacobian / sqrt( length_squared[6] * length_squared[7] * length_squared[11] );
}
if( current_scaled_jacobian < scaled_jacobian ) shear = scaled_jacobian = current_scaled_jacobian;
}
if( metrics_request_flag & V_HEX_CONDITION )
{
current_condition = condition_comp( -edges[6], edges[7], -edges[11] );
if( current_condition > condition )
{
condition = current_condition;
}
}
if( metrics_request_flag & V_HEX_ODDY )
{
current_oddy = oddy_comp( -edges[6], edges[7], -edges[11] );
if( current_oddy > oddy )
{
oddy = current_oddy;
}
}
if( metrics_request_flag & ( V_HEX_SHAPE | V_HEX_SHAPE_AND_SIZE ) )
{
if( current_jacobian > VERDICT_DBL_MIN )
current_shape = 3 * pow( current_jacobian, two_thirds ) /
( length_squared[6] + length_squared[7] + length_squared[11] );
else
current_shape = 0;
if( current_shape < shape )
{
shape = current_shape;
}
}
if( metrics_request_flag & ( V_HEX_RELATIVE_SIZE_SQUARED | V_HEX_SHAPE_AND_SIZE | V_HEX_SHEAR_AND_SIZE ) )
{
if( det_sum > VERDICT_DBL_MIN && rel_size_error != VERDICT_TRUE )
{
double tau = det_sum / ( 8 * detw );
metric_vals->relative_size_squared = (double)VERDICT_MIN( tau * tau, 1.0 / tau / tau );
}
else
metric_vals->relative_size_squared = 0.0;
}
// set values from above calculations
if( metrics_request_flag & V_HEX_JACOBIAN ) metric_vals->jacobian = (double)jacobian;
if( metrics_request_flag & V_HEX_SCALED_JACOBIAN ) metric_vals->scaled_jacobian = (double)scaled_jacobian;
if( metrics_request_flag & V_HEX_CONDITION ) metric_vals->condition = (double)( condition / 3.0 );
if( metrics_request_flag & V_HEX_SHEAR )
{
if( shear < VERDICT_DBL_MIN ) // shear has range 0 to +1
shear = 0;
metric_vals->shear = (double)shear;
}
if( metrics_request_flag & V_HEX_SHAPE ) metric_vals->shape = (double)shape;
if( metrics_request_flag & V_HEX_SHAPE_AND_SIZE )
metric_vals->shape_and_size = (double)( shape * metric_vals->relative_size_squared );
if( metrics_request_flag & V_HEX_SHEAR_AND_SIZE )
metric_vals->shear_and_size = (double)( shear * metric_vals->relative_size_squared );
if( metrics_request_flag & V_HEX_ODDY ) metric_vals->oddy = (double)oddy;
if( metrics_request_flag & V_HEX_STRETCH )
{
static const double HEX_STRETCH_SCALE_FACTOR = sqrt( 3.0 );
double min_edge = length_squared[0];
for( int j = 1; j < 12; j++ )
min_edge = VERDICT_MIN( min_edge, length_squared[j] );
double max_diag = diag_length( 1, coordinates );
metric_vals->stretch = (double)( HEX_STRETCH_SCALE_FACTOR * ( safe_ratio( sqrt( min_edge ), max_diag ) ) );
}
}
if( metrics_request_flag & V_HEX_DIAGONAL ) metric_vals->diagonal = v_hex_diagonal( num_nodes, coordinates );
if( metrics_request_flag & V_HEX_DIMENSION ) metric_vals->dimension = v_hex_dimension( num_nodes, coordinates );
if( metrics_request_flag & V_HEX_DISTORTION ) metric_vals->distortion = v_hex_distortion( num_nodes, coordinates );
// take care of any overflow problems
// max_edge_ratio
if( metric_vals->max_edge_ratio > 0 )
metric_vals->max_edge_ratio = (double)VERDICT_MIN( metric_vals->max_edge_ratio, VERDICT_DBL_MAX );
else
metric_vals->max_edge_ratio = (double)VERDICT_MAX( metric_vals->max_edge_ratio, -VERDICT_DBL_MAX );
// skew
if( metric_vals->skew > 0 )
metric_vals->skew = (double)VERDICT_MIN( metric_vals->skew, VERDICT_DBL_MAX );
else
metric_vals->skew = (double)VERDICT_MAX( metric_vals->skew, -VERDICT_DBL_MAX );
// taper
if( metric_vals->taper > 0 )
metric_vals->taper = (double)VERDICT_MIN( metric_vals->taper, VERDICT_DBL_MAX );
else
metric_vals->taper = (double)VERDICT_MAX( metric_vals->taper, -VERDICT_DBL_MAX );
// volume
if( metric_vals->volume > 0 )
metric_vals->volume = (double)VERDICT_MIN( metric_vals->volume, VERDICT_DBL_MAX );
else
metric_vals->volume = (double)VERDICT_MAX( metric_vals->volume, -VERDICT_DBL_MAX );
// stretch
if( metric_vals->stretch > 0 )
metric_vals->stretch = (double)VERDICT_MIN( metric_vals->stretch, VERDICT_DBL_MAX );
else
metric_vals->stretch = (double)VERDICT_MAX( metric_vals->stretch, -VERDICT_DBL_MAX );
// diagonal
if( metric_vals->diagonal > 0 )
metric_vals->diagonal = (double)VERDICT_MIN( metric_vals->diagonal, VERDICT_DBL_MAX );
else
metric_vals->diagonal = (double)VERDICT_MAX( metric_vals->diagonal, -VERDICT_DBL_MAX );
// dimension
if( metric_vals->dimension > 0 )
metric_vals->dimension = (double)VERDICT_MIN( metric_vals->dimension, VERDICT_DBL_MAX );
else
metric_vals->dimension = (double)VERDICT_MAX( metric_vals->dimension, -VERDICT_DBL_MAX );
// oddy
if( metric_vals->oddy > 0 )
metric_vals->oddy = (double)VERDICT_MIN( metric_vals->oddy, VERDICT_DBL_MAX );
else
metric_vals->oddy = (double)VERDICT_MAX( metric_vals->oddy, -VERDICT_DBL_MAX );
// condition
if( metric_vals->condition > 0 )
metric_vals->condition = (double)VERDICT_MIN( metric_vals->condition, VERDICT_DBL_MAX );
else
metric_vals->condition = (double)VERDICT_MAX( metric_vals->condition, -VERDICT_DBL_MAX );
// jacobian
if( metric_vals->jacobian > 0 )
metric_vals->jacobian = (double)VERDICT_MIN( metric_vals->jacobian, VERDICT_DBL_MAX );
else
metric_vals->jacobian = (double)VERDICT_MAX( metric_vals->jacobian, -VERDICT_DBL_MAX );
// scaled_jacobian
if( metric_vals->scaled_jacobian > 0 )
metric_vals->scaled_jacobian = (double)VERDICT_MIN( metric_vals->scaled_jacobian, VERDICT_DBL_MAX );
else
metric_vals->scaled_jacobian = (double)VERDICT_MAX( metric_vals->scaled_jacobian, -VERDICT_DBL_MAX );
// shear
if( metric_vals->shear > 0 )
metric_vals->shear = (double)VERDICT_MIN( metric_vals->shear, VERDICT_DBL_MAX );
else
metric_vals->shear = (double)VERDICT_MAX( metric_vals->shear, -VERDICT_DBL_MAX );
// shape
if( metric_vals->shape > 0 )
metric_vals->shape = (double)VERDICT_MIN( metric_vals->shape, VERDICT_DBL_MAX );
else
metric_vals->shape = (double)VERDICT_MAX( metric_vals->shape, -VERDICT_DBL_MAX );
// relative_size_squared
if( metric_vals->relative_size_squared > 0 )
metric_vals->relative_size_squared = (double)VERDICT_MIN( metric_vals->relative_size_squared, VERDICT_DBL_MAX );
else
metric_vals->relative_size_squared =
(double)VERDICT_MAX( metric_vals->relative_size_squared, -VERDICT_DBL_MAX );
// shape_and_size
if( metric_vals->shape_and_size > 0 )
metric_vals->shape_and_size = (double)VERDICT_MIN( metric_vals->shape_and_size, VERDICT_DBL_MAX );
else
metric_vals->shape_and_size = (double)VERDICT_MAX( metric_vals->shape_and_size, -VERDICT_DBL_MAX );
// shear_and_size
if( metric_vals->shear_and_size > 0 )
metric_vals->shear_and_size = (double)VERDICT_MIN( metric_vals->shear_and_size, VERDICT_DBL_MAX );
else
metric_vals->shear_and_size = (double)VERDICT_MAX( metric_vals->shear_and_size, -VERDICT_DBL_MAX );
// distortion
if( metric_vals->distortion > 0 )
metric_vals->distortion = (double)VERDICT_MIN( metric_vals->distortion, VERDICT_DBL_MAX );
else
metric_vals->distortion = (double)VERDICT_MAX( metric_vals->distortion, -VERDICT_DBL_MAX );
if( metrics_request_flag & V_HEX_MED_ASPECT_FROBENIUS )
metric_vals->med_aspect_frobenius = v_hex_med_aspect_frobenius( 8, coordinates );
if( metrics_request_flag & V_HEX_EDGE_RATIO ) metric_vals->edge_ratio = v_hex_edge_ratio( 8, coordinates );
}
|