IGMBaseLib 1.0
|
正二十面体格子点およびコントロールボリュームの座標データ等の計算を行う手続きを提供するモジュール. [詳細]
関数/サブルーチン | |
subroutine, public | construct_STDGrid (icgrid) |
STD-grid を構築する. | |
subroutine, public | calc_STDGrid_GP_vertex (icgrid) |
STD-Grid の格子点座標を計算する. | |
subroutine, public | calc_CV_vertex (icgrid) |
正二十面格子の各格子点に付随するコントロールボリュームの頂点座標を計算する. | |
subroutine, public | calc_CV_Area (icgrid) |
正二十面格子の各格子点に付随するコントロールボリュームの面積を計算する. | |
recursive subroutine | split_rectangle_region (rc_AGrid, ic_radius, row1, row2, col1, col2) |
ある矩形領域を 4 分割し, 新しく生成された 5 個の格子点の座標を計算する. | |
subroutine, public | project_on_sphere (radius, vec) |
矩形領域の平面上にある一点を原点に対して球面上へ射影する. | |
real(DP), dimension(3) | calc_each_CV_vertex (rc_AGrid, GP3_index, ic_radius, idMin) |
コントールボリュームの複数ある頂点のうちの一個の頂点座標を計算する. |
正二十面体格子点およびコントロールボリュームの座標データ等の計算を行う手続きを提供するモジュール.
Copyright (C) GFD Dennou Club, 2011-2012. All rights reserved.
license ??
subroutine,public STDGrid_Builder::calc_CV_Area | ( | type(IcGrid2D_FVM),intent(inout) | icgrid | ) |
正二十面格子の各格子点に付随するコントロールボリュームの面積を計算する.
求めたコントロールボリュームの面積は, 引数で指定された構造体 IcGrid2D_FVM の成分 rcs_CV に保存される.
[in,out] | icgrid | 構造体 IcGrid2D_FVM の変数. |
STDGrid_Builder.f90 の 314 行で定義されています。
subroutine,public STDGrid_Builder::calc_CV_vertex | ( | type(IcGrid2D_FVM),intent(inout) | icgrid | ) |
正二十面格子の各格子点に付随するコントロールボリュームの頂点座標を計算する.
求めたコントロールボリュームの頂点座標は, 引数で指定された構造体 IcGrid2D_FVM の成分 rcs_CV に保存される.
[in,out] | icgrid | 構造体 IcGrid2D_FVM の変数. |
STDGrid_Builder.f90 の 217 行で定義されています。
real(DP),dimension(3) STDGrid_Builder::calc_each_CV_vertex | ( | real(DP),dimension(idmin:,idmin:,:),intent(in) | rc_AGrid, |
integer,dimension(:,:),intent(in) | GP3_index, | ||
real(DP),intent(in) | ic_radius, | ||
integer,intent(in) | idMin | ||
) | [private] |
コントールボリュームの複数ある頂点のうちの一個の頂点座標を計算する.
[in] | idMin | 配列 rc_AGrid の第 1 or 2 添字の下限インデックス. |
[in] | rc_AGrid | 対象となるコントロールボリュームを含む矩形領域内の格子点の位置ベクトルを保持する配列. |
[in] | GP3_index | コントロールボリュームの頂点座標の計算に必要な格子点と対応する配列添字の情報を保持する配列. |
[in] | ic_radius | 正二十面体格子を内包する球の半径. |
STDGrid_Builder.f90 の 499 行で定義されています。
subroutine,public STDGrid_Builder::calc_STDGrid_GP_vertex | ( | type(IcGrid2D_FVM),intent(inout) | icgrid | ) |
STD-Grid の格子点座標を計算する.
求めた格子点座標は, 引数で指定された構造体 IcGrid2D_FVM の成分 rcs_AGrid に保存される.
[in,out] | icgrid | 構造型 IcGrid2D_FVM の変数. |
STDGrid_Builder.f90 の 133 行で定義されています。
subroutine,public STDGrid_Builder::construct_STDGrid | ( | type(IcGrid2D_FVM),intent(inout) | icgrid | ) |
subroutine,public STDGrid_Builder::project_on_sphere | ( | real(8),intent(in) | radius, |
real(8),dimension(3),intent(inout) | vec | ||
) |
矩形領域の平面上にある一点を原点に対して球面上へ射影する.
[in] | radius | 正二十面体格子を内包する球の半径. |
[in,out] | vec | 球面に射影する点の位置ベクトル. その結果は vec に代入される. |
STDGrid_Builder.f90 の 463 行で定義されています。
recursive subroutine STDGrid_Builder::split_rectangle_region | ( | real(DP),dimension(:,:,:),intent(inout) | rc_AGrid, |
real(DP),intent(in) | ic_radius, | ||
integer,intent(in) | row1, | ||
integer,intent(in) | row2, | ||
integer,intent(in) | col1, | ||
integer,intent(in) | col2 | ||
) | [private] |
ある矩形領域を 4 分割し, 新しく生成された 5 個の格子点の座標を計算する.
[in,out] | rc_AGrid | 対象となる矩形領域内の格子点の位置ベクトルを保持する配列. |
[in] | ic_radius | 正二十面体格子を内包する球の半径. |
[in] | row1 | |
[in] | row2 | |
[in] | col1 | |
[in] | col2 |
STDGrid_Builder.f90 の 398 行で定義されています。