20 #include "libmesh/libmesh_common.h"
22 #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
26 #include "libmesh/dof_map.h"
27 #include "libmesh/libmesh_logging.h"
28 #include "libmesh/numeric_vector.h"
29 #include "libmesh/nlopt_optimization_solver.h"
30 #include "libmesh/sparse_matrix.h"
31 #include "libmesh/int_range.h"
32 #include "libmesh/utility.h"
42 LOG_SCOPE(
"objective()",
"NloptOptimizationSolver");
53 sys.solution->set(i, x[i]);
56 sys.solution->close();
59 sys.get_dof_map().enforce_constraints_exactly(sys);
72 libmesh_error_msg(
"Objective function not defined in __libmesh_nlopt_objective");
86 std::vector<double> grad;
87 sys.rhs->localize_to_one(grad);
88 for (
unsigned int i = 0; i < n; ++i)
89 gradient[i] = grad[i];
92 libmesh_error_msg(
"Gradient function not defined in __libmesh_nlopt_objective");
114 LOG_SCOPE(
"equality_constraints()",
"NloptOptimizationSolver");
126 if (sys.solution->size() != n)
127 libmesh_error_msg(
"Error: Input vector x has different length than sys.solution!");
130 sys.solution->set(i, x[i]);
131 sys.solution->close();
134 sys.get_dof_map().enforce_constraints_exactly(sys);
153 for (
unsigned int i = 0; i < m; ++i)
154 result[i] = (*sys.C_eq)(i);
170 sys.C_eq_jac->close();
174 for (
const auto & dof_index : sys.eq_constraint_jac_sparsity[i])
175 gradient[n*i+dof_index] = (*sys.C_eq_jac)(i,dof_index);
178 libmesh_error_msg(
"Jacobian function not defined in __libmesh_nlopt_equality_constraints");
183 libmesh_error_msg(
"Constraints function not defined in __libmesh_nlopt_equality_constraints");
194 LOG_SCOPE(
"inequality_constraints()",
"NloptOptimizationSolver");
206 if (sys.solution->size() != n)
207 libmesh_error_msg(
"Error: Input vector x has different length than sys.solution!");
210 sys.solution->set(i, x[i]);
211 sys.solution->close();
214 sys.get_dof_map().enforce_constraints_exactly(sys);
233 for (
unsigned int i = 0; i < m; ++i)
234 result[i] = (*sys.C_ineq)(i);
250 sys.C_ineq_jac->close();
254 for (
const auto & dof_index : sys.ineq_constraint_jac_sparsity[i])
255 gradient[n*i+dof_index] = (*sys.C_ineq_jac)(i,dof_index);
258 libmesh_error_msg(
"Jacobian function not defined in __libmesh_nlopt_inequality_constraints");
263 libmesh_error_msg(
"Constraints function not defined in __libmesh_nlopt_inequality_constraints");
272 template <
typename T>
277 _result(NLOPT_SUCCESS),
279 _constraints_tolerance(1.e-8)
287 libmesh_error_msg(
"The NloptOptimizationSolver should not be used with dof_id_type != unsigned int.");
292 template <
typename T>
300 template <
typename T>
313 template <
typename T>
322 std::string nlopt_algorithm_name =
"LD_SLSQP";
326 nlopt_algorithm_name);
329 _opt = nlopt_create(libmesh_map_find(_nlopt_algorithms, nlopt_algorithm_name),
330 this->system().solution->size());
336 template <
typename T>
339 LOG_SCOPE(
"solve()",
"NloptOptimizationSolver");
343 unsigned int nlopt_size = this->system().solution->size();
351 nlopt_set_min_objective(_opt,
355 libmesh_error_msg(
"NLopt failed to set min objective: " <<
ierr);
358 if (this->lower_and_upper_bounds_object)
361 this->lower_and_upper_bounds_object->lower_and_upper_bounds(this->system());
363 std::vector<Real> nlopt_lb(nlopt_size);
364 std::vector<Real> nlopt_ub(nlopt_size);
365 for (
unsigned int i = 0; i < nlopt_size; ++i)
367 nlopt_lb[i] = this->system().get_vector(
"lower_bounds")(i);
368 nlopt_ub[i] = this->system().get_vector(
"upper_bounds")(i);
371 nlopt_set_lower_bounds(_opt, nlopt_lb.data());
372 nlopt_set_upper_bounds(_opt, nlopt_ub.data());
376 if (this->equality_constraints_object)
381 std::vector<double> equality_constraints_tolerances(this->system().C_eq->size(),
382 _constraints_tolerance);
388 nlopt_add_equality_mconstraint(_opt,
389 equality_constraints_tolerances.size(),
392 equality_constraints_tolerances.data());
395 libmesh_error_msg(
"NLopt failed to add equality constraint: " <<
ierr);
399 if (this->inequality_constraints_object)
402 std::vector<double> inequality_constraints_tolerances(this->system().C_ineq->size(),
403 _constraints_tolerance);
405 nlopt_add_inequality_mconstraint(_opt,
406 inequality_constraints_tolerances.size(),
409 inequality_constraints_tolerances.data());
413 nlopt_set_ftol_rel(_opt, this->objective_function_relative_tolerance);
416 nlopt_set_maxeval(_opt, this->max_objective_function_evaluations);
419 this->_iteration_count = 0;
422 std::vector<Real> x(nlopt_size);
424 _result = nlopt_optimize(_opt, x.data(), &min_val);
432 << this->get_iteration_count()
438 template <
typename T>
441 libMesh::out <<
"NLopt optimization solver convergence/divergence reason: "
442 << this->get_converged_reason() << std::endl;
447 template <
typename T>
450 return static_cast<int>(_result);
462 #endif // #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)