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

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

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


Revision 1.21 - (show annotations) (download)
Sun Aug 29 19:48:06 2004 UTC (19 years, 9 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint55c_post, checkpoint54e_post, checkpoint55d_pre, checkpoint55j_post, checkpoint56b_post, checkpoint55h_post, checkpoint55b_post, checkpoint56c_post, checkpoint55, checkpoint57a_post, checkpoint54f_post, checkpoint55g_post, checkpoint55f_post, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, checkpoint55e_post, checkpoint55a_post, checkpoint56a_post, checkpoint55d_post
Changes since 1.20: +2 -2 lines
Fix cloud fraction diagnostic counter - increment when bi=bj=1

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_swrad.F,v 1.20 2004/08/26 14:41:56 molod Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 subroutine swrio (nymd,nhms,bi,bj,ndswr,myid,istrip,npcs,
6 . low_level,mid_level,im,jm,lm,
7 . pz,plz,plze,dpres,pkht,pkz,tz,qz,oz,co2,
8 . albvisdr,albvisdf,albnirdr,albnirdf,
9 . dtradsw,dtswclr,radswg,swgclr,
10 . fdifpar,fdirpar,osr,osrclr,
11 . ptop,nswcld,cldsw,cswmo,nswlz,swlz,
12 . lpnt,imstturb,qliqave,fccave,landtype,xlats,xlons)
13
14 implicit none
15 #ifdef ALLOW_DIAGNOSTICS
16 #include "SIZE.h"
17 #include "diagnostics_SIZE.h"
18 #include "diagnostics.h"
19 #endif
20
21 c Input Variables
22 c ---------------
23 integer nymd,nhms,bi,bj,ndswr,myid,istrip,npcs
24 integer mid_level,low_level
25 integer im,jm,lm
26 _RL ptop
27 _RL pz(im,jm),plz(im,jm,lm),plze(im,jm,lm+1),dpres(im,jm,lm)
28 _RL pkht(im,jm,lm+1),pkz(im,jm,lm)
29 _RL tz(im,jm,lm),qz(im,jm,lm)
30 _RL oz(im,jm,lm)
31 _RL co2
32 _RL albvisdr(im,jm),albvisdf(im,jm),albnirdr(im,jm)
33 _RL albnirdf(im,jm)
34 _RL radswg(im,jm),swgclr(im,jm),fdifpar(im,jm),fdirpar(im,jm)
35 _RL osr(im,jm),osrclr(im,jm),dtradsw(im,jm,lm)
36 _RL dtswclr(im,jm,lm)
37 integer nswcld,nswlz
38 _RL cldsw(im,jm,lm),cswmo(im,jm,lm),swlz(im,jm,lm)
39 logical lpnt
40 integer imstturb
41 _RL qliqave(im,jm,lm),fccave(im,jm,lm)
42 integer landtype(im,jm)
43 _RL xlats(im,jm),xlons(im,jm)
44
45 c Local Variables
46 c ---------------
47 integer i,j,L,nn,nsecf
48 integer ntmstp,nymd2,nhms2
49 _RL getcon,grav,cp,undef
50 _RL ra,alf,reffw,reffi,tminv
51
52 parameter ( reffw = 10.0 )
53 parameter ( reffi = 65.0 )
54
55 _RL tdry(im,jm,lm)
56 _RL TEMP1(im,jm)
57 _RL TEMP2(im,jm)
58 _RL zenith (im,jm)
59 _RL cldtot (im,jm,lm)
60 _RL cldmxo (im,jm,lm)
61 _RL totcld (im,jm)
62 _RL cldlow (im,jm)
63 _RL cldmid (im,jm)
64 _RL cldhi (im,jm)
65 _RL taulow (im,jm)
66 _RL taumid (im,jm)
67 _RL tauhi (im,jm)
68 _RL tautype(im,jm,lm,3)
69 _RL tau(im,jm,lm)
70 _RL albedo(im,jm)
71
72 _RL PK(ISTRIP,lm)
73 _RL qzl(ISTRIP,lm),CLRO(ISTRIP,lm)
74 _RL TZL(ISTRIP,lm)
75 _RL OZL(ISTRIP,lm)
76 _RL PLE(ISTRIP,lm+1)
77 _RL COSZ(ISTRIP)
78 _RL dpstrip(ISTRIP,lm)
79
80 _RL albuvdr(ISTRIP),albuvdf(ISTRIP)
81 _RL albirdr(ISTRIP),albirdf(ISTRIP)
82 _RL difpar (ISTRIP),dirpar (ISTRIP)
83
84 _RL fdirir(istrip),fdifir(istrip)
85 _RL fdiruv(istrip),fdifuv(istrip)
86
87 _RL flux(istrip,lm+1)
88 _RL fluxclr(istrip,lm+1)
89 _RL dtsw(istrip,lm)
90 _RL dtswc(istrip,lm)
91
92 _RL taul (istrip,lm)
93 _RL reff (istrip,lm,2)
94 _RL tauc (istrip,lm,2)
95 _RL taua (istrip,lm)
96 _RL tstrip (istrip)
97
98 logical first
99 data first /.true./
100
101 C **********************************************************************
102 C **** INITIALIZATION ****
103 C **********************************************************************
104
105 grav = getcon('GRAVITY')
106 cp = getcon('CP')
107 undef = getcon('UNDEF')
108
109 NTMSTP = nsecf(NDSWR)
110 TMINV = 1./float(ntmstp)
111
112 C Compute Temperature from Theta
113 C ------------------------------
114 do L=1,lm
115 do j=1,jm
116 do i=1,im
117 tdry(i,j,L) = tz(i,j,L)*pkz(i,j,L)
118 enddo
119 enddo
120 enddo
121
122 if (first .and. myid.eq.1 ) then
123 print *
124 print *,'Low-Level Clouds are Grouped between levels: ',
125 . lm,' and ',low_level
126 print *,'Mid-Level Clouds are Grouped between levels: ',
127 . low_level-1,' and ',mid_level
128 print *
129 first = .false.
130 endif
131
132 C **********************************************************************
133 C **** CALCULATE COSINE OF THE ZENITH ANGLE ****
134 C **********************************************************************
135
136 CALL ASTRO ( NYMD, NHMS, XLATS,XLONS, im*jm, TEMP1,RA )
137 NYMD2 = NYMD
138 NHMS2 = NHMS
139 CALL TICK ( NYMD2, NHMS2, NTMSTP )
140 CALL ASTRO ( NYMD2, NHMS2, XLATS,XLONS, im*jm, TEMP2,RA )
141
142 do j = 1,jm
143 do i = 1,im
144 zenith(I,j) = TEMP1(I,j) + TEMP2(I,j)
145 IF( zenith(I,j) .GT. 0.0 ) THEN
146 zenith(I,j) = (2./3.)*( zenith(I,j)-TEMP2(I,j)*TEMP1(I,j)
147 . / zenith(I,j) )
148 ENDIF
149 ENDDO
150 ENDDO
151
152 C **********************************************************************
153 c **** Compute Two-Dimension Total Cloud Fraction (0-1) ****
154 C **********************************************************************
155
156 c Initialize Clear Sky Probability for Low, Mid, and High Clouds
157 c --------------------------------------------------------------
158 do j =1,jm
159 do i =1,im
160 cldlow(i,j) = 0.0
161 cldmid(i,j) = 0.0
162 cldhi (i,j) = 0.0
163 enddo
164 enddo
165
166 c Adjust cloud fractions and cloud liquid water due to moist turbulence
167 c ---------------------------------------------------------------------
168 if(imstturb.ne.0) then
169 do L =1,lm
170 do j =1,jm
171 do i =1,im
172 cldtot(i,j,L)=min(1.0,max(cldsw(i,j,L),fccave(i,j,L)/imstturb))
173 cldmxo(i,j,L)=min(1.0,cswmo(i,j,L))
174 swlz(i,j,L)=swlz(i,j,L)+qliqave(i,j,L)/imstturb
175 enddo
176 enddo
177 enddo
178 else
179 do L =1,lm
180 do j =1,jm
181 do i =1,im
182 cldtot(i,j,L) = min( 1.0,cldsw(i,j,L) )
183 cldmxo(i,j,L) = min( 1.0,cswmo(i,j,L) )
184 enddo
185 enddo
186 enddo
187 endif
188
189 c Compute 3-D Cloud Fractions
190 c ---------------------------
191 if( nswcld.ne.0 ) then
192 do L = 1,lm
193 do j = 1,jm
194 do i = 1,im
195 c Compute Low-Mid-High Maximum Overlap Cloud Fractions
196 c ----------------------------------------------------
197 if( L.lt.mid_level ) then
198 cldhi (i,j) = max( cldtot(i,j,L),cldhi (i,j) )
199 elseif( L.ge.mid_level .and.
200 . L.lt.low_level ) then
201 cldmid(i,j) = max( cldtot(i,j,L),cldmid(i,j) )
202 elseif( L.ge.low_level ) then
203 cldlow(i,j) = max( cldtot(i,j,L),cldlow(i,j) )
204 endif
205
206 enddo
207 enddo
208 enddo
209 endif
210
211 c Totcld => Product of Clear Sky Probabilities for Low, Mid, and High Clouds
212 c --------------------------------------------------------------------------
213 do j = 1,jm
214 do i = 1,im
215 totcld(i,j) = 1.0 - (1.-cldhi (i,j))
216 . * (1.-cldmid(i,j))
217 . * (1.-cldlow(i,j))
218 enddo
219 enddo
220
221 c Compute Cloud Diagnostics
222 c -------------------------
223 if(icldfrc.gt.0) then
224 do j=1,jm
225 do i=1,im
226 qdiag(i,j,icldfrc,bi,bj) = qdiag(i,j,icldfrc,bi,bj) + totcld(i,j)
227 enddo
228 enddo
229 if ( (bi.eq.1) .and. (bj.eq.1) ) ncldfrc = ncldfrc + 1
230 endif
231
232 if( icldras.gt.0 ) then
233 do L=1,lm
234 do j=1,jm
235 do i=1,im
236 qdiag(i,j,icldras+L-1,bi,bj) = qdiag(i,j,icldras+L-1,bi,bj) +
237 . cswmo(i,j,L)
238 enddo
239 enddo
240 enddo
241 if ( (bi.eq.1) .and. (bj.eq.1) ) ncldras = ncldras + 1
242 endif
243
244 if( icldtot.gt.0 ) then
245 do L=1,lm
246 do j=1,jm
247 do i=1,im
248 qdiag(i,j,icldtot+L-1,bi,bj) = qdiag(i,j,icldtot+L-1,bi,bj) +
249 . cldtot(i,j,L)
250 enddo
251 enddo
252 enddo
253 if ( (bi.eq.1) .and. (bj.eq.1) ) ncldtot = ncldtot + 1
254 endif
255
256 if( icldlow.gt.0 ) then
257 do j=1,jm
258 do i=1,im
259 qdiag(i,j,icldlow,bi,bj) = qdiag(i,j,icldlow,bi,bj) + cldlow(i,j)
260 enddo
261 enddo
262 if ( (bi.eq.1) .and. (bj.eq.1) ) ncldlow = ncldlow + 1
263 endif
264
265 if( icldmid.gt.0 ) then
266 do j=1,jm
267 do i=1,im
268 qdiag(i,j,icldmid,bi,bj) = qdiag(i,j,icldmid,bi,bj) + cldmid(i,j)
269 enddo
270 enddo
271 if ( (bi.eq.1) .and. (bj.eq.1) ) ncldmid = ncldmid + 1
272 endif
273
274 if( icldhi.gt.0 ) then
275 do j=1,jm
276 do i=1,im
277 qdiag(i,j,icldhi,bi,bj) = qdiag(i,j,icldhi,bi,bj) + cldhi(i,j)
278 enddo
279 enddo
280 if ( (bi.eq.1) .and. (bj.eq.1) ) ncldhi = ncldhi + 1
281 endif
282
283 if( ilzrad.gt.0 ) then
284 do L=1,lm
285 do j=1,jm
286 do i=1,im
287 qdiag(i,j,ilzrad+L-1,bi,bj) = qdiag(i,j,ilzrad+L-1,bi,bj) +
288 . swlz(i,j,L)*1.0e6
289 enddo
290 enddo
291 enddo
292 if ( (bi.eq.1) .and. (bj.eq.1) ) nlzrad = nlzrad + 1
293 endif
294
295 c Albedo Diagnostics
296 c ------------------
297 if( ialbvisdr.gt.0 ) then
298 do j=1,jm
299 do i=1,im
300 qdiag(i,j,ialbvisdr,bi,bj) = qdiag(i,j,ialbvisdr,bi,bj) +
301 . albvisdr(i,j)
302 enddo
303 enddo
304 if ( (bi.eq.1) .and. (bj.eq.1) ) nalbvisdr = nalbvisdr + 1
305 endif
306
307 if( ialbvisdf.gt.0 ) then
308 do j=1,jm
309 do i=1,im
310 qdiag(i,j,ialbvisdf,bi,bj) = qdiag(i,j,ialbvisdf,bi,bj) +
311 . albvisdf(i,j)
312 enddo
313 enddo
314 if ( (bi.eq.1) .and. (bj.eq.1) ) nalbvisdf = nalbvisdf + 1
315 endif
316
317 if( ialbnirdr.gt.0 ) then
318 do j=1,jm
319 do i=1,im
320 qdiag(i,j,ialbnirdr,bi,bj) = qdiag(i,j,ialbnirdr,bi,bj) +
321 . albnirdr(i,j)
322 enddo
323 enddo
324 if ( (bi.eq.1) .and. (bj.eq.1) ) nalbnirdr = nalbnirdr + 1
325 endif
326
327 if( ialbnirdf.gt.0 ) then
328 do j=1,jm
329 do i=1,im
330 qdiag(i,j,ialbnirdf,bi,bj) = qdiag(i,j,ialbnirdf,bi,bj) +
331 . albnirdf(i,j)
332 enddo
333 enddo
334 if ( (bi.eq.1) .and. (bj.eq.1) ) nalbnirdf = nalbnirdf + 1
335 endif
336
337 C Compute Optical Thicknesses and Diagnostics
338 C -------------------------------------------
339 call opthk(tdry,plz,plze,swlz,cldtot,cldmxo,landtype,im,jm,lm,
340 . tautype)
341
342 do L = 1,lm
343 do j = 1,jm
344 do i = 1,im
345 tau(i,j,L) = tautype(i,j,L,1)+tautype(i,j,L,2)+tautype(i,j,L,3)
346 enddo
347 enddo
348 enddo
349
350 if( itauave.gt.0 ) then
351 do L=1,lm
352 do j=1,jm
353 do i=1,im
354 qdiag(i,j,itauave+L-1,bi,bj) = qdiag(i,j,itauave+L-1,bi,bj) +
355 . tau(i,j,L)*100/(plze(i,j,L+1)-plze(i,j,L))
356 enddo
357 enddo
358 enddo
359 if ( (bi.eq.1) .and. (bj.eq.1) ) ntauave = ntauave + 1
360 endif
361
362 if( itaucld.gt.0 ) then
363 do L=1,lm
364 do j=1,jm
365 do i=1,im
366 if( cldtot(i,j,L).ne.0.0 ) then
367 qdiag(i,j,itaucld +L-1,bi,bj) = qdiag(i,j,itaucld +L-1,bi,bj) +
368 . tau(i,j,L)*100/(plze(i,j,L+1)-plze(i,j,L))
369 qdiag(i,j,itaucldc+L-1,bi,bj) =
370 . qdiag(i,j,itaucldc+L-1,bi,bj) + 1.0
371 endif
372 enddo
373 enddo
374 enddo
375 endif
376
377 c Compute Low, Mid, and High Cloud Optical Depth Diagnostics
378 c ----------------------------------------------------------
379 if( itaulow.ne.0 ) then
380 do j = 1,jm
381 do i = 1,im
382 if( cldlow(i,j).ne.0.0 ) then
383 taulow(i,j) = 0.0
384 do L = low_level,lm
385 taulow(i,j) = taulow(i,j) + tau(i,j,L)
386 enddo
387 qdiag(i,j,itaulow,bi,bj ) = qdiag(i,j,itaulow,bi,bj ) +
388 . taulow(i,j)
389 qdiag(i,j,itaulowc,bi,bj) = qdiag(i,j,itaulowc,bi,bj) + 1.0
390 endif
391 enddo
392 enddo
393 endif
394
395 if( itaumid.ne.0 ) then
396 do j = 1,jm
397 do i = 1,im
398 if( cldmid(i,j).ne.0.0 ) then
399 taumid(i,j) = 0.0
400 do L = mid_level,low_level+1
401 taumid(i,j) = taumid(i,j) + tau(i,j,L)
402 enddo
403 qdiag(i,j,itaumid,bi,bj ) = qdiag(i,j,itaumid,bi,bj ) +
404 . taumid(i,j)
405 qdiag(i,j,itaumidc,bi,bj) = qdiag(i,j,itaumidc,bi,bj) + 1.0
406 endif
407 enddo
408 enddo
409 endif
410
411 if( itauhi.ne.0 ) then
412 do j = 1,jm
413 do i = 1,im
414 if( cldhi(i,j).ne.0.0 ) then
415 tauhi(i,j) = 0.0
416 do L = 1,mid_level+1
417 tauhi(i,j) = tauhi(i,j) + tau(i,j,L)
418 enddo
419 qdiag(i,j,itauhi,bi,bj ) = qdiag(i,j,itauhi,bi,bj ) +
420 . tauhi(i,j)
421 qdiag(i,j,itauhic,bi,bj) = qdiag(i,j,itauhic,bi,bj) + 1.0
422 endif
423 enddo
424 enddo
425 endif
426
427 C***********************************************************************
428 C **** LOOP OVER GLOBAL REGIONS ****
429 C **********************************************************************
430
431 do 1000 nn = 1,npcs
432
433 C **********************************************************************
434 C **** VARIABLE INITIALIZATION ****
435 C **********************************************************************
436
437 CALL STRIP ( zenith,COSZ,im*jm,ISTRIP,1,NN )
438
439 CALL STRIP ( plze, ple ,im*jm,ISTRIP,lm+1,NN)
440 CALL STRIP ( pkz , pk ,im*jm,ISTRIP,lm ,NN)
441 CALL STRIP ( dpres,dpstrip,im*jm,ISTRIP,lm ,NN)
442 CALL STRIP ( tdry, tzl ,im*jm,ISTRIP,lm ,NN)
443 CALL STRIP ( qz , qzl ,im*jm,ISTRIP,lm ,NN)
444 CALL STRIP ( oz , ozl ,im*jm,ISTRIP,lm ,NN)
445 CALL STRIP ( tau , taul ,im*jm,ISTRIP,lm ,NN)
446
447 CALL STRIP ( albvisdr,albuvdr,im*jm,ISTRIP,1,NN )
448 CALL STRIP ( albvisdf,albuvdf,im*jm,ISTRIP,1,NN )
449 CALL STRIP ( albnirdr,albirdr,im*jm,ISTRIP,1,NN )
450 CALL STRIP ( albnirdf,albirdf,im*jm,ISTRIP,1,NN )
451
452 call strip ( cldtot,clro,im*jm,istrip,lm,nn )
453
454 c Partition Tau between Water and Ice Particles
455 c ---------------------------------------------
456 do L= 1,lm
457 do i= 1,istrip
458 alf = min( max((tzl(i,l)-253.15)/20.,0.) ,1.)
459 taua(i,L) = 0.
460
461 if( alf.ne.0.0 .and. alf.ne.1.0 ) then
462 tauc(i,L,1) = taul(i,L)/(1.+alf/(1-alf) * (reffi/reffw*0.80) )
463 tauc(i,L,2) = taul(i,L)/(1.+(1-alf)/alf * (reffw/reffi*1.25) )
464 endif
465
466 if( alf.eq.0.0 ) then
467 tauc(i,L,1) = taul(i,L)
468 tauc(i,L,2) = 0.0
469 endif
470
471 if( alf.eq.1.0 ) then
472 tauc(i,L,1) = 0.0
473 tauc(i,L,2) = taul(i,L)
474 endif
475
476 reff(i,L,1) = reffi
477 reff(i,L,2) = reffw
478 enddo
479 enddo
480
481 call sorad ( istrip,1,1,lm,ple,tzl,qzl,ozl,co2,
482 . tauc,reff,clro,mid_level,low_level,taua,
483 . albirdr,albirdf,albuvdr,albuvdf,cosz,
484 . flux,fluxclr,fdirir,fdifir,dirpar,difpar,
485 . fdiruv,fdifuv )
486
487 C **********************************************************************
488 C **** Compute Mass-Weighted Theta Tendencies from Fluxes ****
489 C **********************************************************************
490
491 do l=1,lm
492 do i=1,istrip
493 alf = grav*(ple(i,Lm+1)-ptop)/(cp*dpstrip(i,L)*100)
494 dtsw (i,L) = alf*( flux (i,L)-flux (i,L+1) )/pk(i,L)
495 dtswc(i,L) = alf*( fluxclr(i,L)-fluxclr(i,L+1) )/pk(i,L)
496 enddo
497 enddo
498
499 call paste ( dtsw , dtradsw ,istrip,im*jm,lm,nn )
500 call paste ( dtswc, dtswclr ,istrip,im*jm,lm,nn )
501
502 call paste ( flux (1,1),osr ,istrip,im*jm,1,nn )
503 call paste ( fluxclr(1,1),osrclr,istrip,im*jm,1,nn )
504
505 call paste ( flux (1,lm+1),radswg,istrip,im*jm,1,nn )
506 call paste ( fluxclr(1,lm+1),swgclr,istrip,im*jm,1,nn )
507
508 call paste ( dirpar,fdirpar,istrip,im*jm,1,nn )
509 call paste ( difpar,fdifpar,istrip,im*jm,1,nn )
510
511 c Calculate Mean Albedo
512 c ---------------------
513 do i=1,istrip
514 if( cosz(i).gt.0.0 ) then
515 tstrip(i) = 1.0 - flux(i,lm+1)/
516 . ( fdirir(i)+fdifir(i)+dirpar(i)+difpar(i) + fdiruv(i)+fdifuv(i) )
517 if( tstrip(i).lt.0.0 ) tstrip(i) = undef
518 else
519 tstrip(i) = undef
520 endif
521 enddo
522 call paste ( tstrip,albedo,istrip,im*jm,1,nn )
523
524 1000 CONTINUE
525
526 c Mean Albedo Diagnostic
527 c ----------------------
528 if( ialbedo.gt.0 ) then
529 do j=1,jm
530 do i=1,im
531 if( albedo(i,j).ne.undef ) then
532 qdiag(i,j,ialbedo,bi,bj ) = qdiag(i,j,ialbedo,bi,bj )+albedo(i,j)
533 qdiag(i,j,ialbedoc,bi,bj) = qdiag(i,j,ialbedoc,bi,bj) + 1.0
534 endif
535 enddo
536 enddo
537 endif
538
539 C **********************************************************************
540 C **** ZERO-OUT OR BUMP COUNTERS ****
541 C **********************************************************************
542
543 nswlz = 0
544 nswcld = 0
545 imstturb = 0
546
547 do L = 1,lm
548 do j = 1,jm
549 do i = 1,im
550 fccave(i,j,L) = 0.
551 qliqave(i,j,L) = 0.
552 enddo
553 enddo
554 enddo
555
556 return
557 end
558 subroutine opthk ( tl,pl,ple,lz,cf,cfm,lwi,im,jm,lm,tau )
559 C***********************************************************************
560 C
561 C PURPOSE:
562 C ========
563 C Compute cloud optical thickness using temperature and
564 C cloud fractions.
565 C
566 C INPUT:
567 C ======
568 C tl ......... Temperature at Model Layers (K)
569 C pl ......... Model Layer Pressures (mb)
570 C ple ........ Model Edge Pressures (mb)
571 C lz ......... Diagnosed Convective and Large-Scale Cloud Water (g/g)
572 C cf ......... Total Cloud Fraction (Random + Maximum Overlap)
573 C cfm ........ Maximum Overlap Cloud Fraction
574 C lwi ........ Surface Land Type
575 C im ......... Longitudinal dimension
576 C jm ......... Latitudinal dimension
577 C lm ......... Vertical dimension
578 C
579 C OUTPUT:
580 C =======
581 C tau ........ Cloud Optical Thickness (non-dimensional)
582 C tau(im,jm,lm,1): Suspended Ice
583 C tau(im,jm,lm,2): Suspended Water
584 C tau(im,jm,lm,3): Raindrops
585 C
586 C***********************************************************************
587
588 implicit none
589
590 integer im,jm,lm,i,j,L
591
592 _RL tl(im,jm,lm)
593 _RL pl(im,jm,lm)
594 _RL ple(im,jm,lm+1)
595 _RL lz(im,jm,lm)
596 _RL cf(im,jm,lm)
597 _RL cfm(im,jm,lm)
598 _RL tau(im,jm,lm,3)
599 integer lwi(im,jm)
600
601 _RL dp, alf, fracls, fraccu
602 _RL tauice, tauh2o, tauras
603
604 c Compute Cloud Optical Depths
605 c ----------------------------
606 do L=1,lm
607 do j=1,jm
608 do i=1,im
609 alf = min( max(( tl(i,j,L)-233.15)/20.,0.) ,1.)
610 dp = ple(i,j,L+1)-ple(i,j,L)
611 tau(i,j,L,1) = 0.0
612 tau(i,j,L,2) = 0.0
613 tau(i,j,L,3) = 0.0
614 if( cf(i,j,L).ne.0.0 ) then
615
616 c Determine fraction of large-scale and cumulus clouds
617 c ----------------------------------------------------
618 fracls = ( cf(i,j,L)-cfm(i,j,L) )/cf(i,j,L)
619 fraccu = 1.0-fracls
620
621 c Define tau for large-scale ice and water clouds
622 c Note: tauice is scaled between (0.02 & .2) for: 0 < lz < 2 mg/kg over Land
623 c Note: tauh2o is scaled between (0.20 & 20) for: 0 < lz < 5 mg/kg over Land
624 c Note: tauh2o is scaled between (0.20 & 12) for: 0 < lz < 50 mg/kg over Ocean
625 c -------------------------------------------------------------------------------
626
627 c Large-Scale Ice
628 c ---------------
629 tauice = max( 0.0002, 0.002*min(500*lz(i,j,L)*1000,1.0) )
630 tau(i,j,L,1) = fracls*(1-alf)*tauice*dp
631
632 c Large-Scale Water
633 c -----------------
634 C Over Land
635 if( lwi(i,j).le.10 ) then
636 tauh2o = max( 0.0020, 0.200*min(200*lz(i,j,L)*1000,1.0) )
637 tau(i,j,L,3) = fracls*alf*tauh2o*dp
638 else
639 C Non-Precipitation Clouds Over Ocean
640 if( lz(i,j,L).eq.0.0 ) then
641 tauh2o = .12
642 tau(i,j,L,2) = fracls*alf*tauh2o*dp
643 else
644 C Over Ocean
645 tauh2o = max( 0.0020, 0.120*min( 20*lz(i,j,L)*1000,1.0) )
646 tau(i,j,L,3) = fracls*alf*tauh2o*dp
647 endif
648 endif
649
650 c Sub-Grid Convective
651 c -------------------
652 if( tl(i,j,L).gt.210.0 ) then
653 tauras = .16
654 tau(i,j,L,3) = tau(i,j,L,3) + fraccu*tauras*dp
655 else
656 tauras = tauice
657 tau(i,j,L,1) = tau(i,j,L,1) + fraccu*tauras*dp
658 endif
659
660 c Define tau for large-scale and cumulus clouds
661 c ---------------------------------------------
662 ccc tau(i,j,L) = ( fracls*( alf*tauh2o + (1-alf)*tauice )
663 ccc . + fraccu*tauras )*dp
664
665 endif
666
667 enddo
668 enddo
669 enddo
670
671 return
672 end
673 subroutine sorad(m,n,ndim,np,pl,ta,wa,oa,co2,
674 * taucld,reff,fcld,ict,icb,
675 * taual,rsirbm,rsirdf,rsuvbm,rsuvdf,cosz,
676 * flx,flc,fdirir,fdifir,fdirpar,fdifpar,
677 * fdiruv,fdifuv)
678 c********************************************************************
679 c
680 c This routine computes solar fluxes due to the absoption by water
681 c vapor, ozone, co2, o2, clouds, and aerosols and due to the
682 c scattering by clouds, aerosols, and gases.
683 c
684 c This is a vectorized code. It computes the fluxes simultaneous for
685 c (m x n) soundings, which is a subset of the (m x ndim) soundings.
686 c In a global climate model, m and ndim correspond to the numbers of
687 c grid boxes in the zonal and meridional directions, respectively.
688 c
689 c Ice and liquid cloud particles are allowed to co-exist in any of the
690 c np layers. Two sets of cloud parameters are required as inputs, one
691 c for ice paticles and the other for liquid particles. These parameters
692 c are optical thickness (taucld) and effective particle size (reff).
693 c
694 c If no information is available for reff, a default value of
695 c 10 micron for liquid water and 75 micron for ice can be used.
696 c
697 c Clouds are grouped into high, middle, and low clouds separated by the
698 c level indices ict and icb. For detail, see the subroutine cldscale.
699 c
700 c----- Input parameters:
701 c units size
702 c number of soundings in zonal direction (m) n/d 1
703 c number of soundings in meridional direction (n) n/d 1
704 c maximum number of soundings in n/d 1
705 c meridional direction (ndim)
706 c number of atmospheric layers (np) n/d 1
707 c level pressure (pl) mb m*ndim*(np+1)
708 c layer temperature (ta) k m*ndim*np
709 c layer specific humidity (wa) gm/gm m*ndim*np
710 c layer ozone concentration (oa) gm/gm m*ndim*np
711 c co2 mixing ratio by volumn (co2) parts/part 1
712 c cloud optical thickness (taucld) n/d m*ndim*np*2
713 c index 1 for ice particles
714 c index 2 for liquid drops
715 c effective cloud-particle size (reff) micrometer m*ndim*np*2
716 c index 1 for ice particles
717 c index 2 for liquid drops
718 c cloud amount (fcld) fraction m*ndim*np
719 c level index separating high and middle n/d 1
720 c clouds (ict)
721 c level index separating middle and low clouds n/d 1
722 c clouds (icb)
723 c aerosol optical thickness (taual) n/d m*ndim*np
724 c solar ir surface albedo for beam fraction m*ndim
725 c radiation (rsirbm)
726 c solar ir surface albedo for diffuse fraction m*ndim
727 c radiation (rsirdf)
728 c uv + par surface albedo for beam fraction m*ndim
729 c radiation (rsuvbm)
730 c uv + par surface albedo for diffuse fraction m*ndim
731 c radiation (rsuvdf)
732 c cosine of solar zenith angle (cosz) n/d m*ndim
733 c
734 c----- Output parameters
735 c
736 c all-sky flux (downward minus upward) (flx) fraction m*ndim*(np+1)
737 c clear-sky flux (downward minus upward) (flc) fraction m*ndim*(np+1)
738 c all-sky direct downward ir (0.7-10 micron)
739 c flux at the surface (fdirir) fraction m*ndim
740 c all-sky diffuse downward ir flux at
741 c the surface (fdifir) fraction m*ndim
742 c all-sky direct downward par (0.4-0.7 micron)
743 c flux at the surface (fdirpar) fraction m*ndim
744 c all-sky diffuse downward par flux at
745 c the surface (fdifpar) fraction m*ndim
746 *
747 c----- Notes:
748 c
749 c (1) The unit of flux is fraction of the incoming solar radiation
750 c at the top of the atmosphere. Therefore, fluxes should
751 c be equal to flux multiplied by the extra-terrestrial solar
752 c flux and the cosine of solar zenith angle.
753 c (2) Clouds and aerosols can be included in any layers by specifying
754 c fcld(i,j,k), taucld(i,j,k,*) and taual(i,j,k), k=1,np.
755 c For an atmosphere without clouds and aerosols,
756 c set fcld(i,j,k)=taucld(i,j,k,*)=taual(i,j,k)=0.0.
757 c (3) Aerosol single scattering albedos and asymmetry
758 c factors are specified in the subroutines solir and soluv.
759 c (4) pl(i,j,1) is the pressure at the top of the model, and
760 c pl(i,j,np+1) is the surface pressure.
761 c (5) the pressure levels ict and icb correspond approximately
762 c to 400 and 700 mb.
763 c
764 c**************************************************************************
765
766 implicit none
767
768 c-----Explicit Inline Directives
769
770 #ifdef CRAY
771 #ifdef f77
772 cfpp$ expand (expmn)
773 #endif
774 #endif
775 _RL expmn
776
777 c-----input parameters
778
779 integer m,n,ndim,np,ict,icb
780 _RL pl(m,ndim,np+1),ta(m,ndim,np),wa(m,ndim,np),oa(m,ndim,np)
781 _RL taucld(m,ndim,np,2),reff(m,ndim,np,2)
782 _RL fcld(m,ndim,np),taual(m,ndim,np)
783 _RL rsirbm(m,ndim),rsirdf(m,ndim),
784 * rsuvbm(m,ndim),rsuvdf(m,ndim),cosz(m,ndim),co2
785
786 c-----output parameters
787
788 _RL flx(m,ndim,np+1),flc(m,ndim,np+1)
789 _RL fdirir(m,ndim),fdifir(m,ndim)
790 _RL fdirpar(m,ndim),fdifpar(m,ndim)
791 _RL fdiruv(m,ndim),fdifuv(m,ndim)
792
793 c-----temporary array
794
795 integer i,j,k
796 _RL cc(m,n,3),tauclb(m,n,np),tauclf(m,n,np)
797 _RL dp(m,n,np),wh(m,n,np),oh(m,n,np),scal(m,n,np)
798 _RL swh(m,n,np+1),so2(m,n,np+1),df(m,n,np+1)
799 _RL sdf(m,n),sclr(m,n),csm(m,n),x
800
801 c-----------------------------------------------------------------
802
803 do j= 1, n
804 do i= 1, m
805
806 swh(i,j,1)=0.
807 so2(i,j,1)=0.
808
809 c-----csm is the effective secant of the solar zenith angle
810 c see equation (12) of Lacis and Hansen (1974, JAS)
811
812 csm(i,j)=35./sqrt(1224.*cosz(i,j)*cosz(i,j)+1.)
813
814 enddo
815 enddo
816
817
818 do k= 1, np
819 do j= 1, n
820 do i= 1, m
821
822 c-----compute layer thickness and pressure-scaling function.
823 c indices for the surface level and surface layer
824 c are np+1 and np, respectively.
825
826 dp(i,j,k)=pl(i,j,k+1)-pl(i,j,k)
827 scal(i,j,k)=dp(i,j,k)*(.5*(pl(i,j,k)+pl(i,j,k+1))/300.)**.8
828
829 c-----compute scaled water vapor amount, unit is g/cm**2
830
831 wh(i,j,k)=1.02*wa(i,j,k)*scal(i,j,k)*
832 * (1.+0.00135*(ta(i,j,k)-240.))
833 swh(i,j,k+1)=swh(i,j,k)+wh(i,j,k)
834
835 c-----compute ozone amount, unit is (cm-atm)stp.
836
837 oh(i,j,k)=1.02*oa(i,j,k)*dp(i,j,k)*466.7
838
839 enddo
840 enddo
841 enddo
842
843
844 c-----scale cloud optical thickness in each layer from taucld (with
845 c cloud amount fcld) to tauclb and tauclf (with cloud amount cc).
846 c tauclb is the scaled optical thickness for beam radiation and
847 c tauclf is for diffuse radiation.
848
849 call cldscale(m,n,ndim,np,cosz,fcld,taucld,ict,icb,
850 * cc,tauclb,tauclf)
851
852 c-----initialize fluxes for all-sky (flx), clear-sky (flc), and
853 c flux reduction (df)
854
855 do k=1, np+1
856 do j=1, n
857 do i=1, m
858 flx(i,j,k)=0.
859 flc(i,j,k)=0.
860 df(i,j,k)=0.
861 enddo
862 enddo
863 enddo
864
865 c-----compute solar ir fluxes
866
867 call solir (m,n,ndim,np,wh,taucld,tauclb,tauclf,reff,ict,icb
868 * ,fcld,cc,taual,csm,rsirbm,rsirdf,flx,flc,fdirir,fdifir)
869
870 c-----compute and update uv and par fluxes
871
872 call soluv (m,n,ndim,np,oh,dp,taucld,tauclb,tauclf,reff,ict,icb
873 * ,fcld,cc,taual,csm,rsuvbm,rsuvdf,flx,flc
874 * ,fdirpar,fdifpar,fdiruv,fdifuv)
875
876 c-----compute scaled amount of o2 (so2), unit is (cm-atm)stp.
877
878 do k= 1, np
879 do j= 1, n
880 do i= 1, m
881 so2(i,j,k+1)=so2(i,j,k)+165.22*scal(i,j,k)
882 enddo
883 enddo
884 enddo
885
886 c-----compute flux reduction due to oxygen following
887 c chou (J. climate, 1990). The fraction 0.0287 is the
888 c extraterrestrial solar flux in the o2 bands.
889
890 do k= 2, np+1
891 do j= 1, n
892 do i= 1, m
893 x=so2(i,j,k)*csm(i,j)
894 df(i,j,k)=df(i,j,k)+0.0287*(1.-expmn(-0.00027*sqrt(x)))
895 enddo
896 enddo
897 enddo
898
899 c-----compute scaled amounts for co2 (so2). unit is (cm-atm)stp.
900
901 do k= 1, np
902 do j= 1, n
903 do i= 1, m
904 so2(i,j,k+1)=so2(i,j,k)+co2*789.*scal(i,j,k)
905 enddo
906 enddo
907 enddo
908
909 c-----compute and update flux reduction due to co2 following
910 c chou (J. Climate, 1990)
911
912 call flxco2(m,n,np,so2,swh,csm,df)
913
914 c-----adjust for the effect of o2 cnd co2 on clear-sky fluxes.
915
916 do k= 2, np+1
917 do j= 1, n
918 do i= 1, m
919 flc(i,j,k)=flc(i,j,k)-df(i,j,k)
920 enddo
921 enddo
922 enddo
923
924 c-----adjust for the all-sky fluxes due to o2 and co2. It is
925 c assumed that o2 and co2 have no effects on solar radiation
926 c below clouds. clouds are assumed randomly overlapped.
927
928 do j=1,n
929 do i=1,m
930 sdf(i,j)=0.0
931 sclr(i,j)=1.0
932 enddo
933 enddo
934
935 do k=1,np
936 do j=1,n
937 do i=1,m
938
939 c-----sclr is the fraction of clear sky.
940 c sdf is the flux reduction below clouds.
941
942 if(fcld(i,j,k).gt.0.01) then
943 sdf(i,j)=sdf(i,j)+df(i,j,k)*sclr(i,j)*fcld(i,j,k)
944 sclr(i,j)=sclr(i,j)*(1.-fcld(i,j,k))
945 endif
946 flx(i,j,k+1)=flx(i,j,k+1)-sdf(i,j)-df(i,j,k+1)*sclr(i,j)
947
948 enddo
949 enddo
950 enddo
951
952 c-----adjust for the direct downward ir flux.
953 do j= 1, n
954 do i= 1, m
955 x = (1.-rsirdf(i,j))*fdifir(i,j) + (1.-rsirbm(i,j))*fdirir(i,j)
956 x = (-sdf(i,j)-df(i,j,np+1)*sclr(i,j))/x
957 fdirir(i,j)=fdirir(i,j)*(1.+x)
958 fdifir(i,j)=fdifir(i,j)*(1.+x)
959 enddo
960 enddo
961
962 return
963 end
964
965 c********************************************************************
966
967 subroutine cldscale (m,n,ndim,np,cosz,fcld,taucld,ict,icb,
968 * cc,tauclb,tauclf)
969
970 c********************************************************************
971 c
972 c This subroutine computes the covers of high, middle, and
973 c low clouds and scales the cloud optical thickness.
974 c
975 c To simplify calculations in a cloudy atmosphere, clouds are
976 c grouped into high, middle and low clouds separated by the levels
977 c ict and icb (level 1 is the top of the atmosphere).
978 c
979 c Within each of the three groups, clouds are assumed maximally
980 c overlapped, and the cloud cover (cc) of a group is the maximum
981 c cloud cover of all the layers in the group. The optical thickness
982 c (taucld) of a given layer is then scaled to new values (tauclb and
983 c tauclf) so that the layer reflectance corresponding to the cloud
984 c cover cc is the same as the original reflectance with optical
985 c thickness taucld and cloud cover fcld.
986 c
987 c---input parameters
988 c
989 c number of grid intervals in zonal direction (m)
990 c number of grid intervals in meridional direction (n)
991 c maximum number of grid intervals in meridional direction (ndim)
992 c number of atmospheric layers (np)
993 c cosine of the solar zenith angle (cosz)
994 c fractional cloud cover (fcld)
995 c cloud optical thickness (taucld)
996 c index separating high and middle clouds (ict)
997 c index separating middle and low clouds (icb)
998 c
999 c---output parameters
1000 c
1001 c fractional cover of high, middle, and low clouds (cc)
1002 c scaled cloud optical thickness for beam radiation (tauclb)
1003 c scaled cloud optical thickness for diffuse radiation (tauclf)
1004 c
1005 c********************************************************************
1006
1007 implicit none
1008
1009 c-----input parameters
1010
1011 integer m,n,ndim,np,ict,icb
1012 _RL cosz(m,ndim),fcld(m,ndim,np),taucld(m,ndim,np,2)
1013
1014 c-----output parameters
1015
1016 _RL cc(m,n,3),tauclb(m,n,np),tauclf(m,n,np)
1017
1018 c-----temporary variables
1019
1020 integer i,j,k,im,it,ia,kk
1021 _RL fm,ft,fa,xai,taux
1022
1023 c-----pre-computed table
1024
1025 integer nm,nt,na
1026 parameter (nm=11,nt=9,na=11)
1027 _RL dm,dt,da,t1,caib(nm,nt,na),caif(nt,na)
1028 parameter (dm=0.1,dt=0.30103,da=0.1,t1=-0.9031)
1029
1030 c-----include the pre-computed table for cai
1031
1032 #include "cai-dat.h"
1033 save caib,caif
1034
1035
1036 c-----clouds within each of the high, middle, and low clouds are
1037 c assumed maximally overlapped, and the cloud cover (cc)
1038 c for a group is the maximum cloud cover of all the layers
1039 c in the group
1040
1041 do j=1,n
1042 do i=1,m
1043 cc(i,j,1)=0.0
1044 cc(i,j,2)=0.0
1045 cc(i,j,3)=0.0
1046 enddo
1047 enddo
1048
1049 do k=1,ict-1
1050 do j=1,n
1051 do i=1,m
1052 cc(i,j,1)=max(cc(i,j,1),fcld(i,j,k))
1053 enddo
1054 enddo
1055 enddo
1056
1057 do k=ict,icb-1
1058 do j=1,n
1059 do i=1,m
1060 cc(i,j,2)=max(cc(i,j,2),fcld(i,j,k))
1061 enddo
1062 enddo
1063 enddo
1064
1065 do k=icb,np
1066 do j=1,n
1067 do i=1,m
1068 cc(i,j,3)=max(cc(i,j,3),fcld(i,j,k))
1069 enddo
1070 enddo
1071 enddo
1072
1073 c-----scale the cloud optical thickness.
1074 c taucld(i,j,k,1) is the optical thickness for ice particles, and
1075 c taucld(i,j,k,2) is the optical thickness for liquid particles.
1076
1077 do k=1,np
1078
1079 if(k.lt.ict) then
1080 kk=1
1081 elseif(k.ge.ict .and. k.lt.icb) then
1082 kk=2
1083 else
1084 kk=3
1085 endif
1086
1087 do j=1,n
1088 do i=1,m
1089
1090 tauclb(i,j,k) = 0.0
1091 tauclf(i,j,k) = 0.0
1092
1093 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1094 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1095
1096 c-----normalize cloud cover
1097
1098 fa=fcld(i,j,k)/cc(i,j,kk)
1099
1100 c-----table look-up
1101
1102 taux=min(taux,32.)
1103
1104 fm=cosz(i,j)/dm
1105 ft=(log10(taux)-t1)/dt
1106 fa=fa/da
1107
1108 im=int(fm+1.5)
1109 it=int(ft+1.5)
1110 ia=int(fa+1.5)
1111
1112 im=max(im,2)
1113 it=max(it,2)
1114 ia=max(ia,2)
1115
1116 im=min(im,nm-1)
1117 it=min(it,nt-1)
1118 ia=min(ia,na-1)
1119
1120 fm=fm-float(im-1)
1121 ft=ft-float(it-1)
1122 fa=fa-float(ia-1)
1123
1124 c-----scale cloud optical thickness for beam radiation.
1125 c the scaling factor, xai, is a function of the solar zenith
1126 c angle, optical thickness, and cloud cover.
1127
1128 xai= (-caib(im-1,it,ia)*(1.-fm)+
1129 * caib(im+1,it,ia)*(1.+fm))*fm*.5+caib(im,it,ia)*(1.-fm*fm)
1130
1131 xai=xai+(-caib(im,it-1,ia)*(1.-ft)+
1132 * caib(im,it+1,ia)*(1.+ft))*ft*.5+caib(im,it,ia)*(1.-ft*ft)
1133
1134 xai=xai+(-caib(im,it,ia-1)*(1.-fa)+
1135 * caib(im,it,ia+1)*(1.+fa))*fa*.5+caib(im,it,ia)*(1.-fa*fa)
1136
1137 xai= xai-2.*caib(im,it,ia)
1138 xai=max(xai,0.0)
1139
1140 tauclb(i,j,k) = taux*xai
1141
1142 c-----scale cloud optical thickness for diffuse radiation.
1143 c the scaling factor, xai, is a function of the cloud optical
1144 c thickness and cover but not the solar zenith angle.
1145
1146 xai= (-caif(it-1,ia)*(1.-ft)+
1147 * caif(it+1,ia)*(1.+ft))*ft*.5+caif(it,ia)*(1.-ft*ft)
1148
1149 xai=xai+(-caif(it,ia-1)*(1.-fa)+
1150 * caif(it,ia+1)*(1.+fa))*fa*.5+caif(it,ia)*(1.-fa*fa)
1151
1152 xai= xai-caif(it,ia)
1153 xai=max(xai,0.0)
1154
1155 tauclf(i,j,k) = taux*xai
1156
1157 endif
1158
1159 enddo
1160 enddo
1161 enddo
1162
1163 return
1164 end
1165 c***********************************************************************
1166
1167 subroutine solir (m,n,ndim,np,wh,taucld,tauclb,tauclf,reff,
1168 * ict,icb,fcld,cc,taual,csm,rsirbm,rsirdf,flx,flc,fdirir,fdifir)
1169
1170 c************************************************************************
1171 c compute solar flux in the infrared region. The spectrum is divided
1172 c into three bands:
1173 c
1174 c band wavenumber(/cm) wavelength (micron)
1175 c 1 1000-4400 2.27-10.0
1176 c 2 4400-8200 1.22-2.27
1177 c 3 8200-14300 0.70-1.22
1178 c
1179 c----- Input parameters: units size
1180 c
1181 c number of soundings in zonal direction (m) n/d 1
1182 c number of soundings in meridional direction (n) n/d 1
1183 c maximum number of soundings in n/d 1
1184 c meridional direction (ndim)
1185 c number of atmospheric layers (np) n/d 1
1186 c layer water vapor content (wh) gm/cm^2 m*n*np
1187 c cloud optical thickness (taucld) n/d m*ndim*np*2
1188 c index 1 for ice paticles
1189 c index 2 for liquid particles
1190 c scaled cloud optical thickness n/d m*n*np
1191 c for beam radiation (tauclb)
1192 c scaled cloud optical thickness n/d m*n*np
1193 c for diffuse radiation (tauclf)
1194 c effective cloud-particle size (reff) micrometer m*ndim*np*2
1195 c index 1 for ice paticles
1196 c index 2 for liquid particles
1197 c level index separating high and n/d m*n
1198 c middle clouds (ict)
1199 c level index separating middle and n/d m*n
1200 c low clouds (icb)
1201 c cloud amount (fcld) fraction m*ndim*np
1202 c cloud amount of high, middle, and n/d m*n*3
1203 c low clouds (cc)
1204 c aerosol optical thickness (taual) n/d m*ndim*np
1205 c cosecant of the solar zenith angle (csm) n/d m*n
1206 c near ir surface albedo for beam fraction m*ndim
1207 c radiation (rsirbm)
1208 c near ir surface albedo for diffuse fraction m*ndim
1209 c radiation (rsirdf)
1210 c
1211 c----- output (updated) parameters:
1212 c
1213 c all-sky flux (downward-upward) (flx) fraction m*ndim*(np+1)
1214 c clear-sky flux (downward-upward) (flc) fraction m*ndim*(np+1)
1215 c all-sky direct downward ir flux at
1216 c the surface (fdirir) fraction m*ndim
1217 c all-sky diffuse downward ir flux at
1218 c the surface (fdifir) fraction m*ndim
1219 c
1220 c----- note: the following parameters must be specified by users:
1221 c aerosol single scattering albedo (ssaal) n/d nband
1222 c aerosol asymmetry factor (asyal) n/d nband
1223 c
1224 c*************************************************************************
1225
1226 implicit none
1227
1228 c-----Explicit Inline Directives
1229
1230 #ifdef CRAY
1231 #ifdef f77
1232 cfpp$ expand (deledd)
1233 cfpp$ expand (sagpol)
1234 cfpp$ expand (expmn)
1235 #endif
1236 #endif
1237 _RL expmn
1238
1239 c-----input parameters
1240
1241 integer m,n,ndim,np,ict,icb
1242 _RL taucld(m,ndim,np,2),reff(m,ndim,np,2),fcld(m,ndim,np)
1243 _RL tauclb(m,n,np),tauclf(m,n,np),cc(m,n,3)
1244 _RL rsirbm(m,ndim),rsirdf(m,ndim)
1245 _RL wh(m,n,np),taual(m,ndim,np),csm(m,n)
1246
1247 c-----output (updated) parameters
1248
1249 _RL flx(m,ndim,np+1),flc(m,ndim,np+1)
1250 _RL fdirir(m,ndim),fdifir(m,ndim)
1251
1252 c-----static parameters
1253
1254 integer nk,nband
1255 parameter (nk=10,nband=3)
1256 _RL xk(nk),hk(nband,nk),ssaal(nband),asyal(nband)
1257 _RL aia(nband,3),awa(nband,3),aig(nband,3),awg(nband,3)
1258
1259 c-----temporary array
1260
1261 integer ib,ik,i,j,k
1262 _RL ssacl(m,n,np),asycl(m,n,np)
1263 _RL rr(m,n,np+1,2),tt(m,n,np+1,2),td(m,n,np+1,2),
1264 * rs(m,n,np+1,2),ts(m,n,np+1,2)
1265 _RL fall(m,n,np+1),fclr(m,n,np+1)
1266 _RL fsdir(m,n),fsdif(m,n)
1267
1268 _RL tauwv,tausto,ssatau,asysto,tauto,ssato,asyto
1269 _RL taux,reff1,reff2,w1,w2,g1,g2
1270 _RL ssaclt(m,n),asyclt(m,n)
1271 _RL rr1t(m,n),tt1t(m,n),td1t(m,n),rs1t(m,n),ts1t(m,n)
1272 _RL rr2t(m,n),tt2t(m,n),td2t(m,n),rs2t(m,n),ts2t(m,n)
1273
1274 c-----water vapor absorption coefficient for 10 k-intervals.
1275 c unit: cm^2/gm
1276
1277 data xk/
1278 1 0.0010, 0.0133, 0.0422, 0.1334, 0.4217,
1279 2 1.334, 5.623, 31.62, 177.8, 1000.0/
1280
1281 c-----water vapor k-distribution function,
1282 c the sum of hk is 0.52926. unit: fraction
1283
1284 data hk/
1285 1 .01074,.08236,.20673, .00360,.01157,.03497,
1286 2 .00411,.01133,.03011, .00421,.01143,.02260,
1287 3 .00389,.01240,.01336, .00326,.01258,.00696,
1288 4 .00499,.01381,.00441, .00465,.00650,.00115,
1289 5 .00245,.00244,.00026, .00145,.00094,.00000/
1290
1291 c-----aerosol single-scattering albedo and asymmetry factor
1292
1293 data ssaal,asyal/nband*0.999,nband*0.850/
1294
1295 c-----coefficients for computing the single scattering albedo of
1296 c ice clouds from ssa=1-(aia(*,1)+aia(*,2)*reff+aia(*,3)*reff**2)
1297
1298 data aia/
1299 1 .08938331, .00215346,-.00000260,
1300 2 .00299387, .00073709, .00000746,
1301 3 -.00001038,-.00000134, .00000000/
1302
1303 c-----coefficients for computing the single scattering albedo of
1304 c liquid clouds from ssa=1-(awa(*,1)+awa(*,2)*reff+awa(*,3)*reff**2)
1305
1306 data awa/
1307 1 .01209318,-.00019934, .00000007,
1308 2 .01784739, .00088757, .00000845,
1309 3 -.00036910,-.00000650,-.00000004/
1310
1311 c-----coefficients for computing the asymmetry factor of ice clouds
1312 c from asycl=aig(*,1)+aig(*,2)*reff+aig(*,3)*reff**2
1313
1314 data aig/
1315 1 .84090400, .76098937, .74935228,
1316 2 .00126222, .00141864, .00119715,
1317 3 -.00000385,-.00000396,-.00000367/
1318
1319 c-----coefficients for computing the asymmetry factor of liquid clouds
1320 c from asycl=awg(*,1)+awg(*,2)*reff+awg(*,3)*reff**2
1321
1322 data awg/
1323 1 .83530748, .74513197, .79375035,
1324 2 .00257181, .01370071, .00832441,
1325 3 .00005519,-.00038203,-.00023263/
1326
1327 c-----initialize surface fluxes, reflectances, and transmittances
1328
1329 do j= 1, n
1330 do i= 1, m
1331 fdirir(i,j)=0.0
1332 fdifir(i,j)=0.0
1333 rr(i,j,np+1,1)=rsirbm(i,j)
1334 rr(i,j,np+1,2)=rsirbm(i,j)
1335 rs(i,j,np+1,1)=rsirdf(i,j)
1336 rs(i,j,np+1,2)=rsirdf(i,j)
1337 td(i,j,np+1,1)=0.0
1338 td(i,j,np+1,2)=0.0
1339 tt(i,j,np+1,1)=0.0
1340 tt(i,j,np+1,2)=0.0
1341 ts(i,j,np+1,1)=0.0
1342 ts(i,j,np+1,2)=0.0
1343 enddo
1344 enddo
1345
1346 c-----integration over spectral bands
1347
1348 do 100 ib=1,nband
1349
1350 c-----compute cloud single scattering albedo and asymmetry factor
1351 c for a mixture of ice and liquid particles.
1352
1353 do k= 1, np
1354
1355 do j= 1, n
1356 do i= 1, m
1357
1358 ssaclt(i,j)=1.0
1359 asyclt(i,j)=1.0
1360
1361 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1362 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1363
1364 reff1=min(reff(i,j,k,1),130.)
1365 reff2=min(reff(i,j,k,2),20.0)
1366
1367 w1=(1.-(aia(ib,1)+(aia(ib,2)+
1368 * aia(ib,3)*reff1)*reff1))*taucld(i,j,k,1)
1369 w2=(1.-(awa(ib,1)+(awa(ib,2)+
1370 * awa(ib,3)*reff2)*reff2))*taucld(i,j,k,2)
1371 ssaclt(i,j)=(w1+w2)/taux
1372
1373 g1=(aig(ib,1)+(aig(ib,2)+aig(ib,3)*reff1)*reff1)*w1
1374 g2=(awg(ib,1)+(awg(ib,2)+awg(ib,3)*reff2)*reff2)*w2
1375 asyclt(i,j)=(g1+g2)/(w1+w2)
1376
1377 endif
1378
1379 enddo
1380 enddo
1381
1382 do j=1,n
1383 do i=1,m
1384 ssacl(i,j,k)=ssaclt(i,j)
1385 enddo
1386 enddo
1387 do j=1,n
1388 do i=1,m
1389 asycl(i,j,k)=asyclt(i,j)
1390 enddo
1391 enddo
1392
1393 enddo
1394
1395 c-----integration over the k-distribution function
1396
1397 do 200 ik=1,nk
1398
1399 do 300 k= 1, np
1400
1401 do j= 1, n
1402 do i= 1, m
1403
1404 tauwv=xk(ik)*wh(i,j,k)
1405
1406 c-----compute total optical thickness, single scattering albedo,
1407 c and asymmetry factor.
1408
1409 tausto=tauwv+taual(i,j,k)+1.0e-8
1410 ssatau=ssaal(ib)*taual(i,j,k)
1411 asysto=asyal(ib)*ssaal(ib)*taual(i,j,k)
1412
1413 c-----compute reflectance and transmittance for cloudless layers
1414
1415 tauto=tausto
1416 ssato=ssatau/tauto+1.0e-8
1417
1418 c if (ssato .gt. 0.001) then
1419
1420 c ssato=min(ssato,0.999999)
1421 c asyto=asysto/(ssato*tauto)
1422
1423 c call deledd(tauto,ssato,asyto,csm(i,j),
1424 c * rr1t(i,j),tt1t(i,j),td1t(i,j))
1425
1426 c call sagpol (tauto,ssato,asyto,rs1t(i,j),ts1t(i,j))
1427
1428 c else
1429
1430 td1t(i,j)=expmn(-tauto*csm(i,j))
1431 ts1t(i,j)=expmn(-1.66*tauto)
1432 tt1t(i,j)=0.0
1433 rr1t(i,j)=0.0
1434 rs1t(i,j)=0.0
1435
1436 c endif
1437
1438 c-----compute reflectance and transmittance for cloud layers
1439
1440 if (tauclb(i,j,k) .lt. 0.01) then
1441
1442 rr2t(i,j)=rr1t(i,j)
1443 tt2t(i,j)=tt1t(i,j)
1444 td2t(i,j)=td1t(i,j)
1445 rs2t(i,j)=rs1t(i,j)
1446 ts2t(i,j)=ts1t(i,j)
1447
1448 else
1449
1450 tauto=tausto+tauclb(i,j,k)
1451 ssato=(ssatau+ssacl(i,j,k)*tauclb(i,j,k))/tauto+1.0e-8
1452 ssato=min(ssato,0.999999)
1453 asyto=(asysto+asycl(i,j,k)*ssacl(i,j,k)*tauclb(i,j,k))/
1454 * (ssato*tauto)
1455
1456 call deledd(tauto,ssato,asyto,csm(i,j),
1457 * rr2t(i,j),tt2t(i,j),td2t(i,j))
1458
1459 tauto=tausto+tauclf(i,j,k)
1460 ssato=(ssatau+ssacl(i,j,k)*tauclf(i,j,k))/tauto+1.0e-8
1461 ssato=min(ssato,0.999999)
1462 asyto=(asysto+asycl(i,j,k)*ssacl(i,j,k)*tauclf(i,j,k))/
1463 * (ssato*tauto)
1464
1465 call sagpol (tauto,ssato,asyto,rs2t(i,j),ts2t(i,j))
1466
1467 endif
1468
1469 enddo
1470 enddo
1471
1472 c-----the following code appears to be awkward, but it is efficient
1473 c in certain parallel processors.
1474
1475 do j=1,n
1476 do i=1,m
1477 rr(i,j,k,1)=rr1t(i,j)
1478 enddo
1479 enddo
1480 do j=1,n
1481 do i=1,m
1482 tt(i,j,k,1)=tt1t(i,j)
1483 enddo
1484 enddo
1485 do j=1,n
1486 do i=1,m
1487 td(i,j,k,1)=td1t(i,j)
1488 enddo
1489 enddo
1490 do j=1,n
1491 do i=1,m
1492 rs(i,j,k,1)=rs1t(i,j)
1493 enddo
1494 enddo
1495 do j=1,n
1496 do i=1,m
1497 ts(i,j,k,1)=ts1t(i,j)
1498 enddo
1499 enddo
1500
1501 do j=1,n
1502 do i=1,m
1503 rr(i,j,k,2)=rr2t(i,j)
1504 enddo
1505 enddo
1506 do j=1,n
1507 do i=1,m
1508 tt(i,j,k,2)=tt2t(i,j)
1509 enddo
1510 enddo
1511 do j=1,n
1512 do i=1,m
1513 td(i,j,k,2)=td2t(i,j)
1514 enddo
1515 enddo
1516 do j=1,n
1517 do i=1,m
1518 rs(i,j,k,2)=rs2t(i,j)
1519 enddo
1520 enddo
1521 do j=1,n
1522 do i=1,m
1523 ts(i,j,k,2)=ts2t(i,j)
1524 enddo
1525 enddo
1526
1527 300 continue
1528
1529 c-----flux calculations
1530
1531 do k= 1, np+1
1532 do j= 1, n
1533 do i= 1, m
1534 fclr(i,j,k) = 0.
1535 fall(i,j,k) = 0.
1536 enddo
1537 enddo
1538 enddo
1539 do j= 1, n
1540 do i= 1, m
1541 fsdir(i,j) = 0.
1542 fsdif(i,j) = 0.
1543 enddo
1544 enddo
1545
1546 call cldflx (m,n,np,ict,icb,cc,rr,tt,td,rs,ts,
1547 * fclr,fall,fsdir,fsdif)
1548
1549 do k= 1, np+1
1550 do j= 1, n
1551 do i= 1, m
1552 flx(i,j,k) = flx(i,j,k)+fall(i,j,k)*hk(ib,ik)
1553 enddo
1554 enddo
1555 do j= 1, n
1556 do i= 1, m
1557 flc(i,j,k) = flc(i,j,k)+fclr(i,j,k)*hk(ib,ik)
1558 enddo
1559 enddo
1560 enddo
1561
1562 do j= 1, n
1563 do i= 1, m
1564 fdirir(i,j) = fdirir(i,j)+fsdir(i,j)*hk(ib,ik)
1565 fdifir(i,j) = fdifir(i,j)+fsdif(i,j)*hk(ib,ik)
1566 enddo
1567 enddo
1568
1569 200 continue
1570 100 continue
1571
1572 return
1573 end
1574
1575 c************************************************************************
1576
1577 subroutine soluv (m,n,ndim,np,oh,dp,taucld,tauclb,tauclf,reff,
1578 * ict,icb,fcld,cc,taual,csm,rsuvbm,rsuvdf,flx,flc
1579 * ,fdirpar,fdifpar,fdiruv,fdifuv)
1580
1581 c************************************************************************
1582 c compute solar fluxes in the uv+visible region. the spectrum is
1583 c grouped into 8 bands:
1584 c
1585 c Band Micrometer
1586 c
1587 c UV-C 1. .175 - .225
1588 c 2. .225 - .245
1589 c .260 - .280
1590 c 3. .245 - .260
1591 c
1592 c UV-B 4. .280 - .295
1593 c 5. .295 - .310
1594 c 6. .310 - .320
1595 c
1596 c UV-A 7. .320 - .400
1597 c
1598 c PAR 8. .400 - .700
1599 c
1600 c----- Input parameters: units size
1601 c
1602 c number of soundings in zonal direction (m) n/d 1
1603 c number of soundings in meridional direction (n) n/d 1
1604 c maximum number of soundings in n/d 1
1605 c meridional direction (ndim)
1606 c number of atmospheric layers (np) n/d 1
1607 c layer ozone content (oh) (cm-atm)stp m*n*np
1608 c layer pressure thickness (dp) mb m*n*np
1609 c cloud optical thickness (taucld) n/d m*ndim*np*2
1610 c index 1 for ice paticles
1611 c index 2 for liquid particles
1612 c scaled cloud optical thickness n/d m*n*np
1613 c for beam radiation (tauclb)
1614 c scaled cloud optical thickness n/d m*n*np
1615 c for diffuse radiation (tauclf)
1616 c effective cloud-particle size (reff) micrometer m*ndim*np*2
1617 c index 1 for ice paticles
1618 c index 2 for liquid particles
1619 c level indiex separating high and n/d m*n
1620 c middle clouds (ict)
1621 c level indiex separating middle and n/d m*n
1622 c low clouds (icb)
1623 c cloud amount (fcld) fraction m*ndim*np
1624 c cloud amounts of high, middle, and n/d m*n*3
1625 c low clouds (cc)
1626 c aerosol optical thickness (taual) n/d m*ndim*np
1627 c cosecant of the solar zenith angle (csm) n/d m*n
1628 c uv+par surface albedo for beam fraction m*ndim
1629 c radiation (rsuvbm)
1630 c uv+par surface albedo for diffuse fraction m*ndim
1631 c radiation (rsuvdf)
1632 c
1633 c----- output (updated) parameters:
1634 c
1635 c all-sky net downward flux (flx) fraction m*ndim*(np+1)
1636 c clear-sky net downward flux (flc) fraction m*ndim*(np+1)
1637 c all-sky direct downward par flux at
1638 c the surface (fdirpar) fraction m*ndim
1639 c all-sky diffuse downward par flux at
1640 c the surface (fdifpar) fraction m*ndim
1641 c
1642 c----- note: the following parameters must be specified by users:
1643 c
1644 c aerosol single scattering albedo (ssaal) n/d 1
1645 c aerosol asymmetry factor (asyal) n/d 1
1646 c
1647 *
1648 c***********************************************************************
1649
1650 implicit none
1651
1652 c-----Explicit Inline Directives
1653
1654 #ifdef CRAY
1655 #ifdef f77
1656 cfpp$ expand (deledd)
1657 cfpp$ expand (sagpol)
1658 #endif
1659 #endif
1660
1661 c-----input parameters
1662
1663 integer m,n,ndim,np,ict,icb
1664 _RL taucld(m,ndim,np,2),reff(m,ndim,np,2),fcld(m,ndim,np)
1665 _RL tauclb(m,n,np),tauclf(m,n,np),cc(m,n,3)
1666 _RL oh(m,n,np),dp(m,n,np),taual(m,ndim,np)
1667 _RL rsuvbm(m,ndim),rsuvdf(m,ndim),csm(m,n)
1668
1669 c-----output (updated) parameter
1670
1671 _RL flx(m,ndim,np+1),flc(m,ndim,np+1)
1672 _RL fdirpar(m,ndim),fdifpar(m,ndim)
1673 _RL fdiruv(m,ndim),fdifuv(m,ndim)
1674
1675 c-----static parameters
1676
1677 integer nband
1678 parameter (nband=8)
1679 _RL hk(nband),xk(nband),ry(nband)
1680 _RL asyal(nband),ssaal(nband),aig(3),awg(3)
1681
1682 c-----temporary array
1683
1684 integer i,j,k,ib
1685 _RL taurs,tauoz,tausto,ssatau,asysto,tauto,ssato,asyto
1686 _RL taux,reff1,reff2,g1,g2,asycl(m,n,np)
1687 _RL td(m,n,np+1,2),rr(m,n,np+1,2),tt(m,n,np+1,2),
1688 * rs(m,n,np+1,2),ts(m,n,np+1,2)
1689 _RL fall(m,n,np+1),fclr(m,n,np+1),fsdir(m,n),fsdif(m,n)
1690 _RL asyclt(m,n)
1691 _RL rr1t(m,n),tt1t(m,n),td1t(m,n),rs1t(m,n),ts1t(m,n)
1692 _RL rr2t(m,n),tt2t(m,n),td2t(m,n),rs2t(m,n),ts2t(m,n)
1693
1694 c-----hk is the fractional extra-terrestrial solar flux.
1695 c the sum of hk is 0.47074.
1696
1697 data hk/.00057, .00367, .00083, .00417,
1698 * .00600, .00556, .05913, .39081/
1699
1700 c-----xk is the ozone absorption coefficient. unit: /(cm-atm)stp
1701
1702 data xk /30.47, 187.2, 301.9, 42.83,
1703 * 7.09, 1.25, 0.0345, 0.0539/
1704
1705 c-----ry is the extinction coefficient for Rayleigh scattering.
1706 c unit: /mb.
1707
1708 data ry /.00604, .00170, .00222, .00132,
1709 * .00107, .00091, .00055, .00012/
1710
1711 c-----aerosol single-scattering albedo and asymmetry factor
1712
1713 data ssaal,asyal/nband*0.999,nband*0.850/
1714
1715 c-----coefficients for computing the asymmetry factor of ice clouds
1716 c from asycl=aig(*,1)+aig(*,2)*reff+aig(*,3)*reff**2
1717
1718 data aig/.74625000,.00105410,-.00000264/
1719
1720 c-----coefficients for computing the asymmetry factor of liquid
1721 c clouds from asycl=awg(*,1)+awg(*,2)*reff+awg(*,3)*reff**2
1722
1723 data awg/.82562000,.00529000,-.00014866/
1724
1725 c-----initialize surface reflectances and transmittances
1726
1727 do j= 1, n
1728 do i= 1, m
1729 rr(i,j,np+1,1)=rsuvbm(i,j)
1730 rr(i,j,np+1,2)=rsuvbm(i,j)
1731 rs(i,j,np+1,1)=rsuvdf(i,j)
1732 rs(i,j,np+1,2)=rsuvdf(i,j)
1733 td(i,j,np+1,1)=0.0
1734 td(i,j,np+1,2)=0.0
1735 tt(i,j,np+1,1)=0.0
1736 tt(i,j,np+1,2)=0.0
1737 ts(i,j,np+1,1)=0.0
1738 ts(i,j,np+1,2)=0.0
1739 enddo
1740 enddo
1741
1742 c-----compute cloud asymmetry factor for a mixture of
1743 c liquid and ice particles. unit of reff is micrometers.
1744
1745 do k= 1, np
1746
1747 do j= 1, n
1748 do i= 1, m
1749
1750 asyclt(i,j)=1.0
1751
1752 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1753 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1754
1755 reff1=min(reff(i,j,k,1),130.)
1756 reff2=min(reff(i,j,k,2),20.0)
1757
1758 g1=(aig(1)+(aig(2)+aig(3)*reff1)*reff1)*taucld(i,j,k,1)
1759 g2=(awg(1)+(awg(2)+awg(3)*reff2)*reff2)*taucld(i,j,k,2)
1760 asyclt(i,j)=(g1+g2)/taux
1761
1762 endif
1763
1764 enddo
1765 enddo
1766
1767 do j=1,n
1768 do i=1,m
1769 asycl(i,j,k)=asyclt(i,j)
1770 enddo
1771 enddo
1772
1773 enddo
1774
1775 do j=1,n
1776 do i=1,m
1777 fdiruv(i,j)=0.0
1778 fdifuv(i,j)=0.0
1779 enddo
1780 enddo
1781
1782 c-----integration over spectral bands
1783
1784 do 100 ib=1,nband
1785
1786 do 300 k= 1, np
1787
1788 do j= 1, n
1789 do i= 1, m
1790
1791 c-----compute ozone and rayleigh optical thicknesses
1792
1793 taurs=ry(ib)*dp(i,j,k)
1794 tauoz=xk(ib)*oh(i,j,k)
1795
1796 c-----compute total optical thickness, single scattering albedo,
1797 c and asymmetry factor
1798
1799 tausto=taurs+tauoz+taual(i,j,k)+1.0e-8
1800 ssatau=ssaal(ib)*taual(i,j,k)+taurs
1801 asysto=asyal(ib)*ssaal(ib)*taual(i,j,k)
1802
1803 c-----compute reflectance and transmittance for cloudless layers
1804
1805 tauto=tausto
1806 ssato=ssatau/tauto+1.0e-8
1807 ssato=min(ssato,0.999999)
1808 asyto=asysto/(ssato*tauto)
1809
1810 call deledd(tauto,ssato,asyto,csm(i,j),
1811 * rr1t(i,j),tt1t(i,j),td1t(i,j))
1812
1813 call sagpol (tauto,ssato,asyto,rs1t(i,j),ts1t(i,j))
1814
1815 c-----compute reflectance and transmittance for cloud layers
1816
1817 if (tauclb(i,j,k) .lt. 0.01) then
1818
1819 rr2t(i,j)=rr1t(i,j)
1820 tt2t(i,j)=tt1t(i,j)
1821 td2t(i,j)=td1t(i,j)
1822 rs2t(i,j)=rs1t(i,j)
1823 ts2t(i,j)=ts1t(i,j)
1824
1825 else
1826
1827 tauto=tausto+tauclb(i,j,k)
1828 ssato=(ssatau+tauclb(i,j,k))/tauto+1.0e-8
1829 ssato=min(ssato,0.999999)
1830 asyto=(asysto+asycl(i,j,k)*tauclb(i,j,k))/(ssato*tauto)
1831
1832 call deledd(tauto,ssato,asyto,csm(i,j),
1833 * rr2t(i,j),tt2t(i,j),td2t(i,j))
1834
1835 tauto=tausto+tauclf(i,j,k)
1836 ssato=(ssatau+tauclf(i,j,k))/tauto+1.0e-8
1837 ssato=min(ssato,0.999999)
1838 asyto=(asysto+asycl(i,j,k)*tauclf(i,j,k))/(ssato*tauto)
1839
1840 call sagpol (tauto,ssato,asyto,rs2t(i,j),ts2t(i,j))
1841
1842 endif
1843
1844 enddo
1845 enddo
1846
1847 do j=1,n
1848 do i=1,m
1849 rr(i,j,k,1)=rr1t(i,j)
1850 enddo
1851 enddo
1852 do j=1,n
1853 do i=1,m
1854 tt(i,j,k,1)=tt1t(i,j)
1855 enddo
1856 enddo
1857 do j=1,n
1858 do i=1,m
1859 td(i,j,k,1)=td1t(i,j)
1860 enddo
1861 enddo
1862 do j=1,n
1863 do i=1,m
1864 rs(i,j,k,1)=rs1t(i,j)
1865 enddo
1866 enddo
1867 do j=1,n
1868 do i=1,m
1869 ts(i,j,k,1)=ts1t(i,j)
1870 enddo
1871 enddo
1872
1873 do j=1,n
1874 do i=1,m
1875 rr(i,j,k,2)=rr2t(i,j)
1876 enddo
1877 enddo
1878 do j=1,n
1879 do i=1,m
1880 tt(i,j,k,2)=tt2t(i,j)
1881 enddo
1882 enddo
1883 do j=1,n
1884 do i=1,m
1885 td(i,j,k,2)=td2t(i,j)
1886 enddo
1887 enddo
1888 do j=1,n
1889 do i=1,m
1890 rs(i,j,k,2)=rs2t(i,j)
1891 enddo
1892 enddo
1893 do j=1,n
1894 do i=1,m
1895 ts(i,j,k,2)=ts2t(i,j)
1896 enddo
1897 enddo
1898
1899 300 continue
1900
1901 c-----flux calculations
1902
1903 do k= 1, np+1
1904 do j= 1, n
1905 do i= 1, m
1906 fclr(i,j,k) = 0.
1907 fall(i,j,k) = 0.
1908 enddo
1909 enddo
1910 enddo
1911 do j= 1, n
1912 do i= 1, m
1913 fsdir(i,j) = 0.
1914 fsdif(i,j) = 0.
1915 enddo
1916 enddo
1917 call cldflx (m,n,np,ict,icb,cc,rr,tt,td,rs,ts,
1918 * fclr,fall,fsdir,fsdif)
1919
1920 do k= 1, np+1
1921 do j= 1, n
1922 do i= 1, m
1923 flx(i,j,k)=flx(i,j,k)+fall(i,j,k)*hk(ib)
1924 enddo
1925 enddo
1926 do j= 1, n
1927 do i= 1, m
1928 flc(i,j,k)=flc(i,j,k)+fclr(i,j,k)*hk(ib)
1929 enddo
1930 enddo
1931 enddo
1932
1933 if(ib.eq.nband) then
1934 do j=1,n
1935 do i=1,m
1936 fdirpar(i,j)=fsdir(i,j)*hk(ib)
1937 fdifpar(i,j)=fsdif(i,j)*hk(ib)
1938 enddo
1939 enddo
1940 else
1941 do j=1,n
1942 do i=1,m
1943 fdiruv(i,j)=fdiruv(i,j)+fsdir(i,j)*hk(ib)
1944 fdifuv(i,j)=fdifuv(i,j)+fsdif(i,j)*hk(ib)
1945 enddo
1946 enddo
1947 endif
1948
1949 100 continue
1950
1951 return
1952 end
1953
1954 c*********************************************************************
1955
1956 subroutine deledd(tau,ssc,g0,csm,rr,tt,td)
1957
1958 c*********************************************************************
1959 c
1960 c-----uses the delta-eddington approximation to compute the
1961 c bulk scattering properties of a single layer
1962 c coded following King and Harshvardhan (JAS, 1986)
1963 c
1964 c inputs:
1965 c
1966 c tau: the effective optical thickness
1967 c ssc: the effective single scattering albedo
1968 c g0: the effective asymmetry factor
1969 c csm: the effective secant of the zenith angle
1970 c
1971 c outputs:
1972 c
1973 c rr: the layer reflection of the direct beam
1974 c tt: the layer diffuse transmission of the direct beam
1975 c td: the layer direct transmission of the direct beam
1976
1977 c*********************************************************************
1978
1979 implicit none
1980
1981 c-----Explicit Inline Directives
1982
1983 #ifdef CRAY
1984 #ifdef f77
1985 cfpp$ expand (expmn)
1986 #endif
1987 #endif
1988 _RL expmn
1989
1990 _RL zero,one,two,three,four,fourth,seven,tumin
1991 parameter (one=1., three=3.)
1992 parameter (seven=7., two=2.)
1993 parameter (four=4., fourth=.25)
1994 parameter (zero=0., tumin=1.e-20)
1995
1996 c-----input parameters
1997 _RL tau,ssc,g0,csm
1998
1999 c-----output parameters
2000 _RL rr,tt,td
2001
2002 c-----temporary parameters
2003
2004 _RL zth,ff,xx,taup,sscp,gp,gm1,gm2,gm3,akk,alf1,alf2,
2005 * all1,bll,st7,st8,cll,dll,fll,ell,st1,st2,st3,st4
2006 c
2007 zth = one / csm
2008
2009 c delta-eddington scaling of single scattering albedo,
2010 c optical thickness, and asymmetry factor,
2011 c K & H eqs(27-29)
2012
2013 ff = g0*g0
2014 xx = one-ff*ssc
2015 taup= tau*xx
2016 sscp= ssc*(one-ff)/xx
2017 gp = g0/(one+g0)
2018
2019 c extinction of the direct beam transmission
2020
2021 td = expmn(-taup*csm)
2022
2023 c gamma1, gamma2, gamma3 and gamma4. see table 2 and eq(26) K & H
2024 c ssc and gp are the d-s single scattering
2025 c albedo and asymmetry factor.
2026
2027 xx = three*gp
2028 gm1 = (seven - sscp*(four+xx))*fourth
2029 gm2 = -(one - sscp*(four-xx))*fourth
2030 gm3 = (two - zth*xx )*fourth
2031
2032 c akk is k as defined in eq(25) of K & H
2033
2034 xx = gm1-gm2
2035 akk = sqrt((gm1+gm2)*xx)
2036
2037 c alf1 and alf2 are alpha1 and alpha2 from eqs (23) & (24) of K & H
2038
2039 alf1 = gm1 - gm3 * xx
2040 alf2 = gm2 + gm3 * xx
2041
2042 c all1 is last term in eq(21) of K & H
2043 c bll is last term in eq(22) of K & H
2044
2045 xx = akk * two
2046 all1 = (gm3 - alf2 * zth )*xx*td
2047 bll = (one - gm3 + alf1*zth)*xx
2048
2049 xx = akk * zth
2050 st7 = one - xx
2051 st8 = one + xx
2052
2053 xx = akk * gm3
2054 cll = (alf2 + xx) * st7
2055 dll = (alf2 - xx) * st8
2056
2057 xx = akk * (one-gm3)
2058 fll = (alf1 + xx) * st8
2059 ell = (alf1 - xx) * st7
2060
2061 st3 = max(abs(st7*st8),tumin)
2062 st3 = sign (st3,st7)
2063
2064 st2 = expmn(-akk*taup)
2065 st4 = st2 * st2
2066
2067 st1 = sscp / ((akk+gm1 + (akk-gm1)*st4) * st3)
2068
2069 c rr is r-hat of eq(21) of K & H
2070 c tt is diffuse part of t-hat of eq(22) of K & H
2071
2072 rr = ( cll-dll*st4 -all1*st2)*st1
2073 tt = - ((fll-ell*st4)*td-bll*st2)*st1
2074
2075 rr = max(rr,zero)
2076 tt = max(tt,zero)
2077
2078 return
2079 end
2080
2081 c*********************************************************************
2082
2083 subroutine sagpol(tau,ssc,g0,rll,tll)
2084
2085 c*********************************************************************
2086 c-----transmittance (tll) and reflectance (rll) of diffuse radiation
2087 c follows Sagan and Pollock (JGR, 1967).
2088 c also, eq.(31) of Lacis and Hansen (JAS, 1974).
2089 c
2090 c-----input parameters:
2091 c
2092 c tau: the effective optical thickness
2093 c ssc: the effective single scattering albedo
2094 c g0: the effective asymmetry factor
2095 c
2096 c-----output parameters:
2097 c
2098 c rll: the layer reflection of diffuse radiation
2099 c tll: the layer transmission of diffuse radiation
2100 c
2101 c*********************************************************************
2102
2103 implicit none
2104
2105 c-----Explicit Inline Directives
2106
2107 #ifdef CRAY
2108 #ifdef f77
2109 cfpp$ expand (expmn)
2110 #endif
2111 #endif
2112 _RL expmn
2113
2114 _RL one,three,four
2115 parameter (one=1., three=3., four=4.)
2116
2117 c-----output parameters:
2118
2119 _RL tau,ssc,g0
2120
2121 c-----output parameters:
2122
2123 _RL rll,tll
2124
2125 c-----temporary arrays
2126
2127 _RL xx,uuu,ttt,emt,up1,um1,st1
2128
2129 xx = one-ssc*g0
2130 uuu = sqrt( xx/(one-ssc))
2131 ttt = sqrt( xx*(one-ssc)*three )*tau
2132 emt = expmn(-ttt)
2133 up1 = uuu + one
2134 um1 = uuu - one
2135 xx = um1*emt
2136 st1 = one / ((up1+xx) * (up1-xx))
2137 rll = up1*um1*(one-emt*emt)*st1
2138 tll = uuu*four*emt *st1
2139
2140 return
2141 end
2142
2143 c*******************************************************************
2144
2145 function expmn(fin)
2146
2147 c*******************************************************************
2148 c compute exponential for arguments in the range 0> fin > -10.
2149 c*******************************************************************
2150 implicit none
2151 _RL fin,expmn
2152
2153 _RL one,expmin,e1,e2,e3,e4
2154 parameter (one=1.0, expmin=-10.0)
2155 parameter (e1=1.0, e2=-2.507213e-1)
2156 parameter (e3=2.92732e-2, e4=-3.827800e-3)
2157
2158 if (fin .lt. expmin) fin = expmin
2159 expmn = ((e4*fin + e3)*fin+e2)*fin+e1
2160 expmn = expmn * expmn
2161 expmn = one / (expmn * expmn)
2162
2163 return
2164 end
2165
2166 c*******************************************************************
2167
2168 subroutine cldflx (m,n,np,ict,icb,cc,rr,tt,td,rs,ts,
2169 * fclr,fall,fsdir,fsdif)
2170
2171 c*******************************************************************
2172 c compute upward and downward fluxes using a two-stream adding method
2173 c following equations (3)-(5) of Chou (1992, JAS).
2174 c
2175 c clouds are grouped into high, middle, and low clouds which are
2176 c assumed randomly overlapped. It involves eight sets of calculations.
2177 c In each set of calculations, each atmospheric layer is homogeneous,
2178 c either with clouds or without clouds.
2179
2180 c input parameters:
2181 c m: number of soundings in zonal direction
2182 c n: number of soundings in meridional direction
2183 c np: number of atmospheric layers
2184 c ict: the level separating high and middle clouds
2185 c icb: the level separating middle and low clouds
2186 c cc: effective cloud covers for high, middle and low clouds
2187 c tt: diffuse transmission of a layer illuminated by beam radiation
2188 c td: direct beam tranmssion
2189 c ts: transmission of a layer illuminated by diffuse radiation
2190 c rr: reflection of a layer illuminated by beam radiation
2191 c rs: reflection of a layer illuminated by diffuse radiation
2192 c
2193 c output parameters:
2194 c fclr: clear-sky flux (downward minus upward)
2195 c fall: all-sky flux (downward minus upward)
2196 c fdndir: surface direct downward flux
2197 c fdndif: surface diffuse downward flux
2198 c*********************************************************************c
2199
2200 implicit none
2201
2202 c-----input parameters
2203
2204 integer m,n,np,ict,icb
2205
2206 _RL rr(m,n,np+1,2),tt(m,n,np+1,2),td(m,n,np+1,2)
2207 _RL rs(m,n,np+1,2),ts(m,n,np+1,2)
2208 _RL cc(m,n,3)
2209
2210 c-----temporary array
2211
2212 integer i,j,k,ih,im,is
2213 _RL rra(m,n,np+1,2,2),tta(m,n,np+1,2,2),tda(m,n,np+1,2,2)
2214 _RL rsa(m,n,np+1,2,2),rxa(m,n,np+1,2,2)
2215 _RL ch(m,n),cm(m,n),ct(m,n),flxdn(m,n,np+1)
2216 _RL fdndir(m,n),fdndif(m,n),fupdif
2217 _RL denm,xx
2218
2219 c-----output parameters
2220
2221 _RL fclr(m,n,np+1),fall(m,n,np+1)
2222 _RL fsdir(m,n),fsdif(m,n)
2223
2224 c-----initialize all-sky flux (fall) and surface downward fluxes
2225
2226 do k=1,np+1
2227 do j=1,n
2228 do i=1,m
2229 fall(i,j,k)=0.0
2230 enddo
2231 enddo
2232 enddo
2233
2234 do j=1,n
2235 do i=1,m
2236 fsdir(i,j)=0.0
2237 fsdif(i,j)=0.0
2238 enddo
2239 enddo
2240
2241 c-----compute transmittances and reflectances for a composite of
2242 c layers. layers are added one at a time, going down from the top.
2243 c tda is the composite transmittance illuminated by beam radiation
2244 c tta is the composite diffuse transmittance illuminated by
2245 c beam radiation
2246 c rsa is the composite reflectance illuminated from below
2247 c by diffuse radiation
2248 c tta and rsa are computed from eqs. (4b) and (3b) of Chou
2249
2250 c-----for high clouds. indices 1 and 2 denote clear and cloudy
2251 c situations, respectively.
2252
2253 do 10 ih=1,2
2254
2255 do j= 1, n
2256 do i= 1, m
2257 tda(i,j,1,ih,1)=td(i,j,1,ih)
2258 tta(i,j,1,ih,1)=tt(i,j,1,ih)
2259 rsa(i,j,1,ih,1)=rs(i,j,1,ih)
2260 tda(i,j,1,ih,2)=td(i,j,1,ih)
2261 tta(i,j,1,ih,2)=tt(i,j,1,ih)
2262 rsa(i,j,1,ih,2)=rs(i,j,1,ih)
2263 enddo
2264 enddo
2265
2266 do k= 2, ict-1
2267 do j= 1, n
2268 do i= 1, m
2269 denm = ts(i,j,k,ih)/( 1.-rsa(i,j,k-1,ih,1)*rs(i,j,k,ih))
2270 tda(i,j,k,ih,1)= tda(i,j,k-1,ih,1)*td(i,j,k,ih)
2271 tta(i,j,k,ih,1)= tda(i,j,k-1,ih,1)*tt(i,j,k,ih)
2272 * +(tda(i,j,k-1,ih,1)*rr(i,j,k,ih)
2273 * *rsa(i,j,k-1,ih,1)+tta(i,j,k-1,ih,1))*denm
2274 rsa(i,j,k,ih,1)= rs(i,j,k,ih)+ts(i,j,k,ih)
2275 * *rsa(i,j,k-1,ih,1)*denm
2276 tda(i,j,k,ih,2)= tda(i,j,k,ih,1)
2277 tta(i,j,k,ih,2)= tta(i,j,k,ih,1)
2278 rsa(i,j,k,ih,2)= rsa(i,j,k,ih,1)
2279 enddo
2280 enddo
2281 enddo
2282
2283 c-----for middle clouds
2284
2285 do 10 im=1,2
2286
2287 do k= ict, icb-1
2288 do j= 1, n
2289 do i= 1, m
2290 denm = ts(i,j,k,im)/( 1.-rsa(i,j,k-1,ih,im)*rs(i,j,k,im))
2291 tda(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*td(i,j,k,im)
2292 tta(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*tt(i,j,k,im)
2293 * +(tda(i,j,k-1,ih,im)*rr(i,j,k,im)
2294 * *rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2295 rsa(i,j,k,ih,im)= rs(i,j,k,im)+ts(i,j,k,im)
2296 * *rsa(i,j,k-1,ih,im)*denm
2297 enddo
2298 enddo
2299 enddo
2300
2301 10 continue
2302
2303 c-----layers are added one at a time, going up from the surface.
2304 c rra is the composite reflectance illuminated by beam radiation
2305 c rxa is the composite reflectance illuminated from above
2306 c by diffuse radiation
2307 c rra and rxa are computed from eqs. (4a) and (3a) of Chou
2308
2309 c-----for the low clouds
2310
2311 do 20 is=1,2
2312
2313 do j= 1, n
2314 do i= 1, m
2315 rra(i,j,np+1,1,is)=rr(i,j,np+1,is)
2316 rxa(i,j,np+1,1,is)=rs(i,j,np+1,is)
2317 rra(i,j,np+1,2,is)=rr(i,j,np+1,is)
2318 rxa(i,j,np+1,2,is)=rs(i,j,np+1,is)
2319 enddo
2320 enddo
2321
2322 do k=np,icb,-1
2323 do j= 1, n
2324 do i= 1, m
2325 denm=ts(i,j,k,is)/( 1.-rs(i,j,k,is)*rxa(i,j,k+1,1,is) )
2326 rra(i,j,k,1,is)=rr(i,j,k,is)+(td(i,j,k,is)
2327 * *rra(i,j,k+1,1,is)+tt(i,j,k,is)*rxa(i,j,k+1,1,is))*denm
2328 rxa(i,j,k,1,is)= rs(i,j,k,is)+ts(i,j,k,is)
2329 * *rxa(i,j,k+1,1,is)*denm
2330 rra(i,j,k,2,is)=rra(i,j,k,1,is)
2331 rxa(i,j,k,2,is)=rxa(i,j,k,1,is)
2332 enddo
2333 enddo
2334 enddo
2335
2336 c-----for middle clouds
2337
2338 do 20 im=1,2
2339
2340 do k= icb-1,ict,-1
2341 do j= 1, n
2342 do i= 1, m
2343 denm=ts(i,j,k,im)/( 1.-rs(i,j,k,im)*rxa(i,j,k+1,im,is) )
2344 rra(i,j,k,im,is)= rr(i,j,k,im)+(td(i,j,k,im)
2345 * *rra(i,j,k+1,im,is)+tt(i,j,k,im)*rxa(i,j,k+1,im,is))*denm
2346 rxa(i,j,k,im,is)= rs(i,j,k,im)+ts(i,j,k,im)
2347 * *rxa(i,j,k+1,im,is)*denm
2348 enddo
2349 enddo
2350 enddo
2351
2352 20 continue
2353
2354 c-----integration over eight sky situations.
2355 c ih, im, is denotes high, middle and low cloud groups.
2356
2357 do 100 ih=1,2
2358
2359 c-----clear portion
2360
2361 if(ih.eq.1) then
2362 do j=1,n
2363 do i=1,m
2364 ch(i,j)=1.0-cc(i,j,1)
2365 enddo
2366 enddo
2367
2368 else
2369
2370 c-----cloudy portion
2371
2372 do j=1,n
2373 do i=1,m
2374 ch(i,j)=cc(i,j,1)
2375 enddo
2376 enddo
2377
2378 endif
2379
2380 do 100 im=1,2
2381
2382 c-----clear portion
2383
2384 if(im.eq.1) then
2385
2386 do j=1,n
2387 do i=1,m
2388 cm(i,j)=ch(i,j)*(1.0-cc(i,j,2))
2389 enddo
2390 enddo
2391
2392 else
2393
2394 c-----cloudy portion
2395
2396 do j=1,n
2397 do i=1,m
2398 cm(i,j)=ch(i,j)*cc(i,j,2)
2399 enddo
2400 enddo
2401
2402 endif
2403
2404 do 100 is=1,2
2405
2406 c-----clear portion
2407
2408 if(is.eq.1) then
2409
2410 do j=1,n
2411 do i=1,m
2412 ct(i,j)=cm(i,j)*(1.0-cc(i,j,3))
2413 enddo
2414 enddo
2415
2416 else
2417
2418 c-----cloudy portion
2419
2420 do j=1,n
2421 do i=1,m
2422 ct(i,j)=cm(i,j)*cc(i,j,3)
2423 enddo
2424 enddo
2425
2426 endif
2427
2428 c-----add one layer at a time, going down.
2429
2430 do k= icb, np
2431 do j= 1, n
2432 do i= 1, m
2433 denm = ts(i,j,k,is)/( 1.-rsa(i,j,k-1,ih,im)*rs(i,j,k,is) )
2434 tda(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*td(i,j,k,is)
2435 tta(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*tt(i,j,k,is)
2436 * +(tda(i,j,k-1,ih,im)*rr(i,j,k,is)
2437 * *rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2438 rsa(i,j,k,ih,im)= rs(i,j,k,is)+ts(i,j,k,is)
2439 * *rsa(i,j,k-1,ih,im)*denm
2440 enddo
2441 enddo
2442 enddo
2443
2444 c-----add one layer at a time, going up.
2445
2446 do k= ict-1,1,-1
2447 do j= 1, n
2448 do i= 1, m
2449 denm =ts(i,j,k,ih)/(1.-rs(i,j,k,ih)*rxa(i,j,k+1,im,is))
2450 rra(i,j,k,im,is)= rr(i,j,k,ih)+(td(i,j,k,ih)
2451 * *rra(i,j,k+1,im,is)+tt(i,j,k,ih)*rxa(i,j,k+1,im,is))*denm
2452 rxa(i,j,k,im,is)= rs(i,j,k,ih)+ts(i,j,k,ih)
2453 * *rxa(i,j,k+1,im,is)*denm
2454 enddo
2455 enddo
2456 enddo
2457
2458 c-----compute fluxes following eq (5) of Chou (1992)
2459
2460 c fdndir is the direct downward flux
2461 c fdndif is the diffuse downward flux
2462 c fupdif is the diffuse upward flux
2463
2464 do k=2,np+1
2465 do j=1, n
2466 do i=1, m
2467 denm= 1./(1.- rxa(i,j,k,im,is)*rsa(i,j,k-1,ih,im))
2468 fdndir(i,j)= tda(i,j,k-1,ih,im)
2469 xx = tda(i,j,k-1,ih,im)*rra(i,j,k,im,is)
2470 fdndif(i,j)= (xx*rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2471 fupdif= (xx+tta(i,j,k-1,ih,im)*rxa(i,j,k,im,is))*denm
2472 flxdn(i,j,k)=fdndir(i,j)+fdndif(i,j)-fupdif
2473 enddo
2474 enddo
2475 enddo
2476
2477 do j=1, n
2478 do i=1, m
2479 flxdn(i,j,1)=1.0-rra(i,j,1,im,is)
2480 enddo
2481 enddo
2482
2483 c-----summation of fluxes over all (eight) sky situations.
2484
2485 do k=1,np+1
2486 do j=1,n
2487 do i=1,m
2488 if(ih.eq.1 .and. im.eq.1 .and. is.eq.1) then
2489 fclr(i,j,k)=flxdn(i,j,k)
2490 endif
2491 fall(i,j,k)=fall(i,j,k)+flxdn(i,j,k)*ct(i,j)
2492 enddo
2493 enddo
2494 enddo
2495
2496 do j=1,n
2497 do i=1,m
2498 fsdir(i,j)=fsdir(i,j)+fdndir(i,j)*ct(i,j)
2499 fsdif(i,j)=fsdif(i,j)+fdndif(i,j)*ct(i,j)
2500 enddo
2501 enddo
2502
2503 100 continue
2504
2505 return
2506 end
2507
2508 c*****************************************************************
2509
2510 subroutine flxco2(m,n,np,swc,swh,csm,df)
2511
2512 c*****************************************************************
2513
2514 c-----compute the reduction of clear-sky downward solar flux
2515 c due to co2 absorption.
2516
2517 implicit none
2518
2519 c-----input parameters
2520
2521 integer m,n,np
2522 _RL csm(m,n),swc(m,n,np+1),swh(m,n,np+1),cah(22,19)
2523
2524 c-----output (undated) parameter
2525
2526 _RL df(m,n,np+1)
2527
2528 c-----temporary array
2529
2530 integer i,j,k,ic,iw
2531 _RL xx,clog1,wlog,dc,dw,x1,x2,y2
2532
2533 c********************************************************************
2534 c-----include co2 look-up table
2535
2536 #include "cah-dat.h"
2537 c save cah
2538
2539 c********************************************************************
2540 c-----table look-up for the reduction of clear-sky solar
2541 c radiation due to co2. The fraction 0.0343 is the
2542 c extraterrestrial solar flux in the co2 bands.
2543
2544 do k= 2, np+1
2545 do j= 1, n
2546 do i= 1, m
2547 xx=1./.3
2548 clog1=log10(swc(i,j,k)*csm(i,j))
2549 wlog=log10(swh(i,j,k)*csm(i,j))
2550 ic=int( (clog1+3.15)*xx+1.)
2551 iw=int( (wlog+4.15)*xx+1.)
2552 if(ic.lt.2)ic=2
2553 if(iw.lt.2)iw=2
2554 if(ic.gt.22)ic=22
2555 if(iw.gt.19)iw=19
2556 dc=clog1-float(ic-2)*.3+3.
2557 dw=wlog-float(iw-2)*.3+4.
2558 x1=cah(1,iw-1)+(cah(1,iw)-cah(1,iw-1))*xx*dw
2559 x2=cah(ic-1,iw-1)+(cah(ic-1,iw)-cah(ic-1,iw-1))*xx*dw
2560 y2=x2+(cah(ic,iw-1)-cah(ic-1,iw-1))*xx*dc
2561 if (x1.lt.y2) x1=y2
2562 df(i,j,k)=df(i,j,k)+0.0343*(x1-y2)
2563 enddo
2564 enddo
2565 enddo
2566
2567 return
2568 end

  ViewVC Help
Powered by ViewVC 1.1.22