cherk(3f) - [BLAS:COMPLEX_BLAS_LEVEL3] performs one of the hermitian rank k operations
C:=alphaATRANSPOSE(A)+beta*C, C hermitian.
subroutine cherk(uplo,trans,n,k,alpha,a,lda,beta,c,ldc)
.. Scalar Arguments ..
real,intent(in) :: alpha,beta
integer,intent(in) :: k,lda,ldc,n
character,intent(in) :: trans,uplo
..
.. Array Arguments ..
complex,intent(in) :: a(lda,*)
complex,intent(inout) :: c(ldc,*)
..
CHERK performs one of the hermitian rank k operations
C := alpha*A*A**H + beta*C,
or
C := alpha*A**H*A + beta*C,
where alpha and beta are real scalars, C is an n by n hermitian 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**H + beta*C.
TRANS = 'C' or 'c' C := alpha*A**H*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 = 'C' or 'c', K specifies the number of rows of the
matrix A. K must be at least zero.
ALPHA
ALPHA is REAL
On entry, ALPHA specifies the scalar alpha.
A
A is COMPLEX 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 REAL
On entry, BETA specifies the scalar beta.
C
C is COMPLEX 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 hermitian 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 hermitian 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.
Note that the imaginary parts of the diagonal elements need
not be set, they are assumed to be zero, and on exit they
are set to zero.
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.
– Modified 8-Nov-93 to set C(J,J) to REAL( C(J,J) ) when BETA = 1. Ed Anderson, Cray Research Inc.
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 | |||
real, | intent(in) | :: | alpha | |||
complex, | intent(in) | :: | a(lda,*) | |||
integer, | intent(in) | :: | lda | |||
real, | intent(in) | :: | beta | |||
complex, | 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 | ||||
real, | public, | parameter | :: | one | = | 1.0e+0 | |
real, | public | :: | rtemp | ||||
complex, | public | :: | temp | ||||
logical, | public | :: | upper | ||||
real, | public, | parameter | :: | zero | = | 0.0e+0 |
subroutine cherk(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 ..
real,intent(in) :: alpha,beta
integer,intent(in) :: k,lda,ldc,n
character,intent(in) :: trans,uplo
! ..
! .. Array Arguments ..
complex,intent(in) :: a(lda,*)
complex,intent(inout) :: c(ldc,*)
! ..
!
! =====================================================================
!
! .. External Functions ..
! ..
! .. External Subroutines ..
! ..
! .. Intrinsic Functions ..
intrinsic cmplx,conjg,max,real
! ..
! .. Local Scalars ..
complex temp
real rtemp
integer i,info,j,l,nrowa
logical upper
! ..
! .. Parameters ..
real one,zero
parameter (one=1.0e+0,zero=0.0e+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,'C'))) 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('CHERK ',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
do i = 1,j
c(i,j) = zero
enddo
enddo
else
do j = 1,n
do i = 1,j - 1
c(i,j) = beta*c(i,j)
enddo
c(j,j) = beta*real(c(j,j))
enddo
endif
else
if (beta.eq.zero) then
do j = 1,n
do i = j,n
c(i,j) = zero
enddo
enddo
else
do j = 1,n
c(j,j) = beta*real(c(j,j))
do i = j + 1,n
c(i,j) = beta*c(i,j)
enddo
enddo
endif
endif
return
endif
!
! Start the operations.
!
if (lsame(trans,'N')) then
!
! Form C := alpha*A*A**H + beta*C.
!
if (upper) then
do j = 1,n
if (beta.eq.zero) then
c(1:j,j) = zero
elseif (beta.ne.one) then
do i = 1,j - 1
c(i,j) = beta*c(i,j)
enddo
c(j,j) = beta*real(c(j,j))
else
c(j,j) = real(c(j,j))
endif
do l = 1,k
if (a(j,l).ne.cmplx(zero)) then
temp = alpha*conjg(a(j,l))
do i = 1,j - 1
c(i,j) = c(i,j) + temp*a(i,l)
enddo
c(j,j) = real(c(j,j)) + real(temp*a(i,l))
endif
enddo
enddo
else
do j = 1,n
if (beta.eq.zero) then
do i = j,n
c(i,j) = zero
enddo
elseif (beta.ne.one) then
c(j,j) = beta*real(c(j,j))
do i = j + 1,n
c(i,j) = beta*c(i,j)
enddo
else
c(j,j) = real(c(j,j))
endif
do l = 1,k
if (a(j,l).ne.cmplx(zero)) then
temp = alpha*conjg(a(j,l))
c(j,j) = real(c(j,j)) + real(temp*a(j,l))
do i = j + 1,n
c(i,j) = c(i,j) + temp*a(i,l)
enddo
endif
enddo
enddo
endif
else
!
! Form C := alpha*A**H*A + beta*C.
!
if (upper) then
do j = 1,n
do i = 1,j - 1
temp = zero
do l = 1,k
temp = temp + conjg(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
rtemp = zero
do l = 1,k
rtemp = rtemp + conjg(a(l,j))*a(l,j)
enddo
if (beta.eq.zero) then
c(j,j) = alpha*rtemp
else
c(j,j) = alpha*rtemp + beta*real(c(j,j))
endif
enddo
else
do j = 1,n
rtemp = zero
do l = 1,k
rtemp = rtemp + conjg(a(l,j))*a(l,j)
enddo
if (beta.eq.zero) then
c(j,j) = alpha*rtemp
else
c(j,j) = alpha*rtemp + beta*real(c(j,j))
endif
do i = j + 1,n
temp = zero
do l = 1,k
temp = temp + conjg(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 cherk