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,
78  double* lambda)
79 {
80  // Here, we assume we only have starting values for x
81  assert(init_x == true);
82  assert(init_z == false);
83  assert(init_lambda == false);
84 
85  VectorXd x_all = nlp_->GetVariableValues();
86  Eigen::Map<VectorXd>(&x[0], x_all.rows()) = x_all;
87 
88  return true;
89 }
90 
91 bool IpoptAdapter::eval_f(Index n, const double* x, bool new_x, 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, double* grad_f)
98 {
99  Eigen::VectorXd grad = nlp_->EvaluateCostFunctionGradient(x, finite_diff_);
100  Eigen::Map<Eigen::MatrixXd>(grad_f,n,1) = grad;
101  return true;
102 }
103 
104 bool IpoptAdapter::eval_g(Index n, const double* x, bool new_x, Index m, double* g)
105 {
106  VectorXd g_eig = nlp_->EvaluateConstraints(x);
107  Eigen::Map<VectorXd>(g,m) = g_eig;
108  return true;
109 }
110 
111 bool IpoptAdapter::eval_jac_g(Index n, const double* x, bool new_x,
112  Index m, Index nele_jac, Index* iRow, Index *jCol,
113  double* values)
114 {
115  // defines the positions of the nonzero elements of the jacobian
116  if (values == NULL) {
117  // If "jacobian_approximation" option is set as "finite-difference-values", the Jacobian is dense!
118  Index nele=0;
119  if (finite_diff_) { // dense jacobian
120  for (Index row = 0; row < m; row++) {
121  for (Index col = 0; col < n; col++) {
122  iRow[nele] = row;
123  jCol[nele] = col;
124  nele++;
125  }
126  }
127  }
128  else { // sparse jacobian
129  auto jac = nlp_->GetJacobianOfConstraints();
130  for (int k=0; k<jac.outerSize(); ++k) {
131  for (Jacobian::InnerIterator it(jac,k); it; ++it) {
132  iRow[nele] = it.row();
133  jCol[nele] = it.col();
134  nele++;
135  }
136  }
137  }
138 
139  assert(nele == nele_jac); // initial sparsity structure is never allowed to change
140  }
141  else {
142  // only gets used if "jacobian_approximation finite-difference-values" is not set
143  nlp_->EvalNonzerosOfJacobian(x, values);
144  }
145 
146  return true;
147 }
148 
149 bool IpoptAdapter::intermediate_callback(AlgorithmMode mode,
150  Index iter, double obj_value,
151  double inf_pr, double inf_du,
152  double mu, double d_norm,
153  double regularization_size,
154  double alpha_du, double alpha_pr,
155  Index ls_trials,
156  const IpoptData* ip_data,
157  IpoptCalculatedQuantities* ip_cq)
158 {
159  nlp_->SaveCurrent();
160  return true;
161 }
162 
163 void IpoptAdapter::finalize_solution(SolverReturn status,
164  Index n, const double* x, const double* z_L, const double* z_U,
165  Index m, const double* g, const double* lambda,
166  double obj_value,
167  const IpoptData* ip_data,
168  IpoptCalculatedQuantities* ip_cq)
169 {
170  nlp_->SetVariables(x);
171  nlp_->SaveCurrent();
172 }
173 
174 } // namespace Ipopt
A generic optimization problem with variables, costs and constraints.
Definition: problem.h:97
virtual bool eval_grad_f(Index n, const double *x, bool new_x, double *grad_f)
VecBound GetBoundsOnConstraints() const
The upper and lower bound of each individual constraint.
Definition: problem.cc:130
double EvaluateCostFunction(const double *x)
The scalar cost for current optimization variables x.
Definition: problem.cc:89
IpoptAdapter(Problem &nlp, bool finite_diff=false)
Creates an IpoptAdapter wrapping the nlp.
virtual bool get_bounds_info(Index n, double *x_l, double *x_u, Index m, double *g_l, double *g_u)
VectorXd GetVariableValues() const
The current value of the optimization variables.
Definition: problem.cc:77
virtual bool get_nlp_info(Index &n, Index &m, Index &nnz_jac_g, Index &nnz_h_lag, IndexStyleEnum &index_style)
namespace defined by the Ipopt solver.
Definition: ipopt_adapter.h:42
void SaveCurrent()
Saves the current values of the optimization variables in x_prev.
Definition: problem.cc:177
VecBound GetBoundsOnOptimizationVariables() const
The maximum and minimum value each optimization variable is allowed to have.
Definition: problem.cc:71
void SetVariables(const double *x)
Updates the variables with the values of the raw pointer x.
Definition: problem.cc:83
Jacobian GetJacobianOfConstraints() const
The sparse-matrix representation of Jacobian of the constraints.
Definition: problem.cc:165
virtual bool eval_f(Index n, const double *x, bool new_x, double &obj_value)
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)
Problem * nlp_
The solver independent problem definition.
Definition: ipopt_adapter.h:70
virtual bool eval_jac_g(Index n, const double *x, bool new_x, Index m, Index nele_jac, Index *iRow, Index *jCol, double *values)
Problem::VectorXd VectorXd
Definition: ipopt_adapter.h:57
void EvalNonzerosOfJacobian(const double *x, double *values)
Extracts those entries from constraint Jacobian that are not zero.
Definition: problem.cc:155
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)
int GetNumberOfConstraints() const
The number of individual constraints.
Definition: problem.cc:136
virtual bool eval_g(Index n, const double *x, bool new_x, Index m, double *g)
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)
VectorXd EvaluateConstraints(const double *x)
Each constraint value g(x) for current optimization variables x.
Definition: problem.cc:142
bool finite_diff_
Flag that indicates the "finite-difference-values" option is set.
Definition: ipopt_adapter.h:71
VectorXd EvaluateCostFunctionGradient(const double *x, bool use_finite_difference_approximation=false)
The column-vector of derivatives of the cost w.r.t. each variable.
Definition: problem.cc:100
int GetNumberOfOptimizationVariables() const
The number of optimization variables.
Definition: problem.cc:65


ifopt
Author(s): Alexander W. Winkler
autogenerated on Fri Jan 22 2021 03:47:32