Class | Algebra |
In: |
algebra.f90
|
代数演算を主に行うモジュール
Subroutine : | |||
nx : | integer, intent(in)
| ||
a(nx,nx) : | real, intent(inout)
| ||
b(nx) : | real, intent(inout)
| ||
eps : | real, intent(in)
| ||
x(nx) : | real, intent(inout)
|
ガウスザイデル法による連立 1 次方程式ソルバ
subroutine Gau_Sei(nx, a, b, eps, x) ! ガウスザイデル法による連立 1 次方程式ソルバ implicit none integer, intent(in) :: nx ! 第 1 成分の要素数 real, intent(inout) :: a(nx,nx) ! 係数行列 (第 1 要素が行列の行成分を表す) real, intent(inout) :: b(nx) ! ax=b のベクトル real, intent(in) :: eps ! 収束条件 real, intent(inout) :: x(nx) ! 解く解 integer :: i, j, k, l, m, n ! イテレーション用添字 real :: xn ! 更新した x(i) のテンプ領域 real :: err, err_max ! 誤差 !-- 初期値は 0.0 からスタートする --- x=0.0 !-- ピボッティング call Pivot_part( nx, a, b ) !-- 以下, while を使用するため, 1 回目のイテレートは単独で行う --- err_max=0.0 do i=1,nx xn=0.0 do j=1,nx if(j/=i)then xn=xn+a(i,j)*x(j) end if end do xn=(b(i)-xn)/a(i,i) err=errata(x(i),xn,1) write(*,*) "err_max", x(i), nx, err_max,err if(err_max<=err)then err_max=err end if x(i)=xn end do !-- 以下より, 収束条件を満たすまでループする --- do while(err_max>=eps) err_max=0.0 do i=1,nx xn=0.0 do j=1,nx if(j/=i)then xn=xn+a(i,j)*x(j) end if end do xn=(b(i)-xn)/a(i,i) err=errata(x(i),xn,1) if(err_max<=err)then err_max=err end if x(i)=xn end do end do end subroutine Gau_Sei
Subroutine : | |||
nx : | integer, intent(in)
| ||
a(nx,nx) : | real, intent(inout)
| ||
b(nx) : | real, intent(inout)
| ||
eps : | real, intent(in)
| ||
x(nx) : | real, intent(inout)
|
ヤコビ法による連立 1 次方程式ソルバ
subroutine Jacobi_algebra(nx, a, b, eps, x) ! ヤコビ法による連立 1 次方程式ソルバ implicit none integer, intent(in) :: nx ! 第 1 成分の要素数 real, intent(inout) :: a(nx,nx) ! 係数行列 (第 1 要素が行列の行を表す) real, intent(inout) :: b(nx) ! ax=b のベクトル real, intent(in) :: eps ! 収束条件 real, intent(inout) :: x(nx) ! 解く解 real :: y(nx) ! ヤコビ法で使用する一時格納用配列, この配列で一斉更新する.(xn の代わり) integer :: i, j, k, l, m, n ! イテレーション用添字 real :: err, err_max ! 誤差 !-- 初期値は 0,0 からスタートする --- x=0.0 y=0.0 !-- ピボッティング call Pivot_part( nx, a, b ) !-- 以下, 実際のソルバ(while を使用するため, 1 回目のイテレートは単独で行う) --- err_max=0.0 do i=1,nx y(i)=0.0 do j=1,nx if(j/=i)then y(i)=y(i)+a(i,j)*x(j) end if end do y(i)=(b(i)-y(i))/a(i,i) err=errata(x(i),y(i),1) write(*,*) "err_max", x(i), nx, err_max,err if(err_max<=err)then err_max=err end if end do do i=1,nx ! データの一斉更新 x(i)=y(i) end do !-- 以下より, 収束条件を満たすまでループする --- do while(err_max>=eps) err_max=0.0 do i=1,nx y(i)=0.0 do j=1,nx if(j/=i)then y(i)=y(i)+a(i,j)*x(j) end if end do y(i)=(b(i)-y(i))/a(i,i) err=errata(x(i),y(i),1) if(err_max<=err)then err_max=err end if end do do i=1,nx ! データの一斉更新 x(i)=y(i) end do end do end subroutine Jacobi_algebra
Subroutine : | |||
nmax : | integer, intent(in)
| ||
a(nmax,nmax) : | real, intent(inout)
| ||
b(nmax) : | real, intent(inout)
| ||
x(nmax) : | real, intent(inout)
| ||
itermax : | integer |
— LU 分解を計算するサブルーチン —
subroutine LU_devs( nmax, a, b, x, itermax ) !-- LU 分解を計算するサブルーチン --- implicit none integer, intent(in) :: nmax ! 行列要素数 real, intent(inout) :: a(nmax,nmax) ! 係数行列 (第 1 要素が行を表現) real, intent(inout) :: b(nmax) ! 右辺のベクトル real, intent(inout) :: x(nmax) ! 解 real :: d(nmax,nmax), r(nmax), y(nmax) integer :: ip(nmax) real :: scale(nmax), dx(nmax) real :: s, t, pivot, xnorm, dxnorm real :: s1, s2, s3, s4, s5, t0, t1, t2, t3, t4, eps integer :: iter, itermax integer :: p, itemp, i, j, k !-- 反復改良での精度の設定 --- t4=6.0 !-- 配列 x(i) の初期化 --- do i=1,nmax x(i)=0.0 end do do i=1,nmax do j=1,nmax d(i,j)=a(i,j) end do ip(i)=i !-- 最大値を計算するループ --- s=d(i,1) do j=2,nmax if(d(i,j).gt.s)then s=d(i,j) end if end do scale(i)=1.0/s end do do k=1,nmax t=d(ip(k),k)*scale(ip(k)) p=k do i=k,nmax t0=d(ip(i),k)*scale(ip(i)) if(t0.gt.t)then t=t0 p=i end if end do !-- ip(p) と ip(k) の入れ替え --- if(p.ne.k)then itemp=ip(p) ip(p)=ip(k) ip(k)=itemp end if pivot=d(ip(k),k) do i=k+1,nmax d(ip(i),k)=d(ip(i),k)/pivot do j=k+1,nmax d(ip(i),j)=d(ip(i),j)-d(ip(i),k)*d(ip(k),j) end do end do if(k.ge.nmax-1)then exit end if end do !-- 前進消去 --- y(1)=b(ip(1)) do i=2,nmax s1=0.0 do j=1,i-1 s1=s1+d(ip(i),j)*y(j) end do y(i)=b(ip(i))-s1 end do !-- 後退代入 --- x(nmax)=y(nmax)/d(ip(nmax),nmax) do i=nmax-1,1,-1 s2=0.0 do j=i+1,nmax s2=s2+d(ip(i),j)*y(j) end do x(i)=(y(i)-s2)/d(ip(i),i) end do t1=x(1) xnorm=x(1) do i=2,nmax if(x(i).gt.t1)then t1=x(i) xnorm=x(i) end if end do !-- 反復改良 --- eps=10**(-t4) ! 標準精度を 10 進数の t4 桁とする do iter=1,itermax if(xnorm==0.0)then exit end if !-- 残差の計算 --- do i=1,nmax s3=0.0 do j=1,nmax s3=s3+a(i,j)*x(j) end do r(i)=b(i)-s3 end do !-- 前進消去 --- y(1)=r(ip(1)) do i=2,nmax s4=0.0 do j=1,i-1 s4=s4+d(ip(i),j)*y(j) end do y(i)=r(ip(i))-s4 end do !-- 後退代入 --- dx(nmax)=y(nmax)/d(ip(nmax),nmax) do i=nmax-1,1,-1 s5=0.0 do j=i+1,nmax s5=s5+d(ip(i),j)*y(j) end do dx(i)=(y(i)-s5)/d(ip(i),i) end do do i=1,nmax x(i)=x(i)+dx(i) end do t3=dx(1) dxnorm=dx(1) do i=1,nmax if(dx(i).gt.t3)then t3=dx(i) dxnorm=dx(i) end if end do if(dxnorm/xnorm.le.eps)then exit end if end do end subroutine LU_devs
Subroutine : | |||
nx : | integer, intent(in)
| ||
a(nx,nx) : | real, intent(inout)
| ||
b(nx) : | real, intent(inout)
| ||
eps : | real, intent(in)
| ||
accel : | real, intent(in)
| ||
x(nx) : | real, intent(inout)
|
ガウスザイデル法かつ, SOR で加速による連立 1 次方程式ソルバ
subroutine SOR_Gau_Sei(nx, a, b, eps, accel, x) ! ガウスザイデル法かつ, SOR で加速による連立 1 次方程式ソルバ implicit none integer, intent(in) :: nx ! 第 1 成分の要素数 real, intent(inout) :: a(nx,nx) ! 係数行列 (第 1 要素が行列の行を表す) real, intent(inout) :: b(nx) ! ax=b のベクトル real, intent(in) :: eps ! 収束条件 real, intent(in) :: accel ! 加速係数. ただし, 数学的に accel >= 2 では発散するので, ! この値以上が設定されるとエラーで止める. real, intent(inout) :: x(nx) ! 解く解 integer :: i, j, k, l, m, n ! イテレーション用添字 real :: xn ! 更新した x(i) のテンプ領域 real :: err, err_max ! 誤差 !-- 加速パラメータの確認 if(accel>=2.0)then write(*,*) "***** ERROR *****" write(*,*) "accel parameter must be less than 2.0. STOP." stop end if !-- 初期値は 0.0 からスタートする --- x=0.0 !-- ピボッティング call Pivot_part( nx, a, b ) !-- 以下, while を使用するため, 1 回目のイテレートは単独で行う --- err_max=0.0 do i=1,nx xn=0.0 do j=1,nx if(j/=i)then xn=xn+a(i,j)*x(j) end if end do xn=(b(i)-xn)/a(i,i) xn=x(i)+accel*(xn-x(i)) err=errata(x(i),xn,1) write(*,*) "err_max", x(i), nx, err_max,err if(err_max<=err)then err_max=err end if x(i)=xn end do !-- 以下より, 収束条件を満たすまでループする --- do while(err_max>=eps) err_max=0.0 do i=1,nx xn=0.0 do j=1,nx if(j/=i)then xn=xn+a(i,j)*x(j) end if end do xn=(b(i)-xn)/a(i,i) xn=x(i)+accel*(xn-x(i)) err=errata(x(i),xn,1) if(err_max<=err)then err_max=err end if x(i)=xn end do end do end subroutine SOR_Gau_Sei
Subroutine : | |||
nx : | integer, intent(in)
| ||
a(nx,nx) : | real, intent(inout)
| ||
b(nx) : | real, intent(inout)
| ||
eps : | real, intent(in)
| ||
accel : | real, intent(in)
| ||
x(nx) : | real, intent(inout)
|
ヤコビ法かつ SOR 加速による連立 1 次方程式ソルバ
subroutine SOR_Jacobi_algebra(nx, a, b, eps, accel, x) ! ヤコビ法かつ SOR 加速による連立 1 次方程式ソルバ implicit none integer, intent(in) :: nx ! 第 1 成分の要素数 real, intent(inout) :: a(nx,nx) ! 係数行列 (第 1 要素が行列の行を表す) real, intent(inout) :: b(nx) ! ax=b のベクトル real, intent(in) :: eps ! 収束条件 real, intent(in) :: accel ! 加速係数. ただし, 数学的に accel >= 2 では発散するので, ! この値以上が設定されるとエラーで止める. real, intent(inout) :: x(nx) ! 解く解 real :: y(nx) ! ヤコビ法で使用する一時格納用配列, この配列で一斉更新する.(xn の代わり) integer :: i, j, k, l, m, n ! イテレーション用添字 real :: err, err_max ! 誤差 !-- 加速パラメータの確認 if(accel>=2.0)then write(*,*) "***** ERROR *****" write(*,*) "accel parameter must be less than 2.0. STOP." stop end if !-- 初期値は 0,0 からスタートする --- x=0.0 y=0.0 !-- ピボッティング call Pivot_part( nx, a, b ) !-- 以下, 実際のソルバ(while を使用するため, 1 回目のイテレートは単独で行う) --- err_max=0.0 do i=1,nx y(i)=0.0 do j=1,nx if(j/=i)then y(i)=y(i)+a(i,j)*x(j) end if end do y(i)=(b(i)-y(i))/a(i,i) err=errata(x(i),y(i),1) write(*,*) "err_max", x(i), nx, err_max,err if(err_max<=err)then err_max=err end if end do do i=1,nx ! データの一斉更新 x(i)=x(i)+accel*(y(i)-x(i)) end do !-- 以下より, 収束条件を満たすまでループする --- do while(err_max>=eps) err_max=0.0 do i=1,nx y(i)=0.0 do j=1,nx if(j/=i)then y(i)=y(i)+a(i,j)*x(j) end if end do y(i)=(b(i)-y(i))/a(i,i) y(i)=x(i)+accel*(y(i)-x(i)) err=errata(x(i),y(i),1) if(err_max<=err)then err_max=err end if end do do i=1,nx ! データの一斉更新 x(i)=y(i) end do end do end subroutine SOR_Jacobi_algebra
Subroutine : | |||
nx : | integer, intent(in)
| ||
a(nx,nx) : | real, intent(in)
| ||
eps : | real, intent(in)
| ||
eigenval : | real, intent(inout)
| ||
eigenvec(nx) : | real, intent(inout)
|
べき乗法を用いて行列の最大固有値とその固有値に対応する固有ベクトルを求める.
subroutine eigenvalue_power( nx, a, eps, eigenval, eigenvec ) ! べき乗法を用いて行列の最大固有値とその固有値に対応する固有ベクトルを求める. implicit none integer, intent(in) :: nx ! 行列の要素数 real, intent(in) :: a(nx,nx) ! 固有値を求める行列 (第 1 要素が行を表す) real, intent(in) :: eps ! 収束判定条件 real, intent(inout) :: eigenval ! 行列 a の最大固有値 real, intent(inout) :: eigenvec(nx) ! 固有ベクトル integer :: i, j, k real, dimension(nx) :: x, y real :: tmp1, tmp2, err_max, err do i=1,nx x(i)=1.0 ! 反復法の初期値として非ゼロのベクトルを定義 end do tmp1=sqrt(vec_dot( x, x )) ! 初期ベクトルのノルムを計算 do i=1,nx x(i)=x(i)/tmp1 ! 初期ベクトルを規格化 end do err_max=eps ! while 文に入れるための処理 !-- 反復法開始 do while(err_max>=eps) err_max=0.0 do i=1,nx y(i)=0.0 ! 配列の初期化 do j=1,nx y(i)=y(i)+a(i,j)*x(j) end do end do tmp1=sqrt(vec_dot( y, y )) ! ベクトルのノルムを計算 do i=1,nx x(i)=y(i)/tmp1 ! x(i) の更新 end do !-- 固有値計算 tmp2=vec_dot( x, y ) ! 上で計算した Ay に y^t (つまり, 固有ベクトルの転置) をかける. err_max=errata( eigenval, tmp2, 1 ) ! 過去の x(i) と更新した y(i) の誤差比較 eigenval=tmp2 end do !-- 反復法終了 do i=1,nx eigenvec(i)=x(i) ! 固有ベクトルの変数へ代入 end do end subroutine eigenvalue_power
Subroutine : | |||
nmax : | integer, intent(in)
| ||
a(nmax,nmax) : | real, intent(inout)
| ||
b(nmax) : | real, intent(inout) | ||
x(nmax) : | real, intent(inout) |
部分ピボット付きガウスの消去法
subroutine gausss( nmax, a, b, x ) ! 部分ピボット付きガウスの消去法 implicit none integer, intent(in) :: nmax ! 配列の要素数 real, intent(inout) :: a(nmax,nmax) ! 係数行列 (第 1 要素が行を表す) real, intent(inout) :: b(nmax) real, intent(inout) :: x(nmax) real :: s, pivotb real :: pivot(nmax+1) integer :: piv, i, j, k !-- 前進消去 --- !-- A(I,J) の前進消去 --- do k=1,nmax-1 !-- PIVOT の選択 --- !-- まず, I 成分の最大値を決定する --- piv=k do i=k+1,nmax if(abs(a(i,k)).gt.abs(a(piv,k)))then piv=i end if end do !-- ここまでで, 最大値が決定される --- !-- 以下で, 最大値をとる成分の行を入れ替える --- do j=k,nmax pivot(j)=a(k,j) a(k,j)=a(piv,j) a(piv,j)=pivot(j) end do pivotb=b(k) b(k)=b(piv) b(piv)=pivotb !-- PIVOT 処理ここまで --- do i=k+1,nmax a(k,i)=a(k,i)/a(k,k) end do b(k)=b(k)/a(k,k) a(k,k)=1.0 do j=k+1,nmax do i=k+1,nmax a(j,i)=a(j,i)-a(k,i)*a(j,k) end do b(j)=b(j)-b(k)*a(j,k) a(j,k)=0.0 end do end do b(nmax)=b(nmax)/a(nmax,nmax) a(nmax,nmax)=1.0 !-- X(I) の後進代入 x(nmax)=b(nmax) do i=nmax-1,1,-1 s=b(i) do j=i+1,nmax s=s-a(i,j)*x(j) end do x(i)=s end do end subroutine gausss
Subroutine : | |||
nx : | integer, intent(in)
| ||
ax(nx,nx) : | real, intent(in)
| ||
xx(nx,nx) : | real, intent(inout)
|
ガウスの消去法を拡張して, 逆行列を計算する. 具体的には, 右辺の b に単位ベクトルを 1 つずつ代入し,消去した結果のベクトルを並べて 逆行列にする.
subroutine invert_mat( nx, ax, xx ) ! ガウスの消去法を拡張して, 逆行列を計算する. ! 具体的には, 右辺の b に単位ベクトルを 1 つずつ代入し,消去した結果のベクトルを並べて ! 逆行列にする. implicit none integer, intent(in) :: nx ! 配列の要素数 real, intent(in) :: ax(nx,nx) ! 逆行列を求める行列 real, intent(inout) :: xx(nx,nx) ! 求めた逆行列 integer :: i, j, k real :: c(nx,nx) ! 単位行列 real :: d(nx,nx) ! a(i,j) をガウスルーチンに渡すと, 結果が変化するのでこれに一時退避 c=0.0 do i=1,nx c(i,i)=1.0 end do do i=1,nx do j=1,nx do k=1,nx d(k,j)=ax(k,j) ! ダミー変数へ代入 end do end do call gausss( nx, d, c(:,i), xx(:,i) ) ! 配列第 2 成分が列を表すので, この順番で OK. end do end subroutine invert_mat
Subroutine : | |||
x(:) : | real, intent(in)
| ||
y(size(x)) : | real, intent(in)
| ||
bot : | real, intent(in)
| ||
top : | real, intent(in)
| ||
res : | real, intent(inout)
|
1 次元台形積分 不等間隔でも計算可能であるが, 精度は保証しない.
subroutine rectangle_int( x, y, bot, top, res ) ! 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 ! 台形積分の積分値 integer :: i, j, k, nx real :: dx nx=size(x) res=0.0 do i=1,nx-1 if(x(i)>=bot.and.x(i)<=top)then ! 積分開始位置を決定 res=res+0.5*(x(i+1)-x(i))*(y(i+1)+y(i)) end if end do end subroutine rectangle_int
Subroutine : | |||
u(:,:) : | real, intent(in)
| ||
v(size(u,1),size(u,2)) : | real, intent(inout)
|
シュミットの直交化法を用いて, nx 次元ベクトルを正規直交化する. 引数の配列は第 1 要素がベクトルの各成分, 第 2 成分がベクトル群の 1 ベクトルを表す. つまり, u(i,j) は j 番目ベクトルの第 i 成分ということを表す. 行列で表現するなら, 縦ベクトルを横に並べた形と等しい.
subroutine schumit_norm( u, v ) ! シュミットの直交化法を用いて, nx 次元ベクトルを正規直交化する. ! 引数の配列は第 1 要素がベクトルの各成分, 第 2 成分がベクトル群の 1 ベクトルを表す. ! つまり, u(i,j) は j 番目ベクトルの第 i 成分ということを表す. ! 行列で表現するなら, 縦ベクトルを横に並べた形と等しい. implicit none real, intent(in) :: u(:,:) ! 直交化する前のベクトル real, intent(inout) :: v(size(u,1),size(u,2)) ! 直交化後のベクトル integer :: i, j, k, nx, ny real :: tmp real :: tmpn(size(u,2)) ! 正規化を行うときの, ノルムの値を格納する. real :: tmps(size(u,1)) ! 総和演算の際の一時格納に使用. nx=size(u,1) ny=size(u,2) tmpn(1)=sqrt(vec_dot( u(:,1), u(:,1) )) do i=1,nx ! 1 本目のベクトルを基底の基準とする. v(i,1)=u(i,1)/tmpn(1) end do do j=2,ny do i=1,nx tmps(i)=0.0 end do do k=1,j-1 do i=1,nx tmps(i)=tmps(i)+vec_dot( u(:,j), v(:,k) )*u(i,j) end do end do do i=1,nx v(i,j)=u(i,j)-tmps(i) end do ! 以下で正規化を行う tmpn(j)=sqrt(vec_dot( v(:,j), v(:,j) )) do i=1,nx v(i,j)=v(i,j)/tmpn(j) end do end do end subroutine schumit_norm