/[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.2 - (hide annotations) (download)
Tue Jun 15 16:06:03 2004 UTC (20 years ago) by molod
Branch: MAIN
Changes since 1.1: +6 -5 lines
Include cpp options and diagnostics common if needed

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

  ViewVC Help
Powered by ViewVC 1.1.22