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

729  ! -- modules
731  ! -- dummy
732  class(GweLkeType) :: this
733  !
734  ! -- Allocate scalars in TspAptType
735  call this%TspAptType%allocate_scalars()
736  !
737  ! -- Allocate
738  call mem_allocate(this%idxbudrain, 'IDXBUDRAIN', this%memoryPath)
739  call mem_allocate(this%idxbudevap, 'IDXBUDEVAP', this%memoryPath)
740  call mem_allocate(this%idxbudroff, 'IDXBUDROFF', this%memoryPath)
741  call mem_allocate(this%idxbudiflw, 'IDXBUDIFLW', this%memoryPath)
742  call mem_allocate(this%idxbudwdrl, 'IDXBUDWDRL', this%memoryPath)
743  call mem_allocate(this%idxbudoutf, 'IDXBUDOUTF', this%memoryPath)
744  call mem_allocate(this%idxbudlbcd, 'IDXBUDLBCD', this%memoryPath)
745  !
746  ! -- Initialize
747  this%idxbudrain = 0
748  this%idxbudevap = 0
749  this%idxbudroff = 0
750  this%idxbudiflw = 0
751  this%idxbudwdrl = 0
752  this%idxbudoutf = 0
753  this%idxbudlbcd = 0
754  !
755  ! -- Return
756  return

◆ find_lke_package()

subroutine gwelkemodule::find_lke_package ( class(gwelketype this)

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

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

◆ lke_allocate_arrays()

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

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

763  ! -- modules
765  ! -- dummy
766  class(GweLkeType), intent(inout) :: this
767  ! -- local
768  integer(I4B) :: n
769  !
770  ! -- Time series
771  call mem_allocate(this%temprain, this%ncv, 'TEMPRAIN', this%memoryPath)
772  call mem_allocate(this%tempevap, this%ncv, 'TEMPEVAP', this%memoryPath)
773  call mem_allocate(this%temproff, this%ncv, 'TEMPROFF', this%memoryPath)
774  call mem_allocate(this%tempiflw, this%ncv, 'TEMPIFLW', this%memoryPath)
775  !
776  ! -- Call standard TspAptType allocate arrays
777  call this%TspAptType%apt_allocate_arrays()
778  !
779  ! -- Initialize
780  do n = 1, this%ncv
781  this%temprain(n) = dzero
782  this%tempevap(n) = dzero
783  this%temproff(n) = dzero
784  this%tempiflw(n) = dzero
785  end do
786  !
787  !
788  ! -- Return
789  return

◆ 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 1117 of file gwe-lke.f90.

1118  ! -- dummy
1119  class(GweLkeType), intent(inout) :: this
1120  character(len=*), intent(in) :: obstypeid
1121  real(DP), intent(inout) :: v
1122  integer(I4B), intent(in) :: jj
1123  logical, intent(inout) :: found
1124  ! -- local
1125  integer(I4B) :: n1, n2
1126  !
1127  found = .true.
1128  select case (obstypeid)
1129  case ('RAINFALL')
1130  if (this%iboundpak(jj) /= 0) then
1131  call this%lke_rain_term(jj, n1, n2, v)
1132  end if
1133  case ('EVAPORATION')
1134  if (this%iboundpak(jj) /= 0) then
1135  call this%lke_evap_term(jj, n1, n2, v)
1136  end if
1137  case ('RUNOFF')
1138  if (this%iboundpak(jj) /= 0) then
1139  call this%lke_roff_term(jj, n1, n2, v)
1140  end if
1141  case ('EXT-INFLOW')
1142  if (this%iboundpak(jj) /= 0) then
1143  call this%lke_iflw_term(jj, n1, n2, v)
1144  end if
1145  case ('WITHDRAWAL')
1146  if (this%iboundpak(jj) /= 0) then
1147  call this%lke_wdrl_term(jj, n1, n2, v)
1148  end if
1149  case ('EXT-OUTFLOW')
1150  if (this%iboundpak(jj) /= 0) then
1151  call this%lke_outf_term(jj, n1, n2, v)
1152  end if
1153  case default
1154  found = .false.
1155  end select
1156  !
1157  ! -- Return
1158  return

◆ 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
158  !
159  ! -- Return
160  return
Here is the caller graph for this function:

◆ lke_da()

subroutine gwelkemodule::lke_da ( class(gwelketype this)

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

795  ! -- modules
797  ! -- dummy
798  class(GweLkeType) :: this
799  !
800  ! -- Deallocate scalars
801  call mem_deallocate(this%idxbudrain)
802  call mem_deallocate(this%idxbudevap)
803  call mem_deallocate(this%idxbudroff)
804  call mem_deallocate(this%idxbudiflw)
805  call mem_deallocate(this%idxbudwdrl)
806  call mem_deallocate(this%idxbudoutf)
807  call mem_deallocate(this%idxbudlbcd)
808  !
809  ! -- Deallocate time series
810  call mem_deallocate(this%temprain)
811  call mem_deallocate(this%tempevap)
812  call mem_deallocate(this%temproff)
813  call mem_deallocate(this%tempiflw)
814  !
815  ! -- Deallocate scalars in TspAptType
816  call this%TspAptType%bnd_da()
817  !
818  ! -- Return
819  return

◆ 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 1005 of file gwe-lke.f90.

1006  ! -- dummy
1007  class(GweLkeType) :: this
1008  ! -- local
1009  integer(I4B) :: indx
1010  !
1011  ! -- Store obs type and assign procedure pointer
1012  ! for temperature observation type.
1013  call this%obs%StoreObsType('temperature', .false., indx)
1014  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1015  !
1016  ! -- Store obs type and assign procedure pointer
1017  ! for flow between features, such as lake to lake.
1018  call this%obs%StoreObsType('flow-ja-face', .true., indx)
1019  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid12
1020  !
1021  ! -- Store obs type and assign procedure pointer
1022  ! for from-mvr observation type.
1023  call this%obs%StoreObsType('from-mvr', .true., indx)
1024  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1025  !
1026  ! -- Store obs type and assign procedure pointer
1027  ! for to-mvr observation type.
1028  call this%obs%StoreObsType('to-mvr', .true., indx)
1029  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1030  !
1031  ! -- Store obs type and assign procedure pointer
1032  ! for storage observation type.
1033  call this%obs%StoreObsType('storage', .true., indx)
1034  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1035  !
1036  ! -- Store obs type and assign procedure pointer
1037  ! for constant observation type.
1038  call this%obs%StoreObsType('constant', .true., indx)
1039  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1040  !
1041  ! -- Store obs type and assign procedure pointer
1042  ! for observation type: lke
1043  call this%obs%StoreObsType('lke', .true., indx)
1044  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1045  !
1046  ! -- Store obs type and assign procedure pointer
1047  ! for rainfall observation type.
1048  call this%obs%StoreObsType('rainfall', .true., indx)
1049  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1050  !
1051  ! -- Store obs type and assign procedure pointer
1052  ! for evaporation observation type.
1053  call this%obs%StoreObsType('evaporation', .true., indx)
1054  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1055  !
1056  ! -- Store obs type and assign procedure pointer
1057  ! for runoff observation type.
1058  call this%obs%StoreObsType('runoff', .true., indx)
1059  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1060  !
1061  ! -- Store obs type and assign procedure pointer
1062  ! for inflow observation type.
1063  call this%obs%StoreObsType('ext-inflow', .true., indx)
1064  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1065  !
1066  ! -- Store obs type and assign procedure pointer
1067  ! for withdrawal observation type.
1068  call this%obs%StoreObsType('withdrawal', .true., indx)
1069  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1070  !
1071  ! -- Store obs type and assign procedure pointer
1072  ! for ext-outflow observation type.
1073  call this%obs%StoreObsType('ext-outflow', .true., indx)
1074  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
1075  !
1076  ! -- Return
1077  return
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 852 of file gwe-lke.f90.

854  ! -- dummy
855  class(GweLkeType) :: this
856  integer(I4B), intent(in) :: ientry
857  integer(I4B), intent(inout) :: n1
858  integer(I4B), intent(inout) :: n2
859  real(DP), intent(inout), optional :: rrate
860  real(DP), intent(inout), optional :: rhsval
861  real(DP), intent(inout), optional :: hcofval
862  ! -- local
863  real(DP) :: qbnd
864  real(DP) :: heatlat
865  !
866  n1 = this%flowbudptr%budterm(this%idxbudevap)%id1(ientry)
867  n2 = this%flowbudptr%budterm(this%idxbudevap)%id2(ientry)
868  ! -- Note that qbnd is negative for evap
869  qbnd = this%flowbudptr%budterm(this%idxbudevap)%flow(ientry)
870  heatlat = this%gwecommon%gwerhow * this%gwecommon%gwelatheatvap
871  if (present(rrate)) rrate = qbnd * heatlat
872  if (present(rhsval)) rhsval = -rrate
873  if (present(hcofval)) hcofval = dzero
874  !
875  ! -- Return
876  return

◆ 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 291 of file gwe-lke.f90.

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

◆ 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 616 of file gwe-lke.f90.

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

◆ lke_get_nbudterms()

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

This overrides a function in the parent class.

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

477  ! -- dummy
478  class(GweLkeType) :: this
479  ! -- Return
480  integer(I4B) :: nbudterms
481  !
482  ! -- Number of budget terms is 7
483  ! 1) rainfall
484  ! 2) evap
485  ! 3) runoff
486  ! 4) ext-inflow
487  ! 5) withdrawal
488  ! 6) ext-outflow
489  ! 7) lakebed-cond
490  !
491  nbudterms = 7
492  !
493  ! -- Return
494  return

◆ 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 912 of file gwe-lke.f90.

914  ! -- dummy
915  class(GweLkeType) :: this
916  integer(I4B), intent(in) :: ientry
917  integer(I4B), intent(inout) :: n1
918  integer(I4B), intent(inout) :: n2
919  real(DP), intent(inout), optional :: rrate
920  real(DP), intent(inout), optional :: rhsval
921  real(DP), intent(inout), optional :: hcofval
922  ! -- local
923  real(DP) :: qbnd
924  real(DP) :: ctmp
925  !
926  n1 = this%flowbudptr%budterm(this%idxbudiflw)%id1(ientry)
927  n2 = this%flowbudptr%budterm(this%idxbudiflw)%id2(ientry)
928  qbnd = this%flowbudptr%budterm(this%idxbudiflw)%flow(ientry)
929  ctmp = this%tempiflw(n1)
930  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
931  if (present(rhsval)) rhsval = -rrate
932  if (present(hcofval)) hcofval = dzero
933  !
934  ! -- Return
935  return

◆ 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 974 of file gwe-lke.f90.

976  ! -- dummy
977  class(GweLkeType) :: this
978  integer(I4B), intent(in) :: ientry
979  integer(I4B), intent(inout) :: n1
980  integer(I4B), intent(inout) :: n2
981  real(DP), intent(inout), optional :: rrate
982  real(DP), intent(inout), optional :: rhsval
983  real(DP), intent(inout), optional :: hcofval
984  ! -- local
985  real(DP) :: qbnd
986  real(DP) :: ctmp
987  !
988  n1 = this%flowbudptr%budterm(this%idxbudoutf)%id1(ientry)
989  n2 = this%flowbudptr%budterm(this%idxbudoutf)%id2(ientry)
990  qbnd = this%flowbudptr%budterm(this%idxbudoutf)%flow(ientry)
991  ctmp = this%xnewpak(n1)
992  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
993  if (present(rhsval)) rhsval = dzero
994  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac
995  !
996  ! -- Return
997  return

◆ 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 824 of file gwe-lke.f90.

826  ! -- dummy
827  class(GweLkeType) :: this
828  integer(I4B), intent(in) :: ientry
829  integer(I4B), intent(inout) :: n1
830  integer(I4B), intent(inout) :: n2
831  real(DP), intent(inout), optional :: rrate
832  real(DP), intent(inout), optional :: rhsval
833  real(DP), intent(inout), optional :: hcofval
834  ! -- local
835  real(DP) :: qbnd
836  real(DP) :: ctmp
837  !
838  n1 = this%flowbudptr%budterm(this%idxbudrain)%id1(ientry)
839  n2 = this%flowbudptr%budterm(this%idxbudrain)%id2(ientry)
840  qbnd = this%flowbudptr%budterm(this%idxbudrain)%flow(ientry)
841  ctmp = this%temprain(n1)
842  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
843  if (present(rhsval)) rhsval = -rrate
844  if (present(hcofval)) hcofval = dzero
845  !
846  ! -- Return
847  return

◆ 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 881 of file gwe-lke.f90.

883  ! -- dummy
884  class(GweLkeType) :: this
885  integer(I4B), intent(in) :: ientry
886  integer(I4B), intent(inout) :: n1
887  integer(I4B), intent(inout) :: n2
888  real(DP), intent(inout), optional :: rrate
889  real(DP), intent(inout), optional :: rhsval
890  real(DP), intent(inout), optional :: hcofval
891  ! -- local
892  real(DP) :: qbnd
893  real(DP) :: ctmp
894  !
895  n1 = this%flowbudptr%budterm(this%idxbudroff)%id1(ientry)
896  n2 = this%flowbudptr%budterm(this%idxbudroff)%id2(ientry)
897  qbnd = this%flowbudptr%budterm(this%idxbudroff)%flow(ientry)
898  ctmp = this%temproff(n1)
899  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
900  if (present(rhsval)) rhsval = -rrate
901  if (present(hcofval)) hcofval = dzero
902  !
903  ! -- Return
904  return

◆ 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 1084 of file gwe-lke.f90.

1085  ! -- dummy
1086  class(GweLkeType), intent(inout) :: this !< package class
1087  type(ObserveType), intent(inout) :: obsrv !< observation object
1088  logical, intent(inout) :: found !< indicate whether observation was found
1089  !
1090  found = .true.
1091  select case (obsrv%ObsTypeId)
1092  case ('RAINFALL')
1093  call this%rp_obs_byfeature(obsrv)
1094  case ('EVAPORATION')
1095  call this%rp_obs_byfeature(obsrv)
1096  case ('RUNOFF')
1097  call this%rp_obs_byfeature(obsrv)
1098  case ('EXT-INFLOW')
1099  call this%rp_obs_byfeature(obsrv)
1100  case ('WITHDRAWAL')
1101  call this%rp_obs_byfeature(obsrv)
1102  case ('EXT-OUTFLOW')
1103  call this%rp_obs_byfeature(obsrv)
1104  case ('TO-MVR')
1105  call this%rp_obs_budterm(obsrv, &
1106  this%flowbudptr%budterm(this%idxbudtmvr))
1107  case default
1108  found = .false.
1109  end select
1110  !
1111  ! -- Return
1112  return

◆ 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 1163 of file gwe-lke.f90.

1164  ! -- modules
1166  ! -- dummy
1167  class(GweLkeType), intent(inout) :: this
1168  integer(I4B), intent(in) :: itemno
1169  character(len=*), intent(in) :: keyword
1170  logical, intent(inout) :: found
1171  ! -- local
1172  character(len=LINELENGTH) :: text
1173  integer(I4B) :: ierr
1174  integer(I4B) :: jj
1175  real(DP), pointer :: bndElem => null()
1176  !
1177  ! RAINFALL <rainfall>
1178  ! EVAPORATION <evaporation>
1179  ! RUNOFF <runoff>
1180  ! EXT-INFLOW <inflow>
1181  ! WITHDRAWAL <withdrawal>
1182  !
1183  found = .true.
1184  select case (keyword)
1185  case ('RAINFALL')
1186  ierr = this%apt_check_valid(itemno)
1187  if (ierr /= 0) then
1188  goto 999
1189  end if
1190  call this%parser%GetString(text)
1191  jj = 1
1192  bndelem => this%temprain(itemno)
1193  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1194  this%packName, 'BND', this%tsManager, &
1195  this%iprpak, 'RAINFALL')
1196  case ('EVAPORATION')
1197  ierr = this%apt_check_valid(itemno)
1198  if (ierr /= 0) then
1199  goto 999
1200  end if
1201  call this%parser%GetString(text)
1202  jj = 1
1203  bndelem => this%tempevap(itemno)
1204  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1205  this%packName, 'BND', this%tsManager, &
1206  this%iprpak, 'EVAPORATION')
1207  case ('RUNOFF')
1208  ierr = this%apt_check_valid(itemno)
1209  if (ierr /= 0) then
1210  goto 999
1211  end if
1212  call this%parser%GetString(text)
1213  jj = 1
1214  bndelem => this%temproff(itemno)
1215  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1216  this%packName, 'BND', this%tsManager, &
1217  this%iprpak, 'RUNOFF')
1218  case ('EXT-INFLOW')
1219  ierr = this%apt_check_valid(itemno)
1220  if (ierr /= 0) then
1221  goto 999
1222  end if
1223  call this%parser%GetString(text)
1224  jj = 1
1225  bndelem => this%tempiflw(itemno)
1226  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1227  this%packName, 'BND', this%tsManager, &
1228  this%iprpak, 'EXT-INFLOW')
1229  case default
1230  !
1231  ! -- Keyword not recognized so return to caller with found = .false.
1232  found = .false.
1233  end select
1234  !
1235 999 continue
1236  !
1237  ! -- Return
1238  return
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 499 of file gwe-lke.f90.

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

◆ lke_solve()

subroutine gwelkemodule::lke_solve ( class(gwelketype this)

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

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

◆ 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 943 of file gwe-lke.f90.

945  ! -- dummy
946  class(GweLkeType) :: this
947  integer(I4B), intent(in) :: ientry
948  integer(I4B), intent(inout) :: n1
949  integer(I4B), intent(inout) :: n2
950  real(DP), intent(inout), optional :: rrate
951  real(DP), intent(inout), optional :: rhsval
952  real(DP), intent(inout), optional :: hcofval
953  ! -- local
954  real(DP) :: qbnd
955  real(DP) :: ctmp
956  !
957  n1 = this%flowbudptr%budterm(this%idxbudwdrl)%id1(ientry)
958  n2 = this%flowbudptr%budterm(this%idxbudwdrl)%id2(ientry)
959  qbnd = this%flowbudptr%budterm(this%idxbudwdrl)%flow(ientry)
960  ctmp = this%xnewpak(n1)
961  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
962  if (present(rhsval)) rhsval = dzero
963  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac
964  !
965  ! -- Return
966  return

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'