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

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

  ViewVC Help
Powered by ViewVC 1.1.22