tao_optimization_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_TAO_OPTIMIZATION_SOLVER_H
21 #define LIBMESH_TAO_OPTIMIZATION_SOLVER_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 // Petsc include files.
26 #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
27 
28 // Local includes
29 #include "libmesh/petsc_macro.h"
31 
32 // Include header for the Tao optimization library
33 #include <petsctao.h>
34 
35 namespace libMesh
36 {
37 
38 // Allow users access to these functions in case they want to reuse them. Users shouldn't
39 // need access to these most of the time as they are used internally by this object.
40 extern "C"
41 {
42  PetscErrorCode __libmesh_tao_objective (Tao tao, Vec x, PetscReal * objective, void * ctx);
43  PetscErrorCode __libmesh_tao_gradient(Tao tao, Vec x, Vec g, void * ctx);
44  PetscErrorCode __libmesh_tao_hessian(Tao tao, Vec x, Mat h, Mat pc, void * ctx);
45  PetscErrorCode __libmesh_tao_equality_constraints(Tao tao, Vec x, Vec ce, void * ctx);
46  PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx);
47  PetscErrorCode __libmesh_tao_inequality_constraints(Tao tao, Vec x, Vec cineq, void * ctx);
48  PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx);
49 }
50 
57 template <typename T>
59 {
60 public:
61 
66 
70  explicit
72 
77 
81  virtual void clear () override;
82 
86  virtual void init () override;
87 
91  Tao tao() { this->init(); return _tao; }
92 
96  virtual void solve () override;
97 
103  virtual void get_dual_variables() override;
104 
109  virtual void print_converged_reason() override;
110 
117  virtual int get_converged_reason() override;
118 
119 protected:
120 
124  Tao _tao;
125 
136  TaoConvergedReason _reason;
137 
138 private:
139 
140  friend PetscErrorCode __libmesh_tao_objective (Tao tao, Vec x, PetscReal * objective, void * ctx);
141  friend PetscErrorCode __libmesh_tao_gradient(Tao tao, Vec x, Vec g, void * ctx);
142  friend PetscErrorCode __libmesh_tao_hessian(Tao tao, Vec x, Mat h, Mat pc, void * ctx);
143  friend PetscErrorCode __libmesh_tao_equality_constraints(Tao tao, Vec x, Vec ce, void * ctx);
144  friend PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx);
145  friend PetscErrorCode __libmesh_tao_inequality_constraints(Tao tao, Vec x, Vec cineq, void * ctx);
146  friend PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx);
147 };
148 
149 
150 
151 } // namespace libMesh
152 
153 
154 #endif // #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
155 #endif // LIBMESH_TAO_OPTIMIZATION_SOLVER_H
PetscErrorCode __libmesh_tao_objective(Tao, Vec x, PetscReal *objective, void *ctx)
friend PetscErrorCode __libmesh_tao_inequality_constraints(Tao tao, Vec x, Vec cineq, void *ctx)
PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao, Vec x, Mat J, Mat Jpre, void *ctx)
friend PetscErrorCode __libmesh_tao_hessian(Tao tao, Vec x, Mat h, Mat pc, void *ctx)
friend PetscErrorCode __libmesh_tao_equality_constraints(Tao tao, Vec x, Vec ce, void *ctx)
virtual int get_converged_reason() override
PetscErrorCode __libmesh_tao_gradient(Tao, Vec x, Vec g, void *ctx)
virtual void print_converged_reason() override
friend PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void *ctx)
friend PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void *ctx)
friend PetscErrorCode __libmesh_tao_objective(Tao tao, Vec x, PetscReal *objective, void *ctx)
friend PetscErrorCode __libmesh_tao_gradient(Tao tao, Vec x, Vec g, void *ctx)
PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao, Vec x, Mat J, Mat Jpre, void *ctx)
const sys_type & system() const
PetscErrorCode __libmesh_tao_hessian(Tao, Vec x, Mat h, Mat pc, void *ctx)
PetscErrorCode __libmesh_tao_equality_constraints(Tao, Vec x, Vec ce, void *ctx)
PetscErrorCode __libmesh_tao_inequality_constraints(Tao, Vec x, Vec cineq, void *ctx)
virtual void get_dual_variables() override