Class | cloud_simple |
In: |
radiation/cloud_simple.f90
|
Note that Japanese and English are described in parallel.
簡単雲モデルによる雲の計算.
In this module, the amount of cloud is calculated by use of a simple cloud model.
!$ ! RadiationFluxDennouAGCM : | 放射フラックスの計算 |
!$ ! ———— : | ———— |
!$ ! RadiationFluxDennouAGCM : | Calculate radiation flux |
Subroutine : | |
xyr_Press( 0:imax-1, 1:jmax, 0:kmax ) : | real(DP), intent(in ) |
xyz_Temp( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_DQH2OLiqDtCum( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_DQH2OLiqDtLSC( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_QH2OLiq( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(inout) |
xyz_DTempDtCond( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(inout) |
xy_RainCum( 0:imax-1, 1:jmax ) : | real(DP), intent(out ) |
xy_SnowCum( 0:imax-1, 1:jmax ) : | real(DP), intent(out ) |
xy_RainLsc( 0:imax-1, 1:jmax ) : | real(DP), intent(out ) |
xy_SnowLsc( 0:imax-1, 1:jmax ) : | real(DP), intent(out ) |
subroutine CloudSimple( xyr_Press, xyz_Temp, xyz_DQH2OLiqDtCum, xyz_DQH2OLiqDtLSC, xyz_QH2OLiq, xyz_DTempDtCond, xy_RainCum, xy_SnowCum, xy_RainLsc, xy_SnowLsc ) ! USE statements ! ! 時刻管理 ! Time control ! use timeset, only: DelTime ! $ \Delta t $ [s] ! 雲関系ルーチン ! Cloud-related routines ! use cloud_utils, only : CloudUtilsCalcPRCPKeyLLTemp3D real(DP), intent(in ) :: xyr_Press ( 0:imax-1, 1:jmax, 0:kmax ) real(DP), intent(in ) :: xyz_Temp ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(in ) :: xyz_DQH2OLiqDtCum( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(in ) :: xyz_DQH2OLiqDtLSC( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(inout) :: xyz_QH2OLiq ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(inout) :: xyz_DTempDtCond ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(out ) :: xy_RainCum ( 0:imax-1, 1:jmax ) real(DP), intent(out ) :: xy_SnowCum ( 0:imax-1, 1:jmax ) real(DP), intent(out ) :: xy_RainLsc ( 0:imax-1, 1:jmax ) real(DP), intent(out ) :: xy_SnowLsc ( 0:imax-1, 1:jmax ) real(DP) :: xyz_QH2OLiqB( 0:imax-1, 1:jmax, 1:kmax ) real(DP) :: xyz_DQRainDt( 0:imax-1, 1:jmax, 1:kmax ) !!$ real(DP) :: xyz_DTempDtPrcpPCCum( 0:imax-1, 1:jmax, 1:kmax ) !!$ real(DP) :: xyz_DTempDtPrcpPCLsc( 0:imax-1, 1:jmax, 1:kmax ) ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. cloud_simple_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if ! Cloud optical depth ! ! Numerical solution !!$ xyz_DQCloudWaterDt = xyz_DQCloudWaterDt & !!$ & - xyz_QCloudWater / ( CloudLifeTime + 1.0d-100 ) !!$ ( X_{t+1} - X_{t-1} ) / ( 2 \Delta t ) = Q - X_{t+1} / \tau !!$ !!$ X_{t+1} / ( 2 \Delta t ) + X_{t+1} / \tau = X_{t-1} / ( 2 \Delta t ) + Q !!$ ( 1 / ( 2 \Delta t ) + 1 / \tau ) X_{t+1} = X_{t-1} / ( 2 \Delta t ) + Q !!$ X_{t+1} = ( X_{t-1} / ( 2 \Delta t ) + Q ) / ( 1 / ( 2 \Delta t ) + 1 / \tau ) xyz_QH2OLiq = ( xyz_QH2OLiq / ( 2.0_DP * DelTime ) + xyz_DQH2OLiqDtCum + xyz_DQH2OLiqDtLSC ) / ( 1.0_DP / ( 2.0_DP * DelTime ) + 1.0_DP / ( CloudLifeTime + 1.0d-100 ) ) call CloudUtilsCalcPRCPKeyLLTemp3D( xyr_Press, xyz_Temp, xyz_DQH2OLiqDtCum, xy_RainCum, xy_SnowCum ) call CloudUtilsCalcPRCPKeyLLTemp3D( xyr_Press, xyz_Temp, xyz_DQH2OLiqDtLsc, xy_RainLsc, xy_SnowLsc ) !!$ ! Analytical solution !!$ !!$ ! save cloud water before adjustment !!$ xyz_QH2OLiqB = xyz_QH2OLiq !!$ !!$ xyz_QH2OLiq = & !!$ & xyz_QH2OLiq * exp( - 2.0_DP * DelTime / ( CloudLifeTime + 1.0d-100 ) ) & !!$ & + ( xyz_DQH2OLiqDtCum + xyz_DQH2OLiqDtLSC ) * CloudLifeTime & !!$ & * ( 1.0_DP - exp( - 2.0_DP * DelTime / ( CloudLifeTime + 1.0d-100 ) ) ) !!$ !!$ xyz_DQRainDt = & !!$ & xyz_QH2OLiqB & !!$ & + ( xyz_DQH2OLiqDtCum + xyz_DQH2OLiqDtLSC ) * 2.0_DP * DelTime & !!$ & - xyz_QH2OLiq !!$ xyz_DQRainDt = xyz_DQRainDt / ( 2.0_DP * DelTime ) !!$ !!$ call CloudUtilsCalcPRCPKeyLLTemp( & !!$ & xyr_Press, xyz_Temp, xyz_DQRainDt, & ! (in ) !!$ & xy_RainCum, xy_SnowCum & ! (out) !!$ & ) !!$ xy_RainLsc = 0.0_DP !!$ xy_SnowLsc = 0.0_DP !!$ call CloudSimpleCalcPRCPWithPC( & !!$ & xyr_Press, xyz_TempA, xyz_DQH2OLiqDtCum, & ! (in ) !!$ & xy_RainCumulus, xy_SnowCumulus, & ! (out) !!$ & xyz_DTempDtPrcpPCCum & ! (out) !!$ & ) !!$ call CloudSimpleCalcPRCPWithPC( & !!$ & xyr_Press, xyz_TempA, xyz_DQH2OLiqDtLsc, & ! (in ) !!$ & xy_RainLsc, xy_SnowLsc, & ! (out) !!$ & xyz_DTempDtPrcpPCLsc & ! (out) !!$ & ) !!$ xyz_TempA = xyz_TempA & !!$ & + ( xyz_DTempDtPrcpPCCum + xyz_DTempDtPrcpPCLsc ) & !!$ & * 2.0_DP * DelTime !!$ xyz_DTempDtCond = xyz_DTempDtCond & !!$ & + xyz_DTempDtPrcpPCCum + xyz_DTempDtPrcpPCLsc end subroutine CloudSimple
Subroutine : | |
xyz_Press( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_Temp( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_QH2OTot( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_CloudCover( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(out) |
subroutine CloudSimpleCalcCloudCover( xyz_Press, xyz_Temp, xyz_QH2OTot, xyz_CloudCover ) ! USE statements ! ! 飽和比湿の算出 ! Evaluate saturation specific humidity ! use saturate, only: xyz_CalcQVapSat real(DP), intent(in ) :: xyz_Press ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(in ) :: xyz_Temp ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(in ) :: xyz_QH2OTot ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(out) :: xyz_CloudCover( 0:imax-1, 1:jmax, 1:kmax ) real(DP) :: xyz_RH(0:imax-1, 1:jmax, 1:kmax) ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. cloud_simple_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if select case ( IDCloudCoverMethod ) case ( IDCloudCoverMethodConst ) xyz_CloudCover = CloudCover case ( IDCloudCoverMethodRH ) ! see Sundqvist et al. (1989), Del Genio et al. (1996) xyz_RH = xyz_QH2OTot / xyz_CalcQVapSat( xyz_Temp, xyz_Press ) xyz_RH = min( xyz_RH, 1.0_DP ) xyz_CloudCover = 1.0d0 - sqrt( ( 1.0d0 - xyz_RH ) / ( 1.0d0 - RHCrtl ) ) xyz_CloudCover = max( xyz_CloudCover, 0.0_DP ) xyz_CloudCover = min( xyz_CloudCover, 1.0_DP ) end select end subroutine CloudSimpleCalcCloudCover
Subroutine : | |
xyr_Press( 0:imax-1, 1:jmax, 0:kmax ) : | real(DP), intent(in ) |
xyz_Temp( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xyz_DQH2OLiqDt( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) |
xy_SurfRainFlux( 0:imax-1, 1:jmax ) : | real(DP), intent(out) |
xy_SurfSnowFlux( 0:imax-1, 1:jmax ) : | real(DP), intent(out) |
xyz_DTempDt( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(out) |
subroutine CloudSimpleCalcPRCPWithPC( xyr_Press, xyz_Temp, xyz_DQH2OLiqDt, xy_SurfRainFlux, xy_SurfSnowFlux, xyz_DTempDt ) ! 時刻管理 ! Time control ! use timeset, only: DelTime ! $ \Delta t $ [s] ! 物理定数設定 ! Physical constants settings ! use constants, only: CpDry, Grav, LatentHeatFusion ! $ L $ [J kg-1] . ! 融解の潜熱. ! Latent heat of fusion ! 雲関系ルーチン ! Cloud-related routines ! use cloud_utils, only : CloudUtilsWatFraction real(DP), intent(in ) :: xyr_Press ( 0:imax-1, 1:jmax, 0:kmax ) real(DP), intent(in ) :: xyz_Temp ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(in ) :: xyz_DQH2OLiqDt ( 0:imax-1, 1:jmax, 1:kmax ) real(DP), intent(out) :: xy_SurfRainFlux( 0:imax-1, 1:jmax ) real(DP), intent(out) :: xy_SurfSnowFlux( 0:imax-1, 1:jmax ) real(DP), intent(out) :: xyz_DTempDt ( 0:imax-1, 1:jmax, 1:kmax ) ! 作業変数 ! Work variables ! real(DP) :: xyz_WatFrac( 0:imax-1, 1:jmax, 1:kmax ) real(DP) :: xy_Rain ( 0:imax-1, 1:jmax ) real(DP) :: xy_Snow ( 0:imax-1, 1:jmax ) real(DP) :: xy_Tot ( 0:imax-1, 1:jmax ) integer:: k ! 初期化確認 ! Initialization check ! if ( .not. cloud_simple_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if call CloudUtilsWatFraction( xyz_Temp, xyz_WatFrac ) xy_SurfRainFlux = 0.0_DP xy_SurfSnowFlux = 0.0_DP do k = kmax, 1, -1 xy_Rain = xy_SurfRainFlux + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav xy_Snow = xy_SurfSnowFlux xy_Tot = xy_Rain + xy_Snow xyz_DTempDt(:,:,k) = + ( xy_Rain - xyz_WatFrac(:,:,k) * xy_Tot ) * 2.0_DP * DelTime / ( ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav ) * LatentHeatFusion / CpDry / ( 2.0_DP * DelTime ) xy_SurfRainFlux = xyz_WatFrac(:,:,k) * xy_Tot xy_SurfSnowFlux = ( 1.0_DP - xyz_WatFrac(:,:,k) ) * xy_Tot end do end subroutine CloudSimpleCalcPRCPWithPC
Subroutine : | |
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in ) |
xyz_QH2OWatAndIce(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in ) |
xyz_QH2OWat(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out) |
xyz_QH2OIce(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out) |
subroutine CloudSimpleDivideWatAndIce( xyz_Temp, xyz_QH2OWatAndIce, xyz_QH2OWat, xyz_QH2OIce ) ! USE statements ! ! 雲関系ルーチン ! Cloud-related routines ! use cloud_utils, only : CloudUtilsWatFraction real(DP), intent(in ) :: xyz_Temp (0:imax-1, 1:jmax, 1:kmax) real(DP), intent(in ) :: xyz_QH2OWatAndIce(0:imax-1, 1:jmax, 1:kmax) real(DP), intent(out) :: xyz_QH2OWat (0:imax-1, 1:jmax, 1:kmax) real(DP), intent(out) :: xyz_QH2OIce (0:imax-1, 1:jmax, 1:kmax) real(DP) :: xyz_WatFrac(0:imax-1, 1:jmax, 1:kmax) ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. cloud_simple_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if call CloudUtilsWatFraction( xyz_Temp, xyz_WatFrac ) xyz_QH2OWat = xyz_QH2OWatAndIce * xyz_WatFrac xyz_QH2OIce = xyz_QH2OWatAndIce * ( 1.0_DP - xyz_WatFrac ) end subroutine CloudSimpleDivideWatAndIce
Subroutine : | |
ArgFlagSnow : | logical, intent(in) |
This procedure input/output NAMELIST#cloud_simple_nml .
subroutine CloudSimpleInit( ArgFlagSnow ) ! ファイル入出力補助 ! File I/O support ! use dc_iounit, only: FileOpen ! NAMELIST ファイル入力に関するユーティリティ ! Utilities for NAMELIST file input ! use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid ! ヒストリデータ出力 ! History data output ! use gtool_historyauto, only: HistoryAutoAddVariable ! 大規模凝結 (非対流性凝結) ! Large scale condensation (non-convective condensation) ! use lscond, only : LScaleCondInit ! 宣言文 ; Declaration statements ! logical, intent(in) :: ArgFlagSnow character(STRING) :: CloudCoverMethod integer:: unit_nml ! NAMELIST ファイルオープン用装置番号. ! Unit number for NAMELIST file open integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT. ! IOSTAT of NAMELIST read ! NAMELIST 変数群 ! NAMELIST group name ! namelist /cloud_simple_nml/ CloudLifeTime, CloudCoverMethod, RHCrtl, CloudCover ! ! デフォルト値については初期化手続 "cloud_simple#CloudSimpleInit" ! のソースコードを参照のこと. ! ! Refer to source codes in the initialization procedure ! "cloud_simple#CloudSimpleInit" for the default values. ! ! 実行文 ; Executable statement ! if ( cloud_simple_inited ) return FlagSnow = ArgFlagSnow ! デフォルト値の設定 ! Default values settings ! CloudLifeTime = 3600.0_DP CloudCoverMethod = 'Const' !!$ CloudCoverMethod = 'RH' RHCrtl = 0.8_DP CloudCover = 1.0_DP ! 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 = cloud_simple_nml, iostat = iostat_nml ) ! (out) close( unit_nml ) call NmlutilMsg( iostat_nml, module_name ) ! (in) end if select case ( CloudCoverMethod ) case ( 'Const' ) IDCloudCoverMethod = IDCloudCoverMethodConst case ( 'RH' ) IDCloudCoverMethod = IDCloudCoverMethodRH case default call MessageNotify( 'E', module_name, 'CloudCoverMethod=<%c> is not supported.', c1 = trim(CloudCoverMethod) ) end select ! Initialization of modules used in this module ! ! 大規模凝結 (非対流性凝結) (Manabe, 1965) ! Large scale condensation (non-convective condensation) (Le Treut and Li, 1991) ! call LScaleCondInit ! ヒストリデータ出力のためのへの変数登録 ! Register of variables for history data output ! !!$ call HistoryAutoAddVariable( 'EffCloudCover', & !!$ & (/ 'lon ', 'lat ', 'time' /), & !!$ & 'effective cloud cover', '1' ) ! 印字 ; Print ! call MessageNotify( 'M', module_name, '----- Initialization Messages -----' ) call MessageNotify( 'M', module_name, 'CloudLifeTime = %f', d = (/ CloudLifeTime /) ) call MessageNotify( 'M', module_name, 'CloudCoverMethod = %c', c1 = trim(CloudCoverMethod) ) call MessageNotify( 'M', module_name, 'RHCrtl = %f', d = (/ RHCrtl /) ) call MessageNotify( 'M', module_name, 'CloudCover = %f', d = (/ CloudCover /) ) call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) ) cloud_simple_inited = .true. end subroutine CloudSimpleInit
Variable : | |||
cloud_simple_inited = .false. : | logical, save
|
Constant : | |||
module_name = ‘cloud_simple‘ : | character(*), parameter
|
Constant : | |||
version = ’$Name: dcpam5-20120921 $’ // ’$Id: cloud_simple.f90,v 1.1 2012-09-08 15:16:39 yot Exp $’ : | character(*), parameter
|