forked from cpmech/gosl
-
Notifications
You must be signed in to change notification settings - Fork 0
/
matrix.go
615 lines (561 loc) · 14.9 KB
/
matrix.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
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
// Copyright 2016 The Gosl 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 la implements functions and structure for Linear Algebra computations. It defines a
// Vector and Matrix types for computations with dense data and also a Triplet and CCMatrix for
// sparse data.
package la
import (
"math"
"strings"
"github.com/cpmech/gosl/chk"
"github.com/cpmech/gosl/io"
"github.com/cpmech/gosl/la/oblas"
"github.com/cpmech/gosl/utl"
)
// Matrix implements a column-major representation of a matrix by using a linear array that can be passed to Fortran code
//
// NOTE: the functions related to Matrix do not check for the limits of indices and dimensions.
// Panic may occur then.
//
// Example:
// _ _
// | 0 3 |
// A = | 1 4 |
// |_ 2 5 _|(m x n)
//
// data[i+j*m] = A[i][j]
//
type Matrix struct {
M, N int // dimensions
Data []float64 // data array. column-major => Fortran
}
// NewMatrix allocates a new (empty) Matrix with given (m,n) (row/col sizes)
func NewMatrix(m, n int) (o *Matrix) {
o = new(Matrix)
o.M, o.N = m, n
o.Data = make([]float64, m*n)
return
}
// NewMatrixDeep2 allocates a new Matrix from given (Deep2) nested slice.
// NOTE: make sure to have at least 1x1 item
func NewMatrixDeep2(a [][]float64) (o *Matrix) {
o = new(Matrix)
o.M, o.N = len(a), len(a[0])
o.Data = make([]float64, o.M*o.N)
o.SetFromDeep2(a)
return
}
// NewMatrixRaw creates a new Matrix using given raw data
// Input:
// rawdata -- data organized as column-major; e.g. Fortran format
// NOTE:
// (1) rawdata is not copied!
// (2) the external slice rawdata should not be changed or deleted
func NewMatrixRaw(m, n int, rawdata []float64) (o *Matrix) {
o = new(Matrix)
o.M, o.N = m, n
o.Data = rawdata
return
}
// SetFromDeep2 sets matrix with data from a nested slice (Deep2) structure
func (o *Matrix) SetFromDeep2(a [][]float64) {
k := 0
for j := 0; j < o.N; j++ {
for i := 0; i < o.M; i++ {
o.Data[k] = a[i][j]
k++
}
}
}
// SetDiag sets diagonal matrix with diagonal components equal to val
func (o *Matrix) SetDiag(val float64) {
for i := 0; i < o.M; i++ {
for j := 0; j < o.N; j++ {
if i == j {
o.Data[i+j*o.M] = val
} else {
o.Data[i+j*o.M] = 0
}
}
}
}
// Set sets value
func (o *Matrix) Set(i, j int, val float64) {
o.Data[i+j*o.M] = val // col-major
}
// Get gets value
func (o *Matrix) Get(i, j int) float64 {
return o.Data[i+j*o.M] // col-major
}
// GetDeep2 returns nested slice representation
func (o *Matrix) GetDeep2() (M [][]float64) {
M = make([][]float64, o.M)
for i := 0; i < o.M; i++ {
M[i] = make([]float64, o.N)
for j := 0; j < o.N; j++ {
M[i][j] = o.Data[i+j*o.M]
}
}
return
}
// GetCopy returns a copy of this matrix
func (o *Matrix) GetCopy() (clone *Matrix) {
clone = NewMatrix(o.M, o.N)
copy(clone.Data, o.Data)
return
}
// GetTranspose returns the tranpose matrix
func (o *Matrix) GetTranspose() (tran *Matrix) {
tran = NewMatrix(o.N, o.M)
for i := 0; i < o.N; i++ {
for j := 0; j < o.M; j++ {
tran.Set(i, j, o.Get(j, i))
}
}
return
}
// GetComplex returns a complex version of this matrix
func (o *Matrix) GetComplex() (b *MatrixC) {
b = NewMatrixC(o.M, o.N)
for i := 0; i < o.M; i++ {
for j := 0; j < o.N; j++ {
b.Set(i, j, complex(o.Get(i, j), 0))
}
}
return
}
// CopyInto copies the scaled components of this matrix into another one (result)
// result := α * this ⇒ result[ij] := α * this[ij]
func (o *Matrix) CopyInto(result *Matrix, α float64) {
for k := 0; k < o.M*o.N; k++ {
result.Data[k] = α * o.Data[k]
}
}
// Add adds value to (i,j) location
func (o *Matrix) Add(i, j int, val float64) {
o.Data[i+j*o.M] += val // col-major
}
// Fill fills this matrix with a single number val
// aij = val
func (o *Matrix) Fill(val float64) {
for k := 0; k < o.M*o.N; k++ {
o.Data[k] = val
}
}
// ClearRC clear rows and columns and set diagonal components
// _ _ _ _
// Example: | 1 2 3 4 | | 1 2 3 4 |
// A = | 5 6 7 8 | ⇒ clear([1,2], [], 1.0) ⇒ A = | 0 1 0 0 |
// |_ 4 3 2 1 _| |_ 0 0 1 0 _|
//
func (o *Matrix) ClearRC(rows, cols []int, diag float64) {
for _, r := range rows {
for j := 0; j < o.N; j++ {
if r == j {
o.Set(r, j, diag)
} else {
o.Set(r, j, 0.0)
}
}
}
for _, c := range cols {
for i := 0; i < o.M; i++ {
if i == c {
o.Set(i, c, diag)
} else {
o.Set(i, c, 0.0)
}
}
}
}
// ClearBry clears boundaries
// _ _ _ _
// Example: | 1 2 3 | | 1 0 0 |
// A = | 4 5 6 | ⇒ clear(1.0) ⇒ A = | 0 5 0 |
// |_ 7 8 9 _| |_ 0 0 1 _|
//
func (o *Matrix) ClearBry(diag float64) {
o.ClearRC([]int{0, o.M - 1}, []int{0, o.N - 1}, diag)
}
// MaxDiff returns the maximum difference between the components of this and another matrix
func (o *Matrix) MaxDiff(another *Matrix) (maxdiff float64) {
maxdiff = math.Abs(o.Data[0] - another.Data[0])
for k := 1; k < o.M*o.N; k++ {
diff := math.Abs(o.Data[k] - another.Data[k])
if diff > maxdiff {
maxdiff = diff
}
}
return
}
// Largest returns the largest component |a[ij]| of this matrix, normalised by den
// largest := |a[ij]| / den
func (o *Matrix) Largest(den float64) (largest float64) {
largest = math.Abs(o.Data[0])
for k := 1; k < o.M*o.N; k++ {
tmp := math.Abs(o.Data[k])
if tmp > largest {
largest = tmp
}
}
return largest / den
}
// Col access column j of this matrix. No copies are made since the internal data are in
// col-major format already.
// NOTE: this method can be used to modify the columns; e.g. with o.Col(0)[0] = 123
func (o *Matrix) Col(j int) Vector {
return o.Data[j*o.M : (j+1)*o.M]
}
// GetRow returns row i of this matrix
func (o *Matrix) GetRow(i int) (row Vector) {
row = make([]float64, o.N)
for j := 0; j < o.N; j++ {
row[j] = o.Data[i+j*o.M]
}
return
}
// GetCol returns column j of this matrix
func (o *Matrix) GetCol(j int) (col Vector) {
col = make([]float64, o.M)
copy(col, o.Data[j*o.M:(j+1)*o.M])
return
}
// ExtractCols returns columns from j=start to j=endp1-1
// start -- first column
// endp1 -- "end-plus-one", the number of the last requested column + 1
func (o *Matrix) ExtractCols(start, endp1 int) (reduced *Matrix) {
if endp1 <= start {
chk.Panic("endp1 'end-plus-one' must be greater than start. start=%d, endp1=%d invalid\n", start, endp1)
}
ncolNew := endp1 - start
reduced = NewMatrix(o.M, ncolNew)
copy(reduced.Data, o.Data[start*o.M:endp1*o.M])
return
}
// SetCol sets the values of a column j with a single value
func (o *Matrix) SetCol(j int, value float64) {
for k := j * o.M; k < (j+1)*o.M; k++ {
o.Data[k] = value
}
}
// NormFrob returns the Frobenious norm of this matrix
// nrm := ‖a‖_F = sqrt(Σ_i Σ_j a[ij]⋅a[ij]) = ‖a‖_2
func (o *Matrix) NormFrob() (nrm float64) {
for k := 0; k < o.M*o.N; k++ {
nrm += o.Data[k] * o.Data[k]
}
return math.Sqrt(nrm)
}
// NormInf returns the infinite norm of this matrix
// nrm := ‖a‖_∞ = max_i ( Σ_j a[ij] )
func (o *Matrix) NormInf() (nrm float64) {
for j := 0; j < o.N; j++ { // sum first row
nrm += math.Abs(o.Data[j*o.M])
}
var sumrow float64
for i := 1; i < o.M; i++ {
sumrow = 0.0
for j := 0; j < o.N; j++ { // sum the other rows
sumrow += math.Abs(o.Data[i+j*o.M])
if sumrow > nrm {
nrm = sumrow
}
}
}
return
}
// Apply sets this matrix with the scaled components of another matrix
// this := α * another ⇒ this[i] := α * another[i]
// NOTE: "another" may be "this"
func (o Matrix) Apply(α float64, another *Matrix) {
for k := 0; k < o.M*o.N; k++ {
o.Data[k] = α * another.Data[k]
}
}
// Det computes the determinant of matrix using the LU factorization
// NOTE: this method may fail due to overflow...
func (o *Matrix) Det() (det float64) {
if o.M != o.N {
chk.Panic("matrix must be square to compute determinant. %d × %d is invalid\n", o.M, o.N)
}
ai := make([]float64, len(o.Data))
copy(ai, o.Data)
ipiv := make([]int32, utl.Imin(o.M, o.N))
oblas.Dgetrf(o.M, o.N, ai, o.M, ipiv) // NOTE: ipiv are 1-based indices
det = 1.0
for i := 0; i < o.M; i++ {
if ipiv[i]-1 == int32(i) { // NOTE: ipiv are 1-based indices
det = +det * ai[i+i*o.M]
} else {
det = -det * ai[i+i*o.M]
}
}
return
}
// Print prints matrix (without commas or brackets)
func (o *Matrix) Print(nfmt string) (l string) {
if nfmt == "" {
nfmt = "%g "
}
for i := 0; i < o.M; i++ {
if i > 0 {
l += "\n"
}
for j := 0; j < o.N; j++ {
l += io.Sf(nfmt, o.Get(i, j))
}
}
return
}
// PrintGo prints matrix in Go format
func (o *Matrix) PrintGo(nfmt string) (l string) {
if nfmt == "" {
nfmt = "%10g"
}
l = "[][]float64{\n"
for i := 0; i < o.M; i++ {
l += " {"
for j := 0; j < o.N; j++ {
if j > 0 {
l += ","
}
l += io.Sf(nfmt, o.Get(i, j))
}
l += "},\n"
}
l += "}"
return
}
// PrintPy prints matrix in Python format
func (o *Matrix) PrintPy(nfmt string) (l string) {
if nfmt == "" {
nfmt = "%10g"
}
l = "np.matrix([\n"
for i := 0; i < o.M; i++ {
l += " ["
for j := 0; j < o.N; j++ {
if j > 0 {
l += ","
}
l += io.Sf(nfmt, o.Get(i, j))
}
l += "],\n"
}
l += "], dtype=float)"
return
}
// complex ///////////////////////////////////////////////////////////////////////////////////////
// MatrixC implements a column-major representation of a matrix of complex numbers by using a linear
// array that can be passed to Fortran code.
//
// NOTE: the functions related to MatrixC do not check for the limits of indices and dimensions.
// Panic may occur then.
//
// Example:
// _ _
// | 0+0i 3+3i |
// A = | 1+1i 4+4i |
// |_ 2+2i 5+5i _|(m x n)
//
// data[i+j*m] = A[i][j]
//
type MatrixC struct {
M, N int // dimensions
Data []complex128 // data array. column-major => Fortran
}
// NewMatrixC allocates a new (empty) MatrixC with given (m,n) (row/col sizes)
func NewMatrixC(m, n int) (o *MatrixC) {
o = new(MatrixC)
o.M, o.N = m, n
o.Data = make([]complex128, m*n)
return
}
// NewMatrixDeep2c allocates a new MatrixC from given (Deep2c) nested slice.
// NOTE: make sure to have at least 1x1 items
func NewMatrixDeep2c(a [][]complex128) (o *MatrixC) {
o = new(MatrixC)
o.M, o.N = len(a), len(a[0])
o.Data = make([]complex128, o.M*o.N)
o.SetFromDeep2c(a)
return
}
// SetFromDeep2c sets matrix with data from a nested slice (Deep2c) structure
func (o *MatrixC) SetFromDeep2c(a [][]complex128) {
k := 0
for j := 0; j < o.N; j++ {
for i := 0; i < o.M; i++ {
o.Data[k] = a[i][j]
k++
}
}
}
// Set sets value
func (o *MatrixC) Set(i, j int, val complex128) {
o.Data[i+j*o.M] = val // col-major
}
// Get gets value
func (o *MatrixC) Get(i, j int) complex128 {
return o.Data[i+j*o.M] // col-major
}
// GetDeep2 returns nested slice representation
func (o *MatrixC) GetDeep2() (M [][]complex128) {
M = make([][]complex128, o.M)
for i := 0; i < o.M; i++ {
M[i] = make([]complex128, o.N)
for j := 0; j < o.N; j++ {
M[i][j] = o.Data[i+j*o.M]
}
}
return
}
// GetCopy returns a copy of this matrix
func (o *MatrixC) GetCopy() (clone *MatrixC) {
clone = NewMatrixC(o.M, o.N)
copy(clone.Data, o.Data)
return
}
// GetTranspose returns the tranpose matrix
func (o *MatrixC) GetTranspose() (tran *MatrixC) {
tran = NewMatrixC(o.N, o.M)
for i := 0; i < o.N; i++ {
for j := 0; j < o.M; j++ {
tran.Set(i, j, o.Get(j, i))
}
}
return
}
// Add adds value to (i,j) location
func (o *MatrixC) Add(i, j int, val complex128) {
o.Data[i+j*o.M] += val // col-major
}
// Fill fills this matrix with a single number val
// aij = val
func (o *MatrixC) Fill(val complex128) {
for k := 0; k < o.M*o.N; k++ {
o.Data[k] = val
}
}
// Col access column j of this matrix. No copies are made since the internal data are in
// col-major format already.
// NOTE: this method can be used to modify the columns; e.g. with o.Col(0)[0] = 123
func (o *MatrixC) Col(j int) VectorC {
return o.Data[j*o.M : (j+1)*o.M]
}
// GetRow returns row i of this matrix
func (o *MatrixC) GetRow(i int) (row VectorC) {
row = make([]complex128, o.N)
for j := 0; j < o.N; j++ {
row[j] = o.Data[i+j*o.M]
}
return
}
// GetCol returns column j of this matrix
func (o *MatrixC) GetCol(j int) (col VectorC) {
col = make([]complex128, o.M)
copy(col, o.Data[j*o.M:(j+1)*o.M])
return
}
// GetColReal returns column j of this matrix considering that the imaginary part is zero
func (o *MatrixC) GetColReal(j int, checkZeroImag bool) (col Vector) {
col = make([]float64, o.M)
if checkZeroImag {
for i := 0; i < o.M; i++ {
if math.Abs(imag(o.Data[i+j*o.M])) > 0.0 {
chk.Panic("imaginary part is not zero\n")
}
col[i] = real(o.Data[i+j*o.M])
}
return
}
for i := 0; i < o.M; i++ {
col[i] = real(o.Data[i+j*o.M])
}
return
}
// Apply sets this matrix with the scaled components of another matrix
// this := α * another ⇒ this[i] := α * another[i]
// NOTE: "another" may be "this"
func (o MatrixC) Apply(α complex128, another *MatrixC) {
for k := 0; k < o.M*o.N; k++ {
o.Data[k] = α * another.Data[k]
}
}
// Print prints matrix (without commas or brackets).
// NOTE: if non-empty, nfmtI must have '+' e.g. %+g
func (o *MatrixC) Print(nfmtR, nfmtI string) (l string) {
if nfmtR == "" {
nfmtR = "%g"
}
if nfmtI == "" {
nfmtI = "%+g"
}
if !strings.ContainsRune(nfmtI, '+') {
nfmtI = strings.Replace(nfmtI, "%", "%+", -1)
}
for i := 0; i < o.M; i++ {
if i > 0 {
l += "\n"
}
for j := 0; j < o.N; j++ {
if j > 0 {
l += " "
}
l += io.Sf(nfmtR, real(o.Get(i, j))) + io.Sf(nfmtI, imag(o.Get(i, j))) + "i"
}
}
return
}
// PrintGo prints matrix in Go format
// NOTE: if non-empty, nfmtI must have '+' e.g. %+g
func (o *MatrixC) PrintGo(nfmtR, nfmtI string) (l string) {
if nfmtR == "" {
nfmtR = "%g"
}
if nfmtI == "" {
nfmtI = "%+g"
}
if !strings.ContainsRune(nfmtI, '+') {
nfmtI = strings.Replace(nfmtI, "%", "%+", -1)
}
l = "[][]complex128{\n"
for i := 0; i < o.M; i++ {
l += " {"
for j := 0; j < o.N; j++ {
if j > 0 {
l += ","
}
l += io.Sf(nfmtR, real(o.Get(i, j))) + io.Sf(nfmtI, imag(o.Get(i, j))) + "i"
}
l += "},\n"
}
l += "}"
return
}
// PrintPy prints matrix in Python format
// NOTE: if non-empty, nfmtI must have '+' e.g. %+g
func (o *MatrixC) PrintPy(nfmtR, nfmtI string) (l string) {
if nfmtR == "" {
nfmtR = "%g"
}
if nfmtI == "" {
nfmtI = "%+g"
}
if !strings.ContainsRune(nfmtI, '+') {
nfmtI = strings.Replace(nfmtI, "%", "%+", -1)
}
l = "np.matrix([\n"
for i := 0; i < o.M; i++ {
l += " ["
for j := 0; j < o.N; j++ {
if j > 0 {
l += ","
}
l += io.Sf(nfmtR, real(o.Get(i, j))) + io.Sf(nfmtI, imag(o.Get(i, j))) + "j"
}
l += "],\n"
}
l += "], dtype=complex)"
return
}