LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zla_gbrcond_c.f
Go to the documentation of this file.
1*> \brief \b ZLA_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*> Download ZLA_GBRCOND_C + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zla_gbrcond_c.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zla_gbrcond_c.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zla_gbrcond_c.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* DOUBLE PRECISION FUNCTION ZLA_GBRCOND_C( TRANS, N, KL, KU, AB,
20* LDAB, AFB, LDAFB, IPIV,
21* C, CAPPLY, INFO, WORK,
22* RWORK )
23*
24* .. Scalar Arguments ..
25* CHARACTER TRANS
26* LOGICAL CAPPLY
27* INTEGER N, KL, KU, KD, KE, LDAB, LDAFB, INFO
28* ..
29* .. Array Arguments ..
30* INTEGER IPIV( * )
31* COMPLEX*16 AB( LDAB, * ), AFB( LDAFB, * ), WORK( * )
32* DOUBLE PRECISION C( * ), RWORK( * )
33*
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> ZLA_GBRCOND_C Computes the infinity norm condition number of
42*> op(A) * inv(diag(C)) where C is a DOUBLE PRECISION vector.
43*> \endverbatim
44*
45* Arguments:
46* ==========
47*
48*> \param[in] TRANS
49*> \verbatim
50*> TRANS is CHARACTER*1
51*> Specifies the form of the system of equations:
52*> = 'N': A * X = B (No transpose)
53*> = 'T': A**T * X = B (Transpose)
54*> = 'C': A**H * X = B (Conjugate Transpose = Transpose)
55*> \endverbatim
56*>
57*> \param[in] N
58*> \verbatim
59*> N is INTEGER
60*> The number of linear equations, i.e., the order of the
61*> matrix A. N >= 0.
62*> \endverbatim
63*>
64*> \param[in] KL
65*> \verbatim
66*> KL is INTEGER
67*> The number of subdiagonals within the band of A. KL >= 0.
68*> \endverbatim
69*>
70*> \param[in] KU
71*> \verbatim
72*> KU is INTEGER
73*> The number of superdiagonals within the band of A. KU >= 0.
74*> \endverbatim
75*>
76*> \param[in] AB
77*> \verbatim
78*> AB is COMPLEX*16 array, dimension (LDAB,N)
79*> On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
80*> The j-th column of A is stored in the j-th column of the
81*> array AB as follows:
82*> AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
83*> \endverbatim
84*>
85*> \param[in] LDAB
86*> \verbatim
87*> LDAB is INTEGER
88*> The leading dimension of the array AB. LDAB >= KL+KU+1.
89*> \endverbatim
90*>
91*> \param[in] AFB
92*> \verbatim
93*> AFB is COMPLEX*16 array, dimension (LDAFB,N)
94*> Details of the LU factorization of the band matrix A, as
95*> computed by ZGBTRF. U is stored as an upper triangular
96*> band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
97*> and the multipliers used during the factorization are stored
98*> in rows KL+KU+2 to 2*KL+KU+1.
99*> \endverbatim
100*>
101*> \param[in] LDAFB
102*> \verbatim
103*> LDAFB is INTEGER
104*> The leading dimension of the array AFB. LDAFB >= 2*KL+KU+1.
105*> \endverbatim
106*>
107*> \param[in] IPIV
108*> \verbatim
109*> IPIV is INTEGER array, dimension (N)
110*> The pivot indices from the factorization A = P*L*U
111*> as computed by ZGBTRF; row i of the matrix was interchanged
112*> with row IPIV(i).
113*> \endverbatim
114*>
115*> \param[in] C
116*> \verbatim
117*> C is DOUBLE PRECISION array, dimension (N)
118*> The vector C in the formula op(A) * inv(diag(C)).
119*> \endverbatim
120*>
121*> \param[in] CAPPLY
122*> \verbatim
123*> CAPPLY is LOGICAL
124*> If .TRUE. then access the vector C in the formula above.
125*> \endverbatim
126*>
127*> \param[out] INFO
128*> \verbatim
129*> INFO is INTEGER
130*> = 0: Successful exit.
131*> i > 0: The ith argument is invalid.
132*> \endverbatim
133*>
134*> \param[out] WORK
135*> \verbatim
136*> WORK is COMPLEX*16 array, dimension (2*N).
137*> Workspace.
138*> \endverbatim
139*>
140*> \param[out] RWORK
141*> \verbatim
142*> RWORK is DOUBLE PRECISION array, dimension (N).
143*> Workspace.
144*> \endverbatim
145*
146* Authors:
147* ========
148*
149*> \author Univ. of Tennessee
150*> \author Univ. of California Berkeley
151*> \author Univ. of Colorado Denver
152*> \author NAG Ltd.
153*
154*> \ingroup la_gbrcond
155*
156* =====================================================================
157 DOUBLE PRECISION FUNCTION zla_gbrcond_c( TRANS, N, KL, KU, AB,
158 $ LDAB, AFB, LDAFB, IPIV,
159 $ 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*16 ab( ldab, * ), afb( ldafb, * ), work( * )
174 DOUBLE PRECISION c( * ), rwork( * )
175*
176*
177* =====================================================================
178*
179* .. Local Scalars ..
180 LOGICAL notrans
181 INTEGER kase, i, j
182 DOUBLE PRECISION ainvnm, anorm, tmp
183 COMPLEX*16 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 zlacn2, zgbtrs, xerbla
194* ..
195* .. Intrinsic Functions ..
196 INTRINSIC abs, max
197* ..
198* .. Statement Functions ..
199 DOUBLE PRECISION cabs1
200* ..
201* .. Statement Function Definitions ..
202 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
203* ..
204* .. Executable Statements ..
205 zla_gbrcond_c = 0.0d+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( 'ZLA_GBRCOND_C', -info )
225 RETURN
226 END IF
227*
228* Compute norm of op(A)*op2(C).
229*
230 anorm = 0.0d+0
231 kd = ku + 1
232 ke = kl + 1
233 IF ( notrans ) THEN
234 DO i = 1, n
235 tmp = 0.0d+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.0d+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 zla_gbrcond_c = 1.0d+0
269 RETURN
270 ELSE IF( anorm .EQ. 0.0d+0 ) THEN
271 RETURN
272 END IF
273*
274* Estimate the norm of inv(op(A)).
275*
276 ainvnm = 0.0d+0
277*
278 kase = 0
279 10 CONTINUE
280 CALL zlacn2( 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 zgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
292 $ ipiv, work, n, info )
293 ELSE
294 CALL zgbtrs( '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 zgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
317 $ ldafb, ipiv, work, n, info )
318 ELSE
319 CALL zgbtrs( '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.0d+0 )
335 $ zla_gbrcond_c = 1.0d+0 / ainvnm
336*
337 RETURN
338*
339* End of ZLA_GBRCOND_C
340*
341 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
ZGBTRS
Definition zgbtrs.f:137
double precision function zla_gbrcond_c(trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, c, capply, info, work, rwork)
ZLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general banded ma...
subroutine zlacn2(n, v, x, est, kase, isave)
ZLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition zlacn2.f:131
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48