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 "OptimizationReporter.h" 11 : 12 : #include "libmesh/int_range.h" 13 : 14 : registerMooseObjectDeprecated("OptimizationApp", OptimizationReporter, "12/31/2024 24:00"); 15 : 16 : InputParameters 17 184 : OptimizationReporter::validParams() 18 : { 19 184 : InputParameters params = OptimizationDataTempl<OptimizationReporterBase>::validParams(); 20 184 : params.addClassDescription("Computes objective function, gradient and contains reporters for " 21 : "communicating between optimizeSolve and subapps"); 22 368 : params.addRequiredParam<std::vector<dof_id_type>>( 23 : "num_values", 24 : "Number of parameter values associated with each parameter group in 'parameter_names'."); 25 368 : 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 368 : 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 368 : 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 184 : return params; 41 0 : } 42 : 43 92 : OptimizationReporter::OptimizationReporter(const InputParameters & parameters) 44 92 : : OptimizationDataTempl<OptimizationReporterBase>(parameters) 45 : { 46 92 : mooseDeprecated( 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 92 : } 51 : void 52 184 : OptimizationReporter::setICsandBounds() 53 : { 54 552 : _nvalues = getParam<std::vector<dof_id_type>>("num_values"); 55 184 : _ndof = std::accumulate(_nvalues.begin(), _nvalues.end(), 0); 56 : 57 : // size checks 58 184 : if (_parameter_names.size() != _nvalues.size()) 59 0 : paramError( 60 : "num_parameters", 61 : "There should be a number in \'num_parameters\' for each name in \'parameter_names\'."); 62 : 63 444 : for (const auto & param_id : make_range(_nparams)) 64 : { 65 260 : _gradients[param_id]->resize(_nvalues[param_id]); 66 : 67 260 : std::vector<Real> ic(parseInputData("initial_condition", 0, param_id)); 68 : std::vector<Real> lb( 69 260 : parseInputData("lower_bounds", std::numeric_limits<Real>::lowest(), param_id)); 70 : std::vector<Real> ub( 71 260 : parseInputData("upper_bounds", std::numeric_limits<Real>::max(), param_id)); 72 : 73 260 : _lower_bounds.insert(_lower_bounds.end(), lb.begin(), lb.end()); 74 260 : _upper_bounds.insert(_upper_bounds.end(), ub.begin(), ub.end()); 75 : 76 260 : _parameters[param_id]->assign(ic.begin(), ic.end()); 77 260 : } 78 184 : } 79 : 80 : Real 81 0 : OptimizationReporter::computeObjective() 82 : { 83 : // This will only be executed if measurement_values are available on the main app 84 0 : for (const auto i : index_range(_measurement_values)) 85 0 : _misfit_values[i] = _simulation_values[i] - _measurement_values[i]; 86 : 87 : Real val = 0.0; 88 0 : for (auto & misfit : _misfit_values) 89 0 : val += misfit * misfit; 90 : 91 0 : if (_tikhonov_coeff > 0.0) 92 : { 93 : Real param_norm_sqr = 0; 94 0 : for (const auto & data : _parameters) 95 0 : for (const auto & val : *data) 96 0 : param_norm_sqr += val * val; 97 : 98 0 : val += _tikhonov_coeff * param_norm_sqr; 99 : } 100 : 101 0 : return val * 0.5; 102 : } 103 : 104 : void 105 0 : OptimizationReporter::setMisfitToSimulatedValues() 106 : { 107 0 : _misfit_values = _simulation_values; 108 0 : }