8 #ifndef IMPISD_GAUSSIAN_PROCESS_INTERPOLATION_H
9 #define IMPISD_GAUSSIAN_PROCESS_INTERPOLATION_H
11 #include <IMP/isd/isd_config.h>
13 #include <boost/scoped_ptr.hpp>
17 #include <IMP/algebra/eigen3/Eigen/Dense>
18 #include <IMP/algebra/eigen3/Eigen/Cholesky>
20 IMPISD_BEGIN_NAMESPACE
22 class GaussianProcessInterpolationRestraint;
23 class GaussianProcessInterpolationScoreState;
51 Floats sample_std,
unsigned n_obs,
55 double sparse_cutoff = 1e-7);
80 double get_posterior_mean(
Floats x)
const;
81 double get_posterior_covariance(
Floats x1,
Floats x2)
const;
84 IMP_Eigen::MatrixXd get_posterior_covariance_matrix(
FloatsList x)
const;
96 IMP_Eigen::VectorXd get_posterior_covariance_derivative(
Floats x)
const;
99 Floats get_posterior_covariance_derivative(
Floats x,
bool)
const;
108 IMP_Eigen::MatrixXd get_posterior_covariance_hessian(
Floats x)
const;
118 void force_mean_update();
119 void force_covariance_update();
123 unsigned get_number_of_m_particles()
const;
127 bool get_m_particle_is_optimized(
unsigned i)
const;
131 unsigned get_number_of_Omega_particles()
const;
135 bool get_Omega_particle_is_optimized(
unsigned i)
const;
139 Floats get_data_mean()
const;
143 friend class GaussianProcessInterpolationScoreState;
149 IMP_Eigen::VectorXd get_I()
const {
return I_; }
151 IMP_Eigen::VectorXd get_m()
const;
153 IMP_Eigen::VectorXd get_m_derivative(
unsigned particle)
const;
155 IMP_Eigen::VectorXd get_m_second_derivative(
unsigned particle1,
156 unsigned particle2)
const;
158 void add_to_m_particle_derivative(
unsigned particle,
double value,
161 IMP_Eigen::VectorXd get_wx_vector(
Floats xval)
const;
163 IMP_Eigen::DiagonalMatrix<double, IMP_Eigen::Dynamic> get_S()
const {
167 IMP_Eigen::MatrixXd get_W()
const;
169 IMP_Eigen::MatrixXd get_Omega()
const;
171 IMP_Eigen::MatrixXd get_Omega_derivative(
unsigned particle)
const;
173 IMP_Eigen::MatrixXd get_Omega_second_derivative(
unsigned particle1,
174 unsigned particle2)
const;
176 void add_to_Omega_particle_derivative(
unsigned particle,
double value,
179 IMP_Eigen::LDLT<IMP_Eigen::MatrixXd, IMP_Eigen::Upper> get_ldlt()
const;
181 IMP_Eigen::MatrixXd get_Omi()
const;
183 IMP_Eigen::VectorXd get_OmiIm()
const;
189 void update_flags_mean();
190 void update_flags_covariance();
195 void compute_Omega();
201 void compute_OmiIm();
204 IMP_Eigen::VectorXd get_wx_vector_derivative(
Floats q,
unsigned i)
const;
206 IMP_Eigen::VectorXd get_wx_vector_second_derivative(
Floats q,
unsigned i,
210 IMP_Eigen::VectorXd get_dcov_dwq(
Floats q)
const;
212 IMP_Eigen::MatrixXd get_dcov_dOm(
Floats q)
const;
214 IMP_Eigen::MatrixXd get_d2cov_dwq_dwq()
const;
216 IMP_Eigen::MatrixXd get_d2cov_dwq_dOm(
Floats q,
unsigned m)
const;
218 IMP_Eigen::MatrixXd get_d2cov_dOm_dOm(
Floats q,
unsigned m,
unsigned n)
const;
221 void compute_I(
Floats mean);
223 void compute_S(
Floats std);
236 IMP_Eigen::VectorXd I_, m_;
237 IMP_Eigen::MatrixXd W_, Omega_, Omi_;
238 IMP_Eigen::DiagonalMatrix<double, IMP_Eigen::Dynamic> S_;
239 IMP_Eigen::VectorXd OmiIm_;
240 bool flag_m_, flag_m_gpir_, flag_Omi_, flag_OmiIm_, flag_W_, flag_Omega_,
241 flag_Omega_gpir_, flag_ldlt_;
245 IMP_Eigen::LDLT<IMP_Eigen::MatrixXd, IMP_Eigen::Upper> ldlt_;
Base class for functions of one variable.
#define IMP_OBJECT_METHODS(Name)
Define the basic things needed by any Object.
A decorator for scale parameters particles.
Classes for general functions.
Various general useful macros for IMP.
Common base class for heavy weight IMP objects.
gaussian process restraint
A smart pointer to a ref-counted Object that is a class member.
Base class for functions of two variables.
GaussianProcessInterpolation.
Class to handle individual particles of a Model object.
Class for adding derivatives from restraints to the model.
Classes for general functions.