IGMBaseLib 1.0
関数/サブルーチン | 変数

SPRGCGrid_Builderモジュール

STD-Grid をバネ力学を用いて改良した SPR-GC-grid を作成する手続きを提供するモジュール. [詳細]

関数/サブルーチン

subroutine, public construct_SPRGCGrid (icgrid, tstep_limit)
 SPR-GC-grid を構築する.
subroutine, public set_iteration_limiter (limiter_flag)
 時間積分に伴う逐次計算ループの反復回数を制限するか否かのフラグを設定する.
subroutine modify_STDGrid_with_spring_dyn (icgrid)
 STD-grid を修正することによって得られる, SPR-GC-grid を生成する.
subroutine move_GP_with_spring_dyn (icgrid, idMin, idMax, EMin, EMax)
 全球領域の各格子点に対してバネ力学を適用し, 仮想的なつり合い位置を求める.
subroutine check_global_convergence (icgrid, conv_total_ave, conv_max, EMin, EMax)
 バネでつながれた格子点の振動が, 全球で収束したかを調べる.
real(8) check_convergence (r0, pts, pt_num)
 収束度合いを, Tomita etal, 2001 の式(22) を用いてチェックする.
integer generate_near_GP6index (rcID, GP0_i, GP0_j, EMin, EMax)
 格子点(rcID, GP0_i, GP0_j) の近傍にある格子点 6 個に対応する, 配列のインデックスを生成する.
integer generate_near_GP5index (icgrid, rcID, GP0_i, GP0_j)
 格子点(rcID, GP0_i, GP0_j) の近傍にある格子点 5 個に対応する, 配列のインデックスを生成する.
subroutine integ_ode_1step_spring_dynamics (r0, w0, GP, pt_num)
 固定された 6 点とバネでつながれた質点の位置の時間発展を計算する.
real(DP), dimension(3) f1 (r0, w0)
 以下の常微分方程式の( Tomita etal, 2001 の式(21) ) の右辺の値を計算する.

\[ \DD{\Dvect{r_0}}{t} = f_1(\Dvect{r}_0, \Dvect{w}_0) = \Dvect{w}_0 \]


real(8), dimension(1:3) f2 (r0, w0, pts, pt_num, alpha, k, eq_d, radius)
 以下の常微分方程式の( Tomita etal, 2001 の式(20) ) の右辺の値を計算する.

\[ M \DD{\Dvect{w_0}}{t} = f_2(\Dvect{r}_0, \Dvect{w}_0) \ = \sum_{i=1}^6 k(d_i - \bar{d}) \Dvect{e}_i - \alpha \Dvect{w}_0 \]


real(8), dimension(3) calc_ei (r0, ri)
 点 G_0(r0) から 点 Gi(ri) を結ぶベクトルの単位ベクトルを取得する.

変数

real(DP) eq_d
 力が課されていない時のバネの長さ.
real(DP) radius
 正二十面体を内包する球の半径.
real(DP) beta = 0.4d0
 チューニングパラメータ.
real(DP) alpha = 2000.0d0
 摩擦定数.
real(DP) k = 1.0d06
 バネ定数.
real(DP) dt = 0.001d0
 時間刻み幅.
integer end_tstep = 20000
 時間積分のステップ数.
logical tstep_limiter = .false.
 時間積分に伴う逐次計算ループの反復回数を制限するか否かのフラグ. デフォルトは, ループ回数を制限しない.
integer conv_check_step = 100
 収束判定を行うステップ間隔
integer, allocatable GP_IJ
 ある格子点を取り囲む格子点(6 個 or 5 個)のインデックスを保持する配列.

説明

STD-Grid をバネ力学を用いて改良した SPR-GC-grid を作成する手続きを提供するモジュール.



Copyright (C) GFD Dennou Club, 2011-2012. All rights reserved.
license ??

作者:
Yuta Kawai

関数/サブルーチン

real(8),dimension(3) SPRGCGrid_Builder::calc_ei ( real(8),dimension(3),intent(in)  r0,
real(8),dimension(3),intent(in)  ri 
) [private]

点 G_0(r0) から 点 Gi(ri) を結ぶベクトルの単位ベクトルを取得する.

引数:
[in]r0点 r0 の位置ベクトル.
[in]ri点 ri の位置ベクトル.
戻り値:
$\Dvect{r}_i-\Dvect{r}_0$の単位ベクトル.

SPRGCGrid_Builder.f90894 行で定義されています。

real(8) SPRGCGrid_Builder::check_convergence ( real(8),dimension(3),intent(in)  r0,
real(8),dimension(pt_num, 3),intent(in)  pts,
integer,intent(in)  pt_num 
) [private]

収束度合いを, Tomita etal, 2001 の式(22) を用いてチェックする.

戻りがゼロに近ければ近いほど, 質点はつり合いの位置に存在することになる.

引数:
[in]r0点 P0 の位置ベクトル.
[in]pts点 P0 を取り囲む点の位置ベクトル.
[in]pt_num点 P0 を取り囲む近傍の格子点数.
戻り値:
Tomita etal, 2001 の式(22)の左辺の値.

SPRGCGrid_Builder.f90560 行で定義されています。

subroutine SPRGCGrid_Builder::check_global_convergence ( type(IcGrid2D_FVM),intent(in)  icgrid,
real(DP),intent(out)  conv_total_ave,
real(DP),intent(out)  conv_max,
integer,intent(in)  EMin,
integer,intent(in)  EMax 
) [private]

バネでつながれた格子点の振動が, 全球で収束したかを調べる.

引数:
[in]conv_toal_ave
[in]conv_max
[in]EMin正二十面格子の座標データを保持する配列の袖領域も含めないインデックスの最小値.
[in]EMax正二十面格子の座標データを保持する配列の袖領域も含めないインデックスの最大値.

SPRGCGrid_Builder.f90482 行で定義されています。

subroutine,public SPRGCGrid_Builder::construct_SPRGCGrid ( type(IcGrid2D_FVM),intent(inout)  icgrid,
integer,intent(in),optional  tstep_limit 
)

SPR-GC-grid を構築する.

引数:
[in,out]icgrid構造型 IcGrid2D_FVM の変数.

SPRGCGrid_Builder.f90156 行で定義されています。

real(DP),dimension(3) SPRGCGrid_Builder::f1 ( real(DP),dimension(3),intent(in)  r0,
real(DP),dimension(3),intent(in)  w0 
) [private]

以下の常微分方程式の( Tomita etal, 2001 の式(21) ) の右辺の値を計算する.

\[ \DD{\Dvect{r_0}}{t} = f_1(\Dvect{r}_0, \Dvect{w}_0) = \Dvect{w}_0 \]

引数:
[in]r0点 r0 の位置ベクトル.
[in]w0点 r0 の速度ベクトル.
戻り値:
関数 f1 の値.

SPRGCGrid_Builder.f90793 行で定義されています。

real(8),dimension(1:3) SPRGCGrid_Builder::f2 ( real(DP),dimension(3),intent(in)  r0,
real(DP),dimension(3),intent(in)  w0,
real(DP),dimension(pt_num,3),intent(in)  pts,
integer,intent(in)  pt_num,
real(DP),intent(in)  alpha,
real(DP),intent(in)  k,
real(DP),intent(in)  eq_d,
real(DP),intent(in)  radius 
) [private]

以下の常微分方程式の( Tomita etal, 2001 の式(20) ) の右辺の値を計算する.

\[ M \DD{\Dvect{w_0}}{t} = f_2(\Dvect{r}_0, \Dvect{w}_0) \ = \sum_{i=1}^6 k(d_i - \bar{d}) \Dvect{e}_i - \alpha \Dvect{w}_0 \]

引数:
[in]r0点 r0 の位置ベクトル.
[in]w0点 r0 の速度ベクトル.
[in]pts
[in]pt_num
[in]alpha
[in]k
[in]eq_d
[in]radius正二十面体を内包する球の半径.
戻り値:
 関数 f2 の値.

SPRGCGrid_Builder.f90842 行で定義されています。

integer SPRGCGrid_Builder::generate_near_GP5index ( type(IcGrid2D_FVM),intent(inout)  icgrid,
integer,intent(in)  rcID,
integer,intent(in)  GP0_i,
integer,intent(in)  GP0_j 
) [private]

格子点(rcID, GP0_i, GP0_j) の近傍にある格子点 5 個に対応する, 配列のインデックスを生成する.

引数:
[in,out]icgrid構造体 IcGrid2D_FVM の変数.
[in]rcID
[in]GP0_i
[in]GP0_j

SPRGCGrid_Builder.f90667 行で定義されています。

integer SPRGCGrid_Builder::generate_near_GP6index ( integer,intent(in)  rcID,
integer,intent(in)  GP0_i,
integer,intent(in)  GP0_j,
integer  EMin,
integer  EMax 
) [private]

格子点(rcID, GP0_i, GP0_j) の近傍にある格子点 6 個に対応する, 配列のインデックスを生成する.

引数:
[in]rcID
[in]GP0_i
[in]GP0_j
[in]EMin
[in]EMax

SPRGCGrid_Builder.f90609 行で定義されています。

subroutine SPRGCGrid_Builder::integ_ode_1step_spring_dynamics ( real(DP),dimension(3),intent(inout)  r0,
real(DP),dimension(3),intent(inout)  w0,
real(DP),dimension(pt_num,3),intent(in)  GP,
integer,intent(in)  pt_num 
) [private]

固定された 6 点とバネでつながれた質点の位置の時間発展を計算する.

常微分方程式系の時間積分の方法は, 3 次の Runge=Kutta 法である.

引数:
[in,out]r0
[in,out]w0
[in]GP
[in]pt_num
[in]dt
[in]alpha
[in]k
[in]eq_d
[in]ic_radius正二十面体を内包する球の半径.

SPRGCGrid_Builder.f90724 行で定義されています。

subroutine SPRGCGrid_Builder::modify_STDGrid_with_spring_dyn ( type(IcGrid2D_FVM),intent(inout)  icgrid) [private]

STD-grid を修正することによって得られる, SPR-GC-grid を生成する.

具体的には, バネ力学を用いて STD-grid の配置を修正する. その後, STD-GC-grid と同様の方法で, コントロールボリュームを定義して, コントロールボリュームの重心に 格子点を移動させる.

このバネ力学を用いた修正により, コントロールボリュームの面積・歪みが, 全球領域で単調に変化するようになる. STD-grid の場合, コントロールボリュームの面積・歪みがフラクタル的な分布をしているために, 数値的に微分値を評価する際, 組織的な数値誤差が発生してしまうが, この修正により数値誤差は改善される.

引数:
[in,out]icgrid構造体 IcGrid2D_FVM の変数.

SPRGCGrid_Builder.f90268 行で定義されています。

subroutine SPRGCGrid_Builder::move_GP_with_spring_dyn ( type(IcGrid2D_FVM),intent(inout)  icgrid,
integer,intent(in)  idMin,
integer,intent(in)  idMax,
integer,intent(in)  EMin,
integer,intent(in)  EMax 
) [private]

全球領域の各格子点に対してバネ力学を適用し, 仮想的なつり合い位置を求める.

各格子点を仮想的な質点とみなし, それぞれ近傍の質点 6 個とを仮想的なバネでつなぐ. そして, 全質点がつり合いの状態にある位置を求める.

引数:
[in]idMin正二十面格子の座標データを保持する配列の袖領域も含めたインデックスの最小値.
[in]idMax正二十面格子の座標データを保持する配列の袖領域も含めたインデックスの最大値.
[in]EMin正二十面格子の座標データを保持する配列の袖領域も含めないインデックスの最小値.
[in]EMax正二十面格子の座標データを保持する配列の袖領域も含めないインデックスの最大値.

SPRGCGrid_Builder.f90348 行で定義されています。

subroutine,public SPRGCGrid_Builder::set_iteration_limiter ( logical,intent(in)  limiter_flag)

時間積分に伴う逐次計算ループの反復回数を制限するか否かのフラグを設定する.

デフォルトではフラグ tstep_limiter は真であり, 最大 2 万回で逐次計算ループを抜ける. このフラグを偽に設定した場合には, check_convergence による収束判定が真を返すまで, 逐次計算ループを続行される.

引数:
[in]limiter_flag反復回数を制限するか否かのフラグ.

SPRGCGrid_Builder.f90221 行で定義されています。


変数

real(DP) SPRGCGrid_Builder::alpha = 2000.0d0

摩擦定数.

SPRGCGrid_Builder.f90112 行で定義されています。

real(DP) SPRGCGrid_Builder::beta = 0.4d0

チューニングパラメータ.

SPRGCGrid_Builder.f90105 行で定義されています。

収束判定を行うステップ間隔

SPRGCGrid_Builder.f90137 行で定義されています。

real(DP) SPRGCGrid_Builder::dt = 0.001d0

時間刻み幅.

SPRGCGrid_Builder.f90120 行で定義されています。

時間積分のステップ数.

SPRGCGrid_Builder.f90124 行で定義されています。

力が課されていない時のバネの長さ.

SPRGCGrid_Builder.f9097 行で定義されています。

integer,allocatable SPRGCGrid_Builder::GP_IJ

ある格子点を取り囲む格子点(6 個 or 5 個)のインデックスを保持する配列.

SPRGCGrid_Builder.f90141 行で定義されています。

real(DP) SPRGCGrid_Builder::k = 1.0d06

バネ定数.

SPRGCGrid_Builder.f90116 行で定義されています。

正二十面体を内包する球の半径.

SPRGCGrid_Builder.f90101 行で定義されています。

時間積分に伴う逐次計算ループの反復回数を制限するか否かのフラグ. デフォルトは, ループ回数を制限しない.

SPRGCGrid_Builder.f90130 行で定義されています。

 全て クラス ネームスペース ファイル 関数 変数