1 *> \brief <b> DPPSV computes the solution to system of linear equations A * X = B for OTHER matrices</b>
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DPPSV + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dppsv.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dppsv.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dppsv.f">
21 * SUBROUTINE DPPSV( UPLO, N, NRHS, AP, B, LDB, INFO )
23 * .. Scalar Arguments ..
25 * INTEGER INFO, LDB, N, NRHS
27 * .. Array Arguments ..
28 * DOUBLE PRECISION AP( * ), B( LDB, * )
37 *> DPPSV computes the solution to a real system of linear equations
39 *> where A is an N-by-N symmetric positive definite matrix stored in
40 *> packed format and X and B are N-by-NRHS matrices.
42 *> The Cholesky decomposition is used to factor A as
43 *> A = U**T* U, if UPLO = 'U', or
44 *> A = L * L**T, if UPLO = 'L',
45 *> where U is an upper triangular matrix and L is a lower triangular
46 *> matrix. The factored form of A is then used to solve the system of
47 *> equations A * X = B.
55 *> UPLO is CHARACTER*1
56 *> = 'U': Upper triangle of A is stored;
57 *> = 'L': Lower triangle of A is stored.
63 *> The number of linear equations, i.e., the order of the
70 *> The number of right hand sides, i.e., the number of columns
71 *> of the matrix B. NRHS >= 0.
76 *> AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
77 *> On entry, the upper or lower triangle of the symmetric matrix
78 *> A, packed columnwise in a linear array. The j-th column of A
79 *> is stored in the array AP as follows:
80 *> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
81 *> if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
82 *> See below for further details.
84 *> On exit, if INFO = 0, the factor U or L from the Cholesky
85 *> factorization A = U**T*U or A = L*L**T, in the same storage
91 *> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
92 *> On entry, the N-by-NRHS right hand side matrix B.
93 *> On exit, if INFO = 0, the N-by-NRHS solution matrix X.
99 *> The leading dimension of the array B. LDB >= max(1,N).
105 *> = 0: successful exit
106 *> < 0: if INFO = -i, the i-th argument had an illegal value
107 *> > 0: if INFO = i, the leading minor of order i of A is not
108 *> positive definite, so the factorization could not be
109 *> completed, and the solution has not been computed.
115 *> \author Univ. of Tennessee
116 *> \author Univ. of California Berkeley
117 *> \author Univ. of Colorado Denver
120 *> \date November 2011
122 *> \ingroup doubleOTHERsolve
124 *> \par Further Details:
125 * =====================
129 *> The packed storage scheme is illustrated by the following example
130 *> when N = 4, UPLO = 'U':
132 *> Two-dimensional storage of the symmetric matrix A:
136 *> a33 a34 (aij = conjg(aji))
139 *> Packed storage of the upper triangle of A:
141 *> AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]
144 * =====================================================================
145 SUBROUTINE DPPSV( UPLO, N, NRHS, AP, B, LDB, INFO )
147 * -- LAPACK driver routine (version 3.4.0) --
148 * -- LAPACK is a software package provided by Univ. of Tennessee, --
149 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
152 * .. Scalar Arguments ..
154 INTEGER INFO, LDB, N, NRHS
156 * .. Array Arguments ..
157 DOUBLE PRECISION AP( * ), B( LDB, * )
160 * =====================================================================
162 * .. External Functions ..
166 * .. External Subroutines ..
167 EXTERNAL DPPTRF, DPPTRS, XERBLA
169 * .. Intrinsic Functions ..
172 * .. Executable Statements ..
174 * Test the input parameters.
177 IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
179 ELSE IF( N.LT.0 ) THEN
181 ELSE IF( NRHS.LT.0 ) THEN
183 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
187 CALL XERBLA( 'DPPSV ', -INFO )
191 * Compute the Cholesky factorization A = U**T*U or A = L*L**T.
193 CALL DPPTRF( UPLO, N, AP, INFO )
196 * Solve the system A*X = B, overwriting B with X.
198 CALL DPPTRS( UPLO, N, NRHS, AP, B, LDB, INFO )