/[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.2 - (show annotations) (download)
Tue Jun 15 16:06:03 2004 UTC (20 years ago) by molod
Branch: MAIN
Changes since 1.1: +85 -96 lines
Include cpp options and diagnostics common if needed

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

  ViewVC Help
Powered by ViewVC 1.1.22