Class Statistics
In: statistics.f90

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

Methods

Public Instance methods

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

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

[Source]

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

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

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

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

[Source]

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

  if(present(error))then
     call Mean_2d( nx, ny, x, ave, error )
     do j=1,ny
        do i=1,nx
           if(x(i,j)==error)then
              anor(i,j)=0.0
           else
              anor(i,j)=x(i,j)-ave
           end if
        end do
     end do
  else
     call Mean_2d( nx, ny, 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 :
nx :integer, intent(in)
: データの要素数 1
ny :integer, intent(in)
: データの要素数 2
nz :integer, intent(in)
: データの要素数 3
x(nx,ny,nz) :real, intent(in)
: データ
anor(nx,ny,nz) :real, intent(inout)
: 各 x(i,j,k) に対応する偏差 anor(i,j,k)
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

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

[Source]

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

  if(present(error))then
     call Mean_3d( nx, ny, nz, 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)=0.0
              else
                 anor(i,j,k)=x(i,j,k)-ave
              end if
           end do
        end do
     end do
  else
     call Mean_3d( nx, ny, nz, 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
Subroutine :
nx :integer, intent(in)
: データ個数
x(nx) :real, intent(in)
: データ要素 1
y(nx) :real, intent(in)
: データ要素 2
cc :real, intent(inout)
: 相関係数
error :real, intent(in), optional
: 欠損値

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

[Source]

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

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

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

end subroutine Cor_Coe
Subroutine :
nx :integer, intent(in)
: データ数
nx :integer, intent(in)
x(nx) :real, intent(in)
: データ要素 1
y(nx) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片

最小二乗法による傾きと切片計算

[Source]

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

  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

  call summ(nx,v,a)
  call summ(nx,x,b)
  call summ(nx,y,c)
  call summ(nx,u,d)

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

  contains
  subroutine summ(nx,z,add)
    implicit none
    integer, intent(in) :: nx
    real, intent(in) :: z(nx)
    real, intent(inout) :: add
    integer :: i

    add=0.0
    do i=1,nx
       add=add+z(i)
    end do

  end subroutine summ

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

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

[Source]

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

  summ=0.0
  nt=0

  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 :
nx :integer, intent(in)
: データの要素数 1
ny :integer, intent(in)
: データの要素数 2
x(nx,ny) :real, intent(in)
: データ
ave :real, intent(inout)
: 計算する平均値
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

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

[Source]

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

  summ=0.0
  nt=0

  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 :
nx :integer, intent(in)
: データの要素数 1
ny :integer, intent(in)
: データの要素数 2
nz :integer, intent(in)
: データの要素数 2
x(nx,ny,nz) :real, intent(in)
: データ
ave :real, intent(inout)
: 計算する平均値
error :real, intent(in), optional
: 欠損値が存在するデータセットの場合の欠損値

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

[Source]

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

  summ=0.0
  nt=0

  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 :
nx :integer, intent(in)
: データ数
x(nx) :real, intent(in)
: データ要素 1
y(nx) :real, intent(in)
: データ要素 2
slope :real, intent(inout)
: 最適な傾き
intercept :real, intent(inout)
: 最適な切片

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

[Source]

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

  call Anomaly_1d( nx, x, u )
  call Anomaly_1d( nx, y, v )
  call LSM( nx, u, v, slope, intercept )

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

  cov=0.0

  if(present(error))then
     call Anomaly_1d( nx, x, an1, error )
     call Anomaly_1d( nx, 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( nx, x, an1 )
     call Anomaly_1d( nx, y, an2 )
     do i=1,nx
        cov=cov+an1(i)*an2(i)
     end do
  end if

end subroutine covariance
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=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(x(j)<=point)then
        i=j
     else
        exit
     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=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=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 :
tmin :real, intent(in)
: 内挿点の左端
tmax :real, intent(in)
: 内挿点の右端
xmin :real, intent(in)
: tmin での値
xmax :real, intent(in)
: tmax での値
point :real, intent(in)
: 内挿点
val :real, intent(inout)
: 内挿点での値

1 次の線形内挿ルーチン

[Source]

subroutine interpolation_1d( tmin, tmax, xmin, xmax, point, val )
  ! 1 次の線形内挿ルーチン
  implicit none
  real, intent(in) :: tmin  ! 内挿点の左端
  real, intent(in) :: tmax  ! 内挿点の右端
  real, intent(in) :: xmin  ! tmin での値
  real, intent(in) :: xmax  ! tmax での値
  real, intent(in) :: point  ! 内挿点
  real, intent(inout) :: val  ! 内挿点での値
  real :: fd, dt

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

  anor=0.0

  if(present(error))then
     call Anomaly_1d( nx, 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( nx, x, an )
     do i=1,nx
        anor=anor+an(i)**2
     end do
  end if

end subroutine stand_vari