/[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.26 - (show annotations) (download)
Wed Apr 1 19:54:17 2009 UTC (15 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint62, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.25: +138 -130 lines
- perpetual spring equinox with "#define FIZHI_USE_FIXED_DAY"

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

  ViewVC Help
Powered by ViewVC 1.1.22