/[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.4 - (show annotations) (download)
Wed Jul 7 19:33:48 2004 UTC (20 years ago) by molod
Branch: MAIN
CVS Tags: checkpoint54b_post
Changes since 1.3: +3 -13 lines
Almost last code to get rid of references to sigma coordinate

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

  ViewVC Help
Powered by ViewVC 1.1.22