1 /++ 2 $(H1 Lubeck 2 - `@nogc` Linear Algebra) 3 +/ 4 module kaleidic.lubeck2; 5 6 import mir.algorithm.iteration: equal, each; 7 import mir.blas; 8 import mir.exception; 9 import mir.internal.utility: isComplex, realType; 10 import mir.lapack; 11 import mir.math.common: sqrt, approxEqual, pow, fabs; 12 import mir.ndslice; 13 import mir.rc.array; 14 import mir.utility: min, max; 15 import std.traits: isFloatingPoint, Unqual; 16 import std.typecons: Flag, Yes, No; 17 import mir.complex: Complex; 18 19 /++ 20 Identity matrix. 21 22 Params: 23 n = number of columns 24 m = optional number of rows, default n 25 Results: 26 Matrix which is 1 on the diagonal and 0 elsewhere 27 +/ 28 @safe pure nothrow @nogc 29 Slice!(RCI!T, 2) eye(T = double)( 30 size_t n, 31 size_t m = 0 32 ) 33 if (isFloatingPoint!T || isComplex!T) 34 in 35 { 36 assert(n>0); 37 assert(m>=0); 38 } 39 out (i) 40 { 41 assert(i.length!0==n); 42 assert(i.length!1== (m==0 ? n : m)); 43 } 44 do 45 { 46 auto c = rcslice!T([n, (m==0?n:m)], cast(T)0); 47 c.diagonal[] = cast(T)1; 48 return c; 49 } 50 51 /// Real numbers 52 @safe pure nothrow 53 unittest 54 { 55 import mir.ndslice; 56 import mir.math; 57 58 assert(eye(1)== [ 59 [1]]); 60 assert(eye(2)== [ 61 [1, 0], 62 [0, 1]]); 63 assert(eye(3)== [ 64 [1, 0, 0], 65 [0, 1, 0], 66 [0, 0, 1]]); 67 assert(eye(1,2) == [ 68 [1,0]]); 69 assert(eye(2,1) == [ 70 [1], 71 [0]]); 72 } 73 74 /++ 75 General matrix-matrix multiplication. Allocates result to using Mir refcounted arrays. 76 Params: 77 a = m(rows) x k(cols) matrix 78 b = k(rows) x n(cols) matrix 79 Result: 80 m(rows) x n(cols) 81 +/ 82 @safe pure nothrow @nogc 83 Slice!(RCI!T, 2) mtimes(T, SliceKind kindA, SliceKind kindB)( 84 Slice!(const(T)*, 2, kindA) a, 85 Slice!(const(T)*, 2, kindB) b 86 ) 87 if (isFloatingPoint!T || isComplex!T) 88 in 89 { 90 assert(a.length!1 == b.length!0); 91 } 92 out (c) 93 { 94 assert(c.length!0 == a.length!0); 95 assert(c.length!1 == b.length!1); 96 } 97 do 98 { 99 // optimisations for spcecial cases can be added in the future 100 auto c = mininitRcslice!T(a.length!0, b.length!1); 101 gemm(cast(T)1, a, b, cast(T)0, c.lightScope); 102 return c; 103 } 104 105 /// ditto 106 @safe pure nothrow @nogc 107 Slice!(RCI!(Unqual!A), 2) mtimes(A, B, SliceKind kindA, SliceKind kindB)( 108 auto ref const Slice!(RCI!A, 2, kindA) a, 109 auto ref const Slice!(RCI!B, 2, kindB) b 110 ) 111 if (is(Unqual!A == Unqual!B)) 112 in 113 { 114 assert(a.length!1 == b.length!0); 115 } 116 do 117 { 118 auto scopeA = a.lightScope.lightConst; 119 auto scopeB = b.lightScope.lightConst; 120 return .mtimes(scopeA, scopeB); 121 } 122 123 @safe pure nothrow @nogc 124 Slice!(RCI!(Unqual!A), 2) mtimes(A, B, SliceKind kindA, SliceKind kindB)( 125 auto ref const Slice!(RCI!A, 2, kindA) a, 126 Slice!(const(B)*, 2, kindB) b 127 ) 128 if (is(Unqual!A == Unqual!B)) 129 in 130 { 131 assert(a.length!1 == b.length!0); 132 } 133 do 134 { 135 auto scopeA = a.lightScope.lightConst; 136 return .mtimes(scopeA, b); 137 } 138 139 @safe pure nothrow @nogc 140 Slice!(RCI!(Unqual!A), 2) mtimes(A, B, SliceKind kindA, SliceKind kindB)( 141 Slice!(const(A)*, 2, kindA) a, 142 auto ref const Slice!(RCI!B, 2, kindB) b 143 ) 144 if (is(Unqual!A == Unqual!B)) 145 in 146 { 147 assert(a.length!1 == b.length!0); 148 } 149 do 150 { 151 auto scopeB = b.lightScope.lightConst; 152 return .mtimes(a, scopeB); 153 } 154 /// Real numbers 155 @safe pure nothrow 156 unittest 157 { 158 import mir.ndslice; 159 import mir.math; 160 161 auto a = mininitRcslice!double(3, 5); 162 auto b = mininitRcslice!double(5, 4); 163 164 a[] = 165 [[-5, 1, 7, 7, -4], 166 [-1, -5, 6, 3, -3], 167 [-5, -2, -3, 6, 0]]; 168 169 b[] = 170 [[-5, -3, 3, 1], 171 [ 4, 3, 6, 4], 172 [-4, -2, -2, 2], 173 [-1, 9, 4, 8], 174 [ 9, 8, 3, -2]]; 175 176 assert(mtimes!(double, double)(a, b) == 177 [[-42, 35, -7, 77], 178 [-69, -21, -42, 21], 179 [ 23, 69, 3, 29]]); 180 } 181 182 /// Complex numbers 183 @safe pure nothrow 184 unittest 185 { 186 import mir.complex; 187 import mir.ndslice; 188 import mir.math; 189 190 auto a = mininitRcslice!(Complex!double)(3, 5); 191 auto b = mininitRcslice!(Complex!double)(5, 4); 192 193 a[] = 194 [[-5, 1, 7, 7, -4], 195 [-1, -5, 6, 3, -3], 196 [-5, -2, -3, 6, 0]]; 197 198 b[] = 199 [[-5, -3, 3, 1], 200 [ 4, 3, 6, 4], 201 [-4, -2, -2, 2], 202 [-1, 9, 4, 8], 203 [ 9, 8, 3, -2]]; 204 205 assert(mtimes!((Complex!double), (Complex!double))(a, b) == 206 [[-42, 35, -7, 77], 207 [-69, -21, -42, 21], 208 [ 23, 69, 3, 29]]); 209 } 210 211 /++ 212 Solve systems of linear equations AX = B for X. 213 Computes minimum-norm solution to a linear least squares problem 214 if A is not a square matrix. 215 +/ 216 @safe pure @nogc 217 Slice!(RCI!T, 2) mldivide (T, SliceKind kindA, SliceKind kindB)( 218 Slice!(const(T)*, 2, kindA) a, 219 Slice!(const(T)*, 2, kindB) b, 220 ) 221 if (isFloatingPoint!T || isComplex!T) 222 { 223 enforce!"mldivide: parameter shapes mismatch"(a.length!0 == b.length!0); 224 225 auto rcat = a.transposed.as!T.rcslice; 226 auto at = rcat.lightScope.canonical; 227 auto rcbt = b.transposed.as!T.rcslice; 228 auto bt = rcbt.lightScope.canonical; 229 size_t info; 230 if (a.length!0 == a.length!1) 231 { 232 auto rcipiv = at.length.mininitRcslice!lapackint; 233 auto ipiv = rcipiv.lightScope; 234 foreach(i; 0 .. ipiv.length) 235 ipiv[i] = 0; 236 info = gesv!T(at, ipiv, bt); 237 //info > 0 means some diagonla elem of U is 0 so no solution 238 } 239 else 240 { 241 static if(!isComplex!T) 242 { 243 size_t liwork = void; 244 auto lwork = gelsd_wq(at, bt, liwork); 245 auto rcs = min(at.length!0, at.length!1).mininitRcslice!T; 246 auto s = rcs.lightScope; 247 auto rcwork = lwork.rcslice!T; 248 auto work = rcwork.lightScope; 249 auto rciwork = liwork.rcslice!lapackint; 250 auto iwork = rciwork.lightScope; 251 size_t rank = void; 252 T rcond = -1; 253 254 info = gelsd!T(at, bt, s, rcond, rank, work, iwork); 255 //info > 0 means that many components failed to converge 256 } 257 else 258 { 259 size_t liwork = void; 260 size_t lrwork = void; 261 auto lwork = gelsd_wq(at, bt, lrwork, liwork); 262 auto rcs = min(at.length!0, at.length!1).mininitRcslice!(realType!T); 263 auto s = rcs.lightScope; 264 auto rcwork = lwork.rcslice!T; 265 auto work = rcwork.lightScope; 266 auto rciwork = liwork.rcslice!lapackint; 267 auto iwork = rciwork.lightScope; 268 auto rcrwork = lrwork.rcslice!(realType!T); 269 auto rwork = rcrwork.lightScope; 270 size_t rank = void; 271 realType!T rcond = -1; 272 273 info = gelsd!T(at, bt, s, rcond, rank, work, rwork, iwork); 274 //info > 0 means that many components failed to converge 275 } 276 bt = bt[0 .. $, 0 .. at.length!0]; 277 } 278 enforce!"mldivide: some off-diagonal elements of an intermediate bidiagonal form did not converge to zero."(!info); 279 return bt.transposed.as!T.rcslice; 280 } 281 282 /// ditto 283 @safe pure @nogc 284 Slice!(RCI!(Unqual!A), 2) 285 mldivide 286 (A, B, SliceKind kindA, SliceKind kindB)( 287 auto ref const Slice!(RCI!A, 2, kindA) a, 288 auto ref const Slice!(RCI!B, 2, kindB) b 289 ) 290 do 291 { 292 auto al = a.lightScope.lightConst; 293 auto bl = b.lightScope.lightConst; 294 return mldivide(al, bl); 295 } 296 297 /// ditto 298 @safe pure @nogc 299 Slice!(RCI!(Unqual!A), 1) 300 mldivide 301 (A, B, SliceKind kindA, SliceKind kindB)( 302 auto ref const Slice!(RCI!A, 2, kindA) a, 303 auto ref const Slice!(RCI!B, 1, kindB) b 304 ) 305 do 306 { 307 auto al = a.lightScope.lightConst; 308 auto bl = b.lightScope.lightConst; 309 return mldivide(al, bl); 310 } 311 312 /// ditto 313 @safe pure @nogc 314 Slice!(RCI!T, 1) mldivide (T, SliceKind kindA, SliceKind kindB)( 315 Slice!(const(T)*, 2, kindA) a, 316 Slice!(const(T)*, 1, kindB) b, 317 ) 318 if (isFloatingPoint!T || isComplex!T) 319 { 320 import mir.ndslice.topology: flattened; 321 return mldivide(a, b.sliced(b.length, 1)).flattened; 322 } 323 324 pure unittest 325 { 326 import mir.complex; 327 auto a = mininitRcslice!double(2, 2); 328 a[] = [[2,3], 329 [1, 4]]; 330 auto res = mldivide(eye(2), a); 331 assert(equal!approxEqual(res, a)); 332 auto b = mininitRcslice!(Complex!double)(2, 2); 333 b[] = [[Complex!double(2, 1),Complex!double(3, 2)], 334 [Complex!double(1, 3), Complex!double(4, 4)]]; 335 auto cres = mldivide(eye!(Complex!double)(2), b); 336 assert(cres == b); 337 auto c = mininitRcslice!double(2, 2); 338 c[] = [[5,3], 339 [2,6]]; 340 auto d = mininitRcslice!double(2,1); 341 d[] = [[4], 342 [1]]; 343 auto e = mininitRcslice!double(2,1); 344 e[] = [[23], 345 [14]]; 346 res = mldivide(c, e); 347 assert(equal!approxEqual(res, d)); 348 } 349 350 pure unittest 351 { 352 import mir.complex; 353 import mir.ndslice; 354 import mir.math; 355 356 auto a = mininitRcslice!double(6, 4); 357 a[] = [ 358 -0.57, -1.28, -0.39, 0.25, 359 -1.93, 1.08, -0.31, -2.14, 360 2.30, 0.24, 0.40, -0.35, 361 -1.93, 0.64, -0.66, 0.08, 362 0.15, 0.30, 0.15, -2.13, 363 -0.02, 1.03, -1.43, 0.50, 364 ].sliced(6, 4); 365 366 auto b = mininitRcslice!double(6,1); 367 b[] = [ 368 -2.67, 369 -0.55, 370 3.34, 371 -0.77, 372 0.48, 373 4.10, 374 ].sliced(6,1); 375 376 auto x = mininitRcslice!double(4,1); 377 x[] = [ 378 1.5339, 379 1.8707, 380 -1.5241, 381 0.0392 382 ].sliced(4,1); 383 384 auto res = mldivide(a, b); 385 assert(equal!((a, b) => fabs(a - b) < 5e-5)(res, x)); 386 387 auto ca = mininitRcslice!(Complex!double)(6, 4); 388 ca[] = [ 389 -0.57, -1.28, -0.39, 0.25, 390 -1.93, 1.08, -0.31, -2.14, 391 2.30, 0.24, 0.40, -0.35, 392 -1.93, 0.64, -0.66, 0.08, 393 0.15, 0.30, 0.15, -2.13, 394 -0.02, 1.03, -1.43, 0.50, 395 ].sliced(6, 4); 396 397 auto cb = mininitRcslice!(Complex!double)(6,1); 398 cb[] = [ 399 -2.67, 400 -0.55, 401 3.34, 402 -0.77, 403 0.48, 404 4.10, 405 ].sliced(6,1); 406 407 auto cx = mininitRcslice!(Complex!double)(4,1); 408 cx[] = [ 409 1.5339, 410 1.8707, 411 -1.5241, 412 0.0392 413 ].sliced(4,1); 414 415 auto cres = mldivide(ca, cb); 416 assert(equal!((a, b) => fabs(a - b) < 5e-5)(cres, x)); 417 } 418 419 /++ 420 Solve systems of linear equations AX = I for X, where I is the identity. 421 X is the right inverse of A if it exists, it's also a (Moore-Penrose) Pseudoinverse if A is invertible then X is the inverse. 422 Computes minimum-norm solution to a linear least squares problem 423 if A is not a square matrix. 424 +/ 425 @safe pure @nogc Slice!(RCI!A, 2) mlinverse(A, SliceKind kindA)( 426 auto ref Slice!(RCI!A, 2, kindA) a 427 ) 428 { 429 auto aScope = a.lightScope.lightConst; 430 return mlinverse!A(aScope); 431 } 432 433 @safe pure @nogc Slice!(RCI!A, 2) mlinverse(A, SliceKind kindA)( 434 Slice!(const(A)*, 2, kindA) a 435 ) 436 { 437 auto a_i = a.as!A.rcslice; 438 auto a_i_light = a_i.lightScope.canonical; 439 auto rcipiv = min(a_i.length!0, a_i.length!1).mininitRcslice!lapackint; 440 auto ipiv = rcipiv.lightScope; 441 auto info = getrf!A(a_i_light, ipiv); 442 if (info == 0) 443 { 444 auto rcwork = getri_wq!A(a_i_light).mininitRcslice!A; 445 auto work = rcwork.lightScope; 446 info = getri!A(a_i_light, ipiv, work); 447 } 448 enforce!"Matrix is not invertible as has zero determinant"(!info); 449 return a_i; 450 } 451 452 pure unittest 453 { 454 import mir.ndslice; 455 import mir.math; 456 457 auto a = mininitRcslice!double(2, 2); 458 a[] = [[1,0], 459 [0,-1]]; 460 auto ans = mlinverse!double(a); 461 assert(equal!approxEqual(ans, a)); 462 } 463 464 pure 465 unittest 466 { 467 import mir.ndslice; 468 import mir.math; 469 470 auto a = mininitRcslice!double(2, 2); 471 a[] = [[ 0, 1], 472 [-1, 0]]; 473 auto aInv = mininitRcslice!double(2, 2); 474 aInv[] = [[0, -1], 475 [1, 0]]; 476 auto ans = a.mlinverse; 477 assert(equal!approxEqual(ans, aInv)); 478 } 479 480 ///Singuar value decomposition result 481 struct SvdResult(T) 482 { 483 /// 484 Slice!(RCI!T, 2) u; 485 ///Singular Values 486 Slice!(RCI!T) sigma; 487 /// 488 Slice!(RCI!T, 2) vt; 489 } 490 491 /++ 492 Computes the singular value decomposition. 493 Params: 494 a = input `M x N` matrix 495 slim = If true the first `min(M,N)` columns of `u` and the first 496 `min(M,N)` rows of `vt` are returned in the ndslices `u` and `vt`. 497 out result = $(LREF SvdResult). ] 498 Returns: error code from CBlas 499 +/ 500 @safe pure @nogc SvdResult!T svd 501 ( 502 T, 503 string algorithm = "gesvd", 504 SliceKind kind, 505 )( 506 Slice!(const(T)*, 2, kind) a, 507 Flag!"slim" slim = No.slim, 508 ) 509 if (algorithm == "gesvd" || algorithm == "gesdd") 510 { 511 auto m = cast(lapackint)a.length!1; 512 auto n = cast(lapackint)a.length!0; 513 514 auto s = mininitRcslice!T(min(m, n)); 515 auto u = mininitRcslice!T(slim ? s.length : m, m); 516 auto vt = mininitRcslice!T(n, slim ? s.length : n); 517 518 if (m == 0 || n == 0) 519 { 520 u.lightScope[] = 0; 521 u.lightScope.diagonal[] = 1; 522 vt.lightScope[] = 0; 523 vt.lightScope.diagonal[] = 1; 524 } 525 else 526 { 527 static if (algorithm == "gesvd") 528 { 529 auto jobu = slim ? 'S' : 'A'; 530 auto jobvt = slim ? 'S' : 'A'; 531 auto rca_sliced = a.as!T.rcslice; 532 auto rca = rca_sliced.lightScope.canonical; 533 auto rcwork = gesvd_wq(jobu, jobvt, rca, u.lightScope.canonical, vt.lightScope.canonical).mininitRcslice!T; 534 auto work = rcwork.lightScope; 535 auto info = gesvd(jobu, jobvt, rca, s.lightScope, u.lightScope.canonical, vt.lightScope.canonical, work); 536 } 537 else // gesdd 538 { 539 auto rciwork = mininitRcslice!lapackint(s.length * 8); 540 auto iwork = rciwork.lightScope; 541 auto jobz = slim ? 'S' : 'A'; 542 auto rca_sliced = a.as!T.rcslice; 543 auto rca = rca_sliced.lightScope.canonical; 544 auto rcwork = gesdd_wq(jobz, rca, u.lightScope, vt.lightScope).minitRcslice!T; 545 auto work = rcwork.lightScope; 546 auto info = gesdd(jobz, rca, s.lightScope, u.lightScope, vt.lightScope, work, iwork); 547 } 548 enum msg = (algorithm == "gesvd" ? "svd: DBDSDC did not converge, updating process failed" : "svd: DBDSQR did not converge"); 549 enforce!("svd: " ~ msg)(!info); 550 } 551 return SvdResult!T(vt, s, u); //transposed 552 } 553 554 @safe pure SvdResult!T svd 555 ( 556 T, 557 string algorithm = "gesvd", 558 SliceKind kind, 559 )( 560 auto ref scope const Slice!(RCI!T,2,kind) matrix, 561 Flag!"slim" slim = No.slim 562 ) 563 { 564 auto matrixScope = matrix.lightScope.lightConst; 565 return svd!(T, algorithm)(matrixScope, slim); 566 } 567 568 pure unittest 569 { 570 import mir.ndslice; 571 import mir.math; 572 573 auto a = mininitRcslice!double(6, 4); 574 a[] = [[7.52, -1.10, -7.95, 1.08], 575 [-0.76, 0.62, 9.34, -7.10], 576 [5.13, 6.62, -5.66, 0.87], 577 [-4.75, 8.52, 5.75, 5.30], 578 [1.33, 4.91, -5.49, -3.52], 579 [-2.40, -6.77, 2.34, 3.95]]; 580 581 auto r1 = a.svd; 582 583 auto sigma1 = rcslice!double(a.shape, 0); 584 sigma1.diagonal[] = r1.sigma; 585 auto m1 = (r1.u).mtimes(sigma1).mtimes(r1.vt); 586 assert(equal!((a, b) => fabs(a-b)< 1e-8)(a, m1)); 587 588 auto r2 = a.svd; 589 590 auto sigma2 = rcslice!double(a.shape, 0.0); 591 sigma2.diagonal[] = r2.sigma; 592 auto m2 = r2.u.mtimes(sigma2).mtimes(r2.vt); 593 594 assert(equal!((a, b) => fabs(a-b)< 1e-8)(a, m2)); 595 596 auto r = a.svd(Yes.slim); 597 assert(r.u.shape == [6, 4]); 598 assert(r.vt.shape == [4, 4]); 599 } 600 601 pure unittest 602 { 603 import mir.ndslice; 604 import mir.math; 605 606 auto a = mininitRcslice!double(6, 4); 607 a[] = [[7.52, -1.10, -7.95, 1.08], 608 [-0.76, 0.62, 9.34, -7.10], 609 [5.13, 6.62, -5.66, 0.87], 610 [-4.75, 8.52, 5.75, 5.30], 611 [1.33, 4.91, -5.49, -3.52], 612 [-2.40, -6.77, 2.34, 3.95]]; 613 614 auto r1 = a.svd(No.slim); 615 616 auto sigma1 = rcslice!double(a.shape, 0.0); 617 sigma1.diagonal[] = r1.sigma; 618 auto m1 = r1.u.mtimes(sigma1).mtimes(r1.vt); 619 620 assert(equal!((a, b) => fabs(a-b)< 1e-8)(a, m1)); 621 } 622 623 unittest 624 { 625 import mir.algorithm.iteration: all; 626 627 // empty matrix as input means that u or vt is identity matrix 628 auto identity = slice!double([4, 4], 0); 629 identity.diagonal[] = 1; 630 631 auto a = slice!double(0, 4); 632 auto res = a.svd; 633 634 import mir.conv: to; 635 636 assert(res.u.shape == [0, 0]); 637 assert(res.vt.shape == [4, 4]); 638 assert(res.vt.all!approxEqual(identity), res.vt.to!string); 639 640 auto b = slice!double(4, 0); 641 res = b.svd; 642 643 assert(res.u.shape == [4, 4]); 644 assert(res.vt.shape == [0, 0]); 645 646 assert(res.u.all!approxEqual(identity), res.u.to!string); 647 } 648 649 struct EigenResult(T) 650 { 651 Slice!(RCI!(complexType!T), 2) eigenvectors; 652 Slice!(RCI!(complexType!T)) eigenvalues; 653 } 654 655 struct QRResult(T) 656 { 657 Slice!(RCI!T,2) Q; 658 Slice!(RCI!T,2) R; 659 660 @safe this 661 ( 662 SliceKind kindA, 663 SliceKind kindTau 664 )( 665 Slice!(RCI!T, 2, kindA) a, 666 Slice!(RCI!T, 1, kindTau) tau 667 ) 668 { 669 auto aScope = a.lightScope.lightConst; 670 auto tauScope = tau.lightScope.lightConst; 671 this(aScope, tauScope); 672 } 673 674 @safe this 675 ( 676 SliceKind kindA, 677 SliceKind kindTau 678 )( 679 Slice!(const(T)*, 2, kindA) a, 680 Slice!(const(T)*, 1, kindTau) tau 681 ) 682 { 683 R = mininitRcslice!T(a.shape); 684 foreach (i; 0 .. R.length!0) 685 { 686 foreach (j; 0 .. R.length!1) 687 { 688 if (i >= j) 689 { 690 R[j, i] = a[i, j]; 691 } 692 else 693 { 694 R[j ,i] = cast(T)0; 695 } 696 } 697 } 698 auto rcwork = mininitRcslice!T(a.length!0); 699 auto work = rcwork.lightScope; 700 auto aSliced = a.as!T.rcslice; 701 auto aScope = aSliced.lightScope.canonical; 702 auto tauSliced = tau.as!T.rcslice; 703 auto tauScope = tauSliced.lightScope; 704 705 static if(is(T == double) || is(T == float)) 706 orgqr!T(aScope, tauScope, work); 707 else 708 ungqr!T(aScope, tauScope, work); 709 710 Q = aScope.transposed.as!T.rcslice; 711 } 712 } 713 714 @safe pure QRResult!T qr(T, SliceKind kind)( 715 auto ref const Slice!(RCI!T, 2, kind) matrix 716 ) 717 { 718 auto matrixScope = matrix.lightScope.lightConst; 719 return qr!(T, kind)(matrixScope); 720 } 721 722 @safe pure QRResult!T qr(T, SliceKind kind)( 723 auto ref const Slice!(const(T)*, 2, kind) matrix 724 ) 725 { 726 auto a = matrix.transposed.as!T.rcslice; 727 auto tau = mininitRcslice!T(a.length!0); 728 auto rcwork = mininitRcslice!T(a.length!0); 729 auto work = rcwork.lightScope; 730 auto aScope = a.lightScope.canonical; 731 auto tauScope = tau.lightScope; 732 geqrf!T(aScope, tauScope, work); 733 return QRResult!T(aScope, tauScope); 734 } 735 736 pure nothrow 737 unittest 738 { 739 import mir.ndslice; 740 import mir.math; 741 742 auto data = mininitRcslice!double(3, 3); 743 data[] = [[12, -51, 4], 744 [ 6, 167, -68], 745 [-4, 24, -41]]; 746 747 auto res = qr(data); 748 auto q = mininitRcslice!double(3, 3); 749 q[] = [[-6.0/7.0, 69.0/175.0, 58.0/175.0], 750 [-3.0/7.0, -158.0/175.0, -6.0/175.0], 751 [ 2.0/7.0, -6.0/35.0 , 33.0/35.0 ]]; 752 auto aE = function (double x, double y) => approxEqual(x, y, 0.00005, 0.00005); 753 auto r = mininitRcslice!double(3, 3); 754 r[] = [[-14, -21, 14], 755 [ 0, -175, 70], 756 [ 0, 0, -35]]; 757 assert(equal!approxEqual(mtimes(res.Q, res.R), data)); 758 } 759 760 pure nothrow 761 unittest 762 { 763 import mir.ndslice; 764 import mir.math; 765 766 auto data = mininitRcslice!(Complex!double)(3, 3); 767 data[] = [[12, -51, 4], 768 [ 6, 167, -68], 769 [-4, 24, -41]]; 770 771 auto res = qr(data); 772 auto q = mininitRcslice!(Complex!double)(3, 3); 773 q[] = [[-6.0/7.0, 69.0/175.0, 58.0/175.0], 774 [-3.0/7.0, -158.0/175.0, -6.0/175.0], 775 [ 2.0/7.0, -6.0/35.0 , 33.0/35.0 ]]; 776 auto aE = function (Complex!double x, Complex!double y) => approxEqual(x, y, 0.00005, 0.00005); 777 auto r = mininitRcslice!(Complex!double)(3, 3); 778 r[] = [[-14, -21, 14], 779 [ 0, -175, 70], 780 [ 0, 0, -35]]; 781 assert(equal!approxEqual(mtimes(res.Q, res.R), data)); 782 } 783 784 785 @safe pure @nogc 786 EigenResult!(realType!T) eigen(T, SliceKind kind)( 787 auto ref Slice!(const(T)*,2, kind) a, 788 ) 789 { 790 enforce!"eigen: input matrix must be square"(a.length!0 == a.length!1); 791 const n = a.length; 792 auto rcw = n.mininitRcslice!(complexType!T); 793 auto w = rcw.lightScope; 794 auto rcwork = mininitRcslice!T(16 * n); 795 auto work = rcwork.lightScope; 796 auto z = [n, n].mininitRcslice!(complexType!T);//vl (left eigenvectors) 797 auto rca = a.transposed.as!T.rcslice; 798 auto as = rca.lightScope; 799 static if (isComplex!T) 800 { 801 auto rcrwork = [2 * n].mininitRcslice!(realType!T); 802 auto rwork = rcrwork.lightScope; 803 auto info = geev!(T, realType!T)('N', 'V', as.canonical, w, Slice!(T*, 2, Canonical).init, z.lightScope.canonical, work, rwork); 804 enforce!"eigen failed"(!info); 805 } 806 else 807 { 808 alias C = complexType!T; 809 auto wr = sliced((cast(T[]) w.field)[0 .. n]); 810 auto wi = sliced((cast(T[]) w.field)[n .. n * 2]); 811 auto rczr = [n, n].mininitRcslice!T; 812 auto zr = rczr.lightScope; 813 auto info = geev!T('N', 'V', as.canonical, wr, wi, Slice!(T*, 2, Canonical).init, zr.canonical, work); 814 enforce!"eigen failed"(!info); 815 work[0 .. n] = wr; 816 work[n .. n * 2] = wi; 817 foreach (i, ref e; w.field) 818 { 819 e = C(work[i], work[n + i]); 820 auto zi = z.lightScope[i]; 821 if (e.im > 0) 822 zi[] = zip(zr[i], zr[i + 1]).map!((a, b) => C(a, b)); 823 else 824 if (e.im < 0) 825 zi[] = zip(zr[i - 1], zr[i]).map!((a, b) => C(a, -b)); 826 else 827 zi[] = zr[i].as!C; 828 } 829 } 830 return typeof(return)(z, rcw); 831 } 832 833 //@safe pure @nogc 834 EigenResult!(realType!T) eigen(T, SliceKind kind)( 835 auto ref Slice!(RCI!T,2, kind) a 836 ) 837 { 838 auto as = a.lightScope.lightConst; 839 return eigen(as); 840 } 841 842 /// 843 // pure 844 unittest 845 { 846 import mir.blas; 847 import mir.complex; 848 import mir.ndslice; 849 import mir.math; 850 851 auto data = 852 [[ 0, 1], 853 [-1, 0]].fuse.as!double.rcslice; 854 auto eigenvalues = [Complex!double(0, 1), Complex!double(0, -1)].sliced; 855 auto eigenvectors = 856 [[Complex!double(0, -1), Complex!double(1, 0)], 857 [Complex!double(0, 1), Complex!double(1, 0)]].fuse; 858 859 auto res = data.eigen; 860 861 assert(res.eigenvalues.equal!approxEqual(eigenvalues)); 862 foreach (i; 0 .. eigenvectors.length) 863 assert((res.eigenvectors.lightScope[i] / eigenvectors[i]).diff.slice.nrm2.approxEqual(0)); 864 } 865 866 867 /// 868 @safe pure 869 unittest 870 { 871 import mir.complex; 872 import mir.ndslice; 873 import mir.math; 874 import mir.blas; 875 876 auto data = 877 [[0, 1, 0], 878 [0, 0, 1], 879 [1, 0, 0]].fuse.as!double.rcslice; 880 auto c = 3.0.sqrt; 881 auto eigenvalues = [Complex!double(-1, c) / 2, Complex!double(-1, - c) / 2, Complex!double(1, 0)]; 882 883 auto eigenvectors = 884 [[Complex!double(-1, +c), Complex!double(-1, -c) , Complex!double(2, 0)], 885 [Complex!double(-1, -c), Complex!double(-1, +c) , Complex!double(2, 0)], 886 [Complex!double(1, 0), Complex!double(1, 0), Complex!double(1, 0)]].fuse; 887 888 auto res = data.eigen; 889 890 assert(res.eigenvalues.equal!approxEqual(eigenvalues)); 891 foreach (i; 0 .. eigenvectors.length) 892 assert((res.eigenvectors.lightScope[i] / eigenvectors[i]).diff.slice.nrm2.approxEqual(0)); 893 894 auto cdata = data.lightScope.as!(Complex!double).rcslice; 895 res = cdata.eigen; 896 897 assert(res.eigenvalues.equal!approxEqual(eigenvalues)); 898 foreach (i; 0 .. eigenvectors.length) 899 assert((res.eigenvectors.lightScope[i] / eigenvectors[i]).diff.slice.nrm2.approxEqual(0)); 900 } 901 902 903 /// Principal component analysis result. 904 struct PcaResult(T) 905 { 906 /// Eigenvectors (Eigenvectors[i] is the ith eigenvector) 907 Slice!(RCI!T,2) eigenvectors; 908 /// Principal component scores. (Input matrix rotated to basis of eigenvectors) 909 Slice!(RCI!T, 2) scores; 910 /// Principal component variances. (Eigenvalues) 911 Slice!(RCI!T) eigenvalues; 912 /// The means of each data column (0 if not centred) 913 Slice!(RCI!T) mean; 914 /// The standard deviations of each column (1 if not normalized) 915 Slice!(RCI!T) stdDev; 916 /// Principal component Loadings vectors (Eigenvectors times sqrt(Eigenvalue)) 917 Slice!(RCI!T, 2) loadings() 918 { 919 auto result = mininitRcslice!T(eigenvectors.shape); 920 for (size_t i = 0; i < eigenvectors.length!0 && i < eigenvalues.length; i++){ 921 if(eigenvalues[i] != 0) 922 foreach (j; 0 .. eigenvectors.length!1) 923 result[i, j] = eigenvectors[i, j] * sqrt(eigenvalues[i]); 924 else 925 result[i][] = eigenvectors[i][]; 926 } 927 return result; 928 } 929 930 Slice!(RCI!T, 2) loadingsScores() {// normalized data in basis of {sqrt(eval)*evect} 931 //if row i is a vector p, the original data point is mean[] + p[] * stdDev[] 932 auto result = mininitRcslice!T(scores.shape); 933 foreach (i; 0 .. scores.length!0){ 934 for (size_t j=0; j < scores.length!1 && j < eigenvalues.length; j++) 935 { 936 if(eigenvalues[j] != 0) 937 result[i, j] = scores[i, j] / sqrt(eigenvalues[j]); 938 else 939 result[i, j] = scores[i, j]; 940 } 941 } 942 return result; 943 } 944 945 Slice!(RCI!T) explainedVariance() { 946 import mir.math.sum: sum; 947 //auto totalVariance = eigenvalues.sum!"kb2"; 948 auto totalVariance = 0.0; 949 foreach (val; eigenvalues) 950 totalVariance += val; 951 auto result = mininitRcslice!double(eigenvalues.shape); 952 foreach (i; 0 .. result.length!0) 953 result[i] = eigenvalues[i] / totalVariance; 954 return result; 955 } 956 957 Slice!(RCI!T,2) q() 958 { 959 return eigenvectors; 960 } 961 962 //returns matrix to transform into basis of 'n' most significatn eigenvectors 963 Slice!(RCI!(T),2) q(size_t n) 964 in { 965 assert(n <= eigenvectors.length!0); 966 } 967 do { 968 auto res = mininitRcslice!T(eigenvectors.shape); 969 res[] = eigenvectors; 970 for ( ; n < eigenvectors.length!0; n++) 971 res[n][] = 0; 972 return res; 973 } 974 //returns matrix to transform into basis of eigenvectors with value larger than delta 975 Slice!(RCI!(T),2) q(T delta) 976 do 977 { 978 auto res = mininitRcslice!T(eigenvectors.shape); 979 res[] = eigenvectors; 980 foreach (i; 0 .. eigenvectors.length!0) 981 if (fabs(eigenvalues[i]) <= delta) 982 res[i][] = 0; 983 return res; 984 } 985 986 //transforms data into eigenbasis 987 Slice!(RCI!(T),2) transform(Slice!(RCI!(T),2) data) 988 { 989 return q.mlinverse.mtimes(data); 990 } 991 992 //transforms data into eigenbasis 993 Slice!(RCI!(T),2) transform(Slice!(RCI!(T),2) data, size_t n) 994 { 995 return q(n).mlinverse.mtimes(data); 996 } 997 //transforms data into eigenbasis 998 Slice!(RCI!(T),2) transform(Slice!(RCI!(T),2) data, T delta) 999 { 1000 return q(delta).mlinverse.mtimes(data); 1001 } 1002 } 1003 1004 1005 /++ 1006 Principal component analysis of raw data. 1007 Template: 1008 correlation = Flag to use correlation matrix instead of covariance 1009 Params: 1010 data = input `M x N` matrix, where 'M (rows)>= N(cols)' 1011 devEst = 1012 meanEst = 1013 fixEigenvectorDirections = 1014 Returns: $(LREF PcaResult) 1015 +/ 1016 @safe pure @nogc 1017 PcaResult!T pca( 1018 SliceKind kind, 1019 T 1020 )( 1021 Slice!(const(T)*, 2, kind) data, 1022 DeviationEstimator devEst = DeviationEstimator.sample, 1023 MeanEstimator meanEst = MeanEstimator.average, 1024 Flag!"fixEigenvectorDirections" fixEigenvectorDirections = Yes.fixEigenvectorDirections, 1025 ) 1026 in 1027 { 1028 assert(data.length!0 >= data.length!1); 1029 } 1030 do 1031 { 1032 real n = (data.length!0 <= 1 ? 1 : data.length!0 -1 );//num observations 1033 auto mean = rcslice!(T,1)([data.length!1], cast(T)0); 1034 auto stdDev = rcslice!(T,1)([data.length!1], cast(T)1); 1035 auto centeredData = centerColumns!T(data, mean, meanEst); 1036 //this part gets the eigenvectors of the sample covariance without explicitly calculating the covariance matrix 1037 //to implement a minimum covariance deteriminant this block would need to be redone 1038 //firstly one would calculate the MCD then call eigen to get it's eigenvectors 1039 auto processedData = normalizeColumns!T(centeredData, stdDev, devEst); 1040 auto svdResult = processedData.svd(Yes.slim); 1041 with (svdResult) 1042 { 1043 //u[i][] is the ith eigenvector 1044 foreach (i; 0 .. u.length!0){ 1045 foreach (j; 0 .. u.length!1){ 1046 u[i, j] *= sigma[j]; 1047 } 1048 } 1049 auto eigenvalues = mininitRcslice!double(sigma.shape); 1050 for (size_t i = 0; i < sigma.length && i < eigenvalues.length; i++){ 1051 eigenvalues[i] = sigma[i] * sigma[i] / n; //square singular values to get eigenvalues 1052 } 1053 if (fixEigenvectorDirections) 1054 { 1055 foreach(size_t i; 0 .. sigma.length) 1056 { 1057 //these are directed so the 0th component is +ve 1058 if (vt[i, 0] < 0) 1059 { 1060 vt[i][] *= -1; 1061 u[0 .. $, i] *= -1; 1062 } 1063 } 1064 } 1065 PcaResult!T result; 1066 result.scores = u; 1067 result.eigenvectors = vt; 1068 result.eigenvalues = eigenvalues; 1069 result.mean = mean; 1070 result.stdDev = stdDev; 1071 return result; 1072 } 1073 } 1074 1075 @safe pure @nogc 1076 PcaResult!T pca(T, SliceKind kind) 1077 ( 1078 auto ref const Slice!(RCI!T, 2, kind) data, 1079 DeviationEstimator devEst = DeviationEstimator.sample, 1080 MeanEstimator meanEst = MeanEstimator.average, 1081 Flag!"fixEigenvectorDirections" fixEigenvectorDirections = Yes.fixEigenvectorDirections, 1082 ) 1083 do 1084 { 1085 auto d = data.lightScope.lightConst; 1086 return d.pca(devEst, meanEst, fixEigenvectorDirections); 1087 } 1088 1089 pure 1090 unittest 1091 { 1092 import mir.ndslice; 1093 import mir.math; 1094 1095 auto data = mininitRcslice!double(3, 2); 1096 data[] = [[ 1, -1], 1097 [ 0, 1], 1098 [-1, 0]]; 1099 //cov =0.5 * [[ 2, -1], 1100 // [-1, 2]] 1101 const auto const_data = data; 1102 auto mean = mininitRcslice!double(2); 1103 assert(data == centerColumns(const_data, mean)); 1104 assert(mean == [0,0]); 1105 PcaResult!double res = const_data.pca; 1106 1107 auto evs = mininitRcslice!double(2, 2); 1108 evs[] = [[1, -1], 1109 [1, 1]]; 1110 evs[] /= sqrt(2.0); 1111 assert(equal!approxEqual(res.eigenvectors, evs)); 1112 1113 auto score = mininitRcslice!double(3, 2); 1114 score[] = [[ 1.0, 0.0], 1115 [-0.5, 0.5], 1116 [-0.5, -0.5]]; 1117 score[] *= sqrt(2.0); 1118 assert(equal!approxEqual(res.scores, score)); 1119 1120 auto evals = mininitRcslice!double(2); 1121 evals[] = [1.5, 0.5]; 1122 assert(equal!approxEqual(res.eigenvalues, evals)); 1123 } 1124 1125 pure 1126 unittest 1127 { 1128 import mir.ndslice; 1129 import mir.math; 1130 1131 auto data = mininitRcslice!double(10, 3); 1132 data[] = [[7, 4, 3], 1133 [4, 1, 8], 1134 [6, 3, 5], 1135 [8, 6, 1], 1136 [8, 5, 7], 1137 [7, 2, 9], 1138 [5, 3, 3], 1139 [9, 5, 8], 1140 [7, 4, 5], 1141 [8, 2, 2]]; 1142 1143 auto m1 = 69.0/10.0; 1144 auto m2 = 35.0/10.0; 1145 auto m3 = 51.0/10.0; 1146 1147 auto centeredData = mininitRcslice!double(10, 3); 1148 centeredData[] = [[7-m1, 4-m2, 3-m3], 1149 [4-m1, 1-m2, 8-m3], 1150 [6-m1, 3-m2, 5-m3], 1151 [8-m1, 6-m2, 1-m3], 1152 [8-m1, 5-m2, 7-m3], 1153 [7-m1, 2-m2, 9-m3], 1154 [5-m1, 3-m2, 3-m3], 1155 [9-m1, 5-m2, 8-m3], 1156 [7-m1, 4-m2, 5-m3], 1157 [8-m1, 2-m2, 2-m3]]; 1158 auto mean = mininitRcslice!double(3); 1159 auto cenRes = centerColumns(data, mean); 1160 1161 assert(equal!approxEqual(centeredData, cenRes)); 1162 assert(equal!approxEqual(mean, [m1,m2,m3])); 1163 1164 auto res = data.pca; 1165 auto coeff = mininitRcslice!double(3, 3); 1166 coeff[] = [[0.6420046 , 0.6863616 , -0.3416692 ], 1167 [0.38467229, 0.09713033, 0.91792861], 1168 [0.6632174 , -0.7207450 , -0.2016662 ]]; 1169 1170 auto score = mininitRcslice!double(10, 3); 1171 score[] = [[ 0.5148128, -0.63083556, -0.03351152], 1172 [-2.6600105, 0.06280922, -0.33089322], 1173 [-0.5840389, -0.29060575, -0.15658900], 1174 [ 2.0477577, -0.90963475, -0.36627323], 1175 [ 0.8832739, 0.99120250, -0.34153847], 1176 [-1.0837642, 1.20857108, 0.44706241], 1177 [-0.7618703, -1.19712391, -0.44810271], 1178 [ 1.1828371, 1.57067601, 0.02182598], 1179 [ 0.2713493, 0.02325373, -0.17721300], 1180 [ 0.1896531, -0.82831257, 1.38523276]]; 1181 1182 auto stdDev = mininitRcslice!double(3); 1183 stdDev[] = [1.3299527, 0.9628478, 0.5514979]; 1184 1185 auto eigenvalues = mininitRcslice!double(3); 1186 eigenvalues[] = [1.768774, 0.9270759, 0.3041499]; 1187 1188 assert(equal!approxEqual(res.eigenvectors, coeff)); 1189 assert(equal!approxEqual(res.scores, score)); 1190 assert(equal!approxEqual(res.eigenvalues, eigenvalues)); 1191 } 1192 1193 pure 1194 unittest 1195 { 1196 import mir.ndslice; 1197 import mir.math; 1198 1199 auto data = mininitRcslice!double(13, 4); 1200 data[] =[[ 7, 26, 6, 60], 1201 [ 1, 29, 15, 52], 1202 [11, 56, 8, 20], 1203 [11, 31, 8, 47], 1204 [ 7, 52, 6, 33], 1205 [11, 55, 9, 22], 1206 [ 3, 71, 17, 6], 1207 [ 1, 31, 22, 44], 1208 [ 2, 54, 18, 22], 1209 [21, 47, 4, 26], 1210 [ 1, 40, 23, 34], 1211 [11, 66, 9, 12], 1212 [10, 68, 8, 12]]; 1213 1214 auto m1 = 97.0/13.0; 1215 auto m2 = 626.0/13.0; 1216 auto m3 = 153.0/13.0; 1217 auto m4 = 390.0/13.0; 1218 1219 auto centeredData = mininitRcslice!double(13, 4); 1220 centeredData[] = [[ 7-m1, 26-m2, 6-m3, 60-m4], 1221 [ 1-m1, 29-m2, 15-m3, 52-m4], 1222 [11-m1, 56-m2, 8-m3, 20-m4], 1223 [11-m1, 31-m2, 8-m3, 47-m4], 1224 [ 7-m1, 52-m2, 6-m3, 33-m4], 1225 [11-m1, 55-m2, 9-m3, 22-m4], 1226 [ 3-m1, 71-m2, 17-m3, 6-m4], 1227 [ 1-m1, 31-m2, 22-m3, 44-m4], 1228 [ 2-m1, 54-m2, 18-m3, 22-m4], 1229 [21-m1, 47-m2, 4-m3, 26-m4], 1230 [ 1-m1, 40-m2, 23-m3, 34-m4], 1231 [11-m1, 66-m2, 9-m3, 12-m4], 1232 [10-m1, 68-m2 ,8-m3, 12-m4]]; 1233 auto mean = mininitRcslice!double(4); 1234 auto cenRes = centerColumns(data, mean); 1235 1236 assert(equal!approxEqual(centeredData, cenRes)); 1237 assert(mean == [m1,m2,m3,m4]); 1238 1239 auto res = data.pca(DeviationEstimator.none); 1240 auto coeff = mininitRcslice!double(4, 4); 1241 coeff[] = [[0.067799985695474, 0.678516235418647, -0.029020832106229, -0.730873909451461], 1242 [0.646018286568728, 0.019993340484099, -0.755309622491133, 0.108480477171676], 1243 [0.567314540990512, -0.543969276583817, 0.403553469172668, -0.468397518388289], 1244 [0.506179559977705, 0.493268092159297, 0.515567418476836, 0.484416225289198]]; 1245 1246 auto score = mininitRcslice!double(13, 4); 1247 score[] = [[-36.821825999449700, 6.870878154227367, -4.590944457629745, 0.396652582713912], 1248 [-29.607273420710964, -4.610881963526308, -2.247578163663940, -0.395843536696492], 1249 [ 12.981775719737618, 4.204913183175938, 0.902243082694698, -1.126100587210615], 1250 [-23.714725720918022, 6.634052554708721, 1.854742000806314, -0.378564808384691], 1251 [ 0.553191676624597, 4.461732123178686, -6.087412652325177, 0.142384896047281], 1252 [ 10.812490833309816, 3.646571174544059, 0.912970791674604, -0.134968810314680], 1253 [ 32.588166608817929, -8.979846284936063, -1.606265913996588, 0.081763927599947], 1254 [-22.606395499005586, -10.725906457369449, 3.236537714483416, 0.324334774646368], 1255 [ 9.262587237675838, -8.985373347478788, -0.016909578102172, -0.543746175981799], 1256 [ 3.283969329640680, 14.157277337500918, 7.046512994833761, 0.340509860960606], 1257 [ -9.220031117829379, -12.386080787220454, 3.428342878284624, 0.435152769664895], 1258 [ 25.584908517429557, 2.781693148152386, -0.386716066864491, 0.446817950545605], 1259 [ 26.903161834677597, 2.930971165042989, -2.445522630195304, 0.411607156409658]]; 1260 1261 auto eigenvalues = mininitRcslice!double(4); 1262 1263 1264 eigenvalues[] = [517.7968780739053, 67.4964360487231, 12.4054300480810, 0.2371532651878]; 1265 1266 assert(equal!approxEqual(res.eigenvectors, coeff)); 1267 assert(equal!approxEqual(res.scores, score)); 1268 assert(equal!approxEqual(res.eigenvalues, eigenvalues)); 1269 } 1270 1271 ///complex extensions 1272 1273 private T conj(T)( 1274 const T z 1275 ) 1276 if (isComplex!T) 1277 { 1278 return z.re - (1fi* z.im); 1279 } 1280 1281 private template complexType(C) 1282 { 1283 import mir.complex: Complex; 1284 static if (isComplex!C) 1285 alias complexType = Unqual!C; 1286 else static if (is(Unqual!C == double)) 1287 alias complexType = Complex!(Unqual!C); 1288 } 1289 1290 /// 1291 enum MeanEstimator 1292 { 1293 /// 1294 none, 1295 /// 1296 average, 1297 /// 1298 median 1299 } 1300 1301 /// 1302 @safe pure nothrow @nogc 1303 T median(T)(auto ref Slice!(RCI!T) data) 1304 { 1305 auto dataScope = data.lightScope.lightConst; 1306 return median!T(dataScope); 1307 } 1308 1309 /// ditto 1310 @safe pure nothrow @nogc 1311 T median(T)(Slice!(const(T)*) data) 1312 { 1313 import mir.ndslice.sorting: sort; 1314 size_t len = cast(int) data.length; 1315 size_t n = len / 2; 1316 auto temp = data.as!T.rcslice; 1317 temp.lightScope.sort(); 1318 return len % 2 ? temp[n] : 0.5f * (temp[n - 1] + temp[n]); 1319 } 1320 1321 /// 1322 @safe pure 1323 unittest 1324 { 1325 import mir.ndslice; 1326 import mir.math; 1327 1328 auto a = mininitRcslice!double(3); 1329 a[] = [3, 1, 7]; 1330 auto med = median!double(a.flattened); 1331 assert(med == 3.0); 1332 assert(a == [3, 1, 7]);//add in stddev out param 1333 double aDev; 1334 auto aCenter = centerColumns(a, aDev, MeanEstimator.median); 1335 assert(aCenter == [0.0, -2.0, 4.0]); 1336 assert(aDev == 3.0); 1337 auto b = mininitRcslice!double(4); 1338 b[] = [4,2,5,1]; 1339 auto medB = median!double(b.flattened); 1340 assert(medB == 3.0); 1341 assert(b == [4,2,5,1]); 1342 double bDev; 1343 auto bCenter = centerColumns(b, bDev, MeanEstimator.median); 1344 assert(bCenter == [1.0, -1.0, 2.0, -2.0]); 1345 assert(bDev == 3.0); 1346 } 1347 1348 /++ 1349 Mean Centring of raw data. 1350 Params: 1351 matrix = input `M x N` matrix 1352 mean = column means 1353 est = mean estimation method 1354 Returns: 1355 `M x N` matrix with each column translated by the column mean 1356 +/ 1357 @safe pure nothrow @nogc 1358 Slice!(RCI!T,2) centerColumns(T, SliceKind kind) 1359 ( 1360 Slice!(const(T)*, 2, kind) matrix, 1361 out Slice!(RCI!T) mean, 1362 MeanEstimator est = MeanEstimator.average, 1363 ) 1364 { 1365 mean = rcslice!T([matrix.length!1], cast(T)0); 1366 if (est == MeanEstimator.none) 1367 { 1368 return matrix.as!T.rcslice; 1369 } 1370 auto at = matrix.transposed.as!T.rcslice; 1371 auto len = at.length!1; 1372 foreach (i; 0 .. at.length!0) 1373 { 1374 if (est == MeanEstimator.average) 1375 { 1376 foreach(j; 0 .. at.length!1) 1377 mean[i] += (at[i][j]/len); 1378 } 1379 else // (est == MeanEstimator.median) 1380 { 1381 mean[i] = median(at[i].flattened); 1382 } 1383 at[i][] -= mean[i]; 1384 } 1385 auto atSliced = at.transposed.as!T.rcslice; 1386 return atSliced; 1387 } 1388 1389 /// ditto 1390 @safe pure nothrow @nogc 1391 Slice!(RCI!T) centerColumns(T) 1392 ( 1393 Slice!(const(T)*) col, 1394 out T mean, 1395 MeanEstimator est = MeanEstimator.average, 1396 ) 1397 { 1398 mean = cast(T)0; 1399 if (est == MeanEstimator.none) 1400 { 1401 return col.as!T.rcslice; 1402 } 1403 auto len = col.length; 1404 if (est == MeanEstimator.average) 1405 { 1406 foreach(j; 0 .. len) 1407 mean += (col[j]/len); 1408 } 1409 else // (est == MeanEstimator.median) 1410 { 1411 mean = median(col); 1412 } 1413 auto result = mininitRcslice!T(len); 1414 foreach (j; 0 .. len) 1415 result[j] = col[j] - mean; 1416 return result; 1417 } 1418 1419 /// ditto 1420 @safe pure nothrow @nogc 1421 Slice!(RCI!T) centerColumns(T) 1422 ( 1423 auto ref const Slice!(RCI!T) col, 1424 out T mean, 1425 MeanEstimator est = MeanEstimator.average, 1426 ) 1427 { 1428 auto colScope = col.lightScope; 1429 return centerColumns!(T)(colScope, mean, est); 1430 } 1431 1432 /// ditto 1433 @safe pure nothrow @nogc 1434 Slice!(RCI!T,2) centerColumns(T, SliceKind kind) 1435 ( 1436 auto ref const Slice!(RCI!T,2,kind) matrix, 1437 out Slice!(RCI!T) mean, 1438 MeanEstimator est = MeanEstimator.average, 1439 ) 1440 { 1441 auto matrixScope = matrix.lightScope; 1442 return centerColumns(matrixScope, mean, est); 1443 } 1444 1445 /// 1446 @safe pure nothrow 1447 unittest 1448 { 1449 import mir.ndslice; 1450 import mir.math; 1451 1452 auto data = mininitRcslice!double(2,1); 1453 data[] = [[1], 1454 [3]]; 1455 auto mean = mininitRcslice!double(1); 1456 auto res = centerColumns(data, mean, MeanEstimator.average); 1457 assert(mean[0] == 2); 1458 assert(res == [[-1],[1]]); 1459 } 1460 1461 /// 1462 enum DeviationEstimator 1463 { 1464 /// 1465 none, 1466 /// 1467 sample, 1468 /// median absolute deviation 1469 mad 1470 } 1471 1472 /++ 1473 Normalization of raw data. 1474 Params: 1475 matrix = input `M x N` matrix, each row an observation and each column mean centred 1476 stdDev = column standard deviation 1477 devEst = estimation method 1478 Returns: 1479 `M x N` matrix with each column divided by it's standard deviation 1480 +/ 1481 @safe pure nothrow @nogc Slice!(RCI!T,2) normalizeColumns(T, SliceKind kind)( 1482 auto ref const Slice!(const(T)*,2,kind) matrix, 1483 out Slice!(RCI!T) stdDev, 1484 DeviationEstimator devEst = DeviationEstimator.sample 1485 ) 1486 { 1487 stdDev = rcslice!T([matrix.length!1], cast(T)0); 1488 if (devEst == DeviationEstimator.none) 1489 { 1490 auto matrixSliced = matrix.as!T.rcslice; 1491 return matrixSliced; 1492 } 1493 else 1494 { 1495 import mir.math.sum: sum; 1496 auto mTSliced = matrix.transposed.as!T.rcslice; 1497 auto mT = mTSliced.lightScope.canonical; 1498 foreach (i; 0 .. mT.length!0) 1499 { 1500 auto processedRow = mininitRcslice!T(mT.length!1); 1501 if (devEst == DeviationEstimator.sample) 1502 { 1503 foreach (j; 0 .. mT.length!1) 1504 processedRow[j] = mT[i, j] * mT[i, j]; 1505 stdDev[i] = sqrt(processedRow.sum!"kb2" / (mT[i].length - 1)); 1506 } 1507 else if (devEst == DeviationEstimator.mad) 1508 { 1509 foreach (j; 0 .. mT.length!1) 1510 processedRow[j] = fabs(mT[i,j]); 1511 stdDev[i] = median!T(processedRow); 1512 } 1513 mT[i][] /= stdDev[i]; 1514 } 1515 auto mSliced = mT.transposed.rcslice; 1516 return mSliced; 1517 } 1518 } 1519 1520 /// ditto 1521 @safe pure nothrow @nogc Slice!(RCI!T,2) normalizeColumns(T, SliceKind kind)( 1522 auto ref const Slice!(RCI!T,2,kind) matrix, 1523 out Slice!(RCI!T) stdDev, 1524 DeviationEstimator devEst = DeviationEstimator.sample, 1525 ) 1526 { 1527 auto matrixScope = matrix.lightScope.lightConst; 1528 return normalizeColumns!(T, kind)(matrixScope, stdDev, devEst); 1529 } 1530 1531 /// 1532 @safe pure nothrow unittest 1533 { 1534 import mir.ndslice; 1535 import mir.math; 1536 1537 auto data = mininitRcslice!double(2,2); 1538 data[] = [[ 2, -1], 1539 [-2, 1]]; 1540 //sd1 = 2 * sqrt(2.0); 1541 //sd2 = sqrt(2.0); 1542 auto x = 1.0 / sqrt(2.0); 1543 auto scaled = mininitRcslice!double(2,2); 1544 scaled[] = [[ x, -x], 1545 [-x, x]]; 1546 auto stdDev = mininitRcslice!double(2); 1547 assert(normalizeColumns(data, stdDev) == scaled); 1548 assert(stdDev == [2*sqrt(2.0), sqrt(2.0)]); 1549 }