MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
gwe-esl.f90
Go to the documentation of this file.
2  !
3  use kindmodule, only: dp, i4b
5  use bndextmodule, only: bndexttype
9  !
10  implicit none
11  !
12  private
13  public :: esl_create
14  !
15  character(len=LENFTYPE) :: ftype = 'ESL'
16  character(len=16) :: text = ' ESL'
17  !
18  type, extends(bndexttype) :: gweesltype
19 
20  type(gweinputdatatype), pointer :: gwecommon => null() !< pointer to shared gwe data used by multiple packages but set in mst
21  real(dp), dimension(:), pointer, contiguous :: senerrate => null() !< energy source loading rate
22 
23  contains
24 
25  procedure :: allocate_scalars => esl_allocate_scalars
26  procedure :: allocate_arrays => esl_allocate_arrays
27  procedure :: bnd_cf => esl_cf
28  procedure :: bnd_fc => esl_fc
29  procedure :: bnd_da => esl_da
30  procedure :: define_listlabel
31  procedure :: bound_value => esl_bound_value
32  ! -- methods for observations
33  procedure, public :: bnd_obs_supported => esl_obs_supported
34  procedure, public :: bnd_df_obs => esl_df_obs
35 
36  end type gweesltype
37 
38 contains
39 
40  !> @brief Create an energy source loading package
41  !!
42  !! This subroutine points bndobj to the newly created package
43  !<
44  subroutine esl_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
45  gwecommon, input_mempath)
46  ! -- modules
47  use bndmodule, only: bndtype
48  ! -- dummy
49  class(bndtype), pointer :: packobj
50  integer(I4B), intent(in) :: id
51  integer(I4B), intent(in) :: ibcnum
52  integer(I4B), intent(in) :: inunit
53  integer(I4B), intent(in) :: iout
54  character(len=*), intent(in) :: namemodel
55  character(len=*), intent(in) :: pakname
56  type(gweinputdatatype), intent(in), target :: gwecommon !< shared data container for use by multiple GWE packages
57  character(len=*), intent(in) :: input_mempath
58  ! -- local
59  type(gweesltype), pointer :: eslobj
60  !
61  ! -- Allocate the object and assign values to object variables
62  allocate (eslobj)
63  packobj => eslobj
64  !
65  ! -- Create name and memory path
66  call packobj%set_names(ibcnum, namemodel, pakname, ftype, input_mempath)
67  packobj%text = text
68  !
69  ! -- Allocate scalars
70  call eslobj%allocate_scalars()
71  !
72  ! -- Initialize package
73  call packobj%pack_initialize()
74  !
75  packobj%inunit = inunit
76  packobj%iout = iout
77  packobj%id = id
78  packobj%ibcnum = ibcnum
79  packobj%ncolbnd = 1
80  packobj%iscloc = 1
81  !
82  ! -- Store pointer to shared data module for accessing cpw, rhow
83  ! for the budget calculations, and for accessing the latent heat of
84  ! vaporization for evaporative cooling.
85  eslobj%gwecommon => gwecommon
86  end subroutine esl_create
87 
88  !> @brief Deallocate memory
89  !<
90  subroutine esl_da(this)
91  ! -- modules
93  ! -- dummy
94  class(gweesltype) :: this
95  !
96  ! -- Deallocate parent package
97  call this%BndExtType%bnd_da()
98  end subroutine esl_da
99 
100  !> @brief Allocate scalars
101  !!
102  !! Allocate scalars specific to this energy source loading package
103  !<
104  subroutine esl_allocate_scalars(this)
105  ! -- modules
107  ! -- dummy
108  class(gweesltype) :: this
109  !
110  ! -- base class allocate scalars
111  call this%BndExtType%allocate_scalars()
112  !
113  ! -- allocate the object and assign values to object variables
114  !
115  ! -- Set values
116  end subroutine esl_allocate_scalars
117 
118  !> @brief Allocate arrays
119  !<
120  subroutine esl_allocate_arrays(this, nodelist, auxvar)
122  ! -- dummy
123  class(gweesltype) :: this
124  ! -- local
125  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist !< package nodelist
126  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar !< package aux variable array
127 
128  ! -- base class allocate arrays
129  call this%BndExtType%allocate_arrays(nodelist, auxvar)
130 
131  ! -- set input context pointers
132  call mem_setptr(this%senerrate, 'SENERRATE', this%input_mempath)
133  !
134  ! -- checkin input context pointers
135  call mem_checkin(this%senerrate, 'SENERRATE', this%memoryPath, &
136  'SENERRATE', this%input_mempath)
137  end subroutine esl_allocate_arrays
138 
139  !> @brief Formulate the HCOF and RHS terms
140  !!
141  !! This subroutine:
142  !! - calculates hcof and rhs terms
143  !! - skip if no sources
144  !<
145  subroutine esl_cf(this)
146  ! -- dummy
147  class(gweesltype) :: this
148  ! -- local
149  integer(I4B) :: i, node
150  real(DP) :: q
151  !
152  ! -- Return if no sources
153  if (this%nbound == 0) return
154  !
155  ! -- Calculate hcof and rhs for each source entry
156  do i = 1, this%nbound
157  node = this%nodelist(i)
158  this%hcof(i) = dzero
159  if (this%ibound(node) <= 0) then
160  this%rhs(i) = dzero
161  cycle
162  end if
163  q = this%bound_value(1, i)
164  this%rhs(i) = -q
165  end do
166  end subroutine esl_cf
167 
168  !> @brief Add matrix terms related to specified energy source loading
169  !!
170  !! Copy rhs and hcof into solution rhs and amat
171  !<
172  subroutine esl_fc(this, rhs, ia, idxglo, matrix_sln)
173  ! -- dummy
174  class(gweesltype) :: this
175  real(DP), dimension(:), intent(inout) :: rhs
176  integer(I4B), dimension(:), intent(in) :: ia
177  integer(I4B), dimension(:), intent(in) :: idxglo
178  class(matrixbasetype), pointer :: matrix_sln
179  ! -- local
180  integer(I4B) :: i, n, ipos
181  !
182  ! -- pakmvrobj fc
183  if (this%imover == 1) then
184  call this%pakmvrobj%fc()
185  end if
186  !
187  ! -- Copy package rhs and hcof into solution rhs and amat
188  do i = 1, this%nbound
189  n = this%nodelist(i)
190  rhs(n) = rhs(n) + this%rhs(i)
191  ipos = ia(n)
192  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
193  !
194  ! -- If mover is active and mass is being withdrawn,
195  ! store available mass (as positive value).
196  if (this%imover == 1 .and. this%rhs(i) > dzero) then
197  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
198  end if
199  end do
200  end subroutine esl_fc
201 
202  !> @brief Define list labels
203  !!
204  !! Define the list heading that is written to iout when
205  !! PRINT_INPUT option is used.
206  !<
207  subroutine define_listlabel(this)
208  ! -- dummy
209  class(gweesltype), intent(inout) :: this
210  !
211  ! -- Create the header list label
212  this%listlabel = trim(this%filtyp)//' NO.'
213  if (this%dis%ndim == 3) then
214  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
215  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
216  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
217  elseif (this%dis%ndim == 2) then
218  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
219  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
220  else
221  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
222  end if
223  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
224  if (this%inamedbound == 1) then
225  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
226  end if
227  end subroutine define_listlabel
228 
229  ! -- Procedures related to observations
230 
231  !> @brief Support function for specified energy source loading observations
232  !!
233  !! This function:
234  !! - returns true because ESL package supports observations.
235  !! - overrides BndType%bnd_obs_supported()
236  !<
237  logical function esl_obs_supported(this)
238  implicit none
239  ! -- dummy
240  class(gweesltype) :: this
241  !
242  esl_obs_supported = .true.
243  end function esl_obs_supported
244 
245  !> @brief Define observations
246  !!
247  !! This subroutine:
248  !! - stores observation types supported by ESL package.
249  !! - overrides BndType%bnd_df_obs
250  !<
251  subroutine esl_df_obs(this)
252  implicit none
253  ! -- dummy
254  class(gweesltype) :: this
255  ! -- local
256  integer(I4B) :: indx
257  !
258  call this%obs%StoreObsType('esl', .true., indx)
259  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
260  !
261  ! -- Store obs type and assign procedure pointer
262  ! for to-mvr observation type.
263  call this%obs%StoreObsType('to-mvr', .true., indx)
264  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
265  end subroutine esl_df_obs
266 
267  !> @ brief Return a bound value
268  !!
269  !! Return a bound value associated with an ncolbnd index
270  !! and row.
271  !<
272  function esl_bound_value(this, col, row) result(bndval)
273  ! -- modules
274  use constantsmodule, only: dzero
275  ! -- dummy variables
276  class(gweesltype), intent(inout) :: this
277  integer(I4B), intent(in) :: col
278  integer(I4B), intent(in) :: row
279  ! -- result
280  real(dp) :: bndval
281  !
282  select case (col)
283  case (1)
284  bndval = this%senerrate(row)
285  case default
286  end select
287  end function esl_bound_value
288 
289 end module gweeslmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:103
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public esl_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, gwecommon, input_mempath)
Create an energy source loading package.
Definition: gwe-esl.f90:46
real(dp) function esl_bound_value(this, col, row)
@ brief Return a bound value
Definition: gwe-esl.f90:273
subroutine esl_da(this)
Deallocate memory.
Definition: gwe-esl.f90:91
subroutine esl_allocate_arrays(this, nodelist, auxvar)
Allocate arrays.
Definition: gwe-esl.f90:121
subroutine esl_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwe-esl.f90:146
character(len=lenftype) ftype
Definition: gwe-esl.f90:15
subroutine define_listlabel(this)
Define list labels.
Definition: gwe-esl.f90:208
subroutine esl_df_obs(this)
Define observations.
Definition: gwe-esl.f90:252
subroutine esl_fc(this, rhs, ia, idxglo, matrix_sln)
Add matrix terms related to specified energy source loading.
Definition: gwe-esl.f90:173
logical function esl_obs_supported(this)
Support function for specified energy source loading observations.
Definition: gwe-esl.f90:238
character(len=16) text
Definition: gwe-esl.f90:16
subroutine esl_allocate_scalars(this)
Allocate scalars.
Definition: gwe-esl.f90:105
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
@ brief BndType
Data for sharing among multiple packages. Originally read in from.