/[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.2 - (show annotations) (download)
Fri Jul 16 19:52:21 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.1: +19 -16 lines
Fizhi land model include files

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

  ViewVC Help
Powered by ViewVC 1.1.22