LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zhemv()

subroutine zhemv ( character uplo,
integer n,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(*) x,
integer incx,
complex*16 beta,
complex*16, dimension(*) y,
integer incy )

ZHEMV

Purpose:
!>
!> ZHEMV  performs the matrix-vector  operation
!>
!>    y := alpha*A*x + beta*y,
!>
!> where alpha and beta are scalars, x and y are n element vectors and
!> A is an n by n hermitian matrix.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On entry, UPLO specifies whether the upper or lower
!>           triangular part of the array A is to be referenced as
!>           follows:
!>
!>              UPLO = 'U' or 'u'   Only the upper triangular part of A
!>                                  is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the lower triangular part of A
!>                                  is to be referenced.
!> 
[in]N
!>          N is INTEGER
!>           On entry, N specifies the order of the matrix A.
!>           N must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, N )
!>           Before entry with  UPLO = 'U' or 'u', the leading n by n
!>           upper triangular part of the array A must contain the upper
!>           triangular part of the hermitian matrix and the strictly
!>           lower triangular part of A is not referenced.
!>           Before entry with UPLO = 'L' or 'l', the leading n by n
!>           lower triangular part of the array A must contain the lower
!>           triangular part of the hermitian matrix and the strictly
!>           upper triangular part of A is not referenced.
!>           Note that the imaginary parts of the diagonal elements need
!>           not be set and are assumed to be zero.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in the calling (sub) program. LDA must be at least
!>           max( 1, n ).
!> 
[in]X
!>          X is COMPLEX*16 array, dimension at least
!>           ( 1 + ( n - 1 )*abs( INCX ) ).
!>           Before entry, the incremented array X must contain the n
!>           element vector x.
!> 
[in]INCX
!>          INCX is INTEGER
!>           On entry, INCX specifies the increment for the elements of
!>           X. INCX must not be zero.
!> 
[in]BETA
!>          BETA is COMPLEX*16
!>           On entry, BETA specifies the scalar beta. When BETA is
!>           supplied as zero then Y need not be set on input.
!> 
[in,out]Y
!>          Y is COMPLEX*16 array, dimension at least
!>           ( 1 + ( n - 1 )*abs( INCY ) ).
!>           Before entry, the incremented array Y must contain the n
!>           element vector y. On exit, Y is overwritten by the updated
!>           vector y.
!> 
[in]INCY
!>          INCY is INTEGER
!>           On entry, INCY specifies the increment for the elements of
!>           Y. INCY must not be zero.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 2 Blas routine.
!>  The vector and matrix arguments are not referenced when N = 0, or M = 0
!>
!>  -- Written on 22-October-1986.
!>     Jack Dongarra, Argonne National Lab.
!>     Jeremy Du Croz, Nag Central Office.
!>     Sven Hammarling, Nag Central Office.
!>     Richard Hanson, Sandia National Labs.
!> 

Definition at line 153 of file zhemv.f.

154*
155* -- Reference BLAS level2 routine --
156* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
157* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
158*
159* .. Scalar Arguments ..
160 COMPLEX*16 ALPHA,BETA
161 INTEGER INCX,INCY,LDA,N
162 CHARACTER UPLO
163* ..
164* .. Array Arguments ..
165 COMPLEX*16 A(LDA,*),X(*),Y(*)
166* ..
167*
168* =====================================================================
169*
170* .. Parameters ..
171 COMPLEX*16 ONE
172 parameter(one= (1.0d+0,0.0d+0))
173 COMPLEX*16 ZERO
174 parameter(zero= (0.0d+0,0.0d+0))
175* ..
176* .. Local Scalars ..
177 COMPLEX*16 TEMP1,TEMP2
178 INTEGER I,INFO,IX,IY,J,JX,JY,KX,KY
179* ..
180* .. External Functions ..
181 LOGICAL LSAME
182 EXTERNAL lsame
183* ..
184* .. External Subroutines ..
185 EXTERNAL xerbla
186* ..
187* .. Intrinsic Functions ..
188 INTRINSIC dble,dconjg,max
189* ..
190*
191* Test the input parameters.
192*
193 info = 0
194 IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
195 info = 1
196 ELSE IF (n.LT.0) THEN
197 info = 2
198 ELSE IF (lda.LT.max(1,n)) THEN
199 info = 5
200 ELSE IF (incx.EQ.0) THEN
201 info = 7
202 ELSE IF (incy.EQ.0) THEN
203 info = 10
204 END IF
205 IF (info.NE.0) THEN
206 CALL xerbla('ZHEMV ',info)
207 RETURN
208 END IF
209*
210* Quick return if possible.
211*
212 IF ((n.EQ.0) .OR. ((alpha.EQ.zero).AND. (beta.EQ.one))) RETURN
213*
214* Set up the start points in X and Y.
215*
216 IF (incx.GT.0) THEN
217 kx = 1
218 ELSE
219 kx = 1 - (n-1)*incx
220 END IF
221 IF (incy.GT.0) THEN
222 ky = 1
223 ELSE
224 ky = 1 - (n-1)*incy
225 END IF
226*
227* Start the operations. In this version the elements of A are
228* accessed sequentially with one pass through the triangular part
229* of A.
230*
231* First form y := beta*y.
232*
233 IF (beta.NE.one) THEN
234 IF (incy.EQ.1) THEN
235 IF (beta.EQ.zero) THEN
236 DO 10 i = 1,n
237 y(i) = zero
238 10 CONTINUE
239 ELSE
240 DO 20 i = 1,n
241 y(i) = beta*y(i)
242 20 CONTINUE
243 END IF
244 ELSE
245 iy = ky
246 IF (beta.EQ.zero) THEN
247 DO 30 i = 1,n
248 y(iy) = zero
249 iy = iy + incy
250 30 CONTINUE
251 ELSE
252 DO 40 i = 1,n
253 y(iy) = beta*y(iy)
254 iy = iy + incy
255 40 CONTINUE
256 END IF
257 END IF
258 END IF
259 IF (alpha.EQ.zero) RETURN
260 IF (lsame(uplo,'U')) THEN
261*
262* Form y when A is stored in upper triangle.
263*
264 IF ((incx.EQ.1) .AND. (incy.EQ.1)) THEN
265 DO 60 j = 1,n
266 temp1 = alpha*x(j)
267 temp2 = zero
268 DO 50 i = 1,j - 1
269 y(i) = y(i) + temp1*a(i,j)
270 temp2 = temp2 + dconjg(a(i,j))*x(i)
271 50 CONTINUE
272 y(j) = y(j) + temp1*dble(a(j,j)) + alpha*temp2
273 60 CONTINUE
274 ELSE
275 jx = kx
276 jy = ky
277 DO 80 j = 1,n
278 temp1 = alpha*x(jx)
279 temp2 = zero
280 ix = kx
281 iy = ky
282 DO 70 i = 1,j - 1
283 y(iy) = y(iy) + temp1*a(i,j)
284 temp2 = temp2 + dconjg(a(i,j))*x(ix)
285 ix = ix + incx
286 iy = iy + incy
287 70 CONTINUE
288 y(jy) = y(jy) + temp1*dble(a(j,j)) + alpha*temp2
289 jx = jx + incx
290 jy = jy + incy
291 80 CONTINUE
292 END IF
293 ELSE
294*
295* Form y when A is stored in lower triangle.
296*
297 IF ((incx.EQ.1) .AND. (incy.EQ.1)) THEN
298 DO 100 j = 1,n
299 temp1 = alpha*x(j)
300 temp2 = zero
301 y(j) = y(j) + temp1*dble(a(j,j))
302 DO 90 i = j + 1,n
303 y(i) = y(i) + temp1*a(i,j)
304 temp2 = temp2 + dconjg(a(i,j))*x(i)
305 90 CONTINUE
306 y(j) = y(j) + alpha*temp2
307 100 CONTINUE
308 ELSE
309 jx = kx
310 jy = ky
311 DO 120 j = 1,n
312 temp1 = alpha*x(jx)
313 temp2 = zero
314 y(jy) = y(jy) + temp1*dble(a(j,j))
315 ix = jx
316 iy = jy
317 DO 110 i = j + 1,n
318 ix = ix + incx
319 iy = iy + incy
320 y(iy) = y(iy) + temp1*a(i,j)
321 temp2 = temp2 + dconjg(a(i,j))*x(ix)
322 110 CONTINUE
323 y(jy) = y(jy) + alpha*temp2
324 jx = jx + incx
325 jy = jy + incy
326 120 CONTINUE
327 END IF
328 END IF
329*
330 RETURN
331*
332* End of ZHEMV
333*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: