https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | List of all members
OptimizeSolve Class Reference

solveObject to interface with Petsc Tao More...

#include <OptimizeSolve.h>

Inheritance diagram for OptimizeSolve:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 OptimizeSolve (Executioner &ex)
 
virtual bool solve () override
 
const OptimizationReporterBasegetOptimizationReporter () const
 
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. More...
 
virtual void initialSetup ()
 
virtual void setInnerSolve (SolveObject &solve)
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
PerfGraphperfGraph ()
 
bool isDefaultPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessor (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessorByName (const PostprocessorName &name) const
 
std::size_t coupledPostprocessors (const std::string &param_name) const
 
const PostprocessorName & getPostprocessorName (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 

Public Attributes

const ConsoleStream _console
 

Protected Member Functions

virtual PetscErrorCode variableBounds (Tao tao)
 Bounds routine. More...
 
virtual Real objectiveFunction ()
 Objective routine. More...
 
virtual void gradientFunction (libMesh::PetscVector< Number > &gradient)
 Gradient routine. More...
 
virtual PetscErrorCode applyHessian (libMesh::PetscVector< Number > &s, libMesh::PetscVector< Number > &Hs)
 Hessian application routine. More...
 
OptimizationReporterBasegetObjFunction ()
 function to get the objective reporter More...
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 
std::string timedSectionName (const std::string &section_name) const
 
virtual void addPostprocessorDependencyHelper (const PostprocessorName &) const
 

Protected Attributes

const libMesh::Parallel::Communicator _my_comm
 Communicator used for operations. More...
 
const ExecFlagEnum_solve_on
 List of execute flags for when to solve the system. More...
 
OptimizationReporterBase_obj_function = nullptr
 objective function defining objective, gradient, and hessian More...
 
Tao _tao
 Tao optimization object. More...
 
Executioner_executioner
 
FEProblemBase_problem
 
DisplacedProblem_displaced_problem
 
MooseMesh_mesh
 
MooseMesh_displaced_mesh
 
SystemBase_solver_sys
 
AuxiliarySystem_aux
 
SolveObject_inner_solve
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
MooseApp_pg_moose_app
 
const std::string _prefix
 
const Parallel::Communicator & _communicator
 

Private Types

enum  TaoSolverEnum {
  TaoSolverEnum::NEWTON_TRUST_REGION, TaoSolverEnum::BOUNDED_NEWTON_TRUST_REGION, TaoSolverEnum::BOUNDED_CONJUGATE_GRADIENT, TaoSolverEnum::NEWTON_LINE_SEARCH,
  TaoSolverEnum::BOUNDED_NEWTON_LINE_SEARCH, TaoSolverEnum::BOUNDED_QUASI_NEWTON_TRUST_REGION, TaoSolverEnum::NEWTON_TRUST_LINE, TaoSolverEnum::BOUNDED_NEWTON_TRUST_LINE,
  TaoSolverEnum::QUASI_NEWTON, TaoSolverEnum::BOUNDED_QUASI_NEWTON, TaoSolverEnum::NELDER_MEAD, TaoSolverEnum::BOUNDED_QUASI_NEWTON_LINE_SEARCH,
  TaoSolverEnum::ORTHANT_QUASI_NEWTON, TaoSolverEnum::GRADIENT_PROJECTION_CONJUGATE_GRADIENT, TaoSolverEnum::BUNDLE_RISK_MIN, TaoSolverEnum::AUGMENTED_LAGRANGIAN_MULTIPLIER_METHOD
}
 Enum of tao solver types. More...
 

Private Member Functions

PetscErrorCode taoSolve ()
 Here is where we call tao and solve. More...
 
void setTaoSolutionStatus (double f, int its, double gnorm, double cnorm, double xdiff)
 output optimization iteration solve data More...
 
PetscErrorCode taoALCreate ()
 Used for creating petsc structures when using the ALMM algorithm. More...
 
PetscErrorCode taoALDestroy ()
 Used for destroying petsc structures when using the ALMM algorithm. More...
 

Static Private Member Functions

static PetscErrorCode objectiveFunctionWrapper (Tao tao, Vec x, Real *objective, void *ctx)
 
static PetscErrorCode hessianFunctionWrapper (Tao tao, Vec x, Mat hessian, Mat pc, void *ctx)
 
static PetscErrorCode applyHessianWrapper (Mat H, Vec s, Vec Hs)
 
static PetscErrorCode objectiveAndGradientFunctionWrapper (Tao tao, Vec x, Real *objective, Vec gradient, void *ctx)
 
static PetscErrorCode variableBoundsWrapper (Tao, Vec xl, Vec xu, void *ctx)
 
static PetscErrorCode monitor (Tao tao, void *ctx)
 
static PetscErrorCode equalityFunctionWrapper (Tao tao, Vec x, Vec ce, void *ctx)
 
static PetscErrorCode equalityGradientFunctionWrapper (Tao tao, Vec x, Mat gradient_e, Mat gradient_epre, void *ctx)
 
static PetscErrorCode inequalityFunctionWrapper (Tao tao, Vec x, Vec ci, void *ctx)
 
static PetscErrorCode inequalityGradientFunctionWrapper (Tao tao, Vec x, Mat gradient_i, Mat gradient_ipre, void *ctx)
 

Private Attributes

bool _verbose
 control optimization executioner output More...
 
bool _output_opt_iters
 Use time step as the iteration counter for purposes of outputting. More...
 
std::vector< int_total_iterate_vec
 total solves per iteration More...
 
std::vector< int_obj_iterate_vec
 number of objective solves per iteration More...
 
std::vector< int_grad_iterate_vec
 gradient solves per iteration More...
 
std::vector< int_hess_iterate_vec
 Hessian solves per iteration. More...
 
std::vector< int_function_solve_vec
 total solves per iteration More...
 
std::vector< double > _f_vec
 objective value per iteration More...
 
std::vector< double > _gnorm_vec
 gradient norm per iteration More...
 
std::vector< double > _cnorm_vec
 infeasibility norm per iteration More...
 
std::vector< double > _xdiff_vec
 step length per iteration More...
 
enum OptimizeSolve::TaoSolverEnum _tao_solver_enum
 
dof_id_type _ndof
 Number of parameters being optimized. More...
 
std::unique_ptr< libMesh::PetscVector< Number > > _parameters
 Parameters (solution) given to TAO. More...
 
Mat _hessian
 Hessian (matrix) - usually a matrix-free representation. More...
 
Vec _ce
 Equality constraint vector. More...
 
Vec _ci
 Inequality constraint vector. More...
 
Mat _gradient_e
 Equality constraint gradient. More...
 
Mat _gradient_i
 Inequality constraint gradient. More...
 
int _obj_iterate = 0
 
int _grad_iterate = 0
 
int _hess_iterate = 0
 
Moose::PetscSupport::PetscOptions _petsc_options
 
SolverParams _solver_params
 

Detailed Description

solveObject to interface with Petsc Tao

Definition at line 28 of file OptimizeSolve.h.

Member Enumeration Documentation

◆ TaoSolverEnum

enum OptimizeSolve::TaoSolverEnum
strongprivate

Enum of tao solver types.

Enumerator
NEWTON_TRUST_REGION 
BOUNDED_NEWTON_TRUST_REGION 
BOUNDED_CONJUGATE_GRADIENT 
NEWTON_LINE_SEARCH 
BOUNDED_NEWTON_LINE_SEARCH 
BOUNDED_QUASI_NEWTON_TRUST_REGION 
NEWTON_TRUST_LINE 
BOUNDED_NEWTON_TRUST_LINE 
QUASI_NEWTON 
BOUNDED_QUASI_NEWTON 
NELDER_MEAD 
BOUNDED_QUASI_NEWTON_LINE_SEARCH 
ORTHANT_QUASI_NEWTON 
GRADIENT_PROJECTION_CONJUGATE_GRADIENT 
BUNDLE_RISK_MIN 
AUGMENTED_LAGRANGIAN_MULTIPLIER_METHOD 

Definition at line 154 of file OptimizeSolve.h.

154  {
155  NEWTON_TRUST_REGION,
156  BOUNDED_NEWTON_TRUST_REGION,
157  BOUNDED_CONJUGATE_GRADIENT,
158  NEWTON_LINE_SEARCH,
159  BOUNDED_NEWTON_LINE_SEARCH,
160  BOUNDED_QUASI_NEWTON_TRUST_REGION,
161  NEWTON_TRUST_LINE,
162  BOUNDED_NEWTON_TRUST_LINE,
163  QUASI_NEWTON,
164  BOUNDED_QUASI_NEWTON,
165  NELDER_MEAD,
166  BOUNDED_QUASI_NEWTON_LINE_SEARCH,
167  ORTHANT_QUASI_NEWTON,
168  GRADIENT_PROJECTION_CONJUGATE_GRADIENT,
169  BUNDLE_RISK_MIN,
170  AUGMENTED_LAGRANGIAN_MULTIPLIER_METHOD
enum OptimizeSolve::TaoSolverEnum _tao_solver_enum

Constructor & Destructor Documentation

◆ OptimizeSolve()

OptimizeSolve::OptimizeSolve ( Executioner ex)

Definition at line 43 of file OptimizeSolve.C.

44  : SolveObject(ex),
45  _my_comm(MPI_COMM_SELF),
46  _solve_on(getParam<ExecFlagEnum>("solve_on")),
47  _verbose(getParam<bool>("verbose")),
48  _output_opt_iters(getParam<bool>("output_optimization_iterations")),
49  _tao_solver_enum(getParam<MooseEnum>("tao_solver").getEnum<TaoSolverEnum>()),
51 {
52  if (libMesh::n_threads() > 1)
53  mooseError("OptimizeSolve does not currently support threaded execution");
54 
57  "moose", 27225, "Outputting for transient executioners has not been implemented.");
58 }
FEProblemBase & _problem
unsigned int n_threads()
const libMesh::Parallel::Communicator _my_comm
Communicator used for operations.
Definition: OptimizeSolve.h:75
bool _verbose
control optimization executioner output
Definition: OptimizeSolve.h:91
void mooseDocumentedError(const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
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
SolveObject(Executioner &ex)
void mooseError(Args &&... args) const
enum OptimizeSolve::TaoSolverEnum _tao_solver_enum
virtual bool isTransient() const override
bool _output_opt_iters
Use time step as the iteration counter for purposes of outputting.
Definition: OptimizeSolve.h:94

Member Function Documentation

◆ applyHessian()

PetscErrorCode OptimizeSolve::applyHessian ( libMesh::PetscVector< Number > &  s,
libMesh::PetscVector< Number > &  Hs 
)
protectedvirtual

Hessian application routine.

Definition at line 435 of file OptimizeSolve.C.

436 {
438  TIME_SECTION("applyHessian", 2, "Hessian forward/adjoint solve");
439  // What happens for material inversion when the Hessian
440  // is dependent on the parameters? Deal with it later???
441  // see notes on how this needs to change for Material inversion
442  if (_problem.hasMultiApps() &&
444  mooseError("Hessian based optimization algorithms require a sub-app with:\n"
445  " execute_on = HOMOGENEOUS_FORWARD");
447 
452  mooseError("Homogeneous forward solve multiapp failed!");
454  _inner_solve->solve();
455 
457 
462  mooseError("Adjoint solve multiapp failed!");
464  _inner_solve->solve();
465 
467  _hess_iterate++;
468  PetscFunctionReturn(PETSC_SUCCESS);
469 }
FEProblemBase & _problem
virtual void setMisfitToSimulatedValues()
Function to override misfit values with the simulated values from the matrix free hessian forward sol...
Moose::PetscSupport::PetscOptions _petsc_options
PetscFunctionBegin
const ExecFlagType EXEC_HOMOGENEOUS_FORWARD
virtual void execute(const ExecFlagType &exec_type)
SolveObject * _inner_solve
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
bool isValueSet(const std::string &value) const
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
const ExecFlagType EXEC_ADJOINT
virtual void computeGradient(libMesh::PetscVector< Number > &gradient) const
Function to compute gradient.
void restoreMultiApps(ExecFlagType type, bool force=false)
virtual void updateParameters(const libMesh::PetscVector< Number > &x)
Function to set parameters.
void mooseError(Args &&... args) const
bool execMultiApps(ExecFlagType type, bool auto_advance=true)
bool hasMultiApps() const
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
SolverParams _solver_params
virtual bool solve()=0

◆ applyHessianWrapper()

PetscErrorCode OptimizeSolve::applyHessianWrapper ( Mat  H,
Vec  s,
Vec  Hs 
)
staticprivate

Definition at line 370 of file OptimizeSolve.C.

Referenced by hessianFunctionWrapper().

371 {
372  void * ctx;
373 
375  LibmeshPetscCallQ(MatShellGetContext(H, &ctx));
376 
377  auto * solver = static_cast<OptimizeSolve *>(ctx);
378  libMesh::PetscVector<Number> sbar(s, solver->_my_comm);
379  libMesh::PetscVector<Number> Hsbar(Hs, solver->_my_comm);
380  return solver->applyHessian(sbar, Hsbar);
381 }
PetscFunctionBegin
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
void * ctx
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ equalityFunctionWrapper()

PetscErrorCode OptimizeSolve::equalityFunctionWrapper ( Tao  tao,
Vec  x,
Vec  ce,
void ctx 
)
staticprivate

Definition at line 492 of file OptimizeSolve.C.

Referenced by taoALCreate().

493 {
495  // grab the solver
496  auto * solver = static_cast<OptimizeSolve *>(ctx);
497  libMesh::PetscVector<Number> eq_con(ce, solver->_my_comm);
498  // use the OptimizationReporterBase class to actually compute equality constraints
499  OptimizationReporterBase * obj_func = solver->getObjFunction();
500  obj_func->computeEqualityConstraints(eq_con);
501  PetscFunctionReturn(PETSC_SUCCESS);
502 }
PetscFunctionBegin
Base class for optimization objects, implements routines for calculating misfit.
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
virtual void computeEqualityConstraints(libMesh::PetscVector< Number > &eqs_constraints) const
Function to compute the equality constraints.
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ equalityGradientFunctionWrapper()

PetscErrorCode OptimizeSolve::equalityGradientFunctionWrapper ( Tao  tao,
Vec  x,
Mat  gradient_e,
Mat  gradient_epre,
void ctx 
)
staticprivate

Definition at line 505 of file OptimizeSolve.C.

Referenced by taoALCreate().

507 {
509  // grab the solver
510  auto * solver = static_cast<OptimizeSolve *>(ctx);
511  libMesh::PetscMatrix<Number> grad_eq(gradient_e, solver->_my_comm);
512  // use the OptimizationReporterBase class to actually compute equality
513  // constraints gradient
514  OptimizationReporterBase * obj_func = solver->getObjFunction();
515  obj_func->computeEqualityGradient(grad_eq);
516  PetscFunctionReturn(PETSC_SUCCESS);
517 }
PetscFunctionBegin
virtual void computeEqualityGradient(libMesh::PetscMatrix< Number > &gradient) const
Function to compute the gradient of the equality constraints/ This is the last call of the equality c...
Base class for optimization objects, implements routines for calculating misfit.
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ getObjFunction()

OptimizationReporterBase* OptimizeSolve::getObjFunction ( )
inlineprotected

function to get the objective reporter

Definition at line 84 of file OptimizeSolve.h.

84 { return _obj_function; }
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81

◆ getOptimizationReporter()

const OptimizationReporterBase& OptimizeSolve::getOptimizationReporter ( ) const
inline

Definition at line 36 of file OptimizeSolve.h.

36 { return *_obj_function; }
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81

◆ getTaoSolutionStatus()

void OptimizeSolve::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.

Parameters
tot_iterstotal solves per iteration
gnormgradient norm per iteration
obj_itersnumber of objective solves per iteration
cnorminfeasibility norm per iteration
grad_itersgradient solves per iteration
xdiffstep length per iteration
hess_itersHessian solves per iteration
fobjective value per iteration
tot_solvestotal solves per iteration

Definition at line 233 of file OptimizeSolve.C.

Referenced by OptimizationInfo::execute().

242 {
243  const auto num = _total_iterate_vec.size();
244  tot_iters.resize(num);
245  obj_iters.resize(num);
246  grad_iters.resize(num);
247  hess_iters.resize(num);
248  tot_solves.resize(num);
249  f.resize(num);
250  gnorm.resize(num);
251  cnorm.resize(num);
252  xdiff.resize(num);
253 
254  for (const auto i : make_range(num))
255  {
256  tot_iters[i] = _total_iterate_vec[i];
257  obj_iters[i] = _obj_iterate_vec[i];
258  grad_iters[i] = _grad_iterate_vec[i];
259  hess_iters[i] = _hess_iterate_vec[i];
260  tot_solves[i] = _function_solve_vec[i];
261  f[i] = _f_vec[i];
262  gnorm[i] = _gnorm_vec[i];
263  cnorm[i] = _cnorm_vec[i];
264  xdiff[i] = _xdiff_vec[i];
265  }
266 }
std::vector< int > _hess_iterate_vec
Hessian solves per iteration.
Real f(Real x)
Test function for Brents method.
std::vector< int > _grad_iterate_vec
gradient solves per iteration
std::vector< double > _f_vec
objective value per iteration
std::vector< int > _function_solve_vec
total solves per iteration
std::vector< double > _cnorm_vec
infeasibility norm per iteration
IntRange< T > make_range(T beg, T end)
std::vector< double > _xdiff_vec
step length per iteration
std::vector< double > _gnorm_vec
gradient norm per iteration
std::vector< int > _total_iterate_vec
total solves per iteration
std::vector< int > _obj_iterate_vec
number of objective solves per iteration

◆ gradientFunction()

void OptimizeSolve::gradientFunction ( libMesh::PetscVector< Number > &  gradient)
protectedvirtual

Gradient routine.

Definition at line 417 of file OptimizeSolve.C.

418 {
419  TIME_SECTION("gradientFunction", 2, "Gradient adjoint solve");
421 
426  mooseError("Adjoint solve multiapp failed!");
428  _inner_solve->solve();
429 
430  _grad_iterate++;
431  _obj_function->computeGradient(gradient);
432 }
FEProblemBase & _problem
Moose::PetscSupport::PetscOptions _petsc_options
virtual void execute(const ExecFlagType &exec_type)
std::unique_ptr< libMesh::PetscVector< Number > > _parameters
Parameters (solution) given to TAO.
SolveObject * _inner_solve
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
bool isValueSet(const std::string &value) const
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
const ExecFlagType EXEC_ADJOINT
virtual void computeGradient(libMesh::PetscVector< Number > &gradient) const
Function to compute gradient.
void restoreMultiApps(ExecFlagType type, bool force=false)
virtual void updateParameters(const libMesh::PetscVector< Number > &x)
Function to set parameters.
void mooseError(Args &&... args) const
bool execMultiApps(ExecFlagType type, bool auto_advance=true)
SolverParams _solver_params
virtual bool solve()=0

◆ hessianFunctionWrapper()

PetscErrorCode OptimizeSolve::hessianFunctionWrapper ( Tao  tao,
Vec  x,
Mat  hessian,
Mat  pc,
void ctx 
)
staticprivate

Definition at line 358 of file OptimizeSolve.C.

Referenced by taoSolve().

360 {
362  // Define Hessian-vector multiplication routine
363  auto * solver = static_cast<OptimizeSolve *>(ctx);
364  LibmeshPetscCallQ(MatShellSetOperation(
365  solver->_hessian, MATOP_MULT, (void (*)(void))OptimizeSolve::applyHessianWrapper));
366  PetscFunctionReturn(PETSC_SUCCESS);
367 }
PetscFunctionBegin
static PetscErrorCode applyHessianWrapper(Mat H, Vec s, Vec Hs)
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
class infix_ostream_iterator if void
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ inequalityFunctionWrapper()

PetscErrorCode OptimizeSolve::inequalityFunctionWrapper ( Tao  tao,
Vec  x,
Vec  ci,
void ctx 
)
staticprivate

Definition at line 520 of file OptimizeSolve.C.

Referenced by taoALCreate().

521 {
523  // grab the solver
524  auto * solver = static_cast<OptimizeSolve *>(ctx);
525  libMesh::PetscVector<Number> ineq_con(ci, solver->_my_comm);
526  // use the OptimizationReporterBase class to actually compute equality constraints
527  OptimizationReporterBase * obj_func = solver->getObjFunction();
528  obj_func->computeInequalityConstraints(ineq_con);
529  PetscFunctionReturn(PETSC_SUCCESS);
530 }
virtual void computeInequalityConstraints(libMesh::PetscVector< Number > &ineqs_constraints) const
Function to compute the inequality constraints.
PetscFunctionBegin
Base class for optimization objects, implements routines for calculating misfit.
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ inequalityGradientFunctionWrapper()

PetscErrorCode OptimizeSolve::inequalityGradientFunctionWrapper ( Tao  tao,
Vec  x,
Mat  gradient_i,
Mat  gradient_ipre,
void ctx 
)
staticprivate

Definition at line 533 of file OptimizeSolve.C.

Referenced by taoALCreate().

535 {
537  // grab the solver
538  auto * solver = static_cast<OptimizeSolve *>(ctx);
539  libMesh::PetscMatrix<Number> grad_ineq(gradient_i, solver->_my_comm);
540  // use the OptimizationReporterBase class to actually compute equality
541  // constraints gradient
542  OptimizationReporterBase * obj_func = solver->getObjFunction();
543  obj_func->computeInequalityGradient(grad_ineq);
544  PetscFunctionReturn(PETSC_SUCCESS);
545 }
PetscFunctionBegin
Base class for optimization objects, implements routines for calculating misfit.
virtual void computeInequalityGradient(libMesh::PetscMatrix< Number > &gradient) const
Function to compute the gradient of the inequality constraints/ This is the last call of the inequali...
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ monitor()

PetscErrorCode OptimizeSolve::monitor ( Tao  tao,
void ctx 
)
staticprivate

Definition at line 311 of file OptimizeSolve.C.

Referenced by taoSolve().

312 {
313  TaoConvergedReason reason;
314  PetscInt its;
315  PetscReal f, gnorm, cnorm, xdiff;
316 
318  LibmeshPetscCallQ(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
319 
320  auto * solver = static_cast<OptimizeSolve *>(ctx);
321  solver->setTaoSolutionStatus((double)f, (int)its, (double)gnorm, (double)cnorm, (double)xdiff);
322 
323  PetscFunctionReturn(PETSC_SUCCESS);
324 }
PetscFunctionBegin
Real f(Real x)
Test function for Brents method.
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28
void setTaoSolutionStatus(double f, int its, double gnorm, double cnorm, double xdiff)
output optimization iteration solve data

◆ objectiveAndGradientFunctionWrapper()

PetscErrorCode OptimizeSolve::objectiveAndGradientFunctionWrapper ( Tao  tao,
Vec  x,
Real objective,
Vec  gradient,
void ctx 
)
staticprivate

Definition at line 341 of file OptimizeSolve.C.

Referenced by taoSolve().

343 {
345  auto * solver = static_cast<OptimizeSolve *>(ctx);
346 
347  libMesh::PetscVector<Number> param(x, solver->_my_comm);
348  solver->_parameters->swap(param);
349 
350  (*objective) = solver->objectiveFunction();
351  libMesh::PetscVector<Number> grad(gradient, solver->_my_comm);
352  solver->gradientFunction(grad);
353  solver->_parameters->swap(param);
354  PetscFunctionReturn(PETSC_SUCCESS);
355 }
PetscFunctionBegin
const std::vector< double > x
std::string grad(const std::string &var)
Definition: NS.h:91
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ objectiveFunction()

Real OptimizeSolve::objectiveFunction ( )
protectedvirtual

Objective routine.

Definition at line 394 of file OptimizeSolve.C.

395 {
396  TIME_SECTION("objectiveFunction", 2, "Objective forward solve");
398 
401 
404  {
405  // We do this so we can output for failed solves.
407  mooseError("Forward solve multiapp failed!");
408  }
410  _inner_solve->solve();
411 
412  _obj_iterate++;
414 }
FEProblemBase & _problem
Moose::PetscSupport::PetscOptions _petsc_options
const ExecFlagType EXEC_FORWARD
virtual void execute(const ExecFlagType &exec_type)
std::unique_ptr< libMesh::PetscVector< Number > > _parameters
Parameters (solution) given to TAO.
SolveObject * _inner_solve
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
virtual Real computeObjective()=0
Function to compute objective.
bool isValueSet(const std::string &value) const
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
void restoreMultiApps(ExecFlagType type, bool force=false)
virtual void updateParameters(const libMesh::PetscVector< Number > &x)
Function to set parameters.
void mooseError(Args &&... args) const
bool execMultiApps(ExecFlagType type, bool auto_advance=true)
SolverParams _solver_params
virtual bool solve()=0
virtual void outputStep(ExecFlagType type)

◆ objectiveFunctionWrapper()

PetscErrorCode OptimizeSolve::objectiveFunctionWrapper ( Tao  tao,
Vec  x,
Real objective,
void ctx 
)
staticprivate

Function wrappers for tao

Definition at line 327 of file OptimizeSolve.C.

Referenced by taoSolve().

328 {
330  auto * solver = static_cast<OptimizeSolve *>(ctx);
331 
332  libMesh::PetscVector<Number> param(x, solver->_my_comm);
333  solver->_parameters->swap(param);
334 
335  (*objective) = solver->objectiveFunction();
336  solver->_parameters->swap(param);
337  PetscFunctionReturn(PETSC_SUCCESS);
338 }
PetscFunctionBegin
const std::vector< double > x
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

◆ setTaoSolutionStatus()

void OptimizeSolve::setTaoSolutionStatus ( double  f,
int  its,
double  gnorm,
double  cnorm,
double  xdiff 
)
private

output optimization iteration solve data

Definition at line 269 of file OptimizeSolve.C.

Referenced by monitor().

270 {
271  // set data from TAO
272  _total_iterate_vec.push_back(its);
273  _f_vec.push_back(f);
274  _gnorm_vec.push_back(gnorm);
275  _cnorm_vec.push_back(cnorm);
276  _xdiff_vec.push_back(xdiff);
277  // set data we collect on this optimization iteration and then reset for next iteration
278  _obj_iterate_vec.push_back(_obj_iterate);
279  _grad_iterate_vec.push_back(_grad_iterate);
280  _hess_iterate_vec.push_back(_hess_iterate);
281  // count total number of FE solves
282  int solves = _obj_iterate + _grad_iterate + 2 * _hess_iterate;
283  _function_solve_vec.push_back(solves);
284  _obj_iterate = 0;
285  _grad_iterate = 0;
286  _hess_iterate = 0;
287 
288  // Pass down the iteration number if the subapp is of the Steady/SteadyAndAdjoint type.
289  // This enables exodus per-iteration output.
290  for (auto & sub_app : _app.getExecutioner()->feProblem().getMultiAppWarehouse().getObjects())
291  {
292  if (auto steady = dynamic_cast<Steady *>(sub_app->getExecutioner(0)))
293  steady->setIterationNumberOutput((unsigned int)its);
294  }
295 
296  // Output the converged iteration outputs
298 
299  // Increment timestep. In steady problems timestep = time for outputting.
300  // See Output.C
301  if (_output_opt_iters)
302  _problem.timeStep() += 1;
303 
304  // print verbose per iteration output
305  if (_verbose)
306  _console << "TAO SOLVER: iteration=" << its << "\tf=" << f << "\tgnorm=" << gnorm
307  << "\tcnorm=" << cnorm << "\txdiff=" << xdiff << std::endl;
308 }
FEProblemBase & _problem
FEProblemBase & feProblem()
const ExecFlagType EXEC_FORWARD
ExecuteMooseObjectWarehouse< MultiApp > & getMultiAppWarehouse()
std::vector< int > _hess_iterate_vec
Hessian solves per iteration.
bool _verbose
control optimization executioner output
Definition: OptimizeSolve.h:91
const std::vector< std::shared_ptr< MultiApp > > & getObjects(THREAD_ID tid=0) const
Real f(Real x)
Test function for Brents method.
Executioner * getExecutioner() const
virtual int & timeStep() const
std::vector< int > _grad_iterate_vec
gradient solves per iteration
std::vector< double > _f_vec
objective value per iteration
std::vector< int > _function_solve_vec
total solves per iteration
MooseApp & _app
std::vector< double > _cnorm_vec
infeasibility norm per iteration
std::vector< double > _xdiff_vec
step length per iteration
std::vector< double > _gnorm_vec
gradient norm per iteration
const ConsoleStream _console
std::vector< int > _total_iterate_vec
total solves per iteration
bool _output_opt_iters
Use time step as the iteration counter for purposes of outputting.
Definition: OptimizeSolve.h:94
std::vector< int > _obj_iterate_vec
number of objective solves per iteration
virtual void outputStep(ExecFlagType type)

◆ solve()

bool OptimizeSolve::solve ( )
overridevirtual

Implements SolveObject.

Definition at line 61 of file OptimizeSolve.C.

Referenced by Optimize::execute().

62 {
63  TIME_SECTION("optimizeSolve", 1, "Optimization Solve");
64  // Initial solve
66 
67  // Grab objective function
68  if (!_problem.hasUserObject("OptimizationReporter"))
69  mooseError("No OptimizationReporter object found.");
71 
72  // Initialize solution and matrix
74  _ndof = _parameters->size();
75 
76  // time step defaults 1, we want to start at 0 for first iteration to be
77  // consistent with TAO iterations.
79  _problem.timeStep() = 0;
80  bool solveInfo = (taoSolve() == 0);
81  return solveInfo;
82 }
FEProblemBase & _problem
T & getUserObject(const std::string &name, unsigned int tid=0) const
PetscErrorCode taoSolve()
Here is where we call tao and solve.
Definition: OptimizeSolve.C:85
bool hasUserObject(const std::string &name) const
Base class for optimization objects, implements routines for calculating misfit.
void setInitialCondition(libMesh::PetscVector< Number > &param)
Function to initialize petsc vectors from vpp data.
std::unique_ptr< libMesh::PetscVector< Number > > _parameters
Parameters (solution) given to TAO.
SolveObject * _inner_solve
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81
virtual int & timeStep() const
void mooseError(Args &&... args) const
bool _output_opt_iters
Use time step as the iteration counter for purposes of outputting.
Definition: OptimizeSolve.h:94
virtual bool solve()=0
dof_id_type _ndof
Number of parameters being optimized.

◆ taoALCreate()

PetscErrorCode OptimizeSolve::taoALCreate ( )
private

Used for creating petsc structures when using the ALMM algorithm.

Definition at line 548 of file OptimizeSolve.C.

Referenced by taoSolve().

549 {
552  {
553  // Create equality vector
554  LibmeshPetscCallQ(VecCreate(_my_comm.get(), &_ce));
557  LibmeshPetscCallQ(VecSetFromOptions(_ce));
558  LibmeshPetscCallQ(VecSetUp(_ce));
559 
560  // Set equality jacobian matrix
561  LibmeshPetscCallQ(MatCreate(_my_comm.get(), &_gradient_e));
562  LibmeshPetscCallQ(MatSetSizes(
564  LibmeshPetscCallQ(MatSetFromOptions(_gradient_e));
565  LibmeshPetscCallQ(MatSetUp(_gradient_e));
566 
567  // Set the Equality Constraints
568  LibmeshPetscCallQ(TaoSetEqualityConstraintsRoutine(_tao, _ce, equalityFunctionWrapper, this));
569 
570  // Set the Equality Constraints Jacobian
571  LibmeshPetscCallQ(TaoSetJacobianEqualityRoutine(
573  }
574 
576  {
577  // Create inequality vector
578  LibmeshPetscCallQ(VecCreate(_my_comm.get(), &_ci));
581  LibmeshPetscCallQ(VecSetFromOptions(_ci));
582  LibmeshPetscCallQ(VecSetUp(_ci));
583 
584  // Set inequality jacobian matrix
585  LibmeshPetscCallQ(MatCreate(_my_comm.get(), &_gradient_i));
586  LibmeshPetscCallQ(MatSetSizes(_gradient_i,
588  _ndof,
590  _ndof));
591  LibmeshPetscCallQ(MatSetFromOptions(_gradient_i));
592  LibmeshPetscCallQ(MatSetUp(_gradient_i));
593 
594  // Set the Inequality constraints
596  TaoSetInequalityConstraintsRoutine(_tao, _ci, inequalityFunctionWrapper, this));
597 
598  // Set the Inequality constraints Jacobian
599  LibmeshPetscCallQ(TaoSetJacobianInequalityRoutine(
601  }
602  PetscFunctionReturn(PETSC_SUCCESS);
603 }
Mat _gradient_e
Equality constraint gradient.
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)
dof_id_type getNumEqCons() const
Function to get the total number of equalities.
PetscFunctionBegin
Vec _ce
Equality constraint vector.
Vec _ci
Inequality constraint vector.
Mat _gradient_i
Inequality constraint gradient.
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81
static PetscErrorCode inequalityGradientFunctionWrapper(Tao tao, Vec x, Mat gradient_i, Mat gradient_ipre, void *ctx)
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
static PetscErrorCode equalityGradientFunctionWrapper(Tao tao, Vec x, Mat gradient_e, Mat gradient_epre, void *ctx)
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
static PetscErrorCode equalityFunctionWrapper(Tao tao, Vec x, Vec ce, void *ctx)
dof_id_type getNumInEqCons() const
Function to get the total number of inequalities.
dof_id_type _ndof
Number of parameters being optimized.

◆ taoALDestroy()

PetscErrorCode OptimizeSolve::taoALDestroy ( )
private

Used for destroying petsc structures when using the ALMM algorithm.

Definition at line 606 of file OptimizeSolve.C.

Referenced by taoSolve().

607 {
610  {
611  LibmeshPetscCallQ(VecDestroy(&_ce));
612  LibmeshPetscCallQ(MatDestroy(&_gradient_e));
613  }
615  {
616  LibmeshPetscCallQ(VecDestroy(&_ci));
617  LibmeshPetscCallQ(MatDestroy(&_gradient_i));
618  }
619 
620  PetscFunctionReturn(PETSC_SUCCESS);
621 }
Mat _gradient_e
Equality constraint gradient.
dof_id_type getNumEqCons() const
Function to get the total number of equalities.
PetscFunctionBegin
Vec _ce
Equality constraint vector.
Vec _ci
Inequality constraint vector.
Mat _gradient_i
Inequality constraint gradient.
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
dof_id_type getNumInEqCons() const
Function to get the total number of inequalities.

◆ taoSolve()

PetscErrorCode OptimizeSolve::taoSolve ( )
private

Here is where we call tao and solve.

Definition at line 85 of file OptimizeSolve.C.

Referenced by solve().

86 {
88  // Initialize tao object
89  LibmeshPetscCallQ(TaoCreate(_my_comm.get(), &_tao));
90 
91 #if PETSC_RELEASE_LESS_THAN(3, 21, 0)
92  LibmeshPetscCallQ(TaoSetMonitor(_tao, monitor, this, nullptr));
93 #else
94  LibmeshPetscCallQ(TaoMonitorSet(_tao, monitor, this, nullptr));
95 #endif
96 
97  switch (_tao_solver_enum)
98  {
100  LibmeshPetscCallQ(TaoSetType(_tao, TAONTR));
101  break;
103  LibmeshPetscCallQ(TaoSetType(_tao, TAOBNTR));
104  break;
106  LibmeshPetscCallQ(TaoSetType(_tao, TAOBNCG));
107  break;
109  LibmeshPetscCallQ(TaoSetType(_tao, TAONLS));
110  break;
112  LibmeshPetscCallQ(TaoSetType(_tao, TAOBNLS));
113  break;
115  LibmeshPetscCallQ(TaoSetType(_tao, TAOBQNKTR));
116  break;
118  LibmeshPetscCallQ(TaoSetType(_tao, TAONTL));
119  break;
121  LibmeshPetscCallQ(TaoSetType(_tao, TAOBNTL));
122  break;
124  LibmeshPetscCallQ(TaoSetType(_tao, TAOLMVM));
125  break;
127  LibmeshPetscCallQ(TaoSetType(_tao, TAOBLMVM));
128  break;
129 
131  LibmeshPetscCallQ(TaoSetType(_tao, TAONM));
132  break;
133 
135  LibmeshPetscCallQ(TaoSetType(_tao, TAOBQNLS));
136  break;
138  LibmeshPetscCallQ(TaoSetType(_tao, TAOOWLQN));
139  break;
141  LibmeshPetscCallQ(TaoSetType(_tao, TAOGPCG));
142  break;
144  LibmeshPetscCallQ(TaoSetType(_tao, TAOBMRM));
145  break;
147 #if !PETSC_VERSION_LESS_THAN(3, 15, 0)
148  LibmeshPetscCallQ(TaoSetType(_tao, TAOALMM));
149  // Need to cancel monitors for ALMM, if not there is a segfault at MOOSE destruction. Setup
150  // default constraint monitor.
151 #if PETSC_RELEASE_GREATER_EQUALS(3, 21, 0)
152  LibmeshPetscCallQ(TaoMonitorCancel(_tao));
153 #else
154  LibmeshPetscCallQ(TaoCancelMonitors(_tao));
155 #endif
156  LibmeshPetscCallQ(PetscOptionsSetValue(NULL, "-tao_cmonitor", NULL));
157  break;
158 #else
159  mooseError("ALMM is only compatible with PETSc versions above 3.14. ");
160 #endif
161 
162  default:
163  mooseError("Invalid Tao solve type");
164  }
165 
166  // Set objective and gradient functions
167 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
168  LibmeshPetscCallQ(TaoSetObjective(_tao, objectiveFunctionWrapper, this));
169 #else
170  LibmeshPetscCallQ(TaoSetObjectiveRoutine(_tao, objectiveFunctionWrapper, this));
171 #endif
172 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
174  TaoSetObjectiveAndGradient(_tao, NULL, objectiveAndGradientFunctionWrapper, this));
175 #else
177  TaoSetObjectiveAndGradientRoutine(_tao, objectiveAndGradientFunctionWrapper, this));
178 #endif
179 
180  // Set matrix-free version of the Hessian function
181  LibmeshPetscCallQ(MatCreateShell(_my_comm.get(), _ndof, _ndof, _ndof, _ndof, this, &_hessian));
182  // Link matrix-free Hessian to Tao
183 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
185 #else
186  LibmeshPetscCallQ(TaoSetHessianRoutine(_tao, _hessian, _hessian, hessianFunctionWrapper, this));
187 #endif
188 
189  // Set initial guess
190 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
191  LibmeshPetscCallQ(TaoSetSolution(_tao, _parameters->vec()));
192 #else
193  LibmeshPetscCallQ(TaoSetInitialVector(_tao, _parameters->vec()));
194 #endif
195 
196  // Set TAO petsc options
197  LibmeshPetscCallQ(TaoSetFromOptions(_tao));
198 
199  // save nonTAO PETSC options to reset before every call to execute()
201  // We only use a single system solve at this point
203 
204  // Set bounds for bounded optimization
205  LibmeshPetscCallQ(TaoSetVariableBoundsRoutine(_tao, variableBoundsWrapper, this));
206 
209 
210  // Backup multiapps so transient problems start with the same initial condition
214 
215  // Solve optimization
216  LibmeshPetscCallQ(TaoSolve(_tao));
217 
218  // Print solve statistics
219  if (getParam<bool>("verbose"))
220  LibmeshPetscCallQ(TaoView(_tao, PETSC_VIEWER_STDOUT_WORLD));
221 
222  LibmeshPetscCallQ(TaoDestroy(&_tao));
223 
224  LibmeshPetscCallQ(MatDestroy(&_hessian));
225 
228 
229  PetscFunctionReturn(PETSC_SUCCESS);
230 }
FEProblemBase & _problem
Moose::PetscSupport::PetscOptions & getPetscOptions()
Moose::PetscSupport::PetscOptions _petsc_options
const libMesh::Parallel::Communicator _my_comm
Communicator used for operations.
Definition: OptimizeSolve.h:75
Tao _tao
Tao optimization object.
Definition: OptimizeSolve.h:87
PetscFunctionBegin
const ExecFlagType EXEC_HOMOGENEOUS_FORWARD
const ExecFlagType EXEC_FORWARD
PetscErrorCode taoALDestroy()
Used for destroying petsc structures when using the ALMM algorithm.
static PetscErrorCode variableBoundsWrapper(Tao, Vec xl, Vec xu, void *ctx)
static PetscErrorCode objectiveAndGradientFunctionWrapper(Tao tao, Vec x, Real *objective, Vec gradient, void *ctx)
std::unique_ptr< libMesh::PetscVector< Number > > _parameters
Parameters (solution) given to TAO.
static PetscErrorCode monitor(Tao tao, void *ctx)
void backupMultiApps(ExecFlagType type)
const ExecFlagType EXEC_ADJOINT
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
Mat _hessian
Hessian (matrix) - usually a matrix-free representation.
void mooseError(Args &&... args) const
SolverParams & solverParams(unsigned int solver_sys_num=0)
enum OptimizeSolve::TaoSolverEnum _tao_solver_enum
static PetscErrorCode hessianFunctionWrapper(Tao tao, Vec x, Mat hessian, Mat pc, void *ctx)
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
SolverParams _solver_params
PetscErrorCode taoALCreate()
Used for creating petsc structures when using the ALMM algorithm.
static PetscErrorCode objectiveFunctionWrapper(Tao tao, Vec x, Real *objective, void *ctx)
dof_id_type _ndof
Number of parameters being optimized.

◆ validParams()

InputParameters OptimizeSolve::validParams ( )
static

Definition at line 18 of file OptimizeSolve.C.

Referenced by Optimize::validParams().

19 {
21  MooseEnum tao_solver_enum(
22  "taontr taobntr taobncg taonls taobnls taobqnktr taontl taobntl taolmvm "
23  "taoblmvm taonm taobqnls taoowlqn taogpcg taobmrm taoalmm");
25  "tao_solver", tao_solver_enum, "Tao solver to use for optimization.");
26  ExecFlagEnum exec_enum = ExecFlagEnum();
27  exec_enum.addAvailableFlags(EXEC_NONE,
34  params.addParam<ExecFlagEnum>(
35  "solve_on", exec_enum, "List of flags indicating when inner system solve should occur.");
36  params.addParam<bool>(
37  "output_optimization_iterations",
38  false,
39  "Use the time step as the current iteration for outputting optimization history.");
40  return params;
41 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const ExecFlagType EXEC_HOMOGENEOUS_FORWARD
const ExecFlagType EXEC_FORWARD
void addRequiredParam(const std::string &name, const std::string &doc_string)
InputParameters emptyInputParameters()
const ExecFlagType EXEC_ADJOINT

◆ variableBounds()

PetscErrorCode OptimizeSolve::variableBounds ( Tao  tao)
protectedvirtual

Bounds routine.

Definition at line 472 of file OptimizeSolve.C.

473 {
475  unsigned int sz = _obj_function->getNumParams();
476 
479 
480  // copy values from upper and lower bounds to xl and xu
481  for (const auto i : make_range(sz))
482  {
483  xl.set(i, _obj_function->getLowerBound(i));
484  xu.set(i, _obj_function->getUpperBound(i));
485  }
486  // set upper and lower bounds in tao solver
487  LibmeshPetscCallQ(TaoSetVariableBounds(tao, xl.vec(), xu.vec()));
488  PetscFunctionReturn(PETSC_SUCCESS);
489 }
Real getLowerBound(dof_id_type i) const
const libMesh::Parallel::Communicator _my_comm
Communicator used for operations.
Definition: OptimizeSolve.h:75
PetscFunctionBegin
virtual dof_id_type getNumParams() const
Function to get the total number of parameters.
OptimizationReporterBase * _obj_function
objective function defining objective, gradient, and hessian
Definition: OptimizeSolve.h:81
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
IntRange< T > make_range(T beg, T end)
Real getUpperBound(dof_id_type i) const
Upper and lower bounds for each parameter being controlled.
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)

◆ variableBoundsWrapper()

PetscErrorCode OptimizeSolve::variableBoundsWrapper ( Tao  tao,
Vec  xl,
Vec  xu,
void ctx 
)
staticprivate

Definition at line 384 of file OptimizeSolve.C.

Referenced by taoSolve().

385 {
387  auto * solver = static_cast<OptimizeSolve *>(ctx);
388 
389  LibmeshPetscCallQ(solver->variableBounds(tao));
390  PetscFunctionReturn(PETSC_SUCCESS);
391 }
PetscFunctionBegin
LibmeshPetscCallQ(DMShellGetContext(dm, &ctx))
void * ctx
PetscFunctionReturn(LIBMESH_PETSC_SUCCESS)
solveObject to interface with Petsc Tao
Definition: OptimizeSolve.h:28

Member Data Documentation

◆ _ce

Vec OptimizeSolve::_ce
private

Equality constraint vector.

Definition at line 183 of file OptimizeSolve.h.

Referenced by taoALCreate(), and taoALDestroy().

◆ _ci

Vec OptimizeSolve::_ci
private

Inequality constraint vector.

Definition at line 186 of file OptimizeSolve.h.

Referenced by taoALCreate(), and taoALDestroy().

◆ _cnorm_vec

std::vector<double> OptimizeSolve::_cnorm_vec
private

infeasibility norm per iteration

Definition at line 117 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _f_vec

std::vector<double> OptimizeSolve::_f_vec
private

objective value per iteration

Definition at line 113 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _function_solve_vec

std::vector<int> OptimizeSolve::_function_solve_vec
private

total solves per iteration

Definition at line 111 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _gnorm_vec

std::vector<double> OptimizeSolve::_gnorm_vec
private

gradient norm per iteration

Definition at line 115 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _grad_iterate

int OptimizeSolve::_grad_iterate = 0
private

Definition at line 99 of file OptimizeSolve.h.

Referenced by gradientFunction(), and setTaoSolutionStatus().

◆ _grad_iterate_vec

std::vector<int> OptimizeSolve::_grad_iterate_vec
private

gradient solves per iteration

Definition at line 107 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _gradient_e

Mat OptimizeSolve::_gradient_e
private

Equality constraint gradient.

Definition at line 189 of file OptimizeSolve.h.

Referenced by taoALCreate(), and taoALDestroy().

◆ _gradient_i

Mat OptimizeSolve::_gradient_i
private

Inequality constraint gradient.

Definition at line 192 of file OptimizeSolve.h.

Referenced by taoALCreate(), and taoALDestroy().

◆ _hess_iterate

int OptimizeSolve::_hess_iterate = 0
private

Definition at line 100 of file OptimizeSolve.h.

Referenced by applyHessian(), and setTaoSolutionStatus().

◆ _hess_iterate_vec

std::vector<int> OptimizeSolve::_hess_iterate_vec
private

Hessian solves per iteration.

Definition at line 109 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _hessian

Mat OptimizeSolve::_hessian
private

Hessian (matrix) - usually a matrix-free representation.

Definition at line 180 of file OptimizeSolve.h.

Referenced by taoSolve().

◆ _my_comm

const libMesh::Parallel::Communicator OptimizeSolve::_my_comm
protected

Communicator used for operations.

Definition at line 75 of file OptimizeSolve.h.

Referenced by taoALCreate(), taoSolve(), and variableBounds().

◆ _ndof

dof_id_type OptimizeSolve::_ndof
private

Number of parameters being optimized.

Definition at line 174 of file OptimizeSolve.h.

Referenced by solve(), taoALCreate(), and taoSolve().

◆ _obj_function

OptimizationReporterBase* OptimizeSolve::_obj_function = nullptr
protected

objective function defining objective, gradient, and hessian

Definition at line 81 of file OptimizeSolve.h.

Referenced by applyHessian(), getObjFunction(), getOptimizationReporter(), gradientFunction(), objectiveFunction(), solve(), taoALCreate(), taoALDestroy(), and variableBounds().

◆ _obj_iterate

int OptimizeSolve::_obj_iterate = 0
private

count individual solves for output

Definition at line 98 of file OptimizeSolve.h.

Referenced by objectiveFunction(), and setTaoSolutionStatus().

◆ _obj_iterate_vec

std::vector<int> OptimizeSolve::_obj_iterate_vec
private

number of objective solves per iteration

Definition at line 105 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _output_opt_iters

bool OptimizeSolve::_output_opt_iters
private

Use time step as the iteration counter for purposes of outputting.

Definition at line 94 of file OptimizeSolve.h.

Referenced by OptimizeSolve(), setTaoSolutionStatus(), and solve().

◆ _parameters

std::unique_ptr<libMesh::PetscVector<Number> > OptimizeSolve::_parameters
private

Parameters (solution) given to TAO.

Definition at line 177 of file OptimizeSolve.h.

Referenced by gradientFunction(), objectiveFunction(), solve(), and taoSolve().

◆ _petsc_options

Moose::PetscSupport::PetscOptions OptimizeSolve::_petsc_options
private

These are needed to reset the petsc options for the optimization solve using Moose::PetscSupport::petscSetOptions This only sets the finite difference options, the other optimizeSolve options are set-up in TAO using TaoSetFromOptions()

Definition at line 126 of file OptimizeSolve.h.

Referenced by applyHessian(), gradientFunction(), objectiveFunction(), and taoSolve().

◆ _solve_on

const ExecFlagEnum& OptimizeSolve::_solve_on
protected

List of execute flags for when to solve the system.

Definition at line 78 of file OptimizeSolve.h.

Referenced by applyHessian(), gradientFunction(), and objectiveFunction().

◆ _solver_params

SolverParams OptimizeSolve::_solver_params
private

Definition at line 127 of file OptimizeSolve.h.

Referenced by applyHessian(), gradientFunction(), objectiveFunction(), and taoSolve().

◆ _tao

Tao OptimizeSolve::_tao
protected

Tao optimization object.

Definition at line 87 of file OptimizeSolve.h.

Referenced by taoALCreate(), and taoSolve().

◆ _tao_solver_enum

enum OptimizeSolve::TaoSolverEnum OptimizeSolve::_tao_solver_enum
private

Referenced by taoSolve().

◆ _total_iterate_vec

std::vector<int> OptimizeSolve::_total_iterate_vec
private

total solves per iteration

Definition at line 103 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().

◆ _verbose

bool OptimizeSolve::_verbose
private

control optimization executioner output

Definition at line 91 of file OptimizeSolve.h.

Referenced by setTaoSolutionStatus().

◆ _xdiff_vec

std::vector<double> OptimizeSolve::_xdiff_vec
private

step length per iteration

Definition at line 119 of file OptimizeSolve.h.

Referenced by getTaoSolutionStatus(), and setTaoSolutionStatus().


The documentation for this class was generated from the following files: