Source

pflotran-dev / src / pflotran / condition_control.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
module Condition_Control_module
  
  ! This module store routines that operate on conditions from a level above 
  ! that of the realization_module.  This is necessary to access capability
  ! such as HDF5 which is unavailable from within the realization object
  ! and below.  Routines in this module will loop over realization, levels,
  ! and patches without calling underlying level/patch versions of the
  ! subroutines, which is common in realization.F90 - GEH
 
  implicit none

  private

#include "definitions.h"

  public :: CondControlAssignFlowInitCond, &
            CondControlAssignTranInitCond, &
#ifdef SURFACE_FLOW
            CondControlAssignFlowInitCondSurface, &
#endif
            CondControlScaleSourceSink
 
contains

! ************************************************************************** !
!
! CondControlAssignFlowInitCond: Assigns flow initial conditions to model
! author: Glenn Hammond
! date: 11/02/07, 10/18/11
!
! ************************************************************************** !
subroutine CondControlAssignFlowInitCond(realization)

  use Realization_class
  use Discretization_module
  use Region_module
  use Option_module
  use Field_module
  use Coupler_module
  use Condition_module
  use Dataset_Aux_module
  use Grid_module
  use Level_module
  use Patch_module
  use Water_EOS_module
#ifdef DASVYAT
  use MFD_module, only : MFDInitializeMassMatrices
#endif

  use Global_Aux_module
  use General_Aux_module
  
  implicit none

#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"
  
  type(realization_type) :: realization
  
  PetscInt :: icell, iconn, idof, iface
  PetscInt :: local_id, ghosted_id, iend, ibegin
  PetscReal, pointer :: xx_p(:), iphase_loc_p(:), xx_faces_p(:)
  PetscErrorCode :: ierr
  
  character(len=MAXSTRINGLENGTH) :: string
  
  type(option_type), pointer :: option
  type(field_type), pointer :: field  
  type(patch_type), pointer :: patch
  type(grid_type), pointer :: grid
  type(discretization_type), pointer :: discretization
  type(coupler_type), pointer :: initial_condition
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  type(flow_general_condition_type), pointer :: general
  type(dataset_type), pointer :: dataset
  type(global_auxvar_type) :: global_aux
  type(general_auxvar_type) :: general_aux
  PetscBool :: use_dataset
  PetscBool :: dataset_flag(realization%option%nflowdof)
  PetscInt :: num_connections
  PetscInt, pointer :: conn_id_ptr(:)
  PetscInt :: ghosted_offset
  PetscReal :: x(realization%option%nflowdof)
  PetscReal :: temperature, p_sat

  option => realization%option
  discretization => realization%discretization
  field => realization%field
  patch => realization%patch

  if (option%iflowmode == G_MODE) then
    call GlobalAuxVarInit(global_aux,option)
    call GeneralAuxVarInit(general_aux,option)
  endif
  
  cur_level => realization%level_list%first
  do 
    if (.not.associated(cur_level)) exit
    cur_patch => cur_level%patch_list%first
    do
      if (.not.associated(cur_patch)) exit

      grid => cur_patch%grid

      select case(option%iflowmode)
      
        case(G_MODE) ! general phase mode

          call GridVecGetArrayF90(grid,field%flow_xx,xx_p, ierr); CHKERRQ(ierr)
          call GridVecGetArrayF90(grid,field%iphas_loc,iphase_loc_p,ierr)
      
          xx_p = -999.d0
      
          initial_condition => cur_patch%initial_conditions%first
          do
      
            if (.not.associated(initial_condition)) exit

            if (.not.associated(initial_condition%flow_aux_real_var)) then
              if (.not.associated(initial_condition%flow_condition)) then
                option%io_buffer = 'Flow condition is NULL in initial condition'
                call printErrMsg(option)
              endif
              
              general => initial_condition%flow_condition%general
              
              string = 'in flow condition "' // &
                trim(initial_condition%flow_condition%name) // &
                '" within initial condition "' // &
                trim(initial_condition%flow_condition%name) // &
                '" must be of type Dirichlet or Hydrostatic'
              ! error checking.  the data must match the state
              select case(initial_condition%flow_condition%iphase)
                case(TWO_PHASE_STATE)  
                  if (.not. &
                      (general%gas_pressure%itype == DIRICHLET_BC .or. &
                       general%gas_pressure%itype == HYDROSTATIC_BC)) then
                    option%io_buffer = 'Gas pressure ' // trim(string)
                    call printErrMsg(option)
                  endif
                  if (.not. &
                      (general%gas_saturation%itype == DIRICHLET_BC .or. &
                       general%gas_saturation%itype == HYDROSTATIC_BC)) then
                    option%io_buffer = 'Gas saturation ' // trim(string)
                    call printErrMsg(option)
                  endif
                case(LIQUID_STATE)
                  if (.not. &
                      (general%liquid_pressure%itype == DIRICHLET_BC .or. &
                       general%liquid_pressure%itype == HYDROSTATIC_BC)) then
                    option%io_buffer = 'Liquid pressure ' // trim(string)
                    call printErrMsg(option)
                  endif
                  if (.not. &
                      (general%mole_fraction%itype == DIRICHLET_BC .or. &
                       general%mole_fraction%itype == HYDROSTATIC_BC)) then
                    option%io_buffer = 'Mole fraction ' // trim(string)
                    call printErrMsg(option)
                  endif
                case(GAS_STATE)
                  if (.not. &
                      (general%gas_pressure%itype == DIRICHLET_BC .or. &
                       general%gas_pressure%itype == HYDROSTATIC_BC)) then
                    option%io_buffer = 'Gas pressure ' // trim(string)
                    call printErrMsg(option)
                  endif
                  if (.not. &
                      (general%mole_fraction%itype == DIRICHLET_BC .or. &
                       general%mole_fraction%itype == HYDROSTATIC_BC)) then
                    option%io_buffer = 'Gas saturation ' // trim(string)
                    call printErrMsg(option)
                  endif
              end select
              if (.not. &
                  (general%temperature%itype == DIRICHLET_BC .or. &
                   general%temperature%itype == HYDROSTATIC_BC)) then
                option%io_buffer = 'Temperature ' // trim(string)
                call printErrMsg(option)
              endif                              
              
              
              do icell=1,initial_condition%region%num_cells
                local_id = initial_condition%region%cell_ids(icell)
                ghosted_id = grid%nL2G(local_id)
                iend = local_id*option%nflowdof
                ibegin = iend-option%nflowdof+1
                if (cur_patch%imat(ghosted_id) <= 0) then
                  xx_p(ibegin:iend) = 0.d0
                  iphase_loc_p(ghosted_id) = 0
                  cycle
                endif
                ! decrement ibegin to give a local offset of 0
                ibegin = ibegin - 1
                select case(initial_condition%flow_condition%iphase)
                  case(TWO_PHASE_STATE)
                    xx_p(ibegin+GENERAL_GAS_PRESSURE_DOF) = &
                      general%gas_pressure%flow_dataset%time_series%cur_value(1)
                    xx_p(ibegin+GENERAL_GAS_SATURATION_DOF) = &
                      general%gas_saturation%flow_dataset%time_series%cur_value(1)
                    temperature = general%temperature%flow_dataset%time_series%cur_value(1)
                    call psat(temperature,p_sat,ierr)
                    ! p_a = p_g - p_s(T)
                    xx_p(ibegin+GENERAL_AIR_PRESSURE_DOF) = &
                      general%gas_pressure%flow_dataset%time_series%cur_value(1) - &
                      p_sat
                  case(LIQUID_STATE)
                    xx_p(ibegin+GENERAL_LIQUID_PRESSURE_DOF) = &
                      general%liquid_pressure%flow_dataset%time_series%cur_value(1)
                    xx_p(ibegin+GENERAL_LIQUID_STATE_MOLE_FRACTION_DOF) = &
                      general%mole_fraction%flow_dataset%time_series%cur_value(1)
                    xx_p(ibegin+GENERAL_LIQUID_STATE_TEMPERATURE_DOF) = &
                      general%temperature%flow_dataset%time_series%cur_value(1)
                  case(GAS_STATE)
                    xx_p(ibegin+GENERAL_GAS_PRESSURE_DOF) = &
                      general%gas_pressure%flow_dataset%time_series%cur_value(1)
                    xx_p(ibegin+GENERAL_AIR_PRESSURE_DOF) = &
                      general%gas_pressure%flow_dataset%time_series%cur_value(1) * &
                      general%mole_fraction%flow_dataset%time_series%cur_value(1)
                    xx_p(ibegin+GENERAL_GAS_STATE_TEMPERATURE_DOF) = &
                      general%temperature%flow_dataset%time_series%cur_value(1)
                end select
                iphase_loc_p(ghosted_id) = initial_condition%flow_condition%iphase
                cur_patch%aux%Global%aux_vars(ghosted_id)%istate = &
                  initial_condition%flow_condition%iphase
              enddo
            else
              do iconn=1,initial_condition%connection_set%num_connections
                local_id = initial_condition%connection_set%id_dn(iconn)
                ghosted_id = grid%nL2G(local_id)
                iend = local_id*option%nflowdof
                ibegin = iend-option%nflowdof+1
                if (cur_patch%imat(ghosted_id) <= 0) then
                  xx_p(ibegin:iend) = 0.d0
                  iphase_loc_p(ghosted_id) = 0
                  cycle
                endif
                xx_p(ibegin:iend) = &
                  initial_condition%flow_aux_real_var(1:option%nflowdof,iconn)
                iphase_loc_p(ghosted_id) = initial_condition%flow_condition%iphase
                cur_patch%aux%Global%aux_vars(ghosted_id)%istate = &
                  initial_condition%flow_condition%iphase
              enddo
            endif
            initial_condition => initial_condition%next
          enddo
     
          call GridVecRestoreArrayF90(grid,field%flow_xx,xx_p, ierr)
          call GridVecRestoreArrayF90(grid,field%iphas_loc,iphase_loc_p,ierr)
              
        case default
          ! assign initial conditions values to domain
          if (discretization%itype == STRUCTURED_GRID_MIMETIC.or. &
              discretization%itype == UNSTRUCTURED_GRID_MIMETIC) then
            call GridVecGetArrayF90(grid,field%flow_xx, xx_p, ierr); CHKERRQ(ierr)
            call VecGetArrayF90(field%flow_xx_faces, xx_faces_p, ierr); CHKERRQ(ierr)        
          else
            call GridVecGetArrayF90(grid,field%flow_xx,xx_p, ierr); CHKERRQ(ierr)
          end if
          call GridVecGetArrayF90(grid,field%iphas_loc,iphase_loc_p,ierr)
      
          xx_p = -999.d0
      
          initial_condition => cur_patch%initial_conditions%first
          do
      
            if (.not.associated(initial_condition)) exit

            if (discretization%itype == STRUCTURED_GRID_MIMETIC.or. &
                discretization%itype == UNSTRUCTURED_GRID_MIMETIC) then
#ifdef DASVYAT
              if (.not.associated(initial_condition%flow_aux_real_var)) then
                do icell=1,initial_condition%region%num_cells
                  local_id = initial_condition%region%cell_ids(icell)
                  ghosted_id = grid%nL2G(local_id)
                  if (cur_patch%imat(ghosted_id) <= 0) then
                    iphase_loc_p(ghosted_id) = 0
                    cycle
                  endif
                  iphase_loc_p(ghosted_id)=initial_condition%flow_condition%iphase
                 enddo
              else
                do iface=1,initial_condition%numfaces_set
                  ghosted_id = initial_condition%faces_set(iface)
                  local_id = grid%fG2L(ghosted_id)
                  if (local_id > 0) then
                    iend = local_id*option%nflowdof
                    ibegin = iend-option%nflowdof+1
                    xx_faces_p(ibegin:iend) = &
                    initial_condition%flow_aux_real_var(1:option%nflowdof,iface)
                  endif
                enddo
                do iconn=1,initial_condition%connection_set%num_connections
                  local_id = initial_condition%region%cell_ids(iconn)
                  ghosted_id = grid%nL2G(local_id)
                  iend = local_id*option%nflowdof
                  ibegin = iend-option%nflowdof+1
                  if (cur_patch%imat(ghosted_id) <= 0) then
                    xx_p(ibegin:iend) = 0.d0
                    iphase_loc_p(ghosted_id) = 0
                    cycle
                  endif
                  xx_p(ibegin:iend) = &
                    initial_condition%flow_aux_real_var(1:option%nflowdof, &
                                                        iconn + &
                                                        initial_condition%numfaces_set)
                  xx_faces_p(grid%nlmax_faces + &
                            ibegin:grid%nlmax_faces + iend) = &
                                xx_p(ibegin:iend) ! for LP -formulation
                  iphase_loc_p(ghosted_id) = &
                    initial_condition%flow_aux_int_var(1,iconn + &
                                                      initial_condition%numfaces_set)
                enddo
              endif
#endif
            else 
              use_dataset = PETSC_FALSE
              dataset_flag = PETSC_FALSE
              do idof = 1, option%nflowdof
                dataset =>  initial_condition%flow_condition% &
                                 sub_condition_ptr(idof)%ptr% &
                                 flow_dataset%dataset
                if (associated(dataset)) then
                  if (dataset%is_cell_indexed) then
                    use_dataset = PETSC_TRUE
                    dataset_flag(idof) = PETSC_TRUE
                    call ConditionControlMapDatasetToVec(realization, &
                           initial_condition%flow_condition% &
                             sub_condition_ptr(idof)%ptr% &
                             flow_dataset%dataset,idof,field%flow_xx,GLOBAL)
                  endif
                endif
              enddo            
              if (.not.associated(initial_condition%flow_aux_real_var) .and. &
                  .not.associated(initial_condition%flow_condition)) then
                option%io_buffer = 'Flow condition is NULL in initial condition'
                call printErrMsg(option)
              endif
              if (associated(initial_condition%flow_aux_real_var)) then
                num_connections = &
                  initial_condition%connection_set%num_connections
                conn_id_ptr => initial_condition%connection_set%id_dn
              else
                num_connections = initial_condition%region%num_cells
                conn_id_ptr => initial_condition%region%cell_ids
              endif
              do iconn=1, num_connections
                local_id = conn_id_ptr(iconn)
                ghosted_id = grid%nL2G(local_id)
                iend = local_id*option%nflowdof
                ibegin = iend-option%nflowdof+1
                if (cur_patch%imat(ghosted_id) <= 0) then
                  xx_p(ibegin:iend) = 0.d0
                  iphase_loc_p(ghosted_id) = 0
                  cycle
                endif
                if (associated(initial_condition%flow_aux_real_var)) then
                  do idof = 1, option%nflowdof
                    if (.not.dataset_flag(idof)) then
                      xx_p(ibegin+idof-1) =  &
                        initial_condition%flow_aux_real_var(idof,iconn)
                    endif
                  enddo
                else
                  do idof = 1, option%nflowdof
                    if (.not.dataset_flag(idof)) then
                      xx_p(ibegin+idof-1) = &
                        initial_condition%flow_condition% &
                          sub_condition_ptr(idof)%ptr%flow_dataset% &
                          time_series%cur_value(1)
                    endif
                  enddo
                endif
                iphase_loc_p(ghosted_id) = &
                  initial_condition%flow_condition%iphase
                if (option%iflowmode == G_MODE) then
                  cur_patch%aux%Global%aux_vars(ghosted_id)%istate = &
                    int(iphase_loc_p(ghosted_id))
                endif
              enddo
            end if
            initial_condition => initial_condition%next
          enddo
     
          call GridVecRestoreArrayF90(grid,field%flow_xx,xx_p, ierr)

      end select 
   
      cur_patch => cur_patch%next
    enddo
    cur_level => cur_level%next
  enddo
  
  if (option%iflowmode == G_MODE) then
    call GlobalAuxVarStrip(global_aux)
    call GeneralAuxVarStrip(general_aux)
  endif  
   
  ! update dependent vectors
  call DiscretizationGlobalToLocal(discretization,field%flow_xx,field%flow_xx_loc,NFLOWDOF)  
  

  call VecCopy(field%flow_xx, field%flow_yy, ierr)
  call DiscretizationLocalToLocal(discretization,field%iphas_loc,field%iphas_loc,ONEDOF)  
  call DiscretizationLocalToLocal(discretization,field%iphas_loc,field%iphas_old_loc,ONEDOF)

#ifdef DASVYAT
  if (discretization%itype == STRUCTURED_GRID_MIMETIC.or. &
      discretization%itype == UNSTRUCTURED_GRID_MIMETIC) then

    call VecRestoreArrayF90(field%flow_xx_faces,xx_faces_p, ierr)
    call RealizationSetUpBC4Faces(realization)

    !call DiscretizationGlobalToLocalFaces(discretization, field%flow_xx_faces, field%flow_xx_loc_faces, NFLOWDOF)
    call DiscretizationGlobalToLocalLP(discretization, field%flow_xx_faces, field%flow_xx_loc_faces, NFLOWDOF)
    call VecCopy(field%flow_xx_faces, field%flow_yy_faces, ierr)
    call MFDInitializeMassMatrices(realization%discretization%grid,&
                                  realization%field, &
                                  realization%discretization%MFD, realization%option)
    patch%aux%Richards%aux_vars_cell_pressures_up_to_date = PETSC_TRUE

  endif
#endif
!  stop 
end subroutine CondControlAssignFlowInitCond

! ************************************************************************** !
!
! CondControlAssignTranInitCond: Assigns transport initial conditions to model
! author: Glenn Hammond
! date: 11/02/07, 10/18/11
!
! ************************************************************************** !
subroutine CondControlAssignTranInitCond(realization)

  use Realization_class
  use Discretization_module
  use Region_module
  use Option_module
  use Field_module
  use Coupler_module
  use Condition_module
  use Constraint_module
  use Grid_module
  use Dataset_Aux_module
  use Level_module
  use Patch_module
  use Reactive_Transport_Aux_module
  use Reaction_Aux_module
  use Global_Aux_module
  use Reaction_module
  use HDF5_module
  
  implicit none

#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"
  
  type(realization_type) :: realization
  
  PetscInt :: icell, iconn, idof, isub_condition, temp_int, iimmobile
  PetscInt :: local_id, ghosted_id, iend, ibegin
  PetscInt :: irxn, isite, imnrl, ikinrxn
  PetscReal, pointer :: xx_p(:), xx_loc_p(:), porosity_loc(:), vec_p(:)
  PetscErrorCode :: ierr
  
  type(option_type), pointer :: option
  type(field_type), pointer :: field  
  type(patch_type), pointer :: patch
  type(grid_type), pointer :: grid
  type(discretization_type), pointer :: discretization
  type(coupler_type), pointer :: initial_condition
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  type(reaction_type), pointer :: reaction
  type(reactive_transport_auxvar_type), pointer :: rt_aux_vars(:)
  type(global_auxvar_type), pointer :: global_aux_vars(:)
  type(tran_constraint_coupler_type), pointer :: constraint_coupler

  PetscInt :: iphase
  PetscInt :: offset
  PetscBool :: re_equilibrate_at_each_cell
  character(len=MAXSTRINGLENGTH) :: string, string2
  type(dataset_type), pointer :: dataset
  PetscInt :: aq_dataset_to_idof(realization%reaction%naqcomp)
  PetscInt :: iaqdataset, num_aq_datasets
  PetscBool :: use_aq_dataset
  PetscReal :: ave_num_iterations
  PetscReal :: tempreal
  PetscLogDouble :: tstart, tend
  
  option => realization%option
  discretization => realization%discretization
  field => realization%field
  patch => realization%patch
  reaction => realization%reaction
  
  iphase = 1
  
  cur_level => realization%level_list%first
  do 
    if (.not.associated(cur_level)) exit
    cur_patch => cur_level%patch_list%first
    do
      if (.not.associated(cur_patch)) exit

      grid => cur_patch%grid
      rt_aux_vars => cur_patch%aux%RT%aux_vars
      global_aux_vars => cur_patch%aux%Global%aux_vars

      ! assign initial conditions values to domain
      call GridVecGetArrayF90(grid,field%tran_xx,xx_p,ierr)
      call GridVecGetArrayF90(grid,field%porosity_loc,porosity_loc,ierr)
      
      xx_p = -999.d0
      
      initial_condition => cur_patch%initial_conditions%first
      do
      
        if (.not.associated(initial_condition)) exit
        
        constraint_coupler => initial_condition%tran_condition%cur_constraint_coupler

        re_equilibrate_at_each_cell = PETSC_FALSE
        use_aq_dataset = PETSC_FALSE
        num_aq_datasets = 0
        aq_dataset_to_idof = 0
        do idof = 1, reaction%naqcomp ! primary aqueous concentrations
          if (constraint_coupler%aqueous_species%external_dataset(idof)) then
            num_aq_datasets = num_aq_datasets + 1
            aq_dataset_to_idof(num_aq_datasets) = idof
            re_equilibrate_at_each_cell = PETSC_TRUE
            use_aq_dataset = PETSC_TRUE
            string = 'constraint ' // trim(constraint_coupler%constraint_name)
            dataset => DatasetGetPointer(realization%datasets, &
                         constraint_coupler%aqueous_species%constraint_aux_string(idof), &
                         string,option)
            call ConditionControlMapDatasetToVec(realization,dataset,idof, &
                                                 field%tran_xx_loc,LOCAL)
          endif
        enddo

        ! read in heterogeneous mineral volume fractions
        if (associated(constraint_coupler%minerals)) then
          do imnrl = 1, reaction%mineral%nkinmnrl
            if (constraint_coupler%minerals%external_dataset(imnrl)) then
              re_equilibrate_at_each_cell = PETSC_TRUE
              string = 'constraint ' // trim(constraint_coupler%constraint_name)
              dataset => DatasetGetPointer(realization%datasets, &
                           constraint_coupler%minerals%constraint_aux_string(imnrl), &
                           string,option)
              string = '' ! group name
              string2 = dataset%h5_dataset_name ! dataset name
              call HDF5ReadCellIndexedRealArray(realization,field%work, &
                                                dataset%filename, &
                                                string,string2, &
                                                dataset%realization_dependent)
              call DiscretizationGlobalToLocal(discretization,field%work,field%work_loc,ONEDOF)
              call GridVecGetArrayF90(grid,field%work_loc,vec_p,ierr)
              do icell=1,initial_condition%region%num_cells
                local_id = initial_condition%region%cell_ids(icell)
                ghosted_id = grid%nL2G(local_id)
                rt_aux_vars(ghosted_id)%mnrl_volfrac0(imnrl) = vec_p(ghosted_id)
                rt_aux_vars(ghosted_id)%mnrl_volfrac(imnrl) = vec_p(ghosted_id)
              enddo
              call GridVecRestoreArrayF90(grid,field%work_loc,vec_p,ierr)
            endif
          enddo
        endif
          
        ! read in heterogeneous immobile
        if (associated(constraint_coupler%immobile_species)) then
          do iimmobile = 1, reaction%immobile%nimmobile
            if (constraint_coupler%immobile_species%external_dataset(iimmobile)) then
              ! no need to requilibrate at each cell
              string = 'constraint ' // trim(constraint_coupler%constraint_name)
              dataset => DatasetGetPointer(realization%datasets, &
                  constraint_coupler%immobile_species%constraint_aux_string(iimmobile), &
                  string,option)
              string = '' ! group name
              string2 = dataset%h5_dataset_name ! dataset name
              call HDF5ReadCellIndexedRealArray(realization,field%work, &
                                                dataset%filename, &
                                                string,string2, &
                                                dataset%realization_dependent)
              call DiscretizationGlobalToLocal(discretization,field%work,field%work_loc,ONEDOF)
              call GridVecGetArrayF90(grid,field%work_loc,vec_p,ierr)
              do icell=1,initial_condition%region%num_cells
                local_id = initial_condition%region%cell_ids(icell)
                ghosted_id = grid%nL2G(local_id)
                rt_aux_vars(ghosted_id)%immobile(iimmobile) = vec_p(ghosted_id)
              enddo
              call GridVecRestoreArrayF90(grid,field%work_loc,vec_p,ierr)
            endif
          enddo
        endif
          
        if (.not.option%use_isothermal) then
          re_equilibrate_at_each_cell = PETSC_TRUE
        endif
        
        if (use_aq_dataset) then
          call GridVecGetArrayF90(grid,field%tran_xx_loc,xx_loc_p,ierr); CHKERRQ(ierr)
          call PetscTime(tstart,ierr) 
        endif
        
        ave_num_iterations = 0.d0
        do icell=1,initial_condition%region%num_cells
          local_id = initial_condition%region%cell_ids(icell)
          ghosted_id = grid%nL2G(local_id)
          iend = local_id*option%ntrandof
          ibegin = iend-option%ntrandof+1
          if (cur_patch%imat(ghosted_id) <= 0) then
            xx_p(ibegin:iend) = 1.d-200
            cycle
          endif
          if (re_equilibrate_at_each_cell) then
            if (use_aq_dataset) then
              offset = (ghosted_id-1)*option%ntrandof
              do iaqdataset = 1, num_aq_datasets
                ! remember that xx_loc_p holds the data set values that were read in
                temp_int = aq_dataset_to_idof(iaqdataset)
                constraint_coupler%aqueous_species%constraint_conc(temp_int) = &
                  xx_loc_p(offset+temp_int)
              enddo
            endif
            option%iflag = grid%nG2A(grid%nL2G(local_id))
            if (icell == 1) then
              call ReactionEquilibrateConstraint(rt_aux_vars(ghosted_id), &
                global_aux_vars(ghosted_id),reaction, &
                constraint_coupler%constraint_name, &
                constraint_coupler%aqueous_species, &
                constraint_coupler%minerals, &
                constraint_coupler%surface_complexes, &
                constraint_coupler%colloids, &
                constraint_coupler%immobile_species, &
                porosity_loc(ghosted_id), &
                constraint_coupler%num_iterations, &
                PETSC_FALSE,option)
            else
!geh              call RTAuxVarCopy(rt_aux_vars(ghosted_id), &
!geh                rt_aux_vars(grid%nL2G(initial_condition%region%cell_ids(icell-1))), &
!geh                option)
              rt_aux_vars(ghosted_id)%pri_molal = &
                rt_aux_vars(grid%nL2G(initial_condition%region%cell_ids(icell-1)))%pri_molal
              call ReactionEquilibrateConstraint(rt_aux_vars(ghosted_id), &
                global_aux_vars(ghosted_id),reaction, &
                constraint_coupler%constraint_name, &
                constraint_coupler%aqueous_species, &
                constraint_coupler%minerals, &
                constraint_coupler%surface_complexes, &
                constraint_coupler%colloids, &
                constraint_coupler%immobile_species, &
                porosity_loc(ghosted_id), &
                constraint_coupler%num_iterations, &
                PETSC_TRUE,option)
            endif
            option%iflag = 0
            ave_num_iterations = ave_num_iterations + &
              constraint_coupler%num_iterations
          endif
          ! ibegin is the local non-ghosted offset: (local_id-1)*option%ntrandof+1
          offset = ibegin + reaction%offset_aqueous - 1
          ! primary aqueous concentrations
          do idof = 1, reaction%naqcomp 
            xx_p(offset+idof) = &
              constraint_coupler%aqueous_species%basis_molarity(idof) / &
              global_aux_vars(ghosted_id)%den_kg(iphase)*1000.d0 ! convert molarity -> molality
          enddo
          ! mineral volume fractions
          if (associated(constraint_coupler%minerals)) then
            do imnrl = 1, reaction%mineral%nkinmnrl
              ! if read from a dataset, the vol frac was set above.  Don't want to
              ! overwrite
              if (.not.constraint_coupler%minerals%external_dataset(imnrl)) then
                rt_aux_vars(ghosted_id)%mnrl_volfrac0(imnrl) = &
                  constraint_coupler%minerals%constraint_vol_frac(imnrl)
                rt_aux_vars(ghosted_id)%mnrl_volfrac(imnrl) = &
                  constraint_coupler%minerals%constraint_vol_frac(imnrl)
              endif
              rt_aux_vars(ghosted_id)%mnrl_area0(imnrl) = &
                constraint_coupler%minerals%constraint_area(imnrl)
              rt_aux_vars(ghosted_id)%mnrl_area(imnrl) = &
                constraint_coupler%minerals%constraint_area(imnrl)
            enddo
          endif
          ! kinetic surface complexes
          if (associated(constraint_coupler%surface_complexes)) then
            do idof = 1, reaction%surface_complexation%nkinsrfcplx
              rt_aux_vars(ghosted_id)%kinsrfcplx_conc(idof,-1) = & !geh: to catch bug
                constraint_coupler%surface_complexes%constraint_conc(idof)
            enddo
            do ikinrxn = 1, reaction%surface_complexation%nkinsrfcplxrxn
              irxn = reaction%surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
              isite = reaction%surface_complexation%srfcplxrxn_to_surf(irxn)
              rt_aux_vars(ghosted_id)%kinsrfcplx_free_site_conc(isite) = &
                constraint_coupler%surface_complexes%basis_free_site_conc(isite)
            enddo
          endif
          ! this is for the multi-rate surface complexation model
          if (reaction%surface_complexation%nkinmrsrfcplxrxn > 0 .and. &
            ! geh: if we re-equilibrate at each grid cell, we do not want to
            ! overwrite the reequilibrated values with those from the constraint
              .not. re_equilibrate_at_each_cell) then
            ! copy over total sorbed concentration
            rt_aux_vars(ghosted_id)%kinmr_total_sorb = &
              constraint_coupler%rt_auxvar%kinmr_total_sorb
            ! copy over free site concentration
            rt_aux_vars(ghosted_id)%srfcplxrxn_free_site_conc = &
              constraint_coupler%rt_auxvar%srfcplxrxn_free_site_conc
          endif
          ! colloids fractions
          if (associated(constraint_coupler%colloids)) then
            offset = ibegin + reaction%offset_colloid - 1
            do idof = 1, reaction%ncoll ! primary aqueous concentrations
              xx_p(offset+idof) = &
                constraint_coupler%colloids%basis_conc_mob(idof) / &
                global_aux_vars(ghosted_id)%den_kg(iphase)*1000.d0 ! convert molarity -> molality
              rt_aux_vars(ghosted_id)%colloid%conc_imb(idof) = &
                constraint_coupler%colloids%basis_conc_imb(idof)
            enddo
          endif
          ! immobile
          if (associated(constraint_coupler%immobile_species)) then
            offset = ibegin + reaction%offset_immobile - 1
            do iimmobile = 1, reaction%immobile%nimmobile
              if (constraint_coupler%immobile_species%external_dataset(iimmobile)) then
                ! already read into rt_aux_vars above.
                xx_p(offset+iimmobile) = &
                  rt_aux_vars(ghosted_id)%immobile(iimmobile)
              else
                xx_p(offset+iimmobile) = &
                  constraint_coupler%immobile_species%constraint_conc(iimmobile)
                rt_aux_vars(ghosted_id)%immobile(iimmobile) = &
                  constraint_coupler%immobile_species%constraint_conc(iimmobile)
              endif
            enddo
          endif
        enddo ! icell=1,initial_condition%region%num_cells
        if (use_aq_dataset) then
          call PetscTime(tend,ierr) 
          call GridVecRestoreArrayF90(grid,field%tran_xx_loc,xx_loc_p,ierr)
          ave_num_iterations = ave_num_iterations / &
            initial_condition%region%num_cells
          write(option%io_buffer,&
                '("Average number of iterations in ReactionEquilibrateConstraint():", &
                & f5.1)') ave_num_iterations
          call printMsg(option)
          write(option%io_buffer,'(f10.2," Seconds to equilibrate constraints")') &
            tend-tstart
          call printMsg(option)
        endif
        initial_condition => initial_condition%next
      enddo
      
      call GridVecRestoreArrayF90(grid,field%tran_xx,xx_p, ierr)
      call GridVecRestoreArrayF90(grid,field%porosity_loc,porosity_loc,ierr)

      cur_patch => cur_patch%next
    enddo
    cur_level => cur_level%next
  enddo
  
  ! check to ensure that minimum concentration is not less than or equal
  ! to zero
  call VecMin(field%tran_xx,PETSC_NULL_INTEGER,tempreal,ierr)
  if (tempreal <= 0.d0) then
    option%io_buffer = 'ERROR: Zero concentrations found in initial ' // &
      'transport solution.'
    call printMsg(option)
    ! now figure out which species have zero concentrations
    do idof = 1, option%ntrandof
      call VecStrideMin(field%tran_xx,idof-1,offset,tempreal,ierr)
      if (tempreal <= 0.d0) then
        write(string,*) tempreal
        if (idof <= reaction%naqcomp) then
          string2 = '  Aqueous species "' // &
            trim(reaction%primary_species_names(idof))
        else
          string2 = '  Immobile species "' // &
            trim(reaction%immobile%names(idof-reaction%offset_immobile))
        endif
          string2 = trim(string2) // &
            '" has zero concentration (' // &
            trim(adjustl(string)) // ').'
        call printMsg(option)
      endif
    enddo
    option%io_buffer = 'Free ion concentations must be positive.  Try ' // &
      'using a small value such as 1.e-20 or 1.e-40 instead of zero.'
    call printErrMsg(option)
  endif
  
  ! update dependent vectors
  call DiscretizationGlobalToLocal(discretization,field%tran_xx, &
                                   field%tran_xx_loc,NTRANDOF)  
  call VecCopy(field%tran_xx, field%tran_yy, ierr)

end subroutine CondControlAssignTranInitCond

! ************************************************************************** !
!
! ConditionControlMapDatasetToVec: maps an external dataset to a PETSc vec
!                                  representing values at each grid cell
! author: Glenn Hammond
! date: 03/23/12
!
! ************************************************************************** !
subroutine ConditionControlMapDatasetToVec(realization,dataset,idof, &
                                           mdof_vec,vec_type)
  use Realization_class
  use Option_module
  use Field_module
  use Dataset_Aux_module
  use HDF5_module
  use Discretization_module

  implicit none
  
#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"  
  
  type(realization_type) :: realization
  type(dataset_type), pointer :: dataset
  PetscInt :: idof
  Vec :: mdof_vec
  PetscInt :: vec_type

  type(field_type), pointer :: field
  type(option_type), pointer :: option
  character(len=MAXSTRINGLENGTH) :: string, string2
  PetscErrorCode :: ierr

  field => realization%field
  option => realization%option
  
  if (associated(dataset)) then
    string = '' ! group name
    ! have to copy to string2 due to mismatch in string size
    string2 = dataset%h5_dataset_name
    call HDF5ReadCellIndexedRealArray(realization,field%work, &
                                      dataset%filename, &
                                      string,string2, &
                                      dataset%realization_dependent)
    if (vec_type == GLOBAL) then
      call VecStrideScatter(field%work,idof-1,mdof_vec, &
                            INSERT_VALUES,ierr)    
    else
      call DiscretizationGlobalToLocal(realization%discretization, &
                                       field%work, &
                                       field%work_loc,ONEDOF)
      call VecStrideScatter(field%work_loc,idof-1,mdof_vec, &
                            INSERT_VALUES,ierr)    
    endif
  endif

end subroutine ConditionControlMapDatasetToVec

! ************************************************************************** !
!
! CondControlScaleSourceSink: Scales select source/sinks based on perms
! author: Glenn Hammond
! date: 09/03/08, 10/18/11
!
! ************************************************************************** !
subroutine CondControlScaleSourceSink(realization)

  use Realization_class
  use Discretization_module
  use Region_module
  use Option_module
  use Field_module
  use Coupler_module
  use Connection_module
  use Condition_module
  use Grid_module
  use Level_module
  use Patch_module
  
  implicit none

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

  
  type(realization_type) :: realization
  
  PetscErrorCode :: ierr
  
  type(option_type), pointer :: option
  type(field_type), pointer :: field  
  type(patch_type), pointer :: patch
  type(grid_type), pointer :: grid
  type(discretization_type), pointer :: discretization
  type(coupler_type), pointer :: cur_source_sink
  type(connection_set_type), pointer :: cur_connection_set
  
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  PetscReal, pointer :: vec_ptr(:)
  PetscReal, pointer :: perm_loc_ptr(:)
  PetscInt :: local_id
  PetscInt :: ghosted_id, neighbor_ghosted_id
  PetscInt :: iconn
  PetscReal :: scale, sum
  PetscInt :: icount
  PetscInt :: x_count, y_count, z_count
  PetscInt, parameter :: x_width = 1, y_width = 1, z_width = 0
  
  PetscInt :: ghosted_neighbors(0:27)
  
  option => realization%option
  discretization => realization%discretization
  field => realization%field
  patch => realization%patch

  ! GB: grid was uninitialized
  grid => patch%grid
  call GridVecGetArrayF90(grid,field%perm_xx_loc,perm_loc_ptr,ierr)

  cur_level => realization%level_list%first
  do 
    if (.not.associated(cur_level)) exit
    cur_patch => cur_level%patch_list%first
    do
      if (.not.associated(cur_patch)) exit
      ! BIG-TIME warning here.  I assume that all source/sink cells are within 
      ! a single patch - geh

      grid => cur_patch%grid

      cur_source_sink => cur_patch%source_sinks%first
      do
        if (.not.associated(cur_source_sink)) exit

        call VecZeroEntries(field%work,ierr)
        call GridVecGetArrayF90(grid,field%work,vec_ptr,ierr)

        cur_connection_set => cur_source_sink%connection_set
    
        do iconn = 1, cur_connection_set%num_connections
          local_id = cur_connection_set%id_dn(iconn)
          ghosted_id = grid%nL2G(local_id)

          select case(option%iflowmode)
            case(RICHARDS_MODE,G_MODE)
               call GridGetGhostedNeighbors(grid,ghosted_id,DMDA_STENCIL_STAR, &
                                            x_width,y_width,z_width, &
                                            x_count,y_count,z_count, &
                                            ghosted_neighbors,option)
               ! ghosted neighbors is ordered first in x, then, y, then z
               icount = 0
               sum = 0.d0
               ! x-direction
               do while (icount < x_count)
                 icount = icount + 1
                 neighbor_ghosted_id = ghosted_neighbors(icount)
                 sum = sum + perm_loc_ptr(neighbor_ghosted_id)* &
                             grid%structured_grid%dy(neighbor_ghosted_id)* &
                             grid%structured_grid%dz(neighbor_ghosted_id)
                 
               enddo
               ! y-direction
               do while (icount < x_count + y_count)
                 icount = icount + 1
                 neighbor_ghosted_id = ghosted_neighbors(icount)                 
                 sum = sum + perm_loc_ptr(neighbor_ghosted_id)* &
                             grid%structured_grid%dx(neighbor_ghosted_id)* &
                             grid%structured_grid%dz(neighbor_ghosted_id)
                 
               enddo
               ! z-direction
               do while (icount < x_count + y_count + z_count)
                 icount = icount + 1
                 neighbor_ghosted_id = ghosted_neighbors(icount)                 
                 sum = sum + perm_loc_ptr(neighbor_ghosted_id)* &
                             grid%structured_grid%dx(neighbor_ghosted_id)* &
                             grid%structured_grid%dy(neighbor_ghosted_id)
               enddo
               vec_ptr(local_id) = vec_ptr(local_id) + sum
            case(TH_MODE)
            case(THC_MODE)
            case(THMC_MODE)
            case(MPH_MODE)
            case(IMS_MODE)
            case(MIS_MODE)
            case(FLASH2_MODE)
          end select

        enddo
        
        call GridVecRestoreArrayF90(grid,field%work,vec_ptr,ierr)
        call VecNorm(field%work,NORM_1,scale,ierr)
        scale = 1.d0/scale
        call VecScale(field%work,scale,ierr)

        call GridVecGetArrayF90(grid,field%work,vec_ptr, ierr)
        do iconn = 1, cur_connection_set%num_connections      
          local_id = cur_connection_set%id_dn(iconn)
          select case(option%iflowmode)
            case(RICHARDS_MODE,G_MODE)
              cur_source_sink%flow_aux_real_var(ONE_INTEGER,iconn) = &
                vec_ptr(local_id)
            case(TH_MODE)
            case(THC_MODE)
            case(THMC_MODE)
            case(MPH_MODE)
            case(IMS_MODE)
            case(MIS_MODE)
            case(FLASH2_MODE)
          end select

        enddo
        call GridVecRestoreArrayF90(grid,field%work,vec_ptr,ierr)
        
        cur_source_sink => cur_source_sink%next
      enddo
      cur_patch => cur_patch%next
    enddo
    cur_level => cur_level%next
  enddo

  call GridVecRestoreArrayF90(grid,field%perm_xx_loc,perm_loc_ptr, ierr)
   
end subroutine CondControlScaleSourceSink

! ************************************************************************** !
#ifdef SURFACE_FLOW
subroutine CondControlAssignFlowInitCondSurface(surf_realization)

  use Surface_Realization_class
  use Discretization_module
  use Region_module
  use Option_module
  use Surface_Field_module
  use Coupler_module
  use Condition_module
  use Grid_module
  use Level_module
  use Patch_module
  use Water_EOS_module
  
  implicit none

#include "finclude/petscvec.h"
#include "finclude/petscvec.h90"
  
  type(surface_realization_type) :: surf_realization
  
  PetscInt :: icell, iconn, idof, iface
  PetscInt :: local_id, ghosted_id, iend, ibegin
  PetscReal, pointer :: xx_p(:)!, iphase_loc_p(:), xx_faces_p(:)
  PetscErrorCode :: ierr
  
  PetscReal :: temperature, p_sat
  character(len=MAXSTRINGLENGTH) :: string
  
  type(option_type), pointer :: option
  type(surface_field_type), pointer :: surf_field  
  type(patch_type), pointer :: patch
  type(grid_type), pointer :: grid
  type(discretization_type), pointer :: discretization
  type(coupler_type), pointer :: initial_condition
  type(level_type), pointer :: cur_level
  type(patch_type), pointer :: cur_patch
  type(flow_general_condition_type), pointer :: general

  option => surf_realization%option
  discretization => surf_realization%discretization
  surf_field => surf_realization%surf_field
  patch => surf_realization%patch


  cur_level => surf_realization%level_list%first
  do 
    if (.not.associated(cur_level)) exit
    cur_patch => cur_level%patch_list%first
    do
      if (.not.associated(cur_patch)) exit

      grid => cur_patch%grid

      select case(option%iflowmode)
      
        case(G_MODE) ! general phase mode

        case (RICHARDS_MODE,TH_MODE)
          ! assign initial conditions values to domain
          call GridVecGetArrayF90(grid,surf_field%flow_xx,xx_p, ierr); CHKERRQ(ierr)
    
          xx_p = -999.d0
      
          initial_condition => cur_patch%initial_conditions%first
          do
      
            if (.not.associated(initial_condition)) exit

              if (.not.associated(initial_condition%flow_aux_real_var)) then
                if (.not.associated(initial_condition%flow_condition)) then
                  option%io_buffer = 'Flow condition is NULL in initial condition'
                  call printErrMsg(option)
                endif
                do icell=1,initial_condition%region%num_cells
                  local_id = initial_condition%region%cell_ids(icell)
                  ghosted_id = grid%nL2G(local_id)
                  iend = local_id*option%nflowdof
                  ibegin = iend-option%nflowdof+1
                  if (cur_patch%imat(ghosted_id) <= 0) then
                    xx_p(ibegin:iend) = 0.d0
                    cycle
                  endif
                  do idof = 1, option%nflowdof
                    xx_p(ibegin+idof-1) = &
                          initial_condition%flow_condition% &
                          sub_condition_ptr(idof)%ptr%flow_dataset%time_series%cur_value(1)
                    !TODO(GB): Correct the initialization of surface flow condition
                    if (idof == 1) xx_p(ibegin+idof-1) = 0.d0
                    !if (idof == 1.and.option%iflowmode==TH_MODE) then
                    !  xx_p(ibegin+idof-1) = 0.d0+option%reference_pressure
                    !endif
                  enddo
                enddo
              else
                do iconn=1,initial_condition%connection_set%num_connections
                  local_id = initial_condition%connection_set%id_dn(iconn)
                  ghosted_id = grid%nL2G(local_id)
                  iend = local_id*option%nflowdof
                  ibegin = iend-option%nflowdof+1
                  if (cur_patch%imat(ghosted_id) <= 0) then
                    xx_p(ibegin:iend) = 0.d0
                   cycle
                  endif
                  xx_p(ibegin:iend) = &
                        initial_condition%flow_aux_real_var(1:option%nflowdof,iconn)
                  !TODO(gb): Correct the initialization of surface flow condition
                  xx_p(ibegin) = 0.d0
                  !if (idof == 1.and.option%iflowmode==TH_MODE) then
                  !  xx_p(ibegin) = 0.d0 + option%reference_pressure
                  !endif
                enddo
              endif
            initial_condition => initial_condition%next
          enddo
     
          call GridVecRestoreArrayF90(grid,surf_field%flow_xx,xx_p, ierr)
        case default
          option%io_buffer = 'CondControlAssignFlowInitCondSurface not ' // &
            'for this mode'
          call printErrMsg(option)
      end select 
   
      cur_patch => cur_patch%next
    enddo
    cur_level => cur_level%next
  enddo
   
  ! update dependent vectors
  call DiscretizationGlobalToLocal(discretization,surf_field%flow_xx,surf_field%flow_xx_loc,NFLOWDOF)

end subroutine CondControlAssignFlowInitCondSurface
#endif
! SURFACE_FLOW

end module Condition_Control_module