Class | at_module |
In: |
libsrc/at_module/at_module.f90
|
Authors: | Shin-ichi Takehiro, Youhei SASAKI |
Version: | $Id: at_module.f90 598 2013-08-20 03:23:44Z takepiro $ |
Copyright&License: | See COPYRIGHT |
spml/at_module モジュールは 1 次元有限領域の下での流体運動を チェビシェフ変換によりスペクトル数値計算するための Fortran90 関数を 提供する.
2 次元データの 1 次元に関して同時にスペクトル計算を実行するための 関数も提供しており, 2, 3 次元領域での計算のベースも提供する.
このモジュールは内部で ISPACK/FTPACK の Fortran77 サブルーチンを 呼んでいる. スペクトルデータおよび格子点データの格納方法については ISPACK/FTPACK のマニュアルを参照されたい.
t_ : | チェビシェフデータ |
g_ : | 1 次元格子点データ |
at_ : | 1 次元チェビシェフデータが複数並んだ 2 次元データ |
ag_ : | 1 次元格子点データが複数並んだ 2 次元データ. |
_t : | チェビシェフデータ |
_g : | 1 次元格子点データ |
_at : | 1 次元チェビシェフデータが複数並んだ 2 次元データ |
_ag : | 1 次元格子点データが複数並んだ 2 次元データ |
at_Initial : | チェビシェフ変換の格子点数, 波数, 領域の大きさの設定 |
g_X : | 格子点座標(X)を格納した 1 次元配列 |
g_X_Weight : | 重み座標を格納した 1 次元配列 |
g_t, ag_at : | チェビシェフデータから格子データへの変換 |
t_g, at_ag : | 格子データからチェビシェフデータへの変換 |
t_Dx_t, at_Dx_at : | チェビシェフデータに X 微分を作用させる |
Interpolate_t, a_Interpolate_at :: チェビシェフデータから任意の点での値を求める
at_Boundaries_DD, at_Boundaries_DN, at_Boundaries_ND, at_Boundaries_NN :: ディリクレ,ノイマン境界条件の適用
at_BoundariesTau_DD, at_BoundariesTau_DN, at_BoundariesTau_ND, at_BoundariesTau_NN :: ディリクレ, ノイマン境界条件の適用(タウ法)
at_BoundariesGrid_DD, at_BoundariesGrid_DN, at_BoundariesGrid_ND, at_BoundariesGrid_NN :: ディリクレ, ノイマン境界条件の適用(選点法)
a_Int_ag, a_Avr_ag : | 1 次元格子点データの並んだ 2 次元配列の積分および平均 |
Int_g, Avr_g : | 1 次元格子点データの積分および平均 |
Function : | |||
Interpolate_t : | real(8)
| ||
t_data : | real(8), dimension(0:), intent(in)
| ||
xval : | real(8), intent(in)
|
function Interpolate_t(t_data,xval) real(8), dimension(0:), intent(in) :: t_data !(in) 入力チェビシェフデータ real(8), intent(in) :: xval ! 補間する点の座標 real(8) :: Interpolate_t ! 補間した結果の値 integer :: kmax ! 入力配列の最大次数 real(8) :: y2, y1, y0, x ! Crenshow's reccurence formula 計算用変数 integer :: k ! DO 文変数 kmax = size(t_data)-1 x =(xval -(g_X(0)+g_X(im))/2 )/(g_X(0)-g_X(im))*2 y2 = 0 y1 = 0 do k=kmax,1,-1 y0 = 2*x*y1 - y2 + t_data(k) y2 = y1 y1 = y0 enddo Interpolate_t = - y2 + x*y1 + t_data(0)/2 if ( kmax == im ) then Interpolate_t = Interpolate_t -t_data(kmax)/2*cos(kmax*acos(x)) endif end function Interpolate_t
Function : | |||
a_Avr_ag : | real(8), dimension(size(ag,1))
| ||
ag : | real(8), dimension(:,0:), intent(in)
|
1 次元格子点データが並んだ 2 次元配列の平均
function a_Avr_ag(ag) ! ! 1 次元格子点データが並んだ 2 次元配列の平均 ! real(8), dimension(:,0:), intent(in) :: ag !(in)入力格子点データ real(8), dimension(size(ag,1)) :: a_Avr_ag !(out) 平均したデータ a_Avr_ag = a_Int_ag(ag)/sum(g_X_Weight) end function a_Avr_ag
Function : | |||
a_Int_ag : | real(8), dimension(size(ag,1))
| ||
ag : | real(8), dimension(:,0:), intent(in)
|
1 次元格子点データが並んだ 2 次元配列の積分
function a_Int_ag(ag) ! ! 1 次元格子点データが並んだ 2 次元配列の積分 ! real(8), dimension(:,0:), intent(in) :: ag !(in)入力格子点データ real(8), dimension(size(ag,1)) :: a_Int_ag !(out) 積分したデータ integer :: i if ( size(ag,2) < im+1 ) then call MessageNotify('E','ae_ag', 'The Grid points of input data too small.') elseif ( size(ag,2) > im+1 ) then call MessageNotify('W','ae_ag', 'The Grid points of input data too large.') endif a_Int_ag = 0.0d0 do i=0,im a_Int_ag(:) = a_Int_ag(:) + ag(:,i)*g_X_Weight(i) enddo end function a_Int_ag
Function : | |||
a_Interpolate_at : | real(8), dimension(size(at_data,1))
| ||
at_data : | real(8), dimension(:,0:), intent(in)
| ||
xval : | real(8), intent(in)
|
function a_Interpolate_at(at_data,xval) real(8), dimension(:,0:), intent(in) :: at_data !(in) 入力チェビシェフデータ real(8), intent(in) :: xval ! 補間する点の座標 real(8), dimension(size(at_data,1)) :: a_Interpolate_at ! 補間した結果の値 integer :: kmax ! 入力配列の最大次数 real(8), dimension(size(at_data,1)) :: y2, y1, y0 real(8) :: x ! Crenshow's reccurence formula 計算用変数 integer :: k ! DO 文変数 kmax = size(at_data,2)-1 x =(xval -(g_X(0)+g_X(im))/2 )/(g_X(0)-g_X(im))*2 y2 = 0 y1 = 0 do k=kmax,1,-1 y0 = 2*x*y1 - y2 + at_data(:,k) y2 = y1 y1 = y0 enddo a_Interpolate_at = - y2 + x*y1 + at_data(:,0)/2 if ( kmax == im ) then a_Interpolate_at = a_Interpolate_at -at_data(:,kmax)/2*cos(kmax*acos(x)) endif end function a_Interpolate_at
Function : | |||
ag_at : | real(8), dimension(size(at_data,1),0:im)
| ||
at_data : | real(8), dimension(:,:), intent(in)
|
チェビシェフデータから格子データへ変換する(2 次元配列用).
function ag_at(at_data) ! ! チェビシェフデータから格子データへ変換する(2 次元配列用). ! real(8), dimension(:,:), intent(in) :: at_data !(in) チェビシェフデータ real(8), dimension(size(at_data,1),0:im) :: ag_at !(out) 格子点データ real(8), dimension(size(at_data,1)*im) :: y ! 作業用配列 integer :: m m = size(at_data,1) if ( size(at_data,2)-1 < km ) then call MessageNotify('E','ag_at', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','ag_at', 'The Chebyshev dimension of input data too large.') endif ag_at = 0.0D0 ag_at(:,0:km)=at_data call fttctb(m,im,ag_at,y,it,t) end function ag_at
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet 型境界条件の適用(実空間での評価, 2 次元配列用) 両境界での値を与える.
subroutine at_BoundariesGrid_DD_2d(at_data,values) ! ! Dirichlet 型境界条件の適用(実空間での評価, 2 次元配列用) ! 両境界での値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(size(at_data,1),0:im) :: ag_data real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( im /= km ) then call MessageNotify('E','at_BoundariesGrid_DD', 'Chebyshev truncation and number of grid points should be same.') endif if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesGrid_DD', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesGrid_DD', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:im,0:km),kp(0:im)) tt_data = 0 do k=0,km tt_data(k,k)=1.0 enddo tg_data = ag_at(tt_data) alu = transpose(tg_data) ! alu(km-1,:) = tg_data(:,0) ! alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif ag_data = ag_at(at_data) ag_data(:,0) = value1 ag_data(:,im) = value2 at_data = lusolve(alu,kp,ag_data) end subroutine at_BoundariesGrid_DD_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Dirichlet 型境界条件の適用(実空間での評価, 1 次元配列用) 両境界での値を与える.
subroutine at_BoundariesGrid_DD_1d(t_data,values) ! ! Dirichlet 型境界条件の適用(実空間での評価, 1 次元配列用) ! 両境界での値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesGrid_DD_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesGrid_DD_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(実空間での評価, 2 次元配列用) i=0 で値, i=im で勾配の値を与える.
subroutine at_BoundariesGrid_DN_2d(at_data,values) ! ! Dirichlet/Neumann 型境界条件の適用(実空間での評価, 2 次元配列用) ! i=0 で値, i=im で勾配の値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(size(at_data,1),0:im) :: ag_data real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( im /= km ) then call MessageNotify('E','at_BoundariesGrid_DN', 'Chebyshev truncation and number of grid points should be same.') endif if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesGrid_DN', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesGrid_DN', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:im,0:km),kp(0:im)) tt_data = 0 do k=0,km tt_data(k,k)=1.0 enddo tg_data = ag_at(tt_data) alu = transpose(tg_data) tg_data = ag_at(at_dx_at(tt_data)) alu(im,:) = tg_data(:,im) call ludecomp(alu,kp) endif ag_data = ag_at(at_data) ag_data(:,0) = value1 ag_data(:,im) = value2 at_data = lusolve(alu,kp,ag_data) end subroutine at_BoundariesGrid_DN_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(実空間での評価, 1 次元配列用) i=0 で値, i=im で勾配の値を与える.
subroutine at_BoundariesGrid_DN_1d(t_data,values) ! ! Dirichlet/Neumann 型境界条件の適用(実空間での評価, 1 次元配列用) ! i=0 で値, i=im で勾配の値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesGrid_DN_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesGrid_DN_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(実空間での評価, 2 次元配列用) i=0 で勾配の値, i=im で値を与える.
subroutine at_BoundariesGrid_ND_2d(at_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(実空間での評価, 2 次元配列用) ! i=0 で勾配の値, i=im で値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(size(at_data,1),0:im) :: ag_data real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( im /= km ) then call MessageNotify('E','at_BoundariesGrid_ND', 'Chebyshev truncation and number of grid points should be same.') endif if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesGrid_ND', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesGrid_DD', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:im,0:km),kp(0:im)) tt_data = 0 do k=0,km tt_data(k,k)=1.0 enddo tg_data = ag_at(tt_data) alu = transpose(tg_data) tg_data = ag_at(at_dx_at(tt_data)) alu(0,:) = tg_data(:,0) call ludecomp(alu,kp) endif ag_data = ag_at(at_data) ag_data(:,0) = value1 ag_data(:,im) = value2 at_data = lusolve(alu,kp,ag_data) end subroutine at_BoundariesGrid_ND_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Neumann 型境界条件の適用(実空間での評価, 1 次元配列用) i=0 で勾配の値, i=im で値を与える.
subroutine at_BoundariesGrid_ND_1d(t_data,values) ! ! Neumann 型境界条件の適用(実空間での評価, 1 次元配列用) ! i=0 で勾配の値, i=im で値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesGrid_ND_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesGrid_ND_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Neumann 型境界条件の適用(実空間での評価, 2 次元配列用) 両境界で勾配の値を与える.
subroutine at_BoundariesGrid_NN_2d(at_data,values) ! ! Neumann 型境界条件の適用(実空間での評価, 2 次元配列用) ! 両境界で勾配の値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(size(at_data,1),0:im) :: ag_data real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( im /= km ) then call MessageNotify('E','at_BoundariesGrid_NN', 'Chebyshev truncation and number of grid points should be same.') endif if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesGrid_NN', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesGrid_NN', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:im,0:km),kp(0:im)) tt_data = 0 do k=0,km tt_data(k,k)=1.0 enddo tg_data = ag_at(tt_data) alu = transpose(tg_data) tg_data = ag_at(at_dx_at(tt_data)) alu(0,:) = tg_data(:,0) alu(im,:) = tg_data(:,im) call ludecomp(alu,kp) endif ag_data = ag_at(at_data) ag_data(:,0) = value1 ag_data(:,im) = value2 at_data = lusolve(alu,kp,ag_data) end subroutine at_BoundariesGrid_NN_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Neumann 型境界条件の適用(実空間での評価, 1 次元配列用) 両境界で勾配の値を与える.
subroutine at_BoundariesGrid_NN_1d(t_data,values) ! ! Neumann 型境界条件の適用(実空間での評価, 1 次元配列用) ! 両境界で勾配の値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesGrid_NN_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesGrid_NN_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) 両境界での値を与える.
subroutine at_BoundariesTau_DD_2d(at_data,values) ! ! Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) ! 両境界での値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data ! データ(m,0:km) !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_DD', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_DD', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data=0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(tt_data) alu(km-1,:) = tg_data(:,0) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_DD_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) 両境界での値を与える.
subroutine at_BoundariesTau_DD_1d(t_data,values) ! ! Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) ! 両境界での値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_DD_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_DD_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(タウ法, 2 次元配列用) i=0 で値, i=im で勾配の値を与える.
subroutine at_BoundariesTau_DN_2d(at_data,values) ! ! Dirichlet/Neumann 型境界条件の適用(タウ法, 2 次元配列用) ! i=0 で値, i=im で勾配の値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_DN', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_DN', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data = 0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(tt_data) alu(km-1,:) = tg_data(:,0) tg_data = ag_at(at_Dx_at(tt_data)) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_DN_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(タウ法, 1 次元配列用) i=0 で値, i=im で勾配の値を与える.
subroutine at_BoundariesTau_DN_1d(t_data,values) ! ! Dirichlet/Neumann 型境界条件の適用(タウ法, 1 次元配列用) ! i=0 で値, i=im で勾配の値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_DN_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_DN_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) i=0 で勾配の値, i=im で値を与える.
subroutine at_BoundariesTau_ND_2d(at_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) ! i=0 で勾配の値, i=im で値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_ND', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_ND', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data = 0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(at_Dx_at(tt_data)) alu(km-1,:) = tg_data(:,0) tg_data = ag_at(tt_data) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_ND_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) i=0 で勾配の値, i=im で値を与える.
subroutine at_BoundariesTau_ND_1d(t_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) ! i=0 で勾配の値, i=im で値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_ND_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_ND_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) 両境界で勾配の値を与える.
subroutine at_BoundariesTau_NN_2d(at_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) ! 両境界で勾配の値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_NN', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_NN', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data = 0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(at_Dx_at(tt_data)) alu(km-1,:) = tg_data(:,0) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_NN_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) 両境界で勾配の値を与える.
このサブルーチンを直接使うことを勧めない. 共通インターフェース at_Boundaries_NN を用いること.
subroutine at_BoundariesTau_NN_1d(t_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) ! 両境界で勾配の値を与える. ! ! このサブルーチンを直接使うことを勧めない. ! 共通インターフェース at_Boundaries_NN を用いること. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_NN_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_NN_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) 両境界での値を与える.
subroutine at_BoundariesTau_DD_2d(at_data,values) ! ! Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) ! 両境界での値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data ! データ(m,0:km) !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_DD', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_DD', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data=0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(tt_data) alu(km-1,:) = tg_data(:,0) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_DD_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) 両境界での値を与える.
subroutine at_BoundariesTau_DD_1d(t_data,values) ! ! Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) ! 両境界での値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_DD_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_DD_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(タウ法, 2 次元配列用) i=0 で値, i=im で勾配の値を与える.
subroutine at_BoundariesTau_DN_2d(at_data,values) ! ! Dirichlet/Neumann 型境界条件の適用(タウ法, 2 次元配列用) ! i=0 で値, i=im で勾配の値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_DN', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_DN', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data = 0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(tt_data) alu(km-1,:) = tg_data(:,0) tg_data = ag_at(at_Dx_at(tt_data)) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_DN_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(タウ法, 1 次元配列用) i=0 で値, i=im で勾配の値を与える.
subroutine at_BoundariesTau_DN_1d(t_data,values) ! ! Dirichlet/Neumann 型境界条件の適用(タウ法, 1 次元配列用) ! i=0 で値, i=im で勾配の値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_DN_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_DN_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) i=0 で勾配の値, i=im で値を与える.
subroutine at_BoundariesTau_ND_2d(at_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) ! i=0 で勾配の値, i=im で値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_ND', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_ND', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data = 0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(at_Dx_at(tt_data)) alu(km-1,:) = tg_data(:,0) tg_data = ag_at(tt_data) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_ND_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) i=0 で勾配の値, i=im で値を与える.
subroutine at_BoundariesTau_ND_1d(t_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) ! i=0 で勾配の値, i=im で値を与える. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_ND_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_ND_1d
Subroutine : | |||
at_data : | real(8), dimension(:,0:),intent(inout)
| ||
values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) 両境界で勾配の値を与える.
subroutine at_BoundariesTau_NN_2d(at_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) ! 両境界で勾配の値を与える. ! real(8), dimension(:,0:),intent(inout) :: at_data !(inout) 境界条件を適用するチェビシェフデータ(m,0:km) real(8), dimension(:,:), intent(in), optional :: values !(in) 境界値(m,2) real(8), dimension(:,:), allocatable :: alu integer, dimension(:), allocatable :: kp real(8), dimension(0:km,0:km) :: tt_data real(8), dimension(0:km,0:im) :: tg_data real(8), dimension(size(at_data,1)) :: value1, value2 ! 境界値 logical :: first = .true. integer :: k save :: alu, kp, first if ( size(at_data,2)-1 < km ) then call MessageNotify('E','at_BoundariesTau_NN', 'The Chebyshev dimension of input data too small.') elseif ( size(at_data,2)-1 > km ) then call MessageNotify('W','at_BoundariesTau_NN', 'The Chebyshev dimension of input data too large.') endif if (.not. present(values)) then value1=0 value2=0 else value1 = values(:,1) value2 = values(:,2) endif if ( first ) then first = .false. allocate(alu(0:km,0:km),kp(0:km)) tt_data = 0 do k=0,km tt_data(k,k)=1 enddo alu = tt_data tg_data = ag_at(at_Dx_at(tt_data)) alu(km-1,:) = tg_data(:,0) alu(km,:) = tg_data(:,im) call ludecomp(alu,kp) endif at_data(:,km-1) = value1 at_data(:,km) = value2 at_data = lusolve(alu,kp,at_data) end subroutine at_BoundariesTau_NN_2d
Subroutine : | |||
t_data : | real(8), dimension(0:km),intent(inout)
| ||
values : | real(8), dimension(2), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) 両境界で勾配の値を与える.
このサブルーチンを直接使うことを勧めない. 共通インターフェース at_Boundaries_NN を用いること.
subroutine at_BoundariesTau_NN_1d(t_data,values) ! ! Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) ! 両境界で勾配の値を与える. ! ! このサブルーチンを直接使うことを勧めない. ! 共通インターフェース at_Boundaries_NN を用いること. ! real(8), dimension(0:km),intent(inout) :: t_data !(inout) 境界条件を適用するチェビシェフデータ(0:km) real(8), dimension(2), intent(in), optional :: values !(in) 境界値 real(8), dimension(1,0:km) :: at_work real(8), dimension(1,2) :: vwork ! 境界値 if (.not. present(values)) then vwork(1,1)=0 vwork(1,2)=0 else vwork(1,:) = values endif at_work(1,:)=t_data call at_BoundariesTau_NN_2d(at_work,vwork) t_data=at_work(1,:) end subroutine at_BoundariesTau_NN_1d
Function : | |||
at_Dx_at : | real(8), dimension(size(at_data,1),0:size(at_data,2)-1)
| ||
at_data : | real(8), dimension(:,0:), intent(in)
|
入力チェビシェフデータに X 微分を作用する(2 次元配列用).
チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.
function at_Dx_at(at_data) ! ! 入力チェビシェフデータに X 微分を作用する(2 次元配列用). ! ! チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を ! 作用させたデータのチェビシェフ変換のことである. ! ! real(8), dimension(:,0:), intent(in) :: at_data !(in) 入力チェビシェフデータ real(8), dimension(size(at_data,1),0:size(at_data,2)-1) :: at_Dx_at !(out) チェビシェフデータの X 微分 integer :: m, k integer :: nm, kmax nm=size(at_data,1) kmax=size(at_data,2)-1 if ( kmax < km ) then call MessageNotify('W','at_Dx_at', 'The Chebyshev dimension of input data too small.') elseif ( kmax > km ) then call MessageNotify('E','at_Dx_at', 'The Chebyshev dimension of input data too large.') endif if ( kmax == im ) then do m=1,nm at_Dx_at(m,kmax) = 0. at_Dx_at(m,kmax-1) = 2 * km * at_data(m,kmax) /2 enddo else do m=1,nm at_Dx_at(m,kmax) = 0. at_Dx_at(m,kmax-1) = 2 * km * at_data(m,kmax) ! スタートはグリッド対応最大波数未満. Factor 1/2 不要 enddo endif do k=kmax-2,0,-1 do m=1,nm at_Dx_at(m,k) = at_Dx_at(m,k+2) + 2*(k+1)*at_data(m,k+1) enddo enddo do k=0,kmax do m=1,nm at_Dx_at(m,k) = 2/xl * at_Dx_at(m,k) enddo enddo end function at_Dx_at
Subroutine : | |||
i : | integer,intent(in)
| ||
k : | integer,intent(in)
| ||
xmin : | real(8),intent(in)
| ||
xmax : | real(8),intent(in)
|
チェビシェフ変換の格子点数, 波数, 領域の大きさを設定する.
他の関数や変数を呼ぶ前に, 最初にこのサブルーチンを呼んで 初期設定をしなければならない.
subroutine at_Initial(i,k,xmin,xmax) ! ! チェビシェフ変換の格子点数, 波数, 領域の大きさを設定する. ! ! 他の関数や変数を呼ぶ前に, 最初にこのサブルーチンを呼んで ! 初期設定をしなければならない. ! integer,intent(in) :: i !(in) 格子点数 integer,intent(in) :: k !(in) 切断波数 real(8),intent(in) :: xmin, xmax !(in) 座標の範囲 integer :: ii,kk im=i km=k xl = xmax-xmin if ( im <= 0 .or. km <= 0 ) then call MessageNotify('E','at_initial', 'Number of grid points and waves should be positive') elseif ( mod(im,2) /= 0 ) then call MessageNotify('E','at_initial','Number of grid points should be even') elseif ( km > im ) then call MessageNotify('E','at_initial','KM shoud be less equal IM') endif if ( allocated(t) ) deallocate(t) allocate(t(3*im)) call fttcti(im,it,t) if ( allocated(g_X) ) deallocate(g_X) allocate(g_X(0:im)) do ii=0,im g_X(ii) = (xmax+xmin)/2 + xl/2 * cos(pi*ii/im) enddo if ( allocated(g_X_Weight) ) deallocate(g_X_Weight) allocate(g_X_Weight(0:im)) do ii=0,im g_X_Weight(ii) = 1.0 do kk=2,km,2 g_X_Weight(ii) = g_X_Weight(ii) + 2/(1D0-kk**2) * cos(kk*ii*pi/im) enddo if ( (km == im) .and. (mod(im,2)==0) ) then ! 最後の和は factor 1/2. g_X_Weight(ii) = g_X_Weight(ii) - 1/(1D0-km**2)* cos(km*ii*pi/im) endif g_X_Weight(ii) = 2D0/im * g_X_Weight(ii) * xl/2 enddo g_X_Weight(0) = g_X_Weight(0) / 2 g_X_Weight(im) = g_X_Weight(im) / 2 call MessageNotify('M','at_initial','at_module (2009/07/31) is initialized') end subroutine at_Initial
Function : | |||
at_ag : | real(8), dimension(size(ag_data,1),0:km)
| ||
ag_data : | real(8), dimension(:,:), intent(in)
|
格子データからチェビシェフデータへ変換する(2 次元配列用).
function at_ag(ag_data) ! ! 格子データからチェビシェフデータへ変換する(2 次元配列用). ! real(8), dimension(:,:), intent(in) :: ag_data !(in) 格子点データ real(8), dimension(size(ag_data,1),0:km) :: at_ag !(out) チェビシェフデータ real(8), dimension(size(ag_data,1)*im) :: y real(8), dimension(size(ag_data,1),0:im) :: ag_work integer :: m m = size(ag_data,1) if ( size(ag_data,2)-1 < im ) then call MessageNotify('E','at_ag', 'The Grid points of input data too small.') elseif ( size(ag_data,2)-1 > im ) then call MessageNotify('W','at_ag', 'The Grid points of input data too large.') endif ag_work = ag_data call fttctf(m,im,ag_work,y,it,t) at_ag = ag_work(:,0:km) end function at_ag
Function : | |||
g_t : | real(8), dimension(0:im)
| ||
t_data : | real(8), dimension(:), intent(in)
|
チェビシェフデータから格子データへ変換する(1 次元配列用).
function g_t(t_data) ! ! チェビシェフデータから格子データへ変換する(1 次元配列用). ! real(8), dimension(:), intent(in) :: t_data !(in) チェビシェフデータ real(8), dimension(0:im) :: g_t !(out) 格子点データ real(8), dimension(1,size(t_data)) :: t_work ! 作業用配列 real(8), dimension(1,0:im) :: g_work ! 作業用配列 t_work(1,:) = t_data g_work = ag_at(t_work) g_t = g_work(1,:) end function g_t
Function : | |||
t_Dx_t : | real(8), dimension(size(t_data))
| ||
t_data : | real(8), dimension(:), intent(in)
|
入力チェビシェフデータに X 微分を作用する(1 次元配列用).
チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.
function t_Dx_t(t_data) ! ! 入力チェビシェフデータに X 微分を作用する(1 次元配列用). ! ! チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を ! 作用させたデータのチェビシェフ変換のことである. ! ! real(8), dimension(:), intent(in) :: t_data !(in) 入力チェビシェフデータ real(8), dimension(size(t_data)) :: t_Dx_t !(out) チェビシェフデータの X 微分 real(8), dimension(1,size(t_data)) :: at_work ! 作業用配列 at_work(1,:) = t_data at_work = at_Dx_at(at_work) t_Dx_t = at_work(1,:) end function t_Dx_t
Function : | |||
t_g : | real(8), dimension(0:km)
| ||
g_data : | real(8), dimension(:), intent(in)
|
台形格子 -> スペクトル
格子データからチェビシェフデータへ変換する(1 次元配列用).
function t_g(g_data) ! 台形格子 -> スペクトル ! ! 格子データからチェビシェフデータへ変換する(1 次元配列用). ! real(8), dimension(:), intent(in) :: g_data !(in) 格子点データ real(8), dimension(0:km) :: t_g !(out) チェビシェフデータ real(8), dimension(1,size(g_data)) :: ag_work real(8), dimension(1,0:km) :: at_work ag_work(1,:) = g_data at_work = at_ag(ag_work) t_g = at_work(1,:) end function t_g