MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
GweInputData.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: i4b, dp
5 
6  implicit none
7  private
8  public :: gweinputdatatype
9  public :: gweshared_dat_cr
10  public :: set_gwe_dat_ptrs
11 
12  !> Data for sharing among multiple packages. Originally read in from
13  !< the EST package
14 
16 
17  ! dim
18  integer(I4B) :: nnodes !< number of cells
19 
20  ! strings
21  character(len=LENMEMPATH) :: memorypath = '' !< the location in the memory manager where the variables are stored
22 
23  ! est data to be share across multiple packages
24  real(dp), pointer :: gwerhow => null() !< Density of water (for GWE purposes, a constant scalar)
25  real(dp), pointer :: gwecpw => null() !< Heat capacity of water (non-spatially varying)
26  real(dp), pointer :: gwelatheatvap => null() !< latent heat of vaporization
27  real(dp), dimension(:), pointer, contiguous :: gwerhos => null() !< Density of the aquifer material
28  real(dp), dimension(:), pointer, contiguous :: gwecps => null() !< Heat capacity of solids (spatially varying)
29 
30  contains
31 
32  ! public
33  procedure, public :: set_gwe_dat_ptrs
34  procedure, public :: gweshared_dat_da
35  ! private
36  procedure, private :: set_gwe_scalar_ptrs
37  procedure, private :: set_gwe_array_ptrs
38 
39  end type gweinputdatatype
40 
41 contains
42 
43  !> @brief Allocate the shared data
44  !<
45  subroutine gweshared_dat_cr(gwe_dat)
46  ! modules
47  ! dummy
48  type(gweinputdatatype), pointer :: gwe_dat !< the input data block
49 
50  ! create the object
51  allocate (gwe_dat)
52 
53  end subroutine gweshared_dat_cr
54 
55  !> @brief Allocate and read data from EST
56  !!
57  !! EST data, including heat capacity of water (cpw), density of water
58  !! (rhow), latent heat of vaporization (latheatvap), heat capacity of
59  !! the aquifer material (cps), and density of the aquifer material
60  !! (rhow) is used among other packages and is therefore stored in a
61  !! separate class
62  !<
63  subroutine set_gwe_dat_ptrs(this, gwerhow, gwecpw, gwelatheatvap, &
64  gwerhos, gwecps)
65  ! dummy
66  class(gweinputdatatype) :: this !< the input data block
67  real(dp), intent(in), pointer :: gwerhow !< ptr to density of water specified in EST
68  real(dp), intent(in), pointer :: gwecpw !< ptr to heat capacity of water specified in EST
69  real(dp), intent(in), pointer :: gwelatheatvap !< ptr to latent heat of vaporization specified in EST
70  real(dp), dimension(:), pointer, contiguous :: gwerhos !< ptr to sptially-variably density of aquifer material specified in EST
71  real(dp), dimension(:), pointer, contiguous :: gwecps !< ptr to sptially-variably heat capacity of aquifer material specified in EST
72 
73  ! allocate scalars
74  call this%set_gwe_scalar_ptrs(gwerhow, gwecpw, gwelatheatvap)
75 
76  ! allocate arrays
77  call this%set_gwe_array_ptrs(gwerhos, gwecps)
78 
79  end subroutine set_gwe_dat_ptrs
80 
81  !> @brief Set pointers to scalars read by the EST package
82  !! for use by other packages
83  !!
84  !! Set pointers to GWE-related scalars and arrays for use
85  !! by multiple packages. For example, a package capable of
86  !! simulating evaporation will need access to latent heat of
87  !! of vaporization.
88  !<
89  subroutine set_gwe_scalar_ptrs(this, gwerhow, gwecpw, gwelatheatvap)
90  ! dummy
91  class(gweinputdatatype) :: this !< GweInputDataType object
92  real(DP), pointer, intent(in) :: gwerhow !< density of water
93  real(DP), pointer, intent(in) :: gwecpw !< mass-based heat capacity of water
94  real(DP), pointer, intent(in), optional :: gwelatheatvap !< latent heat of vaporization
95 
96  ! set the pointers
97  ! fixed density of water to be used by GWE (vs GWT-based density)
98  this%gwerhow => gwerhow
99  ! spatially constant heat capacity of water ! kluge note: "specific heat" (which is heat capacity per unit mass) is probably the more correct term
100  this%gwecpw => gwecpw
101  ! latent heat of vaporization
102  if (present(gwelatheatvap)) then
103  this%gwelatheatvap => gwelatheatvap
104  end if
105 
106  end subroutine set_gwe_scalar_ptrs
107 
108  !> @brief Set pointers to data arrays read by the EST package
109  !! for use by other packages
110  !!
111  !! Set pointers to GWE-related arrays for use by multiple packages
112  !<
113  subroutine set_gwe_array_ptrs(this, gwerhos, gwecps)
114  ! dummy
115  class(gweinputdatatype) :: this !< GweInputDataType object
116  real(DP), dimension(:), pointer, contiguous, intent(in) :: gwerhos
117  real(DP), dimension(:), pointer, contiguous, intent(in) :: gwecps
118 
119  ! set the pointers
120  ! spatially-variable density of aquifer solid material
121  this%gwerhos => gwerhos
122  ! spatially-variable heat capacity of aquifer solid material
123  this%gwecps => gwecps
124 
125  end subroutine set_gwe_array_ptrs
126 
127  !> @ brief Deallocate memory
128  !!
129  !! Set pointers to null
130  !<
131  subroutine gweshared_dat_da(this)
132  ! dummy
133  class(gweinputdatatype) :: this !< the input data block
134 
135  ! scalars
136  this%gwelatheatvap => null()
137  this%gwerhow => null()
138  this%gwecpw => null()
139 
140  ! arrays
141  this%gwerhos => null()
142  this%gwecps => null()
143 
144  end subroutine gweshared_dat_da
145 
146 end module gweinputdatamodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dep3
real constant 1000
Definition: Constants.f90:88
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
subroutine, public set_gwe_dat_ptrs(this, gwerhow, gwecpw, gwelatheatvap, gwerhos, gwecps)
Allocate and read data from EST.
subroutine, public gweshared_dat_cr(gwe_dat)
Allocate the shared data.
subroutine set_gwe_array_ptrs(this, gwerhos, gwecps)
Set pointers to data arrays read by the EST package for use by other packages.
subroutine set_gwe_scalar_ptrs(this, gwerhow, gwecpw, gwelatheatvap)
Set pointers to scalars read by the EST package for use by other packages.
subroutine gweshared_dat_da(this)
@ brief Deallocate memory
This module defines variable data types.
Definition: kind.f90:8
Data for sharing among multiple packages. Originally read in from.