condition_control.F90       coverage:  66.67 %func     50.88 %block


     1) module Condition_Control_module
     2)   
     3)   ! This module store routines that operate on conditions from a level above 
     4)   ! that of the realization_module.  This is necessary to access capability
     5)   ! such as HDF5 which is unavailable from within the realization object
     6)   ! and below.  Routines in this module will loop over realization, levels,
     7)   ! and patches without calling underlying level/patch versions of the
     8)   ! subroutines, which is common in realization.F90 - GEH
     9)  
    10)   use PFLOTRAN_Constants_module
    11) 
    12)   implicit none
    13) 
    14)   private
    15) 
    16) #include "petsc/finclude/petscsys.h"
    17) #include "petsc/finclude/petscvec.h"
    18) #include "petsc/finclude/petscvec.h90"
    19) 
    20)   public :: CondControlAssignFlowInitCond, &
    21)             CondControlAssignTranInitCond, &
    22)             CondControlAssignFlowInitCondSurface, &
    23)             CondControlScaleSourceSink
    24)  
    25) contains
    26) 
    27) ! ************************************************************************** !
    28) 
    29) subroutine CondControlAssignFlowInitCond(realization)
    30)   ! 
    31)   ! Assigns flow initial conditions to model
    32)   ! 
    33)   ! Author: Glenn Hammond
    34)   ! Date: 11/02/07, 10/18/11
    35)   ! 
    36) 
    37)   use Realization_Subsurface_class
    38)   use Discretization_module
    39)   use Region_module
    40)   use Option_module
    41)   use Field_module
    42)   use Coupler_module
    43)   use Condition_module
    44)   use Dataset_Base_class
    45)   use Dataset_Gridded_HDF5_class
    46)   use Dataset_Common_HDF5_class
    47)   use Grid_module
    48)   use Patch_module
    49)   use EOS_Water_module
    50) 
    51)   use Global_module
    52)   use Global_Aux_module
    53)   use General_Aux_module
    54)   use TOilIms_Aux_module
    55)   
    56)   implicit none
    57) 
    58) #include "petsc/finclude/petscvec.h"
    59) #include "petsc/finclude/petscvec.h90"
    60)   
    61)   class(realization_subsurface_type) :: realization
    62)   
    63)   PetscInt :: icell, iconn, idof, iface
    64)   PetscInt :: local_id, ghosted_id, iend, ibegin
    65)   PetscReal, pointer :: xx_p(:), iphase_loc_p(:)
    66)   PetscErrorCode :: ierr
    67)   
    68)   character(len=MAXSTRINGLENGTH) :: string
    69)   
    70)   type(option_type), pointer :: option
    71)   type(field_type), pointer :: field  
    72)   type(patch_type), pointer :: patch
    73)   type(grid_type), pointer :: grid
    74)   type(discretization_type), pointer :: discretization
    75)   type(coupler_type), pointer :: initial_condition
    76)   type(patch_type), pointer :: cur_patch
    77)   type(flow_general_condition_type), pointer :: general
    78)   type(flow_toil_ims_condition_type), pointer :: toil_ims
    79)   class(dataset_base_type), pointer :: dataset
    80)   type(global_auxvar_type) :: global_aux
    81)   type(general_auxvar_type) :: general_aux
    82)   PetscBool :: use_dataset
    83)   PetscBool :: dataset_flag(realization%option%nflowdof)
    84)   PetscInt :: num_connections
    85)   PetscInt, pointer :: conn_id_ptr(:)
    86)   PetscInt :: offset, istate
    87)   PetscReal :: x(realization%option%nflowdof)
    88)   PetscReal :: temperature, p_sat
    89)   PetscReal :: tempreal
    90) 
    91)   option => realization%option
    92)   discretization => realization%discretization
    93)   field => realization%field
    94)   patch => realization%patch
    95) 
    96)   ! to catch uninitialized grid cells.  see VecMin check at bottom.
    97)   call VecGetArrayF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
    98)   iphase_loc_p = UNINITIALIZED_DOUBLE
    99)   call VecRestoreArrayF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
   100) 
   101)   if (option%iflowmode == G_MODE) then
   102)     call GlobalAuxVarInit(global_aux,option)
   103)     call GeneralAuxVarInit(general_aux,PETSC_FALSE,option)
   104)   endif
   105)   
   106)   cur_patch => realization%patch_list%first
   107)   do
   108)     if (.not.associated(cur_patch)) exit
   109) 
   110)     grid => cur_patch%grid
   111) 
   112)     select case(option%iflowmode)
   113)       
   114)       case(G_MODE) ! general phase mode
   115) 
   116)         call VecGetArrayF90(field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
   117)         call VecGetArrayF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
   118)       
   119)         xx_p = UNINITIALIZED_DOUBLE
   120)       
   121)         initial_condition => cur_patch%initial_condition_list%first
   122)         do
   123)       
   124)           if (.not.associated(initial_condition)) exit
   125) 
   126)           if (.not.associated(initial_condition%flow_aux_real_var)) then
   127)             if (.not.associated(initial_condition%flow_condition)) then
   128)               option%io_buffer = 'Flow condition is NULL in initial condition'
   129)               call printErrMsg(option)
   130)             endif
   131)               
   132)             general => initial_condition%flow_condition%general
   133)               
   134)             string = 'in flow condition "' // &
   135)               trim(initial_condition%flow_condition%name) // &
   136)               '" within initial condition "' // &
   137)               trim(initial_condition%flow_condition%name) // &
   138)               '" must be of type Dirichlet or Hydrostatic'
   139)             ! error checking.  the data must match the state
   140)             select case(initial_condition%flow_condition%iphase)
   141)               case(TWO_PHASE_STATE)  
   142)                 if (.not. &
   143)                     (general%gas_pressure%itype == DIRICHLET_BC .or. &
   144)                       general%gas_pressure%itype == HYDROSTATIC_BC)) then
   145)                   option%io_buffer = 'Gas pressure ' // trim(string)
   146)                   call printErrMsg(option)
   147)                 endif
   148)                 if (.not. &
   149)                     (general%gas_saturation%itype == DIRICHLET_BC .or. &
   150)                       general%gas_saturation%itype == HYDROSTATIC_BC)) then
   151)                   option%io_buffer = 'Gas saturation ' // trim(string)
   152)                   call printErrMsg(option)
   153)                 endif
   154)               case(LIQUID_STATE)
   155)                 if (.not. &
   156)                     (general%liquid_pressure%itype == DIRICHLET_BC .or. &
   157)                       general%liquid_pressure%itype == HYDROSTATIC_BC)) then
   158)                   option%io_buffer = 'Liquid pressure ' // trim(string)
   159)                   call printErrMsg(option)
   160)                 endif
   161)                 if (.not. &
   162)                     (general%mole_fraction%itype == DIRICHLET_BC .or. &
   163)                       general%mole_fraction%itype == HYDROSTATIC_BC)) then
   164)                   option%io_buffer = 'Mole fraction ' // trim(string)
   165)                   call printErrMsg(option)
   166)                 endif
   167)               case(GAS_STATE)
   168)                 if (.not. &
   169)                     (general%gas_pressure%itype == DIRICHLET_BC .or. &
   170)                       general%gas_pressure%itype == HYDROSTATIC_BC)) then
   171)                   option%io_buffer = 'Gas pressure ' // trim(string)
   172)                   call printErrMsg(option)
   173)                 endif
   174)                 if (.not. &
   175)                     (general%mole_fraction%itype == DIRICHLET_BC .or. &
   176)                       general%mole_fraction%itype == HYDROSTATIC_BC)) then
   177)                   option%io_buffer = 'Gas saturation ' // trim(string)
   178)                   call printErrMsg(option)
   179)                 endif
   180)             end select
   181)             if (.not. &
   182)                 (general%temperature%itype == DIRICHLET_BC .or. &
   183)                   general%temperature%itype == HYDROSTATIC_BC)) then
   184)               option%io_buffer = 'Temperature ' // trim(string)
   185)               call printErrMsg(option)
   186)             endif                              
   187)               
   188)               
   189)             do icell=1,initial_condition%region%num_cells
   190)               local_id = initial_condition%region%cell_ids(icell)
   191)               ghosted_id = grid%nL2G(local_id)
   192)               iend = local_id*option%nflowdof
   193)               ibegin = iend-option%nflowdof+1
   194)               if (cur_patch%imat(ghosted_id) <= 0) then
   195)                 xx_p(ibegin:iend) = 0.d0
   196)                 iphase_loc_p(ghosted_id) = 0
   197)                 cycle
   198)               endif
   199)               ! decrement ibegin to give a local offset of 0
   200)               ibegin = ibegin - 1
   201)               select case(initial_condition%flow_condition%iphase)
   202)                 case(TWO_PHASE_STATE)
   203)                   xx_p(ibegin+GENERAL_GAS_PRESSURE_DOF) = &
   204)                     general%gas_pressure%dataset%rarray(1)
   205)                   xx_p(ibegin+GENERAL_GAS_SATURATION_DOF) = &
   206)                     general%gas_saturation%dataset%rarray(1)
   207)                   temperature = general%temperature%dataset%rarray(1)
   208)                   if (general_2ph_energy_dof == GENERAL_TEMPERATURE_INDEX) then
   209)                     xx_p(ibegin+GENERAL_ENERGY_DOF) = temperature
   210)                   else
   211)                     call EOSWaterSaturationPressure(temperature,p_sat,ierr)
   212)                     ! p_a = p_g - p_s(T)
   213)                     xx_p(ibegin+GENERAL_2PH_STATE_AIR_PRESSURE_DOF) = &
   214)                       general%gas_pressure%dataset%rarray(1) - &
   215)                       p_sat
   216)                   endif
   217)                 case(LIQUID_STATE)
   218)                   xx_p(ibegin+GENERAL_LIQUID_PRESSURE_DOF) = &
   219)                     general%liquid_pressure%dataset%rarray(1)
   220)                   xx_p(ibegin+GENERAL_LIQUID_STATE_X_MOLE_DOF) = &
   221)                     general%mole_fraction%dataset%rarray(1)
   222)                   xx_p(ibegin+GENERAL_ENERGY_DOF) = &
   223)                     general%temperature%dataset%rarray(1)
   224)                 case(GAS_STATE)
   225)                   xx_p(ibegin+GENERAL_GAS_PRESSURE_DOF) = &
   226)                     general%gas_pressure%dataset%rarray(1)
   227)                   xx_p(ibegin+GENERAL_GAS_STATE_AIR_PRESSURE_DOF) = &
   228)                     general%gas_pressure%dataset%rarray(1) * &
   229)                     general%mole_fraction%dataset%rarray(1)
   230)                   xx_p(ibegin+GENERAL_ENERGY_DOF) = &
   231)                     general%temperature%dataset%rarray(1)
   232)               end select
   233)               iphase_loc_p(ghosted_id) = initial_condition%flow_condition%iphase
   234)               cur_patch%aux%Global%auxvars(ghosted_id)%istate = &
   235)                 initial_condition%flow_condition%iphase
   236)             enddo
   237)           else
   238)             do iconn=1,initial_condition%connection_set%num_connections
   239)               local_id = initial_condition%connection_set%id_dn(iconn)
   240)               ghosted_id = grid%nL2G(local_id)
   241)               if (cur_patch%imat(ghosted_id) <= 0) then
   242)                 iend = local_id*option%nflowdof
   243)                 ibegin = iend-option%nflowdof+1
   244)                 xx_p(ibegin:iend) = 0.d0
   245)                 iphase_loc_p(ghosted_id) = 0
   246)                 cycle
   247)               endif
   248)               offset = (local_id-1)*option%nflowdof
   249)               istate = initial_condition%flow_aux_int_var(1,iconn)
   250)               do idof = 1, option%nflowdof
   251)                 xx_p(offset+idof) = &
   252)                   initial_condition%flow_aux_real_var( &
   253)                     initial_condition%flow_aux_mapping( &
   254)                       dof_to_primary_variable(idof,istate)),iconn)
   255)               enddo
   256)               iphase_loc_p(ghosted_id) = istate
   257)               cur_patch%aux%Global%auxvars(ghosted_id)%istate = istate
   258)             enddo
   259)           endif
   260)           initial_condition => initial_condition%next
   261)         enddo
   262)      
   263)         call VecRestoreArrayF90(field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
   264)         call VecRestoreArrayF90(field%iphas_loc,iphase_loc_p, &
   265)                                 ierr);CHKERRQ(ierr)
   266) 
   267)       case(TOIL_IMS_MODE)
   268) 
   269)         call VecGetArrayF90(field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
   270)         call VecGetArrayF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
   271)       
   272)         xx_p = UNINITIALIZED_DOUBLE
   273)       
   274)         initial_condition => cur_patch%initial_condition_list%first
   275) 
   276)         do
   277)       
   278)           if (.not.associated(initial_condition)) exit
   279) 
   280)           if (.not.associated(initial_condition%flow_aux_real_var)) then
   281)             if (.not.associated(initial_condition%flow_condition)) then
   282)               option%io_buffer = 'Flow condition is NULL in initial condition'
   283)               call printErrMsg(option)
   284)             endif
   285)               
   286)             toil_ims => initial_condition%flow_condition%toil_ims
   287)               
   288)             string = 'in flow condition "' // &
   289)               trim(initial_condition%flow_condition%name) // &
   290)               '" within initial condition "' // &
   291)               trim(initial_condition%flow_condition%name) // &
   292)               '" must be of type Dirichlet or Hydrostatic'
   293)             ! check pressure condition  
   294)             if (.not. &
   295)                (toil_ims%pressure%itype == DIRICHLET_BC .or. &
   296)                  toil_ims%pressure%itype == HYDROSTATIC_BC)) then
   297)                  option%io_buffer = 'Oil pressure ' // trim(string)
   298)                  call printErrMsg(option)
   299)             endif
   300)             ! check saturation condition
   301)             if (.not. &
   302)                (toil_ims%saturation%itype == DIRICHLET_BC .or. &
   303)                  toil_ims%saturation%itype == HYDROSTATIC_BC)) then
   304)                  option%io_buffer = 'Oil saturation ' // trim(string)
   305)                 call printErrMsg(option)
   306)             endif
   307)             ! check temperature condition 
   308)             if (.not. &
   309)                 (toil_ims%temperature%itype == DIRICHLET_BC .or. &
   310)                   toil_ims%temperature%itype == HYDROSTATIC_BC)) then
   311)               option%io_buffer = 'Temperature ' // trim(string)
   312)               call printErrMsg(option)
   313)             endif                              
   314)             ! error checking.
   315)             do icell=1,initial_condition%region%num_cells
   316)               local_id = initial_condition%region%cell_ids(icell)
   317)               ghosted_id = grid%nL2G(local_id)
   318)               iend = local_id*option%nflowdof
   319)               ibegin = iend-option%nflowdof+1
   320)               if (cur_patch%imat(ghosted_id) <= 0) then
   321)                 xx_p(ibegin:iend) = 0.d0
   322)                 iphase_loc_p(ghosted_id) = 0
   323)                 cycle
   324)               endif
   325)               ! decrement ibegin to give a local offset of 0
   326)               ibegin = ibegin - 1
   327)               ! assign initial conditions
   328)               xx_p(ibegin+TOIL_IMS_PRESSURE_DOF) = &
   329)                  toil_ims%pressure%dataset%rarray(1)
   330)               xx_p(ibegin+TOIL_IMS_SATURATION_DOF) = &
   331)                  toil_ims%saturation%dataset%rarray(1)
   332)               xx_p(ibegin+TOIL_IMS_ENERGY_DOF) = & 
   333)                  toil_ims%temperature%dataset%rarray(1)
   334)               ! iphase not required  
   335)               ! assign 0 to pass internal routine test at the end
   336)               iphase_loc_p(ghosted_id) = 0
   337)               !iphase_loc_p(ghosted_id) = initial_condition%flow_condition%iphase
   338)               !cur_patch%aux%Global%auxvars(ghosted_id)%istate = &
   339)               !  initial_condition%flow_condition%iphase
   340)             enddo
   341)           else ! if initial condition values defined in flow_aux_real_var
   342)             do iconn=1,initial_condition%connection_set%num_connections
   343)               local_id = initial_condition%connection_set%id_dn(iconn)
   344)               ghosted_id = grid%nL2G(local_id)
   345)               if (cur_patch%imat(ghosted_id) <= 0) then
   346)                 iend = local_id*option%nflowdof
   347)                 ibegin = iend-option%nflowdof+1
   348)                 xx_p(ibegin:iend) = 0.d0
   349)                 ! iphase not required 
   350)                 ! assign 0 to pass internal routine test at the end
   351)                 iphase_loc_p(ghosted_id) = 0
   352)                 cycle
   353)               endif
   354)               offset = (local_id-1)*option%nflowdof
   355)               !istate = initial_condition%flow_aux_int_var(1,iconn)
   356)               do idof = 1, option%nflowdof
   357)                 xx_p(offset+idof) = &
   358)                   initial_condition%flow_aux_real_var( &
   359)                     initial_condition%flow_aux_mapping( &
   360)                       toil_ims_dof_to_primary_vars(idof)),iconn)
   361)                       !dof_to_primary_variable(idof,istate)),iconn)
   362)                       !toil_ims_dof_to_primary_vars(3)
   363)               enddo
   364)               ! iphase not required 
   365)               ! assign 0 to pass internal routine test at the end
   366)               iphase_loc_p(ghosted_id) = 0 
   367)               !iphase_loc_p(ghosted_id) = istate
   368)               !cur_patch%aux%Global%auxvars(ghosted_id)%istate = istate
   369)             enddo
   370)           endif
   371)           initial_condition => initial_condition%next
   372)         enddo
   373)      
   374)         call VecRestoreArrayF90(field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
   375)         call VecRestoreArrayF90(field%iphas_loc,iphase_loc_p, &
   376)                                 ierr);CHKERRQ(ierr)
   377)               
   378)       case default
   379)         ! assign initial conditions values to domain
   380)         call VecGetArrayF90(field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
   381)         call VecGetArrayF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
   382)       
   383)         xx_p = UNINITIALIZED_DOUBLE
   384)       
   385)         initial_condition => cur_patch%initial_condition_list%first
   386)         do
   387)       
   388)           if (.not.associated(initial_condition)) exit
   389) 
   390)           use_dataset = PETSC_FALSE
   391)           dataset_flag = PETSC_FALSE
   392)           do idof = 1, option%nflowdof
   393)             dataset =>  initial_condition%flow_condition% &
   394)                               sub_condition_ptr(idof)%ptr%dataset
   395)             select type(dataset_ptr => dataset)
   396)               class is(dataset_gridded_hdf5_type)
   397)                 ! already mapped to flow_aux_real_var
   398)               class is(dataset_common_hdf5_type)
   399)                 use_dataset = PETSC_TRUE
   400)                 dataset_flag(idof) = PETSC_TRUE
   401)                 call ConditionControlMapDatasetToVec(realization, &
   402)                         initial_condition%flow_condition% &
   403)                           sub_condition_ptr(idof)%ptr%dataset,idof, &
   404)                         field%flow_xx,GLOBAL)
   405)               class default
   406)             end select
   407)           enddo            
   408)           if (.not.associated(initial_condition%flow_aux_real_var) .and. &
   409)               .not.associated(initial_condition%flow_condition)) then
   410)             option%io_buffer = 'Flow condition is NULL in initial condition'
   411)             call printErrMsg(option)
   412)           endif
   413)           if (associated(initial_condition%flow_aux_real_var)) then
   414)             num_connections = &
   415)               initial_condition%connection_set%num_connections
   416)             conn_id_ptr => initial_condition%connection_set%id_dn
   417)           else
   418)             num_connections = initial_condition%region%num_cells
   419)             conn_id_ptr => initial_condition%region%cell_ids
   420)           endif
   421)           do iconn=1, num_connections
   422)             local_id = conn_id_ptr(iconn)
   423)             ghosted_id = grid%nL2G(local_id)
   424)             iend = local_id*option%nflowdof
   425)             ibegin = iend-option%nflowdof+1
   426)             if (cur_patch%imat(ghosted_id) <= 0) then
   427)               xx_p(ibegin:iend) = 0.d0
   428)               iphase_loc_p(ghosted_id) = 0
   429)               cycle
   430)             endif
   431)             if (associated(initial_condition%flow_aux_real_var)) then
   432)               do idof = 1, option%nflowdof
   433)                 if (.not.dataset_flag(idof)) then
   434)                   xx_p(ibegin+idof-1) =  &
   435)                     initial_condition%flow_aux_real_var(idof,iconn)
   436)                 endif
   437)               enddo
   438)             else
   439)               do idof = 1, option%nflowdof
   440)                 if (.not.dataset_flag(idof)) then
   441)                   xx_p(ibegin+idof-1) = &
   442)                     initial_condition%flow_condition% &
   443)                       sub_condition_ptr(idof)%ptr%dataset%rarray(1)
   444)                 endif
   445)               enddo
   446)             endif
   447)             ! TODO(geh): phase out field%iphas_loc
   448)             iphase_loc_p(ghosted_id) = &
   449)               initial_condition%flow_condition%iphase
   450)             if (option%iflowmode == G_MODE) then
   451)               cur_patch%aux%Global%auxvars(ghosted_id)%istate = &
   452)                 int(iphase_loc_p(ghosted_id))
   453)             endif
   454)           enddo
   455)           initial_condition => initial_condition%next
   456)         enddo
   457)         call VecRestoreArrayF90(field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
   458) 
   459)     end select 
   460)    
   461)     cur_patch => cur_patch%next
   462)   enddo
   463)   
   464)   if (option%iflowmode == G_MODE) then
   465)     call GlobalUpdateState(realization)
   466)     call GlobalAuxVarStrip(global_aux)
   467)     call GeneralAuxVarStrip(general_aux)
   468)   endif  
   469)    
   470)   ! update dependent vectors
   471)   call DiscretizationGlobalToLocal(discretization,field%flow_xx, &
   472)                                    field%flow_xx_loc,NFLOWDOF)  
   473) 
   474)   call VecCopy(field%flow_xx, field%flow_yy, ierr);CHKERRQ(ierr)
   475)   call DiscretizationLocalToLocal(discretization,field%iphas_loc, &
   476)                                   field%iphas_loc,ONEDOF)  
   477)   call DiscretizationLocalToLocal(discretization,field%iphas_loc, &
   478)                                   field%iphas_old_loc,ONEDOF)
   479) 
   480)   ! cannot perform VecMin on local vector as the ghosted corner values are not
   481)   ! updated during the local to local update.
   482)   call DiscretizationLocalToGlobal(discretization,field%iphas_loc,field%work, &
   483)                                    ONEDOF)
   484)   call VecMin(field%work,PETSC_NULL_INTEGER,tempreal,ierr);CHKERRQ(ierr)
   485)   if (tempreal < 0.d0) then
   486) !    print *, tempreal
   487)     option%io_buffer = 'Uninitialized cells in domain.'
   488)     call printErrMsg(option)
   489)   endif
   490) 
   491) end subroutine CondControlAssignFlowInitCond
   492) 
   493) ! ************************************************************************** !
   494) 
   495) subroutine CondControlAssignTranInitCond(realization)
   496)   ! 
   497)   ! Assigns transport initial conditions to model
   498)   ! 
   499)   ! Author: Glenn Hammond
   500)   ! Date: 11/02/07, 10/18/11
   501)   ! 
   502) 
   503)   use Realization_Subsurface_class
   504)   use Discretization_module
   505)   use Region_module
   506)   use Option_module
   507)   use Field_module
   508)   use Coupler_module
   509)   use Condition_module
   510)   use Transport_Constraint_module
   511)   use Grid_module
   512)   use Dataset_Base_class
   513)   use Patch_module
   514)   use Reactive_Transport_module, only : RTUpdateAuxVars
   515)   use Reactive_Transport_Aux_module
   516)   use Reaction_Aux_module
   517)   use Global_Aux_module
   518)   use Material_Aux_class
   519)   use Reaction_module
   520)   use HDF5_module
   521)   
   522)   implicit none
   523) 
   524) #include "petsc/finclude/petscvec.h"
   525) #include "petsc/finclude/petscvec.h90"
   526)   
   527)   class(realization_subsurface_type) :: realization
   528)   
   529)   PetscInt :: icell, iconn, idof, isub_condition, temp_int, iimmobile
   530)   PetscInt :: local_id, ghosted_id, iend, ibegin
   531)   PetscInt :: irxn, isite, imnrl, ikinrxn
   532)   PetscReal, pointer :: xx_p(:), xx_loc_p(:), vec_p(:)
   533)   PetscErrorCode :: ierr
   534)   
   535)   type(option_type), pointer :: option
   536)   type(field_type), pointer :: field  
   537)   type(patch_type), pointer :: patch
   538)   type(grid_type), pointer :: grid
   539)   type(discretization_type), pointer :: discretization
   540)   type(coupler_type), pointer :: initial_condition
   541)   type(patch_type), pointer :: cur_patch
   542)   type(reaction_type), pointer :: reaction
   543)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:)
   544)   type(global_auxvar_type), pointer :: global_auxvars(:)
   545)   type(tran_constraint_coupler_type), pointer :: constraint_coupler
   546)   class(material_auxvar_type), pointer :: material_auxvars(:)
   547) 
   548)   PetscInt :: iphase
   549)   PetscInt :: offset
   550)   PetscBool :: re_equilibrate_at_each_cell
   551)   character(len=MAXSTRINGLENGTH) :: string, string2
   552)   class(dataset_base_type), pointer :: dataset
   553)   PetscInt :: aq_dataset_to_idof(realization%reaction%naqcomp)
   554)   PetscInt :: iaqdataset, num_aq_datasets
   555)   PetscBool :: use_aq_dataset
   556)   PetscReal :: ave_num_iterations
   557)   PetscReal :: tempreal
   558)   PetscInt :: prev_equilibrated_ghosted_id
   559)   PetscReal, pointer :: iphase_loc_p(:)
   560)   PetscReal, pointer :: flow_xx_p(:)
   561)   PetscLogDouble :: tstart, tend
   562)   
   563)   option => realization%option
   564)   discretization => realization%discretization
   565)   field => realization%field
   566)   patch => realization%patch
   567)   reaction => realization%reaction
   568)   
   569)   iphase = 1
   570)   
   571)   cur_patch => realization%patch_list%first
   572)   do
   573)     if (.not.associated(cur_patch)) exit
   574) 
   575)     grid => cur_patch%grid
   576)     rt_auxvars => cur_patch%aux%RT%auxvars
   577)     global_auxvars => cur_patch%aux%Global%auxvars
   578)     material_auxvars => cur_patch%aux%Material%auxvars
   579) 
   580)     ! assign initial conditions values to domain
   581)     call VecGetArrayF90(field%tran_xx,xx_p,ierr);CHKERRQ(ierr)
   582)     select case(option%iflowmode)
   583)       case(MPH_MODE,FLASH2_MODE)
   584)         call VecGetArrayReadF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
   585)         call VecGetArrayF90(field%flow_xx,flow_xx_p, ierr);CHKERRQ(ierr)
   586)     end select
   587)       
   588)     xx_p = UNINITIALIZED_DOUBLE
   589)       
   590)     initial_condition => cur_patch%initial_condition_list%first
   591)     do
   592)       
   593)       if (.not.associated(initial_condition)) exit
   594)         
   595)       constraint_coupler => initial_condition%tran_condition%cur_constraint_coupler
   596) 
   597)       re_equilibrate_at_each_cell = PETSC_FALSE
   598)       use_aq_dataset = PETSC_FALSE
   599)       num_aq_datasets = 0
   600)       aq_dataset_to_idof = 0
   601)       do idof = 1, reaction%naqcomp ! primary aqueous concentrations
   602)         if (constraint_coupler%aqueous_species%external_dataset(idof)) then
   603)           num_aq_datasets = num_aq_datasets + 1
   604)           aq_dataset_to_idof(num_aq_datasets) = idof
   605)           re_equilibrate_at_each_cell = PETSC_TRUE
   606)           use_aq_dataset = PETSC_TRUE
   607)           string = 'constraint ' // trim(constraint_coupler%constraint_name)
   608)           dataset => DatasetBaseGetPointer(realization%datasets, &
   609)                         constraint_coupler%aqueous_species%constraint_aux_string(idof), &
   610)                         string,option)
   611)           call ConditionControlMapDatasetToVec(realization,dataset,idof, &
   612)                                                 field%tran_xx_loc,LOCAL)
   613)         endif
   614)       enddo
   615) 
   616)       ! read in heterogeneous mineral volume fractions
   617)       if (associated(constraint_coupler%minerals)) then
   618)         do imnrl = 1, reaction%mineral%nkinmnrl
   619)           if (constraint_coupler%minerals%external_vol_frac_dataset(imnrl)) then
   620)             re_equilibrate_at_each_cell = PETSC_TRUE
   621)             string = 'constraint ' // trim(constraint_coupler%constraint_name)
   622)             dataset => DatasetBaseGetPointer(realization%datasets, &
   623)                           constraint_coupler%minerals% &
   624)                             constraint_vol_frac_string(imnrl), &
   625)                           string,option)
   626)             idof = ONE_INTEGER
   627)             call ConditionControlMapDatasetToVec(realization,dataset,idof, &
   628)                                                   field%work_loc,LOCAL)
   629)             call VecGetArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   630)             do icell=1,initial_condition%region%num_cells
   631)               local_id = initial_condition%region%cell_ids(icell)
   632)               ghosted_id = grid%nL2G(local_id)
   633)               rt_auxvars(ghosted_id)%mnrl_volfrac0(imnrl) = vec_p(ghosted_id)
   634)               rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl) = vec_p(ghosted_id)
   635)             enddo
   636)             call VecRestoreArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   637)           endif
   638)         enddo
   639)       endif
   640)           
   641)       ! read in heterogeneous mineral surface area
   642)       if (associated(constraint_coupler%minerals)) then
   643)         do imnrl = 1, reaction%mineral%nkinmnrl
   644)           if (constraint_coupler%minerals%external_area_dataset(imnrl)) then
   645)             re_equilibrate_at_each_cell = PETSC_TRUE
   646)             string = 'constraint ' // trim(constraint_coupler%constraint_name)
   647)             dataset => DatasetBaseGetPointer(realization%datasets, &
   648)                           constraint_coupler%minerals% &
   649)                           constraint_area_string(imnrl), &
   650)                           string,option)
   651)             idof = ONE_INTEGER
   652)             call ConditionControlMapDatasetToVec(realization,dataset,idof, &
   653)                                                   field%work_loc,LOCAL)
   654)             call VecGetArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   655)             do icell=1,initial_condition%region%num_cells
   656)               local_id = initial_condition%region%cell_ids(icell)
   657)               ghosted_id = grid%nL2G(local_id)
   658)               rt_auxvars(ghosted_id)%mnrl_area0(imnrl) = vec_p(ghosted_id)
   659)               rt_auxvars(ghosted_id)%mnrl_area(imnrl) = vec_p(ghosted_id)
   660)             enddo
   661)             call VecRestoreArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   662)           endif
   663)         enddo
   664)       endif
   665)           
   666)       ! read in heterogeneous immobile
   667)       if (associated(constraint_coupler%immobile_species)) then
   668)         do iimmobile = 1, reaction%immobile%nimmobile
   669)           if (constraint_coupler%immobile_species%external_dataset(iimmobile)) then
   670)             ! no need to requilibrate at each cell
   671)             string = 'constraint ' // trim(constraint_coupler%constraint_name)
   672)             dataset => DatasetBaseGetPointer(realization%datasets, &
   673)                 constraint_coupler%immobile_species%constraint_aux_string(iimmobile), &
   674)                 string,option)
   675)             idof = ONE_INTEGER
   676)             call ConditionControlMapDatasetToVec(realization,dataset,idof, &
   677)                                                   field%work_loc,LOCAL)
   678)             call VecGetArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   679)             do icell=1,initial_condition%region%num_cells
   680)               local_id = initial_condition%region%cell_ids(icell)
   681)               ghosted_id = grid%nL2G(local_id)
   682)               rt_auxvars(ghosted_id)%immobile(iimmobile) = vec_p(ghosted_id)
   683)             enddo
   684)             call VecRestoreArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   685)           endif
   686)         enddo
   687)       endif
   688)           
   689)       if (.not.option%use_isothermal) then
   690)         re_equilibrate_at_each_cell = PETSC_TRUE
   691)       endif
   692)         
   693)       if (use_aq_dataset) then
   694)         call VecGetArrayF90(field%tran_xx_loc,xx_loc_p,ierr);CHKERRQ(ierr)
   695)         call PetscTime(tstart,ierr);CHKERRQ(ierr)
   696)       endif
   697)         
   698)       ave_num_iterations = 0.d0
   699)       prev_equilibrated_ghosted_id = 0
   700)       do icell=1,initial_condition%region%num_cells
   701)         local_id = initial_condition%region%cell_ids(icell)
   702)         ghosted_id = grid%nL2G(local_id)
   703)         iend = local_id*option%ntrandof
   704)         ibegin = iend-option%ntrandof+1
   705)         if (cur_patch%imat(ghosted_id) <= 0) then
   706)           xx_p(ibegin:iend) = 1.d-200
   707)           cycle
   708)         endif
   709)         if (re_equilibrate_at_each_cell) then
   710)           if (use_aq_dataset) then
   711)             offset = (ghosted_id-1)*option%ntrandof
   712)             do iaqdataset = 1, num_aq_datasets
   713)               ! remember that xx_loc_p holds the data set values that were read in
   714)               temp_int = aq_dataset_to_idof(iaqdataset)
   715)               constraint_coupler%aqueous_species%constraint_conc(temp_int) = &
   716)                 xx_loc_p(offset+temp_int)
   717)             enddo
   718)           endif
   719)           option%iflag = grid%nG2A(grid%nL2G(local_id))
   720)           if (prev_equilibrated_ghosted_id == 0) then
   721)             call ReactionEquilibrateConstraint(rt_auxvars(ghosted_id), &
   722)               global_auxvars(ghosted_id),material_auxvars(ghosted_id), &
   723)               reaction, &
   724)               constraint_coupler%constraint_name, &
   725)               constraint_coupler%aqueous_species, &
   726)               constraint_coupler%free_ion_guess, &
   727)               constraint_coupler%minerals, &
   728)               constraint_coupler%surface_complexes, &
   729)               constraint_coupler%colloids, &
   730)               constraint_coupler%immobile_species, &
   731)               constraint_coupler%num_iterations, &
   732)               PETSC_FALSE,option)
   733)           else
   734)             ! copy molalities from previous equilibrated auxvar as initial guess
   735)             rt_auxvars(ghosted_id)%pri_molal = &
   736)               rt_auxvars(prev_equilibrated_ghosted_id)%pri_molal
   737)             call ReactionEquilibrateConstraint(rt_auxvars(ghosted_id), &
   738)               global_auxvars(ghosted_id),material_auxvars(ghosted_id), &
   739)               reaction, &
   740)               constraint_coupler%constraint_name, &
   741)               constraint_coupler%aqueous_species, &
   742)               constraint_coupler%free_ion_guess, &
   743)               constraint_coupler%minerals, &
   744)               constraint_coupler%surface_complexes, &
   745)               constraint_coupler%colloids, &
   746)               constraint_coupler%immobile_species, &
   747)               constraint_coupler%num_iterations, &
   748)               PETSC_TRUE,option)
   749)           endif
   750)           option%iflag = 0
   751)           ave_num_iterations = ave_num_iterations + &
   752)             constraint_coupler%num_iterations
   753)           ! update CO2 mole fraction for CO2 modes
   754) #if 0
   755)           ! TODO(geh): ideally, the intermingling of the flow process model
   756)           ! with transport is not ideal.  Peter should be looking into whether
   757)           ! we can remove this code in favor of a slighly less accurate
   758)           ! solution.
   759)           select case(option%iflowmode)
   760)             case(MPH_MODE,FLASH2_MODE)
   761)               if (int(iphase_loc_p(ghosted_id)) == 1) then
   762)                 tempreal = &
   763)                   RCO2MoleFraction(rt_auxvars(ghosted_id), &
   764)                                    global_auxvars(ghosted_id),reaction,option)
   765)                 ! concentration dof in flow solution vector
   766)                 flow_xx_p(local_id*option%nflowdof) = tempreal
   767)               endif
   768)           end select
   769) #endif
   770)           ! prev_eq_ghosted_id is only updated to the prev active cell
   771)           prev_equilibrated_ghosted_id = ghosted_id
   772)         endif
   773)         ! ibegin is the local non-ghosted offset: (local_id-1)*option%ntrandof+1
   774)         offset = ibegin + reaction%offset_aqueous - 1
   775)         ! primary aqueous concentrations
   776)         do idof = 1, reaction%naqcomp 
   777)           xx_p(offset+idof) = &
   778)             constraint_coupler%aqueous_species%basis_molarity(idof) / &
   779)             global_auxvars(ghosted_id)%den_kg(iphase)*1000.d0 ! convert molarity -> molality
   780)         enddo
   781)         ! mineral volume fractions
   782)         if (associated(constraint_coupler%minerals)) then
   783)           do imnrl = 1, reaction%mineral%nkinmnrl
   784)             ! if read from a dataset, the vol frac was set above.  Don't want to
   785)             ! overwrite
   786)             if (.not.constraint_coupler%minerals% &
   787)                   external_vol_frac_dataset(imnrl)) then
   788)               rt_auxvars(ghosted_id)%mnrl_volfrac0(imnrl) = &
   789)                 constraint_coupler%minerals%constraint_vol_frac(imnrl)
   790)               rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl) = &
   791)                 constraint_coupler%minerals%constraint_vol_frac(imnrl)
   792)             endif
   793)             if (.not.constraint_coupler%minerals% &
   794)                   external_area_dataset(imnrl)) then
   795)               rt_auxvars(ghosted_id)%mnrl_area0(imnrl) = &
   796)                 constraint_coupler%minerals%constraint_area(imnrl)
   797)               rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
   798)                 constraint_coupler%minerals%constraint_area(imnrl)
   799)             endif
   800)           enddo
   801)         endif
   802)         ! kinetic surface complexes
   803)         if (associated(constraint_coupler%surface_complexes)) then
   804)           do idof = 1, reaction%surface_complexation%nkinsrfcplx
   805)             rt_auxvars(ghosted_id)%kinsrfcplx_conc(idof,-1) = & !geh: to catch bug
   806)               constraint_coupler%surface_complexes%constraint_conc(idof)
   807)           enddo
   808)           do ikinrxn = 1, reaction%surface_complexation%nkinsrfcplxrxn
   809)             irxn = reaction%surface_complexation%kinsrfcplxrxn_to_srfcplxrxn(ikinrxn)
   810)             isite = reaction%surface_complexation%srfcplxrxn_to_surf(irxn)
   811)             rt_auxvars(ghosted_id)%kinsrfcplx_free_site_conc(isite) = &
   812)               constraint_coupler%surface_complexes%basis_free_site_conc(isite)
   813)           enddo
   814)         endif
   815)         ! this is for the multi-rate surface complexation model
   816)         if (reaction%surface_complexation%nkinmrsrfcplxrxn > 0 .and. &
   817)           ! geh: if we re-equilibrate at each grid cell, we do not want to
   818)           ! overwrite the reequilibrated values with those from the constraint
   819)             .not. re_equilibrate_at_each_cell) then
   820)           ! copy over total sorbed concentration
   821)           rt_auxvars(ghosted_id)%kinmr_total_sorb = &
   822)             constraint_coupler%rt_auxvar%kinmr_total_sorb
   823)           ! copy over free site concentration
   824)           rt_auxvars(ghosted_id)%srfcplxrxn_free_site_conc = &
   825)             constraint_coupler%rt_auxvar%srfcplxrxn_free_site_conc
   826)         endif
   827)         ! colloids fractions
   828)         if (associated(constraint_coupler%colloids)) then
   829)           offset = ibegin + reaction%offset_colloid - 1
   830)           do idof = 1, reaction%ncoll ! primary aqueous concentrations
   831)             xx_p(offset+idof) = &
   832)               constraint_coupler%colloids%basis_conc_mob(idof) / &
   833)               global_auxvars(ghosted_id)%den_kg(iphase)*1000.d0 ! convert molarity -> molality
   834)             rt_auxvars(ghosted_id)%colloid%conc_imb(idof) = &
   835)               constraint_coupler%colloids%basis_conc_imb(idof)
   836)           enddo
   837)         endif
   838)         ! immobile
   839)         if (associated(constraint_coupler%immobile_species)) then
   840)           offset = ibegin + reaction%offset_immobile - 1
   841)           do iimmobile = 1, reaction%immobile%nimmobile
   842)             if (constraint_coupler%immobile_species%external_dataset(iimmobile)) then
   843)               ! already read into rt_auxvars above.
   844)               xx_p(offset+iimmobile) = &
   845)                 rt_auxvars(ghosted_id)%immobile(iimmobile)
   846)             else
   847)               xx_p(offset+iimmobile) = &
   848)                 constraint_coupler%immobile_species%constraint_conc(iimmobile)
   849)               rt_auxvars(ghosted_id)%immobile(iimmobile) = &
   850)                 constraint_coupler%immobile_species%constraint_conc(iimmobile)
   851)             endif
   852)           enddo
   853)         endif
   854)       enddo ! icell=1,initial_condition%region%num_cells
   855)       if (use_aq_dataset) then
   856)         call PetscTime(tend,ierr);CHKERRQ(ierr)
   857)         call VecRestoreArrayF90(field%tran_xx_loc,xx_loc_p,ierr);CHKERRQ(ierr)
   858)         ave_num_iterations = ave_num_iterations / &
   859)           initial_condition%region%num_cells
   860)         write(option%io_buffer,&
   861)               '("Average number of iterations in ReactionEquilibrateConstraint():", &
   862)               & f5.1)') ave_num_iterations
   863)         call printMsg(option)
   864)         write(option%io_buffer,'(f10.2," Seconds to equilibrate constraints")') &
   865)           tend-tstart
   866)         call printMsg(option)
   867)       endif
   868)       initial_condition => initial_condition%next
   869)     enddo
   870)       
   871)     call VecRestoreArrayF90(field%tran_xx,xx_p, ierr);CHKERRQ(ierr)
   872)     select case(option%iflowmode)
   873)       case(MPH_MODE,FLASH2_MODE)
   874)         call VecRestoreArrayReadF90(field%iphas_loc,iphase_loc_p,ierr);CHKERRQ(ierr)
   875)         call VecRestoreArrayF90(field%flow_xx,flow_xx_p, ierr);CHKERRQ(ierr)
   876)     end select
   877) 
   878)     cur_patch => cur_patch%next
   879)   enddo
   880)   
   881)   ! check to ensure that minimum concentration is not less than or equal
   882)   ! to zero
   883)   call VecMin(field%tran_xx,PETSC_NULL_INTEGER,tempreal,ierr);CHKERRQ(ierr)
   884)   if (tempreal <= 0.d0) then
   885)     option%io_buffer = 'ERROR: Zero concentrations found in initial ' // &
   886)       'transport solution.'
   887)     call printMsg(option)
   888)     ! now figure out which species have zero concentrations
   889)     do idof = 1, option%ntrandof
   890)       call VecStrideMin(field%tran_xx,idof-1,offset,tempreal, &
   891)                         ierr);CHKERRQ(ierr)
   892)       if (tempreal <= 0.d0) then
   893)         write(string,*) tempreal
   894)         if (idof <= reaction%naqcomp) then
   895)           string2 = '  Aqueous species "' // &
   896)             trim(reaction%primary_species_names(idof))
   897)         else
   898)           string2 = '  Immobile species "' // &
   899)             trim(reaction%immobile%names(idof-reaction%offset_immobile))
   900)         endif
   901)           string2 = trim(string2) // &
   902)             '" has zero concentration (' // &
   903)             trim(adjustl(string)) // ').'
   904)         call printMsg(option,string2)
   905)       endif
   906)     enddo
   907)     option%io_buffer = ''
   908)     call printMsg(option)
   909)     option%io_buffer = '*** Begin Note'
   910)     call printMsg(option)
   911)     option%io_buffer = 'If concentrations = -999., they have not ' // &
   912)               'been initialized properly.'
   913)     call printMsg(option)
   914)     option%io_buffer = '*** End Note'
   915)     call printMsg(option)
   916)     option%io_buffer = 'Free ion concentations must be positive.  Try ' // &
   917)       'using a small value such as 1.e-20 or 1.e-40 instead of zero.'
   918)     call printErrMsg(option)
   919)   endif
   920)   
   921)   ! update dependent vectors
   922)   call DiscretizationGlobalToLocal(discretization,field%tran_xx, &
   923)                                    field%tran_xx_loc,NTRANDOF)  
   924)   call VecCopy(field%tran_xx, field%tran_yy, ierr);CHKERRQ(ierr)
   925) 
   926)   ! override initial conditions if they are to be read from a file
   927)   if (len_trim(option%initialize_transport_filename) > 1) then
   928)     call CondControlReadTransportIC(realization, &
   929)                                     option%initialize_transport_filename)
   930)   endif
   931)   ! PETSC_FALSE = no activity coefficients
   932)   call RTUpdateAuxVars(realization,PETSC_TRUE,PETSC_FALSE,PETSC_FALSE)
   933)   ! at this point the auxvars have been computed with activity coef = 1.d0
   934)   ! to use intitial condition with activity coefs /= 1.d0, must update
   935)   ! activity coefs and recompute auxvars
   936)   if (realization%reaction%act_coef_update_frequency /= &
   937)       ACT_COEF_FREQUENCY_OFF) then
   938)     call RTUpdateAuxVars(realization,PETSC_TRUE,PETSC_FALSE,PETSC_TRUE)
   939)     !geh: you may ask, why call this twice....  We need to iterate at least
   940)     !     once to ensure that the activity coefficients are more accurate.
   941)     !     Otherwise, the total component concentrations can be quite
   942)     !     different from what is defined in the input file.
   943)     call RTUpdateAuxVars(realization,PETSC_TRUE,PETSC_FALSE,PETSC_TRUE)
   944)   endif
   945) 
   946) end subroutine CondControlAssignTranInitCond
   947) 
   948) ! ************************************************************************** !
   949) 
   950) subroutine ConditionControlMapDatasetToVec(realization,dataset,idof, &
   951)                                            mdof_vec,vec_type)
   952)   ! 
   953)   ! maps an external dataset to a PETSc vec
   954)   ! representing values at each grid cell
   955)   ! 
   956)   ! Author: Glenn Hammond
   957)   ! Date: 03/23/12
   958)   ! 
   959)   use Realization_Subsurface_class
   960)   use Option_module
   961)   use Field_module
   962)   use Dataset_Common_HDF5_class
   963)   use Dataset_Base_class
   964)   use HDF5_module
   965)   use Discretization_module
   966) 
   967)   implicit none
   968)   
   969) #include "petsc/finclude/petscvec.h"
   970) #include "petsc/finclude/petscvec.h90"  
   971)   
   972)   class(realization_subsurface_type) :: realization
   973)   class(dataset_base_type), pointer :: dataset
   974)   PetscInt :: idof
   975)   Vec :: mdof_vec
   976)   PetscInt :: vec_type
   977) 
   978)   type(field_type), pointer :: field
   979)   type(option_type), pointer :: option
   980)   character(len=MAXSTRINGLENGTH) :: string, string2
   981)   PetscErrorCode :: ierr
   982) 
   983)   field => realization%field
   984)   option => realization%option
   985)   
   986)   if (associated(dataset)) then
   987)     select type(dataset)
   988)       class is (dataset_common_hdf5_type)
   989)         string = '' ! group name
   990)         ! have to copy to string2 due to mismatch in string size
   991)         string2 = dataset%hdf5_dataset_name
   992)         call HDF5ReadCellIndexedRealArray(realization,field%work, &
   993)                                           dataset%filename, &
   994)                                           string,string2, &
   995)                                           dataset%realization_dependent)
   996)         if (vec_type == GLOBAL) then
   997)           call VecStrideScatter(field%work,idof-1,mdof_vec, &
   998)                                 INSERT_VALUES,ierr);CHKERRQ(ierr)
   999)         else
  1000)           call DiscretizationGlobalToLocal(realization%discretization, &
  1001)                                            field%work, &
  1002)                                            field%work_loc,ONEDOF)
  1003)           call VecStrideScatter(field%work_loc,idof-1,mdof_vec, &
  1004)                                 INSERT_VALUES,ierr);CHKERRQ(ierr)
  1005)         endif
  1006)       class default
  1007)         option%io_buffer = 'Dataset "' // trim(dataset%name) // &
  1008)           '" not supported in ConditionControlMapDatasetToVec.'
  1009)         call printErrMsg(option)
  1010)     end select
  1011)   endif
  1012) 
  1013) end subroutine ConditionControlMapDatasetToVec
  1014) 
  1015) ! ************************************************************************** !
  1016) 
  1017) subroutine CondControlScaleSourceSink(realization)
  1018)   ! 
  1019)   ! Scales select source/sinks based on perms
  1020)   ! 
  1021)   ! Author: Glenn Hammond
  1022)   ! Date: 09/03/08, 10/18/11
  1023)   ! 
  1024) 
  1025)   use Realization_Subsurface_class
  1026)   use Discretization_module
  1027)   use Region_module
  1028)   use Option_module
  1029)   use Field_module
  1030)   use Coupler_module
  1031)   use Connection_module
  1032)   use Condition_module
  1033)   use Grid_module
  1034)   use Patch_module
  1035)   use Material_Aux_class
  1036)   use Variables_module, only : PERMEABILITY_X
  1037) 
  1038)   implicit none
  1039) 
  1040) #include "petsc/finclude/petscvec.h"
  1041) #include "petsc/finclude/petscvec.h90"
  1042) #include "petsc/finclude/petscdmda.h"
  1043) 
  1044)   
  1045)   class(realization_subsurface_type) :: realization
  1046)   
  1047)   PetscErrorCode :: ierr
  1048)   
  1049)   type(option_type), pointer :: option
  1050)   type(field_type), pointer :: field  
  1051)   type(patch_type), pointer :: patch
  1052)   type(grid_type), pointer :: grid
  1053)   type(discretization_type), pointer :: discretization
  1054)   type(coupler_type), pointer :: cur_source_sink
  1055)   type(connection_set_type), pointer :: cur_connection_set
  1056)   class(material_auxvar_type), pointer :: material_auxvars(:)
  1057)   type(patch_type), pointer :: cur_patch
  1058)   PetscReal, pointer :: vec_ptr(:)
  1059)   PetscInt :: local_id
  1060)   PetscInt :: ghosted_id, neighbor_ghosted_id
  1061)   PetscInt :: iconn
  1062)   PetscReal :: scale, sum
  1063)   PetscInt :: icount
  1064)   PetscInt :: x_count, y_count, z_count
  1065)   PetscInt, parameter :: x_width = 1, y_width = 1, z_width = 0
  1066)   
  1067)   PetscInt :: ghosted_neighbors(0:27)
  1068)   
  1069)   option => realization%option
  1070)   discretization => realization%discretization
  1071)   field => realization%field
  1072)   patch => realization%patch
  1073)   material_auxvars => realization%patch%aux%Material%auxvars
  1074)   
  1075)   ! GB: grid was uninitialized
  1076)   grid => patch%grid
  1077) 
  1078)   cur_patch => realization%patch_list%first
  1079)   do
  1080)     if (.not.associated(cur_patch)) exit
  1081)     ! BIG-TIME warning here.  I assume that all source/sink cells are within 
  1082)     ! a single patch - geh
  1083) 
  1084)     grid => cur_patch%grid
  1085) 
  1086)     cur_source_sink => cur_patch%source_sink_list%first
  1087)     do
  1088)       if (.not.associated(cur_source_sink)) exit
  1089) 
  1090)       call VecZeroEntries(field%work,ierr);CHKERRQ(ierr)
  1091)       call VecGetArrayF90(field%work,vec_ptr,ierr);CHKERRQ(ierr)
  1092) 
  1093)       cur_connection_set => cur_source_sink%connection_set
  1094)     
  1095)       do iconn = 1, cur_connection_set%num_connections
  1096)         local_id = cur_connection_set%id_dn(iconn)
  1097)         ghosted_id = grid%nL2G(local_id)
  1098) 
  1099)         select case(option%iflowmode)
  1100)           case(RICHARDS_MODE,G_MODE)
  1101)               call GridGetGhostedNeighbors(grid,ghosted_id,DMDA_STENCIL_STAR, &
  1102)                                           x_width,y_width,z_width, &
  1103)                                           x_count,y_count,z_count, &
  1104)                                           ghosted_neighbors,option)
  1105)               ! ghosted neighbors is ordered first in x, then, y, then z
  1106)               icount = 0
  1107)               sum = 0.d0
  1108)               ! x-direction
  1109)               do while (icount < x_count)
  1110)                 icount = icount + 1
  1111)                 neighbor_ghosted_id = ghosted_neighbors(icount)
  1112)                 sum = sum + MaterialAuxVarGetValue(material_auxvars( &
  1113)                               neighbor_ghosted_id),PERMEABILITY_X) * &
  1114)                             grid%structured_grid%dy(neighbor_ghosted_id)* &
  1115)                             grid%structured_grid%dz(neighbor_ghosted_id)
  1116)                  
  1117)               enddo
  1118)               ! y-direction
  1119)               do while (icount < x_count + y_count)
  1120)                 icount = icount + 1
  1121)                 neighbor_ghosted_id = ghosted_neighbors(icount)                 
  1122)                 sum = sum + MaterialAuxVarGetValue(material_auxvars( &
  1123)                               neighbor_ghosted_id),PERMEABILITY_X) * &
  1124)                             grid%structured_grid%dx(neighbor_ghosted_id)* &
  1125)                             grid%structured_grid%dz(neighbor_ghosted_id)
  1126)                  
  1127)               enddo
  1128)               ! z-direction
  1129)               do while (icount < x_count + y_count + z_count)
  1130)                 icount = icount + 1
  1131)                 neighbor_ghosted_id = ghosted_neighbors(icount)                 
  1132)                 sum = sum + MaterialAuxVarGetValue(material_auxvars( &
  1133)                               neighbor_ghosted_id),PERMEABILITY_X) * &
  1134)                             grid%structured_grid%dx(neighbor_ghosted_id)* &
  1135)                             grid%structured_grid%dy(neighbor_ghosted_id)
  1136)               enddo
  1137)               vec_ptr(local_id) = vec_ptr(local_id) + sum
  1138)           case(TH_MODE)
  1139)           case(MPH_MODE)
  1140)           case(IMS_MODE)
  1141)           case(MIS_MODE)
  1142)           case(FLASH2_MODE)
  1143)         end select
  1144) 
  1145)       enddo
  1146)         
  1147)       call VecRestoreArrayF90(field%work,vec_ptr,ierr);CHKERRQ(ierr)
  1148)       call VecNorm(field%work,NORM_1,scale,ierr);CHKERRQ(ierr)
  1149)       scale = 1.d0/scale
  1150)       call VecScale(field%work,scale,ierr);CHKERRQ(ierr)
  1151) 
  1152)       call VecGetArrayF90(field%work,vec_ptr, ierr);CHKERRQ(ierr)
  1153)       do iconn = 1, cur_connection_set%num_connections      
  1154)         local_id = cur_connection_set%id_dn(iconn)
  1155)         select case(option%iflowmode)
  1156)           case(RICHARDS_MODE,G_MODE)
  1157)             cur_source_sink%flow_aux_real_var(ONE_INTEGER,iconn) = &
  1158)               vec_ptr(local_id)
  1159)           case(TH_MODE)
  1160)           case(MPH_MODE)
  1161)           case(IMS_MODE)
  1162)           case(MIS_MODE)
  1163)           case(FLASH2_MODE)
  1164)         end select
  1165) 
  1166)       enddo
  1167)       call VecRestoreArrayF90(field%work,vec_ptr,ierr);CHKERRQ(ierr)
  1168)         
  1169)       cur_source_sink => cur_source_sink%next
  1170)     enddo
  1171)     cur_patch => cur_patch%next
  1172)   enddo
  1173) 
  1174) end subroutine CondControlScaleSourceSink
  1175) 
  1176) ! ************************************************************************** !
  1177) 
  1178) subroutine CondControlReadTransportIC(realization,filename)
  1179)   ! 
  1180)   ! Assigns transport initial condition from
  1181)   ! HDF5 file
  1182)   ! 
  1183)   ! Author: Glenn Hammond
  1184)   ! Date: 03/05/10
  1185)   ! 
  1186) 
  1187)   use Realization_Subsurface_class
  1188)   use Option_module
  1189)   use Field_module
  1190)   use Grid_module
  1191)   use Patch_module
  1192)   use Reactive_Transport_module
  1193)   use Reaction_Aux_module
  1194)   use Discretization_module
  1195)   use HDF5_module
  1196)   
  1197)   implicit none
  1198) 
  1199) #include "petsc/finclude/petscvec.h"
  1200) #include "petsc/finclude/petscvec.h90"
  1201)   
  1202)   class(realization_subsurface_type) :: realization
  1203)   character(len=MAXSTRINGLENGTH) :: filename
  1204)   
  1205)   PetscInt :: local_id, idx, offset, idof
  1206)   PetscReal, pointer :: xx_p(:)
  1207)   character(len=MAXSTRINGLENGTH) :: group_name
  1208)   character(len=MAXSTRINGLENGTH) :: dataset_name
  1209)   PetscReal, pointer :: vec_p(:)  
  1210)   PetscErrorCode :: ierr
  1211)   
  1212)   type(option_type), pointer :: option
  1213)   type(field_type), pointer :: field  
  1214)   type(patch_type), pointer :: patch
  1215)   type(grid_type), pointer :: grid
  1216)   type(discretization_type), pointer :: discretization
  1217)   type(patch_type), pointer :: cur_patch
  1218)   type(reaction_type), pointer :: reaction
  1219) 
  1220)   option => realization%option
  1221)   discretization => realization%discretization
  1222)   field => realization%field
  1223)   patch => realization%patch
  1224)   reaction => realization%reaction
  1225) 
  1226)   cur_patch => realization%patch_list%first
  1227)   do
  1228)     if (.not.associated(cur_patch)) exit
  1229) 
  1230)     grid => cur_patch%grid
  1231) 
  1232)       ! assign initial conditions values to domain
  1233)     call VecGetArrayF90(field%tran_xx,xx_p, ierr);CHKERRQ(ierr)
  1234) 
  1235)     ! Primary species concentrations for all modes 
  1236)     do idof = 1, option%ntrandof ! primary aqueous concentrations
  1237)       offset = idof
  1238)       group_name = ''
  1239)       dataset_name = reaction%primary_species_names(idof)
  1240)       call HDF5ReadCellIndexedRealArray(realization,field%work, &
  1241)                                         filename,group_name, &
  1242)                                         dataset_name,option%id>0)
  1243)       call VecGetArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1244)       do local_id=1, grid%nlmax
  1245)         if (cur_patch%imat(grid%nL2G(local_id)) <= 0) cycle
  1246)         if (vec_p(local_id) < 1.d-40) then
  1247)           print *,  option%myrank, grid%nG2A(grid%nL2G(local_id)), &
  1248)             ': Zero free-ion concentration in Initial Condition read from file.'
  1249)         endif
  1250)         idx = (local_id-1)*option%ntrandof + offset
  1251)         xx_p(idx) = vec_p(local_id)
  1252)       enddo
  1253)       call VecRestoreArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1254)      
  1255)     enddo     
  1256) 
  1257)     call VecRestoreArrayF90(field%tran_xx,xx_p, ierr);CHKERRQ(ierr)
  1258)         
  1259)     cur_patch => cur_patch%next
  1260)   enddo
  1261)    
  1262)   ! update dependent vectors
  1263)   call DiscretizationGlobalToLocal(discretization,field%tran_xx, &
  1264)                                    field%tran_xx_loc,NTRANDOF)  
  1265)   call VecCopy(field%tran_xx, field%tran_yy, ierr);CHKERRQ(ierr)
  1266)   
  1267) end subroutine CondControlReadTransportIC
  1268) 
  1269) ! ************************************************************************** !
  1270) 
  1271) subroutine CondControlAssignFlowInitCondSurface(surf_realization)
  1272) 
  1273)   use Realization_Surface_class
  1274)   use Discretization_module
  1275)   use Region_module
  1276)   use Option_module
  1277)   use Surface_Field_module
  1278)   use Coupler_module
  1279)   use Condition_module
  1280)   use Grid_module
  1281)   use Patch_module
  1282)   use EOS_Water_module
  1283)   use Surface_TH_Aux_module
  1284)   use Surface_Global_Aux_module
  1285)   
  1286)   implicit none
  1287) 
  1288) #include "petsc/finclude/petscvec.h"
  1289) #include "petsc/finclude/petscvec.h90"
  1290)   
  1291)   class(realization_surface_type) :: surf_realization
  1292)   
  1293)   PetscInt :: icell, iconn, idof, iface
  1294)   PetscInt :: local_id, ghosted_id, iend, ibegin
  1295)   PetscReal, pointer :: xx_p(:)!, iphase_loc_p(:)
  1296)   PetscErrorCode :: ierr
  1297)   
  1298)   PetscReal :: temperature, p_sat
  1299)   PetscReal :: pw, dw_kg, dw_mol, hw
  1300)   PetscReal :: temp
  1301)   PetscReal :: dpsat_dt
  1302)   character(len=MAXSTRINGLENGTH) :: string
  1303)   
  1304)   type(option_type), pointer :: option
  1305)   type(surface_field_type), pointer :: surf_field  
  1306)   type(patch_type), pointer :: patch
  1307)   type(grid_type), pointer :: grid
  1308)   type(discretization_type), pointer :: discretization
  1309)   type(coupler_type), pointer :: initial_condition
  1310)   type(patch_type), pointer :: cur_patch
  1311)   type(flow_general_condition_type), pointer :: general
  1312)   type(Surface_TH_auxvar_type), pointer :: surf_th_auxvars(:)
  1313)   type(surface_global_auxvar_type), pointer :: surf_global_auxvars(:)
  1314) 
  1315)   option => surf_realization%option
  1316)   discretization => surf_realization%discretization
  1317)   surf_field => surf_realization%surf_field
  1318)   patch => surf_realization%patch
  1319) 
  1320)   if (option%iflowmode == TH_MODE) then
  1321)     surf_th_auxvars => patch%surf_aux%SurfaceTH%auxvars
  1322)     surf_global_auxvars => patch%surf_aux%SurfaceGlobal%auxvars
  1323)   endif
  1324) 
  1325)   cur_patch => surf_realization%patch_list%first
  1326)   do
  1327)     if (.not.associated(cur_patch)) exit
  1328) 
  1329)     grid => cur_patch%grid
  1330) 
  1331)     select case(option%iflowmode)
  1332)       
  1333)       case(G_MODE) ! general phase mode
  1334) 
  1335)       case (RICHARDS_MODE,TH_MODE)
  1336)         ! assign initial conditions values to domain
  1337)         call VecGetArrayF90(surf_field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
  1338)     
  1339)         xx_p = UNINITIALIZED_DOUBLE
  1340)       
  1341)         initial_condition => cur_patch%initial_condition_list%first
  1342)         do
  1343)       
  1344)           if (.not.associated(initial_condition)) exit
  1345) 
  1346)             if (.not.associated(initial_condition%flow_aux_real_var)) then
  1347)               if (.not.associated(initial_condition%flow_condition)) then
  1348)                 option%io_buffer = 'Flow condition is NULL in initial condition'
  1349)                 call printErrMsg(option)
  1350)               endif
  1351)               do icell=1,initial_condition%region%num_cells
  1352)                 local_id = initial_condition%region%cell_ids(icell)
  1353)                 ghosted_id = grid%nL2G(local_id)
  1354)                 iend = local_id*option%nflowdof
  1355)                 ibegin = iend-option%nflowdof+1
  1356)                 if (cur_patch%imat(ghosted_id) <= 0) then
  1357)                   xx_p(ibegin:iend) = 0.d0
  1358)                   cycle
  1359)                 endif
  1360)                 do idof = 1, option%nflowdof
  1361)                   select case (idof)
  1362)                     case (ONE_INTEGER)
  1363)                       xx_p(ibegin+idof-1) = &
  1364)                         initial_condition%flow_condition% &
  1365)                         sub_condition_ptr(idof)%ptr%dataset%rarray(1)
  1366)                     case (TWO_INTEGER)
  1367)                       temp = &
  1368)                         initial_condition%flow_condition% &
  1369)                         sub_condition_ptr(idof)%ptr%dataset%rarray(1)
  1370)                       pw = option%reference_pressure
  1371)                         
  1372)                       call EOSWaterDensity(temp,pw,dw_kg,dw_mol,ierr)
  1373)                       ! [rho*h*T*Cw]
  1374)                       xx_p(ibegin+idof-1) = dw_kg*xx_p(ibegin)* &
  1375)                                             (temp + 273.15d0)* &
  1376)                                             surf_th_auxvars(ghosted_id)%Cw
  1377)                       surf_global_auxvars(ghosted_id)%den_kg(1) = dw_kg
  1378)                       surf_global_auxvars(ghosted_id)%temp = temp
  1379)                   end select
  1380)                 enddo
  1381)                 xx_p(ibegin:iend) = 0.0d0
  1382)               enddo
  1383)             else
  1384)               do iconn=1,initial_condition%connection_set%num_connections
  1385)                 local_id = initial_condition%connection_set%id_dn(iconn)
  1386)                 ghosted_id = grid%nL2G(local_id)
  1387)                 iend = local_id*option%nflowdof
  1388)                 ibegin = iend-option%nflowdof+1
  1389)                 if (cur_patch%imat(ghosted_id) <= 0) then
  1390)                   xx_p(ibegin:iend) = 0.d0
  1391)                   cycle
  1392)                 endif
  1393)                 xx_p(ibegin:iend) = &
  1394)                       initial_condition%flow_aux_real_var(1:option%nflowdof,iconn)
  1395)                 xx_p(ibegin:iend) = 0.0d0
  1396)               enddo
  1397)             endif
  1398)           initial_condition => initial_condition%next
  1399)         enddo
  1400)      
  1401)         call VecRestoreArrayF90(surf_field%flow_xx,xx_p, ierr);CHKERRQ(ierr)
  1402)       case default
  1403)         option%io_buffer = 'CondControlAssignFlowInitCondSurface not ' // &
  1404)           'for this mode'
  1405)         call printErrMsg(option)
  1406)     end select 
  1407)    
  1408)     cur_patch => cur_patch%next
  1409)   enddo
  1410)    
  1411)   ! update dependent vectors
  1412)   call DiscretizationGlobalToLocal(discretization, surf_field%flow_xx, &
  1413)                                    surf_field%flow_xx_loc, NFLOWDOF)
  1414) 
  1415) end subroutine CondControlAssignFlowInitCondSurface
  1416) 
  1417) end module Condition_Control_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody