[go: up one dir, main page]

Menu

[r484]: / trunk / src / engine.c  Maximize  Restore  History

Download this file

2466 lines (2000 with data), 81.8 kB

   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
/*******************************************************************************
* This file is part of mdcore.
* Coypright (c) 2010 Pedro Gonnet (pedro.gonnet@durham.ac.uk)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*
******************************************************************************/
/* include some standard header files */
#include <stdlib.h>
#include <stdio.h>
#include <pthread.h>
#include <math.h>
#include <float.h>
#include <string.h>
#ifdef CELL
#include <libspe2.h>
#endif
/* Include conditional headers. */
#include "../config.h"
#ifdef WITH_MPI
#include <mpi.h>
#endif
#ifdef HAVE_OPENMP
#include <omp.h>
#endif
#ifdef WITH_METIS
#include <metis.h>
#endif
/* include local headers */
#include "cycle.h"
#include "errs.h"
#include "fptype.h"
#include "lock.h"
#include "part.h"
#include "cell.h"
#include "task.h"
#include "queue.h"
#include "space.h"
#include "potential.h"
#include "runner.h"
#include "bond.h"
#include "rigid.h"
#include "angle.h"
#include "dihedral.h"
#include "exclusion.h"
#include "reader.h"
#include "engine.h"
/** ID of the last error. */
int engine_err = engine_err_ok;
/* the error macro. */
#define error(id) ( engine_err = errs_register( id , engine_err_msg[-(id)] , __LINE__ , __FUNCTION__ , __FILE__ ) )
/* list of error messages. */
char *engine_err_msg[30] = {
"Nothing bad happened.",
"An unexpected NULL pointer was encountered.",
"A call to malloc failed, probably due to insufficient memory.",
"An error occured when calling a space function.",
"A call to a pthread routine failed.",
"An error occured when calling a runner function.",
"One or more values were outside of the allowed range.",
"An error occured while calling a cell function.",
"The computational domain is too small for the requested operation.",
"mdcore was not compiled with MPI.",
"An error occured while calling an MPI function.",
"An error occured when calling a bond function.",
"An error occured when calling an angle function.",
"An error occured when calling a reader function.",
"An error occured while interpreting the PSF file.",
"An error occured while interpreting the PDB file.",
"An error occured while interpreting the CPF file.",
"An error occured when calling a potential function.",
"An error occured when calling an exclusion function.",
"An error occured while computing the bonded sets.",
"An error occured when calling a dihedral funtion.",
"An error occured when calling a CUDA funtion.",
"mdcore was not compiled with CUDA support.",
"CUDA support is only available in single-precision.",
"Max. number of parts per cell exceeded.",
"An error occured when calling a queue funtion.",
"An error occured when evaluating a rigid constraint.",
"Cell cutoff size doesn't work with METIS",
"METIS library undefined",
"Not yet implemented."
};
/**
* @brief Re-shuffle the particles in the engine.
*
* @param e The #engine on which to run.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_shuffle ( struct engine *e ) {
int cid, k;
struct cell *c;
struct space *s = &e->s;
/* Flush the ghost cells (to avoid overlapping particles) */
#pragma omp parallel for schedule(static), private(cid)
for ( cid = 0 ; cid < s->nr_ghost ; cid++ )
cell_flush( &(s->cells[s->cid_ghost[cid]]) , s->partlist , s->celllist );
/* Shuffle the domain. */
if ( space_shuffle_local( s ) < 0 )
return error(engine_err_space);
#ifdef WITH_MPI
/* Get the incomming particle from other procs if needed. */
if ( e->flags & engine_flag_mpi )
if ( engine_exchange_incomming( e ) < 0 )
return error(engine_err);
#endif
/* Welcome the new particles in each cell, unhook the old ones. */
#pragma omp parallel for schedule(static), private(cid,c,k)
for ( cid = 0 ; cid < s->nr_marked ; cid++ ) {
c = &(s->cells[s->cid_marked[cid]]);
if ( !(c->flags & cell_flag_ghost) )
cell_welcome( c , s->partlist );
else {
for ( k = 0 ; k < c->incomming_count ; k++ )
e->s.partlist[ c->incomming[k].id ] = NULL;
c->incomming_count = 0;
}
}
/* return quietly */
return engine_err_ok;
}
/**
* @brief Set all the engine timers to 0.
*
* @param e The #engine.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_timers_reset ( struct engine *e ) {
int k;
/* Check input nonsense. */
if ( e == NULL )
return error(engine_err_null);
/* Run through the timers and set them to 0. */
for ( k = 0 ; k < engine_timer_last ; k++ )
e->timers[k] = 0;
/* What, that's it? */
return engine_err_ok;
}
/**
* @brief Check if the Verlet-list needs to be updated.
*
* @param e The #engine.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_verlet_update ( struct engine *e ) {
int cid, pid, k;
double dx, w, maxdx = 0.0, skin;
struct cell *c;
struct part *p;
struct space *s = &e->s;
ticks tic;
#ifdef HAVE_OPENMP
int step;
double lmaxdx;
#endif
/* Do we really need to do this? */
if ( !(e->flags & engine_flag_verlet) )
return engine_err_ok;
/* Get the skin width. */
skin = fmin( s->h[0] , fmin( s->h[1] , s->h[2] ) ) - s->cutoff;
/* Get the maximum particle movement. */
if ( !s->verlet_rebuild ) {
#ifdef HAVE_OPENMP
#pragma omp parallel private(c,cid,pid,p,dx,k,w,step,lmaxdx)
{
lmaxdx = 0.0; step = omp_get_num_threads();
for ( cid = omp_get_thread_num() ; cid < s->nr_real ; cid += step ) {
c = &(s->cells[s->cid_real[cid]]);
for ( pid = 0 ; pid < c->count ; pid++ ) {
p = &(c->parts[pid]);
for ( dx = 0.0 , k = 0 ; k < 3 ; k++ ) {
w = p->x[k] - c->oldx[ 4*pid + k ];
dx += w*w;
}
lmaxdx = fmax( dx , lmaxdx );
}
}
#pragma omp critical
maxdx = fmax( lmaxdx , maxdx );
}
#else
for ( cid = 0 ; cid < s->nr_real ; cid++ ) {
c = &(s->cells[s->cid_real[cid]]);
for ( pid = 0 ; pid < c->count ; pid++ ) {
p = &(c->parts[pid]);
for ( dx = 0.0 , k = 0 ; k < 3 ; k++ ) {
w = p->x[k] - c->oldx[ 4*pid + k ];
dx += w*w;
}
maxdx = fmax( dx , maxdx );
}
}
#endif
#ifdef WITH_MPI
/* Collect the maximum displacement from other nodes. */
if ( ( e->flags & engine_flag_mpi ) && ( e->nr_nodes > 1 ) ) {
/* Do not use in-place as it is buggy when async is going on in the background. */
if ( MPI_Allreduce( MPI_IN_PLACE , &maxdx , 1 , MPI_DOUBLE , MPI_MAX , e->comm ) != MPI_SUCCESS )
return error(engine_err_mpi);
}
#endif
/* Are we still in the green? */
maxdx = sqrt(maxdx);
s->verlet_rebuild = ( 2.0*maxdx > skin );
}
/* Do we have to rebuild the Verlet list? */
if ( s->verlet_rebuild ) {
/* printf("engine_verlet_update: re-building verlet lists next step...\n");
printf("engine_verlet_update: maxdx=%e, skin=%e.\n",maxdx,skin); */
/* Wait for any unterminated exchange. */
tic = getticks();
#ifdef WITH_MPI
if ( e->flags & engine_flag_async )
if ( engine_exchange_wait( e ) < 0 )
return error(engine_err);
#endif
tic = getticks() - tic;
e->timers[engine_timer_exchange1] += tic;
e->timers[engine_timer_verlet] -= tic;
/* Move the particles to the respecitve cells. */
if ( engine_shuffle( e ) < 0 )
return error(engine_err);
/* Store the current positions as a reference. */
#pragma omp parallel for schedule(static), private(cid,c,pid,p,k)
for ( cid = 0 ; cid < s->nr_real ; cid++ ) {
c = &(s->cells[s->cid_real[cid]]);
if ( c->oldx == NULL || c->oldx_size < c->count ) {
free(c->oldx);
c->oldx_size = c->size + 20;
c->oldx = (FPTYPE *)malloc( sizeof(FPTYPE) * 4 * c->oldx_size );
}
for ( pid = 0 ; pid < c->count ; pid++ ) {
p = &(c->parts[pid]);
for ( k = 0 ; k < 3 ; k++ )
c->oldx[ 4*pid + k ] = p->x[k];
}
}
/* Set the maximum displacement to zero. */
s->maxdx = 0;
}
/* Otherwise, just store the maximum displacement. */
else
s->maxdx = maxdx;
/* All done! */
return engine_err_ok;
}
/**
* @brief Set-up the engine for distributed-memory parallel operation.
*
* @param e The #engine to set-up.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* This function assumes that #engine_split_bisect or some similar
* function has already been called and that #nodeID, #nr_nodes as
* well as the #cell @c nodeIDs have been set.
*/
int engine_split ( struct engine *e ) {
int i, k, cid, cjd;
struct cell *ci, *cj, *ct;
struct space *s = &(e->s);
/* Check for nonsense inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Start by allocating and initializing the send/recv lists. */
if ( ( e->send = (struct engine_comm *)malloc( sizeof(struct engine_comm) * e->nr_nodes ) ) == NULL ||
( e->recv = (struct engine_comm *)malloc( sizeof(struct engine_comm) * e->nr_nodes ) ) == NULL )
return error(engine_err_malloc);
for ( k = 0 ; k < e->nr_nodes ; k++ ) {
if ( ( e->send[k].cellid = (int *)malloc( sizeof(int) * 100 ) ) == NULL )
return error(engine_err_malloc);
e->send[k].size = 100;
e->send[k].count = 0;
if ( ( e->recv[k].cellid = (int *)malloc( sizeof(int) * 100 ) ) == NULL )
return error(engine_err_malloc);
e->recv[k].size = 100;
e->recv[k].count = 0;
}
/* Un-mark all cells. */
for ( cid = 0 ; cid < s->nr_cells ; cid++ )
s->cells[cid].flags &= ~cell_flag_marked;
/* Loop over each cell pair... */
for ( i = 0 ; i < s->nr_tasks ; i++ ) {
/* Is this task a pair? */
if ( s->tasks[i].type != task_type_pair )
continue;
/* Get the cells in this pair. */
cid = s->tasks[i].i;
cjd = s->tasks[i].j;
ci = &( s->cells[ cid ] );
cj = &( s->cells[ cjd ] );
/* If it is a ghost-ghost pair, skip it. */
if ( (ci->flags & cell_flag_ghost) && (cj->flags & cell_flag_ghost) )
continue;
/* Mark the cells. */
ci->flags |= cell_flag_marked;
cj->flags |= cell_flag_marked;
/* Make cj the ghost cell and bail if both are real. */
if ( ci->flags & cell_flag_ghost ) {
ct = ci; ci = cj; cj = ct;
k = cid; cid = cjd; cjd = k;
}
else if ( !( cj->flags & cell_flag_ghost ) )
continue;
/* Store the communication between cid and cjd. */
/* Store the send, if not already there... */
for ( k = 0 ; k < e->send[cj->nodeID].count && e->send[cj->nodeID].cellid[k] != cid ; k++ );
if ( k == e->send[cj->nodeID].count ) {
if ( e->send[cj->nodeID].count == e->send[cj->nodeID].size ) {
e->send[cj->nodeID].size += 100;
if ( ( e->send[cj->nodeID].cellid = (int *)realloc( e->send[cj->nodeID].cellid , sizeof(int) * e->send[cj->nodeID].size ) ) == NULL )
return error(engine_err_malloc);
}
e->send[cj->nodeID].cellid[ e->send[cj->nodeID].count++ ] = cid;
}
/* Store the recv, if not already there... */
for ( k = 0 ; k < e->recv[cj->nodeID].count && e->recv[cj->nodeID].cellid[k] != cjd ; k++ );
if ( k == e->recv[cj->nodeID].count ) {
if ( e->recv[cj->nodeID].count == e->recv[cj->nodeID].size ) {
e->recv[cj->nodeID].size += 100;
if ( ( e->recv[cj->nodeID].cellid = (int *)realloc( e->recv[cj->nodeID].cellid , sizeof(int) * e->recv[cj->nodeID].size ) ) == NULL )
return error(engine_err_malloc);
}
e->recv[cj->nodeID].cellid[ e->recv[cj->nodeID].count++ ] = cjd;
}
}
/* Nuke all ghost-ghost tasks. */
i = 0;
while ( i < s->nr_tasks ) {
/* Pair? */
if ( s->tasks[i].type == task_type_pair ) {
/* Get the cells in this pair. */
ci = &( s->cells[ s->tasks[i].i ] );
cj = &( s->cells[ s->tasks[i].j ] );
/* If it is a ghost-ghost pair, skip it. */
if ( (ci->flags & cell_flag_ghost) && (cj->flags & cell_flag_ghost) )
s->tasks[i] = s->tasks[ --(s->nr_tasks) ];
else
i += 1;
}
/* Self? */
else if ( s->tasks[i].type == task_type_self ) {
/* Get the cells in this pair. */
ci = &( s->cells[ s->tasks[i].i ] );
/* If it is a ghost-ghost pair, skip it. */
if ( ci->flags & cell_flag_ghost )
s->tasks[i] = s->tasks[ --(s->nr_tasks) ];
else
i += 1;
}
/* Sort? */
else if ( s->tasks[i].type == task_type_sort ) {
/* Get the cells in this pair. */
ci = &( s->cells[ s->tasks[i].i ] );
/* If it is a ghost-ghost pair, skip it. */
if ( !(ci->flags & cell_flag_marked) )
s->tasks[i] = s->tasks[ --(s->nr_tasks) ];
else
i += 1;
}
/* Bonded? */
else if ( s->tasks[i].type == task_type_bonded ) {
/* Clean out the ghost cells. */
for ( k = 0 ; k < e->sets[ s->tasks[i].i ].nr_cells ; k++ )
if ( !( s->cells[ e->sets[ s->tasks[i].i ].cells[k] ].flags & cell_flag_marked ) )
e->sets[ s->tasks[i].i ].cells[k--] = e->sets[ s->tasks[i].i ].cells[ --(e->sets[ s->tasks[i].i ].nr_cells) ];
/* If there are no cells left, skip it. */
if ( e->sets[ s->tasks[i].i ].nr_cells == 0 )
s->tasks[i] = s->tasks[ --(s->nr_tasks) ];
else
i += 1;
}
/* Otherwise? */
else
return error(engine_err_nyi);
}
/* Clear all task dependencies and re-link each sort task with its cell. */
for ( i = 0 ; i < s->nr_tasks ; i++ ) {
s->tasks[i].nr_unlock = 0;
if ( s->tasks[i].type == task_type_sort ) {
s->cells[ s->tasks[i].i ].sort = &s->tasks[i];
s->tasks[i].flags = 0;
}
}
/* Run through the tasks and make each pair depend on the sorts.
Also set the flags for each sort. */
for ( k = 0 ; k < s->nr_tasks ; k++ )
if ( s->tasks[k].type == task_type_pair ) {
if ( task_addunlock( s->cells[ s->tasks[k].i ].sort , &s->tasks[k] ) != 0 ||
task_addunlock( s->cells[ s->tasks[k].j ].sort , &s->tasks[k] ) != 0 )
return error(space_err_task);
s->cells[ s->tasks[k].i ].sort->flags |= 1 << s->tasks[k].flags;
s->cells[ s->tasks[k].j ].sort->flags |= 1 << s->tasks[k].flags;
}
/* Empty unmarked cells. */
for ( k = 0 ; k < s->nr_cells ; k++ )
if ( !( s->cells[k].flags & cell_flag_marked ) )
cell_flush( &s->cells[k] , s->partlist , s->celllist );
/* Set ghost markings on particles. */
for ( cid = 0 ; cid < s->nr_cells ; cid++ )
if ( s->cells[cid].flags & cell_flag_ghost )
for ( k = 0 ; k < s->cells[cid].count ; k++ )
s->cells[cid].parts[k].flags |= part_flag_ghost;
/* Fill the cid lists with marked, local and ghost cells. */
s->nr_real = 0; s->nr_ghost = 0; s->nr_marked = 0;
for ( cid = 0 ; cid < s->nr_cells ; cid++ )
if ( s->cells[cid].flags & cell_flag_marked ) {
s->cid_marked[ s->nr_marked++ ] = cid;
if ( s->cells[cid].flags & cell_flag_ghost ) {
s->cells[cid].id = -s->nr_cells;
s->cid_ghost[ s->nr_ghost++ ] = cid;
}
else {
s->cells[cid].id = s->nr_real;
s->cid_real[ s->nr_real++ ] = cid;
}
}
/* Done deal. */
return engine_err_ok;
}
/**
* @brief Split the computation domain over a number of nodes using METIS graph partitioning.
*
*@param e The #engine to split up.
*@param N The number of computational nodes.
*@param flags Flag telling whether to split the space for MPI or for GPUs.
*
*@return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_split_METIS ( struct engine *e, int N, int flags ){
#ifdef WITH_METIS
//printf("Using METIS algorithm to split the space\n");
int currentIndex, i,j,shiftDim,neighbor;
idx_t vw; //Temporary vertex Weight store
idx_t ew; //Temporary edge Weight store
//Do single GPU version ie. N = 1
if(N==1)
{
if( flags == engine_split_MPI )
{
for(i = 0; i < e->s.nr_cells; i++ )
{
e->s.cells[i].nodeID = 0;
}
}else if ( flags == engine_split_GPU )
{
for( i = 0 ; i < e->s.nr_cells ; i++ )
{
e->s.cells[i].GPUID = 0;
}
}
return engine_err_ok;
}
//Values to adjust weighting of edges dependent on spatial share.
//These values are taken from random simulations of vector distances as described in GONNET 2007.
float FACE = 0.50004f;
float EDGE = 0.16176f;
float CORNER = 0.036213f;
int nr_pairs = 0;
for( i = 0; i < e->s.nr_tasks; i++ )
if(e->s.tasks[i].type == task_type_pair)
nr_pairs++;
nr_pairs *= 2;
/* Check inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Check cell size >= cutoff distance*/
if( e->s.h[0] < e->s.cutoff || e->s.h[1] < e->s.cutoff || e->s.h[2] < e->s.cutoff)
return error(engine_err_cutoff);
//Need to include #include <metis.h>
/*Allocate memory required for METIS*/
/*Number of cells = number of nodes. */
idx_t *xadj = (idx_t*) malloc((e->s.nr_cells+1) * sizeof(idx_t));
if(xadj == NULL) return error(engine_err_malloc);
/*Number of edges = number of cellpairs */
idx_t *adjncy = (idx_t*) malloc (nr_pairs *2 * sizeof(idx_t));
if(adjncy == NULL) return error(engine_err_malloc);
/*Vertex Weights */
idx_t *vwgt = (idx_t*) malloc(e->s.nr_cells * sizeof(idx_t));
if(vwgt == NULL) return error(engine_err_malloc);
/*Edge Weights */
idx_t *adjwgt = (idx_t*) malloc(nr_pairs * 2 * sizeof(idx_t));
if(adjwgt == NULL) return error(engine_err_malloc);
/*Number vertices */
idx_t *nvtxs = (idx_t*) malloc(sizeof(idx_t));
if(nvtxs==NULL)return 1;
*nvtxs = e->s.nr_cells;
//Number constraints
idx_t *ncon = (idx_t*) malloc(sizeof(idx_t));
if(ncon==NULL)return 1;
*ncon = 1;
//Number partitions
idx_t *nparts = (idx_t*) malloc(sizeof(idx_t));
if(nparts==NULL)return 1;
*nparts = N;
/*results*/
idx_t *objval = (idx_t*) malloc(sizeof(idx_t));
if(objval==NULL)return 1;
idx_t *part = (idx_t*) malloc(e->s.nr_cells * (sizeof(idx_t)));
if(part==NULL)return 1;
//Loop over cell pairs and add to array if valid. Needs to be double loop.
currentIndex=0;
for(j=0; j<e->s.nr_cells; j++)
{
vw=0;
ew=0;
xadj[j]=currentIndex;
for(i=0; i<e->s.nr_tasks; i++)
{
if(e->s.tasks[i].type == task_type_sort)
continue;
shiftDim=0;
//If the pair involves cell j.
if(e->s.tasks[i].i==j || e->s.tasks[i].j==j)
{
//If type = self this is a self interaction. All particles interact (n^2-n interactions)
//Increment vertex weight by e->s.cell[j].count ^ 2 - e->s.cell[j].count
if(e->s.tasks[i].type == task_type_self)
{
vw+= e->s.cells[j].count*e->s.cells[j].count - e->s.cells[j].count;
}else{
//Find neighbor cell index.
if(e->s.tasks[i].i==j)
neighbor = e->s.tasks[i].j;
else
neighbor = e->s.tasks[i].i;
//Check how many shift==0.
if(e->s.cells[e->s.tasks[i].i].loc[0] == e->s.cells[e->s.tasks[i].j].loc[0])
shiftDim++;
if(e->s.cells[e->s.tasks[i].i].loc[1] == e->s.cells[e->s.tasks[i].j].loc[1])
shiftDim++;
if(e->s.cells[e->s.tasks[i].i].loc[2] == e->s.cells[e->s.tasks[i].j].loc[2])
shiftDim++;
//If shiftDim = 2 this is a face interaction. Add edge from j to i
//Edge has weight e->s.cell[j].count * e->s.cell[neighbor].count * FACE (Estimated number of distance calculations)
if(shiftDim==2)
{
ew = e->s.cells[j].count * e->s.cells[neighbor].count * FACE;
//Vertex Weights are the sum of edge weights + self interaction
vw+= ew;
//Add an edge and the edge weight to the graph
adjncy[currentIndex] = neighbor;
adjwgt[currentIndex] = ew;
currentIndex++;
}
//If shiftDim = 1 this is an edge interaction. Add edge from j to i
//Edge has weight e->s.cell[j].count * e->s.cell[neightbor.count * EDGE (Estimated number of distance calculations)
if(shiftDim==1)
{
ew = e->s.cells[j].count * e->s.cells[neighbor].count * EDGE;
//Vertex Weights are the sum of edge weights + self interaction
vw+=ew;
//Add an edge and the edge weight to the graph
adjncy[currentIndex] = neighbor;
adjwgt[currentIndex] = ew;
currentIndex++;
}
//If shiftDim = 0 this is an corner interaction. Add edge from j to i
//Edge has weight e->s.cell[j].count * e->s.cell[neightbor.count * CORNER (Estimated number of distance calculations)
if(shiftDim==0)
{
ew = e->s.cells[j].count * e->s.cells[neighbor].count * CORNER;
//Vertex Weights are the sum of edge weights + self interaction
vw+=ew;
//Add an edge and the edge weight to the graph
adjncy[currentIndex] = neighbor;
adjwgt[currentIndex] = ew;
currentIndex++;
}
}
}
}
//We now know vertex weight.
vwgt[j] = vw;
}
//Need to add the final thing METIS needs to xadj
xadj[e->s.nr_cells]=currentIndex;
//Setup METIS options
idx_t *options = (idx_t*) malloc(METIS_NOPTIONS * sizeof(idx_t));
if(options==NULL)return 1;
METIS_SetDefaultOptions(options);
options[METIS_OPTION_PTYPE]=METIS_PTYPE_KWAY;
options[METIS_OPTION_OBJTYPE]=METIS_OBJTYPE_CUT;
options[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM;
options[METIS_OPTION_IPTYPE]=METIS_IPTYPE_GROW;
options[METIS_OPTION_NCUTS] = 1;
options[METIS_OPTION_NSEPS] = 1;
options[METIS_OPTION_NUMBERING] = 0; //C-style numbering :)
options[METIS_OPTION_NITER] = 10;
options[METIS_OPTION_SEED] = 1;
options[METIS_OPTION_CONTIG] = 1;
options[METIS_OPTION_UFACTOR] = 30;
// options[METIS_OPTION_DBGLVL] = METIS_DBG_INFO;
//Run METIS to partition the graph
METIS_PartGraphKway( nvtxs , ncon , xadj , adjncy , vwgt , NULL , adjwgt , nparts , NULL , NULL , options , objval , part );
if( flags == engine_split_MPI )
{
for(i = 0; i < e->s.nr_cells; i++ )
{
e->s.cells[i].nodeID = part[i];
//Not my cell? Mark as ghost.
if(part[i] != e->nodeID)
e->s.cells[i].flags |= cell_flag_ghost;
e->nr_nodes = N;
}
}else if ( flags == engine_split_GPU )
{
for( i = 0 ; i < e->s.nr_cells ; i++ )
{
e->s.cells[i].GPUID = part[i];
}
}
/* Store the number of nodes. */
/*Free memory used by METIS*/
free(xadj);
free(adjncy);
free(vwgt);
free(adjwgt);
free(nvtxs);
free(ncon);
free(nparts);
free(objval);
free(part);
free(options);
printf("Successfully split the space\n");
/* Call it a day. */
return engine_err_ok;
#else
/* Bisect recursively */
/* Interior, recursive function that actually does the split. */
int engine_split_bisect_rec( int N_min , int N_max , int x_min , int x_max , int y_min , int y_max , int z_min , int z_max , int flags) {
int i, j, k, m, Nm;
int hx, hy, hz;
unsigned int flag = 0;
struct cell *c;
/* Check inputs. */
if ( x_max < x_min || y_max < y_min || z_max < z_min )
return error(engine_err_domain);
/* Is there nothing left to split? */
if ( N_min == N_max ) {
/* Flag as ghost or not? */
if( flags == engine_split_MPI )
{
if ( N_min != e->nodeID )
flag = cell_flag_ghost;
/* printf("engine_split_bisect: marking range [ %i..%i , %i..%i , %i..%i ] with flag %i.\n",
x_min, x_max, y_min, y_max, z_min, z_max, flag ); */
/* Run through the cells. */
for ( i = x_min ; i < x_max ; i++ )
for ( j = y_min ; j < y_max ; j++ )
for ( k = z_min ; k < z_max ; k++ ) {
c = &( e->s.cells[ space_cellid(&(e->s),i,j,k) ] );
c->flags |= flag;
c->nodeID = N_min;
}
}else{
for ( i = x_min ; i < x_max ; i++ )
for ( j = y_min ; j < y_max ; j++ )
for ( k = z_min ; k < z_max ; k++ ) {
c = &( e->s.cells[ space_cellid(&(e->s),i,j,k) ] );
c->GPUID = N_min;
}
}
}
/* Otherwise, bisect. */
else {
hx = x_max - x_min;
hy = y_max - y_min;
hz = z_max - z_min;
Nm = (N_min + N_max) / 2;
/* Is the x-axis the largest? */
if ( hx > hy && hx > hz ) {
m = (x_min + x_max) / 2;
if ( engine_split_bisect_rec( N_min , Nm , x_min , m , y_min , y_max , z_min , z_max , flags) < 0 ||
engine_split_bisect_rec( Nm+1 , N_max , m , x_max , y_min , y_max , z_min , z_max , flags) < 0 )
return error(engine_err);
}
/* Nope, maybe the y-axis? */
else if ( hy > hz ) {
m = (y_min + y_max) / 2;
if ( engine_split_bisect_rec( N_min , Nm , x_min , x_max , y_min , m , z_min , z_max , flags ) < 0 ||
engine_split_bisect_rec( Nm+1 , N_max , x_min , x_max , m , y_max , z_min , z_max , flags) < 0 )
return error(engine_err);
}
/* Then it has to be the z-axis. */
else {
m = (z_min + z_max) / 2;
if ( engine_split_bisect_rec( N_min , Nm , x_min , x_max , y_min , y_max , z_min , m , flags) < 0 ||
engine_split_bisect_rec( Nm+1 , N_max , x_min , x_max , y_min , y_max , m , z_max , flags) < 0 )
return error(engine_err);
}
}
/* So far, so good! */
return engine_err_ok;
}
/* Check inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Call the recursive bisection. */
if ( engine_split_bisect_rec( 0 , N-1 , 0 , e->s.cdim[0] , 0 , e->s.cdim[1] , 0 , e->s.cdim[2] , flags) < 0 )
return error(engine_err);
/* Store the number of nodes. */
e->nr_nodes = N;
/* Call it a day. */
return engine_err_ok;
#endif
}
/**
* @brief Split the computational domain over a number of nodes using
* bisection.
*
* @param e The #engine to split up.
* @param N The number of computational nodes.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_split_bisect ( struct engine *e , int N ) {
/* Interior, recursive function that actually does the split. */
int engine_split_bisect_rec( int N_min , int N_max , int x_min , int x_max , int y_min , int y_max , int z_min , int z_max ) {
int i, j, k, m, Nm;
int hx, hy, hz;
unsigned int flag = 0;
struct cell *c;
/* Check inputs. */
if ( x_max < x_min || y_max < y_min || z_max < z_min )
return error(engine_err_domain);
/* Is there nothing left to split? */
if ( N_min == N_max ) {
/* Flag as ghost or not? */
if ( N_min != e->nodeID )
flag = cell_flag_ghost;
/* printf("engine_split_bisect: marking range [ %i..%i , %i..%i , %i..%i ] with flag %i.\n",
x_min, x_max, y_min, y_max, z_min, z_max, flag ); */
/* Run through the cells. */
for ( i = x_min ; i < x_max ; i++ )
for ( j = y_min ; j < y_max ; j++ )
for ( k = z_min ; k < z_max ; k++ ) {
c = &( e->s.cells[ space_cellid(&(e->s),i,j,k) ] );
c->flags |= flag;
c->nodeID = N_min;
}
}
/* Otherwise, bisect. */
else {
hx = x_max - x_min;
hy = y_max - y_min;
hz = z_max - z_min;
Nm = (N_min + N_max) / 2;
/* Is the x-axis the largest? */
if ( hx > hy && hx > hz ) {
m = (x_min + x_max) / 2;
if ( engine_split_bisect_rec( N_min , Nm , x_min , m , y_min , y_max , z_min , z_max ) < 0 ||
engine_split_bisect_rec( Nm+1 , N_max , m , x_max , y_min , y_max , z_min , z_max ) < 0 )
return error(engine_err);
}
/* Nope, maybe the y-axis? */
else if ( hy > hz ) {
m = (y_min + y_max) / 2;
if ( engine_split_bisect_rec( N_min , Nm , x_min , x_max , y_min , m , z_min , z_max ) < 0 ||
engine_split_bisect_rec( Nm+1 , N_max , x_min , x_max , m , y_max , z_min , z_max ) < 0 )
return error(engine_err);
}
/* Then it has to be the z-axis. */
else {
m = (z_min + z_max) / 2;
if ( engine_split_bisect_rec( N_min , Nm , x_min , x_max , y_min , y_max , z_min , m ) < 0 ||
engine_split_bisect_rec( Nm+1 , N_max , x_min , x_max , y_min , y_max , m , z_max ) < 0 )
return error(engine_err);
}
}
/* So far, so good! */
return engine_err_ok;
}
/* Check inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Call the recursive bisection. */
if ( engine_split_bisect_rec( 0 , N-1 , 0 , e->s.cdim[0] , 0 , e->s.cdim[1] , 0 , e->s.cdim[2] ) < 0 )
return error(engine_err);
/* Store the number of nodes. */
e->nr_nodes = N;
/* Call it a day. */
return engine_err_ok;
}
/**
* @brief Clear all particles from this #engine.
*
* @param e The #engine to flush.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_flush ( struct engine *e ) {
/* check input. */
if ( e == NULL )
return error(engine_err_null);
/* Clear the space. */
if ( space_flush( &e->s ) < 0 )
return error(engine_err_space);
/* done for now. */
return engine_err_ok;
}
/**
* @brief Clear all particles from this #engine's ghost cells.
*
* @param e The #engine to flush.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_flush_ghosts ( struct engine *e ) {
/* check input. */
if ( e == NULL )
return error(engine_err_null);
/* Clear the space. */
if ( space_flush_ghosts( &e->s ) < 0 )
return error(engine_err_space);
/* done for now. */
return engine_err_ok;
}
/**
* @brief Set the explicit electrostatic potential.
*
* @param e The #engine.
* @param ep The electrostatic #potential.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* If @c ep is not @c NULL, the flag #engine_flag_explepot is set,
* otherwise it is cleared.
*/
int engine_setexplepot ( struct engine *e , struct potential *ep ) {
/* check inputs. */
if ( e == NULL )
return error(engine_err_null);
/* was a potential supplied? */
if ( ep != NULL ) {
/* set the flag. */
e->flags |= engine_flag_explepot;
/* set the potential. */
e->ep = ep;
}
/* otherwise, just clear the flag. */
else
e->flags &= ~engine_flag_explepot;
/* done for now. */
return engine_err_ok;
}
/**
* @brief Unload a set of particle data from the #engine.
*
* @param e The #engine.
* @param x An @c N times 3 array of the particle positions.
* @param v An @c N times 3 array of the particle velocities.
* @param type A vector of length @c N of the particle type IDs.
* @param pid A vector of length @c N of the particle IDs.
* @param vid A vector of length @c N of the particle virtual IDs.
* @param q A vector of length @c N of the individual particle charges.
* @param flags A vector of length @c N of the particle flags.
* @param epot A pointer to a #double in which to store the total potential energy.
* @param N the maximum number of particles.
*
* @return The number of particles unloaded or < 0 on
* error (see #engine_err).
*
* The fields @c x, @c v, @c type, @c pid, @c vid, @c q, @c epot and/or @c flags may be NULL.
*/
int engine_unload ( struct engine *e , double *x , double *v , int *type , int *pid , int *vid , double *q , unsigned int *flags , double *epot , int N ) {
struct part *p;
struct cell *c;
int j, k, cid, count = 0, *ind;
double epot_acc = 0.0;
/* check the inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Allocate and fill the indices. */
if ( ( ind = (int *)alloca( sizeof(int) * (e->s.nr_cells + 1) ) ) == NULL )
return error(engine_err_malloc);
ind[0] = 0;
for ( k = 0 ; k < e->s.nr_cells ; k++ )
ind[k+1] = ind[k] + e->s.cells[k].count;
if ( ind[e->s.nr_cells] > N )
return error(engine_err_range);
/* Loop over each cell. */
#pragma omp parallel for schedule(static), private(cid,count,c,k,p,j), reduction(+:epot_acc)
for ( cid = 0 ; cid < e->s.nr_cells ; cid++ ) {
/* Get a hold of the cell. */
c = &( e->s.cells[cid] );
count = ind[cid];
/* Collect the potential energy if requested. */
epot_acc += c->epot;
/* Loop over the parts in this cell. */
for ( k = 0 ; k < c->count ; k++ ) {
/* Get a hold of the particle. */
p = &( c->parts[k] );
/* get this particle's data, where requested. */
if ( x != NULL )
for ( j = 0 ; j < 3 ; j++ )
x[count*3+j] = c->origin[j] + p->x[j];
if ( v != NULL)
for ( j = 0 ; j < 3 ; j++ )
v[count*3+j] = p->v[j];
if ( type != NULL )
type[count] = p->type;
if ( pid != NULL )
pid[count] = p->id;
if ( vid != NULL )
vid[count] = p->vid;
if ( q != NULL )
q[count] = p->q;
if ( flags != NULL )
flags[count] = p->flags;
/* Step-up the counter. */
count += 1;
}
}
/* Write back the potential energy, if requested. */
if ( epot != NULL )
*epot += epot_acc;
/* to the pub! */
return ind[e->s.nr_cells];
}
/**
* @brief Unload a set of particle data from the marked cells of an #engine
*
* @param e The #engine.
* @param x An @c N times 3 array of the particle positions.
* @param v An @c N times 3 array of the particle velocities.
* @param type A vector of length @c N of the particle type IDs.
* @param pid A vector of length @c N of the particle IDs.
* @param vid A vector of length @c N of the particle virtual IDs.
* @param q A vector of length @c N of the individual particle charges.
* @param flags A vector of length @c N of the particle flags.
* @param epot A pointer to a #double in which to store the total potential energy.
* @param N the maximum number of particles.
*
* @return The number of particles unloaded or < 0 on
* error (see #engine_err).
*
* The fields @c x, @c v, @c type, @c pid, @c vid, @c q, @c epot and/or @c flags may be NULL.
*/
int engine_unload_marked ( struct engine *e , double *x , double *v , int *type , int *pid , int *vid , double *q , unsigned int *flags , double *epot , int N ) {
struct part *p;
struct cell *c;
int j, k, cid, count = 0, *ind;
double epot_acc = 0.0;
/* check the inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Allocate and fill the indices. */
if ( ( ind = (int *)alloca( sizeof(int) * (e->s.nr_cells + 1) ) ) == NULL )
return error(engine_err_malloc);
ind[0] = 0;
for ( k = 0 ; k < e->s.nr_cells ; k++ )
if ( e->s.cells[k].flags & cell_flag_marked )
ind[k+1] = ind[k] + e->s.cells[k].count;
else
ind[k+1] = ind[k];
if ( ind[e->s.nr_cells] > N )
return error(engine_err_range);
/* Loop over each cell. */
#pragma omp parallel for schedule(static), private(cid,count,c,k,p,j), reduction(+:epot_acc)
for ( cid = 0 ; cid < e->s.nr_marked ; cid++ ) {
/* Get a hold of the cell. */
c = &( e->s.cells[e->s.cid_marked[cid]] );
count = ind[e->s.cid_marked[cid]];
/* Collect the potential energy if requested. */
epot_acc += c->epot;
/* Loop over the parts in this cell. */
for ( k = 0 ; k < c->count ; k++ ) {
/* Get a hold of the particle. */
p = &( c->parts[k] );
/* get this particle's data, where requested. */
if ( x != NULL )
for ( j = 0 ; j < 3 ; j++ )
x[count*3+j] = c->origin[j] + p->x[j];
if ( v != NULL)
for ( j = 0 ; j < 3 ; j++ )
v[count*3+j] = p->v[j];
if ( type != NULL )
type[count] = p->type;
if ( pid != NULL )
pid[count] = p->id;
if ( vid != NULL )
vid[count] = p->vid;
if ( q != NULL )
q[count] = p->q;
if ( flags != NULL )
flags[count] = p->flags;
/* Step-up the counter. */
count += 1;
}
}
/* Write back the potential energy, if requested. */
if ( epot != NULL )
*epot += epot_acc;
/* to the pub! */
return ind[e->s.nr_cells];
}
/**
* @brief Unload real particles that may have wandered into a ghost cell.
*
* @param e The #engine.
* @param x An @c N times 3 array of the particle positions.
* @param v An @c N times 3 array of the particle velocities.
* @param type A vector of length @c N of the particle type IDs.
* @param pid A vector of length @c N of the particle IDs.
* @param vid A vector of length @c N of the particle virtual IDs.
* @param q A vector of length @c N of the individual particle charges.
* @param flags A vector of length @c N of the particle flags.
* @param epot A pointer to a #double in which to store the total potential energy.
* @param N the maximum number of particles.
*
* @return The number of particles unloaded or < 0 on
* error (see #engine_err).
*
* The fields @c x, @c v, @c type, @c vid, @c pid, @c q, @c epot and/or @c flags may be NULL.
*/
int engine_unload_strays ( struct engine *e , double *x , double *v , int *type , int *pid , int *vid , double *q , unsigned int *flags , double *epot , int N ) {
struct part *p;
struct cell *c;
int j, k, cid, count = 0;
double epot_acc = 0.0;
/* check the inputs. */
if ( e == NULL )
return error(engine_err_null);
/* Loop over each cell. */
for ( cid = 0 ; cid < e->s.nr_real ; cid++ ) {
/* Get a hold of the cell. */
c = &( e->s.cells[e->s.cid_real[cid]] );
/* Collect the potential energy if requested. */
epot_acc += c->epot;
/* Loop over the parts in this cell. */
for ( k = c->count-1 ; k >= 0 && !(c->parts[k].flags & part_flag_ghost) ; k-- ) {
/* Get a hold of the particle. */
p = &( c->parts[k] );
if ( p->flags & part_flag_ghost )
continue;
/* get this particle's data, where requested. */
if ( x != NULL )
for ( j = 0 ; j < 3 ; j++ )
x[count*3+j] = c->origin[j] + p->x[j];
if ( v != NULL)
for ( j = 0 ; j < 3 ; j++ )
v[count*3+j] = p->v[j];
if ( type != NULL )
type[count] = p->type;
if ( pid != NULL )
pid[count] = p->id;
if ( vid != NULL )
vid[count] = p->vid;
if ( q != NULL )
q[count] = p->q;
if ( flags != NULL )
flags[count] = p->flags;
/* increase the counter. */
count += 1;
}
}
/* Write back the potential energy, if requested. */
if ( epot != NULL )
*epot += epot_acc;
/* to the pub! */
return count;
}
/**
* @brief Load a set of particle data.
*
* @param e The #engine.
* @param x An @c N times 3 array of the particle positions.
* @param v An @c N times 3 array of the particle velocities.
* @param type A vector of length @c N of the particle type IDs.
* @param pid A vector of length @c N of the particle IDs.
* @param vid A vector of length @c N of the particle virtual IDs.
* @param q A vector of length @c N of the individual particle charges.
* @param flags A vector of length @c N of the particle flags.
* @param N the number of particles to load.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* If the parameters @c v, @c flags, @c vid or @c q are @c NULL, then
* these values are set to zero.
*/
int engine_load ( struct engine *e , double *x , double *v , int *type , int *pid , int *vid , double *q , unsigned int *flags , int N ) {
struct part p;
struct space *s;
int j, k;
/* check the inputs. */
if ( e == NULL || x == NULL || type == NULL )
return error(engine_err_null);
/* Get a handle on the space. */
s = &(e->s);
/* init the velocity and charge in case not specified. */
p.v[0] = 0.0; p.v[1] = 0.0; p.v[2] = 0.0;
p.f[0] = 0.0; p.f[1] = 0.0; p.f[2] = 0.0;
p.q = 0.0;
p.flags = part_flag_none;
/* loop over the entries. */
for ( j = 0 ; j < N ; j++ ) {
/* set the particle data. */
p.type = type[j];
if ( pid != NULL )
p.id = pid[j];
else
p.id = j;
if ( vid != NULL )
p.vid = vid[j];
if ( flags != NULL )
p.flags = flags[j];
if ( v != NULL )
for ( k = 0 ; k < 3 ; k++ )
p.v[k] = v[j*3+k];
if ( q != 0 )
p.q = q[j];
/* add the part to the space. */
if ( space_addpart( s , &p , &x[3*j] ) < 0 )
return error(engine_err_space);
}
/* to the pub! */
return engine_err_ok;
}
/**
* @brief Load a set of particle data as ghosts
*
* @param e The #engine.
* @param x An @c N times 3 array of the particle positions.
* @param v An @c N times 3 array of the particle velocities.
* @param type A vector of length @c N of the particle type IDs.
* @param pid A vector of length @c N of the particle IDs.
* @param vid A vector of length @c N of the particle virtual IDs.
* @param q A vector of length @c N of the individual particle charges.
* @param flags A vector of length @c N of the particle flags.
* @param N the number of particles to load.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* If the parameters @c v, @c flags, @c vid or @c q are @c NULL, then
* these values are set to zero.
*/
int engine_load_ghosts ( struct engine *e , double *x , double *v , int *type , int *pid , int *vid , double *q , unsigned int *flags , int N ) {
struct part p;
struct space *s;
int j, k;
/* check the inputs. */
if ( e == NULL || x == NULL || type == NULL )
return error(engine_err_null);
/* Get a handle on the space. */
s = &(e->s);
/* init the velocity and charge in case not specified. */
p.v[0] = 0.0; p.v[1] = 0.0; p.v[2] = 0.0;
p.f[0] = 0.0; p.f[1] = 0.0; p.f[2] = 0.0;
p.q = 0.0;
p.flags = part_flag_ghost;
/* loop over the entries. */
for ( j = 0 ; j < N ; j++ ) {
/* set the particle data. */
p.type = type[j];
if ( pid != NULL )
p.id = pid[j];
else
p.id = j;
if ( vid != NULL )
p.vid = vid[j];
if ( flags != NULL )
p.flags = flags[j] | part_flag_ghost;
if ( v != NULL )
for ( k = 0 ; k < 3 ; k++ )
p.v[k] = v[j*3+k];
if ( q != 0 )
p.q = q[j];
/* add the part to the space. */
if ( space_addpart( s , &p , &x[3*j] ) < 0 )
return error(engine_err_space);
}
/* to the pub! */
return engine_err_ok;
}
/**
* @brief Look for a given type by name.
*
* @param e The #engine.
* @param name The type name.
*
* @return The type ID or < 0 on error (see #engine_err).
*/
int engine_gettype ( struct engine *e , char *name ) {
int k;
/* check for nonsense. */
if ( e == NULL || name == NULL )
return error(engine_err_null);
/* Loop over the types... */
for ( k = 0 ; k < e->nr_types ; k++ ) {
/* Compare the name. */
if ( strcmp( e->types[k].name , name ) == 0 )
return k;
}
/* Otherwise, nothing found... */
return engine_err_range;
}
/**
* @brief Look for a given type by its second name.
*
* @param e The #engine.
* @param name2 The type name2.
*
* @return The type ID or < 0 on error (see #engine_err).
*/
int engine_gettype2 ( struct engine *e , char *name2 ) {
int k;
/* check for nonsense. */
if ( e == NULL || name2 == NULL )
return error(engine_err_null);
/* Loop over the types... */
for ( k = 0 ; k < e->nr_types ; k++ ) {
/* Compare the name. */
if ( strcmp( e->types[k].name2 , name2 ) == 0 )
return k;
}
/* Otherwise, nothing found... */
return engine_err_range;
}
/**
* @brief Add a type definition.
*
* @param e The #engine.
* @param mass The particle type mass.
* @param charge The particle type charge.
* @param name Particle name, can be @c NULL.
* @param name2 Particle second name, can be @c NULL.
*
* @return The type ID or < 0 on error (see #engine_err).
*
* The particle type ID must be an integer greater or equal to 0
* and less than the value @c max_type specified in #engine_init.
*/
int engine_addtype ( struct engine *e , double mass , double charge , char *name , char *name2 ) {
/* check for nonsense. */
if ( e == NULL )
return error(engine_err_null);
if ( e->nr_types >= e->max_type )
return error(engine_err_range);
/* set the type. */
e->types[e->nr_types].mass = mass;
e->types[e->nr_types].imass = 1.0 / mass;
e->types[e->nr_types].charge = charge;
if ( name != NULL )
strcpy( e->types[e->nr_types].name , name );
else
strcpy( e->types[e->nr_types].name , "X" );
if ( name2 != NULL )
strcpy( e->types[e->nr_types].name2 , name2 );
else
strcpy( e->types[e->nr_types].name2 , "X" );
/* bring good tidings. */
return e->nr_types++;
}
/**
* @brief Add an interaction potential.
*
* @param e The #engine.
* @param p The #potential to add to the #engine.
* @param i ID of particle type for this interaction.
* @param j ID of second particle type for this interaction.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* Adds the given potential for pairs of particles of type @c i and @c j,
* where @c i and @c j may be the same type ID.
*/
int engine_addpot ( struct engine *e , struct potential *p , int i , int j ) {
/* check for nonsense. */
if ( e == NULL )
return error(engine_err_null);
if ( i < 0 || i >= e->max_type || j < 0 || j >= e->max_type )
return error(engine_err_range);
/* store the potential. */
e->p[ i * e->max_type + j ] = p;
if ( i != j )
e->p[ j * e->max_type + i ] = p;
/* end on a good note. */
return engine_err_ok;
}
/**
* @brief Start the runners in the given #engine.
*
* @param e The #engine to start.
* @param nr_runners The number of runners start.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* Allocates and starts the specified number of #runner. Also initializes
* the Verlet lists.
*/
int engine_start ( struct engine *e , int nr_runners , int nr_queues ) {
int cid, pid, k, i;
struct cell *c;
struct part *p;
struct space *s = &e->s;
/* Is MPI really needed? */
if ( e->flags & engine_flag_mpi && e->nr_nodes == 1 )
e->flags &= ~( engine_flag_mpi | engine_flag_async );
#ifdef WITH_MPI
/* Set up async communication? */
if ( e->flags & engine_flag_async ) {
/* Init the mutex and condition variable for the asynchronous communication. */
if ( pthread_mutex_init( &e->xchg_mutex , NULL ) != 0 ||
pthread_cond_init( &e->xchg_cond , NULL ) != 0 ||
pthread_mutex_init( &e->xchg2_mutex , NULL ) != 0 ||
pthread_cond_init( &e->xchg2_cond , NULL ) != 0 )
return error(engine_err_pthread);
/* Set the exchange flags. */
e->xchg_started = 0;
e->xchg_running = 0;
e->xchg2_started = 0;
e->xchg2_running = 0;
/* Start a thread with the async exchange. */
if ( pthread_create( &e->thread_exchg , NULL , (void *(*)(void *))engine_exchange_async_run , e ) != 0 )
return error(engine_err_pthread);
if ( pthread_create( &e->thread_exchg2 , NULL , (void *(*)(void *))engine_exchange_rigid_async_run , e ) != 0 )
return error(engine_err_pthread);
}
#endif
/* Fill-in the Verlet lists if needed. */
if ( e->flags & engine_flag_verlet ) {
/* Shuffle the domain. */
if ( engine_shuffle( e ) < 0 )
return error(engine_err);
/* Store the current positions as a reference. */
#pragma omp parallel for schedule(static), private(cid,c,pid,p,k)
for ( cid = 0 ; cid < s->nr_real ; cid++ ) {
c = &(s->cells[s->cid_real[cid]]);
if ( c->oldx == NULL || c->oldx_size < c->count ) {
free(c->oldx);
c->oldx_size = c->size + 20;
c->oldx = (FPTYPE *)malloc( sizeof(FPTYPE) * 4 * c->oldx_size );
}
for ( pid = 0 ; pid < c->count ; pid++ ) {
p = &(c->parts[pid]);
for ( k = 0 ; k < 3 ; k++ )
c->oldx[ 4*pid + k ] = p->x[k];
}
}
/* Re-set the Verlet rebuild flag. */
s->verlet_rebuild = 1;
}
/* Set-up the bonded tasks if needed. */
if ( e->flags & engine_flag_parbonded ) {
int grid[3] = { e->s.cdim[0] / 2 , e->s.cdim[1] / 2 , e->s.cdim[2] / 2 };
if ( engine_bonded_makesets( e , grid ) < 0 )
return error(engine_err);
}
/* Do we even need runners? */
if ( e->flags & engine_flag_cuda ) {
/* Set the number of runners. */
e->nr_runners = nr_runners;
#if defined(HAVE_CUDA) && defined(WITH_CUDA)
/* Load the potentials and pairs to the CUDA device. */
if ( engine_cuda_load( e ) < 0 )
return error(engine_err);
#else
/* Was not compiled with CUDA support. */
return error(engine_err_nocuda);
#endif
}
else {
/* Allocate the queues */
if ( ( e->queues = (struct queue *)malloc( sizeof(struct queue) * nr_queues )) == NULL )
return error(engine_err_malloc);
e->nr_queues = nr_queues;
/* Initialize and fill the queues. */
for ( i = 0 ; i < e->nr_queues ; i++ )
if ( queue_init( &e->queues[i] , 2*s->nr_tasks/e->nr_queues , e , s->tasks ) != queue_err_ok )
return error(engine_err_queue);
for ( i = 0 ; i < s->nr_tasks ; i++ )
if ( queue_insert( &e->queues[ i % e->nr_queues ] , &s->tasks[i] ) < 0 )
return error(engine_err_queue);
/* (Allocate the runners */
if ( ( e->runners = (struct runner *)malloc( sizeof(struct runner) * nr_runners )) == NULL )
return error(engine_err_malloc);
e->nr_runners = nr_runners;
/* initialize the runners. */
for ( i = 0 ; i < nr_runners ; i++ )
if ( runner_init( &e->runners[ i ] , e , i ) < 0 )
return error(engine_err_runner);
/* wait for the runners to be in place */
while (e->barrier_count != e->nr_runners)
if (pthread_cond_wait(&e->done_cond,&e->barrier_mutex) != 0)
return error(engine_err_pthread);
}
/* Set the number of runners. */
e->nr_runners = nr_runners;
/* all is well... */
return engine_err_ok;
}
/**
* @brief Start the SPU-associated runners in the given #engine.
*
* @param e The #engine to start.
* @param nr_runners The number of runners start.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* Allocates and starts the specified number of #runner.
*/
int engine_start_SPU ( struct engine *e , int nr_runners ) {
int i;
struct runner *temp;
#ifdef WITH_MPI
/* Set up async communication? */
if ( e->flags & engine_flag_async ) {
/* Init the mutex and condition variable for the asynchronous communication. */
if ( pthread_mutex_init( &e->xchg_mutex , NULL ) != 0 ||
pthread_cond_init( &e->xchg_cond , NULL ) != 0 )
return error(engine_err_pthread);
/* Set the exchange flags. */
e->xchg_started = 0;
e->xchg_running = 0;
/* Start a thread with the async exchange. */
if ( pthread_create( &e->thread_exchg , NULL , (void *(*)(void *))engine_exchange_async_run , e ) != 0 )
return error(engine_err_pthread);
}
#endif
/* (re)allocate the runners */
if ( e->nr_runners == 0 ) {
if ( ( e->runners = (struct runner *)malloc( sizeof(struct runner) * nr_runners )) == NULL )
return error(engine_err_malloc);
}
else {
if ( ( temp = (struct runner *)malloc( sizeof(struct runner) * (e->nr_runners + nr_runners) )) == NULL )
return error(engine_err_malloc);
memcpy( temp , e->runners , sizeof(struct runner) * e->nr_runners );
free( e->runners );
e->runners = temp;
}
/* initialize the runners. */
for ( i = 0 ; i < nr_runners ; i++ )
if ( runner_init_SPU(&e->runners[e->nr_runners + i],e,e->nr_runners + i) < 0 )
return error(engine_err_runner);
e->nr_runners += nr_runners;
/* wait for the runners to be in place */
while (e->barrier_count != e->nr_runners)
if (pthread_cond_wait(&e->done_cond,&e->barrier_mutex) != 0)
return error(engine_err_pthread);
/* all is well... */
return engine_err_ok;
}
/**
* @brief Compute the nonbonded interactions in the current step.
*
* @param e The #engine on which to run.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* This routine advances the timestep counter by one, prepares the #space
* for a timestep, releases the #runner's associated with the #engine
* and waits for them to finnish.
*/
int engine_nonbond_eval ( struct engine *e ) {
int k;
/* Re-set the queues. */
for ( k = 0 ; k < e->nr_queues ; k++ )
e->queues[k].next = 0;
/* open the door for the runners */
e->barrier_count = -e->barrier_count;
if ( e->nr_runners == 1 ) {
if (pthread_cond_signal(&e->barrier_cond) != 0)
return error(engine_err_pthread);
}
else {
if (pthread_cond_broadcast(&e->barrier_cond) != 0)
return error(engine_err_pthread);
}
/* wait for the runners to come home */
while (e->barrier_count < e->nr_runners)
if (pthread_cond_wait(&e->done_cond,&e->barrier_mutex) != 0)
return error(engine_err_pthread);
/* All in a days work. */
return engine_err_ok;
}
/**
* @brief Update the particle velocities and positions, re-shuffle if
* appropriate.
* @param e The #engine on which to run.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_advance ( struct engine *e ) {
int cid, pid, k, delta[3], step;
struct cell *c, *c_dest;
struct part *p;
struct space *s;
FPTYPE dt, w, h[3];
double epot = 0.0, epot_local;
/* Get a grip on the space. */
s = &(e->s);
dt = e->dt;
for ( k = 0 ; k < 3 ; k++ )
h[k] = s->h[k];
/* update the particle velocities and positions */
if ( e->flags & engine_flag_verlet || e->flags & engine_flag_mpi ) {
/* Collect potential energy from ghosts. */
for ( cid = 0 ; cid < s->nr_ghost ; cid++ )
epot += s->cells[ s->cid_ghost[cid] ].epot;
#pragma omp parallel private(cid,c,pid,p,w,k,epot_local)
{
step = omp_get_num_threads(); epot_local = 0.0;
for ( cid = omp_get_thread_num() ; cid < s->nr_real ; cid += step ) {
c = &(s->cells[ s->cid_real[cid] ]);
epot_local += c->epot;
for ( pid = 0 ; pid < c->count ; pid++ ) {
p = &( c->parts[pid] );
w = dt * e->types[p->type].imass;
for ( k = 0 ; k < 3 ; k++ ) {
p->v[k] += p->f[k] * w;
p->x[k] += dt * p->v[k];
}
}
}
#pragma omp atomic
epot += epot_local;
}
}
else {
/* Collect potential energy from ghosts. */
for ( cid = 0 ; cid < s->nr_ghost ; cid++ )
epot += s->cells[ s->cid_ghost[cid] ].epot;
#pragma omp parallel private(cid,c,pid,p,w,k,delta,c_dest,epot_local)
{
step = omp_get_num_threads(); epot_local = 0.0;
for ( cid = omp_get_thread_num() ; cid < s->nr_real ; cid += step ) {
c = &(s->cells[ s->cid_real[cid] ]);
epot_local += c->epot;
pid = 0;
while ( pid < c->count ) {
p = &( c->parts[pid] );
w = dt * e->types[p->type].imass;
for ( k = 0 ; k < 3 ; k++ ) {
p->v[k] += p->f[k] * w;
p->x[k] += dt * p->v[k];
delta[k] = __builtin_isgreaterequal( p->x[k] , h[k] ) - __builtin_isless( p->x[k] , 0.0 );
}
/* do we have to move this particle? */
if ( ( delta[0] != 0 ) || ( delta[1] != 0 ) || ( delta[2] != 0 ) ) {
for ( k = 0 ; k < 3 ; k++ )
p->x[k] -= delta[k] * h[k];
c_dest = &( s->cells[ space_cellid( s ,
(c->loc[0] + delta[0] + s->cdim[0]) % s->cdim[0] ,
(c->loc[1] + delta[1] + s->cdim[1]) % s->cdim[1] ,
(c->loc[2] + delta[2] + s->cdim[2]) % s->cdim[2] ) ] );
pthread_mutex_lock(&c_dest->cell_mutex);
cell_add_incomming( c_dest , p );
pthread_mutex_unlock(&c_dest->cell_mutex);
s->celllist[ p->id ] = c_dest;
c->count -= 1;
if ( pid < c->count ) {
c->parts[pid] = c->parts[c->count];
s->partlist[ c->parts[pid].id ] = &( c->parts[pid] );
}
}
else
pid += 1;
}
}
#pragma omp atomic
epot += epot_local;
}
/* Welcome the new particles in each cell. */
#pragma omp parallel for schedule(static)
for ( cid = 0 ; cid < s->nr_marked ; cid++ )
cell_welcome( &(s->cells[ s->cid_marked[cid] ]) , s->partlist );
}
/* Store the accumulated potential energy. */
s->epot_nonbond += epot;
s->epot += epot;
/* return quietly */
return engine_err_ok;
}
/**
* @brief Run the engine for a single time step.
*
* @param e The #engine on which to run.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* This routine advances the timestep counter by one, prepares the #space
* for a timestep, releases the #runner's associated with the #engine
* and waits for them to finnish.
*
* Once all the #runner's are done, the particle velocities and positions
* are updated and the particles are re-sorted in the #space.
*/
int engine_step ( struct engine *e ) {
int k;
ticks tic, tic_step = getticks();
/* increase the time stepper */
e->time += 1;
/* prepare the space */
tic = getticks();
if ( space_prepare( &e->s ) != space_err_ok )
return error(engine_err_space);
e->timers[engine_timer_prepare] += getticks() - tic;
/* Clear the runner timers. */
for ( k = 0 ; k < runner_timer_count ; k++ )
runner_timers[k] = 0;
/* Make sure the verlet lists are up to date. */
if ( e->flags & engine_flag_verlet ) {
/* Start the clock. */
tic = getticks();
/* Check particle movement and update cells if necessary. */
if ( engine_verlet_update( e ) < 0 )
return error(engine_err);
/* Store the timing. */
e->timers[engine_timer_verlet] += getticks() - tic;
}
/* Otherwise, if async MPI, move the particles accross the
node boundaries. */
else { // if ( e->flags & engine_flag_async ) {
tic = getticks();
if ( engine_shuffle( e ) < 0 )
return error(engine_err_space);
e->timers[engine_timer_advance] += getticks() - tic;
}
#ifdef WITH_MPI
/* Re-distribute the particles to the processors. */
if ( e->flags & engine_flag_mpi ) {
/* Start the clock. */
tic = getticks();
if ( e->flags & engine_flag_async ) {
if ( engine_exchange_async( e ) < 0 )
return error(engine_err);
}
else {
if ( engine_exchange( e ) < 0 )
return error(engine_err);
}
/* Store the timing. */
e->timers[engine_timer_exchange1] += getticks() - tic;
}
#endif
/* Compute the non-bonded interactions. */
tic = getticks();
#if defined(HAVE_CUDA) && defined(WITH_CUDA)
if ( e->flags & engine_flag_cuda ) {
if ( engine_nonbond_cuda( e ) < 0 )
return error(engine_err);
}
else
#endif
if ( engine_nonbond_eval( e ) < 0 )
return error(engine_err);
e->timers[engine_timer_nonbond] += getticks() - tic;
/* Clear the verlet-rebuild flag if it was set. */
if ( e->flags & engine_flag_verlet && e->s.verlet_rebuild )
e->s.verlet_rebuild = 0;
/* Do bonded interactions. */
if ( !( e->flags & engine_flag_parbonded ) ) {
tic = getticks();
if ( engine_bonded_eval( e ) < 0 )
return error(engine_err);
e->timers[engine_timer_bonded] += getticks() - tic;
}
else
e->timers[engine_timer_bonded] += runner_timers[ runner_timer_bonded ];
/* update the particle velocities and positions. */
tic = getticks();
if ( engine_advance( e ) < 0 )
return error(engine_err);
e->timers[engine_timer_advance] += getticks() - tic;
/* Shake the particle positions? */
if ( e->nr_rigids > 0 ) {
#ifdef WITH_MPI
/* If we have to do some communication first... */
if ( e->flags & engine_flag_mpi ) {
/* Sort the constraints. */
tic = getticks();
if ( engine_rigid_sort( e ) != 0 )
return error(engine_err);
e->timers[engine_timer_rigid] += getticks() - tic;
/* Start the clock. */
tic = getticks();
if ( e->flags & engine_flag_async ) {
if ( engine_exchange_rigid_async( e ) != 0 )
return error(engine_err);
}
else {
if ( engine_exchange_rigid( e ) != 0 )
return error(engine_err);
}
/* Store the timing. */
e->timers[engine_timer_exchange2] += getticks() - tic;
}
#endif
/* Resolve the constraints. */
tic = getticks();
if ( engine_rigid_eval( e ) != 0 )
return error(engine_err);
e->timers[engine_timer_rigid] += getticks() - tic;
}
/* Stop the clock. */
e->timers[engine_timer_step] += getticks() - tic_step;
/* return quietly */
return engine_err_ok;
}
/**
* @brief Barrier routine to hold the @c runners back.
*
* @param e The #engine to wait on.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*
* After being initialized, and after every timestep, every #runner
* calls this routine which blocks until all the runners have returned
* and the #engine signals the next timestep.
*/
int engine_barrier ( struct engine *e ) {
/* lock the barrier mutex */
if (pthread_mutex_lock(&e->barrier_mutex) != 0)
return error(engine_err_pthread);
/* wait for the barrier to close */
while (e->barrier_count < 0)
if (pthread_cond_wait(&e->barrier_cond,&e->barrier_mutex) != 0)
return error(engine_err_pthread);
/* if i'm the last thread in, signal that the barrier is full */
if (++e->barrier_count == e->nr_runners) {
if (pthread_cond_signal(&e->done_cond) != 0)
return error(engine_err_pthread);
}
/* wait for the barrier to re-open */
while (e->barrier_count > 0)
if (pthread_cond_wait(&e->barrier_cond,&e->barrier_mutex) != 0)
return error(engine_err_pthread);
/* if i'm the last thread out, signal to those waiting to get back in */
if (++e->barrier_count == 0)
if (pthread_cond_broadcast(&e->barrier_cond) != 0)
return error(engine_err_pthread);
/* free the barrier mutex */
if (pthread_mutex_unlock(&e->barrier_mutex) != 0)
return error(engine_err_pthread);
/* all is well... */
return engine_err_ok;
}
/**
* @brief Initialize an #engine with the given data and MPI enabled.
*
* @param e The #engine to initialize.
* @param origin An array of three doubles containing the cartesian origin
* of the space.
* @param dim An array of three doubles containing the size of the space.
* @param L The minimum cell edge length, should be at least @c cutoff.
* @param cutoff The maximum interaction cutoff to use.
* @param period A bitmask describing the periodicity of the domain
* (see #space_periodic_full).
* @param max_type The maximum number of particle types that will be used
* by this engine.
* @param flags Bit-mask containing the flags for this engine.
* @param comm The MPI comm to use.
* @param rank The ID of this node.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
#ifdef WITH_MPI
int engine_init_mpi ( struct engine *e , const double *origin , const double *dim , double *L , double cutoff , unsigned int period , int max_type , unsigned int flags , MPI_Comm comm , int rank ) {
/* Init the engine. */
if ( engine_init( e , origin , dim , L , cutoff , period , max_type , flags | engine_flag_mpi ) < 0 )
return error(engine_err);
/* Store the MPI Comm and rank. */
e->comm = comm;
e->nodeID = rank;
/* Bail. */
return engine_err_ok;
}
#endif
/**
* @brief Kill all runners and de-allocate the data of an engine.
*
* @param e the #engine to finalize.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_finalize ( struct engine *e ) {
int j, k;
/* make sure the inputs are ok */
if ( e == NULL )
return error(engine_err_null);
/* Shut down the runners, if they were started. */
if ( e->runners != NULL ) {
for ( k = 0 ; k < e->nr_runners ; k++ )
if ( pthread_cancel( e->runners[k].thread ) != 0 )
return error(engine_err_pthread);
free( e->runners );
free( e->queues );
}
/* Finalize the space. */
// if ( space_finalize( &e->s ) < 0 )
// return error(engine_err_space);
/* Free-up the types. */
free( e->types );
/* Free the potentials. */
if ( e->p != NULL ) {
for ( j = 0 ; j < e->nr_types ; j++ )
for ( k = j ; k < e->nr_types ; k++ ) {
if ( e->p[ j*e->max_type + k ] != NULL )
potential_clear( e->p[ j*e->max_type + k ] );
if ( e->p_bond[ j*e->max_type + k ] != NULL )
potential_clear( e->p_bond[ j*e->max_type + k ] );
}
for ( k = 0 ; k < e->nr_anglepots ; k++ )
potential_clear( e->p_angle[k] );
for ( k = 0 ; k < e->nr_dihedralpots ; k++ )
potential_clear( e->p_dihedral[k] );
free( e->p );
}
if ( e->p_bond != NULL )
free( e->p_bond );
if ( e->p_angle != NULL )
free( e->p_angle );
if ( e->p_dihedral != NULL )
free( e->p_dihedral );
/* Free the communicators, if needed. */
if ( e->flags & engine_flag_mpi ) {
for ( k = 0 ; k < e->nr_nodes ; k++ ) {
free( e->send[k].cellid );
free( e->recv[k].cellid );
}
free( e->send );
free( e->recv );
}
/* Free the bonded interactions. */
free( e->bonds );
free( e->angles );
free( e->dihedrals );
free( e->exclusions );
free( e->rigids );
free( e->part2rigid );
/* If we have bonded sets, kill them. */
for ( k = 0 ; k < e->nr_sets ; k++ ) {
free( e->sets[k].bonds );
free( e->sets[k].angles );
free( e->sets[k].dihedrals );
free( e->sets[k].exclusions );
free( e->sets[k].cells );
}
free( e->sets );
/* Clear all the counts and what not. */
bzero( e , sizeof( struct engine ) );
/* Happy and I know it... */
return engine_err_ok;
}
/**
* @brief Initialize an #engine with the given data.
*
* @param e The #engine to initialize.
* @param origin An array of three doubles containing the cartesian origin
* of the space.
* @param dim An array of three doubles containing the size of the space.
* @param L The minimum cell edge length in each dimension.
* @param cutoff The maximum interaction cutoff to use.
* @param period A bitmask describing the periodicity of the domain
* (see #space_periodic_full).
* @param max_type The maximum number of particle types that will be used
* by this engine.
* @param flags Bit-mask containing the flags for this engine.
*
* @return #engine_err_ok or < 0 on error (see #engine_err).
*/
int engine_init ( struct engine *e , const double *origin , const double *dim , double *L , double cutoff , unsigned int period , int max_type , unsigned int flags ) {
int cid;
/* make sure the inputs are ok */
if ( e == NULL || origin == NULL || dim == NULL || L == NULL )
return error(engine_err_null);
/* Check for bad flags. */
#ifdef FPTYPE_DOUBLE
if ( e->flags & engine_flag_cuda )
return error(engine_err_cudasp);
#endif
/* init the space with the given parameters */
if ( space_init( &(e->s) , origin , dim , L , cutoff , period ) < 0 )
return error(engine_err_space);
/* Set some flag implications. */
if ( flags & engine_flag_verlet_pseudo )
flags |= engine_flag_verlet_pairwise;
if ( flags & engine_flag_verlet_pairwise )
flags |= engine_flag_verlet;
if ( flags & engine_flag_cuda )
flags |= engine_flag_nullpart;
/* Set the flags. */
e->flags = flags;
/* By default there is only one node. */
e->nr_nodes = 1;
/* Init the timers. */
if ( engine_timers_reset( e ) < 0 )
return error(engine_err);
/* Init the runners to 0. */
e->runners = NULL;
e->nr_runners = 0;
/* Start with no queues. */
e->queues = NULL;
e->nr_queues = 0;
/* Init the bonds array. */
e->bonds_size = 100;
if ( ( e->bonds = (struct bond *)malloc( sizeof( struct bond ) * e->bonds_size ) ) == NULL )
return error(engine_err_malloc);
e->nr_bonds = 0;
/* Init the exclusions array. */
e->exclusions_size = 100;
if ( ( e->exclusions = (struct exclusion *)malloc( sizeof( struct exclusion ) * e->exclusions_size ) ) == NULL )
return error(engine_err_malloc);
e->nr_exclusions = 0;
/* Init the rigids array. */
e->rigids_size = 100;
if ( ( e->rigids = (struct rigid *)malloc( sizeof( struct rigid ) * e->rigids_size ) ) == NULL )
return error(engine_err_malloc);
e->nr_rigids = 0;
e->tol_rigid = 1e-6;
e->nr_constr = 0;
e->part2rigid = NULL;
/* Init the angles array. */
e->angles_size = 100;
if ( ( e->angles = (struct angle *)malloc( sizeof( struct angle ) * e->angles_size ) ) == NULL )
return error(engine_err_malloc);
e->nr_angles = 0;
/* Init the dihedrals array. */
e->dihedrals_size = 100;
if ( ( e->dihedrals = (struct dihedral *)malloc( sizeof( struct dihedral ) * e->dihedrals_size ) ) == NULL )
return error(engine_err_malloc);
e->nr_dihedrals = 0;
/* set the maximum nr of types */
if ( flags & engine_flag_nullpart )
max_type += 1;
e->max_type = max_type;
e->nr_types = 0;
if ( ( e->types = (struct part_type *)malloc( sizeof(struct part_type) * max_type ) ) == NULL )
return error(engine_err_malloc);
if ( flags & engine_flag_nullpart ) {
e->types[0].id = 0;
e->types[0].mass = 0.0;
e->types[0].imass = 0.0;
e->types[0].charge = 0.0;
e->types[0].eps = 0.0;
e->types[0].rmin = 0.0;
strcpy( e->types[0].name , "NULL" );
strcpy( e->types[0].name2 , "NULL" );
e->nr_types = 1;
}
/* Init the sets. */
e->sets = NULL;
e->nr_sets = 0;
/* allocate the interaction matrices */
if ( ( e->p = (struct potential **)malloc( sizeof(struct potential *) * max_type * max_type ) ) == NULL )
return error(engine_err_malloc);
bzero( e->p , sizeof(struct potential *) * max_type * max_type );
if ( (e->p_bond = (struct potential **)malloc( sizeof(struct potential *) * max_type * max_type )) == NULL)
return error(engine_err_malloc);
bzero( e->p_bond , sizeof(struct potential *) * max_type * max_type );
e->anglepots_size = 100;
if ( (e->p_angle = (struct potential **)malloc( sizeof(struct potential *) * e->anglepots_size )) == NULL)
return error(engine_err_malloc);
bzero( e->p_angle , sizeof(struct potential *) * e->anglepots_size );
e->nr_anglepots = 0;
e->dihedralpots_size = 100;
if ( (e->p_dihedral = (struct potential **)malloc( sizeof(struct potential *) * e->dihedralpots_size )) == NULL)
return error(engine_err_malloc);
bzero( e->p_dihedral , sizeof(struct potential *) * e->dihedralpots_size );
e->nr_dihedralpots = 0;
/* Make sortlists? */
if ( flags & engine_flag_verlet_pseudo ) {
for ( cid = 0 ; cid < e->s.nr_cells ; cid++ )
if ( e->s.cells[cid].flags & cell_flag_marked )
if ( ( e->s.cells[cid].sortlist = (unsigned int *)malloc( sizeof(unsigned int) * 13 * e->s.cells[cid].size ) ) == NULL )
return error(engine_err_malloc);
}
/* init the barrier variables */
e->barrier_count = 0;
if ( pthread_mutex_init( &e->barrier_mutex , NULL ) != 0 ||
pthread_cond_init( &e->barrier_cond , NULL ) != 0 ||
pthread_cond_init( &e->done_cond , NULL ) != 0)
return error(engine_err_pthread);
/* init the barrier */
if (pthread_mutex_lock(&e->barrier_mutex) != 0)
return error(engine_err_pthread);
e->barrier_count = 0;
/* Init the comm arrays. */
e->send = NULL;
e->recv = NULL;
/* all is well... */
return engine_err_ok;
}