Class | phy_implicit_utils |
In: |
phy_implicit/phy_implicit_utils.f90
|
Note that Japanese and English are described in parallel.
PhyImplEvalRadLFluxA : | 長波フラックス補正 |
———— : | ———— |
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 ! ! 初期化確認 ! Initialization check ! if ( .not. phy_implicit_utils_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if ! 計算時間計測開始 ! Start measurement of computation time ! call TimesetClockStart( module_name ) ! $ 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. * DelTime end do ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine PhyImplEvalRadLFluxA
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
Subroutine : |
phy_implicit_utils モジュールの初期化を行います. NAMELIST#phy_implicit_nml の読み込みはこの手続きで行われます.
"phy_implicit_utils" module is initialized. "NAMELIST#phy_implicit_utils_nml" is loaded in this procedure.
subroutine PhyImplUtilsInit ! ! phy_implicit_utils モジュールの初期化を行います. ! NAMELIST#phy_implicit_nml の読み込みはこの手続きで行われます. ! ! "phy_implicit_utils" module is initialized. ! "NAMELIST#phy_implicit_utils_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/ !!$ ! !!$ ! デフォルト値については初期化手続 "phy_implicit#PhyImplInit" !!$ ! のソースコードを参照のこと. !!$ ! !!$ ! Refer to source codes in the initialization procedure !!$ ! "phy_implicit#PhyImplInit" for the default values. !!$ ! ! 実行文 ; Executable statement ! if ( phy_implicit_utils_inited ) return ! デフォルト値の設定 ! 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 ! 印字 ; Print ! call MessageNotify( 'M', module_name, '----- Initialization Messages -----' ) call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) ) phy_implicit_utils_inited = .true. end subroutine PhyImplUtilsInit
Constant : | |||
module_name = ‘phy_implicit_utils‘ : | character(*), parameter
|
Variable : | |||
phy_implicit_utils_inited = .false. : | logical, save
|
Constant : | |||
version = ’$Name: dcpam5-20140204 $’ // ’$Id: phy_implicit_utils.f90,v 1.2 2011-07-08 03:55:19 yot Exp $’ : | character(*), parameter
|