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 : #pragma once 11 : 12 : #include "GeneralReporter.h" 13 : #include "libmesh/petsc_matrix.h" 14 : 15 : namespace libMesh 16 : { 17 : template <typename Number> 18 : class PetscVector; 19 : } 20 : 21 : /** 22 : * Base class for optimization objects, implements routines for calculating misfit. Derived classes 23 : * are responsible for parameter members and gradient computation. 24 : */ 25 : class OptimizationReporterBase : public GeneralReporter 26 : { 27 : public: 28 : static InputParameters validParams(); 29 : OptimizationReporterBase(const InputParameters & parameters); 30 : 31 12104 : virtual void initialize() override {} 32 12104 : virtual void execute() override {} 33 12104 : virtual void finalize() override {} 34 : /** 35 : * Function to compute objective. 36 : * This is the last function called in objective routine 37 : */ 38 : virtual Real computeObjective() = 0; 39 : 40 : /** 41 : * Function to compute gradient. 42 : * This is the last call of the gradient routine. 43 : */ 44 : virtual void computeGradient(libMesh::PetscVector<Number> & gradient) const; 45 : 46 : /** 47 : * Function to initialize petsc vectors from vpp data 48 : */ 49 : void setInitialCondition(libMesh::PetscVector<Number> & param); 50 : 51 : /** 52 : * Function to override misfit values with the simulated values from the matrix free hessian 53 : * forward solve. 54 : */ 55 117 : virtual void setMisfitToSimulatedValues(){}; 56 : 57 : /** 58 : * Upper and lower bounds for each parameter being controlled 59 : * 60 : * @param i Parameter index 61 : * @return The upper/lower bound for parameter i 62 : */ 63 : Real getUpperBound(dof_id_type i) const; 64 : Real getLowerBound(dof_id_type i) const; 65 : 66 : /** 67 : * Function to get the total number of parameters 68 : * @return total number of parameters 69 : */ 70 157 : virtual dof_id_type getNumParams() const { return _ndof; } 71 : 72 : /** 73 : * Function to compute the equality constraints. 74 : * This is the last call of the equality function routine. 75 : */ 76 : virtual void computeEqualityConstraints(libMesh::PetscVector<Number> & eqs_constraints) const; 77 : 78 : /** 79 : * Function to compute the inequality constraints. 80 : * This is the last call of the inequality function routine. 81 : */ 82 : virtual void computeInequalityConstraints(libMesh::PetscVector<Number> & ineqs_constraints) const; 83 : 84 : /** 85 : * Function to compute the gradient of the equality constraints/ 86 : * This is the last call of the equality constraint gradient routine. 87 : */ 88 : virtual void computeEqualityGradient(libMesh::PetscMatrix<Number> & gradient) const; 89 : 90 : /** 91 : * Function to compute the gradient of the inequality constraints/ 92 : * This is the last call of the inequality constraint gradient routine. 93 : */ 94 : virtual void computeInequalityGradient(libMesh::PetscMatrix<Number> & gradient) const; 95 : 96 : /** 97 : * Function to get the total number of equalities 98 : * @return total number of equality constraints 99 : */ 100 94 : dof_id_type getNumEqCons() const { return _n_eq_cons; } 101 : 102 : /** 103 : * Function to get the total number of inequalities 104 : * @return total number of inequalities constraints 105 : */ 106 74 : dof_id_type getNumInEqCons() const { return _n_ineq_cons; } 107 : 108 : protected: 109 : /** 110 : * Function to set parameters. 111 : * This is the first function called in objective/gradient/hessian routine 112 : */ 113 : virtual void updateParameters(const libMesh::PetscVector<Number> & x); 114 : 115 : /** 116 : * Function to to parse bounds and initial conditions from input file 117 : */ 118 : std::vector<Real> 119 : parseInputData(std::string type, Real default_value, unsigned int param_id) const; 120 : 121 : /// Sets the initial conditions and bounds right before it is needed. 122 0 : virtual void setICsandBounds(){}; 123 : 124 : /// Parameter names 125 : const std::vector<ReporterValueName> & _parameter_names; 126 : /// Number of parameter vectors 127 : const unsigned int _nparams; 128 : 129 : /// Parameter values declared as reporter data 130 : std::vector<std::vector<Real> *> _parameters; 131 : /// Gradient values declared as reporter data 132 : std::vector<std::vector<Real> *> _gradients; 133 : 134 : /// Tikhonov Coefficient for regularization 135 : const Real _tikhonov_coeff; 136 : 137 : /// Equality constraint names 138 : const std::vector<ReporterValueName> * _equality_names; 139 : /// Number of equality constraint names 140 : const unsigned int _n_eq_cons; 141 : /// Equality values declared as reporter data 142 : std::vector<std::vector<Real> *> _eq_constraints; 143 : /// Gradient values declared as reporter data 144 : std::vector<std::vector<Real> *> _eq_gradients; 145 : 146 : /// Inequality constraint names 147 : const std::vector<ReporterValueName> * _inequality_names; 148 : /// Number of inequality constraint names 149 : const unsigned int _n_ineq_cons; 150 : /// Inequality values declared as reporter data 151 : std::vector<std::vector<Real> *> _ineq_constraints; 152 : /// Gradient values declared as reporter data 153 : std::vector<std::vector<Real> *> _ineq_gradients; 154 : 155 : /// Bounds of the parameters 156 : std::vector<Real> _lower_bounds; 157 : std::vector<Real> _upper_bounds; 158 : 159 : /// Number of values for each parameter 160 : std::vector<dof_id_type> _nvalues; 161 : /// Total number of parameters 162 : dof_id_type _ndof; 163 : 164 : private: 165 : friend class OptimizeSolve; 166 : friend class OptimizationReporterTest; 167 : };