Source

pflotran-dev / src / pflotran / realization.F90

   1
   2
   3
   4
   5
   6
   7
   8
   9
  10
  11
  12
  13
  14
  15
  16
  17
  18
  19
  20
  21
  22
  23
  24
  25
  26
  27
  28
  29
  30
  31
  32
  33
  34
  35
  36
  37
  38
  39
  40
  41
  42
  43
  44
  45
  46
  47
  48
  49
  50
  51
  52
  53
  54
  55
  56
  57
  58
  59
  60
  61
  62
  63
  64
  65
  66
  67
  68
  69
  70
  71
  72
  73
  74
  75
  76
  77
  78
  79
  80
  81
  82
  83
  84
  85
  86
  87
  88
  89
  90
  91
  92
  93
  94
  95
  96
  97
  98
  99
 100
 101
 102
 103
 104
 105
 106
 107
 108
 109
 110
 111
 112
 113
 114
 115
 116
 117
 118
 119
 120
 121
 122
 123
 124
 125
 126
 127
 128
 129
 130
 131
 132
 133
 134
 135
 136
 137
 138
 139
 140
 141
 142
 143
 144
 145
 146
 147
 148
 149
 150
 151
 152
 153
 154
 155
 156
 157
 158
 159
 160
 161
 162
 163
 164
 165
 166
 167
 168
 169
 170
 171
 172
 173
 174
 175
 176
 177
 178
 179
 180
 181
 182
 183
 184
 185
 186
 187
 188
 189
 190
 191
 192
 193
 194
 195
 196
 197
 198
 199
 200
 201
 202
 203
 204
 205
 206
 207
 208
 209
 210
 211
 212
 213
 214
 215
 216
 217
 218
 219
 220
 221
 222
 223
 224
 225
 226
 227
 228
 229
 230
 231
 232
 233
 234
 235
 236
 237
 238
 239
 240
 241
 242
 243
 244
 245
 246
 247
 248
 249
 250
 251
 252
 253
 254
 255
 256
 257
 258
 259
 260
 261
 262
 263
 264
 265
 266
 267
 268
 269
 270
 271
 272
 273
 274
 275
 276
 277
 278
 279
 280
 281
 282
 283
 284
 285
 286
 287
 288
 289
 290
 291
 292
 293
 294
 295
 296
 297
 298
 299
 300
 301
 302
 303
 304
 305
 306
 307
 308
 309
 310
 311
 312
 313
 314
 315
 316
 317
 318
 319
 320
 321
 322
 323
 324
 325
 326
 327
 328
 329
 330
 331
 332
 333
 334
 335
 336
 337
 338
 339
 340
 341
 342
 343
 344
 345
 346
 347
 348
 349
 350
 351
 352
 353
 354
 355
 356
 357
 358
 359
 360
 361
 362
 363
 364
 365
 366
 367
 368
 369
 370
 371
 372
 373
 374
 375
 376
 377
 378
 379
 380
 381
 382
 383
 384
 385
 386
 387
 388
 389
 390
 391
 392
 393
 394
 395
 396
 397
 398
 399
 400
 401
 402
 403
 404
 405
 406
 407
 408
 409
 410
 411
 412
 413
 414
 415
 416
 417
 418
 419
 420
 421
 422
 423
 424
 425
 426
 427
 428
 429
 430
 431
 432
 433
 434
 435
 436
 437
 438
 439
 440
 441
 442
 443
 444
 445
 446
 447
 448
 449
 450
 451
 452
 453
 454
 455
 456
 457
 458
 459
 460
 461
 462
 463
 464
 465
 466
 467
 468
 469
 470
 471
 472
 473
 474
 475
 476
 477
 478
 479
 480
 481
 482
 483
 484
 485
 486
 487
 488
 489
 490
 491
 492
 493
 494
 495
 496
 497
 498
 499
 500
 501
 502
 503
 504
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 546
 547
 548
 549
 550
 551
 552
 553
 554
 555
 556
 557
 558
 559
 560
 561
 562
 563
 564
 565
 566
 567
 568
 569
 570
 571
 572
 573
 574
 575
 576
 577
 578
 579
 580
 581
 582
 583
 584
 585
 586
 587
 588
 589
 590
 591
 592
 593
 594
 595
 596
 597
 598
 599
 600
 601
 602
 603
 604
 605
 606
 607
 608
 609
 610
 611
 612
 613
 614
 615
 616
 617
 618
 619
 620
 621
 622
 623
 624
 625
 626
 627
 628
 629
 630
 631
 632
 633
 634
 635
 636
 637
 638
 639
 640
 641
 642
 643
 644
 645
 646
 647
 648
 649
 650
 651
 652
 653
 654
 655
 656
 657
 658
 659
 660
 661
 662
 663
 664
 665
 666
 667
 668
 669
 670
 671
 672
 673
 674
 675
 676
 677
 678
 679
 680
 681
 682
 683
 684
 685
 686
 687
 688
 689
 690
 691
 692
 693
 694
 695
 696
 697
 698
 699
 700
 701
 702
 703
 704
 705
 706
 707
 708
 709
 710
 711
 712
 713
 714
 715
 716
 717
 718
 719
 720
 721
 722
 723
 724
 725
 726
 727
 728
 729
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
module Realization_class
  
  use Realization_Base_class

  use Option_module
  use Output_Aux_module
  use Input_module
  use Region_module
  use Condition_module
  use Constraint_module
  use Material_module
  use Saturation_Function_module
  use Dataset_Aux_module
  use Fluid_module
  use Discretization_module
  use Field_module
  use Debug_module
  use Uniform_Velocity_module
  use Waypoint_module
  use Output_Aux_module
  
  use Reaction_Aux_module
  
  use Level_module
  use Patch_module
  
  implicit none

private

#include "definitions.h"
  type, public, extends(realization_base_type) :: realization_type

    type(region_list_type), pointer :: regions
    type(condition_list_type), pointer :: flow_conditions
    type(tran_condition_list_type), pointer :: transport_conditions
    type(tran_constraint_list_type), pointer :: transport_constraints
    
    type(tran_constraint_type), pointer :: sec_transport_constraint
    type(material_property_type), pointer :: material_properties
    type(material_property_ptr_type), pointer :: material_property_array(:)
    type(fluid_property_type), pointer :: fluid_properties
    type(fluid_property_type), pointer :: fluid_property_array(:)
    type(saturation_function_type), pointer :: saturation_functions
    type(dataset_type), pointer :: datasets
    type(saturation_function_ptr_type), pointer :: saturation_function_array(:)
    
    type(uniform_velocity_dataset_type), pointer :: uniform_velocity_dataset
    character(len=MAXSTRINGLENGTH) :: nonuniform_velocity_filename
    
    type(waypoint_list_type), pointer :: waypoints
    
  end type realization_type

  interface RealizationCreate
    module procedure RealizationCreate1
    module procedure RealizationCreate2
  end interface
  
  public :: RealizationCreate, &
            RealizationDestroy, &
            RealizationProcessCouplers, &
            RealizationInitAllCouplerAuxVars, &
            RealizationProcessConditions, &
            RealizationUpdate, &
            RealizationAddWaypointsToList, &
            RealizationCreateDiscretization, &
            RealizationLocalizeRegions, &
            RealizationAddCoupler, &
            RealizationAddStrata, &
            RealizationAddObservation, &
            RealizUpdateUniformVelocity, &
            RealizationRevertFlowParameters, &
!            RealizationGetDataset, &
!            RealizGetDatasetValueAtCell, &
!            RealizationSetDataset, &
            RealizationPrintCouplers, &
            RealizationInitConstraints, &
            RealProcessMatPropAndSatFunc, &
            RealProcessFluidProperties, &
            RealizationUpdateProperties, &
            RealizationCountCells, &
            RealizationPrintGridStatistics, &
            RealizationSetUpBC4Faces, &
            RealizatonPassPtrsToPatches, &
            RealLocalToLocalWithArray, &
            RealizationCalculateCFL1Timestep, &
            RealizationNonInitializedData

  !TODO(intel)
  ! public from Realization_Base_class
  !public :: RealizationGetDataset

contains
  
! ************************************************************************** !
!
! RealizationCreate1: Allocates and initializes a new Realization object
! author: Glenn Hammond
! date: 10/25/07
!
! ************************************************************************** !
function RealizationCreate1()

  implicit none
  
  type(realization_type), pointer :: RealizationCreate1
  
  type(realization_type), pointer :: realization
  type(option_type), pointer :: option
  
  nullify(option)
  RealizationCreate1 => RealizationCreate2(option)
  
end function RealizationCreate1  

! ************************************************************************** !
!
! RealizationCreate2: Allocates and initializes a new Realization object
! author: Glenn Hammond
! date: 10/25/07
!
! ************************************************************************** !
function RealizationCreate2(option)

  implicit none
  
  type(option_type), pointer :: option
  
  type(realization_type), pointer :: RealizationCreate2
  
  type(realization_type), pointer :: realization
  
  allocate(realization)
  call RealizationBaseInit(realization,option)

  allocate(realization%regions)
  call RegionInitList(realization%regions)

  allocate(realization%flow_conditions)
  call FlowConditionInitList(realization%flow_conditions)
  allocate(realization%transport_conditions)
  call TranConditionInitList(realization%transport_conditions)
  allocate(realization%transport_constraints)
  call TranConstraintInitList(realization%transport_constraints)

  nullify(realization%material_properties)
  nullify(realization%material_property_array)
  nullify(realization%fluid_properties)
  nullify(realization%fluid_property_array)
  nullify(realization%saturation_functions)
  nullify(realization%saturation_function_array)
  nullify(realization%datasets)
  nullify(realization%uniform_velocity_dataset)
  nullify(realization%sec_transport_constraint)
  realization%nonuniform_velocity_filename = ''

  nullify(realization%waypoints)

  RealizationCreate2 => realization
  
end function RealizationCreate2 

! ************************************************************************** !
!
! RealizationCreateDiscretization: Creates grid
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationCreateDiscretization(realization)

  use Grid_module
  use Unstructured_Grid_Aux_module
  use Unstructured_Grid_module, only : UGridEnsureRightHandRule
  use Structured_Grid_module, only : StructGridCreateTVDGhosts
  use MFD_module
  use Coupler_module
  use Discretization_module
  use Unstructured_Cell_module
  
  implicit none
  
#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"

  type(realization_type) :: realization
  
  type(discretization_type), pointer :: discretization
  type(grid_type), pointer :: grid
  type(field_type), pointer :: field
  type(option_type), pointer :: option
  type(coupler_type), pointer :: boundary_condition
  PetscErrorCode :: ierr
  PetscInt, allocatable :: int_tmp(:)
  PetscInt :: test,j, num_LP_dof
  PetscOffset :: i_da
  PetscReal, pointer :: real_tmp(:)
  type(dm_ptr_type), pointer :: dm_ptr
  Vec :: is_bnd_vec
  PetscInt :: ivar

  option => realization%option
  field => realization%field
 
  discretization => realization%discretization
  
  call DiscretizationCreateDMs(discretization,option)

  ! 1 degree of freedom, global
  call DiscretizationCreateVector(discretization,ONEDOF,field%porosity0, &
                                  GLOBAL,option)
 
  call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                     field%tortuosity0)
  call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                     field%volume)

  call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                     field%work)
  
  ! 1 degree of freedom, local
  call DiscretizationCreateVector(discretization,ONEDOF,field%porosity_loc, &
                                  LOCAL,option)
  call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                     field%tortuosity_loc)

  call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                     field%work_loc)
  
  if (option%nflowdof > 0) then

    ! 1-dof global  
    call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                       field%perm0_xx)
    call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                       field%perm0_yy)
    call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                       field%perm0_zz)
    if (discretization%itype == STRUCTURED_GRID_MIMETIC.or. &
        discretization%itype == UNSTRUCTURED_GRID_MIMETIC) then
      call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                         field%perm0_xz)
      call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                         field%perm0_xy)
      call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                         field%perm0_yz)
    endif
    call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                       field%perm_pow)

    ! 1-dof local
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%ithrm_loc)
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%icap_loc)
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%iphas_loc)
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%iphas_old_loc)
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%perm_xx_loc)
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%perm_yy_loc)
    call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                       field%perm_zz_loc)
    if (discretization%itype == STRUCTURED_GRID_MIMETIC.or. &
        discretization%itype == UNSTRUCTURED_GRID_MIMETIC) then
      call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                         field%perm_xz_loc)
      call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                         field%perm_xy_loc)
      call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                         field%perm_yz_loc)
    endif

    ! ndof degrees of freedom, global
    call DiscretizationCreateVector(discretization,NFLOWDOF,field%flow_xx, &
                                    GLOBAL,option)
    call DiscretizationDuplicateVector(discretization,field%flow_xx, &
                                       field%flow_yy)
    call DiscretizationDuplicateVector(discretization,field%flow_xx, &
                                       field%flow_dxx)
    call DiscretizationDuplicateVector(discretization,field%flow_xx, &
                                       field%flow_r)
    call DiscretizationDuplicateVector(discretization,field%flow_xx, &
                                       field%flow_accum)

    ! ndof degrees of freedom, local
    call DiscretizationCreateVector(discretization,NFLOWDOF,field%flow_xx_loc, &
                                    LOCAL,option)
  endif

  if (option%ntrandof > 0) then
    if (option%reactive_transport_coupling == GLOBAL_IMPLICIT) then
      ! ndof degrees of freedom, global
      call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx, &
                                      GLOBAL,option)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_yy)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_dxx)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_r)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_accum)

      ! ndof degrees of freedom, local
      call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx_loc, &
                                      LOCAL,option)
                                      
      if (realization%reaction%use_log_formulation) then
        call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                           field%tran_log_xx)
        call DiscretizationDuplicateVector(discretization,field%tran_xx_loc, &
                                           field%tran_work_loc)
      endif
 
    else ! operator splitting
      ! ndof degrees of freedom, global
      ! create the 1 dof vector for solving the individual linear systems
      call DiscretizationCreateVector(discretization,ONEDOF,field%tran_rhs_coef, &
                                      GLOBAL,option)
      ! create the ntran dof vector for storage of the solution
      call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx, &
                                      GLOBAL,option)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_yy)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_dxx)
      call DiscretizationDuplicateVector(discretization,field%tran_xx, &
                                         field%tran_rhs)

      ! ndof degrees of freedom, local
      ! again, just for storage of the current colution
      call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx_loc, &
                                      LOCAL,option)

    endif
    
  endif

  select case(discretization%itype)
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
      grid => discretization%grid
      ! set up nG2L, nL2G, etc.
      call GridMapIndices(grid, discretization%dm_1dof%dm, &
                          discretization%stencil_type,&
                          discretization%lsm_flux_method, &
                          option)
      if (option%itranmode == EXPLICIT_ADVECTION) then
        call StructGridCreateTVDGhosts(grid%structured_grid, &
                                       realization%reaction%naqcomp, &
                                       field%tran_xx, &
                                       discretization%dm_1dof%dm, &
                                       field%tvd_ghosts, &
                                       discretization%tvd_ghost_scatter, &
                                       option)
      endif
      call GridComputeSpacing(grid,option)
      call GridComputeCoordinates(grid,discretization%origin,option)
      call GridComputeVolumes(grid,field%volume,option)
      ! set up internal connectivity, distance, etc.
      call GridComputeInternalConnect(grid,option)
      if (discretization%itype == STRUCTURED_GRID_MIMETIC) then
          call GridComputeCell2FaceConnectivity(grid, discretization%MFD, option)
      end if
      if (discretization%lsm_flux_method) then
        call DiscretizationDuplicateVector(discretization,field%porosity_loc, &
                                          is_bnd_vec)
        call GridComputeNeighbors(grid,field%work_loc,option)
        call DiscretizationLocalToLocal(discretization,field%work_loc,is_bnd_vec,ONEDOF)
        call GridSaveBoundaryCellInfo(discretization%grid,is_bnd_vec,option)
        call VecDestroy(is_bnd_vec,ierr)
      endif
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      grid => discretization%grid
      ! set up nG2L, NL2G, etc.
      call UGridMapIndices(grid%unstructured_grid, &
                           discretization%dm_1dof%ugdm, &
                           grid%nG2L,grid%nL2G,grid%nG2A,grid%nG2P,option)
      call GridComputeCoordinates(grid,discretization%origin,option, & 
                                    discretization%dm_1dof%ugdm) 
      if (grid%itype == IMPLICIT_UNSTRUCTURED_GRID) then
        call UGridEnsureRightHandRule(grid%unstructured_grid,grid%x, &
                                      grid%y,grid%z,grid%nG2A,grid%nL2G, &
                                      option)
      endif
      ! set up internal connectivity, distance, etc.
      call GridComputeInternalConnect(grid,option, &
                                      discretization%dm_1dof%ugdm) 
      call GridComputeVolumes(grid,field%volume,option)
#ifdef MFD_UGRID
      call GridComputeCell2FaceConnectivity(discretization%grid,discretization%MFD,option)
#endif
  end select 
 
  ! Vectors with face degrees of freedom
#ifdef DASVYAT
  if (discretization%itype == STRUCTURED_GRID_MIMETIC .or. &
        discretization%itype == UNSTRUCTURED_GRID_MIMETIC) then

    if (option%nflowdof > 0) then
      num_LP_dof = (grid%nlmax_faces + grid%nlmax)*option%nflowdof
      call VecCreateMPI(option%mycomm, num_LP_dof, &
                  PETSC_DETERMINE,field%flow_xx_faces,ierr)
      call VecSetBlockSize(field%flow_xx_faces,option%nflowdof,ierr)

      call DiscretizationDuplicateVector(discretization, field%flow_xx_faces, &
                                        field%flow_r_faces)
      call DiscretizationDuplicateVector(discretization, field%flow_xx_faces, &
                                        field%flow_dxx_faces)
      call DiscretizationDuplicateVector(discretization, field%flow_xx_faces, &
                                        field%flow_yy_faces)

      call VecCreateSeq(PETSC_COMM_SELF, (grid%ngmax_faces + grid%ngmax)*option%nflowdof, &
                                              field%flow_xx_loc_faces, ierr)
      call VecSetBlockSize(field%flow_xx_loc_faces,option%nflowdof,ierr)

      call DiscretizationDuplicateVector(discretization, field%flow_xx_loc_faces, &
                                          field%flow_r_loc_faces)
      call DiscretizationDuplicateVector(discretization, field%flow_xx_loc_faces, &
                                          field%flow_bc_loc_faces)
      call DiscretizationDuplicateVector(discretization, field%flow_xx_loc_faces, &
                                          field%work_loc_faces)
     endif

    call RealizationCreatenG2LP(realization)

    dm_ptr => DiscretizationGetDMPtrFromIndex(discretization, NFLOWDOF)
    call GridComputeGlobalCell2FaceConnectivity(grid, discretization%MFD, &
                                                  dm_ptr%dm, NFLOWDOF, option)
   endif
#endif
 
  ! initialize to -999.d0 for check later that verifies all values 
  ! have been set
  call VecSet(field%porosity0,-999.d0,ierr)

  ! Allocate vectors to hold temporally average output quantites
  if(realization%output_option%aveg_output_variable_list%nvars>0) then

    field%nvars = realization%output_option%aveg_output_variable_list%nvars
    allocate(field%avg_vars_vec(field%nvars))

    do ivar=1,field%nvars
      call DiscretizationDuplicateVector(discretization,field%porosity0, &
                                         field%avg_vars_vec(ivar))
      call VecSet(field%avg_vars_vec(ivar),0.d0,ierr)
    enddo
  endif
       
  ! Allocate vectors to hold flowrate quantities
  if(realization%output_option%print_hdf5_mass_flowrate.or. &
     realization%output_option%print_hdf5_energy_flowrate.or. &
     realization%output_option%print_hdf5_aveg_mass_flowrate.or. &
     realization%output_option%print_hdf5_aveg_energy_flowrate) then

    call VecCreateMPI(option%mycomm, &
         (option%nflowdof*MAX_FACE_PER_CELL+1)*realization%patch%grid%nlmax, &
          PETSC_DETERMINE,field%flowrate_inst,ierr)
    call VecSet(field%flowrate_inst,0.d0,ierr)

    ! If average flowrate has to be saved, create a vector for it
    if(realization%output_option%print_hdf5_aveg_mass_flowrate.or. &
       realization%output_option%print_hdf5_aveg_energy_flowrate) then
      call VecCreateMPI(option%mycomm, &
          (option%nflowdof*MAX_FACE_PER_CELL+1)*realization%patch%grid%nlmax, &
          PETSC_DETERMINE,field%flowrate_aveg,ierr)
    call VecSet(field%flowrate_aveg,0.d0,ierr)
    endif
  endif

end subroutine RealizationCreateDiscretization

! ************************************************************************** !
!
! RealizationLocalizeRegions: Localizes regions within each patch
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationLocalizeRegions(realization)

  use Option_module
  use String_module
  use Grid_module

  implicit none
  
  type(realization_type) :: realization
  
  type (region_type), pointer :: cur_region, cur_region2
  type(option_type), pointer :: option

  option => realization%option

  ! check to ensure that region names are not duplicated
  cur_region => realization%regions%first
  do
    if (.not.associated(cur_region)) exit
    cur_region2 => cur_region%next
    do
      if (.not.associated(cur_region2)) exit
      if (StringCompare(cur_region%name,cur_region2%name,MAXWORDLENGTH)) then
        option%io_buffer = 'Duplicate region names: ' // trim(cur_region%name)
        call printErrMsg(option)
      endif
      cur_region2 => cur_region2%next
    enddo
    cur_region => cur_region%next
  enddo

  call PatchLocalizeRegions(realization%patch,realization%regions, &
                            realization%option)

end subroutine RealizationLocalizeRegions

! ************************************************************************** !
!
! RealizatonPassPtrsToPatches: Sets patch%field => realization%field
! author: Glenn Hammond
! date: 01/12/11
!
! ************************************************************************** !
subroutine RealizatonPassPtrsToPatches(realization)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  
  realization%patch%field => realization%field
  realization%patch%datasets => realization%datasets
  realization%patch%reaction => realization%reaction
  
end subroutine RealizatonPassPtrsToPatches

! ************************************************************************** !
!
! RealizationAddCoupler: Adds a copy of a coupler to a list
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationAddCoupler(realization,coupler)

  use Coupler_module

  implicit none
  
  type(realization_type) :: realization
  type(coupler_type), pointer :: coupler
  
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: patch
  
  type(coupler_type), pointer :: new_coupler
  
  patch => realization%patch
  
  ! only add to flow list for now, since they will be split out later
  new_coupler => CouplerCreate(coupler)
  select case(coupler%itype)
    case(BOUNDARY_COUPLER_TYPE)
      call CouplerAddToList(new_coupler,patch%boundary_conditions)
    case(INITIAL_COUPLER_TYPE)
      call CouplerAddToList(new_coupler,patch%initial_conditions)
    case(SRC_SINK_COUPLER_TYPE)
      call CouplerAddToList(new_coupler,patch%source_sinks)
  end select
  nullify(new_coupler)

  call CouplerDestroy(coupler)
 
end subroutine RealizationAddCoupler

! ************************************************************************** !
!> This routine sets up nG2LP() mapping for MIMETIC discretization.
!! nG2LP: For a given ghosted cell ID, return the index within the PETSc
!!        solution vector that contains solution at cell centers + cell faces.
!!        The index returned is in PETSc order (0-based).
!!
!> @author
!! ???
!!
!! date: ???
! ************************************************************************** !
subroutine RealizationCreatenG2LP(realization)

  use Grid_module

  implicit none
#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"
#include "finclude/petscmat.h"
#include "finclude/petscmat.h90"
#include "finclude/petscdm.h"
#include "finclude/petscdm.h90"
#include "finclude/petscis.h"
#include "finclude/petscis.h90"
#include "finclude/petscviewer.h"
#include "finclude/petscsnes.h"
#include "finclude/petscpc.h"
#include "definitions.h"

  type(realization_type) :: realization

  type(option_type), pointer :: option
  type(discretization_type), pointer :: discretization
  type(grid_type), pointer :: grid
  PetscInt :: global_offset
  PetscInt :: ghosted_id
  PetscInt :: local_id
  PetscInt :: num_ghosted
  PetscErrorCode :: ierr

  Vec :: vec_LP_cell_id
  Vec :: vec_LP_cell_id_loc

  IS :: is_ghosted, is_global
  VecScatter :: VC_global2ghosted

  PetscScalar, pointer :: lp_cell_ids(:), lp_cell_ids_loc(:)
  PetscInt, pointer :: int_tmp_gh(:), int_tmp_gl(:)

  option => realization%option
  discretization => realization%discretization
  grid => discretization%grid
  
  global_offset = 0
  grid%global_faces_offset = 0
  grid%global_cell_offset = 0

  allocate(grid%nG2LP(grid%ngmax))

  call MPI_Exscan(grid%nlmax_faces, grid%global_faces_offset, &
                  ONE_INTEGER,MPI_INTEGER,MPI_SUM,option%mycomm,ierr)
  call MPI_Exscan(grid%nlmax, grid%global_cell_offset, &
                  ONE_INTEGER,MPI_INTEGER,MPI_SUM,option%mycomm,ierr)

  global_offset = grid%global_faces_offset + grid%global_cell_offset

  call DiscretizationCreateVector(discretization,ONEDOF,vec_LP_cell_id, &
                                  GLOBAL,option)
  call DiscretizationCreateVector(discretization,ONEDOF,vec_LP_cell_id_loc, &
                                  LOCAL,option)
  call VecGetArrayF90(vec_LP_cell_id,lp_cell_ids,ierr)
  do local_id=1,grid%nlmax
    grid%nG2LP(grid%nL2G(local_id))=global_offset+grid%nlmax_faces+local_id-1
    lp_cell_ids(local_id)=global_offset+grid%nlmax_faces+local_id
  enddo
  call VecRestoreArrayF90(vec_LP_cell_id,lp_cell_ids,ierr)

  allocate(int_tmp_gh(grid%ngmax-grid%nlmax))
  allocate(int_tmp_gl(grid%ngmax-grid%nlmax))

  num_ghosted = 1
  do ghosted_id = 1,grid%ngmax
    if (grid%nG2L(ghosted_id) < 1) then
      int_tmp_gh(num_ghosted) = ghosted_id - 1
      int_tmp_gl(num_ghosted) = grid%nG2P(ghosted_id)
      num_ghosted=num_ghosted+1
    endif
  enddo

  call ISCreateBlock(option%mycomm, ONEDOF, grid%ngmax - grid%nlmax, &
                     int_tmp_gh, PETSC_COPY_VALUES, is_ghosted, ierr)
  call ISCreateBlock(option%mycomm, ONEDOF, grid%ngmax - grid%nlmax, &
                     int_tmp_gl, PETSC_COPY_VALUES, is_global, ierr)
  call VecScatterCreate(vec_LP_cell_id, is_global, vec_LP_cell_id_loc, &
                        is_ghosted, VC_global2ghosted, ierr)
  deallocate(int_tmp_gh)
  deallocate(int_tmp_gl)

  call VecScatterBegin(VC_global2ghosted, vec_LP_cell_id, vec_LP_cell_id_loc, &
                      INSERT_VALUES,SCATTER_FORWARD,ierr)
  call VecScatterEnd(VC_global2ghosted, vec_LP_cell_id, vec_LP_cell_id_loc, &
                      INSERT_VALUES,SCATTER_FORWARD,ierr)

  call VecGetArrayF90(vec_LP_cell_id_loc,lp_cell_ids_loc,ierr)
  do ghosted_id=1,grid%ngmax
    if (grid%nG2L(ghosted_id)<1) then
      grid%nG2LP(ghosted_id)=int(lp_cell_ids_loc(ghosted_id))-1
    endif
  end do
  call VecRestoreArrayF90(vec_LP_cell_id_loc,lp_cell_ids_loc,ierr)

  call VecDestroy(vec_LP_cell_id, ierr)
  call VecDestroy(vec_LP_cell_id_loc, ierr)
  call VecScatterDestroy(VC_global2ghosted , ierr)
  call ISDestroy(is_ghosted, ierr)
  call ISDestroy(is_global, ierr)

end subroutine RealizationCreatenG2LP

! ************************************************************************** !
!
! RealizationAddStrata: Adds a copy of a strata to a list
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationAddStrata(realization,strata)

  use Strata_module

  implicit none
  
  type(realization_type) :: realization
  type(strata_type), pointer :: strata
  
  type(strata_type), pointer :: new_strata
  
  new_strata => StrataCreate(strata)
  call StrataAddToList(new_strata,realization%patch%strata)
  nullify(new_strata)
  
  call StrataDestroy(strata)
 
end subroutine RealizationAddStrata

! ************************************************************************** !
!
! RealizationAddObservation: Adds a copy of a observation object to a list
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationAddObservation(realization,observation)

  use Observation_module

  implicit none
  
  type(realization_type) :: realization
  type(observation_type), pointer :: observation
  
  type(observation_type), pointer :: new_observation
  
  new_observation => ObservationCreate(observation)
  call ObservationAddToList(new_observation, &
                            realization%patch%observation)
  nullify(new_observation)

  call ObservationDestroy(observation)
 
end subroutine RealizationAddObservation

! ************************************************************************** !
!
! RealizationProcessCouplers: Sets connectivity and pointers for couplers
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationProcessCouplers(realization)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  
  call PatchProcessCouplers( realization%patch,realization%flow_conditions, &
                             realization%transport_conditions, &
                             realization%option)
  
end subroutine RealizationProcessCouplers

! ************************************************************************** !
!
! RealizationProcessConditions: Sets up auxiliary data associated with 
!                               conditions
! author: Glenn Hammond
! date: 10/14/08
!
! ************************************************************************** !
subroutine RealizationProcessConditions(realization)

  use Dataset_module
  
  implicit none
  
  type(realization_type) :: realization

  call DatasetProcessDatasets(realization%datasets,realization%option)
  
  if (realization%option%nflowdof > 0) then
    call RealProcessFlowConditions(realization)
  endif
  if (realization%option%ntrandof > 0) then
    call RealProcessTranConditions(realization)
  endif
 
end subroutine RealizationProcessConditions

! ************************************************************************** !
!
! RealProcessMatPropAndSatFunc: Sets up linkeage between material properties
!                               and saturation function, auxiliary arrays
!                               and datasets
! author: Glenn Hammond
! date: 01/21/09, 01/12/11
!
! ************************************************************************** !
subroutine RealProcessMatPropAndSatFunc(realization)

  use String_module
  
  implicit none
  
  type(realization_type) :: realization
  
  PetscBool :: found
  PetscInt :: i
  type(option_type), pointer :: option
  type(material_property_type), pointer :: cur_material_property
  type(patch_type), pointer :: patch
  character(len=MAXSTRINGLENGTH) :: string

  option => realization%option
  patch => realization%patch
  
  ! organize lists
  call MaterialPropConvertListToArray(realization%material_properties, &
                                      realization%material_property_array, &
                                      option)
  call SaturatFuncConvertListToArray(realization%saturation_functions, &
                                      realization%saturation_function_array, &
                                      option)

  ! set up mirrored pointer arrays within patches to saturation functions
  ! and material properties
  patch%material_properties => realization%material_properties
  call MaterialPropConvertListToArray(patch%material_properties, &
                                      patch%material_property_array, &
                                      option)
  patch%saturation_functions => realization%saturation_functions
  call SaturatFuncConvertListToArray(patch%saturation_functions, &
                                      patch%saturation_function_array, &
                                      option)
    
  cur_material_property => realization%material_properties                            
  do                                      
    if (.not.associated(cur_material_property)) exit

    ! obtain saturation function id
    if (option%iflowmode /= NULL_MODE) then
      cur_material_property%saturation_function_id = &
        SaturationFunctionGetID(realization%saturation_functions, &
                                cur_material_property%saturation_function_name, &
                                cur_material_property%name,option)
    endif
    
    ! if named, link dataset to property
    if (.not.StringNull(cur_material_property%porosity_dataset_name)) then
      string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
               '),POROSITY'
      cur_material_property%porosity_dataset => &
        DatasetGetPointer(realization%datasets, &
                          cur_material_property%porosity_dataset_name, &
                          string,option)
    endif
    if (.not.StringNull(cur_material_property%permeability_dataset_name)) then
      string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
               '),PERMEABILITY'
      cur_material_property%permeability_dataset => &
        DatasetGetPointer(realization%datasets, &
                          cur_material_property%permeability_dataset_name, &
                          string,option)
    endif
    
    cur_material_property => cur_material_property%next
  enddo
  
  
end subroutine RealProcessMatPropAndSatFunc

! ************************************************************************** !
!
! RealProcessFluidProperties: Sets up linkeage with fluid properties
! author: Glenn Hammond
! date: 01/21/09
!
! ************************************************************************** !
subroutine RealProcessFluidProperties(realization)
  
  implicit none
  
  type(realization_type) :: realization
  
  PetscBool :: found
  type(option_type), pointer :: option
  type(fluid_property_type), pointer :: cur_fluid_property
  
  option => realization%option
  
  found = PETSC_FALSE
  cur_fluid_property => realization%fluid_properties                            
  do                                      
    if (.not.associated(cur_fluid_property)) exit
    found = PETSC_TRUE
    select case(trim(cur_fluid_property%phase_name))
      case('LIQUID_PHASE')
        cur_fluid_property%phase_id = LIQUID_PHASE
      case('GAS_PHASE')
        cur_fluid_property%phase_id = GAS_PHASE
      case default
        cur_fluid_property%phase_id = LIQUID_PHASE
    end select
    cur_fluid_property => cur_fluid_property%next
  enddo
  
  if (option%ntrandof > 0 .and. .not.found) then
    option%io_buffer = 'A fluid property must be present in input file' // &
                       ' for solute transport'
  endif
  
end subroutine RealProcessFluidProperties

! ************************************************************************** !
!
! RealProcessFlowConditions: Sets linkage of flow conditions to dataset
! author: Glenn Hammond
! date: 10/26/11
!
! ************************************************************************** !
subroutine RealProcessFlowConditions(realization)

  use Dataset_module

  implicit none

  type(realization_type) :: realization
  
  type(flow_condition_type), pointer :: cur_flow_condition
  type(flow_sub_condition_type), pointer :: cur_flow_sub_condition
  type(option_type), pointer :: option
  character(len=MAXSTRINGLENGTH) :: string
  character(len=MAXWORDLENGTH) :: dataset_name
  PetscInt :: i
  type(dataset_type), pointer :: dataset
  
  option => realization%option
  
  ! loop over flow conditions looking for linkage to datasets
  cur_flow_condition => realization%flow_conditions%first
  do
    if (.not.associated(cur_flow_condition)) exit
    !TODO(geh): could destroy the time_series here if dataset allocated
    select case(option%iflowmode)
      case(G_MODE)
      case(RICHARDS_MODE,MIS_MODE,TH_MODE)
        do i = 1, size(cur_flow_condition%sub_condition_ptr)
          ! check for dataset in flow_dataset
          if (associated(cur_flow_condition%sub_condition_ptr(i)%ptr% &
                          flow_dataset%dataset)) then
            dataset_name = cur_flow_condition%sub_condition_ptr(i)%ptr% &
                          flow_dataset%dataset%name
            ! delete the dataset since it is solely a placeholder
            call DatasetDestroy(cur_flow_condition%sub_condition_ptr(i)%ptr% &
                                flow_dataset%dataset)
            ! get dataset from list
            string = 'flow_condition ' // trim(cur_flow_condition%name)
            dataset => &
              DatasetGetPointer(realization%datasets,dataset_name,string,option)
            cur_flow_condition%sub_condition_ptr(i)%ptr%flow_dataset%dataset => &
              dataset
            nullify(dataset)
          endif
          if (associated(cur_flow_condition%sub_condition_ptr(i)%ptr% &
                          datum%dataset)) then
            dataset_name = cur_flow_condition%sub_condition_ptr(i)%ptr% &
                          datum%dataset%name
            ! delete the dataset since it is solely a placeholder
            call DatasetDestroy(cur_flow_condition%sub_condition_ptr(i)%ptr% &
                                datum%dataset)
            ! get dataset from list
            string = 'flow_condition ' // trim(cur_flow_condition%name)
            dataset => &
              DatasetGetPointer(realization%datasets,dataset_name,string,option)
            cur_flow_condition%sub_condition_ptr(i)%ptr%datum%dataset => &
              dataset
            nullify(dataset)
          endif
          if (associated(cur_flow_condition%sub_condition_ptr(i)%ptr% &
                          gradient%dataset)) then
            dataset_name = cur_flow_condition%sub_condition_ptr(i)%ptr% &
                          gradient%dataset%name
            ! delete the dataset since it is solely a placeholder
            call DatasetDestroy(cur_flow_condition%sub_condition_ptr(i)%ptr% &
                                gradient%dataset)
            ! get dataset from list
            string = 'flow_condition ' // trim(cur_flow_condition%name)
            dataset => &
              DatasetGetPointer(realization%datasets,dataset_name,string,option)
            cur_flow_condition%sub_condition_ptr(i)%ptr%gradient%dataset => &
              dataset
            nullify(dataset)
          endif
        enddo
    end select
    cur_flow_condition => cur_flow_condition%next
  enddo

end subroutine RealProcessFlowConditions

! ************************************************************************** !
!
! RealProcessTranConditions: Sets up auxiliary data associated with 
!                            transport conditions
! author: Glenn Hammond
! date: 10/14/08
!
! ************************************************************************** !
subroutine RealProcessTranConditions(realization)

  use String_module
  use Reaction_module
  use Constraint_module
  
  implicit none
  
  type(realization_type) :: realization
  
  
  PetscBool :: found
  type(option_type), pointer :: option
  type(tran_condition_type), pointer :: cur_condition
  type(tran_constraint_coupler_type), pointer :: cur_constraint_coupler
  type(tran_constraint_type), pointer :: cur_constraint, another_constraint
  
  option => realization%option
  
  ! check for duplicate constraint names
  cur_constraint => realization%transport_constraints%first
  do
    if (.not.associated(cur_constraint)) exit
      another_constraint => cur_constraint%next
      ! now compare names
      found = PETSC_FALSE
      do
        if (.not.associated(another_constraint)) exit
        if (StringCompare(cur_constraint%name,another_constraint%name, &
            MAXWORDLENGTH)) then
          found = PETSC_TRUE
        endif
        another_constraint => another_constraint%next
      enddo
      if (found) then
        option%io_buffer = 'Duplicate transport constraints named "' // &
                 trim(cur_constraint%name) // '"'
        call printErrMsg(realization%option)
      endif
    cur_constraint => cur_constraint%next
  enddo
  
  ! initialize constraints
  cur_constraint => realization%transport_constraints%first
  do
    if (.not.associated(cur_constraint)) exit
    call ReactionProcessConstraint(realization%reaction, &
                                   cur_constraint%name, &
                                   cur_constraint%aqueous_species, &
                                   cur_constraint%minerals, &
                                   cur_constraint%surface_complexes, &
                                   cur_constraint%colloids, &
                                   cur_constraint%immobile_species, &
                                   realization%option)
    cur_constraint => cur_constraint%next
  enddo
  
  if (option%use_mc) then
    call ReactionProcessConstraint(realization%reaction, &
                                   realization%sec_transport_constraint%name, &
                                   realization%sec_transport_constraint%aqueous_species, &
                                   realization%sec_transport_constraint%minerals, &
                                   realization%sec_transport_constraint%surface_complexes, &
                                   realization%sec_transport_constraint%colloids, &
                                   realization%sec_transport_constraint%immobile_species, &
                                   realization%option)
  endif
  
  ! tie constraints to couplers, if not already associated
  cur_condition => realization%transport_conditions%first
  do

    if (.not.associated(cur_condition)) exit
    cur_constraint_coupler => cur_condition%constraint_coupler_list
    do
      if (.not.associated(cur_constraint_coupler)) exit
      if (.not.associated(cur_constraint_coupler%aqueous_species)) then
        cur_constraint => realization%transport_constraints%first
        do
          if (.not.associated(cur_constraint)) exit
          if (StringCompare(cur_constraint%name, &
                             cur_constraint_coupler%constraint_name, &
                             MAXWORDLENGTH)) then
            cur_constraint_coupler%aqueous_species => cur_constraint%aqueous_species
            cur_constraint_coupler%minerals => cur_constraint%minerals
            cur_constraint_coupler%surface_complexes => cur_constraint%surface_complexes
            cur_constraint_coupler%colloids => cur_constraint%colloids
            cur_constraint_coupler%immobile_species => cur_constraint%immobile_species
            exit
          endif
          cur_constraint => cur_constraint%next
        enddo
        if (.not.associated(cur_constraint_coupler%aqueous_species)) then
          option%io_buffer = 'Transport constraint "' // &
                   trim(cur_constraint_coupler%constraint_name) // &
                   '" not found in input file constraints.'
          call printErrMsg(realization%option)
        endif
      endif
      cur_constraint_coupler => cur_constraint_coupler%next
    enddo
    if (associated(cur_condition%constraint_coupler_list%next)) then ! more than one
      cur_condition%is_transient = PETSC_TRUE
    else
      cur_condition%is_transient = PETSC_FALSE
    endif
    cur_condition => cur_condition%next
  enddo
 
  ! final details for setup
  cur_condition => realization%transport_conditions%first
  do
    if (.not.associated(cur_condition)) exit
    ! is the condition transient?
    if (associated(cur_condition%constraint_coupler_list%next)) then ! more than one
      cur_condition%is_transient = PETSC_TRUE
    else
      cur_condition%is_transient = PETSC_FALSE
    endif
    ! set pointer to first constraint coupler
    cur_condition%cur_constraint_coupler => cur_condition%constraint_coupler_list
    
    cur_condition => cur_condition%next
  enddo

end subroutine RealProcessTranConditions

! ************************************************************************** !
!
! RealizationInitConstraints: Initializes constraint concentrations
! author: Glenn Hammond
! date: 12/04/08
!
! ************************************************************************** !
subroutine RealizationInitConstraints(realization)

  implicit none

  type(realization_type) :: realization
  
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  
  cur_level => realization%level_list%first
  do 
    if (.not.associated(cur_level)) exit
    cur_patch => cur_level%patch_list%first
    do
      if (.not.associated(cur_patch)) exit
      call PatchInitConstraints(cur_patch,realization%reaction, &
                                realization%option)
      cur_patch => cur_patch%next
    enddo
    cur_level => cur_level%next
  enddo            
 
end subroutine RealizationInitConstraints

! ************************************************************************** !
!
! RealizationPrintCouplers: Print boundary and initial condition data
! author: Glenn Hammond
! date: 10/28/08
!
! ************************************************************************** !
subroutine RealizationPrintCouplers(realization)

  use Coupler_module
  
  implicit none
  
  type(realization_type) :: realization
  
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  type(coupler_type), pointer :: cur_coupler
  type(option_type), pointer :: option
  type(reaction_type), pointer :: reaction
 
  option => realization%option
  reaction => realization%reaction
 
  if (.not.OptionPrintToFile(option)) return
  
  cur_level => realization%level_list%first
  do 
    if (.not.associated(cur_level)) exit
    cur_patch => cur_level%patch_list%first
    do
      if (.not.associated(cur_patch)) exit

      cur_coupler => cur_patch%initial_conditions%first
      do
        if (.not.associated(cur_coupler)) exit
        call RealizationPrintCoupler(cur_coupler,reaction,option)    
        cur_coupler => cur_coupler%next
      enddo
     
      cur_coupler => cur_patch%boundary_conditions%first
      do
        if (.not.associated(cur_coupler)) exit
        call RealizationPrintCoupler(cur_coupler,reaction,option)    
        cur_coupler => cur_coupler%next
      enddo
     
      cur_coupler => cur_patch%source_sinks%first
      do
        if (.not.associated(cur_coupler)) exit
        call RealizationPrintCoupler(cur_coupler,reaction,option)    
        cur_coupler => cur_coupler%next
      enddo

      cur_patch => cur_patch%next
    enddo
    cur_level => cur_level%next
  enddo            
 
end subroutine RealizationPrintCouplers

! ************************************************************************** !
!
! RealizationPrintCoupler: Prints boundary and initial condition coupler 
! author: Glenn Hammond
! date: 10/28/08
!
! ************************************************************************** !
subroutine RealizationPrintCoupler(coupler,reaction,option)

  use Coupler_module
  use Reaction_module
  
  implicit none
  
  type(coupler_type) :: coupler
  type(option_type) :: option
  type(reaction_type), pointer :: reaction
  
  character(len=MAXSTRINGLENGTH) :: string
  
  type(flow_condition_type), pointer :: flow_condition
  type(tran_condition_type), pointer :: tran_condition
  type(region_type), pointer :: region
  type(tran_constraint_coupler_type), pointer :: constraint_coupler
   
98 format(40('=+'))
99 format(80('-'))
100 format(a)
  
  flow_condition => coupler%flow_condition
  tran_condition => coupler%tran_condition
  region => coupler%region

  write(option%fid_out,*)
  write(option%fid_out,98)


  select case(coupler%itype)
    case(INITIAL_COUPLER_TYPE)
      string = 'Initial Condition'
    case(BOUNDARY_COUPLER_TYPE)
      string = 'Boundary Condition'
    case(SRC_SINK_COUPLER_TYPE)
      string = 'Source Sink'
  end select
  write(option%fid_out,'(/,2x,a,/)') trim(string)

  write(option%fid_out,99)
101 format(5x,'     Flow Condition: ',2x,a)
  if (associated(flow_condition)) &
    write(option%fid_out,101) trim(flow_condition%name)
102 format(5x,'Transport Condition: ',2x,a)
  if (associated(tran_condition)) &
    write(option%fid_out,102) trim(tran_condition%name)
103 format(5x,'             Region: ',2x,a)
  if (associated(region)) &
    write(option%fid_out,103) trim(region%name)
  write(option%fid_out,99)
  
  if (associated(flow_condition)) then
    call FlowConditionPrint(flow_condition,option)
  endif
  if (associated(tran_condition)) then
    constraint_coupler => tran_condition%cur_constraint_coupler
    write(option%fid_out,'(/,2x,''Transport Condition: '',a)') &
      trim(tran_condition%name)
    if (associated(reaction)) then
      call ReactionPrintConstraint(constraint_coupler,reaction,option)
      write(option%fid_out,'(/)')
      write(option%fid_out,99)
    endif
  endif
 
end subroutine RealizationPrintCoupler

! ************************************************************************** !
!
! RealizationInitCouplerAuxVars: Initializes coupler auxillary variables 
!                                within list
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizationInitAllCouplerAuxVars(realization)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  
  !geh: Must update conditions prior to initializing the aux vars.  
  !     Otherwise, datasets will not have been read for routines such as
  !     hydrostatic and auxvars will be initialized to garbage.
  call FlowConditionUpdate(realization%flow_conditions,realization%option, &
                           realization%option%time)
  call TranConditionUpdate(realization%transport_conditions, &
                           realization%option, &
                           realization%option%time)  
  call PatchInitAllCouplerAuxVars(realization%patch,realization%option)
   
end subroutine RealizationInitAllCouplerAuxVars

! ************************************************************************** !
!
! RealizUpdateAllCouplerAuxVars: Updates auxiliary variables associated 
!                                  with couplers in lis
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizUpdateAllCouplerAuxVars(realization,force_update_flag)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  PetscBool :: force_update_flag

  call PatchUpdateAllCouplerAuxVars(realization%patch,force_update_flag, &
                                    realization%option)

end subroutine RealizUpdateAllCouplerAuxVars

! ************************************************************************** !
!
! RealizationUpdate: Update parameters in realization (e.g. conditions, bcs, srcs)
! author: Glenn Hammond
! date: 11/09/07
!
! ************************************************************************** !
subroutine RealizationUpdate(realization)

  implicit none
  
  type(realization_type) :: realization
  
  PetscBool :: force_update_flag = PETSC_FALSE
  
  ! must update conditions first
  call FlowConditionUpdate(realization%flow_conditions,realization%option, &
                           realization%option%time)
  call TranConditionUpdate(realization%transport_conditions, &
                           realization%option, &
                           realization%option%time)
  call RealizUpdateAllCouplerAuxVars(realization,force_update_flag)
  if (associated(realization%uniform_velocity_dataset)) then
    call RealizUpdateUniformVelocity(realization)
  endif
! currently don't use aux_vars, just condition for src/sinks
!  call RealizationUpdateSrcSinks(realization)

end subroutine RealizationUpdate

! ************************************************************************** !
!
! RealizationRevertFlowParameters: Assigns initial porosity/perms to vecs
! author: Glenn Hammond
! date: 05/09/08
!
! ************************************************************************** !
subroutine RealizationRevertFlowParameters(realization)

  use Option_module
  use Field_module
  use Discretization_module

  implicit none
  
  type(realization_type) :: realization
  
  type(field_type), pointer :: field
  type(option_type), pointer :: option
  type(discretization_type), pointer :: discretization
  
  option => realization%option
  field => realization%field
  discretization => realization%discretization

  if (option%nflowdof > 0) then
    call DiscretizationGlobalToLocal(discretization,field%perm0_xx, &
                           field%perm_xx_loc,ONEDOF)  
    call DiscretizationGlobalToLocal(discretization,field%perm0_yy, &
                           field%perm_yy_loc,ONEDOF)  
    call DiscretizationGlobalToLocal(discretization,field%perm0_zz, &
                           field%perm_zz_loc,ONEDOF)   
  endif   
  call DiscretizationGlobalToLocal(discretization,field%porosity0, &
                                   field%porosity_loc,ONEDOF)
  call DiscretizationGlobalToLocal(discretization,field%tortuosity0, &
                                   field%tortuosity_loc,ONEDOF)
                           
end subroutine RealizationRevertFlowParameters

! ************************************************************************** !
!
! RealizUpdateUniformVelocity: Assigns uniform velocity for transport
! author: Glenn Hammond
! date: 02/22/08
!
! ************************************************************************** !
subroutine RealizUpdateUniformVelocity(realization)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  
  call UniformVelocityDatasetUpdate(realization%option, &
                                    realization%option%time, &
                                    realization%uniform_velocity_dataset)
  call PatchUpdateUniformVelocity(realization%patch, &
                            realization%uniform_velocity_dataset%cur_value, &
                            realization%option)
 
end subroutine RealizUpdateUniformVelocity

! ************************************************************************** !
!
! RealizationAddWaypointsToList: Creates waypoints associated with source/sinks
!                             boundary conditions, etc. and add to list
! author: Glenn Hammond
! date: 11/01/07
!
! ************************************************************************** !
subroutine RealizationAddWaypointsToList(realization)

  use Option_module
  use Waypoint_module

  implicit none
  
  type(realization_type) :: realization
  
  type(waypoint_list_type), pointer :: waypoint_list
  type(flow_condition_type), pointer :: cur_flow_condition
  type(tran_condition_type), pointer :: cur_tran_condition
  type(flow_sub_condition_type), pointer :: sub_condition
  type(tran_constraint_coupler_type), pointer :: cur_constraint_coupler
  type(waypoint_type), pointer :: waypoint, cur_waypoint
  type(option_type), pointer :: option
  PetscInt :: itime, isub_condition
  PetscReal :: temp_real, final_time
  PetscReal, pointer :: times(:)

  option => realization%option
  waypoint_list => realization%waypoints
  nullify(times)
  
  ! set flag for final output
  cur_waypoint => waypoint_list%first
  do
    if (.not.associated(cur_waypoint)) exit
    if (cur_waypoint%final) then
      cur_waypoint%print_output = realization%output_option%print_final
      exit
    endif
    cur_waypoint => cur_waypoint%next
  enddo
  ! use final time in conditional below
  if (associated(cur_waypoint)) then
    final_time = cur_waypoint%time
  else
    option%io_buffer = 'Final time not found in RealizationAddWaypointsToList'
    call printErrMsg(option)
  endif

  ! add update of flow conditions
  cur_flow_condition => realization%flow_conditions%first
  do
    if (.not.associated(cur_flow_condition)) exit
    if (cur_flow_condition%sync_time_with_update) then
      do isub_condition = 1, cur_flow_condition%num_sub_conditions
        sub_condition => cur_flow_condition%sub_condition_ptr(isub_condition)%ptr
        !TODO(geh): check if this updated more than simply the flow_dataset (i.e. datum and gradient)
        !geh: followup - no, datum/gradient are not considered.  Should they be considered?
        call FlowConditionDatasetGetTimes(option,sub_condition,final_time, &
                                          times)
        if (size(times) > 1000) then
          option%io_buffer = 'For flow condition "' // &
            trim(cur_flow_condition%name) // &
            '" dataset "' // trim(sub_condition%name) // &
            '", the number of times is excessive for synchronization ' // &
            'with waypoints.'
          call printErrMsg(option)
        endif
        do itime = 1, size(times)
          waypoint => WaypointCreate()
          waypoint%time = times(itime)
          waypoint%update_conditions = PETSC_TRUE
          call WaypointInsertInList(waypoint,waypoint_list)
        enddo
        deallocate(times)
        nullify(times)
      enddo
    endif
    cur_flow_condition => cur_flow_condition%next
  enddo
      
  ! add update of transport conditions
  cur_tran_condition => realization%transport_conditions%first
  do
    if (.not.associated(cur_tran_condition)) exit
    if (cur_tran_condition%is_transient) then
      cur_constraint_coupler => cur_tran_condition%constraint_coupler_list
      do
        if (.not.associated(cur_constraint_coupler)) exit
        if (cur_constraint_coupler%time > 1.d-40) then
          waypoint => WaypointCreate()
          waypoint%time = cur_constraint_coupler%time
          waypoint%update_conditions = PETSC_TRUE
          call WaypointInsertInList(waypoint,waypoint_list)
        endif
        cur_constraint_coupler => cur_constraint_coupler%next
      enddo
    endif
    cur_tran_condition => cur_tran_condition%next
  enddo

  ! add update of velocity fields
  if (associated(realization%uniform_velocity_dataset)) then
    if (realization%uniform_velocity_dataset%times(1) > 1.d-40 .or. &
        size(realization%uniform_velocity_dataset%times) > 1) then
      do itime = 1, size(realization%uniform_velocity_dataset%times)
        waypoint => WaypointCreate()
        waypoint%time = realization%uniform_velocity_dataset%times(itime)
        waypoint%update_conditions = PETSC_TRUE
        call WaypointInsertInList(waypoint,waypoint_list)
      enddo
    endif
  endif

  ! add waypoints for periodic output
  if (realization%output_option%periodic_output_time_incr > 0.d0 .or. &
      realization%output_option%periodic_tr_output_time_incr > 0.d0) then

    if (realization%output_option%periodic_output_time_incr > 0.d0) then
      ! standard output
      temp_real = 0.d0
      do
        temp_real = temp_real + realization%output_option%periodic_output_time_incr
        if (temp_real > final_time) exit
        waypoint => WaypointCreate()
        waypoint%time = temp_real
        waypoint%print_output = PETSC_TRUE
        call WaypointInsertInList(waypoint,realization%waypoints)
      enddo
    endif
    
    if (realization%output_option%periodic_tr_output_time_incr > 0.d0) then
      ! transient observation output
      temp_real = 0.d0
      do
        temp_real = temp_real + realization%output_option%periodic_tr_output_time_incr
        if (temp_real > final_time) exit
        waypoint => WaypointCreate()
        waypoint%time = temp_real
        waypoint%print_tr_output = PETSC_TRUE 
        call WaypointInsertInList(waypoint,realization%waypoints)
      enddo
    endif

  endif

end subroutine RealizationAddWaypointsToList

#if 0
! ************************************************************************** !
!
! RealizationGetDataset: Extracts variables indexed by ivar and isubvar from a 
!                        realization
! author: Glenn Hammond
! date: 09/12/08
!
! ************************************************************************** !
subroutine RealizationGetDataset(realization,vec,ivar,isubvar,isubvar1)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  Vec :: vec
  PetscInt :: ivar
  PetscInt :: isubvar
  PetscInt, optional :: isubvar1
  
  call PatchGetDataset(realization%patch,realization%field, &
                       realization%reaction,realization%option, &
                       realization%output_option,vec,ivar,isubvar,isubvar1)

end subroutine RealizationGetDataset

! ************************************************************************** !
!
! RealizGetDatasetValueAtCell: Extracts variables indexed by ivar and isubvar
!                              from a realization
! author: Glenn Hammond
! date: 09/12/08
!
! ************************************************************************** !
function RealizGetDatasetValueAtCell(realization,ivar,isubvar,ghosted_id, &
                                     isubvar1)

  use Option_module

  implicit none
  
  PetscReal :: RealizGetDatasetValueAtCell
  type(realization_type) :: realization
  PetscInt :: ivar
  PetscInt :: isubvar
  PetscInt, optional :: isubvar1
  PetscInt :: ghosted_id
  
  PetscReal :: value
  
  value = PatchGetDatasetValueAtCell(realization%patch,realization%field, &
                                     realization%reaction, &
                                     realization%option, &
                                     realization%output_option, &
                                     ivar,isubvar,ghosted_id,isubvar1)
  RealizGetDatasetValueAtCell = value

end function RealizGetDatasetValueAtCell

! ************************************************************************** !
!
! RealizationSetDataset: Sets variables indexed by ivar and isubvar in a 
!                        realization
! author: Glenn Hammond
! date: 09/12/08
!
! ************************************************************************** !
subroutine RealizationSetDataset(realization,vec,vec_format,ivar,isubvar)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  Vec :: vec
  PetscInt :: vec_format
  PetscInt :: ivar
  PetscInt :: isubvar

  call PatchSetDataset(realization%patch,realization%field, &
                       realization%option, &
                       vec,vec_format,ivar,isubvar)

end subroutine RealizationSetDataset

#endif

! ************************************************************************** !
!
! RealizationUpdateProperties: Updates coupled properties at each grid cell
! author: Glenn Hammond
! date: 08/05/09
!
! ************************************************************************** !
subroutine RealizationUpdateProperties(realization)

  implicit none

  type(realization_type) :: realization
  
  type(option_type), pointer :: option  
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  PetscReal :: min_value  
  PetscInt :: ivalue
  PetscErrorCode :: ierr
  
  option => realization%option
    
  call RealizationUpdatePropertiesPatch(realization)
  
  ! perform check to ensure that porosity is bounded between 0 and 1
  ! since it is calculated as 1.d-sum_volfrac, it cannot be > 1
  call VecMin(realization%field%porosity_loc,ivalue,min_value,ierr)
  if (min_value < 0.d0) then
    write(option%io_buffer,*) 'Sum of mineral volume fractions has ' // &
      'exceeded 1.d0 at cell (note PETSc numbering): ', ivalue
    call printErrMsg(option)
  endif
   
end subroutine RealizationUpdateProperties

! ************************************************************************** !
!
! RealizationUpdatePropertiesPatch: Updates coupled properties at each grid cell 
! author: Glenn Hammond
! date: 08/05/09
!
! ************************************************************************** !
subroutine RealizationUpdatePropertiesPatch(realization)

  use Grid_module
  use Reactive_Transport_Aux_module
 
  implicit none
  
  type(realization_type) :: realization

  type(option_type), pointer :: option
  type(patch_type), pointer :: patch
  type(field_type), pointer :: field
  type(reaction_type), pointer :: reaction
  type(grid_type), pointer :: grid
  type(material_property_ptr_type), pointer :: material_property_array(:)
  type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:) 
  type(discretization_type), pointer :: discretization

  PetscInt :: local_id, ghosted_id
  PetscInt :: imnrl, imnrl1, imnrl_armor
  PetscReal :: sum_volfrac
  PetscReal :: scale, porosity_scale, volfrac_scale
  PetscBool :: porosity_updated
  PetscReal, pointer :: vec_p(:)
  PetscReal, pointer :: porosity_loc_p(:), porosity0_p(:)
  PetscReal, pointer :: tortuosity_loc_p(:), tortuosity0_p(:)
  PetscReal, pointer :: perm0_xx_p(:), perm0_yy_p(:), perm0_zz_p(:)
  PetscReal, pointer :: perm_xx_loc_p(:), perm_yy_loc_p(:), perm_zz_loc_p(:)
  PetscErrorCode :: ierr

  option => realization%option
  discretization => realization%discretization
  patch => realization%patch
  field => realization%field
  reaction => realization%reaction
  grid => patch%grid
  material_property_array => realization%material_property_array
  rt_auxvars => patch%aux%RT%aux_vars

  if (.not.associated(patch%imat)) then
    option%io_buffer = 'Materials IDs not present in run.  Material ' // &
      ' properties cannot be updated without material ids ask Glenn'
    call printErrMsg(option)
  endif

  porosity_updated = PETSC_FALSE
  if (reaction%update_porosity) then
    porosity_updated = PETSC_TRUE
  
    call GridVecGetArrayF90(grid,field%porosity_loc,porosity_loc_p,ierr)

    if (reaction%mineral%nkinmnrl > 0) then
      do local_id = 1, grid%nlmax
        ghosted_id = grid%nL2G(local_id)

        ! Go ahead and compute for inactive cells since their porosity does
        ! not matter (avoid check on active/inactive)
        sum_volfrac = 0.d0
        do imnrl = 1, reaction%mineral%nkinmnrl
          sum_volfrac = sum_volfrac + &
                        rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl)
        enddo 
        porosity_loc_p(ghosted_id) = max(1.d0-sum_volfrac, &
                                         reaction%minimum_porosity)
      enddo
    endif

    call GridVecRestoreArrayF90(grid,field%porosity_loc,porosity_loc_p,ierr)
    
  endif
  
  if ((porosity_updated .and. &
       (reaction%update_tortuosity .or. &
        reaction%update_permeability)) .or. &
      ! if porosity ratio is used in mineral surface area update, we must
      ! recalculate it every time.
      (reaction%update_mineral_surface_area .and. &
       reaction%update_mnrl_surf_with_porosity)) then
    call GridVecGetArrayF90(grid,field%porosity0,porosity0_p,ierr)
    call GridVecGetArrayF90(grid,field%porosity_loc,porosity_loc_p,ierr)
    call GridVecGetArrayF90(grid,field%work,vec_p,ierr)
    do local_id = 1, grid%nlmax
      ghosted_id = grid%nL2G(local_id)
      vec_p(local_id) = porosity_loc_p(ghosted_id)/porosity0_p(local_id)
    enddo
    call GridVecRestoreArrayF90(grid,field%porosity0,porosity0_p,ierr)
    call GridVecRestoreArrayF90(grid,field%porosity_loc,porosity_loc_p,ierr)
    call GridVecRestoreArrayF90(grid,field%work,vec_p,ierr)
  endif      

  if (reaction%update_mineral_surface_area) then

    if (reaction%update_mnrl_surf_with_porosity) then
      ! placing the get/restore array calls within the condition will
      ! avoid improper access.
      call GridVecGetArrayF90(grid,field%work,vec_p,ierr)
    endif

    do local_id = 1, grid%nlmax
      ghosted_id = grid%nL2G(local_id)
      do imnrl = 1, reaction%mineral%nkinmnrl

        porosity_scale = 1.d0
        if (reaction%update_mnrl_surf_with_porosity) then
          porosity_scale = vec_p(local_id)** &
             reaction%mineral%kinmnrl_surf_area_porosity_pwr(imnrl)
!       geh: srf_area_vol_frac_pwr must be defined on a per mineral basis, not
!       solely material type.
!       material_property_array(patch%imat(ghosted_id))%ptr%mnrl_surf_area_porosity_pwr
        endif

        volfrac_scale = 1.d0
        if (rt_auxvars(ghosted_id)%mnrl_volfrac0(imnrl) > 0.d0) then
          volfrac_scale = (rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl)/ &
                         rt_auxvars(ghosted_id)%mnrl_volfrac0(imnrl))** &
             reaction%mineral%kinmnrl_surf_area_vol_frac_pwr(imnrl)
!       geh: srf_area_vol_frac_pwr must be defined on a per mineral basis, not
!       solely material type.
!       material_property_array(patch%imat(ghosted_id))%ptr%mnrl_surf_area_volfrac_pwr
!         rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
!           rt_auxvars(ghosted_id)%mnrl_area0(imnrl)*porosity_scale*volfrac_scale
!       else
!         rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
!           rt_auxvars(ghosted_id)%mnrl_area0(imnrl)
        endif

        rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
            rt_auxvars(ghosted_id)%mnrl_area0(imnrl)*porosity_scale*volfrac_scale

        if (reaction%update_armor_mineral_surface .and. &
            reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl) > 0.d0) then
          imnrl_armor = imnrl
          do imnrl1 = 1, reaction%mineral%nkinmnrl
            if (reaction%mineral%kinmnrl_armor_min_names(imnrl) == &
                reaction%mineral%kinmnrl_names(imnrl1)) then
              imnrl_armor = imnrl1
              exit
            endif
          enddo

!         print *,'update-armor: ',imnrl,imnrl_armor, &
!         reaction%mineral%kinmnrl_armor_min_names(imnrl_armor)

!       check for negative surface area armoring correction
          if (reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl) > &
              rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl_armor)) then

            if (reaction%update_armor_mineral_surface_flag == 0) then ! surface unarmored
              rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
                rt_auxvars(ghosted_id)%mnrl_area(imnrl) * &
                ((reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl) &
                - rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl_armor))/ &
                reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl))** &
                reaction%mineral%kinmnrl_surf_area_vol_frac_pwr(imnrl)
            else
              rt_auxvars(ghosted_id)%mnrl_area(imnrl) = rt_auxvars(ghosted_id)%mnrl_area0(imnrl)
              reaction%update_armor_mineral_surface_flag = 0
            endif
          else
            rt_auxvars(ghosted_id)%mnrl_area(imnrl) = 0.d0
            reaction%update_armor_mineral_surface_flag = 1 ! surface armored
          endif
        endif

!       print *,'update min srf: ',imnrl,local_id,reaction%mineral%kinmnrl_names(imnrl), &
!       reaction%mineral%kinmnrl_armor_min_names(imnrl), &
!       reaction%update_armor_mineral_surface, &
!       rt_auxvars(ghosted_id)%mnrl_area(imnrl), &
!       reaction%mineral%kinmnrl_armor_pwr(imnrl), &
!       reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl), &
!       rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl_armor), &
!       rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl)
      enddo
    enddo

    if (reaction%update_mnrl_surf_with_porosity) then
      call GridVecRestoreArrayF90(grid,field%work,vec_p,ierr)
    endif

    call DiscretizationLocalToLocal(discretization,field%tortuosity_loc, &
                                     field%tortuosity_loc,ONEDOF)
  endif
      
  if (reaction%update_tortuosity) then
    call GridVecGetArrayF90(grid,field%tortuosity_loc,tortuosity_loc_p,ierr)  
    call GridVecGetArrayF90(grid,field%tortuosity0,tortuosity0_p,ierr)  
    call GridVecGetArrayF90(grid,field%work,vec_p,ierr)
    do local_id = 1, grid%nlmax
      ghosted_id = grid%nL2G(local_id)
      scale = vec_p(local_id)** &
        material_property_array(patch%imat(ghosted_id))%ptr%tortuosity_pwr
      tortuosity_loc_p(ghosted_id) = tortuosity0_p(local_id)*scale
    enddo
    call GridVecRestoreArrayF90(grid,field%tortuosity_loc,tortuosity_loc_p,ierr)  
    call GridVecRestoreArrayF90(grid,field%tortuosity0,tortuosity0_p,ierr)  
    call GridVecRestoreArrayF90(grid,field%work,vec_p,ierr)

    call DiscretizationLocalToLocal(discretization,field%tortuosity_loc, &
                                     field%tortuosity_loc,ONEDOF)
  endif
      
  if (reaction%update_permeability) then
    call GridVecGetArrayF90(grid,field%porosity0,porosity0_p,ierr)
    call GridVecGetArrayF90(grid,field%porosity_loc,porosity_loc_p,ierr)
    call GridVecGetArrayF90(grid,field%perm0_xx,perm0_xx_p,ierr)
    call GridVecGetArrayF90(grid,field%perm0_zz,perm0_zz_p,ierr)
    call GridVecGetArrayF90(grid,field%perm0_yy,perm0_yy_p,ierr)
    call GridVecGetArrayF90(grid,field%perm_xx_loc,perm_xx_loc_p,ierr)
    call GridVecGetArrayF90(grid,field%perm_zz_loc,perm_zz_loc_p,ierr)
    call GridVecGetArrayF90(grid,field%perm_yy_loc,perm_yy_loc_p,ierr)
!   call GridVecGetArrayF90(grid,field%work,vec_p,ierr)
    do local_id = 1, grid%nlmax
      ghosted_id = grid%nL2G(local_id)
      if (porosity_loc_p(ghosted_id) >= material_property_array(patch%imat(ghosted_id))%ptr%permeability_crit_por) then
        scale = ((porosity_loc_p(ghosted_id)-material_property_array(patch%imat(ghosted_id))%ptr%permeability_crit_por) &
        /(porosity0_p(local_id)-material_property_array(patch%imat(ghosted_id))%ptr%permeability_crit_por))** &
        material_property_array(patch%imat(ghosted_id))%ptr%permeability_pwr

#ifdef PERM
        scale = scale*((1.001-porosity0_p(local_id)**2)/(1.001-porosity_loc_p(ghosted_id)**2))
#endif

        if (scale < material_property_array(patch%imat(ghosted_id))%ptr%permeability_min_scale_fac) &
          scale = material_property_array(patch%imat(ghosted_id))%ptr%permeability_min_scale_fac
      else
        scale = material_property_array(patch%imat(ghosted_id))%ptr%permeability_min_scale_fac
      endif
!     scale = vec_p(local_id)** &
!             material_property_array(patch%imat(ghosted_id))%ptr%permeability_pwr
      perm_xx_loc_p(ghosted_id) = perm0_xx_p(local_id)*scale
      perm_yy_loc_p(ghosted_id) = perm0_yy_p(local_id)*scale
      perm_zz_loc_p(ghosted_id) = perm0_zz_p(local_id)*scale
    enddo
    call GridVecRestoreArrayF90(grid,field%porosity0,porosity0_p,ierr)
    call GridVecRestoreArrayF90(grid,field%porosity_loc,porosity_loc_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm0_xx,perm0_xx_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm0_zz,perm0_zz_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm0_yy,perm0_yy_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm_xx_loc,perm_xx_loc_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm_zz_loc,perm_zz_loc_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm_yy_loc,perm_yy_loc_p,ierr)
!   call GridVecRestoreArrayF90(grid,field%work,vec_p,ierr)

    call DiscretizationLocalToLocal(discretization,field%perm_xx_loc, &
                                    field%perm_xx_loc,ONEDOF)
    call DiscretizationLocalToLocal(discretization,field%perm_yy_loc, &
                                    field%perm_yy_loc,ONEDOF)
    call DiscretizationLocalToLocal(discretization,field%perm_zz_loc, &
                                    field%perm_zz_loc,ONEDOF)
  endif  
  
end subroutine RealizationUpdatePropertiesPatch

! ************************************************************************** !
!
! RealLocalToLocalWithArray: Takes an F90 array that is ghosted
!                            and updates the ghosted values
! author: Glenn Hammond
! date: 06/09/11
!
! ************************************************************************** !
subroutine RealLocalToLocalWithArray(realization,array_id)

  use Grid_module

  implicit none

  type(realization_type) :: realization
  PetscInt :: array_id
  
  type(patch_type), pointer :: patch
  type(grid_type), pointer :: grid
  type(field_type), pointer :: field

  field => realization%field
  patch => realization%patch

  grid => patch%grid
  select case(array_id)
    case(MATERIAL_ID_ARRAY)
      call GridCopyIntegerArrayToVec(grid,patch%imat,field%work_loc, &
                                     grid%ngmax)
    case(SATURATION_FUNCTION_ID_ARRAY)
      call GridCopyIntegerArrayToVec(grid,patch%sat_func_id, &
                                     field%work_loc, grid%ngmax)
  end select

  call DiscretizationLocalToLocal(realization%discretization,field%work_loc, &
                                  field%work_loc,ONEDOF)

  select case(array_id)
    case(MATERIAL_ID_ARRAY)
      call GridCopyVecToIntegerArray(grid,patch%imat,field%work_loc, &
                                      grid%ngmax)
    case(SATURATION_FUNCTION_ID_ARRAY)
      call GridCopyVecToIntegerArray(grid,patch%sat_func_id, &
                                      field%work_loc, grid%ngmax)
  end select

end subroutine RealLocalToLocalWithArray

! ************************************************************************** !
!
! RealizationCountCells: Counts # of active and inactive grid cells 
! author: Glenn Hammond
! date: 06/01/10
!
! ************************************************************************** !
subroutine RealizationCountCells(realization,global_total_count, &
                                 global_active_count,total_count,active_count)

  use Option_module

  implicit none
  
  type(realization_type) :: realization
  PetscInt :: global_total_count
  PetscInt :: global_active_count
  PetscInt :: total_count
  PetscInt :: active_count
  
  PetscInt :: patch_total_count
  PetscInt :: patch_active_count
  PetscInt :: temp_int_in(2), temp_int_out(2)
  PetscErrorCode :: ierr
  
  type(patch_type), pointer :: patch
  
  total_count = 0
  active_count = 0
    
  patch => realization%patch
  call PatchCountCells(patch,patch_total_count,patch_active_count)
  total_count = total_count + patch_total_count
  active_count = active_count + patch_active_count
  
  temp_int_in(1) = total_count
  temp_int_in(2) = active_count
  call MPI_Allreduce(temp_int_in,temp_int_out,TWO_INTEGER_MPI,MPIU_INTEGER, &
                     MPI_SUM,realization%option%mycomm,ierr)
  global_total_count = temp_int_out(1)
  global_active_count = temp_int_out(2)

end subroutine RealizationCountCells

! ************************************************************************** !
! ************************************************************************** !
subroutine RealizationSetUpBC4Faces(realization)

  use Connection_module
  use Coupler_module
  use Patch_module
  use Grid_module
  use Field_module
  use MFD_Aux_module
  

#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"

  type(realization_type) :: realization

#ifdef DASVYAT
  type(grid_type), pointer :: grid
  type(patch_type), pointer :: patch
  type(field_type), pointer :: field
  

  type(mfd_auxvar_type), pointer :: aux_var
  type(connection_set_type), pointer :: conn
  type(coupler_type), pointer ::  boundary_condition

  PetscReal, pointer :: bc_faces_p(:), xx_faces_p(:)
  PetscInt :: iconn, sum_connection, bc_type, bound_id
  PetscInt :: local_id, ghosted_id, ghost_face_id, j, jface, local_face_id
  PetscErrorCode :: ierr

  patch => realization%patch
  grid => patch%grid
  field => realization%field

  call VecGetArrayF90(field%flow_bc_loc_faces, bc_faces_p, ierr)
  call VecGetArrayF90(field%flow_xx_faces, xx_faces_p, ierr)

  boundary_condition => patch%boundary_conditions%first
  sum_connection = 0
  do
    if (.not.associated(boundary_condition)) exit
    bc_type = boundary_condition%flow_condition%itype(RICHARDS_PRESSURE_DOF)

    do iconn = 1, boundary_condition%numfaces_set
      sum_connection = sum_connection + 1

      local_id = boundary_condition%region%cell_ids(iconn)
      ghosted_id = grid%nL2G(local_id)

      aux_var => grid%MFD%aux_vars(local_id)
      do j = 1, aux_var%numfaces
        ghost_face_id = aux_var%face_id_gh(j)
        local_face_id = grid%fG2L(ghost_face_id)
        conn => grid%faces(ghost_face_id)%conn_set_ptr
        jface = grid%faces(ghost_face_id)%id
        if (boundary_condition%faces_set(iconn) == ghost_face_id) then
          if ((bc_type == DIRICHLET_BC).or.(bc_type == HYDROSTATIC_BC)  &
              .or.(bc_type == SEEPAGE_BC).or.(bc_type == CONDUCTANCE_BC) ) then
            bc_faces_p(ghost_face_id) = boundary_condition%flow_aux_real_var(1,iconn)*conn%area(jface)
            xx_faces_p(local_face_id) = boundary_condition%flow_aux_real_var(1,iconn)
          else if ((bc_type == NEUMANN_BC)) then
            bc_faces_p(ghost_face_id) = boundary_condition%flow_aux_real_var(1,iconn)*conn%area(jface)
            bound_id = grid%fL2B(local_face_id)
            if (bound_id>0) then
              patch%boundary_velocities(realization%option%nphase, bound_id) = &
                boundary_condition%flow_aux_real_var(1,iconn)
            endif
          endif
        endif
      enddo
    enddo
    boundary_condition => boundary_condition%next
  enddo

  call VecRestoreArrayF90(field%flow_xx_faces, xx_faces_p, ierr)
  call VecRestoreArrayF90(field%flow_bc_loc_faces, bc_faces_p, ierr)

#endif

end subroutine RealizationSetUpBC4Faces

! ************************************************************************** !
!
! RealizationPrintGridStatistics: Prints statistics regarding the numerical
!                                 discretization 
! author: Glenn Hammond
! date: 06/01/10
!
! ************************************************************************** !
subroutine RealizationPrintGridStatistics(realization)

  use Grid_module

  implicit none

  type(realization_type) :: realization
  
  type(option_type), pointer :: option
  type(grid_type), pointer :: grid

  PetscInt :: i1, i2, i3
  PetscReal :: r1, r2, r3
  PetscInt :: global_total_count, global_active_count
  PetscInt :: total_count, active_count
  PetscReal :: total_min, total_max, total_mean, total_variance
  PetscReal :: active_min, active_max, active_mean, active_variance
  PetscInt :: inactive_histogram(12), temp_int_out(12)
  PetscReal :: inactive_percentages(12)
  PetscErrorCode :: ierr

  option => realization%option
  grid => realization%patch%grid

  ! print # of active and inactive grid cells
  call RealizationCountCells(realization,global_total_count, &
                             global_active_count,total_count,active_count)
  r1 = dble(total_count)
  call OptionMaxMinMeanVariance(r1,total_max, &
                                total_min,total_mean, &
                                total_variance,PETSC_TRUE,option)
  r1 = dble(active_count)
  call OptionMaxMinMeanVariance(r1,active_max, &
                                active_min,active_mean, &
                                active_variance,PETSC_TRUE,option)
                  
  r1 = dble(active_count) / dble(total_count)    
  inactive_histogram = 0                          
  if (r1 >= (1.d0-1.d-8)) then
    inactive_histogram(12) = 1
  else if (r1 >= .9d0 .and. r1 < (1.d0-1.d-8)) then
    inactive_histogram(11) = 1
  else if (r1 >= .8d0 .and. r1 < .9d0) then
    inactive_histogram(10) = 1
  else if (r1 >= .7d0 .and. r1 < .8d0) then
    inactive_histogram(9) = 1
  else if (r1 >= .6d0 .and. r1 < .7d0) then
    inactive_histogram(8) = 1
  else if (r1 >= .5d0 .and. r1 < .6d0) then
    inactive_histogram(7) = 1
  else if (r1 >= .4d0 .and. r1 < .5d0) then
    inactive_histogram(6) = 1
  else if (r1 >= .3d0 .and. r1 < .4d0) then
    inactive_histogram(5) = 1
  else if (r1 >= .2d0 .and. r1 < .3d0) then
    inactive_histogram(4) = 1
  else if (r1 >= .1d0 .and. r1 < .2d0) then
    inactive_histogram(3) = 1
  else if (r1 > 1.d-20 .and. r1 < .1d0) then
    inactive_histogram(2) = 1
  else if (r1 < 1.d-20) then
    inactive_histogram(1) = 1
  endif
  
  call MPI_Allreduce(inactive_histogram,temp_int_out,TWELVE_INTEGER_MPI, &
                     MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)

  ! why I cannot use *100, I do not know....geh
  inactive_percentages = dble(temp_int_out)/dble(option%mycommsize)*10.d0
  inactive_percentages = inactive_percentages+1.d-8

  r1 = 0.d0
  do i1 = 1, 12
    r1 = r1 + inactive_percentages(i1)
  enddo
                                
  i1 = -999
  i2 = -999
  i3 = -999
  if (associated(grid%structured_grid)) then
    i1 = grid%structured_grid%npx_final
    i2 = grid%structured_grid%npy_final
    i3 = grid%structured_grid%npz_final
  endif
  if (OptionPrintToScreen(option)) then
    write(*,'(/," Grid Stats:",/, &
              & "                       Global # cells: ",i12,/, &
              & "                Global # active cells: ",i12,/, &
              & "                              # cores: ",i12,/, &
              & "         Processor core decomposition: ",3i6,/, &
              & "               Maximum # cells / core: ",i12,/, &
              & "               Minimum # cells / core: ",i12,/, &
              & "               Average # cells / core: ",1pe12.4,/, &
              & "               Std Dev # cells / core: ",1pe12.4,/, &
              & "        Maximum # active cells / core: ",i12,/, &
              & "        Minimum # active cells / core: ",i12,/, &
              & "        Average # active cells / core: ",1pe12.4,/, &
              & "        Std Dev # active cells / core: ",1pe12.4,/,/, &
              & "        % cores with % active cells =       0%: ",1f7.2,/, &
              & "        % cores with % active cells =  0.1-10%: ",1f7.2,/, &
              & "        % cores with % active cells =   10-20%: ",1f7.2,/, &
              & "        % cores with % active cells =   20-30%: ",1f7.2,/, &
              & "        % cores with % active cells =   30-40%: ",1f7.2,/, &
              & "        % cores with % active cells =   40-50%: ",1f7.2,/, &
              & "        % cores with % active cells =   50-60%: ",1f7.2,/, &
              & "        % cores with % active cells =   60-70%: ",1f7.2,/, &
              & "        % cores with % active cells =   70-80%: ",1f7.2,/, &
              & "        % cores with % active cells =   80-90%: ",1f7.2,/, &
              & "        % cores with % active cells = 90-99.9%: ",1f7.2,/, &
              & "        % cores with % active cells =     100%: ",1f7.2,/, &
              & "                                        Check : ",1f7.2,/)') &
           global_total_count, &
           global_active_count, &
           option%mycommsize, &
           i1,i2,i3, &
           int(total_max+1.d-4), &
           int(total_min+1.d-4), &
           total_mean, sqrt(total_variance), &
           int(active_max+1.d-4), &
           int(active_min+1.d-4), &
           active_mean, sqrt(active_variance), &
           inactive_percentages(1), &
           inactive_percentages(2), &
           inactive_percentages(3), &
           inactive_percentages(4), &
           inactive_percentages(5), &
           inactive_percentages(6), &
           inactive_percentages(7), &
           inactive_percentages(8), &
           inactive_percentages(9), &
           inactive_percentages(10), &
           inactive_percentages(11), &
           inactive_percentages(12), &
           r1
  endif
  if (OptionPrintToFile(option)) then
    write(option%fid_out,'(/," Grid Stats:",/, &
               & "                       Global # cells: ",i12,/, &
               & "                Global # active cells: ",i12,/, &
               & "                              # cores: ",i12,/, &
               & "         Processor core decomposition: ",3i6,/, &
               & "               Maximum # cells / core: ",i12,/, &
               & "               Minimum # cells / core: ",i12,/, &
               & "               Average # cells / core: ",1pe12.4,/, &
               & "               Std Dev # cells / core: ",1pe12.4,/, &
               & "        Maximum # active cells / core: ",i12,/, &
               & "        Minimum # active cells / core: ",i12,/, &
               & "        Average # active cells / core: ",1pe12.4,/, &
               & "        Std Dev # active cells / core: ",1pe12.4,/,/, &
               & "        % cores with % active cells =       0%: ",1f7.2,/, &
               & "        % cores with % active cells =  0.1-10%: ",1f7.2,/, &
               & "        % cores with % active cells =   10-20%: ",1f7.2,/, &
               & "        % cores with % active cells =   20-30%: ",1f7.2,/, &
               & "        % cores with % active cells =   30-40%: ",1f7.2,/, &
               & "        % cores with % active cells =   40-50%: ",1f7.2,/, &
               & "        % cores with % active cells =   50-60%: ",1f7.2,/, &
               & "        % cores with % active cells =   60-70%: ",1f7.2,/, &
               & "        % cores with % active cells =   70-80%: ",1f7.2,/, &
               & "        % cores with % active cells =   80-90%: ",1f7.2,/, &
               & "        % cores with % active cells = 90-99.9%: ",1f7.2,/, &
               & "        % cores with % active cells =     100%: ",1f7.2,/, &
               & "                                        Check : ",1f7.2,/)') &
           global_total_count, &
           global_active_count, &
           option%mycommsize, &
           i1,i2,i3, &
           int(total_max+1.d-4), &
           int(total_min+1.d-4), &
           total_mean, sqrt(total_variance), &
           int(active_max+1.d-4), &
           int(active_min+1.d-4), &
           active_mean, sqrt(active_variance), &
           inactive_percentages(1), &
           inactive_percentages(2), &
           inactive_percentages(3), &
           inactive_percentages(4), &
           inactive_percentages(5), &
           inactive_percentages(6), &
           inactive_percentages(7), &
           inactive_percentages(8), &
           inactive_percentages(9), &
           inactive_percentages(10), &
           inactive_percentages(11), &
           inactive_percentages(12), &
           r1
  endif

end subroutine RealizationPrintGridStatistics

! ************************************************************************** !
!
! RealizationCalculateCFL1Timestep: Calculates largest time step that  
!                                   preserves a CFL # of 1 in a realization
! author: Glenn Hammond
! date: 10/07/11
!
! ************************************************************************** !
subroutine RealizationCalculateCFL1Timestep(realization,max_dt_cfl_1)

  implicit none

  type(realization_type) realization
  PetscReal :: max_dt_cfl_1
  
  type(patch_type), pointer :: patch
  PetscReal :: max_dt_cfl_1_patch
  PetscReal :: tempreal
  PetscErrorCode :: ierr
  
  max_dt_cfl_1 = 1.d20
  patch => realization%patch
  call PatchCalculateCFL1Timestep(patch,realization%option, &
                                  max_dt_cfl_1_patch)
  max_dt_cfl_1 = min(max_dt_cfl_1,max_dt_cfl_1_patch)

  ! get the minimum across all cores
  call MPI_Allreduce(max_dt_cfl_1,tempreal,ONE_INTEGER_MPI, &
                     MPI_DOUBLE_PRECISION,MPI_MIN, &
                     realization%option%mycomm,ierr)
  max_dt_cfl_1 = tempreal

end subroutine RealizationCalculateCFL1Timestep

! ************************************************************************** !
!
! RealizationNonInitializedData: Checks for non-initialized data sets
!                                i.e. porosity, permeability
! author: Glenn Hammond
! date: 02/08/13
!
! ************************************************************************** !
subroutine RealizationNonInitializedData(realization)

  use Grid_module
  use Patch_module
  use Option_module
  use Field_module

  implicit none
  
  type(realization_type) :: realization
  
  type(patch_type), pointer :: patch
  type(grid_type), pointer :: grid
  type(option_type), pointer :: option
  type(field_type), pointer :: field  
  PetscReal, pointer :: vec_p(:), vecy_p(:), vecz_p(:)
  PetscReal :: min_value, global_value
  PetscInt :: local_id
  PetscErrorCode :: ierr
  
  patch => realization%patch
  grid => patch%grid
  option => realization%option
  field => realization%field
  
  ! cannot use VecMin as there may be inactive cells without data assigned.
  
  min_value = 1.d20
  ! porosity
  call GridVecGetArrayF90(grid,field%porosity0,vec_p,ierr)
  do local_id = 1, grid%nlmax
    if (patch%imat(grid%nL2G(local_id)) <= 0) cycle
    min_value = min(min_value,vec_p(local_id)) 
  enddo
  call GridVecRestoreArrayF90(grid,field%porosity0,vec_p,ierr)
  call MPI_Allreduce(min_value,global_value,ONE_INTEGER_MPI, &
                     MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
  
  if (global_value < -998.d0) then
    option%io_buffer = 'Porosity not initialized at all cells.  ' // &
                       'Ensure that REGIONS cover entire domain!!!'
    call printErrMsg(option)
  endif
  
  if (option%iflowmode /= NULL_MODE) then
    call GridVecGetArrayF90(grid,field%perm0_xx,vec_p,ierr)
    call GridVecGetArrayF90(grid,field%perm0_yy,vecy_p,ierr)
    call GridVecGetArrayF90(grid,field%perm0_zz,vecz_p,ierr)
    min_value = 1.d20
    do local_id = 1, grid%nlmax
      if (patch%imat(grid%nL2G(local_id)) <= 0) cycle
      min_value = min(min_value,vec_p(local_id),vecy_p(local_id), &
                      vecz_p(local_id)) 
    enddo        
    call GridVecRestoreArrayF90(grid,field%perm0_xx,vec_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm0_yy,vecy_p,ierr)
    call GridVecRestoreArrayF90(grid,field%perm0_zz,vecz_p,ierr)
    call MPI_Allreduce(min_value,global_value,ONE_INTEGER_MPI, &
                       MPI_DOUBLE_PRECISION,MPI_MIN,option%mycomm,ierr)
    if (global_value < 1.d-60) then
      option%io_buffer = &
        'A positive non-zero permeability must be defined throughout ' // &
        'domain in X, Y and Z.'
      call printErrMsg(option)
    endif
  endif

end subroutine RealizationNonInitializedData

! ************************************************************************** !
!
! RealizationDestroy: Deallocates a realization
! author: Glenn Hammond
! date: 11/01/07
!
! ************************************************************************** !
subroutine RealizationDestroy(realization)

  implicit none
  
  type(realization_type), pointer :: realization
  
  if (.not.associated(realization)) return
    
  call FieldDestroy(realization%field)

!  call OptionDestroy(realization%option) !geh it will be destroy externally
  call OutputOptionDestroy(realization%output_option)
  call RegionDestroyList(realization%regions)
  
  call FlowConditionDestroyList(realization%flow_conditions)
  call TranConditionDestroyList(realization%transport_conditions)
  call TranConstraintDestroyList(realization%transport_constraints)

  call LevelDestroyList(realization%level_list)

  if (associated(realization%debug)) deallocate(realization%debug)
  nullify(realization%debug)
  
  if (associated(realization%fluid_property_array)) &
    deallocate(realization%fluid_property_array)
  nullify(realization%fluid_property_array)
  call FluidPropertyDestroy(realization%fluid_properties)
  
  if (associated(realization%material_property_array)) &
    deallocate(realization%material_property_array)
  nullify(realization%material_property_array)
  call MaterialPropertyDestroy(realization%material_properties)

  if (associated(realization%saturation_function_array)) &
    deallocate(realization%saturation_function_array)
  nullify(realization%saturation_function_array)
  call SaturationFunctionDestroy(realization%saturation_functions)

  call DatasetDestroy(realization%datasets)
  
  call UniformVelocityDatasetDestroy(realization%uniform_velocity_dataset)
  
  call DiscretizationDestroy(realization%discretization)
  
  call ReactionDestroy(realization%reaction)
  
  call TranConstraintDestroy(realization%sec_transport_constraint)  
  
end subroutine RealizationDestroy

end module Realization_class