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