realization_subsurface.F90       coverage:  82.86 %func     63.04 %block


     1) module Realization_Subsurface_class
     2)   
     3)   use Realization_Base_class
     4) 
     5)   use Option_module
     6)   use Output_Aux_module
     7)   use Input_Aux_module
     8)   use Region_module
     9)   use Condition_module
    10)   use Transport_Constraint_module
    11)   use Material_module
    12)   use Saturation_Function_module
    13)   use Characteristic_Curves_module
    14)   use Dataset_Base_class
    15)   use Fluid_module
    16)   use Discretization_module
    17)   use Field_module
    18)   use Debug_module
    19)   use Uniform_Velocity_module
    20)   use Output_Aux_module
    21)   
    22)   use Reaction_Aux_module
    23)   
    24)   use Patch_module
    25)   
    26)   use PFLOTRAN_Constants_module
    27) 
    28)   implicit none
    29) 
    30) private
    31) 
    32) #include "petsc/finclude/petscsys.h"
    33) #include "petsc/finclude/petscvec.h"
    34) #include "petsc/finclude/petscvec.h90"
    35)   type, public, extends(realization_base_type) :: realization_subsurface_type
    36) 
    37)     type(region_list_type), pointer :: region_list
    38)     type(condition_list_type), pointer :: flow_conditions
    39)     type(tran_condition_list_type), pointer :: transport_conditions
    40)     type(tran_constraint_list_type), pointer :: transport_constraints
    41)     
    42)     type(tran_constraint_type), pointer :: sec_transport_constraint
    43)     type(material_property_type), pointer :: material_properties
    44)     type(fluid_property_type), pointer :: fluid_properties
    45)     type(fluid_property_type), pointer :: fluid_property_array(:)
    46)     type(saturation_function_type), pointer :: saturation_functions
    47)     class(characteristic_curves_type), pointer :: characteristic_curves
    48)     class(dataset_base_type), pointer :: datasets
    49)     
    50)     type(uniform_velocity_dataset_type), pointer :: uniform_velocity_dataset
    51)     character(len=MAXSTRINGLENGTH) :: nonuniform_velocity_filename
    52)     
    53)   end type realization_subsurface_type
    54) 
    55)   interface RealizationCreate
    56)     module procedure RealizationCreate1
    57)     module procedure RealizationCreate2
    58)   end interface
    59)   
    60)   public :: RealizationCreate, &
    61)             RealizationStrip, &
    62)             RealizationDestroyLegacy, &
    63)             RealizationProcessCouplers, &
    64)             RealizationInitAllCouplerAuxVars, &
    65)             RealizationProcessConditions, &
    66)             RealizationProcessDatasets, &
    67)             RealizationAddWaypointsToList, &
    68)             RealizationCreateDiscretization, &
    69)             RealizationLocalizeRegions, &
    70)             RealizationAddCoupler, &
    71)             RealizationAddStrata, &
    72)             RealizUpdateUniformVelocity, &
    73)             RealizationRevertFlowParameters, &
    74) !            RealizationGetVariable, &
    75) !            RealizGetVariableValueAtCell, &
    76) !            RealizationSetVariable, &
    77)             RealizationPrintCouplers, &
    78)             RealizationInitConstraints, &
    79)             RealProcessMatPropAndSatFunc, &
    80)             RealProcessFluidProperties, &
    81)             RealizationUpdatePropertiesTS, &
    82)             RealizationUpdatePropertiesNI, &
    83)             RealizationCalcMineralPorosity, &
    84)             RealizationCountCells, &
    85)             RealizationPrintGridStatistics, &
    86)             RealizationPassPtrsToPatches, &
    87)             RealLocalToLocalWithArray, &
    88)             RealizationCalculateCFL1Timestep, &
    89)             RealizUpdateAllCouplerAuxVars, &
    90)             RealizUnInitializedVarsFlow, &
    91)             RealizUnInitializedVarsTran, &
    92)             RealizSetSoilReferencePressure
    93) 
    94)   !TODO(intel)
    95)   ! public from Realization_Base_class
    96)   !public :: RealizationGetVariable
    97) 
    98) contains
    99) 
   100) ! ************************************************************************** !
   101) 
   102) function RealizationCreate1()
   103)   ! 
   104)   ! Allocates and initializes a new Realization object
   105)   ! 
   106)   ! Author: Glenn Hammond
   107)   ! Date: 10/25/07
   108)   ! 
   109) 
   110)   implicit none
   111)   
   112)   class(realization_subsurface_type), pointer :: RealizationCreate1
   113)   
   114)   class(realization_subsurface_type), pointer :: realization
   115)   type(option_type), pointer :: option
   116)   
   117)   nullify(option)
   118)   RealizationCreate1 => RealizationCreate2(option)
   119)   
   120) end function RealizationCreate1  
   121) 
   122) ! ************************************************************************** !
   123) 
   124) function RealizationCreate2(option)
   125)   ! 
   126)   ! Allocates and initializes a new Realization object
   127)   ! 
   128)   ! Author: Glenn Hammond
   129)   ! Date: 10/25/07
   130)   ! 
   131) 
   132)   implicit none
   133)   
   134)   type(option_type), pointer :: option
   135)   
   136)   class(realization_subsurface_type), pointer :: RealizationCreate2
   137)   
   138)   class(realization_subsurface_type), pointer :: realization
   139)   
   140)   allocate(realization)
   141)   call RealizationBaseInit(realization,option)
   142) 
   143)   allocate(realization%region_list)
   144)   call RegionInitList(realization%region_list)
   145) 
   146)   allocate(realization%flow_conditions)
   147)   call FlowConditionInitList(realization%flow_conditions)
   148)   allocate(realization%transport_conditions)
   149)   call TranConditionInitList(realization%transport_conditions)
   150)   allocate(realization%transport_constraints)
   151)   call TranConstraintInitList(realization%transport_constraints)
   152) 
   153)   nullify(realization%material_properties)
   154)   nullify(realization%fluid_properties)
   155)   nullify(realization%fluid_property_array)
   156)   nullify(realization%saturation_functions)
   157)   nullify(realization%characteristic_curves)
   158)   nullify(realization%datasets)
   159)   nullify(realization%uniform_velocity_dataset)
   160)   nullify(realization%sec_transport_constraint)
   161)   realization%nonuniform_velocity_filename = ''
   162) 
   163)   RealizationCreate2 => realization
   164)   
   165) end function RealizationCreate2 
   166) 
   167) ! ************************************************************************** !
   168) 
   169) subroutine RealizationCreateDiscretization(realization)
   170)   ! 
   171)   ! Creates grid
   172)   ! 
   173)   ! Author: Glenn Hammond
   174)   ! Date: 02/22/08
   175)   ! 
   176) 
   177)   use Grid_module
   178)   use Grid_Unstructured_Aux_module
   179)   use Grid_Unstructured_module, only : UGridEnsureRightHandRule
   180)   use Grid_Structured_module, only : StructGridCreateTVDGhosts
   181)   use Coupler_module
   182)   use Discretization_module
   183)   use Grid_Unstructured_Cell_module
   184)   use DM_Kludge_module
   185)   use Variables_module, only : VOLUME
   186)   use Communicator_Structured_class, only : StructuredCommunicatorCreate
   187)   use Communicator_Unstructured_class, only : UnstructuredCommunicatorCreate
   188)   
   189)   implicit none
   190)   
   191) #include "petsc/finclude/petscvec.h"
   192) #include "petsc/finclude/petscvec.h90"
   193) 
   194)   class(realization_subsurface_type) :: realization
   195)   
   196)   type(discretization_type), pointer :: discretization
   197)   type(grid_type), pointer :: grid
   198)   type(field_type), pointer :: field
   199)   type(option_type), pointer :: option
   200)   type(coupler_type), pointer :: boundary_condition
   201)   PetscErrorCode :: ierr
   202)   PetscInt :: ivar
   203) 
   204)   option => realization%option
   205)   field => realization%field
   206)  
   207)   discretization => realization%discretization
   208)   
   209)   call DiscretizationCreateDMs(discretization, option%nflowdof, &
   210)                                option%ntrandof, option%nphase, &
   211)                                option%ngeomechdof, option%n_stress_strain_dof, &
   212)                                option)
   213) 
   214)   ! 1 degree of freedom, global
   215)   call DiscretizationCreateVector(discretization,ONEDOF,field%work, &
   216)                                   GLOBAL,option)
   217)   call DiscretizationDuplicateVector(discretization,field%work, &
   218)                                      field%porosity0)
   219)   call DiscretizationDuplicateVector(discretization,field%work, &
   220)                                      field%tortuosity0)
   221)   call DiscretizationDuplicateVector(discretization,field%work, &
   222)                                      field%volume0)
   223) !geh: this is now allocated in 
   224) !     init_subsurface.F90:SubsurfAllocMatPropDataStructs()
   225) !  call DiscretizationDuplicateVector(discretization,field%work, &
   226) !                                     field%compressibility0)
   227)   if (option%flow%transient_porosity) then
   228)     call DiscretizationDuplicateVector(discretization,field%work, &
   229)                                        field%porosity_base_store)
   230)     call DiscretizationDuplicateVector(discretization,field%work, &
   231)                                        field%porosity_t)
   232)     call DiscretizationDuplicateVector(discretization,field%work, &
   233)                                        field%porosity_tpdt)
   234)   endif
   235) 
   236)   ! 1 degree of freedom, local
   237)   call DiscretizationCreateVector(discretization,ONEDOF,field%work_loc, &
   238)                                   LOCAL,option)
   239)   
   240)   if (option%nflowdof > 0) then
   241) 
   242)     ! 1-dof global  
   243)     call DiscretizationDuplicateVector(discretization,field%work, &
   244)                                        field%perm0_xx)
   245)     call DiscretizationDuplicateVector(discretization,field%work, &
   246)                                        field%perm0_yy)
   247)     call DiscretizationDuplicateVector(discretization,field%work, &
   248)                                        field%perm0_zz)
   249) 
   250)     ! 1-dof local
   251)     call DiscretizationDuplicateVector(discretization,field%work_loc, &
   252)                                        field%ithrm_loc)
   253)     call DiscretizationDuplicateVector(discretization,field%work_loc, &
   254)                                        field%icap_loc)
   255)     call DiscretizationDuplicateVector(discretization,field%work_loc, &
   256)                                        field%iphas_loc)
   257)     call DiscretizationDuplicateVector(discretization,field%work_loc, &
   258)                                        field%iphas_old_loc)
   259)     
   260)     ! ndof degrees of freedom, global
   261)     call DiscretizationCreateVector(discretization,NFLOWDOF,field%flow_xx, &
   262)                                     GLOBAL,option)
   263)     call DiscretizationDuplicateVector(discretization,field%flow_xx, &
   264)                                        field%flow_yy)
   265)     call DiscretizationDuplicateVector(discretization,field%flow_xx, &
   266)                                        field%flow_dxx)
   267)     call DiscretizationDuplicateVector(discretization,field%flow_xx, &
   268)                                        field%flow_r)
   269)     call DiscretizationDuplicateVector(discretization,field%flow_xx, &
   270)                                        field%flow_accum)
   271)     call DiscretizationDuplicateVector(discretization,field%flow_xx, &
   272)                                        field%flow_accum2)
   273) 
   274)     ! ndof degrees of freedom, local
   275)     call DiscretizationCreateVector(discretization,NFLOWDOF,field%flow_xx_loc, &
   276)                                     LOCAL,option)
   277)   endif
   278) 
   279)   if (option%ntrandof > 0) then
   280)     if (option%transport%reactive_transport_coupling == GLOBAL_IMPLICIT) then
   281)       ! ndof degrees of freedom, global
   282)       call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx, &
   283)                                       GLOBAL,option)
   284)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   285)                                          field%tran_yy)
   286)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   287)                                          field%tran_dxx)
   288)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   289)                                          field%tran_r)
   290)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   291)                                          field%tran_accum)
   292) 
   293)       ! ndof degrees of freedom, local
   294)       call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx_loc, &
   295)                                       LOCAL,option)
   296)                                       
   297)       if (realization%reaction%use_log_formulation) then
   298)         call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   299)                                            field%tran_log_xx)
   300)         call DiscretizationDuplicateVector(discretization,field%tran_xx_loc, &
   301)                                            field%tran_work_loc)
   302)       endif
   303)  
   304)     else ! operator splitting
   305)       ! ndof degrees of freedom, global
   306)       ! create the 1 dof vector for solving the individual linear systems
   307)       call DiscretizationCreateVector(discretization,ONEDOF,field%tran_rhs_coef, &
   308)                                       GLOBAL,option)
   309)       ! create the ntran dof vector for storage of the solution
   310)       call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx, &
   311)                                       GLOBAL,option)
   312)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   313)                                          field%tran_yy)
   314)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   315)                                          field%tran_dxx)
   316)       call DiscretizationDuplicateVector(discretization,field%tran_xx, &
   317)                                          field%tran_rhs)
   318) 
   319)       ! ndof degrees of freedom, local
   320)       ! again, just for storage of the current colution
   321)       call DiscretizationCreateVector(discretization,NTRANDOF,field%tran_xx_loc, &
   322)                                       LOCAL,option)
   323) 
   324)     endif
   325)     
   326)   endif
   327) 
   328)   select case(discretization%itype)
   329)     case(STRUCTURED_GRID)
   330)       grid => discretization%grid
   331)       ! set up nG2L, nL2G, etc.
   332)       call GridMapIndices(grid, &
   333)                           discretization%dm_1dof, &
   334)                           discretization%stencil_type,&
   335)                           option)
   336)       if (option%itranmode == EXPLICIT_ADVECTION) then
   337)         call StructGridCreateTVDGhosts(grid%structured_grid, &
   338)                                        realization%reaction%naqcomp, &
   339)                                        field%tran_xx, &
   340)                                        discretization%dm_1dof%dm, &
   341)                                        field%tvd_ghosts, &
   342)                                        discretization%tvd_ghost_scatter, &
   343)                                        option)
   344)       endif
   345)       call GridComputeSpacing(grid,discretization%origin_global,option)
   346)       call GridComputeCoordinates(grid,discretization%origin_global,option)
   347)       call GridComputeVolumes(grid,field%volume0,option)
   348)       ! set up internal connectivity, distance, etc.
   349)       call GridComputeInternalConnect(grid,option)
   350)     case(UNSTRUCTURED_GRID)
   351)       grid => discretization%grid
   352)       ! set up nG2L, NL2G, etc.
   353)       call GridMapIndices(grid, &
   354)                           discretization%dm_1dof, &
   355)                           discretization%stencil_type,&
   356)                           option)
   357)       call GridComputeCoordinates(grid,discretization%origin_global,option, & 
   358)                                     discretization%dm_1dof%ugdm) 
   359)       if (grid%itype == IMPLICIT_UNSTRUCTURED_GRID) then
   360)         call UGridEnsureRightHandRule(grid%unstructured_grid,grid%x, &
   361)                                       grid%y,grid%z,grid%nG2A,grid%nL2G, &
   362)                                       option)
   363)       endif
   364)       ! set up internal connectivity, distance, etc.
   365)       call GridComputeInternalConnect(grid,option, &
   366)                                       discretization%dm_1dof%ugdm) 
   367)       call GridComputeVolumes(grid,field%volume0,option)
   368)   end select
   369)  
   370)   ! initialize to UNINITIALIZED_DOUBLE for check later that verifies all values 
   371)   ! have been set
   372)   call VecSet(field%porosity0,UNINITIALIZED_DOUBLE,ierr);CHKERRQ(ierr)
   373) 
   374)   ! Allocate vectors to hold temporally average output quantites
   375)   if (realization%output_option%aveg_output_variable_list%nvars>0) then
   376) 
   377)     field%nvars = realization%output_option%aveg_output_variable_list%nvars
   378)     allocate(field%avg_vars_vec(field%nvars))
   379) 
   380)     do ivar=1,field%nvars
   381)       call DiscretizationDuplicateVector(discretization,field%work, &
   382)                                          field%avg_vars_vec(ivar))
   383)       call VecSet(field%avg_vars_vec(ivar),0.d0,ierr);CHKERRQ(ierr)
   384)     enddo
   385)   endif
   386)        
   387)   ! Allocate vectors to hold flowrate quantities
   388)   if (realization%output_option%print_hdf5_mass_flowrate.or. &
   389)      realization%output_option%print_hdf5_energy_flowrate.or. &
   390)      realization%output_option%print_hdf5_aveg_mass_flowrate.or. &
   391)      realization%output_option%print_hdf5_aveg_energy_flowrate) then
   392)     call VecCreateMPI(option%mycomm, &
   393)         (option%nflowdof*MAX_FACE_PER_CELL+1)*realization%patch%grid%nlmax, &
   394)         PETSC_DETERMINE,field%flowrate_inst,ierr);CHKERRQ(ierr)
   395)     call VecSet(field%flowrate_inst,0.d0,ierr);CHKERRQ(ierr)
   396)   endif
   397) 
   398)   ! Allocate vectors to hold velocity at face
   399)   if (realization%output_option%print_hdf5_vel_face) then
   400) 
   401)     ! vx
   402)     call VecCreateMPI(option%mycomm, &
   403)         (option%nflowspec*MAX_FACE_PER_CELL+1)*realization%patch%grid%nlmax, &
   404)         PETSC_DETERMINE,field%vx_face_inst,ierr);CHKERRQ(ierr)
   405)     call VecSet(field%vx_face_inst,0.d0,ierr);CHKERRQ(ierr)
   406) 
   407)     ! vy and vz
   408)     call VecDuplicate(field%vx_face_inst,field%vy_face_inst, &
   409)                       ierr);CHKERRQ(ierr)
   410)     call VecDuplicate(field%vx_face_inst,field%vz_face_inst, &
   411)                       ierr);CHKERRQ(ierr)
   412)   endif
   413) 
   414)   if (realization%output_option%print_explicit_flowrate) then
   415)     call VecCreateMPI(option%mycomm, &
   416)          size(grid%unstructured_grid%explicit_grid%connections,2), &
   417)          PETSC_DETERMINE,field%flowrate_inst,ierr);CHKERRQ(ierr)
   418)     call VecSet(field%flowrate_inst,0.d0,ierr);CHKERRQ(ierr)
   419)   endif
   420)     
   421)   ! If average flowrate has to be saved, create a vector for it
   422)   if (realization%output_option%print_hdf5_aveg_mass_flowrate.or. &
   423)       realization%output_option%print_hdf5_aveg_energy_flowrate) then
   424)     call VecCreateMPI(option%mycomm, &
   425)         (option%nflowdof*MAX_FACE_PER_CELL+1)*realization%patch%grid%nlmax, &
   426)         PETSC_DETERMINE,field%flowrate_aveg,ierr);CHKERRQ(ierr)
   427)     call VecSet(field%flowrate_aveg,0.d0,ierr);CHKERRQ(ierr)
   428)   endif
   429) 
   430)   select case(realization%discretization%itype)
   431)     case(STRUCTURED_GRID)
   432)       realization%comm1 => StructuredCommunicatorCreate()
   433)     case(UNSTRUCTURED_GRID)
   434)       realization%comm1 => UnstructuredCommunicatorCreate()
   435)   end select
   436)   call realization%comm1%SetDM(discretization%dm_1dof)
   437) 
   438)   if (option%flow%quasi_3d) then
   439)     call RealizCreateFlowMassTransferVec(realization)
   440)   endif
   441) 
   442) end subroutine RealizationCreateDiscretization
   443) 
   444) ! ************************************************************************** !
   445) 
   446) subroutine RealizationLocalizeRegions(realization)
   447)   ! 
   448)   ! Localizes regions within each patch
   449)   ! 
   450)   ! Author: Glenn Hammond
   451)   ! Date: 02/22/08
   452)   ! 
   453) 
   454)   use Option_module
   455)   use String_module
   456)   use Grid_module
   457) 
   458)   implicit none
   459)   
   460)   class(realization_subsurface_type) :: realization
   461)   
   462)   type (region_type), pointer :: cur_region, cur_region2
   463)   type(option_type), pointer :: option
   464) 
   465)   option => realization%option
   466) 
   467)   ! check to ensure that region names are not duplicated
   468)   cur_region => realization%region_list%first
   469)   do
   470)     if (.not.associated(cur_region)) exit
   471)     cur_region2 => cur_region%next
   472)     do
   473)       if (.not.associated(cur_region2)) exit
   474)       if (StringCompare(cur_region%name,cur_region2%name,MAXWORDLENGTH)) then
   475)         option%io_buffer = 'Duplicate region names: ' // trim(cur_region%name)
   476)         call printErrMsg(option)
   477)       endif
   478)       cur_region2 => cur_region2%next
   479)     enddo
   480)     cur_region => cur_region%next
   481)   enddo
   482) 
   483)   call PatchLocalizeRegions(realization%patch,realization%region_list, &
   484)                             realization%option)
   485) 
   486) end subroutine RealizationLocalizeRegions
   487) 
   488) ! ************************************************************************** !
   489) 
   490) subroutine RealizationPassPtrsToPatches(realization)
   491)   ! 
   492)   ! Sets patch%field => realization%field
   493)   ! 
   494)   ! Author: Glenn Hammond
   495)   ! Date: 01/12/11
   496)   ! 
   497) 
   498)   use Option_module
   499) 
   500)   implicit none
   501)   
   502)   class(realization_subsurface_type) :: realization
   503)   
   504)   realization%patch%field => realization%field
   505)   realization%patch%datasets => realization%datasets
   506)   realization%patch%reaction => realization%reaction
   507)   
   508) end subroutine RealizationPassPtrsToPatches
   509) 
   510) ! ************************************************************************** !
   511) 
   512) subroutine RealizationAddCoupler(realization,coupler)
   513)   ! 
   514)   ! Adds a copy of a coupler to a list
   515)   ! 
   516)   ! Author: Glenn Hammond
   517)   ! Date: 02/22/08
   518)   ! 
   519) 
   520)   use Coupler_module
   521) 
   522)   implicit none
   523)   
   524)   class(realization_subsurface_type) :: realization
   525)   type(coupler_type), pointer :: coupler
   526)   
   527)   type(patch_type), pointer :: patch
   528)   
   529)   type(coupler_type), pointer :: new_coupler
   530)   
   531)   patch => realization%patch
   532)   
   533)   ! only add to flow list for now, since they will be split out later
   534)   new_coupler => CouplerCreate(coupler)
   535)   select case(coupler%itype)
   536)     case(BOUNDARY_COUPLER_TYPE)
   537)       call CouplerAddToList(new_coupler,patch%boundary_condition_list)
   538)     case(INITIAL_COUPLER_TYPE)
   539)       call CouplerAddToList(new_coupler,patch%initial_condition_list)
   540)     case(SRC_SINK_COUPLER_TYPE)
   541)       call CouplerAddToList(new_coupler,patch%source_sink_list)
   542)   end select
   543)   nullify(new_coupler)
   544) 
   545)   call CouplerDestroy(coupler)
   546)  
   547) end subroutine RealizationAddCoupler
   548) 
   549) ! ************************************************************************** !
   550) 
   551) subroutine RealizationAddStrata(realization,strata)
   552)   ! 
   553)   ! Adds a copy of a strata to a list
   554)   ! 
   555)   ! Author: Glenn Hammond
   556)   ! Date: 02/22/08
   557)   ! 
   558) 
   559)   use Strata_module
   560) 
   561)   implicit none
   562)   
   563)   class(realization_subsurface_type) :: realization
   564)   type(strata_type), pointer :: strata
   565)   
   566)   type(strata_type), pointer :: new_strata
   567)   
   568)   new_strata => StrataCreate(strata)
   569)   call StrataAddToList(new_strata,realization%patch%strata_list)
   570)   nullify(new_strata)
   571)   
   572)   call StrataDestroy(strata)
   573)  
   574) end subroutine RealizationAddStrata
   575) 
   576) 
   577) ! ************************************************************************** !
   578) 
   579) subroutine RealizationProcessCouplers(realization)
   580)   ! 
   581)   ! Sets connectivity and pointers for couplers
   582)   ! 
   583)   ! Author: Glenn Hammond
   584)   ! Date: 02/22/08
   585)   ! 
   586) 
   587)   use Option_module
   588) 
   589)   implicit none
   590)   
   591)   class(realization_subsurface_type) :: realization
   592)   
   593)   call PatchProcessCouplers( realization%patch,realization%flow_conditions, &
   594)                              realization%transport_conditions, &
   595)                              realization%option)
   596)   
   597) end subroutine RealizationProcessCouplers
   598) 
   599) ! ************************************************************************** !
   600) 
   601) subroutine RealizationProcessDatasets(realization)
   602)   ! 
   603)   ! Processes datasets before they are linked to anything else
   604)   ! 
   605)   ! Author: Glenn Hammond
   606)   ! Date: 01/20/16
   607)   ! 
   608)   use Dataset_module
   609)   
   610)   implicit none
   611)   
   612)   class(realization_subsurface_type) :: realization
   613) 
   614)   call DatasetScreenForNonCellIndexed(realization%datasets,realization%option)
   615)   
   616) end subroutine RealizationProcessDatasets
   617) 
   618) ! ************************************************************************** !
   619) 
   620) subroutine RealizationProcessConditions(realization)
   621)   ! 
   622)   ! Sets up auxiliary data associated with
   623)   ! conditions
   624)   ! 
   625)   ! Author: Glenn Hammond
   626)   ! Date: 10/14/08
   627)   ! 
   628)   use Data_Mediator_Base_class
   629)   use Data_Mediator_Dataset_class
   630)   
   631)   implicit none
   632)   
   633)   class(realization_subsurface_type) :: realization
   634)   class(data_mediator_base_type), pointer :: cur_data_mediator
   635) 
   636)   if (realization%option%nflowdof > 0) then
   637)     call RealProcessFlowConditions(realization)
   638)   endif
   639)   if (realization%option%ntrandof > 0) then
   640)     call RealProcessTranConditions(realization)
   641)   endif
   642)   
   643)   ! update data mediators
   644)   cur_data_mediator => realization%flow_data_mediator_list
   645)   do
   646)     if (.not.associated(cur_data_mediator)) exit
   647)     call RealizCreateFlowMassTransferVec(realization)
   648)     select type(cur_data_mediator)
   649)       class is(data_mediator_dataset_type)
   650)         call DataMediatorDatasetInit(cur_data_mediator, &
   651)                                      realization%discretization, &
   652)                                      realization%datasets, &
   653)                                      realization%option)
   654)         call cur_data_mediator%Update(realization%field%flow_mass_transfer, &
   655)                                       realization%option)
   656)       class default
   657)     end select
   658)     cur_data_mediator => cur_data_mediator%next
   659)   enddo
   660) 
   661)   cur_data_mediator => realization%tran_data_mediator_list
   662)   do
   663)     if (.not.associated(cur_data_mediator)) exit
   664)     call RealizCreateTranMassTransferVec(realization)
   665)     select type(cur_data_mediator)
   666)       class is(data_mediator_dataset_type)
   667)         call DataMediatorDatasetInit(cur_data_mediator, &
   668)                                      realization%discretization, &
   669)                                      realization%datasets, &
   670)                                      realization%option)
   671)         call cur_data_mediator%Update(realization%field%tran_mass_transfer, &
   672)                                       realization%option)
   673)       class default
   674)     end select
   675)     cur_data_mediator => cur_data_mediator%next
   676)   enddo
   677) 
   678) end subroutine RealizationProcessConditions
   679) 
   680) ! ************************************************************************** !
   681) 
   682) subroutine RealProcessMatPropAndSatFunc(realization)
   683)   ! 
   684)   ! Sets up linkeage between material properties
   685)   ! and saturation function, auxiliary arrays
   686)   ! and datasets
   687)   ! 
   688)   ! Author: Glenn Hammond
   689)   ! Date: 01/21/09, 01/12/11
   690)   ! 
   691) 
   692)   use String_module
   693)   use Dataset_Common_HDF5_class
   694)   use Dataset_module
   695)   
   696)   implicit none
   697)   
   698)   class(realization_subsurface_type) :: realization
   699)   
   700)   PetscBool :: found
   701)   PetscInt :: i
   702)   type(option_type), pointer :: option
   703)   type(material_property_type), pointer :: cur_material_property
   704)   type(patch_type), pointer :: patch
   705)   character(len=MAXSTRINGLENGTH) :: string
   706)   class(dataset_base_type), pointer :: dataset
   707) 
   708)   option => realization%option
   709)   patch => realization%patch
   710)   
   711)   ! set up mirrored pointer arrays within patches to saturation functions
   712)   ! and material properties
   713)   patch%material_properties => realization%material_properties
   714)   call MaterialPropConvertListToArray(patch%material_properties, &
   715)                                       patch%material_property_array, &
   716)                                       option)
   717)   if (associated(realization%saturation_functions)) then
   718)     patch%saturation_functions => realization%saturation_functions
   719)     call SaturatFuncConvertListToArray(patch%saturation_functions, &
   720)                                        patch%saturation_function_array, &
   721)                                        option)
   722)   endif
   723)   if (associated(realization%characteristic_curves)) then
   724)     patch%characteristic_curves => realization%characteristic_curves
   725)     call CharCurvesConvertListToArray(patch%characteristic_curves, &
   726)                                       patch%characteristic_curves_array, &
   727)                                       option)
   728)   endif
   729)                                       
   730)   ! create mapping of internal to external material id
   731)   call MaterialCreateIntToExtMapping(patch%material_property_array, &
   732)                                      patch%imat_internal_to_external)
   733)     
   734)   cur_material_property => realization%material_properties                            
   735)   do                                      
   736)     if (.not.associated(cur_material_property)) exit
   737) 
   738)     ! obtain saturation function id
   739)     if (option%iflowmode /= NULL_MODE) then
   740)       if (associated(patch%saturation_function_array)) then
   741)         cur_material_property%saturation_function_id = &
   742)           SaturationFunctionGetID(patch%saturation_functions, &
   743)                              cur_material_property%saturation_function_name, &
   744)                              cur_material_property%name,option)
   745)       endif
   746)       if (associated(patch%characteristic_curves_array)) then
   747)         cur_material_property%saturation_function_id = &
   748)           CharacteristicCurvesGetID(patch%characteristic_curves_array, &
   749)                              cur_material_property%saturation_function_name, &
   750)                              cur_material_property%name,option)
   751)       endif
   752)       if (cur_material_property%saturation_function_id == 0) then
   753)         option%io_buffer = 'Characteristic curve "' // &
   754)           trim(cur_material_property%saturation_function_name) // &
   755)           '" not found.'
   756)         call printErrMsg(option)
   757)       endif
   758)     endif
   759)     
   760)     ! if named, link dataset to property
   761)     if (associated(cur_material_property%porosity_dataset)) then
   762)       string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
   763)                '),POROSITY'
   764)       dataset => &
   765)         DatasetBaseGetPointer(realization%datasets, &
   766)                               cur_material_property%porosity_dataset%name, &
   767)                               string,option)
   768)       call DatasetDestroy(cur_material_property%porosity_dataset)
   769)       select type(dataset)
   770)         class is (dataset_common_hdf5_type)
   771)           cur_material_property%porosity_dataset => dataset
   772)         class default
   773)           option%io_buffer = 'Incorrect dataset type for porosity.'
   774)           call printErrMsg(option)
   775)       end select
   776)     endif
   777)     if (associated(cur_material_property%tortuosity_dataset)) then
   778)       string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
   779)                '),TORTUOSITY'
   780)       dataset => &
   781)         DatasetBaseGetPointer(realization%datasets, &
   782)                               cur_material_property%tortuosity_dataset%name, &
   783)                               string,option)
   784)       call DatasetDestroy(cur_material_property%tortuosity_dataset)
   785)       select type(dataset)
   786)         class is (dataset_common_hdf5_type)
   787)           cur_material_property%tortuosity_dataset => dataset
   788)         class default
   789)           option%io_buffer = 'Incorrect dataset type for tortuosity.'
   790)           call printErrMsg(option)
   791)       end select
   792)     endif
   793)     if (associated(cur_material_property%permeability_dataset)) then
   794)       string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
   795)                '),PERMEABILITY or PERMEABILITY X'
   796)       dataset => &
   797)         DatasetBaseGetPointer(realization%datasets, &
   798)                             cur_material_property%permeability_dataset%name, &
   799)                             string,option)
   800)       call DatasetDestroy(cur_material_property%permeability_dataset)
   801)       select type(dataset)
   802)         class is (dataset_common_hdf5_type)
   803)           cur_material_property%permeability_dataset => dataset
   804)         class default
   805)           option%io_buffer = 'Incorrect dataset type for permeability or &
   806)                              &permeability X.'
   807)           call printErrMsg(option)
   808)       end select      
   809)     endif
   810)     if (associated(cur_material_property%permeability_dataset_y)) then
   811)       string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
   812)                '),PERMEABILITY Y'
   813)       dataset => &
   814)         DatasetBaseGetPointer(realization%datasets, &
   815)                             cur_material_property%permeability_dataset_y%name, &
   816)                             string,option)
   817)       call DatasetDestroy(cur_material_property%permeability_dataset_y)
   818)       select type(dataset)
   819)         class is (dataset_common_hdf5_type)
   820)           cur_material_property%permeability_dataset_y => dataset
   821)         class default
   822)           option%io_buffer = 'Incorrect dataset type for permeability Y.'
   823)           call printErrMsg(option)
   824)       end select      
   825)     endif
   826)     if (associated(cur_material_property%permeability_dataset_z)) then
   827)       string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
   828)                '),PERMEABILITY Z'
   829)       dataset => &
   830)         DatasetBaseGetPointer(realization%datasets, &
   831)                             cur_material_property%permeability_dataset_z%name, &
   832)                             string,option)
   833)       call DatasetDestroy(cur_material_property%permeability_dataset_z)
   834)       select type(dataset)
   835)         class is (dataset_common_hdf5_type)
   836)           cur_material_property%permeability_dataset_z => dataset
   837)         class default
   838)           option%io_buffer = 'Incorrect dataset type for permeability Z.'
   839)           call printErrMsg(option)
   840)       end select      
   841)     endif
   842)     if (associated(cur_material_property%compressibility_dataset)) then
   843)       string = 'MATERIAL_PROPERTY(' // trim(cur_material_property%name) // &
   844)                '),SOIL_COMPRESSIBILITY'
   845)       dataset => &
   846)         DatasetBaseGetPointer(realization%datasets, &
   847)                          cur_material_property%compressibility_dataset%name, &
   848)                          string,option)
   849)       call DatasetDestroy(cur_material_property%compressibility_dataset)
   850)       select type(dataset)
   851)         class is (dataset_common_hdf5_type)
   852)           cur_material_property%compressibility_dataset => dataset
   853)         class default
   854)           option%io_buffer = 'Incorrect dataset type for soil_compressibility.'
   855)           call printErrMsg(option)
   856)       end select      
   857)     endif
   858)     
   859)     cur_material_property => cur_material_property%next
   860)   enddo
   861)   
   862)   
   863) end subroutine RealProcessMatPropAndSatFunc
   864) 
   865) ! ************************************************************************** !
   866) 
   867) subroutine RealProcessFluidProperties(realization)
   868)   ! 
   869)   ! Sets up linkeage with fluid properties
   870)   ! 
   871)   ! Author: Glenn Hammond
   872)   ! Date: 01/21/09
   873)   ! 
   874)   
   875)   implicit none
   876)   
   877)   class(realization_subsurface_type) :: realization
   878)   
   879)   PetscBool :: found
   880)   type(option_type), pointer :: option
   881)   type(fluid_property_type), pointer :: cur_fluid_property
   882)   
   883)   option => realization%option
   884)   
   885)   found = PETSC_FALSE
   886)   cur_fluid_property => realization%fluid_properties                            
   887)   do                                      
   888)     if (.not.associated(cur_fluid_property)) exit
   889)     found = PETSC_TRUE
   890)     select case(trim(cur_fluid_property%phase_name))
   891)       case('LIQUID')
   892)         cur_fluid_property%phase_id = LIQUID_PHASE
   893)       case('GAS')
   894)         cur_fluid_property%phase_id = GAS_PHASE
   895)       case default
   896)         cur_fluid_property%phase_id = LIQUID_PHASE
   897)     end select
   898)     cur_fluid_property => cur_fluid_property%next
   899)   enddo
   900)   
   901)   if (option%ntrandof > 0 .and. .not.found) then
   902)     option%io_buffer = 'A fluid property must be present in input file' // &
   903)                        ' for solute transport'
   904)   endif
   905)   
   906) end subroutine RealProcessFluidProperties
   907) 
   908) ! ************************************************************************** !
   909) 
   910) subroutine RealProcessFlowConditions(realization)
   911)   ! 
   912)   ! Sets linkage of flow conditions to dataset
   913)   ! 
   914)   ! Author: Glenn Hammond
   915)   ! Date: 10/26/11
   916)   ! 
   917) 
   918)   use Dataset_Base_class
   919)   use Dataset_module
   920) 
   921)   implicit none
   922) 
   923)   class(realization_subsurface_type) :: realization
   924)   
   925)   type(flow_condition_type), pointer :: cur_flow_condition
   926)   type(flow_sub_condition_type), pointer :: cur_flow_sub_condition
   927)   type(option_type), pointer :: option
   928)   character(len=MAXSTRINGLENGTH) :: string
   929)   PetscInt :: i
   930)   
   931)   option => realization%option
   932)   
   933)   ! loop over flow conditions looking for linkage to datasets
   934)   cur_flow_condition => realization%flow_conditions%first
   935)   do
   936)     if (.not.associated(cur_flow_condition)) exit
   937)     string = 'flow_condition ' // trim(cur_flow_condition%name)
   938)     ! find datum dataset
   939)     call DatasetFindInList(realization%datasets,cur_flow_condition%datum, &
   940)                            cur_flow_condition%default_time_storage, &
   941)                            string,option)
   942)     select case(option%iflowmode)
   943)       case default
   944)         do i = 1, size(cur_flow_condition%sub_condition_ptr)
   945)           ! find dataset
   946)           call DatasetFindInList(realization%datasets, &
   947)                  cur_flow_condition%sub_condition_ptr(i)%ptr%dataset, &
   948)                  cur_flow_condition%default_time_storage, &
   949)                  string,option)
   950)           ! find gradient dataset
   951)           call DatasetFindInList(realization%datasets, &
   952)                  cur_flow_condition%sub_condition_ptr(i)%ptr%gradient, &
   953)                  cur_flow_condition%default_time_storage, &
   954)                  string,option)
   955)         enddo
   956)     end select
   957)     cur_flow_condition => cur_flow_condition%next
   958)   enddo
   959) 
   960) end subroutine RealProcessFlowConditions
   961) 
   962) ! ************************************************************************** !
   963) 
   964) subroutine RealProcessTranConditions(realization)
   965)   ! 
   966)   ! Sets up auxiliary data associated with
   967)   ! transport conditions
   968)   ! 
   969)   ! Author: Glenn Hammond
   970)   ! Date: 10/14/08
   971)   ! 
   972) 
   973)   use String_module
   974)   use Reaction_module
   975)   use Transport_Constraint_module
   976)   
   977)   implicit none
   978)   
   979)   class(realization_subsurface_type) :: realization
   980)   
   981)   
   982)   PetscBool :: found
   983)   type(option_type), pointer :: option
   984)   type(tran_condition_type), pointer :: cur_condition
   985)   type(tran_constraint_coupler_type), pointer :: cur_constraint_coupler
   986)   type(tran_constraint_type), pointer :: cur_constraint, another_constraint
   987)   
   988)   option => realization%option
   989)   
   990)   ! check for duplicate constraint names
   991)   cur_constraint => realization%transport_constraints%first
   992)   do
   993)     if (.not.associated(cur_constraint)) exit
   994)       another_constraint => cur_constraint%next
   995)       ! now compare names
   996)       found = PETSC_FALSE
   997)       do
   998)         if (.not.associated(another_constraint)) exit
   999)         if (StringCompare(cur_constraint%name,another_constraint%name, &
  1000)             MAXWORDLENGTH)) then
  1001)           found = PETSC_TRUE
  1002)         endif
  1003)         another_constraint => another_constraint%next
  1004)       enddo
  1005)       if (found) then
  1006)         option%io_buffer = 'Duplicate transport constraints named "' // &
  1007)                  trim(cur_constraint%name) // '"'
  1008)         call printErrMsg(realization%option)
  1009)       endif
  1010)     cur_constraint => cur_constraint%next
  1011)   enddo
  1012)   
  1013)   ! initialize constraints
  1014)   cur_constraint => realization%transport_constraints%first
  1015)   do
  1016)     if (.not.associated(cur_constraint)) exit
  1017)     call ReactionProcessConstraint(realization%reaction, &
  1018)                                    cur_constraint%name, &
  1019)                                    cur_constraint%aqueous_species, &
  1020)                                    cur_constraint%free_ion_guess, &
  1021)                                    cur_constraint%minerals, &
  1022)                                    cur_constraint%surface_complexes, &
  1023)                                    cur_constraint%colloids, &
  1024)                                    cur_constraint%immobile_species, &
  1025)                                    realization%option)
  1026)     cur_constraint => cur_constraint%next
  1027)   enddo
  1028)   
  1029)   if (option%use_mc) then
  1030)     call ReactionProcessConstraint(realization%reaction, &
  1031)                                    realization%sec_transport_constraint%name, &
  1032)                                    realization%sec_transport_constraint%aqueous_species, &
  1033)                                    realization%sec_transport_constraint%free_ion_guess, &
  1034)                                    realization%sec_transport_constraint%minerals, &
  1035)                                    realization%sec_transport_constraint%surface_complexes, &
  1036)                                    realization%sec_transport_constraint%colloids, &
  1037)                                    realization%sec_transport_constraint%immobile_species, &
  1038)                                    realization%option)
  1039)   endif
  1040)   
  1041)   ! tie constraints to couplers, if not already associated
  1042)   cur_condition => realization%transport_conditions%first
  1043)   do
  1044) 
  1045)     if (.not.associated(cur_condition)) exit
  1046)     cur_constraint_coupler => cur_condition%constraint_coupler_list
  1047)     do
  1048)       if (.not.associated(cur_constraint_coupler)) exit
  1049)       if (.not.associated(cur_constraint_coupler%aqueous_species)) then
  1050)         cur_constraint => realization%transport_constraints%first
  1051)         do
  1052)           if (.not.associated(cur_constraint)) exit
  1053)           if (StringCompare(cur_constraint%name, &
  1054)                              cur_constraint_coupler%constraint_name, &
  1055)                              MAXWORDLENGTH)) then
  1056)             cur_constraint_coupler%aqueous_species => cur_constraint%aqueous_species
  1057)             cur_constraint_coupler%free_ion_guess => cur_constraint%free_ion_guess
  1058)             cur_constraint_coupler%minerals => cur_constraint%minerals
  1059)             cur_constraint_coupler%surface_complexes => cur_constraint%surface_complexes
  1060)             cur_constraint_coupler%colloids => cur_constraint%colloids
  1061)             cur_constraint_coupler%immobile_species => cur_constraint%immobile_species
  1062)             exit
  1063)           endif
  1064)           cur_constraint => cur_constraint%next
  1065)         enddo
  1066)         if (.not.associated(cur_constraint_coupler%aqueous_species)) then
  1067)           option%io_buffer = 'Transport constraint "' // &
  1068)                    trim(cur_constraint_coupler%constraint_name) // &
  1069)                    '" not found in input file constraints.'
  1070)           call printErrMsg(realization%option)
  1071)         endif
  1072)       endif
  1073)       cur_constraint_coupler => cur_constraint_coupler%next
  1074)     enddo
  1075)     if (associated(cur_condition%constraint_coupler_list%next)) then ! more than one
  1076)       cur_condition%is_transient = PETSC_TRUE
  1077)     else
  1078)       cur_condition%is_transient = PETSC_FALSE
  1079)     endif
  1080)     cur_condition => cur_condition%next
  1081)   enddo
  1082)  
  1083)   ! final details for setup
  1084)   cur_condition => realization%transport_conditions%first
  1085)   do
  1086)     if (.not.associated(cur_condition)) exit
  1087)     ! is the condition transient?
  1088)     if (associated(cur_condition%constraint_coupler_list%next)) then ! more than one
  1089)       cur_condition%is_transient = PETSC_TRUE
  1090)     else
  1091)       cur_condition%is_transient = PETSC_FALSE
  1092)     endif
  1093)     ! set pointer to first constraint coupler
  1094)     cur_condition%cur_constraint_coupler => cur_condition%constraint_coupler_list
  1095)     
  1096)     cur_condition => cur_condition%next
  1097)   enddo
  1098) 
  1099) end subroutine RealProcessTranConditions
  1100) 
  1101) ! ************************************************************************** !
  1102) 
  1103) subroutine RealizationInitConstraints(realization)
  1104)   ! 
  1105)   ! Initializes constraint concentrations
  1106)   ! 
  1107)   ! Author: Glenn Hammond
  1108)   ! Date: 12/04/08
  1109)   ! 
  1110) 
  1111)   implicit none
  1112) 
  1113)   class(realization_subsurface_type) :: realization
  1114)   
  1115)   type(patch_type), pointer :: cur_patch
  1116)   
  1117)   cur_patch => realization%patch_list%first
  1118)   do
  1119)     if (.not.associated(cur_patch)) exit
  1120)     call PatchInitConstraints(cur_patch,realization%reaction, &
  1121)                               realization%option)
  1122)     cur_patch => cur_patch%next
  1123)   enddo
  1124)  
  1125) end subroutine RealizationInitConstraints
  1126) 
  1127) ! ************************************************************************** !
  1128) 
  1129) subroutine RealizationPrintCouplers(realization)
  1130)   ! 
  1131)   ! Print boundary and initial condition data
  1132)   ! 
  1133)   ! Author: Glenn Hammond
  1134)   ! Date: 10/28/08
  1135)   ! 
  1136) 
  1137)   use Coupler_module
  1138)   
  1139)   implicit none
  1140)   
  1141)   class(realization_subsurface_type) :: realization
  1142)   
  1143)   type(patch_type), pointer :: cur_patch
  1144)   type(coupler_type), pointer :: cur_coupler
  1145)   type(option_type), pointer :: option
  1146)   type(reaction_type), pointer :: reaction
  1147)  
  1148)   option => realization%option
  1149)   reaction => realization%reaction
  1150)  
  1151)   if (.not.OptionPrintToFile(option)) return
  1152)   
  1153)   cur_patch => realization%patch_list%first
  1154)   do
  1155)     if (.not.associated(cur_patch)) exit
  1156) 
  1157)     cur_coupler => cur_patch%initial_condition_list%first
  1158)     do
  1159)       if (.not.associated(cur_coupler)) exit
  1160)       call RealizationPrintCoupler(cur_coupler,reaction,option)    
  1161)       cur_coupler => cur_coupler%next
  1162)     enddo
  1163)      
  1164)     cur_coupler => cur_patch%boundary_condition_list%first
  1165)     do
  1166)       if (.not.associated(cur_coupler)) exit
  1167)       call RealizationPrintCoupler(cur_coupler,reaction,option)    
  1168)       cur_coupler => cur_coupler%next
  1169)     enddo
  1170)      
  1171)     cur_coupler => cur_patch%source_sink_list%first
  1172)     do
  1173)       if (.not.associated(cur_coupler)) exit
  1174)       call RealizationPrintCoupler(cur_coupler,reaction,option)    
  1175)       cur_coupler => cur_coupler%next
  1176)     enddo
  1177) 
  1178)     cur_patch => cur_patch%next
  1179)   enddo
  1180)     
  1181) end subroutine RealizationPrintCouplers
  1182) 
  1183) ! ************************************************************************** !
  1184) 
  1185) subroutine RealizationPrintCoupler(coupler,reaction,option)
  1186)   ! 
  1187)   ! Prints boundary and initial condition coupler
  1188)   ! 
  1189)   ! Author: Glenn Hammond
  1190)   ! Date: 10/28/08
  1191)   ! 
  1192) 
  1193)   use Coupler_module
  1194)   use Reaction_module
  1195)   
  1196)   implicit none
  1197)   
  1198)   type(coupler_type) :: coupler
  1199)   type(option_type) :: option
  1200)   type(reaction_type), pointer :: reaction
  1201)   
  1202)   character(len=MAXSTRINGLENGTH) :: string
  1203)   
  1204)   type(flow_condition_type), pointer :: flow_condition
  1205)   type(tran_condition_type), pointer :: tran_condition
  1206)   type(region_type), pointer :: region
  1207)   type(tran_constraint_coupler_type), pointer :: constraint_coupler
  1208)    
  1209) 98 format(40('=+'))
  1210) 99 format(80('-'))
  1211)   
  1212)   flow_condition => coupler%flow_condition
  1213)   tran_condition => coupler%tran_condition
  1214)   region => coupler%region
  1215) 
  1216)   write(option%fid_out,*)
  1217)   write(option%fid_out,98)
  1218) 
  1219) 
  1220)   select case(coupler%itype)
  1221)     case(INITIAL_COUPLER_TYPE)
  1222)       string = 'Initial Condition'
  1223)     case(BOUNDARY_COUPLER_TYPE)
  1224)       string = 'Boundary Condition'
  1225)     case(SRC_SINK_COUPLER_TYPE)
  1226)       string = 'Source Sink'
  1227)   end select
  1228)   write(option%fid_out,'(/,2x,a,/)') trim(string)
  1229) 
  1230)   write(option%fid_out,99)
  1231) 101 format(5x,'     Flow Condition: ',2x,a)
  1232)   if (associated(flow_condition)) &
  1233)     write(option%fid_out,101) trim(flow_condition%name)
  1234) 102 format(5x,'Transport Condition: ',2x,a)
  1235)   if (associated(tran_condition)) &
  1236)     write(option%fid_out,102) trim(tran_condition%name)
  1237) 103 format(5x,'             Region: ',2x,a)
  1238)   if (associated(region)) &
  1239)     write(option%fid_out,103) trim(region%name)
  1240)   write(option%fid_out,99)
  1241)   
  1242)   if (associated(flow_condition)) then
  1243)     call FlowConditionPrint(flow_condition,option)
  1244)   endif
  1245)   if (associated(tran_condition)) then
  1246)     constraint_coupler => tran_condition%cur_constraint_coupler
  1247)     write(option%fid_out,'(/,2x,''Transport Condition: '',a)') &
  1248)       trim(tran_condition%name)
  1249)     if (associated(reaction)) then
  1250)       call ReactionPrintConstraint(constraint_coupler,reaction,option)
  1251)       write(option%fid_out,'(/)')
  1252)       write(option%fid_out,99)
  1253)     endif
  1254)   endif
  1255)  
  1256) end subroutine RealizationPrintCoupler
  1257) 
  1258) ! ************************************************************************** !
  1259) 
  1260) subroutine RealizationInitAllCouplerAuxVars(realization)
  1261)   ! 
  1262)   ! RealizationInitCouplerAuxVars: Initializes coupler auxillary variables
  1263)   ! within list
  1264)   ! 
  1265)   ! Author: Glenn Hammond
  1266)   ! Date: 02/22/08
  1267)   ! 
  1268) 
  1269)   use Option_module
  1270) 
  1271)   implicit none
  1272)   
  1273)   class(realization_subsurface_type) :: realization
  1274)   
  1275)   !geh: Must update conditions prior to initializing the aux vars.  
  1276)   !     Otherwise, datasets will not have been read for routines such as
  1277)   !     hydrostatic and auxvars will be initialized to garbage.
  1278)   call FlowConditionUpdate(realization%flow_conditions,realization%option, &
  1279)                            realization%option%time)
  1280)   call TranConditionUpdate(realization%transport_conditions, &
  1281)                            realization%option, &
  1282)                            realization%option%time)  
  1283)   call PatchInitAllCouplerAuxVars(realization%patch,realization%option)
  1284)    
  1285) end subroutine RealizationInitAllCouplerAuxVars
  1286) 
  1287) ! ************************************************************************** !
  1288) 
  1289) subroutine RealizUpdateAllCouplerAuxVars(realization,force_update_flag)
  1290)   ! 
  1291)   ! Updates auxiliary variables associated
  1292)   ! with couplers in lis
  1293)   ! 
  1294)   ! Author: Glenn Hammond
  1295)   ! Date: 02/22/08
  1296)   ! 
  1297) 
  1298)   use Option_module
  1299) 
  1300)   implicit none
  1301)   
  1302)   class(realization_subsurface_type) :: realization
  1303)   PetscBool :: force_update_flag
  1304) 
  1305)   !TODO(geh): separate flow from transport in these calls
  1306)   call PatchUpdateAllCouplerAuxVars(realization%patch,force_update_flag, &
  1307)                                     realization%option)
  1308) 
  1309) end subroutine RealizUpdateAllCouplerAuxVars
  1310) 
  1311) ! ************************************************************************** !
  1312) 
  1313) subroutine RealizationRevertFlowParameters(realization)
  1314)   ! 
  1315)   ! Assigns initial porosity/perms to vecs
  1316)   ! 
  1317)   ! Author: Glenn Hammond
  1318)   ! Date: 05/09/08
  1319)   ! 
  1320) 
  1321)   use Option_module
  1322)   use Field_module
  1323)   use Discretization_module
  1324)   use Material_Aux_class
  1325)   use Variables_module
  1326) 
  1327)   implicit none
  1328)   
  1329)   class(realization_subsurface_type) :: realization
  1330)   
  1331)   type(field_type), pointer :: field
  1332)   type(option_type), pointer :: option
  1333)   type(discretization_type), pointer :: discretization
  1334)   type(material_type), pointer :: Material
  1335)   
  1336)   option => realization%option
  1337)   field => realization%field
  1338)   discretization => realization%discretization
  1339)   Material => realization%patch%aux%Material
  1340) 
  1341)   if (option%nflowdof > 0) then
  1342)     call DiscretizationGlobalToLocal(discretization,field%perm0_xx, &
  1343)                                      field%work_loc,ONEDOF)  
  1344)     call MaterialSetAuxVarVecLoc(Material,field%work_loc,PERMEABILITY_X, &
  1345)                                  ZERO_INTEGER)
  1346)     call DiscretizationGlobalToLocal(discretization,field%perm0_yy, &
  1347)                                      field%work_loc,ONEDOF)  
  1348)     call MaterialSetAuxVarVecLoc(Material,field%work_loc,PERMEABILITY_Y, &
  1349)                                  ZERO_INTEGER)
  1350)     call DiscretizationGlobalToLocal(discretization,field%perm0_zz, &
  1351)                                      field%work_loc,ONEDOF)  
  1352)     call MaterialSetAuxVarVecLoc(Material,field%work_loc,PERMEABILITY_Z, &
  1353)                                  ZERO_INTEGER)
  1354)   endif   
  1355)   call DiscretizationGlobalToLocal(discretization,field%porosity0, &
  1356)                                    field%work_loc,ONEDOF)  
  1357)   call MaterialSetAuxVarVecLoc(Material,field%work_loc,POROSITY, &
  1358)                                ZERO_INTEGER)
  1359)   call DiscretizationGlobalToLocal(discretization,field%tortuosity0, &
  1360)                                    field%work_loc,ONEDOF)  
  1361)   call MaterialSetAuxVarVecLoc(Material,field%work_loc,TORTUOSITY, &
  1362)                                ZERO_INTEGER)
  1363) 
  1364) end subroutine RealizationRevertFlowParameters
  1365) 
  1366) ! ************************************************************************** !
  1367) 
  1368) subroutine RealizUpdateUniformVelocity(realization)
  1369)   ! 
  1370)   ! Assigns uniform velocity for transport
  1371)   ! 
  1372)   ! Author: Glenn Hammond
  1373)   ! Date: 02/22/08
  1374)   ! 
  1375) 
  1376)   use Option_module
  1377) 
  1378)   implicit none
  1379)   
  1380)   class(realization_subsurface_type) :: realization
  1381)   
  1382)   call UniformVelocityDatasetUpdate(realization%option, &
  1383)                                     realization%option%time, &
  1384)                                     realization%uniform_velocity_dataset)
  1385)   call PatchUpdateUniformVelocity(realization%patch, &
  1386)                             realization%uniform_velocity_dataset%cur_value, &
  1387)                             realization%option)
  1388)  
  1389) end subroutine RealizUpdateUniformVelocity
  1390) 
  1391) ! ************************************************************************** !
  1392) 
  1393) subroutine RealizationAddWaypointsToList(realization,waypoint_list)
  1394)   ! 
  1395)   ! Creates waypoints associated with source/sinks
  1396)   ! boundary conditions, etc. and add to list
  1397)   ! 
  1398)   ! Author: Glenn Hammond
  1399)   ! Date: 11/01/07
  1400)   ! 
  1401) 
  1402)   use Option_module
  1403)   use Waypoint_module
  1404)   use Time_Storage_module
  1405)   use Data_Mediator_Base_class
  1406)   use Data_Mediator_Dataset_class
  1407)   use Strata_module
  1408) 
  1409)   implicit none
  1410)   
  1411)   class(realization_subsurface_type) :: realization
  1412)   type(waypoint_list_type) :: waypoint_list
  1413) 
  1414)   type(flow_condition_type), pointer :: cur_flow_condition
  1415)   type(tran_condition_type), pointer :: cur_tran_condition
  1416)   type(flow_sub_condition_type), pointer :: sub_condition
  1417)   type(tran_constraint_coupler_type), pointer :: cur_constraint_coupler
  1418)   class(data_mediator_base_type), pointer :: cur_data_mediator
  1419)   type(waypoint_type), pointer :: waypoint, cur_waypoint
  1420)   type(option_type), pointer :: option
  1421)   type(strata_type), pointer :: cur_strata
  1422)   PetscInt :: itime, isub_condition
  1423)   PetscReal :: temp_real, final_time
  1424)   PetscReal, pointer :: times(:)
  1425) 
  1426)   option => realization%option
  1427)   nullify(times)
  1428)   
  1429)   ! set flag for final output
  1430)   cur_waypoint => waypoint_list%first
  1431)   do
  1432)     if (.not.associated(cur_waypoint)) exit
  1433)     if (cur_waypoint%final) then
  1434)       cur_waypoint%print_snap_output = &
  1435)         realization%output_option%print_final_snap
  1436)       exit
  1437)     endif
  1438)     cur_waypoint => cur_waypoint%next
  1439)   enddo
  1440)   ! use final time in conditional below
  1441)   if (associated(cur_waypoint)) then
  1442)     final_time = cur_waypoint%time
  1443)   else
  1444)     option%io_buffer = 'Final time not found in RealizationAddWaypointsToList'
  1445)     call printErrMsg(option)
  1446)   endif
  1447) 
  1448)   ! add update of flow conditions
  1449)   cur_flow_condition => realization%flow_conditions%first
  1450)   do
  1451)     if (.not.associated(cur_flow_condition)) exit
  1452)     if (cur_flow_condition%sync_time_with_update) then
  1453)       do isub_condition = 1, cur_flow_condition%num_sub_conditions
  1454)         sub_condition => cur_flow_condition%sub_condition_ptr(isub_condition)%ptr
  1455)         !TODO(geh): check if this updated more than simply the flow_dataset (i.e. datum and gradient)
  1456)         !geh: followup - no, datum/gradient are not considered.  Should they be considered?
  1457)         call TimeStorageGetTimes(sub_condition%dataset%time_storage, option, &
  1458)                                 final_time, times)
  1459)         if (associated(times)) then
  1460)           if (size(times) > 1000) then
  1461)             option%io_buffer = 'For flow condition "' // &
  1462)               trim(cur_flow_condition%name) // &
  1463)               '" dataset "' // trim(sub_condition%name) // &
  1464)               '", the number of times is excessive for synchronization ' // &
  1465)               'with waypoints.'
  1466)             call printErrMsg(option)
  1467)           endif
  1468)           do itime = 1, size(times)
  1469)             waypoint => WaypointCreate()
  1470)             waypoint%time = times(itime)
  1471)             waypoint%update_conditions = PETSC_TRUE
  1472)             call WaypointInsertInList(waypoint,waypoint_list)
  1473)           enddo
  1474)           deallocate(times)
  1475)           nullify(times)
  1476)         endif
  1477)       enddo
  1478)     endif
  1479)     cur_flow_condition => cur_flow_condition%next
  1480)   enddo
  1481)       
  1482)   ! add update of transport conditions
  1483)   cur_tran_condition => realization%transport_conditions%first
  1484)   do
  1485)     if (.not.associated(cur_tran_condition)) exit
  1486)     if (cur_tran_condition%is_transient) then
  1487)       cur_constraint_coupler => cur_tran_condition%constraint_coupler_list
  1488)       do
  1489)         if (.not.associated(cur_constraint_coupler)) exit
  1490)         if (cur_constraint_coupler%time > 1.d-40) then
  1491)           waypoint => WaypointCreate()
  1492)           waypoint%time = cur_constraint_coupler%time
  1493)           waypoint%update_conditions = PETSC_TRUE
  1494)           call WaypointInsertInList(waypoint,waypoint_list)
  1495)         endif
  1496)         cur_constraint_coupler => cur_constraint_coupler%next
  1497)       enddo
  1498)     endif
  1499)     cur_tran_condition => cur_tran_condition%next
  1500)   enddo
  1501) 
  1502)   ! add update of velocity fields
  1503)   if (associated(realization%uniform_velocity_dataset)) then
  1504)     if (realization%uniform_velocity_dataset%times(1) > 1.d-40 .or. &
  1505)         size(realization%uniform_velocity_dataset%times) > 1) then
  1506)       do itime = 1, size(realization%uniform_velocity_dataset%times)
  1507)         waypoint => WaypointCreate()
  1508)         waypoint%time = realization%uniform_velocity_dataset%times(itime)
  1509)         waypoint%update_conditions = PETSC_TRUE
  1510)         call WaypointInsertInList(waypoint,waypoint_list)
  1511)       enddo
  1512)     endif
  1513)   endif
  1514)   
  1515)   ! add waypoints for flow mass transfer
  1516)   if (associated(realization%flow_data_mediator_list)) then
  1517)     cur_data_mediator => realization%flow_data_mediator_list
  1518)     do
  1519)       if (.not.associated(cur_data_mediator)) exit
  1520)       select type(cur_data_mediator)
  1521)         class is(data_mediator_dataset_type)
  1522)           if (associated(cur_data_mediator%dataset%time_storage)) then
  1523)             do itime = 1, cur_data_mediator%dataset%time_storage%max_time_index
  1524)               waypoint => WaypointCreate()
  1525)               waypoint%time = &
  1526)                 cur_data_mediator%dataset%time_storage%times(itime)
  1527)               waypoint%update_conditions = PETSC_TRUE
  1528)               call WaypointInsertInList(waypoint,waypoint_list)
  1529)             enddo
  1530)           endif
  1531)         class default
  1532)       end select 
  1533)       cur_data_mediator => cur_data_mediator%next
  1534)     enddo
  1535)   endif  
  1536) 
  1537)   ! add waypoints for rt mass transfer
  1538)   if (associated(realization%tran_data_mediator_list)) then
  1539)     cur_data_mediator => realization%tran_data_mediator_list
  1540)     do
  1541)       if (.not.associated(cur_data_mediator)) exit
  1542)       select type(cur_data_mediator)
  1543)         class is(data_mediator_dataset_type)
  1544)           if (associated(cur_data_mediator%dataset%time_storage)) then
  1545)             do itime = 1, cur_data_mediator%dataset%time_storage%max_time_index
  1546)               waypoint => WaypointCreate()
  1547)               waypoint%time = &
  1548)                 cur_data_mediator%dataset%time_storage%times(itime)
  1549)               waypoint%update_conditions = PETSC_TRUE
  1550)               call WaypointInsertInList(waypoint,waypoint_list)
  1551)             enddo
  1552)           endif
  1553)         class default
  1554)       end select           
  1555)       cur_data_mediator => cur_data_mediator%next
  1556)     enddo
  1557)   endif 
  1558) 
  1559)   ! add in strata that change over time
  1560)   cur_strata => realization%patch%strata_list%first
  1561)   do
  1562)     if (.not.associated(cur_strata)) exit
  1563)     if (Initialized(cur_strata%start_time)) then
  1564)       waypoint => WaypointCreate()
  1565)       waypoint%time = cur_strata%start_time
  1566)       waypoint%sync = PETSC_TRUE
  1567)       call WaypointInsertInList(waypoint,waypoint_list)
  1568)     endif
  1569)     if (Initialized(cur_strata%final_time)) then
  1570)       waypoint => WaypointCreate()
  1571)       waypoint%time = cur_strata%final_time
  1572)       waypoint%sync = PETSC_TRUE
  1573)       call WaypointInsertInList(waypoint,waypoint_list)
  1574)     endif
  1575)     cur_strata => cur_strata%next
  1576)   enddo
  1577) 
  1578) end subroutine RealizationAddWaypointsToList
  1579) 
  1580) ! ************************************************************************** !
  1581) 
  1582) subroutine RealizationUpdatePropertiesTS(realization)
  1583)   ! 
  1584)   ! Updates coupled properties at each grid cell
  1585)   ! 
  1586)   ! Author: Glenn Hammond
  1587)   ! Date: 08/05/09
  1588)   ! 
  1589) 
  1590)   use Grid_module
  1591)   use Reactive_Transport_Aux_module
  1592)   use Material_Aux_class
  1593)   use Variables_module, only : POROSITY, TORTUOSITY, PERMEABILITY_X, &
  1594)                                PERMEABILITY_Y, PERMEABILITY_Z
  1595)  
  1596)   implicit none
  1597)   
  1598)   class(realization_subsurface_type) :: realization
  1599) 
  1600)   type(option_type), pointer :: option
  1601)   type(patch_type), pointer :: patch
  1602)   type(field_type), pointer :: field
  1603)   type(reaction_type), pointer :: reaction
  1604)   type(grid_type), pointer :: grid
  1605)   type(material_property_ptr_type), pointer :: material_property_array(:)
  1606)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:) 
  1607)   type(discretization_type), pointer :: discretization
  1608)   class(material_auxvar_type), pointer :: material_auxvars(:)
  1609) 
  1610)   PetscInt :: local_id, ghosted_id
  1611)   PetscInt :: imnrl, imnrl1, imnrl_armor, imat
  1612)   PetscReal :: sum_volfrac
  1613)   PetscReal :: scale, porosity_scale, volfrac_scale
  1614)   PetscBool :: porosity_updated
  1615)   PetscReal, pointer :: vec_p(:)
  1616)   PetscReal, pointer :: porosity0_p(:)
  1617)   PetscReal, pointer :: tortuosity0_p(:)
  1618)   PetscReal, pointer :: perm0_xx_p(:), perm0_yy_p(:), perm0_zz_p(:)
  1619)   PetscReal, pointer :: perm_ptr(:)
  1620)   PetscReal :: min_value
  1621)   PetscReal :: critical_porosity
  1622)   PetscReal :: porosity_base
  1623)   PetscInt :: ivalue
  1624)   PetscErrorCode :: ierr
  1625) 
  1626)   option => realization%option
  1627)   discretization => realization%discretization
  1628)   patch => realization%patch
  1629)   field => realization%field
  1630)   reaction => realization%reaction
  1631)   grid => patch%grid
  1632)   material_property_array => patch%material_property_array
  1633)   rt_auxvars => patch%aux%RT%auxvars
  1634)   material_auxvars => patch%aux%Material%auxvars
  1635) 
  1636)   porosity_updated = PETSC_FALSE
  1637)   if (reaction%update_porosity) then
  1638)     porosity_updated = PETSC_TRUE
  1639)     call RealizationCalcMineralPorosity(realization)
  1640)   endif
  1641)   
  1642)   if ((porosity_updated .and. &
  1643)        (reaction%update_tortuosity .or. &
  1644)         reaction%update_permeability)) .or. &
  1645)       ! if porosity ratio is used in mineral surface area update, we must
  1646)       ! recalculate it every time.
  1647)       (reaction%update_mineral_surface_area .and. &
  1648)        reaction%update_mnrl_surf_with_porosity)) then
  1649)     call VecGetArrayF90(field%porosity0,porosity0_p,ierr);CHKERRQ(ierr)
  1650)     call VecGetArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1651)     do local_id = 1, grid%nlmax
  1652)       ghosted_id = grid%nL2G(local_id)
  1653)       vec_p(local_id) = material_auxvars(ghosted_id)%porosity_base / &
  1654)                         porosity0_p(local_id)
  1655)     enddo
  1656)     call VecRestoreArrayF90(field%porosity0,porosity0_p,ierr);CHKERRQ(ierr)
  1657)     call VecRestoreArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1658)   endif      
  1659) 
  1660)   if (reaction%update_mineral_surface_area) then
  1661) 
  1662)     if (reaction%update_mnrl_surf_with_porosity) then
  1663)       ! placing the get/restore array calls within the condition will
  1664)       ! avoid improper access.
  1665)       call VecGetArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1666)     endif
  1667) 
  1668)     do local_id = 1, grid%nlmax
  1669)       ghosted_id = grid%nL2G(local_id)
  1670)       do imnrl = 1, reaction%mineral%nkinmnrl
  1671) 
  1672)         porosity_scale = 1.d0
  1673)         if (reaction%update_mnrl_surf_with_porosity) then
  1674)           porosity_scale = vec_p(local_id)** &
  1675)              reaction%mineral%kinmnrl_surf_area_porosity_pwr(imnrl)
  1676) !       geh: srf_area_vol_frac_pwr must be defined on a per mineral basis, not
  1677) !       solely material type.
  1678) !       material_property_array(patch%imat(ghosted_id))%ptr%mnrl_surf_area_porosity_pwr
  1679)         endif
  1680) 
  1681)         volfrac_scale = 1.d0
  1682)         if (rt_auxvars(ghosted_id)%mnrl_volfrac0(imnrl) > 0.d0) then
  1683)           volfrac_scale = (rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl)/ &
  1684)                          rt_auxvars(ghosted_id)%mnrl_volfrac0(imnrl))** &
  1685)              reaction%mineral%kinmnrl_surf_area_vol_frac_pwr(imnrl)
  1686) !       geh: srf_area_vol_frac_pwr must be defined on a per mineral basis, not
  1687) !       solely material type.
  1688) !       material_property_array(patch%imat(ghosted_id))%ptr%mnrl_surf_area_volfrac_pwr
  1689) !         rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
  1690) !           rt_auxvars(ghosted_id)%mnrl_area0(imnrl)*porosity_scale*volfrac_scale
  1691) !       else
  1692) !         rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
  1693) !           rt_auxvars(ghosted_id)%mnrl_area0(imnrl)
  1694)         endif
  1695) 
  1696)         rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
  1697)             rt_auxvars(ghosted_id)%mnrl_area0(imnrl)*porosity_scale*volfrac_scale
  1698) 
  1699)         if (reaction%update_armor_mineral_surface .and. &
  1700)             reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl) > 0.d0) then
  1701)           imnrl_armor = imnrl
  1702)           do imnrl1 = 1, reaction%mineral%nkinmnrl
  1703)             if (reaction%mineral%kinmnrl_armor_min_names(imnrl) == &
  1704)                 reaction%mineral%kinmnrl_names(imnrl1)) then
  1705)               imnrl_armor = imnrl1
  1706)               exit
  1707)             endif
  1708)           enddo
  1709) 
  1710) !         print *,'update-armor: ',imnrl,imnrl_armor, &
  1711) !         reaction%mineral%kinmnrl_armor_min_names(imnrl_armor)
  1712) 
  1713) !       check for negative surface area armoring correction
  1714)           if (reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl) > &
  1715)               rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl_armor)) then
  1716) 
  1717)             if (reaction%update_armor_mineral_surface_flag == 0) then ! surface unarmored
  1718)               rt_auxvars(ghosted_id)%mnrl_area(imnrl) = &
  1719)                 rt_auxvars(ghosted_id)%mnrl_area(imnrl) * &
  1720)                 ((reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl) &
  1721)                 - rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl_armor))/ &
  1722)                 reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl))** &
  1723)                 reaction%mineral%kinmnrl_surf_area_vol_frac_pwr(imnrl)
  1724)             else
  1725)               rt_auxvars(ghosted_id)%mnrl_area(imnrl) = rt_auxvars(ghosted_id)%mnrl_area0(imnrl)
  1726)               reaction%update_armor_mineral_surface_flag = 0
  1727)             endif
  1728)           else
  1729)             rt_auxvars(ghosted_id)%mnrl_area(imnrl) = 0.d0
  1730)             reaction%update_armor_mineral_surface_flag = 1 ! surface armored
  1731)           endif
  1732)         endif
  1733) 
  1734) !       print *,'update min srf: ',imnrl,local_id,reaction%mineral%kinmnrl_names(imnrl), &
  1735) !       reaction%mineral%kinmnrl_armor_min_names(imnrl), &
  1736) !       reaction%update_armor_mineral_surface, &
  1737) !       rt_auxvars(ghosted_id)%mnrl_area(imnrl), &
  1738) !       reaction%mineral%kinmnrl_armor_pwr(imnrl), &
  1739) !       reaction%mineral%kinmnrl_armor_crit_vol_frac(imnrl), &
  1740) !       rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl_armor), &
  1741) !       rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl)
  1742)       enddo
  1743)     enddo
  1744) 
  1745)     if (reaction%update_mnrl_surf_with_porosity) then
  1746)       call VecRestoreArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1747)     endif
  1748) !geh:remove
  1749)     call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1750)                                  TORTUOSITY,ZERO_INTEGER)
  1751)     call DiscretizationLocalToLocal(discretization,field%work_loc, &
  1752)                                     field%work_loc,ONEDOF)
  1753)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1754)                                  TORTUOSITY,ZERO_INTEGER)
  1755)   endif
  1756)       
  1757)   if (reaction%update_tortuosity) then
  1758)     call VecGetArrayF90(field%tortuosity0,tortuosity0_p,ierr);CHKERRQ(ierr)
  1759)     call VecGetArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1760)     do local_id = 1, grid%nlmax
  1761)       ghosted_id = grid%nL2G(local_id)
  1762)       scale = vec_p(local_id)** &
  1763)         material_property_array(patch%imat(ghosted_id))%ptr%tortuosity_pwr
  1764)       material_auxvars(ghosted_id)%tortuosity = &
  1765)         tortuosity0_p(local_id)*scale
  1766)     enddo
  1767)     call VecRestoreArrayF90(field%tortuosity0,tortuosity0_p, &
  1768)                             ierr);CHKERRQ(ierr)
  1769)     call VecRestoreArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  1770)     call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1771)                                  TORTUOSITY,ZERO_INTEGER)
  1772)     call DiscretizationLocalToLocal(discretization,field%work_loc, &
  1773)                                     field%work_loc,ONEDOF)
  1774)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1775)                                  TORTUOSITY,ZERO_INTEGER)
  1776)   endif
  1777)       
  1778)   if (reaction%update_permeability) then
  1779)     call VecGetArrayF90(field%perm0_xx,perm0_xx_p,ierr);CHKERRQ(ierr)
  1780)     call VecGetArrayF90(field%perm0_zz,perm0_zz_p,ierr);CHKERRQ(ierr)
  1781)     call VecGetArrayF90(field%perm0_yy,perm0_yy_p,ierr);CHKERRQ(ierr)
  1782)     call VecGetArrayF90(field%porosity0,porosity0_p,ierr);CHKERRQ(ierr)
  1783)     do local_id = 1, grid%nlmax
  1784)       ghosted_id = grid%nL2G(local_id)
  1785)       imat = patch%imat(ghosted_id)
  1786)       critical_porosity = material_property_array(imat)%ptr% &
  1787)                             permeability_crit_por
  1788)       porosity_base = material_auxvars(ghosted_id)%porosity_base
  1789)       scale = 0.d0
  1790)       if (porosity_base > critical_porosity .and. &
  1791)           porosity0_p(local_id) > critical_porosity) then
  1792)         scale = ((porosity_base - critical_porosity) / &
  1793)                  (porosity0_p(local_id) - critical_porosity)) ** &
  1794)                 material_property_array(imat)%ptr%permeability_pwr
  1795)       endif
  1796)       scale = max(material_property_array(imat)%ptr% &
  1797)                     permeability_min_scale_fac,scale)
  1798)       !geh: this is a kludge for gfortran.  the code reports errors when 
  1799)       !     material_auxvars(ghosted_id)%permeability is used.
  1800)       ! This is not an issue with Intel
  1801) #if 1
  1802)       perm_ptr => material_auxvars(ghosted_id)%permeability
  1803)       perm_ptr(perm_xx_index) = perm0_xx_p(local_id)*scale
  1804)       perm_ptr(perm_yy_index) = perm0_yy_p(local_id)*scale
  1805)       perm_ptr(perm_zz_index) = perm0_zz_p(local_id)*scale
  1806) #else
  1807)       material_auxvars(ghosted_id)%permeability(perm_xx_index) = &
  1808)         perm0_xx_p(local_id)*scale
  1809)       material_auxvars(ghosted_id)%permeability(perm_yy_index) = &
  1810)         perm0_yy_p(local_id)*scale
  1811)       material_auxvars(ghosted_id)%permeability(perm_zz_index) = &
  1812)         perm0_zz_p(local_id)*scale
  1813) #endif
  1814)     enddo
  1815)     call VecRestoreArrayF90(field%perm0_xx,perm0_xx_p,ierr);CHKERRQ(ierr)
  1816)     call VecRestoreArrayF90(field%perm0_zz,perm0_zz_p,ierr);CHKERRQ(ierr)
  1817)     call VecRestoreArrayF90(field%perm0_yy,perm0_yy_p,ierr);CHKERRQ(ierr)
  1818) 
  1819)     call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1820)                                  PERMEABILITY_X,ZERO_INTEGER)
  1821)     call DiscretizationLocalToLocal(discretization,field%work_loc, &
  1822)                                     field%work_loc,ONEDOF)
  1823)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1824)                                  PERMEABILITY_X,ZERO_INTEGER)
  1825)     call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1826)                                  PERMEABILITY_Y,ZERO_INTEGER)
  1827)     call DiscretizationLocalToLocal(discretization,field%work_loc, &
  1828)                                     field%work_loc,ONEDOF)
  1829)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1830)                                  PERMEABILITY_Y,ZERO_INTEGER)
  1831)     call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1832)                                  PERMEABILITY_Z,ZERO_INTEGER)
  1833)     call DiscretizationLocalToLocal(discretization,field%work_loc, &
  1834)                                     field%work_loc,ONEDOF)
  1835)     call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1836)                                  PERMEABILITY_Z,ZERO_INTEGER)
  1837)   endif  
  1838)   
  1839)   ! perform check to ensure that porosity is bounded between 0 and 1
  1840)   ! since it is calculated as 1.d-sum_volfrac, it cannot be > 1
  1841)   call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1842)                                POROSITY,POROSITY_MINERAL)
  1843)   call DiscretizationLocalToGlobal(discretization,field%work_loc, &
  1844)                                   field%work,ONEDOF)
  1845)   call VecMin(field%work,ivalue,min_value,ierr);CHKERRQ(ierr)
  1846)   if (min_value < 0.d0) then
  1847)     write(option%io_buffer,*) 'Sum of mineral volume fractions has ' // &
  1848)       'exceeded 1.d0 at cell (note PETSc numbering): ', ivalue
  1849)     call printErrMsg(option)
  1850)   endif
  1851)    
  1852) end subroutine RealizationUpdatePropertiesTS
  1853) 
  1854) ! ************************************************************************** !
  1855) 
  1856) subroutine RealizationUpdatePropertiesNI(realization)
  1857)   ! 
  1858)   ! Updates coupled properties at each grid cell
  1859)   ! 
  1860)   ! Author: Glenn Hammond
  1861)   ! Date: 08/05/09
  1862)   ! 
  1863) 
  1864)   use Grid_module
  1865)   use Reactive_Transport_Aux_module
  1866)   use Material_Aux_class
  1867)   use Variables_module, only : POROSITY, TORTUOSITY, PERMEABILITY_X, &
  1868)                                PERMEABILITY_Y, PERMEABILITY_Z
  1869)  
  1870)   implicit none
  1871)   
  1872)   class(realization_subsurface_type) :: realization
  1873) 
  1874) #if 0
  1875)   type(option_type), pointer :: option
  1876)   type(patch_type), pointer :: patch
  1877)   type(field_type), pointer :: field
  1878)   type(reaction_type), pointer :: reaction
  1879)   type(grid_type), pointer :: grid
  1880)   type(material_property_ptr_type), pointer :: material_property_array(:)
  1881)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:) 
  1882)   type(discretization_type), pointer :: discretization
  1883)   class(material_auxvar_type), pointer :: material_auxvars(:)
  1884) 
  1885)   PetscInt :: local_id, ghosted_id
  1886)   PetscInt :: imnrl, imnrl1, imnrl_armor, imat
  1887)   PetscReal :: sum_volfrac
  1888)   PetscReal :: scale, porosity_scale, volfrac_scale
  1889)   PetscBool :: porosity_updated
  1890)   PetscReal, pointer :: vec_p(:)
  1891)   PetscReal, pointer :: porosity0_p(:)
  1892)   PetscReal, pointer :: porosity_mnrl_loc_p(:)
  1893)   PetscReal, pointer :: tortuosity0_p(:)
  1894)   PetscReal, pointer :: perm0_xx_p(:), perm0_yy_p(:), perm0_zz_p(:)
  1895)   PetscReal :: min_value  
  1896)   PetscInt :: ivalue
  1897)   PetscErrorCode :: ierr
  1898) 
  1899)   option => realization%option
  1900)   discretization => realization%discretization
  1901)   patch => realization%patch
  1902)   field => realization%field
  1903)   reaction => realization%reaction
  1904)   grid => patch%grid
  1905)   material_property_array => patch%material_property_array
  1906)   rt_auxvars => patch%aux%RT%auxvars
  1907)   material_auxvars => patch%aux%Material%auxvars
  1908) #endif
  1909) 
  1910) end subroutine RealizationUpdatePropertiesNI
  1911) 
  1912) ! ************************************************************************** !
  1913) 
  1914) subroutine RealizationCalcMineralPorosity(realization)
  1915)   ! 
  1916)   ! Calculates porosity based on the sum of mineral volume fractions
  1917)   ! 
  1918)   ! Author: Glenn Hammond
  1919)   ! Date: 11/03/14
  1920)   !
  1921) 
  1922)   use Grid_module
  1923)   use Reactive_Transport_Aux_module
  1924)   use Material_Aux_class
  1925)   use Variables_module, only : POROSITY
  1926)  
  1927)   implicit none
  1928)   
  1929)   class(realization_subsurface_type) :: realization
  1930) 
  1931)   type(option_type), pointer :: option
  1932)   type(patch_type), pointer :: patch
  1933)   type(field_type), pointer :: field
  1934)   type(reaction_type), pointer :: reaction
  1935)   type(grid_type), pointer :: grid
  1936)   type(reactive_transport_auxvar_type), pointer :: rt_auxvars(:) 
  1937)   type(discretization_type), pointer :: discretization
  1938)   class(material_auxvar_type), pointer :: material_auxvars(:)
  1939) 
  1940)   PetscInt :: local_id, ghosted_id
  1941)   PetscInt :: imnrl 
  1942)   PetscReal :: sum_volfrac
  1943)   PetscErrorCode :: ierr
  1944) 
  1945)   option => realization%option
  1946)   discretization => realization%discretization
  1947)   patch => realization%patch
  1948)   field => realization%field
  1949)   reaction => realization%reaction
  1950)   grid => patch%grid
  1951)   rt_auxvars => patch%aux%RT%auxvars
  1952)   material_auxvars => patch%aux%Material%auxvars
  1953) 
  1954)   if (reaction%mineral%nkinmnrl > 0) then
  1955)     do local_id = 1, grid%nlmax
  1956)       ghosted_id = grid%nL2G(local_id)
  1957)       ! Go ahead and compute for inactive cells since their porosity does
  1958)       ! not matter (avoid check on active/inactive)
  1959)       sum_volfrac = 0.d0
  1960)       do imnrl = 1, reaction%mineral%nkinmnrl
  1961)         sum_volfrac = sum_volfrac + &
  1962)                       rt_auxvars(ghosted_id)%mnrl_volfrac(imnrl)
  1963)       enddo 
  1964)       ! the adjusted porosity becomes:
  1965)       ! 1 - sum(mineral volume fractions), but is truncated.
  1966)       material_auxvars(ghosted_id)%porosity_base = &
  1967)         max(1.d0-sum_volfrac,reaction%minimum_porosity)
  1968)     enddo
  1969)   endif
  1970)   ! update ghosted porosities
  1971)   call MaterialGetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1972)                                POROSITY,POROSITY_MINERAL)
  1973)   call DiscretizationLocalToLocal(discretization,field%work_loc, &
  1974)                                   field%work_loc,ONEDOF)
  1975)   call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1976)                                POROSITY,POROSITY_CURRENT)
  1977)   call MaterialSetAuxVarVecLoc(patch%aux%Material,field%work_loc, &
  1978)                                POROSITY,POROSITY_MINERAL)
  1979) 
  1980) end subroutine RealizationCalcMineralPorosity
  1981) 
  1982) ! ************************************************************************** !
  1983) 
  1984) subroutine RealLocalToLocalWithArray(realization,array_id)
  1985)   ! 
  1986)   ! Takes an F90 array that is ghosted
  1987)   ! and updates the ghosted values
  1988)   ! 
  1989)   ! Author: Glenn Hammond
  1990)   ! Date: 06/09/11
  1991)   ! 
  1992) 
  1993)   use Grid_module
  1994) 
  1995)   implicit none
  1996) 
  1997)   class(realization_subsurface_type) :: realization
  1998)   PetscInt :: array_id
  1999)   
  2000)   type(patch_type), pointer :: patch
  2001)   type(grid_type), pointer :: grid
  2002)   type(field_type), pointer :: field
  2003) 
  2004)   field => realization%field
  2005)   patch => realization%patch
  2006) 
  2007)   grid => patch%grid
  2008)   select case(array_id)
  2009)     case(MATERIAL_ID_ARRAY)
  2010)       call GridCopyIntegerArrayToVec(grid,patch%imat,field%work_loc, &
  2011)                                      grid%ngmax)
  2012)     case(SATURATION_FUNCTION_ID_ARRAY)
  2013)       call GridCopyIntegerArrayToVec(grid,patch%sat_func_id, &
  2014)                                      field%work_loc, grid%ngmax)
  2015)   end select
  2016) 
  2017)   call DiscretizationLocalToLocal(realization%discretization,field%work_loc, &
  2018)                                   field%work_loc,ONEDOF)
  2019) 
  2020)   select case(array_id)
  2021)     case(MATERIAL_ID_ARRAY)
  2022)       call GridCopyVecToIntegerArray(grid,patch%imat,field%work_loc, &
  2023)                                       grid%ngmax)
  2024)     case(SATURATION_FUNCTION_ID_ARRAY)
  2025)       call GridCopyVecToIntegerArray(grid,patch%sat_func_id, &
  2026)                                       field%work_loc, grid%ngmax)
  2027)   end select
  2028) 
  2029) end subroutine RealLocalToLocalWithArray
  2030) 
  2031) ! ************************************************************************** !
  2032) 
  2033) subroutine RealizationCountCells(realization,global_total_count, &
  2034)                                  global_active_count,total_count,active_count)
  2035)   ! 
  2036)   ! Counts # of active and inactive grid cells
  2037)   ! 
  2038)   ! Author: Glenn Hammond
  2039)   ! Date: 06/01/10
  2040)   ! 
  2041) 
  2042)   use Option_module
  2043) 
  2044)   implicit none
  2045)   
  2046)   class(realization_subsurface_type) :: realization
  2047)   PetscInt :: global_total_count
  2048)   PetscInt :: global_active_count
  2049)   PetscInt :: total_count
  2050)   PetscInt :: active_count
  2051)   
  2052)   PetscInt :: patch_total_count
  2053)   PetscInt :: patch_active_count
  2054)   PetscInt :: temp_int_in(2), temp_int_out(2)
  2055)   PetscErrorCode :: ierr
  2056)   
  2057)   type(patch_type), pointer :: patch
  2058)   
  2059)   total_count = 0
  2060)   active_count = 0
  2061)     
  2062)   patch => realization%patch
  2063)   call PatchCountCells(patch,patch_total_count,patch_active_count)
  2064)   total_count = total_count + patch_total_count
  2065)   active_count = active_count + patch_active_count
  2066)   
  2067)   temp_int_in(1) = total_count
  2068)   temp_int_in(2) = active_count
  2069)   call MPI_Allreduce(temp_int_in,temp_int_out,TWO_INTEGER_MPI,MPIU_INTEGER, &
  2070)                      MPI_SUM,realization%option%mycomm,ierr)
  2071)   global_total_count = temp_int_out(1)
  2072)   global_active_count = temp_int_out(2)
  2073) 
  2074) end subroutine RealizationCountCells
  2075) 
  2076) ! ************************************************************************** !
  2077) 
  2078) subroutine RealizationPrintGridStatistics(realization)
  2079)   ! 
  2080)   ! Prints statistics regarding the numerical
  2081)   ! discretization
  2082)   ! 
  2083)   ! Author: Glenn Hammond
  2084)   ! Date: 06/01/10
  2085)   ! 
  2086) 
  2087)   use Grid_module
  2088) 
  2089)   implicit none
  2090) 
  2091)   class(realization_subsurface_type) :: realization
  2092)   
  2093)   type(option_type), pointer :: option
  2094)   type(grid_type), pointer :: grid
  2095) 
  2096)   PetscInt :: i1, i2, i3
  2097)   PetscReal :: r1, r2, r3
  2098)   PetscInt :: global_total_count, global_active_count
  2099)   PetscInt :: total_count, active_count
  2100)   PetscReal :: total_min, total_max, total_mean, total_variance
  2101)   PetscReal :: active_min, active_max, active_mean, active_variance
  2102)   PetscInt :: inactive_histogram(12), temp_int_out(12)
  2103)   PetscReal :: inactive_percentages(12)
  2104)   PetscErrorCode :: ierr
  2105) 
  2106)   option => realization%option
  2107)   grid => realization%patch%grid
  2108) 
  2109)   ! print # of active and inactive grid cells
  2110)   call RealizationCountCells(realization,global_total_count, &
  2111)                              global_active_count,total_count,active_count)
  2112)   r1 = dble(total_count)
  2113)   call OptionMaxMinMeanVariance(r1,total_max, &
  2114)                                 total_min,total_mean, &
  2115)                                 total_variance,PETSC_TRUE,option)
  2116)   r1 = dble(active_count)
  2117)   call OptionMaxMinMeanVariance(r1,active_max, &
  2118)                                 active_min,active_mean, &
  2119)                                 active_variance,PETSC_TRUE,option)
  2120)                   
  2121)   r1 = dble(active_count) / dble(total_count)    
  2122)   inactive_histogram = 0                          
  2123)   if (r1 >= (1.d0-1.d-8)) then
  2124)     inactive_histogram(12) = 1
  2125)   else if (r1 >= .9d0 .and. r1 < (1.d0-1.d-8)) then
  2126)     inactive_histogram(11) = 1
  2127)   else if (r1 >= .8d0 .and. r1 < .9d0) then
  2128)     inactive_histogram(10) = 1
  2129)   else if (r1 >= .7d0 .and. r1 < .8d0) then
  2130)     inactive_histogram(9) = 1
  2131)   else if (r1 >= .6d0 .and. r1 < .7d0) then
  2132)     inactive_histogram(8) = 1
  2133)   else if (r1 >= .5d0 .and. r1 < .6d0) then
  2134)     inactive_histogram(7) = 1
  2135)   else if (r1 >= .4d0 .and. r1 < .5d0) then
  2136)     inactive_histogram(6) = 1
  2137)   else if (r1 >= .3d0 .and. r1 < .4d0) then
  2138)     inactive_histogram(5) = 1
  2139)   else if (r1 >= .2d0 .and. r1 < .3d0) then
  2140)     inactive_histogram(4) = 1
  2141)   else if (r1 >= .1d0 .and. r1 < .2d0) then
  2142)     inactive_histogram(3) = 1
  2143)   else if (r1 > 1.d-20 .and. r1 < .1d0) then
  2144)     inactive_histogram(2) = 1
  2145)   else if (r1 < 1.d-20) then
  2146)     inactive_histogram(1) = 1
  2147)   endif
  2148)   
  2149)   call MPI_Allreduce(inactive_histogram,temp_int_out,TWELVE_INTEGER_MPI, &
  2150)                      MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
  2151) 
  2152)   ! why I cannot use *100, I do not know....geh
  2153)   inactive_percentages = dble(temp_int_out)/dble(option%mycommsize)*10.d0
  2154)   inactive_percentages = inactive_percentages+1.d-8
  2155) 
  2156)   r1 = 0.d0
  2157)   do i1 = 1, 12
  2158)     r1 = r1 + inactive_percentages(i1)
  2159)   enddo
  2160)                                 
  2161)   i1 = UNINITIALIZED_INTEGER
  2162)   i2 = UNINITIALIZED_INTEGER
  2163)   i3 = UNINITIALIZED_INTEGER
  2164)   if (associated(grid%structured_grid)) then
  2165)     i1 = grid%structured_grid%npx_final
  2166)     i2 = grid%structured_grid%npy_final
  2167)     i3 = grid%structured_grid%npz_final
  2168)   endif
  2169)   if (OptionPrintToScreen(option)) then
  2170)     write(*,'(/," Grid Stats:",/, &
  2171)               & "                       Global # cells: ",i12,/, &
  2172)               & "                Global # active cells: ",i12,/, &
  2173)               & "                              # cores: ",i12,/, &
  2174)               & "         Processor core decomposition: ",3i6,/, &
  2175)               & "               Maximum # cells / core: ",i12,/, &
  2176)               & "               Minimum # cells / core: ",i12,/, &
  2177)               & "               Average # cells / core: ",1pe12.4,/, &
  2178)               & "               Std Dev # cells / core: ",1pe12.4,/, &
  2179)               & "        Maximum # active cells / core: ",i12,/, &
  2180)               & "        Minimum # active cells / core: ",i12,/, &
  2181)               & "        Average # active cells / core: ",1pe12.4,/, &
  2182)               & "        Std Dev # active cells / core: ",1pe12.4,/,/, &
  2183)               & "        % cores with % active cells =       0%: ",1f7.2,/, &
  2184)               & "        % cores with % active cells =  0.1-10%: ",1f7.2,/, &
  2185)               & "        % cores with % active cells =   10-20%: ",1f7.2,/, &
  2186)               & "        % cores with % active cells =   20-30%: ",1f7.2,/, &
  2187)               & "        % cores with % active cells =   30-40%: ",1f7.2,/, &
  2188)               & "        % cores with % active cells =   40-50%: ",1f7.2,/, &
  2189)               & "        % cores with % active cells =   50-60%: ",1f7.2,/, &
  2190)               & "        % cores with % active cells =   60-70%: ",1f7.2,/, &
  2191)               & "        % cores with % active cells =   70-80%: ",1f7.2,/, &
  2192)               & "        % cores with % active cells =   80-90%: ",1f7.2,/, &
  2193)               & "        % cores with % active cells = 90-99.9%: ",1f7.2,/, &
  2194)               & "        % cores with % active cells =     100%: ",1f7.2,/, &
  2195)               & "                                        Check : ",1f7.2,/)') &
  2196)            global_total_count, &
  2197)            global_active_count, &
  2198)            option%mycommsize, &
  2199)            i1,i2,i3, &
  2200)            int(total_max+1.d-4), &
  2201)            int(total_min+1.d-4), &
  2202)            total_mean, sqrt(total_variance), &
  2203)            int(active_max+1.d-4), &
  2204)            int(active_min+1.d-4), &
  2205)            active_mean, sqrt(active_variance), &
  2206)            inactive_percentages(1), &
  2207)            inactive_percentages(2), &
  2208)            inactive_percentages(3), &
  2209)            inactive_percentages(4), &
  2210)            inactive_percentages(5), &
  2211)            inactive_percentages(6), &
  2212)            inactive_percentages(7), &
  2213)            inactive_percentages(8), &
  2214)            inactive_percentages(9), &
  2215)            inactive_percentages(10), &
  2216)            inactive_percentages(11), &
  2217)            inactive_percentages(12), &
  2218)            r1
  2219)   endif
  2220)   if (OptionPrintToFile(option)) then
  2221)     write(option%fid_out,'(/," Grid Stats:",/, &
  2222)                & "                       Global # cells: ",i12,/, &
  2223)                & "                Global # active cells: ",i12,/, &
  2224)                & "                              # cores: ",i12,/, &
  2225)                & "         Processor core decomposition: ",3i6,/, &
  2226)                & "               Maximum # cells / core: ",i12,/, &
  2227)                & "               Minimum # cells / core: ",i12,/, &
  2228)                & "               Average # cells / core: ",1pe12.4,/, &
  2229)                & "               Std Dev # cells / core: ",1pe12.4,/, &
  2230)                & "        Maximum # active cells / core: ",i12,/, &
  2231)                & "        Minimum # active cells / core: ",i12,/, &
  2232)                & "        Average # active cells / core: ",1pe12.4,/, &
  2233)                & "        Std Dev # active cells / core: ",1pe12.4,/,/, &
  2234)                & "        % cores with % active cells =       0%: ",1f7.2,/, &
  2235)                & "        % cores with % active cells =  0.1-10%: ",1f7.2,/, &
  2236)                & "        % cores with % active cells =   10-20%: ",1f7.2,/, &
  2237)                & "        % cores with % active cells =   20-30%: ",1f7.2,/, &
  2238)                & "        % cores with % active cells =   30-40%: ",1f7.2,/, &
  2239)                & "        % cores with % active cells =   40-50%: ",1f7.2,/, &
  2240)                & "        % cores with % active cells =   50-60%: ",1f7.2,/, &
  2241)                & "        % cores with % active cells =   60-70%: ",1f7.2,/, &
  2242)                & "        % cores with % active cells =   70-80%: ",1f7.2,/, &
  2243)                & "        % cores with % active cells =   80-90%: ",1f7.2,/, &
  2244)                & "        % cores with % active cells = 90-99.9%: ",1f7.2,/, &
  2245)                & "        % cores with % active cells =     100%: ",1f7.2,/, &
  2246)                & "                                        Check : ",1f7.2,/)') &
  2247)            global_total_count, &
  2248)            global_active_count, &
  2249)            option%mycommsize, &
  2250)            i1,i2,i3, &
  2251)            int(total_max+1.d-4), &
  2252)            int(total_min+1.d-4), &
  2253)            total_mean, sqrt(total_variance), &
  2254)            int(active_max+1.d-4), &
  2255)            int(active_min+1.d-4), &
  2256)            active_mean, sqrt(active_variance), &
  2257)            inactive_percentages(1), &
  2258)            inactive_percentages(2), &
  2259)            inactive_percentages(3), &
  2260)            inactive_percentages(4), &
  2261)            inactive_percentages(5), &
  2262)            inactive_percentages(6), &
  2263)            inactive_percentages(7), &
  2264)            inactive_percentages(8), &
  2265)            inactive_percentages(9), &
  2266)            inactive_percentages(10), &
  2267)            inactive_percentages(11), &
  2268)            inactive_percentages(12), &
  2269)            r1
  2270)   endif
  2271) 
  2272) end subroutine RealizationPrintGridStatistics
  2273) 
  2274) ! ************************************************************************** !
  2275) 
  2276) subroutine RealizationCalculateCFL1Timestep(realization,max_dt_cfl_1)
  2277)   ! 
  2278)   ! Calculates largest time step that
  2279)   ! preserves a CFL # of 1 in a realization
  2280)   ! 
  2281)   ! Author: Glenn Hammond
  2282)   ! Date: 10/07/11
  2283)   ! 
  2284) 
  2285)   implicit none
  2286) 
  2287)   class(realization_subsurface_type) realization
  2288)   PetscReal :: max_dt_cfl_1
  2289)   
  2290)   type(patch_type), pointer :: patch
  2291)   PetscReal :: max_dt_cfl_1_patch
  2292)   PetscReal :: tempreal
  2293)   PetscErrorCode :: ierr
  2294)   
  2295)   max_dt_cfl_1 = 1.d20
  2296)   patch => realization%patch
  2297)   call PatchCalculateCFL1Timestep(patch,realization%option, &
  2298)                                   max_dt_cfl_1_patch)
  2299)   max_dt_cfl_1 = min(max_dt_cfl_1,max_dt_cfl_1_patch)
  2300) 
  2301)   ! get the minimum across all cores
  2302)   call MPI_Allreduce(max_dt_cfl_1,tempreal,ONE_INTEGER_MPI, &
  2303)                      MPI_DOUBLE_PRECISION,MPI_MIN, &
  2304)                      realization%option%mycomm,ierr)
  2305)   max_dt_cfl_1 = tempreal
  2306) 
  2307) end subroutine RealizationCalculateCFL1Timestep
  2308) 
  2309) ! ************************************************************************** !
  2310) 
  2311) subroutine RealizUnInitializedVarsFlow(realization)
  2312)   ! 
  2313)   ! Checks for uninitialized flow variables
  2314)   ! 
  2315)   ! Author: Glenn Hammond
  2316)   ! Date: 07/06/16
  2317)   ! 
  2318)   use Option_module
  2319)   use Material_Aux_class
  2320)   use Variables_module, only : VOLUME, MINERAL_POROSITY, PERMEABILITY_X, &
  2321)                                PERMEABILITY_Y, PERMEABILITY_Z
  2322) 
  2323)   implicit none
  2324)   
  2325)   class(realization_subsurface_type) :: realization
  2326) 
  2327)   character(len=MAXWORDLENGTH) :: var_name
  2328)   PetscInt :: i
  2329) 
  2330)   call RealizUnInitializedVar1(realization,VOLUME,'volume')
  2331)   ! mineral porosity is the base, unmodified porosity
  2332)   call RealizUnInitializedVar1(realization,MINERAL_POROSITY,'porosity')
  2333)   call RealizUnInitializedVar1(realization,PERMEABILITY_X,'permeability X')
  2334)   call RealizUnInitializedVar1(realization,PERMEABILITY_Y,'permeability Y')
  2335)   call RealizUnInitializedVar1(realization,PERMEABILITY_Z,'permeability Z')
  2336)   do i = 1, max_material_index
  2337)     var_name = MaterialAuxIndexToPropertyName(i)
  2338)     call RealizUnInitializedVar1(realization,i,var_name)
  2339)   enddo
  2340) 
  2341) end subroutine RealizUnInitializedVarsFlow
  2342) 
  2343) ! ************************************************************************** !
  2344) 
  2345) subroutine RealizUnInitializedVarsTran(realization)
  2346)   ! 
  2347)   ! Checks for uninitialized transport variables
  2348)   ! 
  2349)   ! Author: Glenn Hammond
  2350)   ! Date: 07/06/16
  2351)   ! 
  2352) 
  2353)   use Grid_module
  2354)   use Patch_module
  2355)   use Option_module
  2356)   use Material_module
  2357)   use Material_Aux_class
  2358)   use Variables_module, only : VOLUME, MINERAL_POROSITY, TORTUOSITY
  2359) 
  2360)   implicit none
  2361)   
  2362)   class(realization_subsurface_type) :: realization
  2363) 
  2364)   call RealizUnInitializedVar1(realization,VOLUME,'volume')
  2365)   ! mineral porosity is the base, unmodified porosity
  2366)   call RealizUnInitializedVar1(realization,MINERAL_POROSITY,'porosity')
  2367)   call RealizUnInitializedVar1(realization,TORTUOSITY,'tortuosity')
  2368) 
  2369) end subroutine RealizUnInitializedVarsTran
  2370) 
  2371) ! ************************************************************************** !
  2372) 
  2373) subroutine RealizUnInitializedVar1(realization,ivar,var_name)
  2374)   ! 
  2375)   ! Checks whether a variable is initialized at all active grid cells
  2376)   ! 
  2377)   ! Author: Glenn Hammond
  2378)   ! Date: 07/06/16
  2379)   ! 
  2380) 
  2381)   use Option_module
  2382)   use Field_module
  2383)   use Patch_module
  2384)   use Grid_module
  2385) 
  2386)   implicit none
  2387)   
  2388)   class(realization_subsurface_type) :: realization
  2389)   PetscInt :: ivar
  2390)   character(len=*) :: var_name
  2391)   
  2392)   type(option_type), pointer :: option
  2393)   type(field_type), pointer :: field
  2394)   type(grid_type), pointer :: grid
  2395)   type(patch_type), pointer :: patch
  2396)   PetscReal, pointer :: vec_p(:)
  2397)   PetscInt :: local_id
  2398)   PetscInt :: imin
  2399)   PetscReal :: rmin
  2400)   character(len=MAXWORDLENGTH) :: word
  2401)   PetscErrorCode :: ierr
  2402) 
  2403)   option => realization%option
  2404)   field => realization%field
  2405)   patch => realization%patch
  2406)   grid => patch%grid
  2407) 
  2408)   call RealizationGetVariable(realization,realization%field%work, &
  2409)                               ivar,ZERO_INTEGER)
  2410)   ! apply mask to filter inactive cells
  2411)   call VecGetArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  2412)   do local_id = 1, grid%nlmax
  2413)     ! if inactive, set to 1.d-40
  2414)     if (patch%imat(grid%nL2G(local_id)) <= 0) vec_p(local_id) = 1.d-40
  2415)   enddo
  2416)   call VecRestoreArrayF90(field%work,vec_p,ierr);CHKERRQ(ierr)
  2417)   call VecMin(field%work,imin,rmin,ierr);CHKERRQ(ierr)
  2418)   if (Uninitialized(rmin)) then
  2419)     write(word,*) imin+1 ! zero to one based indexing
  2420)     option%io_buffer = 'Incorrect assignment of variable (' &
  2421)       // trim(var_name) // ',cell=' // trim(adjustl(word)) // &
  2422)       '). Please send this error message and your input file to &
  2423)       &pflotran-dev@googlegroups.com.'
  2424)     call printErrMsg(option)
  2425)   endif
  2426) 
  2427) end subroutine RealizUnInitializedVar1
  2428) 
  2429) ! ************************************************************************** !
  2430) 
  2431) subroutine RealizSetSoilReferencePressure(realization)
  2432)   ! 
  2433)   ! Deallocates a realization
  2434)   ! 
  2435)   ! Author: Glenn Hammond
  2436)   ! Date: 07/06/16
  2437)   ! 
  2438)   use Patch_module
  2439)   use Grid_module
  2440)   use Material_Aux_class
  2441)   use Fracture_module
  2442)   use Variables_module, only : MAXIMUM_PRESSURE, SOIL_REFERENCE_PRESSURE
  2443) 
  2444)   implicit none
  2445) 
  2446)   type(realization_subsurface_type) :: realization
  2447) 
  2448)   type(patch_type), pointer :: patch
  2449)   type(grid_type), pointer :: grid
  2450)   class(material_auxvar_type), pointer :: material_auxvars(:)
  2451)   type(material_type), pointer :: Material
  2452)   type(material_property_ptr_type), pointer :: material_property_array(:)
  2453)   PetscReal, pointer :: vec_loc_p(:)
  2454) 
  2455)   PetscInt :: ghosted_id
  2456)   PetscInt :: imat
  2457)   PetscErrorCode :: ierr
  2458) 
  2459)   patch => realization%patch
  2460)   grid => patch%grid
  2461)   material_property_array => patch%material_property_array
  2462)   material_auxvars => patch%aux%Material%auxvars
  2463) 
  2464)   call RealizationGetVariable(realization,realization%field%work, &
  2465)                               MAXIMUM_PRESSURE,ZERO_INTEGER)
  2466)   call DiscretizationGlobalToLocal(realization%discretization, &
  2467)                                    realization%field%work, &
  2468)                                    realization%field%work_loc, &
  2469)                                    ONEDOF)
  2470)   call VecGetArrayReadF90(realization%field%work_loc,vec_loc_p, &
  2471)                           ierr); CHKERRQ(ierr)
  2472) 
  2473)   do ghosted_id = 1, grid%ngmax
  2474)     imat = patch%imat(ghosted_id)
  2475)     if (imat <= 0) cycle
  2476)     if (associated(material_auxvars(ghosted_id)%fracture)) then
  2477)       call FractureSetInitialPressure(material_auxvars(ghosted_id)%fracture, &
  2478)                                       vec_loc_p(ghosted_id))
  2479)     endif
  2480)     if (material_property_array(imat)%ptr%soil_reference_pressure_initial) then
  2481)       call MaterialAuxVarSetValue(material_auxvars(ghosted_id), &
  2482)                                   SOIL_REFERENCE_PRESSURE, &
  2483)                                   vec_loc_p(ghosted_id))
  2484)     endif
  2485)   enddo
  2486) 
  2487)   call VecRestoreArrayReadF90(realization%field%work_loc,vec_loc_p, &
  2488)                               ierr); CHKERRQ(ierr)
  2489) 
  2490) end subroutine RealizSetSoilReferencePressure
  2491) 
  2492) ! ************************************************************************** !
  2493) 
  2494) subroutine RealizationDestroyLegacy(realization)
  2495)   ! 
  2496)   ! Deallocates a realization
  2497)   ! 
  2498)   ! Author: Glenn Hammond
  2499)   ! Date: 11/01/07
  2500)   ! 
  2501) 
  2502)   use Dataset_module
  2503) 
  2504)   implicit none
  2505)   
  2506)   class(realization_subsurface_type), pointer :: realization
  2507)   
  2508)   if (.not.associated(realization)) return
  2509)     
  2510)   call FieldDestroy(realization%field)
  2511) 
  2512) !  call OptionDestroy(realization%option) !geh it will be destroy externally
  2513)   call OutputOptionDestroy(realization%output_option)
  2514)   call RegionDestroyList(realization%region_list)
  2515)   
  2516)   call FlowConditionDestroyList(realization%flow_conditions)
  2517)   call TranConditionDestroyList(realization%transport_conditions)
  2518)   call TranConstraintDestroyList(realization%transport_constraints)
  2519) 
  2520)   call PatchDestroyList(realization%patch_list)
  2521) 
  2522)   if (associated(realization%debug)) deallocate(realization%debug)
  2523)   nullify(realization%debug)
  2524)   
  2525)   if (associated(realization%fluid_property_array)) &
  2526)     deallocate(realization%fluid_property_array)
  2527)   nullify(realization%fluid_property_array)
  2528)   call FluidPropertyDestroy(realization%fluid_properties)
  2529)   
  2530)   call MaterialPropertyDestroy(realization%material_properties)
  2531) 
  2532)   call SaturationFunctionDestroy(realization%saturation_functions)
  2533)   print *, 'RealizationDestroyLegacy cannot be removed.'
  2534)   stop
  2535)   call CharacteristicCurvesDestroy(realization%characteristic_curves)
  2536) 
  2537)   call DatasetDestroy(realization%datasets)
  2538)   
  2539)   call UniformVelocityDatasetDestroy(realization%uniform_velocity_dataset)
  2540)   
  2541)   call DiscretizationDestroy(realization%discretization)
  2542)   
  2543)   call ReactionDestroy(realization%reaction,realization%option)
  2544)   
  2545)   call TranConstraintDestroy(realization%sec_transport_constraint)
  2546)   
  2547)   deallocate(realization)
  2548)   nullify(realization)
  2549)   
  2550) end subroutine RealizationDestroyLegacy
  2551) 
  2552) ! ************************************************************************** !
  2553) 
  2554) subroutine RealizationStrip(this)
  2555)   ! 
  2556)   ! Deallocates a realization
  2557)   ! 
  2558)   ! Author: Glenn Hammond
  2559)   ! Date: 11/01/07
  2560)   ! 
  2561) 
  2562)   use Dataset_module
  2563) 
  2564)   implicit none
  2565)   
  2566)   class(realization_subsurface_type) :: this
  2567)   
  2568)   call RealizationBaseStrip(this)
  2569)   call RegionDestroyList(this%region_list)
  2570)   
  2571)   call FlowConditionDestroyList(this%flow_conditions)
  2572)   call TranConditionDestroyList(this%transport_conditions)
  2573)   call TranConstraintDestroyList(this%transport_constraints)
  2574) 
  2575)   if (associated(this%fluid_property_array)) &
  2576)     deallocate(this%fluid_property_array)
  2577)   nullify(this%fluid_property_array)
  2578)   call FluidPropertyDestroy(this%fluid_properties)
  2579)   
  2580)   call MaterialPropertyDestroy(this%material_properties)
  2581) 
  2582)   call SaturationFunctionDestroy(this%saturation_functions)
  2583)   call CharacteristicCurvesDestroy(this%characteristic_curves)  
  2584) 
  2585)   call DatasetDestroy(this%datasets)
  2586)   
  2587)   call UniformVelocityDatasetDestroy(this%uniform_velocity_dataset)
  2588)   
  2589)   call ReactionDestroy(this%reaction,this%option)
  2590)   
  2591)   call TranConstraintDestroy(this%sec_transport_constraint)
  2592)   
  2593) end subroutine RealizationStrip
  2594) 
  2595) end module Realization_Subsurface_class

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