LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgbmv.f
Go to the documentation of this file.
1*> \brief \b DGBMV
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 DGBMV(TRANS,M,N,KL,KU,ALPHA,A,LDA,X,INCX,BETA,Y,INCY)
12*
13* .. Scalar Arguments ..
14* DOUBLE PRECISION ALPHA,BETA
15* INTEGER INCX,INCY,KL,KU,LDA,M,N
16* CHARACTER TRANS
17* ..
18* .. Array Arguments ..
19* DOUBLE PRECISION A(LDA,*),X(*),Y(*)
20* ..
21*
22*
23*> \par Purpose:
24* =============
25*>
26*> \verbatim
27*>
28*> DGBMV performs one of the matrix-vector operations
29*>
30*> y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y,
31*>
32*> where alpha and beta are scalars, x and y are vectors and A is an
33*> m by n band matrix, with kl sub-diagonals and ku super-diagonals.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] TRANS
40*> \verbatim
41*> TRANS is CHARACTER*1
42*> On entry, TRANS specifies the operation to be performed as
43*> follows:
44*>
45*> TRANS = 'N' or 'n' y := alpha*A*x + beta*y.
46*>
47*> TRANS = 'T' or 't' y := alpha*A**T*x + beta*y.
48*>
49*> TRANS = 'C' or 'c' y := alpha*A**T*x + beta*y.
50*> \endverbatim
51*>
52*> \param[in] M
53*> \verbatim
54*> M is INTEGER
55*> On entry, M specifies the number of rows of the matrix A.
56*> M must be at least zero.
57*> \endverbatim
58*>
59*> \param[in] N
60*> \verbatim
61*> N is INTEGER
62*> On entry, N specifies the number of columns of the matrix A.
63*> N must be at least zero.
64*> \endverbatim
65*>
66*> \param[in] KL
67*> \verbatim
68*> KL is INTEGER
69*> On entry, KL specifies the number of sub-diagonals of the
70*> matrix A. KL must satisfy 0 .le. KL.
71*> \endverbatim
72*>
73*> \param[in] KU
74*> \verbatim
75*> KU is INTEGER
76*> On entry, KU specifies the number of super-diagonals of the
77*> matrix A. KU must satisfy 0 .le. KU.
78*> \endverbatim
79*>
80*> \param[in] ALPHA
81*> \verbatim
82*> ALPHA is DOUBLE PRECISION.
83*> On entry, ALPHA specifies the scalar alpha.
84*> \endverbatim
85*>
86*> \param[in] A
87*> \verbatim
88*> A is DOUBLE PRECISION array, dimension ( LDA, N )
89*> Before entry, the leading ( kl + ku + 1 ) by n part of the
90*> array A must contain the matrix of coefficients, supplied
91*> column by column, with the leading diagonal of the matrix in
92*> row ( ku + 1 ) of the array, the first super-diagonal
93*> starting at position 2 in row ku, the first sub-diagonal
94*> starting at position 1 in row ( ku + 2 ), and so on.
95*> Elements in the array A that do not correspond to elements
96*> in the band matrix (such as the top left ku by ku triangle)
97*> are not referenced.
98*> The following program segment will transfer a band matrix
99*> from conventional full matrix storage to band storage:
100*>
101*> DO 20, J = 1, N
102*> K = KU + 1 - J
103*> DO 10, I = MAX( 1, J - KU ), MIN( M, J + KL )
104*> A( K + I, J ) = matrix( I, J )
105*> 10 CONTINUE
106*> 20 CONTINUE
107*> \endverbatim
108*>
109*> \param[in] LDA
110*> \verbatim
111*> LDA is INTEGER
112*> On entry, LDA specifies the first dimension of A as declared
113*> in the calling (sub) program. LDA must be at least
114*> ( kl + ku + 1 ).
115*> \endverbatim
116*>
117*> \param[in] X
118*> \verbatim
119*> X is DOUBLE PRECISION array, dimension at least
120*> ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n'
121*> and at least
122*> ( 1 + ( m - 1 )*abs( INCX ) ) otherwise.
123*> Before entry, the incremented array X must contain the
124*> vector x.
125*> \endverbatim
126*>
127*> \param[in] INCX
128*> \verbatim
129*> INCX is INTEGER
130*> On entry, INCX specifies the increment for the elements of
131*> X. INCX must not be zero.
132*> \endverbatim
133*>
134*> \param[in] BETA
135*> \verbatim
136*> BETA is DOUBLE PRECISION.
137*> On entry, BETA specifies the scalar beta. When BETA is
138*> supplied as zero then Y need not be set on input.
139*> \endverbatim
140*>
141*> \param[in,out] Y
142*> \verbatim
143*> Y is DOUBLE PRECISION array, dimension at least
144*> ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n'
145*> and at least
146*> ( 1 + ( n - 1 )*abs( INCY ) ) otherwise.
147*> Before entry, the incremented array Y must contain the
148*> vector y. On exit, Y is overwritten by the updated vector y.
149*> \endverbatim
150*>
151*> \param[in] INCY
152*> \verbatim
153*> INCY is INTEGER
154*> On entry, INCY specifies the increment for the elements of
155*> Y. INCY must not be zero.
156*> \endverbatim
157*
158* Authors:
159* ========
160*
161*> \author Univ. of Tennessee
162*> \author Univ. of California Berkeley
163*> \author Univ. of Colorado Denver
164*> \author NAG Ltd.
165*
166*> \ingroup double_blas_level2
167*
168*> \par Further Details:
169* =====================
170*>
171*> \verbatim
172*>
173*> Level 2 Blas routine.
174*> The vector and matrix arguments are not referenced when N = 0, or M = 0
175*>
176*> -- Written on 22-October-1986.
177*> Jack Dongarra, Argonne National Lab.
178*> Jeremy Du Croz, Nag Central Office.
179*> Sven Hammarling, Nag Central Office.
180*> Richard Hanson, Sandia National Labs.
181*> \endverbatim
182*>
183* =====================================================================
184 SUBROUTINE dgbmv(TRANS,M,N,KL,KU,ALPHA,A,LDA,X,INCX,BETA,Y,INCY)
185*
186* -- Reference BLAS level2 routine --
187* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
188* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189*
190* .. Scalar Arguments ..
191 DOUBLE PRECISION ALPHA,BETA
192 INTEGER INCX,INCY,KL,KU,LDA,M,N
193 CHARACTER TRANS
194* ..
195* .. Array Arguments ..
196 DOUBLE PRECISION A(LDA,*),X(*),Y(*)
197* ..
198*
199* =====================================================================
200*
201* .. Parameters ..
202 DOUBLE PRECISION ONE,ZERO
203 parameter(one=1.0d+0,zero=0.0d+0)
204* ..
205* .. Local Scalars ..
206 DOUBLE PRECISION TEMP
207 INTEGER I,INFO,IX,IY,J,JX,JY,K,KUP1,KX,KY,LENX,LENY
208* ..
209* .. External Functions ..
210 LOGICAL LSAME
211 EXTERNAL lsame
212* ..
213* .. External Subroutines ..
214 EXTERNAL xerbla
215* ..
216* .. Intrinsic Functions ..
217 INTRINSIC max,min
218* ..
219*
220* Test the input parameters.
221*
222 info = 0
223 IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
224 + .NOT.lsame(trans,'C')) THEN
225 info = 1
226 ELSE IF (m.LT.0) THEN
227 info = 2
228 ELSE IF (n.LT.0) THEN
229 info = 3
230 ELSE IF (kl.LT.0) THEN
231 info = 4
232 ELSE IF (ku.LT.0) THEN
233 info = 5
234 ELSE IF (lda.LT. (kl+ku+1)) THEN
235 info = 8
236 ELSE IF (incx.EQ.0) THEN
237 info = 10
238 ELSE IF (incy.EQ.0) THEN
239 info = 13
240 END IF
241 IF (info.NE.0) THEN
242 CALL xerbla('DGBMV ',info)
243 RETURN
244 END IF
245*
246* Quick return if possible.
247*
248 IF ((m.EQ.0) .OR. (n.EQ.0) .OR.
249 + ((alpha.EQ.zero).AND. (beta.EQ.one))) RETURN
250*
251* Set LENX and LENY, the lengths of the vectors x and y, and set
252* up the start points in X and Y.
253*
254 IF (lsame(trans,'N')) THEN
255 lenx = n
256 leny = m
257 ELSE
258 lenx = m
259 leny = n
260 END IF
261 IF (incx.GT.0) THEN
262 kx = 1
263 ELSE
264 kx = 1 - (lenx-1)*incx
265 END IF
266 IF (incy.GT.0) THEN
267 ky = 1
268 ELSE
269 ky = 1 - (leny-1)*incy
270 END IF
271*
272* Start the operations. In this version the elements of A are
273* accessed sequentially with one pass through the band part of A.
274*
275* First form y := beta*y.
276*
277 IF (beta.NE.one) THEN
278 IF (incy.EQ.1) THEN
279 IF (beta.EQ.zero) THEN
280 DO 10 i = 1,leny
281 y(i) = zero
282 10 CONTINUE
283 ELSE
284 DO 20 i = 1,leny
285 y(i) = beta*y(i)
286 20 CONTINUE
287 END IF
288 ELSE
289 iy = ky
290 IF (beta.EQ.zero) THEN
291 DO 30 i = 1,leny
292 y(iy) = zero
293 iy = iy + incy
294 30 CONTINUE
295 ELSE
296 DO 40 i = 1,leny
297 y(iy) = beta*y(iy)
298 iy = iy + incy
299 40 CONTINUE
300 END IF
301 END IF
302 END IF
303 IF (alpha.EQ.zero) RETURN
304 kup1 = ku + 1
305 IF (lsame(trans,'N')) THEN
306*
307* Form y := alpha*A*x + y.
308*
309 jx = kx
310 IF (incy.EQ.1) THEN
311 DO 60 j = 1,n
312 temp = alpha*x(jx)
313 k = kup1 - j
314 DO 50 i = max(1,j-ku),min(m,j+kl)
315 y(i) = y(i) + temp*a(k+i,j)
316 50 CONTINUE
317 jx = jx + incx
318 60 CONTINUE
319 ELSE
320 DO 80 j = 1,n
321 temp = alpha*x(jx)
322 iy = ky
323 k = kup1 - j
324 DO 70 i = max(1,j-ku),min(m,j+kl)
325 y(iy) = y(iy) + temp*a(k+i,j)
326 iy = iy + incy
327 70 CONTINUE
328 jx = jx + incx
329 IF (j.GT.ku) ky = ky + incy
330 80 CONTINUE
331 END IF
332 ELSE
333*
334* Form y := alpha*A**T*x + y.
335*
336 jy = ky
337 IF (incx.EQ.1) THEN
338 DO 100 j = 1,n
339 temp = zero
340 k = kup1 - j
341 DO 90 i = max(1,j-ku),min(m,j+kl)
342 temp = temp + a(k+i,j)*x(i)
343 90 CONTINUE
344 y(jy) = y(jy) + alpha*temp
345 jy = jy + incy
346 100 CONTINUE
347 ELSE
348 DO 120 j = 1,n
349 temp = zero
350 ix = kx
351 k = kup1 - j
352 DO 110 i = max(1,j-ku),min(m,j+kl)
353 temp = temp + a(k+i,j)*x(ix)
354 ix = ix + incx
355 110 CONTINUE
356 y(jy) = y(jy) + alpha*temp
357 jy = jy + incy
358 IF (j.GT.ku) kx = kx + incx
359 120 CONTINUE
360 END IF
361 END IF
362*
363 RETURN
364*
365* End of DGBMV
366*
367 END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine dgbmv(TRANS, M, N, KL, KU, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGBMV
Definition: dgbmv.f:185