csyr2k Subroutine

public subroutine csyr2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)

NAME

csyr2k(3f) - [BLAS:COMPLEX_BLAS_LEVEL3]

C:=alphaATRANSPOSE(B)+alphaBTRANSPOSE(A)+beta*C, C symmetric.

SYNOPSIS

 subroutine csyr2k(uplo,trans,n,k,alpha,a,lda,b,ldb,beta,c,ldc)

   .. Scalar Arguments ..
   complex,intent(in)    :: alpha,beta
   integer,intent(in)    :: k,lda,ldb,ldc,n
   character,intent(in)  :: trans,uplo
   ..
   .. Array Arguments ..
   complex,intent(in)    :: a(lda,*),b(ldb,*)
   complex,intent(inout) :: c(ldc,*)
   ..

DEFINITION

CSYR2K performs one of the symmetric rank 2k operations

 C := alpha*A*B**T + alpha*B*A**T + beta*C,

or

 C := alpha*A**T*B + alpha*B**T*A + beta*C,

where alpha and beta are scalars, C is an n by n symmetric matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.

OPTIONS

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*B**T + alpha*B*A**T +
                                      beta*C.

           TRANS = 'T' or 't'    C := alpha*A**T*B + alpha*B**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 matrices A and B, and on entry with
        TRANS = 'T' or 't', K specifies the number of rows of the
        matrices A and B. K must be at least zero.

ALPHA

       ALPHA is COMPLEX
        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 ).

B

       B is COMPLEX array, dimension ( LDB, kb ), where kb 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 B must contain the matrix B, otherwise
        the leading k by n part of the array B must contain the
        matrix B.

LDB

       LDB is INTEGER
        On entry, LDB specifies the first dimension of B as declared
        in the calling (sub) program. When TRANS = 'N' or 'n'
        then LDB must be at least max( 1, n ), otherwise LDB must
        be at least max( 1, k ).

BETA

       BETA is COMPLEX
        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 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 ).

AUTHORS

  • Univ. of Tennessee
  • Univ. of California Berkeley
  • Univ. of Colorado Denver
  • NAG Ltd.

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.

SEE ALSO

Online html documentation available at
http://www.netlib.org/lapack/explore-html/

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
integer, intent(in) :: n
integer, intent(in) :: k
complex, intent(in) :: alpha
complex, intent(in) :: a(lda,*)
integer, intent(in) :: lda
complex, intent(in) :: b(ldb,*)
integer, intent(in) :: ldb
complex, intent(in) :: beta
complex, intent(inout) :: c(ldc,*)
integer, intent(in) :: ldc

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer, public :: i
integer, public :: info
integer, public :: j
integer, public :: l
integer, public :: nrowa
complex, public, parameter :: one = (1.0e+0,0.0e+0)
complex, public :: temp1
complex, public :: temp2
logical, public :: upper
complex, public, parameter :: zero = (0.0e+0,0.0e+0)

Source Code

       subroutine csyr2k(uplo,trans,n,k,alpha,a,lda,b,ldb,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,intent(in)    :: alpha,beta
      integer,intent(in)    :: k,lda,ldb,ldc,n
      character,intent(in)  :: trans,uplo
!     ..
!     .. Array Arguments ..
      complex,intent(in)    :: a(lda,*),b(ldb,*)
      complex,intent(inout) :: c(ldc,*)
!     ..
!
!  =====================================================================
!
!     .. External Functions ..
!     ..
!     .. External Subroutines ..
!     ..
!     .. Intrinsic Functions ..
      intrinsic max
!     ..
!     .. Local Scalars ..
      complex temp1,temp2
      integer i,info,j,l,nrowa
      logical upper
!     ..
!     .. Parameters ..
      complex one
      parameter (one= (1.0e+0,0.0e+0))
      complex zero
      parameter (zero= (0.0e+0,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,'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 (ldb.lt.max(1,nrowa)) then
          info = 9
      elseif (ldc.lt.max(1,n)) then
          info = 12
      endif
      if (info.ne.0) then
          call xerbla('CSYR2K',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*B**T + alpha*B*A**T + 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) .or. (b(j,l).ne.zero)) then
                          temp1 = alpha*b(j,l)
                          temp2 = alpha*a(j,l)
                          do i = 1,j
                              c(i,j) = c(i,j) + a(i,l)*temp1 + b(i,l)*temp2
                          enddo
                      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) .or. (b(j,l).ne.zero)) then
                          temp1 = alpha*b(j,l)
                          temp2 = alpha*a(j,l)
                          do i = j,n
                              c(i,j) = c(i,j) + a(i,l)*temp1 + b(i,l)*temp2
                          enddo
                      endif
                  enddo
              enddo
          endif
      else
!
!        Form  C := alpha*A**T*B + alpha*B**T*A + C.
!
          if (upper) then
              do j = 1,n
                  do i = 1,j
                      temp1 = zero
                      temp2 = zero
                      do l = 1,k
                          temp1 = temp1 + a(l,i)*b(l,j)
                          temp2 = temp2 + b(l,i)*a(l,j)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp1 + alpha*temp2
                      else
                          c(i,j) = beta*c(i,j) + alpha*temp1 + alpha*temp2
                      endif
                  enddo
              enddo
          else
              do j = 1,n
                  do i = j,n
                      temp1 = zero
                      temp2 = zero
                      do l = 1,k
                          temp1 = temp1 + a(l,i)*b(l,j)
                          temp2 = temp2 + b(l,i)*a(l,j)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp1 + alpha*temp2
                      else
                          c(i,j) = beta*c(i,j) + alpha*temp1 + alpha*temp2
                      endif
                  enddo
              enddo
          endif
      endif
!
!     End of CSYR2K.
!
      end subroutine csyr2k