/[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.1 - (hide annotations) (download)
Tue Jun 15 14:47:23 2004 UTC (20 years ago) by molod
Branch: MAIN
Add CVS header and name stuff (and rename)

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

  ViewVC Help
Powered by ViewVC 1.1.22