LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dposvxx.f
Go to the documentation of this file.
1*> \brief <b> DPOSVXX computes the solution to system of linear equations A * X = B for PO matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DPOSVXX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dposvxx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dposvxx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dposvxx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
20* S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
21* N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
22* NPARAMS, PARAMS, WORK, IWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER EQUED, FACT, UPLO
26* INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
27* $ N_ERR_BNDS
28* DOUBLE PRECISION RCOND, RPVGRW
29* ..
30* .. Array Arguments ..
31* INTEGER IWORK( * )
32* DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
33* $ X( LDX, * ), WORK( * )
34* DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ),
35* $ ERR_BNDS_NORM( NRHS, * ),
36* $ ERR_BNDS_COMP( NRHS, * )
37* ..
38*
39*
40*> \par Purpose:
41* =============
42*>
43*> \verbatim
44*>
45*> DPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
46*> to compute the solution to a double precision system of linear equations
47*> A * X = B, where A is an N-by-N symmetric positive definite matrix
48*> and X and B are N-by-NRHS matrices.
49*>
50*> If requested, both normwise and maximum componentwise error bounds
51*> are returned. DPOSVXX will return a solution with a tiny
52*> guaranteed error (O(eps) where eps is the working machine
53*> precision) unless the matrix is very ill-conditioned, in which
54*> case a warning is returned. Relevant condition numbers also are
55*> calculated and returned.
56*>
57*> DPOSVXX accepts user-provided factorizations and equilibration
58*> factors; see the definitions of the FACT and EQUED options.
59*> Solving with refinement and using a factorization from a previous
60*> DPOSVXX call will also produce a solution with either O(eps)
61*> errors or warnings, but we cannot make that claim for general
62*> user-provided factorizations and equilibration factors if they
63*> differ from what DPOSVXX would itself produce.
64*> \endverbatim
65*
66*> \par Description:
67* =================
68*>
69*> \verbatim
70*>
71*> The following steps are performed:
72*>
73*> 1. If FACT = 'E', double precision scaling factors are computed to equilibrate
74*> the system:
75*>
76*> diag(S)*A*diag(S) *inv(diag(S))*X = diag(S)*B
77*>
78*> Whether or not the system will be equilibrated depends on the
79*> scaling of the matrix A, but if equilibration is used, A is
80*> overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
81*>
82*> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
83*> factor the matrix A (after equilibration if FACT = 'E') as
84*> A = U**T* U, if UPLO = 'U', or
85*> A = L * L**T, if UPLO = 'L',
86*> where U is an upper triangular matrix and L is a lower triangular
87*> matrix.
88*>
89*> 3. If the leading principal minor of order i is not positive,
90*> then the routine returns with INFO = i. Otherwise, the factored
91*> form of A is used to estimate the condition number of the matrix
92*> A (see argument RCOND). If the reciprocal of the condition number
93*> is less than machine precision, the routine still goes on to solve
94*> for X and compute error bounds as described below.
95*>
96*> 4. The system of equations is solved for X using the factored form
97*> of A.
98*>
99*> 5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
100*> the routine will use iterative refinement to try to get a small
101*> error and error bounds. Refinement calculates the residual to at
102*> least twice the working precision.
103*>
104*> 6. If equilibration was used, the matrix X is premultiplied by
105*> diag(S) so that it solves the original system before
106*> equilibration.
107*> \endverbatim
108*
109* Arguments:
110* ==========
111*
112*> \verbatim
113*> Some optional parameters are bundled in the PARAMS array. These
114*> settings determine how refinement is performed, but often the
115*> defaults are acceptable. If the defaults are acceptable, users
116*> can pass NPARAMS = 0 which prevents the source code from accessing
117*> the PARAMS argument.
118*> \endverbatim
119*>
120*> \param[in] FACT
121*> \verbatim
122*> FACT is CHARACTER*1
123*> Specifies whether or not the factored form of the matrix A is
124*> supplied on entry, and if not, whether the matrix A should be
125*> equilibrated before it is factored.
126*> = 'F': On entry, AF contains the factored form of A.
127*> If EQUED is not 'N', the matrix A has been
128*> equilibrated with scaling factors given by S.
129*> A and AF are not modified.
130*> = 'N': The matrix A will be copied to AF and factored.
131*> = 'E': The matrix A will be equilibrated if necessary, then
132*> copied to AF and factored.
133*> \endverbatim
134*>
135*> \param[in] UPLO
136*> \verbatim
137*> UPLO is CHARACTER*1
138*> = 'U': Upper triangle of A is stored;
139*> = 'L': Lower triangle of A is stored.
140*> \endverbatim
141*>
142*> \param[in] N
143*> \verbatim
144*> N is INTEGER
145*> The number of linear equations, i.e., the order of the
146*> matrix A. N >= 0.
147*> \endverbatim
148*>
149*> \param[in] NRHS
150*> \verbatim
151*> NRHS is INTEGER
152*> The number of right hand sides, i.e., the number of columns
153*> of the matrices B and X. NRHS >= 0.
154*> \endverbatim
155*>
156*> \param[in,out] A
157*> \verbatim
158*> A is DOUBLE PRECISION array, dimension (LDA,N)
159*> On entry, the symmetric matrix A, except if FACT = 'F' and EQUED =
160*> 'Y', then A must contain the equilibrated matrix
161*> diag(S)*A*diag(S). If UPLO = 'U', the leading N-by-N upper
162*> triangular part of A contains the upper triangular part of the
163*> matrix A, and the strictly lower triangular part of A is not
164*> referenced. If UPLO = 'L', the leading N-by-N lower triangular
165*> part of A contains the lower triangular part of the matrix A, and
166*> the strictly upper triangular part of A is not referenced. A is
167*> not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
168*> 'N' on exit.
169*>
170*> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
171*> diag(S)*A*diag(S).
172*> \endverbatim
173*>
174*> \param[in] LDA
175*> \verbatim
176*> LDA is INTEGER
177*> The leading dimension of the array A. LDA >= max(1,N).
178*> \endverbatim
179*>
180*> \param[in,out] AF
181*> \verbatim
182*> AF is DOUBLE PRECISION array, dimension (LDAF,N)
183*> If FACT = 'F', then AF is an input argument and on entry
184*> contains the triangular factor U or L from the Cholesky
185*> factorization A = U**T*U or A = L*L**T, in the same storage
186*> format as A. If EQUED .ne. 'N', then AF is the factored
187*> form of the equilibrated matrix diag(S)*A*diag(S).
188*>
189*> If FACT = 'N', then AF is an output argument and on exit
190*> returns the triangular factor U or L from the Cholesky
191*> factorization A = U**T*U or A = L*L**T of the original
192*> matrix A.
193*>
194*> If FACT = 'E', then AF is an output argument and on exit
195*> returns the triangular factor U or L from the Cholesky
196*> factorization A = U**T*U or A = L*L**T of the equilibrated
197*> matrix A (see the description of A for the form of the
198*> equilibrated matrix).
199*> \endverbatim
200*>
201*> \param[in] LDAF
202*> \verbatim
203*> LDAF is INTEGER
204*> The leading dimension of the array AF. LDAF >= max(1,N).
205*> \endverbatim
206*>
207*> \param[in,out] EQUED
208*> \verbatim
209*> EQUED is CHARACTER*1
210*> Specifies the form of equilibration that was done.
211*> = 'N': No equilibration (always true if FACT = 'N').
212*> = 'Y': Both row and column equilibration, i.e., A has been
213*> replaced by diag(S) * A * diag(S).
214*> EQUED is an input argument if FACT = 'F'; otherwise, it is an
215*> output argument.
216*> \endverbatim
217*>
218*> \param[in,out] S
219*> \verbatim
220*> S is DOUBLE PRECISION array, dimension (N)
221*> The row scale factors for A. If EQUED = 'Y', A is multiplied on
222*> the left and right by diag(S). S is an input argument if FACT =
223*> 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED
224*> = 'Y', each element of S must be positive. If S is output, each
225*> element of S is a power of the radix. If S is input, each element
226*> of S should be a power of the radix to ensure a reliable solution
227*> and error estimates. Scaling by powers of the radix does not cause
228*> rounding errors unless the result underflows or overflows.
229*> Rounding errors during scaling lead to refining with a matrix that
230*> is not equivalent to the input matrix, producing error estimates
231*> that may not be reliable.
232*> \endverbatim
233*>
234*> \param[in,out] B
235*> \verbatim
236*> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
237*> On entry, the N-by-NRHS right hand side matrix B.
238*> On exit,
239*> if EQUED = 'N', B is not modified;
240*> if EQUED = 'Y', B is overwritten by diag(S)*B;
241*> \endverbatim
242*>
243*> \param[in] LDB
244*> \verbatim
245*> LDB is INTEGER
246*> The leading dimension of the array B. LDB >= max(1,N).
247*> \endverbatim
248*>
249*> \param[out] X
250*> \verbatim
251*> X is DOUBLE PRECISION array, dimension (LDX,NRHS)
252*> If INFO = 0, the N-by-NRHS solution matrix X to the original
253*> system of equations. Note that A and B are modified on exit if
254*> EQUED .ne. 'N', and the solution to the equilibrated system is
255*> inv(diag(S))*X.
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*> Reciprocal scaled condition number. This is an estimate of the
268*> reciprocal Skeel condition number of the matrix A after
269*> equilibration (if done). If this is less than the machine
270*> precision (in particular, if it is zero), the matrix is singular
271*> to working precision. Note that the error may still be small even
272*> if this number is very small and the matrix appears ill-
273*> conditioned.
274*> \endverbatim
275*>
276*> \param[out] RPVGRW
277*> \verbatim
278*> RPVGRW is DOUBLE PRECISION
279*> Reciprocal pivot growth. On exit, this contains the reciprocal
280*> pivot growth factor norm(A)/norm(U). The "max absolute element"
281*> norm is used. If this is much less than 1, then the stability of
282*> the LU factorization of the (equilibrated) matrix A could be poor.
283*> This also means that the solution X, estimated condition numbers,
284*> and error bounds could be unreliable. If factorization fails with
285*> 0<INFO<=N, then this contains the reciprocal pivot growth factor
286*> for the leading INFO columns of A.
287*> \endverbatim
288*>
289*> \param[out] BERR
290*> \verbatim
291*> BERR is DOUBLE PRECISION array, dimension (NRHS)
292*> Componentwise relative backward error. This is the
293*> componentwise relative backward error of each solution vector X(j)
294*> (i.e., the smallest relative change in any element of A or B that
295*> makes X(j) an exact solution).
296*> \endverbatim
297*>
298*> \param[in] N_ERR_BNDS
299*> \verbatim
300*> N_ERR_BNDS is INTEGER
301*> Number of error bounds to return for each right hand side
302*> and each type (normwise or componentwise). See ERR_BNDS_NORM and
303*> ERR_BNDS_COMP below.
304*> \endverbatim
305*>
306*> \param[out] ERR_BNDS_NORM
307*> \verbatim
308*> ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
309*> For each right-hand side, this array contains information about
310*> various error bounds and condition numbers corresponding to the
311*> normwise relative error, which is defined as follows:
312*>
313*> Normwise relative error in the ith solution vector:
314*> max_j (abs(XTRUE(j,i) - X(j,i)))
315*> ------------------------------
316*> max_j abs(X(j,i))
317*>
318*> The array is indexed by the type of error information as described
319*> below. There currently are up to three pieces of information
320*> returned.
321*>
322*> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
323*> right-hand side.
324*>
325*> The second index in ERR_BNDS_NORM(:,err) contains the following
326*> three fields:
327*> err = 1 "Trust/don't trust" boolean. Trust the answer if the
328*> reciprocal condition number is less than the threshold
329*> sqrt(n) * dlamch('Epsilon').
330*>
331*> err = 2 "Guaranteed" error bound: The estimated forward error,
332*> almost certainly within a factor of 10 of the true error
333*> so long as the next entry is greater than the threshold
334*> sqrt(n) * dlamch('Epsilon'). This error bound should only
335*> be trusted if the previous boolean is true.
336*>
337*> err = 3 Reciprocal condition number: Estimated normwise
338*> reciprocal condition number. Compared with the threshold
339*> sqrt(n) * dlamch('Epsilon') to determine if the error
340*> estimate is "guaranteed". These reciprocal condition
341*> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
342*> appropriately scaled matrix Z.
343*> Let Z = S*A, where S scales each row by a power of the
344*> radix so all absolute row sums of Z are approximately 1.
345*>
346*> See Lapack Working Note 165 for further details and extra
347*> cautions.
348*> \endverbatim
349*>
350*> \param[out] ERR_BNDS_COMP
351*> \verbatim
352*> ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
353*> For each right-hand side, this array contains information about
354*> various error bounds and condition numbers corresponding to the
355*> componentwise relative error, which is defined as follows:
356*>
357*> Componentwise relative error in the ith solution vector:
358*> abs(XTRUE(j,i) - X(j,i))
359*> max_j ----------------------
360*> abs(X(j,i))
361*>
362*> The array is indexed by the right-hand side i (on which the
363*> componentwise relative error depends), and the type of error
364*> information as described below. There currently are up to three
365*> pieces of information returned for each right-hand side. If
366*> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
367*> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS < 3, then at most
368*> the first (:,N_ERR_BNDS) entries are returned.
369*>
370*> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
371*> right-hand side.
372*>
373*> The second index in ERR_BNDS_COMP(:,err) contains the following
374*> three fields:
375*> err = 1 "Trust/don't trust" boolean. Trust the answer if the
376*> reciprocal condition number is less than the threshold
377*> sqrt(n) * dlamch('Epsilon').
378*>
379*> err = 2 "Guaranteed" error bound: The estimated forward error,
380*> almost certainly within a factor of 10 of the true error
381*> so long as the next entry is greater than the threshold
382*> sqrt(n) * dlamch('Epsilon'). This error bound should only
383*> be trusted if the previous boolean is true.
384*>
385*> err = 3 Reciprocal condition number: Estimated componentwise
386*> reciprocal condition number. Compared with the threshold
387*> sqrt(n) * dlamch('Epsilon') to determine if the error
388*> estimate is "guaranteed". These reciprocal condition
389*> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
390*> appropriately scaled matrix Z.
391*> Let Z = S*(A*diag(x)), where x is the solution for the
392*> current right-hand side and S scales each row of
393*> A*diag(x) by a power of the radix so all absolute row
394*> sums of Z are approximately 1.
395*>
396*> See Lapack Working Note 165 for further details and extra
397*> cautions.
398*> \endverbatim
399*>
400*> \param[in] NPARAMS
401*> \verbatim
402*> NPARAMS is INTEGER
403*> Specifies the number of parameters set in PARAMS. If <= 0, the
404*> PARAMS array is never referenced and default values are used.
405*> \endverbatim
406*>
407*> \param[in,out] PARAMS
408*> \verbatim
409*> PARAMS is DOUBLE PRECISION array, dimension NPARAMS
410*> Specifies algorithm parameters. If an entry is < 0.0, then
411*> that entry will be filled with default value used for that
412*> parameter. Only positions up to NPARAMS are accessed; defaults
413*> are used for higher-numbered parameters.
414*>
415*> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
416*> refinement or not.
417*> Default: 1.0D+0
418*> = 0.0: No refinement is performed, and no error bounds are
419*> computed.
420*> = 1.0: Use the extra-precise refinement algorithm.
421*> (other values are reserved for future use)
422*>
423*> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
424*> computations allowed for refinement.
425*> Default: 10
426*> Aggressive: Set to 100 to permit convergence using approximate
427*> factorizations or factorizations other than LU. If
428*> the factorization uses a technique other than
429*> Gaussian elimination, the guarantees in
430*> err_bnds_norm and err_bnds_comp may no longer be
431*> trustworthy.
432*>
433*> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
434*> will attempt to find a solution with small componentwise
435*> relative error in the double-precision algorithm. Positive
436*> is true, 0.0 is false.
437*> Default: 1.0 (attempt componentwise convergence)
438*> \endverbatim
439*>
440*> \param[out] WORK
441*> \verbatim
442*> WORK is DOUBLE PRECISION array, dimension (4*N)
443*> \endverbatim
444*>
445*> \param[out] IWORK
446*> \verbatim
447*> IWORK is INTEGER array, dimension (N)
448*> \endverbatim
449*>
450*> \param[out] INFO
451*> \verbatim
452*> INFO is INTEGER
453*> = 0: Successful exit. The solution to every right-hand side is
454*> guaranteed.
455*> < 0: If INFO = -i, the i-th argument had an illegal value
456*> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
457*> has been completed, but the factor U is exactly singular, so
458*> the solution and error bounds could not be computed. RCOND = 0
459*> is returned.
460*> = N+J: The solution corresponding to the Jth right-hand side is
461*> not guaranteed. The solutions corresponding to other right-
462*> hand sides K with K > J may not be guaranteed as well, but
463*> only the first such right-hand side is reported. If a small
464*> componentwise error is not requested (PARAMS(3) = 0.0) then
465*> the Jth right-hand side is the first with a normwise error
466*> bound that is not guaranteed (the smallest J such
467*> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
468*> the Jth right-hand side is the first with either a normwise or
469*> componentwise error bound that is not guaranteed (the smallest
470*> J such that either ERR_BNDS_NORM(J,1) = 0.0 or
471*> ERR_BNDS_COMP(J,1) = 0.0). See the definition of
472*> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
473*> about all of the right-hand sides check ERR_BNDS_NORM or
474*> ERR_BNDS_COMP.
475*> \endverbatim
476*
477* Authors:
478* ========
479*
480*> \author Univ. of Tennessee
481*> \author Univ. of California Berkeley
482*> \author Univ. of Colorado Denver
483*> \author NAG Ltd.
484*
485*> \ingroup posvxx
486*
487* =====================================================================
488 SUBROUTINE dposvxx( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF,
489 $ EQUED,
490 $ S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
491 $ N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
492 $ NPARAMS, PARAMS, WORK, IWORK, INFO )
493*
494* -- LAPACK driver routine --
495* -- LAPACK is a software package provided by Univ. of Tennessee, --
496* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
497*
498* .. Scalar Arguments ..
499 CHARACTER EQUED, FACT, UPLO
500 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
501 $ N_ERR_BNDS
502 DOUBLE PRECISION RCOND, RPVGRW
503* ..
504* .. Array Arguments ..
505 INTEGER IWORK( * )
506 DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
507 $ X( LDX, * ), WORK( * )
508 DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ),
509 $ err_bnds_norm( nrhs, * ),
510 $ err_bnds_comp( nrhs, * )
511* ..
512*
513* ==================================================================
514*
515* .. Parameters ..
516 DOUBLE PRECISION ZERO, ONE
517 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
518 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
519 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
520 INTEGER CMP_ERR_I, PIV_GROWTH_I
521 parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
522 $ berr_i = 3 )
523 parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
524 parameter( cmp_rcond_i = 7, cmp_err_i = 8,
525 $ piv_growth_i = 9 )
526* ..
527* .. Local Scalars ..
528 LOGICAL EQUIL, NOFACT, RCEQU
529 INTEGER INFEQU, J
530 DOUBLE PRECISION AMAX, BIGNUM, SMIN, SMAX,
531 $ SCOND, SMLNUM
532* ..
533* .. External Functions ..
534 EXTERNAL lsame, dlamch, dla_porpvgrw
535 LOGICAL LSAME
536 DOUBLE PRECISION DLAMCH, DLA_PORPVGRW
537* ..
538* .. External Subroutines ..
539 EXTERNAL dpoequb, dpotrf, dpotrs, dlacpy,
540 $ dlaqsy,
542* ..
543* .. Intrinsic Functions ..
544 INTRINSIC max, min
545* ..
546* .. Executable Statements ..
547*
548 info = 0
549 nofact = lsame( fact, 'N' )
550 equil = lsame( fact, 'E' )
551 smlnum = dlamch( 'Safe minimum' )
552 bignum = one / smlnum
553 IF( nofact .OR. equil ) THEN
554 equed = 'N'
555 rcequ = .false.
556 ELSE
557 rcequ = lsame( equed, 'Y' )
558 ENDIF
559*
560* Default is failure. If an input parameter is wrong or
561* factorization fails, make everything look horrible. Only the
562* pivot growth is set here, the rest is initialized in DPORFSX.
563*
564 rpvgrw = zero
565*
566* Test the input parameters. PARAMS is not tested until DPORFSX.
567*
568 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.
569 $ lsame( fact, 'F' ) ) THEN
570 info = -1
571 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
572 $ .NOT.lsame( uplo, 'L' ) ) THEN
573 info = -2
574 ELSE IF( n.LT.0 ) THEN
575 info = -3
576 ELSE IF( nrhs.LT.0 ) THEN
577 info = -4
578 ELSE IF( lda.LT.max( 1, n ) ) THEN
579 info = -6
580 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
581 info = -8
582 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
583 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
584 info = -9
585 ELSE
586 IF ( rcequ ) THEN
587 smin = bignum
588 smax = zero
589 DO 10 j = 1, n
590 smin = min( smin, s( j ) )
591 smax = max( smax, s( j ) )
592 10 CONTINUE
593 IF( smin.LE.zero ) THEN
594 info = -10
595 ELSE IF( n.GT.0 ) THEN
596 scond = max( smin, smlnum ) / min( smax, bignum )
597 ELSE
598 scond = one
599 END IF
600 END IF
601 IF( info.EQ.0 ) THEN
602 IF( ldb.LT.max( 1, n ) ) THEN
603 info = -12
604 ELSE IF( ldx.LT.max( 1, n ) ) THEN
605 info = -14
606 END IF
607 END IF
608 END IF
609*
610 IF( info.NE.0 ) THEN
611 CALL xerbla( 'DPOSVXX', -info )
612 RETURN
613 END IF
614*
615 IF( equil ) THEN
616*
617* Compute row and column scalings to equilibrate the matrix A.
618*
619 CALL dpoequb( n, a, lda, s, scond, amax, infequ )
620 IF( infequ.EQ.0 ) THEN
621*
622* Equilibrate the matrix.
623*
624 CALL dlaqsy( uplo, n, a, lda, s, scond, amax, equed )
625 rcequ = lsame( equed, 'Y' )
626 END IF
627 END IF
628*
629* Scale the right-hand side.
630*
631 IF( rcequ ) CALL dlascl2( n, nrhs, s, b, ldb )
632*
633 IF( nofact .OR. equil ) THEN
634*
635* Compute the Cholesky factorization of A.
636*
637 CALL dlacpy( uplo, n, n, a, lda, af, ldaf )
638 CALL dpotrf( uplo, n, af, ldaf, info )
639*
640* Return if INFO is non-zero.
641*
642 IF( info.NE.0 ) THEN
643*
644* Pivot in column INFO is exactly 0
645* Compute the reciprocal pivot growth factor of the
646* leading rank-deficient INFO columns of A.
647*
648 rpvgrw = dla_porpvgrw( uplo, info, a, lda, af, ldaf,
649 $ work )
650 RETURN
651 ENDIF
652 END IF
653*
654* Compute the reciprocal growth factor RPVGRW.
655*
656 rpvgrw = dla_porpvgrw( uplo, n, a, lda, af, ldaf, work )
657*
658* Compute the solution matrix X.
659*
660 CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
661 CALL dpotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
662*
663* Use iterative refinement to improve the computed solution and
664* compute error bounds and backward error estimates for it.
665*
666 CALL dporfsx( uplo, equed, n, nrhs, a, lda, af, ldaf,
667 $ s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm,
668 $ err_bnds_comp, nparams, params, work, iwork, info )
669
670*
671* Scale solutions.
672*
673 IF ( rcequ ) THEN
674 CALL dlascl2 ( n, nrhs, s, x, ldx )
675 END IF
676*
677 RETURN
678*
679* End of DPOSVXX
680*
681 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
double precision function dla_porpvgrw(uplo, ncols, a, lda, af, ldaf, work)
DLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian...
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine dlaqsy(uplo, n, a, lda, s, scond, amax, equed)
DLAQSY scales a symmetric/Hermitian matrix, using scaling factors computed by spoequ.
Definition dlaqsy.f:131
subroutine dlascl2(m, n, d, x, ldx)
DLASCL2 performs diagonal scaling on a matrix.
Definition dlascl2.f:88
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpoequb(n, a, lda, s, scond, amax, info)
DPOEQUB
Definition dpoequb.f:116
subroutine dporfsx(uplo, equed, n, nrhs, a, lda, af, ldaf, s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
DPORFSX
Definition dporfsx.f:393
subroutine dposvxx(fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
DPOSVXX computes the solution to system of linear equations A * X = B for PO matrices
Definition dposvxx.f:493
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:105
subroutine dpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
DPOTRS
Definition dpotrs.f:108