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