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

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

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


Revision 1.5 - (hide annotations) (download)
Wed Jul 7 19:33:48 2004 UTC (19 years, 10 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint54b_post
Changes since 1.4: +74 -129 lines
Almost last code to get rid of references to sigma coordinate

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

  ViewVC Help
Powered by ViewVC 1.1.22