github.com/gonum/matrix@v0.0.0-20181209220409-c518dec07be9/mat64/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 mat64 6 7 import ( 8 "math" 9 "math/rand" 10 "testing" 11 12 "github.com/gonum/blas/testblas" 13 "github.com/gonum/matrix" 14 ) 15 16 func TestCholesky(t *testing.T) { 17 for _, test := range []struct { 18 a *SymDense 19 20 cond float64 21 want *TriDense 22 posdef bool 23 }{ 24 { 25 a: NewSymDense(3, []float64{ 26 4, 1, 1, 27 0, 2, 3, 28 0, 0, 6, 29 }), 30 cond: 37, 31 want: NewTriDense(3, true, []float64{ 32 2, 0.5, 0.5, 33 0, 1.3228756555322954, 2.0788046015507495, 34 0, 0, 1.195228609334394, 35 }), 36 posdef: true, 37 }, 38 } { 39 _, n := test.a.Dims() 40 for _, chol := range []*Cholesky{ 41 {}, 42 {chol: NewTriDense(n-1, true, nil)}, 43 {chol: NewTriDense(n, true, nil)}, 44 {chol: NewTriDense(n+1, true, nil)}, 45 } { 46 ok := chol.Factorize(test.a) 47 if ok != test.posdef { 48 t.Errorf("unexpected return from Cholesky factorization: got: ok=%t want: ok=%t", ok, test.posdef) 49 } 50 fc := DenseCopyOf(chol.chol) 51 if !Equal(fc, test.want) { 52 t.Error("incorrect Cholesky factorization") 53 } 54 if math.Abs(test.cond-chol.cond) > 1e-13 { 55 t.Errorf("Condition number mismatch: Want %v, got %v", test.cond, chol.cond) 56 } 57 var U TriDense 58 U.UFromCholesky(chol) 59 aCopy := DenseCopyOf(test.a) 60 var a Dense 61 a.Mul(U.TTri(), &U) 62 if !EqualApprox(&a, aCopy, 1e-14) { 63 t.Error("unexpected Cholesky factor product") 64 } 65 66 var L TriDense 67 L.LFromCholesky(chol) 68 a.Mul(&L, L.TTri()) 69 if !EqualApprox(&a, aCopy, 1e-14) { 70 t.Error("unexpected Cholesky factor product") 71 } 72 } 73 } 74 } 75 76 func TestCholeskySolve(t *testing.T) { 77 for _, test := range []struct { 78 a *SymDense 79 b *Dense 80 ans *Dense 81 }{ 82 { 83 a: NewSymDense(2, []float64{ 84 1, 0, 85 0, 1, 86 }), 87 b: NewDense(2, 1, []float64{5, 6}), 88 ans: NewDense(2, 1, []float64{5, 6}), 89 }, 90 { 91 a: NewSymDense(3, []float64{ 92 53, 59, 37, 93 0, 83, 71, 94 37, 71, 101, 95 }), 96 b: NewDense(3, 1, []float64{5, 6, 7}), 97 ans: NewDense(3, 1, []float64{0.20745069393718094, -0.17421475529583694, 0.11577794010226464}), 98 }, 99 } { 100 var chol Cholesky 101 ok := chol.Factorize(test.a) 102 if !ok { 103 t.Fatal("unexpected Cholesky factorization failure: not positive definite") 104 } 105 106 var x Dense 107 x.SolveCholesky(&chol, test.b) 108 if !EqualApprox(&x, test.ans, 1e-12) { 109 t.Error("incorrect Cholesky solve solution") 110 } 111 112 var ans Dense 113 ans.Mul(test.a, &x) 114 if !EqualApprox(&ans, test.b, 1e-12) { 115 t.Error("incorrect Cholesky solve solution product") 116 } 117 } 118 } 119 120 func TestSolveTwoChol(t *testing.T) { 121 for _, test := range []struct { 122 a, b *SymDense 123 }{ 124 { 125 a: NewSymDense(2, []float64{ 126 1, 0, 127 0, 1, 128 }), 129 b: NewSymDense(2, []float64{ 130 1, 0, 131 0, 1, 132 }), 133 }, 134 { 135 a: NewSymDense(2, []float64{ 136 1, 0, 137 0, 1, 138 }), 139 b: NewSymDense(2, []float64{ 140 2, 0, 141 0, 2, 142 }), 143 }, 144 { 145 a: NewSymDense(3, []float64{ 146 53, 59, 37, 147 59, 83, 71, 148 37, 71, 101, 149 }), 150 b: NewSymDense(3, []float64{ 151 2, -1, 0, 152 -1, 2, -1, 153 0, -1, 2, 154 }), 155 }, 156 } { 157 var chola, cholb Cholesky 158 ok := chola.Factorize(test.a) 159 if !ok { 160 t.Fatal("unexpected Cholesky factorization failure for a: not positive definite") 161 } 162 ok = cholb.Factorize(test.b) 163 if !ok { 164 t.Fatal("unexpected Cholesky factorization failure for b: not positive definite") 165 } 166 167 var x Dense 168 x.solveTwoChol(&chola, &cholb) 169 170 var ans Dense 171 ans.Mul(test.a, &x) 172 if !EqualApprox(&ans, test.b, 1e-12) { 173 var y Dense 174 y.Solve(test.a, test.b) 175 t.Errorf("incorrect Cholesky solve solution product\ngot solution:\n%.4v\nwant solution\n%.4v", 176 Formatted(&x), Formatted(&y)) 177 } 178 } 179 } 180 181 func TestCholeskySolveVec(t *testing.T) { 182 for _, test := range []struct { 183 a *SymDense 184 b *Vector 185 ans *Vector 186 }{ 187 { 188 a: NewSymDense(2, []float64{ 189 1, 0, 190 0, 1, 191 }), 192 b: NewVector(2, []float64{5, 6}), 193 ans: NewVector(2, []float64{5, 6}), 194 }, 195 { 196 a: NewSymDense(3, []float64{ 197 53, 59, 37, 198 0, 83, 71, 199 0, 0, 101, 200 }), 201 b: NewVector(3, []float64{5, 6, 7}), 202 ans: NewVector(3, []float64{0.20745069393718094, -0.17421475529583694, 0.11577794010226464}), 203 }, 204 } { 205 var chol Cholesky 206 ok := chol.Factorize(test.a) 207 if !ok { 208 t.Fatal("unexpected Cholesky factorization failure: not positive definite") 209 } 210 211 var x Vector 212 x.SolveCholeskyVec(&chol, test.b) 213 if !EqualApprox(&x, test.ans, 1e-12) { 214 t.Error("incorrect Cholesky solve solution") 215 } 216 217 var ans Vector 218 ans.MulVec(test.a, &x) 219 if !EqualApprox(&ans, test.b, 1e-12) { 220 t.Error("incorrect Cholesky solve solution product") 221 } 222 } 223 } 224 225 func TestFromCholesky(t *testing.T) { 226 for _, test := range []*SymDense{ 227 NewSymDense(3, []float64{ 228 53, 59, 37, 229 0, 83, 71, 230 0, 0, 101, 231 }), 232 } { 233 var chol Cholesky 234 ok := chol.Factorize(test) 235 if !ok { 236 t.Fatal("unexpected Cholesky factorization failure: not positive definite") 237 } 238 var s SymDense 239 s.FromCholesky(&chol) 240 241 if !EqualApprox(&s, test, 1e-12) { 242 t.Errorf("Cholesky reconstruction not equal to original matrix.\nWant:\n% v\nGot:\n% v\n", Formatted(test), Formatted(&s)) 243 } 244 } 245 } 246 247 func TestCloneCholesky(t *testing.T) { 248 for _, test := range []*SymDense{ 249 NewSymDense(3, []float64{ 250 53, 59, 37, 251 0, 83, 71, 252 0, 0, 101, 253 }), 254 } { 255 var chol Cholesky 256 ok := chol.Factorize(test) 257 if !ok { 258 panic("bad test") 259 } 260 var chol2 Cholesky 261 chol2.Clone(&chol) 262 263 if chol.cond != chol2.cond { 264 t.Errorf("condition number mismatch from zero") 265 } 266 if !Equal(chol.chol, chol2.chol) { 267 t.Errorf("chol mismatch from zero") 268 } 269 270 // Corrupt chol2 and try again 271 chol2.cond = math.NaN() 272 chol2.chol = NewTriDense(2, matrix.Upper, nil) 273 chol2.Clone(&chol) 274 if chol.cond != chol2.cond { 275 t.Errorf("condition number mismatch from non-zero") 276 } 277 if !Equal(chol.chol, chol2.chol) { 278 t.Errorf("chol mismatch from non-zero") 279 } 280 } 281 } 282 283 func TestInverseCholesky(t *testing.T) { 284 for _, n := range []int{1, 3, 5, 9} { 285 data := make([]float64, n*n) 286 for i := range data { 287 data[i] = rand.NormFloat64() 288 } 289 var s SymDense 290 s.SymOuterK(1, NewDense(n, n, data)) 291 292 var chol Cholesky 293 ok := chol.Factorize(&s) 294 if !ok { 295 t.Errorf("Bad test, cholesky decomposition failed") 296 } 297 298 var sInv SymDense 299 sInv.InverseCholesky(&chol) 300 301 var ans Dense 302 ans.Mul(&sInv, &s) 303 if !equalApprox(eye(n), &ans, 1e-8, false) { 304 var diff Dense 305 diff.Sub(eye(n), &ans) 306 t.Errorf("SymDense times Cholesky inverse not identity. Norm diff = %v", Norm(&diff, 2)) 307 } 308 } 309 } 310 311 func TestCholeskySymRankOne(t *testing.T) { 312 rand.Seed(1) 313 for _, n := range []int{1, 2, 3, 4, 5, 7, 10, 20, 50, 100} { 314 for k := 0; k < 10; k++ { 315 data := make([]float64, n*n) 316 for i := range data { 317 data[i] = rand.NormFloat64() 318 } 319 320 var a SymDense 321 a.SymOuterK(1, NewDense(n, n, data)) 322 323 xdata := make([]float64, n) 324 for i := range xdata { 325 xdata[i] = rand.NormFloat64() 326 } 327 x := NewVector(n, xdata) 328 329 var chol Cholesky 330 ok := chol.Factorize(&a) 331 if !ok { 332 t.Errorf("Bad random test, Cholesky factorization failed") 333 continue 334 } 335 336 alpha := rand.Float64() 337 ok = chol.SymRankOne(&chol, alpha, x) 338 if !ok { 339 t.Errorf("n=%v, alpha=%v: unexpected failure", n, alpha) 340 continue 341 } 342 a.SymRankOne(&a, alpha, x) 343 344 var achol SymDense 345 achol.FromCholesky(&chol) 346 if !EqualApprox(&achol, &a, 1e-13) { 347 t.Errorf("n=%v, alpha=%v: mismatch between updated matrix and from Cholesky:\nupdated:\n%v\nfrom Cholesky:\n%v", 348 n, alpha, Formatted(&a), Formatted(&achol)) 349 } 350 } 351 } 352 353 for i, test := range []struct { 354 a *SymDense 355 alpha float64 356 x []float64 357 358 wantOk bool 359 }{ 360 { 361 // Update (to positive definite matrix). 362 a: NewSymDense(4, []float64{ 363 1, 1, 1, 1, 364 0, 2, 3, 4, 365 0, 0, 6, 10, 366 0, 0, 0, 20, 367 }), 368 alpha: 1, 369 x: []float64{0, 0, 0, 1}, 370 wantOk: true, 371 }, 372 { 373 // Downdate to singular matrix. 374 a: NewSymDense(4, []float64{ 375 1, 1, 1, 1, 376 0, 2, 3, 4, 377 0, 0, 6, 10, 378 0, 0, 0, 20, 379 }), 380 alpha: -1, 381 x: []float64{0, 0, 0, 1}, 382 wantOk: false, 383 }, 384 { 385 // Downdate to positive definite matrix. 386 a: NewSymDense(4, []float64{ 387 1, 1, 1, 1, 388 0, 2, 3, 4, 389 0, 0, 6, 10, 390 0, 0, 0, 20, 391 }), 392 alpha: -1 / 2, 393 x: []float64{0, 0, 0, 1}, 394 wantOk: true, 395 }, 396 } { 397 var chol Cholesky 398 ok := chol.Factorize(test.a) 399 if !ok { 400 t.Errorf("Case %v: bad test, Cholesky factorization failed", i) 401 continue 402 } 403 404 x := NewVector(len(test.x), test.x) 405 ok = chol.SymRankOne(&chol, test.alpha, x) 406 if !ok { 407 if test.wantOk { 408 t.Errorf("Case %v: unexpected failure from SymRankOne", i) 409 } 410 continue 411 } 412 if ok && !test.wantOk { 413 t.Errorf("Case %v: expected a failure from SymRankOne", i) 414 } 415 416 a := test.a 417 a.SymRankOne(a, test.alpha, x) 418 419 var achol SymDense 420 achol.FromCholesky(&chol) 421 if !EqualApprox(&achol, a, 1e-13) { 422 t.Errorf("Case %v: mismatch between updated matrix and from Cholesky:\nupdated:\n%v\nfrom Cholesky:\n%v", 423 i, Formatted(a), Formatted(&achol)) 424 } 425 } 426 } 427 428 func BenchmarkCholeskySmall(b *testing.B) { 429 benchmarkCholesky(b, 2) 430 } 431 432 func BenchmarkCholeskyMedium(b *testing.B) { 433 benchmarkCholesky(b, testblas.MediumMat) 434 } 435 436 func BenchmarkCholeskyLarge(b *testing.B) { 437 benchmarkCholesky(b, testblas.LargeMat) 438 } 439 440 func benchmarkCholesky(b *testing.B, n int) { 441 base := make([]float64, n*n) 442 for i := range base { 443 base[i] = rand.Float64() 444 } 445 bm := NewDense(n, n, base) 446 bm.Mul(bm.T(), bm) 447 am := NewSymDense(n, bm.mat.Data) 448 449 var chol Cholesky 450 b.ResetTimer() 451 for i := 0; i < b.N; i++ { 452 ok := chol.Factorize(am) 453 if !ok { 454 panic("not pos def") 455 } 456 } 457 }