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 decay ET function from mf-2005. More...
 
real(dp) function, public etfunc_nlin (efflndsrf, extdp, resid_pet, deriv_et, trhs, thcof, hgwf)
 Calculate gwf ET using a square decay ET function with smoothing at the specified extinction depth. 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 107 of file UzfEtUtil.f90.

108  ! dummy
109  real(DP), intent(in) :: hgwf !< groundwater head
110  real(DP), intent(in) :: lndsrf !< effective land surface (after applying surfdep to land surface)
111  real(DP), intent(in) :: extdp !< extinction depth
112  ! return
113  real(DP) :: sclfac
114  !
115  sclfac = (hgwf - (lndsrf - extdp)) / extdp
116  ! the calculated scaling factor cannot exceed 1.0
117  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 subtracting simulated ET in the 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 16 of file UzfEtUtil.f90.

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

◆ etfunc_nlin()

real(dp) function, public uzfetutilmodule::etfunc_nlin ( 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 
)
Parameters
[in]efflndsrfeffective land surface elevation after subtracting off 0.5*surfdep
[in]extdpextinction depth
[in]resid_petresidual pET remaining after subtracting simulated ET in the 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

Definition at line 78 of file UzfEtUtil.f90.

79  ! -- return
80  real(DP) :: etfunc_nlin
81  ! -- dummy
82  real(DP), intent(in) :: efflndsrf !< effective land surface elevation after subtracting off 0.5*surfdep
83  real(DP), intent(in) :: extdp !< extinction depth
84  real(DP), intent(in) :: resid_pet !< residual pET remaining after subtracting simulated ET in the unsaturated zone
85  real(DP), intent(inout) :: deriv_et !< derivative of gw ET for Newton addition to equations in _fn()
86  real(DP), intent(inout) :: trhs !< total uzf rhs contribution to GWF model
87  real(DP), intent(inout) :: thcof !< total uzf hcof contribution to GWF model
88  real(DP), intent(in) :: hgwf !< calculated groundwater head
89  ! -- local
90  real(DP) :: etgw
91  real(DP) :: range
92  real(DP) :: depth, scale
93  !
94  depth = hgwf - (efflndsrf - extdp)
95  if (depth < dzero) depth = dzero
96  etgw = resid_pet
97  range = dem3 * extdp
98  call scubic(depth, range, deriv_et, scale)
99  etgw = etgw * scale
100  trhs = etgw
101  deriv_et = -deriv_et * etgw
102  etfunc_nlin = etgw
Here is the call graph for this function:
Here is the caller graph for this function: