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