![]() |
Stan
2.14.0
probability, sampling & optimization
|
Variational family approximation with mean-field (diagonal covariance) multivariate normal distribution. More...
#include <normal_meanfield.hpp>
Public Member Functions | |
normal_meanfield (size_t dimension) | |
Construct a variational distribution of the specified dimensionality with a zero mean and zero log standard deviation (unit standard deviation). More... | |
normal_meanfield (const Eigen::VectorXd &cont_params) | |
Construct a variational distribution with the specified mean vector and zero log standard deviation (unit standard deviation). More... | |
normal_meanfield (const Eigen::VectorXd &mu, const Eigen::VectorXd &omega) | |
Construct a variational distribution with the specified mean and log standard deviation vectors. More... | |
int | dimension () const |
Return the dimensionality of the approximation. More... | |
const Eigen::VectorXd & | mu () const |
Return the mean vector. More... | |
const Eigen::VectorXd & | omega () const |
Return the log standard deviation vector. More... | |
void | set_mu (const Eigen::VectorXd &mu) |
Set the mean vector to the specified value. More... | |
void | set_omega (const Eigen::VectorXd &omega) |
Set the log standard deviation vector to the specified value. More... | |
void | set_to_zero () |
Sets the mean and log standard deviation vector for this approximation to zero. More... | |
normal_meanfield | square () const |
Return a new mean field approximation resulting from squaring the entries in the mean and log standard deviation. More... | |
normal_meanfield | sqrt () const |
Return a new mean field approximation resulting from taking the square root of the entries in the mean and log standard deviation. More... | |
normal_meanfield & | operator= (const normal_meanfield &rhs) |
Return this approximation after setting its mean vector and Cholesky factor for covariance to the values given by the specified approximation. More... | |
normal_meanfield & | operator+= (const normal_meanfield &rhs) |
Add the mean and Cholesky factor of the covariance matrix of the specified approximation to this approximation. More... | |
normal_meanfield & | operator/= (const normal_meanfield &rhs) |
Return this approximation after elementwise division by the specified approximation's mean and log standard deviation vectors. More... | |
normal_meanfield & | operator+= (double scalar) |
Return this approximation after adding the specified scalar to each entry in the mean and log standard deviation vectors. More... | |
normal_meanfield & | operator*= (double scalar) |
Return this approximation after multiplying by the specified scalar to each entry in the mean and log standard deviation vectors. More... | |
const Eigen::VectorXd & | mean () const |
Returns the mean vector for this approximation. More... | |
double | entropy () const |
Return the entropy of the approximation. More... | |
Eigen::VectorXd | transform (const Eigen::VectorXd &eta) const |
Return the transform of the sepcified vector using the Cholesky factor and mean vector. More... | |
template<class BaseRNG > | |
void | sample (BaseRNG &rng, Eigen::VectorXd &eta) const |
Assign a draw from this mean field approximation to the specified vector using the specified random number generator. More... | |
template<class M , class BaseRNG > | |
void | calc_grad (normal_meanfield &elbo_grad, M &m, Eigen::VectorXd &cont_params, int n_monte_carlo_grad, BaseRNG &rng, interface_callbacks::writer::base_writer &message_writer) const |
Calculates the "blackbox" gradient with respect to both the location vector (mu) and the log-std vector (omega) in parallel. More... | |
![]() | |
base_family () | |
base_family | square () const |
base_family | sqrt () const |
base_family | operator= (const base_family &rhs) |
base_family | operator+= (const base_family &rhs) |
base_family | operator/= (const base_family &rhs) |
base_family | operator+= (double scalar) |
base_family | operator*= (double scalar) |
const Eigen::VectorXd & | mean () const |
double | entropy () const |
Eigen::VectorXd | transform (const Eigen::VectorXd &eta) const |
template<class BaseRNG > | |
void | sample (BaseRNG &rng, Eigen::VectorXd &eta) const |
template<class M , class BaseRNG > | |
void | calc_grad (base_family &elbo_grad, M &m, Eigen::VectorXd &cont_params, int n_monte_carlo_grad, BaseRNG &rng, interface_callbacks::writer::base_writer &message_writer) const |
Additional Inherited Members | |
![]() | |
void | write_error_msg_ (std::ostream *error_msgs, const std::exception &e) const |
Variational family approximation with mean-field (diagonal covariance) multivariate normal distribution.
Definition at line 20 of file normal_meanfield.hpp.
|
inlineexplicit |
Construct a variational distribution of the specified dimensionality with a zero mean and zero log standard deviation (unit standard deviation).
[in] | dimension | Dimensionality of distribution. |
Definition at line 45 of file normal_meanfield.hpp.
|
inlineexplicit |
Construct a variational distribution with the specified mean vector and zero log standard deviation (unit standard deviation).
[in] | cont_params | Mean vector. |
Definition at line 58 of file normal_meanfield.hpp.
|
inline |
Construct a variational distribution with the specified mean and log standard deviation vectors.
[in] | mu | Mean vector. |
[in] | omega | Log standard deviation vector. |
std::domain_error | If the sizes of mean and log standard deviation vectors are different, or if either contains a not-a-number value. |
Definition at line 74 of file normal_meanfield.hpp.
|
inline |
Calculates the "blackbox" gradient with respect to both the location vector (mu) and the log-std vector (omega) in parallel.
It uses the same gradient computed from a set of Monte Carlo samples.
M | Model class. |
BaseRNG | Class of base random number generator. |
[in] | elbo_grad | Parameters to store "blackbox" gradient |
[in] | m | Model. |
[in] | cont_params | Continuous parameters. |
[in] | n_monte_carlo_grad | Number of samples for gradient computation. |
[in,out] | rng | Random number generator. |
[in,out] | message_writer | writer for messages |
std::domain_error | If the number of divergent iterations exceeds its specified bounds. |
Definition at line 363 of file normal_meanfield.hpp.
|
inline |
Return the dimensionality of the approximation.
Definition at line 89 of file normal_meanfield.hpp.
|
inline |
Return the entropy of the approximation.
The entropy is defined by 0.5 * dim * (1+log2pi) + 0.5 * log det diag(sigma^2) = 0.5 * dim * (1+log2pi) + sum(log(sigma)) = 0.5 * dim * (1+log2pi) + sum(omega)
Definition at line 298 of file normal_meanfield.hpp.
|
inline |
Returns the mean vector for this approximation.
See: mu()
.
Definition at line 284 of file normal_meanfield.hpp.
|
inline |
Return the mean vector.
Definition at line 94 of file normal_meanfield.hpp.
|
inline |
Return the log standard deviation vector.
Definition at line 99 of file normal_meanfield.hpp.
|
inline |
Return this approximation after multiplying by the specified scalar to each entry in the mean and log standard deviation vectors.
Warning: No finiteness check is made on the scalar, so it may introduce NaNs.
[in] | scalar | Scalar to add. |
Definition at line 271 of file normal_meanfield.hpp.
|
inline |
Add the mean and Cholesky factor of the covariance matrix of the specified approximation to this approximation.
[in] | rhs | Approximation from which to gather the mean and log standard deviation vectors. |
std::domain_error | If the size of the specified approximation does not match the size of this approximation. |
Definition at line 208 of file normal_meanfield.hpp.
|
inline |
Return this approximation after adding the specified scalar to each entry in the mean and log standard deviation vectors.
Warning: No finiteness check is made on the scalar, so it may introduce NaNs.
[in] | scalar | Scalar to add. |
Definition at line 253 of file normal_meanfield.hpp.
|
inline |
Return this approximation after elementwise division by the specified approximation's mean and log standard deviation vectors.
[in] | rhs | Approximation from which to gather the mean and log standard deviation vectors. |
std::domain_error | If the dimensionality of the specified approximation does not match this approximation's dimensionality. |
Definition at line 231 of file normal_meanfield.hpp.
|
inline |
Return this approximation after setting its mean vector and Cholesky factor for covariance to the values given by the specified approximation.
[in] | rhs | Approximation from which to gather the mean and log standard deviation vectors. |
std::domain_error | If the dimensionality of the specified approximation does not match this approximation's dimensionality. |
Definition at line 186 of file normal_meanfield.hpp.
|
inline |
Assign a draw from this mean field approximation to the specified vector using the specified random number generator.
BaseRNG | Class of random number generator. |
[in] | rng | Base random number generator. |
[in,out] | eta | Vector to which the draw is assigned; must already be properly sized. |
std::range_error | If the index is out of range. |
Definition at line 337 of file normal_meanfield.hpp.
|
inline |
Set the mean vector to the specified value.
[in] | mu | Mean vector. |
std::domain_error | If the mean vector's size does not match this approximation's dimensionality, or if it contains not-a-number values. |
Definition at line 109 of file normal_meanfield.hpp.
|
inline |
Set the log standard deviation vector to the specified value.
[in] | omega | Log standard deviation vector. |
std::domain_error | If the log standard deviation vector's size does not match this approximation's dimensionality, or if it contains not-a-number values. |
Definition at line 129 of file normal_meanfield.hpp.
|
inline |
Sets the mean and log standard deviation vector for this approximation to zero.
Definition at line 144 of file normal_meanfield.hpp.
|
inline |
Return a new mean field approximation resulting from taking the square root of the entries in the mean and log standard deviation.
The new approximation does not hold any references to this approximation.
Warning: No checks are carried out to ensure the entries are non-negative before taking square roots, so not-a-number values may result.
Definition at line 170 of file normal_meanfield.hpp.
|
inline |
Return a new mean field approximation resulting from squaring the entries in the mean and log standard deviation.
The new approximation does not hold any references to this approximation.
Definition at line 155 of file normal_meanfield.hpp.
|
inline |
Return the transform of the sepcified vector using the Cholesky factor and mean vector.
The transform is defined by S^{-1}(eta) = sigma * eta + mu = exp(omega) * eta + mu.
[in] | eta | Vector to transform. |
std::domain_error | If the specified vector's size does not match the dimensionality of this approximation. |
Definition at line 315 of file normal_meanfield.hpp.