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