LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dtbmv.f
Go to the documentation of this file.
1*> \brief \b DTBMV
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 DTBMV(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* DOUBLE PRECISION A(LDA,*),X(*)
19* ..
20*
21*
22*> \par Purpose:
23* =============
24*>
25*> \verbatim
26*>
27*> DTBMV performs one of the matrix-vector operations
28*>
29*> x := A*x, or x := A**T*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**T*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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 dtbmv(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 DOUBLE PRECISION A(LDA,*),X(*)
197* ..
198*
199* =====================================================================
200*
201* .. Parameters ..
202 DOUBLE PRECISION ZERO
203 parameter(zero=0.0d+0)
204* ..
205* .. Local Scalars ..
206 DOUBLE PRECISION TEMP
207 INTEGER I,INFO,IX,J,JX,KPLUS1,KX,L
208 LOGICAL NOUNIT
209* ..
210* .. External Functions ..
211 LOGICAL LSAME
212 EXTERNAL lsame
213* ..
214* .. External Subroutines ..
215 EXTERNAL xerbla
216* ..
217* .. Intrinsic Functions ..
218 INTRINSIC 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('DTBMV ',info)
244 RETURN
245 END IF
246*
247* Quick return if possible.
248*
249 IF (n.EQ.0) RETURN
250*
251 nounit = lsame(diag,'N')
252*
253* Set up the start point in X if the increment is not unity. This
254* will be ( N - 1 )*INCX too small for descending loops.
255*
256 IF (incx.LE.0) THEN
257 kx = 1 - (n-1)*incx
258 ELSE IF (incx.NE.1) THEN
259 kx = 1
260 END IF
261*
262* Start the operations. In this version the elements of A are
263* accessed sequentially with one pass through A.
264*
265 IF (lsame(trans,'N')) THEN
266*
267* Form x := A*x.
268*
269 IF (lsame(uplo,'U')) THEN
270 kplus1 = k + 1
271 IF (incx.EQ.1) THEN
272 DO 20 j = 1,n
273 IF (x(j).NE.zero) THEN
274 temp = x(j)
275 l = kplus1 - j
276 DO 10 i = max(1,j-k),j - 1
277 x(i) = x(i) + temp*a(l+i,j)
278 10 CONTINUE
279 IF (nounit) x(j) = x(j)*a(kplus1,j)
280 END IF
281 20 CONTINUE
282 ELSE
283 jx = kx
284 DO 40 j = 1,n
285 IF (x(jx).NE.zero) THEN
286 temp = x(jx)
287 ix = kx
288 l = kplus1 - j
289 DO 30 i = max(1,j-k),j - 1
290 x(ix) = x(ix) + temp*a(l+i,j)
291 ix = ix + incx
292 30 CONTINUE
293 IF (nounit) x(jx) = x(jx)*a(kplus1,j)
294 END IF
295 jx = jx + incx
296 IF (j.GT.k) kx = kx + incx
297 40 CONTINUE
298 END IF
299 ELSE
300 IF (incx.EQ.1) THEN
301 DO 60 j = n,1,-1
302 IF (x(j).NE.zero) THEN
303 temp = x(j)
304 l = 1 - j
305 DO 50 i = min(n,j+k),j + 1,-1
306 x(i) = x(i) + temp*a(l+i,j)
307 50 CONTINUE
308 IF (nounit) x(j) = x(j)*a(1,j)
309 END IF
310 60 CONTINUE
311 ELSE
312 kx = kx + (n-1)*incx
313 jx = kx
314 DO 80 j = n,1,-1
315 IF (x(jx).NE.zero) THEN
316 temp = x(jx)
317 ix = kx
318 l = 1 - j
319 DO 70 i = min(n,j+k),j + 1,-1
320 x(ix) = x(ix) + temp*a(l+i,j)
321 ix = ix - incx
322 70 CONTINUE
323 IF (nounit) x(jx) = x(jx)*a(1,j)
324 END IF
325 jx = jx - incx
326 IF ((n-j).GE.k) kx = kx - incx
327 80 CONTINUE
328 END IF
329 END IF
330 ELSE
331*
332* Form x := A**T*x.
333*
334 IF (lsame(uplo,'U')) THEN
335 kplus1 = k + 1
336 IF (incx.EQ.1) THEN
337 DO 100 j = n,1,-1
338 temp = x(j)
339 l = kplus1 - j
340 IF (nounit) temp = temp*a(kplus1,j)
341 DO 90 i = j - 1,max(1,j-k),-1
342 temp = temp + a(l+i,j)*x(i)
343 90 CONTINUE
344 x(j) = temp
345 100 CONTINUE
346 ELSE
347 kx = kx + (n-1)*incx
348 jx = kx
349 DO 120 j = n,1,-1
350 temp = x(jx)
351 kx = kx - incx
352 ix = kx
353 l = kplus1 - j
354 IF (nounit) temp = temp*a(kplus1,j)
355 DO 110 i = j - 1,max(1,j-k),-1
356 temp = temp + a(l+i,j)*x(ix)
357 ix = ix - incx
358 110 CONTINUE
359 x(jx) = temp
360 jx = jx - incx
361 120 CONTINUE
362 END IF
363 ELSE
364 IF (incx.EQ.1) THEN
365 DO 140 j = 1,n
366 temp = x(j)
367 l = 1 - j
368 IF (nounit) temp = temp*a(1,j)
369 DO 130 i = j + 1,min(n,j+k)
370 temp = temp + a(l+i,j)*x(i)
371 130 CONTINUE
372 x(j) = temp
373 140 CONTINUE
374 ELSE
375 jx = kx
376 DO 160 j = 1,n
377 temp = x(jx)
378 kx = kx + incx
379 ix = kx
380 l = 1 - j
381 IF (nounit) temp = temp*a(1,j)
382 DO 150 i = j + 1,min(n,j+k)
383 temp = temp + a(l+i,j)*x(ix)
384 ix = ix + incx
385 150 CONTINUE
386 x(jx) = temp
387 jx = jx + incx
388 160 CONTINUE
389 END IF
390 END IF
391 END IF
392*
393 RETURN
394*
395* End of DTBMV
396*
397 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dtbmv(uplo, trans, diag, n, k, a, lda, x, incx)
DTBMV
Definition dtbmv.f:186