-
Notifications
You must be signed in to change notification settings - Fork 27
/
lapack.go
2982 lines (2866 loc) · 102 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
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
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// 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 (
"math"
"github.com/gonum/blas"
"github.com/gonum/lapack"
"github.com/gonum/lapack/cgo/lapacke"
)
// Copied from lapack/native. Keep in sync.
const (
absIncNotOne = "lapack: increment not one or negative one"
badAlpha = "lapack: bad alpha length"
badAuxv = "lapack: auxv has insufficient length"
badBeta = "lapack: bad beta length"
badD = "lapack: d has insufficient length"
badDecompUpdate = "lapack: bad decomp update"
badDiag = "lapack: bad diag"
badDims = "lapack: bad input dimensions"
badDirect = "lapack: bad direct"
badE = "lapack: e has insufficient length"
badEVComp = "lapack: bad EVComp"
badEVJob = "lapack: bad EVJob"
badEVSide = "lapack: bad EVSide"
badGSVDJob = "lapack: bad GSVDJob"
badHowMany = "lapack: bad HowMany"
badIlo = "lapack: ilo out of range"
badIhi = "lapack: ihi out of range"
badIpiv = "lapack: bad permutation length"
badJob = "lapack: bad Job"
badK1 = "lapack: k1 out of range"
badK2 = "lapack: k2 out of range"
badKperm = "lapack: incorrect permutation length"
badLdA = "lapack: index of a out of range"
badNb = "lapack: nb out of range"
badNorm = "lapack: bad norm"
badPivot = "lapack: bad pivot"
badS = "lapack: s has insufficient length"
badShifts = "lapack: bad shifts"
badSide = "lapack: bad side"
badSlice = "lapack: bad input slice length"
badSort = "lapack: bad Sort"
badStore = "lapack: bad store"
badTau = "lapack: tau has insufficient length"
badTauQ = "lapack: tauQ has insufficient length"
badTauP = "lapack: tauP has insufficient length"
badTrans = "lapack: bad trans"
badVn1 = "lapack: vn1 has insufficient length"
badVn2 = "lapack: vn2 has insufficient length"
badUplo = "lapack: illegal triangle"
badWork = "lapack: insufficient working memory"
badWorkStride = "lapack: insufficient working array stride"
badZ = "lapack: insufficient z length"
kGTM = "lapack: k > m"
kGTN = "lapack: k > n"
kLT0 = "lapack: k < 0"
mLT0 = "lapack: m < 0"
mLTN = "lapack: m < n"
nanScale = "lapack: NaN scale factor"
negDimension = "lapack: negative matrix dimension"
negZ = "lapack: negative z value"
nLT0 = "lapack: n < 0"
nLTM = "lapack: n < m"
offsetGTM = "lapack: offset > m"
shortWork = "lapack: working array shorter than declared"
zeroDiv = "lapack: zero divisor"
)
func min(m, n int) int {
if m < n {
return m
}
return n
}
func max(m, n int) int {
if m < n {
return n
}
return m
}
// 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 n < 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")
}
}
// checkVector verifies the parameters of a vector input.
// Copied from lapack/native. Keep in sync.
func checkVector(n int, v []float64, inc int) {
if n < 0 {
panic("lapack: negative vector length")
}
if (inc > 0 && (n-1)*inc >= len(v)) || (inc < 0 && (1-n)*inc >= len(v)) {
panic("lapack: insufficient vector slice length")
}
}
// Implementation is the cgo-based C implementation of LAPACK routines.
type Implementation struct{}
var _ lapack.Float64 = Implementation{}
// Dgeqp3 computes a QR factorization with column pivoting of the
// m×n matrix A: A*P = Q*R using Level 3 BLAS.
//
// The matrix Q is represented as a product of elementary reflectors
// Q = H_0 H_1 . . . H_{k-1}, where k = min(m,n).
// Each H_i has the form
// H_i = I - tau * v * v^T
// where tau and v are real vectors with v[0:i-1] = 0 and v[i] = 1;
// v[i:m] is stored on exit in A[i:m, i], and tau in tau[i].
//
// jpvt specifies a column pivot to be applied to A. If
// jpvt[j] is at least zero, the jth column of A is permuted
// to the front of A*P (a leading column), if jpvt[j] is -1
// the jth column of A is a free column. If jpvt[j] < -1, Dgeqp3
// will panic. On return, jpvt holds the permutation that was
// applied; the jth column of A*P was the jpvt[j] column of A.
// jpvt must have length n or Dgeqp3 will panic.
//
// tau holds the scalar factors of the elementary reflectors.
// It must have length min(m, n), otherwise Dgeqp3 will panic.
//
// work must have length at least max(1,lwork), and lwork must be at least
// 3*n+1, otherwise Dgeqp3 will panic. For optimal performance lwork must
// be at least 2*n+(n+1)*nb, where nb is the optimal blocksize. On return,
// work[0] will contain the optimal value of lwork.
//
// If lwork == -1, instead of performing Dgeqp3, only the optimal value of lwork
// will be stored in work[0].
//
// Dgeqp3 is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dgeqp3(m, n int, a []float64, lda int, jpvt []int, tau, work []float64, lwork int) {
checkMatrix(m, n, a, lda)
if len(jpvt) != n {
panic(badIpiv)
}
if len(tau) != min(m, n) {
panic(badTau)
}
if len(work) < max(1, lwork) {
panic(badWork)
}
// Don't update jpvt if querying lwkopt.
if lwork == -1 {
lapacke.Dgeqp3(m, n, a, lda, nil, nil, work, -1)
return
}
jpvt32 := make([]int32, len(jpvt))
for i, v := range jpvt {
v++
if v != int(int32(v)) || v < 0 || n < v {
panic("lapack: jpvt element out of range")
}
jpvt32[i] = int32(v)
}
lapacke.Dgeqp3(m, n, a, lda, jpvt32, tau, work, lwork)
for i, v := range jpvt32 {
jpvt[i] = int(v - 1)
}
}
// Dgerqf computes an RQ factorization of the m×n matrix A,
// A = R * Q.
// On exit, if m <= n, the upper triangle of the subarray
// A[0:m, n-m:n] contains the m×m upper triangular matrix R.
// If m >= n, the elements on and above the (m-n)-th subdiagonal
// contain the m×n upper trapezoidal matrix R.
// The remaining elements, with tau, represent the
// orthogonal matrix Q as a product of min(m,n) elementary
// reflectors.
//
// The matrix Q is represented as a product of elementary reflectors
// Q = H_0 H_1 . . . H_{min(m,n)-1}.
// Each H(i) has the form
// H_i = I - tau_i * v * v^T
// where v is a vector with v[0:n-k+i-1] stored in A[m-k+i, 0:n-k+i-1],
// v[n-k+i:n] = 0 and v[n-k+i] = 1.
//
// tau must have length min(m,n), work must have length max(1, lwork),
// and lwork must be -1 or at least max(1, m), otherwise Dgerqf will panic.
// On exit, work[0] will contain the optimal length for work.
//
// Dgerqf is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dgerqf(m, n int, a []float64, lda int, tau, work []float64, lwork int) {
checkMatrix(m, n, a, lda)
if len(work) < max(1, lwork) {
panic(shortWork)
}
if lwork != -1 && lwork < max(1, m) {
panic(badWork)
}
k := min(m, n)
if len(tau) != k {
panic(badTau)
}
lapacke.Dgerqf(m, n, a, lda, tau, work, lwork)
}
// Dlacn2 estimates the 1-norm of an n×n matrix A using sequential updates with
// matrix-vector products provided externally.
//
// Dlacn2 is called sequentially and it returns the value of est and kase to be
// used on the next call.
// On the initial call, kase must be 0.
// In between calls, x must be overwritten by
// A * X if kase was returned as 1,
// A^T * X if kase was returned as 2,
// and all other parameters must not be changed.
// On the final return, kase is returned as 0, v contains A*W where W is a
// vector, and est = norm(V)/norm(W) is a lower bound for 1-norm of A.
//
// v, x, and isgn must all have length n and n must be at least 1, otherwise
// Dlacn2 will panic. isave is used for temporary storage.
//
// Dlacn2 is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dlacn2(n int, v, x []float64, isgn []int, est float64, kase int, isave *[3]int) (float64, int) {
if n < 1 {
panic("lapack: non-positive n")
}
checkVector(n, x, 1)
checkVector(n, v, 1)
if len(isgn) < n {
panic("lapack: insufficient isgn length")
}
if isave[0] < 0 || isave[0] > 5 {
panic("lapack: bad isave value")
}
if isave[0] == 0 && kase != 0 {
panic("lapack: bad isave value")
}
isgn32 := make([]int32, n)
for i, v := range isgn {
isgn32[i] = int32(v)
}
pest := []float64{est}
// Save one allocation by putting isave and kase into the same slice.
isavekase := []int32{int32(isave[0]), int32(isave[1]), int32(isave[2]), int32(kase)}
lapacke.Dlacn2(n, v, x, isgn32, pest, isavekase[3:], isavekase[:3])
for i, v := range isgn32 {
isgn[i] = int(v)
}
isave[0] = int(isavekase[0])
isave[1] = int(isavekase[1])
isave[2] = int(isavekase[2])
return pest[0], int(isavekase[3])
}
// Dlacpy copies the elements of A specified by uplo into B. Uplo can specify
// a triangular portion with blas.Upper or blas.Lower, or can specify all of the
// elemest with blas.All.
func (impl Implementation) Dlacpy(uplo blas.Uplo, m, n int, a []float64, lda int, b []float64, ldb int) {
checkMatrix(m, n, a, lda)
checkMatrix(m, n, b, ldb)
lapacke.Dlacpy(uplo, m, n, a, lda, b, ldb)
}
// Dlapmt rearranges the columns of the m×n matrix X as specified by the
// permutation k_0, k_1, ..., k_n-1 of the integers 0, ..., n-1.
//
// If forward is true a forward permutation is performed:
//
// X[0:m, k[j]] is moved to X[0:m, j] for j = 0, 1, ..., n-1.
//
// otherwise a backward permutation is performed:
//
// X[0:m, j] is moved to X[0:m, k[j]] for j = 0, 1, ..., n-1.
//
// k must have length n, otherwise Dlapmt will panic. k is zero-indexed.
func (impl Implementation) Dlapmt(forward bool, m, n int, x []float64, ldx int, k []int) {
checkMatrix(m, n, x, ldx)
if len(k) != n {
panic(badKperm)
}
if n <= 1 {
return
}
var forwrd int32
if forward {
forwrd = 1
}
k32 := make([]int32, len(k))
for i, v := range k {
v++ // Convert to 1-based indexing.
if v != int(int32(v)) {
panic("lapack: k element out of range")
}
k32[i] = int32(v)
}
lapacke.Dlapmt(forwrd, m, n, x, ldx, k32)
}
// Dlapy2 is the LAPACK version of math.Hypot.
//
// Dlapy2 is an internal routine. It is exported for testing purposes.
func (Implementation) Dlapy2(x, y float64) float64 {
return lapacke.Dlapy2(x, y)
}
// Dlarfb applies a block reflector to a matrix.
//
// In the call to Dlarfb, the mxn c is multiplied by the implicitly defined matrix h as follows:
// c = h * c if side == Left and trans == NoTrans
// c = c * h if side == Right and trans == NoTrans
// c = h^T * c if side == Left and trans == Trans
// c = c * h^T if side == Right and trans == Trans
// h is a product of elementary reflectors. direct sets the direction of multiplication
// h = h_1 * h_2 * ... * h_k if direct == Forward
// h = h_k * h_k-1 * ... * h_1 if direct == Backward
// The combination of direct and store defines the orientation of the elementary
// reflectors. In all cases the ones on the diagonal are implicitly represented.
//
// If direct == lapack.Forward and store == lapack.ColumnWise
// V = [ 1 ]
// [v1 1 ]
// [v1 v2 1]
// [v1 v2 v3]
// [v1 v2 v3]
// If direct == lapack.Forward and store == lapack.RowWise
// V = [ 1 v1 v1 v1 v1]
// [ 1 v2 v2 v2]
// [ 1 v3 v3]
// If direct == lapack.Backward and store == lapack.ColumnWise
// V = [v1 v2 v3]
// [v1 v2 v3]
// [ 1 v2 v3]
// [ 1 v3]
// [ 1]
// If direct == lapack.Backward and store == lapack.RowWise
// V = [v1 v1 1 ]
// [v2 v2 v2 1 ]
// [v3 v3 v3 v3 1]
// An elementary reflector can be explicitly constructed by extracting the
// corresponding elements of v, placing a 1 where the diagonal would be, and
// placing zeros in the remaining elements.
//
// t is a k×k matrix containing the block reflector, and this function will panic
// if t is not of sufficient size. See Dlarft for more information.
//
// work is a temporary storage matrix with stride ldwork.
// work must be of size at least n×k side == Left and m×k if side == Right, and
// this function will panic if this size is not met.
//
// Dlarfb is an internal routine. It is exported for testing purposes.
func (Implementation) Dlarfb(side blas.Side, trans blas.Transpose, direct lapack.Direct, store lapack.StoreV, m, n, k int, v []float64, ldv int, t []float64, ldt int, c []float64, ldc int, work []float64, ldwork int) {
if side != blas.Left && side != blas.Right {
panic(badSide)
}
if trans != blas.Trans && trans != blas.NoTrans {
panic(badTrans)
}
if direct != lapack.Forward && direct != lapack.Backward {
panic(badDirect)
}
if store != lapack.ColumnWise && store != lapack.RowWise {
panic(badStore)
}
checkMatrix(m, n, c, ldc)
if k < 0 {
panic(kLT0)
}
checkMatrix(k, k, t, ldt)
nv := m
nw := n
if side == blas.Right {
nv = n
nw = m
}
if store == lapack.ColumnWise {
checkMatrix(nv, k, v, ldv)
} else {
checkMatrix(k, nv, v, ldv)
}
// TODO(vladimir-ch): Replace the following two lines with
// checkMatrix(nw, k, work, ldwork)
// if and when the issue
// https://github.com/Reference-LAPACK/lapack/issues/37
// has been resolved.
ldwork = nw
work = make([]float64, ldwork*k)
lapacke.Dlarfb(side, trans, byte(direct), byte(store), m, n, k, v, ldv, t, ldt, c, ldc, work, ldwork)
}
// Dlarfg generates an elementary reflector for a Householder matrix. It creates
// a real elementary reflector of order n such that
// H * (alpha) = (beta)
// ( x) ( 0)
// H^T * H = I
// H is represented in the form
// H = 1 - tau * (1; v) * (1 v^T)
// where tau is a real scalar.
//
// On entry, x contains the vector x, on exit it contains v.
//
// Dlarfg is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dlarfg(n int, alpha float64, x []float64, incX int) (beta, tau float64) {
if n < 0 {
panic(nLT0)
}
if n <= 1 {
return alpha, 0
}
checkVector(n-1, x, incX)
_alpha := []float64{alpha}
_tau := []float64{0}
lapacke.Dlarfg(n, _alpha, x, incX, _tau)
return _alpha[0], _tau[0]
}
// Dlarft forms the triangular factor T of a block reflector H, storing the answer
// in t.
// H = I - V * T * V^T if store == lapack.ColumnWise
// H = I - V^T * T * V if store == lapack.RowWise
// H is defined by a product of the elementary reflectors where
// H = H_0 * H_1 * ... * H_{k-1} if direct == lapack.Forward
// H = H_{k-1} * ... * H_1 * H_0 if direct == lapack.Backward
//
// t is a k×k triangular matrix. t is upper triangular if direct = lapack.Forward
// and lower triangular otherwise. This function will panic if t is not of
// sufficient size.
//
// store describes the storage of the elementary reflectors in v. Please see
// Dlarfb for a description of layout.
//
// tau contains the scalar factors of the elementary reflectors H_i.
//
// Dlarft is an internal routine. It is exported for testing purposes.
func (Implementation) Dlarft(direct lapack.Direct, store lapack.StoreV, n, k int,
v []float64, ldv int, tau []float64, t []float64, ldt int) {
if n == 0 {
return
}
if n < 0 || k < 0 {
panic(negDimension)
}
if direct != lapack.Forward && direct != lapack.Backward {
panic(badDirect)
}
if store != lapack.RowWise && store != lapack.ColumnWise {
panic(badStore)
}
if len(tau) < k {
panic(badTau)
}
checkMatrix(k, k, t, ldt)
lapacke.Dlarft(byte(direct), byte(store), n, k, v, ldv, tau, t, ldt)
}
// Dlange computes the matrix norm of the general m×n matrix a. The input norm
// specifies the norm computed.
// lapack.MaxAbs: the maximum absolute value of an element.
// lapack.MaxColumnSum: the maximum column sum of the absolute values of the entries.
// lapack.MaxRowSum: the maximum row sum of the absolute values of the entries.
// lapack.NormFrob: the square root of the sum of the squares of the entries.
// If norm == lapack.MaxColumnSum, work must be of length n, and this function will panic otherwise.
// There are no restrictions on work for the other matrix norms.
func (impl Implementation) Dlange(norm lapack.MatrixNorm, m, n int, a []float64, lda int, work []float64) float64 {
checkMatrix(m, n, a, lda)
switch norm {
case lapack.MaxRowSum, lapack.MaxColumnSum, lapack.NormFrob, lapack.MaxAbs:
default:
panic(badNorm)
}
if norm == lapack.MaxColumnSum && len(work) < n {
panic(badWork)
}
return lapacke.Dlange(byte(norm), m, n, a, lda, work)
}
// Dlansy computes the specified norm of an n×n symmetric matrix. If
// norm == lapack.MaxColumnSum or norm == lapackMaxRowSum work must have length
// at least n, otherwise work is unused.
func (impl Implementation) Dlansy(norm lapack.MatrixNorm, uplo blas.Uplo, n int, a []float64, lda int, work []float64) float64 {
checkMatrix(n, n, a, lda)
switch norm {
case lapack.MaxRowSum, lapack.MaxColumnSum, lapack.NormFrob, lapack.MaxAbs:
default:
panic(badNorm)
}
if (norm == lapack.MaxColumnSum || norm == lapack.MaxRowSum) && len(work) < n {
panic(badWork)
}
if uplo != blas.Upper && uplo != blas.Lower {
panic(badUplo)
}
return lapacke.Dlansy(byte(norm), uplo, n, a, lda, work)
}
// Dlantr computes the specified norm of an m×n trapezoidal matrix A. If
// norm == lapack.MaxColumnSum work must have length at least n, otherwise work
// is unused.
func (impl Implementation) Dlantr(norm lapack.MatrixNorm, uplo blas.Uplo, diag blas.Diag, m, n int, a []float64, lda int, work []float64) float64 {
checkMatrix(m, n, a, lda)
switch norm {
case lapack.MaxRowSum, lapack.MaxColumnSum, lapack.NormFrob, lapack.MaxAbs:
default:
panic(badNorm)
}
if uplo != blas.Upper && uplo != blas.Lower {
panic(badUplo)
}
if diag != blas.Unit && diag != blas.NonUnit {
panic(badDiag)
}
if norm == lapack.MaxColumnSum && len(work) < n {
panic(badWork)
}
return lapacke.Dlantr(byte(norm), uplo, diag, m, n, a, lda, work)
}
// Dlarfx applies an elementary reflector H to a real m×n matrix C, from either
// the left or the right, with loop unrolling when the reflector has order less
// than 11.
//
// H is represented in the form
// H = I - tau * v * v^T,
// where tau is a real scalar and v is a real vector. If tau = 0, then H is
// taken to be the identity matrix.
//
// v must have length equal to m if side == blas.Left, and equal to n if side ==
// blas.Right, otherwise Dlarfx will panic.
//
// c and ldc represent the m×n matrix C. On return, C is overwritten by the
// matrix H * C if side == blas.Left, or C * H if side == blas.Right.
//
// work must have length at least n if side == blas.Left, and at least m if side
// == blas.Right, otherwise Dlarfx will panic. work is not referenced if H has
// order < 11.
func (impl Implementation) Dlarfx(side blas.Side, m, n int, v []float64, tau float64, c []float64, ldc int, work []float64) {
checkMatrix(m, n, c, ldc)
switch side {
case blas.Left:
checkVector(m, v, 1)
if len(work) < n && m > 10 {
panic(badWork)
}
case blas.Right:
checkVector(n, v, 1)
if len(work) < m && n > 10 {
panic(badWork)
}
default:
panic(badSide)
}
lapacke.Dlarfx(side, m, n, v, tau, c, ldc, work)
}
// Dlascl multiplies an m×n matrix by the scalar cto/cfrom.
//
// cfrom must not be zero, and cto and cfrom must not be NaN, otherwise Dlascl
// will panic.
//
// Dlascl is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dlascl(kind lapack.MatrixType, kl, ku int, cfrom, cto float64, m, n int, a []float64, lda int) {
checkMatrix(m, n, a, lda)
if cfrom == 0 {
panic(zeroDiv)
}
if math.IsNaN(cfrom) || math.IsNaN(cto) {
panic(nanScale)
}
lapacke.Dlascl(byte(kind), kl, ku, cfrom, cto, m, n, a, lda)
}
// Dlaset sets the off-diagonal elements of A to alpha, and the diagonal
// elements to beta. If uplo == blas.Upper, only the elements in the upper
// triangular part are set. If uplo == blas.Lower, only the elements in the
// lower triangular part are set. If uplo is otherwise, all of the elements of A
// are set.
//
// Dlaset is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dlaset(uplo blas.Uplo, m, n int, alpha, beta float64, a []float64, lda int) {
checkMatrix(m, n, a, lda)
lapacke.Dlaset(uplo, m, n, alpha, beta, a, lda)
}
// Dlasrt sorts the numbers in the input slice d. If s == lapack.SortIncreasing,
// the elements are sorted in increasing order. If s == lapack.SortDecreasing,
// the elements are sorted in decreasing order. For other values of s Dlasrt
// will panic.
//
// Dlasrt is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dlasrt(s lapack.Sort, n int, d []float64) {
checkVector(n, d, 1)
switch s {
default:
panic(badSort)
case lapack.SortIncreasing, lapack.SortDecreasing:
}
lapacke.Dlasrt(byte(s), n, d[:n])
}
// Dlaswp swaps the rows k1 to k2 of a rectangular matrix A according to the
// indices in ipiv so that row k is swapped with ipiv[k].
//
// n is the number of columns of A and incX is the increment for ipiv. If incX
// is 1, the swaps are applied from k1 to k2. If incX is -1, the swaps are
// applied in reverse order from k2 to k1. For other values of incX Dlaswp will
// panic. ipiv must have length k2+1, otherwise Dlaswp will panic.
//
// The indices k1, k2, and the elements of ipiv are zero-based.
//
// Dlaswp is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dlaswp(n int, a []float64, lda, k1, k2 int, ipiv []int, incX int) {
switch {
case n < 0:
panic(nLT0)
case k2 < 0:
panic(badK2)
case k1 < 0 || k2 < k1:
panic(badK1)
case len(ipiv) != k2+1:
panic(badIpiv)
case incX != 1 && incX != -1:
panic(absIncNotOne)
}
ipiv32 := make([]int32, len(ipiv))
for i, v := range ipiv {
ipiv32[i] = int32(v + 1)
}
lapacke.Dlaswp(n, a, lda, k1+1, k2+1, ipiv32, incX)
}
// 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 lapacke.Dpotrf.
if n < 0 {
panic(nLT0)
}
if lda < n {
panic(badLdA)
}
if n == 0 {
return true
}
return lapacke.Dpotrf(ul, n, a, lda)
}
// Dgebal balances an n×n matrix A. Balancing consists of two stages, permuting
// and scaling. Both steps are optional and depend on the value of job.
//
// Permuting consists of applying a permutation matrix P such that the matrix
// that results from P^T*A*P takes the upper block triangular form
// [ T1 X Y ]
// P^T A P = [ 0 B Z ],
// [ 0 0 T2 ]
// where T1 and T2 are upper triangular matrices and B contains at least one
// nonzero off-diagonal element in each row and column. The indices ilo and ihi
// mark the starting and ending columns of the submatrix B. The eigenvalues of A
// isolated in the first 0 to ilo-1 and last ihi+1 to n-1 elements on the
// diagonal can be read off without any roundoff error.
//
// Scaling consists of applying a diagonal similarity transformation D such that
// D^{-1}*B*D has the 1-norm of each row and its corresponding column nearly
// equal. The output matrix is
// [ T1 X*D Y ]
// [ 0 inv(D)*B*D inv(D)*Z ].
// [ 0 0 T2 ]
// Scaling may reduce the 1-norm of the matrix, and improve the accuracy of
// the computed eigenvalues and/or eigenvectors.
//
// job specifies the operations that will be performed on A.
// If job is lapack.None, Dgebal sets scale[i] = 1 for all i and returns ilo=0, ihi=n-1.
// If job is lapack.Permute, only permuting will be done.
// If job is lapack.Scale, only scaling will be done.
// If job is lapack.PermuteScale, both permuting and scaling will be done.
//
// On return, if job is lapack.Permute or lapack.PermuteScale, it will hold that
// A[i,j] == 0, for i > j and j ∈ {0, ..., ilo-1, ihi+1, ..., n-1}.
// If job is lapack.None or lapack.Scale, or if n == 0, it will hold that
// ilo == 0 and ihi == n-1.
//
// On return, scale will contain information about the permutations and scaling
// factors applied to A. If π(j) denotes the index of the column interchanged
// with column j, and D[j,j] denotes the scaling factor applied to column j,
// then
// scale[j] == π(j), for j ∈ {0, ..., ilo-1, ihi+1, ..., n-1},
// == D[j,j], for j ∈ {ilo, ..., ihi}.
// scale must have length equal to n, otherwise Dgebal will panic.
//
// Dgebal is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dgebal(job lapack.Job, n int, a []float64, lda int, scale []float64) (ilo, ihi int) {
switch job {
default:
panic(badJob)
case lapack.None, lapack.Permute, lapack.Scale, lapack.PermuteScale:
}
checkMatrix(n, n, a, lda)
if len(scale) != n {
panic("lapack: bad length of scale")
}
ilo32 := make([]int32, 1)
ihi32 := make([]int32, 1)
lapacke.Dgebal(job, n, a, lda, ilo32, ihi32, scale)
ilo = int(ilo32[0]) - 1
ihi = int(ihi32[0]) - 1
for j := 0; j < ilo; j++ {
scale[j]--
}
for j := ihi + 1; j < n; j++ {
scale[j]--
}
return ilo, ihi
}
// Dgebak transforms an n×m matrix V as
// V = P D V, if side == blas.Right,
// V = P D^{-1} V, if side == blas.Left,
// where P and D are n×n permutation and scaling matrices, respectively,
// implicitly represented by job, scale, ilo and ihi as returned by Dgebal.
//
// Typically, columns of the matrix V contain the right or left (determined by
// side) eigenvectors of the balanced matrix output by Dgebal, and Dgebak forms
// the eigenvectors of the original matrix.
//
// Dgebak is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dgebak(job lapack.Job, side lapack.EVSide, n, ilo, ihi int, scale []float64, m int, v []float64, ldv int) {
switch job {
default:
panic(badJob)
case lapack.None, lapack.Permute, lapack.Scale, lapack.PermuteScale:
}
var bside blas.Side
switch side {
default:
panic(badSide)
case lapack.LeftEV:
bside = blas.Left
case lapack.RightEV:
bside = blas.Right
}
checkMatrix(n, m, v, ldv)
switch {
case ilo < 0 || max(0, n-1) < ilo:
panic(badIlo)
case ihi < min(ilo, n-1) || n <= ihi:
panic(badIhi)
}
// Convert permutation indices to 1-based.
for j := 0; j < ilo; j++ {
scale[j]++
}
for j := ihi + 1; j < n; j++ {
scale[j]++
}
lapacke.Dgebak(job, bside, n, ilo+1, ihi+1, scale, m, v, ldv)
// Convert permutation indices back to 0-based.
for j := 0; j < ilo; j++ {
scale[j]--
}
for j := ihi + 1; j < n; j++ {
scale[j]--
}
}
// Dbdsqr performs a singular value decomposition of a real n×n bidiagonal matrix.
//
// The SVD of the bidiagonal matrix B is
// B = Q * S * P^T
// where S is a diagonal matrix of singular values, Q is an orthogonal matrix of
// left singular vectors, and P is an orthogonal matrix of right singular vectors.
//
// Q and P are only computed if requested. If left singular vectors are requested,
// this routine returns U * Q instead of Q, and if right singular vectors are
// requested P^T * VT is returned instead of P^T.
//
// Frequently Dbdsqr is used in conjunction with Dgebrd which reduces a general
// matrix A into bidiagonal form. In this case, the SVD of A is
// A = (U * Q) * S * (P^T * VT)
//
// This routine may also compute Q^T * C.
//
// d and e contain the elements of the bidiagonal matrix b. d must have length at
// least n, and e must have length at least n-1. Dbdsqr will panic if there is
// insufficient length. On exit, D contains the singular values of B in decreasing
// order.
//
// VT is a matrix of size n×ncvt whose elements are stored in vt. The elements
// of vt are modified to contain P^T * VT on exit. VT is not used if ncvt == 0.
//
// U is a matrix of size nru×n whose elements are stored in u. The elements
// of u are modified to contain U * Q on exit. U is not used if nru == 0.
//
// C is a matrix of size n×ncc whose elements are stored in c. The elements
// of c are modified to contain Q^T * C on exit. C is not used if ncc == 0.
//
// work contains temporary storage and must have length at least 4*n. Dbdsqr
// will panic if there is insufficient working memory.
//
// Dbdsqr returns whether the decomposition was successful.
func (impl Implementation) Dbdsqr(uplo blas.Uplo, n, ncvt, nru, ncc int, d, e, vt []float64, ldvt int, u []float64, ldu int, c []float64, ldc int, work []float64) (ok bool) {
if uplo != blas.Upper && uplo != blas.Lower {
panic(badUplo)
}
if ncvt != 0 {
checkMatrix(n, ncvt, vt, ldvt)
}
if nru != 0 {
checkMatrix(nru, n, u, ldu)
}
if ncc != 0 {
checkMatrix(n, ncc, c, ldc)
}
if len(d) < n {
panic(badD)
}
if len(e) < n-1 {
panic(badE)
}
if len(work) < 4*n {
panic(badWork)
}
// An address must be passed to cgo. If lengths are zero, allocate a slice.
if len(vt) == 0 {
vt = make([]float64, 1)
}
if len(u) == 0 {
vt = make([]float64, 1)
}
if len(c) == 0 {
c = make([]float64, 1)
}
return lapacke.Dbdsqr(uplo, n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c, ldc, work)
}
// Dgebrd reduces a general m×n matrix A to upper or lower bidiagonal form B by
// an orthogonal transformation:
// Q^T * A * P = B.
// The diagonal elements of B are stored in d and the off-diagonal elements are
// stored in e. These are additionally stored along the diagonal of A and the
// off-diagonal of A. If m >= n B is an upper-bidiagonal matrix, and if m < n B
// is a lower-bidiagonal matrix.
//
// The remaining elements of A store the data needed to construct Q and P.
// The matrices Q and P are products of elementary reflectors
// if m >= n, Q = H_0 * H_1 * ... * H_{n-1},
// P = G_0 * G_1 * ... * G_{n-2},
// if m < n, Q = H_0 * H_1 * ... * H_{m-2},
// P = G_0 * G_1 * ... * G_{m-1},
// where
// H_i = I - tauQ[i] * v_i * v_i^T,
// G_i = I - tauP[i] * u_i * u_i^T.
//
// As an example, on exit the entries of A when m = 6, and n = 5
// [ d e u1 u1 u1]
// [v1 d e u2 u2]
// [v1 v2 d e u3]
// [v1 v2 v3 d e]
// [v1 v2 v3 v4 d]
// [v1 v2 v3 v4 v5]
// and when m = 5, n = 6
// [ d u1 u1 u1 u1 u1]
// [ e d u2 u2 u2 u2]
// [v1 e d u3 u3 u3]
// [v1 v2 e d u4 u4]
// [v1 v2 v3 e d u5]
//
// d, tauQ, and tauP must all have length at least min(m,n), and e must have
// length min(m,n) - 1, unless lwork is -1 when there is no check except for
// work which must have a length of at least one.
//
// work is temporary storage, and lwork specifies the usable memory length.
// At minimum, lwork >= max(1,m,n) or be -1 and this function will panic otherwise.
// Dgebrd is blocked decomposition, but the block size is limited
// by the temporary space available. If lwork == -1, instead of performing Dgebrd,
// the optimal work length will be stored into work[0].
//
// Dgebrd is an internal routine. It is exported for testing purposes.
func (impl Implementation) Dgebrd(m, n int, a []float64, lda int, d, e, tauQ, tauP, work []float64, lwork int) {
checkMatrix(m, n, a, lda)
minmn := min(m, n)
if len(d) < minmn {
panic(badD)
}
if len(e) < minmn-1 {
panic(badE)
}
if len(tauQ) < minmn {
panic(badTauQ)
}
if len(tauP) < minmn {
panic(badTauP)
}
if lwork != -1 && lwork < max(1, max(m, n)) {
panic(badWork)
}
if len(work) < max(1, lwork) {
panic(badWork)
}
lapacke.Dgebrd(m, n, a, lda, d, e, tauQ, tauP, work, lwork)
}
// Dgecon estimates the reciprocal of the condition number of the n×n matrix A
// given the LU decomposition of the matrix. The condition number computed may
// be based on the 1-norm or the ∞-norm.
//
// The slice a contains the result of the LU decomposition of A as computed by Dgetrf.
//
// anorm is the corresponding 1-norm or ∞-norm of the original matrix A.
//
// work is a temporary data slice of length at least 4*n and Dgecon will panic otherwise.
//
// iwork is a temporary data slice of length at least n and Dgecon will panic otherwise.
// Elements of iwork must fit within the int32 type or Dgecon will panic.
func (impl Implementation) Dgecon(norm lapack.MatrixNorm, n int, a []float64, lda int, anorm float64, work []float64, iwork []int) float64 {
checkMatrix(n, n, a, lda)
if norm != lapack.MaxColumnSum && norm != lapack.MaxRowSum {
panic("bad norm")
}
if len(work) < 4*n {
panic(badWork)
}
if len(iwork) < n {
panic(badWork)
}
rcond := make([]float64, 1)
_iwork := make([]int32, len(iwork))
for i, v := range iwork {
if v != int(int32(v)) {
panic("lapack: iwork element out of range")
}
_iwork[i] = int32(v)
}
lapacke.Dgecon(byte(norm), n, a, lda, anorm, rcond, work, _iwork)
for i, v := range _iwork {
iwork[i] = int(v)
}
return rcond[0]
}
// 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
// orthonormal 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-1} * ... * H_1 * H_0,
// where k = min(m,n).
//
// 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)
if len(tau) < min(m, n) {
panic(badTau)
}
if len(work) < m {
panic(badWork)
}
lapacke.Dgelq2(m, n, a, lda, tau, work)
}