1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
|
SUBROUTINE SLA_GBAMV( TRANS, M, N, KL, KU, ALPHA, AB, LDAB, X,
$ INCX, BETA, Y, INCY )
*
* -- LAPACK routine (version 3.2) --
* -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and --
* -- Jason Riedy of Univ. of California Berkeley. --
* -- November 2008 --
*
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley and NAG Ltd. --
*
IMPLICIT NONE
* ..
* .. Scalar Arguments ..
REAL ALPHA, BETA
INTEGER INCX, INCY, LDAB, M, N, KL, KU, TRANS
* ..
* .. Array Arguments ..
REAL AB( LDAB, * ), X( * ), Y( * )
* ..
*
* Purpose
* =======
*
* SLA_GEAMV performs one of the matrix-vector operations
*
* y := alpha*abs(A)*abs(x) + beta*abs(y),
* or y := alpha*abs(A)'*abs(x) + beta*abs(y),
*
* where alpha and beta are scalars, x and y are vectors and A is an
* m by n matrix.
*
* This function is primarily used in calculating error bounds.
* To protect against underflow during evaluation, components in
* the resulting vector are perturbed away from zero by (N+1)
* times the underflow threshold. To prevent unnecessarily large
* errors for block-structure embedded in general matrices,
* "symbolically" zero components are not perturbed. A zero
* entry is considered "symbolic" if all multiplications involved
* in computing that entry have at least one zero multiplicand.
*
* Parameters
* ==========
*
* TRANS - INTEGER
* On entry, TRANS specifies the operation to be performed as
* follows:
*
* BLAS_NO_TRANS y := alpha*abs(A)*abs(x) + beta*abs(y)
* BLAS_TRANS y := alpha*abs(A')*abs(x) + beta*abs(y)
* BLAS_CONJ_TRANS y := alpha*abs(A')*abs(x) + beta*abs(y)
*
* Unchanged on exit.
*
* M - INTEGER
* On entry, M specifies the number of rows of the matrix A.
* M must be at least zero.
* Unchanged on exit.
*
* N - INTEGER
* On entry, N specifies the number of columns of the matrix A.
* N must be at least zero.
* Unchanged on exit.
*
* KL - INTEGER
* The number of subdiagonals within the band of A. KL >= 0.
*
* KU - INTEGER
* The number of superdiagonals within the band of A. KU >= 0.
*
* ALPHA - REAL
* On entry, ALPHA specifies the scalar alpha.
* Unchanged on exit.
*
* A - REAL array of DIMENSION ( LDA, n )
* Before entry, the leading m by n part of the array A must
* contain the matrix of coefficients.
* Unchanged on exit.
*
* LDA - INTEGER
* On entry, LDA specifies the first dimension of A as declared
* in the calling (sub) program. LDA must be at least
* max( 1, m ).
* Unchanged on exit.
*
* X - REAL array of DIMENSION at least
* ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n'
* and at least
* ( 1 + ( m - 1 )*abs( INCX ) ) otherwise.
* Before entry, the incremented array X must contain the
* vector x.
* Unchanged on exit.
*
* INCX - INTEGER
* On entry, INCX specifies the increment for the elements of
* X. INCX must not be zero.
* Unchanged on exit.
*
* BETA - REAL
* On entry, BETA specifies the scalar beta. When BETA is
* supplied as zero then Y need not be set on input.
* Unchanged on exit.
*
* Y - REAL array of DIMENSION at least
* ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n'
* and at least
* ( 1 + ( n - 1 )*abs( INCY ) ) otherwise.
* Before entry with BETA non-zero, the incremented array Y
* must contain the vector y. On exit, Y is overwritten by the
* updated vector y.
*
* INCY - INTEGER
* On entry, INCY specifies the increment for the elements of
* Y. INCY must not be zero.
* Unchanged on exit.
*
*
* Level 2 Blas routine.
* ..
* .. Parameters ..
REAL ONE, ZERO
PARAMETER ( ONE = 1.0E+0, ZERO = 0.0E+0 )
* ..
* .. Local Scalars ..
LOGICAL SYMB_ZERO
REAL TEMP, SAFE1
INTEGER I, INFO, IY, J, JX, KX, KY, LENX, LENY, KD
* ..
* .. External Subroutines ..
EXTERNAL XERBLA, SLAMCH
REAL SLAMCH
* ..
* .. External Functions ..
EXTERNAL ILATRANS
INTEGER ILATRANS
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, ABS, SIGN
* ..
* .. Executable Statements ..
*
* Test the input parameters.
*
INFO = 0
IF ( .NOT.( ( TRANS.EQ.ILATRANS( 'N' ) )
$ .OR. ( TRANS.EQ.ILATRANS( 'T' ) )
$ .OR. ( TRANS.EQ.ILATRANS( 'C' ) ) ) ) THEN
INFO = 1
ELSE IF( M.LT.0 )THEN
INFO = 2
ELSE IF( N.LT.0 )THEN
INFO = 3
ELSE IF( KL.LT.0 ) THEN
INFO = 4
ELSE IF( KU.LT.0 ) THEN
INFO = 5
ELSE IF( LDAB.LT.KL+KU+1 )THEN
INFO = 6
ELSE IF( INCX.EQ.0 )THEN
INFO = 8
ELSE IF( INCY.EQ.0 )THEN
INFO = 11
END IF
IF( INFO.NE.0 )THEN
CALL XERBLA( 'SLA_GBAMV ', INFO )
RETURN
END IF
*
* Quick return if possible.
*
IF( ( M.EQ.0 ).OR.( N.EQ.0 ).OR.
$ ( ( ALPHA.EQ.ZERO ).AND.( BETA.EQ.ONE ) ) )
$ RETURN
*
* Set LENX and LENY, the lengths of the vectors x and y, and set
* up the start points in X and Y.
*
IF( TRANS.EQ.ILATRANS( 'N' ) )THEN
LENX = N
LENY = M
ELSE
LENX = M
LENY = N
END IF
IF( INCX.GT.0 )THEN
KX = 1
ELSE
KX = 1 - ( LENX - 1 )*INCX
END IF
IF( INCY.GT.0 )THEN
KY = 1
ELSE
KY = 1 - ( LENY - 1 )*INCY
END IF
*
* Set SAFE1 essentially to be the underflow threshold times the
* number of additions in each row.
*
SAFE1 = SLAMCH( 'Safe minimum' )
SAFE1 = (N+1)*SAFE1
*
* Form y := alpha*abs(A)*abs(x) + beta*abs(y).
*
* The O(M*N) SYMB_ZERO tests could be replaced by O(N) queries to
* the inexact flag. Still doesn't help change the iteration order
* to per-column.
*
KD = KU + 1
IY = KY
IF ( INCX.EQ.1 ) THEN
DO I = 1, LENY
IF ( BETA .EQ. ZERO ) THEN
SYMB_ZERO = .TRUE.
Y( IY ) = 0.0
ELSE IF ( Y( IY ) .EQ. ZERO ) THEN
SYMB_ZERO = .TRUE.
ELSE
SYMB_ZERO = .FALSE.
Y( IY ) = BETA * ABS( Y( IY ) )
END IF
IF ( ALPHA .NE. ZERO ) THEN
DO J = MAX( I-KU, 1 ), MIN( I+KL, LENX )
IF( TRANS.EQ.ILATRANS( 'N' ) )THEN
TEMP = ABS( AB( KD+I-J, J ) )
ELSE
TEMP = ABS( AB( J, KD+I-J ) )
END IF
SYMB_ZERO = SYMB_ZERO .AND.
$ ( X( J ) .EQ. ZERO .OR. TEMP .EQ. ZERO )
Y( IY ) = Y( IY ) + ALPHA*ABS( X( J ) )*TEMP
END DO
END IF
IF ( .NOT.SYMB_ZERO )
$ Y( IY ) = Y( IY ) + SIGN( SAFE1, Y( IY ) )
IY = IY + INCY
END DO
ELSE
DO I = 1, LENY
IF ( BETA .EQ. ZERO ) THEN
SYMB_ZERO = .TRUE.
Y( IY ) = 0.0
ELSE IF ( Y( IY ) .EQ. ZERO ) THEN
SYMB_ZERO = .TRUE.
ELSE
SYMB_ZERO = .FALSE.
Y( IY ) = BETA * ABS( Y( IY ) )
END IF
IF ( ALPHA .NE. ZERO ) THEN
JX = KX
DO J = MAX( I-KU, 1 ), MIN( I+KL, LENX )
IF( TRANS.EQ.ILATRANS( 'N' ) )THEN
TEMP = ABS( AB( KD+I-J, J ) )
ELSE
TEMP = ABS( AB( J, KD+I-J ) )
END IF
SYMB_ZERO = SYMB_ZERO .AND.
$ ( X( JX ) .EQ. ZERO .OR. TEMP .EQ. ZERO )
Y( IY ) = Y( IY ) + ALPHA*ABS( X( JX ) )*TEMP
JX = JX + INCX
END DO
END IF
IF ( .NOT.SYMB_ZERO )
$ Y( IY ) = Y( IY ) + SIGN( SAFE1, Y( IY ) )
IY = IY + INCY
END DO
END IF
*
RETURN
*
* End of SLA_GBAMV
*
END
|