nlopt_optimization_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 #include "libmesh/libmesh_common.h"
21 
22 #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
23 
24 
25 // C++ includes
26 
27 // Local Includes
28 #include "libmesh/dof_map.h"
29 #include "libmesh/numeric_vector.h"
31 #include "libmesh/sparse_matrix.h"
32 #include "libmesh/int_range.h"
33 
34 namespace libMesh
35 {
36 
37 double __libmesh_nlopt_objective(unsigned n,
38  const double * x,
39  double * gradient,
40  void * data)
41 {
42  LOG_SCOPE("objective()", "NloptOptimizationSolver");
43 
44  // ctx should be a pointer to the solver (it was passed in as void *)
46  static_cast<NloptOptimizationSolver<Number> *> (data);
47 
48  OptimizationSystem & sys = solver->system();
49 
50  // We'll use current_local_solution below, so let's ensure that it's consistent
51  // with the vector x that was passed in.
52  for (auto i : index_range(*sys.solution))
53  sys.solution->set(i, x[i]);
54 
55  // Make sure the solution vector is parallel-consistent
56  sys.solution->close();
57 
58  // Impose constraints on X
59  sys.get_dof_map().enforce_constraints_exactly(sys);
60 
61  // Update sys.current_local_solution based on X
62  sys.update();
63 
64  Real objective;
65  if (solver->objective_object != nullptr)
66  {
67  objective =
68  solver->objective_object->objective(*(sys.current_local_solution), sys);
69  }
70  else
71  {
72  libmesh_error_msg("Objective function not defined in __libmesh_nlopt_objective");
73  }
74 
75  // If the gradient has been requested, fill it in
76  if (gradient)
77  {
78  if (solver->gradient_object != nullptr)
79  {
80  solver->gradient_object->gradient(*(sys.current_local_solution), *(sys.rhs), sys);
81 
82  // we've filled up sys.rhs with the gradient data, now copy it
83  // to the nlopt data structure
84  libmesh_assert(sys.rhs->size() == n);
85 
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];
90  }
91  else
92  libmesh_error_msg("Gradient function not defined in __libmesh_nlopt_objective");
93  }
94 
95  // Increment the iteration count.
96  solver->get_iteration_count()++;
97 
98  // Possibly print the current value of the objective function
99  if (solver->verbose)
100  libMesh::out << objective << std::endl;
101 
102  return objective;
103 }
104 
105 
106 
108  double * result,
109  unsigned n,
110  const double * x,
111  double * gradient,
112  void * data)
113 {
114  LOG_SCOPE("equality_constraints()", "NloptOptimizationSolver");
115 
116  libmesh_assert(data);
117 
118  // data should be a pointer to the solver (it was passed in as void *)
120  static_cast<NloptOptimizationSolver<Number> *> (data);
121 
122  OptimizationSystem & sys = solver->system();
123 
124  // We'll use current_local_solution below, so let's ensure that it's consistent
125  // with the vector x that was passed in.
126  if (sys.solution->size() != n)
127  libmesh_error_msg("Error: Input vector x has different length than sys.solution!");
128 
129  for (auto i : index_range(*sys.solution))
130  sys.solution->set(i, x[i]);
131  sys.solution->close();
132 
133  // Impose constraints on the solution vector
134  sys.get_dof_map().enforce_constraints_exactly(sys);
135 
136  // Update sys.current_local_solution based on the solution vector
137  sys.update();
138 
139  // Call the user's equality constraints function if there is one.
141  if (eco)
142  {
143  eco->equality_constraints(*sys.current_local_solution,
144  *sys.C_eq,
145  sys);
146 
147  sys.C_eq->close();
148 
149  // Copy the values out of eq_constraints into 'result'.
150  // TODO: Even better would be if we could use 'result' directly
151  // as the storage of eq_constraints. Perhaps a serial-only
152  // NumericVector variant which supports this option?
153  for (unsigned int i = 0; i < m; ++i)
154  result[i] = (*sys.C_eq)(i);
155 
156  // If gradient != nullptr, then the Jacobian matrix of the equality
157  // constraints has been requested. The incoming 'gradient'
158  // array is of length m*n and d(c_i)/d(x_j) = gradient[n*i+j].
159  if (gradient)
160  {
163 
164  if (eco_jac)
165  {
166  eco_jac->equality_constraints_jacobian(*sys.current_local_solution,
167  *sys.C_eq_jac,
168  sys);
169 
170  sys.C_eq_jac->close();
171 
172  // copy the Jacobian data to the gradient array
173  for (numeric_index_type i=0; i<m; i++)
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);
176  }
177  else
178  libmesh_error_msg("Jacobian function not defined in __libmesh_nlopt_equality_constraints");
179  }
180 
181  }
182  else
183  libmesh_error_msg("Constraints function not defined in __libmesh_nlopt_equality_constraints");
184 }
185 
186 
188  double * result,
189  unsigned n,
190  const double * x,
191  double * gradient,
192  void * data)
193 {
194  LOG_SCOPE("inequality_constraints()", "NloptOptimizationSolver");
195 
196  libmesh_assert(data);
197 
198  // data should be a pointer to the solver (it was passed in as void *)
200  static_cast<NloptOptimizationSolver<Number> *> (data);
201 
202  OptimizationSystem & sys = solver->system();
203 
204  // We'll use current_local_solution below, so let's ensure that it's consistent
205  // with the vector x that was passed in.
206  if (sys.solution->size() != n)
207  libmesh_error_msg("Error: Input vector x has different length than sys.solution!");
208 
209  for (auto i : index_range(*sys.solution))
210  sys.solution->set(i, x[i]);
211  sys.solution->close();
212 
213  // Impose constraints on the solution vector
214  sys.get_dof_map().enforce_constraints_exactly(sys);
215 
216  // Update sys.current_local_solution based on the solution vector
217  sys.update();
218 
219  // Call the user's inequality constraints function if there is one.
221  if (ineco)
222  {
223  ineco->inequality_constraints(*sys.current_local_solution,
224  *sys.C_ineq,
225  sys);
226 
227  sys.C_ineq->close();
228 
229  // Copy the values out of ineq_constraints into 'result'.
230  // TODO: Even better would be if we could use 'result' directly
231  // as the storage of ineq_constraints. Perhaps a serial-only
232  // NumericVector variant which supports this option?
233  for (unsigned int i = 0; i < m; ++i)
234  result[i] = (*sys.C_ineq)(i);
235 
236  // If gradient != nullptr, then the Jacobian matrix of the equality
237  // constraints has been requested. The incoming 'gradient'
238  // array is of length m*n and d(c_i)/d(x_j) = gradient[n*i+j].
239  if (gradient)
240  {
243 
244  if (ineco_jac)
245  {
246  ineco_jac->inequality_constraints_jacobian(*sys.current_local_solution,
247  *sys.C_ineq_jac,
248  sys);
249 
250  sys.C_ineq_jac->close();
251 
252  // copy the Jacobian data to the gradient array
253  for (numeric_index_type i=0; i<m; i++)
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);
256  }
257  else
258  libmesh_error_msg("Jacobian function not defined in __libmesh_nlopt_inequality_constraints");
259  }
260 
261  }
262  else
263  libmesh_error_msg("Constraints function not defined in __libmesh_nlopt_inequality_constraints");
264 }
265 
266 //---------------------------------------------------------------------
267 
268 
269 
270 //---------------------------------------------------------------------
271 // NloptOptimizationSolver<> methods
272 template <typename T>
274  :
275  OptimizationSolver<T>(system_in),
276  _opt(nullptr),
277  _result(NLOPT_SUCCESS),
278  _iteration_count(0),
279  _constraints_tolerance(1.e-8)
280 {
281  // The nlopt interfaces all use unsigned int as their index type, so
282  // don't risk using the NloptOptimizationSolver with a libmesh that
283  // is configured to use 64-bit indices. We can detect this at
284  // configure time, so it should not be possible to actually reach
285  // this error message... it's here just in case.
286  if (sizeof(dof_id_type) != sizeof(unsigned int))
287  libmesh_error_msg("The NloptOptimizationSolver should not be used with dof_id_type != unsigned int.");
288 }
289 
290 
291 
292 template <typename T>
294 {
295  this->clear ();
296 }
297 
298 
299 
300 template <typename T>
302 {
303  if (this->initialized())
304  {
305  this->_is_initialized = false;
306 
307  nlopt_destroy(_opt);
308  }
309 }
310 
311 
312 
313 template <typename T>
315 {
316  // Initialize the data structures if not done so already.
317  if (!this->initialized())
318  {
319  this->_is_initialized = true;
320 
321  // By default, use the LD_SLSQP solver
322  std::string nlopt_algorithm_name = "LD_SLSQP";
323 
324  if (libMesh::on_command_line("--nlopt-algorithm"))
325  nlopt_algorithm_name = libMesh::command_line_next ("--nlopt-algorithm",
326  nlopt_algorithm_name);
327 
328  // Convert string to an nlopt algorithm type
329  auto it = _nlopt_algorithms.find(nlopt_algorithm_name);
330 
331  if (it == _nlopt_algorithms.end())
332  libmesh_error_msg("Invalid nlopt algorithm requested on command line: " \
333  << nlopt_algorithm_name);
334 
335  _opt = nlopt_create(it->second, this->system().solution->size());
336  }
337 }
338 
339 
340 
341 template <typename T>
343 {
344  LOG_SCOPE("solve()", "NloptOptimizationSolver");
345 
346  this->init ();
347 
348  unsigned int nlopt_size = this->system().solution->size();
349 
350  // We have to have an objective function
351  libmesh_assert( this->objective_object );
352 
353  // Set routine for objective and (optionally) gradient evaluation
354  {
355  nlopt_result ierr =
356  nlopt_set_min_objective(_opt,
358  this);
359  if (ierr < 0)
360  libmesh_error_msg("NLopt failed to set min objective: " << ierr);
361  }
362 
363  if (this->lower_and_upper_bounds_object)
364  {
365  // Need to actually compute the bounds vectors first
366  this->lower_and_upper_bounds_object->lower_and_upper_bounds(this->system());
367 
368  std::vector<Real> nlopt_lb(nlopt_size);
369  std::vector<Real> nlopt_ub(nlopt_size);
370  for (unsigned int i = 0; i < nlopt_size; ++i)
371  {
372  nlopt_lb[i] = this->system().get_vector("lower_bounds")(i);
373  nlopt_ub[i] = this->system().get_vector("upper_bounds")(i);
374  }
375 
376  nlopt_set_lower_bounds(_opt, nlopt_lb.data());
377  nlopt_set_upper_bounds(_opt, nlopt_ub.data());
378  }
379 
380  // If we have an equality constraints object, tell NLopt about it.
381  if (this->equality_constraints_object)
382  {
383  // NLopt requires a vector to specify the tolerance for each constraint.
384  // NLopt makes a copy of this vector internally, so it's safe for us to
385  // let it go out of scope.
386  std::vector<double> equality_constraints_tolerances(this->system().C_eq->size(),
387  _constraints_tolerance);
388 
389  // It would be nice to call the C interface directly, at least it should return an error
390  // code we could parse... unfortunately, there does not seem to be a way to extract
391  // the underlying nlopt_opt object from the nlopt::opt class!
392  nlopt_result ierr =
393  nlopt_add_equality_mconstraint(_opt,
394  equality_constraints_tolerances.size(),
396  this,
397  equality_constraints_tolerances.data());
398 
399  if (ierr < 0)
400  libmesh_error_msg("NLopt failed to add equality constraint: " << ierr);
401  }
402 
403  // If we have an inequality constraints object, tell NLopt about it.
404  if (this->inequality_constraints_object)
405  {
406  // NLopt requires a vector to specify the tolerance for each constraint
407  std::vector<double> inequality_constraints_tolerances(this->system().C_ineq->size(),
408  _constraints_tolerance);
409 
410  nlopt_add_inequality_mconstraint(_opt,
411  inequality_constraints_tolerances.size(),
413  this,
414  inequality_constraints_tolerances.data());
415  }
416 
417  // Set a relative tolerance on the optimization parameters
418  nlopt_set_ftol_rel(_opt, this->objective_function_relative_tolerance);
419 
420  // Set the maximum number of allowed objective function evaluations
421  nlopt_set_maxeval(_opt, this->max_objective_function_evaluations);
422 
423  // Reset internal iteration counter
424  this->_iteration_count = 0;
425 
426  // Perform the optimization
427  std::vector<Real> x(nlopt_size);
428  Real min_val = 0.;
429  _result = nlopt_optimize(_opt, x.data(), &min_val);
430 
431  if (_result < 0)
432  libMesh::out << "NLopt failed!" << std::endl;
433  else
434  libMesh::out << "NLopt obtained optimal value: "
435  << min_val
436  << " in "
437  << this->get_iteration_count()
438  << " iterations."
439  << std::endl;
440 }
441 
442 
443 template <typename T>
445 {
446  libMesh::out << "NLopt optimization solver convergence/divergence reason: "
447  << this->get_converged_reason() << std::endl;
448 }
449 
450 
451 
452 template <typename T>
454 {
455  return static_cast<int>(_result);
456 }
457 
458 
459 //------------------------------------------------------------------
460 // Explicit instantiations
461 template class NloptOptimizationSolver<Number>;
462 
463 } // namespace libMesh
464 
465 
466 
467 #endif // #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
OptimizationSystem::ComputeObjective * objective_object
virtual void gradient(const NumericVector< Number > &X, NumericVector< Number > &grad_f, sys_type &S)=0
virtual void inequality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_ineq_jac, sys_type &S)=0
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Definition: int_range.h:104
OptimizationSystem::ComputeGradient * gradient_object
OptimizationSystem::ComputeInequalityConstraints * inequality_constraints_object
OptimizationSystem::ComputeEqualityConstraintsJacobian * equality_constraints_jacobian_object
dof_id_type numeric_index_type
Definition: id_types.h:92
void init(triangulateio &t)
T command_line_next(std::string name, T value)
Definition: libmesh.C:941
virtual void equality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_eq, sys_type &S)=0
void __libmesh_nlopt_inequality_constraints(unsigned m, double *result, unsigned n, const double *x, double *gradient, void *data)
const sys_type & system() const
virtual Number objective(const NumericVector< Number > &X, sys_type &S)=0
virtual void print_converged_reason() override
virtual void equality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_eq_jac, sys_type &S)=0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
PetscErrorCode ierr
double __libmesh_nlopt_objective(unsigned n, const double *x, double *gradient, void *data)
bool initialized()
Definition: libmesh.C:258
bool on_command_line(std::string arg)
Definition: libmesh.C:876
IterBase * data
OStreamProxy out(std::cout)
virtual void inequality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_ineq, sys_type &S)=0
OptimizationSystem::ComputeInequalityConstraintsJacobian * inequality_constraints_jacobian_object
uint8_t dof_id_type
Definition: id_types.h:64
OptimizationSystem::ComputeEqualityConstraints * equality_constraints_object
void __libmesh_nlopt_equality_constraints(unsigned m, double *result, unsigned n, const double *x, double *gradient, void *data)