gonum.org/v1/gonum@v0.14.0/lapack/gonum/dpotrs.go (about)

     1  // Copyright ©2018 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  // Dpotrs solves a system of n linear equations A*X = B where A is an n×n
    13  // symmetric positive definite matrix and B is an n×nrhs matrix. The matrix A is
    14  // represented by its Cholesky factorization
    15  //
    16  //	A = Uᵀ*U  if uplo == blas.Upper
    17  //	A = L*Lᵀ  if uplo == blas.Lower
    18  //
    19  // as computed by Dpotrf. On entry, B contains the right-hand side matrix B, on
    20  // return it contains the solution matrix X.
    21  func (Implementation) Dpotrs(uplo blas.Uplo, n, nrhs int, a []float64, lda int, b []float64, ldb int) {
    22  	switch {
    23  	case uplo != blas.Upper && uplo != blas.Lower:
    24  		panic(badUplo)
    25  	case n < 0:
    26  		panic(nLT0)
    27  	case nrhs < 0:
    28  		panic(nrhsLT0)
    29  	case lda < max(1, n):
    30  		panic(badLdA)
    31  	case ldb < max(1, nrhs):
    32  		panic(badLdB)
    33  	}
    34  
    35  	// Quick return if possible.
    36  	if n == 0 || nrhs == 0 {
    37  		return
    38  	}
    39  
    40  	switch {
    41  	case len(a) < (n-1)*lda+n:
    42  		panic(shortA)
    43  	case len(b) < (n-1)*ldb+nrhs:
    44  		panic(shortB)
    45  	}
    46  
    47  	bi := blas64.Implementation()
    48  
    49  	if uplo == blas.Upper {
    50  		// Solve Uᵀ * U * X = B where U is stored in the upper triangle of A.
    51  
    52  		// Solve Uᵀ * X = B, overwriting B with X.
    53  		bi.Dtrsm(blas.Left, blas.Upper, blas.Trans, blas.NonUnit, n, nrhs, 1, a, lda, b, ldb)
    54  		// Solve U * X = B, overwriting B with X.
    55  		bi.Dtrsm(blas.Left, blas.Upper, blas.NoTrans, blas.NonUnit, n, nrhs, 1, a, lda, b, ldb)
    56  	} else {
    57  		// Solve L * Lᵀ * X = B where L is stored in the lower triangle of A.
    58  
    59  		// Solve L * X = B, overwriting B with X.
    60  		bi.Dtrsm(blas.Left, blas.Lower, blas.NoTrans, blas.NonUnit, n, nrhs, 1, a, lda, b, ldb)
    61  		// Solve Lᵀ * X = B, overwriting B with X.
    62  		bi.Dtrsm(blas.Left, blas.Lower, blas.Trans, blas.NonUnit, n, nrhs, 1, a, lda, b, ldb)
    63  	}
    64  }