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