mirror of
https://github.com/gonum/gonum.git
synced 2025-10-05 23:26:52 +08:00
lapack/lapack64: add Potri to Float64 interface
This commit is contained in:

committed by
Vladimír Chalupecký

parent
ab2339bae3
commit
ce5163176b
@@ -29,6 +29,7 @@ type Float64 interface {
|
|||||||
Dormlq(side blas.Side, trans blas.Transpose, m, n, k int, a []float64, lda int, tau, c []float64, ldc int, work []float64, lwork int)
|
Dormlq(side blas.Side, trans blas.Transpose, m, n, k int, a []float64, lda int, tau, c []float64, ldc int, work []float64, lwork int)
|
||||||
Dpocon(uplo blas.Uplo, n int, a []float64, lda int, anorm float64, work []float64, iwork []int) float64
|
Dpocon(uplo blas.Uplo, n int, a []float64, lda int, anorm float64, work []float64, iwork []int) float64
|
||||||
Dpotrf(ul blas.Uplo, n int, a []float64, lda int) (ok bool)
|
Dpotrf(ul blas.Uplo, n int, a []float64, lda int) (ok bool)
|
||||||
|
Dpotri(ul blas.Uplo, n int, a []float64, lda int) (ok bool)
|
||||||
Dpotrs(ul blas.Uplo, n, nrhs int, a []float64, lda int, b []float64, ldb int)
|
Dpotrs(ul blas.Uplo, n, nrhs int, a []float64, lda int, b []float64, ldb int)
|
||||||
Dsyev(jobz EVJob, uplo blas.Uplo, n int, a []float64, lda int, w, work []float64, lwork int) (ok bool)
|
Dsyev(jobz EVJob, uplo blas.Uplo, n int, a []float64, lda int, w, work []float64, lwork int) (ok bool)
|
||||||
Dtrcon(norm MatrixNorm, uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int, work []float64, iwork []int) float64
|
Dtrcon(norm MatrixNorm, uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int, work []float64, iwork []int) float64
|
||||||
|
@@ -37,6 +37,26 @@ func Potrf(a blas64.Symmetric) (t blas64.Triangular, ok bool) {
|
|||||||
return
|
return
|
||||||
}
|
}
|
||||||
|
|
||||||
|
// Potri computes the inverse of a real symmetric positive definite matrix A
|
||||||
|
// using its Cholesky factorization.
|
||||||
|
//
|
||||||
|
// On entry, t contains the triangular factor U or L from the Cholesky
|
||||||
|
// factorization A = U^T*U or A = L*L^T, as computed by Potrf.
|
||||||
|
//
|
||||||
|
// On return, the upper or lower triangle of the (symmetric) inverse of A is
|
||||||
|
// stored in t, overwriting the input factor U or L, and also returned in a. The
|
||||||
|
// underlying data between a and t is shared.
|
||||||
|
//
|
||||||
|
// The returned bool indicates whether the inverse was computed successfully.
|
||||||
|
func Potri(t blas64.Triangular) (a blas64.Symmetric, ok bool) {
|
||||||
|
ok = lapack64.Dpotri(t.Uplo, t.N, t.Data, t.Stride)
|
||||||
|
a.Uplo = t.Uplo
|
||||||
|
a.N = t.N
|
||||||
|
a.Data = t.Data
|
||||||
|
a.Stride = t.Stride
|
||||||
|
return
|
||||||
|
}
|
||||||
|
|
||||||
// Potrs solves a system of n linear equations A*X = B where A is an n×n
|
// Potrs solves a system of n linear equations A*X = B where A is an n×n
|
||||||
// symmetric positive definite matrix and B is an n×nrhs matrix, using the
|
// symmetric positive definite matrix and B is an n×nrhs matrix, using the
|
||||||
// Cholesky factorization A = U^T*U or A = L*L^T. t contains the corresponding
|
// Cholesky factorization A = U^T*U or A = L*L^T. t contains the corresponding
|
||||||
|
Reference in New Issue
Block a user