1use crate::errors::{GpError, Result};
2use crate::mean_models::*;
3use crate::optimization::{CobylaParams, optimize_params, prepare_multistart};
4use crate::parameters::{GpParams, GpValidParams};
5use crate::utils::{DistanceMatrix, NormalizedData, pairwise_differences};
6use crate::{ThetaTuning, correlation_models::*};
7
8use linfa::dataset::{WithLapack, WithoutLapack};
9use linfa::prelude::{Dataset, DatasetBase, Fit, Float, PredictInplace};
10
11#[cfg(not(feature = "blas"))]
12use linfa_linalg::{cholesky::*, eigh::*, qr::*, svd::*, triangular::*};
13#[cfg(feature = "blas")]
14use log::warn;
15#[cfg(feature = "blas")]
16use ndarray_linalg::{cholesky::*, eigh::*, qr::*, svd::*, triangular::*};
17
18use linfa_pls::PlsRegression;
19use ndarray::{Array, Array1, Array2, ArrayBase, Axis, Data, Ix1, Ix2, Zip};
20
21use ndarray_rand::RandomExt;
22use ndarray_rand::rand_distr::Normal;
23use ndarray_stats::QuantileExt;
24
25use log::debug;
26use rayon::prelude::*;
27#[cfg(feature = "serializable")]
28use serde::{Deserialize, Serialize};
29use std::fmt;
30use std::time::Instant;
31
32pub const GP_OPTIM_N_START: usize = 10;
34pub const GP_COBYLA_MIN_EVAL: usize = 25;
36pub const GP_COBYLA_MAX_EVAL: usize = 1000;
38
39#[derive(Default, Debug)]
42#[cfg_attr(
43 feature = "serializable",
44 derive(Serialize, Deserialize),
45 serde(bound(deserialize = "F: Deserialize<'de>"))
46)]
47pub(crate) struct GpInnerParams<F: Float> {
48 sigma2: F,
50 beta: Array2<F>,
52 gamma: Array2<F>,
54 r_chol: Array2<F>,
56 ft: Array2<F>,
58 ft_qr_r: Array2<F>,
60}
61
62impl<F: Float> Clone for GpInnerParams<F> {
63 fn clone(&self) -> Self {
64 Self {
65 sigma2: self.sigma2.to_owned(),
66 beta: self.beta.to_owned(),
67 gamma: self.gamma.to_owned(),
68 r_chol: self.r_chol.to_owned(),
69 ft: self.ft.to_owned(),
70 ft_qr_r: self.ft_qr_r.to_owned(),
71 }
72 }
73}
74
75#[derive(Debug)]
166#[cfg_attr(
167 feature = "serializable",
168 derive(Serialize, Deserialize),
169 serde(bound(
170 serialize = "F: Serialize, Mean: Serialize, Corr: Serialize",
171 deserialize = "F: Deserialize<'de>, Mean: Deserialize<'de>, Corr: Deserialize<'de>"
172 ))
173)]
174pub struct GaussianProcess<F: Float, Mean: RegressionModel<F>, Corr: CorrelationModel<F>> {
175 theta: Array1<F>,
177 likelihood: F,
180 inner_params: GpInnerParams<F>,
182 w_star: Array2<F>,
184 xt_norm: NormalizedData<F>,
186 yt_norm: NormalizedData<F>,
188 pub(crate) training_data: (Array2<F>, Array1<F>),
190 pub(crate) params: GpValidParams<F, Mean, Corr>,
192}
193
194pub(crate) enum GpSamplingMethod {
195 Cholesky,
196 EigenValues,
197}
198
199pub type Kriging<F> = GpParams<F, ConstantMean, SquaredExponentialCorr>;
201
202impl<F: Float> Kriging<F> {
203 pub fn params() -> GpParams<F, ConstantMean, SquaredExponentialCorr> {
205 GpParams::new(ConstantMean(), SquaredExponentialCorr())
206 }
207}
208
209impl<F: Float, Mean: RegressionModel<F>, Corr: CorrelationModel<F>> Clone
210 for GaussianProcess<F, Mean, Corr>
211{
212 fn clone(&self) -> Self {
213 Self {
214 theta: self.theta.to_owned(),
215 likelihood: self.likelihood,
216 inner_params: self.inner_params.clone(),
217 w_star: self.w_star.to_owned(),
218 xt_norm: self.xt_norm.clone(),
219 yt_norm: self.yt_norm.clone(),
220 training_data: self.training_data.clone(),
221 params: self.params.clone(),
222 }
223 }
224}
225
226impl<F: Float, Mean: RegressionModel<F>, Corr: CorrelationModel<F>> fmt::Display
227 for GaussianProcess<F, Mean, Corr>
228{
229 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
230 write!(
231 f,
232 "GP(mean={}, corr={}, theta={}, variance={}, likelihood={})",
233 self.params.mean,
234 self.params.corr,
235 self.theta,
236 self.inner_params.sigma2,
237 self.likelihood,
238 )
239 }
240}
241
242impl<F: Float, Mean: RegressionModel<F>, Corr: CorrelationModel<F>> GaussianProcess<F, Mean, Corr> {
243 pub fn params<NewMean: RegressionModel<F>, NewCorr: CorrelationModel<F>>(
245 mean: NewMean,
246 corr: NewCorr,
247 ) -> GpParams<F, NewMean, NewCorr> {
248 GpParams::new(mean, corr)
249 }
250
251 pub fn predict(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>) -> Result<Array1<F>> {
254 let xnorm = (x - &self.xt_norm.mean) / &self.xt_norm.std;
255 let f = self.params.mean.value(&xnorm);
257 let corr = self._compute_correlation(&xnorm);
259 let y_ = &f.dot(&self.inner_params.beta) + &corr.dot(&self.inner_params.gamma);
261 Ok((&y_ * &self.yt_norm.std + &self.yt_norm.mean).remove_axis(Axis(1)))
263 }
264
265 pub fn predict_var(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>) -> Result<Array1<F>> {
268 let (rt, u, _) = self._compute_rt_u(x);
269
270 let mut mse = Array::ones(rt.ncols()) - rt.mapv(|v| v * v).sum_axis(Axis(0))
271 + u.mapv(|v: F| v * v).sum_axis(Axis(0));
272 mse.mapv_inplace(|v| self.inner_params.sigma2 * v);
273
274 Ok(mse.mapv(|v| if v < F::zero() { F::zero() } else { F::cast(v) }))
277 }
278
279 fn _compute_covariance(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>) -> Array2<F> {
281 let (rt, u, xnorm) = self._compute_rt_u(x);
282
283 let cross_dx = pairwise_differences(&xnorm, &xnorm);
284 let k = self.params.corr.value(&cross_dx, &self.theta, &self.w_star);
285 let k = k
286 .into_shape_with_order((xnorm.nrows(), xnorm.nrows()))
287 .unwrap();
288
289 let mut cov_matrix = k - rt.t().to_owned().dot(&rt) + u.t().dot(&u);
292 cov_matrix.mapv_inplace(|v| self.inner_params.sigma2 * v);
293 cov_matrix
294 }
295
296 fn _compute_rt_u(
299 &self,
300 x: &ArrayBase<impl Data<Elem = F>, Ix2>,
301 ) -> (Array2<F>, Array2<F>, Array2<F>) {
302 let xnorm = (x - &self.xt_norm.mean) / &self.xt_norm.std;
303 let corr = self._compute_correlation(&xnorm);
304 let inners = &self.inner_params;
305
306 let corr_t = corr.t().to_owned();
307 #[cfg(feature = "blas")]
308 let rt = inners
309 .r_chol
310 .to_owned()
311 .with_lapack()
312 .solve_triangular(UPLO::Lower, Diag::NonUnit, &corr_t.with_lapack())
313 .unwrap()
314 .without_lapack();
315 #[cfg(not(feature = "blas"))]
316 let rt = inners
317 .r_chol
318 .solve_triangular(&corr_t, UPLO::Lower)
319 .unwrap();
320
321 let rhs = inners.ft.t().dot(&rt) - self.params.mean.value(&xnorm).t();
322 #[cfg(feature = "blas")]
323 let u = inners
324 .ft_qr_r
325 .to_owned()
326 .t()
327 .with_lapack()
328 .solve_triangular(UPLO::Upper, Diag::NonUnit, &rhs.with_lapack())
329 .unwrap()
330 .without_lapack();
331 #[cfg(not(feature = "blas"))]
332 let u = inners
333 .ft_qr_r
334 .t()
335 .solve_triangular(&rhs, UPLO::Lower)
336 .unwrap();
337 (rt, u, xnorm)
338 }
339
340 fn _compute_correlation(&self, xnorm: &ArrayBase<impl Data<Elem = F>, Ix2>) -> Array2<F> {
342 let dx = pairwise_differences(xnorm, &self.xt_norm.data);
344 let r = self.params.corr.value(&dx, &self.theta, &self.w_star);
346 let n_obs = xnorm.nrows();
347 let nt = self.xt_norm.data.nrows();
348 r.into_shape_with_order((n_obs, nt)).unwrap().to_owned()
349 }
350
351 pub fn sample_chol(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>, n_traj: usize) -> Array2<F> {
353 self._sample(x, n_traj, GpSamplingMethod::Cholesky)
354 }
355
356 pub fn sample_eig(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>, n_traj: usize) -> Array2<F> {
358 self._sample(x, n_traj, GpSamplingMethod::EigenValues)
359 }
360
361 pub fn sample(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>, n_traj: usize) -> Array2<F> {
363 self.sample_eig(x, n_traj)
364 }
365
366 fn _sample(
371 &self,
372 x: &ArrayBase<impl Data<Elem = F>, Ix2>,
373 n_traj: usize,
374 method: GpSamplingMethod,
375 ) -> Array2<F> {
376 let mean = self.predict(x).unwrap();
377 let cov = self._compute_covariance(x);
378 sample(x, mean.insert_axis(Axis(1)), cov, n_traj, method)
379 }
380
381 pub fn theta(&self) -> &Array1<F> {
383 &self.theta
384 }
385
386 pub fn variance(&self) -> F {
388 self.inner_params.sigma2
389 }
390
391 pub fn likelihood(&self) -> F {
393 self.likelihood
394 }
395
396 pub fn kpls_dim(&self) -> Option<usize> {
398 if self.w_star.ncols() < self.xt_norm.ncols() {
399 Some(self.w_star.ncols())
400 } else {
401 None
402 }
403 }
404
405 pub fn dims(&self) -> (usize, usize) {
407 (self.xt_norm.ncols(), self.yt_norm.ncols())
408 }
409
410 pub fn predict_kth_derivatives(
413 &self,
414 x: &ArrayBase<impl Data<Elem = F>, Ix2>,
415 kx: usize,
416 ) -> Array1<F> {
417 let xnorm = (x - &self.xt_norm.mean) / &self.xt_norm.std;
418 let corr = self._compute_correlation(&xnorm);
419
420 let beta = &self.inner_params.beta;
421 let gamma = &self.inner_params.gamma;
422
423 let df_dx_kx = if self.inner_params.beta.nrows() <= 1 + self.xt_norm.data.ncols() {
424 let df = self.params.mean.jacobian(&x.row(0));
426 let df_dx = df.t().row(kx).dot(beta);
427 df_dx.broadcast((x.nrows(), 1)).unwrap().to_owned()
428 } else {
429 let mut dfdx = Array2::zeros((x.nrows(), 1));
431 Zip::from(dfdx.rows_mut())
432 .and(xnorm.rows())
433 .for_each(|mut dfxi, xi| {
434 let df = self.params.mean.jacobian(&xi);
435 let df_dx = (df.t().row(kx)).dot(beta);
436 dfxi.assign(&df_dx);
437 });
438 dfdx
439 };
440
441 let nr = x.nrows();
442 let nc = self.xt_norm.data.nrows();
443 let d_dx_1 = &xnorm
444 .column(kx)
445 .to_owned()
446 .into_shape_with_order((nr, 1))
447 .unwrap()
448 .broadcast((nr, nc))
449 .unwrap()
450 .to_owned();
451
452 let d_dx_2 = self
453 .xt_norm
454 .data
455 .column(kx)
456 .to_owned()
457 .as_standard_layout()
458 .into_shape_with_order((1, nc))
459 .unwrap()
460 .to_owned();
461
462 let d_dx = d_dx_1 - d_dx_2;
463
464 let theta = &self.theta.to_owned();
466 let d_dx_corr = d_dx * corr;
467
468 let res = (df_dx_kx - d_dx_corr.dot(gamma).map(|v| F::cast(2.) * theta[kx] * *v))
472 * self.yt_norm.std[0]
473 / self.xt_norm.std[kx];
474 res.column(0).to_owned()
475 }
476
477 pub fn predict_gradients(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>) -> Array2<F> {
480 let mut drv = Array2::<F>::zeros((x.nrows(), self.xt_norm.data.ncols()));
481 Zip::from(drv.rows_mut())
482 .and(x.rows())
483 .for_each(|mut row, xi| {
484 let pred = self.predict_jacobian(&xi);
485 row.assign(&pred.column(0));
486 });
487 drv
488 }
489
490 pub fn predict_jacobian(&self, x: &ArrayBase<impl Data<Elem = F>, Ix1>) -> Array2<F> {
493 let xx = x.to_owned().insert_axis(Axis(0));
494 let mut jac = Array2::zeros((xx.ncols(), 1));
495
496 let xnorm = (xx - &self.xt_norm.mean) / &self.xt_norm.std;
497
498 let beta = &self.inner_params.beta;
499 let gamma = &self.inner_params.gamma;
500
501 let df = self.params.mean.jacobian(&xnorm.row(0));
502 let df_dx = df.t().dot(beta);
503
504 let dr =
505 self.params
506 .corr
507 .jacobian(&xnorm.row(0), &self.xt_norm.data, &self.theta, &self.w_star);
508
509 let dr_dx = df_dx + dr.t().dot(gamma);
510 Zip::from(jac.rows_mut())
511 .and(dr_dx.rows())
512 .and(&self.xt_norm.std)
513 .for_each(|mut jc, dr_i, std_i| {
514 let jc_i = dr_i.map(|v| *v * self.yt_norm.std[0] / *std_i);
515 jc.assign(&jc_i)
516 });
517
518 jac
519 }
520
521 #[cfg(not(feature = "blas"))]
524 pub fn predict_var_gradients_single(
525 &self,
526 x: &ArrayBase<impl Data<Elem = F>, Ix1>,
527 ) -> Array1<F> {
528 let x = &(x.to_owned().insert_axis(Axis(0)));
529 let xnorm = (x - &self.xt_norm.mean) / &self.xt_norm.std;
530 let dx = pairwise_differences(&xnorm, &self.xt_norm.data);
531 let sigma2 = self.inner_params.sigma2;
532 let r_chol = &self.inner_params.r_chol;
533
534 let r = self.params.corr.value(&dx, &self.theta, &self.w_star);
535 let dr =
536 self.params
537 .corr
538 .jacobian(&xnorm.row(0), &self.xt_norm.data, &self.theta, &self.w_star);
539
540 let rho1 = r_chol.solve_triangular(&r, UPLO::Lower).unwrap();
542
543 let inv_kr = r_chol.t().solve_triangular(&rho1, UPLO::Upper).unwrap();
545
546 let p2 = inv_kr.t().dot(&dr);
551
552 let f_x = self.params.mean.value(&xnorm).t().to_owned();
553 let f_mean = self.params.mean.value(&self.xt_norm.data);
554
555 let rho2 = r_chol.solve_triangular(&f_mean, UPLO::Lower).unwrap();
557 let inv_kf = r_chol.t().solve_triangular(&rho2, UPLO::Upper).unwrap();
559
560 let a_mat = f_x.t().to_owned() - r.t().dot(&inv_kf);
562
563 let b_mat = f_mean.t().dot(&inv_kf);
565 let rho3 = b_mat.cholesky().unwrap();
567 let inv_bat = rho3.solve_triangular(&a_mat.t(), UPLO::Lower).unwrap();
569 let d_mat = rho3.t().solve_triangular(&inv_bat, UPLO::Upper).unwrap();
571
572 let df = self.params.mean.jacobian(&xnorm.row(0));
573
574 let d_a = df.t().to_owned() - dr.t().dot(&inv_kf);
576
577 let p4 = d_mat.t().dot(&d_a.t());
582 let two = F::cast(2.);
583 let prime = (p4 - p2).mapv(|v| two * v);
584
585 let x_std = &self.xt_norm.std;
586 let dvar = (prime / x_std).mapv(|v| v * sigma2);
587 dvar.row(0).into_owned()
588 }
589
590 #[cfg(feature = "blas")]
592 pub fn predict_var_gradients_single(
593 &self,
594 x: &ArrayBase<impl Data<Elem = F>, Ix1>,
595 ) -> Array1<F> {
596 let x = &(x.to_owned().insert_axis(Axis(0)));
597 let xnorm = (x - &self.xt_norm.mean) / &self.xt_norm.std;
598
599 let dx = pairwise_differences(&xnorm, &self.xt_norm.data);
600
601 let sigma2 = self.inner_params.sigma2;
602 let r_chol = &self.inner_params.r_chol.to_owned().with_lapack();
603
604 let r = self
605 .params
606 .corr
607 .value(&dx, &self.theta, &self.w_star)
608 .with_lapack();
609 let dr = self
610 .params
611 .corr
612 .jacobian(&xnorm.row(0), &self.xt_norm.data, &self.theta, &self.w_star)
613 .with_lapack();
614
615 let rho1 = r_chol
616 .solve_triangular(UPLO::Lower, Diag::NonUnit, &r)
617 .unwrap();
618 let inv_kr = r_chol
619 .t()
620 .solve_triangular(UPLO::Upper, Diag::NonUnit, &rho1)
621 .unwrap();
622
623 let p2 = inv_kr.t().dot(&dr);
626
627 let f_x = self.params.mean.value(x).t().to_owned();
628 let f_mean = self.params.mean.value(&self.xt_norm.data).with_lapack();
629
630 let rho2 = r_chol
631 .solve_triangular(UPLO::Lower, Diag::NonUnit, &f_mean)
632 .unwrap();
633 let inv_kf = r_chol
634 .t()
635 .solve_triangular(UPLO::Upper, Diag::NonUnit, &rho2)
636 .unwrap();
637
638 let a_mat = f_x.t().to_owned().with_lapack() - r.t().dot(&inv_kf);
639
640 let b_mat = f_mean.t().dot(&inv_kf);
641
642 let d_mat = match b_mat.cholesky(UPLO::Lower) {
643 Ok(rho3) => {
644 let inv_bat = rho3
645 .solve_triangular(UPLO::Upper, Diag::NonUnit, &a_mat.t().to_owned())
646 .unwrap();
647 rho3.t()
648 .solve_triangular(UPLO::Upper, Diag::NonUnit, &inv_bat)
649 .unwrap()
650 }
651 Err(_) => {
652 warn!("Cholesky decomposition error during variance dervivatives computation");
653 Array2::zeros((b_mat.nrows(), b_mat.ncols()))
654 }
655 };
656
657 let df = self.params.mean.jacobian(&xnorm.row(0)).with_lapack();
658
659 let d_a = df.t().to_owned() - dr.t().dot(&inv_kf);
660 let p4 = d_mat.t().dot(&d_a.t());
662
663 let two = F::cast(2.);
664 let prime_t = (p4 - p2).without_lapack().mapv(|v| two * v);
665
666 let x_std = &self.xt_norm.std;
667 let dvar = (prime_t / x_std).mapv(|v| v * sigma2);
668 dvar.row(0).into_owned()
669 }
670
671 pub fn predict_var_gradients(&self, x: &ArrayBase<impl Data<Elem = F>, Ix2>) -> Array2<F> {
674 let mut derivs = Array::zeros((x.nrows(), x.ncols()));
675 Zip::from(derivs.rows_mut())
676 .and(x.rows())
677 .for_each(|mut der, x| der.assign(&self.predict_var_gradients_single(&x)));
678 derivs
679 }
680}
681
682impl<F, D, Mean, Corr> PredictInplace<ArrayBase<D, Ix2>, Array1<F>>
683 for GaussianProcess<F, Mean, Corr>
684where
685 F: Float,
686 D: Data<Elem = F>,
687 Mean: RegressionModel<F>,
688 Corr: CorrelationModel<F>,
689{
690 fn predict_inplace(&self, x: &ArrayBase<D, Ix2>, y: &mut Array1<F>) {
691 assert_eq!(
692 x.nrows(),
693 y.len(),
694 "The number of data points must match the number of output targets."
695 );
696
697 let values = self.predict(x).expect("GP Prediction");
698 *y = values;
699 }
700
701 fn default_target(&self, x: &ArrayBase<D, Ix2>) -> Array1<F> {
702 Array1::zeros((x.nrows(),))
703 }
704}
705
706#[allow(dead_code)]
708pub struct GpVariancePredictor<'a, F, Mean, Corr>(&'a GaussianProcess<F, Mean, Corr>)
709where
710 F: Float,
711 Mean: RegressionModel<F>,
712 Corr: CorrelationModel<F>;
713
714impl<F, D, Mean, Corr> PredictInplace<ArrayBase<D, Ix2>, Array1<F>>
715 for GpVariancePredictor<'_, F, Mean, Corr>
716where
717 F: Float,
718 D: Data<Elem = F>,
719 Mean: RegressionModel<F>,
720 Corr: CorrelationModel<F>,
721{
722 fn predict_inplace(&self, x: &ArrayBase<D, Ix2>, y: &mut Array1<F>) {
723 assert_eq!(
724 x.nrows(),
725 y.len(),
726 "The number of data points must match the number of output targets."
727 );
728
729 let values = self.0.predict_var(x).expect("GP Prediction");
730 *y = values;
731 }
732
733 fn default_target(&self, x: &ArrayBase<D, Ix2>) -> Array1<F> {
734 Array1::zeros(x.nrows())
735 }
736}
737
738impl<F: Float, Mean: RegressionModel<F>, Corr: CorrelationModel<F>, D: Data<Elem = F>>
739 Fit<ArrayBase<D, Ix2>, ArrayBase<D, Ix1>, GpError> for GpValidParams<F, Mean, Corr>
740{
741 type Object = GaussianProcess<F, Mean, Corr>;
742
743 fn fit(
745 &self,
746 dataset: &DatasetBase<ArrayBase<D, Ix2>, ArrayBase<D, Ix1>>,
747 ) -> Result<Self::Object> {
748 let x = dataset.records();
749 let y = dataset.targets().to_owned().insert_axis(Axis(1));
750
751 if let Some(d) = self.kpls_dim()
752 && *d > x.ncols()
753 {
754 return Err(GpError::InvalidValueError(format!(
755 "Dimension reduction {} should be smaller than actual \
756 training input dimensions {}",
757 d,
758 x.ncols()
759 )));
760 }
761
762 let dim = if let Some(n_components) = self.kpls_dim() {
763 *n_components
764 } else {
765 x.ncols()
766 };
767
768 let (x, y, active, init) = match self.theta_tuning() {
769 ThetaTuning::Fixed(init) | ThetaTuning::Full { init, bounds: _ } => (
770 x.to_owned(),
771 y.to_owned(),
772 (0..dim).collect::<Vec<_>>(),
773 init,
774 ),
775 ThetaTuning::Partial {
776 init,
777 bounds: _,
778 active,
779 } => (x.to_owned(), y.to_owned(), active.to_vec(), init),
780 };
781 let theta0_dim = init.len();
783 let theta0 = if theta0_dim == 1 {
784 Array1::from_elem(dim, init[0])
785 } else if theta0_dim == dim {
786 init.to_owned()
787 } else {
788 panic!(
789 "Initial guess for theta should be either 1-dim or dim of xtrain (w_star.ncols()), got {theta0_dim}"
790 )
791 };
792
793 let xtrain = NormalizedData::new(&x);
794 let ytrain = NormalizedData::new(&y);
795
796 let mut w_star = Array2::eye(x.ncols());
797 if let Some(n_components) = self.kpls_dim() {
798 let ds = Dataset::new(x.to_owned(), y.to_owned());
799 w_star = PlsRegression::params(*n_components).fit(&ds).map_or_else(
800 |e| match e {
801 linfa_pls::PlsError::PowerMethodConstantResidualError() => {
802 Ok(Array2::zeros((x.ncols(), *n_components)))
803 }
804 err => Err(err),
805 },
806 |v| Ok(v.rotations().0.to_owned()),
807 )?;
808 };
809 let x_distances = DistanceMatrix::new(&xtrain.data);
810 let sums = x_distances
811 .d
812 .mapv(|v| num_traits::float::Float::abs(v))
813 .sum_axis(Axis(1));
814 if *sums.min().unwrap() == F::zero() {
815 println!(
816 "Warning: multiple x input features have the same value (at least same row twice)."
817 );
818 }
819 let fx = self.mean().value(&xtrain.data);
820
821 let opt_params = match self.theta_tuning() {
822 ThetaTuning::Fixed(init) => {
823 init.to_owned()
825 }
826 ThetaTuning::Full { init: _, bounds }
827 | ThetaTuning::Partial {
828 init: _,
829 bounds,
830 active: _,
831 } => {
832 let base: f64 = 10.;
833 let objfn = |x: &[f64], _gradient: Option<&mut [f64]>, _params: &mut ()| -> f64 {
834 let mut theta = theta0.to_owned();
835 let xarr = x.iter().map(|v| base.powf(*v)).collect::<Vec<_>>();
836 std::iter::zip(active.clone(), xarr).for_each(|(i, xi)| theta[i] = F::cast(xi));
837
838 for v in theta.iter() {
839 if v.is_nan() {
841 return f64::INFINITY;
843 }
844 }
845 let rxx = self.corr().value(&x_distances.d, &theta, &w_star);
846 match reduced_likelihood(&fx, rxx, &x_distances, &ytrain, self.nugget()) {
847 Ok(r) => unsafe { -(*(&r.0 as *const F as *const f64)) },
848 Err(_) => f64::INFINITY,
849 }
850 };
851
852 let bounds_dim = bounds.len();
855 let bounds = if bounds_dim == 1 {
856 vec![bounds[0]; w_star.ncols()]
857 } else if bounds_dim == w_star.ncols() {
858 bounds.to_vec()
859 } else {
860 panic!(
861 "Bounds for theta should be either 1-dim or dim of xtrain ({}), got {}",
862 w_star.ncols(),
863 bounds_dim
864 )
865 };
866
867 let active_bounds = bounds
869 .iter()
870 .enumerate()
871 .filter(|(i, _)| active.contains(i))
872 .map(|(_, &b)| b)
873 .collect::<Vec<_>>();
874 let (theta_inits, bounds) = prepare_multistart(
875 self.n_start(),
876 &theta0.select(Axis(0), &active),
877 &active_bounds,
878 );
879 debug!("Optimize with multistart theta = {theta_inits:?} and bounds = {bounds:?}");
880 let now = Instant::now();
881 let opt_params = (0..theta_inits.nrows())
882 .into_par_iter()
883 .map(|i| {
884 optimize_params(
885 objfn,
886 &theta_inits.row(i).to_owned(),
887 &bounds,
888 CobylaParams {
889 maxeval: (10 * theta_inits.ncols())
890 .clamp(GP_COBYLA_MIN_EVAL, self.max_eval()),
891 ..CobylaParams::default()
892 },
893 )
894 })
895 .reduce(
896 || (f64::INFINITY, Array::ones((theta_inits.ncols(),))),
897 |a, b| if b.0 < a.0 { b } else { a },
898 );
899 debug!("elapsed optim = {:?}", now.elapsed().as_millis());
900 opt_params.1.mapv(|v| F::cast(base.powf(v)))
901 }
902 };
903
904 let opt_params = match self.theta_tuning() {
906 ThetaTuning::Fixed(_) | ThetaTuning::Full { init: _, bounds: _ } => opt_params,
907 ThetaTuning::Partial {
908 init,
909 bounds: _,
910 active,
911 } => {
912 let mut opt_theta = init.to_owned();
913 std::iter::zip(active.clone(), opt_params)
914 .for_each(|(i, xi)| opt_theta[i] = F::cast(xi));
915 opt_theta
916 }
917 };
918
919 let rxx = self.corr().value(&x_distances.d, &opt_params, &w_star);
920 let (lkh, inner_params) =
921 reduced_likelihood(&fx, rxx, &x_distances, &ytrain, self.nugget())?;
922 Ok(GaussianProcess {
923 theta: opt_params,
924 likelihood: lkh,
925 inner_params,
926 w_star,
927 xt_norm: xtrain,
928 yt_norm: ytrain,
929 training_data: (x.to_owned(), y.to_owned().remove_axis(Axis(1))),
930 params: self.clone(),
931 })
932 }
933}
934
935#[cfg(not(feature = "blas"))]
942fn reduced_likelihood<F: Float>(
943 fx: &ArrayBase<impl Data<Elem = F>, Ix2>,
944 rxx: ArrayBase<impl Data<Elem = F>, Ix2>,
945 x_distances: &DistanceMatrix<F>,
946 ytrain: &NormalizedData<F>,
947 nugget: F,
948) -> Result<(F, GpInnerParams<F>)> {
949 let mut r_mx: Array2<F> = Array2::<F>::eye(x_distances.n_obs).mapv(|v| v + v * nugget);
951 for (i, ij) in x_distances.d_indices.outer_iter().enumerate() {
952 r_mx[[ij[0], ij[1]]] = rxx[[i, 0]];
953 r_mx[[ij[1], ij[0]]] = rxx[[i, 0]];
954 }
955 let fxl = fx;
956 let r_chol = r_mx.cholesky()?;
958 let ft = r_chol.solve_triangular(fxl, UPLO::Lower)?;
960 let (ft_qr_q, ft_qr_r) = ft.qr().unwrap().into_decomp();
961
962 let (_, sv_qr_r, _) = ft_qr_r.svd(false, false).unwrap();
964 let cond_ft = sv_qr_r[sv_qr_r.len() - 1] / sv_qr_r[0];
965 if F::cast(cond_ft) < F::cast(1e-10) {
966 let (_, sv_f, _) = &fxl.svd(false, false).unwrap();
967 let cond_fx = sv_f[0] / sv_f[sv_f.len() - 1];
968 if F::cast(cond_fx) > F::cast(1e15) {
969 return Err(GpError::LikelihoodComputationError(
970 "F is too ill conditioned. Poor combination \
971 of regression model and observations."
972 .to_string(),
973 ));
974 } else {
975 return Err(GpError::LikelihoodComputationError(
977 "ft is too ill conditioned, try another theta again".to_string(),
978 ));
979 }
980 }
981 let yt = r_chol.solve_triangular(&ytrain.data, UPLO::Lower)?;
982
983 let beta = ft_qr_r.solve_triangular_into(ft_qr_q.t().dot(&yt), UPLO::Upper)?;
984 let rho = yt - ft.dot(&beta);
985 let rho_sqr = rho.mapv(|v| v * v).sum_axis(Axis(0));
986
987 let gamma = r_chol.t().solve_triangular_into(rho, UPLO::Upper)?;
988 let n_obs: F = F::cast(x_distances.n_obs);
991
992 let logdet = r_chol.diag().mapv(|v: F| v.log10()).sum() * F::cast(2.) / n_obs;
993
994 let sigma2 = rho_sqr / n_obs;
996 let reduced_likelihood = -n_obs * (sigma2.sum().log10() + logdet);
997
998 Ok((
999 reduced_likelihood,
1000 GpInnerParams {
1001 sigma2: sigma2[0] * ytrain.std[0] * ytrain.std[0],
1002 beta,
1003 gamma,
1004 r_chol,
1005 ft,
1006 ft_qr_r,
1007 },
1008 ))
1009}
1010
1011#[cfg(feature = "blas")]
1013fn reduced_likelihood<F: Float>(
1014 fx: &ArrayBase<impl Data<Elem = F>, Ix2>,
1015 rxx: ArrayBase<impl Data<Elem = F>, Ix2>,
1016 x_distances: &DistanceMatrix<F>,
1017 ytrain: &NormalizedData<F>,
1018 nugget: F,
1019) -> Result<(F, GpInnerParams<F>)> {
1020 let mut r_mx: Array2<F> = Array2::<F>::eye(x_distances.n_obs).mapv(|v| v + v * nugget);
1022 for (i, ij) in x_distances.d_indices.outer_iter().enumerate() {
1023 r_mx[[ij[0], ij[1]]] = rxx[[i, 0]];
1024 r_mx[[ij[1], ij[0]]] = rxx[[i, 0]];
1025 }
1026
1027 let fxl = fx.to_owned().with_lapack();
1028
1029 let r_chol = r_mx.with_lapack().cholesky(UPLO::Lower)?;
1031
1032 let ft = r_chol.solve_triangular(UPLO::Lower, Diag::NonUnit, &fxl)?;
1034 let (ft_qr_q, ft_qr_r) = ft.qr().unwrap();
1035
1036 let (_, sv_qr_r, _) = ft_qr_r.svd(false, false).unwrap();
1038 let cond_ft = sv_qr_r[sv_qr_r.len() - 1] / sv_qr_r[0];
1039 if F::cast(cond_ft) < F::cast(1e-10) {
1040 let (_, sv_f, _) = &fxl.svd(false, false).unwrap();
1041 let cond_fx = sv_f[0] / sv_f[sv_f.len() - 1];
1042 if F::cast(cond_fx) > F::cast(1e15) {
1043 return Err(GpError::LikelihoodComputationError(
1044 "F is too ill conditioned. Poor combination \
1045 of regression model and observations."
1046 .to_string(),
1047 ));
1048 } else {
1049 return Err(GpError::LikelihoodComputationError(
1051 "ft is too ill conditioned, try another theta again".to_string(),
1052 ));
1053 }
1054 }
1055
1056 let yt = r_chol.solve_triangular(
1057 UPLO::Lower,
1058 Diag::NonUnit,
1059 &ytrain.data.to_owned().with_lapack(),
1060 )?;
1061
1062 let beta = ft_qr_r.solve_triangular_into(UPLO::Upper, Diag::NonUnit, ft_qr_q.t().dot(&yt))?;
1063
1064 let rho = yt - ft.dot(&beta);
1065 let rho_sqr = rho.mapv(|v| v * v).sum_axis(Axis(0));
1066 let rho_sqr = rho_sqr.without_lapack();
1067
1068 let gamma = r_chol
1069 .t()
1070 .solve_triangular_into(UPLO::Upper, Diag::NonUnit, rho)?;
1071
1072 let n_obs: F = F::cast(x_distances.n_obs);
1075
1076 let logdet = r_chol
1077 .to_owned()
1078 .without_lapack()
1079 .diag()
1080 .mapv(|v: F| v.log10())
1081 .sum()
1082 * F::cast(2.)
1083 / n_obs;
1084
1085 let sigma2: Array1<F> = rho_sqr / n_obs;
1087 let reduced_likelihood = -n_obs * (sigma2.sum().log10() + logdet);
1088 Ok((
1089 reduced_likelihood,
1090 GpInnerParams {
1091 sigma2: sigma2[0] * ytrain.std[0] * ytrain.std[0],
1092 beta: beta.without_lapack(),
1093 gamma: gamma.without_lapack(),
1094 r_chol: r_chol.without_lapack(),
1095 ft: ft.without_lapack(),
1096 ft_qr_r: ft_qr_r.without_lapack(),
1097 },
1098 ))
1099}
1100
1101pub(crate) fn sample<F: Float>(
1107 x: &ArrayBase<impl Data<Elem = F>, Ix2>,
1108 mean_x: Array2<F>,
1109 cov_x: Array2<F>,
1110 n_traj: usize,
1111 method: GpSamplingMethod,
1112) -> Array2<F> {
1113 let n_eval = x.nrows();
1114 let c = match method {
1115 GpSamplingMethod::Cholesky => {
1116 #[cfg(not(feature = "blas"))]
1117 let c = cov_x.with_lapack().cholesky().unwrap();
1118 #[cfg(feature = "blas")]
1119 let c = cov_x.with_lapack().cholesky(UPLO::Lower).unwrap();
1120 c
1121 }
1122 GpSamplingMethod::EigenValues => {
1123 #[cfg(feature = "blas")]
1124 let (v, w) = cov_x.with_lapack().eigh(UPLO::Lower).unwrap();
1125 #[cfg(not(feature = "blas"))]
1126 let (v, w) = cov_x.with_lapack().eigh_into().unwrap();
1127 let v = v.mapv(F::cast);
1128 let v = v.mapv(|x| {
1129 if x < F::cast(1e-9) {
1131 return F::zero();
1132 }
1133 x.sqrt()
1134 });
1135 let d = Array2::from_diag(&v).with_lapack();
1136 #[cfg(feature = "blas")]
1137 let c = w.dot(&d);
1138 #[cfg(not(feature = "blas"))]
1139 let c = w.dot(&d);
1140 c
1141 }
1142 }
1143 .without_lapack();
1144 let normal = Normal::new(0., 1.).unwrap();
1145 let ary = Array::random((n_eval, n_traj), normal).mapv(|v| F::cast(v));
1146 mean_x.to_owned() + c.dot(&ary)
1147}
1148
1149#[cfg(test)]
1150mod tests {
1151 use super::*;
1152 use approx::{assert_abs_diff_eq, assert_abs_diff_ne};
1153 use argmin_testfunctions::rosenbrock;
1154 use egobox_doe::{Lhs, LhsKind, SamplingMethod};
1155 use linfa::prelude::Predict;
1156 #[cfg(not(feature = "blas"))]
1157 use linfa_linalg::norm::Norm;
1158 use ndarray::{Array, Zip, arr1, arr2, array};
1159 #[cfg(feature = "blas")]
1160 use ndarray_linalg::Norm;
1161 use ndarray_npy::write_npy;
1162 use ndarray_rand::RandomExt;
1163 use ndarray_rand::rand::SeedableRng;
1164 use ndarray_rand::rand_distr::Uniform;
1165 use ndarray_stats::DeviationExt;
1166 use paste::paste;
1167 use rand_xoshiro::Xoshiro256Plus;
1168
1169 #[test]
1170 fn test_constant_function() {
1171 let dim = 3;
1172 let lim = array![[0., 1.]];
1173 let xlimits = lim.broadcast((dim, 2)).unwrap();
1174 let rng = Xoshiro256Plus::seed_from_u64(42);
1175 let nt = 5;
1176 let xt = Lhs::new(&xlimits).with_rng(rng).sample(nt);
1177 let yt = Array::from_vec(vec![3.1; nt]);
1178 let gp = GaussianProcess::<f64, ConstantMean, SquaredExponentialCorr>::params(
1179 ConstantMean::default(),
1180 SquaredExponentialCorr::default(),
1181 )
1182 .theta_init(array![0.1])
1183 .kpls_dim(Some(1))
1184 .fit(&Dataset::new(xt, yt))
1185 .expect("GP fit error");
1186 let rng = Xoshiro256Plus::seed_from_u64(43);
1187 let xtest = Lhs::new(&xlimits).with_rng(rng).sample(nt);
1188 let ytest = gp.predict(&xtest).expect("prediction error");
1189 assert_abs_diff_eq!(Array::from_elem((nt,), 3.1), ytest, epsilon = 1e-6);
1190 }
1191
1192 macro_rules! test_gp {
1193 ($regr:ident, $corr:ident) => {
1194 paste! {
1195
1196 #[test]
1197 fn [<test_gp_ $regr:snake _ $corr:snake >]() {
1198 let xt = array![[0.0], [1.0], [2.0], [3.0], [4.0]];
1199 let xplot = Array::linspace(0., 4., 100).insert_axis(Axis(1));
1200 let yt = array![0.0, 1.0, 1.5, 0.9, 1.0];
1201 let gp = GaussianProcess::<f64, [<$regr Mean>], [<$corr Corr>] >::params(
1202 [<$regr Mean>]::default(),
1203 [<$corr Corr>]::default(),
1204 )
1205 .theta_init(array![0.1])
1206 .fit(&Dataset::new(xt, yt))
1207 .expect("GP fit error");
1208 let yvals = gp
1209 .predict(&arr2(&[[1.0], [3.5]]))
1210 .expect("prediction error");
1211 let expected_y = arr1(&[1.0, 0.9]);
1212 assert_abs_diff_eq!(expected_y, yvals, epsilon = 0.5);
1213
1214 let gpr_vals = gp.predict(&xplot).unwrap();
1215
1216 let yvars = gp
1217 .predict_var(&arr2(&[[1.0], [3.5]]))
1218 .expect("prediction error");
1219 let expected_vars = arr1(&[0., 0.1]);
1220 assert_abs_diff_eq!(expected_vars, yvars, epsilon = 0.5);
1221
1222 let gpr_vars = gp.predict_var(&xplot).unwrap();
1223
1224 let test_dir = "target/tests";
1225 std::fs::create_dir_all(test_dir).ok();
1226
1227 let xplot_file = stringify!([<gp_x_ $regr:snake _ $corr:snake >]);
1228 let file_path = format!("{}/{}.npy", test_dir, xplot_file);
1229 write_npy(file_path, &xplot).expect("x saved");
1230
1231 let gp_vals_file = stringify!([<gp_vals_ $regr:snake _ $corr:snake >]);
1232 let file_path = format!("{}/{}.npy", test_dir, gp_vals_file);
1233 write_npy(file_path, &gpr_vals).expect("gp vals saved");
1234
1235 let gp_vars_file = stringify!([<gp_vars_ $regr:snake _ $corr:snake >]);
1236 let file_path = format!("{}/{}.npy", test_dir, gp_vars_file);
1237 write_npy(file_path, &gpr_vars).expect("gp vars saved");
1238 }
1239 }
1240 };
1241 }
1242
1243 test_gp!(Constant, SquaredExponential);
1244 test_gp!(Constant, AbsoluteExponential);
1245 test_gp!(Constant, Matern32);
1246 test_gp!(Constant, Matern52);
1247
1248 test_gp!(Linear, SquaredExponential);
1249 test_gp!(Linear, AbsoluteExponential);
1250 test_gp!(Linear, Matern32);
1251 test_gp!(Linear, Matern52);
1252
1253 test_gp!(Quadratic, SquaredExponential);
1254 test_gp!(Quadratic, AbsoluteExponential);
1255 test_gp!(Quadratic, Matern32);
1256 test_gp!(Quadratic, Matern52);
1257
1258 fn griewank(x: &Array2<f64>) -> Array1<f64> {
1259 let dim = x.ncols();
1260 let d = Array1::linspace(1., dim as f64, dim).mapv(|v| v.sqrt());
1261 let mut y = Array1::zeros((x.nrows(),));
1262 Zip::from(&mut y).and(x.rows()).for_each(|y, x| {
1263 let s = x.mapv(|v| v * v).sum() / 4000.;
1264 let p = (x.to_owned() / &d)
1265 .mapv(|v| v.cos())
1266 .fold(1., |acc, x| acc * x);
1267 *y = s - p + 1.;
1268 });
1269 y
1270 }
1271
1272 #[test]
1273 fn test_griewank() {
1274 let x = array![[1., 1., 1., 1., 1.], [2., 2., 2., 2., 2.]];
1275 assert_abs_diff_eq!(array![0.72890641, 1.01387135], griewank(&x), epsilon = 1e-8);
1276 }
1277
1278 #[test]
1279 fn test_kpls_griewank() {
1280 let dims = [5]; let nts = [100]; let lim = array![[-600., 600.]];
1283
1284 let test_dir = "target/tests";
1285 std::fs::create_dir_all(test_dir).ok();
1286
1287 (0..dims.len()).for_each(|i| {
1288 let dim = dims[i];
1289 let nt = nts[i];
1290 let xlimits = lim.broadcast((dim, 2)).unwrap();
1291
1292 let prefix = "griewank";
1293 let xfilename = format!("{test_dir}/{prefix}_xt_{nt}x{dim}.npy");
1294 let yfilename = format!("{test_dir}/{prefix}_yt_{nt}x1.npy");
1295
1296 let rng = Xoshiro256Plus::seed_from_u64(42);
1297 let xt = Lhs::new(&xlimits).with_rng(rng).sample(nt);
1298 write_npy(xfilename, &xt).expect("cannot save xt");
1299 let yt = griewank(&xt);
1300 write_npy(yfilename, &yt).expect("cannot save yt");
1301
1302 let gp = GaussianProcess::<f64, ConstantMean, SquaredExponentialCorr>::params(
1303 ConstantMean::default(),
1304 SquaredExponentialCorr::default(),
1305 )
1306 .kpls_dim(Some(3))
1307 .fit(&Dataset::new(xt, yt))
1308 .expect("GP fit error");
1309
1310 let rng = Xoshiro256Plus::seed_from_u64(0);
1311 let xtest = Lhs::new(&xlimits).with_rng(rng).sample(100);
1312 let ytest = gp.predict(&xtest).expect("prediction error");
1314 let ytrue = griewank(&xtest);
1315
1316 let nrmse = (ytrue.to_owned() - &ytest).norm_l2() / ytrue.norm_l2();
1317 println!(
1318 "diff={} ytrue={} nrsme={}",
1319 (ytrue.to_owned() - &ytest).norm_l2(),
1320 ytrue.norm_l2(),
1321 nrmse
1322 );
1323 assert_abs_diff_eq!(nrmse, 0., epsilon = 1e-2);
1324 });
1325 }
1326
1327 fn tensor_product_exp(x: &ArrayBase<impl Data<Elem = f64>, Ix2>) -> Array1<f64> {
1328 x.mapv(|v| v.exp()).map_axis(Axis(1), |row| row.product())
1329 }
1330
1331 #[test]
1332 fn test_kpls_tp_exp() {
1333 let dim = 3;
1334 let nt = 300;
1335 let lim = array![[-1., 1.]];
1336 let xlimits = lim.broadcast((dim, 2)).unwrap();
1337 let rng = Xoshiro256Plus::seed_from_u64(42);
1338 let xt = Lhs::new(&xlimits).with_rng(rng).sample(nt);
1339 let yt = tensor_product_exp(&xt);
1340
1341 let gp = GaussianProcess::<f64, ConstantMean, SquaredExponentialCorr>::params(
1342 ConstantMean::default(),
1343 SquaredExponentialCorr::default(),
1344 )
1345 .kpls_dim(Some(1))
1346 .fit(&Dataset::new(xt, yt))
1347 .expect("GP training");
1348
1349 let xv = Lhs::new(&xlimits).sample(100);
1350 let yv = tensor_product_exp(&xv);
1351
1352 let ytest = gp.predict(&xv).unwrap();
1353 let err = ytest.l2_dist(&yv).unwrap() / yv.norm_l2();
1354 assert_abs_diff_eq!(err, 0., epsilon = 2e-2);
1355 }
1356
1357 fn rosenb(x: &ArrayBase<impl Data<Elem = f64>, Ix2>) -> Array1<f64> {
1358 let mut y: Array1<f64> = Array1::zeros((x.nrows(),));
1359 Zip::from(&mut y).and(x.rows()).par_for_each(|yi, xi| {
1360 *yi = rosenbrock(&xi.to_vec());
1361 });
1362 y
1363 }
1364
1365 #[test]
1366 fn test_kpls_rosenb() {
1367 let dim = 20;
1368 let nt = 30;
1369 let lim = array![[-1., 1.]];
1370 let xlimits = lim.broadcast((dim, 2)).unwrap();
1371 let rng = Xoshiro256Plus::seed_from_u64(42);
1372 let xt = Lhs::new(&xlimits).with_rng(rng).sample(nt);
1373 let yt = rosenb(&xt);
1374
1375 let gp = GaussianProcess::<f64, ConstantMean, Matern32Corr>::params(
1376 ConstantMean::default(),
1377 Matern52Corr::default(),
1378 )
1379 .kpls_dim(Some(1))
1380 .fit(&Dataset::new(xt.to_owned(), yt))
1381 .expect("GP training");
1382
1383 let rng2 = Xoshiro256Plus::seed_from_u64(41);
1384 let xv = Lhs::new(&xlimits).with_rng(rng2).sample(300);
1385 let yv = rosenb(&xv);
1386
1387 let ytest = gp.predict(&xv).expect("GP prediction");
1388 let err = ytest.l2_dist(&yv).unwrap() / yv.norm_l2();
1389 assert_abs_diff_eq!(err, 0., epsilon = 4e-1);
1390
1391 let var = GpVariancePredictor(&gp).predict(&xt);
1392 assert_abs_diff_eq!(var, Array1::zeros(nt), epsilon = 2e-1);
1393 }
1394
1395 fn sphere(x: &Array2<f64>) -> Array1<f64> {
1396 (x * x).sum_axis(Axis(1))
1397 }
1398
1399 fn dsphere(x: &Array2<f64>) -> Array2<f64> {
1400 x.mapv(|v| 2. * v)
1401 }
1402
1403 fn norm1(x: &Array2<f64>) -> Array1<f64> {
1404 x.mapv(|v| v.abs()).sum_axis(Axis(1)).to_owned()
1405 }
1406
1407 fn dnorm1(x: &Array2<f64>) -> Array2<f64> {
1408 x.mapv(|v| if v > 0. { 1. } else { -1. })
1409 }
1410
1411 macro_rules! test_gp_derivatives {
1412 ($regr:ident, $corr:ident, $func:ident, $limit:expr_2021, $nt:expr_2021) => {
1413 paste! {
1414
1415 #[test]
1416 fn [<test_gp_derivatives_ $regr:snake _ $corr:snake>]() {
1417 let mut rng = Xoshiro256Plus::seed_from_u64(42);
1418 let xt = egobox_doe::Lhs::new(&array![[-$limit, $limit], [-$limit, $limit]])
1419 .kind(egobox_doe::LhsKind::CenteredMaximin)
1420 .with_rng(rng.clone())
1421 .sample($nt);
1422
1423 let yt = [<$func>](&xt);
1424 let gp = GaussianProcess::<f64, [<$regr Mean>], [<$corr Corr>] >::params(
1425 [<$regr Mean>]::default(),
1426 [<$corr Corr>]::default(),
1427 )
1428 .fit(&Dataset::new(xt, yt))
1429 .expect("GP fitting");
1430
1431 let x = Array::random_using((2,), Uniform::new(-$limit, $limit), &mut rng);
1432 let xa: f64 = x[0];
1434 let xb: f64 = x[1];
1435 let e = 1e-5;
1436
1437 let x = array![
1438 [xa, xb],
1439 [xa + e, xb],
1440 [xa - e, xb],
1441 [xa, xb + e],
1442 [xa, xb - e]
1443 ];
1444
1445 let y_pred = gp.predict(&x).unwrap();
1446 println!("value at [{},{}] = {}", xa, xb, y_pred);
1447 let y_deriv = gp.predict_gradients(&x);
1448 println!("deriv at [{},{}] = {}", xa, xb, y_deriv);
1449 let true_deriv = [<d $func>](&array![[xa, xb]]);
1450 println!("true deriv at [{},{}] = {}", xa, xb, true_deriv);
1451 println!("jacob = at [{},{}] = {}", xa, xb, gp.predict_jacobian(&array![xa, xb]));
1452
1453 let diff_g = (y_pred[1] - y_pred[2]) / (2. * e);
1454 let diff_d = (y_pred[3] - y_pred[4]) / (2. * e);
1455
1456 if (diff_g-true_deriv[[0, 0]]).abs() < 10. {
1458 assert_rel_or_abs_error(y_deriv[[0, 0]], diff_g);
1459 }
1460 if (diff_d-true_deriv[[0, 1]]).abs() < 10. {
1461 assert_rel_or_abs_error(y_deriv[[0, 1]], diff_d);
1462 }
1463 }
1464 }
1465 };
1466 }
1467
1468 test_gp_derivatives!(Constant, SquaredExponential, sphere, 10., 10);
1469 test_gp_derivatives!(Linear, SquaredExponential, sphere, 10., 10);
1470 test_gp_derivatives!(Quadratic, SquaredExponential, sphere, 10., 10);
1471 test_gp_derivatives!(Constant, AbsoluteExponential, sphere, 10., 10);
1472 test_gp_derivatives!(Linear, AbsoluteExponential, norm1, 10., 16);
1473 test_gp_derivatives!(Quadratic, AbsoluteExponential, norm1, 10., 16);
1474 test_gp_derivatives!(Constant, Matern32, norm1, 10., 16);
1475 test_gp_derivatives!(Linear, Matern32, norm1, 10., 16);
1476 test_gp_derivatives!(Quadratic, Matern32, sphere, 10., 16);
1477 test_gp_derivatives!(Constant, Matern52, norm1, 10., 16);
1478 test_gp_derivatives!(Linear, Matern52, norm1, 10., 16);
1479 test_gp_derivatives!(Quadratic, Matern52, sphere, 10., 10);
1480
1481 #[allow(unused_macros)]
1482 macro_rules! test_gp_variance_derivatives {
1483 ($regr:ident, $corr:ident, $func:ident, $limit:expr_2021, $nt:expr_2021) => {
1484 paste! {
1485
1486 #[test]
1487 fn [<test_gp_variance_derivatives_ $regr:snake _ $corr:snake _ $func:snake>]() {
1488 let mut rng = Xoshiro256Plus::seed_from_u64(42);
1489 let xt = egobox_doe::Lhs::new(&array![[-$limit, $limit], [-$limit, $limit]]).with_rng(rng.clone()).sample($nt);
1490 let yt = [<$func>](&xt);
1491 println!(stringify!(<$func>));
1492
1493 let gp = GaussianProcess::<f64, [<$regr Mean>], [<$corr Corr>] >::params(
1494 [<$regr Mean>]::default(),
1495 [<$corr Corr>]::default(),
1496 )
1497 .fit(&Dataset::new(xt, yt))
1498 .expect("GP fitting");
1499
1500 for _ in 0..10 {
1501 let x = Array::random_using((2,), Uniform::new(-$limit, $limit), &mut rng);
1502 let xa: f64 = x[0];
1503 let xb: f64 = x[1];
1504 let e = 1e-5;
1505
1506 let x = array![
1507 [xa, xb],
1508 [xa + e, xb],
1509 [xa - e, xb],
1510 [xa, xb + e],
1511 [xa, xb - e]
1512 ];
1513 println!("****************************************");
1514 let y_pred = gp.predict(&x).unwrap();
1515 println!("value at [{},{}] = {}", xa, xb, y_pred);
1516 let y_deriv = gp.predict_gradients(&x);
1517 println!("deriv at [{},{}] = {}", xa, xb, y_deriv);
1518 let y_pred = gp.predict_var(&x).unwrap();
1519 println!("variance at [{},{}] = {}", xa, xb, y_pred);
1520 let y_deriv = gp.predict_var_gradients(&x);
1521 println!("variance deriv at [{},{}] = {}", xa, xb, y_deriv);
1522
1523 let diff_g = (y_pred[1] - y_pred[2]) / (2. * e);
1524 let diff_d = (y_pred[3] - y_pred[4]) / (2. * e);
1525
1526 assert_rel_or_abs_error(y_deriv[[0, 0]], diff_g);
1527 assert_rel_or_abs_error(y_deriv[[0, 1]], diff_d);
1528 }
1529 }
1530 }
1531 };
1532 }
1533
1534 test_gp_variance_derivatives!(Constant, SquaredExponential, sphere, 10., 100);
1535 test_gp_variance_derivatives!(Linear, SquaredExponential, sphere, 10., 100);
1536 test_gp_variance_derivatives!(Quadratic, SquaredExponential, sphere, 10., 100);
1537 #[cfg(not(feature = "nlopt"))]
1539 test_gp_variance_derivatives!(Constant, AbsoluteExponential, norm1, 10., 100);
1540 test_gp_variance_derivatives!(Linear, AbsoluteExponential, norm1, 1., 50);
1541 test_gp_variance_derivatives!(Quadratic, AbsoluteExponential, sphere, 10., 100);
1542 test_gp_variance_derivatives!(Constant, Matern32, sphere, 10., 100);
1543 test_gp_variance_derivatives!(Linear, Matern32, norm1, 1., 50);
1544 test_gp_variance_derivatives!(Quadratic, Matern32, sphere, 10., 100);
1545 test_gp_variance_derivatives!(Constant, Matern52, sphere, 10., 100);
1546 test_gp_variance_derivatives!(Linear, Matern52, norm1, 1., 50);
1547 test_gp_variance_derivatives!(Quadratic, Matern52, sphere, 10., 100);
1548
1549 #[test]
1550 fn test_variance_derivatives() {
1551 let xt = egobox_doe::FullFactorial::new(&array![[-10., 10.], [-10., 10.]]).sample(10);
1552 let yt = sphere(&xt);
1553
1554 let gp = GaussianProcess::<f64, ConstantMean, SquaredExponentialCorr>::params(
1555 ConstantMean::default(),
1556 SquaredExponentialCorr::default(),
1557 )
1558 .fit(&Dataset::new(xt, yt))
1559 .expect("GP fitting");
1560
1561 for _ in 0..20 {
1562 let mut rng = Xoshiro256Plus::seed_from_u64(42);
1563 let x = Array::random_using((2,), Uniform::new(-10., 10.), &mut rng);
1564 let xa: f64 = x[0];
1565 let xb: f64 = x[1];
1566 let e = 1e-5;
1567
1568 let x = array![
1569 [xa, xb],
1570 [xa + e, xb],
1571 [xa - e, xb],
1572 [xa, xb + e],
1573 [xa, xb - e]
1574 ];
1575 let y_pred = gp.predict_var(&x).unwrap();
1576 println!("variance at [{xa},{xb}] = {y_pred}");
1577 let y_deriv = gp.predict_var_gradients(&x);
1578 println!("variance deriv at [{xa},{xb}] = {y_deriv}");
1579
1580 let diff_g = (y_pred[1] - y_pred[2]) / (2. * e);
1581 let diff_d = (y_pred[3] - y_pred[4]) / (2. * e);
1582
1583 if y_pred[0].abs() > 1e-1 && y_pred[0].abs() > 1e-1 {
1584 assert_rel_or_abs_error(y_deriv[[0, 0]], diff_g);
1586 }
1587 if y_pred[0].abs() > 1e-1 && y_pred[0].abs() > 1e-1 {
1588 assert_rel_or_abs_error(y_deriv[[0, 1]], diff_d);
1590 }
1591 }
1592 }
1593
1594 #[test]
1595 fn test_fixed_theta() {
1596 let xt = array![[0.0], [1.0], [2.0], [3.0], [4.0]];
1597 let yt = array![0.0, 1.0, 1.5, 0.9, 1.0];
1598 let gp = Kriging::params()
1599 .fit(&Dataset::new(xt.clone(), yt.clone()))
1600 .expect("GP fit error");
1601 let default = ThetaTuning::default();
1602 assert_abs_diff_ne!(*gp.theta(), default.init());
1603 let expected = gp.theta();
1604
1605 let gp = Kriging::params()
1606 .theta_tuning(ThetaTuning::Fixed(expected.clone()))
1607 .fit(&Dataset::new(xt, yt))
1608 .expect("GP fit error");
1609 assert_abs_diff_eq!(*gp.theta(), expected);
1610 }
1611
1612 fn x2sinx(x: &Array2<f64>) -> Array1<f64> {
1613 ((x * x) * (x).mapv(|v| v.sin())).remove_axis(Axis(1))
1614 }
1615
1616 #[test]
1617 fn test_sampling() {
1618 let xdoe = array![[-8.5], [-4.0], [-3.0], [-1.0], [4.0], [7.5]];
1619 let ydoe = x2sinx(&xdoe);
1620 let krg = Kriging::<f64>::params()
1621 .fit(&Dataset::new(xdoe, ydoe))
1622 .expect("Kriging training");
1623 let n_plot = 35;
1624 let n_traj = 10;
1625 let (x_min, x_max) = (-10., 10.);
1626 let x = Array::linspace(x_min, x_max, n_plot)
1627 .into_shape_with_order((n_plot, 1))
1628 .unwrap();
1629 let trajs = krg.sample(&x, n_traj);
1630 assert_eq!(&[n_plot, n_traj], trajs.shape())
1631 }
1632
1633 #[test]
1634 fn test_sampling_eigen() {
1635 let xdoe = array![[-8.5], [-4.0], [-3.0], [-1.0], [4.0], [7.5]];
1636 let ydoe = x2sinx(&xdoe);
1637 let krg = Kriging::<f64>::params()
1638 .fit(&Dataset::new(xdoe, ydoe))
1639 .expect("Kriging training");
1640 let n_plot = 500;
1641 let n_traj = 10;
1642 let (x_min, x_max) = (-10., 10.);
1643 let x = Array::linspace(x_min, x_max, n_plot)
1644 .into_shape_with_order((n_plot, 1))
1645 .unwrap();
1646 let trajs = krg.sample_eig(&x, n_traj);
1647 assert_eq!(&[n_plot, n_traj], trajs.shape());
1648 assert!(!trajs.fold(false, |acc, v| acc || v.is_nan())); }
1650
1651 fn assert_rel_or_abs_error(y_deriv: f64, fdiff: f64) {
1652 println!("analytic deriv = {y_deriv}, fdiff = {fdiff}");
1653 if fdiff.abs() < 1. {
1654 let atol = 1.;
1655 println!("Check absolute error: abs({y_deriv}) should be < {atol}");
1656 assert_abs_diff_eq!(y_deriv, 0.0, epsilon = atol); } else {
1658 let rtol = 6e-1;
1659 let rel_error = (y_deriv - fdiff).abs() / fdiff.abs(); println!("Check relative error: {rel_error} should be < {rtol}");
1661 assert_abs_diff_eq!(rel_error, 0.0, epsilon = rtol);
1662 }
1663 }
1664
1665 fn sin_linear(x: &Array2<f64>) -> Array2<f64> {
1666 let x1 = x.column(0).to_owned().mapv(|v| v.sin());
1668 let x2 = x.column(0).mapv(|v| 2. * v) + x.column(1).mapv(|v| 5. * v);
1669 (x1 + x2)
1670 .mapv(|v| v + 10.)
1671 .into_shape_with_order((x.nrows(), 1))
1672 .unwrap()
1673 }
1674
1675 #[test]
1676 fn test_bug_var_derivatives() {
1677 let _xt = egobox_doe::Lhs::new(&array![[-5., 10.], [-5., 10.]])
1678 .kind(LhsKind::Centered)
1679 .sample(12);
1680 let _yt = sin_linear(&_xt);
1681
1682 let xt = array![
1683 [6.875, -4.375],
1684 [-3.125, 1.875],
1685 [1.875, -1.875],
1686 [-4.375, 3.125],
1687 [8.125, 9.375],
1688 [4.375, 4.375],
1689 [0.625, 0.625],
1690 [9.375, 6.875],
1691 [5.625, 8.125],
1692 [-0.625, -3.125],
1693 [3.125, 5.625],
1694 [-1.875, -0.625]
1695 ];
1696 let yt = array![
1697 2.43286801,
1698 13.10840811,
1699 5.32908578,
1700 17.81862219,
1701 74.08849877,
1702 39.68137781,
1703 14.96009727,
1704 63.17475741,
1705 61.26331775,
1706 -7.46009727,
1707 44.39159189,
1708 2.17091422,
1709 ];
1710
1711 let gp = GaussianProcess::<f64, ConstantMean, SquaredExponentialCorr>::params(
1712 ConstantMean::default(),
1713 SquaredExponentialCorr::default(),
1714 )
1715 .theta_tuning(ThetaTuning::Fixed(array![
1716 f64::sqrt(2. * 0.0437386),
1717 f64::sqrt(2. * 0.00697978)
1718 ]))
1719 .fit(&Dataset::new(xt, yt))
1720 .expect("GP fitting");
1721
1722 let e = 5e-6;
1723 let xa = -1.3;
1724 let xb = 2.5;
1725 let x = array![
1726 [xa, xb],
1727 [xa + e, xb],
1728 [xa - e, xb],
1729 [xa, xb + e],
1730 [xa, xb - e]
1731 ];
1732 let y_pred = gp.predict_var(&x).unwrap();
1733 let y_deriv = gp.predict_var_gradients(&array![[xa, xb]]);
1734 let diff_g = (y_pred[1] - y_pred[2]) / (2. * e);
1735 let diff_d = (y_pred[3] - y_pred[4]) / (2. * e);
1736
1737 assert_abs_diff_eq!(y_deriv[[0, 0]], diff_g, epsilon = 1e-5);
1738 assert_abs_diff_eq!(y_deriv[[0, 1]], diff_d, epsilon = 1e-5);
1739 }
1740}