table of contents
hfrk(3) | Library Functions Manual | hfrk(3) |
NAME¶
hfrk - hfrk: Hermitian rank-k update, RFP format
SYNOPSIS¶
Functions¶
subroutine CHFRK (transr, uplo, trans, n, k, alpha, a, lda,
beta, c)
CHFRK performs a Hermitian rank-k operation for matrix in RFP format.
subroutine DSFRK (transr, uplo, trans, n, k, alpha, a, lda, beta, c)
DSFRK performs a symmetric rank-k operation for matrix in RFP format.
subroutine SSFRK (transr, uplo, trans, n, k, alpha, a, lda, beta, c)
SSFRK performs a symmetric rank-k operation for matrix in RFP format.
subroutine ZHFRK (transr, uplo, trans, n, k, alpha, a, lda, beta, c)
ZHFRK performs a Hermitian rank-k operation for matrix in RFP format.
Detailed Description¶
Function Documentation¶
subroutine CHFRK (character transr, character uplo, character trans, integer n, integer k, real alpha, complex, dimension( lda, * ) a, integer lda, real beta, complex, dimension( * ) c)¶
CHFRK performs a Hermitian rank-k operation for matrix in RFP format.
Purpose:
!> !> Level 3 BLAS like routine for C in RFP Format. !> !> CHFRK 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. !>
Parameters
!> TRANSR is CHARACTER*1 !> = 'N': The Normal Form of RFP A is stored; !> = 'C': The Conjugate-transpose Form of RFP A is stored. !>
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. !> !> Unchanged on exit. !>
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. !> !> Unchanged on exit. !>
N
!> N is INTEGER !> On entry, N specifies the order of the matrix C. N must be !> at least zero. !> Unchanged on exit. !>
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. !> Unchanged on exit. !>
ALPHA
!> ALPHA is REAL !> On entry, ALPHA specifies the scalar alpha. !> Unchanged on exit. !>
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. !> Unchanged on exit. !>
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 ). !> Unchanged on exit. !>
BETA
!> BETA is REAL !> On entry, BETA specifies the scalar beta. !> Unchanged on exit. !>
C
!> C is COMPLEX array, dimension (N*(N+1)/2) !> On entry, the matrix A in RFP Format. RFP Format is !> described by TRANSR, UPLO and N. 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. !>
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Definition at line 166 of file chfrk.f.
subroutine DSFRK (character transr, character uplo, character trans, integer n, integer k, double precision alpha, double precision, dimension( lda, * ) a, integer lda, double precision beta, double precision, dimension( * ) c)¶
DSFRK performs a symmetric rank-k operation for matrix in RFP format.
Purpose:
!> !> Level 3 BLAS like routine for C in RFP Format. !> !> DSFRK 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 real 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. !>
Parameters
!> TRANSR is CHARACTER*1 !> = 'N': The Normal Form of RFP A is stored; !> = 'T': The Transpose Form of RFP A is stored. !>
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. !> !> Unchanged on exit. !>
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. !> !> Unchanged on exit. !>
N
!> N is INTEGER !> On entry, N specifies the order of the matrix C. N must be !> at least zero. !> Unchanged on exit. !>
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. !> Unchanged on exit. !>
ALPHA
!> ALPHA is DOUBLE PRECISION !> On entry, ALPHA specifies the scalar alpha. !> Unchanged on exit. !>
A
!> A is DOUBLE PRECISION 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. !> Unchanged on exit. !>
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 ). !> Unchanged on exit. !>
BETA
!> BETA is DOUBLE PRECISION !> On entry, BETA specifies the scalar beta. !> Unchanged on exit. !>
C
!> C is DOUBLE PRECISION array, dimension (NT) !> NT = N*(N+1)/2. On entry, the symmetric matrix C in RFP !> Format. RFP Format is described by TRANSR, UPLO and N. !>
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Definition at line 164 of file dsfrk.f.
subroutine SSFRK (character transr, character uplo, character trans, integer n, integer k, real alpha, real, dimension( lda, * ) a, integer lda, real beta, real, dimension( * ) c)¶
SSFRK performs a symmetric rank-k operation for matrix in RFP format.
Purpose:
!> !> Level 3 BLAS like routine for C in RFP Format. !> !> SSFRK 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 real 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. !>
Parameters
!> TRANSR is CHARACTER*1 !> = 'N': The Normal Form of RFP A is stored; !> = 'T': The Transpose Form of RFP A is stored. !>
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. !> !> Unchanged on exit. !>
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. !> !> Unchanged on exit. !>
N
!> N is INTEGER !> On entry, N specifies the order of the matrix C. N must be !> at least zero. !> Unchanged on exit. !>
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. !> Unchanged on exit. !>
ALPHA
!> ALPHA is REAL !> On entry, ALPHA specifies the scalar alpha. !> Unchanged on exit. !>
A
!> A is REAL 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. !> Unchanged on exit. !>
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 ). !> Unchanged on exit. !>
BETA
!> BETA is REAL !> On entry, BETA specifies the scalar beta. !> Unchanged on exit. !>
C
!> C is REAL array, dimension (NT) !> NT = N*(N+1)/2. On entry, the symmetric matrix C in RFP !> Format. RFP Format is described by TRANSR, UPLO and N. !>
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Definition at line 164 of file ssfrk.f.
subroutine ZHFRK (character transr, character uplo, character trans, integer n, integer k, double precision alpha, complex*16, dimension( lda, * ) a, integer lda, double precision beta, complex*16, dimension( * ) c)¶
ZHFRK performs a Hermitian rank-k operation for matrix in RFP format.
Purpose:
!> !> Level 3 BLAS like routine for C in RFP Format. !> !> ZHFRK 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. !>
Parameters
!> TRANSR is CHARACTER*1 !> = 'N': The Normal Form of RFP A is stored; !> = 'C': The Conjugate-transpose Form of RFP A is stored. !>
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. !> !> Unchanged on exit. !>
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. !> !> Unchanged on exit. !>
N
!> N is INTEGER !> On entry, N specifies the order of the matrix C. N must be !> at least zero. !> Unchanged on exit. !>
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. !> Unchanged on exit. !>
ALPHA
!> ALPHA is DOUBLE PRECISION !> On entry, ALPHA specifies the scalar alpha. !> Unchanged on exit. !>
A
!> A is COMPLEX*16 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. !> Unchanged on exit. !>
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 ). !> Unchanged on exit. !>
BETA
!> BETA is DOUBLE PRECISION !> On entry, BETA specifies the scalar beta. !> Unchanged on exit. !>
C
!> C is COMPLEX*16 array, dimension (N*(N+1)/2) !> On entry, the matrix A in RFP Format. RFP Format is !> described by TRANSR, UPLO and N. 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. !>
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Definition at line 166 of file zhfrk.f.
Author¶
Generated automatically by Doxygen for LAPACK from the source code.
Version 3.12.0 | LAPACK |