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

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_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)
: 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_DTempDt, xyz_DQVapDt, 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(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)
                              !
                              ! 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_DTempDt = xyz_DTempDt + xyz_DTempDtLsc
    xyz_DQVapDt = xyz_DQVapDt + xyz_DQVapDtLsc


    xyz_DQH2OLiqDt = - xyz_DQVapDtLsc


    ! calculation for output
    xy_RainLsc     = 0.0d0
    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, 'RainLsc',    xy_RainLsc     )
      call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
      call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
    else
      if ( FlagOutput ) then
        call HistoryAutoPut( TimeN, 'RainLsc',    xy_RainLsc     )
        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 :

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
    !
    ! 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

    ! 宣言文 ; 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 /lscond_nml/ CrtlRH, ItrtMax, FlagSublimation
          !
          ! デフォルト値については初期化手続 "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
    !
    CrtlRH          = 1.0_DP
    ItrtMax         = 3
    FlagSublimation = .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

    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
    call HistoryAutoAddVariable( 'RainLsc', (/ 'lon ', 'lat ', 'time' /), 'precipitation by large scale condensation', 'kg m-2 s-1' )
    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, '  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, '-- version = %c', c1 = trim(version) )

    lscond_inited = .true.

  end subroutine LScaleCondInit

Private Instance methods

CrtlRH
Variable :
CrtlRH :real(DP), save
: 臨界相対湿度. Critical relative humidity
FlagSublimation
Variable :
FlagSublimation :logical, save
: flag for treating sublimation
ItrtMax
Variable :
ItrtMax :integer, save
: イテレーション回数. Number of iteration
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-20140204-5 $’ // ’$Id: lscond.f90,v 1.16 2013-09-30 03:02:14 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version