Class phy_implicit
In: phy_implicit/phy_implicit.F90

陰解法による時間積分

Time integration with implicit scheme

Note that Japanese and English are described in parallel.

Procedures List

PhyImplTendency :時間変化率の計算
PhyImplEvalRadLFluxA :長波フラックス補正
———— :————
PhyImplTendency :Calculate tendency
PhyImplEvalRadLFluxA :Longwave flux correction

Methods

Included Modules

gridset composition dc_types dc_message constants timeset gtool_historyauto saturate_nha1992 saturate_t1930 Bucket_Model namelist_util dc_iounit dc_string axesset

Public Instance methods

Subroutine :
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 長波フラックス. Longwave flux
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ DP{T}{t} $ . 温度変化. Temperature tendency
xy_DSurfTempDt(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度変化率. Surface temperature tendency
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(in)
: 長波地表温度変化. Surface temperature tendency with longwave
xyr_RadLFluxA(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: $ t-\Delta t $ における変化率を元に 算出された $ t+Delta t $ における 長波フラックス.

Longwave flux at $ t+Delta t $ calculated from the tendency at $ t-\Delta t $ .

$ t-\Delta t $ における変化率を元に, $ t+Delta t $ の長波フラックス (xyr_RadLFluxA) を算出します.

Evaluate longwave flux at $ t+Delta t $ (xyr_RadLFluxA) from the tendency at $ t-\Delta t $ .

[Source]

  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. * 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)
: 東西風速フラックス. Eastward wind flux
xyr_VFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 南北風速フラックス. Northward wind flux
xyr_TempFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 温度フラックス. Temperature flux
xyrf_QMixFlux(0:imax-1, 1:jmax, 0:kmax, 1:ncmax) :real(DP), intent(in)
: 質量フラックス. Mass flux of constituents
xyr_RadSFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 短波 (日射) フラックス. Shortwave (insolation) flux
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 長波フラックス. Longwave flux
xy_GroundTempFlux(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地中熱フラックス. Ground temperature flux
xy_SurfTemp(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度. Surface temperature
xy_SurfHumidCoef(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表湿潤度. Surface humidity coefficient
xy_SurfCond(0:imax-1, 1:jmax) :integer, intent(in)
: 地表状態. Surface condition
xy_SurfHeatCapacity(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表熱容量. Surface heat capacity
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(in)
: 長波地表温度変化. Surface temperature tendency with longwave
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
xyz_Exner(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: Exner 関数 (整数レベル). Exner function (full level)
xyr_Exner(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: Exner 関数 (半整数レベル). Exner function (half level)
xyr_VelTransCoef(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 輸送係数:運動量. Transfer coefficient: velocity
xyr_TempTransCoef(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 輸送係数:温度. Transfer coefficient: temperature
xyr_QMixTransCoef(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 輸送係数:質量. Transfer coefficient: mass of constituents
xy_SurfVelTransCoef(0:imax-1, 1:jmax) :real(DP), intent(in)
: 輸送係数:運動量. Diffusion coefficient: velocity
xy_SurfTempTransCoef(0:imax-1, 1:jmax) :real(DP), intent(in)
: 輸送係数:温度. Transfer coefficient: temperature
xy_SurfQVapTransCoef(0:imax-1, 1:jmax) :real(DP), intent(in)
: 輸送係数:比湿. Transfer coefficient: specific humidity
xyz_DUDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ DP{u}{t} $ . 東西風速変化. Eastward wind tendency
xyz_DVDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ DP{v}{t} $ . 南北風速変化. Northward wind tendency
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ DP{T}{t} $ . 温度変化. Temperature tendency
xyzf_DQMixDt(0:imax-1, 1:jmax, 1:kmax, 1:ncmax) :real(DP), intent(out)
: $ DP{q}{t} $ . 比湿変化. Temperature tendency
xy_DSurfTempDt(0:imax-1, 1:jmax) :real(DP), intent(out)
: 地表面温度変化率. Surface temperature tendency

時間変化率の計算を行います.

Calculate tendencies.

[Source]

  subroutine PhyImplTendency( xyr_UFlux, xyr_VFlux, xyr_TempFlux, xyrf_QMixFlux, xyr_RadSFlux, xyr_RadLFlux, xy_GroundTempFlux, xy_SurfTemp, xy_SurfHumidCoef, xy_SurfCond, xy_SurfHeatCapacity, xyra_DelRadLFlux, xyr_Press, xyz_Exner, xyr_Exner, xyr_VelTransCoef, xyr_TempTransCoef, xyr_QMixTransCoef, xy_SurfVelTransCoef, xy_SurfTempTransCoef, xy_SurfQVapTransCoef, xyz_DUDt, xyz_DVDt, xyz_DTempDt, xyzf_DQMixDt, 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

    ! バケツモデル
    ! bucket model
    !
    use Bucket_Model, only: FlagBucketModel, FlagBucketModelSnow


    ! 宣言文 ; 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):: xyrf_QMixFlux (0:imax-1, 1:jmax, 0:kmax, 1:ncmax)
                              ! 質量フラックス. 
                              ! Mass flux of constituents

    real(DP), intent(in):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 短波 (日射) フラックス. 
                              ! Shortwave (insolation) flux
    real(DP), intent(in):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Longwave flux

    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_SurfHumidCoef (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):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(in):: xyz_Exner (0:imax-1, 1:jmax, 1:kmax)
                              ! Exner 関数 (整数レベル). 
                              ! Exner function (full level)
    real(DP), intent(in):: xyr_Exner (0:imax-1, 1:jmax, 0:kmax)
                              ! Exner 関数 (半整数レベル). 
                              ! Exner function (half level)

    real(DP), intent(in):: xyr_VelTransCoef (0:imax-1, 1:jmax, 0:kmax)
                              ! 輸送係数:運動量. 
                              ! Transfer coefficient: velocity
    real(DP), intent(in):: xyr_TempTransCoef (0:imax-1, 1:jmax, 0:kmax)
                              ! 輸送係数:温度. 
                              ! Transfer coefficient: temperature
    real(DP), intent(in):: xyr_QMixTransCoef (0:imax-1, 1:jmax, 0:kmax)
                              ! 輸送係数:質量. 
                              ! Transfer coefficient: mass of constituents

    real(DP), intent(in):: xy_SurfVelTransCoef (0:imax-1, 1:jmax)
                              ! 輸送係数:運動量. 
                              ! Diffusion coefficient: velocity
    real(DP), intent(in):: xy_SurfTempTransCoef (0:imax-1, 1:jmax)
                              ! 輸送係数:温度. 
                              ! Transfer coefficient: temperature
    real(DP), intent(in):: xy_SurfQVapTransCoef (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):: xyzf_DQMixDt(0:imax-1, 1:jmax, 1:kmax, 1:ncmax)
                              ! $ \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):: xyza_QMixMtx (0:imax-1, 1:jmax, 1:kmax, -1:1)
                              ! 質量混合比陰解行列. 
                              ! Implicit matrix about mass mixing ratio
    real(DP):: xyaa_SurfMtx (0:imax-1, 1:jmax, 0:0, -1:1)
                              ! 惑星表面エネルギー収支用陰解行列
                              ! Implicit matrix for surface energy balance
    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):: xyza_QMixLUMtx(0:imax-1, 1:jmax, 1:kmax,-1:1)
                              ! LU 行列. 
                              ! LU matrix
    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
    integer:: n               ! 組成方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in dimension of constituents

    ! 実行文 ; Executable statement
    !

    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )

    ! 初期化
    ! Initialization
    !
    if ( .not. phy_implicit_inited ) call PhyImplInit


    ! バケツモデルの扱いの確認
    ! Check use of bucket model
    !
    ! Bucket model with this routine is not supported fully, although some parts of 
    ! calculation is implemented. 
    !
    if ( FlagBucketModel ) then
      call MessageNotify( 'E', module_name, 'Bucket model cannot be used with this routine.' )
    end if

    ! 雪の扱いの確認
    ! Check about treatment of snow
    !
    if ( FlagBucketModelSnow ) then
      call MessageNotify( 'E', module_name, 'Snow melt is not treated in this routine.' )
    end if


    ! 陰解法のための行列作成
    ! Create matrices for implicit scheme
    !

    ! 鉛直拡散スキームの輸送係数から陰解行列の計算 (速度)
    ! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (velocity)
    !
    k = 1
    xyza_UVMtx  (:,:,k,-1) = 0.0d0
    xyza_UVMtx  (:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xy_SurfVelTransCoef(:,:) + xyr_VelTransCoef(:,:,k  )
    xyza_UVMtx  (:,:,k, 1) = - xyr_VelTransCoef(:,:,k)

    do k = 2, kmax-1
      xyza_UVMtx  (:,:,k,-1) = - xyr_VelTransCoef(:,:,k-1)
      xyza_UVMtx  (:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xyr_VelTransCoef(:,:,k-1) + xyr_VelTransCoef(:,:,k  )
      xyza_UVMtx  (:,:,k, 1) = - xyr_VelTransCoef(:,:,k)
    end do

    k = kmax
    xyza_UVMtx  (:,:,k,-1) = - xyr_VelTransCoef(:,:,k-1)
    xyza_UVMtx  (:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xyr_VelTransCoef(:,:,k-1)
    xyza_UVMtx  (:,:,k, 1) = 0.0d0


    ! 鉛直拡散スキームの輸送係数から陰解行列の計算 (温度)
    ! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (temperature)
    !
    k = 1
    xyra_TempMtx(:,:,k,-1) = - CpDry * xy_SurfTempTransCoef(:,:)
    xyra_TempMtx(:,:,k, 0) = - CpDry * ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k  ) * xy_SurfTempTransCoef(:,:) + CpDry * xyr_Exner(:,:,k  ) / xyz_Exner(:,:,k  ) * xyr_TempTransCoef(:,:,k  )
    xyra_TempMtx(:,:,k, 1) = - CpDry * xyr_Exner(:,:,k  ) / xyz_Exner(:,:,k+1) * xyr_TempTransCoef(:,:,k  )

    do k = 2, kmax-1
      xyra_TempMtx(:,:,k,-1) = - CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k-1) * xyr_TempTransCoef(:,:,k-1)
      xyra_TempMtx(:,:,k, 0) = - CpDry * ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k  ) * xyr_TempTransCoef(:,:,k-1) + CpDry * xyr_Exner(:,:,k  ) / xyz_Exner(:,:,k  ) * xyr_TempTransCoef(:,:,k  )
      xyra_TempMtx(:,:,k, 1) = - CpDry * xyr_Exner(:,:,k  ) / xyz_Exner(:,:,k+1) * xyr_TempTransCoef(:,:,k  )
    end do

    k = kmax
    xyra_TempMtx(:,:,k,-1) = - CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k-1) * xyr_TempTransCoef(:,:,k-1)
    xyra_TempMtx(:,:,k, 0) = - CpDry * ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k  ) * xyr_TempTransCoef(:,:,k-1)
    xyra_TempMtx(:,:,k, 1) = 0.0d0




    ! 鉛直拡散スキームの輸送係数から陰解行列の計算 (比湿)
    ! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (specific humidity)
    !
    ! 飽和比湿の計算
    ! 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

    k = 1
    if ( FlagBucketModel ) then
      xyza_QVapMtx(:,:,k,-1) = 0.0d0
      xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xyr_QMixTransCoef(:,:,k  )
    else
      xyza_QVapMtx(:,:,k,-1) = - xy_SurfHumidCoef(:,:) * xy_SurfQVapTransCoef(:,:) * xy_SurfDQVapSatDTemp(:,:)
      xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xy_SurfHumidCoef(:,:) * xy_SurfQVapTransCoef(:,:) + xyr_QMixTransCoef(:,:,k  )
    end if
    xyza_QVapMtx(:,:,k, 1) = - xyr_QMixTransCoef(:,:,k  )

    do k = 2, kmax-1
      xyza_QVapMtx(:,:,k,-1) = - xyr_QMixTransCoef(:,:,k-1)
      xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xyr_QMixTransCoef(:,:,k-1) + xyr_QMixTransCoef(:,:,k  )
      xyza_QVapMtx(:,:,k, 1) = - xyr_QMixTransCoef(:,:,k  )
    end do

    k = kmax
    xyza_QVapMtx(:,:,k,-1) = - xyr_QMixTransCoef(:,:,k-1)
    xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xyr_QMixTransCoef(:,:,k-1)
    xyza_QVapMtx(:,:,k, 1) = 0.0d0


    ! 鉛直拡散スキームの輸送係数から陰解行列の計算 (比湿を除く混合比)
    ! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (mixing ratio except for specific humidity)
    !
    k = 1
    xyza_QMixMtx(:,:,k,-1) = 0.0d0
    xyza_QMixMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2. * DelTime ) + xyr_QMixTransCoef(:,:,k  )
    xyza_QMixMtx(:,:,k, 1) = - xyr_QMixTransCoef(:,:,k  )

    do l = -1, 1
      do k = 2, kmax
        xyza_QMixMtx(:,:,k,l) = xyza_QVapMtx(:,:,k,l)
      end do
    end do


    ! 地表面過程の輸送係数から陰解行列の計算
    ! Calculate implicit matrices from transfer coefficient of surface process
    !
    do i = 0, imax-1
      do j = 1, jmax
        if ( xy_SurfCond(i,j) >= 1 ) then

          if ( FlagBucketModel ) then
            xyaa_SurfMtx(i,j,0,-1) = 0.0d0
            xyaa_SurfMtx(i,j,0, 0) = xy_SurfHeatCapacity(i,j) / ( 2. * DelTime ) + CpDry * xy_SurfTempTransCoef(i,j) + xyra_DelRadLFlux(i,j,0,0)
          else
            xyaa_SurfMtx(i,j,0,-1) = - LatentHeat * xy_SurfHumidCoef(i,j) * xy_SurfQVapTransCoef(i,j)
            xyaa_SurfMtx(i,j,0, 0) = xy_SurfHeatCapacity(i,j) / ( 2. * DelTime ) + CpDry * xy_SurfTempTransCoef(i,j) + xyra_DelRadLFlux(i,j,0,0) + LatentHeat * xy_SurfHumidCoef(i,j) * xy_SurfQVapTransCoef(i,j) * xy_SurfDQVapSatDTemp(i,j)
          end if
          xyaa_SurfMtx(i,j,0, 1) = - CpDry * xyr_Exner(i,j,0) / xyz_Exner(i,j,1) * xy_SurfTempTransCoef(i,j) + xyra_DelRadLFlux(i,j,0,1)
        else
          xyaa_SurfMtx(i,j,0,-1) = 0.0d0
          xyaa_SurfMtx(i,j,0, 0) = 1.0d0
          xyaa_SurfMtx(i,j,0, 1) = 0.0d0
        end if
      end do
    end do


    ! 東西風速, 南北風速の計算
    ! Calculate eastward and northward wind
    !
    xyza_UVLUMtx = xyza_UVMtx

    call PhyImplLUDecomp3( xyza_UVLUMtx, imax * jmax, kmax ) ! (in)

    do k = 1, kmax
      xyz_DUDt(:,:,k) = - ( xyr_UFlux(:,:,k) - xyr_UFlux(:,:,k-1) )
      xyz_DVDt(:,:,k) = - ( xyr_VFlux(:,:,k) - xyr_VFlux(:,:,k-1) )
    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) = xyza_QVapMtx(:,:,k,l)
      end do

      k = 0
      xyza_TempQVapLUMtx(:,:, k, l) = xyaa_SurfMtx(:,:,0,l)

      do k =  1,  kmax
        xyza_TempQVapLUMtx(:,:, k, l) = xyra_TempMtx(:,:,k,l)
      end do

    end do

    call PhyImplLUDecomp3( xyza_TempQVapLUMtx, imax * jmax, (2 * kmax) + 1 ) ! (in)

    do k = 1, kmax
      xyz_DelTempQVap(:,:,-k) = - ( xyrf_QMixFlux(:,:,k,IndexH2OVap) - xyrf_QMixFlux(:,:,k-1,IndexH2OVap) )
    end do

    k = 0
    do j = 1, jmax
      do i = 0, imax-1
        if ( xy_SurfCond(i,j) >= 1 ) then
          xyz_DelTempQVap(i,j,k) = - xyr_RadSFlux(i,j,0) - xyr_RadLFlux(i,j,0) - xyr_TempFlux(i,j,0) - LatentHeat * xyrf_QMixFlux(i,j,0,IndexH2OVap) + xy_GroundTempFlux(i,j)
        else
          xyz_DelTempQVap(i,j,k) = 0.0d0
        end if
      end do
    end do

    do k = 1, kmax
      xyz_DelTempQVap(:,:,k) = - ( xyr_TempFlux(:,:, k) - xyr_TempFlux(:,:,  k-1 ) )
    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. * DelTime )
      xyz_DVDt(:,:,k)                 = xyz_DVDt(:,:,k)         / ( 2. * DelTime )
      xyz_DTempDt(:,:,k)              = xyz_DelTempQVap(:,:, k) / ( 2. * DelTime )
      xyzf_DQMixDt(:,:,k,IndexH2OVap) = xyz_DelTempQVap(:,:,-k) / ( 2. * DelTime )
    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. * DelTime )
        else
          xy_DSurfTempDt(i,j) = 0.
        end if
      end do
    end do


    ! 比湿を除く質量混合比
    ! Calculate mass mixing ratio except for specific humidity
    !
    xyza_QMixLUMtx = xyza_QMixMtx

    call PhyImplLUDecomp3( xyza_QMixLUMtx, imax * jmax, kmax )

    do n = 1, ncmax
      if ( n == IndexH2OVap ) cycle

      do k = 1, kmax
        xyzf_DQMixDt(:,:,k,n) = - ( xyrf_QMixFlux(:,:,k,n) - xyrf_QMixFlux(:,:,k-1,n) )
      end do

      call PhyImplLUSolve3( xyzf_DQMixDt(:,:,:,n), xyza_QMixLUMtx, 1, imax * jmax, kmax )
    end do



!!$    !#########################################################
!!$    ! code for debug, this will be removed, (Y. O. Takahashi, 2009/04/07)
!!$    i = 1
!!$    j = jmax / 2
!!$    write( 6, * ) &
!!$      & - xyr_RadSFlux(i,j,0),                                                      &
!!$      & - ( xyr_RadLFlux(i,j,0)                                                     &
!!$      &   + xyra_DelRadLFlux(i,j,0,0) * xy_DSurfTempDt(i,j) * ( 2.0d0 * DelTime )   &
!!$      &   + xyra_DelRadLFlux(i,j,0,1) * xyz_DTempDt(i,j,1) * ( 2.0d0 * DelTime ) ), &
!!$      & - ( xyr_TempFlux(i,j,0)                                                     &
!!$      &   - CpDry * xyr_Exner(i,j,0) * xy_SurfTempTransCoef(i,j)                    &
!!$      &     * ( xyz_DTempDt(i,j,1) / xyz_Exner(i,j,1)                               &
!!$      &       - xy_DSurfTempDt(i,j) / xyr_Exner(i,j,0) ) * ( 2.0d0 * DelTime ) ),   &
!!$      & - LatentHeat                                                                &
!!$      &   * ( xyr_QVapFlux(i,j,0)                                                   &
!!$      &     - xy_SurfQVapTransCoef(i,j)                                             &
!!$      &       * ( xyz_DQVapDt(i,j,1)                                                &
!!$      &         - xy_SurfDQVapSatDTemp(i,j) * xy_DSurfTempDt(i,j) )                 &
!!$      &       * ( 2.0d0 * DelTime ) ) !, &
!!$!      & + xy_GroundTempFlux(i,j)
!!$
!!$    xy_SurfQVapSat(i,j) =                                                           &
!!$      & - xyr_RadSFlux(i,j,0)                                                       &
!!$      & - ( xyr_RadLFlux(i,j,0)                                                     &
!!$      &   + xyra_DelRadLFlux(i,j,0,0) * xy_DSurfTempDt(i,j) * ( 2.0d0 * DelTime )   &
!!$      &   + xyra_DelRadLFlux(i,j,0,1) * xyz_DTempDt(i,j,1) * ( 2.0d0 * DelTime ) )  &
!!$      & - ( xyr_TempFlux(i,j,0)                                                     &
!!$      &   - CpDry * xyr_Exner(i,j,0) * xy_SurfTempTransCoef(i,j)                    &
!!$      &     * ( xyz_DTempDt(i,j,1) / xyz_Exner(i,j,1)                               &
!!$      &       - xy_DSurfTempDt(i,j) / xyr_Exner(i,j,0) ) * ( 2.0d0 * DelTime ) )    &
!!$      & - LatentHeat                                                                &
!!$!      &   * ( xyr_QVapFlux(i,j,0)                                                   &
!!$!      &     - xy_SurfQVapTransCoef(i,j)                                             &
!!$!      &       * ( xyz_DQVapDt(i,j,1)                                                &
!!$!      &         - xy_SurfDQVapSatDTemp(i,j) * xy_DSurfTempDt(i,j) )                 &
!!$!      &       * ( 2.0d0 * DelTime ) )
!!$      &   * xyr_QVapFlux(i,j,0)
!!$    write( 6, * ) xy_SurfQVapSat(i,j)
!!$    !#########################################################



    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )

  end subroutine PhyImplTendency
phy_implicit_inited
Variable :
phy_implicit_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  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.

[Source]

  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/
!!$          !
!!$          ! デフォルト値については初期化手続 "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

    ! 印字 ; 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)
: LU 行列. LU matrix
JDim :integer, intent(in)
NDim :integer, intent(in)

3 重対角行列の LU 分解を行います.

LU decomposition of triple diagonal matrix.

[Source]

  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)
: 右辺ベクトル / 解. Right-hand side vector / solution
jna_LUMtx(JDim, NDim, -1:1) :real(DP), intent(in)
: LU 行列. LU matrix
IDim :integer, intent(in)
JDim :integer, intent(in)
NDim :integer, intent(in)

LU 分解による解の計算 (3重対角行列用) を行います.

Solve with LU decomposition (For triple diagonal matrix).

[Source]

  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
module_name
Constant :
module_name = ‘phy_implicit :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20100424 $’ // ’$Id: phy_implicit.F90,v 1.18 2009-10-05 14:25:32 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version