Scippy

SCIP

Solving Constraint Integer Programs

nlpi_ipopt.cpp
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2 /* */
3 /* This file is part of the program and library */
4 /* SCIP --- Solving Constraint Integer Programs */
5 /* */
6 /* Copyright (C) 2002-2018 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not email to scip@zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file nlpi_ipopt.cpp
17  * @ingroup NLPIS
18  * @brief Ipopt NLP interface
19  * @author Stefan Vigerske
20  * @author Benjamin Müller
21  *
22  * @todo warm starts
23  * @todo use new_x: Ipopt sets new_x = false if any function has been evaluated for the current x already, while oracle allows new_x to be false only if the current function has been evaluated for the current x before
24  *
25  * This file can only be compiled if Ipopt is available.
26  * Otherwise, to resolve public functions, use nlpi_ipopt_dummy.c.
27  * Since the dummy code is C instead of C++, it has been moved into a separate file.
28  */
29 
30 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
31 
32 #include "nlpi/nlpi_ipopt.h"
33 
34 #include "nlpi/nlpi.h"
35 #include "nlpi/nlpioracle.h"
36 #include "nlpi/exprinterpret.h"
37 #include "scip/interrupt.h"
38 #include "scip/pub_misc.h"
39 #include "scip/misc.h"
40 
41 #include <new> /* for std::bad_alloc */
42 #include <sstream>
43 
44 #ifdef __GNUC__
45 #pragma GCC diagnostic ignored "-Wshadow"
46 #endif
47 #include "IpoptConfig.h"
48 #include "IpIpoptApplication.hpp"
49 #include "IpIpoptCalculatedQuantities.hpp"
50 #include "IpSolveStatistics.hpp"
51 #include "IpJournalist.hpp"
52 #include "IpIpoptData.hpp"
53 #include "IpTNLPAdapter.hpp"
54 #include "IpOrigIpoptNLP.hpp"
55 #include "IpLapack.hpp"
56 #ifdef __GNUC__
57 #pragma GCC diagnostic warning "-Wshadow"
58 #endif
59 
60 #if (IPOPT_VERSION_MAJOR < 3 || (IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 12))
61 #error "The Ipopt interface requires at least 3.12."
62 #endif
63 
64 using namespace Ipopt;
65 
66 #define NLPI_NAME "ipopt" /**< short concise name of solver */
67 #define NLPI_DESC "Ipopt interface" /**< description of solver */
68 #define NLPI_PRIORITY 1000 /**< priority */
69 
70 #ifdef SCIP_DEBUG
71 #define DEFAULT_PRINTLEVEL J_WARNING /**< default print level of Ipopt */
72 #else
73 #define DEFAULT_PRINTLEVEL J_STRONGWARNING /**< default print level of Ipopt */
74 #endif
75 #define DEFAULT_MAXITER 3000 /**< default iteration limit for Ipopt */
76 
77 #define MAXPERTURB 0.01 /**< maximal perturbation of bounds in starting point heuristic */
78 #define FEASTOLFACTOR 0.05 /**< factor for user-given feasibility tolerance to get feasibility tolerance that is actually passed to Ipopt */
79 
80 #define DEFAULT_RANDSEED 71 /**< initial random seed */
81 
82 /* Convergence check (see ScipNLP::intermediate_callback)
83  *
84  * If the fastfail option is enabled, then we stop Ipopt if the reduction in
85  * primal infeasibility is not sufficient for a consecutive number of iterations.
86  * With the parameters as given below, we require Ipopt to
87  * - not increase the primal infeasibility after 5 iterations
88  * - reduce the primal infeasibility by at least 50% within 10 iterations
89  * - reduce the primal infeasibility by at least 90% within 30 iterations
90  * The targets are updated once they are reached and the limit on allowed iterations to reach the new target is reset.
91  *
92  * In certain situations, it is allowed to exceed an iteration limit:
93  * - If we are in the first 10 (convcheck_startiter) iterations.
94  * - If we are within 10 (convcheck_startiter) iterations after the restoration phase ended.
95  * The reason for this is that during feasibility restoration phase Ipopt aims completely on
96  * reducing constraint violation, completely forgetting the objective function.
97  * When returning from feasibility restoration and considering the original objective again,
98  * it is unlikely that Ipopt will continue to decrease primal infeasibility, since it may now target on
99  * more on optimality again. Thus, we do not check convergence for a number of iterations.
100  * - If the target on dual infeasibility reduction has been achieved, we are below twice the iteration limit, and
101  * we are not in restoration mode.
102  * The reason for this is that if Ipopt makes good progress towards optimality,
103  * we want to allow some more iterations where primal infeasibility is not reduced.
104  * However, in restoration mode, dual infeasibility does not correspond to the original problem and
105  * the complete aim is to restore primal infeasibility.
106  */
107 static const int convcheck_nchecks = 3; /**< number of convergence checks */
108 static const int convcheck_startiter = 10; /**< iteration where to start convergence checking */
109 static const int convcheck_maxiter[convcheck_nchecks] = { 5, 15, 30 }; /**< maximal number of iterations to achieve each convergence check */
110 static const SCIP_Real convcheck_minred[convcheck_nchecks] = { 1.0, 0.5, 0.1 }; /**< minimal required infeasibility reduction in each convergence check */
111 
112 class ScipNLP;
113 
114 struct SCIP_NlpiData
115 {
116 public:
117  BMS_BLKMEM* blkmem; /**< block memory */
118  SCIP_MESSAGEHDLR* messagehdlr; /**< message handler */
119  SCIP_Real infinity; /**< initial value for infinity */
120  std::string defoptions; /**< modified default options for Ipopt */
121 
122  /** constructor */
123  explicit SCIP_NlpiData(
124  BMS_BLKMEM* blkmem_ /**< block memory */
125  )
126  : blkmem(blkmem_),
127  messagehdlr(NULL),
129  { }
130 };
131 
132 struct SCIP_NlpiProblem
133 {
134 public:
135  SCIP_NLPIORACLE* oracle; /**< Oracle-helper to store and evaluate NLP */
136 
137  SmartPtr<IpoptApplication> ipopt; /**< Ipopt application */
138  SmartPtr<ScipNLP> nlp; /**< NLP in Ipopt form */
139  std::string optfile; /**< name of options file */
140  bool storeintermediate;/**< whether to store intermediate solutions */
141  bool fastfail; /**< whether to stop Ipopt if convergence seems slow */
142 
143  SCIP_Bool firstrun; /**< whether the next NLP solve will be the first one (with the current problem structure) */
144  SCIP_Real* initguess; /**< initial values for primal variables, or NULL if not known */
145 
146  SCIP_NLPSOLSTAT lastsolstat; /**< solution status from last run */
147  SCIP_NLPTERMSTAT lasttermstat; /**< termination status from last run */
148  SCIP_Real* lastsolprimals; /**< primal solution values from last run, if available */
149  SCIP_Real* lastsoldualcons; /**< dual solution values of constraints from last run, if available */
150  SCIP_Real* lastsoldualvarlb; /**< dual solution values of variable lower bounds from last run, if available */
151  SCIP_Real* lastsoldualvarub; /**< dual solution values of variable upper bounds from last run, if available */
152  SCIP_Real lastsolinfeas;/**< infeasibility (constraint violation) of solution stored in lastsolprimals */
153  int lastniter; /**< number of iterations in last run */
154  SCIP_Real lasttime; /**< time spend in last run */
155 
156  /** constructor */
158  : oracle(NULL),
159  storeintermediate(false), fastfail(false),
160  firstrun(TRUE), initguess(NULL),
161  lastsolstat(SCIP_NLPSOLSTAT_UNKNOWN), lasttermstat(SCIP_NLPTERMSTAT_OTHER),
162  lastsolprimals(NULL), lastsoldualcons(NULL), lastsoldualvarlb(NULL), lastsoldualvarub(NULL),
163  lastniter(-1), lasttime(-1.0)
164  { }
165 };
166 
167 /** TNLP implementation for SCIPs NLP */
168 class ScipNLP : public TNLP
169 {
170 private:
171  SCIP_NLPIPROBLEM* nlpiproblem; /**< NLPI problem data */
172  SCIP_RANDNUMGEN* randnumgen; /**< random number generator */
173  BMS_BLKMEM* blkmem; /**< block memory */
174 
175  SCIP_Real conv_prtarget[convcheck_nchecks]; /**< target primal infeasibility for each convergence check */
176  SCIP_Real conv_dutarget[convcheck_nchecks]; /**< target dual infeasibility for each convergence check */
177  int conv_iterlim[convcheck_nchecks]; /**< iteration number where target primal infeasibility should to be achieved */
178  int conv_lastrestoiter; /**< last iteration number in restoration mode, or -1 if none */
179 
180 public:
181  bool approxhessian; /**< do we tell Ipopt to approximate the hessian? (may also be false if user set to approx. hessian via option file) */
182 
183  // cppcheck-suppress uninitMemberVar
184  /** constructor */
185  ScipNLP(
186  SCIP_NLPIPROBLEM* nlpiproblem_ = NULL,/**< NLPI problem data */
187  BMS_BLKMEM* blkmem_ = NULL /**< block memory */
188  )
189  : nlpiproblem(nlpiproblem_), randnumgen(NULL), blkmem(blkmem_), conv_lastrestoiter(-1), approxhessian(false)
190  {
191  assert(blkmem != NULL);
193  }
194 
195  /** destructor */
196  ~ScipNLP()
197  {
198  assert(randnumgen != NULL);
199  SCIPrandomFree(&randnumgen, blkmem);
200  }
201 
202  /** sets NLPI data structure */
203  void setNLPIPROBLEM(SCIP_NLPIPROBLEM* nlpiproblem_)
204  {
205  assert(nlpiproblem_ != NULL);
206  nlpiproblem = nlpiproblem_;
207  }
208 
209  /** Method to return some info about the nlp */
210  bool get_nlp_info(
211  Index& n, /**< place to store number of variables */
212  Index& m, /**< place to store number of constraints */
213  Index& nnz_jac_g, /**< place to store number of nonzeros in jacobian */
214  Index& nnz_h_lag, /**< place to store number of nonzeros in hessian */
215  IndexStyleEnum& index_style /**< place to store used index style (0-based or 1-based) */
216  );
217 
218  /** Method to return the bounds for my problem */
219  bool get_bounds_info(
220  Index n, /**< number of variables */
221  Number* x_l, /**< buffer to store lower bounds on variables */
222  Number* x_u, /**< buffer to store upper bounds on variables */
223  Index m, /**< number of constraints */
224  Number* g_l, /**< buffer to store lower bounds on constraints */
225  Number* g_u /**< buffer to store lower bounds on constraints */
226  );
227 
228  /** Method to return the starting point for the algorithm */
229  bool get_starting_point(
230  Index n, /**< number of variables */
231  bool init_x, /**< whether initial values for primal values are requested */
232  Number* x, /**< buffer to store initial primal values */
233  bool init_z, /**< whether initial values for dual values of variable bounds are requested */
234  Number* z_L, /**< buffer to store dual values for variable lower bounds */
235  Number* z_U, /**< buffer to store dual values for variable upper bounds */
236  Index m, /**< number of constraints */
237  bool init_lambda, /**< whether initial values for dual values of constraints are required */
238  Number* lambda /**< buffer to store dual values of constraints */
239  );
240 
241  /** Method to return the variables linearity. */
242  bool get_variables_linearity(
243  Index n, /**< number of variables */
244  LinearityType* var_types /**< buffer to store linearity types of variables */
245  );
246 
247  /** Method to return the constraint linearity. */
248  bool get_constraints_linearity(
249  Index m, /**< number of constraints */
250  LinearityType* const_types /**< buffer to store linearity types of constraints */
251  );
252 
253  /** Method to return the number of nonlinear variables. */
254  Index get_number_of_nonlinear_variables();
255 
256  /** Method to return the indices of the nonlinear variables */
257  bool get_list_of_nonlinear_variables(
258  Index num_nonlin_vars, /**< number of nonlinear variables */
259  Index* pos_nonlin_vars /**< array to fill with indices of nonlinear variables */
260  );
261 
262  /** Method to return metadata about variables and constraints */
263  bool get_var_con_metadata(
264  Index n, /**< number of variables */
265  StringMetaDataMapType& var_string_md, /**< variable meta data of string type */
266  IntegerMetaDataMapType& var_integer_md,/**< variable meta data of integer type */
267  NumericMetaDataMapType& var_numeric_md,/**< variable meta data of numeric type */
268  Index m, /**< number of constraints */
269  StringMetaDataMapType& con_string_md, /**< constraint meta data of string type */
270  IntegerMetaDataMapType& con_integer_md,/**< constraint meta data of integer type */
271  NumericMetaDataMapType& con_numeric_md /**< constraint meta data of numeric type */
272  );
273 
274  /** Method to return the objective value */
275  bool eval_f(
276  Index n, /**< number of variables */
277  const Number* x, /**< point to evaluate */
278  bool new_x, /**< whether some function evaluation method has been called for this point before */
279  Number& obj_value /**< place to store objective function value */
280  );
281 
282  /** Method to return the gradient of the objective */
283  bool eval_grad_f(
284  Index n, /**< number of variables */
285  const Number* x, /**< point to evaluate */
286  bool new_x, /**< whether some function evaluation method has been called for this point before */
287  Number* grad_f /**< buffer to store objective gradient */
288  );
289 
290  /** Method to return the constraint residuals */
291  bool eval_g(
292  Index n, /**< number of variables */
293  const Number* x, /**< point to evaluate */
294  bool new_x, /**< whether some function evaluation method has been called for this point before */
295  Index m, /**< number of constraints */
296  Number* g /**< buffer to store constraint function values */
297  );
298 
299  /** Method to return:
300  * 1) The structure of the jacobian (if "values" is NULL)
301  * 2) The values of the jacobian (if "values" is not NULL)
302  */
303  bool eval_jac_g(
304  Index n, /**< number of variables */
305  const Number* x, /**< point to evaluate */
306  bool new_x, /**< whether some function evaluation method has been called for this point before */
307  Index m, /**< number of constraints */
308  Index nele_jac, /**< number of nonzero entries in jacobian */
309  Index* iRow, /**< buffer to store row indices of nonzero jacobian entries, or NULL if values
310  * are requested */
311  Index* jCol, /**< buffer to store column indices of nonzero jacobian entries, or NULL if values
312  * are requested */
313  Number* values /**< buffer to store values of nonzero jacobian entries, or NULL if structure is
314  * requested */
315  );
316 
317  /** Method to return:
318  * 1) The structure of the hessian of the lagrangian (if "values" is NULL)
319  * 2) The values of the hessian of the lagrangian (if "values" is not NULL)
320  */
321  bool eval_h(
322  Index n, /**< number of variables */
323  const Number* x, /**< point to evaluate */
324  bool new_x, /**< whether some function evaluation method has been called for this point before */
325  Number obj_factor, /**< weight for objective function */
326  Index m, /**< number of constraints */
327  const Number* lambda, /**< weights for constraint functions */
328  bool new_lambda, /**< whether the hessian has been evaluated for these values of lambda before */
329  Index nele_hess, /**< number of nonzero entries in hessian */
330  Index* iRow, /**< buffer to store row indices of nonzero hessian entries, or NULL if values
331  * are requested */
332  Index* jCol, /**< buffer to store column indices of nonzero hessian entries, or NULL if values
333  * are requested */
334  Number* values /**< buffer to store values of nonzero hessian entries, or NULL if structure is requested */
335  );
336 
337  /** Method called by the solver at each iteration.
338  *
339  * Checks whether Ctrl-C was hit.
340  */
341  bool intermediate_callback(
342  AlgorithmMode mode, /**< current mode of algorithm */
343  Index iter, /**< current iteration number */
344  Number obj_value, /**< current objective value */
345  Number inf_pr, /**< current primal infeasibility */
346  Number inf_du, /**< current dual infeasibility */
347  Number mu, /**< current barrier parameter */
348  Number d_norm, /**< current gradient norm */
349  Number regularization_size,/**< current size of regularization */
350  Number alpha_du, /**< current dual alpha */
351  Number alpha_pr, /**< current primal alpha */
352  Index ls_trials, /**< current number of linesearch trials */
353  const IpoptData* ip_data, /**< pointer to Ipopt Data */
354  IpoptCalculatedQuantities* ip_cq /**< pointer to current calculated quantities */
355  );
356 
357  /** This method is called when the algorithm is complete so the TNLP can store/write the solution. */
358  void finalize_solution(
359  SolverReturn status, /**< solve and solution status */
360  Index n, /**< number of variables */
361  const Number* x, /**< primal solution values */
362  const Number* z_L, /**< dual values of variable lower bounds */
363  const Number* z_U, /**< dual values of variable upper bounds */
364  Index m, /**< number of constraints */
365  const Number* g, /**< values of constraints */
366  const Number* lambda, /**< dual values of constraints */
367  Number obj_value, /**< objective function value */
368  const IpoptData* data, /**< pointer to Ipopt Data */
369  IpoptCalculatedQuantities* cq /**< pointer to calculated quantities */
370  );
371 };
372 
373 /** A particular Ipopt::Journal implementation that uses the SCIP message routines for output. */
374 class ScipJournal : public Ipopt::Journal {
375 private:
376  /** reference to message handler pointer in NLPI data */
377  SCIP_MESSAGEHDLR*& messagehdlr;
378 
379 public:
380  ScipJournal(
381  const char* name, /**< name of journal */
382  Ipopt::EJournalLevel default_level, /**< default verbosity level */
383  SCIP_MESSAGEHDLR*& messagehdlr_ /**< pointer where to get message handler from */
384  )
385  : Ipopt::Journal(name, default_level),
386  messagehdlr(messagehdlr_)
387  { }
388 
389  ~ScipJournal() { }
390 
391 protected:
392  void PrintImpl(
393  Ipopt::EJournalCategory category, /**< category of message */
394  Ipopt::EJournalLevel level, /**< verbosity level of message */
395  const char* str /**< message to print */
396  )
397  {
398  if( level == J_ERROR )
399  {
401  }
402  else
403  {
404  SCIPmessagePrintInfo(messagehdlr, str);
405  }
406  }
407 
408  void PrintfImpl(
409  Ipopt::EJournalCategory category, /**< category of message */
410  Ipopt::EJournalLevel level, /**< verbosity level of message */
411  const char* pformat, /**< message printing format */
412  va_list ap /**< arguments of message */
413  )
414  {
415  if( level == J_ERROR )
416  {
417  SCIPmessageVPrintError(pformat, ap);
418  }
419  else
420  {
421  SCIPmessageVPrintInfo(messagehdlr, pformat, ap);
422  }
423  }
424 
425  void FlushBufferImpl() { }
426 };
427 
428 /** clears the last solution arrays and sets the solstat and termstat to unknown and other, resp. */
429 static
431  SCIP_NLPIPROBLEM* problem /**< data structure of problem */
432  )
433 {
434  assert(problem != NULL);
435 
442  problem->lastsolinfeas = SCIP_INVALID;
443 }
444 
445 /** sets feasibility tolerance parameters in Ipopt
446  *
447  * Sets tol and constr_viol_tol to FEASTOLFACTOR*feastol and acceptable_tol and acceptable_viol_tol to feastol.
448  * Since the users and Ipopts conception of feasibility may differ, we let Ipopt try to compute solutions
449  * that are more accurate (w.r.t. constraint violation) than requested by the user.
450  * Only if Ipopt has problems to achieve this accuracy, we also accept solutions that are accurate w.r.t. feastol only.
451  * The additional effort for computing a more accurate solution should be small if one can assume fast convergence when close to a local minimizer.
452  */
453 static
455  SCIP_NLPIPROBLEM* nlpiproblem,
456  SCIP_Real feastol
457  )
458 {
459  assert(nlpiproblem != NULL);
460 
461  nlpiproblem->ipopt->Options()->SetNumericValue("tol", FEASTOLFACTOR * feastol);
462  nlpiproblem->ipopt->Options()->SetNumericValue("constr_viol_tol", FEASTOLFACTOR * feastol);
463 
464  nlpiproblem->ipopt->Options()->SetNumericValue("acceptable_tol", feastol);
465  nlpiproblem->ipopt->Options()->SetNumericValue("acceptable_constr_viol_tol", feastol);
466 
467  /* It seem to be better to let Ipopt relax bounds a bit to ensure that a relative interior exists.
468  * However, if we relax the bounds too much, then the solutions tend to be slightly infeasible.
469  * If the user wants to set a tight feasibility tolerance, then (s)he has probably difficulties to compute accurate enough solutions.
470  * Thus, we turn off the bound_relax_factor completely if it would be below its default value of 1e-8.
471  */
472  nlpiproblem->ipopt->Options()->SetNumericValue("bound_relax_factor", feastol < 1e-8/FEASTOLFACTOR ? 0.0 : FEASTOLFACTOR * feastol);
473 }
474 
475 /** sets optimality tolerance parameters in Ipopt
476  *
477  * Sets dual_inf_tol and compl_inf_tol to opttol.
478  * We leave acceptable_dual_inf_tol and acceptable_compl_inf_tol untouched, which means that if Ipopt has convergence problems, then
479  * it can stop with a solution that is still feasible (see setFeastol), but essentially without a proof of local optimality.
480  * Note, that we report the solution as feasible only if Ipopt stopped on an "acceptable point" (see ScipNLP::finalize_solution).
481  *
482  * Note, that parameters tol and acceptable_tol (set in setFeastol) also apply.
483  */
484 static
486  SCIP_NLPIPROBLEM* nlpiproblem,
487  SCIP_Real opttol
488  )
489 {
490  assert(nlpiproblem != NULL);
491 
492  nlpiproblem->ipopt->Options()->SetNumericValue("dual_inf_tol", opttol);
493  nlpiproblem->ipopt->Options()->SetNumericValue("compl_inf_tol", opttol);
494 }
495 
496 /** copy method of NLP interface (called when SCIP copies plugins)
497  *
498  * input:
499  * - blkmem block memory of target SCIP
500  * - sourcenlpi the NLP interface to copy
501  * - targetnlpi buffer to store pointer to copy of NLP interface
502  */
503 static
504 SCIP_DECL_NLPICOPY(nlpiCopyIpopt)
505 {
506  SCIP_NLPIDATA* sourcedata;
507  SCIP_NLPIDATA* targetdata;
508 
509  assert(sourcenlpi != NULL);
510  assert(targetnlpi != NULL);
511 
512  sourcedata = SCIPnlpiGetData(sourcenlpi);
513  assert(sourcedata != NULL);
514 
515  SCIP_CALL( SCIPcreateNlpSolverIpopt(blkmem, targetnlpi) );
516  assert(*targetnlpi != NULL);
517 
518  SCIP_CALL( SCIPnlpiSetRealPar(*targetnlpi, NULL, SCIP_NLPPAR_INFINITY, sourcedata->infinity) );
519  SCIP_CALL( SCIPnlpiSetMessageHdlr(*targetnlpi, sourcedata->messagehdlr) );
520 
521  targetdata = SCIPnlpiGetData(*targetnlpi);
522  assert(targetdata != NULL);
523 
524  targetdata->defoptions = sourcedata->defoptions;
525 
526  return SCIP_OKAY;
527 }
528 
529 /** destructor of NLP interface to free nlpi data
530  *
531  * input:
532  * - nlpi datastructure for solver interface
533  */
534 static
535 SCIP_DECL_NLPIFREE(nlpiFreeIpopt)
536 {
537  SCIP_NLPIDATA* data;
538 
539  assert(nlpi != NULL);
540 
541  data = SCIPnlpiGetData(nlpi);
542  assert(data != NULL);
543 
544  delete data;
545 
546  return SCIP_OKAY;
547 }
548 
549 /** gets pointer for NLP solver to do dirty stuff
550  *
551  * input:
552  * - nlpi datastructure for solver interface
553  *
554  * return: void pointer to solver
555  */
556 static
557 SCIP_DECL_NLPIGETSOLVERPOINTER(nlpiGetSolverPointerIpopt)
558 {
559  assert(nlpi != NULL);
560 
561  return NULL;
562 }
563 
564 /** creates a problem instance
565  *
566  * input:
567  * - nlpi datastructure for solver interface
568  * - problem pointer to store the problem data
569  * - name name of problem, can be NULL
570  */
571 static
572 SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemIpopt)
573 {
574  SCIP_NLPIDATA* data;
575 
576  assert(nlpi != NULL);
577  assert(problem != NULL);
578 
579  data = SCIPnlpiGetData(nlpi);
580  assert(data != NULL);
581 
582  *problem = new SCIP_NLPIPROBLEM;
583  if( *problem == NULL )
584  return SCIP_NOMEMORY;
585 
586  SCIP_CALL( SCIPnlpiOracleCreate(data->blkmem, &(*problem)->oracle) );
587  SCIP_CALL( SCIPnlpiOracleSetInfinity((*problem)->oracle, data->infinity) );
588  SCIP_CALL( SCIPnlpiOracleSetProblemName((*problem)->oracle, name) );
589 
590  try
591  {
592  /* initialize IPOPT without default journal */
593  (*problem)->ipopt = new IpoptApplication(false);
594  if( IsNull((*problem)->ipopt) )
595  throw std::bad_alloc();
596 
597  /* plugin our journal to get output through SCIP message handler */
598  SmartPtr<Journal> jrnl = new ScipJournal("console", J_ITERSUMMARY, data->messagehdlr);
599  if( IsNull(jrnl) )
600  throw std::bad_alloc();
601  jrnl->SetPrintLevel(J_DBG, J_NONE);
602  if( !(*problem)->ipopt->Jnlst()->AddJournal(jrnl) )
603  {
604  SCIPerrorMessage("Failed to register ScipJournal for IPOPT output.");
605  }
606 
607  /* initialize Ipopt/SCIP NLP interface */
608  (*problem)->nlp = new ScipNLP(*problem, data->blkmem);
609  if( IsNull((*problem)->nlp) )
610  throw std::bad_alloc();
611  }
612  catch( std::bad_alloc )
613  {
614  SCIPerrorMessage("Not enough memory to initialize Ipopt.\n");
615  return SCIP_NOMEMORY;
616  }
617 
618  /* modify Ipopt's default settings to what we believe is appropriate */
619  (*problem)->ipopt->RegOptions()->AddStringOption2("store_intermediate", "whether to store the most feasible intermediate solutions", "no", "yes", "", "no", "", "useful when Ipopt looses a once found feasible solution and then terminates with an infeasible point");
620  (*problem)->ipopt->Options()->SetIntegerValue("print_level", DEFAULT_PRINTLEVEL);
621  /* (*problem)->ipopt->Options()->SetStringValue("print_timing_statistics", "yes"); */
622  (*problem)->ipopt->Options()->SetStringValue("mu_strategy", "adaptive");
623  (*problem)->ipopt->Options()->SetStringValue("expect_infeasible_problem", "yes");
624  (*problem)->ipopt->Options()->SetIntegerValue("max_iter", DEFAULT_MAXITER);
625  (*problem)->ipopt->Options()->SetNumericValue("nlp_lower_bound_inf", -data->infinity, false);
626  (*problem)->ipopt->Options()->SetNumericValue("nlp_upper_bound_inf", data->infinity, false);
627  (*problem)->ipopt->Options()->SetNumericValue("diverging_iterates_tol", data->infinity, false);
628  /* (*problem)->ipopt->Options()->SetStringValue("dependency_detector", "ma28"); */
629  /* if the expression interpreter does not give hessians, tell Ipopt to approximate hessian */
630 #ifdef SCIP_DEBUG
631  (*problem)->ipopt->Options()->SetStringValue("derivative_test", "second-order");
632 #endif
633  setFeastol(*problem, SCIP_DEFAULT_FEASTOL);
635 
636  /* apply user's given modifications to Ipopt's default settings */
637  if( data->defoptions.length() > 0 )
638  {
639  std::istringstream is(data->defoptions);
640 
641 #if (IPOPT_VERSION_MAJOR > 3) || (IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR > 12) || (IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR == 12 && IPOPT_VERSION_RELEASE >= 5)
642  if( !(*problem)->ipopt->Options()->ReadFromStream(*(*problem)->ipopt->Jnlst(), is, true) )
643 #else
644  if( !(*problem)->ipopt->Options()->ReadFromStream(*(*problem)->ipopt->Jnlst(), is) )
645 #endif
646  {
647  SCIPerrorMessage("Error when modifiying Ipopt options using options string\n%s\n", data->defoptions.c_str());
648  return SCIP_ERROR;
649  }
650  }
651 
652  /* apply user's given options file (this one is NLPI problem specific) */
653  if( (*problem)->ipopt->Initialize((*problem)->optfile) != Solve_Succeeded )
654  {
655  SCIPerrorMessage("Error during initialization of Ipopt using optionfile \"%s\"\n", (*problem)->optfile.c_str());
656  return SCIP_ERROR;
657  }
658 
659 
660  return SCIP_OKAY;
661 }
662 
663 /** free a problem instance
664  *
665  * input:
666  * - nlpi datastructure for solver interface
667  * - problem pointer where problem data is stored
668  */
669 static
670 SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemIpopt)
671 {
672  assert(nlpi != NULL);
673  assert(problem != NULL);
674  assert(*problem != NULL);
675 
676  if( (*problem)->oracle != NULL )
677  {
678  SCIP_CALL( SCIPnlpiOracleFree(&(*problem)->oracle) );
679  }
680 
681  BMSfreeMemoryArrayNull(&(*problem)->initguess);
682  BMSfreeMemoryArrayNull(&(*problem)->lastsolprimals);
683  BMSfreeMemoryArrayNull(&(*problem)->lastsoldualcons);
684  BMSfreeMemoryArrayNull(&(*problem)->lastsoldualvarlb);
685  BMSfreeMemoryArrayNull(&(*problem)->lastsoldualvarub);
686 
687  delete *problem;
688  *problem = NULL;
689 
690  return SCIP_OKAY;
691 }
692 
693 /** gets pointer to solver-internal problem instance to do dirty stuff
694  *
695  * input:
696  * - nlpi datastructure for solver interface
697  * - problem datastructure for problem instance
698  *
699  * return: void pointer to problem instance
700  */
701 static
702 SCIP_DECL_NLPIGETPROBLEMPOINTER(nlpiGetProblemPointerIpopt)
703 {
704  assert(nlpi != NULL);
705  assert(problem != NULL);
706 
707  return GetRawPtr(problem->nlp);
708 }
709 
710 /** add variables
711  *
712  * input:
713  * - nlpi datastructure for solver interface
714  * - problem datastructure for problem instance
715  * - nvars number of variables
716  * - lbs lower bounds of variables, can be NULL if -infinity
717  * - ubs upper bounds of variables, can be NULL if +infinity
718  * - varnames names of variables, can be NULL
719  */
720 static
721 SCIP_DECL_NLPIADDVARS(nlpiAddVarsIpopt)
722 {
723  assert(nlpi != NULL);
724  assert(problem != NULL);
725  assert(problem->oracle != NULL);
726 
727  SCIP_CALL( SCIPnlpiOracleAddVars(problem->oracle, nvars, lbs, ubs, varnames) );
728 
729  problem->firstrun = TRUE;
730  BMSfreeMemoryArrayNull(&problem->initguess);
731  invalidateSolution(problem);
732 
733  return SCIP_OKAY;
734 }
735 
736 /** add constraints
737  *
738  * quadratic coefficiens: row oriented matrix for each constraint
739  *
740  * input:
741  * - nlpi datastructure for solver interface
742  * - problem datastructure for problem instance
743  * - ncons number of added constraints
744  * - lhss left hand sides of constraints
745  * - rhss right hand sides of constraints
746  * - linoffsets start index of each constraints linear coefficients in lininds and linvals
747  * length: ncons + 1, linoffsets[ncons] gives length of lininds and linvals
748  * may be NULL in case of no linear part
749  * - lininds variable indices
750  * may be NULL in case of no linear part
751  * - linvals coefficient values
752  * may be NULL in case of no linear part
753  * - nquadelems number of quadratic elements for each constraint
754  * may be NULL in case of no quadratic part
755  * - quadelems quadratic elements for each constraint
756  * may be NULL in case of no quadratic part
757  * - exprvaridxs indices of variables in expression tree, maps variable indices in expression tree to indices in nlp
758  * entry of array may be NULL in case of no expression tree
759  * may be NULL in case of no expression tree in any constraint
760  * - exprtrees expression tree for nonquadratic part of constraints
761  * entry of array may be NULL in case of no nonquadratic part
762  * may be NULL in case of no nonquadratic part in any constraint
763  * - names of constraints, may be NULL or entries may be NULL
764  */
765 static
766 SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsIpopt)
767 {
768  assert(nlpi != NULL);
769  assert(problem != NULL);
770  assert(problem->oracle != NULL);
771 
772  SCIP_CALL( SCIPnlpiOracleAddConstraints(problem->oracle,
773  ncons, lhss, rhss,
774  nlininds, lininds, linvals,
775  nquadelems, quadelems,
776  exprvaridxs, exprtrees, names) );
777 
778  problem->firstrun = TRUE;
779  invalidateSolution(problem);
780 
781  return SCIP_OKAY;
782 }
783 
784 /** sets or overwrites objective, a minimization problem is expected
785  *
786  * May change sparsity pattern.
787  *
788  * input:
789  * - nlpi datastructure for solver interface
790  * - problem datastructure for problem instance
791  * - nlins number of linear variables
792  * - lininds variable indices
793  * may be NULL in case of no linear part
794  * - linvals coefficient values
795  * may be NULL in case of no linear part
796  * - nquadelems number of elements in matrix of quadratic part
797  * - quadelems elements of quadratic part
798  * may be NULL in case of no quadratic part
799  * - exprvaridxs indices of variables in expression tree, maps variable indices in expression tree to indices in nlp
800  * may be NULL in case of no expression tree
801  * - exprtree expression tree for nonquadratic part of objective function
802  * may be NULL in case of no nonquadratic part
803  * - constant objective value offset
804  */
805 static
806 SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveIpopt)
807 {
808  assert(nlpi != NULL);
809  assert(problem != NULL);
810  assert(problem->oracle != NULL);
811 
812  /* We pass the objective gradient in dense form to Ipopt, so if the sparsity of that gradient changes, we do not need to reset Ipopt (firstrun=TRUE).
813  * However, if the sparsity of the Hessian matrix of the objective changes, then the sparsity pattern of the Hessian of the Lagrangian may change.
814  * Thus, reset Ipopt if the objective was and/or becomes nonlinear, but leave firstrun untouched if it was and stays linear.
815  */
816  if( nquadelems > 0 || exprtree != NULL || SCIPnlpiOracleGetConstraintDegree(problem->oracle, -1) > 1 )
817  problem->firstrun = TRUE;
818 
819  SCIP_CALL( SCIPnlpiOracleSetObjective(problem->oracle,
820  constant, nlins, lininds, linvals,
821  nquadelems, quadelems,
822  exprvaridxs, exprtree) );
823 
824  invalidateSolution(problem);
825 
826  return SCIP_OKAY;
827 }
828 
829 /** change variable bounds
830  *
831  * input:
832  * - nlpi datastructure for solver interface
833  * - problem datastructure for problem instance
834  * - nvars number of variables to change bounds
835  * - indices indices of variables to change bounds
836  * - lbs new lower bounds
837  * - ubs new upper bounds
838  */
839 static
840 SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsIpopt)
841 {
842  int i;
843 
844  assert(nlpi != NULL);
845  assert(problem != NULL);
846  assert(problem->oracle != NULL);
847 
848  /* If some variable is fixed or unfixed, then better don't reoptimize the NLP in the next solve.
849  * Calling Optimize instead of ReOptimize should remove fixed variables from the problem that is solved by Ipopt.
850  * This way, the variable fixing is satisfied exactly in a solution, see also #1254.
851  */
852  for( i = 0; i < nvars && !problem->firstrun; ++i )
853  if( (SCIPnlpiOracleGetVarLbs(problem->oracle)[indices[i]] == SCIPnlpiOracleGetVarUbs(problem->oracle)[indices[i]]) != (lbs[i] == ubs[i]) )
854  problem->firstrun = TRUE;
855 
856  SCIP_CALL( SCIPnlpiOracleChgVarBounds(problem->oracle, nvars, indices, lbs, ubs) );
857 
858  invalidateSolution(problem);
859 
860  return SCIP_OKAY;
861 }
862 
863 /** change constraint bounds
864  *
865  * input:
866  * - nlpi datastructure for solver interface
867  * - problem datastructure for problem instance
868  * - nconss number of constraints to change sides
869  * - indices indices of constraints to change sides
870  * - lhss new left hand sides
871  * - rhss new right hand sides
872  */
873 static
874 SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesIpopt)
875 {
876  assert(nlpi != NULL);
877  assert(problem != NULL);
878  assert(problem->oracle != NULL);
879 
880  SCIP_CALL( SCIPnlpiOracleChgConsSides(problem->oracle, nconss, indices, lhss, rhss) );
881 
882  invalidateSolution(problem);
883 
884  return SCIP_OKAY;
885 }
886 
887 /** delete a set of variables
888  *
889  * input:
890  * - nlpi datastructure for solver interface
891  * - problem datastructure for problem instance
892  * - nlpi datastructure for solver interface
893  * - dstats deletion status of vars; 1 if var should be deleted, 0 if not
894  *
895  * output:
896  * - dstats new position of var, -1 if var was deleted
897  */
898 static
899 SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetIpopt)
900 {
901  assert(nlpi != NULL);
902  assert(problem != NULL);
903  assert(problem->oracle != NULL);
904 
905  SCIP_CALL( SCIPnlpiOracleDelVarSet(problem->oracle, dstats) );
906 
907  problem->firstrun = TRUE;
908  BMSfreeMemoryArrayNull(&problem->initguess); // @TODO keep initguess for remaining variables
909 
910  invalidateSolution(problem);
911 
912  return SCIP_OKAY;
913 }
914 
915 /** delete a set of constraints
916  *
917  * input:
918  * - nlpi datastructure for solver interface
919  * - problem datastructure for problem instance
920  * - dstats deletion status of rows; 1 if row should be deleted, 0 if not
921  *
922  * output:
923  * - dstats new position of row, -1 if row was deleted
924  */
925 static
926 SCIP_DECL_NLPIDELCONSSET(nlpiDelConstraintSetIpopt)
927 {
928  assert(nlpi != NULL);
929  assert(problem != NULL);
930  assert(problem->oracle != NULL);
931 
932  SCIP_CALL( SCIPnlpiOracleDelConsSet(problem->oracle, dstats) );
933 
934  problem->firstrun = TRUE;
935 
936  invalidateSolution(problem);
937 
938  return SCIP_OKAY;
939 }
940 
941 /** change one linear coefficient in a constraint or objective
942  *
943  * input:
944  * - nlpi datastructure for solver interface
945  * - problem datastructure for problem instance
946  * - idx index of constraint or -1 for objective
947  * - nvals number of values in linear constraint
948  * - varidxs indices of variable
949  * - vals new values for coefficient
950  *
951  * return: Error if coefficient did not exist before
952  */
953 static
954 SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsIpopt)
955 {
956  assert(nlpi != NULL);
957  assert(problem != NULL);
958  assert(problem->oracle != NULL);
959 
960  SCIP_CALL( SCIPnlpiOracleChgLinearCoefs(problem->oracle, idx, nvals, varidxs, vals) );
961  invalidateSolution(problem);
962 
963  return SCIP_OKAY;
964 }
965 
966 /** change one coefficient in the quadratic part of a constraint or objective
967  *
968  * input:
969  * - nlpi datastructure for solver interface
970  * - problem datastructure for problem instance
971  * - idx index of constraint or -1 for objective
972  * - nquadelems number of entries in quadratic matrix to change
973  * - quadelems new elements in quadratic matrix (replacing already existing ones or adding new ones)
974  *
975  * return: Error if coefficient did not exist before
976  */
977 static
978 SCIP_DECL_NLPICHGQUADCOEFS(nlpiChgQuadraticCoefsIpopt)
979 {
980  assert(nlpi != NULL);
981  assert(problem != NULL);
982  assert(problem->oracle != NULL);
983 
984  SCIP_CALL( SCIPnlpiOracleChgQuadCoefs(problem->oracle, idx, nquadelems, quadelems) );
985  invalidateSolution(problem);
986 
987  return SCIP_OKAY;
988 }
989 
990 /** replaces the expression tree of a constraint or objective
991  *
992  * input:
993  * - nlpi datastructure for solver interface
994  * - problem datastructure for problem instance
995  * - idxcons index of constraint or -1 for objective
996  * - exprtree new expression tree for constraint or objective, or NULL to only remove previous tree
997  */
998 static
999 SCIP_DECL_NLPICHGEXPRTREE(nlpiChgExprtreeIpopt)
1000 {
1001  assert(nlpi != NULL);
1002  assert(problem != NULL);
1003  assert(problem->oracle != NULL);
1004 
1005  SCIP_CALL( SCIPnlpiOracleChgExprtree(problem->oracle, idxcons, exprvaridxs, exprtree) );
1006  invalidateSolution(problem);
1007 
1008  return SCIP_OKAY;
1009 }
1010 
1011 /** change the value of one parameter in the nonlinear part
1012  *
1013  * input:
1014  * - nlpi datastructure for solver interface
1015  * - problem datastructure for problem instance
1016  * - idxcons index of constraint or -1 for objective
1017  * - idxparam index of parameter
1018  * - value new value for nonlinear parameter
1019  *
1020  * return: Error if parameter does not exist
1021  */
1022 static
1023 SCIP_DECL_NLPICHGNONLINCOEF(nlpiChgNonlinCoefIpopt)
1024 {
1025  assert(nlpi != NULL);
1026  assert(problem != NULL);
1027  assert(problem->oracle != NULL);
1028 
1029  SCIP_CALL( SCIPnlpiOracleChgExprParam(problem->oracle, idxcons, idxparam, value) );
1030  invalidateSolution(problem);
1031 
1032  return SCIP_OKAY;
1033 }
1034 
1035 /** change the constant offset in the objective
1036  *
1037  * input:
1038  * - nlpi datastructure for solver interface
1039  * - problem datastructure for problem instance
1040  * - objconstant new value for objective constant
1041  */
1042 static
1043 SCIP_DECL_NLPICHGOBJCONSTANT( nlpiChgObjConstantIpopt )
1044 {
1045  assert(nlpi != NULL);
1046  assert(problem != NULL);
1047  assert(problem->oracle != NULL);
1048 
1049  SCIP_CALL( SCIPnlpiOracleChgObjConstant(problem->oracle, objconstant) );
1050 
1051  return SCIP_OKAY;
1052 }
1053 
1054 /** sets initial guess for primal variables
1055  *
1056  * input:
1057  * - nlpi datastructure for solver interface
1058  * - problem datastructure for problem instance
1059  * - primalvalues initial primal values for variables, or NULL to clear previous values
1060  * - consdualvalues initial dual values for constraints, or NULL to clear previous values
1061  * - varlbdualvalues initial dual values for variable lower bounds, or NULL to clear previous values
1062  * - varubdualvalues initial dual values for variable upper bounds, or NULL to clear previous values
1063  */
1064 static
1065 SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessIpopt)
1066 {
1067  assert(nlpi != NULL);
1068  assert(problem != NULL);
1069  assert(problem->oracle != NULL);
1070 
1071  if( primalvalues != NULL )
1072  {
1073  if( !problem->initguess )
1074  {
1075  if( BMSduplicateMemoryArray(&problem->initguess, primalvalues, SCIPnlpiOracleGetNVars(problem->oracle)) == NULL )
1076  return SCIP_NOMEMORY;
1077  }
1078  else
1079  {
1080  BMScopyMemoryArray(problem->initguess, primalvalues, SCIPnlpiOracleGetNVars(problem->oracle));
1081  }
1082  }
1083  else
1084  {
1085  BMSfreeMemoryArrayNull(&problem->initguess);
1086  }
1087 
1088  return SCIP_OKAY;
1089 }
1090 
1091 /** tries to solve NLP
1092  *
1093  * input:
1094  * - nlpi datastructure for solver interface
1095  * - problem datastructure for problem instance
1096  */
1097 static
1098 SCIP_DECL_NLPISOLVE(nlpiSolveIpopt)
1099 {
1100  ApplicationReturnStatus status;
1101 
1102  assert(nlpi != NULL);
1103  assert(problem != NULL);
1104  assert(problem->oracle != NULL);
1105 
1106  assert(IsValid(problem->ipopt));
1107  assert(IsValid(problem->nlp));
1108 
1109  problem->nlp->setNLPIPROBLEM(problem);
1110 
1111  problem->lastniter = -1;
1112  problem->lasttime = -1.0;
1113  problem->lastsolinfeas = SCIP_INVALID;
1114 
1115  try
1116  {
1117  SmartPtr<SolveStatistics> stats;
1118 
1119  if( problem->firstrun )
1120  {
1122 
1123  cap = SCIPexprintGetCapability() & SCIPnlpiOracleGetEvalCapability(problem->oracle);
1124 
1125  /* if the expression interpreter or some user expression do not support function values and gradients and Hessians, and the problem is not at most quadratic,
1126  * change NLP parameters or give an error
1127  */
1129  && SCIPnlpiOracleGetMaxDegree(problem->oracle) > 2 )
1130  {
1131  /* @todo could enable Jacobian approximation in Ipopt */
1134  {
1135  SCIPerrorMessage("Do not have expression interpreter that can compute function values and gradients. Cannot solve NLP with Ipopt.\n");
1136  problem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
1137  problem->lasttermstat = SCIP_NLPTERMSTAT_OTHER;
1138  return SCIP_OKAY;
1139  }
1140 
1141  /* enable Hessian approximation if we are nonquadratic and the expression interpreter or user expression do not support Hessians */
1142  if( !(cap & SCIP_EXPRINTCAPABILITY_HESSIAN) )
1143  {
1144  problem->ipopt->Options()->SetStringValue("hessian_approximation", "limited-memory");
1145  problem->nlp->approxhessian = true;
1146  }
1147  else
1148  problem->nlp->approxhessian = false;
1149  }
1150 
1151  status = problem->ipopt->OptimizeTNLP(GetRawPtr(problem->nlp));
1152  }
1153  else
1154  {
1155  status = problem->ipopt->ReOptimizeTNLP(GetRawPtr(problem->nlp));
1156  }
1157 
1158  // catch the very bad status codes
1159  switch( status ) {
1160  case Invalid_Problem_Definition:
1161  case Invalid_Option:
1162  case Unrecoverable_Exception:
1163  case NonIpopt_Exception_Thrown:
1164  case Internal_Error:
1165  SCIPerrorMessage("Ipopt returned with application return status %d\n", status);
1166  return SCIP_ERROR;
1167  case Insufficient_Memory:
1168  SCIPerrorMessage("Ipopt returned with status \"Insufficient Memory\"\n");
1169  return SCIP_NOMEMORY;
1170  case Invalid_Number_Detected:
1171  SCIPdebugMessage("Ipopt failed because of an invalid number in function or derivative value\n");
1172  invalidateSolution(problem);
1173  problem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
1174  problem->lasttermstat = SCIP_NLPTERMSTAT_EVALERR;
1175  default: ;
1176  }
1177 
1178  stats = problem->ipopt->Statistics();
1179  if( IsValid(stats) )
1180  {
1181  problem->lastniter = stats->IterationCount();
1182  problem->lasttime = stats->TotalCPUTime();
1183  }
1184  else
1185  {
1186  /* Ipopt does not provide access to the statistics when all variables have been fixed */
1187  problem->lastniter = 0;
1188  problem->lasttime = 0.0;
1189  }
1190  }
1191  catch( IpoptException& except )
1192  {
1193  SCIPerrorMessage("Ipopt returned with exception: %s\n", except.Message().c_str());
1194  return SCIP_ERROR;
1195  }
1196 
1197  problem->firstrun = FALSE;
1198 
1199  return SCIP_OKAY;
1200 }
1201 
1202 /** gives solution status
1203  *
1204  * input:
1205  * - nlpi datastructure for solver interface
1206  * - problem datastructure for problem instance
1207  *
1208  * return: Solution Status
1209  */
1210 static
1211 SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatIpopt)
1212 {
1213  assert(nlpi != NULL);
1214  assert(problem != NULL);
1215 
1216  return problem->lastsolstat;
1217 }
1218 
1219 /** gives termination reason
1220  *
1221  * input:
1222  * - nlpi datastructure for solver interface
1223  * - problem datastructure for problem instance
1224  *
1225  * return: Termination Status
1226  */
1227 static
1228 SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatIpopt)
1229 {
1230  assert(nlpi != NULL);
1231  assert(problem != NULL);
1232 
1233  return problem->lasttermstat;
1234 }
1235 
1236 /** gives primal and dual solution values
1237  *
1238  * input:
1239  * - nlpi datastructure for solver interface
1240  * - problem datastructure for problem instance
1241  * - primalvalues buffer to store pointer to array to primal values, or NULL if not needed
1242  * - consdualvalues buffer to store pointer to array to dual values of constraints, or NULL if not needed
1243  * - varlbdualvalues buffer to store pointer to array to dual values of variable lower bounds, or NULL if not needed
1244  * - varubdualvalues buffer to store pointer to array to dual values of variable lower bounds, or NULL if not needed
1245  * - objval buffer store the objective value, or NULL if not needed
1246  */
1247 static
1248 SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionIpopt)
1249 {
1250  assert(nlpi != NULL);
1251  assert(problem != NULL);
1252 
1253  if( primalvalues != NULL )
1254  *primalvalues = problem->lastsolprimals;
1255 
1256  if( consdualvalues != NULL )
1257  *consdualvalues = problem->lastsoldualcons;
1258 
1259  if( varlbdualvalues != NULL )
1260  *varlbdualvalues = problem->lastsoldualvarlb;
1261 
1262  if( varubdualvalues != NULL )
1263  *varubdualvalues = problem->lastsoldualvarub;
1264 
1265  if( objval != NULL )
1266  {
1267  if( problem->lastsolprimals != NULL )
1268  {
1269  /* TODO store last solution value instead of reevaluating the objective function */
1270  SCIP_CALL( SCIPnlpiOracleEvalObjectiveValue(problem->oracle, problem->lastsolprimals, objval) );
1271  }
1272  else
1273  *objval = SCIP_INVALID;
1274  }
1275 
1276  return SCIP_OKAY;
1277 }
1278 
1279 /** gives solve statistics
1280  *
1281  * input:
1282  * - nlpi datastructure for solver interface
1283  * - problem datastructure for problem instance
1284  * - statistics pointer to store statistics
1285  *
1286  * output:
1287  * - statistics solve statistics
1288  */
1289 static
1290 SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsIpopt)
1291 {
1292  assert(nlpi != NULL);
1293  assert(problem != NULL);
1294 
1295  SCIPnlpStatisticsSetNIterations(statistics, problem->lastniter);
1296  SCIPnlpStatisticsSetTotalTime (statistics, problem->lasttime);
1297 
1298  return SCIP_OKAY;
1299 }
1300 
1301 /** gives required size of a buffer to store a warmstart object
1302  *
1303  * input:
1304  * - nlpi datastructure for solver interface
1305  * - problem datastructure for problem instance
1306  * - size pointer to store required size for warmstart buffer
1307  *
1308  * output:
1309  * - size required size for warmstart buffer
1310  */
1311 static
1312 SCIP_DECL_NLPIGETWARMSTARTSIZE(nlpiGetWarmstartSizeIpopt)
1313 {
1314  SCIPerrorMessage("method of Ipopt nonlinear solver is not implemented\n");
1315  return SCIP_ERROR;
1316 }
1317 
1318 /** stores warmstart information in buffer
1319  *
1320  * Required size of buffer should have been obtained by SCIPnlpiGetWarmstartSize before.
1321  *
1322  * input:
1323  * - nlpi datastructure for solver interface
1324  * - problem datastructure for problem instance
1325  * - buffer memory to store warmstart information
1326  *
1327  * output:
1328  * - buffer warmstart information in solver specific data structure
1329  */
1330 static
1331 SCIP_DECL_NLPIGETWARMSTARTMEMO(nlpiGetWarmstartMemoIpopt)
1332 {
1333  SCIPerrorMessage("method of Ipopt nonlinear solver is not implemented\n");
1334  return SCIP_ERROR;
1335 }
1336 
1337 /** sets warmstart information in solver
1338  *
1339  * Write warmstart to buffer.
1340  *
1341  * input:
1342  * - nlpi datastructure for solver interface
1343  * - problem datastructure for problem instance
1344  * - buffer warmstart information
1345  */
1346 static
1347 SCIP_DECL_NLPISETWARMSTARTMEMO(nlpiSetWarmstartMemoIpopt)
1348 {
1349  SCIPerrorMessage("method of Ipopt nonlinear solver is not implemented\n");
1350  SCIPABORT();
1351  return SCIP_OKAY;
1352 }
1353 
1354 /** gets integer parameter of NLP
1355  *
1356  * input:
1357  * - nlpi datastructure for solver interface
1358  * - problem datastructure for problem instance
1359  * - type parameter number
1360  * - ival pointer to store the parameter value
1361  *
1362  * output:
1363  * - ival parameter value
1364  */
1365 static
1366 SCIP_DECL_NLPIGETINTPAR(nlpiGetIntParIpopt)
1367 {
1368  assert(nlpi != NULL);
1369  assert(ival != NULL);
1370  assert(problem != NULL);
1371  assert(IsValid(problem->ipopt));
1372 
1373  //@TODO try-catch block for Ipopt exceptions
1374  switch( type )
1375  {
1377  {
1378  *ival = 1;
1379  break;
1380  }
1381 
1382  case SCIP_NLPPAR_VERBLEVEL:
1383  {
1384  int printlevel;
1385  problem->ipopt->Options()->GetIntegerValue("print_level", printlevel, "");
1386  if( printlevel <= J_STRONGWARNING )
1387  *ival = 0;
1388  else if( printlevel >= J_DETAILED )
1389  *ival = printlevel - J_ITERSUMMARY + 1;
1390  else /* J_SUMMARY or J_WARNING or J_ITERSUMMARY */
1391  *ival = 1;
1392  break;
1393  }
1394 
1395  case SCIP_NLPPAR_FEASTOL:
1396  {
1397  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
1398  return SCIP_PARAMETERWRONGTYPE;
1399  }
1400 
1401  case SCIP_NLPPAR_RELOBJTOL:
1402  {
1403  SCIPerrorMessage("relative objective tolerance parameter is of type real.\n");
1404  return SCIP_PARAMETERWRONGTYPE;
1405  }
1406 
1407  case SCIP_NLPPAR_LOBJLIM:
1408  {
1409  SCIPerrorMessage("objective limit parameter is of type real.\n");
1410  return SCIP_PARAMETERWRONGTYPE;
1411  }
1412 
1413  case SCIP_NLPPAR_INFINITY:
1414  {
1415  SCIPerrorMessage("infinity parameter is of type real.\n");
1416  return SCIP_PARAMETERWRONGTYPE;
1417  }
1418 
1419  case SCIP_NLPPAR_ITLIM:
1420  {
1421  problem->ipopt->Options()->GetIntegerValue("max_iter", *ival, "");
1422  break;
1423  }
1424 
1425  case SCIP_NLPPAR_TILIM:
1426  {
1427  SCIPerrorMessage("time limit parameter is of type real.\n");
1428  return SCIP_PARAMETERWRONGTYPE;
1429  }
1430 
1431  case SCIP_NLPPAR_OPTFILE:
1432  {
1433  SCIPerrorMessage("optfile parameter is of type string.\n");
1434  return SCIP_PARAMETERWRONGTYPE;
1435  }
1436 
1437  case SCIP_NLPPAR_FASTFAIL:
1438  {
1439  *ival = problem->fastfail ? 1 : 0;
1440  break;
1441  }
1442 
1443  default:
1444  {
1445  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
1446  return SCIP_PARAMETERUNKNOWN;
1447  }
1448  }
1449 
1450  return SCIP_OKAY;
1451 }
1452 
1453 /** sets integer parameter of NLP
1454  *
1455  * input:
1456  * - nlpi datastructure for solver interface
1457  * - problem datastructure for problem instance
1458  * - type parameter number
1459  * - ival parameter value
1460  */
1461 static
1462 SCIP_DECL_NLPISETINTPAR(nlpiSetIntParIpopt)
1463 {
1464  assert(nlpi != NULL);
1465  assert(problem != NULL);
1466  assert(IsValid(problem->ipopt));
1467 
1468  switch( type )
1469  {
1471  {
1472  if( ival == 0 || ival == 1 )
1473  {
1474  SCIP_NLPIDATA* data;
1475 
1476  data = SCIPnlpiGetData(nlpi);
1477  assert(data != NULL);
1478 
1479  SCIPmessagePrintWarning(data->messagehdlr, "from scratch parameter not supported by Ipopt interface yet. Ignored.\n");
1480  }
1481  else
1482  {
1483  SCIPerrorMessage("Value %d for parameter from scratch out of range {0, 1}\n", ival);
1484  return SCIP_PARAMETERWRONGVAL;
1485  }
1486  break;
1487  }
1488 
1489  case SCIP_NLPPAR_VERBLEVEL:
1490  {
1491  switch( ival )
1492  {
1493  case 0:
1494  problem->ipopt->Options()->SetIntegerValue("print_level", J_STRONGWARNING);
1495  break;
1496  case 1:
1497  problem->ipopt->Options()->SetIntegerValue("print_level", J_ITERSUMMARY);
1498  break;
1499  case 2:
1500  problem->ipopt->Options()->SetIntegerValue("print_level", J_DETAILED);
1501  break;
1502  default:
1503  if( ival > 2 )
1504  {
1505  problem->ipopt->Options()->SetIntegerValue("print_level", MIN(J_ITERSUMMARY + (ival-1), J_ALL));
1506  break;
1507  }
1508  else
1509  {
1510  SCIPerrorMessage("Value %d for parameter from verbosity level out of range {0, 1, 2}\n", ival);
1511  return SCIP_PARAMETERWRONGVAL;
1512  }
1513  }
1514  break;
1515  }
1516 
1517  case SCIP_NLPPAR_FEASTOL:
1518  {
1519  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
1520  return SCIP_PARAMETERWRONGTYPE;
1521  }
1522 
1523  case SCIP_NLPPAR_RELOBJTOL:
1524  {
1525  SCIPerrorMessage("relative objective tolerance parameter is of type real.\n");
1526  return SCIP_PARAMETERWRONGTYPE;
1527  }
1528 
1529  case SCIP_NLPPAR_LOBJLIM:
1530  {
1531  SCIPerrorMessage("objective limit parameter is of type real.\n");
1532  return SCIP_PARAMETERWRONGTYPE;
1533  }
1534 
1535  case SCIP_NLPPAR_INFINITY:
1536  {
1537  SCIPerrorMessage("infinity parameter is of type real.\n");
1538  return SCIP_PARAMETERWRONGTYPE;
1539  }
1540 
1541  case SCIP_NLPPAR_ITLIM:
1542  {
1543  if( ival >= 0 )
1544  {
1545  problem->ipopt->Options()->SetIntegerValue("max_iter", ival);
1546  }
1547  else
1548  {
1549  SCIPerrorMessage("Value %d for parameter iteration limit is negative\n", ival);
1550  return SCIP_PARAMETERWRONGVAL;
1551  }
1552  break;
1553  }
1554 
1555  case SCIP_NLPPAR_TILIM:
1556  {
1557  SCIPerrorMessage("time limit parameter is of type real.\n");
1558  return SCIP_PARAMETERWRONGTYPE;
1559  }
1560 
1561  case SCIP_NLPPAR_OPTFILE:
1562  {
1563  SCIPerrorMessage("optfile parameter is of type string.\n");
1564  return SCIP_PARAMETERWRONGTYPE;
1565  }
1566 
1567  case SCIP_NLPPAR_FASTFAIL:
1568  {
1569  if( ival == 0 || ival == 1 )
1570  {
1571  problem->fastfail = (bool)ival;
1572  problem->storeintermediate = (bool)ival;
1573  }
1574  else
1575  {
1576  SCIPerrorMessage("Value %d for parameter fastfail out of range {0, 1}\n", ival);
1577  return SCIP_PARAMETERWRONGVAL;
1578  }
1579  break;
1580  }
1581 
1582  default:
1583  {
1584  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
1585  return SCIP_PARAMETERUNKNOWN;
1586  }
1587  }
1588 
1589  return SCIP_OKAY;
1590 }
1591 
1592 /** gets floating point parameter of NLP
1593  *
1594  * input:
1595  * - nlpi datastructure for solver interface
1596  * - problem datastructure for problem instance, can be NULL only if type == SCIP_NLPPAR_INFINITY
1597  * - type parameter number
1598  * - dval pointer to store the parameter value
1599  *
1600  * output:
1601  * - dval parameter value
1602  */
1603 static
1604 SCIP_DECL_NLPIGETREALPAR(nlpiGetRealParIpopt)
1605 {
1606  assert(nlpi != NULL);
1607  assert(dval != NULL);
1608  assert(type == SCIP_NLPPAR_INFINITY || problem != NULL);
1609  assert(type == SCIP_NLPPAR_INFINITY || IsValid(problem->ipopt));
1610 
1611  switch( type )
1612  {
1614  {
1615  SCIPerrorMessage("from scratch parameter is of type int.\n");
1616  return SCIP_PARAMETERWRONGTYPE;
1617  }
1618 
1619  case SCIP_NLPPAR_VERBLEVEL:
1620  {
1621  SCIPerrorMessage("verbosity level parameter is of type int.\n");
1622  return SCIP_PARAMETERWRONGTYPE;
1623  }
1624 
1625  case SCIP_NLPPAR_FEASTOL:
1626  {
1627  problem->ipopt->Options()->GetNumericValue("acceptable_constr_viol_tol", *dval, "");
1628  break;
1629  }
1630 
1631  case SCIP_NLPPAR_RELOBJTOL:
1632  {
1633  problem->ipopt->Options()->GetNumericValue("dual_inf_tol", *dval, "");
1634  break;
1635  }
1636 
1637  case SCIP_NLPPAR_LOBJLIM:
1638  {
1639  *dval = -SCIPnlpiOracleGetInfinity(problem->oracle);
1640  break;
1641  }
1642 
1643  case SCIP_NLPPAR_INFINITY:
1644  {
1645  if( problem )
1646  {
1647  *dval = SCIPnlpiOracleGetInfinity(problem->oracle);
1648  }
1649  else
1650  {
1651  SCIP_NLPIDATA* data = SCIPnlpiGetData(nlpi);
1652  assert(data != NULL);
1653  *dval = data->infinity;
1654  }
1655  break;
1656  }
1657 
1658  case SCIP_NLPPAR_ITLIM:
1659  {
1660  SCIPerrorMessage("iteration limit parameter is of type int.\n");
1661  return SCIP_PARAMETERWRONGTYPE;
1662  }
1663 
1664  case SCIP_NLPPAR_TILIM:
1665  {
1666  problem->ipopt->Options()->GetNumericValue("max_cpu_time", *dval, "");
1667  break;
1668  }
1669 
1670  case SCIP_NLPPAR_OPTFILE:
1671  {
1672  SCIPerrorMessage("option file parameter is of type string.\n");
1673  return SCIP_PARAMETERWRONGTYPE;
1674  }
1675 
1676  case SCIP_NLPPAR_FASTFAIL:
1677  {
1678  SCIPerrorMessage("fastfail parameter is of type int.\n");
1679  return SCIP_PARAMETERWRONGTYPE;
1680  }
1681 
1682  default:
1683  {
1684  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
1685  return SCIP_PARAMETERUNKNOWN;
1686  }
1687  }
1688 
1689  return SCIP_OKAY;
1690 }
1691 
1692 /** sets floating point parameter of NLP
1693  *
1694  * input:
1695  * - nlpi datastructure for solver interface
1696  * - problem datastructure for problem instance, can be NULL only if type == SCIP_NLPPAR_INFINITY
1697  * - type parameter number
1698  * - dval parameter value
1699  */
1700 static
1701 SCIP_DECL_NLPISETREALPAR(nlpiSetRealParIpopt)
1702 {
1703  assert(nlpi != NULL);
1704  assert(type == SCIP_NLPPAR_INFINITY || problem != NULL);
1705  assert(type == SCIP_NLPPAR_INFINITY || IsValid(problem->ipopt));
1706 
1707  switch( type )
1708  {
1710  {
1711  SCIPerrorMessage("from scratch parameter is of type int.\n");
1712  return SCIP_PARAMETERWRONGTYPE;
1713  }
1714 
1715  case SCIP_NLPPAR_VERBLEVEL:
1716  {
1717  SCIPerrorMessage("verbosity level parameter is of type int.\n");
1718  return SCIP_PARAMETERWRONGTYPE;
1719  }
1720 
1721  case SCIP_NLPPAR_FEASTOL:
1722  {
1723  if( dval >= 0 )
1724  {
1725  setFeastol(problem, dval);
1726  }
1727  else
1728  {
1729  SCIPerrorMessage("Value %g for parameter feasibility tolerance is negative\n", dval);
1730  return SCIP_PARAMETERWRONGVAL;
1731  }
1732  break;
1733  }
1734 
1735  case SCIP_NLPPAR_RELOBJTOL:
1736  {
1737  if( dval >= 0 )
1738  {
1739  setOpttol(problem, dval);
1740  }
1741  else
1742  {
1743  SCIPerrorMessage("Value %g for parameter relative objective tolerance is negative\n", dval);
1744  return SCIP_PARAMETERWRONGVAL;
1745  }
1746  break;
1747  }
1748 
1749  case SCIP_NLPPAR_LOBJLIM:
1750  {
1751  SCIP_NLPIDATA* data;
1752 
1753  data = SCIPnlpiGetData(nlpi);
1754  assert(data != NULL);
1755 
1756  SCIPmessagePrintWarning(data->messagehdlr, "Parameter lower objective limit not supported by Ipopt interface yet. Ignored.\n");
1757  break;
1758  }
1759 
1760  case SCIP_NLPPAR_INFINITY:
1761  {
1762  if( dval < 0.0 )
1763  return SCIP_PARAMETERWRONGVAL;
1764  if( problem )
1765  {
1766  problem->ipopt->Options()->SetNumericValue("diverging_iterates_tol", dval);
1767  problem->ipopt->Options()->SetNumericValue("nlp_lower_bound_inf", -dval);
1768  problem->ipopt->Options()->SetNumericValue("nlp_upper_bound_inf", dval);
1769  SCIPnlpiOracleSetInfinity(problem->oracle, dval);
1770  }
1771  else
1772  {
1773  SCIP_NLPIDATA* data = SCIPnlpiGetData(nlpi);
1774  assert(data != NULL);
1775  data->infinity = dval;
1776  }
1777  break;
1778  }
1779 
1780  case SCIP_NLPPAR_ITLIM:
1781  {
1782  SCIPerrorMessage("iteration limit parameter is of type int.\n");
1783  return SCIP_PARAMETERWRONGTYPE;
1784  }
1785 
1786  case SCIP_NLPPAR_TILIM:
1787  {
1788  if( dval >= 0 )
1789  {
1790  problem->ipopt->Options()->SetNumericValue("max_cpu_time", dval);
1791  }
1792  else
1793  {
1794  SCIPerrorMessage("Value %g for parameter time limit is negative\n", dval);
1795  return SCIP_PARAMETERWRONGVAL;
1796  }
1797  break;
1798  }
1799 
1800  case SCIP_NLPPAR_OPTFILE:
1801  {
1802  SCIPerrorMessage("option file parameter is of type string.\n");
1803  return SCIP_PARAMETERWRONGTYPE;
1804  }
1805 
1806  case SCIP_NLPPAR_FASTFAIL:
1807  {
1808  SCIPerrorMessage("fastfail parameter is of type int.\n");
1809  return SCIP_PARAMETERWRONGTYPE;
1810  }
1811 
1812  default:
1813  {
1814  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
1815  return SCIP_PARAMETERUNKNOWN;
1816  }
1817  }
1818 
1819  return SCIP_OKAY;
1820 }
1821 
1822 /** gets string parameter of NLP
1823  *
1824  * input:
1825  * - nlpi NLP interface structure
1826  * - problem datastructure for problem instance
1827  * - type parameter number
1828  * - sval pointer to store the string value, the user must not modify the string
1829  *
1830  * output:
1831  * - sval parameter value
1832  */
1833 static
1834 SCIP_DECL_NLPIGETSTRINGPAR( nlpiGetStringParIpopt )
1835 {
1836  assert(nlpi != NULL);
1837  assert(problem != NULL);
1838 
1839  switch( type )
1840  {
1842  {
1843  SCIPerrorMessage("from scratch parameter is of type int.\n");
1844  return SCIP_PARAMETERWRONGTYPE;
1845  }
1846 
1847  case SCIP_NLPPAR_VERBLEVEL:
1848  {
1849  SCIPerrorMessage("verbosity level parameter is of type int.\n");
1850  return SCIP_PARAMETERWRONGTYPE;
1851  }
1852 
1853  case SCIP_NLPPAR_FEASTOL:
1854  {
1855  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
1856  return SCIP_PARAMETERWRONGTYPE;
1857  }
1858 
1859  case SCIP_NLPPAR_RELOBJTOL:
1860  {
1861  SCIPerrorMessage("objective tolerance parameter is of type real.\n");
1862  return SCIP_PARAMETERWRONGTYPE;
1863  }
1864 
1865  case SCIP_NLPPAR_LOBJLIM:
1866  {
1867  SCIPerrorMessage("objective limit parameter is of type real.\n");
1868  return SCIP_PARAMETERWRONGTYPE;
1869  }
1870 
1871  case SCIP_NLPPAR_INFINITY:
1872  {
1873  SCIPerrorMessage("infinity parameter is of type real.\n");
1874  return SCIP_PARAMETERWRONGTYPE;
1875  }
1876 
1877  case SCIP_NLPPAR_ITLIM:
1878  {
1879  SCIPerrorMessage("iteration limit parameter is of type int.\n");
1880  return SCIP_PARAMETERWRONGTYPE;
1881  }
1882 
1883  case SCIP_NLPPAR_TILIM:
1884  {
1885  SCIPerrorMessage("time limit parameter is of type real.\n");
1886  return SCIP_PARAMETERWRONGTYPE;
1887  }
1888 
1889  case SCIP_NLPPAR_OPTFILE:
1890  {
1891  if( !problem->optfile.empty() )
1892  *sval = problem->optfile.c_str();
1893  else
1894  *sval = NULL;
1895  return SCIP_OKAY;
1896  }
1897 
1898  case SCIP_NLPPAR_FASTFAIL:
1899  {
1900  SCIPerrorMessage("fastfail parameter is of type int.\n");
1901  return SCIP_PARAMETERWRONGTYPE;
1902  }
1903 
1904  default:
1905  {
1906  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
1907  return SCIP_PARAMETERUNKNOWN;
1908  }
1909  }
1910 
1911  return SCIP_OKAY;
1912 }
1913 
1914 /** sets string parameter of NLP
1915  *
1916  * input:
1917  * - nlpi NLP interface structure
1918  * - problem datastructure for problem instance
1919  * - type parameter number
1920  * - sval parameter value
1921  */
1922 static
1923 SCIP_DECL_NLPISETSTRINGPAR( nlpiSetStringParIpopt )
1924 {
1925  assert(nlpi != NULL);
1926  assert(problem != NULL);
1927 
1928  switch( type )
1929  {
1931  {
1932  SCIPerrorMessage("from scratch parameter is of type int.\n");
1933  return SCIP_PARAMETERWRONGTYPE;
1934  }
1935 
1936  case SCIP_NLPPAR_VERBLEVEL:
1937  {
1938  SCIPerrorMessage("verbosity level parameter is of type int.\n");
1939  return SCIP_PARAMETERWRONGTYPE;
1940  }
1941 
1942  case SCIP_NLPPAR_FEASTOL:
1943  {
1944  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
1945  return SCIP_PARAMETERWRONGTYPE;
1946  }
1947 
1948  case SCIP_NLPPAR_RELOBJTOL:
1949  {
1950  SCIPerrorMessage("objective tolerance parameter is of type real.\n");
1951  return SCIP_PARAMETERWRONGTYPE;
1952  }
1953 
1954  case SCIP_NLPPAR_LOBJLIM:
1955  {
1956  SCIPerrorMessage("objective limit parameter is of type real.\n");
1957  return SCIP_PARAMETERWRONGTYPE;
1958  }
1959 
1960  case SCIP_NLPPAR_INFINITY:
1961  {
1962  SCIPerrorMessage("infinity parameter is of type real.\n");
1963  return SCIP_PARAMETERWRONGTYPE;
1964  }
1965 
1966  case SCIP_NLPPAR_ITLIM:
1967  {
1968  SCIPerrorMessage("iteration limit parameter is of type int.\n");
1969  return SCIP_PARAMETERWRONGTYPE;
1970  }
1971 
1972  case SCIP_NLPPAR_TILIM:
1973  {
1974  SCIPerrorMessage("time limit parameter is of type real.\n");
1975  return SCIP_PARAMETERWRONGTYPE;
1976  }
1977 
1978  case SCIP_NLPPAR_OPTFILE:
1979  {
1980  if( sval != NULL )
1981  problem->optfile = sval;
1982  else
1983  problem->optfile.clear();
1984 
1985  if( problem->ipopt->Initialize(problem->optfile) != Solve_Succeeded )
1986  {
1987  SCIPerrorMessage("Error initializing Ipopt using optionfile \"%s\"\n", problem->optfile.c_str());
1988  return SCIP_ERROR;
1989  }
1990  problem->ipopt->Options()->GetBoolValue("store_intermediate", problem->storeintermediate, "");
1991  problem->firstrun = TRUE;
1992 
1993  return SCIP_OKAY;
1994  }
1995 
1996  case SCIP_NLPPAR_FASTFAIL:
1997  {
1998  SCIPerrorMessage("fastfail parameter is of type int.\n");
1999  return SCIP_PARAMETERWRONGTYPE;
2000  }
2001 
2002  default:
2003  {
2004  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
2005  return SCIP_PARAMETERUNKNOWN;
2006  }
2007  }
2008 
2009  return SCIP_OKAY;
2010 }
2011 
2012 /** sets message handler for message output
2013  *
2014  * input:
2015  * - nlpi NLP interface structure
2016  * - messagehdlr SCIP message handler, or NULL to suppress all output
2017  */
2018 static
2019 SCIP_DECL_NLPISETMESSAGEHDLR( nlpiSetMessageHdlrIpopt )
2020 {
2021  SCIP_NLPIDATA* nlpidata;
2022 
2023  assert(nlpi != NULL);
2024 
2025  nlpidata = SCIPnlpiGetData(nlpi);
2026  assert(nlpidata != NULL);
2027 
2028  nlpidata->messagehdlr = messagehdlr;
2029 
2030  return SCIP_OKAY; /*lint !e527*/
2031 } /*lint !e715*/
2032 
2033 /** create solver interface for Ipopt solver */
2035  BMS_BLKMEM* blkmem, /**< block memory data structure */
2036  SCIP_NLPI** nlpi /**< pointer to buffer for nlpi address */
2037  )
2038 {
2039  SCIP_NLPIDATA* nlpidata;
2040 
2041  assert(blkmem != NULL);
2042  assert(nlpi != NULL);
2043 
2044  SCIP_ALLOC( nlpidata = new SCIP_NLPIDATA(blkmem) );
2045 
2046  SCIP_CALL( SCIPnlpiCreate(nlpi,
2048  nlpiCopyIpopt, nlpiFreeIpopt, nlpiGetSolverPointerIpopt,
2049  nlpiCreateProblemIpopt, nlpiFreeProblemIpopt, nlpiGetProblemPointerIpopt,
2050  nlpiAddVarsIpopt, nlpiAddConstraintsIpopt, nlpiSetObjectiveIpopt,
2051  nlpiChgVarBoundsIpopt, nlpiChgConsSidesIpopt, nlpiDelVarSetIpopt, nlpiDelConstraintSetIpopt,
2052  nlpiChgLinearCoefsIpopt, nlpiChgQuadraticCoefsIpopt, nlpiChgExprtreeIpopt, nlpiChgNonlinCoefIpopt,
2053  nlpiChgObjConstantIpopt, nlpiSetInitialGuessIpopt, nlpiSolveIpopt, nlpiGetSolstatIpopt, nlpiGetTermstatIpopt,
2054  nlpiGetSolutionIpopt, nlpiGetStatisticsIpopt,
2055  nlpiGetWarmstartSizeIpopt, nlpiGetWarmstartMemoIpopt, nlpiSetWarmstartMemoIpopt,
2056  nlpiGetIntParIpopt, nlpiSetIntParIpopt, nlpiGetRealParIpopt, nlpiSetRealParIpopt,
2057  nlpiGetStringParIpopt, nlpiSetStringParIpopt, nlpiSetMessageHdlrIpopt,
2058  nlpidata) );
2059 
2060  return SCIP_OKAY;
2061 }
2062 
2063 /** gets string that identifies Ipopt (version number) */
2064 const char* SCIPgetSolverNameIpopt(void)
2065 {
2066  return "Ipopt " IPOPT_VERSION;
2067 }
2068 
2069 /** gets string that describes Ipopt */
2070 const char* SCIPgetSolverDescIpopt(void)
2071 {
2072  return "Interior Point Optimizer developed by A. Waechter et.al. (www.coin-or.org/Ipopt)";
2073 }
2074 
2075 /** returns whether Ipopt is available, i.e., whether it has been linked in */
2077 {
2078  return TRUE;
2079 }
2080 
2081 /** gives a pointer to the IpoptApplication object stored in Ipopt-NLPI's NLPI problem data structure */
2083  SCIP_NLPIPROBLEM* nlpiproblem /**< NLP problem of Ipopt-NLPI */
2084  )
2085 {
2086  assert(nlpiproblem != NULL);
2087 
2088  return (void*)GetRawPtr(nlpiproblem->ipopt);
2089 }
2090 
2091 /** gives a pointer to the NLPIORACLE object stored in Ipopt-NLPI's NLPI problem data structure */
2093  SCIP_NLPIPROBLEM* nlpiproblem /**< NLP problem of Ipopt-NLPI */
2094  )
2095 {
2096  assert(nlpiproblem != NULL);
2097 
2098  return nlpiproblem->oracle;
2099 }
2100 
2101 /** sets modified default settings that are used when setting up an Ipopt problem
2102  *
2103  * Do not forget to add a newline after the last option in optionsstring.
2104  */
2106  SCIP_NLPI* nlpi, /**< Ipopt NLP interface */
2107  const char* optionsstring /**< string with options as in Ipopt options file */
2108  )
2109 {
2110  SCIP_NLPIDATA* data;
2111 
2112  assert(nlpi != NULL);
2113 
2114  data = SCIPnlpiGetData(nlpi);
2115  assert(data != NULL);
2116 
2117  data->defoptions = optionsstring;
2118 }
2119 
2120 /** Method to return some info about the nlp */
2121 bool ScipNLP::get_nlp_info(
2122  Index& n, /**< place to store number of variables */
2123  Index& m, /**< place to store number of constraints */
2124  Index& nnz_jac_g, /**< place to store number of nonzeros in jacobian */
2125  Index& nnz_h_lag, /**< place to store number of nonzeros in hessian */
2126  IndexStyleEnum& index_style /**< place to store used index style (0-based or 1-based) */
2127  )
2128 {
2129  const int* offset;
2130  SCIP_RETCODE retcode;
2131 
2132  assert(nlpiproblem != NULL);
2133  assert(nlpiproblem->oracle != NULL);
2134 
2135  n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
2136  m = SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle);
2137 
2138  retcode = SCIPnlpiOracleGetJacobianSparsity(nlpiproblem->oracle, &offset, NULL);
2139  if( retcode != SCIP_OKAY )
2140  return false;
2141  assert(offset != NULL);
2142  nnz_jac_g = offset[m];
2143 
2144  if( !approxhessian )
2145  {
2146  retcode = SCIPnlpiOracleGetHessianLagSparsity(nlpiproblem->oracle, &offset, NULL);
2147  if( retcode != SCIP_OKAY )
2148  return false;
2149  assert(offset != NULL);
2150  nnz_h_lag = offset[n];
2151  }
2152  else
2153  {
2154  nnz_h_lag = 0;
2155  }
2156 
2157  index_style = TNLP::C_STYLE;
2158 
2159  return true;
2160 }
2161 
2162 /** Method to return the bounds for my problem */
2163 bool ScipNLP::get_bounds_info(
2164  Index n, /**< number of variables */
2165  Number* x_l, /**< buffer to store lower bounds on variables */
2166  Number* x_u, /**< buffer to store upper bounds on variables */
2167  Index m, /**< number of constraints */
2168  Number* g_l, /**< buffer to store lower bounds on constraints */
2169  Number* g_u /**< buffer to store lower bounds on constraints */
2170  )
2171 {
2172  assert(nlpiproblem != NULL);
2173  assert(nlpiproblem->oracle != NULL);
2174 
2175  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2176  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2177 
2178  assert(SCIPnlpiOracleGetVarLbs(nlpiproblem->oracle) != NULL);
2179  assert(SCIPnlpiOracleGetVarUbs(nlpiproblem->oracle) != NULL);
2180 
2181  BMScopyMemoryArray(x_l, SCIPnlpiOracleGetVarLbs(nlpiproblem->oracle), n);
2182  BMScopyMemoryArray(x_u, SCIPnlpiOracleGetVarUbs(nlpiproblem->oracle), n);
2183 #ifndef NDEBUG
2184  for( int i = 0; i < n; ++i )
2185  assert(x_l[i] <= x_u[i]);
2186 #endif
2187 
2188  for( int i = 0; i < m; ++i )
2189  {
2190  g_l[i] = SCIPnlpiOracleGetConstraintLhs(nlpiproblem->oracle, i);
2191  g_u[i] = SCIPnlpiOracleGetConstraintRhs(nlpiproblem->oracle, i);
2192  assert(g_l[i] <= g_u[i]);
2193  }
2194 
2195  return true;
2196 }
2197 
2198 /** Method to return the starting point for the algorithm */
2199 bool ScipNLP::get_starting_point(
2200  Index n, /**< number of variables */
2201  bool init_x, /**< whether initial values for primal values are requested */
2202  Number* x, /**< buffer to store initial primal values */
2203  bool init_z, /**< whether initial values for dual values of variable bounds are requested */
2204  Number* z_L, /**< buffer to store dual values for variable lower bounds */
2205  Number* z_U, /**< buffer to store dual values for variable upper bounds */
2206  Index m, /**< number of constraints */
2207  bool init_lambda, /**< whether initial values for dual values of constraints are required */
2208  Number* lambda /**< buffer to store dual values of constraints */
2209  )
2210 {
2211  assert(nlpiproblem != NULL);
2212  assert(nlpiproblem->oracle != NULL);
2213 
2214  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2215  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2216 
2217  if( init_x )
2218  {
2219  if( nlpiproblem->initguess )
2220  {
2221  BMScopyMemoryArray(x, nlpiproblem->initguess, n);
2222  }
2223  else
2224  {
2225  SCIP_Real lb, ub;
2226 
2227  SCIPdebugMessage("Ipopt started without intial primal values; make up starting guess by projecting 0 onto variable bounds\n");
2228 
2229  for( int i = 0; i < n; ++i )
2230  {
2231  lb = SCIPnlpiOracleGetVarLbs(nlpiproblem->oracle)[i];
2232  ub = SCIPnlpiOracleGetVarUbs(nlpiproblem->oracle)[i];
2233  if( lb > 0.0 )
2234  x[i] = SCIPrandomGetReal(randnumgen, lb, lb + MAXPERTURB*MIN(1.0, ub-lb));
2235  else if( ub < 0.0 )
2236  x[i] = SCIPrandomGetReal(randnumgen, ub - MAXPERTURB*MIN(1.0, ub-lb), ub);
2237  else
2238  x[i] = SCIPrandomGetReal(randnumgen, MAX(lb, -MAXPERTURB*MIN(1.0, ub-lb)), MIN(ub, MAXPERTURB*MIN(1.0, ub-lb)));
2239  }
2240  }
2241  }
2242  if( init_z || init_lambda )
2243  return false;
2244 
2245  return true;
2246 }
2247 
2248 /** Method to return the variables linearity. */
2249 bool ScipNLP::get_variables_linearity(
2250  Index n, /**< number of variables */
2251  LinearityType* var_types /**< buffer to store linearity types of variables */
2252  )
2253 {
2254  assert(nlpiproblem != NULL);
2255  assert(nlpiproblem->oracle != NULL);
2256 
2257  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2258 
2259  for( int i = 0; i < n; ++i )
2260  var_types[i] = (SCIPnlpiOracleGetVarDegree(nlpiproblem->oracle, i) <= 1 ? LINEAR : NON_LINEAR);
2261 
2262  return true;
2263 }
2264 
2265 /** Method to return the constraint linearity. */
2266 bool ScipNLP::get_constraints_linearity(
2267  Index m, /**< number of constraints */
2268  LinearityType* const_types /**< buffer to store linearity types of constraints */
2269  )
2270 {
2271  int i;
2272 
2273  assert(nlpiproblem != NULL);
2274  assert(nlpiproblem->oracle != NULL);
2275 
2276  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2277 
2278  for( i = 0; i < m; ++i )
2279  const_types[i] = (SCIPnlpiOracleGetConstraintDegree(nlpiproblem->oracle, i) <= 1 ? LINEAR : NON_LINEAR);
2280 
2281  return true;
2282 }
2283 
2284 /** Method to return the number of nonlinear variables. */
2285 Index ScipNLP::get_number_of_nonlinear_variables()
2286 {
2287  int count;
2288  int n;
2289 
2290  assert(nlpiproblem != NULL);
2291  assert(nlpiproblem->oracle != NULL);
2292 
2293  n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
2294 
2295  count = 0;
2296  for( int i = 0; i < n; ++i )
2297  if (SCIPnlpiOracleGetVarDegree(nlpiproblem->oracle, i) > 1)
2298  ++count;
2299 
2300  return count;
2301 }
2302 
2303 /** Method to return the indices of the nonlinear variables */
2304 bool ScipNLP::get_list_of_nonlinear_variables(
2305  Index num_nonlin_vars, /**< number of nonlinear variables */
2306  Index* pos_nonlin_vars /**< array to fill with indices of nonlinear variables */
2307  )
2308 {
2309  int count;
2310  int n;
2311 
2312  assert(nlpiproblem != NULL);
2313  assert(nlpiproblem->oracle != NULL);
2314 
2315  n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
2316 
2317  count = 0;
2318  for( int i = 0; i < n; ++i )
2319  if (SCIPnlpiOracleGetVarDegree(nlpiproblem->oracle, i) > 1)
2320  {
2321  assert(count < num_nonlin_vars);
2322  pos_nonlin_vars[count++] = i;
2323  }
2324 
2325  assert(count == num_nonlin_vars);
2326 
2327  return true;
2328 }
2329 
2330 /** Method to return metadata about variables and constraints */
2331 bool ScipNLP::get_var_con_metadata(
2332  Index n, /**< number of variables */
2333  StringMetaDataMapType& var_string_md, /**< variable meta data of string type */
2334  IntegerMetaDataMapType& var_integer_md,/**< variable meta data of integer type */
2335  NumericMetaDataMapType& var_numeric_md,/**< variable meta data of numeric type */
2336  Index m, /**< number of constraints */
2337  StringMetaDataMapType& con_string_md, /**< constraint meta data of string type */
2338  IntegerMetaDataMapType& con_integer_md,/**< constraint meta data of integer type */
2339  NumericMetaDataMapType& con_numeric_md /**< constraint meta data of numeric type */
2340  )
2341 {
2342  assert(nlpiproblem != NULL);
2343  assert(nlpiproblem->oracle != NULL);
2344  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2345  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2346 
2347  char** varnames = SCIPnlpiOracleGetVarNames(nlpiproblem->oracle);
2348  if( varnames != NULL )
2349  {
2350  std::vector<std::string>& varnamesvec(var_string_md["idx_names"]);
2351  varnamesvec.reserve(n);
2352  for( int i = 0; i < n; ++i )
2353  {
2354  if( varnames[i] != NULL )
2355  {
2356  varnamesvec.push_back(varnames[i]);
2357  }
2358  else
2359  {
2360  char buffer[20];
2361  sprintf(buffer, "nlpivar%8d", i);
2362  varnamesvec.push_back(buffer);
2363  }
2364  }
2365  }
2366 
2367  std::vector<std::string>& consnamesvec(con_string_md["idx_names"]);
2368  consnamesvec.reserve(m);
2369  for( int i = 0; i < m; ++i )
2370  {
2371  if( SCIPnlpiOracleGetConstraintName(nlpiproblem->oracle, i) != NULL )
2372  {
2373  consnamesvec.push_back(SCIPnlpiOracleGetConstraintName(nlpiproblem->oracle, i));
2374  }
2375  else
2376  {
2377  char buffer[20];
2378  sprintf(buffer, "nlpicons%8d", i);
2379  consnamesvec.push_back(buffer);
2380  }
2381  }
2382 
2383  return true;
2384 }
2385 
2386 /** Method to return the objective value */
2387 bool ScipNLP::eval_f(
2388  Index n, /**< number of variables */
2389  const Number* x, /**< point to evaluate */
2390  bool new_x, /**< whether some function evaluation method has been called for this point before */
2391  Number& obj_value /**< place to store objective function value */
2392  )
2393 {
2394  assert(nlpiproblem != NULL);
2395  assert(nlpiproblem->oracle != NULL);
2396 
2397  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2398 
2399  return SCIPnlpiOracleEvalObjectiveValue(nlpiproblem->oracle, x, &obj_value) == SCIP_OKAY;
2400 }
2401 
2402 /** Method to return the gradient of the objective */
2403 bool ScipNLP::eval_grad_f(
2404  Index n, /**< number of variables */
2405  const Number* x, /**< point to evaluate */
2406  bool new_x, /**< whether some function evaluation method has been called for this point before */
2407  Number* grad_f /**< buffer to store objective gradient */
2408  )
2409 {
2410  SCIP_Real dummy;
2411 
2412  assert(nlpiproblem != NULL);
2413  assert(nlpiproblem->oracle != NULL);
2414 
2415  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2416 
2417  return SCIPnlpiOracleEvalObjectiveGradient(nlpiproblem->oracle, x, TRUE, &dummy, grad_f) == SCIP_OKAY;
2418 }
2419 
2420 /** Method to return the constraint residuals */
2421 bool ScipNLP::eval_g(
2422  Index n, /**< number of variables */
2423  const Number* x, /**< point to evaluate */
2424  bool new_x, /**< whether some function evaluation method has been called for this point before */
2425  Index m, /**< number of constraints */
2426  Number* g /**< buffer to store constraint function values */
2427  )
2428 {
2429  assert(nlpiproblem != NULL);
2430  assert(nlpiproblem->oracle != NULL);
2431 
2432  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2433 
2434  return SCIPnlpiOracleEvalConstraintValues(nlpiproblem->oracle, x, g) == SCIP_OKAY;
2435 }
2436 
2437 /** Method to return:
2438  * 1) The structure of the jacobian (if "values" is NULL)
2439  * 2) The values of the jacobian (if "values" is not NULL)
2440  */
2441 bool ScipNLP::eval_jac_g(
2442  Index n, /**< number of variables */
2443  const Number* x, /**< point to evaluate */
2444  bool new_x, /**< whether some function evaluation method has been called for this point before */
2445  Index m, /**< number of constraints */
2446  Index nele_jac, /**< number of nonzero entries in jacobian */
2447  Index* iRow, /**< buffer to store row indices of nonzero jacobian entries, or NULL if values are requested */
2448  Index* jCol, /**< buffer to store column indices of nonzero jacobian entries, or NULL if values are requested */
2449  Number* values /**< buffer to store values of nonzero jacobian entries, or NULL if structure is requested */
2450  )
2451 {
2452  assert(nlpiproblem != NULL);
2453  assert(nlpiproblem->oracle != NULL);
2454 
2455  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2456  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2457 
2458  if( values == NULL )
2459  { /* Ipopt wants to know sparsity structure */
2460  const int* jacoffset;
2461  const int* jaccol;
2462  int j;
2463  int i;
2464 
2465  assert(iRow != NULL);
2466  assert(jCol != NULL);
2467 
2468  if( SCIPnlpiOracleGetJacobianSparsity(nlpiproblem->oracle, &jacoffset, &jaccol) != SCIP_OKAY )
2469  return false;
2470 
2471  assert(jacoffset[0] == 0);
2472  assert(jacoffset[m] == nele_jac);
2473  j = jacoffset[0];
2474  for( i = 0; i < m; ++i )
2475  for( ; j < jacoffset[i+1]; ++j )
2476  iRow[j] = i;
2477 
2478  BMScopyMemoryArray(jCol, jaccol, nele_jac);
2479  }
2480  else
2481  {
2482  if( SCIPnlpiOracleEvalJacobian(nlpiproblem->oracle, x, TRUE, NULL, values) != SCIP_OKAY )
2483  return false;
2484  }
2485 
2486  return true;
2487 }
2488 
2489 /** Method to return:
2490  * 1) The structure of the hessian of the lagrangian (if "values" is NULL)
2491  * 2) The values of the hessian of the lagrangian (if "values" is not NULL)
2492  */
2493 bool ScipNLP::eval_h(
2494  Index n, /**< number of variables */
2495  const Number* x, /**< point to evaluate */
2496  bool new_x, /**< whether some function evaluation method has been called for this point before */
2497  Number obj_factor, /**< weight for objective function */
2498  Index m, /**< number of constraints */
2499  const Number* lambda, /**< weights for constraint functions */
2500  bool new_lambda, /**< whether the hessian has been evaluated for these values of lambda before */
2501  Index nele_hess, /**< number of nonzero entries in hessian */
2502  Index* iRow, /**< buffer to store row indices of nonzero hessian entries, or NULL if values are requested */
2503  Index* jCol, /**< buffer to store column indices of nonzero hessian entries, or NULL if values are requested */
2504  Number* values /**< buffer to store values of nonzero hessian entries, or NULL if structure is requested */
2505  )
2506 {
2507  assert(nlpiproblem != NULL);
2508  assert(nlpiproblem->oracle != NULL);
2509 
2510  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2511  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2512 
2513  if( values == NULL )
2514  { /* Ipopt wants to know sparsity structure */
2515  const int* heslagoffset;
2516  const int* heslagcol;
2517  int j;
2518  int i;
2519 
2520  assert(iRow != NULL);
2521  assert(jCol != NULL);
2522 
2523  if( SCIPnlpiOracleGetHessianLagSparsity(nlpiproblem->oracle, &heslagoffset, &heslagcol) != SCIP_OKAY )
2524  return false;
2525 
2526  assert(heslagoffset[0] == 0);
2527  assert(heslagoffset[n] == nele_hess);
2528  j = heslagoffset[0];
2529  for( i = 0; i < n; ++i )
2530  for( ; j < heslagoffset[i+1]; ++j )
2531  iRow[j] = i;
2532 
2533  BMScopyMemoryArray(jCol, heslagcol, nele_hess);
2534  }
2535  else
2536  {
2537  if( SCIPnlpiOracleEvalHessianLag(nlpiproblem->oracle, x, TRUE, obj_factor, lambda, values) != SCIP_OKAY )
2538  return false;
2539  }
2540 
2541  return true;
2542 }
2543 
2544 /** Method called by the solver at each iteration.
2545  *
2546  * Checks whether Ctrl-C was hit.
2547  */
2548 bool ScipNLP::intermediate_callback(
2549  AlgorithmMode mode, /**< current mode of algorithm */
2550  Index iter, /**< current iteration number */
2551  Number obj_value, /**< current objective value */
2552  Number inf_pr, /**< current primal infeasibility */
2553  Number inf_du, /**< current dual infeasibility */
2554  Number mu, /**< current barrier parameter */
2555  Number d_norm, /**< current gradient norm */
2556  Number regularization_size,/**< current size of regularization */
2557  Number alpha_du, /**< current dual alpha */
2558  Number alpha_pr, /**< current primal alpha */
2559  Index ls_trials, /**< current number of linesearch trials */
2560  const IpoptData* ip_data, /**< pointer to Ipopt Data */
2561  IpoptCalculatedQuantities* ip_cq /**< pointer to current calculated quantities */
2562  )
2563 {
2564  if( nlpiproblem->storeintermediate && mode == RegularMode && inf_pr < nlpiproblem->lastsolinfeas )
2565  {
2566  Ipopt::TNLPAdapter* tnlp_adapter;
2567 
2568  tnlp_adapter = NULL;
2569  if( ip_cq != NULL )
2570  {
2571  Ipopt::OrigIpoptNLP* orignlp;
2572 
2573  orignlp = dynamic_cast<OrigIpoptNLP*>(GetRawPtr(ip_cq->GetIpoptNLP()));
2574  if( orignlp != NULL )
2575  tnlp_adapter = dynamic_cast<TNLPAdapter*>(GetRawPtr(orignlp->nlp()));
2576  }
2577 
2578  if( tnlp_adapter != NULL && ip_data != NULL && IsValid(ip_data->curr()) )
2579  {
2580  SCIPdebugMessage("update lastsol: inf_pr old = %g -> new = %g\n", nlpiproblem->lastsolinfeas, inf_pr);
2581 
2582  if( nlpiproblem->lastsolprimals == NULL )
2583  {
2584  assert(nlpiproblem->lastsoldualcons == NULL);
2585  assert(nlpiproblem->lastsoldualvarlb == NULL);
2586  assert(nlpiproblem->lastsoldualvarub == NULL);
2587  if( BMSallocMemoryArray(&nlpiproblem->lastsolprimals, SCIPnlpiOracleGetNVars(nlpiproblem->oracle)) == NULL ||
2588  BMSallocMemoryArray(&nlpiproblem->lastsoldualcons, SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle)) == NULL ||
2589  BMSallocMemoryArray(&nlpiproblem->lastsoldualvarlb, SCIPnlpiOracleGetNVars(nlpiproblem->oracle)) == NULL ||
2590  BMSallocMemoryArray(&nlpiproblem->lastsoldualvarub, SCIPnlpiOracleGetNVars(nlpiproblem->oracle)) == NULL )
2591  {
2592  SCIPerrorMessage("out-of-memory in ScipNLP::intermediate_callback()\n");
2593  return TRUE;
2594  }
2595  }
2596 
2597  assert(IsValid(ip_data->curr()->x()));
2598  tnlp_adapter->ResortX(*ip_data->curr()->x(), nlpiproblem->lastsolprimals);
2599  nlpiproblem->lastsolinfeas = inf_pr;
2600 
2601  assert(IsValid(ip_data->curr()->y_c()));
2602  assert(IsValid(ip_data->curr()->y_d()));
2603  tnlp_adapter->ResortG(*ip_data->curr()->y_c(), *ip_data->curr()->y_d(), nlpiproblem->lastsoldualcons);
2604 
2605  // need to clear arrays first because ResortBnds only sets values for non-fixed variables
2606  BMSclearMemoryArray(nlpiproblem->lastsoldualvarlb, SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2607  BMSclearMemoryArray(nlpiproblem->lastsoldualvarub, SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2608  assert(IsValid(ip_data->curr()->z_L()));
2609  assert(IsValid(ip_data->curr()->z_U()));
2610  tnlp_adapter->ResortBnds(*ip_data->curr()->z_L(), nlpiproblem->lastsoldualvarlb, *ip_data->curr()->z_U(), nlpiproblem->lastsoldualvarub);
2611 
2612  }
2613  }
2614 
2615  /* do convergence test if fastfail is enabled */
2616  if( nlpiproblem->fastfail )
2617  {
2618  int i;
2619 
2620  if( iter == 0 )
2621  {
2622  conv_lastrestoiter = -1;
2623  }
2624  else if( mode == RestorationPhaseMode )
2625  {
2626  conv_lastrestoiter = iter;
2627  }
2628  else if( conv_lastrestoiter == iter-1 )
2629  {
2630  /* just switched back from restoration mode, reset dual reduction targets */
2631  for( i = 0; i < convcheck_nchecks; ++i )
2632  conv_dutarget[i] = convcheck_minred[i] * inf_du;
2633  }
2634 
2635  if( iter == convcheck_startiter )
2636  {
2637  /* define initial targets and iteration limits */
2638  for( i = 0; i < convcheck_nchecks; ++i )
2639  {
2640  conv_prtarget[i] = convcheck_minred[i] * inf_pr;
2641  conv_dutarget[i] = convcheck_minred[i] * inf_du;
2642  conv_iterlim[i] = iter + convcheck_maxiter[i];
2643  }
2644  }
2645  else if( iter > convcheck_startiter )
2646  {
2647  /* check if we should stop */
2648  for( i = 0; i < convcheck_nchecks; ++i )
2649  {
2650  if( inf_pr <= conv_prtarget[i] )
2651  {
2652  /* sufficient reduction w.r.t. primal infeasibility target
2653  * reset target w.r.t. current infeasibilities
2654  */
2655  conv_prtarget[i] = convcheck_minred[i] * inf_pr;
2656  conv_dutarget[i] = convcheck_minred[i] * inf_du;
2657  conv_iterlim[i] = iter + convcheck_maxiter[i];
2658  }
2659  else if( iter >= conv_iterlim[i] )
2660  {
2661  /* we hit a limit, should we really stop? */
2662  SCIPdebugMessage("convcheck %d: inf_pr = %e > target %e; inf_du = %e target %e: ",
2663  i, inf_pr, conv_prtarget[i], inf_du, conv_dutarget[i]);
2664  if( mode == RegularMode && iter <= conv_lastrestoiter + convcheck_startiter )
2665  {
2666  /* if we returned from feasibility restoration recently, we allow some more iterations,
2667  * because Ipopt may go for optimality for some iterations, at the costs of infeasibility
2668  */
2669  SCIPdebugPrintf("continue, because restoration phase only %d iters ago\n", iter - conv_lastrestoiter);
2670  }
2671  else if( mode == RegularMode && inf_du <= conv_dutarget[i] && iter < conv_iterlim[i] + convcheck_maxiter[i] )
2672  {
2673  /* if dual reduction is sufficient, we allow for twice the number of iterations to reach primal infeas reduction */
2674  SCIPdebugPrintf("continue, because dual infeas. red. sufficient and only %d iters above limit\n", iter - conv_iterlim[i]);
2675  }
2676  else
2677  {
2678  SCIPdebugPrintf("abort\n");
2679  return false;
2680  }
2681  }
2682  }
2683  }
2684  }
2685 
2686  return (SCIPinterrupted() == FALSE);
2687 }
2688 
2689 /** This method is called when the algorithm is complete so the TNLP can store/write the solution. */
2690 void ScipNLP::finalize_solution(
2691  SolverReturn status, /**< solve and solution status */
2692  Index n, /**< number of variables */
2693  const Number* x, /**< primal solution values */
2694  const Number* z_L, /**< dual values of variable lower bounds */
2695  const Number* z_U, /**< dual values of variable upper bounds */
2696  Index m, /**< number of constraints */
2697  const Number* g, /**< values of constraints */
2698  const Number* lambda, /**< dual values of constraints */
2699  Number obj_value, /**< objective function value */
2700  const IpoptData* data, /**< pointer to Ipopt Data */
2701  IpoptCalculatedQuantities* cq /**< pointer to calculated quantities */
2702  )
2703 {
2704  assert(nlpiproblem != NULL);
2705  assert(nlpiproblem->oracle != NULL);
2706 
2707  assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2708  assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2709 
2710  bool check_feasibility = false; // whether we should check x for feasibility, if not NULL
2711  switch( status )
2712  {
2713  case SUCCESS:
2714  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_LOCOPT;
2715  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OKAY;
2716  assert(x != NULL);
2717  break;
2718 
2719  case STOP_AT_ACCEPTABLE_POINT:
2720  /* if stop at acceptable point, then dual infeasibility can be arbitrary large, so claim only feasibility */
2721  case FEASIBLE_POINT_FOUND:
2722  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_FEASIBLE;
2723  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OKAY;
2724  assert(x != NULL);
2725  break;
2726 
2727  case MAXITER_EXCEEDED:
2728  check_feasibility = true;
2729  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2730  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_ITLIM;
2731  break;
2732 
2733  case CPUTIME_EXCEEDED:
2734  check_feasibility = true;
2735  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2736  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_TILIM;
2737  break;
2738 
2739  case STOP_AT_TINY_STEP:
2740  case RESTORATION_FAILURE:
2741  case ERROR_IN_STEP_COMPUTATION:
2742  check_feasibility = true;
2743  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2744  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_NUMERR;
2745  break;
2746 
2747  case LOCAL_INFEASIBILITY:
2748  /* still check feasibility, since we let Ipopt solve with higher tolerance than actually required */
2749  check_feasibility = true;
2750  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_LOCINFEASIBLE;
2751  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OKAY;
2752  break;
2753 
2754  case DIVERGING_ITERATES:
2755  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNBOUNDED;
2756  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OKAY;
2757  break;
2758 
2759  case INVALID_NUMBER_DETECTED:
2760  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2761  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_EVALERR;
2762  break;
2763 
2764  case USER_REQUESTED_STOP:
2765  check_feasibility = true;
2766  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2767  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OKAY;
2768  break;
2769 
2770  case TOO_FEW_DEGREES_OF_FREEDOM:
2771  case INTERNAL_ERROR:
2772  case INVALID_OPTION:
2773  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2774  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OTHER;
2775  break;
2776 
2777  case OUT_OF_MEMORY:
2778  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2779  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_MEMERR;
2780  break;
2781 
2782  default:
2783  SCIPerrorMessage("Ipopt returned with unknown solution status %d\n", status);
2784  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2785  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_OTHER;
2786  break;
2787  }
2788 
2789  /* if Ipopt reports its solution as locally infeasible or we don't know feasibility, then report the intermediate point with lowest constraint violation, if available */
2790  if( (x == NULL || nlpiproblem->lastsolstat == SCIP_NLPSOLSTAT_LOCINFEASIBLE || nlpiproblem->lastsolstat == SCIP_NLPSOLSTAT_UNKNOWN) && nlpiproblem->lastsolinfeas != SCIP_INVALID )
2791  {
2792  /* if infeasibility of lastsol is not invalid, then lastsol values should exist */
2793  assert(nlpiproblem->lastsolprimals != NULL);
2794  assert(nlpiproblem->lastsoldualcons != NULL);
2795  assert(nlpiproblem->lastsoldualvarlb != NULL);
2796  assert(nlpiproblem->lastsoldualvarub != NULL);
2797 
2798  /* check if lastsol is feasible */
2799  Number constrvioltol;
2800  nlpiproblem->ipopt->Options()->GetNumericValue("acceptable_constr_viol_tol", constrvioltol, "");
2801  if( nlpiproblem->lastsolinfeas <= constrvioltol )
2802  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_FEASIBLE;
2803  else
2804  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2805 
2806  SCIPdebugMessage("drop Ipopt's final point and report intermediate locally %sfeasible solution with infeas %g instead (acceptable: %g)\n",
2807  nlpiproblem->lastsolstat == SCIP_NLPSOLSTAT_LOCINFEASIBLE ? "in" : "", nlpiproblem->lastsolinfeas, constrvioltol);
2808  }
2809  else
2810  {
2811  assert(x != NULL);
2812  assert(lambda != NULL);
2813  assert(z_L != NULL);
2814  assert(z_U != NULL);
2815 
2816  if( nlpiproblem->lastsolprimals == NULL )
2817  {
2818  assert(nlpiproblem->lastsoldualcons == NULL);
2819  assert(nlpiproblem->lastsoldualvarlb == NULL);
2820  assert(nlpiproblem->lastsoldualvarub == NULL);
2821  BMSallocMemoryArray(&nlpiproblem->lastsolprimals, n);
2822  BMSallocMemoryArray(&nlpiproblem->lastsoldualcons, m);
2823  BMSallocMemoryArray(&nlpiproblem->lastsoldualvarlb, n);
2824  BMSallocMemoryArray(&nlpiproblem->lastsoldualvarub, n);
2825 
2826  if( nlpiproblem->lastsolprimals == NULL || nlpiproblem->lastsoldualcons == NULL ||
2827  nlpiproblem->lastsoldualvarlb == NULL || nlpiproblem->lastsoldualvarub == NULL )
2828  {
2829  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2830  nlpiproblem->lasttermstat = SCIP_NLPTERMSTAT_MEMERR;
2831  return;
2832  }
2833  }
2834 
2835  BMScopyMemoryArray(nlpiproblem->lastsolprimals, x, n);
2836  BMScopyMemoryArray(nlpiproblem->lastsoldualcons, lambda, m);
2837  BMScopyMemoryArray(nlpiproblem->lastsoldualvarlb, z_L, n);
2838  BMScopyMemoryArray(nlpiproblem->lastsoldualvarub, z_U, n);
2839 
2840  if( check_feasibility && cq != NULL )
2841  {
2842  Number constrviol;
2843  Number constrvioltol;
2844 
2845  constrviol = cq->curr_constraint_violation();
2846 
2847  nlpiproblem->ipopt->Options()->GetNumericValue("acceptable_constr_viol_tol", constrvioltol, "");
2848  if( constrviol <= constrvioltol )
2849  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_FEASIBLE;
2850  else
2851  nlpiproblem->lastsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
2852  }
2853  }
2854 }
2855 
2856 /** Calls Lapacks Dsyev routine to compute eigenvalues and eigenvectors of a dense matrix.
2857  *
2858  * It's here, because we use Ipopt's interface to Lapack.
2859  */
2861  SCIP_Bool computeeigenvectors,/**< should also eigenvectors should be computed ? */
2862  int N, /**< dimension */
2863  SCIP_Real* a, /**< matrix data on input (size N*N); eigenvectors on output if computeeigenvectors == TRUE */
2864  SCIP_Real* w /**< buffer to store eigenvalues (size N) */
2865  )
2866 {
2867  int info;
2868 
2869  IpLapackDsyev(computeeigenvectors, N, a, N, w, info);
2870 
2871  if( info != 0 )
2872  {
2873  SCIPerrorMessage("There was an error when calling DSYEV. INFO = %d\n", info);
2874  return SCIP_ERROR;
2875  }
2876 
2877  return SCIP_OKAY;
2878 }
2879 
2880 /** solves a linear problem of the form Ax = b for a regular matrix 3*3 A */
2881 static
2883  SCIP_Real* A, /**< matrix data on input (size 3*3); filled column-wise */
2884  SCIP_Real* b, /**< right hand side vector (size 3) */
2885  SCIP_Real* x, /**< buffer to store solution (size 3) */
2886  SCIP_Bool* success /**< pointer to store if the solving routine was successful */
2887  )
2888 {
2889  SCIP_Real Acopy[9];
2890  SCIP_Real bcopy[3];
2891  int pivotcopy[3];
2892  const int N = 3;
2893  int info;
2894 
2895  assert(A != NULL);
2896  assert(b != NULL);
2897  assert(x != NULL);
2898  assert(success != NULL);
2899 
2900  BMScopyMemoryArray(Acopy, A, N*N);
2901  BMScopyMemoryArray(bcopy, b, N);
2902 
2903  /* compute the LU factorization */
2904  IpLapackDgetrf(N, Acopy, pivotcopy, N, info);
2905 
2906  if( info != 0 )
2907  {
2908  SCIPerrorMessage("There was an error when calling Dgetrf. INFO = %d\n", info);
2909  *success = FALSE;
2910  }
2911  else
2912  {
2913  *success = TRUE;
2914 
2915  /* solve linear problem */
2916  IpLapackDgetrs(N, 1, Acopy, N, pivotcopy, bcopy, N);
2917 
2918  /* copy the solution */
2919  BMScopyMemoryArray(x, bcopy, N);
2920  }
2921 
2922  return SCIP_OKAY;
2923 }
2924 
2925 /** solves a linear problem of the form Ax = b for a regular matrix A
2926  *
2927  * Calls Lapacks IpLapackDgetrf routine to calculate a LU factorization and uses this factorization to solve
2928  * the linear problem Ax = b.
2929  * It's here, because Ipopt is linked against Lapack.
2930  */
2932  int N, /**< dimension */
2933  SCIP_Real* A, /**< matrix data on input (size N*N); filled column-wise */
2934  SCIP_Real* b, /**< right hand side vector (size N) */
2935  SCIP_Real* x, /**< buffer to store solution (size N) */
2936  SCIP_Bool* success /**< pointer to store if the solving routine was successful */
2937  )
2938 {
2939  SCIP_Real* Acopy;
2940  SCIP_Real* bcopy;
2941  int* pivotcopy;
2942  int info;
2943 
2944  assert(N > 0);
2945  assert(A != NULL);
2946  assert(b != NULL);
2947  assert(x != NULL);
2948  assert(success != NULL);
2949 
2950  /* call SCIPsolveLinearProb3() for performance reasons */
2951  if( N == 3 )
2952  {
2953  SCIP_CALL( SCIPsolveLinearProb3(A, b, x, success) );
2954  return SCIP_OKAY;
2955  }
2956 
2957  Acopy = NULL;
2958  bcopy = NULL;
2959  pivotcopy = NULL;
2960 
2961  SCIP_ALLOC( BMSduplicateMemoryArray(&Acopy, A, N*N) );
2962  SCIP_ALLOC( BMSduplicateMemoryArray(&bcopy, b, N) );
2963  SCIP_ALLOC( BMSallocMemoryArray(&pivotcopy, N) );
2964 
2965  /* compute the LU factorization */
2966  IpLapackDgetrf(N, Acopy, pivotcopy, N, info);
2967 
2968  if( info != 0 )
2969  {
2970  SCIPerrorMessage("There was an error when calling Dgetrf. INFO = %d\n", info);
2971  *success = FALSE;
2972  }
2973  else
2974  {
2975  *success = TRUE;
2976 
2977  /* solve linear problem */
2978  IpLapackDgetrs(N, 1, Acopy, N, pivotcopy, bcopy, N);
2979 
2980  /* copy the solution */
2981  BMScopyMemoryArray(x, bcopy, N);
2982  }
2983 
2984  BMSfreeMemoryArray(&pivotcopy);
2985  BMSfreeMemoryArray(&bcopy);
2986  BMSfreeMemoryArray(&Acopy);
2987 
2988  return SCIP_OKAY;
2989 }
SCIP_Bool SCIPisIpoptAvailableIpopt(void)
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValues(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *convals)
Definition: nlpioracle.c:2435
enum SCIP_NlpTermStat SCIP_NLPTERMSTAT
Definition: type_nlpi.h:84
void * SCIPgetNlpiOracleIpopt(SCIP_NLPIPROBLEM *nlpiproblem)
SCIP_Real * lastsoldualvarlb
Definition: nlpi_ipopt.cpp:150
const SCIP_Real * SCIPnlpiOracleGetVarUbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2225
#define BMSfreeMemoryArrayNull(ptr)
Definition: memory.h:130
methods to interpret (evaluate) an expression tree "fast"
SCIP_RETCODE SCIPnlpiOracleGetJacobianSparsity(SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2510
SCIP_EXPRINTCAPABILITY SCIPnlpiOracleGetEvalCapability(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2375
SmartPtr< ScipNLP > nlp
Definition: nlpi_ipopt.cpp:138
void SCIPsetModifiedDefaultSettingsIpopt(SCIP_NLPI *nlpi, const char *optionsstring)
SCIP_RETCODE SCIPnlpiOracleGetHessianLagSparsity(SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2790
SCIP_NLPSOLSTAT lastsolstat
Definition: nlpi_ipopt.cpp:146
static SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemIpopt)
Definition: nlpi_ipopt.cpp:572
int SCIPnlpiOracleGetNVars(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2195
#define infinity
Definition: gastrans.c:71
static SCIP_DECL_NLPIGETWARMSTARTMEMO(nlpiGetWarmstartMemoIpopt)
static SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemIpopt)
Definition: nlpi_ipopt.cpp:670
#define NLPI_NAME
Definition: nlpi_ipopt.cpp:66
struct SCIP_NlpiProblem SCIP_NLPIPROBLEM
Definition: type_nlpi.h:39
internal methods for NLPI solver interfaces
SCIP_RETCODE SCIPnlpiOracleSetObjective(SCIP_NLPIORACLE *oracle, const SCIP_Real constant, int nlin, const int *lininds, const SCIP_Real *linvals, int nquadelems, const SCIP_QUADELEM *quadelems, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree)
Definition: nlpioracle.c:1607
methods to store an NLP and request function, gradient, and hessian values
SCIP_Real lasttime
Definition: nlpi_ipopt.cpp:154
SCIP_NLPIORACLE * oracle
#define FALSE
Definition: def.h:64
SCIP_Real SCIPnlpiOracleGetConstraintLhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2276
static SCIP_DECL_NLPIFREE(nlpiFreeIpopt)
Definition: nlpi_ipopt.cpp:535
static SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionIpopt)
#define NLPI_DESC
Definition: nlpi_ipopt.cpp:67
#define FEASTOLFACTOR
Definition: nlpi_ipopt.cpp:78
#define TRUE
Definition: def.h:63
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:53
#define DEFAULT_PRINTLEVEL
Definition: nlpi_ipopt.cpp:73
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveValue(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *objval)
Definition: nlpioracle.c:2397
static SCIP_DECL_NLPIGETINTPAR(nlpiGetIntParIpopt)
#define BMSallocMemoryArray(ptr, num)
Definition: memory.h:105
static SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessIpopt)
void SCIPmessageVPrintInfo(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr, va_list ap)
Definition: message.c:598
#define SCIP_CALL_ABORT_QUIET(x)
Definition: def.h:324
static SCIP_DECL_NLPISETINTPAR(nlpiSetIntParIpopt)
#define SCIPdebugMessage
Definition: pub_message.h:77
SCIP_RETCODE SCIPnlpiSetMessageHdlr(SCIP_NLPI *nlpi, SCIP_MESSAGEHDLR *messagehdlr)
Definition: nlpi.c:720
unsigned int SCIP_EXPRINTCAPABILITY
static SCIP_DECL_NLPIGETPROBLEMPOINTER(nlpiGetProblemPointerIpopt)
Definition: nlpi_ipopt.cpp:702
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveGradient(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *objval, SCIP_Real *objgrad)
Definition: nlpioracle.c:2461
static SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveIpopt)
Definition: nlpi_ipopt.cpp:806
SCIP_RETCODE SCIPnlpiOracleDelConsSet(SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1824
static void setOpttol(SCIP_NLPIPROBLEM *nlpiproblem, SCIP_Real opttol)
Definition: nlpi_ipopt.cpp:485
SCIP_RETCODE SCIPnlpiOracleChgObjConstant(SCIP_NLPIORACLE *oracle, SCIP_Real objconstant)
Definition: nlpioracle.c:2179
static SCIP_DECL_NLPISETSTRINGPAR(nlpiSetStringParIpopt)
int SCIPnlpiOracleGetNConstraints(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2205
std::string optfile
Definition: nlpi_ipopt.cpp:139
SCIP_Real * lastsoldualvarub
Definition: nlpi_ipopt.cpp:151
#define DEFAULT_RANDSEED
Definition: nlpi_ipopt.cpp:80
SCIP_Real lastsolinfeas
Definition: nlpi_ipopt.cpp:152
#define DEFAULT_MAXITER
Definition: nlpi_ipopt.cpp:75
void SCIPmessageVPrintError(const char *formatstr, va_list ap)
Definition: message.c:794
static SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsIpopt)
int SCIPnlpiOracleGetConstraintDegree(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2317
const char * SCIPgetSolverNameIpopt(void)
#define BMSfreeMemoryArray(ptr)
Definition: memory.h:129
SCIP_RETCODE SCIPnlpiOracleChgExprtree(SCIP_NLPIORACLE *oracle, int considx, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree)
Definition: nlpioracle.c:2097
static SCIP_DECL_NLPICHGOBJCONSTANT(nlpiChgObjConstantIpopt)
SCIP_RETCODE SCIPnlpiOracleCreate(BMS_BLKMEM *blkmem, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1325
SCIP_Bool firstrun
Definition: nlpi_ipopt.cpp:143
void SCIPmessagePrintError(const char *formatstr,...)
Definition: message.c:781
#define SCIPerrorMessage
Definition: pub_message.h:45
#define NLPI_PRIORITY
Definition: nlpi_ipopt.cpp:68
#define SCIPdebugPrintf
Definition: pub_message.h:80
const SCIP_Real * SCIPnlpiOracleGetVarLbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2215
SCIP_RETCODE SCIPnlpiOracleChgVarBounds(SCIP_NLPIORACLE *oracle, int nvars, const int *indices, const SCIP_Real *lbs, const SCIP_Real *ubs)
Definition: nlpioracle.c:1644
struct SCIP_NlpiData SCIP_NLPIDATA
Definition: type_nlpi.h:38
enum SCIP_NlpSolStat SCIP_NLPSOLSTAT
Definition: type_nlpi.h:69
static SCIP_DECL_NLPICHGQUADCOEFS(nlpiChgQuadraticCoefsIpopt)
Definition: nlpi_ipopt.cpp:978
static SCIP_DECL_NLPICHGNONLINCOEF(nlpiChgNonlinCoefIpopt)
SCIP_RETCODE SCIPnlpiOracleAddVars(SCIP_NLPIORACLE *oracle, int nvars, const SCIP_Real *lbs, const SCIP_Real *ubs, const char **varnames)
Definition: nlpioracle.c:1444
static SCIP_DECL_NLPIGETREALPAR(nlpiGetRealParIpopt)
void SCIPmessagePrintWarning(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr,...)
Definition: message.c:417
internal miscellaneous methods
void SCIPrandomFree(SCIP_RANDNUMGEN **randnumgen, BMS_BLKMEM *blkmem)
Definition: misc.c:9350
static void invalidateSolution(SCIP_NLPIPROBLEM *problem)
Definition: nlpi_ipopt.cpp:430
#define SCIP_CALL(x)
Definition: def.h:350
SCIP_Real * lastsolprimals
Definition: nlpi_ipopt.cpp:148
#define SCIP_DEFAULT_FEASTOL
Definition: def.h:157
static SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsIpopt)
Definition: nlpi_ipopt.cpp:840
void SCIPmessagePrintInfo(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr,...)
Definition: message.c:584
SCIP_EXPRINTCAPABILITY SCIPexprintGetCapability(void)
static SCIP_DECL_NLPICOPY(nlpiCopyIpopt)
Definition: nlpi_ipopt.cpp:504
SCIP_RETCODE LapackDsyev(SCIP_Bool computeeigenvectors, int N, SCIP_Real *a, SCIP_Real *w)
SCIP_Real SCIPnlpiOracleGetConstraintRhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2289
SCIP_Real SCIPnlpiOracleGetInfinity(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1397
#define BMSduplicateMemoryArray(ptr, source, num)
Definition: memory.h:125
SCIP_RETCODE SCIPnlpiOracleFree(SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1352
static const int convcheck_nchecks
Definition: nlpi_ipopt.cpp:107
methods for catching the user CTRL-C interrupt
Ipopt NLP interface.
SmartPtr< IpoptApplication > ipopt
Definition: nlpi_ipopt.cpp:137
SCIP_RETCODE SCIPnlpiCreate(SCIP_NLPI **nlpi, const char *name, const char *description, int priority, SCIP_DECL_NLPICOPY((*nlpicopy)), SCIP_DECL_NLPIFREE((*nlpifree)), SCIP_DECL_NLPIGETSOLVERPOINTER((*nlpigetsolverpointer)), SCIP_DECL_NLPICREATEPROBLEM((*nlpicreateproblem)), SCIP_DECL_NLPIFREEPROBLEM((*nlpifreeproblem)), SCIP_DECL_NLPIGETPROBLEMPOINTER((*nlpigetproblempointer)), SCIP_DECL_NLPIADDVARS((*nlpiaddvars)), SCIP_DECL_NLPIADDCONSTRAINTS((*nlpiaddconstraints)), SCIP_DECL_NLPISETOBJECTIVE((*nlpisetobjective)), SCIP_DECL_NLPICHGVARBOUNDS((*nlpichgvarbounds)), SCIP_DECL_NLPICHGCONSSIDES((*nlpichgconssides)), SCIP_DECL_NLPIDELVARSET((*nlpidelvarset)), SCIP_DECL_NLPIDELCONSSET((*nlpidelconsset)), SCIP_DECL_NLPICHGLINEARCOEFS((*nlpichglinearcoefs)), SCIP_DECL_NLPICHGQUADCOEFS((*nlpichgquadcoefs)), SCIP_DECL_NLPICHGEXPRTREE((*nlpichgexprtree)), SCIP_DECL_NLPICHGNONLINCOEF((*nlpichgnonlincoef)), SCIP_DECL_NLPICHGOBJCONSTANT((*nlpichgobjconstant)), SCIP_DECL_NLPISETINITIALGUESS((*nlpisetinitialguess)), SCIP_DECL_NLPISOLVE((*nlpisolve)), SCIP_DECL_NLPIGETSOLSTAT((*nlpigetsolstat)), SCIP_DECL_NLPIGETTERMSTAT((*nlpigettermstat)), SCIP_DECL_NLPIGETSOLUTION((*nlpigetsolution)), SCIP_DECL_NLPIGETSTATISTICS((*nlpigetstatistics)), SCIP_DECL_NLPIGETWARMSTARTSIZE((*nlpigetwarmstartsize)), SCIP_DECL_NLPIGETWARMSTARTMEMO((*nlpigetwarmstartmemo)), SCIP_DECL_NLPISETWARMSTARTMEMO((*nlpisetwarmstartmemo)), SCIP_DECL_NLPIGETINTPAR((*nlpigetintpar)), SCIP_DECL_NLPISETINTPAR((*nlpisetintpar)), SCIP_DECL_NLPIGETREALPAR((*nlpigetrealpar)), SCIP_DECL_NLPISETREALPAR((*nlpisetrealpar)), SCIP_DECL_NLPIGETSTRINGPAR((*nlpigetstringpar)), SCIP_DECL_NLPISETSTRINGPAR((*nlpisetstringpar)), SCIP_DECL_NLPISETMESSAGEHDLR((*nlpisetmessagehdlr)), SCIP_NLPIDATA *nlpidata)
Definition: nlpi.c:40
public data structures and miscellaneous methods
#define SCIP_EXPRINTCAPABILITY_HESSIAN
static SCIP_DECL_NLPISETWARMSTARTMEMO(nlpiSetWarmstartMemoIpopt)
#define SCIP_Bool
Definition: def.h:61
SCIP_Bool SCIPinterrupted(void)
Definition: interrupt.c:147
#define MAXPERTURB
Definition: nlpi_ipopt.cpp:77
SCIP_RETCODE SCIPnlpiOracleSetProblemName(SCIP_NLPIORACLE *oracle, const char *name)
Definition: nlpioracle.c:1409
static void setFeastol(SCIP_NLPIPROBLEM *nlpiproblem, SCIP_Real feastol)
Definition: nlpi_ipopt.cpp:454
SCIP_RETCODE SCIPnlpiOracleAddConstraints(SCIP_NLPIORACLE *oracle, int nconss, const SCIP_Real *lhss, const SCIP_Real *rhss, const int *nlininds, int *const *lininds, SCIP_Real *const *linvals, const int *nquadelems, SCIP_QUADELEM *const *quadelems, int *const *exprvaridxs, SCIP_EXPRTREE *const *exprtrees, const char **consnames)
Definition: nlpioracle.c:1529
static SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsIpopt)
Definition: nlpi_ipopt.cpp:954
void SCIPnlpStatisticsSetNIterations(SCIP_NLPSTATISTICS *statistics, int niterations)
Definition: nlpi.c:836
static const int convcheck_startiter
Definition: nlpi_ipopt.cpp:108
static SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatIpopt)
#define MAX(x, y)
Definition: tclique_def.h:75
void SCIPnlpStatisticsSetTotalTime(SCIP_NLPSTATISTICS *statistics, SCIP_Real totaltime)
Definition: nlpi.c:846
static SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatIpopt)
SCIP_RETCODE SCIPnlpiOracleChgLinearCoefs(SCIP_NLPIORACLE *oracle, int considx, int nentries, const int *varidxs, const SCIP_Real *newcoefs)
Definition: nlpioracle.c:1902
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:116
static SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsIpopt)
Definition: nlpi_ipopt.cpp:766
SCIP_RETCODE SCIPnlpiOracleChgQuadCoefs(SCIP_NLPIORACLE *oracle, int considx, int nquadelems, const SCIP_QUADELEM *quadelems)
Definition: nlpioracle.c:1999
SCIP_RETCODE SCIPcreateNlpSolverIpopt(BMS_BLKMEM *blkmem, SCIP_NLPI **nlpi)
SCIP_RETCODE SCIPnlpiOracleChgConsSides(SCIP_NLPIORACLE *oracle, int nconss, const int *indices, const SCIP_Real *lhss, const SCIP_Real *rhss)
Definition: nlpioracle.c:1680
SCIP_Real SCIPrandomGetReal(SCIP_RANDNUMGEN *randnumgen, SCIP_Real minrandval, SCIP_Real maxrandval)
Definition: misc.c:9388
SCIP_RETCODE SCIPnlpiOracleDelVarSet(SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1714
static SCIP_DECL_NLPIGETSOLVERPOINTER(nlpiGetSolverPointerIpopt)
Definition: nlpi_ipopt.cpp:557
SCIP_RETCODE SCIPsolveLinearProb(int N, SCIP_Real *A, SCIP_Real *b, SCIP_Real *x, SCIP_Bool *success)
char * SCIPnlpiOracleGetConstraintName(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2302
#define SCIP_DEFAULT_INFINITY
Definition: def.h:154
#define SCIP_EXPRINTCAPABILITY_FUNCVALUE
#define SCIP_DEFAULT_DUALFEASTOL
Definition: def.h:160
void * SCIPgetIpoptApplicationPointerIpopt(SCIP_NLPIPROBLEM *nlpiproblem)
SCIP_RETCODE SCIPnlpiOracleChgExprParam(SCIP_NLPIORACLE *oracle, int considx, int paramidx, SCIP_Real paramval)
Definition: nlpioracle.c:2154
int SCIPnlpiOracleGetVarDegree(SCIP_NLPIORACLE *oracle, int varidx)
Definition: nlpioracle.c:2247
SCIP_NLPIDATA * SCIPnlpiGetData(SCIP_NLPI *nlpi)
Definition: nlpi.c:733
static SCIP_DECL_NLPIGETWARMSTARTSIZE(nlpiGetWarmstartSizeIpopt)
SCIP_RETCODE SCIPnlpiOracleEvalHessianLag(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real objfactor, const SCIP_Real *lambda, SCIP_Real *hessian)
Definition: nlpioracle.c:2889
const char * SCIPgetSolverDescIpopt(void)
static SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetIpopt)
Definition: nlpi_ipopt.cpp:899
#define SCIP_Real
Definition: def.h:149
SCIP_NLPTERMSTAT lasttermstat
Definition: nlpi_ipopt.cpp:147
static SCIP_DECL_NLPISETMESSAGEHDLR(nlpiSetMessageHdlrIpopt)
SCIP_Real * lastsoldualcons
Definition: nlpi_ipopt.cpp:149
SCIP_RETCODE SCIPrandomCreate(SCIP_RANDNUMGEN **randnumgen, BMS_BLKMEM *blkmem, unsigned int initialseed)
Definition: misc.c:9334
#define SCIP_INVALID
Definition: def.h:169
static SCIP_DECL_NLPISOLVE(nlpiSolveIpopt)
SCIP_RETCODE SCIPnlpiOracleSetInfinity(SCIP_NLPIORACLE *oracle, SCIP_Real infinity)
Definition: nlpioracle.c:1381
static SCIP_DECL_NLPIGETSTRINGPAR(nlpiGetStringParIpopt)
static SCIP_DECL_NLPIDELCONSSET(nlpiDelConstraintSetIpopt)
Definition: nlpi_ipopt.cpp:926
static SCIP_DECL_NLPISETREALPAR(nlpiSetRealParIpopt)
SCIP_RETCODE SCIPnlpiOracleEvalJacobian(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *convals, SCIP_Real *jacobi)
Definition: nlpioracle.c:2651
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:112
static const SCIP_Real convcheck_minred[convcheck_nchecks]
Definition: nlpi_ipopt.cpp:110
struct BMS_BlkMem BMS_BLKMEM
Definition: memory.h:419
char ** SCIPnlpiOracleGetVarNames(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2235
static SCIP_DECL_NLPICHGEXPRTREE(nlpiChgExprtreeIpopt)
Definition: nlpi_ipopt.cpp:999
#define SCIP_ALLOC(x)
Definition: def.h:361
#define SCIPABORT()
Definition: def.h:322
static SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesIpopt)
Definition: nlpi_ipopt.cpp:874
static SCIP_RETCODE SCIPsolveLinearProb3(SCIP_Real *A, SCIP_Real *b, SCIP_Real *x, SCIP_Bool *success)
static SCIP_DECL_NLPIADDVARS(nlpiAddVarsIpopt)
Definition: nlpi_ipopt.cpp:721
#define SCIP_EXPRINTCAPABILITY_GRADIENT
int SCIPnlpiOracleGetMaxDegree(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2349
SCIP_RETCODE SCIPnlpiSetRealPar(SCIP_NLPI *nlpi, SCIP_NLPIPROBLEM *problem, SCIP_NLPPARAM type, SCIP_Real dval)
Definition: nlpi.c:671
static const int convcheck_maxiter[convcheck_nchecks]
Definition: nlpi_ipopt.cpp:109