AutoPas  3.0.0
Loading...
Searching...
No Matches
GaussianProcess.h
Go to the documentation of this file.
1
7#pragma once
8
9#include <Eigen/Core>
10#include <utility>
11
12#include "AcquisitionFunction.h"
16#include "autopas/utils/Math.h"
19
20namespace autopas {
21
32
33 using Vector = Eigen::VectorXd;
34
35 // number of samples to find optimal hyperparameters
36 static constexpr size_t hp_sample_size = 10000;
37 // number of hyperparameters
38 static constexpr size_t hp_size = 100;
39
40 public:
47 GaussianProcess(size_t dims, double sigma, Random &rngRef);
48
49 virtual ~GaussianProcess();
50
55 void setDimension(size_t dims);
56
60 void clear();
61
66 [[nodiscard]] size_t numEvidence() const;
67
72 [[nodiscard]] std::pair<const std::vector<Vector> &, const Vector &> getEvidence() const;
73
81 void addEvidence(const Vector &input, double output, bool tuneHypers);
82
87 [[nodiscard]] const Vector &getEvidenceMax() const;
88
95 [[nodiscard]] double predictMean(const Vector &input) const;
96
101 [[nodiscard]] double getDefaultVar() const;
102
108 [[nodiscard]] double predictVar(const Vector &input) const;
109
116 [[nodiscard]] double predictOutputPDF(const Vector &input, double output) const;
117
125 [[nodiscard]] double predictOutputScaledPDF(const Vector &input, double output) const;
126
134 [[nodiscard]] double calcAcquisition(AcquisitionFunctionOption af, const Vector &input) const;
135
143 [[nodiscard]] Vector sampleAquisitionMax(AcquisitionFunctionOption af, const std::vector<Vector> &samples) const;
144
155 [[nodiscard]] static std::tuple<std::vector<double>, std::vector<double>,
156 std::vector<autopas::GaussianProcess::Vector>>
157 generateHyperparameterSamples(size_t sampleSize, Random &rng, size_t dims, double sigma, double evidenceMinValue,
158 double evidenceMaxValue);
159
164 [[nodiscard]] std::vector<GaussianHyperparameters> &getHyperparameters();
165
174 void setHyperparameters(const std::vector<double> &sample_means, const std::vector<double> &sample_thetas,
175 const std::vector<autopas::GaussianProcess::Vector> &sample_dimScales);
176
181
182 private:
189 void tuneHyperparameters();
190
200 [[nodiscard]] static double kernel(const Vector &input1, const Vector &input2, double theta,
201 const autopas::GaussianProcess::Vector &dimScale);
202
210 [[nodiscard]] autopas::GaussianProcess::Vector kernelVector(const Vector &input, double theta,
211 const autopas::GaussianProcess::Vector &dimScale) const;
212
213 std::vector<Vector> _inputs;
214 Vector _outputs;
215
219 size_t _dims;
220
224 double _evidenceMinValue;
228 double _evidenceMaxValue;
232 Vector _evidenceMaxVector;
233
237 const double _sigma;
238
242 std::vector<GaussianHyperparameters> _hypers;
243
244 Random &_rng;
245};
246} // namespace autopas
Hyperparameters of Gaussian models and derived matrices used for prediction.
Definition: GaussianHyperparameters.h:15
Gaussian process is a stochastical model.
Definition: GaussianProcess.h:30
double predictMean(const Vector &input) const
Try to predict f(x) using the evidence provided so far.
Definition: GaussianProcess.cpp:84
void setDimension(size_t dims)
Change input dimension.
Definition: GaussianProcess.cpp:25
double predictOutputScaledPDF(const Vector &input, double output) const
Calculate the scaled probability density of provided output given provided input.
Definition: GaussianProcess.cpp:139
const Vector & getEvidenceMax() const
Get the evidence with the highest output value.
Definition: GaussianProcess.cpp:76
Vector sampleAquisitionMax(AcquisitionFunctionOption af, const std::vector< Vector > &samples) const
Find the input in samples which maximizes given aquisition function.
Definition: GaussianProcess.cpp:151
double predictVar(const Vector &input) const
The variance of the predicted f(x) from predictMean().
Definition: GaussianProcess.cpp:113
std::vector< GaussianHyperparameters > & getHyperparameters()
Get current hyperparameters.
Definition: GaussianProcess.cpp:218
double getDefaultVar() const
Get the variance if evidence are ignored.
Definition: GaussianProcess.cpp:105
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.
Definition: GaussianProcess.cpp:220
void clear()
Discard all evidence.
Definition: GaussianProcess.cpp:30
double calcAcquisition(AcquisitionFunctionOption af, const Vector &input) const
Calculates the acquisition function for given input.
Definition: GaussianProcess.cpp:146
std::pair< const std::vector< Vector > &, const Vector & > getEvidence() const
Get all currently stored evidence.
Definition: GaussianProcess.cpp:39
size_t numEvidence() const
Get the number of evidence provided.
Definition: GaussianProcess.cpp:36
void addEvidence(const Vector &input, double output, bool tuneHypers)
Provide a input-output pair as evidence.
Definition: GaussianProcess.cpp:43
void normalizeHyperparameters()
Normalize weights of hyperparameters and truncate lowest weights.
Definition: GaussianProcess.cpp:240
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.
Definition: GaussianProcess.cpp:170
double predictOutputPDF(const Vector &input, double output) const
Calculate the probability density of provided output given provided input.
Definition: GaussianProcess.cpp:133
Class for random algorithms.
Definition: Random.h:22
This is the main namespace of AutoPas.
Definition: AutoPasDecl.h:32