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

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

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


Revision 1.5 - (show annotations) (download)
Mon Jul 12 21:33:37 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.4: +51 -82 lines
Corrections for argument list match and proper declarations

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

  ViewVC Help
Powered by ViewVC 1.1.22