dpofa Subroutine

public subroutine dpofa(a, lda, n, info)

dpofa factors a real symmetric positive definite matrix.

History

  • linpack. this version dated 08/14/78 . cleve moler, university of new mexico, argonne national lab.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout) :: a(lda,*)

Dimension (lda, n):

  • On entry: the symmetric matrix to be factored. only the diagonal and upper triangle are used.
  • On return: an upper triangular matrix r so that a = trans(r)*r where trans(r) is the transpose. the strict lower triangle is unaltered. if info /= 0, the factorization is not complete.
integer, intent(in) :: lda

the leading dimension of the array a.

integer, intent(in) :: n

the order of the matrix a.

integer, intent(out) :: info
  • info = 0 for normal return.
  • info = k signals an error condition. the leading minor of order k is not positive definite.

Calls

proc~~dpofa~~CallsGraph proc~dpofa lbfgsb_linpack_module::dpofa proc~ddot lbfgsb_blas_module::ddot proc~dpofa->proc~ddot

Called by

proc~~dpofa~~CalledByGraph proc~dpofa lbfgsb_linpack_module::dpofa proc~formk lbfgsb_module::formk proc~formk->proc~dpofa proc~formt lbfgsb_module::formt proc~formt->proc~dpofa proc~mainlb lbfgsb_module::mainlb proc~mainlb->proc~formk proc~mainlb->proc~formt proc~setulb lbfgsb_module::setulb proc~setulb->proc~mainlb

Source Code

   subroutine dpofa(a,lda,n,info)

   integer,intent(in) :: lda !! the leading dimension of the array `a`.
   integer,intent(in) :: n !! the order of the matrix a.
   integer,intent(out) :: info !!  * `info = 0` for normal return.
                               !!  * `info = k` signals an error condition.  the leading minor
                               !!     of order `k` is not positive definite.
   real(wp),intent(inout) :: a(lda,*) !! Dimension `(lda, n)`:
                                      !!
                                      !!  * On entry: the symmetric matrix to be factored.  only the
                                      !!    diagonal and upper triangle are used.
                                      !!  * On return: an upper triangular matrix `r` so that `a = trans(r)*r`
                                      !!    where `trans(r)` is the transpose.
                                      !!    the strict lower triangle is unaltered.
                                      !!    if `info /= 0`, the factorization is not complete.

   real(wp) :: t
   real(wp) :: s
   integer :: j,jm1,k

      do j = 1, n
         info = j
         s = 0.0_wp
         jm1 = j - 1
         if (jm1 >= 1) then
            do k = 1, jm1
               t = a(k,j) - ddot(k-1,a(1,k),1,a(1,j),1)
               t = t/a(k,k)
               a(k,j) = t
               s = s + t*t
            end do
         end if
         s = a(j,j) - s
         if (s <= 0.0_wp) return
         a(j,j) = sqrt(s)
      end do
      info = 0
   end subroutine dpofa