Radiation Scheme in CCPP
radlw_main.f
Go to the documentation of this file.
1 !!!!! ============================================================== !!!!!
2 !!!!! lw-rrtm3 radiation package description !!!!!
3 !!!!! ============================================================== !!!!!
4 ! !
5 ! this package includes ncep's modifications of the rrtm-lw radiation !
6 ! code from aer inc. !
7 ! !
8 ! the lw-rrtm3 package includes these parts: !
9 ! !
10 ! 'radlw_rrtm3_param.f' !
11 ! 'radlw_rrtm3_datatb.f' !
12 ! 'radlw_rrtm3_main.f' !
13 ! !
14 ! the 'radlw_rrtm3_param.f' contains: !
15 ! !
16 ! 'module_radlw_parameters' -- band parameters set up !
17 ! !
18 ! the 'radlw_rrtm3_datatb.f' contains: !
19 ! !
20 ! 'module_radlw_avplank' -- plank flux data !
21 ! 'module_radlw_ref' -- reference temperature and pressure !
22 ! 'module_radlw_cldprlw' -- cloud property coefficients !
23 ! 'module_radlw_kgbnn' -- absorption coeffients for 16 !
24 ! bands, where nn = 01-16 !
25 ! !
26 ! the 'radlw_rrtm3_main.f' contains: !
27 ! !
28 ! 'module_radlw_main' -- main lw radiation transfer !
29 ! !
30 ! in the main module 'module_radlw_main' there are only two !
31 ! externally callable subroutines: !
32 ! !
33 ! !
34 ! 'lwrad' -- main lw radiation routine !
35 ! inputs: !
36 ! (plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, !
37 ! clouds,icseed,aerosols,sfemis,sfgtmp, !
38 ! npts, nlay, nlp1, lprnt, !
39 ! outputs: !
40 ! hlwc,topflx,sfcflx, !
41 !! optional outputs: !
42 ! HLW0,HLWB,FLXPRF) !
43 ! !
44 ! 'rlwinit' -- initialization routine !
45 ! inputs: !
46 ! ( me ) !
47 ! outputs: !
48 ! (none) !
49 ! !
50 ! all the lw radiation subprograms become contained subprograms !
51 ! in module 'module_radlw_main' and many of them are not directly !
52 ! accessable from places outside the module. !
53 ! !
54 ! derived data type constructs used: !
55 ! !
56 ! 1. radiation flux at toa: (from module 'module_radlw_parameters') !
57 ! topflw_type - derived data type for toa rad fluxes !
58 ! upfxc total sky upward flux at toa !
59 ! upfx0 clear sky upward flux at toa !
60 ! !
61 ! 2. radiation flux at sfc: (from module 'module_radlw_parameters') !
62 ! sfcflw_type - derived data type for sfc rad fluxes !
63 ! upfxc total sky upward flux at sfc !
64 ! upfx0 clear sky upward flux at sfc !
65 ! dnfxc total sky downward flux at sfc !
66 ! dnfx0 clear sky downward flux at sfc !
67 ! !
68 ! 3. radiation flux profiles(from module 'module_radlw_parameters') !
69 ! proflw_type - derived data type for rad vertical prof !
70 ! upfxc level upward flux for total sky !
71 ! dnfxc level downward flux for total sky !
72 ! upfx0 level upward flux for clear sky !
73 ! dnfx0 level downward flux for clear sky !
74 ! !
75 ! external modules referenced: !
76 ! !
77 ! 'module physparam' !
78 ! 'module physcons' !
79 ! 'mersenne_twister' !
80 ! !
81 ! compilation sequence is: !
82 ! !
83 ! 'radlw_rrtm3_param.f' !
84 ! 'radlw_rrtm3_datatb.f' !
85 ! 'radlw_rrtm3_main.f' !
86 ! !
87 ! and all should be put in front of routines that use lw modules !
88 ! !
89 !==========================================================================!
90 ! !
91 ! the original aer's program declarations: !
92 ! !
93 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
94 ! |
95 ! Copyright 2002-2007, Atmospheric & Environmental Research, Inc. (AER). |
96 ! This software may be used, copied, or redistributed as long as it is |
97 ! not sold and this copyright notice is reproduced on each copy made. |
98 ! This model is provided as is without any express or implied warranties. |
99 ! (http://www.rtweb.aer.com/) |
100 ! |
101 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
102 ! !
103 ! ************************************************************************ !
104 ! !
105 ! rrtmg_lw !
106 ! !
107 ! !
108 ! a rapid radiative transfer model !
109 ! for the longwave region !
110 ! for application to general circulation models !
111 ! !
112 ! !
113 ! atmospheric and environmental research, inc. !
114 ! 131 hartwell avenue !
115 ! lexington, ma 02421 !
116 ! !
117 ! eli j. mlawer !
118 ! jennifer s. delamere !
119 ! michael j. iacono !
120 ! shepard a. clough !
121 ! !
122 ! !
123 ! email: miacono@aer.com !
124 ! email: emlawer@aer.com !
125 ! email: jdelamer@aer.com !
126 ! !
127 ! the authors wish to acknowledge the contributions of the !
128 ! following people: steven j. taubman, karen cady-pereira, !
129 ! patrick d. brown, ronald e. farren, luke chen, robert bergstrom. !
130 ! !
131 ! ************************************************************************ !
132 ! !
133 ! references: !
134 ! (rrtm_lw/rrtmg_lw): !
135 ! clough, s.A., m.w. shephard, e.j. mlawer, j.s. delamere, !
136 ! m.j. iacono, k. cady-pereira, s. boukabara, and p.d. brown: !
137 ! atmospheric radiative transfer modeling: a summary of the aer !
138 ! codes, j. quant. spectrosc. radiat. transfer, 91, 233-244, 2005. !
139 ! !
140 ! mlawer, e.j., s.j. taubman, p.d. brown, m.j. iacono, and s.a. !
141 ! clough: radiative transfer for inhomogeneous atmospheres: rrtm, !
142 ! a validated correlated-k model for the longwave. j. geophys. res., !
143 ! 102, 16663-16682, 1997. !
144 ! !
145 ! (mcica): !
146 ! pincus, r., h. w. barker, and j.-j. morcrette: a fast, flexible, !
147 ! approximation technique for computing radiative transfer in !
148 ! inhomogeneous cloud fields, j. geophys. res., 108(d13), 4376, !
149 ! doi:10.1029/2002JD003322, 2003. !
150 ! !
151 ! ************************************************************************ !
152 ! !
153 ! aer's revision history: !
154 ! this version of rrtmg_lw has been modified from rrtm_lw to use a !
155 ! reduced set of g-points for application to gcms. !
156 ! !
157 ! -- original version (derived from rrtm_lw), reduction of g-points, !
158 ! other revisions for use with gcms. !
159 ! 1999: m. j. iacono, aer, inc. !
160 ! -- adapted for use with ncar/cam3. !
161 ! may 2004: m. j. iacono, aer, inc. !
162 ! -- revised to add mcica capability. !
163 ! nov 2005: m. j. iacono, aer, inc. !
164 ! -- conversion to f90 formatting for consistency with rrtmg_sw. !
165 ! feb 2007: m. j. iacono, aer, inc. !
166 ! -- modifications to formatting to use assumed-shape arrays. !
167 ! aug 2007: m. j. iacono, aer, inc. !
168 ! !
169 ! ************************************************************************ !
170 ! !
171 ! ncep modifications history log: !
172 ! !
173 ! nov 1999, ken campana -- received the original code from !
174 ! aer (1998 ncar ccm version), updated to link up with !
175 ! ncep mrf model !
176 ! jun 2000, ken campana -- added option to switch random and !
177 ! maximum/random cloud overlap !
178 ! 2001, shrinivas moorthi -- further updates for mrf model !
179 ! may 2001, yu-tai hou -- updated on trace gases and cloud !
180 ! property based on rrtm_v3.0 codes. !
181 ! dec 2001, yu-tai hou -- rewritten code into fortran 90 std !
182 ! set ncep radiation structure standard that contains !
183 ! three plug-in compatable fortran program files: !
184 ! 'radlw_param.f', 'radlw_datatb.f', 'radlw_main.f' !
185 ! fixed bugs in subprograms taugb14, taugb2, etc. added !
186 ! out-of-bounds protections. (a detailed note of !
187 ! up_to_date modifications/corrections by ncep was sent !
188 ! to aer in 2002) !
189 ! jun 2004, yu-tai hou -- added mike iacono's apr 2004 !
190 ! modification of variable diffusivity angles. !
191 ! apr 2005, yu-tai hou -- minor modifications on module !
192 ! structures include rain/snow effect (this version of !
193 ! code was given back to aer in jun 2006) !
194 ! mar 2007, yu-tai hou -- added aerosol effect for ncep !
195 ! models using the generallized aerosol optical property!
196 ! scheme for gfs model. !
197 ! apr 2007, yu-tai hou -- added spectral band heating as an !
198 ! optional output to support the 500 km gfs model's !
199 ! upper stratospheric radiation calculations. and !
200 ! restructure optional outputs for easy access by !
201 ! different models. !
202 ! oct 2008, yu-tai hou -- modified to include new features !
203 ! from aer's newer release v4.4-v4.7, including the !
204 ! mcica sub-grid cloud option. add rain/snow optical !
205 ! properties support to cloudy sky calculations. !
206 ! correct errors in mcica cloud optical properties for !
207 ! ebert & curry scheme (ilwcice=1) that needs band !
208 ! index conversion. simplified and unified sw and lw !
209 ! sub-column cloud subroutines into one module by using !
210 ! optional parameters. !
211 ! mar 2009, yu-tai hou -- replaced the original random number!
212 ! generator coming from the original code with ncep w3 !
213 ! library to simplify the program and moved sub-column !
214 ! cloud subroutines inside the main module. added !
215 ! option of user provided permutation seeds that could !
216 ! be randomly generated from forecast time stamp. !
217 ! oct 2009, yu-tai hou -- modified subrtines "cldprop" and !
218 ! "rlwinit" according updats from aer's rrtmg_lw v4.8. !
219 ! nov 2009, yu-tai hou -- modified subrtine "taumol" according
220 ! updats from aer's rrtmg_lw version 4.82. notice the !
221 ! cloud ice/liquid are assumed as in-cloud quantities, !
222 ! not as grid averaged quantities. !
223 ! jun 2010, yu-tai hou -- optimized code to improve efficiency
224 ! apr 2012, b. ferrier and y. hou -- added conversion factor to fu's!
225 ! cloud-snow optical property scheme. !
226 ! nov 2012, yu-tai hou -- modified control parameters thru !
227 ! module 'physparam'. !
228 ! !
229 !!!!! ============================================================== !!!!!
230 !!!!! end descriptions !!!!!
231 !!!!! ============================================================== !!!!!
232 
235 !========================================!
237 !........................................!
238 !
239  use physparam, only : ilwrate, ilwrgas, ilwcliq, ilwcice, &
240  & isubclw, icldflg, iovrlw, ivflip, &
241  & kind_phys
242  use physcons, only : con_g, con_cp, con_avgd, con_amd, &
243  & con_amw, con_amo3
244  use mersenne_twister, only : random_setseed, random_number, &
245  & random_stat
246 
248 !
249  use module_radlw_avplank, only : totplnk
250  use module_radlw_ref, only : preflog, tref, chi_mls
251 !
252  implicit none
253 !
254  private
255 !
256 ! ... version tag and last revision date
257  character(40), parameter :: &
258  & VTAGLW='NCEP LW v5.1 Nov 2012 -RRTMG-LW v4.82 '
259 ! & VTAGLW='NCEP LW v5.0 Aug 2012 -RRTMG-LW v4.82 '
260 ! & VTAGLW='RRTMG-LW v4.82 Nov 2009 '
261 ! & VTAGLW='RRTMG-LW v4.8 Oct 2009 '
262 ! & VTAGLW='RRTMG-LW v4.71 Mar 2009 '
263 ! & VTAGLW='RRTMG-LW v4.4 Oct 2008 '
264 ! & VTAGLW='RRTM-LW v2.3g Mar 2007 '
265 ! & VTAGLW='RRTM-LW v2.3g Apr 2004 '
266 
267 ! --- constant values
268  real (kind=kind_phys), parameter :: eps = 1.0e-6
269  real (kind=kind_phys), parameter :: oneminus= 1.0-eps
270  real (kind=kind_phys), parameter :: cldmin = 1.0e-80
271  real (kind=kind_phys), parameter :: bpade = 1.0/0.278 ! pade approx constant
272  real (kind=kind_phys), parameter :: stpfac = 296.0/1013.0
273  real (kind=kind_phys), parameter :: wtdiff = 0.5 ! weight for radiance to flux conversion
274  real (kind=kind_phys), parameter :: tblint = ntbl ! lookup table conversion factor
275  real (kind=kind_phys), parameter :: f_zero = 0.0
276  real (kind=kind_phys), parameter :: f_one = 1.0
277 
278 ! ... atomic weights for conversion from mass to volume mixing ratios
279  real (kind=kind_phys), parameter :: amdw = con_amd/con_amw
280  real (kind=kind_phys), parameter :: amdo3 = con_amd/con_amo3
281 
282 ! ... band indices
283  integer, dimension(nbands) :: nspa, nspb
284 
285  data nspa / 1, 1, 9, 9, 9, 1, 9, 1, 9, 1, 1, 9, 9, 1, 9, 9 /
286  data nspb / 1, 1, 5, 5, 5, 0, 1, 1, 1, 1, 1, 0, 0, 1, 0, 0 /
287 
288 ! ... band wavenumber intervals
289 ! real (kind=kind_phys) :: wavenum1(nbands), wavenum2(nbands)
290 ! data wavenum1/ &
291 ! & 10., 350., 500., 630., 700., 820., 980., 1080., &
292 !err & 1180., 1390., 1480., 1800., 2080., 2250., 2390., 2600. /
293 ! & 1180., 1390., 1480., 1800., 2080., 2250., 2380., 2600. /
294 ! data wavenum2/ &
295 ! & 350., 500., 630., 700., 820., 980., 1080., 1180., &
296 !err & 1390., 1480., 1800., 2080., 2250., 2390., 2600., 3250. /
297 ! & 1390., 1480., 1800., 2080., 2250., 2380., 2600., 3250. /
298 ! real (kind=kind_phys) :: delwave(nbands)
299 ! data delwave / 340., 150., 130., 70., 120., 160., 100., 100., &
300 ! & 210., 90., 320., 280., 170., 130., 220., 650. /
301 
302 ! --- reset diffusivity angle for Bands 2-3 and 5-9 to vary (between 1.50
303 ! and 1.80) as a function of total column water vapor. the function
304 ! has been defined to minimize flux and cooling rate errors in these bands
305 ! over a wide range of precipitable water values.
306  real (kind=kind_phys), dimension(nbands) :: a0, a1, a2
307 
308  data a0 / 1.66, 1.55, 1.58, 1.66, 1.54, 1.454, 1.89, 1.33, &
309  & 1.668, 1.66, 1.66, 1.66, 1.66, 1.66, 1.66, 1.66 /
310  data a1 / 0.00, 0.25, 0.22, 0.00, 0.13, 0.446, -0.10, 0.40, &
311  & -0.006, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00 /
312  data a2 / 0.00, -12.0, -11.7, 0.00, -0.72,-0.243, 0.19,-0.062, &
313  & 0.414, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00 /
314 
315 !! --- logical flags for optional output fields
316 
317  logical :: lhlwb = .false.
318  logical :: lhlw0 = .false.
319  logical :: lflxprf= .false.
320 
321 ! --- those data will be set up only once by "rlwinit"
322 
323 ! ... fluxfac, heatfac are factors for fluxes (in w/m**2) and heating
324 ! rates (in k/day, or k/sec set by subroutine 'rlwinit')
325 ! semiss0 are default surface emissivity for each bands
326 
327  real (kind=kind_phys) :: fluxfac, heatfac, semiss0(nbands)
328  data semiss0(:) / nbands*1.0 /
329 
330  real (kind=kind_phys) :: tau_tbl(0:ntbl) !clr-sky opt dep (for cldy transfer)
331  real (kind=kind_phys) :: exp_tbl(0:ntbl) !transmittance lookup table
332  real (kind=kind_phys) :: tfn_tbl(0:ntbl) !tau transition function; i.e. the
333  !transition of planck func from mean lyr
334  !temp to lyr boundary temp as a func of
335  !opt dep. "linear in tau" method is used.
336 
337 ! --- the following variables are used for sub-column cloud scheme
338 
339  integer, parameter :: ipsdlw0 = ngptlw ! initial permutation seed
340 
341 ! --- public accessable subprograms
342 
343  public lwrad, rlwinit
344 
345 
346 ! ================
347  contains
348 ! ================
349 
412 ! --------------------------------
413  subroutine lwrad &
414 ! --------------------------------
416 ! --- inputs:
417  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
418  & clouds,icseed,aerosols,sfemis,sfgtmp, &
419  & npts, nlay, nlp1, lprnt, &
420 ! --- outputs:
421  & hlwc,topflx,sfcflx &
422 !! --- optional:
423  &, hlw0,hlwb,flxprf &
424  & )
425 
426 ! ==================== defination of variables ==================== !
427 ! !
428 ! input variables: !
429 ! plyr (npts,nlay) : layer mean pressures (mb) !
430 ! plvl (npts,nlp1) : interface pressures (mb) !
431 ! tlyr (npts,nlay) : layer mean temperature (k) !
432 ! tlvl (npts,nlp1) : interface temperatures (k) !
433 ! qlyr (npts,nlay) : layer specific humidity (gm/gm) *see inside !
434 ! olyr (npts,nlay) : layer ozone concentration (gm/gm) *see inside !
435 ! gasvmr(npts,nlay,:): atmospheric gases amount: !
436 ! (check module_radiation_gases for definition) !
437 ! gasvmr(:,:,1) - co2 volume mixing ratio !
438 ! gasvmr(:,:,2) - n2o volume mixing ratio !
439 ! gasvmr(:,:,3) - ch4 volume mixing ratio !
440 ! gasvmr(:,:,4) - o2 volume mixing ratio !
441 ! gasvmr(:,:,5) - co volume mixing ratio !
442 ! gasvmr(:,:,6) - cfc11 volume mixing ratio !
443 ! gasvmr(:,:,7) - cfc12 volume mixing ratio !
444 ! gasvmr(:,:,8) - cfc22 volume mixing ratio !
445 ! gasvmr(:,:,9) - ccl4 volume mixing ratio !
446 ! clouds(npts,nlay,:): layer cloud profiles: !
447 ! (check module_radiation_clouds for definition) !
448 ! --- for ilwcliq > 0 --- !
449 ! clouds(:,:,1) - layer total cloud fraction !
450 ! clouds(:,:,2) - layer in-cloud liq water path (g/m**2) !
451 ! clouds(:,:,3) - mean eff radius for liq cloud (micron) !
452 ! clouds(:,:,4) - layer in-cloud ice water path (g/m**2) !
453 ! clouds(:,:,5) - mean eff radius for ice cloud (micron) !
454 ! clouds(:,:,6) - layer rain drop water path (g/m**2) !
455 ! clouds(:,:,7) - mean eff radius for rain drop (micron) !
456 ! clouds(:,:,8) - layer snow flake water path (g/m**2) !
457 ! clouds(:,:,9) - mean eff radius for snow flake (micron) !
458 ! --- for ilwcliq = 0 --- !
459 ! clouds(:,:,1) - layer total cloud fraction !
460 ! clouds(:,:,2) - layer cloud optical depth !
461 ! clouds(:,:,3) - layer cloud single scattering albedo !
462 ! clouds(:,:,4) - layer cloud asymmetry factor !
463 ! icseed(npts) : auxiliary special cloud related array !
464 ! when module variable isubclw=2, it provides !
465 ! permutation seed for each column profile that !
466 ! are used for generating random numbers. !
467 ! when isubclw /=2, it will not be used. !
468 ! aerosols(npts,nlay,nbands,:) : aerosol optical properties !
469 ! (check module_radiation_aerosols for definition)!
470 ! (:,:,:,1) - optical depth !
471 ! (:,:,:,2) - single scattering albedo !
472 ! (:,:,:,3) - asymmetry parameter !
473 ! sfemis (npts) : surface emissivity !
474 ! sfgtmp (npts) : surface ground temperature (k) !
475 ! npts : total number of horizontal points !
476 ! nlay, nlp1 : total number of vertical layers, levels !
477 ! lprnt : cntl flag for diagnostic print out !
478 ! !
479 ! output variables: !
480 ! hlwc (npts,nlay): total sky heating rate (k/day or k/sec) !
481 ! topflx(npts) : radiation fluxes at top, component: !
482 ! (check module_radlw_paramters for definition) !
483 ! upfxc - total sky upward flux at top (w/m2) !
484 ! upfx0 - clear sky upward flux at top (w/m2) !
485 ! sfcflx(npts) : radiation fluxes at sfc, component: !
486 ! (check module_radlw_paramters for definition) !
487 ! upfxc - total sky upward flux at sfc (w/m2) !
488 ! upfx0 - clear sky upward flux at sfc (w/m2) !
489 ! dnfxc - total sky downward flux at sfc (w/m2) !
490 ! dnfx0 - clear sky downward flux at sfc (w/m2) !
491 ! !
492 !! optional output variables: !
493 ! hlwb(npts,nlay,nbands): spectral band total sky heating rates !
494 ! hlw0 (npts,nlay): clear sky heating rate (k/day or k/sec) !
495 ! flxprf(npts,nlp1): level radiative fluxes (w/m2), components: !
496 ! (check module_radlw_paramters for definition) !
497 ! upfxc - total sky upward flux !
498 ! dnfxc - total sky dnward flux !
499 ! upfx0 - clear sky upward flux !
500 ! dnfx0 - clear sky dnward flux !
501 ! !
502 ! external module variables: (in physparam) !
503 ! ilwrgas - control flag for rare gases (ch4,n2o,o2,cfcs, etc.) !
504 ! =0: do not include rare gases !
505 ! >0: include all rare gases !
506 ! ilwcliq - control flag for liq-cloud optical properties !
507 ! =0: input cloud optical depth, ignor ilwcice !
508 ! =1: input cld liqp & reliq, hu & stamnes (1993) !
509 ! =2: not used !
510 ! ilwcice - control flag for ice-cloud optical properties !
511 ! *** if ilwcliq==0, ilwcice is ignored !
512 ! =1: input cld icep & reice, ebert & curry (1997) !
513 ! =2: input cld icep & reice, streamer (1996) !
514 ! =3: input cld icep & reice, fu (1998) !
515 ! isubclw - sub-column cloud approximation control flag !
516 ! =0: no sub-col cld treatment, use grid-mean cld quantities !
517 ! =1: mcica sub-col, prescribed seeds to get random numbers !
518 ! =2: mcica sub-col, providing array icseed for random numbers!
519 ! iovrlw - cloud overlapping control flag !
520 ! =0: random overlapping clouds !
521 ! =1: maximum/random overlapping clouds !
522 ! =2: maximum overlap cloud (used for isubclw>0 only) !
523 ! ivflip - control flag for vertical index direction !
524 ! =0: vertical index from toa to surface !
525 ! =1: vertical index from surface to toa !
526 ! !
527 ! module parameters, control variables: !
528 ! nbands - number of longwave spectral bands !
529 ! maxgas - maximum number of absorbing gaseous !
530 ! maxxsec - maximum number of cross-sections !
531 ! ngptlw - total number of g-point subintervals !
532 ! ng## - number of g-points in band (##=1-16) !
533 ! ngb(ngptlw) - band indices for each g-point !
534 ! bpade - pade approximation constant (1/0.278) !
535 ! nspa,nspb(nbands)- number of lower/upper ref atm's per band !
536 ! delwave(nbands) - longwave band width (wavenumbers) !
537 ! ipsdlw0 - permutation seed for mcica sub-col clds !
538 ! !
539 ! major local variables: !
540 ! pavel (nlay) - layer pressures (mb) !
541 ! delp (nlay) - layer pressure thickness (mb) !
542 ! tavel (nlay) - layer temperatures (k) !
543 ! tz (0:nlay) - level (interface) temperatures (k) !
544 ! semiss (nbands) - surface emissivity for each band !
545 ! wx (nlay,maxxsec) - cross-section molecules concentration !
546 ! coldry (nlay) - dry air column amount !
547 ! (1.e-20*molecules/cm**2) !
548 ! cldfrc (0:nlp1) - layer cloud fraction !
549 ! taucld (nbands,nlay) - layer cloud optical depth for each band !
550 ! cldfmc (ngptlw,nlay) - layer cloud fraction for each g-point !
551 ! tauaer (nbands,nlay) - aerosol optical depths !
552 ! fracs (ngptlw,nlay) - planck fractions !
553 ! tautot (ngptlw,nlay) - total optical depths (gaseous+aerosols) !
554 ! colamt (nlay,maxgas) - column amounts of absorbing gases !
555 ! 1-maxgas are for watervapor, carbon !
556 ! dioxide, ozone, nitrous oxide, methane, !
557 ! oxigen, carbon monoxide, respectively !
558 ! (molecules/cm**2) !
559 ! pwvcm - column precipitable water vapor (cm) !
560 ! secdiff(nbands) - variable diffusivity angle defined as !
561 ! an exponential function of the column !
562 ! water amount in bands 2-3 and 5-9. !
563 ! this reduces the bias of several w/m2 in !
564 ! downward surface flux in high water !
565 ! profiles caused by using the constant !
566 ! diffusivity angle of 1.66. (mji) !
567 ! facij (nlay) - indicator of interpolation factors !
568 ! =0/1: indicate lower/higher temp & height !
569 ! selffac(nlay) - scale factor for self-continuum, equals !
570 ! (w.v. density)/(atm density at 296K,1013 mb) !
571 ! selffrac(nlay) - factor for temp interpolation of ref !
572 ! self-continuum data !
573 ! indself(nlay) - index of the lower two appropriate ref !
574 ! temp for the self-continuum interpolation !
575 ! forfac (nlay) - scale factor for w.v. foreign-continuum !
576 ! forfrac(nlay) - factor for temp interpolation of ref !
577 ! w.v. foreign-continuum data !
578 ! indfor (nlay) - index of the lower two appropriate ref !
579 ! temp for the foreign-continuum interp !
580 ! laytrop - tropopause layer index at which switch is !
581 ! made from one conbination kew species to !
582 ! another. !
583 ! jp(nlay),jt(nlay),jt1(nlay) !
584 ! - lookup table indexes !
585 ! totuflux(0:nlay) - total-sky upward longwave flux (w/m2) !
586 ! totdflux(0:nlay) - total-sky downward longwave flux (w/m2) !
587 ! htr(nlay) - total-sky heating rate (k/day or k/sec) !
588 ! totuclfl(0:nlay) - clear-sky upward longwave flux (w/m2) !
589 ! totdclfl(0:nlay) - clear-sky downward longwave flux (w/m2) !
590 ! htrcl(nlay) - clear-sky heating rate (k/day or k/sec) !
591 ! fnet (0:nlay) - net longwave flux (w/m2) !
592 ! fnetc (0:nlay) - clear-sky net longwave flux (w/m2) !
593 ! !
594 ! !
595 ! ====================== end of definitions =================== !
596 
597 ! --- inputs:
598  integer, intent(in) :: npts, nlay, nlp1
599  integer, intent(in) :: icseed(npts)
600 
601  logical, intent(in) :: lprnt
602 
603  real (kind=kind_phys), dimension(npts,nlp1), intent(in) :: plvl, &
604  & tlvl
605  real (kind=kind_phys), dimension(npts,nlay), intent(in) :: plyr, &
606  & tlyr, qlyr, olyr
607 
608  real (kind=kind_phys), dimension(npts,nlay,9),intent(in):: gasvmr
609  real (kind=kind_phys), dimension(npts,nlay,9),intent(in):: clouds
610 
611  real (kind=kind_phys), dimension(npts), intent(in) :: sfemis, &
612  & sfgtmp
613 
614  real (kind=kind_phys), dimension(npts,nlay,nbands,3),intent(in):: &
615  & aerosols
616 
617 ! --- outputs:
618  real (kind=kind_phys), dimension(npts,nlay), intent(out) :: hlwc
619 
620  type(topflw_type), dimension(npts), intent(out) :: topflx
621  type(sfcflw_type), dimension(npts), intent(out) :: sfcflx
622 
623 !! --- optional outputs:
624  real (kind=kind_phys), dimension(npts,nlay,nbands),optional, &
625  & intent(out) :: hlwb
626  real (kind=kind_phys), dimension(npts,nlay), optional, &
627  & intent(out) :: hlw0
628  type (proflw_type), dimension(npts,nlp1), optional, &
629  & intent(out) :: flxprf
630 
631 ! --- locals:
632  real (kind=kind_phys), dimension(0:nlp1) :: cldfrc
633 
634  real (kind=kind_phys), dimension(0:nlay) :: totuflux, totdflux, &
635  & totuclfl, totdclfl, tz
636 
637  real (kind=kind_phys), dimension(nlay) :: htr, htrcl
638 
639  real (kind=kind_phys), dimension(nlay) :: pavel, tavel, delp, &
640  & clwp, ciwp, relw, reiw, cda1, cda2, cda3, cda4, &
641  & coldry, colbrd, h2ovmr, o3vmr, fac00, fac01, fac10, fac11, &
642  & selffac, selffrac, forfac, forfrac, minorfrac, scaleminor, &
643  & scaleminorn2, temcol
644 
645  real (kind=kind_phys), dimension(nbands,0:nlay) :: pklev, pklay
646 
647  real (kind=kind_phys), dimension(nlay,nbands) :: htrb
648  real (kind=kind_phys), dimension(nbands,nlay) :: taucld, tauaer
649  real (kind=kind_phys), dimension(ngptlw,nlay) :: fracs, tautot, &
650  & cldfmc
651 
652  real (kind=kind_phys), dimension(nbands) :: semiss, secdiff
653 
654 ! --- column amount of absorbing gases:
655 ! (:,m) m = 1-h2o, 2-co2, 3-o3, 4-n2o, 5-ch4, 6-o2, 7-co
656  real (kind=kind_phys) :: colamt(nlay,maxgas)
657 
658 ! --- column cfc cross-section amounts:
659 ! (:,m) m = 1-ccl4, 2-cfc11, 3-cfc12, 4-cfc22
660  real (kind=kind_phys) :: wx(nlay,maxxsec)
661 
662 ! --- reference ratios of binary species parameter in lower atmosphere:
663 ! (:,m,:) m = 1-h2o/co2, 2-h2o/o3, 3-h2o/n2o, 4-h2o/ch4, 5-n2o/co2, 6-o3/co2
664  real (kind=kind_phys) :: rfrate(nlay,nrates,2)
665 
666  real (kind=kind_phys) :: tem0, tem1, tem2, pwvcm, summol, stemp
667 
668  integer, dimension(npts) :: ipseed
669  integer, dimension(nlay) :: jp, jt, jt1, indself, indfor, indminor
670  integer :: laytrop, iplon, i, j, k, k1
671  logical :: lcf1
672 
673 !
674 !===> ... begin here
675 !
676 
677 ! --- ... initialization
678 
679  lhlwb = present ( hlwb )
680  lhlw0 = present ( hlw0 )
681  lflxprf= present ( flxprf )
682 
683  colamt(:,:) = f_zero
684 
686 ! --- ... change random number seed value for each radiation invocation
687 
688  if ( isubclw == 1 ) then ! advance prescribed permutation seed
689  do i = 1, npts
690  ipseed(i) = ipsdlw0 + i
691  enddo
692  elseif ( isubclw == 2 ) then ! use input array of permutaion seeds
693  do i = 1, npts
694  ipseed(i) = icseed(i)
695  enddo
696  endif
697 
698 ! if ( lprnt ) then
699 ! print *,' In radlw, isubclw, ipsdlw0,ipseed =', &
700 ! & isubclw, ipsdlw0, ipseed
701 ! endif
702 
703 ! --- ... loop over horizontal npts profiles
704 
705  lab_do_iplon : do iplon = 1, npts
706 
708  if (sfemis(iplon) > eps .and. sfemis(iplon) <= 1.0) then ! input surface emissivity
709  do j = 1, nbands
710  semiss(j) = sfemis(iplon)
711  enddo
712  else ! use default values
713  do j = 1, nbands
714  semiss(j) = semiss0(j)
715  enddo
716  endif
717 
718  stemp = sfgtmp(iplon) ! surface ground temp
720 ! --- ... prepare atmospheric profile for use in rrtm
721 ! the vertical index of internal array is from surface to top
722 
723 ! --- ... molecular amounts are input or converted to volume mixing ratio
724 ! and later then converted to molecular amount (molec/cm2) by the
725 ! dry air column coldry (in molec/cm2) which is calculated from the
726 ! layer pressure thickness (in mb), based on the hydrostatic equation
727 ! --- ... and includes a correction to account for h2o in the layer.
728 
729  if (ivflip == 0) then ! input from toa to sfc
730 
731  tem1 = 100.0 * con_g
732  tem2 = 1.0e-20 * 1.0e3 * con_avgd
733  tz(0) = tlvl(iplon,nlp1)
734 
735  do k = 1, nlay
736  k1 = nlp1 - k
737  pavel(k)= plyr(iplon,k1)
738  delp(k) = plvl(iplon,k1+1) - plvl(iplon,k1)
739  tavel(k)= tlyr(iplon,k1)
740  tz(k) = tlvl(iplon,k1)
741 
743 ! --- ... set absorber amount
744 !test use
745 ! h2ovmr(k)= max(f_zero,qlyr(iplon,k1)*amdw) ! input mass mixing ratio
746 ! h2ovmr(k)= max(f_zero,qlyr(iplon,k1)) ! input vol mixing ratio
747 ! o3vmr (k)= max(f_zero,olyr(iplon,k1)) ! input vol mixing ratio
748 !ncep model use
749  h2ovmr(k)= max(f_zero,qlyr(iplon,k1) &
750  & *amdw/(f_one-qlyr(iplon,k1))) ! input specific humidity
751  o3vmr(k)= max(f_zero,olyr(iplon,k1)*amdo3) ! input mass mixing ratio
752 
753 ! --- ... tem0 is the molecular weight of moist air
754  tem0 = (f_one - h2ovmr(k))*con_amd + h2ovmr(k)*con_amw
755  coldry(k) = tem2*delp(k) / (tem1*tem0*(f_one+h2ovmr(k)))
756  temcol(k) = 1.0e-12 * coldry(k)
757 
758  colamt(k,1) = max(f_zero, coldry(k)*h2ovmr(k)) ! h2o
759  colamt(k,2) = max(temcol(k), coldry(k)*gasvmr(iplon,k1,1)) ! co2
760  colamt(k,3) = max(temcol(k), coldry(k)*o3vmr(k)) ! o3
761  enddo
762 
765 ! --- ... set up col amount for rare gases, convert from volume mixing ratio
766 ! to molec/cm2 based on coldry (scaled to 1.0e-20)
767 
768  if (ilwrgas > 0) then
769  do k = 1, nlay
770  k1 = nlp1 - k
771  colamt(k,4)=max(temcol(k), coldry(k)*gasvmr(iplon,k1,2)) ! n2o
772  colamt(k,5)=max(temcol(k), coldry(k)*gasvmr(iplon,k1,3)) ! ch4
773  colamt(k,6)=max(f_zero, coldry(k)*gasvmr(iplon,k1,4)) ! o2
774  colamt(k,7)=max(f_zero, coldry(k)*gasvmr(iplon,k1,5)) ! co
775 
776  wx(k,1) = max( f_zero, coldry(k)*gasvmr(iplon,k1,9) ) ! ccl4
777  wx(k,2) = max( f_zero, coldry(k)*gasvmr(iplon,k1,6) ) ! cf11
778  wx(k,3) = max( f_zero, coldry(k)*gasvmr(iplon,k1,7) ) ! cf12
779  wx(k,4) = max( f_zero, coldry(k)*gasvmr(iplon,k1,8) ) ! cf22
780  enddo
781  else
782  do k = 1, nlay
783  colamt(k,4) = f_zero ! n2o
784  colamt(k,5) = f_zero ! ch4
785  colamt(k,6) = f_zero ! o2
786  colamt(k,7) = f_zero ! co
787 
788  wx(k,1) = f_zero
789  wx(k,2) = f_zero
790  wx(k,3) = f_zero
791  wx(k,4) = f_zero
792  enddo
793  endif
794 
796 ! --- ... set aerosol optical properties
797 
798  do k = 1, nlay
799  k1 = nlp1 - k
800  do j = 1, nbands
801  tauaer(j,k) = aerosols(iplon,k1,j,1) &
802  & * (f_one - aerosols(iplon,k1,j,2))
803  enddo
804  enddo
805 
807  if (ilwcliq > 0) then ! use prognostic cloud method
808  do k = 1, nlay
809  k1 = nlp1 - k
810  cldfrc(k)= clouds(iplon,k1,1)
811  clwp(k) = clouds(iplon,k1,2)
812  relw(k) = clouds(iplon,k1,3)
813  ciwp(k) = clouds(iplon,k1,4)
814  reiw(k) = clouds(iplon,k1,5)
815  cda1(k) = clouds(iplon,k1,6)
816  cda2(k) = clouds(iplon,k1,7)
817  cda3(k) = clouds(iplon,k1,8)
818  cda4(k) = clouds(iplon,k1,9)
819  enddo
820  else ! use diagnostic cloud method
821  do k = 1, nlay
822  k1 = nlp1 - k
823  cldfrc(k)= clouds(iplon,k1,1)
824  cda1(k) = clouds(iplon,k1,2)
825  enddo
826  endif ! end if_ilwcliq
827 
828  cldfrc(0) = f_one ! padding value only
829  cldfrc(nlp1) = f_zero ! padding value only
830 
832 ! --- ... compute precipitable water vapor for diffusivity angle adjustments
833 
834  tem1 = f_zero
835  tem2 = f_zero
836  do k = 1, nlay
837  tem1 = tem1 + coldry(k) + colamt(k,1)
838  tem2 = tem2 + colamt(k,1)
839  enddo
840 
841  tem0 = 10.0 * tem2 / (amdw * tem1 * con_g)
842  pwvcm = tem0 * plvl(iplon,nlp1)
843 
844  else ! input from sfc to toa
845 
846  tem1 = 100.0 * con_g
847  tem2 = 1.0e-20 * 1.0e3 * con_avgd
848  tz(0) = tlvl(iplon,1)
849 
850  do k = 1, nlay
851  pavel(k)= plyr(iplon,k)
852  delp(k) = plvl(iplon,k) - plvl(iplon,k+1)
853  tavel(k)= tlyr(iplon,k)
854  tz(k) = tlvl(iplon,k+1)
855 
856 ! --- ... set absorber amount
857 !test use
858 ! h2ovmr(k)= max(f_zero,qlyr(iplon,k)*amdw) ! input mass mixing ratio
859 ! h2ovmr(k)= max(f_zero,qlyr(iplon,k)) ! input vol mixing ratio
860 ! o3vmr (k)= max(f_zero,olyr(iplon,k)) ! input vol mixing ratio
861 !ncep model use
862  h2ovmr(k)= max(f_zero,qlyr(iplon,k) &
863  & *amdw/(f_one-qlyr(iplon,k))) ! input specific humidity
864  o3vmr(k)= max(f_zero,olyr(iplon,k)*amdo3) ! input mass mixing ratio
865 
866 ! --- ... tem0 is the molecular weight of moist air
867  tem0 = (f_one - h2ovmr(k))*con_amd + h2ovmr(k)*con_amw
868  coldry(k) = tem2*delp(k) / (tem1*tem0*(f_one+h2ovmr(k)))
869  temcol(k) = 1.0e-12 * coldry(k)
870 
871  colamt(k,1) = max(f_zero, coldry(k)*h2ovmr(k)) ! h2o
872  colamt(k,2) = max(temcol(k), coldry(k)*gasvmr(iplon,k,1)) ! co2
873  colamt(k,3) = max(temcol(k), coldry(k)*o3vmr(k)) ! o3
874  enddo
875 
876 ! --- ... set up col amount for rare gases, convert from volume mixing ratio
877 ! to molec/cm2 based on coldry (scaled to 1.0e-20)
878 
879  if (ilwrgas > 0) then
880  do k = 1, nlay
881  colamt(k,4)=max(temcol(k), coldry(k)*gasvmr(iplon,k,2)) ! n2o
882  colamt(k,5)=max(temcol(k), coldry(k)*gasvmr(iplon,k,3)) ! ch4
883  colamt(k,6)=max(f_zero, coldry(k)*gasvmr(iplon,k,4)) ! o2
884  colamt(k,7)=max(f_zero, coldry(k)*gasvmr(iplon,k,5)) ! co
885 
886  wx(k,1) = max( f_zero, coldry(k)*gasvmr(iplon,k,9) ) ! ccl4
887  wx(k,2) = max( f_zero, coldry(k)*gasvmr(iplon,k,6) ) ! cf11
888  wx(k,3) = max( f_zero, coldry(k)*gasvmr(iplon,k,7) ) ! cf12
889  wx(k,4) = max( f_zero, coldry(k)*gasvmr(iplon,k,8) ) ! cf22
890  enddo
891  else
892  do k = 1, nlay
893  colamt(k,4) = f_zero ! n2o
894  colamt(k,5) = f_zero ! ch4
895  colamt(k,6) = f_zero ! o2
896  colamt(k,7) = f_zero ! co
897 
898  wx(k,1) = f_zero
899  wx(k,2) = f_zero
900  wx(k,3) = f_zero
901  wx(k,4) = f_zero
902  enddo
903  endif
904 
905 ! --- ... set aerosol optical properties
906 
907  do j = 1, nbands
908  do k = 1, nlay
909  tauaer(j,k) = aerosols(iplon,k,j,1) &
910  & * (f_one - aerosols(iplon,k,j,2))
911  enddo
912  enddo
913 
914  if (ilwcliq > 0) then ! use prognostic cloud method
915  do k = 1, nlay
916  cldfrc(k)= clouds(iplon,k,1)
917  clwp(k) = clouds(iplon,k,2)
918  relw(k) = clouds(iplon,k,3)
919  ciwp(k) = clouds(iplon,k,4)
920  reiw(k) = clouds(iplon,k,5)
921  cda1(k) = clouds(iplon,k,6)
922  cda2(k) = clouds(iplon,k,7)
923  cda3(k) = clouds(iplon,k,8)
924  cda4(k) = clouds(iplon,k,9)
925  enddo
926  else ! use diagnostic cloud method
927  do k = 1, nlay
928  cldfrc(k)= clouds(iplon,k,1)
929  cda1(k) = clouds(iplon,k,2)
930  enddo
931  endif ! end if_ilwcliq
932 
933  cldfrc(0) = f_one ! padding value only
934  cldfrc(nlp1) = f_zero ! padding value only
935 
936 ! --- ... compute precipitable water vapor for diffusivity angle adjustments
937 
938  tem1 = f_zero
939  tem2 = f_zero
940  do k = 1, nlay
941  tem1 = tem1 + coldry(k) + colamt(k,1)
942  tem2 = tem2 + colamt(k,1)
943  enddo
944 
945  tem0 = 10.0 * tem2 / (amdw * tem1 * con_g)
946  pwvcm = tem0 * plvl(iplon,1)
947 
948  endif ! if_ivflip
949 
951 ! --- ... compute column amount for broadening gases
952 
953  do k = 1, nlay
954  summol = f_zero
955  do i = 2, maxgas
956  summol = summol + colamt(k,i)
957  enddo
958  colbrd(k) = coldry(k) - summol
959  enddo
960 
962 ! --- ... compute diffusivity angle adjustments
963 
964  tem1 = 1.80
965  tem2 = 1.50
966  do j = 1, nbands
967  if (j==1 .or. j==4 .or. j==10) then
968  secdiff(j) = 1.66
969  else
970  secdiff(j) = min( tem1, max( tem2, &
971  & a0(j)+a1(j)*exp(a2(j)*pwvcm) ))
972  endif
973  enddo
974 
975 ! if (lprnt) then
976 ! print *,' coldry',coldry
977 ! print *,' wx(*,1) ',(wx(k,1),k=1,NLAY)
978 ! print *,' wx(*,2) ',(wx(k,2),k=1,NLAY)
979 ! print *,' wx(*,3) ',(wx(k,3),k=1,NLAY)
980 ! print *,' wx(*,4) ',(wx(k,4),k=1,NLAY)
981 ! print *,' iplon ',iplon
982 ! print *,' pavel ',pavel
983 ! print *,' delp ',delp
984 ! print *,' tavel ',tavel
985 ! print *,' tz ',tz
986 ! print *,' h2ovmr ',h2ovmr
987 ! print *,' o3vmr ',o3vmr
988 ! endif
989 
991 ! --- ... for cloudy atmosphere, use cldprop to set cloud optical properties
992 
993  lcf1 = .false.
994  lab_do_k0 : do k = 1, nlay
995  if ( cldfrc(k) > eps ) then
996  lcf1 = .true.
997  exit lab_do_k0
998  endif
999  enddo lab_do_k0
1000 
1001  if ( lcf1 ) then
1002 
1003  call cldprop &
1004 ! --- inputs:
1005  & ( cldfrc,clwp,relw,ciwp,reiw,cda1,cda2,cda3,cda4, &
1006  & nlay, nlp1, ipseed(iplon), &
1007 ! --- outputs:
1008  & cldfmc, taucld &
1009  & )
1010 
1011  else
1012  cldfmc = f_zero
1013  taucld = f_zero
1014  endif
1015 
1016 ! if (lprnt) then
1017 ! print *,' after cldprop'
1018 ! print *,' clwp',clwp
1019 ! print *,' ciwp',ciwp
1020 ! print *,' relw',relw
1021 ! print *,' reiw',reiw
1022 ! print *,' taucl',cda1
1023 ! print *,' cldfrac',cldfrc
1024 ! endif
1025 
1028  call setcoef &
1029 ! --- inputs:
1030  & ( pavel,tavel,tz,stemp,h2ovmr,colamt,coldry,colbrd, &
1031  & nlay, nlp1, &
1032 ! --- outputs:
1033  & laytrop,pklay,pklev,jp,jt,jt1, &
1034  & rfrate,fac00,fac01,fac10,fac11, &
1035  & selffac,selffrac,indself,forfac,forfrac,indfor, &
1036  & minorfrac,scaleminor,scaleminorn2,indminor &
1037  & )
1038 
1039 ! if (lprnt) then
1040 ! print *,'laytrop',laytrop
1041 ! print *,'colh2o',(colamt(k,1),k=1,NLAY)
1042 ! print *,'colco2',(colamt(k,2),k=1,NLAY)
1043 ! print *,'colo3', (colamt(k,3),k=1,NLAY)
1044 ! print *,'coln2o',(colamt(k,4),k=1,NLAY)
1045 ! print *,'colch4',(colamt(k,5),k=1,NLAY)
1046 ! print *,'fac00',fac00
1047 ! print *,'fac01',fac01
1048 ! print *,'fac10',fac10
1049 ! print *,'fac11',fac11
1050 ! print *,'jp',jp
1051 ! print *,'jt',jt
1052 ! print *,'jt1',jt1
1053 ! print *,'selffac',selffac
1054 ! print *,'selffrac',selffrac
1055 ! print *,'indself',indself
1056 ! print *,'forfac',forfac
1057 ! print *,'forfrac',forfrac
1058 ! print *,'indfor',indfor
1059 ! endif
1060 
1062 ! --- ... calculate the gaseous optical depths and Planck fractions for
1063 ! each longwave spectral band.
1064 
1065  call taumol &
1066 ! --- inputs:
1067  & ( laytrop,pavel,coldry,colamt,colbrd,wx,tauaer, &
1068  & rfrate,fac00,fac01,fac10,fac11,jp,jt,jt1, &
1069  & selffac,selffrac,indself,forfac,forfrac,indfor, &
1070  & minorfrac,scaleminor,scaleminorn2,indminor, &
1071  & nlay, &
1072 ! --- outputs:
1073  & fracs, tautot &
1074  & )
1075 
1076 ! if (lprnt) then
1077 ! print *,' after taumol'
1078 ! do k = 1, nlay
1079 ! write(6,121) k
1080 !121 format(' k =',i3,5x,'FRACS')
1081 ! write(6,122) (fracs(j,k),j=1,ngptlw)
1082 !122 format(10e14.7)
1083 ! write(6,123) k
1084 !123 format(' k =',i3,5x,'TAUTOT')
1085 ! write(6,122) (tautot(j,k),j=1,ngptlw)
1086 ! enddo
1087 ! endif
1088 
1093 ! --- ... call the radiative transfer routine based on cloud scheme
1094 ! selection. clear sky calculation is done at the same time.
1095 
1096  if (isubclw <= 0) then
1097 
1098  if (iovrlw <= 0) then
1099 
1100  call rtrn &
1101 ! --- inputs:
1102  & ( semiss,delp,cldfrc,taucld,tautot,pklay,pklev, &
1103  & fracs,secdiff,nlay,nlp1, &
1104 ! --- outputs:
1105  & totuflux,totdflux,htr, totuclfl,totdclfl,htrcl, htrb &
1106  & )
1107 
1108  else
1109 
1110  call rtrnmr &
1111 ! --- inputs:
1112  & ( semiss,delp,cldfrc,taucld,tautot,pklay,pklev, &
1113  & fracs,secdiff,nlay,nlp1, &
1114 ! --- outputs:
1115  & totuflux,totdflux,htr, totuclfl,totdclfl,htrcl, htrb &
1116  & )
1117 
1118  endif ! end if_iovrlw_block
1119 
1120  else
1121 
1122  call rtrnmc &
1123 ! --- inputs:
1124  & ( semiss,delp,cldfmc,taucld,tautot,pklay,pklev, &
1125  & fracs,secdiff,nlay,nlp1, &
1126 ! --- outputs:
1127  & totuflux,totdflux,htr, totuclfl,totdclfl,htrcl, htrb &
1128  & )
1129 
1130  endif ! end if_isubclw_block
1132 ! --- ... output total-sky and clear-sky fluxes and heating rates
1133 
1134  topflx(iplon)%upfxc = totuflux(nlay)
1135  topflx(iplon)%upfx0 = totuclfl(nlay)
1136 
1137  sfcflx(iplon)%upfxc = totuflux(0)
1138  sfcflx(iplon)%upfx0 = totuclfl(0)
1139  sfcflx(iplon)%dnfxc = totdflux(0)
1140  sfcflx(iplon)%dnfx0 = totdclfl(0)
1141 
1142  if (ivflip == 0) then ! output from toa to sfc
1143 
1144 !! --- ... optional fluxes
1145  if ( lflxprf ) then
1146  do k = 0, nlay
1147  k1 = nlp1 - k
1148  flxprf(iplon,k1)%upfxc = totuflux(k)
1149  flxprf(iplon,k1)%dnfxc = totdflux(k)
1150  flxprf(iplon,k1)%upfx0 = totuclfl(k)
1151  flxprf(iplon,k1)%dnfx0 = totdclfl(k)
1152  enddo
1153  endif
1154 
1155  do k = 1, nlay
1156  k1 = nlp1 - k
1157  hlwc(iplon,k1) = htr(k)
1158  enddo
1159 
1160 !! --- ... optional clear sky heating rate
1161  if ( lhlw0 ) then
1162  do k = 1, nlay
1163  k1 = nlp1 - k
1164  hlw0(iplon,k1) = htrcl(k)
1165  enddo
1166  endif
1167 
1168 !! --- ... optional spectral band heating rate
1169  if ( lhlwb ) then
1170  do j = 1, nbands
1171  do k = 1, nlay
1172  k1 = nlp1 - k
1173  hlwb(iplon,k1,j) = htrb(k,j)
1174  enddo
1175  enddo
1176  endif
1177 
1178  else ! output from sfc to toa
1179 
1180 !! --- ... optional fluxes
1181  if ( lflxprf ) then
1182  do k = 0, nlay
1183  flxprf(iplon,k+1)%upfxc = totuflux(k)
1184  flxprf(iplon,k+1)%dnfxc = totdflux(k)
1185  flxprf(iplon,k+1)%upfx0 = totuclfl(k)
1186  flxprf(iplon,k+1)%dnfx0 = totdclfl(k)
1187  enddo
1188  endif
1189 
1190  do k = 1, nlay
1191  hlwc(iplon,k) = htr(k)
1192  enddo
1193 
1194 !! --- ... optional clear sky heating rate
1195  if ( lhlw0 ) then
1196  do k = 1, nlay
1197  hlw0(iplon,k) = htrcl(k)
1198  enddo
1199  endif
1200 
1201 !! --- ... optional spectral band heating rate
1202  if ( lhlwb ) then
1203  do j = 1, nbands
1204  do k = 1, nlay
1205  hlwb(iplon,k,j) = htrb(k,j)
1206  enddo
1207  enddo
1208  endif
1209 
1210  endif ! if_ivflip
1211 
1212  enddo lab_do_iplon
1213 
1214 !...................................
1215  end subroutine lwrad
1216 !-----------------------------------
1218 
1246 !-----------------------------------
1247  subroutine rlwinit &
1248 !...................................
1249 ! --- inputs:
1250  & ( me )
1251 ! --- outputs: (none)
1252 
1253 ! =================== program usage description =================== !
1254 ! !
1255 ! purpose: initialize non-varying module variables, conversion factors,!
1256 ! and look-up tables. !
1257 ! !
1258 ! subprograms called: none !
1259 ! !
1260 ! ==================== defination of variables ==================== !
1261 ! !
1262 ! inputs: !
1263 ! me - print control for parallel process !
1264 ! !
1265 ! outputs: (none) !
1266 ! !
1267 ! external module variables: (in physparam) !
1268 ! ilwrate - heating rate unit selections !
1269 ! =1: output in k/day !
1270 ! =2: output in k/second !
1271 ! ilwrgas - control flag for rare gases (ch4,n2o,o2,cfcs, etc.) !
1272 ! =0: do not include rare gases !
1273 ! >0: include all rare gases !
1274 ! ilwcliq - liquid cloud optical properties contrl flag !
1275 ! =0: input cloud opt depth from diagnostic scheme !
1276 ! >0: input cwp,rew, and other cloud content parameters !
1277 ! isubclw - sub-column cloud approximation control flag !
1278 ! =0: no sub-col cld treatment, use grid-mean cld quantities !
1279 ! =1: mcica sub-col, prescribed seeds to get random numbers !
1280 ! =2: mcica sub-col, providing array icseed for random numbers!
1281 ! icldflg - cloud scheme control flag !
1282 ! =0: diagnostic scheme gives cloud tau, omiga, and g. !
1283 ! =1: prognostic scheme gives cloud liq/ice path, etc. !
1284 ! iovrlw - clouds vertical overlapping control flag !
1285 ! =0: random overlapping clouds !
1286 ! =1: maximum/random overlapping clouds !
1287 ! =2: maximum overlap cloud (isubcol>0 only) !
1288 ! !
1289 ! ******************************************************************* !
1290 ! original code description !
1291 ! !
1292 ! original version: michael j. iacono; july, 1998 !
1293 ! first revision for ncar ccm: september, 1998 !
1294 ! second revision for rrtm_v3.0: september, 2002 !
1295 ! !
1296 ! this subroutine performs calculations necessary for the initialization
1297 ! of the longwave model. lookup tables are computed for use in the lw !
1298 ! radiative transfer, and input absorption coefficient data for each !
1299 ! spectral band are reduced from 256 g-point intervals to 140. !
1300 ! !
1301 ! ******************************************************************* !
1302 ! !
1303 ! definitions: !
1304 ! arrays for 10000-point look-up tables: !
1305 ! tau_tbl - clear-sky optical depth (used in cloudy radiative transfer!
1306 ! exp_tbl - exponential lookup table for tansmittance !
1307 ! tfn_tbl - tau transition function; i.e. the transition of the Planck!
1308 ! function from that for the mean layer temperature to that !
1309 ! for the layer boundary temperature as a function of optical
1310 ! depth. the "linear in tau" method is used to make the table
1311 ! !
1312 ! ******************************************************************* !
1313 ! !
1314 ! ====================== end of description block ================= !
1315 
1316 ! --- inputs:
1317  integer, intent(in) :: me
1318 
1319 ! --- outputs: none
1320 
1321 ! --- locals:
1322  real (kind=kind_phys), parameter :: expeps = 1.e-20
1323 
1324  real (kind=kind_phys) :: tfn, pival, explimit
1325 
1326  integer :: i
1327 
1328 !
1329 !===> ... begin here
1330 !
1331  if ( iovrlw<0 .or. iovrlw>2 ) then
1332  print *,' *** Error in specification of cloud overlap flag', &
1333  & ' IOVRLW=',iovrlw,' in RLWINIT !!'
1334  stop
1335  elseif ( iovrlw==2 .and. isubclw==0 ) then
1336  if (me == 0) then
1337  print *,' *** IOVRLW=2 - maximum cloud overlap, is not yet', &
1338  & ' available for ISUBCLW=0 setting!!'
1339  print *,' The program uses maximum/random overlap', &
1340  & ' instead.'
1341  endif
1342 
1343  iovrlw = 1
1344  endif
1345 
1346  if (me == 0) then
1347  print *,' - Using AER Longwave Radiation, Version: ', vtaglw
1348 
1349  if (ilwrgas > 0) then
1350  print *,' --- Include rare gases N2O, CH4, O2, CFCs ', &
1351  & 'absorptions in LW'
1352  else
1353  print *,' --- Rare gases effect is NOT included in LW'
1354  endif
1355 
1356  if ( isubclw == 0 ) then
1357  print *,' --- Using standard grid average clouds, no ', &
1358  & 'sub-column clouds approximation applied'
1359  elseif ( isubclw == 1 ) then
1360  print *,' --- Using MCICA sub-colum clouds approximation ', &
1361  & 'with a prescribed sequence of permutaion seeds'
1362  elseif ( isubclw == 2 ) then
1363  print *,' --- Using MCICA sub-colum clouds approximation ', &
1364  & 'with provided input array of permutation seeds'
1365  else
1366  print *,' *** Error in specification of sub-column cloud ', &
1367  & ' control flag isubclw =',isubclw,' !!'
1368  stop
1369  endif
1370  endif
1371 
1372 ! --- ... check cloud flags for consistency
1373 
1374  if ((icldflg == 0 .and. ilwcliq /= 0) .or. &
1375  & (icldflg == 1 .and. ilwcliq == 0)) then
1376  print *,' *** Model cloud scheme inconsistent with LW', &
1377  & ' radiation cloud radiative property setup !!'
1378  stop
1379  endif
1380 
1381 ! --- ... setup default surface emissivity for each band here
1382 
1383  semiss0(:) = f_one
1384 
1385 ! --- ... setup constant factors for flux and heating rate
1386 ! the 1.0e-2 is to convert pressure from mb to N/m**2
1387 
1388  pival = 2.0 * asin(f_one)
1389  fluxfac = pival * 2.0d4
1390 ! fluxfac = 62831.85307179586 ! = 2 * pi * 1.0e4
1391 
1392  if (ilwrate == 1) then
1393 ! heatfac = 8.4391
1394 ! heatfac = con_g * 86400. * 1.0e-2 / con_cp ! (in k/day)
1395  heatfac = con_g * 864.0 / con_cp ! (in k/day)
1396  else
1397  heatfac = con_g * 1.0e-2 / con_cp ! (in k/second)
1398  endif
1399 
1400 ! --- ... compute lookup tables for transmittance, tau transition
1401 ! function, and clear sky tau (for the cloudy sky radiative
1402 ! transfer). tau is computed as a function of the tau
1403 ! transition function, transmittance is calculated as a
1404 ! function of tau, and the tau transition function is
1405 ! calculated using the linear in tau formulation at values of
1406 ! tau above 0.01. tf is approximated as tau/6 for tau < 0.01.
1407 ! all tables are computed at intervals of 0.001. the inverse
1408 ! of the constant used in the pade approximation to the tau
1409 ! transition function is set to b.
1410 
1411  tau_tbl(0) = f_zero
1412  exp_tbl(0) = f_one
1413  tfn_tbl(0) = f_zero
1414 
1415  tau_tbl(ntbl) = 1.e10
1416  exp_tbl(ntbl) = expeps
1417  tfn_tbl(ntbl) = f_one
1418 
1419  explimit = aint( -log(tiny(exp_tbl(0))) )
1420 
1421  do i = 1, ntbl-1
1422 !org tfn = float(i) / float(ntbl)
1423 !org tau_tbl(i) = bpade * tfn / (f_one - tfn)
1424  tfn = real(i, kind_phys) / real(ntbl-i, kind_phys)
1425  tau_tbl(i) = bpade * tfn
1426  if (tau_tbl(i) >= explimit) then
1427  exp_tbl(i) = expeps
1428  else
1429  exp_tbl(i) = exp( -tau_tbl(i) )
1430  endif
1431 
1432  if (tau_tbl(i) < 0.06) then
1433  tfn_tbl(i) = tau_tbl(i) / 6.0
1434  else
1435  tfn_tbl(i) = f_one - 2.0*( (f_one / tau_tbl(i)) &
1436  & - ( exp_tbl(i) / (f_one - exp_tbl(i)) ) )
1437  endif
1438  enddo
1439 
1440 !...................................
1441  end subroutine rlwinit
1442 !-----------------------------------
1443 
1446 ! ----------------------------
1447  subroutine cldprop &
1448 ! ............................
1449 ! --- inputs:
1450  & ( cfrac,cliqp,reliq,cicep,reice,cdat1,cdat2,cdat3,cdat4, &
1451  & nlay, nlp1, ipseed, &
1452 ! --- outputs:
1453  & cldfmc, taucld &
1454  & )
1455 
1456 ! =================== program usage description =================== !
1457 ! !
1458 ! purpose: compute the cloud optical depth(s) for each cloudy layer !
1459 ! and g-point interval. !
1460 ! !
1461 ! subprograms called: none !
1462 ! !
1463 ! ==================== defination of variables ==================== !
1464 ! !
1465 ! inputs: -size- !
1466 ! cfrac - real, layer cloud fraction 0:nlp1 !
1467 ! ..... for ilwcliq > 0 (prognostic cloud sckeme) - - - !
1468 ! cliqp - real, layer in-cloud liq water path (g/m**2) nlay !
1469 ! reliq - real, mean eff radius for liq cloud (micron) nlay !
1470 ! cicep - real, layer in-cloud ice water path (g/m**2) nlay !
1471 ! reice - real, mean eff radius for ice cloud (micron) nlay !
1472 ! cdat1 - real, layer rain drop water path (g/m**2) nlay !
1473 ! cdat2 - real, effective radius for rain drop (microm) nlay !
1474 ! cdat3 - real, layer snow flake water path (g/m**2) nlay !
1475 ! cdat4 - real, effective radius for snow flakes (micron) nlay !
1476 ! ..... for ilwcliq = 0 (diagnostic cloud sckeme) - - - !
1477 ! cdat1 - real, input cloud optical depth nlay !
1478 ! cdat2 - real, layer cloud single scattering albedo nlay !
1479 ! cdat3 - real, layer cloud asymmetry factor nlay !
1480 ! cdat4 - real, optional use nlay !
1481 ! cliqp - not used nlay !
1482 ! reliq - not used nlay !
1483 ! cicep - not used nlay !
1484 ! reice - not used nlay !
1485 ! !
1486 ! nlay - integer, number of vertical layers 1 !
1487 ! nlp1 - integer, number of vertical levels 1 !
1488 ! ipseed- permutation seed for generating random numbers (isubclw>0) !
1489 ! !
1490 ! outputs: !
1491 ! cldfmc - real, cloud fraction for each sub-column ngptlw*nlay!
1492 ! taucld - real, cld opt depth for bands (non-mcica) nbands*nlay!
1493 ! !
1494 ! explanation of the method for each value of ilwcliq, and ilwcice. !
1495 ! set up in module "module_radlw_cntr_para" !
1496 ! !
1497 ! ilwcliq=0 : input cloud optical property (tau, ssa, asy). !
1498 ! (used for diagnostic cloud method) !
1499 ! ilwcliq>0 : input cloud liq/ice path and effective radius, also !
1500 ! require the user of 'ilwcice' to specify the method !
1501 ! used to compute aborption due to water/ice parts. !
1502 ! ................................................................... !
1503 ! !
1504 ! ilwcliq=1: the water droplet effective radius (microns) is input!
1505 ! and the opt depths due to water clouds are computed !
1506 ! as in hu and stamnes, j., clim., 6, 728-742, (1993). !
1507 ! the values for absorption coefficients appropriate for
1508 ! the spectral bands in rrtm have been obtained for a !
1509 ! range of effective radii by an averaging procedure !
1510 ! based on the work of j. pinto (private communication).
1511 ! linear interpolation is used to get the absorption !
1512 ! coefficients for the input effective radius. !
1513 ! !
1514 ! ilwcice=1: the cloud ice path (g/m2) and ice effective radius !
1515 ! (microns) are input and the optical depths due to ice!
1516 ! clouds are computed as in ebert and curry, jgr, 97, !
1517 ! 3831-3836 (1992). the spectral regions in this work !
1518 ! have been matched with the spectral bands in rrtm to !
1519 ! as great an extent as possible: !
1520 ! e&c 1 ib = 5 rrtm bands 9-16 !
1521 ! e&c 2 ib = 4 rrtm bands 6-8 !
1522 ! e&c 3 ib = 3 rrtm bands 3-5 !
1523 ! e&c 4 ib = 2 rrtm band 2 !
1524 ! e&c 5 ib = 1 rrtm band 1 !
1525 ! ilwcice=2: the cloud ice path (g/m2) and ice effective radius !
1526 ! (microns) are input and the optical depths due to ice!
1527 ! clouds are computed as in rt code, streamer v3.0 !
1528 ! (ref: key j., streamer user's guide, cooperative !
1529 ! institute for meteorological satellite studies, 2001,!
1530 ! 96 pp.) valid range of values for re are between 5.0 !
1531 ! and 131.0 micron. !
1532 ! ilwcice=3: the ice generalized effective size (dge) is input and!
1533 ! the optical properties, are calculated as in q. fu, !
1534 ! j. climate, (1998). q. fu provided high resolution !
1535 ! tales which were appropriately averaged for the bands!
1536 ! in rrtm_lw. linear interpolation is used to get the !
1537 ! coeff from the stored tables. valid range of values !
1538 ! for deg are between 5.0 and 140.0 micron. !
1539 ! !
1540 ! other cloud control module variables: !
1541 ! isubclw =0: standard cloud scheme, no sub-col cloud approximation !
1542 ! >0: mcica sub-col cloud scheme using ipseed as permutation!
1543 ! seed for generating rundom numbers !
1544 ! !
1545 ! ====================== end of description block ================= !
1546 !
1548 
1549 ! --- inputs:
1550  integer, intent(in) :: nlay, nlp1, ipseed
1551 
1552  real (kind=kind_phys), dimension(0:nlp1), intent(in) :: cfrac
1553  real (kind=kind_phys), dimension(nlay), intent(in) :: cliqp, &
1554  & reliq, cicep, reice, cdat1, cdat2, cdat3, cdat4
1555 
1556 ! --- outputs:
1557  real (kind=kind_phys), dimension(ngptlw,nlay),intent(out):: cldfmc
1558  real (kind=kind_phys), dimension(nbands,nlay),intent(out):: taucld
1559 
1560 ! --- locals:
1561  real (kind=kind_phys), dimension(nbands) :: tauliq, tauice
1562  real (kind=kind_phys), dimension(nlay) :: cldf
1563 
1564  real (kind=kind_phys) :: dgeice, factor, fint, tauran, tausnw, &
1565  & cldliq, refliq, cldice, refice
1566 
1567  logical :: lcloudy(ngptlw,nlay)
1568  integer :: ia, ib, ig, k, index
1569 
1570 !
1571 !===> ... begin here
1572 !
1573  do k = 1, nlay
1574  do ib = 1, nbands
1575  taucld(ib,k) = f_zero
1576  enddo
1577  enddo
1578 
1579  do k = 1, nlay
1580  do ig = 1, ngptlw
1581  cldfmc(ig,k) = f_zero
1582  enddo
1583  enddo
1584 
1585 ! --- ... compute cloud radiative properties for a cloudy column
1586 
1587  lab_if_ilwcliq : if (ilwcliq > 0) then
1588 
1589  lab_do_k : do k = 1, nlay
1590  lab_if_cld : if (cfrac(k) > cldmin) then
1591 
1592  tauran = absrain * cdat1(k) ! ncar formula
1593 !! tausnw = abssnow1 * cdat3(k) ! ncar formula
1594 ! --- if use fu's formula it needs to be normalized by snow density
1595 ! !not use snow density = 0.1 g/cm**3 = 0.1 g/(mu * m**2)
1596 ! use ice density = 0.9167 g/cm**3 = 0.9167 g/(mu * m**2)
1597 ! factor 1.5396=8/(3*sqrt(3)) converts reff to generalized ice particle size
1598 ! use newer factor value 1.0315
1599 ! 1/(0.9167*1.0315) = 1.05756
1600  if (cdat3(k)>f_zero .and. cdat4(k)>10.0_kind_phys) then
1601  tausnw = abssnow0*1.05756*cdat3(k)/cdat4(k) ! fu's formula
1602  else
1603  tausnw = f_zero
1604  endif
1605 
1606  cldliq = cliqp(k)
1607  cldice = cicep(k)
1608 ! refliq = max(2.5e0, min(60.0e0, reliq(k) ))
1609 ! refice = max(5.0e0, reice(k) )
1610  refliq = reliq(k)
1611  refice = reice(k)
1612 
1613 ! --- ... calculation of absorption coefficients due to water clouds.
1614 
1615  if ( cldliq <= f_zero ) then
1616  do ib = 1, nbands
1617  tauliq(ib) = f_zero
1618  enddo
1619  else
1620  if ( ilwcliq == 1 ) then
1621 
1622  factor = refliq - 1.5
1623  index = max( 1, min( 57, int( factor ) ))
1624  fint = factor - float(index)
1625 
1626  do ib = 1, nbands
1627  tauliq(ib) = max(f_zero, cldliq*(absliq1(index,ib) &
1628  & + fint*(absliq1(index+1,ib)-absliq1(index,ib)) ))
1629  enddo
1630  endif ! end if_ilwcliq_block
1631  endif ! end if_cldliq_block
1632 
1633 ! --- ... calculation of absorption coefficients due to ice clouds.
1634 
1635  if ( cldice <= f_zero ) then
1636  do ib = 1, nbands
1637  tauice(ib) = f_zero
1638  enddo
1639  else
1640 
1641 ! --- ... ebert and curry approach for all particle sizes though somewhat
1642 ! unjustified for large ice particles
1643 
1644  if ( ilwcice == 1 ) then
1645  refice = min(130.0, max(13.0, real(refice) ))
1646 
1647  do ib = 1, nbands
1648  ia = ipat(ib) ! eb_&_c band index for ice cloud coeff
1649  tauice(ib) = max(f_zero, cldice*(absice1(1,ia) &
1650  & + absice1(2,ia)/refice) )
1651  enddo
1652 
1653 ! --- ... streamer approach for ice effective radius between 5.0 and 131.0 microns
1654 ! and ebert and curry approach for ice eff radius greater than 131.0 microns.
1655 ! no smoothing between the transition of the two methods.
1656 
1657  elseif ( ilwcice == 2 ) then
1658 
1659  factor = (refice - 2.0) / 3.0
1660  index = max( 1, min( 42, int( factor ) ))
1661  fint = factor - float(index)
1662 
1663  do ib = 1, nbands
1664  tauice(ib) = max(f_zero, cldice*(absice2(index,ib) &
1665  & + fint*(absice2(index+1,ib) - absice2(index,ib)) ))
1666  enddo
1667 
1668 ! --- ... fu's approach for ice effective radius between 4.8 and 135 microns
1669 ! (generalized effective size from 5 to 140 microns)
1670 
1671  elseif ( ilwcice == 3 ) then
1672 
1673 ! dgeice = max(5.0, 1.5396*refice) ! v4.4 value
1674  dgeice = max(5.0, 1.0315*refice) ! v4.71 value
1675  factor = (dgeice - 2.0) / 3.0
1676  index = max( 1, min( 45, int( factor ) ))
1677  fint = factor - float(index)
1678 
1679  do ib = 1, nbands
1680  tauice(ib) = max(f_zero, cldice*(absice3(index,ib) &
1681  & + fint*(absice3(index+1,ib) - absice3(index,ib)) ))
1682  enddo
1683 
1684  endif ! end if_ilwcice_block
1685  endif ! end if_cldice_block
1686 
1687  do ib = 1, nbands
1688  taucld(ib,k) = tauice(ib) + tauliq(ib) + tauran + tausnw
1689  enddo
1690 
1691  endif lab_if_cld
1692  enddo lab_do_k
1693 
1694  else lab_if_ilwcliq
1695 
1696  do k = 1, nlay
1697  if (cfrac(k) > cldmin) then
1698  do ib = 1, nbands
1699  taucld(ib,k) = cdat1(k)
1700  enddo
1701  endif
1702  enddo
1703 
1704  endif lab_if_ilwcliq
1705 
1706 ! --- distribute cloud properties to each g-point
1707 
1708  if ( isubclw > 0 ) then ! mcica sub-col clouds approx
1709  do k = 1, nlay
1710  if ( cfrac(k) < cldmin ) then
1711  cldf(k) = f_zero
1712  else
1713  cldf(k) = cfrac(k)
1714  endif
1715  enddo
1716 
1717 ! --- ... call sub-column cloud generator
1718 
1719  call mcica_subcol &
1720 ! --- inputs:
1721  & ( cldf, nlay, ipseed, &
1722 ! --- output:
1723  & lcloudy &
1724  & )
1725 
1726  do k = 1, nlay
1727  do ig = 1, ngptlw
1728  if ( lcloudy(ig,k) ) then
1729  cldfmc(ig,k) = f_one
1730  else
1731  cldfmc(ig,k) = f_zero
1732  endif
1733  enddo
1734  enddo
1735 
1736  endif ! end if_isubclw_block
1737 
1738  return
1739 ! ..................................
1740  end subroutine cldprop
1741 ! ----------------------------------
1742 
1751 ! ----------------------------------
1752  subroutine mcica_subcol &
1753 ! ..................................
1754 ! --- inputs:
1755  & ( cldf, nlay, ipseed, &
1756 ! --- outputs:
1757  & lcloudy &
1758  & )
1759 
1760 ! ==================== defination of variables ==================== !
1761 ! !
1762 ! input variables: size !
1763 ! cldf - real, layer cloud fraction nlay !
1764 ! nlay - integer, number of model vertical layers 1 !
1765 ! ipseed - integer, permute seed for random num generator 1 !
1766 ! ** note : if the cloud generator is called multiple times, need !
1767 ! to permute the seed between each call; if between calls !
1768 ! for lw and sw, use values differ by the number of g-pts. !
1769 ! !
1770 ! output variables: !
1771 ! lcloudy - logical, sub-colum cloud profile flag array ngptlw*nlay!
1772 ! !
1773 ! other control flags from module variables: !
1774 ! iovrlw : control flag for cloud overlapping method !
1775 ! =0:random; =1:maximum/random: =2:maximum !
1776 ! !
1777 ! ===================== end of definitions ==================== !
1778 
1779  implicit none
1780 
1781 ! --- inputs:
1782  integer, intent(in) :: nlay, ipseed
1783 
1784  real (kind=kind_phys), dimension(nlay), intent(in) :: cldf
1785 
1786 ! --- outputs:
1787  logical, dimension(ngptlw,nlay), intent(out) :: lcloudy
1788 
1789 ! --- locals:
1790  real (kind=kind_phys) :: cdfunc(ngptlw,nlay), rand1d(ngptlw), &
1791  & rand2d(nlay*ngptlw), tem1
1792 
1793  type(random_stat) :: stat ! for thread safe random generator
1794 
1795  integer :: k, n, k1
1796 !
1797 !===> ... begin here
1798 !
1799 ! --- ... advance randum number generator by ipseed values
1800 
1801  call random_setseed &
1802 ! --- inputs:
1803  & ( ipseed, &
1804 ! --- outputs:
1805  & stat &
1806  & )
1807 
1808 ! --- ... sub-column set up according to overlapping assumption
1809 
1810  select case ( iovrlw )
1811 
1812  case( 0 ) ! random overlap, pick a random value at every level
1813 
1814  call random_number &
1815 ! --- inputs: ( none )
1816 ! --- outputs:
1817  & ( rand2d, stat )
1818 
1819  k1 = 0
1820  do n = 1, ngptlw
1821  do k = 1, nlay
1822  k1 = k1 + 1
1823  cdfunc(n,k) = rand2d(k1)
1824  enddo
1825  enddo
1826 
1827  case( 1 ) ! max-ran overlap
1828 
1829  call random_number &
1830 ! --- inputs: ( none )
1831 ! --- outputs:
1832  & ( rand2d, stat )
1833 
1834  k1 = 0
1835  do n = 1, ngptlw
1836  do k = 1, nlay
1837  k1 = k1 + 1
1838  cdfunc(n,k) = rand2d(k1)
1839  enddo
1840  enddo
1841 
1842 ! --- first pick a random number for bottom (or top) layer.
1843 ! then walk up the column: (aer's code)
1844 ! if layer below is cloudy, use the same rand num in the layer below
1845 ! if layer below is clear, use a new random number
1846 
1847 ! --- from bottom up
1848  do k = 2, nlay
1849  k1 = k - 1
1850  tem1 = f_one - cldf(k1)
1851 
1852  do n = 1, ngptlw
1853  if ( cdfunc(n,k1) > tem1 ) then
1854  cdfunc(n,k) = cdfunc(n,k1)
1855  else
1856  cdfunc(n,k) = cdfunc(n,k) * tem1
1857  endif
1858  enddo
1859  enddo
1860 
1861 ! --- or walk down the column: (if use original author's method)
1862 ! if layer above is cloudy, use the same rand num in the layer above
1863 ! if layer above is clear, use a new random number
1864 
1865 ! --- from top down
1866 ! do k = nlay-1, 1, -1
1867 ! k1 = k + 1
1868 ! tem1 = f_one - cldf(k1)
1869 
1870 ! do n = 1, ngptlw
1871 ! if ( cdfunc(n,k1) > tem1 ) then
1872 ! cdfunc(n,k) = cdfunc(n,k1)
1873 ! else
1874 ! cdfunc(n,k) = cdfunc(n,k) * tem1
1875 ! endif
1876 ! enddo
1877 ! enddo
1878 
1879  case( 2 ) ! maximum overlap, pick same random numebr at every level
1880 
1881  call random_number &
1882 ! --- inputs: ( none )
1883 ! --- outputs:
1884  & ( rand1d, stat )
1885 
1886  do n = 1, ngptlw
1887  tem1 = rand1d(n)
1888 
1889  do k = 1, nlay
1890  cdfunc(n,k) = tem1
1891  enddo
1892  enddo
1893 
1894  end select
1895 
1896 ! --- ... generate subcolumns for homogeneous clouds
1897 
1898  do k = 1, nlay
1899  tem1 = f_one - cldf(k)
1900 
1901  do n = 1, ngptlw
1902  lcloudy(n,k) = cdfunc(n,k) >= tem1
1903  enddo
1904  enddo
1905 
1906  return
1907 ! ..................................
1908  end subroutine mcica_subcol
1909 ! ----------------------------------
1910 
1913 ! ----------------------------------
1914  subroutine setcoef &
1915 ! ..................................
1916 ! --- inputs:
1917  & ( pavel,tavel,tz,stemp,h2ovmr,colamt,coldry,colbrd, &
1918  & nlay, nlp1, &
1919 ! --- outputs:
1920  & laytrop,pklay,pklev,jp,jt,jt1, &
1921  & rfrate,fac00,fac01,fac10,fac11, &
1922  & selffac,selffrac,indself,forfac,forfrac,indfor, &
1923  & minorfrac,scaleminor,scaleminorn2,indminor &
1924  & )
1925 
1926 ! =================== program usage description =================== !
1927 ! !
1928 ! purpose: compute various coefficients needed in radiative transfer !
1929 ! calculations. !
1930 ! !
1931 ! subprograms called: none !
1932 ! !
1933 ! ==================== defination of variables ==================== !
1934 ! !
1935 ! inputs: -size- !
1936 ! pavel - real, layer pressures (mb) nlay !
1937 ! tavel - real, layer temperatures (k) nlay !
1938 ! tz - real, level (interface) temperatures (k) 0:nlay !
1939 ! stemp - real, surface ground temperature (k) 1 !
1940 ! h2ovmr - real, layer w.v. volum mixing ratio (kg/kg) nlay !
1941 ! colamt - real, column amounts of absorbing gases nlay*maxgas!
1942 ! 2nd indices range: 1-maxgas, for watervapor, !
1943 ! carbon dioxide, ozone, nitrous oxide, methane, !
1944 ! oxigen, carbon monoxide,etc. (molecules/cm**2) !
1945 ! coldry - real, dry air column amount nlay !
1946 ! colbrd - real, column amount of broadening gases nlay !
1947 ! nlay/nlp1 - integer, total number of vertical layers, levels 1 !
1948 ! !
1949 ! outputs: !
1950 ! laytrop - integer, tropopause layer index (unitless) 1 !
1951 ! pklay - real, integrated planck func at lay temp nbands*0:nlay!
1952 ! pklev - real, integrated planck func at lev temp nbands*0:nlay!
1953 ! jp - real, indices of lower reference pressure nlay !
1954 ! jt, jt1 - real, indices of lower reference temperatures nlay !
1955 ! rfrate - real, ref ratios of binary species param nlay*nrates*2!
1956 ! (:,m,:)m=1-h2o/co2,2-h2o/o3,3-h2o/n2o,4-h2o/ch4,5-n2o/co2,6-o3/co2!
1957 ! (:,:,n)n=1,2: the rates of ref press at the 2 sides of the layer !
1958 ! facij - real, factors multiply the reference ks, nlay !
1959 ! i,j=0/1 for lower/higher of the 2 appropriate !
1960 ! temperatures and altitudes. !
1961 ! selffac - real, scale factor for w. v. self-continuum nlay !
1962 ! equals (w. v. density)/(atmospheric density !
1963 ! at 296k and 1013 mb) !
1964 ! selffrac - real, factor for temperature interpolation of nlay !
1965 ! reference w. v. self-continuum data !
1966 ! indself - integer, index of lower ref temp for selffac nlay !
1967 ! forfac - real, scale factor for w. v. foreign-continuum nlay !
1968 ! forfrac - real, factor for temperature interpolation of nlay !
1969 ! reference w.v. foreign-continuum data !
1970 ! indfor - integer, index of lower ref temp for forfac nlay !
1971 ! minorfrac - real, factor for minor gases nlay !
1972 ! scaleminor,scaleminorn2 !
1973 ! - real, scale factors for minor gases nlay !
1974 ! indminor - integer, index of lower ref temp for minor gases nlay !
1975 ! !
1976 ! ====================== end of definitions =================== !
1977 
1978 ! --- inputs:
1979  integer, intent(in) :: nlay, nlp1
1980 
1981  real (kind=kind_phys), dimension(nlay,maxgas),intent(in):: colamt
1982  real (kind=kind_phys), dimension(0:nlay), intent(in):: tz
1983 
1984  real (kind=kind_phys), dimension(nlay), intent(in) :: pavel, &
1985  & tavel, h2ovmr, coldry, colbrd
1986 
1987  real (kind=kind_phys), intent(in) :: stemp
1988 
1989 ! --- outputs:
1990  integer, dimension(nlay), intent(out) :: jp, jt, jt1, indself, &
1991  & indfor, indminor
1992 
1993  integer, intent(out) :: laytrop
1994 
1995  real (kind=kind_phys), dimension(nlay,nrates,2), intent(out) :: &
1996  & rfrate
1997  real (kind=kind_phys), dimension(nbands,0:nlay), intent(out) :: &
1998  & pklev, pklay
1999 
2000  real (kind=kind_phys), dimension(nlay), intent(out) :: &
2001  & fac00, fac01, fac10, fac11, selffac, selffrac, forfac, &
2002  & forfrac, minorfrac, scaleminor, scaleminorn2
2003 
2004 ! --- locals:
2005  real (kind=kind_phys) :: tlvlfr, tlyrfr, plog, fp, ft, ft1, &
2006  & tem1, tem2
2007 
2008  integer :: i, k, jp1, indlev, indlay
2009 !
2010 !===> ... begin here
2011 !
2012 ! --- ... calculate information needed by the radiative transfer routine
2013 ! that is specific to this atmosphere, especially some of the
2014 ! coefficients and indices needed to compute the optical depths
2015 ! by interpolating data from stored reference atmospheres.
2016 
2017  indlay = min(180, max(1, int(stemp-159.0) ))
2018  indlev = min(180, max(1, int(tz(0)-159.0) ))
2019  tlyrfr = stemp - int(stemp)
2020  tlvlfr = tz(0) - int(tz(0))
2021  do i = 1, nbands
2022  tem1 = totplnk(indlay+1,i) - totplnk(indlay,i)
2023  tem2 = totplnk(indlev+1,i) - totplnk(indlev,i)
2024  pklay(i,0) = delwave(i) * (totplnk(indlay,i) + tlyrfr*tem1)
2025  pklev(i,0) = delwave(i) * (totplnk(indlev,i) + tlvlfr*tem2)
2026  enddo
2027 
2028 ! --- ... begin layer loop
2029 ! calculate the integrated Planck functions for each band at the
2030 ! surface, level, and layer temperatures.
2031 
2032  laytrop = 0
2033 
2034  do k = 1, nlay
2035 
2036  indlay = min(180, max(1, int(tavel(k)-159.0) ))
2037  tlyrfr = tavel(k) - int(tavel(k))
2038 
2039  indlev = min(180, max(1, int(tz(k)-159.0) ))
2040  tlvlfr = tz(k) - int(tz(k))
2041 
2042 ! --- ... begin spectral band loop
2043 
2044  do i = 1, nbands
2045  pklay(i,k) = delwave(i) * (totplnk(indlay,i) + tlyrfr &
2046  & * (totplnk(indlay+1,i) - totplnk(indlay,i)) )
2047  pklev(i,k) = delwave(i) * (totplnk(indlev,i) + tlvlfr &
2048  & * (totplnk(indlev+1,i) - totplnk(indlev,i)) )
2049  enddo
2050 
2051 ! --- ... find the two reference pressures on either side of the
2052 ! layer pressure. store them in jp and jp1. store in fp the
2053 ! fraction of the difference (in ln(pressure)) between these
2054 ! two values that the layer pressure lies.
2055 
2056  plog = log(pavel(k))
2057  jp(k)= max(1, min(58, int(36.0 - 5.0*(plog+0.04)) ))
2058  jp1 = jp(k) + 1
2059 ! --- ... limit pressure extrapolation at the top
2060  fp = max(f_zero, min(f_one, 5.0*(preflog(jp(k))-plog) ))
2061 !org fp = 5.0 * (preflog(jp(k)) - plog)
2062 
2063 ! --- ... determine, for each reference pressure (jp and jp1), which
2064 ! reference temperature (these are different for each
2065 ! reference pressure) is nearest the layer temperature but does
2066 ! not exceed it. store these indices in jt and jt1, resp.
2067 ! store in ft (resp. ft1) the fraction of the way between jt
2068 ! (jt1) and the next highest reference temperature that the
2069 ! layer temperature falls.
2070 
2071  tem1 = (tavel(k)-tref(jp(k))) / 15.0
2072  tem2 = (tavel(k)-tref(jp1 )) / 15.0
2073  jt(k) = max(1, min(4, int(3.0 + tem1) ))
2074  jt1(k) = max(1, min(4, int(3.0 + tem2) ))
2075 ! --- ... restrict extrapolation ranges by limiting abs(det t) < 37.5 deg
2076  ft = max(-0.5, min(1.5, tem1 - float(jt(k) - 3) ))
2077  ft1 = max(-0.5, min(1.5, tem2 - float(jt1(k) - 3) ))
2078 !org ft = tem1 - float(jt (k) - 3)
2079 !org ft1 = tem2 - float(jt1(k) - 3)
2080 
2081 ! --- ... we have now isolated the layer ln pressure and temperature,
2082 ! between two reference pressures and two reference temperatures
2083 ! (for each reference pressure). we multiply the pressure
2084 ! fraction fp with the appropriate temperature fractions to get
2085 ! the factors that will be needed for the interpolation that yields
2086 ! the optical depths (performed in routines taugbn for band n)
2087 
2088  tem1 = f_one - fp
2089  fac10(k) = tem1 * ft
2090  fac00(k) = tem1 * (f_one - ft)
2091  fac11(k) = fp * ft1
2092  fac01(k) = fp * (f_one - ft1)
2093 
2094  forfac(k) = pavel(k)*stpfac / (tavel(k)*(1.0 + h2ovmr(k)))
2095  selffac(k) = h2ovmr(k) * forfac(k)
2096 
2097 ! --- ... set up factors needed to separately include the minor gases
2098 ! in the calculation of absorption coefficient
2099 
2100  scaleminor(k) = pavel(k) / tavel(k)
2101  scaleminorn2(k) = (pavel(k) / tavel(k)) &
2102  & * (colbrd(k)/(coldry(k) + colamt(k,1)))
2103  tem1 = (tavel(k) - 180.8) / 7.2
2104  indminor(k) = min(18, max(1, int(tem1)))
2105  minorfrac(k) = tem1 - float(indminor(k))
2106 
2107 ! --- ... if the pressure is less than ~100mb, perform a different
2108 ! set of species interpolations.
2109 
2110  if (plog > 4.56) then
2111 
2112  laytrop = laytrop + 1
2113 
2114  tem1 = (332.0 - tavel(k)) / 36.0
2115  indfor(k) = min(2, max(1, int(tem1)))
2116  forfrac(k) = tem1 - float(indfor(k))
2117 
2118 ! --- ... set up factors needed to separately include the water vapor
2119 ! self-continuum in the calculation of absorption coefficient.
2120 
2121  tem1 = (tavel(k) - 188.0) / 7.2
2122  indself(k) = min(9, max(1, int(tem1)-7))
2123  selffrac(k) = tem1 - float(indself(k) + 7)
2124 
2125 ! --- ... setup reference ratio to be used in calculation of binary
2126 ! species parameter in lower atmosphere.
2127 
2128  rfrate(k,1,1) = chi_mls(1,jp(k)) / chi_mls(2,jp(k))
2129  rfrate(k,1,2) = chi_mls(1,jp(k)+1) / chi_mls(2,jp(k)+1)
2130 
2131  rfrate(k,2,1) = chi_mls(1,jp(k)) / chi_mls(3,jp(k))
2132  rfrate(k,2,2) = chi_mls(1,jp(k)+1) / chi_mls(3,jp(k)+1)
2133 
2134  rfrate(k,3,1) = chi_mls(1,jp(k)) / chi_mls(4,jp(k))
2135  rfrate(k,3,2) = chi_mls(1,jp(k)+1) / chi_mls(4,jp(k)+1)
2136 
2137  rfrate(k,4,1) = chi_mls(1,jp(k)) / chi_mls(6,jp(k))
2138  rfrate(k,4,2) = chi_mls(1,jp(k)+1) / chi_mls(6,jp(k)+1)
2139 
2140  rfrate(k,5,1) = chi_mls(4,jp(k)) / chi_mls(2,jp(k))
2141  rfrate(k,5,2) = chi_mls(4,jp(k)+1) / chi_mls(2,jp(k)+1)
2142 
2143  else
2144 
2145  tem1 = (tavel(k) - 188.0) / 36.0
2146  indfor(k) = 3
2147  forfrac(k) = tem1 - f_one
2148 
2149  indself(k) = 0
2150  selffrac(k) = f_zero
2151 
2152 ! --- ... setup reference ratio to be used in calculation of binary
2153 ! species parameter in upper atmosphere.
2154 
2155  rfrate(k,1,1) = chi_mls(1,jp(k)) / chi_mls(2,jp(k))
2156  rfrate(k,1,2) = chi_mls(1,jp(k)+1) / chi_mls(2,jp(k)+1)
2157 
2158  rfrate(k,6,1) = chi_mls(3,jp(k)) / chi_mls(2,jp(k))
2159  rfrate(k,6,2) = chi_mls(3,jp(k)+1) / chi_mls(2,jp(k)+1)
2160 
2161  endif
2162 
2163 ! --- ... rescale selffac and forfac for use in taumol
2164 
2165  selffac(k) = colamt(k,1) * selffac(k)
2166  forfac(k) = colamt(k,1) * forfac(k)
2167 
2168  enddo ! end do_k layer loop
2169 
2170  return
2171 ! ..................................
2172  end subroutine setcoef
2173 ! ----------------------------------
2174 
2203 ! ----------------------------------
2204  subroutine rtrn &
2205 ! ..................................
2206 ! --- inputs:
2207  & ( semiss,delp,cldfrc,taucld,tautot,pklay,pklev, &
2208  & fracs,secdif, nlay,nlp1, &
2209 ! --- outputs:
2210  & totuflux,totdflux,htr, totuclfl,totdclfl,htrcl, htrb &
2211  & )
2212 
2213 ! =================== program usage description =================== !
2214 ! !
2215 ! purpose: compute the upward/downward radiative fluxes, and heating !
2216 ! rates for both clear or cloudy atmosphere. clouds are assumed as !
2217 ! randomly overlaping in a vertical colum. !
2218 ! !
2219 ! subprograms called: none !
2220 ! !
2221 ! ==================== defination of variables ==================== !
2222 ! !
2223 ! inputs: -size- !
2224 ! semiss - real, lw surface emissivity nbands!
2225 ! delp - real, layer pressure thickness (mb) nlay !
2226 ! cldfrc - real, layer cloud fraction 0:nlp1 !
2227 ! taucld - real, layer cloud opt depth nbands,nlay!
2228 ! tautot - real, total optical depth (gas+aerosols) ngptlw,nlay!
2229 ! pklay - real, integrated planck func at lay temp nbands*0:nlay!
2230 ! pklev - real, integrated planck func at lev temp nbands*0:nlay!
2231 ! fracs - real, planck fractions ngptlw,nlay!
2232 ! secdif - real, secant of diffusivity angle nbands!
2233 ! nlay - integer, number of vertical layers 1 !
2234 ! nlp1 - integer, number of vertical levels (interfaces) 1 !
2235 ! !
2236 ! outputs: !
2237 ! totuflux- real, total sky upward flux (w/m2) 0:nlay !
2238 ! totdflux- real, total sky downward flux (w/m2) 0:nlay !
2239 ! htr - real, total sky heating rate (k/sec or k/day) nlay !
2240 ! totuclfl- real, clear sky upward flux (w/m2) 0:nlay !
2241 ! totdclfl- real, clear sky downward flux (w/m2) 0:nlay !
2242 ! htrcl - real, clear sky heating rate (k/sec or k/day) nlay !
2243 ! htrb - real, spectral band lw heating rate (k/day) nlay*nbands!
2244 ! !
2245 ! module veriables: !
2246 ! ngb - integer, band index for each g-value ngptlw!
2247 ! fluxfac - real, conversion factor for fluxes (pi*2.e4) 1 !
2248 ! heatfac - real, conversion factor for heating rates (g/cp*1e-2) 1 !
2249 ! tblint - real, conversion factor for look-up tbl (float(ntbl) 1 !
2250 ! bpade - real, pade approx constant (1/0.278) 1 !
2251 ! wtdiff - real, weight for radiance to flux conversion 1 !
2252 ! ntbl - integer, dimension of look-up tables 1 !
2253 ! tau_tbl - real, clr-sky opt dep lookup table 0:ntbl !
2254 ! exp_tbl - real, transmittance lookup table 0:ntbl !
2255 ! tfn_tbl - real, tau transition function 0:ntbl !
2256 ! !
2257 ! local variables: !
2258 ! itgas - integer, index for gases contribution look-up table 1 !
2259 ! ittot - integer, index for gases plus clouds look-up table 1 !
2260 ! reflct - real, surface reflectance 1 !
2261 ! atrgas - real, gaseous absorptivity 1 !
2262 ! atrtot - real, gaseous and cloud absorptivity 1 !
2263 ! odcld - real, cloud optical depth 1 !
2264 ! efclrfr- real, effective clear sky fraction (1-efcldfr) nlay !
2265 ! odepth - real, optical depth of gaseous only 1 !
2266 ! odtot - real, optical depth of gas and cloud 1 !
2267 ! gasfac - real, gas-only pade factor, used for planck fn 1 !
2268 ! totfac - real, gas+cld pade factor, used for planck fn 1 !
2269 ! bbdgas - real, gas-only planck function for downward rt 1 !
2270 ! bbugas - real, gas-only planck function for upward rt 1 !
2271 ! bbdtot - real, gas and cloud planck function for downward rt 1 !
2272 ! bbutot - real, gas and cloud planck function for upward rt 1 !
2273 ! gassrcu- real, upwd source radiance due to gas only nlay!
2274 ! totsrcu- real, upwd source radiance due to gas+cld nlay!
2275 ! gassrcd- real, dnwd source radiance due to gas only 1 !
2276 ! totsrcd- real, dnwd source radiance due to gas+cld 1 !
2277 ! radtotu- real, spectrally summed total sky upwd radiance 1 !
2278 ! radclru- real, spectrally summed clear sky upwd radiance 1 !
2279 ! radtotd- real, spectrally summed total sky dnwd radiance 1 !
2280 ! radclrd- real, spectrally summed clear sky dnwd radiance 1 !
2281 ! toturad- real, total sky upward radiance by layer 0:nlay*nbands!
2282 ! clrurad- real, clear sky upward radiance by layer 0:nlay*nbands!
2283 ! totdrad- real, total sky downward radiance by layer 0:nlay*nbands!
2284 ! clrdrad- real, clear sky downward radiance by layer 0:nlay*nbands!
2285 ! fnet - real, net longwave flux (w/m2) 0:nlay !
2286 ! fnetc - real, clear sky net longwave flux (w/m2) 0:nlay !
2287 ! !
2288 ! !
2289 ! ******************************************************************* !
2290 ! original code description !
2291 ! !
2292 ! original version: e. j. mlawer, et al. rrtm_v3.0 !
2293 ! revision for gcms: michael j. iacono; october, 2002 !
2294 ! revision for f90: michael j. iacono; june, 2006 !
2295 ! !
2296 ! this program calculates the upward fluxes, downward fluxes, and !
2297 ! heating rates for an arbitrary clear or cloudy atmosphere. the input !
2298 ! to this program is the atmospheric profile, all Planck function !
2299 ! information, and the cloud fraction by layer. a variable diffusivity!
2300 ! angle (secdif) is used for the angle integration. bands 2-3 and 5-9 !
2301 ! use a value for secdif that varies from 1.50 to 1.80 as a function !
2302 ! of the column water vapor, and other bands use a value of 1.66. the !
2303 ! gaussian weight appropriate to this angle (wtdiff=0.5) is applied !
2304 ! here. note that use of the emissivity angle for the flux integration!
2305 ! can cause errors of 1 to 4 W/m2 within cloudy layers. !
2306 ! clouds are treated with a random cloud overlap method. !
2307 ! !
2308 ! ******************************************************************* !
2309 ! ====================== end of description block ================= !
2310 
2311 ! --- inputs:
2312  integer, intent(in) :: nlay, nlp1
2313 
2314  real (kind=kind_phys), dimension(0:nlp1), intent(in) :: cldfrc
2315  real (kind=kind_phys), dimension(nbands), intent(in) :: semiss, &
2316  & secdif
2317  real (kind=kind_phys), dimension(nlay), intent(in) :: delp
2318 
2319  real (kind=kind_phys), dimension(nbands,nlay),intent(in):: taucld
2320  real (kind=kind_phys), dimension(ngptlw,nlay),intent(in):: fracs, &
2321  & tautot
2322 
2323  real (kind=kind_phys), dimension(nbands,0:nlay), intent(in) :: &
2324  & pklev, pklay
2325 
2326 ! --- outputs:
2327  real (kind=kind_phys), dimension(nlay), intent(out) :: htr, htrcl
2328 
2329  real (kind=kind_phys), dimension(nlay,nbands),intent(out) :: htrb
2330 
2331  real (kind=kind_phys), dimension(0:nlay), intent(out) :: &
2332  & totuflux, totdflux, totuclfl, totdclfl
2333 
2334 ! --- locals:
2335  real (kind=kind_phys), parameter :: rec_6 = 0.166667
2336 
2337  real (kind=kind_phys), dimension(0:nlay,nbands) :: clrurad, &
2338  & clrdrad, toturad, totdrad
2339 
2340  real (kind=kind_phys), dimension(nlay) :: gassrcu, totsrcu, &
2341  & trngas, efclrfr, rfdelp
2342  real (kind=kind_phys), dimension(0:nlay) :: fnet, fnetc
2343 
2344  real (kind=kind_phys) :: totsrcd, gassrcd, tblind, odepth, odtot, &
2345  & odcld, atrtot, atrgas, reflct, totfac, gasfac, flxfac, &
2346  & plfrac, blay, bbdgas, bbdtot, bbugas, bbutot, dplnku, &
2347  & dplnkd, radtotu, radclru, radtotd, radclrd, rad0, &
2348  & clfr, trng, gasu
2349 
2350  integer :: ittot, itgas, ib, ig, k
2351 !
2352 !===> ... begin here
2353 !
2354  do ib = 1, nbands
2355  do k = 0, nlay
2356  toturad(k,ib) = f_zero
2357  totdrad(k,ib) = f_zero
2358  clrurad(k,ib) = f_zero
2359  clrdrad(k,ib) = f_zero
2360  enddo
2361  enddo
2362 
2363  do k = 0, nlay
2364  totuflux(k) = f_zero
2365  totdflux(k) = f_zero
2366  totuclfl(k) = f_zero
2367  totdclfl(k) = f_zero
2368  enddo
2369 
2370 ! --- ... loop over all g-points
2371 
2372  do ig = 1, ngptlw
2373  ib = ngb(ig)
2374 
2375  radtotd = f_zero
2376  radclrd = f_zero
2377 
2378 ! --- ... downward radiative transfer loop.
2379 
2380  do k = nlay, 1, -1
2381 
2382 ! --- ... clear sky, gases contribution
2383 
2384  odepth = max( f_zero, secdif(ib)*tautot(ig,k) )
2385  if (odepth <= 0.06) then
2386  atrgas = odepth - 0.5*odepth*odepth
2387  trng = f_one - atrgas
2388  gasfac = rec_6 * odepth
2389  else
2390  tblind = odepth / (bpade + odepth)
2391  itgas = tblint*tblind + 0.5
2392  trng = exp_tbl(itgas)
2393  atrgas = f_one - trng
2394  gasfac = tfn_tbl(itgas)
2395  odepth = tau_tbl(itgas)
2396  endif
2397 
2398  plfrac = fracs(ig,k)
2399  blay = pklay(ib,k)
2400 
2401  dplnku = pklev(ib,k ) - blay
2402  dplnkd = pklev(ib,k-1) - blay
2403  bbdgas = plfrac * (blay + dplnkd*gasfac)
2404  bbugas = plfrac * (blay + dplnku*gasfac)
2405  gassrcd= bbdgas * atrgas
2406  gassrcu(k)= bbugas * atrgas
2407  trngas(k) = trng
2408 
2409 ! --- ... total sky, gases+clouds contribution
2410 
2411  clfr = cldfrc(k)
2412  if (clfr >= eps) then
2413 ! --- ... cloudy layer
2414 
2415  odcld = secdif(ib) * taucld(ib,k)
2416  efclrfr(k) = f_one-(f_one - exp(-odcld))*clfr
2417  odtot = odepth + odcld
2418  if (odtot < 0.06) then
2419  totfac = rec_6 * odtot
2420  atrtot = odtot - 0.5*odtot*odtot
2421  else
2422  tblind = odtot / (bpade + odtot)
2423  ittot = tblint*tblind + 0.5
2424  totfac = tfn_tbl(ittot)
2425  atrtot = f_one - exp_tbl(ittot)
2426  endif
2427 
2428  bbdtot = plfrac * (blay + dplnkd*totfac)
2429  bbutot = plfrac * (blay + dplnku*totfac)
2430  totsrcd= bbdtot * atrtot
2431  totsrcu(k)= bbutot * atrtot
2432 
2433 ! --- ... total sky radiance
2434  radtotd = radtotd*trng*efclrfr(k) + gassrcd &
2435  & + clfr*(totsrcd - gassrcd)
2436  totdrad(k-1,ib) = totdrad(k-1,ib) + radtotd
2437 
2438 ! --- ... clear sky radiance
2439  radclrd = radclrd*trng + gassrcd
2440  clrdrad(k-1,ib) = clrdrad(k-1,ib) + radclrd
2441 
2442  else
2443 ! --- ... clear layer
2444 
2445 ! --- ... total sky radiance
2446  radtotd = radtotd*trng + gassrcd
2447  totdrad(k-1,ib) = totdrad(k-1,ib) + radtotd
2448 
2449 ! --- ... clear sky radiance
2450  radclrd = radclrd*trng + gassrcd
2451  clrdrad(k-1,ib) = clrdrad(k-1,ib) + radclrd
2452 
2453  endif ! end if_clfr_block
2454 
2455  enddo ! end do_k_loop
2456 
2457 ! --- ... spectral emissivity & reflectance
2458 ! include the contribution of spectrally varying longwave emissivity
2459 ! and reflection from the surface to the upward radiative transfer.
2460 ! note: spectral and Lambertian reflection are identical for the
2461 ! diffusivity angle flux integration used here.
2462 
2463  reflct = f_one - semiss(ib)
2464  rad0 = semiss(ib) * fracs(ig,1) * pklay(ib,0)
2465 
2466 ! --- ... total sky radiance
2467  radtotu = rad0 + reflct*radtotd
2468  toturad(0,ib) = toturad(0,ib) + radtotu
2469 
2470 ! --- ... clear sky radiance
2471  radclru = rad0 + reflct*radclrd
2472  clrurad(0,ib) = clrurad(0,ib) + radclru
2473 
2474 ! --- ... upward radiative transfer
2475 
2476  do k = 1, nlay
2477  clfr = cldfrc(k)
2478  trng = trngas(k)
2479  gasu = gassrcu(k)
2480 
2481  if (clfr >= eps) then
2482 ! --- ... cloudy layer
2483 
2484 ! --- ... total sky radiance
2485  radtotu = radtotu*trng*efclrfr(k) + gasu &
2486  & + clfr*(totsrcu(k) - gasu)
2487  toturad(k,ib) = toturad(k,ib) + radtotu
2488 
2489 ! --- ... clear sky radiance
2490  radclru = radclru*trng + gasu
2491  clrurad(k,ib) = clrurad(k,ib) + radclru
2492 
2493  else
2494 ! --- ... clear layer
2495 
2496 ! --- ... total sky radiance
2497  radtotu = radtotu*trng + gasu
2498  toturad(k,ib) = toturad(k,ib) + radtotu
2499 
2500 ! --- ... clear sky radiance
2501  radclru = radclru*trng + gasu
2502  clrurad(k,ib) = clrurad(k,ib) + radclru
2503 
2504  endif ! end if_clfr_block
2505 
2506  enddo ! end do_k_loop
2507 
2508  enddo ! end do_ig_loop
2509 
2510 ! --- ... process longwave output from band for total and clear streams.
2511 ! calculate upward, downward, and net flux.
2512 
2513  flxfac = wtdiff * fluxfac
2514 
2515  do k = 0, nlay
2516  do ib = 1, nbands
2517  totuflux(k) = totuflux(k) + toturad(k,ib)
2518  totdflux(k) = totdflux(k) + totdrad(k,ib)
2519  totuclfl(k) = totuclfl(k) + clrurad(k,ib)
2520  totdclfl(k) = totdclfl(k) + clrdrad(k,ib)
2521  enddo
2522 
2523  totuflux(k) = totuflux(k) * flxfac
2524  totdflux(k) = totdflux(k) * flxfac
2525  totuclfl(k) = totuclfl(k) * flxfac
2526  totdclfl(k) = totdclfl(k) * flxfac
2527  enddo
2528 
2529 ! --- ... calculate net fluxes and heating rates
2530  fnet(0) = totuflux(0) - totdflux(0)
2531 
2532  do k = 1, nlay
2533  rfdelp(k) = heatfac / delp(k)
2534  fnet(k) = totuflux(k) - totdflux(k)
2535  htr(k) = (fnet(k-1) - fnet(k)) * rfdelp(k)
2536  enddo
2537 
2538 !! --- ... optional clear sky heating rates
2539  if ( lhlw0 ) then
2540  fnetc(0) = totuclfl(0) - totdclfl(0)
2541 
2542  do k = 1, nlay
2543  fnetc(k) = totuclfl(k) - totdclfl(k)
2544  htrcl(k) = (fnetc(k-1) - fnetc(k)) * rfdelp(k)
2545  enddo
2546  endif
2547 
2548 !! --- ... optional spectral band heating rates
2549  if ( lhlwb ) then
2550  do ib = 1, nbands
2551  fnet(0) = (toturad(0,ib) - totdrad(0,ib)) * flxfac
2552 
2553  do k = 1, nlay
2554  fnet(k) = (toturad(k,ib) - totdrad(k,ib)) * flxfac
2555  htrb(k,ib) = (fnet(k-1) - fnet(k)) * rfdelp(k)
2556  enddo
2557  enddo
2558  endif
2559 
2560 ! ..................................
2561  end subroutine rtrn
2562 ! ----------------------------------
2563 
2564 
2565 ! ----------------------------------
2566  subroutine rtrnmr &
2567 ! ..................................
2568 ! --- inputs:
2569  & ( semiss,delp,cldfrc,taucld,tautot,pklay,pklev, &
2570  & fracs,secdif, nlay,nlp1, &
2571 ! --- outputs:
2572  & totuflux,totdflux,htr, totuclfl,totdclfl,htrcl, htrb &
2573  & )
2574 
2575 ! =================== program usage description =================== !
2576 ! !
2577 ! purpose: compute the upward/downward radiative fluxes, and heating !
2578 ! rates for both clear or cloudy atmosphere. clouds are assumed as in !
2579 ! maximum-randomly overlaping in a vertical colum. !
2580 ! !
2581 ! subprograms called: none !
2582 ! !
2583 ! ==================== defination of variables ==================== !
2584 ! !
2585 ! inputs: -size- !
2586 ! semiss - real, lw surface emissivity nbands!
2587 ! delp - real, layer pressure thickness (mb) nlay !
2588 ! cldfrc - real, layer cloud fraction 0:nlp1 !
2589 ! taucld - real, layer cloud opt depth nbands,nlay!
2590 ! tautot - real, total optical depth (gas+aerosols) ngptlw,nlay!
2591 ! pklay - real, integrated planck func at lay temp nbands*0:nlay!
2592 ! pklev - real, integrated planck func at lev temp nbands*0:nlay!
2593 ! fracs - real, planck fractions ngptlw,nlay!
2594 ! secdif - real, secant of diffusivity angle nbands!
2595 ! nlay - integer, number of vertical layers 1 !
2596 ! nlp1 - integer, number of vertical levels (interfaces) 1 !
2597 ! !
2598 ! outputs: !
2599 ! totuflux- real, total sky upward flux (w/m2) 0:nlay !
2600 ! totdflux- real, total sky downward flux (w/m2) 0:nlay !
2601 ! htr - real, total sky heating rate (k/sec or k/day) nlay !
2602 ! totuclfl- real, clear sky upward flux (w/m2) 0:nlay !
2603 ! totdclfl- real, clear sky downward flux (w/m2) 0:nlay !
2604 ! htrcl - real, clear sky heating rate (k/sec or k/day) nlay !
2605 ! htrb - real, spectral band lw heating rate (k/day) nlay*nbands!
2606 ! !
2607 ! module veriables: !
2608 ! ngb - integer, band index for each g-value ngptlw!
2609 ! fluxfac - real, conversion factor for fluxes (pi*2.e4) 1 !
2610 ! heatfac - real, conversion factor for heating rates (g/cp*1e-2) 1 !
2611 ! tblint - real, conversion factor for look-up tbl (float(ntbl) 1 !
2612 ! bpade - real, pade approx constant (1/0.278) 1 !
2613 ! wtdiff - real, weight for radiance to flux conversion 1 !
2614 ! ntbl - integer, dimension of look-up tables 1 !
2615 ! tau_tbl - real, clr-sky opt dep lookup table 0:ntbl !
2616 ! exp_tbl - real, transmittance lookup table 0:ntbl !
2617 ! tfn_tbl - real, tau transition function 0:ntbl !
2618 ! !
2619 ! local variables: !
2620 ! itgas - integer, index for gases contribution look-up table 1 !
2621 ! ittot - integer, index for gases plus clouds look-up table 1 !
2622 ! reflct - real, surface reflectance 1 !
2623 ! atrgas - real, gaseous absorptivity 1 !
2624 ! atrtot - real, gaseous and cloud absorptivity 1 !
2625 ! odcld - real, cloud optical depth 1 !
2626 ! odepth - real, optical depth of gaseous only 1 !
2627 ! odtot - real, optical depth of gas and cloud 1 !
2628 ! gasfac - real, gas-only pade factor, used for planck fn 1 !
2629 ! totfac - real, gas+cld pade factor, used for planck fn 1 !
2630 ! bbdgas - real, gas-only planck function for downward rt 1 !
2631 ! bbugas - real, gas-only planck function for upward rt 1 !
2632 ! bbdtot - real, gas and cloud planck function for downward rt 1 !
2633 ! bbutot - real, gas and cloud planck function for upward rt 1 !
2634 ! gassrcu- real, upwd source radiance due to gas only nlay!
2635 ! totsrcu- real, upwd source radiance due to gas + cld nlay!
2636 ! gassrcd- real, dnwd source radiance due to gas only 1 !
2637 ! totsrcd- real, dnwd source radiance due to gas + cld 1 !
2638 ! radtotu- real, spectrally summed total sky upwd radiance 1 !
2639 ! radclru- real, spectrally summed clear sky upwd radiance 1 !
2640 ! radtotd- real, spectrally summed total sky dnwd radiance 1 !
2641 ! radclrd- real, spectrally summed clear sky dnwd radiance 1 !
2642 ! toturad- real, total sky upward radiance by layer 0:nlay*nbands!
2643 ! clrurad- real, clear sky upward radiance by layer 0:nlay*nbands!
2644 ! totdrad- real, total sky downward radiance by layer 0:nlay*nbands!
2645 ! clrdrad- real, clear sky downward radiance by layer 0:nlay*nbands!
2646 ! fnet - real, net longwave flux (w/m2) 0:nlay !
2647 ! fnetc - real, clear sky net longwave flux (w/m2) 0:nlay !
2648 ! !
2649 ! !
2650 ! ******************************************************************* !
2651 ! original code description !
2652 ! !
2653 ! original version: e. j. mlawer, et al. rrtm_v3.0 !
2654 ! revision for gcms: michael j. iacono; october, 2002 !
2655 ! revision for f90: michael j. iacono; june, 2006 !
2656 ! !
2657 ! this program calculates the upward fluxes, downward fluxes, and !
2658 ! heating rates for an arbitrary clear or cloudy atmosphere. the input !
2659 ! to this program is the atmospheric profile, all Planck function !
2660 ! information, and the cloud fraction by layer. a variable diffusivity!
2661 ! angle (secdif) is used for the angle integration. bands 2-3 and 5-9 !
2662 ! use a value for secdif that varies from 1.50 to 1.80 as a function !
2663 ! of the column water vapor, and other bands use a value of 1.66. the !
2664 ! gaussian weight appropriate to this angle (wtdiff=0.5) is applied !
2665 ! here. note that use of the emissivity angle for the flux integration!
2666 ! can cause errors of 1 to 4 W/m2 within cloudy layers. !
2667 ! clouds are treated with a maximum-random cloud overlap method. !
2668 ! !
2669 ! ******************************************************************* !
2670 ! ====================== end of description block ================= !
2671 
2672 ! --- inputs:
2673  integer, intent(in) :: nlay, nlp1
2674 
2675  real (kind=kind_phys), dimension(0:nlp1), intent(in) :: cldfrc
2676  real (kind=kind_phys), dimension(nbands), intent(in) :: semiss, &
2677  & secdif
2678  real (kind=kind_phys), dimension(nlay), intent(in) :: delp
2679 
2680  real (kind=kind_phys), dimension(nbands,nlay),intent(in):: taucld
2681  real (kind=kind_phys), dimension(ngptlw,nlay),intent(in):: fracs, &
2682  & tautot
2683 
2684  real (kind=kind_phys), dimension(nbands,0:nlay), intent(in) :: &
2685  & pklev, pklay
2686 
2687 ! --- outputs:
2688  real (kind=kind_phys), dimension(nlay), intent(out) :: htr, htrcl
2689 
2690  real (kind=kind_phys), dimension(nlay,nbands),intent(out) :: htrb
2691 
2692  real (kind=kind_phys), dimension(0:nlay), intent(out) :: &
2693  & totuflux, totdflux, totuclfl, totdclfl
2694 
2695 ! --- locals:
2696  real (kind=kind_phys), parameter :: rec_6 = 0.166667
2697 
2698  real (kind=kind_phys), dimension(0:nlay,nbands) :: clrurad, &
2699  & clrdrad, toturad, totdrad
2700 
2701  real (kind=kind_phys), dimension(nlay) :: gassrcu, totsrcu, &
2702  & trngas, trntot, rfdelp
2703  real (kind=kind_phys), dimension(0:nlay) :: fnet, fnetc
2704 
2705  real (kind=kind_phys) :: totsrcd, gassrcd, tblind, odepth, odtot, &
2706  & odcld, atrtot, atrgas, reflct, totfac, gasfac, flxfac, &
2707  & plfrac, blay, bbdgas, bbdtot, bbugas, bbutot, dplnku, &
2708  & dplnkd, radtotu, radclru, radtotd, radclrd, rad0, rad, &
2709  & totradd, clrradd, totradu, clrradu, fmax, fmin, rat1, rat2,&
2710  & radmod, clfr, trng, trnt, gasu, totu
2711 
2712  integer :: ittot, itgas, ib, ig, k
2713 
2714 ! dimensions for cloud overlap adjustment
2715  real (kind=kind_phys), dimension(nlp1) :: faccld1u, faccld2u, &
2716  & facclr1u, facclr2u, faccmb1u, faccmb2u
2717  real (kind=kind_phys), dimension(0:nlay) :: faccld1d, faccld2d, &
2718  & facclr1d, facclr2d, faccmb1d, faccmb2d
2719 
2720  logical :: lstcldu(nlay), lstcldd(nlay)
2721 !
2722 !===> ... begin here
2723 !
2724  do k = 1, nlp1
2725  faccld1u(k) = f_zero
2726  faccld2u(k) = f_zero
2727  facclr1u(k) = f_zero
2728  facclr2u(k) = f_zero
2729  faccmb1u(k) = f_zero
2730  faccmb2u(k) = f_zero
2731  enddo
2732 
2733  lstcldu(1) = cldfrc(1) > eps
2734  rat1 = f_zero
2735  rat2 = f_zero
2736 
2737  do k = 1, nlay-1
2738 
2739  lstcldu(k+1) = cldfrc(k+1)>eps .and. cldfrc(k)<=eps
2740 
2741  if (cldfrc(k) > eps) then
2742 ! --- ... maximum/random cloud overlap
2743 
2744  if (cldfrc(k+1) >= cldfrc(k)) then
2745  if (lstcldu(k)) then
2746  if (cldfrc(k) < f_one) then
2747  facclr2u(k+1) = (cldfrc(k+1) - cldfrc(k)) &
2748  & / (f_one - cldfrc(k))
2749  endif
2750  facclr2u(k) = f_zero
2751  faccld2u(k) = f_zero
2752  else
2753  fmax = max(cldfrc(k), cldfrc(k-1))
2754  if (cldfrc(k+1) > fmax) then
2755  facclr1u(k+1) = rat2
2756  facclr2u(k+1) = (cldfrc(k+1) - fmax)/(f_one - fmax)
2757  elseif (cldfrc(k+1) < fmax) then
2758  facclr1u(k+1) = (cldfrc(k+1) - cldfrc(k)) &
2759  & / (cldfrc(k-1) - cldfrc(k))
2760  else
2761  facclr1u(k+1) = rat2
2762  endif
2763  endif
2764 
2765  if (facclr1u(k+1)>f_zero .or. facclr2u(k+1)>f_zero) then
2766  rat1 = f_one
2767  rat2 = f_zero
2768  else
2769  rat1 = f_zero
2770  rat2 = f_zero
2771  endif
2772  else
2773  if (lstcldu(k)) then
2774  faccld2u(k+1) = (cldfrc(k) - cldfrc(k+1)) / cldfrc(k)
2775  facclr2u(k) = f_zero
2776  faccld2u(k) = f_zero
2777  else
2778  fmin = min(cldfrc(k), cldfrc(k-1))
2779  if (cldfrc(k+1) <= fmin) then
2780  faccld1u(k+1) = rat1
2781  faccld2u(k+1) = (fmin - cldfrc(k+1)) / fmin
2782  else
2783  faccld1u(k+1) = (cldfrc(k) - cldfrc(k+1)) &
2784  & / (cldfrc(k) - fmin)
2785  endif
2786  endif
2787 
2788  if (faccld1u(k+1)>f_zero .or. faccld2u(k+1)>f_zero) then
2789  rat1 = f_zero
2790  rat2 = f_one
2791  else
2792  rat1 = f_zero
2793  rat2 = f_zero
2794  endif
2795  endif
2796 
2797  faccmb1u(k+1) = facclr1u(k+1) * faccld2u(k) * cldfrc(k-1)
2798  faccmb2u(k+1) = faccld1u(k+1) * facclr2u(k) &
2799  & * (f_one - cldfrc(k-1))
2800  endif
2801 
2802  enddo
2803 
2804  do k = 0, nlay
2805  faccld1d(k) = f_zero
2806  faccld2d(k) = f_zero
2807  facclr1d(k) = f_zero
2808  facclr2d(k) = f_zero
2809  faccmb1d(k) = f_zero
2810  faccmb2d(k) = f_zero
2811  enddo
2812 
2813  lstcldd(nlay) = cldfrc(nlay) > eps
2814  rat1 = f_zero
2815  rat2 = f_zero
2816 
2817  do k = nlay, 2, -1
2818 
2819  lstcldd(k-1) = cldfrc(k-1) > eps .and. cldfrc(k)<=eps
2820 
2821  if (cldfrc(k) > eps) then
2822 
2823  if (cldfrc(k-1) >= cldfrc(k)) then
2824  if (lstcldd(k)) then
2825  if (cldfrc(k) < f_one) then
2826  facclr2d(k-1) = (cldfrc(k-1) - cldfrc(k)) &
2827  & / (f_one - cldfrc(k))
2828  endif
2829 
2830  facclr2d(k) = f_zero
2831  faccld2d(k) = f_zero
2832  else
2833  fmax = max(cldfrc(k), cldfrc(k+1))
2834 
2835  if (cldfrc(k-1) > fmax) then
2836  facclr1d(k-1) = rat2
2837  facclr2d(k-1) = (cldfrc(k-1) - fmax) / (f_one - fmax)
2838  elseif (cldfrc(k-1) < fmax) then
2839  facclr1d(k-1) = (cldfrc(k-1) - cldfrc(k)) &
2840  & / (cldfrc(k+1) - cldfrc(k))
2841  else
2842  facclr1d(k-1) = rat2
2843  endif
2844  endif
2845 
2846  if (facclr1d(k-1)>f_zero .or. facclr2d(k-1)>f_zero) then
2847  rat1 = f_one
2848  rat2 = f_zero
2849  else
2850  rat1 = f_zero
2851  rat2 = f_zero
2852  endif
2853  else
2854  if (lstcldd(k)) then
2855  faccld2d(k-1) = (cldfrc(k) - cldfrc(k-1)) / cldfrc(k)
2856  facclr2d(k) = f_zero
2857  faccld2d(k) = f_zero
2858  else
2859  fmin = min(cldfrc(k), cldfrc(k+1))
2860 
2861  if (cldfrc(k-1) <= fmin) then
2862  faccld1d(k-1) = rat1
2863  faccld2d(k-1) = (fmin - cldfrc(k-1)) / fmin
2864  else
2865  faccld1d(k-1) = (cldfrc(k) - cldfrc(k-1)) &
2866  & / (cldfrc(k) - fmin)
2867  endif
2868  endif
2869 
2870  if (faccld1d(k-1)>f_zero .or. faccld2d(k-1)>f_zero) then
2871  rat1 = f_zero
2872  rat2 = f_one
2873  else
2874  rat1 = f_zero
2875  rat2 = f_zero
2876  endif
2877  endif
2878 
2879  faccmb1d(k-1) = facclr1d(k-1) * faccld2d(k) * cldfrc(k+1)
2880  faccmb2d(k-1) = faccld1d(k-1) * facclr2d(k) &
2881  & * (f_one - cldfrc(k+1))
2882  endif
2883 
2884  enddo
2885 
2886 ! --- ... initialize for radiative transfer.
2887 
2888  do ib = 1, nbands
2889  do k = 0, nlay
2890  toturad(k,ib) = f_zero
2891  totdrad(k,ib) = f_zero
2892  clrurad(k,ib) = f_zero
2893  clrdrad(k,ib) = f_zero
2894  enddo
2895  enddo
2896 
2897  do k = 0, nlay
2898  totuflux(k) = f_zero
2899  totdflux(k) = f_zero
2900  totuclfl(k) = f_zero
2901  totdclfl(k) = f_zero
2902  enddo
2903 
2904 ! --- ... loop over all g-points
2905 
2906  do ig = 1, ngptlw
2907  ib = ngb(ig)
2908 
2909  radtotd = f_zero
2910  radclrd = f_zero
2911 
2912 ! --- ... downward radiative transfer loop.
2913 
2914  do k = nlay, 1, -1
2915 
2916 ! --- ... clear sky, gases contribution
2917 
2918  odepth = max( f_zero, secdif(ib)*tautot(ig,k) )
2919  if (odepth <= 0.06) then
2920  atrgas = odepth - 0.5*odepth*odepth
2921  trng = f_one - atrgas
2922  gasfac = rec_6 * odepth
2923  else
2924  tblind = odepth / (bpade + odepth)
2925  itgas = tblint*tblind + 0.5
2926  trng = exp_tbl(itgas)
2927  atrgas = f_one - trng
2928  gasfac = tfn_tbl(itgas)
2929  odepth = tau_tbl(itgas)
2930  endif
2931 
2932  plfrac = fracs(ig,k)
2933  blay = pklay(ib,k)
2934 
2935  dplnku = pklev(ib,k ) - blay
2936  dplnkd = pklev(ib,k-1) - blay
2937  bbdgas = plfrac * (blay + dplnkd*gasfac)
2938  bbugas = plfrac * (blay + dplnku*gasfac)
2939  gassrcd = bbdgas * atrgas
2940  gassrcu(k)= bbugas * atrgas
2941  trngas(k) = trng
2942 
2943 ! --- ... total sky, gases+clouds contribution
2944 
2945  clfr = cldfrc(k)
2946  if (lstcldd(k)) then
2947  totradd = clfr * radtotd
2948  clrradd = radtotd - totradd
2949  rad = f_zero
2950  endif
2951 
2952  if (clfr >= eps) then
2953 ! --- ... cloudy layer
2954 
2955  odcld = secdif(ib) * taucld(ib,k)
2956  odtot = odepth + odcld
2957  if (odtot < 0.06) then
2958  totfac = rec_6 * odtot
2959  atrtot = odtot - 0.5*odtot*odtot
2960  trnt = f_one - atrtot
2961  else
2962  tblind = odtot / (bpade + odtot)
2963  ittot = tblint*tblind + 0.5
2964  totfac = tfn_tbl(ittot)
2965  trnt = exp_tbl(ittot)
2966  atrtot = f_one - trnt
2967  endif
2968 
2969  bbdtot = plfrac * (blay + dplnkd*totfac)
2970  bbutot = plfrac * (blay + dplnku*totfac)
2971  totsrcd = bbdtot * atrtot
2972  totsrcu(k)= bbutot * atrtot
2973  trntot(k) = trnt
2974 
2975  totradd = totradd*trnt + clfr*totsrcd
2976  clrradd = clrradd*trng + (f_one - clfr)*gassrcd
2977 
2978 ! --- ... total sky radiance
2979  radtotd = totradd + clrradd
2980  totdrad(k-1,ib) = totdrad(k-1,ib) + radtotd
2981 
2982 ! --- ... clear sky radiance
2983  radclrd = radclrd*trng + gassrcd
2984  clrdrad(k-1,ib) = clrdrad(k-1,ib) + radclrd
2985 
2986  radmod = rad*(facclr1d(k-1)*trng + faccld1d(k-1)*trnt) &
2987  & - faccmb1d(k-1)*gassrcd + faccmb2d(k-1)*totsrcd
2988 
2989  rad = -radmod + facclr2d(k-1)*(clrradd + radmod) &
2990  & - faccld2d(k-1)*(totradd - radmod)
2991  totradd = totradd + rad
2992  clrradd = clrradd - rad
2993 
2994  else
2995 ! --- ... clear layer
2996 
2997 ! --- ... total sky radiance
2998  radtotd = radtotd*trng + gassrcd
2999  totdrad(k-1,ib) = totdrad(k-1,ib) + radtotd
3000 
3001 ! --- ... clear sky radiance
3002  radclrd = radclrd*trng + gassrcd
3003  clrdrad(k-1,ib) = clrdrad(k-1,ib) + radclrd
3004 
3005  endif ! end if_clfr_block
3006 
3007  enddo ! end do_k_loop
3008 
3009 ! --- ... spectral emissivity & reflectance
3010 ! include the contribution of spectrally varying longwave emissivity
3011 ! and reflection from the surface to the upward radiative transfer.
3012 ! note: spectral and Lambertian reflection are identical for the
3013 ! diffusivity angle flux integration used here.
3014 
3015  reflct = f_one - semiss(ib)
3016  rad0 = semiss(ib) * fracs(ig,1) * pklay(ib,0)
3017 
3018 ! --- ... total sky radiance
3019  radtotu = rad0 + reflct*radtotd
3020  toturad(0,ib) = toturad(0,ib) + radtotu
3021 
3022 ! --- ... clear sky radiance
3023  radclru = rad0 + reflct*radclrd
3024  clrurad(0,ib) = clrurad(0,ib) + radclru
3025 
3026 ! --- ... upward radiative transfer loop.
3027 
3028  do k = 1, nlay
3029 
3030  clfr = cldfrc(k)
3031  trng = trngas(k)
3032  gasu = gassrcu(k)
3033 
3034  if (lstcldu(k)) then
3035  totradu = clfr * radtotu
3036  clrradu = radtotu - totradu
3037  rad = f_zero
3038  endif
3039 
3040  if (clfr >= eps) then
3041 ! --- ... cloudy layer
3042 
3043  trnt = trntot(k)
3044  totu = totsrcu(k)
3045  totradu = totradu*trnt + clfr*totu
3046  clrradu = clrradu*trng + (f_one - clfr)*gasu
3047 
3048 ! --- ... total sky radiance
3049  radtotu = totradu + clrradu
3050  toturad(k,ib) = toturad(k,ib) + radtotu
3051 
3052 ! --- ... clear sky radiance
3053  radclru = radclru*trng + gasu
3054  clrurad(k,ib) = clrurad(k,ib) + radclru
3055 
3056  radmod = rad*(facclr1u(k+1)*trng + faccld1u(k+1)*trnt) &
3057  & - faccmb1u(k+1)*gasu + faccmb2u(k+1)*totu
3058  rad = -radmod + facclr2u(k+1)*(clrradu + radmod) &
3059  & - faccld2u(k+1)*(totradu - radmod)
3060  totradu = totradu + rad
3061  clrradu = clrradu - rad
3062 
3063  else
3064 ! --- ... clear layer
3065 
3066 ! --- ... total sky radiance
3067  radtotu = radtotu*trng + gasu
3068  toturad(k,ib) = toturad(k,ib) + radtotu
3069 
3070 ! --- ... clear sky radiance
3071  radclru = radclru*trng + gasu
3072  clrurad(k,ib) = clrurad(k,ib) + radclru
3073 
3074  endif ! end if_clfr_block
3075 
3076  enddo ! end do_k_loop
3077 
3078  enddo ! end do_ig_loop
3079 
3080 ! --- ... process longwave output from band for total and clear streams.
3081 ! calculate upward, downward, and net flux.
3082 
3083  flxfac = wtdiff * fluxfac
3084 
3085  do k = 0, nlay
3086  do ib = 1, nbands
3087  totuflux(k) = totuflux(k) + toturad(k,ib)
3088  totdflux(k) = totdflux(k) + totdrad(k,ib)
3089  totuclfl(k) = totuclfl(k) + clrurad(k,ib)
3090  totdclfl(k) = totdclfl(k) + clrdrad(k,ib)
3091  enddo
3092 
3093  totuflux(k) = totuflux(k) * flxfac
3094  totdflux(k) = totdflux(k) * flxfac
3095  totuclfl(k) = totuclfl(k) * flxfac
3096  totdclfl(k) = totdclfl(k) * flxfac
3097  enddo
3098 
3099 ! --- ... calculate net fluxes and heating rates
3100  fnet(0) = totuflux(0) - totdflux(0)
3101 
3102  do k = 1, nlay
3103  rfdelp(k) = heatfac / delp(k)
3104  fnet(k) = totuflux(k) - totdflux(k)
3105  htr(k) = (fnet(k-1) - fnet(k)) * rfdelp(k)
3106  enddo
3107 
3108 !! --- ... optional clear sky heating rates
3109  if ( lhlw0 ) then
3110  fnetc(0) = totuclfl(0) - totdclfl(0)
3111 
3112  do k = 1, nlay
3113  fnetc(k) = totuclfl(k) - totdclfl(k)
3114  htrcl(k) = (fnetc(k-1) - fnetc(k)) * rfdelp(k)
3115  enddo
3116  endif
3117 
3118 !! --- ... optional spectral band heating rates
3119  if ( lhlwb ) then
3120  do ib = 1, nbands
3121  fnet(0) = (toturad(0,ib) - totdrad(0,ib)) * flxfac
3122 
3123  do k = 1, nlay
3124  fnet(k) = (toturad(k,ib) - totdrad(k,ib)) * flxfac
3125  htrb(k,ib) = (fnet(k-1) - fnet(k)) * rfdelp(k)
3126  enddo
3127  enddo
3128  endif
3129 
3130 ! .................................
3131  end subroutine rtrnmr
3132 ! ---------------------------------
3133 
3134 
3135 ! ---------------------------------
3136  subroutine rtrnmc &
3137 ! .................................
3138 ! --- inputs:
3139  & ( semiss,delp,cldfmc,taucld,tautot,pklay,pklev, &
3140  & fracs,secdif, nlay,nlp1, &
3141 ! --- outputs:
3142  & totuflux,totdflux,htr, totuclfl,totdclfl,htrcl, htrb &
3143  & )
3144 
3145 ! =================== program usage description =================== !
3146 ! !
3147 ! purpose: compute the upward/downward radiative fluxes, and heating !
3148 ! rates for both clear or cloudy atmosphere. clouds are treated with !
3149 ! the mcica stochastic approach. !
3150 ! !
3151 ! subprograms called: none !
3152 ! !
3153 ! ==================== defination of variables ==================== !
3154 ! !
3155 ! inputs: -size- !
3156 ! semiss - real, lw surface emissivity nbands!
3157 ! delp - real, layer pressure thickness (mb) nlay !
3158 ! cldfmc - real, layer cloud fraction (sub-column) ngptlw*nlay!
3159 ! taucld - real, layer cloud opt depth nbands*nlay!
3160 ! tautot - real, total optical depth (gas+aerosols) ngptlw*nlay!
3161 ! pklay - real, integrated planck func at lay temp nbands*0:nlay!
3162 ! pklev - real, integrated planck func at lev temp nbands*0:nlay!
3163 ! fracs - real, planck fractions ngptlw*nlay!
3164 ! secdif - real, secant of diffusivity angle nbands!
3165 ! nlay - integer, number of vertical layers 1 !
3166 ! nlp1 - integer, number of vertical levels (interfaces) 1 !
3167 ! !
3168 ! outputs: !
3169 ! totuflux- real, total sky upward flux (w/m2) 0:nlay !
3170 ! totdflux- real, total sky downward flux (w/m2) 0:nlay !
3171 ! htr - real, total sky heating rate (k/sec or k/day) nlay !
3172 ! totuclfl- real, clear sky upward flux (w/m2) 0:nlay !
3173 ! totdclfl- real, clear sky downward flux (w/m2) 0:nlay !
3174 ! htrcl - real, clear sky heating rate (k/sec or k/day) nlay !
3175 ! htrb - real, spectral band lw heating rate (k/day) nlay*nbands!
3176 ! !
3177 ! module veriables: !
3178 ! ngb - integer, band index for each g-value ngptlw!
3179 ! fluxfac - real, conversion factor for fluxes (pi*2.e4) 1 !
3180 ! heatfac - real, conversion factor for heating rates (g/cp*1e-2) 1 !
3181 ! tblint - real, conversion factor for look-up tbl (float(ntbl) 1 !
3182 ! bpade - real, pade approx constant (1/0.278) 1 !
3183 ! wtdiff - real, weight for radiance to flux conversion 1 !
3184 ! ntbl - integer, dimension of look-up tables 1 !
3185 ! tau_tbl - real, clr-sky opt dep lookup table 0:ntbl !
3186 ! exp_tbl - real, transmittance lookup table 0:ntbl !
3187 ! tfn_tbl - real, tau transition function 0:ntbl !
3188 ! !
3189 ! local variables: !
3190 ! itgas - integer, index for gases contribution look-up table 1 !
3191 ! ittot - integer, index for gases plus clouds look-up table 1 !
3192 ! reflct - real, surface reflectance 1 !
3193 ! atrgas - real, gaseous absorptivity 1 !
3194 ! atrtot - real, gaseous and cloud absorptivity 1 !
3195 ! odcld - real, cloud optical depth 1 !
3196 ! efclrfr- real, effective clear sky fraction (1-efcldfr) nlay!
3197 ! odepth - real, optical depth of gaseous only 1 !
3198 ! odtot - real, optical depth of gas and cloud 1 !
3199 ! gasfac - real, gas-only pade factor, used for planck function 1 !
3200 ! totfac - real, gas and cloud pade factor, used for planck fn 1 !
3201 ! bbdgas - real, gas-only planck function for downward rt 1 !
3202 ! bbugas - real, gas-only planck function for upward rt 1 !
3203 ! bbdtot - real, gas and cloud planck function for downward rt 1 !
3204 ! bbutot - real, gas and cloud planck function for upward rt 1 !
3205 ! gassrcu- real, upwd source radiance due to gas nlay!
3206 ! totsrcu- real, upwd source radiance due to gas+cld nlay!
3207 ! gassrcd- real, dnwd source radiance due to gas 1 !
3208 ! totsrcd- real, dnwd source radiance due to gas+cld 1 !
3209 ! radtotu- real, spectrally summed total sky upwd radiance 1 !
3210 ! radclru- real, spectrally summed clear sky upwd radiance 1 !
3211 ! radtotd- real, spectrally summed total sky dnwd radiance 1 !
3212 ! radclrd- real, spectrally summed clear sky dnwd radiance 1 !
3213 ! toturad- real, total sky upward radiance by layer 0:nlay*nbands!
3214 ! clrurad- real, clear sky upward radiance by layer 0:nlay*nbands!
3215 ! totdrad- real, total sky downward radiance by layer 0:nlay*nbands!
3216 ! clrdrad- real, clear sky downward radiance by layer 0:nlay*nbands!
3217 ! fnet - real, net longwave flux (w/m2) 0:nlay !
3218 ! fnetc - real, clear sky net longwave flux (w/m2) 0:nlay !
3219 ! !
3220 ! !
3221 ! ******************************************************************* !
3222 ! original code description !
3223 ! !
3224 ! original version: e. j. mlawer, et al. rrtm_v3.0 !
3225 ! revision for gcms: michael j. iacono; october, 2002 !
3226 ! revision for f90: michael j. iacono; june, 2006 !
3227 ! !
3228 ! this program calculates the upward fluxes, downward fluxes, and !
3229 ! heating rates for an arbitrary clear or cloudy atmosphere. the input !
3230 ! to this program is the atmospheric profile, all Planck function !
3231 ! information, and the cloud fraction by layer. a variable diffusivity!
3232 ! angle (secdif) is used for the angle integration. bands 2-3 and 5-9 !
3233 ! use a value for secdif that varies from 1.50 to 1.80 as a function !
3234 ! of the column water vapor, and other bands use a value of 1.66. the !
3235 ! gaussian weight appropriate to this angle (wtdiff=0.5) is applied !
3236 ! here. note that use of the emissivity angle for the flux integration!
3237 ! can cause errors of 1 to 4 W/m2 within cloudy layers. !
3238 ! clouds are treated with the mcica stochastic approach and !
3239 ! maximum-random cloud overlap. !
3240 ! !
3241 ! ******************************************************************* !
3242 ! ====================== end of description block ================= !
3243 
3244 ! --- inputs:
3245  integer, intent(in) :: nlay, nlp1
3246 
3247  real (kind=kind_phys), dimension(nbands), intent(in) :: semiss, &
3248  & secdif
3249  real (kind=kind_phys), dimension(nlay), intent(in) :: delp
3250 
3251  real (kind=kind_phys), dimension(nbands,nlay),intent(in):: taucld
3252  real (kind=kind_phys), dimension(ngptlw,nlay),intent(in):: fracs, &
3253  & tautot, cldfmc
3254 
3255  real (kind=kind_phys), dimension(nbands,0:nlay), intent(in) :: &
3256  & pklev, pklay
3257 
3258 ! --- outputs:
3259  real (kind=kind_phys), dimension(nlay), intent(out) :: htr, htrcl
3260 
3261  real (kind=kind_phys), dimension(nlay,nbands),intent(out) :: htrb
3262 
3263  real (kind=kind_phys), dimension(0:nlay), intent(out) :: &
3264  & totuflux, totdflux, totuclfl, totdclfl
3265 
3266 ! --- locals:
3267  real (kind=kind_phys), parameter :: rec_6 = 0.166667
3268 
3269  real (kind=kind_phys), dimension(0:nlay,nbands) :: clrurad, &
3270  & clrdrad, toturad, totdrad
3271 
3272  real (kind=kind_phys), dimension(nlay) :: gassrcu, totsrcu, &
3273  & trngas, efclrfr, rfdelp
3274  real (kind=kind_phys), dimension(0:nlay) :: fnet, fnetc
3275 
3276  real (kind=kind_phys) :: totsrcd, gassrcd, tblind, odepth, odtot, &
3277  & odcld, atrtot, atrgas, reflct, totfac, gasfac, flxfac, &
3278  & plfrac, blay, bbdgas, bbdtot, bbugas, bbutot, dplnku, &
3279  & dplnkd, radtotu, radclru, radtotd, radclrd, rad0, &
3280  & clfm, trng, gasu
3281 
3282  integer :: ittot, itgas, ib, ig, k
3283 !
3284 !===> ... begin here
3285 !
3286  do ib = 1, nbands
3287  do k = 0, nlay
3288  toturad(k,ib) = f_zero
3289  totdrad(k,ib) = f_zero
3290  clrurad(k,ib) = f_zero
3291  clrdrad(k,ib) = f_zero
3292  enddo
3293  enddo
3294 
3295  do k = 0, nlay
3296  totuflux(k) = f_zero
3297  totdflux(k) = f_zero
3298  totuclfl(k) = f_zero
3299  totdclfl(k) = f_zero
3300  enddo
3301 
3302 ! --- ... loop over all g-points
3303 
3304  do ig = 1, ngptlw
3305  ib = ngb(ig)
3306 
3307  radtotd = f_zero
3308  radclrd = f_zero
3309 
3310 ! --- ... downward radiative transfer loop.
3311 
3312  do k = nlay, 1, -1
3313 
3314 ! --- ... clear sky, gases contribution
3315 
3316  odepth = max( f_zero, secdif(ib)*tautot(ig,k) )
3317  if (odepth <= 0.06) then
3318  atrgas = odepth - 0.5*odepth*odepth
3319  trng = f_one - atrgas
3320  gasfac = rec_6 * odepth
3321  else
3322  tblind = odepth / (bpade + odepth)
3323  itgas = tblint*tblind + 0.5
3324  trng = exp_tbl(itgas)
3325  atrgas = f_one - trng
3326  gasfac = tfn_tbl(itgas)
3327  odepth = tau_tbl(itgas)
3328  endif
3329 
3330  plfrac = fracs(ig,k)
3331  blay = pklay(ib,k)
3332 
3333  dplnku = pklev(ib,k ) - blay
3334  dplnkd = pklev(ib,k-1) - blay
3335  bbdgas = plfrac * (blay + dplnkd*gasfac)
3336  bbugas = plfrac * (blay + dplnku*gasfac)
3337  gassrcd= bbdgas * atrgas
3338  gassrcu(k)= bbugas * atrgas
3339  trngas(k) = trng
3340 
3341 ! --- ... total sky, gases+clouds contribution
3342 
3343  clfm = cldfmc(ig,k)
3344  if (clfm >= eps) then
3345 ! --- ... cloudy layer
3346 
3347  odcld = secdif(ib) * taucld(ib,k)
3348  efclrfr(k) = f_one - (f_one - exp(-odcld))*clfm
3349  odtot = odepth + odcld
3350  if (odtot < 0.06) then
3351  totfac = rec_6 * odtot
3352  atrtot = odtot - 0.5*odtot*odtot
3353  else
3354  tblind = odtot / (bpade + odtot)
3355  ittot = tblint*tblind + 0.5
3356  totfac = tfn_tbl(ittot)
3357  atrtot = f_one - exp_tbl(ittot)
3358  endif
3359 
3360  bbdtot = plfrac * (blay + dplnkd*totfac)
3361  bbutot = plfrac * (blay + dplnku*totfac)
3362  totsrcd= bbdtot * atrtot
3363  totsrcu(k)= bbutot * atrtot
3364 
3365 ! --- ... total sky radiance
3366  radtotd = radtotd*trng*efclrfr(k) + gassrcd &
3367  & + clfm*(totsrcd - gassrcd)
3368  totdrad(k-1,ib) = totdrad(k-1,ib) + radtotd
3369 
3370 ! --- ... clear sky radiance
3371  radclrd = radclrd*trng + gassrcd
3372  clrdrad(k-1,ib) = clrdrad(k-1,ib) + radclrd
3373 
3374  else
3375 ! --- ... clear layer
3376 
3377 ! --- ... total sky radiance
3378  radtotd = radtotd*trng + gassrcd
3379  totdrad(k-1,ib) = totdrad(k-1,ib) + radtotd
3380 
3381 ! --- ... clear sky radiance
3382  radclrd = radclrd*trng + gassrcd
3383  clrdrad(k-1,ib) = clrdrad(k-1,ib) + radclrd
3384 
3385  endif ! end if_clfm_block
3386 
3387  enddo ! end do_k_loop
3388 
3389 ! --- ... spectral emissivity & reflectance
3390 ! include the contribution of spectrally varying longwave emissivity
3391 ! and reflection from the surface to the upward radiative transfer.
3392 ! note: spectral and Lambertian reflection are identical for the
3393 ! diffusivity angle flux integration used here.
3394 
3395  reflct = f_one - semiss(ib)
3396  rad0 = semiss(ib) * fracs(ig,1) * pklay(ib,0)
3397 
3398 ! --- ... total sky radiance
3399  radtotu = rad0 + reflct*radtotd
3400  toturad(0,ib) = toturad(0,ib) + radtotu
3401 
3402 ! --- ... clear sky radiance
3403  radclru = rad0 + reflct*radclrd
3404  clrurad(0,ib) = clrurad(0,ib) + radclru
3405 
3406 ! --- ... upward radiative transfer
3407 ! toturad holds summed radiance for total sky stream
3408 ! clrurad holds summed radiance for clear sky stream
3409 
3410  do k = 1, nlay
3411  clfm = cldfmc(ig,k)
3412  trng = trngas(k)
3413  gasu = gassrcu(k)
3414 
3415  if (clfm > eps) then
3416 ! --- ... cloudy layer
3417 
3418 ! --- ... total sky radiance
3419  radtotu = radtotu*trng*efclrfr(k) + gasu &
3420  & + clfm*(totsrcu(k) - gasu)
3421  toturad(k,ib) = toturad(k,ib) + radtotu
3422 
3423 ! --- ... clear sky radiance
3424  radclru = radclru*trng + gasu
3425  clrurad(k,ib) = clrurad(k,ib) + radclru
3426 
3427  else
3428 ! --- ... clear layer
3429 
3430 ! --- ... total sky radiance
3431  radtotu = radtotu*trng + gasu
3432  toturad(k,ib) = toturad(k,ib) + radtotu
3433 
3434 ! --- ... clear sky radiance
3435  radclru = radclru*trng + gasu
3436  clrurad(k,ib) = clrurad(k,ib) + radclru
3437 
3438  endif ! end if_clfm_block
3439 
3440  enddo ! end do_k_loop
3441 
3442  enddo ! end do_ig_loop
3443 
3444 ! --- ... process longwave output from band for total and clear streams.
3445 ! calculate upward, downward, and net flux.
3446 
3447  flxfac = wtdiff * fluxfac
3448 
3449  do k = 0, nlay
3450  do ib = 1, nbands
3451  totuflux(k) = totuflux(k) + toturad(k,ib)
3452  totdflux(k) = totdflux(k) + totdrad(k,ib)
3453  totuclfl(k) = totuclfl(k) + clrurad(k,ib)
3454  totdclfl(k) = totdclfl(k) + clrdrad(k,ib)
3455  enddo
3456 
3457  totuflux(k) = totuflux(k) * flxfac
3458  totdflux(k) = totdflux(k) * flxfac
3459  totuclfl(k) = totuclfl(k) * flxfac
3460  totdclfl(k) = totdclfl(k) * flxfac
3461  enddo
3462 
3463 ! --- ... calculate net fluxes and heating rates
3464  fnet(0) = totuflux(0) - totdflux(0)
3465 
3466  do k = 1, nlay
3467  rfdelp(k) = heatfac / delp(k)
3468  fnet(k) = totuflux(k) - totdflux(k)
3469  htr(k) = (fnet(k-1) - fnet(k)) * rfdelp(k)
3470  enddo
3471 
3472 !! --- ... optional clear sky heating rates
3473  if ( lhlw0 ) then
3474  fnetc(0) = totuclfl(0) - totdclfl(0)
3475 
3476  do k = 1, nlay
3477  fnetc(k) = totuclfl(k) - totdclfl(k)
3478  htrcl(k) = (fnetc(k-1) - fnetc(k)) * rfdelp(k)
3479  enddo
3480  endif
3481 
3482 !! --- ... optional spectral band heating rates
3483  if ( lhlwb ) then
3484  do ib = 1, nbands
3485  fnet(0) = (toturad(0,ib) - totdrad(0,ib)) * flxfac
3486 
3487  do k = 1, nlay
3488  fnet(k) = (toturad(k,ib) - totdrad(k,ib)) * flxfac
3489  htrb(k,ib) = (fnet(k-1) - fnet(k)) * rfdelp(k)
3490  enddo
3491  enddo
3492  endif
3493 
3494 ! ..................................
3495  end subroutine rtrnmc
3496 ! ----------------------------------
3497 
3498 
3499 ! ----------------------------------
3500  subroutine taumol &
3501 ! ..................................
3502 ! --- inputs:
3503  & ( laytrop,pavel,coldry,colamt,colbrd,wx,tauaer, &
3504  & rfrate,fac00,fac01,fac10,fac11,jp,jt,jt1, &
3505  & selffac,selffrac,indself,forfac,forfrac,indfor, &
3506  & minorfrac,scaleminor,scaleminorn2,indminor, &
3507  & nlay, &
3508 ! --- outputs:
3509  & fracs, tautot &
3510  & )
3511 
3512 ! ************ original subprogram description *************** !
3513 ! !
3514 ! optical depths developed for the !
3515 ! !
3516 ! rapid radiative transfer model (rrtm) !
3517 ! !
3518 ! atmospheric and environmental research, inc. !
3519 ! 131 hartwell avenue !
3520 ! lexington, ma 02421 !
3521 ! !
3522 ! eli j. mlawer !
3523 ! jennifer delamere !
3524 ! steven j. taubman !
3525 ! shepard a. clough !
3526 ! !
3527 ! email: mlawer@aer.com !
3528 ! email: jdelamer@aer.com !
3529 ! !
3530 ! the authors wish to acknowledge the contributions of the !
3531 ! following people: karen cady-pereira, patrick d. brown, !
3532 ! michael j. iacono, ronald e. farren, luke chen, !
3533 ! robert bergstrom. !
3534 ! !
3535 ! revision for g-point reduction: michael j. iacono; aer, inc. !
3536 ! !
3537 ! taumol !
3538 ! !
3539 ! this file contains the subroutines taugbn (where n goes from !
3540 ! 1 to 16). taugbn calculates the optical depths and planck !
3541 ! fractions per g-value and layer for band n. !
3542 ! !
3543 ! ******************************************************************* !
3544 ! ================== program usage description ================== !
3545 ! !
3546 ! call taumol !
3547 ! inputs: !
3548 ! ( laytrop,pavel,coldry,colamt,colbrd,wx,tauaer, !
3549 ! rfrate,fac00,fac01,fac10,fac11,jp,jt,jt1, !
3550 ! selffac,selffrac,indself,forfac,forfrac,indfor, !
3551 ! minorfrac,scaleminor,scaleminorn2,indminor, !
3552 ! nlay, !
3553 ! outputs: !
3554 ! fracs, tautot ) !
3555 ! !
3556 ! subprograms called: taugb## (## = 01 -16) !
3557 ! !
3558 ! !
3559 ! ==================== defination of variables ==================== !
3560 ! !
3561 ! inputs: size !
3562 ! laytrop - integer, tropopause layer index (unitless) 1 !
3563 ! layer at which switch is made for key species !
3564 ! pavel - real, layer pressures (mb) nlay !
3565 ! coldry - real, column amount for dry air (mol/cm2) nlay !
3566 ! colamt - real, column amounts of h2o, co2, o3, n2o, ch4, !
3567 ! o2, co (mol/cm**2) nlay*maxgas!
3568 ! colbrd - real, column amount of broadening gases nlay !
3569 ! wx - real, cross-section amounts(mol/cm2) nlay*maxxsec!
3570 ! tauaer - real, aerosol optical depth nbands*nlay !
3571 ! rfrate - real, reference ratios of binary species parameter !
3572 ! (:,m,:)m=1-h2o/co2,2-h2o/o3,3-h2o/n2o,4-h2o/ch4,5-n2o/co2,6-o3/co2!
3573 ! (:,:,n)n=1,2: the rates of ref press at the 2 sides of the layer !
3574 ! nlay*nrates*2!
3575 ! facij - real, factors multiply the reference ks, i,j of 0/1 !
3576 ! for lower/higher of the 2 appropriate temperatures !
3577 ! and altitudes nlay !
3578 ! jp - real, index of lower reference pressure nlay !
3579 ! jt, jt1 - real, indices of lower reference temperatures nlay !
3580 ! for pressure levels jp and jp+1, respectively !
3581 ! selffac - real, scale factor for water vapor self-continuum !
3582 ! equals (water vapor density)/(atmospheric density !
3583 ! at 296k and 1013 mb) nlay !
3584 ! selffrac - real, factor for temperature interpolation of !
3585 ! reference water vapor self-continuum data nlay !
3586 ! indself - integer, index of lower reference temperature for !
3587 ! the self-continuum interpolation nlay !
3588 ! forfac - real, scale factor for w. v. foreign-continuum nlay !
3589 ! forfrac - real, factor for temperature interpolation of !
3590 ! reference w.v. foreign-continuum data nlay !
3591 ! indfor - integer, index of lower reference temperature for !
3592 ! the foreign-continuum interpolation nlay !
3593 ! minorfrac - real, factor for minor gases nlay !
3594 ! scaleminor,scaleminorn2 !
3595 ! - real, scale factors for minor gases nlay !
3596 ! indminor - integer, index of lower reference temperature for !
3597 ! minor gases nlay !
3598 ! nlay - integer, total number of layers 1 !
3599 ! !
3600 ! outputs: !
3601 ! fracs - real, planck fractions ngptlw,nlay!
3602 ! tautot - real, total optical depth (gas+aerosols) ngptlw,nlay!
3603 ! !
3604 ! internal variables: !
3605 ! ng## - integer, number of g-values in band ## (##=01-16) 1 !
3606 ! nspa - integer, for lower atmosphere, the number of ref !
3607 ! atmos, each has different relative amounts of the !
3608 ! key species for the band nbands!
3609 ! nspb - integer, same but for upper atmosphere nbands!
3610 ! absa - real, k-values for lower ref atmospheres (no w.v. !
3611 ! self-continuum) (cm**2/molecule) nspa(##)*5*13*ng##!
3612 ! absb - real, k-values for high ref atmospheres (all sources) !
3613 ! (cm**2/molecule) nspb(##)*5*13:59*ng##!
3614 ! ka_m'mgas'- real, k-values for low ref atmospheres minor species !
3615 ! (cm**2/molecule) mmn##*ng##!
3616 ! kb_m'mgas'- real, k-values for high ref atmospheres minor species !
3617 ! (cm**2/molecule) mmn##*ng##!
3618 ! selfref - real, k-values for w.v. self-continuum for ref atmos !
3619 ! used below laytrop (cm**2/mol) 10*ng##!
3620 ! forref - real, k-values for w.v. foreign-continuum for ref atmos
3621 ! used below/above laytrop (cm**2/mol) 4*ng##!
3622 ! !
3623 ! ****************************************************************** !
3624 
3625 ! --- inputs:
3626  integer, intent(in) :: nlay, laytrop
3627 
3628  integer, dimension(nlay), intent(in) :: jp, jt, jt1, indself, &
3629  & indfor, indminor
3630 
3631  real (kind=kind_phys), dimension(nlay), intent(in) :: pavel, &
3632  & coldry, colbrd, fac00, fac01, fac10, fac11, selffac, &
3633  & selffrac, forfac, forfrac, minorfrac, scaleminor, &
3634  & scaleminorn2
3635 
3636  real (kind=kind_phys), dimension(nlay,maxgas), intent(in):: colamt
3637  real (kind=kind_phys), dimension(nlay,maxxsec),intent(in):: wx
3638 
3639  real (kind=kind_phys), dimension(nbands,nlay), intent(in):: tauaer
3640 
3641  real (kind=kind_phys), dimension(nlay,nrates,2), intent(in) :: &
3642  & rfrate
3643 
3644 ! --- outputs:
3645  real (kind=kind_phys), dimension(ngptlw,nlay), intent(out) :: &
3646  & fracs, tautot
3647 
3648 ! --- locals
3649  real (kind=kind_phys), dimension(ngptlw,nlay) :: taug
3650 
3651  integer :: ib, ig, k
3652 !
3653 !===> ... begin here
3654 !
3655  call taugb01
3656  call taugb02
3657  call taugb03
3658  call taugb04
3659  call taugb05
3660  call taugb06
3661  call taugb07
3662  call taugb08
3663  call taugb09
3664  call taugb10
3665  call taugb11
3666  call taugb12
3667  call taugb13
3668  call taugb14
3669  call taugb15
3670  call taugb16
3671 
3672 ! --- combine gaseous and aerosol optical depths
3673 
3674  do ig = 1, ngptlw
3675  ib = ngb(ig)
3676 
3677  do k = 1, nlay
3678  tautot(ig,k) = taug(ig,k) + tauaer(ib,k)
3679  enddo
3680  enddo
3681 
3682 ! =================
3683  contains
3684 ! =================
3685 
3686 ! ----------------------------------
3687  subroutine taugb01
3688 ! ..................................
3689 
3690 ! ------------------------------------------------------------------ !
3691 ! written by eli j. mlawer, atmospheric & environmental research. !
3692 ! revised by michael j. iacono, atmospheric & environmental research. !
3693 ! !
3694 ! band 1: 10-350 cm-1 (low key - h2o; low minor - n2) !
3695 ! (high key - h2o; high minor - n2) !
3696 ! !
3697 ! compute the optical depth by interpolating in ln(pressure) and !
3698 ! temperature. below laytrop, the water vapor self-continuum and !
3699 ! foreign continuum is interpolated (in temperature) separately. !
3700 ! ------------------------------------------------------------------ !
3701 
3702  use module_radlw_kgb01
3703 
3704 ! --- locals:
3705  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
3706  & indm, indmp, ig
3707 
3708  real (kind=kind_phys) :: pp, corradj, scalen2, tauself, taufor, &
3709  & taun2
3710 !
3711 !===> ... begin here
3712 !
3713 ! --- minor gas mapping levels:
3714 ! lower - n2, p = 142.5490 mbar, t = 215.70 k
3715 ! upper - n2, p = 142.5490 mbar, t = 215.70 k
3716 
3717 ! --- ... lower atmosphere loop
3718 
3719  do k = 1, laytrop
3720  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(1) + 1
3721  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(1) + 1
3722  inds = indself(k)
3723  indf = indfor(k)
3724  indm = indminor(k)
3725 
3726  ind0p = ind0 + 1
3727  ind1p = ind1 + 1
3728  indsp = inds + 1
3729  indfp = indf + 1
3730  indmp = indm + 1
3731 
3732  pp = pavel(k)
3733  scalen2 = colbrd(k) * scaleminorn2(k)
3734  if (pp < 250.0) then
3735  corradj = f_one - 0.15 * (250.0-pp) / 154.4
3736  else
3737  corradj = f_one
3738  endif
3739 
3740  do ig = 1, ng01
3741  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
3742  & * (selfref(ig,indsp) - selfref(ig,inds)))
3743  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
3744  & * (forref(ig,indfp) - forref(ig,indf)))
3745  taun2 = scalen2 * (ka_mn2(ig,indm) + minorfrac(k) &
3746  & * (ka_mn2(ig,indmp) - ka_mn2(ig,indm)))
3747 
3748  taug(ig,k) = corradj * (colamt(k,1) &
3749  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
3750  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
3751  & + tauself + taufor + taun2)
3752 
3753  fracs(ig,k) = fracrefa(ig)
3754  enddo
3755  enddo
3756 
3757 ! --- ... upper atmosphere loop
3758 
3759  do k = laytrop+1, nlay
3760  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(1) + 1
3761  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(1) + 1
3762  indf = indfor(k)
3763  indm = indminor(k)
3764 
3765  ind0p = ind0 + 1
3766  ind1p = ind1 + 1
3767  indfp = indf + 1
3768  indmp = indm + 1
3769 
3770  scalen2 = colbrd(k) * scaleminorn2(k)
3771  corradj = f_one - 0.15 * (pavel(k) / 95.6)
3772 
3773  do ig = 1, ng01
3774  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
3775  & * (forref(ig,indfp) - forref(ig,indf)))
3776  taun2 = scalen2 * (kb_mn2(ig,indm) + minorfrac(k) &
3777  & * (kb_mn2(ig,indmp) - kb_mn2(ig,indm)))
3778 
3779  taug(ig,k) = corradj * (colamt(k,1) &
3780  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
3781  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
3782  & + taufor + taun2)
3783 
3784  fracs(ig,k) = fracrefb(ig)
3785  enddo
3786  enddo
3787 
3788 ! ..................................
3789  end subroutine taugb01
3790 ! ----------------------------------
3791 
3792 ! ----------------------------------
3793  subroutine taugb02
3794 ! ..................................
3795 
3796 ! ------------------------------------------------------------------ !
3797 ! band 2: 350-500 cm-1 (low key - h2o; high key - h2o) !
3798 ! ------------------------------------------------------------------ !
3799 
3800  use module_radlw_kgb02
3801 
3802 ! --- locals:
3803  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
3804  & ig
3805 
3806  real (kind=kind_phys) :: corradj, tauself, taufor
3807 !
3808 !===> ... begin here
3809 !
3810 ! --- ... lower atmosphere loop
3811 
3812  do k = 1, laytrop
3813  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(2) + 1
3814  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(2) + 1
3815  inds = indself(k)
3816  indf = indfor(k)
3817 
3818  ind0p = ind0 + 1
3819  ind1p = ind1 + 1
3820  indsp = inds + 1
3821  indfp = indf + 1
3822 
3823  corradj = f_one - 0.05 * (pavel(k) - 100.0) / 900.0
3824 
3825  do ig = 1, ng02
3826  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
3827  & * (selfref(ig,indsp) - selfref(ig,inds)))
3828  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
3829  & * (forref(ig,indfp) - forref(ig,indf)))
3830 
3831  taug(ns02+ig,k) = corradj * (colamt(k,1) &
3832  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
3833  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
3834  & + tauself + taufor)
3835 
3836  fracs(ns02+ig,k) = fracrefa(ig)
3837  enddo
3838  enddo
3839 
3840 ! --- ... upper atmosphere loop
3841 
3842  do k = laytrop+1, nlay
3843  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(2) + 1
3844  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(2) + 1
3845  indf = indfor(k)
3846 
3847  ind0p = ind0 + 1
3848  ind1p = ind1 + 1
3849  indfp = indf + 1
3850 
3851  do ig = 1, ng02
3852  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
3853  & * (forref(ig,indfp) - forref(ig,indf)))
3854 
3855  taug(ns02+ig,k) = colamt(k,1) &
3856  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
3857  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
3858  & + taufor
3859 
3860  fracs(ns02+ig,k) = fracrefb(ig)
3861  enddo
3862  enddo
3863 
3864 ! ..................................
3865  end subroutine taugb02
3866 ! ----------------------------------
3867 
3868 ! ----------------------------------
3869  subroutine taugb03
3870 ! ..................................
3871 
3872 ! ------------------------------------------------------------------ !
3873 ! band 3: 500-630 cm-1 (low key - h2o,co2; low minor - n2o) !
3874 ! (high key - h2o,co2; high minor - n2o) !
3875 ! ------------------------------------------------------------------ !
3876 
3877  use module_radlw_kgb03
3878 
3879 ! --- locals:
3880  integer :: k, ind0, ind1, inds, indsp, indf, indfp, indm, indmp, &
3881  & id000, id010, id100, id110, id200, id210, jmn2o, jmn2op, &
3882  & id001, id011, id101, id111, id201, id211, jpl, jplp, &
3883  & ig, js, js1
3884 
3885  real (kind=kind_phys) :: absn2o, ratn2o, adjfac, adjcoln2o, &
3886  & speccomb, specparm, specmult, fs, &
3887  & speccomb1, specparm1, specmult1, fs1, &
3888  & speccomb_mn2o, specparm_mn2o, specmult_mn2o, fmn2o, &
3889  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
3890  & refrat_planck_a, refrat_planck_b, refrat_m_a, refrat_m_b, &
3891  & fac000, fac100, fac200, fac010, fac110, fac210, &
3892  & fac001, fac101, fac201, fac011, fac111, fac211, &
3893  & tau_major, tau_major1, tauself, taufor, n2om1, n2om2, &
3894  & p, p4, fk0, fk1, fk2
3895 !
3896 !===> ... begin here
3897 !
3898 ! --- ... minor gas mapping levels:
3899 ! lower - n2o, p = 706.272 mbar, t = 278.94 k
3900 ! upper - n2o, p = 95.58 mbar, t = 215.7 k
3901 
3902  refrat_planck_a = chi_mls(1,9)/chi_mls(2,9) ! P = 212.725 mb
3903  refrat_planck_b = chi_mls(1,13)/chi_mls(2,13) ! P = 95.58 mb
3904  refrat_m_a = chi_mls(1,3)/chi_mls(2,3) ! P = 706.270 mb
3905  refrat_m_b = chi_mls(1,13)/chi_mls(2,13) ! P = 95.58 mb
3906 
3907 ! --- ... lower atmosphere loop
3908 
3909  do k = 1, laytrop
3910  speccomb = colamt(k,1) + rfrate(k,1,1)*colamt(k,2)
3911  specparm = colamt(k,1) / speccomb
3912  specmult = 8.0 * min(specparm, oneminus)
3913  js = 1 + int(specmult)
3914  fs = mod(specmult, f_one)
3915  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(3) + js
3916 
3917  speccomb1 = colamt(k,1) + rfrate(k,1,2)*colamt(k,2)
3918  specparm1 = colamt(k,1) / speccomb1
3919  specmult1 = 8.0 * min(specparm1, oneminus)
3920  js1 = 1 + int(specmult1)
3921  fs1 = mod(specmult1, f_one)
3922  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(3) + js1
3923 
3924  speccomb_mn2o = colamt(k,1) + refrat_m_a*colamt(k,2)
3925  specparm_mn2o = colamt(k,1) / speccomb_mn2o
3926  specmult_mn2o = 8.0 * min(specparm_mn2o, oneminus)
3927  jmn2o = 1 + int(specmult_mn2o)
3928  fmn2o = mod(specmult_mn2o, f_one)
3929 
3930  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,2)
3931  specparm_planck = colamt(k,1) / speccomb_planck
3932  specmult_planck = 8.0 * min(specparm_planck, oneminus)
3933  jpl = 1 + int(specmult_planck)
3934  fpl = mod(specmult_planck, f_one)
3935 
3936  inds = indself(k)
3937  indf = indfor(k)
3938  indm = indminor(k)
3939  indsp = inds + 1
3940  indfp = indf + 1
3941  indmp = indm + 1
3942  jmn2op= jmn2o+ 1
3943  jplp = jpl + 1
3944 
3945 ! --- ... in atmospheres where the amount of n2O is too great to be considered
3946 ! a minor species, adjust the column amount of n2O by an empirical factor
3947 ! to obtain the proper contribution.
3948 
3949  p = coldry(k) * chi_mls(4,jp(k)+1)
3950  ratn2o = colamt(k,4) / p
3951  if (ratn2o > 1.5) then
3952  adjfac = 0.5 + (ratn2o - 0.5)**0.65
3953  adjcoln2o = adjfac * p
3954  else
3955  adjcoln2o = colamt(k,4)
3956  endif
3957 
3958  if (specparm < 0.125) then
3959  p = fs - f_one
3960  p4 = p**4
3961  fk0 = p4
3962  fk1 = f_one - p - 2.0*p4
3963  fk2 = p + p4
3964  id000 = ind0
3965  id010 = ind0 + 9
3966  id100 = ind0 + 1
3967  id110 = ind0 +10
3968  id200 = ind0 + 2
3969  id210 = ind0 +11
3970  else if (specparm > 0.875) then
3971  p = -fs
3972  p4 = p**4
3973  fk0 = p4
3974  fk1 = f_one - p - 2.0*p4
3975  fk2 = p + p4
3976  id000 = ind0 + 1
3977  id010 = ind0 +10
3978  id100 = ind0
3979  id110 = ind0 + 9
3980  id200 = ind0 - 1
3981  id210 = ind0 + 8
3982  else
3983  fk0 = f_one - fs
3984  fk1 = fs
3985  fk2 = f_zero
3986  id000 = ind0
3987  id010 = ind0 + 9
3988  id100 = ind0 + 1
3989  id110 = ind0 +10
3990  id200 = ind0
3991  id210 = ind0
3992  endif
3993 
3994  fac000 = fk0*fac00(k)
3995  fac100 = fk1*fac00(k)
3996  fac200 = fk2*fac00(k)
3997  fac010 = fk0*fac10(k)
3998  fac110 = fk1*fac10(k)
3999  fac210 = fk2*fac10(k)
4000 
4001  if (specparm1 < 0.125) then
4002  p = fs1 - f_one
4003  p4 = p**4
4004  fk0 = p4
4005  fk1 = f_one - p - 2.0*p4
4006  fk2 = p + p4
4007  id001 = ind1
4008  id011 = ind1 + 9
4009  id101 = ind1 + 1
4010  id111 = ind1 +10
4011  id201 = ind1 + 2
4012  id211 = ind1 +11
4013  elseif (specparm1 > 0.875) then
4014  p = -fs1
4015  p4 = p**4
4016  fk0 = p4
4017  fk1 = f_one - p - 2.0*p4
4018  fk2 = p + p4
4019  id001 = ind1 + 1
4020  id011 = ind1 +10
4021  id101 = ind1
4022  id111 = ind1 + 9
4023  id201 = ind1 - 1
4024  id211 = ind1 + 8
4025  else
4026  fk0 = f_one - fs1
4027  fk1 = fs1
4028  fk2 = f_zero
4029  id001 = ind1
4030  id011 = ind1 + 9
4031  id101 = ind1 + 1
4032  id111 = ind1 +10
4033  id201 = ind1
4034  id211 = ind1
4035  endif
4036 
4037  fac001 = fk0*fac01(k)
4038  fac101 = fk1*fac01(k)
4039  fac201 = fk2*fac01(k)
4040  fac011 = fk0*fac11(k)
4041  fac111 = fk1*fac11(k)
4042  fac211 = fk2*fac11(k)
4043 
4044  do ig = 1, ng03
4045  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
4046  & * (selfref(ig,indsp) - selfref(ig,inds)))
4047  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
4048  & * (forref(ig,indfp) - forref(ig,indf)))
4049  n2om1 = ka_mn2o(ig,jmn2o,indm) + fmn2o &
4050  & * (ka_mn2o(ig,jmn2op,indm) - ka_mn2o(ig,jmn2o,indm))
4051  n2om2 = ka_mn2o(ig,jmn2o,indmp) + fmn2o &
4052  & * (ka_mn2o(ig,jmn2op,indmp) - ka_mn2o(ig,jmn2o,indmp))
4053  absn2o = n2om1 + minorfrac(k) * (n2om2 - n2om1)
4054 
4055  tau_major = speccomb &
4056  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
4057  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
4058  & + fac200*absa(ig,id200) + fac210*absa(ig,id210))
4059 
4060  tau_major1 = speccomb1 &
4061  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
4062  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
4063  & + fac201*absa(ig,id201) + fac211*absa(ig,id211))
4064 
4065  taug(ns03+ig,k) = tau_major + tau_major1 &
4066  & + tauself + taufor + adjcoln2o*absn2o
4067 
4068  fracs(ns03+ig,k) = fracrefa(ig,jpl) + fpl &
4069  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
4070  enddo ! end do_k_loop
4071  enddo ! end do_ig_loop
4072 
4073 ! --- ... upper atmosphere loop
4074 
4075  do k = laytrop+1, nlay
4076  speccomb = colamt(k,1) + rfrate(k,1,1)*colamt(k,2)
4077  specparm = colamt(k,1) / speccomb
4078  specmult = 4.0 * min(specparm, oneminus)
4079  js = 1 + int(specmult)
4080  fs = mod(specmult, f_one)
4081  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(3) + js
4082 
4083  speccomb1 = colamt(k,1) + rfrate(k,1,2)*colamt(k,2)
4084  specparm1 = colamt(k,1) / speccomb1
4085  specmult1 = 4.0 * min(specparm1, oneminus)
4086  js1 = 1 + int(specmult1)
4087  fs1 = mod(specmult1, f_one)
4088  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(3) + js1
4089 
4090  speccomb_mn2o = colamt(k,1) + refrat_m_b*colamt(k,2)
4091  specparm_mn2o = colamt(k,1) / speccomb_mn2o
4092  specmult_mn2o = 4.0 * min(specparm_mn2o, oneminus)
4093  jmn2o = 1 + int(specmult_mn2o)
4094  fmn2o = mod(specmult_mn2o, f_one)
4095 
4096  speccomb_planck = colamt(k,1) + refrat_planck_b*colamt(k,2)
4097  specparm_planck = colamt(k,1) / speccomb_planck
4098  specmult_planck = 4.0 * min(specparm_planck, oneminus)
4099  jpl = 1 + int(specmult_planck)
4100  fpl = mod(specmult_planck, f_one)
4101 
4102  indf = indfor(k)
4103  indm = indminor(k)
4104  indfp = indf + 1
4105  indmp = indm + 1
4106  jmn2op= jmn2o+ 1
4107  jplp = jpl + 1
4108 
4109  id000 = ind0
4110  id010 = ind0 + 5
4111  id100 = ind0 + 1
4112  id110 = ind0 + 6
4113  id001 = ind1
4114  id011 = ind1 + 5
4115  id101 = ind1 + 1
4116  id111 = ind1 + 6
4117 
4118 ! --- ... in atmospheres where the amount of n2o is too great to be considered
4119 ! a minor species, adjust the column amount of N2O by an empirical factor
4120 ! to obtain the proper contribution.
4121 
4122  p = coldry(k) * chi_mls(4,jp(k)+1)
4123  ratn2o = colamt(k,4) / p
4124  if (ratn2o > 1.5) then
4125  adjfac = 0.5 + (ratn2o - 0.5)**0.65
4126  adjcoln2o = adjfac * p
4127  else
4128  adjcoln2o = colamt(k,4)
4129  endif
4130 
4131  fk0 = f_one - fs
4132  fk1 = fs
4133  fac000 = fk0*fac00(k)
4134  fac010 = fk0*fac10(k)
4135  fac100 = fk1*fac00(k)
4136  fac110 = fk1*fac10(k)
4137 
4138  fk0 = f_one - fs1
4139  fk1 = fs1
4140  fac001 = fk0*fac01(k)
4141  fac011 = fk0*fac11(k)
4142  fac101 = fk1*fac01(k)
4143  fac111 = fk1*fac11(k)
4144 
4145  do ig = 1, ng03
4146  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
4147  & * (forref(ig,indfp) - forref(ig,indf)))
4148  n2om1 = kb_mn2o(ig,jmn2o,indm) + fmn2o &
4149  & * (kb_mn2o(ig,jmn2op,indm) - kb_mn2o(ig,jmn2o,indm))
4150  n2om2 = kb_mn2o(ig,jmn2o,indmp) + fmn2o &
4151  & * (kb_mn2o(ig,jmn2op,indmp) - kb_mn2o(ig,jmn2o,indmp))
4152  absn2o = n2om1 + minorfrac(k) * (n2om2 - n2om1)
4153 
4154  tau_major = speccomb &
4155  & * (fac000*absb(ig,id000) + fac010*absb(ig,id010) &
4156  & + fac100*absb(ig,id100) + fac110*absb(ig,id110))
4157 
4158  tau_major1 = speccomb1 &
4159  & * (fac001*absb(ig,id001) + fac011*absb(ig,id011) &
4160  & + fac101*absb(ig,id101) + fac111*absb(ig,id111))
4161 
4162  taug(ns03+ig,k) = tau_major + tau_major1 &
4163  & + taufor + adjcoln2o*absn2o
4164 
4165  fracs(ns03+ig,k) = fracrefb(ig,jpl) + fpl &
4166  & * (fracrefb(ig,jplp) - fracrefb(ig,jpl))
4167  enddo
4168  enddo
4169 
4170 ! ..................................
4171  end subroutine taugb03
4172 ! ----------------------------------
4173 
4174 ! ----------------------------------
4175  subroutine taugb04
4176 ! ..................................
4177 
4178 ! ------------------------------------------------------------------ !
4179 ! band 4: 630-700 cm-1 (low key - h2o,co2; high key - o3,co2) !
4180 ! ------------------------------------------------------------------ !
4181 
4182  use module_radlw_kgb04
4183 
4184 ! --- locals:
4185  integer :: k, ind0, ind1, inds, indsp, indf, indfp, jpl, jplp, &
4186  & id000, id010, id100, id110, id200, id210, ig, js, js1, &
4187  & id001, id011, id101, id111, id201, id211
4188 
4189  real (kind=kind_phys) :: tauself, taufor, p, p4, fk0, fk1, fk2, &
4190  & speccomb, specparm, specmult, fs, &
4191  & speccomb1, specparm1, specmult1, fs1, &
4192  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
4193  & fac000, fac100, fac200, fac010, fac110, fac210, &
4194  & fac001, fac101, fac201, fac011, fac111, fac211, &
4195  & refrat_planck_a, refrat_planck_b, tau_major, tau_major1
4196 !
4197 !===> ... begin here
4198 !
4199  refrat_planck_a = chi_mls(1,11)/chi_mls(2,11) ! P = 142.5940 mb
4200  refrat_planck_b = chi_mls(3,13)/chi_mls(2,13) ! P = 95.58350 mb
4201 
4202 ! --- ... lower atmosphere loop
4203 
4204  do k = 1, laytrop
4205  speccomb = colamt(k,1) + rfrate(k,1,1)*colamt(k,2)
4206  specparm = colamt(k,1) / speccomb
4207  specmult = 8.0 * min(specparm, oneminus)
4208  js = 1 + int(specmult)
4209  fs = mod(specmult, f_one)
4210  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(4) + js
4211 
4212  speccomb1 = colamt(k,1) + rfrate(k,1,2)*colamt(k,2)
4213  specparm1 = colamt(k,1) / speccomb1
4214  specmult1 = 8.0 * min(specparm1, oneminus)
4215  js1 = 1 + int(specmult1)
4216  fs1 = mod(specmult1, f_one)
4217  ind1 = ( jp(k)*5 + (jt1(k)-1)) * nspa(4) + js1
4218 
4219  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,2)
4220  specparm_planck = colamt(k,1) / speccomb_planck
4221  specmult_planck = 8.0 * min(specparm_planck, oneminus)
4222  jpl = 1 + int(specmult_planck)
4223  fpl = mod(specmult_planck, 1.0)
4224 
4225  inds = indself(k)
4226  indf = indfor(k)
4227  indsp = inds + 1
4228  indfp = indf + 1
4229  jplp = jpl + 1
4230 
4231  if (specparm < 0.125) then
4232  p = fs - f_one
4233  p4 = p**4
4234  fk0 = p4
4235  fk1 = f_one - p - 2.0*p4
4236  fk2 = p + p4
4237  id000 = ind0
4238  id010 = ind0 + 9
4239  id100 = ind0 + 1
4240  id110 = ind0 +10
4241  id200 = ind0 + 2
4242  id210 = ind0 +11
4243  elseif (specparm > 0.875) then
4244  p = -fs
4245  p4 = p**4
4246  fk0 = p4
4247  fk1 = f_one - p - 2.0*p4
4248  fk2 = p + p4
4249  id000 = ind0 + 1
4250  id010 = ind0 +10
4251  id100 = ind0
4252  id110 = ind0 + 9
4253  id200 = ind0 - 1
4254  id210 = ind0 + 8
4255  else
4256  fk0 = f_one - fs
4257  fk1 = fs
4258  fk2 = f_zero
4259  id000 = ind0
4260  id010 = ind0 + 9
4261  id100 = ind0 + 1
4262  id110 = ind0 +10
4263  id200 = ind0
4264  id210 = ind0
4265  endif
4266 
4267  fac000 = fk0*fac00(k)
4268  fac100 = fk1*fac00(k)
4269  fac200 = fk2*fac00(k)
4270  fac010 = fk0*fac10(k)
4271  fac110 = fk1*fac10(k)
4272  fac210 = fk2*fac10(k)
4273 
4274  if (specparm1 < 0.125) then
4275  p = fs1 - f_one
4276  p4 = p**4
4277  fk0 = p4
4278  fk1 = f_one - p - 2.0*p4
4279  fk2 = p + p4
4280  id001 = ind1
4281  id011 = ind1 + 9
4282  id101 = ind1 + 1
4283  id111 = ind1 +10
4284  id201 = ind1 + 2
4285  id211 = ind1 +11
4286  elseif (specparm1 > 0.875) then
4287  p = -fs1
4288  p4 = p**4
4289  fk0 = p4
4290  fk1 = f_one - p - 2.0*p4
4291  fk2 = p + p4
4292  id001 = ind1 + 1
4293  id011 = ind1 +10
4294  id101 = ind1
4295  id111 = ind1 + 9
4296  id201 = ind1 - 1
4297  id211 = ind1 + 8
4298  else
4299  fk0 = f_one - fs1
4300  fk1 = fs1
4301  fk2 = f_zero
4302  id001 = ind1
4303  id011 = ind1 + 9
4304  id101 = ind1 + 1
4305  id111 = ind1 +10
4306  id201 = ind1
4307  id211 = ind1
4308  endif
4309 
4310  fac001 = fk0*fac01(k)
4311  fac101 = fk1*fac01(k)
4312  fac201 = fk2*fac01(k)
4313  fac011 = fk0*fac11(k)
4314  fac111 = fk1*fac11(k)
4315  fac211 = fk2*fac11(k)
4316 
4317  do ig = 1, ng04
4318  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
4319  & * (selfref(ig,indsp) - selfref(ig,inds)))
4320  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
4321  & * (forref(ig,indfp) - forref(ig,indf)))
4322 
4323  tau_major = speccomb &
4324  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
4325  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
4326  & + fac200*absa(ig,id200) + fac210*absa(ig,id210))
4327 
4328  tau_major1 = speccomb1 &
4329  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
4330  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
4331  & + fac201*absa(ig,id201) + fac211*absa(ig,id211))
4332 
4333  taug(ns04+ig,k) = tau_major + tau_major1 + tauself + taufor
4334 
4335  fracs(ns04+ig,k) = fracrefa(ig,jpl) + fpl &
4336  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
4337  enddo ! end do_k_loop
4338  enddo ! end do_ig_loop
4339 
4340 ! --- ... upper atmosphere loop
4341 
4342  do k = laytrop+1, nlay
4343  speccomb = colamt(k,3) + rfrate(k,6,1)*colamt(k,2)
4344  specparm = colamt(k,3) / speccomb
4345  specmult = 4.0 * min(specparm, oneminus)
4346  js = 1 + int(specmult)
4347  fs = mod(specmult, f_one)
4348  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(4) + js
4349 
4350  speccomb1 = colamt(k,3) + rfrate(k,6,2)*colamt(k,2)
4351  specparm1 = colamt(k,3) / speccomb1
4352  specmult1 = 4.0 * min(specparm1, oneminus)
4353  js1 = 1 + int(specmult1)
4354  fs1 = mod(specmult1, f_one)
4355  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(4) + js1
4356 
4357  speccomb_planck = colamt(k,3) + refrat_planck_b*colamt(k,2)
4358  specparm_planck = colamt(k,3) / speccomb_planck
4359  specmult_planck = 4.0 * min(specparm_planck, oneminus)
4360  jpl = 1 + int(specmult_planck)
4361  fpl = mod(specmult_planck, f_one)
4362  jplp = jpl + 1
4363 
4364  id000 = ind0
4365  id010 = ind0 + 5
4366  id100 = ind0 + 1
4367  id110 = ind0 + 6
4368  id001 = ind1
4369  id011 = ind1 + 5
4370  id101 = ind1 + 1
4371  id111 = ind1 + 6
4372 
4373  fk0 = f_one - fs
4374  fk1 = fs
4375  fac000 = fk0*fac00(k)
4376  fac010 = fk0*fac10(k)
4377  fac100 = fk1*fac00(k)
4378  fac110 = fk1*fac10(k)
4379 
4380  fk0 = f_one - fs1
4381  fk1 = fs1
4382  fac001 = fk0*fac01(k)
4383  fac011 = fk0*fac11(k)
4384  fac101 = fk1*fac01(k)
4385  fac111 = fk1*fac11(k)
4386 
4387  do ig = 1, ng04
4388  tau_major = speccomb &
4389  & * (fac000*absb(ig,id000) + fac010*absb(ig,id010) &
4390  & + fac100*absb(ig,id100) + fac110*absb(ig,id110))
4391  tau_major1 = speccomb1 &
4392  & * (fac001*absb(ig,id001) + fac011*absb(ig,id011) &
4393  & + fac101*absb(ig,id101) + fac111*absb(ig,id111))
4394 
4395  taug(ns04+ig,k) = tau_major + tau_major1
4396 
4397  fracs(ns04+ig,k) = fracrefb(ig,jpl) + fpl &
4398  & * (fracrefb(ig,jplp) - fracrefb(ig,jpl))
4399  enddo
4400 
4401 ! --- ... empirical modification to code to improve stratospheric cooling rates
4402 ! for co2. revised to apply weighting for g-point reduction in this band.
4403 
4404  taug(ns04+ 8,k) = taug(ns04+ 8,k) * 0.92
4405  taug(ns04+ 9,k) = taug(ns04+ 9,k) * 0.88
4406  taug(ns04+10,k) = taug(ns04+10,k) * 1.07
4407  taug(ns04+11,k) = taug(ns04+11,k) * 1.1
4408  taug(ns04+12,k) = taug(ns04+12,k) * 0.99
4409  taug(ns04+13,k) = taug(ns04+13,k) * 0.88
4410  taug(ns04+14,k) = taug(ns04+14,k) * 0.943
4411  enddo
4412 
4413 ! ..................................
4414  end subroutine taugb04
4415 ! ----------------------------------
4416 
4417 ! ----------------------------------
4418  subroutine taugb05
4419 ! ..................................
4420 
4421 ! ------------------------------------------------------------------ !
4422 ! band 5: 700-820 cm-1 (low key - h2o,co2; low minor - o3, ccl4) !
4423 ! (high key - o3,co2) !
4424 ! ------------------------------------------------------------------ !
4425 
4426  use module_radlw_kgb05
4427 
4428 ! --- locals:
4429  integer :: k, ind0, ind1, inds, indsp, indf, indfp, indm, indmp, &
4430  & id000, id010, id100, id110, id200, id210, jmo3, jmo3p, &
4431  & id001, id011, id101, id111, id201, id211, jpl, jplp, &
4432  & ig, js, js1
4433 
4434  real (kind=kind_phys) :: tauself, taufor, o3m1, o3m2, abso3, &
4435  & speccomb, specparm, specmult, fs, &
4436  & speccomb1, specparm1, specmult1, fs1, &
4437  & speccomb_mo3, specparm_mo3, specmult_mo3, fmo3, &
4438  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
4439  & refrat_planck_a, refrat_planck_b, refrat_m_a, &
4440  & fac000, fac100, fac200, fac010, fac110, fac210, &
4441  & fac001, fac101, fac201, fac011, fac111, fac211, &
4442  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21
4443 !
4444 !===> ... begin here
4445 !
4446 ! --- ... minor gas mapping level :
4447 ! lower - o3, p = 317.34 mbar, t = 240.77 k
4448 ! lower - ccl4
4449 
4450 ! --- ... calculate reference ratio to be used in calculation of Planck
4451 ! fraction in lower/upper atmosphere.
4452 
4453  refrat_planck_a = chi_mls(1,5)/chi_mls(2,5) ! P = 473.420 mb
4454  refrat_planck_b = chi_mls(3,43)/chi_mls(2,43) ! P = 0.2369 mb
4455  refrat_m_a = chi_mls(1,7)/chi_mls(2,7) ! P = 317.348 mb
4456 
4457 ! --- ... lower atmosphere loop
4458 
4459  do k = 1, laytrop
4460  speccomb = colamt(k,1) + rfrate(k,1,1)*colamt(k,2)
4461  specparm = colamt(k,1) / speccomb
4462  specmult = 8.0 * min(specparm, oneminus)
4463  js = 1 + int(specmult)
4464  fs = mod(specmult, f_one)
4465  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(5) + js
4466 
4467  speccomb1 = colamt(k,1) + rfrate(k,1,2)*colamt(k,2)
4468  specparm1 = colamt(k,1) / speccomb1
4469  specmult1 = 8.0 * min(specparm1, oneminus)
4470  js1 = 1 + int(specmult1)
4471  fs1 = mod(specmult1, f_one)
4472  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(5) + js1
4473 
4474  speccomb_mo3 = colamt(k,1) + refrat_m_a*colamt(k,2)
4475  specparm_mo3 = colamt(k,1) / speccomb_mo3
4476  specmult_mo3 = 8.0 * min(specparm_mo3, oneminus)
4477  jmo3 = 1 + int(specmult_mo3)
4478  fmo3 = mod(specmult_mo3, f_one)
4479 
4480  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,2)
4481  specparm_planck = colamt(k,1) / speccomb_planck
4482  specmult_planck = 8.0 * min(specparm_planck, oneminus)
4483  jpl = 1 + int(specmult_planck)
4484  fpl = mod(specmult_planck, f_one)
4485 
4486  inds = indself(k)
4487  indf = indfor(k)
4488  indm = indminor(k)
4489  indsp = inds + 1
4490  indfp = indf + 1
4491  indmp = indm + 1
4492  jplp = jpl + 1
4493  jmo3p = jmo3 + 1
4494 
4495  if (specparm < 0.125 .and. specparm1 < 0.125) then
4496  p0 = fs - f_one
4497  p40 = p0**4
4498  fk00 = p40
4499  fk10 = f_one - p0 - 2.0*p40
4500  fk20 = p0 + p40
4501 
4502  p1 = fs1 - f_one
4503  p41 = p1**4
4504  fk01 = p41
4505  fk11 = f_one - p1 - 2.0*p41
4506  fk21 = p1 + p41
4507 
4508  id000 = ind0
4509  id010 = ind0 + 9
4510  id100 = ind0 + 1
4511  id110 = ind0 +10
4512  id200 = ind0 + 2
4513  id210 = ind0 +11
4514 
4515  id001 = ind1
4516  id011 = ind1 + 9
4517  id101 = ind1 + 1
4518  id111 = ind1 +10
4519  id201 = ind1 + 2
4520  id211 = ind1 +11
4521  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
4522  p0 = -fs
4523  p40 = p0**4
4524  fk00 = p40
4525  fk10 = f_one - p0 - 2.0*p40
4526  fk20 = p0 + p40
4527 
4528  p1 = -fs1
4529  p41 = p1**4
4530  fk01 = p41
4531  fk11 = f_one - p1 - 2.0*p41
4532  fk21 = p1 + p41
4533 
4534  id000 = ind0 + 1
4535  id010 = ind0 +10
4536  id100 = ind0
4537  id110 = ind0 + 9
4538  id200 = ind0 - 1
4539  id210 = ind0 + 8
4540 
4541  id001 = ind1 + 1
4542  id011 = ind1 +10
4543  id101 = ind1
4544  id111 = ind1 + 9
4545  id201 = ind1 - 1
4546  id211 = ind1 + 8
4547  else
4548  fk00 = f_one - fs
4549  fk10 = fs
4550  fk20 = f_zero
4551 
4552  fk01 = f_one - fs1
4553  fk11 = fs1
4554  fk21 = f_zero
4555 
4556  id000 = ind0
4557  id010 = ind0 + 9
4558  id100 = ind0 + 1
4559  id110 = ind0 +10
4560  id200 = ind0
4561  id210 = ind0
4562 
4563  id001 = ind1
4564  id011 = ind1 + 9
4565  id101 = ind1 + 1
4566  id111 = ind1 +10
4567  id201 = ind1
4568  id211 = ind1
4569  endif
4570 
4571  fac000 = fk00 * fac00(k)
4572  fac100 = fk10 * fac00(k)
4573  fac200 = fk20 * fac00(k)
4574  fac010 = fk00 * fac10(k)
4575  fac110 = fk10 * fac10(k)
4576  fac210 = fk20 * fac10(k)
4577 
4578  fac001 = fk01 * fac01(k)
4579  fac101 = fk11 * fac01(k)
4580  fac201 = fk21 * fac01(k)
4581  fac011 = fk01 * fac11(k)
4582  fac111 = fk11 * fac11(k)
4583  fac211 = fk21 * fac11(k)
4584 
4585  do ig = 1, ng05
4586  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
4587  & * (selfref(ig,indsp) - selfref(ig,inds)))
4588  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
4589  & * (forref(ig,indfp) - forref(ig,indf)))
4590  o3m1 = ka_mo3(ig,jmo3,indm) + fmo3 &
4591  & * (ka_mo3(ig,jmo3p,indm) - ka_mo3(ig,jmo3,indm))
4592  o3m2 = ka_mo3(ig,jmo3,indmp) + fmo3 &
4593  & * (ka_mo3(ig,jmo3p,indmp) - ka_mo3(ig,jmo3,indmp))
4594  abso3 = o3m1 + minorfrac(k)*(o3m2 - o3m1)
4595 
4596  taug(ns05+ig,k) = speccomb &
4597  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
4598  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
4599  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
4600  & + speccomb1 &
4601  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
4602  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
4603  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
4604  & + tauself + taufor+abso3*colamt(k,3)+wx(k,1)*ccl4(ig)
4605 
4606  fracs(ns05+ig,k) = fracrefa(ig,jpl) + fpl &
4607  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
4608  enddo
4609  enddo
4610 
4611 ! --- ... upper atmosphere loop
4612 
4613  do k = laytrop+1, nlay
4614  speccomb = colamt(k,3) + rfrate(k,6,1)*colamt(k,2)
4615  specparm = colamt(k,3) / speccomb
4616  specmult = 4.0 * min(specparm, oneminus)
4617  js = 1 + int(specmult)
4618  fs = mod(specmult, f_one)
4619  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(5) + js
4620 
4621  speccomb1 = colamt(k,3) + rfrate(k,6,2)*colamt(k,2)
4622  specparm1 = colamt(k,3) / speccomb1
4623  specmult1 = 4.0 * min(specparm1, oneminus)
4624  js1 = 1 + int(specmult1)
4625  fs1 = mod(specmult1, f_one)
4626  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(5) + js1
4627 
4628  speccomb_planck = colamt(k,3) + refrat_planck_b*colamt(k,2)
4629  specparm_planck = colamt(k,3) / speccomb_planck
4630  specmult_planck = 4.0 * min(specparm_planck, oneminus)
4631  jpl = 1 + int(specmult_planck)
4632  fpl = mod(specmult_planck, f_one)
4633  jplp= jpl + 1
4634 
4635  id000 = ind0
4636  id010 = ind0 + 5
4637  id100 = ind0 + 1
4638  id110 = ind0 + 6
4639  id001 = ind1
4640  id011 = ind1 + 5
4641  id101 = ind1 + 1
4642  id111 = ind1 + 6
4643 
4644  fk00 = f_one - fs
4645  fk10 = fs
4646 
4647  fk01 = f_one - fs1
4648  fk11 = fs1
4649 
4650  fac000 = fk00 * fac00(k)
4651  fac010 = fk00 * fac10(k)
4652  fac100 = fk10 * fac00(k)
4653  fac110 = fk10 * fac10(k)
4654 
4655  fac001 = fk01 * fac01(k)
4656  fac011 = fk01 * fac11(k)
4657  fac101 = fk11 * fac01(k)
4658  fac111 = fk11 * fac11(k)
4659 
4660  do ig = 1, ng05
4661  taug(ns05+ig,k) = speccomb &
4662  & * (fac000*absb(ig,id000) + fac010*absb(ig,id010) &
4663  & + fac100*absb(ig,id100) + fac110*absb(ig,id110)) &
4664  & + speccomb1 &
4665  & * (fac001*absb(ig,id001) + fac011*absb(ig,id011) &
4666  & + fac101*absb(ig,id101) + fac111*absb(ig,id111)) &
4667  & + wx(k,1) * ccl4(ig)
4668 
4669  fracs(ns05+ig,k) = fracrefb(ig,jpl) + fpl &
4670  & * (fracrefb(ig,jplp) - fracrefb(ig,jpl))
4671  enddo
4672  enddo
4673 
4674 ! ..................................
4675  end subroutine taugb05
4676 ! ----------------------------------
4677 
4678 ! ----------------------------------
4679  subroutine taugb06
4680 ! ..................................
4681 
4682 ! ------------------------------------------------------------------ !
4683 ! band 6: 820-980 cm-1 (low key - h2o; low minor - co2) !
4684 ! (high key - none; high minor - cfc11, cfc12)
4685 ! ------------------------------------------------------------------ !
4686 
4687  use module_radlw_kgb06
4688 
4689 ! --- locals:
4690  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
4691  & indm, indmp, ig
4692 
4693  real (kind=kind_phys) :: ratco2, adjfac, adjcolco2, tauself, &
4694  & taufor, absco2, temp
4695 !
4696 !===> ... begin here
4697 !
4698 ! --- ... minor gas mapping level:
4699 ! lower - co2, p = 706.2720 mb, t = 294.2 k
4700 ! upper - cfc11, cfc12
4701 
4702 ! --- ... lower atmosphere loop
4703 
4704  do k = 1, laytrop
4705  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(6) + 1
4706  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(6) + 1
4707 
4708  inds = indself(k)
4709  indf = indfor(k)
4710  indm = indminor(k)
4711  indsp = inds + 1
4712  indfp = indf + 1
4713  indmp = indm + 1
4714  ind0p = ind0 + 1
4715  ind1p = ind1 + 1
4716 
4717 ! --- ... in atmospheres where the amount of co2 is too great to be considered
4718 ! a minor species, adjust the column amount of co2 by an empirical factor
4719 ! to obtain the proper contribution.
4720 
4721  temp = coldry(k) * chi_mls(2,jp(k)+1)
4722  ratco2 = colamt(k,2) / temp
4723  if (ratco2 > 3.0) then
4724  adjfac = 2.0 + (ratco2-2.0)**0.77
4725  adjcolco2 = adjfac * temp
4726  else
4727  adjcolco2 = colamt(k,2)
4728  endif
4729 
4730  do ig = 1, ng06
4731  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
4732  & * (selfref(ig,indsp) - selfref(ig,inds)))
4733  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
4734  & * (forref(ig,indfp) - forref(ig,indf)))
4735  absco2 = ka_mco2(ig,indm) + minorfrac(k) &
4736  & * (ka_mco2(ig,indmp) - ka_mco2(ig,indm))
4737 
4738  taug(ns06+ig,k) = colamt(k,1) &
4739  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
4740  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
4741  & + tauself + taufor + adjcolco2*absco2 &
4742  & + wx(k,2)*cfc11adj(ig) + wx(k,3)*cfc12(ig)
4743 
4744  fracs(ns06+ig,k) = fracrefa(ig)
4745  enddo
4746  enddo
4747 
4748 ! --- ... upper atmosphere loop
4749 ! nothing important goes on above laytrop in this band.
4750 
4751  do k = laytrop+1, nlay
4752  do ig = 1, ng06
4753  taug(ns06+ig,k) = wx(k,2)*cfc11adj(ig) + wx(k,3)*cfc12(ig)
4754 
4755  fracs(ns06+ig,k) = fracrefa(ig)
4756  enddo
4757  enddo
4758 
4759 ! ..................................
4760  end subroutine taugb06
4761 ! ----------------------------------
4762 
4763 ! ----------------------------------
4764  subroutine taugb07
4765 ! ..................................
4766 
4767 ! ------------------------------------------------------------------ !
4768 ! band 7: 980-1080 cm-1 (low key - h2o,o3; low minor - co2) !
4769 ! (high key - o3; high minor - co2) !
4770 ! ------------------------------------------------------------------ !
4771 
4772  use module_radlw_kgb07
4773 
4774 ! --- locals:
4775  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
4776  & id000, id010, id100, id110, id200, id210, indm, indmp, &
4777  & id001, id011, id101, id111, id201, id211, jmco2, jmco2p, &
4778  & jpl, jplp, ig, js, js1
4779 
4780  real (kind=kind_phys) :: tauself, taufor, co2m1, co2m2, absco2, &
4781  & speccomb, specparm, specmult, fs, &
4782  & speccomb1, specparm1, specmult1, fs1, &
4783  & speccomb_mco2, specparm_mco2, specmult_mco2, fmco2, &
4784  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
4785  & refrat_planck_a, refrat_m_a, ratco2, adjfac, adjcolco2, &
4786  & fac000, fac100, fac200, fac010, fac110, fac210, &
4787  & fac001, fac101, fac201, fac011, fac111, fac211, &
4788  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21, temp
4789 !
4790 !===> ... begin here
4791 !
4792 ! --- ... minor gas mapping level :
4793 ! lower - co2, p = 706.2620 mbar, t= 278.94 k
4794 ! upper - co2, p = 12.9350 mbar, t = 234.01 k
4795 
4796 ! --- ... calculate reference ratio to be used in calculation of Planck
4797 ! fraction in lower atmosphere.
4798 
4799  refrat_planck_a = chi_mls(1,3)/chi_mls(3,3) ! P = 706.2620 mb
4800  refrat_m_a = chi_mls(1,3)/chi_mls(3,3) ! P = 706.2720 mb
4801 
4802 ! --- ... lower atmosphere loop
4803 
4804  do k = 1, laytrop
4805  speccomb = colamt(k,1) + rfrate(k,2,1)*colamt(k,3)
4806  specparm = colamt(k,1) / speccomb
4807  specmult = 8.0 * min(specparm, oneminus)
4808  js = 1 + int(specmult)
4809  fs = mod(specmult, f_one)
4810  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(7) + js
4811 
4812  speccomb1 = colamt(k,1) + rfrate(k,2,2)*colamt(k,3)
4813  specparm1 = colamt(k,1) / speccomb1
4814  specmult1 = 8.0 * min(specparm1, oneminus)
4815  js1 = 1 + int(specmult1)
4816  fs1 = mod(specmult1, f_one)
4817  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(7) + js1
4818 
4819  speccomb_mco2 = colamt(k,1) + refrat_m_a*colamt(k,3)
4820  specparm_mco2 = colamt(k,1) / speccomb_mco2
4821  specmult_mco2 = 8.0 * min(specparm_mco2, oneminus)
4822  jmco2 = 1 + int(specmult_mco2)
4823  fmco2 = mod(specmult_mco2, f_one)
4824 
4825  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,3)
4826  specparm_planck = colamt(k,1) / speccomb_planck
4827  specmult_planck = 8.0 * min(specparm_planck, oneminus)
4828  jpl = 1 + int(specmult_planck)
4829  fpl = mod(specmult_planck, f_one)
4830 
4831  inds = indself(k)
4832  indf = indfor(k)
4833  indm = indminor(k)
4834  indsp = inds + 1
4835  indfp = indf + 1
4836  indmp = indm + 1
4837  jplp = jpl + 1
4838  jmco2p= jmco2+ 1
4839  ind0p = ind0 + 1
4840  ind1p = ind1 + 1
4841 
4842 ! --- ... in atmospheres where the amount of CO2 is too great to be considered
4843 ! a minor species, adjust the column amount of CO2 by an empirical factor
4844 ! to obtain the proper contribution.
4845 
4846  temp = coldry(k) * chi_mls(2,jp(k)+1)
4847  ratco2 = colamt(k,2) / temp
4848  if (ratco2 > 3.0) then
4849  adjfac = 3.0 + (ratco2-3.0)**0.79
4850  adjcolco2 = adjfac * temp
4851  else
4852  adjcolco2 = colamt(k,2)
4853  endif
4854 
4855  if (specparm < 0.125 .and. specparm1 < 0.125) then
4856  p0 = fs - f_one
4857  p40 = p0**4
4858  fk00 = p40
4859  fk10 = f_one - p0 - 2.0*p40
4860  fk20 = p0 + p40
4861 
4862  p1 = fs1 - f_one
4863  p41 = p1**4
4864  fk01 = p41
4865  fk11 = f_one - p1 - 2.0*p41
4866  fk21 = p1 + p41
4867 
4868  id000 = ind0
4869  id010 = ind0 + 9
4870  id100 = ind0 + 1
4871  id110 = ind0 +10
4872  id200 = ind0 + 2
4873  id210 = ind0 +11
4874 
4875  id001 = ind1
4876  id011 = ind1 + 9
4877  id101 = ind1 + 1
4878  id111 = ind1 +10
4879  id201 = ind1 + 2
4880  id211 = ind1 +11
4881  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
4882  p0 = -fs
4883  p40 = p0**4
4884  fk00 = p40
4885  fk10 = f_one - p0 - 2.0*p40
4886  fk20 = p0 + p40
4887 
4888  p1 = -fs1
4889  p41 = p1**4
4890  fk01 = p41
4891  fk11 = f_one - p1 - 2.0*p41
4892  fk21 = p1 + p41
4893 
4894  id000 = ind0 + 1
4895  id010 = ind0 +10
4896  id100 = ind0
4897  id110 = ind0 + 9
4898  id200 = ind0 - 1
4899  id210 = ind0 + 8
4900 
4901  id001 = ind1 + 1
4902  id011 = ind1 +10
4903  id101 = ind1
4904  id111 = ind1 + 9
4905  id201 = ind1 - 1
4906  id211 = ind1 + 8
4907  else
4908  fk00 = f_one - fs
4909  fk10 = fs
4910  fk20 = f_zero
4911 
4912  fk01 = f_one - fs1
4913  fk11 = fs1
4914  fk21 = f_zero
4915 
4916  id000 = ind0
4917  id010 = ind0 + 9
4918  id100 = ind0 + 1
4919  id110 = ind0 +10
4920  id200 = ind0
4921  id210 = ind0
4922 
4923  id001 = ind1
4924  id011 = ind1 + 9
4925  id101 = ind1 + 1
4926  id111 = ind1 +10
4927  id201 = ind1
4928  id211 = ind1
4929  endif
4930 
4931  fac000 = fk00 * fac00(k)
4932  fac100 = fk10 * fac00(k)
4933  fac200 = fk20 * fac00(k)
4934  fac010 = fk00 * fac10(k)
4935  fac110 = fk10 * fac10(k)
4936  fac210 = fk20 * fac10(k)
4937 
4938  fac001 = fk01 * fac01(k)
4939  fac101 = fk11 * fac01(k)
4940  fac201 = fk21 * fac01(k)
4941  fac011 = fk01 * fac11(k)
4942  fac111 = fk11 * fac11(k)
4943  fac211 = fk21 * fac11(k)
4944 
4945  do ig = 1, ng07
4946  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
4947  & * (selfref(ig,indsp) - selfref(ig,inds)))
4948  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
4949  & * (forref(ig,indfp) - forref(ig,indf)))
4950  co2m1 = ka_mco2(ig,jmco2,indm) + fmco2 &
4951  & * (ka_mco2(ig,jmco2p,indm) - ka_mco2(ig,jmco2,indm))
4952  co2m2 = ka_mco2(ig,jmco2,indmp) + fmco2 &
4953  & * (ka_mco2(ig,jmco2p,indmp) - ka_mco2(ig,jmco2,indmp))
4954  absco2 = co2m1 + minorfrac(k) * (co2m2 - co2m1)
4955 
4956  taug(ns07+ig,k) = speccomb &
4957  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
4958  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
4959  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
4960  & + speccomb1 &
4961  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
4962  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
4963  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
4964  & + tauself + taufor + adjcolco2*absco2
4965 
4966  fracs(ns07+ig,k) = fracrefa(ig,jpl) + fpl &
4967  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
4968  enddo
4969  enddo
4970 
4971 ! --- ... upper atmosphere loop
4972 
4973 ! --- ... in atmospheres where the amount of co2 is too great to be considered
4974 ! a minor species, adjust the column amount of co2 by an empirical factor
4975 ! to obtain the proper contribution.
4976 
4977  do k = laytrop+1, nlay
4978  temp = coldry(k) * chi_mls(2,jp(k)+1)
4979  ratco2 = colamt(k,2) / temp
4980  if (ratco2 > 3.0) then
4981  adjfac = 2.0 + (ratco2-2.0)**0.79
4982  adjcolco2 = adjfac * temp
4983  else
4984  adjcolco2 = colamt(k,2)
4985  endif
4986 
4987  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(7) + 1
4988  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(7) + 1
4989 
4990  indm = indminor(k)
4991  indmp = indm + 1
4992  ind0p = ind0 + 1
4993  ind1p = ind1 + 1
4994 
4995  do ig = 1, ng07
4996  absco2 = kb_mco2(ig,indm) + minorfrac(k) &
4997  & * (kb_mco2(ig,indmp) - kb_mco2(ig,indm))
4998 
4999  taug(ns07+ig,k) = colamt(k,3) &
5000  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
5001  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
5002  & + adjcolco2 * absco2
5003 
5004  fracs(ns07+ig,k) = fracrefb(ig)
5005  enddo
5006 
5007 ! --- ... empirical modification to code to improve stratospheric cooling rates
5008 ! for o3. revised to apply weighting for g-point reduction in this band.
5009 
5010  taug(ns07+ 6,k) = taug(ns07+ 6,k) * 0.92
5011  taug(ns07+ 7,k) = taug(ns07+ 7,k) * 0.88
5012  taug(ns07+ 8,k) = taug(ns07+ 8,k) * 1.07
5013  taug(ns07+ 9,k) = taug(ns07+ 9,k) * 1.1
5014  taug(ns07+10,k) = taug(ns07+10,k) * 0.99
5015  taug(ns07+11,k) = taug(ns07+11,k) * 0.855
5016  enddo
5017 
5018 ! ..................................
5019  end subroutine taugb07
5020 ! ----------------------------------
5021 
5022 ! ----------------------------------
5023  subroutine taugb08
5024 ! ..................................
5025 
5026 ! ------------------------------------------------------------------ !
5027 ! band 8: 1080-1180 cm-1 (low key - h2o; low minor - co2,o3,n2o) !
5028 ! (high key - o3; high minor - co2, n2o) !
5029 ! ------------------------------------------------------------------ !
5030 
5031  use module_radlw_kgb08
5032 
5033 ! --- locals:
5034  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
5035  & indm, indmp, ig
5036 
5037  real (kind=kind_phys) :: tauself, taufor, absco2, abso3, absn2o, &
5038  & ratco2, adjfac, adjcolco2, temp
5039 !
5040 !===> ... begin here
5041 !
5042 ! --- ... minor gas mapping level:
5043 ! lower - co2, p = 1053.63 mb, t = 294.2 k
5044 ! lower - o3, p = 317.348 mb, t = 240.77 k
5045 ! lower - n2o, p = 706.2720 mb, t= 278.94 k
5046 ! lower - cfc12,cfc11
5047 ! upper - co2, p = 35.1632 mb, t = 223.28 k
5048 ! upper - n2o, p = 8.716e-2 mb, t = 226.03 k
5049 
5050 ! --- ... lower atmosphere loop
5051 
5052  do k = 1, laytrop
5053  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(8) + 1
5054  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(8) + 1
5055 
5056  inds = indself(k)
5057  indf = indfor(k)
5058  indm = indminor(k)
5059  ind0p = ind0 + 1
5060  ind1p = ind1 + 1
5061  indsp = inds + 1
5062  indfp = indf + 1
5063  indmp = indm + 1
5064 
5065 ! --- ... in atmospheres where the amount of co2 is too great to be considered
5066 ! a minor species, adjust the column amount of co2 by an empirical factor
5067 ! to obtain the proper contribution.
5068 
5069  temp = coldry(k) * chi_mls(2,jp(k)+1)
5070  ratco2 = colamt(k,2) / temp
5071  if (ratco2 > 3.0) then
5072  adjfac = 2.0 + (ratco2-2.0)**0.65
5073  adjcolco2 = adjfac * temp
5074  else
5075  adjcolco2 = colamt(k,2)
5076  endif
5077 
5078  do ig = 1, ng08
5079  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
5080  & * (selfref(ig,indsp) - selfref(ig,inds)))
5081  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5082  & * (forref(ig,indfp) - forref(ig,indf)))
5083  absco2 = (ka_mco2(ig,indm) + minorfrac(k) &
5084  & * (ka_mco2(ig,indmp) - ka_mco2(ig,indm)))
5085  abso3 = (ka_mo3(ig,indm) + minorfrac(k) &
5086  & * (ka_mo3(ig,indmp) - ka_mo3(ig,indm)))
5087  absn2o = (ka_mn2o(ig,indm) + minorfrac(k) &
5088  & * (ka_mn2o(ig,indmp) - ka_mn2o(ig,indm)))
5089 
5090  taug(ns08+ig,k) = colamt(k,1) &
5091  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
5092  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
5093  & + tauself+taufor + adjcolco2*absco2 &
5094  & + colamt(k,3)*abso3 + colamt(k,4)*absn2o &
5095  & + wx(k,3)*cfc12(ig) + wx(k,4)*cfc22adj(ig)
5096 
5097  fracs(ns08+ig,k) = fracrefa(ig)
5098  enddo
5099  enddo
5100 
5101 ! --- ... upper atmosphere loop
5102 
5103  do k = laytrop+1, nlay
5104  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(8) + 1
5105  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(8) + 1
5106 
5107  indm = indminor(k)
5108  ind0p = ind0 + 1
5109  ind1p = ind1 + 1
5110  indmp = indm + 1
5111 
5112 ! --- ... in atmospheres where the amount of co2 is too great to be considered
5113 ! a minor species, adjust the column amount of co2 by an empirical factor
5114 ! to obtain the proper contribution.
5115 
5116  temp = coldry(k) * chi_mls(2,jp(k)+1)
5117  ratco2 = colamt(k,2) / temp
5118  if (ratco2 > 3.0) then
5119  adjfac = 2.0 + (ratco2-2.0)**0.65
5120  adjcolco2 = adjfac * temp
5121  else
5122  adjcolco2 = colamt(k,2)
5123  endif
5124 
5125  do ig = 1, ng08
5126  absco2 = (kb_mco2(ig,indm) + minorfrac(k) &
5127  & * (kb_mco2(ig,indmp) - kb_mco2(ig,indm)))
5128  absn2o = (kb_mn2o(ig,indm) + minorfrac(k) &
5129  & * (kb_mn2o(ig,indmp) - kb_mn2o(ig,indm)))
5130 
5131  taug(ns08+ig,k) = colamt(k,3) &
5132  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
5133  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
5134  & + adjcolco2*absco2 + colamt(k,4)*absn2o &
5135  & + wx(k,3)*cfc12(ig) + wx(k,4)*cfc22adj(ig)
5136 
5137  fracs(ns08+ig,k) = fracrefb(ig)
5138  enddo
5139  enddo
5140 
5141 ! ..................................
5142  end subroutine taugb08
5143 ! ----------------------------------
5144 
5145 ! ----------------------------------
5146  subroutine taugb09
5147 ! ..................................
5148 
5149 ! ------------------------------------------------------------------ !
5150 ! band 9: 1180-1390 cm-1 (low key - h2o,ch4; low minor - n2o) !
5151 ! (high key - ch4; high minor - n2o) !
5152 ! ------------------------------------------------------------------ !
5153 
5154  use module_radlw_kgb09
5155 
5156 ! --- locals:
5157  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
5158  & id000, id010, id100, id110, id200, id210, indm, indmp, &
5159  & id001, id011, id101, id111, id201, id211, jmn2o, jmn2op, &
5160  & jpl, jplp, ig, js, js1
5161 
5162  real (kind=kind_phys) :: tauself, taufor, n2om1, n2om2, absn2o, &
5163  & speccomb, specparm, specmult, fs, &
5164  & speccomb1, specparm1, specmult1, fs1, &
5165  & speccomb_mn2o, specparm_mn2o, specmult_mn2o, fmn2o, &
5166  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
5167  & refrat_planck_a, refrat_m_a, ratn2o, adjfac, adjcoln2o, &
5168  & fac000, fac100, fac200, fac010, fac110, fac210, &
5169  & fac001, fac101, fac201, fac011, fac111, fac211, &
5170  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21, temp
5171 !
5172 !===> ... begin here
5173 !
5174 ! --- ... minor gas mapping level :
5175 ! lower - n2o, p = 706.272 mbar, t = 278.94 k
5176 ! upper - n2o, p = 95.58 mbar, t = 215.7 k
5177 
5178 ! --- ... calculate reference ratio to be used in calculation of Planck
5179 ! fraction in lower/upper atmosphere.
5180 
5181  refrat_planck_a = chi_mls(1,9)/chi_mls(6,9) ! P = 212 mb
5182  refrat_m_a = chi_mls(1,3)/chi_mls(6,3) ! P = 706.272 mb
5183 
5184 ! --- ... lower atmosphere loop
5185 
5186  do k = 1, laytrop
5187  speccomb = colamt(k,1) + rfrate(k,4,1)*colamt(k,5)
5188  specparm = colamt(k,1) / speccomb
5189  specmult = 8.0 * min(specparm, oneminus)
5190  js = 1 + int(specmult)
5191  fs = mod(specmult, f_one)
5192  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(9) + js
5193 
5194  speccomb1 = colamt(k,1) + rfrate(k,4,2)*colamt(k,5)
5195  specparm1 = colamt(k,1) / speccomb1
5196  specmult1 = 8.0 * min(specparm1, oneminus)
5197  js1 = 1 + int(specmult1)
5198  fs1 = mod(specmult1, f_one)
5199  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(9) + js1
5200 
5201  speccomb_mn2o = colamt(k,1) + refrat_m_a*colamt(k,5)
5202  specparm_mn2o = colamt(k,1) / speccomb_mn2o
5203  specmult_mn2o = 8.0 * min(specparm_mn2o, oneminus)
5204  jmn2o = 1 + int(specmult_mn2o)
5205  fmn2o = mod(specmult_mn2o, f_one)
5206 
5207  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,5)
5208  specparm_planck = colamt(k,1) / speccomb_planck
5209  specmult_planck = 8.0 * min(specparm_planck, oneminus)
5210  jpl = 1 + int(specmult_planck)
5211  fpl = mod(specmult_planck, f_one)
5212 
5213  inds = indself(k)
5214  indf = indfor(k)
5215  indm = indminor(k)
5216  indsp = inds + 1
5217  indfp = indf + 1
5218  indmp = indm + 1
5219  jplp = jpl + 1
5220  jmn2op= jmn2o+ 1
5221 
5222 ! --- ... in atmospheres where the amount of n2o is too great to be considered
5223 ! a minor species, adjust the column amount of n2o by an empirical factor
5224 ! to obtain the proper contribution.
5225 
5226  temp = coldry(k) * chi_mls(4,jp(k)+1)
5227  ratn2o = colamt(k,4) / temp
5228  if (ratn2o > 1.5) then
5229  adjfac = 0.5 + (ratn2o-0.5)**0.65
5230  adjcoln2o = adjfac * temp
5231  else
5232  adjcoln2o = colamt(k,4)
5233  endif
5234 
5235  if (specparm < 0.125 .and. specparm1 < 0.125) then
5236  p0 = fs - f_one
5237  p40 = p0**4
5238  fk00 = p40
5239  fk10 = f_one - p0 - 2.0*p40
5240  fk20 = p0 + p40
5241 
5242  p1 = fs1 - f_one
5243  p41 = p1**4
5244  fk01 = p41
5245  fk11 = f_one - p1 - 2.0*p41
5246  fk21 = p1 + p41
5247 
5248  id000 = ind0
5249  id010 = ind0 + 9
5250  id100 = ind0 + 1
5251  id110 = ind0 +10
5252  id200 = ind0 + 2
5253  id210 = ind0 +11
5254 
5255  id001 = ind1
5256  id011 = ind1 + 9
5257  id101 = ind1 + 1
5258  id111 = ind1 +10
5259  id201 = ind1 + 2
5260  id211 = ind1 +11
5261 
5262  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
5263  p0 = -fs
5264  p40 = p0**4
5265  fk00 = p40
5266  fk10 = f_one - p0 - 2.0*p40
5267  fk20 = p0 + p40
5268 
5269  p1 = -fs1
5270  p41 = p1**4
5271  fk01 = p41
5272  fk11 = f_one - p1 - 2.0*p41
5273  fk21 = p1 + p41
5274 
5275  id000 = ind0 + 1
5276  id010 = ind0 +10
5277  id100 = ind0
5278  id110 = ind0 + 9
5279  id200 = ind0 - 1
5280  id210 = ind0 + 8
5281 
5282  id001 = ind1 + 1
5283  id011 = ind1 +10
5284  id101 = ind1
5285  id111 = ind1 + 9
5286  id201 = ind1 - 1
5287  id211 = ind1 + 8
5288  else
5289  fk00 = f_one - fs
5290  fk10 = fs
5291  fk20 = f_zero
5292 
5293  fk01 = f_one - fs1
5294  fk11 = fs1
5295  fk21 = f_zero
5296 
5297  id000 = ind0
5298  id010 = ind0 + 9
5299  id100 = ind0 + 1
5300  id110 = ind0 +10
5301  id200 = ind0
5302  id210 = ind0
5303 
5304  id001 = ind1
5305  id011 = ind1 + 9
5306  id101 = ind1 + 1
5307  id111 = ind1 +10
5308  id201 = ind1
5309  id211 = ind1
5310  endif
5311 
5312  fac000 = fk00 * fac00(k)
5313  fac100 = fk10 * fac00(k)
5314  fac200 = fk20 * fac00(k)
5315  fac010 = fk00 * fac10(k)
5316  fac110 = fk10 * fac10(k)
5317  fac210 = fk20 * fac10(k)
5318 
5319  fac001 = fk01 * fac01(k)
5320  fac101 = fk11 * fac01(k)
5321  fac201 = fk21 * fac01(k)
5322  fac011 = fk01 * fac11(k)
5323  fac111 = fk11 * fac11(k)
5324  fac211 = fk21 * fac11(k)
5325 
5326  do ig = 1, ng09
5327  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
5328  & * (selfref(ig,indsp) - selfref(ig,inds)))
5329  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5330  & * (forref(ig,indfp) - forref(ig,indf)))
5331  n2om1 = ka_mn2o(ig,jmn2o,indm) + fmn2o &
5332  & * (ka_mn2o(ig,jmn2op,indm) - ka_mn2o(ig,jmn2o,indm))
5333  n2om2 = ka_mn2o(ig,jmn2o,indmp) + fmn2o &
5334  & * (ka_mn2o(ig,jmn2op,indmp) - ka_mn2o(ig,jmn2o,indmp))
5335  absn2o = n2om1 + minorfrac(k) * (n2om2 - n2om1)
5336 
5337  taug(ns09+ig,k) = speccomb &
5338  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
5339  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
5340  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
5341  & + speccomb1 &
5342  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
5343  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
5344  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
5345  & + tauself + taufor + adjcoln2o*absn2o
5346 
5347  fracs(ns09+ig,k) = fracrefa(ig,jpl) + fpl &
5348  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
5349  enddo
5350  enddo
5351 
5352 ! --- ... upper atmosphere loop
5353 
5354  do k = laytrop+1, nlay
5355  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(9) + 1
5356  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(9) + 1
5357 
5358  indm = indminor(k)
5359  ind0p = ind0 + 1
5360  ind1p = ind1 + 1
5361  indmp = indm + 1
5362 
5363 ! --- ... in atmospheres where the amount of n2o is too great to be considered
5364 ! a minor species, adjust the column amount of n2o by an empirical factor
5365 ! to obtain the proper contribution.
5366 
5367  temp = coldry(k) * chi_mls(4,jp(k)+1)
5368  ratn2o = colamt(k,4) / temp
5369  if (ratn2o > 1.5) then
5370  adjfac = 0.5 + (ratn2o - 0.5)**0.65
5371  adjcoln2o = adjfac * temp
5372  else
5373  adjcoln2o = colamt(k,4)
5374  endif
5375 
5376  do ig = 1, ng09
5377  absn2o = kb_mn2o(ig,indm) + minorfrac(k) &
5378  & * (kb_mn2o(ig,indmp) - kb_mn2o(ig,indm))
5379 
5380  taug(ns09+ig,k) = colamt(k,5) &
5381  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
5382  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
5383  & + adjcoln2o*absn2o
5384 
5385  fracs(ns09+ig,k) = fracrefb(ig)
5386  enddo
5387  enddo
5388 
5389 ! ..................................
5390  end subroutine taugb09
5391 ! ----------------------------------
5392 
5393 ! ----------------------------------
5394  subroutine taugb10
5395 ! ..................................
5396 
5397 ! ------------------------------------------------------------------ !
5398 ! band 10: 1390-1480 cm-1 (low key - h2o; high key - h2o) !
5399 ! ------------------------------------------------------------------ !
5400 
5401  use module_radlw_kgb10
5402 
5403 ! --- locals:
5404  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
5405  & ig
5406 
5407  real (kind=kind_phys) :: tauself, taufor
5408 !
5409 !===> ... begin here
5410 !
5411 ! --- ... lower atmosphere loop
5412 
5413  do k = 1, laytrop
5414  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(10) + 1
5415  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(10) + 1
5416 
5417  inds = indself(k)
5418  indf = indfor(k)
5419  ind0p = ind0 + 1
5420  ind1p = ind1 + 1
5421  indsp = inds + 1
5422  indfp = indf + 1
5423 
5424  do ig = 1, ng10
5425  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
5426  & * (selfref(ig,indsp) - selfref(ig,inds)))
5427  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5428  & * (forref(ig,indfp) - forref(ig,indf)))
5429 
5430  taug(ns10+ig,k) = colamt(k,1) &
5431  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
5432  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
5433  & + tauself + taufor
5434 
5435  fracs(ns10+ig,k) = fracrefa(ig)
5436  enddo
5437  enddo
5438 
5439 ! --- ... upper atmosphere loop
5440 
5441  do k = laytrop+1, nlay
5442  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(10) + 1
5443  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(10) + 1
5444 
5445  indf = indfor(k)
5446  ind0p = ind0 + 1
5447  ind1p = ind1 + 1
5448  indfp = indf + 1
5449 
5450  do ig = 1, ng10
5451  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5452  & * (forref(ig,indfp) - forref(ig,indf)))
5453 
5454  taug(ns10+ig,k) = colamt(k,1) &
5455  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
5456  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
5457  & + taufor
5458 
5459  fracs(ns10+ig,k) = fracrefb(ig)
5460  enddo
5461  enddo
5462 
5463 ! ..................................
5464  end subroutine taugb10
5465 ! ----------------------------------
5466 
5467 ! ----------------------------------
5468  subroutine taugb11
5469 ! ..................................
5470 
5471 ! ------------------------------------------------------------------ !
5472 ! band 11: 1480-1800 cm-1 (low - h2o; low minor - o2) !
5473 ! (high key - h2o; high minor - o2) !
5474 ! ------------------------------------------------------------------ !
5475 
5476  use module_radlw_kgb11
5477 
5478 ! --- locals:
5479  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
5480  & indm, indmp, ig
5481 
5482  real (kind=kind_phys) :: scaleo2, tauself, taufor, tauo2
5483 !
5484 !===> ... begin here
5485 !
5486 ! --- ... minor gas mapping level :
5487 ! lower - o2, p = 706.2720 mbar, t = 278.94 k
5488 ! upper - o2, p = 4.758820 mbarm t = 250.85 k
5489 
5490 ! --- ... lower atmosphere loop
5491 
5492  do k = 1, laytrop
5493  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(11) + 1
5494  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(11) + 1
5495 
5496  inds = indself(k)
5497  indf = indfor(k)
5498  indm = indminor(k)
5499  ind0p = ind0 + 1
5500  ind1p = ind1 + 1
5501  indsp = inds + 1
5502  indfp = indf + 1
5503  indmp = indm + 1
5504 
5505  scaleo2 = colamt(k,6) * scaleminor(k)
5506 
5507  do ig = 1, ng11
5508  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
5509  & * (selfref(ig,indsp) - selfref(ig,inds)))
5510  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5511  & * (forref(ig,indfp) - forref(ig,indf)))
5512  tauo2 = scaleo2 * (ka_mo2(ig,indm) + minorfrac(k) &
5513  & * (ka_mo2(ig,indmp) - ka_mo2(ig,indm)))
5514 
5515  taug(ns11+ig,k) = colamt(k,1) &
5516  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
5517  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
5518  & + tauself + taufor + tauo2
5519 
5520  fracs(ns11+ig,k) = fracrefa(ig)
5521  enddo
5522  enddo
5523 
5524 ! --- ... upper atmosphere loop
5525 
5526  do k = laytrop+1, nlay
5527  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(11) + 1
5528  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(11) + 1
5529 
5530  indf = indfor(k)
5531  indm = indminor(k)
5532  ind0p = ind0 + 1
5533  ind1p = ind1 + 1
5534  indfp = indf + 1
5535  indmp = indm + 1
5536 
5537  scaleo2 = colamt(k,6) * scaleminor(k)
5538 
5539  do ig = 1, ng11
5540  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5541  & * (forref(ig,indfp) - forref(ig,indf)))
5542  tauo2 = scaleo2 * (kb_mo2(ig,indm) + minorfrac(k) &
5543  & * (kb_mo2(ig,indmp) - kb_mo2(ig,indm)))
5544 
5545  taug(ns11+ig,k) = colamt(k,1) &
5546  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
5547  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p)) &
5548  & + taufor + tauo2
5549 
5550  fracs(ns11+ig,k) = fracrefb(ig)
5551  enddo
5552  enddo
5553 
5554 ! ..................................
5555  end subroutine taugb11
5556 ! ----------------------------------
5557 
5558 ! ----------------------------------
5559  subroutine taugb12
5560 ! ..................................
5561 
5562 ! ------------------------------------------------------------------ !
5563 ! band 12: 1800-2080 cm-1 (low - h2o,co2; high - nothing) !
5564 ! ------------------------------------------------------------------ !
5565 
5566  use module_radlw_kgb12
5567 
5568 ! --- locals:
5569  integer :: k, ind0, ind1, inds, indsp, indf, indfp, jpl, jplp, &
5570  & id000, id010, id100, id110, id200, id210, ig, js, js1, &
5571  & id001, id011, id101, id111, id201, id211
5572 
5573  real (kind=kind_phys) :: tauself, taufor, refrat_planck_a, &
5574  & speccomb, specparm, specmult, fs, &
5575  & speccomb1, specparm1, specmult1, fs1, &
5576  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
5577  & fac000, fac100, fac200, fac010, fac110, fac210, &
5578  & fac001, fac101, fac201, fac011, fac111, fac211, &
5579  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21
5580 !
5581 !===> ... begin here
5582 !
5583 ! --- ... calculate reference ratio to be used in calculation of Planck
5584 ! fraction in lower/upper atmosphere.
5585 
5586  refrat_planck_a = chi_mls(1,10)/chi_mls(2,10) ! P = 174.164 mb
5587 
5588 ! --- ... lower atmosphere loop
5589 
5590  do k = 1, laytrop
5591  speccomb = colamt(k,1) + rfrate(k,1,1)*colamt(k,2)
5592  specparm = colamt(k,1) / speccomb
5593  specmult = 8.0 * min(specparm, oneminus)
5594  js = 1 + int(specmult)
5595  fs = mod(specmult, f_one)
5596  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(12) + js
5597 
5598  speccomb1 = colamt(k,1) + rfrate(k,1,2)*colamt(k,2)
5599  specparm1 = colamt(k,1) / speccomb1
5600  specmult1 = 8.0 * min(specparm1, oneminus)
5601  js1 = 1 + int(specmult1)
5602  fs1 = mod(specmult1, f_one)
5603  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(12) + js1
5604 
5605  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,2)
5606  specparm_planck = colamt(k,1) / speccomb_planck
5607  if (specparm_planck >= oneminus) specparm_planck=oneminus
5608  specmult_planck = 8.0 * specparm_planck
5609  jpl = 1 + int(specmult_planck)
5610  fpl = mod(specmult_planck, f_one)
5611 
5612  inds = indself(k)
5613  indf = indfor(k)
5614  indsp = inds + 1
5615  indfp = indf + 1
5616  jplp = jpl + 1
5617 
5618  if (specparm < 0.125 .and. specparm1 < 0.125) then
5619  p0 = fs - f_one
5620  p40 = p0**4
5621  fk00 = p40
5622  fk10 = f_one - p0 - 2.0*p40
5623  fk20 = p0 + p40
5624 
5625  p1 = fs1 - f_one
5626  p41 = p1**4
5627  fk01 = p41
5628  fk11 = f_one - p1 - 2.0*p41
5629  fk21 = p1 + p41
5630 
5631  id000 = ind0
5632  id010 = ind0 + 9
5633  id100 = ind0 + 1
5634  id110 = ind0 +10
5635  id200 = ind0 + 2
5636  id210 = ind0 +11
5637 
5638  id001 = ind1
5639  id011 = ind1 + 9
5640  id101 = ind1 + 1
5641  id111 = ind1 +10
5642  id201 = ind1 + 2
5643  id211 = ind1 +11
5644  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
5645  p0 = -fs
5646  p40 = p0**4
5647  fk00 = p40
5648  fk10 = f_one - p0 - 2.0*p40
5649  fk20 = p0 + p40
5650 
5651  p1 = -fs1
5652  p41 = p1**4
5653  fk01 = p41
5654  fk11 = f_one - p1 - 2.0*p41
5655  fk21 = p1 + p41
5656 
5657  id000 = ind0 + 1
5658  id010 = ind0 +10
5659  id100 = ind0
5660  id110 = ind0 + 9
5661  id200 = ind0 - 1
5662  id210 = ind0 + 8
5663 
5664  id001 = ind1 + 1
5665  id011 = ind1 +10
5666  id101 = ind1
5667  id111 = ind1 + 9
5668  id201 = ind1 - 1
5669  id211 = ind1 + 8
5670  else
5671  fk00 = f_one - fs
5672  fk10 = fs
5673  fk20 = f_zero
5674 
5675  fk01 = f_one - fs1
5676  fk11 = fs1
5677  fk21 = f_zero
5678 
5679  id000 = ind0
5680  id010 = ind0 + 9
5681  id100 = ind0 + 1
5682  id110 = ind0 +10
5683  id200 = ind0
5684  id210 = ind0
5685 
5686  id001 = ind1
5687  id011 = ind1 + 9
5688  id101 = ind1 + 1
5689  id111 = ind1 +10
5690  id201 = ind1
5691  id211 = ind1
5692  endif
5693 
5694  fac000 = fk00 * fac00(k)
5695  fac100 = fk10 * fac00(k)
5696  fac200 = fk20 * fac00(k)
5697  fac010 = fk00 * fac10(k)
5698  fac110 = fk10 * fac10(k)
5699  fac210 = fk20 * fac10(k)
5700 
5701  fac001 = fk01 * fac01(k)
5702  fac101 = fk11 * fac01(k)
5703  fac201 = fk21 * fac01(k)
5704  fac011 = fk01 * fac11(k)
5705  fac111 = fk11 * fac11(k)
5706  fac211 = fk21 * fac11(k)
5707 
5708  do ig = 1, ng12
5709  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
5710  & * (selfref(ig,indsp) - selfref(ig,inds)))
5711  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5712  & * (forref(ig,indfp) - forref(ig,indf)))
5713 
5714  taug(ns12+ig,k) = speccomb &
5715  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
5716  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
5717  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
5718  & + speccomb1 &
5719  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
5720  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
5721  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
5722  & + tauself + taufor
5723 
5724  fracs(ns12+ig,k) = fracrefa(ig,jpl) + fpl &
5725  & *(fracrefa(ig,jplp) - fracrefa(ig,jpl))
5726  enddo
5727  enddo
5728 
5729 ! --- ... upper atmosphere loop
5730 
5731  do k = laytrop+1, nlay
5732  do ig = 1, ng12
5733  taug(ns12+ig,k) = f_zero
5734  fracs(ns12+ig,k) = f_zero
5735  enddo
5736  enddo
5737 
5738 ! ..................................
5739  end subroutine taugb12
5740 ! ----------------------------------
5741 
5742 ! ----------------------------------
5743  subroutine taugb13
5744 ! ..................................
5745 
5746 ! ------------------------------------------------------------------ !
5747 ! band 13: 2080-2250 cm-1 (low key-h2o,n2o; high minor-o3 minor) !
5748 ! ------------------------------------------------------------------ !
5749 
5750  use module_radlw_kgb13
5751 
5752 ! --- locals:
5753  integer :: k, ind0, ind1, inds, indsp, indf, indfp, indm, indmp, &
5754  & id000, id010, id100, id110, id200, id210, jmco2, jpl, &
5755  & id001, id011, id101, id111, id201, id211, jmco2p, jplp, &
5756  & jmco, jmcop, ig, js, js1
5757 
5758  real (kind=kind_phys) :: tauself, taufor, co2m1, co2m2, absco2, &
5759  & speccomb, specparm, specmult, fs, &
5760  & speccomb1, specparm1, specmult1, fs1, &
5761  & speccomb_mco2, specparm_mco2, specmult_mco2, fmco2, &
5762  & speccomb_mco, specparm_mco, specmult_mco, fmco, &
5763  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
5764  & refrat_planck_a, refrat_m_a, refrat_m_a3, ratco2, &
5765  & adjfac, adjcolco2, com1, com2, absco, abso3, &
5766  & fac000, fac100, fac200, fac010, fac110, fac210, &
5767  & fac001, fac101, fac201, fac011, fac111, fac211, &
5768  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21, temp
5769 !
5770 !===> ... begin here
5771 !
5772 ! --- ... minor gas mapping levels :
5773 ! lower - co2, p = 1053.63 mb, t = 294.2 k
5774 ! lower - co, p = 706 mb, t = 278.94 k
5775 ! upper - o3, p = 95.5835 mb, t = 215.7 k
5776 
5777 ! --- ... calculate reference ratio to be used in calculation of Planck
5778 ! fraction in lower/upper atmosphere.
5779 
5780  refrat_planck_a = chi_mls(1,5)/chi_mls(4,5) ! P = 473.420 mb (Level 5)
5781  refrat_m_a = chi_mls(1,1)/chi_mls(4,1) ! P = 1053. (Level 1)
5782  refrat_m_a3 = chi_mls(1,3)/chi_mls(4,3) ! P = 706. (Level 3)
5783 
5784 ! --- ... lower atmosphere loop
5785 
5786  do k = 1, laytrop
5787  speccomb = colamt(k,1) + rfrate(k,3,1)*colamt(k,4)
5788  specparm = colamt(k,1) / speccomb
5789  specmult = 8.0 * min(specparm, oneminus)
5790  js = 1 + int(specmult)
5791  fs = mod(specmult, f_one)
5792  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(13) + js
5793 
5794  speccomb1 = colamt(k,1) + rfrate(k,3,2)*colamt(k,4)
5795  specparm1 = colamt(k,1) / speccomb1
5796  specmult1 = 8.0 * min(specparm1, oneminus)
5797  js1 = 1 + int(specmult1)
5798  fs1 = mod(specmult1, f_one)
5799  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(13) + js1
5800 
5801  speccomb_mco2 = colamt(k,1) + refrat_m_a*colamt(k,4)
5802  specparm_mco2 = colamt(k,1) / speccomb_mco2
5803  specmult_mco2 = 8.0 * min(specparm_mco2, oneminus)
5804  jmco2 = 1 + int(specmult_mco2)
5805  fmco2 = mod(specmult_mco2, f_one)
5806 
5807 ! --- ... in atmospheres where the amount of co2 is too great to be considered
5808 ! a minor species, adjust the column amount of co2 by an empirical factor
5809 ! to obtain the proper contribution.
5810 
5811  speccomb_mco = colamt(k,1) + refrat_m_a3*colamt(k,4)
5812  specparm_mco = colamt(k,1) / speccomb_mco
5813  specmult_mco = 8.0 * min(specparm_mco, oneminus)
5814  jmco = 1 + int(specmult_mco)
5815  fmco = mod(specmult_mco, f_one)
5816 
5817  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,4)
5818  specparm_planck = colamt(k,1) / speccomb_planck
5819  specmult_planck = 8.0 * min(specparm_planck, oneminus)
5820  jpl = 1 + int(specmult_planck)
5821  fpl = mod(specmult_planck, f_one)
5822 
5823  inds = indself(k)
5824  indf = indfor(k)
5825  indm = indminor(k)
5826  indsp = inds + 1
5827  indfp = indf + 1
5828  indmp = indm + 1
5829  jplp = jpl + 1
5830  jmco2p= jmco2+ 1
5831  jmcop = jmco + 1
5832 
5833 ! --- ... in atmospheres where the amount of co2 is too great to be considered
5834 ! a minor species, adjust the column amount of co2 by an empirical factor
5835 ! to obtain the proper contribution.
5836 
5837  temp = coldry(k) * 3.55e-4
5838  ratco2 = colamt(k,2) / temp
5839  if (ratco2 > 3.0) then
5840  adjfac = 2.0 + (ratco2-2.0)**0.68
5841  adjcolco2 = adjfac * temp
5842  else
5843  adjcolco2 = colamt(k,2)
5844  endif
5845 
5846  if (specparm < 0.125 .and. specparm1 < 0.125) then
5847  p0 = fs - f_one
5848  p40 = p0**4
5849  fk00 = p40
5850  fk10 = f_one - p0 - 2.0*p40
5851  fk20 = p0 + p40
5852 
5853  p1 = fs1 - f_one
5854  p41 = p1**4
5855  fk01 = p41
5856  fk11 = f_one - p1 - 2.0*p41
5857  fk21 = p1 + p41
5858 
5859  id000 = ind0
5860  id010 = ind0 + 9
5861  id100 = ind0 + 1
5862  id110 = ind0 +10
5863  id200 = ind0 + 2
5864  id210 = ind0 +11
5865 
5866  id001 = ind1
5867  id011 = ind1 + 9
5868  id101 = ind1 + 1
5869  id111 = ind1 +10
5870  id201 = ind1 + 2
5871  id211 = ind1 +11
5872  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
5873  p0 = -fs
5874  p40 = p0**4
5875  fk00 = p40
5876  fk10 = f_one - p0 - 2.0*p40
5877  fk20 = p0 + p40
5878 
5879  p1 = -fs1
5880  p41 = p1**4
5881  fk01 = p41
5882  fk11 = f_one - p1 - 2.0*p41
5883  fk21 = p1 + p41
5884 
5885  id000 = ind0 + 1
5886  id010 = ind0 +10
5887  id100 = ind0
5888  id110 = ind0 + 9
5889  id200 = ind0 - 1
5890  id210 = ind0 + 8
5891 
5892  id001 = ind1 + 1
5893  id011 = ind1 +10
5894  id101 = ind1
5895  id111 = ind1 + 9
5896  id201 = ind1 - 1
5897  id211 = ind1 + 8
5898  else
5899  fk00 = f_one - fs
5900  fk10 = fs
5901  fk20 = f_zero
5902 
5903  fk01 = f_one - fs1
5904  fk11 = fs1
5905  fk21 = f_zero
5906 
5907  id000 = ind0
5908  id010 = ind0 + 9
5909  id100 = ind0 + 1
5910  id110 = ind0 +10
5911  id200 = ind0
5912  id210 = ind0
5913 
5914  id001 = ind1
5915  id011 = ind1 + 9
5916  id101 = ind1 + 1
5917  id111 = ind1 +10
5918  id201 = ind1
5919  id211 = ind1
5920  endif
5921 
5922  fac000 = fk00 * fac00(k)
5923  fac100 = fk10 * fac00(k)
5924  fac200 = fk20 * fac00(k)
5925  fac010 = fk00 * fac10(k)
5926  fac110 = fk10 * fac10(k)
5927  fac210 = fk20 * fac10(k)
5928 
5929  fac001 = fk01 * fac01(k)
5930  fac101 = fk11 * fac01(k)
5931  fac201 = fk21 * fac01(k)
5932  fac011 = fk01 * fac11(k)
5933  fac111 = fk11 * fac11(k)
5934  fac211 = fk21 * fac11(k)
5935 
5936  do ig = 1, ng13
5937  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
5938  & * (selfref(ig,indsp) - selfref(ig,inds)))
5939  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
5940  & * (forref(ig,indfp) - forref(ig,indf)))
5941  co2m1 = ka_mco2(ig,jmco2,indm) + fmco2 &
5942  & * (ka_mco2(ig,jmco2p,indm) - ka_mco2(ig,jmco2,indm))
5943  co2m2 = ka_mco2(ig,jmco2,indmp) + fmco2 &
5944  & * (ka_mco2(ig,jmco2p,indmp) - ka_mco2(ig,jmco2,indmp))
5945  absco2 = co2m1 + minorfrac(k) * (co2m2 - co2m1)
5946  com1 = ka_mco(ig,jmco,indm) + fmco &
5947  & * (ka_mco(ig,jmcop,indm) - ka_mco(ig,jmco,indm))
5948  com2 = ka_mco(ig,jmco,indmp) + fmco &
5949  & * (ka_mco(ig,jmcop,indmp) - ka_mco(ig,jmco,indmp))
5950  absco = com1 + minorfrac(k) * (com2 - com1)
5951 
5952  taug(ns13+ig,k) = speccomb &
5953  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
5954  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
5955  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
5956  & + speccomb1 &
5957  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
5958  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
5959  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
5960  & + tauself + taufor + adjcolco2*absco2 &
5961  & + colamt(k,7)*absco
5962 
5963  fracs(ns13+ig,k) = fracrefa(ig,jpl) + fpl &
5964  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
5965  enddo
5966  enddo
5967 
5968 ! --- ... upper atmosphere loop
5969 
5970  do k = laytrop+1, nlay
5971  indm = indminor(k)
5972  indmp = indm + 1
5973 
5974  do ig = 1, ng13
5975  abso3 = kb_mo3(ig,indm) + minorfrac(k) &
5976  & * (kb_mo3(ig,indmp) - kb_mo3(ig,indm))
5977 
5978  taug(ns13+ig,k) = colamt(k,3)*abso3
5979 
5980  fracs(ns13+ig,k) = fracrefb(ig)
5981  enddo
5982  enddo
5983 
5984 ! ..................................
5985  end subroutine taugb13
5986 ! ----------------------------------
5987 
5988 ! ----------------------------------
5989  subroutine taugb14
5990 ! ..................................
5991 
5992 ! ------------------------------------------------------------------ !
5993 ! band 14: 2250-2380 cm-1 (low - co2; high - co2) !
5994 ! ------------------------------------------------------------------ !
5995 
5996  use module_radlw_kgb14
5997 
5998 ! --- locals:
5999  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
6000  & ig
6001 
6002  real (kind=kind_phys) :: tauself, taufor
6003 !
6004 !===> ... begin here
6005 !
6006 ! --- ... lower atmosphere loop
6007 
6008  do k = 1, laytrop
6009  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(14) + 1
6010  ind1 = ( jp(k) *5 + (jt1(k)-1)) * nspa(14) + 1
6011 
6012  inds = indself(k)
6013  indf = indfor(k)
6014  ind0p = ind0 + 1
6015  ind1p = ind1 + 1
6016  indsp = inds + 1
6017  indfp = indf + 1
6018 
6019  do ig = 1, ng14
6020  tauself = selffac(k) * (selfref(ig,inds) + selffrac(k) &
6021  & * (selfref(ig,indsp) - selfref(ig,inds)))
6022  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
6023  & * (forref(ig,indfp) - forref(ig,indf)))
6024 
6025  taug(ns14+ig,k) = colamt(k,2) &
6026  & * (fac00(k)*absa(ig,ind0) + fac10(k)*absa(ig,ind0p) &
6027  & + fac01(k)*absa(ig,ind1) + fac11(k)*absa(ig,ind1p)) &
6028  & + tauself + taufor
6029 
6030  fracs(ns14+ig,k) = fracrefa(ig)
6031  enddo
6032  enddo
6033 
6034 ! --- ... upper atmosphere loop
6035 
6036  do k = laytrop+1, nlay
6037  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(14) + 1
6038  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(14) + 1
6039 
6040  ind0p = ind0 + 1
6041  ind1p = ind1 + 1
6042 
6043  do ig = 1, ng14
6044  taug(ns14+ig,k) = colamt(k,2) &
6045  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
6046  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p))
6047 
6048  fracs(ns14+ig,k) = fracrefb(ig)
6049  enddo
6050  enddo
6051 
6052 ! ..................................
6053  end subroutine taugb14
6054 ! ----------------------------------
6055 
6056 ! ----------------------------------
6057  subroutine taugb15
6058 ! ..................................
6059 
6060 ! ------------------------------------------------------------------ !
6061 ! band 15: 2380-2600 cm-1 (low - n2o,co2; low minor - n2) !
6062 ! (high - nothing) !
6063 ! ------------------------------------------------------------------ !
6064 
6065  use module_radlw_kgb15
6066 
6067 ! --- locals:
6068  integer :: k, ind0, ind1, inds, indsp, indf, indfp, indm, indmp, &
6069  & id000, id010, id100, id110, id200, id210, jpl, jplp, &
6070  & id001, id011, id101, id111, id201, id211, jmn2, jmn2p, &
6071  & ig, js, js1
6072 
6073  real (kind=kind_phys) :: scalen2, tauself, taufor, &
6074  & speccomb, specparm, specmult, fs, &
6075  & speccomb1, specparm1, specmult1, fs1, &
6076  & speccomb_mn2, specparm_mn2, specmult_mn2, fmn2, &
6077  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
6078  & refrat_planck_a, refrat_m_a, n2m1, n2m2, taun2, &
6079  & fac000, fac100, fac200, fac010, fac110, fac210, &
6080  & fac001, fac101, fac201, fac011, fac111, fac211, &
6081  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21
6082 !
6083 !===> ... begin here
6084 !
6085 ! --- ... minor gas mapping level :
6086 ! lower - nitrogen continuum, P = 1053., T = 294.
6087 
6088 ! --- ... calculate reference ratio to be used in calculation of Planck
6089 ! fraction in lower atmosphere.
6090 
6091  refrat_planck_a = chi_mls(4,1)/chi_mls(2,1) ! P = 1053. mb (Level 1)
6092  refrat_m_a = chi_mls(4,1)/chi_mls(2,1) ! P = 1053. mb
6093 
6094 ! --- ... lower atmosphere loop
6095 
6096  do k = 1, laytrop
6097  speccomb = colamt(k,4) + rfrate(k,5,1)*colamt(k,2)
6098  specparm = colamt(k,4) / speccomb
6099  specmult = 8.0 * min(specparm, oneminus)
6100  js = 1 + int(specmult)
6101  fs = mod(specmult, f_one)
6102  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(15) + js
6103 
6104  speccomb1 = colamt(k,4) + rfrate(k,5,2)*colamt(k,2)
6105  specparm1 = colamt(k,4) / speccomb1
6106  specmult1 = 8.0 * min(specparm1, oneminus)
6107  js1 = 1 + int(specmult1)
6108  fs1 = mod(specmult1, f_one)
6109  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(15) + js1
6110 
6111  speccomb_mn2 = colamt(k,4) + refrat_m_a*colamt(k,2)
6112  specparm_mn2 = colamt(k,4) / speccomb_mn2
6113  specmult_mn2 = 8.0 * min(specparm_mn2, oneminus)
6114  jmn2 = 1 + int(specmult_mn2)
6115  fmn2 = mod(specmult_mn2, f_one)
6116 
6117  speccomb_planck = colamt(k,4) + refrat_planck_a*colamt(k,2)
6118  specparm_planck = colamt(k,4) / speccomb_planck
6119  specmult_planck = 8.0 * min(specparm_planck, oneminus)
6120  jpl = 1 + int(specmult_planck)
6121  fpl = mod(specmult_planck, f_one)
6122 
6123  scalen2 = colbrd(k) * scaleminor(k)
6124 
6125  inds = indself(k)
6126  indf = indfor(k)
6127  indm = indminor(k)
6128  indsp = inds + 1
6129  indfp = indf + 1
6130  indmp = indm + 1
6131  jplp = jpl + 1
6132  jmn2p = jmn2 + 1
6133 
6134  if (specparm < 0.125 .and. specparm1 < 0.125) then
6135  p0 = fs - f_one
6136  p40 = p0**4
6137  fk00 = p40
6138  fk10 = f_one - p0 - 2.0*p40
6139  fk20 = p0 + p40
6140 
6141  p1 = fs1 - f_one
6142  p41 = p1**4
6143  fk01 = p41
6144  fk11 = f_one - p1 - 2.0*p41
6145  fk21 = p1 + p41
6146 
6147  id000 = ind0
6148  id010 = ind0 + 9
6149  id100 = ind0 + 1
6150  id110 = ind0 +10
6151  id200 = ind0 + 2
6152  id210 = ind0 +11
6153 
6154  id001 = ind1
6155  id011 = ind1 + 9
6156  id101 = ind1 + 1
6157  id111 = ind1 +10
6158  id201 = ind1 + 2
6159  id211 = ind1 +11
6160  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
6161  p0 = -fs
6162  p40 = p0**4
6163  fk00 = p40
6164  fk10 = f_one - p0 - 2.0*p40
6165  fk20 = p0 + p40
6166 
6167  p1 = -fs1
6168  p41 = p1**4
6169  fk01 = p41
6170  fk11 = f_one - p1 - 2.0*p41
6171  fk21 = p1 + p41
6172 
6173  id000 = ind0 + 1
6174  id010 = ind0 +10
6175  id100 = ind0
6176  id110 = ind0 + 9
6177  id200 = ind0 - 1
6178  id210 = ind0 + 8
6179 
6180  id001 = ind1 + 1
6181  id011 = ind1 +10
6182  id101 = ind1
6183  id111 = ind1 + 9
6184  id201 = ind1 - 1
6185  id211 = ind1 + 8
6186  else
6187  fk00 = f_one - fs
6188  fk10 = fs
6189  fk20 = f_zero
6190 
6191  fk01 = f_one - fs1
6192  fk11 = fs1
6193  fk21 = f_zero
6194 
6195  id000 = ind0
6196  id010 = ind0 + 9
6197  id100 = ind0 + 1
6198  id110 = ind0 +10
6199  id200 = ind0
6200  id210 = ind0
6201 
6202  id001 = ind1
6203  id011 = ind1 + 9
6204  id101 = ind1 + 1
6205  id111 = ind1 +10
6206  id201 = ind1
6207  id211 = ind1
6208  endif
6209 
6210  fac000 = fk00 * fac00(k)
6211  fac100 = fk10 * fac00(k)
6212  fac200 = fk20 * fac00(k)
6213  fac010 = fk00 * fac10(k)
6214  fac110 = fk10 * fac10(k)
6215  fac210 = fk20 * fac10(k)
6216 
6217  fac001 = fk01 * fac01(k)
6218  fac101 = fk11 * fac01(k)
6219  fac201 = fk21 * fac01(k)
6220  fac011 = fk01 * fac11(k)
6221  fac111 = fk11 * fac11(k)
6222  fac211 = fk21 * fac11(k)
6223 
6224  do ig = 1, ng15
6225  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
6226  & * (selfref(ig,indsp) - selfref(ig,inds)))
6227  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
6228  & * (forref(ig,indfp) - forref(ig,indf)))
6229  n2m1 = ka_mn2(ig,jmn2,indm) + fmn2 &
6230  & * (ka_mn2(ig,jmn2p,indm) - ka_mn2(ig,jmn2,indm))
6231  n2m2 = ka_mn2(ig,jmn2,indmp) + fmn2 &
6232  & * (ka_mn2(ig,jmn2p,indmp) - ka_mn2(ig,jmn2,indmp))
6233  taun2 = scalen2 * (n2m1 + minorfrac(k) * (n2m2 - n2m1))
6234 
6235  taug(ns15+ig,k) = speccomb &
6236  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
6237  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
6238  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
6239  & + speccomb1 &
6240  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
6241  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
6242  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
6243  & + tauself + taufor + taun2
6244 
6245  fracs(ns15+ig,k) = fracrefa(ig,jpl) + fpl &
6246  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
6247  enddo
6248  enddo
6249 
6250 ! --- ... upper atmosphere loop
6251 
6252  do k = laytrop+1, nlay
6253  do ig = 1, ng15
6254  taug(ns15+ig,k) = f_zero
6255 
6256  fracs(ns15+ig,k) = f_zero
6257  enddo
6258  enddo
6259 
6260 ! ..................................
6261  end subroutine taugb15
6262 ! ----------------------------------
6263 
6264 ! ----------------------------------
6265  subroutine taugb16
6266 ! ..................................
6267 
6268 ! ------------------------------------------------------------------ !
6269 ! band 16: 2600-3250 cm-1 (low key- h2o,ch4; high key - ch4) !
6270 ! ------------------------------------------------------------------ !
6271 
6272  use module_radlw_kgb16
6273 
6274 ! --- locals:
6275  integer :: k, ind0, ind0p, ind1, ind1p, inds, indsp, indf, indfp, &
6276  & id000, id010, id100, id110, id200, id210, jpl, jplp, &
6277  & id001, id011, id101, id111, id201, id211, ig, js, js1
6278 
6279  real (kind=kind_phys) :: tauself, taufor, refrat_planck_a, &
6280  & speccomb, specparm, specmult, fs, &
6281  & speccomb1, specparm1, specmult1, fs1, &
6282  & speccomb_planck,specparm_planck,specmult_planck,fpl, &
6283  & fac000, fac100, fac200, fac010, fac110, fac210, &
6284  & fac001, fac101, fac201, fac011, fac111, fac211, &
6285  & p0, p40, fk00, fk10, fk20, p1, p41, fk01, fk11, fk21
6286 !
6287 !===> ... begin here
6288 !
6289 ! --- ... calculate reference ratio to be used in calculation of Planck
6290 ! fraction in lower atmosphere.
6291 
6292  refrat_planck_a = chi_mls(1,6)/chi_mls(6,6) ! P = 387. mb (Level 6)
6293 
6294 ! --- ... lower atmosphere loop
6295 
6296  do k = 1, laytrop
6297  speccomb = colamt(k,1) + rfrate(k,4,1)*colamt(k,5)
6298  specparm = colamt(k,1) / speccomb
6299  specmult = 8.0 * min(specparm, oneminus)
6300  js = 1 + int(specmult)
6301  fs = mod(specmult, f_one)
6302  ind0 = ((jp(k)-1)*5 + (jt(k)-1)) * nspa(16) + js
6303 
6304  speccomb1 = colamt(k,1) + rfrate(k,4,2)*colamt(k,5)
6305  specparm1 = colamt(k,1) / speccomb1
6306  specmult1 = 8.0 * min(specparm1, oneminus)
6307  js1 = 1 + int(specmult1)
6308  fs1 = mod(specmult1, f_one)
6309  ind1 = (jp(k)*5 + (jt1(k)-1)) * nspa(16) + js1
6310 
6311  speccomb_planck = colamt(k,1) + refrat_planck_a*colamt(k,5)
6312  specparm_planck = colamt(k,1) / speccomb_planck
6313  specmult_planck = 8.0 * min(specparm_planck, oneminus)
6314  jpl = 1 + int(specmult_planck)
6315  fpl = mod(specmult_planck, f_one)
6316 
6317  inds = indself(k)
6318  indf = indfor(k)
6319  indsp = inds + 1
6320  indfp = indf + 1
6321  jplp = jpl + 1
6322 
6323  if (specparm < 0.125 .and. specparm1 < 0.125) then
6324  p0 = fs - f_one
6325  p40 = p0**4
6326  fk00 = p40
6327  fk10 = f_one - p0 - 2.0*p40
6328  fk20 = p0 + p40
6329 
6330  p1 = fs1 - f_one
6331  p41 = p1**4
6332  fk01 = p41
6333  fk11 = f_one - p1 - 2.0*p41
6334  fk21 = p1 + p41
6335 
6336  id000 = ind0
6337  id010 = ind0 + 9
6338  id100 = ind0 + 1
6339  id110 = ind0 +10
6340  id200 = ind0 + 2
6341  id210 = ind0 +11
6342 
6343  id001 = ind1
6344  id011 = ind1 + 9
6345  id101 = ind1 + 1
6346  id111 = ind1 +10
6347  id201 = ind1 + 2
6348  id211 = ind1 +11
6349  elseif (specparm > 0.875 .and. specparm1 > 0.875) then
6350  p0 = -fs
6351  p40 = p0**4
6352  fk00 = p40
6353  fk10 = f_one - p0 - 2.0*p40
6354  fk20 = p0 + p40
6355 
6356  p1 = -fs1
6357  p41 = p1**4
6358  fk01 = p41
6359  fk11 = f_one - p1 - 2.0*p41
6360  fk21 = p1 + p41
6361 
6362  id000 = ind0 + 1
6363  id010 = ind0 +10
6364  id100 = ind0
6365  id110 = ind0 + 9
6366  id200 = ind0 - 1
6367  id210 = ind0 + 8
6368 
6369  id001 = ind1 + 1
6370  id011 = ind1 +10
6371  id101 = ind1
6372  id111 = ind1 + 9
6373  id201 = ind1 - 1
6374  id211 = ind1 + 8
6375  else
6376  fk00 = f_one - fs
6377  fk10 = fs
6378  fk20 = f_zero
6379 
6380  fk01 = f_one - fs1
6381  fk11 = fs1
6382  fk21 = f_zero
6383 
6384  id000 = ind0
6385  id010 = ind0 + 9
6386  id100 = ind0 + 1
6387  id110 = ind0 +10
6388  id200 = ind0
6389  id210 = ind0
6390 
6391  id001 = ind1
6392  id011 = ind1 + 9
6393  id101 = ind1 + 1
6394  id111 = ind1 +10
6395  id201 = ind1
6396  id211 = ind1
6397  endif
6398 
6399  fac000 = fk00 * fac00(k)
6400  fac100 = fk10 * fac00(k)
6401  fac200 = fk20 * fac00(k)
6402  fac010 = fk00 * fac10(k)
6403  fac110 = fk10 * fac10(k)
6404  fac210 = fk20 * fac10(k)
6405 
6406  fac001 = fk01 * fac01(k)
6407  fac101 = fk11 * fac01(k)
6408  fac201 = fk21 * fac01(k)
6409  fac011 = fk01 * fac11(k)
6410  fac111 = fk11 * fac11(k)
6411  fac211 = fk21 * fac11(k)
6412 
6413  do ig = 1, ng16
6414  tauself = selffac(k)* (selfref(ig,inds) + selffrac(k) &
6415  & * (selfref(ig,indsp) - selfref(ig,inds)))
6416  taufor = forfac(k) * (forref(ig,indf) + forfrac(k) &
6417  & * (forref(ig,indfp) - forref(ig,indf)))
6418 
6419  taug(ns16+ig,k) = speccomb &
6420  & * (fac000*absa(ig,id000) + fac010*absa(ig,id010) &
6421  & + fac100*absa(ig,id100) + fac110*absa(ig,id110) &
6422  & + fac200*absa(ig,id200) + fac210*absa(ig,id210)) &
6423  & + speccomb1 &
6424  & * (fac001*absa(ig,id001) + fac011*absa(ig,id011) &
6425  & + fac101*absa(ig,id101) + fac111*absa(ig,id111) &
6426  & + fac201*absa(ig,id201) + fac211*absa(ig,id211)) &
6427  & + tauself + taufor
6428 
6429  fracs(ns16+ig,k) = fracrefa(ig,jpl) + fpl &
6430  & * (fracrefa(ig,jplp) - fracrefa(ig,jpl))
6431  enddo
6432  enddo
6433 
6434 ! --- ... upper atmosphere loop
6435 
6436  do k = laytrop+1, nlay
6437  ind0 = ((jp(k)-13)*5 + (jt(k)-1)) * nspb(16) + 1
6438  ind1 = ((jp(k)-12)*5 + (jt1(k)-1)) * nspb(16) + 1
6439 
6440  ind0p = ind0 + 1
6441  ind1p = ind1 + 1
6442 
6443  do ig = 1, ng16
6444  taug(ns16+ig,k) = colamt(k,5) &
6445  & * (fac00(k)*absb(ig,ind0) + fac10(k)*absb(ig,ind0p) &
6446  & + fac01(k)*absb(ig,ind1) + fac11(k)*absb(ig,ind1p))
6447 
6448  fracs(ns16+ig,k) = fracrefb(ig)
6449  enddo
6450  enddo
6451 
6452 ! ..................................
6453  end subroutine taugb16
6454 ! ----------------------------------
6455 
6456 ! ..................................
6457  end subroutine taumol
6458 !-----------------------------------
6459 
6460 
6461 !
6462 !........................................!
6463  end module module_radlw_main !
6464 !========================================!
6465 
real(kind=kind_phys), parameter con_amo3
molecular wght of o3 (g/mol)
Definition: physcons.f:134
real(kind=kind_phys), dimension(ng09, mfr09), public forref
This module sets up absorption coefficients for band 12: 1800-2080 cm-1 (low - h2o, co2; high - /)
real(kind=kind_phys), dimension(ng10, mfr10), public forref
integer, parameter ngptlw
num of total g-points
Definition: radlw_param.f:110
real(kind=kind_phys), dimension(ng02, msa02), public absa
real(kind=kind_phys), dimension(ng14), public fracrefa
real(kind=kind_phys), dimension(ng08, msb08), public absb
real(kind=kind_phys), dimension(ng01), public fracrefa
real(kind=kind_phys), dimension(nbands) a0
Definition: radlw_main.f:306
subroutine taugb10
Definition: radlw_main.f:5395
real(kind=kind_phys), parameter cldmin
Definition: radlw_main.f:270
real(kind=kind_phys), parameter bpade
Definition: radlw_main.f:271
real(kind=kind_phys), dimension(ng13, mfr13), public forref
real(kind=kind_phys), dimension(ng13, msf13), public selfref
real(kind=kind_phys), dimension(ng01, mfr01), public forref
real(kind=kind_phys), parameter wtdiff
Definition: radlw_main.f:273
real(kind=kind_phys), dimension(ng08, mmc08), public ka_mco2
This module sets up absorption coefficients for band 15: 2380-2600 cm-1 (low - n2o, co2; high - /)
real(kind=kind_phys), dimension(nbands) delwave
Definition: radlw_param.f:152
real(kind=kind_phys), dimension(ng04, maf04), public fracrefa
real(kind=kind_phys), dimension(58, nbands) absliq1
Hu and Stamnes method. the liquid water absorption coefficients are listed for a range of effective r...
Definition: radlw_datatb.f:952
real(kind=kind_phys), dimension(ng04, msf04), public selfref
real(kind=kind_phys), dimension(ng03, maf03), public fracrefa
real(kind=kind_phys), dimension(ng03, maf03, mmn03), public ka_mn2o
real(kind=kind_phys), dimension(ng12, maf12), public fracrefa
real(kind=kind_phys), dimension(ng14, mfr14), public forref
integer, parameter ipsdlw0
Definition: radlw_main.f:339
real(kind=kind_phys), dimension(ng05, mbf05), public fracrefb
real(kind=kind_phys), dimension(ng03, msa03), public absa
This module sets up absorption coefficients for band 06: 820-980 cm-1 (low - h2o; high - /) ...
real(kind=kind_phys), dimension(ng08, msa08), public absa
real(kind=kind_phys), dimension(ng04, msa04), public absa
integer, save ilwcice
lw optical property for ice clouds (only ilwcliq>0) =0:not defined yet =1:input cip...
Definition: physparam.f:116
real(kind=kind_phys), dimension(ng11, mfr11), public forref
real(kind=kind_phys), dimension(ng16, msf16), public selfref
real(kind=kind_phys), dimension(ng06, msa06), public absa
subroutine taugb15
Definition: radlw_main.f:6058
real(kind=kind_phys), dimension(nbands) a2
Definition: radlw_main.f:306
integer, parameter ilwrgas
lw rare gases effect control flag (ch4,n2o,o2,cfcs,...) =0:no; =1:yes.
Definition: physparam.f:101
real(kind=kind_phys), dimension(43, nbands) absice2
for iflagice =2, absice2 are the ice water absorption coefficients used for streamer method...
real(kind=kind_phys), dimension(ng14, msb14), public absb
real(kind=kind_phys), dimension(59) preflog
Definition: radlw_datatb.f:743
real(kind=kind_phys), dimension(ng06), public cfc11adj
integer, parameter maxxsec
num of halocarbon gasees
Definition: radlw_param.f:116
This module sets up absorption coefficients for band 09: 1180-1390 cm-1 (low - h2o, ch4; high - ch4)
real(kind=kind_phys), dimension(ng13, maf13, mmo13), public ka_mco2
This module sets up absorption coefficients for band 11: 1480-1800 cm-1 (low - h2o; high - h2o) ...
real(kind=kind_phys), dimension(ng02), public fracrefb
real(kind=kind_phys), dimension(ng11, msa11), public absa
real(kind=kind_phys), dimension(ng11, msb11), public absb
define type construct for optional radiation flux profiles
Definition: radlw_param.f:94
real(kind=kind_phys), parameter f_one
Definition: radlw_main.f:276
This module sets up absorption coefficients for band 02: 250-500 cm-1 (low - h2o; high - h2o) ...
real(kind=kind_phys), dimension(ng07, maf07), public fracrefa
real(kind=kind_phys), dimension(ng10), public fracrefb
subroutine taugb16
Definition: radlw_main.f:6266
real(kind=kind_phys), parameter absrain
absrain is the rain drop absorption coefficient (m2/g)
Definition: radlw_datatb.f:928
real(kind=kind_phys), dimension(ng04, msb04), public absb
real(kind=kind_phys), dimension(ng07, msb07), public absb
real(kind=kind_phys), dimension(ng13, maf13), public fracrefa
real(kind=kind_phys), dimension(ng16, msb16), public absb
real(kind=kind_phys), dimension(ng01, mmn01), public ka_mn2
integer, dimension(nbands) ipat
ipat is bands index for ebert&curry ice cloud (for iflagice=1)
Definition: radlw_datatb.f:921
real(kind=kind_phys), dimension(ng15, maf15), public fracrefa
real(kind=kind_phys), dimension(ng05, msa05), public absa
This module contains some the most frequently used math and physics constants for gcm models...
Definition: physcons.f:29
real(kind=kind_phys), parameter oneminus
Definition: radlw_main.f:269
real(kind=kind_phys), dimension(ng08, msf08), public selfref
real(kind=kind_phys), dimension(ng16, mfr16), public forref
real(kind=kind_phys), dimension(nbands) a1
Definition: radlw_main.f:306
This module sets up absorption coefficients for band 05: 700-820 cm-1 (low - h2o, co2; high - co2...
real(kind=kind_phys), dimension(ng01, msb01), public absb
subroutine taugb13
Definition: radlw_main.f:5744
real(kind=kind_phys), dimension(ng02, msf02), public selfref
real(kind=kind_phys), dimension(ng15, msa15), public absa
real(kind=kind_phys), parameter amdw
Definition: radlw_main.f:279
real(kind=kind_phys), dimension(7, 59) chi_mls
Definition: radlw_datatb.f:791
real(kind=kind_phys), parameter abssnow0
abssnow0 is the snow flake absorption coefficient (micron), fu coeff
Definition: radlw_datatb.f:932
real(kind=kind_phys), dimension(ng07, mfr07), public forref
subroutine mcica_subcol
This suroutine computes sub-colum cloud profile flag array.
Definition: radlw_main.f:1754
real(kind=kind_phys), parameter amdo3
Definition: radlw_main.f:280
real(kind=kind_phys), dimension(ng06, mfr06), public forref
This module contains reference temperature and pressure.
Definition: radlw_datatb.f:733
real(kind=kind_phys), dimension(ng10, msf10), public selfref
real(kind=kind_phys), dimension(ng02), public fracrefa
subroutine taugb08
Definition: radlw_main.f:5024
integer, save iovrlw
cloud overlapping control flag for lw
Definition: physparam.f:199
real(kind=kind_phys), dimension(nplnk, nbands), public totplnk
Definition: radlw_datatb.f:68
real(kind=kind_phys), dimension(ng14, msf14), public selfref
This module sets up absorption coefficients for band 14: 2250-2380 cm-1 (low - co2; high - co2) ...
This module defines commonly used control variables/parameters in physics related programs...
Definition: physparam.f:27
real(kind=kind_phys), dimension(ng12, msf12), public selfref
real(kind=kind_phys), dimension(46, nbands) absice3
for iflagice = 3, absice3 are the ice water absorption coefficients used for fu parameterization. particle size 5 - 140 micron in increments of 3 microns. units = m2/g. hexagonal ice particle parameterization absorption units (abs coef/iwc): [(m^-1)/(g m^-3)]
subroutine taugb12
Definition: radlw_main.f:5560
real(kind=kind_phys), dimension(nbands) semiss0
Definition: radlw_main.f:327
real(kind=kind_phys), dimension(ng15, msf15), public selfref
This module sets up absorption coefficients for band 01: 10-250 cm-1 (low - h2o; high - h2o) ...
real(kind=kind_phys), dimension(ng09, msa09), public absa
real(kind=kind_phys), dimension(ng12, mfr12), public forref
real(kind=kind_phys), dimension(ng09, mmn09), public kb_mn2o
real(kind=kind_phys), dimension(ng11), public fracrefb
real(kind=kind_phys), dimension(ng07, msf07), public selfref
real(kind=kind_phys), dimension(0:ntbl) tau_tbl
Definition: radlw_main.f:330
real(kind=kind_phys), dimension(ng03, mbf03), public fracrefb
integer, parameter ntbl
lookup table dimension
Definition: radlw_param.f:112
integer, dimension(nbands) nspa
Definition: radlw_main.f:283
integer, save isubclw
sub-column cloud approx flag in lw radiation
Definition: physparam.f:228
integer, dimension(nbands) nspb
Definition: radlw_main.f:283
real(kind=kind_phys), dimension(ng02, mfr02), public forref
real(kind=kind_phys), dimension(ng11, msf11), public selfref
subroutine rtrn
This subroutine computes the upward/downward radiative fluxes, and heating rates for both clear or cl...
Definition: radlw_main.f:2206
real(kind=kind_phys), dimension(ng01, msf01), public selfref
integer, save ivflip
vertical profile indexing flag
Definition: physparam.f:224
real(kind=kind_phys), dimension(ng05, msf05), public selfref
This module sets up absorption coefficients for band 08: 1080-1180 cm-1 (low - h2o; high - o3) ...
real(kind=kind_phys), dimension(ng09), public fracrefb
real(kind=kind_phys), dimension(ng13, mmo13), public kb_mo3
real(kind=kind_phys), dimension(59) tref
Definition: radlw_datatb.f:743
real(kind=kind_phys), dimension(ng13), public fracrefb
real(kind=kind_phys), dimension(ng10, msb10), public absb
This module sets up absorption coefficients for band 10: 1390-1480 cm-1 (low - h2o; high - h2o) ...
real(kind=kind_phys), dimension(ng11, mmo11), public ka_mo2
integer, save ilwcliq
lw optical property for liquid clouds =0:input cld opt depth, ignoring ilwcice setting =1:input c...
Definition: physparam.f:107
define type construct for radiation fluxes at surface
Definition: radlw_param.f:80
real(kind=kind_phys), dimension(ng07), public fracrefb
real(kind=kind_phys), dimension(ng05, maf05, mmo05), public ka_mo3
real(kind=kind_phys), dimension(ng04, mfr04), public forref
real(kind=kind_phys), dimension(ng08), public fracrefb
real(kind=kind_phys), parameter eps
Definition: radlw_main.f:268
subroutine taugb04
Definition: radlw_main.f:4176
real(kind=kind_phys), dimension(ng01), public fracrefb
subroutine setcoef
This subroutine computes various coefficients needed in radiative transfer calculations.
Definition: radlw_main.f:1916
real(kind=kind_phys), dimension(ng05, mfr05), public forref
real(kind=kind_phys), dimension(ng14, msa14), public absa
This module contains LW band parameters set up.
Definition: radlw_param.f:58
This module includes ncep&#39;s modifications of the rrtm-lw radiation ! code from aer inc...
Definition: radlw_main.f:236
real(kind=kind_phys), dimension(ng16, maf16), public fracrefa
real(kind=kind_phys), dimension(ng09, maf09), public fracrefa
real(kind=kind_phys), parameter con_amw
molecular wght of water vapor (g/mol)
Definition: physcons.f:132
integer, parameter ilwrate
lw heating rate unit =1:k/day; =2:k/second.
Definition: physparam.f:97
This module sets up absorption coefficients for band 13: 2080-2250 cm-1 (low - h2o, n2o; high - /)
real(kind=kind_phys), dimension(ng01, mmn01), public kb_mn2
subroutine taugb07
Definition: radlw_main.f:4765
real(kind=kind_phys), dimension(0:ntbl) tfn_tbl
Definition: radlw_main.f:332
real(kind=kind_phys) fluxfac
Definition: radlw_main.f:327
real(kind=kind_phys), dimension(ng07, maf07, mmc07), public ka_mco2
real(kind=kind_phys), dimension(ng08, mmc08), public kb_mn2o
integer, parameter nrates
num of ref rates of binary species
Definition: radlw_param.f:118
real(kind=kind_phys), parameter con_amd
molecular wght of dry air (g/mol)
Definition: physcons.f:130
real(kind=kind_phys), dimension(ng16), public fracrefb
character(40), parameter vtaglw
Definition: radlw_main.f:257
real(kind=kind_phys), parameter f_zero
Definition: radlw_main.f:275
real(kind=kind_phys), dimension(0:ntbl) exp_tbl
Definition: radlw_main.f:331
This module sets up absorption coefficients for band 04: 630-700 cm-1 (low - h2o, co2; high - co2...
real(kind=kind_phys), dimension(ng05), public ccl4
real(kind=kind_phys), dimension(ng03, mbf03, mmn03), public kb_mn2o
real(kind=kind_phys), dimension(2, 5) absice1
for iflagice = 1, absice1 are the ice water absorption coefficients used for ebert and curry method ...
real(kind=kind_phys), dimension(ng08, mmc08), public kb_mco2
subroutine cldprop
This subroutine computes the cloud optical depth(s) for each cloudy layer and g-point interval...
Definition: radlw_main.f:1449
subroutine rtrnmr
Definition: radlw_main.f:2568
real(kind=kind_phys), dimension(ng06, msf06), public selfref
real(kind=kind_phys), parameter con_g
gravity (m/s2)
Definition: physcons.f:52
subroutine rtrnmc
Definition: radlw_main.f:3138
subroutine taugb09
Definition: radlw_main.f:5147
real(kind=kind_phys), dimension(ng09, msb09), public absb
real(kind=kind_phys), dimension(ng15, maf15, mmn15), public ka_mn2
real(kind=kind_phys), dimension(ng10), public fracrefa
subroutine taugb11
Definition: radlw_main.f:5469
subroutine taugb02
Definition: radlw_main.f:3794
This module sets up absorption coefficients for band 07: 980-1080 cm-1 (low - h2o, o3; high - o3)
subroutine, public lwrad
This subroutine is the main lw radiation routine.
Definition: radlw_main.f:415
real(kind=kind_phys), dimension(ng07, msa07), public absa
real(kind=kind_phys), dimension(ng07, mmc07), public kb_mco2
real(kind=kind_phys), dimension(ng08, mmc08), public ka_mo3
subroutine taugb03
Definition: radlw_main.f:3870
real(kind=kind_phys), dimension(ng06), public cfc12
real(kind=kind_phys), dimension(ng15, mfr15), public forref
subroutine taugb14
Definition: radlw_main.f:5990
real(kind=kind_phys), dimension(ng11), public fracrefa
subroutine taugb05
Definition: radlw_main.f:4419
real(kind=kind_phys), dimension(ng13, msa13), public absa
define type construct for radiation fluxes at toa
Definition: radlw_param.f:70
real(kind=kind_phys), dimension(ng05, msb05), public absb
subroutine taumol
Definition: radlw_main.f:3502
real(kind=kind_phys), dimension(ng16, msa16), public absa
real(kind=kind_phys), dimension(ng06, mmc06), public ka_mco2
real(kind=kind_phys), dimension(ng04, mbf04), public fracrefb
real(kind=kind_phys), dimension(ng13, maf13, mmo13), public ka_mco
real(kind=kind_phys) heatfac
Definition: radlw_main.f:327
real(kind=kind_phys), dimension(ng11, mmo11), public kb_mo2
real(kind=kind_phys), parameter con_cp
spec heat air at p (J/kg/K)
Definition: physcons.f:73
This module sets up absorption coefficients for band 16: 2600-3000 cm-1 (low - h2o, ch4; high - /)
real(kind=kind_phys), dimension(ng09, maf09, mmn09), public ka_mn2o
real(kind=kind_phys), dimension(ng03, mfr03), public forref
real(kind=kind_phys), dimension(ng08, mfr08), public forref
real(kind=kind_phys), dimension(ng06), public fracrefa
subroutine taugb06
Definition: radlw_main.f:4680
real(kind=kind_phys), dimension(ng08), public cfc12
real(kind=kind_phys), dimension(ng10, msa10), public absa
This module contains cloud property coefficients.
Definition: radlw_datatb.f:910
real(kind=kind_phys), parameter con_avgd
avogadro constant (1/mol)
Definition: physcons.f:125
integer, parameter maxgas
max num of absorbing gases
Definition: radlw_param.f:114
subroutine, public rlwinit
This subroutine performs calculations necessary for the initialization of the longwave model...
Definition: radlw_main.f:1249
This module sets up absorption coefficients for band 03: 500-630 cm-1 (low - h2o, co2; high - h2o...
real(kind=kind_phys), dimension(ng05, maf05), public fracrefa
real(kind=kind_phys), dimension(ng08), public fracrefa
real(kind=kind_phys), dimension(ng01, msa01), public absa
subroutine taugb01
Definition: radlw_main.f:3688
real(kind=kind_phys), parameter tblint
Definition: radlw_main.f:274
real(kind=kind_phys), dimension(ng03, msf03), public selfref
real(kind=kind_phys), dimension(ng12, msa12), public absa
real(kind=kind_phys), dimension(ng03, msb03), public absb
real(kind=kind_phys), dimension(ng02, msb02), public absb
This module contains plank flux data.
Definition: radlw_datatb.f:58
real(kind=kind_phys), dimension(ng14), public fracrefb
real(kind=kind_phys), dimension(ng09, msf09), public selfref
integer, parameter nbands
num of total spectral bands
Definition: radlw_param.f:108
real(kind=kind_phys), parameter stpfac
Definition: radlw_main.f:272
real(kind=kind_phys), dimension(ng08, mmc08), public ka_mn2o
integer, dimension(ngptlw) ngb
Definition: radlw_param.f:139
real(kind=kind_phys), dimension(ng08), public cfc22adj
integer, save icldflg
cloud optical property scheme control flag
Definition: physparam.f:193