1 次元有限領域: チェビシェフ変換 More...
Functions/Subroutines | |
subroutine, public | at_initial (i, k, xmin, xmax) |
チェビシェフ変換の格子点数, 波数, 領域の大きさを設定する. More... | |
real(dp) function, dimension(size(at_data, 1), 0:im), public | ag_at (at_data) |
チェビシェフデータから格子データへ変換する(2 次元配列用). More... | |
real(dp) function, dimension(0:im), public | g_t (t_data) |
チェビシェフデータから格子データへ変換する(1 次元配列用). More... | |
real(dp) function, dimension(size(ag_data, 1), 0:km), public | at_ag (ag_data) |
格子データからチェビシェフデータへ変換する(2 次元配列用). More... | |
real(dp) function, dimension(0:km), public | t_g (g_data) |
格子データからチェビシェフデータへ変換する(1 次元配列用). More... | |
real(dp) function, dimension(size(at_data, 1), 0:size(at_data, 2) -1), public | at_dx_at (at_data) |
入力チェビシェフデータに X 微分を作用する(2 次元配列用). More... | |
real(dp) function, dimension(size(t_data)), public | t_dx_t (t_data) |
入力チェビシェフデータに X 微分を作用する(1 次元配列用). More... | |
real(dp) function, dimension(size(ag, 1)), public | a_int_ag (ag) |
1 次元格子点データが並んだ 2 次元配列の積分 More... | |
real(dp) function, public | int_g (g) |
1 次元格子点データの積分 More... | |
real(dp) function, dimension(size(ag, 1)), public | a_avr_ag (ag) |
1 次元格子点データが並んだ 2 次元配列の平均 More... | |
real(dp) function, public | avr_g (g) |
1 次元格子点データの平均 More... | |
real(dp) function, public | interpolate_t (t_data, xval) |
1 次元データの補間 More... | |
real(dp) function, dimension(size(at_data, 1)), public | a_interpolate_at (at_data, xval) |
1 次元データの並んだ 2次元配列の補間 More... | |
subroutine, public | at_boundaries (at_data, cond, values) |
境界条件の適用(タウ法, 2 次元配列用) More... | |
subroutine, public | t_boundaries (t_data, cond, values) |
境界条件の適用(タウ法, 1 次元配列用) 境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0 More... | |
subroutine, public | at_boundariesgrid (at_data, cond, values) |
境界条件の適用(実空間での評価, 2 次元配列用) 境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0 More... | |
subroutine, public | t_boundariesgrid (t_data, cond, values) |
境界条件の適用(タウ法, 1 次元配列用) 境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0 More... | |
subroutine, public | at_finalize |
モジュールの終了処理(割り付け配列の解放)をおこなう. More... | |
Variables | |
real(dp), dimension(:), allocatable, save, public | g_x |
格子点座標(X)を格納した 1 次元配列 More... | |
real(dp), dimension(:), allocatable, save, public | g_x_weight |
重み座標を格納した 1 次元配列 More... | |
1 次元有限領域: チェビシェフ変換
1 次元有限領域の下での流体運動をチェビシェフ変換によりスペクトル数値 計算するための Fortran90 関数を提供する.
2 次元データの 1 次元に関して同時にスペクトル計算を実行するための 関数も提供しており, 2, 3 次元領域での計算のベースも提供する.
このモジュールは内部で ISPACK3/FXPACK の Fortran77 サブルーチンを 呼んでいる. スペクトルデータおよび格子点データの格納方法については ISPACK3/FXPACK のマニュアルを参照されたい.
t_
, g_
, at_
, ag_
) は, 返す値の形を示している.t_
:: チェビシェフデータg_
:: 1 次元格子点データat_
:: 1 次元チェビシェフデータが複数並んだ 2 次元データag_
:: 1 次元格子点データが複数並んだ 2 次元データ.Dx
)は, その関数の作用を表している._e
,_at
,_g
,_ag
) は, 入力変数の形がチェビシェフデータおよび格子点データであることを示している._t
:: チェビシェフデータ_g
:: 1 次元格子点データ_at
:: 1 次元チェビシェフデータが複数並んだ 2 次元データ_ag
:: 1 次元格子点データが複数並んだ 2 次元データDP = kind(1.0D)
であるg
: 1 次元格子点データ.real(DP), dimension(0:im)
.im
は X 座標の格子点数であり, サブルーチン at_initial() にて あらかじめ設定しておく.t
: チェビシェフデータ.real(DP), dimension(0:km)
.km
は X 方向の最大波数であり, サブルーチン at_initial() にて あらかじめ設定しておく.ag
: 1 次元(X)格子点データの並んだ 2 次元データ.real(DP), dimension(:,0:im)
. 第 2 次元が X 方向を表す.at
: 1 次元チェビシェフデータの並んだ 2 次元データ.real(DP), dimension(:,0:km)
. 第 2 次元がスペクトルを表す.g_
で始まる関数が返す値は 1 次元格子点データに同じ.t_
で始まる関数が返す値はチェビシェフデータに同じ.ag_
で始まる関数が返す値は 1 次元格子点データの並んだ 2 次元データに同じ.at_
で始まる関数が返す値は 1 次元チェビシェフデータの並んだ 2 次元データに同じ.real(dp) function, dimension(size(ag,1)), public at_module::a_avr_ag | ( | real(dp), dimension(:,0:), intent(in) | ag | ) |
1 次元格子点データが並んだ 2 次元配列の平均
[in] | ag | 入力格子点データ |
Definition at line 517 of file at_module.f90.
References a_int_ag(), and g_x_weight.
real(dp) function, dimension(size(ag,1)), public at_module::a_int_ag | ( | real(dp), dimension(:,0:), intent(in) | ag | ) |
1 次元格子点データが並んだ 2 次元配列の積分
[in] | ag | 入力格子点データ |
Definition at line 483 of file at_module.f90.
References g_x_weight.
real(dp) function, dimension(size(at_data,1)), public at_module::a_interpolate_at | ( | real(dp), dimension(:,0:), intent(in) | at_data, |
real(dp), intent(in) | xval | ||
) |
1 次元データの並んだ 2次元配列の補間
[in] | at_data | 入力チェビシェフデータ |
[in] | xval | 補間する点の座標 |
Definition at line 576 of file at_module.f90.
References g_x.
real(dp) function, dimension(size(at_data,1),0:im), public at_module::ag_at | ( | real(dp), dimension(:,0:), intent(in) | at_data | ) |
チェビシェフデータから格子データへ変換する(2 次元配列用).
[in] | at_data | チェビシェフデータ |
Definition at line 328 of file at_module.f90.
real(dp) function, dimension(size(ag_data,1),0:km), public at_module::at_ag | ( | real(dp), dimension(:,0:), intent(in) | ag_data | ) |
格子データからチェビシェフデータへ変換する(2 次元配列用).
[in] | ag_data | 格子点データ |
Definition at line 374 of file at_module.f90.
subroutine, public at_module::at_boundaries | ( | real(dp), dimension(:,:), intent(inout) | at_data, |
character(len=2), intent(in), optional | cond, | ||
real(dp), dimension(:,:), intent(in), optional | values | ||
) |
境界条件の適用(タウ法, 2 次元配列用)
境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0
[in,out] | at_data | 境界条件を適用するチェビシェフデータ(m,0:km) |
[in] | cond | 境界条件.
|
[in] | values | 境界値(m,2) |
Definition at line 619 of file at_module.f90.
subroutine, public at_module::at_boundariesgrid | ( | real(dp), dimension(:,:), intent(inout) | at_data, |
character(len=2), intent(in), optional | cond, | ||
real(dp), dimension(:,:), intent(in), optional | values | ||
) |
境界条件の適用(実空間での評価, 2 次元配列用) 境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0
[in,out] | at_data | 境界条件を適用するチェビシェフデータ(m,0:km) |
[in] | cond | 境界条件.
|
[in] | values | 境界値(m,2) |
Definition at line 711 of file at_module.f90.
real(dp) function, dimension(size(at_data,1),0:size(at_data,2)-1), public at_module::at_dx_at | ( | real(dp), dimension(:,0:), intent(in) | at_data | ) |
入力チェビシェフデータに X 微分を作用する(2 次元配列用).
チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.
[in] | at_data | 入力チェビシェフデータ |
Definition at line 420 of file at_module.f90.
subroutine, public at_module::at_finalize | ( | ) |
モジュールの終了処理(割り付け配列の解放)をおこなう.
Definition at line 796 of file at_module.f90.
References ag_at(), at_dx_at(), g_x, and g_x_weight.
subroutine, public at_module::at_initial | ( | integer, intent(in) | i, |
integer, intent(in) | k, | ||
real(dp), intent(in) | xmin, | ||
real(dp), intent(in) | xmax | ||
) |
チェビシェフ変換の格子点数, 波数, 領域の大きさを設定する.
他の関数や変数を呼ぶ前に, 最初にこのサブルーチンを呼んで 初期設定をしなければならない.
[in] | i | 格子点数 |
[in] | k | 切断波数 |
[in] | xmin | 座標の範囲 |
[in] | xmax | 座標の範囲 |
Definition at line 268 of file at_module.f90.
References g_x, and g_x_weight.
real(dp) function, public at_module::avr_g | ( | real(dp), dimension(0:im), intent(in) | g | ) |
1 次元格子点データの平均
[in] | g | 格子点データ |
Definition at line 528 of file at_module.f90.
References g_x_weight, and int_g().
real(dp) function, dimension(0:im), public at_module::g_t | ( | real(dp), dimension(:), intent(in) | t_data | ) |
チェビシェフデータから格子データへ変換する(1 次元配列用).
[in] | t_data | チェビシェフデータ |
Definition at line 357 of file at_module.f90.
References ag_at().
real(dp) function, public at_module::int_g | ( | real(dp), dimension(0:im), intent(in) | g | ) |
1 次元格子点データの積分
[in] | g | 格子点データ |
Definition at line 507 of file at_module.f90.
References g_x_weight.
real(dp) function, public at_module::interpolate_t | ( | real(dp), dimension(0:), intent(in) | t_data, |
real(dp), intent(in) | xval | ||
) |
1 次元データの補間
[in] | t_data | 入力チェビシェフデータ |
[in] | xval | 補間する点の座標 |
Definition at line 539 of file at_module.f90.
References g_x.
subroutine, public at_module::t_boundaries | ( | real(dp), dimension(:), intent(inout) | t_data, |
character(len=2), intent(in), optional | cond, | ||
real(dp), dimension(2), intent(in), optional | values | ||
) |
境界条件の適用(タウ法, 1 次元配列用) 境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0
[in,out] | t_data | 境界条件を適用するチェビシェフデータ(0:km) |
[in] | cond | 境界条件.
|
[in] | values | 境界値(2) |
Definition at line 664 of file at_module.f90.
subroutine, public at_module::t_boundariesgrid | ( | real(dp), dimension(:), intent(inout) | t_data, |
character(len=2), intent(in), optional | cond, | ||
real(dp), dimension(2), intent(in), optional | values | ||
) |
境界条件の適用(タウ法, 1 次元配列用) 境界条件と両境界での値をオプションで与える. デフォルトは DD, values = 0
[in,out] | t_data | 境界条件を適用するチェビシェフデータ(0:km) |
[in] | cond | 境界条件.
|
[in] | values | 境界値(2) |
Definition at line 755 of file at_module.f90.
real(dp) function, dimension(size(t_data)), public at_module::t_dx_t | ( | real(dp), dimension(:), intent(in) | t_data | ) |
入力チェビシェフデータに X 微分を作用する(1 次元配列用).
チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.
[in] | t_data | 入力チェビシェフデータ |
Definition at line 467 of file at_module.f90.
References at_dx_at().
real(dp) function, dimension(0:km), public at_module::t_g | ( | real(dp), dimension(:), intent(in) | g_data | ) |
格子データからチェビシェフデータへ変換する(1 次元配列用).
[in] | g_data | 格子点データ |
Definition at line 399 of file at_module.f90.
References at_ag().
real(dp), dimension(:), allocatable, save, public at_module::g_x |
格子点座標(X)を格納した 1 次元配列
Definition at line 218 of file at_module.f90.
real(dp), dimension(:), allocatable, save, public at_module::g_x_weight |
重み座標を格納した 1 次元配列
Definition at line 220 of file at_module.f90.