https://mooseframework.inl.gov
Classes | Public Member Functions | Protected Attributes | List of all members
StochasticTools::GaussianProcess Class Reference

Utility class dedicated to hold structures and functions commont to Gaussian Processes. More...

#include <GaussianProcess.h>

Classes

struct  GPOptimizerOptions
 Structure containing the optimization options for hyperparameter-tuning. More...
 

Public Member Functions

 GaussianProcess ()
 
void initialize (CovarianceFunctionBase *covariance_function, const std::vector< std::string > &params_to_tune, const std::vector< Real > &min=std::vector< Real >(), const std::vector< Real > &max=std::vector< Real >())
 Initializes the most important structures in the Gaussian Process: the covariance function and a tuning map which is used if the user requires parameter tuning. More...
 
void setupCovarianceMatrix (const RealEigenMatrix &training_params, const RealEigenMatrix &training_data, const GPOptimizerOptions &opts)
 Sets up the covariance matrix given data and optimization options. More...
 
void setupStoredMatrices (const RealEigenMatrix &input)
 Sets up the Cholesky decomposition and inverse action of the covariance matrix. More...
 
void linkCovarianceFunction (CovarianceFunctionBase *covariance_function)
 Finds and links the covariance function to this object. More...
 
void generateTuningMap (const std::vector< std::string > &params_to_tune, const std::vector< Real > &min=std::vector< Real >(), const std::vector< Real > &max=std::vector< Real >())
 Sets up the tuning map which is used if the user requires parameter tuning. More...
 
void standardizeParameters (RealEigenMatrix &parameters, bool keep_moments=false)
 Standardizes the vector of input parameters (x values). More...
 
void standardizeData (RealEigenMatrix &data, bool keep_moments=false)
 Standardizes the vector of responses (y values). More...
 
void tuneHyperParamsAdam (const RealEigenMatrix &training_params, const RealEigenMatrix &training_data, const GPOptimizerOptions &opts)
 
Real getLoss (RealEigenMatrix &inputs, RealEigenMatrix &outputs)
 
std::vector< RealgetGradient (RealEigenMatrix &inputs) const
 
void mapToVec (const std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real >> &tuning_data, const std::unordered_map< std::string, Real > &scalar_map, const std::unordered_map< std::string, std::vector< Real >> &vector_map, std::vector< Real > &vec) const
 Function used to convert the hyperparameter maps in this object to vectors. More...
 
void vecToMap (const std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real >> &tuning_data, std::unordered_map< std::string, Real > &scalar_map, std::unordered_map< std::string, std::vector< Real >> &vector_map, const std::vector< Real > &vec) const
 Function used to convert the vectors back to hyperparameter maps. More...
 
const StochasticTools::StandardizergetParamStandardizer () const
 Get constant reference to the contained structures. More...
 
const StochasticTools::StandardizergetDataStandardizer () const
 
const RealEigenMatrixgetK () const
 
const RealEigenMatrixgetKResultsSolve () const
 
const Eigen::LLT< RealEigenMatrix > & getKCholeskyDecomp () const
 
const CovarianceFunctionBasegetCovarFunction () const
 
const CovarianceFunctionBasegetCovarFunctionPtr () const
 
const std::string & getCovarType () const
 
const std::string & getCovarName () const
 
const std::vector< UserObjectName > & getDependentCovarNames () const
 
const std::map< UserObjectName, std::string > & getDependentCovarTypes () const
 
const unsigned intgetCovarNumOutputs () const
 
const unsigned intgetNumTunableParams () const
 
const std::unordered_map< std::string, Real > & getHyperParamMap () const
 
const std::unordered_map< std::string, std::vector< Real > > & getHyperParamVectorMap () const
 
StochasticTools::StandardizerparamStandardizer ()
 Get non-constant reference to the contained structures (if they need to be modified from the utside) More...
 
StochasticTools::StandardizerdataStandardizer ()
 
RealEigenMatrixK ()
 
RealEigenMatrixKResultsSolve ()
 
Eigen::LLT< RealEigenMatrix > & KCholeskyDecomp ()
 
CovarianceFunctionBasecovarFunctionPtr ()
 
CovarianceFunctionBasecovarFunction ()
 
std::string & covarType ()
 
std::string & covarName ()
 
std::map< UserObjectName, std::string > & dependentCovarTypes ()
 
std::vector< UserObjectName > & dependentCovarNames ()
 
unsigned intcovarNumOutputs ()
 
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > & tuningData ()
 
std::unordered_map< std::string, Real > & hyperparamMap ()
 
std::unordered_map< std::string, std::vector< Real > > & hyperparamVectorMap ()
 

Protected Attributes

CovarianceFunctionBase_covariance_function = nullptr
 Covariance function object. More...
 
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > _tuning_data
 Contains tuning inforation. Index of hyperparam, size, and min/max bounds. More...
 
unsigned int _num_tunable
 Number of tunable hyperparameters. More...
 
std::string _covar_type
 Type of covariance function used for this GP. More...
 
std::string _covar_name
 The name of the covariance function used in this GP. More...
 
std::vector< UserObjectName > _dependent_covar_names
 The names of the covariance functions the used covariance function depends on. More...
 
std::map< UserObjectName, std::string > _dependent_covar_types
 The types of the covariance functions the used covariance function depends on. More...
 
unsigned int _num_outputs
 The number of outputs of the GP. More...
 
std::unordered_map< std::string, Real_hyperparam_map
 Scalar hyperparameters. Stored for use in surrogate. More...
 
std::unordered_map< std::string, std::vector< Real > > _hyperparam_vec_map
 Vector hyperparameters. Stored for use in surrogate. More...
 
StochasticTools::Standardizer _param_standardizer
 Standardizer for use with params (x) More...
 
StochasticTools::Standardizer _data_standardizer
 Standardizer for use with data (y) More...
 
RealEigenMatrix _K
 An _n_sample by _n_sample covariance matrix constructed from the selected kernel function. More...
 
RealEigenMatrix _K_results_solve
 A solve of Ax=b via Cholesky. More...
 
Eigen::LLT< RealEigenMatrix_K_cho_decomp
 Cholesky decomposition Eigen object. More...
 
const RealEigenMatrix_training_params
 Paramaters (x) used for training, along with statistics. More...
 
const RealEigenMatrix_training_data
 Data (y) used for training. More...
 
unsigned int _batch_size
 The batch size for Adam optimization. More...
 

Detailed Description

Utility class dedicated to hold structures and functions commont to Gaussian Processes.

It can be used to standardize parameters, manipulate covariance data and compute additional stored matrices.

Definition at line 25 of file GaussianProcess.h.

Constructor & Destructor Documentation

◆ GaussianProcess()

StochasticTools::GaussianProcess::GaussianProcess ( )

Definition at line 46 of file GaussianProcess.C.

46 {}

Member Function Documentation

◆ covarFunction()

CovarianceFunctionBase& StochasticTools::GaussianProcess::covarFunction ( )
inline

Definition at line 209 of file GaussianProcess.h.

209 { return *_covariance_function; }
CovarianceFunctionBase * _covariance_function
Covariance function object.

◆ covarFunctionPtr()

CovarianceFunctionBase* StochasticTools::GaussianProcess::covarFunctionPtr ( )
inline

Definition at line 208 of file GaussianProcess.h.

208 { return _covariance_function; }
CovarianceFunctionBase * _covariance_function
Covariance function object.

◆ covarName()

std::string& StochasticTools::GaussianProcess::covarName ( )
inline

Definition at line 211 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

211 { return _covar_name; }
std::string _covar_name
The name of the covariance function used in this GP.

◆ covarNumOutputs()

unsigned int& StochasticTools::GaussianProcess::covarNumOutputs ( )
inline

Definition at line 214 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

214 { return _num_outputs; }
unsigned int _num_outputs
The number of outputs of the GP.

◆ covarType()

std::string& StochasticTools::GaussianProcess::covarType ( )
inline

Definition at line 210 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

210 { return _covar_type; }
std::string _covar_type
Type of covariance function used for this GP.

◆ dataStandardizer()

StochasticTools::Standardizer& StochasticTools::GaussianProcess::dataStandardizer ( )
inline

Definition at line 204 of file GaussianProcess.h.

Referenced by dataLoad(), dataStore(), GaussianProcessTrainer::postTrain(), and ActiveLearningGaussianProcess::reTrain().

204 { return _data_standardizer; }
StochasticTools::Standardizer _data_standardizer
Standardizer for use with data (y)

◆ dependentCovarNames()

std::vector<UserObjectName>& StochasticTools::GaussianProcess::dependentCovarNames ( )
inline

Definition at line 213 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

213 { return _dependent_covar_names; }
std::vector< UserObjectName > _dependent_covar_names
The names of the covariance functions the used covariance function depends on.

◆ dependentCovarTypes()

std::map<UserObjectName, std::string>& StochasticTools::GaussianProcess::dependentCovarTypes ( )
inline

Definition at line 212 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

212 { return _dependent_covar_types; }
std::map< UserObjectName, std::string > _dependent_covar_types
The types of the covariance functions the used covariance function depends on.

◆ generateTuningMap()

void StochasticTools::GaussianProcess::generateTuningMap ( const std::vector< std::string > &  params_to_tune,
const std::vector< Real > &  min = std::vector<Real>(),
const std::vector< Real > &  max = std::vector<Real>() 
)

Sets up the tuning map which is used if the user requires parameter tuning.

Parameters
params_to_tuneList of parameters which need to be tuned.
minList of lower bounds for the parameter tuning.
maxList of upper bounds for parameter tuning.

Definition at line 102 of file GaussianProcess.C.

Referenced by initialize().

105 {
106  _num_tunable = 0;
107 
108  const bool upper_bounds_specified = min_vector.size();
109  const bool lower_bounds_specified = max_vector.size();
110 
111  for (const auto param_i : index_range(params_to_tune))
112  {
113  const auto & hp = params_to_tune[param_i];
115  {
116  unsigned int size;
117  Real min;
118  Real max;
119  // Get size and default min/max
120  const bool found = _covariance_function->getTuningData(hp, size, min, max);
121 
122  if (!found)
123  ::mooseError("The covariance parameter ", hp, " could not be found!");
124 
125  // Check for overridden min/max
126  min = lower_bounds_specified ? min_vector[param_i] : min;
127  max = upper_bounds_specified ? max_vector[param_i] : max;
128  // Save data in tuple
129  _tuning_data[hp] = std::make_tuple(_num_tunable, size, min, max);
130  _num_tunable += size;
131  }
132  }
133 }
unsigned int _num_tunable
Number of tunable hyperparameters.
CovarianceFunctionBase * _covariance_function
Covariance function object.
auto max(const L &left, const R &right)
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > _tuning_data
Contains tuning inforation. Index of hyperparam, size, and min/max bounds.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
auto min(const L &left, const R &right)
virtual bool getTuningData(const std::string &name, unsigned int &size, Real &min, Real &max) const
Get the default minimum and maximum and size of a hyperparameter.
auto index_range(const T &sizable)
virtual bool isTunable(const std::string &name) const
Check if a given parameter is tunable.

◆ getCovarFunction()

const CovarianceFunctionBase& StochasticTools::GaussianProcess::getCovarFunction ( ) const
inline

Definition at line 177 of file GaussianProcess.h.

Referenced by GaussianProcessSurrogate::evaluate(), and GaussianProcessTrainer::GaussianProcessTrainer().

177 { return *_covariance_function; }
CovarianceFunctionBase * _covariance_function
Covariance function object.

◆ getCovarFunctionPtr()

const CovarianceFunctionBase* StochasticTools::GaussianProcess::getCovarFunctionPtr ( ) const
inline

Definition at line 178 of file GaussianProcess.h.

Referenced by GaussianProcessSurrogate::setupCovariance().

178 { return _covariance_function; }
CovarianceFunctionBase * _covariance_function
Covariance function object.

◆ getCovarName()

const std::string& StochasticTools::GaussianProcess::getCovarName ( ) const
inline

Definition at line 180 of file GaussianProcess.h.

180 { return _covar_name; }
std::string _covar_name
The name of the covariance function used in this GP.

◆ getCovarNumOutputs()

const unsigned int& StochasticTools::GaussianProcess::getCovarNumOutputs ( ) const
inline

Definition at line 189 of file GaussianProcess.h.

189 { return _num_outputs; }
unsigned int _num_outputs
The number of outputs of the GP.

◆ getCovarType()

const std::string& StochasticTools::GaussianProcess::getCovarType ( ) const
inline

Definition at line 179 of file GaussianProcess.h.

179 { return _covar_type; }
std::string _covar_type
Type of covariance function used for this GP.

◆ getDataStandardizer()

const StochasticTools::Standardizer& StochasticTools::GaussianProcess::getDataStandardizer ( ) const
inline

Definition at line 173 of file GaussianProcess.h.

Referenced by GaussianProcessSurrogate::evaluate().

173 { return _data_standardizer; }
StochasticTools::Standardizer _data_standardizer
Standardizer for use with data (y)

◆ getDependentCovarNames()

const std::vector<UserObjectName>& StochasticTools::GaussianProcess::getDependentCovarNames ( ) const
inline

Definition at line 181 of file GaussianProcess.h.

182  {
183  return _dependent_covar_names;
184  }
std::vector< UserObjectName > _dependent_covar_names
The names of the covariance functions the used covariance function depends on.

◆ getDependentCovarTypes()

const std::map<UserObjectName, std::string>& StochasticTools::GaussianProcess::getDependentCovarTypes ( ) const
inline

Definition at line 185 of file GaussianProcess.h.

186  {
187  return _dependent_covar_types;
188  }
std::map< UserObjectName, std::string > _dependent_covar_types
The types of the covariance functions the used covariance function depends on.

◆ getGradient()

std::vector< Real > StochasticTools::GaussianProcess::getGradient ( RealEigenMatrix inputs) const

Definition at line 251 of file GaussianProcess.C.

Referenced by tuneHyperParamsAdam().

252 {
255  std::vector<Real> grad_vec;
256  grad_vec.resize(_num_tunable);
257  for (auto iter = _tuning_data.begin(); iter != _tuning_data.end(); ++iter)
258  {
259  std::string hyper_param_name = iter->first;
260  const auto first_index = std::get<0>(iter->second);
261  const auto num_entries = std::get<1>(iter->second);
262  for (unsigned int ii = 0; ii < num_entries; ++ii)
263  {
264  const auto global_index = first_index + ii;
265  _covariance_function->computedKdhyper(dKdhp, inputs, hyper_param_name, ii);
266  RealEigenMatrix tmp = alpha * dKdhp - _K_cho_decomp.solve(dKdhp);
267  grad_vec[global_index] = -tmp.trace() / 2.0;
268  }
269  }
270  return grad_vec;
271 }
unsigned int _num_tunable
Number of tunable hyperparameters.
RealEigenMatrix _K_results_solve
A solve of Ax=b via Cholesky.
CovarianceFunctionBase * _covariance_function
Covariance function object.
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > _tuning_data
Contains tuning inforation. Index of hyperparam, size, and min/max bounds.
Eigen::LLT< RealEigenMatrix > _K_cho_decomp
Cholesky decomposition Eigen object.
unsigned int _batch_size
The batch size for Adam optimization.
static const std::string alpha
Definition: NS.h:134
virtual bool computedKdhyper(RealEigenMatrix &dKdhp, const RealEigenMatrix &x, const std::string &hyper_param_name, unsigned int ind) const
Redirect dK/dhp for hyperparameter "hp".

◆ getHyperParamMap()

const std::unordered_map<std::string, Real>& StochasticTools::GaussianProcess::getHyperParamMap ( ) const
inline

Definition at line 191 of file GaussianProcess.h.

Referenced by GaussianProcessData::initialize().

191 { return _hyperparam_map; }
std::unordered_map< std::string, Real > _hyperparam_map
Scalar hyperparameters. Stored for use in surrogate.

◆ getHyperParamVectorMap()

const std::unordered_map<std::string, std::vector<Real> >& StochasticTools::GaussianProcess::getHyperParamVectorMap ( ) const
inline

Definition at line 192 of file GaussianProcess.h.

Referenced by GaussianProcessData::initialize().

193  {
194  return _hyperparam_vec_map;
195  }
std::unordered_map< std::string, std::vector< Real > > _hyperparam_vec_map
Vector hyperparameters. Stored for use in surrogate.

◆ getK()

const RealEigenMatrix& StochasticTools::GaussianProcess::getK ( ) const
inline

Definition at line 174 of file GaussianProcess.h.

174 { return _K; }
RealEigenMatrix _K
An _n_sample by _n_sample covariance matrix constructed from the selected kernel function.

◆ getKCholeskyDecomp()

const Eigen::LLT<RealEigenMatrix>& StochasticTools::GaussianProcess::getKCholeskyDecomp ( ) const
inline

Definition at line 176 of file GaussianProcess.h.

Referenced by GaussianProcessSurrogate::evaluate().

176 { return _K_cho_decomp; }
Eigen::LLT< RealEigenMatrix > _K_cho_decomp
Cholesky decomposition Eigen object.

◆ getKResultsSolve()

const RealEigenMatrix& StochasticTools::GaussianProcess::getKResultsSolve ( ) const
inline

Definition at line 175 of file GaussianProcess.h.

Referenced by GaussianProcessSurrogate::evaluate().

175 { return _K_results_solve; }
RealEigenMatrix _K_results_solve
A solve of Ax=b via Cholesky.

◆ getLoss()

Real StochasticTools::GaussianProcess::getLoss ( RealEigenMatrix inputs,
RealEigenMatrix outputs 
)

Definition at line 234 of file GaussianProcess.C.

Referenced by tuneHyperParamsAdam().

235 {
236  _covariance_function->computeCovarianceMatrix(_K, inputs, inputs, true);
237 
238  RealEigenMatrix flattened_data = outputs.reshaped(outputs.rows() * outputs.cols(), 1);
239 
240  setupStoredMatrices(flattened_data);
241 
242  Real log_likelihood = 0;
243  log_likelihood += -(flattened_data.transpose() * _K_results_solve)(0, 0);
244  log_likelihood += -std::log(_K.determinant());
245  log_likelihood -= _batch_size * std::log(2 * M_PI);
246  log_likelihood = -log_likelihood / 2;
247  return log_likelihood;
248 }
RealEigenMatrix _K_results_solve
A solve of Ax=b via Cholesky.
CovarianceFunctionBase * _covariance_function
Covariance function object.
RealEigenMatrix _K
An _n_sample by _n_sample covariance matrix constructed from the selected kernel function.
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
unsigned int _batch_size
The batch size for Adam optimization.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void computeCovarianceMatrix(RealEigenMatrix &K, const RealEigenMatrix &x, const RealEigenMatrix &xp, const bool is_self_covariance) const =0
Generates the Covariance Matrix given two sets of points in the parameter space.
void setupStoredMatrices(const RealEigenMatrix &input)
Sets up the Cholesky decomposition and inverse action of the covariance matrix.

◆ getNumTunableParams()

const unsigned int& StochasticTools::GaussianProcess::getNumTunableParams ( ) const
inline

Definition at line 190 of file GaussianProcess.h.

190 { return _num_tunable; }
unsigned int _num_tunable
Number of tunable hyperparameters.

◆ getParamStandardizer()

const StochasticTools::Standardizer& StochasticTools::GaussianProcess::getParamStandardizer ( ) const
inline

Get constant reference to the contained structures.

Definition at line 172 of file GaussianProcess.h.

Referenced by GaussianProcessSurrogate::evaluate().

172 { return _param_standardizer; }
StochasticTools::Standardizer _param_standardizer
Standardizer for use with params (x)

◆ hyperparamMap()

std::unordered_map<std::string, Real>& StochasticTools::GaussianProcess::hyperparamMap ( )
inline

Definition at line 219 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

219 { return _hyperparam_map; }
std::unordered_map< std::string, Real > _hyperparam_map
Scalar hyperparameters. Stored for use in surrogate.

◆ hyperparamVectorMap()

std::unordered_map<std::string, std::vector<Real> >& StochasticTools::GaussianProcess::hyperparamVectorMap ( )
inline

Definition at line 220 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

221  {
222  return _hyperparam_vec_map;
223  }
std::unordered_map< std::string, std::vector< Real > > _hyperparam_vec_map
Vector hyperparameters. Stored for use in surrogate.

◆ initialize()

void StochasticTools::GaussianProcess::initialize ( CovarianceFunctionBase covariance_function,
const std::vector< std::string > &  params_to_tune,
const std::vector< Real > &  min = std::vector<Real>(),
const std::vector< Real > &  max = std::vector<Real>() 
)

Initializes the most important structures in the Gaussian Process: the covariance function and a tuning map which is used if the user requires parameter tuning.

Parameters
covariance_functionPointer to the covariance function that needs to be used for the Gaussian Process.
params_to_tuneList of parameters which need to be tuned.
minList of lower bounds for the parameter tuning.
maxList of upper bounds for parameter tuning.

Definition at line 49 of file GaussianProcess.C.

Referenced by ActiveLearningGaussianProcess::ActiveLearningGaussianProcess(), and GaussianProcessTrainer::GaussianProcessTrainer().

53 {
54  linkCovarianceFunction(covariance_function);
55  generateTuningMap(params_to_tune, min, max);
56 }
void linkCovarianceFunction(CovarianceFunctionBase *covariance_function)
Finds and links the covariance function to this object.
void generateTuningMap(const std::vector< std::string > &params_to_tune, const std::vector< Real > &min=std::vector< Real >(), const std::vector< Real > &max=std::vector< Real >())
Sets up the tuning map which is used if the user requires parameter tuning.

◆ K()

RealEigenMatrix& StochasticTools::GaussianProcess::K ( )
inline

Definition at line 205 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

205 { return _K; }
RealEigenMatrix _K
An _n_sample by _n_sample covariance matrix constructed from the selected kernel function.

◆ KCholeskyDecomp()

Eigen::LLT<RealEigenMatrix>& StochasticTools::GaussianProcess::KCholeskyDecomp ( )
inline

Definition at line 207 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

207 { return _K_cho_decomp; }
Eigen::LLT< RealEigenMatrix > _K_cho_decomp
Cholesky decomposition Eigen object.

◆ KResultsSolve()

RealEigenMatrix& StochasticTools::GaussianProcess::KResultsSolve ( )
inline

Definition at line 206 of file GaussianProcess.h.

Referenced by dataLoad(), and dataStore().

206 { return _K_results_solve; }
RealEigenMatrix _K_results_solve
A solve of Ax=b via Cholesky.

◆ linkCovarianceFunction()

void StochasticTools::GaussianProcess::linkCovarianceFunction ( CovarianceFunctionBase covariance_function)

Finds and links the covariance function to this object.

Used mainly in the covariance data action.

Parameters
covariance_functionPointer to the covariance function that needs to be used for the Gaussian Process.

Definition at line 59 of file GaussianProcess.C.

Referenced by initialize(), and GaussianProcessSurrogate::setupCovariance().

60 {
61  _covariance_function = covariance_function;
67 }
CovarianceFunctionBase * _covariance_function
Covariance function object.
std::string _covar_name
The name of the covariance function used in this GP.
std::map< UserObjectName, std::string > _dependent_covar_types
The types of the covariance functions the used covariance function depends on.
virtual const std::string & name() const
unsigned int _num_outputs
The number of outputs of the GP.
std::vector< UserObjectName > _dependent_covar_names
The names of the covariance functions the used covariance function depends on.
const std::string & type() const
const std::vector< UserObjectName > & dependentCovarianceNames() const
Get the names of the dependent covariances.
void dependentCovarianceTypes(std::map< UserObjectName, std::string > &name_type_map) const
Populate a map with the names and types of the dependent covariance functions.
std::string _covar_type
Type of covariance function used for this GP.
unsigned int numOutputs() const
Return the number of outputs assumed for this covariance function.

◆ mapToVec()

void StochasticTools::GaussianProcess::mapToVec ( const std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real >> &  tuning_data,
const std::unordered_map< std::string, Real > &  scalar_map,
const std::unordered_map< std::string, std::vector< Real >> &  vector_map,
std::vector< Real > &  vec 
) const

Function used to convert the hyperparameter maps in this object to vectors.

Definition at line 274 of file GaussianProcess.C.

Referenced by tuneHyperParamsAdam().

280 {
281  for (auto iter : tuning_data)
282  {
283  const std::string & param_name = iter.first;
284  const auto scalar_it = scalar_map.find(param_name);
285  if (scalar_it != scalar_map.end())
286  vec[std::get<0>(iter.second)] = scalar_it->second;
287  else
288  {
289  const auto vector_it = vector_map.find(param_name);
290  if (vector_it != vector_map.end())
291  for (unsigned int ii = 0; ii < std::get<1>(iter.second); ++ii)
292  vec[std::get<0>(iter.second) + ii] = (vector_it->second)[ii];
293  }
294  }
295 }

◆ paramStandardizer()

StochasticTools::Standardizer& StochasticTools::GaussianProcess::paramStandardizer ( )
inline

Get non-constant reference to the contained structures (if they need to be modified from the utside)

Definition at line 203 of file GaussianProcess.h.

Referenced by dataLoad(), dataStore(), GaussianProcessTrainer::postTrain(), and ActiveLearningGaussianProcess::reTrain().

203 { return _param_standardizer; }
StochasticTools::Standardizer _param_standardizer
Standardizer for use with params (x)

◆ setupCovarianceMatrix()

void StochasticTools::GaussianProcess::setupCovarianceMatrix ( const RealEigenMatrix training_params,
const RealEigenMatrix training_data,
const GPOptimizerOptions opts 
)

Sets up the covariance matrix given data and optimization options.

Parameters
training_paramsThe training parameter values (x values) for the covariance matrix.
training_dataThe training data (y values) for the inversion of the covariance matrix.
optsThe optimizer options.

Definition at line 70 of file GaussianProcess.C.

Referenced by GaussianProcessTrainer::postTrain(), and ActiveLearningGaussianProcess::reTrain().

73 {
74  const bool batch_decision = opts.batch_size > 0 && (opts.batch_size <= training_params.rows());
75  _batch_size = batch_decision ? opts.batch_size : training_params.rows();
77 
78  if (_tuning_data.size())
79  tuneHyperParamsAdam(training_params, training_data, opts);
80 
81  _K.resize(training_params.rows() * training_data.cols(),
82  training_params.rows() * training_data.cols());
83  _covariance_function->computeCovarianceMatrix(_K, training_params, training_params, true);
84 
85  RealEigenMatrix flattened_data =
86  training_data.reshaped(training_params.rows() * training_data.cols(), 1);
87 
88  // Compute the Cholesky decomposition and inverse action of the covariance matrix
89  setupStoredMatrices(flattened_data);
90 
92 }
CovarianceFunctionBase * _covariance_function
Covariance function object.
std::unordered_map< std::string, std::vector< Real > > _hyperparam_vec_map
Vector hyperparameters. Stored for use in surrogate.
std::unordered_map< std::string, Real > _hyperparam_map
Scalar hyperparameters. Stored for use in surrogate.
void tuneHyperParamsAdam(const RealEigenMatrix &training_params, const RealEigenMatrix &training_data, const GPOptimizerOptions &opts)
unsigned int _num_outputs
The number of outputs of the GP.
void buildHyperParamMap(std::unordered_map< std::string, Real > &map, std::unordered_map< std::string, std::vector< Real >> &vec_map) const
Populates the input maps with the owned hyperparameters.
RealEigenMatrix _K
An _n_sample by _n_sample covariance matrix constructed from the selected kernel function.
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > _tuning_data
Contains tuning inforation. Index of hyperparam, size, and min/max bounds.
unsigned int _batch_size
The batch size for Adam optimization.
virtual void computeCovarianceMatrix(RealEigenMatrix &K, const RealEigenMatrix &x, const RealEigenMatrix &xp, const bool is_self_covariance) const =0
Generates the Covariance Matrix given two sets of points in the parameter space.
void setupStoredMatrices(const RealEigenMatrix &input)
Sets up the Cholesky decomposition and inverse action of the covariance matrix.

◆ setupStoredMatrices()

void StochasticTools::GaussianProcess::setupStoredMatrices ( const RealEigenMatrix input)

Sets up the Cholesky decomposition and inverse action of the covariance matrix.

Parameters
inputThe vector/matrix which right multiples the inverse of the covariance matrix.

Definition at line 95 of file GaussianProcess.C.

Referenced by getLoss(), and setupCovarianceMatrix().

96 {
97  _K_cho_decomp = _K.llt();
98  _K_results_solve = _K_cho_decomp.solve(input);
99 }
RealEigenMatrix _K_results_solve
A solve of Ax=b via Cholesky.
RealEigenMatrix _K
An _n_sample by _n_sample covariance matrix constructed from the selected kernel function.
Eigen::LLT< RealEigenMatrix > _K_cho_decomp
Cholesky decomposition Eigen object.

◆ standardizeData()

void StochasticTools::GaussianProcess::standardizeData ( RealEigenMatrix data,
bool  keep_moments = false 
)

Standardizes the vector of responses (y values).

Parameters
dataThe vector/matrix of input data.
keep_momentsIf previously computed or new moments are to be used.

Definition at line 144 of file GaussianProcess.C.

Referenced by GaussianProcessTrainer::postTrain(), and ActiveLearningGaussianProcess::reTrain().

145 {
146  if (!keep_moments)
149 }
StochasticTools::Standardizer _data_standardizer
Standardizer for use with data (y)
void getStandardized(RealEigenMatrix &input) const
Returns the standardized (centered and scaled) of the provided input.
Definition: Standardizer.C:80
void computeSet(const RealEigenMatrix &input)
Methods for computing and setting mean and standard deviation.
Definition: Standardizer.C:58

◆ standardizeParameters()

void StochasticTools::GaussianProcess::standardizeParameters ( RealEigenMatrix parameters,
bool  keep_moments = false 
)

Standardizes the vector of input parameters (x values).

Parameters
parametersThe vector/matrix of input data.
keep_momentsIf previously computed or new moments are to be used.

Definition at line 136 of file GaussianProcess.C.

Referenced by GaussianProcessTrainer::postTrain(), and ActiveLearningGaussianProcess::reTrain().

137 {
138  if (!keep_moments)
141 }
void getStandardized(RealEigenMatrix &input) const
Returns the standardized (centered and scaled) of the provided input.
Definition: Standardizer.C:80
StochasticTools::Standardizer _param_standardizer
Standardizer for use with params (x)
void computeSet(const RealEigenMatrix &input)
Methods for computing and setting mean and standard deviation.
Definition: Standardizer.C:58

◆ tuneHyperParamsAdam()

void StochasticTools::GaussianProcess::tuneHyperParamsAdam ( const RealEigenMatrix training_params,
const RealEigenMatrix training_data,
const GPOptimizerOptions opts 
)

Definition at line 152 of file GaussianProcess.C.

Referenced by setupCovarianceMatrix().

155 {
156  std::vector<Real> theta(_num_tunable, 0.0);
158 
160 
161  // Internal params for Adam; set to the recommended values in the paper
162  Real b1 = opts.b1;
163  Real b2 = opts.b2;
164  Real eps = opts.eps;
165 
166  std::vector<Real> m0(_num_tunable, 0.0);
167  std::vector<Real> v0(_num_tunable, 0.0);
168 
169  Real new_val;
170  Real m_hat;
171  Real v_hat;
172  Real store_loss = 0.0;
173  std::vector<Real> grad1;
174 
175  // Initialize randomizer
176  std::vector<unsigned int> v_sequence(training_params.rows());
177  std::iota(std::begin(v_sequence), std::end(v_sequence), 0);
178  RealEigenMatrix inputs(_batch_size, training_params.cols());
179  RealEigenMatrix outputs(_batch_size, training_data.cols());
180  if (opts.show_every_nth_iteration)
181  Moose::out << "OPTIMIZING GP HYPER-PARAMETERS USING Adam" << std::endl;
182  for (unsigned int ss = 0; ss < opts.num_iter; ++ss)
183  {
184  // Shuffle data
185  MooseRandom generator;
186  generator.seed(0, 1980);
187  generator.saveState();
188  MooseUtils::shuffle<unsigned int>(v_sequence, generator, 0);
189  for (unsigned int ii = 0; ii < _batch_size; ++ii)
190  {
191  for (unsigned int jj = 0; jj < training_params.cols(); ++jj)
192  inputs(ii, jj) = training_params(v_sequence[ii], jj);
193 
194  for (unsigned int jj = 0; jj < training_data.cols(); ++jj)
195  outputs(ii, jj) = training_data(v_sequence[ii], jj);
196  }
197 
198  store_loss = getLoss(inputs, outputs);
199  if (opts.show_every_nth_iteration && ((ss + 1) % opts.show_every_nth_iteration == 0))
200  Moose::out << "Iteration: " << ss + 1 << " LOSS: " << store_loss << std::endl;
201  grad1 = getGradient(inputs);
202  for (auto iter = _tuning_data.begin(); iter != _tuning_data.end(); ++iter)
203  {
204  const auto first_index = std::get<0>(iter->second);
205  const auto num_entries = std::get<1>(iter->second);
206  for (unsigned int ii = 0; ii < num_entries; ++ii)
207  {
208  const auto global_index = first_index + ii;
209  m0[global_index] = b1 * m0[global_index] + (1 - b1) * grad1[global_index];
210  v0[global_index] =
211  b2 * v0[global_index] + (1 - b2) * grad1[global_index] * grad1[global_index];
212  m_hat = m0[global_index] / (1 - std::pow(b1, (ss + 1)));
213  v_hat = v0[global_index] / (1 - std::pow(b2, (ss + 1)));
214  new_val = theta[global_index] - opts.learning_rate * m_hat / (std::sqrt(v_hat) + eps);
215 
216  const auto min_value = std::get<2>(iter->second);
217  const auto max_value = std::get<3>(iter->second);
218 
219  theta[global_index] = std::min(std::max(new_val, min_value), max_value);
220  }
221  }
224  }
225  if (opts.show_every_nth_iteration)
226  {
227  Moose::out << "OPTIMIZED GP HYPER-PARAMETERS:" << std::endl;
228  Moose::out << Moose::stringify(theta) << std::endl;
229  Moose::out << "FINAL LOSS: " << store_loss << std::endl;
230  }
231 }
unsigned int _num_tunable
Number of tunable hyperparameters.
void loadHyperParamMap(const std::unordered_map< std::string, Real > &map, const std::unordered_map< std::string, std::vector< Real >> &vec_map)
Load some hyperparameters into the local maps contained in this object.
CovarianceFunctionBase * _covariance_function
Covariance function object.
std::unordered_map< std::string, std::vector< Real > > _hyperparam_vec_map
Vector hyperparameters. Stored for use in surrogate.
const Real eps
void saveState()
std::unordered_map< std::string, Real > _hyperparam_map
Scalar hyperparameters. Stored for use in surrogate.
void seed(std::size_t i, unsigned int seed)
Real getLoss(RealEigenMatrix &inputs, RealEigenMatrix &outputs)
void buildHyperParamMap(std::unordered_map< std::string, Real > &map, std::unordered_map< std::string, std::vector< Real >> &vec_map) const
Populates the input maps with the owned hyperparameters.
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
std::string stringify(const T &t)
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > _tuning_data
Contains tuning inforation. Index of hyperparam, size, and min/max bounds.
unsigned int _batch_size
The batch size for Adam optimization.
std::vector< Real > getGradient(RealEigenMatrix &inputs) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void vecToMap(const std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real >> &tuning_data, std::unordered_map< std::string, Real > &scalar_map, std::unordered_map< std::string, std::vector< Real >> &vector_map, const std::vector< Real > &vec) const
Function used to convert the vectors back to hyperparameter maps.
void mapToVec(const std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real >> &tuning_data, const std::unordered_map< std::string, Real > &scalar_map, const std::unordered_map< std::string, std::vector< Real >> &vector_map, std::vector< Real > &vec) const
Function used to convert the hyperparameter maps in this object to vectors.
MooseUnits pow(const MooseUnits &, int)

◆ tuningData()

std::unordered_map<std::string, std::tuple<unsigned int, unsigned int, Real, Real> >& StochasticTools::GaussianProcess::tuningData ( )
inline

Definition at line 215 of file GaussianProcess.h.

216  {
217  return _tuning_data;
218  }
std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real > > _tuning_data
Contains tuning inforation. Index of hyperparam, size, and min/max bounds.

◆ vecToMap()

void StochasticTools::GaussianProcess::vecToMap ( const std::unordered_map< std::string, std::tuple< unsigned int, unsigned int, Real, Real >> &  tuning_data,
std::unordered_map< std::string, Real > &  scalar_map,
std::unordered_map< std::string, std::vector< Real >> &  vector_map,
const std::vector< Real > &  vec 
) const

Function used to convert the vectors back to hyperparameter maps.

Definition at line 298 of file GaussianProcess.C.

Referenced by tuneHyperParamsAdam().

304 {
305  for (auto iter : tuning_data)
306  {
307  const std::string & param_name = iter.first;
308  if (scalar_map.find(param_name) != scalar_map.end())
309  scalar_map[param_name] = vec[std::get<0>(iter.second)];
310  else if (vector_map.find(param_name) != vector_map.end())
311  for (unsigned int ii = 0; ii < std::get<1>(iter.second); ++ii)
312  vector_map[param_name][ii] = vec[std::get<0>(iter.second) + ii];
313  }
314 }

Member Data Documentation

◆ _batch_size

unsigned int StochasticTools::GaussianProcess::_batch_size
protected

The batch size for Adam optimization.

Definition at line 279 of file GaussianProcess.h.

Referenced by getGradient(), getLoss(), setupCovarianceMatrix(), and tuneHyperParamsAdam().

◆ _covar_name

std::string StochasticTools::GaussianProcess::_covar_name
protected

The name of the covariance function used in this GP.

Definition at line 240 of file GaussianProcess.h.

Referenced by covarName(), getCovarName(), and linkCovarianceFunction().

◆ _covar_type

std::string StochasticTools::GaussianProcess::_covar_type
protected

Type of covariance function used for this GP.

Definition at line 237 of file GaussianProcess.h.

Referenced by covarType(), getCovarType(), and linkCovarianceFunction().

◆ _covariance_function

CovarianceFunctionBase* StochasticTools::GaussianProcess::_covariance_function = nullptr
protected

◆ _data_standardizer

StochasticTools::Standardizer StochasticTools::GaussianProcess::_data_standardizer
protected

Standardizer for use with data (y)

Definition at line 261 of file GaussianProcess.h.

Referenced by dataStandardizer(), getDataStandardizer(), and standardizeData().

◆ _dependent_covar_names

std::vector<UserObjectName> StochasticTools::GaussianProcess::_dependent_covar_names
protected

The names of the covariance functions the used covariance function depends on.

Definition at line 243 of file GaussianProcess.h.

Referenced by dependentCovarNames(), getDependentCovarNames(), and linkCovarianceFunction().

◆ _dependent_covar_types

std::map<UserObjectName, std::string> StochasticTools::GaussianProcess::_dependent_covar_types
protected

The types of the covariance functions the used covariance function depends on.

Definition at line 246 of file GaussianProcess.h.

Referenced by dependentCovarTypes(), getDependentCovarTypes(), and linkCovarianceFunction().

◆ _hyperparam_map

std::unordered_map<std::string, Real> StochasticTools::GaussianProcess::_hyperparam_map
protected

Scalar hyperparameters. Stored for use in surrogate.

Definition at line 252 of file GaussianProcess.h.

Referenced by getHyperParamMap(), hyperparamMap(), setupCovarianceMatrix(), and tuneHyperParamsAdam().

◆ _hyperparam_vec_map

std::unordered_map<std::string, std::vector<Real> > StochasticTools::GaussianProcess::_hyperparam_vec_map
protected

Vector hyperparameters. Stored for use in surrogate.

Definition at line 255 of file GaussianProcess.h.

Referenced by getHyperParamVectorMap(), hyperparamVectorMap(), setupCovarianceMatrix(), and tuneHyperParamsAdam().

◆ _K

RealEigenMatrix StochasticTools::GaussianProcess::_K
protected

An _n_sample by _n_sample covariance matrix constructed from the selected kernel function.

Definition at line 264 of file GaussianProcess.h.

Referenced by getK(), getLoss(), K(), setupCovarianceMatrix(), and setupStoredMatrices().

◆ _K_cho_decomp

Eigen::LLT<RealEigenMatrix> StochasticTools::GaussianProcess::_K_cho_decomp
protected

Cholesky decomposition Eigen object.

Definition at line 270 of file GaussianProcess.h.

Referenced by getGradient(), getKCholeskyDecomp(), KCholeskyDecomp(), and setupStoredMatrices().

◆ _K_results_solve

RealEigenMatrix StochasticTools::GaussianProcess::_K_results_solve
protected

A solve of Ax=b via Cholesky.

Definition at line 267 of file GaussianProcess.h.

Referenced by getGradient(), getKResultsSolve(), getLoss(), KResultsSolve(), and setupStoredMatrices().

◆ _num_outputs

unsigned int StochasticTools::GaussianProcess::_num_outputs
protected

The number of outputs of the GP.

Definition at line 249 of file GaussianProcess.h.

Referenced by covarNumOutputs(), getCovarNumOutputs(), linkCovarianceFunction(), and setupCovarianceMatrix().

◆ _num_tunable

unsigned int StochasticTools::GaussianProcess::_num_tunable
protected

Number of tunable hyperparameters.

Definition at line 234 of file GaussianProcess.h.

Referenced by generateTuningMap(), getGradient(), getNumTunableParams(), and tuneHyperParamsAdam().

◆ _param_standardizer

StochasticTools::Standardizer StochasticTools::GaussianProcess::_param_standardizer
protected

Standardizer for use with params (x)

Definition at line 258 of file GaussianProcess.h.

Referenced by getParamStandardizer(), paramStandardizer(), and standardizeParameters().

◆ _training_data

const RealEigenMatrix* StochasticTools::GaussianProcess::_training_data
protected

Data (y) used for training.

Definition at line 276 of file GaussianProcess.h.

◆ _training_params

const RealEigenMatrix* StochasticTools::GaussianProcess::_training_params
protected

Paramaters (x) used for training, along with statistics.

Definition at line 273 of file GaussianProcess.h.

◆ _tuning_data

std::unordered_map<std::string, std::tuple<unsigned int, unsigned int, Real, Real> > StochasticTools::GaussianProcess::_tuning_data
protected

Contains tuning inforation. Index of hyperparam, size, and min/max bounds.

Definition at line 231 of file GaussianProcess.h.

Referenced by generateTuningMap(), getGradient(), setupCovarianceMatrix(), tuneHyperParamsAdam(), and tuningData().


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