MODFLOW 6  version 6.7.0.dev0
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)
 Create an energy source loading package. More...
 
subroutine esl_da (this)
 Deallocate memory. More...
 
subroutine esl_allocate_scalars (this)
 Allocate scalars. 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...
 
subroutine esl_rp_ts (this)
 Procedure related to time series. 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 185 of file gwe-esl.f90.

186  ! -- dummy
187  class(GweEslType), intent(inout) :: this
188  !
189  ! -- Create the header list label
190  this%listlabel = trim(this%filtyp)//' NO.'
191  if (this%dis%ndim == 3) then
192  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
193  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
194  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
195  elseif (this%dis%ndim == 2) then
196  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
197  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
198  else
199  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
200  end if
201  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
202  if (this%inamedbound == 1) then
203  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
204  end if

◆ esl_allocate_scalars()

subroutine gweeslmodule::esl_allocate_scalars ( class(gweesltype this)

Allocate scalars specific to this energy source loading package

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

104  ! -- modules
106  ! -- dummy
107  class(GweEslType) :: this
108  !
109  ! -- call standard BndType allocate scalars
110  call this%BndType%allocate_scalars()
111  !
112  ! -- allocate the object and assign values to object variables
113  !
114  ! -- Set values

◆ esl_cf()

subroutine gweeslmodule::esl_cf ( class(gweesltype this)

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

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

124  ! -- dummy
125  class(GweEslType) :: this
126  ! -- local
127  integer(I4B) :: i, node
128  real(DP) :: q
129  !
130  ! -- Return if no sources
131  if (this%nbound == 0) return
132  !
133  ! -- Calculate hcof and rhs for each source entry
134  do i = 1, this%nbound
135  node = this%nodelist(i)
136  this%hcof(i) = dzero
137  if (this%ibound(node) <= 0) then
138  this%rhs(i) = dzero
139  cycle
140  end if
141  q = this%bound(1, i)
142  this%rhs(i) = -q
143  end do

◆ 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 
)

This subroutine points bndobj to the newly created package

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

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

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  ! -- local
58  type(GweEslType), pointer :: eslobj
59  !
60  ! -- Allocate the object and assign values to object variables
61  allocate (eslobj)
62  packobj => eslobj
63  !
64  ! -- Create name and memory path
65  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
66  packobj%text = text
67  !
68  ! -- Allocate scalars
69  call eslobj%allocate_scalars()
70  !
71  ! -- Initialize package
72  call packobj%pack_initialize()
73  !
74  packobj%inunit = inunit
75  packobj%iout = iout
76  packobj%id = id
77  packobj%ibcnum = ibcnum
78  packobj%ncolbnd = 1
79  packobj%iscloc = 1
80  !
81  ! -- Store pointer to shared data module for accessing cpw, rhow
82  ! for the budget calculations, and for accessing the latent heat of
83  ! vaporization for evaporative cooling.
84  eslobj%gwecommon => gwecommon
Here is the caller graph for this function:

◆ esl_da()

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

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

90  ! -- modules
92  ! -- dummy
93  class(GweEslType) :: this
94  !
95  ! -- Deallocate parent package
96  call this%BndType%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 229 of file gwe-esl.f90.

230  implicit none
231  ! -- dummy
232  class(GweEslType) :: this
233  ! -- local
234  integer(I4B) :: indx
235  !
236  call this%obs%StoreObsType('esl', .true., indx)
237  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
238  !
239  ! -- Store obs type and assign procedure pointer
240  ! for to-mvr observation type.
241  call this%obs%StoreObsType('to-mvr', .true., indx)
242  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 150 of file gwe-esl.f90.

151  ! -- dummy
152  class(GweEslType) :: this
153  real(DP), dimension(:), intent(inout) :: rhs
154  integer(I4B), dimension(:), intent(in) :: ia
155  integer(I4B), dimension(:), intent(in) :: idxglo
156  class(MatrixBaseType), pointer :: matrix_sln
157  ! -- local
158  integer(I4B) :: i, n, ipos
159  !
160  ! -- pakmvrobj fc
161  if (this%imover == 1) then
162  call this%pakmvrobj%fc()
163  end if
164  !
165  ! -- Copy package rhs and hcof into solution rhs and amat
166  do i = 1, this%nbound
167  n = this%nodelist(i)
168  rhs(n) = rhs(n) + this%rhs(i)
169  ipos = ia(n)
170  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
171  !
172  ! -- If mover is active and mass is being withdrawn,
173  ! store available mass (as positive value).
174  if (this%imover == 1 .and. this%rhs(i) > dzero) then
175  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
176  end if
177  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 215 of file gwe-esl.f90.

216  implicit none
217  ! -- dummy
218  class(GweEslType) :: this
219  !
220  esl_obs_supported = .true.

◆ esl_rp_ts()

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

Assign tsLinkText appropriately for all time series in use by package. In the ESL package only the SENERRATE variable can be controlled by time series.

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

252  ! -- dummy
253  class(GweEslType), intent(inout) :: this
254  ! -- local
255  integer(I4B) :: i, nlinks
256  type(TimeSeriesLinkType), pointer :: tslink => null()
257  !
258  nlinks = this%TsManager%boundtslinks%Count()
259  do i = 1, nlinks
260  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
261  if (associated(tslink)) then
262  if (tslink%JCol == 1) then
263  tslink%Text = 'SMASSRATE'
264  end if
265  end if
266  end do
Here is the call graph for this function:

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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