Class Special_Function
In: special_function.f90

— 特殊関数を計算するモジュール —

Methods

Public Instance methods

Full_Ellip1_Func( k ) result(Full_Ellip1_Func_d)
Function :
Full_Ellip1_Func_d :double precision
k :double precision, intent(in)
: 関数の引数

第 1 種完全楕円関数計算

Alias for Full_Ellip1_Func_d

Full_Ellip1_Func( k ) result(Full_Ellip1_Func_f)
Function :
Full_Ellip1_Func_f :real
k :real, intent(in)
: 関数の引数

第 1 種完全楕円関数計算

Alias for Full_Ellip1_Func_f

Function :
Full_Ellip1_Func_d :double precision
k :double precision, intent(in)
: 関数の引数

第 1 種完全楕円関数計算

[Source]

double precision function Full_Ellip1_Func_d(k)  ! 第 1 種完全楕円関数計算
  implicit none
  double precision, intent(in) :: k  ! 関数の引数
  double precision :: pi, m, dt, t, tmin, tmax
  integer :: i, j
  integer, parameter :: nmax=1000
  double precision :: f, x

  f(m,x) = 1.0d0/dsqrt(1.0d0-(m*dsin(x))**2)

  if(k.ge.1.0d0)then
     write(*,*) "(error ! : k must 0=<k<1.)"
     return
  end if

  pi = 3.1415926535898d0

  tmin = 0.0d0
  tmax = pi/2.0d0
  dt = (tmax-tmin)/dble(nmax-1)

  Full_Ellip1_Func_d = 0.5d0*dt*(f(k,tmin)+f(k,tmax))
  do i=1,nmax-2
     t = tmin+dt*dble(i)
     Full_Ellip1_Func_d = Full_Ellip1_Func_d+dt*f(k,t)
  end do

  return
end function
Function :
Full_Ellip1_Func_f :real
k :real, intent(in)
: 関数の引数

第 1 種完全楕円関数計算

[Source]

real function Full_Ellip1_Func_f(k)  ! 第 1 種完全楕円関数計算
  implicit none
  real, intent(in) :: k  ! 関数の引数
  real :: pi, m, dt, t, tmin, tmax
  integer :: i, j
  integer, parameter :: nmax=1000
  real :: f, x

  f(m,x) = 1.0/sqrt(1.0-(m*sin(x))**2)

  if(k.ge.1.0)then
     write(*,*) "(error ! : k must 0=<k<1.)"
     return
  end if

  pi = 3.14159265

  tmin = 0.0
  tmax = pi/2.0
  dt = (tmax-tmin)/(nmax-1)

  Full_Ellip1_Func_f = 0.5*dt*(f(k,tmin)+f(k,tmax))
  do i=1,nmax-2
     t = tmin+dt*i
     Full_Ellip1_Func_f = Full_Ellip1_Func_f+dt*f(k,t)
  end do

  return
end function
Full_Ellip2_Func( k ) result(Full_Ellip2_Func_d)
Function :
Full_Ellip2_Func_d :double precision
k :double precision, intent(in)
: 関数の引数

第二種完全楕円関数

Alias for Full_Ellip2_Func_d

Full_Ellip2_Func( k ) result(Full_Ellip2_Func_f)
Function :
Full_Ellip2_Func_f :real
k :real, intent(in)
: 関数の引数

第二種完全楕円関数

Alias for Full_Ellip2_Func_f

Function :
Full_Ellip2_Func_d :double precision
k :double precision, intent(in)
: 関数の引数

第二種完全楕円関数

[Source]

double precision function Full_Ellip2_Func_d(k)  ! 第二種完全楕円関数
  implicit none
  double precision, intent(in) :: k  ! 関数の引数
  double precision :: pi, m, dt, t, tmin, tmax
  integer :: i, j
  integer, parameter :: nmax=1000
  double precision :: f, x

  f(m,x) = dsqrt(1.0d0-(m*dsin(x))**2)

  pi = 3.1415926535898d0

  if(k.gt.1.0d0)then
     write(*,*) "(error) ! : k must 0=<k=<1."
     return
  end if

  tmin = 0.0d0
  tmax = pi/2.0d0
  dt = (tmax-tmin)/dble(nmax-1)

  Full_Ellip2_Func_d = 0.5d0*dt*(f(k,tmin)+f(k,tmax))
  do i=1,nmax-2
     t = tmin+dt*dble(i)
     Full_Ellip2_Func_d = Full_Ellip2_Func_d+dt*f(k,t)
  end do

  return
end function
Function :
Full_Ellip2_Func_f :real
k :real, intent(in)
: 関数の引数

第二種完全楕円関数

[Source]

real function Full_Ellip2_Func_f(k)  ! 第二種完全楕円関数
  implicit none
  real, intent(in) :: k  ! 関数の引数
  real :: pi, m, dt, t, tmin, tmax
  integer :: i, j
  integer, parameter :: nmax=1000
  real :: f, x

  f(m,x) = sqrt(1.0-(m*sin(x))**2)

  pi = 3.14159265

  if(k.gt.1.0)then
     write(*,*) "(error) ! : k must 0=<k=<1."
     return
  end if

  tmin = 0.0
  tmax = pi/2.0
  dt = (tmax-tmin)/(nmax-1)

  Full_Ellip2_Func_f = 0.5*dt*(f(k,tmin)+f(k,tmax))
  do i=1,nmax-2
     t = tmin+dt*i
     Full_Ellip2_Func_f = Full_Ellip2_Func_f+dt*f(k,t)
  end do

  return
end function
Subroutine :
nmax :integer, intent(in)
: ベッセル関数のゼロ点の最大個数
mmax :integer, intent(in)
: ベッセル関数の最大次数
k(0:nmax,mmax) :double precision, intent(inout)
: mmax 次までの nmax+1 個のゼロ点を格納する

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

 ベッセル関数のゼロ点を計算する *

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

[Source]

subroutine besdzero(nmax,mmax,k)
!**********************************
!  ベッセル関数のゼロ点を計算する *
!**********************************
  implicit none
  integer, intent(in) :: nmax          ! ベッセル関数のゼロ点の最大個数
  integer, intent(in) :: mmax          ! ベッセル関数の最大次数
  double precision, intent(inout) :: k(0:nmax,mmax)  ! mmax 次までの nmax+1 個のゼロ点を格納する
  double precision :: bessj_d, a, b, c, d, e, f, g, lim, dx
  integer :: i, j, n

!-- 二分法の解と近似する条件 ---
  lim=1.0d-6    ! 収束条件

!-- 二分法の二点を決定するための, 刻み幅 ---
!-- ベッセル関数のゼロ点の間隔はおよそ 3 ごとであるので,
!-- 0.5 ずつ刻めば, まあいいか
  dx=0.5d0
!-- 配列の初期化 ---
  do i=0,nmax
     do j=1,mmax
        k(i,j)=0.0d0
     end do
  end do

!-- 0 次計算 ---
  k(0,1)=0.0d0
  d=k(0,1)

  do 10 i=1,mmax

     if(i.gt.1)then
        d=k(0,i-1)+dx
     end if

     do while (k(0,mmax).eq.0.0d0)
        a=d
        e=bessj_d(0,a)
        b=a+dx
        f=bessj_d(0,b)
        d=d+dx

        do while (e*f.lt.0.0d0)
           c=0.5d0*(a+b)
           g=bessj_d(0,c)
           if(e*g.lt.0.0d0)then
              b=c
           else
              a=c
           end if

           if(abs(g).lt.lim)then
              k(0,i)=c
              go to 10
           end if

        end do
     end do 
  10 continue

  if(nmax > 0)then
!-- 1 次以上計算 ---
     do n=1,nmax
     do 21 i=1,mmax
        d=k(n-1,i)+dx
        do while (k(n,mmax).eq.0.0d0)
           a=d
           e=bessj_d(n,a)
           b=a+dx
           f=bessj_d(n,b)
           d=d+dx
           do while (e*f.lt.0.0d0)
              c=0.5d0*(a+b)
              g=bessj_d(n,c)
              if(e*g.lt.0.0d0)then
                 b=c
              else
                 a=c
              end if
              if(abs(g).lt.lim)then
                 k(n,i)=c
                 go to 21
              end if
           end do
        end do
  21  continue
     end do
  end if

end subroutine
Subroutine :
nmax :integer, intent(in)
: ベッセル関数のゼロ点の最大個数
mmax :integer, intent(in)
: ベッセル関数の最大次数
k(0:nmax,mmax) :real, intent(inout)
: mmax 次までの nmax+1 個のゼロ点を格納する

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

 ベッセル関数のゼロ点を計算する *

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

[Source]

subroutine besfzero(nmax,mmax,k)
!**********************************
!  ベッセル関数のゼロ点を計算する *
!**********************************
  implicit none
  integer, intent(in) :: nmax          ! ベッセル関数のゼロ点の最大個数
  integer, intent(in) :: mmax          ! ベッセル関数の最大次数
  real, intent(inout) :: k(0:nmax,mmax)  ! mmax 次までの nmax+1 個のゼロ点を格納する
  real :: bessj_f, a, b, c, d, e, f, g, lim, dx
  integer :: i, j, n

!-- 二分法の解と近似する条件 ---
  lim=1.0e-6    ! 収束条件

!-- 二分法の二点を決定するための, 刻み幅 ---
!-- ベッセル関数のゼロ点の間隔はおよそ 3 ごとであるので,
!-- 0.5 ずつ刻めば, まあいいか
!-- (注意)実際使用の際は, bessj_f 関数が参照されているかを確認のこと.
!-- バグ検証中
  dx=0.5
!-- 配列の初期化 ---
  do i=0,nmax
     do j=1,mmax
        k(i,j)=0.0
     end do
  end do

!-- 0 次計算 ---
  k(0,1)=0.0
  d=k(0,1)

  do 10 i=1,mmax

     if(i.gt.1)then
        d=k(0,i-1)+dx
     end if

     do while (k(0,mmax).eq.0.0)
        a=d
        e=bessj_f(0,a)
        b=a+dx
        f=bessj_f(0,b)
        d=d+dx

        do while (e*f.lt.0.0)
           c=0.5*(a+b)
           g=bessj_f(0,c)
           if(e*g.lt.0.0)then
              b=c
           else
              a=c
           end if

           if(abs(g).lt.lim)then
              k(0,i)=c
              go to 10
           end if

        end do
     end do 
  10 continue

  if(nmax > 0)then
!-- 1 次以上計算 ---
     do n=1,nmax
     do 21 i=1,mmax
        d=k(n-1,i)+dx
        do while (k(n,mmax).eq.0.0)
           a=d
           e=bessj_f(n,a)
           b=a+dx
           f=bessj_f(n,b)
           d=d+dx
           do while (e*f.lt.0.0)
              c=0.5*(a+b)
              g=bessj_f(n,c)
              if(e*g.lt.0.0)then
                 b=c
              else
                 a=c
              end if
              if(abs(g).lt.lim)then
                 k(n,i)=c
                 go to 21
              end if
           end do
        end do
  21  continue
     end do
  end if

end subroutine
bessj( m, t ) result(bessj_d)
Function :
bessj_d :double precision
m :integer, intent(in)
: 計算する次数
t :double precision, intent(in)
: 引数

第 I 種ベッセル関数を計算する

Alias for bessj_d

bessj( m, t ) result(bessj_f)
Function :
bessj_f :real
m :integer, intent(in)
: 計算する次数
t :real, intent(in)
: 引数

第 I 種ベッセル関数を計算する

Alias for bessj_f

Function :
bessj_d :double precision
m :integer, intent(in)
: 計算する次数
t :double precision, intent(in)
: 引数

第 I 種ベッセル関数を計算する

[Source]

double precision function bessj_d(m,t)  ! 第 I 種ベッセル関数を計算する
  implicit none
  integer, intent(in) :: m  ! 計算する次数
  double precision, intent(in) :: t  ! 引数
  integer :: istep, n
  double precision :: x
  integer, parameter :: mmax = 100 ! 数値積分用の配列
  double precision, parameter :: pis=3.14159265
  double precision, parameter :: xmin = 0.0d0, xmax = 2.0d0*pis
  double precision, parameter :: dx = (xmax-xmin)/(mmax-1)

!-- 負の次数であった場合を分ける ---
  if(m < 0)then
     n=-m
  else
     n=m
  end if

!-- ベッセル関数の積分計算 ---
  bessj_d=0.0d0

  do istep=2,mmax-1
     x=xmin+dx*dble(istep-1)
     bessj_d=bessj_d+dx*(dcos(t*dsin(x)-dble(n)*x))
  end do

  bessj_d=bessj_d+0.5d0*dx*(dcos(t*dsin(xmin)-dble(n)*xmin) +dcos(t*dsin(xmax)-dble(n)*xmax))
  bessj_d=bessj_d/(2.0d0*pis)

!-- 負の次数であった場合を分ける ---
  if(m.lt.0)then
     bessj_d=((-1.0d0)**n)*bessj_d
  end if

  return
end function
Function :
bessj_f :real
m :integer, intent(in)
: 計算する次数
t :real, intent(in)
: 引数

第 I 種ベッセル関数を計算する

[Source]

real function bessj_f(m,t)  ! 第 I 種ベッセル関数を計算する
  implicit none
  integer, intent(in) :: m  ! 計算する次数
  real, intent(in) :: t  ! 引数
  integer :: istep, n
  real :: x
  integer, parameter :: mmax = 100 ! 数値積分用の配列
  real, parameter :: pis=3.14159265
  real, parameter :: xmin = 0.0d0, xmax = 2.0d0*pis
  real, parameter :: dx = (xmax-xmin)/(mmax-1)

!-- 負の次数であった場合を分ける ---
  if(m < 0)then
     n=-m
  else
     n=m
  end if

!-- ベッセル関数の積分計算 ---
  bessj_f=0.0

  do istep=2,mmax-1
     x=xmin+dx*(istep-1)
     bessj_f=bessj_f+dx*(cos(t*sin(x)-real(n)*x))
  end do

  bessj_f=bessj_f+0.5*dx*(cos(t*sin(xmin)-real(n)*xmin) +cos(t*sin(xmax)-real(n)*xmax))
  bessj_f=bessj_f/(2.0*pis)

!-- 負の次数であった場合を分ける ---
  if(m.lt.0)then
     bessj_f=((-1.0)**n)*bessj_f
  end if

  return
end function
beszero( nmax, mmax, k )
Subroutine :
nmax :integer, intent(in)
: ベッセル関数のゼロ点の最大個数
mmax :integer, intent(in)
: ベッセル関数の最大次数
k(0:nmax,mmax) :double precision, intent(inout)
: mmax 次までの nmax+1 個のゼロ点を格納する

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

 ベッセル関数のゼロ点を計算する *

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

Alias for besdzero

beszero( nmax, mmax, k )
Subroutine :
nmax :integer, intent(in)
: ベッセル関数のゼロ点の最大個数
mmax :integer, intent(in)
: ベッセル関数の最大次数
k(0:nmax,mmax) :real, intent(inout)
: mmax 次までの nmax+1 個のゼロ点を格納する

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

 ベッセル関数のゼロ点を計算する *

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

Alias for besfzero

Function :
delta :real
t :integer, intent(in)
: 行成分
u :integer, intent(in)
: 列成分

クロネッカーのデルタを計算するサブルーチン

[Source]

real function delta(t,u)  ! クロネッカーのデルタを計算するサブルーチン
  implicit none
  integer, intent(in) :: t  ! 行成分
  integer, intent(in) :: u  ! 列成分

  if(t==u)then
     delta=1.0
  else
     delta=0.0
  end if

  return
end function
Function :
digamma :real
k :integer, intent(in)
: (k+1) 項目までの計算

— ダイガンマ関数を計算するサブルーチン — — 使い方 — — 関数名は "digamma(n)" で, 引数は必ず整数でなければならない

[Source]

real function digamma(k)
  !-- ダイガンマ関数を計算するサブルーチン ---
  !-- 使い方 ---
  !-- 関数名は "digamma(n)" で, 引数は必ず整数でなければならない
  implicit none
  integer, intent(in) :: k  ! (k+1) 項目までの計算
  integer :: j

  if (k.gt.1) then
     digamma=0.0
     do j=1,k
        digamma=digamma+1.0/j
     end do
  else
     if (k.eq.1) then
        digamma=1.0
     else
        digamma=0.0
     end if
  end if
  return
end function
Function :
epsilon :real
i :integer, intent(in)
: 第 1 成分
j :integer, intent(in)
: 第 1 成分
k :integer, intent(in)
: 第 1 成分

— エディントンのイプシロンを計算するサブルーチン — — F77 版では利用できなかった CASE 文で振り分けを行う — — i,j,k は 1..3 の 3 種類しか存在しないという仮定のもとの関数であるので, — 相対性理論でのテンソルには適用できない. —

[Source]

real function epsilon(i,j,k)
!-- エディントンのイプシロンを計算するサブルーチン ---
!-- F77 版では利用できなかった CASE 文で振り分けを行う ---
!-- i,j,k は 1..3 の 3 種類しか存在しないという仮定のもとの関数であるので,
!-- 相対性理論でのテンソルには適用できない. ---
  implicit none
  integer, intent(in) :: i  ! 第 1 成分
  integer, intent(in) :: j  ! 第 1 成分
  integer, intent(in) :: k  ! 第 1 成分

  select case (i)

  case (1)

     select case (j)

     case (1)
        epsilon=0.0

     case (2)

        select case (k)

        case (1)
           epsilon=0.0

        case (2)
           epsilon=0.0

        case (3)
           epsilon=1.0

        end select

     case (3)

        select case (k)

        case (1)
           epsilon=0.0

        case (2)
           epsilon=-1.0

        case (3)
           epsilon=0.0

        end select
     end select

  case (2)

     select case (j)

     case (1)

        select case (k)

        case (1)
           epsilon=0.0

        case (2)
           epsilon=0.0

        case (3)
           epsilon=-1.0

        end select

     case (2)
        epsilon=0.0

     case (3)

        select case (k)

        case (1)
           epsilon=1.0

        case (2)
           epsilon=0.0

        case (3)
           epsilon=0.0

        end select
     end select

  case (3)

     select case (j)

     case (1)

        select case (k)

        case (1)
           epsilon=0.0

        case (2)
           epsilon=1.0

        case (3)
           epsilon=0.0

        end select

     case (2)

        select case (k)

        case (1)
           epsilon=-1.0

        case (2)
           epsilon=0.0

        case (3)
           epsilon=0.0

        end select

     case (3)
        epsilon=0.0

     end select
  end select

  return
end function
Function :
kaijo :real
k :integer, intent(in)

— 階乗関数を計算するサブルーチン — — 使い方 — — 関数名は "kaijo(k)" で, 引数 "k" には整数のみ入れること

[Source]

real function kaijo(k)
  !-- 階乗関数を計算するサブルーチン ---
  !-- 使い方 ---
  !-- 関数名は "kaijo(k)" で, 引数 "k" には整数のみ入れること
  implicit none
  integer, intent(in) :: k
  integer :: j

  if (k.lt.2) then
     kaijo=1.0
  else
     kaijo=1.0
     do j=1,k
        kaijo=j*kaijo
     end do
  end if

  return
end function