MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
gweeslmodule Module Reference

Data Types

type  gweesltype
 

Functions/Subroutines

subroutine, public esl_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, gwecommon, input_mempath)
 Create an energy source loading package. More...
 
subroutine esl_da (this)
 Deallocate memory. More...
 
subroutine esl_allocate_scalars (this)
 Allocate scalars. More...
 
subroutine esl_allocate_arrays (this, nodelist, auxvar)
 Allocate arrays. More...
 
subroutine esl_ck (this)
 Check energy source loading boundary condition data. More...
 
subroutine esl_cf (this)
 Formulate the HCOF and RHS terms. More...
 
subroutine esl_fc (this, rhs, ia, idxglo, matrix_sln)
 Add matrix terms related to specified energy source loading. More...
 
subroutine define_listlabel (this)
 Define list labels. More...
 
logical function esl_obs_supported (this)
 Support function for specified energy source loading observations. More...
 
subroutine esl_df_obs (this)
 Define observations. More...
 
real(dp) function esl_bound_value (this, col, row)
 @ brief Return a bound value More...
 
real(dp) function ener_mult (this, row)
 Return a value that applies a multiplier. More...
 

Variables

character(len=lenftype) ftype = 'ESL'
 
character(len=16) text = ' ESL'
 

Function/Subroutine Documentation

◆ define_listlabel()

subroutine gweeslmodule::define_listlabel ( class(gweesltype), intent(inout)  this)
private

Define the list heading that is written to iout when PRINT_INPUT option is used.

Definition at line 244 of file gwe-esl.f90.

245  ! -- dummy
246  class(GweEslType), intent(inout) :: this
247  !
248  ! -- Create the header list label
249  this%listlabel = trim(this%filtyp)//' NO.'
250  if (this%dis%ndim == 3) then
251  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
252  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
253  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
254  elseif (this%dis%ndim == 2) then
255  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
256  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
257  else
258  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
259  end if
260  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
261  if (this%inamedbound == 1) then
262  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
263  end if

◆ ener_mult()

real(dp) function gweeslmodule::ener_mult ( class(gweesltype), intent(inout)  this,
integer(i4b), intent(in)  row 
)

Apply multiplier to specified energy load depending on user-selected option

Parameters
[in,out]thisBndExtType object

Definition at line 331 of file gwe-esl.f90.

332  ! -- modules
333  use constantsmodule, only: dzero
334  ! -- dummy variables
335  class(GweEslType), intent(inout) :: this !< BndExtType object
336  integer(I4B), intent(in) :: row
337  ! -- result
338  real(DP) :: ener
339  !
340  if (this%iauxmultcol > 0) then
341  ener = this%senerrate(row) * this%auxvar(this%iauxmultcol, row)
342  else
343  ener = this%senerrate(row)
344  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ esl_allocate_arrays()

subroutine gweeslmodule::esl_allocate_arrays ( class(gweesltype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
Parameters
nodelistpackage nodelist
auxvarpackage aux variable array

Definition at line 124 of file gwe-esl.f90.

126  ! -- dummy
127  class(GweEslType) :: this
128  ! -- local
129  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist !< package nodelist
130  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar !< package aux variable array
131 
132  ! -- base class allocate arrays
133  call this%BndExtType%allocate_arrays(nodelist, auxvar)
134 
135  ! -- set input context pointers
136  call mem_setptr(this%senerrate, 'SENERRATE', this%input_mempath)
137  !
138  ! -- checkin input context pointers
139  call mem_checkin(this%senerrate, 'SENERRATE', this%memoryPath, &
140  'SENERRATE', this%input_mempath)

◆ esl_allocate_scalars()

subroutine gweeslmodule::esl_allocate_scalars ( class(gweesltype this)

Allocate scalars specific to this energy source loading package

Definition at line 108 of file gwe-esl.f90.

109  ! -- modules
111  ! -- dummy
112  class(GweEslType) :: this
113  !
114  ! -- base class allocate scalars
115  call this%BndExtType%allocate_scalars()
116  !
117  ! -- allocate the object and assign values to object variables
118  !
119  ! -- Set values

◆ esl_bound_value()

real(dp) function gweeslmodule::esl_bound_value ( class(gweesltype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private

Return a bound value associated with an ncolbnd index and row.

Definition at line 309 of file gwe-esl.f90.

310  ! -- modules
311  use constantsmodule, only: dzero
312  ! -- dummy variables
313  class(GweEslType), intent(inout) :: this
314  integer(I4B), intent(in) :: col
315  integer(I4B), intent(in) :: row
316  ! -- result
317  real(DP) :: bndval
318  !
319  select case (col)
320  case (1)
321  bndval = this%senerrate(row)
322  case default
323  end select

◆ esl_cf()

subroutine gweeslmodule::esl_cf ( class(gweesltype this)
private

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

Definition at line 179 of file gwe-esl.f90.

180  ! -- dummy
181  class(GweEslType) :: this
182  ! -- local
183  integer(I4B) :: i, node
184  real(DP) :: q
185  !
186  ! -- Return if no sources
187  if (this%nbound == 0) return
188  !
189  ! -- Calculate hcof and rhs for each source entry
190  do i = 1, this%nbound
191  node = this%nodelist(i)
192  this%hcof(i) = dzero
193  if (this%ibound(node) <= 0) then
194  this%rhs(i) = dzero
195  cycle
196  end if
197  !
198  ! -- set energy loading rate accounting for multiplier
199  q = this%ener_mult(i)
200  !
201  this%rhs(i) = -q
202  end do

◆ esl_ck()

subroutine gweeslmodule::esl_ck ( class(gweesltype), intent(inout)  this)

Definition at line 145 of file gwe-esl.f90.

146  ! -- dummy
147  class(GweEslType), intent(inout) :: this
148  ! -- local
149  integer(I4B) :: i
150  integer(I4B) :: node
151  ! -- formats
152  character(len=*), parameter :: fmtenermulterr = &
153  "('ESL BOUNDARY (',i0,') ESL MULTIPLIER (',g10.3,') IS &
154  &LESS THAN ZERO THEREBY REVERSING THE ORIGINAL SIGN ON THE &
155  &AMOUNT OF ENERGY ENTERING OR EXITING THE MODEL.')"
156  !
157  ! -- check stress period data
158  do i = 1, this%nbound
159  node = this%nodelist(i)
160  !
161  ! -- accumulate warnings
162  if (this%iauxmultcol > 0) then
163  if (this%auxvar(this%iauxmultcol, i) < dzero) then
164  write (warnmsg, fmt=fmtenermulterr) &
165  i, this%auxvar(this%iauxmultcol, i)
166  call store_warning(warnmsg)
167  write (this%iout, '(/1x,a)') 'WARNING: '//trim(warnmsg)
168  end if
169  end if
170  end do
Here is the call graph for this function:

◆ esl_create()

subroutine, public gweeslmodule::esl_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=*), intent(in)  pakname,
type(gweinputdatatype), intent(in), target  gwecommon,
character(len=*), intent(in)  input_mempath 
)

This subroutine points bndobj to the newly created package

Parameters
[in]gwecommonshared data container for use by multiple GWE packages

Definition at line 48 of file gwe-esl.f90.

50  ! -- modules
51  use bndmodule, only: bndtype
52  ! -- dummy
53  class(BndType), pointer :: packobj
54  integer(I4B), intent(in) :: id
55  integer(I4B), intent(in) :: ibcnum
56  integer(I4B), intent(in) :: inunit
57  integer(I4B), intent(in) :: iout
58  character(len=*), intent(in) :: namemodel
59  character(len=*), intent(in) :: pakname
60  type(GweInputDataType), intent(in), target :: gwecommon !< shared data container for use by multiple GWE packages
61  character(len=*), intent(in) :: input_mempath
62  ! -- local
63  type(GweEslType), pointer :: eslobj
64  !
65  ! -- Allocate the object and assign values to object variables
66  allocate (eslobj)
67  packobj => eslobj
68  !
69  ! -- Create name and memory path
70  call packobj%set_names(ibcnum, namemodel, pakname, ftype, input_mempath)
71  packobj%text = text
72  !
73  ! -- Allocate scalars
74  call eslobj%allocate_scalars()
75  !
76  ! -- Initialize package
77  call packobj%pack_initialize()
78  !
79  packobj%inunit = inunit
80  packobj%iout = iout
81  packobj%id = id
82  packobj%ibcnum = ibcnum
83  packobj%ncolbnd = 1
84  packobj%iscloc = 1
85  !
86  ! -- Store pointer to shared data module for accessing cpw, rhow
87  ! for the budget calculations, and for accessing the latent heat of
88  ! vaporization for evaporative cooling.
89  eslobj%gwecommon => gwecommon
This module contains the base boundary package.
@ brief BndType
Here is the caller graph for this function:

◆ esl_da()

subroutine gweeslmodule::esl_da ( class(gweesltype this)

Definition at line 94 of file gwe-esl.f90.

95  ! -- modules
97  ! -- dummy
98  class(GweEslType) :: this
99  !
100  ! -- Deallocate parent package
101  call this%BndExtType%bnd_da()

◆ esl_df_obs()

subroutine gweeslmodule::esl_df_obs ( class(gweesltype this)
private

This subroutine:

  • stores observation types supported by ESL package.
  • overrides BndTypebnd_df_obs

Definition at line 288 of file gwe-esl.f90.

289  implicit none
290  ! -- dummy
291  class(GweEslType) :: this
292  ! -- local
293  integer(I4B) :: indx
294  !
295  call this%obs%StoreObsType('esl', .true., indx)
296  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
297  !
298  ! -- Store obs type and assign procedure pointer
299  ! for to-mvr observation type.
300  call this%obs%StoreObsType('to-mvr', .true., indx)
301  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ esl_fc()

subroutine gweeslmodule::esl_fc ( class(gweesltype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)
private

Copy rhs and hcof into solution rhs and amat

Definition at line 209 of file gwe-esl.f90.

210  ! -- dummy
211  class(GweEslType) :: this
212  real(DP), dimension(:), intent(inout) :: rhs
213  integer(I4B), dimension(:), intent(in) :: ia
214  integer(I4B), dimension(:), intent(in) :: idxglo
215  class(MatrixBaseType), pointer :: matrix_sln
216  ! -- local
217  integer(I4B) :: i, n, ipos
218  !
219  ! -- pakmvrobj fc
220  if (this%imover == 1) then
221  call this%pakmvrobj%fc()
222  end if
223  !
224  ! -- Copy package rhs and hcof into solution rhs and amat
225  do i = 1, this%nbound
226  n = this%nodelist(i)
227  rhs(n) = rhs(n) + this%rhs(i)
228  ipos = ia(n)
229  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
230  !
231  ! -- If mover is active and mass is being withdrawn,
232  ! store available mass (as positive value).
233  if (this%imover == 1 .and. this%rhs(i) > dzero) then
234  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
235  end if
236  end do

◆ esl_obs_supported()

logical function gweeslmodule::esl_obs_supported ( class(gweesltype this)
private

This function:

  • returns true because ESL package supports observations.
  • overrides BndTypebnd_obs_supported()

Definition at line 274 of file gwe-esl.f90.

275  implicit none
276  ! -- dummy
277  class(GweEslType) :: this
278  !
279  esl_obs_supported = .true.

Variable Documentation

◆ ftype

character(len=lenftype) gweeslmodule::ftype = 'ESL'
private

Definition at line 17 of file gwe-esl.f90.

17  character(len=LENFTYPE) :: ftype = 'ESL'

◆ text

character(len=16) gweeslmodule::text = ' ESL'
private

Definition at line 18 of file gwe-esl.f90.

18  character(len=16) :: text = ' ESL'