Radiation Scheme in CCPP
module_radlw_main Module Reference

This module includes ncep's modifications of the rrtm-lw radiation ! code from aer inc.

Functions/Subroutines

subroutine, public rlwinit
 This subroutine performs calculations necessary for the initialization of the longwave model. lookup tables are computed for use in the lw radiative transfer, and input absorption coefficient data for each spectral band are reduced from 256 g-point intervals to 140. More...
 
subroutine cldprop
 This subroutine computes the cloud optical depth(s) for each cloudy layer and g-point interval. More...
 
subroutine mcica_subcol
 This suroutine computes sub-colum cloud profile flag array. More...
 
subroutine setcoef
 This subroutine computes various coefficients needed in radiative transfer calculations. More...
 
subroutine rtrn
 This subroutine computes the upward/downward radiative fluxes, and heating rates for both clear or cloudy atmosphere. Clouds assumed as randomly overlaping in a vertical column. More...
 
subroutine rtrnmr
 
subroutine rtrnmc
 
subroutine taumol
 
subroutine, public lwrad
 This subroutine is the main lw radiation routine. More...
 

Variables

character(40), parameter vtaglw ='NCEP LW v5.1 Nov 2012 -RRTMG-LW v4.82 '
 
real(kind=kind_phys), parameter eps = 1.0e-6
 
real(kind=kind_phys), parameter oneminus = 1.0-eps
 
real(kind=kind_phys), parameter cldmin = 1.0e-80
 
real(kind=kind_phys), parameter bpade = 1.0/0.278
 
real(kind=kind_phys), parameter stpfac = 296.0/1013.0
 
real(kind=kind_phys), parameter wtdiff = 0.5
 
real(kind=kind_phys), parameter tblint = ntbl
 
real(kind=kind_phys), parameter f_zero = 0.0
 
real(kind=kind_phys), parameter f_one = 1.0
 
real(kind=kind_phys), parameter amdw = con_amd/con_amw
 
real(kind=kind_phys), parameter amdo3 = con_amd/con_amo3
 
integer, dimension(nbands) nspa
 
integer, dimension(nbands) nspb
 
real(kind=kind_phys), dimension(nbands) a0
 
real(kind=kind_phys), dimension(nbands) a1
 
real(kind=kind_phys), dimension(nbands) a2
 
logical lhlwb = .false.
 
logical lhlw0 = .false.
 
logical lflxprf = .false.
 
real(kind=kind_phys) fluxfac
 
real(kind=kind_phys) heatfac
 
real(kind=kind_phys), dimension(nbands) semiss0
 
real(kind=kind_phys), dimension(0:ntbl) tau_tbl
 
real(kind=kind_phys), dimension(0:ntbl) exp_tbl
 
real(kind=kind_phys), dimension(0:ntbl) tfn_tbl
 
integer, parameter ipsdlw0 = ngptlw
 

Function/Subroutine Documentation

subroutine module_radlw_main::cldprop ( )
private

Definition at line 1449 of file radlw_main.f.

References module_radlw_cldprlw::absice1, module_radlw_cldprlw::absice2, module_radlw_cldprlw::absice3, module_radlw_cldprlw::absliq1, module_radlw_cldprlw::absrain, module_radlw_cldprlw::abssnow0, cldmin, f_one, f_zero, physparam::ilwcice, physparam::ilwcliq, module_radlw_cldprlw::ipat, physparam::isubclw, mcica_subcol(), module_radlw_parameters::nbands, and module_radlw_parameters::ngptlw.

Referenced by lwrad().

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 ! ..................................
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(43, nbands) absice2
for iflagice =2, absice2 are the ice water absorption coefficients used for streamer method...
real(kind=kind_phys), parameter absrain
absrain is the rain drop absorption coefficient (m2/g)
Definition: radlw_datatb.f:928
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), parameter abssnow0
abssnow0 is the snow flake absorption coefficient (micron), fu coeff
Definition: radlw_datatb.f:932
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)]
real(kind=kind_phys), dimension(2, 5) absice1
for iflagice = 1, absice1 are the ice water absorption coefficients used for ebert and curry method ...
This module contains cloud property coefficients.
Definition: radlw_datatb.f:910

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine, public module_radlw_main::lwrad ( )
Parameters
[in]plyr(npts,nlay),model layer mean pressure in mb
[in]plvl(npts,nlp1),model interface pressure in mb
[in]tlyr(npts,nlay),model layer mean temperature in K
[in]tlvl(npts,nlp1),model interface temperature in K
[in]qlyr(npts,nlay),layer specific humidity in gm/gm *see inside
[in]olyr(npts,nlay),layer ozone concentration in gm/gm *see inside
[in]gasvmr(npts,nlay,:),atmospheric gases amount:(check module_radiation_gases for definition)
(:,:,1) - co2 volume mixing ratio
(:,:,2) - n2o volume mixing ratio
(:,:,3) - ch4 volume mixing ratio
(:,:,4) - o2 volume mixing ratio
(:,:,5) - co volume mixing ratio
(:,:,6) - cfc11 volume mixing ratio
(:,:,7) - cfc12 volume mixing ratio
(:,:,8) - cfc22 volume mixing ratio
(:,:,9) - ccl4 volume mixing ratio
[in]clouds(npts,nlay,:),layer cloud profile(check module_radiation_clouds for definition)
— for ilwcliq > 0 —
(:,:,1) - layer total cloud fraction
(:,:,2) - layer in-cloud liq water path ( \( g/m^2 \))
(:,:,3) - mean eff radius for liq cloud (micron)
(:,:,4) - layer in-cloud ice water path ( \( g/m^2 \))
(:,:,5) - mean eff radius for ice cloud (micron)
(:,:,6) - layer rain drop water path ( \( g/m^2 \))
(:,:,7) - mean eff radius for rain drop (micron)
(:,:,8) - layer snow flake water path ( \( g/m^2 \))
(:,:,9) - mean eff radius for snow flake (micron)
— for ilwcliq = 0 —
(:,:,1) - layer total cloud fraction
(:,:,2) - layer cloud optical depth
(:,:,3) - layer cloud single scattering albedo
(:,:,4) - layer cloud asymmetry factor
[in]icseed(npts),auxiliary special cloud related array.when module variable isubclw=2, it provides permutation seed for each column profile that are used for generating random numbers.when isubclw /=2, it will not be used.
[in]aerosols(npts,nlay,nbdsw,:),aerosol optical properties (check module_radiation_aerosols for definition)
(:,:,:,1) - optical depth
(:,:,:,2) - single scattering albedo
(:,:,:,3) - asymmetry parameter
[in]sfemis(npts),surface emissivity
[in]sfgtmp(npts),surface ground temperature in K
[in]nptstotal number of horizontal points
[in]nlay,nlp1total number of vertical layers, levels
[in]lprntcntl flag for diagnostic print out
[out]hlwc(npts,nlay),total sky heating rate in k/day or k/sec
[out]topflx(npts),radiation fluxes at top, components (check module_radlw_parameters for definition)
upfxc - total sky upward flux at top ( \( w/m^2 \))
upfx0 - clear sky upward flux at top ( \( w/m^2 \))
[out]sfcflx(npts),radiation fluxes at sfc, components (check module_radlw_parameters for definition)
upfxc - total sky upward flux at sfc ( \( w/m^2 \))
dnfxc - total sky downward flux at sfc ( \( w/m^2 \))
upfx0 - clear sky upward flux at sfc ( \( w/m^2 \))
dnfx0 - clear sky downward flux at sfc ( \( w/m^2 \))
Optional:
[out]hlwb(npts,nlay,nbands),spectral band total sky heating rates
[out]hlw0(npts,nlay),clear sky heating rates (k/sec or k/day)
[out]flxprf(npts,nlp1),level radiation fluxes ( \( w/m^2 \)), components (check module_radlw_parameters for definition)
dnfxc - total sky downward flux
upfxc - total sky upward flux
dnfx0 - clear sky downward flux
upfx0 - clear sky upward flux

General Algorithm

  1. Change random number seed value for each radiation invocation (isubclw =1 or 2)
  2. Read surface emissivity
  3. Prepare atmospheric profile for use in rrtm
  4. Set absorber amount for h2o, co2, and o3
  5. Set up column amount for rare gases n2o,ch4,o2,co,ccl4,cf11,cf12,cf22, convert from volume mixing ratio to molec/cm2 based on coldry (scaled to 1.0e-20).
  6. Set aerosol optical properties
  7. Read cloud optical properties
  8. Compute precipitable water vapor for diffusivity angle adjustments
  9. Compute column amount for broadening gases
  10. Compute diffusivity angle adjustments
  11. For cloudy atmosphere, call cldprop (Cloud Optical Property in LW) to set cloud optical properties
  12. Calling setcoef (LW Coefficients calculation) to compute various coefficients needed in radiative transfer calculations.
  13. Call taumol to calculte the gaseous optical depths and Plank fractions for each longwave spectral band.
  14. Call the radiative transfer routine based on cloud scheme selection. compute the upward/downward radiative fluxes, and heating rates for both clear or cloudy atmosphere.
    - call rtrn(): clouds are assumed as randomly overlaping in a vertical column
    - call rtrnmr (The Longwave Ratiative Transfer Routine): clouds are assumed as in maximum-randomly overlaping in a vertical column;
    - call rtrnmc: clouds are treated with the mcica stochastic approach.
  15. Save outputs

Definition at line 415 of file radlw_main.f.

References a0, a1, a2, amdo3, amdw, cldprop(), physcons::con_amd, physcons::con_amw, physcons::con_avgd, physcons::con_g, eps, f_one, f_zero, physparam::ilwcliq, physparam::ilwrgas, physparam::iovrlw, ipsdlw0, physparam::isubclw, physparam::ivflip, lflxprf, lhlw0, lhlwb, module_radlw_parameters::maxgas, module_radlw_parameters::nbands, rtrn(), rtrnmc(), rtrnmr(), semiss0, setcoef(), and taumol().

Referenced by module_radiation_driver::grrad().

415 
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 !...................................

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine module_radlw_main::mcica_subcol ( )
Parameters
[in]cldfreal, (nlay), layer cloud fraction
[in]nlayinteger, 1, number of model vertical layers
[in]ipseedinteger, 1, permute seed for random num generator
[out]lcloudylogical, (ngptlw,nlay),sub-colum cloud profile flag array

External Module Variables


physparam::iovrlw - control flag for cloud overlapping method
=0:random; =1:maximum/random: =2:maximum

Definition at line 1754 of file radlw_main.f.

References f_one, physparam::iovrlw, and module_radlw_parameters::ngptlw.

Referenced by cldprop().

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 ! ..................................

Here is the caller graph for this function:

subroutine, public module_radlw_main::rlwinit ( )
Parameters
[in]meinteger, print control for parallel process
[out]NONE

External Module Variables


physparam::ilwrate - heating rate unit selections
=1: output in k/day
=2: output in k/second
physparam::ilwrgas - control flag for rare gases (ch4,n2o,o2,cfcs, etc.)
=0: do not include rare gases
>0: include all rare gases
physparam::ilwcliq - liquid cloud optical properties contrl flag
=0: input cloud opt depth from diagnostic scheme
>0: input cwp,rew, and other cloud content parameters
physparam::isubclw - sub-column cloud approximation control flag
=0: no sub-col cld treatment, use grid-mean cld quantities
=1: mcica sub-col, prescribed seeds to get random numbers
=2: mcica sub-col, providing array icseed for random numbers
physparam::icldflg - cloud scheme control flag
=0: diagnostic scheme gives cloud tau, omiga, and g.
=1: prognostic scheme gives cloud liq/ice path, etc.
physparam::iovrlw - clouds vertical overlapping control flag
=0: random overlapping clouds
=1: maximum/random overlapping clouds
=2: maximum overlap cloud (isubcol>0 only)

Definition at line 1249 of file radlw_main.f.

References bpade, physcons::con_cp, physcons::con_g, exp_tbl, f_one, f_zero, fluxfac, heatfac, physparam::icldflg, physparam::ilwcliq, physparam::ilwrate, physparam::ilwrgas, physparam::iovrlw, physparam::isubclw, module_radlw_parameters::ntbl, semiss0, tau_tbl, tfn_tbl, and vtaglw.

Referenced by module_radiation_driver::radinit().

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 !...................................

Here is the caller graph for this function:

subroutine module_radlw_main::rtrn ( )
private

Original Code Description: this program calculates the upward fluxes, downward fluxes, and heating rates for an arbitrary clear or cloudy atmosphere. The input to this program is the atmospheric profile, all Planck function information, and the cloud fraction by layer. A variable diffusivity angle (secdif) is used for the angle integration. Bands 2-3 and 5-9 use a value for secdif that varies from 1.50 to 1.80 as a function of the column water vapor, and other bands use a value of 1.66. The gaussian weight appropriate to this angle (wtdiff =0.5) is applied here. Note that use of the emissivity angle for the flux integration can cause errors of 1 to 4 \(W/m^2\) within cloudy layers. Clouds are treated with a random cloud overlap method.

Parameters
[in]semissreal,(nbands), lw surface emissivity
[in]delpreal, (nlay), layer pressure thickness (mb)
[in]cldfrcreal, (0:nlp1), layer cloud fraction
[in]taucldreal, (nbands,nlay), layer cloud opt depth
[in]tautotreal, (ngptlw,nlay), total optical depth (gas+aerosols)
[in]pklayreal, (nbands,0:nlay), integrated planck func at lay temp
[in]pklevreal, (nbands,0:nlay), integrated planck func at lev temp
[in]fracsreal, (ngptlw,nlay), planck fractions
[in]secdifreal, (nbands), secant of diffusivity angle
[in]nlayinteger, 1, number of vertical layers
[in]nlp1integer, 1, number of vertical levels (interfaces)
[out]totufluxreal, (0:nlay), total sky upward flux \((w/m^2)\)
[out]totdfluxreal, (0:nlay), total sky downward flux \((w/m^2)\)
[out]htrreal, (nlay), total sky heating rate (k/sec or k/day)
[out]totuclflreal, (0:nlay), clear sky upward flux \((w/m^2)\)
[out]totdclflreal, (0:nlay), clear sky downward flux \((w/m^2)\)
[out]htrclreal, (nlay), clear sky heating rate (k/sec or k/day)
[out]htrbreal, (nlay,nbands), spectral band lw heating rate (k/day)

Definition at line 2206 of file radlw_main.f.

References bpade, eps, exp_tbl, f_one, f_zero, fluxfac, heatfac, lhlw0, lhlwb, module_radlw_parameters::nbands, module_radlw_parameters::ngb, module_radlw_parameters::ngptlw, tau_tbl, tblint, tfn_tbl, and wtdiff.

Referenced by lwrad().

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 ! ..................................

Here is the caller graph for this function:

subroutine module_radlw_main::rtrnmc ( )
private

Definition at line 3138 of file radlw_main.f.

References bpade, eps, exp_tbl, f_one, f_zero, fluxfac, heatfac, lhlw0, lhlwb, module_radlw_parameters::nbands, module_radlw_parameters::ngb, module_radlw_parameters::ngptlw, tau_tbl, tblint, tfn_tbl, and wtdiff.

Referenced by lwrad().

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 ! ..................................

Here is the caller graph for this function:

subroutine module_radlw_main::rtrnmr ( )
private

Definition at line 2568 of file radlw_main.f.

References bpade, eps, exp_tbl, f_one, f_zero, fluxfac, heatfac, lhlw0, lhlwb, module_radlw_parameters::nbands, module_radlw_parameters::ngb, module_radlw_parameters::ngptlw, tau_tbl, tblint, tfn_tbl, and wtdiff.

Referenced by lwrad().

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 ! .................................

Here is the caller graph for this function:

subroutine module_radlw_main::setcoef ( )
private

Definition at line 1916 of file radlw_main.f.

References module_radlw_ref::chi_mls, module_radlw_parameters::delwave, f_one, f_zero, module_radlw_parameters::nbands, module_radlw_ref::preflog, stpfac, module_radlw_avplank::totplnk, and module_radlw_ref::tref.

Referenced by lwrad().

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 ! ..................................

Here is the caller graph for this function:

subroutine module_radlw_main::taumol ( )
private

Definition at line 3502 of file radlw_main.f.

References module_radlw_parameters::ngb, module_radlw_parameters::ngptlw, taugb01(), taugb02(), taugb03(), taugb04(), taugb05(), taugb06(), taugb07(), taugb08(), taugb09(), taugb10(), taugb11(), taugb12(), taugb13(), taugb14(), taugb15(), and taugb16().

Referenced by lwrad().

Here is the call graph for this function:

Here is the caller graph for this function:

Variable Documentation

real (kind=kind_phys), dimension(nbands) module_radlw_main::a0
private

Definition at line 306 of file radlw_main.f.

Referenced by lwrad().

306  real (kind=kind_phys), dimension(nbands) :: a0, a1, a2
real (kind=kind_phys), dimension(nbands) module_radlw_main::a1
private

Definition at line 306 of file radlw_main.f.

Referenced by lwrad().

real (kind=kind_phys), dimension(nbands) module_radlw_main::a2
private

Definition at line 306 of file radlw_main.f.

Referenced by lwrad().

real (kind=kind_phys), parameter module_radlw_main::amdo3 = con_amd/con_amo3
private

Definition at line 280 of file radlw_main.f.

Referenced by lwrad().

280  real (kind=kind_phys), parameter :: amdo3 = con_amd/con_amo3
real (kind=kind_phys), parameter module_radlw_main::amdw = con_amd/con_amw
private

Definition at line 279 of file radlw_main.f.

Referenced by lwrad().

279  real (kind=kind_phys), parameter :: amdw = con_amd/con_amw
real (kind=kind_phys), parameter module_radlw_main::bpade = 1.0/0.278
private

Definition at line 271 of file radlw_main.f.

Referenced by rlwinit(), rtrn(), rtrnmc(), and rtrnmr().

271  real (kind=kind_phys), parameter :: bpade = 1.0/0.278 ! pade approx constant
real (kind=kind_phys), parameter module_radlw_main::cldmin = 1.0e-80
private

Definition at line 270 of file radlw_main.f.

Referenced by cldprop().

270  real (kind=kind_phys), parameter :: cldmin = 1.0e-80
real (kind=kind_phys), parameter module_radlw_main::eps = 1.0e-6
private

Definition at line 268 of file radlw_main.f.

Referenced by lwrad(), rtrn(), rtrnmc(), and rtrnmr().

268  real (kind=kind_phys), parameter :: eps = 1.0e-6
real (kind=kind_phys), dimension(0:ntbl) module_radlw_main::exp_tbl
private

Definition at line 331 of file radlw_main.f.

Referenced by rlwinit(), rtrn(), rtrnmc(), and rtrnmr().

331  real (kind=kind_phys) :: exp_tbl(0:ntbl) !transmittance lookup table
real (kind=kind_phys), parameter module_radlw_main::f_one = 1.0
private

Definition at line 276 of file radlw_main.f.

Referenced by cldprop(), lwrad(), mcica_subcol(), rlwinit(), rtrn(), rtrnmc(), rtrnmr(), setcoef(), taugb01(), taugb02(), taugb03(), taugb04(), taugb05(), taugb07(), taugb09(), taugb12(), taugb13(), taugb15(), and taugb16().

276  real (kind=kind_phys), parameter :: f_one = 1.0
real (kind=kind_phys), parameter module_radlw_main::f_zero = 0.0
private

Definition at line 275 of file radlw_main.f.

Referenced by cldprop(), lwrad(), rlwinit(), rtrn(), rtrnmc(), rtrnmr(), setcoef(), taugb03(), taugb04(), taugb05(), taugb07(), taugb09(), taugb12(), taugb13(), taugb15(), and taugb16().

275  real (kind=kind_phys), parameter :: f_zero = 0.0
real (kind=kind_phys) module_radlw_main::fluxfac
private

Definition at line 327 of file radlw_main.f.

Referenced by rlwinit(), rtrn(), rtrnmc(), and rtrnmr().

327  real (kind=kind_phys) :: fluxfac, heatfac, semiss0(nbands)
real (kind=kind_phys) module_radlw_main::heatfac
private

Definition at line 327 of file radlw_main.f.

Referenced by rlwinit(), rtrn(), rtrnmc(), and rtrnmr().

integer, parameter module_radlw_main::ipsdlw0 = ngptlw
private

Definition at line 339 of file radlw_main.f.

Referenced by lwrad().

339  integer, parameter :: ipsdlw0 = ngptlw ! initial permutation seed
logical module_radlw_main::lflxprf = .false.
private

Definition at line 319 of file radlw_main.f.

Referenced by lwrad().

319  logical :: lflxprf= .false.
logical module_radlw_main::lhlw0 = .false.
private

Definition at line 318 of file radlw_main.f.

Referenced by lwrad(), rtrn(), rtrnmc(), and rtrnmr().

318  logical :: lhlw0 = .false.
logical module_radlw_main::lhlwb = .false.
private

Definition at line 317 of file radlw_main.f.

Referenced by lwrad(), rtrn(), rtrnmc(), and rtrnmr().

317  logical :: lhlwb = .false.
integer, dimension(nbands) module_radlw_main::nspa
private

Definition at line 283 of file radlw_main.f.

Referenced by taugb01(), taugb02(), taugb03(), taugb04(), taugb05(), taugb06(), taugb07(), taugb08(), taugb09(), taugb10(), taugb11(), taugb12(), taugb13(), taugb14(), taugb15(), and taugb16().

283  integer, dimension(nbands) :: nspa, nspb
integer, dimension(nbands) module_radlw_main::nspb
private
real (kind=kind_phys), parameter module_radlw_main::oneminus = 1.0-eps
private

Definition at line 269 of file radlw_main.f.

Referenced by taugb03(), taugb04(), taugb05(), taugb07(), taugb09(), taugb12(), taugb13(), taugb15(), and taugb16().

269  real (kind=kind_phys), parameter :: oneminus= 1.0-eps
real (kind=kind_phys), dimension(nbands) module_radlw_main::semiss0
private

Definition at line 327 of file radlw_main.f.

Referenced by lwrad(), and rlwinit().

real (kind=kind_phys), parameter module_radlw_main::stpfac = 296.0/1013.0
private

Definition at line 272 of file radlw_main.f.

Referenced by setcoef().

272  real (kind=kind_phys), parameter :: stpfac = 296.0/1013.0
real (kind=kind_phys), dimension(0:ntbl) module_radlw_main::tau_tbl
private

Definition at line 330 of file radlw_main.f.

Referenced by rlwinit(), rtrn(), rtrnmc(), and rtrnmr().

330  real (kind=kind_phys) :: tau_tbl(0:ntbl) !clr-sky opt dep (for cldy transfer)
real (kind=kind_phys), parameter module_radlw_main::tblint = ntbl
private

Definition at line 274 of file radlw_main.f.

Referenced by rtrn(), rtrnmc(), and rtrnmr().

274  real (kind=kind_phys), parameter :: tblint = ntbl ! lookup table conversion factor
real (kind=kind_phys), dimension(0:ntbl) module_radlw_main::tfn_tbl
private

Definition at line 332 of file radlw_main.f.

Referenced by rlwinit(), rtrn(), rtrnmc(), and rtrnmr().

332  real (kind=kind_phys) :: tfn_tbl(0:ntbl) !tau transition function; i.e. the
character(40), parameter module_radlw_main::vtaglw ='NCEP LW v5.1 Nov 2012 -RRTMG-LW v4.82 '
private

Definition at line 257 of file radlw_main.f.

Referenced by rlwinit().

257  character(40), parameter :: &
258  & VTAGLW='NCEP LW v5.1 Nov 2012 -RRTMG-LW v4.82 '
real (kind=kind_phys), parameter module_radlw_main::wtdiff = 0.5
private

Definition at line 273 of file radlw_main.f.

Referenced by rtrn(), rtrnmc(), and rtrnmr().

273  real (kind=kind_phys), parameter :: wtdiff = 0.5 ! weight for radiance to flux conversion