Public Member Functions | Public Attributes | Protected Attributes | List of all members
o2scl::prob_dens_mdim_bound_gaussian< vec_t, mat_t > Class Template Reference

Gaussian distribution bounded by a hypercube. More...

#include <prob_dens_func.h>

Inheritance diagram for o2scl::prob_dens_mdim_bound_gaussian< vec_t, mat_t >:
o2scl::prob_dens_mdim_gaussian< vec_t, mat_t > o2scl::prob_dens_mdim< vec_t >

Detailed Description

template<class vec_t = boost::numeric::ublas::vector<double>, class mat_t = boost::numeric::ublas::matrix<double>>
class o2scl::prob_dens_mdim_bound_gaussian< vec_t, mat_t >

Note
This class naively resamples the Gaussian until a sample is within bounds. This is a temporary hack and can be very slow depending on the size of the volume excluded.
Warning
The PDF is not yet properly normalized

Definition at line 1314 of file prob_dens_func.h.

Public Member Functions

 prob_dens_mdim_bound_gaussian ()
 Create an empty distribution.
 
 prob_dens_mdim_bound_gaussian (size_t p_ndim, vec_t &p_peak, mat_t &covar, vec_t &p_low, vec_t &p_high)
 Create a distribution with the specified peak, covariance matrix, lower limits, and upper limits.
 
void set (size_t p_ndim, vec_t &p_peak, mat_t &covar, vec_t &p_low, vec_t &p_high)
 Set the peak, covariance matrix, lower limits, and upper limits. More...
 
virtual double pdf (const vec_t &x) const
 Compute the probability density function (arbitrary normalization)
 
virtual double log_pdf (const vec_t &x) const
 Compute the natural log of the probability density function (arbitrary normalization)
 
virtual void operator() (vec_t &x) const
 Sample the distribution.
 
- Public Member Functions inherited from o2scl::prob_dens_mdim_gaussian< vec_t, mat_t >
const mat_t & get_chol ()
 Get the Cholesky decomposition.
 
const mat_t & get_covar_inv ()
 Get the inverse of the covariance matrix.
 
const vec_t & get_peak ()
 Get the peak location.
 
const double & get_norm ()
 Get the normalization.
 
virtual size_t dim () const
 The dimensionality.
 
 prob_dens_mdim_gaussian ()
 Create an empty distribution.
 
 prob_dens_mdim_gaussian (const prob_dens_mdim_gaussian &pdmg)
 Copy constructor.
 
prob_dens_mdim_gaussianoperator= (const prob_dens_mdim_gaussian &pdmg)
 Copy constructor with operator=.
 
template<class mat2_t , class vec2_t , class mat2_col_t = const_matrix_column_gen<mat2_t>>
 prob_dens_mdim_gaussian (size_t p_mdim, size_t n_pts, const mat2_t &pts, const vec2_t &vals)
 Create a distribution from a set of samples from a multidimensional Gaussian.
 
 prob_dens_mdim_gaussian (size_t p_ndim, vec_t &p_peak, mat_t &covar)
 Create a distribution from the covariance matrix.
 
void set (size_t p_ndim, vec_t &p_peak, mat_t &covar)
 Set the peak and covariance matrix for the distribution. More...
 
void set_alt (size_t p_ndim, vec_t &p_peak, mat_t &p_chol, mat_t &p_covar_inv, double p_norm)
 Alternate set function for use when covariance matrix has already been decomposed and inverted.
 
template<class vec_vec_t , class mat_col_t , class func_t >
void set_gproc (size_t n_dim, size_t n_init, vec_vec_t &x, vec_t &y, func_t &fcovar)
 Given a data set and a covariance function, construct probability distribution based on a Gaussian process which includes noise. More...
 

Public Attributes

size_t samp_max
 Maximum number of samples.
 
- Public Attributes inherited from o2scl::prob_dens_mdim_gaussian< vec_t, mat_t >
o2scl::prob_dens_gaussian pdg
 Standard normal.
 

Protected Attributes

vec_t low
 Lower limits.
 
vec_t high
 Upper limits.
 
- Protected Attributes inherited from o2scl::prob_dens_mdim_gaussian< vec_t, mat_t >
mat_t chol
 Cholesky decomposition.
 
mat_t covar_inv
 Inverse of the covariance matrix.
 
vec_t peak
 Location of the peak.
 
double norm
 Normalization factor, $ \det ( 2 \pi \Sigma)^{-1/2} $.
 
size_t ndim
 Number of dimensions.
 
vec_t q
 Temporary storage 1.
 
vec_t vtmp
 Temporary storage 2.
 

Member Function Documentation

◆ set()

template<class vec_t = boost::numeric::ublas::vector<double>, class mat_t = boost::numeric::ublas::matrix<double>>
void o2scl::prob_dens_mdim_bound_gaussian< vec_t, mat_t >::set ( size_t  p_ndim,
vec_t &  p_peak,
mat_t &  covar,
vec_t &  p_low,
vec_t &  p_high 
)
inline
Note
This function is called in constructors and thus should not be virtual.

Definition at line 1354 of file prob_dens_func.h.


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

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).