Class Statistics
In: statistics.f90

統計解析関係のルーチン集

Methods

Included Modules

Matrix_Calc

Public Instance methods

Subroutine :
x(:) :real, intent(in)
: データ
anor(size(x)) :real, intent(inout)
: 各 x(i) に対応する偏差 anor(i)
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

1 次元データ配列の偏差を返す

[Source]

subroutine Anomaly_1d( x, anor, error )  ! 1 次元データ配列の偏差を返す
  implicit none
  real, intent(in) :: x(:)  ! データ
  real, intent(inout) :: anor(size(x))  ! 各 x(i) に対応する偏差 anor(i)
  real, intent(in), optional :: error  ! 欠損値が存在するデータセットの場合の欠損値
  integer :: i
  integer :: nx  ! データの要素数
  real :: ave

  nx=size(x)

  if(present(error))then
     call Mean_1d( x, ave, error )
     do i=1,nx
        if(x(i)==error)then
           anor(i)=error
        else
           anor(i)=x(i)-ave
        end if
     end do
  else
     call Mean_1d( x, ave )
     do i=1,nx
        anor(i)=x(i)-ave
     end do
  end if

end subroutine Anomaly_1d
Subroutine :
x(:,:) :real, intent(in)
: データ
anor(size(x,1),size(x,2)) :real, intent(inout)
: 各 x(i,j) に対応する偏差 anor(i,j)
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

2 次元データ配列の偏差を返す

[Source]

subroutine Anomaly_2d( x, anor, error )  ! 2 次元データ配列の偏差を返す
  implicit none
  real, intent(in) :: x(:,:)  ! データ
  real, intent(inout) :: anor(size(x,1),size(x,2))  ! 各 x(i,j) に対応する偏差 anor(i,j)
  real, intent(in), optional :: error  ! 欠損値が存在するデータセットの場合の欠損値
  integer :: i, j
  integer :: nx  ! データの要素数 1
  integer :: ny  ! データの要素数 2
  real :: ave

  nx=size(x,1)
  ny=size(x,2)

  if(present(error))then
     call Mean_2d( x, ave, error )
     do j=1,ny
        do i=1,nx
           if(x(i,j)==error)then
              anor(i,j)=error
           else
              anor(i,j)=x(i,j)-ave
           end if
        end do
     end do
  else
     call Mean_2d( x, ave, error )
     do j=1,ny
        do i=1,nx
           anor(i,j)=x(i,j)-ave
        end do
     end do
  end if

end subroutine Anomaly_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ
anor(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 各 x(i,j,k) に対応する偏差 anor(i,j,k)
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

3 次元データ配列の偏差を返す

[Source]

subroutine Anomaly_3d( x, anor, error )  ! 3 次元データ配列の偏差を返す
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ
  real, intent(inout) :: anor(size(x,1),size(x,2),size(x,3))  ! 各 x(i,j,k) に対応する偏差 anor(i,j,k)
  real, intent(in), optional :: error  ! 欠損値が存在するデータセットの場合の欠損値
  integer :: i, j, k
  integer :: nx  ! データの要素数 1
  integer :: ny  ! データの要素数 2
  integer :: nz  ! データの要素数 3
  real :: ave

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  if(present(error))then
     call Mean_3d( x, ave, error )
     do k=1,nz
        do j=1,ny
           do i=1,nx
              if(x(i,j,k)==error)then
                 anor(i,j,k)=error
              else
                 anor(i,j,k)=x(i,j,k)-ave
              end if
           end do
        end do
     end do
  else
     call Mean_3d( x, ave, error )
     do k=1,nz
        do j=1,ny
           do i=1,nx
              anor(i,j,k)=x(i,j,k)-ave
           end do
        end do
     end do
  end if

end subroutine Anomaly_3d
Bubble_Sort( a, b, sig )
Subroutine :
a(:) :integer, intent(in)
: ソートする配列
b(size(a)) :integer, intent(inout)
: ソートした結果を格納する配列
sig :character(1), intent(in)
: ソートの順番 ‘i’ = 要素番号の若いものに小さい値が入る ‘r’ = 要素番号の若いものに大きい値が入る

バブルソートを用いて数値データを sig の方向にソートする.

Alias for Bubble_Sort_i

Bubble_Sort( a, b, sig )
Subroutine :
a(:) :real, intent(in)
: ソートする配列
b(size(a)) :real, intent(inout)
: ソートした結果を格納する配列
sig :character(1), intent(in)
: ソートの順番 ‘i’ = 要素番号の若いものに小さい値が入る ‘r’ = 要素番号の若いものに大きい値が入る

バブルソートを用いて数値データを sig の方向にソートする.

Alias for Bubble_Sort_f

Subroutine :
a(:) :real, intent(in)
: ソートする配列
b(size(a)) :real, intent(inout)
: ソートした結果を格納する配列
sig :character(1), intent(in)
: ソートの順番 ‘i’ = 要素番号の若いものに小さい値が入る ‘r’ = 要素番号の若いものに大きい値が入る

バブルソートを用いて数値データを sig の方向にソートする.

[Source]

subroutine Bubble_Sort_f( a, b, sig )
! バブルソートを用いて数値データを sig の方向にソートする.
  implicit none
  real, intent(in) :: a(:)  ! ソートする配列
  real, intent(inout) :: b(size(a))  ! ソートした結果を格納する配列
  character(1), intent(in) :: sig  ! ソートの順番
                                   ! 'i' = 要素番号の若いものに小さい値が入る
                                   ! 'r' = 要素番号の若いものに大きい値が入る
  integer :: i, j, n
  real :: tmp

  n=size(a)

  if(sig/='i'.and.sig/='r')then
     write(*,*) "### ERROR ###"
     write(*,*) "sig flag is 'r' .or. 'i', STOP."
     stop
  end if

  do i=1,n
     b(i)=a(i)
  end do

  if(sig=='i')then  ! 昇べきソート
     do i=1,n
        do j=1,n-1
           if(b(j)>b(j+1))then
              tmp=b(j+1)
              b(j+1)=b(j)
              b(j)=tmp
           end if
        end do
     end do
  else
     do i=1,n
        do j=1,n-1
           if(b(j)<b(j+1))then
              tmp=b(j+1)
              b(j+1)=b(j)
              b(j)=tmp
           end if
        end do
     end do
  end if

end subroutine Bubble_Sort_f
Subroutine :
a(:) :integer, intent(in)
: ソートする配列
b(size(a)) :integer, intent(inout)
: ソートした結果を格納する配列
sig :character(1), intent(in)
: ソートの順番 ‘i’ = 要素番号の若いものに小さい値が入る ‘r’ = 要素番号の若いものに大きい値が入る

バブルソートを用いて数値データを sig の方向にソートする.

[Source]

subroutine Bubble_Sort_i( a, b, sig )
! バブルソートを用いて数値データを sig の方向にソートする.
  implicit none
  integer, intent(in) :: a(:)  ! ソートする配列
  integer, intent(inout) :: b(size(a))  ! ソートした結果を格納する配列
  character(1), intent(in) :: sig  ! ソートの順番
                                   ! 'i' = 要素番号の若いものに小さい値が入る
                                   ! 'r' = 要素番号の若いものに大きい値が入る
  integer :: i, j, n
  integer :: tmp

  n=size(a)

  if(sig/='i'.and.sig/='r')then
     write(*,*) "### ERROR ###"
     write(*,*) "sig flag is 'r' .or. 'i', STOP."
     stop
  end if

  do i=1,n
     b(i)=a(i)
  end do

  if(sig=='i')then  ! 昇べきソート
     do i=1,n
        do j=1,n-1
           if(b(j)>b(j+1))then
              tmp=b(j+1)
              b(j+1)=b(j)
              b(j)=tmp
           end if
        end do
     end do
  else
     do i=1,n
        do j=1,n-1
           if(b(j)<b(j+1))then
              tmp=b(j+1)
              b(j+1)=b(j)
              b(j)=tmp
           end if
        end do
     end do
  end if

end subroutine Bubble_Sort_i
Cor_Coe( x, y, cc, [error] )
Subroutine :
x(:) :real, intent(in)
: データ要素 1
y(size(x)) :real, intent(in)
: データ要素 2
cc :real, intent(inout)
: 相関係数
error :real, intent(in), optional
: 欠損値

2 データの相関係数を計算するルーチン

Alias for Cor_Coe_1d

Subroutine :
x(:) :real, intent(in)
: データ要素 1
y(size(x)) :real, intent(in)
: データ要素 2
cc :real, intent(inout)
: 相関係数
error :real, intent(in), optional
: 欠損値

2 データの相関係数を計算するルーチン

[Source]

subroutine Cor_Coe_1d( x, y ,cc, error )  ! 2 データの相関係数を計算するルーチン
  implicit none
  real, intent(in) :: x(:)  ! データ要素 1
  real, intent(in) :: y(size(x))  ! データ要素 2
  real, intent(inout) :: cc  ! 相関係数
  real, intent(in), optional :: error  ! 欠損値
  integer :: nx  ! データ個数
  real :: cov, anor1, anor2

  nx=size(x)

  if(present(error))then
     call covariance( x, y, cov, error )
     call stand_vari( x, anor1, error )
     call stand_vari( y, anor2, error )
  else
     call covariance( x, y, cov )
     call stand_vari( x, anor1 )
     call stand_vari( y, anor2 )
  end if

  cc=cov/(sqrt(anor1)*sqrt(anor2))

end subroutine Cor_Coe_1d
Subroutine :
x(:,:) :real, intent(in)
: データ要素 1
y(size(x,1),size(x,2)) :real, intent(in)
: データ要素 2
cc :real, intent(inout)
: 相関係数
error :real, intent(in), optional
: 欠損値

2 データの相関係数を計算するルーチン (2 次元版)

[Source]

subroutine Cor_Coe_2d( x, y ,cc, error )  ! 2 データの相関係数を計算するルーチン (2 次元版)
  implicit none
  real, intent(in) :: x(:,:)  ! データ要素 1
  real, intent(in) :: y(size(x,1),size(x,2))  ! データ要素 2
  real, intent(inout) :: cc  ! 相関係数
  real, intent(in), optional :: error  ! 欠損値
  integer :: i, j, counter
  integer :: nx  ! データ個数 1
  integer :: ny  ! データ個数 2
  real, dimension(size(x,1)*size(x,2)) :: val1, val2

  nx=size(x,1)
  ny=size(x,2)

  counter=0

  do j=1,ny
     do i=1,nx
        counter=counter+1
        val1(counter)=x(i,j)
        val2(counter)=y(i,j)
     end do
  end do

  if(present(error))then
     call Cor_Coe_1d( val1, val2, cc, error )
  else
     call Cor_Coe_1d( val1, val2, cc )
  end if

end subroutine Cor_Coe_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ要素 1
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: データ要素 2
cc :real, intent(inout)
: 相関係数
error :real, intent(in), optional
: 欠損値

2 データの相関係数を計算するルーチン (3 次元版)

[Source]

subroutine Cor_Coe_3d( x, y ,cc, error )  ! 2 データの相関係数を計算するルーチン (3 次元版)
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ要素 1
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! データ要素 2
  real, intent(inout) :: cc  ! 相関係数
  real, intent(in), optional :: error  ! 欠損値
  integer :: i, j, k, counter
  integer :: nx  ! データ個数 1
  integer :: ny  ! データ個数 2
  integer :: nz  ! データ個数 2
  real, dimension(size(x,1)*size(x,2)*size(x,3)) :: val1, val2

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  counter=0

  do k=1,nz
     do j=1,ny
        do i=1,nx
           counter=counter+1
           val1(counter)=x(i,j,k)
           val2(counter)=y(i,j,k)
        end do
     end do
  end do

  if(present(error))then
     call Cor_Coe_1d( val1, val2, cc, error )
  else
     call Cor_Coe_1d( val1, val2, cc )
  end if

end subroutine Cor_Coe_3d
LSM( x, y, slope, intercept, [undef] )
Subroutine :
x(:) :real, intent(in)
: データ要素 1
y(size(x)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
undef :real, intent(in), optional
: undef

最小二乗法による傾きと切片計算 (1 次元データ版)

Alias for LSM_1d

Subroutine :
x(:) :real, intent(in)
: データ要素 1
y(size(x)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
undef :real, intent(in), optional
: undef

最小二乗法による傾きと切片計算 (1 次元データ版)

[Source]

subroutine LSM_1d( x, y, slope, intercept, undef )  ! 最小二乗法による傾きと切片計算 (1 次元データ版)
  implicit none
  real, intent(in) :: x(:)  ! データ要素 1
  real, intent(in) :: y(size(x))  ! データ要素 2
  real, intent(inout) :: slope  ! 最適な傾き
  real, intent(inout) :: intercept  ! 最適な切片
  real, intent(in), optional :: undef ! undef 
  real :: u(size(x)), v(size(x))
  integer :: i
  integer :: nx  ! データ数
  real :: a, b, c, d

  nx=size(x)
  a=0.0
  b=0.0
  c=0.0
  d=0.0

!$omp parallel do shared(u, v, x, y) private(i)
  do i=1,nx
     u(i)=x(i)*x(i)
     v(i)=x(i)*y(i)
  end do
!$omp end parallel do

  if(present(undef))then
     call summ(v,a,undef)
     call summ(x,b,undef)
     call summ(y,c,undef)
     call summ(u,d,undef)
  else
     call summ(v,a)
     call summ(x,b)
     call summ(y,c)
     call summ(u,d)
  end if

  slope=(nx*a-b*c)/(nx*d-b**2)
  intercept=(c*d-a*b)/(nx*d-b**2)

end subroutine LSM_1d
Subroutine :
x(:,:) :real, intent(in)
: データ要素 1
y(size(x,1),size(x,2)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
undef :real, intent(in), optional
: undef

最小二乗法による傾きと切片計算 (2 次元データ版)

[Source]

subroutine LSM_2d( x, y, slope, intercept, undef )  ! 最小二乗法による傾きと切片計算 (2 次元データ版)
  implicit none
  real, intent(in) :: x(:,:)  ! データ要素 1
  real, intent(in) :: y(size(x,1),size(x,2))  ! データ要素 2
  real, intent(inout) :: slope  ! 最適な傾き
  real, intent(inout) :: intercept  ! 最適な切片
  real, intent(in), optional :: undef ! undef 
  real :: u(size(x,1)*size(x,2)), v(size(x,1)*size(x,2))
  integer :: i, j, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2

  nx=size(x,1)
  ny=size(x,2)

  counter=0

  do j=1,ny
     do i=1,nx
        counter=counter+1
        u(counter)=x(i,j)
        v(counter)=y(i,j)
     end do
  end do

  if(present(undef))then
     call LSM_1d( u, v, slope, intercept, undef )
  else
     call LSM_1d( u, v, slope, intercept )
  end if

end subroutine LSM_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ要素 1
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
undef :real, intent(in), optional
: undef

最小二乗法による傾きと切片計算 (3 次元データ版)

[Source]

subroutine LSM_3d( x, y, slope, intercept, undef )  ! 最小二乗法による傾きと切片計算 (3 次元データ版)
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ要素 1
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! データ要素 2
  real, intent(inout) :: slope  ! 最適な傾き
  real, intent(inout) :: intercept  ! 最適な切片
  real, intent(in), optional :: undef ! undef 
  real :: u(size(x,1)*size(x,2)*size(x,3)), v(size(x,1)*size(x,2)*size(x,3))
  integer :: i, j, k, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2
  integer :: nz  ! データ数 3

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  counter=0

  do k=1,nz
     do j=1,ny
        do i=1,nx
           counter=counter+1
           u(counter)=x(i,j,k)
           v(counter)=y(i,j,k)
        end do
     end do
  end do

  if(present(undef))then
     call LSM_1d( u, v, slope, intercept, undef )
  else
     call LSM_1d( u, v, slope, intercept )
  end if

end subroutine LSM_3d
LSM_poly( x, y, a, intercept, [undef] )
Subroutine :
x(:) :real, intent(in)
: データ要素配列 1
y(size(x)) :real, intent(in)
: データ要素配列 2
a(:) :real, intent(inout)
: 多項式の係数
intercept :real, intent(inout)
: y 切片. a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が あり, 紛らわしいと判断したため, a_0 である y 切片を 独立で引数として渡すことにした.
undef :real, intent(in), optional
: 未定義値.

LSM の多項式近似バージョン. LSM では, F(x)=a_0+a_1x の直線近似を行っていたが, LSM_poly では, F(x)=sum^{N}_{n=0}{a_nx^n} の任意次数の多項式曲線近似を行うことが可能. アルゴリズムは最小二乗法を用いており, 係数のソルバには gausss ルーチンを使用.

Alias for LSM_poly_1d

Subroutine :
x(:) :real, intent(in)
: データ要素配列 1
y(size(x)) :real, intent(in)
: データ要素配列 2
a(:) :real, intent(inout)
: 多項式の係数
intercept :real, intent(inout)
: y 切片. a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が あり, 紛らわしいと判断したため, a_0 である y 切片を 独立で引数として渡すことにした.
undef :real, intent(in), optional
: 未定義値.

LSM の多項式近似バージョン. LSM では, F(x)=a_0+a_1x の直線近似を行っていたが, LSM_poly では, F(x)=sum^{N}_{n=0}{a_nx^n} の任意次数の多項式曲線近似を行うことが可能. アルゴリズムは最小二乗法を用いており, 係数のソルバには gausss ルーチンを使用.

[Source]

subroutine LSM_poly_1d( x, y, a, intercept, undef )
! LSM の多項式近似バージョン.
! LSM では, F(x)=a_0+a_1x の直線近似を行っていたが,
! LSM_poly では, F(x)=\sum^{N}_{n=0}{a_nx^n}
! の任意次数の多項式曲線近似を行うことが可能.
! アルゴリズムは最小二乗法を用いており, 係数のソルバには gausss ルーチンを使用.
  use Matrix_Calc
  implicit none
  real, intent(in) :: x(:)  ! データ要素配列 1
  real, intent(in) :: y(size(x))  ! データ要素配列 2
  real, intent(inout) :: a(:)  ! 多項式の係数
  real, intent(inout) :: intercept  ! y 切片. 
                         ! a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が
                         ! あり, 紛らわしいと判断したため, a_0 である y 切片を
                         ! 独立で引数として渡すことにした.
  real, intent(in), optional :: undef  ! 未定義値.
  integer :: i, j, k
  integer :: nx  ! データの個数
  integer :: poly_n  ! 近似する曲線の最高次数. 1 なら, LSM と同じ.
  real :: coe(0:size(a)), tmpa_coe(0:size(a),0:size(a)), tmpb_coe(0:size(a))
          ! coe は a_n が入る. tmp_coe はデータの総和が入る.
          ! [注意] : 第一要素が行. 第二要素が列.
  real :: tmp(size(x))  ! べき乗計算の一時配列

  nx=size(x)
  poly_n=size(a)

!-- gausss に渡しやすいように, 用意した配列に引数を代入.
  if(present(undef))then
     do k=0,poly_n  ! 列成分の計算
        do j=0,poly_n  ! 行成分の計算. 行成分の計算が先に回ることに注意.
           if(j >= k)then  ! 行成分(j)より列成分(k)の要素数が小さい場合, 値を
                           ! まじめに計算する.
              do i=1,nx
                 if(x(i)/=undef)then
                    tmp(i)=x(i)**(j+k)
                 else
                    tmp(i)=undef
                 end if
              end do
              call summ( tmp, tmpa_coe(j,k), undef )
           else  ! 行成分(j)より列成分(k)の要素数が大きい場合, 解く係数行列が
                 ! 対称行列であることから, 値の参照代入のみ行う.
              tmpa_coe(j,k)=tmpa_coe(k,j)  ! 対称成分の代入(すでに計算済み)
           end if
        end do
     end do
     do j=0,poly_n
        do i=1,nx
           if(x(i)/=undef)then
              tmp(i)=y(i)*(x(i)**j)
           else
              tmp(i)=undef
           end if
        end do
        call summ( tmp, tmpb_coe(j), undef )
     end do
  else  ! undef 処理がないとき.
     do k=0,poly_n  ! 列成分の計算
        do j=0,poly_n  ! 行成分の計算. 行成分の計算が先に回ることに注意.
           if(j >= k)then  ! 行成分(j)より列成分(k)の要素数が小さい場合, 値を
                           ! まじめに計算する.
              do i=1,nx
                 tmp(i)=x(i)**(j+k)
              end do
              call summ( tmp, tmpa_coe(j,k), undef )
           else  ! 行成分(j)より列成分(k)の要素数が大きい場合, 解く係数行列が
                 ! 対称行列であることから, 値の参照代入のみ行う.
              tmpa_coe(j,k)=tmpa_coe(k,j)  ! 対称成分の代入(すでに計算済み)
           end if
        end do
     end do
     do j=0,poly_n
        do i=1,nx
           tmp(i)=y(i)*(x(i)**j)
        end do
        call summ( tmp, tmpb_coe(j), undef )
     end do
  end if

!  以上で係数行列に値が入った.

  call gausss( tmpa_coe(0:poly_n,0:poly_n), tmpb_coe(0:poly_n), coe(0:poly_n) )

  do i=1,poly_n
     a(i)=coe(i)
  end do
  intercept=coe(0)

end subroutine LSM_poly_1d
Subroutine :
x(:,:) :real, intent(in)
: データ要素配列 1
y(size(x,1),size(x,2)) :real, intent(in)
: データ要素配列 2
a(:) :real, intent(inout)
: 多項式の係数
intercept :real, intent(inout)
: y 切片. a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が あり, 紛らわしいと判断したため, a_0 である y 切片を 独立で引数として渡すことにした.
undef :real, intent(in), optional
: 未定義値.

LSM の多項式近似バージョン. (2 次元データ版)

[Source]

subroutine LSM_poly_2d( x, y, a, intercept, undef )
! LSM の多項式近似バージョン. (2 次元データ版)
  use Matrix_Calc
  implicit none
  real, intent(in) :: x(:,:)  ! データ要素配列 1
  real, intent(in) :: y(size(x,1),size(x,2))  ! データ要素配列 2
  real, intent(inout) :: a(:)  ! 多項式の係数
  real, intent(inout) :: intercept  ! y 切片. 
                         ! a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が
                         ! あり, 紛らわしいと判断したため, a_0 である y 切片を
                         ! 独立で引数として渡すことにした.
  real, intent(in), optional :: undef  ! 未定義値.
  integer :: i, j, counter
  integer :: nx  ! データの個数 1
  integer :: ny  ! データの個数 2
  real, dimension(size(x,1)*size(x,2)) :: val1, val2

  nx=size(x,1)
  ny=size(x,2)

  counter=0

  do j=1,ny
     do i=1,nx
        counter=counter+1
        val1(counter)=x(i,j)
        val2(counter)=y(i,j)
     end do
  end do

  if(present(undef))then
     call LSM_poly_1d( val1, val2, a, intercept, undef )
  else
     call LSM_poly_1d( val1, val2, a, intercept )
  end if

end subroutine LSM_poly_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ要素配列 1
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: データ要素配列 2
a(:) :real, intent(inout)
: 多項式の係数
intercept :real, intent(inout)
: y 切片. a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が あり, 紛らわしいと判断したため, a_0 である y 切片を 独立で引数として渡すことにした.
undef :real, intent(in), optional
: 未定義値.

LSM の多項式近似バージョン. (3 次元データ版)

[Source]

subroutine LSM_poly_3d( x, y, a, intercept, undef )
! LSM の多項式近似バージョン. (3 次元データ版)
  use Matrix_Calc
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ要素配列 1
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! データ要素配列 2
  real, intent(inout) :: a(:)  ! 多項式の係数
  real, intent(inout) :: intercept  ! y 切片. 
                         ! a に組み込むと引数を渡すとき, poly_n+1 で渡す必要が
                         ! あり, 紛らわしいと判断したため, a_0 である y 切片を
                         ! 独立で引数として渡すことにした.
  real, intent(in), optional :: undef  ! 未定義値.
  integer :: i, j, k, counter
  integer :: nx  ! データの個数 1
  integer :: ny  ! データの個数 2
  integer :: nz  ! データの個数 3
  real, dimension(size(x,1)*size(x,2)*size(x,3)) :: val1, val2

  nx=size(x,1)
  ny=size(x,2)
  ny=size(x,3)

  counter=0

  do k=1,nz
     do j=1,ny
        do i=1,nx
           counter=counter+1
           val1(counter)=x(i,j,k)
           val2(counter)=y(i,j,k)
        end do
     end do
  end do

  if(present(undef))then
     call LSM_poly_1d( val1, val2, a, intercept, undef )
  else
     call LSM_poly_1d( val1, val2, a, intercept )
  end if

end subroutine LSM_poly_3d
Subroutine :
x(:) :real, intent(in)
: データ
ave :real, intent(inout)
: 計算する平均値
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

1 次元配列平均値計算ルーチン

[Source]

subroutine Mean_1d( x, ave, error )  ! 1 次元配列平均値計算ルーチン
  implicit none
  real, intent(in) :: x(:)  ! データ
  real, intent(inout) :: ave  ! 計算する平均値
  real, intent(in), optional :: error  ! 欠損値が存在するデータセットの場合の欠損値
  integer :: i, nt
  integer :: nx  ! データの要素数
  real :: summ

  summ=0.0
  nt=0
  nx=size(x)

  if(present(error))then
     do i=1,nx
        if(x(i)/=error)then
           summ=summ+x(i)
           nt=1+nt
        end if
     end do

     if(nt/=0)then
        ave=summ/nt
     else
        ave=0.0
     end if

  else

     do i=1,nx
        summ=summ+x(i)
     end do

     ave=summ/nx

  end if

end subroutine Mean_1d
Subroutine :
x(:,:) :real, intent(in)
: データ
ave :real, intent(inout)
: 計算する平均値
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

2 次元配列平均値計算ルーチン

[Source]

subroutine Mean_2d( x, ave, error )  ! 2 次元配列平均値計算ルーチン
  implicit none
  real, intent(in) :: x(:,:)  ! データ
  real, intent(inout) :: ave  ! 計算する平均値
  real, intent(in), optional :: error  ! 欠損値が存在するデータセットの場合の欠損値
  integer :: i, j, nt
  integer :: nx  ! データの要素数 1
  integer :: ny  ! データの要素数 2
  real :: summ

  summ=0.0
  nt=0
  nx=size(x,1)
  ny=size(x,2)

  if(present(error))then
     do j=1,ny
        do i=1,nx
           if(x(i,j)/=error)then
              summ=summ+x(i,j)
              nt=1+nt
           end if
        end do
     end do

     if(nt/=0)then
        ave=summ/nt
     else
        ave=0.0
     end if

  else

     do j=1,ny
        do i=1,nx
           summ=summ+x(i,j)
        end do
     end do

     ave=summ/(nx*ny)

  end if

end subroutine Mean_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ
ave :real, intent(inout)
: 計算する平均値
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

3 次元配列平均値計算ルーチン

[Source]

subroutine Mean_3d( x, ave, error )  ! 3 次元配列平均値計算ルーチン
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ
  real, intent(inout) :: ave  ! 計算する平均値
  real, intent(in), optional :: error  ! 欠損値が存在するデータセットの場合の欠損値
  integer :: i, j, k, nt
  integer :: nx  ! データの要素数 1
  integer :: ny  ! データの要素数 2
  integer :: nz  ! データの要素数 2
  real :: summ

  summ=0.0
  nt=0
  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  if(present(error))then
     do k=1,nz
        do j=1,ny
           do i=1,nx
              if(x(i,j,k)/=error)then
                 summ=summ+x(i,j,k)
                 nt=1+nt
              end if
           end do
        end do
     end do

     if(nt/=0)then
        ave=summ/nt
     else
        ave=0.0
     end if

  else

     do k=1,nz
        do j=1,ny
           do i=1,nx
              summ=summ+x(i,j,k)
           end do
        end do
     end do

     ave=summ/(nx*ny*nz)

  end if

end subroutine Mean_3d
Subroutine :
x(:) :real, intent(in)
: データ
n :integer, intent(in)
: 平均をとる数
y(size(x)) :real, intent(inout)
: 平均化した後のデータ. 実際は, y(1:n-1) までの配列にはゼロが入る.
error :real, intent(in), optional
: 欠損値

移動平均を計算するルーチン

[Source]

subroutine Move_ave( x, n, y, error )
! 移動平均を計算するルーチン
  implicit none
  real, intent(in) :: x(:)  ! データ
  integer, intent(in) :: n  ! 平均をとる数
  real, intent(inout) :: y(size(x))  ! 平均化した後のデータ.
                      ! 実際は, y(1:n-1) までの配列にはゼロが入る.
  real, intent(in), optional :: error  ! 欠損値
  integer :: nx, i
  real :: tmp

  nx=size(x)

  if(nx<n.or.n<2)then
     write(*,*) "### ERROR ### (Move_ave)"
     write(*,*) "x(nx) : nx must be more than n or n must be more than 2."
     write(*,*) "STOP"
     stop
  end if

  if(present(error))then
     call Mean_1d( x(1:n), tmp, error )
     y(1:n-1)=0.0
     y(n)=tmp

     do i=n+1,nx
        y(i)=y(i-1)+(x(i)-x(i-n))/real(n)
     end do
  else
     call Mean_1d( x(1:n), tmp )
     y(1:n-1)=0.0
     y(n)=tmp

     do i=n+1,nx
        y(i)=y(i-1)+(x(i)-x(i-n))/real(n)
     end do
  end if

end subroutine Move_ave
Reg_Line( x, y, slope, intercept, [undef] )
Subroutine :
x(:) :real, intent(in)
: データ要素 1
y(size(x)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
undef :real, intent(in), optional
: 未定義値

LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン

Alias for Reg_Line_1d

Subroutine :
x(:) :real, intent(in)
: データ要素 1
y(size(x)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
undef :real, intent(in), optional
: 未定義値

LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン

[Source]

subroutine Reg_Line_1d( x, y, slope, intercept, undef )
  ! LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン
  implicit none
  real, intent(in) :: x(:)  ! データ要素 1
  real, intent(in) :: y(size(x))  ! データ要素 2
  real, intent(inout) :: slope  ! 最適な傾き
  real, intent(inout) :: intercept  ! 最適な切片
  real, intent(in), optional :: undef  ! 未定義値
  real :: u(size(x)), v(size(x))
  integer :: nx  ! データ数

  nx=size(x)

  if(present(undef))then
     call Anomaly_1d( x, u, undef )
     call Anomaly_1d( y, v, undef )
     call LSM( u, v, slope, intercept, undef )
  else
     call Anomaly_1d( x, u )
     call Anomaly_1d( y, v )
     call LSM( u, v, slope, intercept )
  end if

end subroutine Reg_Line_1d
Subroutine :
x(:,:) :real, intent(in)
: データ要素 1
y(size(x,1),size(x,2)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
error :real, intent(in), optional

LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン (2 次元版)

[Source]

subroutine Reg_Line_2d( x, y, slope, intercept, error )
  ! LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン (2 次元版)
  implicit none
  real, intent(in) :: x(:,:)  ! データ要素 1
  real, intent(in) :: y(size(x,1),size(x,2))  ! データ要素 2
  real, intent(inout) :: slope  ! 最適な傾き
  real, intent(inout) :: intercept  ! 最適な切片
  real, intent(in), optional :: error
  real, dimension(size(x,1)*size(x,2)) :: u, v
  integer :: i, j, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2

  nx=size(x,1)
  ny=size(x,2)

  counter=0

  do j=1,ny
     do i=1,nx
        counter=counter+1
        u(counter)=x(i,j)
        v(counter)=y(i,j)
     end do
  end do

  if(present(error))then
     call Reg_Line_1d( u, v, slope, intercept, error )
  else
     call Reg_Line_1d( u, v, slope, intercept )
  end if

end subroutine Reg_Line_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ要素 1
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片
error :real, intent(in), optional

LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン (3 次元版)

[Source]

subroutine Reg_Line_3d( x, y, slope, intercept, error )
  ! LSM を用いて回帰直線の傾き slope と切片 intercept を計算するルーチン (3 次元版)
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ要素 1
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! データ要素 2
  real, intent(inout) :: slope  ! 最適な傾き
  real, intent(inout) :: intercept  ! 最適な切片
  real, intent(in), optional :: error
  real, dimension(size(x,1)*size(x,2)*size(x,3)) :: u, v
  integer :: i, j, k, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2
  integer :: nz  ! データ数 3

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  counter=0

  do k=1,nz
     do j=1,ny
        do i=1,nx
           counter=counter+1
           u(counter)=x(i,j,k)
           v(counter)=y(i,j,k)
        end do
     end do
  end do

  if(present(error))then
     call Reg_Line_1d( u, v, slope, intercept, error )
  else
     call Reg_Line_1d( u, v, slope, intercept )
  end if

end subroutine Reg_Line_3d
Subroutine :
x(:) :real, intent(in)
: 元座標
r(:) :real, intent(in)
: 内挿座標
u(size(x)) :real, intent(in)
: 元データ
v(size(r)) :real, intent(inout)
: 内挿したデータ
undef :integer, intent(in), optional
: 未定義値

座標 x で定義されているデータ u を 座標 r で定義されるデータ v に自動で内挿する.

[Source]

subroutine auto_interpolation_1d( x, r, u, v, undef )
  ! 座標 x で定義されているデータ u を
  ! 座標 r で定義されるデータ v に自動で内挿する.
  implicit none
  real, intent(in) :: x(:)  ! 元座標
  real, intent(in) :: r(:)  ! 内挿座標
  real, intent(in) :: u(size(x))  ! 元データ
  real, intent(inout) :: v(size(r))  ! 内挿したデータ
  integer, intent(in), optional :: undef  ! 未定義値
  integer :: i, nx, nr, ix, ir
  integer :: defun

  if(present(undef))then
     defun=undef
  else
     defun=0
  end if

  nx=size(x)
  nr=size(r)

  do i=1, nr
     call interpo_search_1d( x, r(i), ir, defun )

     if(ir/=defun.and.ir<nx)then
        call interpolation_1d( x(ir:ir+1), u(ir:ir+1), r(i), v(i) )
     else
        v(i)=real(defun)
     end if
  end do

end subroutine auto_interpolation_1d
Subroutine :
x(:) :real, intent(in)
: 元座標 1
y(:) :real, intent(in)
: 元座標 2
r(:) :real, intent(in)
: 内挿座標 1
q(:) :real, intent(in)
: 内挿座標 2
u(size(x),size(y)) :real, intent(in)
: 元データ
v(size(r),size(q)) :real, intent(inout)
: 内挿したデータ
undef :integer, intent(in), optional
: 未定義値

座標 x, y で定義されているデータ u を 座標 r, q で定義されるデータ v に自動で内挿する.

[Source]

subroutine auto_interpolation_2d( x, y, r, q, u, v, undef )
  ! 座標 x, y で定義されているデータ u を
  ! 座標 r, q で定義されるデータ v に自動で内挿する.
  implicit none
  real, intent(in) :: x(:)  ! 元座標 1
  real, intent(in) :: y(:)  ! 元座標 2
  real, intent(in) :: r(:)  ! 内挿座標 1
  real, intent(in) :: q(:)  ! 内挿座標 2
  real, intent(in) :: u(size(x),size(y))  ! 元データ
  real, intent(inout) :: v(size(r),size(q))  ! 内挿したデータ
  integer, intent(in), optional :: undef  ! 未定義値
  integer :: i, j, nx, ny, nr, nq, ix, iy, ir, iq
  integer :: defun

  if(present(undef))then
     defun=undef
  else
     defun=0
  end if

  nx=size(x)
  ny=size(y)
  nr=size(r)
  nq=size(q)

  do j=1, nq
     do i=1, nr
        call interpo_search_2d( x, y, r(i), q(j), ir, iq, defun )

        if(ir/=defun.and.iq/=defun.and.ir<nx.and.iq<ny)then
           call interpolation_2d( x(ir:ir+1), y(iq:iq+1), u(ir:ir+1,iq:iq+1), (/r(i), q(j)/), v(i,j) )
        else
           v(i,j)=real(defun)
        end if
     end do
  end do

end subroutine auto_interpolation_2d
Subroutine :
x(:) :real, intent(in)
: 元座標 1
y(:) :real, intent(in)
: 元座標 2
z(:) :real, intent(in)
: 元座標 3
r(:) :real, intent(in)
: 内挿座標 1
q(:) :real, intent(in)
: 内挿座標 2
p(:) :real, intent(in)
: 内挿座標 3
u(size(x),size(y),size(z)) :real, intent(in)
: 元データ
v(size(r),size(q),size(p)) :real, intent(inout)
: 内挿したデータ
undef :integer, intent(in), optional
: 未定義値

座標 x, y, z で定義されているデータ u を 座標 r, q, p で定義されるデータ v に自動で内挿する.

[Source]

subroutine auto_interpolation_3d( x, y, z, r, q, p, u, v, undef )
  ! 座標 x, y, z で定義されているデータ u を
  ! 座標 r, q, p で定義されるデータ v に自動で内挿する.
  implicit none
  real, intent(in) :: x(:)  ! 元座標 1
  real, intent(in) :: y(:)  ! 元座標 2
  real, intent(in) :: z(:)  ! 元座標 3
  real, intent(in) :: r(:)  ! 内挿座標 1
  real, intent(in) :: q(:)  ! 内挿座標 2
  real, intent(in) :: p(:)  ! 内挿座標 3
  real, intent(in) :: u(size(x),size(y),size(z))  ! 元データ
  real, intent(inout) :: v(size(r),size(q),size(p))  ! 内挿したデータ
  integer, intent(in), optional :: undef  ! 未定義値
  integer :: i, j, k, nx, ny, nz, nr, nq, np, ix, iy, iz, ir, iq, ip
  integer :: defun

  if(present(undef))then
     defun=undef
  else
     defun=0
  end if

  nx=size(x)
  ny=size(y)
  nz=size(z)
  nr=size(r)
  nq=size(q)
  np=size(p)

  do k=1, np
     do j=1, nq
        do i=1, nr
           call interpo_search_3d( x, y, z, r(i), q(j), p(k), ir, iq, ip, defun )

           if(ir/=defun.and.iq/=defun.and.ip/=defun.and. ir<nx.and.iq<ny.and.ip<nz)then
              call interpolation_3d( x(ir:ir+1), y(iq:iq+1), z(ip:ip+1), u(ir:ir+1,iq:iq+1,ip:ip+1), (/r(i), q(j), p(k)/), v(i,j,k) )
           else
              v(i,j,k)=real(defun)
           end if
        end do
     end do
  end do

end subroutine auto_interpolation_3d
covariance( x, y, cov, [error] )
Subroutine :
x(:) :real, intent(in)
: データ 1
y(size(x)) :real, intent(in)
: データ 2
cov :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

2 つの 1 次元データの共分散を計算 共分散$sigma $の定義は, $$sigma =sum^{nx}_{i=1}{(x-\bar{x})(y-\bar{y})} $$

Alias for covariance_1d

Subroutine :
x(:) :real, intent(in)
: データ 1
y(size(x)) :real, intent(in)
: データ 2
cov :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

2 つの 1 次元データの共分散を計算 共分散$sigma $の定義は, $$sigma =sum^{nx}_{i=1}{(x-\bar{x})(y-\bar{y})} $$

[Source]

subroutine covariance_1d( x, y, cov, error )  ! 2 つの 1 次元データの共分散を計算
  ! 共分散$\sigma $の定義は,
  ! $$\sigma =\sum^{nx}_{i=1}{(x-\bar{x})(y-\bar{y})} $$
  implicit none
  real, intent(in) :: x(:)  ! データ 1
  real, intent(in) :: y(size(x))  ! データ 2
  real, intent(inout) :: cov  ! 標準偏差
  real, intent(in), optional :: error  ! 欠損値
  integer :: i
  integer :: nx  ! データ数
  real :: an1(size(x)), an2(size(x))

  nx=size(x)
  cov=0.0

  if(present(error))then
     call Anomaly_1d( x, an1, error )
     call Anomaly_1d( y, an2, error )
     do i=1,nx
        if(x(i)/=error)then
           cov=cov+an1(i)*an2(i)
        end if
     end do
  else
     call Anomaly_1d( x, an1 )
     call Anomaly_1d( y, an2 )
     do i=1,nx
        cov=cov+an1(i)*an2(i)
     end do
  end if

end subroutine covariance_1d
Subroutine :
x(:,:) :real, intent(in)
: データ 1
y(size(x,1),size(x,2)) :real, intent(in)
: データ 2
cov :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

2 つの 2 次元データの共分散を計算

[Source]

subroutine covariance_2d( x, y, cov, error )  ! 2 つの 2 次元データの共分散を計算
  implicit none
  real, intent(in) :: x(:,:)  ! データ 1
  real, intent(in) :: y(size(x,1),size(x,2))  ! データ 2
  real, intent(inout) :: cov  ! 標準偏差
  real, intent(in), optional :: error  ! 欠損値
  integer :: i, j, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2
  real :: val1(size(x,1)*size(x,2)), val2(size(x,1)*size(x,2))

  nx=size(x,1)
  ny=size(x,2)

  counter=0

  do j=1,ny
     do i=1,nx
        counter=counter+1
        val1(counter)=x(i,j)
        val2(counter)=y(i,j)
     end do
  end do

  cov=0.0

  if(present(error))then
     call covariance_1d( val1, val2, cov, error )
  else
     call covariance_1d( val1, val2, cov )
  end if

end subroutine covariance_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ 1
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: データ 2
cov :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

2 つの 3 次元データの共分散を計算

[Source]

subroutine covariance_3d( x, y, cov, error )  ! 2 つの 3 次元データの共分散を計算
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ 1
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! データ 2
  real, intent(inout) :: cov  ! 標準偏差
  real, intent(in), optional :: error  ! 欠損値
  integer :: i, j, k, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2
  integer :: nz  ! データ数 3
  real :: val1(size(x,1)*size(x,2)*size(x,3)), val2(size(x,1)*size(x,2)*size(x,3))

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  counter=0

  do k=1,nz
     do j=1,ny
        do i=1,nx
           counter=counter+1
           val1(counter)=x(i,j,k)
           val2(counter)=y(i,j,k)
        end do
     end do
  end do

  cov=0.0

  if(present(error))then
     call covariance_1d( val1, val2, cov, error )
  else
     call covariance_1d( val1, val2, cov )
  end if

end subroutine covariance_3d
Subroutine :
x(:) :real, intent(in)
: 漸増配列
point :real, intent(in)
: この点
i :integer, intent(inout)
: point の値を越えない最大の値をもつ要素番号
undeff :integer, intent(in), optional
: 探索範囲の配列要素より小さい値を探索しようとした際, undef を返すが, その undef 値を設定する. default では 0.

漸増配列(要素数が増えるごとに値が大きくなる配列)のなかで, point の前に来る要素番号を出力する.

[Source]

subroutine interpo_search_1d( x, point, i, undeff )
  ! 漸増配列(要素数が増えるごとに値が大きくなる配列)のなかで,
  ! point の前に来る要素番号を出力する.
  implicit none
  real, intent(in) :: x(:)  ! 漸増配列
  real, intent(in) :: point  ! この点
  integer, intent(inout) :: i  ! point の値を越えない最大の値をもつ要素番号
  integer, intent(in), optional :: undeff  ! 探索範囲の配列要素より小さい値を探索しようとした際, undef を返すが, その undef 値を設定する. default では 0.
  integer :: nx, j
  integer :: just

  nx=size(x)
  if(present(undeff))then
     just=int(undeff)
  else
     just=0
  end if

  do j=1,nx
     if(x(1)>point)then
        write(*,*) "****** WARNING ******"
        write(*,*) "searching point was not found."
        write(*,*) "Abort. Exit.!!!"
        i=just
        exit
     end if

     if(present(undeff))then
        if(x(j)/=undeff)then
           if(x(j)<=point)then
              i=j
           else
              exit
           end if
        end if
     else
        if(x(j)<=point)then
           i=j
        else
           exit
        end if
     end if
  end do

end subroutine interpo_search_1d
Subroutine :
x(:) :real, intent(in)
: 漸増配列 x
y(:) :real, intent(in)
: 漸増配列 y
pointx :real, intent(in)
: この点 x
pointy :real, intent(in)
: この点 y
i :integer, intent(inout)
: pointx の値を越えない最大の値をもつ要素番号
j :integer, intent(inout)
: pointy の値を越えない最大の値をもつ要素番号
undeff :integer, intent(in), optional
: 探索範囲の配列要素より小さい値を探索しようとした際, undef を返すが, その undef 値を設定する. default では 0.

漸増配列(要素数が増えるごとに値が大きくなる配列)のなかで, point の前に来る要素番号を出力する.

[Source]

subroutine interpo_search_2d( x, y, pointx, pointy, i, j, undeff )
  ! 漸増配列(要素数が増えるごとに値が大きくなる配列)のなかで,
  ! point の前に来る要素番号を出力する.
  implicit none
  real, intent(in) :: x(:)  ! 漸増配列 x
  real, intent(in) :: y(:)  ! 漸増配列 y
  real, intent(in) :: pointx  ! この点 x
  real, intent(in) :: pointy  ! この点 y
  integer, intent(inout) :: i  ! pointx の値を越えない最大の値をもつ要素番号
  integer, intent(inout) :: j  ! pointy の値を越えない最大の値をもつ要素番号
  integer, intent(in), optional :: undeff  ! 探索範囲の配列要素より小さい値を探索しようとした際, undef を返すが, その undef 値を設定する. default では 0.
  integer :: just

  if(present(undeff))then
     just=int(undeff)
     call interpo_search_1d( x, pointx, i, just )
     call interpo_search_1d( y, pointy, j, just )
  else
     call interpo_search_1d( x, pointx, i )
     call interpo_search_1d( y, pointy, j )
  end if

end subroutine interpo_search_2d
Subroutine :
x(:) :real, intent(in)
: 漸増配列 x
y(:) :real, intent(in)
: 漸増配列 y
z(:) :real, intent(in)
: 漸増配列 z
pointx :real, intent(in)
: この点 x
pointy :real, intent(in)
: この点 y
pointz :real, intent(in)
: この点 z
i :integer, intent(inout)
: pointx の値を越えない最大の値をもつ要素番号
j :integer, intent(inout)
: pointy の値を越えない最大の値をもつ要素番号
k :integer, intent(inout)
: pointz の値を越えない最大の値をもつ要素番号
undeff :integer, intent(in), optional
: 探索範囲の配列要素より小さい値を探索しようとした際, undef を返すが, その undef 値を設定する. default では 0.

漸増配列(要素数が増えるごとに値が大きくなる配列)のなかで, point の前に来る要素番号を出力する.

[Source]

subroutine interpo_search_3d( x, y, z, pointx, pointy, pointz, i, j, k, undeff )
  ! 漸増配列(要素数が増えるごとに値が大きくなる配列)のなかで,
  ! point の前に来る要素番号を出力する.
  implicit none
  real, intent(in) :: x(:)  ! 漸増配列 x
  real, intent(in) :: y(:)  ! 漸増配列 y
  real, intent(in) :: z(:)  ! 漸増配列 z
  real, intent(in) :: pointx  ! この点 x
  real, intent(in) :: pointy  ! この点 y
  real, intent(in) :: pointz  ! この点 z
  integer, intent(inout) :: i  ! pointx の値を越えない最大の値をもつ要素番号
  integer, intent(inout) :: j  ! pointy の値を越えない最大の値をもつ要素番号
  integer, intent(inout) :: k  ! pointz の値を越えない最大の値をもつ要素番号
  integer, intent(in), optional :: undeff  ! 探索範囲の配列要素より小さい値を探索しようとした際, undef を返すが, その undef 値を設定する. default では 0.
  integer :: just

  if(present(undeff))then
     just=int(undeff)
     call interpo_search_1d( x, pointx, i, just )
     call interpo_search_1d( y, pointy, j, just )
     call interpo_search_1d( z, pointz, k, just )
  else
     call interpo_search_1d( x, pointx, i )
     call interpo_search_1d( y, pointy, j )
     call interpo_search_1d( z, pointz, k )
  end if

end subroutine interpo_search_3d
Subroutine :
x(2) :real, intent(in)
: 内挿点の左右端
y(2) :real, intent(in)
: x の点で定義されている値
point :real, intent(in)
: 内挿点
val :real, intent(inout)
: 内挿点での値

1 次の線形内挿ルーチン

[Source]

subroutine interpolation_1d( x, y, point, val )
  ! 1 次の線形内挿ルーチン
  implicit none
  real, intent(in) :: x(2)  ! 内挿点の左右端
  real, intent(in) :: y(2)  ! x の点で定義されている値
  real, intent(in) :: point  ! 内挿点
  real, intent(inout) :: val  ! 内挿点での値
  real :: fd, dt
  real :: tmin
  real :: tmax
  real :: xmin
  real :: xmax

  tmin=x(1)
  tmax=x(2)
  xmin=y(1)
  xmax=y(2)

  dt=point-tmin

  fd=(xmax-xmin)/(tmax-tmin)

  val=xmin+dt*fd
end subroutine interpolation_1d
Subroutine :
x(2) :real, intent(in)
: 内挿の空間点 x 方向の左右端
y(2) :real, intent(in)
: 内挿の空間点 y 方向の左右端
z(2,2) :real, intent(in)
: x, y での各点での値, (i,j) について, i<=x, j<=y
point(2) :real, intent(in)
: 内挿点 point(1)<=x 座標, point(2)<=y 座標
val :real, intent(inout)
: 内挿点での値

2 次の重線形内挿ルーチン 本ルーチンは直線直交座標空間でのみ使用可能.

[Source]

subroutine interpolation_2d( x, y, z, point, val )
  ! 2 次の重線形内挿ルーチン
  ! 本ルーチンは直線直交座標空間でのみ使用可能.
  implicit none
  real, intent(in) :: x(2)  ! 内挿の空間点 x 方向の左右端
  real, intent(in) :: y(2)  ! 内挿の空間点 y 方向の左右端
  real, intent(in) :: z(2,2)  ! x, y での各点での値, (i,j) について, i<=x, j<=y
  real, intent(in) :: point(2)  ! 内挿点 point(1)<=x 座標, point(2)<=y 座標
  real, intent(inout) :: val  ! 内挿点での値
  real :: valx(2)

  ! y(1) での x 方向の内挿点での値
  call interpolation_1d( (/x(1), x(2)/), (/z(1,1), z(2,1)/), point(1), valx(1) )

  ! y(2) での x 方向の内挿点での値
  call interpolation_1d( (/x(1), x(2)/), (/z(1,2), z(2,2)/), point(1), valx(2) )

  ! x の内挿点からの y 方向の内挿点での値(これが求める内挿点)
  call interpolation_1d( (/y(1), y(2)/), (/valx(1), valx(2)/), point(2), val )

end subroutine interpolation_2d
Subroutine :
x(2) :real, intent(in)
: 内挿の空間点 x 方向の左右端
y(2) :real, intent(in)
: 内挿の空間点 y 方向の左右端
z(2) :real, intent(in)
: 内挿の空間点 z 方向の左右端
u(2,2,2) :real, intent(in)
: x, y, z での各点での値, (i,j,k) について, i<=x, j<=y, k<=z
point(3) :real, intent(in)
: 内挿点 point(1)<=x 座標, point(2)<=y 座標, point(3)<=z 座標
val :real, intent(inout)
: 内挿点での値

3 次の重線形内挿ルーチン 本ルーチンは直線直交座標空間でのみ使用可能.

[Source]

subroutine interpolation_3d( x, y, z, u, point, val )
  ! 3 次の重線形内挿ルーチン
  ! 本ルーチンは直線直交座標空間でのみ使用可能.
  implicit none
  real, intent(in) :: x(2)  ! 内挿の空間点 x 方向の左右端
  real, intent(in) :: y(2)  ! 内挿の空間点 y 方向の左右端
  real, intent(in) :: z(2)  ! 内挿の空間点 z 方向の左右端
  real, intent(in) :: u(2,2,2)  ! x, y, z での各点での値, (i,j,k) について, i<=x, j<=y, k<=z
  real, intent(in) :: point(3)  ! 内挿点 point(1)<=x 座標, point(2)<=y 座標, point(3)<=z 座標
  real, intent(inout) :: val  ! 内挿点での値
  real :: valx(2)

  ! z(1) での x-y 平面での重線形内挿の値
  call interpolation_2d( x, y, u(:,:,1), point(1:2), valx(1) )

  ! z(2) での x 方向の内挿点での値
  call interpolation_2d( x, y, u(:,:,2), point(1:2), valx(2) )

  ! z(1) の内挿点からの z 方向の内挿点での値(これが求める内挿点)
  call interpolation_1d( (/z(1), z(2)/), (/valx(1), valx(2)/), point(3), val )

end subroutine interpolation_3d
Subroutine :
x(:) :real, intent(in)
: 漸増配列
point :real, intent(in)
: この点
i :integer, intent(inout)
: point の最近傍地点の要素番号

1 次元最近傍探索ルーチン interpo_search_1d から値を求め, その値と +1 した値の距離を比較して 距離の短い方を選択する.

[Source]

subroutine nearest_search_1d( x, point, i )
  ! 1 次元最近傍探索ルーチン
  ! interpo_search_1d から値を求め, その値と +1 した値の距離を比較して
  ! 距離の短い方を選択する.
  implicit none
  real, intent(in) :: x(:)  ! 漸増配列
  real, intent(in) :: point  ! この点
  integer, intent(inout) :: i  ! point の最近傍地点の要素番号
  real :: tmp1, tmp2
  integer :: j, nx

  nx=size(x)

  call interpo_search_1d( x, point, j )

  if(j==0)then  ! i=1 にしたいので, tmp1 にx(1), tmp2 に x(2) を入れれば, 後の if 文
           ! でうまく処理される.
     tmp1=x(j+1)
     tmp2=x(j+2)
  else
     if(j==nx)then  ! i=nx にしたいので, tmp2 に x(nx), tmp1 に x(nx-1) を入れれば,
            ! 後の if 文でうまく処理される.
        tmp1=x(j)
        tmp2=x(j-1)
     else
        tmp1=x(j)
        tmp2=x(j+1)
     end if
  end if

  if(abs(point-tmp1)>abs(tmp2-point))then
     i=j+1
  else
     i=j
  end if

end subroutine nearest_search_1d
Subroutine :
x(:) :real, intent(in)
: 漸増配列 x
y(:) :real, intent(in)
: 漸増配列 y
pointx :real, intent(in)
: この点 x
pointy :real, intent(in)
: この点 y
i :integer, intent(inout)
: pointx の最近要素番号
j :integer, intent(inout)
: pointy の最近要素番号

2 次元最近傍探索ルーチン nearest_search_1d から値を求める. 本来, 2 次元であるため, 周囲 4 点の最近を計算する必要があるが, ここでは直交直線座標を考えているので, 各軸方向独立で最近点を計算し, どちらも最近の点が求めたい 2 次元の最近点となる.

[Source]

subroutine nearest_search_2d( x, y, pointx, pointy, i, j )
  ! 2 次元最近傍探索ルーチン
  ! nearest_search_1d から値を求める.
  ! 本来, 2 次元であるため, 周囲 4 点の最近を計算する必要があるが,
  ! ここでは直交直線座標を考えているので, 各軸方向独立で最近点を計算し,
  ! どちらも最近の点が求めたい 2 次元の最近点となる.
  implicit none
  real, intent(in) :: x(:)  ! 漸増配列 x
  real, intent(in) :: y(:)  ! 漸増配列 y
  real, intent(in) :: pointx  ! この点 x
  real, intent(in) :: pointy  ! この点 y
  integer, intent(inout) :: i  ! pointx の最近要素番号
  integer, intent(inout) :: j  ! pointy の最近要素番号

  call nearest_search_1d( x, pointx, i )
  call nearest_search_1d( y, pointy, j )

end subroutine nearest_search_2d
Subroutine :
x(:) :real, intent(in)
: 漸増配列 x
y(:) :real, intent(in)
: 漸増配列 y
z(:) :real, intent(in)
: 漸増配列 z
pointx :real, intent(in)
: この点 x
pointy :real, intent(in)
: この点 y
pointz :real, intent(in)
: この点 z
i :integer, intent(inout)
: pointx の最近要素番号
j :integer, intent(inout)
: pointy の最近要素番号
k :integer, intent(inout)
: pointz の最近要素番号

2 次元最近傍探索ルーチン nearest_search_1d から値を求める. 本来, 2 次元であるため, 周囲 4 点の最近を計算する必要があるが, ここでは直交直線座標を考えているので, 各軸方向独立で最近点を計算し, どちらも最近の点が求めたい 2 次元の最近点となる.

[Source]

subroutine nearest_search_3d( x, y, z, pointx, pointy, pointz, i, j, k )
  ! 2 次元最近傍探索ルーチン
  ! nearest_search_1d から値を求める.
  ! 本来, 2 次元であるため, 周囲 4 点の最近を計算する必要があるが,
  ! ここでは直交直線座標を考えているので, 各軸方向独立で最近点を計算し,
  ! どちらも最近の点が求めたい 2 次元の最近点となる.
  implicit none
  real, intent(in) :: x(:)  ! 漸増配列 x
  real, intent(in) :: y(:)  ! 漸増配列 y
  real, intent(in) :: z(:)  ! 漸増配列 z
  real, intent(in) :: pointx  ! この点 x
  real, intent(in) :: pointy  ! この点 y
  real, intent(in) :: pointz  ! この点 z
  integer, intent(inout) :: i  ! pointx の最近要素番号
  integer, intent(inout) :: j  ! pointy の最近要素番号
  integer, intent(inout) :: k  ! pointz の最近要素番号

  call nearest_search_1d( x, pointx, i )
  call nearest_search_1d( y, pointy, j )
  call nearest_search_1d( z, pointz, k )

end subroutine nearest_search_3d
stand_vari( x, anor, [error] )
Subroutine :
x(:) :real, intent(in)
: データ
anor :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

1 次元データの標準偏差を計算 標準偏差$sigma $の定義は, $$sigma =sum^{nx}_{i=1}{epsilon ^2} $$ ただし, $epsilon $は平均値からのずれ$x-\bar{x}$である.

Alias for stand_vari_1d

Subroutine :
x(:) :real, intent(in)
: データ
anor :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

1 次元データの標準偏差を計算 標準偏差$sigma $の定義は, $$sigma =sum^{nx}_{i=1}{epsilon ^2} $$ ただし, $epsilon $は平均値からのずれ$x-\bar{x}$である.

[Source]

subroutine stand_vari_1d( x, anor, error )  ! 1 次元データの標準偏差を計算
  ! 標準偏差$\sigma $の定義は,
  ! $$\sigma =\sum^{nx}_{i=1}{epsilon ^2} $$
  ! ただし, $\epsilon $は平均値からのずれ$x-\bar{x}$である.
  implicit none
  real, intent(in) :: x(:)  ! データ
  real, intent(inout) :: anor  ! 標準偏差
  real, intent(in), optional :: error  ! 欠損値
  integer :: i
  integer :: nx  ! データ数
  real :: an(size(x))

  nx=size(x)
  anor=0.0

  if(present(error))then
     call Anomaly_1d( x, an, error )
     do i=1,nx
        if(x(i)/=error)then
           anor=anor+an(i)**2
        end if
     end do
  else
     call Anomaly_1d( x, an )
     do i=1,nx
        anor=anor+an(i)**2
     end do
  end if

end subroutine stand_vari_1d
Subroutine :
x(:,:) :real, intent(in)
: データ
anor :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

2 次元データの標準偏差を計算

[Source]

subroutine stand_vari_2d( x, anor, error )  ! 2 次元データの標準偏差を計算
  implicit none
  real, intent(in) :: x(:,:)  ! データ
  real, intent(inout) :: anor  ! 標準偏差
  real, intent(in), optional :: error  ! 欠損値
  integer :: i, j, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2
  real :: val(size(x,1)*size(x,2))

  nx=size(x,1)
  ny=size(x,2)

  counter=0
  do j=1,ny
     do i=1,nx
        counter=counter+1
        val(counter)=x(i,j)
     end do
  end do

  anor=0.0

  if(present(error))then
     call stand_vari_1d( val, anor, error )
  else
     call stand_vari_1d( val, anor )
  end if

end subroutine stand_vari_2d
Subroutine :
x(:,:,:) :real, intent(in)
: データ
anor :real, intent(inout)
: 標準偏差
error :real, intent(in), optional
: 欠損値

3 次元データの標準偏差を計算

[Source]

subroutine stand_vari_3d( x, anor, error )  ! 3 次元データの標準偏差を計算
  implicit none
  real, intent(in) :: x(:,:,:)  ! データ
  real, intent(inout) :: anor  ! 標準偏差
  real, intent(in), optional :: error  ! 欠損値
  integer :: i, j, k, counter
  integer :: nx  ! データ数 1
  integer :: ny  ! データ数 2
  integer :: nz  ! データ数 3
  real :: val(size(x,1)*size(x,2)*size(x,3))

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  counter=0
  do k=1,nz
     do j=1,ny
        do i=1,nx
           counter=counter+1
           val(counter)=x(i,j,k)
        end do
     end do
  end do

  anor=0.0

  if(present(error))then
     call stand_vari_1d( val, anor, error )
  else
     call stand_vari_1d( val, anor )
  end if

end subroutine stand_vari_3d