22 #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS) 51 LOG_SCOPE(
"objective()",
"TaoOptimizationSolver");
53 PetscErrorCode
ierr = 0;
56 libmesh_assert(objective);
82 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
87 libmesh_error_msg(
"Objective function not defined in __libmesh_tao_objective");
99 LOG_SCOPE(
"gradient()",
"TaoOptimizationSolver");
101 PetscErrorCode
ierr = 0;
133 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
138 libmesh_error_msg(
"Gradient function not defined in __libmesh_tao_gradient");
150 LOG_SCOPE(
"hessian()",
"TaoOptimizationSolver");
152 PetscErrorCode
ierr = 0;
181 hessian.attach_dof_map(sys.get_dof_map());
184 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
193 libmesh_error_msg(
"Hessian function not defined in __libmesh_tao_hessian");
207 LOG_SCOPE(
"equality_constraints()",
"TaoOptimizationSolver");
209 PetscErrorCode
ierr = 0;
238 eq_constraints.
zero();
241 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
246 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_equality_constraints");
248 eq_constraints.close();
259 LOG_SCOPE(
"equality_constraints_jacobian()",
"TaoOptimizationSolver");
261 PetscErrorCode
ierr = 0;
265 libmesh_assert(Jpre);
290 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
295 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_equality_constraints_jacobian");
308 LOG_SCOPE(
"inequality_constraints()",
"TaoOptimizationSolver");
310 PetscErrorCode
ierr = 0;
313 libmesh_assert(cineq);
339 ineq_constraints.
zero();
342 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
347 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_inequality_constraints");
349 ineq_constraints.close();
360 LOG_SCOPE(
"inequality_constraints_jacobian()",
"TaoOptimizationSolver");
362 PetscErrorCode
ierr = 0;
366 libmesh_assert(Jpre);
391 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
396 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_inequality_constraints_jacobian");
411 template <
typename T>
414 _reason(TAO_CONVERGED_USER)
420 template <
typename T>
428 template <
typename T>
435 PetscErrorCode
ierr=0;
437 ierr = TaoDestroy(&_tao);
438 LIBMESH_CHKERR(
ierr);
444 template <
typename T>
452 PetscErrorCode
ierr=0;
454 ierr = TaoCreate(this->comm().
get(),&_tao);
455 LIBMESH_CHKERR(
ierr);
459 template <
typename T>
462 LOG_SCOPE(
"solve()",
"TaoOptimizationSolver");
466 this->system().solution->zero();
468 PetscMatrix<T> * hessian = cast_ptr<PetscMatrix<T> *>(this->system().matrix);
470 PetscVector<T> * x = cast_ptr<PetscVector<T> *>(this->system().solution.get());
471 PetscVector<T> * ceq = cast_ptr<PetscVector<T> *>(this->system().C_eq.get());
472 PetscMatrix<T> * ceq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_eq_jac.get());
473 PetscVector<T> * cineq = cast_ptr<PetscVector<T> *>(this->system().C_ineq.get());
474 PetscMatrix<T> * cineq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_ineq_jac.get());
475 PetscVector<T> * lb = cast_ptr<PetscVector<T> *>(&this->system().get_vector(
"lower_bounds"));
476 PetscVector<T> * ub = cast_ptr<PetscVector<T> *>(&this->system().get_vector(
"upper_bounds"));
481 PetscErrorCode
ierr = 0;
488 ierr = TaoSetFromOptions(_tao);
489 LIBMESH_CHKERR(
ierr);
498 ierr = TaoSetTolerances(_tao,
499 #
if PETSC_RELEASE_LESS_THAN(3,7,0)
506 this->objective_function_relative_tolerance,
508 LIBMESH_CHKERR(
ierr);
512 ierr = TaoSetMaximumFunctionEvaluations(_tao, this->max_objective_function_evaluations);
513 LIBMESH_CHKERR(
ierr);
522 ierr = TaoSetInitialVector(_tao, x->vec());
523 LIBMESH_CHKERR(
ierr);
526 libmesh_assert( this->objective_object );
530 LIBMESH_CHKERR(
ierr);
532 if (this->gradient_object)
535 LIBMESH_CHKERR(
ierr);
538 if (this->hessian_object)
541 LIBMESH_CHKERR(
ierr);
544 if (this->lower_and_upper_bounds_object)
547 this->lower_and_upper_bounds_object->lower_and_upper_bounds(this->system());
549 ierr = TaoSetVariableBounds(_tao,
552 LIBMESH_CHKERR(
ierr);
555 if (this->equality_constraints_object)
558 LIBMESH_CHKERR(
ierr);
561 if (this->equality_constraints_jacobian_object)
563 ierr = TaoSetJacobianEqualityRoutine(_tao,
568 LIBMESH_CHKERR(
ierr);
572 if (this->inequality_constraints_object)
575 LIBMESH_CHKERR(
ierr);
579 if (this->inequality_constraints_jacobian_object)
581 ierr = TaoSetJacobianInequalityRoutine(_tao,
586 LIBMESH_CHKERR(
ierr);
590 ierr = TaoSetFromOptions(_tao);
591 LIBMESH_CHKERR(
ierr);
594 ierr = TaoSolve(_tao);
595 LIBMESH_CHKERR(
ierr);
601 this->system().get_dof_map().enforce_constraints_exactly(this->system());
604 ierr = TaoGetConvergedReason(_tao, &_reason);
605 LIBMESH_CHKERR(
ierr);
609 template <
typename T>
612 LOG_SCOPE(
"get_dual_variables()",
"TaoOptimizationSolver");
615 cast_ptr<PetscVector<T> *>(this->system().lambda_eq.get());
617 cast_ptr<PetscVector<T> *>(this->system().lambda_ineq.get());
619 Vec lambda_eq_petsc_vec = lambda_eq_petsc->
vec();
620 Vec lambda_ineq_petsc_vec = lambda_ineq_petsc->vec();
622 PetscErrorCode
ierr = 0;
623 ierr = TaoGetDualVariables(_tao,
624 &lambda_eq_petsc_vec,
625 &lambda_ineq_petsc_vec);
626 LIBMESH_CHKERR(
ierr);
630 template <
typename T>
633 libMesh::out <<
"Tao optimization solver convergence/divergence reason: " 634 << TaoConvergedReasons[this->get_converged_reason()] << std::endl;
639 template <
typename T>
642 PetscErrorCode
ierr=0;
646 ierr = TaoGetConvergedReason(_tao, &_reason);
647 LIBMESH_CHKERR(
ierr);
650 return static_cast<int>(_reason);
662 #endif // #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS) PetscErrorCode __libmesh_tao_objective(Tao, Vec x, PetscReal *objective, void *ctx)
OptimizationSystem::ComputeObjective * objective_object
PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao, Vec x, Mat J, Mat Jpre, void *ctx)
NumericVector interface to PETSc Vec.
virtual void gradient(const NumericVector< Number > &X, NumericVector< Number > &grad_f, sys_type &S)=0
virtual void zero() override
virtual void solve() override
virtual void hessian(const NumericVector< Number > &X, SparseMatrix< Number > &H_f, sys_type &S)=0
virtual void inequality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_ineq_jac, sys_type &S)=0
const Parallel::Communicator & comm() const
OptimizationSystem::ComputeGradient * gradient_object
virtual int get_converged_reason() override
PetscErrorCode __libmesh_tao_gradient(Tao, Vec x, Vec g, void *ctx)
virtual void print_converged_reason() override
OptimizationSystem::ComputeHessian * hessian_object
virtual void swap(NumericVector< T > &v) override
OptimizationSystem::ComputeInequalityConstraints * inequality_constraints_object
OptimizationSystem::ComputeEqualityConstraintsJacobian * equality_constraints_jacobian_object
TaoOptimizationSolver(sys_type &system)
void init(triangulateio &t)
PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao, Vec x, Mat J, Mat Jpre, void *ctx)
virtual void equality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_eq, sys_type &S)=0
const sys_type & system() const
void attach_dof_map(const DofMap &dof_map)
virtual Number objective(const NumericVector< Number > &X, sys_type &S)=0
virtual void equality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_eq_jac, sys_type &S)=0
SparseMatrix interface to PETSc Mat.
virtual void clear() override
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)
OStreamProxy out(std::cout)
virtual void init() override
PetscErrorCode __libmesh_tao_inequality_constraints(Tao, Vec x, Vec cineq, void *ctx)
virtual void get_dual_variables() override
virtual void inequality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_ineq, sys_type &S)=0
OptimizationSystem::ComputeInequalityConstraintsJacobian * inequality_constraints_jacobian_object
OptimizationSystem::ComputeEqualityConstraints * equality_constraints_object