Class | DynImpFunc |
In: |
dynamic/dynimpfunc.f90
|
Copyright (C) GFD Dennou Club, 2004, 2006. All rights reserved.
* Developer: SUGIYAMA Ko-ichiro * Version: $Id: dynimpfunc.f90,v 1.1.1.1 2006/04/25 03:43:58 deepconv Exp $ * Tag Name: $Name: $ * Change History:
陰解法を用いたエクスナー関数の計算. deepconv/arare では時間積分として He-VI 法を利用しているので, エクスナー関数は陰解法で解く.
* 離散化する際, 上下境界条件として鉛直速度が零を与えている. * 空間方向に 2 次精度の離散化を陽に利用しているため, differentiate_center4 モジュールを指定することはできないので注意.
subroutine DynImpFunc_init () ! !エクスナー関数を陰解法で解く際に必要となる, 係数行列の要素を決め, !LU 分解を行う. ! !暗黙の型宣言禁止 implicit none !配列の割り付け allocate( & & A(RegZMin+1:RegZMax), & & B(RegZMin+2:RegZMax), & & C(RegZMin+1:RegZMax-1), & & B2(RegZMin+3:RegZMax), & & xz_F1BasicZ(DimXMin:DimXMax, DimZMin:DimZMax), & & xr_F2BasicZ(DimXMin:DimXMax, DimZMin:DimZMax), & & xz_VPotTempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax), & & IPIV(RegZMin+1:RegZMax) ) !仮温位の定義 xz_VPotTempBasicZ = xz_PotTempBasicZ / xz_EffMolWtBasicZ !係数行列の計算 xz_F1BasicZ = ( xz_VelSoundBasicZ ** 2.0d0 ) & & / (CpDry * xz_DensBasicZ *(xz_VPotTempBasicZ ** 2.0d0)) xr_F2BasicZ = xr_avr_xz( & & CpDry * xz_DensBasicZ * (xz_VPotTempBasicZ ** 2.0d0) ) A(RegZMin+2: RegZMax-1) = & & 1.0d0 + (beta ** 2.0d0) & & * xz_F1BasicZ(RegXMax, RegZMin+2: RegZMax-1) & & * ( xr_F2BasicZ(RegXMax, RegZMin+2: RegZMax-1) & & + xr_F2BasicZ(RegXMax, RegZMin+1: RegZMax-2) ) & & * (DelTimeShort ** 2.0d0) / ( DelZ ** 2.0d0) A(RegZMin+1) = 1.0d0 & & + (beta ** 2.0d0) & & * xz_F1BasicZ(RegXMax, RegZMin+1) & & * xr_F2BasicZ(RegXMax, RegZMin+1) & & * (DelTimeShort ** 2.0d0) / ( DelZ ** 2.0d0) A(RegZMax) = 1.0d0 & & + (beta ** 2.0d0) & & * xz_F1BasicZ(RegXMax, RegZMax) & & * xr_F2BasicZ(RegXMax, RegZMax-1) & & * (DelTimeShort ** 2.0d0) / ( DelZ ** 2.0d0) B(RegZMin+2: RegZMax) = & & - ( beta ** 2.0d0 ) & & * xz_F1BasicZ(RegXMax, RegZMin+1: RegZMax-1) & & * xr_F2BasicZ(RegXMax, RegZMin+1: RegZMax-1) & & * (DelTimeShort ** 2.0d0) / ( DelZ ** 2.0d0 ) C(RegZMin+1: RegZMax-1) = & & - ( beta ** 2.0d0 ) & & * xz_F1BasicZ(RegXMax, RegZMin+2: RegZMax) & & * xr_F2BasicZ(RegXMax, RegZMin+1: RegZMax-1) & & * (DelTimeShort ** 2.0d0) / ( DelZ ** 2.0d0 ) !配列の大きさを保管 LN = RegZMax - RegZMin !係数行列を LU 分解 (LAPACK) call LinSolvLU(A, B, C, B2, IPIV) end subroutine
subroutine LinSolv (D) ! !LU 分解されている実 3 項行列の連立 1 次方程式(倍精度). LAPACK 利用 ! !暗黙の型宣言禁止 implicit none !変数定義 real(8), intent(inout) :: D(1,LN) !定数/解行列 integer :: NRHS integer :: INFO integer :: LDB character(1),parameter :: TRANS = 'N' !変数の初期化 NRHS = 0 INFO = 0 LDB = LN !解行列の計算. LAPACK を使用. call DGTTRS(TRANS, LN, NRHS, C, A, B, B2, IPIV, D, LDB, INFO) ! !解のコンディションをチェック. ! if (INFO /= 0) then ! call MessageNotify("Error", "lapack_linear", "INFO is not 0") ! stop ! end if end subroutine
subroutine LinSolvLU (A, B, C, B2, IPIV) ! !実 3 項行列の LU 分解(倍精度). LAPACK 利用 ! !暗黙の型宣言禁止 implicit none !変数定義 real(8), intent(inout) :: A(LN) !係数行列 real(8), intent(inout) :: B(LN-1) !係数行列 real(8), intent(inout) :: C(LN-1) !係数行列 real(8), intent(out) :: B2(LN-2) !係数行列 real(8), intent(out) :: IPIV(LN) !部分ピボット交換の情報を格納 integer :: INFO !変数の初期化 INFO = 0 !解行列の計算. LAPACK を使用. call DGTTRF(LN, C, A, B, B2, IPIV, INFO) ! !解のコンディションをチェック. ! if (INFO /= 0) then ! call MessageNotify("Error", "lapack_linear", "INFO is not 0") ! stop ! end if end subroutine