/[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.24 - (show annotations) (download)
Thu Jun 16 16:46:12 2005 UTC (19 years ago) by ce107
Branch: MAIN
Changes since 1.23: +9 -9 lines
Fix fizhi constants to _d 0 form as the IBM XL compiler complains on mixed
precision calls to intrinsics like max, min. Only problematic cases have
been altered - consider compiling with -qdpc=e to fix the precision of the
rest.

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

  ViewVC Help
Powered by ViewVC 1.1.22