MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
gwelkemodule Module Reference

Data Types

type  gwelketype
 

Functions/Subroutines

subroutine, public lke_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
 Create a new lke package. More...
 
subroutine find_lke_package (this)
 Find corresponding lke package. More...
 
subroutine lke_fc_expanded (this, rhs, ia, idxglo, matrix_sln)
 Add matrix terms related to LKE. More...
 
subroutine lke_solve (this)
 Add terms specific to lakes to the explicit lake solve. More...
 
integer(i4b) function lke_get_nbudterms (this)
 Function to return the number of budget terms just for this package. More...
 
subroutine lke_setup_budobj (this, idx)
 Set up the budget object that stores all the lake flows. More...
 
subroutine lke_fill_budobj (this, idx, x, flowja, ccratin, ccratout)
 Copy flow terms into thisbudobj. More...
 
subroutine allocate_scalars (this)
 Allocate scalars specific to the lake energy transport (LKE) package. More...
 
subroutine lke_allocate_arrays (this)
 Allocate arrays specific to the lake energy transport (LKE) package. More...
 
subroutine lke_da (this)
 Deallocate memory. More...
 
subroutine lke_rain_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Rain term. More...
 
subroutine lke_evap_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Evaporative term. More...
 
subroutine lke_roff_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Runoff term. More...
 
subroutine lke_iflw_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Inflow Term. More...
 
subroutine lke_wdrl_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Specified withdrawal term. More...
 
subroutine lke_outf_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Outflow term. More...
 
subroutine lke_df_obs (this)
 Defined observation types. More...
 
subroutine lke_rp_obs (this, obsrv, found)
 Process package specific obs. More...
 
subroutine lke_bd_obs (this, obstypeid, jj, v, found)
 Calculate observation value and pass it back to APT. More...
 
subroutine lke_set_stressperiod (this, itemno, keyword, found)
 Sets the stress period attributes for keyword use. More...
 

Variables

character(len= *), parameter ftype = 'LKE'
 
character(len= *), parameter flowtype = 'LAK'
 
character(len=16) text = ' LKE'
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine gwelkemodule::allocate_scalars ( class(gwelketype this)

Definition at line 707 of file gwe-lke.f90.

708  ! -- modules
710  ! -- dummy
711  class(GweLkeType) :: this
712  !
713  ! -- Allocate scalars in TspAptType
714  call this%TspAptType%allocate_scalars()
715  !
716  ! -- Allocate
717  call mem_allocate(this%idxbudrain, 'IDXBUDRAIN', this%memoryPath)
718  call mem_allocate(this%idxbudevap, 'IDXBUDEVAP', this%memoryPath)
719  call mem_allocate(this%idxbudroff, 'IDXBUDROFF', this%memoryPath)
720  call mem_allocate(this%idxbudiflw, 'IDXBUDIFLW', this%memoryPath)
721  call mem_allocate(this%idxbudwdrl, 'IDXBUDWDRL', this%memoryPath)
722  call mem_allocate(this%idxbudoutf, 'IDXBUDOUTF', this%memoryPath)
723  call mem_allocate(this%idxbudlbcd, 'IDXBUDLBCD', this%memoryPath)
724  !
725  ! -- Initialize
726  this%idxbudrain = 0
727  this%idxbudevap = 0
728  this%idxbudroff = 0
729  this%idxbudiflw = 0
730  this%idxbudwdrl = 0
731  this%idxbudoutf = 0
732  this%idxbudlbcd = 0

◆ find_lke_package()

subroutine gwelkemodule::find_lke_package ( class(gwelketype this)

Definition at line 162 of file gwe-lke.f90.

163  ! -- modules
165  ! -- dummy
166  class(GweLkeType) :: this
167  ! -- local
168  character(len=LINELENGTH) :: errmsg
169  class(BndType), pointer :: packobj
170  integer(I4B) :: ip, icount
171  integer(I4B) :: nbudterm
172  logical :: found
173  !
174  ! -- Initialize found to false, and error later if flow package cannot
175  ! be found
176  found = .false.
177  !
178  ! -- If user is specifying flows in a binary budget file, then set up
179  ! the budget file reader, otherwise set a pointer to the flow package
180  ! budobj
181  if (this%fmi%flows_from_file) then
182  call this%fmi%set_aptbudobj_pointer(this%flowpackagename, this%flowbudptr)
183  if (associated(this%flowbudptr)) found = .true.
184  !
185  else
186  if (associated(this%fmi%gwfbndlist)) then
187  ! -- Look through gwfbndlist for a flow package with the same name as
188  ! this transport package name
189  do ip = 1, this%fmi%gwfbndlist%Count()
190  packobj => getbndfromlist(this%fmi%gwfbndlist, ip)
191  if (packobj%packName == this%flowpackagename) then
192  found = .true.
193  !
194  ! -- Store BndType pointer to packobj, and then
195  ! use the select type to point to the budobj in flow package
196  this%flowpackagebnd => packobj
197  select type (packobj)
198  type is (laktype)
199  this%flowbudptr => packobj%budobj
200  end select
201  end if
202  if (found) exit
203  end do
204  end if
205  end if
206  !
207  ! -- Error if flow package not found
208  if (.not. found) then
209  write (errmsg, '(a)') 'Could not find flow package with name '&
210  &//trim(adjustl(this%flowpackagename))//'.'
211  call store_error(errmsg)
212  call this%parser%StoreErrorUnit()
213  end if
214  !
215  ! -- Allocate space for idxbudssm, which indicates whether this is a
216  ! special budget term or one that is a general source and sink
217  nbudterm = this%flowbudptr%nbudterm
218  call mem_allocate(this%idxbudssm, nbudterm, 'IDXBUDSSM', this%memoryPath)
219  !
220  ! -- Process budget terms and identify special budget terms
221  write (this%iout, '(/, a, a)') &
222  'PROCESSING '//ftype//' INFORMATION FOR ', this%packName
223  write (this%iout, '(a)') ' IDENTIFYING FLOW TERMS IN '//flowtype//' PACKAGE'
224  write (this%iout, '(a, i0)') &
225  ' NUMBER OF '//flowtype//' = ', this%flowbudptr%ncv
226  icount = 1
227  do ip = 1, this%flowbudptr%nbudterm
228  select case (trim(adjustl(this%flowbudptr%budterm(ip)%flowtype)))
229  case ('FLOW-JA-FACE')
230  this%idxbudfjf = ip
231  this%idxbudssm(ip) = 0
232  case ('GWF')
233  this%idxbudgwf = ip
234  this%idxbudlbcd = ip
235  this%idxbudssm(ip) = 0
236  case ('STORAGE')
237  this%idxbudsto = ip
238  this%idxbudssm(ip) = 0
239  case ('RAINFALL')
240  this%idxbudrain = ip
241  this%idxbudssm(ip) = 0
242  case ('EVAPORATION')
243  this%idxbudevap = ip
244  this%idxbudssm(ip) = 0
245  case ('RUNOFF')
246  this%idxbudroff = ip
247  this%idxbudssm(ip) = 0
248  case ('EXT-INFLOW')
249  this%idxbudiflw = ip
250  this%idxbudssm(ip) = 0
251  case ('WITHDRAWAL')
252  this%idxbudwdrl = ip
253  this%idxbudssm(ip) = 0
254  case ('EXT-OUTFLOW')
255  this%idxbudoutf = ip
256  this%idxbudssm(ip) = 0
257  case ('TO-MVR')
258  this%idxbudtmvr = ip
259  this%idxbudssm(ip) = 0
260  case ('FROM-MVR')
261  this%idxbudfmvr = ip
262  this%idxbudssm(ip) = 0
263  case ('AUXILIARY')
264  this%idxbudaux = ip
265  this%idxbudssm(ip) = 0
266  case default
267  !
268  ! -- Set idxbudssm equal to a column index for where the temperatures
269  ! are stored in the tempbud(nbudssm, ncv) array
270  this%idxbudssm(ip) = icount
271  icount = icount + 1
272  end select
273  write (this%iout, '(a, i0, " = ", a,/, a, i0)') &
274  ' TERM ', ip, trim(adjustl(this%flowbudptr%budterm(ip)%flowtype)), &
275  ' MAX NO. OF ENTRIES = ', this%flowbudptr%budterm(ip)%maxlist
276  end do
277  write (this%iout, '(a, //)') 'DONE PROCESSING '//ftype//' INFORMATION'
Here is the call graph for this function:

◆ lke_allocate_arrays()

subroutine gwelkemodule::lke_allocate_arrays ( class(gwelketype), intent(inout)  this)

Definition at line 738 of file gwe-lke.f90.

739  ! -- modules
741  ! -- dummy
742  class(GweLkeType), intent(inout) :: this
743  ! -- local
744  integer(I4B) :: n
745  !
746  ! -- Time series
747  call mem_allocate(this%temprain, this%ncv, 'TEMPRAIN', this%memoryPath)
748  call mem_allocate(this%tempevap, this%ncv, 'TEMPEVAP', this%memoryPath)
749  call mem_allocate(this%temproff, this%ncv, 'TEMPROFF', this%memoryPath)
750  call mem_allocate(this%tempiflw, this%ncv, 'TEMPIFLW', this%memoryPath)
751  !
752  ! -- Call standard TspAptType allocate arrays
753  call this%TspAptType%apt_allocate_arrays()
754  !
755  ! -- Initialize
756  do n = 1, this%ncv
757  this%temprain(n) = dzero
758  this%tempevap(n) = dzero
759  this%temproff(n) = dzero
760  this%tempiflw(n) = dzero
761  end do
762  !

◆ lke_bd_obs()

subroutine gwelkemodule::lke_bd_obs ( class(gwelketype), intent(inout)  this,
character(len=*), intent(in)  obstypeid,
integer(i4b), intent(in)  jj,
real(dp), intent(inout)  v,
logical, intent(inout)  found 
)

Definition at line 1063 of file gwe-lke.f90.

1064  ! -- dummy
1065  class(GweLkeType), intent(inout) :: this
1066  character(len=*), intent(in) :: obstypeid
1067  real(DP), intent(inout) :: v
1068  integer(I4B), intent(in) :: jj
1069  logical, intent(inout) :: found
1070  ! -- local
1071  integer(I4B) :: n1, n2
1072  !
1073  found = .true.
1074  select case (obstypeid)
1075  case ('RAINFALL')
1076  if (this%iboundpak(jj) /= 0) then
1077  call this%lke_rain_term(jj, n1, n2, v)
1078  end if
1079  case ('EVAPORATION')
1080  if (this%iboundpak(jj) /= 0) then
1081  call this%lke_evap_term(jj, n1, n2, v)
1082  end if
1083  case ('RUNOFF')
1084  if (this%iboundpak(jj) /= 0) then
1085  call this%lke_roff_term(jj, n1, n2, v)
1086  end if
1087  case ('EXT-INFLOW')
1088  if (this%iboundpak(jj) /= 0) then
1089  call this%lke_iflw_term(jj, n1, n2, v)
1090  end if
1091  case ('WITHDRAWAL')
1092  if (this%iboundpak(jj) /= 0) then
1093  call this%lke_wdrl_term(jj, n1, n2, v)
1094  end if
1095  case ('EXT-OUTFLOW')
1096  if (this%iboundpak(jj) /= 0) then
1097  call this%lke_outf_term(jj, n1, n2, v)
1098  end if
1099  case default
1100  found = .false.
1101  end select

◆ lke_create()

subroutine, public gwelkemodule::lke_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(tspfmitype), pointer  fmi,
real(dp), intent(in), pointer  eqnsclfac,
type(gweinputdatatype), intent(in), target  gwecommon,
character(len=*), intent(in)  dvt,
character(len=*), intent(in)  dvu,
character(len=*), intent(in)  dvua 
)
Parameters
[in]eqnsclfacgoverning equation scale factor
[in]gwecommonshared data container for use by multiple GWE packages
[in]dvtFor GWE, set to "TEMPERATURE" in TspAptType
[in]dvuFor GWE, set to "energy" in TspAptType
[in]dvuaFor GWE, set to "E" in TspAptType

Definition at line 101 of file gwe-lke.f90.

103  ! -- dummy
104  class(BndType), pointer :: packobj
105  integer(I4B), intent(in) :: id
106  integer(I4B), intent(in) :: ibcnum
107  integer(I4B), intent(in) :: inunit
108  integer(I4B), intent(in) :: iout
109  character(len=*), intent(in) :: namemodel
110  character(len=*), intent(in) :: pakname
111  type(TspFmiType), pointer :: fmi
112  real(DP), intent(in), pointer :: eqnsclfac !< governing equation scale factor
113  type(GweInputDataType), intent(in), target :: gwecommon !< shared data container for use by multiple GWE packages
114  character(len=*), intent(in) :: dvt !< For GWE, set to "TEMPERATURE" in TspAptType
115  character(len=*), intent(in) :: dvu !< For GWE, set to "energy" in TspAptType
116  character(len=*), intent(in) :: dvua !< For GWE, set to "E" in TspAptType
117  ! -- local
118  type(GweLkeType), pointer :: lkeobj
119  !
120  ! -- Allocate the object and assign values to object variables
121  allocate (lkeobj)
122  packobj => lkeobj
123  !
124  ! -- Create name and memory path
125  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
126  packobj%text = text
127  !
128  ! -- Allocate scalars
129  call lkeobj%allocate_scalars()
130  !
131  ! -- Initialize package
132  call packobj%pack_initialize()
133  !
134  packobj%inunit = inunit
135  packobj%iout = iout
136  packobj%id = id
137  packobj%ibcnum = ibcnum
138  packobj%ncolbnd = 1
139  packobj%iscloc = 1
140  !
141  ! -- Store pointer to flow model interface. When the GwfGwe exchange is
142  ! created, it sets fmi%bndlist so that the GWE model has access to all
143  ! the flow packages
144  lkeobj%fmi => fmi
145  !
146  ! -- Store pointer to governing equation scale factor
147  lkeobj%eqnsclfac => eqnsclfac
148  !
149  ! -- Store pointer to shared data module for accessing cpw, rhow
150  ! for the budget calculations, and for accessing the latent heat of
151  ! vaporization for evaporative cooling.
152  lkeobj%gwecommon => gwecommon
153  !
154  ! -- Set labels that will be used in generalized APT class
155  lkeobj%depvartype = dvt
156  lkeobj%depvarunit = dvu
157  lkeobj%depvarunitabbrev = dvua
Here is the caller graph for this function:

◆ lke_da()

subroutine gwelkemodule::lke_da ( class(gwelketype this)

Definition at line 767 of file gwe-lke.f90.

768  ! -- modules
770  ! -- dummy
771  class(GweLkeType) :: this
772  !
773  ! -- Deallocate scalars
774  call mem_deallocate(this%idxbudrain)
775  call mem_deallocate(this%idxbudevap)
776  call mem_deallocate(this%idxbudroff)
777  call mem_deallocate(this%idxbudiflw)
778  call mem_deallocate(this%idxbudwdrl)
779  call mem_deallocate(this%idxbudoutf)
780  call mem_deallocate(this%idxbudlbcd)
781  !
782  ! -- Deallocate time series
783  call mem_deallocate(this%temprain)
784  call mem_deallocate(this%tempevap)
785  call mem_deallocate(this%temproff)
786  call mem_deallocate(this%tempiflw)
787  !
788  ! -- Deallocate scalars in TspAptType
789  call this%TspAptType%bnd_da()

◆ lke_df_obs()

subroutine gwelkemodule::lke_df_obs ( class(gwelketype this)

Store the observation type supported by the APT package and override BndTypebnd_df_obs

Definition at line 957 of file gwe-lke.f90.

958  ! -- dummy
959  class(GweLkeType) :: this
960  ! -- local
961  integer(I4B) :: indx
962  !
963  ! -- Store obs type and assign procedure pointer
964  ! for temperature observation type.
965  call this%obs%StoreObsType('temperature', .false., indx)
966  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
967  !
968  ! -- Store obs type and assign procedure pointer
969  ! for flow between features, such as lake to lake.
970  call this%obs%StoreObsType('flow-ja-face', .true., indx)
971  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid12
972  !
973  ! -- Store obs type and assign procedure pointer
974  ! for from-mvr observation type.
975  call this%obs%StoreObsType('from-mvr', .true., indx)
976  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
977  !
978  ! -- Store obs type and assign procedure pointer
979  ! for to-mvr observation type.
980  call this%obs%StoreObsType('to-mvr', .true., indx)
981  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
982  !
983  ! -- Store obs type and assign procedure pointer
984  ! for storage observation type.
985  call this%obs%StoreObsType('storage', .true., indx)
986  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
987  !
988  ! -- Store obs type and assign procedure pointer
989  ! for constant observation type.
990  call this%obs%StoreObsType('constant', .true., indx)
991  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
992  !
993  ! -- Store obs type and assign procedure pointer
994  ! for observation type: lke
995  call this%obs%StoreObsType('lke', .true., indx)
996  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
997  !
998  ! -- Store obs type and assign procedure pointer
999  ! for rainfall observation type.
1000  call this%obs%StoreObsType('rainfall', .true., indx)
1001  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1002  !
1003  ! -- Store obs type and assign procedure pointer
1004  ! for evaporation observation type.
1005  call this%obs%StoreObsType('evaporation', .true., indx)
1006  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1007  !
1008  ! -- Store obs type and assign procedure pointer
1009  ! for runoff observation type.
1010  call this%obs%StoreObsType('runoff', .true., indx)
1011  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1012  !
1013  ! -- Store obs type and assign procedure pointer
1014  ! for inflow observation type.
1015  call this%obs%StoreObsType('ext-inflow', .true., indx)
1016  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1017  !
1018  ! -- Store obs type and assign procedure pointer
1019  ! for withdrawal observation type.
1020  call this%obs%StoreObsType('withdrawal', .true., indx)
1021  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1022  !
1023  ! -- Store obs type and assign procedure pointer
1024  ! for ext-outflow observation type.
1025  call this%obs%StoreObsType('ext-outflow', .true., indx)
1026  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
Here is the call graph for this function:

◆ lke_evap_term()

subroutine gwelkemodule::lke_evap_term ( class(gwelketype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  n2,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)

Definition at line 819 of file gwe-lke.f90.

821  ! -- dummy
822  class(GweLkeType) :: this
823  integer(I4B), intent(in) :: ientry
824  integer(I4B), intent(inout) :: n1
825  integer(I4B), intent(inout) :: n2
826  real(DP), intent(inout), optional :: rrate
827  real(DP), intent(inout), optional :: rhsval
828  real(DP), intent(inout), optional :: hcofval
829  ! -- local
830  real(DP) :: qbnd
831  real(DP) :: heatlat
832  !
833  n1 = this%flowbudptr%budterm(this%idxbudevap)%id1(ientry)
834  n2 = this%flowbudptr%budterm(this%idxbudevap)%id2(ientry)
835  ! -- Note that qbnd is negative for evap
836  qbnd = this%flowbudptr%budterm(this%idxbudevap)%flow(ientry)
837  heatlat = this%gwecommon%gwerhow * this%gwecommon%gwelatheatvap
838  if (present(rrate)) rrate = qbnd * heatlat
839  if (present(rhsval)) rhsval = -rrate
840  if (present(hcofval)) hcofval = dzero

◆ lke_fc_expanded()

subroutine gwelkemodule::lke_fc_expanded ( class(gwelketype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)

This will be called from TspAptTypeapt_fc_expanded() in order to add matrix terms specifically for LKE

Definition at line 285 of file gwe-lke.f90.

286  ! -- dummy
287  class(GweLkeType) :: this
288  real(DP), dimension(:), intent(inout) :: rhs
289  integer(I4B), dimension(:), intent(in) :: ia
290  integer(I4B), dimension(:), intent(in) :: idxglo
291  class(MatrixBaseType), pointer :: matrix_sln
292  ! -- local
293  integer(I4B) :: j, n, n1, n2
294  integer(I4B) :: iloc
295  integer(I4B) :: iposd, iposoffd
296  integer(I4B) :: ipossymd, ipossymoffd
297  integer(I4B) :: auxpos
298  real(DP) :: rrate
299  real(DP) :: rhsval
300  real(DP) :: hcofval
301  real(DP) :: ctherm !< thermal conductance
302  real(DP) :: wa !< wetted area
303  real(DP) :: ktf !< thermal conductivity of streambed material
304  real(DP) :: s !< thickness of conductive streambed material
305  !
306  ! -- Add rainfall contribution
307  if (this%idxbudrain /= 0) then
308  do j = 1, this%flowbudptr%budterm(this%idxbudrain)%nlist
309  call this%lke_rain_term(j, n1, n2, rrate, rhsval, hcofval)
310  iloc = this%idxlocnode(n1)
311  iposd = this%idxpakdiag(n1)
312  call matrix_sln%add_value_pos(iposd, hcofval)
313  rhs(iloc) = rhs(iloc) + rhsval
314  end do
315  end if
316  !
317  ! -- Add evaporation contribution
318  if (this%idxbudevap /= 0) then
319  do j = 1, this%flowbudptr%budterm(this%idxbudevap)%nlist
320  call this%lke_evap_term(j, n1, n2, rrate, rhsval, hcofval)
321  iloc = this%idxlocnode(n1)
322  iposd = this%idxpakdiag(n1)
323  call matrix_sln%add_value_pos(iposd, hcofval)
324  rhs(iloc) = rhs(iloc) + rhsval
325  end do
326  end if
327  !
328  ! -- Add runoff contribution
329  if (this%idxbudroff /= 0) then
330  do j = 1, this%flowbudptr%budterm(this%idxbudroff)%nlist
331  call this%lke_roff_term(j, n1, n2, rrate, rhsval, hcofval)
332  iloc = this%idxlocnode(n1)
333  iposd = this%idxpakdiag(n1)
334  call matrix_sln%add_value_pos(iposd, hcofval)
335  rhs(iloc) = rhs(iloc) + rhsval
336  end do
337  end if
338  !
339  ! -- Add inflow contribution
340  if (this%idxbudiflw /= 0) then
341  do j = 1, this%flowbudptr%budterm(this%idxbudiflw)%nlist
342  call this%lke_iflw_term(j, n1, n2, rrate, rhsval, hcofval)
343  iloc = this%idxlocnode(n1)
344  iposd = this%idxpakdiag(n1)
345  call matrix_sln%add_value_pos(iposd, hcofval)
346  rhs(iloc) = rhs(iloc) + rhsval
347  end do
348  end if
349  !
350  ! -- Add withdrawal contribution
351  if (this%idxbudwdrl /= 0) then
352  do j = 1, this%flowbudptr%budterm(this%idxbudwdrl)%nlist
353  call this%lke_wdrl_term(j, n1, n2, rrate, rhsval, hcofval)
354  iloc = this%idxlocnode(n1)
355  iposd = this%idxpakdiag(n1)
356  call matrix_sln%add_value_pos(iposd, hcofval)
357  rhs(iloc) = rhs(iloc) + rhsval
358  end do
359  end if
360  !
361  ! -- Add outflow contribution
362  if (this%idxbudoutf /= 0) then
363  do j = 1, this%flowbudptr%budterm(this%idxbudoutf)%nlist
364  call this%lke_outf_term(j, n1, n2, rrate, rhsval, hcofval)
365  iloc = this%idxlocnode(n1)
366  iposd = this%idxpakdiag(n1)
367  call matrix_sln%add_value_pos(iposd, hcofval)
368  rhs(iloc) = rhs(iloc) + rhsval
369  end do
370  end if
371  !
372  ! -- Add lakebed conduction contribution
373  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
374  !
375  ! -- Set n to feature number and process if active feature
376  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
377  if (this%iboundpak(n) /= 0) then
378  !
379  ! -- Set acoef and rhs to negative so they are relative to sfe and not gwe
380  auxpos = this%flowbudptr%budterm(this%idxbudgwf)%naux
381  wa = this%flowbudptr%budterm(this%idxbudgwf)%auxvar(auxpos, j)
382  ktf = this%ktf(n)
383  s = this%rfeatthk(n)
384  ctherm = ktf * wa / s
385  !
386  ! -- Add to sfe row
387  iposd = this%idxdglo(j)
388  iposoffd = this%idxoffdglo(j)
389  call matrix_sln%add_value_pos(iposd, -ctherm) ! kluge note: make sure the signs on ctherm are correct here and below
390  call matrix_sln%add_value_pos(iposoffd, ctherm)
391  !
392  ! -- Add to gwe row for sfe connection
393  ipossymd = this%idxsymdglo(j)
394  ipossymoffd = this%idxsymoffdglo(j)
395  call matrix_sln%add_value_pos(ipossymd, -ctherm)
396  call matrix_sln%add_value_pos(ipossymoffd, ctherm)
397  end if
398  end do

◆ lke_fill_budobj()

subroutine gwelkemodule::lke_fill_budobj ( class(gwelketype this,
integer(i4b), intent(inout)  idx,
real(dp), dimension(:), intent(in)  x,
real(dp), dimension(:), intent(inout), contiguous  flowja,
real(dp), intent(inout)  ccratin,
real(dp), intent(inout)  ccratout 
)

Definition at line 598 of file gwe-lke.f90.

599  ! -- dummy
600  class(GweLkeType) :: this
601  integer(I4B), intent(inout) :: idx
602  real(DP), dimension(:), intent(in) :: x
603  real(DP), dimension(:), contiguous, intent(inout) :: flowja
604  real(DP), intent(inout) :: ccratin
605  real(DP), intent(inout) :: ccratout
606  ! -- local
607  integer(I4B) :: j, n1, n2
608  integer(I4B) :: nlist
609  integer(I4B) :: igwfnode
610  integer(I4B) :: idiag
611  integer(I4B) :: auxpos
612  real(DP) :: q
613  real(DP) :: ctherm !< thermal conductance
614  real(DP) :: wa !< wetted area
615  real(DP) :: ktf !< thermal conductivity of streambed material
616  real(DP) :: s !< thickness of conductive streambed materia
617  !
618  ! -- Rain
619  idx = idx + 1
620  nlist = this%flowbudptr%budterm(this%idxbudrain)%nlist
621  call this%budobj%budterm(idx)%reset(nlist)
622  do j = 1, nlist
623  call this%lke_rain_term(j, n1, n2, q)
624  call this%budobj%budterm(idx)%update_term(n1, n2, q)
625  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
626  end do
627  !
628  ! -- Evaporation
629  idx = idx + 1
630  nlist = this%flowbudptr%budterm(this%idxbudevap)%nlist
631  call this%budobj%budterm(idx)%reset(nlist)
632  do j = 1, nlist
633  call this%lke_evap_term(j, n1, n2, q)
634  call this%budobj%budterm(idx)%update_term(n1, n2, q)
635  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
636  end do
637  !
638  ! -- Runoff
639  idx = idx + 1
640  nlist = this%flowbudptr%budterm(this%idxbudroff)%nlist
641  call this%budobj%budterm(idx)%reset(nlist)
642  do j = 1, nlist
643  call this%lke_roff_term(j, n1, n2, q)
644  call this%budobj%budterm(idx)%update_term(n1, n2, q)
645  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
646  end do
647  !
648  ! -- Est-Inflow
649  idx = idx + 1
650  nlist = this%flowbudptr%budterm(this%idxbudiflw)%nlist
651  call this%budobj%budterm(idx)%reset(nlist)
652  do j = 1, nlist
653  call this%lke_iflw_term(j, n1, n2, q)
654  call this%budobj%budterm(idx)%update_term(n1, n2, q)
655  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
656  end do
657  !
658  ! -- Withdrawal
659  idx = idx + 1
660  nlist = this%flowbudptr%budterm(this%idxbudwdrl)%nlist
661  call this%budobj%budterm(idx)%reset(nlist)
662  do j = 1, nlist
663  call this%lke_wdrl_term(j, n1, n2, q)
664  call this%budobj%budterm(idx)%update_term(n1, n2, q)
665  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
666  end do
667  !
668  ! -- Ext-Outflow
669  idx = idx + 1
670  nlist = this%flowbudptr%budterm(this%idxbudoutf)%nlist
671  call this%budobj%budterm(idx)%reset(nlist)
672  do j = 1, nlist
673  call this%lke_outf_term(j, n1, n2, q)
674  call this%budobj%budterm(idx)%update_term(n1, n2, q)
675  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
676  end do
677  !
678  ! -- Lakebed-Cond
679  idx = idx + 1
680  call this%budobj%budterm(idx)%reset(this%maxbound)
681  do j = 1, this%flowbudptr%budterm(this%idxbudlbcd)%nlist
682  q = dzero
683  n1 = this%flowbudptr%budterm(this%idxbudlbcd)%id1(j)
684  if (this%iboundpak(n1) /= 0) then
685  igwfnode = this%flowbudptr%budterm(this%idxbudlbcd)%id2(j)
686  auxpos = this%flowbudptr%budterm(this%idxbudgwf)%naux ! for now there is only 1 aux variable under 'GWF'
687  wa = this%flowbudptr%budterm(this%idxbudgwf)%auxvar(auxpos, j)
688  ktf = this%ktf(n1)
689  s = this%rfeatthk(n1)
690  ctherm = ktf * wa / s
691  q = ctherm * (x(igwfnode) - this%xnewpak(n1))
692  end if
693  call this%budobj%budterm(idx)%update_term(n1, igwfnode, q)
694  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
695  if (this%iboundpak(n1) /= 0) then
696  ! -- Contribution to gwe cell budget
697  this%simvals(n1) = this%simvals(n1) - q
698  idiag = this%dis%con%ia(igwfnode)
699  flowja(idiag) = flowja(idiag) - q
700  end if
701  end do

◆ lke_get_nbudterms()

integer(i4b) function gwelkemodule::lke_get_nbudterms ( class(gwelketype this)

This overrides a function in the parent class.

Definition at line 464 of file gwe-lke.f90.

465  ! -- dummy
466  class(GweLkeType) :: this
467  ! -- Return
468  integer(I4B) :: nbudterms
469  !
470  ! -- Number of budget terms is 7
471  ! 1) rainfall
472  ! 2) evap
473  ! 3) runoff
474  ! 4) ext-inflow
475  ! 5) withdrawal
476  ! 6) ext-outflow
477  ! 7) lakebed-cond
478  !
479  nbudterms = 7

◆ lke_iflw_term()

subroutine gwelkemodule::lke_iflw_term ( class(gwelketype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  n2,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)

Accounts for energy flowing into a lake from a connected stream, for example.

Definition at line 873 of file gwe-lke.f90.

875  ! -- dummy
876  class(GweLkeType) :: this
877  integer(I4B), intent(in) :: ientry
878  integer(I4B), intent(inout) :: n1
879  integer(I4B), intent(inout) :: n2
880  real(DP), intent(inout), optional :: rrate
881  real(DP), intent(inout), optional :: rhsval
882  real(DP), intent(inout), optional :: hcofval
883  ! -- local
884  real(DP) :: qbnd
885  real(DP) :: ctmp
886  !
887  n1 = this%flowbudptr%budterm(this%idxbudiflw)%id1(ientry)
888  n2 = this%flowbudptr%budterm(this%idxbudiflw)%id2(ientry)
889  qbnd = this%flowbudptr%budterm(this%idxbudiflw)%flow(ientry)
890  ctmp = this%tempiflw(n1)
891  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
892  if (present(rhsval)) rhsval = -rrate
893  if (present(hcofval)) hcofval = dzero

◆ lke_outf_term()

subroutine gwelkemodule::lke_outf_term ( class(gwelketype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  n2,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)

Accounts for the energy leaving a lake, for example, energy exiting a lake via a flow into a draining stream channel.

Definition at line 929 of file gwe-lke.f90.

931  ! -- dummy
932  class(GweLkeType) :: this
933  integer(I4B), intent(in) :: ientry
934  integer(I4B), intent(inout) :: n1
935  integer(I4B), intent(inout) :: n2
936  real(DP), intent(inout), optional :: rrate
937  real(DP), intent(inout), optional :: rhsval
938  real(DP), intent(inout), optional :: hcofval
939  ! -- local
940  real(DP) :: qbnd
941  real(DP) :: ctmp
942  !
943  n1 = this%flowbudptr%budterm(this%idxbudoutf)%id1(ientry)
944  n2 = this%flowbudptr%budterm(this%idxbudoutf)%id2(ientry)
945  qbnd = this%flowbudptr%budterm(this%idxbudoutf)%flow(ientry)
946  ctmp = this%xnewpak(n1)
947  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
948  if (present(rhsval)) rhsval = dzero
949  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac

◆ lke_rain_term()

subroutine gwelkemodule::lke_rain_term ( class(gwelketype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  n2,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)

Definition at line 794 of file gwe-lke.f90.

796  ! -- dummy
797  class(GweLkeType) :: this
798  integer(I4B), intent(in) :: ientry
799  integer(I4B), intent(inout) :: n1
800  integer(I4B), intent(inout) :: n2
801  real(DP), intent(inout), optional :: rrate
802  real(DP), intent(inout), optional :: rhsval
803  real(DP), intent(inout), optional :: hcofval
804  ! -- local
805  real(DP) :: qbnd
806  real(DP) :: ctmp
807  !
808  n1 = this%flowbudptr%budterm(this%idxbudrain)%id1(ientry)
809  n2 = this%flowbudptr%budterm(this%idxbudrain)%id2(ientry)
810  qbnd = this%flowbudptr%budterm(this%idxbudrain)%flow(ientry)
811  ctmp = this%temprain(n1)
812  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
813  if (present(rhsval)) rhsval = -rrate
814  if (present(hcofval)) hcofval = dzero

◆ lke_roff_term()

subroutine gwelkemodule::lke_roff_term ( class(gwelketype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  n2,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)

Definition at line 845 of file gwe-lke.f90.

847  ! -- dummy
848  class(GweLkeType) :: this
849  integer(I4B), intent(in) :: ientry
850  integer(I4B), intent(inout) :: n1
851  integer(I4B), intent(inout) :: n2
852  real(DP), intent(inout), optional :: rrate
853  real(DP), intent(inout), optional :: rhsval
854  real(DP), intent(inout), optional :: hcofval
855  ! -- local
856  real(DP) :: qbnd
857  real(DP) :: ctmp
858  !
859  n1 = this%flowbudptr%budterm(this%idxbudroff)%id1(ientry)
860  n2 = this%flowbudptr%budterm(this%idxbudroff)%id2(ientry)
861  qbnd = this%flowbudptr%budterm(this%idxbudroff)%flow(ientry)
862  ctmp = this%temproff(n1)
863  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
864  if (present(rhsval)) rhsval = -rrate
865  if (present(hcofval)) hcofval = dzero

◆ lke_rp_obs()

subroutine gwelkemodule::lke_rp_obs ( class(gwelketype), intent(inout)  this,
type(observetype), intent(inout)  obsrv,
logical, intent(inout)  found 
)

Method to process specific observations for this package.

Parameters
[in,out]thispackage class
[in,out]obsrvobservation object
[in,out]foundindicate whether observation was found

Definition at line 1033 of file gwe-lke.f90.

1034  ! -- dummy
1035  class(GweLkeType), intent(inout) :: this !< package class
1036  type(ObserveType), intent(inout) :: obsrv !< observation object
1037  logical, intent(inout) :: found !< indicate whether observation was found
1038  !
1039  found = .true.
1040  select case (obsrv%ObsTypeId)
1041  case ('RAINFALL')
1042  call this%rp_obs_byfeature(obsrv)
1043  case ('EVAPORATION')
1044  call this%rp_obs_byfeature(obsrv)
1045  case ('RUNOFF')
1046  call this%rp_obs_byfeature(obsrv)
1047  case ('EXT-INFLOW')
1048  call this%rp_obs_byfeature(obsrv)
1049  case ('WITHDRAWAL')
1050  call this%rp_obs_byfeature(obsrv)
1051  case ('EXT-OUTFLOW')
1052  call this%rp_obs_byfeature(obsrv)
1053  case ('TO-MVR')
1054  call this%rp_obs_budterm(obsrv, &
1055  this%flowbudptr%budterm(this%idxbudtmvr))
1056  case default
1057  found = .false.
1058  end select

◆ lke_set_stressperiod()

subroutine gwelkemodule::lke_set_stressperiod ( class(gwelketype), intent(inout)  this,
integer(i4b), intent(in)  itemno,
character(len=*), intent(in)  keyword,
logical, intent(inout)  found 
)

Definition at line 1106 of file gwe-lke.f90.

1107  ! -- modules
1109  ! -- dummy
1110  class(GweLkeType), intent(inout) :: this
1111  integer(I4B), intent(in) :: itemno
1112  character(len=*), intent(in) :: keyword
1113  logical, intent(inout) :: found
1114  ! -- local
1115  character(len=LINELENGTH) :: text
1116  integer(I4B) :: ierr
1117  integer(I4B) :: jj
1118  real(DP), pointer :: bndElem => null()
1119  !
1120  ! RAINFALL <rainfall>
1121  ! EVAPORATION <evaporation>
1122  ! RUNOFF <runoff>
1123  ! EXT-INFLOW <inflow>
1124  ! WITHDRAWAL <withdrawal>
1125  !
1126  found = .true.
1127  select case (keyword)
1128  case ('RAINFALL')
1129  ierr = this%apt_check_valid(itemno)
1130  if (ierr /= 0) then
1131  goto 999
1132  end if
1133  call this%parser%GetString(text)
1134  jj = 1
1135  bndelem => this%temprain(itemno)
1136  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1137  this%packName, 'BND', this%tsManager, &
1138  this%iprpak, 'RAINFALL')
1139  case ('EVAPORATION')
1140  ierr = this%apt_check_valid(itemno)
1141  if (ierr /= 0) then
1142  goto 999
1143  end if
1144  call this%parser%GetString(text)
1145  jj = 1
1146  bndelem => this%tempevap(itemno)
1147  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1148  this%packName, 'BND', this%tsManager, &
1149  this%iprpak, 'EVAPORATION')
1150  case ('RUNOFF')
1151  ierr = this%apt_check_valid(itemno)
1152  if (ierr /= 0) then
1153  goto 999
1154  end if
1155  call this%parser%GetString(text)
1156  jj = 1
1157  bndelem => this%temproff(itemno)
1158  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1159  this%packName, 'BND', this%tsManager, &
1160  this%iprpak, 'RUNOFF')
1161  case ('EXT-INFLOW')
1162  ierr = this%apt_check_valid(itemno)
1163  if (ierr /= 0) then
1164  goto 999
1165  end if
1166  call this%parser%GetString(text)
1167  jj = 1
1168  bndelem => this%tempiflw(itemno)
1169  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1170  this%packName, 'BND', this%tsManager, &
1171  this%iprpak, 'EXT-INFLOW')
1172  case default
1173  !
1174  ! -- Keyword not recognized so return to caller with found = .false.
1175  found = .false.
1176  end select
1177  !
1178 999 continue
subroutine, public read_value_or_time_series_adv(textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, varName)
Call this subroutine from advanced packages to define timeseries link for a variable (varName).
Here is the call graph for this function:

◆ lke_setup_budobj()

subroutine gwelkemodule::lke_setup_budobj ( class(gwelketype this,
integer(i4b), intent(inout)  idx 
)

Definition at line 484 of file gwe-lke.f90.

485  ! -- modules
486  use constantsmodule, only: lenbudtxt
487  ! -- dummy
488  class(GweLkeType) :: this
489  integer(I4B), intent(inout) :: idx
490  ! -- local
491  integer(I4B) :: n, n1, n2
492  integer(I4B) :: maxlist, naux
493  real(DP) :: q
494  character(len=LENBUDTXT) :: text
495  !
496  ! -- Addition of heat associated with rainfall directly on the lake surface
497  text = ' RAINFALL'
498  idx = idx + 1
499  maxlist = this%flowbudptr%budterm(this%idxbudrain)%maxlist
500  naux = 0
501  call this%budobj%budterm(idx)%initialize(text, &
502  this%name_model, &
503  this%packName, &
504  this%name_model, &
505  this%packName, &
506  maxlist, .false., .false., &
507  naux)
508  !
509  ! -- Evaporative cooling from lake surface
510  text = ' EVAPORATION'
511  idx = idx + 1
512  maxlist = this%flowbudptr%budterm(this%idxbudevap)%maxlist
513  naux = 0
514  call this%budobj%budterm(idx)%initialize(text, &
515  this%name_model, &
516  this%packName, &
517  this%name_model, &
518  this%packName, &
519  maxlist, .false., .false., &
520  naux)
521  !
522  ! -- Addition of heat associated with runoff that flows to the lake
523  text = ' RUNOFF'
524  idx = idx + 1
525  maxlist = this%flowbudptr%budterm(this%idxbudroff)%maxlist
526  naux = 0
527  call this%budobj%budterm(idx)%initialize(text, &
528  this%name_model, &
529  this%packName, &
530  this%name_model, &
531  this%packName, &
532  maxlist, .false., .false., &
533  naux)
534  !
535  ! -- Addition of heat associated with user-specified inflow to the lake
536  text = ' EXT-INFLOW'
537  idx = idx + 1
538  maxlist = this%flowbudptr%budterm(this%idxbudiflw)%maxlist
539  naux = 0
540  call this%budobj%budterm(idx)%initialize(text, &
541  this%name_model, &
542  this%packName, &
543  this%name_model, &
544  this%packName, &
545  maxlist, .false., .false., &
546  naux)
547  !
548  ! -- Removal of heat associated with user-specified withdrawal from lake
549  text = ' WITHDRAWAL'
550  idx = idx + 1
551  maxlist = this%flowbudptr%budterm(this%idxbudwdrl)%maxlist
552  naux = 0
553  call this%budobj%budterm(idx)%initialize(text, &
554  this%name_model, &
555  this%packName, &
556  this%name_model, &
557  this%packName, &
558  maxlist, .false., .false., &
559  naux)
560  !
561  ! -- Removal of heat associated with outflow from lake that leaves
562  ! model domain
563  text = ' EXT-OUTFLOW'
564  idx = idx + 1
565  maxlist = this%flowbudptr%budterm(this%idxbudoutf)%maxlist
566  naux = 0
567  call this%budobj%budterm(idx)%initialize(text, &
568  this%name_model, &
569  this%packName, &
570  this%name_model, &
571  this%packName, &
572  maxlist, .false., .false., &
573  naux)
574  !
575  ! -- Conductive exchange of heat through the wetted lakebed
576  text = ' LAKEBED-COND'
577  idx = idx + 1
578  maxlist = this%flowbudptr%budterm(this%idxbudlbcd)%maxlist
579  naux = 0
580  call this%budobj%budterm(idx)%initialize(text, &
581  this%name_model, &
582  this%packName, &
583  this%name_model, &
584  this%packName, &
585  maxlist, .false., .false., &
586  naux)
587  call this%budobj%budterm(idx)%reset(maxlist)
588  q = dzero
589  do n = 1, maxlist
590  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(n)
591  n2 = this%flowbudptr%budterm(this%idxbudgwf)%id2(n)
592  call this%budobj%budterm(idx)%update_term(n1, n2, q)
593  end do
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37

◆ lke_solve()

subroutine gwelkemodule::lke_solve ( class(gwelketype this)

Definition at line 403 of file gwe-lke.f90.

404  ! -- dummy
405  class(GweLkeType) :: this
406  ! -- local
407  integer(I4B) :: j
408  integer(I4B) :: n1, n2
409  real(DP) :: rrate
410  !
411  ! -- Add rainfall contribution
412  if (this%idxbudrain /= 0) then
413  do j = 1, this%flowbudptr%budterm(this%idxbudrain)%nlist
414  call this%lke_rain_term(j, n1, n2, rrate)
415  this%dbuff(n1) = this%dbuff(n1) + rrate
416  end do
417  end if
418  !
419  ! -- Add evaporation contribution
420  if (this%idxbudevap /= 0) then
421  do j = 1, this%flowbudptr%budterm(this%idxbudevap)%nlist
422  call this%lke_evap_term(j, n1, n2, rrate)
423  this%dbuff(n1) = this%dbuff(n1) + rrate
424  end do
425  end if
426  !
427  ! -- Add runoff contribution
428  if (this%idxbudroff /= 0) then
429  do j = 1, this%flowbudptr%budterm(this%idxbudroff)%nlist
430  call this%lke_roff_term(j, n1, n2, rrate)
431  this%dbuff(n1) = this%dbuff(n1) + rrate
432  end do
433  end if
434  !
435  ! -- Add inflow contribution
436  if (this%idxbudiflw /= 0) then
437  do j = 1, this%flowbudptr%budterm(this%idxbudiflw)%nlist
438  call this%lke_iflw_term(j, n1, n2, rrate)
439  this%dbuff(n1) = this%dbuff(n1) + rrate
440  end do
441  end if
442  !
443  ! -- Add withdrawal contribution
444  if (this%idxbudwdrl /= 0) then
445  do j = 1, this%flowbudptr%budterm(this%idxbudwdrl)%nlist
446  call this%lke_wdrl_term(j, n1, n2, rrate)
447  this%dbuff(n1) = this%dbuff(n1) + rrate
448  end do
449  end if
450  !
451  ! -- Add outflow contribution
452  if (this%idxbudoutf /= 0) then
453  do j = 1, this%flowbudptr%budterm(this%idxbudoutf)%nlist
454  call this%lke_outf_term(j, n1, n2, rrate)
455  this%dbuff(n1) = this%dbuff(n1) + rrate
456  end do
457  end if

◆ lke_wdrl_term()

subroutine gwelkemodule::lke_wdrl_term ( class(gwelketype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  n2,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)

Accounts for energy associated with a withdrawal of water from a lake or group of lakes.

Definition at line 901 of file gwe-lke.f90.

903  ! -- dummy
904  class(GweLkeType) :: this
905  integer(I4B), intent(in) :: ientry
906  integer(I4B), intent(inout) :: n1
907  integer(I4B), intent(inout) :: n2
908  real(DP), intent(inout), optional :: rrate
909  real(DP), intent(inout), optional :: rhsval
910  real(DP), intent(inout), optional :: hcofval
911  ! -- local
912  real(DP) :: qbnd
913  real(DP) :: ctmp
914  !
915  n1 = this%flowbudptr%budterm(this%idxbudwdrl)%id1(ientry)
916  n2 = this%flowbudptr%budterm(this%idxbudwdrl)%id2(ientry)
917  qbnd = this%flowbudptr%budterm(this%idxbudwdrl)%flow(ientry)
918  ctmp = this%xnewpak(n1)
919  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
920  if (present(rhsval)) rhsval = dzero
921  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac

Variable Documentation

◆ flowtype

character(len=*), parameter gwelkemodule::flowtype = 'LAK'

Definition at line 53 of file gwe-lke.f90.

53  character(len=*), parameter :: flowtype = 'LAK'

◆ ftype

character(len=*), parameter gwelkemodule::ftype = 'LKE'

Definition at line 52 of file gwe-lke.f90.

52  character(len=*), parameter :: ftype = 'LKE'

◆ text

character(len=16) gwelkemodule::text = ' LKE'

Definition at line 54 of file gwe-lke.f90.

54  character(len=16) :: text = ' LKE'