| Path: | src/lumatrix.f90 |
| Last Update: | Thu Jan 20 17:02:01 JST 2005 |
Copyright (c) 2002-2005 Shin-ichi Takehiro. All rights reserved.
表題 lumatrix : 行列の LU 分解による線形連立方程式の解
履歴 2002/01/20 竹広真一
2002/06/10 竹広真一 ベクトル長問題回避のため lusol2 を使用
2005/01/10 竹広真一 msgdmp -> MessageNotify に変更
| alu(:,: ) : | real(8), intent(inout)
| ||
| kp(size(alu,1)) : | integer, intent(out)
|
subroutine ludecomp21(alu,kp)
! ALU(NDIM,NDIM), KP(NDIM)
! NDIM x NDIM の行列を LU 分解.
! LU行列は 入力行列に上書きされる.
real(8), intent(inout) :: alu(:,: ) ! 入力/LU行列
integer, intent(out) :: kp(size(alu,1)) ! ピボット
if ( size(alu,1) > size(alu,2) ) then
call MessageNotify('E','ludecomp', 'The third dimension is less than the second')
elseif( size(alu,1) < size(alu,2) ) then
call MessageNotify('W','ludecomp', 'The third dimension is grater than the second')
endif
!" 行列のLU分解(部分ピボット選択)
call LUMAKE( alu, kp, 1, size(alu,1) )
end subroutine ludecomp21
| alu(:,:,:) : | real(8), intent(inout)
| ||
| kp(size(alu,1),size(alu,2)) : | integer, intent(out)
|
subroutine ludecomp32(alu,kp)
! ALU(JDIM,NDIM,NDIM), KP(JDIM,NDIM)
! NDIM x NDIM の行列 JDIM 個を一度に LU 分解.
! LU行列は 入力行列に上書きされる.
real(8), intent(inout) :: alu(:,:,:) ! 入力/LU行列
integer, intent(out) :: kp(size(alu,1),size(alu,2)) ! ピボット
if ( size(alu,2) > size(alu,3) ) then
call MessageNotify('E','ludecomp', 'The third dimension is less than the second')
elseif( size(alu,2) < size(alu,3) ) then
call MessageNotify('W','ludecomp', 'The third dimension is grater than the second')
endif
!" 行列のLU分解(部分ピボット選択)
call LUMAKE( alu, kp, size(alu,1), size(alu,2) )
end subroutine ludecomp32
| lusolve211(size(b)) : | real(8)
| ||
| alu(:,:) : | real(8), intent(in)
| ||
| kp(:) : | integer, intent(in)
| ||
| b(:) : | real(8), intent(in)
|
function lusolve211(alu,kp,b)
! ALU(NDIM,NDIM), KP(NDIM), B(NDIM)
! NDIM x NDIM 型行列の連立方程式
! A X = B を 1 個の B に対して計算する.
real(8), intent(in) :: alu(:,:) ! 入力/LU行列
integer, intent(in) :: kp(:) ! ピボット
real(8), intent(in) :: b(:) ! 右辺ベクトル
real(8) :: lusolve211(size(b)) ! 解
lusolve211 = b
call LUSOL2( lusolve211, alu , kp, 1, size(b) )
end function lusolve211
| lusolve212(size(b,1),size(b,2)) : | real(8)
| ||
| alu(:,:) : | real(8), intent(in)
| ||
| kp(:) : | integer, intent(in)
| ||
| b(:,:) : | real(8), intent(in)
|
function lusolve212(alu,kp,b)
! ALU(NDIM,NDIM), KP(NDIM), B(JDIM,NDIM)
! NDIM x NDIM 型行列の連立方程式
! A X = B を JDIM 個の B に対して計算する.
real(8), intent(in) :: alu(:,:) ! 入力/LU行列
integer, intent(in) :: kp(:) ! ピボット
real(8), intent(in) :: b(:,:) ! 右辺ベクトル
real(8) :: lusolve212(size(b,1),size(b,2)) ! 解
lusolve212 = b
call LUSOL2( lusolve212, alu , kp, size(b,1), size(b,2) )
end function lusolve212
| lusolve322(size(b,1),size(b,2)) : | real(8)
| ||
| alu(:,:,:) : | real(8), intent(in)
| ||
| kp(:,:) : | integer, intent(in)
| ||
| b(:,:) : | real(8), intent(in)
|
function lusolve322(alu,kp,b)
! ALU(JDIM,NDIM,NDIM), KP(JDIM,NDIM), B(JDIM,NDIM)
! NDIM x NDIM 型行列を JDIM 個並べた連立方程式
! A X = B をひとつの B の並びに対して計算する.
real(8), intent(in) :: alu(:,:,:) ! 入力/LU行列
integer, intent(in) :: kp(:,:) ! ピボット
real(8), intent(in) :: b(:,:) ! 右辺ベクトル
real(8) :: lusolve322(size(b,1),size(b,2)) ! 解
lusolve322 = b
call LUSOLV( lusolve322, alu , kp, 1, size(b,1), size(b,2) )
end function lusolve322
| lusolve323(size(b,1),size(b,2),size(b,3)) : | real(8)
| ||
| alu(:,:,:) : | real(8), intent(in)
| ||
| kp(:,:) : | integer, intent(in)
| ||
| b(:,:,:) : | real(8), intent(in)
|
function lusolve323(alu,kp,b)
! ALU(JDIM,NDIM,NDIM), KP(JDIM,NDIM), B(IDIM,JDIM,NDIM)
! NDIM x NDIM 型行列を JDIM 個並べた連立方程式
! A X = B を IDIM 個の B に対して計算する.
real(8), intent(in) :: alu(:,:,:) ! 入力/LU行列
integer, intent(in) :: kp(:,:) ! ピボット
real(8), intent(in) :: b(:,:,:) ! 右辺ベクトル
real(8) :: lusolve323(size(b,1),size(b,2),size(b,3)) ! 解
lusolve323 = b
call LUSOLV( lusolve323, alu , kp, size(b,1), size(b,2), size(b,3) )
end function lusolve323