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

Functions/Subroutines

real(dp) function, public etfunc_lin (efflndsrf, extdp, resid_pet, deriv_et, trhs, thcof, hgwf, celtop, celbot)
 Calculate gwf et using linear ET function from mf-2005. More...
 
pure real(dp) function, public calc_lin_scaling_fac (hgwf, lndsrf, extdp)
 Calculate the linear scaling factor. More...
 

Function/Subroutine Documentation

◆ calc_lin_scaling_fac()

pure real(dp) function, public uzfetutilmodule::calc_lin_scaling_fac ( real(dp), intent(in)  hgwf,
real(dp), intent(in)  lndsrf,
real(dp), intent(in)  extdp 
)
Parameters
[in]hgwfgroundwater head
[in]lndsrfeffective land surface (after applying surfdep to land surface)
[in]extdpextinction depth

Definition at line 75 of file UzfEtUtil.f90.

76  ! dummy
77  real(DP), intent(in) :: hgwf !< groundwater head
78  real(DP), intent(in) :: lndsrf !< effective land surface (after applying surfdep to land surface)
79  real(DP), intent(in) :: extdp !< extinction depth
80  ! return
81  real(DP) :: sclfac
82  !
83  sclfac = (hgwf - (lndsrf - extdp)) / extdp
84  ! the calculated scaling factor cannot exceed 1.0
85  if (sclfac > done) sclfac = done
Here is the caller graph for this function:

◆ etfunc_lin()

real(dp) function, public uzfetutilmodule::etfunc_lin ( real(dp), intent(in)  efflndsrf,
real(dp), intent(in)  extdp,
real(dp), intent(in)  resid_pet,
real(dp), intent(inout)  deriv_et,
real(dp), intent(inout)  trhs,
real(dp), intent(inout)  thcof,
real(dp), intent(in)  hgwf,
real(dp), intent(in)  celtop,
real(dp), intent(in)  celbot 
)
Parameters
[in]efflndsrfeffective land surface elevation after subtracting off 0.5*surfdep
[in]extdpextinction depth
[in]resid_petresidual pET remaining after applying actual ET from unsaturated zone
[in,out]deriv_etderivative of gw ET for Newton addition to equations in _fn()
[in,out]trhstotal uzf rhs contribution to GWF model
[in,out]thcoftotal uzf hcof contribution to GWF model
[in]hgwfcalculated groundwater head
[in]celtopelevation of the top of the cell
[in]celbotelevation of the bottom of the cell

Definition at line 14 of file UzfEtUtil.f90.

16  ! return
17  real(DP) :: etfunc_lin
18  ! dummy
19  real(DP), intent(in) :: efflndsrf !< effective land surface elevation after subtracting off 0.5*surfdep
20  real(DP), intent(in) :: extdp !< extinction depth
21  real(DP), intent(in) :: resid_pet !< residual pET remaining after applying actual ET from unsaturated zone
22  real(DP), intent(inout) :: deriv_et !< derivative of gw ET for Newton addition to equations in _fn()
23  real(DP), intent(inout) :: trhs !< total uzf rhs contribution to GWF model
24  real(DP), intent(inout) :: thcof !< total uzf hcof contribution to GWF model
25  real(DP), intent(in) :: hgwf !< calculated groundwater head
26  real(DP), intent(in) :: celtop !< elevation of the top of the cell
27  real(DP), intent(in) :: celbot !< elevation of the bottom of the cell
28  ! local
29  real(DP) :: etgw
30  real(DP) :: range
31  real(DP) :: depth, scale, thick, lin_scaling_fac
32  !
33  ! initialize
34  etgw = dzero
35  !
36  ! if water table between ET surface and extinction depth
37  ! (extdp starts at the bottom of the effective land surface elevation,
38  ! where the effective land surface elevation accounts for surfdep)
39  if (hgwf > (efflndsrf - extdp) .and. hgwf < efflndsrf) THEN
40  lin_scaling_fac = calc_lin_scaling_fac(hgwf, efflndsrf, extdp)
41  etgw = resid_pet * lin_scaling_fac
42  !
43  trhs = resid_pet - resid_pet * efflndsrf / extdp
44  thcof = -resid_pet / extdp
45  etgw = trhs - (thcof * hgwf)
46  !
47  ! water table above land surface
48  else if (hgwf >= efflndsrf) then
49  trhs = resid_pet
50  etgw = resid_pet
51  !
52  ! water table below extinction depth
53  else
54  etfunc_lin = dzero
55  return
56  end if
57  !
58  ! calculate rate
59  depth = hgwf - (efflndsrf - extdp)
60  thick = celtop - celbot
61  if (depth > thick) depth = thick
62  if (depth < dzero) depth = dzero
63  range = dem4 * extdp
64  call scubic(depth, range, deriv_et, scale)
65  trhs = scale * trhs
66  thcof = scale * thcof
67  etgw = trhs - (thcof * hgwf)
68  deriv_et = -deriv_et * etgw
69  etfunc_lin = etgw
70  !
Here is the call graph for this function:
Here is the caller graph for this function: