| Class | LAPACK_linear |
| In: |
src/util/lapack_linear.f90
|
| A(N) : | real(8), intent(inout)
| ||
| B(N-1) : | real(8), intent(inout)
| ||
| C(N-1) : | real(8), intent(inout)
| ||
| D(1,N) : | real(8), intent(inout)
|
実 3 項行列の連立 1 次方程式(倍精度)
subroutine lapack_dgtsv(A, B, C, D)
!=== Dependency
!=end
!=== 暗黙の型宣言禁止
implicit none
!=begin
!=== In/Out
real(8), intent(inout) :: A(N) !係数行列
real(8), intent(inout) :: B(N-1) !係数行列
real(8), intent(inout) :: C(N-1) !係数行列
real(8), intent(inout) :: D(1,N) !定数/解行列
!=end
!=== Work
integer :: NRHS
integer :: INFO
integer :: LDB
!=== 変数の初期化
NRHS = 1
INFO = 0
LDB = N
!=== 解行列の計算. LAPACK を使用.
call DGTSV(N, NRHS, C, A, B, D, LDB, INFO)
!=== 解のコンディションをチェック.
if (INFO /= 0) then
call MessageNotify("Error", "lapack_linear", "INFO is not 0")
stop
end if
end subroutine lapack_dgtsv
| M : | integer, intent(in)
|
subroutine lapack_dgtsv_init(M)
!=== Input
integer, intent(in) :: M !配列の大きさ
!=end
N = M
end subroutine lapack_dgtsv_init
| A(N) : | real(8), intent(inout)
| ||
| B(N-1) : | real(8), intent(inout)
| ||
| C(N-1) : | real(8), intent(inout)
| ||
| B2(N-2) : | real(8), intent(out)
| ||
| IPIV(N) : | real(8), intent(out)
|
実 3 項行列の LU 分解(倍精度)
subroutine lapack_dgttrf(A, B, C, B2, IPIV)
!=== Dependency
!=end
!=== 暗黙の型宣言禁止
implicit none
!=begin
!=== In/Out
real(8), intent(inout) :: A(N) !係数行列
real(8), intent(inout) :: B(N-1) !係数行列
real(8), intent(inout) :: C(N-1) !係数行列
real(8), intent(out) :: B2(N-2) !係数行列
real(8), intent(out) :: IPIV(N) !部分ピボット交換の情報を格納
!=end
!=== Work
integer :: INFO
!=== 変数の初期化
INFO = 0
!=== 解行列の計算. LAPACK を使用.
call DGTTRF(N, C, A, B, B2, IPIV, INFO)
!=== 解のコンディションをチェック.
if (INFO /= 0) then
call MessageNotify("Error", "lapack_linear", "INFO is not 0")
stop
end if
end subroutine lapack_dgttrf
| A(N) : | real(8), intent(in)
| ||
| B(N-1) : | real(8), intent(in)
| ||
| C(N-1) : | real(8), intent(in)
| ||
| D(1,N) : | real(8), intent(inout)
| ||
| B2(N-2) : | real(8), intent(in)
| ||
| IPIV(N) : | real(8), intent(in)
|
subroutine lapack_dgttrs(A, B, C, D, B2, IPIV)
!=== Dependency
!=end
!=== 暗黙の型宣言禁止
implicit none
!=begin
!=== In/Out
real(8), intent(in) :: A(N) !係数行列
real(8), intent(in) :: B(N-1) !係数行列
real(8), intent(in) :: C(N-1) !係数行列
real(8), intent(in) :: B2(N-2) !係数行列
real(8), intent(inout) :: D(1,N) !定数/解行列
real(8), intent(in) :: IPIV(N) !部分ピボット交換の情報を格納
!=end
!=== Work
integer :: NRHS
integer :: INFO
integer :: LDB
character(1),parameter :: TRANS = 'N'
!=== 変数の初期化
NRHS = 0
INFO = 0
LDB = N
!=== 解行列の計算. LAPACK を使用.
call DGTTRS(TRANS, N, NRHS, C, A, B, B2, IPIV, D, LDB, INFO)
!=== 解のコンディションをチェック.
if (INFO /= 0) then
call MessageNotify("Error", "lapack_linear", "INFO is not 0")
stop
end if
end subroutine lapack_dgttrs