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

source

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
source

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));
source

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.

source

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.

source

pub fn standard(dims: usize) -> Result<Self, MvGaussianError>

Create a standard Gaussian distribution with zero mean and identity covariance matrix.

source

pub fn ndims(&self) -> usize

Get the number of dimensions

Example
let mvg = MvGaussian::standard(4).unwrap();
assert_eq!(mvg.ndims(), 4);
source

pub fn mu(&self) -> &DVector<f64>

Get a reference to the mean

source

pub fn cov(&self) -> &DMatrix<f64>

Get a reference to the covariance

source

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);
source

pub fn set_mu_unchecked(&mut self, mu: DVector<f64>)

source

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);
source

pub fn set_cov_unchecked(&mut self, cov: DMatrix<f64>)

Set the covariance matrix without input validation

Trait Implementationsยง

sourceยง

impl Clone for MvGaussian

sourceยง

fn clone(&self) -> MvGaussian

Returns a copy of the value. Read more
1.0.0 ยท sourceยง

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
sourceยง

impl ConjugatePrior<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>, MvGaussian> for NormalInvWishart

ยง

type Posterior = NormalInvWishart

Type of the posterior distribution
ยง

type LnMCache = f64

Type of the ln_m cache
ยง

type LnPpCache = (NormalInvWishart, f64)

Type of the ln_pp cache
sourceยง

fn posterior( &self, x: &DataOrSuffStat<'_, DVector<f64>, MvGaussian> ) -> NormalInvWishart

Computes the posterior distribution from the data
sourceยง

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

Log marginal likelihood with supplied cache.
sourceยง

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

Log posterior predictive of y given x with supplied ln(norm)
sourceยง

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

Log posterior predictive of y given x
sourceยง

fn m(&self, x: &DataOrSuffStat<'_, X, Fx>) -> f64

Marginal likelihood of x
sourceยง

fn pp(&self, y: &X, x: &DataOrSuffStat<'_, X, Fx>) -> f64

Posterior Predictive distribution
sourceยง

impl ContinuousDistr<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian

sourceยง

fn pdf(&self, x: &X) -> f64

The value of the Probability Density Function (PDF) at x Read more
sourceยง

fn ln_pdf(&self, x: &X) -> f64

The value of the log Probability Density Function (PDF) at x Read more
sourceยง

impl ContinuousDistr<MvGaussian> for NormalInvWishart

sourceยง

fn pdf(&self, x: &X) -> f64

The value of the Probability Density Function (PDF) at x Read more
sourceยง

fn ln_pdf(&self, x: &X) -> f64

The value of the log Probability Density Function (PDF) at x Read more
sourceยง

impl Debug for MvGaussian

sourceยง

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
sourceยง

impl<'de> Deserialize<'de> for MvGaussian

sourceยง

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

sourceยง

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
sourceยง

impl Entropy for MvGaussian

sourceยง

fn entropy(&self) -> f64

The entropy, H(X)
sourceยง

impl From<&MvGaussian> for String

sourceยง

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

ยง

type Stat = MvGaussianSuffStat

sourceยง

fn empty_suffstat(&self) -> Self::Stat

sourceยง

fn ln_f_stat(&self, stat: &Self::Stat) -> f64

Return the log likelihood for the data represented by the sufficient statistic.
sourceยง

impl Mean<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian

sourceยง

fn mean(&self) -> Option<DVector<f64>>

Returns None if the mean is undefined
sourceยง

impl Mode<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian

sourceยง

fn mode(&self) -> Option<DVector<f64>>

Returns None if the mode is undefined or is not a single value
sourceยง

impl PartialEq<MvGaussian> for MvGaussian

sourceยง

fn eq(&self, other: &MvGaussian) -> bool

This method tests for self and other values to be equal, and is used by ==.
1.0.0 ยท sourceยง

fn ne(&self, other: &Rhs) -> bool

This method tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
sourceยง

impl Rv<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian

sourceยง

fn ln_f(&self, x: &DVector<f64>) -> f64

Probability function Read more
sourceยง

fn draw<R: Rng>(&self, rng: &mut R) -> DVector<f64>

Single draw from the Rv Read more
sourceยง

fn f(&self, x: &X) -> f64

Probability function Read more
sourceยง

fn sample<R: Rng>(&self, n: usize, rng: &mut R) -> Vec<X>

Multiple draws of the Rv Read more
sourceยง

fn sample_stream<'r, R: Rng>( &'r self, rng: &'r mut R ) -> Box<dyn Iterator<Item = X> + 'r>

Create a never-ending iterator of samples Read more
sourceยง

impl Rv<MvGaussian> for NormalInvWishart

sourceยง

fn ln_f(&self, x: &MvGaussian) -> f64

Probability function Read more
sourceยง

fn draw<R: Rng>(&self, rng: &mut R) -> MvGaussian

Single draw from the Rv Read more
sourceยง

fn f(&self, x: &X) -> f64

Probability function Read more
sourceยง

fn sample<R: Rng>(&self, n: usize, rng: &mut R) -> Vec<X>

Multiple draws of the Rv Read more
sourceยง

fn sample_stream<'r, R: Rng>( &'r self, rng: &'r mut R ) -> Box<dyn Iterator<Item = X> + 'r>

Create a never-ending iterator of samples Read more
sourceยง

impl Serialize for MvGaussian

sourceยง

fn serialize<__S>(&self, __serializer: __S) -> Result<__S::Ok, __S::Error>where __S: Serializer,

Serialize this value into the given Serde serializer. Read more
sourceยง

impl Support<Matrix<f64, Dyn, Const<1>, VecStorage<f64, Dyn, Const<1>>>> for MvGaussian

sourceยง

fn supports(&self, x: &DVector<f64>) -> bool

Returns true if x is in the support of the Rv Read more
sourceยง

impl Support<MvGaussian> for NormalInvWishart

sourceยง

fn supports(&self, x: &MvGaussian) -> bool

Returns true if x is in the support of the Rv Read more
sourceยง

impl Variance<Matrix<f64, Dyn, Dyn, VecStorage<f64, Dyn, Dyn>>> for MvGaussian

sourceยง

fn variance(&self) -> Option<DMatrix<f64>>

Returns None if the variance is undefined

Auto Trait Implementationsยง

Blanket Implementationsยง

sourceยง

impl<T> Any for Twhere T: 'static + ?Sized,

sourceยง

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
sourceยง

impl<T> Borrow<T> for Twhere T: ?Sized,

sourceยง

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
sourceยง

impl<T> BorrowMut<T> for Twhere T: ?Sized,

sourceยง

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
sourceยง

impl<T> From<T> for T

sourceยง

fn from(t: T) -> T

Returns the argument unchanged.

sourceยง

impl<T, U> Into<U> for Twhere U: From<T>,

sourceยง

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

sourceยง

impl<Fx> Rv<Datum> for Fxwhere Fx: RvDatum,

sourceยง

fn ln_f(&self, x: &Datum) -> f64

Probability function Read more
sourceยง

fn draw<R>(&self, rng: &mut R) -> Datumwhere R: Rng,

Single draw from the Rv Read more
sourceยง

fn sample<R>(&self, n: usize, rng: &mut R) -> Vec<Datum, Global>where R: Rng,

Multiple draws of the Rv Read more
sourceยง

fn sample_stream<R, 'r>( &'r self, rng: &'r mut R ) -> Box<dyn Iterator<Item = Datum> + 'r, Global>where R: Rng,

Create a never-ending iterator of samples Read more
sourceยง

fn f(&self, x: &X) -> f64

Probability function Read more
sourceยง

impl<T> Same<T> for T

ยง

type Output = T

Should always be Self
ยง

impl<SS, SP> SupersetOf<SS> for SPwhere SS: SubsetOf<SP>,

ยง

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

Checks if self is actually part of its subset T (and can be converted to it).
ยง

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

The inclusion map: converts self to the equivalent element of its superset.
sourceยง

impl<T> ToOwned for Twhere T: Clone,

ยง

type Owned = T

The resulting type after obtaining ownership.
sourceยง

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
sourceยง

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
sourceยง

impl<T> ToString for Twhere T: Display + ?Sized,

sourceยง

default fn to_string(&self) -> String

Converts the given value to a String. Read more
sourceยง

impl<T, U> TryFrom<U> for Twhere U: Into<T>,

ยง

type Error = Infallible

The type returned in the event of a conversion error.
sourceยง

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
sourceยง

impl<T, U> TryInto<U> for Twhere U: TryFrom<T>,

ยง

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
sourceยง

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
ยง

impl<V, T> VZip<V> for Twhere V: MultiLane<T>,

ยง

fn vzip(self) -> V

sourceยง

impl<T> DeserializeOwned for Twhere T: for<'de> Deserialize<'de>,

sourceยง

impl<T> DeserializeOwnedAlias for Twhere T: DeserializeOwned,

sourceยง

impl<T> Scalar for Twhere T: 'static + Clone + PartialEq<T> + Debug,

sourceยง

impl<T> SendAlias for T

ยง

impl<T> SendSyncUnwindSafe for Twhere T: Send + Sync + UnwindSafe + ?Sized,

sourceยง

impl<T> SerializeAlias for Twhere T: Serialize,

sourceยง

impl<T> SyncAlias for T