LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssysvx.f
Go to the documentation of this file.
1*> \brief <b> SSYSVX computes the solution to system of linear equations A * X = B for SY matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SSYSVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssysvx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssysvx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssysvx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SSYSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV, B,
20* LDB, X, LDX, RCOND, FERR, BERR, WORK, LWORK,
21* IWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER FACT, UPLO
25* INTEGER INFO, LDA, LDAF, LDB, LDX, LWORK, N, NRHS
26* REAL RCOND
27* ..
28* .. Array Arguments ..
29* INTEGER IPIV( * ), IWORK( * )
30* REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
31* $ BERR( * ), FERR( * ), WORK( * ), X( LDX, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> SSYSVX uses the diagonal pivoting factorization to compute the
41*> solution to a real system of linear equations A * X = B,
42*> where A is an N-by-N symmetric matrix and X and B are N-by-NRHS
43*> matrices.
44*>
45*> Error bounds on the solution and a condition estimate are also
46*> provided.
47*> \endverbatim
48*
49*> \par Description:
50* =================
51*>
52*> \verbatim
53*>
54*> The following steps are performed:
55*>
56*> 1. If FACT = 'N', the diagonal pivoting method is used to factor A.
57*> The form of the factorization is
58*> A = U * D * U**T, if UPLO = 'U', or
59*> A = L * D * L**T, if UPLO = 'L',
60*> where U (or L) is a product of permutation and unit upper (lower)
61*> triangular matrices, and D is symmetric and block diagonal with
62*> 1-by-1 and 2-by-2 diagonal blocks.
63*>
64*> 2. If some D(i,i)=0, so that D is exactly singular, then the routine
65*> returns with INFO = i. Otherwise, the factored form of A is used
66*> to estimate the condition number of the matrix A. If the
67*> reciprocal of the condition number is less than machine precision,
68*> INFO = N+1 is returned as a warning, but the routine still goes on
69*> to solve for X and compute error bounds as described below.
70*>
71*> 3. The system of equations is solved for X using the factored form
72*> of A.
73*>
74*> 4. Iterative refinement is applied to improve the computed solution
75*> matrix and calculate error bounds and backward error estimates
76*> for it.
77*> \endverbatim
78*
79* Arguments:
80* ==========
81*
82*> \param[in] FACT
83*> \verbatim
84*> FACT is CHARACTER*1
85*> Specifies whether or not the factored form of A has been
86*> supplied on entry.
87*> = 'F': On entry, AF and IPIV contain the factored form of
88*> A. AF and IPIV will not be modified.
89*> = 'N': The matrix A will be copied to AF and factored.
90*> \endverbatim
91*>
92*> \param[in] UPLO
93*> \verbatim
94*> UPLO is CHARACTER*1
95*> = 'U': Upper triangle of A is stored;
96*> = 'L': Lower triangle of A is stored.
97*> \endverbatim
98*>
99*> \param[in] N
100*> \verbatim
101*> N is INTEGER
102*> The number of linear equations, i.e., the order of the
103*> matrix A. N >= 0.
104*> \endverbatim
105*>
106*> \param[in] NRHS
107*> \verbatim
108*> NRHS is INTEGER
109*> The number of right hand sides, i.e., the number of columns
110*> of the matrices B and X. NRHS >= 0.
111*> \endverbatim
112*>
113*> \param[in] A
114*> \verbatim
115*> A is REAL array, dimension (LDA,N)
116*> The symmetric matrix A. If UPLO = 'U', the leading N-by-N
117*> upper triangular part of A contains the upper triangular part
118*> of the matrix A, and the strictly lower triangular part of A
119*> is not referenced. If UPLO = 'L', the leading N-by-N lower
120*> triangular part of A contains the lower triangular part of
121*> the matrix A, and the strictly upper triangular part of A is
122*> not referenced.
123*> \endverbatim
124*>
125*> \param[in] LDA
126*> \verbatim
127*> LDA is INTEGER
128*> The leading dimension of the array A. LDA >= max(1,N).
129*> \endverbatim
130*>
131*> \param[in,out] AF
132*> \verbatim
133*> AF is REAL array, dimension (LDAF,N)
134*> If FACT = 'F', then AF is an input argument and on entry
135*> contains the block diagonal matrix D and the multipliers used
136*> to obtain the factor U or L from the factorization
137*> A = U*D*U**T or A = L*D*L**T as computed by SSYTRF.
138*>
139*> If FACT = 'N', then AF is an output argument and on exit
140*> returns the block diagonal matrix D and the multipliers used
141*> to obtain the factor U or L from the factorization
142*> A = U*D*U**T or A = L*D*L**T.
143*> \endverbatim
144*>
145*> \param[in] LDAF
146*> \verbatim
147*> LDAF is INTEGER
148*> The leading dimension of the array AF. LDAF >= max(1,N).
149*> \endverbatim
150*>
151*> \param[in,out] IPIV
152*> \verbatim
153*> IPIV is INTEGER array, dimension (N)
154*> If FACT = 'F', then IPIV is an input argument and on entry
155*> contains details of the interchanges and the block structure
156*> of D, as determined by SSYTRF.
157*> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
158*> interchanged and D(k,k) is a 1-by-1 diagonal block.
159*> If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
160*> columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
161*> is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) =
162*> IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
163*> interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
164*>
165*> If FACT = 'N', then IPIV is an output argument and on exit
166*> contains details of the interchanges and the block structure
167*> of D, as determined by SSYTRF.
168*> \endverbatim
169*>
170*> \param[in] B
171*> \verbatim
172*> B is REAL array, dimension (LDB,NRHS)
173*> The N-by-NRHS right hand side matrix B.
174*> \endverbatim
175*>
176*> \param[in] LDB
177*> \verbatim
178*> LDB is INTEGER
179*> The leading dimension of the array B. LDB >= max(1,N).
180*> \endverbatim
181*>
182*> \param[out] X
183*> \verbatim
184*> X is REAL array, dimension (LDX,NRHS)
185*> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X.
186*> \endverbatim
187*>
188*> \param[in] LDX
189*> \verbatim
190*> LDX is INTEGER
191*> The leading dimension of the array X. LDX >= max(1,N).
192*> \endverbatim
193*>
194*> \param[out] RCOND
195*> \verbatim
196*> RCOND is REAL
197*> The estimate of the reciprocal condition number of the matrix
198*> A. If RCOND is less than the machine precision (in
199*> particular, if RCOND = 0), the matrix is singular to working
200*> precision. This condition is indicated by a return code of
201*> INFO > 0.
202*> \endverbatim
203*>
204*> \param[out] FERR
205*> \verbatim
206*> FERR is REAL array, dimension (NRHS)
207*> The estimated forward error bound for each solution vector
208*> X(j) (the j-th column of the solution matrix X).
209*> If XTRUE is the true solution corresponding to X(j), FERR(j)
210*> is an estimated upper bound for the magnitude of the largest
211*> element in (X(j) - XTRUE) divided by the magnitude of the
212*> largest element in X(j). The estimate is as reliable as
213*> the estimate for RCOND, and is almost always a slight
214*> overestimate of the true error.
215*> \endverbatim
216*>
217*> \param[out] BERR
218*> \verbatim
219*> BERR is REAL array, dimension (NRHS)
220*> The componentwise relative backward error of each solution
221*> vector X(j) (i.e., the smallest relative change in
222*> any element of A or B that makes X(j) an exact solution).
223*> \endverbatim
224*>
225*> \param[out] WORK
226*> \verbatim
227*> WORK is REAL array, dimension (MAX(1,LWORK))
228*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
229*> \endverbatim
230*>
231*> \param[in] LWORK
232*> \verbatim
233*> LWORK is INTEGER
234*> The length of WORK. LWORK >= max(1,3*N), and for best
235*> performance, when FACT = 'N', LWORK >= max(1,3*N,N*NB), where
236*> NB is the optimal blocksize for SSYTRF.
237*>
238*> If LWORK = -1, then a workspace query is assumed; the routine
239*> only calculates the optimal size of the WORK array, returns
240*> this value as the first entry of the WORK array, and no error
241*> message related to LWORK is issued by XERBLA.
242*> \endverbatim
243*>
244*> \param[out] IWORK
245*> \verbatim
246*> IWORK is INTEGER array, dimension (N)
247*> \endverbatim
248*>
249*> \param[out] INFO
250*> \verbatim
251*> INFO is INTEGER
252*> = 0: successful exit
253*> < 0: if INFO = -i, the i-th argument had an illegal value
254*> > 0: if INFO = i, and i is
255*> <= N: D(i,i) is exactly zero. The factorization
256*> has been completed but the factor D is exactly
257*> singular, so the solution and error bounds could
258*> not be computed. RCOND = 0 is returned.
259*> = N+1: D is nonsingular, but RCOND is less than machine
260*> precision, meaning that the matrix is singular
261*> to working precision. Nevertheless, the
262*> solution and error bounds are computed because
263*> there are a number of situations where the
264*> computed solution can be more accurate than the
265*> value of RCOND would suggest.
266*> \endverbatim
267*
268* Authors:
269* ========
270*
271*> \author Univ. of Tennessee
272*> \author Univ. of California Berkeley
273*> \author Univ. of Colorado Denver
274*> \author NAG Ltd.
275*
276*> \ingroup hesvx
277*
278* =====================================================================
279 SUBROUTINE ssysvx( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV,
280 $ B,
281 $ LDB, X, LDX, RCOND, FERR, BERR, WORK, LWORK,
282 $ IWORK, INFO )
283*
284* -- LAPACK driver routine --
285* -- LAPACK is a software package provided by Univ. of Tennessee, --
286* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
287*
288* .. Scalar Arguments ..
289 CHARACTER FACT, UPLO
290 INTEGER INFO, LDA, LDAF, LDB, LDX, LWORK, N, NRHS
291 REAL RCOND
292* ..
293* .. Array Arguments ..
294 INTEGER IPIV( * ), IWORK( * )
295 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
296 $ BERR( * ), FERR( * ), WORK( * ), X( LDX, * )
297* ..
298*
299* =====================================================================
300*
301* .. Parameters ..
302 REAL ZERO
303 PARAMETER ( ZERO = 0.0e+0 )
304* ..
305* .. Local Scalars ..
306 LOGICAL LQUERY, NOFACT
307 INTEGER LWKMIN, LWKOPT, NB
308 REAL ANORM
309* ..
310* .. External Functions ..
311 LOGICAL LSAME
312 INTEGER ILAENV
313 REAL SLAMCH, SLANSY, SROUNDUP_LWORK
314 EXTERNAL ilaenv, lsame, slamch,
315 $ slansy, sroundup_lwork
316* ..
317* .. External Subroutines ..
318 EXTERNAL slacpy, ssycon, ssyrfs, ssytrf, ssytrs,
319 $ xerbla
320* ..
321* .. Intrinsic Functions ..
322 INTRINSIC max
323* ..
324* .. Executable Statements ..
325*
326* Test the input parameters.
327*
328 info = 0
329 nofact = lsame( fact, 'N' )
330 lquery = ( lwork.EQ.-1 )
331 lwkmin = max( 1, 3*n )
332 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
333 info = -1
334 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
335 $ .NOT.lsame( uplo, 'L' ) )
336 $ THEN
337 info = -2
338 ELSE IF( n.LT.0 ) THEN
339 info = -3
340 ELSE IF( nrhs.LT.0 ) THEN
341 info = -4
342 ELSE IF( lda.LT.max( 1, n ) ) THEN
343 info = -6
344 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
345 info = -8
346 ELSE IF( ldb.LT.max( 1, n ) ) THEN
347 info = -11
348 ELSE IF( ldx.LT.max( 1, n ) ) THEN
349 info = -13
350 ELSE IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
351 info = -18
352 END IF
353*
354 IF( info.EQ.0 ) THEN
355 lwkopt = lwkmin
356 IF( nofact ) THEN
357 nb = ilaenv( 1, 'SSYTRF', uplo, n, -1, -1, -1 )
358 lwkopt = max( lwkopt, n*nb )
359 END IF
360 work( 1 ) = sroundup_lwork(lwkopt)
361 END IF
362*
363 IF( info.NE.0 ) THEN
364 CALL xerbla( 'SSYSVX', -info )
365 RETURN
366 ELSE IF( lquery ) THEN
367 RETURN
368 END IF
369*
370 IF( nofact ) THEN
371*
372* Compute the factorization A = U*D*U**T or A = L*D*L**T.
373*
374 CALL slacpy( uplo, n, n, a, lda, af, ldaf )
375 CALL ssytrf( uplo, n, af, ldaf, ipiv, work, lwork, info )
376*
377* Return if INFO is non-zero.
378*
379 IF( info.GT.0 )THEN
380 rcond = zero
381 RETURN
382 END IF
383 END IF
384*
385* Compute the norm of the matrix A.
386*
387 anorm = slansy( 'I', uplo, n, a, lda, work )
388*
389* Compute the reciprocal of the condition number of A.
390*
391 CALL ssycon( uplo, n, af, ldaf, ipiv, anorm, rcond, work,
392 $ iwork,
393 $ info )
394*
395* Compute the solution vectors X.
396*
397 CALL slacpy( 'Full', n, nrhs, b, ldb, x, ldx )
398 CALL ssytrs( uplo, n, nrhs, af, ldaf, ipiv, x, ldx, info )
399*
400* Use iterative refinement to improve the computed solutions and
401* compute error bounds and backward error estimates for them.
402*
403 CALL ssyrfs( uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x,
404 $ ldx, ferr, berr, work, iwork, info )
405*
406* Set INFO = N+1 if the matrix is singular to working precision.
407*
408 IF( rcond.LT.slamch( 'Epsilon' ) )
409 $ info = n + 1
410*
411 work( 1 ) = sroundup_lwork(lwkopt)
412*
413 RETURN
414*
415* End of SSYSVX
416*
417 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssycon(uplo, n, a, lda, ipiv, anorm, rcond, work, iwork, info)
SSYCON
Definition ssycon.f:128
subroutine ssyrfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
SSYRFS
Definition ssyrfs.f:190
subroutine ssysvx(fact, uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, lwork, iwork, info)
SSYSVX computes the solution to system of linear equations A * X = B for SY matrices
Definition ssysvx.f:283
subroutine ssytrf(uplo, n, a, lda, ipiv, work, lwork, info)
SSYTRF
Definition ssytrf.f:180
subroutine ssytrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
SSYTRS
Definition ssytrs.f:118
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101