1 /++
2 $(H1 Lubeck - Linear Algebra)
3 
4 Authors: Ilya Yaroshenko, Lars Tandle Kyllingstad (SciD author)
5 +/
6 module kaleidic.lubeck;
7 
8 import cblas : Diag;
9 import mir.blas;
10 import mir.internal.utility : realType, isComplex;
11 import mir.lapack;
12 import mir.math.common;
13 import mir.ndslice.allocation;
14 import mir.ndslice.dynamic: transposed;
15 import mir.ndslice.slice;
16 import mir.ndslice.topology;
17 import mir.ndslice.traits : isMatrix;
18 import mir.utility;
19 import std.meta;
20 import std.traits;
21 import std.typecons: Flag, Yes, No;
22 import mir.complex: Complex;
23 public import mir.lapack: lapackint;
24 
25 template CommonType(A)
26 {
27     alias CommonType = A;
28 }
29 
30 template CommonType(A, B)
31 {
32     static if (isComplex!A || isComplex!B)
33         alias CommonType = Complex!(CommonType!(realType!A, realType!B));
34     else
35         alias CommonType = typeof(A.init + B.init);
36 }
37 
38 version(LDC)
39     import ldc.attributes: fastmath;
40 else
41     enum { fastmath };
42 
43 private template IterationType(Iterator)
44 {
45     alias T = Unqual!(typeof(Iterator.init[0]));
46     static if (isIntegral!T || is(T == real))
47         alias IterationType = double;
48     else
49     static if (is(T == Complex!real))
50         alias IterationType = Complex!double;
51     else
52     {
53         static assert(
54             is(T == double) ||
55             is(T == float) ||
56             is(T == Complex!double) ||
57             is(T == Complex!float));
58         alias IterationType = T;
59     }
60 }
61 
62 /++
63 Gets the type that can be used with Blas routines that all types can be implicitly converted to. 
64 
65 +/
66 alias BlasType(Iterators...) =
67     CommonType!(staticMap!(IterationType, Iterators));
68 
69 /++
70 General matrix-matrix multiplication. Allocates result to an uninitialized slice using GC.
71 Params:
72     a = m(rows) x k(cols) matrix
73     b = k(rows) x n(cols) matrix
74 Result: 
75     m(rows) x n(cols)
76 +/
77 Slice!(BlasType!(IteratorA, IteratorB)*, 2)
78     mtimes(IteratorA, SliceKind kindA, IteratorB, SliceKind kindB)(
79         Slice!(IteratorA, 2, kindA) a,
80         Slice!(IteratorB, 2, kindB) b)
81 {
82     assert(a.length!1 == b.length!0);
83 
84     // reallocate data if required
85     alias A = BlasType!IteratorA;
86     alias B = BlasType!IteratorB;
87     alias C = CommonType!(A, B);
88     static if (!is(Unqual!IteratorA == C*))
89         return .mtimes(a.as!C.slice, b);
90     else
91     static if (!is(Unqual!IteratorB == C*))
92         return .mtimes(a, b.as!C.slice);
93     else
94     {
95         static if (kindA != Contiguous)
96             if (a._stride!0 != 1 && a._stride!1 != 1
97                 || a._stride!0 <= 0
98                 || a._stride!1 <= 0)
99                 return .mtimes(a.slice, b);
100 
101         static if (kindB != Contiguous)
102             if (b._stride!0 != 1 && b._stride!1 != 1
103                 || b._stride!0 <= 0
104                 || b._stride!1 <= 0)
105                 return .mtimes(a, b.slice);
106 
107         auto c = uninitSlice!C(a.length!0, b.length!1);
108 
109         if (a.length!1 == 1 && b.length!0 == 1)
110         {
111             c[] = cast(C) 0;
112             ger(cast(C)1, a.front!1, b.front, c);
113         }
114         else
115         {
116             gemm(cast(C)1, a, b, cast(C)0, c);
117         }
118         return c;
119     }
120 }
121 
122 ///
123 unittest
124 {
125     import mir.ndslice;
126 
127     auto a =
128         [-5,  1,  7, 7, -4,
129          -1, -5,  6, 3, -3,
130          -5, -2, -3, 6,  0].sliced(3, 5);
131 
132     auto b = slice!double(5, 4);
133     b[] =
134         [[-5, -3,  3,  1],
135          [ 4,  3,  6,  4],
136          [-4, -2, -2,  2],
137          [-1,  9,  4,  8],
138          [ 9, 8,  3, -2]];
139 
140     assert(mtimes(a, b) ==
141         [[-42,  35,  -7, 77],
142          [-69, -21, -42, 21],
143          [ 23,  69,   3, 29]]
144         );
145 }
146 
147 /// ger specialized case in mtimes
148 unittest
149 {
150     import mir.ndslice;
151 
152     // from https://github.com/kaleidicassociates/lubeck/issues/8
153     {
154         auto a = [1.0f, 2.0f].sliced(2, 1);
155         auto b = [1.0f, 2.0f].sliced(2, 1);
156         assert(mtimes(a, b.transposed) == [[1, 2], [2, 4]]);
157     }
158     {
159         auto a = [1.0, 2.0].sliced(1, 2);
160         auto b = [1.0, 2.0].sliced(1, 2);
161         assert(mtimes(a.transposed, b) == [[1, 2], [2, 4]]);
162     }
163 }
164 
165 ///
166 unittest
167 {
168     import mir.ndslice;
169 
170     // from https://github.com/kaleidicassociates/lubeck/issues/3
171     Slice!(float*, 2) a = slice!float(1, 1);
172     Slice!(float*, 2, Universal) b1 = slice!float(16, 1).transposed;
173     Slice!(float*, 2) b2 = slice!float(1, 16);
174 
175     a[] = 3;
176     b1[] = 4;
177     b2[] = 4;
178 
179     // Confirm that this message does not appear
180     // Outputs: ** On entry to SGEMM  parameter number  8 had an illegal value
181     assert(a.mtimes(b1) == a.mtimes(b2));
182 }
183 
184 /++
185 General matrix-matrix multiplication. Allocates result to an uninitialized slice using GC.
186 Params:
187     a = m(rows) x k(cols) matrix
188     b = k(rows) x 1(cols) vector
189 Result:
190     m(rows) x 1(cols)
191 +/
192 Slice!(BlasType!(IteratorA, IteratorB)*)
193     mtimes(IteratorA, SliceKind kindA, IteratorB, SliceKind kindB)(
194         Slice!(IteratorA, 2, kindA) a,
195         Slice!(IteratorB, 1, kindB) b)
196 {
197     assert(a.length!1 == b.length!0);
198 
199     // reallocate data if required
200     alias A = BlasType!IteratorA;
201     alias B = BlasType!IteratorB;
202     alias C = CommonType!(A, B);
203     static if (!is(Unqual!IteratorA == C*))
204         return .mtimes(a.as!C.slice, b);
205     else
206     static if (!is(Unqual!IteratorB == C*))
207         return .mtimes(a, b.as!C.slice);
208     else
209     {
210         static if (kindA != Contiguous)
211             if (a._stride!0 != 1 && a._stride!1 != 1
212                 || a._stride!0 <= 0
213                 || a._stride!1 <= 0)
214                 return .mtimes(a.slice, b);
215 
216         static if (kindB != Contiguous)
217             if (b._stride!1 <= 0)
218                 return .mtimes(a, b.slice);
219 
220         auto c = uninitSlice!C(a.length!0);
221         gemv(cast(C)1, a, b, cast(C)0, c);
222         return c;
223     }
224 }
225 
226 /++
227 General matrix-matrix multiplication.
228 Params:
229     a = 1(rows) x k(cols) vector
230     b = k(rows) x n(cols) matrix
231 Result:
232     1(rows) x n(cols)
233 +/
234 Slice!(BlasType!(IteratorA, IteratorB)*)
235     mtimes(IteratorA, SliceKind kindA, IteratorB, SliceKind kindB)(
236         Slice!(IteratorB, 1, kindB) a,
237         Slice!(IteratorA, 2, kindA) b,
238         )
239 {
240     return .mtimes(b.universal.transposed, a);
241 }
242 
243 ///
244 unittest
245 {
246     import mir.ndslice;
247 
248     auto a =
249         [-5,  1,  7, 7, -4,
250          -1, -5,  6, 3, -3,
251          -5, -2, -3, 6,  0]
252             .sliced(3, 5)
253             .universal
254             .transposed;
255 
256     auto b = slice!double(5);
257     b[] = [-5, 4,-4,-1, 9];
258 
259     assert(mtimes(b, a) == [-42, -69, 23]);
260 }
261 
262 /++
263 Vector-vector multiplication (dot product).
264 Params:
265     a = 1(rows) x k(cols) vector
266     b = k(rows) x 1(cols) matrix
267 Result:
268     scalar
269 +/
270 CommonType!(BlasType!IteratorA, BlasType!IteratorB)
271     mtimes(IteratorA, SliceKind kindA, IteratorB, SliceKind kindB)(
272         Slice!(IteratorB, 1, kindB) a,
273         Slice!(IteratorA, 1, kindA) b,
274         )
275 {
276     alias A = BlasType!IteratorA;
277     alias B = BlasType!IteratorB;
278     alias C = CommonType!(A, B);
279     static if (is(IteratorB == C*) && is(IteratorA == C*))
280     {
281         return dot(a, b);
282     }
283     else
284     {
285         auto c = cast(typeof(return)) 0;
286         import mir.algorithm.iteration: reduce;
287         return c.reduce!"a + b * c"(a.as!(typeof(return)), b.as!(typeof(return)));
288     }
289 }
290 
291 ///
292 unittest
293 {
294     import mir.ndslice;
295 
296     auto a = [1, 2, 4].sliced;
297     auto b = [3, 4, 2].sliced;
298     assert(a.mtimes(b) == 19);
299 }
300 
301 /++
302 Calculates the inverse of a matrix.
303 +/
304 auto inv(Iterator, SliceKind kind)(Slice!(Iterator, 2, kind) a)
305 in
306 {
307     assert (a.length!0 == a.length!1, "matrix must be square");
308 }
309 do
310 {
311     alias T = BlasType!Iterator;
312 
313     auto m = a.as!T.slice.canonical;
314     auto ipiv = m.length.uninitSlice!lapackint;
315 
316     auto info = getrf!T(m, ipiv);
317     if (info == 0)
318     {
319         info = getri(m, ipiv, m.getri_wq.uninitSlice!T);
320     }
321 
322     import std.exception: enforce;
323     enforce(info == 0, "inv: matrix is singular");
324     return m;
325 }
326 
327 ///
328 unittest
329 {
330     import mir.complex;
331     import mir.ndslice;
332 
333     auto a =  [
334         1, 0, 2,
335         2, 2, 0,
336         0, 1, 1]
337         .sliced(3, 3);
338     
339     enum : double { _13 = 1.0/3.0, _16 = 1.0/6.0, _23 = 2.0/3.0 }
340     auto ans = [
341         _13, _13, -_23,
342         -_13,_16, _23,
343         _13, -_16, _13]
344         .sliced(a.shape);
345 
346     import mir.algorithm.iteration: equal;
347     import mir.math.common: approxEqual;
348     assert(equal!((a, b) => a.approxEqual(b, 1e-10L, 1e-10L))(a.inv, ans));
349     assert(equal!((a, b) => a.approxEqual(b, 1e-10L, 1e-10L))(a.map!(a => Complex!double(a, 0)).inv.member!"re", ans));
350 }
351 
352 ///
353 unittest
354 {
355     import mir.ndslice.topology: iota;
356 
357     try
358     {
359         auto m = [3, 3].iota.inv;
360         assert (false, "Matrix should be detected as singular");
361     }
362     catch (Exception e)
363     {
364         assert (true);
365     }
366 }
367 
368 ///
369 struct SvdResult(T)
370 {
371     ///
372     Slice!(T*, 2) u;
373     ///
374     Slice!(T*) sigma;
375     ///
376     Slice!(T*, 2) vt;
377 }
378 
379 /++
380 Computes the singular value decomposition.
381 
382 Params:
383     matrix = input `M x N` matrix
384     slim = If true the first `min(M,N)` columns of `u` and the first
385         `min(M,N)` rows of `vt` are returned in the ndslices `u` and `vt`.
386 Returns: $(LREF SvdResult). Results are allocated by the GC.
387 +/
388 auto svd(
389         Flag!"allowDestroy" allowDestroy = No.allowDestroy,
390         string algorithm = "gesvd",
391         SliceKind kind,
392         Iterator
393     )(
394         Slice!(Iterator, 2, kind) matrix,
395         Flag!"slim" slim = No.slim,
396     )
397     if (algorithm == "gesvd" || algorithm == "gesdd")
398 {
399     import lapack;
400     alias T = BlasType!Iterator;
401     static if (allowDestroy && kind != Universal && is(Iterstor == T*))
402         alias a = matrix.canonical;
403     else
404         auto a = matrix.as!T.slice.canonical;
405 
406     auto m = cast(lapackint)a.length!1;
407     auto n = cast(lapackint)a.length!0;
408 
409     auto s = uninitSlice!T(min(m, n));
410     auto u = uninitSlice!T(slim ? s.length : m, m);
411     auto vt = uninitSlice!T(n, slim ? s.length : n);
412 
413     if (m == 0 || n == 0)
414     {
415         u[] = 0;
416         u.diagonal[] = 1;
417         vt[] = 0;
418         vt.diagonal[] = 1;
419     }
420     else
421     {
422         static if (algorithm == "gesvd")
423         {
424             auto jobu = slim ? 'S' : 'A';
425             auto jobvt = slim ? 'S' : 'A';
426             auto work = gesvd_wq(jobu, jobvt, a, u.canonical, vt.canonical).uninitSlice!T;
427             auto info = gesvd(jobu, jobvt, a, s, u.canonical, vt.canonical, work);
428             enum msg = "svd: DBDSQR did not converge";
429         }
430         else // gesdd
431         {
432             auto iwork = uninitSlice!lapackint(s.length * 8);
433             auto jobz = slim ? 'S' : 'A';
434             auto work = gesdd_wq(jobz, a, u.canonical, vt.canonical).uninitSlice!T;
435             auto info = gesdd(jobz, a, s, u.canonical, vt.canonical, work, iwork);
436             enum msg = "svd: DBDSDC did not converge, updating process failed";
437         }
438         import std.exception: enforce;
439         enforce(info == 0, msg);
440     }
441     return SvdResult!T(vt, s, u); //transposed
442 }
443 
444 ///
445 unittest
446 {
447     import mir.ndslice;
448 
449     auto a =  [
450          7.52,  -1.10,  -7.95,   1.08,
451         -0.76,   0.62,   9.34,  -7.10,
452          5.13,   6.62,  -5.66,   0.87,
453         -4.75,   8.52,   5.75,   5.30,
454          1.33,   4.91,  -5.49,  -3.52,
455         -2.40,  -6.77,   2.34,   3.95]
456         .sliced(6, 4);
457 
458     auto r = a.svd;
459 
460     auto sigma = slice!double(a.shape, 0);
461     sigma.diagonal[] = r.sigma;
462     auto m = r.u.mtimes(sigma).mtimes(r.vt);
463 
464     import mir.algorithm.iteration: equal;
465     import mir.math.common: approxEqual;
466     assert(equal!((a, b) => a.approxEqual(b, 1e-8, 1e-8))(a, m));
467 }
468 
469 ///
470 unittest
471 {
472     import std.typecons: Yes;
473     import mir.ndslice;
474 
475     auto a =  [
476          7.52,  -1.10,  -7.95,   1.08,
477         -0.76,   0.62,   9.34,  -7.10,
478          5.13,   6.62,  -5.66,   0.87,
479         -4.75,   8.52,   5.75,   5.30,
480          1.33,   4.91,  -5.49,  -3.52,
481         -2.40,  -6.77,   2.34,   3.95]
482         .sliced(6, 4);
483 
484     auto r = a.svd(Yes.slim);
485     assert(r.u.shape == [6, 4]);
486     assert(r.vt.shape == [4, 4]);
487 }
488 
489 unittest
490 {
491     import mir.algorithm.iteration: all;
492 
493     // empty matrix as input means that u or vt is identity matrix
494     auto identity = slice!double([4, 4], 0);
495     identity.diagonal[] = 1;
496 
497     auto a = slice!double(0, 4);
498     auto res = a.svd;
499 
500     import mir.conv: to;
501 
502     assert(res.u.shape == [0, 0]);
503     assert(res.vt.shape == [4, 4]);
504     assert(res.vt.all!approxEqual(identity), res.vt.to!string);
505 
506     auto b = slice!double(4, 0);
507     res = b.svd;
508 
509     assert(res.u.shape == [4, 4]);
510     assert(res.vt.shape == [0, 0]);
511 
512     assert(res.u.all!approxEqual(identity), res.u.to!string);
513 }
514 
515 /++
516 Solve systems of linear equations AX = B for X.
517 Computes minimum-norm solution to a linear least squares problem
518 if A is not a square matrix.
519 +/
520 Slice!(BlasType!(IteratorA, IteratorB)*, 2)
521     mldivide
522     (IteratorA, SliceKind kindA, IteratorB, SliceKind kindB)(
523         Slice!(IteratorA, 2, kindA) a,
524         Slice!(IteratorB, 2, kindB) b)
525 {
526     import std.exception: enforce;
527     import std.conv: to;
528 
529     assert(a.length!0 == b.length!0);
530 
531     alias A = BlasType!IteratorA;
532     alias B = BlasType!IteratorB;
533     alias C = CommonType!(A, B);
534 
535     auto a_ = a.universal.transposed.as!C.slice.canonical;
536     auto b_ = b.universal.transposed.as!C.slice.canonical;
537 
538     if (a.length!0 == a.length!1)
539     {
540         auto ipiv = a_.length.uninitSlice!lapackint;
541         auto info = gesv(a_, ipiv, b_);
542     }
543     else
544     {
545         static if(!isComplex!C)
546         {
547             size_t liwork = void;
548             auto lwork = gelsd_wq(a_, b_, liwork);
549             auto s = min(a_.length!0, a_.length!1).uninitSlice!C;
550             auto work = lwork.uninitSlice!C;
551             auto iwork = liwork.uninitSlice!lapackint;
552             size_t rank;
553             C rcond = -1;
554 
555             auto info = gelsd(a_, b_, s, rcond, rank, work, iwork);
556         }
557         else
558         {
559             size_t liwork = void;
560             size_t lrwork = void;
561             auto lwork = gelsd_wq(a_, b_, lrwork, liwork);
562             auto s = min(a_.length!0, a_.length!1).uninitSlice!(realType!C);
563             auto work = lwork.uninitSlice!C;
564             auto iwork = liwork.uninitSlice!lapackint;
565             auto rwork = lrwork.uninitSlice!(realType!C);
566             size_t rank;
567             realType!C rcond = -1;
568 
569             auto info = gelsd!C(a_, b_, s, rcond, rank, work, rwork, iwork);
570         }
571 
572         enforce(info == 0, to!string(info) ~ " off-diagonal elements of an intermediate bidiagonal form did not converge to zero.");
573         b_ = b_[0 .. $, 0 .. a_.length!0];
574     }
575 
576     return b_.universal.transposed.slice;
577 }
578 
579 /// ditto
580 Slice!(BlasType!(IteratorA, IteratorB)*)
581     mldivide
582     (IteratorA, SliceKind kindA, IteratorB, SliceKind kindB)(
583         Slice!(IteratorA, 2, kindA) a,
584         Slice!(IteratorB, 1, kindB) b)
585 {
586     return a.mldivide(b.repeat(1).unpack.universal.transposed).front!1.assumeContiguous;
587 }
588 
589 /// AX=B
590 unittest
591 {
592     import mir.complex;
593     import std.meta: AliasSeq;
594     import mir.ndslice;
595 
596     foreach(C; AliasSeq!(double, Complex!double))
597     {
598         static if(is(C == Complex!double))
599             alias transform = a => C(a, 0);
600         else
601             enum transform = "a";
602 
603         auto a = [
604             1, -1,  1,
605             2,  2, -4,
606             -1,  5,  0].sliced(3, 3).map!transform;
607         auto b = [
608             2.0,  0,
609             -6  , -6,
610             9  ,  1].sliced(3, 2).map!transform;
611         auto t = [
612             1.0, -1,
613             2  ,  0,
614             3  ,  1].sliced(3, 2).map!transform;
615 
616         auto x = mldivide(a, b);
617         assert(x == t);
618     }
619 }
620 
621 /// Ax=B
622 unittest
623 {
624     import mir.complex;
625     import std.meta: AliasSeq;
626     import mir.ndslice;
627 
628     foreach(C; AliasSeq!(double, Complex!double))
629     {
630         static if(is(C == Complex!double))
631             alias transform = a => C(a, 0);
632         else
633             enum transform = "a";
634 
635         auto a = [
636             1, -1,  1,
637             2,  2, -4,
638             -1,  5,  0].sliced(3, 3).map!transform;
639         auto b = [
640             2.0,
641             -6  ,
642             9  ].sliced(3).map!transform;
643         auto t = [
644             1.0,
645             2  ,
646             3  ].sliced(3).map!transform;
647 
648         auto x = mldivide(a, b);
649         assert(x == t);
650     }
651 }
652 
653 /// Least-Squares Solution of Underdetermined System
654 unittest
655 {
656     import mir.complex;
657     import std.meta: AliasSeq;
658     import mir.ndslice;
659 
660     foreach(C; AliasSeq!(double, )) //Complex!double fails for DMD>=2085
661     {
662         static if(is(C == Complex!double))
663             alias transform = a => C(a, 0);
664         else
665             enum transform = "a";
666 
667         auto a = [
668             -0.57,  -1.28,  -0.39,   0.25,
669             -1.93,   1.08,  -0.31,  -2.14,
670             2.30,   0.24,   0.40,  -0.35,
671             -1.93,   0.64,  -0.66,   0.08,
672             0.15,   0.30,   0.15,  -2.13,
673             -0.02,   1.03,  -1.43,   0.50,
674             ].sliced(6, 4).map!transform;
675 
676         auto b = [
677         -2.67,
678         -0.55,
679         3.34,
680         -0.77,
681         0.48,
682         4.10,
683         ].sliced.map!transform;
684 
685         auto x = [
686             1.5339,
687             1.8707,
688             -1.5241,
689             0.0392].sliced.map!transform;
690 
691         import mir.math.common: approxEqual;
692         import mir.algorithm.iteration: all;
693         alias appr = all!((a, b) => approxEqual(a, b, 1e-3, 1e-3));
694         assert(appr(a.mldivide(b), x));
695     }
696 }
697 
698 /// Principal component analises result.
699 struct PcaResult(T)
700 {
701     /// Principal component coefficients, also known as loadings.
702     Slice!(T*, 2) coeff;
703     /// Principal component scores.
704     Slice!(T*, 2) score;
705     /// Principal component variances.
706     Slice!(T*) latent;
707 }
708 
709 /++
710 Principal component analysis of raw data.
711 
712 Params:
713     matrix = input `M x N` matrix, where 'M (rows)>= N(cols)'
714     cc = Flag to centern columns. True by default.
715 Returns: $(LREF PcaResult)
716 +/
717 auto pca(Flag!"allowDestroy" allowDestroy = No.allowDestroy, Iterator, SliceKind kind)(Slice!(Iterator, 2, kind) matrix, in Flag!"centerColumns" cc = Yes.centerColumns)
718 in
719 {
720     assert(matrix.length!0 >= matrix.length!1);
721 }
722 do
723 {
724     import mir.math.sum: sum;
725     import mir.algorithm.iteration: maxIndex, eachUploPair;
726     import mir.utility: swap;
727 
728     alias T = BlasType!Iterator;
729     SvdResult!T svdResult;
730     if (cc)
731     {
732         static if (allowDestroy && kind != Universal && is(Iterstor == T*))
733             alias m = matrix;
734         else
735             auto m = matrix.as!T.slice;
736         foreach (col; m.universal.transposed)
737             col[] -= col.sum!"kb2" / col.length;
738         svdResult = m.svd!(Yes.allowDestroy)(Yes.slim);
739     }
740     else
741     {
742         svdResult = matrix.svd!(allowDestroy)(Yes.slim);
743     }
744     with (svdResult)
745     {
746         foreach (row; u)
747             row[] *= sigma;
748         T c = max(0, ptrdiff_t(matrix.length) - cast(bool) cc);
749         foreach (ref s; sigma)
750             s = s * s / c;
751         foreach(size_t i; 0 .. sigma.length)
752         {
753             auto col = vt[i];
754             if (col[col.map!fabs.maxIndex] < 0)
755             {
756                 col[] *= -1;
757                 u[0 .. $, i] *= -1;
758             }
759         }
760         vt.eachUploPair!swap;
761         return PcaResult!T(vt, u, sigma);
762     }
763 }
764 
765 ///
766 unittest
767 {
768     import mir.ndslice;
769 
770     import mir.math.common: approxEqual;
771     import mir.algorithm.iteration: equal;
772 
773     auto ingedients = [
774          7,  26,   6,  60,
775          1,  29,  15,  52,
776         11,  56,   8,  20,
777         11,  31,   8,  47,
778          7,  52,   6,  33,
779         11,  55,   9,  22,
780          3,  71,  17,   6,
781          1,  31,  22,  44,
782          2,  54,  18,  22,
783         21,  47,   4,  26,
784          1,  40,  23,  34,
785         11,  66,   9,  12,
786         10,  68,   8,  12].sliced(13, 4);
787 
788     auto res = ingedients.pca;
789 
790     auto coeff = [
791         -0.067799985695474,  -0.646018286568728,   0.567314540990512,   0.506179559977705,
792         -0.678516235418647,  -0.019993340484099,  -0.543969276583817,   0.493268092159297,
793          0.029020832106229,   0.755309622491133,   0.403553469172668,   0.515567418476836,
794          0.730873909451461,  -0.108480477171676,  -0.468397518388289,   0.484416225289198,
795     ].sliced(4, 4);
796 
797     auto score = [
798          36.821825999449700,  -6.870878154227367,  -4.590944457629745,   0.396652582713912,
799          29.607273420710964,   4.610881963526308,  -2.247578163663940,  -0.395843536696492,
800         -12.981775719737618,  -4.204913183175938,   0.902243082694698,  -1.126100587210615,
801          23.714725720918022,  -6.634052554708721,   1.854742000806314,  -0.378564808384691,
802          -0.553191676624597,  -4.461732123178686,  -6.087412652325177,   0.142384896047281,
803         -10.812490833309816,  -3.646571174544059,   0.912970791674604,  -0.134968810314680,
804         -32.588166608817929,   8.979846284936063,  -1.606265913996588,   0.081763927599947,
805          22.606395499005586,  10.725906457369449,   3.236537714483416,   0.324334774646368,
806          -9.262587237675838,   8.985373347478788,  -0.016909578102172,  -0.543746175981799,
807          -3.283969329640680, -14.157277337500918,   7.046512994833761,   0.340509860960606,
808           9.220031117829379,  12.386080787220454,   3.428342878284624,   0.435152769664895,
809         -25.584908517429557,  -2.781693148152386,  -0.386716066864491,   0.446817950545605,
810         -26.903161834677597,  -2.930971165042989,  -2.445522630195304,   0.411607156409658,
811     ].sliced(13, 4);
812 
813     auto latent = [5.177968780739053, 0.674964360487231, 0.124054300480810, 0.002371532651878].sliced;
814     latent[] *= 100;
815 
816     assert(equal!approxEqual(res.coeff, coeff));
817     assert(equal!approxEqual(res.score, score));
818     assert(equal!approxEqual(res.latent, latent));
819 }
820 
821 /++
822 Computes Moore-Penrose pseudoinverse of matrix.
823 
824 Params:
825     matrix = Input `M x N` matrix.
826     tolerance = The computation is based on SVD and any singular values less than tolerance are treated as zero.
827 Returns: Moore-Penrose pseudoinverse matrix
828 +/
829 Slice!(BlasType!Iterator*, 2)
830     pinv(Flag!"allowDestroy" allowDestroy = No.allowDestroy, Iterator, SliceKind kind)(Slice!(Iterator, 2, kind) matrix, double tolerance = double.nan)
831 {
832     import mir.algorithm.iteration: find, each;
833     import std.math: nextUp;
834 
835     auto svd = matrix.svd!allowDestroy(Yes.slim);
836     if (tolerance != tolerance)
837     {
838         auto n = svd.sigma.front;
839         auto eps = n.nextUp - n;
840         tolerance = max(matrix.length!0, matrix.length!1) * eps;
841     }
842     auto st = svd.sigma.find!(a => !(a >= tolerance));
843     static if (is(typeof(st) : sizediff_t))
844         alias si = st;
845     else
846         auto si = st[0];
847     auto s = svd.sigma[0 .. $ - si];
848     s.each!"a = 1 / a";
849     svd.vt[0 .. s.length].pack!1.map!"a".zip(s).each!"a.a[] *= a.b";
850     auto v = svd.vt[0 .. s.length].universal.transposed;
851     auto ut = svd.u.universal.transposed[0 .. s.length];
852     return v.mtimes(ut);
853 }
854 
855 ///
856 unittest
857 {
858     import mir.ndslice;
859 
860     auto a = [
861         64,  2,  3, 61, 60,  6,
862          9, 55, 54, 12, 13, 51,
863         17, 47, 46, 20, 21, 43,
864         40, 26, 27, 37, 36, 30,
865         32, 34, 35, 29, 28, 38,
866         41, 23, 22, 44, 45, 19,
867         49, 15, 14, 52, 53, 11,
868          8, 58, 59,  5,  4, 62].sliced(8, 6);
869 
870     auto b = a.pinv;
871 
872     auto result = [
873         0.0177, -0.0165, -0.0164,  0.0174,  0.0173, -0.0161, -0.0160,  0.0170,
874        -0.0121,  0.0132,  0.0130, -0.0114, -0.0112,  0.0124,  0.0122, -0.0106,
875        -0.0055,  0.0064,  0.0060, -0.0043, -0.0040,  0.0049,  0.0045, -0.0028,
876        -0.0020,  0.0039,  0.0046, -0.0038, -0.0044,  0.0064,  0.0070, -0.0063,
877        -0.0086,  0.0108,  0.0115, -0.0109, -0.0117,  0.0139,  0.0147, -0.0141,
878         0.0142, -0.0140, -0.0149,  0.0169,  0.0178, -0.0176, -0.0185,  0.0205].sliced(6, 8);
879 
880     import mir.math.common: approxEqual;
881     import mir.algorithm.iteration: all;
882 
883     assert(b.all!((a, b) => approxEqual(a, b, 1e-2, 1e-2))(result));
884 }
885 
886 /++
887 Covariance matrix.
888 
889 Params:
890     matrix = matrix whose rows represent observations and whose columns represent random variables.
891 Returns:
892     Normalized by `N-1` covariance matrix.
893 +/
894 Slice!(BlasType!Iterator*, 2)
895     cov(Iterator, SliceKind kind)(Slice!(Iterator, 2, kind) matrix)
896 {
897     import mir.math.sum: sum;
898     import mir.algorithm.iteration: each, eachUploPair;
899     alias A = BlasType!Iterator;
900     static if (kind == Contiguous)
901         auto mc = matrix.canonical;
902     else
903         alias mc = matrix;
904     auto m = mc.shape.uninitSlice!A.canonical;
905     auto s = m;
906 
907     auto factor = 1 / A(m.length!0 - 1).sqrt;
908     while(m.length!1)
909     {
910         auto shift = - mc.front!1.sum!A / m.length!0;
911         m.front!1.each!((ref a, b) { a = (b + shift) * factor; })(mc.front!1);
912         mc.popFront!1;
913         m.popFront!1;
914     }
915 
916     auto alpha = cast(A) 1;
917     auto beta = cast(A) 0;
918     auto c = [s.length!1, s.length!1].uninitSlice!A;
919     syrk(Uplo.Upper, alpha, s.universal.transposed, beta, c);
920 
921     c.eachUploPair!"b = a";
922     return c;
923 }
924 
925 ///
926 unittest
927 {
928     import mir.ndslice;
929 
930     import std.stdio;
931     import mir.ndslice;
932 
933     auto c = 8.magic[0..$-1].cov;
934 
935     auto result = [
936          350.0000, -340.6667, -331.3333,  322.0000,  312.6667, -303.3333, -294.0000,  284.6667,
937         -340.6667,  332.4762,  324.2857, -316.0952, -307.9048,  299.7143,  291.5238, -283.3333,
938         -331.3333,  324.2857,  317.2381, -310.1905, -303.1429,  296.0952,  289.0476, -282.0000,
939          322.0000, -316.0952, -310.1905,  304.2857,  298.3810, -292.4762, -286.5714,  280.6667,
940          312.6667, -307.9048, -303.1429,  298.3810,  293.6190, -288.8571, -284.0952,  279.3333,
941         -303.3333,  299.7143,  296.0952, -292.4762, -288.8571,  285.2381,  281.6190, -278.0000,
942         -294.0000,  291.5238,  289.0476, -286.5714, -284.0952,  281.6190,  279.1429, -276.6667,
943          284.6667, -283.3333, -282.0000,  280.6667,  279.3333, -278.0000, -276.6667,  275.3333].sliced(8, 8);
944     import mir.math.common: approxEqual;
945     import mir.algorithm.iteration: all;
946     assert(c.all!((a, b) => approxEqual(a, b, 1e-5, 1e-5))(result));
947 }
948 
949 /++
950 Matrix determinant.
951 +/
952 auto detSymmetric(Iterator, SliceKind kind)(char store, Slice!(Iterator, 2, kind) a)
953 in
954 {
955     assert(store == 'U' || store == 'L');
956     assert (a.length!0 == a.length!1, "matrix must be square");
957     assert (a.length!0, "matrix must not be empty");
958 }
959 do
960 {
961     import mir.algorithm.iteration: each;
962     import mir.ndslice.topology: diagonal;
963     import mir.math.numeric: ProdAccumulator;
964 
965     alias T = BlasType!Iterator;
966 
967     auto packed = uninitSlice!T(a.length * (a.length + 1) / 2);
968     auto ipiv = a.length.uninitSlice!lapackint;
969     int sign;
970     ProdAccumulator!T prod;
971     if (store == 'L')
972     {
973         auto pck = packed.stairs!"+"(a.length);
974         auto gen = a.stairs!"+";
975         pck.each!"a[] = b"(gen);
976         auto info = sptrf(pck, ipiv);
977         if (info > 0)
978             return cast(T) 0;
979         for (size_t j; j < ipiv.length; j++)
980         {
981             auto i = ipiv[j];
982             // 1x1 block at m[k,k]
983             if (i > 0)
984             {
985                 prod.put(pck[j].back);
986                 sign ^= i != j + 1; // i.e. row interchanged with another
987             }
988             else
989             {
990                 i = -i;
991                 auto offDiag = pck[j + 1][$ - 2];
992                 auto blockDet = pck[j].back * pck[j + 1].back - offDiag * offDiag;
993                 prod.put(blockDet);
994                 sign ^= i != j + 1 && i != j + 2; // row interchanged with other
995                 j++;
996             }
997         }
998     }
999     else
1000     {
1001         auto pck = packed.stairs!"-"(a.length);
1002         auto gen = a.stairs!"-";
1003         pck.each!"a[] = b"(gen);
1004         auto info = sptrf(pck, ipiv);
1005         if (info > 0)
1006             return cast(T) 0;
1007         for (size_t j; j < ipiv.length; j++)
1008         {
1009             auto i = ipiv[j];
1010             // 1x1 block at m[k,k]
1011             if (i > 0)
1012             {
1013                 prod.put(pck[j].front);
1014                 sign ^= i != j + 1; // i.e. row interchanged with another
1015             }
1016             else
1017             {
1018                 i = -i;
1019                 auto offDiag = pck[j][1];
1020                 auto blockDet = pck[j].front * pck[j + 1].front - offDiag * offDiag;
1021                 prod.put(blockDet);
1022                 sign ^= i != j + 1 && i != j + 2; // row interchanged with other
1023                 j++;
1024             }
1025         }
1026     }
1027     if(sign & 1)
1028         prod.x = -prod.x;
1029     return prod.prod;
1030 }
1031 
1032 /// ditto
1033 auto det(Iterator, SliceKind kind)(Slice!(Iterator, 2, kind) a)
1034 in
1035 {
1036     assert (a.length!0 == a.length!1, "matrix must be square");
1037 }
1038 do
1039 {
1040     import mir.ndslice.topology: diagonal, zip, iota;
1041     import mir.math.numeric: ProdAccumulator;
1042 
1043     alias T = BlasType!Iterator;
1044 
1045     auto m = a.as!T.slice.canonical;
1046     auto ipiv = a.length.uninitSlice!lapackint;
1047 
1048     // LU factorization
1049     auto info = m.getrf(ipiv);
1050 
1051     // If matrix is singular, determinant is zero.
1052     if (info > 0)
1053     {
1054         return cast(T) 0;
1055     }
1056 
1057     // The determinant is the product of the diagonal entries
1058     // of the upper triangular matrix. The array ipiv contains
1059     // the pivots.
1060     int sign;
1061     ProdAccumulator!T prod;
1062     foreach (tup; m.diagonal.zip(ipiv, [ipiv.length].iota(1)))
1063     {
1064         prod.put(tup.a);
1065         sign ^= tup.b != tup.c; // i.e. row interchanged with another
1066     }
1067     if(sign & 1)
1068         prod.x = -prod.x;
1069     return prod.prod;
1070 }
1071 
1072 ///
1073 unittest
1074 {
1075     import mir.ndslice;
1076     import mir.math;
1077 
1078     // Check for zero-determinant shortcut.
1079     auto ssing = [4, 2, 2, 1].sliced(2, 2);
1080     auto ssingd = det(ssing);
1081     assert (det(ssing) == 0);
1082     assert (detSymmetric('L', ssing) == 0);
1083 
1084     // check determinant of empty matrix
1085     assert(slice!double(0, 0).det == 1);
1086     // check determinant of zero matrix
1087     assert(repeat(0, 9).sliced(3, 3).det == 0);
1088 
1089     // General dense matrix.
1090     int dn = 101;
1091     auto d = uninitSlice!double(dn, dn);
1092     foreach (k; 0 .. dn)
1093     foreach (l; 0 .. dn)
1094         d[k,l] = 0.5 * (k == l ? (k + 1) * (k + 1) + 1 : 2 * (k + 1) * (l + 1));
1095 
1096     auto dd = det(d);
1097     import mir.math.common: approxEqual;
1098     assert (approxEqual(dd, 3.539152633479803e289, double.epsilon.sqrt));
1099 
1100     // Symmetric packed matrix
1101     auto spa = [ 1.0, -2, 3, 4, 5, -6, -7, -8, -9, 10].sliced.stairs!"+"(4);
1102     auto sp = [spa.length, spa.length].uninitSlice!double;
1103     import mir.algorithm.iteration: each;
1104     sp.stairs!"+".each!"a[] = b"(spa);
1105     assert (detSymmetric('L', sp).approxEqual(5874.0, double.epsilon.sqrt));
1106     assert (detSymmetric('U', sp.universal.transposed).approxEqual(5874.0, double.epsilon.sqrt));
1107 }
1108 
1109 /++
1110 Eigenvalues and eigenvectors POD.
1111 
1112 See_also: $(LREF eigSymmetric).
1113 +/
1114 struct EigSymmetricResult(T)
1115 {
1116     /// Eigenvalues
1117     Slice!(T*) values;
1118     /// Eigenvectors stored in rows
1119     Slice!(T*, 2) vectors;
1120 }
1121 
1122 /++
1123 Eigenvalues and eigenvectors of symmetric matrix.
1124 
1125 Returns:
1126     $(LREF EigSymmetricResult)
1127 +/
1128 auto eigSymmetric(Flag!"computeVectors" cv = Yes.computeVectors, Iterator, SliceKind kind)(char store, Slice!(Iterator, 2, kind) a)
1129 in
1130 {
1131     assert(store == 'U' || store == 'L');
1132     assert (a.length!0 == a.length!1, "matrix must be square");
1133     assert (a.length!0, "matrix must not be empty");
1134 }
1135 do
1136 {
1137     import mir.algorithm.iteration: each;
1138     import mir.ndslice.topology: diagonal;
1139     import mir.math.numeric: ProdAccumulator;
1140 
1141     alias T = BlasType!Iterator;
1142 
1143     auto packed = [a.length * (a.length + 1) / 2].uninitSlice!T;
1144     auto w = [a.length].uninitSlice!T;
1145     static if (cv)
1146     {
1147         auto z = [a.length, a.length].uninitSlice!T;
1148     }
1149     else
1150     {
1151         T[1] _vData = void;
1152         auto z = _vData[].sliced(1, 1);
1153     }
1154     auto work = [a.length * 3].uninitSlice!T;
1155     size_t info = void;
1156     auto jobz = cv ? 'V' : 'N';
1157     if (store == 'L')
1158     {
1159         auto pck = packed.stairs!"+"(a.length);
1160         auto gen = a.stairs!"+";
1161         pck.each!"a[] = b"(gen);
1162         info = spev!T(jobz, pck, w, z.canonical, work);
1163     }
1164     else
1165     {
1166         auto pck = packed.stairs!"-"(a.length);
1167         auto gen = a.stairs!"-";
1168         pck.each!"a[] = b"(gen);
1169         info = spev!T(jobz, pck, w, z.canonical, work);
1170     }
1171     import std.exception: enforce;
1172     import std.format: format;
1173     enforce (info == 0, format("The algorithm failed to converge." ~
1174         "%s off-diagonal elements of an intermediate tridiagonal form did not converge to zero.", info));
1175     static if (cv)
1176     {
1177         return EigSymmetricResult!T(w, z);
1178     }
1179     else
1180     {
1181         return EigSymmetricResult!T(w);
1182     }
1183 }
1184 
1185 ///
1186 unittest
1187 {
1188     import mir.ndslice;
1189 
1190     import mir.ndslice.slice: sliced;
1191     import mir.ndslice.topology: universal, map;
1192     import mir.ndslice.dynamic: transposed;
1193     import mir.math.common: approxEqual;
1194     import mir.algorithm.iteration: all;
1195 
1196     auto a = [
1197         1.0000, 0.5000, 0.3333, 0.2500,
1198         0.5000, 1.0000, 0.6667, 0.5000,
1199         0.3333, 0.6667, 1.0000, 0.7500,
1200         0.2500, 0.5000, 0.7500, 1.0000].sliced(4, 4);
1201 
1202     auto eigr = eigSymmetric('L', a);
1203 
1204     alias appr = all!((a, b) => approxEqual(a, b, 1e-3, 1e-3));
1205 
1206     assert(appr(eigr.values, [0.2078,0.4078,0.8482,2.5362]));
1207 
1208     auto test = [
1209          0.0693, -0.4422, -0.8105, 0.3778,
1210         -0.3618,  0.7420, -0.1877, 0.5322,
1211          0.7694,  0.0486,  0.3010, 0.5614,
1212         -0.5219, -0.5014,  0.4662, 0.5088].sliced(4, 4).transposed;
1213 
1214 
1215     foreach (i; 0 .. 4)
1216         assert(appr(eigr.vectors[i], test[i]) || appr(eigr.vectors[i].map!"-a", test[i]));
1217 }
1218 
1219 version (unittest)
1220 {
1221 /++
1222 Swaps rows of input matrix
1223 Params
1224     ipiv = pivot points
1225 Returns:
1226     a = shifted matrix
1227 +/
1228     private void moveRows(Iterator, SliceKind kind)
1229                         (Slice!(Iterator, 2, kind) a,
1230                         Slice!(lapackint*) ipiv)
1231     {
1232         import mir.algorithm.iteration: each;
1233         foreach_reverse(i;0..ipiv.length)
1234         {
1235             if(ipiv[i] == i + 1)
1236                 continue;
1237             each!swap(a[i], a[ipiv[i] - 1]);
1238         }
1239     }
1240 }
1241 
1242 unittest
1243 {
1244     import mir.ndslice;
1245 
1246     auto A = 
1247            [ 9,  9,  9,
1248              8,  8,  8,
1249              7,  7,  7,
1250              6,  6,  6,
1251              5,  5,  5,
1252              4,  4,  4,
1253              3,  3,  3,
1254              2,  2,  2,
1255              1,  1,  1,
1256              0,  0,  0 ]
1257              .sliced(10, 3)
1258              .as!double.slice;
1259     auto ipiv = [ lapackint(10), 9, 8, 7, 6, 6, 7, 8, 9, 10 ].sliced(10);
1260     moveRows(A, ipiv);
1261 
1262     auto B = 
1263            [ 0,  0,  0,
1264              1,  1,  1,
1265              2,  2,  2,
1266              3,  3,  3,
1267              4,  4,  4,
1268              5,  5,  5,
1269              6,  6,  6,
1270              7,  7,  7,
1271              8,  8,  8,
1272              9,  9,  9 ]
1273              .sliced(10, 3)
1274              .as!double.slice;
1275 
1276     import mir.algorithm.iteration: equal;
1277     assert(equal!((a, b) => fabs(a - b) < 1e-12)(B, A));
1278 }
1279 
1280 unittest
1281 {
1282     auto A = 
1283            [ 1,  1,  1,
1284              2,  2,  2,
1285              3,  3,  3,
1286              4,  4,  4,
1287              5,  5,  5,
1288              6,  6,  6,
1289              7,  7,  7,
1290              8,  8,  8,
1291              9,  9,  9,
1292              0,  0,  0 ]
1293              .sliced(10, 3)
1294              .as!double.slice;
1295     auto ipiv = [ lapackint(2), 3, 4, 5, 6, 7, 8, 9, 10, 10 ].sliced(10);
1296     moveRows(A, ipiv);
1297     
1298     auto B = 
1299            [ 0,  0,  0,
1300              1,  1,  1,
1301              2,  2,  2,
1302              3,  3,  3,
1303              4,  4,  4,
1304              5,  5,  5,
1305              6,  6,  6,
1306              7,  7,  7,
1307              8,  8,  8,
1308              9,  9,  9 ]
1309              .sliced(10, 3)
1310              .as!double.slice;
1311 
1312     import mir.algorithm.iteration: equal;
1313     assert(equal!((a, b) => fabs(a - b) < 1e-12)(B, A));
1314 }
1315 
1316 ///LUResult consist lu factorization.
1317 struct LUResult(T)
1318 {
1319     /++
1320     Matrix in witch lower triangular is L part of factorization
1321     (diagonal elements of L are not stored), upper triangular
1322     is U part of factorization.
1323     +/
1324     Slice!(T*, 2, Canonical) lut;
1325     /++
1326     The pivot indices, for 1 <= i <= min(M,N), row i of the matrix
1327     was interchanged with row ipiv(i).
1328     +/
1329     Slice!(lapackint*) ipiv;
1330     ///L part of the factorization.
1331     auto l() @property
1332     {
1333         import mir.algorithm.iteration: eachUpper;
1334         auto l = lut.transposed[0..lut.length!1, 0..min(lut.length!0, lut.length!1)].slice.canonical;
1335         l.eachUpper!"a = 0";
1336         l.diagonal[] = 1;
1337         return l;
1338     }
1339     ///U part of the factorization.
1340     auto u() @property
1341     {
1342         import mir.algorithm.iteration: eachLower;
1343         auto u = lut.transposed[0..min(lut.length!0, lut.length!1), 0..lut.length!0].slice.canonical;
1344         u.eachLower!"a = 0";
1345         return u;
1346     }
1347 
1348     ///
1349     auto solve(Flag!"allowDestroy" allowDestroy = No.allowDestroy, Iterator, size_t N, SliceKind kind)(
1350         char trans,
1351         Slice!(Iterator, N, kind) b)
1352     {
1353         return luSolve!(allowDestroy)(trans, lut, ipiv, b);
1354     }
1355 }
1356 
1357 /++
1358 Computes LU factorization of a general 'M x N' matrix 'A' using partial
1359 pivoting with row interchanges.
1360 The factorization has the form:
1361     \A = P * L * U
1362 Where P is a permutation matrix, L is lower triangular with unit
1363 diagonal elements (lower trapezoidal if m > n), and U is upper
1364 triangular (upper trapezoidal if m < n).
1365 Params:
1366     allowDestroy = flag to delete the source matrix.
1367     a = input 'M x N' matrix for factorization.
1368 Returns: $(LREF LUResalt)
1369 +/
1370 auto luDecomp(Flag!"allowDestroy" allowDestroy = No.allowDestroy,
1371               Iterator, SliceKind kind)
1372              (Slice!(Iterator, 2, kind) a)
1373 {
1374     alias T = BlasType!Iterator;
1375     auto ipiv = uninitSlice!lapackint(min(a.length!0, a.length!1));
1376     auto b = a.transposed;
1377     auto m = (allowDestroy && b._stride!1 == 1) ? b.assumeCanonical : a.transposed.as!T.slice.canonical;
1378     
1379     getrf(m, ipiv);
1380     return LUResult!T(m, ipiv);
1381 }
1382 
1383 /++
1384 Solves a system of linear equations
1385     \A * X = B, or
1386     \A**T * X = B
1387 with a general 'N x N' matrix 'A' using the LU factorization computed by luDecomp.
1388 Params:
1389     allowDestroy = flag to delete the source matrix.
1390     lut = factorization of matrix 'A', A = P * L * U.
1391     ipiv = the pivot indices from luDecomp.
1392     b = the right hand side matrix B.
1393     trans = specifies the form of the system of equations:
1394           'N': A * X = B (No transpose)
1395           'T': A**T * X = B (Transpose)
1396           'C': A**T * X = B (Conjugate transpose = Transpose)
1397 Returns:
1398     Return solve of the system linear equations.
1399 +/
1400 auto luSolve(Flag!"allowDestroy" allowDestroy = No.allowDestroy, SliceKind kindB, size_t N, IteratorB, IteratorLU)(
1401     char trans,
1402     Slice!(IteratorLU, 2, Canonical) lut,
1403     Slice!(lapackint*) ipiv,
1404     Slice!(IteratorB, N, kindB) b,
1405     )
1406 in
1407 {
1408     assert(lut.length!0 == lut.length!1, "matrix must be squared");
1409     assert(ipiv.length == lut.length, "size of ipiv must be equal to the number of rows a");
1410     assert(lut.length!1 == b.length!0, "number of columns a should be equal to the number of rows b");
1411 }
1412 do
1413 {
1414     alias LU = BlasType!IteratorLU;
1415     alias B = BlasType!IteratorB;
1416     alias T = CommonType!(LU, B);
1417     static if(is(T* == IteratorLU))
1418         auto lut_ = lut;
1419     else
1420         auto lut_ = lut.as!T.slice.canonical;
1421     
1422     //convect vector to matrix.
1423     static if(N == 1)
1424         auto k = b.sliced(1, b.length);
1425     else
1426         auto k = b.transposed;
1427 
1428     static if(is(IteratorB == T*))
1429         auto m = (allowDestroy && k._stride!1 == 1) ? k.assumeCanonical : k.as!T.slice.canonical;
1430     else
1431         auto m = k.as!T.slice.canonical;
1432     getrs!T(trans, lut_, m, ipiv);
1433     return m.transposed;
1434 }
1435 
1436 unittest
1437 {
1438     auto A =
1439         [ 1,  4, -3,  5,  6,
1440          -2,  8,  5,  7,  8,
1441           3,  4,  7,  9,  1,
1442           2,  4,  6,  3,  2,
1443           6,  8,  3,  5,  2 ]
1444             .sliced(5, 5)
1445             .as!double.slice
1446             .canonical;
1447     auto B =
1448         [ 1,  3,  4,  8,  0,  0,  0,
1449           2,  1,  7,  1,  0,  0,  0,
1450           3,  5,  7,  7,  0,  0,  0,
1451           4,  4,  9,  8,  0,  0,  0,
1452           5,  5,  8,  1,  0,  0,  0 ]
1453             .sliced(5, 7)
1454             .as!double.slice
1455             .universal;
1456 
1457     auto B_ = B[0..$, 0..4];
1458     auto LU = A.luDecomp();
1459     auto m = luSolve('N', LU.lut, LU.ipiv, B_);
1460 
1461     import mir.math.common: approxEqual;
1462     import mir.algorithm.iteration: equal;
1463     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1464     assert(appr(mtimes(A, m), B_));
1465 }
1466 
1467 ///
1468 unittest
1469 {
1470     import mir.math;
1471     import mir.ndslice;
1472 
1473     auto A =
1474         [ 1,  4, -3,  5,  6,
1475          -2,  8,  5,  7,  8,
1476           3,  4,  7,  9,  1,
1477           2,  4,  6,  3,  2,
1478           6,  8,  3,  5,  2 ]
1479             .sliced(5, 5)
1480             .as!double.slice
1481             .canonical;
1482     
1483     import mir.random.variable;
1484     import mir.random.algorithm;
1485     auto B = randomSlice(uniformVar(-100, 100), 5, 100);
1486     
1487     auto LU = A.luDecomp();
1488     auto X = LU.solve('N', B);
1489 
1490     import mir.algorithm.iteration: equal;
1491     assert(equal!((a, b) => fabs(a - b) < 1e-12)(mtimes(A, X), B));
1492 }
1493 
1494 unittest
1495 {
1496     auto A =
1497         [ 1,  4, -3,  5,  6,
1498          -2,  8,  5,  7,  8,
1499           3,  4,  7,  9,  1,
1500           2,  4,  6,  3,  2,
1501           6,  8,  3,  5,  2 ]
1502             .sliced(5, 5)
1503             .as!double.slice
1504             .canonical;
1505     auto B =
1506         [ 2,  8,  3,  5,  8,
1507           8,  1,  4,  9,  86,
1508           1,  6,  7,  1,  67,
1509           6,  1,  5,  4,  45,
1510           1,  2,  3,  1,  11 ]
1511             .sliced(5, 5)
1512             .as!double.slice
1513             .universal;
1514     auto C = B.slice;
1515 
1516     auto LU = A.luDecomp();
1517     auto m = luSolve!(Yes.allowDestroy)('N', LU.lut, LU.ipiv, B.transposed);
1518     auto m2 = LU.solve('N', C);
1519 
1520     import mir.math.common: approxEqual;
1521     import mir.algorithm.iteration: equal;
1522     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1523     assert(appr(mtimes(A, m), C.transposed));
1524     assert(equal!approxEqual(mtimes(A, m2), C));
1525 }
1526 
1527 unittest
1528 {
1529     auto A =
1530         [ 1,  4, -3,  5,  6,
1531          -2,  8,  5,  7,  8,
1532           3,  4,  7,  9,  1,
1533           2,  4,  6,  3,  2,
1534           6,  8,  3,  5,  2 ]
1535             .sliced(5, 5)
1536             .as!float.slice
1537             .canonical;
1538     auto B = [ 1,  2,  3,  4,  5 ].sliced(5).as!double.slice;
1539     auto C = B.slice.sliced(5, 1);
1540 
1541     auto LU = A.luDecomp();
1542     auto m = luSolve!(Yes.allowDestroy)('N', LU.lut, LU.ipiv, B);
1543 
1544     import mir.math.common: approxEqual;
1545     import mir.algorithm.iteration: equal;
1546     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1547     assert(appr(mtimes(A, m), C));
1548 }
1549 
1550 unittest
1551 {
1552     auto A =
1553         [ 1,  4, -3,  5,  6,
1554          -2,  8,  5,  7,  8,
1555           3,  4,  7,  9,  1,
1556           2,  4,  6,  3,  2,
1557           6,  8,  3,  5,  2 ]
1558             .sliced(5, 5)
1559             .as!double.slice
1560             .canonical;
1561     auto B = [ 1,  15,  4,  5,  8,
1562                3,  20,  1,  9,  11 ].sliced(5, 2).as!float.slice;
1563 
1564     auto LU = A.luDecomp();
1565     auto m = luSolve('N', LU.lut, LU.ipiv, B);
1566 
1567     import mir.math.common: approxEqual;
1568     import mir.algorithm.iteration: equal;
1569     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1570     assert(appr(mtimes(A, m), B));
1571 }
1572 
1573 unittest
1574 {
1575     auto A =
1576         [ 11,  14, -31,  53,  62,
1577          -92,  83,  52,  74,  83,
1578           31,  45,  73,  96,  17,
1579           23,  14,  65,  35,  26,
1580           62,  28,  34,  51,  25 ]
1581             .sliced(5, 5)
1582             .as!float.slice
1583             .universal;
1584     auto B =
1585         [ 6,  1,  3,  1,  11,
1586          12,  5,  7,  6,  78,
1587           8,  4,  1,  5,  54,
1588           3,  1,  8,  1,  45,
1589           1,  6,  8,  6,  312 ]
1590             .sliced(5, 5)
1591             .as!double.slice;
1592     auto B2 = B.slice;
1593     auto C = B.slice;
1594 
1595     auto LU = luDecomp(A.transposed);
1596     auto m = luSolve!(Yes.allowDestroy)('T', LU.lut, LU.ipiv, B);
1597     auto m2 = luSolve!(Yes.allowDestroy)('N', LU.lut, LU.ipiv, B2);
1598 
1599     import mir.math.common: approxEqual;
1600     import mir.algorithm.iteration: equal;
1601     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1602 
1603     assert(appr(mtimes(A, m), C));
1604     assert(appr(mtimes(A.transposed, m2), C));
1605 }
1606 
1607 unittest
1608 {
1609     auto A =
1610         [ 54,  93,  14,  44,  33,
1611           51,  85,  28,  81,  75,
1612           89,  17,  15,  44,  58,
1613           75,  80,  18,  35,  14,
1614           21,  48,  72,  21,  88 ]
1615             .sliced(5, 5)
1616             .as!double.slice
1617             .universal;
1618     auto B =
1619         [ 5,  7,  8,  3,  78,
1620           1,  2,  5,  4,  5,
1621           2,  4,  1,  5,  15,
1622           1,  1,  4,  1,  154,
1623           1,  3,  1,  8,  17 ]
1624             .sliced(5, 5)
1625             .as!float.slice
1626             .canonical;
1627 
1628     auto LU = A.luDecomp();
1629     auto m = luSolve('N', LU.lut, LU.ipiv, B);
1630 
1631     import mir.math.common: approxEqual;
1632     import mir.algorithm.iteration: equal;
1633     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1634     assert(appr(mtimes(A, m), B));
1635 }
1636 
1637 unittest
1638 {
1639     
1640     auto B =
1641         [ 1,  4, -3,
1642          -2,  8,  5,
1643           3,  4,  7,
1644           2,  4,  6 ]
1645             .sliced(4, 3)
1646             .as!double.slice;
1647 
1648     auto LU = B.luDecomp!(Yes.allowDestroy)();
1649     LU.l; LU.u;
1650     auto res = mtimes(LU.l, LU.u);
1651     moveRows(res, LU.ipiv);
1652 
1653     import mir.algorithm.iteration: equal;
1654     import mir.math.common: approxEqual;
1655     assert(res.equal!approxEqual(B));
1656 }
1657 
1658 unittest
1659 {
1660     import mir.ndslice;
1661 
1662     auto B =
1663         [ 3, -7, -2,  2,
1664          -3,  5,  1,  0,
1665           6, -4,  0, -5,
1666          -9,  5, -5, 12 ]
1667             .sliced(4, 4)
1668             .as!double.slice
1669             .canonical;
1670     auto C = B.transposed.slice;
1671 
1672     auto LU = B.transposed.luDecomp!(Yes.allowDestroy)();
1673     auto res = mtimes(LU.l, LU.u);
1674     moveRows(res, LU.ipiv);
1675 
1676     import mir.algorithm.iteration: equal;
1677     import mir.math.common: approxEqual;
1678     assert(res.equal!approxEqual(C));
1679 }
1680 
1681 ///Consist LDL factorization;
1682 struct LDLResult(T)
1683 {
1684     /++
1685     uplo = 'U': Upper triangle is stored;
1686          'L': lower triangle is stored.
1687     +/
1688     char uplo;
1689     /++
1690     Matrix in witch lower triangular matrix is 'L' part of
1691     factorization, diagonal is 'D' part.
1692     +/
1693     Slice!(T*, 2, Canonical) matrix;
1694     /++
1695     The pivot indices.
1696     If ipiv(k) > 0, then rows and columns k and ipiv(k) were
1697     interchanged and D(k, k) is a '1 x 1' diagonal block.
1698     If ipiv(k) = ipiv(k + 1) < 0, then rows and columns k+1 and
1699     -ipiv(k) were interchanged and D(k:k+1, k:k+1) is a '2 x 2'
1700     diagonal block.
1701     +/
1702     Slice!(lapackint*) ipiv;
1703     /++
1704     Return solves a system of linear equations
1705         \A * X = B,
1706     using LDL factorization.
1707     +/
1708     auto solve(Flag!"allowDestroy" allowDestroy = No.allowDestroy, SliceKind kindB, size_t N, IteratorB)(
1709         Slice!(IteratorB, N, kindB) b)
1710     {
1711         return ldlSolve!(allowDestroy)(uplo, matrix, ipiv, b);
1712     }
1713 }
1714 
1715 /++
1716 Computes the factorization of a real symmetric matrix A using the
1717 Bunch-Kaufman diagonal pivoting method.
1718 The for of the factorization is:
1719     \A = L*D*L**T
1720 Where L is product if permutation and unit lower triangular matrices,
1721 and D is symmetric and block diagonal with '1 x 1' and '2 x 2'
1722 diagonal blocks.
1723 Params:
1724     allowDestroy = flag to delete the source matrix.
1725     a = input symmetric 'n x n' matrix for factorization.
1726     uplo = 'U': Upper triangle is stored;
1727          'L': lower triangle is stored.
1728 Returns:$(LREF LDLResult)
1729 +/
1730 auto ldlDecomp(Flag!"allowDestroy" allowDestroy = No.allowDestroy, Iterator, SliceKind kind)(
1731     char uplo,
1732     Slice!(Iterator, 2, kind) a)
1733 in
1734 {
1735     assert(a.length!0 == a.length!1, "matrix must be squared");
1736 }
1737 do
1738 {
1739     alias T = BlasType!Iterator;
1740     auto work = [T.sizeof * a.length].uninitSlice!T;
1741     auto ipiv = a.length.uninitSlice!lapackint;
1742     auto m = (allowDestroy && a._stride!1 == 1) ? a.assumeCanonical : a.transposed.as!T.slice.canonical;
1743 
1744     sytrf!T(uplo, m, ipiv, work);
1745     return LDLResult!T(uplo, m, ipiv);
1746 }
1747 
1748 /++
1749 Solves a system of linear equations \A * X = B with symmetric matrix 'A' using the
1750 factorization
1751 \A = U * D * U**T, or
1752 \A = L * D * L**T
1753 computed by ldlDecomp.
1754 Params:
1755     allowDestroy = flag to delete the source matrix.
1756     a = 'LD' or 'UD' matrix computed by ldlDecomp.
1757     ipiv = details of the interchanges and the block structure of D as determined by ldlDecomp.
1758     b = the right hand side matrix.
1759     uplo = specifies whether the details of the factorization are stored as an upper or
1760            lower triangular matrix:
1761          'U': Upper triangular, form is \A = U * D * U**T;
1762          'L': Lower triangular, form is \A = L * D * L**T.
1763 Returns:
1764     The solution matrix.
1765 +/
1766 auto ldlSolve(Flag!"allowDestroy" allowDestroy = No.allowDestroy, SliceKind kindB, size_t N, IteratorB, IteratorA)(
1767     char uplo,
1768     Slice!(IteratorA, 2, Canonical) a,
1769     Slice!(lapackint*) ipiv,
1770     Slice!(IteratorB, N, kindB) b)
1771 in
1772 {
1773     assert(a.length!0 == a.length!1, "matrix must be squared");
1774     assert(ipiv.length == a.length, "size of ipiv must be equal to the number of rows a");
1775     assert(a.length!1 == b.length!0, "number of columns a should be equal to the number of rows b");
1776 }
1777 do
1778 {
1779     alias A = BlasType!IteratorA;
1780     alias B = BlasType!IteratorB;
1781     alias T = CommonType!(A, B);
1782     static if(is(T* == IteratorA))
1783         auto a_ = a;
1784     else
1785         auto a_ = a.as!T.slice.canonical;
1786     
1787     //convect vector to matrix.
1788     static if(N == 1)
1789         auto k = b.sliced(1, b.length);
1790     else
1791         auto k = b.transposed;
1792 
1793     auto work = [T.sizeof * a.length].uninitSlice!T;
1794     static if(is(IteratorB == T*))
1795         auto m = (allowDestroy && k._stride!1 == 1) ? k.assumeCanonical : k.as!T.slice.canonical;
1796     else
1797         auto m = k.as!T.slice.canonical;
1798     sytrs2!T(a_, m, ipiv, work, uplo);
1799     return m.transposed;
1800 }
1801 
1802 ///
1803 unittest
1804 {
1805     import mir.ndslice;
1806 
1807     auto A =
1808         [ 2.07,  3.87,  4.20, -1.15,
1809           3.87, -0.21,  1.87,  0.63,
1810           4.20,  1.87,  1.15,  2.06,
1811          -1.15,  0.63,  2.06, -1.81 ]
1812             .sliced(4, 4)
1813             .as!double.slice
1814             .canonical;
1815 
1816     import mir.random.variable;
1817     import mir.random.algorithm;
1818     auto B = randomSlice(uniformVar(-100, 100), 4, 100);
1819 
1820     auto LDL = ldlDecomp('L', A);
1821     auto X = LDL.solve(B);
1822 
1823     import mir.math.common: approxEqual;
1824     import mir.algorithm.iteration: equal;
1825     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1826     assert(appr(mtimes(A, X), B));
1827 }
1828 
1829 unittest
1830 {
1831     auto A =
1832         [ 9, -1,  2,
1833          -1,  8, -5,
1834           2, -5,  7 ]
1835             .sliced(3, 3)
1836             .as!float.slice
1837             .canonical;
1838     auto A_ = A.slice;
1839 
1840     auto B =
1841         [ 5,  7,  1,
1842           1,  8,  5,
1843           9,  3,  2 ]
1844           .sliced(3, 3)
1845           .as!double.slice
1846           .canonical;
1847     auto B_ = B.slice;
1848 
1849     auto LDL = ldlDecomp!(Yes.allowDestroy)('L', A);
1850     auto X = ldlSolve!(Yes.allowDestroy)(LDL.uplo, A, LDL.ipiv, B.transposed);
1851 
1852     import mir.math.common: approxEqual;
1853     import mir.algorithm.iteration: equal;
1854     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1855     assert(appr(mtimes(A_, X), B_.transposed));
1856 }
1857 
1858 unittest
1859 {
1860     auto A =
1861         [10, 20, 30,
1862          20, 45, 80,
1863          30, 80, 171 ]
1864             .sliced(3, 3)
1865             .as!double.slice
1866             .canonical;
1867     auto B = [ 1, 4, 7 ].sliced(3).as!float.slice.canonical;
1868     auto B_ = B.sliced(3, 1);
1869 
1870     auto LDL = ldlDecomp('L', A);
1871     auto X = LDL.solve(B);
1872     
1873     import mir.math.common: approxEqual;
1874     import mir.algorithm.iteration: equal;
1875     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
1876     assert(appr(mtimes(A, X), B_));
1877 }
1878 
1879 struct choleskyResult(T)
1880 {
1881     /++
1882     If uplo = 'L': lower triangle of 'matrix' is stored.
1883     If uplo = 'U': upper triangle of 'matrix' is stored.
1884     +/
1885     char uplo;
1886     /++
1887     if uplo = 'L', the leading 'N x N' lower triangular part of A
1888     contains the lower triangular part of the matrix A, and the
1889     strictly upper triangular part of A is zeroed.
1890     if uplo = 'U', the leading 'N x N' upper triangular part of A
1891     contains the upper triangular part of the matrix A, and the
1892     strictly lower triangular part of A is zeroed.
1893     +/
1894     Slice!(T*, 2, Canonical) matrix;
1895     /++
1896     Return solves a system of linear equations
1897         \A * X = B,
1898     using Cholesky factorization.
1899     +/
1900     auto solve(Flag!"allowDestroy" allowDestroy = No.allowDestroy,
1901                Iterator, size_t N, SliceKind kind)
1902               (Slice!(Iterator, N, kind) b)
1903     {
1904         return choleskySolve!(allowDestroy)(uplo, matrix, b);
1905     }
1906 }
1907 /++
1908 Computs Cholesky decomposition of symmetric positive definite matrix 'A'.
1909 The factorization has the form:
1910     \A = U**T * U, if UPLO = 'U', or
1911     \A = L * L**T, if UPLO = 'L'
1912 Where U is an upper triangular matrix and L is lower triangular.
1913 Params:
1914     allowDestroy = flag to delete the source matrix.
1915     a = symmetric 'N x N' matrix.
1916     uplo = if uplo is Upper, then upper triangle of A is stored, else
1917     lower.
1918 Returns: $(LREF choleskyResult)
1919 +/
1920 auto choleskyDecomp( Flag!"allowDestroy" allowDestroy = No.allowDestroy, Iterator, SliceKind kind)(
1921     char uplo,
1922     Slice!(Iterator, 2, kind) a)
1923 in
1924 {
1925     assert(uplo == 'L' || uplo == 'U');
1926     assert(a.length!0 == a.length!1, "matrix must be squared");
1927 }
1928 do
1929 {
1930     import mir.exception: MirException;
1931     alias T = BlasType!Iterator;
1932     static if(is(Iterator == T*))
1933         auto m = (allowDestroy && a._stride!1 == 1) ? a.assumeCanonical : a.as!T.slice.canonical;
1934     else
1935         auto m = a.as!T.slice.canonical;
1936     if (auto info = potrf!T(uplo == 'U' ? 'L' : 'U', m))
1937         throw new MirException("Leading minor of order ", info, " is not positive definite, and the factorization could not be completed.");
1938     import mir.algorithm.iteration: eachUploPair;
1939     auto d = m.universal;
1940     if (uplo == 'U')
1941         d = d.transposed;
1942     d.eachUploPair!("a = 0");
1943     return choleskyResult!T(uplo, m);
1944 }
1945 
1946 ///
1947 unittest
1948 {
1949     import mir.algorithm.iteration: equal, eachUploPair;
1950     import mir.ndslice;
1951     import mir.random.algorithm;
1952     import mir.random.variable;
1953     import mir.math.common: approxEqual;
1954 
1955     auto A =
1956            [ 25, double.nan, double.nan,
1957              15, 18,  double.nan,
1958              -5,  0, 11 ]
1959              .sliced(3, 3);
1960     
1961     auto B = randomSlice(uniformVar(-100, 100), 3, 100);
1962 
1963     auto C = choleskyDecomp('L', A);
1964     auto X = C.solve(B);
1965 
1966     A.eachUploPair!"a = b";
1967     assert(equal!approxEqual(mtimes(A, X), B));
1968 }
1969 
1970 /++
1971     Solves a system of linear equations A * X = B with a symmetric matrix A using the
1972     Cholesky factorization:
1973     \A = U**T * U or
1974     \A = L * L**T
1975     computed by choleskyDecomp.
1976 Params:
1977     allowDestroy = flag to delete the source matrix.
1978     c = the triangular factor 'U' or 'L' from the Cholesky factorization
1979         \A = U**T * U or
1980         \A = L * L**T,
1981     as computed by choleskyDecomp.
1982     b = the right hand side matrix.
1983     uplo = 'U': Upper triangle of A is stored;
1984          'L': Lower triangle of A is stored.
1985 Returns:
1986     The solution matrix X.
1987 +/
1988 auto choleskySolve(Flag!"allowDestroy" allowDestroy = No.allowDestroy, SliceKind kindB, size_t N, IteratorB, IteratorC)(
1989     char uplo,
1990     Slice!(IteratorC, 2, Canonical) c,
1991     Slice!(IteratorB, N, kindB) b)
1992 in
1993 {
1994     assert(uplo == 'L' || uplo == 'U');
1995     assert(c.length!0 == c.length!1, "matrix must be squared");
1996     assert(c.length!1 == b.length!0, "number of columns a should be equal to the number of rows b");
1997 }
1998 do
1999 {
2000     uplo = uplo == 'U' ? 'L' : 'U';
2001     alias B = BlasType!IteratorB;
2002     alias C = BlasType!IteratorC;
2003     alias T = CommonType!(B, C);
2004     static if(is(T* == IteratorC))
2005         auto c_ = c;
2006     else
2007         auto c_ = c.as!T.slice.canonical;
2008 
2009     //convect vector to matrix.
2010     static if(N == 1)
2011         auto k = b.sliced(1, b.length);
2012     else
2013         auto k = b.transposed;
2014 
2015     static if(is(IteratorB == T*))
2016         auto m = (allowDestroy && k._stride!1 == 1) ? k.assumeCanonical : k.as!T.slice.canonical;
2017     else
2018         auto m = k.as!T.slice.canonical;
2019     potrs!T(uplo, c_, m);
2020     return m.transposed;
2021 }
2022 
2023 ///
2024 unittest
2025 {
2026     import mir.ndslice.slice: sliced;
2027     import mir.ndslice.topology: as;
2028     import std.typecons: Flag, Yes;
2029 
2030     auto A =
2031             [ 1.0,  1,  3,
2032               1  ,  5,  5,
2033               3  ,  5, 19 ].sliced(3, 3);
2034 
2035     auto B = [ 10.0,  157,  80 ].sliced;
2036     auto C_ = B.dup.sliced(3, 1);
2037 
2038     auto C = choleskyDecomp('U', A);
2039     auto X = choleskySolve!(Yes.allowDestroy)(C.uplo, C.matrix, B);
2040 
2041     import mir.math.common: approxEqual;
2042     import mir.algorithm.iteration: equal;
2043     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
2044     assert(appr(mtimes(A, X), C_));
2045 }
2046 
2047 unittest
2048 {
2049     auto A =
2050             [6.0f,  15,  55,
2051              15,  55, 225,
2052              55, 225, 979 ].sliced(3, 3).canonical;
2053     auto B =
2054             [ 7.0,  3,
2055               2,  1,
2056               1,  8 ].sliced(3, 2).universal;
2057 
2058     auto C = choleskyDecomp('L', A);
2059     auto X = choleskySolve(C.uplo, C.matrix, B);
2060 
2061     import mir.math.common: approxEqual;
2062     import mir.algorithm.iteration: equal;
2063     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
2064 
2065     assert(appr(mtimes(A, X), B));
2066 }
2067 
2068 
2069 ///
2070 struct QRResult(T)
2071 {
2072     /++
2073     Matrix in witch the elements on and above the diagonal of the array contain the min(M, N) x N
2074     upper trapezoidal matrix 'R' (R is upper triangular if m >= n). The elements below the
2075     diagonal, with the array tau, represent the orthogonal matrix 'Q' as product of min(m, n).
2076     +/
2077     Slice!(T*, 2, Canonical) matrix;
2078     ///The scalar factors of the elementary reflectors
2079     Slice!(T*) tau;
2080     /++
2081     Solve the least squares problem:
2082         \min ||A * X - B||
2083     Using the QR factorization:
2084         \A = Q * R
2085     computed by qrDecomp.
2086     +/
2087     auto solve(Flag!"allowDestroy" allowDestroy = No.allowDestroy,
2088                Iterator, size_t N, SliceKind kind)
2089               (Slice!(Iterator, N, kind) b)
2090     {
2091         return qrSolve!(allowDestroy)(matrix, tau, b);
2092     }
2093     
2094     /++
2095     Extract the Q matrix
2096     +/
2097     auto Q(Flag!"allowDestroy" allowDestroy = No.allowDestroy)
2098     {
2099         auto work = [matrix.length].uninitSlice!T;
2100 
2101         auto m = (allowDestroy && matrix._stride!1 == 1) ? matrix.assumeCanonical : matrix.as!T.slice.canonical;
2102 
2103         static if(is(T == double) || is(T == float))
2104             orgqr!T(m, tau, work);
2105         else
2106             ungqr!T(m, tau, work);
2107         return m.transposed;
2108     }
2109     
2110     /++
2111     Extract the R matrix
2112     +/
2113     auto R()
2114     {
2115         import mir.algorithm.iteration: eachLower;
2116 
2117         auto r = [tau.length, tau.length].uninitSlice!T;
2118         if (matrix.shape[0] == matrix.shape[1]) {
2119             r[] = matrix.transposed.slice;
2120         } else {
2121             r[] = matrix[0..tau.length, 0..tau.length].transposed.slice;
2122         }
2123         r.eachLower!((ref a) {a = cast(T)0;});
2124         return r.universal;
2125     }
2126 
2127     /++
2128     Reconstruct the original matrix given a QR decomposition
2129     +/
2130     auto reconstruct()
2131     {
2132         auto q = Q();
2133         auto r = R();
2134         return reconstruct(q, r);
2135     }
2136 
2137     /++
2138     Reconstruct the original matrix given a QR decomposition
2139     +/
2140     auto reconstruct(T, U)(T q, U r)
2141         if (isMatrix!T && isMatrix!U)
2142     {
2143         return mtimes(q, r).universal;
2144     }
2145 }
2146 
2147 ///
2148 unittest
2149 {
2150     import mir.ndslice;
2151 
2152     auto A =
2153             [ 1,  1,  0,
2154               1,  0,  1,
2155               0,  1,  1 ]
2156               .sliced(3, 3)
2157               .as!double.slice;
2158 
2159     auto Q_test =
2160             [ -0.7071068,  0.4082483,  -0.5773503,
2161               -0.7071068, -0.4082483,   0.5773503,
2162                        0,  0.8164966,   0.5773503]
2163               .sliced(3, 3)
2164               .as!double.slice;
2165 
2166     auto R_test =
2167             [ -1.414214,  -0.7071068,   -0.7071068,
2168                       0,   1.2247449,    0.4082483,
2169                       0,           0,    1.1547005]
2170               .sliced(3, 3)
2171               .as!double.slice;
2172 
2173     auto val = qrDecomp(A);
2174 
2175     //saving these values to doublecheck they don't change later
2176     auto val_matrix = val.matrix.slice;
2177     auto val_tau = val.tau.slice;
2178 
2179     import mir.math.common: approxEqual;
2180     import mir.ndslice : equal;
2181     
2182     auto r = val.R;
2183     assert(equal!approxEqual(val.R, R_test));
2184 
2185     auto q = val.Q;
2186     assert(equal!approxEqual(val.Q, Q_test));
2187 
2188     //double-checking values do not change
2189     assert(equal!approxEqual(val_matrix, val.matrix));
2190     assert(equal!approxEqual(val_tau, val.tau));
2191 
2192     auto a = val.reconstruct;
2193     assert(equal!approxEqual(A, a));
2194 }
2195 
2196 unittest
2197 {
2198     auto A =
2199             [  3,  -6,
2200                4,  -8,
2201                0,   1]
2202               .sliced(3, 2)
2203               .as!double.slice;
2204 
2205     auto Q_check =
2206             [ -0.6,  0,
2207               -0.8,  0,
2208                0.0, -1]
2209               .sliced(3, 2)
2210               .as!double.slice;
2211 
2212     auto R_check =
2213             [ -5,  10,
2214                0,  -1]
2215               .sliced(2, 2)
2216               .as!double.slice;
2217 
2218     auto C = qrDecomp(A);
2219     auto q = C.Q;
2220     auto r = C.R;
2221     auto A_reconstructed = C.reconstruct(q, r);
2222 
2223     import mir.math.common: approxEqual;
2224     import mir.algorithm.iteration: equal;
2225     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
2226     assert(appr(q, Q_check));
2227     assert(equal!approxEqual(r, R_check));
2228     assert(equal!approxEqual(A_reconstructed, A));
2229 }
2230 
2231 /++
2232 Computes a QR factorization of matrix 'a'.
2233 Params:
2234     allowDestroy = flag to delete the source matrix.
2235     a = initial matrix
2236 Returns: $(LREF QRResult)
2237 +/
2238 auto qrDecomp(Flag!"allowDestroy" allowDestroy = No.allowDestroy,
2239               Iterator, SliceKind kind)
2240              (Slice!(Iterator, 2, kind) a)
2241 {
2242     alias T = BlasType!Iterator;
2243     auto work = [T.sizeof * a.length].uninitSlice!T;
2244     auto tau = (cast(int) min(a.length!0, a.length!1)).uninitSlice!T;
2245     auto m = (allowDestroy && a._stride!1 == 1) ? a.assumeCanonical : a.transposed.as!T.slice.canonical;
2246 
2247     geqrf!T(m, tau, work);
2248     return QRResult!T(m, tau);
2249 }
2250 
2251 /++
2252 Solve the least squares problem:
2253     \min ||A * X - B||
2254 Using the QR factorization:
2255     \A = Q * R
2256 computed by qrDecomp.
2257 Params:
2258     allowDestroy = flag to delete the source matrix.
2259     a = detalis of the QR factorization of the original matrix as returned by qrDecomp.
2260     tau = details of the orhtogonal matrix Q.
2261     b = right hand side matrix.
2262 Returns: solution matrix.
2263 +/
2264 auto qrSolve(Flag!"allowDestroy" allowDestroy = No.allowDestroy,
2265              SliceKind kindB, size_t N, IteratorB, IteratorA, IteratorT)
2266             (Slice!(IteratorA, 2, Canonical) a,
2267              Slice!(IteratorT) tau,
2268              Slice!(IteratorB, N, kindB) b
2269             )
2270 in
2271 {
2272     assert(a.length!1 == b.length!0, "number of columns a should be equal to the number of rows b");
2273 }
2274 do
2275 {
2276     alias A = BlasType!IteratorA;
2277     alias B = BlasType!IteratorB;
2278     alias T = CommonType!(A, B);
2279     static if(is(T* == IteratorA))
2280         auto a_ = a;
2281     else
2282         auto a_ = a.as!T.slice.canonical;
2283     static if(is(T* == IteratorT))
2284         auto tau_ = tau;
2285     else
2286         auto tau_ = tau.as!T.slice.canonical;
2287 
2288     //convect vector to matrix.
2289     static if(N == 1)
2290         auto k = b.sliced(1, b.length);
2291     else
2292         auto k = b.transposed;
2293 
2294     static if(is(IteratorB == T*))
2295         auto m = (allowDestroy && k._stride!1 == 1) ? k.assumeCanonical : k.as!T.slice.canonical;
2296     else
2297         auto m = k.as!T.slice.canonical;
2298     auto work = [m.length!0].uninitSlice!T;
2299     static if(is(T == double) || is(T == float))
2300         ormqr!T('L', 'T', a_, tau_, m, work);
2301     else
2302         unmqr!T('L', 'C', a_, tau_, m, work);
2303 
2304     if (a_.length!0 != a_.length!1) {
2305         a_ = a_[0..tau.length, 0..tau.length];
2306         m = m.selectFront!1(tau.length);
2307     }
2308     trsm!T(Side.Right, Uplo.Lower, Diag.NonUnit, cast(T) 1.0, a_, m);
2309 
2310     return m.transposed;
2311 }
2312 
2313 ///
2314 unittest
2315 {
2316     import mir.ndslice;
2317 
2318     auto A =
2319             [ 3,  1, -1,  2,
2320              -5,  1,  3, -4,
2321               2,  0,  1, -1,
2322               1, -5,  3, -3 ]
2323               .sliced(4, 4)
2324               .as!double.slice;
2325 
2326     import mir.random.variable;
2327     import mir.random.algorithm;
2328     auto B = randomSlice(uniformVar(-100, 100), 4, 100);
2329 
2330     auto C = qrDecomp(A);
2331     auto X = C.solve(B);
2332 
2333     import mir.math.common: approxEqual;
2334     import mir.algorithm.iteration: equal;
2335     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
2336     assert(appr(mtimes(A, X), B));
2337 }
2338 
2339 unittest
2340 {
2341     auto A =
2342             [ 3,  1, -1,  2,
2343              -5,  1,  3, -4,
2344               2,  0,  1, -1,
2345               1, -5,  3, -3 ]
2346               .sliced(4, 4)
2347               .as!float.slice;
2348     auto B = [ 6, -12, 1, 3 ].sliced(4).as!double.slice.canonical;
2349     auto B_ = B.slice.sliced(B.length, 1);
2350     auto C = qrDecomp(A);
2351     auto X = qrSolve(C.matrix, C.tau, B);
2352 
2353     import mir.math.common: approxEqual;
2354     import mir.algorithm.iteration: equal;
2355     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
2356     assert(appr(mtimes(A, X), B_));
2357 }
2358 
2359 unittest
2360 {
2361     auto A =
2362             [ 1,  1,  0,
2363               1,  0,  1,
2364               0,  1,  1 ]
2365               .sliced(3, 3)
2366               .as!double.slice;
2367 
2368     auto B =
2369             [ 7,  6,  98,
2370               4,  8,  17,
2371               5,  3,  24 ]
2372               .sliced(3, 3)
2373               .as!float.slice;
2374     auto C = qrDecomp(A);
2375     auto X = qrSolve(C.matrix, C.tau, B);
2376 
2377     import mir.math.common: approxEqual;
2378     import mir.algorithm.iteration: equal;
2379     
2380     assert(equal!approxEqual(mtimes(A, X), B));
2381 }
2382 
2383 unittest
2384 {
2385     import mir.complex;
2386     auto A =
2387             [ 1,  1,  0,
2388               1,  0,  1,
2389               0,  1,  1 ]
2390               .sliced(3, 3)
2391               .as!(Complex!double).slice;
2392 
2393     auto B =
2394             [ 15,  78,  11,
2395               21,  47,  71,
2396               81,  11,  81 ]
2397               .sliced(3, 3)
2398               .as!(Complex!float).slice;
2399     auto C = qrDecomp(A);
2400     auto X = qrSolve(C.matrix, C.tau, B);
2401     auto res = mtimes(A, X);
2402 
2403     import mir.algorithm.iteration: equal;
2404     assert(equal!((a, b) => fabs(a - b) < 1e-12)(res, B));
2405 }
2406 
2407 unittest
2408 {
2409     auto A =
2410             [  3,  -6,
2411                4,  -8,
2412                0,   1]
2413               .sliced(3, 2)
2414               .as!double.slice;
2415 
2416     auto B = [-1, 7, 2]
2417                .sliced(3)
2418                .as!double.slice.canonical;
2419 
2420     auto X_check = [5, 2]
2421                     .sliced(2, 1)
2422                     .as!double.slice;
2423 
2424     auto C = qrDecomp(A);
2425     auto X = qrSolve(C.matrix, C.tau, B);
2426 
2427     import mir.math.common: approxEqual;
2428     import mir.algorithm.iteration: equal;
2429     alias appr = equal!((a, b) => approxEqual(a, b, 1e-5, 1e-5));
2430     assert(appr(X, X_check));
2431 }