zsyrk(3f) - [BLAS:COMPLEX16_BLAS_LEVEL3]
subroutine zsyrk(uplo,trans,n,k,alpha,a,lda,beta,c,ldc)
.. Scalar Arguments ..
complex(kind=real64),intent(inout) :: alpha,beta
integer,intent(in) :: k,lda,ldc,n
character,intent(in) :: trans,uplo
..
.. Array Arguments ..
complex(kind=real64),intent(in) :: a(lda,*)
complex(kind=real64),intent(inout) :: c(ldc,*)
..
ZSYRK performs one of the symmetric rank k operations
C := alpha*A*A**T + beta*C,
or
C := alpha*A**T*A + beta*C,
where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case and a k by n matrix in the second case.
UPLO
UPLO is CHARACTER*1
On entry, UPLO specifies whether the upper or lower
triangular part of the array C is to be referenced as
follows:
UPLO = 'U' or 'u' Only the upper triangular part of C
is to be referenced.
UPLO = 'L' or 'l' Only the lower triangular part of C
is to be referenced.
TRANS
TRANS is CHARACTER*1
On entry, TRANS specifies the operation to be performed as
follows:
TRANS = 'N' or 'n' C := alpha*A*A**T + beta*C.
TRANS = 'T' or 't' C := alpha*A**T*A + beta*C.
N
N is INTEGER
On entry, N specifies the order of the matrix C. N must be
at least zero.
K
K is INTEGER
On entry with TRANS = 'N' or 'n', K specifies the number
of columns of the matrix A, and on entry with
TRANS = 'T' or 't', K specifies the number of rows of the
matrix A. K must be at least zero.
ALPHA
ALPHA is complex(kind=real64)
On entry, ALPHA specifies the scalar alpha.
A
A is complex(kind=real64) array, dimension ( LDA, ka ), where ka is
k when TRANS = 'N' or 'n', and is n otherwise.
Before entry with TRANS = 'N' or 'n', the leading n by k
part of the array A must contain the matrix A, otherwise
the leading k by n part of the array A must contain the
matrix A.
LDA
LDA is INTEGER
On entry, LDA specifies the first dimension of A as declared
in the calling (sub) program. When TRANS = 'N' or 'n'
then LDA must be at least max( 1, n ), otherwise LDA must
be at least max( 1, k ).
BETA
BETA is complex(kind=real64)
On entry, BETA specifies the scalar beta.
C
C is complex(kind=real64) array, dimension ( LDC, N )
Before entry with UPLO = 'U' or 'u', the leading n by n
upper triangular part of the array C must contain the upper
triangular part of the symmetric matrix and the strictly
lower triangular part of C is not referenced. On exit, the
upper triangular part of the array C is overwritten by the
upper triangular part of the updated matrix.
Before entry with UPLO = 'L' or 'l', the leading n by n
lower triangular part of the array C must contain the lower
triangular part of the symmetric matrix and the strictly
upper triangular part of C is not referenced. On exit, the
lower triangular part of the array C is overwritten by the
lower triangular part of the updated matrix.
LDC
LDC is INTEGER
On entry, LDC specifies the first dimension of C as declared
in the calling (sub) program. LDC must be at least
max( 1, n ).
date:December 2016
FURTHER DETAILS
Level 3 Blas routine.
– Written on 8-February-1989. Jack Dongarra, Argonne National Laboratory. Iain Duff, AERE Harwell. Jeremy Du Croz, Numerical Algorithms Group Ltd. Sven Hammarling, Numerical Algorithms Group Ltd.
Online html documentation available at
http://www.netlib.org/lapack/explore-html/
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=1), | intent(in) | :: | uplo | |||
character(len=1), | intent(in) | :: | trans | |||
integer, | intent(in) | :: | n | |||
integer, | intent(in) | :: | k | |||
complex(kind=real64), | intent(inout) | :: | alpha | |||
complex(kind=real64), | intent(in) | :: | a(lda,*) | |||
integer, | intent(in) | :: | lda | |||
complex(kind=real64), | intent(inout) | :: | beta | |||
complex(kind=real64), | intent(inout) | :: | c(ldc,*) | |||
integer, | intent(in) | :: | ldc |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | i | ||||
integer, | public | :: | info | ||||
integer, | public | :: | j | ||||
integer, | public | :: | l | ||||
integer, | public | :: | nrowa | ||||
complex(kind=real64), | public, | parameter | :: | one | = | (1.0d+0,0.0d+0) | |
complex(kind=real64), | public | :: | temp | ||||
logical, | public | :: | upper | ||||
complex(kind=real64), | public, | parameter | :: | zero | = | (0.0d+0,0.0d+0) |
subroutine zsyrk(uplo,trans,n,k,alpha,a,lda,beta,c,ldc)
implicit none
!
! -- Reference BLAS level3 routine (version 3.7.0) --
! -- Reference BLAS is a software package provided by Univ. of Tennessee, --
! -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
! December 2016
!
! .. Scalar Arguments ..
complex(kind=real64),intent(inout) :: alpha,beta
integer,intent(in) :: k,lda,ldc,n
character,intent(in) :: trans,uplo
! ..
! .. Array Arguments ..
complex(kind=real64),intent(in) :: a(lda,*)
complex(kind=real64),intent(inout) :: c(ldc,*)
! ..
!
! =====================================================================
!
! .. External Functions ..
! LOGICAL LSAME
! EXTERNAL LSAME
! ..
! .. External Subroutines ..
! EXTERNAL XERBLA
! ..
! .. Intrinsic Functions ..
intrinsic max
! ..
! .. Local Scalars ..
complex(kind=real64) :: temp
integer i,info,j,l,nrowa
logical upper
! ..
! .. Parameters ..
complex(kind=real64) :: one
parameter (one= (1.0d+0,0.0d+0))
complex(kind=real64) :: zero
parameter (zero= (0.0d+0,0.0d+0))
! ..
!
! Test the input parameters.
!
if (lsame(trans,'N')) then
nrowa = n
else
nrowa = k
endif
upper = lsame(uplo,'U')
!
info = 0
if ((.not.upper) .and. (.not.lsame(uplo,'L'))) then
info = 1
elseif ((.not.lsame(trans,'N')) .and. (.not.lsame(trans,'T'))) then
info = 2
elseif (n.lt.0) then
info = 3
elseif (k.lt.0) then
info = 4
elseif (lda.lt.max(1,nrowa)) then
info = 7
elseif (ldc.lt.max(1,n)) then
info = 10
endif
if (info.ne.0) then
call xerbla('ZSYRK ',info)
return
endif
!
! Quick return if possible.
!
if ((n.eq.0) .or. (((alpha.eq.zero).or. (k.eq.0)).and. (beta.eq.one))) return
!
! And when alpha.eq.zero.
!
if (alpha.eq.zero) then
if (upper) then
if (beta.eq.zero) then
do j = 1,n
c(1:j,j) = zero
enddo
else
do j = 1,n
c(1:j,j) = beta*c(1:j,j)
enddo
endif
else
if (beta.eq.zero) then
do j = 1,n
c(j:n,j) = zero
enddo
else
do j = 1,n
c(j:n,j) = beta*c(j:n,j)
enddo
endif
endif
return
endif
!
! Start the operations.
!
if (lsame(trans,'N')) then
!
! Form C := alpha*A*A**T + beta*C.
!
if (upper) then
do j = 1,n
if (beta.eq.zero) then
c(1:j,j) = zero
elseif (beta.ne.one) then
c(1:j,j) = beta*c(1:j,j)
endif
do l = 1,k
if (a(j,l).ne.zero) then
temp = alpha*a(j,l)
c(1:j,j) = c(1:j,j) + temp*a(1:j,l)
endif
enddo
enddo
else
do j = 1,n
if (beta.eq.zero) then
c(j:n,j) = zero
elseif (beta.ne.one) then
c(j:n,j) = beta*c(j:n,j)
endif
do l = 1,k
if (a(j,l).ne.zero) then
temp = alpha*a(j,l)
c(j:n,j) = c(j:n,j) + temp*a(j:n,l)
endif
enddo
enddo
endif
else
!
! Form C := alpha*A**T*A + beta*C.
!
if (upper) then
do j = 1,n
do i = 1,j
temp = zero
do l = 1,k
temp = temp + a(l,i)*a(l,j)
enddo
if (beta.eq.zero) then
c(i,j) = alpha*temp
else
c(i,j) = alpha*temp + beta*c(i,j)
endif
enddo
enddo
else
do j = 1,n
do i = j,n
temp = zero
do l = 1,k
temp = temp + a(l,i)*a(l,j)
enddo
if (beta.eq.zero) then
c(i,j) = alpha*temp
else
c(i,j) = alpha*temp + beta*c(i,j)
endif
enddo
enddo
endif
endif
end subroutine zsyrk