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 <Eigen/Dense>
18 #include <Eigen/Cholesky>
20 IMPISD_BEGIN_NAMESPACE
22 using Eigen::MatrixXd;
23 using Eigen::VectorXd;
24 using Eigen::RowVectorXd;
27 class GaussianProcessInterpolationRestraint;
28 class GaussianProcessInterpolationScoreState;
63 double sparse_cutoff=1e-7);
87 double get_posterior_mean(
Floats x)
const;
88 double get_posterior_covariance(
Floats x1,
92 MatrixXd get_posterior_covariance_matrix(
FloatsList x)
const;
104 VectorXd get_posterior_covariance_derivative(
Floats x)
const;
107 Floats get_posterior_covariance_derivative(
Floats x,
bool)
const;
116 MatrixXd get_posterior_covariance_hessian(
Floats x)
const;
127 void force_mean_update();
128 void force_covariance_update();
132 unsigned get_number_of_m_particles()
const;
136 bool get_m_particle_is_optimized(
unsigned i)
const;
140 unsigned get_number_of_Omega_particles()
const;
144 bool get_Omega_particle_is_optimized(
unsigned i)
const;
148 Floats get_data_mean()
const;
152 friend class GaussianProcessInterpolationScoreState;
158 VectorXd get_I()
const {
return I_;}
160 VectorXd get_m()
const;
162 VectorXd get_m_derivative(
unsigned particle)
const;
164 VectorXd get_m_second_derivative(
unsigned particle1,
unsigned particle2)
167 void add_to_m_particle_derivative(
unsigned particle,
double value,
170 VectorXd get_wx_vector(
Floats xval)
const;
172 Eigen::DiagonalMatrix<double, Eigen::Dynamic> get_S()
const
177 MatrixXd get_W()
const;
179 MatrixXd get_Omega()
const;
181 MatrixXd get_Omega_derivative(
unsigned particle)
const;
183 MatrixXd get_Omega_second_derivative(
unsigned particle1,
unsigned particle2)
186 void add_to_Omega_particle_derivative(
unsigned particle,
double value,
189 MatrixXd get_Omi()
const;
191 VectorXd get_OmiIm()
const;
198 void update_flags_mean();
199 void update_flags_covariance();
204 void compute_Omega();
208 void compute_OmiIm();
211 VectorXd get_wx_vector_derivative(
Floats q,
unsigned i)
const;
213 VectorXd get_wx_vector_second_derivative(
Floats q,
unsigned i,
unsigned j)
217 VectorXd get_dcov_dwq(
Floats q)
const;
219 MatrixXd get_dcov_dOm(
Floats q)
const;
221 MatrixXd get_d2cov_dwq_dwq()
const;
223 MatrixXd get_d2cov_dwq_dOm(
Floats q,
unsigned m)
const;
225 MatrixXd get_d2cov_dOm_dOm(
Floats q,
unsigned m,
unsigned n)
const;
228 void compute_I(
Floats mean);
230 void compute_S(
Floats std);
240 IMP::internal::OwnerPointer<UnivariateFunction> mean_function_;
242 IMP::internal::OwnerPointer<BivariateFunction> covariance_function_;
244 MatrixXd W_,Omega_,Omi_;
245 Eigen::DiagonalMatrix<double, Eigen::Dynamic> S_;
247 bool flag_m_, flag_m_gpir_, flag_Omi_, flag_OmiIm_, flag_W_,
248 flag_Omega_, flag_Omega_gpir_;
249 Pointer<Particle> sigma_;