Struct rv::dist::MvGaussian
source ยท pub struct MvGaussian { /* private fields */ }
Expand description
Multivariate Gaussian/Normal Distribution, ๐ฉ(ฮผ, ฮฃ).
Example
Generate a Wishart random 3x3 matrix ฮฃ ~ Wฮฝ(S)
use nalgebra::{DMatrix, DVector};
use rv::prelude::*;
let k = 3; // number of dimensions
let df = 6; // The degrees of freedom is an unsigned int > k
// The scale matrix of the wishar distribution
let scale_mat: DMatrix<f64> = DMatrix::identity(k, k);
// The draw procedure outlined in the appendices of "Bayesian Data
// Analysis" by Andrew Gelman and colleagues.
let mut rng = rand::thread_rng();
// 1. Create a multivariate normal with zero mean and covariance matrix S
let mvg = MvGaussian::new(DVector::zeros(k), scale_mat).unwrap();
// 2. Draw ฮฝ (df) vectors {x_1, ..., x_ฮฝ}
let xs = mvg.sample(df, &mut rng);
// 3. Compute the sum ฮฃ xx'
let mat = xs
.iter()
.fold(DMatrix::<f64>::zeros(k, k), |acc, x: &DVector<f64>| {
acc +x*x.transpose()
});
// Check that the matrix is square and has the right size
assert_eq!(mat.nrows(), k);
assert_eq!(mat.ncols(), k);
// check that the matrix is positive definite by attempting a Cholesky
// decomposition.
assert!(mat.cholesky().is_some());
Implementationsยง
sourceยงimpl MvGaussian
impl MvGaussian
sourcepub fn new(mu: DVector<f64>, cov: DMatrix<f64>) -> Result<Self, MvGaussianError>
pub fn new(mu: DVector<f64>, cov: DMatrix<f64>) -> Result<Self, MvGaussianError>
Create a new multivariate Gaussian distribution
Arguments
- mu: k-length mean vector
- cov: k-by-k positive-definite covariance matrix
sourcepub fn new_cholesky(
mu: DVector<f64>,
cov_chol: Cholesky<f64, Dyn>
) -> Result<Self, MvGaussianError>
pub fn new_cholesky( mu: DVector<f64>, cov_chol: Cholesky<f64, Dyn> ) -> Result<Self, MvGaussianError>
Create a new multivariate Gaussian distribution from Cholesky factorized Dov
Arguments
- mu: k-length mean vector
- cov_chol: Choleksy decomposition of k-by-k positive-definite covariance matrix
use nalgebra::{DMatrix, DVector};
use rv::prelude::*;
let mu = DVector::zeros(3);
let cov = DMatrix::from_row_slice(3, 3, &[
2.0, 1.0, 0.0,
1.0, 2.0, 0.0,
0.0, 0.0, 1.0
]);
let chol = cov.clone().cholesky().unwrap();
let mvg_r = MvGaussian::new_cholesky(mu, chol);
assert!(mvg_r.is_ok());
let mvg = mvg_r.unwrap();
assert!(cov.relative_eq(mvg.cov(), 1E-8, 1E-8));
sourcepub fn new_unchecked(mu: DVector<f64>, cov: DMatrix<f64>) -> Self
pub fn new_unchecked(mu: DVector<f64>, cov: DMatrix<f64>) -> Self
Creates a new MvGaussian from mean and covariance without checking whether the parameters are valid.
sourcepub fn new_cholesky_unchecked(
mu: DVector<f64>,
cov_chol: Cholesky<f64, Dyn>
) -> Self
pub fn new_cholesky_unchecked( mu: DVector<f64>, cov_chol: Cholesky<f64, Dyn> ) -> Self
Creates a new MvGaussian from mean and covarianceโs Cholesky factorization without checking whether the parameters are valid.
sourcepub fn standard(dims: usize) -> Result<Self, MvGaussianError>
pub fn standard(dims: usize) -> Result<Self, MvGaussianError>
Create a standard Gaussian distribution with zero mean and identity covariance matrix.
sourcepub fn ndims(&self) -> usize
pub fn ndims(&self) -> usize
Get the number of dimensions
Example
let mvg = MvGaussian::standard(4).unwrap();
assert_eq!(mvg.ndims(), 4);
sourcepub fn set_mu(&mut self, mu: DVector<f64>) -> Result<(), MvGaussianError>
pub fn set_mu(&mut self, mu: DVector<f64>) -> Result<(), MvGaussianError>
Set the mean
Example
let mut mvg = MvGaussian::standard(3).unwrap();
let x = DVector::<f64>::zeros(3);
assert::close(mvg.ln_f(&x), -2.756815599614018, 1E-8);
let cov_vals = vec![
1.01742788,
0.36586652,
-0.65620486,
0.36586652,
1.00564553,
-0.42597261,
-0.65620486,
-0.42597261,
1.27247972,
];
let cov: DMatrix<f64> = DMatrix::from_row_slice(3, 3, &cov_vals);
let mu = DVector::<f64>::from_column_slice(&vec![0.5, 3.1, -6.2]);
mvg.set_mu(mu).unwrap();
mvg.set_cov(cov).unwrap();
assert::close(mvg.ln_f(&x), -24.602370253215661, 1E-8);
pub fn set_mu_unchecked(&mut self, mu: DVector<f64>)
sourcepub fn set_cov(&mut self, cov: DMatrix<f64>) -> Result<(), MvGaussianError>
pub fn set_cov(&mut self, cov: DMatrix<f64>) -> Result<(), MvGaussianError>
Set the covariance matrix
Example
let mut mvg = MvGaussian::standard(3).unwrap();
let x = DVector::<f64>::zeros(3);
assert::close(mvg.ln_f(&x), -2.756815599614018, 1E-8);
let cov_vals = vec![
1.01742788,
0.36586652,
-0.65620486,
0.36586652,
1.00564553,
-0.42597261,
-0.65620486,
-0.42597261,
1.27247972,
];
let cov: DMatrix<f64> = DMatrix::from_row_slice(3, 3, &cov_vals);
let mu = DVector::<f64>::from_column_slice(&vec![0.5, 3.1, -6.2]);
mvg.set_mu(mu).unwrap();
mvg.set_cov(cov).unwrap();
assert::close(mvg.ln_f(&x), -24.602370253215661, 1E-8);
sourcepub fn set_cov_unchecked(&mut self, cov: DMatrix<f64>)
pub fn set_cov_unchecked(&mut self, cov: DMatrix<f64>)
Set the covariance matrix without input validation
Trait Implementationsยง
sourceยงimpl Clone for MvGaussian
impl Clone for MvGaussian
sourceยงfn clone(&self) -> MvGaussian
fn clone(&self) -> MvGaussian
Returns a copy of the value. Read more
1.0.0 ยท sourceยงfn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
Performs copy-assignment from
source
. Read moresourceยงimpl ConjugatePrior<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>, MvGaussian> for NormalInvWishart
impl ConjugatePrior<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>, MvGaussian> for NormalInvWishart
ยงtype Posterior = NormalInvWishart
type Posterior = NormalInvWishart
Type of the posterior distribution
ยงtype LnPpCache = (NormalInvWishart, f64)
type LnPpCache = (NormalInvWishart, f64)
Type of the
ln_pp
cachesourceยงfn posterior(
&self,
x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian>
) -> NormalInvWishart
fn posterior( &self, x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian> ) -> NormalInvWishart
Computes the posterior distribution from the data
sourceยงfn ln_m_cache(&self) -> f64
fn ln_m_cache(&self) -> f64
Compute the cache for the log marginal likelihood.
sourceยงfn ln_m_with_cache(
&self,
cache: &Self::LnMCache,
x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian>
) -> f64
fn ln_m_with_cache( &self, cache: &Self::LnMCache, x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian> ) -> f64
Log marginal likelihood with supplied cache.
sourceยงfn ln_pp_cache(
&self,
x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian>
) -> Self::LnPpCache
fn ln_pp_cache( &self, x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian> ) -> Self::LnPpCache
Compute the cache for the Log posterior predictive of y given x. Read more
sourceยงfn ln_pp_with_cache(&self, cache: &Self::LnPpCache, y: &DVector<f64>) -> f64
fn ln_pp_with_cache(&self, cache: &Self::LnPpCache, y: &DVector<f64>) -> f64
Log posterior predictive of y given x with supplied ln(norm)
sourceยงfn ln_m(&self, x: &DataOrSuffStat<'_, X, Fx>) -> f64
fn ln_m(&self, x: &DataOrSuffStat<'_, X, Fx>) -> f64
The log marginal likelihood
sourceยงfn ln_pp(&self, y: &X, x: &DataOrSuffStat<'_, X, Fx>) -> f64
fn ln_pp(&self, y: &X, x: &DataOrSuffStat<'_, X, Fx>) -> f64
Log posterior predictive of y given x
sourceยงfn m(&self, x: &DataOrSuffStat<'_, X, Fx>) -> f64
fn m(&self, x: &DataOrSuffStat<'_, X, Fx>) -> f64
Marginal likelihood of x
sourceยงimpl ContinuousDistr<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian
impl ContinuousDistr<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian
sourceยงimpl Debug for MvGaussian
impl Debug for MvGaussian
sourceยงimpl<'de> Deserialize<'de> for MvGaussian
impl<'de> Deserialize<'de> for MvGaussian
sourceยงfn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where __D: Deserializer<'de>,
Deserialize this value from the given Serde deserializer. Read more
sourceยงimpl Display for MvGaussian
impl Display for MvGaussian
sourceยงimpl From<&MvGaussian> for String
impl From<&MvGaussian> for String
sourceยงfn from(mvg: &MvGaussian) -> String
fn from(mvg: &MvGaussian) -> String
Converts to this type from the input type.
sourceยงimpl HasSuffStat<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian
impl HasSuffStat<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian
sourceยงimpl PartialEq<MvGaussian> for MvGaussian
impl PartialEq<MvGaussian> for MvGaussian
sourceยงfn eq(&self, other: &MvGaussian) -> bool
fn eq(&self, other: &MvGaussian) -> bool
This method tests for
self
and other
values to be equal, and is used
by ==
.sourceยงimpl Rv<MvGaussian> for NormalInvWishart
impl Rv<MvGaussian> for NormalInvWishart
sourceยงimpl Serialize for MvGaussian
impl Serialize for MvGaussian
sourceยงimpl Support<MvGaussian> for NormalInvWishart
impl Support<MvGaussian> for NormalInvWishart
Auto Trait Implementationsยง
impl RefUnwindSafe for MvGaussian
impl Send for MvGaussian
impl Sync for MvGaussian
impl Unpin for MvGaussian
impl UnwindSafe for MvGaussian
Blanket Implementationsยง
sourceยงimpl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere T: ?Sized,
sourceยงfn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
sourceยงimpl<Fx> Rv<Datum> for Fxwhere
Fx: RvDatum,
impl<Fx> Rv<Datum> for Fxwhere Fx: RvDatum,
ยงimpl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere SS: SubsetOf<SP>,
ยงfn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
The inverse inclusion map: attempts to construct
self
from the equivalent element of its
superset. Read moreยงfn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
Checks if
self
is actually part of its subset T
(and can be converted to it).ยงfn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
Use with care! Same as
self.to_subset
but without any property checks. Always succeeds.ยงfn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
The inclusion map: converts
self
to the equivalent element of its superset.