Dakota  Version
Explore and Predict with Confidence
Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
GaussianProcess Class Reference

The GaussianProcess constructs a Gaussian Process regressor surrogate given a matrix of data. More...

Inheritance diagram for GaussianProcess:
Surrogate

Public Member Functions

 GaussianProcess ()
 Constructor that uses defaultConfigOptions and does not build.
 
 GaussianProcess (const ParameterList &param_list)
 Constructor that sets configOptions but does not build. More...
 
 GaussianProcess (const std::string &param_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 &param_list)
 Constructor for the GaussianProcess that sets configOptions and builds the GP. More...
 
 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. 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, const int qoi) override
 Evaluate the Gaussian Process at a set of prediction points for a single qoi. More...
 
VectorXd value (const MatrixXd &eval_points)
 Evaluate the Gaussian Process at a set of prediction points for QoI index 0. More...
 
MatrixXd gradient (const MatrixXd &eval_points, const int qoi) override
 Evaluate the gradient of the Gaussian process at a set of prediction points for a single QoI. More...
 
MatrixXd gradient (const MatrixXd &eval_points)
 Evaluate the gradient of the Gaussian process at a set of prediction points for QoI index 0. More...
 
MatrixXd hessian (const MatrixXd &eval_point, const int qoi) override
 Evaluate the Hessian of the Gaussian process at a single point for a single QoI. More...
 
MatrixXd hessian (const MatrixXd &eval_point)
 Evaluate the Hessian of the Gaussian process at a single point for QoI index 0. 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< Surrogateclone () 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 &param_list)
 Constructor that sets configOptions but does not build. More...
 
 Surrogate (const MatrixXd &samples, const MatrixXd &response, const ParameterList &param_list)
 Constructor for the Surrogate that sets configOptions and builds the surrogate (does nothing in the base class). More...
 
virtual ~Surrogate ()
 Default destructor.
 
VectorXd value (const MatrixXd &eval_points)
 Evaluate the Surrogate at a set of prediction points for QoI index 0. More...
 
MatrixXd gradient (const MatrixXd &eval_points)
 Evaluate the gradient of the Surrogate at a set of prediction points for QoI index 0. More...
 
MatrixXd hessian (const MatrixXd &eval_point)
 Evaluate the Hessian of the Surrogate at a single point for QoI index 0. 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.
 
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< MatrixXdGramMatrixDerivs
 Derivatives of the Gram matrix w.r.t. the hyperparameters.
 
std::vector< MatrixXdcwiseDists2
 Squared component-wise distances between points in the surrogate dataset.
 
std::vector< MatrixXdcwiseMixedDists
 Component-wise distances between prediction and build points.
 
std::vector< MatrixXdcwiseMixedDists2
 Squared component-wise distances between prediction and build points.
 
std::vector< MatrixXdcwisePredDists2
 Component-wise distances between prediction points.
 
Eigen::LDLT< MatrixXdCholFact
 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< PolynomialRegressionpolyRegression
 PolynomialRegression for trend function.
 
std::string kernel_type
 Kernel type.
 
std::shared_ptr< Kernelkernel
 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< Surrogateload (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.
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ GaussianProcess() [1/4]

GaussianProcess ( const ParameterList param_list)

Constructor that sets configOptions but does not build.

Parameters
[in]param_listList that overrides entries in defaultConfigOptions.

References Surrogate::configOptions, GaussianProcess::default_options(), and Surrogate::defaultConfigOptions.

◆ GaussianProcess() [2/4]

GaussianProcess ( const std::string &  param_list_yaml_filename)

Constructor for the GaussianProcess that sets configOptions but does not build the GP.

Parameters
[in]param_list_yaml_filenameA 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() [3/4]

GaussianProcess ( const MatrixXd samples,
const MatrixXd response,
const ParameterList param_list 
)

Constructor for the GaussianProcess that sets configOptions and builds the GP.

Parameters
[in]samplesMatrix of data for surrogate construction - (num_samples by num_features)
[in]responseVector of targets for surrogate construction - (num_samples by num_qoi = 1; only 1 response is supported currently).
[in]param_listList that overrides entries in defaultConfigOptions

References GaussianProcess::build(), Surrogate::configOptions, and GaussianProcess::default_options().

◆ GaussianProcess() [4/4]

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.

Parameters
[in]samplesMatrix of data for surrogate construction - (num_samples by num_features)
[in]responseVector of targets for surrogate construction - (num_samples by num_qoi = 1; only 1 response is supported currently).
[in]param_list_yaml_filenameA 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().

Member Function Documentation

◆ build()

void build ( const MatrixXd eval_points,
const MatrixXd response 
)
overridevirtual

Build the GP using specified build data.

Parameters
[in]eval_pointsMatrix of data for surrogate construction - (num_samples by num_features)
[in]responseVector 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().

◆ value() [1/2]

VectorXd value ( const MatrixXd eval_points,
const int  qoi 
)
overridevirtual

◆ value() [2/2]

VectorXd value ( const MatrixXd eval_points)
inline

Evaluate the Gaussian Process at a set of prediction points for QoI index 0.

Parameters
[in]eval_pointsMatrix for prediction points - (num_points by num_features).
Returns
Mean of the Gaussian process at the prediction points.

References Surrogate::value().

◆ gradient() [1/2]

MatrixXd gradient ( const MatrixXd eval_points,
const int  qoi 
)
overridevirtual

Evaluate the gradient of the Gaussian process at a set of prediction points for a single QoI.

Parameters
[in]eval_pointsCoordinates of the prediction points - (num_pts by num_features).
[in]qoiIndex of response/QoI for which to compute derivatives.
Returns
Matrix of gradient vectors at 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(), dakota::silence_unused_args(), GaussianProcess::targetValues, and GaussianProcess::thetaValues.

◆ gradient() [2/2]

MatrixXd gradient ( const MatrixXd eval_points)
inline

Evaluate the gradient of the Gaussian process at a set of prediction points for QoI index 0.

Parameters
[in]eval_pointsCoordinates of the prediction points - (num_pts by num_features).
Returns
Matrix of gradient vectors at the prediction points - (num_pts by num_features).

References Surrogate::gradient().

◆ hessian() [1/2]

MatrixXd hessian ( const MatrixXd eval_point,
const int  qoi 
)
overridevirtual

◆ hessian() [2/2]

MatrixXd hessian ( const MatrixXd eval_point)
inline

Evaluate the Hessian of the Gaussian process at a single point for QoI index 0.

Parameters
[in]eval_pointCoordinates of the prediction point - (1 by num_features).
Returns
Hessian matrix at the prediction point - (num_features by num_features).

References Surrogate::hessian().

◆ covariance() [1/2]

MatrixXd covariance ( const MatrixXd eval_points,
const int  qoi 
)

◆ covariance() [2/2]

MatrixXd covariance ( const MatrixXd eval_points)
inline

Evaluate the covariance matrix for the Gaussian Process at a set of prediction points for QoI index 0.

Parameters
[in]eval_pointsMatrix for the prediction points - (num_points by num_features).
Returns
[out] Covariance of the Gaussian process at the prediction points.

References GaussianProcess::covariance().

◆ variance() [1/2]

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.

Parameters
[in]eval_pointsMatrix for the prediction points - (num_points by num_features).
[in]qoiIndex of response/QoI for which to compute derivatives
Returns
[out] Variance of the Gaussian process at the prediction points.

References GaussianProcess::covariance(), and dakota::silence_unused_args().

Referenced by PYBIND11_MODULE(), and GaussianProcess::variance().

◆ variance() [2/2]

VectorXd variance ( const MatrixXd eval_points)
inline

Evaluate the variance of the Gaussian Process at a set of prediction points for QoI index 0.

Parameters
[in]eval_pointsMatrix for the prediction points - (num_points by num_features).
Returns
[out] Variance of the Gaussian process at the prediction points.

References GaussianProcess::variance().

◆ negative_marginal_log_likelihood()

void negative_marginal_log_likelihood ( bool  compute_grad,
bool  compute_gram,
double &  obj_value,
VectorXd obj_gradient 
)

◆ setup_hyperparameter_bounds()

void setup_hyperparameter_bounds ( VectorXd sigma_bounds,
MatrixXd length_scale_bounds,
VectorXd nugget_bounds 
)

Initialize the hyperparameter bounds for MLE from values in configOptions.

Parameters
[out]sigma_boundsBounds for the sigma (i.e. scale) hyperparameter.
[out]length_scale_boundsBounds for the anisotropic length-scale hyperparameters.
[out]nugget_boundsBounds for the estimated nugget hyperparameter.

References Surrogate::configOptions, GaussianProcess::estimateNugget, GaussianProcess::numNuggetTerms, and Surrogate::numVariables.

Referenced by GaussianProcess::build().

◆ get_num_opt_variables()

int get_num_opt_variables ( )

Get the number of optimization variables.

Returns
Number of total optimization variables (hyperparameters + trend coefficients + nugget)

References GaussianProcess::numNuggetTerms, GaussianProcess::numPolyTerms, and Surrogate::numVariables.

Referenced by GP_Objective::GP_Objective().

◆ get_num_variables()

int get_num_variables ( ) const

Get the dimension of the feature space.

Returns
numVariables The dimension of the feature space.

References Surrogate::numVariables.

◆ get_objective_function_history()

VectorXd get_objective_function_history ( )
inline

Get the history of objective function values from MLE with restarts.

Returns
objectiveFunctionHistory Vector of final objective function values.

References GaussianProcess::objectiveFunctionHistory.

◆ get_objective_gradient_history()

MatrixXd get_objective_gradient_history ( )
inline

Get the history of objective function gradients from MLE with restarts.

Returns
objectiveGradientHistory Matrix of final objective funciton values
  • (num_restarts, num_hyperparameters).

References GaussianProcess::objectiveGradientHistory.

◆ get_theta_history()

MatrixXd get_theta_history ( )
inline

Get the history of hyperparameter values from MLE with restarts.

Returns
thetaHistory Vector of final hyperparameter (theta) values
  • (num_restarts, num_hyperparameters).

References GaussianProcess::thetaHistory.

Referenced by PYBIND11_MODULE().

◆ set_opt_params()

void set_opt_params ( const std::vector< double > &  opt_params)

◆ compute_pred_dists()

void compute_pred_dists ( const MatrixXd scaled_pred_pts)
private

Compute distances between build and prediction points. This includes build-prediction and prediction-prediction distance matrices.

Parameters
[in]scaled_pred_ptsMatrix 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().

◆ compute_gram()

void compute_gram ( const std::vector< MatrixXd > &  dists2,
bool  add_nugget,
bool  compute_derivs,
MatrixXd gram 
)
private

Compute a Gram matrix given a vector of squared distances and optionally compute its derivatives and/or adds nugget terms.

Parameters
[in]dists2Vector of squared distance matrices.
[in]add_nuggetBool for whether or add nugget terms.
[in]compute_derivsBool for whether or not to compute the derivatives of the Gram matrix.
[out]gramGram 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().

◆ generate_initial_guesses()

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 
)
private

Randomly generate initial guesses for the optimization routine.

Parameters
[in]sigma_boundsBounds for the scaling hyperparameter (sigma).
[in]length_scale_boundsBounds for the length scales for each feature (l).
[in]nugget_boundsBounds for the nugget term.
[in]num_restartsNumber of restarts for the optimizer.
[in]seedSeed for the random number generator.
[out]initial_guessesMatrix 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().

◆ setup_default_optimization_params()

void setup_default_optimization_params ( Teuchos::RCP< ParameterList rol_params)
private

Set the default optimization parameters for ROL for GP hyperparameter estimation.

Parameters
[in]rol_paramsRCP to a Teuchos::ParameterList of ROL's options.

Referenced by GaussianProcess::build().


The documentation for this class was generated from the following files: