Scippy

SCIP

Solving Constraint Integer Programs

nlpioracle.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-2022 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not visit scipopt.org. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file nlpioracle.c
17  * @ingroup OTHER_CFILES
18  * @brief implementation of NLPI oracle
19  * @author Stefan Vigerske
20  *
21  * @todo jacobi evaluation should be sparse
22  */
23 
24 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
25 
26 #include "scip/scip.h"
27 #include "scip/nlpioracle.h"
28 #include "scip/exprinterpret.h"
29 #include "scip/expr_pow.h"
30 #include "scip/expr_varidx.h"
31 
32 #include <string.h> /* for strlen */
33 
34 /**@name NLPI Oracle data structures */
35 /**@{ */
36 
38 {
39  SCIP_Real lhs; /**< left hand side (for constraint) or constant (for objective) */
40  SCIP_Real rhs; /**< right hand side (for constraint) or constant (for objective) */
41 
42  int linsize; /**< length of linidxs and lincoefs arrays */
43  int nlinidxs; /**< number of linear variable indices and coefficients */
44  int* linidxs; /**< variable indices in linear part, or NULL if none */
45  SCIP_Real* lincoefs; /**< variable coefficients in linear part, of NULL if none */
46 
47  SCIP_EXPR* expr; /**< expression for nonlinear part, or NULL if none */
48  SCIP_EXPRINTDATA* exprintdata; /**< expression interpret data for expression, or NULL if no expr or not compiled yet */
49 
50  char* name; /**< name of constraint */
51 };
53 
55 {
56  char* name; /**< name of problem */
57 
58  int varssize; /**< length of variables related arrays */
59  int nvars; /**< number of variables */
60  SCIP_Real* varlbs; /**< array with variable lower bounds */
61  SCIP_Real* varubs; /**< array with variable upper bounds */
62  char** varnames; /**< array with variable names */
63  int* varlincount; /**< array with number of appearances of variable in linear part of objective or constraints */
64  int* varnlcount; /**< array with number of appearances of variable in nonlinear part of objective or constraints */
65 
66  int consssize; /**< length of constraints related arrays */
67  int nconss; /**< number of constraints */
68  SCIP_NLPIORACLECONS** conss; /**< constraints, or NULL if none */
69 
70  SCIP_NLPIORACLECONS* objective; /**< objective */
71 
72  int* jacoffsets; /**< rowwise jacobi sparsity pattern: constraint offsets in jaccols */
73  int* jaccols; /**< rowwise jacobi sparsity pattern: indices of variables appearing in constraints */
74 
75  int* heslagoffsets; /**< rowwise sparsity pattern of hessian matrix of Lagrangian: row offsets in heslagcol */
76  int* heslagcols; /**< rowwise sparsity pattern of hessian matrix of Lagrangian: column indices; sorted for each row */
77 
78  SCIP_EXPRINT* exprinterpreter; /**< interpreter for expressions: evaluation and derivatives */
79  SCIP_CLOCK* evalclock; /**< clock measuring evaluation time */
80 };
81 
82 /**@} */
83 
84 /*lint -e440*/
85 /*lint -e441*/
86 /*lint -e866*/
87 
88 /**@name Local functions */
89 /**@{ */
90 
91 /** ensures that those arrays in oracle that store information on variables have at least a given length */
92 static
94  SCIP* scip, /**< SCIP data structure */
95  SCIP_NLPIORACLE* oracle, /**< NLPIORACLE data structure */
96  int minsize /**< minimal required size */
97  )
98 {
99  assert(oracle != NULL);
100 
101  if( minsize > oracle->varssize )
102  {
103  int newsize;
104 
105  newsize = SCIPcalcMemGrowSize(scip, minsize);
106  assert(newsize >= minsize);
107 
108  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &oracle->varlbs, oracle->varssize, newsize) );
109  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &oracle->varubs, oracle->varssize, newsize) );
110  if( oracle->varnames != NULL )
111  {
112  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &oracle->varnames, oracle->varssize, newsize) );
113  }
114  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &oracle->varlincount, oracle->varssize, newsize) );
115  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &oracle->varnlcount, oracle->varssize, newsize) );
116 
117  oracle->varssize = newsize;
118  }
119  assert(oracle->varssize >= minsize);
120 
121  return SCIP_OKAY;
122 }
123 
124 /** ensures that constraints array in oracle has at least a given length */
125 static
127  SCIP* scip, /**< SCIP data structure */
128  SCIP_NLPIORACLE* oracle, /**< NLPIORACLE data structure */
129  int minsize /**< minimal required size */
130  )
131 {
132  assert(oracle != NULL);
133 
134  SCIP_CALL( SCIPensureBlockMemoryArray(scip, &oracle->conss, &oracle->consssize, minsize) );
135  assert(oracle->consssize >= minsize);
136 
137  return SCIP_OKAY;
138 }
139 
140 /** ensures that arrays for linear part in a oracle constraints have at least a given length */
141 static
143  SCIP* scip, /**< SCIP data structure */
144  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
145  int minsize /**< minimal required size */
146  )
147 {
148  assert(cons != NULL);
149 
150  if( minsize > cons->linsize )
151  {
152  int newsize;
153 
154  newsize = SCIPcalcMemGrowSize(scip, minsize);
155  assert(newsize >= minsize);
156 
157  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &cons->linidxs, cons->linsize, newsize) );
158  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &cons->lincoefs, cons->linsize, newsize) );
159  cons->linsize = newsize;
160  }
161  assert(cons->linsize >= minsize);
162 
163  return SCIP_OKAY;
164 }
165 
166 /** ensures that a given array of integers has at least a given length */
167 static
169  SCIP* scip, /**< SCIP data structure */
170  int** intarray, /**< array of integers */
171  int* len, /**< length of array (modified if reallocated) */
172  int minsize /**< minimal required array length */
173  )
174 {
175  assert(intarray != NULL);
176  assert(len != NULL);
177 
178  SCIP_CALL( SCIPensureBlockMemoryArray(scip, intarray, len, minsize) );
179  assert(*len >= minsize);
180 
181  return SCIP_OKAY;
182 }
183 
184 /** Invalidates the sparsity pattern of the Jacobian.
185  * Should be called when constraints are added or deleted.
186  */
187 static
189  SCIP* scip, /**< SCIP data structure */
190  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
191  )
192 {
193  assert(oracle != NULL);
194 
195  SCIPdebugMessage("%p invalidate jacobian sparsity\n", (void*)oracle);
196 
197  if( oracle->jacoffsets == NULL )
198  { /* nothing to do */
199  assert(oracle->jaccols == NULL);
200  return;
201  }
202 
203  assert(oracle->jaccols != NULL);
204  SCIPfreeBlockMemoryArray(scip, &oracle->jaccols, oracle->jacoffsets[oracle->nconss]);
205  SCIPfreeBlockMemoryArray(scip, &oracle->jacoffsets, oracle->nconss + 1);
206 }
207 
208 /** Invalidates the sparsity pattern of the Hessian of the Lagragian.
209  * Should be called when the objective is set or constraints are added or deleted.
210  */
211 static
213  SCIP* scip, /**< SCIP data structure */
214  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
215  )
216 {
217  assert(oracle != NULL);
218 
219  SCIPdebugMessage("%p invalidate hessian lag sparsity\n", (void*)oracle);
220 
221  if( oracle->heslagoffsets == NULL )
222  { /* nothing to do */
223  assert(oracle->heslagcols == NULL);
224  return;
225  }
226 
227  assert(oracle->heslagcols != NULL);
228  SCIPfreeBlockMemoryArray(scip, &oracle->heslagcols, oracle->heslagoffsets[oracle->nvars]);
229  SCIPfreeBlockMemoryArray(scip, &oracle->heslagoffsets, oracle->nvars + 1);
230 }
231 
232 /** increases or decreases variable counts in oracle w.r.t. linear and nonlinear appearance */
233 static
235  SCIP* scip, /**< SCIP data structure */
236  SCIP_NLPIORACLE* oracle, /**< oracle data structure */
237  int factor, /**< whether to add (factor=1) or remove (factor=-1) variable counts */
238  int nlinidxs, /**< number of linear indices */
239  const int* linidxs, /**< indices of variables in linear part */
240  SCIP_EXPR* expr /**< expression */
241  )
242 {
243  int j;
244 
245  assert(oracle != NULL);
246  assert(oracle->varlincount != NULL || (nlinidxs == 0 && expr == NULL));
247  assert(oracle->varnlcount != NULL || (nlinidxs == 0 && expr == NULL));
248  assert(factor == 1 || factor == -1);
249  assert(nlinidxs == 0 || linidxs != NULL);
250 
251  for( j = 0; j < nlinidxs; ++j )
252  {
253  oracle->varlincount[linidxs[j]] += factor;
254  assert(oracle->varlincount[linidxs[j]] >= 0);
255  }
256 
257  if( expr != NULL )
258  {
259  SCIP_EXPRITER* it;
260 
261  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
263 
264  for( ; !SCIPexpriterIsEnd(it); expr = SCIPexpriterGetNext(it) )
265  if( SCIPisExprVaridx(scip, expr) )
266  {
267  oracle->varnlcount[SCIPgetIndexExprVaridx(expr)] += factor;
268  assert(oracle->varnlcount[SCIPgetIndexExprVaridx(expr)] >= 0);
269  }
270 
271  SCIPfreeExpriter(&it);
272  }
273 
274  return SCIP_OKAY;
275 }
276 
277 /** sorts a linear term, merges duplicate entries and removes entries with coefficient 0.0 */
278 static
280  int* nidxs, /**< number of variables */
281  int* idxs, /**< indices of variables */
282  SCIP_Real* coefs /**< coefficients of variables */
283  )
284 {
285  int offset;
286  int j;
287 
288  assert(nidxs != NULL);
289  assert(idxs != NULL || *nidxs == 0);
290  assert(coefs != NULL || *nidxs == 0);
291 
292  if( *nidxs == 0 )
293  return;
294 
295  SCIPsortIntReal(idxs, coefs, *nidxs);
296 
297  offset = 0;
298  j = 0;
299  while( j+offset < *nidxs )
300  {
301  assert(idxs[j] >= 0); /*lint !e613*/
302 
303  /* move j+offset to j, if different */
304  if( offset > 0 )
305  {
306  idxs[j] = idxs[j+offset]; /*lint !e613*/
307  coefs[j] = coefs[j+offset]; /*lint !e613*/
308  }
309 
310  /* add up coefs for j+offset+1... as long as they have the same index */
311  while( j+offset+1 < *nidxs && idxs[j] == idxs[j+offset+1] ) /*lint !e613*/
312  {
313  coefs[j] += coefs[j+offset+1]; /*lint !e613*/
314  ++offset;
315  }
316 
317  /* if j'th element is 0, increase offset, otherwise increase j */
318  if( coefs[j] == 0.0 ) /*lint !e613*/
319  ++offset;
320  else
321  ++j;
322  }
323  *nidxs -= offset;
324 }
325 
326 /** creates a NLPI constraint from given constraint data */
327 static
329  SCIP* scip, /**< SCIP data structure */
330  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
331  SCIP_NLPIORACLECONS** cons, /**< buffer where to store pointer to constraint */
332  int nlinidxs, /**< length of linear part */
333  const int* linidxs, /**< indices of linear part, or NULL if nlinidxs == 0 */
334  const SCIP_Real* lincoefs, /**< coefficients of linear part, or NULL if nlinidxs == 0 */
335  SCIP_EXPR* expr, /**< expression, or NULL */
336  SCIP_Real lhs, /**< left-hand-side of constraint */
337  SCIP_Real rhs, /**< right-hand-side of constraint */
338  const char* name /**< name of constraint, or NULL */
339  )
340 {
341  assert(cons != NULL);
342  assert(nlinidxs >= 0);
343  assert(linidxs != NULL || nlinidxs == 0);
344  assert(lincoefs != NULL || nlinidxs == 0);
345  assert(EPSLE(lhs, rhs, SCIP_DEFAULT_EPSILON));
346 
347  SCIP_CALL( SCIPallocClearBlockMemory(scip, cons) );
348  assert(*cons != NULL);
349 
350  if( nlinidxs > 0 )
351  {
352  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(*cons)->linidxs, linidxs, nlinidxs) );
353  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(*cons)->lincoefs, lincoefs, nlinidxs) );
354  (*cons)->linsize = nlinidxs;
355  (*cons)->nlinidxs = nlinidxs;
356 
357  /* sort, merge duplicates, remove zero's */
358  sortLinearCoefficients(&(*cons)->nlinidxs, (*cons)->linidxs, (*cons)->lincoefs);
359  assert((*cons)->linidxs[0] >= 0);
360  }
361 
362  if( expr != NULL )
363  {
364  (*cons)->expr = expr;
365  SCIPcaptureExpr(expr);
366 
367  SCIP_CALL( SCIPexprintCompile(scip, oracle->exprinterpreter, (*cons)->expr, &(*cons)->exprintdata) );
368  }
369 
370  if( lhs > rhs )
371  {
372  assert(EPSEQ(lhs, rhs, SCIP_DEFAULT_EPSILON));
373  lhs = rhs;
374  }
375  (*cons)->lhs = lhs;
376  (*cons)->rhs = rhs;
377 
378  if( name != NULL )
379  {
380  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(*cons)->name, name, strlen(name)+1) );
381  }
382 
383  /* add variable counts */
384  SCIP_CALL( updateVariableCounts(scip, oracle, 1, (*cons)->nlinidxs, (*cons)->linidxs, (*cons)->expr) );
385 
386  return SCIP_OKAY;
387 }
388 
389 /** frees a constraint */
390 static
392  SCIP* scip, /**< SCIP data structure */
393  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
394  SCIP_NLPIORACLECONS** cons, /**< pointer to constraint that should be freed */
395  SCIP_Bool updatevarcount /**< whether the update variable counts (typically TRUE) */
396  )
397 {
398  assert(oracle != NULL);
399  assert(cons != NULL);
400  assert(*cons != NULL);
401 
402  SCIPdebugMessage("free constraint %p\n", (void*)*cons);
403 
404  /* remove variable counts */
405  if( updatevarcount )
406  {
407  SCIP_CALL( updateVariableCounts(scip, oracle, -1, (*cons)->nlinidxs, (*cons)->linidxs, (*cons)->expr) );
408  }
409 
410  SCIPfreeBlockMemoryArrayNull(scip, &(*cons)->linidxs, (*cons)->linsize);
411  SCIPfreeBlockMemoryArrayNull(scip, &(*cons)->lincoefs, (*cons)->linsize);
412 
413  if( (*cons)->expr != NULL )
414  {
415  SCIP_CALL( SCIPexprintFreeData(scip, oracle->exprinterpreter, (*cons)->expr, &(*cons)->exprintdata) );
416  SCIP_CALL( SCIPreleaseExpr(scip, &(*cons)->expr) );
417  }
418 
419  if( (*cons)->name != NULL )
420  {
421  SCIPfreeBlockMemoryArrayNull(scip, &(*cons)->name, strlen((*cons)->name)+1);
422  }
423 
424  SCIPfreeBlockMemory(scip, cons);
425  assert(*cons == NULL);
426 
427  return SCIP_OKAY;
428 }
429 
430 /** frees all constraints
431  *
432  * \attention This omits updating the variable counts in the oracle.
433  */
434 static
436  SCIP* scip, /**< SCIP data structure */
437  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
438  )
439 {
440  int i;
441 
442  assert(oracle != NULL);
443 
444  SCIPdebugMessage("%p free constraints\n", (void*)oracle);
445 
446  for( i = 0; i < oracle->nconss; ++i )
447  {
448  SCIP_CALL( freeConstraint(scip, oracle, &oracle->conss[i], FALSE) );
449  assert(oracle->conss[i] == NULL);
450  }
451  oracle->nconss = 0;
452 
453  SCIPfreeBlockMemoryArrayNull(scip, &oracle->conss, oracle->consssize);
454  oracle->consssize = 0;
455 
456  return SCIP_OKAY;
457 }
458 
459 /** moves one variable
460  * The place where it moves to need to be empty (all NULL) but allocated.
461  * Note that this function does not update the variable indices in the constraints!
462  */
463 static
465  SCIP* scip, /**< SCIP data structure */
466  SCIP_NLPIORACLE* oracle, /**< pointer to store NLPIORACLE data structure */
467  int fromidx, /**< index of variable to move */
468  int toidx /**< index of place where to move variable to */
469  )
470 {
471  assert(oracle != NULL);
472 
473  SCIPdebugMessage("%p move variable\n", (void*)oracle);
474 
475  assert(0 <= fromidx);
476  assert(0 <= toidx);
477  assert(fromidx < oracle->nvars);
478  assert(toidx < oracle->nvars);
479 
480  assert(oracle->varnames == NULL || oracle->varnames[toidx] == NULL);
481 
482  oracle->varlbs[toidx] = oracle->varlbs[fromidx];
483  oracle->varubs[toidx] = oracle->varubs[fromidx];
484  oracle->varlbs[fromidx] = -SCIPinfinity(scip);
485  oracle->varubs[fromidx] = SCIPinfinity(scip);
486 
487  oracle->varlincount[toidx] = oracle->varlincount[fromidx];
488  oracle->varnlcount[toidx] = oracle->varnlcount[fromidx];
489  oracle->varlincount[fromidx] = 0;
490  oracle->varnlcount[fromidx] = 0;
491 
492  if( oracle->varnames != NULL )
493  {
494  oracle->varnames[toidx] = oracle->varnames[fromidx];
495  oracle->varnames[fromidx] = NULL;
496  }
497 
498  return SCIP_OKAY;
499 }
500 
501 /** frees all variables */
502 static
504  SCIP* scip, /**< SCIP data structure */
505  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
506  )
507 {
508  int i;
509 
510  assert(oracle != NULL);
511 
512  SCIPdebugMessage("%p free variables\n", (void*)oracle);
513 
514  if( oracle->varnames != NULL )
515  {
516  for( i = 0; i < oracle->nvars; ++i )
517  {
518  if( oracle->varnames[i] != NULL )
519  {
520  SCIPfreeBlockMemoryArray(scip, &oracle->varnames[i], strlen(oracle->varnames[i])+1); /*lint !e866*/
521  }
522  }
523  SCIPfreeBlockMemoryArrayNull(scip, &oracle->varnames, oracle->varssize);
524  }
525  oracle->nvars = 0;
526 
527  SCIPfreeBlockMemoryArrayNull(scip, &oracle->varlbs, oracle->varssize);
528  SCIPfreeBlockMemoryArrayNull(scip, &oracle->varubs, oracle->varssize);
529  SCIPfreeBlockMemoryArrayNull(scip, &oracle->varlincount, oracle->varssize);
530  SCIPfreeBlockMemoryArrayNull(scip, &oracle->varnlcount, oracle->varssize);
531 
532  oracle->varssize = 0;
533 }
534 
535 /** applies a mapping of indices to one array of indices */
536 static
538  int* indexmap, /**< mapping from old variable indices to new indices */
539  int nindices, /**< number of indices in indices1 and indices2 */
540  int* indices /**< array of indices to adjust */
541  )
542 {
543  assert(indexmap != NULL);
544  assert(nindices == 0 || indices != NULL);
545 
546  for( ; nindices ; --nindices, ++indices )
547  {
548  assert(indexmap[*indices] >= 0);
549  *indices = indexmap[*indices];
550  }
551 }
552 
553 /** removes entries with index -1 (marked as deleted) from array of linear elements
554  * assumes that array is sorted by index, i.e., all -1 are at the beginning
555  */
556 static
558  int** linidxs, /**< variable indices */
559  SCIP_Real** coefs, /**< variable coefficients */
560  int* nidxs /**< number of indices */
561  )
562 {
563  int i;
564  int offset;
565 
566  SCIPdebugMessage("clear deleted linear elements\n");
567 
568  assert(linidxs != NULL);
569  assert(*linidxs != NULL);
570  assert(coefs != NULL);
571  assert(*coefs != NULL);
572  assert(nidxs != NULL);
573  assert(*nidxs > 0);
574 
575  /* search for beginning of non-delete entries @todo binary search? */
576  for( offset = 0; offset < *nidxs; ++offset )
577  if( (*linidxs)[offset] >= 0 )
578  break;
579 
580  /* nothing was deleted */
581  if( offset == 0 )
582  return;
583 
584  /* some or all elements were deleted -> move remaining ones front */
585  for( i = 0; i < *nidxs - offset; ++i )
586  {
587  (*linidxs)[i] = (*linidxs)[i+offset];
588  (*coefs)[i] = (*coefs) [i+offset];
589  }
590  *nidxs -= offset;
591 }
592 
593 /** computes the value of a function */
594 static
596  SCIP* scip, /**< SCIP data structure */
597  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
598  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
599  const SCIP_Real* x, /**< the point where to evaluate */
600  SCIP_Real* val /**< pointer to store function value */
601  )
602 { /*lint --e{715}*/
603  assert(oracle != NULL);
604  assert(cons != NULL);
605  assert(x != NULL || oracle->nvars == 0);
606  assert(val != NULL);
607 
608  SCIPdebugMessage("%p eval function value\n", (void*)oracle);
609 
610  *val = 0.0;
611 
612  if( cons->nlinidxs > 0 )
613  {
614  int* linidxs;
616  int nlin;
617 
618  nlin = cons->nlinidxs;
619  linidxs = cons->linidxs;
620  lincoefs = cons->lincoefs;
621  assert(linidxs != NULL);
622  assert(lincoefs != NULL);
623  assert(x != NULL);
624 
625  for( ; nlin > 0; --nlin, ++linidxs, ++lincoefs )
626  *val += *lincoefs * x[*linidxs];
627  }
628 
629  if( cons->expr != NULL )
630  {
631  SCIP_Real nlval;
632 
633  SCIP_CALL( SCIPexprintEval(scip, oracle->exprinterpreter, cons->expr, cons->exprintdata, (SCIP_Real*)x, &nlval) );
634  if( !SCIPisFinite(nlval) || SCIPisInfinity(scip, ABS(nlval)) )
635  *val = nlval;
636  else
637  *val += nlval;
638  }
639 
640  return SCIP_OKAY;
641 }
642 
643 /** computes the value and gradient of a function
644  *
645  * @return SCIP_INVALIDDATA, if the function or its gradient could not be evaluated (domain error, etc.)
646  */
647 static
649  SCIP* scip, /**< SCIP data structure */
650  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
651  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
652  const SCIP_Real* x, /**< the point where to evaluate */
653  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
654  SCIP_Real* RESTRICT val, /**< pointer to store function value */
655  SCIP_Real* RESTRICT grad /**< pointer to store function gradient */
656  )
657 { /*lint --e{715}*/
658  assert(oracle != NULL);
659  assert(x != NULL || oracle->nvars == 0);
660  assert(val != NULL);
661  assert(grad != NULL);
662 
663  SCIPdebugMessage("%p eval function gradient\n", (void*)oracle);
664 
665  *val = 0.0;
666  BMSclearMemoryArray(grad, oracle->nvars);
667 
668  if( cons->expr != NULL )
669  {
670  SCIP_Real nlval;
671  int i;
672 
673  SCIPdebugMsg(scip, "eval gradient of ");
674  SCIPdebug( if( isnewx ) {printf("\nx ="); for( i = 0; i < oracle->nvars; ++i) printf(" %g", x[i]); printf("\n");} )
675 
676  SCIP_CALL( SCIPexprintGrad(scip, oracle->exprinterpreter, cons->expr, cons->exprintdata, (SCIP_Real*)x, isnewx, &nlval, grad) );
677 
678  SCIPdebug( printf("g ="); for( i = 0; i < oracle->nvars; ++i) printf(" %g", grad[i]); printf("\n"); )
679 
680  /* check for eval error */
681  if( !SCIPisFinite(nlval) || SCIPisInfinity(scip, ABS(nlval)) )
682  {
683  SCIPdebugMessage("gradient evaluation yield invalid function value %g\n", nlval);
684  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
685  }
686  for( i = 0; i < oracle->nvars; ++i )
687  if( !SCIPisFinite(grad[i]) )
688  {
689  SCIPdebugMessage("gradient evaluation yield invalid gradient value %g\n", grad[i]);
690  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
691  }
692 
693  *val += nlval;
694  }
695 
696  if( cons->nlinidxs > 0 )
697  {
698  int* linidxs;
700  int nlin;
701 
702  nlin = cons->nlinidxs;
703  linidxs = cons->linidxs;
704  lincoefs = cons->lincoefs;
705  assert(linidxs != NULL);
706  assert(lincoefs != NULL);
707  assert(x != NULL);
708 
709  for( ; nlin > 0; --nlin, ++linidxs, ++lincoefs )
710  {
711  *val += *lincoefs * x[*linidxs];
712  grad[*linidxs] += *lincoefs;
713  }
714  }
715 
716  return SCIP_OKAY;
717 }
718 
719 /** collects indices of nonzero entries in the lower-left part of the hessian matrix of an expression
720  * adds the indices to a given set of indices, avoiding duplicates */
721 static
723  SCIP* scip, /**< SCIP data structure */
724  SCIP_NLPIORACLE* oracle, /**< NLPI oracle */
725  int** colnz, /**< indices of nonzero entries for each column */
726  int* collen, /**< space allocated to store indices of nonzeros for each column */
727  int* colnnz, /**< number of nonzero entries for each column */
728  int* nzcount, /**< counter for total number of nonzeros; should be increased when nzflag is set to 1 the first time */
729  SCIP_EXPR* expr, /**< expression */
730  SCIP_EXPRINTDATA* exprintdata, /**< expression interpreter data for expression */
731  int dim /**< dimension of matrix */
732  )
733 {
734  SCIP_Real* x;
735  int* rowidxs;
736  int* colidxs;
737  int nnz;
738  int row;
739  int col;
740  int pos;
741  int i;
742 
743  assert(oracle != NULL);
744  assert(colnz != NULL);
745  assert(collen != NULL);
746  assert(colnnz != NULL);
747  assert(nzcount != NULL);
748  assert(expr != NULL);
749  assert(dim >= 0);
750 
751  SCIPdebugMessage("%p hess lag sparsity set nzflag for expr\n", (void*)oracle);
752 
753  SCIP_CALL( SCIPallocBufferArray(scip, &x, oracle->nvars) );
754  for( i = 0; i < oracle->nvars; ++i )
755  x[i] = 2.0; /* hope that this value does not make much trouble for the evaluation routines */
756 
757  SCIP_CALL( SCIPexprintHessianSparsity(scip, oracle->exprinterpreter, expr, exprintdata, x, &rowidxs, &colidxs, &nnz) );
758 
759  for( i = 0; i < nnz; ++i )
760  {
761  row = rowidxs[i];
762  col = colidxs[i];
763 
764  assert(row < oracle->nvars);
765  assert(col <= row);
766 
767  if( colnz[row] == NULL || !SCIPsortedvecFindInt(colnz[row], col, colnnz[row], &pos) )
768  {
769  SCIP_CALL( ensureIntArraySize(scip, &colnz[row], &collen[row], colnnz[row]+1) );
770  SCIPsortedvecInsertInt(colnz[row], col, &colnnz[row], NULL);
771  ++*nzcount;
772  }
773  }
774 
775  SCIPfreeBufferArray(scip, &x);
776 
777  return SCIP_OKAY;
778 }
779 
780 /** adds hessian of an expression into hessian structure */
781 static
783  SCIP* scip, /**< SCIP data structure */
784  SCIP_NLPIORACLE* oracle, /**< oracle */
785  SCIP_Real weight, /**< weight of quadratic part */
786  const SCIP_Real* x, /**< point for which hessian should be returned */
787  SCIP_Bool new_x, /**< whether point has been evaluated before */
788  SCIP_EXPR* expr, /**< expression */
789  SCIP_EXPRINTDATA* exprintdata, /**< expression interpreter data for expression */
790  int* hesoffset, /**< row offsets in sparse matrix that is to be filled */
791  int* hescol, /**< column indices in sparse matrix that is to be filled */
792  SCIP_Real* values /**< buffer for values of sparse matrix that is to be filled */
793  )
794 {
795  SCIP_Real val;
796  SCIP_Real* h;
797  int* rowidxs;
798  int* colidxs;
799  int nnz;
800  int row;
801  int col;
802  int pos;
803  int i;
804 
805  SCIPdebugMessage("%p hess lag add expr\n", (void*)oracle);
806 
807  assert(oracle != NULL);
808  assert(x != NULL || new_x == FALSE);
809  assert(expr != NULL);
810  assert(hesoffset != NULL);
811  assert(hescol != NULL);
812  assert(values != NULL);
813 
814  SCIP_CALL( SCIPexprintHessian(scip, oracle->exprinterpreter, expr, exprintdata, (SCIP_Real*)x, new_x, &val, &rowidxs, &colidxs, &h, &nnz) );
815  if( !SCIPisFinite(val) )
816  {
817  SCIPdebugMessage("hessian evaluation yield invalid function value %g\n", val);
818  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
819  }
820 
821  for( i = 0; i < nnz; ++i )
822  {
823  if( !SCIPisFinite(h[i]) )
824  {
825  SCIPdebugMessage("hessian evaluation yield invalid hessian value %g\n", *h);
826  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
827  }
828 
829  if( h[i] == 0.0 )
830  continue;
831 
832  row = rowidxs[i];
833  col = colidxs[i];
834 
835  if( !SCIPsortedvecFindInt(&hescol[hesoffset[row]], col, hesoffset[row+1] - hesoffset[row], &pos) )
836  {
837  SCIPerrorMessage("Could not find entry (%d, %d) in hessian sparsity\n", row, col);
838  return SCIP_ERROR;
839  }
840 
841  values[hesoffset[row] + pos] += weight * h[i];
842  }
843 
844  return SCIP_OKAY;
845 }
846 
847 /** prints a name, if available, makes sure it has not more than 64 characters, and adds a unique prefix if the longnames flag is set */
848 static
850  char* buffer, /**< buffer to print to, has to be not NULL and should be at least 65 bytes */
851  char* name, /**< name, or NULL */
852  int idx, /**< index of var or cons which the name corresponds to */
853  char prefix, /**< a letter (typically 'x' or 'e') to distinguish variable and equation names, if names[idx] is not available */
854  const char* suffix, /**< a suffer to add to the name, or NULL */
855  SCIP_Bool longnames /**< whether prefixes for long names should be added */
856  )
857 {
858  assert(idx >= 0 && idx < 100000); /* to ensure that we do not exceed the size of the buffer */
859 
860  if( longnames )
861  {
862  if( name != NULL )
863  (void) SCIPsnprintf(buffer, 64, "%c%05d%.*s%s", prefix, idx, suffix != NULL ? (int)(57-strlen(suffix)) : 57, name, suffix ? suffix : "");
864  else
865  (void) SCIPsnprintf(buffer, 64, "%c%05d", prefix, idx);
866  }
867  else
868  {
869  if( name != NULL )
870  {
871  assert(strlen(name) + (suffix != NULL ? strlen(suffix) : 0) <= 64);
872  (void) SCIPsnprintf(buffer, 64, "%s%s", name, suffix != NULL ? suffix : "");
873  }
874  else
875  {
876  assert(1 + 5 + (suffix != NULL ? strlen(suffix) : 0) <= 64);
877  (void) SCIPsnprintf(buffer, 64, "%c%d%s", prefix, idx, suffix != NULL ? suffix : "");
878  }
879  }
880 }
881 
882 /** prints a function */
883 static
885  SCIP* scip, /**< SCIP data structure */
886  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
887  FILE* file, /**< file to print to, has to be not NULL */
888  SCIP_NLPIORACLECONS* cons, /**< constraint which function to print */
889  SCIP_Bool longvarnames /**< whether variable names need to be shorten to 64 characters */
890  )
891 { /*lint --e{715}*/
892  int i;
893  char namebuf[70];
894 
895  SCIPdebugMessage("%p print function\n", (void*)oracle);
896 
897  assert(oracle != NULL);
898  assert(file != NULL);
899  assert(cons != NULL);
900 
901  for( i = 0; i < cons->nlinidxs; ++i )
902  {
903  printName(namebuf, oracle->varnames != NULL ? oracle->varnames[cons->linidxs[i]] : NULL, cons->linidxs[i], 'x', NULL, longvarnames);
904  SCIPinfoMessage(scip, file, "%+.15g*%s", cons->lincoefs[i], namebuf);
905  if( i % 10 == 9 )
906  SCIPinfoMessage(scip, file, "\n");
907  }
908 
909  if( cons->expr != NULL )
910  {
911  /* TODO SCIPprintExpr does not use the variable names in oracle->varnames, probably that should be changed */
912  SCIPinfoMessage(scip, file, " +");
913  SCIP_CALL( SCIPprintExpr(scip, cons->expr, file) );
914  }
915 
916  return SCIP_OKAY;
917 }
918 
919 /** returns whether an expression contains nonsmooth operands (min, max, abs, ...) */
920 static
922  SCIP* scip, /**< SCIP data structure */
923  SCIP_EXPR* expr, /**< expression */
924  SCIP_Bool* nonsmooth /**< buffer to store whether expression seems nonsmooth */
925  )
926 {
927  SCIP_EXPRITER* it;
928 
929  assert(expr != NULL);
930  assert(nonsmooth != NULL);
931 
932  *nonsmooth = FALSE;
933 
934  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
936 
937  for( ; !SCIPexpriterIsEnd(it); expr = SCIPexpriterGetNext(it) )
938  {
939  const char* hdlrname;
940  if( SCIPisExprSignpower(scip, expr) )
941  {
942  *nonsmooth = TRUE;
943  break;
944  }
945  hdlrname = SCIPexprhdlrGetName(SCIPexprGetHdlr(expr));
946  if( strcmp(hdlrname, "abs") == 0 )
947  {
948  *nonsmooth = TRUE;
949  break;
950  }
951  if( strcmp(hdlrname, "min") == 0 )
952  {
953  *nonsmooth = TRUE;
954  break;
955  }
956  if( strcmp(hdlrname, "max") == 0 )
957  {
958  *nonsmooth = TRUE;
959  break;
960  }
961  }
962 
963  SCIPfreeExpriter(&it);
964 
965  return SCIP_OKAY;
966 }
967 
968 /**@} */
969 
970 /**@name public function */
971 /**@{ */
972 
973 /** creates an NLPIORACLE data structure */
975  SCIP* scip, /**< SCIP data structure */
976  SCIP_NLPIORACLE** oracle /**< pointer to store NLPIORACLE data structure */
977  )
978 {
979  SCIP_Bool nlpieval;
980 
981  assert(oracle != NULL);
982 
983  SCIPdebugMessage("%p oracle create\n", (void*)oracle);
984 
985  SCIP_CALL( SCIPallocMemory(scip, oracle) );
986  BMSclearMemory(*oracle);
987 
988  SCIPdebugMessage("Oracle initializes expression interpreter %s\n", SCIPexprintGetName());
989  SCIP_CALL( SCIPexprintCreate(scip, &(*oracle)->exprinterpreter) );
990 
991  SCIP_CALL( SCIPcreateClock(scip, &(*oracle)->evalclock) );
992 
993  SCIP_CALL( SCIPgetBoolParam(scip, "timing/nlpieval", &nlpieval) );
994  if( !nlpieval )
995  SCIPsetClockEnabled((*oracle)->evalclock, FALSE);
996 
997  /* create zero objective function */
998  SCIP_CALL( createConstraint(scip, *oracle, &(*oracle)->objective, 0, NULL, NULL, NULL, 0.0, 0.0, NULL) );
999 
1000  return SCIP_OKAY;
1001 }
1002 
1003 /** frees an NLPIORACLE data structure */
1005  SCIP* scip, /**< SCIP data structure */
1006  SCIP_NLPIORACLE** oracle /**< pointer to NLPIORACLE data structure */
1007  )
1008 {
1009  assert(oracle != NULL);
1010  assert(*oracle != NULL);
1011 
1012  SCIPdebugMessage("%p oracle free\n", (void*)oracle);
1013 
1014  invalidateJacobiSparsity(scip, *oracle);
1015  invalidateHessianLagSparsity(scip, *oracle);
1016 
1017  SCIP_CALL( freeConstraint(scip, *oracle, &(*oracle)->objective, FALSE) );
1018  SCIP_CALL( freeConstraints(scip, *oracle) );
1019  freeVariables(scip, *oracle);
1020 
1021  SCIP_CALL( SCIPfreeClock(scip, &(*oracle)->evalclock) );
1022 
1023  SCIP_CALL( SCIPexprintFree(scip, &(*oracle)->exprinterpreter) );
1024 
1025  if( (*oracle)->name != NULL )
1026  {
1027  SCIP_CALL( SCIPnlpiOracleSetProblemName(scip, *oracle, NULL) );
1028  }
1029 
1030  BMSfreeMemory(oracle);
1031 
1032  return SCIP_OKAY;
1033 }
1034 
1035 /** sets the problem name (used for printing) */
1037  SCIP* scip, /**< SCIP data structure */
1038  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1039  const char* name /**< name of problem */
1040  )
1041 {
1042  assert(oracle != NULL);
1043 
1044  SCIPdebugMessage("%p set problem name\n", (void*)oracle);
1045 
1046  if( oracle->name != NULL )
1047  {
1048  SCIPfreeBlockMemoryArray(scip, &oracle->name, strlen(oracle->name)+1);
1049  }
1050 
1051  if( name != NULL )
1052  {
1053  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &oracle->name, name, strlen(name)+1) );
1054  }
1055 
1056  return SCIP_OKAY;
1057 }
1058 
1059 /** gets the problem name, or NULL if none set */
1061  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1062  )
1063 {
1064  assert(oracle != NULL);
1065 
1066  SCIPdebugMessage("%p get problem name\n", (void*)oracle);
1067 
1068  return oracle->name;
1069 }
1070 
1071 /** adds variables */
1073  SCIP* scip, /**< SCIP data structure */
1074  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1075  int nvars, /**< number of variables to add */
1076  const SCIP_Real* lbs, /**< array with lower bounds of new variables, or NULL if all -infinity */
1077  const SCIP_Real* ubs, /**< array with upper bounds of new variables, or NULL if all +infinity */
1078  const char** varnames /**< array with names of new variables, or NULL if no names should be stored */
1079  )
1080 {
1081  int i;
1082 
1083  assert(oracle != NULL);
1084 
1085  SCIPdebugMessage("%p add vars\n", (void*)oracle);
1086 
1087  if( nvars == 0 )
1088  return SCIP_OKAY;
1089 
1090  assert(nvars > 0);
1091 
1092  SCIP_CALL( ensureVarsSize(scip, oracle, oracle->nvars + nvars) );
1093 
1094  if( lbs != NULL )
1095  {
1096  BMScopyMemoryArray(&oracle->varlbs[oracle->nvars], lbs, nvars);
1097  }
1098  else
1099  for( i = 0; i < nvars; ++i )
1100  oracle->varlbs[oracle->nvars+i] = -SCIPinfinity(scip);
1101 
1102  if( ubs != NULL )
1103  {
1104  BMScopyMemoryArray(&oracle->varubs[oracle->nvars], ubs, nvars);
1105 
1106  /* ensure variable bounds are consistent */
1107  for( i = oracle->nvars; i < oracle->nvars + nvars; ++i )
1108  {
1109  if( oracle->varlbs[i] > oracle->varubs[i] )
1110  {
1111  assert(EPSEQ(oracle->varlbs[i], oracle->varubs[i], SCIP_DEFAULT_EPSILON));
1112  oracle->varlbs[i] = oracle->varubs[i];
1113  }
1114  }
1115  }
1116  else
1117  for( i = 0; i < nvars; ++i )
1118  oracle->varubs[oracle->nvars+i] = SCIPinfinity(scip);
1119 
1120  if( varnames != NULL )
1121  {
1122  if( oracle->varnames == NULL )
1123  {
1124  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &oracle->varnames, oracle->varssize) );
1125  }
1126 
1127  for( i = 0; i < nvars; ++i )
1128  {
1129  if( varnames[i] != NULL )
1130  {
1131  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &oracle->varnames[oracle->nvars+i], varnames[i], strlen(varnames[i])+1) );
1132  }
1133  else
1134  oracle->varnames[oracle->nvars+i] = NULL;
1135  }
1136  }
1137  else if( oracle->varnames != NULL )
1138  {
1139  BMSclearMemoryArray(&oracle->varnames[oracle->nvars], nvars);
1140  }
1141 
1142  BMSclearMemoryArray(&oracle->varlincount[oracle->nvars], nvars);
1143  BMSclearMemoryArray(&oracle->varnlcount[oracle->nvars], nvars);
1144 
1145  /* @TODO update sparsity pattern by extending heslagoffsets */
1146  invalidateHessianLagSparsity(scip, oracle);
1147 
1148  oracle->nvars += nvars;
1149 
1150  return SCIP_OKAY;
1151 }
1152 
1153 /** adds constraints
1154  *
1155  * linear coefficients: row(=constraint) oriented matrix;
1156  * quadratic coefficients: row oriented matrix for each constraint
1157  */
1159  SCIP* scip, /**< SCIP data structure */
1160  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1161  int nconss, /**< number of constraints to add */
1162  const SCIP_Real* lhss, /**< array with left-hand sides of constraints, or NULL if all -infinity */
1163  const SCIP_Real* rhss, /**< array with right-hand sides of constraints, or NULL if all +infinity */
1164  const int* nlininds, /**< number of linear coefficients for each constraint, may be NULL in case of no linear part */
1165  int* const* lininds, /**< indices of variables for linear coefficients for each constraint, may be NULL in case of no linear part */
1166  SCIP_Real* const* linvals, /**< values of linear coefficient for each constraint, may be NULL in case of no linear part */
1167  SCIP_EXPR** exprs, /**< NULL if no nonlinear parts, otherwise exprs[.] gives nonlinear part,
1168  * or NULL if no nonlinear part in this constraint */
1169  const char** consnames /**< names of new constraints, or NULL if no names should be stored */
1170  )
1171 { /*lint --e{715}*/
1172  SCIP_NLPIORACLECONS* cons;
1173  SCIP_Bool addednlcon; /* whether a nonlinear constraint was added */
1174  int c;
1175 
1176  assert(oracle != NULL);
1177 
1178  SCIPdebugMessage("%p add constraints\n", (void*)oracle);
1179 
1180  if( nconss == 0 )
1181  return SCIP_OKAY;
1182 
1183  assert(nconss > 0);
1184 
1185  addednlcon = FALSE;
1186 
1187  invalidateJacobiSparsity(scip, oracle); /* @TODO we could also update (extend) the sparsity pattern */
1188 
1189  SCIP_CALL( ensureConssSize(scip, oracle, oracle->nconss + nconss) );
1190  for( c = 0; c < nconss; ++c )
1191  {
1192  SCIP_CALL( createConstraint(scip, oracle, &cons,
1193  nlininds != NULL ? nlininds[c] : 0,
1194  lininds != NULL ? lininds[c] : NULL,
1195  linvals != NULL ? linvals[c] : NULL,
1196  exprs != NULL ? exprs[c] : NULL,
1197  lhss != NULL ? lhss[c] : -SCIPinfinity(scip),
1198  rhss != NULL ? rhss[c] : SCIPinfinity(scip),
1199  consnames != NULL ? consnames[c] : NULL
1200  ) );
1201 
1202  if( cons->expr != NULL )
1203  addednlcon = TRUE;
1204 
1205  oracle->conss[oracle->nconss+c] = cons;
1206  }
1207  oracle->nconss += nconss;
1208 
1209  if( addednlcon == TRUE )
1210  invalidateHessianLagSparsity(scip, oracle);
1211 
1212  return SCIP_OKAY;
1213 }
1214 
1215 /** sets or overwrites objective, a minimization problem is expected
1216  *
1217  * May change sparsity pattern.
1218  */
1220  SCIP* scip, /**< SCIP data structure */
1221  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1222  const SCIP_Real constant, /**< constant part of objective */
1223  int nlin, /**< number of linear variable coefficients */
1224  const int* lininds, /**< indices of linear variables, or NULL if no linear part */
1225  const SCIP_Real* linvals, /**< coefficients of linear variables, or NULL if no linear part */
1226  SCIP_EXPR* expr /**< expression of nonlinear part, or NULL if no nonlinear part */
1227  )
1228 { /*lint --e{715}*/
1229  assert(oracle != NULL);
1230  assert(!SCIPisInfinity(scip, REALABS(constant)));
1231 
1232  SCIPdebugMessage("%p set objective\n", (void*)oracle);
1233 
1234  if( expr != NULL || oracle->objective->expr != NULL )
1235  invalidateHessianLagSparsity(scip, oracle);
1236 
1237  /* clear previous objective */
1238  SCIP_CALL( freeConstraint(scip, oracle, &oracle->objective, TRUE) );
1239 
1240  /* create new objective */
1241  SCIP_CALL( createConstraint(scip, oracle, &oracle->objective,
1242  nlin, lininds, linvals, expr, constant, constant, NULL) );
1243 
1244  return SCIP_OKAY;
1245 }
1246 
1247 /** change variable bounds */
1249  SCIP* scip, /**< SCIP data structure */
1250  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1251  int nvars, /**< number of variables to change bounds */
1252  const int* indices, /**< indices of variables to change bounds */
1253  const SCIP_Real* lbs, /**< new lower bounds, or NULL if all should be -infty */
1254  const SCIP_Real* ubs /**< new upper bounds, or NULL if all should be +infty */
1255  )
1256 {
1257  int i;
1258 
1259  assert(oracle != NULL);
1260  assert(indices != NULL || nvars == 0);
1261 
1262  SCIPdebugMessage("%p chg var bounds\n", (void*)oracle);
1263 
1264  for( i = 0; i < nvars; ++i )
1265  {
1266  assert(indices != NULL);
1267  assert(indices[i] >= 0);
1268  assert(indices[i] < oracle->nvars);
1269 
1270  oracle->varlbs[indices[i]] = (lbs != NULL ? lbs[i] : -SCIPinfinity(scip));
1271  oracle->varubs[indices[i]] = (ubs != NULL ? ubs[i] : SCIPinfinity(scip));
1272 
1273  if( oracle->varlbs[indices[i]] > oracle->varubs[indices[i]] )
1274  {
1275  /* inconsistent bounds; let's assume it's due to rounding and make them equal */
1276  assert(EPSEQ(oracle->varlbs[indices[i]], oracle->varubs[indices[i]], SCIP_DEFAULT_EPSILON));
1277  oracle->varlbs[indices[i]] = oracle->varubs[indices[i]];
1278  }
1279  }
1280 
1281  return SCIP_OKAY;
1282 }
1283 
1284 /** change constraint sides */
1286  SCIP* scip, /**< SCIP data structure */
1287  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1288  int nconss, /**< number of constraints to change bounds */
1289  const int* indices, /**< indices of constraints to change bounds */
1290  const SCIP_Real* lhss, /**< new left-hand sides, or NULL if all should be -infty */
1291  const SCIP_Real* rhss /**< new right-hand sides, or NULL if all should be +infty */
1292  )
1293 {
1294  int i;
1295 
1296  assert(oracle != NULL);
1297  assert(indices != NULL || nconss == 0);
1298 
1299  SCIPdebugMessage("%p chg cons sides\n", (void*)oracle);
1300 
1301  for( i = 0; i < nconss; ++i )
1302  {
1303  assert(indices != NULL);
1304  assert(indices[i] >= 0);
1305  assert(indices[i] < oracle->nconss);
1306 
1307  oracle->conss[indices[i]]->lhs = (lhss != NULL ? lhss[i] : -SCIPinfinity(scip));
1308  oracle->conss[indices[i]]->rhs = (rhss != NULL ? rhss[i] : SCIPinfinity(scip));
1309  if( oracle->conss[indices[i]]->lhs > oracle->conss[indices[i]]->rhs )
1310  {
1311  assert(EPSEQ(oracle->conss[indices[i]]->lhs, oracle->conss[indices[i]]->rhs, SCIP_DEFAULT_EPSILON));
1312  oracle->conss[indices[i]]->lhs = oracle->conss[indices[i]]->rhs;
1313  }
1314  }
1315 
1316  return SCIP_OKAY;
1317 }
1318 
1319 /** deletes a set of variables */
1321  SCIP* scip, /**< SCIP data structure */
1322  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1323  int* delstats /**< deletion status of vars in input (1 if var should be deleted, 0 if not);
1324  * new position of var in output (-1 if var was deleted) */
1325  )
1326 { /*lint --e{715}*/
1327  int c;
1328  int lastgood; /* index of the last variable that should be kept */
1329  SCIP_NLPIORACLECONS* cons;
1330  SCIP_EXPRITER* it;
1331 
1332  assert(oracle != NULL);
1333 
1334  SCIPdebugMessage("%p del var set\n", (void*)oracle);
1335 
1336  invalidateJacobiSparsity(scip, oracle);
1337  invalidateHessianLagSparsity(scip, oracle);
1338 
1339  lastgood = oracle->nvars - 1;
1340  while( lastgood >= 0 && delstats[lastgood] == 1 )
1341  --lastgood;
1342  if( lastgood < 0 )
1343  {
1344  /* all variables should be deleted */
1345  assert(oracle->nconss == 0); /* we could relax this by checking that all constraints are constant */
1346  oracle->objective->nlinidxs = 0;
1347  for( c = 0; c < oracle->nvars; ++c )
1348  delstats[c] = -1;
1349  freeVariables(scip, oracle);
1350  return SCIP_OKAY;
1351  }
1352 
1353  /* delete variables at the end */
1354  for( c = oracle->nvars - 1; c > lastgood; --c )
1355  {
1356  if( oracle->varnames && oracle->varnames[c] != NULL )
1357  {
1358  SCIPfreeBlockMemoryArray(scip, &oracle->varnames[c], strlen(oracle->varnames[c])+1);
1359  }
1360  delstats[c] = -1;
1361  }
1362 
1363  /* go through variables from the beginning on
1364  * if variable should be deleted, free it and move lastgood variable to this position
1365  * then update lastgood */
1366  for( c = 0; c <= lastgood; ++c )
1367  {
1368  if( delstats[c] == 0 )
1369  { /* variable should not be deleted and is kept on position c */
1370  delstats[c] = c;
1371  continue;
1372  }
1373  assert(delstats[c] == 1); /* variable should be deleted */
1374 
1375  if( oracle->varnames != NULL && oracle->varnames[c] != NULL )
1376  {
1377  SCIPfreeBlockMemoryArray(scip, &oracle->varnames[c], strlen(oracle->varnames[c])+1);
1378  }
1379  delstats[c] = -1;
1380 
1381  /* move variable at position lastgood to position c */
1382  SCIP_CALL( moveVariable(scip, oracle, lastgood, c) );
1383  delstats[lastgood] = c; /* mark that lastgood variable is now at position c */
1384 
1385  /* move lastgood forward, delete variables on the way */
1386  --lastgood;
1387  while( lastgood > c && delstats[lastgood] == 1)
1388  {
1389  if( oracle->varnames && oracle->varnames[lastgood] != NULL )
1390  {
1391  SCIPfreeBlockMemoryArray(scip, &oracle->varnames[lastgood], strlen(oracle->varnames[lastgood])+1);
1392  }
1393  delstats[lastgood] = -1;
1394  --lastgood;
1395  }
1396  }
1397  assert(c == lastgood);
1398 
1399  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
1400 
1401  for( c = -1; c < oracle->nconss; ++c )
1402  {
1403  cons = c < 0 ? oracle->objective : oracle->conss[c];
1404  assert(cons != NULL);
1405 
1406  /* update indices in linear part, sort indices, and then clear elements that are marked as deleted */
1407  mapIndices(delstats, cons->nlinidxs, cons->linidxs);
1408  SCIPsortIntReal(cons->linidxs, cons->lincoefs, cons->nlinidxs);
1409  clearDeletedLinearElements(&cons->linidxs, &cons->lincoefs, &cons->nlinidxs);
1410 
1411  if( cons->expr != NULL )
1412  {
1413  /* update variable indices in varidx expressions */
1414  SCIP_EXPR* expr;
1415  SCIP_Bool keptvar = FALSE; /* whether any of the variables in expr was not deleted */
1416 #ifndef NDEBUG
1417  SCIP_Bool delvar = FALSE; /* whether any of the variables in expr was deleted */
1418 #endif
1419 
1421  for( expr = cons->expr; !SCIPexpriterIsEnd(it); expr = SCIPexpriterGetNext(it) )
1422  {
1423  if( !SCIPisExprVaridx(scip, expr) )
1424  continue;
1425 
1426  if( delstats[SCIPgetIndexExprVaridx(expr)] >= 0 )
1427  {
1428  /* if variable is not deleted, then set its new index */
1429  keptvar = TRUE;
1430  SCIPsetIndexExprVaridx(expr, delstats[SCIPgetIndexExprVaridx(expr)]);
1431 
1432  /* if variable is kept, then there must not have been any variable that was deleted */
1433  assert(!delvar);
1434  }
1435  else
1436  {
1437 #ifndef NDEBUG
1438  delvar = TRUE;
1439 #endif
1440  /* if variable is deleted, then there must not have been any variable that was kept
1441  * (either all variables are deleted, which removes the expr, or none)
1442  */
1443  assert(!keptvar);
1444  }
1445  }
1446  if( !keptvar )
1447  {
1448  SCIP_CALL( SCIPexprintFreeData(scip, oracle->exprinterpreter, cons->expr, &cons->exprintdata) );
1449  SCIP_CALL( SCIPreleaseExpr(scip, &cons->expr) );
1450  }
1451  }
1452  }
1453 
1454  SCIPfreeExpriter(&it);
1455 
1456  oracle->nvars = lastgood+1;
1457 
1458  return SCIP_OKAY;
1459 }
1460 
1461 /** deletes a set of constraints */
1463  SCIP* scip, /**< SCIP data structure */
1464  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1465  int* delstats /**< array with deletion status of rows in input (1 if row should be deleted, 0 if not);
1466  * new position of row in output (-1 if row was deleted) */
1467  )
1468 { /*lint --e{715}*/
1469  int c;
1470  int lastgood; /* index of the last constraint that should be kept */
1471 
1472  assert(oracle != NULL);
1473 
1474  SCIPdebugMessage("%p del cons set\n", (void*)oracle);
1475 
1476  invalidateJacobiSparsity(scip, oracle);
1477  invalidateHessianLagSparsity(scip, oracle);
1478 
1479  lastgood = oracle->nconss - 1;
1480  while( lastgood >= 0 && delstats[lastgood] == 1)
1481  --lastgood;
1482  if( lastgood < 0 )
1483  {
1484  /* all constraints should be deleted */
1485  for( c = 0; c < oracle->nconss; ++c )
1486  delstats[c] = -1;
1487  SCIP_CALL( freeConstraints(scip, oracle) );
1488 
1489  /* the previous call did not keep variable counts uptodate
1490  * since we only have an objective function left, we reset the counts to the ones of the objective
1491  */
1492  BMSclearMemoryArray(oracle->varlincount, oracle->nvars);
1493  BMSclearMemoryArray(oracle->varnlcount, oracle->nvars);
1494  SCIP_CALL( updateVariableCounts(scip, oracle, 1, oracle->objective->nlinidxs, oracle->objective->linidxs, oracle->objective->expr) );
1495 
1496  return SCIP_OKAY;
1497  }
1498 
1499  /* delete constraints at the end */
1500  for( c = oracle->nconss - 1; c > lastgood; --c )
1501  {
1502  SCIP_CALL( freeConstraint(scip, oracle, &oracle->conss[c], TRUE) );
1503  assert(oracle->conss[c] == NULL);
1504  delstats[c] = -1;
1505  }
1506 
1507  /* go through constraint from the beginning on
1508  * if constraint should be deleted, free it and move lastgood constraint to this position
1509  * then update lastgood */
1510  for( c = 0; c <= lastgood; ++c )
1511  {
1512  if( delstats[c] == 0 )
1513  {
1514  /* constraint should not be deleted and is kept on position c */
1515  delstats[c] = c;
1516  continue;
1517  }
1518  assert(delstats[c] == 1); /* constraint should be deleted */
1519 
1520  SCIP_CALL( freeConstraint(scip, oracle, &oracle->conss[c], TRUE) );
1521  assert(oracle->conss[c] == NULL);
1522  delstats[c] = -1;
1523 
1524  /* move constraint at position lastgood to position c */
1525  oracle->conss[c] = oracle->conss[lastgood];
1526  assert(oracle->conss[c] != NULL);
1527  delstats[lastgood] = c; /* mark that lastgood constraint is now at position c */
1528  oracle->conss[lastgood] = NULL;
1529  --lastgood;
1530 
1531  /* move lastgood forward, delete constraints on the way */
1532  while( lastgood > c && delstats[lastgood] == 1)
1533  {
1534  SCIP_CALL( freeConstraint(scip, oracle, &oracle->conss[lastgood], TRUE) );
1535  assert(oracle->conss[lastgood] == NULL);
1536  delstats[lastgood] = -1;
1537  --lastgood;
1538  }
1539  }
1540  assert(c == lastgood+1);
1541 
1542  oracle->nconss = lastgood+1;
1543 
1544  return SCIP_OKAY;
1545 }
1546 
1547 /** changes (or adds) linear coefficients in one constraint or objective */
1549  SCIP* scip, /**< SCIP data structure */
1550  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1551  int considx, /**< index of constraint where linear coefficients should be changed, or -1 for objective */
1552  int nentries, /**< number of coefficients to change */
1553  const int* varidxs, /**< array with indices of variables which coefficients should be changed */
1554  const SCIP_Real* newcoefs /**< array with new coefficients of variables */
1555  )
1556 { /*lint --e{715}*/
1557  SCIP_NLPIORACLECONS* cons;
1558  SCIP_Bool needsort;
1559  int i;
1560 
1561  SCIPdebugMessage("%p chg linear coefs\n", (void*)oracle);
1562 
1563  assert(oracle != NULL);
1564  assert(varidxs != NULL || nentries == 0);
1565  assert(newcoefs != NULL || nentries == 0);
1566  assert(considx >= -1);
1567  assert(considx < oracle->nconss);
1568 
1569  if( nentries == 0 )
1570  return SCIP_OKAY;
1571 
1572  SCIPdebugMessage("change %d linear coefficients in cons %d\n", nentries, considx);
1573 
1574  needsort = FALSE;
1575 
1576  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
1577 
1578  if( cons->linsize == 0 )
1579  {
1580  /* first time we have linear coefficients in this constraint (or objective) */
1581  assert(cons->linidxs == NULL);
1582  assert(cons->lincoefs == NULL);
1583 
1584  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &cons->linidxs, varidxs, nentries) );
1585  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &cons->lincoefs, newcoefs, nentries) );
1586  cons->linsize = nentries;
1587  cons->nlinidxs = nentries;
1588 
1589  SCIP_CALL( updateVariableCounts(scip, oracle, 1, nentries, varidxs, NULL) );
1590 
1591  needsort = TRUE;
1592  }
1593  else
1594  {
1595  int pos;
1596 
1597  for( i = 0; i < nentries; ++i )
1598  {
1599  assert(varidxs[i] >= 0); /*lint !e613*/
1600  assert(varidxs[i] < oracle->nvars); /*lint !e613*/
1601 
1602  if( SCIPsortedvecFindInt(cons->linidxs, varidxs[i], cons->nlinidxs, &pos) ) /*lint !e613*/
1603  {
1604  SCIPdebugMessage("replace coefficient of var %d at pos %d by %g\n", varidxs[i], pos, newcoefs[i]); /*lint !e613*/
1605 
1606  cons->lincoefs[pos] = newcoefs[i]; /*lint !e613*/
1607 
1608  /* remember that we need to sort/merge/squeeze array if coefficient became zero here */
1609  needsort |= (newcoefs[i] == 0.0); /*lint !e613 !e514*/
1610 
1611  if( newcoefs[i] == 0.0 )
1612  {
1613  --oracle->varlincount[varidxs[i]];
1614  assert(oracle->varlincount[varidxs[i]] >= 0);
1615  }
1616  }
1617  else if( newcoefs[i] != 0.0 ) /*lint !e613*/
1618  {
1619  /* append new entry */
1620  SCIPdebugMessage("add coefficient of var %d at pos %d, value %g\n", varidxs[i], cons->nlinidxs, newcoefs[i]); /*lint !e613*/
1621 
1622  SCIP_CALL( ensureConsLinSize(scip, cons, cons->nlinidxs + (nentries-i)) );
1623  cons->linidxs[cons->nlinidxs] = varidxs[i]; /*lint !e613*/
1624  cons->lincoefs[cons->nlinidxs] = newcoefs[i]; /*lint !e613*/
1625  ++cons->nlinidxs;
1626 
1627  ++oracle->varlincount[varidxs[i]];
1628 
1629  needsort = TRUE;
1630  }
1631  }
1632  }
1633 
1634  if( needsort )
1635  {
1636  invalidateJacobiSparsity(scip, oracle);
1637  sortLinearCoefficients(&cons->nlinidxs, cons->linidxs, cons->lincoefs);
1638  }
1639 
1640  return SCIP_OKAY;
1641 }
1642 
1643 /** replaces expression of one constraint or objective */
1645  SCIP* scip, /**< SCIP data structure */
1646  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1647  int considx, /**< index of constraint where expression should be changed, or -1 for objective */
1648  SCIP_EXPR* expr /**< new expression, or NULL */
1649  )
1650 {
1651  SCIP_NLPIORACLECONS* cons;
1652 
1653  SCIPdebugMessage("%p chg expr\n", (void*)oracle);
1654 
1655  assert(oracle != NULL);
1656  assert(considx >= -1);
1657  assert(considx < oracle->nconss);
1658 
1659  invalidateHessianLagSparsity(scip, oracle);
1660  invalidateJacobiSparsity(scip, oracle);
1661 
1662  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
1663 
1664  /* free previous expression */
1665  if( cons->expr != NULL )
1666  {
1667  SCIP_CALL( updateVariableCounts(scip, oracle, -1, 0, NULL, cons->expr) );
1668  SCIP_CALL( SCIPexprintFreeData(scip, oracle->exprinterpreter, cons->expr, &cons->exprintdata) );
1669  SCIP_CALL( SCIPreleaseExpr(scip, &cons->expr) );
1670  }
1671 
1672  /* if user did not want to set new expr, then we are done */
1673  if( expr == NULL )
1674  return SCIP_OKAY;
1675 
1676  assert(oracle->exprinterpreter != NULL);
1677 
1678  /* install new expression */
1679  cons->expr = expr;
1680  SCIPcaptureExpr(cons->expr);
1681  SCIP_CALL( SCIPexprintCompile(scip, oracle->exprinterpreter, cons->expr, &cons->exprintdata) );
1682 
1683  /* keep variable counts up to date */
1684  SCIP_CALL( updateVariableCounts(scip, oracle, 1, 0, NULL, cons->expr) );
1685 
1686  return SCIP_OKAY;
1687 }
1688 
1689 /** changes the constant value in the objective function */ /*lint -e{715}*/
1691  SCIP* scip, /**< SCIP data structure */
1692  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1693  SCIP_Real objconstant /**< new value for objective constant */
1694  )
1695 { /*lint --e{715}*/
1696  assert(oracle != NULL);
1697 
1698  SCIPdebugMessage("%p chg obj constant\n", (void*)oracle);
1699 
1700  oracle->objective->lhs = objconstant;
1701  oracle->objective->rhs = objconstant;
1702 
1703  return SCIP_OKAY;
1704 }
1705 
1706 /** gives the current number of variables */
1708  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1709  )
1710 {
1711  assert(oracle != NULL);
1712 
1713  return oracle->nvars;
1714 }
1715 
1716 /** gives the current number of constraints */
1718  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1719  )
1720 {
1721  assert(oracle != NULL);
1722 
1723  return oracle->nconss;
1724 }
1725 
1726 /** gives the variables lower bounds */
1728  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1729  )
1730 {
1731  assert(oracle != NULL);
1732 
1733  return oracle->varlbs;
1734 }
1735 
1736 /** gives the variables upper bounds */
1738  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1739  )
1740 {
1741  assert(oracle != NULL);
1742 
1743  return oracle->varubs;
1744 }
1745 
1746 /** gives the variables names, or NULL if not set */
1748  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1749  )
1750 {
1751  assert(oracle != NULL);
1752 
1753  return oracle->varnames;
1754 }
1755 
1756 /** indicates whether variable appears nonlinear in any objective or constraint */ /*lint --e{715}*/
1758  SCIP* scip, /**< SCIP data structure */
1759  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1760  int varidx /**< the variable to check */
1761  )
1762 {
1763  assert(oracle != NULL);
1764  assert(varidx >= 0);
1765  assert(varidx < oracle->nvars);
1766  assert(oracle->varnlcount != NULL);
1767 
1768  return oracle->varnlcount[varidx] > 0;
1769 }
1770 
1771 /** returns number of linear and nonlinear appearances of variables in objective and constraints */ /*lint --e{715}*/
1773  SCIP* scip, /**< SCIP data structure */
1774  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1775  const int** lincounts, /**< buffer to return pointer to array of counts of linear appearances */
1776  const int** nlcounts /**< buffer to return pointer to array of counts of nonlinear appearances */
1777  )
1778 {
1779  assert(oracle != NULL);
1780  assert(lincounts != NULL);
1781  assert(nlcounts != NULL);
1782 
1783  *lincounts = oracle->varlincount;
1784  *nlcounts = oracle->varnlcount;
1785 }
1786 
1787 /** gives constant term of objective */
1789  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1790  )
1791 {
1792  assert(oracle != NULL);
1793  assert(oracle->objective->lhs == oracle->objective->rhs); /*lint !e777*/
1794 
1795  return oracle->objective->lhs;
1796 }
1797 
1798 /** gives left-hand side of a constraint */
1800  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1801  int considx /**< constraint index */
1802  )
1803 {
1804  assert(oracle != NULL);
1805  assert(considx >= 0);
1806  assert(considx < oracle->nconss);
1807 
1808  return oracle->conss[considx]->lhs;
1809 }
1810 
1811 /** gives right-hand side of a constraint */
1813  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1814  int considx /**< constraint index */
1815  )
1816 {
1817  assert(oracle != NULL);
1818  assert(considx >= 0);
1819  assert(considx < oracle->nconss);
1820 
1821  return oracle->conss[considx]->rhs;
1822 }
1823 
1824 /** gives name of a constraint, may be NULL */
1826  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1827  int considx /**< constraint index */
1828  )
1829 {
1830  assert(oracle != NULL);
1831  assert(considx >= 0);
1832  assert(considx < oracle->nconss);
1833 
1834  return oracle->conss[considx]->name;
1835 }
1836 
1837 /** indicates whether constraint is nonlinear */
1839  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1840  int considx /**< index of constraint for which nonlinearity status is returned, or -1 for objective */
1841  )
1842 {
1843  SCIP_NLPIORACLECONS* cons;
1844 
1845  assert(oracle != NULL);
1846  assert(considx >= -1);
1847  assert(considx < oracle->nconss);
1848 
1849  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
1850 
1851  return cons->expr != NULL;
1852 }
1853 
1854 /** gives the evaluation capabilities that are shared among all expressions in the problem */
1856  SCIP* scip, /**< SCIP data structure */
1857  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1858  )
1859 {
1860  int c;
1861  SCIP_EXPRINTCAPABILITY evalcapability;
1862 
1863  assert(oracle != NULL);
1864 
1865  if( oracle->objective->expr != NULL )
1866  evalcapability = SCIPexprintGetExprCapability(scip, oracle->exprinterpreter, oracle->objective->expr, oracle->objective->exprintdata);
1867  else
1868  evalcapability = SCIP_EXPRINTCAPABILITY_ALL;
1869 
1870  for( c = 0; c < oracle->nconss; ++c )
1871  if( oracle->conss[c]->expr != NULL )
1872  evalcapability &= SCIPexprintGetExprCapability(scip, oracle->exprinterpreter, oracle->conss[c]->expr, oracle->conss[c]->exprintdata);
1873 
1874  return evalcapability;
1875 }
1876 
1877 /** evaluates the objective function in a given point */
1879  SCIP* scip, /**< SCIP data structure */
1880  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1881  const SCIP_Real* x, /**< point where to evaluate */
1882  SCIP_Real* objval /**< pointer to store objective value */
1883  )
1884 {
1885  SCIP_RETCODE retcode;
1886 
1887  assert(oracle != NULL);
1888 
1889  SCIPdebugMessage("%p eval obj value\n", (void*)oracle);
1890 
1891  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
1892  retcode = evalFunctionValue(scip, oracle, oracle->objective, x, objval);
1893  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
1894 
1895  assert(oracle->objective->lhs == oracle->objective->rhs); /*lint !e777*/
1896  if( retcode == SCIP_OKAY )
1897  *objval += oracle->objective->lhs;
1898 
1899  return retcode;
1900 }
1901 
1902 /** evaluates one constraint function in a given point */
1904  SCIP* scip, /**< SCIP data structure */
1905  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1906  int considx, /**< index of constraint to evaluate */
1907  const SCIP_Real* x, /**< point where to evaluate */
1908  SCIP_Real* conval /**< pointer to store constraint value */
1909  )
1910 {
1911  SCIP_RETCODE retcode;
1912 
1913  assert(oracle != NULL);
1914  assert(x != NULL || oracle->nvars == 0);
1915  assert(conval != NULL);
1916 
1917  SCIPdebugMessage("%p eval cons value\n", (void*)oracle);
1918 
1919  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
1920  retcode = evalFunctionValue(scip, oracle, oracle->conss[considx], x, conval);
1921  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
1922 
1923  return retcode;
1924 }
1925 
1926 /** evaluates all constraint functions in a given point */
1928  SCIP* scip, /**< SCIP data structure */
1929  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1930  const SCIP_Real* x, /**< point where to evaluate */
1931  SCIP_Real* convals /**< buffer to store constraint values */
1932  )
1933 {
1934  SCIP_RETCODE retcode = SCIP_OKAY;
1935  int i;
1936 
1937  SCIPdebugMessage("%p eval cons values\n", (void*)oracle);
1938 
1939  assert(oracle != NULL);
1940  assert(x != NULL || oracle->nvars == 0);
1941  assert(convals != NULL);
1942 
1943  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
1944  for( i = 0; i < oracle->nconss; ++i )
1945  {
1946  retcode = evalFunctionValue(scip, oracle, oracle->conss[i], x, &convals[i]);
1947  if( retcode != SCIP_OKAY )
1948  break;
1949  }
1950  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
1951 
1952  return retcode;
1953 }
1954 
1955 /** computes the objective gradient in a given point
1956  *
1957  * @return SCIP_INVALIDDATA, if the function or its gradient could not be evaluated (domain error, etc.)
1958  */
1960  SCIP* scip, /**< SCIP data structure */
1961  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1962  const SCIP_Real* x, /**< point where to evaluate */
1963  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
1964  SCIP_Real* objval, /**< pointer to store objective value */
1965  SCIP_Real* objgrad /**< pointer to store (dense) objective gradient */
1966  )
1967 {
1968  SCIP_RETCODE retcode;
1969  assert(oracle != NULL);
1970 
1971  SCIPdebugMessage("%p eval obj grad\n", (void*)oracle);
1972 
1973  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
1974  retcode = evalFunctionGradient(scip, oracle, oracle->objective, x, isnewx, objval, objgrad);
1975  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
1976 
1977  assert(oracle->objective->lhs == oracle->objective->rhs); /*lint !e777*/
1978  if( retcode == SCIP_OKAY )
1979  *objval += oracle->objective->lhs;
1980 
1981  return retcode;
1982 }
1983 
1984 /** computes a constraints gradient in a given point
1985  *
1986  * @return SCIP_INVALIDDATA, if the function or its gradient could not be evaluated (domain error, etc.)
1987  */
1989  SCIP* scip, /**< SCIP data structure */
1990  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1991  const int considx, /**< index of constraint to compute gradient for */
1992  const SCIP_Real* x, /**< point where to evaluate */
1993  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
1994  SCIP_Real* conval, /**< pointer to store constraint value */
1995  SCIP_Real* congrad /**< pointer to store (dense) constraint gradient */
1996  )
1997 {
1998  SCIP_RETCODE retcode;
1999 
2000  assert(oracle != NULL);
2001  assert(x != NULL || oracle->nvars == 0);
2002  assert(conval != NULL);
2003 
2004  SCIPdebugMessage("%p eval cons grad\n", (void*)oracle);
2005 
2006  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
2007  retcode = evalFunctionGradient(scip, oracle, oracle->conss[considx], x, isnewx, conval, congrad);
2008  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
2009 
2010  return retcode;
2011 }
2012 
2013 /** gets sparsity pattern (rowwise) of Jacobian matrix
2014  *
2015  * Note that internal data is returned in *offset and *col, thus the user does not need to allocate memory there.
2016  * Adding or deleting constraints destroys the sparsity structure and make another call to this function necessary.
2017  */
2019  SCIP* scip, /**< SCIP data structure */
2020  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2021  const int** offset, /**< pointer to store pointer that stores the offsets to each rows sparsity pattern in col, can be NULL */
2022  const int** col /**< pointer to store pointer that stores the indices of variables that appear in each row, offset[nconss] gives length of col, can be NULL */
2023  )
2024 {
2025  SCIP_NLPIORACLECONS* cons;
2026  SCIP_EXPRITER* it;
2027  SCIP_EXPR* expr;
2028  SCIP_Bool* nzflag;
2029  int nnz;
2030  int maxnnz;
2031  int i;
2032  int j;
2033 
2034  assert(oracle != NULL);
2035 
2036  SCIPdebugMessage("%p get jacobian sparsity\n", (void*)oracle);
2037 
2038  if( oracle->jacoffsets != NULL )
2039  {
2040  assert(oracle->jaccols != NULL);
2041  if( offset != NULL )
2042  *offset = oracle->jacoffsets;
2043  if( col != NULL )
2044  *col = oracle->jaccols;
2045  return SCIP_OKAY;
2046  }
2047 
2048  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
2049 
2050  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &oracle->jacoffsets, oracle->nconss + 1) );
2051 
2052  maxnnz = MIN(oracle->nvars, 10) * oracle->nconss; /* initial guess */
2053  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &oracle->jaccols, maxnnz) );
2054 
2055  if( maxnnz == 0 )
2056  {
2057  /* no variables */
2058  BMSclearMemoryArray(oracle->jacoffsets, oracle->nconss + 1);
2059  if( offset != NULL )
2060  *offset = oracle->jacoffsets;
2061  if( col != NULL )
2062  *col = oracle->jaccols;
2063 
2064  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
2065 
2066  return SCIP_OKAY;
2067  }
2068  nnz = 0;
2069 
2070  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &nzflag, oracle->nvars) );
2071 
2072  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
2074 
2075  for( i = 0; i < oracle->nconss; ++i )
2076  {
2077  oracle->jacoffsets[i] = nnz;
2078 
2079  cons = oracle->conss[i];
2080  assert(cons != NULL);
2081 
2082  if( cons->expr == NULL )
2083  {
2084  /* for a linear constraint, we can just copy the linear indices from the constraint into the sparsity pattern */
2085  if( cons->nlinidxs > 0 )
2086  {
2087  SCIP_CALL( ensureIntArraySize(scip, &oracle->jaccols, &maxnnz, nnz + cons->nlinidxs) );
2088  BMScopyMemoryArray(&oracle->jaccols[nnz], cons->linidxs, cons->nlinidxs);
2089  nnz += cons->nlinidxs;
2090  }
2091  continue;
2092  }
2093 
2094  /* check which variables appear in constraint i
2095  * @todo this could be done faster for very sparse constraint by assembling all appearing variables, sorting, and removing duplicates
2096  */
2097  BMSclearMemoryArray(nzflag, oracle->nvars);
2098 
2099  for( j = 0; j < cons->nlinidxs; ++j )
2100  nzflag[cons->linidxs[j]] = TRUE;
2101 
2102  for( expr = SCIPexpriterRestartDFS(it, cons->expr); !SCIPexpriterIsEnd(it); expr = SCIPexpriterGetNext(it) )
2103  if( SCIPisExprVaridx(scip, expr) )
2104  {
2105  assert(SCIPgetIndexExprVaridx(expr) < oracle->nvars);
2106  nzflag[SCIPgetIndexExprVaridx(expr)] = TRUE;
2107  }
2108 
2109  /* store variables indices in jaccols */
2110  for( j = 0; j < oracle->nvars; ++j )
2111  {
2112  if( nzflag[j] == FALSE )
2113  continue;
2114 
2115  SCIP_CALL( ensureIntArraySize(scip, &oracle->jaccols, &maxnnz, nnz + 1) );
2116  oracle->jaccols[nnz] = j;
2117  ++nnz;
2118  }
2119  }
2120 
2121  SCIPfreeExpriter(&it);
2122 
2123  oracle->jacoffsets[oracle->nconss] = nnz;
2124 
2125  /* shrink jaccols array to nnz */
2126  if( nnz < maxnnz )
2127  {
2128  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &oracle->jaccols, maxnnz, nnz) );
2129  }
2130 
2131  SCIPfreeBlockMemoryArray(scip, &nzflag, oracle->nvars);
2132 
2133  if( offset != NULL )
2134  *offset = oracle->jacoffsets;
2135  if( col != NULL )
2136  *col = oracle->jaccols;
2137 
2138  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
2139 
2140  return SCIP_OKAY;
2141 }
2142 
2143 /** evaluates the Jacobian matrix in a given point
2144  *
2145  * The values in the Jacobian matrix are returned in the same order as specified by the offset and col arrays obtained by SCIPnlpiOracleGetJacobianSparsity().
2146  * The user need to call SCIPnlpiOracleGetJacobianSparsity() at least ones before using this function.
2147  *
2148  * @return SCIP_INVALIDDATA, if the Jacobian could not be evaluated (domain error, etc.)
2149  */
2151  SCIP* scip, /**< SCIP data structure */
2152  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2153  const SCIP_Real* x, /**< point where to evaluate */
2154  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
2155  SCIP_Real* convals, /**< pointer to store constraint values, can be NULL */
2156  SCIP_Real* jacobi /**< pointer to store sparse jacobian values */
2157  )
2158 {
2159  SCIP_NLPIORACLECONS* cons;
2160  SCIP_RETCODE retcode;
2161  SCIP_Real* grad;
2162  SCIP_Real nlval;
2163  int i;
2164  int j;
2165  int k;
2166  int l;
2167 
2168  SCIPdebugMessage("%p eval jacobian\n", (void*)oracle);
2169 
2170  assert(oracle != NULL);
2171  assert(jacobi != NULL);
2172 
2173  assert(oracle->jacoffsets != NULL);
2174  assert(oracle->jaccols != NULL);
2175 
2176  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
2177 
2178  SCIP_CALL( SCIPallocCleanBufferArray(scip, &grad, oracle->nvars) );
2179 
2180  retcode = SCIP_OKAY;
2181 
2182  j = oracle->jacoffsets[0]; /* TODO isn't oracle->jacoffsets[0] == 0 and thus always j == k ? */
2183  k = 0;
2184  for( i = 0; i < oracle->nconss; ++i )
2185  {
2186  cons = oracle->conss[i];
2187  assert(cons != NULL);
2188 
2189  if( cons->expr == NULL )
2190  {
2191  if( convals != NULL )
2192  convals[i] = 0.0;
2193 
2194  /* for a linear constraint, we can just copy the linear coefs from the constraint into the jacobian */
2195  if( cons->nlinidxs > 0 )
2196  {
2197  BMScopyMemoryArray(&jacobi[k], cons->lincoefs, cons->nlinidxs);
2198  j += cons->nlinidxs;
2199  k += cons->nlinidxs;
2200  if( convals != NULL )
2201  for( l = 0; l < cons->nlinidxs; ++l )
2202  convals[i] += cons->lincoefs[l] * x[cons->linidxs[l]];
2203  }
2204  assert(j == oracle->jacoffsets[i+1]);
2205  continue;
2206  }
2207 
2208  /* eval grad for nonlinear and add to jacobi */
2209  SCIPdebugMsg(scip, "eval gradient of ");
2210  SCIPdebug( if( isnewx ) {printf("\nx ="); for( l = 0; l < oracle->nvars; ++l) printf(" %g", x[l]); printf("\n");} )
2211 
2212  SCIP_CALL( SCIPexprintGrad(scip, oracle->exprinterpreter, cons->expr, cons->exprintdata, (SCIP_Real*)x, isnewx, &nlval, grad) );
2213 
2214  SCIPdebug( printf("g ="); for( l = oracle->jacoffsets[i]; l < oracle->jacoffsets[i+1]; ++l) printf(" %g", grad[oracle->jaccols[l]]); printf("\n"); )
2215 
2216  if( !SCIPisFinite(nlval) || SCIPisInfinity(scip, ABS(nlval)) )
2217  {
2218  SCIPdebugMessage("gradient evaluation yield invalid function value %g\n", nlval);
2219  retcode = SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
2220  goto TERMINATE;
2221  }
2222  if( convals != NULL )
2223  convals[i] = nlval;
2224 
2225  /* add linear part to grad */
2226  for( l = 0; l < cons->nlinidxs; ++l )
2227  {
2228  if( convals != NULL )
2229  convals[i] += cons->lincoefs[l] * x[cons->linidxs[l]];
2230  /* if grad[cons->linidxs[l]] is not finite, then adding a finite value doesn't change that, so don't check that here */
2231  grad[cons->linidxs[l]] += cons->lincoefs[l];
2232  }
2233 
2234  /* store complete gradient (linear + nonlinear) in jacobi
2235  * use the already evaluated sparsity pattern to pick only elements from grad that could have been set
2236  */
2237  assert(j == oracle->jacoffsets[i]);
2238  for( ; j < oracle->jacoffsets[i+1]; ++j )
2239  {
2240  if( !SCIPisFinite(grad[oracle->jaccols[j]]) )
2241  {
2242  SCIPdebugMessage("gradient evaluation yield invalid gradient value %g\n", grad[l]);
2243  retcode = SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
2244  goto TERMINATE;
2245  }
2246  jacobi[k++] = grad[oracle->jaccols[j]];
2247  /* reset to 0 for next constraint */
2248  grad[oracle->jaccols[j]] = 0.0;
2249  }
2250 
2251 #ifndef NDEBUG
2252  /* check that exprint really wrote only into expected elements of grad
2253  * TODO remove after some testing for better performance of debug runs */
2254  for( l = 0; l < oracle->nvars; ++l )
2255  assert(grad[l] == 0.0);
2256 #endif
2257  }
2258 
2259 TERMINATE:
2260  /* if there was an eval error, then we may have interrupted before cleaning up the grad buffer */
2261  if( retcode == SCIP_INVALIDDATA )
2262  BMSclearMemoryArray(grad, oracle->nvars);
2263 
2264  SCIPfreeCleanBufferArray(scip, &grad);
2265 
2266  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
2267 
2268  return retcode;
2269 }
2270 
2271 /** gets sparsity pattern of the Hessian matrix of the Lagrangian
2272  *
2273  * Note that internal data is returned in *offset and *col, thus the user must not to allocate memory there.
2274  * Adding or deleting variables, objective, or constraints may destroy the sparsity structure and make another call to this function necessary.
2275  * Only elements of the lower left triangle and the diagonal are counted.
2276  */
2278  SCIP* scip, /**< SCIP data structure */
2279  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2280  const int** offset, /**< pointer to store pointer that stores the offsets to each rows sparsity pattern in col, can be NULL */
2281  const int** col /**< pointer to store pointer that stores the indices of variables that appear in each row, offset[nconss] gives length of col, can be NULL */
2282  )
2283 {
2284  int** colnz; /* nonzeros in Hessian corresponding to one column */
2285  int* collen; /* collen[i] is length of array colnz[i] */
2286  int* colnnz; /* colnnz[i] is number of entries in colnz[i] (<= collen[i]) */
2287  int nnz;
2288  int i;
2289  int j;
2290  int cnt;
2291 
2292  assert(oracle != NULL);
2293 
2294  SCIPdebugMessage("%p get hessian lag sparsity\n", (void*)oracle);
2295 
2296  if( oracle->heslagoffsets != NULL )
2297  {
2298  assert(oracle->heslagcols != NULL);
2299  if( offset != NULL )
2300  *offset = oracle->heslagoffsets;
2301  if( col != NULL )
2302  *col = oracle->heslagcols;
2303  return SCIP_OKAY;
2304  }
2305 
2306  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
2307 
2308  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &oracle->heslagoffsets, oracle->nvars + 1) );
2309 
2310  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &colnz, oracle->nvars) );
2311  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &collen, oracle->nvars) );
2312  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &colnnz, oracle->nvars) );
2313  BMSclearMemoryArray(colnz, oracle->nvars);
2314  BMSclearMemoryArray(collen, oracle->nvars);
2315  BMSclearMemoryArray(colnnz, oracle->nvars);
2316  nnz = 0;
2317 
2318  if( oracle->objective->expr != NULL )
2319  {
2320  SCIP_CALL( hessLagSparsitySetNzFlagForExpr(scip, oracle, colnz, collen, colnnz, &nnz, oracle->objective->expr, oracle->objective->exprintdata, oracle->nvars) );
2321  }
2322 
2323  for( i = 0; i < oracle->nconss; ++i )
2324  {
2325  if( oracle->conss[i]->expr != NULL )
2326  {
2327  SCIP_CALL( hessLagSparsitySetNzFlagForExpr(scip, oracle, colnz, collen, colnnz, &nnz, oracle->conss[i]->expr, oracle->conss[i]->exprintdata, oracle->nvars) );
2328  }
2329  }
2330 
2331  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &oracle->heslagcols, nnz) );
2332 
2333  /* set hessian sparsity from colnz, colnnz */
2334  cnt = 0;
2335  for( i = 0; i < oracle->nvars; ++i )
2336  {
2337  oracle->heslagoffsets[i] = cnt;
2338  for( j = 0; j < colnnz[i]; ++j )
2339  {
2340  assert(cnt < nnz);
2341  oracle->heslagcols[cnt++] = colnz[i][j];
2342  }
2343  SCIPfreeBlockMemoryArrayNull(scip, &colnz[i], collen[i]);
2344  collen[i] = 0;
2345  }
2346  oracle->heslagoffsets[oracle->nvars] = cnt;
2347  assert(cnt == nnz);
2348 
2349  SCIPfreeBlockMemoryArray(scip, &colnz, oracle->nvars);
2350  SCIPfreeBlockMemoryArray(scip, &colnnz, oracle->nvars);
2351  SCIPfreeBlockMemoryArray(scip, &collen, oracle->nvars);
2352 
2353  if( offset != NULL )
2354  *offset = oracle->heslagoffsets;
2355  if( col != NULL )
2356  *col = oracle->heslagcols;
2357 
2358  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
2359 
2360  return SCIP_OKAY;
2361 }
2362 
2363 /** evaluates the Hessian matrix of the Lagrangian in a given point
2364  *
2365  * The values in the Hessian matrix are returned in the same order as specified by the offset and col arrays obtained by SCIPnlpiOracleGetHessianLagSparsity().
2366  * The user must call SCIPnlpiOracleGetHessianLagSparsity() at least ones before using this function.
2367  * Only elements of the lower left triangle and the diagonal are computed.
2368  *
2369  * @return SCIP_INVALIDDATA, if the Hessian could not be evaluated (domain error, etc.)
2370  */
2372  SCIP* scip, /**< SCIP data structure */
2373  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2374  const SCIP_Real* x, /**< point where to evaluate */
2375  SCIP_Bool isnewx_obj, /**< has the point x changed since the last call to an objective evaluation function? */
2376  SCIP_Bool isnewx_cons, /**< has the point x changed since the last call to the constraint evaluation function? */
2377  SCIP_Real objfactor, /**< weight for objective function */
2378  const SCIP_Real* lambda, /**< weights (Lagrangian multipliers) for the constraints */
2379  SCIP_Real* hessian /**< pointer to store sparse hessian values */
2380  )
2381 { /*lint --e{715}*/
2382  SCIP_RETCODE retcode = SCIP_OKAY;
2383  int i;
2384 
2385  assert(oracle != NULL);
2386  assert(x != NULL);
2387  assert(lambda != NULL || oracle->nconss == 0);
2388  assert(hessian != NULL);
2389 
2390  assert(oracle->heslagoffsets != NULL);
2391  assert(oracle->heslagcols != NULL);
2392 
2393  SCIPdebugMessage("%p eval hessian lag\n", (void*)oracle);
2394 
2395  SCIP_CALL( SCIPstartClock(scip, oracle->evalclock) );
2396 
2397  BMSclearMemoryArray(hessian, oracle->heslagoffsets[oracle->nvars]);
2398 
2399  if( objfactor != 0.0 && oracle->objective->expr != NULL )
2400  {
2401  retcode = hessLagAddExpr(scip, oracle, objfactor, x, isnewx_obj, oracle->objective->expr, oracle->objective->exprintdata, oracle->heslagoffsets, oracle->heslagcols, hessian);
2402  }
2403 
2404  for( i = 0; i < oracle->nconss && retcode == SCIP_OKAY; ++i )
2405  {
2406  assert( lambda != NULL ); /* for lint */
2407  if( lambda[i] == 0.0 || oracle->conss[i]->expr == NULL )
2408  continue;
2409  retcode = hessLagAddExpr(scip, oracle, lambda[i], x, isnewx_cons, oracle->conss[i]->expr, oracle->conss[i]->exprintdata, oracle->heslagoffsets, oracle->heslagcols, hessian);
2410  }
2411 
2412  SCIP_CALL( SCIPstopClock(scip, oracle->evalclock) );
2413 
2414  return retcode;
2415 }
2416 
2417 /** resets clock that measures evaluation time */
2419  SCIP* scip, /**< SCIP data structure */
2420  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2421  )
2422 {
2423  assert(oracle != NULL);
2424 
2425  SCIP_CALL( SCIPresetClock(scip, oracle->evalclock) );
2426 
2427  return SCIP_OKAY;
2428 }
2429 
2430 /** gives time spend in evaluation since last reset of clock
2431  *
2432  * Gives 0 if the eval clock is disabled.
2433  */
2435  SCIP* scip, /**< SCIP data structure */
2436  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2437  )
2438 {
2439  assert(oracle != NULL);
2440 
2441  return SCIPgetClockTime(scip, oracle->evalclock);
2442 }
2443 
2444 /** prints the problem to a file. */
2446  SCIP* scip, /**< SCIP data structure */
2447  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2448  FILE* file /**< file to print to, or NULL for standard output */
2449  )
2450 { /*lint --e{777} */
2451  int i;
2452  SCIP_Real lhs;
2453  SCIP_Real rhs;
2454 
2455  assert(oracle != NULL);
2456 
2457  SCIPdebugMessage("%p print problem\n", (void*)oracle);
2458 
2459  if( file == NULL )
2460  file = stdout;
2461 
2462  SCIPinfoMessage(scip, file, "NLPI Oracle %s: %d variables and %d constraints\n", oracle->name ? oracle->name : "", oracle->nvars, oracle->nconss);
2463  for( i = 0; i < oracle->nvars; ++i )
2464  {
2465  if( oracle->varnames != NULL && oracle->varnames[i] != NULL )
2466  SCIPinfoMessage(scip, file, "%10s (x%d)", oracle->varnames[i], i); /* give also name x%d as it will be by expression-print (printFunction) */
2467  else
2468  SCIPinfoMessage(scip, file, "x%09d", i);
2469  SCIPinfoMessage(scip, file, ": [%8g, %8g]", oracle->varlbs[i], oracle->varubs[i]);
2470  SCIPinfoMessage(scip, file, "\t #linear: %d #nonlinear: %d\n", oracle->varlincount[i], oracle->varnlcount[i]);
2471  }
2472 
2473  SCIPinfoMessage(scip, file, "objective: ");
2474  SCIP_CALL( printFunction(scip, oracle, file, oracle->objective, FALSE) );
2475  if( oracle->objective->lhs != 0.0 )
2476  SCIPinfoMessage(scip, file, "%+.15g", oracle->objective->lhs);
2477  SCIPinfoMessage(scip, file, "\n");
2478 
2479  for( i = 0; i < oracle->nconss; ++i )
2480  {
2481  if( oracle->conss[i]->name != NULL )
2482  SCIPinfoMessage(scip, file, "%10s", oracle->conss[i]->name);
2483  else
2484  SCIPinfoMessage(scip, file, "con%07d", i);
2485 
2486  lhs = oracle->conss[i]->lhs;
2487  rhs = oracle->conss[i]->rhs;
2488  SCIPinfoMessage(scip, file, ": ");
2489  if( !SCIPisInfinity(scip, -lhs) && !SCIPisInfinity(scip, rhs) && lhs != rhs )
2490  SCIPinfoMessage(scip, file, "%.15g <= ", lhs);
2491 
2492  SCIP_CALL( printFunction(scip, oracle, file, oracle->conss[i], FALSE) );
2493 
2494  if( lhs == rhs )
2495  SCIPinfoMessage(scip, file, " = %.15g", rhs);
2496  else if( !SCIPisInfinity(scip, rhs) )
2497  SCIPinfoMessage(scip, file, " <= %.15g", rhs);
2498  else if( !SCIPisInfinity(scip, -lhs) )
2499  SCIPinfoMessage(scip, file, " >= %.15g", lhs);
2500 
2501  SCIPinfoMessage(scip, file, "\n");
2502  }
2503 
2504  return SCIP_OKAY;
2505 }
2506 
2507 /** prints the problem to a file in GAMS format
2508  *
2509  * If there are variable (equation, resp.) names with more than 9 characters, then variable (equation, resp.) names are prefixed with an unique identifier.
2510  * This is to make it easier to identify variables solution output in the listing file.
2511  * Names with more than 64 characters are shorten to 64 letters due to GAMS limits.
2512  */
2514  SCIP* scip, /**< SCIP data structure */
2515  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2516  SCIP_Real* initval, /**< starting point values for variables or NULL */
2517  FILE* file /**< file to print to, or NULL for standard output */
2518  )
2519 { /*lint --e{777} */
2520  int i;
2521  int nllevel; /* level of nonlinearity of problem: linear = 0, quadratic, smooth nonlinear, nonsmooth */
2522  static const char* nllevelname[4] = { "LP", "QCP", "NLP", "DNLP" };
2523  char problemname[SCIP_MAXSTRLEN];
2524  char namebuf[70];
2525  SCIP_Bool havelongvarnames;
2526  SCIP_Bool havelongequnames;
2527 
2528  SCIPdebugMessage("%p print problem gams\n", (void*)oracle);
2529 
2530  assert(oracle != NULL);
2531 
2532  if( file == NULL )
2533  file = stdout;
2534 
2535  nllevel = 0;
2536 
2537  havelongvarnames = FALSE;
2538  for( i = 0; i < oracle->nvars; ++i )
2539  if( oracle->varnames != NULL && oracle->varnames[i] != NULL && strlen(oracle->varnames[i]) > 9 )
2540  {
2541  havelongvarnames = TRUE;
2542  break;
2543  }
2544 
2545  havelongequnames = FALSE;
2546  for( i = 0; i < oracle->nconss; ++i )
2547  if( oracle->conss[i]->name && strlen(oracle->conss[i]->name) > 9 )
2548  {
2549  havelongequnames = TRUE;
2550  break;
2551  }
2552 
2553  SCIPinfoMessage(scip, file, "$offlisting\n");
2554  SCIPinfoMessage(scip, file, "$offdigit\n");
2555  SCIPinfoMessage(scip, file, "* NLPI Oracle Problem %s\n", oracle->name ? oracle->name : "");
2556  SCIPinfoMessage(scip, file, "Variables ");
2557  for( i = 0; i < oracle->nvars; ++i )
2558  {
2559  printName(namebuf, oracle->varnames != NULL ? oracle->varnames[i] : NULL, i, 'x', NULL, havelongvarnames);
2560  SCIPinfoMessage(scip, file, "%s, ", namebuf);
2561  if( i % 10 == 9 )
2562  SCIPinfoMessage(scip, file, "\n");
2563  }
2564  SCIPinfoMessage(scip, file, "NLPIORACLEOBJVAR;\n\n");
2565  for( i = 0; i < oracle->nvars; ++i )
2566  {
2567  char* name;
2568  name = oracle->varnames != NULL ? oracle->varnames[i] : NULL;
2569  if( oracle->varlbs[i] == oracle->varubs[i] )
2570  {
2571  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
2572  SCIPinfoMessage(scip, file, "%s.fx = %.15g;\t", namebuf, oracle->varlbs[i]);
2573  }
2574  else
2575  {
2576  if( !SCIPisInfinity(scip, -oracle->varlbs[i]) )
2577  {
2578  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
2579  SCIPinfoMessage(scip, file, "%s.lo = %.15g;\t", namebuf, oracle->varlbs[i]);
2580  }
2581  if( !SCIPisInfinity(scip, oracle->varubs[i]) )
2582  {
2583  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
2584  SCIPinfoMessage(scip, file, "%s.up = %.15g;\t", namebuf, oracle->varubs[i]);
2585  }
2586  }
2587  if( initval != NULL )
2588  {
2589  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
2590  SCIPinfoMessage(scip, file, "%s.l = %.15g;\t", namebuf, initval[i]);
2591  }
2592  SCIPinfoMessage(scip, file, "\n");
2593  }
2594  SCIPinfoMessage(scip, file, "\n");
2595 
2596  SCIPinfoMessage(scip, file, "Equations ");
2597  for( i = 0; i < oracle->nconss; ++i )
2598  {
2599  printName(namebuf, oracle->conss[i]->name, i, 'e', NULL, havelongequnames);
2600  SCIPinfoMessage(scip, file, "%s, ", namebuf);
2601 
2602  if( !SCIPisInfinity(scip, -oracle->conss[i]->lhs) && !SCIPisInfinity(scip, oracle->conss[i]->rhs) && oracle->conss[i]->lhs != oracle->conss[i]->rhs )
2603  {
2604  /* ranged row: add second constraint */
2605  printName(namebuf, oracle->conss[i]->name, i, 'e', "_RNG", havelongequnames);
2606  SCIPinfoMessage(scip, file, "%s, ", namebuf);
2607  }
2608  if( i % 10 == 9 )
2609  SCIPinfoMessage(scip, file, "\n");
2610  }
2611  SCIPinfoMessage(scip, file, "NLPIORACLEOBJ;\n\n");
2612 
2613  SCIPinfoMessage(scip, file, "NLPIORACLEOBJ.. NLPIORACLEOBJVAR =E= ");
2614  SCIP_CALL( printFunction(scip, oracle, file, oracle->objective, havelongvarnames) );
2615  if( oracle->objective->lhs != 0.0 )
2616  SCIPinfoMessage(scip, file, "%+.15g", oracle->objective->lhs);
2617  SCIPinfoMessage(scip, file, ";\n");
2618 
2619  for( i = 0; i < oracle->nconss; ++i )
2620  {
2621  SCIP_Real lhs;
2622  SCIP_Real rhs;
2623 
2624  printName(namebuf, oracle->conss[i]->name, i, 'e', NULL, havelongequnames);
2625  SCIPinfoMessage(scip, file, "%s.. ", namebuf);
2626 
2627  SCIP_CALL( printFunction(scip, oracle, file, oracle->conss[i], havelongvarnames) );
2628 
2629  lhs = oracle->conss[i]->lhs;
2630  rhs = oracle->conss[i]->rhs;
2631 
2632  if( lhs == rhs )
2633  SCIPinfoMessage(scip, file, " =E= %.15g", rhs);
2634  else if( !SCIPisInfinity(scip, rhs) )
2635  SCIPinfoMessage(scip, file, " =L= %.15g", rhs);
2636  else if( !SCIPisInfinity(scip, -lhs) )
2637  SCIPinfoMessage(scip, file, " =G= %.15g", lhs);
2638  else
2639  SCIPinfoMessage(scip, file, " =N= 0");
2640  SCIPinfoMessage(scip, file, ";\n");
2641 
2642  if( !SCIPisInfinity(scip, lhs) && !SCIPisInfinity(scip, rhs) && lhs != rhs )
2643  {
2644  printName(namebuf, oracle->conss[i]->name, i, 'e', "_RNG", havelongequnames);
2645  SCIPinfoMessage(scip, file, "%s.. ", namebuf);
2646 
2647  SCIP_CALL( printFunction(scip, oracle, file, oracle->conss[i], havelongvarnames) );
2648 
2649  SCIPinfoMessage(scip, file, " =G= %.15g;\n", lhs);
2650  }
2651 
2652  if( nllevel <= 1 && oracle->conss[i]->expr != NULL )
2653  nllevel = 2;
2654  if( nllevel <= 2 && oracle->conss[i]->expr != NULL )
2655  {
2656  SCIP_Bool nonsmooth;
2657  SCIP_CALL( exprIsNonSmooth(scip, oracle->conss[i]->expr, &nonsmooth) );
2658  if( nonsmooth )
2659  nllevel = 3;
2660  }
2661  }
2662 
2663  (void) SCIPsnprintf(problemname, SCIP_MAXSTRLEN, "%s", oracle->name ? oracle->name : "m");
2664 
2665  SCIPinfoMessage(scip, file, "Model %s / all /;\n", problemname);
2666  SCIPinfoMessage(scip, file, "option limrow = 0;\n");
2667  SCIPinfoMessage(scip, file, "option limcol = 0;\n");
2668  SCIPinfoMessage(scip, file, "Solve %s minimizing NLPIORACLEOBJVAR using %s;\n", problemname, nllevelname[nllevel]);
2669 
2670  return SCIP_OKAY;
2671 }
2672 
2673 /**@} */
#define SCIPfreeBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:101
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveGradient(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *objval, SCIP_Real *objgrad)
Definition: nlpioracle.c:1959
#define SCIPreallocBlockMemoryArray(scip, ptr, oldnum, newnum)
Definition: scip_mem.h:90
SCIP_RETCODE SCIPexpriterInit(SCIP_EXPRITER *iterator, SCIP_EXPR *expr, SCIP_EXPRITER_TYPE type, SCIP_Bool allowrevisit)
Definition: expriter.c:491
const char * SCIPexprintGetName(void)
SCIP_RETCODE SCIPnlpiOracleChgLinearCoefs(SCIP *scip, SCIP_NLPIORACLE *oracle, int considx, int nentries, const int *varidxs, const SCIP_Real *newcoefs)
Definition: nlpioracle.c:1548
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:84
static SCIP_RETCODE ensureConsLinSize(SCIP *scip, SCIP_NLPIORACLECONS *cons, int minsize)
Definition: nlpioracle.c:142
methods to interpret (evaluate) an expression "fast"
SCIP_RETCODE SCIPprintExpr(SCIP *scip, SCIP_EXPR *expr, FILE *file)
Definition: scip_expr.c:1476
SCIP_RETCODE SCIPexprintCompile(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA **exprintdata)
#define SCIP_MAXSTRLEN
Definition: def.h:293
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValue(SCIP *scip, SCIP_NLPIORACLE *oracle, int considx, const SCIP_Real *x, SCIP_Real *conval)
Definition: nlpioracle.c:1903
const char * SCIPexprhdlrGetName(SCIP_EXPRHDLR *exprhdlr)
Definition: expr.c:525
#define SCIPallocClearBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:88
int SCIPcalcMemGrowSize(SCIP *scip, int num)
Definition: scip_mem.c:130
SCIP_RETCODE SCIPexprintFreeData(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA **exprintdata)
static SCIP_RETCODE freeConstraint(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS **cons, SCIP_Bool updatevarcount)
Definition: nlpioracle.c:391
SCIP_Real * lincoefs
Definition: nlpioracle.c:45
methods to store an NLP and request function, gradient, and Hessian values
SCIP_RETCODE SCIPstopClock(SCIP *scip, SCIP_CLOCK *clck)
Definition: scip_timing.c:169
SCIP_RETCODE SCIPnlpiOracleEvalHessianLag(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx_obj, SCIP_Bool isnewx_cons, SCIP_Real objfactor, const SCIP_Real *lambda, SCIP_Real *hessian)
Definition: nlpioracle.c:2371
SCIP_RETCODE SCIPexprintHessianSparsity(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata, SCIP_Real *varvals, int **rowidxs, int **colidxs, int *nnz)
#define FALSE
Definition: def.h:87
#define EPSEQ(x, y, eps)
Definition: def.h:202
static void printName(char *buffer, char *name, int idx, char prefix, const char *suffix, SCIP_Bool longnames)
Definition: nlpioracle.c:849
SCIP_RETCODE SCIPnlpiOracleAddConstraints(SCIP *scip, SCIP_NLPIORACLE *oracle, int nconss, const SCIP_Real *lhss, const SCIP_Real *rhss, const int *nlininds, int *const *lininds, SCIP_Real *const *linvals, SCIP_EXPR **exprs, const char **consnames)
Definition: nlpioracle.c:1158
void SCIPsetClockEnabled(SCIP_CLOCK *clck, SCIP_Bool enable)
Definition: scip_timing.c:182
SCIP_Real SCIPinfinity(SCIP *scip)
int SCIPsnprintf(char *t, int len, const char *s,...)
Definition: misc.c:10755
#define TRUE
Definition: def.h:86
#define SCIPdebug(x)
Definition: pub_message.h:84
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:54
static void freeVariables(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:503
SCIP_RETCODE SCIPnlpiOracleSetProblemName(SCIP *scip, SCIP_NLPIORACLE *oracle, const char *name)
Definition: nlpioracle.c:1036
static void clearDeletedLinearElements(int **linidxs, SCIP_Real **coefs, int *nidxs)
Definition: nlpioracle.c:557
SCIP_Real SCIPnlpiOracleGetEvalTime(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2434
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:99
#define SCIPdebugMessage
Definition: pub_message.h:87
SCIP_RETCODE SCIPnlpiOracleSetObjective(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real constant, int nlin, const int *lininds, const SCIP_Real *linvals, SCIP_EXPR *expr)
Definition: nlpioracle.c:1219
SCIP_RETCODE SCIPfreeClock(SCIP *scip, SCIP_CLOCK **clck)
Definition: scip_timing.c:118
SCIP_RETCODE SCIPexprintGrad(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata, SCIP_Real *varvals, SCIP_Bool new_varvals, SCIP_Real *val, SCIP_Real *gradient)
static void mapIndices(int *indexmap, int nindices, int *indices)
Definition: nlpioracle.c:537
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:127
static SCIP_RETCODE moveVariable(SCIP *scip, SCIP_NLPIORACLE *oracle, int fromidx, int toidx)
Definition: nlpioracle.c:464
unsigned int SCIP_EXPRINTCAPABILITY
#define BMSfreeMemory(ptr)
Definition: memory.h:138
void SCIPcaptureExpr(SCIP_EXPR *expr)
Definition: scip_expr.c:1399
SCIP_RETCODE SCIPnlpiOraclePrintProblem(SCIP *scip, SCIP_NLPIORACLE *oracle, FILE *file)
Definition: nlpioracle.c:2445
int * varnlcount
Definition: nlpioracle.c:64
SCIP_Real SCIPnlpiOracleGetConstraintRhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:1812
static void invalidateHessianLagSparsity(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:212
int * heslagoffsets
Definition: nlpioracle.c:75
#define SCIPdebugMsg
Definition: scip_message.h:69
SCIP_VAR ** x
Definition: circlepacking.c:54
void SCIPinfoMessage(SCIP *scip, FILE *file, const char *formatstr,...)
Definition: scip_message.c:199
static SCIP_RETCODE exprIsNonSmooth(SCIP *scip, SCIP_EXPR *expr, SCIP_Bool *nonsmooth)
Definition: nlpioracle.c:921
#define SCIP_DEFAULT_EPSILON
Definition: def.h:183
#define SCIPallocCleanBufferArray(scip, ptr, num)
Definition: scip_mem.h:133
SCIP_CLOCK * evalclock
Definition: nlpioracle.c:79
#define SCIPduplicateBlockMemoryArray(scip, ptr, source, num)
Definition: scip_mem.h:96
int * jacoffsets
Definition: nlpioracle.c:72
SCIP_RETCODE SCIPcreateClock(SCIP *scip, SCIP_CLOCK **clck)
Definition: scip_timing.c:67
static SCIP_RETCODE ensureIntArraySize(SCIP *scip, int **intarray, int *len, int minsize)
Definition: nlpioracle.c:168
handler for variable index expressions
static void sortLinearCoefficients(int *nidxs, int *idxs, SCIP_Real *coefs)
Definition: nlpioracle.c:279
#define SCIPerrorMessage
Definition: pub_message.h:55
int SCIPnlpiOracleGetNConstraints(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1717
int SCIPnlpiOracleGetNVars(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1707
SCIP_RETCODE SCIPnlpiOracleResetEvalTime(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2418
char * SCIPnlpiOracleGetConstraintName(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:1825
SCIP_RETCODE SCIPnlpiOracleChgVarBounds(SCIP *scip, SCIP_NLPIORACLE *oracle, int nvars, const int *indices, const SCIP_Real *lbs, const SCIP_Real *ubs)
Definition: nlpioracle.c:1248
static SCIP_RETCODE hessLagAddExpr(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_Real weight, const SCIP_Real *x, SCIP_Bool new_x, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata, int *hesoffset, int *hescol, SCIP_Real *values)
Definition: nlpioracle.c:782
SCIP_Real SCIPnlpiOracleGetObjectiveConstant(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1788
SCIP_RETCODE SCIPgetBoolParam(SCIP *scip, const char *name, SCIP_Bool *value)
Definition: scip_param.c:241
#define NULL
Definition: lpi_spx1.cpp:155
power and signed power expression handlers
#define REALABS(x)
Definition: def.h:201
SCIP_RETCODE SCIPnlpiOracleGetJacobianSparsity(SCIP *scip, SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2018
SCIP_RETCODE SCIPexprintCreate(SCIP *scip, SCIP_EXPRINT **exprint)
SCIP_RETCODE SCIPnlpiOracleAddVars(SCIP *scip, SCIP_NLPIORACLE *oracle, int nvars, const SCIP_Real *lbs, const SCIP_Real *ubs, const char **varnames)
Definition: nlpioracle.c:1072
#define SCIP_CALL(x)
Definition: def.h:384
#define SCIPensureBlockMemoryArray(scip, ptr, arraysizeptr, minsize)
Definition: scip_mem.h:98
SCIP_VAR * h
Definition: circlepacking.c:59
SCIP_RETCODE SCIPnlpiOracleChgExpr(SCIP *scip, SCIP_NLPIORACLE *oracle, int considx, SCIP_EXPR *expr)
Definition: nlpioracle.c:1644
SCIP_Real SCIPnlpiOracleGetConstraintLhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:1799
SCIP_Bool SCIPsortedvecFindInt(int *intarray, int val, int len, int *pos)
SCIP_RETCODE SCIPnlpiOracleEvalConstraintGradient(SCIP *scip, SCIP_NLPIORACLE *oracle, const int considx, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *conval, SCIP_Real *congrad)
Definition: nlpioracle.c:1988
SCIP_RETCODE SCIPcreateExpriter(SCIP *scip, SCIP_EXPRITER **iterator)
Definition: scip_expr.c:2300
static SCIP_RETCODE createConstraint(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS **cons, int nlinidxs, const int *linidxs, const SCIP_Real *lincoefs, SCIP_EXPR *expr, SCIP_Real lhs, SCIP_Real rhs, const char *name)
Definition: nlpioracle.c:328
SCIP_RETCODE SCIPexprintEval(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata, SCIP_Real *varvals, SCIP_Real *val)
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:115
SCIP_EXPR * expr
Definition: nlpioracle.c:47
#define SCIP_EXPRINTCAPABILITY_ALL
SCIP_RETCODE SCIPnlpiOracleEvalJacobian(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *convals, SCIP_Real *jacobi)
Definition: nlpioracle.c:2150
SCIP_EXPRINTCAPABILITY SCIPexprintGetExprCapability(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata)
SCIP_RETCODE SCIPnlpiOracleDelVarSet(SCIP *scip, SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1320
char ** varnames
Definition: nlpioracle.c:62
#define SCIP_Bool
Definition: def.h:84
SCIP_EXPR * SCIPexpriterRestartDFS(SCIP_EXPRITER *iterator, SCIP_EXPR *expr)
Definition: expriter.c:620
SCIP_NLPIORACLECONS * objective
Definition: nlpioracle.c:70
const char * SCIPnlpiOracleGetProblemName(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1060
SCIP_Real SCIPgetClockTime(SCIP *scip, SCIP_CLOCK *clck)
Definition: scip_timing.c:310
SCIP_RETCODE SCIPreleaseExpr(SCIP *scip, SCIP_EXPR **expr)
Definition: scip_expr.c:1407
SCIP_RETCODE SCIPexprintFree(SCIP *scip, SCIP_EXPRINT **exprint)
SCIP_RETCODE SCIPnlpiOracleCreate(SCIP *scip, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:974
SCIP_Real * varubs
Definition: nlpioracle.c:61
void SCIPsetIndexExprVaridx(SCIP_EXPR *expr, int newindex)
Definition: expr_varidx.c:268
struct SCIP_ExprInt SCIP_EXPRINT
SCIP_Real * varlbs
Definition: nlpioracle.c:60
const SCIP_Real * SCIPnlpiOracleGetVarLbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1727
SCIP_EXPR * SCIPexpriterGetNext(SCIP_EXPRITER *iterator)
Definition: expriter.c:848
#define EPSLE(x, y, eps)
Definition: def.h:204
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:127
SCIP_EXPRHDLR * SCIPexprGetHdlr(SCIP_EXPR *expr)
Definition: expr.c:3821
SCIP_Bool SCIPnlpiOracleIsConstraintNonlinear(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:1838
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
SCIP_RETCODE SCIPexprintHessian(SCIP *scip, SCIP_EXPRINT *exprint, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata, SCIP_Real *varvals, SCIP_Bool new_varvals, SCIP_Real *val, int **rowidxs, int **colidxs, SCIP_Real **hessianvals, int *nnz)
#define BMSclearMemory(ptr)
Definition: memory.h:122
const SCIP_Real * SCIPnlpiOracleGetVarUbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1737
SCIP_RETCODE SCIPnlpiOracleFree(SCIP *scip, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1004
static SCIP_RETCODE evalFunctionGradient(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS *cons, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *RESTRICT val, SCIP_Real *RESTRICT grad)
Definition: nlpioracle.c:648
void SCIPfreeExpriter(SCIP_EXPRITER **iterator)
Definition: scip_expr.c:2314
static SCIP_RETCODE updateVariableCounts(SCIP *scip, SCIP_NLPIORACLE *oracle, int factor, int nlinidxs, const int *linidxs, SCIP_EXPR *expr)
Definition: nlpioracle.c:234
SCIP_Bool SCIPnlpiOracleIsVarNonlinear(SCIP *scip, SCIP_NLPIORACLE *oracle, int varidx)
Definition: nlpioracle.c:1757
static SCIP_RETCODE ensureVarsSize(SCIP *scip, SCIP_NLPIORACLE *oracle, int minsize)
Definition: nlpioracle.c:93
static SCIP_RETCODE hessLagSparsitySetNzFlagForExpr(SCIP *scip, SCIP_NLPIORACLE *oracle, int **colnz, int *collen, int *colnnz, int *nzcount, SCIP_EXPR *expr, SCIP_EXPRINTDATA *exprintdata, int dim)
Definition: nlpioracle.c:722
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveValue(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *objval)
Definition: nlpioracle.c:1878
static SCIP_RETCODE freeConstraints(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:435
SCIP_RETCODE SCIPresetClock(SCIP *scip, SCIP_CLOCK *clck)
Definition: scip_timing.c:135
SCIP_RETCODE SCIPnlpiOraclePrintProblemGams(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_Real *initval, FILE *file)
Definition: nlpioracle.c:2513
SCIP_Bool SCIPisExprSignpower(SCIP *scip, SCIP_EXPR *expr)
Definition: expr_pow.c:3215
SCIP_Bool SCIPisExprVaridx(SCIP *scip, SCIP_EXPR *expr)
Definition: expr_varidx.c:242
static void invalidateJacobiSparsity(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:188
static SCIP_RETCODE evalFunctionValue(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS *cons, const SCIP_Real *x, SCIP_Real *val)
Definition: nlpioracle.c:595
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValues(SCIP *scip, SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *convals)
Definition: nlpioracle.c:1927
#define SCIP_Real
Definition: def.h:177
SCIP_RETCODE SCIPnlpiOracleGetHessianLagSparsity(SCIP *scip, SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2277
#define SCIPfreeCleanBufferArray(scip, ptr)
Definition: scip_mem.h:137
SCIP_EXPRINTCAPABILITY SCIPnlpiOracleGetEvalCapability(SCIP *scip, SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1855
SCIP_EXPRINT * exprinterpreter
Definition: nlpioracle.c:78
SCIP_EXPRINTDATA * exprintdata
Definition: nlpioracle.c:48
#define SCIPisFinite(x)
Definition: pub_misc.h:1892
static SCIP_RETCODE ensureConssSize(SCIP *scip, SCIP_NLPIORACLE *oracle, int minsize)
Definition: nlpioracle.c:126
#define SCIPallocMemory(scip, ptr)
Definition: scip_mem.h:51
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:102
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:123
int * varlincount
Definition: nlpioracle.c:63
#define SCIPallocClearBlockMemory(scip, ptr)
Definition: scip_mem.h:82
SCIP_Bool SCIPexpriterIsEnd(SCIP_EXPRITER *iterator)
Definition: expriter.c:959
char ** SCIPnlpiOracleGetVarNames(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1747
int SCIPgetIndexExprVaridx(SCIP_EXPR *expr)
Definition: expr_varidx.c:257
struct SCIP_ExprIntData SCIP_EXPRINTDATA
SCIP_RETCODE SCIPnlpiOracleChgConsSides(SCIP *scip, SCIP_NLPIORACLE *oracle, int nconss, const int *indices, const SCIP_Real *lhss, const SCIP_Real *rhss)
Definition: nlpioracle.c:1285
SCIP_RETCODE SCIPstartClock(SCIP *scip, SCIP_CLOCK *clck)
Definition: scip_timing.c:152
SCIP_RETCODE SCIPnlpiOracleDelConsSet(SCIP *scip, SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1462
static SCIP_RETCODE printFunction(SCIP *scip, SCIP_NLPIORACLE *oracle, FILE *file, SCIP_NLPIORACLECONS *cons, SCIP_Bool longvarnames)
Definition: nlpioracle.c:884
void SCIPsortedvecInsertInt(int *intarray, int keyval, int *len, int *pos)
int * heslagcols
Definition: nlpioracle.c:76
SCIP callable library.
void SCIPnlpiOracleGetVarCounts(SCIP *scip, SCIP_NLPIORACLE *oracle, const int **lincounts, const int **nlcounts)
Definition: nlpioracle.c:1772
SCIP_NLPIORACLECONS ** conss
Definition: nlpioracle.c:68
void SCIPsortIntReal(int *intarray, SCIP_Real *realarray, int len)
SCIP_RETCODE SCIPnlpiOracleChgObjConstant(SCIP *scip, SCIP_NLPIORACLE *oracle, SCIP_Real objconstant)
Definition: nlpioracle.c:1690