init_subsurface.F90       coverage:  100.00 %func     72.88 %block


     1) module Init_Subsurface_module
     2) 
     3)   !geh: there can be no dependencies on simulation object in this file
     4)   use PFLOTRAN_Constants_module
     5) 
     6)   implicit none
     7) 
     8)   private
     9) 
    10) #include "petsc/finclude/petscsys.h"
    11) 
    12)   public :: InitSubsurfAssignMatIDsToRegns, &
    13)             InitSubsurfAssignMatProperties, &
    14)             SubsurfInitMaterialProperties, &
    15)             SubsurfAssignVolsToMatAuxVars, &
    16)             SubsurfSandboxesSetup, &
    17)             InitSubsurfaceSetupZeroArrays
    18)   
    19) contains
    20) 
    21) ! ************************************************************************** !
    22) 
    23) subroutine SubsurfInitMaterialProperties(realization)
    24)   ! 
    25)   ! Initializes material property data structres and assign them to the domain.
    26)   ! 
    27)   ! Author: Glenn Hammond
    28)   ! Date: 10/07/14
    29)   ! 
    30)   use Realization_Subsurface_class
    31)   
    32)   implicit none
    33)   
    34)   class(realization_subsurface_type) :: realization
    35)   
    36)   call SubsurfAllocMatPropDataStructs(realization)
    37)   call InitSubsurfAssignMatIDsToRegns(realization)
    38)   call InitSubsurfAssignMatProperties(realization)
    39)   
    40) end subroutine SubsurfInitMaterialProperties
    41) 
    42) ! ************************************************************************** !
    43) 
    44) subroutine SubsurfAllocMatPropDataStructs(realization)
    45)   ! 
    46)   ! Allocates data structures associated with storage of material properties
    47)   ! 
    48)   ! Author: Glenn Hammond
    49)   ! Date: 10/07/14
    50)   ! 
    51)   use Realization_Subsurface_class
    52)   use Material_module
    53)   use Option_module
    54)   use Discretization_module
    55)   use Grid_module
    56)   use Patch_module
    57)   use Material_Aux_class
    58)   
    59)   implicit none
    60)   
    61)   class(realization_subsurface_type) :: realization
    62)   
    63)   PetscInt :: ghosted_id
    64)   PetscInt :: istart, iend
    65)   PetscInt :: i
    66)   
    67)   type(option_type), pointer :: option
    68)   type(grid_type), pointer :: grid
    69)   type(patch_type), pointer :: cur_patch
    70)   class(material_auxvar_type), pointer :: material_auxvars(:)
    71) 
    72)   option => realization%option
    73) 
    74)   ! initialize material auxiliary indices.  this does not have to be done
    75)   ! for each patch, just once.
    76)   call MaterialInitAuxIndices(realization%patch%material_property_array,option)
    77) 
    78)   ! create mappinging
    79)   ! loop over all patches and allocation material id arrays
    80)   cur_patch => realization%patch_list%first
    81)   do
    82)     if (.not.associated(cur_patch)) exit
    83)     grid => cur_patch%grid
    84)     if (.not.associated(cur_patch%imat)) then
    85)       allocate(cur_patch%imat(grid%ngmax))
    86)       ! initialize to "unset"
    87)       cur_patch%imat = UNINITIALIZED_INTEGER
    88)       ! also allocate saturation function id
    89)       allocate(cur_patch%sat_func_id(grid%ngmax))
    90)       cur_patch%sat_func_id = UNINITIALIZED_INTEGER
    91)     endif
    92)     
    93)     cur_patch%aux%Material => MaterialAuxCreate()
    94)     allocate(material_auxvars(grid%ngmax))
    95)     do ghosted_id = 1, grid%ngmax
    96)       call MaterialAuxVarInit(material_auxvars(ghosted_id),option)
    97)     enddo
    98)     cur_patch%aux%Material%num_aux = grid%ngmax
    99)     cur_patch%aux%Material%auxvars => material_auxvars
   100)     nullify(material_auxvars)
   101)     
   102)     cur_patch => cur_patch%next
   103)   enddo
   104) 
   105)   ! Create Vec that holds compressibility
   106)   if (soil_compressibility_index > 0) then
   107)     call DiscretizationDuplicateVector(realization%discretization, &
   108)                                        realization%field%work, &
   109)                                        realization%field%compressibility0)
   110)   endif
   111) 
   112) end subroutine SubsurfAllocMatPropDataStructs
   113) 
   114) ! ************************************************************************** !
   115) 
   116) subroutine InitSubsurfAssignMatIDsToRegns(realization)
   117)   ! 
   118)   ! Assigns material properties to associated regions in the model
   119)   ! 
   120)   ! Author: Glenn Hammond
   121)   ! Date: 11/02/07
   122)   ! 
   123)   use Realization_Subsurface_class
   124)   use Strata_module
   125)   use Region_module
   126)   use Option_module
   127)   use Grid_module
   128)   use Patch_module
   129)   use Field_module
   130)   use Material_module
   131)   use Material_Aux_class
   132)   
   133)   implicit none
   134)   
   135)   class(realization_subsurface_type) :: realization
   136)   
   137)   PetscInt :: icell, local_id, ghosted_id
   138)   PetscInt :: istart, iend
   139)   PetscInt :: local_min, global_min
   140)   PetscErrorCode :: ierr
   141)   
   142)   type(option_type), pointer :: option
   143)   type(grid_type), pointer :: grid
   144)   type(field_type), pointer :: field
   145)   type(strata_type), pointer :: strata
   146)   type(patch_type), pointer :: cur_patch
   147) 
   148)   type(material_property_type), pointer :: material_property
   149)   type(region_type), pointer :: region
   150)   class(material_auxvar_type), pointer :: material_auxvars(:)
   151)   
   152)   option => realization%option
   153) 
   154)   cur_patch => realization%patch_list%first
   155)   do
   156)     if (.not.associated(cur_patch)) exit
   157)     ! set material ids to uninitialized
   158)     material_auxvars => cur_patch%aux%Material%auxvars
   159)     cur_patch%imat = UNINITIALIZED_INTEGER
   160)     grid => cur_patch%grid
   161)     strata => cur_patch%strata_list%first
   162)     do
   163)       if (.not.associated(strata)) exit
   164)       ! if not within time period specified, skip the strata.
   165)       ! use a one second tolerance on the start time and end time
   166)       if (StrataWithinTimePeriod(strata,option%time)) then
   167)         ! Read in cell by cell material ids if they exist
   168)         if (.not.associated(strata%region) .and. strata%active) then
   169)           call SubsurfReadMaterialIDsFromFile(realization, &
   170)                                               strata%realization_dependent, &
   171)                                               strata%material_property_filename)
   172)         ! Otherwise, set based on region
   173)         else
   174)           region => strata%region
   175)           material_property => strata%material_property
   176)           if (associated(region)) then
   177)             istart = 1
   178)             iend = region%num_cells
   179)           else
   180)             istart = 1
   181)             iend = grid%nlmax
   182)           endif
   183)           do icell=istart, iend
   184)             if (associated(region)) then
   185)               local_id = region%cell_ids(icell)
   186)             else
   187)               local_id = icell
   188)             endif
   189)             ghosted_id = grid%nL2G(local_id)
   190)             if (strata%active) then
   191)               cur_patch%imat(ghosted_id) = material_property%internal_id
   192)             else
   193)               ! if not active, set material id to zero
   194)               cur_patch%imat(ghosted_id) = 0
   195)             endif
   196)           enddo
   197)         endif
   198)       endif
   199)       strata => strata%next
   200)     enddo
   201)     cur_patch => cur_patch%next
   202)   enddo
   203)   
   204)   ! ensure that ghosted values for material ids are up to date
   205)   call RealLocalToLocalWithArray(realization,MATERIAL_ID_ARRAY)
   206)   
   207)   ! set material ids in material auxvar.  this must come after the update of 
   208)   ! ghost values.
   209)   cur_patch => realization%patch_list%first
   210)   do
   211)     if (.not.associated(cur_patch)) exit
   212)     ! set material ids to uninitialized
   213)     material_auxvars => cur_patch%aux%Material%auxvars
   214)     grid => cur_patch%grid
   215)     do ghosted_id = 1, grid%ngmax
   216)       material_auxvars(ghosted_id)%id = cur_patch%imat(ghosted_id)
   217)     enddo
   218)     cur_patch => cur_patch%next
   219)   enddo
   220) 
   221) end subroutine InitSubsurfAssignMatIDsToRegns
   222) 
   223) ! ************************************************************************** !
   224) 
   225) subroutine InitSubsurfAssignMatProperties(realization)
   226)   ! 
   227)   ! Assigns material properties based on material ids
   228)   ! 
   229)   ! Author: Glenn Hammond
   230)   ! Date: 10/07/14
   231)   ! 
   232)   use Realization_Subsurface_class
   233)   use Grid_module
   234)   use Discretization_module
   235)   use Field_module
   236)   use Patch_module
   237)   use Material_Aux_class
   238)   use Material_module
   239)   use Option_module
   240)   use Creep_Closure_module
   241)   use Fracture_module
   242)   use Variables_module, only : PERMEABILITY_X, PERMEABILITY_Y, &
   243)                                PERMEABILITY_Z, PERMEABILITY_XY, &
   244)                                PERMEABILITY_YZ, PERMEABILITY_XZ, &
   245)                                TORTUOSITY, POROSITY, SOIL_COMPRESSIBILITY
   246)   use HDF5_module
   247)   
   248)   implicit none
   249)   
   250) #include "petsc/finclude/petscvec.h"
   251) #include "petsc/finclude/petscvec.h90"
   252) 
   253)   class(realization_subsurface_type) :: realization
   254)   
   255)   PetscReal, pointer :: icap_loc_p(:)
   256)   PetscReal, pointer :: ithrm_loc_p(:)
   257)   PetscReal, pointer :: por0_p(:)
   258)   PetscReal, pointer :: tor0_p(:)
   259)   PetscReal, pointer :: perm_xx_p(:)
   260)   PetscReal, pointer :: perm_yy_p(:)
   261)   PetscReal, pointer :: perm_zz_p(:)
   262)   PetscReal, pointer :: perm_xz_p(:)
   263)   PetscReal, pointer :: perm_xy_p(:)
   264)   PetscReal, pointer :: perm_yz_p(:)
   265)   PetscReal, pointer :: perm_pow_p(:)
   266)   PetscReal, pointer :: vec_p(:)
   267)   PetscReal, pointer :: compress_p(:)
   268)   
   269)   character(len=MAXSTRINGLENGTH) :: string, string2
   270)   type(material_property_type), pointer :: material_property
   271)   type(material_property_type), pointer :: null_material_property
   272)   type(option_type), pointer :: option
   273)   type(discretization_type), pointer :: discretization
   274)   type(grid_type), pointer :: grid
   275)   type(field_type), pointer :: field
   276)   type(patch_type), pointer :: patch
   277)   type(material_type), pointer :: Material
   278)   PetscInt :: local_id, ghosted_id, material_id
   279)   PetscInt :: i
   280)   PetscInt :: tempint
   281)   PetscReal :: tempreal
   282)   PetscErrorCode :: ierr
   283)   
   284)   option => realization%option
   285)   discretization => realization%discretization
   286)   field => realization%field
   287)   patch => realization%patch
   288)   grid => patch%grid
   289)   
   290)   ! set cell by cell material properties
   291)   ! create null material property for inactive cells
   292)   null_material_property => MaterialPropertyCreate()
   293)   if (option%nflowdof > 0) then
   294)     call VecGetArrayF90(field%icap_loc,icap_loc_p,ierr);CHKERRQ(ierr)
   295)     call VecGetArrayF90(field%ithrm_loc,ithrm_loc_p,ierr);CHKERRQ(ierr)
   296)     call VecGetArrayF90(field%perm0_xx,perm_xx_p,ierr);CHKERRQ(ierr)
   297)     call VecGetArrayF90(field%perm0_yy,perm_yy_p,ierr);CHKERRQ(ierr)
   298)     call VecGetArrayF90(field%perm0_zz,perm_zz_p,ierr);CHKERRQ(ierr)
   299)     if (soil_compressibility_index > 0) then
   300)       call VecGetArrayF90(field%compressibility0,compress_p,ierr);CHKERRQ(ierr)
   301)     endif
   302)   endif
   303)   call VecGetArrayF90(field%porosity0,por0_p,ierr);CHKERRQ(ierr)
   304)   call VecGetArrayF90(field%tortuosity0,tor0_p,ierr);CHKERRQ(ierr)
   305)         
   306)   ! have to use Material%auxvars() and not material_auxvars() due to memory
   307)   ! errors in gfortran
   308)   Material => patch%aux%Material
   309) 
   310)   !if material is associated with fracture, then allocate memory.  
   311)   do ghosted_id = 1, grid%ngmax
   312)     material_id = patch%imat(ghosted_id)
   313)     if (material_id > 0) then
   314)       material_property => &
   315)         patch%material_property_array(material_id)%ptr
   316)       call FractureAuxvarInit(material_property%fracture, &
   317)         patch%aux%Material%auxvars(ghosted_id))
   318)     endif
   319)   enddo
   320)   
   321)   do local_id = 1, grid%nlmax
   322)     ghosted_id = grid%nL2G(local_id)
   323)     material_id = patch%imat(ghosted_id)
   324)     if (material_id == 0) then
   325)       material_property => null_material_property
   326)     else if (iabs(material_id) <= &
   327)              size(patch%material_property_array)) then
   328)       if (material_id < 0) then
   329)         material_property => null_material_property
   330)       else
   331)         material_property => &
   332)           patch%material_property_array(material_id)%ptr
   333)         if (.not.associated(material_property)) then
   334)           write(string,*) &
   335)             patch%imat_internal_to_external(material_id)
   336)           option%io_buffer = 'No material property for material id ' // &
   337)                               trim(adjustl(string)) &
   338)                               //  ' defined in input file.'
   339)           call printErrMsgByRank(option)
   340)         endif
   341)       endif
   342)     else if (Uninitialized(material_id)) then 
   343)       write(string,*) grid%nG2A(ghosted_id)
   344)       option%io_buffer = 'Uninitialized material id in patch at cell ' // &
   345)                           trim(adjustl(string))
   346)       call printErrMsgByRank(option)
   347)     else if (material_id > size(patch%material_property_array)) then
   348)       write(option%io_buffer,*) patch%imat_internal_to_external(material_id)
   349)       option%io_buffer = 'Unmatched material id in patch: ' // &
   350)         adjustl(trim(option%io_buffer))
   351)       call printErrMsgByRank(option)
   352)     else
   353)       option%io_buffer = 'Something messed up with material ids. Possibly &
   354)         &material ids not assigned to all grid cells. Contact Glenn!'
   355)       call printErrMsgByRank(option)
   356)     endif
   357)     if (option%nflowdof > 0) then
   358)       patch%sat_func_id(ghosted_id) = &
   359)         material_property%saturation_function_id
   360)       icap_loc_p(ghosted_id) = material_property%saturation_function_id
   361)       ithrm_loc_p(ghosted_id) = iabs(material_property%internal_id)
   362)       perm_xx_p(local_id) = material_property%permeability(1,1)
   363)       perm_yy_p(local_id) = material_property%permeability(2,2)
   364)       perm_zz_p(local_id) = material_property%permeability(3,3)
   365)       if (soil_compressibility_index > 0) then
   366)         compress_p(local_id) = material_property%soil_compressibility
   367)       endif
   368)     endif
   369)     if (associated(Material%auxvars)) then
   370)       call MaterialAssignPropertyToAux(Material%auxvars(ghosted_id), &
   371)                                         material_property,option)
   372)     endif
   373)     por0_p(local_id) = material_property%porosity
   374)     tor0_p(local_id) = material_property%tortuosity
   375)   enddo
   376)   call MaterialPropertyDestroy(null_material_property)
   377) 
   378)   if (option%nflowdof > 0) then
   379)     call VecRestoreArrayF90(field%icap_loc,icap_loc_p,ierr);CHKERRQ(ierr)
   380)     call VecRestoreArrayF90(field%ithrm_loc,ithrm_loc_p,ierr);CHKERRQ(ierr)
   381)     call VecRestoreArrayF90(field%perm0_xx,perm_xx_p,ierr);CHKERRQ(ierr)
   382)     call VecRestoreArrayF90(field%perm0_yy,perm_yy_p,ierr);CHKERRQ(ierr)
   383)     call VecRestoreArrayF90(field%perm0_zz,perm_zz_p,ierr);CHKERRQ(ierr)
   384)     if (soil_compressibility_index > 0) then
   385)       call VecRestoreArrayF90(field%compressibility0,compress_p, &
   386)                               ierr);CHKERRQ(ierr)
   387)     endif
   388)   endif
   389)   call VecRestoreArrayF90(field%porosity0,por0_p,ierr);CHKERRQ(ierr)
   390)   call VecRestoreArrayF90(field%tortuosity0,tor0_p,ierr);CHKERRQ(ierr)
   391)         
   392)   ! read in any user-defined property fields
   393)   do material_id = 1, size(patch%material_property_array)
   394)     material_property => &
   395)             patch%material_property_array(material_id)%ptr
   396)     if (.not.associated(material_property)) cycle
   397)     if (material_property%active) then
   398)       if (associated(material_property%permeability_dataset)) then
   399)         call SubsurfReadPermsFromFile(realization,material_property)
   400)       endif
   401)       if (associated(material_property%compressibility_dataset)) then
   402)         call SubsurfReadDatasetToVecWithMask(realization, &
   403)                material_property%compressibility_dataset, &
   404)                material_property%internal_id,PETSC_FALSE,field%compressibility0)
   405)       endif
   406)       if (associated(material_property%porosity_dataset)) then
   407)         call SubsurfReadDatasetToVecWithMask(realization, &
   408)                material_property%porosity_dataset, &
   409)                material_property%internal_id,PETSC_FALSE,field%porosity0)
   410)       endif
   411)       if (associated(material_property%tortuosity_dataset)) then
   412)         call SubsurfReadDatasetToVecWithMask(realization, &
   413)                material_property%tortuosity_dataset, &
   414)                material_property%internal_id,PETSC_FALSE,field%tortuosity0)
   415)       endif
   416)     endif
   417)   enddo
   418)       
   419)   ! update ghosted values
   420)   if (option%nflowdof > 0) then
   421)     call DiscretizationGlobalToLocal(discretization,field%perm0_xx, &
   422)                                      field%work_loc,ONEDOF)
   423)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   424)                                  PERMEABILITY_X,ZERO_INTEGER)
   425)     call DiscretizationGlobalToLocal(discretization,field%perm0_yy, &
   426)                                      field%work_loc,ONEDOF)  
   427)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   428)                                  PERMEABILITY_Y,ZERO_INTEGER)
   429)     call DiscretizationGlobalToLocal(discretization,field%perm0_zz, &
   430)                                      field%work_loc,ONEDOF)   
   431)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   432)                                  PERMEABILITY_Z,ZERO_INTEGER)
   433)     call DiscretizationLocalToLocal(discretization,field%icap_loc, &
   434)                                     field%icap_loc,ONEDOF)   
   435)     call DiscretizationLocalToLocal(discretization,field%ithrm_loc, &
   436)                                     field%ithrm_loc,ONEDOF)
   437)     call RealLocalToLocalWithArray(realization,SATURATION_FUNCTION_ID_ARRAY)
   438)     
   439)     if (soil_compressibility_index > 0) then
   440)       call DiscretizationGlobalToLocal(discretization,field%compressibility0, &
   441)                                        field%work_loc,ONEDOF)
   442)       call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   443)                                    SOIL_COMPRESSIBILITY,ZERO_INTEGER)
   444)     endif
   445)   endif
   446)   
   447)   call DiscretizationGlobalToLocal(discretization,field%porosity0, &
   448)                                    field%work_loc,ONEDOF)
   449)   call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   450)                                POROSITY,POROSITY_MINERAL)
   451)   call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   452)                                POROSITY,POROSITY_CURRENT)
   453)   call DiscretizationGlobalToLocal(discretization,field%tortuosity0, &
   454)                                     field%work_loc,ONEDOF)
   455)   call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
   456)                                TORTUOSITY,ZERO_INTEGER)
   457) 
   458)   ! copy rock properties to neighboring ghost cells
   459)   do i = 1, max_material_index
   460)     call VecGetArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
   461)     do local_id = 1, patch%grid%nlmax
   462)       vec_p(local_id) = &
   463)           Material%auxvars(patch%grid%nL2G(local_id))%soil_properties(i)
   464)     enddo
   465)     call VecRestoreArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
   466)     call DiscretizationGlobalToLocal(discretization,field%work, &
   467)                                      field%work_loc,ONEDOF)
   468)     call VecGetArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   469)     do ghosted_id = 1, patch%grid%ngmax
   470)       Material%auxvars(ghosted_id)%soil_properties(i) = &
   471)          vec_p(ghosted_id)
   472)     enddo
   473)     call VecRestoreArrayF90(field%work_loc,vec_p,ierr);CHKERRQ(ierr)
   474)   enddo
   475)   
   476)   if (associated(creep_closure)) then
   477)     material_property => &
   478)       MaterialPropGetPtrFromArray(creep_closure%material_name, &
   479)                                   patch%material_property_array)
   480)     if (.not.associated(material_property)) then
   481)       option%io_buffer = 'Creep material "' // &
   482)                         trim(creep_closure%material_name) // &
   483)                         '" not found in material list'
   484)       call printErrMsg(option)
   485)     endif
   486)     creep_closure%imat = material_property%internal_id
   487)   endif
   488)   
   489) end subroutine InitSubsurfAssignMatProperties
   490) 
   491) ! ************************************************************************** !
   492) 
   493) subroutine SubsurfReadMaterialIDsFromFile(realization,realization_dependent, &
   494)                                           filename)
   495)   ! 
   496)   ! Reads in grid cell materials
   497)   ! 
   498)   ! Author: Glenn Hammond
   499)   ! Date: 1/03/08
   500)   ! 
   501) 
   502)   use Realization_Subsurface_class
   503)   use Field_module
   504)   use Grid_module
   505)   use Option_module
   506)   use Patch_module
   507)   use Discretization_module
   508)   use Logging_module
   509)   use Input_Aux_module
   510)   use Material_module
   511) 
   512)   use HDF5_module
   513)   
   514)   implicit none
   515) 
   516) #include "petsc/finclude/petscvec.h"
   517) #include "petsc/finclude/petscvec.h90"
   518)   
   519)   class(realization_subsurface_type) :: realization
   520)   PetscBool :: realization_dependent
   521)   character(len=MAXSTRINGLENGTH) :: filename
   522)   
   523)   type(field_type), pointer :: field
   524)   type(grid_type), pointer :: grid
   525)   type(option_type), pointer :: option
   526)   type(patch_type), pointer :: patch   
   527)   type(input_type), pointer :: input
   528)   type(discretization_type), pointer :: discretization
   529)   character(len=MAXSTRINGLENGTH) :: group_name
   530)   character(len=MAXSTRINGLENGTH) :: dataset_name
   531)   PetscBool :: append_realization_id
   532)   PetscInt :: ghosted_id, natural_id, material_id
   533)   PetscInt :: fid = 86
   534)   PetscInt :: status
   535)   PetscInt, pointer :: external_to_internal_mapping(:)
   536)   Vec :: global_vec
   537)   Vec :: local_vec
   538)   PetscErrorCode :: ierr
   539) 
   540)   field => realization%field
   541)   patch => realization%patch
   542)   grid => patch%grid
   543)   option => realization%option
   544)   discretization => realization%discretization
   545)   
   546)   if (index(filename,'.h5') > 0) then
   547)     group_name = 'Materials'
   548)     dataset_name = 'Material Ids'
   549)     call DiscretizationCreateVector(discretization,ONEDOF,global_vec,GLOBAL, &
   550)                                     option)
   551)     call DiscretizationCreateVector(discretization,ONEDOF,local_vec,LOCAL, &
   552)                                     option)
   553)     call HDF5ReadCellIndexedIntegerArray(realization,global_vec, &
   554)                                          filename,group_name, &
   555)                                          dataset_name,realization_dependent)
   556)     call DiscretizationGlobalToLocal(discretization,global_vec,local_vec,ONEDOF)
   557)     call GridCopyVecToIntegerArray(grid,patch%imat,local_vec,grid%ngmax)
   558)     call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
   559)     call VecDestroy(local_vec,ierr);CHKERRQ(ierr)
   560)   else
   561)     call PetscLogEventBegin(logging%event_hash_map,ierr);CHKERRQ(ierr)
   562)     call GridCreateNaturalToGhostedHash(grid,option)
   563)     input => InputCreate(IUNIT_TEMP,filename,option)
   564)     do
   565)       call InputReadPflotranString(input,option)
   566)       if (InputError(input)) exit
   567)       call InputReadInt(input,option,natural_id)
   568)       call InputErrorMsg(input,option,'natural id','STRATA')
   569)       ! natural ids in hash are zero-based
   570)       ghosted_id = GridGetLocalGhostedIdFromHash(grid,natural_id)
   571)       if (ghosted_id > 0) then
   572)         call InputReadInt(input,option,material_id)
   573)         call InputErrorMsg(input,option,'material id','STRATA')
   574)         patch%imat(ghosted_id) = material_id
   575)       endif
   576)     enddo
   577)     call InputDestroy(input)
   578)     call GridDestroyHashTable(grid)
   579)     call PetscLogEventEnd(logging%event_hash_map,ierr);CHKERRQ(ierr)
   580)   endif
   581)   
   582)   call MaterialCreateExtToIntMapping(patch%material_property_array, &
   583)                                      external_to_internal_mapping)
   584)   call MaterialApplyMapping(external_to_internal_mapping,patch%imat)
   585)   deallocate(external_to_internal_mapping)
   586)   nullify(external_to_internal_mapping)
   587)   
   588) end subroutine SubsurfReadMaterialIDsFromFile
   589) 
   590) ! ************************************************************************** !
   591) 
   592) subroutine SubsurfReadPermsFromFile(realization,material_property)
   593)   ! 
   594)   ! Reads in grid cell permeabilities
   595)   ! 
   596)   ! Author: Glenn Hammond
   597)   ! Date: 01/19/09
   598)   ! 
   599) 
   600)   use Realization_Subsurface_class
   601)   use Field_module
   602)   use Grid_module
   603)   use Option_module
   604)   use Patch_module
   605)   use Discretization_module
   606)   use Logging_module
   607)   use Input_Aux_module
   608)   use Material_module
   609)   use HDF5_module
   610)   use Dataset_Common_HDF5_class
   611)   
   612)   implicit none
   613)   
   614) #include "petsc/finclude/petscvec.h"
   615) #include "petsc/finclude/petscvec.h90"
   616) 
   617)   class(realization_subsurface_type) :: realization
   618)   type(material_property_type) :: material_property
   619) 
   620)   type(field_type), pointer :: field
   621)   type(patch_type), pointer :: patch
   622)   type(grid_type), pointer :: grid
   623)   type(option_type), pointer :: option
   624)   type(input_type), pointer :: input
   625)   type(discretization_type), pointer :: discretization
   626)   character(len=MAXWORDLENGTH) :: word
   627)   class(dataset_common_hdf5_type), pointer :: dataset_common_hdf5_ptr
   628)   PetscInt :: local_id
   629)   PetscInt :: idirection, temp_int
   630)   PetscReal :: ratio, scale
   631)   Vec :: global_vec
   632)   PetscErrorCode :: ierr
   633)   
   634)   PetscReal, pointer :: vec_p(:)
   635)   PetscReal, pointer :: perm_xx_p(:)
   636)   PetscReal, pointer :: perm_yy_p(:)
   637)   PetscReal, pointer :: perm_zz_p(:)
   638) 
   639)   field => realization%field
   640)   patch => realization%patch
   641)   grid => patch%grid
   642)   option => realization%option
   643)   discretization => realization%discretization
   644) 
   645)   call VecGetArrayF90(field%perm0_xx,perm_xx_p,ierr);CHKERRQ(ierr)
   646)   call VecGetArrayF90(field%perm0_yy,perm_yy_p,ierr);CHKERRQ(ierr)
   647)   call VecGetArrayF90(field%perm0_zz,perm_zz_p,ierr);CHKERRQ(ierr)
   648)   
   649)   call DiscretizationCreateVector(discretization,ONEDOF,global_vec,GLOBAL, &
   650)                                   option)
   651)   if (material_property%isotropic_permeability .or. &
   652)       (.not.material_property%isotropic_permeability .and. &
   653)        Initialized(material_property%vertical_anisotropy_ratio))) then
   654)     ! Although the mask of material ID is applied below, we must only read
   655)     ! in the permeabilities that apply to this material so that small, 
   656)     ! localized gridded datasets (that only apply to a subset of the domain)
   657)     ! can be used.
   658)     call VecZeroEntries(global_vec,ierr);CHKERRQ(ierr)
   659)     call SubsurfReadDatasetToVecWithMask(realization, &
   660)                                     material_property%permeability_dataset, &
   661)                                     material_property%internal_id, &
   662)                                     PETSC_FALSE,global_vec)
   663)     call VecGetArrayF90(global_vec,vec_p,ierr);CHKERRQ(ierr)
   664)     ratio = 1.d0
   665)     scale = 1.d0
   666)     !TODO(geh): fix so that ratio and scale work for perms outside
   667)     ! of dataset
   668)     if (Initialized(material_property%vertical_anisotropy_ratio)) then
   669)       ratio = material_property%vertical_anisotropy_ratio
   670)     endif
   671)     if (material_property%permeability_scaling_factor > 0.d0) then
   672)       scale = material_property%permeability_scaling_factor
   673)     endif
   674)     do local_id = 1, grid%nlmax
   675)       if (patch%imat(grid%nL2G(local_id)) == &
   676)           material_property%internal_id) then
   677)         perm_xx_p(local_id) = vec_p(local_id)*scale
   678)         perm_yy_p(local_id) = vec_p(local_id)*scale
   679)         perm_zz_p(local_id) = vec_p(local_id)*ratio*scale
   680)       endif
   681)     enddo
   682)     call VecRestoreArrayF90(global_vec,vec_p,ierr);CHKERRQ(ierr)
   683)   else
   684)     temp_int = Z_DIRECTION
   685)     do idirection = X_DIRECTION,temp_int
   686)       select case(idirection)
   687)         case(X_DIRECTION)
   688)           dataset_common_hdf5_ptr => &
   689)              DatasetCommonHDF5Cast(material_property%permeability_dataset)
   690)         case(Y_DIRECTION)
   691)           dataset_common_hdf5_ptr => &
   692)              DatasetCommonHDF5Cast(material_property%permeability_dataset_y)
   693)         case(Z_DIRECTION)
   694)           dataset_common_hdf5_ptr => &
   695)              DatasetCommonHDF5Cast(material_property%permeability_dataset_z)
   696)       end select
   697)       ! Although the mask of material ID is applied below, we must only read
   698)       ! in the permeabilities that apply to this material so that small, 
   699)       ! localized gridded datasets (that only apply to a subset of the domain)
   700)       ! can be used.
   701)       call VecZeroEntries(global_vec,ierr);CHKERRQ(ierr)
   702)       call SubsurfReadDatasetToVecWithMask(realization, &
   703)                                            dataset_common_hdf5_ptr,&
   704)                                            material_property%internal_id, &
   705)                                            PETSC_FALSE,global_vec)
   706)       call VecGetArrayF90(global_vec,vec_p,ierr);CHKERRQ(ierr)
   707)       select case(idirection)
   708)         case(X_DIRECTION)
   709)           do local_id = 1, grid%nlmax
   710)             if (patch%imat(grid%nL2G(local_id)) == &
   711)                 material_property%internal_id) then
   712)               perm_xx_p(local_id) = vec_p(local_id)
   713)             endif
   714)           enddo
   715)         case(Y_DIRECTION)
   716)           do local_id = 1, grid%nlmax
   717)             if (patch%imat(grid%nL2G(local_id)) == &
   718)                 material_property%internal_id) then
   719)               perm_yy_p(local_id) = vec_p(local_id)
   720)             endif
   721)           enddo
   722)         case(Z_DIRECTION)
   723)           do local_id = 1, grid%nlmax
   724)             if (patch%imat(grid%nL2G(local_id)) == &
   725)                 material_property%internal_id) then
   726)               perm_zz_p(local_id) = vec_p(local_id)
   727)             endif
   728)           enddo
   729)       end select
   730)       call VecRestoreArrayF90(global_vec,vec_p,ierr);CHKERRQ(ierr)
   731)     enddo
   732)   endif
   733)   call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
   734)   
   735)   call VecRestoreArrayF90(field%perm0_xx,perm_xx_p,ierr);CHKERRQ(ierr)
   736)   call VecRestoreArrayF90(field%perm0_yy,perm_yy_p,ierr);CHKERRQ(ierr)
   737)   call VecRestoreArrayF90(field%perm0_zz,perm_zz_p,ierr);CHKERRQ(ierr)
   738)   
   739) end subroutine SubsurfReadPermsFromFile
   740) 
   741) ! ************************************************************************** !
   742) 
   743) subroutine SubsurfReadDatasetToVecWithMask(realization,dataset, &
   744)                                            material_id,read_all_values,vec)
   745)   ! 
   746)   ! Reads a dataset into a PETSc Vec using the material id as a mask
   747)   ! 
   748)   ! Author: Glenn Hammond
   749)   ! Date: 01/19/2016
   750)   ! 
   751) 
   752)   use Realization_Subsurface_class
   753)   use Field_module
   754)   use Grid_module
   755)   use Option_module
   756)   use Patch_module
   757)   use Logging_module
   758)   use Input_Aux_module
   759)   use Material_module
   760)   use HDF5_module
   761)   use Dataset_Base_class
   762)   use Dataset_Common_HDF5_class
   763)   use Dataset_Gridded_HDF5_class
   764)   
   765)   implicit none
   766)   
   767) #include "petsc/finclude/petscvec.h"
   768) #include "petsc/finclude/petscvec.h90"
   769) 
   770)   class(realization_subsurface_type) :: realization
   771)   class(dataset_base_type) :: dataset
   772)   PetscInt :: material_id
   773)   PetscBool :: read_all_values
   774)   Vec :: vec
   775) 
   776)   type(field_type), pointer :: field
   777)   type(patch_type), pointer :: patch
   778)   type(grid_type), pointer :: grid
   779)   type(option_type), pointer :: option
   780)   type(input_type), pointer :: input
   781)   character(len=MAXSTRINGLENGTH) :: group_name
   782)   character(len=MAXSTRINGLENGTH) :: dataset_name
   783)   character(len=MAXSTRINGLENGTH) :: filename
   784)   PetscInt :: local_id, ghosted_id, natural_id
   785)   PetscReal :: tempreal
   786)   PetscErrorCode :: ierr
   787)   
   788)   PetscReal, pointer :: vec_p(:)
   789)   PetscReal, pointer :: work_p(:)
   790)   
   791)   field => realization%field
   792)   patch => realization%patch
   793)   grid => patch%grid
   794)   option => realization%option
   795) 
   796)   call VecGetArrayF90(vec,vec_p,ierr);CHKERRQ(ierr)
   797) 
   798)   if (index(dataset%filename,'.h5') > 0) then
   799)     group_name = ''
   800)     dataset_name = dataset%name
   801)     select type(dataset)
   802)       class is(dataset_gridded_hdf5_type)
   803)         call DatasetGriddedHDF5Load(dataset,option)
   804)         do local_id = 1, grid%nlmax
   805)           ghosted_id = grid%nL2G(local_id)
   806)           if (read_all_values .or. &
   807)               patch%imat(ghosted_id) == material_id) then
   808)             call DatasetGriddedHDF5InterpolateReal(dataset, &
   809)                    grid%x(ghosted_id),grid%y(ghosted_id),grid%z(ghosted_id), &
   810)                    vec_p(local_id),option)
   811)           endif
   812)         enddo
   813)         ! now we strip the dataset to save storage, saving only the name
   814)         ! and filename incase it must be read later
   815)         filename = dataset%filename
   816)         dataset_name = dataset%name
   817)         call DatasetGriddedHDF5Strip(dataset)
   818)         call DatasetGriddedHDF5Init(dataset)
   819)         dataset%filename = filename
   820)         dataset%name = trim(dataset_name)
   821)       class is(dataset_common_hdf5_type)
   822)         dataset_name = dataset%hdf5_dataset_name
   823)         call HDF5ReadCellIndexedRealArray(realization,field%work, &
   824)                                           dataset%filename, &
   825)                                           group_name,dataset_name, &
   826)                                           dataset%realization_dependent)
   827)         call VecGetArrayF90(field%work,work_p,ierr);CHKERRQ(ierr)
   828)         if (read_all_values) then
   829)           vec_p(:) = work_p(:)
   830)         else
   831)           do local_id = 1, grid%nlmax
   832)             if (patch%imat(grid%nL2G(local_id)) == material_id) then
   833)               vec_p(local_id) = work_p(local_id)
   834)             endif
   835)           enddo
   836)         endif
   837)         call VecRestoreArrayF90(field%work,work_p,ierr);CHKERRQ(ierr)
   838)       class default
   839)         option%io_buffer = 'Dataset "' // trim(dataset%name) // '" is of the &
   840)           &wrong type for SubsurfReadDatasetToVecWithMask()'
   841)         call printErrMsg(option)
   842)     end select
   843)   else
   844)     call PetscLogEventBegin(logging%event_hash_map,ierr);CHKERRQ(ierr)
   845)     call GridCreateNaturalToGhostedHash(grid,option)
   846)     input => InputCreate(IUNIT_TEMP,dataset%filename,option)
   847)     do
   848)       call InputReadPflotranString(input,option)
   849)       if (InputError(input)) exit
   850)       call InputReadInt(input,option,natural_id)
   851)       call InputErrorMsg(input,option,'ASCII natural id', &
   852)                          'SubsurfReadDatasetToVecWithMask')
   853)       ghosted_id = GridGetLocalGhostedIdFromHash(grid,natural_id)
   854)       if (ghosted_id > 0) then
   855)         if (read_all_values .or. &
   856)             patch%imat(ghosted_id) == material_id) then
   857)           local_id = grid%nG2L(ghosted_id)
   858)           if (local_id > 0) then
   859)             call InputReadDouble(input,option,tempreal)
   860)             call InputErrorMsg(input,option,'dataset value', &
   861)                                'SubsurfReadDatasetToVecWithMask')
   862)             vec_p(local_id) = tempreal
   863)           endif
   864)         endif
   865)       endif
   866)     enddo
   867)     call InputDestroy(input)
   868)     call GridDestroyHashTable(grid)
   869)     call PetscLogEventEnd(logging%event_hash_map,ierr);CHKERRQ(ierr)
   870)   endif
   871)   
   872)   call VecRestoreArrayF90(vec,vec_p,ierr);CHKERRQ(ierr)
   873)   
   874) end subroutine SubsurfReadDatasetToVecWithMask
   875) 
   876) ! ************************************************************************** !
   877) 
   878) subroutine SubsurfAssignVolsToMatAuxVars(realization)
   879)   ! 
   880)   ! Assigns the cell volumes currently stored in field%volume0 to the 
   881)   ! material auxiliary variable object
   882)   ! 
   883)   ! Author: Glenn Hammond
   884)   ! Date: 01/13/14, 12/04/14
   885)   ! 
   886) 
   887)   use Realization_Subsurface_class
   888)   use Option_module
   889)   use Material_module
   890)   use Discretization_module
   891)   use Field_module
   892)   use Variables_module, only : VOLUME
   893)   
   894)   implicit none
   895)   
   896)   class(realization_subsurface_type) :: realization
   897)   
   898)   type(option_type), pointer :: option
   899)   type(field_type), pointer :: field
   900) 
   901)   option => realization%option
   902)   field => realization%field
   903) 
   904)   call DiscretizationGlobalToLocal(realization%discretization,field%volume0, &
   905)                                    field%work_loc,ONEDOF)
   906)   call MaterialSetAuxVarVecLoc(realization%patch%aux%Material, &
   907)                                field%work_loc,VOLUME,ZERO_INTEGER)
   908) 
   909) end subroutine SubsurfAssignVolsToMatAuxVars
   910) 
   911) ! ************************************************************************** !
   912) 
   913) subroutine SubsurfSandboxesSetup(realization)
   914)   ! 
   915)   ! Initializes sandbox objects.
   916)   ! 
   917)   ! Author: Glenn Hammond
   918)   ! Date: 05/06/14, 12/04/14
   919) 
   920)   use Realization_Subsurface_class
   921)   use SrcSink_Sandbox_module
   922)   
   923)   class(realization_subsurface_type) :: realization
   924)   
   925)   call SSSandboxSetup(realization%patch%grid,realization%option, &
   926)                       realization%output_option)
   927)   
   928) end subroutine SubsurfSandboxesSetup
   929) 
   930) ! ************************************************************************** !
   931) 
   932) subroutine InitSubsurfaceSetupZeroArrays(realization)
   933)   ! 
   934)   ! Initializes sandbox objects.
   935)   ! 
   936)   ! Author: Glenn Hammond
   937)   ! Date: 03/11/16
   938) 
   939)   use Realization_Subsurface_class
   940)   use Option_module
   941) 
   942)   class(realization_subsurface_type) :: realization
   943)   
   944)   type(option_type), pointer :: option
   945)   PetscBool, allocatable :: dof_is_active(:)
   946)   PetscInt :: ndof
   947)   
   948)   option => realization%option
   949)   
   950)   if (option%nflowdof > 0) then
   951)     allocate(dof_is_active(option%nflowdof))
   952)     dof_is_active = PETSC_TRUE
   953) #if defined(ISOTHERMAL)
   954)     select case(option%iflowmode)
   955)       case(TH_MODE)
   956)         ! second equation is energy
   957)         dof_is_active(TWO_INTEGER) = PETSC_FALSE
   958)       case(MPH_MODE,IMS_MODE,MIS_MODE,FLASH2_MODE)
   959)         ! third equation is energy
   960)         dof_is_active(THREE_INTEGER) = PETSC_FALSE
   961)     end select
   962) #endif
   963)     select case(option%iflowmode)
   964)       case(RICHARDS_MODE)
   965)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
   966)                       realization%patch%aux%Richards%zero_rows_local, &
   967)                       realization%patch%aux%Richards%zero_rows_local_ghosted, &
   968)                       realization%patch%aux%Richards%n_zero_rows, &
   969)                       realization%patch%aux%Richards%inactive_cells_exist, &
   970)                       option)
   971)       case(TH_MODE)
   972)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
   973)                       realization%patch%aux%TH%zero_rows_local, &
   974)                       realization%patch%aux%TH%zero_rows_local_ghosted, &
   975)                       realization%patch%aux%TH%n_zero_rows, &
   976)                       realization%patch%aux%TH%inactive_cells_exist, &
   977)                       option)
   978)       case(MPH_MODE)
   979)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
   980)                       realization%patch%aux%Mphase%zero_rows_local, &
   981)                       realization%patch%aux%Mphase%zero_rows_local_ghosted, &
   982)                       realization%patch%aux%Mphase%n_zero_rows, &
   983)                       realization%patch%aux%Mphase%inactive_cells_exist, &
   984)                       option)
   985)       case(G_MODE)
   986)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
   987)                       realization%patch%aux%General%inactive_rows_local, &
   988)                     realization%patch%aux%General%inactive_rows_local_ghosted, &
   989)                       realization%patch%aux%General%n_inactive_rows, &
   990)                       realization%patch%aux%General%inactive_cells_exist, &
   991)                       option)
   992)       case(TOIL_IMS_MODE)
   993)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
   994)                       realization%patch%aux%TOil_ims%inactive_rows_local, &
   995)                    realization%patch%aux%TOil_ims%inactive_rows_local_ghosted, &
   996)                       realization%patch%aux%TOil_ims%n_inactive_rows, &
   997)                       realization%patch%aux%TOil_ims%inactive_cells_exist, &
   998)                       option)
   999)       case(IMS_MODE)
  1000)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
  1001)                       realization%patch%aux%Immis%zero_rows_local, &
  1002)                       realization%patch%aux%Immis%zero_rows_local_ghosted, &
  1003)                       realization%patch%aux%Immis%n_zero_rows, &
  1004)                       realization%patch%aux%Immis%inactive_cells_exist, &
  1005)                       option)
  1006)       case(MIS_MODE)
  1007)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
  1008)                       realization%patch%aux%Miscible%zero_rows_local, &
  1009)                       realization%patch%aux%Miscible%zero_rows_local_ghosted, &
  1010)                       realization%patch%aux%Miscible%n_zero_rows, &
  1011)                       realization%patch%aux%Miscible%inactive_cells_exist, &
  1012)                       option)
  1013)       case(FLASH2_MODE)
  1014)         call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
  1015)                       realization%patch%aux%Flash2%zero_rows_local, &
  1016)                       realization%patch%aux%Flash2%zero_rows_local_ghosted, &
  1017)                       realization%patch%aux%Flash2%n_zero_rows, &
  1018)                       realization%patch%aux%Flash2%inactive_cells_exist, &
  1019)                       option)
  1020)     end select
  1021)     deallocate(dof_is_active)
  1022)   endif
  1023) 
  1024)   if (option%ntrandof > 0) then
  1025)     ! remove ndof above if this is moved
  1026)     if (option%transport%reactive_transport_coupling == GLOBAL_IMPLICIT) then
  1027)       ndof = realization%reaction%ncomp
  1028)     else
  1029)       ndof = 1
  1030)     endif
  1031)     allocate(dof_is_active(ndof))
  1032)     dof_is_active = PETSC_TRUE  
  1033)     call InitSubsurfaceCreateZeroArray(realization%patch,dof_is_active, &
  1034)                   realization%patch%aux%RT%zero_rows_local, &
  1035)                   realization%patch%aux%RT%zero_rows_local_ghosted, &
  1036)                   realization%patch%aux%RT%n_zero_rows, &
  1037)                   realization%patch%aux%RT%inactive_cells_exist, &
  1038)                   option)
  1039)     deallocate(dof_is_active)
  1040)   endif  
  1041) 
  1042) end subroutine InitSubsurfaceSetupZeroArrays
  1043) 
  1044) ! ************************************************************************** !
  1045) 
  1046) subroutine InitSubsurfaceCreateZeroArray(patch,dof_is_active, &
  1047)                                          inactive_rows_local, &
  1048)                                          inactive_rows_local_ghosted, &
  1049)                                          n_inactive_rows, &
  1050)                                          inactive_cells_exist, &
  1051)                                          option)
  1052)   ! 
  1053)   ! Computes the zeroed rows for inactive grid cells
  1054)   ! 
  1055)   ! Author: Glenn Hammond
  1056)   ! Date: 12/13/07, 03/02/16
  1057)   ! 
  1058)   use Realization_Subsurface_class
  1059)   use Patch_module
  1060)   use Grid_module
  1061)   use Option_module
  1062)   use Field_module
  1063)   use Utility_module, only : DeallocateArray
  1064)   
  1065)   implicit none
  1066) 
  1067)   type(patch_type) :: patch
  1068)   PetscBool :: dof_is_active(:)
  1069)   PetscInt, pointer :: inactive_rows_local(:)
  1070)   PetscInt, pointer :: inactive_rows_local_ghosted(:)
  1071)   PetscInt :: n_inactive_rows
  1072)   PetscBool :: inactive_cells_exist
  1073)   type(option_type) :: option
  1074)   
  1075)   PetscInt :: ncount, idof
  1076)   PetscInt :: local_id, ghosted_id
  1077)   PetscInt :: ndof, n_active_dof
  1078) 
  1079)   type(grid_type), pointer :: grid
  1080)   PetscInt :: flag
  1081)   PetscErrorCode :: ierr
  1082)     
  1083)   flag = 0
  1084)   grid => patch%grid
  1085)   ndof = size(dof_is_active)
  1086)   n_active_dof = 0
  1087)   do idof = 1, ndof
  1088)     if (dof_is_active(idof)) n_active_dof = n_active_dof + 1
  1089)   enddo
  1090)   
  1091)   n_inactive_rows = 0
  1092)   inactive_cells_exist = PETSC_FALSE
  1093)   call DeallocateArray(inactive_rows_local)
  1094)   call DeallocateArray(inactive_rows_local_ghosted)
  1095) 
  1096)   do local_id = 1, grid%nlmax
  1097)     ghosted_id = grid%nL2G(local_id)
  1098)     if (patch%imat(ghosted_id) <= 0) then
  1099)       n_inactive_rows = n_inactive_rows + ndof
  1100)     else if (n_active_dof < ndof) then
  1101)       n_inactive_rows = n_inactive_rows + (ndof-n_active_dof)
  1102)     endif
  1103)   enddo
  1104) 
  1105)   allocate(inactive_rows_local(n_inactive_rows))
  1106)   allocate(inactive_rows_local_ghosted(n_inactive_rows))
  1107) 
  1108)   inactive_rows_local = 0
  1109)   inactive_rows_local_ghosted = 0
  1110)   ncount = 0
  1111) 
  1112)   do local_id = 1, grid%nlmax
  1113)     ghosted_id = grid%nL2G(local_id)
  1114)     if (patch%imat(ghosted_id) <= 0) then
  1115)       do idof = 1, ndof
  1116)         ncount = ncount + 1
  1117)         ! 1-based indexing
  1118)         inactive_rows_local(ncount) = (local_id-1)*ndof+idof
  1119)         ! 0-based indexing
  1120)         inactive_rows_local_ghosted(ncount) = (ghosted_id-1)*ndof+idof-1
  1121)       enddo
  1122)     else if (n_active_dof < ndof) then
  1123)       do idof = 1, ndof
  1124)         if (dof_is_active(idof)) cycle
  1125)         ncount = ncount + 1
  1126)         ! 1-based indexing
  1127)         inactive_rows_local(ncount) = (local_id-1)*ndof+idof
  1128)         ! 0-based indexing
  1129)         inactive_rows_local_ghosted(ncount) = (ghosted_id-1)*ndof+idof-1
  1130)       enddo
  1131)     endif
  1132)   enddo
  1133) 
  1134)   call MPI_Allreduce(n_inactive_rows,flag,ONE_INTEGER_MPI,MPIU_INTEGER, &
  1135)                      MPI_MAX,option%mycomm,ierr)
  1136)   if (flag > 0) then
  1137)     inactive_cells_exist = PETSC_TRUE
  1138)   endif
  1139)      
  1140)   if (ncount /= n_inactive_rows) then
  1141)     print *, 'Error:  Mismatch in non-zero row count!', ncount, n_inactive_rows
  1142)     stop
  1143)   endif
  1144) 
  1145) end subroutine InitSubsurfaceCreateZeroArray
  1146) 
  1147) end module Init_Subsurface_module

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