LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cspsvx.f
Go to the documentation of this file.
1*> \brief <b> CSPSVX computes the solution to system of linear equations A * X = B for OTHER matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CSPSVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cspsvx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cspsvx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cspsvx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CSPSVX( FACT, UPLO, N, NRHS, AP, AFP, IPIV, B, LDB, X,
20* LDX, RCOND, FERR, BERR, WORK, RWORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER FACT, UPLO
24* INTEGER INFO, LDB, LDX, N, NRHS
25* REAL RCOND
26* ..
27* .. Array Arguments ..
28* INTEGER IPIV( * )
29* REAL BERR( * ), FERR( * ), RWORK( * )
30* COMPLEX AFP( * ), AP( * ), B( LDB, * ), WORK( * ),
31* $ X( LDX, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> CSPSVX uses the diagonal pivoting factorization A = U*D*U**T or
41*> A = L*D*L**T to compute the solution to a complex system of linear
42*> equations A * X = B, where A is an N-by-N symmetric matrix stored
43*> in packed format and X and B are N-by-NRHS matrices.
44*>
45*> Error bounds on the solution and a condition estimate are also
46*> provided.
47*> \endverbatim
48*
49*> \par Description:
50* =================
51*>
52*> \verbatim
53*>
54*> The following steps are performed:
55*>
56*> 1. If FACT = 'N', the diagonal pivoting method is used to factor A as
57*> A = U * D * U**T, if UPLO = 'U', or
58*> A = L * D * L**T, if UPLO = 'L',
59*> where U (or L) is a product of permutation and unit upper (lower)
60*> triangular matrices and D is symmetric and block diagonal with
61*> 1-by-1 and 2-by-2 diagonal blocks.
62*>
63*> 2. If some D(i,i)=0, so that D is exactly singular, then the routine
64*> returns with INFO = i. Otherwise, the factored form of A is used
65*> to estimate the condition number of the matrix A. If the
66*> reciprocal of the condition number is less than machine precision,
67*> INFO = N+1 is returned as a warning, but the routine still goes on
68*> to solve for X and compute error bounds as described below.
69*>
70*> 3. The system of equations is solved for X using the factored form
71*> of A.
72*>
73*> 4. Iterative refinement is applied to improve the computed solution
74*> matrix and calculate error bounds and backward error estimates
75*> for it.
76*> \endverbatim
77*
78* Arguments:
79* ==========
80*
81*> \param[in] FACT
82*> \verbatim
83*> FACT is CHARACTER*1
84*> Specifies whether or not the factored form of A has been
85*> supplied on entry.
86*> = 'F': On entry, AFP and IPIV contain the factored form
87*> of A. AP, AFP and IPIV will not be modified.
88*> = 'N': The matrix A will be copied to AFP and factored.
89*> \endverbatim
90*>
91*> \param[in] UPLO
92*> \verbatim
93*> UPLO is CHARACTER*1
94*> = 'U': Upper triangle of A is stored;
95*> = 'L': Lower triangle of A is stored.
96*> \endverbatim
97*>
98*> \param[in] N
99*> \verbatim
100*> N is INTEGER
101*> The number of linear equations, i.e., the order of the
102*> matrix A. N >= 0.
103*> \endverbatim
104*>
105*> \param[in] NRHS
106*> \verbatim
107*> NRHS is INTEGER
108*> The number of right hand sides, i.e., the number of columns
109*> of the matrices B and X. NRHS >= 0.
110*> \endverbatim
111*>
112*> \param[in] AP
113*> \verbatim
114*> AP is COMPLEX array, dimension (N*(N+1)/2)
115*> The upper or lower triangle of the symmetric matrix A, packed
116*> columnwise in a linear array. The j-th column of A is stored
117*> in the array AP as follows:
118*> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
119*> if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
120*> See below for further details.
121*> \endverbatim
122*>
123*> \param[in,out] AFP
124*> \verbatim
125*> AFP is COMPLEX array, dimension (N*(N+1)/2)
126*> If FACT = 'F', then AFP is an input argument and on entry
127*> contains the block diagonal matrix D and the multipliers used
128*> to obtain the factor U or L from the factorization
129*> A = U*D*U**T or A = L*D*L**T as computed by CSPTRF, stored as
130*> a packed triangular matrix in the same storage format as A.
131*>
132*> If FACT = 'N', then AFP is an output argument and on exit
133*> contains the block diagonal matrix D and the multipliers used
134*> to obtain the factor U or L from the factorization
135*> A = U*D*U**T or A = L*D*L**T as computed by CSPTRF, stored as
136*> a packed triangular matrix in the same storage format as A.
137*> \endverbatim
138*>
139*> \param[in,out] IPIV
140*> \verbatim
141*> IPIV is INTEGER array, dimension (N)
142*> If FACT = 'F', then IPIV is an input argument and on entry
143*> contains details of the interchanges and the block structure
144*> of D, as determined by CSPTRF.
145*> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
146*> interchanged and D(k,k) is a 1-by-1 diagonal block.
147*> If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
148*> columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
149*> is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) =
150*> IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
151*> interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
152*>
153*> If FACT = 'N', then IPIV is an output argument and on exit
154*> contains details of the interchanges and the block structure
155*> of D, as determined by CSPTRF.
156*> \endverbatim
157*>
158*> \param[in] B
159*> \verbatim
160*> B is COMPLEX array, dimension (LDB,NRHS)
161*> The N-by-NRHS right hand side matrix B.
162*> \endverbatim
163*>
164*> \param[in] LDB
165*> \verbatim
166*> LDB is INTEGER
167*> The leading dimension of the array B. LDB >= max(1,N).
168*> \endverbatim
169*>
170*> \param[out] X
171*> \verbatim
172*> X is COMPLEX array, dimension (LDX,NRHS)
173*> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X.
174*> \endverbatim
175*>
176*> \param[in] LDX
177*> \verbatim
178*> LDX is INTEGER
179*> The leading dimension of the array X. LDX >= max(1,N).
180*> \endverbatim
181*>
182*> \param[out] RCOND
183*> \verbatim
184*> RCOND is REAL
185*> The estimate of the reciprocal condition number of the matrix
186*> A. If RCOND is less than the machine precision (in
187*> particular, if RCOND = 0), the matrix is singular to working
188*> precision. This condition is indicated by a return code of
189*> INFO > 0.
190*> \endverbatim
191*>
192*> \param[out] FERR
193*> \verbatim
194*> FERR is REAL array, dimension (NRHS)
195*> The estimated forward error bound for each solution vector
196*> X(j) (the j-th column of the solution matrix X).
197*> If XTRUE is the true solution corresponding to X(j), FERR(j)
198*> is an estimated upper bound for the magnitude of the largest
199*> element in (X(j) - XTRUE) divided by the magnitude of the
200*> largest element in X(j). The estimate is as reliable as
201*> the estimate for RCOND, and is almost always a slight
202*> overestimate of the true error.
203*> \endverbatim
204*>
205*> \param[out] BERR
206*> \verbatim
207*> BERR is REAL array, dimension (NRHS)
208*> The componentwise relative backward error of each solution
209*> vector X(j) (i.e., the smallest relative change in
210*> any element of A or B that makes X(j) an exact solution).
211*> \endverbatim
212*>
213*> \param[out] WORK
214*> \verbatim
215*> WORK is COMPLEX array, dimension (2*N)
216*> \endverbatim
217*>
218*> \param[out] RWORK
219*> \verbatim
220*> RWORK is REAL array, dimension (N)
221*> \endverbatim
222*>
223*> \param[out] INFO
224*> \verbatim
225*> INFO is INTEGER
226*> = 0: successful exit
227*> < 0: if INFO = -i, the i-th argument had an illegal value
228*> > 0: if INFO = i, and i is
229*> <= N: D(i,i) is exactly zero. The factorization
230*> has been completed but the factor D is exactly
231*> singular, so the solution and error bounds could
232*> not be computed. RCOND = 0 is returned.
233*> = N+1: D is nonsingular, but RCOND is less than machine
234*> precision, meaning that the matrix is singular
235*> to working precision. Nevertheless, the
236*> solution and error bounds are computed because
237*> there are a number of situations where the
238*> computed solution can be more accurate than the
239*> value of RCOND would suggest.
240*> \endverbatim
241*
242* Authors:
243* ========
244*
245*> \author Univ. of Tennessee
246*> \author Univ. of California Berkeley
247*> \author Univ. of Colorado Denver
248*> \author NAG Ltd.
249*
250*> \ingroup hpsvx
251*
252*> \par Further Details:
253* =====================
254*>
255*> \verbatim
256*>
257*> The packed storage scheme is illustrated by the following example
258*> when N = 4, UPLO = 'U':
259*>
260*> Two-dimensional storage of the symmetric matrix A:
261*>
262*> a11 a12 a13 a14
263*> a22 a23 a24
264*> a33 a34 (aij = aji)
265*> a44
266*>
267*> Packed storage of the upper triangle of A:
268*>
269*> AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]
270*> \endverbatim
271*>
272* =====================================================================
273 SUBROUTINE cspsvx( FACT, UPLO, N, NRHS, AP, AFP, IPIV, B, LDB,
274 $ X,
275 $ LDX, RCOND, FERR, BERR, WORK, RWORK, INFO )
276*
277* -- LAPACK driver routine --
278* -- LAPACK is a software package provided by Univ. of Tennessee, --
279* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
280*
281* .. Scalar Arguments ..
282 CHARACTER FACT, UPLO
283 INTEGER INFO, LDB, LDX, N, NRHS
284 REAL RCOND
285* ..
286* .. Array Arguments ..
287 INTEGER IPIV( * )
288 REAL BERR( * ), FERR( * ), RWORK( * )
289 COMPLEX AFP( * ), AP( * ), B( LDB, * ), WORK( * ),
290 $ x( ldx, * )
291* ..
292*
293* =====================================================================
294*
295* .. Parameters ..
296 REAL ZERO
297 PARAMETER ( ZERO = 0.0e+0 )
298* ..
299* .. Local Scalars ..
300 LOGICAL NOFACT
301 REAL ANORM
302* ..
303* .. External Functions ..
304 LOGICAL LSAME
305 REAL CLANSP, SLAMCH
306 EXTERNAL lsame, clansp, slamch
307* ..
308* .. External Subroutines ..
309 EXTERNAL ccopy, clacpy, cspcon, csprfs, csptrf,
310 $ csptrs,
311 $ xerbla
312* ..
313* .. Intrinsic Functions ..
314 INTRINSIC max
315* ..
316* .. Executable Statements ..
317*
318* Test the input parameters.
319*
320 info = 0
321 nofact = lsame( fact, 'N' )
322 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
323 info = -1
324 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
325 $ .NOT.lsame( uplo, 'L' ) )
326 $ THEN
327 info = -2
328 ELSE IF( n.LT.0 ) THEN
329 info = -3
330 ELSE IF( nrhs.LT.0 ) THEN
331 info = -4
332 ELSE IF( ldb.LT.max( 1, n ) ) THEN
333 info = -9
334 ELSE IF( ldx.LT.max( 1, n ) ) THEN
335 info = -11
336 END IF
337 IF( info.NE.0 ) THEN
338 CALL xerbla( 'CSPSVX', -info )
339 RETURN
340 END IF
341*
342 IF( nofact ) THEN
343*
344* Compute the factorization A = U*D*U**T or A = L*D*L**T.
345*
346 CALL ccopy( n*( n+1 ) / 2, ap, 1, afp, 1 )
347 CALL csptrf( uplo, n, afp, ipiv, info )
348*
349* Return if INFO is non-zero.
350*
351 IF( info.GT.0 )THEN
352 rcond = zero
353 RETURN
354 END IF
355 END IF
356*
357* Compute the norm of the matrix A.
358*
359 anorm = clansp( 'I', uplo, n, ap, rwork )
360*
361* Compute the reciprocal of the condition number of A.
362*
363 CALL cspcon( uplo, n, afp, ipiv, anorm, rcond, work, info )
364*
365* Compute the solution vectors X.
366*
367 CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
368 CALL csptrs( uplo, n, nrhs, afp, ipiv, x, ldx, info )
369*
370* Use iterative refinement to improve the computed solutions and
371* compute error bounds and backward error estimates for them.
372*
373 CALL csprfs( uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx,
374 $ ferr,
375 $ berr, work, rwork, info )
376*
377* Set INFO = N+1 if the matrix is singular to working precision.
378*
379 IF( rcond.LT.slamch( 'Epsilon' ) )
380 $ info = n + 1
381*
382 RETURN
383*
384* End of CSPSVX
385*
386 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cspcon(uplo, n, ap, ipiv, anorm, rcond, work, info)
CSPCON
Definition cspcon.f:117
subroutine csprfs(uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CSPRFS
Definition csprfs.f:179
subroutine cspsvx(fact, uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
CSPSVX computes the solution to system of linear equations A * X = B for OTHER matrices
Definition cspsvx.f:276
subroutine csptrf(uplo, n, ap, ipiv, info)
CSPTRF
Definition csptrf.f:156
subroutine csptrs(uplo, n, nrhs, ap, ipiv, b, ldb, info)
CSPTRS
Definition csptrs.f:113
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101