Companion routine to dbnfac. it returns the solution x of the linear system a*x = b in place of b, given the lu-factorization for a in the work array w from dbnfac.
(with , as stored in w), the unit lower triangular system is solved for , and y stored in b. then the upper triangular system is solved for x. the calculations are so arranged that the innermost loops stay within columns.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | dimension(nroww,nrow) | :: | w |
describes the lu-factorization of a banded matrix a of order nrow as constructed in dbnfac. |
|
integer, | intent(in) | :: | nroww |
describes the lu-factorization of a banded matrix a of order nrow as constructed in dbnfac. |
||
integer, | intent(in) | :: | nrow |
describes the lu-factorization of a banded matrix a of order nrow as constructed in dbnfac. |
||
integer, | intent(in) | :: | nbandl |
describes the lu-factorization of a banded matrix a of order nrow as constructed in dbnfac. |
||
integer, | intent(in) | :: | nbandu |
describes the lu-factorization of a banded matrix a of order nrow as constructed in dbnfac. |
||
real(kind=wp), | intent(inout), | dimension(nrow) | :: | b |
in: right side of the system to be solved out: the solution x, of order nrow |
subroutine dbnslv(w,nroww,nrow,nbandl,nbandu,b) integer,intent(in) :: nroww !! describes the lu-factorization of a banded matrix a of order nrow as constructed in [[dbnfac]]. integer,intent(in) :: nrow !! describes the lu-factorization of a banded matrix a of order nrow as constructed in [[dbnfac]]. integer,intent(in) :: nbandl !! describes the lu-factorization of a banded matrix a of order nrow as constructed in [[dbnfac]]. integer,intent(in) :: nbandu !! describes the lu-factorization of a banded matrix a of order nrow as constructed in [[dbnfac]]. real(wp),dimension(nroww,nrow),intent(in) :: w !! describes the lu-factorization of a banded matrix a of order nrow as constructed in [[dbnfac]]. real(wp),dimension(nrow),intent(inout) :: b !! **in**: right side of the system to be solved !! **out**: the solution x, of order nrow integer :: i, j, jmax, middle, nrowm1 middle = nbandu + 1 if (nrow/=1) then nrowm1 = nrow - 1 if (nbandl/=0) then ! forward pass ! for i=1,2,...,nrow-1, subtract right side(i)*(i-th column of l) ! from right side (below i-th row). do i=1,nrowm1 jmax = min(nbandl,nrow-i) do j=1,jmax b(i+j) = b(i+j) - b(i)*w(middle+j,i) enddo enddo endif ! backward pass ! for i=nrow,nrow-1,...,1, divide right side(i) by i-th diagonal ! entry of u, then subtract right side(i)*(i-th column ! of u) from right side (above i-th row). if (nbandu<=0) then ! a is lower triangular. do i=1,nrow b(i) = b(i)/w(1,i) enddo return endif i = nrow do b(i) = b(i)/w(middle,i) jmax = min(nbandu,i-1) do j=1,jmax b(i-j) = b(i-j) - b(i)*w(middle-j,i) enddo i = i - 1 if (i<=1) exit enddo endif b(1) = b(1)/w(middle,1) endsubroutine dbnslv