Skip to content

Commit ce4adcf

Browse files
committed
lapack/lapack64: add Pstrf
1 parent 90a1d53 commit ce4adcf

File tree

2 files changed

+34
-0
lines changed

2 files changed

+34
-0
lines changed

lapack/lapack.go

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -34,6 +34,7 @@ type Float64 interface {
3434
Dpotrf(ul blas.Uplo, n int, a []float64, lda int) (ok bool)
3535
Dpotri(ul blas.Uplo, n int, a []float64, lda int) (ok bool)
3636
Dpotrs(ul blas.Uplo, n, nrhs int, a []float64, lda int, b []float64, ldb int)
37+
Dpstrf(uplo blas.Uplo, n int, a []float64, lda int, piv []int, tol float64, work []float64) (rank int, ok bool)
3738
Dsyev(jobz EVJob, uplo blas.Uplo, n int, a []float64, lda int, w, work []float64, lwork int) (ok bool)
3839
Dtbtrs(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, kd, nrhs int, a []float64, lda int, b []float64, ldb int) (ok bool)
3940
Dtrcon(norm MatrixNorm, uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int, work []float64, iwork []int) float64

lapack/lapack64/lapack64.go

Lines changed: 33 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -128,6 +128,39 @@ func Pbtrs(t blas64.TriangularBand, b blas64.General) {
128128
lapack64.Dpbtrs(t.Uplo, t.N, t.K, b.Cols, t.Data, max(1, t.Stride), b.Data, max(1, b.Stride))
129129
}
130130

131+
// Pstrf computes the Cholesky factorization with complete pivoting of an n×n
132+
// symmetric positive semidefinite matrix A.
133+
//
134+
// The factorization has the form
135+
// Pᵀ * A * P = Uᵀ * U , if a.Uplo = blas.Upper,
136+
// Pᵀ * A * P = L * Lᵀ, if a.Uplo = blas.Lower,
137+
// where U is an upper triangular matrix, L is lower triangular, and P is a
138+
// permutation matrix.
139+
//
140+
// tol is a user-defined tolerance. The algorithm terminates if the pivot is
141+
// less than or equal to tol. If tol is negative, then n*eps*max(A[k,k]) will be
142+
// used instead.
143+
//
144+
// The triangular factor U or L from the Cholesky factorization is returned in t
145+
// and the underlying data between a and t is shared. P is stored on return in
146+
// vector piv such that P[piv[k],k] = 1.
147+
//
148+
// Pstrf also returns the computed rank of A and whether the algorithm completed
149+
// successfully. If ok is false, the matrix A is either rank deficient or is not
150+
// positive semidefinite.
151+
//
152+
// The length of piv must be n and the length of work must be at least 2*n,
153+
// otherwise Pstrf will panic.
154+
func Pstrf(a blas64.Symmetric, piv []int, tol float64, work []float64) (t blas64.Triangular, rank int, ok bool) {
155+
rank, ok = lapack64.Dpstrf(a.Uplo, a.N, a.Data, max(1, a.Stride), piv, tol, work)
156+
t.Uplo = a.Uplo
157+
t.Diag = blas.NonUnit
158+
t.N = a.N
159+
t.Data = a.Data
160+
t.Stride = a.Stride
161+
return t, rank, ok
162+
}
163+
131164
// Gecon estimates the reciprocal of the condition number of the n×n matrix A
132165
// given the LU decomposition of the matrix. The condition number computed may
133166
// be based on the 1-norm or the ∞-norm.

0 commit comments

Comments
 (0)