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

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

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


Revision 1.3 - (hide annotations) (download)
Thu Jun 24 19:57:02 2004 UTC (20 years ago) by molod
Branch: MAIN
CVS Tags: checkpoint54a_pre, checkpoint54a_post, checkpoint54, checkpoint53g_post, checkpoint53f_post
Changes since 1.2: +12 -10 lines
Add bi bj index to fizhi qdiag references

1 molod 1.3 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_lwrad.F,v 1.2 2004/06/15 16:06:03 molod Exp $
2 molod 1.1 C $Name: $
3    
4 molod 1.2 #include "CPP_OPTIONS.h"
5 molod 1.3 subroutine lwrio (nymd,nhms,bi,bj,istrip,npcs,
6 molod 1.1 . pz,tz,qz,plz,plze,pkz,pkht,oz,co2,
7     . cfc11,cfc12,cfc22,
8     . methane,n2o,emissivity,
9     . tgz,radlwg,st4,dst4,
10     . dtradlw,dlwdtg,dtradlwc,lwgclr,
11     . im,jm,lm,ptop,
12     . nlwcld,cldlw,clwmo,nlwlz,lwlz,
13 molod 1.2 . lpnt,imstturb,qliqave,fccave,landtype)
14 molod 1.1
15     implicit none
16 molod 1.2 #ifdef ALLOW_DIAGNOSTICS
17 molod 1.1 #include 'diagnostics.h'
18 molod 1.2 #endif
19 molod 1.1
20     c Input Variables
21     c ---------------
22 molod 1.3 integer nymd,nhms,istrip,npcs,bi,bj
23 molod 1.1
24     integer im,jm,lm
25     real ptop
26    
27     real pz(im,jm)
28     real tz(im,jm,lm)
29     real pkht(im,jm,lm)
30    
31     real co2,cfc11
32     real cfc12,cfc22
33     real methane (lm)
34     real n2o (lm)
35     real oz(im,jm,lm)
36     real qz(im,jm,lm)
37    
38     real radlwg(im,jm)
39     real lwgclr(im,jm)
40     real st4(im,jm)
41     real dst4(im,jm)
42     real dtradlw (im,jm,lm)
43     real dtradlwc(im,jm,lm)
44     real dlwdtg (im,jm,lm)
45    
46     integer nlwcld,nlwlz
47     real cldlw(im,jm,lm)
48     real clwmo(im,jm,lm)
49     real lwlz(im,jm,lm)
50    
51     real emissivity(im,jm,10)
52     real tgz(im,jm)
53     logical lpnt
54     integer imstturb
55     real qliqave(im,jm,lm)
56     real fccave(im,jm,lm)
57    
58     integer landtype(im,jm)
59    
60     c Local Variables
61     c ---------------
62     integer i,j,l,n,nn
63     integer mid_level,low_level
64    
65     real PLZ(im,jm,lm)
66     real PKZ(im,jm,lm)
67     real PLZE(im,jm,lm+1)
68     real cldtot (im,jm,lm)
69     real cldmxo (im,jm,lm)
70    
71     real pl(istrip,lm)
72     real pk(istrip,lm)
73     real pke(istrip,lm)
74     real ple(istrip,lm+1)
75    
76     real ADELPL(ISTRIP,lm)
77     real dtrad(istrip,lm) , dtradc(istrip,lm)
78     real OZL(ISTRIP,lm) , TZL(ISTRIP,lm)
79     real SHZL(ISTRIP,lm) , CLRO(ISTRIP,lm)
80     real CLMO(ISTRIP,lm)
81     real flx(ISTRIP,lm+1) , flxclr(ISTRIP,lm+1)
82     real cldlz(istrip,lm)
83     real dfdts(istrip,lm+1) , dtdtg(istrip,lm)
84    
85     real emiss(istrip,10)
86     real taual(istrip,lm,10)
87     real ssaal(istrip,lm,10)
88     real asyal(istrip,lm,10)
89     real cwc(istrip,lm,3)
90     real reff(istrip,lm,3)
91     real tauc(istrip,lm,3)
92    
93     real SGMT4(ISTRIP)
94     real TSURF(ISTRIP)
95     real dsgmt4(ISTRIP)
96     integer lwi(istrip)
97    
98     real getcon,secday,convrt,pcheck
99    
100     integer koz, kh2o
101     DATA KOZ /20/
102     data kh2o /18/
103     logical high, trace, cldwater
104     data high /.true./
105     data trace /.true./
106     data cldwater /.false./
107    
108    
109     C **********************************************************************
110     C **** INITIALIZATION ****
111     C **********************************************************************
112    
113     SECDAY = GETCON('SDAY')
114     CONVRT = GETCON('GRAVITY') / ( 100.0 * GETCON('CP') )
115    
116     c Determine Level Indices for Low-Mid-High Cloud Regions
117     c ------------------------------------------------------
118     low_level = lm
119     mid_level = lm
120     do L = lm-1,1,-1
121     pcheck = plz(
122    
123     if (pcheck.gt.700.0) low_level = L
124     if (pcheck.gt.400.0) mid_level = L
125     enddo
126    
127     c Adjust cloud fractions and cloud liquid water due to moist turbulence
128     c ---------------------------------------------------------------------
129     if(imstturb.ne.0) then
130     do L =1,lm
131     do j =1,jm
132     do i =1,im
133     cldtot(i,j,L)=min(1.0,max(cldlw(i,j,L),fccave(i,j,L)/imstturb))
134     cldmxo(i,j,L) = min( 1.0 , clwmo(i,j,L) )
135     lwlz(i,j,L) = lwlz(i,j,L) + qliqave(i,j,L)/imstturb
136     enddo
137     enddo
138     enddo
139     else
140     do L =1,lm
141     do j =1,jm
142     do i =1,im
143     cldtot(i,j,L) = min( 1.0,cldlw(i,j,L) )
144     cldmxo(i,j,L) = min( 1.0,clwmo(i,j,L) )
145     enddo
146     enddo
147     enddo
148     endif
149    
150     C***********************************************************************
151     C **** LOOP OVER STRIPS ****
152     C **********************************************************************
153    
154     DO 1000 NN=1,NPCS
155    
156     C **********************************************************************
157     C **** VARIABLE INITIALIZATION ****
158     C **********************************************************************
159    
160     CALL STRIP (OZ,OZL ,im*jm,ISTRIP,lm ,NN)
161     CALL STRIP (PLZE,ple ,im*jm,ISTRIP,lm+1 ,NN)
162     CALL STRIP (PLZ ,pl ,im*jm,ISTRIP,lm ,NN)
163     CALL STRIP (PKZ ,pk ,im*jm,ISTRIP,lm ,NN)
164     CALL STRIP (PKHT,pke ,im*jm,ISTRIP,lm ,NN)
165     CALL STRIP (TZ,TZL ,im*jm,ISTRIP,lm ,NN)
166     CALL STRIP (qz,SHZL ,im*jm,ISTRIP,lm ,NN)
167     CALL STRIP (cldtot,CLRO ,im*jm,ISTRIP,lm ,NN)
168     CALL STRIP (cldmxo,CLMO ,im*jm,ISTRIP,lm ,NN)
169     CALL STRIP (lwlz,cldlz ,im*jm,ISTRIP,lm ,NN)
170     CALL STRIP (tgz,tsurf ,im*jm,ISTRIP,1 ,NN)
171    
172     CALL STRIP (emissivity,emiss,im*jm,ISTRIP,10,NN)
173    
174     call stripitint (landtype,lwi,im*jm,im*jm,istrip,1,nn)
175    
176     DO I = 1,ISTRIP*lm
177     ADELPL(I,1) = convrt / ( ple(I,2)-ple(I,1) )
178     ENDDO
179    
180     C Compute Clouds
181     C --------------
182     IF(NLWCLD.NE.0)THEN
183     DO L = 1,lm
184     DO I = 1,ISTRIP
185     CLRO(I,L) = min( 1.0,clro(i,L) )
186     CLMO(I,L) = min( 1.0,clmo(i,L) )
187     ENDDO
188     ENDDO
189     ENDIF
190    
191     DO L = 1,lm
192     DO I = 1,ISTRIP
193     TZL(I,L) = TZL(I,L)*pk(I,L)
194     ENDDO
195     ENDDO
196     DO I = 1,ISTRIP
197     C To reduce longwave heating in lowest model layer
198     TZL(I,lm) = ( 2*tzl(i,lm)+tsurf(i) )/3.0
199     ENDDO
200    
201     C Compute Optical Thicknesses
202     C ---------------------------
203     call opthk ( tzl,pl,ple,cldlz,clro,clmo,lwi,istrip,1,lm,tauc )
204    
205     do n = 1,3
206     do L = 1,lm
207     do i = 1,istrip
208     tauc(i,L,n) = tauc(i,L,n)*0.75
209     enddo
210     enddo
211     enddo
212    
213     C Set the optical thickness, single scattering albedo and assymetry factor
214     C for aerosols equal to zero to omit the contribution of aerosols
215     c-------------------------------------------------------------------------
216     do n = 1,10
217     do L = 1,lm
218     do i = 1,istrip
219     taual(i,L,n) = 0.
220     ssaal(i,L,n) = 0.
221     asyal(i,L,n) = 0.
222     enddo
223     enddo
224     enddo
225    
226     C Set cwc and reff to zero - when cldwater is .false.
227     c----------------------------------------------------
228     do n = 1,3
229     do L = 1,lm
230     do i = 1,istrip
231     cwc(i,L,n) = 0.
232     reff(i,L,n) = 0.
233     enddo
234     enddo
235     enddo
236    
237     C **********************************************************************
238     C **** CALL M-D Chou RADIATION ****
239     C **********************************************************************
240    
241     call irrad ( istrip,1,1,lm,ple,tzl,shzl,ozl,tsurf,co2,
242     . n2o,methane,cfc11,cfc12,cfc22,emiss,
243     . cldwater,cwc,tauc,reff,clro,mid_level,low_level,
244     . taual,ssaal,asyal,
245     . high,trace,flx,flxclr,dfdts,sgmt4 )
246    
247     C **********************************************************************
248     C **** HEATING RATES ****
249     C **********************************************************************
250    
251     do L = 1,lm
252     do i = 1,istrip
253     dtrad(i,L) = ( flx(i,L)- flx(i,L+1))*adelpl(i,L)
254     dtdtg(i,L) = ( dfdts(i,L)- dfdts(i,L+1))*adelpl(i,L)
255     dtradc(i,L) = (flxclr(i,L)-flxclr(i,L+1))*adelpl(i,L)
256     enddo
257     enddo
258    
259     C **********************************************************************
260     C **** NET UPWARD LONGWAVE FLUX (W/M**2) ****
261     C **********************************************************************
262    
263     DO I = 1,ISTRIP
264     flx (i,1) = -flx (i,1)
265     flxclr(i,1) = -flxclr(i,1)
266     flx (i,lm+1) = -flx (i,lm+1)
267     flxclr(i,lm+1) = -flxclr(i,lm+1)
268     sgmt4(i) = - sgmt4(i)
269     dsgmt4(i) = - dfdts(i,lm+1)
270     ENDDO
271    
272     CALL PASTE ( flx (1,lm+1),RADLWG,ISTRIP,im*jm,1,NN )
273     CALL PASTE ( flxclr(1,lm+1),LWGCLR,ISTRIP,im*jm,1,NN )
274    
275     CALL PASTE ( sgmt4, st4,ISTRIP,im*jm,1,NN )
276     CALL PASTE ( dsgmt4,dst4,ISTRIP,im*jm,1,NN )
277    
278     C **********************************************************************
279     C **** PASTE AND BUMP SOME DIAGNOSTICS ****
280     C **********************************************************************
281    
282 molod 1.3 IF(IOLR.GT.0)CALL PSTBMP(flx(1,1),QDIAG(1,1,IOLR,bi,bj),ISTRIP,
283 molod 1.1 . im*jm, 1,NN)
284 molod 1.3 IF(IOLRCLR.GT.0)CALL PSTBMP(flxclr(1,1),QDIAG(1,1,IOLRCLR,bi,bj),
285     . ISTRIP,im*jm,1,NN)
286     IF(IOZLW.GT.0)CALL PSTBMP(OZL(1,1),QDIAG(1,1,IOZLW,bi,bj),ISTRIP,
287 molod 1.1 . im*jm,lm,NN)
288    
289     C **********************************************************************
290     C **** TENDENCY UPDATES ****
291     C **********************************************************************
292    
293     DO L = 1,lm
294     DO I = 1,ISTRIP
295     DTRAD (I,L) = ( ple(i,lm+1)-PTOP ) * DTRAD (I,L)/pk(I,L)
296     DTRADC(I,L) = ( ple(i,lm+1)-PTOP ) * DTRADC(I,L)/pk(I,L)
297     dtdtg(I,L) = ( ple(i,lm+1)-PTOP ) * dtdtg (I,L)/pk(I,L)
298     ENDDO
299     ENDDO
300     CALL PASTE ( DTRAD ,DTRADLW ,ISTRIP,im*jm,lm,NN )
301     CALL PASTE ( DTRADC,DTRADLWC,ISTRIP,im*jm,lm,NN )
302     CALL PASTE ( dtdtg ,dlwdtg ,ISTRIP,im*jm,lm,NN )
303    
304     1000 CONTINUE
305    
306     C **********************************************************************
307     C **** BUMP DIAGNOSTICS ****
308     C **********************************************************************
309    
310     if(itgrlw.ne.0) then
311     do j = 1,jm
312     do i = 1,im
313 molod 1.3 qdiag(i,j,itgrlw,bi,bj) = qdiag(i,j,itgrlw,bi,bj) + tgz(i,j)
314 molod 1.1 enddo
315     enddo
316     endif
317    
318     if (itlw.ne.0) then
319     do L = 1,lm
320     do j = 1,jm
321     do i = 1,im
322 molod 1.3 qdiag(i,j,itlw+L-1,bi,bj) = qdiag(i,j,itlw+L-1,bi,bj) +
323     . tz(i,j,L)*pkz(i,j,L)
324 molod 1.1 enddo
325     enddo
326     enddo
327     endif
328    
329     if (ishrad.ne.0) then
330     do L = 1,lm
331     do j = 1,jm
332     do i = 1,im
333 molod 1.3 qdiag(i,j,ishrad+L-1,bi,bj) = qdiag(i,j,ishrad+L-1,bi,bj) +
334     . qz(i,j,L)*1000
335 molod 1.1 enddo
336     enddo
337     enddo
338     endif
339    
340    
341     C **********************************************************************
342     C **** Increment Diagnostics Counters and Zero-Out Cloud Info ****
343     C **********************************************************************
344    
345     ntlw = ntlw + 1
346     nshrad = nshrad + 1
347     nozlw = nozlw + 1
348     ntgrlw = ntgrlw + 1
349     nolr = nolr + 1
350     nolrclr = nolrclr + 1
351    
352     nlwlz = 0
353     nlwcld = 0
354     imstturb = 0
355    
356     do L = 1,lm
357     do j = 1,jm
358     do i = 1,im
359     fccave(i,j,L) = 0.
360     qliqave(i,j,L) = 0.
361     enddo
362     enddo
363     enddo
364    
365     return
366     end
367     c********************** November 26, 1997 **************************
368    
369     subroutine irrad (m,n,ndim,np,pl,ta,wa,oa,ts,co2,
370     * n2o,ch4,cfc11,cfc12,cfc22,emiss,
371     * cldwater,cwc,taucl,reff,fcld,ict,icb,
372     * taual,ssaal,asyal,
373     * high,trace,flx,flc,dfdts,st4)
374    
375     c***********************************************************************
376     c
377     c Version IR-5 (September, 1998)
378     c
379     c New features of this version are:
380     c (1) The effect of aerosol scattering on LW fluxes is included.
381     c (2) Absorption and scattering of rain drops are included.
382     c
383     c***********************************************************************
384     c
385     c Version IR-4 (October, 1997)
386     c
387     c New features of this version are:
388     c (1) The surface is treated as non-black. The surface
389     c emissivity, emiss, is an input parameter
390     c (2) The effect of cloud scattering on LW fluxes is included
391     c
392     c*********************************************************************
393     c
394     c This routine computes ir fluxes due to water vapor, co2, o3,
395     c trace gases (n2o, ch4, cfc11, cfc12, cfc22, co2-minor),
396     c clouds, and aerosols.
397     c
398     c This is a vectorized code. It computes fluxes simultaneously for
399     c (m x n) soundings, which is a subset of (m x ndim) soundings.
400     c In a global climate model, m and ndim correspond to the numbers of
401     c grid boxes in the zonal and meridional directions, respectively.
402     c
403     c Some detailed descriptions of the radiation routine are given in
404     c Chou and Suarez (1994).
405     c
406     c Ice and liquid cloud particles are allowed to co-exist in any of the
407     c np layers.
408     c
409     c If no information is available for the effective cloud particle size,
410     c reff, default values of 10 micron for liquid water and 75 micron
411     c for ice can be used.
412     c
413     c The maximum-random assumption is applied for cloud overlapping.
414     c clouds are grouped into high, middle, and low clouds separated by the
415     c level indices ict and icb. Within each of the three groups, clouds
416     c are assumed maximally overlapped, and the cloud cover of a group is
417     c the maximum cloud cover of all the layers in the group. clouds among
418     c the three groups are assumed randomly overlapped. The indices ict and
419     c icb correpond approximately to the 400 mb and 700 mb levels.
420     c
421     c Aerosols are allowed to be in any of the np layers. Aerosol optical
422     c properties can be specified as functions of height and spectral band.
423     c
424     c There are options for computing fluxes:
425     c
426     c If cldwater=.true., taucl is computed from cwc and reff as a
427     c function of height and spectral band.
428     c If cldwater=.false., taucl must be given as input to the radiation
429     c routine. It is independent of spectral band.
430     c
431     c If high = .true., transmission functions in the co2, o3, and the
432     c three water vapor bands with strong absorption are computed using
433     c table look-up. cooling rates are computed accurately from the
434     c surface up to 0.01 mb.
435     c If high = .false., transmission functions are computed using the
436     c k-distribution method with linear pressure scaling for all spectral
437     c bands and gases. cooling rates are not calculated accurately for
438     c pressures less than 20 mb. the computation is faster with
439     c high=.false. than with high=.true.
440    
441     c If trace = .true., absorption due to n2o, ch4, cfcs, and the
442     c two minor co2 bands in the window region is included.
443     c If trace = .false., absorption in minor bands is neglected.
444     c
445     c The IR spectrum is divided into nine bands:
446     c
447     c band wavenumber (/cm) absorber
448     c
449     c 1 0 - 340 h2o
450     c 2 340 - 540 h2o
451     c 3 540 - 800 h2o,cont,co2,n2o
452     c 4 800 - 980 h2o,cont
453     c co2,f11,f12,f22
454     c 5 980 - 1100 h2o,cont,o3
455     c co2,f11
456     c 6 1100 - 1215 h2o,cont
457     c n2o,ch4,f12,f22
458     c 7 1215 - 1380 h2o
459     c n2o,ch4
460     c 8 1380 - 1900 h2o
461     c 9 1900 - 3000 h2o
462     c
463     c In addition, a narrow band in the 15 micrometer region is added to
464     c compute flux reduction due to n2o
465     c
466     c 10 540 - 620 h2o,cont,co2,n2o
467     c
468     c Band 3 (540-800/cm) is further divided into 3 sub-bands :
469     c
470     c subband wavenumber (/cm)
471     c
472     c 1 540 - 620
473     c 2 620 - 720
474     c 3 720 - 800
475     c
476     c---- Input parameters units size
477     c
478     c number of soundings in zonal direction (m) n/d 1
479     c number of soundings in meridional direction (n) n/d 1
480     c maximum number of soundings in
481     c meridional direction (ndim>=n) n/d 1
482     c number of atmospheric layers (np) n/d 1
483     c level pressure (pl) mb m*ndim*(np+1)
484     c layer temperature (ta) k m*ndim*np
485     c layer specific humidity (wa) g/g m*ndim*np
486     c layer ozone mixing ratio by mass (oa) g/g m*ndim*np
487     c surface temperature (ts) k m*ndim
488     c co2 mixing ratio by volumn (co2) pppv 1
489     c n2o mixing ratio by volumn (n2o) pppv np
490     c ch4 mixing ratio by volumn (ch4) pppv np
491     c cfc11 mixing ratio by volumn (cfc11) pppv 1
492     c cfc12 mixing ratio by volumn (cfc12) pppv 1
493     c cfc22 mixing ratio by volumn (cfc22) pppv 1
494     c surface emissivity (emiss) fraction m*ndim*10
495     c input option for cloud optical thickness n/d 1
496     c cldwater="true" if cwc is provided
497     c cldwater="false" if taucl is provided
498     c cloud water mixing ratio (cwc) gm/gm m*ndim*np*3
499     c index 1 for ice particles
500     c index 2 for liquid drops
501     c index 3 for rain drops
502     c cloud optical thickness (taucl) n/d m*ndim*np*3
503     c index 1 for ice particles
504     c index 2 for liquid drops
505     c index 3 for rain drops
506     c effective cloud-particle size (reff) micrometer m*ndim*np*3
507     c index 1 for ice particles
508     c index 2 for liquid drops
509     c index 3 for rain drops
510     c cloud amount (fcld) fraction m*ndim*np
511     c level index separating high and middle n/d 1
512     c clouds (ict)
513     c level index separating middle and low n/d 1
514     c clouds (icb)
515     c aerosol optical thickness (taual) n/d m*ndim*np*10
516     c aerosol single-scattering albedo (ssaal) n/d m*ndim*np*10
517     c aerosol asymmetry factor (asyal) n/d m*ndim*np*10
518     c high (see explanation above) 1
519     c trace (see explanation above) 1
520     c
521     c Data used in table look-up for transmittance calculations:
522     c
523     c c1 , c2, c3: for co2 (band 3)
524     c o1 , o2, o3: for o3 (band 5)
525     c h11,h12,h13: for h2o (band 1)
526     c h21,h22,h23: for h2o (band 2)
527     c h81,h82,h83: for h2o (band 8)
528     c
529     c---- output parameters
530     c
531     c net downward flux, all-sky (flx) w/m**2 m*ndim*(np+1)
532     c net downward flux, clear-sky (flc) w/m**2 m*ndim*(np+1)
533     c sensitivity of net downward flux
534     c to surface temperature (dfdts) w/m**2/k m*ndim*(np+1)
535     c emission by the surface (st4) w/m**2 m*ndim
536     c
537     c Notes:
538     c
539     c (1) Water vapor continuum absorption is included in 540-1380 /cm.
540     c (2) Scattering is parameterized for clouds and aerosols.
541     c (3) Diffuse cloud and aerosol transmissions are computed
542     c from exp(-1.66*tau).
543     c (4) If there are no clouds, flx=flc.
544     c (5) plevel(1) is the pressure at the top of the model atmosphere,
545     c and plevel(np+1) is the surface pressure.
546     c (6) Downward flux is positive and upward flux is negative.
547     c (7) dfdts is negative because upward flux is defined as negative.
548     c (8) For questions and coding errors, contact with Ming-Dah Chou,
549     c Code 913, NASA/Goddard Space Flight Center, Greenbelt, MD 20771.
550     c Phone: 301-286-4012, Fax: 301-286-1759,
551     c e-mail: chou@climate.gsfc.nasa.gov
552     c
553     c-----parameters defining the size of the pre-computed tables for transmittance
554     c calculations using table look-up.
555     c
556     c "nx" is the number of intervals in pressure
557     c "no" is the number of intervals in o3 amount
558     c "nc" is the number of intervals in co2 amount
559     c "nh" is the number of intervals in h2o amount
560     c "nt" is the number of copies to be made from the pre-computed
561     c transmittance tables to reduce "memory-bank conflict"
562     c in parallel machines and, hence, enhancing the speed of
563     c computations using table look-up.
564     c If such advantage does not exist, "nt" can be set to 1.
565     c***************************************************************************
566    
567     cfpp$ expand (h2oexps)
568     cfpp$ expand (conexps)
569     cfpp$ expand (co2exps)
570     cfpp$ expand (n2oexps)
571     cfpp$ expand (ch4exps)
572     cfpp$ expand (comexps)
573     cfpp$ expand (cfcexps)
574     cfpp$ expand (b10exps)
575     cfpp$ expand (tablup)
576     cfpp$ expand (h2okdis)
577     cfpp$ expand (co2kdis)
578     cfpp$ expand (n2okdis)
579     cfpp$ expand (ch4kdis)
580     cfpp$ expand (comkdis)
581     cfpp$ expand (cfckdis)
582     cfpp$ expand (b10kdis)
583    
584     implicit none
585     integer nx,no,nc,nh,nt
586     integer i,j,k,ip,iw,it,ib,ik,iq,isb,k1,k2
587     parameter (nx=26,no=21,nc=24,nh=31,nt=7)
588    
589     c---- input parameters ------
590    
591     integer m,n,ndim,np,ict,icb
592     real pl(m,ndim,np+1),ta(m,ndim,np),wa(m,ndim,np),oa(m,ndim,np),
593     * ts(m,ndim)
594     real co2,n2o(np),ch4(np),cfc11,cfc12,cfc22,emiss(m,ndim,10)
595     real cwc(m,ndim,np,3),taucl(m,ndim,np,3),reff(m,ndim,np,3),
596     * fcld(m,ndim,np)
597     real taual(m,ndim,np,10),ssaal(m,ndim,np,10),asyal(m,ndim,np,10)
598     logical cldwater,high,trace
599    
600     c---- output parameters ------
601    
602     real flx(m,ndim,np+1),flc(m,ndim,np+1),dfdts(m,ndim,np+1),
603     * st4(m,ndim)
604    
605     c---- static data -----
606    
607     real cb(5,10)
608     real xkw(9),aw(9),bw(9),pm(9),fkw(6,9),gkw(6,3),xke(9)
609     real aib(3,10),awb(4,10),aiw(4,10),aww(4,10),aig(4,10),awg(4,10)
610     integer ne(9),mw(9)
611    
612     c---- temporary arrays -----
613    
614     real pa(m,n,np),dt(m,n,np)
615     real sh2o(m,n,np+1),swpre(m,n,np+1),swtem(m,n,np+1)
616     real sco3(m,n,np+1),scopre(m,n,np+1),scotem(m,n,np+1)
617     real dh2o(m,n,np),dcont(m,n,np),dco2(m,n,np),do3(m,n,np)
618     real dn2o(m,n,np),dch4(m,n,np)
619     real df11(m,n,np),df12(m,n,np),df22(m,n,np)
620     real th2o(m,n,6),tcon(m,n,3),tco2(m,n,6,2)
621     real tn2o(m,n,4),tch4(m,n,4),tcom(m,n,2)
622     real tf11(m,n),tf12(m,n),tf22(m,n)
623     real h2oexp(m,n,np,6),conexp(m,n,np,3),co2exp(m,n,np,6,2)
624     real n2oexp(m,n,np,4),ch4exp(m,n,np,4),comexp(m,n,np,2)
625     real f11exp(m,n,np), f12exp(m,n,np), f22exp(m,n,np)
626     real clr(m,n,0:np+1),fclr(m,n)
627     real blayer(m,n,0:np+1),dlayer(m,n,np+1),dbs(m,n)
628     real clrlw(m,n),clrmd(m,n),clrhi(m,n)
629     real cwp(m,n,np,3)
630     real trant(m,n),tranal(m,n),transfc(m,n,np+1),trantcr(m,n,np+1)
631     real flxu(m,n,np+1),flxd(m,n,np+1),flcu(m,n,np+1),flcd(m,n,np+1)
632     real rflx(m,n,np+1),rflc(m,n,np+1)
633    
634     logical oznbnd,co2bnd,h2otbl,conbnd,n2obnd
635     logical ch4bnd,combnd,f11bnd,f12bnd,f22bnd,b10bnd
636    
637     real c1 (nx,nc,nt),c2 (nx,nc,nt),c3 (nx,nc,nt)
638     real o1 (nx,no,nt),o2 (nx,no,nt),o3 (nx,no,nt)
639     real h11(nx,nh,nt),h12(nx,nh,nt),h13(nx,nh,nt)
640     real h21(nx,nh,nt),h22(nx,nh,nt),h23(nx,nh,nt)
641     real h81(nx,nh,nt),h82(nx,nh,nt),h83(nx,nh,nt)
642    
643     real dp,xx,p1,dwe,dpe,a1,b1,fk1,a2,b2,fk2
644     real w1,w2,w3,g1,g2,g3,ww,gg,ff,taux,reff1,reff2
645    
646     c-----the following coefficients (equivalent to table 2 of
647     c chou and suarez, 1995) are for computing spectrally
648     c integrated planck fluxes using eq. (22)
649    
650     data cb/
651     1 -2.6844e-1,-8.8994e-2, 1.5676e-3,-2.9349e-6, 2.2233e-9,
652     2 3.7315e+1,-7.4758e-1, 4.6151e-3,-6.3260e-6, 3.5647e-9,
653     3 3.7187e+1,-3.9085e-1,-6.1072e-4, 1.4534e-5,-1.6863e-8,
654     4 -4.1928e+1, 1.0027e+0,-8.5789e-3, 2.9199e-5,-2.5654e-8,
655     5 -4.9163e+1, 9.8457e-1,-7.0968e-3, 2.0478e-5,-1.5514e-8,
656     6 -4.7107e+1, 8.9130e-1,-5.9735e-3, 1.5596e-5,-9.5911e-9,
657     7 -5.4041e+1, 9.5332e-1,-5.7784e-3, 1.2555e-5,-2.9377e-9,
658     8 -6.9233e+0,-1.5878e-1, 3.9160e-3,-2.4496e-5, 4.9301e-8,
659     9 1.1483e+2,-2.2376e+0, 1.6394e-2,-5.3672e-5, 6.6456e-8,
660     * 1.9668e+1,-3.1161e-1, 1.2886e-3, 3.6835e-7,-1.6212e-9/
661    
662     c-----xkw are the absorption coefficients for the first
663     c k-distribution function due to water vapor line absorption
664     c (tables 4 and 7). units are cm**2/g
665    
666     data xkw / 29.55 , 4.167e-1, 1.328e-2, 5.250e-4,
667     * 5.25e-4, 9.369e-3, 4.719e-2, 1.320e-0, 5.250e-4/
668    
669     c-----xke are the absorption coefficients for the first
670     c k-distribution function due to water vapor continuum absorption
671     c (table 6). units are cm**2/g
672    
673     data xke / 0.00, 0.00, 27.40, 15.8,
674     * 9.40, 7.75, 0.0, 0.0, 0.0/
675    
676     c-----mw are the ratios between neighboring absorption coefficients
677     c for water vapor line absorption (tables 4 and 7).
678    
679     data mw /6,6,8,6,6,8,9,6,16/
680    
681     c-----aw and bw (table 3) are the coefficients for temperature scaling
682     c in eq. (25).
683    
684     data aw/ 0.0021, 0.0140, 0.0167, 0.0302,
685     * 0.0307, 0.0195, 0.0152, 0.0008, 0.0096/
686     data bw/ -1.01e-5, 5.57e-5, 8.54e-5, 2.96e-4,
687     * 2.86e-4, 1.108e-4, 7.608e-5, -3.52e-6, 1.64e-5/
688    
689     data pm/ 1.0, 1.0, 1.0, 1.0, 1.0, 0.77, 0.5, 1.0, 1.0/
690    
691     c-----fkw is the planck-weighted k-distribution function due to h2o
692     c line absorption given in table 4 of Chou and Suarez (1995).
693     c the k-distribution function for the third band, fkw(*,3), is not used
694    
695     data fkw / 0.2747,0.2717,0.2752,0.1177,0.0352,0.0255,
696     2 0.1521,0.3974,0.1778,0.1826,0.0374,0.0527,
697     3 6*1.00,
698     4 0.4654,0.2991,0.1343,0.0646,0.0226,0.0140,
699     5 0.5543,0.2723,0.1131,0.0443,0.0160,0.0000,
700     6 0.5955,0.2693,0.0953,0.0335,0.0064,0.0000,
701     7 0.1958,0.3469,0.3147,0.1013,0.0365,0.0048,
702     8 0.0740,0.1636,0.4174,0.1783,0.1101,0.0566,
703     9 0.1437,0.2197,0.3185,0.2351,0.0647,0.0183/
704    
705     c-----gkw is the planck-weighted k-distribution function due to h2o
706     c line absorption in the 3 subbands (800-720,620-720,540-620 /cm)
707     c of band 3 given in table 7. Note that the order of the sub-bands
708     c is reversed.
709    
710     data gkw/ 0.1782,0.0593,0.0215,0.0068,0.0022,0.0000,
711     2 0.0923,0.1675,0.0923,0.0187,0.0178,0.0000,
712     3 0.0000,0.1083,0.1581,0.0455,0.0274,0.0041/
713    
714     c-----ne is the number of terms used in each band to compute water vapor
715     c continuum transmittance (table 6).
716    
717     data ne /0,0,3,1,1,1,0,0,0/
718     c
719     c-----coefficients for computing the extinction coefficient
720     c for cloud ice particles
721     c polynomial fit: y=a1+a2/x**a3; x is in m**2/gm
722     c
723     data aib / -0.44171, 0.62951, 0.06465,
724     2 -0.13727, 0.61291, 0.28962,
725     3 -0.01878, 1.67680, 0.79080,
726     4 -0.01896, 1.06510, 0.69493,
727     5 -0.04788, 0.88178, 0.54492,
728     6 -0.02265, 1.57390, 0.76161,
729     7 -0.01038, 2.15640, 0.89045,
730     8 -0.00450, 2.51370, 0.95989,
731     9 -0.00044, 3.15050, 1.03750,
732     * -0.02956, 1.44680, 0.71283/
733     c
734     c-----coefficients for computing the extinction coefficient
735     c for cloud liquid drops
736     c polynomial fit: y=a1+a2*x+a3*x**2+a4*x**3; x is in m**2/gm
737     c
738     data awb / 0.08641, 0.01769, -1.5572e-3, 3.4896e-5,
739     2 0.22027, 0.00997, -1.8719e-3, 5.3112e-5,
740     3 0.39252, -0.02817, 7.2931e-4, -3.8151e-6,
741     4 0.39975, -0.03426, 1.2884e-3, -1.7930e-5,
742     5 0.34021, -0.02805, 1.0654e-3, -1.5443e-5,
743     6 0.15587, 0.00371, -7.7705e-4, 2.0547e-5,
744     7 0.05518, 0.04544, -4.2067e-3, 1.0184e-4,
745     8 0.12724, 0.04751, -5.2037e-3, 1.3711e-4,
746     9 0.30390, 0.01656, -3.5271e-3, 1.0828e-4,
747     * 0.63617, -0.06287, 2.2350e-3, -2.3177e-5/
748     c
749     c-----coefficients for computing the single-scattering albedo
750     c for cloud ice particles
751     c polynomial fit: y=a1+a2*x+a3*x**2+a4*x**3; x is in m**2/gm
752     c
753     data aiw/ 0.17201, 1.2229e-2, -1.4837e-4, 5.8020e-7,
754     2 0.81470, -2.7293e-3, 9.7816e-8, 5.7650e-8,
755     3 0.54859, -4.8273e-4, 5.4353e-6, -1.5679e-8,
756     4 0.39218, 4.1717e-3, - 4.8869e-5, 1.9144e-7,
757     5 0.71773, -3.3640e-3, 1.9713e-5, -3.3189e-8,
758     6 0.77345, -5.5228e-3, 4.8379e-5, -1.5151e-7,
759     7 0.74975, -5.6604e-3, 5.6475e-5, -1.9664e-7,
760     8 0.69011, -4.5348e-3, 4.9322e-5, -1.8255e-7,
761     9 0.83963, -6.7253e-3, 6.1900e-5, -2.0862e-7,
762     * 0.64860, -2.8692e-3, 2.7656e-5, -8.9680e-8/
763     c
764     c-----coefficients for computing the single-scattering albedo
765     c for cloud liquid drops
766     c polynomial fit: y=a1+a2*x+a3*x**2+a4*x**3; x is in m**2/gm
767     c
768     data aww/ -7.8566e-2, 8.0875e-2, -4.3403e-3, 8.1341e-5,
769     2 -1.3384e-2, 9.3134e-2, -6.0491e-3, 1.3059e-4,
770     3 3.7096e-2, 7.3211e-2, -4.4211e-3, 9.2448e-5,
771     4 -3.7600e-3, 9.3344e-2, -5.6561e-3, 1.1387e-4,
772     5 0.40212, 7.8083e-2, -5.9583e-3, 1.2883e-4,
773     6 0.57928, 5.9094e-2, -5.4425e-3, 1.2725e-4,
774     7 0.68974, 4.2334e-2, -4.9469e-3, 1.2863e-4,
775     8 0.80122, 9.4578e-3, -2.8508e-3, 9.0078e-5,
776     9 1.02340, -2.6204e-2, 4.2552e-4, 3.2160e-6,
777     * 0.05092, 7.5409e-2, -4.7305e-3, 1.0121e-4/
778     c
779     c-----coefficients for computing the asymmetry factor for cloud ice particles
780     c polynomial fit: y=a1+a2*x+a3*x**2+a4*x**3; x is in m**2/gm
781     c
782     data aig / 0.57867, 1.0135e-2, -1.1142e-4, 4.1537e-7,
783     2 0.72259, 3.1149e-3, -1.9927e-5, 5.6024e-8,
784     3 0.76109, 4.5449e-3, -4.6199e-5, 1.6446e-7,
785     4 0.86934, 2.7474e-3, -3.1301e-5, 1.1959e-7,
786     5 0.89103, 1.8513e-3, -1.6551e-5, 5.5193e-8,
787     6 0.86325, 2.1408e-3, -1.6846e-5, 4.9473e-8,
788     7 0.85064, 2.5028e-3, -2.0812e-5, 6.3427e-8,
789     8 0.86945, 2.4615e-3, -2.3882e-5, 8.2431e-8,
790     9 0.80122, 3.1906e-3, -2.4856e-5, 7.2411e-8,
791     * 0.73290, 4.8034e-3, -4.4425e-5, 1.4839e-7/
792     c
793     c-----coefficients for computing the asymmetry factor for cloud liquid drops
794     c polynomial fit: y=a1+a2*x+a3*x**2+a4*x**3; x is in m**2/gm
795     c
796     data awg / -0.51930, 0.20290, -1.1747e-2, 2.3868e-4,
797     2 -0.22151, 0.19708, -1.2462e-2, 2.6646e-4,
798     3 0.14157, 0.14705, -9.5802e-3, 2.0819e-4,
799     4 0.41590, 0.10482, -6.9118e-3, 1.5115e-4,
800     5 0.55338, 7.7016e-2, -5.2218e-3, 1.1587e-4,
801     6 0.61384, 6.4402e-2, -4.6241e-3, 1.0746e-4,
802     7 0.67891, 4.8698e-2, -3.7021e-3, 9.1966e-5,
803     8 0.78169, 2.0803e-2, -1.4749e-3, 3.9362e-5,
804     9 0.93218, -3.3425e-2, 2.9632e-3, -6.9362e-5,
805     * 0.01649, 0.16561, -1.0723e-2, 2.3220e-4/
806     c
807     c-----include tables used in the table look-up for co2 (band 3),
808     c o3 (band 5), and h2o (bands 1, 2, and 7) transmission functions.
809    
810     logical first
811     data first /.true./
812    
813     include "h2o.tran3"
814     include "co2.tran3"
815     include "o3.tran3"
816    
817     save c1,c2,c3,o1,o2,o3
818     save h11,h12,h13,h21,h22,h23,h81,h82,h83
819    
820     if (first) then
821    
822     c-----tables co2 and h2o are only used with 'high' option
823    
824     if (high) then
825    
826     do iw=1,nh
827     do ip=1,nx
828     h11(ip,iw,1)=1.0-h11(ip,iw,1)
829     h21(ip,iw,1)=1.0-h21(ip,iw,1)
830     h81(ip,iw,1)=1.0-h81(ip,iw,1)
831     enddo
832     enddo
833    
834     do iw=1,nc
835     do ip=1,nx
836     c1(ip,iw,1)=1.0-c1(ip,iw,1)
837     enddo
838     enddo
839    
840     c-----copy tables to enhance the speed of co2 (band 3), o3 (band 5),
841     c and h2o (bands 1, 2, and 8 only) transmission calculations
842     c using table look-up.
843     c-----tables are replicated to avoid memory bank conflicts
844    
845     do it=2,nt
846     do iw=1,nc
847     do ip=1,nx
848     c1 (ip,iw,it)= c1(ip,iw,1)
849     c2 (ip,iw,it)= c2(ip,iw,1)
850     c3 (ip,iw,it)= c3(ip,iw,1)
851     enddo
852     enddo
853     do iw=1,nh
854     do ip=1,nx
855     h11(ip,iw,it)=h11(ip,iw,1)
856     h12(ip,iw,it)=h12(ip,iw,1)
857     h13(ip,iw,it)=h13(ip,iw,1)
858     h21(ip,iw,it)=h21(ip,iw,1)
859     h22(ip,iw,it)=h22(ip,iw,1)
860     h23(ip,iw,it)=h23(ip,iw,1)
861     h81(ip,iw,it)=h81(ip,iw,1)
862     h82(ip,iw,it)=h82(ip,iw,1)
863     h83(ip,iw,it)=h83(ip,iw,1)
864     enddo
865     enddo
866     enddo
867    
868     endif
869    
870     c-----always use table look-up for ozone transmittance
871    
872     do iw=1,no
873     do ip=1,nx
874     o1(ip,iw,1)=1.0-o1(ip,iw,1)
875     enddo
876     enddo
877    
878     do it=2,nt
879     do iw=1,no
880     do ip=1,nx
881     o1 (ip,iw,it)= o1(ip,iw,1)
882     o2 (ip,iw,it)= o2(ip,iw,1)
883     o3 (ip,iw,it)= o3(ip,iw,1)
884     enddo
885     enddo
886     enddo
887    
888     first=.false.
889    
890     endif
891    
892     c-----set the pressure at the top of the model atmosphere
893     c to 1.0e-4 if it is zero
894    
895     do j=1,n
896     do i=1,m
897     if (pl(i,j,1) .eq. 0.0) then
898     pl(i,j,1)=1.0e-4
899     endif
900     enddo
901     enddo
902    
903     c-----compute layer pressure (pa) and layer temperature minus 250K (dt)
904    
905     do k=1,np
906     do j=1,n
907     do i=1,m
908     pa(i,j,k)=0.5*(pl(i,j,k)+pl(i,j,k+1))
909     dt(i,j,k)=ta(i,j,k)-250.0
910     enddo
911     enddo
912     enddo
913    
914     c-----compute layer absorber amount
915    
916     c dh2o : water vapor amount (g/cm**2)
917     c dcont: scaled water vapor amount for continuum absorption
918     c (g/cm**2)
919     c dco2 : co2 amount (cm-atm)stp
920     c do3 : o3 amount (cm-atm)stp
921     c dn2o : n2o amount (cm-atm)stp
922     c dch4 : ch4 amount (cm-atm)stp
923     c df11 : cfc11 amount (cm-atm)stp
924     c df12 : cfc12 amount (cm-atm)stp
925     c df22 : cfc22 amount (cm-atm)stp
926     c the factor 1.02 is equal to 1000/980
927     c factors 789 and 476 are for unit conversion
928     c the factor 0.001618 is equal to 1.02/(.622*1013.25)
929     c the factor 6.081 is equal to 1800/296
930    
931    
932     do k=1,np
933     do j=1,n
934     do i=1,m
935    
936     dp = pl(i,j,k+1)-pl(i,j,k)
937     dh2o (i,j,k) = 1.02*wa(i,j,k)*dp+1.e-10
938     do3 (i,j,k) = 476.*oa(i,j,k)*dp+1.e-10
939     dco2 (i,j,k) = 789.*co2*dp+1.e-10
940     dch4 (i,j,k) = 789.*ch4(k)*dp+1.e-10
941     dn2o (i,j,k) = 789.*n2o(k)*dp+1.e-10
942     df11 (i,j,k) = 789.*cfc11*dp+1.e-10
943     df12 (i,j,k) = 789.*cfc12*dp+1.e-10
944     df22 (i,j,k) = 789.*cfc22*dp+1.e-10
945    
946     c-----compute scaled water vapor amount for h2o continuum absorption
947     c following eq. (43).
948    
949     xx=pa(i,j,k)*0.001618*wa(i,j,k)*wa(i,j,k)*dp
950     dcont(i,j,k) = xx*exp(1800./ta(i,j,k)-6.081)+1.e-10
951    
952     enddo
953     enddo
954     enddo
955    
956     c-----compute column-integrated h2o amoumt, h2o-weighted pressure
957     c and temperature. it follows eqs. (37) and (38).
958    
959     if (high) then
960     call column(m,n,np,pa,dt,dh2o,sh2o,swpre,swtem)
961     endif
962    
963     c-----compute layer cloud water amount (gm/m**2)
964     c index is 1 for ice, 2 for waterdrops and 3 for raindrops.
965     c Rain optical thickness is 0.00307 /(gm/m**2).
966     c It is for a specific drop size distribution provided by Q. Fu.
967    
968     if (cldwater) then
969     do k=1,np
970     do j=1,n
971     do i=1,m
972     dp =pl(i,j,k+1)-pl(i,j,k)
973     cwp(i,j,k,1)=1.02*10000.0*cwc(i,j,k,1)*dp
974     cwp(i,j,k,2)=1.02*10000.0*cwc(i,j,k,2)*dp
975     cwp(i,j,k,3)=1.02*10000.0*cwc(i,j,k,3)*dp
976     taucl(i,j,k,3)=0.00307*cwp(i,j,k,3)
977     enddo
978     enddo
979     enddo
980     endif
981    
982     c-----the surface (np+1) is treated as a layer filled with black clouds.
983     c clr is the equivalent clear fraction of a layer
984     c transfc is the transmttance between the surface and a pressure level.
985     c trantcr is the clear-sky transmttance between the surface and a
986     c pressure level.
987    
988     do j=1,n
989     do i=1,m
990     clr(i,j,0) = 1.0
991     clr(i,j,np+1) = 0.0
992     st4(i,j) = 0.0
993     transfc(i,j,np+1)=1.0
994     trantcr(i,j,np+1)=1.0
995     enddo
996     enddo
997    
998     c-----initialize fluxes
999    
1000     do k=1,np+1
1001     do j=1,n
1002     do i=1,m
1003     flx(i,j,k) = 0.0
1004     flc(i,j,k) = 0.0
1005     dfdts(i,j,k)= 0.0
1006     rflx(i,j,k) = 0.0
1007     rflc(i,j,k) = 0.0
1008     enddo
1009     enddo
1010     enddo
1011    
1012     c-----integration over spectral bands
1013    
1014     do 1000 ib=1,10
1015    
1016     c-----if h2otbl, compute h2o (line) transmittance using table look-up.
1017     c if conbnd, compute h2o (continuum) transmittance in bands 3-6.
1018     c if co2bnd, compute co2 transmittance in band 3.
1019     c if oznbnd, compute o3 transmittance in band 5.
1020     c if n2obnd, compute n2o transmittance in bands 6 and 7.
1021     c if ch4bnd, compute ch4 transmittance in bands 6 and 7.
1022     c if combnd, compute co2-minor transmittance in bands 4 and 5.
1023     c if f11bnd, compute cfc11 transmittance in bands 4 and 5.
1024     c if f12bnd, compute cfc12 transmittance in bands 4 and 6.
1025     c if f22bnd, compute cfc22 transmittance in bands 4 and 6.
1026     c if b10bnd, compute flux reduction due to n2o in band 10.
1027    
1028     h2otbl=high.and.(ib.eq.1.or.ib.eq.2.or.ib.eq.8)
1029     conbnd=ib.ge.3.and.ib.le.6
1030     co2bnd=ib.eq.3
1031     oznbnd=ib.eq.5
1032     n2obnd=ib.eq.6.or.ib.eq.7
1033     ch4bnd=ib.eq.6.or.ib.eq.7
1034     combnd=ib.eq.4.or.ib.eq.5
1035     f11bnd=ib.eq.4.or.ib.eq.5
1036     f12bnd=ib.eq.4.or.ib.eq.6
1037     f22bnd=ib.eq.4.or.ib.eq.6
1038     b10bnd=ib.eq.10
1039    
1040     c-----blayer is the spectrally integrated planck flux of the mean layer
1041     c temperature derived from eq. (22)
1042     c the fitting for the planck flux is valid in the range 160-345 K.
1043    
1044     do k=1,np
1045     do j=1,n
1046     do i=1,m
1047     blayer(i,j,k)=ta(i,j,k)*(ta(i,j,k)*(ta(i,j,k)
1048     * *(ta(i,j,k)*cb(5,ib)+cb(4,ib))+cb(3,ib))
1049     * +cb(2,ib))+cb(1,ib)
1050     enddo
1051     enddo
1052     enddo
1053    
1054     c-----the earth's surface, with an index "np+1", is treated as a layer
1055    
1056     do j=1,n
1057     do i=1,m
1058     blayer(i,j,0) = 0.0
1059     blayer(i,j,np+1)= ( ts(i,j)*(ts(i,j)*(ts(i,j)
1060     * *(ts(i,j)*cb(5,ib)+cb(4,ib))+cb(3,ib))
1061     * +cb(2,ib))+cb(1,ib) )*emiss(i,j,ib)
1062    
1063     c-----dbs is the derivative of the surface emission with respect to
1064     c surface temperature (eq. 59).
1065    
1066     dbs(i,j)=(ts(i,j)*(ts(i,j)*(ts(i,j)*4.*cb(5,ib)
1067     * +3.*cb(4,ib))+2.*cb(3,ib))+cb(2,ib))*emiss(i,j,ib)
1068    
1069     enddo
1070     enddo
1071    
1072     do k=1,np+1
1073     do j=1,n
1074     do i=1,m
1075     dlayer(i,j,k)=blayer(i,j,k-1)-blayer(i,j,k)
1076     enddo
1077     enddo
1078     enddo
1079    
1080     c-----compute column-integrated absorber amoumt, absorber-weighted
1081     c pressure and temperature for co2 (band 3) and o3 (band 5).
1082     c it follows eqs. (37) and (38).
1083    
1084     c-----this is in the band loop to save storage
1085    
1086     if (high .and. co2bnd) then
1087     call column(m,n,np,pa,dt,dco2,sco3,scopre,scotem)
1088     endif
1089    
1090     if (oznbnd) then
1091     call column(m,n,np,pa,dt,do3,sco3,scopre,scotem)
1092     endif
1093    
1094     c-----compute cloud optical thickness
1095    
1096     if (cldwater) then
1097     do k= 1, np
1098     do j= 1, n
1099     do i= 1, m
1100     taucl(i,j,k,1)=cwp(i,j,k,1)*(aib(1,ib)+aib(2,ib)/
1101     * reff(i,j,k,1)**aib(3,ib))
1102     taucl(i,j,k,2)=cwp(i,j,k,2)*(awb(1,ib)+(awb(2,ib)+(awb(3,ib)
1103     * +awb(4,ib)*reff(i,j,k,2))*reff(i,j,k,2))*reff(i,j,k,2))
1104     enddo
1105     enddo
1106     enddo
1107     endif
1108    
1109     c-----compute cloud single-scattering albedo and asymmetry factor for
1110     c a mixture of ice particles, liquid drops, and rain drops
1111     c single-scattering albedo and asymmetry factor of rain are set
1112     c to 0.54 and 0.95, respectively.
1113    
1114     do k= 1, np
1115     do j= 1, n
1116     do i= 1, m
1117    
1118     clr(i,j,k) = 1.0
1119     taux=taucl(i,j,k,1)+taucl(i,j,k,2)+taucl(i,j,k,3)
1120    
1121     if (taux.gt.0.02 .and. fcld(i,j,k).gt.0.01) then
1122    
1123     reff1=min(reff(i,j,k,1),130.)
1124     reff2=min(reff(i,j,k,2),20.0)
1125    
1126     w1=taucl(i,j,k,1)*(aiw(1,ib)+(aiw(2,ib)+(aiw(3,ib)
1127     * +aiw(4,ib)*reff1)*reff1)*reff1)
1128     w2=taucl(i,j,k,2)*(aww(1,ib)+(aww(2,ib)+(aww(3,ib)
1129     * +aww(4,ib)*reff2)*reff2)*reff2)
1130     w3=taucl(i,j,k,3)*0.54
1131     ww=(w1+w2+w3)/taux
1132    
1133     g1=w1*(aig(1,ib)+(aig(2,ib)+(aig(3,ib)
1134     * +aig(4,ib)*reff1)*reff1)*reff1)
1135     g2=w2*(awg(1,ib)+(awg(2,ib)+(awg(3,ib)
1136     * +awg(4,ib)*reff2)*reff2)*reff2)
1137     g3=w3*0.95
1138    
1139     gg=(g1+g2+g3)/(w1+w2+w3)
1140    
1141     c-----parameterization of LW scattering
1142     c compute forward-scattered fraction and scale optical thickness
1143    
1144     ff=0.5+(0.3739+(0.0076+0.1185*gg)*gg)*gg
1145     taux=taux*(1.-ww*ff)
1146    
1147     c-----compute equivalent cloud-free fraction, clr, for each layer
1148     c the cloud diffuse transmittance is approximated by using
1149     c a diffusivity factor of 1.66.
1150    
1151     clr(i,j,k)=1.0-(fcld(i,j,k)*(1.0-exp(-1.66*taux)))
1152    
1153     endif
1154    
1155     enddo
1156     enddo
1157     enddo
1158    
1159     c-----compute the exponential terms (eq. 32) at each layer due to
1160     c water vapor line absorption when k-distribution is used
1161    
1162     if (.not.h2otbl .and. .not.b10bnd) then
1163     call h2oexps(ib,m,n,np,dh2o,pa,dt,xkw,aw,bw,pm,mw,h2oexp)
1164     endif
1165    
1166     c-----compute the exponential terms (eq. 46) at each layer due to
1167     c water vapor continuum absorption
1168    
1169     if (conbnd) then
1170     call conexps(ib,m,n,np,dcont,xke,conexp)
1171     endif
1172    
1173     c-----compute the exponential terms (eq. 32) at each layer due to
1174     c co2 absorption
1175    
1176     if (.not.high .and. co2bnd) then
1177     call co2exps(m,n,np,dco2,pa,dt,co2exp)
1178     endif
1179    
1180     c***** for trace gases *****
1181    
1182     if (trace) then
1183    
1184     c-----compute the exponential terms at each layer due to n2o absorption
1185    
1186     if (n2obnd) then
1187     call n2oexps(ib,m,n,np,dn2o,pa,dt,n2oexp)
1188     endif
1189    
1190     c-----compute the exponential terms at each layer due to ch4 absorption
1191    
1192     if (ch4bnd) then
1193     call ch4exps(ib,m,n,np,dch4,pa,dt,ch4exp)
1194     endif
1195    
1196     c-----compute the exponential terms due to co2 minor absorption
1197    
1198     if (combnd) then
1199     call comexps(ib,m,n,np,dco2,dt,comexp)
1200     endif
1201    
1202     c-----compute the exponential terms due to cfc11 absorption
1203    
1204     if (f11bnd) then
1205     a1 = 1.26610e-3
1206     b1 = 3.55940e-6
1207     fk1 = 1.89736e+1
1208     a2 = 8.19370e-4
1209     b2 = 4.67810e-6
1210     fk2 = 1.01487e+1
1211     call cfcexps(ib,m,n,np,a1,b1,fk1,a2,b2,fk2,df11,dt,f11exp)
1212     endif
1213    
1214     c-----compute the exponential terms due to cfc12 absorption
1215    
1216     if (f12bnd) then
1217     a1 = 8.77370e-4
1218     b1 =-5.88440e-6
1219     fk1 = 1.58104e+1
1220     a2 = 8.62000e-4
1221     b2 =-4.22500e-6
1222     fk2 = 3.70107e+1
1223     call cfcexps(ib,m,n,np,a1,b1,fk1,a2,b2,fk2,df12,dt,f12exp)
1224     endif
1225    
1226     c-----compute the exponential terms due to cfc22 absorption
1227    
1228     if (f22bnd) then
1229     a1 = 9.65130e-4
1230     b1 = 1.31280e-5
1231     fk1 = 6.18536e+0
1232     a2 =-3.00010e-5
1233     b2 = 5.25010e-7
1234     fk2 = 3.27912e+1
1235     call cfcexps(ib,m,n,np,a1,b1,fk1,a2,b2,fk2,df22,dt,f22exp)
1236     endif
1237    
1238     c-----compute the exponential terms at each layer in band 10 due to
1239     c h2o line and continuum, co2, and n2o absorption
1240    
1241     if (b10bnd) then
1242     call b10exps(m,n,np,dh2o,dcont,dco2,dn2o,pa,dt
1243     * ,h2oexp,conexp,co2exp,n2oexp)
1244     endif
1245    
1246     endif
1247    
1248     c-----compute transmittances for regions between levels k1 and k2
1249     c and update the fluxes at the two levels.
1250    
1251    
1252     c-----initialize fluxes
1253    
1254     do k=1,np+1
1255     do j=1,n
1256     do i=1,m
1257     flxu(i,j,k) = 0.0
1258     flxd(i,j,k) = 0.0
1259     flcu(i,j,k) = 0.0
1260     flcd(i,j,k) = 0.0
1261     enddo
1262     enddo
1263     enddo
1264    
1265    
1266     do 2000 k1=1,np
1267    
1268     c-----initialize fclr, th2o, tcon, tco2, and tranal
1269     c fclr is the clear fraction of the line-of-sight
1270     c clrlw, clrmd, and clrhi are the clear-sky fractions of the
1271     c low, middle, and high troposphere, respectively
1272     c tranal is the aerosol transmission function
1273    
1274     do j=1,n
1275     do i=1,m
1276     fclr(i,j) = 1.0
1277     clrlw(i,j) = 1.0
1278     clrmd(i,j) = 1.0
1279     clrhi(i,j) = 1.0
1280     tranal(i,j)= 1.0
1281     enddo
1282     enddo
1283    
1284     c-----for h2o line transmission
1285    
1286     if (.not. h2otbl) then
1287     do ik=1,6
1288     do j=1,n
1289     do i=1,m
1290     th2o(i,j,ik)=1.0
1291     enddo
1292     enddo
1293     enddo
1294     endif
1295    
1296     c-----for h2o continuum transmission
1297    
1298     if (conbnd) then
1299     do iq=1,3
1300     do j=1,n
1301     do i=1,m
1302     tcon(i,j,iq)=1.0
1303     enddo
1304     enddo
1305     enddo
1306     endif
1307    
1308     c-----for co2 transmission using k-distribution method.
1309     c band 3 is divided into 3 sub-bands, but sub-bands 3a and 3c
1310     c are combined in computing the co2 transmittance.
1311    
1312     if (.not.high .and. co2bnd) then
1313     do isb=1,2
1314     do ik=1,6
1315     do j=1,n
1316     do i=1,m
1317     tco2(i,j,ik,isb)=1.0
1318     enddo
1319     enddo
1320     enddo
1321     enddo
1322     endif
1323    
1324     c***** for trace gases *****
1325    
1326     if (trace) then
1327    
1328     c-----for n2o transmission using k-distribution method.
1329    
1330     if (n2obnd) then
1331     do ik=1,4
1332     do j=1,n
1333     do i=1,m
1334     tn2o(i,j,ik)=1.0
1335     enddo
1336     enddo
1337     enddo
1338     endif
1339    
1340     c-----for ch4 transmission using k-distribution method.
1341    
1342     if (ch4bnd) then
1343     do ik=1,4
1344     do j=1,n
1345     do i=1,m
1346     tch4(i,j,ik)=1.0
1347     enddo
1348     enddo
1349     enddo
1350     endif
1351    
1352     c-----for co2-minor transmission using k-distribution method.
1353    
1354     if (combnd) then
1355     do ik=1,2
1356     do j=1,n
1357     do i=1,m
1358     tcom(i,j,ik)=1.0
1359     enddo
1360     enddo
1361     enddo
1362     endif
1363    
1364     c-----for cfc-11 transmission using k-distribution method.
1365    
1366     if (f11bnd) then
1367     do j=1,n
1368     do i=1,m
1369     tf11(i,j)=1.0
1370     enddo
1371     enddo
1372     endif
1373    
1374     c-----for cfc-12 transmission using k-distribution method.
1375    
1376     if (f12bnd) then
1377     do j=1,n
1378     do i=1,m
1379     tf12(i,j)=1.0
1380     enddo
1381     enddo
1382     endif
1383    
1384     c-----for cfc-22 transmission when using k-distribution method.
1385    
1386     if (f22bnd) then
1387     do j=1,n
1388     do i=1,m
1389     tf22(i,j)=1.0
1390     enddo
1391     enddo
1392     endif
1393    
1394     c-----for the transmission in band 10 using k-distribution method.
1395    
1396     if (b10bnd) then
1397     do ik=1,6
1398     do j=1,n
1399     do i=1,m
1400     th2o(i,j,ik)=1.0
1401     tco2(i,j,ik,1)=1.0
1402     enddo
1403     enddo
1404     enddo
1405     do iq=1,3
1406     do j=1,n
1407     do i=1,m
1408     tcon(i,j,iq)=1.0
1409     enddo
1410     enddo
1411     enddo
1412     do ik=1,4
1413     do j=1,n
1414     do i=1,m
1415     tn2o(i,j,ik)=1.0
1416     enddo
1417     enddo
1418     enddo
1419     endif
1420    
1421     endif
1422    
1423     c-----loop over the bottom level of the region (k2)
1424    
1425     do 3000 k2=k1+1,np+1
1426    
1427     c-----initialize total transmission function (trant)
1428    
1429     do j=1,n
1430     do i=1,m
1431     trant(i,j)=1.0
1432     enddo
1433     enddo
1434    
1435     if (h2otbl) then
1436     w1=-8.0
1437     p1=-2.0
1438     dwe=0.3
1439     dpe=0.2
1440    
1441     c-----compute water vapor transmittance using table look-up
1442    
1443     if (ib.eq.1) then
1444     call tablup(k1,k2,m,n,np,nx,nh,nt,sh2o,swpre,swtem,
1445     * w1,p1,dwe,dpe,h11,h12,h13,trant)
1446    
1447     endif
1448     if (ib.eq.2) then
1449     call tablup(k1,k2,m,n,np,nx,nh,nt,sh2o,swpre,swtem,
1450     * w1,p1,dwe,dpe,h21,h22,h23,trant)
1451    
1452     endif
1453     if (ib.eq.8) then
1454     call tablup(k1,k2,m,n,np,nx,nh,nt,sh2o,swpre,swtem,
1455     * w1,p1,dwe,dpe,h81,h82,h83,trant)
1456     endif
1457    
1458     else
1459    
1460     c-----compute water vapor transmittance using k-distribution
1461    
1462     if (.not.b10bnd) then
1463     call h2okdis(ib,m,n,np,k2-1,fkw,gkw,ne,h2oexp,conexp,
1464     * th2o,tcon,trant)
1465     endif
1466    
1467     endif
1468    
1469     if (co2bnd) then
1470    
1471     if (high) then
1472    
1473     c-----compute co2 transmittance using table look-up method
1474    
1475     w1=-4.0
1476     p1=-2.0
1477     dwe=0.3
1478     dpe=0.2
1479     call tablup(k1,k2,m,n,np,nx,nc,nt,sco3,scopre,scotem,
1480     * w1,p1,dwe,dpe,c1,c2,c3,trant)
1481     else
1482    
1483     c-----compute co2 transmittance using k-distribution method
1484     call co2kdis(m,n,np,k2-1,co2exp,tco2,trant)
1485    
1486     endif
1487    
1488     endif
1489    
1490     c-----All use table look-up to compute o3 transmittance.
1491    
1492     if (oznbnd) then
1493     w1=-6.0
1494     p1=-2.0
1495     dwe=0.3
1496     dpe=0.2
1497     call tablup(k1,k2,m,n,np,nx,no,nt,sco3,scopre,scotem,
1498     * w1,p1,dwe,dpe,o1,o2,o3,trant)
1499     endif
1500    
1501     c***** for trace gases *****
1502    
1503     if (trace) then
1504    
1505     c-----compute n2o transmittance using k-distribution method
1506    
1507     if (n2obnd) then
1508     call n2okdis(ib,m,n,np,k2-1,n2oexp,tn2o,trant)
1509     endif
1510    
1511     c-----compute ch4 transmittance using k-distribution method
1512    
1513     if (ch4bnd) then
1514     call ch4kdis(ib,m,n,np,k2-1,ch4exp,tch4,trant)
1515     endif
1516    
1517     c-----compute co2-minor transmittance using k-distribution method
1518    
1519     if (combnd) then
1520     call comkdis(ib,m,n,np,k2-1,comexp,tcom,trant)
1521     endif
1522    
1523     c-----compute cfc11 transmittance using k-distribution method
1524    
1525     if (f11bnd) then
1526     call cfckdis(m,n,np,k2-1,f11exp,tf11,trant)
1527     endif
1528    
1529     c-----compute cfc12 transmittance using k-distribution method
1530    
1531     if (f12bnd) then
1532     call cfckdis(m,n,np,k2-1,f12exp,tf12,trant)
1533     endif
1534    
1535     c-----compute cfc22 transmittance using k-distribution method
1536    
1537     if (f22bnd) then
1538     call cfckdis(m,n,np,k2-1,f22exp,tf22,trant)
1539     endif
1540    
1541     c-----compute transmittance in band 10 using k-distribution method
1542     c here, trant is the change in transmittance due to n2o absorption
1543    
1544     if (b10bnd) then
1545     call b10kdis(m,n,np,k2-1,h2oexp,conexp,co2exp,n2oexp
1546     * ,th2o,tcon,tco2,tn2o,trant)
1547     endif
1548    
1549     endif
1550     c
1551     c-----include aerosol effect
1552     c
1553     do j=1,n
1554     do i=1,m
1555     ff=0.5+(0.3739+(0.0076+0.1185*asyal(i,j,k2-1,ib))
1556     * *asyal(i,j,k2-1,ib))*asyal(i,j,k2-1,ib)
1557     taux=taual(i,j,k2-1,ib)*(1.-ssaal(i,j,k2-1,ib)*ff)
1558     tranal(i,j)=tranal(i,j)*exp(-1.66*taux)
1559     trant (i,j)=trant(i,j) *tranal(i,j)
1560     enddo
1561     enddo
1562    
1563     c-----Identify the clear-sky fractions clrhi, clrmd and clrlw, in the
1564     c high, middle and low cloud groups.
1565     c fclr is the clear-sky fraction between levels k1 and k2 assuming
1566     c the three cloud groups are randomly overlapped.
1567    
1568     do j=1,n
1569     do i=1,m
1570     if( k2 .le. ict ) then
1571     clrhi(i,j)=min(clr(i,j,k2-1),clrhi(i,j))
1572     elseif( k2 .gt. ict .and. k2 .le. icb ) then
1573     clrmd(i,j)=min(clr(i,j,k2-1),clrmd(i,j))
1574     elseif( k2 .gt. icb ) then
1575     clrlw(i,j)=min(clr(i,j,k2-1),clrlw(i,j))
1576     endif
1577     fclr(i,j)=clrlw(i,j)*clrmd(i,j)*clrhi(i,j)
1578    
1579     enddo
1580     enddo
1581    
1582     c-----compute upward and downward fluxes (band 1-9).
1583     c add "boundary" terms to the net downward flux.
1584     c these are the first terms on the right-hand-side of
1585     c eqs. (56a) and (56b). downward fluxes are positive.
1586    
1587     if (.not. b10bnd) then
1588    
1589     if (k2 .eq. k1+1) then
1590    
1591     do j=1,n
1592     do i=1,m
1593    
1594     c-----compute upward and downward fluxes for clear-sky situation
1595    
1596     flcu(i,j,k1)=flcu(i,j,k1)-blayer(i,j,k1)
1597     flcd(i,j,k2)=flcd(i,j,k2)+blayer(i,j,k1)
1598    
1599     c-----compute upward and downward fluxes for all-sky situation
1600    
1601     flxu(i,j,k1)=flxu(i,j,k1)-blayer(i,j,k1)
1602     flxd(i,j,k2)=flxd(i,j,k2)+blayer(i,j,k1)
1603    
1604     enddo
1605     enddo
1606    
1607     endif
1608    
1609     c-----add flux components involving the four layers above and below
1610     c the levels k1 and k2. it follows eqs. (56a) and (56b).
1611    
1612     do j=1,n
1613     do i=1,m
1614     xx=trant(i,j)*dlayer(i,j,k2)
1615     flcu(i,j,k1) =flcu(i,j,k1)+xx
1616     flxu(i,j,k1) =flxu(i,j,k1)+xx*fclr(i,j)
1617     xx=trant(i,j)*dlayer(i,j,k1)
1618     flcd(i,j,k2) =flcd(i,j,k2)+xx
1619     flxd(i,j,k2) =flxd(i,j,k2)+xx*fclr(i,j)
1620     enddo
1621     enddo
1622    
1623     else
1624    
1625     c-----flux reduction due to n2o in band 10
1626    
1627     if (trace) then
1628    
1629     do j=1,n
1630     do i=1,m
1631     rflx(i,j,k1) = rflx(i,j,k1)+trant(i,j)*fclr(i,j)
1632     * *dlayer(i,j,k2)
1633     rflx(i,j,k2) = rflx(i,j,k2)+trant(i,j)*fclr(i,j)
1634     * *dlayer(i,j,k1)
1635     rflc(i,j,k1) = rflc(i,j,k1)+trant(i,j)
1636     * *dlayer(i,j,k2)
1637     rflc(i,j,k2) = rflc(i,j,k2)+trant(i,j)
1638     * *dlayer(i,j,k1)
1639     enddo
1640     enddo
1641    
1642     endif
1643    
1644     endif
1645    
1646     3000 continue
1647    
1648    
1649     c-----transmission between level k1 and the surface
1650    
1651     do j=1,n
1652     do i=1,m
1653     trantcr(i,j,k1) =trant(i,j)
1654     transfc(i,j,k1) =trant(i,j)*fclr(i,j)
1655     enddo
1656     enddo
1657    
1658     c-----compute the partial derivative of fluxes with respect to
1659     c surface temperature (eq. 59)
1660    
1661     if (trace .or. (.not. b10bnd) ) then
1662    
1663     do j=1,n
1664     do i=1,m
1665     dfdts(i,j,k1) =dfdts(i,j,k1)-dbs(i,j)*transfc(i,j,k1)
1666     enddo
1667     enddo
1668    
1669     endif
1670    
1671     2000 continue
1672    
1673     if (.not. b10bnd) then
1674    
1675     c-----add contribution from the surface to the flux terms at the surface
1676    
1677     do j=1,n
1678     do i=1,m
1679     flcu(i,j,np+1)=flcu(i,j,np+1)-blayer(i,j,np+1)
1680     flxu(i,j,np+1)=flxu(i,j,np+1)-blayer(i,j,np+1)
1681     st4(i,j)=st4(i,j)-blayer(i,j,np+1)
1682     dfdts(i,j,np+1)=dfdts(i,j,np+1)-dbs(i,j)
1683     enddo
1684     enddo
1685    
1686     c-----add the flux component reflected by the surface
1687     c note: upward flux is negative
1688    
1689     do k=1, np+1
1690     do j=1,n
1691     do i=1,m
1692     flcu(i,j,k)=flcu(i,j,k)-
1693     * flcd(i,j,np+1)*trantcr(i,j,k)*(1.-emiss(i,j,ib))
1694     flxu(i,j,k)=flxu(i,j,k)-
1695     * flxd(i,j,np+1)*transfc(i,j,k)*(1.-emiss(i,j,ib))
1696     enddo
1697     enddo
1698     enddo
1699    
1700     endif
1701    
1702     c-----summation of fluxes over spectral bands
1703    
1704     do k=1,np+1
1705     do j=1,n
1706     do i=1,m
1707     flc(i,j,k)=flc(i,j,k)+flcd(i,j,k)+flcu(i,j,k)
1708     flx(i,j,k)=flx(i,j,k)+flxd(i,j,k)+flxu(i,j,k)
1709     enddo
1710     enddo
1711     enddo
1712    
1713     1000 continue
1714    
1715     c-----adjust fluxes due to n2o absorption in band 10
1716    
1717     do k=1,np+1
1718     do j=1,n
1719     do i=1,m
1720     flc(i,j,k)=flc(i,j,k)+rflc(i,j,k)
1721     flx(i,j,k)=flx(i,j,k)+rflx(i,j,k)
1722     enddo
1723     enddo
1724     enddo
1725    
1726     return
1727     end
1728     c***********************************************************************
1729     subroutine column (m,n,np,pa,dt,sabs0,sabs,spre,stem)
1730     c***********************************************************************
1731     c-----compute column-integrated (from top of the model atmosphere)
1732     c absorber amount (sabs), absorber-weighted pressure (spre) and
1733     c temperature (stem).
1734     c computations of spre and stem follows eqs. (37) and (38).
1735     c
1736     c--- input parameters
1737     c number of soundings in zonal direction (m)
1738     c number of soundings in meridional direction (n)
1739     c number of atmospheric layers (np)
1740     c layer pressure (pa)
1741     c layer temperature minus 250K (dt)
1742     c layer absorber amount (sabs0)
1743     c
1744     c--- output parameters
1745     c column-integrated absorber amount (sabs)
1746     c column absorber-weighted pressure (spre)
1747     c column absorber-weighted temperature (stem)
1748     c
1749     c--- units of pa and dt are mb and k, respectively.
1750     c units of sabs are g/cm**2 for water vapor and (cm-atm)stp
1751     c for co2 and o3
1752     c***********************************************************************
1753     implicit none
1754     integer m,n,np,i,j,k
1755    
1756     c---- input parameters -----
1757    
1758     real pa(m,n,np),dt(m,n,np),sabs0(m,n,np)
1759    
1760     c---- output parameters -----
1761    
1762     real sabs(m,n,np+1),spre(m,n,np+1),stem(m,n,np+1)
1763    
1764     c*********************************************************************
1765     do j=1,n
1766     do i=1,m
1767     sabs(i,j,1)=0.0
1768     spre(i,j,1)=0.0
1769     stem(i,j,1)=0.0
1770     enddo
1771     enddo
1772    
1773     do k=1,np
1774     do j=1,n
1775     do i=1,m
1776     sabs(i,j,k+1)=sabs(i,j,k)+sabs0(i,j,k)
1777     spre(i,j,k+1)=spre(i,j,k)+pa(i,j,k)*sabs0(i,j,k)
1778     stem(i,j,k+1)=stem(i,j,k)+dt(i,j,k)*sabs0(i,j,k)
1779     enddo
1780     enddo
1781     enddo
1782    
1783     return
1784     end
1785     c**********************************************************************
1786     subroutine h2oexps(ib,m,n,np,dh2o,pa,dt,xkw,aw,bw,pm,mw,h2oexp)
1787     c**********************************************************************
1788     c compute exponentials for water vapor line absorption
1789     c in individual layers.
1790     c
1791     c---- input parameters
1792     c spectral band (ib)
1793     c number of grid intervals in zonal direction (m)
1794     c number of grid intervals in meridional direction (n)
1795     c number of layers (np)
1796     c layer water vapor amount for line absorption (dh2o)
1797     c layer pressure (pa)
1798     c layer temperature minus 250K (dt)
1799     c absorption coefficients for the first k-distribution
1800     c function due to h2o line absorption (xkw)
1801     c coefficients for the temperature and pressure scaling (aw,bw,pm)
1802     c ratios between neighboring absorption coefficients for
1803     c h2o line absorption (mw)
1804     c
1805     c---- output parameters
1806     c 6 exponentials for each layer (h2oexp)
1807     c**********************************************************************
1808     implicit none
1809     integer ib,m,n,np,i,j,k,ik
1810    
1811     c---- input parameters ------
1812    
1813     real dh2o(m,n,np),pa(m,n,np),dt(m,n,np)
1814    
1815     c---- output parameters -----
1816    
1817     real h2oexp(m,n,np,6)
1818    
1819     c---- static data -----
1820    
1821     integer mw(9)
1822     real xkw(9),aw(9),bw(9),pm(9)
1823    
1824     c---- temporary arrays -----
1825    
1826     real xh,xh1
1827    
1828     c**********************************************************************
1829     c note that the 3 sub-bands in band 3 use the same set of xkw, aw,
1830     c and bw, therefore, h2oexp for these sub-bands are identical.
1831     c**********************************************************************
1832    
1833     do k=1,np
1834     do j=1,n
1835     do i=1,m
1836    
1837     c-----xh is the scaled water vapor amount for line absorption
1838     c computed from (27)
1839    
1840     xh = dh2o(i,j,k)*(pa(i,j,k)/500.)**pm(ib)
1841     1 * ( 1.+(aw(ib)+bw(ib)* dt(i,j,k))*dt(i,j,k) )
1842    
1843     c-----h2oexp is the water vapor transmittance of the layer k
1844     c due to line absorption
1845    
1846     h2oexp(i,j,k,1) = exp(-xh*xkw(ib))
1847    
1848     enddo
1849     enddo
1850     enddo
1851    
1852     do ik=2,6
1853    
1854     if (mw(ib).eq.6) then
1855    
1856     do k=1,np
1857     do j=1,n
1858     do i=1,m
1859     xh = h2oexp(i,j,k,ik-1)*h2oexp(i,j,k,ik-1)
1860     h2oexp(i,j,k,ik) = xh*xh*xh
1861     enddo
1862     enddo
1863     enddo
1864    
1865     elseif (mw(ib).eq.8) then
1866    
1867     do k=1,np
1868     do j=1,n
1869     do i=1,m
1870     xh = h2oexp(i,j,k,ik-1)*h2oexp(i,j,k,ik-1)
1871     xh = xh*xh
1872     h2oexp(i,j,k,ik) = xh*xh
1873     enddo
1874     enddo
1875     enddo
1876    
1877     elseif (mw(ib).eq.9) then
1878    
1879     do k=1,np
1880     do j=1,n
1881     do i=1,m
1882     xh=h2oexp(i,j,k,ik-1)*h2oexp(i,j,k,ik-1)*h2oexp(i,j,k,ik-1)
1883     xh1 = xh*xh
1884     h2oexp(i,j,k,ik) = xh*xh1
1885     enddo
1886     enddo
1887     enddo
1888    
1889     else
1890    
1891     do k=1,np
1892     do j=1,n
1893     do i=1,m
1894     xh = h2oexp(i,j,k,ik-1)*h2oexp(i,j,k,ik-1)
1895     xh = xh*xh
1896     xh = xh*xh
1897     h2oexp(i,j,k,ik) = xh*xh
1898     enddo
1899     enddo
1900     enddo
1901    
1902     endif
1903     enddo
1904    
1905     return
1906     end
1907     c**********************************************************************
1908     subroutine conexps(ib,m,n,np,dcont,xke,conexp)
1909     c**********************************************************************
1910     c compute exponentials for continuum absorption in individual layers.
1911     c
1912     c---- input parameters
1913     c spectral band (ib)
1914     c number of grid intervals in zonal direction (m)
1915     c number of grid intervals in meridional direction (n)
1916     c number of layers (np)
1917     c layer scaled water vapor amount for continuum absorption (dcont)
1918     c absorption coefficients for the first k-distribution function
1919     c due to water vapor continuum absorption (xke)
1920     c
1921     c---- output parameters
1922     c 1 or 3 exponentials for each layer (conexp)
1923     c**********************************************************************
1924     implicit none
1925     integer ib,m,n,np,i,j,k,iq
1926    
1927     c---- input parameters ------
1928    
1929     real dcont(m,n,np)
1930    
1931     c---- updated parameters -----
1932    
1933     real conexp(m,n,np,3)
1934    
1935     c---- static data -----
1936    
1937     real xke(9)
1938    
1939     c**********************************************************************
1940    
1941     do k=1,np
1942     do j=1,n
1943     do i=1,m
1944     conexp(i,j,k,1) = exp(-dcont(i,j,k)*xke(ib))
1945     enddo
1946     enddo
1947     enddo
1948    
1949     if (ib .eq. 3) then
1950    
1951     c-----the absorption coefficients for sub-bands 3b (iq=2) and 3a (iq=3)
1952     c are, respectively, two and three times the absorption coefficient
1953     c for sub-band 3c (iq=1) (table 6).
1954    
1955     do iq=2,3
1956     do k=1,np
1957     do j=1,n
1958     do i=1,m
1959     conexp(i,j,k,iq) = conexp(i,j,k,iq-1) *conexp(i,j,k,iq-1)
1960     enddo
1961     enddo
1962     enddo
1963     enddo
1964    
1965     endif
1966    
1967     return
1968     end
1969     c**********************************************************************
1970     subroutine co2exps(m,n,np,dco2,pa,dt,co2exp)
1971     c**********************************************************************
1972     c compute co2 exponentials for individual layers.
1973     c
1974     c---- input parameters
1975     c number of grid intervals in zonal direction (m)
1976     c number of grid intervals in meridional direction (n)
1977     c number of layers (np)
1978     c layer co2 amount (dco2)
1979     c layer pressure (pa)
1980     c layer temperature minus 250K (dt)
1981     c
1982     c---- output parameters
1983     c 6 exponentials for each layer (co2exp)
1984     c**********************************************************************
1985     implicit none
1986     integer m,n,np,i,j,k
1987    
1988     c---- input parameters -----
1989    
1990     real dco2(m,n,np),pa(m,n,np),dt(m,n,np)
1991    
1992     c---- output parameters -----
1993    
1994     real co2exp(m,n,np,6,2)
1995    
1996     c---- temporary arrays -----
1997    
1998     real xc
1999    
2000     c**********************************************************************
2001    
2002     do k=1,np
2003     do j=1,n
2004     do i=1,m
2005    
2006     c-----compute the scaled co2 amount from eq. (27) for band-wings
2007     c (sub-bands 3a and 3c).
2008    
2009     xc = dco2(i,j,k)*(pa(i,j,k)/300.0)**0.5
2010     1 *(1.+(0.0182+1.07e-4*dt(i,j,k))*dt(i,j,k))
2011    
2012     c-----six exponentials by powers of 8 (table 7).
2013    
2014     co2exp(i,j,k,1,1)=exp(-xc*2.656e-5)
2015    
2016     xc=co2exp(i,j,k,1,1)*co2exp(i,j,k,1,1)
2017     xc=xc*xc
2018     co2exp(i,j,k,2,1)=xc*xc
2019    
2020     xc=co2exp(i,j,k,2,1)*co2exp(i,j,k,2,1)
2021     xc=xc*xc
2022     co2exp(i,j,k,3,1)=xc*xc
2023    
2024     xc=co2exp(i,j,k,3,1)*co2exp(i,j,k,3,1)
2025     xc=xc*xc
2026     co2exp(i,j,k,4,1)=xc*xc
2027    
2028     xc=co2exp(i,j,k,4,1)*co2exp(i,j,k,4,1)
2029     xc=xc*xc
2030     co2exp(i,j,k,5,1)=xc*xc
2031    
2032     xc=co2exp(i,j,k,5,1)*co2exp(i,j,k,5,1)
2033     xc=xc*xc
2034     co2exp(i,j,k,6,1)=xc*xc
2035    
2036     c-----compute the scaled co2 amount from eq. (27) for band-center
2037     c region (sub-band 3b).
2038    
2039     xc = dco2(i,j,k)*(pa(i,j,k)/30.0)**0.85
2040     1 *(1.+(0.0042+2.00e-5*dt(i,j,k))*dt(i,j,k))
2041    
2042     co2exp(i,j,k,1,2)=exp(-xc*2.656e-3)
2043    
2044     xc=co2exp(i,j,k,1,2)*co2exp(i,j,k,1,2)
2045     xc=xc*xc
2046     co2exp(i,j,k,2,2)=xc*xc
2047    
2048     xc=co2exp(i,j,k,2,2)*co2exp(i,j,k,2,2)
2049     xc=xc*xc
2050     co2exp(i,j,k,3,2)=xc*xc
2051    
2052     xc=co2exp(i,j,k,3,2)*co2exp(i,j,k,3,2)
2053     xc=xc*xc
2054     co2exp(i,j,k,4,2)=xc*xc
2055    
2056     xc=co2exp(i,j,k,4,2)*co2exp(i,j,k,4,2)
2057     xc=xc*xc
2058     co2exp(i,j,k,5,2)=xc*xc
2059    
2060     xc=co2exp(i,j,k,5,2)*co2exp(i,j,k,5,2)
2061     xc=xc*xc
2062     co2exp(i,j,k,6,2)=xc*xc
2063    
2064     enddo
2065     enddo
2066     enddo
2067    
2068     return
2069     end
2070     c**********************************************************************
2071     subroutine n2oexps(ib,m,n,np,dn2o,pa,dt,n2oexp)
2072     c**********************************************************************
2073     c compute n2o exponentials for individual layers.
2074     c
2075     c---- input parameters
2076     c spectral band (ib)
2077     c number of grid intervals in zonal direction (m)
2078     c number of grid intervals in meridional direction (n)
2079     c number of layers (np)
2080     c layer n2o amount (dn2o)
2081     c layer pressure (pa)
2082     c layer temperature minus 250K (dt)
2083     c
2084     c---- output parameters
2085     c 2 or 4 exponentials for each layer (n2oexp)
2086     c**********************************************************************
2087     implicit none
2088     integer ib,m,n,np,i,j,k
2089    
2090     c---- input parameters -----
2091    
2092     real dn2o(m,n,np),pa(m,n,np),dt(m,n,np)
2093    
2094     c---- output parameters -----
2095    
2096     real n2oexp(m,n,np,4)
2097    
2098     c---- temporary arrays -----
2099    
2100     real xc,xc1,xc2
2101    
2102     c**********************************************************************
2103    
2104     do k=1,np
2105     do j=1,n
2106     do i=1,m
2107    
2108     c-----four exponential by powers of 21 for band 6
2109    
2110     if (ib.eq.6) then
2111    
2112     xc=dn2o(i,j,k)*(1.+(1.9297e-3+4.3750e-6*dt(i,j,k))*dt(i,j,k))
2113     n2oexp(i,j,k,1)=exp(-xc*6.31582e-2)
2114    
2115     xc=n2oexp(i,j,k,1)*n2oexp(i,j,k,1)*n2oexp(i,j,k,1)
2116     xc1=xc*xc
2117     xc2=xc1*xc1
2118     n2oexp(i,j,k,2)=xc*xc1*xc2
2119    
2120     c-----four exponential by powers of 8 for band 7
2121    
2122     else
2123    
2124     xc=dn2o(i,j,k)*(pa(i,j,k)/500.0)**0.48
2125     * *(1.+(1.3804e-3+7.4838e-6*dt(i,j,k))*dt(i,j,k))
2126     n2oexp(i,j,k,1)=exp(-xc*5.35779e-2)
2127    
2128     xc=n2oexp(i,j,k,1)*n2oexp(i,j,k,1)
2129     xc=xc*xc
2130     n2oexp(i,j,k,2)=xc*xc
2131     xc=n2oexp(i,j,k,2)*n2oexp(i,j,k,2)
2132     xc=xc*xc
2133     n2oexp(i,j,k,3)=xc*xc
2134     xc=n2oexp(i,j,k,3)*n2oexp(i,j,k,3)
2135     xc=xc*xc
2136     n2oexp(i,j,k,4)=xc*xc
2137    
2138     endif
2139    
2140     enddo
2141     enddo
2142     enddo
2143    
2144     return
2145     end
2146     c**********************************************************************
2147     subroutine ch4exps(ib,m,n,np,dch4,pa,dt,ch4exp)
2148     c**********************************************************************
2149     c compute ch4 exponentials for individual layers.
2150     c
2151     c---- input parameters
2152     c spectral band (ib)
2153     c number of grid intervals in zonal direction (m)
2154     c number of grid intervals in meridional direction (n)
2155     c number of layers (np)
2156     c layer ch4 amount (dch4)
2157     c layer pressure (pa)
2158     c layer temperature minus 250K (dt)
2159     c
2160     c---- output parameters
2161     c 1 or 4 exponentials for each layer (ch4exp)
2162     c**********************************************************************
2163     implicit none
2164     integer ib,m,n,np,i,j,k
2165    
2166     c---- input parameters -----
2167    
2168     real dch4(m,n,np),pa(m,n,np),dt(m,n,np)
2169    
2170     c---- output parameters -----
2171    
2172     real ch4exp(m,n,np,4)
2173    
2174     c---- temporary arrays -----
2175    
2176     real xc
2177    
2178     c**********************************************************************
2179    
2180     do k=1,np
2181     do j=1,n
2182     do i=1,m
2183    
2184     c-----four exponentials for band 6
2185    
2186     if (ib.eq.6) then
2187    
2188     xc=dch4(i,j,k)*(1.+(1.7007e-2+1.5826e-4*dt(i,j,k))*dt(i,j,k))
2189     ch4exp(i,j,k,1)=exp(-xc*5.80708e-3)
2190    
2191     c-----four exponentials by powers of 12 for band 7
2192    
2193     else
2194    
2195     xc=dch4(i,j,k)*(pa(i,j,k)/500.0)**0.65
2196     * *(1.+(5.9590e-4-2.2931e-6*dt(i,j,k))*dt(i,j,k))
2197     ch4exp(i,j,k,1)=exp(-xc*6.29247e-2)
2198    
2199     xc=ch4exp(i,j,k,1)*ch4exp(i,j,k,1)*ch4exp(i,j,k,1)
2200     xc=xc*xc
2201     ch4exp(i,j,k,2)=xc*xc
2202    
2203     xc=ch4exp(i,j,k,2)*ch4exp(i,j,k,2)*ch4exp(i,j,k,2)
2204     xc=xc*xc
2205     ch4exp(i,j,k,3)=xc*xc
2206    
2207     xc=ch4exp(i,j,k,3)*ch4exp(i,j,k,3)*ch4exp(i,j,k,3)
2208     xc=xc*xc
2209     ch4exp(i,j,k,4)=xc*xc
2210    
2211     endif
2212    
2213     enddo
2214     enddo
2215     enddo
2216    
2217     return
2218     end
2219     c**********************************************************************
2220     subroutine comexps(ib,m,n,np,dcom,dt,comexp)
2221     c**********************************************************************
2222     c compute co2-minor exponentials for individual layers.
2223     c
2224     c---- input parameters
2225     c spectral band (ib)
2226     c number of grid intervals in zonal direction (m)
2227     c number of grid intervals in meridional direction (n)
2228     c number of layers (np)
2229     c layer co2 amount (dcom)
2230     c layer temperature minus 250K (dt)
2231     c
2232     c---- output parameters
2233     c 2 exponentials for each layer (comexp)
2234     c**********************************************************************
2235     implicit none
2236     integer ib,m,n,np,i,j,k
2237    
2238     c---- input parameters -----
2239    
2240     real dcom(m,n,np),dt(m,n,np)
2241    
2242     c---- output parameters -----
2243    
2244     real comexp(m,n,np,2)
2245    
2246     c---- temporary arrays -----
2247    
2248     real xc,xc1,xc2
2249    
2250     c**********************************************************************
2251    
2252     do k=1,np
2253     do j=1,n
2254     do i=1,m
2255    
2256     c-----two exponentials by powers of 60 for band 4
2257    
2258     if (ib.eq.4) then
2259    
2260     xc=dcom(i,j,k)*(1.+(3.5775e-2+4.0447e-4*dt(i,j,k))*dt(i,j,k))
2261     comexp(i,j,k,1)=exp(-xc*2.18947e-5)
2262    
2263     xc=comexp(i,j,k,1)*comexp(i,j,k,1)*comexp(i,j,k,1)
2264     xc=xc*xc
2265     xc1=xc*xc
2266     xc=xc1*xc1
2267     xc=xc*xc
2268     comexp(i,j,k,2)=xc*xc1
2269    
2270     c-----two exponentials by powers of 44 for band 5
2271    
2272     else
2273    
2274     xc=dcom(i,j,k)*(1.+(3.4268e-2+3.7401e-4*dt(i,j,k))*dt(i,j,k))
2275     comexp(i,j,k,1)=exp(-xc*4.74570e-5)
2276    
2277     xc=comexp(i,j,k,1)*comexp(i,j,k,1)
2278     xc1=xc*xc
2279     xc2=xc1*xc1
2280     xc=xc2*xc2
2281     xc=xc*xc
2282     comexp(i,j,k,2)=xc1*xc2*xc
2283    
2284     endif
2285    
2286     enddo
2287     enddo
2288     enddo
2289    
2290     return
2291     end
2292     c**********************************************************************
2293     subroutine cfcexps(ib,m,n,np,a1,b1,fk1,a2,b2,fk2,dcfc,dt,cfcexp)
2294     c**********************************************************************
2295     c compute cfc(-11, -12, -22) exponentials for individual layers.
2296     c
2297     c---- input parameters
2298     c spectral band (ib)
2299     c number of grid intervals in zonal direction (m)
2300     c number of grid intervals in meridional direction (n)
2301     c number of layers (np)
2302     c parameters for computing the scaled cfc amounts
2303     c for temperature scaling (a1,b1,a2,b2)
2304     c the absorption coefficients for the
2305     c first k-distribution function due to cfcs (fk1,fk2)
2306     c layer cfc amounts (dcfc)
2307     c layer temperature minus 250K (dt)
2308     c
2309     c---- output parameters
2310     c 1 exponential for each layer (cfcexp)
2311     c**********************************************************************
2312     implicit none
2313     integer ib,m,n,np,i,j,k
2314    
2315     c---- input parameters -----
2316    
2317     real dcfc(m,n,np),dt(m,n,np)
2318    
2319     c---- output parameters -----
2320    
2321     real cfcexp(m,n,np)
2322    
2323     c---- static data -----
2324    
2325     real a1,b1,fk1,a2,b2,fk2
2326    
2327     c---- temporary arrays -----
2328    
2329     real xf
2330    
2331     c**********************************************************************
2332    
2333     do k=1,np
2334     do j=1,n
2335     do i=1,m
2336    
2337     c-----compute the scaled cfc amount (xf) and exponential (cfcexp)
2338    
2339     if (ib.eq.4) then
2340     xf=dcfc(i,j,k)*(1.+(a1+b1*dt(i,j,k))*dt(i,j,k))
2341     cfcexp(i,j,k)=exp(-xf*fk1)
2342     else
2343     xf=dcfc(i,j,k)*(1.+(a2+b2*dt(i,j,k))*dt(i,j,k))
2344     cfcexp(i,j,k)=exp(-xf*fk2)
2345     endif
2346    
2347     enddo
2348     enddo
2349     enddo
2350    
2351     return
2352     end
2353     c**********************************************************************
2354     subroutine b10exps(m,n,np,dh2o,dcont,dco2,dn2o,pa,dt
2355     * ,h2oexp,conexp,co2exp,n2oexp)
2356     c**********************************************************************
2357     c compute band3a exponentials for individual layers.
2358     c
2359     c---- input parameters
2360     c number of grid intervals in zonal direction (m)
2361     c number of grid intervals in meridional direction (n)
2362     c number of layers (np)
2363     c layer h2o amount for line absorption (dh2o)
2364     c layer h2o amount for continuum absorption (dcont)
2365     c layer co2 amount (dco2)
2366     c layer n2o amount (dn2o)
2367     c layer pressure (pa)
2368     c layer temperature minus 250K (dt)
2369     c
2370     c---- output parameters
2371     c
2372     c exponentials for each layer (h2oexp,conexp,co2exp,n2oexp)
2373     c**********************************************************************
2374     implicit none
2375     integer m,n,np,i,j,k
2376    
2377     c---- input parameters -----
2378    
2379     real dh2o(m,n,np),dcont(m,n,np),dn2o(m,n,np)
2380     real dco2(m,n,np),pa(m,n,np),dt(m,n,np)
2381    
2382     c---- output parameters -----
2383    
2384     real h2oexp(m,n,np,6),conexp(m,n,np,3),co2exp(m,n,np,6,2)
2385     * ,n2oexp(m,n,np,4)
2386    
2387     c---- temporary arrays -----
2388    
2389     real xx,xx1,xx2,xx3
2390    
2391     c**********************************************************************
2392    
2393     do k=1,np
2394     do j=1,n
2395     do i=1,m
2396    
2397     c-----compute the scaled h2o-line amount for the sub-band 3a
2398     c table 3, Chou et al. (JAS, 1993)
2399    
2400     xx=dh2o(i,j,k)*(pa(i,j,k)/500.0)
2401     1 *(1.+(0.0149+6.20e-5*dt(i,j,k))*dt(i,j,k))
2402    
2403     c-----six exponentials by powers of 8
2404     c the constant 0.10624 is equal to 1.66*0.064
2405    
2406     h2oexp(i,j,k,1)=exp(-xx*0.10624)
2407    
2408     xx=h2oexp(i,j,k,1)*h2oexp(i,j,k,1)
2409     xx=xx*xx
2410     h2oexp(i,j,k,2)=xx*xx
2411    
2412     xx=h2oexp(i,j,k,2)*h2oexp(i,j,k,2)
2413     xx=xx*xx
2414     h2oexp(i,j,k,3)=xx*xx
2415    
2416     xx=h2oexp(i,j,k,3)*h2oexp(i,j,k,3)
2417     xx=xx*xx
2418     h2oexp(i,j,k,4)=xx*xx
2419    
2420     xx=h2oexp(i,j,k,4)*h2oexp(i,j,k,4)
2421     xx=xx*xx
2422     h2oexp(i,j,k,5)=xx*xx
2423    
2424     xx=h2oexp(i,j,k,5)*h2oexp(i,j,k,5)
2425     xx=xx*xx
2426     h2oexp(i,j,k,6)=xx*xx
2427    
2428     c-----compute the scaled co2 amount for the sub-band 3a
2429     c table 1, Chou et al. (JAS, 1993)
2430    
2431     xx=dco2(i,j,k)*(pa(i,j,k)/300.0)**0.5
2432     1 *(1.+(0.0179+1.02e-4*dt(i,j,k))*dt(i,j,k))
2433    
2434     c-----six exponentials by powers of 8
2435     c the constant 2.656e-5 is equal to 1.66*1.60e-5
2436    
2437     co2exp(i,j,k,1,1)=exp(-xx*2.656e-5)
2438    
2439     xx=co2exp(i,j,k,1,1)*co2exp(i,j,k,1,1)
2440     xx=xx*xx
2441     co2exp(i,j,k,2,1)=xx*xx
2442    
2443     xx=co2exp(i,j,k,2,1)*co2exp(i,j,k,2,1)
2444     xx=xx*xx
2445     co2exp(i,j,k,3,1)=xx*xx
2446    
2447     xx=co2exp(i,j,k,3,1)*co2exp(i,j,k,3,1)
2448     xx=xx*xx
2449     co2exp(i,j,k,4,1)=xx*xx
2450    
2451     xx=co2exp(i,j,k,4,1)*co2exp(i,j,k,4,1)
2452     xx=xx*xx
2453     co2exp(i,j,k,5,1)=xx*xx
2454    
2455     xx=co2exp(i,j,k,5,1)*co2exp(i,j,k,5,1)
2456     xx=xx*xx
2457     co2exp(i,j,k,6,1)=xx*xx
2458    
2459     c-----one exponential of h2o continuum for sub-band 3a
2460     c tabl 5 of Chou et. al. (JAS, 1993)
2461     c the constant 1.04995e+2 is equal to 1.66*63.25
2462    
2463     conexp(i,j,k,1)=exp(-dcont(i,j,k)*1.04995e+2)
2464    
2465     c-----compute the scaled n2o amount for sub-band 3a
2466    
2467     xx=dn2o(i,j,k)*(1.+(1.4476e-3+3.6656e-6*dt(i,j,k))*dt(i,j,k))
2468    
2469     c-----two exponential2 by powers of 58
2470    
2471     n2oexp(i,j,k,1)=exp(-xx*0.25238)
2472    
2473     xx=n2oexp(i,j,k,1)*n2oexp(i,j,k,1)
2474     xx1=xx*xx
2475     xx1=xx1*xx1
2476     xx2=xx1*xx1
2477     xx3=xx2*xx2
2478     n2oexp(i,j,k,2)=xx*xx1*xx2*xx3
2479    
2480     enddo
2481     enddo
2482     enddo
2483    
2484     return
2485     end
2486     c**********************************************************************
2487     subroutine tablup(k1,k2,m,n,np,nx,nh,nt,sabs,spre,stem,w1,p1,
2488     * dwe,dpe,coef1,coef2,coef3,tran)
2489     c**********************************************************************
2490     c compute water vapor, co2 and o3 transmittances between levels
2491     c k1 and k2 for m x n soundings, using table look-up.
2492     c
2493     c calculations follow Eq. (40) of Chou and Suarez (1994)
2494     c
2495     c---- input ---------------------
2496     c indices for pressure levels (k1 and k2)
2497     c number of grid intervals in zonal direction (m)
2498     c number of grid intervals in meridional direction (n)
2499     c number of atmospheric layers (np)
2500     c number of pressure intervals in the table (nx)
2501     c number of absorber amount intervals in the table (nh)
2502     c number of tables copied (nt)
2503     c column-integrated absorber amount (sabs)
2504     c column absorber amount-weighted pressure (spre)
2505     c column absorber amount-weighted temperature (stem)
2506     c first value of absorber amount (log10) in the table (w1)
2507     c first value of pressure (log10) in the table (p1)
2508     c size of the interval of absorber amount (log10) in the table (dwe)
2509     c size of the interval of pressure (log10) in the table (dpe)
2510     c pre-computed coefficients (coef1, coef2, and coef3)
2511     c
2512     c---- updated ---------------------
2513     c transmittance (tran)
2514     c
2515     c Note:
2516     c (1) units of sabs are g/cm**2 for water vapor and
2517     c (cm-atm)stp for co2 and o3.
2518     c (2) units of spre and stem are, respectively, mb and K.
2519     c (3) there are nt identical copies of the tables (coef1, coef2, and
2520     c coef3). the prupose of using the multiple copies of tables is
2521     c to increase the speed in parallel (vectorized) computations.
2522     C if such advantage does not exist, nt can be set to 1.
2523     c
2524     c**********************************************************************
2525     implicit none
2526     integer k1,k2,m,n,np,nx,nh,nt,i,j
2527    
2528     c---- input parameters -----
2529    
2530     real w1,p1,dwe,dpe
2531     real sabs(m,n,np+1),spre(m,n,np+1),stem(m,n,np+1)
2532     real coef1(nx,nh,nt),coef2(nx,nh,nt),coef3(nx,nh,nt)
2533    
2534     c---- update parameter -----
2535    
2536     real tran(m,n)
2537    
2538     c---- temporary variables -----
2539    
2540     real x1,x2,x3,we,pe,fw,fp,pa,pb,pc,ax,ba,bb,t1,ca,cb,t2
2541     integer iw,ip,nn
2542    
2543     c**********************************************************************
2544    
2545     do j=1,n
2546     do i=1,m
2547    
2548     nn=mod(i,nt)+1
2549    
2550     x1=sabs(i,j,k2)-sabs(i,j,k1)
2551     x2=(spre(i,j,k2)-spre(i,j,k1))/x1
2552     x3=(stem(i,j,k2)-stem(i,j,k1))/x1
2553    
2554     we=(log10(x1)-w1)/dwe
2555     pe=(log10(x2)-p1)/dpe
2556    
2557     we=max(we,w1-2.*dwe)
2558     pe=max(pe,p1)
2559    
2560     iw=int(we+1.5)
2561     ip=int(pe+1.5)
2562    
2563     iw=min(iw,nh-1)
2564     iw=max(iw, 2)
2565    
2566     ip=min(ip,nx-1)
2567     ip=max(ip, 1)
2568    
2569     fw=we-float(iw-1)
2570     fp=pe-float(ip-1)
2571    
2572     c-----linear interpolation in pressure
2573    
2574     pa = coef1(ip,iw-1,nn)*(1.-fp)+coef1(ip+1,iw-1,nn)*fp
2575     pb = coef1(ip,iw, nn)*(1.-fp)+coef1(ip+1,iw, nn)*fp
2576     pc = coef1(ip,iw+1,nn)*(1.-fp)+coef1(ip+1,iw+1,nn)*fp
2577    
2578     c-----quadratic interpolation in absorber amount for coef1
2579    
2580     ax = (-pa*(1.-fw)+pc*(1.+fw)) *fw*0.5 + pb*(1.-fw*fw)
2581    
2582     c-----linear interpolation in absorber amount for coef2 and coef3
2583    
2584     ba = coef2(ip,iw, nn)*(1.-fp)+coef2(ip+1,iw, nn)*fp
2585     bb = coef2(ip,iw+1,nn)*(1.-fp)+coef2(ip+1,iw+1,nn)*fp
2586     t1 = ba*(1.-fw) + bb*fw
2587    
2588     ca = coef3(ip,iw, nn)*(1.-fp)+coef3(ip+1,iw, nn)*fp
2589     cb = coef3(ip,iw+1,nn)*(1.-fp)+coef3(ip+1,iw+1,nn)*fp
2590     t2 = ca*(1.-fw) + cb*fw
2591    
2592     c-----update the total transmittance between levels k1 and k2
2593    
2594     tran(i,j)= (ax + (t1+t2*x3) * x3)*tran(i,j)
2595    
2596     enddo
2597     enddo
2598    
2599     return
2600     end
2601     c**********************************************************************
2602     subroutine h2okdis(ib,m,n,np,k,fkw,gkw,ne,h2oexp,conexp,
2603     * th2o,tcon,tran)
2604     c**********************************************************************
2605     c compute water vapor transmittance between levels k1 and k2 for
2606     c m x n soundings, using the k-distribution method.
2607     c
2608     c computations follow eqs. (34), (46), (50) and (52).
2609     c
2610     c---- input parameters
2611     c spectral band (ib)
2612     c number of grid intervals in zonal direction (m)
2613     c number of grid intervals in meridional direction (n)
2614     c number of levels (np)
2615     c current level (k)
2616     c planck-weighted k-distribution function due to
2617     c h2o line absorption (fkw)
2618     c planck-weighted k-distribution function due to
2619     c h2o continuum absorption (gkw)
2620     c number of terms used in each band to compute water vapor
2621     c continuum transmittance (ne)
2622     c exponentials for line absorption (h2oexp)
2623     c exponentials for continuum absorption (conexp)
2624     c
2625     c---- updated parameters
2626     c transmittance between levels k1 and k2 due to
2627     c water vapor line absorption (th2o)
2628     c transmittance between levels k1 and k2 due to
2629     c water vapor continuum absorption (tcon)
2630     c total transmittance (tran)
2631     c
2632     c**********************************************************************
2633     implicit none
2634     integer ib,m,n,np,k,i,j
2635    
2636     c---- input parameters ------
2637    
2638     real conexp(m,n,np,3),h2oexp(m,n,np,6)
2639     integer ne(9)
2640     real fkw(6,9),gkw(6,3)
2641    
2642     c---- updated parameters -----
2643    
2644     real th2o(m,n,6),tcon(m,n,3),tran(m,n)
2645    
2646     c---- temporary arrays -----
2647    
2648     real trnth2o
2649    
2650     c-----tco2 are the six exp factors between levels k1 and k2
2651     c tran is the updated total transmittance between levels k1 and k2
2652    
2653     c-----th2o is the 6 exp factors between levels k1 and k2 due to
2654     c h2o line absorption.
2655    
2656     c-----tcon is the 3 exp factors between levels k1 and k2 due to
2657     c h2o continuum absorption.
2658    
2659     c-----trnth2o is the total transmittance between levels k1 and k2 due
2660     c to both line and continuum absorption computed from eq. (52).
2661    
2662     do j=1,n
2663     do i=1,m
2664     th2o(i,j,1) = th2o(i,j,1)*h2oexp(i,j,k,1)
2665     th2o(i,j,2) = th2o(i,j,2)*h2oexp(i,j,k,2)
2666     th2o(i,j,3) = th2o(i,j,3)*h2oexp(i,j,k,3)
2667     th2o(i,j,4) = th2o(i,j,4)*h2oexp(i,j,k,4)
2668     th2o(i,j,5) = th2o(i,j,5)*h2oexp(i,j,k,5)
2669     th2o(i,j,6) = th2o(i,j,6)*h2oexp(i,j,k,6)
2670     enddo
2671     enddo
2672    
2673     if (ne(ib).eq.0) then
2674    
2675     do j=1,n
2676     do i=1,m
2677    
2678     trnth2o =(fkw(1,ib)*th2o(i,j,1)
2679     * + fkw(2,ib)*th2o(i,j,2)
2680     * + fkw(3,ib)*th2o(i,j,3)
2681     * + fkw(4,ib)*th2o(i,j,4)
2682     * + fkw(5,ib)*th2o(i,j,5)
2683     * + fkw(6,ib)*th2o(i,j,6))
2684    
2685     tran(i,j)=tran(i,j)*trnth2o
2686    
2687     enddo
2688     enddo
2689    
2690     elseif (ne(ib).eq.1) then
2691    
2692     do j=1,n
2693     do i=1,m
2694    
2695     tcon(i,j,1)= tcon(i,j,1)*conexp(i,j,k,1)
2696    
2697     trnth2o =(fkw(1,ib)*th2o(i,j,1)
2698     * + fkw(2,ib)*th2o(i,j,2)
2699     * + fkw(3,ib)*th2o(i,j,3)
2700     * + fkw(4,ib)*th2o(i,j,4)
2701     * + fkw(5,ib)*th2o(i,j,5)
2702     * + fkw(6,ib)*th2o(i,j,6))*tcon(i,j,1)
2703    
2704     tran(i,j)=tran(i,j)*trnth2o
2705    
2706     enddo
2707     enddo
2708    
2709     else
2710    
2711     do j=1,n
2712     do i=1,m
2713    
2714     tcon(i,j,1)= tcon(i,j,1)*conexp(i,j,k,1)
2715     tcon(i,j,2)= tcon(i,j,2)*conexp(i,j,k,2)
2716     tcon(i,j,3)= tcon(i,j,3)*conexp(i,j,k,3)
2717    
2718    
2719     trnth2o = ( gkw(1,1)*th2o(i,j,1)
2720     * + gkw(2,1)*th2o(i,j,2)
2721     * + gkw(3,1)*th2o(i,j,3)
2722     * + gkw(4,1)*th2o(i,j,4)
2723     * + gkw(5,1)*th2o(i,j,5)
2724     * + gkw(6,1)*th2o(i,j,6) ) * tcon(i,j,1)
2725     * + ( gkw(1,2)*th2o(i,j,1)
2726     * + gkw(2,2)*th2o(i,j,2)
2727     * + gkw(3,2)*th2o(i,j,3)
2728     * + gkw(4,2)*th2o(i,j,4)
2729     * + gkw(5,2)*th2o(i,j,5)
2730     * + gkw(6,2)*th2o(i,j,6) ) * tcon(i,j,2)
2731     * + ( gkw(1,3)*th2o(i,j,1)
2732     * + gkw(2,3)*th2o(i,j,2)
2733     * + gkw(3,3)*th2o(i,j,3)
2734     * + gkw(4,3)*th2o(i,j,4)
2735     * + gkw(5,3)*th2o(i,j,5)
2736     * + gkw(6,3)*th2o(i,j,6) ) * tcon(i,j,3)
2737    
2738     tran(i,j)=tran(i,j)*trnth2o
2739    
2740     enddo
2741     enddo
2742    
2743     endif
2744    
2745     return
2746     end
2747     c**********************************************************************
2748     subroutine co2kdis(m,n,np,k,co2exp,tco2,tran)
2749     c**********************************************************************
2750     c compute co2 transmittances between levels k1 and k2 for
2751     c m x n soundings, using the k-distribution method with linear
2752     c pressure scaling. computations follow eq. (34).
2753     c
2754     c---- input parameters
2755     c number of grid intervals in zonal direction (m)
2756     c number of grid intervals in meridional direction (n)
2757     c number of levels (np)
2758     c current level (k)
2759     c exponentials for co2 absorption (co2exp)
2760     c
2761     c---- updated parameters
2762     c transmittance between levels k1 and k2 due to co2 absorption
2763     c for the various values of the absorption coefficient (tco2)
2764     c total transmittance (tran)
2765     c
2766     c**********************************************************************
2767     implicit none
2768     integer m,n,np,k,i,j
2769    
2770     c---- input parameters -----
2771    
2772     real co2exp(m,n,np,6,2)
2773    
2774     c---- updated parameters -----
2775    
2776     real tco2(m,n,6,2),tran(m,n)
2777    
2778     c---- temporary arrays -----
2779    
2780     real xc
2781    
2782     c-----tco2 is the 6 exp factors between levels k1 and k2.
2783     c xc is the total co2 transmittance given by eq. (53).
2784    
2785     do j=1,n
2786     do i=1,m
2787    
2788     c-----band-wings
2789    
2790     tco2(i,j,1,1)=tco2(i,j,1,1)*co2exp(i,j,k,1,1)
2791     xc= 0.1395 *tco2(i,j,1,1)
2792    
2793     tco2(i,j,2,1)=tco2(i,j,2,1)*co2exp(i,j,k,2,1)
2794     xc=xc+0.1407 *tco2(i,j,2,1)
2795    
2796     tco2(i,j,3,1)=tco2(i,j,3,1)*co2exp(i,j,k,3,1)
2797     xc=xc+0.1549 *tco2(i,j,3,1)
2798    
2799     tco2(i,j,4,1)=tco2(i,j,4,1)*co2exp(i,j,k,4,1)
2800     xc=xc+0.1357 *tco2(i,j,4,1)
2801    
2802     tco2(i,j,5,1)=tco2(i,j,5,1)*co2exp(i,j,k,5,1)
2803     xc=xc+0.0182 *tco2(i,j,5,1)
2804    
2805     tco2(i,j,6,1)=tco2(i,j,6,1)*co2exp(i,j,k,6,1)
2806     xc=xc+0.0220 *tco2(i,j,6,1)
2807    
2808     c-----band-center region
2809    
2810     tco2(i,j,1,2)=tco2(i,j,1,2)*co2exp(i,j,k,1,2)
2811     xc=xc+0.0766 *tco2(i,j,1,2)
2812    
2813     tco2(i,j,2,2)=tco2(i,j,2,2)*co2exp(i,j,k,2,2)
2814     xc=xc+0.1372 *tco2(i,j,2,2)
2815    
2816     tco2(i,j,3,2)=tco2(i,j,3,2)*co2exp(i,j,k,3,2)
2817     xc=xc+0.1189 *tco2(i,j,3,2)
2818    
2819     tco2(i,j,4,2)=tco2(i,j,4,2)*co2exp(i,j,k,4,2)
2820     xc=xc+0.0335 *tco2(i,j,4,2)
2821    
2822     tco2(i,j,5,2)=tco2(i,j,5,2)*co2exp(i,j,k,5,2)
2823     xc=xc+0.0169 *tco2(i,j,5,2)
2824    
2825     tco2(i,j,6,2)=tco2(i,j,6,2)*co2exp(i,j,k,6,2)
2826     xc=xc+0.0059 *tco2(i,j,6,2)
2827    
2828     tran(i,j)=tran(i,j)*xc
2829    
2830     enddo
2831     enddo
2832    
2833     return
2834     end
2835     c**********************************************************************
2836     subroutine n2okdis(ib,m,n,np,k,n2oexp,tn2o,tran)
2837     c**********************************************************************
2838     c compute n2o transmittances between levels k1 and k2 for
2839     c m x n soundings, using the k-distribution method with linear
2840     c pressure scaling.
2841     c
2842     c---- input parameters
2843     c spectral band (ib)
2844     c number of grid intervals in zonal direction (m)
2845     c number of grid intervals in meridional direction (n)
2846     c number of levels (np)
2847     c current level (k)
2848     c exponentials for n2o absorption (n2oexp)
2849     c
2850     c---- updated parameters
2851     c transmittance between levels k1 and k2 due to n2o absorption
2852     c for the various values of the absorption coefficient (tn2o)
2853     c total transmittance (tran)
2854     c
2855     c**********************************************************************
2856     implicit none
2857     integer ib,m,n,np,k,i,j
2858    
2859     c---- input parameters -----
2860    
2861     real n2oexp(m,n,np,4)
2862    
2863     c---- updated parameters -----
2864    
2865     real tn2o(m,n,4),tran(m,n)
2866    
2867     c---- temporary arrays -----
2868    
2869     real xc
2870    
2871     c-----tn2o is the 2 exp factors between levels k1 and k2.
2872     c xc is the total n2o transmittance
2873    
2874     do j=1,n
2875     do i=1,m
2876    
2877     c-----band 6
2878    
2879     if (ib.eq.6) then
2880    
2881     tn2o(i,j,1)=tn2o(i,j,1)*n2oexp(i,j,k,1)
2882     xc= 0.940414*tn2o(i,j,1)
2883    
2884     tn2o(i,j,2)=tn2o(i,j,2)*n2oexp(i,j,k,2)
2885     xc=xc+0.059586*tn2o(i,j,2)
2886    
2887     c-----band 7
2888    
2889     else
2890    
2891     tn2o(i,j,1)=tn2o(i,j,1)*n2oexp(i,j,k,1)
2892     xc= 0.561961*tn2o(i,j,1)
2893    
2894     tn2o(i,j,2)=tn2o(i,j,2)*n2oexp(i,j,k,2)
2895     xc=xc+0.138707*tn2o(i,j,2)
2896    
2897     tn2o(i,j,3)=tn2o(i,j,3)*n2oexp(i,j,k,3)
2898     xc=xc+0.240670*tn2o(i,j,3)
2899    
2900     tn2o(i,j,4)=tn2o(i,j,4)*n2oexp(i,j,k,4)
2901     xc=xc+0.058662*tn2o(i,j,4)
2902    
2903     endif
2904    
2905     tran(i,j)=tran(i,j)*xc
2906    
2907     enddo
2908     enddo
2909    
2910     return
2911     end
2912     c**********************************************************************
2913     subroutine ch4kdis(ib,m,n,np,k,ch4exp,tch4,tran)
2914     c**********************************************************************
2915     c compute ch4 transmittances between levels k1 and k2 for
2916     c m x n soundings, using the k-distribution method with
2917     c linear pressure scaling.
2918     c
2919     c---- input parameters
2920     c spectral band (ib)
2921     c number of grid intervals in zonal direction (m)
2922     c number of grid intervals in meridional direction (n)
2923     c number of levels (np)
2924     c current level (k)
2925     c exponentials for ch4 absorption (ch4exp)
2926     c
2927     c---- updated parameters
2928     c transmittance between levels k1 and k2 due to ch4 absorption
2929     c for the various values of the absorption coefficient (tch4)
2930     c total transmittance (tran)
2931     c
2932     c**********************************************************************
2933     implicit none
2934     integer ib,m,n,np,k,i,j
2935    
2936     c---- input parameters -----
2937    
2938     real ch4exp(m,n,np,4)
2939    
2940     c---- updated parameters -----
2941    
2942     real tch4(m,n,4),tran(m,n)
2943    
2944     c---- temporary arrays -----
2945    
2946     real xc
2947    
2948     c-----tch4 is the 2 exp factors between levels k1 and k2.
2949     c xc is the total ch4 transmittance
2950    
2951     do j=1,n
2952     do i=1,m
2953    
2954     c-----band 6
2955    
2956     if (ib.eq.6) then
2957    
2958     tch4(i,j,1)=tch4(i,j,1)*ch4exp(i,j,k,1)
2959     xc= tch4(i,j,1)
2960    
2961     c-----band 7
2962    
2963     else
2964    
2965     tch4(i,j,1)=tch4(i,j,1)*ch4exp(i,j,k,1)
2966     xc= 0.610650*tch4(i,j,1)
2967    
2968     tch4(i,j,2)=tch4(i,j,2)*ch4exp(i,j,k,2)
2969     xc=xc+0.280212*tch4(i,j,2)
2970    
2971     tch4(i,j,3)=tch4(i,j,3)*ch4exp(i,j,k,3)
2972     xc=xc+0.107349*tch4(i,j,3)
2973    
2974     tch4(i,j,4)=tch4(i,j,4)*ch4exp(i,j,k,4)
2975     xc=xc+0.001789*tch4(i,j,4)
2976    
2977     endif
2978    
2979     tran(i,j)=tran(i,j)*xc
2980    
2981     enddo
2982     enddo
2983    
2984     return
2985     end
2986     c**********************************************************************
2987     subroutine comkdis(ib,m,n,np,k,comexp,tcom,tran)
2988     c**********************************************************************
2989     c compute co2-minor transmittances between levels k1 and k2
2990     c for m x n soundings, using the k-distribution method
2991     c with linear pressure scaling.
2992     c
2993     c---- input parameters
2994     c spectral band (ib)
2995     c number of grid intervals in zonal direction (m)
2996     c number of grid intervals in meridional direction (n)
2997     c number of levels (np)
2998     c current level (k)
2999     c exponentials for co2-minor absorption (comexp)
3000     c
3001     c---- updated parameters
3002     c transmittance between levels k1 and k2 due to co2-minor absorption
3003     c for the various values of the absorption coefficient (tcom)
3004     c total transmittance (tran)
3005     c
3006     c**********************************************************************
3007     implicit none
3008     integer ib,m,n,np,k,i,j
3009    
3010     c---- input parameters -----
3011    
3012     real comexp(m,n,np,2)
3013    
3014     c---- updated parameters -----
3015    
3016     real tcom(m,n,2),tran(m,n)
3017    
3018     c---- temporary arrays -----
3019    
3020     real xc
3021    
3022     c-----tcom is the 2 exp factors between levels k1 and k2.
3023     c xc is the total co2-minor transmittance
3024    
3025     do j=1,n
3026     do i=1,m
3027    
3028     c-----band 4
3029    
3030     if (ib.eq.4) then
3031    
3032     tcom(i,j,1)=tcom(i,j,1)*comexp(i,j,k,1)
3033     xc= 0.972025*tcom(i,j,1)
3034     tcom(i,j,2)=tcom(i,j,2)*comexp(i,j,k,2)
3035     xc=xc+0.027975*tcom(i,j,2)
3036    
3037     c-----band 5
3038    
3039     else
3040    
3041     tcom(i,j,1)=tcom(i,j,1)*comexp(i,j,k,1)
3042     xc= 0.961324*tcom(i,j,1)
3043     tcom(i,j,2)=tcom(i,j,2)*comexp(i,j,k,2)
3044     xc=xc+0.038676*tcom(i,j,2)
3045    
3046     endif
3047    
3048     tran(i,j)=tran(i,j)*xc
3049    
3050     enddo
3051     enddo
3052    
3053     return
3054     end
3055     c**********************************************************************
3056     subroutine cfckdis(m,n,np,k,cfcexp,tcfc,tran)
3057     c**********************************************************************
3058     c compute cfc-(11,12,22) transmittances between levels k1 and k2
3059     c for m x n soundings, using the k-distribution method with
3060     c linear pressure scaling.
3061     c
3062     c---- input parameters
3063     c number of grid intervals in zonal direction (m)
3064     c number of grid intervals in meridional direction (n)
3065     c number of levels (np)
3066     c current level (k)
3067     c exponentials for cfc absorption (cfcexp)
3068     c
3069     c---- updated parameters
3070     c transmittance between levels k1 and k2 due to cfc absorption
3071     c for the various values of the absorption coefficient (tcfc)
3072     c total transmittance (tran)
3073     c
3074     c**********************************************************************
3075     implicit none
3076     integer m,n,np,k,i,j
3077    
3078     c---- input parameters -----
3079    
3080     real cfcexp(m,n,np)
3081    
3082     c---- updated parameters -----
3083    
3084     real tcfc(m,n),tran(m,n)
3085    
3086     c-----tcfc is the exp factors between levels k1 and k2.
3087    
3088     do j=1,n
3089     do i=1,m
3090    
3091     tcfc(i,j)=tcfc(i,j)*cfcexp(i,j,k)
3092     tran(i,j)=tran(i,j)*tcfc(i,j)
3093    
3094     enddo
3095     enddo
3096    
3097     return
3098     end
3099     c**********************************************************************
3100     subroutine b10kdis(m,n,np,k,h2oexp,conexp,co2exp,n2oexp
3101     * ,th2o,tcon,tco2,tn2o,tran)
3102     c**********************************************************************
3103     c
3104     c compute h2o (line and continuum),co2,n2o transmittances between
3105     c levels k1 and k2 for m x n soundings, using the k-distribution
3106     c method with linear pressure scaling.
3107     c
3108     c---- input parameters
3109     c number of grid intervals in zonal direction (m)
3110     c number of grid intervals in meridional direction (n)
3111     c number of levels (np)
3112     c current level (k)
3113     c exponentials for h2o line absorption (h2oexp)
3114     c exponentials for h2o continuum absorption (conexp)
3115     c exponentials for co2 absorption (co2exp)
3116     c exponentials for n2o absorption (n2oexp)
3117     c
3118     c---- updated parameters
3119     c transmittance between levels k1 and k2 due to h2o line absorption
3120     c for the various values of the absorption coefficient (th2o)
3121     c transmittance between levels k1 and k2 due to h2o continuum
3122     c absorption for the various values of the absorption
3123     c coefficient (tcon)
3124     c transmittance between levels k1 and k2 due to co2 absorption
3125     c for the various values of the absorption coefficient (tco2)
3126     c transmittance between levels k1 and k2 due to n2o absorption
3127     c for the various values of the absorption coefficient (tn2o)
3128     c total transmittance (tran)
3129     c
3130     c**********************************************************************
3131     implicit none
3132     integer m,n,np,k,i,j
3133    
3134     c---- input parameters -----
3135    
3136     real h2oexp(m,n,np,6),conexp(m,n,np,3),co2exp(m,n,np,6,2)
3137     * ,n2oexp(m,n,np,4)
3138    
3139     c---- updated parameters -----
3140    
3141     real th2o(m,n,6),tcon(m,n,3),tco2(m,n,6,2),tn2o(m,n,4)
3142     * ,tran(m,n)
3143    
3144     c---- temporary arrays -----
3145    
3146     real xx
3147    
3148     c-----initialize tran
3149    
3150     do j=1,n
3151     do i=1,m
3152     tran(i,j)=1.0
3153     enddo
3154     enddo
3155    
3156     c-----for h2o line
3157    
3158     do j=1,n
3159     do i=1,m
3160    
3161     th2o(i,j,1)=th2o(i,j,1)*h2oexp(i,j,k,1)
3162     xx= 0.3153*th2o(i,j,1)
3163    
3164     th2o(i,j,2)=th2o(i,j,2)*h2oexp(i,j,k,2)
3165     xx=xx+0.4604*th2o(i,j,2)
3166    
3167     th2o(i,j,3)=th2o(i,j,3)*h2oexp(i,j,k,3)
3168     xx=xx+0.1326*th2o(i,j,3)
3169    
3170     th2o(i,j,4)=th2o(i,j,4)*h2oexp(i,j,k,4)
3171     xx=xx+0.0798*th2o(i,j,4)
3172    
3173     th2o(i,j,5)=th2o(i,j,5)*h2oexp(i,j,k,5)
3174     xx=xx+0.0119*th2o(i,j,5)
3175    
3176     tran(i,j)=tran(i,j)*xx
3177    
3178     enddo
3179     enddo
3180    
3181     c-----for h2o continuum
3182    
3183     do j=1,n
3184     do i=1,m
3185    
3186     tcon(i,j,1)=tcon(i,j,1)*conexp(i,j,k,1)
3187     tran(i,j)=tran(i,j)*tcon(i,j,1)
3188    
3189     enddo
3190     enddo
3191    
3192     c-----for co2
3193    
3194     do j=1,n
3195     do i=1,m
3196    
3197     tco2(i,j,1,1)=tco2(i,j,1,1)*co2exp(i,j,k,1,1)
3198     xx= 0.2673*tco2(i,j,1,1)
3199    
3200     tco2(i,j,2,1)=tco2(i,j,2,1)*co2exp(i,j,k,2,1)
3201     xx=xx+ 0.2201*tco2(i,j,2,1)
3202    
3203     tco2(i,j,3,1)=tco2(i,j,3,1)*co2exp(i,j,k,3,1)
3204     xx=xx+ 0.2106*tco2(i,j,3,1)
3205    
3206     tco2(i,j,4,1)=tco2(i,j,4,1)*co2exp(i,j,k,4,1)
3207     xx=xx+ 0.2409*tco2(i,j,4,1)
3208    
3209     tco2(i,j,5,1)=tco2(i,j,5,1)*co2exp(i,j,k,5,1)
3210     xx=xx+ 0.0196*tco2(i,j,5,1)
3211    
3212     tco2(i,j,6,1)=tco2(i,j,6,1)*co2exp(i,j,k,6,1)
3213     xx=xx+ 0.0415*tco2(i,j,6,1)
3214    
3215     tran(i,j)=tran(i,j)*xx
3216    
3217     enddo
3218     enddo
3219    
3220     c-----for n2o
3221    
3222     do j=1,n
3223     do i=1,m
3224    
3225     tn2o(i,j,1)=tn2o(i,j,1)*n2oexp(i,j,k,1)
3226     xx= 0.970831*tn2o(i,j,1)
3227    
3228     tn2o(i,j,2)=tn2o(i,j,2)*n2oexp(i,j,k,2)
3229     xx=xx+0.029169*tn2o(i,j,2)
3230    
3231     tran(i,j)=tran(i,j)*(xx-1.0)
3232    
3233     enddo
3234     enddo
3235    
3236     return
3237     end

  ViewVC Help
Powered by ViewVC 1.1.22