Provided by: liblapack-doc_3.12.0-3build1_all bug

NAME

       ppsv - ppsv: factor and solve

SYNOPSIS

   Functions
       subroutine cppsv (uplo, n, nrhs, ap, b, ldb, info)
            CPPSV computes the solution to system of linear equations A * X = B for OTHER
           matrices
       subroutine dppsv (uplo, n, nrhs, ap, b, ldb, info)
            DPPSV computes the solution to system of linear equations A * X = B for OTHER
           matrices
       subroutine sppsv (uplo, n, nrhs, ap, b, ldb, info)
            SPPSV computes the solution to system of linear equations A * X = B for OTHER
           matrices
       subroutine zppsv (uplo, n, nrhs, ap, b, ldb, info)
            ZPPSV computes the solution to system of linear equations A * X = B for OTHER
           matrices

Detailed Description

Function Documentation

   subroutine cppsv (character uplo, integer n, integer nrhs, complex, dimension( * ) ap,
       complex, dimension( ldb, * ) b, integer ldb, integer info)
        CPPSV computes the solution to system of linear equations A * X = B for OTHER matrices

       Purpose:

            CPPSV computes the solution to a complex system of linear equations
               A * X = B,
            where A is an N-by-N Hermitian positive definite matrix stored in
            packed format and X and B are N-by-NRHS matrices.

            The Cholesky decomposition is used to factor A as
               A = U**H * U,  if UPLO = 'U', or
               A = L * L**H,  if UPLO = 'L',
            where U is an upper triangular matrix and L is a lower triangular
            matrix.  The factored form of A is then used to solve the system of
            equations A * X = B.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     = 'U':  Upper triangle of A is stored;
                     = 'L':  Lower triangle of A is stored.

           N

                     N is INTEGER
                     The number of linear equations, i.e., the order of the
                     matrix A.  N >= 0.

           NRHS

                     NRHS is INTEGER
                     The number of right hand sides, i.e., the number of columns
                     of the matrix B.  NRHS >= 0.

           AP

                     AP is COMPLEX array, dimension (N*(N+1)/2)
                     On entry, the upper or lower triangle of the Hermitian matrix
                     A, packed columnwise in a linear array.  The j-th column of A
                     is stored in the array AP as follows:
                     if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
                     if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
                     See below for further details.

                     On exit, if INFO = 0, the factor U or L from the Cholesky
                     factorization A = U**H*U or A = L*L**H, in the same storage
                     format as A.

           B

                     B is COMPLEX array, dimension (LDB,NRHS)
                     On entry, the N-by-NRHS right hand side matrix B.
                     On exit, if INFO = 0, the N-by-NRHS solution matrix X.

           LDB

                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max(1,N).

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, the leading principal minor of order i
                           of A is not positive, so the factorization could not
                           be completed, and the solution has not been computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Further Details:

             The packed storage scheme is illustrated by the following example
             when N = 4, UPLO = 'U':

             Two-dimensional storage of the Hermitian matrix A:

                a11 a12 a13 a14
                    a22 a23 a24
                        a33 a34     (aij = conjg(aji))
                            a44

             Packed storage of the upper triangle of A:

             AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

   subroutine dppsv (character uplo, integer n, integer nrhs, double precision, dimension( * )
       ap, double precision, dimension( ldb, * ) b, integer ldb, integer info)
        DPPSV computes the solution to system of linear equations A * X = B for OTHER matrices

       Purpose:

            DPPSV computes the solution to a real system of linear equations
               A * X = B,
            where A is an N-by-N symmetric positive definite matrix stored in
            packed format and X and B are N-by-NRHS matrices.

            The Cholesky decomposition is used to factor A as
               A = U**T* U,  if UPLO = 'U', or
               A = L * L**T,  if UPLO = 'L',
            where U is an upper triangular matrix and L is a lower triangular
            matrix.  The factored form of A is then used to solve the system of
            equations A * X = B.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     = 'U':  Upper triangle of A is stored;
                     = 'L':  Lower triangle of A is stored.

           N

                     N is INTEGER
                     The number of linear equations, i.e., the order of the
                     matrix A.  N >= 0.

           NRHS

                     NRHS is INTEGER
                     The number of right hand sides, i.e., the number of columns
                     of the matrix B.  NRHS >= 0.

           AP

                     AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
                     On entry, the upper or lower triangle of the symmetric matrix
                     A, packed columnwise in a linear array.  The j-th column of A
                     is stored in the array AP as follows:
                     if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
                     if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
                     See below for further details.

                     On exit, if INFO = 0, the factor U or L from the Cholesky
                     factorization A = U**T*U or A = L*L**T, in the same storage
                     format as A.

           B

                     B is DOUBLE PRECISION array, dimension (LDB,NRHS)
                     On entry, the N-by-NRHS right hand side matrix B.
                     On exit, if INFO = 0, the N-by-NRHS solution matrix X.

           LDB

                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max(1,N).

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, the leading principal minor of order i
                           of A is not positive, so the factorization could not
                           be completed, and the solution has not been computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Further Details:

             The packed storage scheme is illustrated by the following example
             when N = 4, UPLO = 'U':

             Two-dimensional storage of the symmetric matrix A:

                a11 a12 a13 a14
                    a22 a23 a24
                        a33 a34     (aij = conjg(aji))
                            a44

             Packed storage of the upper triangle of A:

             AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

   subroutine sppsv (character uplo, integer n, integer nrhs, real, dimension( * ) ap, real,
       dimension( ldb, * ) b, integer ldb, integer info)
        SPPSV computes the solution to system of linear equations A * X = B for OTHER matrices

       Purpose:

            SPPSV computes the solution to a real system of linear equations
               A * X = B,
            where A is an N-by-N symmetric positive definite matrix stored in
            packed format and X and B are N-by-NRHS matrices.

            The Cholesky decomposition is used to factor A as
               A = U**T* U,  if UPLO = 'U', or
               A = L * L**T,  if UPLO = 'L',
            where U is an upper triangular matrix and L is a lower triangular
            matrix.  The factored form of A is then used to solve the system of
            equations A * X = B.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     = 'U':  Upper triangle of A is stored;
                     = 'L':  Lower triangle of A is stored.

           N

                     N is INTEGER
                     The number of linear equations, i.e., the order of the
                     matrix A.  N >= 0.

           NRHS

                     NRHS is INTEGER
                     The number of right hand sides, i.e., the number of columns
                     of the matrix B.  NRHS >= 0.

           AP

                     AP is REAL array, dimension (N*(N+1)/2)
                     On entry, the upper or lower triangle of the symmetric matrix
                     A, packed columnwise in a linear array.  The j-th column of A
                     is stored in the array AP as follows:
                     if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
                     if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
                     See below for further details.

                     On exit, if INFO = 0, the factor U or L from the Cholesky
                     factorization A = U**T*U or A = L*L**T, in the same storage
                     format as A.

           B

                     B is REAL array, dimension (LDB,NRHS)
                     On entry, the N-by-NRHS right hand side matrix B.
                     On exit, if INFO = 0, the N-by-NRHS solution matrix X.

           LDB

                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max(1,N).

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, the leading principal minor of order i
                           of A is not positive, so the factorization could not
                           be completed, and the solution has not been computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Further Details:

             The packed storage scheme is illustrated by the following example
             when N = 4, UPLO = 'U':

             Two-dimensional storage of the symmetric matrix A:

                a11 a12 a13 a14
                    a22 a23 a24
                        a33 a34     (aij = conjg(aji))
                            a44

             Packed storage of the upper triangle of A:

             AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

   subroutine zppsv (character uplo, integer n, integer nrhs, complex*16, dimension( * ) ap,
       complex*16, dimension( ldb, * ) b, integer ldb, integer info)
        ZPPSV computes the solution to system of linear equations A * X = B for OTHER matrices

       Purpose:

            ZPPSV computes the solution to a complex system of linear equations
               A * X = B,
            where A is an N-by-N Hermitian positive definite matrix stored in
            packed format and X and B are N-by-NRHS matrices.

            The Cholesky decomposition is used to factor A as
               A = U**H * U,  if UPLO = 'U', or
               A = L * L**H,  if UPLO = 'L',
            where U is an upper triangular matrix and L is a lower triangular
            matrix.  The factored form of A is then used to solve the system of
            equations A * X = B.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     = 'U':  Upper triangle of A is stored;
                     = 'L':  Lower triangle of A is stored.

           N

                     N is INTEGER
                     The number of linear equations, i.e., the order of the
                     matrix A.  N >= 0.

           NRHS

                     NRHS is INTEGER
                     The number of right hand sides, i.e., the number of columns
                     of the matrix B.  NRHS >= 0.

           AP

                     AP is COMPLEX*16 array, dimension (N*(N+1)/2)
                     On entry, the upper or lower triangle of the Hermitian matrix
                     A, packed columnwise in a linear array.  The j-th column of A
                     is stored in the array AP as follows:
                     if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
                     if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
                     See below for further details.

                     On exit, if INFO = 0, the factor U or L from the Cholesky
                     factorization A = U**H*U or A = L*L**H, in the same storage
                     format as A.

           B

                     B is COMPLEX*16 array, dimension (LDB,NRHS)
                     On entry, the N-by-NRHS right hand side matrix B.
                     On exit, if INFO = 0, the N-by-NRHS solution matrix X.

           LDB

                     LDB is INTEGER
                     The leading dimension of the array B.  LDB >= max(1,N).

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, the leading principal minor of order i
                           of A is not positive, so the factorization could not
                           be completed, and the solution has not been computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Further Details:

             The packed storage scheme is illustrated by the following example
             when N = 4, UPLO = 'U':

             Two-dimensional storage of the Hermitian matrix A:

                a11 a12 a13 a14
                    a22 a23 a24
                        a33 a34     (aij = conjg(aji))
                            a44

             Packed storage of the upper triangle of A:

             AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

Author

       Generated automatically by Doxygen for LAPACK from the source code.