LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 of 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 of 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 *> tranformed 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 *> \date November 2011
168 *
169 *> \ingroup double_blas_level2
170 *
171 *> \par Further Details:
172 * =====================
173 *>
174 *> \verbatim
175 *>
176 *> Level 2 Blas routine.
177 *> The vector and matrix arguments are not referenced when N = 0, or M = 0
178 *>
179 *> -- Written on 22-October-1986.
180 *> Jack Dongarra, Argonne National Lab.
181 *> Jeremy Du Croz, Nag Central Office.
182 *> Sven Hammarling, Nag Central Office.
183 *> Richard Hanson, Sandia National Labs.
184 *> \endverbatim
185 *>
186 * =====================================================================
187  SUBROUTINE dtbmv(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
188 *
189 * -- Reference BLAS level2 routine (version 3.4.0) --
190 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
191 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
192 * November 2011
193 *
194 * .. Scalar Arguments ..
195  INTEGER INCX,K,LDA,N
196  CHARACTER DIAG,TRANS,UPLO
197 * ..
198 * .. Array Arguments ..
199  DOUBLE PRECISION A(lda,*),X(*)
200 * ..
201 *
202 * =====================================================================
203 *
204 * .. Parameters ..
205  DOUBLE PRECISION ZERO
206  parameter(zero=0.0d+0)
207 * ..
208 * .. Local Scalars ..
209  DOUBLE PRECISION TEMP
210  INTEGER I,INFO,IX,J,JX,KPLUS1,KX,L
211  LOGICAL NOUNIT
212 * ..
213 * .. External Functions ..
214  LOGICAL LSAME
215  EXTERNAL lsame
216 * ..
217 * .. External Subroutines ..
218  EXTERNAL xerbla
219 * ..
220 * .. Intrinsic Functions ..
221  INTRINSIC max,min
222 * ..
223 *
224 * Test the input parameters.
225 *
226  info = 0
227  IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
228  info = 1
229  ELSE IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
230  + .NOT.lsame(trans,'C')) THEN
231  info = 2
232  ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
233  info = 3
234  ELSE IF (n.LT.0) THEN
235  info = 4
236  ELSE IF (k.LT.0) THEN
237  info = 5
238  ELSE IF (lda.LT. (k+1)) THEN
239  info = 7
240  ELSE IF (incx.EQ.0) THEN
241  info = 9
242  END IF
243  IF (info.NE.0) THEN
244  CALL xerbla('DTBMV ',info)
245  RETURN
246  END IF
247 *
248 * Quick return if possible.
249 *
250  IF (n.EQ.0) RETURN
251 *
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.
334 *
335  IF (lsame(uplo,'U')) THEN
336  kplus1 = k + 1
337  IF (incx.EQ.1) THEN
338  DO 100 j = n,1,-1
339  temp = x(j)
340  l = kplus1 - j
341  IF (nounit) temp = temp*a(kplus1,j)
342  DO 90 i = j - 1,max(1,j-k),-1
343  temp = temp + a(l+i,j)*x(i)
344  90 CONTINUE
345  x(j) = temp
346  100 CONTINUE
347  ELSE
348  kx = kx + (n-1)*incx
349  jx = kx
350  DO 120 j = n,1,-1
351  temp = x(jx)
352  kx = kx - incx
353  ix = kx
354  l = kplus1 - j
355  IF (nounit) temp = temp*a(kplus1,j)
356  DO 110 i = j - 1,max(1,j-k),-1
357  temp = temp + a(l+i,j)*x(ix)
358  ix = ix - incx
359  110 CONTINUE
360  x(jx) = temp
361  jx = jx - incx
362  120 CONTINUE
363  END IF
364  ELSE
365  IF (incx.EQ.1) THEN
366  DO 140 j = 1,n
367  temp = x(j)
368  l = 1 - j
369  IF (nounit) temp = temp*a(1,j)
370  DO 130 i = j + 1,min(n,j+k)
371  temp = temp + a(l+i,j)*x(i)
372  130 CONTINUE
373  x(j) = temp
374  140 CONTINUE
375  ELSE
376  jx = kx
377  DO 160 j = 1,n
378  temp = x(jx)
379  kx = kx + incx
380  ix = kx
381  l = 1 - j
382  IF (nounit) temp = temp*a(1,j)
383  DO 150 i = j + 1,min(n,j+k)
384  temp = temp + a(l+i,j)*x(ix)
385  ix = ix + incx
386  150 CONTINUE
387  x(jx) = temp
388  jx = jx + incx
389  160 CONTINUE
390  END IF
391  END IF
392  END IF
393 *
394  RETURN
395 *
396 * End of DTBMV .
397 *
398  END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dtbmv(UPLO, TRANS, DIAG, N, K, A, LDA, X, INCX)
DTBMV
Definition: dtbmv.f:188