This repository has been archived by the owner on Nov 24, 2018. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 11
/
lapack.go
280 lines (265 loc) · 9.37 KB
/
lapack.go
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
// Copyright ©2015 The gonum Authors. All rights reserved.
// Use of this source code is governed by a BSD-style
// license that can be found in the LICENSE file.
// Package cgo provides an interface to bindings for a C LAPACK library.
package cgo
import (
"github.com/gonum/blas"
"github.com/gonum/lapack"
"github.com/gonum/lapack/cgo/clapack"
)
// Copied from lapack/native. Keep in sync.
const (
absIncNotOne = "lapack: increment not one or negative one"
badDirect = "lapack: bad direct"
badIpiv = "lapack: insufficient permutation length"
badLdA = "lapack: index of a out of range"
badSide = "lapack: bad side"
badStore = "lapack: bad store"
badTau = "lapack: tau has insufficient length"
badTrans = "lapack: bad trans"
badUplo = "lapack: illegal triangle"
badWork = "lapack: insufficient working memory"
badWorkStride = "lapack: insufficient working array stride"
negDimension = "lapack: negative matrix dimension"
nLT0 = "lapack: n < 0"
shortWork = "lapack: working array shorter than declared"
)
func min(m, n int) int {
if m < n {
return m
}
return n
}
// checkMatrix verifies the parameters of a matrix input.
// Copied from lapack/native. Keep in sync.
func checkMatrix(m, n int, a []float64, lda int) {
if m < 0 {
panic("lapack: has negative number of rows")
}
if m < 0 {
panic("lapack: has negative number of columns")
}
if lda < n {
panic("lapack: stride less than number of columns")
}
if len(a) < (m-1)*lda+n {
panic("lapack: insufficient matrix slice length")
}
}
// Implementation is the cgo-based C implementation of LAPACK routines.
type Implementation struct{}
var _ lapack.Float64 = Implementation{}
// Dpotrf computes the cholesky decomposition of the symmetric positive definite
// matrix a. If ul == blas.Upper, then a is stored as an upper-triangular matrix,
// and a = U U^T is stored in place into a. If ul == blas.Lower, then a = L L^T
// is computed and stored in-place into a. If a is not positive definite, false
// is returned. This is the blocked version of the algorithm.
func (impl Implementation) Dpotrf(ul blas.Uplo, n int, a []float64, lda int) (ok bool) {
// ul is checked in clapack.Dpotrf.
if n < 0 {
panic(nLT0)
}
if lda < n {
panic(badLdA)
}
if n == 0 {
return true
}
return clapack.Dpotrf(ul, n, a, lda)
}
// Dgelq2 computes the LQ factorization of the m×n matrix A.
//
// In an LQ factorization, L is a lower triangular m×n matrix, and Q is an n×n
// orthornormal matrix.
//
// a is modified to contain the information to construct L and Q.
// The lower triangle of a contains the matrix L. The upper triangular elements
// (not including the diagonal) contain the elementary reflectors. Tau is modified
// to contain the reflector scales. tau must have length of at least k = min(m,n)
// and this function will panic otherwise.
//
// See Dgeqr2 for a description of the elementary reflectors and orthonormal
// matrix Q. Q is constructed as a product of these elementary reflectors,
// Q = H_k ... H_2*H_1.
//
// Work is temporary storage of length at least m and this function will panic otherwise.
func (impl Implementation) Dgelq2(m, n int, a []float64, lda int, tau, work []float64) {
checkMatrix(m, n, a, lda)
k := min(m, n)
if len(tau) < k {
panic(badTau)
}
if len(work) < m {
panic(badWork)
}
clapack.Dgelq2(m, n, a, lda, tau)
}
// Dgelqf computes the LQ factorization of the m×n matrix A using a blocked
// algorithm. See the documentation for Dgelq2 for a description of the
// parameters at entry and exit.
//
// The C interface does not support providing temporary storage. To provide compatibility
// with native, lwork == -1 will not run Dgeqrf but will instead write the minimum
// work necessary to work[0]. If len(work) < lwork, Dgeqrf will panic.
//
// tau must have length at least min(m,n), and this function will panic otherwise.
func (impl Implementation) Dgelqf(m, n int, a []float64, lda int, tau, work []float64, lwork int) {
if lwork == -1 {
work[0] = float64(m)
return
}
checkMatrix(m, n, a, lda)
if len(work) < lwork {
panic(shortWork)
}
if lwork < m {
panic(badWork)
}
k := min(m, n)
if len(tau) < k {
panic(badTau)
}
clapack.Dgelqf(m, n, a, lda, tau)
}
// Dgeqr2 computes a QR factorization of the m×n matrix A.
//
// In a QR factorization, Q is an m×m orthonormal matrix, and R is an
// upper triangular m×n matrix.
//
// A is modified to contain the information to construct Q and R.
// The upper triangle of a contains the matrix R. The lower triangular elements
// (not including the diagonal) contain the elementary reflectors. Tau is modified
// to contain the reflector scales. tau must have length at least min(m,n), and
// this function will panic otherwise.
//
// The ith elementary reflector can be explicitly constructed by first extracting
// the
// v[j] = 0 j < i
// v[j] = i j == i
// v[j] = a[i*lda+j] j > i
// and computing h_i = I - tau[i] * v * v^T.
//
// The orthonormal matrix Q can be constucted from a product of these elementary
// reflectors, Q = H_1*H_2 ... H_k, where k = min(m,n).
//
// Work is temporary storage of length at least n and this function will panic otherwise.
func (impl Implementation) Dgeqr2(m, n int, a []float64, lda int, tau, work []float64) {
checkMatrix(m, n, a, lda)
if len(work) < n {
panic(badWork)
}
k := min(m, n)
if len(tau) < k {
panic(badTau)
}
clapack.Dgeqr2(m, n, a, lda, tau)
}
// Dgeqrf computes the QR factorization of the m×n matrix A using a blocked
// algorithm. See the documentation for Dgeqr2 for a description of the
// parameters at entry and exit.
//
// The C interface does not support providing temporary storage. To provide compatibility
// with native, lwork == -1 will not run Dgeqrf but will instead write the minimum
// work necessary to work[0]. If len(work) < lwork, Dgeqrf will panic.
//
// tau must have length at least min(m,n), and this function will panic otherwise.
func (impl Implementation) Dgeqrf(m, n int, a []float64, lda int, tau, work []float64, lwork int) {
if lwork == -1 {
work[0] = float64(n)
return
}
checkMatrix(m, n, a, lda)
if len(work) < lwork {
panic(shortWork)
}
if lwork < n {
panic(badWork)
}
k := min(m, n)
if len(tau) < k {
panic(badTau)
}
clapack.Dgeqrf(m, n, a, lda, tau)
}
// Dgetf2 computes the LU decomposition of the m×n matrix A.
// The LU decomposition is a factorization of a into
// A = P * L * U
// where P is a permutation matrix, L is a unit lower triangular matrix, and
// U is a (usually) non-unit upper triangular matrix. On exit, L and U are stored
// in place into a.
//
// ipiv is a permutation vector. It indicates that row i of the matrix was
// changed with ipiv[i]. ipiv must have length at least min(m,n), and will panic
// otherwise. ipiv is zero-indexed.
//
// Dgetf2 returns whether the matrix A is singular. The LU decomposition will
// be computed regardless of the singularity of A, but division by zero
// will occur if the false is returned and the result is used to solve a
// system of equations.
func (Implementation) Dgetf2(m, n int, a []float64, lda int, ipiv []int) (ok bool) {
mn := min(m, n)
checkMatrix(m, n, a, lda)
if len(ipiv) < mn {
panic(badIpiv)
}
ipiv32 := make([]int32, len(ipiv))
ok = clapack.Dgetf2(m, n, a, lda, ipiv32)
for i, v := range ipiv32 {
ipiv[i] = int(v) - 1 // Transform to zero-indexed.
}
return ok
}
// Dgetrf computes the LU decomposition of the m×n matrix A.
// The LU decomposition is a factorization of a into
// A = P * L * U
// where P is a permutation matrix, L is a unit lower triangular matrix, and
// U is a (usually) non-unit upper triangular matrix. On exit, L and U are stored
// in place into a.
//
// ipiv is a permutation vector. It indicates that row i of the matrix was
// changed with ipiv[i]. ipiv must have length at least min(m,n), and will panic
// otherwise. ipiv is zero-indexed.
//
// Dgetrf is the blocked version of the algorithm.
//
// Dgetrf returns whether the matrix A is singular. The LU decomposition will
// be computed regardless of the singularity of A, but division by zero
// will occur if the false is returned and the result is used to solve a
// system of equations.
func (impl Implementation) Dgetrf(m, n int, a []float64, lda int, ipiv []int) (ok bool) {
mn := min(m, n)
checkMatrix(m, n, a, lda)
if len(ipiv) < mn {
panic(badIpiv)
}
ipiv32 := make([]int32, len(ipiv))
ok = clapack.Dgetrf(m, n, a, lda, ipiv32)
for i, v := range ipiv32 {
ipiv[i] = int(v) - 1 // Transform to zero-indexed.
}
return ok
}
// Dgetrs solves a system of equations using an LU factorization.
// The system of equations solved is
// A * X = B if trans == blas.Trans
// A^T * X = B if trans == blas.NoTrans
// A is a general n×n matrix with stride lda. B is a general matrix of size n×nrhs.
//
// On entry b contains the elements of the matrix B. On exit, b contains the
// elements of X, the solution to the system of equations.
//
// a and ipiv contain the LU factorization of A and the permutation indices as
// computed by Dgetrf. ipiv is zero-indexed.
func (impl Implementation) Dgetrs(trans blas.Transpose, n, nrhs int, a []float64, lda int, ipiv []int, b []float64, ldb int) {
checkMatrix(n, n, a, lda)
checkMatrix(n, nrhs, b, ldb)
if len(ipiv) < n {
panic(badIpiv)
}
ipiv32 := make([]int32, len(ipiv))
for i, v := range ipiv {
ipiv32[i] = int32(v) + 1 // Transform to one-indexed.
}
clapack.Dgetrs(trans, n, nrhs, a, lda, ipiv32, b, ldb)
}