LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cla_gbrcond_x.f
Go to the documentation of this file.
1*> \brief \b CLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general banded matrices.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CLA_GBRCOND_X + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cla_gbrcond_x.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cla_gbrcond_x.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cla_gbrcond_x.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* REAL FUNCTION CLA_GBRCOND_X( TRANS, N, KL, KU, AB, LDAB, AFB,
20* LDAFB, IPIV, X, INFO, WORK, RWORK )
21*
22* .. Scalar Arguments ..
23* CHARACTER TRANS
24* INTEGER N, KL, KU, KD, KE, LDAB, LDAFB, INFO
25* ..
26* .. Array Arguments ..
27* INTEGER IPIV( * )
28* COMPLEX AB( LDAB, * ), AFB( LDAFB, * ), WORK( * ),
29* $ X( * )
30* REAL RWORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> CLA_GBRCOND_X Computes the infinity norm condition number of
40*> op(A) * diag(X) where X is a COMPLEX vector.
41*> \endverbatim
42*
43* Arguments:
44* ==========
45*
46*> \param[in] TRANS
47*> \verbatim
48*> TRANS is CHARACTER*1
49*> Specifies the form of the system of equations:
50*> = 'N': A * X = B (No transpose)
51*> = 'T': A**T * X = B (Transpose)
52*> = 'C': A**H * X = B (Conjugate Transpose = Transpose)
53*> \endverbatim
54*>
55*> \param[in] N
56*> \verbatim
57*> N is INTEGER
58*> The number of linear equations, i.e., the order of the
59*> matrix A. N >= 0.
60*> \endverbatim
61*>
62*> \param[in] KL
63*> \verbatim
64*> KL is INTEGER
65*> The number of subdiagonals within the band of A. KL >= 0.
66*> \endverbatim
67*>
68*> \param[in] KU
69*> \verbatim
70*> KU is INTEGER
71*> The number of superdiagonals within the band of A. KU >= 0.
72*> \endverbatim
73*>
74*> \param[in] AB
75*> \verbatim
76*> AB is COMPLEX array, dimension (LDAB,N)
77*> On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
78*> The j-th column of A is stored in the j-th column of the
79*> array AB as follows:
80*> AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
81*> \endverbatim
82*>
83*> \param[in] LDAB
84*> \verbatim
85*> LDAB is INTEGER
86*> The leading dimension of the array AB. LDAB >= KL+KU+1.
87*> \endverbatim
88*>
89*> \param[in] AFB
90*> \verbatim
91*> AFB is COMPLEX array, dimension (LDAFB,N)
92*> Details of the LU factorization of the band matrix A, as
93*> computed by CGBTRF. U is stored as an upper triangular
94*> band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
95*> and the multipliers used during the factorization are stored
96*> in rows KL+KU+2 to 2*KL+KU+1.
97*> \endverbatim
98*>
99*> \param[in] LDAFB
100*> \verbatim
101*> LDAFB is INTEGER
102*> The leading dimension of the array AFB. LDAFB >= 2*KL+KU+1.
103*> \endverbatim
104*>
105*> \param[in] IPIV
106*> \verbatim
107*> IPIV is INTEGER array, dimension (N)
108*> The pivot indices from the factorization A = P*L*U
109*> as computed by CGBTRF; row i of the matrix was interchanged
110*> with row IPIV(i).
111*> \endverbatim
112*>
113*> \param[in] X
114*> \verbatim
115*> X is COMPLEX array, dimension (N)
116*> The vector X in the formula op(A) * diag(X).
117*> \endverbatim
118*>
119*> \param[out] INFO
120*> \verbatim
121*> INFO is INTEGER
122*> = 0: Successful exit.
123*> i > 0: The ith argument is invalid.
124*> \endverbatim
125*>
126*> \param[out] WORK
127*> \verbatim
128*> WORK is COMPLEX array, dimension (2*N).
129*> Workspace.
130*> \endverbatim
131*>
132*> \param[out] RWORK
133*> \verbatim
134*> RWORK is REAL array, dimension (N).
135*> Workspace.
136*> \endverbatim
137*
138* Authors:
139* ========
140*
141*> \author Univ. of Tennessee
142*> \author Univ. of California Berkeley
143*> \author Univ. of Colorado Denver
144*> \author NAG Ltd.
145*
146*> \ingroup la_gbrcond
147*
148* =====================================================================
149 REAL function cla_gbrcond_x( trans, n, kl, ku, ab, ldab, afb,
150 $ ldafb, ipiv, x, info, work, rwork )
151*
152* -- LAPACK computational routine --
153* -- LAPACK is a software package provided by Univ. of Tennessee, --
154* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155*
156* .. Scalar Arguments ..
157 CHARACTER trans
158 INTEGER n, kl, ku, kd, ke, ldab, ldafb, info
159* ..
160* .. Array Arguments ..
161 INTEGER ipiv( * )
162 COMPLEX ab( ldab, * ), afb( ldafb, * ), work( * ),
163 $ x( * )
164 REAL rwork( * )
165* ..
166*
167* =====================================================================
168*
169* .. Local Scalars ..
170 LOGICAL notrans
171 INTEGER kase, i, j
172 REAL ainvnm, anorm, tmp
173 COMPLEX zdum
174* ..
175* .. Local Arrays ..
176 INTEGER isave( 3 )
177* ..
178* .. External Functions ..
179 LOGICAL lsame
180 EXTERNAL lsame
181* ..
182* .. External Subroutines ..
183 EXTERNAL clacn2, cgbtrs, xerbla
184* ..
185* .. Intrinsic Functions ..
186 INTRINSIC abs, max
187* ..
188* .. Statement Functions ..
189 REAL cabs1
190* ..
191* .. Statement Function Definitions ..
192 cabs1( zdum ) = abs( real( zdum ) ) + abs( aimag( zdum ) )
193* ..
194* .. Executable Statements ..
195*
196 cla_gbrcond_x = 0.0e+0
197*
198 info = 0
199 notrans = lsame( trans, 'N' )
200 IF ( .NOT. notrans .AND. .NOT. lsame(trans, 'T') .AND. .NOT.
201 $ lsame( trans, 'C' ) ) THEN
202 info = -1
203 ELSE IF( n.LT.0 ) THEN
204 info = -2
205 ELSE IF( kl.LT.0 .OR. kl.GT.n-1 ) THEN
206 info = -3
207 ELSE IF( ku.LT.0 .OR. ku.GT.n-1 ) THEN
208 info = -4
209 ELSE IF( ldab.LT.kl+ku+1 ) THEN
210 info = -6
211 ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
212 info = -8
213 END IF
214 IF( info.NE.0 ) THEN
215 CALL xerbla( 'CLA_GBRCOND_X', -info )
216 RETURN
217 END IF
218*
219* Compute norm of op(A)*op2(C).
220*
221 kd = ku + 1
222 ke = kl + 1
223 anorm = 0.0
224 IF ( notrans ) THEN
225 DO i = 1, n
226 tmp = 0.0e+0
227 DO j = max( i-kl, 1 ), min( i+ku, n )
228 tmp = tmp + cabs1( ab( kd+i-j, j) * x( j ) )
229 END DO
230 rwork( i ) = tmp
231 anorm = max( anorm, tmp )
232 END DO
233 ELSE
234 DO i = 1, n
235 tmp = 0.0e+0
236 DO j = max( i-kl, 1 ), min( i+ku, n )
237 tmp = tmp + cabs1( ab( ke-i+j, i ) * x( j ) )
238 END DO
239 rwork( i ) = tmp
240 anorm = max( anorm, tmp )
241 END DO
242 END IF
243*
244* Quick return if possible.
245*
246 IF( n.EQ.0 ) THEN
247 cla_gbrcond_x = 1.0e+0
248 RETURN
249 ELSE IF( anorm .EQ. 0.0e+0 ) THEN
250 RETURN
251 END IF
252*
253* Estimate the norm of inv(op(A)).
254*
255 ainvnm = 0.0e+0
256*
257 kase = 0
258 10 CONTINUE
259 CALL clacn2( n, work( n+1 ), work, ainvnm, kase, isave )
260 IF( kase.NE.0 ) THEN
261 IF( kase.EQ.2 ) THEN
262*
263* Multiply by R.
264*
265 DO i = 1, n
266 work( i ) = work( i ) * rwork( i )
267 END DO
268*
269 IF ( notrans ) THEN
270 CALL cgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
271 $ ipiv, work, n, info )
272 ELSE
273 CALL cgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
274 $ ldafb, ipiv, work, n, info )
275 ENDIF
276*
277* Multiply by inv(X).
278*
279 DO i = 1, n
280 work( i ) = work( i ) / x( i )
281 END DO
282 ELSE
283*
284* Multiply by inv(X**H).
285*
286 DO i = 1, n
287 work( i ) = work( i ) / x( i )
288 END DO
289*
290 IF ( notrans ) THEN
291 CALL cgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
292 $ ldafb, ipiv, work, n, info )
293 ELSE
294 CALL cgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
295 $ ipiv, work, n, info )
296 END IF
297*
298* Multiply by R.
299*
300 DO i = 1, n
301 work( i ) = work( i ) * rwork( i )
302 END DO
303 END IF
304 GO TO 10
305 END IF
306*
307* Compute the estimate of the reciprocal condition number.
308*
309 IF( ainvnm .NE. 0.0e+0 )
310 $ cla_gbrcond_x = 1.0e+0 / ainvnm
311*
312 RETURN
313*
314* End of CLA_GBRCOND_X
315*
316 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
CGBTRS
Definition cgbtrs.f:137
real function cla_gbrcond_x(trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, x, info, work, rwork)
CLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general banded matrice...
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48