https://mooseframework.inl.gov
GeneralOptimization.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 "MooseError.h"
11 #include "MooseTypes.h"
12 #include "GeneralOptimization.h"
13 #include "libmesh/id_types.h"
14 #include <numeric>
15 
16 registerMooseObject("OptimizationApp", GeneralOptimization);
17 
20 {
22 
23  params.addRequiredParam<ReporterValueName>(
24  "objective_name", "Preferred name of reporter value defining the objective.");
25  params.addParam<std::vector<dof_id_type>>(
26  "num_values",
27  "Number of parameter values associated with each parameter group in 'parameter_names'.");
28  params.addParam<ReporterValueName>("num_values_name",
29  "Reporter that holds the number of parameter values "
30  "associated with each parameter group in 'parameter_names'.");
31  params.addParam<std::vector<std::vector<Real>>>(
32  "initial_condition",
33  "Initial conditions for each parameter. A vector is given for each parameter group. A "
34  "single value can be given for each group and all parameters in that group will be set to "
35  "that value. The default value is 0.");
36  params.addParam<std::vector<std::vector<Real>>>(
37  "lower_bounds",
38  "Lower bound for each parameter. A vector is given for each parameter group. A single "
39  "value can be given for each group and all parameters in that group will be set to that "
40  "value");
41  params.addParam<std::vector<std::vector<Real>>>(
42  "upper_bounds",
43  "Upper bound for each parameter. A vector is given for each parameter group. A single "
44  "value can be given for each group and all parameters in that group will be set to that "
45  "value");
46 
47  params.addClassDescription("Reporter that provides TAO with the objective, gradient, and "
48  "constraint data, which are supplied by the reporters and "
49  "postprocessors from the forward and adjoint subapps.");
50  return params;
51 }
52 
54  : OptimizationReporterBase(parameters),
55  _objective_val(declareValueByName<Real>(getParam<ReporterValueName>("objective_name"),
57  _num_values_reporter(
58  isParamValid("num_values_name")
59  ? &declareValueByName<std::vector<dof_id_type>>(
60  getParam<ReporterValueName>("num_values_name"), REPORTER_MODE_REPLICATED)
61  : nullptr)
62 {
63 }
64 
65 Real
67 {
68  Real val = 0;
69  if (_tikhonov_coeff > 0.0)
70  {
71  Real param_norm_sqr = 0;
72  for (const auto & data : _parameters)
73  for (const auto & param_val : *data)
74  param_norm_sqr += param_val * param_val;
75  // We multiply by 0.5 to maintain backwards compatibility.
76  val += 0.5 * _tikhonov_coeff * param_norm_sqr;
77  }
78  return _objective_val + val;
79 }
80 
83 {
84  if (_ndof == 0)
85  mooseError(
86  "The number of parameters you have is zero and this shouldn't happen. Make sure you are "
87  "running your forward problem on \'INITIAL\' if you are using a reporter transfer to "
88  "supply "
89  "that information.");
90 
91  return _ndof;
92 }
93 
94 void
96 {
97  // Check that one and only one set of parameters are given.
98  // Set here because some derived reporters use a different method of
99  // determining numbers of dofs
100  if (!(isParamValid("num_values_name") ^ isParamValid("num_values")))
101  paramError("Need to supply one and only one of num_values_name or num_values.");
102 
105  else
106  _nvalues = getParam<std::vector<dof_id_type>>("num_values");
107 
108  _ndof = std::accumulate(_nvalues.begin(), _nvalues.end(), 0);
109 
110  // size checks
111  if (_parameter_names.size() != _nvalues.size())
112  paramError(
113  "num_parameters",
114  "There should be a number in \'num_parameters\' for each name in \'parameter_names\'.");
115 
116  for (const auto & param_id : make_range(_nparams))
117  {
118  _gradients[param_id]->resize(_nvalues[param_id]);
119 
120  std::vector<Real> ic(parseInputData("initial_condition", 0, param_id));
121  std::vector<Real> lb(
122  parseInputData("lower_bounds", std::numeric_limits<Real>::lowest(), param_id));
123  std::vector<Real> ub(
124  parseInputData("upper_bounds", std::numeric_limits<Real>::max(), param_id));
125 
126  _lower_bounds.insert(_lower_bounds.end(), lb.begin(), lb.end());
127  _upper_bounds.insert(_upper_bounds.end(), ub.begin(), ub.end());
128 
129  _parameters[param_id]->assign(ic.begin(), ic.end());
130  }
131 }
std::vector< dof_id_type > * _num_values_reporter
Real & _objective_val
Reporter that will hold the objective value.
const std::vector< ReporterValueName > & _parameter_names
Parameter names.
virtual Real computeObjective() override
Function to compute objective.
virtual dof_id_type getNumParams() const override
Function to get the total number of parameters.
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.
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.
bool isParamValid(const std::string &name) const
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.
GeneralOptimization(const InputParameters &parameters)
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("OptimizationApp", GeneralOptimization)
std::vector< std::vector< Real > * > _gradients
Gradient values declared as reporter data.
IntRange< T > make_range(T beg, T end)
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
const ReporterMode REPORTER_MODE_REPLICATED
static InputParameters validParams()
const Real _tikhonov_coeff
Tikhonov Coefficient for regularization.
virtual void setICsandBounds() override
Sets the initial conditions and bounds right before it is needed.
uint8_t dof_id_type
Optimization reporter that interfaces with TAO.