MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 :: gweshared_dat_df
11  public :: set_gwe_dat_ptrs
12 
13  !> Data for sharing among multiple packages. Originally read in from
14  !< the EST package
15 
17 
18  ! dim
19  integer(I4B) :: nnodes !< number of cells
20 
21  ! strings
22  character(len=LENMEMPATH) :: memorypath = '' !< the location in the memory manager where the variables are stored
23 
24  ! est data to be share across multiple packages
25  real(dp), pointer :: gwerhow => null() !< Density of water (for GWE purposes, a constant scalar)
26  real(dp), pointer :: gwecpw => null() !< Heat capacity of water (non-spatially varying)
27  real(dp), pointer :: gwelatheatvap => null() !< latent heat of vaporization
28  real(dp), dimension(:), pointer, contiguous :: gwerhos => null() !< Density of the aquifer material
29  real(dp), dimension(:), pointer, contiguous :: gwecps => null() !< Heat capacity of solids (spatially varying)
30 
31  contains
32 
33  ! -- public
34  procedure, public :: gweshared_dat_df
35  procedure, public :: set_gwe_dat_ptrs
36  procedure, public :: gweshared_dat_da
37  ! -- private
38  procedure, private :: allocate_shared_vars
39  procedure, private :: set_gwe_shared_scalars
40  procedure, private :: set_gwe_shared_arrays
41 
42  end type gweinputdatatype
43 
44 contains
45 
46 !> @brief Allocate the shared data
47 !<
48  subroutine gweshared_dat_cr(gwe_dat)
49  ! -- modules
50  ! -- dummy
51  type(gweinputdatatype), pointer :: gwe_dat !< the input data block
52  !
53  ! -- Create the object
54  allocate (gwe_dat)
55  !
56  ! -- Return
57  return
58  end subroutine gweshared_dat_cr
59 
60 !> @brief Define the shared data
61 !<
62  subroutine gweshared_dat_df(this, nodes)
63  ! -- dummy
64  class(gweinputdatatype) :: this !< the input data block
65  integer(I4B), intent(in) :: nodes
66  ! -- local
67  !
68  ! -- Allocate variables
69  call this%allocate_shared_vars(nodes)
70  !
71  ! -- Return
72  return
73  end subroutine gweshared_dat_df
74 
75  !> @brief Define the information this object holds
76  !!
77  !! Allocate strings for storing label names
78  !! Intended to be analogous to allocate_scalars()
79  !<
80  subroutine allocate_shared_vars(this, nodes)
81  ! -- dummy
82  class(gweinputdatatype) :: this !< TspLabelsType object
83  integer(I4B), intent(in) :: nodes
84  ! -- local
85  integer(I4B) :: i
86  !
87  allocate (this%gwecpw)
88  allocate (this%gwerhow)
89  allocate (this%gwelatheatvap)
90  allocate (this%gwerhos(nodes))
91  allocate (this%gwecps(nodes))
92  !
93  ! -- Initialize values
94  this%gwecpw = dzero
95  this%gwerhow = dzero
96  this%gwelatheatvap = dzero
97  do i = 1, nodes
98  this%gwecps(i) = dzero
99  this%gwerhos(i) = dzero
100  end do
101  !
102  ! -- Return
103  return
104  end subroutine allocate_shared_vars
105 
106  !> @brief Allocate and read data from EST
107  !!
108  !! EST data, including heat capacity of water (cpw), density of water
109  !! (rhow), latent heat of vaporization (latheatvap), heat capacity of
110  !! the aquifer material (cps), and density of the aquifer material
111  !! (rhow) is used among other packages and is therefore stored in a
112  !! separate class
113  subroutine set_gwe_dat_ptrs(this, gwerhow, gwecpw, gwerhos, gwecps, &
114  gwelatheatvap)
115  ! -- dummy
116  class(gweinputdatatype) :: this !< the input data block
117  real(dp), intent(in) :: gwerhow !< ptr to density of water specified in EST
118  real(dp), intent(in) :: gwecpw !< ptr to heat capacity of water specified in EST
119  real(dp), intent(in) :: gwerhos !< ptr to sptially-variably density of aquifer material specified in EST
120  real(dp), intent(in) :: gwecps !< ptr to sptially-variably heat capacity of aquifer material specified in EST
121  real(dp), intent(in), optional :: gwelatheatvap !< ptr to latent heat of vaporization specified in EST
122  !
123  ! -- Allocate scalars
124  if (present(gwelatheatvap)) then
125  call this%set_gwe_shared_scalars(gwerhow, gwecpw, gwelatheatvap)
126  else
127  call this%set_gwe_shared_scalars(gwerhow, gwecpw)
128  end if
129  !
130  ! -- Allocate arrays
131  call this%set_gwe_shared_arrays(gwerhos, gwecps)
132  !
133  ! -- Return
134  return
135  end subroutine set_gwe_dat_ptrs
136 
137  !> @brief Set pointers to scalars read by the EST package
138  !! for use by other packages
139  !!
140  !! Set pointers to GWE-related scalars and arrays for use
141  !! by multiple packages. For example, a package capable of
142  !! simulating evaporation will need access to latent heat of
143  !! of vaporization.
144  !<
145  subroutine set_gwe_shared_scalars(this, gwerhow, gwecpw, gwelatheatvap)
146  ! -- dummy
147  class(gweinputdatatype) :: this !< GweInputDataType object
148  real(DP), intent(in) :: gwerhow
149  real(DP), intent(in) :: gwecpw
150  real(DP), intent(in), optional :: gwelatheatvap
151  ! -- local
152  !
153  ! -- Set the pointers
154  ! -- Fixed density of water to be used by GWE
155  this%gwerhow = gwerhow
156  ! -- Spatially constant heat capacity of water ! kluge note: "specific heat" (which is heat capacity per unit mass) is probably the more correct term
157  this%gwecpw = gwecpw
158  ! -- Latent heat of vaporization
159  if (present(gwelatheatvap)) then
160  this%gwelatheatvap = gwelatheatvap
161  end if
162  !
163  ! -- Return
164  return
165  end subroutine set_gwe_shared_scalars
166 
167  !> @brief Set pointers to data arrays read by the EST package
168  !! for use by other packages
169  !!
170  !! Set pointers to GWE-related arrays for use
171  !! by multiple packages.
172  !<
173  subroutine set_gwe_shared_arrays(this, gwerhos, gwecps)
174  ! -- dummy
175  class(gweinputdatatype) :: this !< GweInputDataType object
176  real(DP), intent(in) :: gwerhos
177  real(DP), intent(in) :: gwecps
178  ! -- local
179  !
180  ! -- Set the pointers
181  ! -- Spatially-variable density of aquifer solid material
182  this%gwerhos = gwerhos
183  ! -- Spatially-variable heat capacity of aquifer solid material
184  this%gwecps = gwecps
185  !
186  ! -- Return
187  return
188  end subroutine set_gwe_shared_arrays
189 
190  !> @ brief Deallocate memory
191  !!
192  !! Deallocate GWE shared data array memory
193  !<
194  subroutine gweshared_dat_da(this)
195  ! -- dummy
196  class(gweinputdatatype) :: this !< the input data block
197  !
198  ! -- Scalars
199  deallocate (this%gwelatheatvap)
200  deallocate (this%gwerhow)
201  deallocate (this%gwecpw)
202  !
203  ! -- Arrays
204  deallocate (this%gwerhos)
205  deallocate (this%gwecps)
206  !
207  ! -- Return
208  return
209  end subroutine gweshared_dat_da
210 
211 end module gweinputdatamodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
subroutine, public set_gwe_dat_ptrs(this, gwerhow, gwecpw, gwerhos, gwecps, gwelatheatvap)
Allocate and read data from EST.
subroutine set_gwe_shared_scalars(this, gwerhow, gwecpw, gwelatheatvap)
Set pointers to scalars read by the EST package for use by other packages.
subroutine, public gweshared_dat_cr(gwe_dat)
Allocate the shared data.
subroutine gweshared_dat_da(this)
@ brief Deallocate memory
subroutine set_gwe_shared_arrays(this, gwerhos, gwecps)
Set pointers to data arrays read by the EST package for use by other packages.
subroutine, public gweshared_dat_df(this, nodes)
Define the shared data.
subroutine allocate_shared_vars(this, nodes)
Define the information this object holds.
This module defines variable data types.
Definition: kind.f90:8
Data for sharing among multiple packages. Originally read in from.