Class ae_module
In: src/ae_module.f90

Methods

a_Int_ag   ae_Dx_ae   ae_initial   ag_ae  

Included Modules

dc_message

Public Instance methods

a_Int_ag :real(8), dimension(size(ag,1))
ag :real(8), dimension(:,0:), intent(in)

————— 積分計算 ——————

[Source]

    function a_Int_ag(ag)
      real(8), dimension(:,0:), intent(in)     :: ag
      real(8), dimension(size(ag,1))           :: a_Int_ag
      integer :: i

      if ( size(ag,2) < im ) then
         call MessageNotify('E','ae_Int_ag',               'The Grid points of input data too small.')
      elseif ( size(ag,2) > im ) then
         call MessageNotify('W','ae_Int_ag',               'The Grid points of input data too large.')
      endif

      a_Int_ag = 0.0d0
      do i=0,im-1
         a_Int_ag(:) = a_Int_ag(:) + ag(:,i)*g_X_Weight(i)
      enddo
    end function a_Int_ag
ae :real(8), dimension(:,-km:), intent(in)

————— 微分計算 ——————

[Source]

    function ae_Dx_ae(ae)   ! スペクトルに作用する x 微分演算子

      real(8), dimension(:,-km:), intent(in)     :: ae
      real(8), dimension(size(ae,1),-km:km)      :: ae_dx_ae

      integer k

      if ( size(ae,2) < 2*km+1 ) then
         call MessageNotify('W','ae_Dx_ae',               'The Fourier dimension of input data too small.')
      elseif ( size(ae,2) > 2*km+1 ) then
         call MessageNotify('W','ae_Dx_ae',               'The Fourier dimension of input data too large.')
      endif

      do k=-km,km
         ae_Dx_ae(:,k) = -(2*pi*k/xl)*ae(:,-k)
      enddo
    end function ae_dx_ae
i :integer,intent(in)
: 格子点の設定
k :integer,intent(in)
: 切断波数の設定
xmin :real(8),intent(in)
: X 座標範囲
xmax :real(8),intent(in)
: X 座標範囲

————— 初期化 ——————

[Source]

    subroutine ae_initial(i,k,xmin,xmax)

      integer,intent(in) :: i           ! 格子点の設定
      integer,intent(in) :: k           ! 切断波数の設定

      real(8),intent(in) :: xmin, xmax     ! X 座標範囲

      integer :: ii

      im = i
      km = k
      xl = xmax-xmin

      if ( im <= 0 .or. km <= 0 ) then
         call MessageNotify('E','ae_initial',               'Number of grid points and waves should be positive')
      elseif ( mod(im,2) /= 0 ) then
         call MessageNotify('E','ae_initial',               'Number of grid points should be even')
      elseif ( km >= im/2 ) then
         call MessageNotify('E','ae_initial',               'KM shoud be less than IM/2')
      endif

      allocate(ti(im*2))

      call fttrui(im,iti,ti)

      allocate(g_x(0:im-1))
      do ii=0,im
         g_X(ii) = xmin + xl/im*ii
      enddo

      allocate(g_x_weight(0:im-1))
      g_X_Weight = xl/im

    end subroutine ae_initial
ag_ae :real(8), dimension(size(ae,1),0:im-1)
ae :real(8), dimension(:,-km:), intent(in)

————— 基本変換 ——————

[Source]


    function ag_ae(ae) ! スペクトル -> 格子
      real(8), dimension(:,-km:), intent(in)  :: ae
      real(8), dimension(size(ae,1),0:im-1)   :: ag_ae

      real(8), dimension(size(ae,1)*im)       :: y
      integer :: m, k

      m=size(ae,1)
      if ( size(ae,2) < 2*km+1 ) then
         call MessageNotify('E','ag_ae',               'The Fourier dimension of input data too small.')
      elseif ( size(ae,2) > 2*km+1 ) then
         call MessageNotify('W','ag_ae',               'The Fourier dimension of input data too large.')
      endif

      ag_ae = 0.0D0
      ag_ae(:,0)=ae(:,0)
      ag_ae(:,1)=0
      do k=1,km
         ag_ae(:,2*k)=ae(:,k)
         ag_ae(:,2*k+1)=ae(:,-k)
      enddo
      ag_ae(:,2*km+2:im-1)=0

      call fttrub(m,im,ag_ae,y,iti,ti)
    end function ag_ae

[Validate]