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