cholsl Subroutine

public subroutine cholsl(a, n, np, p, b, x, info)

Note

Solves the set of linear equations A x = b, where A is a positive- definite symmetric matrix with physical dimensions np. A and P are are input as the output from choldc. Only the lower triangle of A is accessed. B(1:n) is inout as the right-hand side vector. The solution vector is returned in X(1:n). A, n, np, and P are not modified and can be left in place for successive calls with different right-hand sides B. B is not modified unless you identify B and X in the calling sequence, which is allowed.

(c) after Numerical recipes

Arguments

Type IntentOptional Attributes Name
real(kind=R8), intent(in), dimension(np, np) :: a

system matrix

integer(kind=I4), intent(in) :: n

system size

integer(kind=I4), intent(in) :: np

matrix size

real(kind=R8), intent(in), dimension(np) :: p

diagonal elements

real(kind=R8), intent(inout), dimension(np) :: b

rhs

real(kind=R8), intent(inout), dimension(np) :: x

solution vector

integer(kind=I4), intent(inout) :: info

information ouput


Called by

proc~~cholsl~~CalledByGraph proc~cholsl cholsl proc~moindres_carres moindres_carres proc~moindres_carres->proc~cholsl proc~moindres_carres_lineaire moindres_carres_lineaire proc~moindres_carres_lineaire->proc~cholsl program~test_utils test_utils program~test_utils->proc~cholsl proc~least_squares_tcheby least_squares_tcheby proc~least_squares_tcheby->proc~moindres_carres_lineaire program~test_least test_least program~test_least->proc~moindres_carres program~test_least->proc~moindres_carres_lineaire program~test_tchebychev test_tchebychev program~test_tchebychev->proc~least_squares_tcheby

Source Code

   subroutine cholsl (a, n, np, p, b, x, info)
   !================================================================================================
   !<@note Solves the set of linear equations A x = b, where A is a positive-
   !< definite symmetric matrix with physical dimensions np. A and P are
   !< are input as the output from choldc. Only the lower triangle of A
   !< is accessed. B(1:n) is inout as the right-hand side vector. The
   !< solution vector is returned in X(1:n). A, n, np, and P are not
   !< modified and can be left in place for successive calls with different
   !< right-hand sides B. B is not modified unless you identify B and X in
   !< the calling sequence, which is allowed.
   !<
   !< (c) after Numerical recipes
   !<
   !<@endnote
   !------------------------------------------------------------------------------------------------
   implicit none
   integer(kind=I4), intent(in   )                    :: n     !! *system size*
   integer(kind=I4), intent(in   )                    :: np    !! *matrix size*
   integer(kind=I4), intent(inout)                    :: info  !! *information ouput*
   real   (kind=R8), intent(in   ), dimension(np, np) :: a     !! *system matrix*
   real   (kind=R8), intent(in   ), dimension(np)     :: p     !! *diagonal elements*
   real   (kind=R8), intent(inout), dimension(np)     :: b     !! *rhs*
   real   (kind=R8), intent(inout), dimension(np)     :: x     !! *solution vector*

      integer(kind=I4) :: i, k
      real   (kind=R8) :: ssum

      if (info==1) return
      do i = 1, n
         ssum = b(i)
         do k = i-1, 1, -1
            ssum = ssum -a(i,k)*x(k)
         enddo
         x(i) = ssum/p(i)
      enddo
      do i = n, 1, -1
         ssum = x(i)
         do k = i+1,n
            ssum = ssum -a(k,i)*x(k)
         enddo
         x(i) = ssum/p(i)
      enddo
   return
   endsubroutine cholsl