Source

pflotran-dev / src / pflotran / unstructured_explicit.F90

Full commit
   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
module Unstructured_Explicit_module
  
  use Geometry_module
  use Unstructured_Grid_Aux_module
  
  implicit none

  private 
  
#include "definitions.h"
#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"
#include "finclude/petscis.h"
#include "finclude/petscis.h90"
#if defined(SCORPIO)
  include "scorpiof.h"
#endif

  public :: ExplicitUGridRead, &
            ExplicitUGridDecompose, &
            ExplicitUGridReadInParallel, &
            ExplicitUGridDecomposeNew, &
            ExplicitUGridSetInternConnect, &
            ExplicitUGridSetCellCentroids, &
            ExplicitUGridComputeVolumes, &
            ExplicitUGridSetBoundaryConnect, &
            ExplicitUGridSetConnections

contains

! ************************************************************************** !
!
! ExplicitUGridRead: Reads an explicit unstructured grid
! author: Glenn Hammond
! date: 05/14/12
!
! ************************************************************************** !
subroutine ExplicitUGridRead(explicit_grid,filename,option)

  use Input_module
  use Option_module
  use String_module
  
  implicit none
  
  type(unstructured_explicit_type) :: explicit_grid
  character(len=MAXSTRINGLENGTH) :: filename
  type(option_type) :: option
  
  type(input_type), pointer :: input
  character(len=MAXSTRINGLENGTH) :: string, hint
  character(len=MAXWORDLENGTH) :: word, card
  PetscInt :: fileid, icell, iconn, id_up, id_dn
  
  PetscInt :: num_cells
  PetscInt :: num_connections
  PetscInt :: num_elems
  
  fileid = 86
  input => InputCreate(fileid,filename,option)

! Format of explicit unstructured grid file
! id_, id_up_, id_dn_ = integer
! x_, y_, z_, area_, volume_ = real
! -----------------------------------------------------------------
! CELLS <integer>    integer = # cells (N)
! id_1 x_1 y_1 z_1 volume_1
! id_2 x_2 y_2 z_2 volume_2
! ...
! ...
! id_N x_N y_N z_N volume_N
! CONNECTIONS <integer>   integer = # connections (M)
! id_up_1 id_dn_1 x_1 y_1 z_1 area_1
! id_up_2 id_dn_2 x_2 y_2 z_2 area_2
! ...
! ...
! id_up_M id_dn_M x_M y_M z_M area_M
! -----------------------------------------------------------------


  do
    call InputReadFlotranString(input,option)
    if (InputError(input)) exit

    call InputReadWord(input,option,word,PETSC_FALSE)
    call StringToUpper(word)
    hint = trim(word)
  
    select case(word)
      case('CELLS')
        hint = 'Explicit Unstructured Grid CELLS'
        call InputReadInt(input,option,num_cells)
        call InputErrorMsg(input,option,'number of cells',hint)
        explicit_grid%num_cells_global = num_cells
        allocate(explicit_grid%cell_ids(num_cells))
        explicit_grid%cell_ids = 0
        allocate(explicit_grid%cell_volumes(num_cells))
        explicit_grid%cell_volumes = 0
        allocate(explicit_grid%cell_centroids(num_cells))
        do icell = 1, num_cells
          explicit_grid%cell_centroids(icell)%x = 0.d0
          explicit_grid%cell_centroids(icell)%y = 0.d0
          explicit_grid%cell_centroids(icell)%z = 0.d0
        enddo
        do icell = 1, num_cells
          call InputReadFlotranString(input,option)
          call InputReadStringErrorMsg(input,option,hint)  
          call InputReadInt(input,option,explicit_grid%cell_ids(icell))
          call InputErrorMsg(input,option,'cell id',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%cell_centroids(icell)%x)
          call InputErrorMsg(input,option,'cell x coordinate',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%cell_centroids(icell)%y)
          call InputErrorMsg(input,option,'cell y coordinate',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%cell_centroids(icell)%z)
          call InputErrorMsg(input,option,'cell z coordinate',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%cell_volumes(icell))
          call InputErrorMsg(input,option,'cell volume',hint)
        enddo
      case('CONNECTIONS')
        hint = 'Explicit Unstructured Grid CONNECTIONS'
        call InputReadInt(input,option,num_connections)
        call InputErrorMsg(input,option,'number of connections',hint)
        allocate(explicit_grid%connections(2,num_connections))
        explicit_grid%connections = 0
        allocate(explicit_grid%face_areas(num_connections))
        explicit_grid%face_areas = 0    
        allocate(explicit_grid%face_centroids(num_connections))
        do iconn = 1, num_connections
          explicit_grid%face_centroids(iconn)%x = 0.d0
          explicit_grid%face_centroids(iconn)%y = 0.d0
          explicit_grid%face_centroids(iconn)%z = 0.d0
        enddo
        do iconn = 1, num_connections
          call InputReadFlotranString(input,option)
          call InputReadStringErrorMsg(input,option,hint)  
          call InputReadInt(input,option, &
                            explicit_grid%connections(1,iconn))
          call InputErrorMsg(input,option,'cell id upwind',hint)
          call InputReadInt(input,option, &
                            explicit_grid%connections(2,iconn))
          call InputErrorMsg(input,option,'cell id downwind',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%face_centroids(iconn)%x)
          call InputErrorMsg(input,option,'face x coordinate',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%face_centroids(iconn)%y)
          call InputErrorMsg(input,option,'face y coordinate',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%face_centroids(iconn)%z)
          call InputErrorMsg(input,option,'face z coordinate',hint)
          call InputReadDouble(input,option, &
                               explicit_grid%face_areas(iconn))
          call InputErrorMsg(input,option,'face area',hint)
        enddo
      case('ELEMENTS')
        hint = 'Explicit Unstructured Grid ELEMENTS'
        call InputReadInt(input,option,num_elems)
        call InputErrorMsg(input,option,'number of elements',hint)
        explicit_grid%num_elems = num_elems
        ! Assuming only triangular elements (which is the dual mesh of the voronoi read) in DFN -- Karra 12/17/12
        allocate(explicit_grid%cell_connectivity(3,num_elems)) 
        explicit_grid%cell_connectivity = 0
        do iconn = 1, num_elems
          call InputReadFlotranString(input,option)
          call InputReadStringErrorMsg(input,option,hint)  
          call InputReadInt(input,option, &
                            explicit_grid%cell_connectivity(1,iconn))
          call InputErrorMsg(input,option,'cell vertex 1',hint)
          call InputReadInt(input,option, &
                            explicit_grid%cell_connectivity(2,iconn))
          call InputErrorMsg(input,option,'cell vertex 2',hint)
          call InputReadInt(input,option, &
                            explicit_grid%cell_connectivity(3,iconn))
          call InputErrorMsg(input,option,'cell vertex 3',hint)
        enddo
      case('VERTICES')
        hint = 'Explicit Unstructured Grid VERTICES'     
        allocate(explicit_grid%vertex_coordinates(explicit_grid%num_cells_global))
        do icell = 1, explicit_grid%num_cells_global
          explicit_grid%vertex_coordinates(icell)%x = 0.d0
          explicit_grid%vertex_coordinates(icell)%y = 0.d0
          explicit_grid%vertex_coordinates(icell)%z = 0.d0
        enddo
        do icell = 1, explicit_grid%num_cells_global
          call InputReadFlotranString(input,option)
          call InputReadStringErrorMsg(input,option,hint)  
          call InputReadDouble(input,option, &
                            explicit_grid%vertex_coordinates(icell)%x)
          call InputErrorMsg(input,option,'vertex 1',hint)
          call InputReadDouble(input,option, &
                            explicit_grid%vertex_coordinates(icell)%y)
          call InputErrorMsg(input,option,'vertex 2',hint)
          call InputReadDouble(input,option, &
                            explicit_grid%vertex_coordinates(icell)%z)
          call InputErrorMsg(input,option,'vertex 3',hint)
        enddo
      case default
        option%io_buffer = 'Keyword: ' // trim(word) // &
                           ' not recognized while reading explicit ' // &
                           'unstructured grid.'
        call printErrMsg(option)        
    end select
  enddo

  call InputDestroy(input)

end subroutine ExplicitUGridRead

! ************************************************************************** !
!
! ExplicitUGridReadInParallel: Reads an explicit unstructured grid in parallel
! author: Glenn Hammond
! date: 10/03/12
!
! ************************************************************************** !
subroutine ExplicitUGridReadInParallel(explicit_grid,filename,option)

  use Input_module
  use Option_module
  use String_module
  
  implicit none
  
  type(unstructured_explicit_type) :: explicit_grid
  character(len=MAXSTRINGLENGTH) :: filename
  type(option_type) :: option
  
  type(input_type), pointer :: input
  character(len=MAXSTRINGLENGTH) :: string, hint
  character(len=MAXWORDLENGTH) :: word, card
  PetscInt :: fileid, icell, iconn, irank, remainder, temp_int, num_to_read
  
  PetscInt :: num_cells, num_connections, num_elems
  PetscInt :: num_cells_local, num_cells_local_save
  PetscInt :: num_connections_local, num_connections_local_save
  PetscInt :: num_elems_local, num_elems_local_save
  PetscMPIInt :: status_mpi(MPI_STATUS_SIZE)
  PetscMPIInt :: int_mpi
  PetscErrorCode :: ierr
  PetscReal, allocatable :: temp_real_array(:,:)
  PetscInt, allocatable :: temp_int_array(:,:)
  
! Format of explicit unstructured grid file
! id_, id_up_, id_dn_ = integer
! x_, y_, z_, area_, volume_ = real
! definitions
! id_ = id of grid cell
! id_up_ = id of upwind grid cell in connection
! id_dn_ = id of downwind grid cell in connection
! x_ = x coordinate of cell center
! y_ = y coordinate of cell center
! z_ = z coordinate of cell center
! volume_ = volume of grid cell
! -----------------------------------------------------------------
! CELLS <integer>    integer = # cells (N)
! id_1 x_1 y_1 z_1 volume_1
! id_2 x_2 y_2 z_2 volume_2
! ...
! ...
! id_N x_N y_N z_N volume_N
! CONNECTIONS <integer>   integer = # connections (M)
! id_up_1 id_dn_1 x_1 y_1 z_1 area_1
! id_up_2 id_dn_2 x_2 y_2 z_2 area_2
! ...
! ...
! id_up_M id_dn_M x_M y_M z_M area_M
! -----------------------------------------------------------------

  if (option%myrank == option%io_rank) then
  
    fileid = 86
    input => InputCreate(fileid,filename,option)

    call InputReadFlotranString(input,option)
    ! read CELL card, though we already know the
    call InputReadWord(input,option,card,PETSC_TRUE)
    word = 'CELLS'
    call InputErrorMsg(input,option,word,card)
    if (.not.StringCompare(word,card)) then
      option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
        '" in explicit grid file.'
      call printErrMsgByRank(option)
    endif
  
    hint = 'Explicit Unstructured Grid CELLS'
    call InputReadInt(input,option,temp_int)
    call InputErrorMsg(input,option,'number of cells',hint)
  endif
  
  call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
                 option%mycomm,ierr)
  num_cells = temp_int
  explicit_grid%num_cells_global = num_cells

   ! divide cells across ranks
  num_cells_local = num_cells/option%mycommsize 
  num_cells_local_save = num_cells_local
  remainder = num_cells - &
              num_cells_local*option%mycommsize
  if (option%myrank < remainder) num_cells_local = &
                                 num_cells_local + 1

  allocate(explicit_grid%cell_ids(num_cells_local))
  explicit_grid%cell_ids = 0
  allocate(explicit_grid%cell_volumes(num_cells_local))
  explicit_grid%cell_volumes = 0
  allocate(explicit_grid%cell_centroids(num_cells_local))
  do icell = 1, num_cells_local
    explicit_grid%cell_centroids(icell)%x = 0.d0
    explicit_grid%cell_centroids(icell)%y = 0.d0
    explicit_grid%cell_centroids(icell)%z = 0.d0
  enddo
  
  ! for now, read all cells from ASCII file through io_rank and communicate
  ! to other ranks
  if (option%myrank == option%io_rank) then
    allocate(temp_real_array(5,num_cells_local_save+1))
    ! read for other processors
    do irank = 0, option%mycommsize-1
      temp_real_array = -999.d0
      num_to_read = num_cells_local_save
      if (irank < remainder) num_to_read = num_to_read + 1
      do icell = 1, num_to_read
        call InputReadFlotranString(input,option)
        call InputReadStringErrorMsg(input,option,hint)  
        call InputReadInt(input,option,temp_int)
        call InputErrorMsg(input,option,'cell id',hint)
        temp_real_array(1,icell) = dble(temp_int)
        call InputReadDouble(input,option,temp_real_array(2,icell))
        call InputErrorMsg(input,option,'cell x coordinate',hint)
        call InputReadDouble(input,option,temp_real_array(3,icell))
        call InputErrorMsg(input,option,'cell y coordinate',hint)
        call InputReadDouble(input,option,temp_real_array(4,icell))
        call InputErrorMsg(input,option,'cell z coordinate',hint)
        call InputReadDouble(input,option,temp_real_array(5,icell))
        call InputErrorMsg(input,option,'cell volume',hint)
      enddo

      ! if the cells reside on io_rank
      if (irank == option%io_rank) then
#if UGRID_DEBUG
        write(string,*) num_cells_local
        string = trim(adjustl(string)) // ' cells stored on p0'
        print *, trim(string)
#endif
        do icell = 1, num_cells_local
          explicit_grid%cell_ids(icell) = int(temp_real_array(1,icell))
          explicit_grid%cell_centroids(icell)%x = temp_real_array(2,icell)
          explicit_grid%cell_centroids(icell)%y = temp_real_array(3,icell)
          explicit_grid%cell_centroids(icell)%z = temp_real_array(4,icell)
          explicit_grid%cell_volumes(icell) = temp_real_array(5,icell)
        enddo
      else
        ! otherwise communicate to other ranks
#if UGRID_DEBUG
        write(string,*) num_to_read
        write(word,*) irank
        string = trim(adjustl(string)) // ' cells sent from p0 to p' // &
                 trim(adjustl(word))
        print *, trim(string)
#endif
        int_mpi = num_to_read*5
        call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
                      num_to_read,option%mycomm,ierr)
      endif
    enddo
  else
    ! other ranks post the recv
#if UGRID_DEBUG
    write(string,*) num_cells_local
    write(word,*) option%myrank
    string = trim(adjustl(string)) // ' cells received from p0 at p' // &
              trim(adjustl(word))
    print *, trim(string)
#endif
    allocate(temp_real_array(5,num_cells_local))
    int_mpi = num_cells_local*5
    call MPI_Recv(temp_real_array,int_mpi, &
                  MPI_DOUBLE_PRECISION,option%io_rank, &
                  MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
    do icell = 1, num_cells_local
      explicit_grid%cell_ids(icell) = int(temp_real_array(1,icell))
      explicit_grid%cell_centroids(icell)%x = temp_real_array(2,icell)
      explicit_grid%cell_centroids(icell)%y = temp_real_array(3,icell)
      explicit_grid%cell_centroids(icell)%z = temp_real_array(4,icell)
      explicit_grid%cell_volumes(icell) = temp_real_array(5,icell)
    enddo
    
  endif
  deallocate(temp_real_array)
  
  if (option%myrank == option%io_rank) then
  
 
    call InputReadFlotranString(input,option)
    ! read CONNECTIONS card, though we already know the
    call InputReadWord(input,option,card,PETSC_TRUE)
    word = 'CONNECTIONS'
    call InputErrorMsg(input,option,word,card)
    if (.not.StringCompare(word,card)) then
      option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
        '" in explicit grid file.'
      call printErrMsgByRank(option)
    endif
  
    hint = 'Explicit Unstructured Grid CONNECTIONS'
    call InputReadInt(input,option,temp_int)
    call InputErrorMsg(input,option,'number of connections',hint)
  endif
  
  int_mpi = 1
  call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER,option%io_rank, &
                  option%mycomm,ierr)
  num_connections = temp_int
        
   ! divide cells across ranks
  num_connections_local = num_connections/option%mycommsize 
  num_connections_local_save = num_connections_local
  remainder = num_connections - &
              num_connections_local*option%mycommsize
  if (option%myrank < remainder) num_connections_local = &
                                 num_connections_local + 1

  allocate(explicit_grid%connections(2,num_connections_local))
  explicit_grid%connections = 0
  allocate(explicit_grid%face_areas(num_connections_local))
  explicit_grid%face_areas = 0    
  allocate(explicit_grid%face_centroids(num_connections_local))
  do iconn = 1, num_connections_local
    explicit_grid%face_centroids(iconn)%x = 0.d0
    explicit_grid%face_centroids(iconn)%y = 0.d0
    explicit_grid%face_centroids(iconn)%z = 0.d0
  enddo
        
  ! for now, read all cells from ASCII file through io_rank and communicate
  ! to other ranks
  if (option%myrank == option%io_rank) then
    allocate(temp_real_array(6,num_connections_local_save+1))
    ! read for other processors
    do irank = 0, option%mycommsize-1
      temp_real_array = -999.d0
      num_to_read = num_connections_local_save
      if (irank < remainder) num_to_read = num_to_read + 1
      do iconn = 1, num_to_read
        call InputReadFlotranString(input,option)
        call InputReadStringErrorMsg(input,option,hint)  
        call InputReadInt(input,option,temp_int)
        call InputErrorMsg(input,option,'cell id upwind',hint)
        temp_real_array(1,iconn) = dble(temp_int)
        call InputReadInt(input,option,temp_int)
        call InputErrorMsg(input,option,'cell id downwind',hint)
        temp_real_array(2,iconn) = dble(temp_int)
        call InputReadDouble(input,option,temp_real_array(3,iconn))
        call InputErrorMsg(input,option,'face x coordinate',hint)
        call InputReadDouble(input,option,temp_real_array(4,iconn))
        call InputErrorMsg(input,option,'face y coordinate',hint)
        call InputReadDouble(input,option,temp_real_array(5,iconn))
        call InputErrorMsg(input,option,'face z coordinate',hint)
        call InputReadDouble(input,option,temp_real_array(6,iconn))
        call InputErrorMsg(input,option,'face area',hint)
      enddo

      ! if the cells reside on io_rank
      if (irank == option%io_rank) then
#if UGRID_DEBUG
        write(string,*) num_connections_local
        string = trim(adjustl(string)) // ' connections stored on p0'
        print *, trim(string)
#endif
        do iconn = 1, num_connections_local
          explicit_grid%connections(1,iconn) = int(temp_real_array(1,iconn))
          explicit_grid%connections(2,iconn) = int(temp_real_array(2,iconn))
          explicit_grid%face_centroids(iconn)%x = temp_real_array(3,iconn)
          explicit_grid%face_centroids(iconn)%y = temp_real_array(4,iconn)
          explicit_grid%face_centroids(iconn)%z = temp_real_array(5,iconn)
          explicit_grid%face_areas(iconn) = temp_real_array(6,iconn)
        enddo
      else
        ! otherwise communicate to other ranks
#if UGRID_DEBUG
        write(string,*) num_to_read
        write(word,*) irank
        string = trim(adjustl(string)) // ' connections sent from p0 to p' // &
                 trim(adjustl(word))
        print *, trim(string)
#endif
        int_mpi = num_to_read*6
        call MPI_Send(temp_real_array,int_mpi,MPI_DOUBLE_PRECISION,irank, &
                      num_to_read,option%mycomm,ierr)
      endif
    enddo
  else
    ! other ranks post the recv
#if UGRID_DEBUG
    write(string,*) num_connections_local
    write(word,*) option%myrank
    string = trim(adjustl(string)) // ' connections received from p0 at p' // &
              trim(adjustl(word))
    print *, trim(string)
#endif
    allocate(temp_real_array(6,num_connections_local))
    int_mpi = num_connections_local*6
    call MPI_Recv(temp_real_array,int_mpi, &
                  MPI_DOUBLE_PRECISION,option%io_rank, &
                  MPI_ANY_TAG,option%mycomm,status_mpi,ierr)
    do iconn = 1, num_connections_local
      explicit_grid%connections(1,iconn) = int(temp_real_array(1,iconn))
      explicit_grid%connections(2,iconn) = int(temp_real_array(2,iconn))
      explicit_grid%face_centroids(iconn)%x = temp_real_array(3,iconn)
      explicit_grid%face_centroids(iconn)%y = temp_real_array(4,iconn)
      explicit_grid%face_centroids(iconn)%z = temp_real_array(5,iconn)
      explicit_grid%face_areas(iconn) = temp_real_array(6,iconn)
    enddo
    
  endif
  deallocate(temp_real_array)  
  
  if (option%myrank == option%io_rank) then
    call InputReadFlotranString(input,option)
    ! read ELEMENTS card, we only use this for tecplot output
    ! not used while solving the PDEs
    call InputReadWord(input,option,card,PETSC_TRUE)
    word = 'ELEMENTS'
    if (.not.StringCompare(word,card)) return
    card = 'Explicit Unstruct. Grid ELEMENTS'
    call InputReadInt(input,option,num_elems)
    call InputErrorMsg(input,option,'number of elements',card)
        explicit_grid%num_elems = num_elems
    allocate(explicit_grid%cell_connectivity(3,num_elems)) 
    explicit_grid%cell_connectivity = 0
    do iconn = 1, num_elems
      call InputReadFlotranString(input,option)
      call InputReadStringErrorMsg(input,option,card)  
      call InputReadInt(input,option, &
                        explicit_grid%cell_connectivity(1,iconn))
      call InputErrorMsg(input,option,'cell vertex 1',card)
      call InputReadInt(input,option, &
                        explicit_grid%cell_connectivity(2,iconn))
      call InputErrorMsg(input,option,'cell vertex 2',card)
      call InputReadInt(input,option, &
                        explicit_grid%cell_connectivity(3,iconn))
      call InputErrorMsg(input,option,'cell vertex 3',card)
    enddo
    call InputReadFlotranString(input,option)
    ! read VERTICES card, not used for calcuations, only tecplot output
    call InputReadWord(input,option,card,PETSC_TRUE)
    word = 'VERTICES'
    call InputErrorMsg(input,option,word,card)
    if (.not.StringCompare(word,card)) then
      option%io_buffer = 'Unrecognized keyword "' // trim(card) // &
        '" in explicit grid file.'
      call printErrMsgByRank(option)
    endif
    allocate(explicit_grid%vertex_coordinates(explicit_grid%num_cells_global))
    do icell = 1, explicit_grid%num_cells_global
      explicit_grid%vertex_coordinates(icell)%x = 0.d0
      explicit_grid%vertex_coordinates(icell)%y = 0.d0
      explicit_grid%vertex_coordinates(icell)%z = 0.d0
    enddo
    do icell = 1, explicit_grid%num_cells_global
      call InputReadFlotranString(input,option)
      call InputReadStringErrorMsg(input,option,card)  
      call InputReadDouble(input,option, &
                           explicit_grid%vertex_coordinates(icell)%x)
      call InputErrorMsg(input,option,'vertex 1',card)
      call InputReadDouble(input,option, &
                           explicit_grid%vertex_coordinates(icell)%y)
      call InputErrorMsg(input,option,'vertex 2',card)
      call InputReadDouble(input,option, &
                           explicit_grid%vertex_coordinates(icell)%z)
      call InputErrorMsg(input,option,'vertex 3',card)
    enddo
  endif
  
  if (option%myrank == option%io_rank) then
    call InputDestroy(input)
  endif
      
end subroutine ExplicitUGridReadInParallel

! ************************************************************************** !
!
! ExplicitUGridDecompose: Decomposes an explicit unstructured grid across 
!                         ranks
! author: Glenn Hammond
! date: 05/17/12
!
! ************************************************************************** !
subroutine ExplicitUGridDecompose(explicit_grid, num_ghost_cells, &
                                  global_offset, nmax, nlmax, ngmax, &
                                  cell_ids_natural, cell_ids_petsc, &
                                  ghost_cell_ids_petsc, ao_natural_to_petsc, &
                                  option)
  use Option_module
  use Utility_module, only: reallocateIntArray, SearchOrderedArray
  
  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"
  
  type(unstructured_explicit_type) :: explicit_grid
  PetscInt :: num_ghost_cells
  PetscInt :: global_offset
  PetscInt :: nmax, nlmax, ngmax
  PetscInt, pointer :: cell_ids_natural(:)
  PetscInt, pointer :: cell_ids_petsc(:)
  PetscInt, pointer :: ghost_cell_ids_petsc(:)
  AO :: ao_natural_to_petsc
  type(option_type) :: option
  
  PetscInt :: num_cells_local_new
  PetscInt :: global_offset_new
  PetscInt :: icell
  PetscInt, allocatable :: int_array(:)
  PetscErrorCode :: ierr

  
  num_cells_local_new = size(explicit_grid%cell_ids)
  global_offset_new = 0
  
  allocate(cell_ids_natural(num_cells_local_new))
  cell_ids_natural = explicit_grid%cell_ids

  ! make a list of petsc ids for each local cell (you simply take the global 
  ! offset and add it to the local contiguous cell ids on each processor
  allocate(int_array(num_cells_local_new))
  do icell = 1, num_cells_local_new
    int_array(icell) = icell+global_offset_new
  enddo
  
  ! make the arrays zero-based
  int_array = int_array - 1
  cell_ids_natural = cell_ids_natural - 1
  ! create an application ordering (mapping of natural to petsc ordering)
  call AOCreateBasic(option%mycomm,num_cells_local_new, &
                     cell_ids_natural,int_array, &
                     ao_natural_to_petsc,ierr)
  deallocate(int_array)
  ! make cell_ids_natural 1-based again
  cell_ids_natural = cell_ids_natural + 1
  
  allocate(cell_ids_petsc(num_cells_local_new))
  cell_ids_petsc = cell_ids_natural
  
  nmax = num_cells_local_new
  nlmax = nmax
  ngmax = nmax

  num_ghost_cells = 0
  global_offset = global_offset_new
  
end subroutine ExplicitUGridDecompose

! ************************************************************************** !
!
! ExplicitUGridDecomposeNew: Decomposes an explicit unstructured grid across 
!                         ranks
! author: Glenn Hammond
! date: 05/17/12
!
! ************************************************************************** !
subroutine ExplicitUGridDecomposeNew(ugrid,option)

  use Option_module
  use Utility_module, only: reallocateIntArray, SearchOrderedArray
  
  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"
  
  type(unstructured_grid_type) :: ugrid
  type(option_type) :: option

  type(unstructured_explicit_type), pointer :: explicit_grid
  PetscViewer :: viewer
  
  Mat :: M_mat,M_mat_loc
  Vec :: M_vec
  Mat :: Adj_mat
  Mat :: Dual_mat
  MatPartitioning :: Part
  IS :: is_new
  IS :: is_scatter  
  IS :: is_gather  
  PetscInt :: num_cells_local_new, num_cells_local_old
  Vec :: cells_old, cells_local
  Vec :: connections_old, connections_local
  VecScatter :: vec_scatter
  
  PetscInt :: global_offset_old
  PetscInt :: global_offset_new
  PetscInt :: ghosted_id
  PetscInt, allocatable :: local_connections(:), local_connection_offsets(:)
  PetscInt, allocatable :: int_array(:), int_array2(:), int_array3(:)
  PetscInt, allocatable :: int_array4(:)
  PetscInt, allocatable :: int_array2d(:,:)
  PetscInt :: num_connections_local_old, num_connections_local
  PetscInt :: num_connections_total 
  PetscInt :: num_connections_global, global_connection_offset
  PetscInt :: id_up, id_dn, iconn, icell, count, offset
  PetscInt :: conn_id, dual_id
  PetscBool :: found
  PetscInt :: i, temp_int, idual
  PetscReal :: temp_real
  
  PetscInt :: iflag
  PetscBool :: success
  PetscInt, pointer :: ia_ptr(:), ja_ptr(:)
  PetscInt, pointer :: ia_ptr2(:), ja_ptr2(:)
  PetscReal, pointer :: vec_ptr(:), vec_ptr2(:)
  PetscInt :: num_rows, num_cols, istart, iend, icol
  PetscInt :: cell_stride, dual_offset, connection_offset, connection_stride
  PetscInt :: natural_id_offset
  PetscErrorCode :: ierr
  PetscInt :: icell_up,icell_dn
  
  character(len=MAXSTRINGLENGTH) :: string

  explicit_grid => ugrid%explicit_grid
  
#if UGRID_DEBUG
  call printMsg(option,'Adjacency matrix')
#endif

  num_cells_local_old = size(explicit_grid%cell_ids)
  
  call MPI_Allreduce(num_cells_local_old,ugrid%nmax, &
                     ONE_INTEGER_MPI,MPIU_INTEGER, &
                     MPI_SUM,option%mycomm,ierr)  

  ! determine the global offset from 0 for cells on this rank
  global_offset_old = 0
  call MPI_Exscan(num_cells_local_old,global_offset_old, &
                  ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)  
  
  num_connections_local_old = size(explicit_grid%connections,2)
  
  call MPI_Allreduce(num_connections_local_old,num_connections_global, &
                     ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)

  global_connection_offset = 0
  call MPI_Exscan(num_connections_local_old,global_connection_offset, &
                  ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)

  call VecCreateMPI(option%mycomm,num_cells_local_old,ugrid%nmax, &   
                    M_vec,ierr)
  call VecZeroEntries(M_vec,ierr)
  do iconn = 1, num_connections_local_old
    do i = 1, 2
      icell = explicit_grid%connections(i,iconn)-1
      call VecSetValue(M_vec,icell,1.d0,ADD_VALUES,ierr)
    enddo
  enddo
  call VecAssemblyBegin(M_vec,ierr)
  call VecAssemblyEnd(M_vec,ierr)

#if UGRID_DEBUG
  call PetscViewerASCIIOpen(option%mycomm,'M_vec.out',viewer,ierr)
  call VecView(M_vec,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif

  call VecMax(M_vec,PETSC_NULL_INTEGER,temp_real,ierr)
  call VecDestroy(M_vec,ierr)
  ugrid%max_ndual_per_cell = int(temp_real+0.1d0)
  call MatCreateAIJ(option%mycomm,num_cells_local_old,PETSC_DECIDE, &
                    ugrid%nmax,num_connections_global, &
                    ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
                    ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
                    M_mat,ierr)
  call MatZeroEntries(M_mat,ierr)
  do iconn = 1, num_connections_local_old
    temp_int = iconn+global_connection_offset-1
    do i = 1, 2
      icell = explicit_grid%connections(i,iconn)-1
      call MatSetValue(M_mat,icell,temp_int,1.d0,INSERT_VALUES,ierr)
    enddo
  enddo
  call MatAssemblyBegin(M_mat,MAT_FINAL_ASSEMBLY,ierr)
  call MatAssemblyEnd(M_mat,MAT_FINAL_ASSEMBLY,ierr)

#if UGRID_DEBUG
  call PetscViewerASCIIOpen(option%mycomm,'M_mat.out',viewer,ierr)
  call MatView(M_mat,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif

  ! GB: When MatConvert() is used, the diagonal entries are lost in Adj_mat
  !call MatConvert(M_mat,MATMPIADJ,MAT_INITIAL_MATRIX,Adj_mat,ierr)
  !call MatDestroy(M_mat,ierr)

  ! Alternate method of creating Adj_mat
  if (option%mycommsize>1) then
    call MatMPIAIJGetLocalMat(M_mat,MAT_INITIAL_MATRIX,M_mat_loc,ierr)
    call MatGetRowIJF90(M_mat_loc,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
                        ia_ptr,ja_ptr,success,ierr)
  else
    call MatGetRowIJF90(M_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
                        ia_ptr,ja_ptr,success,ierr)
  endif

  count=0
  do icell = 1,num_rows
    istart = ia_ptr(icell)
    iend = ia_ptr(icell+1)-1
    num_cols = iend-istart+1
    count = count+num_cols
  enddo
  allocate(local_connections(count))
  allocate(local_connection_offsets(num_rows+1))
  local_connection_offsets(1:num_rows+1) = ia_ptr(1:num_rows+1)
  local_connections(1:count)             = ja_ptr(1:count)

  call MPI_Barrier(MPI_COMM_WORLD,ierr)

  call MatCreateMPIAdj(option%mycomm,num_cells_local_old, &
                       num_connections_global, &
                       local_connection_offsets, &
                       local_connections,PETSC_NULL_INTEGER,Adj_mat,ierr)

  if (option%mycommsize>1) then
    call MatRestoreRowIJF90(M_mat_loc,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
                        ia_ptr,ja_ptr,success,ierr)
  else
    call MatRestoreRowIJF90(M_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
                        ia_ptr,ja_ptr,success,ierr)
  endif
  call MatDestroy(M_mat,ierr)

#if UGRID_DEBUG
  call PetscViewerASCIIOpen(option%mycomm,'Adj.out',viewer,ierr)
  call MatView(Adj_mat,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif

  ! Create the Dual matrix.
  call MatCreateAIJ(option%mycomm,num_cells_local_old,PETSC_DECIDE, &
                    ugrid%nmax,ugrid%nmax, &
                    ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
                    ugrid%max_ndual_per_cell,PETSC_NULL_INTEGER, &
                    M_mat,ierr)
  do iconn = 1, num_connections_local_old
    icell_up = explicit_grid%connections(1,iconn)-1
    icell_dn = explicit_grid%connections(2,iconn)-1
    call MatSetValue(M_mat,icell_up,icell_dn,1.d0,INSERT_VALUES,ierr)
    call MatSetValue(M_mat,icell_dn,icell_up,1.d0,INSERT_VALUES,ierr)
  enddo
  
  call MatAssemblyBegin(M_mat,MAT_FINAL_ASSEMBLY,ierr)
  call MatAssemblyEnd(M_mat,MAT_FINAL_ASSEMBLY,ierr)
  
  call MatConvert(M_mat,MATMPIADJ,MAT_INITIAL_MATRIX,Dual_mat,ierr)
  call MatDestroy(M_mat,ierr)

#if UGRID_DEBUG
  call PetscViewerASCIIOpen(option%mycomm,'Dual_mat.out',viewer,ierr)
  call MatView(Dual_mat,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif

  call UGridPartition(ugrid,option,Dual_mat,is_new, &
                      num_cells_local_new)

  ! second argument of ZERO_INTEGER means to use 0-based indexing
  ! MagGetRowIJF90 returns row and column pointers for compressed matrix data
  call MatGetRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
                      ia_ptr,ja_ptr,success,ierr)

  if (.not.success .or. num_rows /= num_cells_local_old) then
    print *, option%myrank, num_rows, success, num_cells_local_old
    option%io_buffer = 'Error getting IJ row indices from dual matrix'
    call printErrMsg(option)
  endif

  call MatRestoreRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
                          num_rows,ia_ptr,ja_ptr,success,ierr)
  
  ! in order to redistributed cell/connection data among ranks, I package it
  ! in a crude way within a strided petsc vec and pass it.  The stride 
  ! determines the size of each cells "packaged" data 
  connection_offset = 6 + 1 ! +1 for -777
  dual_offset = connection_offset + ugrid%max_ndual_per_cell + 1 ! +1 for -888
  cell_stride = dual_offset + ugrid%max_ndual_per_cell + 1 ! +1 for -999999
  natural_id_offset = 2

  ! Information for each cell is packed in a strided petsc vec
  ! The information is ordered within each stride as follows:
  ! -cell_N   ! global cell id (negative indicates 1-based)
  ! natural cell id
  ! cell x coordinate
  ! cell y coordinate
  ! cell z coordinate
  ! cell volume
  ! -777      ! separator between cell info and connection info
  ! conn1     ! connection ids between cell_N and others
  ! conn1
  ! ...
  ! connN     
  ! -888      ! separator between connection info and dual ids
  ! dual1     ! dual ids between cell_N and others
  ! dual2
  ! ...
  ! dualN     
  ! -999999   ! separator indicating end of information for cell_N
  
  ! the purpose of -888, and -999999 is to allow one to use cells of 
  ! various geometry.  
  
  call UGridCreateOldVec(ugrid,option,cells_old, &
                         num_cells_local_old, &
                         is_new,is_scatter,cell_stride)

  ! 0 = 0-based indexing
  ! MagGetRowIJF90 returns row and column pointers for compressed matrix data
  ! pointers to Dual mat
  call MatGetRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,num_rows, &
                      ia_ptr,ja_ptr,success,ierr)
  ! pointers to Adj mat
  call MatGetRowIJF90(Adj_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE,temp_int, &
                      ia_ptr2,ja_ptr2,success,ierr)
  
  if (num_rows /= temp_int) then
    write(string,*) num_rows, temp_int
    option%io_buffer = 'Number of rows in Adj and Dual matrices inconsistent:'
    option%io_buffer = trim(option%io_buffer) // trim(adjustl(string))
    call printErrMsgByRank(option)
  endif

  call VecGetArrayF90(cells_old,vec_ptr,ierr)
  count = 0
  do icell = 1, num_cells_local_old
    count = count + 1
    ! set global cell id
    ! negate to indicate cell id with 1-based numbering (-0 = 0)
    vec_ptr(count) = -(global_offset_old+icell)
    count = count + 1
    vec_ptr(count) = explicit_grid%cell_ids(icell)
    count = count + 1
    vec_ptr(count) = explicit_grid%cell_centroids(icell)%x
    count = count + 1
    vec_ptr(count) = explicit_grid%cell_centroids(icell)%y
    count = count + 1
    vec_ptr(count) = explicit_grid%cell_centroids(icell)%z
    count = count + 1
    vec_ptr(count) = explicit_grid%cell_volumes(icell)
    ! add the separator
    count = count + 1
    vec_ptr(count) = -777  ! help differentiate
 
    ! add the connection ids
    istart = ia_ptr2(icell)
    iend = ia_ptr2(icell+1)-1
    num_cols = iend-istart+1
    if (num_cols > ugrid%max_ndual_per_cell) then
      option%io_buffer = &
        'Number of columns in Adj matrix is larger then max_ndual_per_cell.'
      call printErrMsgByRank(option)
    endif
    do icol = 1, ugrid%max_ndual_per_cell
      count = count + 1
      if (icol <= num_cols) then
        ! increment for 1-based ordering
        vec_ptr(count) = ja_ptr2(icol+istart) + 1
      else
        vec_ptr(count) = 0
      endif
    enddo
    
    ! add the separator
    count = count + 1
    vec_ptr(count) = -888  ! help differentiate

    ! add the dual ids
    istart = ia_ptr(icell)
    iend = ia_ptr(icell+1)-1
    num_cols = iend-istart+1
    if (num_cols > ugrid%max_ndual_per_cell) then
      option%io_buffer = &
        'Number of columns in Dual matrix is larger then max_ndual_per_cell.'
      call printErrMsgByRank(option)
    endif
    do icol = 1, ugrid%max_ndual_per_cell
      count = count + 1
      if (icol <= num_cols) then
        ! increment for 1-based ordering
        vec_ptr(count) = ja_ptr(icol+istart) + 1
      else
        vec_ptr(count) = 0
      endif
    enddo

    ! final separator
    count = count + 1 
    vec_ptr(count) = -999999  ! help differentiate
  enddo
  call VecRestoreArrayF90(cells_old,vec_ptr,ierr)
  
  ! pointers to Dual mat
  call MatRestoreRowIJF90(Dual_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
                          num_rows,ia_ptr,ja_ptr,success,ierr)
  ! pointers to Adj mat
  call MatRestoreRowIJF90(Adj_mat,ZERO_INTEGER,PETSC_FALSE,PETSC_FALSE, &
                          temp_int,ia_ptr2,ja_ptr2,success,ierr)
  call MatDestroy(Dual_mat,ierr)
  call MatDestroy(Adj_mat,ierr)
  
  ! is_scatter is destroyed within UGridNaturalToPetsc
  call UGridNaturalToPetsc(ugrid,option, &
                           cells_old,cells_local, &
                           num_cells_local_new,cell_stride,dual_offset, &
                           natural_id_offset,is_scatter)
  
  call VecDestroy(cells_old,ierr)

  ! set up connections
  connection_stride = 8
  ! create strided vector with the old connection distribution
  call VecCreate(option%mycomm,connections_old,ierr)
  call VecSetSizes(connections_old, &
                   connection_stride*num_connections_local_old, &
                   PETSC_DECIDE,ierr)
  call VecSetFromOptions(connections_old,ierr) 

  call VecGetArrayF90(connections_old,vec_ptr,ierr)
  do iconn = 1, num_connections_local_old
    offset = (iconn-1)*connection_stride
    vec_ptr(offset+1) = explicit_grid%connections(1,iconn)
    vec_ptr(offset+2) = explicit_grid%connections(2,iconn)
    vec_ptr(offset+3) = explicit_grid%face_centroids(iconn)%x
    vec_ptr(offset+4) = explicit_grid%face_centroids(iconn)%y
    vec_ptr(offset+5) = explicit_grid%face_centroids(iconn)%z
    vec_ptr(offset+6) = explicit_grid%face_areas(iconn)
    vec_ptr(offset+7) = 1.d0 ! flag for local connections
    vec_ptr(offset+8) = -888.d0
  enddo
  call VecRestoreArrayF90(connections_old,vec_ptr,ierr)

    
#if UGRID_DEBUG
  call PetscViewerASCIIOpen(option%mycomm,'connections_old.out',viewer,ierr)
  call VecView(connections_old,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif    
  
  ! count the number of cells and their duals  
  call VecGetArrayF90(cells_local,vec_ptr,ierr)
  count = 0
  do ghosted_id=1, ugrid%ngmax
    do iconn = 1, ugrid%max_ndual_per_cell
      conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
      if (conn_id < 1) exit ! here we hit the 0 at the end of last connection
      ! yes, we will be counting them twice
      count = count + 1
    enddo
  enddo   
  num_connections_total = count ! many of these are redundant and will be removed
  ! allocate and fill an array with the natural cell and dual ids
  allocate(int_array(num_connections_total))
  count = 0
  do ghosted_id=1, ugrid%ngmax
    do iconn = 1, ugrid%max_ndual_per_cell
      conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
      if (conn_id < 1) exit ! again we hit the 0 
      count = count + 1
      int_array(count) = conn_id
    enddo
  enddo
  call VecRestoreArrayF90(cells_local,vec_ptr,ierr)
  
  allocate(int_array2(num_connections_total))
  do iconn = 1, num_connections_total
    int_array2(iconn) = iconn
  enddo
  
  ! sort connections - int_array2 will return the reordering while int_array 
  !                    remains the same.
  int_array2 = int_array2 - 1
  call PetscSortIntWithPermutation(num_connections_total,int_array, &
                                   int_array2,ierr)
  int_array2 = int_array2 + 1
  
  ! permute, remove duplicate connections, and renumber to local ordering
  allocate(int_array3(num_connections_total))
  allocate(int_array4(num_connections_total))
  int_array3 = -999
  int_array4 = -999
  int_array3(1) = int_array(int_array2(1))
  count = 1
  do iconn = 1, num_connections_total
    if (int_array(int_array2(iconn)) > int_array3(count)) then
      count = count + 1
      int_array3(count) = int_array(int_array2(iconn))
    endif
    int_array4(int_array2(iconn)) = count
  enddo
  deallocate(int_array)
  deallocate(int_array2)
  
  num_connections_local = count ! new compressed count
  allocate(int_array(num_connections_local))
  int_array = int_array3(1:num_connections_local)
  deallocate(int_array3)
  
  ! replace original connections ids (naturally numbered) with locally 
  ! numbered connection ids (int_array4)
  call VecGetArrayF90(cells_local,vec_ptr,ierr)
  count = 0
  do ghosted_id=1, ugrid%ngmax
    do iconn = 1, ugrid%max_ndual_per_cell
      conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
      if (conn_id < 1) exit ! again we hit the 0 
      count = count + 1
      vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride) = &
        int_array4(count)
    enddo
  enddo
  deallocate(int_array4)
  call VecRestoreArrayF90(cells_local,vec_ptr,ierr)
  
  ! check to ensure that the number before/after are consistent
  if (count /= num_connections_total) then
    write(string,'(2i6)') count, num_connections_total
    option%io_buffer = 'Inconsistent values for num_connections_total: ' // &
      trim(adjustl(string))
    call printErrMsgByRank(option)
  endif
  num_connections_total = -999 ! set to uninitialized value to catch bugs
  
  call VecCreate(PETSC_COMM_SELF,connections_local,ierr)
  call VecSetSizes(connections_local,num_connections_local*connection_stride, &
                   PETSC_DECIDE,ierr)
  call VecSetBlockSize(connections_local,connection_stride,ierr)
  call VecSetFromOptions(connections_local,ierr)

  int_array = int_array-1
  call ISCreateBlock(option%mycomm,connection_stride,num_connections_local, &
                     int_array,PETSC_COPY_VALUES,is_scatter,ierr)
  do iconn = 1, num_connections_local
    int_array(iconn) = iconn-1
  enddo
  call ISCreateBlock(option%mycomm,connection_stride,num_connections_local, &
                     int_array,PETSC_COPY_VALUES,is_gather,ierr)
  deallocate(int_array)

#if UGRID_DEBUG
  call PetscViewerASCIIOpen(option%mycomm,'is_scatter_conn_old_to_local.out',viewer,ierr)
  call ISView(is_scatter,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
  call PetscViewerASCIIOpen(option%mycomm,'is_gather_conn_old_to_local.out',viewer,ierr)
  call ISView(is_gather,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
  call printMsg(option,'Scatter/gathering local connection info')
#endif  
  
  ! scatter all the connection data from the old to local
  call VecScatterCreate(connections_old,is_scatter,connections_local, &
                        is_gather,vec_scatter,ierr)
  call ISDestroy(is_gather,ierr)
  call ISDestroy(is_scatter,ierr)
  call VecScatterBegin(vec_scatter,connections_old,connections_local, &
                       INSERT_VALUES,SCATTER_FORWARD,ierr)
  call VecScatterEnd(vec_scatter,connections_old,connections_local, &
                     INSERT_VALUES,SCATTER_FORWARD,ierr)
  call VecScatterDestroy(vec_scatter,ierr)  

    
#if UGRID_DEBUG
  write(string,*) option%myrank
  string = 'connections_local_nat' // trim(adjustl(string)) // '.out'
  call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer,ierr)
  call VecView(connections_local,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif   
  
  ! loop over cells and change the natural ids in the duals to local ids
  allocate(int_array2d(2,num_connections_local))
  int_array2d = -999
  call VecGetArrayF90(cells_local,vec_ptr,ierr)
  do ghosted_id=1, ugrid%ngmax
    do iconn = 1, ugrid%max_ndual_per_cell
      ! this connection id is now local
      conn_id = int(vec_ptr(iconn + connection_offset + (ghosted_id-1)*cell_stride))
      if (conn_id < 1) exit ! again we hit the 0
      do i = 1, 2
        if (int_array2d(i,conn_id) <= -999) then
          int_array2d(i,conn_id) = ghosted_id
          exit
        endif
        if (i > 2) then
          write(string,'(2i5)') ghosted_id, conn_id
          option%io_buffer = 'Too many local cells match connection: ' // &
            trim(adjustl(string))
          call printErrMsgByRank(option)
        endif
      enddo
    enddo
  enddo
  call VecRestoreArrayF90(cells_local,vec_ptr,ierr)

  ! map natural ids in connections to local ids
  ! negate connection ids as a flag
  int_array2d = -1*int_array2d
  call VecGetArrayF90(connections_local,vec_ptr,ierr)
  do iconn = 1, num_connections_local
    offset = connection_stride*(iconn-1)
    ! all values should be negative at this point, unless uninitialized
    if (maxval(int_array2d(:,iconn)) >= 999) then
      ! connection is between two ghosted cells
      vec_ptr(offset+7) = 0.d0
      cycle
    endif
    id_up = int(vec_ptr(offset+1)) ! this is the natural id
    id_dn = int(vec_ptr(offset+2))
    count = 0
    found = PETSC_FALSE
    do i = 1, 2
      if (ugrid%cell_ids_natural(abs(int_array2d(i,iconn))) == id_up) then
        int_array2d(i,iconn) = abs(int_array2d(i,iconn))
        found = PETSC_TRUE
        count = count + 1
      endif
      if (ugrid%cell_ids_natural(abs(int_array2d(i,iconn))) == id_dn) then
        int_array2d(i,iconn) = abs(int_array2d(i,iconn))
        found = PETSC_TRUE
        count = count - 1
      endif
    enddo
    ! count should be zero, meaning it found the upwind and downwind cell
    ! ids
    if (count /= 0 .or. .not.found) then
      write(string,*) iconn, id_up, id_dn
      if (.not.found) then
        option%io_buffer = 'upwind/downwind cells not found: '
      else if (count < 0) then
        option%io_buffer = 'upwind cell not found: '
      else
        option%io_buffer = 'downwind cell not found: '
      endif
      option%io_buffer = trim(option%io_buffer) // trim(adjustl(string))
      call printErrMsgByRank(option)
    endif
    id_up = int_array2d(1,iconn)
    id_dn = int_array2d(2,iconn)
    if (id_up < id_dn) then
      vec_ptr(offset+1) = id_up  ! now local ids
      vec_ptr(offset+2) = id_dn 
    else
      vec_ptr(offset+1) = id_dn
      vec_ptr(offset+2) = id_up 
    endif
    if (id_up > ugrid%nlmax .and. id_dn > ugrid%nlmax) then
      ! connection is between two ghosted cells
      vec_ptr(offset+7) = 0.d0
    endif
  enddo
  call VecRestoreArrayF90(connections_local,vec_ptr,ierr)
  deallocate(int_array2d)

#if UGRID_DEBUG
  write(string,*) option%myrank
  string = 'connections_local_loc' // trim(adjustl(string)) // '.out'
  call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer,ierr)
  call VecView(connections_local,viewer,ierr)
  call PetscViewerDestroy(viewer,ierr)
#endif   

  ! deallocate/allocate grid cell info locally
  deallocate(explicit_grid%cell_ids)
  deallocate(explicit_grid%cell_volumes)
  deallocate(explicit_grid%cell_centroids)

  allocate(explicit_grid%cell_ids(ugrid%ngmax))
  explicit_grid%cell_ids = -999
  allocate(explicit_grid%cell_volumes(ugrid%ngmax))
  explicit_grid%cell_volumes = -999.d0
  allocate(explicit_grid%cell_centroids(ugrid%ngmax))
  do icell = 1, ugrid%ngmax
    explicit_grid%cell_centroids(icell)%x = -999.d0
    explicit_grid%cell_centroids(icell)%y = -999.d0
    explicit_grid%cell_centroids(icell)%z = -999.d0
  enddo

  call VecGetArrayF90(cells_local,vec_ptr,ierr)
  do ghosted_id=1, ugrid%ngmax
    offset = cell_stride*(ghosted_id-1)
    explicit_grid%cell_ids(ghosted_id) = int(vec_ptr(offset + 2))
    explicit_grid%cell_centroids(ghosted_id)%x = vec_ptr(offset + 3)
    explicit_grid%cell_centroids(ghosted_id)%y = vec_ptr(offset + 4)
    explicit_grid%cell_centroids(ghosted_id)%z = vec_ptr(offset + 5)
    explicit_grid%cell_volumes(ghosted_id) = vec_ptr(offset + 6)
  enddo
  call VecRestoreArrayF90(cells_local,vec_ptr,ierr)

#if UGRID_DEBUG
  write(string,*) option%myrank
  string = 'cells_local_raw' // trim(adjustl(string)) // '.out'
  open(unit=86,file=trim(string))
  do ghosted_id = 1, ugrid%ngmax
    write(86,'(i5,4f10.3)') explicit_grid%cell_ids(ghosted_id), &
                explicit_grid%cell_centroids(ghosted_id)%x, &
                explicit_grid%cell_centroids(ghosted_id)%y, &
                explicit_grid%cell_centroids(ghosted_id)%z, &
                explicit_grid%cell_volumes(ghosted_id)
  enddo
  close(86)
#endif     
  
  ! deallocate/allocate connection info locally
  deallocate(explicit_grid%connections)
  deallocate(explicit_grid%face_areas)
  deallocate(explicit_grid%face_centroids)

  count = 0
  call VecGetArrayF90(connections_local,vec_ptr,ierr)
  do iconn = 1, num_connections_local
    offset = connection_stride*(iconn-1)
    if (vec_ptr(offset+7) > 0.1d0) count = count + 1
  enddo
  call VecRestoreArrayF90(connections_local,vec_ptr,ierr)

  allocate(explicit_grid%connections(2,count))
  explicit_grid%connections = 0
  allocate(explicit_grid%face_areas(count))
  explicit_grid%face_areas = 0    
  allocate(explicit_grid%face_centroids(count))
  do iconn = 1, count
    explicit_grid%face_centroids(iconn)%x = 0.d0
    explicit_grid%face_centroids(iconn)%y = 0.d0
    explicit_grid%face_centroids(iconn)%z = 0.d0
  enddo  
  call VecGetArrayF90(connections_local,vec_ptr,ierr)
  count = 0
  do iconn = 1, num_connections_local
    offset = connection_stride*(iconn-1)
    if (vec_ptr(offset+7) > 0.1d0) then
      count = count + 1
      explicit_grid%connections(1,count) = int(vec_ptr(offset+1))
      explicit_grid%connections(2,count) = int(vec_ptr(offset+2))
      explicit_grid%face_centroids(count)%x = vec_ptr(offset+3)
      explicit_grid%face_centroids(count)%y = vec_ptr(offset+4)
      explicit_grid%face_centroids(count)%z = vec_ptr(offset+5)
      explicit_grid%face_areas(count) = vec_ptr(offset+6)
    endif
  enddo
  call VecRestoreArrayF90(connections_local,vec_ptr,ierr)
  num_connections_local = count

#if UGRID_DEBUG
  write(string,*) option%myrank
  string = 'connections_local_raw' // trim(adjustl(string)) // '.out'
  open(unit=86,file=trim(string))
  do iconn = 1, num_connections_local
    write(86,'(2i5,4f7.3)') explicit_grid%connections(1,iconn), &
                explicit_grid%connections(2,iconn), &
                explicit_grid%face_centroids(iconn)%x, &
                explicit_grid%face_centroids(iconn)%y, &
                explicit_grid%face_centroids(iconn)%z, &
                explicit_grid%face_areas(iconn)
  enddo
  close(86)
#endif     
  
  call VecDestroy(connections_old,ierr)
  call VecDestroy(connections_local,ierr)
  call VecDestroy(cells_local,ierr)
  
end subroutine ExplicitUGridDecomposeNew

! ************************************************************************** !
!
! ExplicitUGridSetCellCentroids: Sets the centroid of each grid cell
! author: Glenn Hammond
! date: 05/17/12
!
! ************************************************************************** !
subroutine ExplicitUGridSetCellCentroids(explicit_grid,x,y,z, &
                                         x_min,x_max,y_min,y_max,z_min,z_max)

  use Option_module

  implicit none
  
  type(unstructured_explicit_type) :: explicit_grid
  PetscReal :: x(:), y(:), z(:)
  PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max

  PetscInt :: icell
  
  do icell = 1, size(explicit_grid%cell_centroids)
    x(icell) = explicit_grid%cell_centroids(icell)%x
    y(icell) = explicit_grid%cell_centroids(icell)%y
    z(icell) = explicit_grid%cell_centroids(icell)%z
  enddo
  
  x_min = minval(x)
  x_max = maxval(x)
  y_min = minval(y)
  y_max = maxval(y)
  z_min = minval(z)
  z_max = maxval(z)
      
end subroutine ExplicitUGridSetCellCentroids
      
! ************************************************************************** !
!
! ExplicitUGridSetInternConnect: Sets up the internal connectivity within  
!                                the connectivity object
! author: Glenn Hammond
! date: 05/17/12
!
! ************************************************************************** !
function ExplicitUGridSetInternConnect(explicit_grid,option)

  use Utility_module
  use Connection_module
  use Option_module

  implicit none
  
  type(connection_set_type), pointer :: ExplicitUGridSetInternConnect
  
  type(unstructured_explicit_type) :: explicit_grid
  type(option_type) :: option
  
  type(connection_set_type), pointer :: connections
  PetscInt :: num_connections
  PetscInt :: iconn
  PetscInt :: id_up, id_dn
  PetscReal :: v(3), v_up(3), v_dn(3)
  PetscReal :: distance
  
  
  num_connections = size(explicit_grid%connections,2)
  connections => ConnectionCreate(num_connections,INTERNAL_CONNECTION_TYPE)
  
  do iconn = 1, num_connections
    id_up = explicit_grid%connections(1,iconn)
    id_dn = explicit_grid%connections(2,iconn)
    connections%id_up(iconn) = id_up
    connections%id_dn(iconn) = id_dn
    
    v_up(1) = explicit_grid%face_centroids(iconn)%x - &
              explicit_grid%cell_centroids(id_up)%x
    v_up(2) = explicit_grid%face_centroids(iconn)%y - &
              explicit_grid%cell_centroids(id_up)%y
    v_up(3) = explicit_grid%face_centroids(iconn)%z - &
              explicit_grid%cell_centroids(id_up)%z

    v_dn(1) = explicit_grid%cell_centroids(id_dn)%x - &
              explicit_grid%face_centroids(iconn)%x
    v_dn(2) = explicit_grid%cell_centroids(id_dn)%y - &
              explicit_grid%face_centroids(iconn)%y
    v_dn(3) = explicit_grid%cell_centroids(id_dn)%z - &
              explicit_grid%face_centroids(iconn)%z

    v = v_up + v_dn
    distance = sqrt(DotProduct(v,v))
    connections%dist(-1,iconn) = sqrt(DotProduct(v_up,v_up))/distance
    connections%dist(0,iconn) = distance
    connections%dist(1:3,iconn) = v/distance
    connections%area(iconn) = explicit_grid%face_areas(iconn)
  enddo
  
  ExplicitUGridSetInternConnect => connections

end function ExplicitUGridSetInternConnect

! ************************************************************************** !
!
! ExplicitUGridComputeVolumes: Sets the volume of each grid cell
! author: Glenn Hammond
! date: 05/17/12
!
! ************************************************************************** !
subroutine ExplicitUGridComputeVolumes(ugrid,option,volume)

  use Option_module

  implicit none
  
  type(unstructured_grid_type) :: ugrid
  type(option_type) :: option
  Vec :: volume

  type(unstructured_explicit_type), pointer :: explicit_grid
  
  PetscInt :: icell
  PetscReal, pointer :: vec_ptr(:)
  PetscErrorCode :: ierr

  explicit_grid => ugrid%explicit_grid

  call VecGetArrayF90(volume,vec_ptr,ierr)
  do icell = 1, ugrid%nlmax
    vec_ptr(icell) = explicit_grid%cell_volumes(icell)
  enddo
  call VecRestoreArrayF90(volume,vec_ptr,ierr)
  
end subroutine ExplicitUGridComputeVolumes

! ************************************************************************** !
!
! ExplicitUGridSetBoundaryConnect: Sets up the boundary connectivity within  
!                                  the connectivity object
! author: Glenn Hammond
! date: 05/18/12
!
! ************************************************************************** !
function ExplicitUGridSetBoundaryConnect(explicit_grid,cell_ids, &
                                         face_centroids,face_areas,option)

  use Utility_module
  use Connection_module
  use Option_module

  implicit none
  
  type(connection_set_type), pointer :: ExplicitUGridSetBoundaryConnect

  type(unstructured_explicit_type) :: explicit_grid
  PetscInt :: cell_ids(:)
  type(point3d_type) :: face_centroids(:)
  PetscReal :: face_areas(:)
  type(option_type) :: option
  
  type(connection_set_type), pointer :: connections
  PetscInt :: num_connections
  PetscInt :: iconn
  PetscInt :: id
  PetscReal :: v(3)
  PetscReal :: distance
  
  
  num_connections = size(cell_ids)
  connections => ConnectionCreate(num_connections,BOUNDARY_CONNECTION_TYPE)
  
  do iconn = 1, num_connections
    id = cell_ids(iconn)
    connections%id_dn(iconn) = id
    
    v(1) = explicit_grid%cell_centroids(id)%x - &
           face_centroids(iconn)%x
    v(2) = explicit_grid%cell_centroids(id)%y - &
           face_centroids(iconn)%y
    v(3) = explicit_grid%cell_centroids(id)%z - &
           face_centroids(iconn)%z

    distance = sqrt(DotProduct(v,v))
    connections%dist(-1,iconn) = 0.d0
    connections%dist(0,iconn) = distance
    connections%dist(1:3,iconn) = v/distance
    connections%area(iconn) = face_areas(iconn)
  enddo
  
  ExplicitUGridSetBoundaryConnect => connections

end function ExplicitUGridSetBoundaryConnect

! ************************************************************************** !
!
! ExplicitUGridSetConnections: Sets up the connectivity for a region
! author: Glenn Hammond
! date: 05/18/12
!
! ************************************************************************** !
function ExplicitUGridSetConnections(explicit_grid,cell_ids,connection_type, &
                                     option)

  use Utility_module
  use Connection_module
  use Option_module

  implicit none
  
  type(connection_set_type), pointer :: ExplicitUGridSetConnections

  type(unstructured_explicit_type) :: explicit_grid
  PetscInt :: cell_ids(:)
  PetscInt :: connection_type
  type(option_type) :: option
  
  type(connection_set_type), pointer :: connections
  PetscInt :: num_connections
  PetscInt :: iconn
  PetscInt :: id
  
  num_connections = size(cell_ids)
  connections => ConnectionCreate(num_connections,connection_type)
  
  do iconn = 1, num_connections
    id = cell_ids(iconn)
    connections%id_dn(iconn) = id
  enddo
  
  ExplicitUGridSetConnections => connections

end function ExplicitUGridSetConnections

end module Unstructured_Explicit_module