Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #pragma once
11 :
12 : #include "Standardizer.h"
13 : #include <Eigen/Dense>
14 :
15 : #include "CovarianceFunctionBase.h"
16 :
17 : namespace StochasticTools
18 : {
19 :
20 : /**
21 : * Utility class dedicated to hold structures and functions commont to
22 : * Gaussian Processes. It can be used to standardize parameters, manipulate
23 : * covariance data and compute additional stored matrices.
24 : */
25 : class GaussianProcess
26 : {
27 : public:
28 : GaussianProcess();
29 :
30 : /**
31 : * Initializes the most important structures in the Gaussian Process: the
32 : * covariance function and a tuning map which is used if the user requires
33 : * parameter tuning.
34 : * @param covariance_function Pointer to the covariance function that
35 : * needs to be used for the Gaussian Process.
36 : * @param params_to_tune List of parameters which need to be tuned.
37 : * @param min List of lower bounds for the parameter tuning.
38 : * @param max List of upper bounds for parameter tuning.
39 : */
40 : void initialize(CovarianceFunctionBase * covariance_function,
41 : const std::vector<std::string> & params_to_tune,
42 : const std::vector<Real> & min = std::vector<Real>(),
43 : const std::vector<Real> & max = std::vector<Real>());
44 :
45 : /// Structure containing the optimization options for
46 : /// hyperparameter-tuning
47 : struct GPOptimizerOptions
48 : {
49 : /// Default constructor
50 : GPOptimizerOptions();
51 : /**
52 : * Construct a new GPOptimizerOptions object using
53 : * input parameters that will control the optimization
54 : * @param show_every_nth_iteration To show the loss value at every n-th iteration, if set to 0,
55 : * nothing is displayed
56 : * @param num_iter The number of iterations we want in the optimization of the GP
57 : * @param batch_size The number of samples in each batch
58 : * @param learning_rate The learning rate for parameter updates
59 : * @param b1 Tuning constant for the Adam algorithm
60 : * @param b2 Tuning constant for the Adam algorithm
61 : * @param eps Tuning constant for the Adam algorithm
62 : * @param lambda Tuning constant for the Adam algorithm
63 : */
64 : GPOptimizerOptions(const bool show_every_nth_iteration = 1,
65 : const unsigned int num_iter = 1000,
66 : const unsigned int batch_size = 0,
67 : const Real learning_rate = 1e-3,
68 : const Real b1 = 0.9,
69 : const Real b2 = 0.999,
70 : const Real eps = 1e-7,
71 : const Real lambda = 0.0);
72 :
73 : /// Switch to enable verbose output for parameter tuning at every n-th iteration
74 : const unsigned int show_every_nth_iteration = false;
75 : /// The number of iterations for Adam optimizer
76 : const unsigned int num_iter = 1000;
77 : /// The batch isize for Adam optimizer
78 : const unsigned int batch_size = 0;
79 : /// The learning rate for Adam optimizer
80 : const Real learning_rate = 1e-3;
81 : /// Tuning parameter from the paper
82 : const Real b1 = 0.9;
83 : /// Tuning parameter from the paper
84 : const Real b2 = 0.999;
85 : /// Tuning parameter from the paper
86 : const Real eps = 1e-7;
87 : /// Tuning parameter from the paper
88 : const Real lambda = 0.0;
89 : };
90 : /**
91 : * Sets up the covariance matrix given data and optimization options.
92 : * @param training_params The training parameter values (x values) for the
93 : * covariance matrix.
94 : * @param training_data The training data (y values) for the inversion of the
95 : * covariance matrix.
96 : * @param opts The optimizer options.
97 : */
98 : void setupCovarianceMatrix(const RealEigenMatrix & training_params,
99 : const RealEigenMatrix & training_data,
100 : const GPOptimizerOptions & opts);
101 :
102 : /**
103 : * Sets up the Cholesky decomposition and inverse action of the covariance matrix.
104 : * @param input The vector/matrix which right multiples the inverse of the covariance matrix.
105 : */
106 : void setupStoredMatrices(const RealEigenMatrix & input);
107 :
108 : /**
109 : * Finds and links the covariance function to this object. Used mainly in the
110 : * covariance data action.
111 : * @param covariance_function Pointer to the covariance function that
112 : * needs to be used for the Gaussian Process.
113 : */
114 : void linkCovarianceFunction(CovarianceFunctionBase * covariance_function);
115 :
116 : /**
117 : * Sets up the tuning map which is used if the user requires parameter tuning.
118 : * @param params_to_tune List of parameters which need to be tuned.
119 : * @param min List of lower bounds for the parameter tuning.
120 : * @param max List of upper bounds for parameter tuning.
121 : */
122 : void generateTuningMap(const std::vector<std::string> & params_to_tune,
123 : const std::vector<Real> & min = std::vector<Real>(),
124 : const std::vector<Real> & max = std::vector<Real>());
125 :
126 : /**
127 : * Standardizes the vector of input parameters (x values).
128 : * @param parameters The vector/matrix of input data.
129 : * @param keep_moments If previously computed or new moments are to be used.
130 : */
131 : void standardizeParameters(RealEigenMatrix & parameters, bool keep_moments = false);
132 :
133 : /**
134 : * Standardizes the vector of responses (y values).
135 : * @param data The vector/matrix of input data.
136 : * @param keep_moments If previously computed or new moments are to be used.
137 : */
138 : void standardizeData(RealEigenMatrix & data, bool keep_moments = false);
139 :
140 : // Tune hyperparameters using Adam
141 : void tuneHyperParamsAdam(const RealEigenMatrix & training_params,
142 : const RealEigenMatrix & training_data,
143 : const GPOptimizerOptions & opts);
144 :
145 : // Computes the loss function
146 : Real getLoss(RealEigenMatrix & inputs, RealEigenMatrix & outputs);
147 :
148 : // Computes Gradient of the loss function
149 : std::vector<Real> getGradient(RealEigenMatrix & inputs) const;
150 :
151 : /// Function used to convert the hyperparameter maps in this object to
152 : /// vectors
153 : void mapToVec(
154 : const std::unordered_map<std::string, std::tuple<unsigned int, unsigned int, Real, Real>> &
155 : tuning_data,
156 : const std::unordered_map<std::string, Real> & scalar_map,
157 : const std::unordered_map<std::string, std::vector<Real>> & vector_map,
158 : std::vector<Real> & vec) const;
159 :
160 : /// Function used to convert the vectors back to hyperparameter maps
161 : void vecToMap(
162 : const std::unordered_map<std::string, std::tuple<unsigned int, unsigned int, Real, Real>> &
163 : tuning_data,
164 : std::unordered_map<std::string, Real> & scalar_map,
165 : std::unordered_map<std::string, std::vector<Real>> & vector_map,
166 : const std::vector<Real> & vec) const;
167 :
168 : /// @{
169 : /**
170 : * Get constant reference to the contained structures
171 : */
172 115964 : const StochasticTools::Standardizer & getParamStandardizer() const { return _param_standardizer; }
173 231928 : const StochasticTools::Standardizer & getDataStandardizer() const { return _data_standardizer; }
174 : const RealEigenMatrix & getK() const { return _K; }
175 115964 : const RealEigenMatrix & getKResultsSolve() const { return _K_results_solve; }
176 115964 : const Eigen::LLT<RealEigenMatrix> & getKCholeskyDecomp() const { return _K_cho_decomp; }
177 348132 : const CovarianceFunctionBase & getCovarFunction() const { return *_covariance_function; }
178 32 : const CovarianceFunctionBase * getCovarFunctionPtr() const { return _covariance_function; }
179 32 : const std::string & getCovarType() const { return _covar_type; }
180 32 : const std::string & getCovarName() const { return _covar_name; }
181 : const std::vector<UserObjectName> & getDependentCovarNames() const
182 : {
183 32 : return _dependent_covar_names;
184 : }
185 : const std::map<UserObjectName, std::string> & getDependentCovarTypes() const
186 : {
187 : return _dependent_covar_types;
188 : }
189 : const unsigned int & getCovarNumOutputs() const { return _num_outputs; }
190 : const unsigned int & getNumTunableParams() const { return _num_tunable; }
191 32 : const std::unordered_map<std::string, Real> & getHyperParamMap() const { return _hyperparam_map; }
192 : const std::unordered_map<std::string, std::vector<Real>> & getHyperParamVectorMap() const
193 : {
194 32 : return _hyperparam_vec_map;
195 : }
196 : ///@}
197 :
198 : /// @{
199 : /**
200 : * Get non-constant reference to the contained structures (if they need to be modified from the
201 : * utside)
202 : */
203 52 : StochasticTools::Standardizer & paramStandardizer() { return _param_standardizer; }
204 52 : StochasticTools::Standardizer & dataStandardizer() { return _data_standardizer; }
205 52 : RealEigenMatrix & K() { return _K; }
206 52 : RealEigenMatrix & KResultsSolve() { return _K_results_solve; }
207 52 : Eigen::LLT<RealEigenMatrix> & KCholeskyDecomp() { return _K_cho_decomp; }
208 : CovarianceFunctionBase * covarFunctionPtr() { return _covariance_function; }
209 : CovarianceFunctionBase & covarFunction() { return *_covariance_function; }
210 52 : std::string & covarType() { return _covar_type; }
211 52 : std::string & covarName() { return _covar_name; }
212 52 : std::map<UserObjectName, std::string> & dependentCovarTypes() { return _dependent_covar_types; }
213 52 : std::vector<UserObjectName> & dependentCovarNames() { return _dependent_covar_names; }
214 52 : unsigned int & covarNumOutputs() { return _num_outputs; }
215 : std::unordered_map<std::string, std::tuple<unsigned int, unsigned int, Real, Real>> & tuningData()
216 : {
217 : return _tuning_data;
218 : }
219 52 : std::unordered_map<std::string, Real> & hyperparamMap() { return _hyperparam_map; }
220 : std::unordered_map<std::string, std::vector<Real>> & hyperparamVectorMap()
221 : {
222 52 : return _hyperparam_vec_map;
223 : }
224 : ///@}
225 :
226 : protected:
227 : /// Covariance function object
228 : CovarianceFunctionBase * _covariance_function = nullptr;
229 :
230 : /// Contains tuning inforation. Index of hyperparam, size, and min/max bounds
231 : std::unordered_map<std::string, std::tuple<unsigned int, unsigned int, Real, Real>> _tuning_data;
232 :
233 : /// Number of tunable hyperparameters
234 : unsigned int _num_tunable;
235 :
236 : /// Type of covariance function used for this GP
237 : std::string _covar_type;
238 :
239 : /// The name of the covariance function used in this GP
240 : std::string _covar_name;
241 :
242 : /// The names of the covariance functions the used covariance function depends on
243 : std::vector<UserObjectName> _dependent_covar_names;
244 :
245 : /// The types of the covariance functions the used covariance function depends on
246 : std::map<UserObjectName, std::string> _dependent_covar_types;
247 :
248 : /// The number of outputs of the GP
249 : unsigned int _num_outputs;
250 :
251 : /// Scalar hyperparameters. Stored for use in surrogate
252 : std::unordered_map<std::string, Real> _hyperparam_map;
253 :
254 : /// Vector hyperparameters. Stored for use in surrogate
255 : std::unordered_map<std::string, std::vector<Real>> _hyperparam_vec_map;
256 :
257 : /// Standardizer for use with params (x)
258 : StochasticTools::Standardizer _param_standardizer;
259 :
260 : /// Standardizer for use with data (y)
261 : StochasticTools::Standardizer _data_standardizer;
262 :
263 : /// An _n_sample by _n_sample covariance matrix constructed from the selected kernel function
264 : RealEigenMatrix _K;
265 :
266 : /// A solve of Ax=b via Cholesky.
267 : RealEigenMatrix _K_results_solve;
268 :
269 : /// Cholesky decomposition Eigen object
270 : Eigen::LLT<RealEigenMatrix> _K_cho_decomp;
271 :
272 : /// Paramaters (x) used for training, along with statistics
273 : const RealEigenMatrix * _training_params;
274 :
275 : /// Data (y) used for training
276 : const RealEigenMatrix * _training_data;
277 :
278 : /// The batch size for Adam optimization
279 : unsigned int _batch_size;
280 : };
281 :
282 : } // StochasticTools namespac
283 :
284 : template <>
285 : void dataStore(std::ostream & stream, Eigen::LLT<RealEigenMatrix> & decomp, void * context);
286 : template <>
287 : void dataLoad(std::istream & stream, Eigen::LLT<RealEigenMatrix> & decomp, void * context);
288 :
289 : template <>
290 : void dataStore(std::ostream & stream, StochasticTools::GaussianProcess & gp_utils, void * context);
291 : template <>
292 : void dataLoad(std::istream & stream, StochasticTools::GaussianProcess & gp_utils, void * context);
|