Class lscond
In: lscond/lscond.f90

大規模凝結 (非対流性凝結)

Large scale condensation (non-convective condensation)

Note that Japanese and English are described in parallel.

大規模凝結過程によって温度と比湿を調節します.

Adjust temperature and specific humidity by a large scale condensation process.

References

 Manabe, S., J. Smagorinsky, R. F. Strickler,
   Simulated climatology of a general circulation model with a hydrologic cycle,
   Mon. Wea. Rev., 93, 769-798, 1965.

Procedures List

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

NAMELIST

NAMELIST#lscond_nml

Methods

Included Modules

gridset dc_types dc_message constants timeset gtool_historyauto saturate namelist_util dc_iounit dc_string cloud_utils

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
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)
: Tendency of H2O liquid mixing ratio
FlagOutput :logical , intent(in ), optional
: Flag for output

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCond( xyz_Temp, xyz_QVap, xyz_Press, xyr_Press, xyz_DQH2OLiqDt, FlagOutput )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

    ! 時刻管理
    ! 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, xyz_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(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)
                              !
                              ! Tendency of H2O liquid mixing ratio
    logical , intent(in ), optional :: FlagOutput
                              !
                              ! Flag for output


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

    real(DP):: xyz_QVapB (0:imax-1, 1:jmax, 1:kmax)
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: xyz_TempB (0:imax-1, 1:jmax, 1:kmax)
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: xyz_QVapSat      (0:imax-1, 1:jmax, 1:kmax)
                              ! 飽和比湿. 
                              ! Saturation specific humidity. 
    real(DP):: xyz_DQVapSatDTemp(0:imax-1, 1:jmax, 1:kmax)
                              ! $ \DD{q_{\rm{sat}}}{T} $
    real(DP):: DelTemp
                              ! 調節による温度変化量. 
                              ! Temperature variation by adjustment

    real(DP):: LatentHeatLocal
                              ! 
                              ! Latent heat used in this routine

    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):: xyz_RainLSC(0:imax-1, 1:jmax, 1:kmax)

    real(DP):: TempBefAdj
    real(DP):: QVapBefAdj

    logical:: xyz_FlagSaturated(0:imax-1, 1:jmax, 1:kmax)


    ! 実行文 ; Executable statement
    !

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


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


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


    ! 調節
    ! Adjustment
    !

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

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

          if ( ( xyz_QVap(i,j,k) / xyz_QVapSat(i,j,k) ) >= CrtlRH ) then
            xyz_FlagSaturated(i,j,k) = .true.
          else
            xyz_FlagSaturated(i,j,k) = .false.
          end if

        end do
      end do
    end do

    ! Set a value for latent heat
    if ( FlagSublimation ) then
      LatentHeatLocal = LatentHeat + LatentHeatFusion
    else
      LatentHeatLocal = LatentHeat
    end if

    do itr = 1, ItrtMax

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


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

            ! 飽和していたら, 温度と比湿の変化を計算
            ! Calculate tendency of temperature and humidity 
            ! if moist is saturation. 
            !
            if ( xyz_FlagSaturated(i,j,k) ) then

              ! 温度の変化分をニュートン法で求める
              ! Calculate variation of temperature
              !
              DelTemp = LatentHeatLocal * ( xyz_QVap(i,j,k) - CrtlRH * xyz_QVapSat(i,j,k) ) / ( CpDry + LatentHeatLocal * CrtlRH * xyz_DQVapSatDTemp(i,j,k) )


              !=========
              ! check routine
              !---------
!!$              TempBefAdj = xyz_Temp(i,j,k)
!!$              QVapBefAdj = xyz_QVap(i,j,k)
              !=========

              ! 温度と比湿の調節
              ! Adjust temperature and specific humidity
              !
              xyz_Temp(i,j,k) = xyz_Temp(i,j,k) + DelTemp
!!$              xyz_QVap(i,j,k) = xyz_QVapSat(i,j,k) + xyz_DQVapSatDTemp(i,j,k) * DelTemp
              xyz_QVap(i,j,k) = CrtlRH * ( xyz_QVapSat(i,j,k) + xyz_DQVapSatDTemp(i,j,k) * DelTemp )

              !=========
              ! check routine
              !---------
!!$              write( 6, * ) '====='
!!$              write( 6, * ) 'Energy difference before and after adjustment and each energy'
!!$              write( 6, * ) &
!!$                &     CpDry * TempBefAdj      + LatentHeatLocal * QVapBefAdj            &
!!$                & - ( CpDry * xyz_Temp(i,j,k) + LatentHeatLocal * xyz_QVap(i,j,k) ),    &
!!$                &     CpDry * TempBefAdj      + LatentHeatLocal * QVapBefAdj,           &
!!$                &   ( CpDry * xyz_Temp(i,j,k) + LatentHeatLocal * xyz_QVap(i,j,k) )
              !=========


            end if

          end do
        end do
      end do

    end do

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

    xyz_DQH2OLiqDt = - xyz_DQVapDtLsc


    ! calculation for output
    xy_RainLsc     = 0.0_DP
    do k = kmax, 1, -1
      xy_RainLsc = xy_RainLsc + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
    end do


    ! ヒストリデータ出力
    ! History data output
    !
    if ( .not. present( FlagOutput ) ) then
      call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
      call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
    else
      if ( FlagOutput ) then
        call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
        call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
      end if
    end if


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

  end subroutine LScaleCond
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
xyz_QH2OLiq(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: Specific liquid water content
xyz_QH2OSol(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: Specific ice content
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)
: Tendency of H2O liquid mixing ratio
xyz_DQH2OSolDt(0:imax-1,1:jmax,1:kmax) :real(DP), intent(out)
: Tendency of H2O ice mixing ratio
FlagOutput :logical , intent(in ), optional
: Flag for output

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCond1D3DWrapper( xyz_Temp, xyz_QVap, xyz_QH2OLiq, xyz_QH2OSol, xyz_Press, xyr_Press, xyz_DQH2OLiqDt, xyz_DQH2OSolDt, FlagOutput )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

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

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

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


    ! 宣言文 ; 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):: xyz_QH2OLiq(0:imax-1, 1:jmax, 1:kmax)
                              ! Specific liquid water content
    real(DP), intent(inout):: xyz_QH2OSol(0:imax-1, 1:jmax, 1:kmax)
                              ! Specific ice content
    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)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP), intent(out) :: xyz_DQH2OSolDt(0:imax-1,1:jmax,1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio
    logical , intent(in ), optional :: FlagOutput
                              !
                              ! Flag for output


    ! 作業変数
    ! Work variables
    !
    real(DP):: xy_RainLsc (0:imax-1, 1:jmax)
                              ! 降水量. 
                              ! Precipitation

    real(DP):: xyz_DTempDtLsc (0:imax-1, 1:jmax, 1:kmax)
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP):: xyz_DQVapDtLsc (0:imax-1, 1:jmax, 1:kmax)
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP) :: z_Temp (1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP) :: z_QVap (1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP) :: z_QH2OLiq(1:kmax)
                              ! Specic liquid water content
    real(DP) :: z_QH2OSol(1:kmax)
                              ! Specic ice content
    real(DP) :: z_Press (1:kmax)
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP) :: r_Press (0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP) :: z_DQH2OLiqDt(1:kmax)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP) :: z_DQH2OSolDt(1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio


    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

    ! 実行文 ; Executable statement
    !

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


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



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

        do k = 1, kmax
          z_Temp   (k) = xyz_Temp   (i,j,k)
          z_QVap   (k) = xyz_QVap   (i,j,k)
          z_QH2OLiq(k) = xyz_QH2OLiq(i,j,k)
          z_QH2OSol(k) = xyz_QH2OSol(i,j,k)
          z_Press  (k) = xyz_Press  (i,j,k)
        end do
        do k = 0, kmax
          r_Press  (k) = xyr_Press  (i,j,k)
        end do

        call LScaleCond1D( z_Temp, z_QVap, z_QH2OLiq, z_QH2OSol, z_Press, r_Press, z_DQH2OLiqDt, z_DQH2OSolDt )

        do k = 1, kmax
          xyz_Temp      (i,j,k) = z_Temp      (k)
          xyz_QVap      (i,j,k) = z_QVap      (k)
          xyz_QH2OLiq   (i,j,k) = z_QH2OLiq   (k)
          xyz_QH2OSol   (i,j,k) = z_QH2OSol   (k)
          xyz_DQH2OLiqDt(i,j,k) = z_DQH2OLiqDt(k)
          xyz_DQH2OSolDt(i,j,k) = z_DQH2OSolDt(k)
        end do

      end do
    end do


    ! calculation for output
    xy_RainLsc     = 0.0_DP
    do k = kmax, 1, -1
      xy_RainLsc = xy_RainLsc + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
    end do


    ! ヒストリデータ出力
    ! History data output
    !
    if ( .not. present( FlagOutput ) ) then
      call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
      call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
    else
      if ( FlagOutput ) then
        call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
        call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
      end if
    end if


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

  end subroutine LScaleCond1D3DWrapper
Subroutine :
Temp :real(DP), intent(inout)
: $ T $ . 温度. Temperature
QVap :real(DP), intent(inout)
: $ q $ . 比湿. Specific humidity
QH2OLiq :real(DP), intent(inout)
: Specific liquid water content
QH2OSol :real(DP), intent(inout)
: Specific ice content
Press :real(DP), intent(in )
: $ p $ . 気圧 (整数レベル). Air pressure (full level)
PressLI :real(DP), intent(in )
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
PressUI :real(DP), intent(in )
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
DQH2OLiqDt :real(DP), intent(out )
: Tendency of H2O liquid mixing ratio
DQH2OSolDt :real(DP), intent(out )
: Tendency of H2O ice mixing ratio

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCond1Grid( Temp, QVap, QH2OLiq, QH2OSol, Press, PressLI, PressUI, DQH2OLiqDt, DQH2OSolDt )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

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

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

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

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


    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(inout):: Temp
                              ! $ T $ .     温度. Temperature
    real(DP), intent(inout):: QVap
                              ! $ q $ .     比湿. Specific humidity
    real(DP), intent(inout):: QH2OLiq
                              ! Specific liquid water content
    real(DP), intent(inout):: QH2OSol
                              ! Specific ice content
    real(DP), intent(in   ):: Press
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP), intent(in   ):: PressLI
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(in   ):: PressUI
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(out  ):: DQH2OLiqDt
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP), intent(out  ):: DQH2OSolDt
                              !
                              ! Tendency of H2O ice mixing ratio


    ! 作業変数
    ! Work variables
    !
    real(DP):: DTempDtLsc
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP):: DQVapDtLsc
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP):: QH2OVap0
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: QH2OLiq0
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: QH2OSol0
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: Temp0
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: QH2OVapB
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: QH2OLiqB
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: QH2OSolB
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: TempB
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: QVapSat
                              ! 飽和比湿. 
                              ! Saturation specific humidity. 
    real(DP):: DQVapSatDTemp
                              ! $ \DD{q_{\rm{sat}}}{T} $
    real(DP):: DelTemp
                              ! 調節による温度変化量. 
                              ! Temperature variation by adjustment

    real(DP) :: WatFrac
    real(DP) :: IceFrac


    integer:: itr             ! イテレーション方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in iteration direction
    logical:: FlagProcess


    ! 実行文 ; Executable statement
    !

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


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


    ! 調節前 "QVap", "Temp" の保存
    ! Store "QVap", "Temp" before adjustment
    !
    QH2OVap0 = QVap
    QH2OLiq0 = QH2OLiq
    QH2OSol0 = QH2OSol
    Temp0    = Temp


    ! 調節
    ! Adjustment
    !

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

    if ( FlagSatAdj ) then
      Temp = Temp + ( LatentHeat * QH2OLiq + ( LatentHeat + LatentHeatFusion ) * QH2OSol ) / CpDry
      QVap    = QVap + QH2OLiq + QH2OSol
      QH2OLiq = 0.0_DP
      QH2OSol = 0.0_DP
    else
      ! Pre-existing cloud water and ice are removed, since those are not 
      ! responsible with condensation (assumption).
      QH2OLiq = 0.0_DP
      QH2OSol = 0.0_DP
    end if
    !
    QVapSat = CalcQVapSat( Temp, Press )
    !
    if ( ( QVap / QVapSat ) >= CrtlRH ) then
      FlagProcess = .true.
    else
      FlagProcess = .false.
    end if


    do itr = 1, ItrtMax

      QH2OVapB = QVap
      QH2OLiqB = QH2OLiq
      QH2OSolB = QH2OSol
      TempB    = Temp

      ! 飽和比湿計算
      ! Calculate saturation specific humidity
      !
      QVapSat       = CalcQVapSat      ( Temp, Press   )
      DQVapSatDTemp = CalcDQVapSatDTemp( Temp, QVapSat )


      ! 飽和していたら, 温度と比湿の変化を計算
      ! Calculate tendency of temperature and humidity 
      ! if moist is saturation. 
      !
      if ( FlagProcess ) then

        ! Liquid water and ice fractions
        call SaturateWatFraction( Temp, WatFrac )
        IceFrac = 1.0_DP - WatFrac

        ! 温度の変化分をニュートン法で求める
        ! Calculate variation of temperature
        !
!!$          DelTemp =                                                 &
!!$            & LatentHeatLocal                                       &
!!$            &   * ( z_QVap(k) - CrtlRH * z_QVapSat(k) )             &
!!$            &   / ( CpDry + LatentHeatLocal * CrtlRH * z_DQVapSatDTemp(k) )
        DelTemp = (   LatentHeat * QVap - LatentHeatFusion * QH2OSol + LatentHeatFusion * IceFrac * ( QVap + QH2OLiq + QH2OSol ) - ( LatentHeat + LatentHeatFusion * IceFrac ) * CrtlRH * QVapSat ) / (   CpDry + ( LatentHeat + LatentHeatFusion * IceFrac ) * CrtlRH * DQVapSatDTemp )


        ! 温度と比湿の調節
        ! Adjust temperature and specific humidity
        !
        Temp = Temp + DelTemp
        QVap = CrtlRH * ( QVapSat + DQVapSatDTemp * DelTemp )
        QH2OLiq = WatFrac * ( QH2OVapB + QH2OLiqB + QH2OSolB - QVap )
        QH2OSol = IceFrac * ( QH2OVapB + QH2OLiqB + QH2OSolB - QVap )

      end if

    end do


    if ( FlagSatAdj ) then
    else
      ! Pre-existing cloud water and ice amount are restored. 
      QH2OLiq = QH2OLiq + QH2OLiq0
      QH2OSol = QH2OSol + QH2OSol0
    end if


    ! 比湿変化率, 温度変化率, 降水量の算出
    ! Calculate specific humidity tendency, temperature tendency, 
    ! precipitation
    !
    DTempDtLsc = ( Temp    - Temp0    ) / ( 2.0_DP * DelTime )
    DQVapDtLsc = ( QVap    - QH2OVap0 ) / ( 2.0_DP * DelTime )
    DQH2OLiqDt = ( QH2OLiq - QH2OLiq0 ) / ( 2.0_DP * DelTime )
    DQH2OSolDt = ( QH2OSol - QH2OSol0 ) / ( 2.0_DP * DelTime )


    call LScaleCond1GridConsChk( Temp0, QH2OVap0, QH2OLiq0, QH2OSol0, Temp , QVap    , QH2OLiq , QH2OSol )


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

  end subroutine LScaleCond1Grid
Subroutine :
FlagSnow :logical, intent(in)

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

"lscond" module is initialized. "NAMELIST#lscond_nml" is loaded in this procedure.

This procedure input/output NAMELIST#lscond_nml .

[Source]

  subroutine LScaleCondInit( FlagSnow )
    !
    ! lscond モジュールの初期化を行います. 
    ! NAMELIST#lscond_nml の読み込みはこの手続きで行われます. 
    !
    ! "lscond" module is initialized. 
    ! "NAMELIST#lscond_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

    ! 雲関系ルーチン
    ! Cloud-related routines
    !
    use cloud_utils, only : CloudUtilsInit

    ! 宣言文 ; Declaration statements
    !
    implicit none

    logical, intent(in) :: FlagSnow


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

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /lscond_nml/ QH2OFlucRatio, CrtlRH, ItrtMax, FlagSublimation, FlagSatAdj
          !
          ! デフォルト値については初期化手続 "lscond#LSCondInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "lscond#LSCondInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( lscond_inited ) return


    ! デフォルト値の設定
    ! Default values settings
    !
    QH2OFlucRatio   = 0.2_DP
    CrtlRH          = 1.0_DP
    ItrtMax         = 3
    FlagSublimation = .false.
    FlagSatAdj      = .false.

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

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

    ! Initialization of modules used in this module
    !

    ! Initialization of modules used in this module
    !
    call SaturateInit

    ! 雲関系ルーチン
    ! Cloud-related routines
    !
    call CloudUtilsInit( FlagSnow )


    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
    call HistoryAutoAddVariable( 'DTempDtLsc', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'large-scale condensation heating', 'K s-1' )
    call HistoryAutoAddVariable( 'DQVapDtLsc', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'large-scale condensation moistening', 'kg kg-1 s-1' )

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, '  QH2OFlucRatio   = %f', d = (/ QH2OFlucRatio /) )
    call MessageNotify( 'M', module_name, '  CrtlRH          = %f', d = (/ CrtlRH /) )
    call MessageNotify( 'M', module_name, '  ItrtMax         = %d', i = (/ ItrtMax /) )
    call MessageNotify( 'M', module_name, '  FlagSublimation = %b', l = (/ FlagSublimation /) )
    call MessageNotify( 'M', module_name, '  FlagSatAdj      = %b', l = (/ FlagSatAdj /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    lscond_inited = .true.

  end subroutine LScaleCondInit
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
xyz_QH2OLiq(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: Specific liquid water content
xyz_QH2OSol(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)
: Specific ice content
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)
: Tendency of H2O liquid mixing ratio
xyz_DQH2OSolDt(0:imax-1,1:jmax,1:kmax) :real(DP), intent(out)
: Tendency of H2O ice mixing ratio
xyz_CloudCover(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out), optional
: Cloud cover
FlagOutput :logical , intent(in ), optional
: Flag for output

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCondLL911D3DWrapper( xyz_Temp, xyz_QVap, xyz_QH2OLiq, xyz_QH2OSol, xyz_Press, xyr_Press, xyz_DQH2OLiqDt, xyz_DQH2OSolDt, xyz_CloudCover, FlagOutput )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

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

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

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


    ! 宣言文 ; 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):: xyz_QH2OLiq(0:imax-1, 1:jmax, 1:kmax)
                              ! Specific liquid water content
    real(DP), intent(inout):: xyz_QH2OSol(0:imax-1, 1:jmax, 1:kmax)
                              ! Specific ice content
    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)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP), intent(out) :: xyz_DQH2OSolDt(0:imax-1,1:jmax,1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio
    real(DP), intent(out), optional :: xyz_CloudCover(0:imax-1, 1:jmax, 1:kmax)
                              !
                              ! Cloud cover
    logical , intent(in ), optional :: FlagOutput
                              !
                              ! Flag for output


    ! 作業変数
    ! Work variables
    !
    real(DP):: xy_RainLsc (0:imax-1, 1:jmax)
                              ! 降水量. 
                              ! Precipitation

    real(DP):: xyz_DTempDtLsc (0:imax-1, 1:jmax, 1:kmax)
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP):: xyz_DQVapDtLsc (0:imax-1, 1:jmax, 1:kmax)
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP) :: z_Temp (1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP) :: z_QVap (1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP) :: z_QH2OLiq(1:kmax)
                              ! Specic liquid water content
    real(DP) :: z_QH2OSol(1:kmax)
                              ! Specic ice content
    real(DP) :: z_Press (1:kmax)
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP) :: r_Press (0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP) :: z_DQH2OLiqDt(1:kmax)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP) :: z_DQH2OSolDt(1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio

    real(DP) :: z_CloudCover(1:kmax)


    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

    ! 実行文 ; Executable statement
    !

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


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



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

        do k = 1, kmax
          z_Temp   (k) = xyz_Temp   (i,j,k)
          z_QVap   (k) = xyz_QVap   (i,j,k)
          z_QH2OLiq(k) = xyz_QH2OLiq(i,j,k)
          z_QH2OSol(k) = xyz_QH2OSol(i,j,k)
          z_Press  (k) = xyz_Press  (i,j,k)
        end do
        do k = 0, kmax
          r_Press  (k) = xyr_Press  (i,j,k)
        end do

        call LScaleCondLL911D( z_Temp, z_QVap, z_QH2OLiq, z_QH2OSol, z_Press, r_Press, z_DQH2OLiqDt, z_DQH2OSolDt, z_CloudCover )

        do k = 1, kmax
          xyz_Temp      (i,j,k) = z_Temp      (k)
          xyz_QVap      (i,j,k) = z_QVap      (k)
          xyz_QH2OLiq   (i,j,k) = z_QH2OLiq   (k)
          xyz_QH2OSol   (i,j,k) = z_QH2OSol   (k)
          xyz_DQH2OLiqDt(i,j,k) = z_DQH2OLiqDt(k)
          xyz_DQH2OSolDt(i,j,k) = z_DQH2OSolDt(k)
        end do
        if ( present( xyz_CloudCover ) ) then
          do k = 1, kmax
            xyz_CloudCover(i,j,k) = z_CloudCover(k)
          end do
        end if

      end do
    end do


    ! calculation for output
    xy_RainLsc     = 0.0_DP
    do k = kmax, 1, -1
      xy_RainLsc = xy_RainLsc + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
    end do


    ! ヒストリデータ出力
    ! History data output
    !
    if ( .not. present( FlagOutput ) ) then
      call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
      call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
    else
      if ( FlagOutput ) then
        call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
        call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
      end if
    end if


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

  end subroutine LScaleCondLL911D3DWrapper

Private Instance methods

CrtlRH
Variable :
CrtlRH :real(DP), save
: 臨界相対湿度. Critical relative humidity
FlagSatAdj
Variable :
FlagSatAdj :logical, save
: flag for saturation adjustment
FlagSublimation
Variable :
FlagSublimation :logical, save
: flag for treating sublimation
ItrtMax
Variable :
ItrtMax :integer, save
: イテレーション回数. Number of iteration
Subroutine :
z_Temp(1:kmax) :real(DP), intent(inout)
: $ T $ . 温度. Temperature
z_QVap(1:kmax) :real(DP), intent(inout)
: $ q $ . 比湿. Specific humidity
z_QH2OLiq(1:kmax) :real(DP), intent(inout)
: Specific liquid water content
z_QH2OSol(1:kmax) :real(DP), intent(inout)
: Specific ice content
z_Press(1:kmax) :real(DP), intent(in )
: $ p $ . 気圧 (整数レベル). Air pressure (full level)
r_Press(0:kmax) :real(DP), intent(in )
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
z_DQH2OLiqDt(1:kmax) :real(DP), intent(out )
: Tendency of H2O liquid mixing ratio
z_DQH2OSolDt(1:kmax) :real(DP), intent(out )
: Tendency of H2O ice mixing ratio

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCond1D( z_Temp, z_QVap, z_QH2OLiq, z_QH2OSol, z_Press, r_Press, z_DQH2OLiqDt, z_DQH2OSolDt )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

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

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

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


    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(inout):: z_Temp (1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(inout):: z_QVap (1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP), intent(inout):: z_QH2OLiq(1:kmax)
                              ! Specific liquid water content
    real(DP), intent(inout):: z_QH2OSol(1:kmax)
                              ! Specific ice content
    real(DP), intent(in   ):: z_Press (1:kmax)
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP), intent(in   ):: r_Press (0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(out  ):: z_DQH2OLiqDt(1:kmax)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP), intent(out  ):: z_DQH2OSolDt(1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio


    ! 作業変数
    ! Work variables
    !
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction


    ! 実行文 ; Executable statement
    !

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


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


    do k = 1, kmax
      call LScaleCond1Grid( z_Temp(k), z_QVap(k), z_QH2OLiq(k), z_QH2OSol(k), z_Press(k), r_Press(k-1), r_Press(k), z_DQH2OLiqDt(k), z_DQH2OSolDt(k) )
    end do


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

  end subroutine LScaleCond1D
Subroutine :
z_TempB(1:kmax) :real(DP), intent(in)
z_QH2OVapB(1:kmax) :real(DP), intent(in)
z_QH2OLiqB(1:kmax) :real(DP), intent(in)
z_QH2OSolB(1:kmax) :real(DP), intent(in)
z_Temp(1:kmax) :real(DP), intent(in)
z_QH2OVap(1:kmax) :real(DP), intent(in)
z_QH2OLiq(1:kmax) :real(DP), intent(in)
z_QH2OSol(1:kmax) :real(DP), intent(in)

[Source]

  subroutine LScaleCond1DConsChk( z_TempB, z_QH2OVapB, z_QH2OLiqB, z_QH2OSolB, z_Temp , z_QH2OVap , z_QH2OLiq , z_QH2OSol )

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

    real(DP), intent(in) :: z_TempB   (1:kmax)
    real(DP), intent(in) :: z_QH2OVapB(1:kmax)
    real(DP), intent(in) :: z_QH2OLiqB(1:kmax)
    real(DP), intent(in) :: z_QH2OSolB(1:kmax)
    real(DP), intent(in) :: z_Temp    (1:kmax)
    real(DP), intent(in) :: z_QH2OVap (1:kmax)
    real(DP), intent(in) :: z_QH2OLiq (1:kmax)
    real(DP), intent(in) :: z_QH2OSol (1:kmax)

    ! Local variables
    !
    real(DP) :: ValB
    real(DP) :: Val
    real(DP) :: Ratio
    integer  :: k


    do k = kmax, 1, -1
      Val =   CpDry * z_TempB(k) + LatentHeat * z_QH2OVapB(k) - LatentHeatFusion * z_QH2OSolB(k)
      ValB = Val

      Val =   CpDry * z_Temp(k) + LatentHeat * z_QH2OVap(k) - LatentHeatFusion * z_QH2OSol(k)

      Ratio = ( Val - ValB ) / ( Val + 1.0d-100 )
      if ( abs( Ratio ) > 1.0d-10 ) then
        call MessageNotify( 'M', module_name, 'Simplified condensate static energy is not conserved, %f.', d = (/ Ratio /) )
      end if
    end do



    do k = kmax, 1, -1
      Val = z_QH2OVapB(k) + z_QH2OLiqB(k) + z_QH2OSolB(k)
      ValB = Val

      Val = z_QH2OVap (k) + z_QH2OLiq (k) + z_QH2OSol (k)

      Ratio = ( Val - ValB ) / ( Val + 1.0d-100 )
      if ( abs( Ratio ) > 1.0d-10 ) then
        call MessageNotify( 'M', module_name, 'H2O mass is not conserved, %f.', d = (/ Ratio /) )
      end if
    end do


  end subroutine LScaleCond1DConsChk
Subroutine :
z_Temp(1:kmax) :real(DP), intent(inout)
: $ T $ . 温度. Temperature
z_QVap(1:kmax) :real(DP), intent(inout)
: $ q $ . 比湿. Specific humidity
z_QH2OLiq(1:kmax) :real(DP), intent(inout)
: Specific liquid water content
z_QH2OSol(1:kmax) :real(DP), intent(inout)
: Specific ice content
z_Press(1:kmax) :real(DP), intent(in )
: $ p $ . 気圧 (整数レベル). Air pressure (full level)
r_Press(0:kmax) :real(DP), intent(in )
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
z_DQH2OLiqDt(1:kmax) :real(DP), intent(out )
: Tendency of H2O liquid mixing ratio
z_DQH2OSolDt(1:kmax) :real(DP), intent(out )
: Tendency of H2O ice mixing ratio

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCond1D_BK( z_Temp, z_QVap, z_QH2OLiq, z_QH2OSol, z_Press, r_Press, z_DQH2OLiqDt, z_DQH2OSolDt )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

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

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

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

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


    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(inout):: z_Temp (1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(inout):: z_QVap (1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP), intent(inout):: z_QH2OLiq(1:kmax)
                              ! Specific liquid water content
    real(DP), intent(inout):: z_QH2OSol(1:kmax)
                              ! Specific ice content
    real(DP), intent(in   ):: z_Press (1:kmax)
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP), intent(in   ):: r_Press (0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(out  ):: z_DQH2OLiqDt(1:kmax)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP), intent(out  ):: z_DQH2OSolDt(1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio


    ! 作業変数
    ! Work variables
    !
    real(DP):: z_DTempDtLsc (1:kmax)
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP):: z_DQVapDtLsc (1:kmax)
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP):: z_QH2OVap0 (1:kmax)
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: z_QH2OLiq0 (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_QH2OSol0 (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_Temp0    (1:kmax)
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: z_QH2OVapB (1:kmax)
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: z_QH2OLiqB (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_QH2OSolB (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_TempB    (1:kmax)
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: z_QVapSat      (1:kmax)
                              ! 飽和比湿. 
                              ! Saturation specific humidity. 
    real(DP):: z_DQVapSatDTemp(1:kmax)
                              ! $ \DD{q_{\rm{sat}}}{T} $
    real(DP):: DelTemp
                              ! 調節による温度変化量. 
                              ! Temperature variation by adjustment

    real(DP) :: WatFrac
    real(DP) :: IceFrac


    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    integer:: itr             ! イテレーション方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in iteration direction
    logical:: z_FlagProcess(1:kmax)


    ! 実行文 ; Executable statement
    !

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


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


    ! 調節前 "QVap", "Temp" の保存
    ! Store "QVap", "Temp" before adjustment
    !
    z_QH2OVap0 = z_QVap
    z_QH2OLiq0 = z_QH2OLiq
    z_QH2OSol0 = z_QH2OSol
    z_Temp0    = z_Temp


    ! 調節
    ! Adjustment
    !

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

    if ( FlagSatAdj ) then
      z_Temp = z_Temp + ( LatentHeat * z_QH2OLiq + ( LatentHeat + LatentHeatFusion ) * z_QH2OSol ) / CpDry
      z_QVap = z_QVap + z_QH2OLiq + z_QH2OSol
      z_QH2OLiq = 0.0_DP
      z_QH2OSol = 0.0_DP
    else
      ! Pre-existing cloud water and ice are removed, since those are not 
      ! responsible with condensation (assumption).
      z_QH2OLiq = 0.0_DP
      z_QH2OSol = 0.0_DP
    end if
    !
    z_QVapSat = a_CalcQVapSat( z_Temp, z_Press )
    !
    do k = 1, kmax
      if ( ( z_QVap(k) / z_QVapSat(k) ) >= CrtlRH ) then
        z_FlagProcess(k) = .true.
      else
        z_FlagProcess(k) = .false.
      end if
    end do


    do itr = 1, ItrtMax

      z_QH2OVapB = z_QVap
      z_QH2OLiqB = z_QH2OLiq
      z_QH2OSolB = z_QH2OSol
      z_TempB    = z_Temp

      ! 飽和比湿計算
      ! Calculate saturation specific humidity
      !
      z_QVapSat       = a_CalcQVapSat      ( z_Temp, z_Press   )
      z_DQVapSatDTemp = a_CalcDQVapSatDTemp( z_Temp, z_QVapSat )


      do k = 1, kmax

        ! 飽和していたら, 温度と比湿の変化を計算
        ! Calculate tendency of temperature and humidity 
        ! if moist is saturation. 
        !
        if ( z_FlagProcess(k) ) then

          ! Liquid water and ice fractions
          call SaturateWatFraction( z_Temp(k), WatFrac )
          IceFrac = 1.0_DP - WatFrac

          ! 温度の変化分をニュートン法で求める
          ! Calculate variation of temperature
          !
          DelTemp = (   LatentHeat * z_QVap(k) - LatentHeatFusion * z_QH2OSol(k) + LatentHeatFusion * IceFrac * ( z_QVap(k) + z_QH2OLiq(k) + z_QH2OSol(k) ) - ( LatentHeat + LatentHeatFusion * IceFrac ) * CrtlRH * z_QVapSat(k) ) / (   CpDry + ( LatentHeat + LatentHeatFusion * IceFrac ) * CrtlRH * z_DQVapSatDTemp(k) )


          ! 温度と比湿の調節
          ! Adjust temperature and specific humidity
          !
          z_Temp(k) = z_Temp(k) + DelTemp
          z_QVap(k) = CrtlRH * ( z_QVapSat(k) + z_DQVapSatDTemp(k) * DelTemp )
          z_QH2OLiq(k) = WatFrac * ( z_QH2OVapB(k) + z_QH2OLiqB(k) + z_QH2OSolB(k) - z_QVap(k) )
          z_QH2OSol(k) = IceFrac * ( z_QH2OVapB(k) + z_QH2OLiqB(k) + z_QH2OSolB(k) - z_QVap(k) )

        end if

      end do

    end do


    if ( FlagSatAdj ) then
    else
      ! Pre-existing cloud water and ice amount are restored. 
      z_QH2OLiq = z_QH2OLiq + z_QH2OLiq0
      z_QH2OSol = z_QH2OSol + z_QH2OSol0
    end if


    ! 比湿変化率, 温度変化率, 降水量の算出
    ! Calculate specific humidity tendency, temperature tendency, 
    ! precipitation
    !
    z_DTempDtLsc = ( z_Temp    - z_Temp0    ) / ( 2.0_DP * DelTime )
    z_DQVapDtLsc = ( z_QVap    - z_QH2OVap0 ) / ( 2.0_DP * DelTime )
    z_DQH2OLiqDt = ( z_QH2OLiq - z_QH2OLiq0 ) / ( 2.0_DP * DelTime )
    z_DQH2OSolDt = ( z_QH2OSol - z_QH2OSol0 ) / ( 2.0_DP * DelTime )


    call LScaleCond1DConsChk( z_Temp0, z_QH2OVap0, z_QH2OLiq0, z_QH2OSol0, z_Temp , z_QVap    , z_QH2OLiq , z_QH2OSol )


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

  end subroutine LScaleCond1D_BK
Subroutine :
TempB :real(DP), intent(in)
QH2OVapB :real(DP), intent(in)
QH2OLiqB :real(DP), intent(in)
QH2OSolB :real(DP), intent(in)
Temp :real(DP), intent(in)
QH2OVap :real(DP), intent(in)
QH2OLiq :real(DP), intent(in)
QH2OSol :real(DP), intent(in)

[Source]

  subroutine LScaleCond1GridConsChk( TempB, QH2OVapB, QH2OLiqB, QH2OSolB, Temp , QH2OVap , QH2OLiq , QH2OSol )

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

    real(DP), intent(in) :: TempB
    real(DP), intent(in) :: QH2OVapB
    real(DP), intent(in) :: QH2OLiqB
    real(DP), intent(in) :: QH2OSolB
    real(DP), intent(in) :: Temp
    real(DP), intent(in) :: QH2OVap
    real(DP), intent(in) :: QH2OLiq
    real(DP), intent(in) :: QH2OSol

    ! Local variables
    !
    real(DP) :: ValB
    real(DP) :: Val
    real(DP) :: Ratio
    integer  :: k


    Val =   CpDry * TempB + LatentHeat * QH2OVapB - LatentHeatFusion * QH2OSolB
    ValB = Val

    Val =   CpDry * Temp + LatentHeat * QH2OVap - LatentHeatFusion * QH2OSol

    Ratio = ( Val - ValB ) / ( Val + 1.0d-100 )
    if ( abs( Ratio ) > 1.0d-10 ) then
      call MessageNotify( 'M', module_name, 'Simplified condensate static energy is not conserved, %f.', d = (/ Ratio /) )
    end if


    Val = QH2OVapB + QH2OLiqB + QH2OSolB
    ValB = Val

    Val = QH2OVap  + QH2OLiq  + QH2OSol

    Ratio = ( Val - ValB ) / ( Val + 1.0d-100 )
    if ( abs( Ratio ) > 1.0d-10 ) then
      call MessageNotify( 'M', module_name, 'H2O mass is not conserved, %f.', d = (/ Ratio /) )
    end if


  end subroutine LScaleCond1GridConsChk
Subroutine :
z_Temp(1:kmax) :real(DP), intent(inout)
: $ T $ . 温度. Temperature
z_QH2OVap(1:kmax) :real(DP), intent(inout)
: $ q $ . 比湿. Specific humidity
z_QH2OLiq(1:kmax) :real(DP), intent(inout)
: Specific liquid water content
z_QH2OSol(1:kmax) :real(DP), intent(inout)
: Specific ice content
z_Press(1:kmax) :real(DP), intent(in )
: $ p $ . 気圧 (整数レベル). Air pressure (full level)
r_Press(0:kmax) :real(DP), intent(in )
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
z_DQH2OLiqDt(1:kmax) :real(DP), intent(out )
: Tendency of H2O liquid mixing ratio
z_DQH2OSolDt(1:kmax) :real(DP), intent(out )
: Tendency of H2O ice mixing ratio
z_CloudCover(1:kmax) :real(DP), intent(out )

大規模凝結スキームにより, 温度と比湿を調節します.

Adjust temperature and specific humidity by large scale condensation scheme.

[Source]

  subroutine LScaleCondLL911D( z_Temp, z_QH2OVap, z_QH2OLiq, z_QH2OSol, z_Press, r_Press, z_DQH2OLiqDt, z_DQH2OSolDt, z_CloudCover )
    !
    ! 大規模凝結スキームにより, 温度と比湿を調節します. 
    !
    ! Adjust temperature and specific humidity by 
    ! large scale condensation scheme. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
                              ! $ L $ [J kg-1] .
                              ! 融解の潜熱.
                              ! Latent heat of fusion

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

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

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

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


    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(inout):: z_Temp (1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(inout):: z_QH2OVap (1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP), intent(inout):: z_QH2OLiq(1:kmax)
                              ! Specific liquid water content
    real(DP), intent(inout):: z_QH2OSol(1:kmax)
                              ! Specific ice content
    real(DP), intent(in   ):: z_Press (1:kmax)
                              ! $ p $ . 気圧 (整数レベル). 
                              ! Air pressure (full level)
    real(DP), intent(in   ):: r_Press (0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(out  ):: z_DQH2OLiqDt(1:kmax)
                              !
                              ! Tendency of H2O liquid mixing ratio
    real(DP), intent(out  ):: z_DQH2OSolDt(1:kmax)
                              !
                              ! Tendency of H2O ice mixing ratio
    real(DP), intent(out  ):: z_CloudCover(1:kmax)


    ! 作業変数
    ! Work variables
    !
    real(DP):: z_DTempDtLsc (1:kmax)
                              ! 温度変化率. 
                              ! Temperature tendency
    real(DP):: z_DQH2OVapDtLsc (1:kmax)
                              ! 比湿変化率. 
                              ! Specific humidity tendency

    real(DP):: z_QH2OTot  (1:kmax)
                              ! 調節前の水量. 
                              ! Specific water content before adjust. 

    real(DP):: z_QH2OVap0 (1:kmax)
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: z_QH2OLiq0 (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_QH2OSol0 (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_Temp0    (1:kmax)
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: z_QH2OVapB (1:kmax)
                              ! 調節前の比湿. 
                              ! Specific humidity before adjust. 
    real(DP):: z_QH2OLiqB (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_QH2OSolB (1:kmax)
                              ! 
                              ! Specific liquid water content before adjust. 
    real(DP):: z_TempB    (1:kmax)
                              ! 調節前の温度. 
                              ! Temperature before adjust. 
                              !
    real(DP):: z_QH2OVapSat      (1:kmax)
                              ! 飽和比湿. 
                              ! Saturation specific humidity. 
    real(DP):: z_DQH2OVapSatDTemp(1:kmax)
                              ! $ \DD{q_{\rm{sat}}}{T} $
    real(DP):: DelTemp
                              ! 調節による温度変化量. 
                              ! Temperature variation by adjustment

    real(DP) :: CoefA
    real(DP) :: CoefB
    real(DP) :: QH2OCond

    real(DP):: z_QH2OTotFluc(1:kmax)

    real(DP) :: WatFrac
    real(DP) :: IceFrac

    real(DP) :: QH2OVapSat

    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    integer:: itr             ! イテレーション方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in iteration direction
    logical:: z_FlagProcess(1:kmax)


    ! 実行文 ; Executable statement
    !

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


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


    ! 調節前 "QH2OVap", "Temp" の保存
    ! Store "QH2OVap", "Temp" before adjustment
    !
    z_QH2OVap0 = z_QH2OVap
    z_QH2OLiq0 = z_QH2OLiq
    z_QH2OSol0 = z_QH2OSol
    z_Temp0    = z_Temp

    ! 調節
    ! Adjustment
    !

    z_QH2OTot     = z_QH2OVap + z_QH2OLiq + z_QH2OSol
    z_QH2OTotFluc = QH2OFlucRatio * z_QH2OTot

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

!!$    z_QH2OVapSat = a_CalcQVapSat( z_Temp, z_Press )
!!$    !
!!$    do k = 1, kmax
!!$      if ( ( z_QH2OTot(k) / z_QH2OVapSat(k) ) >= 1.0_DP ) then
!!$        z_FlagProcess(k) = .true.
!!$      else
!!$        z_FlagProcess(k) = .false.
!!$      end if
!!$    end do

    z_FlagProcess = .true.


    do itr = 1, ItrtMax

      z_QH2OVapB = z_QH2OVap
      z_QH2OLiqB = z_QH2OLiq
      z_QH2OSolB = z_QH2OSol
      z_TempB    = z_Temp

      ! 飽和比湿計算
      ! Calculate saturation specific humidity
      !
      z_QH2OVapSat       = a_CalcQVapSat      ( z_Temp, z_Press      )
      z_DQH2OVapSatDTemp = a_CalcDQVapSatDTemp( z_Temp, z_QH2OVapSat )


      do k = 1, kmax

        ! 飽和していたら, 温度と比湿の変化を計算
        ! Calculate tendency of temperature and humidity 
        ! if moist is saturation. 
        !
        if ( z_FlagProcess(k) ) then

          ! Liquid water and ice fractions
          call SaturateWatFraction( z_Temp(k), WatFrac )
          IceFrac = 1.0_DP - WatFrac


          if ( z_QH2OVapSat(k) < z_QH2OTot(k) - z_QH2OTotFluc(k) ) then
            CoefA =   z_QH2OTot(k) - z_QH2OVapSat(k)
            CoefB = - z_DQH2OVapSatDTemp(k)
          else if ( z_QH2OVapSat(k) < z_QH2OTot(k) + z_QH2OTotFluc(k) ) then
            CoefA =     ( z_QH2OTot(k) + z_QH2OTotFluc(k) - z_QH2OVapSat(k) )**2 / ( 4.0_DP * z_QH2OTotFluc(k) )
            CoefB = -   ( z_QH2OTot(k) + z_QH2OTotFluc(k) - z_QH2OVapSat(k) ) / ( 2.0_DP * z_QH2OTotFluc(k) ) * z_DQH2OVapSatDTemp(k)
          else
            CoefA = 0.0_DP
            CoefB = 0.0_DP
          end if

          ! 温度の変化分をニュートン法で求める
          ! Calculate variation of temperature
          !
          DelTemp = ( - LatentHeat * z_QH2OLiq(k) - ( LatentHeat + LatentHeatFusion ) * z_QH2OSol(k) + ( LatentHeat + LatentHeatFusion * IceFrac ) * CoefA ) / (   CpDry - ( LatentHeat + LatentHeatFusion * IceFrac ) * CoefB )


          ! 温度と比湿の調節
          ! Adjust temperature and specific humidity
          !
          z_Temp(k) = z_Temp(k) + DelTemp
          QH2OCond  = CoefA + CoefB * DelTemp
          z_QH2OVap(k) = z_QH2OTot(k) - QH2OCond
          z_QH2OLiq(k) = WatFrac * QH2OCond
          z_QH2OSol(k) = IceFrac * QH2OCond

          QH2OVapSat = z_QH2OVapSat(k) + z_DQH2OVapSatDTemp(k) * DelTemp

          ! It should be noted that, in general, z_QH2OVapSat(k) is not 
          ! equal to QH2OVapSat.
          if ( z_QH2OVapSat(k) < z_QH2OTot(k) - z_QH2OTotFluc(k) ) then
            z_CloudCover(k) = 1.0_DP
          else if ( z_QH2OVapSat(k) < z_QH2OTot(k) + z_QH2OTotFluc(k) ) then
            z_CloudCover(k) = ( z_QH2OTot(k) + z_QH2OTotFluc(k) - QH2OVapSat ) / ( 2.0_DP * z_QH2OTotFluc(k) )
            z_CloudCover(k) = max( min( z_CloudCover(k), 1.0_DP ), 0.0_DP )
          else
            z_CloudCover(k) = 0.0_DP
          end if
        else
          z_CloudCover(k) = 0.0_DP
        end if

      end do

    end do


    ! 比湿変化率, 温度変化率, 降水量の算出
    ! Calculate specific humidity tendency, temperature tendency, 
    ! precipitation
    !
    z_DTempDtLsc    = ( z_Temp    - z_Temp0    ) / ( 2.0_DP * DelTime )
    z_DQH2OVapDtLsc = ( z_QH2OVap - z_QH2OVap0 ) / ( 2.0_DP * DelTime )
    z_DQH2OLiqDt    = ( z_QH2OLiq - z_QH2OLiq0 ) / ( 2.0_DP * DelTime )
    z_DQH2OSolDt    = ( z_QH2OSol - z_QH2OSol0 ) / ( 2.0_DP * DelTime )


    call LScaleCond1DConsChk( z_Temp0, z_QH2OVap0, z_QH2OLiq0, z_QH2OSol0, z_Temp , z_QH2OVap , z_QH2OLiq , z_QH2OSol )


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

  end subroutine LScaleCondLL911D
QH2OFlucRatio
Variable :
QH2OFlucRatio :real(DP), save
lscond_inited
Variable :
lscond_inited = .false. :logical, save
: 初期設定フラグ. Initialization flag
module_name
Constant :
module_name = ‘lscond :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20150201 $’ // ’$Id: lscond.f90,v 1.21 2015/01/29 12:02:56 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version