Dakota
Version 6.21
Explore and Predict with Confidence
|
The GaussianProcess constructs a Gaussian Process regressor surrogate given a matrix of data. More...
Public Member Functions | |
GaussianProcess () | |
Constructor that uses defaultConfigOptions and does not build. | |
GaussianProcess (const ParameterList ¶m_list) | |
Constructor that sets configOptions but does not build. More... | |
GaussianProcess (const std::string ¶m_list_yaml_filename) | |
Constructor for the GaussianProcess that sets configOptions but does not build the GP. More... | |
GaussianProcess (const MatrixXd &samples, const MatrixXd &response, const ParameterList ¶m_list) | |
Constructor for the GaussianProcess that sets configOptions and builds the GP. More... | |
GaussianProcess (const MatrixXd &samples, const MatrixXd &response, const std::string ¶m_list_yaml_filename) | |
Constructor for the GaussianProcess that sets configOptions and builds the GP. More... | |
~GaussianProcess () | |
Default destructor. | |
void | build (const MatrixXd &eval_points, const MatrixXd &response) override |
Build the GP using specified build data. More... | |
VectorXd | value (const MatrixXd &eval_points) override |
Evaluate the scalar Gaussian Process at a set of prediction points. More... | |
MatrixXd | gradient (const MatrixXd &eval_points) override |
Evaluate the gradient of the scalar Gaussian process at a set of prediction points. More... | |
MatrixXd | hessian (const MatrixXd &eval_point) override |
Evaluate the Hessian of the scalar Gaussian process at a single point. More... | |
MatrixXd | covariance (const MatrixXd &eval_points, const int qoi) |
Evaluate the covariance matrix for the Gaussian Process at a set of prediction points for a single QoI index. More... | |
MatrixXd | covariance (const MatrixXd &eval_points) |
Evaluate the covariance matrix for the Gaussian Process at a set of prediction points for QoI index 0. More... | |
VectorXd | variance (const MatrixXd &eval_points, const int qoi) |
Evaluate the variance of the Gaussian Process at a set of prediction points for a given QoI index. More... | |
VectorXd | variance (const MatrixXd &eval_points) |
Evaluate the variance of the Gaussian Process at a set of prediction points for QoI index 0. More... | |
void | negative_marginal_log_likelihood (bool compute_grad, bool compute_gram, double &obj_value, VectorXd &obj_gradient) |
Evaluate the negative marginal loglikelihood and its gradient. More... | |
void | setup_hyperparameter_bounds (VectorXd &sigma_bounds, MatrixXd &length_scale_bounds, VectorXd &nugget_bounds) |
Initialize the hyperparameter bounds for MLE from values in configOptions. More... | |
int | get_num_opt_variables () |
Get the number of optimization variables. More... | |
int | get_num_variables () const |
Get the dimension of the feature space. More... | |
VectorXd | get_objective_function_history () |
Get the history of objective function values from MLE with restarts. More... | |
MatrixXd | get_objective_gradient_history () |
Get the history of objective function gradients from MLE with restarts. More... | |
MatrixXd | get_theta_history () |
Get the history of hyperparameter values from MLE with restarts. More... | |
void | set_opt_params (const std::vector< double > &opt_params) |
Update the vector of optimization parameters. More... | |
std::shared_ptr< Surrogate > | clone () const override |
clone derived Surrogate class for use in cross-validation | |
Public Member Functions inherited from Surrogate | |
Surrogate () | |
Constructor that uses defaultConfigOptions and does not build. | |
Surrogate (const ParameterList ¶m_list) | |
Constructor that sets configOptions but does not build. More... | |
Surrogate (const MatrixXd &samples, const MatrixXd &response, const ParameterList ¶m_list) | |
Constructor for the Surrogate that sets configOptions and builds the surrogate (does nothing in the base class). More... | |
virtual | ~Surrogate () |
Default destructor. | |
virtual VectorXd | values (const MatrixXd &eval_points) |
Evaluate the Surrogate at a set of prediction points for field QoIs. More... | |
MatrixXd | gradients (const MatrixXd &eval_points) |
Evaluate the gradient of the field Surrogate at a set of prediction points. More... | |
virtual MatrixXd | hessians (const MatrixXd &eval_point) |
Evaluate the Hessian of the field Surrogate at a set of prediction points. More... | |
void | variable_labels (const std::vector< std::string > &var_labels) |
Set the variable/feature names. More... | |
const std::vector< std::string > & | variable_labels () const |
Get the (possibly empty) variable/feature names. More... | |
void | response_labels (const std::vector< std::string > &resp_labels) |
Set the response/QoI names. More... | |
const std::vector< std::string > & | response_labels () const |
Get the (possibly empty) response/QoI names. More... | |
void | set_options (const ParameterList &options) |
Set the Surrogate's configOptions. More... | |
void | get_options (ParameterList &options) |
Get the Surrogate's configOptions. More... | |
void | print_options () |
Print the Surrogate's configOptions. | |
virtual bool | diagnostics_available () |
VectorXd | evaluate_metrics (const StringArray &mnames, const MatrixXd &points, const MatrixXd &ref_values) |
Evalute metrics at specified points (within surrogates) | |
VectorXd | cross_validate (const MatrixXd &samples, const MatrixXd &response, const StringArray &mnames, const int num_folds=5, const int seed=20) |
Perform K-folds cross-validation (within surrogates) | |
template<typename DerivedSurr > | |
void | save (const DerivedSurr &surr_out, const std::string &outfile, const bool binary) |
Serialize a derived (i.e. non-base) surrogate model. More... | |
template<typename DerivedSurr > | |
void | load (const std::string &infile, const bool binary, DerivedSurr &surr_in) |
Load a derived (i.e. non-base) surrogate model. More... | |
Private Member Functions | |
void | default_options () override |
Construct and populate the defaultConfigOptions. | |
void | compute_build_dists () |
Compute squared distances between the scaled build points. | |
void | compute_pred_dists (const MatrixXd &scaled_pred_pts) |
Compute distances between build and prediction points. This includes build-prediction and prediction-prediction distance matrices. More... | |
void | compute_gram (const std::vector< MatrixXd > &dists2, bool add_nugget, bool compute_derivs, MatrixXd &gram) |
Compute a Gram matrix given a vector of squared distances and optionally compute its derivatives and/or adds nugget terms. More... | |
void | generate_initial_guesses (const VectorXd &sigma_bounds, const MatrixXd &length_scale_bounds, const VectorXd &nugget_bounds, const int num_restarts, const int seed, MatrixXd &initial_guesses) |
Randomly generate initial guesses for the optimization routine. More... | |
void | setup_default_optimization_params (Teuchos::RCP< ParameterList > rol_params) |
Set the default optimization parameters for ROL for GP hyperparameter estimation. More... | |
template<class Archive > | |
void | serialize (Archive &archive, const unsigned int version) |
Serializer for save/load. | |
Private Attributes | |
double | fixedNuggetValue |
Small constant added to the diagonal of the Gram matrix to avoid ill-conditioning. | |
MatrixXd | eyeMatrix |
Identity matrix for the build points space. | |
MatrixXd | basisMatrix |
Basis matrix for the sample points in polynomial regression. | |
MatrixXd | targetValues |
Target values for the surrogate dataset. | |
MatrixXd | scaledBuildPoints |
The scaled build points for the surrogate dataset. | |
VectorXd | thetaValues |
Vector of log-space hyperparameters. | |
VectorXd | betaValues |
Vector of polynomial coefficients. | |
double | estimatedNuggetValue |
Estimated nugget term. | |
VectorXd | bestThetaValues |
Vector of best hyperparameters from MLE with restarts. | |
VectorXd | bestBetaValues |
Vector of best polynomial coefficients from MLE with restarts. | |
double | bestEstimatedNuggetValue |
Best estimated nugget value from MLE with restarts. | |
VectorXd | objectiveFunctionHistory |
Final objective function values for each optimization run. | |
MatrixXd | objectiveGradientHistory |
Final objective function gradients for each optimization run. | |
MatrixXd | thetaHistory |
Final hyperparameter values for each optimization run. | |
MatrixXd | GramMatrix |
Gram matrix for the build points. | |
VectorXd | trendTargetResidual |
Difference between target values and trend predictions. | |
VectorXd | GramResidualSolution |
Cholesky solve for Gram matrix with trendTargetResidual rhs. | |
std::vector< MatrixXd > | GramMatrixDerivs |
Derivatives of the Gram matrix w.r.t. the hyperparameters. | |
std::vector< MatrixXd > | cwiseDists2 |
Squared component-wise distances between points in the surrogate dataset. | |
std::vector< MatrixXd > | cwiseMixedDists |
Component-wise distances between prediction and build points. | |
std::vector< MatrixXd > | cwiseMixedDists2 |
Squared component-wise distances between prediction and build points. | |
std::vector< MatrixXd > | cwisePredDists2 |
Component-wise distances between prediction points. | |
Eigen::LDLT< MatrixXd > | CholFact |
Pivoted Cholesky factorization. | |
bool | hasBestCholFact |
Flag for recomputation of the best Cholesky factorization. | |
MatrixXd | predGramMatrix |
Gram matrix for the prediction points. | |
MatrixXd | predMixedGramMatrix |
Gram matrix for the mixed prediction/build points. | |
MatrixXd | predCovariance |
Covariance matrix for the prediction points. | |
MatrixXd | predBasisMatrix |
Polynomial basis matrix for the prediction points. | |
std::shared_ptr< PolynomialRegression > | polyRegression |
PolynomialRegression for trend function. | |
std::string | kernel_type |
Kernel type. | |
std::shared_ptr< Kernel > | kernel |
Kernel. | |
const double | betaBound = 1.0e20 |
Large constant for polynomial coefficients upper/lower bounds. | |
bool | estimateTrend |
Bool for polynomial trend (i.e. semi-parametric GP) estimation. | |
int | numPolyTerms = 0 |
Number of terms in polynomial trend. | |
int | numNuggetTerms = 0 |
Number of terms for the (estimated) nugget parameter. | |
bool | estimateNugget |
Bool for nugget estimation. | |
int | verbosity |
Verbosity level. | |
double | bestObjFunValue = std::numeric_limits<double>::max() |
Final objective function value. | |
const double | PI = 3.14159265358979323846 |
Numerical constant – needed for negative marginal log-likelihood. | |
Friends | |
class | boost::serialization::access |
Allow serializers access to private class data. | |
Additional Inherited Members | |
Static Public Member Functions inherited from Surrogate | |
template<typename SurrHandle > | |
static void | save (const SurrHandle &surr_out, const std::string &outfile, const bool binary) |
serialize Surrogate to file (typically through shared_ptr<Surrogate>, but Derived& or Derived* okay too) | |
template<typename SurrHandle > | |
static void | load (const std::string &infile, const bool binary, SurrHandle &surr_in) |
serialize Surrogate from file (typically through shared_ptr<Surrogate>, but Derived& or Derived* okay too) | |
static std::shared_ptr< Surrogate > | load (const std::string &infile, const bool binary) |
serialize Surrogate from file through pointer to base class (must have been saved via same data type) | |
Public Attributes inherited from Surrogate | |
util::DataScaler | dataScaler |
DataScaler class for a Surrogate's build samples. | |
double | responseOffset = 0. |
Response offset. | |
double | responseScaleFactor = 1. |
Response scale factor. | |
Protected Attributes inherited from Surrogate | |
int | numSamples |
Number of samples in the Surrogate's build samples. | |
int | numVariables |
Number of features/variables in the Surrogate's build samples. | |
std::vector< std::string > | variableLabels |
Names of the variables/features; need not be populated. | |
int | numQOI |
Number of quantities of interest predicted by the surrogate. For scalar-valued surrogates numQOI = 1. | |
std::vector< std::string > | responseLabels |
Names of the responses/QoIs; need not be populated. | |
ParameterList | defaultConfigOptions |
Default Key/value options to configure the surrogate. | |
ParameterList | configOptions |
Key/value options to configure the surrogate - will override defaultConfigOptions. | |
The GaussianProcess constructs a Gaussian Process regressor surrogate given a matrix of data.
The Gaussian Process (GP) uses an anisotropic squared exponential kernel with a constant multiplicative scaling factor. This yields a total of num_features + 1 kernel hyperparameters. These parameters are internally transformed to a log-space vector (theta) for optimization and evaluation of the GP. Polynomial trend and nugget estimation are optional.
The GP's parameters are determined through maximum likelihood estimation (MLE) via minimization of the negative marginal log-likelihood function. ROL's implementation of L-BFGS-B is used to solve the optimization problem, and the algorithm may be run from multiple random initial guesses to increase the chance of finding the global minimum.
Once the GP is constructed its mean, variance, and covariance matrix can be computed for a set of prediction points. Gradients and Hessians are available.
GaussianProcess | ( | const ParameterList & | param_list | ) |
Constructor that sets configOptions but does not build.
[in] | param_list | List that overrides entries in defaultConfigOptions. |
References Surrogate::configOptions, GaussianProcess::default_options(), and Surrogate::defaultConfigOptions.
GaussianProcess | ( | const std::string & | param_list_yaml_filename | ) |
Constructor for the GaussianProcess that sets configOptions but does not build the GP.
[in] | param_list_yaml_filename | A ParameterList file (relative to the location of the Dakota input file) that overrides entries in defaultConfigOptions. |
References Surrogate::configOptions, GaussianProcess::default_options(), and Surrogate::defaultConfigOptions.
GaussianProcess | ( | const MatrixXd & | samples, |
const MatrixXd & | response, | ||
const ParameterList & | param_list | ||
) |
Constructor for the GaussianProcess that sets configOptions and builds the GP.
[in] | samples | Matrix of data for surrogate construction - (num_samples by num_features) |
[in] | response | Vector of targets for surrogate construction - (num_samples by num_qoi = 1; only 1 response is supported currently). |
[in] | param_list | List that overrides entries in defaultConfigOptions |
References GaussianProcess::build(), Surrogate::configOptions, and GaussianProcess::default_options().
GaussianProcess | ( | const MatrixXd & | samples, |
const MatrixXd & | response, | ||
const std::string & | param_list_yaml_filename | ||
) |
Constructor for the GaussianProcess that sets configOptions and builds the GP.
[in] | samples | Matrix of data for surrogate construction - (num_samples by num_features) |
[in] | response | Vector of targets for surrogate construction - (num_samples by num_qoi = 1; only 1 response is supported currently). |
[in] | param_list_yaml_filename | A ParameterList file (relative to the location of the Dakota input file) that overrides entries in defaultConfigOptions. |
References GaussianProcess::build(), Surrogate::configOptions, and GaussianProcess::default_options().
Build the GP using specified build data.
[in] | eval_points | Matrix of data for surrogate construction - (num_samples by num_features) |
[in] | response | Vector of targets for surrogate construction - (num_samples by num_qoi = 1; only 1 response is supported currently). |
Implements Surrogate.
References GaussianProcess::basisMatrix, GaussianProcess::bestBetaValues, GaussianProcess::bestEstimatedNuggetValue, GaussianProcess::bestObjFunValue, GaussianProcess::bestThetaValues, GaussianProcess::betaBound, GaussianProcess::betaValues, GaussianProcess::CholFact, GaussianProcess::compute_build_dists(), GaussianProcess::compute_gram(), Surrogate::configOptions, GaussianProcess::cwiseDists2, Surrogate::dataScaler, Surrogate::defaultConfigOptions, GaussianProcess::estimatedNuggetValue, GaussianProcess::estimateNugget, GaussianProcess::estimateTrend, GaussianProcess::eyeMatrix, GaussianProcess::fixedNuggetValue, GaussianProcess::generate_initial_guesses(), GaussianProcess::GramMatrix, GaussianProcess::GramMatrixDerivs, GaussianProcess::hasBestCholFact, GaussianProcess::kernel, dakota::surrogates::kernel_factory(), GaussianProcess::kernel_type, GaussianProcess::negative_marginal_log_likelihood(), GaussianProcess::numNuggetTerms, GaussianProcess::numPolyTerms, Surrogate::numQOI, Surrogate::numSamples, Surrogate::numVariables, GaussianProcess::objectiveFunctionHistory, GaussianProcess::objectiveGradientHistory, GaussianProcess::polyRegression, Surrogate::responseOffset, Surrogate::responseScaleFactor, DataScaler::scale_samples(), GaussianProcess::scaledBuildPoints, dakota::util::scaler_factory(), DataScaler::scaler_type(), GaussianProcess::setup_default_optimization_params(), GaussianProcess::setup_hyperparameter_bounds(), GaussianProcess::targetValues, GaussianProcess::thetaHistory, GaussianProcess::thetaValues, and GaussianProcess::verbosity.
Referenced by GaussianProcess::GaussianProcess().
Evaluate the scalar Gaussian Process at a set of prediction points.
[in] | eval_points | Matrix for prediction points - (num_points by num_features). |
Reimplemented from Surrogate.
References GaussianProcess::basisMatrix, GaussianProcess::betaValues, GaussianProcess::CholFact, GaussianProcess::compute_gram(), GaussianProcess::compute_pred_dists(), GaussianProcess::cwiseDists2, GaussianProcess::cwiseMixedDists2, Surrogate::dataScaler, GaussianProcess::estimateTrend, GaussianProcess::GramMatrix, GaussianProcess::hasBestCholFact, Surrogate::numVariables, GaussianProcess::polyRegression, GaussianProcess::predBasisMatrix, GaussianProcess::predMixedGramMatrix, Surrogate::responseOffset, Surrogate::responseScaleFactor, DataScaler::scale_samples(), and GaussianProcess::targetValues.
Evaluate the gradient of the scalar Gaussian process at a set of prediction points.
[in] | eval_points | Coordinates of the prediction points - (num_pts by num_features). |
Reimplemented from Surrogate.
References GaussianProcess::basisMatrix, GaussianProcess::betaValues, GaussianProcess::CholFact, GaussianProcess::compute_gram(), GaussianProcess::compute_pred_dists(), GaussianProcess::cwiseDists2, GaussianProcess::cwiseMixedDists, GaussianProcess::cwiseMixedDists2, Surrogate::dataScaler, GaussianProcess::estimateTrend, GaussianProcess::GramMatrix, GaussianProcess::hasBestCholFact, GaussianProcess::kernel, Surrogate::numVariables, GaussianProcess::polyRegression, GaussianProcess::predMixedGramMatrix, Surrogate::responseScaleFactor, DataScaler::scale_samples(), GaussianProcess::targetValues, and GaussianProcess::thetaValues.
Evaluate the Hessian of the scalar Gaussian process at a single point.
[in] | eval_point | Coordinates of the prediction point - (1 by num_features). |
Reimplemented from Surrogate.
References GaussianProcess::basisMatrix, GaussianProcess::betaValues, GaussianProcess::CholFact, GaussianProcess::compute_gram(), GaussianProcess::compute_pred_dists(), GaussianProcess::cwiseDists2, GaussianProcess::cwiseMixedDists, GaussianProcess::cwiseMixedDists2, Surrogate::dataScaler, GaussianProcess::estimateTrend, GaussianProcess::GramMatrix, GaussianProcess::hasBestCholFact, GaussianProcess::kernel, Surrogate::numVariables, GaussianProcess::polyRegression, GaussianProcess::predMixedGramMatrix, Surrogate::responseScaleFactor, DataScaler::scale_samples(), GaussianProcess::targetValues, and GaussianProcess::thetaValues.
Evaluate the covariance matrix for the Gaussian Process at a set of prediction points for a single QoI index.
[in] | eval_points | Matrix for the prediction points - (num_points by num_features). |
[in] | qoi | Index of response/QoI for which to compute derivatives |
References GaussianProcess::basisMatrix, GaussianProcess::betaValues, GaussianProcess::CholFact, GaussianProcess::compute_gram(), GaussianProcess::compute_pred_dists(), GaussianProcess::cwiseDists2, GaussianProcess::cwiseMixedDists2, GaussianProcess::cwisePredDists2, Surrogate::dataScaler, GaussianProcess::estimateTrend, GaussianProcess::GramMatrix, GaussianProcess::hasBestCholFact, Surrogate::numVariables, GaussianProcess::polyRegression, GaussianProcess::predBasisMatrix, GaussianProcess::predCovariance, GaussianProcess::predGramMatrix, GaussianProcess::predMixedGramMatrix, Surrogate::responseScaleFactor, DataScaler::scale_samples(), dakota::silence_unused_args(), and GaussianProcess::targetValues.
Referenced by GaussianProcess::covariance(), and GaussianProcess::variance().
Evaluate the covariance matrix for the Gaussian Process at a set of prediction points for QoI index 0.
[in] | eval_points | Matrix for the prediction points - (num_points by num_features). |
References GaussianProcess::covariance().
Evaluate the variance of the Gaussian Process at a set of prediction points for a given QoI index.
[in] | eval_points | Matrix for the prediction points - (num_points by num_features). |
[in] | qoi | Index of response/QoI for which to compute derivatives |
References GaussianProcess::covariance(), and dakota::silence_unused_args().
Referenced by PYBIND11_MODULE(), and GaussianProcess::variance().
Evaluate the variance of the Gaussian Process at a set of prediction points for QoI index 0.
[in] | eval_points | Matrix for the prediction points - (num_points by num_features). |
References GaussianProcess::variance().
void negative_marginal_log_likelihood | ( | bool | compute_grad, |
bool | compute_gram, | ||
double & | obj_value, | ||
VectorXd & | obj_gradient | ||
) |
Evaluate the negative marginal loglikelihood and its gradient.
[in] | compute_grad | Flag for computation of gradient. |
[in] | compute_gram | Flag for various Gram matrix calculations. |
[out] | obj_value | Value of the objection function. |
[out] | obj_gradient | Gradient of the objective function. |
References GaussianProcess::basisMatrix, GaussianProcess::betaValues, GaussianProcess::CholFact, GaussianProcess::compute_gram(), GaussianProcess::cwiseDists2, GaussianProcess::estimatedNuggetValue, GaussianProcess::estimateNugget, GaussianProcess::estimateTrend, GaussianProcess::eyeMatrix, GaussianProcess::GramMatrix, GaussianProcess::GramMatrixDerivs, GaussianProcess::GramResidualSolution, GaussianProcess::numPolyTerms, Surrogate::numSamples, Surrogate::numVariables, GaussianProcess::PI, GaussianProcess::targetValues, and GaussianProcess::trendTargetResidual.
Referenced by GaussianProcess::build(), GP_Objective::gradient(), and GP_Objective::value().
void setup_hyperparameter_bounds | ( | VectorXd & | sigma_bounds, |
MatrixXd & | length_scale_bounds, | ||
VectorXd & | nugget_bounds | ||
) |
Initialize the hyperparameter bounds for MLE from values in configOptions.
[out] | sigma_bounds | Bounds for the sigma (i.e. scale) hyperparameter. |
[out] | length_scale_bounds | Bounds for the anisotropic length-scale hyperparameters. |
[out] | nugget_bounds | Bounds for the estimated nugget hyperparameter. |
References Surrogate::configOptions, GaussianProcess::estimateNugget, GaussianProcess::numNuggetTerms, and Surrogate::numVariables.
Referenced by GaussianProcess::build().
int get_num_opt_variables | ( | ) |
Get the number of optimization variables.
References GaussianProcess::numNuggetTerms, GaussianProcess::numPolyTerms, and Surrogate::numVariables.
Referenced by GP_Objective::GP_Objective().
int get_num_variables | ( | ) | const |
Get the dimension of the feature space.
References Surrogate::numVariables.
|
inline |
Get the history of objective function values from MLE with restarts.
References GaussianProcess::objectiveFunctionHistory.
|
inline |
Get the history of objective function gradients from MLE with restarts.
References GaussianProcess::objectiveGradientHistory.
|
inline |
Get the history of hyperparameter values from MLE with restarts.
References GaussianProcess::thetaHistory.
Referenced by PYBIND11_MODULE().
void set_opt_params | ( | const std::vector< double > & | opt_params | ) |
Update the vector of optimization parameters.
[in] | opt_params | Vector of optimization parameter values. |
References GaussianProcess::betaValues, GaussianProcess::estimatedNuggetValue, GaussianProcess::estimateNugget, GaussianProcess::estimateTrend, GaussianProcess::numPolyTerms, Surrogate::numVariables, and GaussianProcess::thetaValues.
Referenced by GP_Objective::gradient(), and GP_Objective::value().
|
private |
Compute distances between build and prediction points. This includes build-prediction and prediction-prediction distance matrices.
[in] | scaled_pred_pts | Matrix of scaled prediction points. |
References GaussianProcess::cwiseMixedDists, GaussianProcess::cwiseMixedDists2, GaussianProcess::cwisePredDists2, Surrogate::numSamples, Surrogate::numVariables, and GaussianProcess::scaledBuildPoints.
Referenced by GaussianProcess::covariance(), GaussianProcess::gradient(), GaussianProcess::hessian(), and GaussianProcess::value().
|
private |
Compute a Gram matrix given a vector of squared distances and optionally compute its derivatives and/or adds nugget terms.
[in] | dists2 | Vector of squared distance matrices. |
[in] | add_nugget | Bool for whether or add nugget terms. |
[in] | compute_derivs | Bool for whether or not to compute the derivatives of the Gram matrix. |
[out] | gram | Gram matrix. |
References GaussianProcess::estimatedNuggetValue, GaussianProcess::estimateNugget, GaussianProcess::fixedNuggetValue, GaussianProcess::GramMatrixDerivs, GaussianProcess::kernel, and GaussianProcess::thetaValues.
Referenced by GaussianProcess::build(), GaussianProcess::covariance(), GaussianProcess::gradient(), GaussianProcess::hessian(), GaussianProcess::negative_marginal_log_likelihood(), and GaussianProcess::value().
|
private |
Randomly generate initial guesses for the optimization routine.
[in] | sigma_bounds | Bounds for the scaling hyperparameter (sigma). |
[in] | length_scale_bounds | Bounds for the length scales for each feature (l). |
[in] | nugget_bounds | Bounds for the nugget term. |
[in] | num_restarts | Number of restarts for the optimizer. |
[in] | seed | Seed for the random number generator. |
[out] | initial_guesses | Matrix of initial guesses. |
References dakota::util::create_uniform_random_double_matrix(), GaussianProcess::estimateNugget, GaussianProcess::estimateTrend, GaussianProcess::numNuggetTerms, GaussianProcess::numPolyTerms, and Surrogate::numVariables.
Referenced by GaussianProcess::build().
|
private |
Set the default optimization parameters for ROL for GP hyperparameter estimation.
[in] | rol_params | RCP to a Teuchos::ParameterList of ROL's options. |
Referenced by GaussianProcess::build().