IMP logo
Public Member Functions | Static Public Member Functions | Protected Member Functions
IMP::isd::GaussianProcessInterpolation Class Reference

Detailed Description

GaussianProcessInterpolation.

This class provides methods to perform bayesian interpolation/smoothing of data using a gaussian process prior on the function to be interpolated. It takes a dataset as input (via its sufficient statistics) along with prior mean and covariance functions. It outputs the value of the posterior mean and covariance functions at points requested by the user.

+ Inheritance diagram for IMP::isd::GaussianProcessInterpolation:

List of all members.

Public Member Functions

 GaussianProcessInterpolation (FloatsList x, Floats sample_mean, Floats sample_std, unsigned n_obs, UnivariateFunction *mean_function, BivariateFunction *covariance_function, Particle *sigma, double sparse_cutoff=1e-7)
void force_covariance_update ()
void force_mean_update ()
FloatsList get_data_abscissa () const
Floats get_data_mean () const
FloatsList get_data_variance () const
ContainersTemp get_input_containers () const
ParticlesTemp get_input_particles () const
bool get_m_particle_is_optimized (unsigned i) const
unsigned get_number_of_m_particles () const
unsigned get_number_of_Omega_particles () const
bool get_Omega_particle_is_optimized (unsigned i) const
double get_posterior_covariance (Floats x1, Floats x2) const
VectorXd get_posterior_covariance_derivative (Floats x) const
Floats get_posterior_covariance_derivative (Floats x, bool) const
MatrixXd get_posterior_covariance_hessian (Floats x) const
FloatsList get_posterior_covariance_hessian (Floats x, bool) const
MatrixXd get_posterior_covariance_matrix (FloatsList x) const
FloatsList get_posterior_covariance_matrix (FloatsList x, bool) const
double get_posterior_mean (Floats x) const
- Public Member Functions inherited from IMP::base::Object
std::size_t __hash__ () const
virtual std::string get_type_name () const =0
 Return a string identifying the type of the object.
virtual IMP::base::VersionInfo get_version_info () const =0
 Get information about the module and version of the object.
void set_check_level (CheckLevel l)
void set_log_level (LogLevel l)
 Set the logging level used in this object.
void set_was_used (bool tf) const
void show (std::ostream &out=std::cout) const
const std::string & get_name () const
void set_name (std::string name)

Static Public Member Functions

static
GaussianProcessInterpolation
get_from (IMP::base::Object *o)

Protected Member Functions

void add_to_m_particle_derivative (unsigned particle, double value, DerivativeAccumulator &accum)
void add_to_Omega_particle_derivative (unsigned particle, double value, DerivativeAccumulator &accum)
VectorXd get_I () const
VectorXd get_m () const
VectorXd get_m_derivative (unsigned particle) const
VectorXd get_m_second_derivative (unsigned particle1, unsigned particle2) const
MatrixXd get_Omega () const
MatrixXd get_Omega_derivative (unsigned particle) const
MatrixXd get_Omega_second_derivative (unsigned particle1, unsigned particle2) const
MatrixXd get_Omi () const
VectorXd get_OmiIm () const
Eigen::DiagonalMatrix< double,
Eigen::Dynamic > 
get_S () const
MatrixXd get_W () const
VectorXd get_wx_vector (Floats xval) const
- Protected Member Functions inherited from IMP::base::Object
 Object (std::string name)

Constructor & Destructor Documentation

IMP::isd::GaussianProcessInterpolation::GaussianProcessInterpolation ( FloatsList  x,
Floats  sample_mean,
Floats  sample_std,
unsigned  n_obs,
UnivariateFunction mean_function,
BivariateFunction covariance_function,
Particle sigma,
double  sparse_cutoff = 1e-7 
)

Constructor for the gaussian process

Parameters:
in)x : a list of coordinates in N-dimensional space corresponding to the abscissa of each observation
in)sample_mean $I$ : vector of mean observations at each of the previously defined coordinates
in)sample_std $s$ : vector of sample standard deviations
in)mean_function $m$ : a pointer to the prior mean function to use. Should be compatible with the size of x(i).
in)covariance_function $w$: prior covariance function.
in)sigma : ISD Scale (proportionality factor to S)
in)sparse_cutoff : when to consider that a matrix entry is zero

Computes the necessary matrices and inverses when called.


Member Function Documentation

double IMP::isd::GaussianProcessInterpolation::get_posterior_mean ( Floats  x) const

Get posterior mean and covariance functions, at the points requested Posterior mean is defined as

\[\hat{I}(x) = m(x) + {}^t\mathbf{w}(q) (\mathbf{W}+\mathbf{S})^{-1} (\mathbf{I}-\mathbf{m}) \]

Posterior covariance

\[\hat{\sigma}^2(x,x') = w(x,x') - {}^t\mathbf{w}(x) (\mathbf{W} + \mathbf{S})^{-1} \mathbf{w}(x') \]

where $\mathbf{m}$ is the vector built by evaluating the prior mean function at the observation points; $\mathbf{w}(x)$ is the vector of covariances between each observation point and the current point; $\mathbf{W}$ is the prior covariance matrix built by evaluating the covariance function at each of the observations; $\mathbf{S}$ is the diagonal covariance matrix built from sample_std and n_obs.

Both functions will check if the mean or covariance functions have changed since the last call, and will recompute $(\mathbf{W} + \mathbf{S})^{-1}$ if necessary.


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

Generated on Tue May 22 2012 23:33:35 for IMP by doxygen 1.8.1