/[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.9 - (hide annotations) (download)
Wed Jul 14 17:31:57 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint54c_post
Changes since 1.8: +2 -2 lines
Debugging

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

  ViewVC Help
Powered by ViewVC 1.1.22