/[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.7 - (hide annotations) (download)
Tue Jul 13 21:26:32 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.6: +4 -4 lines
Reflect new transfer coeff file names

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

  ViewVC Help
Powered by ViewVC 1.1.22