gonum.org/v1/gonum@v0.15.1-0.20240517103525-f853624cb1bb/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 an m×n matrix A using partial
    13  // pivoting with row interchanges.
    14  //
    15  // The LU decomposition is a factorization of A into
    16  //
    17  //	A = P * L * U
    18  //
    19  // where P is a permutation matrix, L is a lower triangular with unit diagonal
    20  // elements (lower trapezoidal if m > n), and U is upper triangular (upper
    21  // trapezoidal if m < n).
    22  //
    23  // On entry, a contains the matrix A. On return, L and U are stored in place
    24  // into a, and P is represented by ipiv.
    25  //
    26  // ipiv contains a sequence of row interchanges. It indicates that row i of the
    27  // matrix was interchanged with ipiv[i]. ipiv must have length min(m,n), and
    28  // Dgetrf will panic otherwise. ipiv is zero-indexed.
    29  //
    30  // Dgetrf returns whether the matrix A is nonsingular. The LU decomposition will
    31  // be computed regardless of the singularity of A, but the result should not be
    32  // used to solve a system of equation.
    33  func (impl Implementation) Dgetrf(m, n int, a []float64, lda int, ipiv []int) (ok bool) {
    34  	mn := min(m, n)
    35  	switch {
    36  	case m < 0:
    37  		panic(mLT0)
    38  	case n < 0:
    39  		panic(nLT0)
    40  	case lda < max(1, n):
    41  		panic(badLdA)
    42  	}
    43  
    44  	// Quick return if possible.
    45  	if mn == 0 {
    46  		return true
    47  	}
    48  
    49  	switch {
    50  	case len(a) < (m-1)*lda+n:
    51  		panic(shortA)
    52  	case len(ipiv) != mn:
    53  		panic(badLenIpiv)
    54  	}
    55  
    56  	bi := blas64.Implementation()
    57  
    58  	nb := impl.Ilaenv(1, "DGETRF", " ", m, n, -1, -1)
    59  	if nb <= 1 || mn <= nb {
    60  		// Use the unblocked algorithm.
    61  		return impl.Dgetf2(m, n, a, lda, ipiv)
    62  	}
    63  	ok = true
    64  	for j := 0; j < mn; j += nb {
    65  		jb := min(mn-j, nb)
    66  		blockOk := impl.Dgetf2(m-j, jb, a[j*lda+j:], lda, ipiv[j:j+jb])
    67  		if !blockOk {
    68  			ok = false
    69  		}
    70  		for i := j; i <= min(m-1, j+jb-1); i++ {
    71  			ipiv[i] = j + ipiv[i]
    72  		}
    73  		impl.Dlaswp(j, a, lda, j, j+jb-1, ipiv[:j+jb], 1)
    74  		if j+jb < n {
    75  			impl.Dlaswp(n-j-jb, a[j+jb:], lda, j, j+jb-1, ipiv[:j+jb], 1)
    76  			bi.Dtrsm(blas.Left, blas.Lower, blas.NoTrans, blas.Unit,
    77  				jb, n-j-jb, 1,
    78  				a[j*lda+j:], lda,
    79  				a[j*lda+j+jb:], lda)
    80  			if j+jb < m {
    81  				bi.Dgemm(blas.NoTrans, blas.NoTrans, m-j-jb, n-j-jb, jb, -1,
    82  					a[(j+jb)*lda+j:], lda,
    83  					a[j*lda+j+jb:], lda,
    84  					1, a[(j+jb)*lda+j+jb:], lda)
    85  			}
    86  		}
    87  	}
    88  	return ok
    89  }