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