AutoPas  3.0.0
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | List of all members
autopas::GaussianProcess Class Reference

Gaussian process is a stochastical model. More...

#include <GaussianProcess.h>

Public Member Functions

 GaussianProcess (size_t dims, double sigma, Random &rngRef)
 Constructor.
 
void setDimension (size_t dims)
 Change input dimension.
 
void clear ()
 Discard all evidence.
 
size_t numEvidence () const
 Get the number of evidence provided.
 
std::pair< const std::vector< Vector > &, const Vector & > getEvidence () const
 Get all currently stored evidence.
 
void addEvidence (const Vector &input, double output, bool tuneHypers)
 Provide a input-output pair as evidence.
 
const Vector & getEvidenceMax () const
 Get the evidence with the highest output value.
 
double predictMean (const Vector &input) const
 Try to predict f(x) using the evidence provided so far.
 
double getDefaultVar () const
 Get the variance if evidence are ignored.
 
double predictVar (const Vector &input) const
 The variance of the predicted f(x) from predictMean().
 
double predictOutputPDF (const Vector &input, double output) const
 Calculate the probability density of provided output given provided input.
 
double predictOutputScaledPDF (const Vector &input, double output) const
 Calculate the scaled probability density of provided output given provided input.
 
double calcAcquisition (AcquisitionFunctionOption af, const Vector &input) const
 Calculates the acquisition function for given input.
 
Vector sampleAquisitionMax (AcquisitionFunctionOption af, const std::vector< Vector > &samples) const
 Find the input in samples which maximizes given aquisition function.
 
std::vector< GaussianHyperparameters > & getHyperparameters ()
 Get current hyperparameters.
 
void setHyperparameters (const std::vector< double > &sample_means, const std::vector< double > &sample_thetas, const std::vector< autopas::GaussianProcess::Vector > &sample_dimScales)
 Set the hyperparameters: means, theta, dimScale.
 
void normalizeHyperparameters ()
 Normalize weights of hyperparameters and truncate lowest weights.
 

Static Public Member Functions

static std::tuple< std::vector< double >, std::vector< double >, std::vector< autopas::GaussianProcess::Vector > > generateHyperparameterSamples (size_t sampleSize, Random &rng, size_t dims, double sigma, double evidenceMinValue, double evidenceMaxValue)
 Generate hyperparameter samples.
 

Detailed Description

Gaussian process is a stochastical model.

It predicts the output of a blackbox function f(x) for given input x. To do so, some sample input-output pairs (x,f(x)) should be provided as evidence.

Currently the squared exponential kernel is used. TODO: maybe offer some options.

Constructor & Destructor Documentation

◆ GaussianProcess()

autopas::GaussianProcess::GaussianProcess ( size_t  dims,
double  sigma,
Random rngRef 
)

Constructor.

Parameters
dimsnumber of input dimensions
sigmafixed noise
rngRefreference to rng

Member Function Documentation

◆ addEvidence()

void autopas::GaussianProcess::addEvidence ( const Vector &  input,
double  output,
bool  tuneHypers 
)

Provide a input-output pair as evidence.

Each evidence improve the quality of future predictions.

Parameters
inputx
outputf(x)
tuneHypersif false hyperparemeters need to be set manually

◆ calcAcquisition()

double autopas::GaussianProcess::calcAcquisition ( AcquisitionFunctionOption  af,
const Vector &  input 
) const

Calculates the acquisition function for given input.

Parameters
afacquisition function a:input->double
inputi
Returns
a(i). This value can be compared with values a(x) of other inputs x to weigh which input would give the most gain if its evidence were provided.

◆ generateHyperparameterSamples()

std::tuple< std::vector< double >, std::vector< double >, std::vector< autopas::GaussianProcess::Vector > > autopas::GaussianProcess::generateHyperparameterSamples ( size_t  sampleSize,
autopas::Random rng,
size_t  dims,
double  sigma,
double  evidenceMinValue,
double  evidenceMaxValue 
)
static

Generate hyperparameter samples.

Parameters
sampleSizesize
rngrandom number generator
dimsnumber of input dimension
sigmafixed noise
evidenceMinValuecurrent lowest evidence output
evidenceMaxValuecurrent highest evidence output
Returns
tuple [means, thetas, dimScales]

◆ getDefaultVar()

double autopas::GaussianProcess::getDefaultVar ( ) const

Get the variance if evidence are ignored.

Returns

◆ getEvidence()

std::pair< const std::vector< autopas::GaussianProcess::Vector > &, const autopas::GaussianProcess::Vector & > autopas::GaussianProcess::getEvidence ( ) const

Get all currently stored evidence.

Returns
pair of inputs and outputs

◆ getEvidenceMax()

const autopas::GaussianProcess::Vector & autopas::GaussianProcess::getEvidenceMax ( ) const

Get the evidence with the highest output value.

Returns
input of max

◆ getHyperparameters()

std::vector< autopas::GaussianHyperparameters > & autopas::GaussianProcess::getHyperparameters ( )

Get current hyperparameters.

Returns

◆ numEvidence()

size_t autopas::GaussianProcess::numEvidence ( ) const

Get the number of evidence provided.

Returns

◆ predictMean()

double autopas::GaussianProcess::predictMean ( const Vector &  input) const

Try to predict f(x) using the evidence provided so far.

Parameters
inputx
Returns
expected output of f(x)

◆ predictOutputPDF()

double autopas::GaussianProcess::predictOutputPDF ( const Vector &  input,
double  output 
) const

Calculate the probability density of provided output given provided input.

Parameters
input
output
Returns

◆ predictOutputScaledPDF()

double autopas::GaussianProcess::predictOutputScaledPDF ( const Vector &  input,
double  output 
) const

Calculate the scaled probability density of provided output given provided input.

The probability density is scaled such that the maximum is 1.

Parameters
input
output
Returns

◆ predictVar()

double autopas::GaussianProcess::predictVar ( const Vector &  input) const

The variance of the predicted f(x) from predictMean().

Parameters
inputx
Returns
variance

◆ sampleAquisitionMax()

autopas::GaussianProcess::Vector autopas::GaussianProcess::sampleAquisitionMax ( AcquisitionFunctionOption  af,
const std::vector< Vector > &  samples 
) const

Find the input in samples which maximizes given aquisition function.

TODO: maybe add parameters for hyperparameters of aquisition functions

Parameters
affunction to maximize
samples
Returns

◆ setDimension()

void autopas::GaussianProcess::setDimension ( size_t  dims)

Change input dimension.

Current evidence will be discarded.

Parameters
dims

◆ setHyperparameters()

void autopas::GaussianProcess::setHyperparameters ( const std::vector< double > &  sample_means,
const std::vector< double > &  sample_thetas,
const std::vector< autopas::GaussianProcess::Vector > &  sample_dimScales 
)

Set the hyperparameters: means, theta, dimScale.

The samples are scored equal to the probability that given evidence and hyperparameter-sample generates given output. Hyperparameters weights should be normalized.

Parameters
sample_means
sample_thetas
sample_dimScales

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