https://mooseframework.inl.gov
OptimizationReporterTest.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 
12 #include "libmesh/petsc_vector.h"
13 
14 registerMooseObject("OptimizationTestApp", OptimizationReporterTest);
15 
16 namespace
17 {
18 Real _namespaced_tol = 1e-6;
19 Real
20 differenceBetweenTwoVectors(std::vector<Real> a, std::vector<Real> b)
21 {
22  std::transform(a.cbegin(), a.cend(), b.cbegin(), b.begin(), std::minus<>{});
23  return std::abs(std::accumulate(b.begin(), b.end(), 0.0));
24 }
25 std::string
26 errorCheckMsg(const std::vector<Real> & optReporterBaseBounds,
27  const std::vector<Real> & optReporterTestBounds,
28  std::size_t ndof)
29 {
30  std::string out;
31  if (optReporterBaseBounds.size() != ndof || optReporterTestBounds.size() != ndof)
32  {
33  out += "There should be one bound per parameter read from the OptimizationDataBase object. ";
34  out += "\n Number of parameter values on OptimizationDataBase object: ";
35  out += std::to_string(ndof);
36  out += "\n Number of bounds values read from OptimizationDataBase object: ";
37  out += std::to_string(optReporterBaseBounds.size());
38  out += "\n Number of bounds values read from OptimizationReporterTest object: ";
39  out += std::to_string(optReporterTestBounds.size());
40  }
41  Real diff = differenceBetweenTwoVectors(optReporterBaseBounds, optReporterTestBounds);
42  if (diff > _namespaced_tol)
43  {
44  std::string out("\n OptimizationDataBase object bounds: ");
45  for (const auto & val : optReporterBaseBounds)
46  out += std::to_string(val) + " ";
47  out += "\n OptimizationReporterTest object bounds: ";
48  for (const auto & val : optReporterTestBounds)
49  out += std::to_string(val) + " ";
50  }
51 
52  return out;
53 }
54 }
55 
58 {
60  params.addClassDescription("Tests the requesting of Reporter values.");
61  params.addRequiredParam<std::vector<Real>>("values_to_set_parameters_to",
62  "Testing that OptimizationReporter values can be set");
63  params.addParam<std::vector<Real>>("expected_lower_bounds",
64  "Testing that OptimizationReporter lower bounds "
65  "can be set and used in the objective calculations");
66  params.addParam<std::vector<Real>>("expected_upper_bounds",
67  "Testing that OptimizationReporter upper bounds "
68  "can be set and used in the objective calculations");
69  return params;
70 }
71 
73  : GeneralUserObject(params),
74  _my_comm(MPI_COMM_SELF),
75  _optSolverParameters(std::make_unique<libMesh::PetscVector<Number>>(_my_comm))
76 {
77 }
78 
79 void
81 {
82  if (!_fe_problem.hasUserObject("OptimizationReporter"))
83  mooseError("No optimization reporter object found.");
85  // Set initial conditions
87  std::size_t ndof = _optSolverParameters->size();
88  std::vector<Real> valuesToSetOnOptRepParams(
89  getParam<std::vector<Real>>("values_to_set_parameters_to"));
90  if (ndof != valuesToSetOnOptRepParams.size())
91  mooseError("OptimizationReporter contains ",
92  ndof,
93  " and OptimizationReporterTest contains ",
94  valuesToSetOnOptRepParams.size(),
95  " parameters. They must be the same.");
96 
97  std::vector<Real> lower_bounds(ndof);
98  std::vector<Real> upper_bounds(ndof);
99  for (const auto & i : make_range(ndof))
100  {
101  lower_bounds[i] = _optReporter->getLowerBound(i);
102  upper_bounds[i] = _optReporter->getUpperBound(i);
103  }
104 
105  std::vector<Real> expectedLowerBounds = isParamValid("expected_lower_bounds")
106  ? getParam<std::vector<Real>>("expected_lower_bounds")
107  : std::vector<Real>(ndof, 0.0);
108  std::vector<Real> expectedUpperBounds = isParamValid("expected_upper_bounds")
109  ? getParam<std::vector<Real>>("expected_upper_bounds")
110  : std::vector<Real>(ndof, 0.0);
111 
112  std::string errorCheckLowerBounds(errorCheckMsg(expectedLowerBounds, lower_bounds, ndof));
113  if (!errorCheckLowerBounds.empty())
114  mooseError("Error in lower bounds: ", errorCheckLowerBounds);
115 
116  std::string errorCheckUpperBounds(errorCheckMsg(expectedUpperBounds, upper_bounds, ndof));
117  if (!errorCheckUpperBounds.empty())
118  mooseError("Error in upper bounds: ", errorCheckUpperBounds);
119 }
120 
121 void
123 {
124  // Tests that parameters can be changed correctly
125  std::vector<Real> valuesToSetOnOptRepParams(
126  getParam<std::vector<Real>>("values_to_set_parameters_to"));
127  size_t i = 0;
128  for (auto & val : valuesToSetOnOptRepParams)
129  _optSolverParameters->set(i++, val);
130 
132 }
static InputParameters validParams()
T & getUserObject(const std::string &name, unsigned int tid=0) const
static InputParameters validParams()
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
Real getLowerBound(dof_id_type i) const
bool hasUserObject(const std::string &name) const
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
void addRequiredParam(const std::string &name, const std::string &doc_string)
Base class for optimization objects, implements routines for calculating misfit.
bool isParamValid(const std::string &name) const
std::unique_ptr< libMesh::PetscVector< Number > > _optSolverParameters
void setInitialCondition(libMesh::PetscVector< Number > &param)
Function to initialize petsc vectors from vpp data.
OptimizationReporterBase * _optReporter
const T & getParam(const std::string &name) const
A UserObject that tests the requesting of Reporter values that are actually correct.
OptimizationReporterTest(const InputParameters &params)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("OptimizationTestApp", OptimizationReporterTest)
FEProblemBase & _fe_problem
OStreamProxy out
IntRange< T > make_range(T beg, T end)
virtual void updateParameters(const libMesh::PetscVector< Number > &x)
Function to set parameters.
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
Real getUpperBound(dof_id_type i) const
Upper and lower bounds for each parameter being controlled.
Real Number