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
// *****************************************************************************
/*!
  \file      src/PDE/Limiter.cpp
  \copyright 2012-2015 J. Bakosi,
             2016-2018 Los Alamos National Security, LLC.,
             2019-2021 Triad National Security, LLC.
             All rights reserved. See the LICENSE file for details.
  \brief     Limiters for discontiunous Galerkin methods
  \details   This file contains functions that provide limiter function
    calculations for maintaining monotonicity near solution discontinuities
    for the DG discretization.
*/
// *****************************************************************************

#include <array>
#include <vector>

#include "FaceData.hpp"
#include "Vector.hpp"
#include "Limiter.hpp"
#include "DerivedData.hpp"
#include "Integrate/Quadrature.hpp"
#include "Integrate/Basis.hpp"
#include "Inciter/InputDeck/InputDeck.hpp"
#include "PrefIndicator.hpp"
#include "Reconstruction.hpp"
#include "Integrate/Mass.hpp"
#include "MultiMat/MiscMultiMatFns.hpp"
#include "MultiSpecies/MultiSpeciesIndexing.hpp"

namespace inciter {

extern ctr::InputDeck g_inputdeck;

void
WENO_P1( const std::vector< int >& esuel,
         tk::Fields& U )
// *****************************************************************************
//  Weighted Essentially Non-Oscillatory (WENO) limiter for DGP1
//! \param[in] esuel Elements surrounding elements
//! \param[in,out] U High-order solution vector which gets limited
//! \details This WENO function should be called for transport and compflow
//! \note This limiter function is experimental and untested. Use with caution.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto cweight = inciter::g_inputdeck.get< tag::cweight >();
  auto nelem = esuel.size()/4;
  std::array< std::vector< tk::real >, 3 >
    limU {{ std::vector< tk::real >(nelem),
            std::vector< tk::real >(nelem),
            std::vector< tk::real >(nelem) }};

  std::size_t ncomp = U.nprop()/rdof;

  for (inciter::ncomp_t c=0; c<ncomp; ++c)
  {
    for (std::size_t e=0; e<nelem; ++e)
    {
      WENOLimiting(U, esuel, e, c, rdof, cweight, limU);
    }

    auto mark = c*rdof;

    for (std::size_t e=0; e<nelem; ++e)
    {
      U(e, mark+1) = limU[0][e];
      U(e, mark+2) = limU[1][e];
      U(e, mark+3) = limU[2][e];
    }
  }
}

void
Superbee_P1( const std::vector< int >& esuel,
             const std::vector< std::size_t >& inpoel,
             const std::vector< std::size_t >& ndofel,
             const tk::UnsMesh::Coords& coord,
             tk::Fields& U )
// *****************************************************************************
//  Superbee limiter for DGP1
//! \param[in] esuel Elements surrounding elements
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] coord Array of nodal coordinates
//! \param[in,out] U High-order solution vector which gets limited
//! \details This Superbee function should be called for transport and compflow
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  std::size_t ncomp = U.nprop()/rdof;

  auto beta_lim = 2.0;

  for (std::size_t e=0; e<esuel.size()/4; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1)
    {
      auto phi = SuperbeeLimiting(U, esuel, inpoel, coord, e, ndof, rdof,
                   dof_el, ncomp, beta_lim);

      // apply limiter function
      for (inciter::ncomp_t c=0; c<ncomp; ++c)
      {
        auto mark = c*rdof;
        U(e, mark+1) = phi[c] * U(e, mark+1);
        U(e, mark+2) = phi[c] * U(e, mark+2);
        U(e, mark+3) = phi[c] * U(e, mark+3);
      }
    }
  }
}

void
SuperbeeMultiMat_P1(
  const std::vector< int >& esuel,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  const tk::UnsMesh::Coords& coord,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,
  tk::Fields& P,
  std::size_t nmat )
// *****************************************************************************
//  Superbee limiter for multi-material DGP1
//! \param[in] esuel Elements surrounding elements
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] coord Array of nodal coordinates
//! \param[in] solidx Solid material index indicator
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] P High-order vector of primitives which gets limited
//! \param[in] nmat Number of materials in this PDE system
//! \details This Superbee function should be called for multimat
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  const auto intsharp = inciter::g_inputdeck.get< tag::multimat,
    tag::intsharp >();
  std::size_t ncomp = U.nprop()/rdof;
  std::size_t nprim = P.nprop()/rdof;

  auto beta_lim = 2.0;

  for (std::size_t e=0; e<esuel.size()/4; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1)
    {
      // limit conserved quantities
      auto phic = SuperbeeLimiting(U, esuel, inpoel, coord, e, ndof, rdof,
                    dof_el, ncomp, beta_lim);
      // limit primitive quantities
      auto phip = SuperbeeLimiting(P, esuel, inpoel, coord, e, ndof, rdof,
                    dof_el, nprim, beta_lim);

      std::vector< tk::real > phic_p2;
      std::vector< std::vector< tk::real > > unk, prim;<--- Unused variable: unk<--- Unused variable: prim
      if(ndof > 1)
        BoundPreservingLimiting(nmat, ndof, e, inpoel, coord, U, phic,
          phic_p2);

      // limits under which compression is to be performed
      std::vector< std::size_t > matInt(nmat, 0);
      std::vector< tk::real > alAvg(nmat, 0.0);
      for (std::size_t k=0; k<nmat; ++k)
        alAvg[k] = U(e, volfracDofIdx(nmat,k,rdof,0));
      auto intInd = interfaceIndicator(nmat, alAvg, matInt);
      if ((intsharp > 0) && intInd)
      {
        for (std::size_t k=0; k<nmat; ++k)
        {
          if (matInt[k])
            phic[volfracIdx(nmat,k)] = 1.0;
        }
      }
      else
      {
        if (!g_inputdeck.get< tag::accuracy_test >())
          consistentMultiMatLimiting_P1(nmat, rdof, e, solidx, U, P, phic,
            phic_p2);
      }

      // apply limiter function
      for (inciter::ncomp_t c=0; c<ncomp; ++c)
      {
        auto mark = c*rdof;
        U(e, mark+1) = phic[c] * U(e, mark+1);
        U(e, mark+2) = phic[c] * U(e, mark+2);
        U(e, mark+3) = phic[c] * U(e, mark+3);
      }
      for (inciter::ncomp_t c=0; c<nprim; ++c)
      {
        auto mark = c*rdof;
        P(e, mark+1) = phip[c] * P(e, mark+1);
        P(e, mark+2) = phip[c] * P(e, mark+2);
        P(e, mark+3) = phip[c] * P(e, mark+3);
      }
    }
  }
}

void
VertexBasedTransport_P1(
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  std::size_t nelem,
  const tk::UnsMesh::Coords& coord,
  tk::Fields& U )
// *****************************************************************************
//  Kuzmin's vertex-based limiter for transport DGP1
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] nelem Number of elements
//! \param[in] coord Array of nodal coordinates
//! \param[in,out] U High-order solution vector which gets limited
//! \details This vertex-based limiter function should be called for transport.
//!   For details see: Kuzmin, D. (2010). A vertex-based hierarchical slope
//!   limiter for p-adaptive discontinuous Galerkin methods. Journal of
//!   computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  const auto intsharp = inciter::g_inputdeck.get< tag::transport,
    tag::intsharp >();
  std::size_t ncomp = U.nprop()/rdof;

  for (std::size_t e=0; e<nelem; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1)
    {
      std::vector< tk::real > phi(ncomp, 1.0);
      std::vector< std::size_t > var;
      for (std::size_t c=0; c<ncomp; ++c) var.push_back(c);
      // limit conserved quantities
      VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
        ncomp, phi, var);

      // limits under which compression is to be performed
      std::vector< std::size_t > matInt(ncomp, 0);
      std::vector< tk::real > alAvg(ncomp, 0.0);
      for (std::size_t k=0; k<ncomp; ++k)
        alAvg[k] = U(e,k*rdof);
      auto intInd = interfaceIndicator(ncomp, alAvg, matInt);
      if ((intsharp > 0) && intInd)
      {
        for (std::size_t k=0; k<ncomp; ++k)
        {
          if (matInt[k]) phi[k] = 1.0;
        }
      }

      // apply limiter function
      for (std::size_t c=0; c<ncomp; ++c)
      {
        auto mark = c*rdof;
        U(e, mark+1) = phi[c] * U(e, mark+1);
        U(e, mark+2) = phi[c] * U(e, mark+2);
        U(e, mark+3) = phi[c] * U(e, mark+3);
      }
    }
  }
}

void
VertexBasedCompflow_P1(
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  std::size_t nelem,
  const std::vector< inciter::EOS >& mat_blk,
  const inciter::FaceData& fd,
  const tk::Fields& geoFace,
  const tk::Fields& geoElem,
  const tk::UnsMesh::Coords& coord,
  const tk::FluxFn& flux,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,
  std::vector< std::size_t >& shockmarker )
// *****************************************************************************
//  Kuzmin's vertex-based limiter for single-material DGP1
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] nelem Number of elements
//! \param[in] mat_blk EOS material block
//! \param[in] fd Face connectivity and boundary conditions object
//! \param[in] geoFace Face geometry array
// //! \param[in] geoElem Element geometry array
//! \param[in] coord Array of nodal coordinates
//! \param[in] flux Riemann flux function to use
//! \param[in] solidx Solid material index indicator
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] shockmarker Shock detection marker array
//! \details This vertex-based limiter function should be called for compflow.
//!   For details see: Kuzmin, D. (2010). A vertex-based hierarchical slope
//!   limiter for p-adaptive discontinuous Galerkin methods. Journal of
//!   computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  std::size_t ncomp = U.nprop()/rdof;

  // Null field for MarkShockCells argument
  tk::Fields P;

  if (inciter::g_inputdeck.get< tag::shock_detector_coeff >()
    > 1e-6)
    MarkShockCells(false, nelem, 1, ndof, rdof, mat_blk, ndofel,
      inpoel, coord, fd, geoFace, geoElem, flux, solidx, U, P, shockmarker);

  for (std::size_t e=0; e<nelem; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1 && shockmarker[e])
    {
      std::vector< tk::real > phi(ncomp, 1.0);
      std::vector< std::size_t > var;
      for (std::size_t c=0; c<ncomp; ++c) var.push_back(c);
      // limit conserved quantities
      VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
        ncomp, phi, var);

      // apply limiter function
      for (std::size_t c=0; c<ncomp; ++c)
      {
        auto mark = c*rdof;
        U(e, mark+1) = phi[c] * U(e, mark+1);
        U(e, mark+2) = phi[c] * U(e, mark+2);
        U(e, mark+3) = phi[c] * U(e, mark+3);
      }
    }
  }
}

void
VertexBasedCompflow_P2(
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  std::size_t nelem,
  const std::vector< inciter::EOS >& mat_blk,
  const inciter::FaceData& fd,
  const tk::Fields& geoFace,
  const tk::Fields& geoElem,
  const tk::UnsMesh::Coords& coord,
  [[maybe_unused]] const std::vector< std::size_t >& gid,
  [[maybe_unused]] const std::unordered_map< std::size_t, std::size_t >& bid,
  [[maybe_unused]] const std::vector< std::vector<tk::real> >& uNodalExtrm,
  [[maybe_unused]] const std::vector< std::vector<tk::real> >& mtInv,
  const tk::FluxFn& flux,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,
  std::vector< std::size_t >& shockmarker )
// *****************************************************************************
//  Kuzmin's vertex-based limiter on reference element for single-material DGP2
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] nelem Number of elements
//! \param[in] mat_blk EOS material block
//! \param[in] fd Face connectivity and boundary conditions object
//! \param[in] geoFace Face geometry array
// //! \param[in] geoElem Element geometry array
//! \param[in] coord Array of nodal coordinates
//! \param[in] gid Local->global node id map
//! \param[in] bid Local chare-boundary node ids (value) associated to
//!   global node ids (key)
//! \param[in] uNodalExtrm Chare-boundary nodal extrema for conservative
//!   variables
//! \param[in] mtInv Inverse of Taylor mass matrix
//! \param[in] flux Riemann flux function to use
//! \param[in] solidx Solid material index indicator
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] shockmarker Shock detection marker array
//! \details This vertex-based limiter function should be called for compflow.
//!   For details see: Kuzmin, D. (2010). A vertex-based hierarchical slope
//!   limiter for p-adaptive discontinuous Galerkin methods. Journal of
//!   computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  std::size_t ncomp = U.nprop()/rdof;

  // Null field for MarkShockCells argument
  tk::Fields P;

  if (inciter::g_inputdeck.get< tag::shock_detector_coeff >()
    > 1e-6)
    MarkShockCells(false, nelem, 1, ndof, rdof, mat_blk, ndofel,
      inpoel, coord, fd, geoFace, geoElem, flux, solidx, U, P, shockmarker);

  for (std::size_t e=0; e<nelem; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1 && shockmarker[e])
    {
      // The vector of limiting coefficients for P1 and P2 coefficients
      std::vector< tk::real > phic_p1(ncomp, 1.0);

      // Removing 3rd order DOFs if discontinuity is detected, and applying
      // limiting to the 2nd order/P1 solution
      for (std::size_t c=0; c<ncomp; ++c) {
        for(std::size_t idof = 4; idof < rdof; idof++) {
          auto mark = c * rdof + idof;
          U(e, mark) = 0.0;
        }
      }

      // Obtain limiting coefficient for P1 coefficients
      std::vector< std::size_t > var;
      for (std::size_t c=0; c<ncomp; ++c) var.push_back(c);
      VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
        ncomp, phic_p1, var);

      // apply limiter function to the solution with Taylor basis
      for (std::size_t c=0; c<ncomp; ++c) {
        auto mark = c * rdof;
        for(std::size_t idof=1; idof<4; idof++)
          U(e, mark+idof) = phic_p1[c] * U(e, mark+idof);
      }
    }
  }
}

void
VertexBasedMultiMat_P1(
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  std::size_t nelem,
  const std::vector< inciter::EOS >& mat_blk,
  const inciter::FaceData& fd,
  const tk::Fields& geoFace,
  const tk::Fields& geoElem,
  const tk::UnsMesh::Coords& coord,
  const tk::FluxFn& flux,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,
  tk::Fields& P,
  std::size_t nmat,
  std::vector< std::size_t >& shockmarker )
// *****************************************************************************
//  Kuzmin's vertex-based limiter for multi-material DGP1
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] nelem Number of elements
//! \param[in] mat_blk EOS material block
//! \param[in] fd Face connectivity and boundary conditions object
//! \param[in] geoFace Face geometry array
// //! \param[in] geoElem Element geometry array
//! \param[in] coord Array of nodal coordinates
//! \param[in] flux Riemann flux function to use
//! \param[in] solidx Solid material index indicator
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] P High-order vector of primitives which gets limited
//! \param[in] nmat Number of materials in this PDE system
//! \param[in,out] shockmarker Shock detection marker array
//! \details This vertex-based limiter function should be called for multimat.
//!   For details see: Kuzmin, D. (2010). A vertex-based hierarchical slope
//!   limiter for p-adaptive discontinuous Galerkin methods. Journal of
//!   computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  const auto intsharp = inciter::g_inputdeck.get< tag::multimat,
    tag::intsharp >();
  std::size_t ncomp = U.nprop()/rdof;
  std::size_t nprim = P.nprop()/rdof;

  // Evaluate the interface condition and mark the shock cells
  if (inciter::g_inputdeck.get< tag::shock_detector_coeff >()
    > 1e-6 && ndof > 1)
    MarkShockCells(false, nelem, nmat, ndof, rdof, mat_blk, ndofel,
      inpoel, coord, fd, geoFace, geoElem, flux, solidx, U, P, shockmarker);

  for (std::size_t e=0; e<nelem; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1)
    {
      std::vector< tk::real > phic(ncomp, 1.0);
      std::vector< tk::real > phip(nprim, 1.0);
      if(shockmarker[e]) {
        // When shockmarker is 1, there is discontinuity within the element.
        // Hence, the vertex-based limiter will be applied.

        // limit conserved quantities
        std::vector< std::size_t > varc;
        for (std::size_t c=0; c<ncomp; ++c) varc.push_back(c);
        VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
          ncomp, phic, varc);
        // limit primitive quantities
        std::vector< std::size_t > varp;
        for (std::size_t c=0; c<nprim; ++c) varp.push_back(c);
        VertexBasedLimiting(P, esup, inpoel, coord, e, rdof, dof_el,
          nprim, phip, varp);
      } else {
        // When shockmarker is 0, the volume fraction, density and energy
        // of minor material will still be limited to ensure a stable solution.
        std::vector< std::size_t > vars;
        for (std::size_t k=0; k<nmat; ++k) vars.push_back(volfracIdx(nmat,k));
        VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
          ncomp, phic, vars);

        for(std::size_t k=0; k<nmat; ++k) {
          if(U(e, volfracDofIdx(nmat,k,rdof,0)) < 1e-4) {
            // limit the density and energy of minor materials
            vars.clear();
            vars.push_back(densityIdx(nmat, k));
            vars.push_back(energyIdx(nmat, k));
            if (solidx[k] > 0) {
              for (std::size_t i=0; i<3; ++i)
                for (std::size_t j=0; j<3; ++j)
                  vars.push_back(deformIdx(nmat, solidx[k], i, j));
            }
            VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
              ncomp, phic, vars);

            // limit the pressure of minor materials
            VertexBasedLimiting(P, esup, inpoel, coord, e, rdof, dof_el,
              nprim, phip, std::vector< std::size_t >{pressureIdx(nmat, k)});
          }
        }
      }

      std::vector< tk::real > phic_p2, phip_p2;

      PositivityLimitingMultiMat(nmat, mat_blk, rdof, dof_el, ndofel, e, inpoel,
        coord, fd.Esuel(), U, P, phic, phic_p2, phip, phip_p2);

      // limits under which compression is to be performed
      std::vector< std::size_t > matInt(nmat, 0);
      std::vector< tk::real > alAvg(nmat, 0.0);
      for (std::size_t k=0; k<nmat; ++k)
        alAvg[k] = U(e, volfracDofIdx(nmat,k,rdof,0));
      auto intInd = interfaceIndicator(nmat, alAvg, matInt);
      if ((intsharp > 0) && intInd) {
        for (std::size_t k=0; k<nmat; ++k) {
          if (matInt[k]) {
            phic[volfracIdx(nmat,k)] = 1.0;
          }
        }
      }
      else {
        if(nmat > 1)
          BoundPreservingLimiting(nmat, ndof, e, inpoel, coord, U, phic,
            phic_p2);

        if (!g_inputdeck.get< tag::accuracy_test >())
          consistentMultiMatLimiting_P1(nmat, rdof, e, solidx, U, P, phic,
            phic_p2);
      }

      // apply limiter function
      for (std::size_t c=0; c<ncomp; ++c)
      {
        auto mark = c*rdof;
        U(e, mark+1) = phic[c] * U(e, mark+1);
        U(e, mark+2) = phic[c] * U(e, mark+2);
        U(e, mark+3) = phic[c] * U(e, mark+3);
      }
      for (std::size_t c=0; c<nprim; ++c)
      {
        auto mark = c*rdof;
        P(e, mark+1) = phip[c] * P(e, mark+1);
        P(e, mark+2) = phip[c] * P(e, mark+2);
        P(e, mark+3) = phip[c] * P(e, mark+3);
      }
    }
  }
}

void
VertexBasedMultiMat_P2(
  const bool pref,
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  std::size_t nelem,
  const std::vector< inciter::EOS >& mat_blk,
  const inciter::FaceData& fd,
  const tk::Fields& geoFace,
  const tk::Fields& geoElem,
  const tk::UnsMesh::Coords& coord,
  [[maybe_unused]] const std::vector< std::size_t >& gid,
  [[maybe_unused]] const std::unordered_map< std::size_t, std::size_t >& bid,
  [[maybe_unused]] const std::vector< std::vector<tk::real> >& uNodalExtrm,
  [[maybe_unused]] const std::vector< std::vector<tk::real> >& pNodalExtrm,
  [[maybe_unused]] const std::vector< std::vector<tk::real> >& mtInv,
  const tk::FluxFn& flux,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,
  tk::Fields& P,
  std::size_t nmat,
  std::vector< std::size_t >& shockmarker )
// *****************************************************************************
//  Kuzmin's vertex-based limiter for multi-material DGP2
//! \param[in] pref Indicator for p-adaptive algorithm
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] nelem Number of elements
//! \param[in] mat_blk EOS material block
//! \param[in] fd Face connectivity and boundary conditions object
//! \param[in] geoFace Face geometry array
// //! \param[in] geoElem Element geometry array
//! \param[in] coord Array of nodal coordinates
//! \param[in] gid Local->global node id map
//! \param[in] bid Local chare-boundary node ids (value) associated to
//!   global node ids (key)
//! \param[in] uNodalExtrm Chare-boundary nodal extrema for conservative
//!   variables
//! \param[in] pNodalExtrm Chare-boundary nodal extrema for primitive
//!   variables
//! \param[in] mtInv Inverse of Taylor mass matrix
//! \param[in] flux Riemann flux function to use
//! \param[in] solidx Solid material index indicator
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] P High-order vector of primitives which gets limited
//! \param[in] nmat Number of materials in this PDE system
//! \param[in,out] shockmarker Shock detection marker array
//! \details This vertex-based limiter function should be called for multimat.
//!   For details see: Kuzmin, D. (2010). A vertex-based hierarchical slope
//!   limiter for p-adaptive discontinuous Galerkin methods. Journal of
//!   computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  const auto intsharp = inciter::g_inputdeck.get< tag::multimat,
    tag::intsharp >();
  std::size_t ncomp = U.nprop()/rdof;
  std::size_t nprim = P.nprop()/rdof;

  // Evaluate the interface condition and mark the shock cells
  if (inciter::g_inputdeck.get< tag::shock_detector_coeff >()
    > 1e-6)
    MarkShockCells(pref, nelem, nmat, ndof, rdof, mat_blk, ndofel,
      inpoel, coord, fd, geoFace, geoElem, flux, solidx, U, P, shockmarker);

  for (std::size_t e=0; e<nelem; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    // For multi-material simulation, when dofel = 1, p0p1 is applied and ndof
    // for solution evaluation should be 4
    if(ncomp > 5 && dof_el == 1)
      dof_el = 4;

    if (dof_el > 1)
    {
      // The vector of limiting coefficients for P1
      std::vector< tk::real > phic_p1(ncomp, 1.0), phic_p2(ncomp, 1.0);
      std::vector< tk::real > phip_p1(nprim, 1.0), phip_p2(nprim, 1.0);

      // Only when the cell is marked with discontinuous solution or P0P1 scheme
      // is used, the vertex-based slope limiter will be applied.
      if(shockmarker[e] || dof_el == 4) {
        // Removing 3rd order DOFs if discontinuity is detected, and applying
        // limiting to the 2nd order/P1 solution
        for (std::size_t c=0; c<ncomp; ++c) {
          auto mark = c * rdof;
          for(std::size_t idof = 4; idof < rdof; idof++)
            U(e, mark+idof) = 0.0;
        }
        for (std::size_t c=0; c<nprim; ++c) {
          auto mark = c * rdof;
          for(std::size_t idof = 4; idof < rdof; idof++)
            P(e, mark+idof) = 0.0;
        }

        // Obtain limiter coefficient for P1 conserved quantities
        std::vector< std::size_t > varc;
        for (std::size_t c=0; c<ncomp; ++c) varc.push_back(c);
        VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
          ncomp, phic_p1, varc);
        // Obtain limiter coefficient for P1 primitive quantities
        std::vector< std::size_t > varp;
        for (std::size_t c=0; c<nprim; ++c) varp.push_back(c);
        VertexBasedLimiting(P, esup, inpoel, coord, e, rdof, dof_el,
          nprim, phip_p1, varp);
      } else {
        // When shockmarker is 0, the volume fraction will still be limited to
        // ensure a stable solution. Since the limiting strategy for third order
        // solution will downgrade the accuracy to second order, the density,
        // energy and pressure of minor material will not be limited.
        std::vector< std::size_t > vars;
        for (std::size_t k=0; k<nmat; ++k) vars.push_back(volfracIdx(nmat,k));
        VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
          ncomp, phic_p1, vars);

        //for(std::size_t k=0; k<nmat; ++k) {
        //  if(U(e, volfracDofIdx(nmat,k,rdof,0)) < 1e-4) {
        //    // limit the density of minor materials
        //    VertexBasedLimiting(unk, U, esup, inpoel, coord, e, rdof, dof_el,
        //      ncomp, phic_p1, std::vector< std::size_t >{densityIdx(nmat,k)});

        //    // limit the pressure of minor materials
        //    VertexBasedLimiting(prim, P, esup, inpoel, coord, e, rdof, dof_el,
        //      nprim, phip_p1, std::vector< std::size_t >{pressureIdx(nmat,k)});
        //  }
        //}
      }

      PositivityLimitingMultiMat(nmat, mat_blk, ndof, dof_el, ndofel, e, inpoel,
          coord, fd.Esuel(), U, P, phic_p1, phic_p2, phip_p1, phic_p2);

      // limits under which compression is to be performed
      std::vector< std::size_t > matInt(nmat, 0);
      std::vector< tk::real > alAvg(nmat, 0.0);
      for (std::size_t k=0; k<nmat; ++k)
        alAvg[k] = U(e, volfracDofIdx(nmat,k,rdof,0));
      auto intInd = interfaceIndicator(nmat, alAvg, matInt);
      if ((intsharp > 0) && intInd) {
        for (std::size_t k=0; k<nmat; ++k) {
          if (matInt[k]) {
            phic_p1[volfracIdx(nmat,k)] = 1.0;
          }
        }
      }
      else {
        if(nmat > 1)
          BoundPreservingLimiting(nmat, ndof, e, inpoel, coord, U,
            phic_p1, phic_p2);

        if (!g_inputdeck.get< tag::accuracy_test >())
          consistentMultiMatLimiting_P1(nmat, rdof, e, solidx, U, P, phic_p1,
            phic_p2);
      }

      // apply limiing coefficient
      for (std::size_t c=0; c<ncomp; ++c)
      {
        auto mark = c * rdof;
        for(std::size_t idof=1; idof<4; idof++)
          U(e, mark+idof) = phic_p1[c] * U(e, mark+idof);
        for(std::size_t idof=4; idof<rdof; idof++)
          U(e, mark+idof) = phic_p2[c] * U(e, mark+idof);
      }
      for (std::size_t c=0; c<nprim; ++c)
      {
        auto mark = c * rdof;
        for(std::size_t idof=1; idof<4; idof++)
          P(e, mark+idof) = phip_p1[c] * P(e, mark+idof);
        for(std::size_t idof=4; idof<rdof; idof++)
          P(e, mark+idof) = phip_p2[c] * P(e, mark+idof);
      }
    }
  }
}

void
VertexBasedMultiMat_FV(
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  std::size_t nelem,
  const tk::UnsMesh::Coords& coord,
  const std::vector< int >& srcFlag,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,
  tk::Fields& P,
  std::size_t nmat )
// *****************************************************************************
//  Kuzmin's vertex-based limiter for multi-material FV
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] nelem Number of elements
//! \param[in] coord Array of nodal coordinates
//! \param[in] srcFlag Whether the energy source was added
//! \param[in] solidx Solid material index indicator
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] P High-order vector of primitives which gets limited
//! \param[in] nmat Number of materials in this PDE system
//! \details This vertex-based limiter function should be called for multimat.
//!   For details see: Kuzmin, D. (2010). A vertex-based hierarchical slope
//!   limiter for p-adaptive discontinuous Galerkin methods. Journal of
//!   computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto intsharp = inciter::g_inputdeck.get< tag::multimat,
    tag::intsharp >();
  std::size_t ncomp = U.nprop()/rdof;
  std::size_t nprim = P.nprop()/rdof;

  for (std::size_t e=0; e<nelem; ++e)
  {
    std::vector< tk::real > phic(ncomp, 1.0);
    std::vector< tk::real > phip(nprim, 1.0);
    // limit conserved quantities
    std::vector< std::size_t > var;
    for (std::size_t k=0; k<nmat; ++k) {
      var.push_back(volfracIdx(nmat,k));
      var.push_back(densityIdx(nmat,k));
    }
    VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, rdof, ncomp,
      phic, var);
    // limit primitive quantities
    var.clear();
    for (std::size_t c=0; c<nprim; ++c) var.push_back(c);
    VertexBasedLimiting(P, esup, inpoel, coord, e, rdof, rdof, nprim,
      phip, var);

    // limits under which compression is to be performed
    std::vector< std::size_t > matInt(nmat, 0);
    std::vector< tk::real > alAvg(nmat, 0.0);
    for (std::size_t k=0; k<nmat; ++k)
      alAvg[k] = U(e, volfracDofIdx(nmat,k,rdof,0));
    auto intInd = interfaceIndicator(nmat, alAvg, matInt);
    if ((intsharp > 0) && intInd && srcFlag[e] == 0)
    {
      for (std::size_t k=0; k<nmat; ++k)
      {
        if (matInt[k])
          phic[volfracIdx(nmat,k)] = 1.0;
      }
    }
    else
    {
      if (!g_inputdeck.get< tag::accuracy_test >()) {
        std::vector< tk::real > phic_p2(ncomp, 1.0);
        consistentMultiMatLimiting_P1(nmat, rdof, e, solidx, U, P, phic,
          phic_p2);
      }
    }

    // apply limiter function
    for (std::size_t c=0; c<ncomp; ++c)
    {
      auto mark = c*rdof;
      U(e, mark+1) = phic[c] * U(e, mark+1);
      U(e, mark+2) = phic[c] * U(e, mark+2);
      U(e, mark+3) = phic[c] * U(e, mark+3);
    }
    for (std::size_t c=0; c<nprim; ++c)
    {
      auto mark = c*rdof;
      P(e, mark+1) = phip[c] * P(e, mark+1);
      P(e, mark+2) = phip[c] * P(e, mark+2);
      P(e, mark+3) = phip[c] * P(e, mark+3);
    }
  }
}

void
VertexBasedMultiSpecies_P1(
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const std::vector< std::size_t >& ndofel,
  std::size_t nelem,
  const std::vector< inciter::EOS >& /*mat_blk*/,
  const inciter::FaceData& /*fd*/,
  const tk::Fields& /*geoFace*/,
  const tk::Fields& /*geoElem*/,
  const tk::UnsMesh::Coords& coord,
  const tk::FluxFn& /*flux*/,
  tk::Fields& U,
  std::size_t nspec,
  std::vector< std::size_t >& shockmarker )<--- Parameter 'shockmarker' can be declared with const
// *****************************************************************************
//  Kuzmin's vertex-based limiter for multi-species DGP1
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] ndofel Vector of local number of degrees of freedom
//! \param[in] nelem Number of elements
// //! \param[in] mat_blk EOS material block
// //! \param[in] fd Face connectivity and boundary conditions object
// //! \param[in] geoFace Face geometry array
// //! \param[in] geoElem Element geometry array
//! \param[in] coord Array of nodal coordinates
// //! \param[in] flux Riemann flux function to use
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in] nspec Number of species in this PDE system
//! \param[in,out] shockmarker Shock detection marker array
//! \details This vertex-based limiter function should be called for
//!   multispecies. For details see: Kuzmin, D. (2010). A vertex-based
//!   hierarchical slope limiter for p-adaptive discontinuous Galerkin methods.
//!   Journal of computational and applied mathematics, 233(12), 3077-3085.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ndof = inciter::g_inputdeck.get< tag::ndof >();
  std::size_t ncomp = U.nprop()/rdof;

  //// Evaluate the interface condition and mark the shock cells
  //if (inciter::g_inputdeck.get< tag::shock_detector_coeff >()
  //  > 1e-6 && ndof > 1)
  //  MarkShockCells(false, nelem, nmat, ndof, rdof, mat_blk, ndofel,
  //    inpoel, coord, fd, geoFace, geoElem, flux, solidx, U, P, shockmarker);

  for (std::size_t e=0; e<nelem; ++e)
  {
    // If an rDG method is set up (P0P1), then, currently we compute the P1
    // basis functions and solutions by default. This implies that P0P1 is
    // unsupported in the p-adaptive DG (PDG). This is a workaround until we
    // have rdofel, which is needed to distinguish between ndofs and rdofs per
    // element for pDG.
    std::size_t dof_el;
    if (rdof > ndof)
    {
      dof_el = rdof;
    }
    else
    {
      dof_el = ndofel[e];
    }

    if (dof_el > 1)
    {
      std::vector< std::size_t > vars;
      if(shockmarker[e]) {
        // When shockmarker is 1, there is discontinuity within the element.
        // Hence, the vertex-based limiter will be applied.

        for (std::size_t c=0; c<ncomp; ++c) vars.push_back(c);
      } else {
        // When shockmarker is 0, the density of minor species will still be
        // limited to ensure a stable solution.

        tk::real rhob(0.0);
        for(std::size_t k=0; k<nspec; ++k)
          rhob += U(e, multispecies::densityDofIdx(nspec,k,rdof,0));
        for(std::size_t k=0; k<nspec; ++k) {
          if (U(e, multispecies::densityDofIdx(nspec,k,rdof,0))/rhob < 1e-4) {
            // limit the density of minor species
            vars.push_back(multispecies::densityIdx(nspec, k));
          }
        }
      }

      std::vector< tk::real > phic(ncomp, 1.0);
      VertexBasedLimiting(U, esup, inpoel, coord, e, rdof, dof_el,
        ncomp, phic, vars);

      //std::vector< tk::real > phic_p2, phip_p2;

      //PositivityLimitingMultiMat(nmat, mat_blk, rdof, dof_el, ndofel, e, inpoel,
      //  coord, fd.Esuel(), U, P, phic, phic_p2, phip, phip_p2);

      // TODO: Unit sum of mass fractions is maintained by using common limiter
      // for all species densities. Investigate better approaches.
      if (!g_inputdeck.get< tag::accuracy_test >()) {
        tk::real phi_rhos_p1(1.0);
        for (std::size_t k=0; k<nspec; ++k)
          phi_rhos_p1 = std::min( phi_rhos_p1,
            phic[multispecies::densityIdx(nspec, k)] );
        // same limiter for all densities
        for (std::size_t k=0; k<nspec; ++k)
          phic[multispecies::densityIdx(nspec, k)] = phi_rhos_p1;
      }

      // apply limiter function
      for (std::size_t c=0; c<ncomp; ++c)
      {
        auto mark = c*rdof;
        U(e, mark+1) = phic[c] * U(e, mark+1);
        U(e, mark+2) = phic[c] * U(e, mark+2);
        U(e, mark+3) = phic[c] * U(e, mark+3);
      }
    }
  }
}

void
WENOLimiting( const tk::Fields& U,
              const std::vector< int >& esuel,
              std::size_t e,
              inciter::ncomp_t c,
              std::size_t rdof,
              tk::real cweight,
              std::array< std::vector< tk::real >, 3 >& limU )
// *****************************************************************************
//  WENO limiter function calculation for P1 dofs
//! \param[in] U High-order solution vector which is to be limited
//! \param[in] esuel Elements surrounding elements
//! \param[in] e Id of element whose solution is to be limited
//! \param[in] c Index of component which is to be limited
//! \param[in] rdof Maximum number of reconstructed degrees of freedom
//! \param[in] cweight Weight of the central stencil
//! \param[in,out] limU Limited gradients of component c
// *****************************************************************************
{
  std::array< std::array< tk::real, 3 >, 5 > gradu;
  std::array< tk::real, 5 > wtStencil, osc, wtDof;

  auto mark = c*rdof;

  // reset all stencil values to zero
  for (auto& g : gradu) g.fill(0.0);
  osc.fill(0);
  wtDof.fill(0);
  wtStencil.fill(0);

  // The WENO limiter uses solution data from the neighborhood in the form
  // of stencils to enforce non-oscillatory conditions. The immediate
  // (Von Neumann) neighborhood of a tetrahedral cell consists of the 4
  // cells that share faces with it. These are the 4 neighborhood-stencils
  // for the tetrahedron. The primary stencil is the tet itself. Weights are
  // assigned to these stencils, with the primary stencil usually assigned
  // the highest weight. The lower the primary/central weight, the more
  // dissipative the limiting effect. This central weight is usually problem
  // dependent. It is set higher for relatively weaker discontinuities, and
  // lower for stronger discontinuities.

  // primary stencil
  gradu[0][0] = U(e, mark+1);
  gradu[0][1] = U(e, mark+2);
  gradu[0][2] = U(e, mark+3);
  wtStencil[0] = cweight;

  // stencils from the neighborhood
  for (std::size_t is=1; is<5; ++is)
  {
    auto nel = esuel[ 4*e+(is-1) ];

    // ignore physical domain ghosts
    if (nel == -1)
    {
      gradu[is].fill(0.0);
      wtStencil[is] = 0.0;
      continue;
    }

    std::size_t n = static_cast< std::size_t >( nel );
    gradu[is][0] = U(n, mark+1);
    gradu[is][1] = U(n, mark+2);
    gradu[is][2] = U(n, mark+3);
    wtStencil[is] = 1.0;
  }

  // From these stencils, an oscillation indicator is calculated, which
  // determines the effective weights for the high-order solution DOFs.
  // These effective weights determine the contribution of each of the
  // stencils to the high-order solution DOFs of the current cell which are
  // being limited. If this indicator detects a large oscillation in the
  // solution of the current cell, it reduces the effective weight for the
  // central stencil contribution to its high-order DOFs. This results in
  // a more dissipative and well-behaved solution in the troubled cell.

  // oscillation indicators
  for (std::size_t is=0; is<5; ++is)
    osc[is] = std::sqrt( tk::dot(gradu[is], gradu[is]) );

  tk::real wtotal = 0;

  // effective weights for dofs
  for (std::size_t is=0; is<5; ++is)
  {
    // A small number (1.0e-8) is needed here to avoid dividing by a zero in
    // the case of a constant solution, where osc would be zero. The number
    // is not set to machine zero because it is squared, and a number
    // between 1.0e-8 to 1.0e-6 is needed.
    wtDof[is] = wtStencil[is] * pow( (1.0e-8 + osc[is]), -2 );
    wtotal += wtDof[is];
  }

  for (std::size_t is=0; is<5; ++is)
  {
    wtDof[is] = wtDof[is]/wtotal;
  }

  limU[0][e] = 0.0;
  limU[1][e] = 0.0;
  limU[2][e] = 0.0;

  // limiter function
  for (std::size_t is=0; is<5; ++is)
  {
    limU[0][e] += wtDof[is]*gradu[is][0];
    limU[1][e] += wtDof[is]*gradu[is][1];
    limU[2][e] += wtDof[is]*gradu[is][2];
  }
}

std::vector< tk::real >
SuperbeeLimiting( const tk::Fields& U,
                  const std::vector< int >& esuel,
                  const std::vector< std::size_t >& inpoel,
                  const tk::UnsMesh::Coords& coord,
                  std::size_t e,
                  std::size_t ndof,
                  std::size_t rdof,
                  std::size_t dof_el,
                  inciter:: ncomp_t ncomp,
                  tk::real beta_lim )
// *****************************************************************************
//  Superbee limiter function calculation for P1 dofs
//! \param[in] U High-order solution vector which is to be limited
//! \param[in] esuel Elements surrounding elements
//! \param[in] inpoel Element connectivity
//! \param[in] coord Array of nodal coordinates
//! \param[in] e Id of element whose solution is to be limited
//! \param[in] ndof Maximum number of degrees of freedom
//! \param[in] rdof Maximum number of reconstructed degrees of freedom
//! \param[in] dof_el Local number of degrees of freedom
//! \param[in] ncomp Number of scalar components in this PDE system
//! \param[in] beta_lim Parameter which is equal to 2 for Superbee and 1 for
//!   minmod limiter
//! \return phi Limiter function for solution in element e
// *****************************************************************************
{
  // Superbee is a TVD limiter, which uses min-max bounds that the
  // high-order solution should satisfy, to ensure TVD properties. For a
  // high-order method like DG, this involves the following steps:
  // 1. Find min-max bounds in the immediate neighborhood of cell.
  // 2. Calculate the Superbee function for all the points where solution
  //    needs to be reconstructed to (all quadrature points). From these,
  //    use the minimum value of the limiter function.

  std::vector< tk::real > uMin(ncomp, 0.0), uMax(ncomp, 0.0);

  for (inciter::ncomp_t c=0; c<ncomp; ++c)
  {
    auto mark = c*rdof;
    uMin[c] = U(e, mark);
    uMax[c] = U(e, mark);
  }

  // ----- Step-1: find min/max in the neighborhood
  for (std::size_t is=0; is<4; ++is)
  {
    auto nel = esuel[ 4*e+is ];

    // ignore physical domain ghosts
    if (nel == -1) continue;

    auto n = static_cast< std::size_t >( nel );
    for (inciter::ncomp_t c=0; c<ncomp; ++c)
    {
      auto mark = c*rdof;
      uMin[c] = std::min(uMin[c], U(n, mark));
      uMax[c] = std::max(uMax[c], U(n, mark));
    }
  }

  // ----- Step-2: loop over all quadrature points to get limiter function

  // to loop over all the quadrature points of all faces of element e,
  // coordinates of the quadrature points are needed.
  // Number of quadrature points for face integration
  auto ng = tk::NGfa(ndof);

  // arrays for quadrature points
  std::array< std::vector< tk::real >, 2 > coordgp;
  std::vector< tk::real > wgp;

  coordgp[0].resize( ng );
  coordgp[1].resize( ng );
  wgp.resize( ng );

  // get quadrature point weights and coordinates for triangle
  tk::GaussQuadratureTri( ng, coordgp, wgp );

  const auto& cx = coord[0];
  const auto& cy = coord[1];
  const auto& cz = coord[2];

  // Extract the element coordinates
  std::array< std::array< tk::real, 3>, 4 > coordel {{
    {{ cx[ inpoel[4*e  ] ], cy[ inpoel[4*e  ] ], cz[ inpoel[4*e  ] ] }},
    {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
    {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
    {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }} }};

  // Compute the determinant of Jacobian matrix
  auto detT =
    tk::Jacobian( coordel[0], coordel[1], coordel[2], coordel[3] );

  // initialize limiter function
  std::vector< tk::real > phi(ncomp, 1.0);
  for (std::size_t lf=0; lf<4; ++lf)
  {
    // Extract the face coordinates
    std::array< std::size_t, 3 > inpofa_l {{ inpoel[4*e+tk::lpofa[lf][0]],
                                             inpoel[4*e+tk::lpofa[lf][1]],
                                             inpoel[4*e+tk::lpofa[lf][2]] }};

    std::array< std::array< tk::real, 3>, 3 > coordfa {{
      {{ cx[ inpofa_l[0] ], cy[ inpofa_l[0] ], cz[ inpofa_l[0] ] }},
      {{ cx[ inpofa_l[1] ], cy[ inpofa_l[1] ], cz[ inpofa_l[1] ] }},
      {{ cx[ inpofa_l[2] ], cy[ inpofa_l[2] ], cz[ inpofa_l[2] ] }} }};

    // Gaussian quadrature
    for (std::size_t igp=0; igp<ng; ++igp)
    {
      // Compute the coordinates of quadrature point at physical domain
      auto gp = tk::eval_gp( igp, coordfa, coordgp );

      //Compute the basis functions
      auto B_l = tk::eval_basis( rdof,
            tk::Jacobian( coordel[0], gp, coordel[2], coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], gp, coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], coordel[2], gp ) / detT );

      auto state =
        tk::eval_state(ncomp, rdof, dof_el, e, U, B_l);

      Assert( state.size() == ncomp, "Size mismatch" );

      // compute the limiter function
      for (inciter::ncomp_t c=0; c<ncomp; ++c)
      {
        auto phi_gp = 1.0;
        auto mark = c*rdof;
        auto uNeg = state[c] - U(e, mark);
        if (uNeg > 1.0e-14)
        {
          uNeg = std::max(uNeg, 1.0e-08);
          phi_gp = std::min( 1.0, (uMax[c]-U(e, mark))/(2.0*uNeg) );
        }
        else if (uNeg < -1.0e-14)
        {
          uNeg = std::min(uNeg, -1.0e-08);
          phi_gp = std::min( 1.0, (uMin[c]-U(e, mark))/(2.0*uNeg) );
        }
        else
        {
          phi_gp = 1.0;
        }
        phi_gp = std::max( 0.0,
                           std::max( std::min(beta_lim*phi_gp, 1.0),
                                     std::min(phi_gp, beta_lim) ) );
        phi[c] = std::min( phi[c], phi_gp );
      }
    }
  }

  return phi;
}

void
VertexBasedLimiting(
  const tk::Fields& U,
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  const tk::UnsMesh::Coords& coord,
  std::size_t e,
  std::size_t rdof,
  std::size_t dof_el,
  std::size_t ncomp,
  std::vector< tk::real >& phi,
  const std::vector< std::size_t >& VarList )
// *****************************************************************************
//  Kuzmin's vertex-based limiter function calculation for P1 dofs
//! \param[in] U High-order solution vector which is to be limited
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] coord Array of nodal coordinates
//! \param[in] e Id of element whose solution is to be limited
//! \param[in] rdof Maximum number of reconstructed degrees of freedom
//! \param[in] dof_el Local number of degrees of freedom
//! \param[in] ncomp Number of scalar components in this PDE system
//! \param[in,out] phi Limiter function for solution in element e
//! \param[in] VarList List of variable indices to be limited
// *****************************************************************************
{
  // Kuzmin's vertex-based TVD limiter uses min-max bounds that the
  // high-order solution should satisfy, to ensure TVD properties. For a
  // high-order method like DG, this involves the following steps:
  // 1. Find min-max bounds in the nodal-neighborhood of cell.
  // 2. Calculate the limiter function (Superbee) for all the vertices of cell.
  //    From these, use the minimum value of the limiter function.

  // Prepare for calculating Basis functions
  const auto& cx = coord[0];
  const auto& cy = coord[1];
  const auto& cz = coord[2];

  // Extract the element coordinates
  std::array< std::array< tk::real, 3>, 4 > coordel {{
    {{ cx[ inpoel[4*e  ] ], cy[ inpoel[4*e  ] ], cz[ inpoel[4*e  ] ] }},
    {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
    {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
    {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }} }};

  // Compute the determinant of Jacobian matrix
  auto detT =
    tk::Jacobian( coordel[0], coordel[1], coordel[2], coordel[3] );

  std::vector< tk::real > uMin(VarList.size(), 0.0),
                          uMax(VarList.size(), 0.0);

  // loop over all nodes of the element e
  for (std::size_t lp=0; lp<4; ++lp)
  {
    // reset min/max
    for (std::size_t i=0; i<VarList.size(); ++i)
    {
      auto mark = VarList[i]*rdof;
      uMin[i] = U(e, mark);
      uMax[i] = U(e, mark);
    }
    auto p = inpoel[4*e+lp];
    const auto& pesup = tk::cref_find(esup, p);

    // ----- Step-1: find min/max in the neighborhood of node p
    // loop over all the internal elements surrounding this node p
    for (auto er : pesup)
    {
      for (std::size_t i=0; i<VarList.size(); ++i)
      {
        auto mark = VarList[i]*rdof;
        uMin[i] = std::min(uMin[i], U(er, mark));
        uMax[i] = std::max(uMax[i], U(er, mark));
      }
    }

    // ----- Step-2: compute the limiter function at this node
    // find high-order solution
    std::vector< tk::real > state;
    std::array< tk::real, 3 > gp{cx[p], cy[p], cz[p]};
    auto B_p = tk::eval_basis( rdof,
          tk::Jacobian( coordel[0], gp, coordel[2], coordel[3] ) / detT,
          tk::Jacobian( coordel[0], coordel[1], gp, coordel[3] ) / detT,
          tk::Jacobian( coordel[0], coordel[1], coordel[2], gp ) / detT );
    state = tk::eval_state(ncomp, rdof, dof_el, e, U, B_p);

    Assert( state.size() == ncomp, "Size mismatch" );

    // compute the limiter function
    for (std::size_t i=0; i<VarList.size(); ++i)
    {
      auto c = VarList[i];
      auto phi_gp = 1.0;
      auto mark = c*rdof;
      auto uNeg = state[c] - U(e, mark);
      auto uref = std::max(std::fabs(U(e,mark)), 1e-14);
      if (uNeg > 1.0e-06*uref)
      {
        phi_gp = std::min( 1.0, (uMax[i]-U(e, mark))/uNeg );
      }
      else if (uNeg < -1.0e-06*uref)
      {
        phi_gp = std::min( 1.0, (uMin[i]-U(e, mark))/uNeg );
      }
      else
      {
        phi_gp = 1.0;
      }

    // ----- Step-3: take the minimum of the nodal-limiter functions
      phi[c] = std::min( phi[c], phi_gp );
    }
  }
}

void
VertexBasedLimiting_P2( const std::vector< std::vector< tk::real > >& unk,
  const tk::Fields& U,
  const std::map< std::size_t, std::vector< std::size_t > >& esup,
  const std::vector< std::size_t >& inpoel,
  std::size_t e,
  std::size_t rdof,
  [[maybe_unused]] std::size_t dof_el,
  std::size_t ncomp,
  const std::vector< std::size_t >& gid,
  const std::unordered_map< std::size_t, std::size_t >& bid,
  const std::vector< std::vector<tk::real> >& NodalExtrm,
  const std::vector< std::size_t >& VarList,
  std::vector< tk::real >& phi )
// *****************************************************************************
//  Kuzmin's vertex-based limiter function calculation for P2 dofs
//! \param[in] U High-order solution vector which is to be limited
//! \param[in] esup Elements surrounding points
//! \param[in] inpoel Element connectivity
//! \param[in] e Id of element whose solution is to be limited
//! \param[in] rdof Maximum number of reconstructed degrees of freedom
//! \param[in] dof_el Local number of degrees of freedom
//! \param[in] ncomp Number of scalar components in this PDE system
//! \param[in] gid Local->global node id map
//! \param[in] bid Local chare-boundary node ids (value) associated to
//!   global node ids (key)
//! \param[in] NodalExtrm Chare-boundary nodal extrema
//! \param[in] VarList List of variable indices that need to be limited
//! \param[out] phi Limiter function for solution in element e
//! \details This function limits the P2 dofs of P2 solution in a hierachical
//!   way to P1 dof limiting. Here we treat the first order derivatives the same
//!   way as cell average while second order derivatives represent the gradients
//!   to be limited in the P1 limiting procedure.
// *****************************************************************************
{
  const auto nelem = inpoel.size() / 4;

  std::vector< std::vector< tk::real > > uMin, uMax;
  uMin.resize( VarList.size(), std::vector<tk::real>(3, 0.0) );
  uMax.resize( VarList.size(), std::vector<tk::real>(3, 0.0) );

  // The coordinates of centroid in the reference domain
  std::array< std::vector< tk::real >, 3 > center;
  center[0].resize(1, 0.25);
  center[1].resize(1, 0.25);
  center[2].resize(1, 0.25);

  std::array< std::array< tk::real, 4 >, 3 > cnodes{{
    {{0, 1, 0, 0}},
    {{0, 0, 1, 0}},
    {{0, 0, 0, 1}} }};

  // loop over all nodes of the element e
  for (std::size_t lp=0; lp<4; ++lp)
  {
    // Find the max/min first-order derivatives for internal element
    for (std::size_t i=0; i<VarList.size(); ++i)
    {
      for (std::size_t idir=1; idir < 4; ++idir)
      {
        uMin[i][idir-1] = unk[VarList[i]][idir];
        uMax[i][idir-1] = unk[VarList[i]][idir];
      }
    }

    auto p = inpoel[4*e+lp];
    const auto& pesup = tk::cref_find(esup, p);

    // Step-1: find min/max first order derivative at the centroid in the
    // neighborhood of node p
    for (auto er : pesup)
    {
      if(er < nelem)      // If this is internal element
      {
        // Compute the derivatives of basis function in the reference domain
        auto dBdxi_er = tk::eval_dBdxi(rdof,
          {{center[0][0], center[1][0], center[2][0]}});

        for (std::size_t i=0; i<VarList.size(); ++i)
        {
          auto mark = VarList[i]*rdof;
          for (std::size_t idir = 0; idir < 3; ++idir)
          {
            // The first order derivative at the centroid of element er
            tk::real slope_er(0.0);
            for(std::size_t idof = 1; idof < rdof; idof++)
              slope_er += U(er, mark+idof) * dBdxi_er[idir][idof];

            uMin[i][idir] = std::min(uMin[i][idir], slope_er);
            uMax[i][idir] = std::max(uMax[i][idir], slope_er);

          }
        }
      }
    }
    // If node p is the chare-boundary node, find min/max by comparing with
    // the chare-boundary nodal extrema from vector NodalExtrm
    auto gip = bid.find( gid[p] );
    if(gip != end(bid))
    {
      auto ndof_NodalExtrm = NodalExtrm[0].size() / (ncomp * 2);
      for (std::size_t i=0; i<VarList.size(); ++i)
      {
        for (std::size_t idir = 0; idir < 3; idir++)
        {
          auto max_mark = 2*VarList[i]*ndof_NodalExtrm + 2*idir;
          auto min_mark = max_mark + 1;
          const auto& ex = NodalExtrm[gip->second];
          uMax[i][idir] = std::max(ex[max_mark], uMax[i][idir]);
          uMin[i][idir] = std::min(ex[min_mark], uMin[i][idir]);
        }
      }
    }

    //Step-2: compute the limiter function at this node
    std::array< tk::real, 3 > node{cnodes[0][lp], cnodes[1][lp], cnodes[2][lp]};

    // find high-order solution
    std::vector< std::array< tk::real, 3 > > state;
    state.resize(VarList.size());

    for (std::size_t i=0; i<VarList.size(); ++i)
    {
      auto dx = node[0] - center[0][0];
      auto dy = node[1] - center[1][0];
      auto dz = node[2] - center[2][0];

      auto c = VarList[i];

      state[i][0] = unk[c][1] + unk[c][4]*dx + unk[c][7]*dy + unk[c][8]*dz;
      state[i][1] = unk[c][2] + unk[c][5]*dy + unk[c][7]*dx + unk[c][9]*dz;
      state[i][2] = unk[c][3] + unk[c][6]*dz + unk[c][8]*dx + unk[c][9]*dy;
    }

    // compute the limiter function
    for (std::size_t i=0; i<VarList.size(); ++i)
    {
      auto c = VarList[i];
      tk::real phi_dir(1.0);
      for (std::size_t idir = 1; idir <= 3; ++idir)
      {
        phi_dir = 1.0;<--- phi_dir is assigned
        auto uNeg = state[i][idir-1] - unk[c][idir];
        auto uref = std::max(std::fabs(unk[c][idir]), 1e-14);
        if (uNeg > 1.0e-6*uref)
        {
          phi_dir =<--- phi_dir is overwritten
            std::min( 1.0, ( uMax[i][idir-1] - unk[c][idir])/uNeg );
        }
        else if (uNeg < -1.0e-6*uref)
        {
          phi_dir =
            std::min( 1.0, ( uMin[i][idir-1] - unk[c][idir])/uNeg );
        }
        else
        {
          phi_dir = 1.0;
        }

        phi[c] = std::min( phi[c], phi_dir );
      }
    }
  }
}

void consistentMultiMatLimiting_P1(
  std::size_t nmat,
  std::size_t rdof,
  std::size_t e,
  const std::vector< std::size_t >& solidx,
  tk::Fields& U,<--- Parameter 'U' can be declared with const
  [[maybe_unused]] tk::Fields& P,
  std::vector< tk::real >& phic_p1,
  std::vector< tk::real >& phic_p2 )
// *****************************************************************************
//  Consistent limiter modifications for conservative variables
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] rdof Total number of reconstructed dofs
//! \param[in] e Element being checked for consistency
//! \param[in] solidx Solid material index indicator
//! \param[in] U Vector of conservative variables
//! \param[in] P Vector of primitive variables
//! \param[in,out] phic_p1 Vector of limiter functions for P1 dofs of the
//!   conserved quantities
//! \param[in,out] phip_p2 Vector of limiter functions for P2 dofs of the
//!   conserved quantities
// *****************************************************************************
{
  // find the limiter-function for volume-fractions
  auto phi_al_p1(1.0), phi_al_p2(1.0), almax(0.0), dalmax(0.0);
  //std::size_t nmax(0);
  for (std::size_t k=0; k<nmat; ++k)
  {
    phi_al_p1 = std::min( phi_al_p1, phic_p1[volfracIdx(nmat, k)] );
    if(rdof > 4)
      phi_al_p2 = std::min( phi_al_p2, phic_p2[volfracIdx(nmat, k)] );
    if (almax < U(e,volfracDofIdx(nmat, k, rdof, 0)))
    {
      //nmax = k;
      almax = U(e,volfracDofIdx(nmat, k, rdof, 0));
    }
    tk::real dmax(0.0);
    dmax = std::max(
             std::max(
               std::abs(U(e,volfracDofIdx(nmat, k, rdof, 1))),
               std::abs(U(e,volfracDofIdx(nmat, k, rdof, 2))) ),
               std::abs(U(e,volfracDofIdx(nmat, k, rdof, 3))) );
    dalmax = std::max( dalmax, dmax );
  }

  auto al_band = 1e-4;

  //phi_al = phic[nmax];

  // determine if cell is a material-interface cell based on ad-hoc tolerances.
  // if interface-cell, then modify high-order dofs of conserved unknowns
  // consistently and use same limiter for all equations.
  // Slopes of solution variables \alpha_k \rho_k and \alpha_k \rho_k E_k need
  // to be modified in interface cells, such that slopes in the \rho_k and
  // \rho_k E_k part are ignored and only slopes in \alpha_k are considered.
  // Ideally, we would like to not do this, but this is a necessity to avoid
  // limiter-limiter interactions in multiphase CFD (see "K.-M. Shyue, F. Xiao,
  // An Eulerian interface sharpening algorithm for compressible two-phase flow:
  // the algebraic THINC approach, Journal of Computational Physics 268, 2014,
  // 326–354. doi:10.1016/j.jcp.2014.03.010." and "A. Chiapolino, R. Saurel,
  // B. Nkonga, Sharpening diffuse interfaces with compressible fluids on
  // unstructured meshes, Journal of Computational Physics 340 (2017) 389–417.
  // doi:10.1016/j.jcp.2017.03.042."). This approximation should be applied in
  // as narrow a band of interface-cells as possible. The following if-test
  // defines this band of interface-cells. This tests checks the value of the
  // maximum volume-fraction in the cell (almax) and the maximum change in
  // volume-fraction in the cell (dalmax, calculated from second-order DOFs),
  // to determine the band of interface-cells where the aforementioned fix needs
  // to be applied. This if-test says that, the fix is applied when the change
  // in volume-fraction across a cell is greater than 0.1, *and* the
  // volume-fraction is between 0.1 and 0.9.
  if ( //dalmax > al_band &&
       (almax > al_band && almax < (1.0-al_band)) )
  {
    // 1. consistent high-order dofs
    for (std::size_t k=0; k<nmat; ++k)
    {
      auto alk =
        std::max( 1.0e-14, U(e,volfracDofIdx(nmat, k, rdof, 0)) );
      auto rhok = U(e,densityDofIdx(nmat, k, rdof, 0)) / alk;
      auto rhoE = U(e,energyDofIdx(nmat, k, rdof, 0)) / alk;
      for (std::size_t idof=1; idof<rdof; ++idof)
      {
          U(e,densityDofIdx(nmat, k, rdof, idof)) = rhok *
            U(e,volfracDofIdx(nmat, k, rdof, idof));
          U(e,energyDofIdx(nmat, k, rdof, idof)) = rhoE *
            U(e,volfracDofIdx(nmat, k, rdof, idof));
      }
      if (solidx[k] > 0)
        for (std::size_t i=0; i<3; ++i)
          for (std::size_t j=0; j<3; ++j)
          {
            for (std::size_t idof=1; idof<rdof; ++idof)
              U(e,deformDofIdx(nmat,solidx[k],i,j,rdof,idof)) = 0.0;
          }
    }

    // 2. same limiter for all volume-fractions and densities
    for (std::size_t k=0; k<nmat; ++k)
    {
      phic_p1[volfracIdx(nmat, k)] = phi_al_p1;
      phic_p1[densityIdx(nmat, k)] = phi_al_p1;
      phic_p1[energyIdx(nmat, k)] = phi_al_p1;
      if (solidx[k] > 0)
        for (std::size_t i=0; i<3; ++i)
          for (std::size_t j=0; j<3; ++j)
            phic_p1[deformIdx(nmat,solidx[k],i,j)] = phi_al_p1;
    }
    if(rdof > 4)
    {
      for (std::size_t k=0; k<nmat; ++k)
      {
        phic_p2[volfracIdx(nmat, k)] = phi_al_p2;
        phic_p2[densityIdx(nmat, k)] = phi_al_p2;
        phic_p2[energyIdx(nmat, k)] = phi_al_p2;
        if (solidx[k] > 0)
          for (std::size_t i=0; i<3; ++i)
            for (std::size_t j=0; j<3; ++j)
              phic_p2[deformIdx(nmat,solidx[k],i,j)] = phi_al_p2;
      }
    }
  }
  else
  {
    // same limiter for all volume-fractions
    for (std::size_t k=0; k<nmat; ++k)
      phic_p1[volfracIdx(nmat, k)] = phi_al_p1;
    if(rdof > 4)
      for (std::size_t k=0; k<nmat; ++k)
        phic_p2[volfracIdx(nmat, k)] = phi_al_p2;
  }
}

void BoundPreservingLimiting( std::size_t nmat,
                              std::size_t ndof,
                              std::size_t e,
                              const std::vector< std::size_t >& inpoel,
                              const tk::UnsMesh::Coords& coord,
                              const tk::Fields& U,
                              std::vector< tk::real >& phic_p1,
                              std::vector< tk::real >& phic_p2 )
// *****************************************************************************
//  Bound preserving limiter for volume fractions when MulMat scheme is selected
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] ndof Total number of reconstructed dofs
//! \param[in] e Element being checked for consistency
//! \param[in] inpoel Element connectivity
//! \param[in] coord Array of nodal coordinates
//! \param[in,out] U Second-order solution vector which gets modified near
//!   material interfaces for consistency
//! \param[in] unk Vector of conservative variables based on Taylor basis
//! \param[in,out] phic_p1 Vector of limiter functions for P1 dofs of the
//!   conserved quantities
//! \param[in,out] phic_p2 Vector of limiter functions for P2 dofs of the
//!   conserved quantities
//! \details This bound-preserving limiter is specifically meant to enforce
//!   bounds [0,1], but it does not suppress oscillations like the other 'TVD'
//!   limiters. TVD limiters on the other hand, do not preserve such bounds. A
//!   combination of oscillation-suppressing and bound-preserving limiters can
//!   obtain a non-oscillatory and bounded solution.
// *****************************************************************************
{
  const auto& cx = coord[0];
  const auto& cy = coord[1];
  const auto& cz = coord[2];

  // Extract the element coordinates
  std::array< std::array< tk::real, 3>, 4 > coordel {{
    {{ cx[ inpoel[4*e  ] ], cy[ inpoel[4*e  ] ], cz[ inpoel[4*e  ] ] }},
    {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
    {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
    {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }} }};

  // Compute the determinant of Jacobian matrix
  auto detT =
    tk::Jacobian( coordel[0], coordel[1], coordel[2], coordel[3] );

  std::vector< tk::real > phi_bound(nmat, 1.0);

  // Compute the upper and lower bound for volume fraction
  const tk::real min = 1e-14;
  const tk::real max = 1.0 - min * static_cast<tk::real>(nmat - 1);

  // loop over all faces of the element e
  for (std::size_t lf=0; lf<4; ++lf)
  {
    // Extract the face coordinates
    std::array< std::size_t, 3 > inpofa_l {{ inpoel[4*e+tk::lpofa[lf][0]],
                                             inpoel[4*e+tk::lpofa[lf][1]],
                                             inpoel[4*e+tk::lpofa[lf][2]] }};

    std::array< std::array< tk::real, 3>, 3 > coordfa {{
      {{ cx[ inpofa_l[0] ], cy[ inpofa_l[0] ], cz[ inpofa_l[0] ] }},
      {{ cx[ inpofa_l[1] ], cy[ inpofa_l[1] ], cz[ inpofa_l[1] ] }},
      {{ cx[ inpofa_l[2] ], cy[ inpofa_l[2] ], cz[ inpofa_l[2] ] }} }};

    auto ng = tk::NGfa(ndof);

    // arrays for quadrature points
    std::array< std::vector< tk::real >, 2 > coordgp;
    std::vector< tk::real > wgp;

    coordgp[0].resize( ng );
    coordgp[1].resize( ng );
    wgp.resize( ng );

    // get quadrature point weights and coordinates for triangle
    tk::GaussQuadratureTri( ng, coordgp, wgp );

    // Gaussian quadrature
    for (std::size_t igp=0; igp<ng; ++igp)
    {
      // Compute the coordinates of quadrature point at physical domain
      auto gp = tk::eval_gp( igp, coordfa, coordgp );

      //Compute the basis functions
      auto B = tk::eval_basis( ndof,
            tk::Jacobian( coordel[0], gp, coordel[2], coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], gp, coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], coordel[2], gp ) / detT );

      auto state = eval_state( U.nprop()/ndof, ndof, ndof, e, U, B );

      for(std::size_t imat = 0; imat < nmat; imat++)
      {
        auto phi = BoundPreservingLimitingFunction( min, max,
          state[volfracIdx(nmat, imat)],
          U(e,volfracDofIdx(nmat, imat, ndof, 0)) );
        phi_bound[imat] = std::min( phi_bound[imat], phi );
      }
    }
  }

  // If DG(P2), the bound-preserving limiter should also be applied to the gauss
  // point within the element
  if(ndof > 4)
  {
    auto ng = tk::NGvol(ndof);

    // arrays for quadrature points
    std::array< std::vector< tk::real >, 3 > coordgp;
    std::vector< tk::real > wgp;

    coordgp[0].resize( ng );
    coordgp[1].resize( ng );
    coordgp[2].resize( ng );
    wgp.resize( ng );

    tk::GaussQuadratureTet( ng, coordgp, wgp );

    for (std::size_t igp=0; igp<ng; ++igp)
    {
      // Compute the basis function
      auto B = tk::eval_basis( ndof, coordgp[0][igp], coordgp[1][igp],
        coordgp[2][igp] );

      auto state = tk::eval_state(U.nprop()/ndof, ndof, ndof, e, U, B);

      for(std::size_t imat = 0; imat < nmat; imat++)
      {
        auto phi = BoundPreservingLimitingFunction(min, max,
          state[volfracIdx(nmat, imat)],
          U(e,volfracDofIdx(nmat, imat, ndof, 0)) );
        phi_bound[imat] = std::min( phi_bound[imat], phi );
      }
    }
  }

  for(std::size_t k=0; k<nmat; k++)
    phic_p1[volfracIdx(nmat, k)] = std::min(phi_bound[k],
      phic_p1[volfracIdx(nmat, k)]);

  if(ndof > 4)
    for(std::size_t k=0; k<nmat; k++)
      phic_p2[volfracIdx(nmat, k)] = std::min(phi_bound[k],
        phic_p2[volfracIdx(nmat, k)]);
}

tk::real
BoundPreservingLimitingFunction( const tk::real min,
                                 const tk::real max,
                                 const tk::real al_gp,
                                 const tk::real al_avg )
// *****************************************************************************
//  Bound-preserving limiter function for the volume fractions
//! \param[in] min Minimum bound for volume fraction
//! \param[in] max Maximum bound for volume fraction
//! \param[in] al_gp Volume fraction at the quadrature point
//! \param[in] al_avg Cell-average volume fraction
//! \return The limiting coefficient from the bound-preserving limiter function
// *****************************************************************************
{
  tk::real phi(1.0), al_diff(0.0);
  al_diff = al_gp - al_avg;
  if(al_gp > max && fabs(al_diff) > 1e-15)
    phi = std::fabs( (max - al_avg) / al_diff );
  else if(al_gp < min && fabs(al_diff) > 1e-15)
    phi = std::fabs( (min - al_avg) / al_diff );
  return phi;
}

void PositivityLimitingMultiMat( std::size_t nmat,
                                 const std::vector< inciter::EOS >& mat_blk,
                                 std::size_t rdof,
                                 std::size_t ndof_el,
                                 const std::vector< std::size_t >& ndofel,
                                 std::size_t e,
                                 const std::vector< std::size_t >& inpoel,
                                 const tk::UnsMesh::Coords& coord,
                                 const std::vector< int >& esuel,
                                 const tk::Fields& U,
                                 const tk::Fields& P,
                                 std::vector< tk::real >& phic_p1,
                                 std::vector< tk::real >& phic_p2,
                                 std::vector< tk::real >& phip_p1,
                                 std::vector< tk::real >& phip_p2 )
// *****************************************************************************
//  Positivity preserving limiter for multi-material solver
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] mat_blk EOS material block
//! \param[in] rdof Total number of reconstructed dofs
//! \param[in] ndof_el Number of dofs for element e
//! \param[in] ndofel Vector of local number of degrees of freedome
//! \param[in] e Element being checked for consistency
//! \param[in] inpoel Element connectivity
//! \param[in] coord Array of nodal coordinates
//! \param[in] esuel Elements surrounding elements
//! \param[in] U Vector of conservative variables
//! \param[in] P Vector of primitive variables
//! \param[in,out] phic_p1 Vector of limiter functions for P1 dofs of the
//!   conserved quantities
//! \param[in,out] phic_p2 Vector of limiter functions for P2 dofs of the
//!   conserved quantities
//! \param[in,out] phip_p1 Vector of limiter functions for P1 dofs of the
//!   primitive quantities
//! \param[in,out] phip_p2 Vector of limiter functions for P2 dofs of the
//!   primitive quantities
// *****************************************************************************
{
  const auto ncomp = U.nprop() / rdof;
  const auto nprim = P.nprop() / rdof;

  const auto& cx = coord[0];
  const auto& cy = coord[1];
  const auto& cz = coord[2];

  // Extract the element coordinates
  std::array< std::array< tk::real, 3>, 4 > coordel {{
    {{ cx[ inpoel[4*e  ] ], cy[ inpoel[4*e  ] ], cz[ inpoel[4*e  ] ] }},
    {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
    {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
    {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }} }};

  // Compute the determinant of Jacobian matrix
  auto detT =
    tk::Jacobian( coordel[0], coordel[1], coordel[2], coordel[3] );

  std::vector< tk::real > phic_bound(ncomp, 1.0);
  std::vector< tk::real > phip_bound(nprim, 1.0);

  const tk::real min = 1e-15;

  for (std::size_t lf=0; lf<4; ++lf)
  {
    std::array< std::size_t, 3 > inpofa_l {{ inpoel[4*e+tk::lpofa[lf][0]],
                                             inpoel[4*e+tk::lpofa[lf][1]],
                                             inpoel[4*e+tk::lpofa[lf][2]] }};

    std::array< std::array< tk::real, 3>, 3 > coordfa {{
      {{ cx[ inpofa_l[0] ], cy[ inpofa_l[0] ], cz[ inpofa_l[0] ] }},
      {{ cx[ inpofa_l[1] ], cy[ inpofa_l[1] ], cz[ inpofa_l[1] ] }},
      {{ cx[ inpofa_l[2] ], cy[ inpofa_l[2] ], cz[ inpofa_l[2] ] }} }};

    std::size_t nel;
    if (esuel[4*e+lf] == -1) nel = e;
    else nel = static_cast< std::size_t >(esuel[4*e+lf]);

    auto ng = tk::NGfa(std::max(ndofel[e], ndofel[nel]));

    std::array< std::vector< tk::real >, 2 > coordgp;
    std::vector< tk::real > wgp;

    coordgp[0].resize( ng );
    coordgp[1].resize( ng );
    wgp.resize( ng );

    tk::GaussQuadratureTri( ng, coordgp, wgp );

    for (std::size_t igp=0; igp<ng; ++igp)
    {
      auto gp = tk::eval_gp( igp, coordfa, coordgp );
      auto B = tk::eval_basis( ndof_el,
            tk::Jacobian( coordel[0], gp, coordel[2], coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], gp, coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], coordel[2], gp ) / detT );

      auto state = eval_state(ncomp, rdof, ndof_el, e, U, B);
      auto sprim = eval_state(nprim, rdof, ndof_el, e, P, B);

      for(std::size_t imat = 0; imat < nmat; imat++)
      {
        tk::real phi_rho(1.0), phi_rhoe(1.0), phi_pre(1.0);<--- phi_rho is assigned<--- phi_rhoe is assigned<--- phi_pre is assigned
        // Evaluate the limiting coefficient for material density
        auto rho = state[densityIdx(nmat, imat)];
        auto rho_avg = U(e, densityDofIdx(nmat, imat, rdof, 0));
        phi_rho = PositivityLimiting(min, rho, rho_avg);<--- phi_rho is overwritten
        phic_bound[densityIdx(nmat, imat)] =
          std::min(phic_bound[densityIdx(nmat, imat)], phi_rho);
        // Evaluate the limiting coefficient for material energy
        auto rhoe = state[energyIdx(nmat, imat)];
        auto rhoe_avg = U(e, energyDofIdx(nmat, imat, rdof, 0));
        phi_rhoe = PositivityLimiting(min, rhoe, rhoe_avg);<--- phi_rhoe is overwritten
        phic_bound[energyIdx(nmat, imat)] =
          std::min(phic_bound[energyIdx(nmat, imat)], phi_rhoe);
        // Evaluate the limiting coefficient for material pressure
        auto min_pre = std::max(min, state[volfracIdx(nmat, imat)] *
          mat_blk[imat].compute< EOS::min_eff_pressure >(min, rho,
          state[volfracIdx(nmat, imat)]));
        auto pre = sprim[pressureIdx(nmat, imat)];
        auto pre_avg = P(e, pressureDofIdx(nmat, imat, rdof, 0));
        phi_pre = PositivityLimiting(min_pre, pre, pre_avg);<--- phi_pre is overwritten
        phip_bound[pressureIdx(nmat, imat)] =
          std::min(phip_bound[pressureIdx(nmat, imat)], phi_pre);
      }
    }
  }

  if(ndofel[e] > 4)
  {
    auto ng = tk::NGvol(ndof_el);
    std::array< std::vector< tk::real >, 3 > coordgp;
    std::vector< tk::real > wgp;

    coordgp[0].resize( ng );
    coordgp[1].resize( ng );
    coordgp[2].resize( ng );
    wgp.resize( ng );

    tk::GaussQuadratureTet( ng, coordgp, wgp );

    for (std::size_t igp=0; igp<ng; ++igp)
    {
      auto B = tk::eval_basis( ndof_el, coordgp[0][igp], coordgp[1][igp],
        coordgp[2][igp] );

      auto state = eval_state(ncomp, rdof, ndof_el, e, U, B);
      auto sprim = eval_state(nprim, rdof, ndof_el, e, P, B);

      for(std::size_t imat = 0; imat < nmat; imat++)
      {
        tk::real phi_rho(1.0), phi_rhoe(1.0), phi_pre(1.0);<--- phi_rho is assigned<--- phi_rhoe is assigned<--- phi_pre is assigned
        // Evaluate the limiting coefficient for material density
        auto rho = state[densityIdx(nmat, imat)];
        auto rho_avg = U(e, densityDofIdx(nmat, imat, rdof, 0));
        phi_rho = PositivityLimiting(min, rho, rho_avg);<--- phi_rho is overwritten
        phic_bound[densityIdx(nmat, imat)] =
          std::min(phic_bound[densityIdx(nmat, imat)], phi_rho);
        // Evaluate the limiting coefficient for material energy
        auto rhoe = state[energyIdx(nmat, imat)];
        auto rhoe_avg = U(e, energyDofIdx(nmat, imat, rdof, 0));
        phi_rhoe = PositivityLimiting(min, rhoe, rhoe_avg);<--- phi_rhoe is overwritten
        phic_bound[energyIdx(nmat, imat)] =
          std::min(phic_bound[energyIdx(nmat, imat)], phi_rhoe);
        // Evaluate the limiting coefficient for material pressure
        auto min_pre = std::max(min, state[volfracIdx(nmat, imat)] *
          mat_blk[imat].compute< EOS::min_eff_pressure >(min, rho,
          state[volfracIdx(nmat, imat)]));
        auto pre = sprim[pressureIdx(nmat, imat)];
        auto pre_avg = P(e, pressureDofIdx(nmat, imat, rdof, 0));
        phi_pre = PositivityLimiting(min_pre, pre, pre_avg);<--- phi_pre is overwritten
        phip_bound[pressureIdx(nmat, imat)] =
          std::min(phip_bound[pressureIdx(nmat, imat)], phi_pre);
      }
    }
  }

  // apply new bounds to material quantities
  for (std::size_t k=0; k<nmat; ++k) {
    // mat density
    phic_p1[densityIdx(nmat, k)] = std::min( phic_bound[densityIdx(nmat, k)],
      phic_p1[densityIdx(nmat, k)] );
    // mat energy
    phic_p1[energyIdx(nmat, k)] = std::min( phic_bound[energyIdx(nmat, k)],
      phic_p1[energyIdx(nmat, k)] );
    // mat pressure
    phip_p1[pressureIdx(nmat, k)] = std::min( phip_bound[pressureIdx(nmat, k)],
      phip_p1[pressureIdx(nmat, k)] );

    // for dgp2
    if (ndof_el > 4) {
      // mat density
      phic_p2[densityIdx(nmat, k)] = std::min( phic_bound[densityIdx(nmat, k)],
        phic_p2[densityIdx(nmat, k)] );
      // mat energy
      phic_p2[energyIdx(nmat, k)] = std::min( phic_bound[energyIdx(nmat, k)],
        phic_p2[energyIdx(nmat, k)] );
      // mat pressure
      phip_p2[pressureIdx(nmat, k)] = std::min( phip_bound[pressureIdx(nmat, k)],
        phip_p2[pressureIdx(nmat, k)] );
    }
  }
}

void PositivityPreservingMultiMat_FV(
  const std::vector< std::size_t >& inpoel,
  std::size_t nelem,
  std::size_t nmat,
  const std::vector< inciter::EOS >& mat_blk,
  const tk::UnsMesh::Coords& coord,
  const tk::Fields& /*geoFace*/,
  tk::Fields& U,
  tk::Fields& P )
// *****************************************************************************
//  Positivity preserving limiter for the FV multi-material solver
//! \param[in] inpoel Element connectivity
//! \param[in] nelem Number of elements
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] mat_blk Material EOS block
//! \param[in] coord Array of nodal coordinates
////! \param[in] geoFace Face geometry array
//! \param[in,out] U High-order solution vector which gets limited
//! \param[in,out] P High-order vector of primitives which gets limited
//! \details This positivity preserving limiter function should be called for
//!   FV multimat.
// *****************************************************************************
{
  const auto rdof = inciter::g_inputdeck.get< tag::rdof >();
  const auto ncomp = U.nprop() / rdof;
  const auto nprim = P.nprop() / rdof;

  const auto& cx = coord[0];
  const auto& cy = coord[1];
  const auto& cz = coord[2];

  for (std::size_t e=0; e<nelem; ++e)
  {
    // Extract the element coordinates
    std::array< std::array< tk::real, 3>, 4 > coordel {{
      {{ cx[ inpoel[4*e  ] ], cy[ inpoel[4*e  ] ], cz[ inpoel[4*e  ] ] }},
      {{ cx[ inpoel[4*e+1] ], cy[ inpoel[4*e+1] ], cz[ inpoel[4*e+1] ] }},
      {{ cx[ inpoel[4*e+2] ], cy[ inpoel[4*e+2] ], cz[ inpoel[4*e+2] ] }},
      {{ cx[ inpoel[4*e+3] ], cy[ inpoel[4*e+3] ], cz[ inpoel[4*e+3] ] }} }};

    // Compute the determinant of Jacobian matrix
    auto detT =
      tk::Jacobian( coordel[0], coordel[1], coordel[2], coordel[3] );

    std::vector< tk::real > phic(ncomp, 1.0);
    std::vector< tk::real > phip(nprim, 1.0);

    const tk::real min = 1e-15;

    // 1. Enforce positive density (total energy will be positive if pressure
    //    and density are positive)
    for (std::size_t lf=0; lf<4; ++lf)
    {
      std::array< std::size_t, 3 > inpofa_l {{ inpoel[4*e+tk::lpofa[lf][0]],
                                               inpoel[4*e+tk::lpofa[lf][1]],
                                               inpoel[4*e+tk::lpofa[lf][2]] }};

      // face coordinates
      std::array< std::array< tk::real, 3>, 3 > coordfa {{
        {{ cx[ inpofa_l[0] ], cy[ inpofa_l[0] ], cz[ inpofa_l[0] ] }},
        {{ cx[ inpofa_l[1] ], cy[ inpofa_l[1] ], cz[ inpofa_l[1] ] }},
        {{ cx[ inpofa_l[2] ], cy[ inpofa_l[2] ], cz[ inpofa_l[2] ] }} }};

      // face centroid
      std::array< tk::real, 3 > fc{{
        (coordfa[0][0]+coordfa[1][0]+coordfa[2][0])/3.0 ,
        (coordfa[0][1]+coordfa[1][1]+coordfa[2][1])/3.0 ,
        (coordfa[0][2]+coordfa[1][2]+coordfa[2][2])/3.0 }};

      auto B = tk::eval_basis( rdof,
            tk::Jacobian( coordel[0], fc, coordel[2], coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], fc, coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], coordel[2], fc ) / detT );
      auto state = eval_state(ncomp, rdof, rdof, e, U, B);

      for(std::size_t i=0; i<nmat; i++)
      {
        // Evaluate the limiting coefficient for material density
        auto rho = state[densityIdx(nmat, i)];
        auto rho_avg = U(e, densityDofIdx(nmat, i, rdof, 0));
        auto phi_rho = PositivityLimiting(min, rho, rho_avg);
        phic[densityIdx(nmat, i)] =
          std::min(phic[densityIdx(nmat, i)], phi_rho);
      }
    }
    // apply limiter coefficient
    for(std::size_t i=0; i<nmat; i++)
    {
      U(e, densityDofIdx(nmat,i,rdof,1)) *= phic[densityIdx(nmat,i)];
      U(e, densityDofIdx(nmat,i,rdof,2)) *= phic[densityIdx(nmat,i)];
      U(e, densityDofIdx(nmat,i,rdof,3)) *= phic[densityIdx(nmat,i)];
    }

    // 2. Enforce positive pressure (assuming density is positive)
    for (std::size_t lf=0; lf<4; ++lf)
    {
      std::array< std::size_t, 3 > inpofa_l {{ inpoel[4*e+tk::lpofa[lf][0]],
                                               inpoel[4*e+tk::lpofa[lf][1]],
                                               inpoel[4*e+tk::lpofa[lf][2]] }};

      // face coordinates
      std::array< std::array< tk::real, 3>, 3 > coordfa {{
        {{ cx[ inpofa_l[0] ], cy[ inpofa_l[0] ], cz[ inpofa_l[0] ] }},
        {{ cx[ inpofa_l[1] ], cy[ inpofa_l[1] ], cz[ inpofa_l[1] ] }},
        {{ cx[ inpofa_l[2] ], cy[ inpofa_l[2] ], cz[ inpofa_l[2] ] }} }};

      // face centroid
      std::array< tk::real, 3 > fc{{
        (coordfa[0][0]+coordfa[1][0]+coordfa[2][0])/3.0 ,
        (coordfa[0][1]+coordfa[1][1]+coordfa[2][1])/3.0 ,
        (coordfa[0][2]+coordfa[1][2]+coordfa[2][2])/3.0 }};

      auto B = tk::eval_basis( rdof,
            tk::Jacobian( coordel[0], fc, coordel[2], coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], fc, coordel[3] ) / detT,
            tk::Jacobian( coordel[0], coordel[1], coordel[2], fc ) / detT );
      auto state = eval_state(ncomp, rdof, rdof, e, U, B);
      auto sprim = eval_state(nprim, rdof, rdof, e, P, B);

      for(std::size_t i=0; i<nmat; i++)
      {
        tk::real phi_pre(1.0);<--- phi_pre is assigned
        // Evaluate the limiting coefficient for material pressure
        auto rho = state[densityIdx(nmat, i)];
        auto min_pre = std::max(min, U(e,volfracDofIdx(nmat,i,rdof,0)) *
          mat_blk[i].compute< EOS::min_eff_pressure >(min, rho,
          U(e,volfracDofIdx(nmat,i,rdof,0))));
        auto pre = sprim[pressureIdx(nmat, i)];
        auto pre_avg = P(e, pressureDofIdx(nmat, i, rdof, 0));
        phi_pre = PositivityLimiting(min_pre, pre, pre_avg);<--- phi_pre is overwritten
        phip[pressureIdx(nmat, i)] =
          std::min(phip[pressureIdx(nmat, i)], phi_pre);
      }
    }
    // apply limiter coefficient
    for(std::size_t i=0; i<nmat; i++)
    {
      P(e, pressureDofIdx(nmat,i,rdof,1)) *= phip[pressureIdx(nmat,i)];
      P(e, pressureDofIdx(nmat,i,rdof,2)) *= phip[pressureIdx(nmat,i)];
      P(e, pressureDofIdx(nmat,i,rdof,3)) *= phip[pressureIdx(nmat,i)];
    }
  }
}

tk::real
PositivityLimiting( const tk::real min,
                    const tk::real u_gp,
                    const tk::real u_avg )
// *****************************************************************************
//  Positivity-preserving limiter function
//! \param[in] min Minimum bound for volume fraction
//! \param[in] u_gp Variable quantity at the quadrature point
//! \param[in] u_avg Cell-average variable quantitiy
//! \return The limiting coefficient from the positivity-preserving limiter
//!   function
// *****************************************************************************
{
  tk::real phi(1.0);
  tk::real diff = u_gp - u_avg;
  // Only when u_gp is less than minimum threshold and the high order
  // contribution is not zero, the limiting function will be applied
  if(u_gp < min)
    phi = std::fabs( (min - u_avg) / (diff+std::copysign(1e-15,diff)) );
  return phi;
}

bool
interfaceIndicator( std::size_t nmat,
  const std::vector< tk::real >& al,
  std::vector< std::size_t >& matInt )
// *****************************************************************************
//  Interface indicator function, which checks element for material interface
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] al Cell-averaged volume fractions
//! \param[in] matInt Array indicating which material has an interface
//! \return Boolean which indicates if the element contains a material interface
// *****************************************************************************
{
  bool intInd = false;

  // limits under which compression is to be performed
  auto al_eps = 1e-08;
  auto loLim = 2.0 * al_eps;
  auto hiLim = 1.0 - loLim;

  auto almax = 0.0;
  for (std::size_t k=0; k<nmat; ++k)
  {
    almax = std::max(almax, al[k]);
    matInt[k] = 0;
    if ((al[k] > loLim) && (al[k] < hiLim)) matInt[k] = 1;
  }

  if ((almax > loLim) && (almax < hiLim)) intInd = true;

  return intInd;
}

void MarkShockCells ( const bool pref,
                      const std::size_t nelem,
                      const std::size_t nmat,
                      const std::size_t ndof,
                      const std::size_t rdof,
                      const std::vector< inciter::EOS >& mat_blk,
                      const std::vector< std::size_t >& ndofel,
                      const std::vector< std::size_t >& inpoel,
                      const tk::UnsMesh::Coords& coord,
                      const inciter::FaceData& fd,
                      [[maybe_unused]] const tk::Fields& geoFace,
                      const tk::Fields& geoElem,
                      const tk::FluxFn& flux,
                      const std::vector< std::size_t >& solidx,
                      const tk::Fields& U,
                      const tk::Fields& P,
                      std::vector< std::size_t >& shockmarker )
// *****************************************************************************
//  Mark the cells that contain discontinuity according to the interface
//    condition
//! \param[in] pref Indicator for p-adaptive algorithm
//! \param[in] nelem Number of elements
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] ndof Maximum number of degrees of freedom
//! \param[in] rdof Maximum number of reconstructed degrees of freedom
//! \param[in] mat_blk EOS material block
//! \param[in] ndofel Vector of local number of degrees of freedome
//! \param[in] inpoel Element-node connectivity
//! \param[in] coord Array of nodal coordinates
//! \param[in] fd Face connectivity and boundary conditions object
//! \param[in] geoFace Face geometry array
//! \param[in] geoElem Element geometry array
//! \param[in] flux Flux function to use
//! \param[in] solidx Solid material index indicator
//! \param[in] U Solution vector at recent time step
//! \param[in] P Vector of primitives at recent time step
//! \param[in, out] shockmarker Vector of the shock indicator
//! \details This function computes the discontinuity indicator based on
//!   interface conditon. It is based on the following paper:
//!   Hong L., Gianni A., Robert N. (2021) A moving discontinuous Galerkin
//!   finite element method with interface condition enforcement for
//!   compressible flows. Journal of Computational Physics,
//!   doi: https://doi.org/10.1016/j.jcp.2021.110618
// *****************************************************************************
{
  const auto coeff = g_inputdeck.get< tag::shock_detector_coeff >();

  std::vector< tk::real > IC(U.nunk(), 0.0);
  const auto& esuf = fd.Esuf();
  const auto& inpofa = fd.Inpofa();

  const auto& cx = coord[0];
  const auto& cy = coord[1];
  const auto& cz = coord[2];

  auto ncomp = U.nprop()/rdof;<--- Variable 'ncomp' is assigned a value that is never used.
  auto nprim = P.nprop()/rdof;<--- Variable 'nprim' is assigned a value that is never used.

  // The interface-conservation based indicator will only evaluate the flux jump
  // for the momentum equations
  std::set< std::size_t > vars;
  if(nmat > 1) {          // multi-material flow
    for (std::size_t i=0; i<3; ++i) vars.insert(momentumIdx(nmat, i));
  } else {                // single-material flow
    for (std::size_t i=1; i<=3; ++i) vars.insert(i);
  }

  // Loop over faces
  for (auto f=fd.Nbfac(); f<esuf.size()/2; ++f) {
    Assert( esuf[2*f] > -1 && esuf[2*f+1] > -1, "Interior element detected "
            "as -1" );

    std::size_t el = static_cast< std::size_t >(esuf[2*f]);
    std::size_t er = static_cast< std::size_t >(esuf[2*f+1]);

    // When the number of gauss points for the left and right element are
    // different, choose the larger ng
    auto ng_l = tk::NGfa(ndofel[el]);
    auto ng_r = tk::NGfa(ndofel[er]);

    auto ng = std::max( ng_l, ng_r );

    std::array< std::vector< tk::real >, 2 > coordgp
      { std::vector<tk::real>(ng), std::vector<tk::real>(ng) };
    std::vector< tk::real > wgp( ng );

    tk::GaussQuadratureTri( ng, coordgp, wgp );

    // Extract the element coordinates
    std::array< std::array< tk::real, 3>, 4 > coordel_l {{
      {{ cx[ inpoel[4*el  ] ], cy[ inpoel[4*el  ] ], cz[ inpoel[4*el  ] ] }},
      {{ cx[ inpoel[4*el+1] ], cy[ inpoel[4*el+1] ], cz[ inpoel[4*el+1] ] }},
      {{ cx[ inpoel[4*el+2] ], cy[ inpoel[4*el+2] ], cz[ inpoel[4*el+2] ] }},
      {{ cx[ inpoel[4*el+3] ], cy[ inpoel[4*el+3] ], cz[ inpoel[4*el+3] ] }} }};

    std::array< std::array< tk::real, 3>, 4 > coordel_r {{
      {{ cx[ inpoel[4*er  ] ], cy[ inpoel[4*er  ] ], cz[ inpoel[4*er  ] ] }},
      {{ cx[ inpoel[4*er+1] ], cy[ inpoel[4*er+1] ], cz[ inpoel[4*er+1] ] }},
      {{ cx[ inpoel[4*er+2] ], cy[ inpoel[4*er+2] ], cz[ inpoel[4*er+2] ] }},
      {{ cx[ inpoel[4*er+3] ], cy[ inpoel[4*er+3] ], cz[ inpoel[4*er+3] ] }} }};

    // Compute the determinant of Jacobian matrix
    auto detT_l =
      tk::Jacobian( coordel_l[0], coordel_l[1], coordel_l[2], coordel_l[3] );
    auto detT_r =
      tk::Jacobian( coordel_r[0], coordel_r[1], coordel_r[2], coordel_r[3] );

    std::array< std::array< tk::real, 3>, 3 > coordfa {{
      {{ cx[ inpofa[3*f  ] ], cy[ inpofa[3*f  ] ], cz[ inpofa[3*f  ] ] }},
      {{ cx[ inpofa[3*f+1] ], cy[ inpofa[3*f+1] ], cz[ inpofa[3*f+1] ] }},
      {{ cx[ inpofa[3*f+2] ], cy[ inpofa[3*f+2] ], cz[ inpofa[3*f+2] ] }} }};

    std::array< tk::real, 3 >
      fn{{ geoFace(f,1), geoFace(f,2), geoFace(f,3) }};

    // Numerator and denominator of the shock indicator
    tk::real numer(0.0), denom(0.0);
    std::vector< tk::real > fl_jump, fl_avg;
    fl_jump.resize(3, 0.0);
    fl_avg.resize(3, 0.0);

    for (std::size_t igp=0; igp<ng; ++igp) {
      auto gp = tk::eval_gp( igp, coordfa, coordgp );
      std::size_t dof_el, dof_er;
      if (rdof > ndof)
      {
        dof_el = rdof;
        dof_er = rdof;
      }
      else
      {
        dof_el = ndofel[el];
        dof_er = ndofel[er];
      }
      std::array< tk::real, 3> ref_gp_l{
        tk::Jacobian( coordel_l[0], gp, coordel_l[2], coordel_l[3] ) / detT_l,
        tk::Jacobian( coordel_l[0], coordel_l[1], gp, coordel_l[3] ) / detT_l,
        tk::Jacobian( coordel_l[0], coordel_l[1], coordel_l[2], gp ) / detT_l };
      std::array< tk::real, 3> ref_gp_r{
        tk::Jacobian( coordel_r[0], gp, coordel_r[2], coordel_r[3] ) / detT_r,
        tk::Jacobian( coordel_r[0], coordel_r[1], gp, coordel_r[3] ) / detT_r,
        tk::Jacobian( coordel_r[0], coordel_r[1], coordel_r[2], gp ) / detT_r };
      auto B_l = tk::eval_basis( dof_el, ref_gp_l[0], ref_gp_l[1], ref_gp_l[2] );
      auto B_r = tk::eval_basis( dof_er, ref_gp_r[0], ref_gp_r[1], ref_gp_r[2] );

      std::array< std::vector< tk::real >, 2 > state;

      // Evaluate the high order solution at the qudrature point
      state[0] = tk::evalPolynomialSol(mat_blk, 0, ncomp, nprim, rdof,
        nmat, el, dof_el, inpoel, coord, geoElem, ref_gp_l, B_l, U, P);
      state[1] = tk::evalPolynomialSol(mat_blk, 0, ncomp, nprim, rdof,
        nmat, er, dof_er, inpoel, coord, geoElem, ref_gp_r, B_r, U, P);

      Assert( state[0].size() == ncomp+nprim, "Incorrect size for "
              "appended boundary state vector" );
      Assert( state[1].size() == ncomp+nprim, "Incorrect size for "
              "appended boundary state vector" );

      // Force deformation unknown to first order
      for (std::size_t k=0; k<nmat; ++k)
        if (solidx[k] > 0)
          for (std::size_t i=0; i<3; ++i)
            for (std::size_t j=0; j<3; ++j)
            {
              state[0][deformIdx(nmat, solidx[k], i, j)] = U(el,deformDofIdx(
                nmat, solidx[k], i, j, rdof, 0));
              state[1][deformIdx(nmat, solidx[k], i, j)] = U(er,deformDofIdx(
                nmat, solidx[k], i, j, rdof, 0));
            }

      // Evaluate the flux
      auto fl = flux( ncomp, mat_blk, state[0], {} );
      auto fr = flux( ncomp, mat_blk, state[1], {} );

      std::size_t i(0);
      for (const auto& c : vars) {
        tk::real fn_l(0.0), fn_r(0.0);
        for(std::size_t idir = 0; idir < 3; idir++) {
          fn_l += fl[c][idir] * fn[idir];
          fn_r += fr[c][idir] * fn[idir];
        }
        fl_jump[i] += wgp[igp] * (fn_l - fn_r) * (fn_l - fn_r);
        fl_avg[i]  += wgp[igp] * (fn_l + fn_r) * (fn_l + fn_r) * 0.25;
        ++i;
      }
    }

    // Evaluate the numerator and denominator
    for(std::size_t idir = 0; idir < 3; idir++) {
      numer += std::sqrt(fl_jump[idir]);
      denom += std::sqrt(fl_avg[idir]);
    }

    tk::real Ind(0.0);
    if(denom > 1e-8)
      Ind = numer / denom;
    IC[el] = std::max(IC[el], Ind);
    IC[er] = std::max(IC[er], Ind);
  }

  // Loop over element to mark shock cell
  for (std::size_t e=0; e<nelem; ++e) {
    std::size_t dof_el = pref ? ndofel[e] : rdof;

    tk::real power = 0.0;
    if(dof_el == 10)  power = 1.5;
    else              power = 1.0;

    // Evaluate the threshold
    auto thres = coeff * std::pow(geoElem(e, 4), power);
    if(IC[e] > thres)
      shockmarker[e] = 1;
    else
      shockmarker[e] = 0;
  }
}

void
correctLimConservMultiMat(
  std::size_t nelem,
  const std::vector< EOS >& mat_blk,
  std::size_t nmat,
  const std::vector< std::size_t >& inpoel,
  const tk::UnsMesh::Coords& coord,
  const tk::Fields& geoElem,
  const tk::Fields& prim,
  tk::Fields& unk )
// *****************************************************************************
//  Update the conservative quantities after limiting for multi-material systems
//! \param[in] nelem Number of internal elements
//! \param[in] mat_blk EOS material block
//! \param[in] nmat Number of materials in this PDE system
//! \param[in] inpoel Element-node connectivity
//! \param[in] coord Array of nodal coordinates
//! \param[in] geoElem Element geometry array
//! \param[in] prim Array of primitive variables
//! \param[in,out] unk Array of conservative variables
//! \details This function computes the updated dofs for conservative
//!   quantities based on the limited primitive quantities, to re-instate
//!   consistency between the limited primitive and evolved quantities. For
//!   further details, see Pandare et al. (2023). On the Design of Stable,
//!   Consistent, and Conservative High-Order Methods for Multi-Material
//!   Hydrodynamics. J Comp Phys, 112313.
// *****************************************************************************
{
  const auto rdof = g_inputdeck.get< tag::rdof >();
  std::size_t ncomp = unk.nprop()/rdof;
  std::size_t nprim = prim.nprop()/rdof;
  const auto intsharp = inciter::g_inputdeck.get< tag::multimat,
    tag::intsharp >();

  for (std::size_t e=0; e<nelem; ++e) {
    // Here we pre-compute the right-hand-side vector. The reason that the
    // lhs in DG.cpp is not used is that the size of this vector in this
    // projection procedure should be rdof instead of ndof.
    auto L = tk::massMatrixDubiner(rdof, geoElem(e,0));

    // The right-hand side vector is sized as nprim, i.e. the primitive quantity
    // vector. However, it stores the consistently obtained values of evolved
    // quantities, since nprim is the number of evolved quantities that need to
    // be evaluated consistently. For this reason, accessing R will require
    // the primitive quantity accessors. But this access is intended to give
    // the corresponding evolved quantites, as follows:
    // pressureIdx() - mat. total energy
    // velocityIdx() - bulk momentum components
    // stressIdx() - mat. inverse deformation gradient tensor components
    std::vector< tk::real > R(nprim*rdof, 0.0);

    auto ng = tk::NGvol(rdof);

    // Arrays for quadrature points
    std::array< std::vector< tk::real >, 3 > coordgp;
    std::vector< tk::real > wgp;

    coordgp[0].resize( ng );
    coordgp[1].resize( ng );
    coordgp[2].resize( ng );
    wgp.resize( ng );

    tk::GaussQuadratureTet( ng, coordgp, wgp );

    // Loop over quadrature points in element e
    for (std::size_t igp=0; igp<ng; ++igp) {
      // Compute the basis function
      auto B = tk::eval_basis( rdof, coordgp[0][igp], coordgp[1][igp],
                               coordgp[2][igp] );

      auto w = wgp[igp] * geoElem(e, 0);

      // Evaluate the solution at quadrature point
      auto state = evalPolynomialSol(mat_blk, intsharp, ncomp, nprim,
        rdof, nmat, e, rdof, inpoel, coord, geoElem,
        {{coordgp[0][igp], coordgp[1][igp], coordgp[2][igp]}}, B, unk, prim);

      // Solution vector that stores the material energy and bulk momentum
      std::vector< tk::real > s(nprim, 0.0);

      // Bulk density at quadrature point
      tk::real rhob(0.0);
      for (std::size_t k=0; k<nmat; ++k)
        rhob += state[densityIdx(nmat, k)];

      // Velocity vector at quadrature point
      std::array< tk::real, 3 >
        vel{ state[ncomp+velocityIdx(nmat, 0)],
             state[ncomp+velocityIdx(nmat, 1)],
             state[ncomp+velocityIdx(nmat, 2)] };

      // Compute and store the bulk momentum
      for(std::size_t idir = 0; idir < 3; idir++)
        s[velocityIdx(nmat, idir)] = rhob * vel[idir];

      // Compute and store material energy at quadrature point
      for(std::size_t imat = 0; imat < nmat; imat++) {
        auto alphamat = state[volfracIdx(nmat, imat)];
        auto rhomat = state[densityIdx(nmat, imat)]/alphamat;
        auto premat = state[ncomp+pressureIdx(nmat, imat)]/alphamat;
        auto gmat = getDeformGrad(nmat, imat, state);
        s[pressureIdx(nmat,imat)] = alphamat *
          mat_blk[imat].compute< EOS::totalenergy >( rhomat, vel[0], vel[1],
          vel[2], premat, gmat );
      }

      // Evaluate the righ-hand-side vector
      for(std::size_t k = 0; k < nprim; k++) {
        auto mark = k * rdof;
        for(std::size_t idof = 0; idof < rdof; idof++)
          R[mark+idof] += w * s[k] * B[idof];
      }
    }

    // Update the high order dofs of the material energy
    for(std::size_t imat = 0; imat < nmat; imat++) {
      for(std::size_t idof = 1; idof < rdof; idof++)
        unk(e, energyDofIdx(nmat, imat, rdof, idof)) =
          R[pressureDofIdx(nmat,imat,rdof,idof)] / L[idof];
    }

    // Update the high order dofs of the bulk momentum
    for(std::size_t idir = 0; idir < 3; idir++) {
      for(std::size_t idof = 1; idof < rdof; idof++)
        unk(e, momentumDofIdx(nmat, idir, rdof, idof)) =
          R[velocityDofIdx(nmat,idir,rdof,idof)] / L[idof];
    }
  }
}

tk::real
constrain_pressure( const std::vector< EOS >& mat_blk,
  tk::real apr,
  tk::real arho,
  tk::real alpha=1.0,
  std::size_t imat=0 )
// *****************************************************************************
//  Constrain material partial pressure (alpha_k * p_k)
//! \param[in] apr Material partial pressure (alpha_k * p_k)
//! \param[in] arho Material partial density (alpha_k * rho_k)
//! \param[in] alpha Material volume fraction. Default is 1.0, so that for the
//!   single-material system, this argument can be left unspecified by the
//!   calling code
//! \param[in] imat Material-id who's EoS is required. Default is 0, so that
//!   for the single-material system, this argument can be left unspecified by
//!   the calling code
//! \return Constrained material partial pressure (alpha_k * p_k)
// *****************************************************************************
{
  return std::max(apr, alpha*mat_blk[imat].compute<
    EOS::min_eff_pressure >(1e-12, arho, alpha));
}


} // inciter::