![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_chetrf (magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, magma_int_t *ipiv, magma_int_t *info) |
CHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_chetrf_gpu (magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info) |
CHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_dsidi (magma_uplo_t uplo, double *A, magma_int_t lda, magma_int_t n, magma_int_t *ipiv, double *det, magma_int_t *inert, double *work, magma_int_t job, magma_int_t *info) |
dsidi computes the determinant, inertia and inverse of a double precision symmetric matrix using the factors from dsytrf. | |
magma_int_t | magma_dsytrf (magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, magma_int_t *ipiv, magma_int_t *info) |
DSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_dsytrf_gpu (magma_uplo_t uplo, magma_int_t n, double *dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info) |
DSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_ssidi (magma_uplo_t uplo, float *A, magma_int_t lda, magma_int_t n, magma_int_t *ipiv, float *det, magma_int_t *inert, float *work, magma_int_t job, magma_int_t *info) |
ssidi computes the determinant, inertia and inverse of a float precision symmetric matrix using the factors from ssytrf. | |
magma_int_t | magma_ssytrf (magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, magma_int_t *ipiv, magma_int_t *info) |
SSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_ssytrf_gpu (magma_uplo_t uplo, magma_int_t n, float *dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info) |
SSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_zhetrf (magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, magma_int_t *ipiv, magma_int_t *info) |
ZHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magma_zhetrf_gpu (magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *dA, magma_int_t ldda, magma_int_t *ipiv, magma_int_t *info) |
ZHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method. | |
magma_int_t | magmablas_cdiinertia (magma_int_t n, magmaFloatComplex_const_ptr dA, magma_int_t ldda, int *dneig, magma_queue_t queue) |
magmablas_cdiinertia computes the inertia of a real diagonal matrix. | |
magma_int_t | magmablas_cheinertia (magma_uplo_t uplo, magma_int_t n, magmaFloatComplex_const_ptr dA, magma_int_t ldda, magma_int_t *ipiv, int *dneig, magma_queue_t queue) |
magmablas_cheinertia computes the inertia of a hermitian and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks. | |
magma_int_t | magmablas_ddiinertia (magma_int_t n, magmaDouble_const_ptr dA, magma_int_t ldda, int *dneig, magma_queue_t queue) |
magmablas_ddiinertia computes the inertia of a real diagonal matrix. | |
magma_int_t | magmablas_dsiinertia (magma_uplo_t uplo, magma_int_t n, magmaDouble_const_ptr dA, magma_int_t ldda, magma_int_t *ipiv, int *dneig, magma_queue_t queue) |
magmablas_dsiinertia computes the inertia of a symmetric and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks. | |
magma_int_t | magmablas_sdiinertia (magma_int_t n, magmaFloat_const_ptr dA, magma_int_t ldda, int *dneig, magma_queue_t queue) |
magmablas_sdiinertia computes the inertia of a real diagonal matrix. | |
magma_int_t | magmablas_ssiinertia (magma_uplo_t uplo, magma_int_t n, magmaFloat_const_ptr dA, magma_int_t ldda, magma_int_t *ipiv, int *dneig, magma_queue_t queue) |
magmablas_ssiinertia computes the inertia of a symmetric and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks. | |
magma_int_t | magmablas_zdiinertia (magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, int *dneig, magma_queue_t queue) |
magmablas_zdiinertia computes the inertia of a real diagonal matrix. | |
magma_int_t | magmablas_zheinertia (magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magma_int_t *ipiv, int *dneig, magma_queue_t queue) |
magmablas_zheinertia computes the inertia of a hermitian and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks. | |
magma_int_t magma_chetrf | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaFloatComplex * | A, | ||
magma_int_t | lda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
CHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | A | COMPLEX array, dimension (LDA,N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_chetrf_gpu | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaFloatComplex * | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
CHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | COMPLEX array, dimension (LDA,N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_dsidi | ( | magma_uplo_t | uplo, |
double * | A, | ||
magma_int_t | lda, | ||
magma_int_t | n, | ||
magma_int_t * | ipiv, | ||
double * | det, | ||
magma_int_t * | inert, | ||
double * | work, | ||
magma_int_t | job, | ||
magma_int_t * | info ) |
dsidi computes the determinant, inertia and inverse of a double precision symmetric matrix using the factors from dsytrf.
[in] | uplo | magma_uplo_t
|
[in,out] | A | DOUBLE PRECISION array, dimension (LDA,N) On entry, the output from dsytrf. On exit, the upper triangle of the inverse of the original matrix. The strict lower triangle is never referenced. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in] | ipiv | INTEGER array, dimension (N) The pivot vector from dsytrf. |
[out] | det | DOUBLE PRECISION array, dimension (2) Determinant of the original matrix: Determinant = det(0) * 10^det(1) with 1<=fabs(det(0))<10 or det(0) = 0. |
[out] | inert | INTEGER array, dimension (3) The inertia of the original matrix: inert(0) = number of positive eigenvalues. inert(1) = number of negative eigenvalues. inert(2) = number of zero eigenvalues. |
[in] | work | (workspace) DOUBLE PRECISION array, dimension (N) Work vector. Contents destroyed. |
[in] | job | INTEGER JOB has the decimal expansion abc where if c != 0, the inverse is computed, if b != 0, the determinant is computed, if a != 0, the inertia is computed. |
for example, job = 111 gives all three.
Error condition: a division by zero may occur if the inverse is requested and dsytrf has set info != 0.
This routine is based on the LINPACK dsidi routine (see also LINPACK USERS' guide, page 5.15).
magma_int_t magma_dsytrf | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
double * | A, | ||
magma_int_t | lda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
DSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | A | DOUBLE PRECISION array, dimension (LDA,N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_dsytrf_gpu | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
double * | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
DSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | DOUBLE PRECISION array, dimension (LDA,N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_ssidi | ( | magma_uplo_t | uplo, |
float * | A, | ||
magma_int_t | lda, | ||
magma_int_t | n, | ||
magma_int_t * | ipiv, | ||
float * | det, | ||
magma_int_t * | inert, | ||
float * | work, | ||
magma_int_t | job, | ||
magma_int_t * | info ) |
ssidi computes the determinant, inertia and inverse of a float precision symmetric matrix using the factors from ssytrf.
[in] | uplo | magma_uplo_t
|
[in,out] | A | REAL array, dimension (LDA,N) On entry, the output from ssytrf. On exit, the upper triangle of the inverse of the original matrix. The strict lower triangle is never referenced. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in] | ipiv | INTEGER array, dimension (N) The pivot vector from ssytrf. |
[out] | det | REAL array, dimension (2) Determinant of the original matrix: Determinant = det(0) * 10^det(1) with 1<=fabsf(det(0))<10 or det(0) = 0. |
[out] | inert | INTEGER array, dimension (3) The inertia of the original matrix: inert(0) = number of positive eigenvalues. inert(1) = number of negative eigenvalues. inert(2) = number of zero eigenvalues. |
[in] | work | (workspace) REAL array, dimension (N) Work vector. Contents destroyed. |
[in] | job | INTEGER JOB has the decimal expansion abc where if c != 0, the inverse is computed, if b != 0, the determinant is computed, if a != 0, the inertia is computed. |
for example, job = 111 gives all three.
Error condition: a division by zero may occur if the inverse is requested and ssytrf has set info != 0.
This routine is based on the LINPACK ssidi routine (see also LINPACK USERS' guide, page 5.15).
magma_int_t magma_ssytrf | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
float * | A, | ||
magma_int_t | lda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
SSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | A | REAL array, dimension (LDA,N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_ssytrf_gpu | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
float * | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
SSYTRF computes the factorization of a real symmetric matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | REAL array, dimension (LDA,N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_zhetrf | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaDoubleComplex * | A, | ||
magma_int_t | lda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
ZHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | A | COMPLEX*16 array, dimension (LDA,N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magma_zhetrf_gpu | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaDoubleComplex * | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
magma_int_t * | info ) |
ZHETRF computes the factorization of a complex Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
A = U*D*U^H or A = L*D*L^H
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This is the blocked version of the algorithm, calling Level 3 BLAS.
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrix A. N >= 0. |
[in,out] | dA | COMPLEX*16 array, dimension (LDA,N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, the block diagonal matrix D and the multipliers used to obtain the factor U or L (see below for further details). |
[in] | ldda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[out] | ipiv | INTEGER array, dimension (N) Details of the interchanges and the block structure of D. If IPIV(k) > 0, then rows and columns k and IPIV(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. If UPLO = MagmaUpper and IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. If UPLO = MagmaLower and IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block. |
[out] | info | INTEGER
|
If UPLO = MagmaUpper, then A = U*D*U', where U = P(n)*U(n)* ... P(k)U(k) ..., i.e., U is a product of terms P(k)*U(k), where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s
U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If UPLO = MagmaLower, then A = L*D*L', where L = P(1)*L(1)* ... P(k)*L(k) ..., i.e., L is a product of terms P(k)*L(k), where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by IPIV(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1
L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
magma_int_t magmablas_cdiinertia | ( | magma_int_t | n, |
magmaFloatComplex_const_ptr | dA, | ||
magma_int_t | ldda, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_cdiinertia computes the inertia of a real diagonal matrix.
If matrix entries are complex, magmablas_cdiinertia considers the real part of the diagonal.
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | COMPLEX array of DIMENSION ( LDDA, n ). The input matrix A with diagonal entries for which the inertia is computed. If dA is complex, the computation is done on the real part of the diagonal. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_cheinertia | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaFloatComplex_const_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_cheinertia computes the inertia of a hermitian and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks.
These are matrices comming from the Bunch-Kaufman with diagonal pivoting factorizations (the CHETRF routine).
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | COMPLEX array of DIMENSION ( LDDA, n ). The input matrix A with 1-by-1 and 2-by-2 diagonal block entries for which the inertia is computed. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[in] | ipiv | INTEGER array, dimension (N) The pivot vector from dsytrf. |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_ddiinertia | ( | magma_int_t | n, |
magmaDouble_const_ptr | dA, | ||
magma_int_t | ldda, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_ddiinertia computes the inertia of a real diagonal matrix.
If matrix entries are real, magmablas_ddiinertia considers the real part of the diagonal.
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | DOUBLE PRECISION array of DIMENSION ( LDDA, n ). The input matrix A with diagonal entries for which the inertia is computed. If dA is real, the computation is done on the real part of the diagonal. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_dsiinertia | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaDouble_const_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_dsiinertia computes the inertia of a symmetric and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks.
These are matrices comming from the Bunch-Kaufman with diagonal pivoting factorizations (the DSYTRF routine).
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | DOUBLE PRECISION array of DIMENSION ( LDDA, n ). The input matrix A with 1-by-1 and 2-by-2 diagonal block entries for which the inertia is computed. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[in] | ipiv | INTEGER array, dimension (N) The pivot vector from dsytrf. |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_sdiinertia | ( | magma_int_t | n, |
magmaFloat_const_ptr | dA, | ||
magma_int_t | ldda, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_sdiinertia computes the inertia of a real diagonal matrix.
If matrix entries are real, magmablas_sdiinertia considers the real part of the diagonal.
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | REAL array of DIMENSION ( LDDA, n ). The input matrix A with diagonal entries for which the inertia is computed. If dA is real, the computation is done on the real part of the diagonal. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_ssiinertia | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaFloat_const_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_ssiinertia computes the inertia of a symmetric and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks.
These are matrices comming from the Bunch-Kaufman with diagonal pivoting factorizations (the SSYTRF routine).
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | REAL array of DIMENSION ( LDDA, n ). The input matrix A with 1-by-1 and 2-by-2 diagonal block entries for which the inertia is computed. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[in] | ipiv | INTEGER array, dimension (N) The pivot vector from dsytrf. |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_zdiinertia | ( | magma_int_t | n, |
magmaDoubleComplex_const_ptr | dA, | ||
magma_int_t | ldda, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_zdiinertia computes the inertia of a real diagonal matrix.
If matrix entries are complex, magmablas_zdiinertia considers the real part of the diagonal.
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | COMPLEX_16 array of DIMENSION ( LDDA, n ). The input matrix A with diagonal entries for which the inertia is computed. If dA is complex, the computation is done on the real part of the diagonal. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |
magma_int_t magmablas_zheinertia | ( | magma_uplo_t | uplo, |
magma_int_t | n, | ||
magmaDoubleComplex_const_ptr | dA, | ||
magma_int_t | ldda, | ||
magma_int_t * | ipiv, | ||
int * | dneig, | ||
magma_queue_t | queue ) |
magmablas_zheinertia computes the inertia of a hermitian and block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks.
These are matrices comming from the Bunch-Kaufman with diagonal pivoting factorizations (the ZHETRF routine).
[in] | n | INTEGER. On entry, N specifies the order of the matrix A. N must be at least zero. |
[in] | dA | COMPLEX_16 array of DIMENSION ( LDDA, n ). The input matrix A with 1-by-1 and 2-by-2 diagonal block entries for which the inertia is computed. |
[in] | ldda | INTEGER. On entry, LDDA specifies the leading dimension of A. LDDA must be at least max( 1, n ). |
[in] | ipiv | INTEGER array, dimension (N) The pivot vector from dsytrf. |
[out] | dneig | INTEGER array of DIMENSION 3 on the GPU memory. The number of positive, negative, and zero eigenvalues in this order. |
[in] | queue | magma_queue_t. Queue to execute in. |