20 #include "libmesh/libmesh_common.h"
22 #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
28 #include "libmesh/libmesh_logging.h"
29 #include "libmesh/petsc_vector.h"
30 #include "libmesh/petsc_matrix.h"
31 #include "libmesh/dof_map.h"
32 #include "libmesh/tao_optimization_solver.h"
33 #include "libmesh/equation_systems.h"
52 LOG_SCOPE(
"objective()",
"TaoOptimizationSolver");
54 PetscErrorCode
ierr = 0;
83 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
88 libmesh_error_msg(
"Objective function not defined in __libmesh_tao_objective");
100 LOG_SCOPE(
"gradient()",
"TaoOptimizationSolver");
102 PetscErrorCode
ierr = 0;
134 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
139 libmesh_error_msg(
"Gradient function not defined in __libmesh_tao_gradient");
151 LOG_SCOPE(
"hessian()",
"TaoOptimizationSolver");
153 PetscErrorCode
ierr = 0;
182 hessian.attach_dof_map(sys.get_dof_map());
185 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
194 libmesh_error_msg(
"Hessian function not defined in __libmesh_tao_hessian");
208 LOG_SCOPE(
"equality_constraints()",
"TaoOptimizationSolver");
210 PetscErrorCode
ierr = 0;
239 eq_constraints.
zero();
242 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
247 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_equality_constraints");
249 eq_constraints.close();
260 LOG_SCOPE(
"equality_constraints_jacobian()",
"TaoOptimizationSolver");
262 PetscErrorCode
ierr = 0;
291 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
296 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_equality_constraints_jacobian");
309 LOG_SCOPE(
"inequality_constraints()",
"TaoOptimizationSolver");
311 PetscErrorCode
ierr = 0;
340 ineq_constraints.
zero();
343 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
348 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_inequality_constraints");
350 ineq_constraints.close();
361 LOG_SCOPE(
"inequality_constraints_jacobian()",
"TaoOptimizationSolver");
363 PetscErrorCode
ierr = 0;
392 sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
397 libmesh_error_msg(
"Constraints function not defined in __libmesh_tao_inequality_constraints_jacobian");
412 template <
typename T>
415 _reason(TAO_CONVERGED_USER)
421 template <
typename T>
429 template <
typename T>
436 PetscErrorCode
ierr=0;
438 ierr = TaoDestroy(&_tao);
439 LIBMESH_CHKERR(
ierr);
445 template <
typename T>
453 PetscErrorCode
ierr=0;
455 ierr = TaoCreate(this->comm().
get(),&_tao);
456 LIBMESH_CHKERR(
ierr);
460 template <
typename T>
463 LOG_SCOPE(
"solve()",
"TaoOptimizationSolver");
467 this->system().solution->zero();
469 PetscMatrix<T> * hessian = cast_ptr<PetscMatrix<T> *>(this->system().matrix);
471 PetscVector<T> * x = cast_ptr<PetscVector<T> *>(this->system().solution.get());
472 PetscVector<T> * ceq = cast_ptr<PetscVector<T> *>(this->system().C_eq.get());
473 PetscMatrix<T> * ceq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_eq_jac.get());
474 PetscVector<T> * cineq = cast_ptr<PetscVector<T> *>(this->system().C_ineq.get());
475 PetscMatrix<T> * cineq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_ineq_jac.get());
476 PetscVector<T> * lb = cast_ptr<PetscVector<T> *>(&this->system().get_vector(
"lower_bounds"));
477 PetscVector<T> * ub = cast_ptr<PetscVector<T> *>(&this->system().get_vector(
"upper_bounds"));
482 PetscErrorCode
ierr = 0;
489 ierr = TaoSetFromOptions(_tao);
490 LIBMESH_CHKERR(
ierr);
499 ierr = TaoSetTolerances(_tao,
500 #
if PETSC_RELEASE_LESS_THAN(3,7,0)
507 this->objective_function_relative_tolerance,
509 LIBMESH_CHKERR(
ierr);
513 ierr = TaoSetMaximumFunctionEvaluations(_tao, this->max_objective_function_evaluations);
514 LIBMESH_CHKERR(
ierr);
523 ierr = TaoSetInitialVector(_tao, x->vec());
524 LIBMESH_CHKERR(
ierr);
531 LIBMESH_CHKERR(
ierr);
533 if (this->gradient_object)
536 LIBMESH_CHKERR(
ierr);
539 if (this->hessian_object)
542 LIBMESH_CHKERR(
ierr);
545 if (this->lower_and_upper_bounds_object)
548 this->lower_and_upper_bounds_object->lower_and_upper_bounds(this->system());
550 ierr = TaoSetVariableBounds(_tao,
553 LIBMESH_CHKERR(
ierr);
556 if (this->equality_constraints_object)
559 LIBMESH_CHKERR(
ierr);
562 if (this->equality_constraints_jacobian_object)
564 ierr = TaoSetJacobianEqualityRoutine(_tao,
569 LIBMESH_CHKERR(
ierr);
573 if (this->inequality_constraints_object)
576 LIBMESH_CHKERR(
ierr);
580 if (this->inequality_constraints_jacobian_object)
582 ierr = TaoSetJacobianInequalityRoutine(_tao,
587 LIBMESH_CHKERR(
ierr);
591 ierr = TaoSetFromOptions(_tao);
592 LIBMESH_CHKERR(
ierr);
595 ierr = TaoSolve(_tao);
596 LIBMESH_CHKERR(
ierr);
602 this->system().get_dof_map().enforce_constraints_exactly(this->system());
605 ierr = TaoGetConvergedReason(_tao, &_reason);
606 LIBMESH_CHKERR(
ierr);
610 template <
typename T>
613 LOG_SCOPE(
"get_dual_variables()",
"TaoOptimizationSolver");
616 cast_ptr<PetscVector<T> *>(this->system().lambda_eq.get());
618 cast_ptr<PetscVector<T> *>(this->system().lambda_ineq.get());
620 Vec lambda_eq_petsc_vec = lambda_eq_petsc->
vec();
621 Vec lambda_ineq_petsc_vec = lambda_ineq_petsc->vec();
623 PetscErrorCode
ierr = 0;
624 ierr = TaoGetDualVariables(_tao,
625 &lambda_eq_petsc_vec,
626 &lambda_ineq_petsc_vec);
627 LIBMESH_CHKERR(
ierr);
631 template <
typename T>
634 libMesh::out <<
"Tao optimization solver convergence/divergence reason: "
635 << TaoConvergedReasons[this->get_converged_reason()] << std::endl;
640 template <
typename T>
643 PetscErrorCode
ierr=0;
647 ierr = TaoGetConvergedReason(_tao, &_reason);
648 LIBMESH_CHKERR(
ierr);
651 return static_cast<int>(_reason);
663 #endif // #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)