https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
ExponentialCovariance Class Reference

#include <ExponentialCovariance.h>

Inheritance diagram for ExponentialCovariance:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 ExponentialCovariance (const InputParameters &parameters)
 
void computeCovarianceMatrix (RealEigenMatrix &K, const RealEigenMatrix &x, const RealEigenMatrix &xp, const bool is_self_covariance) const override
 Generates the Covariance Matrix given two points in the parameter space. More...
 
bool computedKdhyper (RealEigenMatrix &dKdhp, const RealEigenMatrix &x, const std::string &hyper_param_name, unsigned int ind) const override
 Redirect dK/dhp for hyperparameter "hp". More...
 
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. More...
 
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. More...
 
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. More...
 
void dependentCovarianceTypes (std::map< UserObjectName, std::string > &name_type_map) const
 Populate a map with the names and types of the dependent covariance functions. More...
 
const std::vector< UserObjectName > & dependentCovarianceNames () const
 Get the names of the dependent covariances. More...
 
virtual bool isTunable (const std::string &name) const
 Check if a given parameter is tunable. More...
 
unsigned int numOutputs () const
 Return the number of outputs assumed for this covariance function. More...
 
std::unordered_map< std::string, Real > & hyperParamMapReal ()
 Get the map of scalar parameters. More...
 
std::unordered_map< std::string, std::vector< Real > > & hyperParamMapVectorReal ()
 Get the map of vector parameters. More...
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 
static void ExponentialFunction (RealEigenMatrix &K, const RealEigenMatrix &x, const RealEigenMatrix &xp, const std::vector< Real > &length_factor, const Real sigma_f_squared, const Real sigma_n_squared, const Real gamma, const bool is_self_covariance)
 
static void computedKdlf (RealEigenMatrix &K, const RealEigenMatrix &x, const std::vector< Real > &length_factor, const Real sigma_f_squared, const Real gamma, const int ind)
 Computes dK/dlf for individual length factors. More...
 

Public Attributes

const ConsoleStream _console
 

Protected Member Functions

const RealaddRealHyperParameter (const std::string &name, const Real value, const bool is_tunable)
 Register a scalar hyperparameter to this covariance function. More...
 
const std::vector< Real > & addVectorRealHyperParameter (const std::string &name, const std::vector< Real > value, const bool is_tunable)
 Register a vector hyperparameter to this covariance function. More...
 
CovarianceFunctionBasegetCovarianceFunctionByName (const UserObjectName &name) const
 Lookup a CovarianceFunction object by name and return pointer. More...
 

Protected Attributes

std::unordered_map< std::string, Real_hp_map_real
 Map of real-valued hyperparameters. More...
 
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
 Map of vector-valued hyperparameters. More...
 
std::unordered_set< std::string > _tunable_hp
 list of tunable hyper-parameters More...
 
const unsigned int _num_outputs
 The number of outputs this covariance function is used to describe. More...
 
const std::vector< UserObjectName > _dependent_covariance_names
 The names of the dependent covariance functions. More...
 
std::vector< std::string > _dependent_covariance_types
 The types of the dependent covariance functions. More...
 
std::vector< CovarianceFunctionBase * > _covariance_functions
 Vector of pointers to the dependent covariance functions. More...
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
const Parallel::Communicator & _communicator
 

Private Attributes

const std::vector< Real > & _length_factor
 lengh factor () for the kernel, in vector form for multiple parameters More...
 
const Real_sigma_f_squared
 signal variance (^2) More...
 
const Real_sigma_n_squared
 noise variance (^2) More...
 
const Real_gamma
 gamma exponential factor for use in kernel More...
 

Detailed Description

Definition at line 14 of file ExponentialCovariance.h.

Constructor & Destructor Documentation

◆ ExponentialCovariance()

ExponentialCovariance::ExponentialCovariance ( const InputParameters parameters)

Definition at line 30 of file ExponentialCovariance.C.

33  "length_factor", getParam<std::vector<Real>>("length_factor"), true)),
35  addRealHyperParameter("signal_variance", getParam<Real>("signal_variance"), true)),
37  addRealHyperParameter("noise_variance", getParam<Real>("noise_variance"), true)),
38  _gamma(addRealHyperParameter("gamma", getParam<Real>("gamma"), false))
39 {
40 }
const Real & _gamma
gamma exponential factor for use in kernel
const std::vector< Real > & _length_factor
lengh factor () for the kernel, in vector form for multiple parameters
CovarianceFunctionBase(const InputParameters &parameters)
const std::vector< Real > & addVectorRealHyperParameter(const std::string &name, const std::vector< Real > value, const bool is_tunable)
Register a vector hyperparameter to this covariance function.
const T & getParam(const std::string &name) const
const Real & _sigma_n_squared
noise variance (^2)
const InputParameters & parameters() const
const Real & _sigma_f_squared
signal variance (^2)
const Real & addRealHyperParameter(const std::string &name, const Real value, const bool is_tunable)
Register a scalar hyperparameter to this covariance function.

Member Function Documentation

◆ addRealHyperParameter()

const Real & CovarianceFunctionBase::addRealHyperParameter ( const std::string &  name,
const Real  value,
const bool  is_tunable 
)
protectedinherited

Register a scalar hyperparameter to this covariance function.

Parameters
nameThe name of the parameter
valueThe initial value of the parameter
is_tunableIf the parameter is tunable during optimization

Definition at line 52 of file CovarianceFunctionBase.C.

55 {
56  const auto prefixed_name = _name + ":" + name;
57  if (is_tunable)
58  _tunable_hp.insert(prefixed_name);
59  return _hp_map_real.emplace(prefixed_name, value).first->second;
60 }
std::unordered_map< std::string, Real > _hp_map_real
Map of real-valued hyperparameters.
std::unordered_set< std::string > _tunable_hp
list of tunable hyper-parameters
virtual const std::string & name() const
const std::string _name

◆ addVectorRealHyperParameter()

const std::vector< Real > & CovarianceFunctionBase::addVectorRealHyperParameter ( const std::string &  name,
const std::vector< Real value,
const bool  is_tunable 
)
protectedinherited

Register a vector hyperparameter to this covariance function.

Parameters
nameThe name of the parameter
valueThe initial value of the parameter
is_tunableIf the parameter is tunable during optimization

Definition at line 63 of file CovarianceFunctionBase.C.

Referenced by LMC::LMC().

66 {
67  const auto prefixed_name = _name + ":" + name;
68  if (is_tunable)
69  _tunable_hp.insert(prefixed_name);
70  return _hp_map_vector_real.emplace(prefixed_name, value).first->second;
71 }
std::unordered_set< std::string > _tunable_hp
list of tunable hyper-parameters
virtual const std::string & name() const
const std::string _name
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
Map of vector-valued hyperparameters.

◆ buildHyperParamMap()

void CovarianceFunctionBase::buildHyperParamMap ( std::unordered_map< std::string, Real > &  map,
std::unordered_map< std::string, std::vector< Real >> &  vec_map 
) const
inherited

Populates the input maps with the owned hyperparameters.

Parameters
mapMap of scalar hyperparameters that should be populated
vec_mapMap of vector hyperparameters that should be populated

Definition at line 115 of file CovarianceFunctionBase.C.

Referenced by StochasticTools::GaussianProcess::setupCovarianceMatrix(), and StochasticTools::GaussianProcess::tuneHyperParamsAdam().

118 {
119  // First, add the hyperparameters of the dependent covariance functions
120  for (const auto dependent_covar : _covariance_functions)
121  dependent_covar->buildHyperParamMap(map, vec_map);
122 
123  // At the end we just append the hyperparameters this object owns
124  for (const auto & iter : _hp_map_real)
125  map[iter.first] = iter.second;
126  for (const auto & iter : _hp_map_vector_real)
127  vec_map[iter.first] = iter.second;
128 }
std::unordered_map< std::string, Real > _hp_map_real
Map of real-valued hyperparameters.
std::vector< CovarianceFunctionBase * > _covariance_functions
Vector of pointers to the dependent covariance functions.
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
Map of vector-valued hyperparameters.

◆ computeCovarianceMatrix()

void ExponentialCovariance::computeCovarianceMatrix ( RealEigenMatrix K,
const RealEigenMatrix x,
const RealEigenMatrix xp,
const bool  is_self_covariance 
) const
overridevirtual

Generates the Covariance Matrix given two points in the parameter space.

Implements CovarianceFunctionBase.

Definition at line 43 of file ExponentialCovariance.C.

47 {
48  if ((unsigned)x.cols() != _length_factor.size())
49  mooseError("length_factor size does not match dimension of trainer input.");
50 
52  K, x, xp, _length_factor, _sigma_f_squared, _sigma_n_squared, _gamma, is_self_covariance);
53 }
const Real & _gamma
gamma exponential factor for use in kernel
const std::vector< Real > & _length_factor
lengh factor () for the kernel, in vector form for multiple parameters
static const std::string K
Definition: NS.h:170
const std::vector< double > x
static void ExponentialFunction(RealEigenMatrix &K, const RealEigenMatrix &x, const RealEigenMatrix &xp, const std::vector< Real > &length_factor, const Real sigma_f_squared, const Real sigma_n_squared, const Real gamma, const bool is_self_covariance)
const Real & _sigma_n_squared
noise variance (^2)
void mooseError(Args &&... args) const
const Real & _sigma_f_squared
signal variance (^2)

◆ computedKdhyper()

bool ExponentialCovariance::computedKdhyper ( RealEigenMatrix dKdhp,
const RealEigenMatrix x,
const std::string &  hyper_param_name,
unsigned int  ind 
) const
overridevirtual

Redirect dK/dhp for hyperparameter "hp".

Reimplemented from CovarianceFunctionBase.

Definition at line 89 of file ExponentialCovariance.C.

93 {
94  if (name().length() + 1 > hyper_param_name.length())
95  return false;
96 
97  const std::string name_without_prefix = hyper_param_name.substr(name().length() + 1);
98 
99  if (name_without_prefix == "noise_variance")
100  {
101  ExponentialFunction(dKdhp, x, x, _length_factor, 0, 1, _gamma, true);
102  return true;
103  }
104 
105  if (name_without_prefix == "signal_variance")
106  {
107  ExponentialFunction(dKdhp, x, x, _length_factor, 1, 0, _gamma, false);
108  return true;
109  }
110 
111  if (name_without_prefix == "length_factor")
112  {
114  return true;
115  }
116 
117  return false;
118 }
const Real & _gamma
gamma exponential factor for use in kernel
const std::vector< Real > & _length_factor
lengh factor () for the kernel, in vector form for multiple parameters
virtual const std::string & name() const
const std::vector< double > x
static void ExponentialFunction(RealEigenMatrix &K, const RealEigenMatrix &x, const RealEigenMatrix &xp, const std::vector< Real > &length_factor, const Real sigma_f_squared, const Real sigma_n_squared, const Real gamma, const bool is_self_covariance)
static void computedKdlf(RealEigenMatrix &K, const RealEigenMatrix &x, const std::vector< Real > &length_factor, const Real sigma_f_squared, const Real gamma, const int ind)
Computes dK/dlf for individual length factors.
const Real & _sigma_f_squared
signal variance (^2)

◆ computedKdlf()

void ExponentialCovariance::computedKdlf ( RealEigenMatrix K,
const RealEigenMatrix x,
const std::vector< Real > &  length_factor,
const Real  sigma_f_squared,
const Real  gamma,
const int  ind 
)
static

Computes dK/dlf for individual length factors.

Definition at line 121 of file ExponentialCovariance.C.

Referenced by computedKdhyper().

127 {
128  unsigned int num_samples_x = x.rows();
129  unsigned int num_params_x = x.cols();
130 
131  mooseAssert(ind < x.cols(), "Incorrect length factor index");
132 
133  for (unsigned int ii = 0; ii < num_samples_x; ++ii)
134  {
135  for (unsigned int jj = 0; jj < num_samples_x; ++jj)
136  {
137  // Compute distance per parameter, scaled by length factor
138  Real r_scaled = 0;
139  for (unsigned int kk = 0; kk < num_params_x; ++kk)
140  r_scaled += pow((x(ii, kk) - x(jj, kk)) / length_factor[kk], 2);
141  r_scaled = sqrt(r_scaled);
142  if (r_scaled != 0)
143  {
144  K(ii, jj) = gamma * std::pow(r_scaled, gamma - 2) * sigma_f_squared *
145  std::exp(-pow(r_scaled, gamma));
146  K(ii, jj) =
147  std::pow(x(ii, ind) - x(jj, ind), 2) / std::pow(length_factor[ind], 3) * K(ii, jj);
148  }
149  else // avoid div by 0. 0/0=0 scenario.
150  K(ii, jj) = 0;
151  }
152  }
153 }
static const std::string K
Definition: NS.h:170
const std::vector< double > x
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
MooseUnits pow(const MooseUnits &, int)

◆ dependentCovarianceNames()

const std::vector<UserObjectName>& CovarianceFunctionBase::dependentCovarianceNames ( ) const
inlineinherited

Get the names of the dependent covariances.

Definition at line 62 of file CovarianceFunctionBase.h.

Referenced by StochasticTools::GaussianProcess::linkCovarianceFunction().

63  {
65  }
const std::vector< UserObjectName > _dependent_covariance_names
The names of the dependent covariance functions.

◆ dependentCovarianceTypes()

void CovarianceFunctionBase::dependentCovarianceTypes ( std::map< UserObjectName, std::string > &  name_type_map) const
inherited

Populate a map with the names and types of the dependent covariance functions.

Parameters
name_type_mapReference to the map which should be populated

Definition at line 163 of file CovarianceFunctionBase.C.

Referenced by StochasticTools::GaussianProcess::linkCovarianceFunction().

165 {
166  for (const auto dependent_covar : _covariance_functions)
167  {
168  dependent_covar->dependentCovarianceTypes(name_type_map);
169  name_type_map.insert(std::make_pair(dependent_covar->name(), dependent_covar->type()));
170  }
171 }
std::vector< CovarianceFunctionBase * > _covariance_functions
Vector of pointers to the dependent covariance functions.

◆ ExponentialFunction()

void ExponentialCovariance::ExponentialFunction ( RealEigenMatrix K,
const RealEigenMatrix x,
const RealEigenMatrix xp,
const std::vector< Real > &  length_factor,
const Real  sigma_f_squared,
const Real  sigma_n_squared,
const Real  gamma,
const bool  is_self_covariance 
)
static

Definition at line 56 of file ExponentialCovariance.C.

Referenced by computeCovarianceMatrix(), and computedKdhyper().

64 {
65  unsigned int num_samples_x = x.rows();
66  unsigned int num_samples_xp = xp.rows();
67  unsigned int num_params_x = x.cols();
68 
69  mooseAssert(num_params_x == xp.cols(),
70  "Number of parameters do not match in covariance kernel calculation");
71 
72  for (unsigned int ii = 0; ii < num_samples_x; ++ii)
73  {
74  for (unsigned int jj = 0; jj < num_samples_xp; ++jj)
75  {
76  // Compute distance per parameter, scaled by length factor
77  Real r_scaled = 0;
78  for (unsigned int kk = 0; kk < num_params_x; ++kk)
79  r_scaled += pow((x(ii, kk) - xp(jj, kk)) / length_factor[kk], 2);
80  r_scaled = sqrt(r_scaled);
81  K(ii, jj) = sigma_f_squared * std::exp(-pow(r_scaled, gamma));
82  }
83  if (is_self_covariance)
84  K(ii, ii) += sigma_n_squared;
85  }
86 }
static const std::string K
Definition: NS.h:170
const std::vector< double > x
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh

◆ getCovarianceFunctionByName()

CovarianceFunctionBase * CovarianceInterface::getCovarianceFunctionByName ( const UserObjectName &  name) const
protectedinherited

Lookup a CovarianceFunction object by name and return pointer.

Definition at line 25 of file CovarianceInterface.C.

Referenced by ActiveLearningGaussianProcess::ActiveLearningGaussianProcess(), CovarianceFunctionBase::CovarianceFunctionBase(), GaussianProcessTrainer::GaussianProcessTrainer(), and GaussianProcessSurrogate::setupCovariance().

26 {
27  std::vector<CovarianceFunctionBase *> models;
29  .query()
30  .condition<AttribName>(name)
31  .condition<AttribSystem>("CovarianceFunction")
32  .queryInto(models);
33  if (models.empty())
34  mooseError("Unable to find a CovarianceFunction object with the name '" + name + "'");
35  return models[0];
36 }
void mooseError(Args &&... args)
TheWarehouse & theWarehouse() const
const std::string name
Definition: Setup.h:20
Query query()
FEProblemBase & _covar_feproblem
Reference to FEProblemBase instance.

◆ getTuningData()

bool CovarianceFunctionBase::getTuningData ( const std::string &  name,
unsigned int size,
Real min,
Real max 
) const
virtualinherited

Get the default minimum and maximum and size of a hyperparameter.

Returns false is the parameter has not been found in this covariance object.

Parameters
nameThe name of the hyperparameter
sizeReference to an unsigned int that will contain the size of the hyperparameter (will be populated with 1 if it is scalar)
minReference to a number which will be populated by the maximum allowed value of the hyperparameter
maxReference to a number which will be populated by the minimum allowed value of the hyperparameter

Definition at line 131 of file CovarianceFunctionBase.C.

Referenced by StochasticTools::GaussianProcess::generateTuningMap().

135 {
136  // First, check the dependent covariances
137  for (const auto dependent_covar : _covariance_functions)
138  if (dependent_covar->getTuningData(name, size, min, max))
139  return true;
140 
141  min = 1e-9;
142  max = 1e9;
143 
144  if (_hp_map_real.find(name) != _hp_map_real.end())
145  {
146  size = 1;
147  return true;
148  }
149  else if (_hp_map_vector_real.find(name) != _hp_map_vector_real.end())
150  {
151  const auto & vector_value = _hp_map_vector_real.find(name);
152  size = vector_value->second.size();
153  return true;
154  }
155  else
156  {
157  size = 0;
158  return false;
159  }
160 }
std::unordered_map< std::string, Real > _hp_map_real
Map of real-valued hyperparameters.
std::vector< CovarianceFunctionBase * > _covariance_functions
Vector of pointers to the dependent covariance functions.
virtual const std::string & name() const
auto max(const L &left, const R &right)
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
Map of vector-valued hyperparameters.
auto min(const L &left, const R &right)

◆ hyperParamMapReal()

std::unordered_map<std::string, Real>& CovarianceFunctionBase::hyperParamMapReal ( )
inlineinherited

Get the map of scalar parameters.

Definition at line 88 of file CovarianceFunctionBase.h.

88 { return _hp_map_real; }
std::unordered_map< std::string, Real > _hp_map_real
Map of real-valued hyperparameters.

◆ hyperParamMapVectorReal()

std::unordered_map<std::string, std::vector<Real> >& CovarianceFunctionBase::hyperParamMapVectorReal ( )
inlineinherited

Get the map of vector parameters.

Definition at line 91 of file CovarianceFunctionBase.h.

92  {
93  return _hp_map_vector_real;
94  }
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
Map of vector-valued hyperparameters.

◆ isTunable()

bool CovarianceFunctionBase::isTunable ( const std::string &  name) const
virtualinherited

Check if a given parameter is tunable.

Parameters
Thename of the hyperparameter

Definition at line 74 of file CovarianceFunctionBase.C.

Referenced by StochasticTools::GaussianProcess::generateTuningMap().

75 {
76  // First, we check if the dependent covariances have the parameter
77  for (const auto dependent_covar : _covariance_functions)
78  if (dependent_covar->isTunable(name))
79  return true;
80 
81  if (_tunable_hp.find(name) != _tunable_hp.end())
82  return true;
83  else if (_hp_map_real.find(name) != _hp_map_real.end() ||
85  mooseError("We found hyperparameter ", name, " but it was not declared tunable!");
86 
87  return false;
88 }
std::unordered_map< std::string, Real > _hp_map_real
Map of real-valued hyperparameters.
std::unordered_set< std::string > _tunable_hp
list of tunable hyper-parameters
std::vector< CovarianceFunctionBase * > _covariance_functions
Vector of pointers to the dependent covariance functions.
virtual const std::string & name() const
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
Map of vector-valued hyperparameters.
void mooseError(Args &&... args) const

◆ loadHyperParamMap()

void CovarianceFunctionBase::loadHyperParamMap ( const std::unordered_map< std::string, Real > &  map,
const std::unordered_map< std::string, std::vector< Real >> &  vec_map 
)
inherited

Load some hyperparameters into the local maps contained in this object.

Parameters
mapInput map of scalar hyperparameters
vec_mapInput map of vector hyperparameters

Definition at line 91 of file CovarianceFunctionBase.C.

Referenced by LoadCovarianceDataAction::load(), and StochasticTools::GaussianProcess::tuneHyperParamsAdam().

94 {
95  // First, load the hyperparameters of the dependent covariance functions
96  for (const auto dependent_covar : _covariance_functions)
97  dependent_covar->loadHyperParamMap(map, vec_map);
98 
99  // Then we load the hyperparameters of this object
100  for (auto & iter : _hp_map_real)
101  {
102  const auto & map_iter = map.find(iter.first);
103  if (map_iter != map.end())
104  iter.second = map_iter->second;
105  }
106  for (auto & iter : _hp_map_vector_real)
107  {
108  const auto & map_iter = vec_map.find(iter.first);
109  if (map_iter != vec_map.end())
110  iter.second = map_iter->second;
111  }
112 }
std::unordered_map< std::string, Real > _hp_map_real
Map of real-valued hyperparameters.
std::vector< CovarianceFunctionBase * > _covariance_functions
Vector of pointers to the dependent covariance functions.
std::unordered_map< std::string, std::vector< Real > > _hp_map_vector_real
Map of vector-valued hyperparameters.

◆ numOutputs()

unsigned int CovarianceFunctionBase::numOutputs ( ) const
inlineinherited

Return the number of outputs assumed for this covariance function.

Definition at line 85 of file CovarianceFunctionBase.h.

Referenced by GaussianProcessSurrogate::evaluate(), GaussianProcessTrainer::GaussianProcessTrainer(), and StochasticTools::GaussianProcess::linkCovarianceFunction().

85 { return _num_outputs; }
const unsigned int _num_outputs
The number of outputs this covariance function is used to describe.

◆ validParams()

InputParameters ExponentialCovariance::validParams ( )
static

Definition at line 16 of file ExponentialCovariance.C.

17 {
19  params.addClassDescription("Exponential covariance function.");
20  params.addRequiredParam<std::vector<Real>>("length_factor",
21  "Length factors to use for Covariance Kernel");
22  params.addRequiredParam<Real>("signal_variance",
23  "Signal Variance ($\\sigma_f^2$) to use for kernel calculation.");
24  params.addParam<Real>(
25  "noise_variance", 0.0, "Noise Variance ($\\sigma_n^2$) to use for kernel calculation.");
26  params.addRequiredParam<Real>("gamma", "Gamma to use for Exponential Covariance Kernel");
27  return params;
28 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
void addRequiredParam(const std::string &name, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _covariance_functions

std::vector<CovarianceFunctionBase *> CovarianceFunctionBase::_covariance_functions
protectedinherited

◆ _dependent_covariance_names

const std::vector<UserObjectName> CovarianceFunctionBase::_dependent_covariance_names
protectedinherited

The names of the dependent covariance functions.

Definition at line 125 of file CovarianceFunctionBase.h.

Referenced by CovarianceFunctionBase::CovarianceFunctionBase(), and CovarianceFunctionBase::dependentCovarianceNames().

◆ _dependent_covariance_types

std::vector<std::string> CovarianceFunctionBase::_dependent_covariance_types
protectedinherited

The types of the dependent covariance functions.

Definition at line 128 of file CovarianceFunctionBase.h.

Referenced by CovarianceFunctionBase::CovarianceFunctionBase().

◆ _gamma

const Real& ExponentialCovariance::_gamma
private

gamma exponential factor for use in kernel

Definition at line 60 of file ExponentialCovariance.h.

Referenced by computeCovarianceMatrix(), and computedKdhyper().

◆ _hp_map_real

std::unordered_map<std::string, Real> CovarianceFunctionBase::_hp_map_real
protectedinherited

◆ _hp_map_vector_real

std::unordered_map<std::string, std::vector<Real> > CovarianceFunctionBase::_hp_map_vector_real
protectedinherited

◆ _length_factor

const std::vector<Real>& ExponentialCovariance::_length_factor
private

lengh factor () for the kernel, in vector form for multiple parameters

Definition at line 51 of file ExponentialCovariance.h.

Referenced by computeCovarianceMatrix(), and computedKdhyper().

◆ _num_outputs

const unsigned int CovarianceFunctionBase::_num_outputs
protectedinherited

The number of outputs this covariance function is used to describe.

Definition at line 122 of file CovarianceFunctionBase.h.

Referenced by LMC::computeAGradient(), LMC::computeBMatrix(), LMC::computeCovarianceMatrix(), LMC::computedKdhyper(), LMC::LMC(), and CovarianceFunctionBase::numOutputs().

◆ _sigma_f_squared

const Real& ExponentialCovariance::_sigma_f_squared
private

signal variance (^2)

Definition at line 54 of file ExponentialCovariance.h.

Referenced by computeCovarianceMatrix(), and computedKdhyper().

◆ _sigma_n_squared

const Real& ExponentialCovariance::_sigma_n_squared
private

noise variance (^2)

Definition at line 57 of file ExponentialCovariance.h.

Referenced by computeCovarianceMatrix().

◆ _tunable_hp

std::unordered_set<std::string> CovarianceFunctionBase::_tunable_hp
protectedinherited

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