Class set_cloud
In: radiation/set_cloud.F90

雲量の設定

Set cloud amount

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#set_cloud_nml

Methods

Included Modules

dc_types gridset constants radiation_two_stream_app gtool_historyauto timeset dc_iounit namelist_util dc_message

Public Instance methods

Subroutine :
xyz_CloudDelOptDep( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(out)

[Source]

  subroutine SetCloudLW( xyz_CloudDelOptDep )

    ! USE statements
    !

    !
    ! Physical constants settings
    !
    use constants, only: Grav, PI       ! $ \pi $ .
                                  ! Circular constant

    use radiation_two_stream_app, only: RadiationTwoStreamApp

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

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

    real(DP), intent(out) :: xyz_CloudDelOptDep( 0:imax-1, 1:jmax, 1:kmax )


    real(DP)            :: xyz_DelLWP     ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP)            :: xyr_CloudOptDep( 0:imax-1, 1:jmax, 0:kmax )

    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: l


    ! 初期化
    ! Initialization
    !
    if ( .not. set_cloud_inited ) call SetCloudInit



    ! Cloud optical depth
    !
    if ( .not. FlagCloud ) then

      xyz_DelLWP         = 0.0d0
      xyz_CloudDelOptDep = 0.0d0

    else

!!$      write( 6, * ) 'Cloud optical depth is calculated by precipitation rate.'
      do k = 1, kmax
        do j = 1, jmax
          do i = 0, imax-1
            if ( ( xyz_DelPRCPCumSave(i,j,k) + xyz_DelPRCPLSCSave(i,j,k) ) > 0.0d0 ) then
              xyz_DelLWP(i,j,k) = ( xyz_DelPRCPCumSave(i,j,k) + xyz_DelPRCPLSCSave(i,j,k) ) * FactorCondToCloudLW
            else
              xyz_DelLWP(i,j,k) = 0.0d0
            end if
          end do
        end do
      end do


      xyz_CloudDelOptDep = xyz_DelLWP

    end if


    ! ヒストリデータ出力
    ! History data output
    !
    xyr_CloudOptDep(:,:,kmax) = 0.0d0
    do k = kmax-1, 0, -1
      xyr_CloudOptDep(:,:,k) = xyr_CloudOptDep(:,:,k+1) + xyz_CloudDelOptDep(:,:,k+1)
    end do
    call HistoryAutoPut( TimeN, 'CODLW'   , xyr_CloudOptDep(:,:,0) )
    call HistoryAutoPut( TimeN, 'DelCODLW', xyz_CloudDelOptDep     )


  end subroutine SetCloudLW
Subroutine :
xyz_DPDt( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in)

[Source]

  subroutine SetCloudRegistDPDt( xyz_DPDt )

    real(DP), intent(in) :: xyz_DPDt( 0:imax-1, 1:jmax, 1:kmax )


    ! 初期化
    ! Initialization
    !
    if ( .not. set_cloud_inited ) call SetCloudInit


    xyz_DPDtSave  = xyz_DPDt


    dpdt_registered = .true.

  end subroutine SetCloudRegistDPDt
Subroutine :
xyz_DelRain( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in)
:
!$ real(DP), intent(in) :xyz_DTempDt( 0:imax-1, 1:jmax, 1:kmax )

[Source]

  subroutine SetCloudRegistPRCPCum( xyz_DelRain )

!!$    real(DP), intent(in) :: xy_Rain    ( 0:imax-1, 1:jmax )
    real(DP), intent(in) :: xyz_DelRain( 0:imax-1, 1:jmax, 1:kmax )
!!$    real(DP), intent(in) :: xyz_DTempDt( 0:imax-1, 1:jmax, 1:kmax )


    ! 初期化
    ! Initialization
    !
    if ( .not. set_cloud_inited ) call SetCloudInit


!!$    xy_PRCPSave         = xy_Rain
    xyz_DelPRCPCumSave  = xyz_DelRain
!!$    xyz_DTempDtCumSave  = xyz_DTempDt


    prcpcum_registered = .true.

  end subroutine SetCloudRegistPRCPCum
Subroutine :
xyz_DelRain( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in)
:
!$ real(DP), intent(in) :xyz_DTempDt( 0:imax-1, 1:jmax, 1:kmax )

[Source]

  subroutine SetCloudRegistPRCPLSC( xyz_DelRain )

!!$    real(DP), intent(in) :: xy_Rain    ( 0:imax-1, 1:jmax )
    real(DP), intent(in) :: xyz_DelRain( 0:imax-1, 1:jmax, 1:kmax )
!!$    real(DP), intent(in) :: xyz_DTempDt( 0:imax-1, 1:jmax, 1:kmax )


    ! 初期化
    ! Initialization
    !
    if ( .not. set_cloud_inited ) call SetCloudInit


!!$    xy_PRCPSave         = xy_Rain
    xyz_DelPRCPLSCSave  = xyz_DelRain
!!$    xyz_DTempDtCumSave  = xyz_DTempDt


    prcplsc_registered = .true.

  end subroutine SetCloudRegistPRCPLSC
Subroutine :
xyz_Press( 0:imax-1, 1:jmax, 1:kmax ) :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_QVap( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
: $ q $ . 混合比. Mass mixing ratio of constituents (1)
xyz_Height( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
xyz_CloudDelOptDep( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(out)

[Source]

  subroutine SetCloudSW( xyz_Press, xyr_Press, xyz_Temp, xyz_QVap, xyz_Height, xyz_CloudDelOptDep )

    ! USE statements
    !

    !
    ! Physical constants settings
    !
    use constants, only: Grav, PI       ! $ \pi $ .
                                  ! Circular constant

    use radiation_two_stream_app, only: RadiationTwoStreamApp

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

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

    real(DP), intent(in ) :: xyz_Press         ( 0:imax-1, 1:jmax, 1:kmax )
    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_QVap          ( 0:imax-1, 1:jmax, 1:kmax )
                              ! $ q $ .   混合比. Mass mixing ratio of constituents (1)
    real(DP), intent(in ) :: xyz_Height        ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(out) :: xyz_CloudDelOptDep( 0:imax-1, 1:jmax, 1:kmax )



    real(DP)            :: xyz_DelLWP     ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP)            :: xyr_CloudOptDep( 0:imax-1, 1:jmax, 0:kmax )

    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: l


    ! 初期化
    ! Initialization
    !
    if ( .not. set_cloud_inited ) call SetCloudInit



    ! Cloud optical depth
    !
    if ( .not. FlagCloud ) then

      xyz_DelLWP         = 0.0d0
      xyz_CloudDelOptDep = 0.0d0

    else

!!$    call RadiationDcpamSWEV1CalcLWP( &
!!$      & xyz_Press, xyr_Press, xyz_Temp, xyz_QVap, xyz_Height, &
!!$      & xyz_DelLWP   &
!!$      & )
!!$
!!$      ! This calculation is based Equation (1) of Slingo (1989) and values of 
!!$      ! ai = 2.7d-2 m2 g-1, bi = 1.3 micron m2 g-1, and r_{eff} = 10 micron.
!!$      !
!!$      xyz_CloudDelOptDep = xyz_DelLWP * ( 2.7d-2 * 1.0d3 + 1.3d0 * 1.0d3 / 10.0d0 )


!!$      write( 6, * ) 'Cloud optical depth is calculated by precipitation rate.'
      do k = 1, kmax
        do j = 1, jmax
          do i = 0, imax-1
            if ( ( xyz_DelPRCPCumSave(i,j,k) + xyz_DelPRCPLSCSave(i,j,k) ) > 0.0d0 ) then
              xyz_DelLWP(i,j,k) = ( xyz_DelPRCPCumSave(i,j,k) + xyz_DelPRCPLSCSave(i,j,k) ) * FactorCondToCloudSW
            else
              xyz_DelLWP(i,j,k) = 0.0d0
            end if
          end do
        end do
      end do


      xyz_CloudDelOptDep = xyz_DelLWP


    end if


    ! ヒストリデータ出力
    ! History data output
    !
    xyr_CloudOptDep(:,:,kmax) = 0.0d0
    do k = kmax-1, 0, -1
      xyr_CloudOptDep(:,:,k) = xyr_CloudOptDep(:,:,k+1) + xyz_CloudDelOptDep(:,:,k+1)
    end do
    call HistoryAutoPut( TimeN, 'CODSW'   , xyr_CloudOptDep(:,:,0) )
    call HistoryAutoPut( TimeN, 'DelCODSW', xyz_CloudDelOptDep     )


  end subroutine SetCloudSW
set_cloud_inited
Variable :
set_cloud_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

FactorCondToCloudLW
Variable :
FactorCondToCloudLW :real(DP), save
: A constant to be used to convert liquid water production rate in a layer to cloud optical depth in long wave radiation.
FactorCondToCloudSW
Variable :
FactorCondToCloudSW :real(DP), save
: A constant to be used to convert liquid water production rate in a layer to cloud optical depth in short wave radiation.
FlagCloud
Variable :
FlagCloud :logical , save
: A flag for cloud set. If FlagCloud is true, the effect of cloud is considered. It should be noticed that in principle FlagCloud = .false. is equivalent to FactorCondToCloudSW = FactorCondToCloudLW = 0. So, this variable is verbose.
Subroutine :
xyz_Press( 0:imax-1, 1:jmax, 1:kmax ) :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_QVap( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
xyz_Height( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(in )
xyz_DelLWP( 0:imax-1, 1:jmax, 1:kmax ) :real(DP), intent(out)

[Source]

  subroutine SetCloudCalcLWP( xyz_Press, xyr_Press, xyz_Temp, xyz_QVap, xyz_Height, xyz_DelLWP )

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: GasRUniv, Grav, GasRDry, CpDry, GasRWet, LatentHeat, EpsV
                              ! $ \epsilon_v $ .
                              ! 水蒸気分子量比.
                              ! Molecular weight of water vapor

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

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


    real(DP), intent(in ) :: xyz_Press ( 0:imax-1, 1:jmax, 1:kmax )
    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_QVap  ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(in ) :: xyz_Height( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(out) :: xyz_DelLWP( 0:imax-1, 1:jmax, 1:kmax )


    integer  :: xy_NumCloudLayer      ( 0:imax-1, 1:jmax )
    real(DP) :: xy_CloudCoverTot      ( 0:imax-1, 1:jmax )
    real(DP) :: xy_CloudCoverEachLayer( 0:imax-1, 1:jmax )
    real(DP) :: xyz_CloudCover        ( 0:imax-1, 1:jmax, 1:kmax )

    real(DP):: xyz_QVapSat            ( 0:imax-1, 1:jmax, 1:kmax )
                              ! 飽和比湿.
                              ! Saturation specific humidity.

    real(DP):: xyz_RH                 ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP):: xyz_RHCor              ( 0:imax-1, 1:jmax, 1:kmax )

    real(DP), parameter :: CloudLiqWatDen0 = 0.18d-3
    real(DP)            :: xy_ColWatVapor           ( 0:imax-1, 1:jmax )
    real(DP)            :: xy_CloudLiqWatScaleHeight( 0:imax-1, 1:jmax )
    real(DP)            :: xyz_CloudLiqWatDen       ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP)            :: xyz_DelCloudLiqWat       ( 0:imax-1, 1:jmax, 1:kmax )

    integer  :: i
    integer  :: j
    integer  :: k


    ! 飽和比湿計算のための文関数定義
    ! Declaration of statement function for
    !   calculation of saturation specific humidity
    !
#ifdef LIB_SATURATE_NHA1992
#include "../saturate/saturate_nha1992_sf.f90"
    EpsVSF     = EpsV
    GasRUnivSF = GasRUniv
#elif LIB_SATURATE_T1930
#include "../saturate/saturate_t1930_sf.f90"
    EpsVSF    = EpsV
    LatHeatSF = LatentHeat
    GasRWetSF = GasRWet
#else
#include "../saturate/saturate_t1930_sf.f90"
    EpsVSF    = EpsV
    LatHeatSF = LatentHeat
    GasRWetSF = GasRWet
#endif


    ! 初期化
    ! Initialization
    !
    if ( .not. set_cloud_inited ) call SetCloudInit


    ! 飽和比湿計算
    ! Calculate saturation specific humidity
    !

    ! Nakajima et al. (1992) を用いた飽和比湿の計算
    ! Calculate saturation specific humidity with Nakajima et al. (1992)
    !
    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          ! CalcQVapSatSF は文関数. (実行文の直前で定義)
          ! "CalcQVapSatSF" is statement function and
          !   is declared just before executable statement.
          !
          xyz_QVapSat(i,j,k) = CalcQVapSatSF( xyz_Temp(i,j,k), xyz_Press(i,j,k) )
        end do
      end do
    end do

    xyz_RH = xyz_QVap / xyz_QVapSat


    xy_NumCloudLayer = 0
    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          if ( xyz_DTempDtCumSave(i,j,k) > 0.0d0 ) then
            xy_NumCloudLayer = xy_NumCloudLayer + 1
          end if
        end do
      end do
    end do

    xy_CloudCoverTot = 0.0d0
    do j = 1, jmax
      do i = 0, imax-1
        if ( xy_PRCPSave(i,j) >= 0.0d0 ) then
          xy_CloudCoverTot(i,j) = 0.20d0 + 0.125d0 * log( 1.0d0 + xy_PRCPSave(i,j) * 86400.0d0 / 1.0d3 * 1.0d3 )
        end if
      end do
    end do
    do j = 1, jmax
      do i = 0, imax-1
        if ( xy_CloudCoverTot(i,j) > 0.8d0 ) then
          xy_CloudCoverTot(i,j) = 0.8d0
        end if
      end do
    end do
    do j = 1, jmax
      do i = 0, imax-1
        if ( xy_NumCloudLayer(i,j) == 0 ) then
          xy_CloudCoverEachLayer(i,j) = 0.0d0
        else
          xy_CloudCoverEachLayer(i,j) = 1.0d0 - ( 1.0d0 - xy_CloudCoverTot(i,j) )**(1.0d0/xy_NumCloudLayer(i,j))
        end if
      end do
    end do

    do k = 1, kmax
      xyz_RHCor(:,:,k) = ( xyz_RH(:,:,k) - xy_CloudCoverEachLayer ) / ( 1.0d0 - xy_CloudCoverEachLayer )
    end do


    xyz_CloudCover = 0.0d0

!    ! Frontal and tropical low-cloud
!    !
!    do k = 1, kmax
!      do j = 1, jmax
!        do i = 0, imax-1
!          if ( xyz_Press(i,j,k) > 750.0d2 .and. xyz_DPDtSave(i,j,k) < 0.0d0 ) then
!            xyz_CloudCover(i,j,k) = ( ( xyz_RHCor(i,j,k) - 0.9d0 ) / 0.1d0 )**2
!          end if
!        end do
!      end do
!    end do


    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          if ( ( xyz_RHCor(i,j,k) > 0.6d0 ) .and. ( xyz_DTempDtCumSave(i,j,k) > 0.0d0 ) ) then
            xyz_CloudCover(i,j,k) = xy_CloudCoverEachLayer(i,j)
          else
            xyz_CloudCover(i,j,k) = 0.0d0
          end if
        end do
      end do
    end do




    xy_ColWatVapor = 0.0d0
    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          xy_ColWatVapor = xy_ColWatVapor + xyz_QVap(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
        end do
      end do
    end do

    xy_CloudLiqWatScaleHeight = 700.0d0 * log( 1.0d0 + 1.0d0 * xy_ColWatVapor ) + 1.0d0

    do k = 1, kmax
      xyz_CloudLiqWatDen(:,:,k) = CloudLiqWatDen0 * exp( - xyz_Height(:,:,k) / xy_CloudLiqWatScaleHeight )
    end do


    do k = 1, kmax
      xyz_DelCloudLiqWat(:,:,k) = xyz_CloudLiqWatDen(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / ( xyz_Press(:,:,k) / ( GasRDry * xyz_Temp(:,:,k) ) * Grav )
    end do

    xyz_DelLWP = xyz_CloudCover * xyz_DelCloudLiqWat


    call HistoryAutoPut( TimeN, 'CloudCover', xyz_CloudCover )


  end subroutine SetCloudCalcLWP
Subroutine :

This procedure input/output NAMELIST#set_cloud_nml .

[Source]

  subroutine SetCloudInit

    ! ファイル入出力補助
    ! 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

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify


    ! 宣言文 ; Declaration statements
    !

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

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /set_cloud_nml/ FactorCondToCloudSW, FactorCondToCloudLW, FlagCloud
          !
          ! デフォルト値については初期化手続 "set_cloud#setCloudInit"
          ! のソースコードを参照のこと.
          !
          ! Refer to source codes in the initialization procedure
          ! "set_cloud#SetCloudInit" for the default values.
          !


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

    FactorCondToCloudSW = 1.0d5
    FactorCondToCloudLW = 1.0d5
    FlagCloud           = .true.


    ! 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 = set_cloud_nml, iostat = iostat_nml )             ! (out)
      close( unit_nml )

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


    allocate( xy_PRCPSave       ( 0:imax-1, 1:jmax ) )
    allocate( xyz_DelPRCPCumSave( 0:imax-1, 1:jmax, 1:kmax ) )
    allocate( xyz_DelPRCPLSCSave( 0:imax-1, 1:jmax, 1:kmax ) )
    allocate( xyz_DTempDtCumSave( 0:imax-1, 1:jmax, 1:kmax ) )
    allocate( xyz_DPDtSave      ( 0:imax-1, 1:jmax, 1:kmax ) )


    xy_PRCPSave        = 0.0d0
    xyz_DelPRCPCumSave = 0.0d0
    xyz_DelPRCPLSCSave = 0.0d0
    xyz_DTempDtCumSave = 0.0d0
    xyz_DPDtSave       = 0.0d0


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

    call HistoryAutoAddVariable( 'DelCODSW', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'cloud optical depth', '1' )

    call HistoryAutoAddVariable( 'CODLW', (/ 'lon ', 'lat ', 'time' /), 'cloud optical depth', '1' )

    call HistoryAutoAddVariable( 'DelCODLW', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'cloud optical depth', '1' )

    call HistoryAutoAddVariable( 'CloudCover', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'cloud cover', '1' )



    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'FactorCondToCloudSW = %f', d = (/ FactorCondToCloudSW /) )
    call MessageNotify( 'M', module_name, 'FactorCondToCloudLW = %f', d = (/ FactorCondToCloudLW /) )
    call MessageNotify( 'M', module_name, 'FlagCloud           = %b', l = (/ FlagCloud /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )


    set_cloud_inited = .true.

  end subroutine SetCloudInit
dpdt_registered
Variable :
dpdt_registered = .false. :logical , save
module_name
Constant :
module_name = ‘set_cloud :character(*), parameter
: モジュールの名称. Module name
prcpcum_registered
Variable :
prcpcum_registered = .false. :logical , save
prcplsc_registered
Variable :
prcplsc_registered = .false. :logical , save
saturate_scheme
Constant :
saturate_scheme = ifdef LIB_SATURATE_NHA1992 elif LIB_SATURATE_T1930 else endif :character(*), parameter
version
Constant :
version = ’$Name: dcpam5-20100224 $’ // ’$Id: set_cloud.F90,v 1.1 2010-01-11 01:28:12 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version
xy_PRCPSave
Variable :
xy_PRCPSave(:,:) :real(DP), allocatable
xyz_DPDtSave
Variable :
xyz_DPDtSave(:,:,:) :real(DP), allocatable
xyz_DTempDtCumSave
Variable :
xyz_DTempDtCumSave(:,:,:) :real(DP), allocatable
xyz_DelPRCPCumSave
Variable :
xyz_DelPRCPCumSave(:,:,:) :real(DP), allocatable
xyz_DelPRCPLSCSave
Variable :
xyz_DelPRCPLSCSave(:,:,:) :real(DP), allocatable

[Validate]