MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 200 of file gwe-esl.f90.

201  ! -- dummy
202  class(GweEslType), intent(inout) :: this
203  !
204  ! -- Create the header list label
205  this%listlabel = trim(this%filtyp)//' NO.'
206  if (this%dis%ndim == 3) then
207  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
208  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
209  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
210  elseif (this%dis%ndim == 2) then
211  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
212  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
213  else
214  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
215  end if
216  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
217  if (this%inamedbound == 1) then
218  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
219  end if
220  !
221  ! -- Return
222  return

◆ esl_allocate_scalars()

subroutine gweeslmodule::esl_allocate_scalars ( class(gweesltype this)

Allocate scalars specific to this energy source loading package

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

110  ! -- modules
112  ! -- dummy
113  class(GweEslType) :: this
114  !
115  ! -- call standard BndType allocate scalars
116  call this%BndType%allocate_scalars()
117  !
118  ! -- allocate the object and assign values to object variables
119  !
120  ! -- Set values
121  !
122  ! -- Return
123  return

◆ esl_cf()

subroutine gweeslmodule::esl_cf ( class(gweesltype this)

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

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

133  ! -- dummy
134  class(GweEslType) :: this
135  ! -- local
136  integer(I4B) :: i, node
137  real(DP) :: q
138  !
139  ! -- Return if no sources
140  if (this%nbound == 0) return
141  !
142  ! -- Calculate hcof and rhs for each source entry
143  do i = 1, this%nbound
144  node = this%nodelist(i)
145  this%hcof(i) = dzero
146  if (this%ibound(node) <= 0) then
147  this%rhs(i) = dzero
148  cycle
149  end if
150  q = this%bound(1, i)
151  this%rhs(i) = -q
152  end do
153  !
154  ! -- Return
155  return

◆ 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
85  !
86  ! -- Return
87  return
Here is the caller graph for this function:

◆ esl_da()

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

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

93  ! -- modules
95  ! -- dummy
96  class(GweEslType) :: this
97  !
98  ! -- Deallocate parent package
99  call this%BndType%bnd_da()
100  !
101  ! -- Return
102  return

◆ 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 250 of file gwe-esl.f90.

251  implicit none
252  ! -- dummy
253  class(GweEslType) :: this
254  ! -- local
255  integer(I4B) :: indx
256  !
257  call this%obs%StoreObsType('esl', .true., indx)
258  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
259  !
260  ! -- Store obs type and assign procedure pointer
261  ! for to-mvr observation type.
262  call this%obs%StoreObsType('to-mvr', .true., indx)
263  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
264  !
265  ! -- Return
266  return
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 162 of file gwe-esl.f90.

163  ! -- dummy
164  class(GweEslType) :: this
165  real(DP), dimension(:), intent(inout) :: rhs
166  integer(I4B), dimension(:), intent(in) :: ia
167  integer(I4B), dimension(:), intent(in) :: idxglo
168  class(MatrixBaseType), pointer :: matrix_sln
169  ! -- local
170  integer(I4B) :: i, n, ipos
171  !
172  ! -- pakmvrobj fc
173  if (this%imover == 1) then
174  call this%pakmvrobj%fc()
175  end if
176  !
177  ! -- Copy package rhs and hcof into solution rhs and amat
178  do i = 1, this%nbound
179  n = this%nodelist(i)
180  rhs(n) = rhs(n) + this%rhs(i)
181  ipos = ia(n)
182  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
183  !
184  ! -- If mover is active and mass is being withdrawn,
185  ! store available mass (as positive value).
186  if (this%imover == 1 .and. this%rhs(i) > dzero) then
187  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
188  end if
189  end do
190  !
191  ! -- Return
192  return

◆ 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 233 of file gwe-esl.f90.

234  implicit none
235  ! -- dummy
236  class(GweEslType) :: this
237  !
238  esl_obs_supported = .true.
239  !
240  ! -- Return
241  return

◆ 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 275 of file gwe-esl.f90.

276  ! -- dummy
277  class(GweEslType), intent(inout) :: this
278  ! -- local
279  integer(I4B) :: i, nlinks
280  type(TimeSeriesLinkType), pointer :: tslink => null()
281  !
282  nlinks = this%TsManager%boundtslinks%Count()
283  do i = 1, nlinks
284  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
285  if (associated(tslink)) then
286  if (tslink%JCol == 1) then
287  tslink%Text = 'SMASSRATE'
288  end if
289  end if
290  end do
291  !
292  ! -- Return
293  return
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'