Class Poly_Function
In: poly_function.f90

直交多項式を計算するサブルーチン集

Methods

Included Modules

special_function Math_Const

Public Instance methods

AS_LEGENDRE( n, m, x, p )
Subroutine :
n :integer, intent(in)
: 計算する legendre 陪関数の最高次数
m :integer, intent(in)
: 計算する legendre 陪関数の随伴次数
x(:) :double precision, intent(in)
: 引数
p(abs(m):n,size(x)) :double precision, intent(inout)
: 計算する legendre 陪関数

***********************************

 Legendre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for AS_LEGENDRE_d

AS_LEGENDRE( n, m, x, p )
Subroutine :
n :integer, intent(in)
: 計算する legendre 陪関数の最高次数
m :integer, intent(in)
: 計算する legendre 陪関数の随伴次数
x(:) :real, intent(in)
: 引数
p(abs(m):n,size(x)) :real, intent(inout)
: 計算する legendre 陪関数

***********************************

 Legendre 陪関数計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) m=随伴次数(n>=m) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for AS_LEGENDRE_f

Subroutine :
n :integer, intent(in)
: 計算する legendre 陪関数の最高次数
m :integer, intent(in)
: 計算する legendre 陪関数の随伴次数
x(:) :double precision, intent(in)
: 引数
p(abs(m):n,size(x)) :double precision, intent(inout)
: 計算する legendre 陪関数

***********************************

 Legendre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine AS_LEGENDRE_d(n, m, x, p)
!***********************************
!  Legendre 多項式計算サブルーチン *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  use special_function
  implicit none
  integer, intent(in) :: n  ! 計算する legendre 陪関数の最高次数
  integer, intent(in) :: m  ! 計算する legendre 陪関数の随伴次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(inout) :: p(abs(m):n,size(x))  ! 計算する legendre 陪関数
  integer :: i, j, abm
  integer :: nmax  ! 引数配列 x の要素数
  double precision :: coe

  nmax=size(x)
  abm=abs(m)

  if(abm>n)then
     write(*,*) "### ERROR : subroutine AS_LEGENDRE ###"
     write(*,*) "P^m_n : m must be less than n, stop."
     stop
  end if

!-- m についての 2 回階乗を計算

  coe=1.0d0

  if(abm>1)then
     do i=2,abm
        coe=(2.0d0*dble(i)-1.0d0)*coe
     end do
  end if

!-- m=0 の場合はルジャンドル多項式と同じなので, それを流用.

  if(m==0)then
     call Legendre_d( n, x, p(:,:) )
  else

  !-- 以下, m/=0の場合の計算.
  !-- 初項の設定 ---

     do i=1,nmax
        p(abm,i)=dsqrt(coe*(1.0d0-x(i)**2)**abm)
     end do

     if(abm<n)then  ! |m|<n の場合の計算.

        do i=1,nmax
           p(abm+1,i)=p(abm,i)*(2.0d0*dble(abm)+1.0d0)*x(i)
        end do

        if(abm+1<n)then  ! |m+1|<n の場合の計算.

        !-- 漸化式の計算 ---
           do j=abm+1,n-1
              do i=1,nmax
                 p(j+1,i)=(1.0d0/dble(j-abm+1))*(p(j,i)*(2.0d0*dble(j)+1.0d0) *x(i)-dble(j+abm)*p(j-1,i))
              end do
           end do

        end if
     end if

     if(m<0)then  ! 負次の場合の追加計算
        coe=dble(kaijo_i(n-m))/dble(kaijo_i(n+m))

        if(mod(m,2)==0)then
           do j=abm,n
              do i=1,nmax
                 p(j,i)=coe*p(j,i)
              end do
           end do
        else
           do j=abm,n
              do i=1,nmax
                 p(j,i)=-coe*p(j,i)
              end do
           end do
        end if

     end if

  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算する legendre 陪関数の最高次数
m :integer, intent(in)
: 計算する legendre 陪関数の随伴次数
x(:) :real, intent(in)
: 引数
p(abs(m):n,size(x)) :real, intent(inout)
: 計算する legendre 陪関数

***********************************

 Legendre 陪関数計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) m=随伴次数(n>=m) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine AS_LEGENDRE_f(n, m, x, p)
!***********************************
!  Legendre 陪関数計算サブルーチン *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! m=随伴次数(n>=m)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  use special_function
  implicit none
  integer, intent(in) :: n  ! 計算する legendre 陪関数の最高次数
  integer, intent(in) :: m  ! 計算する legendre 陪関数の随伴次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(inout) :: p(abs(m):n,size(x))  ! 計算する legendre 陪関数
  integer :: i, j, abm
  integer :: nmax  ! 引数配列 x の要素数
  real :: coe

  nmax=size(x)
  abm=abs(m)

  if(abm>n)then
     write(*,*) "### ERROR : subroutine AS_LEGENDRE ###"
     write(*,*) "P^m_n : m must be less than n, stop."
     stop
  end if

!-- m についての 2 回階乗を計算

  coe=1.0

  if(abm>1)then
     do i=2,abm
        coe=(2.0*real(i)-1.0)*coe
     end do
  end if

!-- m=0 の場合はルジャンドル多項式と同じなので, それを流用.

  if(m==0)then
     call Legendre_f( n, x, p(:,:) )
  else

  !-- 以下, m/=0の場合の計算.
  !-- 初項の設定 ---

     do i=1,nmax
        p(abm,i)=sqrt(coe*(1.0-x(i)**2)**abm)
     end do

     if(abm<n)then  ! |m|<n の場合の計算.

        do i=1,nmax
           p(abm+1,i)=p(abm,i)*(2.0*real(abm)+1.0)*x(i)
        end do

        if(abm+1<n)then  ! |m+1|<n の場合の計算.

        !-- 漸化式の計算 ---
           do j=abm+1,n-1
              do i=1,nmax
                 p(j+1,i)=(1.0/real(j-abm+1))*(p(j,i)*(2.0*real(j)+1.0) *x(i)-real(j+abm)*p(j-1,i))
              end do
           end do

        end if
     end if

     if(m<0)then  ! 負次の場合の追加計算
        coe=real(kaijo_i(n-m))/real(kaijo_i(n+m))

        if(mod(m,2)==0)then
           do j=abm,n
              do i=1,nmax
                 p(j,i)=coe*p(j,i)
              end do
           end do
        else
           do j=abm,n
              do i=1,nmax
                 p(j,i)=-coe*p(j,i)
              end do
           end do
        end if

     end if

  end if

end subroutine
CHEBYSHEV( n, x, che )
Subroutine :
n :integer, intent(in)
: 計算するチェビシェフの最高次数
x(:) :double precision, intent(in)
: チェビシェフ多項式の引数
che(0:n,size(x)) :double precision, intent(inout)

**************************************** *** チェビシェフ漸化式のサブルーチン *** ****************************************

Alias for CHEBYSHEV_d

CHEBYSHEV( n, x, che )
Subroutine :
n :integer, intent(in)
: 計算するチェビシェフの最高次数
x(:) :real, intent(in)
: チェビシェフ多項式の引数
che(0:n,size(x)) :real, intent(inout)
: 計算するチェビシェフ多項式

**************************************** *** チェビシェフ漸化式のサブルーチン *** ****************************************

Alias for CHEBYSHEV_f

Subroutine :
n :integer, intent(in)
: 計算するチェビシェフの最高次数
x(:) :double precision, intent(in)
: チェビシェフ多項式の引数
che(0:n,size(x)) :double precision, intent(inout)

**************************************** *** チェビシェフ漸化式のサブルーチン *** ****************************************

[Source]

subroutine CHEBYSHEV_d(n, x, che)
!****************************************
!*** チェビシェフ漸化式のサブルーチン ***
!****************************************
  implicit none
  integer, intent(in) :: n  ! 計算するチェビシェフの最高次数
  double precision, intent(in) :: x(:)  ! チェビシェフ多項式の引数
  double precision, intent(inout) :: che(0:n,size(x))
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     che(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
        che(1,i)=x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              che(j+1,i)=2.0d0*che(1,i)*che(j,i)-che(j-1,i)
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算するチェビシェフの最高次数
x(:) :real, intent(in)
: チェビシェフ多項式の引数
che(0:n,size(x)) :real, intent(inout)
: 計算するチェビシェフ多項式

**************************************** *** チェビシェフ漸化式のサブルーチン *** ****************************************

[Source]

subroutine CHEBYSHEV_f(n, x, che)
!****************************************
!*** チェビシェフ漸化式のサブルーチン ***
!****************************************
  implicit none
  integer, intent(in) :: n  ! 計算するチェビシェフの最高次数
  real, intent(in) :: x(:)  ! チェビシェフ多項式の引数
  real, intent(inout) :: che(0:n,size(x))  ! 計算するチェビシェフ多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     che(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
        che(1,i)=x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              che(j+1,i)=2.0*che(1,i)*che(j,i)-che(j-1,i)
           end do
        end do
     end if
  end if

end  subroutine
GEGENBAUER( n, x, p, lambda )
Subroutine :
n :integer, intent(in)
: 計算するゲーゲンバウアー多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算するゲーゲンバウアー多項式
lambda :double precision, intent(in)
: ゲーゲンバウアー係数

************************************

  • ゲーゲンバウアー 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

Alias for GEGENBAUER_d

GEGENBAUER( n, x, p, lambda )
Subroutine :
n :integer, intent(in)
: 計算するゲーゲンバウアー多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算するゲーゲンバウアー多項式
lambda :real, intent(in)
: ゲーゲンバウアー係数

************************************

  • ゲーゲンバウアー 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

Alias for GEGENBAUER_f

Subroutine :
n :integer, intent(in)
: 計算するゲーゲンバウアー多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算するゲーゲンバウアー多項式
lambda :double precision, intent(in)
: ゲーゲンバウアー係数

************************************

  • ゲーゲンバウアー 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

[Source]

subroutine GEGENBAUER_d(n, x, p, lambda)
!************************************
!*  ゲーゲンバウアー 多項式計算サブルーチン  *
!************************************
!* 使い方
!* n=次数(0次から指定可能)
!************************************
  implicit none
  integer, intent(in) :: n  ! 計算するゲーゲンバウアー多項式の最高次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(in) :: lambda  ! ゲーゲンバウアー係数
  double precision, intent(inout) :: p(0:n,size(x))  ! 計算するゲーゲンバウアー多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=2.0d0*lambda*x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(1.0d0/dble(j+1))*(2.0d0*(lambda+dble(j))*x(i)*p(j,i) -(2.0d0*lambda+dble(j-1))*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算するゲーゲンバウアー多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算するゲーゲンバウアー多項式
lambda :real, intent(in)
: ゲーゲンバウアー係数

************************************

  • ゲーゲンバウアー 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

[Source]

subroutine GEGENBAUER_f(n, x, p, lambda)
!************************************
!*  ゲーゲンバウアー 多項式計算サブルーチン  *
!************************************
!* 使い方
!* n=次数(0次から指定可能)
!************************************
  implicit none
  integer, intent(in) :: n  ! 計算するゲーゲンバウアー多項式の最高次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(in) :: lambda  ! ゲーゲンバウアー係数
  real, intent(inout) :: p(0:n,size(x))  ! 計算するゲーゲンバウアー多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=2.0*lambda*x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(1.0/real(j+1))*(2.0*(lambda+real(j))*x(i)*p(j,i) -(2.0*lambda+real(j-1))*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
HERMITE( n, x, p )
Subroutine :
n :integer, intent(in)
: 計算する Hermit 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Hermit 多項式

************************************

  • Hermite 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

Alias for HERMITE_d

HERMITE( n, x, p )
Subroutine :
n :integer, intent(in)
: 計算する Hermite 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Hermite 多項式

************************************

  • Hermite 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

Alias for HERMITE_f

Subroutine :
n :integer, intent(in)
: 計算する Hermit 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Hermit 多項式

************************************

  • Hermite 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

[Source]

subroutine HERMITE_d(n, x, p)
!************************************
!*  Hermite 多項式計算サブルーチン  *
!************************************
!* 使い方
!* n=次数(0次から指定可能)
!************************************
  implicit none
  integer, intent(in) :: n  ! 計算する Hermit 多項式の最高次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(inout) :: p(0:n,size(x))  ! 計算する Hermit 多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=2.0d0*x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=2.0d0*(x(i)*p(j,i)-dble(j)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算する Hermite 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Hermite 多項式

************************************

  • Hermite 多項式計算サブルーチン *

************************************

  • 使い方
  • n=次数(0次から指定可能)

************************************

[Source]

subroutine HERMITE_f(n, x, p)
!************************************
!*  Hermite 多項式計算サブルーチン  *
!************************************
!* 使い方
!* n=次数(0次から指定可能)
!************************************
  implicit none
  integer, intent(in) :: n  ! 計算する Hermite 多項式の最高次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(inout) :: p(0:n,size(x))  ! 計算する Hermite 多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=2.0*x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=2.0*(x(i)*p(j,i)-real(j)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
JACOBI_POLY( n, x, p, alpha, beta )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式
alpha :double precision, intent(in)
: 第一引数
beta :double precision, intent(in)
: 第二引数

***********************************

 Jacobi 多項式計算サブルーチン   *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for JACOBI_POLY_d

JACOBI_POLY( n, x, p, alpha, beta )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式
alpha :real, intent(in)
: 第一引数
beta :real, intent(in)
: 第二引数

***********************************

 Jacobi 多項式計算サブルーチン   *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for JACOBI_POLY_f

Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式
alpha :double precision, intent(in)
: 第一引数
beta :double precision, intent(in)
: 第二引数

***********************************

 Jacobi 多項式計算サブルーチン   *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine JACOBI_POLY_d(n, x, p, alpha, beta)
!***********************************
!  Jacobi 多項式計算サブルーチン   *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  double precision, intent(in) :: alpha  ! 第一引数
  double precision, intent(in) :: beta  ! 第二引数
  double precision :: gamma, omega
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 係数の設定 ---
  gamma=alpha+beta
  omega=alpha-beta

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=0.5d0*((gamma+2.0d0)*x(i)+omega)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(0.5d0/(dble(j+1)*dble(j+1+gamma)*dble(2.0*j+gamma))) *((2.0d0*dble(j)+gamma+1.0d0) *(gamma*omega+(2.0d0*dble(j)+gamma) *(2.0d0*dble(j+1)+gamma)*x(i))*p(j,i) -2.0d0*dble(j+alpha)*dble(j+beta) *(2.0d0*dble(j+1)+gamma)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式
alpha :real, intent(in)
: 第一引数
beta :real, intent(in)
: 第二引数

***********************************

 Jacobi 多項式計算サブルーチン   *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine JACOBI_POLY_f(n, x, p, alpha, beta)
!***********************************
!  Jacobi 多項式計算サブルーチン   *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  real, intent(in) :: alpha  ! 第一引数
  real, intent(in) :: beta  ! 第二引数
  real :: gamma, omega
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 係数の設定 ---
  gamma=alpha+beta
  omega=alpha-beta

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=0.5*((gamma+2.0)*x(i)+omega)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(0.5/(real(j+1)*real(j+1+gamma)*real(2.0*j+gamma))) *((2.0*j+gamma+1.0) *(gamma*omega+(2.0*j+gamma)*(2.0*(j+1)+gamma)*x(i)) *p(j,i) -2.0*(j+alpha)*(j+beta)*(2.0*(j+1)+gamma)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
LAGUERRE( n, x, p )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Laguerre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for LAGUERRE_d

LAGUERRE( n, x, p )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Laguerre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for LAGUERRE_f

Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Laguerre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine LAGUERRE_d(n, x, p)
!***********************************
!  Laguerre 多項式計算サブルーチン *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=1.0d0-x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(2.0d0*dble(j)+1.0d0-x(i))*p(j,i) -((dble(j))**2)*p(j-1,i)
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Laguerre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine LAGUERRE_f(n, x, p)
!***********************************
!  Laguerre 多項式計算サブルーチン *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=1.0-x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(2.0*real(j)+1.0-x(i))*p(j,i) -((real(j))**2)*p(j-1,i)
           end do
        end do
     end if
  end if

end subroutine
LEGENDRE( n, x, p )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Legendre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for LEGENDRE_d

LEGENDRE( n, x, p )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Legendre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for LEGENDRE_f

Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Legendre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine LEGENDRE_d(n, x, p)
!***********************************
!  Legendre 多項式計算サブルーチン *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(1.0d0/dble(j+1))*(p(j,i)*(2.0d0*dble(j)+1.0d0) *(p(1,i))-dble(j)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式

***********************************

 Legendre 多項式計算サブルーチン *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine LEGENDRE_f(n, x, p)
!***********************************
!  Legendre 多項式計算サブルーチン *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
        p(1,i)=x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(1.0/real(j+1))*(p(j,i)*(2.0*real(j)+1.0) *(p(1,i))-real(j)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
SONINE( n, x, p, lambda )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式
lambda :double precision, intent(in)
: 第一引数

***********************************

 Sonine 多項式計算サブルーチン  *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for SONINE_d

SONINE( n, x, p, lambda )
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式
lambda :real, intent(in)
: 第一引数

***********************************

 Sonine 多項式計算サブルーチン  *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

Alias for SONINE_f

Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :double precision, intent(in)
: 引数
p(0:n,size(x)) :double precision, intent(inout)
: 計算する Jacobi 多項式
lambda :double precision, intent(in)
: 第一引数

***********************************

 Sonine 多項式計算サブルーチン  *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine SONINE_d(n, x, p, lambda)
!***********************************
!  Sonine 多項式計算サブルーチン  *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  double precision, intent(in) :: x(:)  ! 引数
  double precision, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  double precision, intent(in) :: lambda  ! 第一引数
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0d0
  end do

  if(n > 0)then

     do i=1,nmax
         p(1,i)=lambda+1.0d0-x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(1.0d0/dble(j+1)) *((lambda+2.0d0*dble(j)+1.0d0-x(i))*p(j,i) -(dble(j)+lambda)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: 計算する jacobi 多項式の最高次数
x(:) :real, intent(in)
: 引数
p(0:n,size(x)) :real, intent(inout)
: 計算する Jacobi 多項式
lambda :real, intent(in)
: 第一引数

***********************************

 Sonine 多項式計算サブルーチン  *

*********************************** 使い方 n=次数(0次から指定可能) nmax=空間格子点数 p=p(0:n,nmax) の2次元配列 ***********************************

[Source]

subroutine SONINE_f(n, x, p, lambda)
!***********************************
!  Sonine 多項式計算サブルーチン  *
!***********************************
! 使い方
! n=次数(0次から指定可能)
! nmax=空間格子点数
! p=p(0:n,nmax) の2次元配列
!***********************************
  implicit none
  integer, intent(in) :: n  ! 計算する jacobi 多項式の最高次数
  real, intent(in) :: x(:)  ! 引数
  real, intent(inout) :: p(0:n,size(x))  ! 計算する Jacobi 多項式
  real, intent(in) :: lambda  ! 第一引数
  integer :: i, j
  integer :: nmax  ! 引数配列 x の要素数

  nmax=size(x)

!-- 初項の設定 ---
  do i=1,nmax
     p(0,i)=1.0
  end do

  if(n > 0)then

     do i=1,nmax
         p(1,i)=lambda+1.0-x(i)
     end do

     if(n > 1)then
!-- 漸化式の計算 ---
        do j=1,n-1
           do i=1,nmax
              p(j+1,i)=(1.0/real(j+1)) *((lambda+2.0*j+1.0-x(i))*p(j,i) -(j+lambda)*p(j-1,i))
           end do
        end do
     end if
  end if

end subroutine
Subroutine :
n :integer, intent(in)
: ガウス緯度の分点の数
lat(n) :real, intent(inout)
: ガウス緯度の各点での緯度 [degree]
eps :real, intent(in), optional
: ニュートン法でのゼロ点判定値 デフォルトは 1.0e-6.

ガウス緯度を計算するルーチン

[Source]

subroutine gauss_lat( n, lat, eps )
! ガウス緯度を計算するルーチン
  use Math_Const
  implicit none
  integer, intent(in) :: n  ! ガウス緯度の分点の数
  real, intent(inout) :: lat(n)  ! ガウス緯度の各点での緯度 [degree]
  real, intent(in), optional :: eps  ! ニュートン法でのゼロ点判定値
                                ! デフォルトは 1.0e-6.
  real :: tmp(0:n,1)
  real :: nu(1)
  integer :: k
  real :: diff, eps_max

  if(present(eps))then
     eps_max=eps
  else
     eps_max=1.0e-6
  end if

  do k=1,n  ! k 番目の分点

     nu(1)=cos((pi*(4.0*real(k)-1.0))/(4.0*real(n)+2.0))
     call legendre( n, nu(1:1), tmp(0:n,1:1) )

     do while (eps_max<=tmp(n,1))
        diff=(1.0-nu(1)**2)/(real(n)*(tmp(n-1,1)-nu(1)*tmp(n,1)))
        nu(1)=nu(1)-tmp(n,1)*diff
        call legendre( n, nu(1:1), tmp(0:n,1:1) )
     end do

     lat(k)=asin(nu(1))*180.0/pi

  end do

end subroutine
Subroutine :
n :integer, intent(in)
: 緯度方向の波数
m :integer, intent(in)
: 経度方向の波数
x(:) :real, intent(in)
: 経度変数 [0<=x<=2*pi]
y(:) :real, intent(in)
: 緯度変数 [-1<=y<=1]
p(0:n,-m:m,size(x),size(y)) :complex, intent(inout)
: 球面調和関数

球面調和関数を計算するルーチン

[Source]

subroutine ymn( n, m, x, y, p )
! 球面調和関数を計算するルーチン
  implicit none
  integer, intent(in) :: n  ! 緯度方向の波数
  integer, intent(in) :: m  ! 経度方向の波数
  real, intent(in) :: x(:)  ! 経度変数 [0<=x<=2*pi]
  real, intent(in) :: y(:)  ! 緯度変数 [-1<=y<=1]
  complex, intent(inout) :: p(0:n,-m:m,size(x),size(y))  ! 球面調和関数
  integer :: nmax, mmax
  integer :: i, j, k, l
  complex, parameter :: img=(0.0,1.0)
  real :: pm(0:n,-m:m,size(y))

  nmax=size(x)
  mmax=size(y)

!-- ルジャンドル陪関数を計算.
  do j=-m,m
     call AS_Legendre( n, j, y, pm(abs(j):n,j,:) )
  end do

  do l=-m,m
     do k=abs(l),n
        do j=1,mmax
           do i=1,nmax
              p(k,l,i,j)=pm(k,l,j)*(cos(real(m)*x(i))+img*(sin(real(m)*x(i))))
           end do
        end do
     end do
  end do

end subroutine