https://mooseframework.inl.gov
OptimizationReporter.C
Go to the documentation of this file.
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 #include "OptimizationReporter.h"
11 
12 #include "libmesh/int_range.h"
13 
14 registerMooseObjectDeprecated("OptimizationApp", OptimizationReporter, "12/31/2024 24:00");
15 
18 {
20  params.addClassDescription("Computes objective function, gradient and contains reporters for "
21  "communicating between optimizeSolve and subapps");
22  params.addRequiredParam<std::vector<dof_id_type>>(
23  "num_values",
24  "Number of parameter values associated with each parameter group in 'parameter_names'.");
25  params.addParam<std::vector<std::vector<Real>>>(
26  "initial_condition",
27  "Initial conditions for each parameter. A vector is given for each parameter group. A "
28  "single value can be given for each group and all parameters in that group will be set to "
29  "that value. The default value is 0.");
30  params.addParam<std::vector<std::vector<Real>>>(
31  "lower_bounds",
32  "Lower bound for each parameter. A vector is given for each parameter group. A single "
33  "value can be given for each group and all parameters in that group will be set to that "
34  "value");
35  params.addParam<std::vector<std::vector<Real>>>(
36  "upper_bounds",
37  "Upper bound for each parameter. A vector is given for each parameter group. A single "
38  "value can be given for each group and all parameters in that group will be set to that "
39  "value");
40  return params;
41 }
42 
45 {
47  "The 'OptimizationReporter' is deprecated. Please use 'GeneralOptimization' instead. "
48  "You can achieve the same functionality by using an 'OptimizationData' object in the "
49  "forward application to calculate the objective value, similar to the method used here.");
50 }
51 void
53 {
54  _nvalues = getParam<std::vector<dof_id_type>>("num_values");
55  _ndof = std::accumulate(_nvalues.begin(), _nvalues.end(), 0);
56 
57  // size checks
58  if (_parameter_names.size() != _nvalues.size())
59  paramError(
60  "num_parameters",
61  "There should be a number in \'num_parameters\' for each name in \'parameter_names\'.");
62 
63  for (const auto & param_id : make_range(_nparams))
64  {
65  _gradients[param_id]->resize(_nvalues[param_id]);
66 
67  std::vector<Real> ic(parseInputData("initial_condition", 0, param_id));
68  std::vector<Real> lb(
69  parseInputData("lower_bounds", std::numeric_limits<Real>::lowest(), param_id));
70  std::vector<Real> ub(
71  parseInputData("upper_bounds", std::numeric_limits<Real>::max(), param_id));
72 
73  _lower_bounds.insert(_lower_bounds.end(), lb.begin(), lb.end());
74  _upper_bounds.insert(_upper_bounds.end(), ub.begin(), ub.end());
75 
76  _parameters[param_id]->assign(ic.begin(), ic.end());
77  }
78 }
79 
80 Real
82 {
83  // This will only be executed if measurement_values are available on the main app
84  for (const auto i : index_range(_measurement_values))
86 
87  Real val = 0.0;
88  for (auto & misfit : _misfit_values)
89  val += misfit * misfit;
90 
91  if (_tikhonov_coeff > 0.0)
92  {
93  Real param_norm_sqr = 0;
94  for (const auto & data : _parameters)
95  for (const auto & val : *data)
96  param_norm_sqr += val * val;
97 
98  val += _tikhonov_coeff * param_norm_sqr;
99  }
100 
101  return val * 0.5;
102 }
103 
104 void
106 {
108 }
registerMooseObjectDeprecated("OptimizationApp", OptimizationReporter, "12/31/2024 24:00")
Computes gradient and contains reporters for communicating between optimizeSolve and subapps...
const std::vector< ReporterValueName > & _parameter_names
Parameter names.
std::vector< Real > & _simulation_values
simulated values at measurement xyzt
void mooseDeprecated(Args &&... args) const
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
dof_id_type _ndof
Total number of parameters.
virtual void setMisfitToSimulatedValues() override
Function to override misfit values with the simulated values from the matrix free hessian forward sol...
void addRequiredParam(const std::string &name, const std::string &doc_string)
const unsigned int _nparams
Number of parameter vectors.
Base class for optimization objects, implements routines for calculating misfit.
static InputParameters validParams()
virtual Real computeObjective() override
Function to compute objective.
std::vector< dof_id_type > _nvalues
Number of values for each parameter.
std::vector< Real > parseInputData(std::string type, Real default_value, unsigned int param_id) const
Function to to parse bounds and initial conditions from input file.
std::vector< Real > _lower_bounds
Bounds of the parameters.
void paramError(const std::string &param, Args... args) const
std::vector< std::vector< Real > * > _parameters
Parameter values declared as reporter data.
virtual void setICsandBounds() override
Sets the initial conditions and bounds right before it is needed.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< std::vector< Real > * > _gradients
Gradient values declared as reporter data.
IntRange< T > make_range(T beg, T end)
void addClassDescription(const std::string &doc_string)
const Real _tikhonov_coeff
Tikhonov Coefficient for regularization.
std::vector< Real > & _misfit_values
difference between simulation and measurement values at measurement xyzt
auto index_range(const T &sizable)
OptimizationReporter(const InputParameters &parameters)
static InputParameters validParams()