gonum.org/v1/gonum@v0.14.0/lapack/gonum/dgetrf.go (about) 1 // Copyright ©2015 The Gonum Authors. All rights reserved. 2 // Use of this source code is governed by a BSD-style 3 // license that can be found in the LICENSE file. 4 5 package gonum 6 7 import ( 8 "gonum.org/v1/gonum/blas" 9 "gonum.org/v1/gonum/blas/blas64" 10 ) 11 12 // Dgetrf computes the LU decomposition of the m×n matrix A. 13 // The LU decomposition is a factorization of A into 14 // 15 // A = P * L * U 16 // 17 // where P is a permutation matrix, L is a unit lower triangular matrix, and 18 // U is a (usually) non-unit upper triangular matrix. On exit, L and U are stored 19 // in place into a. 20 // 21 // ipiv is a permutation vector. It indicates that row i of the matrix was 22 // changed with ipiv[i]. ipiv must have length at least min(m,n), and will panic 23 // otherwise. ipiv is zero-indexed. 24 // 25 // Dgetrf is the blocked version of the algorithm. 26 // 27 // Dgetrf returns whether the matrix A is singular. The LU decomposition will 28 // be computed regardless of the singularity of A, but division by zero 29 // will occur if the false is returned and the result is used to solve a 30 // system of equations. 31 func (impl Implementation) Dgetrf(m, n int, a []float64, lda int, ipiv []int) (ok bool) { 32 mn := min(m, n) 33 switch { 34 case m < 0: 35 panic(mLT0) 36 case n < 0: 37 panic(nLT0) 38 case lda < max(1, n): 39 panic(badLdA) 40 } 41 42 // Quick return if possible. 43 if mn == 0 { 44 return true 45 } 46 47 switch { 48 case len(a) < (m-1)*lda+n: 49 panic(shortA) 50 case len(ipiv) != mn: 51 panic(badLenIpiv) 52 } 53 54 bi := blas64.Implementation() 55 56 nb := impl.Ilaenv(1, "DGETRF", " ", m, n, -1, -1) 57 if nb <= 1 || mn <= nb { 58 // Use the unblocked algorithm. 59 return impl.Dgetf2(m, n, a, lda, ipiv) 60 } 61 ok = true 62 for j := 0; j < mn; j += nb { 63 jb := min(mn-j, nb) 64 blockOk := impl.Dgetf2(m-j, jb, a[j*lda+j:], lda, ipiv[j:j+jb]) 65 if !blockOk { 66 ok = false 67 } 68 for i := j; i <= min(m-1, j+jb-1); i++ { 69 ipiv[i] = j + ipiv[i] 70 } 71 impl.Dlaswp(j, a, lda, j, j+jb-1, ipiv[:j+jb], 1) 72 if j+jb < n { 73 impl.Dlaswp(n-j-jb, a[j+jb:], lda, j, j+jb-1, ipiv[:j+jb], 1) 74 bi.Dtrsm(blas.Left, blas.Lower, blas.NoTrans, blas.Unit, 75 jb, n-j-jb, 1, 76 a[j*lda+j:], lda, 77 a[j*lda+j+jb:], lda) 78 if j+jb < m { 79 bi.Dgemm(blas.NoTrans, blas.NoTrans, m-j-jb, n-j-jb, jb, -1, 80 a[(j+jb)*lda+j:], lda, 81 a[j*lda+j+jb:], lda, 82 1, a[(j+jb)*lda+j+jb:], lda) 83 } 84 } 85 } 86 return ok 87 }