LCOV - code coverage report
Current view: top level - physics/rrtmg/aer_src - mcica_subcol_gen_lw.f90 (source / functions) Hit Total Coverage
Test: coverage.info Lines: 71 81 87.7 %
Date: 2025-03-14 01:18:36 Functions: 2 2 100.0 %

          Line data    Source code
       1             : !     path:      $Source: /storm/rc1/cvsroot/rc/rrtmg_lw/src/mcica_subcol_gen_lw.f90,v $
       2             : !     author:    $Author: mike $
       3             : !     revision:  $Revision: 1.3 $
       4             : !     created:   $Date: 2007/08/28 22:38:11 $
       5             : !
       6             : 
       7             : module mcica_subcol_gen_lw
       8             : 
       9             : !  --------------------------------------------------------------------------
      10             : ! |                                                                          |
      11             : ! |  Copyright 2006-2007, Atmospheric & Environmental Research, Inc. (AER).  |
      12             : ! |  This software may be used, copied, or redistributed as long as it is    |
      13             : ! |  not sold and this copyright notice is reproduced on each copy made.     |
      14             : ! |  This model is provided as is without any express or implied warranties. |
      15             : ! |                       (http://www.rtweb.aer.com/)                        |
      16             : ! |                                                                          |
      17             : !  --------------------------------------------------------------------------
      18             : 
      19             : ! Purpose: Create McICA stochastic arrays for cloud physical or optical properties.
      20             : ! Two options are possible:
      21             : ! 1) Input cloud physical properties: cloud fraction, ice and liquid water
      22             : !    paths, ice fraction, and particle sizes.  Output will be stochastic
      23             : !    arrays of these variables.  (inflag = 1)
      24             : ! 2) Input cloud optical properties directly: cloud optical depth, single
      25             : !    scattering albedo and asymmetry parameter.  Output will be stochastic
      26             : !    arrays of these variables.  (inflag = 0; longwave scattering is not
      27             : !    yet available, ssac and asmc are for future expansion)
      28             : 
      29             : ! --------- Modules ----------
      30             : 
      31             : use shr_kind_mod,     only: r8 => shr_kind_r8
      32             : use cam_abortutils,   only: endrun
      33             : 
      34             : use parrrtm, only : nbndlw, ngptlw
      35             : use rrlw_wvn, only: ngb
      36             : 
      37             : implicit none
      38             : private
      39             : 
      40             : public :: mcica_subcol_lw
      41             : 
      42             : !=========================================================================================
      43             : contains
      44             : !=========================================================================================
      45             : 
      46      115200 : subroutine mcica_subcol_lw(lchnk, ncol, nlay, icld, permuteseed, play, &
      47       38400 :                        cldfrac, ciwp, clwp, rei, rel, tauc, cldfmcl, &
      48       76800 :                        ciwpmcl, clwpmcl, reicmcl, relqmcl, taucmcl)
      49             : 
      50             :    ! ----- Input -----
      51             :    ! Control
      52             :    integer, intent(in) :: lchnk           ! chunk identifier
      53             :    integer, intent(in) :: ncol            ! number of columns
      54             :    integer, intent(in) :: nlay            ! number of model layers
      55             :    integer, intent(in) :: icld            ! clear/cloud, cloud overlap flag
      56             :    integer, intent(in) :: permuteseed     ! if the cloud generator is called multiple times, 
      57             :                                                         ! permute the seed between each call.
      58             :                                                         ! between calls for LW and SW, recommended
      59             :                                                         ! permuteseed differes by 'ngpt'
      60             : 
      61             :    ! Atmosphere
      62             :    real(kind=r8), intent(in) :: play(:,:)          ! layer pressures (mb) 
      63             :                                                         !    Dimensions: (ncol,nlay)
      64             : 
      65             :    ! Atmosphere/clouds - cldprop
      66             :    real(kind=r8), intent(in) :: cldfrac(:,:)       ! layer cloud fraction
      67             :                                                      !    Dimensions: (ncol,nlay)
      68             :    real(kind=r8), intent(in) :: tauc(:,:,:)        ! cloud optical depth
      69             :                                                      !    Dimensions: (nbndlw,ncol,nlay)
      70             :    real(kind=r8), intent(in) :: ciwp(:,:)          ! cloud ice water path
      71             :                                                      !    Dimensions: (ncol,nlay)
      72             :    real(kind=r8), intent(in) :: clwp(:,:)          ! cloud liquid water path
      73             :                                                      !    Dimensions: (ncol,nlay)
      74             :    real(kind=r8), intent(in) :: rei(:,:)           ! cloud ice particle size
      75             :                                                      !    Dimensions: (ncol,nlay)
      76             :    real(kind=r8), intent(in) :: rel(:,:)           ! cloud liquid particle size
      77             :                                                      !    Dimensions: (ncol,nlay)
      78             : 
      79             :    ! ----- Output -----
      80             :    ! Atmosphere/clouds - cldprmc [mcica]
      81             :    real(kind=r8), intent(out) :: cldfmcl(:,:,:)    ! cloud fraction [mcica]
      82             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
      83             :    real(kind=r8), intent(out) :: ciwpmcl(:,:,:)    ! cloud ice water path [mcica]
      84             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
      85             :    real(kind=r8), intent(out) :: clwpmcl(:,:,:)    ! cloud liquid water path [mcica]
      86             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
      87             :    real(kind=r8), intent(out) :: relqmcl(:,:)      ! liquid particle size (microns)
      88             :                                                      !    Dimensions: (ncol,nlay)
      89             :    real(kind=r8), intent(out) :: reicmcl(:,:)      ! ice partcle size (microns)
      90             :                                                      !    Dimensions: (ncol,nlay)
      91             :    real(kind=r8), intent(out) :: taucmcl(:,:,:)    ! cloud optical depth [mcica]
      92             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
      93             :    ! ----- Local -----
      94             : 
      95             :    ! Stochastic cloud generator variables [mcica]
      96             :    integer, parameter :: nsubclw = ngptlw ! number of sub-columns (g-point intervals)
      97             :    integer :: km, im, nm                  ! loop indices
      98             : 
      99           0 :    real(kind=r8) :: pmid(ncol, nlay)               ! layer pressures (Pa) 
     100             :    !----------------------------------------------------------------------------
     101             : 
     102             :    ! Return if clear sky; or stop if icld out of range
     103       38400 :    if (icld.eq.0) return
     104       38400 :    if (icld.lt.0.or.icld.gt.3) then 
     105           0 :       call endrun('MCICA_SUBCOL: INVALID ICLD')
     106             :    end if
     107             : 
     108             :    ! NOTE: For GCM mode, permuteseed must be offset between LW and SW by at
     109             :    !       least the number of subcolumns
     110             : 
     111             :    ! Pass particle sizes to new arrays, no subcolumns for these properties yet
     112             :    ! Convert pressures from mb to Pa
     113             : 
     114    33154560 :    reicmcl(:ncol,:nlay) = rei(:ncol,:nlay)
     115    33154560 :    relqmcl(:ncol,:nlay) = rel(:ncol,:nlay)
     116    33154560 :    pmid(:ncol,:nlay)    = play(:ncol,:nlay)*1.e2_r8
     117             : 
     118             :    ! Generate the stochastic subcolumns of cloud optical properties for the longwave;
     119             :    call generate_stochastic_clouds( &
     120             :       ncol, nlay, nsubclw, icld, pmid, &
     121             :       cldfrac, clwp, ciwp, tauc, cldfmcl, &
     122       38400 :       clwpmcl, ciwpmcl, taucmcl, permuteseed)
     123             : 
     124       38400 : end subroutine mcica_subcol_lw
     125             : 
     126             : !=========================================================================================
     127             : 
     128       38400 : subroutine generate_stochastic_clouds( &
     129       76800 :    ncol, nlay, nsubcol, icld, pmid, &
     130       38400 :    cld, clwp, ciwp, tauc, cld_stoch, &
     131       38400 :    clwp_stoch, ciwp_stoch, tauc_stoch, changeSeed) 
     132             : 
     133             :    !----------------------------------------------------------------------------------------------------------------
     134             :    ! ---------------------
     135             :    ! Contact: Cecile Hannay (hannay@ucar.edu)
     136             :    ! 
     137             :    ! Original code: Based on Raisanen et al., QJRMS, 2004.
     138             :    ! 
     139             :    ! Modifications: Generalized for use with RRTMG and added Mersenne Twister as the default
     140             :    !   random number generator, which can be changed to the optional kissvec random number generator
     141             :    !   with flag 'irnd' below. Some extra functionality has been commented or removed.  
     142             :    !   Michael J. Iacono, AER, Inc., February 2007
     143             :    !
     144             :    ! Given a profile of cloud fraction, cloud water and cloud ice, we produce a set of subcolumns.
     145             :    ! Each layer within each subcolumn is homogeneous, with cloud fraction equal to zero or one 
     146             :    ! and uniform cloud liquid and cloud ice concentration.
     147             :    ! The ensemble as a whole reproduces the probability function of cloud liquid and ice within each layer 
     148             :    ! and obeys an overlap assumption in the vertical.   
     149             :    ! 
     150             :    ! Overlap assumption:
     151             :    !  The cloud are consistent with 4 overlap assumptions: random, maximum, maximum-random and exponential. 
     152             :    !  The default option is maximum-random (option 3)
     153             :    !  The options are: 1=random overlap, 2=max/random, 3=maximum overlap, 4=exponential overlap
     154             :    !  This is set with the variable "overlap" 
     155             :    !mji - Exponential overlap option (overlap=4) has been deactivated in this version
     156             :    !  The exponential overlap uses also a length scale, Zo. (real,    parameter  :: Zo = 2500. ) 
     157             :    ! 
     158             :    ! Seed:
     159             :    !  If the stochastic cloud generator is called several times during the same timestep, 
     160             :    !  one should change the seed between the call to insure that the subcolumns are different.
     161             :    !  This is done by changing the argument 'changeSeed'
     162             :    !  For example, if one wants to create a set of columns for the shortwave and another set for the longwave ,
     163             :    !  use 'changeSeed = 1' for the first call and'changeSeed = 2' for the second call 
     164             :    !
     165             :    ! PDF assumption:
     166             :    !  We can use arbitrary complicated PDFS. 
     167             :    !  In the present version, we produce homogeneuous clouds (the simplest case).  
     168             :    !  Future developments include using the PDF scheme of Ben Johnson. 
     169             :    !
     170             :    ! History file:
     171             :    !  Option to add diagnostics variables in the history file. (using FINCL in the namelist)
     172             :    !  nsubcol = number of subcolumns
     173             :    !  overlap = overlap type (1-3)
     174             :    !  Zo = length scale 
     175             :    !  CLOUD_S = mean of the subcolumn cloud fraction ('_S" means Stochastic)
     176             :    !  CLDLIQ_S = mean of the subcolumn cloud water
     177             :    !  CLDICE_S = mean of the subcolumn cloud ice 
     178             :    !
     179             :    ! Note:
     180             :    !   Here: we force that the cloud condensate to be consistent with the cloud fraction 
     181             :    !   i.e we only have cloud condensate when the cell is cloudy. 
     182             :    !   In CAM: The cloud condensate and the cloud fraction are obtained from 2 different equations 
     183             :    !   and the 2 quantities can be inconsistent (i.e. CAM can produce cloud fraction 
     184             :    !   without cloud condensate or the opposite).
     185             :    !---------------------------------------------------------------------------------------------------------------
     186             : 
     187             :    use shr_RandNum_mod, only: ShrIntrinsicRandGen, ShrKissRandGen, &
     188             :                               ShrF95MtRandGen, ShrDsfmtRandGen
     189             : 
     190       38400 :    type(ShrDsfmtRandGen) :: dsfmt_gen
     191       38400 :    type(ShrKissRandGen) :: kiss_gen
     192             : 
     193             :    ! -- Arguments
     194             : 
     195             :    integer, intent(in) :: ncol            ! number of columns
     196             :    integer, intent(in) :: nlay            ! number of layers
     197             :    integer, intent(in) :: icld            ! clear/cloud, cloud overlap flag
     198             :    integer, intent(in) :: nsubcol         ! number of sub-columns (g-point intervals)
     199             :    integer, optional, intent(in) :: changeSeed     ! allows permuting seed
     200             : 
     201             :    real(kind=r8), intent(in) :: pmid(:,:)          ! layer pressure (Pa)
     202             :                                                      !    Dimensions: (ncol,nlay)
     203             :    real(kind=r8), intent(in) :: cld(:,:)           ! cloud fraction 
     204             :                                                      !    Dimensions: (ncol,nlay)
     205             :    real(kind=r8), intent(in) :: clwp(:,:)          ! cloud liquid water path
     206             :                                                      !    Dimensions: (ncol,nlay)
     207             :    real(kind=r8), intent(in) :: ciwp(:,:)          ! cloud ice water path
     208             :                                                      !    Dimensions: (ncol,nlay)
     209             :    real(kind=r8), intent(in) :: tauc(:,:,:)        ! cloud optical depth
     210             :                                                      !    Dimensions: (nbndlw,ncol,nlay)
     211             : 
     212             :    real(kind=r8), intent(out) :: cld_stoch(:,:,:)  ! subcolumn cloud fraction 
     213             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
     214             :    real(kind=r8), intent(out) :: clwp_stoch(:,:,:) ! subcolumn cloud liquid water path
     215             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
     216             :    real(kind=r8), intent(out) :: ciwp_stoch(:,:,:) ! subcolumn cloud ice water path
     217             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
     218             :    real(kind=r8), intent(out) :: tauc_stoch(:,:,:) ! subcolumn cloud optical depth
     219             :                                                      !    Dimensions: (ngptlw,ncol,nlay)
     220             :    ! -- Local variables
     221       76800 :    real(kind=r8) :: cldf(ncol,nlay)                ! cloud fraction 
     222             :     
     223             :    ! Constants (min value for cloud fraction and cloud water and ice)
     224             :    real(kind=r8), parameter :: cldmin = 1.0e-80_r8     ! min cloud fraction
     225             : 
     226             :    ! Variables related to random number and seed 
     227             :    integer :: irnd                        ! flag for random number generator
     228             :                                                         !  0 = kissvec
     229             :                                                         !  1 = Mersenne Twister
     230             : 
     231       76800 :    real(kind=r8), dimension(nsubcol, ncol, nlay) :: CDF, CDF2      ! random numbers
     232       76800 :    integer, dimension(ncol,4) :: kiss_seed
     233       76800 :    real(kind=r8), dimension(ncol,1) :: rand_num_1d ! random number (kissvec)
     234       76800 :    real(kind=r8), dimension(ncol,nlay) :: rand_num ! random number (kissvec)
     235       38400 :    integer, dimension(ncol) :: iseed                       ! seed to create random number (Mersenne Teister)
     236             : 
     237             :    ! Flag to identify cloud fraction in subcolumns
     238       76800 :    logical,  dimension(nsubcol, ncol, nlay) :: iscloudy   ! flag that says whether a gridbox is cloudy
     239             : 
     240             :    ! Indices
     241             :    integer :: ilev, isubcol, i, n         ! indices
     242             :    !----------------------------------------------------------------------------
     243             : 
     244             :    ! Set randum number generator to use (0 = kissvec; 1 = mersennetwister)
     245       38400 :    irnd = 0
     246             : 
     247             :    ! ensure that cloud fractions are in bounds 
     248    33154560 :    cldf(:,:) = cld(:ncol,:nlay)
     249    33154560 :    where (cldf(:,:) < cldmin)
     250             :       cldf(:,:) = 0._r8
     251             :    end where
     252             : 
     253             :    ! ----- Create seed  --------
     254             :    
     255             :    ! Advance randum number generator by changeseed values
     256             :    if (irnd == 0) then   
     257             : 
     258             :       ! For kissvec, create a seed that depends on the state of the columns. Maybe not the best way, but it works.  
     259             :       ! Must use pmid from bottom four layers. 
     260      591360 :       do i=1,ncol
     261      552960 :          if (pmid(i,nlay).lt.pmid(i,nlay-1)) then 
     262           0 :             call endrun('MCICA_SUBCOL: KISSVEC SEED GENERATOR REQUIRES PMID FROM BOTTOM FOUR LAYERS.')
     263             :          end if
     264      552960 :          kiss_seed(i,1) = (pmid(i,nlay)   - int(pmid(i,nlay)))    * 1000000000
     265      552960 :          kiss_seed(i,2) = (pmid(i,nlay-1) - int(pmid(i,nlay-1)))  * 1000000000
     266      552960 :          kiss_seed(i,3) = (pmid(i,nlay-2) - int(pmid(i,nlay-2)))  * 1000000000
     267      591360 :          kiss_seed(i,4) = (pmid(i,nlay-3) - int(pmid(i,nlay-3)))  * 1000000000
     268             :       end do
     269             : 
     270       38400 :       kiss_gen = ShrKissRandGen(kiss_seed)
     271             : 
     272     5798400 :       do i = 1, changeSeed
     273     5798400 :          call kiss_gen%random(rand_num_1d)
     274             :       end do
     275             :    elseif (irnd.eq.1) then
     276             : 
     277             :       do i = 1, ncol
     278             :          if (pmid(i,nlay) < pmid(i,nlay-1)) then
     279             :             call endrun('MCICA_SUBCOL: MT SEED GENERATOR REQUIRES PMID FROM BOTTOM FOUR LAYERS.')
     280             :          end if
     281             :          kiss_seed(i,1) = (pmid(i,nlay)   - int(pmid(i,nlay)))    * 1000000000
     282             :          kiss_seed(i,2) = (pmid(i,nlay-1) - int(pmid(i,nlay-1)))  * 1000000000
     283             :          kiss_seed(i,3) = (pmid(i,nlay-2) - int(pmid(i,nlay-2)))  * 1000000000
     284             :          kiss_seed(i,4) = (pmid(i,nlay-3) - int(pmid(i,nlay-3)))  * 1000000000
     285             :       end do
     286             : 
     287             :       iseed = kiss_seed(:,1) + kiss_seed(:,2) + kiss_seed(:,3) + kiss_seed(:,4)
     288             :       dsfmt_gen =ShrDsfmtRandGen(iseed,1)
     289             : 
     290             :    end if
     291             : 
     292             :    ! ------ Apply overlap assumption --------
     293             : 
     294             :    ! generate the random numbers  
     295             : 
     296       38400 :    select case (icld)
     297             : 
     298             :    case(1) 
     299             :       ! Random overlap
     300             :       ! i) pick a random value at every level
     301             :   
     302             :       if (irnd == 0) then 
     303           0 :          do isubcol = 1,nsubcol
     304           0 :             call kiss_gen%random(rand_num)
     305           0 :             CDF(isubcol,:,:) = rand_num(:,:)
     306             :          end do
     307             :       else if (irnd == 1) then
     308             :          do isubcol = 1, nsubcol
     309             :             call dsfmt_gen%random(rand_num)
     310             :             CDF(isubcol,:,:) = rand_num(:,:) 
     311             :          end do
     312             :       end if
     313             : 
     314             :    case(2) 
     315             :       ! Maximum-Random overlap
     316             :       ! i) pick a random number for top layer.
     317             :       ! ii) walk down the column: 
     318             :       !    - if the layer above is cloudy, we use the same random number than in the layer above
     319             :       !    - if the layer above is clear, we use a new random number 
     320             : 
     321             :       if (irnd == 0) then 
     322     5414400 :          do isubcol = 1, nsubcol
     323     5376000 :             call kiss_gen%random(rand_num)
     324  4641676800 :             CDF(isubcol,:,:) = rand_num(:,:)
     325             :          end do
     326             :       elseif (irnd == 1) then
     327             :          do isubcol = 1, nsubcol
     328             :             call dsfmt_gen%random(rand_num)
     329             :             CDF(isubcol,:,:) = rand_num(:,:)
     330             :          end do
     331             :       end if
     332             : 
     333     2150400 :       do ilev = 2,nlay
     334    32563200 :          do i = 1, ncol
     335  4290316800 :             do isubcol = 1, nsubcol
     336  4288204800 :                if (CDF(isubcol, i, ilev-1) > 1._r8 - cldf(i,ilev-1) ) then
     337   601287544 :                   CDF(isubcol,i,ilev) = CDF(isubcol,i,ilev-1) 
     338             :                else
     339  3656504456 :                   CDF(isubcol,i,ilev) = CDF(isubcol,i,ilev) * (1._r8 - cldf(i,ilev-1))
     340             :                end if
     341             :             end do
     342             :          end do
     343             :       end do
     344             :        
     345             :    case(3) 
     346             :       ! Maximum overlap
     347             :       ! i) pick the same random numebr at every level  
     348             : 
     349       38400 :       if (irnd.eq.0) then 
     350           0 :          do isubcol = 1,nsubcol
     351           0 :             call kiss_gen%random(rand_num_1d)
     352           0 :             do ilev = 1,nlay
     353           0 :                CDF(isubcol,:,ilev) = rand_num_1d(:,1)
     354             :             enddo
     355             :          enddo
     356             :       elseif (irnd.eq.1) then
     357             :          do isubcol = 1, nsubcol
     358             :             call dsfmt_gen%random(rand_num_1d)
     359             :             do ilev = 1, nlay
     360             :                CDF(isubcol,:,ilev) = rand_num_1d(:,1)
     361             :             enddo
     362             :          enddo
     363             :       endif
     364             : 
     365             :    end select
     366             :  
     367             :    ! -- generate subcolumns for homogeneous clouds -----
     368     2188800 :    do ilev = 1,nlay
     369  4368360960 :       iscloudy(:,:,ilev) = (CDF(:,:,ilev) >= 1._r8 - spread(cldf(:,ilev), dim=1, nCopies=nsubcol) )
     370             :    end do
     371             : 
     372             :    ! where the subcolumn is cloudy, the subcolumn cloud fraction is 1;
     373             :    ! where the subcolumn is not cloudy, the subcolumn cloud fraction is 0
     374             : 
     375     2188800 :    do ilev = 1,nlay
     376    33154560 :       do i = 1, ncol
     377  4368322560 :          do isubcol = 1, nsubcol
     378  4366172160 :             if (iscloudy(isubcol,i,ilev) ) then
     379   616669603 :                cld_stoch(isubcol,i,ilev) = 1._r8
     380             :             else
     381  3718536797 :                cld_stoch(isubcol,i,ilev) = 0._r8
     382             :             endif
     383             :          end do
     384             :       end do
     385             :    end do
     386             : 
     387             :    ! where there is a cloud, set the subcolumn cloud properties;
     388             :    ! incoming clwp, ciwp and tauc should be in-cloud quantites and not grid-averaged quantities
     389             : 
     390     2188800 :    do ilev = 1, nlay
     391    33154560 :       do i = 1, ncol
     392  4368322560 :          do isubcol = 1, nsubcol
     393  4366172160 :             if ( iscloudy(isubcol,i,ilev) .and. (cldf(i,ilev) > 0._r8) ) then
     394   616669603 :                clwp_stoch(isubcol,i,ilev) = clwp(i,ilev)
     395   616669603 :                ciwp_stoch(isubcol,i,ilev) = ciwp(i,ilev)
     396             :             else
     397  3718536797 :                clwp_stoch(isubcol,i,ilev) = 0._r8
     398  3718536797 :                ciwp_stoch(isubcol,i,ilev) = 0._r8
     399             :             end if
     400             :          end do
     401             :       end do
     402             :    end do
     403             : 
     404     2188800 :    do ilev = 1, nlay
     405    33154560 :       do i = 1, ncol
     406  4368322560 :          do isubcol = 1,nsubcol
     407  4366172160 :             if ( iscloudy(isubcol,i,ilev) .and. (cldf(i,ilev) > 0._r8) ) then
     408   616669603 :                n = ngb(isubcol)
     409   616669603 :                tauc_stoch(isubcol,i,ilev) = tauc(n,i,ilev)
     410             :             else
     411  3718536797 :                tauc_stoch(isubcol,i,ilev) = 0._r8
     412             :             end if
     413             :          end do
     414             :       end do
     415             :    end do
     416             : 
     417             :    if (irnd == 0) then 
     418       38400 :       call kiss_gen%finalize()
     419             :    else if (irnd == 1) then
     420             :       call dsfmt_gen%finalize()
     421             :    end if
     422             : 
     423       38400 : end subroutine generate_stochastic_clouds
     424             : 
     425             : end module mcica_subcol_gen_lw

Generated by: LCOV version 1.14