gonum.org/v1/gonum@v0.14.0/mat/svd_example_test.go (about) 1 // Copyright ©2020 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_test 6 7 import ( 8 "fmt" 9 "log" 10 11 "gonum.org/v1/gonum/mat" 12 ) 13 14 func ExampleSVD_SolveTo() { 15 // The system described by A is rank deficient. 16 a := mat.NewDense(5, 3, []float64{ 17 -1.7854591879711257, -0.42687285925779594, -0.12730256811265162, 18 -0.5728984211439724, -0.10093393134001777, -0.1181901192353067, 19 1.2484316018707418, 0.5646683943038734, -0.48229492403243485, 20 0.10174927665169475, -0.5805410929482445, 1.3054473231942054, 21 -1.134174808195733, -0.4732430202414438, 0.3528489486370508, 22 }) 23 24 // Perform an SVD retaining all singular vectors. 25 var svd mat.SVD 26 ok := svd.Factorize(a, mat.SVDFull) 27 if !ok { 28 log.Fatal("failed to factorize A") 29 } 30 31 // Determine the rank of the A matrix with a near zero condition threshold. 32 const rcond = 1e-15 33 rank := svd.Rank(rcond) 34 if rank == 0 { 35 log.Fatal("zero rank system") 36 } 37 38 b := mat.NewDense(5, 2, []float64{ 39 -2.318, -4.35, 40 -0.715, 1.451, 41 1.836, -0.119, 42 -0.357, 3.094, 43 -1.636, 0.021, 44 }) 45 46 // Find a least-squares solution using the determined parts of the system. 47 var x mat.Dense 48 svd.SolveTo(&x, b, rank) 49 50 fmt.Printf("singular values = %v\nrank = %d\nx = %.15f", 51 format(svd.Values(nil), 4, rcond), rank, mat.Formatted(&x, mat.Prefix(" "))) 52 53 // Output: 54 // singular values = [2.685 1.526 <1e-15] 55 // rank = 2 56 // x = ⎡ 1.212064313552347 1.507467451093930⎤ 57 // ⎢ 0.415400738264774 -0.624498607705372⎥ 58 // ⎣-0.183184442255280 2.221334193689124⎦ 59 } 60 61 func ExampleSVD_SolveVecTo() { 62 // The system described by A is rank deficient. 63 a := mat.NewDense(5, 3, []float64{ 64 -1.7854591879711257, -0.42687285925779594, -0.12730256811265162, 65 -0.5728984211439724, -0.10093393134001777, -0.1181901192353067, 66 1.2484316018707418, 0.5646683943038734, -0.48229492403243485, 67 0.10174927665169475, -0.5805410929482445, 1.3054473231942054, 68 -1.134174808195733, -0.4732430202414438, 0.3528489486370508, 69 }) 70 71 // Perform an SVD retaining all singular vectors. 72 var svd mat.SVD 73 ok := svd.Factorize(a, mat.SVDFull) 74 if !ok { 75 log.Fatal("failed to factorize A") 76 } 77 78 // Determine the rank of the A matrix with a near zero condition threshold. 79 const rcond = 1e-15 80 rank := svd.Rank(rcond) 81 if rank == 0 { 82 log.Fatal("zero rank system") 83 } 84 85 b := mat.NewVecDense(5, []float64{-2.318, -0.715, 1.836, -0.357, -1.636}) 86 87 // Find a least-squares solution using the determined parts of the system. 88 var x mat.VecDense 89 svd.SolveVecTo(&x, b, rank) 90 91 fmt.Printf("singular values = %v\nrank = %d\nx = %.15f", 92 format(svd.Values(nil), 4, rcond), rank, mat.Formatted(&x, mat.Prefix(" "))) 93 94 // Output: 95 // singular values = [2.685 1.526 <1e-15] 96 // rank = 2 97 // x = ⎡ 1.212064313552347⎤ 98 // ⎢ 0.415400738264774⎥ 99 // ⎣-0.183184442255280⎦ 100 } 101 102 func format(vals []float64, prec int, eps float64) []string { 103 s := make([]string, len(vals)) 104 for i, v := range vals { 105 if v < eps { 106 s[i] = fmt.Sprintf("<%.*g", prec, eps) 107 continue 108 } 109 s[i] = fmt.Sprintf("%.*g", prec, v) 110 } 111 return s 112 }