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 libmesh_error_msg_if(sys.solution->size() != n,
127 "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 libmesh_error_msg_if(sys.solution->size() != n,
"Error: Input vector x has different length than sys.solution!");
209 sys.solution->set(i, x[i]);
210 sys.solution->close();
213 sys.get_dof_map().enforce_constraints_exactly(sys);
232 for (
unsigned int i = 0; i < m; ++i)
233 result[i] = (*sys.C_ineq)(i);
249 sys.C_ineq_jac->close();
253 for (
const auto & dof_index : sys.ineq_constraint_jac_sparsity[i])
254 gradient[n*i+dof_index] = (*sys.C_ineq_jac)(i,dof_index);
257 libmesh_error_msg(
"Jacobian function not defined in __libmesh_nlopt_inequality_constraints");
262 libmesh_error_msg(
"Constraints function not defined in __libmesh_nlopt_inequality_constraints");
271 template <
typename T>
276 _result(NLOPT_SUCCESS),
278 _constraints_tolerance(1.e-8)
285 libmesh_error_msg_if(
sizeof(
dof_id_type) !=
sizeof(
unsigned int),
286 "The NloptOptimizationSolver should not be used with dof_id_type != unsigned int.");
291 template <
typename T>
299 template <
typename T>
312 template <
typename T>
321 std::string nlopt_algorithm_name =
"LD_SLSQP";
325 nlopt_algorithm_name);
328 _opt = nlopt_create(libmesh_map_find(_nlopt_algorithms, nlopt_algorithm_name),
329 this->system().solution->size());
335 template <
typename T>
338 LOG_SCOPE(
"solve()",
"NloptOptimizationSolver");
342 unsigned int nlopt_size = this->system().solution->size();
350 nlopt_set_min_objective(_opt,
353 libmesh_error_msg_if(ierr < 0,
"NLopt failed to set min objective: " << ierr);
356 if (this->lower_and_upper_bounds_object)
359 this->lower_and_upper_bounds_object->lower_and_upper_bounds(this->system());
361 std::vector<Real> nlopt_lb(nlopt_size);
362 std::vector<Real> nlopt_ub(nlopt_size);
363 for (
unsigned int i = 0; i < nlopt_size; ++i)
365 nlopt_lb[i] = this->system().get_vector(
"lower_bounds")(i);
366 nlopt_ub[i] = this->system().get_vector(
"upper_bounds")(i);
369 nlopt_set_lower_bounds(_opt, nlopt_lb.data());
370 nlopt_set_upper_bounds(_opt, nlopt_ub.data());
374 if (this->equality_constraints_object)
379 std::vector<double> equality_constraints_tolerances(this->system().C_eq->size(),
380 _constraints_tolerance);
386 nlopt_add_equality_mconstraint(_opt,
387 equality_constraints_tolerances.size(),
390 equality_constraints_tolerances.data());
392 libmesh_error_msg_if(ierr < 0,
"NLopt failed to add equality constraint: " << ierr);
396 if (this->inequality_constraints_object)
399 std::vector<double> inequality_constraints_tolerances(this->system().C_ineq->size(),
400 _constraints_tolerance);
402 nlopt_add_inequality_mconstraint(_opt,
403 inequality_constraints_tolerances.size(),
406 inequality_constraints_tolerances.data());
410 nlopt_set_ftol_rel(_opt, this->objective_function_relative_tolerance);
413 nlopt_set_maxeval(_opt, this->max_objective_function_evaluations);
416 this->_iteration_count = 0;
419 std::vector<Real> x(nlopt_size);
421 _result = nlopt_optimize(_opt, x.data(), &min_val);
429 << this->get_iteration_count()
435 template <
typename T>
438 libMesh::out <<
"NLopt optimization solver convergence/divergence reason: " 439 << this->get_converged_reason() << std::endl;
444 template <
typename T>
447 return static_cast<int>(_result);
459 #endif // #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS) ~NloptOptimizationSolver()
Destructor.
Abstract base class to be used to calculate the inequality constraints.
OptimizationSystem::ComputeObjective * objective_object
Object that computes the objective function f(X) at the input iterate X.
T command_line_next(std::string name, T default_value)
Use GetPot's search()/next() functions to get following arguments from the command line...
virtual void init() override
Initialize data structures if not done so already.
virtual void gradient(const NumericVector< Number > &X, NumericVector< Number > &grad_f, sys_type &S)=0
This function will be called to compute the gradient of the objective function, and must be implement...
This class provides an interface to the NLopt optimization solvers.
virtual void inequality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_ineq_jac, sys_type &S)=0
This function will be called to evaluate the Jacobian of C_ineq(X).
OptimizationSystem::ComputeGradient * gradient_object
Object that computes the gradient grad_f(X) of the objective function at the input iterate X...
The libMesh namespace provides an interface to certain functionality in the library.
OptimizationSystem::ComputeInequalityConstraints * inequality_constraints_object
Object that computes the inequality constraints vector C_ineq(X).
virtual void solve() override
Call the NLopt solver.
OptimizationSystem::ComputeEqualityConstraintsJacobian * equality_constraints_jacobian_object
Object that computes the Jacobian of C_eq(X).
Abstract base class to be used to calculate the equality constraints.
dof_id_type numeric_index_type
bool _is_initialized
Flag that tells if init() has been called.
This System subclass enables us to assemble an objective function, gradient, Hessian and bounds for o...
NloptOptimizationSolver(sys_type &system)
Constructor.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
virtual void equality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_eq, sys_type &S)=0
This function will be called to evaluate the equality constraints vector C_eq(X). ...
virtual void clear() override
Release all memory and clear data structures.
void __libmesh_nlopt_inequality_constraints(unsigned m, double *result, unsigned n, const double *x, double *gradient, void *data)
const sys_type & system() const
Abstract base class to be used to calculate the Jacobian of the equality constraints.
virtual Number objective(const NumericVector< Number > &X, sys_type &S)=0
This function will be called to compute the objective function to be minimized, and must be implement...
virtual void print_converged_reason() override
Prints a useful message about why the latest optimization solve con(di)verged.
virtual void equality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_eq_jac, sys_type &S)=0
This function will be called to evaluate the Jacobian of C_eq(X).
unsigned & get_iteration_count()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool verbose
Control how much is output from the OptimizationSolver as it's running.
This base class can be inherited from to provide interfaces to optimization solvers from different pa...
double __libmesh_nlopt_objective(unsigned n, const double *x, double *gradient, void *data)
bool initialized()
Checks that library initialization has been done.
bool on_command_line(std::string arg)
virtual int get_converged_reason() override
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
virtual void inequality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_ineq, sys_type &S)=0
This function will be called to evaluate the equality constraints vector C_ineq(X).
OptimizationSystem::ComputeInequalityConstraintsJacobian * inequality_constraints_jacobian_object
Object that computes the Jacobian of C_ineq(X).
OptimizationSystem::ComputeEqualityConstraints * equality_constraints_object
Object that computes the equality constraints vector C_eq(X).
Abstract base class to be used to calculate the Jacobian of the inequality constraints.
void __libmesh_nlopt_equality_constraints(unsigned m, double *result, unsigned n, const double *x, double *gradient, void *data)