/[MITgcm]/MITgcm/pkg/fizhi/fizhi_lsm.F
ViewVC logotype

Contents of /MITgcm/pkg/fizhi/fizhi_lsm.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.9 - (show annotations) (download)
Tue Aug 24 14:36:55 2010 UTC (13 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint62k, checkpoint62j, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, HEAD
Changes since 1.8: +64 -64 lines
remove tabs

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_lsm.F,v 1.8 2010/03/16 00:19:33 jmc Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 SUBROUTINE TILE (
6 I NCH, DTSTEP, ITYP, TRAINL,TRAINC, TSNOW, UM,
7 I ETURB, DEDQA, DEDTC, HSTURB, DHSDQA, DHSDTC,
8 I TM, QM, CD, SUNANG, PARDIR, PARDIF,
9 I SWNET, HLWDWN, PSUR, ZLAI, GREEN, Z2,
10 I SQSCAT, RSOIL1, RSOIL2, RDC, U2FAC,
11 I QSATTC, DQSDTC, ALWRAD, BLWRAD,
12 U TC, TD, QA, SWET1, SWET2, SWET3, CAPAC, SNOW,
13 O EVAP, SHFLUX, RUNOFF, BOMB,
14 O EINT, ESOI, EVEG, ESNO, SMELT, HLATN,
15 O HLWUP,GDRAIN,RUNSRF,FWSOIL,
16 O STRDG1, STRDG2, STRDG3, STRDG4,
17 O STRDG5, STRDG6, STRDG7, STRDG8, STRDG9
18 & )
19
20 C SCCS VERSION @(#)lsm.f 1.2 11/3/92
21
22 C****
23 C**** This subroutine computes the outgoing fluxes of sensible and
24 C**** latent heat from the land surface and updates the surface prognostic
25 C**** variables.
26 C****
27 IMPLICIT NONE
28 C
29 C****
30 #include "sibber.h"
31 C****
32 INTEGER NCH
33 INTEGER ITYP(NCH)
34 _RL DTSTEP, TRAINL(NCH), TRAINC(NCH), TSNOW(NCH), UM(NCH),
35 & ETURB(NCH), DEDQA(NCH), HSTURB(NCH), DHSDTC(NCH),
36 & TM (NCH), CD(NCH), SUNANG(NCH), DHSDQA(NCH),
37 & QM (NCH), PARDIR(NCH), PARDIF(NCH), SWNET(NCH),
38 & HLWDWN(NCH), PSUR(NCH), ZLAI(NCH), GREEN(NCH),
39 & Z2(NCH), SQSCAT(NCH), DEDTC(NCH)
40 _RL RSOIL1(NCH), RSOIL2(NCH), RDC(NCH), U2FAC(NCH),
41 & QSATTC(NCH), DQSDTC(NCH), ALWRAD(NCH), BLWRAD(NCH),
42 & TC(NCH), TD(NCH), QA(NCH), BOMB(NCH),
43 & SWET1(NCH), SWET2(NCH), SWET3(NCH), CAPAC(NCH),
44 & SNOW(NCH), EVAP(NCH), SHFLUX(NCH), RUNOFF(NCH)
45 _RL EINT(NCH), ESOI(NCH), EVEG(NCH), ESNO(NCH),
46 & STRDG1(NCH), STRDG2(NCH), STRDG3(NCH), STRDG4(NCH),
47 & STRDG5(NCH), STRDG6(NCH), STRDG7(NCH), STRDG8(NCH),
48 & STRDG9(NCH),
49 & SMELT(NCH), HLATN(NCH), HLWUP(NCH), GDRAIN(NCH),
50 & RUNSRF(NCH), FWSOIL(NCH)
51 C****
52 INTEGER ChNo
53 _RL SWET(nch,NLAY), VGPSAT(NTYPS), VGCSAT (NTYPS),
54 & VGZDEP(NLAY,NTYPS), VGSLOP(NTYPS), DELTC,
55 & DELEA, VGPH1(NTYPS), VGPH2(NTYPS),
56 & VGRPLN(NTYPS), CSOIL0(NTYPS), WSOI12,
57 & VGBEE(NTYPS), DELZ12(NTYPS),
58 & DELZ23(NTYPS)
59 C****
60
61 _RL PHLAY(nch,NLAY), AKLAY(nch,NLAY), SWET12(nch),
62 & CSOIL(nch),
63 & RCUN(nch), VPDSTR(nch), ESATTX(nch),
64 & VPDSTX(nch), VGBEEX(nch)
65 _RL EMAXRT(nch), VGWMAX(NLAY,NTYPS), FTEMP(nch),
66 & PHR(nch), SOILCO(nch), RC(nch),
67 & EAX(nch), TX(nch), RCX(nch),
68 & DRCDTC(nch), SATCAP(nch), PAR(nch),
69 & PDIR(nch), DUMMY(nch)
70 _RL FTEMPX(nch), DRCDEA(nch), VGPSAX(nch), VGCSAX(nch),
71 & VGZDEX(NLAY,nch), VGSLOX(nch), VGPH1X(nch),
72 & VGPH2X(nch), VGRPLX(nch)
73 _RL DEDEA(nch), DHSDEA(nch), EM(nch), ESATTC(nch),
74 & DESDTC(nch), EA(nch), RA(nch), ALHX(nch),
75 & WETEQ1(nch), WETEQ2(nch),
76 & RX1(nch), RX2(nch), SNWFRC(nch), POTFRC(nch),
77 & ESNFRC(nch), EIRFRC(nch), FCAN(nch), EPFRC,
78 & DEFCIT, EADJST, RTBS
79 _RL cmpbug
80
81 C****
82 DATA VGWMAX /8.4, 621.6, 840.0,
83 2 8.4, 621.6, 840.0,
84 3 8.4, 621.6, 840.0,
85 4 8.4, 197.4, 420.0,
86 5 8.704, 204.5, 435.2,
87 6 8.4, 71.4, 420.0,
88 7 4.000, 4.000, 130.56,
89 8 4.000, 4.000, 130.56,
90 9 8.704, 73.984, 130.56,
91 1 4.000, 4.000, 130.56
92 & /
93 DATA VGBEE / 4., 4., 4., 4., 1.69, 4.,
94 & 4., 1.69, 1.69, 1.69/
95 C DATA VGBEE / 7., 7., 7., 7., 1.69, 7.,
96 C & 7., 1.69, 1.69, 1.69/
97 c DATA VGBEE / 4., 4., 4., 4., 0.95, 4.,
98 c & 4., 0.95, 0.95, 0.95/
99 C DATA VGBEE /7, 7, 7, 7, 2, 7, 7, 4, 2, 4/
100 DATA VGPSAT/ -0.281, -0.281, -0.281, -0.281,
101 5 -0.073, -0.281, -0.281, -0.073,
102 9 -0.073 , -0.073/
103 C DATA VGPSAT/ -0.086, -0.086, -0.086, -0.086,
104 C 5 -0.073, -0.086, -0.086, -0.073,
105 C 9 -0.073 , -0.073/
106 c DATA VGPSAT/ -0.281, -0.281, -0.281, -0.281,
107 c 5 -0.014, -0.281, -0.281, -0.014,
108 c 9 -0.014 , -0.014/
109 c DATA VGCSAT / 0.00002, 0.00002, 0.00002, 0.00002,
110 c 5 0.0000583, 0.00002, 0.00002, 0.0000583,
111 c 9 0.0000583, 0.0000583
112 c & /
113 DATA VGCSAT / 0.0000012, 0.0000012, 0.0000012, 0.0000012,
114 5 0.0000583, 0.0000012, 0.0000012, 0.0000583,
115 9 0.0000583, 0.0000583
116 & /
117 C**** - - - - - - - -
118 C**** THE FOLLOWING 3 VARIABLES MUST MAINTAIN CONSISTENT VALUES. ZDEP12(N) =
119 C**** (VGZDEP(1,N)+VGZDEP(2,N))/2. ; SIMILARLY FOR ZDEP23(N).
120 C****
121 DATA VGZDEP /.02, 1.48, 2.00,
122 2 .02, 1.48, 2.00,
123 3 .02, 1.48, 2.00,
124 4 .02, .47, 1.00,
125 5 .02, .47, 1.00,
126 6 .02, .17, 1.00,
127 7 .0092, .0092, .30,
128 8 .0092, .0092, .30,
129 9 .02, .17, .30,
130 1 .0092, .0092, .30
131 & /
132 DATA DELZ12 / 0.75, 0.75, 0.75, 0.245, 0.245, 0.095, 0.0092,
133 8 0.0092, 0.095, 0.0092 /
134 DATA DELZ23 / 1.74, 1.74, 1.74, 0.735, 0.735, 0.585, 0.155,
135 8 0.155, 0.235, 0.155 /
136 C**** - - - - - - - -
137 DATA VGSLOP / .1736, .1736, .1736, .1736,
138 5 1.000, .1736, .1736, 1.000,
139 & 1.000 , 1.000 /
140 c DATA VGSLOP / .1736, .1736, .1736, .1736,
141 c 5 .0872, .1736, .1736, .0872,
142 c & .0872 , .0872 /
143 c DATA VGSLOP / 1., 1., 1., 1., 1., 1., 1., 1., 1.,1./
144 DATA VGPH1 / -100., -190., -200., -120.,
145 5 -200., -200., -200., -10.,
146 9 -10., -10. /
147 DATA VGPH2 / -500., -250., -250., -230.,
148 5 -400., -400., -250., -100.,
149 9 -100., -100. /
150 DATA VGRPLN /0.245E+09, 0.245E+09, 0.245E+09, 0.250E+09,
151 5 0.250E+09, 0.250E+09, 0.245E+09, 0.1E+09,
152 9 0.1E+09, 0.100E+09 /
153 C****
154 DATA DELTC /0.01/, DELEA /0.001/
155
156 c Note: SNWMID & CSOIL0 parameters modified from (Koster/Suarez)
157 c to obtain improved albedo, radswg, and annual/diurnal cycle (L.Takacs 2/17/00)
158 c ------------------------------------------------------------------------------------
159 DATA CSOIL0 /175000.,175000.,175000.,175000.,175000.,
160 . 175000.,175000.,120000.,175000., 70000./
161 #include "snwmid.h"
162
163 C****
164 C**** ---------------------------------------------------------------------
165 C****
166 CFPP$ EXPAND (TMPFAC, RCANOP, SOIL, WUPDAT, SMFAC)
167 C****
168 C**** Expand data as specified by ITYP into arrays of size NCH.
169 C****
170 DO 50 ChNo = 1, NCH
171 BOMB(CHNO) = 0.
172 VGBEEX(ChNo) = VGBEE(ITYP(ChNo))
173 VGPSAX(ChNo) = VGPSAT(ITYP(ChNo))
174 VGCSAX(ChNo) = VGCSAT(ITYP(ChNo))
175 VGZDEX(SFCLY ,ChNo) = VGZDEP(SFCLY ,ITYP(ChNo))
176 VGZDEX(ROOTLY,ChNo) = VGZDEP(ROOTLY,ITYP(ChNo))
177 VGZDEX(RECHLY,ChNo) = VGZDEP(RECHLY,ITYP(ChNo))
178 VGSLOX(ChNo) = VGSLOP(ITYP(ChNo))
179 VGPH1X(ChNo) = VGPH1(ITYP(ChNo))
180 VGPH2X(ChNo) = VGPH2(ITYP(ChNo))
181 VGRPLX(ChNo) = VGRPLN(ITYP(ChNo))
182 50 CONTINUE
183
184 C****
185 C**** Pre-process input arrays as necessary:
186
187 DO 100 ChNo = 1, NCH
188
189 DEDQA(CHNO) = MAX( DEDQA(CHNO), 500. _d 0/ALHE )
190 DEDTC(CHNO) = MAX( DEDTC(CHNO), 0. _d 0)
191 DHSDQA(CHNO) = MAX( DHSDQA(CHNO), 0. _d 0)
192 DHSDTC(CHNO) = MAX( DHSDTC(CHNO), -10. _d 0)
193
194 EM(CHNO) = QM(CHNO) * PSUR(CHNO) / EPSILON
195 EA(CHNO) = QA(CHNO) * PSUR(CHNO) / EPSILON
196 ESATTC(CHNO) = QSATTC(CHNO) * PSUR(CHNO) / EPSILON
197 DESDTC(CHNO) = DQSDTC(CHNO) * PSUR(CHNO) / EPSILON
198 C DEDEA(CHNO) = DEDQA(CHNO) * EPSILON / ( PSUR(CHNO) * ALHX(CHNO) )
199 DEDEA(CHNO) = DEDQA(CHNO) * EPSILON / PSUR(CHNO)
200 DHSDEA(CHNO) = DHSDQA(CHNO) * EPSILON / PSUR(CHNO)
201 C DEDTC(CHNO) = DEDTC(CHNO) / ALHX(CHNO)
202 C ETURB(CHNO) = ETURB(CHNO) / ALHX(CHNO)
203 PAR(CHNO) = (PARDIR(CHNO) + PARDIF(CHNO) + 1.E-20)
204 PDIR(CHNO) = PARDIR(CHNO) / PAR(CHNO)
205 RA(CHNO) = ONE / ( CD(CHNO) * UM(CHNO) )
206 SATCAP(ChNo) = 0.1 * ZLAI(ChNo)
207 CSOIL(CHNO) = CSOIL0(ITYP(ChNo))
208 SWET(ChNo,SFCLY ) = max( min(SWET1(ChNo),1. _d 0), 0. _d 0)
209 SWET(ChNo,ROOTLY) = max( min(SWET2(ChNo),1. _d 0), 0. _d 0)
210 SWET(ChNo,RECHLY) = max( min(SWET3(ChNo),1. _d 0), 0. _d 0)
211 CAPAC(CHNO) = max( min(CAPAC(ChNo),SATCAP(CHNO)), 0. _d 0)
212 C****
213
214 SNWFRC(CHNO) = SNOW(CHNO) / ( SNOW(CHNO) + SNWMID(ITYP(CHNO)) )
215 FCAN(CHNO) = MIN( 1. _d 0, MAX(0. _d 0,CAPAC(ChNo)/SATCAP(ChNo)) )
216 POTFRC(CHNO)=1.-(1.-SNWFRC(CHNO))*(1.-FCAN(CHNO))
217
218
219 WSOI12=SWET(ChNo,SFCLY ) * VGWMAX(SFCLY ,ITYP(ChNo)) +
220 & SWET(ChNo,ROOTLY) * VGWMAX(ROOTLY,ITYP(ChNo))
221 SWET12(ChNo) = WSOI12 /
222 & (VGWMAX(SFCLY,ITYP(ChNo)) + VGWMAX(ROOTLY,ITYP(ChNo)))
223
224 EMAXRT(CHNO) = ( SNOW(CHNO) + CAPAC(CHNO) + WSOI12 ) / DTSTEP
225
226 C****
227 RUNOFF(CHNO) = 0.
228 RUNSRF(CHNO) = 0.
229 C**** (SMELT is zeroed in FLUXES)
230 C SMELT(CHNO) = 0.
231 FWSOIL(CHNO) = 0.
232 100 CONTINUE
233
234 C****
235 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236 C**** STEP 1: COMPUTE EFFECTIVE RESISTANCE RC FOR ENERGY BALANCE.
237 C**** (RCUNST computes the unstressed resistance,
238 C**** VPDFAC computes the vapor pressure deficit stress term,
239 C**** TMPFAC computes the temperature stress term,
240 C**** SOIL computes soil moisture potentials and conds.,
241 C**** SMFAC computes rc given a leaf water potential stress,
242 C**** RSURFP computes rc given a parallel resist. from the
243 C**** surface,
244 C**** RCANOP computes rc corrected for snow and interception.)
245 C****
246
247 CALL RCUNST (
248 I NCH, ITYP, SUNANG, SQSCAT, PDIR,
249 I PAR, ZLAI, GREEN,
250 O RCUN
251 & )
252
253 CALL VPDFAC (
254 I NCH, ITYP, ESATTC, EA,
255 O VPDSTR
256 & )
257
258 CALL TMPFAC (
259 I NCH, ITYP, TC,
260 O FTEMP
261 & )
262
263 CALL SOIL (
264 I NCH, ITYP, SWET12, VGBEEX, VGPSAX, VGCSAX, DELZ12,
265 O PHR, SOILCO, DUMMY
266 & )
267
268 cmpbug=0.
269
270 CALL SMFAC (
271 I NCH, ITYP, ESATTC, EA, PHR, SOILCO, RCUN,
272 I VPDSTR, FTEMP, TC, PSUR, Z2, RSOIL1, RSOIL2,
273 I VGPH1X, VGPH2X, VGRPLX,
274 O RC
275 & )
276
277 CALL RSURFP (
278 I NCH, UM, U2FAC, Z2, RDC, SWET(1,SFCLY),
279 I ESATTC, EA,
280 U RC,
281 O RX1, RX2
282 & )
283
284 CALL RCANOP (
285 I NCH, CAPAC, SNOW, SATCAP, RA, ETURB, SNWFRC,
286 I POTFRC,
287 U RC
288 & )
289
290 C****
291 C**** - - - - - - - - - - - - - -
292 C****
293 C**** Compute DRC/DT and DRC/DEA using temperature, v.p. perturbations:
294 C****
295
296 DO 120 ChNo = 1, NCH
297 TX(ChNo) = TC(ChNo) + DELTC
298 ESATTX(ChNo) = ESATTC(ChNo) + DESDTC(CHNO) * DELTC
299 EAX(ChNo) = EA(ChNo) + DELEA
300 120 CONTINUE
301 C****
302 C**** temperature:
303 CALL VPDFAC (NCH, ITYP, ESATTX, EA, VPDSTX)
304 CALL TMPFAC (NCH, ITYP, TX, FTEMPX)
305 CALL SMFAC (
306 I NCH, ITYP, ESATTX, EA, PHR, SOILCO, RCUN,
307 I VPDSTX, FTEMPX, TX, PSUR, Z2, RSOIL1, RSOIL2,
308 I VGPH1X, VGPH2X, VGRPLX,
309 O RCX
310 & )
311 CALL RSURFP (NCH, UM, U2FAC, Z2, RDC, SWET(1,SFCLY),
312 I ESATTX, EA,
313 & RCX,
314 O DUMMY,DUMMY
315 & )
316 CALL RCANOP (NCH, CAPAC, SNOW, SATCAP, RA, ETURB, SNWFRC,
317 & POTFRC, RCX)
318 C****
319 DO 125 ChNo = 1, NCH
320 DRCDTC(ChNo) = (RCX(ChNo) - RC(ChNo)) / DELTC
321 125 CONTINUE
322 C****
323 C**** vapor pressure:
324 CALL VPDFAC (NCH, ITYP, ESATTC, EAX, VPDSTX)
325 CALL SMFAC (
326 I NCH, ITYP, ESATTC, EAX, PHR, SOILCO, RCUN,
327 I VPDSTX, FTEMP, TC, PSUR, Z2, RSOIL1, RSOIL2,
328 I VGPH1X, VGPH2X, VGRPLX,
329 O RCX
330 & )
331 CALL RSURFP (NCH, UM, U2FAC, Z2, RDC, SWET(1,SFCLY),
332 I ESATTC, EAX,
333 & RCX,
334 O DUMMY,DUMMY
335 & )
336 CALL RCANOP (NCH, CAPAC, SNOW, SATCAP, RA, ETURB, SNWFRC,
337 & POTFRC, RCX)
338 C****
339 DO 150 ChNo = 1, NCH
340 DRCDEA(ChNo) = (RCX(ChNo) - RC(ChNo)) / DELEA
341 150 CONTINUE
342 C****
343
344 C****
345 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
346 C**** STEP 2: Solve the energy balance at the surface.
347 C****
348 C****
349 C**** Determine effective latent heat of vaporization based on
350 C**** assumed fractional coverage of snow. This requires the
351 C**** determination of the fractions of moisture taken from the various
352 C**** reservoirs:
353 C**** ESNFRC=fraction of total evap. coming from snow
354 C**** EIRFRC=fraction of total evap. coming from interception res.
355
356 DO 200 CHNO=1,NCH
357 RTBS=RX1(CHNO)*RX2(CHNO)/(RX1(CHNO)+RX2(CHNO))
358 EPFRC=POTFRC(CHNO) * ( RA(CHNO) + RTBS ) /
359 & ( RA(CHNO) + POTFRC(CHNO)*RTBS )
360 ESNFRC(CHNO)=EPFRC*SNWFRC(CHNO)/
361 & (SNWFRC(CHNO)+FCAN(CHNO)+1.E-20)
362 EIRFRC(CHNO)=EPFRC*FCAN(CHNO)/(SNWFRC(CHNO)+FCAN(CHNO)+1.E-20)
363 ALHX(CHNO) = (1.-ESNFRC(CHNO))*ALHE + ESNFRC(CHNO)*ALHS
364 200 CONTINUE
365
366 CALL FLUXES (
367 I NCH, ITYP, DTSTEP, ESATTC, DESDTC, ALHX,
368 I ETURB, DEDEA, DEDTC, HSTURB, DHSDEA, DHSDTC,
369 I RC, DRCDEA, DRCDTC,
370 I SWNET, HLWDWN, ALWRAD, BLWRAD, ESNFRC,
371 I TM, EM, CSOIL, PSUR, EMAXRT, VGWMAX,
372 U TC, TD, EA, SWET, SNOW,
373 O RUNOFF, EVAP, SHFLUX, SMELT, HLWUP, BOMB,STRDG1,
374 O STRDG2, STRDG3, STRDG4, STRDG5, STRDG6, STRDG7,
375 O STRDG8, STRDG9
376 & )
377 c****
378 C****
379 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
380 C**** STEP 3: Update the water quantities in the soil and in the
381 C**** interception reservoir.
382 C**** (WUPDAT reduces moisture contents using calculated
383 C**** evap.,
384 C**** GWATER computes darcian flux of water between soil
385 C**** layers.)
386 C****
387
388 CALL WUPDAT (
389 I NCH, ITYP, DTSTEP,
390 I EVAP, SATCAP, VGWMAX, TC, RA, RC,
391 I RX1, RX2,ESNFRC,EIRFRC,
392 U CAPAC, SNOW, SWET, RUNOFF,
393 O EINT, ESOI, EVEG, ESNO, RUNSRF,FWSOIL
394 & )
395
396 C****
397 C**** Correct any energy balance inconsistency.
398 C****
399 DO 300 CHNO=1,NCH
400 IF(EVAP(CHNO) .GT. 0.) THEN
401 EADJST=(ESNO(CHNO)/DTSTEP) - ESNFRC(CHNO)*EVAP(CHNO)
402 SHFLUX(CHNO)=SHFLUX(CHNO)-EADJST*(ALHS-ALHE)
403 ENDIF
404 300 CONTINUE
405
406
407
408 CALL SOIL (
409 I NCH, ITYP, SWET (1, SFCLY), VGBEEX,
410 I VGPSAX, VGCSAX, DELZ12,
411 O PHLAY(1,SFCLY), AKLAY(1,SFCLY), DUMMY
412 & )
413 CALL SOIL (
414 I NCH, ITYP, SWET(1,ROOTLY), VGBEEX,
415 I VGPSAX, VGCSAX, DELZ12,
416 O PHLAY(1,ROOTLY), AKLAY(1,ROOTLY), WETEQ1
417 & )
418 CALL SOIL (
419 I NCH, ITYP, SWET(1,RECHLY), VGBEEX,
420 I VGPSAX, VGCSAX, DELZ23,
421 O PHLAY(1,RECHLY), AKLAY(1,RECHLY), WETEQ2
422 & )
423
424 CALL GWATER (
425 I NCH, ITYP, VGWMAX, PHLAY, AKLAY, TC, DTSTEP,
426 I VGZDEX, VGSLOX, WETEQ1, WETEQ2,
427 I VGPSAX, VGCSAX,
428 U SWET, RUNOFF, GDRAIN
429 & )
430
431 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
432 C****
433 C**** STEP 4: Allow precipitation to fill interception reservoir
434 C**** and top soil layer
435 C****
436
437 CALL SOIL (
438 I NCH, ITYP, SWET(1,ROOTLY), VGBEEX,
439 I VGPSAX, VGCSAX, DELZ12,
440 O PHLAY(1,ROOTLY), AKLAY(1,ROOTLY), WETEQ1
441 & )
442
443 CALL INTERC (
444 I NCH, ITYP, DTSTEP, TRAINL, TRAINC, TSNOW, SATCAP,
445 I VGWMAX, CSOIL, WETEQ1,
446 U TC, CAPAC, SNOW, SWET, RUNOFF, RUNSRF,
447 & SMELT, FWSOIL
448 & )
449
450
451 C****
452 C****
453 C**** Process data for return to GCM:
454
455 DO 2000 ChNo = 1, NCH
456 QA(CHNO) = EA(CHNO) * EPSILON / PSUR(CHNO)
457 C DEDTC(CHNO) = DEDTC(CHNO) * ALHX(CHNO)
458 C ETURB(CHNO) = ETURB(CHNO) * ALHX(CHNO)
459 HLATN (CHNO) = EVAP (CHNO) * ALHX(CHNO)
460 SWET1(ChNo) = SWET(ChNo,SFCLY)
461 SWET2(ChNo) = SWET(ChNo,ROOTLY)
462 SWET3(ChNo) = SWET(ChNo,RECHLY)
463
464 IF(SWET1(ChNo).LT.1.E-10) SWET1(CHNO) = 0.0
465 IF(SWET2(ChNo).LT.1.E-10) SWET2(CHNO) = 0.0
466 IF(SWET3(ChNo).LT.1.E-10) SWET3(CHNO) = 0.0
467 IF(CAPAC(ChNo).LT.1.E-10) CAPAC(CHNO) = 0.0
468 IF(SNOW (ChNo).LT.1.E-10) SNOW (CHNO) = 0.0
469 IF(RUNOFF(ChNo).LT.1.E-10) RUNOFF(CHNO) = 0.0
470
471 EINT(CHNO) = EINT(CHNO) * ALHE / DTSTEP
472 ESOI(CHNO) = ESOI(CHNO) * ALHE / DTSTEP
473 EVEG(CHNO) = EVEG(CHNO) * ALHE / DTSTEP
474 ESNO(CHNO) = ESNO(CHNO) * ALHS / DTSTEP
475
476 DEFCIT = EVAP(CHNO) * ( RC(CHNO) + RA(CHNO) )
477 CCOOMMSTRDG1(CHNO) = DEFCIT / RA(CHNO)
478 CCOOMMSTRDG2(CHNO) = DEFCIT / ( RA(CHNO) + RCUN(CHNO) )
479 CCOOMMSTRDG3(CHNO) = DEFCIT / ( RA(CHNO) + RCUN(CHNO)/VPDSTR(CHNO) )
480 CCOOMMSTRDG4(CHNO) = DEFCIT / ( RA(CHNO) +
481 CCOOMM RCUN(CHNO)/(VPDSTR(CHNO)*FTEMP(CHNO)) )
482
483 2000 CONTINUE
484 C****
485
486 RETURN
487 END
488 C****
489 C**** [ END CHIP ]
490 C****
491 C**** -----------------------------------------------------------------
492 C**** /////////////////////////////////////////////////////////////////
493 C**** -----------------------------------------------------------------
494 C****
495 C**** [ BEGIN INTERC ]
496 C****
497 SUBROUTINE INTERC (
498 I NCH, ITYP, DTSTEP, TRAINL, TRAINC,
499 I TSNOW, SATCAP, WMAX, CSOIL, WETEQ1,
500 U TC, CAPAC, SNOW, SWET1, RUNOFF, RUNSRF,
501 U SMELT, FWSOIL
502 & )
503 C****
504 C**** This routine uses the precipitation forcing to determine
505 C**** changes in interception, snowcover, and soil moisture storage.
506 C****
507 C**** Note: In this formulation, rain that falls on frozen ground
508 C**** runs-off, rather than freezes.
509 C****
510 IMPLICIT NONE
511 C****
512 #include "sibber.h"
513 C****
514 INTEGER NCH
515 INTEGER ITYP(NCH), ChNo
516 _RL TRAINL(NCH), TRAINC(NCH), TSNOW(NCH), SATCAP(NCH),
517 & WMAX(NLAY,NTYPS), TC(NCH), CSOIL(NCH), CAPAC(NCH),
518 & SNOW(NCH), SWET1(NCH), RUNOFF(NCH), SMELT(NCH),
519 & RUNSRF(NCH), FWSOIL(NCH), WETEQ1(NCH), WETINT
520 _RL DTSTEP, SNOWM, WATADD, CAVAIL, THRUC, WRUNC, WRUNL,
521 & TIMFRL, TIMFRC, FWETL, FWETC, THRU1, THRU2, THRUL, XTCORR,
522 & WETFRC
523 C****
524 C DATA FWET0 /0.30/
525 DATA FWETL /1.00/, FWETC /0.20/, TIMFRL/1.00/, TIMFRC/0.125/
526 C****
527 C**** ------------------------------------------------------------------
528 C**** Loop over chips:
529 DO 100 ChNo = 1, NCH
530 C****
531 C**** Add to snow cover. Melt snow if necessary.
532 SNOW(ChNo) = SNOW(ChNo) + TSNOW(ChNo)*DTSTEP
533 SNOWM = 0.
534 IF( SNOW(CHNO).GT.0. .AND. TC(CHNO).GT.TF ) THEN
535 SNOWM = MIN( SNOW(ChNo),
536 & MAX( 0. _d 0, (TC(ChNo)-TF)*CSOIL(ChNo)/ALHM ) )
537 IF( SNOWM .EQ. SNOW(CHNO) ) THEN
538 TC(ChNo) = TC(ChNo) - SNOWM * ALHM / CSOIL(ChNo)
539 SNOW(CHNO)=0.
540 ELSE
541 TC(CHNO)=TF
542 SNOW(ChNo) = SNOW(ChNo) - SNOWM
543 ENDIF
544 SMELT(CHNO)=SMELT(CHNO)+SNOWM/DTSTEP
545 ENDIF
546 C****
547 C**** =======================================================
548 C****
549 C**** Load interception reservoir. STEP 1: Large scale condensation.
550 C****
551 C**** Determine XTCORR, the fraction of a storm that falls on a previously
552 C**** wet surface due to the time correlation of precipitation position.
553 C**** (Time scale TIMFRL for large scale storms set to one for FWETL=1
554 C**** to reflect the effective loss of "position memory" when storm
555 C**** covers entire grid square.)
556
557 XTCORR= (1.-TIMFRL) * MIN( 1. _d 0,(CAPAC(CHNO)/SATCAP(CHNO))/
558 $ FWETL )
559
560 C****
561 C**** Fill interception reservoir with precipitation.
562 C**** THRU1 is first calculated as the amount falling through the
563 C**** canopy under the assumption that all rain falls randomly.
564 C**** only a fraction 1-XTCORR falls randomly, though, so the result
565 C**** is multiplied by 1-XTCORR.
566 C****
567 WATADD = TRAINL(ChNo)*DTSTEP + SMELT(CHNO)*DTSTEP
568 CAVAIL = ( SATCAP(CHNO) - CAPAC(CHNO) ) * FWETL
569 WETINT = CAPAC(CHNO)/SATCAP(CHNO)
570 IF( WATADD*(1.-WETINT) .LT. CAVAIL ) THEN
571 THRU1 = WATADD*WETINT
572 ELSE
573 THRU1 = (WATADD - CAVAIL)
574 ENDIF
575 THRU1=THRU1*(1.-XTCORR)
576
577 C**** THRU2 is the amount that falls immediately through the canopy due
578 C**** to 'position memory'.
579
580 THRU2=XTCORR*WATADD
581
582 THRUL=THRU1+THRU2
583 CAPAC(CHNO)=CAPAC(CHNO)+WATADD-THRU1-THRU2
584 C****
585 C**** ---------------------------------------------------
586 C****
587 C**** STEP 2: Moist convective precipitation.
588 C****
589 C**** Determine XTCORR, the fraction of a storm that falls on a previously
590 C**** wet surface due to the time correlation of precipitation position.
591
592 XTCORR= (1.-TIMFRC) * MIN( 1. _d 0,(CAPAC(CHNO)/SATCAP(CHNO))/
593 $ FWETC )
594
595 C****
596 C**** Fill interception reservoir with precipitation.
597 C**** THRU1 is first calculated as the amount falling through the
598 C**** canopy under the assumption that all rain falls randomly.
599 C**** only a fraction 1-XTCORR falls randomly, though, so the result
600 C**** is multiplied by 1-XTCORR.
601 C****
602 WATADD = TRAINC(ChNo)*DTSTEP
603 CAVAIL = ( SATCAP(CHNO) - CAPAC(CHNO) ) * FWETC
604 WETINT = CAPAC(CHNO)/SATCAP(CHNO)
605 IF( WATADD*(1.-WETINT) .LT. CAVAIL ) THEN
606 THRU1 = WATADD*WETINT
607 ELSE
608 THRU1 = (WATADD - CAVAIL)
609 ENDIF
610 THRU1=THRU1*(1.-XTCORR)
611
612 C**** THRU2 is the amount that falls immediately through the canopy due
613 C**** to 'position memory'.
614
615 THRU2=XTCORR*WATADD
616
617 THRUC=THRU1+THRU2
618 CAPAC(CHNO)=CAPAC(CHNO)+WATADD-THRU1-THRU2
619 C****
620 C***** =================================================================
621 C****
622 C**** Add precipitation moisture to soil, generate runoff. if
623 C**** temperature is below freezing, all precipitation runs off.
624 C****
625 C**** STEP 1: Large scale condensation:
626 C****
627 IF(SWET1(CHNO).GT.WETEQ1(CHNO) .AND. WETEQ1(CHNO).NE.1.) THEN
628 WETFRC=(SWET1(CHNO)-WETEQ1(CHNO))/(1.-WETEQ1(CHNO))
629 ELSE
630 WETFRC=0.
631 ENDIF
632
633 CAVAIL = ( 1.-WETFRC)*FWETL*WMAX(1,ITYP(ChNo))*(1-WETEQ1(CHNO))
634 IF ( THRUL * (1-WETFRC) .LT. CAVAIL ) THEN
635 WRUNL = THRUL * WETFRC
636 SWET1(ChNo) = SWET1(ChNo)
637 * + THRUL * ( 1. - WETFRC) / WMAX(1,ITYP(ChNo))
638 ELSE
639 WRUNL = THRUL - CAVAIL
640 SWET1(CHNO) = SWET1(CHNO) + CAVAIL / WMAX(1,ITYP(ChNo))
641 ENDIF
642 C****
643
644 C**** STEP 2: Moist convective precipitation:
645 C****
646 IF(SWET1(CHNO).GT.WETEQ1(CHNO) .AND. WETEQ1(CHNO).NE.1.) THEN
647 WETFRC=(SWET1(CHNO)-WETEQ1(CHNO))/(1.-WETEQ1(CHNO))
648 ELSE
649 WETFRC=0.
650 ENDIF
651
652 CAVAIL = (1.-WETFRC)*FWETC*WMAX(1,ITYP(ChNo))*(1-WETEQ1(CHNO))
653 IF ( THRUC * (1-WETFRC) .LT. CAVAIL ) THEN
654 WRUNC = THRUC * WETFRC
655 SWET1(ChNo) = SWET1(ChNo)
656 * + THRUC * ( 1. - WETFRC) / WMAX(1,ITYP(ChNo))
657 ELSE
658 WRUNC = THRUC - CAVAIL
659 SWET1(CHNO) = SWET1(CHNO) + CAVAIL / WMAX(1,ITYP(ChNo))
660 ENDIF
661 C****
662 RUNOFF(ChNo) = RUNOFF(ChNo) + (WRUNC+WRUNL)/DTSTEP
663 RUNSRF(ChNo) = RUNSRF(ChNo) + (WRUNC+WRUNL)/DTSTEP
664 FWSOIL(CHNO) = FWSOIL(CHNO) + (THRUC+THRUL-WRUNC-WRUNL)/DTSTEP
665 C****
666 100 CONTINUE
667 C****
668 RETURN
669 END
670 C****
671 C**** [ END INTERC ]
672 C****
673 C**** -----------------------------------------------------------------
674 C**** /////////////////////////////////////////////////////////////////
675 C**** -----------------------------------------------------------------
676 C****
677 C**** [ BEGIN RCUNST ]
678 C****
679 SUBROUTINE RCUNST (
680 I NCH, ITYP, SUNANG, SQSCAT, PDIR,
681 I PAR, ZLAI, GREEN,
682 O RCUN
683 & )
684 C****
685 C**** This subroutine calculates the unstressed canopy resistance.
686 C**** (p. 1353, Sellers 1985.) Extinction coefficients are computed first.
687 C****
688 IMPLICIT NONE
689 C****
690 #include "sibber.h"
691 C****
692 INTEGER NCH
693 INTEGER ITYP(NCH), ChNo
694
695 _RL SUNANG(NCH), PDIR(NCH), PAR(NCH), ZLAI(NCH),
696 & SQSCAT(NCH), GREEN(NCH), RCUN(NCH)
697
698 _RL VGCHIL(NTYPS), VGZMEW(NTYPS),
699 & VGRST1(NTYPS), VGRST2(NTYPS), VGRST3(NTYPS)
700
701 _RL RHO4, EXTK1, EXTK2,
702 & RCINV, GAMMA, EKAT, DUM1,
703 & DUM2, DUM3, AA, BB,
704 & ZK, CC
705
706
707 DATA VGCHIL / 0.1, 0.25, 0.01, -0.3,
708 5 0.01, 0.20, 0.0, 0.0,
709 9 0.0, 0.0 /
710
711 DATA VGZMEW/ 0.9809, 0.9638, 0.9980, 1.0773,
712 5 0.9980, 0.9676, 1.000, 1.000,
713 9 1.000, 1.0000 /
714
715 DATA VGRST1 / 2335.9, 9802.2, 2869.7, 2582.0,
716 5 93989.4, 9802.2, 0.0, 0.0,
717 9 0.0 , 0.0/
718
719 DATA VGRST2 / 0.0, 10.6, 3.7, 1.1,
720 5 0.01, 10.6, 0.0, 0.0,
721 9 0.0 , 0.0/
722
723 DATA VGRST3 / 153.5, 180.0, 233.0, 110.0,
724 5 855.0, 180.0, 1.0, 1.0,
725 9 1.0, 1.0 /
726
727
728
729
730 DO 100 ChNo = 1, NCH
731
732 C**** First compute optical parameters.
733 C**** (Note: CHIL is constrained to be >= 0.01, as in SiB calcs.)
734
735 AA = 0.5 - (0.633 + 0.330*VGCHIL(ITYP(ChNo)))*VGCHIL(ITYP(ChNo))
736 BB = 0.877 * ( ONE - 2.*AA )
737 CC = ( AA + BB*SUNANG(ChNo) ) / SUNANG(ChNo)
738
739 EXTK1 = CC * SQSCAT(ChNo)
740 EXTK2 = (ONE / VGZMEW(ITYP(ChNo))) * SQSCAT(ChNo)
741
742 DUM1 = PDIR(ChNo) * CC
743 DUM2 = (ONE-PDIR(ChNo)) * ( BB*(ONE/3.+PIE/4.) + AA*1.5 )
744
745 C**** Bound extinction coefficient by 50./ZLAI:
746
747 ZK = PDIR(ChNo) *MIN( EXTK1, 50. _d 0/ZLAI(ChNo) ) +
748 & (ONE-PDIR(ChNo))*MIN( EXTK2, 50. _d 0/ZLAI(ChNo) )
749
750 C**** Now compute unstressed canopy resistance:
751
752 GAMMA = VGRST1(ITYP(ChNo)) / VGRST3(ITYP(ChNo)) +
753 & VGRST2(ITYP(ChNo))
754
755 EKAT = EXP( ZK*ZLAI(ChNo) )
756 RHO4 = GAMMA / (PAR(ChNo) * (DUM1 + DUM2))
757
758 DUM1 = (VGRST2(ITYP(ChNo)) - GAMMA) / (GAMMA + 1.E-20)
759 DUM2 = (RHO4 * EKAT + ONE) / (RHO4 + ONE)
760 DUM3 = ZK * VGRST3(ITYP(ChNo))
761 RCINV = ( DUM1*LOG(DUM2) + ZK*ZLAI(ChNo) ) / DUM3
762
763 RCUN(ChNo) = ONE / (RCINV * GREEN(ChNo) + 1.E-10)
764
765 100 CONTINUE
766
767
768 RETURN
769 END
770 C****
771 C**** [ END RCUNST ]
772 C****
773 C**** -----------------------------------------------------------------
774 C**** /////////////////////////////////////////////////////////////////
775 C**** -----------------------------------------------------------------
776 C****
777 C**** [ BEGIN SOIL ]
778 C****
779 SUBROUTINE SOIL (
780 I NCH, ITYP, WET, VGBEEX, VGPSAX, VGCSAX, DELZ,
781 O PHR, SOILCO, WETEQ
782 & )
783 C****
784 C**** This subroutine returns soil moisture potential and conductivity.
785
786 IMPLICIT NONE
787 C****
788 #include "sibber.h"
789 C****
790 INTEGER NCH
791 INTEGER ITYP(NCH), ChNo
792
793 _RL WET(NCH), PHR(NCH), SOILCO(NCH), VGPSAX(NCH),
794 & VGCSAX(NCH), VGBEEX(NCH), DELZ(NTYPS), WETEQ(NCH),
795 & WEXPB, WET0, PHEQ
796
797 DO 100 ChNo = 1, NCH
798
799 WET0 = MAX(WET(CHNO),0.01 _d 0)
800 WEXPB = WET0**VGBEEX(ChNo)
801
802 PHR(ChNo) = VGPSAX(ChNo) / WEXPB
803 SOILCO(ChNo) = VGCSAX(ChNo) * WEXPB * WEXPB * WET0 * WET0 * WET0
804
805 PHEQ = PHR(CHNO) - DELZ(ITYP(CHNO))
806 WETEQ(CHNO) = ( PHEQ/VGPSAX(CHNO) ) ** ( -1/VGBEEX(CHNO) )
807
808 100 CONTINUE
809
810 RETURN
811 END
812 C****
813 C**** [ END SOIL ]
814 C****
815 C**** -----------------------------------------------------------------
816 C**** /////////////////////////////////////////////////////////////////
817 C**** -----------------------------------------------------------------
818 C****
819 C**** [ BEGIN FLUXES ]
820 C****
821 SUBROUTINE FLUXES (
822 I NCH, ITYP, DTSTEP, ESATTC, DESDTC, ALHX,
823 I ETURB, DEDEA, DEDTC, HSTURB, DHSDEA, DHSDTC,
824 I RC, DRCDEA, DRCDTC,
825 I SWNET, HLWDWN, ALWRAD, BLWRAD, ESNFRC,
826 I TM, EM, CSOIL, PSUR, EMAXRT, WMAX,
827 U TC, TD, EA, SWET1, SNOW,
828 O RUNOFF, EVAP, SHFLUX, SMELT, HLWUP, BOMB, STRDG1,
829 O STRDG2, STRDG3, STRDG4, STRDG5, STRDG6, STRDG7,
830 O STRDG8, STRDG9
831 & )
832 C****
833 C**** This subroutine computes the fluxes of latent and sensible heat
834 C**** from the surface through an energy balance calculation.
835 C****
836 IMPLICIT NONE
837 C****
838 #include "sibber.h"
839 C****
840 INTEGER NCH
841 INTEGER ITYP(NCH), ChNo
842
843 _RL DTSTEP, ESATTC(NCH), DESDTC(NCH), ALHX(NCH),
844 & ETURB(NCH), DEDEA(NCH), DEDTC(NCH),
845 & HSTURB(NCH), DHSDEA(NCH), DHSDTC(NCH),
846 & RC(NCH), DRCDEA(NCH), DRCDTC(NCH),
847 & SWNET(NCH), HLWDWN(NCH), ALWRAD(NCH), BLWRAD(NCH),
848 & TM(NCH), EM(NCH), CSOIL(NCH), PSUR(NCH),
849 & EMAXRT(NCH), WMAX(NLAY,NTYPS),
850 & TC(NCH), TD(NCH), EA(NCH), ESNFRC(NCH),
851 & SWET1(NCH,NLAY), SNOW(NCH),
852 & RUNOFF(NCH), EVAP(NCH), SHFLUX(NCH), SMELT(NCH),
853 & HLWUP(NCH), BOMB(NCH)
854 _RL STRDG1(NCH),STRDG2(NCH),STRDG3(NCH),STRDG4(NCH)
855 _RL STRDG5(NCH),STRDG6(NCH),STRDG7(NCH),STRDG8(NCH)
856 _RL STRDG9(NCH)
857
858 _RL HLWTC, CDEEPS, Q0, RHOAIR, CONST, DHLWTC,
859 & EPLANT, A11, A12, A21, A22, F0,
860 & DEA, DTC, SNLEFT, Q0X, Q0SNOW,
861 & EANEW, ESATNW, EHARMN, DETERM, DENOM
862
863 LOGICAL CHOKE
864 _RL deepfac(ntyps)
865 DATA deepfac /1.,1.,1.,1.,1.,1.,1.,1.,1.,1./
866 C****
867 C**** -------------------------------------------------------------------
868
869 DO 200 ChNo = 1, NCH
870 C****
871 HLWTC = ALWRAD(CHNO) + BLWRAD(CHNO) * TC(CHNO)
872 CDEEPS = PIE * CSOIL(ChNo) * 2. / 86400.
873 RHOAIR = PSUR(ChNo) * 100. / (RGAS * TC(ChNo))
874 CONST = RHOAIR * EPSILON / PSUR(ChNo)
875 DHLWTC = BLWRAD(CHNO)
876 C****
877 C**** Compute matrix elements A11, A22, AND Q0 (energy balance equation).
878 C****
879
880 A11 = CSOIL(ChNo)/DTSTEP +
881 & DHLWTC +
882 & DHSDTC(ChNo) +
883 & ALHX(CHNO)*DEDTC(ChNo) +
884 & CDEEPS
885 A12 = DHSDEA(ChNo) + ALHX(CHNO) * DEDEA(ChNo)
886 Q0 = SWNET(ChNo) +
887 & HLWDWN(ChNo) -
888 & HLWTC -
889 & HSTURB(ChNo) -
890 & ALHX(CHNO) * ETURB(ChNo) -
891 & CDEEPS * (TC(ChNo) - TD(ChNo))
892
893 STRDG3(ChNo)=Q0/A11
894 STRDG4(ChNo)=(SWNET(ChNo)+HLWDWN(ChNo)-HLWTC)/A11
895 STRDG5(ChNo)=HSTURB(ChNo)/A11
896 STRDG6(ChNo)=ALHX(CHNO)*ETURB(ChNo)/A11
897 STRDG7(ChNo)=CDEEPS*(TC(ChNo) - TD(ChNo))/A11
898 STRDG9(ChNo)=A11
899
900 C****
901 C**** Compute matrix elements A21, A22, and F0 (canopy water budget
902 C**** equation) and solve for fluxes. Three cases are considered:
903 C****
904 C**** 1. Standard case: RC>0.
905 C**** 2. RC = 0. Can only occur if CIR is full or ETURB is negative.
906 C****
907 CHOKE = .TRUE.
908
909 IF( RC(CHNO) .GT. 0.) THEN
910 EPLANT = CONST * (ESATTC(ChNo) - EA(ChNo)) / RC(ChNo)
911 IF(EPLANT*ETURB(ChNo).GE.0.) THEN
912 EHARMN = 2.*EPLANT*ETURB(CHNO) / (EPLANT + ETURB(ChNo))
913 ELSE
914 EHARMN=0.
915 ENDIF
916 C****
917 C**** Some limitations to A21 and A22 are applied:
918 C**** we assume that the increase in plant evaporation
919 C**** due to an increase in either TC or EA balances
920 C**** or outweighs any decrease due to RC changes.
921 C****
922
923 A21 = -DEDTC(ChNo)*RC(ChNo) +
924 & max(0. _d 0, CONST*DESDTC(ChNo) - EHARMN*DRCDTC(ChNo) )
925 A22 = -( RC(ChNo)*DEDEA(ChNo) +
926 & max( 0. _d 0, CONST + EHARMN*DRCDEA(ChNo) ) )
927
928 F0 = RC(ChNo) * (ETURB(ChNo) - EPLANT)
929 DETERM = MIN( A12*A21/(A11*A22) - 1., -0.1 _d 0)
930 DEA = ( Q0*A21 - A11*F0 ) / ( DETERM * A11*A22 )
931 DTC = ( Q0 - A12*DEA ) / A11
932
933 STRDG1(ChNo)=DTC
934 STRDG2(ChNo)=DEA
935 STRDG8(ChNo)=-A12*DEA/A11
936
937 EVAP(ChNo) = ETURB(ChNo) + DEDEA(ChNo)*DEA + DEDTC(ChNo)*DTC
938 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
939 & DHSDTC(ChNo)*DTC
940 DENOM = DETERM * A11*A22
941 ELSE
942 CHOKE = .FALSE.
943 A21 = -DESDTC(ChNo)
944 A22 = 1.
945 F0 = ESATTC(ChNo) - EA(ChNo)
946 DEA = ( Q0*A21 - A11*F0 ) / ( A12*A21 - A11*A22 )
947 DTC = ( Q0 - A12*DEA ) / A11
948
949 STRDG1(ChNo)=DTC
950 STRDG2(ChNo)=DEA
951 STRDG8(ChNo)=-A12*DEA/A11
952
953 EVAP(ChNo) = ETURB(ChNo) + DEDEA(ChNo)*DEA + DEDTC(ChNo)*DTC
954 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
955 & DHSDTC(ChNo)*DTC
956 DENOM = A12 * A21 - A11*A22
957 ENDIF
958
959
960 C**** - - - - - - - - - - - - - - - - - - - -
961 C**** Account for snowmelt, if necessary.
962 C**** - - - - - - - - - - - - - - - - - - - -
963 SMELT(CHNO)=0.
964 SNLEFT=SNOW(CHNO)-EVAP(CHNO)*DTSTEP*ESNFRC(CHNO)
965 IF(TC(CHNO)+DTC .GT. TF .AND. SNLEFT.GT.0.) THEN
966
967 C**** First re-calculate energy balance under assumption that all
968 C**** snow is melted.
969
970 Q0X=Q0-ALHM*SNLEFT/DTSTEP
971 SMELT(CHNO)=SNLEFT/DTSTEP
972
973 DEA = ( Q0X*A21 - A11*F0 ) / DENOM
974 DTC = ( Q0X - A12*DEA ) / A11
975
976 C**** If TC+DTC is now less than TF, too much snow has melted. Now
977 C**** solve for balance assuming only some of the snow has melted.
978 C**** Set TC to TF.
979
980 IF(TC(CHNO)+DTC .LT. TF) THEN
981 DTC=TF-TC(CHNO)
982 DEA=(F0-A21*DTC)/A22
983 Q0SNOW=A11*DTC+A12*DEA
984 SMELT(CHNO)=(Q0-Q0SNOW)/ALHM
985 ENDIF
986
987 STRDG1(ChNo)=DTC
988 STRDG2(ChNo)=DEA
989 STRDG8(ChNo)=-A12*DEA/A11
990
991 EVAP(ChNo) = ETURB(ChNo) + DEDEA(ChNo)*DEA +
992 & DEDTC(ChNo)*DTC
993 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
994 & DHSDTC(ChNo)*DTC
995
996 ENDIF
997
998 C**** - - - - - - - - - - - - - - - - - - - - - - -
999 C**** Adjustments
1000
1001 C**** 1. Adjust deltas and fluxes if all available water evaporates
1002 C**** during time step:
1003 C****
1004 IF( EVAP(CHNO) .GT. EMAXRT(CHNO) ) THEN
1005 CHOKE = .FALSE.
1006 DEA = EM(CHNO) - EA(CHNO)
1007 DTC =
1008 & (Q0 + ALHX(CHNO)*(ETURB(ChNo)-EMAXRT(CHNO)) - DHSDEA(CHNO)*DEA)
1009 & / ( A11 - ALHX(CHNO)*DEDTC(ChNo) )
1010
1011 STRDG1(ChNo)=DTC
1012 STRDG2(ChNo)=DEA
1013 STRDG8(ChNo)=0.
1014
1015 EVAP(CHNO) = EMAXRT(CHNO)
1016 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
1017 & DHSDTC(ChNo)*DTC
1018 SMELT(CHNO)=0.
1019 ENDIF
1020 C****
1021 C**** Adjust DEA and DTC if solutions were pathological:
1022 C****
1023 ESATNW = ESATTC(CHNO)+DESDTC(CHNO)*DTC
1024 EANEW = EA(CHNO) + DEA
1025
1026
1027
1028 C**** 2. Pathological cases.
1029
1030 C**** Case 1: EVAP is positive in presence of negative gradient.
1031 C**** Case 2: EVAP and ETURB have opposite sign, implying that
1032 C**** "virtual effects" derivatives are meaningless and thus that we
1033 C**** do not know the proper tendency terms.
1034 C**** In both cases, assume zero evaporation for the time step.
1035
1036 C IF( ( EVAP(CHNO) .LT. 0. .AND. EM(CHNO).LT.ESATNW )
1037 C & .OR. ( EVAP(CHNO)*ETURB(CHNO) .LT. 0. ) ) THEN
1038 IF( EVAP(CHNO) .LT. 0. .AND. EM(CHNO).LT.ESATNW ) THEN
1039 CHOKE = .FALSE.
1040 DEA = EM(CHNO) - EA(CHNO)
1041 DTC = ( Q0 + ALHX(CHNO)*ETURB(ChNo) - DHSDEA(CHNO)*DEA ) /
1042 & ( A11 - ALHX(CHNO)*DEDTC(ChNo) )
1043
1044 STRDG1(ChNo)=DTC
1045 STRDG2(ChNo)=DEA
1046 STRDG8(ChNo)=0.
1047
1048 EVAP(CHNO) = 0.
1049 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
1050 & DHSDTC(ChNo)*DTC
1051
1052 IF(TC(CHNO)+DTC .GT. TF .AND. SNOW(CHNO).GT.0.) THEN
1053 Q0X=Q0-ALHM*SNOW(CHNO)/DTSTEP
1054 SMELT(CHNO)=SNOW(CHNO)/DTSTEP
1055 DTC = ( Q0X + ALHX(CHNO)*ETURB(ChNo) - DHSDEA(CHNO)*DEA ) /
1056 & ( A11 - ALHX(CHNO)*DEDTC(ChNo) )
1057 IF(TC(CHNO)+DTC .LT. TF) THEN
1058 DTC=TF-TC(CHNO)
1059 Q0SNOW=A11*DTC+A12*DEA
1060 SMELT(CHNO)=(Q0-Q0SNOW)/ALHM
1061 ENDIF
1062 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
1063 & DHSDTC(ChNo)*DTC
1064 ENDIF
1065 ENDIF
1066
1067
1068
1069 C**** 3. Exceesive dea change: apply "choke".
1070
1071 IF( CHOKE .AND. ABS(DEA) .GT. 0.5*EA(CHNO) ) THEN
1072 DEA = SIGN(.5*EA(CHNO),DEA)
1073 DTC = ( Q0 - A12*DEA ) / A11
1074
1075 STRDG1(ChNo)=DTC
1076 STRDG2(ChNo)=DEA
1077 STRDG8(ChNo)=-A12*DEA/A11
1078
1079 EVAP(ChNo) = ETURB(ChNo) + DEDEA(ChNo)*DEA + DEDTC(ChNo)*DTC
1080 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
1081 & DHSDTC(ChNo)*DTC
1082 IF(TC(CHNO)+DTC .GT. TF .AND. SNOW(CHNO).GT.0.) THEN
1083 SNLEFT=SNOW(CHNO)-EVAP(CHNO)*DTSTEP*ESNFRC(CHNO)
1084 Q0X=Q0-ALHM*SNLEFT/DTSTEP
1085 SMELT(CHNO)=SNLEFT/DTSTEP
1086 DTC = ( Q0X - A12*DEA ) / A11
1087 IF(TC(CHNO)+DTC .LT. TF) THEN
1088 DTC=TF-TC(CHNO)
1089 Q0SNOW=A11*DTC+A12*DEA
1090 SMELT(CHNO)=(Q0-Q0SNOW)/ALHM
1091 ENDIF
1092 EVAP(ChNo) = ETURB(ChNo) + DEDEA(ChNo)*DEA + DEDTC(ChNo)*DTC
1093 SHFLUX(ChNo) = HSTURB(ChNo) + DHSDEA(ChNo)*DEA +
1094 & DHSDTC(ChNo)*DTC
1095 ENDIF
1096 ENDIF
1097
1098 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1099
1100 TC(ChNo) = TC(ChNo) + DTC
1101 EA(ChNo) = EA(ChNo) + DEA
1102 TD(ChNo) = TD(ChNo) +
1103 c CHANGED THIS: deep layer 2 times deeper
1104 & DTSTEP * CDEEPS * (TC(ChNo) - TD(ChNo)) /
1105 . (CSOIL(ChNo)*67.73*deepfac(ityp(ChNo)))
1106 HLWUP(CHNO) = HLWTC + DHLWTC*DTC
1107
1108 SNOW(CHNO)=SNOW(CHNO)-SMELT(CHNO)*DTSTEP
1109
1110 C**** Make sure EA remains positive
1111
1112 EA(CHNO) = MAX(EA(CHNO), 0.0 _d 0)
1113
1114 200 CONTINUE
1115
1116
1117
1118 RETURN
1119 END
1120 C****
1121 C**** [ END FLUXES ]
1122 C****
1123 C**** -----------------------------------------------------------------
1124 C**** /////////////////////////////////////////////////////////////////
1125 C**** -----------------------------------------------------------------
1126 C****
1127 C**** [ BEGIN VPDFAC ]
1128 C****
1129 SUBROUTINE VPDFAC (
1130 I NCH, ITYP, ESATTC, EA,
1131 O VPDSTR
1132 & )
1133 C****
1134 C**** This subroutine computes the vapor pressure deficit stress.
1135 C****
1136 IMPLICIT NONE
1137 C****
1138 #include "sibber.h"
1139 C****
1140 INTEGER NCH
1141 INTEGER ITYP(NCH), ChNo
1142 _RL ESATTC(NCH), EA(NCH), VPDSTR(NCH)
1143 c _RL VGDFAC(NTYPS)
1144 C****
1145 c DATA VGDFAC / .0273, .0357, .0310, .0238,
1146 c 5 .0275, .0275, 0., 0.,
1147 c 9 0., 0. /
1148 C****
1149 C**** -----------------------------------------------------------------
1150
1151 DO 100 ChNo = 1, NCH
1152 C****
1153 c VPDSTR(ChNo) = 1. - (ESATTC(ChNo)-EA(ChNo)) * VGDFAC(ITYP(ChNo))
1154 c VPDSTR (ChNo) = MIN( 1. _d 0, MAX( VPDSTR(ChNo), 1. _d -10 ) )
1155 VPDSTR(CHNO) = 1.
1156 C****
1157 100 CONTINUE
1158 C****
1159 RETURN
1160 END
1161 C****
1162 C**** [ END VPDFAC ]
1163 C****
1164 C**** -----------------------------------------------------------------
1165 C**** /////////////////////////////////////////////////////////////////
1166 C**** -----------------------------------------------------------------
1167 C****
1168 C**** [ BEGIN TMPFAC ]
1169 C****
1170 SUBROUTINE TMPFAC (
1171 I NCH, ITYP, TC,
1172 O FTEMP
1173 & )
1174 C****
1175 C**** Compute temperature stress factor.
1176 C****
1177 IMPLICIT NONE
1178 C****
1179 #include "sibber.h"
1180 C****
1181 INTEGER NCH
1182 INTEGER ITYP(NCH), ChNo, TypPtr
1183 _RL TC(NCH), FTEMP(NCH)
1184 _RL VGTLL(MemFac*NTYPS), VGTU(MemFac*NTYPS),
1185 & VGTCF1(MemFac*NTYPS), VGTCF2(MemFac*NTYPS),
1186 & VGTCF3(MemFac*NTYPS)
1187 C****
1188 DATA VGTLL /MemFac*273., MemFac*273., MemFac*268., MemFac*283.,
1189 5 MemFac*283., MemFac*273., MemFac* 0., MemFac* 0.,
1190 9 MemFac* 0., MemFac* 0. /
1191 DATA VGTU /MemFac*318., MemFac*318., MemFac*313., MemFac*328.,
1192 5 MemFac*323., MemFac*323., MemFac* 0., MemFac* 0.,
1193 9 MemFac* 0. , MemFac* 0./
1194 DATA VGTCF1 / MemFac*-1.43549E-06, MemFac*-6.83584E-07,
1195 3 MemFac* 1.67699E-07, MemFac*-1.43465E-06,
1196 5 MemFac*-2.76097E-06, MemFac*-1.58094E-07,
1197 7 MemFac* 0., MemFac* 0.,
1198 9 MemFac* 0. , MemFac* 0./
1199 DATA VGTCF2 / MemFac* 7.95859E-04, MemFac* 3.72064E-04,
1200 3 MemFac*-7.65944E-05, MemFac* 8.24060E-04,
1201 5 MemFac* 1.57617E-03, MemFac* 8.44847E-05,
1202 7 MemFac* 0., MemFac* 0.,
1203 9 MemFac* 0. , MemFac* 0./
1204 DATA VGTCF3 / MemFac*-1.11575E-01, MemFac*-5.21533E-02,
1205 3 MemFac* 6.14960E-03, MemFac*-1.19602E-01,
1206 5 MemFac*-2.26109E-01, MemFac*-1.27272E-02,
1207 7 MemFac* 0., MemFac* 0.,
1208 9 MemFac* 0. , MemFac* 0./
1209 C****
1210 C**** ----------------------------------------------------------------
1211
1212 DO 100 ChNo = 1, NCH
1213 C****
1214 TypPtr = MOD(ChNo,MemFac) + (ITYP(ChNo)-1)*MemFac + 1
1215 FTEMP(ChNo) = (TC(ChNo) - VGTLL(TypPtr)) *
1216 & (TC(ChNo) - VGTU(TypPtr)) *
1217 & ( VGTCF1(TypPtr)*TC(ChNo)*TC(ChNo) +
1218 & VGTCF2(TypPtr)*TC(ChNo) +
1219 & VGTCF3(TypPtr) )
1220 IF ( TC(ChNo) .LE. VGTLL(TypPtr) .OR. TC(ChNo) .GE. VGTU(TypPtr) )
1221 & FTEMP (ChNo) = 1.E-10
1222 FTEMP(CHNO) = MIN( 1. _d 0, MAX( FTEMP(ChNo), 1. _d -10 ) )
1223 C****
1224 100 CONTINUE
1225 C****
1226 RETURN
1227 END
1228 C****
1229 C**** [ END TMPFAC ]
1230 C****
1231 C**** -----------------------------------------------------------------
1232 C**** /////////////////////////////////////////////////////////////////
1233 C**** -----------------------------------------------------------------
1234 C****
1235 C**** [ BEGIN SMFAC ]
1236 C****
1237 SUBROUTINE SMFAC (
1238 I NCH, ITYP, ESATTC, EA, PHR, SOILCO,
1239 I RCUN, VPDSTR, FTEMP, TC, PSUR, Z2,
1240 I RSOIL1, RSOIL2, VGPH1X, VGPH2X, VGRPLX,
1241 O RC
1242 & )
1243 C****
1244 C**** This subroutine estimates RC after computing the
1245 C**** leaf water potential stress.
1246 C****
1247 IMPLICIT NONE
1248 C****
1249 #include "sibber.h"
1250 C****
1251 INTEGER NCH
1252 INTEGER ITYP(NCH), ChNo
1253 _RL ESATTC(NCH), EA(NCH), PHR(NCH), SOILCO(NCH),
1254 & RCUN(NCH), VPDSTR(NCH), FTEMP(NCH), TC(NCH),
1255 & PSUR(NCH), Z2(NCH), RSOIL1(NCH), RSOIL2(NCH),
1256 & VGPH1X(NCH), VGPH2X(NCH), VGRPLX(NCH), RC(NCH)
1257 _RL RCUNTD, RHOAIR, CONST, DEF, D12, DR2,
1258 & RSOIL, R0, EEST, RSTFAC
1259 C****
1260 C**** -----------------------------------------------------------------
1261
1262 DO 100 ChNo = 1, NCH
1263 C****
1264 RCUNTD = RCUN(ChNo) / ( VPDSTR(ChNo)*FTEMP(ChNo) )
1265 RHOAIR = PSUR(ChNo) * 100. / ( RGAS*TC(ChNo) )
1266
1267 C****
1268 CONST = RHOAIR * EPSILON / PSUR(ChNo)
1269 DEF = ( ESATTC(ChNo) - EA(ChNo) ) * CONST
1270 D12 = VGPH1X(ChNo) - VGPH2X(ChNo)
1271 DR2 = PHR(ChNo) - Z2(ChNo) - VGPH2X(ChNo)
1272 RSOIL = RSOIL1(ChNo) + RSOIL2(ChNo)/SOILCO(ChNo)
1273 R0 = ( VGRPLX(ChNo) + RSOIL ) / RHOW
1274 EEST = DEF*DR2 / ( RCUNTD*D12 + DEF*R0 )
1275 RSTFAC = ( DR2 - R0*EEST ) / D12
1276 RSTFAC = MIN( 1. _d 0, MAX( 0.001 _d 0, RSTFAC ) )
1277 RC(ChNo) = RCUNTD / RSTFAC
1278 C****
1279 100 CONTINUE
1280 C****
1281 RETURN
1282 END
1283 C****
1284 C**** [ END SMFAC ]
1285 C****
1286 C**** -----------------------------------------------------------------
1287 C**** /////////////////////////////////////////////////////////////////
1288 C**** -----------------------------------------------------------------
1289 C****
1290 C**** [ BEGIN RSURFP ]
1291 C****
1292 SUBROUTINE RSURFP (
1293 I NCH, UM, U2FAC, Z2, RDC, WET, ESATTC, EA,
1294 U RC,
1295 O RX1, RX2
1296 & )
1297 C****
1298 IMPLICIT NONE
1299 INTEGER NCH
1300 INTEGER ChNo
1301 _RL UM(NCH), U2FAC(NCH), Z2(NCH), RDC(NCH),
1302 & WET(NCH), ESATTC(NCH), EA(NCH),
1303 & RC(NCH), RX1(NCH), RX2(NCH)
1304 _RL U2, RSURF, HESAT
1305 C****
1306 C**** -----------------------------------------------------------------
1307
1308 DO 100 ChNo = 1, NCH
1309 C****
1310 U2 = UM(ChNo) * U2FAC(ChNo)
1311 c RSURF = RDC(ChNo) / U2 + 30. / (1.E-20 + WET(ChNo))
1312 RSURF = RDC(ChNo) / U2 + 26. + 6. / (1.E-20 + WET(ChNo))**2
1313
1314 C**** Account for subsaturated humidity at soil surface:
1315 C****
1316 HESAT = ESATTC(CHNO) * MIN( 1. _d 0, WET(CHNO)*2. _d 0)
1317 IF( EA(CHNO) .LT. HESAT ) THEN
1318 RSURF=RSURF*( 1. + (ESATTC(CHNO)-HESAT)/(HESAT-EA(CHNO)) )
1319 ELSE
1320 RSURF=1.E10
1321 ENDIF
1322
1323
1324 RX1(CHNO)=RC(CHNO)
1325 RX2(CHNO)=RSURF
1326
1327 RC(ChNo) = RC(CHNO) * RSURF / ( RC(ChNo) + RSURF )
1328 C****
1329 100 CONTINUE
1330 C****
1331 RETURN
1332 END
1333 C****
1334 C**** [ END RSURFP ]
1335 C****
1336 C**** -----------------------------------------------------------------
1337 C**** /////////////////////////////////////////////////////////////////
1338 C**** -----------------------------------------------------------------
1339 C****
1340 C**** [ BEGIN RCANOP ]
1341 C****
1342 SUBROUTINE RCANOP (
1343 I NCH, CAPAC, SNOW, SATCAP, RA, ETURB, SNWFRC,
1344 I POTFRC,
1345 U RC
1346 & )
1347 C****
1348 C**** The effective latent heat resistance RC depends on the quantity
1349 C**** of interception reservoir water and the snow cover. POTFRC
1350 C**** is the fraction of the tile from which potential evaporation
1351 C**** occurs.
1352 C****
1353 IMPLICIT NONE
1354 INTEGER NCH
1355 INTEGER ChNo
1356 _RL CAPAC(NCH), SNOW(NCH), SATCAP(NCH), RA(NCH), ETURB(NCH),
1357 & RC(NCH), SNWFRC(NCH), POTFRC(NCH)
1358 _RL ETCRIT,RAMPFC
1359
1360 C**** (Note: ETCRIT arbitrarily set to ~-5 W/m2, or -2.e-6 mm/sec.)
1361 DATA ETCRIT/ -2.E-6 /
1362 C****
1363 C**** -----------------------------------------------------------------
1364
1365 DO 100 ChNo = 1, NCH
1366
1367 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1368 C**** Case 1: Vegetation present (SATCAP GT .001). Potential evap. from
1369 c**** both interception reservoir and snow.
1370
1371 IF(SATCAP(CHNO).GT..001) THEN
1372 RC(ChNo)=RC(ChNo)*(1.-POTFRC(CHNO))/
1373 & ( 1.+POTFRC(CHNO)*RC(ChNo)/RA(ChNo) )
1374 ENDIF
1375
1376 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1377 C**** Case 2: Vegetation absent (SATCAP LE .001). Potential evap. from
1378 c**** snow only.
1379
1380 IF(SATCAP(CHNO) .LE. .001) THEN
1381 RC(ChNo)=RC(ChNo)*(1.-SNWFRC(CHNO))/
1382 & ( 1.+SNWFRC(CHNO)*RC(ChNo)/RA(ChNo) )
1383 ENDIF
1384
1385 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1386 C**** Assume RC=0 for condensation (dew).
1387 C**** RAMPFC is used to ensure continuity in RC.
1388
1389 RAMPFC=ETURB(CHNO)/ETCRIT
1390 IF ( RAMPFC .GE. 0. ) RC(CHNO) = RC(CHNO)*(1.-RAMPFC)
1391 IF ( RAMPFC .GT. 1. ) RC(ChNo) = 0.
1392 C****
1393 100 CONTINUE
1394 C****
1395 RETURN
1396 END
1397 C****
1398 C**** [ END RCANOP ]
1399 C****
1400 C**** -----------------------------------------------------------------
1401 C**** /////////////////////////////////////////////////////////////////
1402 C**** -----------------------------------------------------------------
1403 C****
1404 C***** [ BEGIN WUPDAT ]
1405 C****
1406 SUBROUTINE WUPDAT (
1407 I NCH, ITYP, DTSTEP,
1408 I EVAP, SATCAP, VGWMAX, TC, RA, RC,
1409 I RX1, RX2,ESNFRC,EIRFRC,
1410 U CAPAC, SNOW, SWET, RUNOFF,
1411 O EINT, ESOI, EVEG, ESNO, RUNSRF, FWSOIL
1412 & )
1413 C****
1414 C**** This subroutine allows evapotranspiration to adjust the water
1415 C**** contents of the interception reservoir and the soil layers.
1416 C****
1417 IMPLICIT NONE
1418 C****
1419 #include "sibber.h"
1420 C****
1421 INTEGER NCH
1422 INTEGER ITYP(NCH), ChNo
1423 _RL EVAP(NCH), SATCAP(NCH), VGWMAX(NLAY,NTYPS),
1424 & TC(NCH), RA(NCH), RC(NCH),
1425 & CAPAC(NCH), SNOW(NCH), SWET(nch,NLAY),
1426 & RUNOFF(NCH), RX1(NCH),
1427 & RX2(NCH), RUNSRF(NCH), FWSOIL(NCH),
1428 & ESNFRC(NCH), EIRFRC(NCH)
1429 _RL EINT(NCH), ESOI(NCH), EVEG(NCH), ESNO(NCH)
1430 _RL DTSTEP, EGRO, FWS, THRU, DEWRUN,
1431 & WTOTAL,WLAY1,WLAY2,ELAY1,ELAY2,EGROI
1432 C****
1433 C**** -----------------------------------------------------------------
1434
1435 DO 100 ChNo = 1, NCH
1436 C****
1437 C**** Partition evap between interception, snow, and ground reservoirs.
1438 C****
1439 WLAY1 = SWET(ChNo,SFCLY) * VGWMAX(SFCLY,ITYP(ChNo))
1440 WLAY2 = SWET(ChNo,ROOTLY) * VGWMAX(ROOTLY,ITYP(ChNo))
1441 WTOTAL = WLAY1 + WLAY2
1442 C****
1443 ESNO(CHNO)=ESNFRC(CHNO)*EVAP(CHNO)*DTSTEP
1444 EINT(CHNO)=EIRFRC(CHNO)*EVAP(CHNO)*DTSTEP
1445 EGRO = EVAP(CHNO)*DTSTEP - ESNO(CHNO) - EINT(CHNO)
1446
1447 C**** Ensure that individual capacities are not exceeded.
1448
1449 IF(ESNO(CHNO) .GT. SNOW(CHNO)) THEN
1450 EINT(CHNO)=EINT(CHNO)+(ESNO(CHNO)-SNOW(CHNO))
1451 ESNO(CHNO)=SNOW(CHNO)
1452 ENDIF
1453 EGROI=EGRO+EINT(CHNO)
1454 IF(EGROI .GT. CAPAC(CHNO)+WTOTAL) THEN
1455 ESNO(CHNO)=ESNO(CHNO)+EGROI-(CAPAC(CHNO)+WTOTAL)
1456 EGROI=CAPAC(CHNO)+WTOTAL
1457 ENDIF
1458
1459 EINT(CHNO)=EGROI-EGRO
1460 IF(EINT(CHNO) .GT. CAPAC(CHNO)) THEN
1461 EGRO=EGRO+EINT(CHNO)-CAPAC(CHNO)
1462 EINT(CHNO)=CAPAC(CHNO)
1463 ENDIF
1464 IF(EGRO .GT. WTOTAL) THEN
1465 EINT(CHNO)=EINT(CHNO)+EGRO-WTOTAL
1466 EGRO=WTOTAL
1467 ENDIF
1468 C****
1469 C**** Separate egro into surface-evaporation/transpiration components:
1470 C****
1471 IF( RX1(CHNO)+RX2(CHNO) .NE. 0. ) THEN
1472 ESOI(CHNO)=EGRO*RX1(CHNO)/(RX1(CHNO)+RX2(CHNO))
1473 EVEG(CHNO)=EGRO - ESOI(CHNO)
1474 ELSE
1475 ESOI(CHNO)=EGRO/2.
1476 EVEG(CHNO)=EGRO/2.
1477 ENDIF
1478 C****
1479 C**** Translate ESOI and EVEG into evaporation fluxes from layers 1 and 2:
1480 C****
1481 IF( WTOTAL .GT. 0. ) THEN
1482 FWS = EVEG(CHNO) / WTOTAL
1483 ELSE
1484 FWS = 1.
1485 ENDIF
1486 ELAY1 = WLAY1*FWS + ESOI(CHNO)
1487 ELAY2 = WLAY2*FWS
1488 C****
1489 C**** Ensure that enough soil water is available in each layer:
1490 C****
1491 IF( ELAY1 .GT. WLAY1 ) THEN
1492 ELAY2 = ELAY2 + (ELAY1 - WLAY1)
1493 ELAY1 = WLAY1
1494 ENDIF
1495 IF( ELAY2 .GT. WLAY2 ) THEN
1496 ELAY1 = ELAY1 + (ELAY2 - WLAY2)
1497 ELAY2 = WLAY2
1498 ENDIF
1499
1500 C**** Special case for condensation:
1501 IF(EVAP(CHNO) .LT. 0.) THEN
1502 EINT(CHNO)=(1.-ESNFRC(CHNO))*EVAP(CHNO)*DTSTEP
1503 ELAY1=0.
1504 ELAY2=0.
1505 ESOI(CHNO)=0.
1506 EVEG(CHNO)=0.
1507 ESNO(CHNO)=ESNFRC(CHNO)*EVAP(CHNO)*DTSTEP
1508 ENDIF
1509
1510 C****
1511 C**** Remove moisture from reservoirs:
1512 C****
1513 SNOW(ChNo) = SNOW(ChNo) - ESNO(CHNO)
1514 CAPAC(ChNo) = CAPAC(ChNo) - EINT(CHNO)
1515 SWET(ChNo,SFCLY) = (WLAY1 - ELAY1) / VGWMAX(SFCLY,ITYP(ChNo))
1516 SWET(ChNo,ROOTLY) = (WLAY2 - ELAY2) / VGWMAX(ROOTLY,ITYP(CHNO))
1517 C****
1518 C**** Ensure against numerical precision problems:
1519 SWET(ChNo,SFCLY) = MIN( 1. _d 0, MAX( 0. _d 0, SWET(ChNo,SFCLY) )
1520 $ )
1521 SWET(ChNo,ROOTLY) = MIN( 1. _d 0, MAX( 0. _d 0, SWET(ChNo,ROOTLY)
1522 $ ) )
1523 C****
1524 C****
1525 C**** -------------------------------------------------
1526 C**** DEWFALL:
1527 C****
1528 C**** If dewfall adds to cir, insure that it does not fill
1529 C**** beyond capacity. If resulting throughfall adds to top soil layer,
1530 C**** insure that it also does not fill beyond capacity.
1531 C****
1532 IF( CAPAC(ChNo) .GT. SATCAP(ChNo) ) THEN
1533 THRU = CAPAC(ChNo) - SATCAP(ChNo)
1534 CAPAC(ChNo) = SATCAP(ChNo)
1535 SWET(ChNo,SFCLY) = SWET(ChNo,SFCLY) +
1536 & THRU / VGWMAX(SFCLY,ITYP(ChNo))
1537 DEWRUN=0.
1538 IF ( SWET(ChNo,SFCLY) .GT. 1. ) THEN
1539 DEWRUN = ( SWET(ChNo,SFCLY) - 1. ) * VGWMAX(SFCLY,ITYP(ChNo))
1540 SWET(ChNo,SFCLY) = 1.
1541 RUNOFF(ChNo) = RUNOFF(ChNo) + DEWRUN/DTSTEP
1542 RUNSRF(ChNo) = RUNSRF(ChNo) + DEWRUN/DTSTEP
1543 ENDIF
1544 FWSOIL(CHNO) = FWSOIL(CHNO) + (THRU-DEWRUN)/DTSTEP
1545 ENDIF
1546 C****
1547 100 CONTINUE
1548 C****
1549 RETURN
1550 END
1551 C****
1552 C**** [ END WUPDAT ]
1553 C****
1554 C**** -----------------------------------------------------------------
1555 C**** /////////////////////////////////////////////////////////////////
1556 C**** -----------------------------------------------------------------
1557 C****
1558 C**** [ BEGIN GWATER ]
1559 C****
1560 SUBROUTINE GWATER (
1561 I NCH, ITYP, WSMAX, PHLAY, AKLAY, TC,
1562 I DTSTEP, VGZDEX, VGSLOX, WETEQ1, WETEQ2,
1563 I PHSAT, AKSAT,
1564 U SWET, RUNOFF, GDRAIN
1565 & )
1566 C****
1567 C**** This subroutine computes diffusion between soil layers.
1568 C****
1569 IMPLICIT NONE
1570 C****
1571 #include "sibber.h"
1572 C****
1573 INTEGER NCH
1574 INTEGER ITYP(NCH), ChNo
1575 _RL VGSLOX(NCH), RUNOFF(NCH), GDRAIN(NCH)
1576 _RL ZDEP12, AKAVE, GWFLUX, ZDEP23, HALFMX, DHDZ,
1577 & FAREA, TFM2, FRAMP
1578
1579 _RL WSMAX(NLAY,NTYPS), PHLAY(nch,NLAY),
1580 & AKLAY(nch,NLAY), TC(NCH),
1581 & DTSTEP, SWET(nch,NLAY),
1582 & VGZDEX(NLAY,nch), WETEQ1(NCH), WETEQ2(NCH),
1583 & PHSAT(NCH), AKSAT(NCH)
1584
1585 C****
1586 C**** ----------------------------------------------------------------
1587
1588 DO 100 ChNo = 1, NCH
1589 C****
1590 C**** Diffusion between layer 1 and 2:
1591 ZDEP12 = VGZDEX(SFCLY,ChNo) + VGZDEX(ROOTLY,ChNo)
1592 IF( SWET(CHNO,SFCLY) .GT. WETEQ1(CHNO) ) THEN
1593 FAREA=(SWET(CHNO,SFCLY) - WETEQ1(CHNO)) / (1. - WETEQ1(CHNO))
1594 DHDZ = 1. + 2.*(PHSAT(CHNO)-PHLAY(CHNO,ROOTLY))/ZDEP12
1595 AKAVE = AKSAT(CHNO)
1596 ELSE
1597 FAREA = 1.
1598 DHDZ = 1. + 2.*(PHLAY(ChNo,SFCLY)-PHLAY(ChNo,ROOTLY))/ZDEP12
1599 AKAVE = AKLAY(ChNo,ROOTLY)
1600 ENDIF
1601 C****
1602 GWFLUX = 1000. * DTSTEP * AKAVE * DHDZ * FAREA
1603 C****
1604 C**** Test for limits on water holding capacity.
1605 C****
1606 HALFMX=0.5*ABS( SWET(CHNO,SFCLY)-WETEQ1(CHNO) )
1607 & * WSMAX(SFCLY,ITYP(CHNO))
1608 IF (GWFLUX .GE. 0.) THEN
1609 GWFLUX = MIN( GWFLUX,
1610 & SWET(ChNo,SFCLY) * WSMAX(SFCLY,ITYP(ChNo)) )
1611 GWFLUX = MIN( GWFLUX, WSMAX(ROOTLY,ITYP(ChNo)) -
1612 & SWET(ChNo,ROOTLY) * WSMAX(ROOTLY,ITYP(ChNo)) )
1613 GWFLUX = MIN( GWFLUX, HALFMX )
1614 ELSE
1615 GWFLUX = -MIN( ABS(GWFLUX),
1616 & SWET(ChNo,ROOTLY) * WSMAX(ROOTLY,ITYP(ChNo)) )
1617 GWFLUX = -MIN( ABS(GWFLUX), WSMAX(SFCLY,ITYP(ChNo)) -
1618 & SWET(ChNo,SFCLY) * WSMAX(SFCLY,ITYP(ChNo)) )
1619 GWFLUX = -MIN( ABS(GWFLUX), HALFMX )
1620 ENDIF
1621 C****
1622 C**** Prevent diffusion when ground is frozen (INCLUDE RAMPING):
1623 TFM2=TF-2.
1624 FRAMP=(TC(CHNO)-TFM2)/2.
1625 FRAMP=MIN(1. _d 0, MAX(0. _d 0,FRAMP) )
1626 GWFLUX=GWFLUX*FRAMP
1627 C****
1628 C**** Update water contents
1629 SWET(ChNo,SFCLY) = SWET(ChNo,SFCLY) -
1630 & GWFLUX / WSMAX(SFCLY,ITYP(ChNo))
1631 SWET(ChNo,ROOTLY) = SWET(ChNo,ROOTLY) +
1632 & GWFLUX / WSMAX(ROOTLY,ITYP(ChNo))
1633 C****
1634 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1635 C**** Diffusion between root and recharge layers:
1636 ZDEP23 = VGZDEX(ROOTLY,ChNo) + VGZDEX(RECHLY,ChNo)
1637 DHDZ=1. + 2. * (PHLAY(ChNo,ROOTLY) - PHLAY(ChNo,RECHLY)) / ZDEP23
1638 IF(DHDZ.GE.0.) THEN
1639 AKAVE = AKLAY(ChNo,ROOTLY)
1640 ELSE
1641 AKAVE = AKLAY(ChNo,RECHLY)
1642 ENDIF
1643 GWFLUX = 1000. * DTSTEP * AKAVE * DHDZ
1644 C****
1645 C**** Test for limits on water holding capacity:
1646 HALFMX=0.5*ABS( SWET(CHNO,ROOTLY)-WETEQ2(CHNO) )
1647 & * WSMAX(ROOTLY,ITYP(CHNO))
1648 IF (GWFLUX .GE. 0.) THEN
1649 GWFLUX = MIN( GWFLUX,
1650 & SWET(ChNo,ROOTLY) * WSMAX(ROOTLY,ITYP(ChNo)) )
1651 GWFLUX = MIN( GWFLUX, WSMAX(RECHLY,ITYP(ChNo)) -
1652 & SWET(ChNo,RECHLY) * WSMAX(RECHLY,ITYP(ChNo)) )
1653 GWFLUX = MIN( GWFLUX, HALFMX )
1654 ELSE
1655 GWFLUX = -MIN( ABS(GWFLUX),
1656 & SWET(ChNo,RECHLY) * WSMAX(RECHLY,ITYP(ChNo)) )
1657 GWFLUX = -MIN( ABS(GWFLUX), WSMAX(ROOTLY,ITYP(ChNo)) -
1658 & SWET(ChNo,ROOTLY) * WSMAX(ROOTLY,ITYP(ChNo)) )
1659 GWFLUX = -MIN( ABS(GWFLUX), HALFMX )
1660 ENDIF
1661 C****
1662 C**** Prevent diffusion when ground is frozen (INCLUDE RAMPING):
1663 FRAMP=(TC(CHNO)-TFM2)/2.
1664 FRAMP=MIN(1. _d 0, MAX(0. _d 0,FRAMP) )
1665 GWFLUX=GWFLUX*FRAMP
1666 C****
1667 C**** Update water contents
1668 SWET(ChNo,ROOTLY) = SWET(ChNo,ROOTLY) -
1669 & GWFLUX / WSMAX(ROOTLY,ITYP(ChNo))
1670 SWET(ChNo,RECHLY) = SWET(ChNo,RECHLY) +
1671 & GWFLUX / WSMAX(RECHLY,ITYP(ChNo))
1672 GDRAIN(CHNO)=GWFLUX/DTSTEP
1673 C****
1674 C**** - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1675 C**** Flux of moisture out of lower soil layer to water table
1676 C**** (approximation to SiB)
1677 C****
1678 GWFLUX = VGSLOX(ChNo) * AKLAY(ChNo,RECHLY) * 1000. * DTSTEP
1679
1680 C**** Prevent diffusion when ground is frozen (INCLUDE RAMPING):
1681 FRAMP=(TC(CHNO)-TFM2)/2.
1682 FRAMP=MIN(1. _d 0, MAX(0. _d 0,FRAMP) )
1683 GWFLUX=GWFLUX*FRAMP
1684
1685 GWFLUX = MIN( GWFLUX,
1686 & SWET(ChNo,RECHLY) * WSMAX(RECHLY,ITYP(ChNo)) )
1687 SWET(ChNo,RECHLY) = SWET(ChNo,RECHLY) -
1688 & GWFLUX / WSMAX(RECHLY,ITYP(ChNo))
1689 C****
1690 RUNOFF (ChNo) = RUNOFF (ChNo) + GWFLUX/DTSTEP
1691 C****
1692 100 CONTINUE
1693 C****
1694 RETURN
1695 END
1696 C****
1697 C**** [ END GWATER ]
1698

  ViewVC Help
Powered by ViewVC 1.1.22