Class Algebra
In: algebra.f90

代数演算を主に行うモジュール 固有値計算, 行列演算, ベクトル演算を主に行う.

Methods

abst   dot_prod   radius   rectangle_int   vec_prod  

Public Instance methods

Subroutine :
x(:,:,:) :real, intent(in)
: x 方向のベクトル成分
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
z(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
dis(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 各点での絶対値ベクトル

3 次元ベクトルの絶対値を計算するルーチン 配列数を調整することにより, 2 次元での計算も可能.

[Source]

subroutine abst(x,y,z,dis)  ! 3 次元ベクトルの絶対値を計算するルーチン
  ! 配列数を調整することにより, 2 次元での計算も可能.
  implicit none
  real, intent(in) :: x(:,:,:)  ! x 方向のベクトル成分
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: z(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(inout) :: dis(size(x,1),size(x,2),size(x,3))  ! 各点での絶対値ベクトル
  integer :: i, j, k, nx, ny, nz

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

!$omp parallel do shared(dis,x,y,z) private(i,j,k)
  do k=1,nz
     do j=1,ny
        do i=1,nx
           dis(i,j,k)=sqrt(x(i,j,k)**2+y(i,j,k)**2+z(i,j,k)**2)
        end do
     end do
  end do
!$omp end parallel do

end subroutine abst
Subroutine :
x(:,:,:) :real, intent(in)
: x 方向のベクトル成分
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
z(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
u(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: x 方向のベクトル成分
v(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
w(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
dot(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 内積
undeff :real, intent(in), optional

2ベクトルの内積計算ルーチン 配列を工夫すると, 2 次元での内積を計算することも可能

[Source]

subroutine dot_prod(x,y,z,u,v,w,dot,undeff)
  ! 2ベクトルの内積計算ルーチン
  ! 配列を工夫すると, 2 次元での内積を計算することも可能
  implicit none
  real, intent(in) :: x(:,:,:)  ! x 方向のベクトル成分
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: z(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(in) :: u(size(x,1),size(x,2),size(x,3))  ! x 方向のベクトル成分
  real, intent(in) :: v(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: w(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(inout) :: dot(size(x,1),size(x,2),size(x,3))  ! 内積
  real, intent(in), optional :: undeff
  integer :: i, j, k, nx, ny, nz

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

  if(present(undeff))then
!$omp parallel do shared(dot,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              if(x(i,j,k)==undeff.or.u(i,j,k)==undeff.or.y(i,j,k)==undeff.or. v(i,j,k)==undeff.or.z(i,j,k)==undeff.or.w(i,j,k)==undeff)then
                 dot(i,j,k)=undeff
              else
                 dot(i,j,k)=x(i,j,k)*u(i,j,k)+y(i,j,k)*v(i,j,k)+z(i,j,k)*w(i,j,k)
              end if
           end do
        end do
     end do
!$omp end parallel do
  else
!$omp parallel do shared(dot,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              dot(i,j,k)=x(i,j,k)*u(i,j,k)+y(i,j,k)*v(i,j,k)+z(i,j,k)*w(i,j,k)
           end do
        end do
     end do
!$omp end parallel do
  end if

end subroutine dot_prod
Subroutine :
xp :real, intent(in)
: 中心位置座標の x 成分
yp :real, intent(in)
: 中心位置座標の y 成分
zp :real, intent(in)
: 中心位置座標の z 成分
x(:) :real, intent(in)
: x 方向の位置座標
y(:) :real, intent(in)
: y 方向の位置座標
z(:) :real, intent(in)
: z 方向の位置座標
rad(size(x),size(y),size(z)) :real, intent(inout)
: 距離

ある位置からの距離を計算するルーチン 配列数を調整することにより, 2 次元での計算も可能.

[Source]

subroutine radius(xp,yp,zp,x,y,z,rad)
  ! ある位置からの距離を計算するルーチン
  ! 配列数を調整することにより, 2 次元での計算も可能.
  implicit none
  real, intent(in) :: xp  ! 中心位置座標の x 成分
  real, intent(in) :: yp  ! 中心位置座標の y 成分
  real, intent(in) :: zp  ! 中心位置座標の z 成分
  real, intent(in) :: x(:)  ! x 方向の位置座標
  real, intent(in) :: y(:)  ! y 方向の位置座標
  real, intent(in) :: z(:)  ! z 方向の位置座標
  real, intent(inout) :: rad(size(x),size(y),size(z))  ! 距離
  integer :: i, j, k, nx, ny, nz

  nx=size(x)
  ny=size(y)
  nz=size(z)

!$omp parallel do shared(rad,x,y,z,xp,yp,zp) private(i,j,k)
  do k=1,nz
     do j=1,ny
        do i=1,nx
           rad(i,j,k)=sqrt((x(i)-xp)**2+(y(j)-yp)**2+(z(k)-zp)**2)
        end do
     end do
  end do
!$omp end parallel do


end subroutine radius
Subroutine :
x(:) :real, intent(in)
: 積分変数
y(size(x)) :real, intent(in)
: 非積分関数
bot :real, intent(in)
: 積分区間左端
top :real, intent(in)
: 積分区間右端
res :real, intent(inout)
: 台形積分の積分値
undeff :real, intent(in), optional

1 次元台形積分 不等間隔でも計算可能であるが, 精度は保証しない.

[Source]

subroutine rectangle_int( x, y, bot, top, res, undeff )  ! 1 次元台形積分
  ! 不等間隔でも計算可能であるが, 精度は保証しない.
  implicit none
  real, intent(in) :: bot  ! 積分区間左端
  real, intent(in) :: top  ! 積分区間右端
  real, intent(in) :: x(:)  ! 積分変数
  real, intent(in) :: y(size(x))  ! 非積分関数
  real, intent(inout) :: res  ! 台形積分の積分値
  real, intent(in), optional :: undeff
  integer :: i, j, k, nx, i_bot, i_top
  real :: dx, y_bot, y_top

  nx=size(x)

  res=0.0

!-- bot < top でなければエラー
  if(bot>top)then
     write(*,*) "#### ERROR ####"
     write(*,*) "integrated interval must be bot < top. STOP"
     stop
  end if

!-- 以下で積分域を決定
!-- 実際には, 積分境界に最近する点
  do i=1,nx
     if(x(i)>bot)then
        if(i==1)then
           write(*,*) "#### WARNING ####"
           write(*,*) "there is NOT the bot in the x(i)."  ! このときは, i_bot=1としておいても, x(1)=bot なので, 余分計算はゼロとなり影響はない.
        end if

        i_bot=i
        exit

     end if
  end do

  do i=1,nx
     if(x(i)<top)then
        if(i==nx)then
           write(*,*) "#### WARNING ####"
           write(*,*) "there is NOT the top in the x(i)."  ! このとき, i_top=nx としても, x(nx)=top なので, 余分計算はゼロとなる.
        end if

        i_top=i

     end if
  end do

!-- 以下で格子に当てはまらない部分を内挿で補完
  y_bot=y(i_bot-1) +((y(i_bot)-y(i_bot-1))/(x(i_bot)-x(i_bot-1)))*(bot-x(i_bot-1))
  y_top=y(i_top) +((y(i_top+1)-y(i_top))/(x(i_top+1)-x(i_top)))*(x(i_top+1)-top)

  if(i_bot<=i_top)then  ! 積分区間内に格子が 1 つ以上あるとき

     if(present(undeff))then
        do i=i_bot,i_top
           if(y(i+1)/=undeff.and.y(i)/=undeff)then
              if(i==i_bot)then
                 if(i_bot<i_top)then  ! 積分区間に 2 格子以上ある場合
                    res=res+0.5*(x(i)-bot)*(y(i)+y_bot) +0.5*(x(i+1)-x(i))*(y(i+1)+y(i))
                        ! 下端の余りと通常の短冊分
                 else
                    if(i_bot==i_top)then
                      ! 積分区間に 1 格子しかない場合, ここで全積分を行って, exit で抜ける.
                       res=res+0.5*(x(i)-bot)*(y(i)+y_bot) +0.5*(top-x(i))*(y_top+y(i))
                       exit
                    end if
                 end if
              else
                 if(i==i_top)then
                    res=res+0.5*(top-x(i))*(y_top+y(i))  ! 上端の余りのみ
                 else
                    res=res+0.5*(x(i+1)-x(i))*(y(i+1)+y(i))  ! 通常の短冊
                 end if
              end if
           end if
        end do
     else
        do i=i_bot,i_top
           if(i==i_bot)then
              if(i_bot<i_top)then  ! 積分区間に 2 格子以上ある場合
                 res=res+0.5*(x(i)-bot)*(y(i)+y_bot) +0.5*(x(i+1)-x(i))*(y(i+1)+y(i))
                     ! 下端の余りと通常の短冊分
              else
                 if(i_bot==i_top)then
                   ! 積分区間に 1 格子しかない場合, ここで全積分を行って, exit で抜ける.
                    res=res+0.5*(x(i)-bot)*(y(i)+y_bot) +0.5*(top-x(i))*(y_top+y(i))
                    exit
                 end if
              end if
           else
              if(i==i_top)then
                 res=res+0.5*(top-x(i))*(y_top+y(i))  ! 上端の余りのみ
              else
                 res=res+0.5*(x(i+1)-x(i))*(y(i+1)+y(i))  ! 通常の短冊
              end if
           end if
        end do
     end if
  else
     if(present(undeff))then
        if(y(i_top)/=undeff.and.y(i_bot)/=undeff)then
           res=0.5*(top-bot)*(y_top+y_bot)
        end if
     else
        res=0.5*(top-bot)*(y_top+y_bot)
     end if
  end if

end subroutine rectangle_int
Subroutine :
x(:,:,:) :real, intent(in)
: x 方向のベクトル成分
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
z(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
u(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: x 方向のベクトル成分
v(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
w(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
vecx(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 外積の x 成分
vecy(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 外積の y 成分
vecz(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 外積の z 成分
undeff :real, intent(in), optional

2ベクトルの外積計算ルーチン 配列の要素数を工夫することで 2 次元外積も計算可能

[Source]

subroutine vec_prod(x,y,z,u,v,w,vecx,vecy,vecz,undeff)
  ! 2ベクトルの外積計算ルーチン
  ! 配列の要素数を工夫することで 2 次元外積も計算可能
  implicit none
  real, intent(in) :: x(:,:,:)  ! x 方向のベクトル成分
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: z(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(in) :: u(size(x,1),size(x,2),size(x,3))  ! x 方向のベクトル成分
  real, intent(in) :: v(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: w(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(inout) :: vecx(size(x,1),size(x,2),size(x,3))  ! 外積の x 成分
  real, intent(inout) :: vecy(size(x,1),size(x,2),size(x,3))  ! 外積の y 成分
  real, intent(inout) :: vecz(size(x,1),size(x,2),size(x,3))  ! 外積の z 成分
  real, intent(in), optional :: undeff
  integer :: i, j, k, nx, ny, nz

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

  if(present(undeff))then
!$omp parallel do shared(vecx,vecy,vecz,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              if(x(i,j,k)==undeff.or.u(i,j,k)==undeff.or.y(i,j,k)==undeff.or. v(i,j,k)==undeff.or.z(i,j,k)==undeff.or.w(i,j,k)==undeff)then
                 vecx(i,j,k)=undeff
                 vecy(i,j,k)=undeff
                 vecz(i,j,k)=undeff
              else
                 vecx(i,j,k)=y(i,j,k)*w(i,j,k)-z(i,j,k)*v(i,j,k)
                 vecy(i,j,k)=z(i,j,k)*u(i,j,k)-x(i,j,k)*w(i,j,k)
                 vecz(i,j,k)=x(i,j,k)*v(i,j,k)-y(i,j,k)*u(i,j,k)
              end if
           end do
        end do
     end do
!$omp end parallel do

  else

!$omp parallel do shared(vecx,vecy,vecz,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              vecx(i,j,k)=y(i,j,k)*w(i,j,k)-z(i,j,k)*v(i,j,k)
              vecy(i,j,k)=z(i,j,k)*u(i,j,k)-x(i,j,k)*w(i,j,k)
              vecz(i,j,k)=x(i,j,k)*v(i,j,k)-y(i,j,k)*u(i,j,k)
           end do
        end do
     end do
!$omp end parallel do
  end if

end subroutine vec_prod