Scippy

SCIP

Solving Constraint Integer Programs

nlpi_filtersqp.c
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-2021 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not visit scipopt.org. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file nlpi_filtersqp.c
17  * @ingroup NLPIS
18  * @brief filterSQP NLP interface
19  * @author Stefan Vigerske
20  *
21  * @todo scaling
22  */
23 
24 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
25 
26 #include <string.h>
27 #if defined(_WIN32)
28 #include <windows.h>
29 #else
30 #include <sys/time.h>
31 #endif
32 
33 #ifndef NPARASCIP
34 #include <pthread.h>
35 #endif
36 
37 #include "scip/misc.h"
38 #include "scip/pub_message.h"
39 #include "nlpi/nlpi_filtersqp.h"
40 #include "nlpi/nlpi.h"
41 #include "nlpi/nlpioracle.h"
42 
43 #define NLPI_NAME "filtersqp" /* short concise name of solver */
44 #define NLPI_DESC "Sequential Quadratic Programming trust region solver by R. Fletcher and S. Leyffer" /* description of solver */
45 #define NLPI_PRIORITY -1000 /* priority of NLP solver */
46 
47 #define RANDSEED 26051979 /**< initial random seed */
48 #define MAXPERTURB 0.01 /**< maximal perturbation of bounds in starting point heuristic */
49 #define MAXNRUNS 3 /**< maximal number of FilterSQP calls per NLP solve (several calls if increasing workspace or decreasing eps) */
50 #define WORKSPACEGROWTHFACTOR 2 /**< factor by which to increase worksapce */
51 #define MINEPS 1e-14 /**< minimal FilterSQP epsilon */
52 #define OPTTOLFACTOR 0.5 /**< factor to apply to optimality tolerance, because FilterSQP do scaling */
53 #define DEFAULT_LOBJLIM (real)(-1e100) /**< default lower objective limit (should mean "unlimited") */
54 #define DEFAULT_FEASOPTTOL 1e-6 /**< default feasibility and optimality tolerance */
55 #define DEFAULT_MAXITER 3000 /**< default iteration limit */
56 
57 /*
58  * Data structures
59  */
60 
61 typedef int fint;
62 typedef double real;
63 typedef long ftnlen;
64 
65 struct SCIP_Time
66 {
67  time_t sec; /**< seconds part of time since epoch */
68  long usec; /**< micro-seconds part of time */
69 };
70 typedef struct SCIP_Time SCIP_TIME;
71 
72 struct SCIP_NlpiData
73 {
74  BMS_BLKMEM* blkmem; /**< block memory */
75  SCIP_MESSAGEHDLR* messagehdlr; /**< message handler */
76  SCIP_Real infinity; /**< initial value for infinity */
77  SCIP_RANDNUMGEN* randnumgen; /**< random number generator, if we have to make up a starting point */
78  SCIP_TIME starttime; /**< time at start of solve */
79 };
80 
81 struct SCIP_NlpiProblem
82 {
83  SCIP_NLPIORACLE* oracle; /**< Oracle-helper to store and evaluate NLP */
84 
85  int varssize; /**< length of variables-related arrays, if allocated */
86  int conssize; /**< length of constraints-related arrays, if allocated */
87 
88  SCIP_Real* initguess; /**< initial values for primal variables, or NULL if not known, size varssize */
89  SCIP_Bool warmstart; /**< whether we could warmstart the next solve */
90 
91  SCIP_Real* primalvalues; /**< primal values of variables in solution, size varssize */
92  SCIP_Real* consdualvalues; /**< dual values of constraints in solution, size conssize */
93  SCIP_Real* varlbdualvalues; /**< dual values of variable lower bounds in solution, size varssize */
94  SCIP_Real* varubdualvalues; /**< dual values of variable upper bounds in solution, size varssize */
95 
96  SCIP_NLPSOLSTAT solstat; /**< solution status from last NLP solve */
97  SCIP_NLPTERMSTAT termstat; /**< termination status from last NLP solve */
98  SCIP_Real solvetime; /**< time spend for last NLP solve */
99  int niterations; /**< number of iterations for last NLP solve */
100 
101  SCIP_Bool fromscratch; /**< value of fromscratch parameter */
102  fint istat[14]; /**< integer solution statistics from last FilterSQP call */
103  real rstat[7]; /**< real solution statistics from last FilterSQP call */
104  real feastol; /**< user-given feasibility tolerance */
105  real opttol; /**< user-given optimality tolerance */
106  real fmin; /**< lower bound on objective value */
107  SCIP_Real maxtime; /**< time limit */
108  fint maxiter; /**< iteration limit */
109  fint iprint; /**< print verbosity level */
110 
111  /* cached FilterSQP data */
112  real* x; /**< variable values, size varssize */
113  real* c; /**< constraint value, size conssize */
114  real* lam; /**< duals, size varssize + conssize */
115  real* bl; /**< variable lower bounds and constraint lhs, size varssize + conssize */
116  real* bu; /**< variable upper bounds and constraint rhs, size varssize + conssize */
117  real* s; /**< scaling factors, size varssize + conssize */
118  char* cstype; /**< constraint linearity, size conssize */
119  real* a; /**< gradients values, size la[0]-1 */
120  fint* la; /**< gradients indices, size lasize */
121  int lasize; /**< length of la array */
122  fint* hessiannz; /**< nonzero information about Hessian, size hessiannzsize */
123  int hessiannzsize;/**< length of hessiannz array */
124  real* ws; /**< real workspace, size mxwk */
125  fint* lws; /**< integer workspace, size mxiwk */
126  fint mxwk; /**< size of real workspace */
127  fint mxiwk; /**< size of integer workspace */
128  real* evalbuffer; /**< buffer to cache evaluation results before passing it to FilterSQP, size evalbufsize */
129  int evalbufsize; /**< size of evaluation buffer */
130 };
131 
132 /*
133  * Local methods
134  */
135 
136 #ifdef FNAME_LCASE_DECOR
137 # define F77_FUNC(name,NAME) name ## _
138 #endif
139 #ifdef FNAME_UCASE_DECOR
140 # define F77_FUNC(name,NAME) NAME ## _
141 #endif
142 #ifdef FNAME_LCASE_NODECOR
143 # define F77_FUNC(name,NAME) name
144 #endif
145 #ifdef FNAME_UCASE_NODECOR
146 # define F77_FUNC(name,NAME) NAME
147 #endif
148 
149 /** FilterSQP main routine.
150  *
151  * Array a has length nnza, which is the number of nonzeros in the gradient of the objective and the Jacobian.
152  * The first entries of a is the objective gradient, next are the gradients of the constraints.
153  *
154  * Array la has length lamax, which is at least nnza+m+2.
155  * la contains the index information of a row-oriented sparse matrix storage. It stores the number of nonzeros, the column indices, and the row starts:
156  * la[0] must be set to nnza+1.
157  * la[1]..la[nnza] are indices of the variables corresponding to the entries in a (colidx).
158  * la[nnza+1]..la[nnza+1+m] contain the index where each row starts in a and la (rowstart).
159  */
160 void F77_FUNC(filtersqp,FILTERSQP)(
161  fint* n, /**< number of variables */
162  fint* m, /**< number of constraints (excluding simple bounds) */
163  fint* kmax, /**< maximum size of null-space (at most n) */
164  fint* maxa, /**< maximum number of nonzero entries allowed in Jacobian */
165  fint* maxf, /**< maximum size of the filter - typically 100 */
166  fint* mlp, /**< maximum level parameter for resolving degeneracy in bqpd - typically 100 */
167  fint* mxwk, /**< maximum size of real workspace ws */
168  fint* mxiwk, /**< maximum size of integer workspace lws */
169  fint* iprint, /**< print flag: 0 is quiet, higher is more */
170  fint* nout, /**< output channel - 6 for screen */
171  fint* ifail, /**< fail flag and warmstart indicator */
172  real* rho, /**< initial trust-region radius - default 10 */
173  real* x, /**< starting point and final solution (array of length n) */
174  real* c, /**< final values of general constraints (array of length m) */
175  real* f, /**< final objective value */
176  real* fmin, /**< lower bound on objective value (as termination criteria) */
177  real* bl, /**< lower bounds of variables and constraints (array of length n+m) */
178  real* bu, /**< upper bounds of variables and constraints (array of length n+m) */
179  real* s, /**< scaling factors (array of length n+m) */
180  real* a, /**< objective gradient (always dense) and Jacobian (sparse or dense) entries */
181  fint* la, /**< column indices and length of rows of entries in a (if sparse) or leading dimension of a (if dense) */
182  real* ws, /**< real workspace (array of length mxwk) */
183  fint* lws, /**< integer workspace (array of length mxiwk) */
184  real* lam, /**< Lagrangian multipliers of simple bounds and general constraints (array of length n+m) */
185  char* cstype, /**< indicator whether constraint is linear ('L') or nonlinear ('N') (array of length m) */
186  real* user, /**< real workspace passed through to user routines */
187  fint* iuser, /**< integer workspace passed through to user routines */
188  fint* maxiter, /**< iteration limit for SQP solver */
189  fint* istat, /**< integer space for solution statistics (array of length 14) */
190  real* rstat, /**< real space for solution statitics (array of length 7) */
191  ftnlen cstype_len /**< 1 ??? */
192  );
193 
194 void F77_FUNC(objfun,OBJFUN)(real *x, fint *n, real *f, real *user, fint *iuser,
195  fint *errflag);
196 
197 void F77_FUNC(confun,CONFUN)(real *x, fint *n, fint *m, real *c, real *a, fint *la,
198  real *user, fint *iuser, fint *errflag);
199 
200 /* TODO what are the arguments of this function and does it need to be implemented?
201  * it's not in the filterSQP manual, but its an undefined symbol in the filterSQP library
202 void F77_FUNC(objgrad,OBJGRAD)(fint *, fint *, fint *, real *, real *, fint *, fint
203  *, real *, fint *, fint *);
204 */
205 void F77_FUNC(objgrad,OBJGRAD)(void);
206 
207 void F77_FUNC(gradient,GRADIENT)(fint *n, fint *m, fint *mxa, real *x, real *a, fint *la,
208  fint *maxa, real *user, fint *iuser, fint *errflag);
209 
210 void F77_FUNC(hessian,HESSIAN)(real *x, fint *n, fint *m, fint *phase, real *lam,
211  real *ws, fint *lws, real *user, fint *iuser,
212  fint *l_hess, fint *li_hess, fint *errflag);
213 
214 /** common block for problemname */
215 /*lint -esym(754,char_l,pname,*::char_l,*::pname) */
216 extern struct
217 {
218  fint char_l;
219  char pname[10];
220 } F77_FUNC(cpname,CPNAME);
221 /*lint -esym(752,cpname_) */
222 
223 /** common block for Hessian storage set to 0, i.e. NO Hessian */
224 /*lint -esym(754,*::phr,*::phc) */
225 extern struct
226 {
228 } F77_FUNC(hessc,HESSC);
229 /*lint -esym(754,phr,phc) */
230 
231 /** common block for upper bound on filter */
232 extern struct
233 {
235 } F77_FUNC(ubdc,UBDC);
236 
237 /** common block for infinity & epsilon */
238 extern struct
239 {
241 } F77_FUNC(nlp_eps_inf,NLP_EPS_INF);
242 
243 /** common block for printing from QP solver */
244 /*lint -esym(754,*::n_bqpd_calls,*::n_bqpd_prfint) */
245 extern struct
246 {
248 } F77_FUNC(bqpd_count,BQPD_COUNT);
249 /*lint -esym(752,bqpd_count_) */
250 /*lint -esym(754,n_bqpd_calls,n_bqpd_prfint) */
251 
252 /** common for scaling: scale_mode = 0 (none), 1 (variables), 2 (vars+cons) */
253 /*lint -esym(754,*::phe) */
254 extern struct
255 {
257 } F77_FUNC(scalec,SCALEC);
258 /*lint -esym(754,phe) */
259 
260 #ifndef NPARASCIP
261 static pthread_mutex_t filtersqpmutex = PTHREAD_MUTEX_INITIALIZER;
262 #endif
263 
264 static
266 {
267  SCIP_TIME t;
268 #ifndef _WIN32
269  struct timeval tp; /*lint !e86*/
270 #endif
271 
272 #ifdef _WIN32
273  t.sec = time(NULL);
274  t.usec = 0;
275 
276 #else
277  (void) gettimeofday(&tp, NULL);
278  t.sec = tp.tv_sec;
279  t.usec = tp.tv_usec;
280 #endif
281 
282  return t;
283 }
284 
285 /* gives time since starttime (in problem) */
286 static
288  SCIP_NLPIDATA* nlpidata /**< NLPI data */
289  )
290 {
291  SCIP_TIME now;
292 
293  assert(nlpidata != NULL);
294 
295  now = gettime();
296 
297  /* now should be after startime */
298  assert(now.sec >= nlpidata->starttime.sec);
299  assert(now.sec > nlpidata->starttime.sec || now.usec >= nlpidata->starttime.usec);
300 
301  return (SCIP_Real)(now.sec - nlpidata->starttime.sec) + 1e-6 * (SCIP_Real)(now.usec - nlpidata->starttime.usec);
302 }
303 
304 static
306  SCIP_NLPIDATA* nlpidata, /**< NLPI data */
307  SCIP_NLPIPROBLEM* nlpiproblem /**< NLPI problem */
308  )
309 {
310  if( nlpiproblem->maxtime == DBL_MAX ) /*lint !e777 */
311  return FALSE;
312 
313  return timeelapsed(nlpidata) >= nlpiproblem->maxtime;
314 }
315 
316 /** Objective function evaluation */ /*lint -e{715} */
317 void F77_FUNC(objfun,OBJFUN)(
318  real* x, /**< value of current variables (array of length n) */
319  fint* n, /**< number of variables */
320  real* f, /**< buffer to store value of objective function */
321  real* user, /**< user real workspace */
322  fint* iuser, /**< user integer workspace */
323  fint* errflag /**< set to 1 if arithmetic exception occurs, otherwise 0 */
324  )
325 { /*lint --e{715} */
326  SCIP_NLPIPROBLEM* problem;
327  real val;
328 
329  problem = (SCIP_NLPIPROBLEM*)(void*)iuser;
330  assert(problem != NULL);
331 
332  if( timelimitreached((SCIP_NLPIDATA*)(void*)user, problem) )
333  {
334  SCIPdebugMessage("timelimit reached, issuing arithmetic exception in objfun\n");
335  *errflag = 1;
336  return;
337  }
338 
339  *errflag = 1;
340  if( SCIPnlpiOracleEvalObjectiveValue(problem->oracle, x, &val) == SCIP_OKAY && SCIPisFinite(val) )
341  {
342  *errflag = 0;
343  *f = val;
344  }
345  else
346  {
347  SCIPdebugMessage("arithmetic exception in objfun\n");
348  }
349 }
350 
351 /** Constraint functions evaluation */ /*lint -e{715} */
352 void F77_FUNC(confun,CONFUN)(
353  real* x, /**< value of current variables (array of length n) */
354  fint* n, /**< number of variables */
355  fint* m, /**< number of constraints */
356  real* c, /**< buffer to store values of constraints (array of length m) */
357  real* a, /**< Jacobian matrix entries */
358  fint* la, /**< Jacobian index information */
359  real* user, /**< user real workspace */
360  fint* iuser, /**< user integer workspace */
361  fint* errflag /**< set to 1 if arithmetic exception occurs, otherwise 0 */
362  )
363 { /*lint --e{715} */
364  SCIP_NLPIPROBLEM* problem;
365  real val;
366  int j;
367 
368  problem = (SCIP_NLPIPROBLEM*)(void*)iuser;
369  assert(problem != NULL);
370 
371  *errflag = 0;
372  for( j = 0; j < *m; ++j )
373  {
374  if( SCIPnlpiOracleEvalConstraintValue(problem->oracle, j, x, &val) != SCIP_OKAY || !SCIPisFinite(val) )
375  {
376  *errflag = 1;
377  SCIPdebugMessage("arithmetic exception in confun for constraint %d\n", j);
378  break;
379  }
380  c[j] = val;
381  }
382 }
383 
384 /** Objective gradient and Jacobian evaluation
385  *
386  * \note If an arithmetic exception occurred, then the gradients must not be modified.
387  */ /*lint -e{715} */
388 void
389 F77_FUNC(gradient,GRADIENT)(
390  fint* n, /**< number of variables */
391  fint* m, /**< number of constraints */
392  fint* mxa, /**< actual number of entries in a */
393  real* x, /**< value of current variables (array of length n) */
394  real* a, /**< Jacobian matrix entries */
395  fint* la, /**< Jacobian index information: column indices and pointers to start of each row */
396  fint* maxa, /**< maximal size of a */
397  real* user, /**< user real workspace */
398  fint* iuser, /**< user integer workspace */
399  fint* errflag /**< set to 1 if arithmetic exception occurs, otherwise 0 */
400  )
401 { /*lint --e{715} */
402  SCIP_NLPIPROBLEM* problem;
403  SCIP_Real dummy;
404 
405  problem = (SCIP_NLPIPROBLEM*)(void*)iuser;
406  assert(problem != NULL);
407  assert(problem->evalbuffer != NULL);
408  assert(problem->evalbufsize >= *maxa);
409 
410  *errflag = 1;
411 
412  if( SCIPnlpiOracleEvalObjectiveGradient(problem->oracle, x, TRUE, &dummy, problem->evalbuffer) == SCIP_OKAY )
413  {
414  if( SCIPnlpiOracleEvalJacobian(problem->oracle, x, TRUE, NULL, problem->evalbuffer+*n) == SCIP_OKAY )
415  {
416  BMScopyMemoryArray(a, problem->evalbuffer, *maxa);
417  *errflag = 0;
418  }
419  else
420  {
421  SCIPdebugMessage("arithmetic exception in gradient for constraints\n");
422  }
423  }
424  else
425  {
426  SCIPdebugMessage("arithmetic exception in gradient for objective\n");
427  }
428 }
429 
430 /* Objective gradient evaluation */
431 /*
432 void F77_FUNC(objgrad,OBJGRAD)(
433  fint*,
434  fint*,
435  fint*,
436  real*,
437  real*,
438  fint*,
439  fint*,
440  real*,
441  fint*,
442  fint*
443  )
444 */
445 void F77_FUNC(objgrad,OBJGRAD)(void)
446 {
447  SCIPerrorMessage("Objgrad not implemented. Should not be called.\n");
448 }
449 
450 /** Hessian of the Lagrangian evaluation
451  *
452  * phase = 1 : Hessian of the Lagrangian without objective Hessian
453  * phase = 2 : Hessian of the Lagrangian (including objective Hessian)
454  *
455  * \note If an arithmetic exception occurred, then the Hessian must not be modified.
456  */ /*lint -e{715} */
457 void
458 F77_FUNC(hessian,HESSIAN)(
459  real* x, /**< value of current variables (array of length n) */
460  fint* n, /**< number of variables */
461  fint* m, /**< number of constraints */
462  fint* phase, /**< indicates what kind of Hessian matrix is required */
463  real* lam, /**< Lagrangian multipliers (array of length n+m) */
464  real* ws, /**< real workspace for Hessian, passed to Wdotd */
465  fint* lws, /**< integer workspace for Hessian, passed to Wdotd */
466  real* user, /**< user real workspace */
467  fint* iuser, /**< user integer workspace */
468  fint* l_hess, /**< space of Hessian real storage ws. On entry: maximal space allowed, on exit: actual amount used */
469  fint* li_hess, /**< space of Hessian integer storage lws. On entry: maximal space allowed, on exit: actual amount used */
470  fint* errflag /**< set to 1 if arithmetic exception occurs, otherwise 0 */
471  )
472 { /*lint --e{715} */
473  SCIP_NLPIPROBLEM* problem;
474  SCIP_Real* lambda;
475  int nnz;
476  int i;
477 
478  problem = (SCIP_NLPIPROBLEM*)(void*)iuser;
479  assert(problem != NULL);
480  assert(problem->evalbuffer != NULL);
481 
482  nnz = problem->hessiannz[0]-1;
483  assert(problem->evalbufsize >= nnz);
484 
485  *errflag = 1;
486 
487  /* initialize lambda to -lam */
488  BMSallocMemoryArray(&lambda, *m);
489  for( i = 0; i < *m; ++i )
490  lambda[i] = -lam[*n+i];
491 
492  if( SCIPnlpiOracleEvalHessianLag(problem->oracle, x, TRUE, (*phase == 1) ? 0.0 : 1.0, lambda, problem->evalbuffer) == SCIP_OKAY )
493  {
494  *l_hess = nnz;
495 
496  BMScopyMemoryArray(ws, problem->evalbuffer, nnz);
497 
498  *errflag = 0;
499 
500  /* copy the complete problem->hessiannz into lws */
501  for( i = 0; i < nnz + *n + 2; ++i )
502  lws[i] = problem->hessiannz[i];
503  *li_hess = nnz + *n + 2;
504  }
505  else
506  {
507  SCIPdebugMessage("arithmetic exception in hessian\n");
508  }
509 
510  BMSfreeMemoryArray(&lambda);
511 }
512 
513 
514 
515 static
517  BMS_BLKMEM* blkmem, /**< block memory */
518  SCIP_NLPIORACLE* oracle, /**< NLPI oracle */
519  fint** la, /**< buffer to store pointer to sparsity structure */
520  int* lasize, /**< buffer to store length of *la array */
521  real** a /**< buffer to store pointer to value buffer */
522  )
523 {
524  const int* offset;
525  const int* col;
526  int nnz; /* number of nonzeros in Jacobian */
527  int nvars;
528  int ncons;
529  int i;
530  int c;
531 
532  assert(la != NULL);
533  assert(lasize != NULL);
534  assert(a != NULL);
535  assert(*la == NULL);
536  assert(*a == NULL);
537 
538  nvars = SCIPnlpiOracleGetNVars(oracle);
539  ncons = SCIPnlpiOracleGetNConstraints(oracle);
540 
541  /* get jacobian sparsity in oracle format: offset are rowstarts in col and col are column indices */
542  SCIP_CALL( SCIPnlpiOracleGetJacobianSparsity(oracle, &offset, &col) );
543  nnz = offset[ncons];
544 
545  /* la stores both column indices (first) and rowstarts (at the end) of the objective gradient and Jacobian
546  * For the objective gradient, always n entries are taken, the Jacobian is stored sparse
547  * la(0) = n+nnz+1 position where rowstarts start in la
548  * la(j) column index of objective gradient or Jacobian row, rowwise
549  * la(la(0)) position of first nonzero element for objective gradient in a()
550  * la(la(0)+i) position of first nonzero element for constraint i gradient in a(), i=1..m
551  * la(la(0)+m+1) = n+nnz first unused position in a
552  * where n = nvars and m = ncons
553  */
554 
555  /* need space for la(0) and column indices and rowstarts (1+ncons+1 for objective, constraints, and end (offsets[ncons])) */
556  *lasize = 1 + (nvars+nnz) + 1+ncons + 1;
557  SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, la, *lasize) );
558 
559  (*la)[0] = nvars+nnz+1;
560 
561  /* the objective gradient is stored in sparse form */
562  for( i = 0; i < nvars; ++i )
563  (*la)[1+i] = 1+i; /* shift by 1 for Fortran */
564  (*la)[(*la)[0]] = 1; /* objective entries start at the beginning in a, shift by 1 for Fortran */
565 
566  /* column indicies are as given by col */
567  for( i = 0; i < nnz; ++i )
568  (*la)[1+nvars+i] = col[i] + 1; /* shift by 1 for Fortran */
569 
570  /* rowstarts are as given by offset, plus extra for objective gradient */
571  for( c = 0; c <= ncons; ++c )
572  (*la)[(*la)[0]+1+c] = offset[c] + nvars + 1; /* shift by nvars for objective, shift by 1 for Fortran */
573 
574  SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, a, nvars + nnz) );
575 
576 #ifdef SCIP_MORE_DEBUG /* enable for debugging Jacobian */
577  for( i = 0; i < 1 + (nvars+nnz) + 1+ncons + 1; ++i )
578  printf("la[%2d] = %2d\n", i, (*la)[i]);
579 #endif
580 
581  return SCIP_OKAY;
582 }
583 
584 static
586  BMS_BLKMEM* blkmem, /**< block memory */
587  SCIP_NLPIORACLE* oracle, /**< NLPI oracle */
588  fint** la, /**< buffer to store pointer to Hessian sparsity structure */
589  int* lasize /**< buffer to store length of *la array */
590  )
591 {
592  const int* offset;
593  const int* col;
594  int nnz; /* number of nonzeros in Jacobian */
595  int nvars;
596  int i;
597  int v;
598 
599  assert(la != NULL);
600  assert(lasize != NULL);
601  assert(*la == NULL);
602 
603  nvars = SCIPnlpiOracleGetNVars(oracle);
604 
605  /* get Hessian sparsity in oracle format: offset are rowstarts in col and col are column indices */
606  SCIP_CALL( SCIPnlpiOracleGetHessianLagSparsity(oracle, &offset, &col) );
607  nnz = offset[nvars];
608 
609  /* la stores both column indices (first) and rowstarts (at the end) of the (sparse) Hessian
610  * la(0) = nnz+1 position where rowstarts start in la
611  * la(j) column index of Hessian row, rowwise
612  * la(la(0)+i) position of first nonzero element for row i, i = 0..n-1
613  * la(la(0)+n) = nnz first unused position in Hessian
614  * where n = nvars
615  */
616 
617  /* need space for la(0) and column indices and rowstarts
618  * 1 for la(0)
619  * nnz for column indices
620  * nvars for rowstarts
621  * 1 for first unused position
622  */
623  *lasize = 1 + nnz + nvars + 1;
624  SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, la, *lasize) );
625 
626  (*la)[0] = nnz+1;
627 
628  /* column indicies are as given by col */
629  for( i = 0; i < nnz; ++i )
630  (*la)[1+i] = col[i] + 1; /* shift by 1 for Fortran */
631 
632  /* rowstarts are as given by offset */
633  for( v = 0; v <= nvars; ++v )
634  (*la)[(*la)[0]+v] = offset[v] + 1; /* shift by 1 for Fortran */
635 
636  F77_FUNC(hessc,HESSC).phl = 1;
637 
638 #ifdef SCIP_MORE_DEBUG /* enable for debugging Hessian */
639  for( i = 0; i < 1 + nnz + nvars + 1; ++i )
640  printf("lw[%2d] = %2d\n", i, (*la)[i]);
641 #endif
642 
643  return SCIP_OKAY;
644 }
645 
646 /** setup starting point for FilterSQP */
647 static
649  SCIP_NLPIDATA* data, /**< NLPI data */
650  SCIP_NLPIPROBLEM* problem, /**< NLPI problem */
651  real* x, /**< array to store initial values */
652  SCIP_Bool* success /**< whether we could setup a point in which functions could be evaluated */
653  )
654 {
655  int i;
656  int n;
657  SCIP_Real val;
658 
659  assert(data != NULL);
660  assert(problem != NULL);
661  assert(x != NULL);
662  assert(success != NULL);
663 
664  n = SCIPnlpiOracleGetNVars(problem->oracle);
665 
666  /* setup starting point */
667  if( problem->initguess != NULL )
668  {
669  for( i = 0; i < n; ++i )
670  x[i] = problem->initguess[i];
671  }
672  else
673  {
674  SCIP_Real lb;
675  SCIP_Real ub;
676 
677  SCIPdebugMessage("FilterSQP started without initial primal values; make up something by projecting 0 onto variable bounds and perturb\n");
678 
679  if( data->randnumgen == NULL )
680  {
681  SCIP_CALL( SCIPrandomCreate(&data->randnumgen, data->blkmem, RANDSEED) );
682  }
683 
684  for( i = 0; i < n; ++i )
685  {
686  lb = SCIPnlpiOracleGetVarLbs(problem->oracle)[i];
687  ub = SCIPnlpiOracleGetVarUbs(problem->oracle)[i];
688  if( lb > 0.0 )
689  x[i] = SCIPrandomGetReal(data->randnumgen, lb, lb + MAXPERTURB*MIN(1.0, ub-lb));
690  else if( ub < 0.0 )
691  x[i] = SCIPrandomGetReal(data->randnumgen, ub - MAXPERTURB*MIN(1.0, ub-lb), ub);
692  else
693  x[i] = SCIPrandomGetReal(data->randnumgen, MAX(lb, -MAXPERTURB*MIN(1.0, ub-lb)), MIN(ub, MAXPERTURB*MIN(1.0, ub-lb)));
694  }
695  }
696 
697  /* check whether objective and constraints can be evaluated and differentiated once in starting point
698  * NOTE: this does not check whether hessian can be computed!
699  */
700  *success = SCIPnlpiOracleEvalObjectiveValue(problem->oracle, x, &val) == SCIP_OKAY && SCIPisFinite(val);
701  *success &= SCIPnlpiOracleEvalObjectiveGradient(problem->oracle, x, FALSE, &val, problem->evalbuffer) == SCIP_OKAY; /*lint !e514*/
702  i = 0;
703  for( ; *success && i < SCIPnlpiOracleGetNConstraints(problem->oracle); ++i )
704  *success = SCIPnlpiOracleEvalConstraintValue(problem->oracle, i, x, &val) == SCIP_OKAY && SCIPisFinite(val);
705  *success &= SCIPnlpiOracleEvalJacobian(problem->oracle, x, FALSE, NULL, problem->evalbuffer) == SCIP_OKAY; /*lint !e514*/
706 
707  if( !*success )
708  {
709  SCIPdebugMessage("could not evaluate or constraint %d in %s starting point or Jacobian not available\n", i-1, problem->initguess != NULL ? "provided" : "made up");
710 
711  if( problem->initguess != NULL )
712  {
713  /* forget given starting point and try to make up our own */
714  BMSfreeBlockMemoryArray(data->blkmem, &problem->initguess, problem->varssize);
715  SCIP_CALL( setupStart(data, problem, x, success) );
716  }
717  }
718 
719  return SCIP_OKAY;
720 }
721 
722 
723 
724 /** sets the solstat and termstat to unknown and other, resp. */
725 static
727  SCIP_NLPIPROBLEM* problem /**< data structure of problem */
728  )
729 {
730  assert(problem != NULL);
731 
732  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
733  problem->termstat = SCIP_NLPTERMSTAT_OTHER;
734 }
735 
736 
737 /** processes results from FilterSQP call */
738 static
740  SCIP_NLPIDATA* nlpidata, /**< NLPI data */
741  SCIP_NLPIPROBLEM* problem, /**< NLPI problem */
742  fint ifail, /**< fail flag from FilterSQP call */
743  real* x, /**< primal solution values from FilterSQP call, or NULL if stopped before filtersqp got called */
744  real* lam /**< dual solution values from FilterSQP call, or NULL if stopped before filtersqp got called */
745  )
746 {
747  int i;
748  int nvars;
749  int ncons;
750 
751  assert(problem != NULL);
752  assert(ifail >= 0);
753  assert((x != NULL) == (lam != NULL));
754 
755  problem->solvetime = timeelapsed(nlpidata);
756 
757  nvars = SCIPnlpiOracleGetNVars(problem->oracle);
758  ncons = SCIPnlpiOracleGetNConstraints(problem->oracle);
759 
760  if( ifail <= 8 && x != NULL )
761  {
762  /* FilterSQP terminated somewhat normally -> store away solution */
763 
764  /* make sure we have memory for solution */
765  if( problem->primalvalues == NULL )
766  {
767  assert(problem->varssize >= nvars); /* ensured in nlpiAddVariables */
768  assert(problem->conssize >= ncons); /* ensured in nlpiAddConstraints */
769  assert(problem->consdualvalues == NULL); /* if primalvalues == NULL, then also consdualvalues should be NULL */
770  assert(problem->varlbdualvalues == NULL); /* if primalvalues == NULL, then also varlbdualvalues should be NULL */
771  assert(problem->varubdualvalues == NULL); /* if primalvalues == NULL, then also varubdualvalues should be NULL */
772 
773  SCIP_ALLOC( BMSallocBlockMemoryArray(nlpidata->blkmem, &problem->primalvalues, problem->varssize) );
774  SCIP_ALLOC( BMSallocBlockMemoryArray(nlpidata->blkmem, &problem->consdualvalues, problem->conssize) );
775  SCIP_ALLOC( BMSallocBlockMemoryArray(nlpidata->blkmem, &problem->varlbdualvalues, problem->varssize) );
776  SCIP_ALLOC( BMSallocBlockMemoryArray(nlpidata->blkmem, &problem->varubdualvalues, problem->varssize) );
777  }
778  else
779  {
780  assert(problem->consdualvalues != NULL);
781  assert(problem->varlbdualvalues != NULL);
782  assert(problem->varubdualvalues != NULL);
783  }
784 
785  for( i = 0; i < nvars; ++i )
786  {
787  problem->primalvalues[i] = x[i];
788 
789  /* if dual from FilterSQP is positive, then it belongs to the lower bound, otherwise to the upper bound */
790  problem->varlbdualvalues[i] = MAX(0.0, lam[i]);
791  problem->varubdualvalues[i] = MAX(0.0, -lam[i]);
792  }
793 
794  /* duals from FilterSQP are negated */
795  for( i = 0; i < ncons; ++i )
796  problem->consdualvalues[i] = -lam[nvars + i];
797  }
798 
799  /* translate ifail to solution and termination status and decide whether we could warmstart next */
800  problem->warmstart = FALSE;
801  switch( ifail )
802  {
803  case 0: /* successful run, solution found */
804  assert(problem->rstat[4] <= problem->feastol); /* should be feasible */
805  problem->solstat = (problem->rstat[0] <= problem->opttol ? SCIP_NLPSOLSTAT_LOCOPT : SCIP_NLPSOLSTAT_FEASIBLE);
806  problem->termstat = SCIP_NLPTERMSTAT_OKAY;
807  problem->warmstart = TRUE;
808  break;
809  case 1: /* unbounded, feasible point with f(x) <= fmin */
810  assert(problem->rstat[4] <= problem->feastol); /* should be feasible */
812  if( problem->fmin == DEFAULT_LOBJLIM ) /*lint !e777*/
813  problem->termstat = SCIP_NLPTERMSTAT_OKAY; /* fmin was not set */
814  else
816  break;
817  case 2: /* linear constraints are inconsistent */
819  problem->termstat = SCIP_NLPTERMSTAT_OKAY;
820  break;
821  case 3: /* (locally) nonlinear infeasible, minimal-infeasible solution found */
822  /* problem->solstat = (problem->rstat[0] <= problem->opttol ? SCIP_NLPSOLSTAT_LOCINFEASIBLE : SCIP_NLPSOLSTAT_UNKNOWN); */
823  problem->solstat = SCIP_NLPSOLSTAT_LOCINFEASIBLE; /* TODO FilterSQP does not set rstat[0] in this case, assuming local infeasibility is valid */
824  problem->termstat = SCIP_NLPTERMSTAT_OKAY;
825  problem->warmstart = TRUE;
826  break;
827  case 4: /* terminate at point with h(x) <= eps (constraint violation below epsilon) but QP infeasible */
828  assert(problem->rstat[4] <= problem->feastol); /* should be feasible */
831  problem->warmstart = TRUE;
832  break;
833  case 5: /* termination with rho < eps (trust region radius below epsilon) */
834  if( problem->rstat[4] <= problem->feastol )
836  else
837  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
839  problem->warmstart = TRUE;
840  break;
841  case 6: /* termination with iter > max_iter */
842  if( problem->rstat[4] <= problem->feastol )
844  else
845  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
846  problem->termstat = SCIP_NLPTERMSTAT_ITLIM;
847  problem->warmstart = TRUE;
848  break;
849  case 7: /* crash in user routine (IEEE error) could not be resolved, or timelimit reached */
850  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
851  if( problem->solvetime >= problem->maxtime )
852  {
853  problem->termstat = SCIP_NLPTERMSTAT_TILIM;
854  problem->warmstart = TRUE;
855  }
856  else
858  break;
859  case 8: /* unexpect ifail from QP solver */
860  if( problem->rstat[4] <= problem->feastol )
862  else
863  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
864  problem->termstat = SCIP_NLPTERMSTAT_OTHER;
865  break;
866  case 9: /* not enough REAL workspace */
867  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
869  break;
870  case 10: /* not enough INTEGER workspace */
871  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
873  break;
874  default:
875  problem->solstat = SCIP_NLPSOLSTAT_UNKNOWN;
876  problem->termstat = SCIP_NLPTERMSTAT_OTHER;
877  break;
878  }
879 
880  return SCIP_OKAY;
881 }
882 
883 
884 /** calculate memory size for dynamically allocated arrays (copied from scip/set.c) */
885 static
887  int num /**< minimum number of entries to store */
888  )
889 {
890  int size;
891 
892  /* calculate the size with this loop, such that the resulting numbers are always the same (-> block memory) */
893  size = 4;
894  while( size < num )
895  size = (int)(1.2 * size + 4);
896 
897  return size;
898 }
899 
900 
901 /*
902  * Callback methods of NLP solver interface
903  */
904 
905 /** copy method of NLP interface (called when SCIP copies plugins)
906  *
907  * input:
908  * - blkmem block memory in target SCIP
909  * - sourcenlpi the NLP interface to copy
910  * - targetnlpi buffer to store pointer to copy of NLP interface
911  */
912 static
913 SCIP_DECL_NLPICOPY( nlpiCopyFilterSQP )
914 {
915  SCIP_NLPIDATA* sourcedata;
916 
917  assert(sourcenlpi != NULL);
918  assert(targetnlpi != NULL);
919 
920  sourcedata = SCIPnlpiGetData(sourcenlpi);
921  assert(sourcedata != NULL);
922 
923  SCIP_CALL( SCIPcreateNlpSolverFilterSQP(blkmem, targetnlpi) );
924  assert(*targetnlpi != NULL);
925 
926  SCIP_CALL( SCIPnlpiSetRealPar(*targetnlpi, NULL, SCIP_NLPPAR_INFINITY, sourcedata->infinity) );
927  SCIP_CALL( SCIPnlpiSetMessageHdlr(*targetnlpi, sourcedata->messagehdlr) );
928 
929  return SCIP_OKAY; /*lint !e527*/
930 } /*lint !e715*/
931 
932 /** destructor of NLP interface to free nlpi data
933  *
934  * input:
935  * - nlpi datastructure for solver interface
936  */
937 static
938 SCIP_DECL_NLPIFREE( nlpiFreeFilterSQP )
939 {
940  SCIP_NLPIDATA* data;
941 
942  assert(nlpi != NULL);
943 
944  data = SCIPnlpiGetData(nlpi);
945  assert(data != NULL);
946 
947  if( data->randnumgen != NULL )
948  {
949  SCIPrandomFree(&data->randnumgen, data->blkmem);
950  }
951 
952  BMSfreeBlockMemory(data->blkmem, &data);
953  assert(data == NULL);
954 
955  return SCIP_OKAY;
956 }
957 
958 /** gets pointer for NLP solver
959  *
960  * to do dirty stuff
961  *
962  * input:
963  * - nlpi datastructure for solver interface
964  *
965  * return: void pointer to solver
966  */
967 static
968 SCIP_DECL_NLPIGETSOLVERPOINTER(nlpiGetSolverPointerFilterSQP)
969 {
970  return NULL; /*lint !e527*/
971 } /*lint !e715*/
972 
973 /** creates a problem instance
974  *
975  * input:
976  * - nlpi datastructure for solver interface
977  * - problem pointer to store the problem data
978  * - name name of problem, can be NULL
979  */
980 static
981 SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemFilterSQP)
982 {
983  SCIP_NLPIDATA* data;
984 
985  assert(nlpi != NULL);
986  assert(problem != NULL);
987 
988  data = SCIPnlpiGetData(nlpi);
989  assert(data != NULL);
990 
991  SCIP_ALLOC( BMSallocBlockMemory(data->blkmem, problem) );
992  BMSclearMemory(*problem);
993 
994  SCIP_CALL( SCIPnlpiOracleCreate(data->blkmem, &(*problem)->oracle) );
995  SCIP_CALL( SCIPnlpiOracleSetInfinity((*problem)->oracle, data->infinity) );
996  SCIP_CALL( SCIPnlpiOracleSetProblemName((*problem)->oracle, name) );
997 
998  (*problem)->feastol = DEFAULT_FEASOPTTOL;
999  (*problem)->opttol = DEFAULT_FEASOPTTOL;
1000  (*problem)->fmin = DEFAULT_LOBJLIM;
1001  (*problem)->maxtime = DBL_MAX;
1002  (*problem)->maxiter = DEFAULT_MAXITER;
1003  (*problem)->iprint = 0;
1004 
1005  invalidateSolution(*problem);
1006 
1007  return SCIP_OKAY; /*lint !e527*/
1008 } /*lint !e715*/
1009 
1010 /** free a problem instance
1011  *
1012  * input:
1013  * - nlpi datastructure for solver interface
1014  * - problem pointer where problem data is stored
1015  */
1016 static
1017 SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemFilterSQP)
1018 {
1019  SCIP_NLPIDATA* data;
1020 
1021  assert(nlpi != NULL);
1022  assert(problem != NULL);
1023  assert(*problem != NULL);
1024 
1025  data = SCIPnlpiGetData(nlpi);
1026  assert(data != NULL);
1027 
1028  if( (*problem)->oracle != NULL )
1029  {
1030  SCIP_CALL( SCIPnlpiOracleFree(&(*problem)->oracle) );
1031  }
1032 
1033  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->initguess, (*problem)->varssize);
1034  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->primalvalues, (*problem)->varssize);
1035  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->consdualvalues, (*problem)->conssize);
1036  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->varlbdualvalues, (*problem)->varssize);
1037  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->varubdualvalues, (*problem)->varssize);
1038  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->cstype, (*problem)->conssize);
1039  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->s, (*problem)->varssize + (*problem)->conssize); /*lint !e776 */
1040  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->bu, (*problem)->varssize + (*problem)->conssize); /*lint !e776 */
1041  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->bl, (*problem)->varssize + (*problem)->conssize); /*lint !e776 */
1042  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->x, (*problem)->varssize);
1043  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->c, (*problem)->conssize);
1044  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->lam, (*problem)->varssize + (*problem)->conssize); /*lint !e776 */
1045  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->a, (*problem)->la != NULL ? (*problem)->la[0]-1 : 0);
1046  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->la, (*problem)->lasize);
1047  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->hessiannz, (*problem)->hessiannzsize);
1048  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->lws, (*problem)->mxiwk);
1049  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->ws, (*problem)->mxwk);
1050  BMSfreeBlockMemoryArrayNull(data->blkmem, &(*problem)->evalbuffer, (*problem)->evalbufsize);
1051 
1052  BMSfreeBlockMemory(data->blkmem, problem);
1053  assert(*problem == NULL);
1054 
1055  return SCIP_OKAY;
1056 } /*lint !e715*/
1057 
1058 /** gets pointer to solver-internal problem instance
1059  *
1060  * to do dirty stuff
1061  *
1062  * input:
1063  * - nlpi datastructure for solver interface
1064  * - problem datastructure for problem instance
1065  *
1066  * return: void pointer to problem instance
1067  */
1068 static
1069 SCIP_DECL_NLPIGETPROBLEMPOINTER(nlpiGetProblemPointerFilterSQP)
1070 {
1071  return NULL; /*lint !e527*/
1072 } /*lint !e715*/
1073 
1074 /** add variables
1075  *
1076  * input:
1077  * - nlpi datastructure for solver interface
1078  * - problem datastructure for problem instance
1079  * - nvars number of variables
1080  * - lbs lower bounds of variables, can be NULL if -infinity
1081  * - ubs upper bounds of variables, can be NULL if +infinity
1082  * - varnames names of variables, can be NULL
1083  */
1084 static
1085 SCIP_DECL_NLPIADDVARS( nlpiAddVarsFilterSQP )
1086 {
1087  SCIP_NLPIDATA* nlpidata;
1088  int oldnvars;
1089 
1090  assert(nlpi != NULL);
1091  assert(problem != NULL);
1092  assert(problem->oracle != NULL);
1093 
1094  nlpidata = SCIPnlpiGetData(nlpi);
1095  assert(nlpidata != NULL);
1096 
1097  oldnvars = SCIPnlpiOracleGetNVars(problem->oracle);
1098 
1099  SCIP_CALL( SCIPnlpiOracleAddVars(problem->oracle, nvars, lbs, ubs, varnames) );
1100 
1101  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->initguess, problem->varssize);
1102  invalidateSolution(problem);
1103  problem->warmstart = FALSE;
1104 
1105  /* increase variables-related arrays in problem, if necessary */
1106  if( problem->varssize < SCIPnlpiOracleGetNVars(problem->oracle) )
1107  {
1108  int newsize = calcGrowSize(SCIPnlpiOracleGetNVars(problem->oracle));
1109 
1110  if( problem->primalvalues != NULL )
1111  {
1112  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->primalvalues, problem->varssize, newsize) );
1113  }
1114 
1115  if( problem->varlbdualvalues != NULL )
1116  {
1117  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->varlbdualvalues, problem->varssize, newsize) );
1118  }
1119 
1120  if( problem->varubdualvalues != NULL )
1121  {
1122  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->varubdualvalues, problem->varssize, newsize) );
1123  }
1124 
1125  if( problem->x != NULL )
1126  {
1127  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->x, problem->varssize, newsize) );
1128  }
1129 
1130  if( problem->lam != NULL )
1131  {
1132  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->lam, problem->varssize + problem->conssize, newsize + problem->conssize) ); /*lint !e776*/
1133  }
1134 
1135  if( problem->bl != NULL )
1136  {
1137  assert(problem->bu != NULL);
1138  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->bl, problem->varssize + problem->conssize, newsize + problem->conssize) ); /*lint !e776*/
1139  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->bu, problem->varssize + problem->conssize, newsize + problem->conssize) ); /*lint !e776*/
1140  }
1141 
1142  if( problem->s != NULL )
1143  {
1144  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->s, problem->varssize + problem->conssize, newsize + problem->conssize) ); /*lint !e776*/
1145  }
1146 
1147  problem->varssize = newsize;
1148  }
1149 
1150  /* update variable bounds in FilterSQP data */
1151  if( problem->bl != NULL )
1152  {
1153  int i;
1154  int nconss;
1155 
1156  nconss = SCIPnlpiOracleGetNConstraints(problem->oracle);
1157 
1158  /* bl and bu have first variable bounds and then constraint sides
1159  * copy the constraint sides to their new position before putting in the new variable bounds
1160  */
1161  for( i = nconss-1; i >= 0; --i )
1162  {
1163  problem->bl[oldnvars+nvars+i] = problem->bl[oldnvars+i];
1164  problem->bu[oldnvars+nvars+i] = problem->bu[oldnvars+i];
1165  }
1166 
1167  /* set bounds of new variable */
1168  for( i = 0; i < nvars; ++i )
1169  {
1170  problem->bl[oldnvars+i] = lbs[i];
1171  problem->bu[oldnvars+i] = ubs[i];
1172  }
1173  }
1174 
1175  /* gradients information is out of date now (objective gradient is stored in dense form) */
1176  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1177  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1178 
1179  /* Hessian information is out of date now (no new entries in Hessian, but also empty cols shows up in sparsity info) */
1180  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1181 
1182  return SCIP_OKAY;
1183 } /*lint !e715*/
1184 
1185 
1186 /** add constraints
1187  * quadratic coefficients: row oriented matrix for each constraint
1188  *
1189  * input:
1190  * - nlpi datastructure for solver interface
1191  * - problem datastructure for problem instance
1192  * - ncons number of added constraints
1193  * - lhss left hand sides of constraints
1194  * - rhss right hand sides of constraints
1195  * - nlininds number of linear coefficients for each constraint
1196  * may be NULL in case of no linear part
1197  * - lininds indices of variables for linear coefficients for each constraint
1198  * may be NULL in case of no linear part
1199  * - linvals values of linear coefficient for each constraint
1200  * may be NULL in case of no linear part
1201  * - nquadrows number of columns in matrix of quadratic part for each constraint
1202  * may be NULL in case of no quadratic part in any constraint
1203  * - quadrowidxs indices of variables for which a quadratic part is specified
1204  * may be NULL in case of no quadratic part in any constraint
1205  * - quadoffsets start index of each rows quadratic coefficients in quadinds[.] and quadvals[.]
1206  * indices are given w.r.t. quadrowidxs., i.e., quadoffsets[.][i] gives the start index of row quadrowidxs[.][i] in quadvals[.]
1207  * quadoffsets[.][nquadrows[.]] gives length of quadinds[.] and quadvals[.]
1208  * entry of array may be NULL in case of no quadratic part
1209  * may be NULL in case of no quadratic part in any constraint
1210  * - quadinds column indices w.r.t. quadrowidxs, i.e., quadrowidxs[quadinds[.][i]] gives the index of the variable corresponding
1211  * to entry i, entry of array may be NULL in case of no quadratic part
1212  * may be NULL in case of no quadratic part in any constraint
1213  * - quadvals coefficient values
1214  * entry of array may be NULL in case of no quadratic part
1215  * may be NULL in case of no quadratic part in any constraint
1216  * - exprvaridxs indices of variables in expression tree, maps variable indices in expression tree to indices in nlp
1217  * entry of array may be NULL in case of no expression tree
1218  * may be NULL in case of no expression tree in any constraint
1219  * - exprtrees expression tree for nonquadratic part of constraints
1220  * entry of array may be NULL in case of no nonquadratic part
1221  * may be NULL in case of no nonquadratic part in any constraint
1222  * - names of constraints, may be NULL or entries may be NULL
1223  */
1224 static
1225 SCIP_DECL_NLPIADDCONSTRAINTS( nlpiAddConstraintsFilterSQP )
1226 {
1227  SCIP_NLPIDATA* nlpidata;
1228  int oldnconss;
1229 
1230  assert(nlpi != NULL);
1231  assert(problem != NULL);
1232  assert(problem->oracle != NULL);
1233 
1234  nlpidata = SCIPnlpiGetData(nlpi);
1235  assert(nlpidata != NULL);
1236 
1237  oldnconss = SCIPnlpiOracleGetNConstraints(problem->oracle);
1238 
1239  SCIP_CALL( SCIPnlpiOracleAddConstraints(problem->oracle,
1240  ncons, lhss, rhss,
1241  nlininds, lininds, linvals,
1242  nquadelems, quadelems,
1243  exprvaridxs, exprtrees, names) );
1244 
1245  invalidateSolution(problem);
1246  problem->warmstart = FALSE;
1247 
1248  /* increase constraints-related arrays in problem, if necessary */
1249  if( SCIPnlpiOracleGetNConstraints(problem->oracle) > problem->conssize )
1250  {
1251  int newsize = calcGrowSize(SCIPnlpiOracleGetNConstraints(problem->oracle));
1252 
1253  if( problem->consdualvalues != NULL )
1254  {
1255  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->consdualvalues, problem->conssize, newsize) );
1256  }
1257 
1258  if( problem->c != NULL )
1259  {
1260  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->c, problem->conssize, newsize) );
1261  }
1262 
1263  if( problem->lam != NULL )
1264  {
1265  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->lam, problem->varssize + problem->conssize, problem->varssize + newsize) ); /*lint !e776*/
1266  }
1267 
1268  if( problem->bl != NULL )
1269  {
1270  assert(problem->bu != NULL);
1271  assert(problem->cstype != NULL);
1272  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->bl, problem->varssize + problem->conssize, problem->varssize + newsize) ); /*lint !e776*/
1273  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->bu, problem->varssize + problem->conssize, problem->varssize + newsize) ); /*lint !e776*/
1274  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->cstype, problem->conssize, newsize) );
1275  }
1276 
1277  if( problem->s != NULL )
1278  {
1279  SCIP_ALLOC( BMSreallocBlockMemoryArray(nlpidata->blkmem, &problem->s, problem->varssize + problem->conssize, problem->varssize + newsize) ); /*lint !e776*/
1280  }
1281 
1282  problem->conssize = newsize;
1283  }
1284 
1285  /* update constraint sides and type in FilterSQP data */
1286  if( problem->bl != NULL )
1287  {
1288  int i;
1289  int nvars;
1290 
1291  nvars = SCIPnlpiOracleGetNVars(problem->oracle);
1292 
1293  for( i = 0; i < ncons; ++i )
1294  {
1295  problem->bl[nvars+oldnconss+i] = lhss[i];
1296  problem->bu[nvars+oldnconss+i] = rhss[i];
1297  problem->cstype[oldnconss+i] = SCIPnlpiOracleGetConstraintDegree(problem->oracle, oldnconss+i) <= 1 ? 'L' : 'N';
1298  }
1299  }
1300 
1301  /* gradients information is out of date now */
1302  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1303  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1304 
1305  /* Hessian information is out of date now */
1306  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1307 
1308  return SCIP_OKAY;
1309 } /*lint !e715*/
1310 
1311 /** sets or overwrites objective, a minimization problem is expected
1312  * May change sparsity pattern.
1313  *
1314  * input:
1315  * - nlpi datastructure for solver interface
1316  * - problem datastructure for problem instance
1317  * - nlins number of linear variables
1318  * - lininds variable indices
1319  * may be NULL in case of no linear part
1320  * - linvals coefficient values
1321  * may be NULL in case of no linear part
1322  * - nquadcols number of columns in matrix of quadratic part
1323  * - quadcols indices of variables for which a quadratic part is specified
1324  * may be NULL in case of no quadratic part
1325  * - quadoffsets start index of each rows quadratic coefficients in quadinds and quadvals
1326  * quadoffsets[.][nquadcols] gives length of quadinds and quadvals
1327  * may be NULL in case of no quadratic part
1328  * - quadinds column indices
1329  * may be NULL in case of no quadratic part
1330  * - quadvals coefficient values
1331  * may be NULL in case of no quadratic part
1332  * - exprvaridxs indices of variables in expression tree, maps variable indices in expression tree to indices in nlp
1333  * may be NULL in case of no expression tree
1334  * - exprtree expression tree for nonquadratic part of objective function
1335  * may be NULL in case of no nonquadratic part
1336  * - constant objective value offset
1337  */
1338 static
1339 SCIP_DECL_NLPISETOBJECTIVE( nlpiSetObjectiveFilterSQP )
1340 {
1341  SCIP_NLPIDATA* nlpidata;
1342 
1343  assert(nlpi != NULL);
1344  assert(problem != NULL);
1345  assert(problem->oracle != NULL);
1346 
1347  nlpidata = SCIPnlpiGetData(nlpi);
1348  assert(nlpidata != NULL);
1349 
1350  SCIP_CALL( SCIPnlpiOracleSetObjective(problem->oracle,
1351  constant, nlins, lininds, linvals,
1352  nquadelems, quadelems,
1353  exprvaridxs, exprtree) );
1354 
1355  invalidateSolution(problem);
1356 
1357  /* gradients info (la,a) should still be ok, as objective gradient is stored in dense form */
1358 
1359  /* Hessian information is out of date now */
1360  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1361 
1362  return SCIP_OKAY;
1363 } /*lint !e715*/
1364 
1365 /** change variable bounds
1366  *
1367  * input:
1368  * - nlpi datastructure for solver interface
1369  * - problem datastructure for problem instance
1370  * - nvars number of variables to change bounds
1371  * - indices indices of variables to change bounds
1372  * - lbs new lower bounds
1373  * - ubs new upper bounds
1374  */
1375 static
1376 SCIP_DECL_NLPICHGVARBOUNDS( nlpiChgVarBoundsFilterSQP )
1377 {
1378  assert(nlpi != NULL);
1379  assert(problem != NULL);
1380  assert(problem->oracle != NULL);
1381 
1382  SCIP_CALL( SCIPnlpiOracleChgVarBounds(problem->oracle, nvars, indices, lbs, ubs) );
1383 
1384  invalidateSolution(problem);
1385 
1386  /* update bounds in FilterSQP data */
1387  if( problem->bl != NULL )
1388  {
1389  int i;
1390 
1391  for( i = 0; i < nvars; ++i )
1392  {
1393  problem->bl[indices[i]] = lbs[i];
1394  problem->bu[indices[i]] = ubs[i];
1395  }
1396  }
1397 
1398  return SCIP_OKAY;
1399 } /*lint !e715*/
1400 
1401 /** change constraint bounds
1402  *
1403  * input:
1404  * - nlpi datastructure for solver interface
1405  * - problem datastructure for problem instance
1406  * - nconss number of constraints to change sides
1407  * - indices indices of constraints to change sides
1408  * - lhss new left hand sides
1409  * - rhss new right hand sides
1410  */
1411 static
1412 SCIP_DECL_NLPICHGCONSSIDES( nlpiChgConsSidesFilterSQP )
1413 {
1414  assert(nlpi != NULL);
1415  assert(problem != NULL);
1416  assert(problem->oracle != NULL);
1417 
1418  SCIP_CALL( SCIPnlpiOracleChgConsSides(problem->oracle, nconss, indices, lhss, rhss) );
1419 
1420  invalidateSolution(problem);
1421 
1422  /* update constraint sense in FilterSQP data */
1423  if( problem->bl != NULL )
1424  {
1425  int i;
1426  int nvars;
1427 
1428  nvars = SCIPnlpiOracleGetNVars(problem->oracle);
1429 
1430  for( i = 0; i < nconss; ++i )
1431  {
1432  problem->bl[nvars+indices[i]] = lhss[i];
1433  problem->bu[nvars+indices[i]] = rhss[i];
1434  }
1435  }
1436 
1437  return SCIP_OKAY;
1438 } /*lint !e715*/
1439 
1440 /** delete a set of variables
1441  *
1442  * input:
1443  * - nlpi datastructure for solver interface
1444  * - problem datastructure for problem instance
1445  * - dstats deletion status of vars; 1 if var should be deleted, 0 if not
1446  *
1447  * output:
1448  * - dstats new position of var, -1 if var was deleted
1449  */
1450 static
1451 SCIP_DECL_NLPIDELVARSET( nlpiDelVarSetFilterSQP )
1452 {
1453  SCIP_NLPIDATA* nlpidata;
1454 
1455  assert(nlpi != NULL);
1456  assert(problem != NULL);
1457  assert(problem->oracle != NULL);
1458 
1459  nlpidata = SCIPnlpiGetData(nlpi);
1460  assert(nlpidata != NULL);
1461 
1462  SCIP_CALL( SCIPnlpiOracleDelVarSet(problem->oracle, dstats) );
1463 
1464  /* @TODO keep initguess and bl, bu for remaining variables? */
1465 
1466  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->initguess, problem->varssize);
1467  invalidateSolution(problem);
1468  problem->warmstart = FALSE;
1469 
1470  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->bl, problem->varssize + problem->conssize); /*lint !e776 */
1471  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->bu, problem->varssize + problem->conssize); /*lint !e776 */
1472  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->cstype, problem->conssize); /* because we assume that cstype is allocated iff bl is allocated */
1473 
1474  /* gradients information is out of date now (objective gradient is stored in dense form) */
1475  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1476  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1477 
1478  /* Hessian information is out of date now */
1479  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1480 
1481  return SCIP_OKAY;
1482 } /*lint !e715*/
1483 
1484 /** delete a set of constraints
1485  *
1486  * input:
1487  * - nlpi datastructure for solver interface
1488  * - problem datastructure for problem instance
1489  * - dstats deletion status of rows; 1 if row should be deleted, 0 if not
1490  *
1491  * output:
1492  * - dstats new position of row, -1 if row was deleted
1493  */
1494 static
1495 SCIP_DECL_NLPIDELCONSSET( nlpiDelConstraintSetFilterSQP )
1496 {
1497  SCIP_NLPIDATA* nlpidata;
1498 
1499  assert(nlpi != NULL);
1500  assert(problem != NULL);
1501  assert(problem->oracle != NULL);
1502 
1503  nlpidata = SCIPnlpiGetData(nlpi);
1504  assert(nlpidata != NULL);
1505 
1506  SCIP_CALL( SCIPnlpiOracleDelConsSet(problem->oracle, dstats) );
1507 
1508  invalidateSolution(problem);
1509  problem->warmstart = FALSE;
1510 
1511  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->bl, problem->varssize + problem->conssize); /*lint !e776 */
1512  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->bu, problem->varssize + problem->conssize); /*lint !e776 */
1513  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->cstype, problem->conssize); /* because we assume that cstype is allocated iff bl is allocated */
1514 
1515  /* gradients information is out of date now */
1516  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1517  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1518 
1519  /* Hessian information is out of date now */
1520  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1521 
1522  return SCIP_OKAY;
1523 } /*lint !e715*/
1524 
1525 /** changes (or adds) linear coefficients in a constraint or objective
1526  *
1527  * input:
1528  * - nlpi datastructure for solver interface
1529  * - problem datastructure for problem instance
1530  * - idx index of constraint or -1 for objective
1531  * - nvals number of values in linear constraint to change
1532  * - varidxs indices of variables which coefficient to change
1533  * - vals new values for coefficients
1534  */
1535 static
1536 SCIP_DECL_NLPICHGLINEARCOEFS( nlpiChgLinearCoefsFilterSQP )
1537 {
1538  SCIP_NLPIDATA* nlpidata;
1539 
1540  assert(nlpi != NULL);
1541  assert(problem != NULL);
1542  assert(problem->oracle != NULL);
1543 
1544  nlpidata = SCIPnlpiGetData(nlpi);
1545  assert(nlpidata != NULL);
1546 
1547  SCIP_CALL( SCIPnlpiOracleChgLinearCoefs(problem->oracle, idx, nvals, varidxs, vals) );
1548 
1549  invalidateSolution(problem);
1550 
1551  /* gradients information (la,a) may have changed if elements were added or removed
1552  * (we only care that sparsity doesn't change, not about actual values in a)
1553  * TODO free only if coefficients were added or removed (SCIPnlpiOracleChgLinearCoefs() could give feedback)
1554  */
1555  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1556  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1557 
1558  /* Hessian information should still be ok */
1559 
1560  return SCIP_OKAY;
1561 } /*lint !e715*/
1562 
1563 /** changes (or adds) coefficients in the quadratic part of a constraint or objective
1564  *
1565  * input:
1566  * - nlpi datastructure for solver interface
1567  * - problem datastructure for problem instance
1568  * - idx index of constraint or -1 for objective
1569  * - nentries number of entries in quadratic matrix to change
1570  * - rows row indices of entries in quadratic matrix where values should be changed
1571  * - cols column indices of entries in quadratic matrix where values should be changed
1572  * - values new values for entries in quadratic matrix
1573  */
1574 static
1575 SCIP_DECL_NLPICHGQUADCOEFS( nlpiChgQuadraticCoefsFilterSQP )
1576 {
1577  SCIP_NLPIDATA* nlpidata;
1578 
1579  assert(nlpi != NULL);
1580  assert(problem != NULL);
1581  assert(problem->oracle != NULL);
1582 
1583  nlpidata = SCIPnlpiGetData(nlpi);
1584  assert(nlpidata != NULL);
1585 
1586  SCIP_CALL( SCIPnlpiOracleChgQuadCoefs(problem->oracle, idx, nquadelems, quadelems) );
1587 
1588  invalidateSolution(problem);
1589 
1590  /* update constraint linearity in FilterSQP data, as we might have changed from linear to nonlinear now */
1591  if( problem->cstype != NULL && idx >= 0 )
1592  problem->cstype[idx] = (SCIPnlpiOracleGetConstraintDegree(problem->oracle, idx) <= 1 ? 'L' : 'N');
1593 
1594  /* gradients information (la,a) may have changed if elements were added or removed
1595  * (we only care that sparsity doesn't change, not about actual values in a)
1596  * TODO free only if coefficients were added or removed (SCIPnlpiOracleChgLinearCoefs() could give feedback)
1597  */
1598  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1599  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1600 
1601  /* Hessian sparsity may have changed if elements were added or removed
1602  * TODO free only if coefficients were added or removed (SCIPnlpiOracleChgLinearCoefs() could give feedback)
1603  */
1604  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1605 
1606  return SCIP_OKAY;
1607 } /*lint !e715*/
1608 
1609 /** replaces the expression tree of a constraint or objective
1610  *
1611  * input:
1612  * - nlpi datastructure for solver interface
1613  * - problem datastructure for problem instance
1614  * - idxcons index of constraint or -1 for objective
1615  * - exprvaridxs indices of variables in expression tree, maps variable indices in expression tree to indices in nlp, or NULL
1616  * - exprtree new expression tree for constraint or objective, or NULL to only remove previous tree
1617  */
1618 static
1619 SCIP_DECL_NLPICHGEXPRTREE( nlpiChgExprtreeFilterSQP )
1620 {
1621  SCIP_NLPIDATA* nlpidata;
1622 
1623  assert(nlpi != NULL);
1624  assert(problem != NULL);
1625  assert(problem->oracle != NULL);
1626 
1627  nlpidata = SCIPnlpiGetData(nlpi);
1628  assert(nlpidata != NULL);
1629 
1630  SCIP_CALL( SCIPnlpiOracleChgExprtree(problem->oracle, idxcons, exprvaridxs, exprtree) );
1631 
1632  invalidateSolution(problem);
1633 
1634  /* update constraint linearity in FilterSQP data, as we might have changed from linear to nonlinear now */
1635  if( problem->cstype != NULL && idxcons >= 0 )
1636  problem->cstype[idxcons] = (SCIPnlpiOracleGetConstraintDegree(problem->oracle, idxcons) <= 1 ? 'L' : 'N');
1637 
1638  /* gradients information (la,a) may have changed */
1639  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1640  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1641 
1642  /* Hessian information may have changed */
1643  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1644 
1645  return SCIP_OKAY;
1646 } /*lint !e715*/
1647 
1648 /** change one coefficient in the nonlinear part
1649  *
1650  * input:
1651  * - nlpi datastructure for solver interface
1652  * - problem datastructure for problem instance
1653  * - idxcons index of constraint or -1 for objective
1654  * - idxparam index of parameter
1655  * - value new value for nonlinear parameter
1656  *
1657  * return: Error if parameter does not exist
1658  */
1659 static
1660 SCIP_DECL_NLPICHGNONLINCOEF( nlpiChgNonlinCoefFilterSQP )
1661 {
1662  SCIP_NLPIDATA* nlpidata;
1663 
1664  assert(nlpi != NULL);
1665  assert(problem != NULL);
1666  assert(problem->oracle != NULL);
1667 
1668  nlpidata = SCIPnlpiGetData(nlpi);
1669  assert(nlpidata != NULL);
1670 
1671  SCIP_CALL( SCIPnlpiOracleChgExprParam(problem->oracle, idxcons, idxparam, value) );
1672 
1673  invalidateSolution(problem);
1674 
1675  /* gradients information (la,a) may have changed (?) */
1676  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->a, problem->la != NULL ? problem->la[0]-1 : 0);
1677  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->la, problem->lasize);
1678 
1679  /* Hessian information may have changed (?) */
1680  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->hessiannz, problem->hessiannzsize);
1681 
1682  return SCIP_OKAY;
1683 } /*lint !e715*/
1684 
1685 /** change the constant offset in the objective
1686  *
1687  * input:
1688  * - nlpi datastructure for solver interface
1689  * - problem datastructure for problem instance
1690  * - objconstant new value for objective constant
1691  */
1692 static
1693 SCIP_DECL_NLPICHGOBJCONSTANT( nlpiChgObjConstantFilterSQP )
1694 {
1695  assert(nlpi != NULL);
1696  assert(problem != NULL);
1697  assert(problem->oracle != NULL);
1698 
1699  SCIP_CALL( SCIPnlpiOracleChgObjConstant(problem->oracle, objconstant) );
1700 
1701  invalidateSolution(problem);
1702 
1703  return SCIP_OKAY;
1704 } /*lint !e715*/
1705 
1706 /** sets initial guess for primal variables
1707  *
1708  * input:
1709  * - nlpi datastructure for solver interface
1710  * - problem datastructure for problem instance
1711  * - primalvalues initial primal values for variables, or NULL to clear previous values
1712  * - consdualvalues initial dual values for constraints, or NULL to clear previous values
1713  * - varlbdualvalues initial dual values for variable lower bounds, or NULL to clear previous values
1714  * - varubdualvalues initial dual values for variable upper bounds, or NULL to clear previous values
1715  */
1716 static
1717 SCIP_DECL_NLPISETINITIALGUESS( nlpiSetInitialGuessFilterSQP )
1718 {
1719  SCIP_NLPIDATA* nlpidata;
1720 
1721  assert(nlpi != NULL);
1722  assert(problem != NULL);
1723  assert(problem->oracle != NULL);
1724 
1725  nlpidata = SCIPnlpiGetData(nlpi);
1726  assert(nlpidata != NULL);
1727 
1728  if( primalvalues != NULL )
1729  {
1730  if( problem->initguess == NULL )
1731  {
1732  SCIP_ALLOC( BMSallocBlockMemoryArray(nlpidata->blkmem, &problem->initguess, problem->varssize) );
1733  }
1734  assert(SCIPnlpiOracleGetNVars(problem->oracle) <= problem->varssize);
1735  BMScopyMemoryArray(problem->initguess, primalvalues, SCIPnlpiOracleGetNVars(problem->oracle));
1736  }
1737  else
1738  {
1739  BMSfreeBlockMemoryArrayNull(nlpidata->blkmem, &problem->initguess, problem->varssize);
1740  }
1741 
1742  return SCIP_OKAY;
1743 } /*lint !e715*/
1744 
1745 /** tries to solve NLP
1746  *
1747  * input:
1748  * - nlpi datastructure for solver interface
1749  * - problem datastructure for problem instance
1750  */
1751 static
1752 SCIP_DECL_NLPISOLVE( nlpiSolveFilterSQP )
1753 {
1754  SCIP_NLPIDATA* data;
1755  SCIP_Bool success;
1756  fint n;
1757  fint m;
1758  fint kmax;
1759  fint maxa;
1760  fint maxf;
1761  fint mlp;
1762  fint lh1;
1763  fint nout;
1764  fint ifail;
1765  fint maxiter;
1766  real rho;
1767  real f;
1768  real* user;
1769  fint* iuser;
1770  ftnlen cstype_len = 1;
1771  fint minmxwk;
1772  fint minmxiwk;
1773  int nruns;
1774  int i;
1775 
1776  data = SCIPnlpiGetData(nlpi);
1777  assert(data != NULL);
1778 
1779  /* start measuring time */
1780  data->starttime = gettime();
1781 
1782  /* if fromscratch parameter is set, then we will not warmstart */
1783  if( problem->fromscratch )
1784  problem->warmstart = FALSE;
1785 
1786  n = SCIPnlpiOracleGetNVars(problem->oracle);
1787  m = SCIPnlpiOracleGetNConstraints(problem->oracle);
1788  kmax = n; /* maximal nullspace dimension */
1789  maxf = 100; /* maximal filter length */
1790  mlp = 100; /* maximum level of degeneracy */
1791 
1792  /* TODO eventually, the output should be redirected to the message handler,
1793  * but even to just redirect to some other file, we would have to open the output-unit in Fortran
1794  */
1795  nout = 6; /* output to stdout for now */
1796  ifail = problem->warmstart ? -1 : 0; /* -1 for warmstart, otherwise 0 */
1797  rho = 10.0; /* initial trust-region radius */
1798 
1799  user = (real*)data;
1800  iuser = (fint*)problem;
1801  if( problem->warmstart ) /* if warmstart, then need to keep istat[0] */
1802  memset(problem->istat+1, 0, sizeof(problem->istat)-sizeof(*problem->istat));
1803  else
1804  memset(problem->istat, 0, sizeof(problem->istat));
1805  memset(problem->rstat, 0, sizeof(problem->rstat));
1806  problem->niterations = 0;
1807 
1808  if( problem->x == NULL )
1809  {
1810  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->x, problem->varssize) );
1811  }
1812  if( problem->c == NULL )
1813  {
1814  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->c, problem->conssize) );
1815  }
1816  if( problem->lam == NULL )
1817  {
1818  SCIP_ALLOC( BMSallocClearBlockMemoryArray(data->blkmem, &problem->lam, problem->varssize + problem->conssize) ); /*lint !e776 */
1819  }
1820  else
1821  {
1822  BMSclearMemoryArray(problem->lam, problem->varssize + problem->conssize);
1823  }
1824  if( problem->s == NULL )
1825  {
1826  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->s, problem->varssize + problem->conssize) );
1827  }
1828 
1829  if( problem->la == NULL )
1830  {
1831  /* allocate la, a and initialize la for Objective Gradient and Jacobian */
1832  SCIP_CALL( setupGradients(data->blkmem, problem->oracle, &problem->la, &problem->lasize, &problem->a) );
1833  }
1834  /* maximal number entries in a = nvars+nnz */
1835  maxa = problem->la[0]-1;
1836 
1837  if( problem->hessiannz == NULL )
1838  {
1839  /* allocate and initialize problem->hessiannz for Hessian */
1840  SCIP_CALL( setupHessian(data->blkmem, problem->oracle, &problem->hessiannz, &problem->hessiannzsize) );
1841  }
1842 
1843  /* setup variable bounds, constraint sides, and constraint types */
1844  if( problem->bl == NULL )
1845  {
1846  assert(problem->bu == NULL);
1847  assert(problem->cstype == NULL);
1848 
1849  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->bl, problem->varssize + problem->conssize) );
1850  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->bu, problem->varssize + problem->conssize) );
1851  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->cstype, problem->conssize) );
1852 
1853  BMScopyMemoryArray(problem->bl, SCIPnlpiOracleGetVarLbs(problem->oracle), n);
1854  BMScopyMemoryArray(problem->bu, SCIPnlpiOracleGetVarUbs(problem->oracle), n);
1855  for( i = 0; i < m; ++i )
1856  {
1857  problem->bl[n+i] = SCIPnlpiOracleGetConstraintLhs(problem->oracle, i);
1858  problem->bu[n+i] = SCIPnlpiOracleGetConstraintRhs(problem->oracle, i);
1859  problem->cstype[i] = SCIPnlpiOracleGetConstraintDegree(problem->oracle, i) <= 1 ? 'L' : 'N';
1860  }
1861  }
1862 
1863  /* buffer for evaluation results (used in setupStart already) */
1864  if( problem->evalbufsize < MAX3(n, problem->hessiannz[0], maxa) )
1865  {
1866  int newsize = calcGrowSize(MAX3(n, problem->hessiannz[0], maxa));
1867  SCIP_ALLOC( BMSreallocBlockMemoryArray(data->blkmem, &problem->evalbuffer, problem->evalbufsize, newsize) );
1868  problem->evalbufsize = newsize;
1869  }
1870 
1871  /* setup starting point */
1872  SCIP_CALL( setupStart(data, problem, problem->x, &success) );
1873  if( !success )
1874  {
1875  /* FilterSQP would crash if starting point cannot be evaluated, so give up */
1876  SCIP_CALL( processSolveOutcome(data, problem, 7, NULL, NULL) );
1877  return SCIP_OKAY;
1878  }
1879 
1880  /* setup workspace */
1881  /* initial guess of real workspace size */
1882  /* FilterSQP manual: mxwk = 21*n + 8*m + mlp + 8*maxf + kmax*(kmax+9)/2 + nprof, with nprof = 20*n as a good guess */
1883  /* Bonmin: mxwk = 21*n + 8*m + mlp + 8*maxf + lh1 + kmax*(kmax+9)/2 + mxwk0,
1884  * with lh1 = nnz_h+8+2*n+m and mxwk0 = 2000000 (parameter) */
1885  lh1 = problem->hessiannz[0]-1 + 8 + 2*n + m;
1886  minmxwk = 21*n + 8*m + mlp + 8*maxf + lh1 + kmax*(kmax+9)/2 + MAX(20*n,2000);
1887  if( problem->ws == NULL )
1888  {
1889  problem->mxwk = minmxwk;
1890  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->ws, problem->mxwk) );
1891  }
1892  else if( problem->mxwk < minmxwk )
1893  {
1894  int newsize = calcGrowSize(minmxwk);
1895  SCIP_ALLOC( BMSreallocBlockMemoryArray(data->blkmem, &problem->ws, problem->mxwk, newsize) );
1896  problem->mxwk = newsize;
1897  }
1898 
1899  /* initial guess of integer workspace size */
1900  /* FilterSQP manual: mxiwk = 13*n + 4*m + mlp + 100 + kmax */
1901  /* Bonmin: mxiwk = 13*n + 4*m + mlp + lh1 + kmax + 113 + mxiwk0, with mxiwk0 = 500000 (parameter) */
1902  minmxiwk = 13*n + 4*m + mlp + lh1 + 100 + kmax + 113 + MAX(5*n,5000);
1903  if( problem->lws == NULL )
1904  {
1905  assert(!problem->warmstart);
1906 
1907  problem->mxiwk = minmxiwk;
1908  SCIP_ALLOC( BMSallocBlockMemoryArray(data->blkmem, &problem->lws, problem->mxiwk) );
1909  }
1910  else if( problem->mxiwk < minmxiwk && !problem->warmstart ) /* if warmstart, then lws should remain untouched (n and m didn't change anyway) */
1911  {
1912  int newsize = calcGrowSize(minmxiwk);
1913  SCIP_ALLOC( BMSreallocBlockMemoryArray(data->blkmem, &problem->lws, problem->mxiwk, newsize) );
1914  problem->mxiwk = newsize;
1915  }
1916  /* in case of some evalerrors, not clearing ws could lead to valgrind warnings about use of uninitialized memory */
1917  memset(problem->ws, 0, (size_t) problem->mxwk * sizeof(real)); /*lint !e571*/
1918 
1919  /* from here on we are not thread-safe: if intended for multithread use, then protect filtersqp call with mutex
1920  * NOTE: we need to make sure that we do not return from nlpiSolve before unlocking the mutex
1921  */
1922 #ifndef NPARASCIP
1923  pthread_mutex_lock(&filtersqpmutex);
1924 #endif
1925 
1926  /* initialize global variables from filtersqp */
1927  /* FilterSQP eps is tolerance for both feasibility and optimality, and also for trust-region radius, etc. */
1928  F77_FUNC(nlp_eps_inf,NLP_EPS_INF).eps = MIN(problem->feastol, problem->opttol * OPTTOLFACTOR);
1929  F77_FUNC(nlp_eps_inf,NLP_EPS_INF).infty = SCIPnlpiOracleGetInfinity(problem->oracle);
1930  F77_FUNC(ubdc,UBDC).ubd = 100.0;
1931  F77_FUNC(ubdc,UBDC).tt = 1.25;
1932  F77_FUNC(scalec,SCALEC).scale_mode = 0;
1933 
1934  for( nruns = 1; ; ++nruns )
1935  {
1936  maxiter = problem->maxiter - problem->niterations;
1937 
1938  F77_FUNC(filtersqp,FILTERSQP)(
1939  &n, &m, &kmax, &maxa,
1940  &maxf, &mlp, &problem->mxwk, &problem->mxiwk,
1941  &problem->iprint, &nout, &ifail, &rho,
1942  problem->x, problem->c, &f, &problem->fmin, problem->bl,
1943  problem->bu, problem->s, problem->a, problem->la, problem->ws,
1944  problem->lws, problem->lam, problem->cstype, user,
1945  iuser, &maxiter, problem->istat,
1946  problem->rstat, cstype_len);
1947 
1948  problem->niterations += problem->istat[1];
1949 
1950  assert(ifail <= 10);
1951  /* if ifail >= 8 (probably the workspace was too small), then retry with larger workspace
1952  * if ifail == 0 (local optimal), but absolute violation of KKT too large, then retry with small eps
1953  */
1954  if( ifail < 8 && (ifail != 0 || problem->rstat[0] <= problem->opttol) )
1955  break;
1956 
1957  if( problem->iprint > 0 )
1958  {
1959  SCIPmessagePrintInfo(data->messagehdlr, "FilterSQP terminated with status %d in run %d, absolute KKT violation is %g\n", ifail, nruns, problem->rstat[0]);
1960  }
1961 
1962  /* if iteration or time limit exceeded, then don't retry */
1963  if( problem->niterations >= problem->maxiter || timelimitreached(data, problem) )
1964  {
1965  if( problem->iprint > 0 )
1966  {
1967  SCIPmessagePrintInfo(data->messagehdlr, "Time or iteration limit reached, not retrying\n");
1968  }
1969  break;
1970  }
1971 
1972  /* if maximal number of runs reached, then stop */
1973  if( nruns >= MAXNRUNS )
1974  {
1975  if( problem->iprint > 0 )
1976  {
1977  SCIPmessagePrintInfo(data->messagehdlr, "Run limit reached, not retrying\n");
1978  }
1979  break;
1980  }
1981 
1982  if( ifail == 0 )
1983  {
1984  SCIP_Real epsfactor;
1985 
1986  if( F77_FUNC(nlp_eps_inf,NLP_EPS_INF).eps <= MINEPS )
1987  {
1988  if( problem->iprint > 0 )
1989  {
1990  SCIPmessagePrintInfo(data->messagehdlr, "Already reached minimal epsilon, not retrying\n");
1991  }
1992  break;
1993  }
1994 
1995  epsfactor = problem->opttol / problem->rstat[0];
1996  assert(epsfactor < 1.0); /* because of the if's above */
1997  epsfactor *= OPTTOLFACTOR;
1998 
1999  F77_FUNC(nlp_eps_inf,NLP_EPS_INF).eps = MAX(MINEPS, F77_FUNC(nlp_eps_inf,NLP_EPS_INF).eps * epsfactor);
2000  if( problem->iprint > 0 )
2001  {
2002  SCIPmessagePrintInfo(data->messagehdlr, "Continue with eps = %g\n", F77_FUNC(nlp_eps_inf,NLP_EPS_INF).eps);
2003  }
2004  ifail = -1; /* do warmstart */
2005 
2006  continue;
2007  }
2008 
2009  /* increase real workspace, if ifail = 9 (real workspace too small) or ifail = 8 (unexpected ifail from QP solver, often also when workspace too small) */
2010  if( ifail == 8 || ifail == 9 )
2011  {
2012  int newsize = calcGrowSize(WORKSPACEGROWTHFACTOR*problem->mxwk);
2013  if( BMSreallocBlockMemoryArray(data->blkmem, &problem->ws, problem->mxwk, newsize) == NULL )
2014  {
2015  /* realloc failed: give up NLP solve */
2016  problem->mxwk = 0;
2017  break;
2018  }
2019  problem->mxwk = newsize;
2020  }
2021 
2022  /* increase integer workspace, if ifail = 10 (integer workspace too small) or ifail = 8 (unexpected ifail from QP solver, often also when workspace too small) */
2023  if( ifail == 8 || ifail == 10 )
2024  {
2025  int newsize = calcGrowSize(WORKSPACEGROWTHFACTOR*problem->mxiwk);
2026  if( BMSreallocBlockMemoryArray(data->blkmem, &problem->lws, problem->mxiwk, newsize) == NULL )
2027  {
2028  /* realloc failed: give up NLP solve */
2029  problem->mxiwk = 0;
2030  break;
2031  }
2032  problem->mxiwk = newsize;
2033 
2034  /* better don't try warmstart for the next trial; warmstart requires that lws is untouched, does extending count as touching? */
2035  ifail = 0;
2036  }
2037 
2038  /* reset starting point, in case it was overwritten by failed solve (return can only be SCIP_OKAY, because randnumgen must exist already)
2039  * NOTE: If initguess is NULL (no user-given starting point), then this will result in a slightly different starting point as in the previous setupStart() call (random numbers)
2040  * However, as no point was given, it shouldn't matter which point we actually start from.
2041  */
2042  (void) setupStart(data, problem, problem->x, &success);
2043  assert(success);
2044  }
2045 
2046 #ifndef NPARASCIP
2047  pthread_mutex_unlock(&filtersqpmutex);
2048 #endif
2049 
2050  SCIP_CALL( processSolveOutcome(data, problem, ifail, problem->x, problem->lam) );
2051 
2052  return SCIP_OKAY; /*lint !e527*/
2053 } /*lint !e715*/
2054 
2055 /** gives solution status
2056  *
2057  * input:
2058  * - nlpi datastructure for solver interface
2059  * - problem datastructure for problem instance
2060  *
2061  * return: Solution Status
2062  */
2063 static
2064 SCIP_DECL_NLPIGETSOLSTAT( nlpiGetSolstatFilterSQP )
2065 {
2066  assert(problem != NULL);
2067 
2068  return problem->solstat;
2069 } /*lint !e715*/
2070 
2071 /** gives termination reason
2072  *
2073  * input:
2074  * - nlpi datastructure for solver interface
2075  * - problem datastructure for problem instance
2076  *
2077  * return: Termination Status
2078  */
2079 static
2080 SCIP_DECL_NLPIGETTERMSTAT( nlpiGetTermstatFilterSQP )
2081 {
2082  assert(problem != NULL);
2083 
2084  return problem->termstat;
2085 } /*lint !e715*/
2086 
2087 /** gives primal and dual solution values
2088  *
2089  * solver can return NULL in dual values if not available
2090  * but if solver provides dual values for one side of variable bounds, then it must also provide those for the other side
2091  *
2092  * for a ranged constraint, the dual variable is positive if the right hand side is active and negative if the left hand side is active
2093  *
2094  * input:
2095  * - nlpi datastructure for solver interface
2096  * - problem datastructure for problem instance
2097  * - primalvalues buffer to store pointer to array to primal values, or NULL if not needed
2098  * - consdualvalues buffer to store pointer to array to dual values of constraints, or NULL if not needed
2099  * - varlbdualvalues buffer to store pointer to array to dual values of variable lower bounds, or NULL if not needed
2100  * - varubdualvalues buffer to store pointer to array to dual values of variable lower bounds, or NULL if not needed
2101  * - objval buffer store the objective value, or NULL if not needed
2102  */
2103 static
2104 SCIP_DECL_NLPIGETSOLUTION( nlpiGetSolutionFilterSQP )
2105 {
2106  assert(problem != NULL);
2107 
2108  if( primalvalues != NULL )
2109  {
2110  assert(problem->primalvalues != NULL);
2111 
2112  *primalvalues = problem->primalvalues;
2113  }
2114 
2115  if( consdualvalues != NULL )
2116  {
2117  assert(problem->consdualvalues != NULL);
2118 
2119  *consdualvalues = problem->consdualvalues;
2120  }
2121 
2122  if( varlbdualvalues != NULL )
2123  {
2124  assert(problem->varlbdualvalues != NULL);
2125 
2126  *varlbdualvalues = problem->varlbdualvalues;
2127  }
2128 
2129  if( varubdualvalues != NULL )
2130  {
2131  assert(problem->varubdualvalues != NULL);
2132 
2133  *varubdualvalues = problem->varubdualvalues;
2134  }
2135 
2136  if( objval != NULL )
2137  {
2138  if( problem->primalvalues != NULL )
2139  {
2140  /* TODO store last solution value instead of reevaluating the objective function */
2141  SCIP_CALL( SCIPnlpiOracleEvalObjectiveValue(problem->oracle, problem->primalvalues, objval) );
2142  }
2143  else
2144  *objval = SCIP_INVALID;
2145  }
2146 
2147  return SCIP_OKAY; /*lint !e527*/
2148 } /*lint !e715*/
2149 
2150 /** gives solve statistics
2151  *
2152  * input:
2153  * - nlpi datastructure for solver interface
2154  * - problem datastructure for problem instance
2155  * - statistics pointer to store statistics
2156  *
2157  * output:
2158  * - statistics solve statistics
2159  */
2160 static
2161 SCIP_DECL_NLPIGETSTATISTICS( nlpiGetStatisticsFilterSQP )
2162 {
2163  assert(problem != NULL);
2164 
2165  SCIPnlpStatisticsSetNIterations(statistics, problem->niterations);
2166  SCIPnlpStatisticsSetTotalTime(statistics, problem->solvetime);
2167 
2168  return SCIP_OKAY; /*lint !e527*/
2169 } /*lint !e715*/
2170 
2171 /** gives required size of a buffer to store a warmstart object
2172  *
2173  * input:
2174  * - nlpi datastructure for solver interface
2175  * - problem datastructure for problem instance
2176  * - size pointer to store required size for warmstart buffer
2177  *
2178  * output:
2179  * - size required size for warmstart buffer
2180  */
2181 static
2182 SCIP_DECL_NLPIGETWARMSTARTSIZE( nlpiGetWarmstartSizeFilterSQP )
2183 {
2184  SCIPerrorMessage("method of filtersqp nonlinear solver is not implemented\n");
2185  SCIPABORT();
2186 
2187  return SCIP_OKAY; /*lint !e527*/
2188 } /*lint !e715*/
2189 
2190 /** stores warmstart information in buffer
2191  *
2192  * required size of buffer should have been obtained by SCIPnlpiGetWarmstartSize before
2193  *
2194  * input:
2195  * - nlpi datastructure for solver interface
2196  * - problem datastructure for problem instance
2197  * - buffer memory to store warmstart information
2198  *
2199  * output:
2200  * - buffer warmstart information in solver specific data structure
2201  */
2202 static
2203 SCIP_DECL_NLPIGETWARMSTARTMEMO( nlpiGetWarmstartMemoFilterSQP )
2204 {
2205  SCIPerrorMessage("method of filtersqp nonlinear solver is not implemented\n");
2206  SCIPABORT();
2207 
2208  return SCIP_OKAY; /*lint !e527*/
2209 } /*lint !e715*/
2210 
2211 /** sets warmstart information in solver
2212  *
2213  * write warmstart to buffer
2214  *
2215  * input:
2216  * - nlpi datastructure for solver interface
2217  * - problem datastructure for problem instance
2218  * - buffer warmstart information
2219  */
2220 static
2221 SCIP_DECL_NLPISETWARMSTARTMEMO( nlpiSetWarmstartMemoFilterSQP )
2222 {
2223  SCIPerrorMessage("method of filtersqp nonlinear solver is not implemented\n");
2224  SCIPABORT();
2225 
2226  return SCIP_OKAY; /*lint !e527*/
2227 } /*lint !e715*/
2228 
2229 /** gets integer parameter of NLP
2230  *
2231  * input:
2232  * - nlpi NLP interface structure
2233  * - problem datastructure for problem instance
2234  * - type parameter number
2235  * - ival pointer to store the parameter value
2236  *
2237  * output:
2238  * - ival parameter value
2239  */
2240 static
2241 SCIP_DECL_NLPIGETINTPAR( nlpiGetIntParFilterSQP )
2242 {
2243  assert(nlpi != NULL);
2244  assert(ival != NULL);
2245  assert(problem != NULL);
2246 
2247  switch( type )
2248  {
2250  {
2251  *ival = problem->fromscratch ? 1 : 0;
2252  break;
2253  }
2254 
2255  case SCIP_NLPPAR_VERBLEVEL:
2256  {
2257  *ival = problem->iprint;
2258  break;
2259  }
2260 
2261  case SCIP_NLPPAR_FEASTOL:
2262  {
2263  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
2264  return SCIP_PARAMETERWRONGTYPE;
2265  }
2266 
2267  case SCIP_NLPPAR_RELOBJTOL:
2268  {
2269  SCIPerrorMessage("relative objective tolerance parameter is of type real.\n");
2270  return SCIP_PARAMETERWRONGTYPE;
2271  }
2272 
2273  case SCIP_NLPPAR_LOBJLIM:
2274  {
2275  SCIPerrorMessage("objective limit parameter is of type real.\n");
2276  return SCIP_PARAMETERWRONGTYPE;
2277  }
2278 
2279  case SCIP_NLPPAR_INFINITY:
2280  {
2281  SCIPerrorMessage("infinity parameter is of type real.\n");
2282  return SCIP_PARAMETERWRONGTYPE;
2283  }
2284 
2285  case SCIP_NLPPAR_ITLIM:
2286  {
2287  *ival = problem->maxiter;
2288  break;
2289  }
2290 
2291  case SCIP_NLPPAR_TILIM:
2292  {
2293  SCIPerrorMessage("time limit parameter is of type real.\n");
2294  return SCIP_PARAMETERWRONGTYPE;
2295  }
2296 
2297  case SCIP_NLPPAR_OPTFILE:
2298  {
2299  SCIPerrorMessage("optfile parameter is of type string.\n");
2300  return SCIP_PARAMETERWRONGTYPE;
2301  }
2302 
2303  case SCIP_NLPPAR_FASTFAIL:
2304  {
2305  *ival = 0;
2306  break;
2307  }
2308 
2309  default:
2310  {
2311  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
2312  return SCIP_PARAMETERUNKNOWN;
2313  }
2314  }
2315 
2316  return SCIP_OKAY;
2317 } /*lint !e715*/
2318 
2319 /** sets integer parameter of NLP
2320  *
2321  * input:
2322  * - nlpi NLP interface structure
2323  * - problem datastructure for problem instance
2324  * - type parameter number
2325  * - ival parameter value
2326  */
2327 static
2328 SCIP_DECL_NLPISETINTPAR( nlpiSetIntParFilterSQP )
2329 {
2330  assert(nlpi != NULL);
2331  assert(problem != NULL);
2332 
2333  switch( type )
2334  {
2336  {
2337  if( ival == 0 || ival == 1 )
2338  {
2339  problem->fromscratch = (SCIP_Bool)ival;
2340  }
2341  else
2342  {
2343  SCIPerrorMessage("Value %d for parameter from scratch out of range {0, 1}\n", ival);
2344  return SCIP_PARAMETERWRONGVAL;
2345  }
2346  break;
2347  }
2348 
2349  case SCIP_NLPPAR_VERBLEVEL:
2350  {
2351  if( ival >= 0 )
2352  {
2353  problem->iprint = ival;
2354  }
2355  else
2356  {
2357  SCIPerrorMessage("Value %d for parameter from verbosity level out of range\n", ival);
2358  return SCIP_PARAMETERWRONGVAL;
2359  }
2360  break;
2361  }
2362 
2363  case SCIP_NLPPAR_FEASTOL:
2364  {
2365  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
2366  return SCIP_PARAMETERWRONGTYPE;
2367  }
2368 
2369  case SCIP_NLPPAR_RELOBJTOL:
2370  {
2371  SCIPerrorMessage("relative objective tolerance parameter is of type real.\n");
2372  return SCIP_PARAMETERWRONGTYPE;
2373  }
2374 
2375  case SCIP_NLPPAR_LOBJLIM:
2376  {
2377  SCIPerrorMessage("objective limit parameter is of type real.\n");
2378  return SCIP_PARAMETERWRONGTYPE;
2379  }
2380 
2381  case SCIP_NLPPAR_INFINITY:
2382  {
2383  SCIPerrorMessage("infinity parameter is of type real.\n");
2384  return SCIP_PARAMETERWRONGTYPE;
2385  }
2386 
2387  case SCIP_NLPPAR_ITLIM:
2388  {
2389  if( ival >= 0 )
2390  {
2391  problem->maxiter = ival;
2392  }
2393  else
2394  {
2395  SCIPerrorMessage("Value %d for parameter iteration limit is negative\n", ival);
2396  return SCIP_PARAMETERWRONGVAL;
2397  }
2398  break;
2399  }
2400 
2401  case SCIP_NLPPAR_TILIM:
2402  {
2403  SCIPerrorMessage("time limit parameter is of type real.\n");
2404  return SCIP_PARAMETERWRONGTYPE;
2405  }
2406 
2407  case SCIP_NLPPAR_OPTFILE:
2408  {
2409  SCIPerrorMessage("optfile parameter is of type string.\n");
2410  return SCIP_PARAMETERWRONGTYPE;
2411  }
2412 
2413  case SCIP_NLPPAR_FASTFAIL:
2414  {
2415  if( ival == 0 || ival == 1 )
2416  {
2417  SCIPdebugMessage("fast fail parameter not supported by FilterSQP interface yet. Ignored.\n");
2418  }
2419  else
2420  {
2421  SCIPerrorMessage("Value %d for parameter fastfail out of range {0, 1}\n", ival);
2422  return SCIP_PARAMETERWRONGVAL;
2423  }
2424  break;
2425  }
2426 
2427  default:
2428  {
2429  SCIPerrorMessage("Parameter %d not known to FilterSQP interface.\n", type);
2430  return SCIP_PARAMETERUNKNOWN;
2431  }
2432  }
2433 
2434  return SCIP_OKAY;
2435 } /*lint !e715*/
2436 
2437 /** gets floating point parameter of NLP
2438  *
2439  * input:
2440  * - nlpi NLP interface structure
2441  * - problem datastructure for problem instance, can be NULL only if type == SCIP_NLPPAR_INFINITY
2442  * - type parameter number
2443  * - dval pointer to store the parameter value
2444  *
2445  * output:
2446  * - dval parameter value
2447  */
2448 static
2449 SCIP_DECL_NLPIGETREALPAR( nlpiGetRealParFilterSQP )
2450 {
2451  assert(nlpi != NULL);
2452  assert(dval != NULL);
2453  assert(type == SCIP_NLPPAR_INFINITY || problem != NULL);
2454 
2455  switch( type )
2456  {
2458  {
2459  SCIPerrorMessage("from scratch parameter is of type int.\n");
2460  return SCIP_PARAMETERWRONGTYPE;
2461  }
2462 
2463  case SCIP_NLPPAR_VERBLEVEL:
2464  {
2465  SCIPerrorMessage("verbosity level parameter is of type int.\n");
2466  return SCIP_PARAMETERWRONGTYPE;
2467  }
2468 
2469  case SCIP_NLPPAR_FEASTOL:
2470  {
2471  *dval = problem->feastol;
2472  break;
2473  }
2474 
2475  case SCIP_NLPPAR_RELOBJTOL:
2476  {
2477  *dval = problem->opttol;
2478  break;
2479  }
2480 
2481  case SCIP_NLPPAR_LOBJLIM:
2482  {
2483  *dval = problem->fmin;
2484  break;
2485  }
2486 
2487  case SCIP_NLPPAR_INFINITY:
2488  {
2489  if( problem )
2490  {
2491  *dval = SCIPnlpiOracleGetInfinity(problem->oracle);
2492  }
2493  else
2494  {
2495  SCIP_NLPIDATA* data = SCIPnlpiGetData(nlpi);
2496  assert(data != NULL);
2497  *dval = data->infinity;
2498  }
2499  break;
2500  }
2501 
2502  case SCIP_NLPPAR_ITLIM:
2503  {
2504  SCIPerrorMessage("iteration limit parameter is of type int.\n");
2505  return SCIP_PARAMETERWRONGTYPE;
2506  }
2507 
2508  case SCIP_NLPPAR_TILIM:
2509  {
2510  *dval = problem->maxtime;
2511  break;
2512  }
2513 
2514  case SCIP_NLPPAR_OPTFILE:
2515  {
2516  SCIPerrorMessage("option file parameter is of type string.\n");
2517  return SCIP_PARAMETERWRONGTYPE;
2518  }
2519 
2520  case SCIP_NLPPAR_FASTFAIL:
2521  {
2522  SCIPerrorMessage("fastfail parameter is of type int.\n");
2523  return SCIP_PARAMETERWRONGTYPE;
2524  }
2525 
2526  default:
2527  {
2528  SCIPerrorMessage("Parameter %d not known to FilterSQP interface.\n", type);
2529  return SCIP_PARAMETERUNKNOWN;
2530  }
2531  }
2532 
2533  return SCIP_OKAY;
2534 } /*lint !e715*/
2535 
2536 /** sets floating point parameter of NLP
2537  *
2538  * input:
2539  * - nlpi NLP interface structure
2540  * - problem datastructure for problem instance, can be NULL only if type == SCIP_NLPPAR_INFINITY
2541  * - type parameter number
2542  * - dval parameter value
2543  */
2544 static
2545 SCIP_DECL_NLPISETREALPAR( nlpiSetRealParFilterSQP )
2546 {
2547  assert(nlpi != NULL);
2548  assert(type == SCIP_NLPPAR_INFINITY || problem != NULL);
2549 
2550  switch( type )
2551  {
2553  {
2554  SCIPerrorMessage("from scratch parameter is of type int.\n");
2555  return SCIP_PARAMETERWRONGTYPE;
2556  }
2557 
2558  case SCIP_NLPPAR_VERBLEVEL:
2559  {
2560  SCIPerrorMessage("verbosity level parameter is of type int.\n");
2561  return SCIP_PARAMETERWRONGTYPE;
2562  }
2563 
2564  case SCIP_NLPPAR_FEASTOL:
2565  {
2566  if( dval >= 0 )
2567  {
2568  problem->feastol = dval;
2569  }
2570  else
2571  {
2572  SCIPerrorMessage("Value %g for parameter feasibility tolerance is negative\n", dval);
2573  return SCIP_PARAMETERWRONGVAL;
2574  }
2575  break;
2576  }
2577 
2578  case SCIP_NLPPAR_RELOBJTOL:
2579  {
2580  if( dval >= 0 )
2581  {
2582  problem->opttol = dval;
2583  }
2584  else
2585  {
2586  SCIPerrorMessage("Value %g for parameter relative objective tolerance is negative\n", dval);
2587  return SCIP_PARAMETERWRONGVAL;
2588  }
2589  break;
2590  }
2591 
2592  case SCIP_NLPPAR_LOBJLIM:
2593  {
2594  problem->fmin = dval;
2595  break;
2596  }
2597 
2598  case SCIP_NLPPAR_INFINITY:
2599  {
2600  if( dval < 0.0 )
2601  return SCIP_PARAMETERWRONGVAL;
2602  if( problem )
2603  {
2604  SCIP_CALL( SCIPnlpiOracleSetInfinity(problem->oracle, dval) );
2605  }
2606  else
2607  {
2608  SCIP_NLPIDATA* data = SCIPnlpiGetData(nlpi);
2609  assert(data != NULL);
2610  data->infinity = dval;
2611  }
2612  break;
2613  }
2614 
2615  case SCIP_NLPPAR_ITLIM:
2616  {
2617  SCIPerrorMessage("iteration limit parameter is of type int.\n");
2618  return SCIP_PARAMETERWRONGTYPE;
2619  }
2620 
2621  case SCIP_NLPPAR_TILIM:
2622  {
2623  if( dval >= 0 )
2624  {
2625  problem->maxtime = dval;
2626  }
2627  else
2628  {
2629  SCIPerrorMessage("Value %g for parameter time limit is negative\n", dval);
2630  return SCIP_PARAMETERWRONGVAL;
2631  }
2632  break;
2633  }
2634 
2635  case SCIP_NLPPAR_OPTFILE:
2636  {
2637  SCIPerrorMessage("option file parameter is of type string.\n");
2638  return SCIP_PARAMETERWRONGTYPE;
2639  }
2640 
2641  case SCIP_NLPPAR_FASTFAIL:
2642  {
2643  SCIPerrorMessage("fastfail parameter is of type int.\n");
2644  return SCIP_PARAMETERWRONGTYPE;
2645  }
2646 
2647  default:
2648  {
2649  SCIPerrorMessage("Parameter %d not known to FilterSQP interface.\n", type);
2650  return SCIP_PARAMETERUNKNOWN;
2651  }
2652  }
2653 
2654  return SCIP_OKAY;
2655 } /*lint !e715*/
2656 
2657 /** gets string parameter of NLP
2658  *
2659  * input:
2660  * - nlpi NLP interface structure
2661  * - problem datastructure for problem instance
2662  * - type parameter number
2663  * - sval pointer to store the string value, the user must not modify the string
2664  *
2665  * output:
2666  * - sval parameter value
2667  */
2668 static
2669 SCIP_DECL_NLPIGETSTRINGPAR( nlpiGetStringParFilterSQP )
2670 {
2671  assert(nlpi != NULL);
2672  assert(problem != NULL);
2673 
2674  switch( type )
2675  {
2677  {
2678  SCIPerrorMessage("from scratch parameter is of type int.\n");
2679  return SCIP_PARAMETERWRONGTYPE;
2680  }
2681 
2682  case SCIP_NLPPAR_VERBLEVEL:
2683  {
2684  SCIPerrorMessage("verbosity level parameter is of type int.\n");
2685  return SCIP_PARAMETERWRONGTYPE;
2686  }
2687 
2688  case SCIP_NLPPAR_FEASTOL:
2689  {
2690  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
2691  return SCIP_PARAMETERWRONGTYPE;
2692  }
2693 
2694  case SCIP_NLPPAR_RELOBJTOL:
2695  {
2696  SCIPerrorMessage("objective tolerance parameter is of type real.\n");
2697  return SCIP_PARAMETERWRONGTYPE;
2698  }
2699 
2700  case SCIP_NLPPAR_LOBJLIM:
2701  {
2702  SCIPerrorMessage("objective limit parameter is of type real.\n");
2703  return SCIP_PARAMETERWRONGTYPE;
2704  }
2705 
2706  case SCIP_NLPPAR_INFINITY:
2707  {
2708  SCIPerrorMessage("infinity parameter is of type real.\n");
2709  return SCIP_PARAMETERWRONGTYPE;
2710  }
2711 
2712  case SCIP_NLPPAR_ITLIM:
2713  {
2714  SCIPerrorMessage("iteration limit parameter is of type int.\n");
2715  return SCIP_PARAMETERWRONGTYPE;
2716  }
2717 
2718  case SCIP_NLPPAR_TILIM:
2719  {
2720  SCIPerrorMessage("time limit parameter is of type real.\n");
2721  return SCIP_PARAMETERWRONGTYPE;
2722  }
2723 
2724  case SCIP_NLPPAR_OPTFILE:
2725  {
2726  *sval = NULL;
2727  return SCIP_OKAY;
2728  }
2729 
2730  case SCIP_NLPPAR_FASTFAIL:
2731  {
2732  SCIPerrorMessage("fastfail parameter is of type int.\n");
2733  return SCIP_PARAMETERWRONGTYPE;
2734  }
2735 
2736  default:
2737  {
2738  SCIPerrorMessage("Parameter %d not known to FilterSQP interface.\n", type);
2739  return SCIP_PARAMETERUNKNOWN;
2740  }
2741  }
2742 } /*lint !e715*/
2743 
2744 /** sets string parameter of NLP
2745  *
2746  * input:
2747  * - nlpi NLP interface structure
2748  * - problem datastructure for problem instance
2749  * - type parameter number
2750  * - sval parameter value
2751  */
2752 static
2753 SCIP_DECL_NLPISETSTRINGPAR( nlpiSetStringParFilterSQP )
2754 {
2755  assert(nlpi != NULL);
2756  assert(problem != NULL);
2757 
2758  switch( type )
2759  {
2761  {
2762  SCIPerrorMessage("from scratch parameter is of type int.\n");
2763  return SCIP_PARAMETERWRONGTYPE;
2764  }
2765 
2766  case SCIP_NLPPAR_VERBLEVEL:
2767  {
2768  SCIPerrorMessage("verbosity level parameter is of type int.\n");
2769  return SCIP_PARAMETERWRONGTYPE;
2770  }
2771 
2772  case SCIP_NLPPAR_FEASTOL:
2773  {
2774  SCIPerrorMessage("feasibility tolerance parameter is of type real.\n");
2775  return SCIP_PARAMETERWRONGTYPE;
2776  }
2777 
2778  case SCIP_NLPPAR_RELOBJTOL:
2779  {
2780  SCIPerrorMessage("objective tolerance parameter is of type real.\n");
2781  return SCIP_PARAMETERWRONGTYPE;
2782  }
2783 
2784  case SCIP_NLPPAR_LOBJLIM:
2785  {
2786  SCIPerrorMessage("objective limit parameter is of type real.\n");
2787  return SCIP_PARAMETERWRONGTYPE;
2788  }
2789 
2790  case SCIP_NLPPAR_INFINITY:
2791  {
2792  SCIPerrorMessage("infinity parameter is of type real.\n");
2793  return SCIP_PARAMETERWRONGTYPE;
2794  }
2795 
2796  case SCIP_NLPPAR_ITLIM:
2797  {
2798  SCIPerrorMessage("iteration limit parameter is of type int.\n");
2799  return SCIP_PARAMETERWRONGTYPE;
2800  }
2801 
2802  case SCIP_NLPPAR_TILIM:
2803  {
2804  SCIPerrorMessage("time limit parameter is of type real.\n");
2805  return SCIP_PARAMETERWRONGTYPE;
2806  }
2807 
2808  case SCIP_NLPPAR_OPTFILE:
2809  {
2810  SCIPmessagePrintWarning(SCIPnlpiGetData(nlpi)->messagehdlr, "Parameter optfile not supported by FilterSQP interface. Ignored.\n");
2811  return SCIP_OKAY;
2812  }
2813 
2814  case SCIP_NLPPAR_FASTFAIL:
2815  {
2816  SCIPerrorMessage("fastfail parameter is of type int.\n");
2817  return SCIP_PARAMETERWRONGTYPE;
2818  }
2819 
2820  default:
2821  {
2822  SCIPerrorMessage("Parameter %d not known to Ipopt interface.\n", type);
2823  return SCIP_PARAMETERUNKNOWN;
2824  }
2825  }
2826 } /*lint !e715*/
2827 
2828 /** sets message handler for message output
2829  *
2830  * input:
2831  * - nlpi NLP interface structure
2832  * - messagehdlr SCIP message handler, or NULL to suppress all output
2833  */
2834 static
2835 SCIP_DECL_NLPISETMESSAGEHDLR( nlpiSetMessageHdlrFilterSQP )
2836 {
2837  SCIP_NLPIDATA* nlpidata;
2838 
2839  assert(nlpi != NULL);
2840 
2841  nlpidata = SCIPnlpiGetData(nlpi);
2842  assert(nlpidata != NULL);
2843 
2844  nlpidata->messagehdlr = messagehdlr;
2845 
2846  return SCIP_OKAY; /*lint !e527*/
2847 } /*lint !e715*/
2848 
2849 /*
2850  * NLP solver interface specific interface methods
2851  */
2852 
2853 /** create solver interface for FilterSQP solver */
2855  BMS_BLKMEM* blkmem, /**< block memory data structure */
2856  SCIP_NLPI** nlpi /**< pointer to buffer for nlpi address */
2857  )
2858 {
2859  SCIP_NLPIDATA* nlpidata;
2860 
2861  assert(blkmem != NULL);
2862  assert(nlpi != NULL);
2863 
2864  /* create filterSQP solver interface data */
2865  SCIP_ALLOC( BMSallocBlockMemory(blkmem, &nlpidata) );
2866 
2867  nlpidata->blkmem = blkmem;
2868  nlpidata->messagehdlr = NULL;
2869  nlpidata->infinity = SCIP_DEFAULT_INFINITY;
2870  nlpidata->randnumgen = NULL;
2871 
2872  /* create solver interface */
2873  SCIP_CALL( SCIPnlpiCreate(nlpi,
2875  nlpiCopyFilterSQP, nlpiFreeFilterSQP, nlpiGetSolverPointerFilterSQP,
2876  nlpiCreateProblemFilterSQP, nlpiFreeProblemFilterSQP, nlpiGetProblemPointerFilterSQP,
2877  nlpiAddVarsFilterSQP, nlpiAddConstraintsFilterSQP, nlpiSetObjectiveFilterSQP,
2878  nlpiChgVarBoundsFilterSQP, nlpiChgConsSidesFilterSQP, nlpiDelVarSetFilterSQP, nlpiDelConstraintSetFilterSQP,
2879  nlpiChgLinearCoefsFilterSQP, nlpiChgQuadraticCoefsFilterSQP, nlpiChgExprtreeFilterSQP, nlpiChgNonlinCoefFilterSQP,
2880  nlpiChgObjConstantFilterSQP, nlpiSetInitialGuessFilterSQP, nlpiSolveFilterSQP, nlpiGetSolstatFilterSQP, nlpiGetTermstatFilterSQP,
2881  nlpiGetSolutionFilterSQP, nlpiGetStatisticsFilterSQP,
2882  nlpiGetWarmstartSizeFilterSQP, nlpiGetWarmstartMemoFilterSQP, nlpiSetWarmstartMemoFilterSQP,
2883  nlpiGetIntParFilterSQP, nlpiSetIntParFilterSQP, nlpiGetRealParFilterSQP, nlpiSetRealParFilterSQP, nlpiGetStringParFilterSQP, nlpiSetStringParFilterSQP,
2884  nlpiSetMessageHdlrFilterSQP,
2885  nlpidata) );
2886 
2887  return SCIP_OKAY;
2888 }
2889 
2890 /** gets string that identifies filterSQP (version number) */
2892  void
2893  )
2894 {
2895  return "filterSQP"; /* TODO version number? */
2896 }
2897 
2898 /** gets string that describes filterSQP */
2900  void
2901  )
2902 {
2903  return NLPI_DESC;
2904 }
2905 
2906 /** returns whether filterSQP is available, i.e., whether it has been linked in */
2908  void
2909  )
2910 {
2911  return TRUE;
2912 }
static SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemFilterSQP)
static SCIP_DECL_NLPIGETWARMSTARTSIZE(nlpiGetWarmstartSizeFilterSQP)
#define BMSfreeBlockMemoryArrayNull(mem, ptr, num)
Definition: memory.h:459
enum SCIP_NlpTermStat SCIP_NLPTERMSTAT
Definition: type_nlpi.h:84
static SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatFilterSQP)
const SCIP_Real * SCIPnlpiOracleGetVarUbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2228
const char * SCIPgetSolverNameFilterSQP(void)
SCIP_RETCODE SCIPnlpiOracleGetJacobianSparsity(SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2513
static SCIP_RETCODE setupGradients(BMS_BLKMEM *blkmem, SCIP_NLPIORACLE *oracle, fint **la, int *lasize, real **a)
SCIP_RETCODE SCIPcreateNlpSolverFilterSQP(BMS_BLKMEM *blkmem, SCIP_NLPI **nlpi)
SCIP_RETCODE SCIPnlpiOracleGetHessianLagSparsity(SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2793
int SCIPnlpiOracleGetNVars(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2198
static SCIP_DECL_NLPIGETWARMSTARTMEMO(nlpiGetWarmstartMemoFilterSQP)
#define infinity
Definition: gastrans.c:71
static SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsFilterSQP)
static SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsFilterSQP)
static SCIP_DECL_NLPISETWARMSTARTMEMO(nlpiSetWarmstartMemoFilterSQP)
#define NLPI_PRIORITY
fint n_bqpd_calls
internal methods for NLPI solver interfaces
SCIP_Bool SCIPisFilterSQPAvailableFilterSQP(void)
SCIP_RETCODE SCIPnlpiOracleSetObjective(SCIP_NLPIORACLE *oracle, const SCIP_Real constant, int nlin, const int *lininds, const SCIP_Real *linvals, int nquadelems, const SCIP_QUADELEM *quadelems, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree)
Definition: nlpioracle.c:1610
SCIP_Real * varubdualvalues
static SCIP_DECL_NLPISOLVE(nlpiSolveFilterSQP)
SCIP_NLPIDATA * SCIPnlpiGetData(SCIP_NLPI *nlpi)
Definition: nlpi.c:735
methods to store an NLP and request function, gradient, and hessian values
SCIP_NLPIORACLE * oracle
#define FALSE
Definition: def.h:73
SCIP_Real SCIPnlpiOracleGetConstraintLhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2279
SCIP_Real SCIPrandomGetReal(SCIP_RANDNUMGEN *randnumgen, SCIP_Real minrandval, SCIP_Real maxrandval)
Definition: misc.c:9981
static SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesFilterSQP)
#define TRUE
Definition: def.h:72
#define NLPI_NAME
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:54
void F77_FUNC(filtersqp, FILTERSQP)
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveValue(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *objval)
Definition: nlpioracle.c:2400
static SCIP_DECL_NLPIADDVARS(nlpiAddVarsFilterSQP)
SCIP_Real * primalvalues
#define WORKSPACEGROWTHFACTOR
static SCIP_TIME gettime(void)
#define BMSallocMemoryArray(ptr, num)
Definition: memory.h:115
static SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionFilterSQP)
#define SCIPdebugMessage
Definition: pub_message.h:87
SCIP_RETCODE SCIPnlpiSetRealPar(SCIP_NLPI *nlpi, SCIP_NLPIPROBLEM *problem, SCIP_NLPPARAM type, SCIP_Real dval)
Definition: nlpi.c:672
#define OPTTOLFACTOR
static SCIP_DECL_NLPIGETPROBLEMPOINTER(nlpiGetProblemPointerFilterSQP)
static pthread_mutex_t filtersqpmutex
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveGradient(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *objval, SCIP_Real *objgrad)
Definition: nlpioracle.c:2464
SCIP_RETCODE SCIPnlpiOracleDelConsSet(SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1827
SCIP_RETCODE SCIPnlpiOracleChgObjConstant(SCIP_NLPIORACLE *oracle, SCIP_Real objconstant)
Definition: nlpioracle.c:2182
SCIP_VAR ** x
Definition: circlepacking.c:54
filterSQP NLP interface
int SCIPnlpiOracleGetNConstraints(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2208
real eps
SCIP_Real * varlbdualvalues
static SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetFilterSQP)
SCIP_RETCODE SCIPnlpiSetMessageHdlr(SCIP_NLPI *nlpi, SCIP_MESSAGEHDLR *messagehdlr)
Definition: nlpi.c:722
static SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatFilterSQP)
long ftnlen
static SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessFilterSQP)
int SCIPnlpiOracleGetConstraintDegree(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2320
#define BMSfreeMemoryArray(ptr)
Definition: memory.h:139
SCIP_RETCODE SCIPnlpiOracleChgExprtree(SCIP_NLPIORACLE *oracle, int considx, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree)
Definition: nlpioracle.c:2100
static SCIP_DECL_NLPICHGQUADCOEFS(nlpiChgQuadraticCoefsFilterSQP)
SCIP_RETCODE SCIPnlpiOracleCreate(BMS_BLKMEM *blkmem, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1328
#define DEFAULT_MAXITER
#define SCIPerrorMessage
Definition: pub_message.h:55
static SCIP_DECL_NLPIGETSOLVERPOINTER(nlpiGetSolverPointerFilterSQP)
const SCIP_Real * SCIPnlpiOracleGetVarLbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2218
SCIP_RETCODE SCIPnlpiOracleChgVarBounds(SCIP_NLPIORACLE *oracle, int nvars, const int *indices, const SCIP_Real *lbs, const SCIP_Real *ubs)
Definition: nlpioracle.c:1647
struct SCIP_NlpiData SCIP_NLPIDATA
Definition: type_nlpi.h:38
enum SCIP_NlpSolStat SCIP_NLPSOLSTAT
Definition: type_nlpi.h:69
SCIP_RETCODE SCIPnlpiOracleAddVars(SCIP_NLPIORACLE *oracle, int nvars, const SCIP_Real *lbs, const SCIP_Real *ubs, const char **varnames)
Definition: nlpioracle.c:1447
SCIP_Bool warmstart
fint scale_mode
void SCIPmessagePrintWarning(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr,...)
Definition: message.c:418
void SCIPnlpStatisticsSetNIterations(SCIP_NLPSTATISTICS *statistics, int niterations)
Definition: nlpi.c:838
#define BMSallocClearBlockMemoryArray(mem, ptr, num)
Definition: memory.h:446
internal miscellaneous methods
#define NULL
Definition: lpi_spx1.cpp:155
int fint
void SCIPrandomFree(SCIP_RANDNUMGEN **randnumgen, BMS_BLKMEM *blkmem)
Definition: misc.c:9943
fint n_bqpd_prfint
static SCIP_DECL_NLPIGETINTPAR(nlpiGetIntParFilterSQP)
#define SCIP_CALL(x)
Definition: def.h:370
static SCIP_DECL_NLPISETINTPAR(nlpiSetIntParFilterSQP)
static SCIP_DECL_NLPICOPY(nlpiCopyFilterSQP)
void SCIPnlpStatisticsSetTotalTime(SCIP_NLPSTATISTICS *statistics, SCIP_Real totaltime)
Definition: nlpi.c:848
const char * SCIPgetSolverDescFilterSQP(void)
void SCIPmessagePrintInfo(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr,...)
Definition: message.c:585
#define MAXNRUNS
static SCIP_DECL_NLPIFREE(nlpiFreeFilterSQP)
static SCIP_DECL_NLPIGETREALPAR(nlpiGetRealParFilterSQP)
SCIP_Real SCIPnlpiOracleGetConstraintRhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2292
SCIP_Real SCIPnlpiOracleGetInfinity(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1400
SCIP_Real solvetime
SCIP_RETCODE SCIPnlpiOracleFree(SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1355
static SCIP_Real timeelapsed(SCIP_NLPIDATA *nlpidata)
real infty
#define BMSfreeBlockMemory(mem, ptr)
Definition: memory.h:456
fint phr
#define MINEPS
static SCIP_Bool timelimitreached(SCIP_NLPIDATA *nlpidata, SCIP_NLPIPROBLEM *nlpiproblem)
#define SCIP_Bool
Definition: def.h:70
#define BMSallocBlockMemoryArray(mem, ptr, num)
Definition: memory.h:445
SCIP_RETCODE SCIPnlpiOracleSetProblemName(SCIP_NLPIORACLE *oracle, const char *name)
Definition: nlpioracle.c:1412
fint phl
#define RANDSEED
SCIP_RETCODE SCIPnlpiOracleAddConstraints(SCIP_NLPIORACLE *oracle, int nconss, const SCIP_Real *lhss, const SCIP_Real *rhss, const int *nlininds, int *const *lininds, SCIP_Real *const *linvals, const int *nquadelems, SCIP_QUADELEM *const *quadelems, int *const *exprvaridxs, SCIP_EXPRTREE *const *exprtrees, const char **consnames)
Definition: nlpioracle.c:1532
real ubd
#define BMSfreeBlockMemoryArray(mem, ptr, num)
Definition: memory.h:458
#define MAX(x, y)
Definition: tclique_def.h:83
static int calcGrowSize(int num)
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValue(SCIP_NLPIORACLE *oracle, int considx, const SCIP_Real *x, SCIP_Real *conval)
Definition: nlpioracle.c:2419
static SCIP_DECL_NLPIGETSTRINGPAR(nlpiGetStringParFilterSQP)
SCIP_RETCODE SCIPnlpiOracleChgLinearCoefs(SCIP_NLPIORACLE *oracle, int considx, int nentries, const int *varidxs, const SCIP_Real *newcoefs)
Definition: nlpioracle.c:1905
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:126
SCIP_RETCODE SCIPnlpiOracleChgQuadCoefs(SCIP_NLPIORACLE *oracle, int considx, int nquadelems, const SCIP_QUADELEM *quadelems)
Definition: nlpioracle.c:2002
#define NLPI_DESC
static SCIP_DECL_NLPISETMESSAGEHDLR(nlpiSetMessageHdlrFilterSQP)
SCIP_Real * initguess
SCIP_RETCODE SCIPnlpiOracleChgConsSides(SCIP_NLPIORACLE *oracle, int nconss, const int *indices, const SCIP_Real *lhss, const SCIP_Real *rhss)
Definition: nlpioracle.c:1683
#define BMSclearMemory(ptr)
Definition: memory.h:121
#define DEFAULT_LOBJLIM
SCIP_Bool fromscratch
static SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsFilterSQP)
static SCIP_DECL_NLPICHGNONLINCOEF(nlpiChgNonlinCoefFilterSQP)
#define DEFAULT_FEASOPTTOL
SCIP_RETCODE SCIPnlpiCreate(SCIP_NLPI **nlpi, const char *name, const char *description, int priority, SCIP_DECL_NLPICOPY((*nlpicopy)), SCIP_DECL_NLPIFREE((*nlpifree)), SCIP_DECL_NLPIGETSOLVERPOINTER((*nlpigetsolverpointer)), SCIP_DECL_NLPICREATEPROBLEM((*nlpicreateproblem)), SCIP_DECL_NLPIFREEPROBLEM((*nlpifreeproblem)), SCIP_DECL_NLPIGETPROBLEMPOINTER((*nlpigetproblempointer)), SCIP_DECL_NLPIADDVARS((*nlpiaddvars)), SCIP_DECL_NLPIADDCONSTRAINTS((*nlpiaddconstraints)), SCIP_DECL_NLPISETOBJECTIVE((*nlpisetobjective)), SCIP_DECL_NLPICHGVARBOUNDS((*nlpichgvarbounds)), SCIP_DECL_NLPICHGCONSSIDES((*nlpichgconssides)), SCIP_DECL_NLPIDELVARSET((*nlpidelvarset)), SCIP_DECL_NLPIDELCONSSET((*nlpidelconsset)), SCIP_DECL_NLPICHGLINEARCOEFS((*nlpichglinearcoefs)), SCIP_DECL_NLPICHGQUADCOEFS((*nlpichgquadcoefs)), SCIP_DECL_NLPICHGEXPRTREE((*nlpichgexprtree)), SCIP_DECL_NLPICHGNONLINCOEF((*nlpichgnonlincoef)), SCIP_DECL_NLPICHGOBJCONSTANT((*nlpichgobjconstant)), SCIP_DECL_NLPISETINITIALGUESS((*nlpisetinitialguess)), SCIP_DECL_NLPISOLVE((*nlpisolve)), SCIP_DECL_NLPIGETSOLSTAT((*nlpigetsolstat)), SCIP_DECL_NLPIGETTERMSTAT((*nlpigettermstat)), SCIP_DECL_NLPIGETSOLUTION((*nlpigetsolution)), SCIP_DECL_NLPIGETSTATISTICS((*nlpigetstatistics)), SCIP_DECL_NLPIGETWARMSTARTSIZE((*nlpigetwarmstartsize)), SCIP_DECL_NLPIGETWARMSTARTMEMO((*nlpigetwarmstartmemo)), SCIP_DECL_NLPISETWARMSTARTMEMO((*nlpisetwarmstartmemo)), SCIP_DECL_NLPIGETINTPAR((*nlpigetintpar)), SCIP_DECL_NLPISETINTPAR((*nlpisetintpar)), SCIP_DECL_NLPIGETREALPAR((*nlpigetrealpar)), SCIP_DECL_NLPISETREALPAR((*nlpisetrealpar)), SCIP_DECL_NLPIGETSTRINGPAR((*nlpigetstringpar)), SCIP_DECL_NLPISETSTRINGPAR((*nlpisetstringpar)), SCIP_DECL_NLPISETMESSAGEHDLR((*nlpisetmessagehdlr)), SCIP_NLPIDATA *nlpidata)
Definition: nlpi.c:41
SCIP_RETCODE SCIPnlpiOracleDelVarSet(SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1717
#define SCIP_DEFAULT_INFINITY
Definition: def.h:168
static SCIP_RETCODE processSolveOutcome(SCIP_NLPIDATA *nlpidata, SCIP_NLPIPROBLEM *problem, fint ifail, real *x, real *lam)
static SCIP_DECL_NLPIDELCONSSET(nlpiDelConstraintSetFilterSQP)
time_t sec
SCIP_RETCODE SCIPnlpiOracleChgExprParam(SCIP_NLPIORACLE *oracle, int considx, int paramidx, SCIP_Real paramval)
Definition: nlpioracle.c:2157
static SCIP_DECL_NLPICHGOBJCONSTANT(nlpiChgObjConstantFilterSQP)
SCIP_RETCODE SCIPnlpiOracleEvalHessianLag(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real objfactor, const SCIP_Real *lambda, SCIP_Real *hessian)
Definition: nlpioracle.c:2892
#define MAXPERTURB
public methods for message output
real tt
SCIP_VAR * a
Definition: circlepacking.c:57
#define SCIP_Real
Definition: def.h:163
static SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemFilterSQP)
SCIP_RETCODE SCIPrandomCreate(SCIP_RANDNUMGEN **randnumgen, BMS_BLKMEM *blkmem, unsigned int initialseed)
Definition: misc.c:9927
static void invalidateSolution(SCIP_NLPIPROBLEM *problem)
static SCIP_RETCODE setupHessian(BMS_BLKMEM *blkmem, SCIP_NLPIORACLE *oracle, fint **la, int *lasize)
#define SCIP_INVALID
Definition: def.h:183
static SCIP_DECL_NLPICHGEXPRTREE(nlpiChgExprtreeFilterSQP)
static SCIP_RETCODE setupStart(SCIP_NLPIDATA *data, SCIP_NLPIPROBLEM *problem, real *x, SCIP_Bool *success)
SCIP_RETCODE SCIPnlpiOracleSetInfinity(SCIP_NLPIORACLE *oracle, SCIP_Real infinity)
Definition: nlpioracle.c:1384
#define SCIPisFinite(x)
Definition: pub_misc.h:1861
SCIP_Real * consdualvalues
SCIP_NLPTERMSTAT termstat
static SCIP_DECL_NLPISETREALPAR(nlpiSetRealParFilterSQP)
SCIP_RETCODE SCIPnlpiOracleEvalJacobian(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *convals, SCIP_Real *jacobi)
Definition: nlpioracle.c:2654
SCIP_NLPSOLSTAT solstat
#define BMSallocBlockMemory(mem, ptr)
Definition: memory.h:443
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:122
static SCIP_DECL_NLPISETSTRINGPAR(nlpiSetStringParFilterSQP)
struct BMS_BlkMem BMS_BLKMEM
Definition: memory.h:429
static SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveFilterSQP)
#define SCIP_ALLOC(x)
Definition: def.h:381
#define SCIPABORT()
Definition: def.h:342
fint phc
fint phe
double real
#define BMSreallocBlockMemoryArray(mem, ptr, oldnum, newnum)
Definition: memory.h:449
static SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsFilterSQP)