IGModel-SW 1.0
|
速度場の時間変化率を計算する(運動方程式を解く)手続きを提供するモジュール. [詳細]
関数/サブルーチン | |
subroutine, public | motion_equation_Init (icgrid) |
motion_equation モジュールの初期化を行う. | |
subroutine, public | calc_motion_eq_dvdt (DVelDtN, xy_VelN, xy_HeightN, xy_Coli, diff_eval) |
運動方程式の半離散式の右辺を評価し, 全格子点に対してタイムレベル n の速度変化率を計算する. | |
real(DP), dimension(3) | local_k (s_pos) |
地理座標系の鉛直方向の単位ベクトルを直交座標系の単位ベクトルを使って表す. | |
変数 | |
type(Field_IcGrid2D), save | xy_Zeta |
相対渦度場データを管理する構造型 Field_IcGrid2D の変数. | |
type(Field_IcGrid2D), save | mecha_energy |
力学的エネルギーの物理場データを管理する構造型 Field_IcGrid2D の変数. | |
type(Field_IcGrid2D), save | grad_mecha_energy |
ジオポテンシャルと運動エネルギーの和の勾配場を管理する構造型 Field_IcGrid2D の変数. |
速度場の時間変化率を計算する(運動方程式を解く)手続きを提供するモジュール.
subroutine,public motion_equation::calc_motion_eq_dvdt | ( | type(Field_IcGrid2D),intent(inout) | DVelDtN, |
type(Field_IcGrid2D),intent(in) | xy_VelN, | ||
type(Field_IcGrid2D),intent(in) | xy_HeightN, | ||
type(Field_IcGrid2D),intent(in) | xy_Coli, | ||
type(Derivate_Field_IcGrid2D),intent(inout) | diff_eval | ||
) |
運動方程式の半離散式の右辺を評価し, 全格子点に対してタイムレベル n の速度変化率を計算する.
具体的には, このサブルーチンでは, 以下の Tomita, etal(2001) の式(6) の半離散式の右辺を評価する.
ここで, 添字の はタイムレベル n の物理場であることを表す.
[in,out] | DVelDtN | タイムレベル n の速度場の時間変化率を管理する構造型 Field_IcGrid2D の変数. |
[in] | xy_VelN | タイムレベル n の速度場を管理する構造型 Field_IcGrid2D の変数. |
[in] | xy_HeightN | タイムレベル n の表面高度場を管理する構造型 Field_IcGrid2D の変数. |
[in] | xy_Coli | |
[in] | diff_eval | 正二十面格子上の物理場の微分演算に対するデータを管理する構造型 Derivate_Field_IcGrid2D の変数. |
motion_equation.f90 の 186 行で定義されています。
real(DP),dimension(3) motion_equation::local_k | ( | real(DP),dimension(3),intent(in) | s_pos | ) | [private] |
地理座標系の鉛直方向の単位ベクトルを直交座標系の単位ベクトルを使って表す.
[in] | s_pos | 地理座標系における位置ベクトル . |
motion_equation.f90 の 305 行で定義されています。
subroutine,public motion_equation::motion_equation_Init | ( | type(IcGrid2D_FVM),intent(inout) | icgrid | ) |
motion_equation モジュールの初期化を行う.
[in,out] | icgrid | 構造型 IcGrid2D_FVM の変数. |
motion_equation.f90 の 122 行で定義されています。
type(Field_IcGrid2D),save motion_equation::grad_mecha_energy |
ジオポテンシャルと運動エネルギーの和の勾配場を管理する構造型 Field_IcGrid2D の変数.
motion_equation.f90 の 107 行で定義されています。
type(Field_IcGrid2D),save motion_equation::mecha_energy |
力学的エネルギーの物理場データを管理する構造型 Field_IcGrid2D の変数.
motion_equation.f90 の 103 行で定義されています。
type(Field_IcGrid2D),save motion_equation::xy_Zeta |
相対渦度場データを管理する構造型 Field_IcGrid2D の変数.
motion_equation.f90 の 99 行で定義されています。