LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cgbcon.f
Go to the documentation of this file.
1*> \brief \b CGBCON
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CGBCON + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgbcon.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgbcon.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgbcon.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CGBCON( NORM, N, KL, KU, AB, LDAB, IPIV, ANORM, RCOND,
20* WORK, RWORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER NORM
24* INTEGER INFO, KL, KU, LDAB, N
25* REAL ANORM, RCOND
26* ..
27* .. Array Arguments ..
28* INTEGER IPIV( * )
29* REAL RWORK( * )
30* COMPLEX AB( LDAB, * ), WORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> CGBCON estimates the reciprocal of the condition number of a complex
40*> general band matrix A, in either the 1-norm or the infinity-norm,
41*> using the LU factorization computed by CGBTRF.
42*>
43*> An estimate is obtained for norm(inv(A)), and the reciprocal of the
44*> condition number is computed as
45*> RCOND = 1 / ( norm(A) * norm(inv(A)) ).
46*> \endverbatim
47*
48* Arguments:
49* ==========
50*
51*> \param[in] NORM
52*> \verbatim
53*> NORM is CHARACTER*1
54*> Specifies whether the 1-norm condition number or the
55*> infinity-norm condition number is required:
56*> = '1' or 'O': 1-norm;
57*> = 'I': Infinity-norm.
58*> \endverbatim
59*>
60*> \param[in] N
61*> \verbatim
62*> N is INTEGER
63*> The order of the matrix A. N >= 0.
64*> \endverbatim
65*>
66*> \param[in] KL
67*> \verbatim
68*> KL is INTEGER
69*> The number of subdiagonals within the band of A. KL >= 0.
70*> \endverbatim
71*>
72*> \param[in] KU
73*> \verbatim
74*> KU is INTEGER
75*> The number of superdiagonals within the band of A. KU >= 0.
76*> \endverbatim
77*>
78*> \param[in] AB
79*> \verbatim
80*> AB is COMPLEX array, dimension (LDAB,N)
81*> Details of the LU factorization of the band matrix A, as
82*> computed by CGBTRF. U is stored as an upper triangular band
83*> matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and
84*> the multipliers used during the factorization are stored in
85*> rows KL+KU+2 to 2*KL+KU+1.
86*> \endverbatim
87*>
88*> \param[in] LDAB
89*> \verbatim
90*> LDAB is INTEGER
91*> The leading dimension of the array AB. LDAB >= 2*KL+KU+1.
92*> \endverbatim
93*>
94*> \param[in] IPIV
95*> \verbatim
96*> IPIV is INTEGER array, dimension (N)
97*> The pivot indices; for 1 <= i <= N, row i of the matrix was
98*> interchanged with row IPIV(i).
99*> \endverbatim
100*>
101*> \param[in] ANORM
102*> \verbatim
103*> ANORM is REAL
104*> If NORM = '1' or 'O', the 1-norm of the original matrix A.
105*> If NORM = 'I', the infinity-norm of the original matrix A.
106*> \endverbatim
107*>
108*> \param[out] RCOND
109*> \verbatim
110*> RCOND is REAL
111*> The reciprocal of the condition number of the matrix A,
112*> computed as RCOND = 1/(norm(A) * norm(inv(A))).
113*> \endverbatim
114*>
115*> \param[out] WORK
116*> \verbatim
117*> WORK is COMPLEX array, dimension (2*N)
118*> \endverbatim
119*>
120*> \param[out] RWORK
121*> \verbatim
122*> RWORK is REAL array, dimension (N)
123*> \endverbatim
124*>
125*> \param[out] INFO
126*> \verbatim
127*> INFO is INTEGER
128*> = 0: successful exit
129*> < 0: if INFO = -i, the i-th argument had an illegal value
130*> \endverbatim
131*
132* Authors:
133* ========
134*
135*> \author Univ. of Tennessee
136*> \author Univ. of California Berkeley
137*> \author Univ. of Colorado Denver
138*> \author NAG Ltd.
139*
140*> \ingroup gbcon
141*
142* =====================================================================
143 SUBROUTINE cgbcon( NORM, N, KL, KU, AB, LDAB, IPIV, ANORM,
144 $ RCOND,
145 $ WORK, RWORK, INFO )
146*
147* -- LAPACK computational routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 CHARACTER NORM
153 INTEGER INFO, KL, KU, LDAB, N
154 REAL ANORM, RCOND
155* ..
156* .. Array Arguments ..
157 INTEGER IPIV( * )
158 REAL RWORK( * )
159 COMPLEX AB( LDAB, * ), WORK( * )
160* ..
161*
162* =====================================================================
163*
164* .. Parameters ..
165 REAL ONE, ZERO
166 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
167* ..
168* .. Local Scalars ..
169 LOGICAL LNOTI, ONENRM
170 CHARACTER NORMIN
171 INTEGER IX, J, JP, KASE, KASE1, KD, LM
172 REAL AINVNM, SCALE, SMLNUM
173 COMPLEX T, ZDUM
174* ..
175* .. Local Arrays ..
176 INTEGER ISAVE( 3 )
177* ..
178* .. External Functions ..
179 LOGICAL LSAME
180 INTEGER ICAMAX
181 REAL SLAMCH
182 COMPLEX CDOTC
183 EXTERNAL lsame, icamax, slamch, cdotc
184* ..
185* .. External Subroutines ..
186 EXTERNAL caxpy, clacn2, clatbs, csrscl,
187 $ xerbla
188* ..
189* .. Intrinsic Functions ..
190 INTRINSIC abs, aimag, min, real
191* ..
192* .. Statement Functions ..
193 REAL CABS1
194* ..
195* .. Statement Function definitions ..
196 cabs1( zdum ) = abs( real( zdum ) ) + abs( aimag( zdum ) )
197* ..
198* .. Executable Statements ..
199*
200* Test the input parameters.
201*
202 info = 0
203 onenrm = norm.EQ.'1' .OR. lsame( norm, 'O' )
204 IF( .NOT.onenrm .AND. .NOT.lsame( norm, 'I' ) ) THEN
205 info = -1
206 ELSE IF( n.LT.0 ) THEN
207 info = -2
208 ELSE IF( kl.LT.0 ) THEN
209 info = -3
210 ELSE IF( ku.LT.0 ) THEN
211 info = -4
212 ELSE IF( ldab.LT.2*kl+ku+1 ) THEN
213 info = -6
214 ELSE IF( anorm.LT.zero ) THEN
215 info = -8
216 END IF
217 IF( info.NE.0 ) THEN
218 CALL xerbla( 'CGBCON', -info )
219 RETURN
220 END IF
221*
222* Quick return if possible
223*
224 rcond = zero
225 IF( n.EQ.0 ) THEN
226 rcond = one
227 RETURN
228 ELSE IF( anorm.EQ.zero ) THEN
229 RETURN
230 END IF
231*
232 smlnum = slamch( 'Safe minimum' )
233*
234* Estimate the norm of inv(A).
235*
236 ainvnm = zero
237 normin = 'N'
238 IF( onenrm ) THEN
239 kase1 = 1
240 ELSE
241 kase1 = 2
242 END IF
243 kd = kl + ku + 1
244 lnoti = kl.GT.0
245 kase = 0
246 10 CONTINUE
247 CALL clacn2( n, work( n+1 ), work, ainvnm, kase, isave )
248 IF( kase.NE.0 ) THEN
249 IF( kase.EQ.kase1 ) THEN
250*
251* Multiply by inv(L).
252*
253 IF( lnoti ) THEN
254 DO 20 j = 1, n - 1
255 lm = min( kl, n-j )
256 jp = ipiv( j )
257 t = work( jp )
258 IF( jp.NE.j ) THEN
259 work( jp ) = work( j )
260 work( j ) = t
261 END IF
262 CALL caxpy( lm, -t, ab( kd+1, j ), 1, work( j+1 ),
263 $ 1 )
264 20 CONTINUE
265 END IF
266*
267* Multiply by inv(U).
268*
269 CALL clatbs( 'Upper', 'No transpose', 'Non-unit', normin,
270 $ n,
271 $ kl+ku, ab, ldab, work, scale, rwork, info )
272 ELSE
273*
274* Multiply by inv(U**H).
275*
276 CALL clatbs( 'Upper', 'Conjugate transpose', 'Non-unit',
277 $ normin, n, kl+ku, ab, ldab, work, scale, rwork,
278 $ info )
279*
280* Multiply by inv(L**H).
281*
282 IF( lnoti ) THEN
283 DO 30 j = n - 1, 1, -1
284 lm = min( kl, n-j )
285 work( j ) = work( j ) - cdotc( lm, ab( kd+1, j ),
286 $ 1,
287 $ work( j+1 ), 1 )
288 jp = ipiv( j )
289 IF( jp.NE.j ) THEN
290 t = work( jp )
291 work( jp ) = work( j )
292 work( j ) = t
293 END IF
294 30 CONTINUE
295 END IF
296 END IF
297*
298* Divide X by 1/SCALE if doing so will not cause overflow.
299*
300 normin = 'Y'
301 IF( scale.NE.one ) THEN
302 ix = icamax( n, work, 1 )
303 IF( scale.LT.cabs1( work( ix ) )*smlnum .OR. scale.EQ.zero )
304 $ GO TO 40
305 CALL csrscl( n, scale, work, 1 )
306 END IF
307 GO TO 10
308 END IF
309*
310* Compute the estimate of the reciprocal condition number.
311*
312 IF( ainvnm.NE.zero )
313 $ rcond = ( one / ainvnm ) / anorm
314*
315 40 CONTINUE
316 RETURN
317*
318* End of CGBCON
319*
320 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine caxpy(n, ca, cx, incx, cy, incy)
CAXPY
Definition caxpy.f:88
subroutine cgbcon(norm, n, kl, ku, ab, ldab, ipiv, anorm, rcond, work, rwork, info)
CGBCON
Definition cgbcon.f:146
subroutine clacn2(n, v, x, est, kase, isave)
CLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition clacn2.f:131
subroutine clatbs(uplo, trans, diag, normin, n, kd, ab, ldab, x, scale, cnorm, info)
CLATBS solves a triangular banded system of equations.
Definition clatbs.f:242
subroutine csrscl(n, sa, sx, incx)
CSRSCL multiplies a vector by the reciprocal of a real scalar.
Definition csrscl.f:82