Class | phy_implicit |
In: |
phy_implicit/phy_implicit.F90
|
Note that Japanese and English are described in parallel.
PhyImplTendency : | 時間変化率の計算 |
PhyImplEvalRadLFluxA : | 長波フラックス補正 |
———— : | ———— |
PhyImplTendency : | Calculate tendency |
PhyImplEvalRadLFluxA : | Longwave flux correction |
Subroutine : | |||
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
xy_DSurfTempDt(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(in)
| ||
xyr_RadLFluxA(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
|
$ t-\Delta t $ における変化率を元に, $ t+Delta t $ の長波フラックス (xyr_RadLFluxA) を算出します.
Evaluate longwave flux at $ t+Delta t $ (xyr_RadLFluxA) from the tendency at $ t-\Delta t $ .
subroutine PhyImplEvalRadLFluxA( xyr_RadLFlux, xyz_DTempDt, xy_DSurfTempDt, xyra_DelRadLFlux, xyr_RadLFluxA ) ! ! $ t-\Delta t $ における変化率を元に, ! $ t+\Delta t $ の長波フラックス (xyr_RadLFluxA) を算出します. ! ! Evaluate longwave flux at $ t+\Delta t $ (xyr_RadLFluxA) ! from the tendency at $ t-\Delta t $ . ! ! モジュール引用 ; USE statements ! ! 時刻管理 ! Time control ! use timeset, only: DelTime, TimesetClockStart, TimesetClockStop ! 宣言文 ; Declaration statements ! implicit none real(DP), intent(in):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Longwave flux real(DP), intent(in):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax) ! $ \DP{T}{t} $ . 温度変化. ! Temperature tendency real(DP), intent(in):: xy_DSurfTempDt (0:imax-1, 1:jmax) ! 地表面温度変化率. ! Surface temperature tendency real(DP), intent(in):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Surface temperature tendency with longwave real(DP), intent(out):: xyr_RadLFluxA (0:imax-1, 1:jmax, 0:kmax) ! $ t-\Delta t $ における変化率を元に ! 算出された $ t+\Delta t $ における ! 長波フラックス. ! ! Longwave flux at $ t+\Delta t $ ! calculated from the tendency at ! $ t-\Delta t $ . ! 作業変数 ! Work variables ! 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 ! $ t+\Delta t $ の長波フラックス (xyr_RadLFluxA) を算出 ! Evaluate longwave flux at $ t+\Delta t $ (xyr_RadLFluxA) ! do k = 0, kmax xyr_RadLFluxA(:,:,k) = xyr_RadLFlux(:,:,k) + ( xy_DSurfTempDt * xyra_DelRadLFlux(:,:,k,0) + xyz_DTempDt(:,:,1) * xyra_DelRadLFlux(:,:,k,1) ) * 2.0_DP * DelTime end do ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine PhyImplEvalRadLFluxA
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)
| ||
xy_SurfTemp(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfHumidCoeff(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfCond(0:imax-1, 1:jmax) : | integer, intent(in)
| ||
xy_SurfHeatCapacity(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | 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)
| ||
xyr_VelTransCoeff(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_TempTransCoeff(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_QVapTransCoeff(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xy_SurfVelTransCoeff(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfTempTransCoeff(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfQVapTransCoeff(0:imax-1, 1:jmax) : | real(DP), 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, xy_SurfTemp, xy_SurfHumidCoeff, xy_SurfCond, xy_SurfHeatCapacity, xyra_DelRadLFlux, xyz_Press, xyr_Press, xyr_VelTransCoeff, xyr_TempTransCoeff, xyr_QVapTransCoeff, xy_SurfVelTransCoeff, xy_SurfTempTransCoeff, xy_SurfQVapTransCoeff, xyz_DUDt, xyz_DVDt, xyz_DTempDt, xyz_DQVapDt, xy_DSurfTempDt ) ! ! 時間変化率の計算を行います. ! ! Calculate tendencies. ! ! モジュール引用 ; USE statements ! ! 物理定数設定 ! Physical constants settings ! use constants, only: Grav, CpDry, LatentHeat, GasRDry ! $ R $ [J kg-1 K-1]. ! 乾燥大気の気体定数. ! Gas constant of air ! 時刻管理 ! 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 ! 飽和比湿計算 ! Evaluate saturation specific humidity ! #ifdef LIB_SATURATE_NHA1992 use saturate_nha1992, only: CalcQVapSat, CalcDQVapSatDTemp #elif LIB_SATURATE_T1930 use saturate_t1930, only: CalcQVapSat, CalcDQVapSatDTemp #else use saturate_t1930, only: CalcQVapSat, CalcDQVapSatDTemp #endif ! 宣言文 ; 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):: xy_SurfTemp (0:imax-1, 1:jmax) ! 地表面温度. ! Surface temperature real(DP), intent(in):: xy_SurfHumidCoeff (0:imax-1, 1:jmax) ! 地表湿潤度. ! Surface humidity coefficient integer, intent(in):: xy_SurfCond (0:imax-1, 1:jmax) ! 地表状態. ! Surface condition real(DP), intent(in):: xy_SurfHeatCapacity (0:imax-1, 1:jmax) ! 地表熱容量. ! Surface heat capacity real(DP), intent(in):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Surface temperature tendency with longwave 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(in):: xyr_VelTransCoeff (0:imax-1, 1:jmax, 0:kmax) ! 輸送係数:運動量. ! Transfer coefficient: velocity real(DP), intent(in):: xyr_TempTransCoeff (0:imax-1, 1:jmax, 0:kmax) ! 輸送係数:温度. ! Transfer coefficient: temperature real(DP), intent(in):: xyr_QVapTransCoeff (0:imax-1, 1:jmax, 0:kmax) ! 輸送係数:比湿. ! Transfer coefficient: specific humidity real(DP), intent(in):: xy_SurfVelTransCoeff (0:imax-1, 1:jmax) ! 輸送係数:運動量. ! Diffusion coefficient: velocity real(DP), intent(in):: xy_SurfTempTransCoeff (0:imax-1, 1:jmax) ! 輸送係数:温度. ! Transfer coefficient: temperature real(DP), intent(in):: xy_SurfQVapTransCoeff (0:imax-1, 1:jmax) ! 輸送係数:比湿. ! Transfer coefficient: specific humidity 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):: xyza_UVMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 速度陰解行列. ! Implicit matrix about velocity real(DP):: xyra_TempMtx (0:imax-1, 1:jmax, 0:kmax, -1:1) ! 温度陰解行列. ! Implicit matrix about temperature real(DP):: xyza_QVapMtx (0:imax-1, 1:jmax, 1:kmax, -1:1) ! 比湿陰解行列. ! Implicit matrix about specific humidity real(DP):: xy_SurfUVMtx (0:imax-1, 1:jmax) ! 速度陰解行列: 地表. ! Implicit matrix about velocity: surface real(DP):: xyaa_SurfTempMtx (0:imax-1, 1:jmax, 0:1, -1:1) ! 温度陰解行列: 地表. ! Implicit matrix about temperature: surface real(DP):: xyaa_SurfQVapMtx (0:imax-1, 1:jmax, 0:1, -1:1) ! 比湿陰解行列: 地表. ! Implicit matrix about specific humidity: surface real(DP):: xya_SurfRadLMtx (0:imax-1, 1:jmax, -1:1) ! $ T $ . 陰解行列: 地表. ! implicit matrix: surface 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 real(DP):: xyz_Exner (0:imax-1, 1:jmax, 1:kmax) ! Exner 関数 (整数レベル). ! Exner function (full level) real(DP):: xyr_Exner (0:imax-1, 1:jmax, 0:kmax) ! Exner 関数 (半整数レベル). ! Exner function (half level) real(DP):: xy_SurfExner (0:imax-1, 1:jmax) ! Exner 関数. ! Exner function real(DP):: xy_SurfQVapSat (0:imax-1, 1:jmax) ! 地表飽和比湿. ! Saturated specific humidity on surface real(DP):: xy_SurfDQVapSatDTemp (0:imax-1, 1:jmax) ! 地表飽和比湿変化. ! Saturated specific humidity tendency on surface 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 long wave implicit matrix ! xya_SurfRadLMtx(:,:,0) = xyra_DelRadLFlux(:,:,0,0) xya_SurfRadLMtx(:,:,1) = xyra_DelRadLFlux(:,:,0,1) xya_SurfRadLMtx(:,:,-1) = 0.0_DP ! Exner 関数算出 ! Calculate Exner functions ! xyz_Exner = ( xyz_Press / RefPress ) ** ( GasRDry / CpDry ) xyr_Exner = ( xyr_Press / RefPress ) ** ( GasRDry / CpDry ) ! 陰解法のための行列作成 ! Create matrices for implicit scheme ! call ImplMatrices( xyr_Press, xy_SurfHeatCapacity, xyza_UVMtx, xyra_TempMtx, xyza_QVapMtx ) ! (out) ! 鉛直拡散スキームの輸送係数から陰解行列の計算 ! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme ! do k = 2, kmax xyza_UVMtx(:,:,k, 0) = xyza_UVMtx(:,:,k,0) + xyr_VelTransCoeff(:,:,k-1) xyza_UVMtx(:,:,k,-1) = - xyr_VelTransCoeff(:,:,k-1) xyra_TempMtx(:,:,k, 0) = xyra_TempMtx(:,:,k,0) + CpDry * xyr_TempTransCoeff(:,:,k-1) * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k) xyra_TempMtx(:,:,k,-1) = - CpDry * xyr_TempTransCoeff(:,:,k-1) * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k-1) xyza_QVapMtx(:,:,k, 0) = xyza_QVapMtx(:,:,k,0) + CpDry * xyr_QVapTransCoeff(:,:,k-1) xyza_QVapMtx(:,:,k,-1) = - CpDry * xyr_QVapTransCoeff(:,:,k-1) end do do k = 1, kmax-1 xyza_UVMtx(:,:,k,0) = xyza_UVMtx(:,:,k,0) + xyr_VelTransCoeff(:,:,k) xyza_UVMtx(:,:,k,1) = - xyr_VelTransCoeff(:,:,k) xyra_TempMtx(:,:,k,0) = xyra_TempMtx(:,:,k,0) + CpDry * xyr_TempTransCoeff(:,:,k) * xyr_Exner(:,:,k) / xyz_Exner(:,:,k) xyra_TempMtx(:,:,k,1) = - CpDry * xyr_TempTransCoeff(:,:,k) * xyr_Exner(:,:,k) / xyz_Exner(:,:,k+1) xyza_QVapMtx(:,:,k,0) = xyza_QVapMtx(:,:,k,0) + CpDry * xyr_QVapTransCoeff(:,:,k) xyza_QVapMtx(:,:,k,1) = - CpDry * xyr_QVapTransCoeff(:,:,k) end do ! Exner 関数算出 ! Calculate Exner functions ! xy_SurfExner = ( xyz_Press(:,:,1) / xyr_Press(:,:,0) )**( GasRDry / CpDry ) ! 飽和比湿の計算 ! Calculate saturated specific humidity ! do i = 0, imax-1 do j = 1, jmax xy_SurfQVapSat(i,j) = CalcQVapSat( xy_SurfTemp(i,j), xyr_Press(i,j,0) ) end do end do do i = 0, imax-1 do j = 1, jmax xy_SurfDQVapSatDTemp(i,j) = CalcDQVapSatDTemp( xy_SurfTemp(i,j), xy_SurfQVapSat(i,j) ) end do end do ! 地表面過程の輸送係数から陰解行列の計算 ! Calculate implicit matrices from transfer coefficient of surface process ! xyaa_SurfTempMtx = 0.0_DP xyaa_SurfQVapMtx = 0.0_DP xy_SurfUVMtx = xy_SurfVelTransCoeff xyaa_SurfTempMtx(:,:,1,0) = CpDry * xy_SurfTempTransCoeff / xy_SurfExner xyaa_SurfTempMtx(:,:,0,1) = - CpDry * xy_SurfTempTransCoeff / xy_SurfExner xyaa_SurfQVapMtx(:,:,1,0) = CpDry * xy_SurfQVapTransCoeff * xy_SurfHumidCoeff xyaa_SurfQVapMtx(:,:,0,1) = - CpDry * xy_SurfQVapTransCoeff * xy_SurfHumidCoeff do i = 0, imax-1 do j = 1, jmax if ( xy_SurfCond(i,j) >= 1 ) then xyaa_SurfTempMtx(i,j,1,-1) = - CpDry * xy_SurfTempTransCoeff(i,j) xyaa_SurfTempMtx(i,j,0,0) = CpDry * xy_SurfTempTransCoeff(i,j) xyaa_SurfQVapMtx(i,j,1,-1) = - LatentHeat * xy_SurfQVapTransCoeff(i,j) * xy_SurfHumidCoeff(i,j) * xy_SurfDQVapSatDTemp(i,j) xyaa_SurfQVapMtx(i,j,0,0) = LatentHeat * xy_SurfQVapTransCoeff(i,j) * xy_SurfHumidCoeff(i,j) * xy_SurfDQVapSatDTemp(i,j) end if end do end do ! 東西風速, 南北風速の計算 ! 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 do j = 1, jmax do i = 0, imax-1 if ( xy_SurfCond(i,j) >= 1 ) then xyza_TempQVapLUMtx(i,j,0, 0) = xyra_TempMtx(i,j,0,0) + xyaa_SurfTempMtx(i,j,0,0) + xyaa_SurfQVapMtx(i,j,0,0) + xya_SurfRadLMtx(i,j,0) xyza_TempQVapLUMtx(i,j,0, 1) = xyaa_SurfTempMtx(i,j,0,1) + xya_SurfRadLMtx(i,j,1) xyza_TempQVapLUMtx(i,j,0,-1) = xyaa_SurfQVapMtx(i,j,0,1) else xyza_TempQVapLUMtx(i,j,0, 0) = 1.0_DP xyza_TempQVapLUMtx(i,j,0, 1) = 0.0_DP xyza_TempQVapLUMtx(i,j,0,-1) = 0.0_DP end if end do end do 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 do j = 1, jmax do i = 0, imax-1 if ( xy_SurfCond(i,j) >= 1 ) then xyz_DelTempQVap(i,j,0) = - xy_SurfRadSFlux(i,j) - xy_SurfRadLFlux(i,j) - xyr_TempFlux(i,j,0) - xyr_QVapFlux(i,j,0) + xy_GroundTempFlux(i,j) else xyz_DelTempQVap(i,j,0) = 0.0_DP end if end do end do 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 ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine PhyImplTendency
Variable : | |||
phy_implicit_inited = .false. : | logical, save, public
|
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)
| ||
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 ImplMatrices( xyr_Press, xy_SurfHeatCapacity, 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 ! $ \Delta t $ [s] ! 宣言文 ; 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 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:: k ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction ! 実行文 ; Executable statement ! ! 質量, 熱容量の項の計算 ! 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 xyra_TempMtx(:,:,0,0) = xy_SurfHeatCapacity(:,:) / ( 2.0_DP * DelTime ) end subroutine ImplMatrices
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.
This procedure input/output NAMELIST#phy_implicit_nml .
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 ! 作業変数 ! Work variables ! 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/ RefPress ! ! デフォルト値については初期化手続 "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 ! RefPress = 1.0e+6_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 = phy_implicit_nml, iostat = iostat_nml ) ! (out) close( unit_nml ) call NmlutilMsg( iostat_nml, module_name ) ! (in) end if ! 印字 ; Print ! call MessageNotify( 'M', module_name, '----- Initialization Messages -----' ) call MessageNotify( 'M', module_name, ' RefPress = %f', d = (/ RefPress /) ) 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-20090225-2 $’ // ’$Id: phy_implicit.F90,v 1.3 2009-02-09 05:46:35 morikawa Exp $’ : | character(*), parameter
|