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;
119 void force_mean_update();
120 void force_covariance_update();
124 unsigned get_number_of_m_particles()
const;
128 bool get_m_particle_is_optimized(
unsigned i)
const;
132 unsigned get_number_of_Omega_particles()
const;
136 bool get_Omega_particle_is_optimized(
unsigned i)
const;
140 Floats get_data_mean()
const;
144 friend class GaussianProcessInterpolationScoreState;
150 IMP_Eigen::VectorXd get_I()
const {
return I_; }
152 IMP_Eigen::VectorXd get_m()
const;
154 IMP_Eigen::VectorXd get_m_derivative(
unsigned particle)
const;
156 IMP_Eigen::VectorXd get_m_second_derivative(
unsigned particle1,
157 unsigned particle2)
const;
159 void add_to_m_particle_derivative(
unsigned particle,
double value,
162 IMP_Eigen::VectorXd get_wx_vector(
Floats xval)
const;
164 IMP_Eigen::DiagonalMatrix<double, IMP_Eigen::Dynamic> get_S()
const {
168 IMP_Eigen::MatrixXd get_W()
const;
170 IMP_Eigen::MatrixXd get_Omega()
const;
172 IMP_Eigen::MatrixXd get_Omega_derivative(
unsigned particle)
const;
174 IMP_Eigen::MatrixXd get_Omega_second_derivative(
unsigned particle1,
175 unsigned particle2)
const;
177 void add_to_Omega_particle_derivative(
unsigned particle,
double value,
180 IMP_Eigen::LDLT<IMP_Eigen::MatrixXd, IMP_Eigen::Upper> get_ldlt()
const;
182 IMP_Eigen::MatrixXd get_Omi()
const;
184 IMP_Eigen::VectorXd get_OmiIm()
const;
190 void update_flags_mean();
191 void update_flags_covariance();
196 void compute_Omega();
202 void compute_OmiIm();
205 IMP_Eigen::VectorXd get_wx_vector_derivative(
Floats q,
unsigned i)
const;
207 IMP_Eigen::VectorXd get_wx_vector_second_derivative(
Floats q,
unsigned i,
211 IMP_Eigen::VectorXd get_dcov_dwq(
Floats q)
const;
213 IMP_Eigen::MatrixXd get_dcov_dOm(
Floats q)
const;
215 IMP_Eigen::MatrixXd get_d2cov_dwq_dwq()
const;
217 IMP_Eigen::MatrixXd get_d2cov_dwq_dOm(
Floats q,
unsigned m)
const;
219 IMP_Eigen::MatrixXd get_d2cov_dOm_dOm(
Floats q,
unsigned m,
unsigned n)
const;
222 void compute_I(
Floats mean);
224 void compute_S(
Floats std);
237 IMP_Eigen::VectorXd I_, m_;
238 IMP_Eigen::MatrixXd W_, Omega_, Omi_;
239 IMP_Eigen::DiagonalMatrix<double, IMP_Eigen::Dynamic> S_;
240 IMP_Eigen::VectorXd OmiIm_;
241 bool flag_m_, flag_m_gpir_, flag_Omi_, flag_OmiIm_, flag_W_, flag_Omega_,
242 flag_Omega_gpir_, flag_ldlt_;
246 IMP_Eigen::LDLT<IMP_Eigen::MatrixXd, IMP_Eigen::Upper> ldlt_;
Base class for functions of one variable.
Class for adding derivatives from restraints to the model.
A smart pointer to a ref-counted Object that is a class memeber.
A decorator for scale parameters particles.
ParticlesTemp get_input_particles(const ModelObjectsTemp &mos)
Classes for general functions.
Import IMP/kernel/macros.h in the namespace.
ContainersTemp get_input_containers(const ModelObjectsTemp &mos)
gaussian process restraint
IMP::base::Vector< IMP::base::WeakPointer< Container > > ContainersTemp
Base class for functions of two variables.
#define IMP_OBJECT_METHODS(Name)
Define the basic things needed by any Object.
Class to handle individual model particles.
Common base class for heavy weight IMP objects.
GaussianProcessInterpolation.
Classes for general functions.