ipopt_adapter.cc
Go to the documentation of this file.
1 /******************************************************************************
2 Copyright (c) 2017, Alexander W. Winkler, ETH Zurich. All rights reserved.
3 
4 Redistribution and use in source and binary forms, with or without modification,
5 are permitted provided that the following conditions are met:
6  * Redistributions of source code must retain the above copyright notice,
7  this list of conditions and the following disclaimer.
8  * Redistributions in binary form must reproduce the above copyright notice,
9  this list of conditions and the following disclaimer in the documentation
10  and/or other materials provided with the distribution.
11  * Neither the name of ETH ZURICH nor the names of its contributors may be
12  used to endorse or promote products derived from this software without
13  specific prior written permission.
14 
15 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
16 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
17 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
18 DISCLAIMED. IN NO EVENT SHALL ETH ZURICH BE LIABLE FOR ANY DIRECT, INDIRECT,
19 INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
20 LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
21 PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
22 LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE
23 OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
24 ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
25 ******************************************************************************/
26 
27 #include <ifopt/ipopt_adapter.h>
28 
29 namespace Ipopt {
30 
31 IpoptAdapter::IpoptAdapter(Problem& nlp, bool finite_diff)
32 {
33  nlp_ = &nlp;
34  finite_diff_ = finite_diff;
35 }
36 
37 bool IpoptAdapter::get_nlp_info(Index& n, Index& m, Index& nnz_jac_g,
38  Index& nnz_h_lag, IndexStyleEnum& index_style)
39 {
42 
43  if (finite_diff_)
44  nnz_jac_g = m * n;
45  else
46  nnz_jac_g = nlp_->GetJacobianOfConstraints().nonZeros();
47 
48  nnz_h_lag = n * n;
49 
50  // start index at 0 for row/col entries
51  index_style = C_STYLE;
52 
53  return true;
54 }
55 
56 bool IpoptAdapter::get_bounds_info(Index n, double* x_lower, double* x_upper,
57  Index m, double* g_l, double* g_u)
58 {
59  auto bounds_x = nlp_->GetBoundsOnOptimizationVariables();
60  for (std::size_t c = 0; c < bounds_x.size(); ++c) {
61  x_lower[c] = bounds_x.at(c).lower_;
62  x_upper[c] = bounds_x.at(c).upper_;
63  }
64 
65  // specific bounds depending on equality and inequality constraints
66  auto bounds_g = nlp_->GetBoundsOnConstraints();
67  for (std::size_t c = 0; c < bounds_g.size(); ++c) {
68  g_l[c] = bounds_g.at(c).lower_;
69  g_u[c] = bounds_g.at(c).upper_;
70  }
71 
72  return true;
73 }
74 
75 bool IpoptAdapter::get_starting_point(Index n, bool init_x, double* x,
76  bool init_z, double* z_L, double* z_U,
77  Index m, bool init_lambda, double* lambda)
78 {
79  // Here, we assume we only have starting values for x
80  assert(init_x == true);
81  assert(init_z == false);
82  assert(init_lambda == false);
83 
84  VectorXd x_all = nlp_->GetVariableValues();
85  Eigen::Map<VectorXd>(&x[0], x_all.rows()) = x_all;
86 
87  return true;
88 }
89 
90 bool IpoptAdapter::eval_f(Index n, const double* x, bool new_x,
91  double& obj_value)
92 {
93  obj_value = nlp_->EvaluateCostFunction(x);
94  return true;
95 }
96 
97 bool IpoptAdapter::eval_grad_f(Index n, const double* x, bool new_x,
98  double* grad_f)
99 {
100  Eigen::VectorXd grad = nlp_->EvaluateCostFunctionGradient(x, finite_diff_);
101  Eigen::Map<Eigen::MatrixXd>(grad_f, n, 1) = grad;
102  return true;
103 }
104 
105 bool IpoptAdapter::eval_g(Index n, const double* x, bool new_x, Index m,
106  double* g)
107 {
108  VectorXd g_eig = nlp_->EvaluateConstraints(x);
109  Eigen::Map<VectorXd>(g, m) = g_eig;
110  return true;
111 }
112 
113 bool IpoptAdapter::eval_jac_g(Index n, const double* x, bool new_x, Index m,
114  Index nele_jac, Index* iRow, Index* jCol,
115  double* values)
116 {
117  // defines the positions of the nonzero elements of the jacobian
118  if (values == NULL) {
119  // If "jacobian_approximation" option is set as "finite-difference-values", the Jacobian is dense!
120  Index nele = 0;
121  if (finite_diff_) { // dense jacobian
122  for (Index row = 0; row < m; row++) {
123  for (Index col = 0; col < n; col++) {
124  iRow[nele] = row;
125  jCol[nele] = col;
126  nele++;
127  }
128  }
129  } else { // sparse jacobian
130  auto jac = nlp_->GetJacobianOfConstraints();
131  for (int k = 0; k < jac.outerSize(); ++k) {
132  for (Jacobian::InnerIterator it(jac, k); it; ++it) {
133  iRow[nele] = it.row();
134  jCol[nele] = it.col();
135  nele++;
136  }
137  }
138  }
139 
140  assert(nele ==
141  nele_jac); // initial sparsity structure is never allowed to change
142  } else {
143  // only gets used if "jacobian_approximation finite-difference-values" is not set
144  nlp_->EvalNonzerosOfJacobian(x, values);
145  }
146 
147  return true;
148 }
149 
151  AlgorithmMode mode, Index iter, double obj_value, double inf_pr,
152  double inf_du, double mu, double d_norm, double regularization_size,
153  double alpha_du, double alpha_pr, Index ls_trials, const IpoptData* ip_data,
154  IpoptCalculatedQuantities* ip_cq)
155 {
156  nlp_->SaveCurrent();
157  return true;
158 }
159 
160 void IpoptAdapter::finalize_solution(SolverReturn status, Index n,
161  const double* x, const double* z_L,
162  const double* z_U, Index m,
163  const double* g, const double* lambda,
164  double obj_value, const IpoptData* ip_data,
165  IpoptCalculatedQuantities* ip_cq)
166 {
167  nlp_->SetVariables(x);
168  nlp_->SaveCurrent();
169 }
170 
171 } // namespace Ipopt
ifopt::Problem::GetNumberOfOptimizationVariables
int GetNumberOfOptimizationVariables() const
The number of optimization variables.
Definition: problem.cc:86
ifopt::Problem::EvaluateCostFunction
double EvaluateCostFunction(const double *x)
The scalar cost for current optimization variables x.
Definition: problem.cc:106
ifopt::Problem::SetVariables
void SetVariables(const double *x)
Updates the variables with the values of the raw pointer x.
Definition: problem.cc:101
ifopt::Problem::EvaluateCostFunctionGradient
VectorXd EvaluateCostFunctionGradient(const double *x, bool use_finite_difference_approximation=false, double epsilon=std::numeric_limits< double >::epsilon())
The column-vector of derivatives of the cost w.r.t. each variable.
Definition: problem.cc:116
ifopt::Problem::GetBoundsOnConstraints
VecBound GetBoundsOnConstraints() const
The upper and lower bound of each individual constraint.
Definition: problem.cc:143
ifopt::Problem::EvalNonzerosOfJacobian
void EvalNonzerosOfJacobian(const double *x, double *values)
Extracts those entries from constraint Jacobian that are not zero.
Definition: problem.cc:164
Ipopt
namespace defined by the Ipopt solver.
Definition: ipopt_adapter.h:42
Ipopt::IpoptAdapter::nlp_
Problem * nlp_
The solver independent problem definition.
Definition: ipopt_adapter.h:118
Ipopt::IpoptAdapter::eval_jac_g
virtual bool eval_jac_g(Index n, const double *x, bool new_x, Index m, Index nele_jac, Index *iRow, Index *jCol, double *values)
Definition: ipopt_adapter.cc:137
ifopt::Problem::GetNumberOfConstraints
int GetNumberOfConstraints() const
The number of individual constraints.
Definition: problem.cc:148
ifopt::Problem::GetJacobianOfConstraints
Jacobian GetJacobianOfConstraints() const
The sparse-matrix representation of Jacobian of the constraints.
Definition: problem.cc:173
ifopt::Problem::EvaluateConstraints
VectorXd EvaluateConstraints(const double *x)
Each constraint value g(x) for current optimization variables x.
Definition: problem.cc:153
ifopt::Problem::GetVariableValues
VectorXd GetVariableValues() const
The current value of the optimization variables.
Definition: problem.cc:96
Ipopt::IpoptAdapter::get_nlp_info
virtual bool get_nlp_info(Index &n, Index &m, Index &nnz_jac_g, Index &nnz_h_lag, IndexStyleEnum &index_style)
Definition: ipopt_adapter.cc:61
ifopt::Problem::GetBoundsOnOptimizationVariables
VecBound GetBoundsOnOptimizationVariables() const
The maximum and minimum value each optimization variable is allowed to have.
Definition: problem.cc:91
Ipopt::IpoptAdapter::eval_f
virtual bool eval_f(Index n, const double *x, bool new_x, double &obj_value)
Definition: ipopt_adapter.cc:114
Ipopt::IpoptAdapter::IpoptAdapter
IpoptAdapter(Problem &nlp, bool finite_diff=false)
Creates an IpoptAdapter wrapping the nlp.
Definition: ipopt_adapter.cc:55
Ipopt::IpoptAdapter::finite_diff_
bool finite_diff_
Flag that indicates the "finite-difference-values" option is set.
Definition: ipopt_adapter.h:120
ifopt::Problem::SaveCurrent
void SaveCurrent()
Saves the current values of the optimization variables in x_prev.
Definition: problem.cc:183
Ipopt::IpoptAdapter::get_bounds_info
virtual bool get_bounds_info(Index n, double *x_l, double *x_u, Index m, double *g_l, double *g_u)
Definition: ipopt_adapter.cc:80
Ipopt::IpoptAdapter::finalize_solution
virtual void finalize_solution(SolverReturn status, Index n, const double *x, const double *z_L, const double *z_U, Index m, const double *g, const double *lambda, double obj_value, const IpoptData *ip_data, IpoptCalculatedQuantities *ip_cq)
Definition: ipopt_adapter.cc:184
Ipopt::IpoptAdapter::get_starting_point
virtual bool get_starting_point(Index n, bool init_x, double *x, bool init_z, double *z_L, double *z_U, Index m, bool init_lambda, double *lambda)
Definition: ipopt_adapter.cc:99
Ipopt::IpoptAdapter::VectorXd
Problem::VectorXd VectorXd
Definition: ipopt_adapter.h:105
ipopt_adapter.h
Ipopt::IpoptAdapter::eval_g
virtual bool eval_g(Index n, const double *x, bool new_x, Index m, double *g)
Definition: ipopt_adapter.cc:129
Ipopt::IpoptAdapter::eval_grad_f
virtual bool eval_grad_f(Index n, const double *x, bool new_x, double *grad_f)
Definition: ipopt_adapter.cc:121
Ipopt::IpoptAdapter::intermediate_callback
virtual bool intermediate_callback(AlgorithmMode mode, Index iter, double obj_value, double inf_pr, double inf_du, double mu, double d_norm, double regularization_size, double alpha_du, double alpha_pr, Index ls_trials, const IpoptData *ip_data, IpoptCalculatedQuantities *ip_cq)
Definition: ipopt_adapter.cc:174


ifopt
Author(s): Alexander W. Winkler
autogenerated on Mon Sep 18 2023 02:14:38