Data Types | Functions/Subroutines | Variables
at_module Module Reference

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

Detailed Description

1 次元有限領域: チェビシェフ変換

Author
Shin-ichi Takehiro, Youhei SASAKI

1 次元有限領域の下での流体運動をチェビシェフ変換によりスペクトル数値 計算するための Fortran90 関数を提供する.

2 次元データの 1 次元に関して同時にスペクトル計算を実行するための 関数も提供しており, 2, 3 次元領域での計算のベースも提供する.

このモジュールは内部で ISPACK3/FXPACK の Fortran77 サブルーチンを 呼んでいる. スペクトルデータおよび格子点データの格納方法については ISPACK3/FXPACK のマニュアルを参照されたい.

命名法

各データの種類の説明

Function/Subroutine Documentation

◆ a_avr_ag()

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

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

Parameters
[in]ag入力格子点データ
Returns
平均したデータ

Definition at line 517 of file at_module.f90.

References a_int_ag(), and g_x_weight.

517 
518  real(DP), dimension(:,0:), intent(in) :: ag
520  real(DP), dimension(size(ag,1)) :: a_avr_ag
521  a_avr_ag = a_int_ag(ag)/sum(g_x_weight)
522 
Here is the call graph for this function:

◆ a_int_ag()

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

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

Parameters
[in]ag入力格子点データ
Returns
積分したデータ

Definition at line 483 of file at_module.f90.

References g_x_weight.

483 
484  real(DP), dimension(:,0:), intent(in) :: ag
486  real(DP), dimension(size(ag,1)) :: a_int_ag
487 
488  integer(8) :: i
489 
490  if ( size(ag,2) < im+1 ) then
491  call messagenotify('E','ae_ag', &
492  'The Grid points of input data too small.')
493  elseif ( size(ag,2) > im+1 ) then
494  call messagenotify('W','ae_ag', &
495  'The Grid points of input data too large.')
496  endif
497 
498  a_int_ag = 0.0d0
499  do i=0,im
500  a_int_ag(:) = a_int_ag(:) + ag(:,i)*g_x_weight(i)
501  enddo

◆ a_interpolate_at()

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次元配列の補間

Parameters
[in]at_data入力チェビシェフデータ
[in]xval補間する点の座標
Returns
補間した結果の値

Definition at line 576 of file at_module.f90.

References g_x.

576 
577  real(DP), dimension(:,0:), intent(in) :: at_data
579  real(DP), intent(in) :: xval
581  real(DP), dimension(size(at_data,1)) :: a_interpolate_at
582 
583 
584  integer :: kmax
585  ! 入力配列の最大次数
586  real(DP), dimension(size(at_data,1)) :: y2, y1, y0
587  real(DP) :: x
588  ! Crenshow's reccurence formula 計算用変数
589  integer :: k
590  ! DO 文変数
591 
592  kmax = size(at_data,2)-1
593 
594  ! x =(xval -(g_X(0)+g_X(im))/2 )/(g_X(0)-g_X(im))*2
595  x =(xval -(g_x(0)+g_x(im))*0.5d0 )/(g_x(0)-g_x(im))*2.0d0
596 
597  y2 = 0.0d0
598  y1 = 0.0d0
599  do k=kmax,1,-1
600  y0 = 2.0d0*x*y1 - y2 + at_data(:,k)
601  y2 = y1
602  y1 = y0
603  enddo
604 
605  ! a_Interpolate_at = - y2 + x*y1 + at_data(:,0)/2
606  a_interpolate_at = - y2 + x*y1 + at_data(:,0)*0.5d0
607  if ( kmax == int(im) ) then
608  a_interpolate_at = a_interpolate_at &
609  & - at_data(:,kmax)/2.0d0*cos(kmax*acos(x))
610  endif

◆ ag_at()

real(dp) function, dimension(size(at_data,1),0:im), public at_module::ag_at ( real(dp), dimension(:,0:), intent(in)  at_data)

チェビシェフデータから格子データへ変換する(2 次元配列用).

Parameters
[in]at_dataチェビシェフデータ
Returns
格子点データ

Definition at line 328 of file at_module.f90.

328 
329  real(DP), dimension(:,0:), intent(in) :: at_data
331  real(DP), dimension(size(at_data,1),0:im) :: ag_at
332 
333  real(DP), dimension(size(at_data,1),0:2*im-1) :: y
334  integer(8) :: nm, k
335 
336  nm=size(at_data,1)
337 
338  y = 0.0d0
339  y(:,0) = at_data(:,0)
340  do k=1,min(km,im-1)
341  y(:,2*k) = at_data(:,k)
342  y(:,2*k+1) = 0.0d0
343  enddo
344  if ( km == im ) then
345  y(:,1) = at_data(:,km)
346  endif
347 
348  call fxrtba(nm,2*im,y,it,t)
349  ! ag_at = y(:,0:im)/2
350  ag_at = y(:,0:im)*0.50d0
351 

◆ at_ag()

real(dp) function, dimension(size(ag_data,1),0:km), public at_module::at_ag ( real(dp), dimension(:,0:), intent(in)  ag_data)

格子データからチェビシェフデータへ変換する(2 次元配列用).

Parameters
[in]ag_data格子点データ
Returns
チェビシェフデータ

Definition at line 374 of file at_module.f90.

374 
375  real(DP), dimension(:,0:), intent(in) :: ag_data
377  real(DP), dimension(size(ag_data,1),0:km) :: at_ag
378 
379  real(DP), dimension(size(ag_data,1),0:2*im-1) :: y
380  integer(8) :: nm, k
381 
382  nm=size(ag_data,1)
383 
384  y(:,0) = ag_data(:,0)
385  do k=1,im-1
386  y(:,2*k) = ag_data(:,k)
387  y(:,2*k+1) = 0.0d0
388  enddo
389  y(:,1) = ag_data(:,im)
390 
391  call fxrtba(nm,2*im,y,it,t)
392  at_ag = y(:,0:km)/im
393 

◆ at_boundaries()

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

Parameters
[in,out]at_data境界条件を適用するチェビシェフデータ(m,0:km)
[in]cond境界条件.
  • DD : 両端ディリクレ
  • DN : 上端ディリクレ, 下端ノイマン
  • ND : 上端ノイマン, 下端ディリクレ
  • NN : 両端ノイマン
[in]values境界値(m,2)

Definition at line 619 of file at_module.f90.

619 
620  real(DP), intent(INOUT) :: at_data(:,:)
626  character(len=2), intent(IN), optional :: cond
628  real(DP), dimension(:,:), intent(in), optional :: values
629 
630  real(DP) :: bcvalues(size(at_data),2)
631  character(len=2) :: bc
632 
633  if ( present(cond) ) then
634  bc = cond
635  else
636  bc = 'DD'
637  endif
638 
639  if ( present(values) ) then
640  bcvalues = values
641  else
642  bcvalues = 0.0d0
643  endif
644  select case(bc)
645  case ('DD')
646  call at_boundaries_dd(at_data,bcvalues)
647  case ('DN')
648  call at_boundaries_dn(at_data,bcvalues)
649  case ('ND')
650  call at_boundaries_nd(at_data,bcvalues)
651  case ('NN')
652  call at_boundaries_nn(at_data,bcvalues)
653  case default
654  call at_boundaries_dd(at_data,bcvalues)
655  end select
656 

◆ at_boundariesgrid()

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

Parameters
[in,out]at_data境界条件を適用するチェビシェフデータ(m,0:km)
[in]cond境界条件.
  • DD : 両端ディリクレ
  • DN : 上端ディリクレ, 下端ノイマン
  • ND : 上端ノイマン, 下端ディリクレ
  • NN : 両端ノイマン
[in]values境界値(m,2)

Definition at line 711 of file at_module.f90.

711 
712  real(DP), intent(INOUT) :: at_data(:,:)
718  character(len=2), intent(IN), optional :: cond
720  real(DP), dimension(:,:), intent(in), optional :: values
721 
722  real(DP) :: bcvalues(size(at_data),2)
723  character(len=2) :: bc
724 
725  if ( present(cond) ) then
726  bc = cond
727  else
728  bc = 'DD'
729  endif
730 
731  if ( present(values) ) then
732  bcvalues = values
733  else
734  bcvalues = 0.0d0
735  endif
736 
737  select case(bc)
738  case ('DD')
739  call at_boundariesgrid_dd(at_data,bcvalues)
740  case ('DN')
741  call at_boundariesgrid_dn(at_data,bcvalues)
742  case ('ND')
743  call at_boundariesgrid_nd(at_data,bcvalues)
744  case ('NN')
745  call at_boundariesgrid_nn(at_data,bcvalues)
746  end select
747 

◆ at_dx_at()

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 微分を 作用させたデータのチェビシェフ変換のことである.

Parameters
[in]at_data入力チェビシェフデータ
Returns
チェビシェフデータの X 微分

Definition at line 420 of file at_module.f90.

420 
421  real(DP), dimension(:,0:), intent(in) :: at_data
423  real(DP), dimension(size(at_data,1),0:size(at_data,2)-1) :: at_dx_at
424 
425  integer(8) :: m, k
426  integer(8) :: nm, kmax
427 
428  nm=size(at_data,1)
429  kmax=size(at_data,2)-1
430 
431  ! TODO: loop 併合できないか?
432  if ( kmax == im ) then
433  do m=1,nm
434  at_dx_at(m,kmax) = 0.0d0
435  ! at_Dx_at(m,kmax-1) = 2.0D0 * km * at_data(m,kmax) /2
436  at_dx_at(m,kmax-1) = 2.0d0 * km * at_data(m,kmax)* 0.5d0
437  enddo
438  else
439  do m=1,nm
440  at_dx_at(m,kmax) = 0.0d0
441  ! スタートはグリッド対応最大波数未満. Factor 1/2 不要
442  at_dx_at(m,kmax-1) = 2.0d0 * km * at_data(m,kmax)
443  enddo
444  endif
445 
446  do k=kmax-2,0,-1
447  do m=1,nm
448  at_dx_at(m,k) = at_dx_at(m,k+2) + 2.0d0*(k+1)*at_data(m,k+1)
449  enddo
450  enddo
451 
452  do k=0,kmax
453  do m=1,nm
454  at_dx_at(m,k) = 2.0d0/xl * at_dx_at(m,k)
455  enddo
456  enddo
457 

◆ at_finalize()

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.

796  if ( .not. at_initialize ) then
797  call messagenotify('W','at_Finalize',&
798  'at_module not initialized yet')
799  return
800  endif
801 
802  deallocate(t) ! 変換用配列
803  deallocate(it) ! 変換用配列
804  deallocate(g_x)
805  deallocate(g_x_weight) ! 格子点座標格納配列
806 
807  at_initialize = .false.
808 
809  at_boundariestau_dd_first = .true.
810  at_boundariestau_dn_first = .true.
811  at_boundariestau_nd_first = .true.
812  at_boundariestau_nn_first = .true.
813  at_boundariesgrid_dd_first = .true.
814  at_boundariesgrid_dn_first = .true.
815  at_boundariesgrid_nd_first = .true.
816  at_boundariesgrid_nn_first = .true.
817 
818  call messagenotify('M','at_Finalize',&
819  & 'at_module (2019/09/29) is finalized')
820 
Here is the call graph for this function:

◆ at_initial()

subroutine, public at_module::at_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]xmin座標の範囲
[in]xmax座標の範囲

Definition at line 268 of file at_module.f90.

References g_x, and g_x_weight.

268 
269  integer,intent(in) :: i
271  integer,intent(in) :: k
273  real(DP),intent(in) :: xmin, xmax
274  integer(8) :: ii,kk
275 
276  im=i
277  km=k
278  xl = xmax-xmin
279 
280  if ( im <= 0 .or. km <= 0 ) then
281  call messagenotify('E','at_initial', &
282  & 'Number of grid points and waves should be positive')
283  elseif ( mod(im,2_8) /= 0 ) then
284  call messagenotify('E','at_initial', &
285  & 'Number of grid points should be even')
286  elseif ( km > im ) then
287  call messagenotify('E','at_initial', &
288  & 'KM shoud be less equal IM')
289  endif
290 
291  allocate(it(im),t(3*im))
292  call fxrini(2*im,it,t)
293 
294  allocate(g_x(0:im))
295  do ii=0,im
296  g_x(ii) = (xmax+xmin)/2 + xl/2 * cos(pi*ii/im)
297  end do
298 
299  allocate(g_x_weight(0:im))
300 
301  do ii=0,im
302  g_x_weight(ii) = 1.0d0
303  do kk=2,km,2
304  g_x_weight(ii) = g_x_weight(ii) &
305  + 2.0d0/(1.0d0 - kk**2) * cos(kk*ii*pi/im)
306  enddo
307  ! 最後の和は factor 1/2.
308  if ( (km == im) .and. (mod(im,2_8) == 0) ) then
309  g_x_weight(ii) = g_x_weight(ii) &
310  - 1.0d0/(1.0d0 - km**2)* cos(km*ii*pi/im)
311  endif
312  g_x_weight(ii) = 2.0d0/im * g_x_weight(ii) * (xl/2.0d0)
313  enddo
314  ! g_X_Weight(0) = g_X_Weight(0) / 2
315  ! g_X_Weight(im) = g_X_Weight(im) / 2
316  g_x_weight(0) = g_x_weight(0) * 0.5d0
317  g_x_weight(im) = g_x_weight(im) * 0.5d0
318 
319  at_initialize = .true.
320 
321  call messagenotify('M','at_initial',&
322  & 'at_module (2019/09/28) is initialized')

◆ avr_g()

real(dp) function, public at_module::avr_g ( real(dp), dimension(0:im), intent(in)  g)

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

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

Definition at line 528 of file at_module.f90.

References g_x_weight, and int_g().

528 
529  real(DP), dimension(0:im), intent(in) :: g
531  real(DP) :: avr_g
532 
533  avr_g = int_g(g)/sum(g_x_weight)
Here is the call graph for this function:

◆ g_t()

real(dp) function, dimension(0:im), public at_module::g_t ( real(dp), dimension(:), intent(in)  t_data)

チェビシェフデータから格子データへ変換する(1 次元配列用).

Parameters
[in]t_dataチェビシェフデータ
Returns
格子点データ

Definition at line 357 of file at_module.f90.

References ag_at().

357 
358  real(DP), dimension(:), intent(in) :: t_data
360  real(DP), dimension(0:im) :: g_t
361 
362  real(DP), dimension(1,size(t_data)) :: t_work
363  real(DP), dimension(1,0:im) :: g_work
364 
365  t_work(1,:) = t_data
366  g_work = ag_at(t_work)
367  g_t = g_work(1,:)
368 
Here is the call graph for this function:

◆ int_g()

real(dp) function, public at_module::int_g ( real(dp), dimension(0:im), intent(in)  g)

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

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

Definition at line 507 of file at_module.f90.

References g_x_weight.

507 
508  real(DP), dimension(0:im), intent(in) :: g
510  real(DP) :: int_g
511  int_g = sum(g*g_x_weight)

◆ interpolate_t()

real(dp) function, public at_module::interpolate_t ( real(dp), dimension(0:), intent(in)  t_data,
real(dp), intent(in)  xval 
)

1 次元データの補間

Parameters
[in]t_data入力チェビシェフデータ
[in]xval補間する点の座標
Returns
補間した結果の値

Definition at line 539 of file at_module.f90.

References g_x.

539 
540  real(DP), dimension(0:), intent(in) :: t_data
542  real(DP), intent(in) :: xval
544  real(DP) :: interpolate_t
545 
546  integer :: kmax
547  ! Crenshow's reccurence formula 計算用変数
548  real(DP) :: y2, y1, y0, x
549  integer :: k
550  ! DO 文変数
551 
552  kmax = size(t_data)-1
553 
554  ! x =(xval -(g_X(0)+g_X(im))/2 )/(g_X(0)-g_X(im))*2
555  x =(xval -(g_x(0)+g_x(im))*0.50d0 )/(g_x(0)-g_x(im))*2.0d0
556 
557  y2 = 0.0d0
558  y1 = 0.0d0
559  do k=kmax,1,-1
560  y0 = 2.0d0*x*y1 - y2 + t_data(k)
561  y2 = y1
562  y1 = y0
563  enddo
564 
565  ! Interpolate_t = - y2 + x*y1 + t_data(0)/2
566  interpolate_t = - y2 + x*y1 + t_data(0)*0.5d0
567  if ( kmax == int(im) ) then
568  interpolate_t = interpolate_t -t_data(kmax)/2.0d0*cos(kmax*acos(x))
569  endif
570 

◆ t_boundaries()

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

Parameters
[in,out]t_data境界条件を適用するチェビシェフデータ(0:km)
[in]cond境界条件.
  • DD : 両端ディリクレ
  • DN : 上端ディリクレ, 下端ノイマン
  • ND : 上端ノイマン, 下端ディリクレ
  • NN : 両端ノイマン
[in]values境界値(2)

Definition at line 664 of file at_module.f90.

664 
665  real(DP), intent(INOUT) :: t_data(:)
671  character(len=2), intent(IN), optional :: cond
673  real(DP), dimension(2), intent(in), optional :: values
674 
675  real(DP) :: bcvalues(2)
676  character(len=2) :: bc
677 
678 
679  if ( present(cond) ) then
680  bc = cond
681  else
682  bc = 'DD'
683  endif
684 
685  if ( present(values) ) then
686  bcvalues = values
687  else
688  bcvalues = 0.0d0
689  endif
690 
691  select case(bc)
692  case ('DD' )
693  call at_boundaries_dd(t_data,bcvalues)
694  case ('DN' )
695  call at_boundaries_dn(t_data,bcvalues)
696  case ('ND' )
697  call at_boundaries_nd(t_data,bcvalues)
698  case ('NN' )
699  call at_boundaries_nn(t_data,bcvalues)
700  case default
701  call at_boundaries_dd(t_data,bcvalues)
702  end select
703 

◆ t_boundariesgrid()

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

Parameters
[in,out]t_data境界条件を適用するチェビシェフデータ(0:km)
[in]cond境界条件.
  • DD : 両端ディリクレ
  • DN : 上端ディリクレ, 下端ノイマン
  • ND : 上端ノイマン, 下端ディリクレ
  • NN : 両端ノイマン
[in]values境界値(2)

Definition at line 755 of file at_module.f90.

755 
756  real(DP), intent(INOUT) :: t_data(:)
762  character(len=2), intent(IN), optional :: cond
764  real(DP), dimension(2), intent(in), optional :: values
765 
766  real(DP) :: bcvalues(2)
767  character(len=2) :: bc
768 
769  if ( present(cond) ) then
770  bc = cond
771  else
772  bc = 'DD'
773  endif
774 
775  if ( present(values) ) then
776  bcvalues = values
777  else
778  bcvalues = 0.0d0
779  endif
780  select case(bc)
781  case ('DD')
782  call at_boundariesgrid_dd(t_data,bcvalues)
783  case ('DN')
784  call at_boundariesgrid_dn(t_data,bcvalues)
785  case ('ND')
786  call at_boundariesgrid_nd(t_data,bcvalues)
787  case ('NN')
788  call at_boundariesgrid_nn(t_data,bcvalues)
789  end select
790 

◆ t_dx_t()

real(dp) function, dimension(size(t_data)), public at_module::t_dx_t ( real(dp), dimension(:), intent(in)  t_data)

入力チェビシェフデータに X 微分を作用する(1 次元配列用).

チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.

Parameters
[in]t_data入力チェビシェフデータ
Returns
チェビシェフデータの X 微分

Definition at line 467 of file at_module.f90.

References at_dx_at().

467 
468  real(DP), dimension(:), intent(in) :: t_data
470  real(DP), dimension(size(t_data)) :: t_dx_t
471 
472  real(DP), dimension(1,size(t_data)) :: at_work
473 
474  at_work(1,:) = t_data
475  at_work = at_dx_at(at_work)
476  t_dx_t = at_work(1,:)
477 
Here is the call graph for this function:

◆ t_g()

real(dp) function, dimension(0:km), public at_module::t_g ( real(dp), dimension(:), intent(in)  g_data)

格子データからチェビシェフデータへ変換する(1 次元配列用).

Parameters
[in]g_data格子点データ
Returns
チェビシェフデータ

Definition at line 399 of file at_module.f90.

References at_ag().

399 
400  real(DP), dimension(:), intent(in) :: g_data
402  real(DP), dimension(0:km) :: t_g
403 
404 
405  real(DP), dimension(1,size(g_data)) :: ag_work
406  real(DP), dimension(1,0:km) :: at_work
407 
408  ag_work(1,:) = g_data
409  at_work = at_ag(ag_work)
410  t_g = at_work(1,:)
411 
Here is the call graph for this function:

Variable Documentation

◆ g_x

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

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

Definition at line 218 of file at_module.f90.

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

◆ g_x_weight

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

重み座標を格納した 1 次元配列

Definition at line 220 of file at_module.f90.

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