github.com/gopherd/gonum@v0.0.4/mat/cholesky_test.go (about) 1 // Copyright ©2013 The Gonum Authors. All rights reserved. 2 // Use of this source code is governed by a BSD-style 3 // license that can be found in the LICENSE file. 4 5 package mat 6 7 import ( 8 "fmt" 9 "math" 10 "strconv" 11 "testing" 12 13 "math/rand" 14 15 "github.com/gopherd/gonum/floats/scalar" 16 ) 17 18 func TestCholesky(t *testing.T) { 19 t.Parallel() 20 for _, test := range []struct { 21 a *SymDense 22 23 cond float64 24 want *TriDense 25 posdef bool 26 }{ 27 { 28 a: NewSymDense(3, []float64{ 29 4, 1, 1, 30 0, 2, 3, 31 0, 0, 6, 32 }), 33 cond: 37, 34 want: NewTriDense(3, true, []float64{ 35 2, 0.5, 0.5, 36 0, 1.3228756555322954, 2.0788046015507495, 37 0, 0, 1.195228609334394, 38 }), 39 posdef: true, 40 }, 41 } { 42 _, n := test.a.Dims() 43 for _, chol := range []*Cholesky{ 44 {}, 45 {chol: NewTriDense(n-1, true, nil)}, 46 {chol: NewTriDense(n, true, nil)}, 47 {chol: NewTriDense(n+1, true, nil)}, 48 } { 49 ok := chol.Factorize(test.a) 50 if ok != test.posdef { 51 t.Errorf("unexpected return from Cholesky factorization: got: ok=%t want: ok=%t", ok, test.posdef) 52 } 53 fc := DenseCopyOf(chol.chol) 54 if !Equal(fc, test.want) { 55 t.Error("incorrect Cholesky factorization") 56 } 57 if math.Abs(test.cond-chol.cond) > 1e-13 { 58 t.Errorf("Condition number mismatch: Want %v, got %v", test.cond, chol.cond) 59 } 60 var U TriDense 61 chol.UTo(&U) 62 aCopy := DenseCopyOf(test.a) 63 var a Dense 64 a.Mul(U.TTri(), &U) 65 if !EqualApprox(&a, aCopy, 1e-14) { 66 t.Error("unexpected Cholesky factor product") 67 } 68 var L TriDense 69 chol.LTo(&L) 70 a.Mul(&L, L.TTri()) 71 if !EqualApprox(&a, aCopy, 1e-14) { 72 t.Error("unexpected Cholesky factor product") 73 } 74 } 75 } 76 } 77 78 func TestCholeskyAt(t *testing.T) { 79 t.Parallel() 80 for _, test := range []*SymDense{ 81 NewSymDense(3, []float64{ 82 53, 59, 37, 83 59, 83, 71, 84 37, 71, 101, 85 }), 86 } { 87 var chol Cholesky 88 ok := chol.Factorize(test) 89 if !ok { 90 t.Fatalf("Matrix not positive definite") 91 } 92 n := test.SymmetricDim() 93 cn := chol.SymmetricDim() 94 if cn != n { 95 t.Errorf("Cholesky size does not match. Got %d, want %d", cn, n) 96 } 97 for i := 0; i < n; i++ { 98 for j := 0; j < n; j++ { 99 got := chol.At(i, j) 100 want := test.At(i, j) 101 if math.Abs(got-want) > 1e-12 { 102 t.Errorf("Cholesky at does not match at %d, %d. Got %v, want %v", i, j, got, want) 103 } 104 } 105 } 106 } 107 } 108 109 func TestCholeskySolveTo(t *testing.T) { 110 t.Parallel() 111 for _, test := range []struct { 112 a *SymDense 113 b *Dense 114 ans *Dense 115 }{ 116 { 117 a: NewSymDense(2, []float64{ 118 1, 0, 119 0, 1, 120 }), 121 b: NewDense(2, 1, []float64{5, 6}), 122 ans: NewDense(2, 1, []float64{5, 6}), 123 }, 124 { 125 a: NewSymDense(3, []float64{ 126 53, 59, 37, 127 0, 83, 71, 128 37, 71, 101, 129 }), 130 b: NewDense(3, 1, []float64{5, 6, 7}), 131 ans: NewDense(3, 1, []float64{0.20745069393718094, -0.17421475529583694, 0.11577794010226464}), 132 }, 133 } { 134 var chol Cholesky 135 ok := chol.Factorize(test.a) 136 if !ok { 137 t.Fatal("unexpected Cholesky factorization failure: not positive definite") 138 } 139 140 var x Dense 141 err := chol.SolveTo(&x, test.b) 142 if err != nil { 143 t.Errorf("unexpected error from Cholesky solve: %v", err) 144 } 145 if !EqualApprox(&x, test.ans, 1e-12) { 146 t.Error("incorrect Cholesky solve solution") 147 } 148 149 var ans Dense 150 ans.Mul(test.a, &x) 151 if !EqualApprox(&ans, test.b, 1e-12) { 152 t.Error("incorrect Cholesky solve solution product") 153 } 154 } 155 } 156 157 func TestCholeskySolveCholTo(t *testing.T) { 158 t.Parallel() 159 for _, test := range []struct { 160 a, b *SymDense 161 }{ 162 { 163 a: NewSymDense(2, []float64{ 164 1, 0, 165 0, 1, 166 }), 167 b: NewSymDense(2, []float64{ 168 1, 0, 169 0, 1, 170 }), 171 }, 172 { 173 a: NewSymDense(2, []float64{ 174 1, 0, 175 0, 1, 176 }), 177 b: NewSymDense(2, []float64{ 178 2, 0, 179 0, 2, 180 }), 181 }, 182 { 183 a: NewSymDense(3, []float64{ 184 53, 59, 37, 185 59, 83, 71, 186 37, 71, 101, 187 }), 188 b: NewSymDense(3, []float64{ 189 2, -1, 0, 190 -1, 2, -1, 191 0, -1, 2, 192 }), 193 }, 194 } { 195 var chola, cholb Cholesky 196 ok := chola.Factorize(test.a) 197 if !ok { 198 t.Fatal("unexpected Cholesky factorization failure for a: not positive definite") 199 } 200 ok = cholb.Factorize(test.b) 201 if !ok { 202 t.Fatal("unexpected Cholesky factorization failure for b: not positive definite") 203 } 204 205 var x Dense 206 err := chola.SolveCholTo(&x, &cholb) 207 if err != nil { 208 t.Errorf("unexpected error from Cholesky solve: %v", err) 209 } 210 211 var ans Dense 212 ans.Mul(test.a, &x) 213 if !EqualApprox(&ans, test.b, 1e-12) { 214 var y Dense 215 err := y.Solve(test.a, test.b) 216 if err != nil { 217 t.Errorf("unexpected error from dense solve: %v", err) 218 } 219 t.Errorf("incorrect Cholesky solve solution product\ngot solution:\n%.4v\nwant solution\n%.4v", 220 Formatted(&x), Formatted(&y)) 221 } 222 } 223 } 224 225 func TestCholeskySolveVecTo(t *testing.T) { 226 t.Parallel() 227 for _, test := range []struct { 228 a *SymDense 229 b *VecDense 230 ans *VecDense 231 }{ 232 { 233 a: NewSymDense(2, []float64{ 234 1, 0, 235 0, 1, 236 }), 237 b: NewVecDense(2, []float64{5, 6}), 238 ans: NewVecDense(2, []float64{5, 6}), 239 }, 240 { 241 a: NewSymDense(3, []float64{ 242 53, 59, 37, 243 0, 83, 71, 244 0, 0, 101, 245 }), 246 b: NewVecDense(3, []float64{5, 6, 7}), 247 ans: NewVecDense(3, []float64{0.20745069393718094, -0.17421475529583694, 0.11577794010226464}), 248 }, 249 } { 250 var chol Cholesky 251 ok := chol.Factorize(test.a) 252 if !ok { 253 t.Fatal("unexpected Cholesky factorization failure: not positive definite") 254 } 255 256 var x VecDense 257 err := chol.SolveVecTo(&x, test.b) 258 if err != nil { 259 t.Errorf("unexpected error from Cholesky solve: %v", err) 260 } 261 if !EqualApprox(&x, test.ans, 1e-12) { 262 t.Error("incorrect Cholesky solve solution") 263 } 264 265 var ans VecDense 266 ans.MulVec(test.a, &x) 267 if !EqualApprox(&ans, test.b, 1e-12) { 268 t.Error("incorrect Cholesky solve solution product") 269 } 270 } 271 } 272 273 func TestCholeskyToSym(t *testing.T) { 274 t.Parallel() 275 for _, test := range []*SymDense{ 276 NewSymDense(3, []float64{ 277 53, 59, 37, 278 0, 83, 71, 279 0, 0, 101, 280 }), 281 } { 282 var chol Cholesky 283 ok := chol.Factorize(test) 284 if !ok { 285 t.Fatal("unexpected Cholesky factorization failure: not positive definite") 286 } 287 var s SymDense 288 chol.ToSym(&s) 289 290 if !EqualApprox(&s, test, 1e-12) { 291 t.Errorf("Cholesky reconstruction not equal to original matrix.\nWant:\n% v\nGot:\n% v\n", Formatted(test), Formatted(&s)) 292 } 293 } 294 } 295 296 func TestCloneCholesky(t *testing.T) { 297 t.Parallel() 298 for _, test := range []*SymDense{ 299 NewSymDense(3, []float64{ 300 53, 59, 37, 301 0, 83, 71, 302 0, 0, 101, 303 }), 304 } { 305 var chol Cholesky 306 ok := chol.Factorize(test) 307 if !ok { 308 panic("bad test") 309 } 310 var chol2 Cholesky 311 chol2.Clone(&chol) 312 313 if chol.cond != chol2.cond { 314 t.Errorf("condition number mismatch from empty") 315 } 316 if !Equal(chol.chol, chol2.chol) { 317 t.Errorf("chol mismatch from empty") 318 } 319 320 // Corrupt chol2 and try again 321 chol2.cond = math.NaN() 322 chol2.chol = NewTriDense(2, Upper, nil) 323 chol2.Clone(&chol) 324 if chol.cond != chol2.cond { 325 t.Errorf("condition number mismatch from non-empty") 326 } 327 if !Equal(chol.chol, chol2.chol) { 328 t.Errorf("chol mismatch from non-empty") 329 } 330 } 331 } 332 333 func TestCholeskyInverseTo(t *testing.T) { 334 t.Parallel() 335 rnd := rand.New(rand.NewSource(1)) 336 for _, n := range []int{1, 3, 5, 9} { 337 data := make([]float64, n*n) 338 for i := range data { 339 data[i] = rnd.NormFloat64() 340 } 341 var s SymDense 342 s.SymOuterK(1, NewDense(n, n, data)) 343 344 var chol Cholesky 345 ok := chol.Factorize(&s) 346 if !ok { 347 t.Errorf("Bad test, cholesky decomposition failed") 348 } 349 350 var sInv SymDense 351 err := chol.InverseTo(&sInv) 352 if err != nil { 353 t.Errorf("unexpected error from Cholesky inverse: %v", err) 354 } 355 356 var ans Dense 357 ans.Mul(&sInv, &s) 358 if !equalApprox(eye(n), &ans, 1e-8, false) { 359 var diff Dense 360 diff.Sub(eye(n), &ans) 361 t.Errorf("SymDense times Cholesky inverse not identity. Norm diff = %v", Norm(&diff, 2)) 362 } 363 } 364 } 365 366 func TestCholeskySymRankOne(t *testing.T) { 367 t.Parallel() 368 rnd := rand.New(rand.NewSource(1)) 369 for _, n := range []int{1, 2, 3, 4, 5, 7, 10, 20, 50, 100} { 370 for k := 0; k < 50; k++ { 371 // Construct a random positive definite matrix. 372 data := make([]float64, n*n) 373 for i := range data { 374 data[i] = rnd.NormFloat64() 375 } 376 var a SymDense 377 a.SymOuterK(1, NewDense(n, n, data)) 378 379 // Construct random data for updating. 380 xdata := make([]float64, n) 381 for i := range xdata { 382 xdata[i] = rnd.NormFloat64() 383 } 384 x := NewVecDense(n, xdata) 385 alpha := rnd.NormFloat64() 386 387 // Compute the updated matrix directly. If alpha > 0, there are no 388 // issues. If alpha < 0, it could be that the final matrix is not 389 // positive definite, so instead switch the two matrices. 390 aUpdate := NewSymDense(n, nil) 391 if alpha > 0 { 392 aUpdate.SymRankOne(&a, alpha, x) 393 } else { 394 aUpdate.CopySym(&a) 395 a.Reset() 396 a.SymRankOne(aUpdate, -alpha, x) 397 } 398 399 // Compare the Cholesky decomposition computed with Cholesky.SymRankOne 400 // with that computed from updating A directly. 401 var chol Cholesky 402 ok := chol.Factorize(&a) 403 if !ok { 404 t.Errorf("Bad random test, Cholesky factorization failed") 405 continue 406 } 407 408 var cholUpdate Cholesky 409 ok = cholUpdate.SymRankOne(&chol, alpha, x) 410 if !ok { 411 t.Errorf("n=%v, alpha=%v: unexpected failure", n, alpha) 412 continue 413 } 414 415 var aCompare SymDense 416 cholUpdate.ToSym(&aCompare) 417 if !EqualApprox(&aCompare, aUpdate, 1e-13) { 418 t.Errorf("n=%v, alpha=%v: mismatch between updated matrix and from Cholesky:\nupdated:\n%v\nfrom Cholesky:\n%v", 419 n, alpha, Formatted(aUpdate), Formatted(&aCompare)) 420 } 421 } 422 } 423 424 for i, test := range []struct { 425 a *SymDense 426 alpha float64 427 x []float64 428 429 wantOk bool 430 }{ 431 { 432 // Update (to positive definite matrix). 433 a: NewSymDense(4, []float64{ 434 1, 1, 1, 1, 435 0, 2, 3, 4, 436 0, 0, 6, 10, 437 0, 0, 0, 20, 438 }), 439 alpha: 1, 440 x: []float64{0, 0, 0, 1}, 441 wantOk: true, 442 }, 443 { 444 // Downdate to singular matrix. 445 a: NewSymDense(4, []float64{ 446 1, 1, 1, 1, 447 0, 2, 3, 4, 448 0, 0, 6, 10, 449 0, 0, 0, 20, 450 }), 451 alpha: -1, 452 x: []float64{0, 0, 0, 1}, 453 wantOk: false, 454 }, 455 { 456 // Downdate to positive definite matrix. 457 a: NewSymDense(4, []float64{ 458 1, 1, 1, 1, 459 0, 2, 3, 4, 460 0, 0, 6, 10, 461 0, 0, 0, 20, 462 }), 463 alpha: -1 / 2, 464 x: []float64{0, 0, 0, 1}, 465 wantOk: true, 466 }, 467 { 468 // Issue #453. 469 a: NewSymDense(1, []float64{1}), 470 alpha: -1, 471 x: []float64{0.25}, 472 wantOk: true, 473 }, 474 } { 475 var chol Cholesky 476 ok := chol.Factorize(test.a) 477 if !ok { 478 t.Errorf("Case %v: bad test, Cholesky factorization failed", i) 479 continue 480 } 481 482 x := NewVecDense(len(test.x), test.x) 483 ok = chol.SymRankOne(&chol, test.alpha, x) 484 if !ok { 485 if test.wantOk { 486 t.Errorf("Case %v: unexpected failure from SymRankOne", i) 487 } 488 continue 489 } 490 if ok && !test.wantOk { 491 t.Errorf("Case %v: expected a failure from SymRankOne", i) 492 } 493 494 a := test.a 495 a.SymRankOne(a, test.alpha, x) 496 497 var achol SymDense 498 chol.ToSym(&achol) 499 if !EqualApprox(&achol, a, 1e-13) { 500 t.Errorf("Case %v: mismatch between updated matrix and from Cholesky:\nupdated:\n%v\nfrom Cholesky:\n%v", 501 i, Formatted(a), Formatted(&achol)) 502 } 503 } 504 } 505 506 func TestCholeskyExtendVecSym(t *testing.T) { 507 t.Parallel() 508 for cas, test := range []struct { 509 a *SymDense 510 }{ 511 { 512 a: NewSymDense(3, []float64{ 513 4, 1, 1, 514 0, 2, 3, 515 0, 0, 6, 516 }), 517 }, 518 } { 519 n := test.a.SymmetricDim() 520 as := test.a.sliceSym(0, n-1) 521 522 // Compute the full factorization to use later (do the full factorization 523 // first to ensure the matrix is positive definite). 524 var cholFull Cholesky 525 ok := cholFull.Factorize(test.a) 526 if !ok { 527 panic("mat: bad test, matrix not positive definite") 528 } 529 530 var chol Cholesky 531 ok = chol.Factorize(as) 532 if !ok { 533 panic("mat: bad test, subset is not positive definite") 534 } 535 row := NewVecDense(n, nil) 536 for i := 0; i < n; i++ { 537 row.SetVec(i, test.a.At(n-1, i)) 538 } 539 540 var cholNew Cholesky 541 ok = cholNew.ExtendVecSym(&chol, row) 542 if !ok { 543 t.Errorf("cas %v: update not positive definite", cas) 544 } 545 var a SymDense 546 cholNew.ToSym(&a) 547 if !EqualApprox(&a, test.a, 1e-12) { 548 t.Errorf("cas %v: mismatch", cas) 549 } 550 551 // test in-place 552 ok = chol.ExtendVecSym(&chol, row) 553 if !ok { 554 t.Errorf("cas %v: in-place update not positive definite", cas) 555 } 556 if !equalChol(&chol, &cholNew) { 557 t.Errorf("cas %v: Cholesky different in-place vs. new", cas) 558 } 559 560 // Test that the factorization is about right compared with the direct 561 // full factorization. Use a high tolerance on the condition number 562 // since the condition number with the updated rule is approximate. 563 if !equalApproxChol(&chol, &cholFull, 1e-12, 0.3) { 564 t.Errorf("cas %v: updated Cholesky does not match full", cas) 565 } 566 } 567 } 568 569 func TestCholeskyScale(t *testing.T) { 570 t.Parallel() 571 for cas, test := range []struct { 572 a *SymDense 573 f float64 574 }{ 575 { 576 a: NewSymDense(3, []float64{ 577 4, 1, 1, 578 0, 2, 3, 579 0, 0, 6, 580 }), 581 f: 0.5, 582 }, 583 } { 584 var chol Cholesky 585 ok := chol.Factorize(test.a) 586 if !ok { 587 t.Errorf("Case %v: bad test, Cholesky factorization failed", cas) 588 continue 589 } 590 591 // Compare the update to a new Cholesky to an update in-place. 592 var cholUpdate Cholesky 593 cholUpdate.Scale(test.f, &chol) 594 chol.Scale(test.f, &chol) 595 if !equalChol(&chol, &cholUpdate) { 596 t.Errorf("Case %d: cholesky mismatch new receiver", cas) 597 } 598 var sym SymDense 599 chol.ToSym(&sym) 600 var comp SymDense 601 comp.ScaleSym(test.f, test.a) 602 if !EqualApprox(&comp, &sym, 1e-14) { 603 t.Errorf("Case %d: cholesky reconstruction doesn't match scaled matrix", cas) 604 } 605 606 var cholTest Cholesky 607 cholTest.Factorize(&comp) 608 if !equalApproxChol(&cholTest, &chol, 1e-12, 1e-12) { 609 t.Errorf("Case %d: cholesky mismatch with scaled matrix. %v, %v", cas, cholTest.cond, chol.cond) 610 } 611 } 612 } 613 614 // equalApproxChol checks that the two Cholesky decompositions are equal. 615 func equalChol(a, b *Cholesky) bool { 616 return Equal(a.chol, b.chol) && a.cond == b.cond 617 } 618 619 // equalApproxChol checks that the two Cholesky decompositions are approximately 620 // the same with the given tolerance on equality for the Triangular component and 621 // condition. 622 func equalApproxChol(a, b *Cholesky, matTol, condTol float64) bool { 623 if !EqualApprox(a.chol, b.chol, matTol) { 624 return false 625 } 626 return scalar.EqualWithinAbsOrRel(a.cond, b.cond, condTol, condTol) 627 } 628 629 func BenchmarkCholeskyFactorize(b *testing.B) { 630 for _, n := range []int{10, 100, 1000} { 631 b.Run("n="+strconv.Itoa(n), func(b *testing.B) { 632 rnd := rand.New(rand.NewSource(1)) 633 634 data := make([]float64, n*n) 635 for i := range data { 636 data[i] = rnd.NormFloat64() 637 } 638 var a SymDense 639 a.SymOuterK(1, NewDense(n, n, data)) 640 641 var chol Cholesky 642 b.ResetTimer() 643 for i := 0; i < b.N; i++ { 644 ok := chol.Factorize(&a) 645 if !ok { 646 panic("not positive definite") 647 } 648 } 649 }) 650 } 651 } 652 653 func BenchmarkCholeskyToSym(b *testing.B) { 654 for _, n := range []int{10, 100, 1000} { 655 b.Run("n="+strconv.Itoa(n), func(b *testing.B) { 656 rnd := rand.New(rand.NewSource(1)) 657 658 data := make([]float64, n*n) 659 for i := range data { 660 data[i] = rnd.NormFloat64() 661 } 662 var a SymDense 663 a.SymOuterK(1, NewDense(n, n, data)) 664 665 var chol Cholesky 666 ok := chol.Factorize(&a) 667 if !ok { 668 panic("not positive definite") 669 } 670 671 dst := NewSymDense(n, nil) 672 673 b.ResetTimer() 674 for i := 0; i < b.N; i++ { 675 chol.ToSym(dst) 676 } 677 }) 678 } 679 } 680 681 func BenchmarkCholeskyInverseTo(b *testing.B) { 682 for _, n := range []int{10, 100, 1000} { 683 b.Run("n="+strconv.Itoa(n), func(b *testing.B) { 684 rnd := rand.New(rand.NewSource(1)) 685 686 data := make([]float64, n*n) 687 for i := range data { 688 data[i] = rnd.NormFloat64() 689 } 690 var a SymDense 691 a.SymOuterK(1, NewDense(n, n, data)) 692 693 var chol Cholesky 694 ok := chol.Factorize(&a) 695 if !ok { 696 panic("not positive definite") 697 } 698 699 dst := NewSymDense(n, nil) 700 701 b.ResetTimer() 702 for i := 0; i < b.N; i++ { 703 err := chol.InverseTo(dst) 704 if err != nil { 705 b.Fatalf("unexpected error from Cholesky inverse: %v", err) 706 } 707 } 708 }) 709 } 710 } 711 712 func TestBandCholeskySolveTo(t *testing.T) { 713 t.Parallel() 714 715 const ( 716 nrhs = 4 717 tol = 1e-14 718 ) 719 rnd := rand.New(rand.NewSource(1)) 720 for _, n := range []int{1, 2, 3, 5, 10} { 721 for _, k := range []int{0, 1, n / 2, n - 1} { 722 k := min(k, n-1) 723 724 a := NewSymBandDense(n, k, nil) 725 for i := 0; i < n; i++ { 726 a.SetSymBand(i, i, rnd.Float64()+float64(n)) 727 for j := i + 1; j < min(i+k+1, n); j++ { 728 a.SetSymBand(i, j, rnd.Float64()) 729 } 730 } 731 732 want := NewDense(n, nrhs, nil) 733 for i := 0; i < n; i++ { 734 for j := 0; j < nrhs; j++ { 735 want.Set(i, j, rnd.NormFloat64()) 736 } 737 } 738 var b Dense 739 b.Mul(a, want) 740 741 for _, typ := range []SymBanded{a, (*basicSymBanded)(a)} { 742 name := fmt.Sprintf("Case n=%d,k=%d,type=%T,nrhs=%d", n, k, typ, nrhs) 743 744 var chol BandCholesky 745 ok := chol.Factorize(typ) 746 if !ok { 747 t.Fatalf("%v: Factorize failed", name) 748 } 749 750 var got Dense 751 err := chol.SolveTo(&got, &b) 752 if err != nil { 753 t.Errorf("%v: unexpected error from SolveTo: %v", name, err) 754 continue 755 } 756 757 var resid Dense 758 resid.Sub(want, &got) 759 diff := Norm(&resid, math.Inf(1)) 760 if diff > tol { 761 t.Errorf("%v: unexpected solution; diff=%v", name, diff) 762 } 763 764 got.Copy(&b) 765 err = chol.SolveTo(&got, &got) 766 if err != nil { 767 t.Errorf("%v: unexpected error from SolveTo when dst==b: %v", name, err) 768 continue 769 } 770 771 resid.Sub(want, &got) 772 diff = Norm(&resid, math.Inf(1)) 773 if diff > tol { 774 t.Errorf("%v: unexpected solution when dst==b; diff=%v", name, diff) 775 } 776 } 777 } 778 } 779 } 780 781 func TestBandCholeskySolveVecTo(t *testing.T) { 782 t.Parallel() 783 784 const tol = 1e-14 785 rnd := rand.New(rand.NewSource(1)) 786 for _, n := range []int{1, 2, 3, 5, 10} { 787 for _, k := range []int{0, 1, n / 2, n - 1} { 788 k := min(k, n-1) 789 790 a := NewSymBandDense(n, k, nil) 791 for i := 0; i < n; i++ { 792 a.SetSymBand(i, i, rnd.Float64()+float64(n)) 793 for j := i + 1; j < min(i+k+1, n); j++ { 794 a.SetSymBand(i, j, rnd.Float64()) 795 } 796 } 797 798 want := NewVecDense(n, nil) 799 for i := 0; i < n; i++ { 800 want.SetVec(i, rnd.NormFloat64()) 801 } 802 var b VecDense 803 b.MulVec(a, want) 804 805 for _, typ := range []SymBanded{a, (*basicSymBanded)(a)} { 806 name := fmt.Sprintf("Case n=%d,k=%d,type=%T", n, k, typ) 807 808 var chol BandCholesky 809 ok := chol.Factorize(typ) 810 if !ok { 811 t.Fatalf("%v: Factorize failed", name) 812 } 813 814 var got VecDense 815 err := chol.SolveVecTo(&got, &b) 816 if err != nil { 817 t.Errorf("%v: unexpected error from SolveVecTo: %v", name, err) 818 continue 819 } 820 821 var resid VecDense 822 resid.SubVec(want, &got) 823 diff := Norm(&resid, math.Inf(1)) 824 if diff > tol { 825 t.Errorf("%v: unexpected solution; diff=%v", name, diff) 826 } 827 828 got.CopyVec(&b) 829 err = chol.SolveVecTo(&got, &got) 830 if err != nil { 831 t.Errorf("%v: unexpected error from SolveVecTo when dst==b: %v", name, err) 832 continue 833 } 834 835 resid.SubVec(want, &got) 836 diff = Norm(&resid, math.Inf(1)) 837 if diff > tol { 838 t.Errorf("%v: unexpected solution when dst==b; diff=%v", name, diff) 839 } 840 } 841 } 842 } 843 } 844 845 func TestBandCholeskyAt(t *testing.T) { 846 t.Parallel() 847 848 const tol = 1e-14 849 rnd := rand.New(rand.NewSource(1)) 850 for _, n := range []int{1, 2, 3, 5, 10} { 851 for _, k := range []int{0, 1, n / 2, n - 1} { 852 k := min(k, n-1) 853 name := fmt.Sprintf("Case n=%d,k=%d", n, k) 854 855 a := NewSymBandDense(n, k, nil) 856 for i := 0; i < n; i++ { 857 a.SetSymBand(i, i, rnd.Float64()+float64(n)) 858 for j := i + 1; j < min(i+k+1, n); j++ { 859 a.SetSymBand(i, j, rnd.Float64()) 860 } 861 } 862 863 var chol BandCholesky 864 ok := chol.Factorize(a) 865 if !ok { 866 t.Fatalf("%v: Factorize failed", name) 867 } 868 869 resid := NewDense(n, n, nil) 870 for i := 0; i < n; i++ { 871 for j := 0; j < n; j++ { 872 resid.Set(i, j, math.Abs(a.At(i, j)-chol.At(i, j))) 873 } 874 } 875 diff := Norm(resid, math.Inf(1)) 876 if diff > tol { 877 t.Errorf("%v: unexpected result; diff=%v, want<=%v", name, diff, tol) 878 } 879 } 880 } 881 } 882 883 func TestBandCholeskyDet(t *testing.T) { 884 t.Parallel() 885 886 const tol = 1e-14 887 rnd := rand.New(rand.NewSource(1)) 888 for _, n := range []int{1, 2, 3, 5, 10} { 889 for _, k := range []int{0, 1, n / 2, n - 1} { 890 k := min(k, n-1) 891 name := fmt.Sprintf("Case n=%d,k=%d", n, k) 892 893 a := NewSymBandDense(n, k, nil) 894 aSym := NewSymDense(n, nil) 895 for i := 0; i < n; i++ { 896 aii := rnd.Float64() + float64(n) 897 a.SetSymBand(i, i, aii) 898 aSym.SetSym(i, i, aii) 899 for j := i + 1; j < min(i+k+1, n); j++ { 900 aij := rnd.Float64() 901 a.SetSymBand(i, j, aij) 902 aSym.SetSym(i, j, aij) 903 } 904 } 905 906 var chol BandCholesky 907 ok := chol.Factorize(a) 908 if !ok { 909 t.Fatalf("%v: Factorize failed", name) 910 } 911 912 var cholDense Cholesky 913 ok = cholDense.Factorize(aSym) 914 if !ok { 915 t.Fatalf("%v: dense Factorize failed", name) 916 } 917 918 want := cholDense.Det() 919 got := chol.Det() 920 diff := math.Abs(got - want) 921 if diff > tol { 922 t.Errorf("%v: unexpected result; got=%v, want=%v (diff=%v)", name, got, want, diff) 923 } 924 } 925 } 926 }