-
Notifications
You must be signed in to change notification settings - Fork 0
/
gw_no3.f
120 lines (102 loc) · 5.34 KB
/
gw_no3.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
subroutine gw_no3
!! ~ ~ ~ PURPOSE ~ ~ ~
!! this subroutine estimates groundwater contribution to
!! streamflow
!! ~ ~ ~ INCOMING VARIABLES ~ ~ ~
!! name |units |definition
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! alpha_bf(:) |1/days |alpha factor for groundwater recession curve
!! alpha_bfe(:)|none |Exp(-alpha_bf(:))
!! deepst(:) |mm H2O |depth of water in deep aquifer
!! ihru |none |HRU number
!! gw_delaye(:)|none |Exp(-1./(delay(:)) where delay(:) is the
!! |groundwater delay (time required for water
!! |leaving the bottom of the root zone to reach
!! |the shallow aquifer; units-days)
!! gw_revap(:) |none |revap coeff: this variable controls the amount
!! |of water moving from the shallow aquifer to
!! |the root zone as a result of soil moisture
!! |depletion
!! gw_spyld(:) |m**3/m**3 |specific yield for shallow aquifer
!! gwht(:) |m |groundwater height
!! gwqmn(:) |mm H2O |threshold depth of water in shallow aquifer
!! |required before groundwater flow will occur
!! pet_day |mm H2O |potential evapotranspiration on current day
!! |in HRU
!! rchrg(:) |mm H2O |amount of water entering shallow aquifer on
!! |previous day in HRU
!! rchrg_dp(:) |none |recharge to deep aquifer: the fraction of
!! |root zone percolation that reaches the deep
!! |aquifer
!! revapmn(:) |mm H2O |threshold depth of water in shallow aquifer
!! |required to allow revap to occur
!! sepbtm |mm H2O |percolation from bottom of soil profile for
!! |the day in HRU
!! shallst(:) |mm H2O |depth of water in shallow aquifer
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! ~ ~ ~ OUTGOING VARIABLES ~ ~ ~
!! name |units |definition
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! deepst(:) |mm H2O |depth of water in deep aquifer
!! gw_q(:) |mm H2O |groundwater contribution to streamflow from
!! |HRU on current day
!! gwht(:) |m |groundwater height
!! gwseep |mm H2O |amount of water recharging deep aquifer on
!! |current day in HRU
!! rchrg(:) |mm H2O |amount of water recharging both aquifers on
!! |current day in HRU
!! revapday |mm H2O |amount of water moving from the shallow
!! |aquifer into the soil profile or being taken
!! |up by plant roots in the shallow aquifer
!! shallst(:) |mm H2O |depth of water in shallow aquifer
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! ~ ~ ~ LOCAL DEFINITIONS ~ ~ ~
!! name |units |definition
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! j |none |HRU number
!! rchrg1 |mm H2O |amount of water entering shallow aquifer on
!! |previous day
!! ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
!! ~ ~ ~ SUBROUTINES/FUNCTIONS CALLED ~ ~ ~
!! Intrinsic: Max
!! ~ ~ ~ ~ ~ ~ END SPECIFICATIONS ~ ~ ~ ~ ~ ~
!! revap is subtracted and rchrg is delayed (johnson, 1977)
use parm
integer :: j
real :: rchrgn1, revapn, gwseepn
j = 0
j = ihru
rchrgn1 = 0.
rchrgn1 = rchrg_n(j)
if (rchrgn1 < 1.e-6) rchrgn1 = 0.0
!! compute nitrate aquifer loading from recharge for current day
rchrg_n(j) = 0.
rchrg_n(j) = (1.- gw_delaye(j)) * percn(j) + gw_delaye(j)
& * rchrgn1
shallst_n(j) = shallst_n(j) + rchrg_n(j)
if (shallst_n(j) < 1.e-6) shallst_n(j) = 0.0
if (shallst(j) < 1.e-6) shallst(j) = 0.0
if (gw_q(j) < 1.e-6) gw_q(j) = 0.0
if (revapday < 1.e-6) revapday = 0.0
if (gwseep < 1.e-6) gwseep = 0.0
!! compute nitrate groundwater contribution to streamflow for day
xx = shallst(j) + gw_q(j) + revapday + gwseep
if (xx > 1.) then
xx = shallst_n(j) / (shallst(j) + gw_q(j) + revapday + gwseep)
else
xx = 0.
end if
if (xx < 1.e-6) xx = 0.0
no3gw(j) = xx * gw_q(j)
revapn = xx * revapday
gwseepn = xx * gwseep
revapn = amax1(1.e-6,revapn)
gwseepn = amax1(1.e-6,gwseepn)
!! subtract nitrate losses from the shallow aquifer
shallst_n(j) = shallst_n(j) - no3gw(j) - revapn - gwseepn
shallst_n(j) = amax1 (0., shallst_n(j))
!! compute nitrate losses in the groundwater
shallst_n(j) = shallst_n(j) * gw_nloss(j)
shallst_n(j) = amax1(0., shallst_n(j))
return
end