| Class | ffts |
| In: |
ffts.f90
|
fft 関連のサブルーチン集
| Subroutine : | |||
| nx : | integer, intent(in)
| ||
| a : | complex, intent(in), dimension(0:nx/2-1)
| ||
| b : | real, intent(inout), dimension(0:nx-1)
| ||
| prim : | character(1), optional, intent(in)
|
1 次元実 FFT 逆変換計算ルーチン
入力配列は複素数数で行い,実数配列を返す.
subroutine c2r_ffttp_1d( nx, a, b, prim )
! 1 次元実 FFT 逆変換計算ルーチン
!
! 入力配列は複素数数で行い,実数配列を返す.
implicit none
integer, intent(in) :: nx ! 1 次元データ要素数
complex, intent(in), dimension(0:nx/2-1) :: a ! 入力複素数データ配列
real, intent(inout), dimension(0:nx-1) :: b ! 出力実数データ配列
character(1), optional, intent(in) :: prim ! 素因数分解をするかどうか
! [o=分解する, x=分解しない] default=分解しない. その場合は, 通常の DFT.
! 素因数分解する場合, nx=2^a*3^b*5^c*7^d までしか分解しないようにする.
! なお, 実 FFT 変換は現在, データ要素数が偶数の変換しか行わないように
! 実装されていることに注意 (奇数データで変換しようとするとエラーを返す).
! また, in 属性の引数は, データ数が半分になっていることに注意.
! 仕様として, 入力引数の n=0 虚部に nx/2 番目の実数データが入っているように
! データを渡す.
! ただし, 'x' のときでも, prim_fact が設定されていれば, そのべきで FFT する.
integer, dimension(4) :: prim_fact ! prim = x のとき設定すると
! そのべきで FFT を行う. prim_fact=(/a,b,c,d/) : 2^a*3^b*5^c*7^d
complex, dimension(0:nx/2-1) :: c, d
complex, parameter :: img=(0.0,1.0)
real, parameter :: pi=3.14159265
integer :: i, j, k
if(mod(nx,2)/=0)then
write(*,*) "*** ERROR ***"
write(*,*) "nx must be even number. Stop."
stop
end if
!-- nx/2 個の複素数データ(独立なフーリエ係数)を同数の複素数データに置き換える.
c(0)=a(0)
do i=1,nx/2-1
c(i)=(a(i)+conjg(a(nx/2-i))) +(img*cos(2.0*pi*i/real(nx))-sin(2.0*pi*i/real(nx))) *(a(i)-conjg(a(nx/2-i)))
end do
!-- FFT 開始
if(present(prim))then
call ffttp_1d( nx/2, c, d, 'i', prim )
else
call ffttp_1d( nx/2, c, d, 'i' )
end if
!-- 変換後の配列を整理
do i=0,nx/2-1
b(2*i)=real(d(i))
b(2*i+1)=aimag(d(i))
end do
end subroutine
| Subroutine : | |||
| nx : | integer, intent(in)
| ||
| a : | complex, intent(in), dimension(0:nx-1)
| ||
| b : | complex, intent(inout), dimension(0:nx-1)
| ||
| csign : | character(1), intent(in)
| ||
| prim : | character(1), optional, intent(in)
| ||
| prim_fact : | integer, dimension(4), optional
|
Temperton's FFT 1d の fft を csign をもとに, 正変換, 逆変換するルーチン
subroutine ffttp_1d( nx, a, b, csign, prim, prim_fact )
! Temperton's FFT
!
! 1d の fft を csign をもとに, 正変換, 逆変換するルーチン
implicit none
integer, intent(in) :: nx ! 入力配列の要素数
complex, intent(in), dimension(0:nx-1) :: a ! 入力配列
complex, intent(inout), dimension(0:nx-1) :: b ! 出力配列
character(1), intent(in) :: csign ! 正逆変換判定 [r=正変換, i=逆変換]
character(1), optional, intent(in) :: prim ! 素因数分解をするかどうか
! [o=分解する, x=分解しない] default=分解しない. その場合は, 通常の DFT.
! 素因数分解する場合, nx=2^a*3^b*5^c*7^d までしか分解しないようにする.
! ただし, 'x' のときでも, prim_fact が設定されていれば, そのべきで FFT する.
integer, dimension(4), optional :: prim_fact ! prim = x のとき設定すると
! そのべきで FFT を行う. prim_fact=(/a,b,c,d/) : 2^a*3^b*5^c*7^d
integer, allocatable, dimension(:) :: l, m, n ! 要素等の作業用配列
real, parameter :: pi=3.14159265
complex, parameter :: img=(0.0,1.0)
complex :: fact, ctmp
integer :: stat, counter, base
integer :: i, j, k, id, jd, kd ! do loop 用配列
integer, parameter, dimension(4) :: prim_dim=(/2, 3, 5, 7/) ! 素因数
integer, dimension(4) :: prim_num ! 各素因数のべき数
complex, dimension(0:1,0:1) :: omega2
complex, dimension(0:2,0:2) :: omega3
complex, dimension(0:4,0:4) :: omega5
complex, dimension(0:6,0:6) :: omega7
complex, allocatable, dimension(:,:) :: omega
complex, dimension(0:nx-1,0:nx-1) :: omegan
complex, dimension(0:nx-1) :: c ! tmp array
base=nx
prim_num=0
counter=0
do i=0,nx-1
b(i)=a(i)
end do
!-- 素因数分解する処理
if(present(prim))then
if(prim=='o')then
do i=1,4
do while(mod(base,prim_dim(i))==0)
base=base/prim_dim(i)
prim_num(i)=prim_num(i)+1
counter=counter+1
end do
end do
if(base==1)then
counter=counter-1
end if
if(counter/=0)then ! prim=='o' であっても, 素因数分解できなければ DFT に送る.
allocate(l(counter+1))
allocate(m(counter+1))
allocate(n(counter+1))
stat=0
do i=1,4
if(prim_num(i)/=0)then
select case(prim_dim(i))
case(2)
n(stat+1:stat+prim_num(i))=2
stat=stat+prim_num(i)
case(3)
n(stat+1:stat+prim_num(i))=3
stat=stat+prim_num(i)
case(5)
n(stat+1:stat+prim_num(i))=5
stat=stat+prim_num(i)
case(7)
n(stat+1:stat+prim_num(i))=7
stat=stat+prim_num(i)
end select
end if
end do
if(base/=1)then
n(counter+1)=base
end if
end if
else
do i=1,4
prim_num(i)=prim_fact(i)
end do
end if
end if
!-- 回転行列を定義
select case(csign)
case('r')
fact=cos(2.0*pi/real(nx))-img*sin(2.0*pi/real(nx))
if(counter/=0)then
omega2=1.0
omega2(1,1)=-1.0
do j=0,2
do i=0,2
omega3(i,j)=cos(2.0*pi*i*j/3.0)-img*sin(2.0*pi*i*j/3.0)
end do
end do
do j=0,4
do i=0,4
omega5(i,j)=cos(2.0*pi*i*j/5.0)-img*sin(2.0*pi*i*j/5.0)
end do
end do
do j=0,6
do i=0,6
omega7(i,j)=cos(2.0*pi*i*j/7.0)-img*sin(2.0*pi*i*j/7.0)
end do
end do
if(base/=1)then
allocate(omega(0:base-1,0:base-1))
do j=0,base-1
do i=0,base-1
omega(i,j)=cos(2.0*pi*i*j/real(base))-img*sin(2.0*pi*i*j/real(base))
end do
end do
end if
else
allocate(omega(0:nx-1,0:nx-1))
do j=0,nx-1
do i=0,nx-1
omega(i,j)=cos(2.0*pi*i*j/real(nx))-img*sin(2.0*pi*i*j/real(nx))
end do
end do
end if
do j=0,nx-1
do i=0,nx-1
omegan(i,j)=cos(2.0*pi*i*j/real(nx))-img*sin(2.0*pi*i*j/real(nx))
end do
end do
case('i')
if(counter/=0)then
fact=exp(2.0*img*pi/real(nx))
fact=cos(2.0*pi/real(nx))+img*sin(2.0*pi/real(nx))
omega2=1.0
omega2(1,1)=-1.0
do j=0,2
do i=0,2
omega3(i,j)=cos(2.0*pi*i*j/3.0)+img*sin(2.0*pi*i*j/3.0)
end do
end do
do j=0,4
do i=0,4
omega5(i,j)=cos(2.0*pi*i*j/5.0)+img*sin(2.0*pi*i*j/5.0)
end do
end do
do j=0,6
do i=0,6
omega7(i,j)=cos(2.0*pi*i*j/7.0)+img*sin(2.0*pi*i*j/7.0)
end do
end do
if(base/=1)then
allocate(omega(0:base-1,0:base-1))
do j=0,base-1
do i=0,base-1
omega(i,j)=cos(2.0*pi*i*j/real(base))+img*sin(2.0*pi*i*j/real(base))
end do
end do
end if
end if
do j=0,nx-1
do i=0,nx-1
omegan(i,j)=cos(2.0*pi*i*j/real(nx))+img*sin(2.0*pi*i*j/real(nx))
end do
end do
case default
write(*,*) "******** ERROR : csign is bad. **********"
write(*,*) "Stop!"
stop
end select
!-- FFT 開始
if(counter/=0)then
!-- 係数行列定義
m(1)=1
l(1)=nx/(n(1)*m(1))
do i=2,counter+1
m(i)=m(i-1)*n(i-1)
l(i)=nx/(n(i)*m(i))
end do
!-- 変換行列 W の定義
do kd=1,counter+1
do jd=0,l(kd)-1
do id=0,n(kd)-1
do k=0,m(kd)-1
ctmp=b(jd*m(kd)+k)
do j=1,n(kd)-1
select case(n(kd))
case(2)
ctmp=ctmp+omega2(id,j)*b(j*l(kd)*m(kd)+jd*m(kd)+k)
case(3)
ctmp=ctmp+omega3(id,j)*b(j*l(kd)*m(kd)+jd*m(kd)+k)
case(5)
ctmp=ctmp+omega5(id,j)*b(j*l(kd)*m(kd)+jd*m(kd)+k)
case(7)
ctmp=ctmp+omega7(id,j)*b(j*l(kd)*m(kd)+jd*m(kd)+k)
case default
ctmp=ctmp+omega(id,j)*b(j*l(kd)*m(kd)+jd*m(kd)+k)
end select
end do
c(jd*n(kd)*m(kd)+id*m(kd)+k)=ctmp*omegan(m(kd),(id*jd))
end do
end do
end do
do id=0,nx-1
b(id)=c(id)
end do
end do
else
do j=0,nx-1
b(j)=a(0)
do i=1,nx-1
b(j)=b(j)+a(i)*omegan(i,j)
end do
end do
end if
if(csign=='r')then
do j=0,nx-1
b(j)=b(j)/real(nx)
end do
end if
end subroutine
| Subroutine : | |||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| a : | complex, intent(in), dimension(0:nx-1,0:ny-1)
| ||
| b : | complex, intent(inout), dimension(0:nx-1,0:ny-1)
| ||
| csign : | character(1), intent(in)
| ||
| prim : | character(1), optional, intent(in)
| ||
| prim_fact : | integer, dimension(4), optional
|
Temperton's FFT (2d ver) 2d の fft を csign をもとに, 正変換, 逆変換するルーチン
subroutine ffttp_2d( nx, ny, a, b, csign, prim, prim_fact )
! Temperton's FFT (2d ver)
!
! 2d の fft を csign をもとに, 正変換, 逆変換するルーチン
implicit none
integer, intent(in) :: nx ! 入力配列の要素数 1
integer, intent(in) :: ny ! 入力配列の要素数 2
complex, intent(in), dimension(0:nx-1,0:ny-1) :: a ! 入力配列
complex, intent(inout), dimension(0:nx-1,0:ny-1) :: b ! 出力配列
character(1), intent(in) :: csign ! 正逆変換判定 [r=正変換, i=逆変換]
character(1), optional, intent(in) :: prim ! 素因数分解をするかどうか
! [o=分解する, x=分解しない] default=分解しない. その場合は, 通常の DFT.
! 素因数分解する場合, nx=2^a*3^b*5^c*7^d までしか分解しないようにする.
! ただし, 'x' のときでも, prim_fact が設定されていれば, そのべきで FFT する.
integer, dimension(4), optional :: prim_fact ! prim = x のとき設定すると
! そのべきで FFT を行う. prim_fact=(/a,b,c,d/) : 2^a*3^b*5^c*7^d
! ここで, prim_fact で設定するべき数は nx 方向のべき数であることに注意.
integer :: i
if(present(prim))then
if(present(prim_fact))then
do i=0,ny-1
call ffttp_1d( nx, a(0:nx-1,i), b(0:nx-1,i), csign, prim, prim_fact )
end do
else
do i=0,ny-1
call ffttp_1d( nx, a(0:nx-1,i), b(0:nx-1,i), csign, prim )
end do
end if
else
do i=0,ny-1
call ffttp_1d( nx, a(0:nx-1,i), b(0:nx-1,i), csign )
end do
end if
end subroutine
| Subroutine : | |||
| nx : | integer, intent(in)
| ||
| a : | real, intent(in), dimension(0:nx-1)
| ||
| b : | complex, intent(inout), dimension(0:nx/2-1)
| ||
| prim : | character(1), optional, intent(in)
|
1 次元実 FFT 計算ルーチン
入力配列は実数で行い, 複素数配列を返す.
subroutine r2c_ffttp_1d( nx, a, b, prim )
! 1 次元実 FFT 計算ルーチン
!
! 入力配列は実数で行い, 複素数配列を返す.
implicit none
integer, intent(in) :: nx ! 1 次元データ要素数
real, intent(in), dimension(0:nx-1) :: a ! 入力実数データ配列
complex, intent(inout), dimension(0:nx/2-1) :: b ! 出力複素数データ配列
character(1), optional, intent(in) :: prim ! 素因数分解をするかどうか
! [o=分解する, x=分解しない] default=分解しない. その場合は, 通常の DFT.
! 素因数分解する場合, nx=2^a*3^b*5^c*7^d までしか分解しないようにする.
! なお, 実 FFT 変換は現在, データ要素数が偶数の変換しか行わないように
! 実装されていることに注意 (奇数データで変換しようとするとエラーを返す).
! また, inout 属性の引数は, データ数が半分になっていることに注意.
! 仕様として, nx/2 番目の実数データは, n=0 番目の虚部に格納されて
! 出力配列に渡されることに注意.
! ただし, 'x' のときでも, prim_fact が設定されていれば, そのべきで FFT する.
integer, dimension(4), optional :: prim_fact ! prim = x のとき設定すると
! そのべきで FFT を行う. prim_fact=(/a,b,c,d/) : 2^a*3^b*5^c*7^d
! ここで入力するべき数は nx/2 のときのべき数である.
complex, dimension(0:nx/2-1) :: c, d
complex, parameter :: img=(0.0,1.0)
real, parameter :: pi=3.14159265
integer :: i, j, k
if(mod(nx,2)/=0)then
write(*,*) "*** ERROR ***"
write(*,*) "nx must be even number. Stop."
stop
end if
!-- nx 個の実数データを nx/2 個の複素数データに置き換える.
do i=0,nx/2-1
c(i)=a(2*i)+img*a(2*i+1)
end do
!-- FFT 開始
if(present(prim))then
call ffttp_1d( nx/2, c, d, 'r', prim )
else
call ffttp_1d( nx/2, c, d, 'r' )
end if
!-- 変換後の配列を整理
!-- b(N) の実部は b(0) の虚部に組み込むことにする.
!-- b(k) の計算で 0.25 をかけるのは, 上の fft で 2/N で規格化しており,
!-- もとの計算では, 1/2N で規格化しなければならないので,
!-- 1/4 をかけることで, 2/N -> 1/2N で規格化したことになる.
!-- b(0) にかかっている係数もその類.
b(0)=0.5*((real(d(0))+aimag(d(0)))+img*(real(d(0))-aimag(d(0))))
do i=1,nx/2-1
b(i)=0.25*((conjg(d(nx/2-i))+d(i)) -(sin(2.0*pi*real(i)/real(nx))+img*cos(2.0*pi*real(i)/real(nx))) *(d(i)-conjg(d(nx/2-i))))
end do
end subroutine