Class | phy_implicit |
In: |
phy_implicit/phy_implicit.f90
|
Note that Japanese and English are described in parallel.
PhyImplGetMatrices : | 陰解行列の作成と取得 |
PhyImplTendency : | 時間変化率の計算 |
PhyImplFluxCorrect : | フラックスの補正 |
———— : | ———— |
PhyImplGetMatrices : | Create and get implicit matrices |
PhyImplTendency : | Calculate tendency |
PhyImplFluxCorrect : | Correct fluxes |
Subroutine : | |||
xyz_DUDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
xyz_DVDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
xyz_DQVapDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
xy_DSurfTempDt(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyza_UVMtx(0:imax-1, 1:jmax, 1:kmax, -1:1) : | real(DP), intent(in)
| ||
xyra_TempMtx(0:imax-1, 1:jmax, 0:kmax, -1:1) : | real(DP), intent(in)
| ||
xyza_QVapMtx(0:imax-1, 1:jmax, 1:kmax, -1:1) : | real(DP), intent(in)
| ||
xy_SurfUVMtx(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyaa_SurfTempMtx(0:imax-1, 1:jmax, 0:1, -1:1) : | real(DP), intent(in) | ||
xyaa_SurfQVapMtx(0:imax-1, 1:jmax, 0:1, -1:1) : | real(DP), intent(in)
| ||
xyr_UFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(inout)
| ||
xyr_VFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(inout)
| ||
xyr_TempFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(inout)
| ||
xyr_QVapFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(inout)
|
フラックスの補正.
Correct fluxes.
subroutine PhyImplFluxCorrect( xyz_DUDt, xyz_DVDt, xyz_DTempDt, xyz_DQVapDt, xy_DSurfTempDt, xyza_UVMtx, xyra_TempMtx, xyza_QVapMtx, xy_SurfUVMtx, xyaa_SurfTempMtx, xyaa_SurfQVapMtx, xyr_UFlux, xyr_VFlux, xyr_TempFlux, xyr_QVapFlux ) ! ! フラックスの補正. ! ! Correct fluxes. ! ! モジュール引用 ; USE statements ! ! 物理定数設定 ! Physical constants settings ! use constants, only: CpDry, LatentHeat ! $ L $ [J kg-1] . ! 凝結の潜熱. ! Latent heat of condensation ! 時刻管理 ! Time control ! use timeset, only: DelTime, TimeN, TimesetClockStart, TimesetClockStop ! ヒストリデータ出力 ! History data output ! use gtool_historyauto, only: HistoryAutoPut ! 宣言文 ; Declaration statements ! implicit none real(DP), intent(in):: xyz_DUDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{u}{t} $ . 東西風速変化. ! Eastward wind tendency real(DP), intent(in):: xyz_DVDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{v}{t} $ . 南北風速変化. ! Northward wind tendency real(DP), intent(in):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{T}{t} $ . 温度変化. ! Temperature tendency real(DP), intent(in):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{q}{t} $ . 比湿変化. ! Temperature tendency real(DP), intent(in):: xy_DSurfTempDt (0:imax-1, 1:jmax) ! 地表面温度変化率. ! Surface temperature tendency real(DP), intent(in):: xyza_UVMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 速度陰解行列. ! Implicit matrix about velocity real(DP), intent(in):: xyra_TempMtx (0:imax-1, 1:jmax, 0:kmax, -1:1) ! 温度陰解行列. ! Implicit matrix about temperature real(DP), intent(in):: xyza_QVapMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 比湿陰解行列. ! Implicit matrix about specific humidity real(DP), intent(in):: xy_SurfUVMtx (0:imax-1, 1:jmax) ! 速度陰解行列: 地表. ! Implicit matrix about velocity: surface real(DP), intent(in):: xyaa_SurfTempMtx (0:imax-1, 1:jmax, 0:1, -1:1) ! 温度陰解行列: 地表. ! Implicit matrix about temperature: surface real(DP), intent(in):: xyaa_SurfQVapMtx (0:imax-1, 1:jmax, 0:1, -1:1) ! 比湿陰解行列: 地表. ! Implicit matrix about specific humidity: surface real(DP), intent(inout):: xyr_UFlux (0:imax-1, 1:jmax, 0:kmax) ! 東西風速フラックス. ! Eastward wind flux real(DP), intent(inout):: xyr_VFlux (0:imax-1, 1:jmax, 0:kmax) ! 南北風速フラックス. ! Northward wind flux real(DP), intent(inout):: xyr_TempFlux (0:imax-1, 1:jmax, 0:kmax) ! 温度フラックス. ! Temperature flux real(DP), intent(inout):: xyr_QVapFlux (0:imax-1, 1:jmax, 0:kmax) ! 比湿フラックス. ! Specific humidity flux ! 作業変数 ! Work variables ! real(DP):: xyr_UFluxB (0:imax-1, 1:jmax, 0:kmax) ! 補正前の東西風速フラックス. ! Eastward wind flux before adjust real(DP):: xyr_VFluxB (0:imax-1, 1:jmax, 0:kmax) ! 補正前の南北風速フラックス. ! Northward wind flux before adjust real(DP):: xyr_TempFluxB (0:imax-1, 1:jmax, 0:kmax) ! 補正前の温度フラックス. ! Temperature flux before adjust real(DP):: xyr_QVapFluxB (0:imax-1, 1:jmax, 0:kmax) ! 補正前の比湿フラックス. ! Specific humidity flux before adjust real(DP):: xyr_DUFluxDt (0:imax-1, 1:jmax, 0:kmax) ! 東西風速フラックスの変化率. ! Eastward wind flux tendency real(DP):: xyr_DVFluxDt (0:imax-1, 1:jmax, 0:kmax) ! 南北風速フラックスの変化率. ! Northward wind flux tendency real(DP):: xyr_DTempFluxDt (0:imax-1, 1:jmax, 0:kmax) ! 温度フラックスの変化率. ! Temperature flux tendency real(DP):: xyr_DQVapFluxDt (0:imax-1, 1:jmax, 0:kmax) ! 比湿フラックスの変化率. ! Specific humidity flux tendency 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 real(DP):: LCp ! $ L / C_p $ [K]. ! 実行文 ; Executable statement ! ! 計算時間計測開始 ! Start measurement of computation time ! call TimesetClockStart( module_name ) ! 初期化 ! Initialization ! if ( .not. phy_implicit_inited ) call PhyImplInit ! 調節前のフラックスの保存 ! Save fluxes before adjustment ! xyr_UFluxB = xyr_UFlux xyr_VFluxB = xyr_VFlux xyr_TempFluxB = xyr_TempFlux xyr_QVapFluxB = xyr_QVapFlux ! フラックス補正 ! Correct fluxes ! LCp = LatentHeat / CpDry do k = 1, kmax-1 do j = 1, jmax do i = 0, imax-1 xyr_UFlux( i,j,k ) = xyr_UFlux( i,j,k ) - ( xyza_UVMtx( i,j,k+1,-1) * xyz_DUDt( i,j,k ) - xyza_UVMtx( i,j,k , 1) * xyz_DUDt( i,j,k+1 ) ) * DelTime xyr_VFlux( i,j,k ) = xyr_VFlux( i,j,k ) - ( xyza_UVMtx( i,j,k+1,-1) * xyz_DVDt( i,j,k ) - xyza_UVMtx( i,j,k , 1) * xyz_DVDt( i,j,k+1 ) ) * DelTime xyr_TempFlux( i,j,k ) = xyr_TempFlux( i,j,k ) - ( xyra_TempMtx( i,j,k+1,-1) * xyz_DTempDt( i,j,k ) - xyra_TempMtx( i,j,k , 1) * xyz_DTempDt( i,j,k+1 ) ) * DelTime xyr_QVapFlux( i,j,k ) = xyr_QVapFlux( i,j,k ) - ( xyza_QVapMtx( i,j,k+1,-1) * xyz_DQVapDt( i,j,k ) - xyza_QVapMtx( i,j,k , 1) * xyz_DQVapDt( i,j,k+1 ) ) * DelTime * LCp end do end do end do do j = 1, jmax do i = 0, imax-1 xyr_UFlux( i,j,0 ) = xyr_UFlux( i,j,0 ) - xy_SurfUVMtx( i,j ) * xyz_DUDt( i,j,1 ) * DelTime xyr_VFlux( i,j,0 ) = xyr_VFlux( i,j,0 ) - xy_SurfUVMtx( i,j ) * xyz_DVDt( i,j,1 ) * DelTime xyr_TempFlux( i,j,0 ) = xyr_TempFlux( i,j,0 ) - ( xyaa_SurfTempMtx( i,j,1,-1 ) * xy_DSurfTempDt( i,j ) + xyaa_SurfTempMtx( i,j,1,0 ) * xyz_DTempDt ( i,j,1 ) ) * DelTime xyr_QVapFlux( i,j,0 ) = xyr_QVapFlux( i,j,0 ) - ( xyaa_SurfQVapMtx( i,j,1,-1 ) * xy_DSurfTempDt( i,j ) + xyaa_SurfQVapMtx( i,j,1,0 ) * xyz_DQVapDt ( i,j,1 ) * LCp ) * DelTime end do end do ! 変化率の算出 ! Calculate tendencies ! xyr_DUFluxDt = ( xyr_UFlux - xyr_UFluxB ) / ( 2.0_DP * DelTime ) xyr_DVFluxDt = ( xyr_VFlux - xyr_VFluxB ) / ( 2.0_DP * DelTime ) xyr_DTempFluxDt = ( xyr_TempFlux - xyr_TempFluxB ) / ( 2.0_DP * DelTime ) xyr_DQVapFluxDt = ( xyr_QVapFlux - xyr_QVapFluxB ) / ( 2.0_DP * DelTime ) ! ヒストリデータ出力 ! History data output ! call HistoryAutoPut( TimeN, 'DUFluxDtPhyImplFluxCor', xyr_DUFluxDt ) call HistoryAutoPut( TimeN, 'DVFluxDtPhyImplFluxCor', xyr_DVFluxDt ) call HistoryAutoPut( TimeN, 'DTempFluxDtPhyImplFluxCor', xyr_DTempFluxDt ) call HistoryAutoPut( TimeN, 'DQVapFluxDtPhyImplFluxCor', xyr_DQVapFluxDt ) ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine PhyImplFluxCorrect
Subroutine : | |||
xyr_Press(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xy_SurfHeatCapacity(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfCond(0:imax-1, 1:jmax) : | integer, intent(in)
| ||
xyr_UFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_VFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_TempFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_QVapFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyza_UVMtx(0:imax-1, 1:jmax, 1:kmax, -1:1) : | real(DP), intent(out)
| ||
xyra_TempMtx(0:imax-1, 1:jmax, 0:kmax, -1:1) : | real(DP), intent(out)
| ||
xyza_QVapMtx(0:imax-1, 1:jmax, 1:kmax, -1:1) : | real(DP), intent(out)
|
陰解行列の作成と取得を行います.
Create and get implicit matrices.
subroutine PhyImplGetMatrices( xyr_Press, xy_SurfHeatCapacity, xy_SurfCond, xyr_UFlux, xyr_VFlux, xyr_TempFlux, xyr_QVapFlux, xyza_UVMtx, xyra_TempMtx, xyza_QVapMtx ) ! ! 陰解行列の作成と取得を行います. ! ! Create and get implicit matrices. ! ! モジュール引用 ; USE statements ! ! 物理定数設定 ! Physical constants settings ! use constants, only: Grav, CpDry ! $ C_p $ [J kg-1 K-1]. ! 乾燥大気の定圧比熱. ! Specific heat of air at constant pressure ! 時刻管理 ! Time control ! use timeset, only: DelTime, TimesetClockStart, TimesetClockStop ! 宣言文 ; Declaration statements ! implicit none real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax) ! $ \hat{p} $ . 気圧 (半整数レベル). ! Air pressure (half level) real(DP), intent(in):: xy_SurfHeatCapacity (0:imax-1, 1:jmax) ! 地表熱容量. ! Surface heat capacity integer, intent(in):: xy_SurfCond (0:imax-1, 1:jmax) ! 地表状態. ! Surface condition real(DP), intent(out):: xyr_UFlux (0:imax-1, 1:jmax, 0:kmax) ! 東西風速フラックス. ! Eastward wind flux real(DP), intent(out):: xyr_VFlux (0:imax-1, 1:jmax, 0:kmax) ! 南北風速フラックス. ! Northward wind flux real(DP), intent(out):: xyr_TempFlux (0:imax-1, 1:jmax, 0:kmax) ! 温度フラックス. ! Temperature flux real(DP), intent(out):: xyr_QVapFlux (0:imax-1, 1:jmax, 0:kmax) ! 比湿フラックス. ! Specific humidity flux real(DP), intent(out):: xyza_UVMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 速度陰解行列. ! Implicit matrix about velocity real(DP), intent(out):: xyra_TempMtx (0:imax-1, 1:jmax, 0:kmax, -1:1) ! 温度陰解行列. ! Implicit matrix about temperature real(DP), intent(out):: xyza_QVapMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 比湿陰解行列. ! Implicit matrix about specific humidity ! 作業変数 ! Work variables ! 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 ! ! 計算時間計測開始 ! Start measurement of computation time ! call TimesetClockStart( module_name ) ! 初期化 ! Initialization ! if ( .not. phy_implicit_inited ) call PhyImplInit ! 質量, 熱容量の項の計算 ! Calculate mass and heat capacity terms ! xyza_UVMtx = 0.0_DP xyra_TempMtx = 0.0_DP xyza_QVapMtx = 0.0_DP do k = 1, kmax xyza_UVMtx(:,:,k,0) = ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav / ( 2.0_DP * DelTime ) xyra_TempMtx(:,:,k,0) = xyza_UVMtx(:,:,k,0) * CpDry xyza_QVapMtx(:,:,k,0) = xyza_UVMtx(:,:,k,0) * CpDry end do do j = 1, jmax do i = 0, imax-1 if ( xy_SurfCond(i,j) >= 1 ) then xyra_TempMtx(i,j,0,0) = xy_SurfHeatCapacity(i,j) / ( 2.0_DP * DelTime ) else xyra_TempMtx(i,j,0,0) = 1.0_DP end if end do end do ! フラックスをリセット ! Reset fluxes ! xyr_UFlux = 0.0_DP xyr_VFlux = 0.0_DP xyr_TempFlux = 0.0_DP xyr_QVapFlux = 0.0_DP ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine PhyImplGetMatrices
Subroutine : | |||
xyr_UFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_VFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_TempFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_QVapFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xy_SurfRadSFlux(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfRadLFlux(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_GroundTempFlux(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyza_UVMtx(0:imax-1, 1:jmax, 1:kmax, -1:1) : | real(DP), intent(in)
| ||
xyra_TempMtx(0:imax-1, 1:jmax, 0:kmax, -1:1) : | real(DP), intent(in)
| ||
xyza_QVapMtx(0:imax-1, 1:jmax, 1:kmax, -1:1) : | real(DP), intent(in)
| ||
xy_SurfUVMtx(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyaa_SurfTempMtx(0:imax-1, 1:jmax, 0:1, -1:1) : | real(DP), intent(in) | ||
xyaa_SurfQVapMtx(0:imax-1, 1:jmax, 0:1, -1:1) : | real(DP), intent(in)
| ||
xya_SurfRadLMtx(0:imax-1, 1:jmax, -1:1) : | real(DP), intent(in)
| ||
xy_SurfCond(0:imax-1, 1:jmax) : | integer, intent(in)
| ||
xyz_DUDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
xyz_DVDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
xyz_DQVapDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
xy_DSurfTempDt(0:imax-1, 1:jmax) : | real(DP), intent(out)
|
時間変化率の計算を行います.
Calculate tendencies.
subroutine PhyImplTendency( xyr_UFlux, xyr_VFlux, xyr_TempFlux, xyr_QVapFlux, xy_SurfRadSFlux, xy_SurfRadLFlux, xy_GroundTempFlux, xyza_UVMtx, xyra_TempMtx, xyza_QVapMtx, xy_SurfUVMtx, xyaa_SurfTempMtx, xyaa_SurfQVapMtx, xya_SurfRadLMtx, xy_SurfCond, xyz_DUDt, xyz_DVDt, xyz_DTempDt, xyz_DQVapDt, xy_DSurfTempDt ) ! ! 時間変化率の計算を行います. ! ! Calculate tendencies. ! ! モジュール引用 ; USE statements ! ! 物理定数設定 ! Physical constants settings ! use constants, only: Grav, CpDry, LatentHeat ! $ L $ [J kg-1] . ! 凝結の潜熱. ! Latent heat of condensation ! 時刻管理 ! Time control ! use timeset, only: DelTime, TimeN, TimesetClockStart, TimesetClockStop ! 時刻管理 ! Time control ! use timeset, only: TimeN ! ステップ $ t $ の時刻. Time of step $ t $. ! ヒストリデータ出力 ! History data output ! use gtool_historyauto, only: HistoryAutoPut ! 宣言文 ; Declaration statements ! implicit none real(DP), intent(in):: xyr_UFlux (0:imax-1, 1:jmax, 0:kmax) ! 東西風速フラックス. ! Eastward wind flux real(DP), intent(in):: xyr_VFlux (0:imax-1, 1:jmax, 0:kmax) ! 南北風速フラックス. ! Northward wind flux real(DP), intent(in):: xyr_TempFlux (0:imax-1, 1:jmax, 0:kmax) ! 温度フラックス. ! Temperature flux real(DP), intent(in):: xyr_QVapFlux (0:imax-1, 1:jmax, 0:kmax) ! 比湿フラックス. ! Specific humidity flux real(DP), intent(in):: xy_SurfRadSFlux (0:imax-1, 1:jmax) ! 短波 (日射) フラックス (地表面). ! Shortwave (insolation) flux on surface real(DP), intent(in):: xy_SurfRadLFlux (0:imax-1, 1:jmax) ! 長波フラックス (地表面). ! Longwave flux on surface real(DP), intent(in):: xy_GroundTempFlux (0:imax-1, 1:jmax) ! 地中熱フラックス. ! Ground temperature flux real(DP), intent(in):: xyza_UVMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 速度陰解行列. ! Implicit matrix about velocity real(DP), intent(in):: xyra_TempMtx (0:imax-1, 1:jmax, 0:kmax, -1:1) ! 温度陰解行列. ! Implicit matrix about temperature real(DP), intent(in):: xyza_QVapMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 比湿陰解行列. ! Implicit matrix about specific humidity real(DP), intent(in):: xy_SurfUVMtx (0:imax-1, 1:jmax) ! 速度陰解行列: 地表. ! Implicit matrix about velocity: surface real(DP), intent(in):: xyaa_SurfTempMtx (0:imax-1, 1:jmax, 0:1, -1:1) ! 温度陰解行列: 地表. ! Implicit matrix about temperature: surface real(DP), intent(in):: xyaa_SurfQVapMtx (0:imax-1, 1:jmax, 0:1, -1:1) ! 比湿陰解行列: 地表. ! Implicit matrix about specific humidity: surface real(DP), intent(in):: xya_SurfRadLMtx (0:imax-1, 1:jmax, -1:1) ! $ T $ . 陰解行列: 地表. ! implicit matrix: surface integer, intent(in):: xy_SurfCond (0:imax-1, 1:jmax) ! 地表状態. ! Surface condition real(DP), intent(out):: xyz_DUDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{u}{t} $ . 東西風速変化. ! Eastward wind tendency real(DP), intent(out):: xyz_DVDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{v}{t} $ . 南北風速変化. ! Northward wind tendency real(DP), intent(out):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{T}{t} $ . 温度変化. ! Temperature tendency real(DP), intent(out):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{q}{t} $ . 比湿変化. ! Temperature tendency real(DP), intent(out):: xy_DSurfTempDt (0:imax-1, 1:jmax) ! 地表面温度変化率. ! Surface temperature tendency ! 作業変数 ! Work variables ! real(DP):: xyz_DelTempQVap (0:imax-1, 1:jmax, -kmax:kmax) ! $ T q $ の時間変化. ! Tendency of $ T q $ real(DP):: xyza_TempQVapLUMtx (0:imax-1, 1:jmax, -kmax:kmax, -1:1) ! LU 行列. ! LU matrix real(DP):: xyza_UVLUMtx (0:imax-1, 1:jmax, 1:kmax,-1:1) ! LU 行列. ! LU matrix 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:: l ! 行列用 DO ループ用作業変数 ! Work variables for DO loop of matrices ! 実行文 ; Executable statement ! ! 計算時間計測開始 ! Start measurement of computation time ! call TimesetClockStart( module_name ) ! 初期化 ! Initialization ! if ( .not. phy_implicit_inited ) call PhyImplInit ! 東西風速, 南北風速の計算 ! Calculate eastward and northward wind ! xyza_UVLUMtx = xyza_UVMtx xyza_UVLUMtx(:,:,1,0) = xyza_UVLUMtx(:,:,1,0) + xy_SurfUVMtx call PhyImplLUDecomp3( xyza_UVLUMtx, imax * jmax, kmax ) ! (in) do k = 1, kmax xyz_DUDt(:,:,k) = xyr_UFlux(:,:,k-1) - xyr_UFlux(:,:,k) xyz_DVDt(:,:,k) = xyr_VFlux(:,:,k-1) - xyr_VFlux(:,:,k) end do call PhyImplLUSolve3( xyz_DUDt, xyza_UVLUMtx, 1, imax * jmax, kmax ) ! (in) call PhyImplLUSolve3( xyz_DVDt, xyza_UVLUMtx, 1, imax * jmax, kmax ) ! (in) ! 温度と比湿の計算 ! Calculate temperature and specific humidity ! do l = -1, 1 do k = 1, kmax xyza_TempQVapLUMtx(:,:,k,l) = xyra_TempMtx(:,:,k,l) xyza_TempQVapLUMtx(:,:,-k,-l) = xyza_QVapMtx(:,:,k,l) end do xyza_TempQVapLUMtx(:,:, 1, l) = xyra_TempMtx(:,:,1,l) + xyaa_SurfTempMtx(:,:,1,l) xyza_TempQVapLUMtx(:,:,-1,-l) = xyza_QVapMtx(:,:,1,l) + xyaa_SurfQVapMtx(:,:,1,l) end do xyza_TempQVapLUMtx(:,:,0, 0) = xyra_TempMtx(:,:,0,0) + xyaa_SurfTempMtx(:,:,0,0) + xyaa_SurfQVapMtx(:,:,0,0) + xya_SurfRadLMtx(:,:,0) xyza_TempQVapLUMtx(:,:,0, 1) = xyaa_SurfTempMtx(:,:,0,1) + xya_SurfRadLMtx(:,:,1) xyza_TempQVapLUMtx(:,:,0,-1) = xyaa_SurfQVapMtx(:,:,0,1) call PhyImplLUDecomp3( xyza_TempQVapLUMtx, imax * jmax, (2 * kmax) + 1 ) ! (in) do k = 1, kmax xyz_DelTempQVap(:,:, k) = xyr_TempFlux(:,:,k-1) - xyr_TempFlux(:,:,k) xyz_DelTempQVap(:,:,-k) = xyr_QVapFlux(:,:,k-1) - xyr_QVapFlux(:,:,k) end do xyz_DelTempQVap(:,:,0) = - xy_SurfRadSFlux - xy_SurfRadLFlux - xyr_TempFlux(:,:,0) - xyr_QVapFlux(:,:,0) + xy_GroundTempFlux call PhyImplLUSolve3( xyz_DelTempQVap, xyza_TempQVapLUMtx, 1, imax * jmax , (2 * kmax) + 1 ) ! (in) ! 時間変化率の計算 ! Calculate tendency ! do k = 1, kmax xyz_DUDt(:,:,k) = xyz_DUDt(:,:,k) / ( 2.0_DP * DelTime ) xyz_DVDt(:,:,k) = xyz_DVDt(:,:,k) / ( 2.0_DP * DelTime ) xyz_DTempDt(:,:,k) = xyz_DelTempQVap(:,:, k) / ( 2.0_DP * DelTime ) xyz_DQVapDt(:,:,k) = xyz_DelTempQVap(:,:,-k) / ( 2.0_DP * DelTime ) / LatentHeat * CpDry end do do j = 1, jmax do i = 0, imax-1 if ( xy_SurfCond(i,j) >= 1 ) then xy_DSurfTempDt(i,j) = xyz_DelTempQVap(i,j,0) / ( 2.0_DP * DelTime ) else xy_DSurfTempDt(i,j) = 0.0_DP end if end do end do ! ヒストリデータ出力 ! History data output ! call HistoryAutoPut( TimeN, 'DUDtVDiff', xyz_DUDt ) call HistoryAutoPut( TimeN, 'DVDtVDiff', xyz_DVDt ) call HistoryAutoPut( TimeN, 'DTempDtVDiff', xyz_DTempDt ) call HistoryAutoPut( TimeN, 'DQVapDtVDiff', xyz_DQVapDt ) ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine PhyImplTendency
Variable : | |||
phy_implicit_inited = .false. : | logical, save, public
|
Subroutine : |
依存モジュールの初期化チェック
Check initialization of dependency modules
subroutine InitCheck ! ! 依存モジュールの初期化チェック ! ! Check initialization of dependency modules ! モジュール引用 ; USE statements ! ! NAMELIST ファイル入力に関するユーティリティ ! Utilities for NAMELIST file input ! use namelist_util, only: namelist_util_inited ! 格子点設定 ! Grid points settings ! use gridset, only: gridset_inited ! 物理定数設定 ! Physical constants settings ! use constants, only: constants_inited ! 座標データ設定 ! Axes data settings ! use axesset, only: axesset_inited ! 時刻管理 ! Time control ! use timeset, only: timeset_inited ! 実行文 ; Executable statement ! if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' ) if ( .not. gridset_inited ) call MessageNotify( 'E', module_name, '"gridset" module is not initialized.' ) if ( .not. constants_inited ) call MessageNotify( 'E', module_name, '"constants" module is not initialized.' ) if ( .not. axesset_inited ) call MessageNotify( 'E', module_name, '"axesset" module is not initialized.' ) if ( .not. timeset_inited ) call MessageNotify( 'E', module_name, '"timeset" module is not initialized.' ) end subroutine InitCheck
Subroutine : |
phy_implicit モジュールの初期化を行います. NAMELIST#phy_implicit_nml の読み込みはこの手続きで行われます.
"phy_implicit" module is initialized. "NAMELIST#phy_implicit_nml" is loaded in this procedure.
subroutine PhyImplInit ! ! phy_implicit モジュールの初期化を行います. ! NAMELIST#phy_implicit_nml の読み込みはこの手続きで行われます. ! ! "phy_implicit" module is initialized. ! "NAMELIST#phy_implicit_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 /phy_implicit_nml/ & ! ! デフォルト値については初期化手続 "phy_implicit#PhyImplInit" ! のソースコードを参照のこと. ! ! Refer to source codes in the initialization procedure ! "phy_implicit#PhyImplInit" for the default values. ! ! 実行文 ; Executable statement ! if ( phy_implicit_inited ) return call InitCheck ! デフォルト値の設定 ! Default values settings ! !!$ ! NAMELIST の読み込み !!$ ! NAMELIST is input !!$ ! !!$ if ( trim(namelist_filename) /= '' ) then !!$ call FileOpen( unit_nml, & ! (out) !!$ & namelist_filename, mode = 'r' ) ! (in) !!$ !!$ rewind( unit_nml ) !!$ read( unit_nml, & ! (in) !!$ & nml = phy_implicit_nml, & ! (out) !!$ & iostat = iostat_nml ) ! (out) !!$ close( unit_nml ) !!$ !!$ call NmlutilMsg( iostat_nml, module_name ) ! (in) !!$ end if ! ヒストリデータ出力のためのへの変数登録 ! Register of variables for history data output ! call HistoryAutoAddVariable( 'DUDtVDiff', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'diffusive force(x)', 'm s-2' ) call HistoryAutoAddVariable( 'DVDtVDiff', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'diffusive force(y)', 'm s-2' ) call HistoryAutoAddVariable( 'DTempDtVDiff', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'diffusive heating', 'K s-1' ) call HistoryAutoAddVariable( 'DQVapDtVDiff', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'diffusive moistening', 'kg kg-1 s-1' ) call HistoryAutoAddVariable( 'DUFluxDtPhyImplFluxCor', (/ 'lon ', 'lat ', 'sigm', 'time' /), 'eastward wind flux tendency', 'N m-2 s-1' ) call HistoryAutoAddVariable( 'DVFluxDtPhyImplFluxCor', (/ 'lon ', 'lat ', 'sigm', 'time' /), 'northward wind flux tendency', 'N m-2 s-1' ) call HistoryAutoAddVariable( 'DTempFluxDtPhyImplFluxCor', (/ 'lon ', 'lat ', 'sigm', 'time' /), 'temperature flux tendency', 'W m-2 s-1' ) call HistoryAutoAddVariable( 'DQVapFluxDtPhyImplFluxCor', (/ 'lon ', 'lat ', 'sigm', 'time' /), 'specific humidity flux tendency', 'W m-2 s-1' ) ! 印字 ; Print ! call MessageNotify( 'M', module_name, '----- Initialization Messages -----' ) call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) ) phy_implicit_inited = .true. end subroutine PhyImplInit
Subroutine : | |||
jna_LUMtx(JDim, NDim, -1:1) : | real(DP), intent(inout)
| ||
JDim : | integer, intent(in) | ||
NDim : | integer, intent(in) |
3 重対角行列の LU 分解を行います.
LU decomposition of triple diagonal matrix.
subroutine PhyImplLUDecomp3( jna_LUMtx, JDim, NDim ) ! ! 3 重対角行列の LU 分解を行います. ! ! LU decomposition of triple diagonal matrix. ! ! 宣言文 ; Declaration statements ! implicit none integer, intent(in):: JDim integer, intent(in):: NDim real(DP), intent(inout):: jna_LUMtx(JDim, NDim, -1:1) ! LU 行列. ! LU matrix ! 作業変数 ! Work variables ! integer:: j, n ! DO ループ用作業変数 ! Work variables for DO loop ! 実行文 ; Executable statement ! ! LU 分解 ! LU decomposition ! do j = 1, JDim jna_LUMtx(j,1,1) = jna_LUMtx(j,1,1) / jna_LUMtx(j,1,0) end do do n = 2, NDim-1 do j = 1, JDim jna_LUMtx(j,n,0) = jna_LUMtx(j,n,0) - jna_LUMtx(j,n,-1) * jna_LUMtx(j,n-1,1) jna_LUMtx(j,n,1) = jna_LUMtx(j,n,1) / jna_LUMtx(j,n,0) end do end do do j = 1, JDim jna_LUMtx(j,NDim,0) = jna_LUMtx(j,NDim, 0) - jna_LUMtx(j,NDim,-1) * jna_LUMtx(j,NDim-1,1) end do end subroutine PhyImplLUDecomp3
Subroutine : | |||
ijn_Vector(IDim, JDim, NDim) : | real(DP), intent(inout)
| ||
jna_LUMtx(JDim, NDim, -1:1) : | real(DP), intent(in)
| ||
IDim : | integer, intent(in) | ||
JDim : | integer, intent(in) | ||
NDim : | integer, intent(in) |
LU 分解による解の計算 (3重対角行列用) を行います.
Solve with LU decomposition (For triple diagonal matrix).
subroutine PhyImplLUSolve3( ijn_Vector, jna_LUMtx, IDim, JDim, NDim ) ! ! LU 分解による解の計算 (3重対角行列用) を行います. ! ! Solve with LU decomposition (For triple diagonal matrix). ! ! 宣言文 ; Declaration statements ! implicit none integer, intent(in):: IDim integer, intent(in):: JDim integer, intent(in):: NDim real(DP), intent(in):: jna_LUMtx(JDim, NDim, -1:1) ! LU 行列. ! LU matrix real(DP), intent(inout):: ijn_Vector(IDim, JDim, NDim) ! 右辺ベクトル / 解. ! Right-hand side vector / solution ! 作業変数 ! Work variables ! integer:: i, j, n ! DO ループ用作業変数 ! Work variables for DO loop ! 実行文 ; Executable statement ! ! 前進代入 ! Forward substitution ! do i = 1, IDim do j = 1, JDim ijn_Vector(i,j,1) = ijn_Vector(i,j,1) / jna_LUMtx(j,1,0) end do end do do n = 2, NDim do i = 1, IDim do j = 1, JDim ijn_Vector(i,j,n) = ( ijn_Vector(i,j,n) - ijn_Vector(i,j,n-1) * jna_LUMtx(j,n,-1) ) / jna_LUMtx(j,n,0) end do end do end do ! 後退代入 ! Backward substitution ! do n = NDim-1, 1, -1 do i = 1, IDim do j = 1, JDim ijn_Vector(i,j,n) = ijn_Vector(i,j,n) - ijn_Vector(i,j,n+1) * jna_LUMtx(j,n,1) end do end do end do end subroutine PhyImplLUSolve3
Constant : | |||
module_name = ‘phy_implicit‘ : | character(*), parameter
|
Constant : | |||
version = ’$Name: dcpam5-20081129 $’ // ’$Id: phy_implicit.f90,v 1.7 2008-11-09 05:46:23 morikawa Exp $’ : | character(*), parameter
|