/[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.6 - (show annotations) (download)
Tue Jul 13 21:18:41 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.5: +2 -1 lines
debugging

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

  ViewVC Help
Powered by ViewVC 1.1.22