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