Functions/Subroutines | Variables
ae_module Module Reference

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...
 

Detailed Description

1 次元周期境界: 実フーリエ変換

Author
Shin-ichi Takehiro, Youhei SASAKI

1 次元周期境界条件の下での流体運動を 実フーリエ変換によるスペクトル法で数値計算するための Fortran90 関数を提供する. 2 次元データの 1 次元に関して同時にスペクトル計算を実行するための 関数も提供しており, 2, 3 次元領域での計算のベースも提供する.

このモジュールは内部で ISPACK3/FXPACK のサブルーチンを呼んでいる. スペクトルデータの格納方法については ISPACK3/FXPACK と異なっているので以下のコメントに注意されたい.

命名法

各データの種類の説明

Function/Subroutine Documentation

◆ a_avr_ag()

real(dp) function, dimension(size(ag,1)), public ae_module::a_avr_ag ( real(dp), dimension(:,0:), intent(in)  ag)

1 次元格子点データが並んだ 2 次元配列の平均

Parameters
[in]ag格子点データ
Returns
平均結果

Definition at line 388 of file ae_module.f90.

References a_int_ag(), and g_x_weight.

388 
389  real(DP), dimension(:,0:), intent(in) :: ag
391  real(DP), dimension(size(ag,1)) :: a_avr_ag
392 
393  a_avr_ag = a_int_ag(ag)/sum(g_x_weight)
394 
Here is the call graph for this function:

◆ a_int_ag()

real(dp) function, dimension(size(ag,1)), public ae_module::a_int_ag ( real(dp), dimension(:,0:), intent(in)  ag)

1 次元格子点データが並んだ 2 次元配列の積分

Parameters
[in]ag格子点データ
Returns
積分結果

Definition at line 352 of file ae_module.f90.

References g_x_weight.

352 
353  real(DP), dimension(:,0:), intent(in) :: ag
355  real(DP), dimension(size(ag,1)) :: a_int_ag
356 
357  integer(8) :: i
358 
359  if ( size(ag,2) < im ) then
360  call messagenotify('E','ae_Int_ag', &
361  'The Grid points of input data too small.')
362  elseif ( size(ag,2) > im ) then
363  call messagenotify('W','ae_Int_ag', &
364  'The Grid points of input data too large.')
365  endif
366 
367  a_int_ag = 0.0d0
368  do i=0,im-1
369  a_int_ag(:) = a_int_ag(:) + ag(:,i)*g_x_weight(i)
370  enddo

◆ a_interpolate_ae()

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次元)

Parameters
[in]ae_data入力スペクトルデータ
[in]xval補間する点の座標
Returns
補間した結果の値

Definition at line 433 of file ae_module.f90.

433 
434  real(DP), dimension(:,-km:), intent(IN) :: ae_data
436  real(DP), intent(in) :: xval
438  real(DP), dimension(size(ae_data,1)) :: a_interpolate_ae
439 
440  integer(8) :: k
441 
442  a_interpolate_ae = ae_data(:,0)
443 
444  do k=1,km
445  a_interpolate_ae = a_interpolate_ae &
446  & + 2.0d0 * ( ae_data(:,k) * cos(2*pi*k/xl*xval) &
447  & - ae_data(:,-k) * sin(2*pi*k/xl*xval) )
448  end do
449 

◆ ae_ag()

real(dp) function, dimension(size(ag,1),-km:km), public ae_module::ae_ag ( real(dp), dimension(:,:), intent(in)  ag)

格子点データからスペクトルデータへ正変換する(2 次元データ用)

スペクトル正変換の定義は以下のとおり:

  • 格子点データ \(g_j \, (j=-0,...,im-1)\)
  • スペクトルデータ \(e_k \, (k=-km,...,km)\)

    \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*}

Parameters
[in]ag格子点データ
Returns
スペクトルデータ

Definition at line 245 of file ae_module.f90.

245 
246  real(DP), dimension(:,:), intent(in) :: ag
248  real(DP), dimension(size(ag,1),-km:km) :: ae_ag
249 
250  real(DP), dimension(size(ag,1),0:im-1) :: ag_work
251  integer(8) :: nm, k
252 
253  nm = size(ag,1)
254  if ( size(ag,2) < im ) then
255  call messagenotify('E','ae_ag', &
256  'The Grid points of input data too small.')
257  elseif ( size(ag,2) > im ) then
258  call messagenotify('W','ae_ag', &
259  'The Grid points of input data too large.')
260  endif
261  ag_work = ag
262 
263  call fxrtfa(nm,im,ag_work,it,t)
264 
265  do k=1,min(km,im/2-1)
266  ae_ag(:, k) = ag_work(:,2*k)
267  ae_ag(:,-k) = ag_work(:,2*k+1)
268  enddo
269  ae_ag(:,0) = ag_work(:,0)
270  if ( km == im/2 ) then
271  ae_ag(:,km) = ag_work(:,1)
272  end if
273 

◆ ae_dx_ae()

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 微分を 作用させたデータのスペクトル変換のことである.

Parameters
[in]ae入力スペクトルデータ
Returns
入力スペクトルデータの X 微分

Definition at line 310 of file ae_module.f90.

310 
311  real(DP), dimension(:,-km:), intent(in) :: ae
313  real(DP), dimension(size(ae,1),-km:km) :: ae_dx_ae
314 
315  integer(8) :: k
316 
317  if ( size(ae,2) < 2*km+1 ) then
318  call messagenotify('W','ae_Dx_ae', &
319  'The Fourier dimension of input data too small.')
320  elseif ( size(ae,2) > 2*km+1 ) then
321  call messagenotify('W','ae_Dx_ae', &
322  'The Fourier dimension of input data too large.')
323  endif
324 
325  do k=-km,km
326  ae_dx_ae(:,k) = -(2*pi*k/xl)*ae(:,-k)
327  enddo

◆ ae_initial()

subroutine, public ae_module::ae_initial ( integer, intent(in)  i,
integer, intent(in)  k,
real(dp), intent(in)  xmin,
real(dp), intent(in)  xmax 
)

スペクトル変換の格子点数, 波数, 領域の大きさを設定する.

他の関数や変数を呼ぶ前に, 最初にこのサブルーチンを呼んで 初期設定をしなければならない.

Parameters
[in]i格子点の数
[in]k切断波数
[in]xminX座標下端
[in]xmaxX座標上端

Definition at line 122 of file ae_module.f90.

References g_x, and g_x_weight.

122 
123  integer, intent(in) :: i
125  integer, intent(in) :: k
127  real(DP),intent(in) :: xmin
129  real(DP),intent(in) :: xmax
130  integer(8) :: ii
131 
132  im = i
133  km = k
134  xl = xmax - xmin
135 
136  if ( im <= 0 .or. km <= 0 ) then
137  call messagenotify('E','ae_Initial', &
138  'Number of grid points and waves should be positive')
139  elseif ( mod(im,2_8) /= 0 ) then
140  call messagenotify('E','ae_Initial', &
141  'Number of grid points should be even')
142  elseif ( km > im/2 ) then
143  call messagenotify('E','ae_Initial', &
144  'KM shoud be less equal to IM/2')
145  endif
146 
147  allocate(it(im/2),t(im*3/2))
148 
149  call fxrini(im,it,t)
150 
151  allocate(g_x(0:im-1))
152  do ii=0,im-1
153  g_x(ii) = xmin + (xl/im)*ii
154  enddo
155 
156  allocate(g_x_weight(0:im-1))
157  g_x_weight = xl/im
158 
159  call messagenotify(&
160  'M','ae_initial','ae_module (2019/09/17) is initialized')
161 

◆ ag_ae()

real(dp) function, dimension(size(ae,1),0:im-1), public ae_module::ag_ae ( real(dp), dimension(:,-km:), intent(in)  ae)

スペクトルデータから格子点データへ逆変換する(2 次元データ用)

スペクトル逆変換の定義は以下のとおり:

  • スペクトルデータ: \(e_k \, (k=-km,...,km)\)
  • 格子点データ: \(g_j \, (j=-0,...,im-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*}

Parameters
[in]aeスペクトルデータ
Returns
格子点データ

Definition at line 176 of file ae_module.f90.

176 
177  real(DP), dimension(:,-km:), intent(in) :: ae
179  real(DP), dimension(size(ae,1),0:im-1) :: ag_ae
180 
181  integer(8) :: nm, k
182  nm=size(ae,1)
183 
184  if ( size(ae,2) < 2*km+1 ) then
185  call messagenotify('E','ag_ae', &
186  'The Fourier dimension of input data too small.')
187  elseif ( size(ae,2) > 2*km+1 ) then
188  call messagenotify('W','ag_ae', &
189  'The Fourier dimension of input data too large.')
190  endif
191 
192  ag_ae = 0.0d0
193  ag_ae(:,0) = ae(:,0)
194  do k=1,min(km,im/2-1)
195  ag_ae(:,2*k) = ae(:,k)
196  ag_ae(:,2*k+1) = ae(:,-k)
197  enddo
198  if ( km == im/2 ) then
199  ag_ae(:,1) = ae(:,km)
200  end if
201 
202  call fxrtba(nm,im,ag_ae,it,t)
203 

◆ avr_g()

real(dp) function, public ae_module::avr_g ( real(dp), dimension(0:im-1), intent(in)  g)

1 次元格子点データの平均

Parameters
[in]g格子点データ
Returns
平均値

Definition at line 400 of file ae_module.f90.

References g_x_weight, and int_g().

400 
401  real(DP), dimension(0:im-1), intent(in) :: g
403  real(DP) :: avr_g
404 
405  avr_g = int_g(g)/sum(g_x_weight)
406 
Here is the call graph for this function:

◆ e_dx_e()

real(dp) function, dimension(-km:km), public ae_module::e_dx_e ( real(dp), dimension(-km:km), intent(in)  e)

入力スペクトルデータに X 微分を作用する(1 次元データ).

スペクトルデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのスペクトル変換のことである.

Parameters
[in]e入力スペクトルデータ
Returns
入力スペクトルデータの X 微分

Definition at line 336 of file ae_module.f90.

References ae_dx_ae().

336 
337  real(DP), dimension(-km:km), intent(in) :: e
339  real(DP), dimension(-km:km) :: e_dx_e
340 
341  real(DP), dimension(1,-km:km) :: ae_work
342 
343  ae_work(1,:) = e
344  ae_work = ae_dx_ae(ae_work)
345  e_dx_e = ae_work(1,:)
346 
Here is the call graph for this function:

◆ e_g()

real(dp) function, dimension(-km:km), public ae_module::e_g ( real(dp), dimension(0:im-1), intent(in)  g)

格子点データからスペクトルデータへ正変換する(1 次元データ用)

スペクトル正変換の定義は以下のとおり:

  • 格子点データ: \(g_j \, (j=-0,...,im-1)\)
  • スペクトルデータ: \(e_k \, (k=-km,...,km)\)

    \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*}

Parameters
[in]g格子点データ
Returns
スペクトルデータ

Definition at line 290 of file ae_module.f90.

References ae_ag().

290 
291  real(DP), dimension(0:im-1), intent(in) :: g
293  real(DP), dimension(-km:km) :: e_g
294 
295  real(DP), dimension(1,size(g)) :: ag_work
296  real(DP), dimension(1,-km:km) :: ae_work
297 
298  ag_work(1,:) = g
299  ae_work = ae_ag(ag_work)
300  e_g = ae_work(1,:)
301 
Here is the call graph for this function:

◆ g_e()

real(dp) function, dimension(0:im-1), public ae_module::g_e ( real(dp), dimension(-km:km), intent(in)  e)

スペクトルデータから格子点データへ逆変換する(1 次元データ用)

  • スペクトルデータ: \(e_k \, (k=-km,...,km)\)
  • 格子点データ: \(g_j \, (j=-0,...,im-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*}

Parameters
[in]eスペクトルデータ
Returns
格子点データ

Definition at line 217 of file ae_module.f90.

References ag_ae().

217 
218  real(DP), dimension(-km:km), intent(in) :: e
220  real(DP), dimension(0:im-1) :: g_e
221 
222  real(DP), dimension(1,size(e)) :: ae_work
223  real(DP), dimension(1,0:im-1) :: ag_work
224 
225  ae_work(1,:) = e
226  ag_work = ag_ae(ae_work)
227  g_e = ag_work(1,:)
228 
Here is the call graph for this function:

◆ int_g()

real(dp) function, public ae_module::int_g ( real(dp), dimension(0:im-1), intent(in)  g)

1 次元格子点データの積分

Parameters
[in]g格子点データ
Returns
積分結果

Definition at line 376 of file ae_module.f90.

References g_x_weight.

376 
377  real(DP), dimension(0:im-1), intent(in) :: g
379  real(DP) :: int_g
380 
381  int_g = sum(g*g_x_weight)
382 

◆ interpolate_e()

real(dp) function, public ae_module::interpolate_e ( real(dp), dimension(-km:km), intent(in)  e_data,
real(dp), intent(in)  xval 
)

補間計算(1次元)

Parameters
[in]e_data入力スペクトルデータ
[in]xval補間する点の座標
Returns
補間した結果の値

Definition at line 412 of file ae_module.f90.

412 
413  real(DP), dimension(-km:km), intent(IN) :: e_data
415  real(DP), intent(in) :: xval
417  real(DP) :: interpolate_e
418 
419  integer(8) :: k
420 
421  interpolate_e = e_data(0)
422  do k=1,km
423  interpolate_e = interpolate_e &
424  & + 2.0d0 * ( e_data(k) * cos(2*pi*k/xl*xval) &
425  & - e_data(-k) * sin(2*pi*k/xl*xval) )
426  end do
427 

Variable Documentation

◆ g_x

real(dp), dimension(:), allocatable, save, public ae_module::g_x

格子点座標(X)を格納した1次元配列.

Definition at line 102 of file ae_module.f90.

102  real(DP), allocatable :: g_x(:)

◆ g_x_weight

real(dp), dimension(:), allocatable, save, public ae_module::g_x_weight

重み座標を格納した 1 次元配列. X方向の格子点の間隔が格納してある.

Definition at line 104 of file ae_module.f90.

104  real(DP), allocatable :: g_x_weight(:)