https://mooseframework.inl.gov
OptimizeSolve.h
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 #pragma once
11 
12 #include "SolveObject.h"
13 
14 #include "SolverParams.h"
15 #include "PetscSupport.h"
16 
17 #include "ExecFlagEnum.h"
18 #include <petsctao.h>
19 
20 #include "libmesh/petsc_vector.h"
21 #include "libmesh/petsc_matrix.h"
22 
24 
28 class OptimizeSolve : public SolveObject
29 {
30 public:
33 
34  virtual bool solve() override;
35 
37 
50  void getTaoSolutionStatus(std::vector<int> & tot_iters,
51  std::vector<double> & gnorm,
52  std::vector<int> & obj_iters,
53  std::vector<double> & cnorm,
54  std::vector<int> & grad_iters,
55  std::vector<double> & xdiff,
56  std::vector<int> & hess_iters,
57  std::vector<double> & f,
58  std::vector<int> & tot_solves) const;
59 
60 protected:
62  virtual PetscErrorCode variableBounds(Tao tao);
63 
65  virtual Real objectiveFunction();
66 
68  virtual void gradientFunction(libMesh::PetscVector<Number> & gradient);
69 
71  virtual PetscErrorCode applyHessian(libMesh::PetscVector<Number> & s,
73 
76 
79 
82 
85 
87  Tao _tao;
88 
89 private:
91  bool _verbose;
92 
95 
98  int _obj_iterate = 0;
99  int _grad_iterate = 0;
100  int _hess_iterate = 0;
103  std::vector<int> _total_iterate_vec;
105  std::vector<int> _obj_iterate_vec;
107  std::vector<int> _grad_iterate_vec;
109  std::vector<int> _hess_iterate_vec;
111  std::vector<int> _function_solve_vec;
113  std::vector<double> _f_vec;
115  std::vector<double> _gnorm_vec;
117  std::vector<double> _cnorm_vec;
119  std::vector<double> _xdiff_vec;
120 
129 
131  PetscErrorCode taoSolve();
132 
134  void setTaoSolutionStatus(double f, int its, double gnorm, double cnorm, double xdiff);
135 
138  static PetscErrorCode objectiveFunctionWrapper(Tao tao, Vec x, Real * objective, void * ctx);
139  static PetscErrorCode hessianFunctionWrapper(Tao tao, Vec x, Mat hessian, Mat pc, void * ctx);
140  static PetscErrorCode applyHessianWrapper(Mat H, Vec s, Vec Hs);
141  static PetscErrorCode
142  objectiveAndGradientFunctionWrapper(Tao tao, Vec x, Real * objective, Vec gradient, void * ctx);
143  static PetscErrorCode variableBoundsWrapper(Tao /*tao*/, Vec xl, Vec xu, void * ctx);
144  static PetscErrorCode monitor(Tao tao, void * ctx);
145  static PetscErrorCode equalityFunctionWrapper(Tao tao, Vec x, Vec ce, void * ctx);
146  static PetscErrorCode
147  equalityGradientFunctionWrapper(Tao tao, Vec x, Mat gradient_e, Mat gradient_epre, void * ctx);
148  static PetscErrorCode inequalityFunctionWrapper(Tao tao, Vec x, Vec ci, void * ctx);
149  static PetscErrorCode
150  inequalityGradientFunctionWrapper(Tao tao, Vec x, Mat gradient_i, Mat gradient_ipre, void * ctx);
152 
154  const enum class TaoSolverEnum {
163  QUASI_NEWTON,
165  NELDER_MEAD,
172 
175 
177  std::unique_ptr<libMesh::PetscVector<Number>> _parameters;
178 
180  Mat _hessian;
181 
183  Vec _ce;
184 
186  Vec _ci;
187 
190 
193 
195  PetscErrorCode taoALCreate();
197  PetscErrorCode taoALDestroy();
198 };
virtual bool solve() override
Definition: OptimizeSolve.C:61
Moose::PetscSupport::PetscOptions _petsc_options
static InputParameters validParams()
Definition: OptimizeSolve.C:18
Mat _gradient_e
Equality constraint gradient.
OptimizationReporterBase * getObjFunction()
function to get the objective reporter
Definition: OptimizeSolve.h:84
const libMesh::Parallel::Communicator _my_comm
Communicator used for operations.
Definition: OptimizeSolve.h:75
Tao _tao
Tao optimization object.
Definition: OptimizeSolve.h:87
static PetscErrorCode inequalityFunctionWrapper(Tao tao, Vec x, Vec ci, void *ctx)
PetscErrorCode taoSolve()
Here is where we call tao and solve.
Definition: OptimizeSolve.C:85
PetscErrorCode taoALDestroy()
Used for destroying petsc structures when using the ALMM algorithm.
Vec _ce
Equality constraint vector.
Vec _ci
Inequality constraint vector.
Mat _gradient_i
Inequality constraint gradient.
std::vector< int > _hess_iterate_vec
Hessian solves per iteration.
static PetscErrorCode applyHessianWrapper(Mat H, Vec s, Vec Hs)
static PetscErrorCode variableBoundsWrapper(Tao, Vec xl, Vec xu, void *ctx)
static PetscErrorCode objectiveAndGradientFunctionWrapper(Tao tao, Vec x, Real *objective, Vec gradient, void *ctx)
OptimizeSolve(Executioner &ex)
Definition: OptimizeSolve.C:43
Base class for optimization objects, implements routines for calculating misfit.
bool _verbose
control optimization executioner output
Definition: OptimizeSolve.h:91
void getTaoSolutionStatus(std::vector< int > &tot_iters, std::vector< double > &gnorm, std::vector< int > &obj_iters, std::vector< double > &cnorm, std::vector< int > &grad_iters, std::vector< double > &xdiff, std::vector< int > &hess_iters, std::vector< double > &f, std::vector< int > &tot_solves) const
Record tao TaoGetSolutionStatus data for output by a reporter.
const OptimizationReporterBase & getOptimizationReporter() const
Definition: OptimizeSolve.h:36
virtual PetscErrorCode applyHessian(libMesh::PetscVector< Number > &s, libMesh::PetscVector< Number > &Hs)
Hessian application routine.
const std::vector< double > x
Real f(Real x)
Test function for Brents method.
std::unique_ptr< libMesh::PetscVector< Number > > _parameters
Parameters (solution) given to TAO.
const ExecFlagEnum & _solve_on
List of execute flags for when to solve the system.
Definition: OptimizeSolve.h:78
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81
static PetscErrorCode monitor(Tao tao, void *ctx)
std::vector< int > _grad_iterate_vec
gradient solves per iteration
std::vector< double > _f_vec
objective value per iteration
static PetscErrorCode inequalityGradientFunctionWrapper(Tao tao, Vec x, Mat gradient_i, Mat gradient_ipre, void *ctx)
TaoSolverEnum
Enum of tao solver types.
std::vector< int > _function_solve_vec
total solves per iteration
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Mat _hessian
Hessian (matrix) - usually a matrix-free representation.
std::vector< double > _cnorm_vec
infeasibility norm per iteration
static PetscErrorCode equalityGradientFunctionWrapper(Tao tao, Vec x, Mat gradient_e, Mat gradient_epre, void *ctx)
std::vector< double > _xdiff_vec
step length per iteration
std::vector< double > _gnorm_vec
gradient norm per iteration
enum OptimizeSolve::TaoSolverEnum _tao_solver_enum
std::vector< int > _total_iterate_vec
total solves per iteration
static PetscErrorCode hessianFunctionWrapper(Tao tao, Vec x, Mat hessian, Mat pc, void *ctx)
SolverParams _solver_params
bool _output_opt_iters
Use time step as the iteration counter for purposes of outputting.
Definition: OptimizeSolve.h:94
virtual void gradientFunction(libMesh::PetscVector< Number > &gradient)
Gradient routine.
static PetscErrorCode equalityFunctionWrapper(Tao tao, Vec x, Vec ce, void *ctx)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28
PetscErrorCode taoALCreate()
Used for creating petsc structures when using the ALMM algorithm.
std::vector< int > _obj_iterate_vec
number of objective solves per iteration
void setTaoSolutionStatus(double f, int its, double gnorm, double cnorm, double xdiff)
output optimization iteration solve data
virtual Real objectiveFunction()
Objective routine.
virtual PetscErrorCode variableBounds(Tao tao)
Bounds routine.
static PetscErrorCode objectiveFunctionWrapper(Tao tao, Vec x, Real *objective, void *ctx)
dof_id_type _ndof
Number of parameters being optimized.
uint8_t dof_id_type