https://mooseframework.inl.gov
ParameterMeshFunction.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 "ParameterMeshFunction.h"
11 
12 #include "AddVariableAction.h"
13 
14 registerMooseObject("OptimizationApp", ParameterMeshFunction);
15 
18 {
20  params.addClassDescription("Optimization function with parameters represented by a mesh and "
21  "finite-element shape functions.");
22 
23  params.addRequiredParam<FileName>("exodus_mesh", "File containing parameter mesh.");
24  params.addParam<MooseEnum>("family",
26  "Family of FE shape functions for parameter.");
27  params.addParam<MooseEnum>("order",
29  "Order of FE shape functions for parameter.");
30 
32  "parameter_name", "Reporter or VectorPostprocessor vector containing parameter values.");
33  params.addParam<ReporterName>("time_name",
34  "Name of vector-postprocessor or reporter vector containing time, "
35  "default assumes time independence.");
36  params.addParam<bool>(
37  "project_points", false, "Whether to find the closest point on parameter mesh.");
38  params.addParam<unsigned int>(
39  "kdtree_candidates",
40  5,
41  "Number of nearest node candidates to consider when projecting points to parameter mesh.");
42 
43  return params;
44 }
45 
47  : OptimizationFunction(parameters),
48  ReporterInterface(this),
49  _parameter_mesh(AddVariableAction::feType(parameters),
50  getParam<FileName>("exodus_mesh"),
51  {},
52  getParam<bool>("project_points"),
53  getParam<unsigned int>("kdtree_candidates")),
54  _values(getReporterValue<std::vector<Real>>("parameter_name")),
55  _coordt(isParamValid("time_name") ? getReporterValue<std::vector<Real>>("time_name")
56  : _empty_vec)
57 {
58 }
59 
60 Real
61 ParameterMeshFunction::value(Real t, const Point & p) const
62 {
63  checkSize();
64 
65  const auto ti = interpolateTime(t);
66  const dof_id_type offset0 = ti[0].first * _parameter_mesh.size();
67  const dof_id_type offset1 = ti[1].first * _parameter_mesh.size();
68 
69  std::vector<dof_id_type> dof_indices;
70  std::vector<Real> weights;
71 
72  _parameter_mesh.getIndexAndWeight(p, dof_indices, weights);
73 
74  Real val = 0;
75  for (const auto & i : index_range(dof_indices))
76  val += (_values[dof_indices[i] + offset0] * ti[0].second +
77  _values[dof_indices[i] + offset1] * ti[1].second) *
78  weights[i];
79  return val;
80 }
81 
83 ParameterMeshFunction::gradient(Real t, const Point & p) const
84 {
85  checkSize();
86 
87  const auto ti = interpolateTime(t);
88  const dof_id_type offset0 = ti[0].first * _parameter_mesh.size();
89  const dof_id_type offset1 = ti[1].first * _parameter_mesh.size();
90 
91  std::vector<dof_id_type> dof_indices;
92  std::vector<RealGradient> weights;
93  _parameter_mesh.getIndexAndWeight(p, dof_indices, weights);
94 
95  RealGradient val(0, 0, 0);
96  for (const auto & i : index_range(dof_indices))
97  val += (_values[dof_indices[i] + offset0] * ti[0].second +
98  _values[dof_indices[i] + offset1] * ti[1].second) *
99  weights[i];
100  return val;
101 }
102 
103 Real
104 ParameterMeshFunction::timeDerivative(Real t, const Point & p) const
105 {
106  checkSize();
107 
108  const auto ti = interpolateTime(t, true);
109  if (ti[0].first == ti[1].first)
110  return 0.0;
111  const dof_id_type offset0 = ti[0].first * _parameter_mesh.size();
112  const dof_id_type offset1 = ti[1].first * _parameter_mesh.size();
113 
114  std::vector<dof_id_type> dof_indices;
115  std::vector<Real> weights;
116  _parameter_mesh.getIndexAndWeight(p, dof_indices, weights);
117 
118  Real val = 0;
119  for (const auto & i : index_range(dof_indices))
120  val += (_values[dof_indices[i] + offset0] * ti[0].second +
121  _values[dof_indices[i] + offset1] * ti[1].second) *
122  weights[i];
123  return val;
124 }
125 
126 std::vector<Real>
127 ParameterMeshFunction::parameterGradient(Real t, const Point & p) const
128 {
129  const auto ti = interpolateTime(t);
130  const dof_id_type offset0 = ti[0].first * _parameter_mesh.size();
131  const dof_id_type offset1 = ti[1].first * _parameter_mesh.size();
132 
133  std::vector<dof_id_type> dof_indices;
134  std::vector<Real> weights;
135  _parameter_mesh.getIndexAndWeight(p, dof_indices, weights);
136 
138  if (!_coordt.empty())
139  sz *= _coordt.size();
140  std::vector<Real> pg(sz, 0.0);
141  for (const auto & i : index_range(dof_indices))
142  {
143  pg[dof_indices[i] + offset0] += weights[i] * ti[0].second;
144  pg[dof_indices[i] + offset1] += weights[i] * ti[1].second;
145  }
146  return pg;
147 }
148 
149 std::array<std::pair<std::size_t, Real>, 2>
150 ParameterMeshFunction::interpolateTime(Real t, bool derivative) const
151 {
152  std::array<std::pair<std::size_t, Real>, 2> ti;
153  if (_coordt.size() <= 1 || MooseUtils::absoluteFuzzyLessEqual(t, _coordt[0]))
154  {
155  ti[0] = {0, 0.5};
156  ti[1] = {0, 0.5};
157  }
159  {
160  ti[0] = {_coordt.size() - 1, 0.5};
161  ti[1] = {_coordt.size() - 1, 0.5};
162  }
163  else
164  for (std::size_t i = 1; i < _coordt.size(); ++i)
166  {
167  const Real dt = _coordt[i] - _coordt[i - 1];
168  if (derivative)
169  {
170  ti[0] = {i - 1, -1.0 / dt};
171  ti[1] = {i, 1.0 / dt};
172  }
173  else
174  {
175  ti[0] = {i - 1, (_coordt[i] - t) / dt};
176  ti[1] = {i, (t - _coordt[i - 1]) / dt};
177  }
178  break;
179  }
180 
181  return ti;
182 }
183 
184 void
186 {
188  if (!_coordt.empty())
189  sz *= _coordt.size();
190  if (sz != _values.size())
191  paramError("parameter_name",
192  "Size of parameter vector (",
193  _values.size(),
194  ") does not match number of degrees of freedom in mesh (",
195  sz,
196  ").");
197 }
std::array< std::pair< std::size_t, Real >, 2 > interpolateTime(Real t, bool derivative=false) const
This function is used to compute the weights for time interpolation.
registerMooseObject("OptimizationApp", ParameterMeshFunction)
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const std::vector< Real > & _values
values from reporter
ParameterMeshFunction(const InputParameters &parameters)
const ParameterMesh _parameter_mesh
Parameter mesh.
virtual RealGradient gradient(Real t, const Point &p) const override
void addRequiredParam(const std::string &name, const std::string &doc_string)
static MooseEnum getNonlinearVariableFamilies()
static MooseEnum getNonlinearVariableOrders()
virtual Real timeDerivative(Real t, const Point &p) const override
void checkSize() const
Used to make sure DoFs in &#39;_parameter_mesh&#39; matches number of values in &#39;_values&#39;.
void paramError(const std::string &param, Args... args) const
dof_id_type size() const
Definition: ParameterMesh.h:54
const std::vector< Real > & _coordt
Time coordinates from reporter.
virtual std::vector< Real > parameterGradient(Real t, const Point &p) const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for functions used in inverse optimization The parameterDerivative function is used in adj...
bool absoluteFuzzyLessEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
static InputParameters validParams()
bool absoluteFuzzyGreaterEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
void addClassDescription(const std::string &doc_string)
static InputParameters validParams()
void getIndexAndWeight(const Point &pt, std::vector< dof_id_type > &dof_indices, std::vector< Real > &weights) const
Interpolate parameters onto the computational mesh getIndexAndWeight is only used by ParameterMeshFun...
auto index_range(const T &sizable)
virtual Real value(Real t, const Point &p) const override
uint8_t dof_id_type