noah 1.14.0
Loading...
Searching...
No Matches
sflx_snippet.f90
Go to the documentation of this file.
1
5 module sflx_snippet
6
7 private
8
9 public frh2o
10
11 contains
12
28 subroutine frh2o &
29! --- inputs:
30 & ( tkelv, smc, sh2o, smcmax, bexp, psis, &
31! --- outputs:
32 & liqwat &
33 & )
34
35! ===================================================================== !
36! description: !
37! !
38! subroutine frh2o calculates amount of supercooled liquid soil water !
39! content if temperature is below 273.15k (t0). requires newton-type !
40! iteration to solve the nonlinear implicit equation given in eqn 17 !
41! of koren et al (1999, jgr, vol 104(d16), 19569-19585). !
42! !
43! new version (june 2001): much faster and more accurate newton !
44! iteration achieved by first taking log of eqn cited above -- less !
45! than 4 (typically 1 or 2) iterations achieves convergence. also, !
46! explicit 1-step solution option for special case of parameter ck=0, !
47! which reduces the original implicit equation to a simpler explicit !
48! form, known as the "flerchinger eqn". improved handling of solution !
49! in the limit of freezing point temperature t0. !
50! !
51! subprogram called: none !
52! !
53! !
54! ==================== defination of variables ==================== !
55! !
56! inputs: size !
57! tkelv - real, temperature (k) 1 !
58! smc - real, total soil moisture content (volumetric) 1 !
59! sh2o - real, liquid soil moisture content (volumetric) 1 !
60! smcmax - real, saturation soil moisture content 1 !
61! bexp - real, soil type "b" parameter 1 !
62! psis - real, saturated soil matric potential 1 !
63! !
64! outputs: !
65! liqwat - real, supercooled liquid water content 1 !
66! !
67! ==================== end of description ===================== !
68!
69! --- constant parameters:
70
71! this block added from physconst.f for snippet
72
73 implicit none
74
75 real, parameter :: gs2 = 9.81
76 real, parameter :: tfreez = 2.7315e+2
77 real, parameter :: lsubf = 3.335e5
78! end block added for snippet
79
80 real, parameter :: ck = 8.0
81 real, parameter :: blim = 5.5
82 real, parameter :: error = 0.005
83
84! --- inputs:
85 real, intent(in) :: tkelv, smc, sh2o, smcmax, bexp, psis
86
87! --- outputs:
88 real, intent(out) :: liqwat
89
90! --- locals:
91 real :: bx, denom, df, dswl, fk, swl, swlk
92
93 integer :: nlog, kcount
94!
95!===> ... begin here
96!
97! --- ... limits on parameter b: b < 5.5 (use parameter blim)
98! simulations showed if b > 5.5 unfrozen water content is
99! non-realistically high at very low temperatures.
100
101 bx = bexp
102 if (bexp > blim) bx = blim
103
104! --- ... initializing iterations counter and iterative solution flag.
105
106 nlog = 0
107 kcount= 0
108
109! --- ... if temperature not significantly below freezing (t0), sh2o = smc
110
111 if (tkelv > (tfreez-1.e-3)) then
112
113 liqwat = smc
114
115 else
116
117 if (ck /= 0.0) then
118
119! --- ... option 1: iterated solution for nonzero ck
120! in koren et al, jgr, 1999, eqn 17
121
122! --- ... initial guess for swl (frozen content)
123
124 swl = smc - sh2o
125
126! --- ... keep within bounds.
127
128 swl = max( min( swl, smc-0.02 ), 0.0 )
129
130! --- ... start of iterations
131
132 do while ( (nlog < 10) .and. (kcount == 0) )
133 nlog = nlog + 1
134
135 df = alog( (psis*gs2/lsubf) * ( (1.0 + ck*swl)**2.0 ) &
136 * (smcmax/(smc-swl))**bx ) - alog(-(tkelv-tfreez)/tkelv)
137
138 denom = 2.0*ck/(1.0 + ck*swl) + bx/(smc - swl)
139 swlk = swl - df/denom
140
141! --- ... bounds useful for mathematical solution.
142
143 swlk = max( min( swlk, smc-0.02 ), 0.0 )
144
145! --- ... mathematical solution bounds applied.
146
147 dswl = abs(swlk - swl)
148 swl = swlk
149
150! --- ... if more than 10 iterations, use explicit method (ck=0 approx.)
151! when dswl less or eq. error, no more iterations required.
152
153 if ( dswl <= error ) then
154 kcount = kcount + 1
155 endif
156 enddo ! end do_while_loop
157
158! --- ... bounds applied within do-block are valid for physical solution.
159
160 liqwat = smc - swl
161
162 endif ! end if_ck_block
163
164! --- ... option 2: explicit solution for flerchinger eq. i.e. ck=0
165! in koren et al., jgr, 1999, eqn 17
166! apply physical bounds to flerchinger solution
167
168 if (kcount == 0) then
169 fk = ( ( (lsubf/(gs2*(-psis))) &
170 * ((tkelv-tfreez)/tkelv) )**(-1/bx) ) * smcmax
171
172 fk = max( fk, 0.02 )
173
174 liqwat = min( fk, smc )
175 endif
176
177 endif ! end if_tkelv_block
178!
179 return
180!...................................
181 end subroutine frh2o
182!-----------------------------------
183 end module sflx_snippet