1 // The libMesh Finite Element Library. 2 // Copyright (C) 2002-2020 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" 30 #include "libmesh/optimization_solver.h" 31 32 // Include header for the Tao optimization library 33 #ifdef I 34 # define LIBMESH_SAW_I 35 #endif 36 #include <petsctao.h> 37 #ifndef LIBMESH_SAW_I 38 # undef I // Avoid complex.h contamination 39 #endif 40 41 namespace libMesh 42 { 43 44 // Allow users access to these functions in case they want to reuse them. Users shouldn't 45 // need access to these most of the time as they are used internally by this object. 46 extern "C" 47 { 48 PetscErrorCode __libmesh_tao_objective (Tao tao, Vec x, PetscReal * objective, void * ctx); 49 PetscErrorCode __libmesh_tao_gradient(Tao tao, Vec x, Vec g, void * ctx); 50 PetscErrorCode __libmesh_tao_hessian(Tao tao, Vec x, Mat h, Mat pc, void * ctx); 51 PetscErrorCode __libmesh_tao_equality_constraints(Tao tao, Vec x, Vec ce, void * ctx); 52 PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx); 53 PetscErrorCode __libmesh_tao_inequality_constraints(Tao tao, Vec x, Vec cineq, void * ctx); 54 PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx); 55 } 56 57 /** 58 * This class provides an interface to the Tao optimization solvers. 59 * 60 * \author David Knezevic 61 * \date 2015 62 */ 63 template <typename T> 64 class TaoOptimizationSolver : public OptimizationSolver<T> 65 { 66 public: 67 68 /** 69 * The type of system that we use in conjunction with this solver. 70 */ 71 typedef OptimizationSystem sys_type; 72 73 /** 74 * Constructor. Initializes Tao data structures. 75 */ 76 explicit 77 TaoOptimizationSolver (sys_type & system); 78 79 /** 80 * Destructor. 81 */ 82 ~TaoOptimizationSolver (); 83 84 /** 85 * Release all memory and clear data structures. 86 */ 87 virtual void clear () override; 88 89 /** 90 * Initialize data structures if not done so already. 91 */ 92 virtual void init () override; 93 94 /** 95 * \returns The raw PETSc Tao context pointer. 96 */ tao()97 Tao tao() { this->init(); return _tao; } 98 99 /** 100 * Call the Tao solver. 101 */ 102 virtual void solve () override; 103 104 /** 105 * Get the current values of dual variables associated with 106 * inequality and equality constraints. The variables will 107 * be stored in _system.lambda_eq and _system.lambda_ineq. 108 */ 109 virtual void get_dual_variables() override; 110 111 /** 112 * Prints a useful message about why the latest optimization solve 113 * con(di)verged. 114 */ 115 virtual void print_converged_reason() override; 116 117 /** 118 * \returns The currently-available (or most recently obtained, if 119 * the Tao object has been destroyed) convergence reason. 120 * 121 * Refer to Tao docs for the meaning of different TaoConvergedReason. 122 */ 123 virtual int get_converged_reason() override; 124 125 protected: 126 127 /** 128 * Optimization solver context 129 */ 130 Tao _tao; 131 132 /** 133 * Store the reason for Tao convergence/divergence for use even 134 * after \p _tao has been cleared. 135 * 136 * \note \p print_converged_reason() will always \e try to get the 137 * current reason with TaoGetConvergedReason(), but if the Tao 138 * object has already been cleared, it will fall back on this stored 139 * value. This value is therefore necessarily \e not cleared by the 140 * clear() function. 141 */ 142 TaoConvergedReason _reason; 143 144 private: 145 146 friend PetscErrorCode __libmesh_tao_objective (Tao tao, Vec x, PetscReal * objective, void * ctx); 147 friend PetscErrorCode __libmesh_tao_gradient(Tao tao, Vec x, Vec g, void * ctx); 148 friend PetscErrorCode __libmesh_tao_hessian(Tao tao, Vec x, Mat h, Mat pc, void * ctx); 149 friend PetscErrorCode __libmesh_tao_equality_constraints(Tao tao, Vec x, Vec ce, void * ctx); 150 friend PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx); 151 friend PetscErrorCode __libmesh_tao_inequality_constraints(Tao tao, Vec x, Vec cineq, void * ctx); 152 friend PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void * ctx); 153 }; 154 155 156 157 } // namespace libMesh 158 159 160 #endif // #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS) 161 #endif // LIBMESH_TAO_OPTIMIZATION_SOLVER_H 162