LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ctbmv.f
Go to the documentation of this file.
1*> \brief \b CTBMV
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE CTBMV(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
12*
13* .. Scalar Arguments ..
14* INTEGER INCX,K,LDA,N
15* CHARACTER DIAG,TRANS,UPLO
16* ..
17* .. Array Arguments ..
18* COMPLEX A(LDA,*),X(*)
19* ..
20*
21*
22*> \par Purpose:
23* =============
24*>
25*> \verbatim
26*>
27*> CTBMV performs one of the matrix-vector operations
28*>
29*> x := A*x, or x := A**T*x, or x := A**H*x,
30*>
31*> where x is an n element vector and A is an n by n unit, or non-unit,
32*> upper or lower triangular band matrix, with ( k + 1 ) diagonals.
33*> \endverbatim
34*
35* Arguments:
36* ==========
37*
38*> \param[in] UPLO
39*> \verbatim
40*> UPLO is CHARACTER*1
41*> On entry, UPLO specifies whether the matrix is an upper or
42*> lower triangular matrix as follows:
43*>
44*> UPLO = 'U' or 'u' A is an upper triangular matrix.
45*>
46*> UPLO = 'L' or 'l' A is a lower triangular matrix.
47*> \endverbatim
48*>
49*> \param[in] TRANS
50*> \verbatim
51*> TRANS is CHARACTER*1
52*> On entry, TRANS specifies the operation to be performed as
53*> follows:
54*>
55*> TRANS = 'N' or 'n' x := A*x.
56*>
57*> TRANS = 'T' or 't' x := A**T*x.
58*>
59*> TRANS = 'C' or 'c' x := A**H*x.
60*> \endverbatim
61*>
62*> \param[in] DIAG
63*> \verbatim
64*> DIAG is CHARACTER*1
65*> On entry, DIAG specifies whether or not A is unit
66*> triangular as follows:
67*>
68*> DIAG = 'U' or 'u' A is assumed to be unit triangular.
69*>
70*> DIAG = 'N' or 'n' A is not assumed to be unit
71*> triangular.
72*> \endverbatim
73*>
74*> \param[in] N
75*> \verbatim
76*> N is INTEGER
77*> On entry, N specifies the order of the matrix A.
78*> N must be at least zero.
79*> \endverbatim
80*>
81*> \param[in] K
82*> \verbatim
83*> K is INTEGER
84*> On entry with UPLO = 'U' or 'u', K specifies the number of
85*> super-diagonals of the matrix A.
86*> On entry with UPLO = 'L' or 'l', K specifies the number of
87*> sub-diagonals of the matrix A.
88*> K must satisfy 0 .le. K.
89*> \endverbatim
90*>
91*> \param[in] A
92*> \verbatim
93*> A is COMPLEX array, dimension ( LDA, N ).
94*> Before entry with UPLO = 'U' or 'u', the leading ( k + 1 )
95*> by n part of the array A must contain the upper triangular
96*> band part of the matrix of coefficients, supplied column by
97*> column, with the leading diagonal of the matrix in row
98*> ( k + 1 ) of the array, the first super-diagonal starting at
99*> position 2 in row k, and so on. The top left k by k triangle
100*> of the array A is not referenced.
101*> The following program segment will transfer an upper
102*> triangular band matrix from conventional full matrix storage
103*> to band storage:
104*>
105*> DO 20, J = 1, N
106*> M = K + 1 - J
107*> DO 10, I = MAX( 1, J - K ), J
108*> A( M + I, J ) = matrix( I, J )
109*> 10 CONTINUE
110*> 20 CONTINUE
111*>
112*> Before entry with UPLO = 'L' or 'l', the leading ( k + 1 )
113*> by n part of the array A must contain the lower triangular
114*> band part of the matrix of coefficients, supplied column by
115*> column, with the leading diagonal of the matrix in row 1 of
116*> the array, the first sub-diagonal starting at position 1 in
117*> row 2, and so on. The bottom right k by k triangle of the
118*> array A is not referenced.
119*> The following program segment will transfer a lower
120*> triangular band matrix from conventional full matrix storage
121*> to band storage:
122*>
123*> DO 20, J = 1, N
124*> M = 1 - J
125*> DO 10, I = J, MIN( N, J + K )
126*> A( M + I, J ) = matrix( I, J )
127*> 10 CONTINUE
128*> 20 CONTINUE
129*>
130*> Note that when DIAG = 'U' or 'u' the elements of the array A
131*> corresponding to the diagonal elements of the matrix are not
132*> referenced, but are assumed to be unity.
133*> \endverbatim
134*>
135*> \param[in] LDA
136*> \verbatim
137*> LDA is INTEGER
138*> On entry, LDA specifies the first dimension of A as declared
139*> in the calling (sub) program. LDA must be at least
140*> ( k + 1 ).
141*> \endverbatim
142*>
143*> \param[in,out] X
144*> \verbatim
145*> X is COMPLEX array, dimension at least
146*> ( 1 + ( n - 1 )*abs( INCX ) ).
147*> Before entry, the incremented array X must contain the n
148*> element vector x. On exit, X is overwritten with the
149*> transformed vector x.
150*> \endverbatim
151*>
152*> \param[in] INCX
153*> \verbatim
154*> INCX is INTEGER
155*> On entry, INCX specifies the increment for the elements of
156*> X. INCX must not be zero.
157*> \endverbatim
158*
159* Authors:
160* ========
161*
162*> \author Univ. of Tennessee
163*> \author Univ. of California Berkeley
164*> \author Univ. of Colorado Denver
165*> \author NAG Ltd.
166*
167*> \ingroup tbmv
168*
169*> \par Further Details:
170* =====================
171*>
172*> \verbatim
173*>
174*> Level 2 Blas routine.
175*> The vector and matrix arguments are not referenced when N = 0, or M = 0
176*>
177*> -- Written on 22-October-1986.
178*> Jack Dongarra, Argonne National Lab.
179*> Jeremy Du Croz, Nag Central Office.
180*> Sven Hammarling, Nag Central Office.
181*> Richard Hanson, Sandia National Labs.
182*> \endverbatim
183*>
184* =====================================================================
185 SUBROUTINE ctbmv(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
186*
187* -- Reference BLAS level2 routine --
188* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
189* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
190*
191* .. Scalar Arguments ..
192 INTEGER INCX,K,LDA,N
193 CHARACTER DIAG,TRANS,UPLO
194* ..
195* .. Array Arguments ..
196 COMPLEX A(LDA,*),X(*)
197* ..
198*
199* =====================================================================
200*
201* .. Parameters ..
202 COMPLEX ZERO
203 parameter(zero= (0.0e+0,0.0e+0))
204* ..
205* .. Local Scalars ..
206 COMPLEX TEMP
207 INTEGER I,INFO,IX,J,JX,KPLUS1,KX,L
208 LOGICAL NOCONJ,NOUNIT
209* ..
210* .. External Functions ..
211 LOGICAL LSAME
212 EXTERNAL lsame
213* ..
214* .. External Subroutines ..
215 EXTERNAL xerbla
216* ..
217* .. Intrinsic Functions ..
218 INTRINSIC conjg,max,min
219* ..
220*
221* Test the input parameters.
222*
223 info = 0
224 IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
225 info = 1
226 ELSE IF (.NOT.lsame(trans,'N') .AND.
227 + .NOT.lsame(trans,'T') .AND.
228 + .NOT.lsame(trans,'C')) THEN
229 info = 2
230 ELSE IF (.NOT.lsame(diag,'U') .AND.
231 + .NOT.lsame(diag,'N')) THEN
232 info = 3
233 ELSE IF (n.LT.0) THEN
234 info = 4
235 ELSE IF (k.LT.0) THEN
236 info = 5
237 ELSE IF (lda.LT. (k+1)) THEN
238 info = 7
239 ELSE IF (incx.EQ.0) THEN
240 info = 9
241 END IF
242 IF (info.NE.0) THEN
243 CALL xerbla('CTBMV ',info)
244 RETURN
245 END IF
246*
247* Quick return if possible.
248*
249 IF (n.EQ.0) RETURN
250*
251 noconj = lsame(trans,'T')
252 nounit = lsame(diag,'N')
253*
254* Set up the start point in X if the increment is not unity. This
255* will be ( N - 1 )*INCX too small for descending loops.
256*
257 IF (incx.LE.0) THEN
258 kx = 1 - (n-1)*incx
259 ELSE IF (incx.NE.1) THEN
260 kx = 1
261 END IF
262*
263* Start the operations. In this version the elements of A are
264* accessed sequentially with one pass through A.
265*
266 IF (lsame(trans,'N')) THEN
267*
268* Form x := A*x.
269*
270 IF (lsame(uplo,'U')) THEN
271 kplus1 = k + 1
272 IF (incx.EQ.1) THEN
273 DO 20 j = 1,n
274 IF (x(j).NE.zero) THEN
275 temp = x(j)
276 l = kplus1 - j
277 DO 10 i = max(1,j-k),j - 1
278 x(i) = x(i) + temp*a(l+i,j)
279 10 CONTINUE
280 IF (nounit) x(j) = x(j)*a(kplus1,j)
281 END IF
282 20 CONTINUE
283 ELSE
284 jx = kx
285 DO 40 j = 1,n
286 IF (x(jx).NE.zero) THEN
287 temp = x(jx)
288 ix = kx
289 l = kplus1 - j
290 DO 30 i = max(1,j-k),j - 1
291 x(ix) = x(ix) + temp*a(l+i,j)
292 ix = ix + incx
293 30 CONTINUE
294 IF (nounit) x(jx) = x(jx)*a(kplus1,j)
295 END IF
296 jx = jx + incx
297 IF (j.GT.k) kx = kx + incx
298 40 CONTINUE
299 END IF
300 ELSE
301 IF (incx.EQ.1) THEN
302 DO 60 j = n,1,-1
303 IF (x(j).NE.zero) THEN
304 temp = x(j)
305 l = 1 - j
306 DO 50 i = min(n,j+k),j + 1,-1
307 x(i) = x(i) + temp*a(l+i,j)
308 50 CONTINUE
309 IF (nounit) x(j) = x(j)*a(1,j)
310 END IF
311 60 CONTINUE
312 ELSE
313 kx = kx + (n-1)*incx
314 jx = kx
315 DO 80 j = n,1,-1
316 IF (x(jx).NE.zero) THEN
317 temp = x(jx)
318 ix = kx
319 l = 1 - j
320 DO 70 i = min(n,j+k),j + 1,-1
321 x(ix) = x(ix) + temp*a(l+i,j)
322 ix = ix - incx
323 70 CONTINUE
324 IF (nounit) x(jx) = x(jx)*a(1,j)
325 END IF
326 jx = jx - incx
327 IF ((n-j).GE.k) kx = kx - incx
328 80 CONTINUE
329 END IF
330 END IF
331 ELSE
332*
333* Form x := A**T*x or x := A**H*x.
334*
335 IF (lsame(uplo,'U')) THEN
336 kplus1 = k + 1
337 IF (incx.EQ.1) THEN
338 DO 110 j = n,1,-1
339 temp = x(j)
340 l = kplus1 - j
341 IF (noconj) THEN
342 IF (nounit) temp = temp*a(kplus1,j)
343 DO 90 i = j - 1,max(1,j-k),-1
344 temp = temp + a(l+i,j)*x(i)
345 90 CONTINUE
346 ELSE
347 IF (nounit) temp = temp*conjg(a(kplus1,j))
348 DO 100 i = j - 1,max(1,j-k),-1
349 temp = temp + conjg(a(l+i,j))*x(i)
350 100 CONTINUE
351 END IF
352 x(j) = temp
353 110 CONTINUE
354 ELSE
355 kx = kx + (n-1)*incx
356 jx = kx
357 DO 140 j = n,1,-1
358 temp = x(jx)
359 kx = kx - incx
360 ix = kx
361 l = kplus1 - j
362 IF (noconj) THEN
363 IF (nounit) temp = temp*a(kplus1,j)
364 DO 120 i = j - 1,max(1,j-k),-1
365 temp = temp + a(l+i,j)*x(ix)
366 ix = ix - incx
367 120 CONTINUE
368 ELSE
369 IF (nounit) temp = temp*conjg(a(kplus1,j))
370 DO 130 i = j - 1,max(1,j-k),-1
371 temp = temp + conjg(a(l+i,j))*x(ix)
372 ix = ix - incx
373 130 CONTINUE
374 END IF
375 x(jx) = temp
376 jx = jx - incx
377 140 CONTINUE
378 END IF
379 ELSE
380 IF (incx.EQ.1) THEN
381 DO 170 j = 1,n
382 temp = x(j)
383 l = 1 - j
384 IF (noconj) THEN
385 IF (nounit) temp = temp*a(1,j)
386 DO 150 i = j + 1,min(n,j+k)
387 temp = temp + a(l+i,j)*x(i)
388 150 CONTINUE
389 ELSE
390 IF (nounit) temp = temp*conjg(a(1,j))
391 DO 160 i = j + 1,min(n,j+k)
392 temp = temp + conjg(a(l+i,j))*x(i)
393 160 CONTINUE
394 END IF
395 x(j) = temp
396 170 CONTINUE
397 ELSE
398 jx = kx
399 DO 200 j = 1,n
400 temp = x(jx)
401 kx = kx + incx
402 ix = kx
403 l = 1 - j
404 IF (noconj) THEN
405 IF (nounit) temp = temp*a(1,j)
406 DO 180 i = j + 1,min(n,j+k)
407 temp = temp + a(l+i,j)*x(ix)
408 ix = ix + incx
409 180 CONTINUE
410 ELSE
411 IF (nounit) temp = temp*conjg(a(1,j))
412 DO 190 i = j + 1,min(n,j+k)
413 temp = temp + conjg(a(l+i,j))*x(ix)
414 ix = ix + incx
415 190 CONTINUE
416 END IF
417 x(jx) = temp
418 jx = jx + incx
419 200 CONTINUE
420 END IF
421 END IF
422 END IF
423*
424 RETURN
425*
426* End of CTBMV
427*
428 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ctbmv(uplo, trans, diag, n, k, a, lda, x, incx)
CTBMV
Definition ctbmv.f:186