Class moist_conv_adjust
In: cumulus/moist_conv_adjust.f90

湿潤対流調節

Moist convective adjustment

Note that Japanese and English are described in parallel.

湿潤対流調節スキームにより, 温度と比湿を調節.

Moist convective adjustment was originally proposed by Manabe et al. (1965). But, the algorithm used in this routine seems to be different from that described by Manabe et al. (1965).

Adjust temperature and specific humidity by convective adjustment scheme.

Procedures List

MoistConvAdjust :温度と比湿の調節
————— :————
MoistConvAdjust :Adjust temperature and specific humidity

NAMELIST

NAMELIST#moist_conv_adjust_nml

Methods

Included Modules

gridset dc_types namelist_util dc_message constants timeset gtool_historyauto saturate dc_iounit dc_string

Public Instance methods

Subroutine :
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: $ T $ . 温度. Temperature
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: $ q $ . 比湿. Specific humidity
!$ real(DP), intent(inout):xy_Rain (0:imax-1, 1:jmax)

!$ ! 降水量. !$ ! Precipitation

xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: 温度変化率. Temperature tendency
xyz_DQVapDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: 比湿変化率. Specific humidity tendency
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ p $ . 気圧 (整数レベル). Air pressure (full level)
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
xyz_DQH2OLiqDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)

湿潤対流調節スキームにより, 温度と比湿を調節.

Adjust temperature and specific humidity by moist convective adjustment

[Source]

  subroutine MoistConvAdjust( xyz_Temp, xyz_QVap, xyz_DTempDt, xyz_DQVapDt, xyz_Press, xyr_Press, xyz_DQH2OLiqDt )
    !
    ! 湿潤対流調節スキームにより, 温度と比湿を調節. 
    !
    ! Adjust temperature and specific humidity by moist convective adjustment
    !

    ! モジュール引用 ; USE statements
    !

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, GasRDry, CpDry, LatentHeat
                              ! $ L $ [J kg-1] . 
                              ! 凝結の潜熱. 
                              ! Latent heat of condensation

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

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

    ! 飽和比湿の算出
    ! Evaluate saturation specific humidity
    !
    use saturate, only: xyz_CalcQVapSat, xy_CalcDQVapSatDTemp


    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(inout):: xyz_Temp (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(inout):: xyz_QVap (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .     比湿. Specific humidity
!!$    real(DP), intent(inout):: xy_Rain (0:imax-1, 1:jmax)
!!$                              ! 降水量. 
!!$                              ! Precipitation
    real(DP), intent(inout):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax)
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP), intent(inout):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax)
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP), intent(in):: xyz_Press (0:imax-1, 1:jmax, 1:kmax)
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(out):: xyz_DQH2OLiqDt(0:imax-1, 1:jmax, 1:kmax)

    ! 作業変数
    ! Work variables
    !
    real(DP):: xy_RainCumulus (0:imax-1, 1:jmax)
                              ! 降水量. 
                              ! Precipitation
    real(DP):: xyz_DTempDtCumulus (0:imax-1, 1:jmax, 1:kmax)
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP):: xyz_DQVapDtCumulus (0:imax-1, 1:jmax, 1:kmax)
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP):: xyz_QVapB (0:imax-1, 1:jmax, 1:kmax)
                              ! 調節前の比湿. 
                              ! Specific humidity before adjustment
    real(DP):: xyz_TempB (0:imax-1, 1:jmax, 1:kmax)
                              ! 調節前の温度. 
                              ! Temperature before adjustment
    logical:: xy_Adjust (0:imax-1, 1:jmax)
                              ! 今回調節されたか否か?. 
                              ! Whether it was adjusted this time or not?
    logical:: xy_AdjustB (0:imax-1, 1:jmax)
                              ! 前回調節されたか否か?. 
                              ! Whether it was adjusted last time or not?
    real(DP):: xyz_DelPress(0:imax-1, 1:jmax, 1:kmax)
                              ! $ \Delta p $
                              !
    real(DP):: xyz_QVapSat (0:imax-1, 1:jmax, 1:kmax)
                              ! 飽和比湿. 
                              ! Saturation specific humidity. 
    real(DP):: xyr_ConvAdjustFactor(0:imax-1, 1:jmax, 0:kmax)
                              ! $ \frac{1}{2} \frac{ R }{Cp} 
                              !   \frac{ p_{k} - p_{k+1} }{ p_{k+1/2} } $

    real(DP):: TempEquivToExcEne
                              ! Temperature equivalent to the excess moist static energy
                              ! (Moist static energy difference devided by specific heat)

    real(DP):: DelQ
    real(DP):: DelTempUppLev
                              ! k+1 番目の層における調節による温度の変化量. 
                              ! Temperature variation by adjustment at k+1 level
    real(DP):: DelTempLowLev
                              ! k 番目の層における調節による温度の変化量. 
                              ! Temperature variation by adjustment at k level
    real(DP):: DQVapSatDTempUppLev
                              ! $ \DD{q^{*}} (k+1)}{T} $
    real(DP):: DQVapSatDTempLowLev
                              ! $ \DD{q^{*}} (k)}{T} $
    real(DP):: GamUppLev
                              ! $ \gamma_{k+1} = \frac{L}{C_p} \DP{q^{*}}{T}_{k+1} $
    real(DP):: GamLowLev
                              ! $ \gamma_{k}   = \frac{L}{C_p} \DP{q^{*}}{T}_{k} $
    logical:: Adjust
                              ! 今回全領域において一度でも調節されたか否か?. 
                              ! Whether it was adjusted even once in global 
                              ! this time or not?

    real(DP):: TempLowLevBefAdj ! Variables for check routine
    real(DP):: TempUppLevBefAdj
    real(DP):: QVapLowLevBefAdj
    real(DP):: QVapUppLevBefAdj

    integer:: i               ! 経度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in longitude
    integer:: j               ! 緯度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in latitude
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    integer:: itr             ! イテレーション方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in iteration direction

    real(DP):: xy_DQVapSatDTempUppLev(0:imax-1, 1:jmax)
    real(DP):: xy_DQVapSatDTempLowLev(0:imax-1, 1:jmax)

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

    real(DP) :: xy_NegDDelLWDt   (0:imax-1, 1:jmax)
    real(DP) :: xyz_DDelLWDtCCPLV(0:imax-1, 1:jmax, 1:kmax)


    ! 実行文 ; Executable statement
    !

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

!!$    if ( .not. FlagUse ) return


    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )


    ! 調節前 "QVap", "Temp" の保存
    ! Store "QVap", "Temp" before adjustment
    !
    xyz_QVapB = xyz_QVap
    xyz_TempB = xyz_Temp

    ! 飽和比湿計算
    ! Calculate saturation specific humidity
    !
    xyz_QVapSat = xyz_CalcQVapSat( xyz_Temp, xyz_Press )

    ! Calculate some values used for moist convective adjustment
    !

    do k = 1, kmax
      xyz_DelPress(:,:,k) = xyr_Press(:,:,k-1) - xyr_Press(:,:,k)
    end do

    ! \frac{1}{2} \frac{ R }{Cp} \frac{ p_{k} - p_{k+1} }{ p_{k+1/2} }
    !
    !   The value at k = 0 is not used.
    k = 0
    xyr_ConvAdjustFactor(:,:,k) = 0.0d0
    !
    do k = 1, kmax-1
      xyr_ConvAdjustFactor(:,:,k) = GasRDry / CpDry * ( xyz_Press(:,:,k) - xyz_Press(:,:,k+1) ) / xyr_Press(:,:,k) / 2.0_DP
    end do
    !   The value at k = kmax is not used.
    k = kmax
    xyr_ConvAdjustFactor(:,:,k) = 0.0d0


    ! 調節
    ! Adjustment
    !
    xy_AdjustB = .true.

    ! 繰り返し
    ! Iteration
    !
    do itr = 1, ItrtMax
      xy_Adjust = .false.

      do k = 1, kmax-1

        xy_DQVapSatDTempUppLev = xy_CalcDQVapSatDTemp( xyz_Temp(:,:,k+1), xyz_QVapSat(:,:,k+1) )
        xy_DQVapSatDTempLowLev = xy_CalcDQVapSatDTemp( xyz_Temp(:,:,k  ), xyz_QVapSat(:,:,k  ) )

        do j = 1, jmax
          do i = 0, imax-1
            if ( xy_AdjustB(i,j) ) then

              ! Temperature equivalent to the excess moist static energy
              ! (Moist static energy difference devided by specific heat)
              !
              TempEquivToExcEne = xyz_Temp(i,j,k) - xyz_Temp(i,j,k+1) + LatentHeat / CpDry * ( xyz_QVapSat(i,j,k) - xyz_QVapSat(i,j,k+1) ) - xyr_ConvAdjustFactor(i,j,k) * ( xyz_Temp(i,j,k) + xyz_Temp(i,j,k+1) )

              ! Check vertical gradient of moist static energy
              !
              if ( TempEquivToExcEne > AdjustCriterion(itr) ) then

                ! Check relative humidity
                !
                if ( ( xyz_QVap(i,j,k+1) / xyz_QVapSat(i,j,k+1) >= CrtlRH ) .and. ( xyz_QVap(i,j,k  ) / xyz_QVapSat(i,j,k  ) >= CrtlRH ) ) then

                  DelQ = xyz_DelPress(i,j,k  ) * ( xyz_QVap(i,j,k  ) - xyz_QVapSat(i,j,k  ) ) + xyz_DelPress(i,j,k+1) * ( xyz_QVap(i,j,k+1) - xyz_QVapSat(i,j,k+1) )

                  DQVapSatDTempUppLev = xy_DQVapSatDTempUppLev(i,j)
                  DQVapSatDTempLowLev = xy_DQVapSatDTempLowLev(i,j)

                  GamUppLev = LatentHeat / CpDry * DQVapSatDTempUppLev
                  GamLowLev = LatentHeat / CpDry * DQVapSatDTempLowLev

                  DelTempUppLev = ( xyz_DelPress(i,j,k) * ( 1.0d0 + GamLowLev ) * TempEquivToExcEne + ( 1.0d0 + GamLowLev - xyr_ConvAdjustFactor(i,j,k) ) * LatentHeat / CpDry * DelQ ) / ( xyr_ConvAdjustFactor(i,j,k) * ( xyz_DelPress(i,j,k  ) * ( 1.0d0 + GamLowLev ) - xyz_DelPress(i,j,k+1) * ( 1.0d0 + GamUppLev ) ) + ( 1.0d0 + GamLowLev ) * ( 1.0d0 + GamUppLev ) * ( xyz_DelPress(i,j,k) + xyz_DelPress(i,j,k+1) ) )

                  DelTempLowLev = ( LatentHeat / CpDry * DelQ - xyz_DelPress(i,j,k+1) * ( 1.0d0 + GamUppLev ) * DelTempUppLev ) / ( ( 1.0d0 + GamLowLev ) * xyz_DelPress(i,j,k) )


                  !=========
                  ! check routine
                  !---------
!!$                  TempLowLevBefAdj = xyz_Temp(i,j,k  )
!!$                  TempUppLevBefAdj = xyz_Temp(i,j,k+1)
!!$                  QVapLowLevBefAdj = xyz_QVap(i,j,k  )
!!$                  QVapUppLevBefAdj = xyz_QVap(i,j,k+1)
                  !=========


                  ! 温度の調節
                  ! Adjust temperature
                  !
                  xyz_Temp(i,j,k  ) = xyz_Temp(i,j,k  ) + DelTempLowLev
                  xyz_Temp(i,j,k+1) = xyz_Temp(i,j,k+1) + DelTempUppLev

                  ! 比湿の調節
                  ! Adjust specific humidity
                  !
                  xyz_QVap(i,j,k  ) = xyz_QVapSat(i,j,k  ) + DQVapSatDTempLowLev * DelTempLowLev
                  xyz_QVap(i,j,k+1) = xyz_QVapSat(i,j,k+1) + DQVapSatDTempUppLev * DelTempUppLev


                  !=========
                  ! check routine
                  !---------
!!$                  write( 6, * ) '====='
!!$                  write( 6, * ) xyz_Temp(i,j,k), xyz_Temp(i,j,k+1), xyz_QVap(i,j,k), xyz_QVap(i,j,k+1)
!!$                  write( 6, * ) DelTempLowLev, DelTempUppLev
!!$                  write( 6, * ) 'Energy difference before and after adjustment and each energy'
!!$                  write( 6, * ) &
!!$                    &   ( CpDry * TempLowLevBefAdj  + LatentHeat * QVapLowLevBefAdj )  &
!!$                    &     * xyz_DelPress(i,j,k  ) / Grav                               &
!!$                    & + ( CpDry * TempUppLevBefAdj  + LatentHeat * QVapUppLevBefAdj )  &
!!$                    &     * xyz_DelPress(i,j,k+1) / Grav                               &
!!$                    & - ( CpDry * xyz_Temp(i,j,k  ) + LatentHeat * xyz_QVap(i,j,k  ) ) &
!!$                    &     * xyz_DelPress(i,j,k  ) / Grav                               &
!!$                    & - ( CpDry * xyz_Temp(i,j,k+1) + LatentHeat * xyz_QVap(i,j,k+1) ) &
!!$                    &     * xyz_DelPress(i,j,k+1) / Grav,                              &
!!$                    &   ( CpDry * TempLowLevBefAdj  + LatentHeat * QVapLowLevBefAdj )  &
!!$                    &     * xyz_DelPress(i,j,k  ) / Grav,                              &
!!$                    &   ( CpDry * TempUppLevBefAdj  + LatentHeat * QVapUppLevBefAdj )  &
!!$                    &     * xyz_DelPress(i,j,k+1) / Grav,                              &
!!$                    &   ( CpDry * xyz_Temp(i,j,k  ) + LatentHeat * xyz_QVap(i,j,k  ) ) &
!!$                    &     * xyz_DelPress(i,j,k  ) / Grav,                              &
!!$                    &   ( CpDry * xyz_Temp(i,j,k+1) + LatentHeat * xyz_QVap(i,j,k+1) ) &
!!$                    &     * xyz_DelPress(i,j,k+1) / Grav
!!$                  write( 6, * ) 'Difference of moist static energy after adjustment'
!!$                  write( 6, * ) &
!!$                    &   ( CpDry * xyz_Temp(i,j,k  ) + LatentHeat * xyz_QVap(i,j,k  ) )  &
!!$                    & - ( CpDry * xyz_Temp(i,j,k+1) + LatentHeat * xyz_QVap(i,j,k+1) )  &
!!$                    & - CpDry * xyr_ConvAdjustFactor(i,j,k)                             &
!!$                    &   * ( xyz_Temp(i,j,k) + xyz_Temp(i,j,k+1) ),                      &
!!$                    &   ( CpDry * xyz_Temp(i,j,k  ) + LatentHeat * xyz_QVap(i,j,k  ) ), &
!!$                    &   ( CpDry * xyz_Temp(i,j,k+1) + LatentHeat * xyz_QVap(i,j,k+1) ), &
!!$                    & - CpDry * xyr_ConvAdjustFactor(i,j,k)                             &
!!$                    &   * ( xyz_Temp(i,j,k) + xyz_Temp(i,j,k+1) )
!!$                  write( 6, * ) 'Difference of water vapor amount before and after adjustment'
!!$                  write( 6, * ) &
!!$                    & - LatentHeat * ( xyz_QVap(i,j,k  ) - QVapLowLevBefAdj ) &
!!$                    & * xyz_DelPress(i,j,k  ) / Grav,                       &
!!$                    & - LatentHeat * ( xyz_QVap(i,j,k+1) - QVapUppLevBefAdj ) &
!!$                    & * xyz_DelPress(i,j,k+1) / Grav
                  !=========


                  xyz_QVapSat(i,j,k  ) = xyz_QVap(i,j,k  )
                  xyz_QVapSat(i,j,k+1) = xyz_QVap(i,j,k+1)

                  ! 調節したか否か?
                  ! Whether it was adjusted or not?
                  !
                  xy_Adjust(i,j) = .true.
                end if

              end if

            end if
          end do
        end do
      end do

      Adjust = .false.
      do i = 0, imax-1
        do j = 1, jmax
          xy_AdjustB(i,j) = xy_Adjust(i,j)
          Adjust          = Adjust .or. xy_Adjust(i,j)
        end do
      end do

      if ( .not. Adjust ) exit

    end do

    ! 比湿変化率, 温度変化率, 降水量の算出
    ! Calculate specific humidity tendency, temperature tendency, precipitation
    !
    xyz_DQVapDtCumulus = ( xyz_QVap - xyz_QVapB ) / ( 2.0_DP * DelTime )

    xyz_DTempDtCumulus = ( xyz_Temp - xyz_TempB ) / ( 2.0_DP * DelTime )



    xyz_DTempDt = xyz_DTempDt + xyz_DTempDtCumulus
    xyz_DQVapDt = xyz_DQVapDt + xyz_DQVapDtCumulus



    ! OLD CODE TO BE DELETED
    !
!!$    xy_RainCumulus = 0.
!!$    do k = 1, kmax
!!$      xy_RainCumulus = xy_RainCumulus &
!!$        & + ( xyz_QVapB(:,:,k) - xyz_QVap(:,:,k) ) &
!!$        &       * xyz_DelPress(:,:,k) / Grav / ( 2.0_DP * DelTime )
!!$    end do
!!$
!!$    xy_RainCumulus = 0.0d0
!!$    do k = kmax, 1, -1
!!$      xy_RainCumulus = xy_RainCumulus                &
!!$        & + xyz_DelPress(:,:,k) / Grav               &
!!$        &   * ( xyz_QVapB(:,:,k) - xyz_QVap(:,:,k) )
!!$    end do
!!$    xy_RainCumulus = xy_RainCumulus / ( 2.0_DP * DelTime )
!!$
!!$    xyz_RainCumulus = ( xyz_QVapB - xyz_QVap ) &
!!$        & * xyz_DelPress / Grav / ( 2.0_DP * DelTime )
!!$    xy_RainCumulus = 0.0d0
!!$    do k = kmax, 1, -1
!!$      xy_RainCumulus = xy_RainCumulus + xyz_RainCumulus(:,:,k)
!!$    end do

!!$    j = jmax/2+1
!!$    do i = 0, imax-1
!!$      if ( xy_RainCumulus(i,j) /= 0.0d0 ) then
!!$        write( 6, * ) xy_RainCumulus(i,j)
!!$      end if
!!$    end do


    xyz_DQH2OLiqDt = ( xyz_QVapB - xyz_QVap ) / ( 2.0_DP * DelTime )

!!$    xyz_RainCumulus = xyz_DQH2OLiqDt * xyz_DelPress / Grav
!!$    xy_RainCumulus = 0.0d0
!!$    do k = kmax, 1, -1
!!$      xy_RainCumulus = xy_RainCumulus + xyz_RainCumulus(:,:,k)
!!$    end do


!!$    j = jmax/2+1
!!$    do i = 0, imax-1
!!$      if ( xy_RainCumulus(i,j) /= 0.0d0 ) then
!!$        write( 6, * ) xy_RainCumulus(i,j)
!!$      end if
!!$    end do
!!$    write( 6, * ) '---'


!!$    xy_Rain     = xy_Rain     + xy_RainCumulus


    ! calculation for output (tentative)
    xyz_RainCumulus = xyz_DQH2OLiqDt * xyz_DelPress / Grav
    xy_RainCumulus = 0.0d0
    do k = kmax, 1, -1
      xy_RainCumulus = xy_RainCumulus + xyz_RainCumulus(:,:,k)
    end do


    ! ヒストリデータ出力
    ! History data output
    !
    call HistoryAutoPut( TimeN, 'RainCumulus',    xy_RainCumulus * LatentHeat )
    call HistoryAutoPut( TimeN, 'DTempDtCumulus', xyz_DTempDtCumulus )
    call HistoryAutoPut( TimeN, 'DQVapDtCumulus', xyz_DQVapDtCumulus )



!!$    if ( present( xyz_DQH2OLiqDt ) ) then
!!$
!!$      xyz_DDelLWDtCCPLV = &
!!$        & + ( xyz_QVapB - xyz_QVap ) &
!!$        &       * xyz_DelPress / Grav / ( 2.0d0 * DelTime )
!!$
!!$      ! Negative cloud production rate is filled with values in lower layers.
!!$      !
!!$      xy_NegDDelLWDt = 0.0d0
!!$      do k = kmax, 1, -1
!!$        do j = 1, jmax
!!$          do i = 0, imax-1
!!$            xyz_DDelLWDtCCPLV(i,j,k) = xyz_DDelLWDtCCPLV(i,j,k) + xy_NegDDelLWDt(i,j)
!!$            if ( xyz_DDelLWDtCCPLV(i,j,k) < 0.0d0 ) then
!!$              xy_NegDDelLWDt(i,j) = xyz_DDelLWDtCCPLV(i,j,k) 
!!$              xyz_DDelLWDtCCPLV(i,j,k) = 0.0d0
!!$            end if
!!$          end do
!!$        end do
!!$      end do
!!$
!!$      xyz_DQH2OLiqDt = xyz_DDelLWDtCCPLV / ( xyz_DelPress / Grav )
!!$
!!$    end if


    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )

  end subroutine MoistConvAdjust
Subroutine :

moist_conv_adjust モジュールの初期化を行います. NAMELIST#moist_conv_adjust_nml の読み込みはこの手続きで行われます.

"moist_conv_adjust" module is initialized. "NAMELIST#moist_conv_adjust_nml" is loaded in this procedure.

This procedure input/output NAMELIST#moist_conv_adjust_nml .

[Source]

  subroutine MoistConvAdjustInit
    !
    ! moist_conv_adjust モジュールの初期化を行います. 
    ! NAMELIST#moist_conv_adjust_nml の読み込みはこの手続きで行われます. 
    !
    ! "moist_conv_adjust" module is initialized. 
    ! "NAMELIST#moist_conv_adjust_nml" is loaded in this procedure. 
    !

    ! モジュール引用 ; USE statements
    !

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

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

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output

    ! 文字列操作
    ! Character handling
    !
    use dc_string, only: StoA

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

    ! 飽和比湿の算出
    ! Evaluate saturation specific humidity
    !
    use saturate, only: SaturateInit


    ! 宣言文 ; Declaration statements
    !
    implicit none

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

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /moist_conv_adjust_nml/ CrtlRH, ItrtMax, AdjustCriterion !, FlagUse

          ! デフォルト値については初期化手続 "moist_conv_adjust#CumAdjInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "moist_conv_adjust#MoistConvAdjustInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( moist_conv_adjust_inited ) return


    ! デフォルト値の設定
    ! Default values settings
    !
    ! default values used in AGCM5
!!$    CrtlRH  = 0.990_DP
!!$    ItrtMax = 10
!!$    AdjustCriterion(1:ItrtMax) = &
!!$      & (/ 0.01_DP, 0.02_DP, 0.02_DP, 0.05_DP, 0.05_DP, &
!!$      &    0.10_DP, 0.10_DP, 0.20_DP, 0.20_DP, 0.40_DP  /)
    !
    CrtlRH  = 1.0d0
    ItrtMax = 10
    AdjustCriterion(1:ItrtMax) = 0.0d0
    !
!!$    FlagUse = .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 = moist_conv_adjust_nml, iostat = iostat_nml )           ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
!!$      if ( iostat_nml == 0 ) write( STDOUT, nml = cumulus_adjust_nml )
    end if

    ! イテレーション回数, 不安定の許容誤差のチェック
    ! Check number of iteration, admissible error of unstability
    !
    call NmlutilAryValid( module_name, AdjustCriterion, 'AdjustCriterion', ItrtMax,    'ItrtMax' )          ! (in)


    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
    call HistoryAutoAddVariable( 'RainCumulus', (/ 'lon ', 'lat ', 'time' /), 'precipitation by cumulus scheme', 'W m-2' )
    call HistoryAutoAddVariable( 'DTempDtCumulus', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'cumulus condensation heating', 'K s-1' )
    call HistoryAutoAddVariable( 'DQVapDtCumulus', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'cumulus condensation moistening', 'kg kg-1 s-1' )


    ! Initialization of modules used in this module
    !
    call SaturateInit


    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
!!$    call MessageNotify( 'M', module_name, '  FlagUse              = %b', l = (/ FlagUse /) )
    call MessageNotify( 'M', module_name, '  CrtlRH               = %f', d = (/ CrtlRH /) )
    call MessageNotify( 'M', module_name, '  ItrtMax              = %d', i = (/ ItrtMax /) )
    call MessageNotify( 'M', module_name, '  AdjustCriterion      = (/ %*r /)', r = real( AdjustCriterion(1:ItrtMax) ), n = (/ ItrtMax /) )
    call MessageNotify( 'M', module_name, '' )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    moist_conv_adjust_inited = .true.

  end subroutine MoistConvAdjustInit

Private Instance methods

AdjustCriterion
Variable :
AdjustCriterion(1:MaxNmlArySize) :real(DP), save
: 調節を行う基準 (湿潤静的エネルギーの差の温度換算値) Criterion of adjustment (temperature difference equivalent to difference of moist static energy)
CrtlRH
Variable :
CrtlRH :real(DP), save
: 臨界相対湿度. Critical relative humidity
ItrtMax
Variable :
ItrtMax :integer, save
: イテレーション回数. Number of iteration
module_name
Constant :
module_name = ‘moist_conv_adjust :character(*), parameter
: モジュールの名称. Module name
moist_conv_adjust_inited
Variable :
moist_conv_adjust_inited = .false. :logical, save
: 初期設定フラグ. Initialization flag
version
Constant :
version = ’$Name: dcpam5-20140314 $’ // ’$Id: moist_conv_adjust.f90,v 1.8 2014-02-14 08:12:22 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version