Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
stan::variational::normal_fullrank Class Reference

#include "stan/variational/families/normal_fullrank.hpp"

Inheritance diagram for stan::variational::normal_fullrank:
stan::variational::base_family

Public Member Functions

 normal_fullrank (size_t dimension)
 
 normal_fullrank (const Eigen::VectorXd &cont_params)
 
 normal_fullrank (const Eigen::VectorXd &mu, const Eigen::MatrixXd &L_chol)
 
int dimension () const
 
const Eigen::VectorXd & mu () const
 
const Eigen::MatrixXd & L_chol () const
 
void set_mu (const Eigen::VectorXd &mu)
 
void set_L_chol (const Eigen::MatrixXd &L_chol)
 
void set_to_zero ()
 
normal_fullrank square () const
 
normal_fullrank sqrt () const
 
normal_fullrankoperator= (const normal_fullrank &rhs)
 
normal_fullrankoperator+= (const normal_fullrank &rhs)
 
normal_fullrankoperator/= (const normal_fullrank &rhs)
 
normal_fullrankoperator+= (double scalar)
 
normal_fullrankoperator*= (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 (normal_fullrank &elbo_grad, M &m, Eigen::VectorXd &cont_params, int n_monte_carlo_grad, BaseRNG &rng, callbacks::logger &logger) const
 
base_family operator+= (const base_family &rhs)
 
base_family operator/= (const base_family &rhs)
 
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, callbacks::logger &logger) const
 

Protected Member Functions

void write_error_msg_ (std::ostream *error_msgs, const std::exception &e) const
 

Private Member Functions

void validate_mean (const char *function, const Eigen::VectorXd &mu)
 
void validate_cholesky_factor (const char *function, const Eigen::MatrixXd &L_chol)
 

Private Attributes

Eigen::VectorXd mu_
 
Eigen::MatrixXd L_chol_
 
const int dimension_
 

Detailed Description

Variational family approximation with full-rank multivariate normal distribution.

Definition at line 20 of file normal_fullrank.hpp.

Constructor & Destructor Documentation

stan::variational::normal_fullrank::normal_fullrank ( size_t  dimension)
inlineexplicit

Construct a variational distribution of the specified dimensionality with a zero mean and Cholesky factor of a zero covariance matrix.

Parameters
[in]dimensionDimensionality of distribution.

Definition at line 88 of file normal_fullrank.hpp.

Referenced by sqrt(), and square().

stan::variational::normal_fullrank::normal_fullrank ( const Eigen::VectorXd &  cont_params)
inlineexplicit

Construct a variational distribution with specified mean vector and Cholesky factor for identity covariance.

Parameters
[in]cont_paramsMean vector.

Definition at line 101 of file normal_fullrank.hpp.

102  : mu_(cont_params),
103  L_chol_(Eigen::MatrixXd::Identity(cont_params.size(),
104  cont_params.size())),
105  dimension_(cont_params.size()) {
106  }
stan::variational::normal_fullrank::normal_fullrank ( const Eigen::VectorXd &  mu,
const Eigen::MatrixXd &  L_chol 
)
inline

Construct a variational distribution with specified mean and Cholesky factor for covariance.

Warning: Positive-definiteness is not enforced for the Cholesky factor.

Parameters
[in]muMean vector.
[in]L_cholCholesky factor of covariance.
Exceptions
std::domain_errorIf the Cholesky factor is not square or not lower triangular, if the mean and Cholesky factor have different dimensionality, or if any of the elements is not-a-number.

Definition at line 122 of file normal_fullrank.hpp.

References validate_cholesky_factor(), and validate_mean().

124  : mu_(mu), L_chol_(L_chol), dimension_(mu.size()) {
125  static const char* function = "stan::variational::normal_fullrank";
126  validate_mean(function, mu);
127  validate_cholesky_factor(function, L_chol);
128  }
void validate_mean(const char *function, const Eigen::VectorXd &mu)
const Eigen::VectorXd & mu() const
void validate_cholesky_factor(const char *function, const Eigen::MatrixXd &L_chol)
const Eigen::MatrixXd & L_chol() const

Member Function Documentation

template<class M , class BaseRNG >
void stan::variational::base_family::calc_grad ( base_family elbo_grad,
M &  m,
Eigen::VectorXd &  cont_params,
int  n_monte_carlo_grad,
BaseRNG &  rng,
callbacks::logger logger 
) const
inherited
template<class M , class BaseRNG >
void stan::variational::normal_fullrank::calc_grad ( normal_fullrank elbo_grad,
M &  m,
Eigen::VectorXd &  cont_params,
int  n_monte_carlo_grad,
BaseRNG &  rng,
callbacks::logger logger 
) const
inline

Calculates the "blackbox" gradient with respect to BOTH the location vector (mu) and the cholesky factor of the scale matrix (L_chol) in parallel. It uses the same gradient computed from a set of Monte Carlo samples

Template Parameters
MModel class.
BaseRNGClass of base random number generator.
Parameters
[in]elbo_gradApproximation to store "blackbox" gradient.
[in]mModel.
[in]cont_paramsContinuous parameters.
[in]n_monte_carlo_gradSample size for gradient computation.
[in,out]rngRandom number generator.
[in,out]loggerlogger for messages
Exceptions
std::domain_errorIf the number of divergent iterations exceeds its specified bounds.

Definition at line 400 of file normal_fullrank.hpp.

References stan::math::check_finite(), stan::math::check_size_match(), d, dimension(), dimension_, stan::math::domain_error(), e, stan::model::gradient(), MECModelEnuComparisons::i, stan::callbacks::logger::info(), stan::math::normal_rng(), set_L_chol(), set_mu(), ss, transform(), submit_syst::y, and Zero().

406  {
407  static const char* function =
408  "stan::variational::normal_fullrank::calc_grad";
410  "Dimension of elbo_grad", elbo_grad.dimension(),
411  "Dimension of variational q", dimension_);
413  "Dimension of variational q", dimension_,
414  "Dimension of variables in model", cont_params.size());
415 
416  Eigen::VectorXd mu_grad = Eigen::VectorXd::Zero(dimension_);
417  Eigen::MatrixXd L_grad = Eigen::MatrixXd::Zero(dimension_, dimension_);
418  double tmp_lp = 0.0;
419  Eigen::VectorXd tmp_mu_grad = Eigen::VectorXd::Zero(dimension_);
420  Eigen::VectorXd eta = Eigen::VectorXd::Zero(dimension_);
421  Eigen::VectorXd zeta = Eigen::VectorXd::Zero(dimension_);
422 
423  // Naive Monte Carlo integration
424  static const int n_retries = 10;
425  for (int i = 0, n_monte_carlo_drop = 0; i < n_monte_carlo_grad; ) {
426  // Draw from standard normal and transform to real-coordinate space
427  for (int d = 0; d < dimension_; ++d) {
428  eta(d) = stan::math::normal_rng(0, 1, rng);
429  }
430  zeta = transform(eta);
431  try {
432  std::stringstream ss;
433  stan::model::gradient(m, zeta, tmp_lp, tmp_mu_grad, &ss);
434  if (ss.str().length() > 0)
435  logger.info(ss);
436  stan::math::check_finite(function, "Gradient of mu", tmp_mu_grad);
437 
438  mu_grad += tmp_mu_grad;
439  for (int ii = 0; ii < dimension_; ++ii) {
440  for (int jj = 0; jj <= ii; ++jj) {
441  L_grad(ii, jj) += tmp_mu_grad(ii) * eta(jj);
442  }
443  }
444  ++i;
445  } catch (const std::exception& e) {
446  ++n_monte_carlo_drop;
447  if (n_monte_carlo_drop >= n_retries * n_monte_carlo_grad) {
448  const char* name = "The number of dropped evaluations";
449  const char* msg1 = "has reached its maximum amount (";
450  int y = n_retries * n_monte_carlo_grad;
451  const char* msg2 = "). Your model may be either severely "
452  "ill-conditioned or misspecified.";
453  stan::math::domain_error(function, name, y, msg1, msg2);
454  }
455  }
456  }
457  mu_grad /= static_cast<double>(n_monte_carlo_grad);
458  L_grad /= static_cast<double>(n_monte_carlo_grad);
459 
460  // Add gradient of entropy term
461  L_grad.diagonal().array() += L_chol_.diagonal().array().inverse();
462 
463  elbo_grad.set_mu(mu_grad);
464  elbo_grad.set_L_chol(L_grad);
465  }
const XML_Char * name
Definition: expat.h:151
void check_finite(const char *function, const char *name, const T_y &y)
Eigen::VectorXd transform(const Eigen::VectorXd &eta) const
Float_t ss
Definition: plot.C:24
::xsd::cxx::tree::exception< char > exception
Definition: Database.h:225
void gradient(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, std::ostream *msgs=0)
Definition: gradient.hpp:15
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
VectorBuilder< true, double, T_loc, T_scale >::type normal_rng(const T_loc &mu, const T_scale &sigma, RNG &rng)
Definition: normal_rng.hpp:35
Float_t d
Definition: plot.C:236
void domain_error(const char *function, const char *name, const T &y, const char *msg1, const char *msg2)
void Zero()
Float_t e
Definition: plot.C:35
int stan::variational::normal_fullrank::dimension ( ) const
inline

Return the dimensionality of the approximation.

Definition at line 133 of file normal_fullrank.hpp.

References dimension_.

Referenced by calc_grad(), operator+=(), operator/=(), operator=(), and TEST().

double stan::variational::normal_fullrank::entropy ( ) const
inline

Return the entropy of this approximation.

The entropy is defined by 0.5 * dim * (1+log2pi) + 0.5 * log det (L^T L) = 0.5 * dim * (1+log2pi) + sum(log(abs(diag(L)))).

Returns
Entropy of this approximation

Definition at line 333 of file normal_fullrank.hpp.

References d, dimension_, stan::math::fabs(), L_chol_, test_ParserArtEvents::log, stan::math::LOG_TWO_PI, fillBadChanDBTables::result, and tmp.

Referenced by TEST().

333  {
334  static double mult = 0.5 * (1.0 + stan::math::LOG_TWO_PI);
335  double result = mult * dimension_;
336  for (int d = 0; d < dimension_; ++d) {
337  double tmp = fabs(L_chol_(d, d));
338  if (tmp != 0.0) result += log(tmp);
339  }
340  return result;
341  }
fvar< T > fabs(const fvar< T > &x)
Definition: fabs.hpp:15
Float_t tmp
Definition: plot.C:36
const double LOG_TWO_PI
Definition: constants.hpp:167
Float_t d
Definition: plot.C:236
const Eigen::MatrixXd& stan::variational::normal_fullrank::L_chol ( ) const
inline

Return the Cholesky factor of the covariance matrix.

Definition at line 143 of file normal_fullrank.hpp.

References L_chol_.

Referenced by operator+=(), operator/=(), operator=(), set_L_chol(), and TEST().

143 { return L_chol_; }
const Eigen::VectorXd& stan::variational::normal_fullrank::mean ( ) const
inline

Returns the mean vector for this approximation.

See: mu().

Returns
Mean vector for this approximation.

Definition at line 320 of file normal_fullrank.hpp.

References mu().

320  {
321  return mu();
322  }
const Eigen::VectorXd & mu() const
const Eigen::VectorXd& stan::variational::normal_fullrank::mu ( ) const
inline

Return the mean vector.

Definition at line 138 of file normal_fullrank.hpp.

References mu_.

Referenced by mean(), operator+=(), operator/=(), operator=(), set_mu(), and TEST().

138 { return mu_; }
normal_fullrank& stan::variational::normal_fullrank::operator*= ( double  scalar)
inline

Return this approximation after multiplying by the specified scalar to each entry in the mean and cholesky factor for covariance.

Warning: No finiteness check is made on the scalar, so it may introduce NaNs.

Parameters
[in]scalarScalar to add.
Returns
This approximation after elementwise addition of the specified scalar.

Definition at line 307 of file normal_fullrank.hpp.

307  {
308  mu_ *= scalar;
309  L_chol_ *= scalar;
310  return *this;
311  }
base_family stan::variational::base_family::operator+= ( const base_family rhs)
inherited
normal_fullrank& stan::variational::normal_fullrank::operator+= ( const normal_fullrank rhs)
inline

Add the mean and Cholesky factor of the covariance matrix of the specified approximation to this approximation.

Parameters
[in]rhsApproximation from which to gather the mean and covariance.
Returns
This approximation after adding the specified approximation.
Exceptions
std::domain_errorIf the dimensionality of the specified approximation does not match this approximation's dimensionality.

Definition at line 241 of file normal_fullrank.hpp.

References stan::math::check_size_match(), dimension(), L_chol(), and mu().

241  {
242  static const char* function =
243  "stan::variational::normal_fullrank::operator+=";
245  "Dimension of lhs", dimension_,
246  "Dimension of rhs", rhs.dimension());
247  mu_ += rhs.mu();
248  L_chol_ += rhs.L_chol();
249  return *this;
250  }
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
normal_fullrank& stan::variational::normal_fullrank::operator+= ( double  scalar)
inline

Return this approximation after adding the specified scalar to each entry in the mean and cholesky factor for covariance.

Warning: No finiteness check is made on the scalar, so it may introduce NaNs.

Parameters
[in]scalarScalar to add.
Returns
This approximation after elementwise addition of the specified scalar.

Definition at line 289 of file normal_fullrank.hpp.

289  {
290  mu_.array() += scalar;
291  L_chol_.array() += scalar;
292  return *this;
293  }
base_family stan::variational::base_family::operator/= ( const base_family rhs)
inherited
normal_fullrank& stan::variational::normal_fullrank::operator/= ( const normal_fullrank rhs)
inline

Return this approximation after elementwise division by the specified approximation's mean and Cholesky factor for covariance.

Parameters
[in]rhsApproximation from which to gather the mean and covariance.
Returns
This approximation after elementwise division by the specified approximation.
Exceptions
std::domain_errorIf the dimensionality of the specified approximation does not match this approximation's dimensionality.

Definition at line 265 of file normal_fullrank.hpp.

References stan::math::check_size_match(), dimension(), L_chol(), and mu().

265  {
266  static const char* function =
267  "stan::variational::normal_fullrank::operator/=";
268 
270  "Dimension of lhs", dimension_,
271  "Dimension of rhs", rhs.dimension());
272 
273  mu_.array() /= rhs.mu().array();
274  L_chol_.array() /= rhs.L_chol().array();
275  return *this;
276  }
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
normal_fullrank& stan::variational::normal_fullrank::operator= ( const normal_fullrank rhs)
inline

Return this approximation after setting its mean vector and Cholesky factor for covariance to the values given by the specified approximation.

Parameters
[in]rhsApproximation from which to gather the mean and covariance.
Returns
This approximation after assignment.
Exceptions
std::domain_errorIf the dimensionality of the specified approximation does not match this approximation's dimensionality.

Definition at line 219 of file normal_fullrank.hpp.

References stan::math::check_size_match(), dimension(), L_chol(), and mu().

219  {
220  static const char* function =
221  "stan::variational::normal_fullrank::operator=";
223  "Dimension of lhs", dimension_,
224  "Dimension of rhs", rhs.dimension());
225  mu_ = rhs.mu();
226  L_chol_ = rhs.L_chol();
227  return *this;
228  }
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
template<class BaseRNG >
void stan::variational::normal_fullrank::sample ( BaseRNG &  rng,
Eigen::VectorXd &  eta 
) const
inline

Set the specified vector to a draw from this variational approximation using the specified random number generator.

Template Parameters
BaseRNGClass of random number generator.
Parameters
[in,out]rngBase random number generator.
[out]etaRandom draw.
Returns
A sample from the variational distribution.

Definition at line 376 of file normal_fullrank.hpp.

References d, dimension_, stan::math::normal_rng(), and transform().

376  {
377  for (int d = 0; d < dimension_; ++d)
378  eta(d) = stan::math::normal_rng(0, 1, rng);
379  eta = transform(eta);
380  }
Eigen::VectorXd transform(const Eigen::VectorXd &eta) const
VectorBuilder< true, double, T_loc, T_scale >::type normal_rng(const T_loc &mu, const T_scale &sigma, RNG &rng)
Definition: normal_rng.hpp:35
Float_t d
Definition: plot.C:236
void stan::variational::normal_fullrank::set_L_chol ( const Eigen::MatrixXd &  L_chol)
inline

Set the Cholesky factor to the specified value.

Parameters
[in]L_cholCholesky factor of covariance matrix.
Exceptions
std::domain_errorIf the specified matrix is not square, is not lower triangular, if its size does not match the dimensionality of this approximation, or if it contains not-a-number values.

Definition at line 167 of file normal_fullrank.hpp.

References L_chol(), and validate_cholesky_factor().

Referenced by calc_grad(), and TEST().

167  {
168  static const char* function = "stan::variational::set_L_chol";
169  validate_cholesky_factor(function, L_chol);
170  L_chol_ = L_chol;
171  }
void validate_cholesky_factor(const char *function, const Eigen::MatrixXd &L_chol)
const Eigen::MatrixXd & L_chol() const
void stan::variational::normal_fullrank::set_mu ( const Eigen::VectorXd &  mu)
inline

Set the mean vector to the specified value.

Parameters
[in]muMean vector.
Exceptions
std::domain_errorIf the size of the specified mean vector does not match the stored dimension of this approximation.

Definition at line 152 of file normal_fullrank.hpp.

References mu(), and validate_mean().

Referenced by calc_grad(), and TEST().

152  {
153  static const char* function = "stan::variational::set_mu";
154  validate_mean(function, mu);
155  mu_ = mu;
156  }
void validate_mean(const char *function, const Eigen::VectorXd &mu)
const Eigen::VectorXd & mu() const
void stan::variational::normal_fullrank::set_to_zero ( )
inline

Set the mean vector and Cholesky factor for the covariance matrix to zero.

Definition at line 177 of file normal_fullrank.hpp.

References Zero().

Referenced by TEST().

normal_fullrank stan::variational::normal_fullrank::sqrt ( ) const
inline

Return a new full rank approximation resulting from taking the square root of the entries in the mean and Cholesky factor for the covariance matrix. 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 203 of file normal_fullrank.hpp.

References normal_fullrank().

203  {
204  return normal_fullrank(Eigen::VectorXd(mu_.array().sqrt()),
205  Eigen::MatrixXd(L_chol_.array().sqrt()));
206  }
normal_fullrank stan::variational::normal_fullrank::square ( ) const
inline

Return a new full rank approximation resulting from squaring the entries in the mean and Cholesky factor for the covariance matrix. The new approximation does not hold any references to this approximation.

Definition at line 188 of file normal_fullrank.hpp.

References normal_fullrank().

188  {
189  return normal_fullrank(Eigen::VectorXd(mu_.array().square()),
190  Eigen::MatrixXd(L_chol_.array().square()));
191  }
Eigen::VectorXd stan::variational::normal_fullrank::transform ( const Eigen::VectorXd &  eta) const
inline

Return the transform of the sepcified vector using the Cholesky factor and mean vector.

The transform is defined by S^{-1}(eta) = L_chol * eta + mu.

Parameters
[in]etaVector to transform.
Exceptions
std::domain_errorIf the specified vector's size does not match the dimensionality of this approximation.
Returns
Transformed vector.

Definition at line 355 of file normal_fullrank.hpp.

References stan::math::check_not_nan(), stan::math::check_size_match(), dimension_, and mu_.

Referenced by calc_grad(), sample(), and TEST().

355  {
356  static const char* function =
357  "stan::variational::normal_fullrank::transform";
359  "Dimension of input vector", eta.size(),
360  "Dimension of mean vector", dimension_);
361  stan::math::check_not_nan(function, "Input vector", eta);
362 
363  return (L_chol_ * eta) + mu_;
364  }
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
void check_not_nan(const char *function, const char *name, const T_y &y)
void stan::variational::normal_fullrank::validate_cholesky_factor ( const char *  function,
const Eigen::MatrixXd &  L_chol 
)
inlineprivate

Raise a domain exception if the specified matrix is not square, not lower triangular, or contains not-a-number values.

Warning: This function does not check that the Cholesky factor is positive definite.

Parameters
[in]L_cholCholesky factor for covariance matrix.
Exceptions
std::domain_errorIf the specified matrix is not square, is not lower triangular, if its size does not match the dimensionality of this approximation, or if it contains not-a-number values.

Definition at line 68 of file normal_fullrank.hpp.

References stan::math::check_lower_triangular(), stan::math::check_not_nan(), stan::math::check_size_match(), and stan::math::check_square().

Referenced by normal_fullrank(), and set_L_chol().

69  {
70  stan::math::check_square(function, "Cholesky factor", L_chol);
72  "Cholesky factor", L_chol);
74  "Dimension of mean vector", dimension_,
75  "Dimension of Cholesky factor", L_chol.rows());
76  stan::math::check_not_nan(function, "Cholesky factor", L_chol);
77  }
void check_lower_triangular(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
void check_not_nan(const char *function, const char *name, const T_y &y)
void check_square(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
const Eigen::MatrixXd & L_chol() const
void stan::variational::normal_fullrank::validate_mean ( const char *  function,
const Eigen::VectorXd &  mu 
)
inlineprivate

Raise a domain exception if the specified vector contains not-a-number values.

Parameters
[in]muMean vector.
Exceptions
std::domain_errorIf the mean vector contains NaN values or does not match this distribution's dimensionality.

Definition at line 46 of file normal_fullrank.hpp.

References stan::math::check_not_nan(), stan::math::check_size_match(), and dimension_.

Referenced by normal_fullrank(), and set_mu().

47  {
48  stan::math::check_not_nan(function, "Mean vector", mu);
50  "Dimension of input vector", mu.size(),
51  "Dimension of current vector", dimension_);
52  }
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
const Eigen::VectorXd & mu() const
void check_not_nan(const char *function, const char *name, const T_y &y)
void stan::variational::base_family::write_error_msg_ ( std::ostream *  error_msgs,
const std::exception &  e 
) const
inlineprotectedinherited

Definition at line 44 of file base_family.hpp.

References allTimeWatchdog::endl, stan::variational::operator*(), stan::variational::operator+(), and stan::variational::operator/().

45  {
46  if (!error_msgs) {
47  return;
48  }
49 
50  *error_msgs
51  << std::endl
52  << "Informational Message: The current gradient evaluation "
53  << "of the ELBO is ignored because of the following issue:"
54  << std::endl
55  << e.what() << std::endl
56  << "If this warning occurs often then your model may be "
57  << "either severely ill-conditioned or misspecified."
58  << std::endl;
59  }
Float_t e
Definition: plot.C:35

Member Data Documentation

const int stan::variational::normal_fullrank::dimension_
private

Dimensionality of distribution.

Definition at line 36 of file normal_fullrank.hpp.

Referenced by calc_grad(), dimension(), entropy(), sample(), transform(), and validate_mean().

Eigen::MatrixXd stan::variational::normal_fullrank::L_chol_
private

Cholesky factor of covariance: Sigma = L_chol * L_chol.transpose()

Definition at line 31 of file normal_fullrank.hpp.

Referenced by entropy(), and L_chol().

Eigen::VectorXd stan::variational::normal_fullrank::mu_
private

Mean vector.

Definition at line 25 of file normal_fullrank.hpp.

Referenced by mu(), and transform().


The documentation for this class was generated from the following file: