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

Data Types

type  gwesfetype
 

Functions/Subroutines

subroutine, public sfe_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
 Create a new sfe package. More...
 
subroutine find_sfe_package (this)
 Find corresponding sfe package. More...
 
subroutine sfe_fc_expanded (this, rhs, ia, idxglo, matrix_sln)
 Add matrix terms related to SFE. More...
 
subroutine sfe_solve (this)
 @ brief Add terms specific to sfr to the explicit sfe solve More...
 
integer(i4b) function sfe_get_nbudterms (this)
 Function to return the number of budget terms just for this package. More...
 
subroutine sfe_setup_budobj (this, idx)
 Set up the budget object that stores all the sfe flows. More...
 
subroutine sfe_fill_budobj (this, idx, x, flowja, ccratin, ccratout)
 Copy flow terms into thisbudobj. More...
 
subroutine allocate_scalars (this)
 Allocate scalars specific to the streamflow energy transport (SFE) package. More...
 
subroutine sfe_allocate_arrays (this)
 Allocate arrays specific to the streamflow energy transport (SFE) package. More...
 
subroutine sfe_ad (this)
 Advance sfe package routine. More...
 
subroutine sfe_da (this)
 Deallocate memory. More...
 
subroutine sfe_rain_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Rain term. More...
 
subroutine sfe_evap_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Evaporative term. More...
 
subroutine sfe_roff_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Runoff term. More...
 
subroutine sfe_iflw_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Inflow Term. More...
 
subroutine sfe_outf_term (this, ientry, n1, n2, rrate, rhsval, hcofval)
 Outflow term. More...
 
subroutine sfe_sbcd_term (this, ientry, n1, igwfnode, rrate, rhsval, hcofval)
 Streambed conduction term. More...
 
subroutine sfe_df_obs (this)
 Observations. More...
 
subroutine sfe_rp_obs (this, obsrv, found)
 Process package specific obs. More...
 
subroutine sfe_bd_obs (this, obstypeid, jj, v, found)
 Calculate observation value and pass it back to APT. More...
 
subroutine sfe_set_stressperiod (this, itemno, keyword, found)
 Sets the stress period attributes for keyword use. More...
 
subroutine sfe_read_cvs (this)
 Read feature information for this advanced package. More...
 
subroutine sfe_get_volumes (this, icv, vnew, vold, delt)
 Return the sfr new volume and old volume. More...
 

Variables

character(len= *), parameter ftype = 'SFE'
 
character(len= *), parameter flowtype = 'SFR'
 
character(len=16) text = ' SFE'
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine gwesfemodule::allocate_scalars ( class(gwesfetype this)
private

Definition at line 650 of file gwe-sfe.f90.

651  ! -- modules
653  ! -- dummy
654  class(GweSfeType) :: this
655  !
656  ! -- Allocate scalars in TspAptType
657  call this%TspAptType%allocate_scalars()
658  !
659  ! -- Allocate
660  call mem_allocate(this%idxbudrain, 'IDXBUDRAIN', this%memoryPath)
661  call mem_allocate(this%idxbudevap, 'IDXBUDEVAP', this%memoryPath)
662  call mem_allocate(this%idxbudroff, 'IDXBUDROFF', this%memoryPath)
663  call mem_allocate(this%idxbudiflw, 'IDXBUDIFLW', this%memoryPath)
664  call mem_allocate(this%idxbudoutf, 'IDXBUDOUTF', this%memoryPath)
665  !
666  ! -- Initialize
667  this%idxbudrain = 0
668  this%idxbudevap = 0
669  this%idxbudroff = 0
670  this%idxbudiflw = 0
671  this%idxbudoutf = 0

◆ find_sfe_package()

subroutine gwesfemodule::find_sfe_package ( class(gwesfetype this)
private

Definition at line 171 of file gwe-sfe.f90.

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

◆ sfe_ad()

subroutine gwesfemodule::sfe_ad ( class(gwesfetype this)

Definition at line 710 of file gwe-sfe.f90.

711  ! -- dummy
712  class(GweSfeType) :: this
713  ! -- local
714  integer(I4B) :: n
715  !
716  ! -- call base bnd_ad
717  call this%TspAptType%bnd_ad()
718  !
719  ! -- update vold
720  do n = 1, this%ncv
721  this%vold(n) = this%vnew(n)
722  end do
723 

◆ sfe_allocate_arrays()

subroutine gwesfemodule::sfe_allocate_arrays ( class(gwesfetype), intent(inout)  this)

Definition at line 677 of file gwe-sfe.f90.

678  ! -- modules
680  ! -- dummy
681  class(GweSfeType), intent(inout) :: this
682  ! -- local
683  integer(I4B) :: n
684  !
685  ! -- Time series
686  call mem_allocate(this%temprain, this%ncv, 'TEMPRAIN', this%memoryPath)
687  call mem_allocate(this%tempevap, this%ncv, 'TEMPEVAP', this%memoryPath)
688  call mem_allocate(this%temproff, this%ncv, 'TEMPROFF', this%memoryPath)
689  call mem_allocate(this%tempiflw, this%ncv, 'TEMPIFLW', this%memoryPath)
690  !
691  call mem_allocate(this%vnew, this%ncv, 'VNEW', this%memoryPath)
692  call mem_allocate(this%vold, this%ncv, 'VOLD', this%memoryPath)
693  !
694  ! -- Call standard TspAptType allocate arrays
695  call this%TspAptType%apt_allocate_arrays()
696  !
697  ! -- Initialize
698  do n = 1, this%ncv
699  this%temprain(n) = dzero
700  this%tempevap(n) = dzero
701  this%temproff(n) = dzero
702  this%tempiflw(n) = dzero
703  this%vnew(n) = dzero
704  this%vold(n) = dzero
705  end do

◆ sfe_bd_obs()

subroutine gwesfemodule::sfe_bd_obs ( class(gwesfetype), intent(inout)  this,
character(len=*), intent(in)  obstypeid,
integer(i4b), intent(in)  jj,
real(dp), intent(inout)  v,
logical, intent(inout)  found 
)
private

Definition at line 1037 of file gwe-sfe.f90.

1038  ! -- dummy
1039  class(GweSfeType), intent(inout) :: this
1040  character(len=*), intent(in) :: obstypeid
1041  real(DP), intent(inout) :: v
1042  integer(I4B), intent(in) :: jj
1043  logical, intent(inout) :: found
1044  ! -- local
1045  integer(I4B) :: n1, n2
1046  !
1047  found = .true.
1048  select case (obstypeid)
1049  case ('RAINFALL')
1050  if (this%iboundpak(jj) /= 0) then
1051  call this%sfe_rain_term(jj, n1, n2, v)
1052  end if
1053  case ('EVAPORATION')
1054  if (this%iboundpak(jj) /= 0) then
1055  call this%sfe_evap_term(jj, n1, n2, v)
1056  end if
1057  case ('RUNOFF')
1058  if (this%iboundpak(jj) /= 0) then
1059  call this%sfe_roff_term(jj, n1, n2, v)
1060  end if
1061  case ('EXT-INFLOW')
1062  if (this%iboundpak(jj) /= 0) then
1063  call this%sfe_iflw_term(jj, n1, n2, v)
1064  end if
1065  case ('EXT-OUTFLOW')
1066  if (this%iboundpak(jj) /= 0) then
1067  call this%sfe_outf_term(jj, n1, n2, v)
1068  end if
1069  case ('STRMBD-COND')
1070  if (this%iboundpak(jj) /= 0) then
1071  call this%sfe_sbcd_term(jj, n1, n2, v)
1072  end if
1073  case default
1074  found = .false.
1075  end select

◆ sfe_create()

subroutine, public gwesfemodule::sfe_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 110 of file gwe-sfe.f90.

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

◆ sfe_da()

subroutine gwesfemodule::sfe_da ( class(gwesfetype this)
private

Definition at line 728 of file gwe-sfe.f90.

729  ! -- modules
731  ! -- dummy
732  class(GweSfeType) :: this
733  !
734  ! -- Deallocate scalars
735  call mem_deallocate(this%idxbudrain)
736  call mem_deallocate(this%idxbudevap)
737  call mem_deallocate(this%idxbudroff)
738  call mem_deallocate(this%idxbudiflw)
739  call mem_deallocate(this%idxbudoutf)
740  !
741  ! -- Deallocate time series
742  call mem_deallocate(this%temprain)
743  call mem_deallocate(this%tempevap)
744  call mem_deallocate(this%temproff)
745  call mem_deallocate(this%tempiflw)
746  !
747  call mem_deallocate(this%vnew)
748  call mem_deallocate(this%vold)
749  !
750  ! -- Deallocate arrays
751  call mem_deallocate(this%ktf)
752  call mem_deallocate(this%rfeatthk)
753  !
754  ! -- Deallocate scalars in TspAptType
755  call this%TspAptType%bnd_da()

◆ sfe_df_obs()

subroutine gwesfemodule::sfe_df_obs ( class(gwesfetype this)
private

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

Definition at line 930 of file gwe-sfe.f90.

931  ! -- modules
932  ! -- dummy
933  class(GweSfeType) :: this
934  ! -- local
935  integer(I4B) :: indx
936  !
937  ! -- Store obs type and assign procedure pointer
938  ! for temperature observation type.
939  call this%obs%StoreObsType('temperature', .false., indx)
940  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
941  !
942  ! -- Store obs type and assign procedure pointer
943  ! for flow between reaches.
944  call this%obs%StoreObsType('flow-ja-face', .true., indx)
945  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid12
946  !
947  ! -- Store obs type and assign procedure pointer
948  ! for from-mvr observation type.
949  call this%obs%StoreObsType('from-mvr', .true., indx)
950  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
951  !
952  ! -- Store obs type and assign procedure pointer
953  ! for to-mvr observation type.
954  call this%obs%StoreObsType('to-mvr', .true., indx)
955  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
956  !
957  ! -- Store obs type and assign procedure pointer
958  ! for storage observation type.
959  call this%obs%StoreObsType('storage', .true., indx)
960  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
961  !
962  ! -- Store obs type and assign procedure pointer
963  ! for constant observation type.
964  call this%obs%StoreObsType('constant', .true., indx)
965  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
966  !
967  ! -- Store obs type and assign procedure pointer
968  ! for observation type: sfe
969  call this%obs%StoreObsType('sfe', .true., indx)
970  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
971  !
972  ! -- Store obs type and assign procedure pointer
973  ! for rainfall observation type.
974  call this%obs%StoreObsType('rainfall', .true., indx)
975  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
976  !
977  ! -- Store obs type and assign procedure pointer
978  ! for evaporation observation type.
979  call this%obs%StoreObsType('evaporation', .true., indx)
980  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
981  !
982  ! -- Store obs type and assign procedure pointer
983  ! for runoff observation type.
984  call this%obs%StoreObsType('runoff', .true., indx)
985  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
986  !
987  ! -- Store obs type and assign procedure pointer
988  ! for inflow observation type.
989  call this%obs%StoreObsType('ext-inflow', .true., indx)
990  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
991  !
992  ! -- Store obs type and assign procedure pointer
993  ! for ext-outflow observation type.
994  call this%obs%StoreObsType('ext-outflow', .true., indx)
995  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
996  !
997  ! -- Store obs type and assign procedure pointer
998  ! for strmbd-cond observation type.
999  call this%obs%StoreObsType('strmbd-cond', .true., indx)
1000  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
Here is the call graph for this function:

◆ sfe_evap_term()

subroutine gwesfemodule::sfe_evap_term ( class(gwesfetype 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 
)
private

Definition at line 784 of file gwe-sfe.f90.

785  ! -- dummy
786  class(GweSfeType) :: this
787  integer(I4B), intent(in) :: ientry
788  integer(I4B), intent(inout) :: n1
789  integer(I4B), intent(inout) :: n2
790  real(DP), intent(inout), optional :: rrate
791  real(DP), intent(inout), optional :: rhsval
792  real(DP), intent(inout), optional :: hcofval
793  ! -- local
794  real(DP) :: qbnd
795  real(DP) :: heatlat
796  !
797  n1 = this%flowbudptr%budterm(this%idxbudevap)%id1(ientry)
798  n2 = this%flowbudptr%budterm(this%idxbudevap)%id2(ientry)
799  ! -- note that qbnd is negative for evap
800  qbnd = this%flowbudptr%budterm(this%idxbudevap)%flow(ientry)
801  heatlat = this%gwecommon%gwerhow * this%gwecommon%gwelatheatvap
802  if (present(rrate)) rrate = qbnd * heatlat
803  if (present(rhsval)) rhsval = -rrate
804  if (present(hcofval)) hcofval = dzero

◆ sfe_fc_expanded()

subroutine gwesfemodule::sfe_fc_expanded ( class(gwesfetype 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 SFE

Definition at line 290 of file gwe-sfe.f90.

291  ! -- dummy
292  class(GweSfeType) :: this
293  real(DP), dimension(:), intent(inout) :: rhs
294  integer(I4B), dimension(:), intent(in) :: ia
295  integer(I4B), dimension(:), intent(in) :: idxglo
296  class(MatrixBaseType), pointer :: matrix_sln
297  ! -- local
298  integer(I4B) :: j, n1, n2
299  integer(I4B) :: iloc
300  integer(I4B) :: iposd, iposoffd
301  integer(I4B) :: ipossymd, ipossymoffd
302  real(DP) :: rrate
303  real(DP) :: rhsval
304  real(DP) :: hcofval
305  !
306  ! -- Add rainfall contribution
307  if (this%idxbudrain /= 0) then
308  do j = 1, this%flowbudptr%budterm(this%idxbudrain)%nlist
309  call this%sfe_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%sfe_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%sfe_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%sfe_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 outflow contribution
351  if (this%idxbudoutf /= 0) then
352  do j = 1, this%flowbudptr%budterm(this%idxbudoutf)%nlist
353  call this%sfe_outf_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 streambed conduction contribution
362  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
363  !
364  ! -- Set n to feature number and process if active feature
365  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
366  if (this%iboundpak(n1) /= 0) then
367  !
368  call this%sfe_sbcd_term(j, n1, n2, rrate, rhsval, hcofval)
369  !
370  ! -- Add to sfe row
371  iposd = this%idxdglo(j)
372  iposoffd = this%idxoffdglo(j)
373  call matrix_sln%add_value_pos(iposd, -hcofval)
374  call matrix_sln%add_value_pos(iposoffd, hcofval)
375  !
376  ! -- Add to gwe row for sfe connection
377  ipossymd = this%idxsymdglo(j)
378  ipossymoffd = this%idxsymoffdglo(j)
379  call matrix_sln%add_value_pos(ipossymd, -hcofval)
380  call matrix_sln%add_value_pos(ipossymoffd, hcofval)
381  end if
382  end do

◆ sfe_fill_budobj()

subroutine gwesfemodule::sfe_fill_budobj ( class(gwesfetype 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 560 of file gwe-sfe.f90.

561  ! -- dummy
562  class(GweSfeType) :: this
563  integer(I4B), intent(inout) :: idx
564  real(DP), dimension(:), intent(in) :: x
565  real(DP), dimension(:), contiguous, intent(inout) :: flowja
566  real(DP), intent(inout) :: ccratin
567  real(DP), intent(inout) :: ccratout
568  ! -- local
569  integer(I4B) :: j, n1, n2
570  integer(I4B) :: igwfnode
571  integer(I4B) :: nlist
572  integer(I4B) :: idiag
573  real(DP) :: q
574  !
575  ! -- Rain
576  idx = idx + 1
577  nlist = this%flowbudptr%budterm(this%idxbudrain)%nlist
578  call this%budobj%budterm(idx)%reset(nlist)
579  do j = 1, nlist
580  call this%sfe_rain_term(j, n1, n2, q)
581  call this%budobj%budterm(idx)%update_term(n1, n2, q)
582  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
583  end do
584  !
585  ! -- Evaporation
586  idx = idx + 1
587  nlist = this%flowbudptr%budterm(this%idxbudevap)%nlist
588  call this%budobj%budterm(idx)%reset(nlist)
589  do j = 1, nlist
590  call this%sfe_evap_term(j, n1, n2, q)
591  call this%budobj%budterm(idx)%update_term(n1, n2, q)
592  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
593  end do
594  !
595  ! -- Runoff
596  idx = idx + 1
597  nlist = this%flowbudptr%budterm(this%idxbudroff)%nlist
598  call this%budobj%budterm(idx)%reset(nlist)
599  do j = 1, nlist
600  call this%sfe_roff_term(j, n1, n2, q)
601  call this%budobj%budterm(idx)%update_term(n1, n2, q)
602  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
603  end do
604  !
605  ! -- Ext-inflow
606  idx = idx + 1
607  nlist = this%flowbudptr%budterm(this%idxbudiflw)%nlist
608  call this%budobj%budterm(idx)%reset(nlist)
609  do j = 1, nlist
610  call this%sfe_iflw_term(j, n1, n2, q)
611  call this%budobj%budterm(idx)%update_term(n1, n2, q)
612  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
613  end do
614  !
615  ! -- Ext-outflow
616  idx = idx + 1
617  nlist = this%flowbudptr%budterm(this%idxbudoutf)%nlist
618  call this%budobj%budterm(idx)%reset(nlist)
619  do j = 1, nlist
620  call this%sfe_outf_term(j, n1, n2, q)
621  call this%budobj%budterm(idx)%update_term(n1, n2, q)
622  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
623  end do
624  !
625  ! -- Strmbd-cond. The idxbudgwf index is used since it contains the
626  ! sfe/cell mapping information
627  idx = idx + 1
628  nlist = this%flowbudptr%budterm(this%idxbudgwf)%nlist
629  call this%budobj%budterm(idx)%reset(nlist)
630  do j = 1, nlist
631  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
632  if (this%iboundpak(n1) /= 0) then
633  ! -- use igwfnode instead of n2 for consistency with usage in apt;
634  ! helps highlight that cell number is sought and used
635  call this%sfe_sbcd_term(j, n1, igwfnode, q)
636  call this%budobj%budterm(idx)%update_term(n1, igwfnode, q)
637  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
638  ! -- Contribution to gwe cell budget, which is different than the
639  ! other terms (i.e., ext-inflow, evap, etc.)
640  this%simvals(n1) = this%simvals(n1) - q
641  idiag = this%dis%con%ia(igwfnode)
642  flowja(idiag) = flowja(idiag) - q
643  end if
644  end do

◆ sfe_get_nbudterms()

integer(i4b) function gwesfemodule::sfe_get_nbudterms ( class(gwesfetype this)
private

This overrides a function in the parent class.

Definition at line 442 of file gwe-sfe.f90.

443  ! -- dummy
444  class(GweSfeType) :: this
445  ! -- return
446  integer(I4B) :: nbudterms
447  !
448  ! -- Number of budget terms is 6:
449  ! 1. rainfall
450  ! 2. evaporation
451  ! 3. runoff
452  ! 4. ext-inflow
453  ! 5. ext-outflow
454  ! 6. strmbd-cond
455  nbudterms = 6

◆ sfe_get_volumes()

subroutine gwesfemodule::sfe_get_volumes ( class(gwesfetype this,
integer(i4b), intent(in)  icv,
real(dp), intent(inout)  vnew,
real(dp), intent(inout)  vold,
real(dp), intent(in)  delt 
)

Definition at line 1329 of file gwe-sfe.f90.

1330  ! -- dummy
1331  class(GweSfeType) :: this
1332  integer(I4B), intent(in) :: icv
1333  real(DP), intent(inout) :: vnew, vold
1334  real(DP), intent(in) :: delt
1335  ! -- local
1336  real(DP) :: qss
1337  !
1338  ! -- get volumes
1339  vold = dzero
1340  vnew = vold
1341  if (this%idxbudsto /= 0) then
1342  qss = this%flowbudptr%budterm(this%idxbudsto)%flow(icv)
1343  vnew = this%flowbudptr%budterm(this%idxbudsto)%auxvar(1, icv)
1344  this%vnew(icv) = vnew
1345  if (qss /= dzero) then
1346  vold = vnew + qss * delt
1347  else
1348  if (vnew == dzero) then
1349  vold = dzero
1350  else
1351  vold = this%vold(icv)
1352  end if
1353  end if
1354  end if

◆ sfe_iflw_term()

subroutine gwesfemodule::sfe_iflw_term ( class(gwesfetype 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 
)
private

Accounts for energy added via streamflow entering into a stream channel; for example, energy entering the model domain via a specified flow in a stream channel.

Definition at line 837 of file gwe-sfe.f90.

838  ! -- dummy
839  class(GweSfeType) :: this
840  integer(I4B), intent(in) :: ientry
841  integer(I4B), intent(inout) :: n1
842  integer(I4B), intent(inout) :: n2
843  real(DP), intent(inout), optional :: rrate
844  real(DP), intent(inout), optional :: rhsval
845  real(DP), intent(inout), optional :: hcofval
846  ! -- local
847  real(DP) :: qbnd
848  real(DP) :: ctmp
849  !
850  n1 = this%flowbudptr%budterm(this%idxbudiflw)%id1(ientry)
851  n2 = this%flowbudptr%budterm(this%idxbudiflw)%id2(ientry)
852  qbnd = this%flowbudptr%budterm(this%idxbudiflw)%flow(ientry)
853  ctmp = this%tempiflw(n1)
854  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
855  if (present(rhsval)) rhsval = -rrate
856  if (present(hcofval)) hcofval = dzero

◆ sfe_outf_term()

subroutine gwesfemodule::sfe_outf_term ( class(gwesfetype 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 
)
private

Accounts for the energy leaving a stream channel, for example, energy exiting the model domain via a flow in a stream channel flowing out of the active domain.

Definition at line 864 of file gwe-sfe.f90.

865  ! -- dummy
866  class(GweSfeType) :: this
867  integer(I4B), intent(in) :: ientry
868  integer(I4B), intent(inout) :: n1
869  integer(I4B), intent(inout) :: n2
870  real(DP), intent(inout), optional :: rrate
871  real(DP), intent(inout), optional :: rhsval
872  real(DP), intent(inout), optional :: hcofval
873  ! -- local
874  real(DP) :: qbnd
875  real(DP) :: ctmp
876  !
877  n1 = this%flowbudptr%budterm(this%idxbudoutf)%id1(ientry)
878  n2 = this%flowbudptr%budterm(this%idxbudoutf)%id2(ientry)
879  qbnd = this%flowbudptr%budterm(this%idxbudoutf)%flow(ientry)
880  ctmp = this%xnewpak(n1)
881  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
882  if (present(rhsval)) rhsval = dzero
883  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac

◆ sfe_rain_term()

subroutine gwesfemodule::sfe_rain_term ( class(gwesfetype 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 760 of file gwe-sfe.f90.

761  ! -- dummy
762  class(GweSfeType) :: this
763  integer(I4B), intent(in) :: ientry
764  integer(I4B), intent(inout) :: n1
765  integer(I4B), intent(inout) :: n2
766  real(DP), intent(inout), optional :: rrate
767  real(DP), intent(inout), optional :: rhsval
768  real(DP), intent(inout), optional :: hcofval
769  ! -- local
770  real(DP) :: qbnd
771  real(DP) :: ctmp
772  !
773  n1 = this%flowbudptr%budterm(this%idxbudrain)%id1(ientry)
774  n2 = this%flowbudptr%budterm(this%idxbudrain)%id2(ientry)
775  qbnd = this%flowbudptr%budterm(this%idxbudrain)%flow(ientry)
776  ctmp = this%temprain(n1)
777  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
778  if (present(rhsval)) rhsval = -rrate
779  if (present(hcofval)) hcofval = dzero

◆ sfe_read_cvs()

subroutine gwesfemodule::sfe_read_cvs ( class(gwesfetype), intent(inout)  this)

Definition at line 1157 of file gwe-sfe.f90.

1158  ! -- modules
1161  ! -- dummy
1162  class(GweSfeType), intent(inout) :: this
1163  ! -- local
1164  character(len=LINELENGTH) :: text
1165  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
1166  character(len=9) :: cno
1167  character(len=50), dimension(:), allocatable :: caux
1168  integer(I4B) :: ierr
1169  logical :: isfound, endOfBlock
1170  integer(I4B) :: n
1171  integer(I4B) :: ii, jj
1172  integer(I4B) :: iaux
1173  integer(I4B) :: itmp
1174  integer(I4B) :: nlak
1175  integer(I4B) :: nconn
1176  integer(I4B), dimension(:), pointer, contiguous :: nboundchk
1177  real(DP), pointer :: bndElem => null()
1178  !
1179  ! -- initialize itmp
1180  itmp = 0
1181  !
1182  ! -- allocate apt data
1183  call mem_allocate(this%strt, this%ncv, 'STRT', this%memoryPath)
1184  call mem_allocate(this%ktf, this%ncv, 'KTF', this%memoryPath)
1185  call mem_allocate(this%rfeatthk, this%ncv, 'RFEATTHK', this%memoryPath)
1186  call mem_allocate(this%lauxvar, this%naux, this%ncv, 'LAUXVAR', &
1187  this%memoryPath)
1188  !
1189  ! -- stream boundary and temperatures
1190  if (this%imatrows == 0) then
1191  call mem_allocate(this%iboundpak, this%ncv, 'IBOUND', this%memoryPath)
1192  call mem_allocate(this%xnewpak, this%ncv, 'XNEWPAK', this%memoryPath)
1193  end if
1194  call mem_allocate(this%xoldpak, this%ncv, 'XOLDPAK', this%memoryPath)
1195  !
1196  ! -- allocate character storage not managed by the memory manager
1197  allocate (this%featname(this%ncv)) ! ditch after boundnames allocated??
1198  !allocate(this%status(this%ncv))
1199  !
1200  do n = 1, this%ncv
1201  this%strt(n) = dep20
1202  this%ktf(n) = dzero
1203  this%rfeatthk(n) = dzero
1204  this%lauxvar(:, n) = dzero
1205  this%xoldpak(n) = dep20
1206  if (this%imatrows == 0) then
1207  this%iboundpak(n) = 1
1208  this%xnewpak(n) = dep20
1209  end if
1210  end do
1211  !
1212  ! -- allocate local storage for aux variables
1213  if (this%naux > 0) then
1214  allocate (caux(this%naux))
1215  end if
1216  !
1217  ! -- allocate and initialize temporary variables
1218  allocate (nboundchk(this%ncv))
1219  do n = 1, this%ncv
1220  nboundchk(n) = 0
1221  end do
1222  !
1223  ! -- get packagedata block
1224  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
1225  supportopenclose=.true.)
1226  !
1227  ! -- parse locations block if detected
1228  if (isfound) then
1229  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text))// &
1230  ' PACKAGEDATA'
1231  nlak = 0
1232  nconn = 0
1233  do
1234  call this%parser%GetNextLine(endofblock)
1235  if (endofblock) exit
1236  n = this%parser%GetInteger()
1237 
1238  if (n < 1 .or. n > this%ncv) then
1239  write (errmsg, '(a,1x,i6)') &
1240  'Itemno must be > 0 and <= ', this%ncv
1241  call store_error(errmsg)
1242  cycle
1243  end if
1244  !
1245  ! -- increment nboundchk
1246  nboundchk(n) = nboundchk(n) + 1
1247  !
1248  ! -- strt
1249  this%strt(n) = this%parser%GetDouble()
1250  !
1251  ! -- read additional thermal conductivity terms
1252  this%ktf(n) = this%parser%GetDouble()
1253  this%rfeatthk(n) = this%parser%GetDouble()
1254  if (this%rfeatthk(n) <= dzero) then
1255  write (errmsg, '(4x,a)') &
1256  '****ERROR. Specified thickness used for thermal &
1257  &conduction MUST BE > 0 else divide by zero error occurs'
1258  call store_error(errmsg)
1259  cycle
1260  end if
1261  !
1262  ! -- get aux data
1263  do iaux = 1, this%naux
1264  call this%parser%GetString(caux(iaux))
1265  end do
1266 
1267  ! -- set default bndName
1268  write (cno, '(i9.9)') n
1269  bndname = 'Feature'//cno
1270 
1271  ! -- featname
1272  if (this%inamedbound /= 0) then
1273  call this%parser%GetStringCaps(bndnametemp)
1274  if (bndnametemp /= '') then
1275  bndname = bndnametemp
1276  end if
1277  end if
1278  this%featname(n) = bndname
1279 
1280  ! -- fill time series aware data
1281  ! -- fill aux data
1282  do jj = 1, this%naux
1283  text = caux(jj)
1284  ii = n
1285  bndelem => this%lauxvar(jj, ii)
1286  call read_value_or_time_series_adv(text, ii, jj, bndelem, &
1287  this%packName, 'AUX', &
1288  this%tsManager, this%iprpak, &
1289  this%auxname(jj))
1290  end do
1291  !
1292  nlak = nlak + 1
1293  end do
1294  !
1295  ! -- check for duplicate or missing lakes
1296  do n = 1, this%ncv
1297  if (nboundchk(n) == 0) then
1298  write (errmsg, '(a,1x,i0)') 'No data specified for feature', n
1299  call store_error(errmsg)
1300  else if (nboundchk(n) > 1) then
1301  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
1302  'Data for feature', n, 'specified', nboundchk(n), 'times'
1303  call store_error(errmsg)
1304  end if
1305  end do
1306  !
1307  write (this%iout, '(1x,a)') &
1308  'END OF '//trim(adjustl(this%text))//' PACKAGEDATA'
1309  else
1310  call store_error('Required packagedata block not found.')
1311  end if
1312  !
1313  ! -- terminate if any errors were detected
1314  if (count_errors() > 0) then
1315  call this%parser%StoreErrorUnit()
1316  end if
1317  !
1318  ! -- deallocate local storage for aux variables
1319  if (this%naux > 0) then
1320  deallocate (caux)
1321  end if
1322  !
1323  ! -- deallocate local storage for nboundchk
1324  deallocate (nboundchk)
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:

◆ sfe_roff_term()

subroutine gwesfemodule::sfe_roff_term ( class(gwesfetype 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 
)
private

Definition at line 809 of file gwe-sfe.f90.

810  ! -- dummy
811  class(GweSfeType) :: this
812  integer(I4B), intent(in) :: ientry
813  integer(I4B), intent(inout) :: n1
814  integer(I4B), intent(inout) :: n2
815  real(DP), intent(inout), optional :: rrate
816  real(DP), intent(inout), optional :: rhsval
817  real(DP), intent(inout), optional :: hcofval
818  ! -- local
819  real(DP) :: qbnd
820  real(DP) :: ctmp
821  !
822  n1 = this%flowbudptr%budterm(this%idxbudroff)%id1(ientry)
823  n2 = this%flowbudptr%budterm(this%idxbudroff)%id2(ientry)
824  qbnd = this%flowbudptr%budterm(this%idxbudroff)%flow(ientry)
825  ctmp = this%temproff(n1)
826  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
827  if (present(rhsval)) rhsval = -rrate
828  if (present(hcofval)) hcofval = dzero

◆ sfe_rp_obs()

subroutine gwesfemodule::sfe_rp_obs ( class(gwesfetype), intent(inout)  this,
type(observetype), intent(inout)  obsrv,
logical, intent(inout)  found 
)
private

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 1007 of file gwe-sfe.f90.

1008  ! -- dummy
1009  class(GweSfeType), intent(inout) :: this !< package class
1010  type(ObserveType), intent(inout) :: obsrv !< observation object
1011  logical, intent(inout) :: found !< indicate whether observation was found
1012  ! -- local
1013  !
1014  found = .true.
1015  select case (obsrv%ObsTypeId)
1016  case ('RAINFALL')
1017  call this%rp_obs_byfeature(obsrv)
1018  case ('EVAPORATION')
1019  call this%rp_obs_byfeature(obsrv)
1020  case ('RUNOFF')
1021  call this%rp_obs_byfeature(obsrv)
1022  case ('EXT-INFLOW')
1023  call this%rp_obs_byfeature(obsrv)
1024  case ('EXT-OUTFLOW')
1025  call this%rp_obs_byfeature(obsrv)
1026  case ('TO-MVR')
1027  call this%rp_obs_byfeature(obsrv)
1028  case ('STRMBD-COND')
1029  call this%rp_obs_byfeature(obsrv)
1030  case default
1031  found = .false.
1032  end select

◆ sfe_sbcd_term()

subroutine gwesfemodule::sfe_sbcd_term ( class(gwesfetype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(inout)  n1,
integer(i4b), intent(inout)  igwfnode,
real(dp), intent(inout), optional  rrate,
real(dp), intent(inout), optional  rhsval,
real(dp), intent(inout), optional  hcofval 
)
private

Accounts for the energy entering or leaving a stream channel by thermal conduction through the streambed.

Definition at line 891 of file gwe-sfe.f90.

892  ! -- dummy
893  class(GweSfeType) :: this
894  integer(I4B), intent(in) :: ientry
895  integer(I4B), intent(inout) :: n1
896  integer(I4B), intent(inout) :: igwfnode
897  real(DP), intent(inout), optional :: rrate
898  real(DP), intent(inout), optional :: rhsval
899  real(DP), intent(inout), optional :: hcofval
900  ! -- local
901  integer(I4B) :: auxpos
902  real(DP) :: wa !< wetted area
903  real(DP) :: ktf !< thermal conductivity of streambed material
904  real(DP) :: s !< thickness of conductive streambed material
905  real(DP) :: ctherm
906  !
907  rrate = dzero
908  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(ientry)
909  ! -- use igwfnode instead of n2 for consistency with usage in apt;
910  ! helps highlight that cell number is sought and used
911  igwfnode = this%flowbudptr%budterm(this%idxbudgwf)%id2(ientry)
912  ! -- For now, there is only 1 aux variable under 'GWF'
913  auxpos = this%flowbudptr%budterm(this%idxbudgwf)%naux
914  wa = this%flowbudptr%budterm(this%idxbudgwf)%auxvar(auxpos, ientry)
915  ktf = this%ktf(n1)
916  s = this%rfeatthk(n1)
917  ctherm = ktf * wa / s
918  ! -- this%xnew available b/c set in parent class (TspAptType) using
919  ! routine set_pointers from the "grandparent" class BndType
920  if (present(rrate)) rrate = ctherm * (this%xnew(igwfnode) - this%xnewpak(n1))
921  if (present(rhsval)) rhsval = dzero
922  if (present(hcofval)) hcofval = ctherm

◆ sfe_set_stressperiod()

subroutine gwesfemodule::sfe_set_stressperiod ( class(gwesfetype), intent(inout)  this,
integer(i4b), intent(in)  itemno,
character(len=*), intent(in)  keyword,
logical, intent(inout)  found 
)
private

Definition at line 1080 of file gwe-sfe.f90.

1081  ! -- modules
1083  ! -- dummy
1084  class(GweSfeType), intent(inout) :: this
1085  integer(I4B), intent(in) :: itemno
1086  character(len=*), intent(in) :: keyword
1087  logical, intent(inout) :: found
1088  ! -- local
1089  character(len=LINELENGTH) :: text
1090  integer(I4B) :: ierr
1091  integer(I4B) :: jj
1092  real(DP), pointer :: bndElem => null()
1093  !
1094  ! RAINFALL <rainfall>
1095  ! EVAPORATION <evaporation>
1096  ! RUNOFF <runoff>
1097  ! INFLOW <inflow>
1098  ! WITHDRAWAL <withdrawal>
1099  !
1100  found = .true.
1101  select case (keyword)
1102  case ('RAINFALL')
1103  ierr = this%apt_check_valid(itemno)
1104  if (ierr /= 0) then
1105  goto 999
1106  end if
1107  call this%parser%GetString(text)
1108  jj = 1
1109  bndelem => this%temprain(itemno)
1110  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1111  this%packName, 'BND', this%tsManager, &
1112  this%iprpak, 'RAINFALL')
1113  case ('EVAPORATION')
1114  ierr = this%apt_check_valid(itemno)
1115  if (ierr /= 0) then
1116  goto 999
1117  end if
1118  call this%parser%GetString(text)
1119  jj = 1
1120  bndelem => this%tempevap(itemno)
1121  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1122  this%packName, 'BND', this%tsManager, &
1123  this%iprpak, 'EVAPORATION')
1124  case ('RUNOFF')
1125  ierr = this%apt_check_valid(itemno)
1126  if (ierr /= 0) then
1127  goto 999
1128  end if
1129  call this%parser%GetString(text)
1130  jj = 1
1131  bndelem => this%temproff(itemno)
1132  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1133  this%packName, 'BND', this%tsManager, &
1134  this%iprpak, 'RUNOFF')
1135  case ('INFLOW')
1136  ierr = this%apt_check_valid(itemno)
1137  if (ierr /= 0) then
1138  goto 999
1139  end if
1140  call this%parser%GetString(text)
1141  jj = 1
1142  bndelem => this%tempiflw(itemno)
1143  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
1144  this%packName, 'BND', this%tsManager, &
1145  this%iprpak, 'INFLOW')
1146  case default
1147  !
1148  ! -- Keyword not recognized so return to caller with found = .false.
1149  found = .false.
1150  end select
1151  !
1152 999 continue
Here is the call graph for this function:

◆ sfe_setup_budobj()

subroutine gwesfemodule::sfe_setup_budobj ( class(gwesfetype this,
integer(i4b), intent(inout)  idx 
)
private

Definition at line 460 of file gwe-sfe.f90.

461  ! -- modules
462  use constantsmodule, only: lenbudtxt
463  ! -- dummy
464  class(GweSfeType) :: this
465  integer(I4B), intent(inout) :: idx
466  ! -- local
467  integer(I4B) :: n, n1, n2
468  integer(I4B) :: maxlist, naux
469  real(DP) :: q
470  character(len=LENBUDTXT) :: text
471  !
472  ! --
473  text = ' RAINFALL'
474  idx = idx + 1
475  maxlist = this%flowbudptr%budterm(this%idxbudrain)%maxlist
476  naux = 0
477  call this%budobj%budterm(idx)%initialize(text, &
478  this%name_model, &
479  this%packName, &
480  this%name_model, &
481  this%packName, &
482  maxlist, .false., .false., &
483  naux)
484  !
485  ! --
486  text = ' EVAPORATION'
487  idx = idx + 1
488  maxlist = this%flowbudptr%budterm(this%idxbudevap)%maxlist
489  naux = 0
490  call this%budobj%budterm(idx)%initialize(text, &
491  this%name_model, &
492  this%packName, &
493  this%name_model, &
494  this%packName, &
495  maxlist, .false., .false., &
496  naux)
497  !
498  ! --
499  text = ' RUNOFF'
500  idx = idx + 1
501  maxlist = this%flowbudptr%budterm(this%idxbudroff)%maxlist
502  naux = 0
503  call this%budobj%budterm(idx)%initialize(text, &
504  this%name_model, &
505  this%packName, &
506  this%name_model, &
507  this%packName, &
508  maxlist, .false., .false., &
509  naux)
510  !
511  ! --
512  text = ' EXT-INFLOW'
513  idx = idx + 1
514  maxlist = this%flowbudptr%budterm(this%idxbudiflw)%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  ! --
525  text = ' EXT-OUTFLOW'
526  idx = idx + 1
527  maxlist = this%flowbudptr%budterm(this%idxbudoutf)%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  ! -- Conduction through the wetted streambed
538  text = ' STRMBD-COND'
539  idx = idx + 1
540  maxlist = this%flowbudptr%budterm(this%idxbudgwf)%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  call this%budobj%budterm(idx)%reset(maxlist)
550  q = dzero
551  do n = 1, maxlist
552  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(n)
553  n2 = this%flowbudptr%budterm(this%idxbudgwf)%id2(n)
554  call this%budobj%budterm(idx)%update_term(n1, n2, q)
555  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

◆ sfe_solve()

subroutine gwesfemodule::sfe_solve ( class(gwesfetype this)
private

Definition at line 387 of file gwe-sfe.f90.

388  ! -- dummy
389  class(GweSfeType) :: this
390  ! -- local
391  integer(I4B) :: j
392  integer(I4B) :: n1, n2
393  real(DP) :: rrate
394  !
395  ! -- Add rainfall contribution
396  if (this%idxbudrain /= 0) then
397  do j = 1, this%flowbudptr%budterm(this%idxbudrain)%nlist
398  call this%sfe_rain_term(j, n1, n2, rrate)
399  this%dbuff(n1) = this%dbuff(n1) + rrate
400  end do
401  end if
402  !
403  ! -- Add evaporation contribution
404  if (this%idxbudevap /= 0) then
405  do j = 1, this%flowbudptr%budterm(this%idxbudevap)%nlist
406  call this%sfe_evap_term(j, n1, n2, rrate)
407  this%dbuff(n1) = this%dbuff(n1) + rrate
408  end do
409  end if
410  !
411  ! -- Add runoff contribution
412  if (this%idxbudroff /= 0) then
413  do j = 1, this%flowbudptr%budterm(this%idxbudroff)%nlist
414  call this%sfe_roff_term(j, n1, n2, rrate)
415  this%dbuff(n1) = this%dbuff(n1) + rrate
416  end do
417  end if
418  !
419  ! -- Add inflow contribution
420  if (this%idxbudiflw /= 0) then
421  do j = 1, this%flowbudptr%budterm(this%idxbudiflw)%nlist
422  call this%sfe_iflw_term(j, n1, n2, rrate)
423  this%dbuff(n1) = this%dbuff(n1) + rrate
424  end do
425  end if
426  !
427  ! -- Add outflow contribution
428  if (this%idxbudoutf /= 0) then
429  do j = 1, this%flowbudptr%budterm(this%idxbudoutf)%nlist
430  call this%sfe_outf_term(j, n1, n2, rrate)
431  this%dbuff(n1) = this%dbuff(n1) + rrate
432  end do
433  end if
434  !
435  ! Note: explicit streambed conduction terms???

Variable Documentation

◆ flowtype

character(len=*), parameter gwesfemodule::flowtype = 'SFR'
private

Definition at line 57 of file gwe-sfe.f90.

57  character(len=*), parameter :: flowtype = 'SFR'

◆ ftype

character(len=*), parameter gwesfemodule::ftype = 'SFE'
private

Definition at line 56 of file gwe-sfe.f90.

56  character(len=*), parameter :: ftype = 'SFE'

◆ text

character(len=16) gwesfemodule::text = ' SFE'
private

Definition at line 58 of file gwe-sfe.f90.

58  character(len=16) :: text = ' SFE'