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 : #include "MooseTypes.h"
11 : #include "OptimizationReporterBase.h"
12 : #include "OptUtils.h"
13 : #include "libmesh/petsc_vector.h"
14 :
15 : InputParameters
16 1308 : OptimizationReporterBase::validParams()
17 : {
18 1308 : InputParameters params = GeneralReporter::validParams();
19 1308 : params.registerBase("OptimizationReporterBase");
20 2616 : params.addRequiredParam<std::vector<ReporterValueName>>(
21 : "parameter_names", "List of parameter names, one for each group of parameters.");
22 3924 : params.addRangeCheckedParam<Real>(
23 2616 : "tikhonov_coeff", 0.0, "tikhonov_coeff >= 0", "Coefficient for Tikhonov Regularization.");
24 1308 : params.addParam<std::vector<ReporterValueName>>(
25 1308 : "equality_names", std::vector<ReporterValueName>(), "List of equality names.");
26 1308 : params.addParam<std::vector<ReporterValueName>>(
27 1308 : "inequality_names", std::vector<ReporterValueName>(), "List of inequality names.");
28 1308 : params.registerBase("OptimizationReporterBase");
29 1308 : return params;
30 0 : }
31 :
32 654 : OptimizationReporterBase::OptimizationReporterBase(const InputParameters & parameters)
33 : : GeneralReporter(parameters),
34 654 : _parameter_names(getParam<std::vector<ReporterValueName>>("parameter_names")),
35 654 : _nparams(_parameter_names.size()),
36 654 : _parameters(_nparams),
37 654 : _gradients(_nparams),
38 1308 : _tikhonov_coeff(getParam<Real>("tikhonov_coeff")),
39 1308 : _equality_names(&getParam<std::vector<ReporterValueName>>("equality_names")),
40 654 : _n_eq_cons(_equality_names->size()),
41 654 : _eq_constraints(_n_eq_cons),
42 654 : _eq_gradients(_n_eq_cons),
43 1308 : _inequality_names(&getParam<std::vector<ReporterValueName>>("inequality_names")),
44 654 : _n_ineq_cons(_inequality_names->size()),
45 654 : _ineq_constraints(_n_ineq_cons),
46 1308 : _ineq_gradients(_n_ineq_cons)
47 : {
48 1560 : for (const auto & i : make_range(_nparams))
49 : {
50 906 : _parameters[i] =
51 1812 : &declareValueByName<std::vector<Real>>(_parameter_names[i], REPORTER_MODE_REPLICATED);
52 3624 : _gradients[i] = &declareValueByName<std::vector<Real>>("grad_" + _parameter_names[i],
53 : REPORTER_MODE_REPLICATED);
54 : }
55 673 : for (const auto & i : make_range(_n_eq_cons))
56 : {
57 19 : _eq_constraints[i] =
58 38 : &declareValueByName<std::vector<Real>>(_equality_names->at(i), REPORTER_MODE_REPLICATED);
59 76 : _eq_gradients[i] = &declareValueByName<std::vector<Real>>("grad_" + _equality_names->at(i),
60 : REPORTER_MODE_REPLICATED);
61 : }
62 663 : for (const auto & i : make_range(_n_ineq_cons))
63 : {
64 9 : _ineq_constraints[i] =
65 18 : &declareValueByName<std::vector<Real>>(_inequality_names->at(i), REPORTER_MODE_REPLICATED);
66 36 : _ineq_gradients[i] = &declareValueByName<std::vector<Real>>("grad_" + _inequality_names->at(i),
67 : REPORTER_MODE_REPLICATED);
68 : }
69 654 : }
70 :
71 : void
72 4172 : OptimizationReporterBase::computeGradient(libMesh::PetscVector<Number> & gradient) const
73 : {
74 8576 : for (const auto & param_group_id : make_range(_nparams))
75 : {
76 4406 : if (_gradients[param_group_id]->size() != _nvalues[param_group_id])
77 2 : mooseError("The gradient for parameter ",
78 2 : _parameter_names[param_group_id],
79 : " has changed, expected ",
80 : _nvalues[param_group_id],
81 : " versus ",
82 2 : _gradients[param_group_id]->size(),
83 : ".");
84 :
85 4404 : if (_tikhonov_coeff > 0.0)
86 : {
87 468 : auto params = _parameters[param_group_id];
88 : auto grads = _gradients[param_group_id];
89 1701 : for (const auto & param_id : make_range(_nvalues[param_group_id]))
90 1233 : (*grads)[param_id] += (*params)[param_id] * _tikhonov_coeff;
91 : }
92 : }
93 :
94 4170 : OptUtils::copyReporterIntoPetscVector(_gradients, gradient);
95 4170 : }
96 :
97 : void
98 654 : OptimizationReporterBase::setInitialCondition(libMesh::PetscVector<Number> & x)
99 : {
100 654 : setICsandBounds();
101 630 : x.init(getNumParams());
102 630 : OptUtils::copyReporterIntoPetscVector(_parameters, x);
103 630 : }
104 :
105 : void
106 15604 : OptimizationReporterBase::updateParameters(const libMesh::PetscVector<Number> & x)
107 : {
108 15604 : OptUtils::copyPetscVectorIntoReporter(x, _parameters);
109 15604 : }
110 :
111 : Real
112 2832 : OptimizationReporterBase::getLowerBound(dof_id_type i) const
113 : {
114 2832 : return _lower_bounds[i];
115 : }
116 :
117 : Real
118 2832 : OptimizationReporterBase::getUpperBound(dof_id_type i) const
119 : {
120 2832 : return _upper_bounds[i];
121 : }
122 :
123 : std::vector<Real>
124 2842 : OptimizationReporterBase::parseInputData(std::string type,
125 : Real default_value,
126 : unsigned int param_id) const
127 : {
128 : // fill with default values
129 2842 : std::vector<Real> parsed_data_id(_nvalues[param_id], default_value);
130 2842 : if (isParamValid(type))
131 : {
132 1541 : std::vector<std::vector<Real>> parsed_data(getParam<std::vector<std::vector<Real>>>(type));
133 1541 : parsed_data_id.assign(parsed_data[param_id].begin(), parsed_data[param_id].end());
134 1541 : if (parsed_data.size() != _nvalues.size())
135 4 : paramError(type,
136 : "There must be a vector of ",
137 : type,
138 : " per parameter group. The ",
139 : type,
140 : " input format is std::vector<std::vector<Real>> so each vector should be "
141 : "seperated by \";\" even if it is a single value per group for a constant ",
142 : type,
143 : ".");
144 : // The case when the initial condition is constant for each parameter group
145 1537 : if (parsed_data[param_id].size() == 1)
146 859 : parsed_data_id.assign(_nvalues[param_id], parsed_data[param_id][0]);
147 678 : else if (parsed_data[param_id].size() != _nvalues[param_id])
148 2 : paramError(type,
149 : "When ",
150 : type,
151 : " are given in input file, there must either be a single value per parameter "
152 : "group or a value for every parameter in the group.");
153 1535 : }
154 :
155 2836 : return parsed_data_id;
156 0 : }
157 :
158 : void
159 0 : OptimizationReporterBase::computeEqualityConstraints(
160 : libMesh::PetscVector<Number> & eqs_constraints) const
161 : {
162 0 : OptUtils::copyReporterIntoPetscVector(_eq_constraints, eqs_constraints);
163 0 : }
164 :
165 : void
166 504 : OptimizationReporterBase::computeInequalityConstraints(
167 : libMesh::PetscVector<Number> & ineqs_constraints) const
168 : {
169 504 : OptUtils::copyReporterIntoPetscVector(_ineq_constraints, ineqs_constraints);
170 504 : }
171 :
172 : void
173 0 : OptimizationReporterBase::computeEqualityGradient(libMesh::PetscMatrix<Number> & jacobian) const
174 : {
175 0 : for (const auto & p : make_range(_n_eq_cons))
176 0 : if (_eq_gradients[p]->size() != _ndof)
177 0 : mooseError("The equality jacobian for parameter ",
178 0 : _parameter_names[p],
179 : " has changed, expected ",
180 0 : _ndof,
181 : " versus ",
182 0 : _eq_gradients[p]->size(),
183 : ".");
184 0 : OptUtils::copyReporterIntoPetscMatrix(_eq_gradients, jacobian);
185 0 : }
186 :
187 : void
188 504 : OptimizationReporterBase::computeInequalityGradient(libMesh::PetscMatrix<Number> & jacobian) const
189 : {
190 1008 : for (const auto & p : make_range(_n_ineq_cons))
191 504 : if (_ineq_gradients[p]->size() != _ndof)
192 0 : mooseError("The inequality jacobian for parameter ",
193 0 : _parameter_names[p],
194 : " has changed, expected ",
195 0 : _ndof,
196 : " versus ",
197 0 : _ineq_gradients[p]->size(),
198 : ".");
199 504 : OptUtils::copyReporterIntoPetscMatrix(_ineq_gradients, jacobian);
200 504 : }
|