Skip to content
This repository was archived by the owner on Nov 24, 2018. It is now read-only.

Commit 06e6694

Browse files
committed
Merge pull request #52 from gonum/adddtrtri
Adddtrtri
2 parents c5ce834 + 2416fc0 commit 06e6694

File tree

7 files changed

+178
-4
lines changed

7 files changed

+178
-4
lines changed

cgo/lapack.go

Lines changed: 17 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -576,6 +576,23 @@ func (impl Implementation) Dtrcon(norm lapack.MatrixNorm, uplo blas.Uplo, diag b
576576
return rcond[0]
577577
}
578578

579+
// Dtrtri computes the inverse of a triangular matrix, storing the result in place
580+
// into a. This is the BLAS level 3 version of the algorithm which builds upon
581+
// Dtrti2 to operate on matrix blocks instead of only individual columns.
582+
//
583+
// Dtrti returns whether the matrix a is singular or whether it's not singular.
584+
// If the matrix is singular the inversion is not performed.
585+
func (impl Implementation) Dtrtri(uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int) (ok bool) {
586+
checkMatrix(n, n, a, lda)
587+
if uplo != blas.Upper && uplo != blas.Lower {
588+
panic(badUplo)
589+
}
590+
if diag != blas.NonUnit && diag != blas.Unit {
591+
panic(badDiag)
592+
}
593+
return clapack.Dtrtri(uplo, diag, n, a, lda)
594+
}
595+
579596
// Dtrtrs solves a triangular system of the form A * X = B or A^T * X = B. Dtrtrs
580597
// returns whether the solve completed successfully. If A is singular, no solve is performed.
581598
func (impl Implementation) Dtrtrs(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, nrhs int, a []float64, lda int, b []float64, ldb int) (ok bool) {

cgo/lapack_test.go

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -94,3 +94,7 @@ func TestDpocon(t *testing.T) {
9494
func TestDtrcon(t *testing.T) {
9595
testlapack.DtrconTest(t, impl)
9696
}
97+
98+
func TestDtrtri(t *testing.T) {
99+
testlapack.DtrtriTest(t, impl)
100+
}

native/dtrti2.go

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -6,7 +6,7 @@ import (
66
)
77

88
// Dtrti2 computes the inverse of a triangular matrix, storing the result in place
9-
// into a.
9+
// into a. This is the BLAS level 2 version of the algorithm.
1010
func (impl Implementation) Dtrti2(uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int) {
1111
checkMatrix(n, n, a, lda)
1212
if uplo != blas.Upper && uplo != blas.Lower {

native/dtrtri.go

Lines changed: 60 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,60 @@
1+
package native
2+
3+
import (
4+
"github.com/gonum/blas"
5+
"github.com/gonum/blas/blas64"
6+
)
7+
8+
// Dtrtri computes the inverse of a triangular matrix, storing the result in place
9+
// into a. This is the BLAS level 3 version of the algorithm which builds upon
10+
// Dtrti2 to operate on matrix blocks instead of only individual columns.
11+
//
12+
// Dtrti returns whether the matrix a is singular or whether it's not singular.
13+
// If the matrix is singular the inversion is not performed.
14+
func (impl Implementation) Dtrtri(uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int) (ok bool) {
15+
checkMatrix(n, n, a, lda)
16+
if uplo != blas.Upper && uplo != blas.Lower {
17+
panic(badUplo)
18+
}
19+
if diag != blas.NonUnit && diag != blas.Unit {
20+
panic(badDiag)
21+
}
22+
if n == 0 {
23+
return
24+
}
25+
nonUnit := diag == blas.NonUnit
26+
if nonUnit {
27+
for i := 0; i < n; i++ {
28+
if a[i*lda+i] == 0 {
29+
return false
30+
}
31+
}
32+
}
33+
34+
bi := blas64.Implementation()
35+
36+
nb := impl.Ilaenv(1, "DTRTRI", "UD", n, -1, -1, -1)
37+
if nb <= 1 || nb > n {
38+
impl.Dtrti2(uplo, diag, n, a, lda)
39+
return true
40+
}
41+
if uplo == blas.Upper {
42+
for j := 0; j < n; j += nb {
43+
jb := min(nb, n-j)
44+
bi.Dtrmm(blas.Left, blas.Upper, blas.NoTrans, diag, j, jb, 1, a, lda, a[j:], lda)
45+
bi.Dtrsm(blas.Right, blas.Upper, blas.NoTrans, diag, j, jb, -1, a[j*lda+j:], lda, a[j:], lda)
46+
impl.Dtrti2(blas.Upper, diag, jb, a[j*lda+j:], lda)
47+
}
48+
return true
49+
}
50+
nn := ((n - 1) / nb) * nb
51+
for j := nn; j >= 0; j -= nb {
52+
jb := min(nb, n-j)
53+
if j+jb <= n-1 {
54+
bi.Dtrmm(blas.Left, blas.Lower, blas.NoTrans, diag, n-j-jb, jb, 1, a[(j+jb)*lda+j+jb:], lda, a[(j+jb)*lda+j:], lda)
55+
bi.Dtrmm(blas.Right, blas.Lower, blas.NoTrans, diag, n-j-jb, jb, -1, a[j*lda+j:], lda, a[(j+jb)*lda+j:], lda)
56+
}
57+
impl.Dtrti2(blas.Lower, diag, jb, a[j*lda+j:], lda)
58+
}
59+
return true
60+
}

native/lapack_test.go

Lines changed: 6 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -116,8 +116,12 @@ func TestDtrcon(t *testing.T) {
116116
testlapack.DtrconTest(t, impl)
117117
}
118118

119-
func TestDtrtri2(t *testing.T) {
120-
testlapack.Dtrtri2Test(t, impl)
119+
func TestDtrti2(t *testing.T) {
120+
testlapack.Dtrti2Test(t, impl)
121+
}
122+
123+
func TestDtrtri(t *testing.T) {
124+
testlapack.DtrtriTest(t, impl)
121125
}
122126

123127
func TestIladlc(t *testing.T) {

testlapack/dtrti2.go

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -14,7 +14,7 @@ type Dtrti2er interface {
1414
Dtrti2(uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int)
1515
}
1616

17-
func Dtrtri2Test(t *testing.T, impl Dtrti2er) {
17+
func Dtrti2Test(t *testing.T, impl Dtrti2er) {
1818
for _, test := range []struct {
1919
a []float64
2020
n int

testlapack/dtrtri.go

Lines changed: 89 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,89 @@
1+
package testlapack
2+
3+
import (
4+
"math"
5+
"math/rand"
6+
"testing"
7+
8+
"github.com/gonum/blas"
9+
"github.com/gonum/blas/blas64"
10+
)
11+
12+
type Dtrtrier interface {
13+
Dtrconer
14+
Dtrtri(uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int) bool
15+
}
16+
17+
func DtrtriTest(t *testing.T, impl Dtrtrier) {
18+
bi := blas64.Implementation()
19+
for _, uplo := range []blas.Uplo{blas.Upper} {
20+
for _, diag := range []blas.Diag{blas.NonUnit, blas.Unit} {
21+
for _, test := range []struct {
22+
n, lda int
23+
}{
24+
{3, 0},
25+
{70, 0},
26+
{200, 0},
27+
{3, 5},
28+
{70, 92},
29+
{200, 205},
30+
} {
31+
n := test.n
32+
lda := test.lda
33+
if lda == 0 {
34+
lda = n
35+
}
36+
a := make([]float64, n*lda)
37+
for i := range a {
38+
a[i] = rand.Float64() + 1 // This keeps the matrices well conditioned.
39+
}
40+
aCopy := make([]float64, len(a))
41+
copy(aCopy, a)
42+
impl.Dtrtri(uplo, diag, n, a, lda)
43+
if uplo == blas.Upper {
44+
for i := 1; i < n; i++ {
45+
for j := 0; j < i; j++ {
46+
aCopy[i*lda+j] = 0
47+
a[i*lda+j] = 0
48+
}
49+
}
50+
} else {
51+
for i := 1; i < n; i++ {
52+
for j := i + 1; j < n; j++ {
53+
aCopy[i*lda+j] = 0
54+
a[i*lda+j] = 0
55+
}
56+
}
57+
}
58+
if diag == blas.Unit {
59+
for i := 0; i < n; i++ {
60+
a[i*lda+i] = 1
61+
aCopy[i*lda+i] = 1
62+
}
63+
}
64+
ans := make([]float64, len(a))
65+
bi.Dgemm(blas.NoTrans, blas.NoTrans, n, n, n, 1, a, lda, aCopy, lda, 0, ans, lda)
66+
iseye := true
67+
for i := 0; i < n; i++ {
68+
for j := 0; j < n; j++ {
69+
if i == j {
70+
if math.Abs(ans[i*lda+i]-1) > 1e-4 {
71+
iseye = false
72+
break
73+
}
74+
} else {
75+
if math.Abs(ans[i*lda+j]) > 1e-4 {
76+
iseye = false
77+
break
78+
}
79+
}
80+
}
81+
}
82+
if !iseye {
83+
t.Errorf("inv(A) * A != I. Upper = %v, unit = %v, n = %v, lda = %v",
84+
uplo == blas.Upper, diag == blas.Unit, n, lda)
85+
}
86+
}
87+
}
88+
}
89+
}

0 commit comments

Comments
 (0)