LCOV - code coverage report
Current view: top level - physics/rrtmg/aer_src - rrtmg_sw_reftra.f90 (source / functions) Hit Total Coverage
Test: coverage.info Lines: 71 75 94.7 %
Date: 2025-03-13 18:55:17 Functions: 1 1 100.0 %

          Line data    Source code
       1             : !     path:      $Source: /storm/rc1/cvsroot/rc/rrtmg_sw/src/rrtmg_sw_reftra.f90,v $
       2             : !     author:    $Author: mike $
       3             : !     revision:  $Revision: 1.2 $
       4             : !     created:   $Date: 2007/08/23 20:40:14 $
       5             : 
       6             :       module rrtmg_sw_reftra
       7             : 
       8             : !  --------------------------------------------------------------------------
       9             : ! |                                                                          |
      10             : ! |  Copyright 2002-2007, Atmospheric & Environmental Research, Inc. (AER).  |
      11             : ! |  This software may be used, copied, or redistributed as long as it is    |
      12             : ! |  not sold and this copyright notice is reproduced on each copy made.     |
      13             : ! |  This model is provided as is without any express or implied warranties. |
      14             : ! |                       (http://www.rtweb.aer.com/)                        |
      15             : ! |                                                                          |
      16             : !  --------------------------------------------------------------------------
      17             : 
      18             : ! ------- Modules -------
      19             : 
      20             :       use shr_kind_mod, only: r8 => shr_kind_r8
      21             : 
      22             :       use rrsw_tbl, only : tblint, bpade, od_lo, exp_tbl
      23             : 
      24             :       implicit none
      25             : 
      26             :       contains
      27             : 
      28             : ! --------------------------------------------------------------------
      29     3617824 :       subroutine reftra_sw(nlayers, ncol, lrtchk, pgg, prmuz, ptau, pw, &
      30     3617824 :                            pref, prefd, ptra, ptrad)
      31             : ! --------------------------------------------------------------------
      32             :   
      33             : ! Purpose: computes the reflectivity and transmissivity of a clear or 
      34             : !   cloudy layer using a choice of various approximations.
      35             : !
      36             : ! Interface:  *rrtmg_sw_reftra* is called by *rrtmg_sw_spcvrt*
      37             : !
      38             : ! Description:
      39             : ! explicit arguments :
      40             : ! --------------------
      41             : ! inputs
      42             : ! ------ 
      43             : !      lrtchk  = .t. for all layers in clear profile
      44             : !      lrtchk  = .t. for cloudy layers in cloud profile 
      45             : !              = .f. for clear layers in cloud profile
      46             : !      pgg     = asymmetry factor
      47             : !      prmuz   = cosine solar zenith angle
      48             : !      ptau    = optical thickness
      49             : !      pw      = single scattering albedo
      50             : !
      51             : ! outputs
      52             : ! -------
      53             : !      pref    : collimated beam reflectivity
      54             : !      prefd   : diffuse beam reflectivity 
      55             : !      ptra    : collimated beam transmissivity
      56             : !      ptrad   : diffuse beam transmissivity
      57             : !
      58             : !
      59             : ! Method:
      60             : ! -------
      61             : !      standard delta-eddington, p.i.f.m., or d.o.m. layer calculations.
      62             : !      kmodts  = 1 eddington (joseph et al., 1976)
      63             : !              = 2 pifm (zdunkowski et al., 1980)
      64             : !              = 3 discrete ordinates (liou, 1973)
      65             : !
      66             : !
      67             : ! Modifications:
      68             : ! --------------
      69             : ! Original: J-JMorcrette, ECMWF, Feb 2003
      70             : ! Revised for F90 reformatting: MJIacono, AER, Jul 2006
      71             : ! Revised to add exponential lookup table: MJIacono, AER, Aug 2007
      72             : !
      73             : ! ------------------------------------------------------------------
      74             : 
      75             : ! ------- Declarations ------
      76             : 
      77             : ! ------- Input -------
      78             : 
      79             :       integer, intent(in) :: nlayers
      80             :       integer, intent(in) :: ncol
      81             :       logical, intent(in) :: lrtchk(ncol,nlayers)                ! Logical flag for reflectivity and
      82             :                                                                  ! and transmissivity calculation; 
      83             :                                                                  !   Dimensions: (nlayers)
      84             : 
      85             :       real(kind=r8), intent(in) :: pgg(ncol,nlayers)             ! asymmetry parameter
      86             :                                                                  !   Dimensions: (nlayers)
      87             :       real(kind=r8), intent(in) :: ptau(ncol,nlayers)            ! optical depth
      88             :                                                                  !   Dimensions: (nlayers)
      89             :       real(kind=r8), intent(in) :: pw(ncol,nlayers)              ! single scattering albedo 
      90             :                                                                  !   Dimensions: (nlayers)
      91             :       real(kind=r8), intent(in) :: prmuz(ncol)                       ! cosine of solar zenith angle
      92             : 
      93             : ! ------- Output -------
      94             : 
      95             :       real(kind=r8), intent(inout) :: pref(ncol,nlayers+1)       ! direct beam reflectivity
      96             :                                                                  !   Dimensions: (nlayers+1)
      97             :       real(kind=r8), intent(inout) :: prefd(ncol,nlayers+1)      ! diffuse beam reflectivity
      98             :                                                                  !   Dimensions: (nlayers+1)
      99             :       real(kind=r8), intent(inout) :: ptra(ncol,nlayers+1)       ! direct beam transmissivity
     100             :                                                                  !   Dimensions: (nlayers+1)
     101             :       real(kind=r8), intent(inout) :: ptrad(ncol,nlayers+1)      ! diffuse beam transmissivity
     102             :                                                                  !   Dimensions: (nlayers+1)
     103             : 
     104             : ! ------- Local -------
     105             : 
     106             :       integer :: jk, icol, kmodts
     107             :       integer :: itind
     108             : 
     109             :       real(kind=r8) :: tblind
     110             :       real(kind=r8) :: za, za1, za2
     111             :       real(kind=r8) :: zbeta, zdend, zdenr, zdent
     112             :       real(kind=r8) :: ze1, ze2, zem1, zem2, zemm, zep1, zep2
     113             :       real(kind=r8) :: zg, zg3, zgamma1, zgamma2, zgamma3, zgamma4, zgt
     114             :       real(kind=r8) :: zr1, zr2, zr3, zr4, zr5
     115             :       real(kind=r8) :: zrk, zrk2, zrkg, zrm1, zrp, zrp1, zrpp
     116             :       real(kind=r8) :: zsr3, zt1, zt2, zt3, zt4, zt5, zto1
     117             :       real(kind=r8) :: zw, zwcrit, zwo
     118             :       real(kind=r8) :: temp1
     119             :       real(kind=r8), parameter :: eps = 1.e-08_r8
     120             : 
     121             : !     ------------------------------------------------------------------
     122             : 
     123             : ! Initialize
     124             : 
     125     3617824 :       zsr3=sqrt(3._r8)
     126     3617824 :       zwcrit=0.9999995_r8
     127     3617824 :       kmodts=2
     128             : 
     129   340075456 :       do jk=1, nlayers
     130  5402251456 :          do icol = 1,ncol
     131  5398633632 :             if (.not.lrtchk(icol,jk)) then
     132  2344836630 :                pref(icol,jk) =0._r8
     133  2344836630 :                ptra(icol,jk) =1._r8
     134  2344836630 :                prefd(icol,jk)=0._r8
     135  2344836630 :                ptrad(icol,jk)=1._r8
     136             :             else
     137  2717339370 :                zto1=ptau(icol,jk)
     138  2717339370 :                zw  =pw(icol,jk)
     139  2717339370 :                zg  =pgg(icol,jk)  
     140             : 
     141             : ! General two-stream expressions
     142             : 
     143  2717339370 :                zg3= 3._r8 * zg
     144             :                if (kmodts == 1) then
     145             :                   zgamma1= (7._r8 - zw * (4._r8 + zg3)) * 0.25_r8
     146             :                   zgamma2=-(1._r8 - zw * (4._r8 - zg3)) * 0.25_r8
     147             :                   zgamma3= (2._r8 - zg3 * prmuz(icol) ) * 0.25_r8
     148             :                else if (kmodts == 2) then  
     149  2717339370 :                   zgamma1= (8._r8 - zw * (5._r8 + zg3)) * 0.25_r8
     150  2717339370 :                   zgamma2=  3._r8 *(zw * (1._r8 - zg )) * 0.25_r8
     151  2717339370 :                   zgamma3= (2._r8 - zg3 * prmuz(icol) ) * 0.25_r8
     152             :                else if (kmodts == 3) then  
     153             :                   zgamma1= zsr3 * (2._r8 - zw * (1._r8 + zg)) * 0.5_r8
     154             :                   zgamma2= zsr3 * zw * (1._r8 - zg ) * 0.5_r8
     155             :                   zgamma3= (1._r8 - zsr3 * zg * prmuz(icol) ) * 0.5_r8
     156             :                end if
     157  2717339370 :                zgamma4= 1._r8 - zgamma3
     158             :     
     159             : ! Recompute original s.s.a. to test for conservative solution
     160             : 
     161  2717339370 :                temp1 = 1._r8 - 2._r8 * zg
     162  2717339370 :                zwo= zw * (temp1 + zg**2)/(temp1 + zg**2 * zw)
     163  2717339370 :                if (zwo >= zwcrit) then
     164             : ! Conservative scattering
     165             : 
     166    28838517 :                   za  = zgamma1 * prmuz(icol) 
     167    28838517 :                   za1 = za - zgamma3
     168    28838517 :                   zgt = zgamma1 * zto1
     169             :         
     170             : ! Homogeneous reflectance and transmittance,
     171             : ! collimated beam
     172             : 
     173    28838517 :                   ze1 = min ( zto1 / prmuz(icol) , 500._r8)
     174    28838517 :                   ze2 = exp( -ze1 )
     175             : 
     176             : ! Use exponential lookup table for transmittance, or expansion of 
     177             : ! exponential for low tau
     178             : !               if (ze1 .le. od_lo) then 
     179             : !                  ze2 = 1._r8 - ze1 + 0.5_r8 * ze1 * ze1
     180             : !               else
     181             : !                  tblind = ze1 / (bpade + ze1)
     182             : !                  itind = tblint * tblind + 0.5_r8
     183             : !                  ze2 = exp_tbl(itind)
     184             : !               endif
     185             : !
     186             : 
     187    28838517 :                   pref(icol,jk) = (zgt - za1 * (1._r8 - ze2)) / ( 1._r8 + zgt)
     188    28838517 :                   ptra(icol,jk) = 1._r8 - pref(icol,jk)
     189             : 
     190             : ! isotropic incidence
     191             : 
     192    28838517 :                   prefd(icol,jk) = zgt / ( 1._r8 + zgt)
     193    28838517 :                   ptrad(icol,jk) = 1._r8 - prefd(icol,jk)        
     194             : 
     195             : ! This is applied for consistency between total (delta-scaled) and direct (unscaled) 
     196             : ! calculations at very low optical depths (tau < 1.e-4) when the exponential lookup
     197             : ! table returns a transmittance of 1.0.
     198    28838517 :                   if (ze2 .eq. 1.0_r8) then 
     199           0 :                      pref(icol,jk) = 0.0_r8
     200           0 :                      ptra(icol,jk) = 1.0_r8
     201           0 :                      prefd(icol,jk) = 0.0_r8
     202           0 :                      ptrad(icol,jk) = 1.0_r8
     203             :                   endif
     204             : 
     205             :                else
     206             : ! Non-conservative scattering
     207             : 
     208  2688500853 :                   za1 = zgamma1 * zgamma4 + zgamma2 * zgamma3
     209  2688500853 :                   za2 = zgamma1 * zgamma3 + zgamma2 * zgamma4
     210  2688500853 :                   zrk = sqrt ( zgamma1**2 - zgamma2**2)
     211  2688500853 :                   zrp = zrk * prmuz(icol)               
     212  2688500853 :                   zrp1 = 1._r8 + zrp
     213  2688500853 :                   zrm1 = 1._r8 - zrp
     214  2688500853 :                   zrk2 = 2._r8 * zrk
     215  2688500853 :                   zrpp = 1._r8 - zrp*zrp
     216  2688500853 :                   zrkg = zrk + zgamma1
     217  2688500853 :                   zr1  = zrm1 * (za2 + zrk * zgamma3)
     218  2688500853 :                   zr2  = zrp1 * (za2 - zrk * zgamma3)
     219  2688500853 :                   zr3  = zrk2 * (zgamma3 - za2 * prmuz(icol) )
     220  2688500853 :                   zr4  = zrpp * zrkg
     221  2688500853 :                   zr5  = zrpp * (zrk - zgamma1)
     222  2688500853 :                   zt1  = zrp1 * (za1 + zrk * zgamma4)
     223  2688500853 :                   zt2  = zrm1 * (za1 - zrk * zgamma4)
     224  2688500853 :                   zt3  = zrk2 * (zgamma4 + za1 * prmuz(icol) )
     225  2688500853 :                   zt4  = zr4
     226  2688500853 :                   zt5  = zr5
     227  2688500853 :                   zbeta = (zgamma1 - zrk) / zrkg !- zr5 / zr4
     228             :         
     229             : ! Homogeneous reflectance and transmittance
     230             : 
     231  2688500853 :                   ze1 = min ( zrk * zto1, 500._r8)
     232  2688500853 :                   ze2 = min ( zto1 / prmuz(icol) , 500._r8)
     233             : !
     234             : ! Original
     235             : !              zep1 = exp( ze1 )
     236             : !              zem1 = exp(-ze1 )
     237             : !              zep2 = exp( ze2 )
     238             : !              zem2 = exp(-ze2 )
     239             : !
     240             : ! Revised original, to reduce exponentials
     241  2688500853 :                   zep1 = exp( ze1 )
     242  2688500853 :                   zem1 = 1._r8 / zep1
     243  2688500853 :                   zep2 = exp( ze2 )
     244  2688500853 :                   zem2 = 1._r8 / zep2
     245             : !
     246             : ! Use exponential lookup table for transmittance, or expansion of 
     247             : ! exponential for low tau
     248             : !               if (ze1 .le. od_lo) then 
     249             : !                  zem1 = 1._r8 - ze1 + 0.5_r8 * ze1 * ze1
     250             : !                  zep1 = 1._r8 / zem1
     251             : !               else
     252             : !                  tblind = ze1 / (bpade + ze1)
     253             : !                  itind = tblint * tblind + 0.5_r8
     254             : !                  zem1 = exp_tbl(itind)
     255             : !                  zep1 = 1._r8 / zem1
     256             : !               endif
     257             : !
     258             : !               if (ze2 .le. od_lo) then 
     259             : !                  zem2 = 1._r8 - ze2 + 0.5_r8 * ze2 * ze2
     260             : !                  zep2 = 1._r8 / zem2
     261             : !               else
     262             : !                  tblind = ze2 / (bpade + ze2)
     263             : !                  itind = tblint * tblind + 0.5_r8
     264             : !                  zem2 = exp_tbl(itind)
     265             : !                  zep2 = 1._r8 / zem2
     266             : !               endif
     267             : 
     268             : ! collimated beam
     269             : 
     270  2688500853 :                   zdenr = zr4*zep1 + zr5*zem1
     271  2688500853 :                   zdent = zt4*zep1 + zt5*zem1
     272  2688500853 :                   if (zdenr .ge. -eps .and. zdenr .le. eps) then
     273           4 :                      pref(icol,jk) = eps
     274           4 :                      ptra(icol,jk) = zem2
     275             :                   else
     276  2688500849 :                      pref(icol,jk) = zw * (zr1*zep1 - zr2*zem1 - zr3*zem2) / zdenr
     277  2688500849 :                      ptra(icol,jk) = zem2 - zem2 * zw * (zt1*zep1 - zt2*zem1 - zt3*zep2) / zdent
     278             :                   endif 
     279             : 
     280             : ! diffuse beam
     281             : 
     282  2688500853 :                   zemm = zem1*zem1
     283  2688500853 :                   zdend = 1._r8 / ( (1._r8 - zbeta*zemm ) * zrkg)
     284  2688500853 :                   prefd(icol,jk) =  zgamma2 * (1._r8 - zemm) * zdend
     285  2688500853 :                   ptrad(icol,jk) =  zrk2*zem1*zdend
     286             : 
     287             :                endif
     288             : 
     289             :             endif         
     290             : 
     291             :          enddo    
     292             :       end do
     293     3617824 :       end subroutine reftra_sw
     294             : 
     295             :       end module rrtmg_sw_reftra
     296             : 

Generated by: LCOV version 1.14