MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
UzfEtUtil.f90
Go to the documentation of this file.
2  use kindmodule, only: dp
3  use constantsmodule, only: dzero, done, dem4
4  use smoothingmodule, only: scubic
5 
6  implicit none
7  private
9 
10 contains
11 
12  !> @brief Calculate gwf et using linear ET function from mf-2005
13  !<
14  function etfunc_lin(efflndsrf, extdp, resid_pet, deriv_et, trhs, thcof, &
15  hgwf, celtop, celbot)
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
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  !
71  end function etfunc_lin
72 
73  !> @brief Calculate the linear scaling factor
74  !<
75  pure function calc_lin_scaling_fac(hgwf, lndsrf, extdp) result(sclfac)
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
86  end function calc_lin_scaling_fac
87 
88 end module uzfetutilmodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dem4
real constant 1e-4
Definition: Constants.f90:107
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module defines variable data types.
Definition: kind.f90:8
subroutine scubic(x, range, dydx, y)
@ brief sCubic
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.
Definition: UzfEtUtil.f90:16
pure real(dp) function, public calc_lin_scaling_fac(hgwf, lndsrf, extdp)
Calculate the linear scaling factor.
Definition: UzfEtUtil.f90:76