$extrastylesheet
Dakota  Version 6.2
Public Member Functions | Protected Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes
GaussProcApproximation Class Reference

Derived approximation class for Gaussian Process implementation. More...

Inheritance diagram for GaussProcApproximation:
Approximation

List of all members.

Public Member Functions

 GaussProcApproximation ()
 default constructor
 GaussProcApproximation (const SharedApproxData &shared_data)
 alternate constructor
 GaussProcApproximation (const ProblemDescDB &problem_db, const SharedApproxData &shared_data)
 standard constructor
 ~GaussProcApproximation ()
 destructor

Protected Member Functions

int min_coefficients () const
 return the minimum number of samples (unknowns) required to build the derived class approximation type in numVars dimensions
int num_constraints () const
 return the number of constraints to be enforced via an anchor point
void build ()
 find the covariance parameters governing the Gaussian process response
Real value (const Variables &vars)
 retrieve the function value for a given parameter set
const RealVector & gradient (const Variables &vars)
 retrieve the function gradient at the predicted value for a given parameter set
Real prediction_variance (const Variables &vars)
 retrieve the variance of the predicted value for a given parameter set

Private Member Functions

void GPmodel_build ()
 Function to compute hyperparameters governing the GP.
void GPmodel_apply (const RealVector &new_x, bool variance_flag, bool gradients_flag)
 Function returns a response value using the GP surface.
void normalize_training_data ()
 Normalizes the initial inputs upon which the GP surface is based.
void get_trend ()
 Gets the trend (basis) functions for the calculation of the mean of the GP If the order = 0, the trend is a constant, if the order = 1, trend is linear, if order = 2, trend is quadratic.
void get_beta_coefficients ()
 Gets the beta coefficients for the calculation of the mean of the GP.
int get_cholesky_factor ()
 Gets the Cholesky factorization of the covariance matrix, with error checking.
void get_process_variance ()
 Gets the estimate of the process variance given the values of beta and the correlation lengthscales.
void get_cov_matrix ()
 calculates the covariance matrix for a given set of input points
void get_cov_vector ()
 calculates the covariance vector between a new point x and the set of inputs upon which the GP is based
void optimize_theta_global ()
 sets up and performs the optimization of the negative log likelihood to determine the optimal values of the covariance parameters using NCSUDirect
void optimize_theta_multipoint ()
 sets up and performs the optimization of the negative log likelihood to determine the optimal values of the covariance parameters using a gradient-based solver and multiple starting points
void predict (bool variance_flag, bool gradients_flag)
 Calculates the predicted new response value for x in normalized space.
Real calc_nll ()
 calculates the negative log likelihood function (based on covariance matrix)
void calc_grad_nll ()
 Gets the gradient of the negative log likelihood function with respect to the correlation lengthscales, theta.
void get_grad_cov_vector ()
 Calculate the derivatives of the covariance vector, with respect to each componeent of x.
void run_point_selection ()
 Runs the point selection algorithm, which will choose a subset of the training set with which to construct the GP model, and estimate the necessary parameters.
void initialize_point_selection ()
 Initializes the point selection routine by choosing a small initial subset of the training points.
void pointsel_get_errors (RealArray &delta)
 Uses the current GP model to compute predictions at all of the training points and find the errors.
int addpoint (int, IntArray &added_index)
 Adds a point to the effective training set. Returns 1 on success.
int pointsel_add_sel (const RealArray &delta)
 Accepts a vector of unsorted prediction errors, determines which points should be added to the effective training set, and adds them.
Real maxval (const RealArray &) const
 Return the maximum value of the elements in a vector.
void pointsel_write_points ()
 Writes out the training set before and after point selection.
void lhood_2d_grid_eval ()
 For problems with 2D input, evaluates the negative log likelihood on a grid.
void writex (const char[])
 Writes out the current training set (in original units) to a specified file.
void writeCovMat (char[])
 Writes out the covariance matrix to a specified file.

Static Private Member Functions

static void negloglik (int mode, int n, const Teuchos::SerialDenseVector< int, double > &X, Real &fx, Teuchos::SerialDenseVector< int, double > &grad_x, int &result_mode)
 static function used by OPT++ as the objective function to optimize the hyperparameters in the covariance of the GP by minimizing the negative log likelihood
static void constraint_eval (int mode, int n, const Teuchos::SerialDenseVector< int, double > &X, Teuchos::SerialDenseVector< int, double > &g, Teuchos::SerialDenseMatrix< int, double > &gradC, int &result_mode)
 static function used by OPT++ as the constraint function in the optimization of the negative log likelihood. Currently this function is empty: it is an unconstrained optimization.
static double negloglikNCSU (const RealVector &x)
 function used by NCSUOptimizer to optimize negloglik objective

Private Attributes

Real approxValue
 value of the approximation returned by value()
Real approxVariance
 value of the approximation returned by prediction_variance()
RealMatrix trainPoints
 A 2-D array (num sample sites = rows, num vars = columns) used to create the Gaussian process.
RealMatrix trainValues
 An array of response values; one response value per sample site.
RealVector trainMeans
 The mean of the input columns of trainPoints.
RealVector trainStdvs
 The standard deviation of the input columns of trainPoints.
RealMatrix normTrainPoints
 Current working set of normalized points upon which the GP is based.
RealMatrix trendFunction
 matrix to hold the trend function
RealMatrix betaCoeffs
 matrix to hold the beta coefficients for the trend function
RealSymMatrix covMatrix
 The covariance matrix where each element (i,j) is the covariance between points Xi and Xj in the initial set of samples.
RealMatrix covVector
 The covariance vector where each element (j,0) is the covariance between a new point X and point Xj from the initial set of samples.
RealMatrix approxPoint
 Point at which a prediction is requested. This is currently a single point, but it could be generalized to be a vector of points.
RealMatrix gradNegLogLikTheta
 matrix to hold the gradient of the negative log likelihood with respect to the theta correlation terms
Teuchos::SerialSpdDenseSolver
< int, Real > 
covSlvr
 The global solver for all computations involving the inverse of the covariance matrix.
RealMatrix gradCovVector
 A matrix, where each column is the derivative of the covVector with respect to a particular componenet of X.
RealMatrix normTrainPointsAll
 Set of all original samples available.
RealMatrix trainValuesAll
 All original samples available.
RealMatrix trendFunctionAll
 Trend function values corresponding to all original samples.
RealMatrix Rinv_YFb
 Matrix for storing inverse of correlation matrix Rinv*(Y-FB)
size_t numObs
 The number of observations on which the GP surface is built.
size_t numObsAll
 The original number of observations.
short trendOrder
 The number of variables in each X variable (number of dimensions of the problem).
RealVector thetaParams
 Theta is the vector of covariance parameters for the GP. We determine the values of theta by optimization Currently, the covariance function is theta[0]*exp(-0.5*sume)+delta*pow(sige,2). sume is the sum squared of weighted distances; it involves a sum of theta[1](Xi(1)-Xj(1))^2 + theta[2](Xi(2)-Xj(2))^2 + ... where Xi(1) is the first dimension value of multi-dimensional variable Xi. delta*pow(sige,2) is a jitter term used to improve matrix computations. delta is zero for the covariance between different points and 1 for the covariance between the same point. sige is the underlying process error.
Real procVar
 The process variance, the multiplier of the correlation matrix.
IntArray pointsAddedIndex
 Used by the point selection algorithm, this vector keeps track all points which have been added.
int cholFlag
 A global indicator for success of the Cholesky factorization.
bool usePointSelection
 a flag to indicate the use of point selection

Static Private Attributes

static GaussProcApproximationGPinstance
 pointer to the active object instance used within the static evaluator

Detailed Description

Derived approximation class for Gaussian Process implementation.

The GaussProcApproximation class provides a global approximation (surrogate) based on a Gaussian process. The Gaussian process is built after normalizing the function values, with zero mean. Opt++ is used to determine the optimal values of the covariance parameters, those which minimize the negative log likelihood function.


Constructor & Destructor Documentation

GaussProcApproximation ( ) [inline]

default constructor

alternate constructor used by EffGlobalOptimization and NonDGlobalReliability that does not use a problem database defaults here are no point selectinn and quadratic trend function.


Member Function Documentation

void GPmodel_apply ( const RealVector &  new_x,
bool  variance_flag,
bool  gradients_flag 
) [private]

Member Data Documentation

short trendOrder [private]

The number of variables in each X variable (number of dimensions of the problem).

The order of the basis function for the mean of the GP If the order = 0, the trend is a constant, if the order = 1, trend is linear, if order = 2, trend is quadratic.

Referenced by GaussProcApproximation::GaussProcApproximation(), GaussProcApproximation::get_beta_coefficients(), GaussProcApproximation::get_trend(), GaussProcApproximation::GPmodel_build(), and GaussProcApproximation::predict().


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