LCOV - code coverage report
Current view: top level - atmos_phys/schemes/cloud_fraction - compute_cloud_fraction.F90 (source / functions) Hit Total Coverage
Test: coverage.info Lines: 130 172 75.6 %
Date: 2025-04-28 18:59:15 Functions: 2 3 66.7 %

          Line data    Source code
       1             : ! Compute cloud fractions using relative humidity (RH) threshold and other methods.
       2             : ! CCPP-ized: Haipeng Lin, February 2025
       3             : module compute_cloud_fraction
       4             :   use ccpp_kinds, only: kind_phys
       5             :   private
       6             :   save
       7             : 
       8             :   public :: compute_cloud_fraction_init
       9             :   public :: compute_cloud_fraction_timestep_init
      10             :   public :: compute_cloud_fraction_run
      11             : 
      12             :   ! Namelist variables
      13             :   logical         :: cldfrc_freeze_dry ! flag for Vavrus correction
      14             :   logical         :: cldfrc_ice        ! flag to compute ice cloud fraction
      15             :   integer         :: iceopt            ! option for ice cloud closure
      16             :                                        ! 1=wang & sassen 2=schiller (iciwc)
      17             :                                        ! 3=wood & field, 4=Wilson (based on smith)
      18             :   real(kind_phys) :: rhminl            ! minimum rh for low stable clouds
      19             :   real(kind_phys) :: rhminl_adj_land   ! rhminl adjustment for snowfree land
      20             :   real(kind_phys) :: rhminh            ! minimum rh for high stable clouds
      21             :   real(kind_phys) :: premit            ! top pressure bound for mid level cloud
      22             :   real(kind_phys) :: premib            ! bottom pressure bound for mid level cloud
      23             :   real(kind_phys) :: icecrit           ! critical RH for ice clouds in Wilson and Ballard closure (smaller = more ice clouds)
      24             : 
      25             :   ! flags
      26             :   logical         :: inversion_cld_off ! turn off stratification-based cloud fraction (inversion_cld)
      27             : 
      28             :   integer         :: k700              ! model level nearest to 700 mb [index]
      29             : 
      30             :   ! tuning constants
      31             :   real(kind_phys), parameter :: pretop = 1.0e2_kind_phys ! pressure bounding high cloud [Pa]
      32             : 
      33             : contains
      34             :   ! Initialize cloud_fraction from namelist parameters
      35             : !> \section arg_table_compute_cloud_fraction_init Argument Table
      36             : !! \htmlinclude compute_cloud_fraction_init.html
      37        1024 :   subroutine compute_cloud_fraction_init( &
      38             :     amIRoot, iulog, &
      39             :     pver, &
      40        1024 :     pref_mid, &
      41             :     inversion_cld_off_in, &
      42             :     cldfrc_freeze_dry_in, cldfrc_ice_in, iceopt_in, &
      43             :     rhminl_in, rhminl_adj_land_in, &
      44             :     rhminh_in, &
      45             :     premit_in, premib_in, &
      46             :     icecrit_in, &
      47        1024 :     errmsg, errflg)
      48             : 
      49             :     ! Input arguments
      50             :     logical,         intent(in)  :: amIRoot
      51             :     integer,         intent(in)  :: iulog
      52             :     integer,         intent(in)  :: pver
      53             :     real(kind_phys), intent(in)  :: pref_mid(:)          ! reference_pressure [Pa]
      54             :     logical,         intent(in)  :: inversion_cld_off_in ! flag to turn off inversion_cld?
      55             :     logical,         intent(in)  :: cldfrc_freeze_dry_in ! flag for Vavrus correction
      56             :     logical,         intent(in)  :: cldfrc_ice_in        ! flag to compute ice cloud fraction
      57             :     integer,         intent(in)  :: iceopt_in            ! option for ice cloud closure
      58             :     real(kind_phys), intent(in)  :: rhminl_in            ! minimum rh for low stable clouds
      59             :     real(kind_phys), intent(in)  :: rhminl_adj_land_in   ! rhminl adjustment for snowfree land
      60             :     real(kind_phys), intent(in)  :: rhminh_in            ! minimum rh for high stable clouds
      61             :     real(kind_phys), intent(in)  :: premit_in            ! top pressure bound for mid level cloud
      62             :     real(kind_phys), intent(in)  :: premib_in            ! bottom pressure bound for mid level cloud
      63             :     real(kind_phys), intent(in)  :: icecrit_in           ! critical RH for ice clouds in Wilson and Ballard closure
      64             : 
      65             :     ! Output arguments
      66             :     character(len=512), intent(out) :: errmsg            ! error message
      67             :     integer,            intent(out) :: errflg            ! error flag
      68             : 
      69        1024 :     errflg = 0
      70        1024 :     errmsg = ''
      71             : 
      72             :     ! Set module variables from input arguments
      73        1024 :     cldfrc_freeze_dry = cldfrc_freeze_dry_in
      74        1024 :     cldfrc_ice        = cldfrc_ice_in
      75        1024 :     iceopt            = iceopt_in
      76        1024 :     rhminl            = rhminl_in
      77        1024 :     rhminl_adj_land   = rhminl_adj_land_in
      78        1024 :     rhminh            = rhminh_in
      79        1024 :     premit            = premit_in
      80        1024 :     premib            = premib_in
      81        1024 :     icecrit           = icecrit_in
      82        1024 :     inversion_cld_off = inversion_cld_off_in
      83             : 
      84        1024 :     if(amIRoot) then
      85           2 :       write(iulog,*) 'tuning parameters compute_cloud_fraction_init: inversion_cld_off ', inversion_cld_off
      86           2 :       write(iulog,*) 'tuning parameters compute_cloud_fraction_init: rhminl ',rhminl,' rhminl_adj_land ',rhminl_adj_land, &
      87           4 :                        ' rhminh ',rhminh,' premit ',premit,' premib ',premib
      88           2 :       write(iulog,*) 'tuning parameters compute_cloud_fraction_init: iceopt ',iceopt,' icecrit ',icecrit
      89             :     endif
      90             : 
      91             :     ! Find vertical level nearest 700 mb.
      92       27648 :     k700 = minloc(abs(pref_mid(:) - 7.e4_kind_phys), 1)
      93             : 
      94        1024 :     if(amIRoot) then
      95           2 :       write(iulog,*) 'compute_cloud_fraction: model level nearest 700 mb is ',k700,' which is ',pref_mid(k700),' pascals '
      96             :     endif
      97             : 
      98        1024 :   end subroutine compute_cloud_fraction_init
      99             : 
     100             :   ! Timestep initialization for cloud fraction
     101             :   ! Specify the perturbation to RH to none by default.
     102             :   ! If a scheme has to perturb RH this quantity should be set by an interstitial
     103             :   ! before calling cloud_fraction a second time.
     104             : !> \section arg_table_compute_cloud_fraction_timestep_init Argument Table
     105             : !! \htmlinclude compute_cloud_fraction_timestep_init.html
     106           0 :   subroutine compute_cloud_fraction_timestep_init(rhpert_flag, errmsg, errflg)
     107             :     logical,            intent(out) :: rhpert_flag
     108             :     character(len=512), intent(out) :: errmsg
     109             :     integer,            intent(out) :: errflg
     110             : 
     111           0 :     errmsg = ''
     112           0 :     errflg = 0
     113             : 
     114           0 :     rhpert_flag = .false.
     115           0 :   end subroutine compute_cloud_fraction_timestep_init
     116             : 
     117             :   ! Compute cloud fraction.
     118             :   !
     119             :   ! Calculates cloud fraction using a relative humidity threshold
     120             :   ! The threshold depends upon pressure, and upon the presence or absence
     121             :   ! of convection as defined by a reasonably large vertical mass flux
     122             :   ! entering that layer from below.
     123             :   !
     124             :   ! Original authors: Many, including Jim McCaa.
     125             : !> \section arg_table_compute_cloud_fraction_run Argument Table
     126             : !! \htmlinclude compute_cloud_fraction_run.html
     127      211176 :   subroutine compute_cloud_fraction_run( &
     128             :     ncol, pver, &
     129             :     cappa, gravit, rair, tmelt, pref, lapse_rate, &
     130             :     top_lev_cloudphys, &
     131      633528 :     pmid, ps, temp, sst, &
     132      422352 :     q, cldice, &
     133      211176 :     phis, &
     134      211176 :     shallowcu, deepcu, concld, & ! convective cloud cover
     135      633528 :     landfrac, ocnfrac, snowh, &
     136             :     rhpert_flag, & ! todo: decide what to do with this
     137      422352 :     cloud, rhcloud, &
     138      211176 :     cldst, &
     139      633528 :     rhu00, icecldf, liqcldf, relhum, &
     140           0 :     errmsg, errflg)
     141             : 
     142             :     ! To-be-ccppized dependencies
     143             :     use wv_saturation, only: qsat, qsat_water, svp_ice_vect
     144             : 
     145             :     ! Arguments
     146             :     integer,         intent(in) :: ncol
     147             :     integer,         intent(in) :: pver
     148             :     real(kind_phys), intent(in) :: cappa
     149             :     real(kind_phys), intent(in) :: gravit
     150             :     real(kind_phys), intent(in) :: rair
     151             :     real(kind_phys), intent(in) :: tmelt
     152             :     real(kind_phys), intent(in) :: pref
     153             :     real(kind_phys), intent(in) :: lapse_rate
     154             :     integer,         intent(in) :: top_lev_cloudphys ! vertical_layer_index_of_cloud_fraction_top [index]
     155             :     real(kind_phys), intent(in) :: pmid(:, :)        ! air_pressure [Pa]
     156             :     real(kind_phys), intent(in) :: ps(:)             ! surface_air_pressure [Pa]
     157             :     real(kind_phys), intent(in) :: temp(:, :)        ! air_temperature [K]
     158             :     real(kind_phys), intent(in) :: sst(:)            ! sea_surface_temperature [K]
     159             :     real(kind_phys), intent(in) :: q(:, :)           ! water_vapor_mixing_ratio_wrt_moist_air_and_condensed_water [kg kg-1]
     160             :     real(kind_phys), intent(in) :: cldice(:, :)      ! cloud_ice_mixing_ratio_wrt_moist_air_and_condensed_water [kg kg-1]
     161             :     real(kind_phys), intent(in) :: phis(:)           ! surface_geopotential [m2 s-2]
     162             :     real(kind_phys), intent(in) :: shallowcu(:, :)   ! shallow convective cloud fraction [fraction]
     163             :     real(kind_phys), intent(in) :: deepcu(:, :)      ! deep convective cloud fraction [fraction]
     164             :     real(kind_phys), intent(in) :: concld(:, :)      ! convective_cloud_area_fraction [fraction]
     165             :     real(kind_phys), intent(in) :: landfrac(:)       ! land_area_fraction [fraction]
     166             :     real(kind_phys), intent(in) :: ocnfrac(:)        ! ocean_area_fraction [fraction]
     167             :     real(kind_phys), intent(in) :: snowh(:)          ! lwe_surface_snow_depth_over_land [m]
     168             :     logical,         intent(in) :: rhpert_flag       ! 0 or 1 to perturb rh
     169             : 
     170             :     ! Output arguments
     171             :     real(kind_phys), intent(out) :: cloud(:, :)      ! cloud_area_fraction [fraction]
     172             :     real(kind_phys), intent(out) :: rhcloud(:, :)    ! cloud fraction
     173             :     real(kind_phys), intent(out) :: cldst(:, :)      ! stratiform_cloud_area_fraction [fraction]
     174             :     real(kind_phys), intent(out) :: rhu00(:, :)      ! RH threshold for cloud [fraction]
     175             :     real(kind_phys), intent(out) :: relhum(:, :)     ! RH for prognostic cldwat [fraction]
     176             :     real(kind_phys), intent(out) :: icecldf(:, :)    ! ice cloud fraction [fraction]
     177             :     real(kind_phys), intent(out) :: liqcldf(:, :)    ! liquid cloud fraction (combined into cloud) [fraction]
     178             :     character(len=512), intent(out) :: errmsg        ! error message
     179             :     integer,            intent(out) :: errflg        ! error flag
     180             : 
     181             :     ! Local variables:
     182             :     real(kind_phys) :: cld                         ! intermediate scratch variable (low cld)
     183      422352 :     real(kind_phys) :: dthdpmn(ncol)               ! most stable lapse rate below 750 mb
     184             :     real(kind_phys) :: dthdp                       ! lapse rate (intermediate variable)
     185      422352 :     real(kind_phys) :: es(ncol, pver)              ! saturation vapor pressure
     186      422352 :     real(kind_phys) :: qs(ncol, pver)              ! saturation specific humidity
     187             :     real(kind_phys) :: rhwght                      ! weighting function for rhlim transition
     188      422352 :     real(kind_phys) :: rh(ncol, pver)              ! relative humidity
     189             :     real(kind_phys) :: rhdif                       ! intermediate scratch variable
     190             :     real(kind_phys) :: strat                       ! intermediate scratch variable
     191      422352 :     real(kind_phys) :: theta(ncol, pver)           ! potential temperature
     192      422352 :     real(kind_phys) :: thetas(ncol)                ! ocean surface potential temperature
     193             :     real(kind_phys) :: rhlim                       ! local rel. humidity threshold estimate
     194             :     real(kind_phys) :: clrsky(ncol)                ! temporary used in random overlap calc
     195      422352 :     real(kind_phys) :: rpdeli(ncol, pver - 1)      ! 1./(pmid(k+1)-pmid(k))
     196             :     real(kind_phys) :: rhpert                      ! the specified perturbation to rh
     197             : 
     198      422352 :     logical  :: cldbnd(ncol)           ! region below high cloud boundary
     199             : 
     200             :     integer  :: i, ierror, k           ! column, level indices
     201             :     integer  :: kp1, ifld
     202      422352 :     integer  :: kdthdp(ncol)
     203             :     integer  :: numkcld                ! number of levels in which to allow clouds
     204             : 
     205             :     !  In Cloud Ice Content variables
     206             :     real(kind_phys) :: a, b, c, as, bs, cs    ! fit parameters
     207             :     real(kind_phys) :: Kc                     ! constant for ice cloud calc (wood & field)
     208             :     real(kind_phys) :: ttmp                   ! limited temperature
     209             :     real(kind_phys) :: icicval                ! empirical iwc value
     210             :     real(kind_phys) :: rho                    ! local air density
     211      422352 :     real(kind_phys) :: esl(ncol, pver)        ! liq sat vapor pressure
     212      211176 :     real(kind_phys) :: esi(ncol, pver)        ! ice sat vapor pressure
     213             :     real(kind_phys) :: ncf, phi               ! Wilson and Ballard parameters
     214             : 
     215             :     ! Statement functions
     216             :     logical land
     217    14524536 :     land(i) = nint(landfrac(i)) == 1
     218             : 
     219      211176 :     errmsg = ''
     220      211176 :     errflg = 0
     221             : 
     222             :     !==================================================================================
     223             :     ! PHILOSOPHY OF PRESENT IMPLEMENTATION
     224             :     ! Modification to philosophy for ice supersaturation
     225             :     ! philosophy below is based on RH water only. This is 'liquid condensation'
     226             :     ! or liquid cloud (even though it will freeze immediately to ice)
     227             :     ! The idea is that the RH limits for condensation are strict only for
     228             :     ! water saturation
     229             :     !
     230             :     ! Ice clouds are formed by explicit parameterization of ice nucleation.
     231             :     ! Closure for ice cloud fraction is done on available cloud ice, such that
     232             :     ! the in-cloud ice content matches an empirical fit
     233             :     ! thus, icecldf = min(cldice/icicval,1) where icicval = f(temp,cldice,numice)
     234             :     ! for a first cut, icicval=f(temp) only.
     235             :     ! Combined cloud fraction is maximum overlap  cloud=max(1,max(icecldf,liqcldf))
     236             :     ! No dA/dt term for ice?
     237             :     !
     238             :     ! There are three co-existing cloud types: convective, inversion related low-level
     239             :     ! stratocumulus, and layered cloud (based on relative humidity).  Layered and
     240             :     ! stratocumulus clouds do not compete with convective cloud for which one creates
     241             :     ! the most cloud.  They contribute collectively to the total grid-box average cloud
     242             :     ! amount.  This is reflected in the way in which the total cloud amount is evaluated
     243             :     ! (a sum as opposed to a logical "or" operation)
     244             :     !
     245             :     !==================================================================================
     246             :     ! set defaults for rhu00 - arbitrary number larger than 1.0 (100%)
     247    85308552 :     rhu00(:, :) = 2.0_kind_phys
     248             :     ! define rh perturbation in order to estimate rhdfda
     249      211176 :     rhpert = 0.01_kind_phys
     250             : 
     251             :     ! Parameters for ice cloud fraction methods
     252             :     ! References for these parameters, and the corresponding methods, are in the compute section
     253             :     ! below.
     254             :     ! iceopt = 1. Wang and Sassen IWC
     255      211176 :     a = 26.87_kind_phys
     256      211176 :     b = 0.569_kind_phys
     257      211176 :     c = 0.002892_kind_phys
     258             :     ! iceopt = 2. Schiller
     259      211176 :     as = -68.4202_kind_phys
     260      211176 :     bs = 0.983917_kind_phys
     261      211176 :     cs = 2.81795_kind_phys
     262             :     ! iceopt = 3. Wood and Field
     263      211176 :     Kc = 75._kind_phys
     264             : 
     265             :     ! Evaluate potential temperature and relative humidity
     266      211176 :     if (cldfrc_ice) then
     267             :       ! If computing ice cloud fraction (CAM5+ MG Morrison and Gettelman microphysics) then water RH
     268           0 :       do k = top_lev_cloudphys, pver
     269           0 :         call qsat_water(temp(1:ncol, k), pmid(1:ncol, k), esl(1:ncol, k), qs(1:ncol, k), ncol)
     270           0 :         call svp_ice_vect(temp(1:ncol, k), esi(1:ncol, k), ncol)
     271             :       end do
     272             :     else
     273             :       ! If not computing ice cloud fraction (CAM4 RK Rasch-Kristjansson microphysics) then hybrid RH,
     274     5701752 :       do k = top_lev_cloudphys, pver
     275     5701752 :         call qsat(temp(1:ncol, k), pmid(1:ncol, k), es(1:ncol, k), qs(1:ncol, k), ncol)
     276             :       end do
     277             :     end if
     278             : 
     279    85308552 :     cloud = 0._kind_phys
     280    85308552 :     icecldf = 0._kind_phys
     281    85308552 :     liqcldf = 0._kind_phys
     282    85308552 :     rhcloud = 0._kind_phys
     283    85308552 :     cldst = 0._kind_phys
     284             : 
     285     5701752 :     do k = top_lev_cloudphys, pver
     286    85097376 :       theta(:ncol, k) = temp(:ncol, k)*(pref/pmid(:ncol, k))**cappa
     287             : 
     288    85308552 :       do i = 1, ncol
     289    79606800 :         if(.not. rhpert_flag) then
     290             :           ! default: no RH perturbation
     291    53071200 :           rh(i, k) = q(i, k)/qs(i, k)
     292             :         else
     293    26535600 :           rh(i, k) = q(i, k)/qs(i, k)*(1.0_kind_phys + rhpert)
     294             :         endif
     295             :         ! record relhum, rh itself will later be modified related with concld
     296    85097376 :         relhum(i, k) = rh(i, k)
     297             :       end do
     298             :     end do
     299             : 
     300             :     ! Initialize other temporary variables
     301      211176 :     ierror = 0
     302     3272976 :     do i = 1, ncol
     303             :       ! Adjust thetas(i) in the presence of non-zero ocean heights.
     304             :       ! This reduces the temperature for positive heights according to a standard lapse rate.
     305     3061800 :       if (ocnfrac(i) > 0.01_kind_phys) then
     306     2217588 :         thetas(i) = (sst(i) - lapse_rate*phis(i)/gravit)*(pref/ps(i))**cappa
     307             :       endif
     308     3272976 :       if (ocnfrac(i) > 0.01_kind_phys .and. sst(i) < 260._kind_phys) then
     309           0 :         ierror = i
     310             :       endif
     311             :     end do
     312             : 
     313      211176 :     if (ierror > 0) then
     314           0 :       write (iulog, *) 'COLDSST: encountered in cldfrc:', ierror, ocnfrac(ierror), sst(ierror)
     315             :     end if
     316             : 
     317     5490576 :     do k = top_lev_cloudphys, pver - 1
     318    82035576 :       rpdeli(:ncol, k) = 1._kind_phys/(pmid(:ncol, k + 1) - pmid(:ncol, k))
     319             :     end do
     320             : 
     321             :     ! shallow and deep convective cloudiness are calculated in the
     322             :     ! convective_cloud_cover CCPP scheme.
     323             : #ifndef PERGRO
     324     5701752 :     do k = top_lev_cloudphys, pver
     325    85308552 :       do i = 1, ncol
     326    85097376 :         rh(i, k) = (rh(i, k) - concld(i, k))/(1.0_kind_phys - concld(i, k))
     327             :       end do
     328             :     end do
     329             : #endif
     330             :     !==================================================================================
     331             :     !
     332             :     !          ****** Compute layer cloudiness ******
     333             :     !
     334             :     !====================================================================
     335             :     ! Begin the evaluation of layered cloud amount based on (modified) RH
     336             :     !====================================================================
     337             :     !
     338      211176 :     numkcld = pver
     339     5490576 :     do k = top_lev_cloudphys + 1, numkcld
     340     5279400 :       kp1 = min(k + 1, pver)
     341    82035576 :       do i = 1, ncol
     342             :         ! This is now designed to apply FOR LIQUID CLOUDS (condensation > RH water)
     343    76545000 :         cldbnd(i) = pmid(i, k) >= pretop
     344             : 
     345    76545000 :         if (pmid(i, k) >= premib) then
     346             :           !==============================================================
     347             :           ! This is the low cloud (below premib) block
     348             :           !==============================================================
     349             :           ! enhance low cloud activation over land with no snow cover
     350    29049072 :           if (land(i) .and. (snowh(i) <= 0.000001_kind_phys)) then
     351     2659221 :             rhlim = rhminl - rhminl_adj_land
     352             :           else
     353    11865315 :             rhlim = rhminl
     354             :           end if
     355             : 
     356    14524536 :           rhdif = (rh(i, k) - rhlim)/(1.0_kind_phys - rhlim)
     357    14524536 :           rhcloud(i, k) = min(0.999_kind_phys, (max(rhdif, 0.0_kind_phys))**2)
     358             : 
     359             :           ! SJV: decrease cloud amount if very low water vapor content
     360             :           ! (thus very cold): "freeze dry"
     361    14524536 :           if (cldfrc_freeze_dry) then
     362    14524536 :             rhcloud(i, k) = rhcloud(i, k)*max(0.15_kind_phys, min(1.0_kind_phys, q(i, k)/0.0030_kind_phys))
     363             :           end if
     364             : 
     365    62020464 :         else if (pmid(i, k) < premit) then
     366             :           !==============================================================
     367             :           ! This is the high cloud (above premit) block
     368             :           !==============================================================
     369             :           !
     370    62020464 :           rhlim = rhminh
     371             :           !
     372    62020464 :           rhdif = (rh(i, k) - rhlim)/(1.0_kind_phys - rhlim)
     373    62020464 :           rhcloud(i, k) = min(0.999_kind_phys, (max(rhdif, 0.0_kind_phys))**2)
     374             :         else
     375             :           !==============================================================
     376             :           ! This is the middle cloud block
     377             :           !==============================================================
     378             :           !
     379             :           !       linear rh threshold transition between thresholds for low & high cloud
     380             :           !
     381           0 :           rhwght = (premib - (max(pmid(i, k), premit)))/(premib - premit)
     382             : 
     383           0 :           if (land(i) .and. (snowh(i) <= 0.000001_kind_phys)) then
     384           0 :             rhlim = rhminh*rhwght + (rhminl - rhminl_adj_land)*(1.0_kind_phys - rhwght)
     385             :           else
     386           0 :             rhlim = rhminh*rhwght + rhminl*(1.0_kind_phys - rhwght)
     387             :           end if
     388           0 :           rhdif = (rh(i, k) - rhlim)/(1.0_kind_phys - rhlim)
     389           0 :           rhcloud(i, k) = min(0.999_kind_phys, (max(rhdif, 0.0_kind_phys))**2)
     390             :         end if
     391             :         !==================================================================================
     392             :         ! WE NEED TO DOCUMENT THE PURPOSE OF THIS TYPE OF CODE (ASSOCIATED WITH 2ND CALL)
     393             :         !==================================================================================
     394             :         !      !
     395             :         !      ! save rhlim to rhu00, it handles well by itself for low/high cloud
     396             :         !      !
     397    76545000 :         rhu00(i, k) = rhlim
     398             :         !==================================================================================
     399             : 
     400    81824400 :         if (cldfrc_ice) then
     401             :           ! Evaluate ice cloud fraction based on in-cloud ice content
     402           0 :           if (iceopt .eq. 1 .or. iceopt .eq. 2) then
     403           0 :             if (iceopt .eq. 1) then
     404             :               ! ICE CLOUD OPTION 1 - Wang & Sassen 2002
     405             :               !         Evaluate desired in-cloud water content
     406             :               !               icicval = f(temp,cldice,numice)
     407             :               !         Start with a function of temperature.
     408             :               !         Wang & Sassen 2002 (JAS), based on ARM site MMCR (midlat cirrus)
     409             :               !         https://doi.org/10.1175/1520-0469(2002)059<2291:CCMPRU>2.0.CO;2
     410             :               !           parameterization valid for 203-253K
     411             :               !           icival > 0 for t>195K
     412           0 :               ttmp = max(195._kind_phys, min(temp(i, k), 253._kind_phys)) - 273.16_kind_phys
     413           0 :               icicval = a + b*ttmp + c*ttmp**2._kind_phys
     414             :               !convert units
     415           0 :               rho = pmid(i, k)/(rair*temp(i, k))
     416           0 :               icicval = icicval*1.e-6_kind_phys/rho
     417             :             else
     418             :               ! ICE CLOUD OPTION 2 - Schiller 2008 (JGR)
     419             :               ! https://doi.org/10.1029/2008JD010342
     420             :               !          Use a curve based on FISH measurements in
     421             :               !          tropics, mid-lats and arctic. Curve is for 180-250K (raise to 273K?)
     422             :               !          use median all flights
     423           0 :               ttmp = max(190._kind_phys, min(temp(i, k), 273.16_kind_phys))
     424           0 :               icicval = 10._kind_phys**(as*bs**ttmp + cs)
     425             :               ! convert units from ppmv to kg kg-1
     426           0 :               icicval = icicval*1.e-6_kind_phys*18._kind_phys/28.97_kind_phys
     427             :             end if
     428             : 
     429             :             ! Set ice cloud fraction for options 1 or 2
     430           0 :             icecldf(i, k) = max(0._kind_phys, min(cldice(i, k)/icicval, 1._kind_phys))
     431           0 :           else if (iceopt .eq. 3) then
     432             :             ! ICE CLOUD OPTION 3 - Wood & Field 2000 (JAS)
     433             :             ! https://doi.org/10.1175/1520-0469(2000)057<1888:RBTWCW>2.0.CO;2
     434             :             ! eq 6: cloud fraction = 1 - exp (-K * qc/qsati)
     435           0 :             icecldf(i, k) = 1._kind_phys - exp(-Kc*cldice(i, k)/(qs(i, k)*(esi(i, k)/esl(i, k))))
     436           0 :             icecldf(i, k) = max(0._kind_phys, min(icecldf(i, k), 1._kind_phys))
     437             :           else
     438             :             ! ICE CLOUD OPTION 4 - Wilson and Ballard 1999
     439             :             ! https://doi.org/10.1002/qj.49712555707
     440             :             ! inversion of smith....
     441             :             !       ncf = cldice / ((1-RHcrit)*qs)
     442             :             ! then a function of ncf....
     443           0 :             ncf = cldice(i, k)/((1._kind_phys - icecrit)*qs(i, k))
     444           0 :             if (ncf <= 0._kind_phys) then
     445           0 :               icecldf(i, k) = 0._kind_phys
     446           0 :             else if (ncf > 0._kind_phys .and. ncf <= 1._kind_phys/6._kind_phys) then
     447           0 :               icecldf(i, k) = 0.5_kind_phys*(6._kind_phys*ncf)**(2._kind_phys/3._kind_phys)
     448           0 :             else if (ncf > 1._kind_phys/6._kind_phys .and. ncf < 1._kind_phys) then
     449           0 :               phi = (acos(3._kind_phys*(1._kind_phys - ncf)/2._kind_phys**(3._kind_phys/2._kind_phys)) + 4._kind_phys*3.1415927_kind_phys)/3._kind_phys
     450           0 :               icecldf(i, k) = (1._kind_phys - 4._kind_phys*cos(phi)*cos(phi))
     451             :             else
     452           0 :               icecldf(i, k) = 1._kind_phys
     453             :             end if
     454           0 :             icecldf(i, k) = max(0._kind_phys, min(icecldf(i, k), 1._kind_phys))
     455             :           end if
     456             : 
     457             :           ! Test code: Combine ice and liquid cloud fraction assuming maximum overlap.
     458             :           ! Combined cloud fraction is maximum overlap
     459             :           !          cloud(i,k)=min(1._kind_phys,max(icecldf(i,k),rhcloud(i,k)))
     460             : 
     461           0 :           liqcldf(i, k) = (1._kind_phys - icecldf(i, k))*rhcloud(i, k)
     462           0 :           cloud(i, k) = liqcldf(i, k) + icecldf(i, k)
     463             :         else ! cldfrc_ice = .false.
     464             :           ! For RK microphysics
     465    76545000 :           cloud(i, k) = rhcloud(i, k)
     466             :         end if
     467             :       end do
     468             :     end do
     469             : 
     470             :     !
     471             :     ! Add in the marine strat
     472             :     ! MARINE STRATUS SHOULD BE A SPECIAL CASE OF LAYERED CLOUD
     473             :     ! CLOUD CURRENTLY CONTAINS LAYERED CLOUD DETERMINED BY RH CRITERIA
     474             :     ! TAKE THE MAXIMUM OF THE DIAGNOSED LAYERED CLOUD OR STRATOCUMULUS
     475             :     !
     476             :     !===================================================================================
     477             :     !
     478             :     !  SOME OBSERVATIONS ABOUT THE FOLLOWING SECTION OF CODE (missed in earlier look)
     479             :     !  K700 IS SET AS A CONSTANT BASED ON HYBRID COORDINATE: IT DOES NOT DEPEND ON
     480             :     !  LOCAL PRESSURE; THERE IS NO PRESSURE RAMP => LOOKS LEVEL DEPENDENT AND
     481             :     !  DISCONTINUOUS IN SPACE (I.E., STRATUS WILL END SUDDENLY WITH NO TRANSITION)
     482             :     !
     483             :     !  IT APPEARS THAT STRAT IS EVALUATED ACCORDING TO KLEIN AND HARTMANN; HOWEVER,
     484             :     !  THE ACTUAL STRATUS AMOUNT (CLDST) APPEARS TO DEPEND DIRECTLY ON THE RH BELOW
     485             :     !  THE STRONGEST PART OF THE LOW LEVEL INVERSION.
     486             :     !PJR answers: 1) the rh limitation is a physical/mathematical limitation
     487             :     !             cant have more cloud than there is RH
     488             :     !             allowed the cloud to exist two layers below the inversion
     489             :     !             because the numerics frequently make 50% relative humidity
     490             :     !             in level below the inversion which would allow no cloud
     491             :     !             2) since  the cloud is only allowed over ocean, it should
     492             :     !             be very insensitive to surface pressure (except due to
     493             :     !             spectral ringing, which also causes so many other problems
     494             :     !             I didnt worry about it.
     495             :     !
     496             :     !==================================================================================
     497      211176 :     if (.not. inversion_cld_off) then
     498             :       !
     499             :       ! Find most stable level below 750 mb for evaluating stratus regimes
     500             :       !
     501     3272976 :       do i = 1, ncol
     502             :         ! Nothing triggers unless a stability greater than this minimum threshold is found
     503     3061800 :         dthdpmn(i) = -0.125_kind_phys
     504     3272976 :         kdthdp(i) = 0
     505             :       end do
     506             : 
     507     5490576 :       do k = top_lev_cloudphys + 1, pver
     508    82035576 :         do i = 1, ncol
     509    81824400 :           if (pmid(i, k) >= premib .and. ocnfrac(i) > 0.01_kind_phys) then
     510             :             ! I think this is done so that dtheta/dp is in units of dg/mb (JJH)
     511    11081184 :             dthdp = 100.0_kind_phys*(theta(i, k) - theta(i, k - 1))*rpdeli(i, k - 1)
     512    11081184 :             if (dthdp < dthdpmn(i)) then
     513      190575 :               dthdpmn(i) = dthdp
     514      190575 :               kdthdp(i) = k     ! index of interface of max inversion
     515             :             end if
     516             :           end if
     517             :         end do
     518             :       end do
     519             : 
     520             :       ! Also check between the bottom layer and the surface
     521             :       ! Only perform this check if the criteria were not met above
     522     3272976 :       do i = 1, ncol
     523     3272976 :         if (kdthdp(i) .eq. 0 .and. ocnfrac(i) > 0.01_kind_phys) then
     524     2030585 :           dthdp = 100.0_kind_phys*(thetas(i) - theta(i, pver))/(ps(i) - pmid(i, pver))
     525     2030585 :           if (dthdp < dthdpmn(i)) then
     526       84191 :             dthdpmn(i) = dthdp
     527       84191 :             kdthdp(i) = pver     ! index of interface of max inversion
     528             :           end if
     529             :         end if
     530             :       end do
     531             : 
     532     3272976 :       do i = 1, ncol
     533     3272976 :         if (kdthdp(i) /= 0) then
     534      271194 :           k = kdthdp(i)
     535      271194 :           kp1 = min(k + 1, pver)
     536             :           ! Note: strat will be zero unless ocnfrac > 0.01
     537      271194 :           strat = min(1._kind_phys, max(0._kind_phys, ocnfrac(i)*((theta(i, k700) - thetas(i))*.057_kind_phys - .5573_kind_phys)))
     538             :           !
     539             :           ! assign the stratus to the layer just below max inversion
     540             :           ! the relative humidity changes so rapidly across the inversion
     541             :           ! that it is not safe to just look immediately below the inversion
     542             :           ! so limit the stratus cloud by rh in both layers below the inversion
     543             :           !
     544      271194 :           cldst(i, k) = min(strat, max(rh(i, k), rh(i, kp1)))
     545             :         end if
     546             :       end do
     547             :     end if  ! .not.inversion_cld_off
     548             : 
     549     5701752 :     do k = top_lev_cloudphys, pver
     550    85308552 :       do i = 1, ncol
     551             :         ! which is greater; standard layered cloud amount or stratocumulus diagnosis
     552    79606800 :         cloud(i, k) = max(rhcloud(i, k), cldst(i, k))
     553             : 
     554             :         ! add in the contributions of convective cloud (determined separately and accounted
     555             :         ! for by modifications to the large-scale relative humidity.
     556    85097376 :         cloud(i, k) = min(cloud(i, k) + concld(i, k), 1.0_kind_phys)
     557             :       end do
     558             :     end do
     559             : 
     560      211176 :   end subroutine compute_cloud_fraction_run
     561             : end module compute_cloud_fraction

Generated by: LCOV version 1.14