LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgesvx.f
Go to the documentation of this file.
1*> \brief <b> ZGESVX computes the solution to system of linear equations A * X = B for GE matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGESVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgesvx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgesvx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgesvx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGESVX( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV,
20* EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR,
21* WORK, RWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER EQUED, FACT, TRANS
25* INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
26* DOUBLE PRECISION RCOND
27* ..
28* .. Array Arguments ..
29* INTEGER IPIV( * )
30* DOUBLE PRECISION BERR( * ), C( * ), FERR( * ), R( * ),
31* $ RWORK( * )
32* COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
33* $ WORK( * ), X( LDX, * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> ZGESVX uses the LU factorization to compute the solution to a complex
43*> system of linear equations
44*> A * X = B,
45*> where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
46*>
47*> Error bounds on the solution and a condition estimate are also
48*> provided.
49*> \endverbatim
50*
51*> \par Description:
52* =================
53*>
54*> \verbatim
55*>
56*> The following steps are performed:
57*>
58*> 1. If FACT = 'E', real scaling factors are computed to equilibrate
59*> the system:
60*> TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
61*> TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
62*> TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
63*> Whether or not the system will be equilibrated depends on the
64*> scaling of the matrix A, but if equilibration is used, A is
65*> overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
66*> or diag(C)*B (if TRANS = 'T' or 'C').
67*>
68*> 2. If FACT = 'N' or 'E', the LU decomposition is used to factor the
69*> matrix A (after equilibration if FACT = 'E') as
70*> A = P * L * U,
71*> where P is a permutation matrix, L is a unit lower triangular
72*> matrix, and U is upper triangular.
73*>
74*> 3. If some U(i,i)=0, so that U is exactly singular, then the routine
75*> returns with INFO = i. Otherwise, the factored form of A is used
76*> to estimate the condition number of the matrix A. If the
77*> reciprocal of the condition number is less than machine precision,
78*> INFO = N+1 is returned as a warning, but the routine still goes on
79*> to solve for X and compute error bounds as described below.
80*>
81*> 4. The system of equations is solved for X using the factored form
82*> of A.
83*>
84*> 5. Iterative refinement is applied to improve the computed solution
85*> matrix and calculate error bounds and backward error estimates
86*> for it.
87*>
88*> 6. If equilibration was used, the matrix X is premultiplied by
89*> diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
90*> that it solves the original system before equilibration.
91*> \endverbatim
92*
93* Arguments:
94* ==========
95*
96*> \param[in] FACT
97*> \verbatim
98*> FACT is CHARACTER*1
99*> Specifies whether or not the factored form of the matrix A is
100*> supplied on entry, and if not, whether the matrix A should be
101*> equilibrated before it is factored.
102*> = 'F': On entry, AF and IPIV contain the factored form of A.
103*> If EQUED is not 'N', the matrix A has been
104*> equilibrated with scaling factors given by R and C.
105*> A, AF, and IPIV are not modified.
106*> = 'N': The matrix A will be copied to AF and factored.
107*> = 'E': The matrix A will be equilibrated if necessary, then
108*> copied to AF and factored.
109*> \endverbatim
110*>
111*> \param[in] TRANS
112*> \verbatim
113*> TRANS is CHARACTER*1
114*> Specifies the form of the system of equations:
115*> = 'N': A * X = B (No transpose)
116*> = 'T': A**T * X = B (Transpose)
117*> = 'C': A**H * X = B (Conjugate transpose)
118*> \endverbatim
119*>
120*> \param[in] N
121*> \verbatim
122*> N is INTEGER
123*> The number of linear equations, i.e., the order of the
124*> matrix A. N >= 0.
125*> \endverbatim
126*>
127*> \param[in] NRHS
128*> \verbatim
129*> NRHS is INTEGER
130*> The number of right hand sides, i.e., the number of columns
131*> of the matrices B and X. NRHS >= 0.
132*> \endverbatim
133*>
134*> \param[in,out] A
135*> \verbatim
136*> A is COMPLEX*16 array, dimension (LDA,N)
137*> On entry, the N-by-N matrix A. If FACT = 'F' and EQUED is
138*> not 'N', then A must have been equilibrated by the scaling
139*> factors in R and/or C. A is not modified if FACT = 'F' or
140*> 'N', or if FACT = 'E' and EQUED = 'N' on exit.
141*>
142*> On exit, if EQUED .ne. 'N', A is scaled as follows:
143*> EQUED = 'R': A := diag(R) * A
144*> EQUED = 'C': A := A * diag(C)
145*> EQUED = 'B': A := diag(R) * A * diag(C).
146*> \endverbatim
147*>
148*> \param[in] LDA
149*> \verbatim
150*> LDA is INTEGER
151*> The leading dimension of the array A. LDA >= max(1,N).
152*> \endverbatim
153*>
154*> \param[in,out] AF
155*> \verbatim
156*> AF is COMPLEX*16 array, dimension (LDAF,N)
157*> If FACT = 'F', then AF is an input argument and on entry
158*> contains the factors L and U from the factorization
159*> A = P*L*U as computed by ZGETRF. If EQUED .ne. 'N', then
160*> AF is the factored form of the equilibrated matrix A.
161*>
162*> If FACT = 'N', then AF is an output argument and on exit
163*> returns the factors L and U from the factorization A = P*L*U
164*> of the original matrix A.
165*>
166*> If FACT = 'E', then AF is an output argument and on exit
167*> returns the factors L and U from the factorization A = P*L*U
168*> of the equilibrated matrix A (see the description of A for
169*> the form of the equilibrated matrix).
170*> \endverbatim
171*>
172*> \param[in] LDAF
173*> \verbatim
174*> LDAF is INTEGER
175*> The leading dimension of the array AF. LDAF >= max(1,N).
176*> \endverbatim
177*>
178*> \param[in,out] IPIV
179*> \verbatim
180*> IPIV is INTEGER array, dimension (N)
181*> If FACT = 'F', then IPIV is an input argument and on entry
182*> contains the pivot indices from the factorization A = P*L*U
183*> as computed by ZGETRF; row i of the matrix was interchanged
184*> with row IPIV(i).
185*>
186*> If FACT = 'N', then IPIV is an output argument and on exit
187*> contains the pivot indices from the factorization A = P*L*U
188*> of the original matrix A.
189*>
190*> If FACT = 'E', then IPIV is an output argument and on exit
191*> contains the pivot indices from the factorization A = P*L*U
192*> of the equilibrated matrix A.
193*> \endverbatim
194*>
195*> \param[in,out] EQUED
196*> \verbatim
197*> EQUED is CHARACTER*1
198*> Specifies the form of equilibration that was done.
199*> = 'N': No equilibration (always true if FACT = 'N').
200*> = 'R': Row equilibration, i.e., A has been premultiplied by
201*> diag(R).
202*> = 'C': Column equilibration, i.e., A has been postmultiplied
203*> by diag(C).
204*> = 'B': Both row and column equilibration, i.e., A has been
205*> replaced by diag(R) * A * diag(C).
206*> EQUED is an input argument if FACT = 'F'; otherwise, it is an
207*> output argument.
208*> \endverbatim
209*>
210*> \param[in,out] R
211*> \verbatim
212*> R is DOUBLE PRECISION array, dimension (N)
213*> The row scale factors for A. If EQUED = 'R' or 'B', A is
214*> multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
215*> is not accessed. R is an input argument if FACT = 'F';
216*> otherwise, R is an output argument. If FACT = 'F' and
217*> EQUED = 'R' or 'B', each element of R must be positive.
218*> \endverbatim
219*>
220*> \param[in,out] C
221*> \verbatim
222*> C is DOUBLE PRECISION array, dimension (N)
223*> The column scale factors for A. If EQUED = 'C' or 'B', A is
224*> multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
225*> is not accessed. C is an input argument if FACT = 'F';
226*> otherwise, C is an output argument. If FACT = 'F' and
227*> EQUED = 'C' or 'B', each element of C must be positive.
228*> \endverbatim
229*>
230*> \param[in,out] B
231*> \verbatim
232*> B is COMPLEX*16 array, dimension (LDB,NRHS)
233*> On entry, the N-by-NRHS right hand side matrix B.
234*> On exit,
235*> if EQUED = 'N', B is not modified;
236*> if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
237*> diag(R)*B;
238*> if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
239*> overwritten by diag(C)*B.
240*> \endverbatim
241*>
242*> \param[in] LDB
243*> \verbatim
244*> LDB is INTEGER
245*> The leading dimension of the array B. LDB >= max(1,N).
246*> \endverbatim
247*>
248*> \param[out] X
249*> \verbatim
250*> X is COMPLEX*16 array, dimension (LDX,NRHS)
251*> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X
252*> to the original system of equations. Note that A and B are
253*> modified on exit if EQUED .ne. 'N', and the solution to the
254*> equilibrated system is inv(diag(C))*X if TRANS = 'N' and
255*> EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C'
256*> and EQUED = 'R' or 'B'.
257*> \endverbatim
258*>
259*> \param[in] LDX
260*> \verbatim
261*> LDX is INTEGER
262*> The leading dimension of the array X. LDX >= max(1,N).
263*> \endverbatim
264*>
265*> \param[out] RCOND
266*> \verbatim
267*> RCOND is DOUBLE PRECISION
268*> The estimate of the reciprocal condition number of the matrix
269*> A after equilibration (if done). If RCOND is less than the
270*> machine precision (in particular, if RCOND = 0), the matrix
271*> is singular to working precision. This condition is
272*> indicated by a return code of INFO > 0.
273*> \endverbatim
274*>
275*> \param[out] FERR
276*> \verbatim
277*> FERR is DOUBLE PRECISION array, dimension (NRHS)
278*> The estimated forward error bound for each solution vector
279*> X(j) (the j-th column of the solution matrix X).
280*> If XTRUE is the true solution corresponding to X(j), FERR(j)
281*> is an estimated upper bound for the magnitude of the largest
282*> element in (X(j) - XTRUE) divided by the magnitude of the
283*> largest element in X(j). The estimate is as reliable as
284*> the estimate for RCOND, and is almost always a slight
285*> overestimate of the true error.
286*> \endverbatim
287*>
288*> \param[out] BERR
289*> \verbatim
290*> BERR is DOUBLE PRECISION array, dimension (NRHS)
291*> The componentwise relative backward error of each solution
292*> vector X(j) (i.e., the smallest relative change in
293*> any element of A or B that makes X(j) an exact solution).
294*> \endverbatim
295*>
296*> \param[out] WORK
297*> \verbatim
298*> WORK is COMPLEX*16 array, dimension (2*N)
299*> \endverbatim
300*>
301*> \param[out] RWORK
302*> \verbatim
303*> RWORK is DOUBLE PRECISION array, dimension (MAX(1,2*N))
304*> On exit, RWORK(1) contains the reciprocal pivot growth
305*> factor norm(A)/norm(U). The "max absolute element" norm is
306*> used. If RWORK(1) is much less than 1, then the stability
307*> of the LU factorization of the (equilibrated) matrix A
308*> could be poor. This also means that the solution X, condition
309*> estimator RCOND, and forward error bound FERR could be
310*> unreliable. If factorization fails with 0<INFO<=N, then
311*> RWORK(1) contains the reciprocal pivot growth factor for the
312*> leading INFO columns of A.
313*> \endverbatim
314*>
315*> \param[out] INFO
316*> \verbatim
317*> INFO is INTEGER
318*> = 0: successful exit
319*> < 0: if INFO = -i, the i-th argument had an illegal value
320*> > 0: if INFO = i, and i is
321*> <= N: U(i,i) is exactly zero. The factorization has
322*> been completed, but the factor U is exactly
323*> singular, so the solution and error bounds
324*> could not be computed. RCOND = 0 is returned.
325*> = N+1: U is nonsingular, but RCOND is less than machine
326*> precision, meaning that the matrix is singular
327*> to working precision. Nevertheless, the
328*> solution and error bounds are computed because
329*> there are a number of situations where the
330*> computed solution can be more accurate than the
331*> value of RCOND would suggest.
332*> \endverbatim
333*
334* Authors:
335* ========
336*
337*> \author Univ. of Tennessee
338*> \author Univ. of California Berkeley
339*> \author Univ. of Colorado Denver
340*> \author NAG Ltd.
341*
342*> \ingroup gesvx
343*
344* =====================================================================
345 SUBROUTINE zgesvx( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF,
346 $ IPIV,
347 $ EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR,
348 $ WORK, RWORK, INFO )
349*
350* -- LAPACK driver routine --
351* -- LAPACK is a software package provided by Univ. of Tennessee, --
352* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
353*
354* .. Scalar Arguments ..
355 CHARACTER EQUED, FACT, TRANS
356 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
357 DOUBLE PRECISION RCOND
358* ..
359* .. Array Arguments ..
360 INTEGER IPIV( * )
361 DOUBLE PRECISION BERR( * ), C( * ), FERR( * ), R( * ),
362 $ RWORK( * )
363 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
364 $ work( * ), x( ldx, * )
365* ..
366*
367* =====================================================================
368*
369* .. Parameters ..
370 DOUBLE PRECISION ZERO, ONE
371 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
372* ..
373* .. Local Scalars ..
374 LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
375 CHARACTER NORM
376 INTEGER I, INFEQU, J
377 DOUBLE PRECISION AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN,
378 $ rowcnd, rpvgrw, smlnum
379* ..
380* .. External Functions ..
381 LOGICAL LSAME
382 DOUBLE PRECISION DLAMCH, ZLANGE, ZLANTR
383 EXTERNAL LSAME, DLAMCH, ZLANGE, ZLANTR
384* ..
385* .. External Subroutines ..
386 EXTERNAL xerbla, zgecon, zgeequ, zgerfs, zgetrf,
387 $ zgetrs,
388 $ zlacpy, zlaqge
389* ..
390* .. Intrinsic Functions ..
391 INTRINSIC max, min
392* ..
393* .. Executable Statements ..
394*
395 info = 0
396 nofact = lsame( fact, 'N' )
397 equil = lsame( fact, 'E' )
398 notran = lsame( trans, 'N' )
399 IF( nofact .OR. equil ) THEN
400 equed = 'N'
401 rowequ = .false.
402 colequ = .false.
403 ELSE
404 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
405 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
406 smlnum = dlamch( 'Safe minimum' )
407 bignum = one / smlnum
408 END IF
409*
410* Test the input parameters.
411*
412 IF( .NOT.nofact .AND.
413 $ .NOT.equil .AND.
414 $ .NOT.lsame( fact, 'F' ) )
415 $ THEN
416 info = -1
417 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
418 $ lsame( trans, 'C' ) ) THEN
419 info = -2
420 ELSE IF( n.LT.0 ) THEN
421 info = -3
422 ELSE IF( nrhs.LT.0 ) THEN
423 info = -4
424 ELSE IF( lda.LT.max( 1, n ) ) THEN
425 info = -6
426 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
427 info = -8
428 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
429 $ ( rowequ .OR. colequ .OR. lsame( equed, 'N' ) ) ) THEN
430 info = -10
431 ELSE
432 IF( rowequ ) THEN
433 rcmin = bignum
434 rcmax = zero
435 DO 10 j = 1, n
436 rcmin = min( rcmin, r( j ) )
437 rcmax = max( rcmax, r( j ) )
438 10 CONTINUE
439 IF( rcmin.LE.zero ) THEN
440 info = -11
441 ELSE IF( n.GT.0 ) THEN
442 rowcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
443 ELSE
444 rowcnd = one
445 END IF
446 END IF
447 IF( colequ .AND. info.EQ.0 ) THEN
448 rcmin = bignum
449 rcmax = zero
450 DO 20 j = 1, n
451 rcmin = min( rcmin, c( j ) )
452 rcmax = max( rcmax, c( j ) )
453 20 CONTINUE
454 IF( rcmin.LE.zero ) THEN
455 info = -12
456 ELSE IF( n.GT.0 ) THEN
457 colcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
458 ELSE
459 colcnd = one
460 END IF
461 END IF
462 IF( info.EQ.0 ) THEN
463 IF( ldb.LT.max( 1, n ) ) THEN
464 info = -14
465 ELSE IF( ldx.LT.max( 1, n ) ) THEN
466 info = -16
467 END IF
468 END IF
469 END IF
470*
471 IF( info.NE.0 ) THEN
472 CALL xerbla( 'ZGESVX', -info )
473 RETURN
474 END IF
475*
476 IF( equil ) THEN
477*
478* Compute row and column scalings to equilibrate the matrix A.
479*
480 CALL zgeequ( n, n, a, lda, r, c, rowcnd, colcnd, amax,
481 $ infequ )
482 IF( infequ.EQ.0 ) THEN
483*
484* Equilibrate the matrix.
485*
486 CALL zlaqge( n, n, a, lda, r, c, rowcnd, colcnd, amax,
487 $ equed )
488 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
489 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
490 END IF
491 END IF
492*
493* Scale the right hand side.
494*
495 IF( notran ) THEN
496 IF( rowequ ) THEN
497 DO 40 j = 1, nrhs
498 DO 30 i = 1, n
499 b( i, j ) = r( i )*b( i, j )
500 30 CONTINUE
501 40 CONTINUE
502 END IF
503 ELSE IF( colequ ) THEN
504 DO 60 j = 1, nrhs
505 DO 50 i = 1, n
506 b( i, j ) = c( i )*b( i, j )
507 50 CONTINUE
508 60 CONTINUE
509 END IF
510*
511 IF( nofact .OR. equil ) THEN
512*
513* Compute the LU factorization of A.
514*
515 CALL zlacpy( 'Full', n, n, a, lda, af, ldaf )
516 CALL zgetrf( n, n, af, ldaf, ipiv, info )
517*
518* Return if INFO is non-zero.
519*
520 IF( info.GT.0 ) THEN
521*
522* Compute the reciprocal pivot growth factor of the
523* leading rank-deficient INFO columns of A.
524*
525 rpvgrw = zlantr( 'M', 'U', 'N', info, info, af, ldaf,
526 $ rwork )
527 IF( rpvgrw.EQ.zero ) THEN
528 rpvgrw = one
529 ELSE
530 rpvgrw = zlange( 'M', n, info, a, lda, rwork ) /
531 $ rpvgrw
532 END IF
533 rwork( 1 ) = rpvgrw
534 rcond = zero
535 RETURN
536 END IF
537 END IF
538*
539* Compute the norm of the matrix A and the
540* reciprocal pivot growth factor RPVGRW.
541*
542 IF( notran ) THEN
543 norm = '1'
544 ELSE
545 norm = 'I'
546 END IF
547 anorm = zlange( norm, n, n, a, lda, rwork )
548 rpvgrw = zlantr( 'M', 'U', 'N', n, n, af, ldaf, rwork )
549 IF( rpvgrw.EQ.zero ) THEN
550 rpvgrw = one
551 ELSE
552 rpvgrw = zlange( 'M', n, n, a, lda, rwork ) / rpvgrw
553 END IF
554*
555* Compute the reciprocal of the condition number of A.
556*
557 CALL zgecon( norm, n, af, ldaf, anorm, rcond, work, rwork,
558 $ info )
559*
560* Compute the solution matrix X.
561*
562 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
563 CALL zgetrs( trans, n, nrhs, af, ldaf, ipiv, x, ldx, info )
564*
565* Use iterative refinement to improve the computed solution and
566* compute error bounds and backward error estimates for it.
567*
568 CALL zgerfs( trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x,
569 $ ldx, ferr, berr, work, rwork, info )
570*
571* Transform the solution matrix X to a solution of the original
572* system.
573*
574 IF( notran ) THEN
575 IF( colequ ) THEN
576 DO 80 j = 1, nrhs
577 DO 70 i = 1, n
578 x( i, j ) = c( i )*x( i, j )
579 70 CONTINUE
580 80 CONTINUE
581 DO 90 j = 1, nrhs
582 ferr( j ) = ferr( j ) / colcnd
583 90 CONTINUE
584 END IF
585 ELSE IF( rowequ ) THEN
586 DO 110 j = 1, nrhs
587 DO 100 i = 1, n
588 x( i, j ) = r( i )*x( i, j )
589 100 CONTINUE
590 110 CONTINUE
591 DO 120 j = 1, nrhs
592 ferr( j ) = ferr( j ) / rowcnd
593 120 CONTINUE
594 END IF
595*
596* Set INFO = N+1 if the matrix is singular to working precision.
597*
598 IF( rcond.LT.dlamch( 'Epsilon' ) )
599 $ info = n + 1
600*
601 rwork( 1 ) = rpvgrw
602 RETURN
603*
604* End of ZGESVX
605*
606 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgecon(norm, n, a, lda, anorm, rcond, work, rwork, info)
ZGECON
Definition zgecon.f:130
subroutine zgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
ZGEEQU
Definition zgeequ.f:138
subroutine zgerfs(trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZGERFS
Definition zgerfs.f:185
subroutine zgesvx(fact, trans, n, nrhs, a, lda, af, ldaf, ipiv, equed, r, c, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZGESVX computes the solution to system of linear equations A * X = B for GE matrices
Definition zgesvx.f:349
subroutine zgetrf(m, n, a, lda, ipiv, info)
ZGETRF
Definition zgetrf.f:106
subroutine zgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
ZGETRS
Definition zgetrs.f:119
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
subroutine zlaqge(m, n, a, lda, r, c, rowcnd, colcnd, amax, equed)
ZLAQGE scales a general rectangular matrix, using row and column scaling factors computed by sgeequ.
Definition zlaqge.f:141