Source

pflotran-dev / src / pflotran / discretization.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
module Discretization_module

  use Grid_module
  use Structured_Grid_module
  use Unstructured_Grid_module
  use Unstructured_Grid_Aux_module
  use Unstructured_Explicit_module
  use MFD_Aux_module
  use MFD_Module

  implicit none

  private
 
#include "definitions.h"

#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/petscdmda.h"
#include "finclude/petscdmshell.h90"

  type, public :: dm_ptr_type
    DM :: dm  ! PETSc DM
    type(ugdm_type), pointer :: ugdm
      ! Unstructured grid "private" dm.  This gets wrapped in a PETSc DM via 
      ! DMShell routines.
  end type dm_ptr_type

  type, public :: discretization_type
    PetscInt :: itype  ! type of discretization (e.g. structured, unstructured, etc.)
    !geh: note that differentiating between implicit and explicit unstructured 
    !     grids is handled within the grid%itype variable, not discritization%itype
    character(len=MAXWORDLENGTH) :: ctype
    PetscReal :: origin(3) ! origin of global domain
    type(grid_type), pointer :: grid  ! pointer to a grid object
    character(len=MAXSTRINGLENGTH) :: filename

    type(dm_ptr_type), pointer :: dmc_nflowdof(:), dmc_ntrandof(:)
      ! Arrays containing hierarchy of coarsened DMs, for use with Galerkin 
      ! multigrid.  Element i of each array is a *finer* DM than element i-1.
    PetscInt :: dm_index_to_ndof(5) ! mapping between a dm_ptr to the number of degrees of freedom
    type(dm_ptr_type), pointer :: dm_1dof
    type(dm_ptr_type), pointer :: dm_nflowdof
    type(dm_ptr_type), pointer :: dm_ntrandof
    type(mfd_type), pointer :: MFD
    VecScatter :: tvd_ghost_scatter
    
    PetscInt :: stencil_width
    PetscInt :: stencil_type
    
    PetscInt :: hydr_flux_method
    PetscInt :: therm_flux_method
    PetscInt :: mech_flux_method
    PetscBool:: lsm_flux_method
    
  end type discretization_type

  public :: DiscretizationCreate, &
            DiscretizationDestroy, &
            DiscretizationReadRequiredCards, &
            DiscretizationRead, &
            DiscretizationCreateVector, &
            DiscretizationDuplicateVector, &         
            DiscretizationCreateJacobian, &
            DiscretizationCreateInterpolation, &
            DiscretizationCreateColoring, &
            DiscretizationGlobalToLocal, &
            DiscretizationGlobalToLocalFaces, &
            DiscretizationLocalToGlobal, &
            DiscretizationLocalToLocal, &
            DiscretizationLocalToLocalFaces, &
            DiscretizationGlobalToNatural, &
            DiscretizationNaturalToGlobal, &
            DiscretizationGlobalToLocalBegin, &
            DiscretizationGlobalToLocalEnd, &
            DiscretizGlobalToLocalFacesBegin, &
            DiscretizGlobalToLocalFacesEnd, &
            DiscretizLocalToLocalFacesBegin, &
            DiscretizLocalToLocalFacesEnd, &
            DiscretizationGlobalToLocalLP, &
            DiscretizationLocalToLocalLP, &
            DiscretizGlobalToLocalLPBegin, &
            DiscretizGlobalToLocalLPEnd, &
            DiscretizLocalToLocalLPBegin, &
            DiscretizLocalToLocalLPEnd, &
            DiscretizationLocalToLocalBegin, &
            DiscretizationLocalToLocalEnd, &
            DiscretizGlobalToNaturalBegin, &
            DiscretizGlobalToNaturalEnd, &
            DiscretizNaturalToGlobalBegin, &
            DiscretizNaturalToGlobalEnd, &
            DiscretizationCreateDMs,&
            DiscretizationGetDMPtrFromIndex, &
            DiscretizationUpdateTVDGhosts, &
            DiscretAOApplicationToPetsc
  
contains

! ************************************************************************** !
!
! DiscretizationCreate: Creates a structured or unstructured discretization
! author: Glenn Hammond
! date: 10/23/07
!
! ************************************************************************** !
function DiscretizationCreate()

  implicit none
  
  type(discretization_type), pointer :: DiscretizationCreate
  
  type(discretization_type), pointer :: discretization
  
  allocate(discretization)
  discretization%ctype = ''
  discretization%itype = 0
  discretization%origin = 0.d0
  discretization%filename = ''

  ! nullify DM pointers
  nullify(discretization%dmc_nflowdof)
  nullify(discretization%dmc_ntrandof)
  allocate(discretization%dm_1dof)
  allocate(discretization%dm_nflowdof)
  allocate(discretization%dm_ntrandof)
  discretization%dm_1dof%dm = 0
  discretization%dm_nflowdof%dm = 0
  discretization%dm_ntrandof%dm = 0
  nullify(discretization%dm_1dof%ugdm)
  nullify(discretization%dm_nflowdof%ugdm)
  nullify(discretization%dm_ntrandof%ugdm)
  
  nullify(discretization%grid)
  nullify(discretization%MFD)
  
  discretization%stencil_width = 1
  discretization%stencil_type = DMDA_STENCIL_STAR
  discretization%hydr_flux_method = TWO_POINT_FLUX
  discretization%therm_flux_method = TWO_POINT_FLUX
  discretization%mech_flux_method = TWO_POINT_FLUX
  discretization%lsm_flux_method = PETSC_FALSE

  discretization%tvd_ghost_scatter = 0
  
  DiscretizationCreate => discretization

end function DiscretizationCreate

! ************************************************************************** !
!
! DiscretizationReadRequiredCards: Reads a discretization from the input file
! author: Glenn Hammond
! date: 11/01/07
!
! ************************************************************************** !
subroutine DiscretizationReadRequiredCards(discretization,input,option)

  use Option_module
  use Input_module
  use String_module

  implicit none

  type(option_type), pointer :: option
  type(input_type), pointer :: input
  type(discretization_type),pointer :: discretization
  character(len=MAXWORDLENGTH) :: word
  type(grid_type), pointer :: grid, grid2
  type(structured_grid_type), pointer :: str_grid
  type(unstructured_grid_type), pointer :: un_str_grid
  character(len=MAXWORDLENGTH) :: structured_grid_ctype
  character(len=MAXWORDLENGTH) :: unstructured_grid_ctype

  character(len=MAXSTRINGLENGTH) :: string

  PetscInt :: structured_grid_itype
  PetscInt :: unstructured_grid_itype
  PetscInt :: nx, ny, nz
  PetscInt :: i
  PetscReal :: tempreal

  nx = 0
  ny = 0
  nz = 0

! we initialize the word to blanks to avoid error reported by valgrind
  word = ''

  do
  
    call InputReadFlotranString(input,option)
    if (input%ierr /= 0) exit

    if (InputCheckExit(input,option)) exit

    call InputReadWord(input,option,word,PETSC_TRUE)
    call InputErrorMsg(input,option,'keyword','GRID')
    call StringToUpper(word)
      
    select case(trim(word))
      case('TYPE')
        call InputReadWord(input,option,discretization%ctype,PETSC_TRUE)
        call InputErrorMsg(input,option,'type','GRID')   
        call StringToLower(discretization%ctype)
        select case(trim(discretization%ctype))
          case('structured')
            discretization%itype = STRUCTURED_GRID
            call InputReadWord(input,option,structured_grid_ctype,PETSC_TRUE)
            call InputDefaultMsg(input,option,'structured_grid_type') 
            call StringToLower(structured_grid_ctype)
            select case(trim(structured_grid_ctype))
              case('cartesian')
                structured_grid_itype = CARTESIAN_GRID
              case('cylindrical')
                structured_grid_itype = CYLINDRICAL_GRID
              case('spherical')
                structured_grid_itype = SPHERICAL_GRID
              case default
                structured_grid_itype = CARTESIAN_GRID
                structured_grid_ctype = 'cartesian'
            end select
          case('unstructured','unstructured_explicit')
            discretization%itype = UNSTRUCTURED_GRID
            word = discretization%ctype
            discretization%ctype = 'unstructured'
            select case(word)
              case('unstructured')
                unstructured_grid_itype = IMPLICIT_UNSTRUCTURED_GRID
                unstructured_grid_ctype = 'implicit unstructured'
              case('unstructured_explicit')
                unstructured_grid_itype = EXPLICIT_UNSTRUCTURED_GRID
                unstructured_grid_ctype = 'explicit unstructured'
            end select
            call InputReadNChars(input,option,discretization%filename,MAXSTRINGLENGTH, &
                                 PETSC_TRUE)
            call InputErrorMsg(input,option,'unstructured filename','GRID')
          case('unstructured_mimetic')
            discretization%itype = UNSTRUCTURED_GRID_MIMETIC
            option%mimetic = PETSC_TRUE
            word = discretization%ctype
            discretization%ctype = 'unstructured_mimetic'
            unstructured_grid_itype = IMPLICIT_UNSTRUCTURED_GRID
            unstructured_grid_ctype = 'implicit unstructured'
            call InputReadNChars(input,option,discretization%filename,MAXSTRINGLENGTH, &
                                 PETSC_TRUE)
            call InputErrorMsg(input,option,'unstructured filename','GRID')
          case('structured_mimetic')
            discretization%itype = STRUCTURED_GRID_MIMETIC
            option%mimetic = PETSC_TRUE
            call InputReadWord(input,option,structured_grid_ctype,PETSC_TRUE)
            call InputDefaultMsg(input,option,'structured_grid_type')   
            call StringToLower(structured_grid_ctype)
            select case(trim(structured_grid_ctype))
              case('cartesian')
                structured_grid_itype = CARTESIAN_GRID
              case('cylindrical')
                structured_grid_itype = CYLINDRICAL_GRID
              case('spherical')
                structured_grid_itype = SPHERICAL_GRID
              case default
                structured_grid_itype = CARTESIAN_GRID
                structured_grid_ctype = 'cartesian'
            end select
          case default
            option%io_buffer = 'Discretization type: ' // &
                               trim(discretization%ctype) // &
                               ' not recognized.'
            call printErrMsg(option)
        end select    
      case('NXYZ')
        call InputReadInt(input,option,nx)
        call InputErrorMsg(input,option,'nx','GRID')
        call InputReadInt(input,option,ny)
        call InputErrorMsg(input,option,'ny','GRID')
        call InputReadInt(input,option,nz)
        call InputErrorMsg(input,option,'nz','GRID')
        if (structured_grid_itype /= CARTESIAN_GRID) then
          ny = 1 ! cylindrical and spherical have 1 cell in Y
          if (structured_grid_itype /= CYLINDRICAL_GRID) nz = 1 ! spherical has 1 cell in Z
        endif
      case('ORIGIN')
        call InputReadDouble(input,option,discretization%origin(X_DIRECTION))
        call InputErrorMsg(input,option,'X direction','Origin')
        call InputReadDouble(input,option,discretization%origin(Y_DIRECTION))
        call InputErrorMsg(input,option,'Y direction','Origin')
        call InputReadDouble(input,option,discretization%origin(Z_DIRECTION))
        call InputErrorMsg(input,option,'Z direction','Origin')        
      case('FILE','GRAVITY','INVERT_Z','MAX_CELLS_SHARING_A_VERTEX',&
           'STENCIL_WIDTH','STENCIL_TYPE','FLUX_METHOD')
      case('DXYZ','BOUNDS')
        call InputSkipToEND(input,option,word) 
      case default
        option%io_buffer = 'Keyword: ' // trim(word) // &
                 ' not recognized in DISCRETIZATION, first read.'
        call printErrMsg(option)          
    end select 
  enddo  

  if (discretization%itype == NULL_GRID) then
    option%io_buffer = 'Discretization type not defined under ' // &
                       'keyword GRID.' 
    call printErrMsg(option)
  endif
  
  grid => GridCreate()
  select case(discretization%itype)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      un_str_grid => UGridCreate()
      select case(unstructured_grid_itype)
        case(IMPLICIT_UNSTRUCTURED_GRID)
          if (index(discretization%filename,'.h5') > 0) then
#if !defined(PETSC_HAVE_HDF5)
            option%io_buffer = 'PFLOTRAN must be built with HDF5 ' // &
              'support to read unstructured grid .h5 files'
            call printErrMsg(option)
#else

#ifdef SCORPIO
            call UGridReadHDF5PIOLib(un_str_grid,discretization%filename,option)
#else
            call UGridReadHDF5(un_str_grid,discretization%filename,option)
#endif
! #ifdef SCORPIO

#endif
!#if !defined(PETSC_HAVE_HDF5)

          else
            call UGridRead(un_str_grid,discretization%filename,option)
          endif
          grid%unstructured_grid => un_str_grid
        case(EXPLICIT_UNSTRUCTURED_GRID)
          un_str_grid%explicit_grid => UGridExplicitCreate()
#if 0          
          call ExplicitUGridRead(un_str_grid%explicit_grid, &
                                 discretization%filename,option)
#else          
          call ExplicitUGridReadInParallel(un_str_grid%explicit_grid, &
                                 discretization%filename,option)
#endif          
          grid%unstructured_grid => un_str_grid
      end select
      grid%itype = unstructured_grid_itype
      grid%ctype = unstructured_grid_ctype
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)      
      if (nx*ny*nz <= 0) &
        call printErrMsg(option,'NXYZ not set correctly for structured grid.')
      str_grid => StructGridCreate()
      str_grid%nx = nx
      str_grid%ny = ny
      str_grid%nz = nz
      str_grid%nxy = str_grid%nx*str_grid%ny
      str_grid%nmax = str_grid%nxy*str_grid%nz
      grid%structured_grid => str_grid
      grid%nmax = str_grid%nmax
      grid%structured_grid%itype = structured_grid_itype
      grid%structured_grid%ctype = structured_grid_ctype
      grid%itype = discretization%itype
      grid%ctype = discretization%ctype
  end select
  grid%discretization_itype=discretization%itype
  discretization%grid => grid

end subroutine DiscretizationReadRequiredCards

! ************************************************************************** !
!
! DiscretizationRead: Reads a discretization from the input file
! author: Glenn Hammond
! date: 11/01/07
!
! ************************************************************************** !
subroutine DiscretizationRead(discretization,input,option)

  use Option_module
  use Input_module
  use String_module

  implicit none

  type(option_type), pointer :: option
  type(input_type), pointer :: input
  type(discretization_type),pointer :: discretization
  character(len=MAXWORDLENGTH) :: word
  type(grid_type), pointer :: grid, grid2
  type(structured_grid_type), pointer :: str_grid
  type(unstructured_grid_type), pointer :: un_str_grid
  character(len=MAXWORDLENGTH) :: structured_grid_ctype
  character(len=MAXSTRINGLENGTH) :: filename
  character(len=MAXSTRINGLENGTH) :: string

  PetscInt :: structured_grid_itype
  PetscInt :: nx, ny, nz
  PetscInt :: i
  PetscReal :: tempreal

  nx = 0
  ny = 0
  nz = 0

! we initialize the word to blanks to avoid error reported by valgrind
  word = ''

  do
  
    call InputReadFlotranString(input,option)
    if (input%ierr /= 0) exit

    if (InputCheckExit(input,option)) exit

    call InputReadWord(input,option,word,PETSC_TRUE)
    call InputErrorMsg(input,option,'keyword','GRID')
    call StringToUpper(word)
      
    select case(trim(word))
      case('TYPE','NXYZ','ORIGIN','FILE')
      case('DXYZ')
        select case(discretization%itype)
          case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
            call StructGridReadDXYZ(discretization%grid%structured_grid,input,option)
          case default
            call printErrMsg(option,'Keyword "DXYZ" not supported for unstructured grid')
        end select
        call InputReadFlotranString(input,option) ! read END card
        call InputReadStringErrorMsg(input,option,'DISCRETIZATION,DXYZ,END')
        if (.not.(InputCheckExit(input,option))) then
          option%io_buffer = 'Card DXYZ should include either 3 entires ' // &
                   '(one for each grid direction or NX+NY+NZ entries)'
          call printErrMsg(option)
        endif
      case('BOUNDS')
        select case(discretization%itype)
          case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
            grid => discretization%grid

            ! read first line and we will split off the legacy approach vs. new
            call InputReadFlotranString(input,option)
            call InputReadStringErrorMsg(input,option,'DISCRETIZATION,BOUNDS,X or Min Coordinate')
            string = input%buf

            do i = 1, 3
              call InputReadDouble(input,option,tempreal)
              if (input%ierr /= 0) exit
            enddo

            input%ierr = 0
            input%buf = string

            if (i == 3) then ! only 2 successfully read
              if (grid%structured_grid%itype == CARTESIAN_GRID .or. &
                  grid%structured_grid%itype == CYLINDRICAL_GRID .or. &
                  grid%structured_grid%itype == SPHERICAL_GRID) then
!geh                  call InputReadFlotranString(input,option) ! x-direction
!geh                  call InputReadStringErrorMsg(input,option,'DISCRETIZATION,BOUNDS,X or R')
                call InputReadDouble(input,option,grid%structured_grid%bounds(X_DIRECTION,LOWER))
                call InputErrorMsg(input,option,'Lower X or R','BOUNDS')
                call InputReadDouble(input,option,grid%structured_grid%bounds(X_DIRECTION,UPPER))
                call InputErrorMsg(input,option,'Upper X or R','BOUNDS')
              endif
              if (grid%structured_grid%itype == CARTESIAN_GRID) then
                call InputReadFlotranString(input,option) ! y-direction
                call InputReadStringErrorMsg(input,option,'DISCRETIZATION,BOUNDS,Y')
                call InputReadDouble(input,option,grid%structured_grid%bounds(Y_DIRECTION,LOWER))
                call InputErrorMsg(input,option,'Lower Y','BOUNDS')
                call InputReadDouble(input,option,grid%structured_grid%bounds(Y_DIRECTION,UPPER))
                call InputErrorMsg(input,option,'Upper Y','BOUNDS')
              else
                grid%structured_grid%bounds(Y_DIRECTION,LOWER) = 0.d0
                grid%structured_grid%bounds(Y_DIRECTION,UPPER) = 1.d0
              endif
              if (grid%structured_grid%itype == CARTESIAN_GRID .or. &
                  grid%structured_grid%itype == CYLINDRICAL_GRID) then
                call InputReadFlotranString(input,option) ! z-direction
                call InputReadStringErrorMsg(input,option,'DISCRETIZATION,BOUNDS,Z')
                call InputReadDouble(input,option,grid%structured_grid%bounds(Z_DIRECTION,LOWER))
                call InputErrorMsg(input,option,'Lower Z','BOUNDS')
                call InputReadDouble(input,option,grid%structured_grid%bounds(Z_DIRECTION,UPPER))
                call InputErrorMsg(input,option,'Upper Z','BOUNDS')
              else
                grid%structured_grid%bounds(Z_DIRECTION,LOWER) = 0.d0
                grid%structured_grid%bounds(Z_DIRECTION,UPPER) = 1.d0
              endif
            else ! new min max coordinate approach
              select case(grid%structured_grid%itype)
                case(CARTESIAN_GRID)
                  i = 3
                case(CYLINDRICAL_GRID)
                  i = 2
                case(SPHERICAL_GRID)
                  i = 1
              end select
              call InputReadNDoubles(input,option, &
                                     grid%structured_grid%bounds(:,LOWER), &
                                     i)
              call InputErrorMsg(input,option,'Minimum Coordinate','BOUNDS')
              call InputReadFlotranString(input,option)
              call InputReadStringErrorMsg(input,option,'DISCRETIZATION,BOUNDS,MAX COORDINATE')
              call InputReadNDoubles(input,option, &
                                     grid%structured_grid%bounds(:,UPPER), &
                                     i)
              call InputErrorMsg(input,option,'Maximum Coordinate','BOUNDS')
              if (grid%structured_grid%itype == CYLINDRICAL_GRID) then
                grid%structured_grid%bounds(Y_DIRECTION,LOWER) = 0.d0
                grid%structured_grid%bounds(Y_DIRECTION,UPPER) = 1.d0
              endif
              if (grid%structured_grid%itype == SPHERICAL_GRID) then
                grid%structured_grid%bounds(Z_DIRECTION,LOWER) = 0.d0
                grid%structured_grid%bounds(Z_DIRECTION,UPPER) = 1.d0
              endif
            endif
            call InputReadFlotranString(input,option)
            call InputReadStringErrorMsg(input,option,'DISCRETIZATION,BOUNDS,END')
            if (.not.(InputCheckExit(input,option))) then
              if (OptionPrintToScreen(option)) then
                if (grid%structured_grid%itype == CARTESIAN_GRID) then
                  print *, 'BOUNDS card for a cartesian structured grid must include ' // &
                           '4 lines.  I.e.'
                  print *, 'BOUNDS'
                  print *, '  x_min  y_min  z_min'
                  print *, '  x_max  y_max  z_max'
                  print *, 'END'
                else if (grid%structured_grid%itype == CYLINDRICAL_GRID) then
                  print *, 'BOUNDS card for a cylindrical structured grid must include ' // &
                           '4 lines.  I.e.'
                  print *, 'BOUNDS'
                  print *, '  r_min  z_min'
                  print *, '  r_max  z_max'
                  print *, 'END'
                else if (grid%structured_grid%itype == SPHERICAL_GRID) then
                  print *, 'BOUNDS card for a spherical structured grid must include ' // &
                           '4 lines.  I.e.'
                  print *, 'BOUNDS'
                  print *, '  r_min'
                  print *, '  r_max'
                  print *, 'END'
                endif
              endif
              stop
            endif            
            discretization%origin(X_DIRECTION) = grid%structured_grid%bounds(X_DIRECTION,LOWER)
            discretization%origin(Y_DIRECTION) = grid%structured_grid%bounds(Y_DIRECTION,LOWER)
            discretization%origin(Z_DIRECTION) = grid%structured_grid%bounds(Z_DIRECTION,LOWER)
        end select
      case ('GRAVITY')
        call InputReadDouble(input,option,option%gravity(X_DIRECTION))
        call InputErrorMsg(input,option,'x-direction','GRAVITY')
        call InputReadDouble(input,option,option%gravity(Y_DIRECTION))
        call InputErrorMsg(input,option,'y-direction','GRAVITY')
        call InputReadDouble(input,option,option%gravity(Z_DIRECTION))
        call InputErrorMsg(input,option,'z-direction','GRAVITY')
        if (option%myrank == option%io_rank .and. &
            option%print_to_screen) &
          write(option%fid_out,'(/," *GRAV",/, &
            & "  gravity    = "," [m/s^2]",3x,1p3e12.4 &
            & )') option%gravity(1:3)
      case ('MAX_CELLS_SHARING_A_VERTEX')
        if (associated(discretization%grid%unstructured_grid)) then
          call InputReadInt(input,option,discretization%grid% &
                            unstructured_grid%max_cells_sharing_a_vertex)
          call InputErrorMsg(input,option,'max_cells_sharing_a_vertex', &
                             'GRID')
        endif          
      case ('INVERT_Z')
      case ('STENCIL_WIDTH')
        call InputReadInt(input,option,discretization%stencil_width)
        call InputErrorMsg(input,option,'stencil_width', &
                           'GRID')
      case ('STENCIL_TYPE')
        call InputReadWord(input,option,word,PETSC_TRUE)
        call InputErrorMsg(input,option,'keyword','GRID')
        call StringToUpper(word)
        select case(trim(word))
          case ('BOX')
            discretization%stencil_type = DMDA_STENCIL_BOX
          case ('STAR')
            discretization%stencil_type = DMDA_STENCIL_STAR
          case default
            option%io_buffer = 'Keyword: ' // trim(word) // &
                 ' not recognized in DISCRETIZATION, second read.'
            call printErrMsg(option)
        end select
      case ('FLUX_METHOD')
        call InputReadWord(input,option,word,PETSC_TRUE)
        call InputErrorMsg(input,option,'keyword','GRID')
        call StringToUpper(word)
        select case(trim(word))
          case ('HYDRAULICS')
            call InputReadWord(input,option,word,PETSC_TRUE)
            call InputErrorMsg(input,option,'keyword','GRID')
            call StringToUpper(word)
            select case(trim(word))
              case ('TWO_POINT_FLUX')
                discretization%hydr_flux_method = TWO_POINT_FLUX
              case ('LSM_FLUX')
                discretization%hydr_flux_method = LSM_FLUX
                discretization%lsm_flux_method = PETSC_TRUE
              case default
                option%io_buffer = 'Keyword: ' // trim(word) // &
                  ' not recognized in DISCRETIZATION, second read.'
                call printErrMsg(option)
            end select
          case ('THERMAL')
            call InputReadWord(input,option,word,PETSC_TRUE)
            call InputErrorMsg(input,option,'keyword','GRID')
            call StringToUpper(word)
            select case(trim(word))
              case ('TWO_POINT_FLUX')
                discretization%therm_flux_method = TWO_POINT_FLUX
              case ('LSM_FLUX')
                discretization%therm_flux_method = LSM_FLUX
                discretization%lsm_flux_method = PETSC_TRUE
              case default
                option%io_buffer = 'Keyword: ' // trim(word) // &
                  ' not recognized in DISCRETIZATION, second read.'
                call printErrMsg(option)
            end select
          case ('MECHANICAL')
            call InputReadWord(input,option,word,PETSC_TRUE)
            call InputErrorMsg(input,option,'keyword','GRID')
            call StringToUpper(word)
            select case(trim(word))
              case ('TWO_POINT_FLUX')
                discretization%mech_flux_method = TWO_POINT_FLUX
              case ('LSM_FLUX')
                discretization%mech_flux_method = LSM_FLUX
                discretization%lsm_flux_method = PETSC_TRUE
              case default
                option%io_buffer = 'Keyword: ' // trim(word) // &
                  ' not recognized in DISCRETIZATION, second read.'
                call printErrMsg(option)
            end select
          case default
            option%io_buffer = 'Keyword: ' // trim(word) // &
            ' not recognized in DISCRETIZATION, second read.'
            call printErrMsg(option)
        end select
      case default
        option%io_buffer = 'Keyword: ' // trim(word) // &
                 ' not recognized in DISCRETIZATION, second read.'
        call printErrMsg(option)          
    end select 
  enddo  

  select case(discretization%itype)
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
      if (discretization%grid%structured_grid%invert_z_axis) then
        option%gravity(Z_DIRECTION) = -1.d0*option%gravity(Z_DIRECTION)
      endif
  end select
  
  if(discretization%lsm_flux_method) then
    if(discretization%stencil_type/=DMDA_STENCIL_BOX) then
      option%io_buffer='For LSM-Flux option, the stencil type needs to be DMDA_STENCIL_BOX'
      call printMsg(option)
      discretization%stencil_type = DMDA_STENCIL_BOX
      discretization%STENCIL_WIDTH = 2
    endif
  endif

end subroutine DiscretizationRead

! ************************************************************************** !
!
! DiscretizationCreateDMs: creates distributed, parallel meshes/grids
! If there are multiple degrees of freedom per grid cell, this will call 
! DiscretizationCreateDM() multiple times to create the DMs corresponding 
! to one degree of freedom grid cell and those corresponding to multiple 
! degrees of freedom per cell.
! author: Glenn Hammond
! date: 02/08/08
!
! ************************************************************************** !
subroutine DiscretizationCreateDMs(discretization,option)
      
  use Option_module    
      
  implicit none
  
  type(discretization_type) :: discretization
  type(option_type) :: option
      
  PetscInt :: ndof
  !PetscInt, parameter :: stencil_width = 1
  PetscErrorCode :: ierr
  PetscInt :: i
  type(unstructured_grid_type), pointer :: ugrid

  select case(discretization%itype)
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
      discretization%dm_index_to_ndof(ONEDOF) = 1
      discretization%dm_index_to_ndof(NPHASEDOF) = option%nphase
      discretization%dm_index_to_ndof(NFLOWDOF) = option%nflowdof
      discretization%dm_index_to_ndof(NTRANDOF) = option%ntrandof
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)

      ! petsc will call parmetis to calculate the graph/dual
#if !defined(PETSC_HAVE_PARMETIS)
      option%io_buffer = &
        'Must compile with Parmetis in order to use unstructured grids.'
      call printErrMsg(option)
#endif
    
      select case(discretization%grid%itype)
        case(IMPLICIT_UNSTRUCTURED_GRID)
          call UGridDecompose(discretization%grid%unstructured_grid, &
                              option)
          if(discretization%lsm_flux_method) then
            call UGridGrowStencilSupport(discretization%grid%unstructured_grid, &
                                         discretization%stencil_width, &
                                         discretization%grid%ghosted_level, &
                                         option)
          endif
        case(EXPLICIT_UNSTRUCTURED_GRID)
          ugrid => discretization%grid%unstructured_grid
#if 0          
          call ExplicitUGridDecompose(ugrid%explicit_grid, &
                                      ugrid%num_ghost_cells, &
                                      ugrid%global_offset, &
                                      ugrid%nmax, &
                                      ugrid%nlmax, &
                                      ugrid%ngmax, &
                                      ugrid%cell_ids_natural, &
                                      ugrid%cell_ids_petsc, &
                                      ugrid%ghost_cell_ids_petsc, &
                                      ugrid%ao_natural_to_petsc, &
                                      option)
#else
          call ExplicitUGridDecomposeNew(ugrid,option)
#endif
      end select
  end select


  !-----------------------------------------------------------------------
  ! Generate the DM objects that will manage communication.
  !-----------------------------------------------------------------------
  ndof = 1
  call DiscretizationCreateDM(discretization,discretization%dm_1dof, &
                              ndof,discretization%stencil_width, &
                              discretization%stencil_type,option)
  
  if (option%nflowdof > 0) then
    ndof = option%nflowdof
    call DiscretizationCreateDM(discretization,discretization%dm_nflowdof, &
                                ndof,discretization%stencil_width, &
                                discretization%stencil_type,option)
  endif
  
  if (option%ntrandof > 0) then
    ndof = option%ntrandof
    call DiscretizationCreateDM(discretization,discretization%dm_ntrandof, &
                                ndof,discretization%stencil_width, &
                                discretization%stencil_type,option)
  endif


  select case(discretization%itype)
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
      ! this function must be called to set up str_grid%lxs, etc.
      call StructGridComputeLocalBounds(discretization%grid%structured_grid, &
                                        discretization%dm_1dof%dm,option)    
      discretization%grid%nlmax = discretization%grid%structured_grid%nlmax
      discretization%grid%ngmax = discretization%grid%structured_grid%ngmax
      discretization%grid%global_offset = &
        discretization%grid%structured_grid%global_offset
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      discretization%grid%nmax = discretization%grid%unstructured_grid%nmax
      discretization%grid%nlmax = discretization%grid%unstructured_grid%nlmax
      discretization%grid%ngmax = discretization%grid%unstructured_grid%ngmax
      discretization%grid%global_offset = &
        discretization%grid%unstructured_grid%global_offset
  end select

end subroutine DiscretizationCreateDMs

! ************************************************************************** !
!
! DiscretizationCreateDM: creates a distributed, parallel mesh/grid
! author: Glenn Hammond
! date: 02/08/08
!
! ************************************************************************** !
subroutine DiscretizationCreateDM(discretization,dm_ptr,ndof,stencil_width, &
                                  stencil_type,option)

  use Option_module
  
  implicit none
  
  type(discretization_type) :: discretization
  type(dm_ptr_type), pointer :: dm_ptr
  PetscInt :: ndof
  PetscInt :: stencil_width,stencil_type
  type(option_type) :: option
  PetscErrorCode :: ierr

  select case(discretization%itype)
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
      call StructGridCreateDM(discretization%grid%structured_grid, &
                                  dm_ptr%dm,ndof,stencil_width,stencil_type,option)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call UGridCreateUGDM(discretization%grid%unstructured_grid, &
                           dm_ptr%ugdm,ndof,option)
      call DMShellCreate(option%mycomm,dm_ptr%dm,ierr)
      call DMShellSetGlobalToLocalVecScatter(dm_ptr%dm,dm_ptr%ugdm%scatter_gtol,ierr)
  end select

end subroutine DiscretizationCreateDM

! ************************************************************************** !
!
! DiscretizationCreateVector: Creates a global PETSc vector
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizationCreateVector(discretization,dm_index,vector, &
                                      vector_type,option)
  use Option_module                                      

  implicit none
  
  type(discretization_type) :: discretization
  PetscInt :: dm_index
  Vec :: vector
  PetscInt :: vector_type
  type(option_type) :: option
  PetscInt :: ndof
  PetscErrorCode :: ierr
  
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)

  select case(discretization%itype)
    case(STRUCTURED_GRID, STRUCTURED_GRID_MIMETIC)
      select case (vector_type)
        case(GLOBAL)
          call DMCreateGlobalVector(dm_ptr%dm,vector,ierr)
        case(LOCAL)
          call DMCreateLocalVector(dm_ptr%dm,vector,ierr)
        case(NATURAL)
          call DMDACreateNaturalVector(dm_ptr%dm,vector,ierr)
      end select
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call UGridDMCreateVector(discretization%grid%unstructured_grid, &
                               dm_ptr%ugdm,vector, &
                               vector_type,option)
  end select
  call VecSet(vector,0.d0,ierr)
  
end subroutine DiscretizationCreateVector

! ************************************************************************** !
!
! DiscretizationDuplicateVector: Creates a global PETSc vector
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizationDuplicateVector(discretization,vector1,vector2)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: vector1
  Vec :: vector2
  
  PetscErrorCode :: ierr
  call VecDuplicate(vector1,vector2,ierr)
  call VecCopy(vector1,vector2,ierr)
  
end subroutine DiscretizationDuplicateVector

! ************************************************************************** !
!
! DiscretizationGetDMPtrFromIndex: Returns the integer pointer for the DM referenced
! author: Glenn Hammond
! date: 02/08/08
!
! ************************************************************************** !
function DiscretizationGetDMPtrFromIndex(discretization,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  PetscInt :: dm_index
  
  type(dm_ptr_type), pointer :: DiscretizationGetDMPtrFromIndex
  
  select case (dm_index)
    case(ONEDOF)
      DiscretizationGetDMPtrFromIndex => discretization%dm_1dof
    case(NFLOWDOF)
      DiscretizationGetDMPtrFromIndex => discretization%dm_nflowdof
    case(NTRANDOF)
      DiscretizationGetDMPtrFromIndex => discretization%dm_ntrandof
  end select  
  
end function DiscretizationGetDMPtrFromIndex

! ************************************************************************** !
! ************************************************************************** !
function DiscretizationGetDMCPtrFromIndex(discretization,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  PetscInt :: dm_index
  
  type(dm_ptr_type), pointer :: DiscretizationGetDMCPtrFromIndex(:)
  
  select case (dm_index)
    case(NFLOWDOF)
      DiscretizationGetDMCPtrFromIndex => discretization%dmc_nflowdof
    case(NTRANDOF)
      DiscretizationGetDMCPtrFromIndex => discretization%dmc_ntrandof
  end select  
  
end function DiscretizationGetDMCPtrFromIndex

! ************************************************************************** !
!
! DiscretizationCreateJacobian: Creates Jacobian matrix associated with discretization
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizationCreateJacobian(discretization,dm_index,mat_type,Jacobian,option)

  use Option_module
  
  implicit none
  
#include "finclude/petscis.h"
#include "finclude/petscis.h90"

  type(discretization_type) :: discretization
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  MatType :: mat_type
  Mat :: Jacobian
  type(option_type) :: option
  PetscInt :: ndof, stencilsize
  PetscInt, pointer :: indices(:)
  PetscInt :: ngmax
  PetscInt :: imax, nlevels, ln, npatches, pn, i
  type(dm_ptr_type), pointer :: dm_ptr
  ISLocalToGlobalMapping :: ptmap
  PetscInt :: islocal

  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)


  select case(discretization%itype)
    case(STRUCTURED_GRID)
#ifndef DMGET
      call DMCreateMatrix(dm_ptr%dm,mat_type,Jacobian,ierr)
#else
      call DMGetMatrix(dm_ptr%dm,mat_type,Jacobian,ierr)
#endif
      call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE,ierr)
      call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr)
    case(UNSTRUCTURED_GRID)
      call UGridDMCreateJacobian(discretization%grid%unstructured_grid, &
                                 dm_ptr%ugdm,mat_type,Jacobian,option)
      call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE,ierr)
      call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr)
    case(UNSTRUCTURED_GRID_MIMETIC)
      select case(dm_index)
        case(NFLOWDOF)
          call MFDCreateJacobianLP(discretization%grid, discretization%MFD, mat_type, Jacobian, option)
          call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE,ierr)
          call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr)
        case(NTRANDOF)
          call UGridDMCreateJacobian(discretization%grid%unstructured_grid, &
                                     dm_ptr%ugdm,mat_type,Jacobian,option)
          call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE,ierr)
          call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr)
      end select
    case(STRUCTURED_GRID_MIMETIC)
#ifdef DASVYAT
      select case(dm_index)
        case(NFLOWDOF)
!          call MFDCreateJacobian(discretization%grid, discretization%MFD, mat_type, Jacobian, option)
          call MFDCreateJacobianLP(discretization%grid, discretization%MFD, mat_type, Jacobian, option)
          call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE,ierr)
          call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr)
        case(NTRANDOF)
#ifndef DMGET
          call DMCreateMatrix(dm_ptr%dm,mat_type,Jacobian,ierr)
#else
          call DMGetMatrix(dm_ptr%dm,mat_type,Jacobian,ierr)
#endif
          call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE,ierr)
          call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr)
      end select
#endif
  end select


end subroutine DiscretizationCreateJacobian

! ************************************************************************** !
!
! DiscretizationCreateInterpolation: Creates interpolation matrix associated 
! with discretization for geometric multigrid.
! author: Richard Mills
! date: 4/25/08.
!
! ************************************************************************** !
subroutine DiscretizationCreateInterpolation(discretization,dm_index, &
                                             interpolation,mg_levels_x, &
                                             mg_levels_y, mg_levels_z, &
                                             option)

  use Option_module
  
  implicit none
  
  type(discretization_type) :: discretization
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  Mat, pointer :: interpolation(:)
  PetscInt :: mg_levels_x, mg_levels_y, mg_levels_z
  type(option_type) :: option

  PetscInt :: mg_levels
  PetscInt :: refine_x, refine_y, refine_z
  type(dm_ptr_type), pointer :: dm_ptr
  type(dm_ptr_type), pointer :: dmc_ptr(:)
  PetscInt :: i
  type(dm_ptr_type), pointer :: dm_fine_ptr
    ! Used to point to finer-grid DM in the loop that constructst the 
    ! interpolation hierarchy.
  
  mg_levels = max(mg_levels_x, mg_levels_y, mg_levels_z)
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
!  dmc_ptr = DiscretizationGetDMCPtrFromIndex(discretization,dm_index)
  select case (dm_index)
    case(NFLOWDOF)
      allocate(discretization%dmc_nflowdof(mg_levels))
      do i=1, mg_levels
        discretization%dmc_nflowdof(i)%dm = 0
        nullify(discretization%dmc_nflowdof(i)%ugdm)
      enddo
      dmc_ptr => discretization%dmc_nflowdof
    case(NTRANDOF)
      allocate(discretization%dmc_ntrandof(mg_levels))
      do i=1, mg_levels
        discretization%dmc_ntrandof(i)%dm = 0
        nullify(discretization%dmc_ntrandof(i)%ugdm)
      enddo
      dmc_ptr => discretization%dmc_ntrandof
  end select  
   
  allocate(interpolation(mg_levels))

  select case(discretization%itype)
    case(STRUCTURED_GRID)
      dm_fine_ptr => dm_ptr
      refine_x = 2; refine_y = 2; refine_z = 2
      do i=mg_levels-1,1,-1
        ! If number of coarsenings performed so far exceeds mg_levels_x-1, 
        ! set refine_x = 1; likewise for y and z.
        if (i <= mg_levels - mg_levels_x ) refine_x = 1
        if (i <= mg_levels - mg_levels_y ) refine_y = 1
        if (i <= mg_levels - mg_levels_z ) refine_z = 1
        call DMDASetRefinementFactor(dm_fine_ptr%dm, refine_x, refine_y, refine_z, &
                                   ierr)
        call DMDASetInterpolationType(dm_fine_ptr%dm, DMDA_Q0, ierr)
        call DMCoarsen(dm_fine_ptr%dm, option%mycomm, dmc_ptr(i)%dm, ierr)
#ifndef DMGET
        call DMCreateInterpolation(dmc_ptr(i)%dm, dm_fine_ptr%dm, &
                                   interpolation(i), PETSC_NULL_OBJECT, ierr)
#else
        call DMGetInterpolation(dmc_ptr(i)%dm, dm_fine_ptr%dm, &
                                interpolation(i), PETSC_NULL_OBJECT, ierr)
#endif
        dm_fine_ptr => dmc_ptr(i)
      enddo
    case(UNSTRUCTURED_GRID)
  end select

end subroutine DiscretizationCreateInterpolation

! ************************************************************************** !
!
! DiscretizationCreateColoring: Creates ISColoring for discretization
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizationCreateColoring(discretization,dm_index,option,coloring)

  use Option_module
  
  implicit none

#include "finclude/petscis.h"
#include "finclude/petscis.h90"
  
  type(discretization_type) :: discretization
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(option_type) :: option
  ISColoring :: coloring

  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
    
  select case(discretization%itype)
    case(STRUCTURED_GRID)
#ifndef DMGET
      call DMCreateColoring(dm_ptr%dm,IS_COLORING_GLOBAL,MATBAIJ,coloring,&
                            ierr)
#else
      call DMGetColoring(dm_ptr%dm,IS_COLORING_GLOBAL,MATBAIJ,coloring,ierr)
#endif
      ! I have set the above to use matrix type MATBAIJ, as that is what we 
      ! usually want (note: for DAs with 1 degree of freedom per grid cell, 
      ! the MATAIJ and MATBAIJ colorings should be equivalent).  What we should 
      ! eventually do here is query the type of the Jacobian matrix, but I'm 
      ! not sure of the best way to do this, as this is currently stashed in 
      ! the 'solver' object. --RTM
    case(UNSTRUCTURED_GRID)
  end select
  
end subroutine DiscretizationCreateColoring

! ************************************************************************** !
!
! DiscretizationGlobalToLocal: Performs global to local communication with DM
! author: Glenn Hammond
! date: 10/24/07
!
! Note that 'dm_index' should correspond to one of the macros defined 
! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
!
! ************************************************************************** !
subroutine DiscretizationGlobalToLocal(discretization,global_vec,local_vec,dm_index)

  implicit none

  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
    
  call DMGlobalToLocalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec,ierr)
  call DMGlobalToLocalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec,ierr)
  
end subroutine DiscretizationGlobalToLocal

! ************************************************************************** !
!
! DiscretizationGlobalToLocalFaces: Performs global to local communication for MFD
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizationGlobalToLocalFaces(discretization,global_vec,local_vec,dm_index)

  implicit none
  

  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
    
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC,UNSTRUCTURED_GRID_MIMETIC)
      call DiscretizGlobalToLocalFacesBegin(discretization,global_vec,local_vec,dm_index)
      call DiscretizGlobalToLocalFacesEnd(discretization,global_vec,local_vec,dm_index)
  end select
  
end subroutine DiscretizationGlobalToLocalFaces


! ************************************************************************** !
!
! DiscretizationGlobalToLocalLP: Performs global to local communication for MFD
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizationGlobalToLocalLP(discretization,global_vec,local_vec,dm_index)

  implicit none
  

  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
    
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC,UNSTRUCTURED_GRID_MIMETIC)
      call DiscretizGlobalToLocalLPBegin(discretization,global_vec,local_vec,dm_index)
      call DiscretizGlobalToLocalLPEnd(discretization,global_vec,local_vec,dm_index)
  end select
  
end subroutine DiscretizationGlobalToLocalLP
  
! ************************************************************************** !
!
! DiscretizationLocalToGlobal: Performs local to global communication with DM
! author: Glenn Hammond
! date: 1/02/08
!
! Note that 'dm_index' should correspond to one of the macros defined 
! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
!
! ************************************************************************** !
subroutine DiscretizationLocalToGlobal(discretization,local_vec,global_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec
  Vec :: global_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMLocalToGlobalBegin(dm_ptr%dm,local_vec,INSERT_VALUES,global_vec,ierr)
      call DMLocalToGlobalEnd(dm_ptr%dm,local_vec,INSERT_VALUES,global_vec,ierr)
   case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_ltog,local_vec,global_vec, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
      call VecScatterEnd(dm_ptr%ugdm%scatter_ltog,local_vec,global_vec, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)
  end select
  
end subroutine DiscretizationLocalToGlobal
  
! ************************************************************************** !
!
! DiscretizationLocalToLocal: Performs local to local communication with DM
! author: Glenn Hammond
! date: 11/14/07
!
! ************************************************************************** !
subroutine DiscretizationLocalToLocal(discretization,local_vec1,local_vec2,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDALocalToLocalBegin(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2,ierr)
      call DMDALocalToLocalEnd(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_ltol,local_vec1,local_vec2, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
      call VecScatterEnd(dm_ptr%ugdm%scatter_ltol,local_vec1,local_vec2, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)    
  end select
  
end subroutine DiscretizationLocalToLocal
  
! ************************************************************************** !
!
! DiscretizationLocalToLocalFaces: Performs local to local communication for face unknowns
! author: Daniil Svyatskiy
! date: 11/14/07
!
! ************************************************************************** !
subroutine DiscretizationLocalToLocalFaces(discretization,local_vec1,local_vec2,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC,UNSTRUCTURED_GRID_MIMETIC)
      call DiscretizLocalToLocalFacesBegin(discretization,local_vec1,local_vec2,dm_index)
      call DiscretizLocalToLocalFacesEnd(discretization,local_vec1,local_vec2,dm_index)
    case(UNSTRUCTURED_GRID)
  end select
  
end subroutine DiscretizationLocalToLocalFaces

! ************************************************************************** !
!
! DiscretizationLocalToLocalLP: Performs local to local communication for face unknowns
! author: Daniil Svyatskiy
! date: 11/14/07
!
! ************************************************************************** !
subroutine DiscretizationLocalToLocalLP(discretization,local_vec1,local_vec2,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC,UNSTRUCTURED_GRID_MIMETIC)
      call DiscretizLocalToLocalLPBegin(discretization,local_vec1,local_vec2,dm_index)
      call DiscretizLocalToLocalLPEnd(discretization,local_vec1,local_vec2,dm_index)
    case(UNSTRUCTURED_GRID)
  end select
  
end subroutine DiscretizationLocalToLocalLP
  
 
! ************************************************************************** !
!
! DiscretizationGlobalToNatural: Performs global to natural communication with DM
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizationGlobalToNatural(discretization,global_vec,natural_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: natural_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDAGlobalToNaturalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec,ierr)
      call DMDAGlobalToNaturalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
      call VecScatterEnd(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)       
  end select
  
end subroutine DiscretizationGlobalToNatural

! ************************************************************************** !
!
! DiscretizationNaturalToGlobal: Performs natural to global communication with DM
! author: Glenn Hammond
! date: 01/12/08
!
! ************************************************************************** !
subroutine DiscretizationNaturalToGlobal(discretization,natural_vec,global_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: natural_vec
  Vec :: global_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDANaturalToGlobalBegin(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec,ierr)
      call DMDANaturalToGlobalEnd(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_ntog,natural_vec,global_vec, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
      call VecScatterEnd(dm_ptr%ugdm%scatter_ntog,natural_vec,global_vec, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)
  end select
  
end subroutine DiscretizationNaturalToGlobal

! ************************************************************************** !
!
! DiscretizationGlobalToLocalBegin: Begins global to local communication with DM
! author: Glenn Hammond
! date: 10/24/07
!
! Note that 'dm_index' should correspond to one of the macros defined 
! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
!
! ************************************************************************** !
subroutine DiscretizationGlobalToLocalBegin(discretization,global_vec,local_vec,dm_index)



  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMGlobalToLocalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_gtol,global_vec,local_vec, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
  end select
  
end subroutine DiscretizationGlobalToLocalBegin
  
! ************************************************************************** !
!
! DiscretizationGlobalToLocalEnd: Ends global to local communication with DM
! author: Glenn Hammond
! date: 10/24/07
!
! Note that 'dm_index' should correspond to one of the macros defined 
! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
!
! ************************************************************************** !
subroutine DiscretizationGlobalToLocalEnd(discretization,global_vec,local_vec,dm_index)

 

 implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMGlobalToLocalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterEnd(dm_ptr%ugdm%scatter_gtol,global_vec,local_vec, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)
  end select
  
end subroutine DiscretizationGlobalToLocalEnd
  
! ************************************************************************** !
!
! DiscretizLocalToLocalFacesBegin: Begins Local to local communication with MFD
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizLocalToLocalFacesBegin(discretization,local_vec1,local_vec2,dm_index)


  use MFD_Aux_module

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC)
      call VecScatterBegin( discretization%MFD%scatter_ltol_faces,local_vec1,local_vec2 , &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizLocalToLocalFacesBegin
  
! ************************************************************************** !
!
! DiscretizLocalToLocalFacesEnd: Ends local  to local communication with DM
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizLocalToLocalFacesEnd(discretization,local_vec1,local_vec2,dm_index)

 
 use MFD_Aux_module

 implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC)
      call VecScatterEnd( discretization%MFD%scatter_ltol_faces,  local_vec1, local_vec2, &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizLocalToLocalFacesEnd
  
! DiscretizGlobalToLocalFacesBegin: Begins global to local communication with MFD
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizGlobalToLocalFacesBegin(discretization,global_vec,local_vec,dm_index)


  use MFD_Aux_module

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC)
      call VecScatterBegin( discretization%MFD%scatter_gtol_faces, global_vec, local_vec, &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizGlobalToLocalFacesBegin
  
! ************************************************************************** !
!
! DiscretizGlobalToLocalFacesEnd: Ends global to local communication with DM
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizGlobalToLocalFacesEnd(discretization,global_vec,local_vec,dm_index)

 
 use MFD_Aux_module

 implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC)
      call VecScatterEnd( discretization%MFD%scatter_gtol_faces, global_vec, local_vec, &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizGlobalToLocalFacesEnd

! ************************************************************************** !
!
! DiscretizLocalToLocalLPBegin: Begins Local to local communication with MFD
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizLocalToLocalLPBegin(discretization,local_vec1,local_vec2,dm_index)


  use MFD_Aux_module

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC)
      call VecScatterBegin( discretization%MFD%scatter_ltol_LP,local_vec1,local_vec2 , &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizLocalToLocalLPBegin
  
! ************************************************************************** !
!
! DiscretizLocalToLocalLPEnd: Ends local  to local communication with DM
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizLocalToLocalLPEnd(discretization,local_vec1,local_vec2,dm_index)

 
 use MFD_Aux_module

 implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC)
      call VecScatterEnd( discretization%MFD%scatter_ltol_LP,  local_vec1, local_vec2, &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizLocalToLocalLPEnd
  
! DiscretizGlobalToLocalLPBegin: Begins global to local communication with MFD
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizGlobalToLocalLPBegin(discretization,global_vec,local_vec,dm_index)


  use MFD_Aux_module

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin( discretization%MFD%scatter_gtol_LP, global_vec, local_vec, &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizGlobalToLocalLPBegin
  
! ************************************************************************** !
!
! DiscretizGlobalToLocalLPEnd: Ends global to local communication with DM
! author: Daniil Svyatskiy
! date: 07/20/10
!
!
! ************************************************************************** !
subroutine DiscretizGlobalToLocalLPEnd(discretization,global_vec,local_vec,dm_index)

 
 use MFD_Aux_module

 implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: local_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  
  select case(discretization%itype)
    case(STRUCTURED_GRID_MIMETIC,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterEnd( discretization%MFD%scatter_gtol_LP, global_vec, local_vec, &
                                INSERT_VALUES,SCATTER_FORWARD, ierr)
  end select
  
end subroutine DiscretizGlobalToLocalLPEnd


  
! ************************************************************************** !
!
! DiscretizationLocalToLocalBegin: Begins local to local communication with DM
! author: Glenn Hammond
! date: 11/14/07
!
! ************************************************************************** !
subroutine DiscretizationLocalToLocalBegin(discretization,local_vec1,local_vec2,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDALocalToLocalBegin(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_ltol,local_vec1,local_vec2, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
  end select

end subroutine DiscretizationLocalToLocalBegin
  
! ************************************************************************** !
!
! DiscretizationLocalToLocalEnd: Ends local to local communication with DM
! author: Glenn Hammond
! date: 11/14/07
!
! ************************************************************************** !
subroutine DiscretizationLocalToLocalEnd(discretization,local_vec1,local_vec2,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: local_vec1
  Vec :: local_vec2
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDALocalToLocalEnd(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterEnd(dm_ptr%ugdm%scatter_ltol,local_vec1,local_vec2, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)    
  end select

end subroutine DiscretizationLocalToLocalEnd
  
! ************************************************************************** !
!
! DiscretizGlobalToNaturalBegin: Begins global to natural communication with DM
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizGlobalToNaturalBegin(discretization,global_vec,natural_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: natural_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDAGlobalToNaturalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterBegin(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
                           INSERT_VALUES,SCATTER_FORWARD,ierr)
  end select
  
end subroutine DiscretizGlobalToNaturalBegin

! ************************************************************************** !
!
! DiscretizGlobalToNaturalEnd: Ends global to natural communication with DM
! author: Glenn Hammond
! date: 10/24/07
!
! ************************************************************************** !
subroutine DiscretizGlobalToNaturalEnd(discretization,global_vec,natural_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: natural_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDAGlobalToNaturalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      call VecScatterEnd(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
                         INSERT_VALUES,SCATTER_FORWARD,ierr)       
  end select
  
end subroutine DiscretizGlobalToNaturalEnd

! ************************************************************************** !
!
! DiscretizNaturalToGlobalBegin: Begins natural to global communication with DM
! author: Glenn Hammond
! date: 01/12/08
!
! ************************************************************************** !
subroutine DiscretizNaturalToGlobalBegin(discretization,natural_vec,global_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: natural_vec
  Vec :: global_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDANaturalToGlobalBegin(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec,ierr)
    case(UNSTRUCTURED_GRID)
  end select
  
end subroutine DiscretizNaturalToGlobalBegin

! ************************************************************************** !
!
! DiscretizNaturalToGlobalEnd: Ends natural to global communication with DM
! author: Glenn Hammond
! date: 01/12/08
!
! ************************************************************************** !
subroutine DiscretizNaturalToGlobalEnd(discretization,natural_vec,global_vec,dm_index)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: natural_vec
  Vec :: global_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  type(dm_ptr_type), pointer :: dm_ptr
  
  dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDANaturalToGlobalEnd(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec,ierr)
    case(UNSTRUCTURED_GRID)
  end select
  
end subroutine DiscretizNaturalToGlobalEnd

! ************************************************************************** !
!
! DiscretizationUpdateTVDGhosts: Updates tvd extended ghost cell values
! author: Glenn Hammond
! date: 02/04/12
!
! ************************************************************************** !
subroutine DiscretizationUpdateTVDGhosts(discretization,global_vec, &
                                         tvd_ghost_vec)

  implicit none
  
  type(discretization_type) :: discretization
  Vec :: global_vec
  Vec :: tvd_ghost_vec
  PetscInt :: dm_index
  PetscErrorCode :: ierr
  
  call VecScatterBegin(discretization%tvd_ghost_scatter,global_vec, &
                       tvd_ghost_vec,INSERT_VALUES,SCATTER_FORWARD,ierr)
  call VecScatterEnd(discretization%tvd_ghost_scatter,global_vec, &
                       tvd_ghost_vec,INSERT_VALUES,SCATTER_FORWARD,ierr)
  
end subroutine DiscretizationUpdateTVDGhosts

! ************************************************************************** !
!
! DiscretAOApplicationToPetsc: Maps application ordering to petsc
! author: Glenn Hammond
! date: 10/12/12
!
! ************************************************************************** !
subroutine DiscretAOApplicationToPetsc(discretization,int_array)

  implicit none
  
#include "finclude/petscao.h"  
  
  type(discretization_type) :: discretization
  PetscInt :: int_array(:)
  PetscErrorCode :: ierr
  
  AO :: ao
  
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      call DMDAGetAO(discretization%dm_1dof,ao,ierr)
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      ao = discretization%grid%unstructured_grid%ao_natural_to_petsc
  end select
  call AOApplicationToPetsc(ao,size(int_array),int_array,ierr)
  
end subroutine DiscretAOApplicationToPetsc

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

  implicit none
  
  type(discretization_type), pointer :: discretization
  
  PetscErrorCode :: ierr
  PetscInt :: i
    
  if (.not.associated(discretization)) return
      
  select case(discretization%itype)
    case(STRUCTURED_GRID,STRUCTURED_GRID_MIMETIC)
      if (discretization%dm_1dof%dm /= 0) &
        call DMDestroy(discretization%dm_1dof%dm,ierr)
      discretization%dm_1dof%dm = 0
      if (discretization%dm_nflowdof%dm /= 0) &
        call DMDestroy(discretization%dm_nflowdof%dm,ierr)
      discretization%dm_nflowdof%dm = 0
      if (discretization%dm_ntrandof%dm /= 0) &
        call DMDestroy(discretization%dm_ntrandof%dm,ierr)
      discretization%dm_ntrandof%dm = 0
      if (associated(discretization%dmc_nflowdof)) then
        do i=1,size(discretization%dmc_nflowdof)
          call DMDestroy(discretization%dmc_nflowdof(i)%dm,ierr)
        enddo
        deallocate(discretization%dmc_nflowdof)
        nullify(discretization%dmc_nflowdof)
      endif
      if (associated(discretization%dmc_ntrandof)) then
        do i=1,size(discretization%dmc_ntrandof)
          call DMDestroy(discretization%dmc_ntrandof(i)%dm,ierr)
        enddo
        deallocate(discretization%dmc_ntrandof)
        nullify(discretization%dmc_ntrandof)
      endif
    case(UNSTRUCTURED_GRID,UNSTRUCTURED_GRID_MIMETIC)
      if (associated(discretization%dm_1dof%ugdm)) &
        call UGridDMDestroy(discretization%dm_1dof%ugdm)
      if (associated(discretization%dm_nflowdof%ugdm)) &
        call UGridDMDestroy(discretization%dm_nflowdof%ugdm)
      if (associated(discretization%dm_ntrandof%ugdm)) &
        call UGridDMDestroy(discretization%dm_ntrandof%ugdm)
  end select
  if (associated(discretization%dm_1dof)) &
    deallocate(discretization%dm_1dof)
  nullify(discretization%dm_1dof)
  if (associated(discretization%dm_nflowdof)) &
    deallocate(discretization%dm_nflowdof)
  nullify(discretization%dm_nflowdof)
  if (associated(discretization%dm_ntrandof)) &
    deallocate(discretization%dm_ntrandof)
  nullify(discretization%dm_ntrandof)


  if (associated(discretization%MFD)) &
        call MFDAuxDestroy(discretization%MFD) 
  nullify(discretization%MFD)

  if (discretization%tvd_ghost_scatter /= 0) &
    call VecScatterDestroy(discretization%tvd_ghost_scatter)
  
end subroutine DiscretizationDestroy
 
end module Discretization_module
Tip: Filter by directory path e.g. /media app.js to search for public/media/app.js.
Tip: Use camelCasing e.g. ProjME to search for ProjectModifiedEvent.java.
Tip: Filter by extension type e.g. /repo .js to search for all .js files in the /repo directory.
Tip: Separate your search with spaces e.g. /ssh pom.xml to search for src/ssh/pom.xml.
Tip: Use ↑ and ↓ arrow keys to navigate and return to view the file.
Tip: You can also navigate files with Ctrl+j (next) and Ctrl+k (previous) and view the file with Ctrl+o.
Tip: You can also navigate files with Alt+j (next) and Alt+k (previous) and view the file with Alt+o.