github.com/gopherd/gonum@v0.0.4/stat/distmv/normal.go (about) 1 // Copyright ©2015 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 distmv 6 7 import ( 8 "math" 9 10 "math/rand" 11 12 "github.com/gopherd/gonum/floats" 13 "github.com/gopherd/gonum/mat" 14 "github.com/gopherd/gonum/stat" 15 "github.com/gopherd/gonum/stat/distuv" 16 ) 17 18 const badInputLength = "distmv: input slice length mismatch" 19 20 // Normal is a multivariate normal distribution (also known as the multivariate 21 // Gaussian distribution). Its pdf in k dimensions is given by 22 // (2 π)^(-k/2) |Σ|^(-1/2) exp(-1/2 (x-μ)'Σ^-1(x-μ)) 23 // where μ is the mean vector and Σ the covariance matrix. Σ must be symmetric 24 // and positive definite. Use NewNormal to construct. 25 type Normal struct { 26 mu []float64 27 28 sigma mat.SymDense 29 30 chol mat.Cholesky 31 logSqrtDet float64 32 dim int 33 34 // If src is altered, rnd must be updated. 35 src rand.Source 36 rnd *rand.Rand 37 } 38 39 // NewNormal creates a new Normal with the given mean and covariance matrix. 40 // NewNormal panics if len(mu) == 0, or if len(mu) != sigma.N. If the covariance 41 // matrix is not positive-definite, the returned boolean is false. 42 func NewNormal(mu []float64, sigma mat.Symmetric, src rand.Source) (*Normal, bool) { 43 if len(mu) == 0 { 44 panic(badZeroDimension) 45 } 46 dim := sigma.SymmetricDim() 47 if dim != len(mu) { 48 panic(badSizeMismatch) 49 } 50 n := &Normal{ 51 src: src, 52 rnd: rand.New(src), 53 dim: dim, 54 mu: make([]float64, dim), 55 } 56 copy(n.mu, mu) 57 ok := n.chol.Factorize(sigma) 58 if !ok { 59 return nil, false 60 } 61 n.sigma = *mat.NewSymDense(dim, nil) 62 n.sigma.CopySym(sigma) 63 n.logSqrtDet = 0.5 * n.chol.LogDet() 64 return n, true 65 } 66 67 // NewNormalChol creates a new Normal distribution with the given mean and 68 // covariance matrix represented by its Cholesky decomposition. NewNormalChol 69 // panics if len(mu) is not equal to chol.Size(). 70 func NewNormalChol(mu []float64, chol *mat.Cholesky, src rand.Source) *Normal { 71 dim := len(mu) 72 if dim != chol.SymmetricDim() { 73 panic(badSizeMismatch) 74 } 75 n := &Normal{ 76 src: src, 77 rnd: rand.New(src), 78 dim: dim, 79 mu: make([]float64, dim), 80 } 81 n.chol.Clone(chol) 82 copy(n.mu, mu) 83 n.logSqrtDet = 0.5 * n.chol.LogDet() 84 return n 85 } 86 87 // NewNormalPrecision creates a new Normal distribution with the given mean and 88 // precision matrix (inverse of the covariance matrix). NewNormalPrecision 89 // panics if len(mu) is not equal to prec.SymmetricDim(). If the precision matrix 90 // is not positive-definite, NewNormalPrecision returns nil for norm and false 91 // for ok. 92 func NewNormalPrecision(mu []float64, prec *mat.SymDense, src rand.Source) (norm *Normal, ok bool) { 93 if len(mu) == 0 { 94 panic(badZeroDimension) 95 } 96 dim := prec.SymmetricDim() 97 if dim != len(mu) { 98 panic(badSizeMismatch) 99 } 100 // TODO(btracey): Computing a matrix inverse is generally numerically instable. 101 // This only has to compute the inverse of a positive definite matrix, which 102 // is much better, but this still loses precision. It is worth considering if 103 // instead the precision matrix should be stored explicitly and used instead 104 // of the Cholesky decomposition of the covariance matrix where appropriate. 105 var chol mat.Cholesky 106 ok = chol.Factorize(prec) 107 if !ok { 108 return nil, false 109 } 110 var sigma mat.SymDense 111 err := chol.InverseTo(&sigma) 112 if err != nil { 113 return nil, false 114 } 115 return NewNormal(mu, &sigma, src) 116 } 117 118 // ConditionNormal returns the Normal distribution that is the receiver conditioned 119 // on the input evidence. The returned multivariate normal has dimension 120 // n - len(observed), where n is the dimension of the original receiver. The updated 121 // mean and covariance are 122 // mu = mu_un + sigma_{ob,un}ᵀ * sigma_{ob,ob}^-1 (v - mu_ob) 123 // sigma = sigma_{un,un} - sigma_{ob,un}ᵀ * sigma_{ob,ob}^-1 * sigma_{ob,un} 124 // where mu_un and mu_ob are the original means of the unobserved and observed 125 // variables respectively, sigma_{un,un} is the unobserved subset of the covariance 126 // matrix, sigma_{ob,ob} is the observed subset of the covariance matrix, and 127 // sigma_{un,ob} are the cross terms. The elements of x_2 have been observed with 128 // values v. The dimension order is preserved during conditioning, so if the value 129 // of dimension 1 is observed, the returned normal represents dimensions {0, 2, ...} 130 // of the original Normal distribution. 131 // 132 // ConditionNormal returns {nil, false} if there is a failure during the update. 133 // Mathematically this is impossible, but can occur with finite precision arithmetic. 134 func (n *Normal) ConditionNormal(observed []int, values []float64, src rand.Source) (*Normal, bool) { 135 if len(observed) == 0 { 136 panic("normal: no observed value") 137 } 138 if len(observed) != len(values) { 139 panic(badInputLength) 140 } 141 for _, v := range observed { 142 if v < 0 || v >= n.Dim() { 143 panic("normal: observed value out of bounds") 144 } 145 } 146 147 _, mu1, sigma11 := studentsTConditional(observed, values, math.Inf(1), n.mu, &n.sigma) 148 if mu1 == nil { 149 return nil, false 150 } 151 return NewNormal(mu1, sigma11, src) 152 } 153 154 // CovarianceMatrix stores the covariance matrix of the distribution in dst. 155 // Upon return, the value at element {i, j} of the covariance matrix is equal 156 // to the covariance of the i^th and j^th variables. 157 // covariance(i, j) = E[(x_i - E[x_i])(x_j - E[x_j])] 158 // If the dst matrix is empty it will be resized to the correct dimensions, 159 // otherwise dst must match the dimension of the receiver or CovarianceMatrix 160 // will panic. 161 func (n *Normal) CovarianceMatrix(dst *mat.SymDense) { 162 if dst.IsEmpty() { 163 *dst = *(dst.GrowSym(n.dim).(*mat.SymDense)) 164 } else if dst.SymmetricDim() != n.dim { 165 panic("normal: input matrix size mismatch") 166 } 167 dst.CopySym(&n.sigma) 168 } 169 170 // Dim returns the dimension of the distribution. 171 func (n *Normal) Dim() int { 172 return n.dim 173 } 174 175 // Entropy returns the differential entropy of the distribution. 176 func (n *Normal) Entropy() float64 { 177 return float64(n.dim)/2*(1+logTwoPi) + n.logSqrtDet 178 } 179 180 // LogProb computes the log of the pdf of the point x. 181 func (n *Normal) LogProb(x []float64) float64 { 182 dim := n.dim 183 if len(x) != dim { 184 panic(badSizeMismatch) 185 } 186 return normalLogProb(x, n.mu, &n.chol, n.logSqrtDet) 187 } 188 189 // NormalLogProb computes the log probability of the location x for a Normal 190 // distribution the given mean and Cholesky decomposition of the covariance matrix. 191 // NormalLogProb panics if len(x) is not equal to len(mu), or if len(mu) != chol.Size(). 192 // 193 // This function saves time and memory if the Cholesky decomposition is already 194 // available. Otherwise, the NewNormal function should be used. 195 func NormalLogProb(x, mu []float64, chol *mat.Cholesky) float64 { 196 dim := len(mu) 197 if len(x) != dim { 198 panic(badSizeMismatch) 199 } 200 if chol.SymmetricDim() != dim { 201 panic(badSizeMismatch) 202 } 203 logSqrtDet := 0.5 * chol.LogDet() 204 return normalLogProb(x, mu, chol, logSqrtDet) 205 } 206 207 // normalLogProb is the same as NormalLogProb, but does not make size checks and 208 // additionally requires log(|Σ|^-0.5) 209 func normalLogProb(x, mu []float64, chol *mat.Cholesky, logSqrtDet float64) float64 { 210 dim := len(mu) 211 c := -0.5*float64(dim)*logTwoPi - logSqrtDet 212 dst := stat.Mahalanobis(mat.NewVecDense(dim, x), mat.NewVecDense(dim, mu), chol) 213 return c - 0.5*dst*dst 214 } 215 216 // MarginalNormal returns the marginal distribution of the given input variables. 217 // That is, MarginalNormal returns 218 // p(x_i) = \int_{x_o} p(x_i | x_o) p(x_o) dx_o 219 // where x_i are the dimensions in the input, and x_o are the remaining dimensions. 220 // See https://en.wikipedia.org/wiki/Marginal_distribution for more information. 221 // 222 // The input src is passed to the call to NewNormal. 223 func (n *Normal) MarginalNormal(vars []int, src rand.Source) (*Normal, bool) { 224 newMean := make([]float64, len(vars)) 225 for i, v := range vars { 226 newMean[i] = n.mu[v] 227 } 228 var s mat.SymDense 229 s.SubsetSym(&n.sigma, vars) 230 return NewNormal(newMean, &s, src) 231 } 232 233 // MarginalNormalSingle returns the marginal of the given input variable. 234 // That is, MarginalNormal returns 235 // p(x_i) = \int_{x_¬i} p(x_i | x_¬i) p(x_¬i) dx_¬i 236 // where i is the input index. 237 // See https://en.wikipedia.org/wiki/Marginal_distribution for more information. 238 // 239 // The input src is passed to the constructed distuv.Normal. 240 func (n *Normal) MarginalNormalSingle(i int, src rand.Source) distuv.Normal { 241 return distuv.Normal{ 242 Mu: n.mu[i], 243 Sigma: math.Sqrt(n.sigma.At(i, i)), 244 Src: src, 245 } 246 } 247 248 // Mean returns the mean of the probability distribution at x. If the 249 // input argument is nil, a new slice will be allocated, otherwise the result 250 // will be put in-place into the receiver. 251 func (n *Normal) Mean(x []float64) []float64 { 252 x = reuseAs(x, n.dim) 253 copy(x, n.mu) 254 return x 255 } 256 257 // Prob computes the value of the probability density function at x. 258 func (n *Normal) Prob(x []float64) float64 { 259 return math.Exp(n.LogProb(x)) 260 } 261 262 // Quantile returns the multi-dimensional inverse cumulative distribution function. 263 // If x is nil, a new slice will be allocated and returned. If x is non-nil, 264 // len(x) must equal len(p) and the quantile will be stored in-place into x. 265 // All of the values of p must be between 0 and 1, inclusive, or Quantile will panic. 266 func (n *Normal) Quantile(x, p []float64) []float64 { 267 dim := n.Dim() 268 if len(p) != dim { 269 panic(badInputLength) 270 } 271 if x == nil { 272 x = make([]float64, dim) 273 } 274 if len(x) != len(p) { 275 panic(badInputLength) 276 } 277 278 // Transform to a standard normal and then transform to a multivariate Gaussian. 279 tmp := make([]float64, len(x)) 280 for i, v := range p { 281 tmp[i] = distuv.UnitNormal.Quantile(v) 282 } 283 n.TransformNormal(x, tmp) 284 return x 285 } 286 287 // Rand generates a random number according to the distributon. 288 // If the input slice is nil, new memory is allocated, otherwise the result is stored 289 // in place. 290 func (n *Normal) Rand(x []float64) []float64 { 291 return NormalRand(x, n.mu, &n.chol, n.src) 292 } 293 294 // NormalRand generates a random number with the given mean and Cholesky 295 // decomposition of the covariance matrix. 296 // If x is nil, new memory is allocated and returned, otherwise the result is stored 297 // in place into x. NormalRand panics if x is non-nil and not equal to len(mu), 298 // or if len(mu) != chol.Size(). 299 // 300 // This function saves time and memory if the Cholesky decomposition is already 301 // available. Otherwise, the NewNormal function should be used. 302 func NormalRand(x, mean []float64, chol *mat.Cholesky, src rand.Source) []float64 { 303 x = reuseAs(x, len(mean)) 304 if len(mean) != chol.SymmetricDim() { 305 panic(badInputLength) 306 } 307 if src == nil { 308 for i := range x { 309 x[i] = rand.NormFloat64() 310 } 311 } else { 312 rnd := rand.New(src) 313 for i := range x { 314 x[i] = rnd.NormFloat64() 315 } 316 } 317 transformNormal(x, x, mean, chol) 318 return x 319 } 320 321 // ScoreInput returns the gradient of the log-probability with respect to the 322 // input x. That is, ScoreInput computes 323 // ∇_x log(p(x)) 324 // If score is nil, a new slice will be allocated and returned. If score is of 325 // length the dimension of Normal, then the result will be put in-place into score. 326 // If neither of these is true, ScoreInput will panic. 327 func (n *Normal) ScoreInput(score, x []float64) []float64 { 328 // Normal log probability is 329 // c - 0.5*(x-μ)' Σ^-1 (x-μ). 330 // So the derivative is just 331 // -Σ^-1 (x-μ). 332 if len(x) != n.Dim() { 333 panic(badInputLength) 334 } 335 if score == nil { 336 score = make([]float64, len(x)) 337 } 338 if len(score) != len(x) { 339 panic(badSizeMismatch) 340 } 341 tmp := make([]float64, len(x)) 342 copy(tmp, x) 343 floats.Sub(tmp, n.mu) 344 345 err := n.chol.SolveVecTo(mat.NewVecDense(len(score), score), mat.NewVecDense(len(tmp), tmp)) 346 if err != nil { 347 panic(err) 348 } 349 floats.Scale(-1, score) 350 return score 351 } 352 353 // SetMean changes the mean of the normal distribution. SetMean panics if len(mu) 354 // does not equal the dimension of the normal distribution. 355 func (n *Normal) SetMean(mu []float64) { 356 if len(mu) != n.Dim() { 357 panic(badSizeMismatch) 358 } 359 copy(n.mu, mu) 360 } 361 362 // TransformNormal transforms the vector, normal, generated from a standard 363 // multidimensional normal into a vector that has been generated under the 364 // distribution of the receiver. 365 // 366 // If dst is non-nil, the result will be stored into dst, otherwise a new slice 367 // will be allocated. TransformNormal will panic if the length of normal is not 368 // the dimension of the receiver, or if dst is non-nil and len(dist) != len(normal). 369 func (n *Normal) TransformNormal(dst, normal []float64) []float64 { 370 if len(normal) != n.dim { 371 panic(badInputLength) 372 } 373 dst = reuseAs(dst, n.dim) 374 if len(dst) != len(normal) { 375 panic(badInputLength) 376 } 377 transformNormal(dst, normal, n.mu, &n.chol) 378 return dst 379 } 380 381 // transformNormal performs the same operation as Normal.TransformNormal except 382 // no safety checks are performed and all memory must be provided. 383 func transformNormal(dst, normal, mu []float64, chol *mat.Cholesky) []float64 { 384 dim := len(mu) 385 dstVec := mat.NewVecDense(dim, dst) 386 srcVec := mat.NewVecDense(dim, normal) 387 // If dst and normal are the same slice, make them the same Vector otherwise 388 // mat complains about being tricky. 389 if &normal[0] == &dst[0] { 390 srcVec = dstVec 391 } 392 dstVec.MulVec(chol.RawU().T(), srcVec) 393 floats.Add(dst, mu) 394 return dst 395 }