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
// *****************************************************************************
/*!
  \file      src/Control/Inciter/InputDeck/Grammar.hpp
  \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     Inciter's input deck grammar definition
  \details   Inciter's input deck grammar definition. We use the Parsing
  Expression Grammar Template Library (PEGTL) to create the grammar and the
  associated parser. Word of advice: read from the bottom up.
*/
// *****************************************************************************
#ifndef InciterInputDeckGrammar_h
#define InciterInputDeckGrammar_h

#include <limits>
#include <cmath>

#include "CommonGrammar.hpp"
#include "CartesianProduct.hpp"
#include "Keywords.hpp"
#include "ContainerUtil.hpp"
#include "Centering.hpp"
#include "PDE/MultiMat/MultiMatIndexing.hpp"

namespace inciter {

extern ctr::InputDeck g_inputdeck_defaults;

//! Inciter input deck facilitating user input for computing shock hydrodynamics
namespace deck {

  //! \brief Specialization of tk::grm::use for Inciter's input deck parser
  template< typename keyword >
  using use = tk::grm::use< keyword, ctr::InputDeck::keywords >;

  // Inciter's InputDeck state

  //! \brief Number of registered equations
  //! \details Counts the number of parsed equation blocks during parsing.
  static tk::TaggedTuple< brigand::list<
             tag::transport,   std::size_t
           , tag::compflow,    std::size_t
           , tag::multimat,    std::size_t
         > > neq;

  //! \brief Parser-lifetime storage for point names
  //! \details Used to track the point names registered so that parsing new ones
  //!    can be required to be unique.
  static std::set< std::string > pointnames;

  //! Parser-lifetime storage of elem or node centering
  static tk::Centering centering = tk::Centering::NODE;

  //! Accepted multimat output variable labels and associated index functions
  //! \details The keys are a list of characters accepted as labels for
  //! denoting (matvar-style) output variables used for multi-material variable
  //! output. We use a case- insesitive comparitor, since when this set is used
  //! we only care about whether the variable is selected or not and not whether
  //! it denotes a full variable (upper case) or a fluctuation (lower case).
  //! This is true when matching these labels for output variables with
  //! instantaenous variables as well terms of products in parsing requested
  //! statistics (for turbulence). The values are associated indexing functions
  //! used to index into the state of the multimaterial system, all must follow
  //! the same signature.
  static std::map< char, tk::MultiMatIdxFn,
                   tk::ctr::CaseInsensitiveCharLess >
    multimatvars{
      { 'd', densityIdx }       // density
    , { 'f', volfracIdx }       // volume fraction
    , { 'm', momentumIdx }      // momentum
    , { 'e', energyIdx }        // specific total energy
    , { 'u', velocityIdx }      // velocity (primitive)
    , { 'p', pressureIdx }      // material pressure (primitive)
  };

} // ::deck
} // ::inciter

namespace tk {
namespace grm {

  using namespace tao;

  // Note that PEGTL action specializations must be in the same namespace as the
  // template being specialized. See http://stackoverflow.com/a/3052604.

  // Inciter's InputDeck actions

  //! Rule used to trigger action
  template< class eq > struct register_inciter_eq : pegtl::success {};
  //! \brief Register differential equation after parsing its block
  //! \details This is used by the error checking functors (check_*) during
  //!    parsing to identify the recently-parsed block.
  template< class eq >
  struct action< register_inciter_eq< eq > > {
    template< typename Input, typename Stack >
    static void apply( const Input&, Stack& ) {
      using inciter::deck::neq;
      ++neq.get< eq >();
    }
  };

  //! Rule used to trigger action
  template< class eq > struct check_mesh : pegtl::success {};
  //! \brief Check mesh ... end block for correctness
  template< class eq >
  struct action< check_mesh< eq > > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::neq;
      auto& mesh = stack.template get< tag::param, eq, tag::mesh >();
      auto& mesh_ref = mesh.template get< tag::reference >();
      // if no mesh reference given by user
      if (mesh_ref.empty() || mesh_ref.size() != neq.get< eq >()) {
        // put in '-', meaning no reference
        mesh_ref.push_back('-');
        auto& location = mesh.template get< tag::location >();
        // if no location, put in the origin
        if (location.size() != neq.get< eq >())
          location.push_back( { 0.0, 0.0, 0.0 } );
        else    // reference was not given, but location was, error out
          Message< Stack, ERROR, MsgKey::LOC_NOMESHREF >( stack, in );
        auto& orientation = mesh.template get< tag::orientation >();
        if (orientation.size() != neq.get< eq >())
          orientation.push_back( { 0.0, 0.0, 0.0 } );
        else    // reference was not given, but orientation was, error out
          Message< Stack, ERROR, MsgKey::ORI_NOMESHREF >( stack, in );
      }

      // Ensure the number of depvars and the number of mesh references equal
      const auto& depvar = stack.template get< tag::param, eq, tag::depvar >();
      Assert( depvar.size() == mesh_ref.size(), "Mesh ref size mismatch" );

      // Ensure mesh ref var is not the same as the current depvar and is
      // defined upstream in input file (by another solver)
      if ( mesh_ref.back() != '-' &&
           (mesh_ref.back() == depvar.back() ||
            depvars.find(mesh_ref.back()) == end(depvars)) )
      {
        Message< Stack, ERROR, MsgKey::DEPVAR_AS_MESHREF >( stack, in );
      }
    }
  };

  //! Rule used to trigger action
  template< class eq > struct check_transport : pegtl::success {};
  //! \brief Set defaults and do error checking on the transport equation block
  //! \details This is error checking that only the transport equation block
  //!   must satisfy. Besides error checking we also set defaults here as
  //!   this block is called when parsing of a transport...end block has
  //!   just finished.
  template< class eq >
  struct action< check_transport< eq > > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::neq;
      using tag::param;

      // Error out if no dependent variable has been selected
      auto& depvar = stack.template get< param, eq, tag::depvar >();
      if (depvar.empty() || depvar.size() != neq.get< eq >())
        Message< Stack, ERROR, MsgKey::NODEPVAR >( stack, in );

      // If no number of components has been selected, default to 1
      auto& ncomp = stack.template get< tag::component, eq >();
      if (ncomp.empty() || ncomp.size() != neq.get< eq >())
        ncomp.push_back( 1 );

      // If physics type is not given, default to 'advection'
      auto& physics = stack.template get< param, eq, tag::physics >();
      if (physics.empty() || physics.size() != neq.get< eq >())
        physics.push_back( inciter::ctr::PhysicsType::ADVECTION );

      // If physics type is advection-diffusion, check for correct number of
      // advection velocity, shear, and diffusion coefficients
      if (physics.back() == inciter::ctr::PhysicsType::ADVDIFF) {
        auto& u0 = stack.template get< param, eq, tag::u0 >();
        if (u0.back().size() != ncomp.back())  // must define 1 component
          Message< Stack, ERROR, MsgKey::WRONGSIZE >( stack, in );
        auto& diff = stack.template get< param, eq, tag::diffusivity >();
        if (diff.back().size() != 3*ncomp.back())  // must define 3 components
          Message< Stack, ERROR, MsgKey::WRONGSIZE >( stack, in );
        auto& lambda = stack.template get< param, eq, tag::lambda >();
        if (lambda.back().size() != 2*ncomp.back()) // must define 2 shear comps
          Message< Stack, ERROR, MsgKey::WRONGSIZE >( stack, in );
      }
      // If problem type is not given, error out
      auto& problem = stack.template get< param, eq, tag::problem >();
      if (problem.empty() || problem.size() != neq.get< eq >())
        Message< Stack, ERROR, MsgKey::NOPROBLEM >( stack, in );
      // Error check Dirichlet boundary condition block for all transport eq
      // configurations
      const auto& bc = stack.template get< param, eq, tag::bc, tag::bcdir >();
      for (const auto& s : bc)
        if (s.empty()) Message< Stack, ERROR, MsgKey::BC_EMPTY >( stack, in );

      // If interface compression is not specified, default to 'false'
      auto& intsharp = stack.template get< param, eq, tag::intsharp >();
      if (intsharp.empty() || intsharp.size() != neq.get< eq >())
        intsharp.push_back( 0 );

      // If interface compression parameter is not specified, default to 1.0
      auto& intsharp_p = stack.template get< param, eq,
                                            tag::intsharp_param >();
      if (intsharp_p.empty() || intsharp_p.size() != neq.get< eq >())
        intsharp_p.push_back( 1.0 );
    }
  };

  //! Rule used to trigger action
  template< class eq > struct check_compflow : pegtl::success {};
  //! \brief Set defaults and do error checking on the compressible flow
  //!   equation block
  //! \details This is error checking that only the compressible flow equation
  //!   block must satisfy. Besides error checking we also set defaults here as
  //!   this block is called when parsing of a compflow...end block has
  //!   just finished.
  template< class eq >
  struct action< check_compflow< eq > > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::neq;
      using tag::param;

      // Error out if no dependent variable has been selected
      auto& depvar = stack.template get< param, eq, tag::depvar >();
      if (depvar.empty() || depvar.size() != neq.get< eq >())
        Message< Stack, ERROR, MsgKey::NODEPVAR >( stack, in );

      // If physics type is not given, default to 'euler'
      auto& physics = stack.template get< param, eq, tag::physics >();
      if (physics.empty() || physics.size() != neq.get< eq >()) {
        physics.push_back( inciter::ctr::PhysicsType::EULER );
      }

      // Set number of components to 5 (mass, 3 x mom, energy)
      stack.template get< tag::component, eq >().push_back( 5 );

      // Set default to sysfct (on/off) if not specified
      auto& sysfct = stack.template get< param, eq, tag::sysfct >();
      if (sysfct.empty() || sysfct.size() != neq.get< eq >())
        sysfct.push_back( 1 );

      // Set default flux to HLLC if not specified
      auto& flux = stack.template get< tag::param, eq, tag::flux >();
      if (flux.empty() || flux.size() != neq.get< eq >())
        flux.push_back( inciter::ctr::FluxType::HLLC );

      // Verify that sysfctvar variables are within bounds (if specified) and
      // defaults if not
      auto& sysfctvar = stack.template get< param, eq, tag::sysfctvar >();
      // If sysfctvar is not specified, use all variables for system FCT
      if (sysfctvar.empty() || sysfctvar.back().empty()) {
        sysfctvar.push_back( {0,1,2,3,4} );
      } else {  // if specified, do error checking on variables
        auto& vars = sysfctvar.back();
        if (vars.size() > 5) {
          Message< Stack, ERROR, MsgKey::SYSFCTVAR >( stack, in );
        }
        for (const auto& i : vars) {
          if (i > 4) Message< Stack, ERROR, MsgKey::SYSFCTVAR >( stack, in );
        }
      }

      // Verify correct number of material properties configured
      auto& matprop = stack.template get< param, eq, tag::material >().back()[0];
      auto& matidxmap = stack.template get< param, eq, tag::matidxmap >();
      matidxmap.template get< tag::eosidx >().resize(1);
      matidxmap.template get< tag::matidx >().resize(1);
      auto& meos = matprop.template get< tag::eos >();
      auto& mat_id = matprop.template get< tag::id >();

      if (mat_id.empty())
        mat_id.push_back(0);
      else if (mat_id.size() != 1)
        Message< Stack, ERROR, MsgKey::NUMMAT >( stack, in );
      else
        mat_id[0] = 0;

      if (meos == inciter::ctr::MaterialType::STIFFENEDGAS) {
        const auto& gamma = matprop.template get< tag::gamma >();
        // If gamma vector is wrong size, error out
        if (gamma.empty() || gamma.size() != 1)
          Message< Stack, ERROR, MsgKey::EOSGAMMA >( stack, in );

        auto& cv = matprop.template get< tag::cv >();
        // As a default, the specific heat of air (717.5 J/Kg-K) is used
        if (cv.empty())
          cv.push_back(717.5);
        // If specific heat vector is wrong size, error out
        if (cv.size() != 1)
          Message< Stack, ERROR, MsgKey::EOSCV >( stack, in );

        auto& pstiff = matprop.template get< tag::pstiff >();
        // As a default, a stiffness coefficient of 0.0 is used
        if (pstiff.empty())
          pstiff.push_back(0.0);
        // If stiffness coefficient vector is wrong size, error out
        if (pstiff.size() != 1)
          Message< Stack, ERROR, MsgKey::EOSPSTIFF >( stack, in );
      }

      // Generate mapping between material index and eos parameter index
      auto& eosmap = matidxmap.template get< tag::eosidx >();
      auto& idxmap = matidxmap.template get< tag::matidx >();
      eosmap[mat_id[0]] = static_cast< std::size_t >(matprop.template get<
        tag::eos >());
      idxmap[mat_id[0]] = 0;

      // If problem type is not given, default to 'user_defined'
      auto& problem = stack.template get< param, eq, tag::problem >();
      if (problem.empty() || problem.size() != neq.get< eq >())
        problem.push_back( inciter::ctr::ProblemType::USER_DEFINED );
      else if (problem.back() == inciter::ctr::ProblemType::VORTICAL_FLOW) {
        const auto& alpha = stack.template get< param, eq, tag::alpha >();
        const auto& beta = stack.template get< param, eq, tag::beta >();
        const auto& p0 = stack.template get< param, eq, tag::p0 >();
        if ( alpha.size() != problem.size() ||
             beta.size() != problem.size() ||
             p0.size() != problem.size() )
          Message< Stack, ERROR, MsgKey::VORTICAL_UNFINISHED >( stack, in );
      }
      else if (problem.back() == inciter::ctr::ProblemType::NL_ENERGY_GROWTH) {
        const auto& alpha = stack.template get< param, eq, tag::alpha >();
        const auto& betax = stack.template get< param, eq, tag::betax >();
        const auto& betay = stack.template get< param, eq, tag::betay >();
        const auto& betaz = stack.template get< param, eq, tag::betaz >();
        const auto& kappa = stack.template get< param, eq, tag::kappa >();
        const auto& r0 = stack.template get< param, eq, tag::r0 >();
        const auto& ce = stack.template get< param, eq, tag::ce >();
        if ( alpha.size() != problem.size() ||
             betax.size() != problem.size() ||
             betay.size() != problem.size() ||
             betaz.size() != problem.size() ||
             kappa.size() != problem.size() ||
             r0.size() != problem.size() ||
             ce.size() != problem.size() )
          Message< Stack, ERROR, MsgKey::ENERGY_UNFINISHED >( stack, in);
      }
      else if (problem.back() == inciter::ctr::ProblemType::RAYLEIGH_TAYLOR) {
        const auto& alpha = stack.template get< param, eq, tag::alpha >();
        const auto& betax = stack.template get< param, eq, tag::betax >();
        const auto& betay = stack.template get< param, eq, tag::betay >();
        const auto& betaz = stack.template get< param, eq, tag::betaz >();
        const auto& kappa = stack.template get< param, eq, tag::kappa >();
        const auto& p0 = stack.template get< param, eq, tag::p0 >();
        const auto& r0 = stack.template get< param, eq, tag::r0 >();
        if ( alpha.size() != problem.size() ||
             betax.size() != problem.size() ||
             betay.size() != problem.size() ||
             betaz.size() != problem.size() ||
             kappa.size() != problem.size() ||
             p0.size() != problem.size() ||
             r0.size() != problem.size() )
          Message< Stack, ERROR, MsgKey::RT_UNFINISHED >( stack, in);
      }

      // Error check on user-defined problem type
      auto& ic = stack.template get< param, eq, tag::ic >();
      auto& bgdensityic = ic.template get< tag::density >();
      auto& bgvelocityic = ic.template get< tag::velocity >();
      auto& bgpressureic = ic.template get< tag::pressure >();
      auto& bgenergyic = ic.template get< tag::energy >();
      auto& bgtemperatureic = ic.template get< tag::temperature >();
      if (problem.back() == inciter::ctr::ProblemType::USER_DEFINED) {
        // must have defined background ICs for user-defined ICs
        auto n = neq.get< eq >();
        if ( bgdensityic.size() != n || bgvelocityic.size() != n ||
             ( bgpressureic.size() != n && bgenergyic.size() != n &&
               bgtemperatureic.size() != n ) )
        {
          Message< Stack, ERROR, MsgKey::BGICMISSING >( stack, in );
        }

        // Error check Dirichlet boundary condition block for all compflow
        // configurations
        const auto& bc = stack.template get< param, eq, tag::bc, tag::bcdir >();
        for (const auto& s : bc)
          if (s.empty()) Message< Stack, ERROR, MsgKey::BC_EMPTY >( stack, in );

        // Error check stagnation BC block
        const auto& stag = stack.template get<tag::param, eq, tag::stag>();
        const auto& spoint = stag.template get< tag::point >();
        const auto& sradius = stag.template get< tag::radius >();
        if ( (!spoint.empty() && !spoint.back().empty() &&
              !sradius.empty() && !sradius.back().empty() &&
              spoint.back().size() != 3*sradius.back().size())
         || (!sradius.empty() && !sradius.back().empty() &&
              !spoint.empty() && !spoint.back().empty() &&
              spoint.back().size() != 3*sradius.back().size())
         || (!spoint.empty() && !spoint.back().empty() &&
              (sradius.empty() || (!sradius.empty() && sradius.back().empty())))
         || (!sradius.empty() && !sradius.back().empty() &&
              (spoint.empty() || (!spoint.empty() && spoint.back().empty()))) )
        {
          Message< Stack, ERROR, MsgKey::STAGBCWRONG >( stack, in );
        }

        // Error check skip BC block
        const auto& skip = stack.template get<tag::param, eq, tag::skip>();
        const auto& kpoint = skip.template get< tag::point >();
        const auto& kradius = skip.template get< tag::radius >();
        if ( (!kpoint.empty() && !kpoint.back().empty() &&
              !kradius.empty() && !kradius.back().empty() &&
              kpoint.back().size() != 3*kradius.back().size())
          || (!kradius.empty() && !kradius.back().empty() &&
              !kpoint.empty() && !kpoint.back().empty() &&
              kpoint.back().size() != 3*kradius.back().size())
          || (!kpoint.empty() && !kpoint.back().empty() &&
              (kradius.empty() || (!kradius.empty() && kradius.back().empty())))
          || (!kradius.empty() && !kradius.back().empty() &&
              (kpoint.empty() || (!kpoint.empty() && kpoint.back().empty()))) )
        {
          Message< Stack, ERROR, MsgKey::SKIPBCWRONG >( stack, in );
        }

        // Error check sponge BC parameter vectors for symmetry BC block
        const auto& sponge =
          stack.template get< tag::param, eq, tag::sponge >();
        const auto& ss = sponge.template get< tag::sideset >();

        const auto& spvel = sponge.template get< tag::velocity >();
        if ( !spvel.empty() && !spvel.back().empty()) {
          if (spvel.back().size() != ss.back().size())
            Message< Stack, ERROR, MsgKey::SPONGEBCWRONG >( stack, in );
          for (const auto& s : spvel.back())
            if ( s < 0.0 || s > 1.0 )
              Message< Stack, ERROR, MsgKey::SPONGEBCWRONG >( stack, in );
        }

        const auto& sppre = sponge.template get< tag::velocity >();
        if ( !sppre.empty() && !sppre.back().empty()) {
          if (sppre.back().size() != ss.back().size())
            Message< Stack, ERROR, MsgKey::SPONGEBCWRONG >( stack, in );
          for (const auto& s : sppre.back())
            if ( s < 0.0 || s > 1.0 )
              Message< Stack, ERROR, MsgKey::SPONGEBCWRONG >( stack, in );
        }

        // Error check user defined time dependent BC for this system
        const auto& tdepbc =
          stack.template get< tag::param, eq, tag::bctimedep >().back();
        // multiple time dependent BCs can be specified on different side sets
        for (const auto& bndry : tdepbc) {
          const auto& s = bndry.template get< tag::sideset >();
          if (s.empty()) Message< Stack, ERROR, MsgKey::BC_EMPTY >( stack, in );
          const auto& f = bndry.template get< tag::fn >();
          if (f.empty() or f.size() % 6 != 0)
            Message< Stack, ERROR, MsgKey::INCOMPLETEUSERFN>( stack, in );
        }
      }
    }
  };

  //! Rule used to trigger action
  template< class eq > struct check_multimat : pegtl::success {};
  //! \brief Set defaults and do error checking on the multimaterial
  //!    compressible flow equation block
  //! \details This is error checking that only the multimaterial compressible
  //!   flow equation block must satisfy. Besides error checking we also set
  //!   defaults here as this block is called when parsing of a
  //!   multimat...end block has just finished.
  template< class eq >
  struct action< check_multimat< eq > > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::neq;
      using tag::param;

      // Error out if no dependent variable has been selected
      auto& depvar = stack.template get< param, eq, tag::depvar >();
      if (depvar.empty() || depvar.size() != neq.get< eq >())
        Message< Stack, ERROR, MsgKey::NODEPVAR >( stack, in );

      // If physics type is not given, default to 'veleq'
      auto& physics = stack.template get< param, eq, tag::physics >();
      if (physics.empty() || physics.size() != neq.get< eq >())
        physics.push_back( inciter::ctr::PhysicsType::VELEQ );

      // Set default flux to AUSM if not specified
      auto& flux = stack.template get< tag::param, eq, tag::flux >();
      if (flux.empty() || flux.size() != neq.get< eq >())
        flux.push_back( inciter::ctr::FluxType::AUSM );

      // Set number of scalar components based on number of materials
      auto& nmat = stack.template get< param, eq, tag::nmat >();
      auto& ncomp = stack.template get< tag::component, eq >();
      if (physics.back() == inciter::ctr::PhysicsType::VELEQ) {
        // physics = veleq: m-material compressible flow
        // scalar components: volfrac:m + mass:m + momentum:3 + energy:m
        // if nmat is unspecified, configure it be 2
        if (nmat.empty() || nmat.size() != neq.get< eq >()) {
          Message< Stack, WARNING, MsgKey::NONMAT >( stack, in );
          nmat.push_back( 2 );
        }
        // set ncomp based on nmat
        auto m = nmat.back();
        ncomp.push_back( m + m + 3 + m );
      }

      // Verify correct number of multi-material properties (gamma, cv, pstiff)
      // have been configured
      auto& matprop = stack.template get< param, eq, tag::material >();
      auto& matidxmap = stack.template get< param, eq, tag::matidxmap >();
      matidxmap.template get< tag::eosidx >().resize(nmat.back());
      matidxmap.template get< tag::matidx >().resize(nmat.back());
      std::size_t tmat(0), i(0);
      std::set< std::size_t > matidset;

      for (auto& mtype : matprop.back()) {
        const auto& meos = mtype.template get< tag::eos >();
        const auto& mat_id = mtype.template get< tag::id >();

        if (meos == inciter::ctr::MaterialType::STIFFENEDGAS) {
          const auto& gamma = mtype.template get< tag::gamma >();
          // If gamma vector is wrong size, error out
          if (gamma.empty() || gamma.size() != mat_id.size())
            Message< Stack, ERROR, MsgKey::EOSGAMMA >( stack, in );

          auto& cv = mtype.template get< tag::cv >();
          // As a default, the specific heat of air (717.5 J/Kg-K) is used
          if (cv.empty()) {
            for (std::size_t k=0; k<mat_id.size(); ++k) {
              cv.push_back(717.5);
            }
          }
          // If specific heat vector is wrong size, error out
          if (cv.size() != mat_id.size())
            Message< Stack, ERROR, MsgKey::EOSCV >( stack, in );

          auto& pstiff = mtype.template get< tag::pstiff >();
          // As a default, a stiffness coefficient of 0.0 is used
          if (pstiff.empty()) {
            for (std::size_t k=0; k<mat_id.size(); ++k) {
              pstiff.push_back(0.0);
            }
          }
          // If stiffness coefficient vector is wrong size, error out
          if (pstiff.size() != mat_id.size())
            Message< Stack, ERROR, MsgKey::EOSPSTIFF >( stack, in );
        }

        // Track total number of materials in multiple material blocks
        tmat += mat_id.size();

        // Check for repeating user specified material ids
        for (auto midx : mat_id) {
          if (!matidset.count(midx))
            matidset.insert(midx);
          else
            Message< Stack, ERROR, MsgKey::REPMATID >( stack, in );
        }

        // Generate mapping between material index and eos parameter index
        auto& eosmap = matidxmap.template get< tag::eosidx >();
        auto& idxmap = matidxmap.template get< tag::matidx >();
        for (auto midx : mat_id) {
          midx -= 1;
          eosmap[midx] = static_cast< std::size_t >(mtype.template get<
            tag::eos >());
          idxmap[midx] = i;
          ++i;
        }
        // end of materials for this eos, thus reset index counter
        i = 0;
      }

      // If total number of materials is incorrect, error out
      if (tmat != nmat.back())
        Message< Stack, ERROR, MsgKey::NUMMAT >( stack, in );

      // Check if material ids are contiguous and 1-based
      if (!matidset.count(1))
        Message< Stack, ERROR, MsgKey::ONEMATID >( stack, in );
      std::size_t icount(1);
      for (auto midx : matidset) {
        if (midx != icount)
          Message< Stack, ERROR, MsgKey::GAPMATID >( stack, in );
        ++icount;
      }

      // If pressure relaxation is not specified, default to 'false'
      auto& prelax = stack.template get< param, eq, tag::prelax >();
      if (prelax.empty() || prelax.size() != neq.get< eq >())
        prelax.push_back( 0 );

      // If pressure relaxation time-scale is not specified, default to 1.0
      auto& prelax_ts = stack.template get< param, eq,
                                            tag::prelax_timescale >();
      if (prelax_ts.empty() || prelax_ts.size() != neq.get< eq >())
        prelax_ts.push_back( 1.0 );

      // If interface compression is not specified, default to 'false'
      auto& intsharp = stack.template get< param, eq, tag::intsharp >();
      if (intsharp.empty() || intsharp.size() != neq.get< eq >())
        intsharp.push_back( 0 );

      // If interface compression parameter is not specified, default to 1.0
      auto& intsharp_p = stack.template get< param, eq,
                                            tag::intsharp_param >();
      if (intsharp_p.empty() || intsharp_p.size() != neq.get< eq >())
        intsharp_p.push_back( 1.0 );

      // If problem type is not given, default to 'user_defined'
      auto& problem = stack.template get< param, eq, tag::problem >();
      if (problem.empty() || problem.size() != neq.get< eq >())
        problem.push_back( inciter::ctr::ProblemType::USER_DEFINED );
      else if (problem.back() == inciter::ctr::ProblemType::VORTICAL_FLOW) {
        const auto& alpha = stack.template get< param, eq, tag::alpha >();
        const auto& beta = stack.template get< param, eq, tag::beta >();
        const auto& p0 = stack.template get< param, eq, tag::p0 >();
        if ( alpha.size() != problem.size() ||
             beta.size() != problem.size() ||
             p0.size() != problem.size() )
          Message< Stack, ERROR, MsgKey::VORTICAL_UNFINISHED >( stack, in );
      }

      // Error check on user-defined problem type
      auto& ic = stack.template get< param, eq, tag::ic >();
      auto& bgmatid = ic.template get< tag::materialid >();
      auto& bgdensityic = ic.template get< tag::density >();
      auto& bgvelocityic = ic.template get< tag::velocity >();
      auto& bgpressureic = ic.template get< tag::pressure >();
      auto& bgenergyic = ic.template get< tag::energy >();
      auto& bgtemperatureic = ic.template get< tag::temperature >();
      if (problem.back() == inciter::ctr::ProblemType::USER_DEFINED) {
        // must have defined background ICs for user-defined ICs
        auto n = neq.get< eq >();
        if (bgmatid.size() != n) {
          Message< Stack, ERROR, MsgKey::BGMATIDMISSING >( stack, in );
        }

        if ( bgdensityic.size() != n || bgvelocityic.size() != n ||
             ( bgpressureic.size() != n && bgenergyic.size() != n &&
               bgtemperatureic.size() != n ) )
        {
          Message< Stack, ERROR, MsgKey::BGICMISSING >( stack, in );
        }

        // each IC box should have material id specified, and it should be
        // within nmat
        auto& icbox = ic.template get< tag::box >();

        if (!icbox.empty()) {
          for (const auto& b : icbox.back()) {   // for all boxes
            auto boxmatid = b.template get< tag::materialid >();
            if (boxmatid == 0) {
              Message< Stack, ERROR, MsgKey::BOXMATIDMISSING >( stack, in );
            }
            else if (boxmatid > nmat.back()) {
              Message< Stack, ERROR, MsgKey::BOXMATIDWRONG >( stack, in );
            }
          }
        }
      }

      // Error check Dirichlet boundary condition block for all multimat
      // configurations
      const auto& bc = stack.template get< param, eq, tag::bc, tag::bcdir >();
      for (const auto& s : bc)
        if (s.empty()) Message< Stack, ERROR, MsgKey::BC_EMPTY >( stack, in );
    }
  };

  //! Rule used to trigger action
  template< class Option, typename...tags >
  struct store_inciter_option : pegtl::success {};
  //! \brief Put option in state at position given by tags
  //! \details This is simply a wrapper around tk::grm::store_option passing the
  //!    stack defaults for inciter.
  template< class Option, typename... tags >
  struct action< store_inciter_option< Option, tags... > > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      store_option< Stack, inciter::deck::use, Option, inciter::ctr::InputDeck,
                    Input, tags... >
                  ( stack, in, inciter::g_inputdeck_defaults );
    }
  };

  //! Function object to ensure disjoint side sets for all boundary conditions
  //! \details This is instantiated using a Cartesian product of all PDE types
  //!    and all BC types at compile time. It goes through all side sets
  //!    configured by the user and triggers an error if a side set is assigned
  //!    a BC more than once (within a solver).
  template< typename Input, typename Stack >
  struct ensure_disjoint {
    const Input& m_input;
    Stack& m_stack;
    explicit ensure_disjoint( const Input& in, Stack& stack ) :
      m_input( in ), m_stack( stack ) {}
    template< typename U > void operator()( brigand::type_<U> ) {
      using Eq = typename brigand::front< U >;
      using BC = typename brigand::back< U >;
      const auto& bc = m_stack.template get< tag::param, Eq, tag::bc, BC >();
      for (const auto& eq : bc) {
        std::unordered_set< int > bcset;
        for (const auto& s : eq) {
          auto id = std::stoi(s);
          if (bcset.find(id) != end(bcset))
            Message< Stack, ERROR, MsgKey::NONDISJOINTBC >( m_stack, m_input );
          else
            bcset.insert( id );
        }
      }
    }
  };

  //! Function object to count the number of meshes assigned to solvers
  //! \details This is instantiated for all PDE types at compile time. It goes
  //!   through all configured solvers (equation system configuration blocks)
  //!   and counts the number of mesh filenames configured.
  template< typename Stack >
  struct count_meshes {
    const Stack& stack;
    std::size_t& count;
    explicit
    count_meshes( const Stack& s, std::size_t& c ) : stack(s), count(c) {}
    template< typename eq > void operator()( brigand::type_<eq> ) {
      count +=
        stack.template get< tag::param, eq, tag::mesh, tag::filename >().size();
    }
  };

  //! Function object to assign mesh ids to solvers
  //! \details This is instantiated for all PDE types at compile time. It goes
  //!   through all configured solvers (equation system configuration blocks)
  //!   and assigns a new mesh id to all solvers configured in the input file.
  template< typename Stack >
  struct assign_meshid {
    Stack& stack;
    std::size_t& meshid;
    explicit assign_meshid( Stack& s, std::size_t& m ) : stack(s), meshid(m) {}
    template< typename eq > void operator()( brigand::type_<eq> ) {
      const auto& eq_mesh_filename =
        stack.template get< tag::param, eq, tag::mesh, tag::filename >();
      auto& id = stack.template get< tag::param, eq, tag::mesh, tag::id >();
      for (std::size_t i=0; i<eq_mesh_filename.size(); ++i)
        id.push_back( meshid++ );
    }
  };

  //! Rule used to trigger action
  struct configure_scheme : pegtl::success {};
  //! Configure scheme selected by user
  //! \details This grammar action configures the number of degrees of freedom
  //! (NDOF) used for DG methods. For finite volume (or DGP0), the DOF are the
  //! cell-averages. This implies ndof=1 for DGP0. Similarly, ndof=4 and 10 for
  //! DGP1 and DGP2 respectively, since they evolve higher (>1) order solution
  //! information (e.g. gradients) as well. "rdof" includes degrees of freedom
  //! that are both, evolved and reconstructed. For rDGPnPm methods (e.g. P0P1
  //! and P1P2), "n" denotes the evolved solution-order and "m" denotes the
  //! reconstructed solution-order; i.e. P0P1 has ndof=1 and rdof=4, whereas
  //! P1P2 has ndof=4 and rdof=10. For a pure DG method without reconstruction
  //! (DGP0, DGP1, DGP2), rdof=ndof. For more information about rDGPnPm methods,
  //! ref. Luo, H. et al. (2013).  A reconstructed discontinuous Galerkin method
  //! based on a hierarchical WENO reconstruction for compressible flows on
  //! tetrahedral grids. Journal of Computational Physics, 236, 477-492.
  template<> struct action< configure_scheme > {
    template< typename Input, typename Stack >
    static void apply( const Input&, Stack& stack ) {
      using inciter::ctr::SchemeType;
      auto& discr = stack.template get< tag::discr >();
      auto& ndof = discr.template get< tag::ndof >();<--- Variable 'ndof' is assigned a value that is never used.
      auto& rdof = discr.template get< tag::rdof >();<--- Variable 'rdof' is assigned a value that is never used.
      auto scheme = discr.template get< tag::scheme >();
      if (scheme == SchemeType::P0P1) {
        ndof = 1; rdof = 4;
      } else if (scheme == SchemeType::DGP1) {
        ndof = rdof = 4;
      } else if (scheme == SchemeType::DGP2) {
        ndof = rdof = 10;
      } else if (scheme == SchemeType::PDG) {
        ndof = rdof = 10;
        stack.template get< tag::pref, tag::pref >() = true;
      }
    }
  };

  //! Function object to do error checking on output time ranges
  template< typename Stack, typename Input >
  struct range_errchk {
    Stack& stack;
    const Input& input;
    explicit range_errchk( Stack& s, const Input& in ) : stack(s), input(in) {}
    template< typename U > void operator()( brigand::type_<U> ) {
      for (const auto& r : stack.template get< tag::output, tag::range, U >())
        if ( r.size() != 3 or r[0] > r[1] or r[2] < 0.0 or r[2] > r[1]-r[0] )
          Message< Stack, ERROR, MsgKey::BADRANGE >( stack, input );
    }
  };

  //! Rule used to trigger action
  struct check_inciter : pegtl::success {};
  //! \brief Do error checking on the inciter block
  template<> struct action< check_inciter > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::neq;
      using inciter::g_inputdeck_defaults;

      // Error out if no dt policy has been selected
      const auto& dt = stack.template get< tag::discr, tag::dt >();
      const auto& cfl = stack.template get< tag::discr, tag::cfl >();
      if ( std::abs(dt - g_inputdeck_defaults.get< tag::discr, tag::dt >()) <
            std::numeric_limits< tk::real >::epsilon() &&
          std::abs(cfl - g_inputdeck_defaults.get< tag::discr, tag::cfl >()) <
            std::numeric_limits< tk::real >::epsilon() )
        Message< Stack, ERROR, MsgKey::NODT >( stack, in );

      // If both dt and cfl are given, warn that dt wins over cfl
      if ( std::abs(dt - g_inputdeck_defaults.get< tag::discr, tag::dt >()) >
            std::numeric_limits< tk::real >::epsilon() &&
          std::abs(cfl - g_inputdeck_defaults.get< tag::discr, tag::cfl >()) >
            std::numeric_limits< tk::real >::epsilon() )
        Message< Stack, WARNING, MsgKey::MULDT >( stack, in );

      // Do error checking on time history points
      const auto& hist = stack.template get< tag::history, tag::point >();
      if (std::any_of( begin(hist), end(hist),
           [](const auto& p){ return p.size() != 3; } ) )
      {
        Message< Stack, ERROR, MsgKey::WRONGSIZE >( stack, in );
      }

      // Do error checking on residual eq system component index
      const auto rc = stack.template get< tag::discr, tag::rescomp >();
      const auto& ncomps = stack.template get< tag::component >();
      if (rc < 1 || rc > ncomps.nprop())
        Message< Stack, ERROR, MsgKey::LARGECOMP >( stack, in );

      // Ensure no different BC types are assigned to the same side set
      using PDETypes = inciter::ctr::parameters::Keys;
      using BCTypes = inciter::ctr::bc::Keys;
      brigand::for_each< tk::cartesian_product< PDETypes, BCTypes > >(
        ensure_disjoint< Input, Stack >( in, stack ) );

      // Do error checking on output time range configuration parameters: they
      // all must be a 3 reals: mintime, maxtime, and dt with maxtime >
      // mintime, and dt<maxtime-mintime.
      brigand::for_each< inciter::ctr::time_range::Keys >
                       ( range_errchk< Stack, Input >( stack, in ) );

      // Do error checking on time history point names (this is a programmer
      // error if triggers, hence assert)
      Assert(
        (stack.template get< tag::history, tag::id >().size() == hist.size()),
        "Number of history points and ids must equal" );

      // If at least a mesh filename is assigned to a solver, all solvers must
      // have a mesh filename assigned
      std::size_t nmesh = 0;
      brigand::for_each< PDETypes >( count_meshes< Stack >( stack, nmesh ) );
      if (nmesh > 0 && nmesh != depvars.size())
        Message< Stack, ERROR, MsgKey::MULTIMESH >( stack, in );

      // Remove duplicate transfer steps
      tk::unique( stack.template get< tag::couple, tag::transfer >() );

      // Now that the inciter ... end block is finished, assign mesh ids to
      // solvers configured
      std::size_t meshid = 0;
      brigand::for_each< PDETypes >( assign_meshid< Stack >( stack, meshid ) );
      Assert( meshid == nmesh, "Not all meshes configured have mesh ids" );
    }
  };

  //! Rule used to trigger action
  struct check_ale : pegtl::success {};
  //! \brief Do error checking on the inciter block
  template<> struct action< check_ale > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::g_inputdeck_defaults;

     // Trigger error if steady state + ALE are both enabled
      auto steady = stack.template get< tag::discr, tag::steady_state >();
      auto ale = stack.template get< tag::ale, tag::ale >();
      if (steady && ale) {
        Message< Stack, ERROR, MsgKey::STEADYALE >( stack, in );
      }

      // Set a sensible default for dvCFL if ALE is enabled and if dvcfl not set
      auto& dvcfl = stack.template get< tag::ale, tag::dvcfl >();
      auto dvcfl_default = g_inputdeck_defaults.get< tag::ale, tag::dvcfl >();
      auto eps = std::numeric_limits< tk::real >::epsilon();
      if (ale && std::abs(dvcfl - dvcfl_default) < eps) dvcfl = 0.01;

      // Set a default of zeros for the mesh force ALE parameters
      auto& meshforce = stack.template get< tag::ale, tag::meshforce >();
      if (ale && meshforce.size() != 4) meshforce = { 0, 0, 0, 0 };

      // Set a default for the ALE mesh motion dimensions
      auto& mesh_motion = stack.template get< tag::ale, tag::mesh_motion >();
      if (ale && mesh_motion.empty()) mesh_motion = { 0, 1, 2 };

      // Error out if mesh motion dimensions are wrong
      if (ale && (mesh_motion.size() > 3 ||
                  std::any_of( begin(mesh_motion), end(mesh_motion),
                     [](auto d){return d > 2;} )) )
      {
        Message< Stack, ERROR, MsgKey::WRONGMESHMOTION >( stack, in );
      }

      // Error checking on user-defined function for ALE's moving sides
      const auto& move = stack.template get< tag::ale, tag::move >();
      for (const auto& s : move) {
        const auto& f = s.template get< tag::fn >();
        if (f.empty() or f.size() % 4 != 0)
          Message< Stack, ERROR, MsgKey::INCOMPLETEUSERFN>( stack, in );
      }
    }
  };

  //! Rule used to trigger action
  template< typename Feature >
  struct enable : pegtl::success {};
  //! Enable adaptive mesh refinement (AMR)
  template< typename Feature >
  struct action< enable< Feature > > {
    template< typename Input, typename Stack >
    static void apply( const Input&, Stack& stack ) {
      stack.template get< Feature, Feature >() = true;
    }
  };

  //! Rule used to trigger action
  struct compute_refvar_idx : pegtl::success {};
  //! Compute indices of refinement variables
  //! \details This functor computes the indices in the unknown vector for all
  //!   refinement variables in the system of systems of dependent variables
  //!   after the refvar...end block has been parsed in the amr...end block.
  //!   After basic error checking, the vector at stack.get<tag::amr,tag::id>()
  //!   is filled.
  template<>
  struct action< compute_refvar_idx > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      // reference variables just parsed by refvar...end block
      const auto& refvar = stack.template get< tag::amr, tag::refvar >();
      // get ncomponents object from this input deck
      const auto& ncomps = stack.template get< tag::component >();
      // compute offset map associating offsets to dependent variables
      auto offsetmap = ncomps.offsetmap( stack );
      // compute number of components associated to dependent variabels
      auto ncompmap = ncomps.ncompmap( stack );
      // reference variable index vector to fill
      auto& refidx = stack.template get< tag::amr, tag::id >();
      // Compute indices for all refvars
      for (const auto& v : refvar) {    // for all reference variables parsed
        // depvar is the first char of a refvar
        auto depvar = v[0];
        // the field ID is optional and is the rest of the depvar string
        std::size_t f = (v.size()>1 ? std::stoul(v.substr(1)) : 1) - 1;
        // field ID must be less than or equal to the number of scalar
        // components configured for the eq system for this dependent variable
        if (f >= tk::cref_find( ncompmap, depvar ))
          Message< Stack, ERROR, MsgKey::NOSUCHCOMPONENT >( stack, in );
        // get offset for depvar
        auto eqsys_offset = tk::cref_find( offsetmap, depvar );
        // the index is the eq offset + field ID
        auto idx = eqsys_offset + f;
        // save refvar index in system of all systems
        refidx.push_back( idx );
      }
    }
  };

  //! Rule used to trigger action
  struct check_amr_errors : pegtl::success {};
  //! Do error checking for the amr...end block
  //! \details This is error checking that only the amr...end block
  //!   must satisfy. Besides error checking this can also set defaults
  //!   as this block is called when parsing of a amr...end block has
  //!   just finished.
  template<>
  struct action< check_amr_errors > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      // Error out if refvar size does not equal refidx size (programmer error)
      Assert( (stack.template get< tag::amr, tag::refvar >().size() ==
               stack.template get< tag::amr, tag::id >().size()),
              "The size of refvar and refidx vectors must equal" );
      const auto& initref = stack.template get< tag::amr, tag::init >();
      const auto& refvar = stack.template get< tag::amr, tag::refvar >();
      const auto& edgelist = stack.template get< tag::amr, tag::edge >();
      // Error out if initref edge list is not divisible by 2 (user error)
      if (edgelist.size() % 2 == 1)
        Message< Stack, ERROR, MsgKey::T0REFODD >( stack, in );
      // Warn if initial AMR will be a no-op
      if ( stack.template get< tag::amr, tag::t0ref >() && initref.empty() )
        Message< Stack, WARNING, MsgKey::T0REFNOOP >( stack, in );
      // Error out if timestepping AMR will be a no-op (user error)
      if ( stack.template get< tag::amr, tag::dtref >() && refvar.empty() )
        Message< Stack, ERROR, MsgKey::DTREFNOOP >( stack, in );
      // Error out if mesh refinement frequency is zero (programmer error)
      Assert( (stack.template get< tag::amr, tag::dtfreq >() > 0),
              "Mesh refinement frequency must be positive" );
    }
  };

  //! Rule used to trigger action
  struct check_pref_errors : pegtl::success {};
  //! Do error checking for the pref...end block
  template<>
  struct action< check_pref_errors > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      auto& tolref = stack.template get< tag::pref, tag::tolref >();<--- Variable 'tolref' can be declared with const
      if (tolref < 0.0 || tolref > 1.0)
        Message< Stack, ERROR, MsgKey::PREFTOL >( stack, in );
    }
  };

  //! Rule used to trigger action
  struct match_pointname : pegtl::success {};
  //! \brief Match PDF name to the registered ones
  //! \details This is used to check the set of PDF names dependent previously
  //!    registered to make sure all are unique.
  template<>
  struct action< match_pointname > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::pointnames;
      // find matched name in set of registered ones
      if (pointnames.find( in.string() ) == pointnames.end()) {
        pointnames.insert( in.string() );
        stack.template get< tag::history, tag::id >().push_back( in.string() );
      }
      else  // error out if name matched var is already registered
        Message< Stack, ERROR, MsgKey::POINTEXISTS >( stack, in );
    }
  };

  //! Rule used to trigger action
  struct push_depvar : pegtl::success {};
  //! Add matched outvar based on depvar into vector of outvars
  //! \details Push outvar based on depvar: use first char of matched token as
  //! OutVar::var, OutVar::name = "" by default. OutVar::name being empty will
  //! be used to differentiate a depvar-based outvar from a human-readable
  //! outvar. Depvar-based outvars can directly access solution arrays using
  //! their field. Human-readable outvars need a mechanism (a function) to read
  //! and compute their variables from solution arrays. The 'getvar' function,
  //! used to compute a physics variable from the numerical solution is assigned
  //! after initial migration and thus not assigned here (during parsing).
  template<>
  struct action< push_depvar > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::centering;
      using inciter::ctr::OutVar;
      auto& vars = stack.template get< tag::cmd, tag::io, tag::outvar >();
      vars.emplace_back(OutVar(in.string()[0], field, centering));
      field = 0;        // reset field
    }
  };

  //! Rule used to trigger action
  struct push_matvar : pegtl::success {};
  //! Add matched outvar based on matvar into vector of outvars
  //! \details Push outvar based on matvar: use depvar char as OutVar::var,
  //! OutVar::name = "" by default. Matvar-based outvars are similar to
  //! depvar-base outvars, in that the OutVar has empty name and the OutVar::var
  //! is a depvar, but instead of having the user try to guess the field id, the
  //! grammar accepts a physics label (accepted multimatvars) and a material
  //! index, which are then converted to a depvar + field index.
  template<>
  struct action< push_matvar > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using tag::param;
      using tag::multimat;
      using inciter::deck::centering;
      using inciter::deck::multimatvars;
      using inciter::ctr::OutVar;
      auto& vars = stack.template get< tag::cmd, tag::io, tag::outvar >();
      auto nmat = stack.template get< param, multimat, tag::nmat >().back();
      auto depvar = stack.template get< param, multimat, tag::depvar >().back();
      // first char of matched token: accepted multimatvar label char
      char v = static_cast<char>( in.string()[0] );
      // Since multimat outvars are configured based on an acceptable character
      // label (see inciter::deck::multimatvars, denoting a physics variable),
      // and a material index (instead of a depvar + a component index),
      // multimat material bounds are checked here. Note that for momentum and
      // velocity, the field id is the spatial direction not the material id.
      // Also note that field (in grammar's state) starts from 0.
      if ( ((v=='u'||v=='U'||v=='m'||v=='M') && field>2) ||
           ((v!='u'&&v!='U'&&v!='m'&&v!='M') && field>=nmat) )
        Message< Stack, ERROR, MsgKey::NOSUCHCOMPONENT >( stack, in );
      // field contains material id, compute multiat component index
      auto comp = tk::cref_find( multimatvars, v )( nmat, field );
      // save depvar + component index based on physics label + material id,
      // also save physics label + material id as matvar
      vars.emplace_back(
        OutVar(depvar, comp, centering, {}, {}, v+std::to_string(field+1)) );
      field = 0;        // reset field
    }
  };

  //! Function object for adding a human-readable output variable
  //! \details Since human-readable outvars do not necessarily have any
  //! reference to the depvar of their system they refer to, nor which system
  //! they refer to, we configure them for all of the systems they are preceded
  //! by. If there is only a single system of the type the outvar is configured,
  //! we simply look up the depvar and use that as OutVar::var. If there are
  //! multiple systems configured upstream to which the outvar could refer to,
  //! we configure an outvar for all systems configured, and postfix the
  //! human-readable OutVar::name with '_' + depvar. Hence this function object
  //! so the code below can be invoked for all equation types.
  template< typename Stack >
  struct AddOutVarHuman {
    Stack& stack;
    const std::string& in_string;
    explicit AddOutVarHuman( Stack& s, const std::string& ins )
      : stack(s), in_string(ins) {}
    template< typename Eq > void operator()( brigand::type_<Eq> ) {
      using inciter::deck::centering;
      using inciter::ctr::OutVar;
      const auto& depvar = stack.template get< tag::param, Eq, tag::depvar >();
      auto& vars = stack.template get< tag::cmd, tag::io, tag::outvar >();
      if (depvar.size() == 1)
        vars.emplace_back( OutVar( depvar[0], 0, centering, in_string ) );
      else
        for (auto d : depvar)
          vars.emplace_back( OutVar( d, 0, centering, in_string + '_' + d ) );<--- Consider using std::transform algorithm instead of a raw loop.
    }
  };

  //! Rule used to trigger action
  struct push_humanvar : pegtl::success {};
  //! Add matched outvar based on depvar into vector of vector of outvars
  //! \details Push outvar based on human readable string for which
  //! OutVar::name = matched token. OutVar::name being not empty will be used to
  //! differentiate a depvar-based outvar from a human-readable outvar.
  //! Depvar-based outvars can directly access solution arrays using their
  //! field. Human-readable outvars need a mechanism (a function) to read and
  //! compute their variables from solution arrays. The 'getvar' function, used
  //! to compute a physics variable from the numerical solution is assigned
  //! after initial migration and thus not assigned here (during parsing).
  template<>
  struct action< push_humanvar > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      brigand::for_each< inciter::ctr::parameters::Keys >
                       ( AddOutVarHuman< Stack >( stack, in.string() ) );
    }
  };

  //! Rule used to trigger action
  struct set_outvar_alias : pegtl::success {};
  //! Set alias of last pushed output variable
  template<>
  struct action< set_outvar_alias > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      // Set alias of last pushed outvar:
      auto& vars = stack.template get< tag::cmd, tag::io, tag::outvar >();
      if (!vars.empty()) vars.back().alias = in.string();
    }
  };

  //! Function object for error checking outvar bounds for each equation type
  template< typename Stack >
  struct OutVarBounds {
    const Stack& stack;
    bool& inbounds;
    explicit OutVarBounds( const Stack& s, bool& i )
      : stack(s), inbounds(i) { inbounds = false; }
    template< typename U > void operator()( brigand::type_<U> ) {
      if (std::is_same_v< U, tag::multimat >) inbounds = true;  // skip multimat
      const auto& depvar = stack.template get< tag::param, U, tag::depvar >();
      const auto& ncomp = stack.template get< tag::component, U >();
      Assert( depvar.size() == ncomp.size(), "Size mismatch" );
      // called after matching each outvar, so only check the last one
      auto& vars = stack.template get< tag::cmd, tag::io, tag::outvar >();
      const auto& last_outvar = vars.back();
      const auto& v = static_cast<char>( std::tolower(last_outvar.var) );
      for (std::size_t e=0; e<depvar.size(); ++e)
        if (v == depvar[e] && last_outvar.field < ncomp[e]) inbounds = true;
    }
  };

  //! Rule used to trigger action
  struct check_outvar : pegtl::success {};
  //! Bounds checking for output variables at the end of a var ... end block
  template<>
  struct action< check_outvar > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      bool inbounds;
      brigand::for_each< inciter::ctr::parameters::Keys >
                       ( OutVarBounds< Stack >( stack, inbounds ) );
      if (!inbounds)
        Message< Stack, ERROR, MsgKey::NOSUCHCOMPONENT >( stack, in );
    }
  };

  //! Rule used to trigger action
  struct set_centering : pegtl::success {};
  //! Set variable centering in parser's state
  template<>
  struct action< set_centering > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& ) {
      inciter::deck::centering =
        (in.string() == "node") ? tk::Centering::NODE : tk::Centering::ELEM;
    }
  };

  //! Rule used to trigger action
  struct match_outvar : pegtl::success {};
  //! Match output variable based on depvar
  template<>
  struct action< match_outvar > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      using inciter::deck::neq;
      using inciter::deck::multimatvars;
      // convert matched string to char
      auto var = stack.template convert< char >( in.string() );
      if (neq.get< tag::multimat >() == 0) {    // if not multimat
        // find matched variable in set of selected ones
        if (depvars.find(var) != end(depvars))
          action< push_depvar >::apply( in, stack );
        else  // error out if matched var is not selected
          Message< Stack, ERROR, MsgKey::NOSUCHOUTVAR >( stack, in );
      } else {    // if multimat
        // find matched variable in set accepted for multimat
        if (multimatvars.find(var) != end(multimatvars))
          action< push_matvar >::apply( in, stack );
        else
          Message< Stack, ERROR, MsgKey::NOSUCHMULTIMATVAR >( stack, in );
      }
    }
  };

  // Store mesh/solver id as a source of a transfer
  template< typename Stack > struct store_transfer_src {
    store_transfer_src( Stack& stack, std::size_t i ) {
      stack.template get< tag::couple, tag::transfer >().emplace_back( i, 0 );
    }
  };

  // Store mesh/solver id as a destination of a transfer
  template< typename Stack > struct store_transfer_dst {
    store_transfer_dst( Stack& stack, std::size_t i ) {
      stack.template get< tag::couple, tag::transfer >().back().dst = i;
    }
  };

  //! Rule used to trigger action
  template< template< class > class StoreTransfer >
  struct push_transfer : pegtl::success {};
  //! Add matched value as a source or destination of solution transfer
  //! \tparam StoreTransfer Type of action to invoke: source or destination
  template< template< class > class StoreTransfer >
  struct action< push_transfer< StoreTransfer > > {
    template< typename Input, typename Stack >
    static void apply( const Input& in, Stack& stack ) {
      // Extract dependent variables for solvers configured
      auto depvar = stack.depvar();
      // Store index of parsed depvar of transfer being configured
      auto c = in.string()[0];
      for (std::size_t i=0; i<depvar.size(); ++i)
        if (depvar[i] == c)
          StoreTransfer< Stack >( stack, i );
    }
  };

} // ::grm
} // ::tk

namespace inciter {

//! Inciter input deck facilitating user input for computing shock hydrodynamics
namespace deck {

  using namespace tao;

  // Inciter's InputDeck grammar

  //! scan and store_back equation keyword and option
  template< typename keyword, class eq >
  struct scan_eq :
         tk::grm::scan< typename keyword::pegtl_string,
                        tk::grm::store_back_option< use,
                                                    ctr::PDE,
                                                    tag::selected,
                                                    tag::pde > > {};

  //! Error checks after an equation...end block has been parsed
  template< class eq, template< class > class eqchecker >
  struct check_errors :
         pegtl::seq<
           // register differential equation block
           tk::grm::register_inciter_eq< eq >,
           // check mesh ... end block
           tk::grm::check_mesh< eq >,
           // do error checking on this block
           eqchecker< eq > > {};

  //! Match discretization option
  template< template< class > class use, class keyword, class Option,
            class Tag >
  struct discroption :
         tk::grm::process< use< keyword >,
                           tk::grm::store_inciter_option<
                             Option, tag::discr, Tag >,
                           pegtl::alpha > {};

  //! Discretization parameters
  struct discretization :
         pegtl::sor<
           tk::grm::discrparam< use, kw::nstep, tag::nstep >,
           tk::grm::discrparam< use, kw::term, tag::term >,
           tk::grm::discrparam< use, kw::t0, tag::t0 >,
           tk::grm::discrparam< use, kw::dt, tag::dt >,
           tk::grm::discrparam< use, kw::cfl, tag::cfl >,
           tk::grm::discrparam< use, kw::residual, tag::residual >,
           tk::grm::discrparam< use, kw::rescomp, tag::rescomp >,
           tk::grm::process< use< kw::fcteps >,
                             tk::grm::Store< tag::discr, tag::fcteps > >,
           tk::grm::process< use< kw::fctclip >,
                             tk::grm::Store< tag::discr, tag::fctclip >,
                             pegtl::alpha >,
           tk::grm::process< use< kw::fct >,
                             tk::grm::Store< tag::discr, tag::fct >,
                             pegtl::alpha >,
           tk::grm::process< use< kw::ctau >,
                             tk::grm::Store< tag::discr, tag::ctau > >,
           tk::grm::process< use< kw::pelocal_reorder >,
                             tk::grm::Store< tag::discr, tag::pelocal_reorder >,
                             pegtl::alpha >,
           tk::grm::process< use< kw::operator_reorder >,
                             tk::grm::Store< tag::discr, tag::operator_reorder >,
                             pegtl::alpha >,
           tk::grm::process< use< kw::steady_state >,
                             tk::grm::Store< tag::discr, tag::steady_state >,
                             pegtl::alpha >,
           tk::grm::interval_iter< use< kw::ttyi >,
                                   tag::output, tag::iter, tag::tty >,
           tk::grm::process_alpha< use< kw::scheme >,
                                   tk::grm::store_inciter_option<
                                     inciter::ctr::Scheme,
                                     tag::discr,
                                     tag::scheme >,
                                   tk::grm::configure_scheme >,
           discroption< use, kw::limiter, inciter::ctr::Limiter, tag::limiter >,
           tk::grm::discrparam< use, kw::cweight, tag::cweight >
         > {};

  //! PDE parameter vector
  template< class keyword, class eq, class param, class... xparams >
  struct pde_parameter_vector :
         tk::grm::parameter_vector< use,
                                    use< keyword >,
                                    tk::grm::Store_back_back,
                                    tk::grm::start_vector,
                                    tk::grm::check_vector,
                                    eq, param, xparams... > {};

   //! Match box parameter
  template< class eq, typename keyword, typename target >
  struct box_parameter :
           pegtl::if_must<
             tk::grm::readkw< typename use<keyword>::pegtl_string >,
             tk::grm::scan<
               pegtl::sor< tk::grm::number,
                           tk::grm::msg< tk::grm::ERROR,
                                         tk::grm::MsgKey::MISSING > >,
               tk::grm::Back_back_store< target,
                 tag::param, eq, tag::ic, tag::box > > > {};

   //! Match box parameter and store deep
  template< class eq, typename keyword, typename target, typename subtarget >
  struct box_deep_parameter :
           pegtl::if_must<
             tk::grm::readkw< typename use<keyword>::pegtl_string >,
             tk::grm::scan<
               pegtl::sor< tk::grm::number,
                           tk::grm::msg< tk::grm::ERROR,
                                         tk::grm::MsgKey::MISSING > >,
               tk::grm::Back_back_deep_store< target, subtarget,
                 tag::param, eq, tag::ic, tag::box > > > {};

   //! Match box parameter vector
  template< class eq, typename keyword, typename target >
  struct box_vector :
         tk::grm::vector< use< keyword >,
                          tk::grm::Back_back_store_back< target,
                            tag::param, eq, tag::ic, tag::box >,
                          use< kw::end > > {};

   //! Match box parameter vector and store deep
  template< class eq, typename keyword, typename target, typename subtarget >
  struct box_deep_vector :
         tk::grm::vector< use< keyword >,
                          tk::grm::Back_back_deep_store_back< target, subtarget,
                            tag::param, eq, tag::ic, tag::box >,
                          use< kw::end > > {};


   //! Match box option
  template< class eq, typename Option, typename keyword, typename target,
            typename subtarget >
  struct box_option :
         tk::grm::process<
           use< keyword >,
           tk::grm::back_back_deep_store_option< target, subtarget, use,
             Option, tag::param, eq, tag::ic, tag::box >,
           pegtl::alpha > {};

   //! Match material option
  template< class eq, typename Option, typename keyword, typename target >
  struct material_option :
         tk::grm::process<
           use< keyword >,
           tk::grm::back_back_store_option< target, use, Option,
             tag::param, eq, tag::material >,
           pegtl::alpha > {};

  //! Match material parameter vector
  template< class eq, typename keyword, typename target >
  struct material_vector :
         tk::grm::vector< use< keyword >,
                          tk::grm::Back_back_store_back< target,
                            tag::param, eq, tag::material >,
                          use< kw::end > > {};

  //! put in PDE parameter for equation matching keyword
  template< typename eq, typename keyword, typename param,
            class kw_type = tk::grm::number >
  struct parameter :
         tk::grm::process< use< keyword >,
                           tk::grm::Store_back< tag::param, eq, param >,
                           kw_type > {};

  //! put in PDE bool parameter for equation matching keyword into vector< int >
  template< typename eq, typename keyword, typename p >
  struct parameter_bool :
         tk::grm::process< use< keyword >,
                           tk::grm::Store_back_bool< tag::param, eq, p >,
                           pegtl::alpha > {};

  //! Boundary conditions block
  template< class keyword, class eq, class param >
  struct bc :
         pegtl::if_must<
           tk::grm::readkw< typename use< keyword >::pegtl_string >,
           tk::grm::block<
             use< kw::end >,
             tk::grm::parameter_vector< use,
                                        use< kw::sideset >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, tag::bc, param > > > {};

  //! Match user-defined function as a discrete list of real numbers
  template< class target, template< class... > class insert, class tag,
    class... tags >
  struct user_fn :
         pegtl::if_must<
           tk::grm::readkw< use< kw::fn >::pegtl_string >,
           tk::grm::block< use< kw::end >,
             tk::grm::scan< tk::grm::number,
               insert< target, tag, tags... > > > > {};

  //! User defined time dependent BC bc_timedep...end block
  template< class eq >
  struct timedep_bc :
         pegtl::if_must<
           tk::grm::readkw< use< kw::bc_timedep >::pegtl_string >,
           tk::grm::start_vector_back< tag::param, eq, tag::bctimedep >,
           tk::grm::block< use< kw::end >,
             user_fn< tag::fn, tk::grm::Back_back_store_back, tag::param, eq,
               tag::bctimedep >,
             pegtl::if_must< tk::grm::vector< use< kw::sideset >,
               tk::grm::Back_back_store_back< tag::sideset, tag::param, eq,
                 tag::bctimedep >,
               use< kw::end > > > > > {};

  //! Stagnation boundary conditions block
  template< class eq, class bc, class kwbc >
  struct bc_spec :
         pegtl::if_must<
           tk::grm::readkw< typename kwbc::pegtl_string >,
           tk::grm::block<
             use< kw::end >,
             tk::grm::parameter_vector< use,
                                        use< kw::radius >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, bc, tag::radius >,
             tk::grm::parameter_vector< use,
                                        use< kw::point >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, bc, tag::point > > > {};

  //! Boundary conditions block
  template< class eq >
  struct sponge :
         pegtl::if_must<
           tk::grm::readkw< typename use< kw::sponge >::pegtl_string >,
           tk::grm::block<
             use< kw::end >,
             tk::grm::parameter_vector< use,
                                        use< kw::velocity >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, tag::sponge, tag::velocity >,
             tk::grm::parameter_vector< use,
                                        use< kw::pressure >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, tag::sponge, tag::pressure >,
             tk::grm::parameter_vector< use,
                                        use< kw::sideset >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, tag::sponge, tag::sideset > > > {};

  //! Farfield boundary conditions block
  template< class keyword, class eq, class param >
  struct farfield_bc :
         pegtl::if_must<
           tk::grm::readkw< typename use< keyword >::pegtl_string >,
           tk::grm::block<
             use< kw::end >,
             parameter< eq, kw::pressure, tag::farfield_pressure >,
             parameter< eq, kw::density, tag::farfield_density >,
             pde_parameter_vector< kw::velocity, eq,
                                   tag::farfield_velocity >,
             tk::grm::parameter_vector< use,
                                        use< kw::sideset >,
                                        tk::grm::Store_back_back,
                                        tk::grm::start_vector,
                                        tk::grm::check_vector,
                                        eq, tag::bc, param > > > {};

  //! edgelist ... end block
  struct edgelist :
         tk::grm::vector< use< kw::amr_edgelist >,
                          tk::grm::Store_back< tag::amr, tag::edge >,
                          use< kw::end >,
                          tk::grm::check_vector< tag::amr, tag::edge > > {};

  //! xminus configuring coordinate-based edge tagging for mesh refinement
  template< typename keyword, typename Tag >
  struct half_world :
         tk::grm::control< use< keyword >, pegtl::digit, tk::grm::Store,
                           tag::amr, Tag > {};

  //! coords ... end block
  struct coords :
           pegtl::if_must<
             tk::grm::readkw< use< kw::amr_coords >::pegtl_string >,
             tk::grm::block< use< kw::end >,
                             half_world< kw::amr_xminus, tag::xminus >,
                             half_world< kw::amr_xplus, tag::xplus >,
                             half_world< kw::amr_yminus, tag::yminus >,
                             half_world< kw::amr_yplus, tag::yplus >,
                             half_world< kw::amr_zminus, tag::zminus >,
                             half_world< kw::amr_zplus, tag::zplus > > > {};

  //! initial conditins box block
  template< class eq >
  struct box :
         pegtl::if_must<
           tk::grm::readkw< use< kw::box >::pegtl_string >,
           tk::grm::start_vector_back< tag::param, eq, tag::ic, tag::box >,
           tk::grm::block< use< kw::end >
             , box_parameter< eq, kw::xmin, tag::xmin >
             , box_parameter< eq, kw::xmax, tag::xmax >
             , box_parameter< eq, kw::ymin, tag::ymin >
             , box_parameter< eq, kw::ymax, tag::ymax >
             , box_parameter< eq, kw::zmin, tag::zmin >
             , box_parameter< eq, kw::zmax, tag::zmax >
             , box_parameter< eq, kw::materialid, tag::materialid >
             , box_parameter< eq, kw::density, tag::density >
             , box_parameter< eq, kw::pressure, tag::pressure >
             , box_parameter< eq, kw::temperature, tag::temperature >
             , box_parameter< eq, kw::energy_content, tag::energy_content >
             , box_parameter< eq, kw::energy, tag::energy >
             , box_parameter< eq, kw::mass, tag::mass >
             , box_vector< eq, kw::velocity, tag::velocity >
             , box_option< eq, ctr::Initiate, kw::initiate, tag::initiate,
                           tag::init >
             , pegtl::if_must<
                 tk::grm::readkw< use< kw::linear >::pegtl_string >,
                 tk::grm::block< use< kw::end >
                   , box_deep_vector< eq, kw::point, tag::initiate, tag::point >
                   , box_deep_parameter< eq, kw::radius, tag::initiate,
                                         tag::radius >
                   , box_deep_parameter< eq, kw::velocity, tag::initiate,
                                         tag::velocity > > >
             > > {};

  //! initial conditions block for compressible flow
  template< class eq >
  struct ic :
         pegtl::if_must<
           tk::grm::readkw< use< kw::ic >::pegtl_string >,
           tk::grm::block< use< kw::end >,
             pegtl::sor<
               pde_parameter_vector< kw::density, eq,
                                     tag::ic, tag::density >,
               pde_parameter_vector< kw::materialid, eq,
                                     tag::ic, tag::materialid >,
               pde_parameter_vector< kw::velocity, eq,
                                     tag::ic, tag::velocity >,
               pde_parameter_vector< kw::pressure, eq,
                                     tag::ic, tag::pressure >,
               pde_parameter_vector< kw::temperature, eq,
                                     tag::ic, tag::temperature >,
               pde_parameter_vector< kw::energy, eq,
                                     tag::ic, tag::energy > >,
               pegtl::seq< box< eq > > > > {};

  //! put in material property for equation matching keyword
  template< typename eq, typename keyword, typename property >
  struct material_property :
         pde_parameter_vector< keyword, eq, property > {};

  //! Material properties block for compressible flow
  template< class eq >
  struct material_properties :
         pegtl::seq<
          pegtl::if_must<
            tk::grm::readkw< use< kw::material >::pegtl_string >,
              tk::grm::start_vector_back< tag::param, eq, tag::material >,
            tk::grm::block< use< kw::end >,
                material_vector< eq, kw::id, tag::id >
              , material_vector< eq, kw::mat_gamma, tag::gamma >
              , material_vector< eq, kw::mat_mu, tag::mu >
              , material_vector< eq, kw::mat_pstiff, tag::pstiff >
              , material_vector< eq, kw::mat_cv, tag::cv >
              , material_vector< eq, kw::mat_k, tag::k >
              , material_option< eq, ctr::Material, kw::eos, tag::eos >
            > > > {};

  //! Mesh ... end block
  template< class eq >
  struct mesh :
         pegtl::if_must<
           tk::grm::readkw< use< kw::mesh >::pegtl_string >,
           tk::grm::block< use< kw::end >,
             tk::grm::filename< use, tag::param, eq, tag::mesh, tag::filename >
           , pde_parameter_vector< kw::location, eq, tag::mesh, tag::location >
           , pde_parameter_vector< kw::orientation, eq,
                                   tag::mesh, tag::orientation >
           , tk::grm::process<
               use< kw::reference >,
               tk::grm::Store_back< tag::param, eq, tag::mesh, tag::reference >,
               pegtl::alpha >
           > > {};

  //! transport equation for scalars
  struct transport :
         pegtl::if_must<
           scan_eq< use< kw::transport >, tag::transport >,
           tk::grm::block< use< kw::end >,
                           tk::grm::policy< use,
                                            use< kw::physics >,
                                            ctr::Physics,
                                            tag::transport,
                                            tag::physics >,
                           tk::grm::policy< use,
                                            use< kw::problem >,
                                            ctr::Problem,
                                            tag::transport,
                                            tag::problem >,
                           tk::grm::depvar< use,
                                            tag::transport,
                                            tag::depvar >,
                           mesh< tag::transport >,
                           tk::grm::component< use< kw::ncomp >,
                                               tag::transport >,
                           pde_parameter_vector< kw::pde_diffusivity,
                                                 tag::transport,
                                                 tag::diffusivity >,
                           pde_parameter_vector< kw::pde_lambda,
                                                 tag::transport,
                                                 tag::lambda >,
                           pde_parameter_vector< kw::pde_u0,
                                                 tag::transport,
                                                 tag::u0 >,
                           bc< kw::bc_dirichlet, tag::transport, tag::bcdir >,
                           bc< kw::bc_sym, tag::transport, tag::bcsym >,
                           bc< kw::bc_inlet, tag::transport, tag::bcinlet >,
                           bc< kw::bc_outlet, tag::transport, tag::bcoutlet >,
                           bc< kw::bc_extrapolate, tag::transport,
                               tag::bcextrapolate >,
                           parameter< tag::transport,
                                      kw::intsharp_param,
                                      tag::intsharp_param >,
                           parameter< tag::transport,
                                      kw::intsharp,
                                      tag::intsharp > >,
           check_errors< tag::transport, tk::grm::check_transport > > {};

  //! compressible flow
  struct compflow :
         pegtl::if_must<
           scan_eq< use< kw::compflow >, tag::compflow >,
           tk::grm::start_vector< tag::param, tag::compflow, tag::ic, tag::box >,
           tk::grm::start_vector< tag::param, tag::compflow, tag::material >,
           tk::grm::start_vector< tag::param, tag::compflow, tag::bctimedep >,
           tk::grm::block< use< kw::end >,
                           tk::grm::policy< use,
                                            use< kw::physics >,
                                            ctr::Physics,
                                            tag::compflow,
                                            tag::physics >,
                           tk::grm::policy< use,
                                            use< kw::problem >,
                                            ctr::Problem,
                                            tag::compflow,
                                            tag::problem >,
                           tk::grm::depvar< use,
                                            tag::compflow,
                                            tag::depvar >,
                           mesh< tag::compflow >,
                           tk::grm::process<
                             use< kw::flux >,
                               tk::grm::store_back_option< use,
                                                           ctr::Flux,
                                                           tag::param,
                                                           tag::compflow,
                                                           tag::flux >,
                             pegtl::alpha >,
                           ic< tag::compflow >,
                           tk::grm::lua< use, tag::param, tag::compflow >,
                           material_properties< tag::compflow >,
                           pde_parameter_vector< kw::sysfctvar,
                                                 tag::compflow,
                                                 tag::sysfctvar >,
                           parameter_bool< tag::compflow,
                                           kw::sysfct,
                                           tag::sysfct >,
                           parameter< tag::compflow, kw::npar,
                                      tag::npar, pegtl::digit >,
                           parameter< tag::compflow, kw::pde_alpha,
                                      tag::alpha >,
                           parameter< tag::compflow, kw::pde_p0,
                                      tag::p0 >,
                           parameter< tag::compflow, kw::pde_betax,
                                      tag::betax >,
                           parameter< tag::compflow, kw::pde_betay,
                                      tag::betay >,
                           parameter< tag::compflow, kw::pde_betaz,
                                      tag::betaz >,
                           parameter< tag::compflow, kw::pde_beta,
                                      tag::beta >,
                           parameter< tag::compflow, kw::pde_r0,
                                      tag::r0 >,
                           parameter< tag::compflow, kw::pde_ce,
                                      tag::ce >,
                           parameter< tag::compflow, kw::pde_kappa,
                                      tag::kappa >,
                           bc< kw::bc_dirichlet, tag::compflow, tag::bcdir >,
                           bc< kw::bc_sym, tag::compflow, tag::bcsym >,
                           bc_spec< tag::compflow, tag::stag, kw::bc_stag >,
                           bc_spec< tag::compflow, tag::skip, kw::bc_skip >,
                           bc< kw::bc_inlet, tag::compflow, tag::bcinlet >,
                           sponge< tag::compflow >,
                           farfield_bc< kw::bc_farfield,
                                        tag::compflow,
                                        tag::bcfarfield >,
                           bc< kw::bc_extrapolate, tag::compflow,
                               tag::bcextrapolate >,
                           timedep_bc< tag::compflow >
                           >,
           check_errors< tag::compflow, tk::grm::check_compflow > > {};

  //! compressible multi-material flow
  struct multimat :
         pegtl::if_must<
           scan_eq< use< kw::multimat >, tag::multimat >,
           tk::grm::start_vector< tag::param, tag::multimat, tag::ic, tag::box >,
           tk::grm::start_vector< tag::param, tag::multimat, tag::material >,
           tk::grm::block< use< kw::end >,
                           tk::grm::policy< use,
                                            use< kw::physics >,
                                            ctr::Physics,
                                            tag::multimat,
                                            tag::physics >,
                           tk::grm::policy< use,
                                            use< kw::problem >,
                                            ctr::Problem,
                                            tag::multimat,
                                            tag::problem >,
                           tk::grm::depvar< use,
                                            tag::multimat,
                                            tag::depvar >,
                           mesh< tag::multimat >,
                           parameter< tag::multimat,
                                      kw::nmat,
                                      tag::nmat >,
                           tk::grm::process<
                             use< kw::flux >,
                               tk::grm::store_back_option< use,
                                                           ctr::Flux,
                                                           tag::param,
                                                           tag::multimat,
                                                           tag::flux >,
                             pegtl::alpha >,
                           ic< tag::multimat >,
                           material_properties< tag::multimat >,
                           parameter< tag::multimat,
                                      kw::pde_alpha,
                                      tag::alpha >,
                           parameter< tag::multimat,
                                      kw::pde_p0,
                                      tag::p0 >,
                           parameter< tag::multimat,
                                      kw::pde_beta,
                                      tag::beta >,
                           bc< kw::bc_dirichlet,
                               tag::multimat,
                               tag::bcdir >,
                           bc< kw::bc_sym,
                               tag::multimat,
                               tag::bcsym >,
                           bc< kw::bc_inlet,
                               tag::multimat,
                               tag::bcinlet >,
                           bc< kw::bc_outlet,
                               tag::multimat,
                               tag::bcoutlet >,
                           bc< kw::bc_extrapolate,
                               tag::multimat,
                               tag::bcextrapolate >,
                           parameter< tag::multimat,
                                      kw::prelax_timescale,
                                      tag::prelax_timescale >,
                           parameter< tag::multimat,
                                      kw::prelax,
                                      tag::prelax >,
                           parameter< tag::multimat,
                                      kw::intsharp_param,
                                      tag::intsharp_param >,
                           parameter< tag::multimat,
                                      kw::intsharp,
                                      tag::intsharp > >,
           check_errors< tag::multimat, tk::grm::check_multimat > > {};

  //! partitioning ... end block
  struct partitioning :
         pegtl::if_must<
           tk::grm::readkw< use< kw::partitioning >::pegtl_string >,
           tk::grm::block< use< kw::end >,
                           tk::grm::process<
                             use< kw::algorithm >,
                             tk::grm::store_inciter_option<
                               tk::ctr::PartitioningAlgorithm,
                               tag::selected,
                               tag::partitioner >,
                             pegtl::alpha > > > {};

  //! equation types
  struct equations :
         pegtl::sor< transport, compflow, multimat > {};

  //! refinement variable(s) (refvar) ... end block
  struct refvars :
         pegtl::if_must<
           tk::grm::vector< use< kw::amr_refvar >,
                            tk::grm::match_depvar<
                              tk::grm::Store_back< tag::amr, tag::refvar > >,
                            use< kw::end >,
                            tk::grm::check_vector< tag::amr, tag::refvar >,
                            tk::grm::fieldvar< pegtl::alpha > >,
           tk::grm::compute_refvar_idx > {};

  //! adaptive mesh refinement (AMR) amr...end block
  struct amr :
         pegtl::if_must<
           tk::grm::readkw< use< kw::amr >::pegtl_string >,
           // enable AMR if amr...end block encountered
           tk::grm::enable< tag::amr >,
           tk::grm::block< use< kw::end >,
                           refvars,
                           edgelist,
                           coords,
                           tk::grm::process<
                             use< kw::amr_initial >,
                             tk::grm::store_back_option< use,
                                                         ctr::AMRInitial,
                                                         tag::amr,
                                                         tag::init >,
                             pegtl::alpha >,
                           tk::grm::process<
                             use< kw::amr_error >,
                             tk::grm::store_inciter_option<
                               ctr::AMRError,
                               tag::amr, tag::error >,
                             pegtl::alpha >,
                           tk::grm::control< use< kw::amr_tolref >,
                                             pegtl::digit,
                                             tk::grm::Store,
                                             tag::amr,
                                             tag::tolref >,
                           tk::grm::control< use< kw::amr_tolderef >,
                                             pegtl::digit,
                                             tk::grm::Store,
                                             tag::amr,
                                             tag::tolderef >,
                           tk::grm::process< use< kw::amr_t0ref >,
                             tk::grm::Store< tag::amr, tag::t0ref >,
                             pegtl::alpha >,
                           tk::grm::process< use< kw::amr_dtref_uniform >,
                             tk::grm::Store< tag::amr, tag::dtref_uniform >,
                             pegtl::alpha >,
                           tk::grm::process< use< kw::amr_dtref >,
                             tk::grm::Store< tag::amr, tag::dtref >,
                             pegtl::alpha >,
                           tk::grm::process< use< kw::amr_dtfreq >,
                             tk::grm::Store< tag::amr, tag::dtfreq >,
                             pegtl::digit > >,
           tk::grm::check_amr_errors > {};


  //! Arbitrary-Lagrangian-Eulerian (ALE) move...end block
  struct moving_sides :
         pegtl::if_must<
           tk::grm::readkw< use< kw::move >::pegtl_string >,
           tk::grm::start_vector< tag::ale, tag::move >,
           tk::grm::block< use< kw::end >,
             tk::grm::process<
                use< kw::fntype >,
                tk::grm::back_store_option< tag::fntype,
                                            use,
                                            tk::ctr::UserTable,
                                            tag::ale, tag::move >,
                pegtl::alpha >,
             user_fn< tag::fn, tk::grm::Back_store_back, tag::ale, tag::move >,
             pegtl::if_must< tk::grm::vector< use< kw::sideset >,
               tk::grm::Back_store_back< tag::sideset, tag::ale, tag::move >,
               use< kw::end > > > > > {};

  //! Arbitrary-Lagrangian-Eulerian (ALE) ale...end block
  struct ale :
         pegtl::if_must<
           tk::grm::readkw< use< kw::ale >::pegtl_string >,
           // enable ALE if ale ...end block encountered
           tk::grm::enable< tag::ale >,
           tk::grm::block< use< kw::end >,
              tk::grm::control< use< kw::dvcfl >,
                                pegtl::digit,
                                tk::grm::Store,
                                tag::ale, tag::dvcfl >,
              tk::grm::control< use< kw::vortmult >,
                                pegtl::digit,
                                tk::grm::Store,
                                tag::ale, tag::vortmult >,
              tk::grm::control< use< kw::meshvel_maxit >,
                                pegtl::digit,
                                tk::grm::Store,
                                tag::ale, tag::maxit >,
              tk::grm::control< use< kw::meshvel_tolerance >,
                                pegtl::digit,
                                tk::grm::Store,
                                tag::ale, tag::tolerance >,
              moving_sides,
              tk::grm::process<
                use< kw::meshvelocity >,
                tk::grm::store_inciter_option< ctr::MeshVelocity,
                                               tag::ale, tag::meshvelocity >,
                pegtl::alpha >,
              tk::grm::process<
                use< kw::smoother >,
                tk::grm::store_inciter_option< ctr::MeshVelocitySmoother,
                                               tag::ale, tag::smoother >,
                pegtl::alpha >,
              pegtl::if_must< tk::grm::dimensions< use< kw::mesh_motion >,
                              tk::grm::Store_back< tag::ale, tag::mesh_motion >,
                              use< kw::end > > >,
              pegtl::if_must< tk::grm::vector< use< kw::meshforce >,
                              tk::grm::Store_back< tag::ale, tag::meshforce >,
                              use< kw::end > > >,
              pegtl::if_must<
                tk::grm::readkw< use< kw::bc_dirichlet >::pegtl_string >,
                tk::grm::block< use< kw::end >,
                  pegtl::if_must< tk::grm::vector< use< kw::sideset >,
                                  tk::grm::Store_back< tag::ale, tag::bcdir >,
                                  use< kw::end > > > > >,
              pegtl::if_must<
                tk::grm::readkw< use< kw::bc_sym >::pegtl_string >,
                tk::grm::block< use< kw::end >,
                  pegtl::if_must< tk::grm::vector< use< kw::sideset >,
                                  tk::grm::Store_back< tag::ale, tag::bcsym >,
                                  use< kw::end > > > > > >,
              tk::grm::check_ale > {};

  //! \brief Match a depvar, defined upstream of control file, coupling a
  //!   solver and store
  template< template< class > class action >
  struct coupled_solver :
         tk::grm::scan_until<
            pegtl::lower,
            tk::grm::match_depvar< tk::grm::push_transfer< action > > > {};

  //! Couple ... end block (used to configure solver coupling)
  struct couple :
         pegtl::if_must<
           tk::grm::readkw< use< kw::couple >::pegtl_string >,
           tk::grm::block< use< kw::end >,
             pegtl::seq<
               coupled_solver< tk::grm::store_transfer_src >,
               pegtl::one<'>'>,
               coupled_solver< tk::grm::store_transfer_dst > > > > {};

  //! p-adaptive refinement (pref) ...end block
  struct pref :
         pegtl::if_must<
           tk::grm::readkw< use< kw::pref >::pegtl_string >,
           tk::grm::block< use< kw::end >,
                           tk::grm::control< use< kw::pref_tolref >,
                                             pegtl::digit,
                                             tk::grm::Store,
                                             tag::pref,
                                             tag::tolref  >,
                           tk::grm::control< use< kw::pref_ndofmax >,
                                             pegtl::digit,
                                             tk::grm::Store,
                                             tag::pref,
                                             tag::ndofmax >,
                           tk::grm::process<
                             use< kw::pref_indicator >,
                             tk::grm::store_inciter_option<
                               ctr::PrefIndicator,
                               tag::pref, tag::indicator >,
                             pegtl::alpha >
                         >,
           tk::grm::check_pref_errors > {};

  //! Match output variable alias
  struct outvar_alias :
         tk::grm::quoted< tk::grm::set_outvar_alias > {};

  //! Match an output variable in a human readable form: var must be a keyword
  template< class var >
  struct outvar_human :
         tk::grm::exact_scan< use< var >, tk::grm::push_humanvar > {};

  //! Match an output variable based on depvar defined upstream of input file
  struct outvar_depvar :
           tk::grm::scan< tk::grm::fieldvar< pegtl::upper >,
             tk::grm::match_outvar, tk::grm::check_outvar > {};

  //! Parse a centering token and if matches, set centering in parser's state
  struct outvar_centering :
         pegtl::sor<
           tk::grm::exact_scan< use< kw::node >, tk::grm::set_centering >,
           tk::grm::exact_scan< use< kw::elem >, tk::grm::set_centering > > {};

  //! outvar ... end block
  struct outvar_block :
         pegtl::if_must<
           tk::grm::readkw< use< kw::outvar >::pegtl_string >,
           tk::grm::block<
             use< kw::end >
           , outvar_centering
           , outvar_depvar
           , outvar_alias
           , outvar_human< kw::outvar_density >
           , outvar_human< kw::outvar_xmomentum >
           , outvar_human< kw::outvar_ymomentum >
           , outvar_human< kw::outvar_zmomentum >
           , outvar_human< kw::outvar_specific_total_energy >
           , outvar_human< kw::outvar_volumetric_total_energy >
           , outvar_human< kw::outvar_xvelocity >
           , outvar_human< kw::outvar_yvelocity >
           , outvar_human< kw::outvar_zvelocity >
           , outvar_human< kw::outvar_pressure >
           , outvar_human< kw::outvar_material_indicator >
           , outvar_human< kw::outvar_analytic >
           > > {};

  //! field_output ... end block
  struct field_output :
         pegtl::if_must<
           tk::grm::readkw< use< kw::field_output >::pegtl_string >,
           tk::grm::block<
             use< kw::end >,
             outvar_block,
             tk::grm::process< use< kw::filetype >,
                               tk::grm::store_inciter_option<
                                 tk::ctr::FieldFile,
                                 tag::selected,
                                 tag::filetype >,
                               pegtl::alpha >,
             tk::grm::interval_iter< use< kw::interval_iter >,
                                     tag::output, tag::iter, tag::field >,
             tk::grm::interval_time< use< kw::interval_time >,
                                     tag::output, tag::time, tag::field >,
             tk::grm::time_range< use, kw::time_range,
                                  tag::output, tag::range, tag::field >,
             tk::grm::process<
               use< kw::refined >,
               tk::grm::Store< tag::cmd, tag::io, tag::refined >,
               pegtl::alpha >,
             pegtl::if_must<
               tk::grm::vector<
                 use< kw::sideset >,
                 tk::grm::Store_back< tag::cmd, tag::io, tag::surface >,
                 use< kw::end > > > > > {};

  //! history_output ... end block
  struct history_output :
         pegtl::if_must<
           tk::grm::readkw< use< kw::history_output >::pegtl_string >,
           tk::grm::block<
             use< kw::end >,
             outvar_block,
             tk::grm::interval_iter< use< kw::interval_iter >,
               tag::output, tag::iter, tag::history >,
             tk::grm::interval_time< use< kw::interval_time >,
               tag::output, tag::time, tag::history >,
             tk::grm::time_range< use, kw::time_range,
                                  tag::output, tag::range, tag::history >,
             tk::grm::precision< use, tag::history >,
             tk::grm::process<
               use< kw::txt_float_format >,
               tk::grm::store_inciter_option< tk::ctr::TxtFloatFormat,
                                              tag::flformat,
                                              tag::history >,
               pegtl::alpha >,
             pegtl::if_must<
               tk::grm::readkw< use< kw::point >::pegtl_string >,
               tk::grm::act< pegtl::identifier, tk::grm::match_pointname >,
               pegtl::seq<
                 tk::grm::start_vector< tag::history, tag::point >,
                 tk::grm::block<
                   use< kw::end >,
                   tk::grm::scan< tk::grm::number,
                     tk::grm::Store_back_back< tag::history, tag::point > > >
               > > > > {};

  //! 'inciter' block
  struct inciter :
         pegtl::if_must<
           tk::grm::readkw< use< kw::inciter >::pegtl_string >,
           pegtl::sor<
             pegtl::seq< tk::grm::block<
                           use< kw::end >,
                           discretization,
                           equations,
                           amr,
                           ale,
                           pref,
                           partitioning,
                           couple,
                           field_output,
                           history_output,
                           tk::grm::diagnostics<
                             use,
                             tk::grm::store_inciter_option > >,
                         tk::grm::check_inciter >,
            tk::grm::msg< tk::grm::MsgType::ERROR,
                          tk::grm::MsgKey::UNFINISHED > > > {};

  //! \brief All keywords
  struct keywords :
         pegtl::sor< tk::grm::title< use >, inciter > {};

  //! \brief Grammar entry point: parse keywords and ignores until eof
  struct read_file :
         tk::grm::read_file< keywords, tk::grm::ignore > {};

} // deck::
} // inciter::

#endif // InciterInputDeckGrammar_h