1 次元周期境界: 実フーリエ変換 More...
Functions/Subroutines | |
subroutine, public | ae_initial (i, k, xmin, xmax) |
スペクトル変換の格子点数, 波数, 領域の大きさを設定する. More... | |
real(dp) function, dimension(size(ae, 1), 0:im-1), public | ag_ae (ae) |
スペクトルデータから格子点データへ逆変換する(2 次元データ用) More... | |
real(dp) function, dimension(0:im-1), public | g_e (e) |
スペクトルデータから格子点データへ逆変換する(1 次元データ用) More... | |
real(dp) function, dimension(size(ag, 1),-km:km), public | ae_ag (ag) |
格子点データからスペクトルデータへ正変換する(2 次元データ用) More... | |
real(dp) function, dimension(-km:km), public | e_g (g) |
格子点データからスペクトルデータへ正変換する(1 次元データ用) More... | |
real(dp) function, dimension(size(ae, 1),-km:km), public | ae_dx_ae (ae) |
入力スペクトルデータに X 微分を作用する(2 次元データ). More... | |
real(dp) function, dimension(-km:km), public | e_dx_e (e) |
入力スペクトルデータに 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_e (e_data, xval) |
補間計算(1次元) More... | |
real(dp) function, dimension(size(ae_data, 1)), public | a_interpolate_ae (ae_data, xval) |
補間計算(2次元) More... | |
Variables | |
real(dp), dimension(:), allocatable, save, public | g_x |
格子点座標(X)を格納した1次元配列. More... | |
real(dp), dimension(:), allocatable, save, public | g_x_weight |
重み座標を格納した 1 次元配列. X方向の格子点の間隔が格納してある. More... | |
1 次元周期境界: 実フーリエ変換
1 次元周期境界条件の下での流体運動を 実フーリエ変換によるスペクトル法で数値計算するための Fortran90 関数を提供する. 2 次元データの 1 次元に関して同時にスペクトル計算を実行するための 関数も提供しており, 2, 3 次元領域での計算のベースも提供する.
このモジュールは内部で ISPACK3/FXPACK のサブルーチンを呼んでいる. スペクトルデータの格納方法については ISPACK3/FXPACK と異なっているので以下のコメントに注意されたい.
e_
, g_
, ae_
, ag_
) は, 返す値の形を示している.e_
:: スペクトルデータ,g_
:: 1 次元格子点データ,ae_
:: 1 次元スペクトルデータが複数並んだ 2 次元データ,ag_
:: 1 次元格子点データが複数並んだ 2 次元データ.Dx
)は, その関数の作用を表している._e
, _ae
, _g
, _ag
) は, 入力変数の形スペクトルデータおよび格子点データであることを示している._e
:: スペクトルデータ_g
:: 1 次元格子点データ_ae
:: 1 次元スペクトルデータが複数並んだ 2 次元データ_ag
:: 1 次元格子点データが複数並んだ 2 次元データDP = kind(1.0D)
であるg
: 1 次元格子点データ.real(DP), dimension(0:im-1)
.im
は X 座標の格子点数であり, サブルーチン ae_initial() にて あらかじめ設定しておく.e
: スペクトルデータ.real(DP), dimension(-km:km)
.km
は X 方向の最大波数であり, サブルーチン ae_initial() にて あらかじめ設定しておく.0:km
までが \(\cos(kx)\) の係数, -km:-1
が \(sin(kx)\) の係数となっている.ag
: 1 次元(X)格子点データの並んだ 2 次元データ.real(DP), dimension(:,0:im-1)
. 第 2 次元が X 方向を表す.ae
: 1 次元スペクトルデータの並んだ 2 次元データ.real(DP), dimension(:,-km:km)
. 第 2 次元がスペクトルを表す.g_
で始まる関数が返す値は 1 次元格子点データに同じ.e_
で始まる関数が返す値はスペクトルデータに同じ.ag_
で始まる関数が返す値は 1 次元格子点データの並んだ 2 次元データに同じ.ae_
で始まる関数が返す値は 1 次元スペクトルデータの並んだ 2 次元データに同じ.real(dp) function, dimension(size(ag,1)), public ae_module::a_avr_ag | ( | real(dp), dimension(:,0:), intent(in) | ag | ) |
1 次元格子点データが並んだ 2 次元配列の平均
[in] | ag | 格子点データ |
Definition at line 388 of file ae_module.f90.
References a_int_ag(), and g_x_weight.
real(dp) function, dimension(size(ag,1)), public ae_module::a_int_ag | ( | real(dp), dimension(:,0:), intent(in) | ag | ) |
1 次元格子点データが並んだ 2 次元配列の積分
[in] | ag | 格子点データ |
Definition at line 352 of file ae_module.f90.
References g_x_weight.
real(dp) function, dimension(size(ae_data,1)), public ae_module::a_interpolate_ae | ( | real(dp), dimension(:,-km:), intent(in) | ae_data, |
real(dp), intent(in) | xval | ||
) |
補間計算(2次元)
[in] | ae_data | 入力スペクトルデータ |
[in] | xval | 補間する点の座標 |
Definition at line 433 of file ae_module.f90.
real(dp) function, dimension(size(ag,1),-km:km), public ae_module::ae_ag | ( | real(dp), dimension(:,:), intent(in) | ag | ) |
格子点データからスペクトルデータへ正変換する(2 次元データ用)
スペクトル正変換の定義は以下のとおり:
\begin{align*} e_0 & = (1/im)\sum_{j=0}^{im-1} g_j \\ e_k & = (1/im)\sum_{j=0}^{im-1} g_j \cos(2\pi jk/im) \quad (k=1,2,...,km) \\ e_{-k} & = - (1/im)\sum_{j=0}^{im-1} g_j \sin(2\pi jk/im) \quad (k=1,2,...,km) \end{align*}
[in] | ag | 格子点データ |
Definition at line 245 of file ae_module.f90.
real(dp) function, dimension(size(ae,1),-km:km), public ae_module::ae_dx_ae | ( | real(dp), dimension(:,-km:), intent(in) | ae | ) |
入力スペクトルデータに X 微分を作用する(2 次元データ).
スペクトルデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのスペクトル変換のことである.
[in] | ae | 入力スペクトルデータ |
Definition at line 310 of file ae_module.f90.
subroutine, public ae_module::ae_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 | X座標下端 |
[in] | xmax | X座標上端 |
Definition at line 122 of file ae_module.f90.
References g_x, and g_x_weight.
real(dp) function, dimension(size(ae,1),0:im-1), public ae_module::ag_ae | ( | real(dp), dimension(:,-km:), intent(in) | ae | ) |
スペクトルデータから格子点データへ逆変換する(2 次元データ用)
スペクトル逆変換の定義は以下のとおり:
\begin{align*} g_j = e_0 + 2\sum_{k=1}^{km}e_k\cos(2\pi jk/im) - e_{-k}\sin(2\pi jk/im)) \quad (j=0,1,...,im-1). \end{align*}
[in] | ae | スペクトルデータ |
Definition at line 176 of file ae_module.f90.
real(dp) function, public ae_module::avr_g | ( | real(dp), dimension(0:im-1), intent(in) | g | ) |
1 次元格子点データの平均
[in] | g | 格子点データ |
Definition at line 400 of file ae_module.f90.
References g_x_weight, and int_g().
real(dp) function, dimension(-km:km), public ae_module::e_dx_e | ( | real(dp), dimension(-km:km), intent(in) | e | ) |
入力スペクトルデータに X 微分を作用する(1 次元データ).
スペクトルデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのスペクトル変換のことである.
[in] | e | 入力スペクトルデータ |
Definition at line 336 of file ae_module.f90.
References ae_dx_ae().
real(dp) function, dimension(-km:km), public ae_module::e_g | ( | real(dp), dimension(0:im-1), intent(in) | g | ) |
格子点データからスペクトルデータへ正変換する(1 次元データ用)
スペクトル正変換の定義は以下のとおり:
\begin{align*} e_0 & = (1/im)\sum_{j=0}^{im-1} g_j \\ e_k & = (1/im)\sum_{j=0}^{im-1} g_j \cos(2\pi jk/im) \quad (k=1,2,...,km) \\ e_{-k} & = - (1/im)\sum_{j=0}^{im-1} g_j \sin(2\pi jk/im) \quad (k=1,2,...,km) \end{align*}
[in] | g | 格子点データ |
Definition at line 290 of file ae_module.f90.
References ae_ag().
real(dp) function, dimension(0:im-1), public ae_module::g_e | ( | real(dp), dimension(-km:km), intent(in) | e | ) |
スペクトルデータから格子点データへ逆変換する(1 次元データ用)
\begin{align*} g_j = e_0 + 2\sum_{k=1}^{km}(e_k\cos(2\pi jk/im) - e_{-k}\sin(2\pi jk/im)) \quad (j=0,1,...,im-1). \end{align*}
[in] | e | スペクトルデータ |
Definition at line 217 of file ae_module.f90.
References ag_ae().
real(dp) function, public ae_module::int_g | ( | real(dp), dimension(0:im-1), intent(in) | g | ) |
1 次元格子点データの積分
[in] | g | 格子点データ |
Definition at line 376 of file ae_module.f90.
References g_x_weight.
real(dp) function, public ae_module::interpolate_e | ( | real(dp), dimension(-km:km), intent(in) | e_data, |
real(dp), intent(in) | xval | ||
) |
補間計算(1次元)
[in] | e_data | 入力スペクトルデータ |
[in] | xval | 補間する点の座標 |
Definition at line 412 of file ae_module.f90.
real(dp), dimension(:), allocatable, save, public ae_module::g_x |
格子点座標(X)を格納した1次元配列.
Definition at line 102 of file ae_module.f90.
real(dp), dimension(:), allocatable, save, public ae_module::g_x_weight |
重み座標を格納した 1 次元配列. X方向の格子点の間隔が格納してある.
Definition at line 104 of file ae_module.f90.