Class cloud_utils
In: radiation/cloud_utils.f90

雲関系ルーチン

Cloud-related routines

Note that Japanese and English are described in parallel.

雲の分布を設定.

In this module, the amount of cloud or cloud optical depth are set. This module is under development and is still a preliminary version.

Procedures List

!$ ! RadiationFluxDennouAGCM :放射フラックスの計算
!$ ! ———— :————
!$ ! RadiationFluxDennouAGCM :Calculate radiation flux

NAMELIST

NAMELIST#cloud_utils_nml

Methods

Included Modules

dc_types dc_message gridset constants_snowseaice constants gtool_historyauto timeset dc_iounit namelist_util

Public Instance methods

Subroutine :
xyz_TransCloudOneLayer(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_CloudCover(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyrr_OverlappedCloudTrans(0:imax-1, 1:jmax, 0:kmax, 0:kmax) :real(DP), intent(out)

[Source]

  subroutine CloudUtilsCalcOverlapCloudTrans( xyz_TransCloudOneLayer, xyz_CloudCover, xyrr_OverlappedCloudTrans )

    ! USE statements
    !

    ! ヒストリデータ出力
    ! History data output
    !
    use gtool_historyauto, only: HistoryAutoPut

    ! 時刻管理
    ! Time control
    !
    use timeset, only: TimeN, EndTime, TimesetClockStart, TimesetClockStop

!!$    use sort, only : SortQuick

    real(DP), intent(in ) :: xyz_TransCloudOneLayer   (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_CloudCover           (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(out) :: xyrr_OverlappedCloudTrans(0:imax-1, 1:jmax, 0:kmax, 0:kmax)


    real(DP) :: xyz_EffCloudCover           (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_CloudCoverSorted        (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_EffCloudCoverSorted     (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_TransCloudOneLayerSorted(0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: CloudCoverSortedCur
    real(DP) :: EffCloudCoverSortedCur
    real(DP) :: TransCloudOneLayerSortedCur
    integer  :: KInsPos
    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: kk
    integer  :: kkk



    ! 実行文 ; Executable statement
    !

    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    ! Cloud optical depth
    !

    select case ( IDCloudOverlapType )
    case ( IDCloudOverlapTypeRandom )

      xyz_EffCloudCover = xyz_CloudCover * ( 1.0_DP - xyz_TransCloudOneLayer )

      do k = 0, kmax
        kk = k
        xyrr_OverlappedCloudTrans(:,:,k,kk) = 1.0_DP
        do kk = k+1, kmax
          xyrr_OverlappedCloudTrans(:,:,k,kk) = xyrr_OverlappedCloudTrans(:,:,k,kk-1) * ( 1.0_DP - xyz_EffCloudCover(:,:,kk) )
        end do
      end do

      do k = 0, kmax
        do kk = 0, k-1
          xyrr_OverlappedCloudTrans(:,:,k,kk) = xyrr_OverlappedCloudTrans(:,:,kk,k)
        end do
      end do

    case ( IDCloudOverlapTypeMaxOverlap )

      ! see Chou et al. (2001)

      xyz_EffCloudCover = xyz_CloudCover * ( 1.0_DP - xyz_TransCloudOneLayer )


      ! Original method (computationally expensive, probably)
      !
!!$        do k = 0, kmax
!!$          kk = k
!!$          xyrr_OverlappedCloudTrans(:,:,k,kk) = 1.0_DP
!!$          do kk = k+1, kmax
!!$
!!$            xyz_CloudCoverSorted         = xyz_CloudCover
!!$            xyz_EffCloudCoverSorted      = xyz_EffCloudCover
!!$            xyz_TransCloudOneLayerSorted = xyz_TransCloudOneLayer
!!$
!!$            call SortQuick( imax, jmax, kk-k,             &
!!$              & xyz_CloudCoverSorted        (:,:,k+1:kk), &
!!$              & xyz_EffCloudCoverSorted     (:,:,k+1:kk), &
!!$              & xyz_TransCloudOneLayerSorted(:,:,k+1:kk)  &
!!$              & )
!!$
!!$            xyrr_OverlappedCloudTrans(:,:,k,kk) = 0.0_DP
!!$            do kkk = k+1, kk
!!$              xyrr_OverlappedCloudTrans(:,:,k,kk) =         &
!!$                & xyz_EffCloudCoverSorted(:,:,kkk)          &
!!$                & + xyrr_OverlappedCloudTrans(:,:,k,kk)     &
!!$                &   * xyz_TransCloudOneLayerSorted(:,:,kkk)
!!$            end do
!!$            xyrr_OverlappedCloudTrans(:,:,k,kk) = &
!!$              & 1.0_DP - xyrr_OverlappedCloudTrans(:,:,k,kk)
!!$
!!$          end do
!!$        end do


      ! Economical method (probably)
      !
      do k = 0, kmax

!!$          do kkk = 1, kmax
!!$              xyz_CloudCoverSorted(:,:,kkk) = real( kmax-kkk ) / real(kmax)
!!$!              xyz_CloudCoverSorted(:,:,kkk) = abs( 0.55d0 - real( kmax-kkk ) / real(kmax) )
!!$          end do
!!$          ! debug output
!!$          if ( k == 0 ) then
!!$            kk = kmax
!!$            do kkk = k+1, kk
!!$              write( 6, * ) kkk, xyz_CloudCoverSorted(0,jmax/2+1,kkk)
!!$            end do
!!$          end if

        xyz_CloudCoverSorted         = xyz_CloudCover
        xyz_EffCloudCoverSorted      = xyz_EffCloudCover
        xyz_TransCloudOneLayerSorted = xyz_TransCloudOneLayer

        kk = k
        xyrr_OverlappedCloudTrans(:,:,k,kk) = 1.0_DP
        do kk = k+1, kmax


          do j = 1, jmax
            do i = 0, imax-1

              ! xyz_CloudCoverSorted(i,j,kk) is inserved in an appropriate position.
              !
              KInsPos = kk
              loop : do kkk = k+1, kk-1

                if ( xyz_CloudCoverSorted(i,j,kk) < xyz_CloudCoverSorted(i,j,kkk) ) then
                  KInsPos = kkk
                  exit loop
                end if

              end do loop

              ! values are saved
              CloudCoverSortedCur         = xyz_CloudCoverSorted        (i,j,kk)
              EffCloudCoverSortedCur      = xyz_EffCloudCoverSorted     (i,j,kk)
              TransCloudOneLayerSortedCur = xyz_TransCloudOneLayerSorted(i,j,kk)

              ! values are shifted upward to empty an array at insert position
              do kkk = kk, KInsPos+1, -1
                xyz_CloudCoverSorted        (i,j,kkk) = xyz_CloudCoverSorted        (i,j,kkk-1)
                xyz_EffCloudCoverSorted     (i,j,kkk) = xyz_EffCloudCoverSorted     (i,j,kkk-1)
                xyz_TransCloudOneLayerSorted(i,j,kkk) = xyz_TransCloudOneLayerSorted(i,j,kkk-1)
              end do
              kkk = KInsPos
              xyz_CloudCoverSorted        (i,j,kkk) = CloudCoverSortedCur
              xyz_EffCloudCoverSorted     (i,j,kkk) = EffCloudCoverSortedCur
              xyz_TransCloudOneLayerSorted(i,j,kkk) = TransCloudOneLayerSortedCur

            end do
          end do


!!$            xyz_CloudCoverSorted         = xyz_CloudCover
!!$            do kkk = 1, kmax
!!$              xyz_CloudCoverSorted(:,:,kkk) = real( kmax-kkk ) / real(kmax)
!!$            end do
!!$            xyz_EffCloudCoverSorted      = xyz_EffCloudCover
!!$            xyz_TransCloudOneLayerSorted = xyz_TransCloudOneLayer
!!$
!!$            call SortQuick( imax, jmax, kk-k,             &
!!$              & xyz_CloudCoverSorted        (:,:,k+1:kk), &
!!$              & xyz_EffCloudCoverSorted     (:,:,k+1:kk), &
!!$              & xyz_TransCloudOneLayerSorted(:,:,k+1:kk)  &
!!$              & )


!!$            ! debug output
!!$            if ( ( k == 0 ) .and. ( kk == kmax-2 ) ) then
!!$              do kkk = k+1, kk
!!$                write( 6, * ) kkk, xyz_CloudCoverSorted(0,jmax/2+1,kkk)
!!$              end do
!!$              write( 6, * ) '-----'
!!$            end if


          xyrr_OverlappedCloudTrans(:,:,k,kk) = 0.0_DP
          do kkk = k+1, kk
            xyrr_OverlappedCloudTrans(:,:,k,kk) = xyz_EffCloudCoverSorted(:,:,kkk) + xyrr_OverlappedCloudTrans(:,:,k,kk) * xyz_TransCloudOneLayerSorted(:,:,kkk)
          end do
          xyrr_OverlappedCloudTrans(:,:,k,kk) = 1.0_DP - xyrr_OverlappedCloudTrans(:,:,k,kk)

        end do
      end do



      do k = 0, kmax
        do kk = 0, k-1
          xyrr_OverlappedCloudTrans(:,:,k,kk) = xyrr_OverlappedCloudTrans(:,:,kk,k)
        end do
      end do


    end select

    ! Output effective cloud cover
    !
!!$    call HistoryAutoPut( TimeN, 'EffCloudCover', &
!!$      & 1.0_DP - xyrr_OverlappedCloudTrans(:,:,0,kmax) )


  end subroutine CloudUtilsCalcOverlapCloudTrans
Subroutine :
xyz_Temp( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
xy_PRCP( 0:imax-1, 1:jmax ) :real(DP), intent(in )
xy_SurfRainFlux( 0:imax-1, 1:jmax ) :real(DP), intent(out)
xy_SurfSnowFlux( 0:imax-1, 1:jmax ) :real(DP), intent(out)

[Source]

  subroutine CloudUtilsCalcPRCPKeyLLTemp( xyz_Temp, xy_PRCP, xy_SurfRainFlux, xy_SurfSnowFlux )


    ! 雪と海氷の定数の設定
    ! Setting constants of snow and sea ice
    !
    use constants_snowseaice, only: TempCondWater


    real(DP), intent(in ) :: xyz_Temp       ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(in ) :: xy_PRCP        ( 0:imax-1, 1:jmax )
    real(DP), intent(out) :: xy_SurfRainFlux( 0:imax-1, 1:jmax )
    real(DP), intent(out) :: xy_SurfSnowFlux( 0:imax-1, 1:jmax )


    ! 作業変数
    ! Work variables
    !
    integer:: i               ! 経度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in longitude
    integer:: j               ! 緯度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in latitude


    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    if ( FlagSnow ) then

      do j = 1, jmax
        do i = 0, imax-1
          if ( xyz_Temp(i,j,1) > TempCondWater ) then
            xy_SurfRainFlux(i,j) = xy_PRCP(i,j)
            xy_SurfSnowFlux(i,j) = 0.0_DP
          else
            xy_SurfRainFlux(i,j) = 0.0_DP
            xy_SurfSnowFlux(i,j) = xy_PRCP(i,j)
          end if
        end do
      end do

    else

      xy_SurfRainFlux = xy_PRCP
      xy_SurfSnowFlux = 0.0_DP

    end if


  end subroutine CloudUtilsCalcPRCPKeyLLTemp
Subroutine :
xyr_Press( 0:imax-1, 1:jmax, 0:kmax ) :real(DP), intent(in )
xyz_Temp( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
xyz_DQH2OLiqDt( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
xy_SurfRainFlux( 0:imax-1, 1:jmax ) :real(DP), intent(out)
xy_SurfSnowFlux( 0:imax-1, 1:jmax ) :real(DP), intent(out)

[Source]

  subroutine CloudUtilsCalcPRCPKeyLLTemp3D( xyr_Press, xyz_Temp, xyz_DQH2OLiqDt, xy_SurfRainFlux, xy_SurfSnowFlux )

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav
                              ! $ g $ [m s-2].
                              ! 重力加速度.
                              ! Gravitational acceleration

    real(DP), intent(in ) :: xyr_Press       ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP), intent(in ) :: xyz_Temp        ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(in ) :: xyz_DQH2OLiqDt  ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(out) :: xy_SurfRainFlux ( 0:imax-1, 1:jmax )
    real(DP), intent(out) :: xy_SurfSnowFlux ( 0:imax-1, 1:jmax )


    ! 作業変数
    ! Work variables
    !
    real(DP) :: xy_PRCP( 0:imax-1, 1:jmax )

    integer  :: k


    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if



    xy_PRCP = 0.0d0
    do k = kmax, 1, -1
      xy_PRCP = xy_PRCP + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
    end do


    call CloudUtilsCalcPRCPKeyLLTemp( xyz_Temp, xy_PRCP, xy_SurfRainFlux, xy_SurfSnowFlux )


  end subroutine CloudUtilsCalcPRCPKeyLLTemp3D
Subroutine :
ArgFlagSnow :logical, intent(in)

This procedure input/output NAMELIST#cloud_utils_nml .

[Source]

  subroutine CloudUtilsInit( ArgFlagSnow )

    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid

    ! ヒストリデータ出力
    ! History data output
    !
    use gtool_historyauto, only: HistoryAutoAddVariable


    ! 宣言文 ; Declaration statements
    !

    logical, intent(in) :: ArgFlagSnow


    character(STRING) :: CloudWatIceFracType
    character(STRING) :: CloudOverlapType

    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号.
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT.
                              ! IOSTAT of NAMELIST read

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /cloud_utils_nml/ CloudWatIceFracType, CloudOverlapType, TempWatLim, TempIceLim
          !
          ! デフォルト値については初期化手続 "cloud_utils#CloudUtilsInit"
          ! のソースコードを参照のこと.
          !
          ! Refer to source codes in the initialization procedure
          ! "cloud_utils#CloudUtilsInit" for the default values.
          !

    ! 実行文 ; Executable statement
    !

    if ( cloud_utils_inited ) return


    FlagSnow = ArgFlagSnow


    ! デフォルト値の設定
    ! Default values settings
    !

    CloudWatIceFracType = "Lin"
    CloudOverlapType    = "Random"
!!$    CloudOverlapType    = "MaxOverlap"

!!$    TempWatLim          = 273.15_DP
!!$    TempIceLim          = 273.15_DP - 40.0_DP
    TempWatLim          = 0.0_DP
    TempIceLim          = 0.0_DP


    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = cloud_utils_nml, iostat = iostat_nml )             ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
    end if


    select case ( CloudOverlapType )
    case ( 'Random' )
      IDCloudOverlapType = IDCloudOverlapTypeRandom
    case ( 'MaxOverlap' )
      IDCloudOverlapType = IDCloudOverlapTypeMaxOverlap
    case default
      call MessageNotify( 'E', module_name, 'CloudOverlapType=<%c> is not supported.', c1 = trim(CloudOverlapType) )
    end select


    select case ( CloudWatIceFracType )
    case ( 'WatOnly' )
      IDCloudWatIceFracType = IDCloudWatIceFracTypeWatOnly
    case ( 'IceOnly' )
      IDCloudWatIceFracType = IDCloudWatIceFracTypeIceOnly
    case ( 'Lin' )
      IDCloudWatIceFracType = IDCloudWatIceFracTypeLin
      if ( TempWatLim < TempIceLim ) then
        call MessageNotify( 'E', module_name, 'TempWatLim must be greater than or equal to TempIceLim.' )
      end if
    case default
      call MessageNotify( 'E', module_name, 'CloudWatIceFracType=<%c> is not supported.', c1 = trim(CloudWatIceFracType) )
    end select


    ! Initialization of modules used in this module
    !


    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
!!$    call HistoryAutoAddVariable( 'EffCloudCover', &
!!$      & (/ 'lon ', 'lat ', 'time' /), &
!!$      & 'effective cloud cover', '1' )



    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'CloudOverlapType    = %c', c1 = trim(CloudOverlapType) )
    call MessageNotify( 'M', module_name, 'CloudWatIceFracType = %c', c1 = trim(CloudWatIceFracType) )
    call MessageNotify( 'M', module_name, 'TempWatLim          = %f', d = (/ TempWatLim /) )
    call MessageNotify( 'M', module_name, 'TempIceLim          = %f', d = (/ TempIceLim /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    cloud_utils_inited = .true.

  end subroutine CloudUtilsInit
Subroutine :
xyz_CloudCover(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_DelCloudOptDep(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)

[Source]

  subroutine CloudUtilsLocalizeCloud( xyz_CloudCover, xyz_DelCloudOptDep )

    ! USE statements
    !

    real(DP), intent(in   ) :: xyz_CloudCover    (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(inout) :: xyz_DelCloudOptDep(0:imax-1, 1:jmax, 1:kmax)


    ! 実行文 ; Executable statement
    !

    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    ! Cloud optical depth is scaled by considering cloud cover less than 1. 

    xyz_DelCloudOptDep = xyz_DelCloudOptDep / max( xyz_CloudCover, 1.0d-3 )


  end subroutine CloudUtilsLocalizeCloud
CloudUtilsWatFraction( Temp, WatFrac )
Subroutine :
Temp :real(DP), intent(in )
WatFrac :real(DP), intent(out)

Alias for CloudUtilsWatFraction0D

CloudUtilsWatFraction( xyz_Temp, xyz_WatFrac )
Subroutine :
xyz_Temp(:,:,:) :real(DP), intent(in )
xyz_WatFrac(:,:,:) :real(DP), intent(out)

Alias for CloudUtilsWatFraction3D

Private Instance methods

Subroutine :
xyz_CloudCover(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_DelCloudOptDep(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)

[Source]

  subroutine CloudUtilsSmearCloudOptDep( xyz_CloudCover, xyz_DelCloudOptDep )

    ! USE statements
    !

    real(DP), intent(in   ) :: xyz_CloudCover    (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(inout) :: xyz_DelCloudOptDep(0:imax-1, 1:jmax, 1:kmax)


    ! 実行文 ; Executable statement
    !

    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    ! Cloud optical depth is scaled by the way of Kiehl et al. (1994).

    xyz_DelCloudOptDep = xyz_DelCloudOptDep * xyz_CloudCover**1.5_DP


  end subroutine CloudUtilsSmearCloudOptDep
Subroutine :
Temp :real(DP), intent(in )
WatFrac :real(DP), intent(out)

[Source]

  subroutine CloudUtilsWatFraction0D( Temp, WatFrac )

    ! USE statements
    !

    real(DP), intent(in ) :: Temp
    real(DP), intent(out) :: WatFrac


    real(DP) :: xyz_Temp   (1,1,1)
    real(DP) :: xyz_WatFrac(1,1,1)

    ! 実行文 ; Executable statement
    !

    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    xyz_Temp = Temp
    call CloudUtilsWatFraction3D( xyz_Temp, xyz_WatFrac )
    WatFrac = xyz_WatFrac(1,1,1)


  end subroutine CloudUtilsWatFraction0D
Subroutine :
xyz_Temp(:,:,:) :real(DP), intent(in )
xyz_WatFrac(:,:,:) :real(DP), intent(out)

[Source]

  subroutine CloudUtilsWatFraction3D( xyz_Temp, xyz_WatFrac )

    ! USE statements
    !

!!$    real(DP), intent(in ) :: xyz_Temp   (0:imax-1, 1:jmax, 1:kmax)
!!$    real(DP), intent(out) :: xyz_WatFrac(0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_Temp   (:,:,:)
    real(DP), intent(out) :: xyz_WatFrac(:,:,:)


    ! 実行文 ; Executable statement
    !

    ! 初期化確認
    ! Initialization check
    !
    if ( .not. cloud_utils_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    if ( FlagSnow ) then

      select case ( IDCloudWatIceFracType )
      case ( IDCloudWatIceFracTypeWatOnly  )
        xyz_WatFrac = 1.0_DP
      case ( IDCloudWatIceFracTypeIceOnly  )
        xyz_WatFrac = 0.0_DP
      case ( IDCloudWatIceFracTypeLin      )
        if ( TempWatLim == TempIceLim ) then
          xyz_WatFrac = ( sign( 1.0_DP, xyz_Temp - TempWatLim ) + 1.0_DP ) / 2.0_DP
        else
          xyz_WatFrac = ( 1.0_DP - 0.0_DP ) / ( TempWatLim - TempIceLim ) * ( xyz_Temp - TempIceLim )
          xyz_WatFrac = min( xyz_WatFrac, 1.0_DP )
          xyz_WatFrac = max( xyz_WatFrac, 0.0_DP )
        end if
      end select

    else

      xyz_WatFrac = 1.0_DP

    end if


  end subroutine CloudUtilsWatFraction3D
FlagSnow
Variable :
FlagSnow :logical , save
: A flag for snow
IDCloudOverlapType
Variable :
IDCloudOverlapType :integer , save
IDCloudOverlapTypeMaxOverlap
Constant :
IDCloudOverlapTypeMaxOverlap = 2 :integer , parameter
IDCloudOverlapTypeRandom
Constant :
IDCloudOverlapTypeRandom = 1 :integer , parameter
IDCloudWatIceFracType
Variable :
IDCloudWatIceFracType :integer , save
IDCloudWatIceFracTypeIceOnly
Constant :
IDCloudWatIceFracTypeIceOnly = 3 :integer , parameter
IDCloudWatIceFracTypeLin
Constant :
IDCloudWatIceFracTypeLin = 1 :integer , parameter
IDCloudWatIceFracTypeWatOnly
Constant :
IDCloudWatIceFracTypeWatOnly = 2 :integer , parameter
TempIceLim
Variable :
TempIceLim :real(DP), save
TempWatLim
Variable :
TempWatLim :real(DP), save
cloud_utils_inited
Variable :
cloud_utils_inited = .false. :logical, save
: 初期設定フラグ. Initialization flag
module_name
Constant :
module_name = ‘cloud_utils :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: $’ // ’$Id: cloud_utils.f90,v 1.5 2014/05/07 10:02:56 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version