| Class | saturation_adjust |
| In: |
lscond/saturation_adjust.f90
|
Note that Japanese and English are described in parallel.
Adjust temperature and specific humidity by the use of saturation adjustment
!$ ! 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.
| SaturationAdjust : | 温度と比湿の調節 |
| —————- : | ———— |
| SaturationAdjust : | Adjust temperature and specific humidity |
| Subroutine : | |||
| xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_QVap(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_QH2OLiq(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
| ||
| xyz_Press(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Press(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xyz_DQH2OLiqDt(0:imax-1,1:jmax,1:kmax) : | real(DP), intent(out) |
Adjust temperature and specific humidity by the use of saturation adjustment
subroutine SaturationAdjust( xyz_Temp, xyz_QVap, xyz_QH2OLiq, xyz_Press, xyr_Press, xyz_DQH2OLiqDt )
!
!
!
! Adjust temperature and specific humidity by the use of saturation
! adjustment
!
! モジュール引用 ; 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
! 大規模凝結 (非対流性凝結) (Manabe, 1965)
! Large scale condensation (non-convective condensation) (Le Treut and Li, 1991)
use lscond, only : LScaleCond
! 宣言文 ; 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_QH2OLiq(0:imax-1, 1:jmax, 1:kmax)
! $ q_w $ . 雲水混合比. Cloud water mixing ratio
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):: xyz_QVapSat (0:imax-1, 1:jmax, 1:kmax)
! 飽和比湿.
! Saturation specific humidity.
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_QH2OLiqB (0:imax-1, 1:jmax, 1:kmax)
! 調節前の雲水混合比.
! Mixing ratio of cloud water before adjust.
real(DP):: xyz_TempB (0:imax-1, 1:jmax, 1:kmax)
! 調節前の温度.
! Temperature before adjust.
!
real(DP):: xyz_RainLSC(0:imax-1, 1:jmax, 1:kmax)
real(DP):: xyz_EvapQH2OLiq(0:imax-1, 1:jmax, 1:kmax)
real(DP):: TempTentative
real(DP):: LatentHeatLocal
!
! Latent heat used in this routine
integer:: i
integer:: j
integer:: k ! 鉛直方向に回る DO ループ用作業変数
! Work variables for DO loop in vertical direction
! 実行文 ; Executable statement
!
! 初期化確認
! Initialization check
!
if ( .not. saturation_adjust_inited ) then
call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
end if
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! Set a value for latent heat
if ( FlagSublimation ) then
LatentHeatLocal = LatentHeat + LatentHeatFusion
else
LatentHeatLocal = LatentHeat
end if
! 調節前 "QVap", "QH2OLiq", "Temp" の保存
! Store "QVap", "QH2OLiq", "Temp" before adjustment
!
xyz_QVapB = xyz_QVap
xyz_QH2OLiqB = xyz_QH2OLiq
xyz_TempB = xyz_Temp
! 飽和比湿計算
! Calculate saturation specific humidity
!
xyz_QVapSat = xyz_CalcQVapSat( xyz_Temp, xyz_Press )
! Evaporate all of cloud liquid water temporarily
!!$ xyz_EvapQH2OLiq = xyz_QH2OLiq
! Evaporate part of cloud liquid water temporarily
do k = 1, kmax
do j = 1, jmax
do i = 0, imax-1
xyz_EvapQH2OLiq(i,j,k) = min( max( xyz_QVapSat(i,j,k), xyz_QVap(i,j,k) ) - xyz_QVap(i,j,k), xyz_QH2OLiq(i,j,k) )
xyz_EvapQH2OLiq(i,j,k) = xyz_EvapQH2OLiq(i,j,k) * ( 1.0_DP - 1.0e-10_DP )
TempTentative = xyz_Temp(i,j,k) - LatentHeatLocal * xyz_EvapQH2OLiq(i,j,k) / CpDry
if ( TempTentative < 1.0_DP ) then
TempTentative = 1.0_DP
xyz_EvapQH2OLiq(i,j,k) = ( xyz_Temp(i,j,k) - TempTentative ) / ( LatentHeatLocal / CpDry )
end if
end do
end do
end do
!
xyz_Temp = xyz_Temp - LatentHeatLocal * xyz_EvapQH2OLiq / CpDry
xyz_QVap = xyz_QVap + xyz_EvapQH2OLiq
! QH2OLiq will be updated in cloud model. Tendency will be updated below.
!!$ xyz_QH2OLiq = 0.0_DP
call LScaleCond( xyz_Temp, xyz_QVap, xyz_Press, xyr_Press, xyz_DQH2OLiqDt, .false. )
! 比湿変化率, 温度変化率, 降水量の算出
! 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_DQH2OLiqDt - xyz_EvapQH2OLiq / ( 2.0_DP * DelTime )
! 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
!
call HistoryAutoPut( TimeN, 'RainLsc', xy_RainLsc )
call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine SaturationAdjust
| Subroutine : | |
| FlagSnow : | logical, intent(in) |
saturation_adjust モジュールの初期化を行います. NAMELIST#saturation_adjust_nml の読み込みはこの手続きで行われます.
"saturation_adjust" module is initialized. "NAMELIST#saturation_adjust_nml" is loaded in this procedure.
This procedure input/output NAMELIST#saturation_adjust_nml .
subroutine SaturationAdjustInit( FlagSnow )
!
! saturation_adjust モジュールの初期化を行います.
! NAMELIST#saturation_adjust_nml の読み込みはこの手続きで行われます.
!
! "saturation_adjust" module is initialized.
! "NAMELIST#saturation_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
! 大規模凝結 (非対流性凝結) (Manabe, 1965)
! Large scale condensation (non-convective condensation) (Le Treut and Li, 1991)
use lscond, only : LScaleCondInit
! 宣言文 ; 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 /saturation_adjust_nml/ FlagSublimation
!
! デフォルト値については初期化手続 "saturation_adjust#SaturationAdjustInit"
! のソースコードを参照のこと.
!
! Refer to source codes in the initialization procedure
! "saturation_adjust#SaturationAdjustInit" for the default values.
!
! 実行文 ; Executable statement
!
if ( saturation_adjust_inited ) return
! デフォルト値の設定
! Default values settings
!
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 = saturation_adjust_nml, iostat = iostat_nml ) ! (out)
close( unit_nml )
call NmlutilMsg( iostat_nml, module_name ) ! (in)
if ( iostat_nml == 0 ) write( STDOUT, nml = saturation_adjust_nml )
end if
! Initialization of modules used in this routine
! 飽和比湿の算出
! Evaluate saturation specific humidity
!
call SaturateInit
! 大規模凝結 (非対流性凝結) (Manabe, 1965)
! Large scale condensation (non-convective condensation) (Le Treut and Li, 1991)
!
call LScaleCondInit( FlagSnow )
! ヒストリデータ出力のためのへの変数登録
! Register of variables for history data output
!
! HistoryAutoAddVariable is called in LScaleCond for variables below.
!
!!$ 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, ' FlagSublimation = %b', l = (/ FlagSublimation /) )
call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
saturation_adjust_inited = .true.
end subroutine SaturationAdjustInit
| Constant : | |||
| module_name = ‘saturation_adjust‘ : | character(*), parameter
|
| Variable : | |||
| saturation_adjust_inited = .false. : | logical, save
|
| Constant : | |||
| version = ’$Name: $’ // ’$Id: saturation_adjust.f90,v 1.4 2015/01/29 12:03:22 yot Exp $’ : | character(*), parameter
|