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-2025 Zuse Institute Berlin (ZIB) */
7/* */
8/* Licensed under the Apache License, Version 2.0 (the "License"); */
9/* you may not use this file except in compliance with the License. */
10/* You may obtain a copy of the License at */
11/* */
12/* http://www.apache.org/licenses/LICENSE-2.0 */
13/* */
14/* Unless required by applicable law or agreed to in writing, software */
15/* distributed under the License is distributed on an "AS IS" BASIS, */
16/* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. */
17/* See the License for the specific language governing permissions and */
18/* limitations under the License. */
19/* */
20/* You should have received a copy of the Apache-2.0 license */
21/* along with SCIP; see the file LICENSE. If not visit scipopt.org. */
22/* */
23/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
24
25/**@file nlpi_ipopt.cpp
26 * @ingroup DEFPLUGINS_NLPI
27 * @brief Ipopt NLP interface
28 * @author Stefan Vigerske
29 * @author Benjamin Müller
30 *
31 * @todo if too few degrees of freedom, solve a slack-minimization problem instead?
32 * @todo automatically switch to Hessian approximation if Hessian is dense or slow? (only do so if opttol/solvertol is large?)
33 *
34 * This file can only be compiled if Ipopt is available.
35 * Otherwise, to resolve public functions, use nlpi_ipopt_dummy.c.
36 * Since the dummy code is C instead of C++, it has been moved into a separate file.
37 */
38
39/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
40
41#include "scip/nlpi_ipopt.h"
42
43#include "scip/nlpioracle.h"
44#include "scip/exprinterpret.h"
45#include "scip/scip_nlpi.h"
47#include "scip/scip_mem.h"
48#include "scip/scip_message.h"
49#include "scip/scip_general.h"
50#include "scip/scip_numerics.h"
51#include "scip/scip_param.h"
52#include "scip/scip_solve.h"
53#include "scip/scip_copy.h"
54#include "scip/pub_misc.h"
55#include "scip/pub_paramset.h"
56#include "scip/pub_message.h"
57
58#include <new> /* for std::bad_alloc */
59#include <sstream>
60#include <cstring>
61
62/* turn off some lint warnings for file */
63/*lint --e{1540,750,3701}*/
64
65#include "IpoptConfig.h"
66
67#if defined(__GNUC__) && IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
68#pragma GCC diagnostic ignored "-Wshadow"
69#endif
70#include "IpIpoptApplication.hpp"
71#include "IpIpoptCalculatedQuantities.hpp"
72#include "IpSolveStatistics.hpp"
73#include "IpJournalist.hpp"
74#include "IpIpoptData.hpp"
75#include "IpTNLPAdapter.hpp"
76#include "IpOrigIpoptNLP.hpp"
77#include "IpLapack.hpp"
78#if defined(__GNUC__) && IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
79#pragma GCC diagnostic warning "-Wshadow"
80#endif
81
82#if IPOPT_VERSION_MAJOR < 3 || (IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 12)
83#error "The Ipopt interface requires at least 3.12.0"
84#endif
85
86/* MUMPS that can be used by Ipopt is not threadsafe
87 * If we want SCIP to be threadsafe (SCIP_THREADSAFE), have std::mutex (C++11 or higher), and use Ipopt before 3.14,
88 * then we protect the call to Ipopt by a mutex if MUMPS is used as linear solver.
89 * Thus, we allow only one Ipopt run at a time.
90 * Ipopt 3.14 has this build-in to its MUMPS interface, so we won't have to take care of this.
91 */
92#if defined(SCIP_THREADSAFE) && __cplusplus >= 201103L && IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
93#define PROTECT_SOLVE_BY_MUTEX
94#include <mutex>
95static std::mutex solve_mutex; /*lint !e1756*/
96#endif
97
98using namespace Ipopt;
99
100#define NLPI_NAME "ipopt" /**< short concise name of solver */
101#define NLPI_DESC "Ipopt interface" /**< description of solver */
102#define NLPI_PRIORITY 1000 /**< priority */
103
104#define MAXPERTURB 0.01 /**< maximal perturbation of bounds in starting point heuristic */
105#define FEASTOLFACTOR 0.9 /**< factor for user-given feasibility tolerance to get feasibility tolerance that is actually passed to Ipopt */
106
107#define DEFAULT_RANDSEED 71 /**< initial random seed */
108
109// enable this to collect statistics on number of iterations and problem characteristics in csv-form in log
110// note that this overwrites given iterlimit
111// see https://git.zib.de/integer/scip/-/snippets/1213 for some script that evaluates the collected data
112// #define COLLECT_SOLVESTATS
113
114/* Convergence check (see ScipNLP::intermediate_callback)
115 *
116 * If the fastfail option is set to aggressive, then we stop Ipopt if the reduction in
117 * primal infeasibility is not sufficient for a consecutive number of iterations.
118 * With the parameters as given below, we require Ipopt to
119 * - not increase the primal infeasibility after 5 iterations
120 * - reduce the primal infeasibility by at least 50% within 10 iterations
121 * - reduce the primal infeasibility by at least 90% within 30 iterations
122 * The targets are updated once they are reached and the limit on allowed iterations to reach the new target is reset.
123 *
124 * In certain situations, it is allowed to exceed an iteration limit:
125 * - If we are in the first 10 (convcheck_startiter) iterations.
126 * - If we are within 10 (convcheck_startiter) iterations after the restoration phase ended.
127 * The reason for this is that during feasibility restoration phase Ipopt aims completely on
128 * reducing constraint violation, completely forgetting the objective function.
129 * When returning from feasibility restoration and considering the original objective again,
130 * it is unlikely that Ipopt will continue to decrease primal infeasibility, since it may now target on
131 * more on optimality again. Thus, we do not check convergence for a number of iterations.
132 * - If the target on dual infeasibility reduction has been achieved, we are below twice the iteration limit, and
133 * we are not in restoration mode.
134 * The reason for this is that if Ipopt makes good progress towards optimality,
135 * we want to allow some more iterations where primal infeasibility is not reduced.
136 * However, in restoration mode, dual infeasibility does not correspond to the original problem and
137 * the complete aim is to restore primal infeasibility.
138 */
139static const int convcheck_nchecks = 3; /**< number of convergence checks */
140static const int convcheck_startiter = 10; /**< iteration where to start convergence checking */
141static const int convcheck_maxiter[convcheck_nchecks] = { 5, 15, 30 }; /**< maximal number of iterations to achieve each convergence check */
142static const SCIP_Real convcheck_minred[convcheck_nchecks] = { 1.0, 0.5, 0.1 }; /**< minimal required infeasibility reduction in each convergence check */
143
144/// integer parameters of Ipopt to make available via SCIP parameters
145static const char* ipopt_int_params[] =
146 { "print_level" }; // print_level must be first
147
148/// string parameters of Ipopt to make available via SCIP parameters
149static const char* ipopt_string_params[] =
150 { "linear_solver",
151 "hsllib",
152 "pardisolib",
153 "linear_system_scaling",
154 "nlp_scaling_method",
155 "mu_strategy",
156 "hessian_approximation"
157 };
158
159class ScipNLP;
160
161struct SCIP_NlpiData
162{
163public:
164 char* optfile; /**< Ipopt options file to read */
165 int print_level; /**< print_level set via nlpi/ipopt/print_level option */
166 SCIP_Real warm_start_push; /**< value to use for Ipopt's warm_start_bound_push/frac options */
167
168 /** constructor */
169 explicit SCIP_NlpiData()
170 : optfile(NULL), print_level(-1), warm_start_push(1e-9)
171 { }
172};
173
174struct SCIP_NlpiProblem
175{
176public:
177 SCIP_NLPIORACLE* oracle; /**< Oracle-helper to store and evaluate NLP */
178 SCIP_RANDNUMGEN* randnumgen; /**< random number generator */
179
180 SmartPtr<IpoptApplication> ipopt; /**< Ipopt application */
181 SmartPtr<ScipNLP> nlp; /**< NLP in Ipopt form */
182
183 bool firstrun; /**< whether the next NLP solve will be the first one */
184 bool samestructure;/**< whether the NLP solved next will still have the same (Ipopt-internal) structure (same number of variables, constraints, bounds, and nonzero pattern) */
185
186 SCIP_NLPSOLSTAT solstat; /**< status of current solution (if any) */
187 SCIP_NLPTERMSTAT termstat; /**< termination status of last solve (if any) */
188 bool solprimalvalid;/**< whether primal solution values are available (solprimals has meaningful values) */
189 bool solprimalgiven;/**< whether primal solution values were set by caller */
190 bool soldualvalid; /**< whether dual solution values are available (soldual* have meaningful values) */
191 bool soldualgiven; /**< whether dual solution values were set by caller */
192 SCIP_Real* solprimals; /**< primal solution values, if available */
193 SCIP_Real* soldualcons; /**< dual solution values of constraints, if available */
194 SCIP_Real* soldualvarlb; /**< dual solution values of variable lower bounds, if available */
195 SCIP_Real* soldualvarub; /**< dual solution values of variable upper bounds, if available */
196 SCIP_Real solobjval; /**< objective function value in solution from last run */
197 SCIP_Real solconsviol; /**< constraint violation of primal solution, if available */
198 SCIP_Real solboundviol; /**< variable bound violation of primal solution, if available */
199 int lastniter; /**< number of iterations in last run */
200 SCIP_Real lasttime; /**< time spend in last run */
201
202 int ncalls; /**< overall number of solver calls */
203 int nsuccess; /**< number of successes (optimal or feasible solution found or proven unbounded) */
204 int nlocinfeas; /**< number of calls resulting in local infeasibility */
205 int nother; /**< number of other calls */
206 int nlimit; /**< number of calls where the solver terminated due to a time or iteration limit */
207
208 /** constructor */
211 firstrun(true), samestructure(true),
213 solprimalvalid(false), solprimalgiven(false), soldualvalid(false), soldualgiven(false),
216 lastniter(-1), lasttime(-1.0), ncalls(0), nsuccess(0), nlocinfeas(0), nother(0), nlimit(0)
217 { }
218};
219
220/** TNLP implementation for SCIPs NLP */
221class ScipNLP : public TNLP
222{
223private:
224 SCIP_NLPIPROBLEM* nlpiproblem; /**< NLPI problem data */
225 SCIP* scip; /**< SCIP data structure */
226 SCIP_NLPPARAM param; /**< NLP solve parameters */
227
228 SCIP_Real conv_prtarget[convcheck_nchecks]; /**< target primal infeasibility for each convergence check */
229 SCIP_Real conv_dutarget[convcheck_nchecks]; /**< target dual infeasibility for each convergence check */
230 int conv_iterlim[convcheck_nchecks]; /**< iteration number where target primal infeasibility should to be achieved */
231 int conv_lastrestoiter; /**< last iteration number in restoration mode, or -1 if none */
232
233 unsigned int current_x; /**< unique number that identifies current iterate (x): incremented when Ipopt calls with new_x=true */
234 unsigned int last_f_eval_x; /**< the number of the iterate for which the objective was last evaluated (eval_f) */
235 unsigned int last_g_eval_x; /**< the number of the iterate for which the constraints were last evaluated (eval_g) */
236
237public:
238 bool approxhessian; /**< do we tell Ipopt to approximate the hessian? (may also be false if user set to approx. hessian via option file) */
239
240 // cppcheck-suppress uninitMemberVar
241 /** constructor */
242 ScipNLP(
243 SCIP_NLPIPROBLEM* nlpiproblem_ = NULL,/**< NLPI problem data */
244 SCIP* scip_ = NULL /**< SCIP data structure */
245 )
246 : nlpiproblem(nlpiproblem_), scip(scip_),
247 conv_lastrestoiter(-1),
248 current_x(1), last_f_eval_x(0), last_g_eval_x(0),
249 approxhessian(false)
250 {
251 assert(scip != NULL);
252 }
253
254 /** destructor */
255 ~ScipNLP()
256 { /*lint --e{1540}*/
257 }
258
259 /** initialize for new solve */
260 void initializeSolve(
261 SCIP_NLPIPROBLEM* nlpiproblem_, /**< NLPI problem */
262 const SCIP_NLPPARAM& nlpparam /**< NLP solve parameters */
263 )
264 {
265 assert(nlpiproblem_ != NULL);
266 nlpiproblem = nlpiproblem_;
267 param = nlpparam;
268
269 // since we are about to start a new solve, use this opportunity to reset the counts on x
270 current_x = 1;
271 last_f_eval_x = 0;
272 last_g_eval_x = 0;
273 }
274
275 /** Method to return some info about the nlp */
276 bool get_nlp_info(
277 Index& n, /**< place to store number of variables */
278 Index& m, /**< place to store number of constraints */
279 Index& nnz_jac_g, /**< place to store number of nonzeros in jacobian */
280 Index& nnz_h_lag, /**< place to store number of nonzeros in hessian */
281 IndexStyleEnum& index_style /**< place to store used index style (0-based or 1-based) */
282 );
283
284 /** Method to return the bounds for my problem */
285 bool get_bounds_info(
286 Index n, /**< number of variables */
287 Number* x_l, /**< buffer to store lower bounds on variables */
288 Number* x_u, /**< buffer to store upper bounds on variables */
289 Index m, /**< number of constraints */
290 Number* g_l, /**< buffer to store lower bounds on constraints */
291 Number* g_u /**< buffer to store lower bounds on constraints */
292 );
293
294 /** Method to return the starting point for the algorithm */
295 bool get_starting_point(
296 Index n, /**< number of variables */
297 bool init_x, /**< whether initial values for primal values are requested */
298 Number* x, /**< buffer to store initial primal values */
299 bool init_z, /**< whether initial values for dual values of variable bounds are requested */
300 Number* z_L, /**< buffer to store dual values for variable lower bounds */
301 Number* z_U, /**< buffer to store dual values for variable upper bounds */
302 Index m, /**< number of constraints */
303 bool init_lambda, /**< whether initial values for dual values of constraints are required */
304 Number* lambda /**< buffer to store dual values of constraints */
305 );
306
307 /** Method to return the number of nonlinear variables. */
308 Index get_number_of_nonlinear_variables();
309
310 /** Method to return the indices of the nonlinear variables */
311 bool get_list_of_nonlinear_variables(
312 Index num_nonlin_vars, /**< number of nonlinear variables */
313 Index* pos_nonlin_vars /**< array to fill with indices of nonlinear variables */
314 );
315
316 /** Method to return metadata about variables and constraints */
317 bool get_var_con_metadata(
318 Index n, /**< number of variables */
319 StringMetaDataMapType& var_string_md, /**< variable meta data of string type */
320 IntegerMetaDataMapType& var_integer_md,/**< variable meta data of integer type */
321 NumericMetaDataMapType& var_numeric_md,/**< variable meta data of numeric type */
322 Index m, /**< number of constraints */
323 StringMetaDataMapType& con_string_md, /**< constraint meta data of string type */
324 IntegerMetaDataMapType& con_integer_md,/**< constraint meta data of integer type */
325 NumericMetaDataMapType& con_numeric_md /**< constraint meta data of numeric type */
326 );
327
328 /** Method to return the objective value */
329 bool eval_f(
330 Index n, /**< number of variables */
331 const Number* x, /**< point to evaluate */
332 bool new_x, /**< whether some function evaluation method has been called for this point before */
333 Number& obj_value /**< place to store objective function value */
334 );
335
336 /** Method to return the gradient of the objective */
337 bool eval_grad_f(
338 Index n, /**< number of variables */
339 const Number* x, /**< point to evaluate */
340 bool new_x, /**< whether some function evaluation method has been called for this point before */
341 Number* grad_f /**< buffer to store objective gradient */
342 );
343
344 /** Method to return the constraint residuals */
345 bool eval_g(
346 Index n, /**< number of variables */
347 const Number* x, /**< point to evaluate */
348 bool new_x, /**< whether some function evaluation method has been called for this point before */
349 Index m, /**< number of constraints */
350 Number* g /**< buffer to store constraint function values */
351 );
352
353 /** Method to return:
354 * 1) The structure of the jacobian (if "values" is NULL)
355 * 2) The values of the jacobian (if "values" is not NULL)
356 */
357 bool eval_jac_g(
358 Index n, /**< number of variables */
359 const Number* x, /**< point to evaluate */
360 bool new_x, /**< whether some function evaluation method has been called for this point before */
361 Index m, /**< number of constraints */
362 Index nele_jac, /**< number of nonzero entries in jacobian */
363 Index* iRow, /**< buffer to store row indices of nonzero jacobian entries, or NULL if values
364 * are requested */
365 Index* jCol, /**< buffer to store column indices of nonzero jacobian entries, or NULL if values
366 * are requested */
367 Number* values /**< buffer to store values of nonzero jacobian entries, or NULL if structure is
368 * requested */
369 );
370
371 /** Method to return:
372 * 1) The structure of the hessian of the lagrangian (if "values" is NULL)
373 * 2) The values of the hessian of the lagrangian (if "values" is not NULL)
374 */
375 bool eval_h(
376 Index n, /**< number of variables */
377 const Number* x, /**< point to evaluate */
378 bool new_x, /**< whether some function evaluation method has been called for this point before */
379 Number obj_factor, /**< weight for objective function */
380 Index m, /**< number of constraints */
381 const Number* lambda, /**< weights for constraint functions */
382 bool new_lambda, /**< whether the hessian has been evaluated for these values of lambda before */
383 Index nele_hess, /**< number of nonzero entries in hessian */
384 Index* iRow, /**< buffer to store row indices of nonzero hessian entries, or NULL if values
385 * are requested */
386 Index* jCol, /**< buffer to store column indices of nonzero hessian entries, or NULL if values
387 * are requested */
388 Number* values /**< buffer to store values of nonzero hessian entries, or NULL if structure is requested */
389 );
390
391 /** Method called by the solver at each iteration.
392 *
393 * Checks whether Ctrl-C was hit.
394 */
395 bool intermediate_callback(
396 AlgorithmMode mode, /**< current mode of algorithm */
397 Index iter, /**< current iteration number */
398 Number obj_value, /**< current objective value */
399 Number inf_pr, /**< current primal infeasibility */
400 Number inf_du, /**< current dual infeasibility */
401 Number mu, /**< current barrier parameter */
402 Number d_norm, /**< current gradient norm */
403 Number regularization_size,/**< current size of regularization */
404 Number alpha_du, /**< current dual alpha */
405 Number alpha_pr, /**< current primal alpha */
406 Index ls_trials, /**< current number of linesearch trials */
407 const IpoptData* ip_data, /**< pointer to Ipopt Data */
408 IpoptCalculatedQuantities* ip_cq /**< pointer to current calculated quantities */
409 );
410
411 /** This method is called when the algorithm is complete so the TNLP can store/write the solution. */
412 void finalize_solution(
413 SolverReturn status, /**< solve and solution status */
414 Index n, /**< number of variables */
415 const Number* x, /**< primal solution values */
416 const Number* z_L, /**< dual values of variable lower bounds */
417 const Number* z_U, /**< dual values of variable upper bounds */
418 Index m, /**< number of constraints */
419 const Number* g, /**< values of constraints */
420 const Number* lambda, /**< dual values of constraints */
421 Number obj_value, /**< objective function value */
422 const IpoptData* data, /**< pointer to Ipopt Data */
423 IpoptCalculatedQuantities* cq /**< pointer to calculated quantities */
424 );
425};
426
427/** A particular Ipopt::Journal implementation that uses the SCIP message routines for output. */
428class ScipJournal : public Ipopt::Journal {
429private:
430 SCIP* scip; /**< SCIP data structure */
431
432public:
433 ScipJournal(
434 const char* name, /**< name of journal */
435 Ipopt::EJournalLevel default_level, /**< default verbosity level */
436 SCIP* scip_ /**< SCIP data structure */
437 )
438 : Ipopt::Journal(name, default_level),
439 scip(scip_)
440 { }
441
442 ~ScipJournal() { }
443
444protected:
445 /*lint -e{715}*/
446 void PrintImpl(
447 Ipopt::EJournalCategory category, /**< category of message */
448 Ipopt::EJournalLevel level, /**< verbosity level of message */
449 const char* str /**< message to print */
450 )
451 { /*lint --e{715} */
452 if( level == J_ERROR )
453 {
454 SCIPerrorMessage("%s", str);
455 }
456 else
457 {
458 SCIPinfoMessage(scip, NULL, "%s", str);
459 }
460 }
461
462 /*lint -e{715}*/
463 void PrintfImpl(
464 Ipopt::EJournalCategory category, /**< category of message */
465 Ipopt::EJournalLevel level, /**< verbosity level of message */
466 const char* pformat, /**< message printing format */
467 va_list ap /**< arguments of message */
468 )
469 { /*lint --e{715} */
470 if( level == J_ERROR )
471 {
472 SCIPmessagePrintErrorHeader(__FILE__, __LINE__);
473 SCIPmessageVPrintError(pformat, ap);
474 }
475 else
476 {
478 }
479 }
480
481 void FlushBufferImpl() { }
482};
483
484/** sets status codes to mark that last NLP solve is no longer valid (usually because the NLP changed) */
485static
487 SCIP_NLPIPROBLEM* problem /**< data structure of problem */
488 )
489{
492 problem->solobjval = SCIP_INVALID;
493 problem->solconsviol = SCIP_INVALID;
494 problem->solboundviol = SCIP_INVALID;
495 problem->lastniter = -1;
496 problem->lasttime = -1.0;
497}
498
499/** sets solution values to be invalid and calls invalidateSolved() */
500static
502 SCIP_NLPIPROBLEM* problem /**< data structure of problem */
503 )
504{
505 assert(problem != NULL);
506
507 problem->solprimalvalid = false;
508 problem->solprimalgiven = false;
509 problem->soldualvalid = false;
510 problem->soldualgiven = false;
511
512 invalidateSolved(problem);
513}
514
515/** makes sure a starting point (initial guess) is available */
516static
518 SCIP* scip, /**< SCIP data structure */
519 SCIP_NLPIPROBLEM* problem, /**< data structure of problem */
520 SCIP_Bool& warmstart /**< whether a warmstart has been requested */
521 )
522{
523 SCIP_Real lb, ub;
524 int n;
525
526 assert(problem != NULL);
527
528 // disable warmstart if no primal or dual solution values are available
529 if( warmstart && (!problem->solprimalvalid || !problem->soldualvalid ))
530 {
531 SCIPdebugMsg(scip, "Disable warmstart as no primal or dual solution available.\n");
532 warmstart = false;
533 }
534
535 // continue below with making up a random primal starting point if
536 // the user did not set a starting point and warmstart is disabled (so the last solution shouldn't be used)
537 // (if warmstart, then due to the checks above we must now have valid primal and dual solution values)
538 if( problem->solprimalgiven || warmstart )
539 {
540 // so we must have a primal solution to start from
541 // if warmstart, then we also need to have a dual solution to start from
542 // if warmstart and primal solution is given by user, then also dual solution should have been given by user
543 assert(problem->solprimalvalid);
544 assert(problem->solprimals != NULL);
545 assert(!warmstart || !problem->solprimalgiven || problem->soldualgiven);
546 assert(!warmstart || problem->soldualcons != NULL);
547 assert(!warmstart || problem->soldualvarlb != NULL);
548 assert(!warmstart || problem->soldualvarub != NULL);
549 SCIPdebugMsg(scip, "Starting solution for %sstart available from %s.\n",
550 warmstart ? "warm" : "cold",
551 problem->solprimalgiven ? "user" : "previous solve");
552 return SCIP_OKAY;
553 }
554
555 SCIPdebugMsg(scip, "Starting solution for coldstart not available. Making up something by projecting 0 onto variable bounds and adding a random perturbation.\n");
556
557 n = SCIPnlpiOracleGetNVars(problem->oracle);
558
559 if( problem->randnumgen == NULL )
560 {
562 }
563
564 if( problem->solprimals == NULL )
565 {
567 }
568
569 for( int i = 0; i < n; ++i )
570 {
571 lb = SCIPnlpiOracleGetVarLbs(problem->oracle)[i];
572 ub = SCIPnlpiOracleGetVarUbs(problem->oracle)[i];
573 if( lb > 0.0 )
574 problem->solprimals[i] = SCIPrandomGetReal(problem->randnumgen, lb, lb + MAXPERTURB*MIN(1.0, ub-lb));
575 else if( ub < 0.0 )
576 problem->solprimals[i] = SCIPrandomGetReal(problem->randnumgen, ub - MAXPERTURB*MIN(1.0, ub-lb), ub);
577 else
578 problem->solprimals[i] = SCIPrandomGetReal(problem->randnumgen, MAX(lb, -MAXPERTURB*MIN(1.0, ub-lb)), MIN(ub, MAXPERTURB*MIN(1.0, ub-lb)));
579 }
580 problem->solprimalvalid = true;
581
582 return SCIP_OKAY;
583}
584
585/// pass NLP solve parameters to Ipopt
586static
588 SCIP* scip, /**< SCIP data structure */
589 SCIP_NLPIDATA* nlpidata, /**< NLPI data */
590 SCIP_NLPIPROBLEM* nlpiproblem, /**< NLP */
591 SCIP_NLPPARAM param /**< solve parameters */
592 )
593{
594 assert(scip != NULL);
595 assert(nlpiproblem != NULL);
596
597 if( nlpidata->print_level < 0 ) // if nlpi/ipopt/print_level param has not been set
598 {
599 switch( param.verblevel )
600 {
601 case 0:
602 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("print_level", J_ERROR);
603 break;
604 case 1:
605 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("print_level", J_SUMMARY);
606 break;
607 case 2:
608 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("print_level", J_ITERSUMMARY);
609 break;
610 case 3:
611 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("print_level", J_DETAILED);
612 break;
613 default:
614 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("print_level", MIN(J_ITERSUMMARY + (param.verblevel-1), J_ALL));
615 break;
616 }
617 }
618
619#ifdef SCIP_DISABLED_CODE
620 if( param.iterlimit < 0 )
621 {
622 if( nlpidata->autoiterlim > 0 )
623 {
624 param.iterlimit = nlpidata->autoiterlim;
625 }
626 else
627 {
628 int nvars = 0; // number of variables, including slacks
629 int nnlvars = 0; // number of nonlinear variables
630 int nvarbnd = 0; // number of finite lower and upper bounds, including slacks
631 int nlincons = 0; // number of linear constraints
632 int nnlcons = 0; // number of nonlinear constraints
633 int jacnnz = 0; // number of nonzeros in Jacobian
634
635 int n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
636 int m = SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle);
637 const SCIP_Real* lbs = SCIPnlpiOracleGetVarLbs(nlpiproblem->oracle);
638 const SCIP_Real* ubs = SCIPnlpiOracleGetVarUbs(nlpiproblem->oracle);
639
640 for( int i = 0; i < n; ++i )
641 {
642 // skip fixed vars
643 if( SCIPisEQ(scip, lbs[i], ubs[i]) )
644 continue;
645
646 ++nvars;
647
648 if( SCIPnlpiOracleIsVarNonlinear(scip, nlpiproblem->oracle, i) )
649 ++nnlvars;
650
651 // every variable lower bound contributes a ln(x-lb) term in the barrier problem
652 if( !SCIPisInfinity(scip, -lbs[i]) )
653 ++nvarbnd;
654
655 // every variable upper bound contributes a ln(ub-x) term in the barrier problem
656 if( !SCIPisInfinity(scip, ubs[i]) )
657 ++nvarbnd;
658 }
659
660 for( int i = 0; i < m; ++i )
661 {
662 if( SCIPnlpiOracleIsConstraintNonlinear(nlpiproblem->oracle, i) )
663 ++nnlcons;
664 else
665 ++nlincons;
666
667 SCIP_Real lhs = SCIPnlpiOracleGetConstraintLhs(nlpiproblem->oracle, i);
668 SCIP_Real rhs = SCIPnlpiOracleGetConstraintRhs(nlpiproblem->oracle, i);
669
670 // every inequality contributes a slack variable
671 if( !SCIPisEQ(scip, lhs, rhs) )
672 ++nvars;
673
674 // every constraint lhs contributes a ln(slack-lhs) term in the barrier problem
675 if( !SCIPisInfinity(scip, -lhs) )
676 ++nvarbnd;
677 // every constraint rhs contributes a ln(rhs-slack) term in the barrier problem
678 if( !SCIPisInfinity(scip, rhs) )
679 ++nvarbnd;
680 }
681
682 const int* offset;
683 int nnlnz;
684 SCIP_CALL( SCIPnlpiOracleGetJacobianRowSparsity(scip, nlpiproblem->oracle, &offset, NULL, NULL, &nnlnz) );
685 jacnnz = offset[m];
686
687 /* fitting data from NLP runs gave the following coefficients (see also !2634):
688 * intercept +40.2756726751
689 * jacnnz -0.0021259769
690 * jacnnz_sqrt +2.0121042012
691 * nlincons -0.0374801925
692 * nlincons_sqrt +2.9562232443
693 * nnlcons -0.0133039200
694 * nnlcons_sqrt -0.0412118434
695 * nnlvars -0.0702890379
696 * nnlvars_sqrt +7.0920920430
697 * nvarbnd +0.0183592749
698 * nvarbnd_sqrt -4.7218258847
699 * nvars +0.0112944627
700 * nvars_sqrt -0.8365873360
701 */
702 param.iterlimit = SCIPfloor(scip, 40.2756726751
703 -0.0021259769 * jacnnz +2.0121042012 * sqrt(jacnnz)
704 -0.0374801925 * nlincons +2.9562232443 * sqrt(nlincons)
705 -0.0133039200 * nnlcons -0.0412118434 * sqrt(nnlcons)
706 -0.0702890379 * nnlvars +7.0920920430 * sqrt(nnlvars)
707 +0.0183592749 * nvarbnd -4.7218258847 * sqrt(nvarbnd)
708 +0.0112944627 * nvars -0.8365873360 * sqrt(nvars));
709 SCIPdebugMsg(scip, "Iteration limit guess: %d\n", param.iterlimit);
710 /* but with all the negative coefficients, let's also ensure some minimal number of iterations */
711 if( param.iterlimit < 50 )
712 param.iterlimit = 50;
713 }
714 if( nlpidata->print_level >= J_SUMMARY || param.verblevel > 0 )
715 SCIPinfoMessage(scip, NULL, "Chosen iteration limit to be %d\n", param.iterlimit);
716 }
717#endif
718 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("max_iter", param.iterlimit);
719
720 (void) nlpiproblem->ipopt->Options()->SetNumericValue("constr_viol_tol", FEASTOLFACTOR * param.feastol);
721 (void) nlpiproblem->ipopt->Options()->SetNumericValue("acceptable_constr_viol_tol", FEASTOLFACTOR * param.feastol);
722
723 /* set optimality tolerance parameters in Ipopt
724 *
725 * Sets dual_inf_tol and compl_inf_tol to opttol and tol to solvertol.
726 * We leave acceptable_dual_inf_tol and acceptable_compl_inf_tol untouched for now, which means that if Ipopt has convergence problems, then
727 * it can stop with a solution that is still feasible, but essentially without a proof of local optimality.
728 * Note, that in this case we report only feasibility and not optimality of the solution (see ScipNLP::finalize_solution).
729 */
730 (void) nlpiproblem->ipopt->Options()->SetNumericValue("dual_inf_tol", param.opttol);
731 (void) nlpiproblem->ipopt->Options()->SetNumericValue("compl_inf_tol", param.opttol);
732 if( param.solvertol > 0.0 )
733 (void) nlpiproblem->ipopt->Options()->SetNumericValue("tol", param.solvertol);
734 else
735#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR > 14 || (IPOPT_VERSION_MINOR == 14 && IPOPT_VERSION_RELEASE >= 2)
736 (void) nlpiproblem->ipopt->Options()->UnsetValue("tol");
737#else
738 (void) nlpiproblem->ipopt->Options()->SetNumericValue("tol", 1e-8); // 1e-8 is the default
739#endif
740
741 /* Ipopt doesn't like a setting of exactly 0 for the max_*_time, so increase as little as possible in that case */
742#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR >= 14
743 (void) nlpiproblem->ipopt->Options()->SetNumericValue("max_wall_time", MAX(param.timelimit, DBL_MIN));
744#else
745 (void) nlpiproblem->ipopt->Options()->SetNumericValue("max_cpu_time", MAX(param.timelimit, DBL_MIN));
746#endif
747
748 // disable acceptable-point heuristic iff fastfail is completely off
749 // by default (fastfail=conservative), it seems useful to have Ipopt stop when it obviously doesn't make progress (like one of the NLPs in the bendersqp ctest)
751 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("acceptable_iter", 0);
752 else
753#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR > 14 || (IPOPT_VERSION_MINOR == 14 && IPOPT_VERSION_RELEASE >= 2)
754 (void) nlpiproblem->ipopt->Options()->UnsetValue("acceptable_iter");
755#else
756 (void) nlpiproblem->ipopt->Options()->SetIntegerValue("acceptable_iter", 15); // 15 is the default
757#endif
758
759 (void) nlpiproblem->ipopt->Options()->SetStringValue("expect_infeasible_problem", param.expectinfeas ? "yes" : "no");
760
761 if( !nlpiproblem->ipopt->Options()->SetStringValue("warm_start_init_point", param.warmstart ? "yes" : "no") && !param.warmstart )
762 {
763 // if we cannot disable warmstarts in Ipopt, then we have a big problem
764 SCIPerrorMessage("Failed to set Ipopt warm_start_init_point option to no.");
765 return SCIP_ERROR;
766 }
767
768 return SCIP_OKAY;
769}
770
771#ifdef COLLECT_SOLVESTATS
772/// writes out some solve status, number of iterations, time, and problem properties
773static
774void collectStatistic(
775 SCIP* scip,
776 ApplicationReturnStatus status,
777 SCIP_NLPIPROBLEM* problem,
778 SmartPtr<SolveStatistics> stats
779 )
780{
781 int nvars = 0; // number of variables, including slacks
782 int nslacks = 0; // number of slack variables
783 int nnlvars = 0; // number of nonlinear variables
784 int nvarlb = 0; // number of finite lower bounds, including slacks
785 int nvarub = 0; // number of finite upper bounds, including slacks
786 int nlincons = 0; // number of linear constraints
787 int nnlcons = 0; // number of nonlinear constraints
788 int objnl = 0; // number of nonlinear objectives
789 int jacnnz = 0; // number of nonzeros in Jacobian
790 int hesnnz = 0; // number of nonzeros in Hessian of Lagrangian
791 int linsys11nz = 0; // number of nonzeros in linear system (11) solved by Ipopt
792 int linsys13nz = 0; // number of nonzeros in linear system (13) solved by Ipopt
793
794 int n = SCIPnlpiOracleGetNVars(problem->oracle);
795 int m = SCIPnlpiOracleGetNConstraints(problem->oracle);
796
797 for( int i = 0; i < n; ++i )
798 {
799 SCIP_Real lb, ub;
800
801 lb = SCIPnlpiOracleGetVarLbs(problem->oracle)[i];
802 ub = SCIPnlpiOracleGetVarUbs(problem->oracle)[i];
803
804 // skip fixed vars
805 if( SCIPisEQ(scip, lb, ub) )
806 continue;
807
808 ++nvars;
809
810 if( SCIPnlpiOracleIsVarNonlinear(scip, problem->oracle, i) )
811 ++nnlvars;
812
813 // every variable lower bound contributes a ln(x-lb) term in the barrier problem
814 if( !SCIPisInfinity(scip, -lb) )
815 ++nvarlb;
816
817 // every variable upper bound contributes a ln(ub-x) term in the barrier problem
818 if( !SCIPisInfinity(scip, ub) )
819 ++nvarub;
820 }
821
822 for( int i = 0; i < m; ++i )
823 {
824 SCIP_Real lhs, rhs;
825
827 ++nnlcons;
828 else
829 ++nlincons;
830
831 lhs = SCIPnlpiOracleGetConstraintLhs(problem->oracle, i);
832 rhs = SCIPnlpiOracleGetConstraintRhs(problem->oracle, i);
833
834 // every inequality contributes a slack variable
835 if( !SCIPisEQ(scip, lhs, rhs) )
836 {
837 ++nvars;
838 ++nslacks;
839 }
840
841 // every constraint lhs contributes a ln(slack-lhs) term in the barrier problem
842 if( !SCIPisInfinity(scip, -lhs) )
843 ++nvarlb;
844 // every constraint rhs contributes a ln(rhs-slack) term in the barrier problem
845 if( !SCIPisInfinity(scip, rhs) )
846 ++nvarub;
847 }
848
849 objnl = SCIPnlpiOracleIsConstraintNonlinear(problem->oracle, -1);
850
851 const int* offset;
852 const int* col;
853 int nnlnz;
855 jacnnz = offset[m];
856
858 hesnnz = offset[n];
859
860 // number of nonzeros of matrix in linear system of barrier problem ((11) in Ipopt paper):
861 // off-diagonal elements of Hessian of Lagrangian (W_k without diagonal)
862 // + number of variables (diagonal of W_k + Sigma_k)
863 // + 2*(elements in Jacobian + number of slacks) (A_k)
864 for( int i = 0; i < n; ++i )
865 for( int j = offset[i]; j < offset[i+1]; ++j )
866 if( col[j] != i )
867 linsys11nz += 2; // off-diagonal element of Lagrangian, once for each triangle of matrix
868 linsys11nz += nvars; // number of variables
869 linsys11nz += 2 * (jacnnz + nslacks); // because each slack var contributes one entry to the Jacobian
870
871 // number of nonzeros of matrix in perturbed linear system of barrier problem ((13) in Ipopt paper):
872 // linsys11nz + number of constraints
873 linsys13nz = linsys11nz + m;
874
875 SCIP_Real linsys11density = (SCIP_Real)linsys11nz / SQR((SCIP_Real)(nvars+m));
876 SCIP_Real linsys13density = (SCIP_Real)linsys13nz / SQR((SCIP_Real)(nvars+m));
877
878 bool expectinfeas;
879 problem->ipopt->Options()->GetBoolValue("expect_infeasible_problem", expectinfeas, "");
880
881 static bool firstwrite = true;
882 if( firstwrite )
883 {
884 printf("IPOPTSTAT status,iter,time,nvars,nnlvars,nvarlb,nvarub,nlincons,nnlcons,objnl,jacnnz,hesnnz,linsys11nz,linsys13nz,linsys11density,linsys13density,expectinfeas\n");
885 firstwrite = false;
886 }
887
888 printf("IPOPTSTAT %d,%d,%g,%d,%d,%d,%d,%d,%d,%d,%d,%d,%d,%d,%f,%f,%d\n",
889 status, stats->IterationCount(), stats->TotalWallclockTime(),
890 nvars, nnlvars, nvarlb, nvarub, nlincons, nnlcons, objnl, jacnnz, hesnnz, linsys11nz, linsys13nz, linsys11density, linsys13density, expectinfeas);
891}
892#endif
893
894/** copy method of NLP interface (called when SCIP copies plugins) */
895static
896SCIP_DECL_NLPICOPY(nlpiCopyIpopt)
897{
899
900 return SCIP_OKAY;
901}
902
903/** destructor of NLP interface to free nlpi data */
904static
905SCIP_DECL_NLPIFREE(nlpiFreeIpopt)
906{
907 assert(nlpidata != NULL);
908
909 delete *nlpidata;
910 *nlpidata = NULL;
911
912 return SCIP_OKAY;
913}
914
915/** gets pointer for NLP solver to do dirty stuff */
916static
917SCIP_DECL_NLPIGETSOLVERPOINTER(nlpiGetSolverPointerIpopt)
918{
919 if( problem == NULL )
920 return NULL;
921
922 return (void*)GetRawPtr(problem->ipopt);
923}
924
925/** creates a problem instance */
926static
927SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemIpopt)
928{
929 SCIP_NLPIDATA* data;
930
931 assert(nlpi != NULL);
932 assert(problem != NULL);
933
934 data = SCIPnlpiGetData(nlpi);
935 assert(data != NULL);
936
937 SCIP_ALLOC( *problem = new SCIP_NLPIPROBLEM ); /*lint !e774*/
938
939 SCIP_CALL( SCIPnlpiOracleCreate(scip, &(*problem)->oracle) );
940 SCIP_CALL( SCIPnlpiOracleSetProblemName(scip, (*problem)->oracle, name) );
941
942 try
943 {
944 /* initialize IPOPT without default journal */
945 (*problem)->ipopt = new IpoptApplication(false);
946
947 /* plugin our journal to get output through SCIP message handler */
948 SmartPtr<Journal> jrnl = new ScipJournal("console", J_ITERSUMMARY, scip);
949 jrnl->SetPrintLevel(J_DBG, J_NONE);
950 if( !(*problem)->ipopt->Jnlst()->AddJournal(jrnl) )
951 {
952 SCIPerrorMessage("Failed to register ScipJournal for IPOPT output.");
953 }
954
955 /* initialize Ipopt/SCIP NLP interface */
956 (*problem)->nlp = new ScipNLP(*problem, scip);
957 }
958 catch( const std::bad_alloc& )
959 {
960 SCIPerrorMessage("Not enough memory to initialize Ipopt.\n");
961 return SCIP_NOMEMORY;
962 }
963
964 for( size_t i = 0; i < sizeof(ipopt_string_params) / sizeof(const char*); ++i )
965 {
966 SCIP_PARAM* param;
968 char* paramval;
969
970 strcpy(paramname, "nlpi/" NLPI_NAME "/");
971 strcat(paramname, ipopt_string_params[i]);
972 param = SCIPgetParam(scip, paramname);
973
974 // skip parameters that we didn't add to SCIP because they didn't exist in this build of Ipopt
975 if( param == NULL )
976 continue;
977
978 // if value wasn't left at the default, then pass to Ipopt and forbid overwriting
979 paramval = SCIPparamGetString(param);
980 assert(paramval != NULL);
981 if( *paramval != '\0' )
982 (void) (*problem)->ipopt->Options()->SetStringValue(ipopt_string_params[i], paramval, false);
983 }
984
985 for( size_t i = 0; i < sizeof(ipopt_int_params) / sizeof(const char*); ++i )
986 {
987 SCIP_PARAM* param;
989 int paramval;
990
991 strcpy(paramname, "nlpi/" NLPI_NAME "/");
992 strcat(paramname, ipopt_int_params[i]);
993 param = SCIPgetParam(scip, paramname);
994
995 // skip parameters that we didn't add to SCIP because they didn't exist in this build of Ipopt
996 if( param == NULL )
997 continue;
998
999 // if value wasn't left at the default, then pass to Ipopt and forbid overwriting
1000 paramval = SCIPparamGetInt(param);
1001 if( paramval != SCIPparamGetIntDefault(param) )
1002 (void) (*problem)->ipopt->Options()->SetIntegerValue(ipopt_int_params[i], paramval, false);
1003 }
1004
1005#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
1006 /* Turn off bound relaxation for older Ipopt, as solutions may be out of bounds by more than constr_viol_tol.
1007 * For Ipopt 3.14, bounds are relaxed by at most constr_viol_tol, so can leave bound_relax_factor at its default.
1008 */
1009 (void) (*problem)->ipopt->Options()->SetNumericValue("bound_relax_factor", 0.0);
1010#endif
1011
1012 /* modify Ipopt's default settings to what we believe is appropriate */
1013 /* (*problem)->ipopt->Options()->SetStringValue("print_timing_statistics", "yes"); */
1014#ifdef SCIP_DEBUG
1015 (void) (*problem)->ipopt->Options()->SetStringValue("print_user_options", "yes");
1016#endif
1017 (void) (*problem)->ipopt->Options()->SetStringValue("sb", "yes");
1018 (void) (*problem)->ipopt->Options()->SetStringValueIfUnset("mu_strategy", "adaptive");
1019 (void) (*problem)->ipopt->Options()->SetIntegerValue("max_iter", INT_MAX);
1020 (void) (*problem)->ipopt->Options()->SetNumericValue("nlp_lower_bound_inf", -SCIPinfinity(scip), false);
1021 (void) (*problem)->ipopt->Options()->SetNumericValue("nlp_upper_bound_inf", SCIPinfinity(scip), false);
1022 (void) (*problem)->ipopt->Options()->SetNumericValue("diverging_iterates_tol", SCIPinfinity(scip), false);
1023
1024 /* when warmstarting, then reduce how much Ipopt modified the starting point */
1025 (void) (*problem)->ipopt->Options()->SetNumericValue("warm_start_bound_push", data->warm_start_push);
1026 (void) (*problem)->ipopt->Options()->SetNumericValue("warm_start_bound_frac", data->warm_start_push);
1027 (void) (*problem)->ipopt->Options()->SetNumericValue("warm_start_slack_bound_push", data->warm_start_push);
1028 (void) (*problem)->ipopt->Options()->SetNumericValue("warm_start_slack_bound_frac", data->warm_start_push);
1029 (void) (*problem)->ipopt->Options()->SetNumericValue("warm_start_mult_bound_push", data->warm_start_push);
1030
1031 /* apply user's given options file */
1032 assert(data->optfile != NULL);
1033 if( (*problem)->ipopt->Initialize(data->optfile) != Solve_Succeeded )
1034 {
1035 SCIPerrorMessage("Error during initialization of Ipopt using optionfile \"%s\"\n", data->optfile);
1036 return SCIP_ERROR;
1037 }
1038
1039 return SCIP_OKAY;
1040}
1041
1042/** free a problem instance */
1043static
1044SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemIpopt)
1045{
1046 int n;
1047 int m;
1048
1049 assert(nlpi != NULL);
1050 assert(problem != NULL);
1051 assert(*problem != NULL);
1052 assert((*problem)->oracle != NULL);
1053
1054 n = SCIPnlpiOracleGetNVars((*problem)->oracle);
1055 m = SCIPnlpiOracleGetNConstraints((*problem)->oracle);
1056
1057 SCIPfreeBlockMemoryArrayNull(scip, &(*problem)->solprimals, n);
1058 SCIPfreeBlockMemoryArrayNull(scip, &(*problem)->soldualcons, m);
1059 SCIPfreeBlockMemoryArrayNull(scip, &(*problem)->soldualvarlb, n);
1060 SCIPfreeBlockMemoryArrayNull(scip, &(*problem)->soldualvarub, n);
1061
1062 SCIP_CALL( SCIPnlpiOracleFree(scip, &(*problem)->oracle) );
1063
1064 if( (*problem)->randnumgen != NULL )
1065 {
1066 SCIPfreeRandom(scip, &(*problem)->randnumgen);
1067 }
1068
1069#ifdef PRINT_NLPSTATS
1070 SCIPinfoMessage(scip, NULL, "\nNLP solver IPOPT stats: ncalls = %d, nsuccess = %d, nlimit = %d, nlocinfeas = %d, nother = %d\n",
1071 (*problem)->ncalls, (*problem)->nsuccess, (*problem)->nlimit, (*problem)->nlocinfeas, (*problem)->nother);
1072#endif
1073
1074 delete *problem;
1075 *problem = NULL;
1076
1077 return SCIP_OKAY;
1078}
1079
1080/** gets pointer to solver-internal problem instance to do dirty stuff */
1081static
1082SCIP_DECL_NLPIGETPROBLEMPOINTER(nlpiGetProblemPointerIpopt)
1083{
1084 assert(nlpi != NULL);
1085 assert(problem != NULL);
1086
1087 return GetRawPtr(problem->nlp);
1088}
1089
1090/** add variables */
1091static
1092SCIP_DECL_NLPIADDVARS(nlpiAddVarsIpopt)
1093{
1094 int oldnvars;
1095
1096 assert(nlpi != NULL);
1097 assert(problem != NULL);
1098 assert(problem->oracle != NULL);
1099
1100 oldnvars = SCIPnlpiOracleGetNVars(problem->oracle);
1101
1102 SCIPfreeBlockMemoryArrayNull(scip, &problem->solprimals, oldnvars);
1103 SCIPfreeBlockMemoryArrayNull(scip, &problem->soldualvarlb, oldnvars);
1104 SCIPfreeBlockMemoryArrayNull(scip, &problem->soldualvarub, oldnvars);
1105 invalidateSolution(problem);
1106
1107 SCIP_CALL( SCIPnlpiOracleAddVars(scip, problem->oracle, nvars, lbs, ubs, varnames) );
1108
1109 problem->samestructure = false;
1110
1111 return SCIP_OKAY;
1112}
1113
1114/** add constraints */
1115static
1116SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsIpopt)
1117{
1118 int oldncons;
1119
1120 assert(nlpi != NULL);
1121 assert(problem != NULL);
1122 assert(problem->oracle != NULL);
1123
1124 oldncons = SCIPnlpiOracleGetNConstraints(problem->oracle);
1125
1126 SCIPfreeBlockMemoryArrayNull(scip, &problem->soldualcons, oldncons);
1127 problem->soldualvalid = false;
1128 problem->soldualgiven = false;
1129
1130 SCIP_CALL( SCIPnlpiOracleAddConstraints(scip, problem->oracle, nconss, lhss, rhss, nlininds, lininds, linvals, exprs, names) );
1131
1132 problem->samestructure = false;
1133
1134 return SCIP_OKAY;
1135}
1136
1137/** sets or overwrites objective, a minimization problem is expected */
1138static
1139SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveIpopt)
1140{
1141 assert(nlpi != NULL);
1142 assert(problem != NULL);
1143 assert(problem->oracle != NULL);
1144
1145 /* We pass the objective gradient in dense form to Ipopt, so if the sparsity of that gradient changes, we do not change the structure of the problem inside Ipopt.
1146 * However, if the sparsity of the Hessian matrix of the objective changes, then the sparsity pattern of the Hessian of the Lagrangian may change.
1147 * Thus, set samestructure=false if the objective was and/or becomes nonlinear, but leave samestructure untouched if it was and stays linear.
1148 */
1149 if( expr != NULL || SCIPnlpiOracleIsConstraintNonlinear(problem->oracle, -1) )
1150 problem->samestructure = false;
1151
1152 SCIP_CALL( SCIPnlpiOracleSetObjective(scip, problem->oracle, constant, nlins, lininds, linvals, expr) );
1153
1154 /* keep solution as valid, but reset solve status and objective value */
1155 invalidateSolved(problem);
1156
1157 return SCIP_OKAY;
1158}
1159
1160/** change variable bounds */
1161static
1162SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsIpopt)
1163{
1164 assert(nlpi != NULL);
1165 assert(problem != NULL);
1166 assert(problem->oracle != NULL);
1167
1168 /* Check whether the structure of the Ipopt internal NLP changes, if problem->samestructure at the moment.
1169 * We need to check whether variables become fixed or unfixed and whether bounds are added or removed.
1170 *
1171 * Project primal solution onto new bounds if currently valid.
1172 */
1173 if( problem->samestructure || problem->solprimalvalid )
1174 {
1175 for( int i = 0; i < nvars; ++i )
1176 {
1177 SCIP_Real oldlb;
1178 SCIP_Real oldub;
1179 oldlb = SCIPnlpiOracleGetVarLbs(problem->oracle)[indices[i]];
1180 oldub = SCIPnlpiOracleGetVarUbs(problem->oracle)[indices[i]];
1181
1182 if( (oldlb == oldub) != (lbs[i] == ubs[i]) ) /*lint !e777*/
1183 problem->samestructure = false;
1184 else if( SCIPisInfinity(scip, -oldlb) != SCIPisInfinity(scip, -lbs[i]) )
1185 problem->samestructure = false;
1186 else if( SCIPisInfinity(scip, oldub) != SCIPisInfinity(scip, ubs[i]) )
1187 problem->samestructure = false;
1188
1189 if( problem->solprimalvalid )
1190 {
1191 assert(problem->solprimals != NULL);
1192 problem->solprimals[i] = MIN(MAX(problem->solprimals[indices[i]], lbs[i]), ubs[i]);
1193 }
1194 }
1195 }
1196
1197 SCIP_CALL( SCIPnlpiOracleChgVarBounds(scip, problem->oracle, nvars, indices, lbs, ubs) );
1198
1199 invalidateSolved(problem);
1200
1201 return SCIP_OKAY;
1202}
1203
1204/** change constraint bounds */
1205static
1206SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesIpopt)
1207{
1208 assert(nlpi != NULL);
1209 assert(problem != NULL);
1210 assert(problem->oracle != NULL);
1211
1212 /* Check whether the structure of the Ipopt internal NLP changes, if problem->samestructure at the moment.
1213 * We need to check whether constraints change from equality to inequality and whether sides are added or removed.
1214 */
1215 for( int i = 0; i < nconss && problem->samestructure; ++i )
1216 {
1217 SCIP_Real oldlhs;
1218 SCIP_Real oldrhs;
1219 oldlhs = SCIPnlpiOracleGetConstraintLhs(problem->oracle, indices[i]);
1220 oldrhs = SCIPnlpiOracleGetConstraintRhs(problem->oracle, indices[i]);
1221
1222 if( (oldlhs == oldrhs) != (lhss[i] == rhss[i]) ) /*lint !e777*/
1223 problem->samestructure = false;
1224 else if( SCIPisInfinity(scip, -oldlhs) != SCIPisInfinity(scip, -lhss[i]) )
1225 problem->samestructure = false;
1226 else if( SCIPisInfinity(scip, oldrhs) != SCIPisInfinity(scip, rhss[i]) )
1227 problem->samestructure = false;
1228 }
1229
1230 SCIP_CALL( SCIPnlpiOracleChgConsSides(scip, problem->oracle, nconss, indices, lhss, rhss) );
1231
1232 invalidateSolved(problem);
1233
1234 return SCIP_OKAY;
1235}
1236
1237/** delete a set of variables */
1238static
1239SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetIpopt)
1240{
1241 int nvars;
1242
1243 assert(nlpi != NULL);
1244 assert(problem != NULL);
1245 assert(problem->oracle != NULL);
1246 assert(SCIPnlpiOracleGetNVars(problem->oracle) == dstatssize);
1247
1248 SCIP_CALL( SCIPnlpiOracleDelVarSet(scip, problem->oracle, dstats) );
1249
1250 nvars = SCIPnlpiOracleGetNVars(problem->oracle);
1251
1252 if( problem->solprimalvalid || problem->soldualvalid )
1253 {
1254 // update existing solution, if valid
1255 assert(!problem->solprimalvalid || problem->solprimals != NULL);
1256 assert(!problem->soldualvalid || problem->soldualvarlb != NULL);
1257 assert(!problem->soldualvalid || problem->soldualvarub != NULL);
1258
1259 int i;
1260 for( i = 0; i < dstatssize; ++i )
1261 {
1262 if( dstats[i] != -1 )
1263 {
1264 assert(dstats[i] >= 0);
1265 assert(dstats[i] < nvars);
1266 if( problem->solprimals != NULL )
1267 problem->solprimals[dstats[i]] = problem->solprimals[i];
1268 if( problem->soldualvarlb != NULL )
1269 {
1270 assert(problem->soldualvarub != NULL);
1271 problem->soldualvarlb[dstats[i]] = problem->soldualvarlb[i];
1272 problem->soldualvarub[dstats[i]] = problem->soldualvarub[i];
1273 }
1274 }
1275 }
1276 }
1277
1278 /* resize solution point arrays */
1279 if( problem->solprimals != NULL )
1280 {
1281 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &problem->solprimals, dstatssize, nvars) );
1282 }
1283 if( problem->soldualvarlb != NULL )
1284 {
1285 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &problem->soldualvarlb, dstatssize, nvars) );
1286 }
1287 if( problem->soldualvarub != NULL )
1288 {
1289 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &problem->soldualvarub, dstatssize, nvars) );
1290 }
1291
1292 problem->samestructure = false;
1293
1294 invalidateSolved(problem);
1295
1296 return SCIP_OKAY;
1297}
1298
1299/** delete a set of constraints */
1300static
1301SCIP_DECL_NLPIDELCONSSET(nlpiDelConstraintSetIpopt)
1302{
1303 int ncons;
1304
1305 assert(nlpi != NULL);
1306 assert(problem != NULL);
1307 assert(problem->oracle != NULL);
1308 assert(SCIPnlpiOracleGetNConstraints(problem->oracle) == dstatssize);
1309
1310 SCIP_CALL( SCIPnlpiOracleDelConsSet(scip, problem->oracle, dstats) );
1311
1312 ncons = SCIPnlpiOracleGetNConstraints(problem->oracle);
1313
1314 if( problem->soldualvalid )
1315 {
1316 // update existing dual solution
1317 assert(problem->soldualcons != NULL);
1318
1319 int i;
1320 for( i = 0; i < dstatssize; ++i )
1321 {
1322 if( dstats[i] != -1 )
1323 {
1324 assert(dstats[i] >= 0);
1325 assert(dstats[i] < ncons);
1326 problem->soldualcons[dstats[i]] = problem->soldualcons[i];
1327 }
1328 }
1329 }
1330
1331 /* resize dual solution point array */
1332 if( problem->soldualcons != NULL )
1333 {
1334 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &problem->soldualcons, dstatssize, ncons) );
1335 }
1336
1337 problem->samestructure = false;
1338
1339 invalidateSolved(problem);
1340
1341 return SCIP_OKAY;
1342}
1343
1344/** change one linear coefficient in a constraint or objective */
1345static
1346SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsIpopt)
1347{
1348 assert(nlpi != NULL);
1349 assert(problem != NULL);
1350 assert(problem->oracle != NULL);
1351
1352 SCIP_CALL( SCIPnlpiOracleChgLinearCoefs(scip, problem->oracle, idx, nvals, varidxs, vals) );
1353
1354 problem->samestructure = false; // nonzero patterns may have changed; TODO SCIPnlpiOracleChgLinearCoefs() should let us know
1355 invalidateSolved(problem);
1356
1357 return SCIP_OKAY;
1358}
1359
1360/** replaces the expression tree of a constraint or objective */
1361static
1362SCIP_DECL_NLPICHGEXPR(nlpiChgExprIpopt)
1363{
1364 assert(nlpi != NULL);
1365 assert(problem != NULL);
1366 assert(problem->oracle != NULL);
1367
1368 SCIP_CALL( SCIPnlpiOracleChgExpr(scip, problem->oracle, idxcons, expr) );
1369
1370 problem->samestructure = false; // nonzero patterns may have changed
1371 invalidateSolved(problem);
1372
1373 return SCIP_OKAY;
1374}
1375
1376/** change the constant offset in the objective */
1377static
1378SCIP_DECL_NLPICHGOBJCONSTANT(nlpiChgObjConstantIpopt)
1379{
1380 SCIP_Real oldconstant;
1381
1382 assert(nlpi != NULL);
1383 assert(problem != NULL);
1384 assert(problem->oracle != NULL);
1385
1386 oldconstant = SCIPnlpiOracleGetObjectiveConstant(problem->oracle);
1387
1388 SCIP_CALL( SCIPnlpiOracleChgObjConstant(scip, problem->oracle, objconstant) );
1389
1390 if( problem->solobjval != SCIP_INVALID )
1391 problem->solobjval += objconstant - oldconstant;
1392
1393 return SCIP_OKAY;
1394}
1395
1396/** sets initial guess for primal variables */
1397static
1398SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessIpopt)
1399{
1400 int nvars;
1401
1402 assert(nlpi != NULL);
1403 assert(problem != NULL);
1404 assert(problem->oracle != NULL);
1405
1406 nvars = SCIPnlpiOracleGetNVars(problem->oracle);
1407
1408 if( primalvalues != NULL )
1409 {
1410 // copy primal solution
1411 SCIPdebugMsg(scip, "set initial guess primal values to user-given\n");
1412 if( problem->solprimals == NULL )
1413 {
1414 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &problem->solprimals, nvars) );
1415 }
1416 BMScopyMemoryArray(problem->solprimals, primalvalues, nvars);
1417 problem->solprimalvalid = true;
1418 problem->solprimalgiven = true;
1419 }
1420 else
1421 {
1422 // invalid current primal solution (if any)
1423 if( problem->solprimalvalid )
1424 {
1425 SCIPdebugMsg(scip, "invalidate initial guess primal values on user-request\n");
1426 }
1427 problem->solprimalvalid = false;
1428 problem->solprimalgiven = false;
1429 }
1430
1431 if( consdualvalues != NULL && varlbdualvalues != NULL && varubdualvalues != NULL )
1432 {
1433 // copy dual solution, if completely given
1434 SCIPdebugMsg(scip, "set initial guess dual values to user-given\n");
1435 int ncons = SCIPnlpiOracleGetNConstraints(problem->oracle);
1436 if( problem->soldualcons == NULL )
1437 {
1439 }
1440 BMScopyMemoryArray(problem->soldualcons, consdualvalues, ncons);
1441
1442 assert((problem->soldualvarlb == NULL) == (problem->soldualvarub == NULL));
1443 if( problem->soldualvarlb == NULL )
1444 {
1447 }
1448 BMScopyMemoryArray(problem->soldualvarlb, varlbdualvalues, nvars);
1449 BMScopyMemoryArray(problem->soldualvarub, varubdualvalues, nvars);
1450 problem->soldualvalid = true;
1451 problem->soldualgiven = true;
1452 }
1453 else
1454 {
1455 // invalid current dual solution (if any)
1456 if( problem->soldualvalid )
1457 {
1458 SCIPdebugMsg(scip, "invalidate initial guess dual values\n");
1459 }
1460 problem->soldualvalid = false;
1461 problem->soldualgiven = false;
1462 }
1463
1464 return SCIP_OKAY;
1465}
1466
1467/** tries to solve NLP */
1468static
1469SCIP_DECL_NLPISOLVE(nlpiSolveIpopt)
1470{
1471 SCIP_NLPIDATA* nlpidata;
1472 ApplicationReturnStatus status;
1473
1474 assert(nlpi != NULL);
1475 assert(problem != NULL);
1476 assert(problem->oracle != NULL);
1477
1478 assert(IsValid(problem->ipopt));
1479 assert(IsValid(problem->nlp));
1480
1481 nlpidata = SCIPnlpiGetData(nlpi);
1482 assert(nlpidata != NULL);
1483
1484 // print parameters if either nlpi/ipopt/print_level has been set high enough or solve called with verblevel>0
1485 if( nlpidata->print_level >= J_SUMMARY || param.verblevel > 0 )
1486 {
1487 SCIPinfoMessage(scip, NULL, "Ipopt solve for problem %s at subSCIP depth %d", SCIPnlpiOracleGetProblemName(problem->oracle), SCIPgetSubscipDepth(scip));
1488 SCIPinfoMessage(scip, NULL, " with parameters " SCIP_NLPPARAM_PRINT(param));
1489 }
1490
1492
1493 if( param.timelimit == 0.0 )
1494 {
1495 /* there is nothing we can do if we are not given any time */
1496 problem->lastniter = 0;
1497 problem->lasttime = 0.0;
1500
1501 return SCIP_OKAY;
1502 }
1503
1504#ifdef COLLECT_SOLVESTATS
1505 // if collecting statistics on how many iterations one may need for a certain Ipopt problem,
1506 // do not use a tight iteration limit (but use some limit to get finished)
1507 param.iterlimit = 1000;
1508#endif
1509
1510 // change status info to unsolved, just in case
1511 invalidateSolved(problem);
1512
1513 // ensure a starting point is available
1514 // also disables param.warmstart if no warmstart available
1515 SCIP_CALL( ensureStartingPoint(scip, problem, param.warmstart) );
1516
1517 // tell NLP that we are about to start a new solve
1518 problem->nlp->initializeSolve(problem, param);
1519
1520 // set Ipopt parameters
1521 SCIP_CALL( handleNlpParam(scip, nlpidata, problem, param) );
1522
1523 try
1524 {
1525#ifdef PROTECT_SOLVE_BY_MUTEX
1526 /* lock solve_mutex if Ipopt is going to use Mumps as linear solver
1527 * unlocking will happen in the destructor of guard, which is called when this block is left
1528 */
1529 std::unique_lock<std::mutex> guard(solve_mutex, std::defer_lock); /*lint !e{728}*/
1530 std::string linsolver;
1531 (void) problem->ipopt->Options()->GetStringValue("linear_solver", linsolver, "");
1532 if( linsolver == "mumps" )
1533 guard.lock();
1534#endif
1535
1536 if( problem->firstrun )
1537 {
1539
1541
1542 /* if the expression interpreter or some user expression do not support function values and gradients and Hessians,
1543 * change NLP parameters or give an error
1544 */
1546 {
1549 {
1550 SCIPerrorMessage("Do not have expression interpreter that can compute function values and gradients. Cannot solve NLP with Ipopt.\n");
1553 return SCIP_OKAY;
1554 }
1555
1556 /* enable Hessian approximation if we are nonquadratic and the expression interpreter or user expression do not support Hessians */
1557 if( !(cap & SCIP_EXPRINTCAPABILITY_HESSIAN) )
1558 {
1559 (void) problem->ipopt->Options()->SetStringValueIfUnset("hessian_approximation", "limited-memory");
1560 problem->nlp->approxhessian = true;
1561 }
1562 else
1563 problem->nlp->approxhessian = false;
1564 }
1565
1566#ifdef SCIP_DEBUG
1567 problem->ipopt->Options()->SetStringValue("derivative_test", problem->nlp->approxhessian ? "first-order" : "second-order");
1568#endif
1569
1570 status = problem->ipopt->OptimizeTNLP(GetRawPtr(problem->nlp));
1571 }
1572 else
1573 {
1574 // TODO to be strict, we should check whether the eval capability has been changed and the Hessian approximation needs to be enabled (in which case we should call OptimizeTNLP instead)
1575 problem->ipopt->Options()->SetStringValue("warm_start_same_structure", problem->samestructure ? "yes" : "no");
1576 status = problem->ipopt->ReOptimizeTNLP(GetRawPtr(problem->nlp));
1577 }
1578
1579 // catch the very bad status codes
1580 switch( status ) {
1581 // everything better than Not_Enough_Degrees_Of_Freedom is a non-serious error
1582 case Solve_Succeeded:
1583 case Solved_To_Acceptable_Level:
1584 case Infeasible_Problem_Detected:
1585 case Search_Direction_Becomes_Too_Small:
1586 case Diverging_Iterates:
1587 case User_Requested_Stop:
1588 case Feasible_Point_Found:
1589 case Maximum_Iterations_Exceeded:
1590 case Restoration_Failed:
1591 case Error_In_Step_Computation:
1592 case Maximum_CpuTime_Exceeded:
1593#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR >= 14
1594 case Maximum_WallTime_Exceeded: // new in Ipopt 3.14
1595 // if Ipopt >= 3.14, finalize_solution should always have been called if we get these status codes
1596 // this should have left us with some solution (unless we ran out of memory in finalize_solution)
1597 assert(problem->solprimalvalid || problem->termstat == SCIP_NLPTERMSTAT_OUTOFMEMORY);
1598 assert(problem->soldualvalid || problem->termstat == SCIP_NLPTERMSTAT_OUTOFMEMORY);
1599#endif
1600 problem->firstrun = false;
1601 problem->samestructure = true;
1602 break;
1603
1604 case Not_Enough_Degrees_Of_Freedom:
1605 assert(problem->termstat == SCIP_NLPTERMSTAT_OTHER);
1606 assert(problem->solstat == SCIP_NLPSOLSTAT_UNKNOWN);
1607 SCIPdebugMsg(scip, "NLP has too few degrees of freedom.\n");
1608 break;
1609
1610 case Invalid_Number_Detected:
1611 SCIPdebugMsg(scip, "Ipopt failed because of an invalid number in function or derivative value\n");
1613 /* Ipopt may or may not have called finalize solution
1614 * if it didn't, then we should still have SCIP_NLPSOLSTAT_UNKNOWN as set in the invalidateSolved() call above
1615 * if it did, then finalize_solution will have set SCIP_NLPSOLSTAT_UNKNOWN or SCIP_NLPSOLSTAT_FEASIBLE
1616 */
1617 assert(problem->solstat == SCIP_NLPSOLSTAT_UNKNOWN || problem->solstat == SCIP_NLPSOLSTAT_FEASIBLE);
1618 break;
1619
1620 case Insufficient_Memory:
1621 assert(problem->termstat == SCIP_NLPTERMSTAT_OTHER);
1622 assert(problem->solstat == SCIP_NLPSOLSTAT_UNKNOWN);
1623 SCIPerrorMessage("Ipopt returned with status \"Insufficient Memory\"\n");
1624 return SCIP_NOMEMORY;
1625
1626 // really bad ones that could be something very unexpected going wrong within Ipopt
1627 case Unrecoverable_Exception:
1628 case Internal_Error:
1629 assert(problem->termstat == SCIP_NLPTERMSTAT_OTHER);
1630 assert(problem->solstat == SCIP_NLPSOLSTAT_UNKNOWN);
1631 SCIPerrorMessage("Ipopt returned with application return status %d\n", status);
1632 break;
1633
1634 // the really bad ones that indicate rather a programming error
1635 case Invalid_Problem_Definition:
1636 case Invalid_Option:
1637 case NonIpopt_Exception_Thrown:
1638 assert(problem->termstat == SCIP_NLPTERMSTAT_OTHER);
1639 assert(problem->solstat == SCIP_NLPSOLSTAT_UNKNOWN);
1640 SCIPerrorMessage("Ipopt returned with application return status %d\n", status);
1641 return SCIP_ERROR;
1642 }
1643
1644#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
1645 SmartPtr<SolveStatistics> stats = problem->ipopt->Statistics();
1646 /* Ipopt does not provide access to the statistics if there was a serious error */
1647 if( IsValid(stats) )
1648 {
1649 problem->lastniter = stats->IterationCount();
1650 problem->lasttime = stats->TotalWallclockTime();
1651
1652#ifdef COLLECT_SOLVESTATS
1653 collectStatistic(scip, status, problem, stats);
1654#endif
1655 }
1656#else
1657 SmartPtr<IpoptData> ip_data = problem->ipopt->IpoptDataObject();
1658 /* I don't think that there is a situation where ip_data is NULL, but check here anyway */
1659 if( IsValid(ip_data) )
1660 {
1661 problem->lastniter = ip_data->iter_count();
1662 problem->lasttime = ip_data->TimingStats().OverallAlgorithm().TotalWallclockTime();
1663 }
1664#endif
1665 else
1666 {
1667 problem->lastniter = 0;
1668 problem->lasttime = 0.0;
1669 }
1670 }
1671 catch( IpoptException& except )
1672 {
1673 SCIPerrorMessage("Ipopt returned with exception: %s\n", except.Message().c_str());
1674 return SCIP_ERROR;
1675 }
1676
1677 return SCIP_OKAY;
1678}
1679
1680/** gives solution status */
1681static
1682SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatIpopt)
1683{
1684 assert(nlpi != NULL);
1685 assert(problem != NULL);
1686
1687 return problem->solstat;
1688}
1689
1690/** gives termination reason */
1691static
1692SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatIpopt)
1693{
1694 assert(nlpi != NULL);
1695 assert(problem != NULL);
1696
1697 return problem->termstat;
1698}
1699
1700/** gives primal and dual solution values */
1701static
1702SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionIpopt)
1703{
1704 assert(nlpi != NULL);
1705 assert(problem != NULL);
1706
1707 if( primalvalues != NULL )
1708 *primalvalues = problem->solprimals;
1709
1710 if( consdualvalues != NULL )
1711 *consdualvalues = problem->soldualcons;
1712
1713 if( varlbdualvalues != NULL )
1714 *varlbdualvalues = problem->soldualvarlb;
1715
1716 if( varubdualvalues != NULL )
1717 *varubdualvalues = problem->soldualvarub;
1718
1719 if( objval != NULL )
1720 *objval = problem->solobjval;
1721
1722 return SCIP_OKAY;
1723}
1724
1725/** gives solve statistics */
1726static
1727SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsIpopt)
1728{
1729 assert(nlpi != NULL);
1730 assert(problem != NULL);
1731 assert(statistics != NULL);
1732
1733 statistics->niterations = problem->lastniter;
1734 statistics->totaltime = problem->lasttime;
1735 statistics->evaltime = SCIPnlpiOracleGetEvalTime(scip, problem->oracle);
1736 statistics->consviol = problem->solconsviol;
1737 statistics->boundviol = problem->solboundviol;
1738
1739 return SCIP_OKAY;
1740}
1741
1742/** create solver interface for Ipopt solver and includes it into SCIP, if Ipopt is available */
1744 SCIP* scip /**< SCIP data structure */
1745 )
1746{
1747 SCIP_NLPIDATA* nlpidata;
1748 SCIP_Bool advanced = FALSE;
1749
1750 assert(scip != NULL);
1751
1752 SCIP_ALLOC( nlpidata = new SCIP_NLPIDATA() ); /*lint !e774*/
1753
1756 nlpiCopyIpopt, nlpiFreeIpopt, nlpiGetSolverPointerIpopt,
1757 nlpiCreateProblemIpopt, nlpiFreeProblemIpopt, nlpiGetProblemPointerIpopt,
1758 nlpiAddVarsIpopt, nlpiAddConstraintsIpopt, nlpiSetObjectiveIpopt,
1759 nlpiChgVarBoundsIpopt, nlpiChgConsSidesIpopt, nlpiDelVarSetIpopt, nlpiDelConstraintSetIpopt,
1760 nlpiChgLinearCoefsIpopt, nlpiChgExprIpopt,
1761 nlpiChgObjConstantIpopt, nlpiSetInitialGuessIpopt, nlpiSolveIpopt, nlpiGetSolstatIpopt, nlpiGetTermstatIpopt,
1762 nlpiGetSolutionIpopt, nlpiGetStatisticsIpopt,
1763 nlpidata) );
1764
1766
1767 SCIP_CALL( SCIPaddStringParam(scip, "nlpi/" NLPI_NAME "/optfile", "name of Ipopt options file",
1768 &nlpidata->optfile, FALSE, "", NULL, NULL) );
1769
1770 SCIP_CALL( SCIPaddRealParam(scip, "nlpi/" NLPI_NAME "/warm_start_push", "amount (relative and absolute) by which starting point is moved away from bounds in warmstarts",
1771 &nlpidata->warm_start_push, FALSE, 1e-9, 0.0, 1.0, NULL, NULL) );
1772
1773 SmartPtr<RegisteredOptions> reg_options = new RegisteredOptions();
1774 IpoptApplication::RegisterAllIpoptOptions(reg_options);
1775
1776 for( size_t i = 0; i < sizeof(ipopt_string_params) / sizeof(const char*); ++i )
1777 {
1778 SmartPtr<const RegisteredOption> option = reg_options->GetOption(ipopt_string_params[i]);
1779
1780 // skip options not available with this build of Ipopt
1781 if( !IsValid(option) )
1782 continue;
1783
1784 assert(option->Type() == OT_String);
1785
1786 // prefix parameter name with nlpi/ipopt
1787 std::string paramname("nlpi/" NLPI_NAME "/");
1788 paramname += option->Name();
1789
1790 // initialize description with short description from Ipopt
1791 std::stringstream descr;
1792 descr << option->ShortDescription();
1793
1794 // add valid values to description, if there are more than one
1795 // the only case where there are less than 2 valid strings should be when anything is valid (in which case there is one valid string with value "*")
1796 std::vector<RegisteredOption::string_entry> validvals = option->GetValidStrings();
1797 if( validvals.size() > 1 )
1798 {
1799 descr << " Valid values if not empty:";
1800 for( std::vector<RegisteredOption::string_entry>::iterator val = validvals.begin(); val != validvals.end(); ++val )
1801 descr << ' ' << val->value_;
1802 }
1803
1804#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR >= 14
1805 // since Ipopt 3.14, Ipopt options have an advanced flag
1806 advanced = option->Advanced();
1807#endif
1808
1809 // we use the empty string as default to recognize later whether the user set has set the option
1810 SCIP_CALL( SCIPaddStringParam(scip, paramname.c_str(), descr.str().c_str(), NULL, advanced, "", NULL, NULL) );
1811 }
1812
1813 for( size_t i = 0; i < sizeof(ipopt_int_params) / sizeof(const char*); ++i )
1814 {
1815 assert(i > 0 || strcmp(ipopt_int_params[0], "print_level") == 0); // we assume print_level at index 0
1816
1817 SmartPtr<const RegisteredOption> option = reg_options->GetOption(ipopt_int_params[i]);
1818
1819 // skip options not available with this build of Ipopt
1820 if( !IsValid(option) )
1821 continue;
1822
1823 assert(option->Type() == OT_Integer);
1824
1825 // prefix parameter name with nlpi/ipopt
1826 std::string paramname("nlpi/" NLPI_NAME "/");
1827 paramname += option->Name();
1828
1829 int lower = option->LowerInteger();
1830 int upper = option->UpperInteger();
1831
1832 // we use value lower-1 as signal that the option was not modified by the user
1833 // for that, we require a finite lower bound
1834 assert(lower > INT_MIN);
1835
1836 // initialize description with short description from Ipopt
1837 std::stringstream descr;
1838 descr << option->ShortDescription();
1839 descr << ' ' << (lower-1) << " to use NLPI or Ipopt default.";
1840
1841#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR >= 14
1842 // since Ipopt 3.14, Ipopt options have an advanced flag
1843 advanced = option->Advanced();
1844#endif
1845
1846 // we use the empty string as default to recognize later whether the user set has set the option
1847 SCIP_CALL( SCIPaddIntParam(scip, paramname.c_str(), descr.str().c_str(),
1848 i == 0 ? &nlpidata->print_level : NULL, advanced,
1849 lower-1, lower-1, upper, NULL, NULL) );
1850 }
1851
1852 return SCIP_OKAY;
1853} /*lint !e429 */
1854
1855/** gets string that identifies Ipopt (version number) */
1856const char* SCIPgetSolverNameIpopt(void)
1857{
1858 return "Ipopt " IPOPT_VERSION;
1859}
1860
1861/** gets string that describes Ipopt */
1862const char* SCIPgetSolverDescIpopt(void)
1863{
1864 return "Interior Point Optimizer developed by A. Waechter et.al. (github.com/coin-or/Ipopt)";
1865}
1866
1867/** returns whether Ipopt is available, i.e., whether it has been linked in */
1869{
1870 return TRUE;
1871}
1872
1873/** gives a pointer to the NLPIORACLE object stored in Ipopt-NLPI's NLPI problem data structure */
1875 SCIP_NLPIPROBLEM* nlpiproblem /**< NLP problem of Ipopt-NLPI */
1876 )
1877{
1878 assert(nlpiproblem != NULL);
1879
1880 return nlpiproblem->oracle;
1881}
1882
1883/** Method to return some info about the nlp */
1884bool ScipNLP::get_nlp_info(
1885 Index& n, /**< place to store number of variables */
1886 Index& m, /**< place to store number of constraints */
1887 Index& nnz_jac_g, /**< place to store number of nonzeros in jacobian */
1888 Index& nnz_h_lag, /**< place to store number of nonzeros in hessian */
1889 IndexStyleEnum& index_style /**< place to store used index style (0-based or 1-based) */
1890 )
1891{
1892 const int* offset;
1893 SCIP_RETCODE retcode;
1894 int nnlnz;
1895
1896 assert(nlpiproblem != NULL);
1897 assert(nlpiproblem->oracle != NULL);
1898
1899 n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
1900 m = SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle);
1901
1902 retcode = SCIPnlpiOracleGetJacobianRowSparsity(scip, nlpiproblem->oracle, &offset, NULL, NULL, &nnlnz);
1903 if( retcode != SCIP_OKAY )
1904 return false;
1905
1906 nnz_jac_g = offset == NULL ? 0 : offset[m];
1907
1908 if( !approxhessian )
1909 {
1910 retcode = SCIPnlpiOracleGetHessianLagSparsity(scip, nlpiproblem->oracle, &offset, NULL, FALSE);
1911 if( retcode != SCIP_OKAY )
1912 return false;
1913 assert(offset != NULL);
1914 nnz_h_lag = offset[n];
1915 }
1916 else
1917 {
1918 nnz_h_lag = 0;
1919 }
1920
1921 index_style = TNLP::C_STYLE;
1922
1923 return true;
1924}
1925
1926/** Method to return the bounds for my problem */
1927bool ScipNLP::get_bounds_info(
1928 Index n, /**< number of variables */
1929 Number* x_l, /**< buffer to store lower bounds on variables */
1930 Number* x_u, /**< buffer to store upper bounds on variables */
1931 Index m, /**< number of constraints */
1932 Number* g_l, /**< buffer to store lower bounds on constraints */
1933 Number* g_u /**< buffer to store lower bounds on constraints */
1934 )
1935{
1936 const int* varlincounts;
1937 const int* varnlcounts;
1938
1939 assert(nlpiproblem != NULL);
1940 assert(nlpiproblem->oracle != NULL);
1941
1942 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
1943 assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
1944
1945 assert(SCIPnlpiOracleGetVarLbs(nlpiproblem->oracle) != NULL || n == 0);
1946 assert(SCIPnlpiOracleGetVarUbs(nlpiproblem->oracle) != NULL || n == 0);
1947
1948 BMScopyMemoryArray(x_l, SCIPnlpiOracleGetVarLbs(nlpiproblem->oracle), n);
1949 BMScopyMemoryArray(x_u, SCIPnlpiOracleGetVarUbs(nlpiproblem->oracle), n);
1950#ifndef NDEBUG
1951 for( int i = 0; i < n; ++i )
1952 assert(x_l[i] <= x_u[i]);
1953#endif
1954
1955 SCIPnlpiOracleGetVarCounts(scip, nlpiproblem->oracle, &varlincounts, &varnlcounts);
1956
1957 /* Ipopt performs better when unused variables do not appear, which we can achieve by fixing them,
1958 * since Ipopts TNLPAdapter will hide them from Ipopts NLP. In the dual solution, bound multipliers (z_L, z_U)
1959 * for these variables should have value 0.0 (they are set to -grad Lagrangian).
1960 */
1961 for( int i = 0; i < n; ++i )
1962 {
1963 if( varlincounts[i] == 0 && varnlcounts[i] == 0 )
1964 {
1965 SCIPdebugMsg(scip, "fix unused variable x%d [%g,%g] to 0.0 or bound\n", i, x_l[i], x_u[i]);
1966 assert(x_l[i] <= x_u[i]);
1967 x_l[i] = x_u[i] = MAX(MIN(x_u[i], 0.0), x_l[i]);
1968 }
1969 }
1970
1971 for( int i = 0; i < m; ++i )
1972 {
1973 g_l[i] = SCIPnlpiOracleGetConstraintLhs(nlpiproblem->oracle, i);
1974 g_u[i] = SCIPnlpiOracleGetConstraintRhs(nlpiproblem->oracle, i);
1975 assert(g_l[i] <= g_u[i]);
1976 }
1977
1978 return true;
1979}
1980
1981/** Method to return the starting point for the algorithm */ /*lint -e{715}*/
1982bool ScipNLP::get_starting_point(
1983 Index n, /**< number of variables */
1984 bool init_x, /**< whether initial values for primal values are requested */
1985 Number* x, /**< buffer to store initial primal values */
1986 bool init_z, /**< whether initial values for dual values of variable bounds are requested */
1987 Number* z_L, /**< buffer to store dual values for variable lower bounds */
1988 Number* z_U, /**< buffer to store dual values for variable upper bounds */
1989 Index m, /**< number of constraints */
1990 bool init_lambda, /**< whether initial values for dual values of constraints are required */
1991 Number* lambda /**< buffer to store dual values of constraints */
1992 )
1993{ /*lint --e{715} */
1994 assert(nlpiproblem != NULL);
1995 assert(nlpiproblem->oracle != NULL);
1996
1997 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
1998 assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
1999
2000 if( init_x )
2001 {
2002 assert(nlpiproblem->solprimalvalid);
2003 assert(nlpiproblem->solprimals != NULL);
2004 BMScopyMemoryArray(x, nlpiproblem->solprimals, n);
2005 }
2006
2007 if( init_z )
2008 {
2009 assert(nlpiproblem->soldualvalid);
2010 assert(nlpiproblem->soldualvarlb != NULL);
2011 assert(nlpiproblem->soldualvarub != NULL);
2012 BMScopyMemoryArray(z_L, nlpiproblem->soldualvarlb, n);
2013 BMScopyMemoryArray(z_U, nlpiproblem->soldualvarub, n);
2014 }
2015
2016 if( init_lambda )
2017 {
2018 assert(nlpiproblem->soldualvalid);
2019 assert(nlpiproblem->soldualcons != NULL);
2020 BMScopyMemoryArray(lambda, nlpiproblem->soldualcons, m);
2021 }
2022
2023 return true;
2024}
2025
2026/** Method to return the number of nonlinear variables. */
2027Index ScipNLP::get_number_of_nonlinear_variables()
2028{
2029 int count;
2030 int n;
2031
2032 assert(nlpiproblem != NULL);
2033 assert(nlpiproblem->oracle != NULL);
2034
2035 n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
2036
2037 count = 0;
2038 for( int i = 0; i < n; ++i )
2039 if( SCIPnlpiOracleIsVarNonlinear(scip, nlpiproblem->oracle, i) )
2040 ++count;
2041
2042 return count;
2043}
2044
2045/** Method to return the indices of the nonlinear variables */
2046bool ScipNLP::get_list_of_nonlinear_variables(
2047 Index num_nonlin_vars, /**< number of nonlinear variables */
2048 Index* pos_nonlin_vars /**< array to fill with indices of nonlinear variables */
2049 )
2050{
2051 int count;
2052 int n;
2053
2054 assert(nlpiproblem != NULL);
2055 assert(nlpiproblem->oracle != NULL);
2056
2057 n = SCIPnlpiOracleGetNVars(nlpiproblem->oracle);
2058
2059 count = 0;
2060 for( int i = 0; i < n; ++i )
2061 {
2062 if( SCIPnlpiOracleIsVarNonlinear(scip, nlpiproblem->oracle, i) )
2063 {
2064 assert(count < num_nonlin_vars);
2065 pos_nonlin_vars[count++] = i;
2066 }
2067 }
2068
2069 assert(count == num_nonlin_vars);
2070
2071 return true;
2072}
2073
2074/** Method to return metadata about variables and constraints */ /*lint -e{715}*/
2075bool ScipNLP::get_var_con_metadata(
2076 Index n, /**< number of variables */
2077 StringMetaDataMapType& var_string_md, /**< variable meta data of string type */
2078 IntegerMetaDataMapType& var_integer_md, /**< variable meta data of integer type */
2079 NumericMetaDataMapType& var_numeric_md, /**< variable meta data of numeric type */
2080 Index m, /**< number of constraints */
2081 StringMetaDataMapType& con_string_md, /**< constraint meta data of string type */
2082 IntegerMetaDataMapType& con_integer_md, /**< constraint meta data of integer type */
2083 NumericMetaDataMapType& con_numeric_md /**< constraint meta data of numeric type */
2084 )
2085{ /*lint --e{715}*/
2086 assert(nlpiproblem != NULL);
2087 assert(nlpiproblem->oracle != NULL);
2088 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2089 assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2090
2091 char** varnames = SCIPnlpiOracleGetVarNames(nlpiproblem->oracle);
2092 if( varnames != NULL )
2093 {
2094 std::vector<std::string>& varnamesvec(var_string_md["idx_names"]);
2095 varnamesvec.reserve((size_t)n);
2096 for( int i = 0; i < n; ++i )
2097 {
2098 if( varnames[i] != NULL )
2099 {
2100 varnamesvec.push_back(varnames[i]); /*lint !e3701*/
2101 }
2102 else
2103 {
2104 char buffer[20];
2105 (void) SCIPsnprintf(buffer, 20, "nlpivar%8d", i);
2106 varnamesvec.push_back(buffer);
2107 }
2108 }
2109 }
2110
2111 std::vector<std::string>& consnamesvec(con_string_md["idx_names"]);
2112 consnamesvec.reserve((size_t)m);
2113 for( int i = 0; i < m; ++i )
2114 {
2115 if( SCIPnlpiOracleGetConstraintName(nlpiproblem->oracle, i) != NULL )
2116 {
2117 consnamesvec.push_back(SCIPnlpiOracleGetConstraintName(nlpiproblem->oracle, i));
2118 }
2119 else
2120 {
2121 char buffer[20];
2122 (void) SCIPsnprintf(buffer, 20, "nlpicons%8d", i);
2123 consnamesvec.push_back(buffer); /*lint !e3701*/
2124 }
2125 }
2126
2127 return true;
2128}
2129
2130/** Method to return the objective value */ /*lint -e{715}*/
2131bool ScipNLP::eval_f(
2132 Index n, /**< number of variables */
2133 const Number* x, /**< point to evaluate */
2134 bool new_x, /**< whether some function evaluation method has been called for this point before */
2135 Number& obj_value /**< place to store objective function value */
2136 )
2137{ /*lint --e{715}*/
2138 assert(nlpiproblem != NULL);
2139 assert(nlpiproblem->oracle != NULL);
2140
2141 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2142
2143 if( new_x )
2144 ++current_x;
2145 last_f_eval_x = current_x;
2146
2147 return SCIPnlpiOracleEvalObjectiveValue(scip, nlpiproblem->oracle, x, &obj_value) == SCIP_OKAY;
2148}
2149
2150/** Method to return the gradient of the objective */ /*lint -e{715}*/
2151bool ScipNLP::eval_grad_f(
2152 Index n, /**< number of variables */
2153 const Number* x, /**< point to evaluate */
2154 bool new_x, /**< whether some function evaluation method has been called for this point before */
2155 Number* grad_f /**< buffer to store objective gradient */
2156 )
2157{ /*lint --e{715}*/
2158 SCIP_Real dummy;
2159
2160 assert(nlpiproblem != NULL);
2161 assert(nlpiproblem->oracle != NULL);
2162
2163 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2164
2165 if( new_x )
2166 ++current_x;
2167 else
2168 {
2169 // pass new_x = TRUE to objective gradient eval iff we have not evaluated the objective function at this point yet
2170 new_x = last_f_eval_x < current_x;
2171 }
2172 // if we evaluate the objective gradient with new_x = true, then this will also evaluate the objective function
2173 // (and if we do with new_x = false, then we already have last_f_eval_x == current_x anyway)
2174 last_f_eval_x = current_x;
2175
2176 return SCIPnlpiOracleEvalObjectiveGradient(scip, nlpiproblem->oracle, x, new_x, &dummy, grad_f) == SCIP_OKAY;
2177}
2178
2179/** Method to return the constraint residuals */ /*lint -e{715}*/
2180bool ScipNLP::eval_g(
2181 Index n, /**< number of variables */
2182 const Number* x, /**< point to evaluate */
2183 bool new_x, /**< whether some function evaluation method has been called for this point before */
2184 Index m, /**< number of constraints */
2185 Number* g /**< buffer to store constraint function values */
2186 )
2187{ /*lint --e{715}*/
2188 assert(nlpiproblem != NULL);
2189 assert(nlpiproblem->oracle != NULL);
2190
2191 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2192
2193 if( new_x )
2194 ++current_x;
2195 last_g_eval_x = current_x;
2196
2197 return SCIPnlpiOracleEvalConstraintValues(scip, nlpiproblem->oracle, x, g) == SCIP_OKAY;
2198}
2199
2200/** Method to return:
2201 * 1) The structure of the jacobian (if "values" is NULL)
2202 * 2) The values of the jacobian (if "values" is not NULL)
2203 */ /*lint -e{715}*/
2204bool ScipNLP::eval_jac_g(
2205 Index n, /**< number of variables */
2206 const Number* x, /**< point to evaluate */
2207 bool new_x, /**< whether some function evaluation method has been called for this point before */
2208 Index m, /**< number of constraints */
2209 Index nele_jac, /**< number of nonzero entries in jacobian */
2210 Index* iRow, /**< buffer to store row indices of nonzero jacobian entries, or NULL if values are requested */
2211 Index* jCol, /**< buffer to store column indices of nonzero jacobian entries, or NULL if values are requested */
2212 Number* values /**< buffer to store values of nonzero jacobian entries, or NULL if structure is requested */
2213 )
2214{ /*lint --e{715}*/
2215 assert(nlpiproblem != NULL);
2216 assert(nlpiproblem->oracle != NULL);
2217
2218 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2219 assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2220
2221 if( values == NULL )
2222 { /* Ipopt wants to know sparsity structure */
2223 const int* jacoffset;
2224 const int* jaccol;
2225 int j;
2226 int i;
2227 int nnlnz;
2228
2229 assert(iRow != NULL);
2230 assert(jCol != NULL);
2231
2232 if( SCIPnlpiOracleGetJacobianRowSparsity(scip, nlpiproblem->oracle, &jacoffset, &jaccol, NULL, &nnlnz) != SCIP_OKAY )
2233 return false;
2234
2235 if( jacoffset != NULL )
2236 {
2237 assert(jacoffset[0] == 0);
2238 assert(jacoffset[m] == nele_jac);
2239 j = jacoffset[0];
2240 for( i = 0; i < m; ++i )
2241 for( ; j < jacoffset[i+1]; ++j )
2242 iRow[j] = i;
2243
2244 BMScopyMemoryArray(jCol, jaccol, nele_jac);
2245 }
2246 }
2247 else
2248 {
2249 if( new_x )
2250 ++current_x;
2251 else
2252 {
2253 // pass new_x = TRUE to Jacobian eval iff we have not evaluated the constraint functions at this point yet
2254 new_x = last_g_eval_x < current_x;
2255 }
2256 // if we evaluate the Jacobian with new_x = true, then this will also evaluate the constraint functions
2257 // (and if we do with new_x = false, then we already have last_g_eval_x == current_x anyway)
2258 last_f_eval_x = current_x;
2259
2260 if( SCIPnlpiOracleEvalJacobian(scip, nlpiproblem->oracle, x, new_x, NULL, values) != SCIP_OKAY )
2261 return false;
2262 }
2263
2264 return true;
2265}
2266
2267/** Method to return:
2268 * 1) The structure of the hessian of the lagrangian (if "values" is NULL)
2269 * 2) The values of the hessian of the lagrangian (if "values" is not NULL)
2270 */ /*lint -e{715}*/
2271bool ScipNLP::eval_h(
2272 Index n, /**< number of variables */
2273 const Number* x, /**< point to evaluate */
2274 bool new_x, /**< whether some function evaluation method has been called for this point before */
2275 Number obj_factor, /**< weight for objective function */
2276 Index m, /**< number of constraints */
2277 const Number* lambda, /**< weights for constraint functions */
2278 bool new_lambda, /**< whether the hessian has been evaluated for these values of lambda before */
2279 Index nele_hess, /**< number of nonzero entries in hessian */
2280 Index* iRow, /**< buffer to store row indices of nonzero hessian entries, or NULL if values are requested */
2281 Index* jCol, /**< buffer to store column indices of nonzero hessian entries, or NULL if values are requested */
2282 Number* values /**< buffer to store values of nonzero hessian entries, or NULL if structure is requested */
2283 )
2284{ /*lint --e{715}*/
2285 assert(nlpiproblem != NULL);
2286 assert(nlpiproblem->oracle != NULL);
2287
2288 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2289 assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2290
2291 if( values == NULL )
2292 { /* Ipopt wants to know sparsity structure */
2293 const int* heslagoffset;
2294 const int* heslagcol;
2295 int j;
2296 int i;
2297
2298 assert(iRow != NULL);
2299 assert(jCol != NULL);
2300
2301 if( SCIPnlpiOracleGetHessianLagSparsity(scip, nlpiproblem->oracle, &heslagoffset, &heslagcol, FALSE) != SCIP_OKAY )
2302 return false;
2303
2304 assert(heslagoffset[0] == 0);
2305 assert(heslagoffset[n] == nele_hess);
2306 j = heslagoffset[0];
2307 for( i = 0; i < n; ++i )
2308 for( ; j < heslagoffset[i+1]; ++j )
2309 iRow[j] = i;
2310
2311 BMScopyMemoryArray(jCol, heslagcol, nele_hess);
2312 }
2313 else
2314 {
2315 bool new_x_obj = new_x;
2316 bool new_x_cons = new_x;
2317 if( new_x )
2318 ++current_x;
2319 else
2320 {
2321 // pass new_x_obj = TRUE iff we have not evaluated the objective function at this point yet
2322 // pass new_x_cons = TRUE iff we have not evaluated the constraint functions at this point yet
2323 new_x_obj = last_f_eval_x < current_x;
2324 new_x_cons = last_g_eval_x < current_x;
2325 }
2326 // evaluating Hessians with new_x will also evaluate the functions itself
2327 last_f_eval_x = current_x;
2328 last_g_eval_x = current_x;
2329
2330 if( SCIPnlpiOracleEvalHessianLag(scip, nlpiproblem->oracle, x, new_x_obj, new_x_cons, obj_factor, lambda, values,
2331 FALSE) != SCIP_OKAY )
2332 return false;
2333 }
2334
2335 return true;
2336}
2337
2338/** Method called by the solver at each iteration.
2339 *
2340 * Checks whether SCIP solve is interrupted, objlimit is reached, or fastfail is triggered.
2341 * Sets solution and termination status accordingly.
2342 */ /*lint -e{715}*/
2343bool ScipNLP::intermediate_callback(
2344 AlgorithmMode mode, /**< current mode of algorithm */
2345 Index iter, /**< current iteration number */
2346 Number obj_value, /**< current objective value */
2347 Number inf_pr, /**< current primal infeasibility */
2348 Number inf_du, /**< current dual infeasibility */
2349 Number mu, /**< current barrier parameter */
2350 Number d_norm, /**< current gradient norm */
2351 Number regularization_size,/**< current size of regularization */
2352 Number alpha_du, /**< current dual alpha */
2353 Number alpha_pr, /**< current primal alpha */
2354 Index ls_trials, /**< current number of linesearch trials */
2355 const IpoptData* ip_data, /**< pointer to Ipopt Data */
2356 IpoptCalculatedQuantities* ip_cq /**< pointer to current calculated quantities */
2357 )
2358{ /*lint --e{715}*/
2360 {
2361 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2362 nlpiproblem->termstat = SCIP_NLPTERMSTAT_INTERRUPT;
2363 return false;
2364 }
2365
2366 /* feasible point with objective value below lower objective limit -> stop */
2367 if( obj_value <= param.lobjlimit && inf_pr <= param.feastol )
2368 {
2369 nlpiproblem->solstat = SCIP_NLPSOLSTAT_FEASIBLE;
2370 nlpiproblem->termstat = SCIP_NLPTERMSTAT_LOBJLIMIT;
2371 return false;
2372 }
2373
2374 /* do convergence test if fastfail is enabled */
2375 if( param.fastfail >= SCIP_NLPPARAM_FASTFAIL_AGGRESSIVE )
2376 {
2377 int i;
2378
2379 if( iter == 0 )
2380 {
2381 conv_lastrestoiter = -1;
2382 }
2383 else if( mode == RestorationPhaseMode )
2384 {
2385 conv_lastrestoiter = iter;
2386 }
2387 else if( conv_lastrestoiter == iter-1 )
2388 {
2389 /* just switched back from restoration mode, reset dual reduction targets */
2390 for( i = 0; i < convcheck_nchecks; ++i )
2391 conv_dutarget[i] = convcheck_minred[i] * inf_du;
2392 }
2393
2394 if( iter == convcheck_startiter )
2395 {
2396 /* define initial targets and iteration limits */
2397 for( i = 0; i < convcheck_nchecks; ++i )
2398 {
2399 conv_prtarget[i] = convcheck_minred[i] * inf_pr;
2400 conv_dutarget[i] = convcheck_minred[i] * inf_du;
2401 conv_iterlim[i] = iter + convcheck_maxiter[i];
2402 }
2403 }
2404 else if( iter > convcheck_startiter )
2405 {
2406 /* check if we should stop */
2407 for( i = 0; i < convcheck_nchecks; ++i )
2408 {
2409 if( inf_pr <= conv_prtarget[i] )
2410 {
2411 /* sufficient reduction w.r.t. primal infeasibility target
2412 * reset target w.r.t. current infeasibilities
2413 */
2414 conv_prtarget[i] = convcheck_minred[i] * inf_pr;
2415 conv_dutarget[i] = convcheck_minred[i] * inf_du;
2416 conv_iterlim[i] = iter + convcheck_maxiter[i];
2417 }
2418 else if( iter >= conv_iterlim[i] )
2419 {
2420 /* we hit a limit, should we really stop? */
2421 SCIPdebugMsg(scip, "convcheck %d: inf_pr = %e > target %e; inf_du = %e target %e: ",
2422 i, inf_pr, conv_prtarget[i], inf_du, conv_dutarget[i]);
2423 if( mode == RegularMode && iter <= conv_lastrestoiter + convcheck_startiter )
2424 {
2425 /* if we returned from feasibility restoration recently, we allow some more iterations,
2426 * because Ipopt may go for optimality for some iterations, at the costs of infeasibility
2427 */
2428 SCIPdebugPrintf("continue, because restoration phase only %d iters ago\n", iter - conv_lastrestoiter);
2429 }
2430 else if( mode == RegularMode && inf_du <= conv_dutarget[i] && iter < conv_iterlim[i] + convcheck_maxiter[i] )
2431 {
2432 /* if dual reduction is sufficient, we allow for twice the number of iterations to reach primal infeas reduction */
2433 SCIPdebugPrintf("continue, because dual infeas. red. sufficient and only %d iters above limit\n", iter - conv_iterlim[i]);
2434 }
2435 else
2436 {
2437 SCIPdebugPrintf("abort solve\n");
2438 if( inf_pr <= param.feastol )
2439 nlpiproblem->solstat = SCIP_NLPSOLSTAT_FEASIBLE;
2440 else
2441 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2442 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OKAY;
2443 return false;
2444 }
2445 }
2446 }
2447 }
2448 }
2449
2450 return true;
2451}
2452
2453/** This method is called when the algorithm is complete so the TNLP can store/write the solution. */ /*lint -e{715}*/
2454void ScipNLP::finalize_solution(
2455 SolverReturn status, /**< solve and solution status */
2456 Index n, /**< number of variables */
2457 const Number* x, /**< primal solution values */
2458 const Number* z_L, /**< dual values of variable lower bounds */
2459 const Number* z_U, /**< dual values of variable upper bounds */
2460 Index m, /**< number of constraints */
2461 const Number* g, /**< values of constraints */
2462 const Number* lambda, /**< dual values of constraints */
2463 Number obj_value, /**< objective function value */
2464 const IpoptData* data, /**< pointer to Ipopt Data */
2465 IpoptCalculatedQuantities* cq /**< pointer to calculated quantities */
2466 )
2467{ /*lint --e{715}*/
2468 assert(nlpiproblem != NULL);
2469 assert(nlpiproblem->oracle != NULL);
2470
2471 assert(n == SCIPnlpiOracleGetNVars(nlpiproblem->oracle));
2472 assert(m == SCIPnlpiOracleGetNConstraints(nlpiproblem->oracle));
2473
2474 bool check_feasibility = false; // whether we should check x for feasibility, if not NULL
2475 ++(nlpiproblem->ncalls);
2476 switch( status )
2477 {
2478 case SUCCESS:
2479 nlpiproblem->solstat = SCIP_NLPSOLSTAT_LOCOPT;
2480 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OKAY;
2481 ++(nlpiproblem->nsuccess);
2482 assert(x != NULL);
2483 break;
2484
2485 case STOP_AT_ACCEPTABLE_POINT:
2486 /* if stop at acceptable point, then dual infeasibility can be arbitrary large, so claim only feasibility */
2487 case FEASIBLE_POINT_FOUND:
2488 nlpiproblem->solstat = SCIP_NLPSOLSTAT_FEASIBLE;
2489 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OKAY;
2490 ++(nlpiproblem->nsuccess);
2491 assert(x != NULL);
2492 break;
2493
2494 case MAXITER_EXCEEDED:
2495 check_feasibility = true;
2496 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2497 nlpiproblem->termstat = SCIP_NLPTERMSTAT_ITERLIMIT;
2498 ++(nlpiproblem->nlimit);
2499 break;
2500
2501 case CPUTIME_EXCEEDED:
2502#if IPOPT_VERSION_MAJOR > 3 || IPOPT_VERSION_MINOR >= 14
2503 case WALLTIME_EXCEEDED:
2504#endif
2505 check_feasibility = true;
2506 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2507 nlpiproblem->termstat = SCIP_NLPTERMSTAT_TIMELIMIT;
2508 ++(nlpiproblem->nlimit);
2509 break;
2510
2511 case STOP_AT_TINY_STEP:
2512 case RESTORATION_FAILURE:
2513 case ERROR_IN_STEP_COMPUTATION:
2514 check_feasibility = true;
2515 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2516 nlpiproblem->termstat = SCIP_NLPTERMSTAT_NUMERICERROR;
2517 ++(nlpiproblem->nother);
2518 break;
2519
2520 case LOCAL_INFEASIBILITY:
2521 nlpiproblem->solstat = SCIP_NLPSOLSTAT_LOCINFEASIBLE;
2522 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OKAY;
2523 ++(nlpiproblem->nlocinfeas);
2524 break;
2525
2526 case USER_REQUESTED_STOP:
2527 // status codes already set in intermediate_callback
2528 break;
2529
2530 case DIVERGING_ITERATES:
2531 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNBOUNDED;
2532 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OKAY;
2533 ++(nlpiproblem->nsuccess);
2534 break;
2535
2536 // for the following status codes, if we get called here at all,
2537 // then Ipopt passes zeros for duals and activities!
2538 // (see https://github.com/coin-or/Ipopt/blob/stable/3.14/src/Interfaces/IpIpoptApplication.cpp#L885-L934)
2539
2540 case INVALID_NUMBER_DETECTED:
2541 // we can get this, if functions can still be evaluated, but are not differentiable
2542 // (so Ipopt couldn't check local optimality)
2543 // so we enable the check below for whether the point is feasible
2544 check_feasibility = true;
2545 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2546 nlpiproblem->termstat = SCIP_NLPTERMSTAT_EVALERROR;
2547 ++(nlpiproblem->nother);
2548 break;
2549
2550 case TOO_FEW_DEGREES_OF_FREEDOM:
2551 case INTERNAL_ERROR:
2552 case INVALID_OPTION:
2553 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2554 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OTHER;
2555 ++(nlpiproblem->nother);
2556 break;
2557
2558 case OUT_OF_MEMORY:
2559 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2560 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OUTOFMEMORY;
2561 ++(nlpiproblem->nother);
2562 break;
2563
2564 default:
2565 SCIPerrorMessage("Ipopt returned with unknown solution status %d\n", status);
2566 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2567 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OTHER;
2568 ++(nlpiproblem->nother);
2569 break;
2570 }
2571
2572 assert(x != NULL);
2573 assert(lambda != NULL);
2574 assert(z_L != NULL);
2575 assert(z_U != NULL);
2576
2577 assert(nlpiproblem->solprimals != NULL);
2578
2579 if( nlpiproblem->soldualcons == NULL )
2580 {
2581 (void) SCIPallocBlockMemoryArray(scip, &nlpiproblem->soldualcons, m);
2582 }
2583 if( nlpiproblem->soldualvarlb == NULL )
2584 {
2585 (void) SCIPallocBlockMemoryArray(scip, &nlpiproblem->soldualvarlb, n);
2586 }
2587 if( nlpiproblem->soldualvarub == NULL )
2588 {
2589 (void) SCIPallocBlockMemoryArray(scip, &nlpiproblem->soldualvarub, n);
2590 }
2591 if( nlpiproblem->soldualcons == NULL || nlpiproblem->soldualvarlb == NULL || nlpiproblem->soldualvarub == NULL )
2592 {
2593 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2594 nlpiproblem->termstat = SCIP_NLPTERMSTAT_OUTOFMEMORY;
2595 return;
2596 }
2597
2598 BMScopyMemoryArray(nlpiproblem->solprimals, x, n);
2599 BMScopyMemoryArray(nlpiproblem->soldualcons, lambda, m);
2600 BMScopyMemoryArray(nlpiproblem->soldualvarlb, z_L, n);
2601 BMScopyMemoryArray(nlpiproblem->soldualvarub, z_U, n);
2602 nlpiproblem->solobjval = obj_value;
2603 nlpiproblem->solprimalvalid = true;
2604 nlpiproblem->solprimalgiven = false;
2605 nlpiproblem->soldualvalid = true;
2606 nlpiproblem->soldualgiven = false;
2607
2608 // get violations, there could be an evaluation error when doing so
2609#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
2610 nlpiproblem->solboundviol = 0.0; // old Ipopt does not calculate bound violations, but for what it's worth, we have set bound_relax_factor=0 then
2611 if( cq == NULL )
2612 {
2613 // with old Ipopt, finalize_solution may be called with cq == NULL if all variables are fixed; we just skip the rest then
2614 nlpiproblem->solconsviol = 0.0;
2615 return;
2616 }
2617#else
2618 assert(cq != NULL);
2619 nlpiproblem->solboundviol = cq->unscaled_curr_orig_bounds_violation(Ipopt::NORM_MAX);
2620#endif
2621 try
2622 {
2623 nlpiproblem->solconsviol = cq->unscaled_curr_nlp_constraint_violation(Ipopt::NORM_MAX);
2624
2625 if( check_feasibility )
2626 {
2627 // we assume that check_feasibility has not been enabled if Ipopt claimed infeasibility, since we should not change solstatus to unknown then
2628 assert(nlpiproblem->solstat != SCIP_NLPSOLSTAT_LOCINFEASIBLE);
2629 if( MAX(nlpiproblem->solconsviol, nlpiproblem->solboundviol) <= param.feastol )
2630 nlpiproblem->solstat = SCIP_NLPSOLSTAT_FEASIBLE;
2631 else
2632 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2633 }
2634 }
2635 catch( const IpoptNLP::Eval_Error& exc )
2636 {
2637 SCIPdebugMsg(scip, "Eval error when checking constraint viol: %s\n", exc.Message().c_str());
2638 assert(status == INVALID_NUMBER_DETECTED);
2639 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2640 nlpiproblem->solconsviol = SCIP_INVALID;
2641 }
2642 catch(...)
2643 {
2644 /* with clang++, an IpoptNLP::Eval_Error wasn't catched by the catch-block above for Ipopt < 3.14.20
2645 * if visibility attributes of the Ipopt exception classes in the Ipopt and SCIP libs differed,
2646 * we we also catch any exception here
2647 */
2648 SCIPdebugMsg(scip, "Unknown exception when checking constraint viol\n");
2649 assert(status == INVALID_NUMBER_DETECTED);
2650 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2651 nlpiproblem->solconsviol = SCIP_INVALID;
2652 }
2653
2654 if( nlpiproblem->solstat == SCIP_NLPSOLSTAT_LOCINFEASIBLE )
2655 {
2656 assert(lambda != NULL);
2657 SCIP_Real tol;
2658 (void) nlpiproblem->ipopt->Options()->GetNumericValue("tol", tol, "");
2659
2660 // Jakobs paper ZR_20-20 says we should have lambda*g(x) + mu*h(x) > 0
2661 // if the NLP is min f(x) s.t. g(x) <= 0, h(x) = 0
2662 // we check this here and change solution status to unknown if the test fails
2663 bool infreasonable = true;
2664 SCIP_Real infproof = 0.0;
2665 for( int i = 0; i < m && infreasonable; ++i )
2666 {
2667 if( fabs(lambda[i]) < tol )
2668 continue;
2669 SCIP_Real side;
2670 if( lambda[i] < 0.0 )
2671 {
2672 // lhs <= g(x) should be active
2673 // in the NLP above, this should be lhs - g(x) <= 0 with negated dual
2674 // so this contributes -lambda*(lhs-g(x)) = lambda*(g(x)-side)
2675 side = SCIPnlpiOracleGetConstraintLhs(nlpiproblem->oracle, i);
2676 if( SCIPisInfinity(scip, -side) )
2677 {
2678 SCIPdebugMessage("inconsistent dual, lambda = %g, but lhs = %g\n", lambda[i], side);
2679 infreasonable = false;
2680 }
2681 }
2682 else
2683 {
2684 // g(x) <= rhs should be active
2685 // in the NLP above, this should be g(x) - rhs <= 0
2686 // so this contributes lambda*(g(x)-rhs)
2687 side = SCIPnlpiOracleGetConstraintRhs(nlpiproblem->oracle, i);
2688 if( SCIPisInfinity(scip, side) )
2689 {
2690 SCIPdebugMessage("inconsistent dual, lambda = %g, but rhs = %g\n", lambda[i], side);
2691 infreasonable = false;
2692 }
2693 }
2694
2695 // g(x) <= 0
2696 infproof += lambda[i] * (g[i] - side);
2697 // SCIPdebugMessage("cons %d lambda %g, slack %g\n", i, lambda[i], g[i] - side);
2698 }
2699 if( infreasonable )
2700 {
2701 SCIPdebugMessage("infproof = %g should be positive to be valid\n", infproof);
2702 if( infproof <= 0.0 )
2703 infreasonable = false;
2704 }
2705
2706 if( !infreasonable )
2707 {
2708 // change status to say we don't know
2709 nlpiproblem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
2710 }
2711 }
2712}
2713
2714/** Calls Lapacks Dsyev routine to compute eigenvalues and eigenvectors of a dense matrix.
2715 *
2716 * This uses Ipopt's interface to Lapack.
2717 */
2719 SCIP_Bool computeeigenvectors,/**< should also eigenvectors should be computed ? */
2720 int N, /**< dimension */
2721 SCIP_Real* a, /**< matrix data on input (size N*N); eigenvectors on output if computeeigenvectors == TRUE */
2722 SCIP_Real* w /**< buffer to store eigenvalues (size N) */
2723 )
2724{
2725 int info;
2726
2727#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
2728 IpLapackDsyev((bool)computeeigenvectors, N, a, N, w, info);
2729#else
2730 IpLapackSyev((bool)computeeigenvectors, N, a, N, w, info);
2731#endif
2732
2733 if( info != 0 )
2734 {
2735 SCIPerrorMessage("There was an error when calling DSYEV. INFO = %d\n", info);
2736 return SCIP_ERROR;
2737 }
2738
2739 return SCIP_OKAY;
2740}
2741
2742/** solves a linear problem of the form Ax = b for a regular matrix 3*3 A */
2743static
2745 SCIP_Real* A, /**< matrix data on input (size 3*3); filled column-wise */
2746 SCIP_Real* b, /**< right hand side vector (size 3) */
2747 SCIP_Real* x, /**< buffer to store solution (size 3) */
2748 SCIP_Bool* success /**< pointer to store if the solving routine was successful */
2749 )
2750{
2751 SCIP_Real Acopy[9];
2752 SCIP_Real bcopy[3];
2753 int pivotcopy[3];
2754 const int N = 3;
2755 int info;
2756
2757 assert(A != NULL);
2758 assert(b != NULL);
2759 assert(x != NULL);
2760 assert(success != NULL);
2761
2762 BMScopyMemoryArray(Acopy, A, N*N);
2763 BMScopyMemoryArray(bcopy, b, N);
2764
2765 /* compute the LU factorization */
2766#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
2767 IpLapackDgetrf(N, Acopy, pivotcopy, N, info);
2768#else
2769 IpLapackGetrf(N, Acopy, pivotcopy, N, info);
2770#endif
2771
2772 if( info != 0 )
2773 {
2774 SCIPdebugMessage("There was an error when calling Dgetrf. INFO = %d\n", info);
2775 *success = FALSE;
2776 }
2777 else
2778 {
2779 *success = TRUE;
2780
2781 /* solve linear problem */
2782#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
2783 IpLapackDgetrs(N, 1, Acopy, N, pivotcopy, bcopy, N);
2784#else
2785 IpLapackGetrs(N, 1, Acopy, N, pivotcopy, bcopy, N);
2786#endif
2787
2788 /* copy the solution */
2789 BMScopyMemoryArray(x, bcopy, N);
2790 }
2791
2792 return SCIP_OKAY;
2793}
2794
2795/** solves a linear problem of the form Ax = b for a regular matrix A
2796 *
2797 * Calls Lapacks DGETRF routine to calculate a LU factorization and uses this factorization to solve
2798 * the linear problem Ax = b.
2799 *
2800 * This uses Ipopt's interface to Lapack.
2801 */
2803 int N, /**< dimension */
2804 SCIP_Real* A, /**< matrix data on input (size N*N); filled column-wise */
2805 SCIP_Real* b, /**< right hand side vector (size N) */
2806 SCIP_Real* x, /**< buffer to store solution (size N) */
2807 SCIP_Bool* success /**< pointer to store if the solving routine was successful */
2808 )
2809{
2810 SCIP_Real* Acopy;
2811 SCIP_Real* bcopy;
2812 int* pivotcopy;
2813 int info;
2814
2815 assert(N > 0);
2816 assert(A != NULL);
2817 assert(b != NULL);
2818 assert(x != NULL);
2819 assert(success != NULL);
2820
2821 /* call solveLinearProb3() for performance reasons */
2822 if( N == 3 )
2823 {
2824 SCIP_CALL( solveLinearProb3(A, b, x, success) );
2825 return SCIP_OKAY;
2826 }
2827
2828 Acopy = NULL;
2829 bcopy = NULL;
2830 pivotcopy = NULL;
2831
2832 SCIP_ALLOC( BMSduplicateMemoryArray(&Acopy, A, N*N) );
2833 SCIP_ALLOC( BMSduplicateMemoryArray(&bcopy, b, N) );
2834 SCIP_ALLOC( BMSallocMemoryArray(&pivotcopy, N) );
2835
2836 /* compute the LU factorization */
2837#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
2838 IpLapackDgetrf(N, Acopy, pivotcopy, N, info);
2839#else
2840 IpLapackGetrf(N, Acopy, pivotcopy, N, info);
2841#endif
2842
2843 if( info != 0 )
2844 {
2845 SCIPdebugMessage("There was an error when calling Dgetrf. INFO = %d\n", info);
2846 *success = FALSE;
2847 }
2848 else
2849 {
2850 *success = TRUE;
2851
2852 /* solve linear problem */
2853#if IPOPT_VERSION_MAJOR == 3 && IPOPT_VERSION_MINOR < 14
2854 IpLapackDgetrs(N, 1, Acopy, N, pivotcopy, bcopy, N);
2855#else
2856 IpLapackGetrs(N, 1, Acopy, N, pivotcopy, bcopy, N);
2857#endif
2858
2859 /* copy the solution */
2860 BMScopyMemoryArray(x, bcopy, N);
2861 }
2862
2863 BMSfreeMemoryArray(&pivotcopy);
2864 BMSfreeMemoryArray(&bcopy);
2865 BMSfreeMemoryArray(&Acopy);
2866
2867 return SCIP_OKAY;
2868}
SCIP_VAR * w
Definition: circlepacking.c:67
SCIP_VAR * a
Definition: circlepacking.c:66
SCIP_VAR ** b
Definition: circlepacking.c:65
SCIP_VAR ** x
Definition: circlepacking.c:63
#define NULL
Definition: def.h:248
#define SCIP_MAXSTRLEN
Definition: def.h:269
#define SCIP_INVALID
Definition: def.h:178
#define SCIP_Bool
Definition: def.h:91
#define MIN(x, y)
Definition: def.h:224
#define SCIP_ALLOC(x)
Definition: def.h:366
#define SCIP_Real
Definition: def.h:156
#define SQR(x)
Definition: def.h:199
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define MAX(x, y)
Definition: def.h:220
#define SCIP_CALL_ABORT(x)
Definition: def.h:334
#define SCIP_CALL(x)
Definition: def.h:355
methods to interpret (evaluate) an expression "fast"
int SCIPgetSubscipDepth(SCIP *scip)
Definition: scip_copy.c:2588
SCIP_EXPRINTCAPABILITY SCIPexprintGetCapability(void)
void SCIPinfoMessage(SCIP *scip, FILE *file, const char *formatstr,...)
Definition: scip_message.c:208
SCIP_MESSAGEHDLR * SCIPgetMessagehdlr(SCIP *scip)
Definition: scip_message.c:88
#define SCIPdebugMsg
Definition: scip_message.h:78
SCIP_RETCODE SCIPincludeNlpSolverIpopt(SCIP *scip)
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveValue(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *objval)
Definition: nlpioracle.c:2124
SCIP_RETCODE SCIPnlpiOracleChgLinearCoefs(SCIP *scip, SCIP_NLPIORACLE *oracle, int considx, int nentries, const int *varidxs, const SCIP_Real *newcoefs)
Definition: nlpioracle.c:1774
SCIP_RETCODE SCIPnlpiOracleGetJacobianRowSparsity(SCIP *scip, SCIP_NLPIORACLE *oracle, const int **rowoffsets, const int **cols, const SCIP_Bool **colnlflags, int *nnlnz)
Definition: nlpioracle.c:2264
SCIP_RETCODE SCIPnlpiOracleGetHessianLagSparsity(SCIP *scip, SCIP_NLPIORACLE *oracle, const int **offset, const int **allnz, SCIP_Bool colwise)
Definition: nlpioracle.c:2670
SCIP_RETCODE SCIPnlpiOracleChgVarBounds(SCIP *scip, SCIP_NLPIORACLE *oracle, int nvars, const int *indices, const SCIP_Real *lbs, const SCIP_Real *ubs)
Definition: nlpioracle.c:1474
SCIP_RETCODE SCIPnlpiOracleAddConstraints(SCIP *scip, SCIP_NLPIORACLE *oracle, int nconss, const SCIP_Real *lhss, const SCIP_Real *rhss, const int *nlininds, int *const *lininds, SCIP_Real *const *linvals, SCIP_EXPR **exprs, const char **consnames)
Definition: nlpioracle.c:1384
SCIP_Bool SCIPnlpiOracleIsConstraintNonlinear(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2084
SCIP_RETCODE SCIPnlpiOracleDelVarSet(SCIP *scip, SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1546
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValues(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *convals)
Definition: nlpioracle.c:2173
SCIP_RETCODE SCIPnlpiOracleCreate(SCIP *scip, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1200
void SCIPnlpiOracleGetVarCounts(SCIP *scip, SCIP_NLPIORACLE *oracle, const int **lincounts, const int **nlcounts)
Definition: nlpioracle.c:1998
char * SCIPnlpiOracleGetConstraintName(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2051
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveGradient(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *objval, SCIP_Real *objgrad)
Definition: nlpioracle.c:2205
SCIP_RETCODE SCIPnlpiOracleResetEvalTime(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2821
SCIP_RETCODE SCIPnlpiOracleSetObjective(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real constant, int nlin, const int *lininds, const SCIP_Real *linvals, SCIP_EXPR *expr)
Definition: nlpioracle.c:1445
SCIP_Real SCIPnlpiOracleGetConstraintRhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2038
SCIP_Real SCIPnlpiOracleGetEvalTime(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2837
SCIP_RETCODE SCIPnlpiOracleChgConsSides(SCIP *scip, SCIP_NLPIORACLE *oracle, int nconss, const int *indices, const SCIP_Real *lhss, const SCIP_Real *rhss)
Definition: nlpioracle.c:1511
SCIP_Real SCIPnlpiOracleGetConstraintLhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2025
SCIP_RETCODE SCIPnlpiOracleAddVars(SCIP *scip, SCIP_NLPIORACLE *oracle, int nvars, const SCIP_Real *lbs, const SCIP_Real *ubs, const char **varnames)
Definition: nlpioracle.c:1298
SCIP_RETCODE SCIPnlpiOracleEvalHessianLag(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx_obj, SCIP_Bool isnewx_cons, SCIP_Real objfactor, const SCIP_Real *lambda, SCIP_Real *hessian, SCIP_Bool colwise)
Definition: nlpioracle.c:2770
int SCIPnlpiOracleGetNVars(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1933
int SCIPnlpiOracleGetNConstraints(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1943
SCIP_EXPRINTCAPABILITY SCIPnlpiOracleGetEvalCapability(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2101
SCIP_Real SCIPnlpiOracleGetObjectiveConstant(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2014
SCIP_Bool SCIPnlpiOracleIsVarNonlinear(SCIP *scip, SCIP_NLPIORACLE *oracle, int varidx)
Definition: nlpioracle.c:1983
SCIP_RETCODE SCIPnlpiOracleEvalJacobian(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *convals, SCIP_Real *jacobi)
Definition: nlpioracle.c:2543
SCIP_RETCODE SCIPnlpiOracleDelConsSet(SCIP *scip, SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1688
SCIP_RETCODE SCIPnlpiOracleSetProblemName(SCIP *scip, SCIP_NLPIORACLE *oracle, const char *name)
Definition: nlpioracle.c:1262
SCIP_RETCODE SCIPnlpiOracleChgObjConstant(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_Real objconstant)
Definition: nlpioracle.c:1916
char ** SCIPnlpiOracleGetVarNames(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1973
const SCIP_Real * SCIPnlpiOracleGetVarLbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1953
const SCIP_Real * SCIPnlpiOracleGetVarUbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1963
SCIP_RETCODE SCIPnlpiOracleFree(SCIP *scip, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1230
const char * SCIPnlpiOracleGetProblemName(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1286
SCIP_RETCODE SCIPnlpiOracleChgExpr(SCIP *scip, SCIP_NLPIORACLE *oracle, int considx, SCIP_EXPR *expr)
Definition: nlpioracle.c:1870
const char * SCIPgetSolverNameIpopt(void)
SCIP_RETCODE SCIPcallLapackDsyevIpopt(SCIP_Bool computeeigenvectors, int N, SCIP_Real *a, SCIP_Real *w)
SCIP_RETCODE SCIPsolveLinearEquationsIpopt(int N, SCIP_Real *A, SCIP_Real *b, SCIP_Real *x, SCIP_Bool *success)
void * SCIPgetNlpiOracleIpopt(SCIP_NLPIPROBLEM *nlpiproblem)
SCIP_Bool SCIPisIpoptAvailableIpopt(void)
const char * SCIPgetSolverDescIpopt(void)
SCIP_RETCODE SCIPaddIntParam(SCIP *scip, const char *name, const char *desc, int *valueptr, SCIP_Bool isadvanced, int defaultvalue, int minvalue, int maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:83
SCIP_PARAM * SCIPgetParam(SCIP *scip, const char *name)
Definition: scip_param.c:234
SCIP_RETCODE SCIPaddStringParam(SCIP *scip, const char *name, const char *desc, char **valueptr, SCIP_Bool isadvanced, const char *defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:194
SCIP_RETCODE SCIPaddRealParam(SCIP *scip, const char *name, const char *desc, SCIP_Real *valueptr, SCIP_Bool isadvanced, SCIP_Real defaultvalue, SCIP_Real minvalue, SCIP_Real maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:139
SCIP_RETCODE SCIPincludeExternalCodeInformation(SCIP *scip, const char *name, const char *description)
Definition: scip_general.c:769
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:93
#define SCIPreallocBlockMemoryArray(scip, ptr, oldnum, newnum)
Definition: scip_mem.h:99
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:111
SCIP_RETCODE SCIPincludeNlpi(SCIP *scip, 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_NLPICHGEXPR((*nlpichgexpr)), 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_NLPIDATA *nlpidata)
Definition: scip_nlpi.c:113
SCIP_NLPIDATA * SCIPnlpiGetData(SCIP_NLPI *nlpi)
Definition: nlpi.c:712
SCIP_Bool SCIPisSolveInterrupted(SCIP *scip)
Definition: scip_solve.c:3580
SCIP_Real SCIPinfinity(SCIP *scip)
SCIP_Real SCIPfloor(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
void SCIPfreeRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen)
SCIP_Real SCIPrandomGetReal(SCIP_RANDNUMGEN *randnumgen, SCIP_Real minrandval, SCIP_Real maxrandval)
Definition: misc.c:10245
SCIP_RETCODE SCIPcreateRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen, unsigned int initialseed, SCIP_Bool useglobalseed)
int SCIPsnprintf(char *t, int len, const char *s,...)
Definition: misc.c:10827
static const char * paramname[]
Definition: lpi_msk.c:5172
#define BMSduplicateMemoryArray(ptr, source, num)
Definition: memory.h:143
#define BMSallocMemoryArray(ptr, num)
Definition: memory.h:123
#define BMSfreeMemoryArray(ptr)
Definition: memory.h:147
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:134
void SCIPmessageVPrintError(const char *formatstr, va_list ap)
Definition: message.c:804
void SCIPmessageVPrintInfo(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr, va_list ap)
Definition: message.c:608
void SCIPmessagePrintErrorHeader(const char *sourcefile, int sourceline)
Definition: message.c:777
static SCIP_DECL_NLPICHGEXPR(nlpiChgExprIpopt)
static SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemIpopt)
Definition: nlpi_ipopt.cpp:927
#define NLPI_PRIORITY
Definition: nlpi_ipopt.cpp:102
#define MAXPERTURB
Definition: nlpi_ipopt.cpp:104
static const SCIP_Real convcheck_minred[convcheck_nchecks]
Definition: nlpi_ipopt.cpp:142
static SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesIpopt)
static const char * ipopt_int_params[]
integer parameters of Ipopt to make available via SCIP parameters
Definition: nlpi_ipopt.cpp:145
static SCIP_DECL_NLPIGETSOLVERPOINTER(nlpiGetSolverPointerIpopt)
Definition: nlpi_ipopt.cpp:917
static SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetIpopt)
static SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsIpopt)
static SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsIpopt)
static SCIP_DECL_NLPIDELCONSSET(nlpiDelConstraintSetIpopt)
static SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatIpopt)
static SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatIpopt)
static SCIP_DECL_NLPICOPY(nlpiCopyIpopt)
Definition: nlpi_ipopt.cpp:896
#define NLPI_NAME
Definition: nlpi_ipopt.cpp:100
static SCIP_DECL_NLPICHGOBJCONSTANT(nlpiChgObjConstantIpopt)
static SCIP_DECL_NLPISOLVE(nlpiSolveIpopt)
static const int convcheck_startiter
Definition: nlpi_ipopt.cpp:140
static const int convcheck_nchecks
Definition: nlpi_ipopt.cpp:139
static SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessIpopt)
static SCIP_DECL_NLPIGETPROBLEMPOINTER(nlpiGetProblemPointerIpopt)
#define DEFAULT_RANDSEED
Definition: nlpi_ipopt.cpp:107
static SCIP_RETCODE solveLinearProb3(SCIP_Real *A, SCIP_Real *b, SCIP_Real *x, SCIP_Bool *success)
#define FEASTOLFACTOR
Definition: nlpi_ipopt.cpp:105
static const char * ipopt_string_params[]
string parameters of Ipopt to make available via SCIP parameters
Definition: nlpi_ipopt.cpp:149
static SCIP_DECL_NLPIFREE(nlpiFreeIpopt)
Definition: nlpi_ipopt.cpp:905
static SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsIpopt)
static SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsIpopt)
static SCIP_RETCODE ensureStartingPoint(SCIP *scip, SCIP_NLPIPROBLEM *problem, SCIP_Bool &warmstart)
Definition: nlpi_ipopt.cpp:517
static SCIP_RETCODE handleNlpParam(SCIP *scip, SCIP_NLPIDATA *nlpidata, SCIP_NLPIPROBLEM *nlpiproblem, SCIP_NLPPARAM param)
pass NLP solve parameters to Ipopt
Definition: nlpi_ipopt.cpp:587
static const int convcheck_maxiter[convcheck_nchecks]
Definition: nlpi_ipopt.cpp:141
static SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionIpopt)
static SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemIpopt)
static SCIP_DECL_NLPIADDVARS(nlpiAddVarsIpopt)
static SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveIpopt)
static void invalidateSolved(SCIP_NLPIPROBLEM *problem)
Definition: nlpi_ipopt.cpp:486
#define NLPI_DESC
Definition: nlpi_ipopt.cpp:101
static void invalidateSolution(SCIP_NLPIPROBLEM *problem)
Definition: nlpi_ipopt.cpp:501
Ipopt NLP interface.
methods to store an NLP and request function, gradient, and Hessian values
char * SCIPparamGetString(SCIP_PARAM *param)
Definition: paramset.c:910
int SCIPparamGetInt(SCIP_PARAM *param)
Definition: paramset.c:733
int SCIPparamGetIntDefault(SCIP_PARAM *param)
Definition: paramset.c:769
public methods for message output
#define SCIPerrorMessage
Definition: pub_message.h:64
#define SCIPdebugMessage
Definition: pub_message.h:96
#define SCIPdebugPrintf
Definition: pub_message.h:99
public data structures and miscellaneous methods
public methods for handling parameter settings
public methods for problem copies
general public methods
public methods for memory management
public methods for message handling
public methods for NLPI solver interfaces
public methods for numerical tolerances
public methods for SCIP parameter handling
public methods for random numbers
public solving methods
SCIP_Real timelimit
Definition: type_nlpi.h:72
SCIP_Real feastol
Definition: type_nlpi.h:69
SCIP_Bool expectinfeas
Definition: type_nlpi.h:76
SCIP_Bool warmstart
Definition: type_nlpi.h:77
SCIP_Real opttol
Definition: type_nlpi.h:70
SCIP_Real solvertol
Definition: type_nlpi.h:71
SCIP_NLPPARAM_FASTFAIL fastfail
Definition: type_nlpi.h:75
unsigned short verblevel
Definition: type_nlpi.h:74
SCIP_Real solconsviol
Definition: nlpi_ipopt.cpp:197
SCIP_Real * soldualcons
Definition: nlpi_ipopt.cpp:193
SCIP_Real solboundviol
Definition: nlpi_ipopt.cpp:198
SmartPtr< IpoptApplication > ipopt
Definition: nlpi_ipopt.cpp:180
SCIP_Bool firstrun
Definition: nlpi_conopt.c:81
SCIP_NLPIORACLE * oracle
Definition: nlpi_conopt.c:78
SCIP_Real lasttime
Definition: nlpi_ipopt.cpp:200
SCIP_RANDNUMGEN * randnumgen
Definition: nlpi_conopt.c:79
SCIP_Real * soldualvarlb
Definition: nlpi_ipopt.cpp:194
SCIP_NLPTERMSTAT termstat
Definition: nlpi_conopt.c:85
SmartPtr< ScipNLP > nlp
Definition: nlpi_ipopt.cpp:181
SCIP_Real solobjval
Definition: nlpi_ipopt.cpp:196
SCIP_NLPSOLSTAT solstat
Definition: nlpi_conopt.c:84
SCIP_Real * soldualvarub
Definition: nlpi_ipopt.cpp:195
SCIP_Real * solprimals
Definition: nlpi_ipopt.cpp:192
#define SCIP_EXPRINTCAPABILITY_GRADIENT
#define SCIP_EXPRINTCAPABILITY_HESSIAN
#define SCIP_EXPRINTCAPABILITY_FUNCVALUE
unsigned int SCIP_EXPRINTCAPABILITY
#define SCIP_NLPPARAM_PRINT(param)
Definition: type_nlpi.h:142
@ SCIP_NLPPARAM_FASTFAIL_OFF
Definition: type_nlpi.h:58
@ SCIP_NLPPARAM_FASTFAIL_AGGRESSIVE
Definition: type_nlpi.h:60
enum SCIP_NlpSolStat SCIP_NLPSOLSTAT
Definition: type_nlpi.h:168
@ SCIP_NLPTERMSTAT_OKAY
Definition: type_nlpi.h:173
@ SCIP_NLPTERMSTAT_TIMELIMIT
Definition: type_nlpi.h:174
@ SCIP_NLPTERMSTAT_NUMERICERROR
Definition: type_nlpi.h:178
@ SCIP_NLPTERMSTAT_OTHER
Definition: type_nlpi.h:182
@ SCIP_NLPTERMSTAT_EVALERROR
Definition: type_nlpi.h:179
@ SCIP_NLPTERMSTAT_LOBJLIMIT
Definition: type_nlpi.h:176
@ SCIP_NLPTERMSTAT_ITERLIMIT
Definition: type_nlpi.h:175
@ SCIP_NLPTERMSTAT_OUTOFMEMORY
Definition: type_nlpi.h:180
@ SCIP_NLPTERMSTAT_INTERRUPT
Definition: type_nlpi.h:177
@ SCIP_NLPSOLSTAT_UNBOUNDED
Definition: type_nlpi.h:165
@ SCIP_NLPSOLSTAT_LOCINFEASIBLE
Definition: type_nlpi.h:163
@ SCIP_NLPSOLSTAT_FEASIBLE
Definition: type_nlpi.h:162
@ SCIP_NLPSOLSTAT_LOCOPT
Definition: type_nlpi.h:161
@ SCIP_NLPSOLSTAT_UNKNOWN
Definition: type_nlpi.h:166
enum SCIP_NlpTermStat SCIP_NLPTERMSTAT
Definition: type_nlpi.h:184
struct SCIP_NlpiData SCIP_NLPIDATA
Definition: type_nlpi.h:52
@ SCIP_NOMEMORY
Definition: type_retcode.h:44
@ SCIP_OKAY
Definition: type_retcode.h:42
@ SCIP_ERROR
Definition: type_retcode.h:43
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63