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-2016 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not email to scip@zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file nlpioracle.c
17  * @brief implementation of NLPI oracle interface
18  * @author Stefan Vigerske
19  *
20  * @todo jacobi evaluation should be sparse
21  */
22 
23 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
24 
25 #include "scip/pub_message.h"
26 #include "scip/pub_misc.h"
27 #include "nlpi/nlpioracle.h"
28 #include "nlpi/pub_expr.h"
29 #include "nlpi/exprinterpret.h"
30 
31 #include <string.h> /* for strlen */
32 
33 /**@name NLPI Oracle data structures */
34 /**@{ */
35 
37 {
38  SCIP_Real lhs; /**< left hand side (for constraint) or constant (for objective) */
39  SCIP_Real rhs; /**< right hand side (for constraint) or constant (for objective) */
40 
41  int linsize; /**< length of linidxs and lincoefs arrays */
42  int nlinidxs; /**< number of linear variable indices and coefficients */
43  int* linidxs; /**< variable indices in linear part, or NULL if none */
44  SCIP_Real* lincoefs; /**< variable coefficients in linear part, of NULL if none */
45 
46  int quadsize; /**< length of quadelems array */
47  int nquadelems; /**< number of quadratic elements */
48  SCIP_QUADELEM* quadelems; /**< quadratic elements, or NULL if none */
49 
50  int* exprvaridxs; /**< indices of variables in expression tree, or NULL if no exprtree */
51  SCIP_EXPRTREE* exprtree; /**< expression tree for nonlinear part, or NULL if none */
52 
53  char* name; /**< name of constraint */
54 };
56 
58 {
59  BMS_BLKMEM* blkmem; /**< block memory */
60  SCIP_Real infinity; /**< value for infinity */
61  char* name; /**< name of problem */
62 
63  int varssize; /**< length of variables related arrays */
64  int nvars; /**< number of variables */
65  SCIP_Real* varlbs; /**< array with variable lower bounds */
66  SCIP_Real* varubs; /**< array with variable upper bounds */
67  char** varnames; /**< array with variable names */
68  int* vardegrees; /**< array with maximal degree of variable over objective and all constraints */
69  SCIP_Bool vardegreesuptodate; /**< whether the variable degrees are up to date */
70 
71  int consssize; /**< length of constraints related arrays */
72  int nconss; /**< number of constraints */
73  SCIP_NLPIORACLECONS** conss; /**< constraints, or NULL if none */
74 
75  SCIP_NLPIORACLECONS* objective; /**< objective */
76 
77  int* jacoffsets; /**< rowwise jacobi sparsity pattern: constraint offsets in jaccols */
78  int* jaccols; /**< rowwise jacobi sparsity pattern: indices of variables appearing in constraints */
79 
80  int* heslagoffsets; /**< rowwise sparsity pattern of hessian matrix of Lagrangian: row offsets in heslagcol */
81  int* heslagcols; /**< rowwise sparsity pattern of hessian matrix of Lagrangian: column indices; sorted for each row */
82 
83 
84  SCIP_EXPRINT* exprinterpreter; /**< interpreter for expression trees: evaluation and derivatives */
85 };
86 
87 /**@} */
88 
89 /**@name Local functions */
90 /**@{ */
91 
92 /** calculate memory size for dynamically allocated arrays (copied from scip/set.c) */
93 static
95  int num /**< minimum number of entries to store */
96  )
97 {
98  int size;
99 
100  /* calculate the size with this loop, such that the resulting numbers are always the same (-> block memory) */
101  size = 4;
102  while( size < num )
103  size = (int)(1.2 * size + 4);
104 
105  return size;
106 }
107 
108 /** ensures that variables related arrays in oracle have at least a given length */
109 static
111  SCIP_NLPIORACLE* oracle, /**< NLPIORACLE data structure */
112  int minsize /**< minimal required size */
113  )
114 {
115  assert(oracle != NULL);
116 
117  if( minsize > oracle->varssize )
118  {
119  int newsize;
120 
121  newsize = calcGrowSize(minsize);
122  assert(newsize >= minsize);
123 
124  SCIP_ALLOC( BMSreallocBlockMemoryArray(oracle->blkmem, &oracle->varlbs, oracle->varssize, newsize) );
125  SCIP_ALLOC( BMSreallocBlockMemoryArray(oracle->blkmem, &oracle->varubs, oracle->varssize, newsize) );
126  if( oracle->varnames != NULL )
127  {
128  SCIP_ALLOC( BMSreallocBlockMemoryArray(oracle->blkmem, &oracle->varnames, oracle->varssize, newsize) );
129  }
130  SCIP_ALLOC( BMSreallocBlockMemoryArray(oracle->blkmem, &oracle->vardegrees, oracle->varssize, newsize) );
131 
132  oracle->varssize = newsize;
133  }
134  assert(oracle->varssize >= minsize);
135 
136  return SCIP_OKAY;
137 }
138 
139 /** ensures that constraints array in oracle has at least a given length */
140 static
142  SCIP_NLPIORACLE* oracle, /**< NLPIORACLE data structure */
143  int minsize /**< minimal required size */
144  )
145 {
146  assert(oracle != NULL);
147 
148  if( minsize > oracle->consssize )
149  {
150  int newsize;
151 
152  newsize = calcGrowSize(minsize);
153  assert(newsize >= minsize);
154 
155  SCIP_ALLOC( BMSreallocBlockMemoryArray(oracle->blkmem, &oracle->conss, oracle->consssize, newsize) );
156  oracle->consssize = newsize;
157  }
158  assert(oracle->consssize >= minsize);
159 
160  return SCIP_OKAY;
161 }
162 
163 /** ensures that arrays for linear part in a oracle constraints have at least a given length */
164 static
166  BMS_BLKMEM* blkmem, /**< block memory */
167  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
168  int minsize /**< minimal required size */
169  )
170 {
171  assert(blkmem != NULL);
172  assert(cons != NULL);
173 
174  if( minsize > cons->linsize )
175  {
176  int newsize;
177 
178  newsize = calcGrowSize(minsize);
179  assert(newsize >= minsize);
180 
181  SCIP_ALLOC( BMSreallocBlockMemoryArray(blkmem, &cons->linidxs, cons->linsize, newsize) );
182  SCIP_ALLOC( BMSreallocBlockMemoryArray(blkmem, &cons->lincoefs, cons->linsize, newsize) );
183  cons->linsize = newsize;
184  }
185  assert(cons->linsize >= minsize);
186 
187  return SCIP_OKAY;
188 }
189 
190 /** ensures that arrays for quadratic part in a oracle constraints have at least a given length */
191 static
193  BMS_BLKMEM* blkmem, /**< block memory */
194  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
195  int minsize /**< minimal required size */
196  )
197 {
198  assert(blkmem != NULL);
199  assert(cons != NULL);
200 
201  if( minsize > cons->quadsize )
202  {
203  int newsize;
204 
205  newsize = calcGrowSize(minsize);
206  assert(newsize >= minsize);
207 
208  SCIP_ALLOC( BMSreallocBlockMemoryArray(blkmem, &cons->quadelems, cons->quadsize, newsize) );
209  cons->quadsize = newsize;
210  }
211  assert(cons->quadsize >= minsize);
212 
213  return SCIP_OKAY;
214 }
215 
216 /** ensures that a given array of integers has at least a given length */
217 static
219  BMS_BLKMEM* blkmem, /**< block memory */
220  int** intarray, /**< array of integers */
221  int* len, /**< length of array (modified if reallocated) */
222  int minsize /**< minimal required array length */
223  )
224 {
225  assert(blkmem != NULL);
226  assert(intarray != NULL);
227  assert(len != NULL);
228 
229  if( minsize > *len )
230  {
231  int newsize;
232 
233  newsize = calcGrowSize(minsize);
234  assert(newsize >= minsize);
235 
236  SCIP_ALLOC( BMSreallocBlockMemoryArray(blkmem, intarray, *len, newsize) );
237  *len = newsize;
238  }
239  assert(*len >= minsize);
240 
241  return SCIP_OKAY;
242 }
243 
244 /** Invalidates the sparsity pattern of the Jacobian.
245  * Should be called when constraints are added or deleted.
246  */
247 static
249  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
250  )
251 {
252  assert(oracle != NULL);
253 
254  SCIPdebugMessage("%p invalidate jacobian sparsity\n", (void*)oracle);
255 
256  if( oracle->jacoffsets == NULL )
257  { /* nothing to do */
258  assert(oracle->jaccols == NULL);
259  return;
260  }
261 
262  assert(oracle->jaccols != NULL);
263  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->jaccols, oracle->jacoffsets[oracle->nconss]);
264  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->jacoffsets, oracle->nconss + 1);
265 }
266 
267 /** Invalidates the sparsity pattern of the Hessian of the Lagragian.
268  * Should be called when the objective is set or constraints are added or deleted.
269  */
270 static
272  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
273  )
274 {
275  assert(oracle != NULL);
276 
277  SCIPdebugMessage("%p invalidate hessian lag sparsity\n", (void*)oracle);
278 
279  if( oracle->heslagoffsets == NULL )
280  { /* nothing to do */
281  assert(oracle->heslagcols == NULL);
282  return;
283  }
284 
285  assert(oracle->heslagcols != NULL);
286  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->heslagcols, oracle->heslagoffsets[oracle->nvars]);
287  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->heslagoffsets, oracle->nvars + 1);
288 }
289 
290 /** sorts a linear term, merges duplicate entries and removes entries with coefficient 0.0 */
291 static
293  int* nidxs, /**< number of variables */
294  int* idxs, /**< indices of variables */
295  SCIP_Real* coefs /**< coefficients of variables */
296  )
297 {
298  int offset;
299  int j;
300 
301  assert(nidxs != NULL);
302  assert(idxs != NULL || *nidxs == 0);
303  assert(coefs != NULL || *nidxs == 0);
304 
305  if( *nidxs == 0 )
306  return;
307 
308  SCIPsortIntReal(idxs, coefs, *nidxs);
309 
310  offset = 0;
311  j = 0;
312  while( j+offset < *nidxs )
313  {
314  assert(idxs[j] >= 0); /*lint !e613*/
315 
316  /* move j+offset to j, if different */
317  if( offset > 0 )
318  {
319  idxs[j] = idxs[j+offset]; /*lint !e613*/
320  coefs[j] = coefs[j+offset]; /*lint !e613*/
321  }
322 
323  /* add up coefs for j+offset+1... as long as they have the same index */
324  while( j+offset+1 < *nidxs && idxs[j] == idxs[j+offset+1] ) /*lint !e613*/
325  {
326  coefs[j] += coefs[j+offset+1]; /*lint !e613*/
327  ++offset;
328  }
329 
330  /* if j'th element is 0, increase offset, otherwise increase j */
331  if( coefs[j] == 0.0 ) /*lint !e613*/
332  ++offset;
333  else
334  ++j;
335  }
336  *nidxs -= offset;
337 }
338 
339 /** creates a NLPI constraint from given constraint data */
340 static
342  BMS_BLKMEM* blkmem, /**< block memory */
343  SCIP_NLPIORACLECONS** cons, /**< buffer where to store pointer to constraint */
344  int nlinidxs, /**< length of linear part */
345  const int* linidxs, /**< indices of linear part, or NULL if nlinidxs == 0 */
346  const SCIP_Real* lincoefs, /**< coefficients of linear part, or NULL if nlinidxs == 0 */
347  int nquadelems, /**< lenght of quadratic part */
348  const SCIP_QUADELEM* quadelems, /**< quadratic elements, or NULL if nquadelems == 0 */
349  const int* exprvaridxs, /**< indicies of variables in expression tree, or NULL if exprtree == NULL */
350  const SCIP_EXPRTREE* exprtree, /**< expression tree, or NULL */
351  SCIP_Real lhs, /**< left-hand-side of constraint */
352  SCIP_Real rhs, /**< right-hand-side of constraint */
353  const char* name /**< name of constraint, or NULL */
354  )
355 {
356  assert(blkmem != NULL);
357  assert(cons != NULL);
358  assert(nlinidxs >= 0);
359  assert(linidxs != NULL || nlinidxs == 0);
360  assert(lincoefs != NULL || nlinidxs == 0);
361  assert(nquadelems >= 0);
362  assert(quadelems != NULL || nquadelems == 0);
363  assert(exprvaridxs != NULL || exprtree == NULL);
364  assert(EPSLE(lhs, rhs, SCIP_DEFAULT_EPSILON));
365 
366  SCIP_ALLOC( BMSallocBlockMemory(blkmem, cons) );
367  assert(*cons != NULL);
368  BMSclearMemory(*cons);
369 
370  if( nlinidxs > 0 )
371  {
372  SCIP_ALLOC( BMSduplicateBlockMemoryArray(blkmem, &(*cons)->linidxs, linidxs, nlinidxs) );
373  SCIP_ALLOC( BMSduplicateBlockMemoryArray(blkmem, &(*cons)->lincoefs, lincoefs, nlinidxs) );
374  (*cons)->linsize = nlinidxs;
375  (*cons)->nlinidxs = nlinidxs;
376 
377  /* sort, merge duplicates, remove zero's */
378  sortLinearCoefficients(&(*cons)->nlinidxs, (*cons)->linidxs, (*cons)->lincoefs);
379  assert((*cons)->linidxs[0] >= 0);
380  }
381 
382  if( nquadelems > 0 )
383  {
384  SCIP_ALLOC( BMSduplicateBlockMemoryArray(blkmem, &(*cons)->quadelems, quadelems, nquadelems) );
385  (*cons)->nquadelems = nquadelems;
386  (*cons)->quadsize = nquadelems;
387 
388  /* sort and squeeze quadratic part */
389  SCIPquadelemSort((*cons)->quadelems, nquadelems);
390  SCIPquadelemSqueeze((*cons)->quadelems, nquadelems, &(*cons)->nquadelems);
391  assert((*cons)->nquadelems == 0 || (*cons)->quadelems[0].idx1 >= 0);
392  assert((*cons)->nquadelems == 0 || (*cons)->quadelems[0].idx2 >= 0);
393  }
394 
395  if( exprtree != NULL )
396  {
397  SCIP_CALL( SCIPexprtreeCopy(blkmem, &(*cons)->exprtree, (SCIP_EXPRTREE*)exprtree) );
398 
399  SCIP_ALLOC( BMSduplicateBlockMemoryArray(blkmem, &(*cons)->exprvaridxs, exprvaridxs, SCIPexprtreeGetNVars((SCIP_EXPRTREE*)exprtree)) );
400  }
401 
402  if( lhs > rhs )
403  {
404  assert(EPSEQ(lhs, rhs, SCIP_DEFAULT_EPSILON));
405  lhs = rhs;
406  }
407  (*cons)->lhs = lhs;
408  (*cons)->rhs = rhs;
409 
410  if( name != NULL )
411  {
412  SCIP_ALLOC( BMSduplicateBlockMemoryArray(blkmem, &(*cons)->name, name, strlen(name)+1) );
413  }
414 
415  return SCIP_OKAY;
416 }
417 
418 /** frees a constraint */
419 static
421  BMS_BLKMEM* blkmem, /**< block memory */
422  SCIP_NLPIORACLECONS** cons /**< pointer to constraint that should be freed */
423  )
424 {
425  assert(blkmem != NULL);
426  assert(cons != NULL);
427  assert(*cons != NULL);
428 
429  SCIPdebugMessage("free constraint %p\n", (void*)*cons);
430 
431  BMSfreeBlockMemoryArrayNull(blkmem, &(*cons)->linidxs, (*cons)->linsize);
432  BMSfreeBlockMemoryArrayNull(blkmem, &(*cons)->lincoefs, (*cons)->linsize);
433 
434  BMSfreeBlockMemoryArrayNull(blkmem, &(*cons)->quadelems, (*cons)->quadsize);
435 
436  if( (*cons)->exprtree != NULL )
437  {
438  BMSfreeBlockMemoryArrayNull(blkmem, &(*cons)->exprvaridxs, SCIPexprtreeGetNVars((*cons)->exprtree));
439  SCIP_CALL_ABORT( SCIPexprtreeFree(&(*cons)->exprtree) );
440  }
441 
442  if( (*cons)->name != NULL )
443  {
444  BMSfreeBlockMemoryArrayNull(blkmem, &(*cons)->name, strlen((*cons)->name)+1);
445  }
446 
447  BMSfreeBlockMemory(blkmem, cons);
448  assert(*cons == NULL);
449 }
450 
451 /** frees all constraints */
452 static
454  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
455  )
456 {
457  int i;
458 
459  assert(oracle != NULL);
460 
461  SCIPdebugMessage("%p free constraints\n", (void*)oracle);
462 
463  for( i = 0; i < oracle->nconss; ++i )
464  {
465  freeConstraint(oracle->blkmem, &oracle->conss[i]);
466  assert(oracle->conss[i] == NULL);
467  }
468  oracle->nconss = 0;
469 
470  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &oracle->conss, oracle->consssize);
471  oracle->consssize = 0;
472 }
473 
474 /** moves one variable
475  * The place where it moves to need to be empty (all NULL) but allocated.
476  * Note that this function does not update the variable indices in the constraints!
477  */
478 static
480  SCIP_NLPIORACLE* oracle, /**< pointer to store NLPIORACLE data structure */
481  int fromidx, /**< index of variable to move */
482  int toidx /**< index of place where to move variable to */
483  )
484 {
485  assert(oracle != NULL);
486 
487  SCIPdebugMessage("%p move variable\n", (void*)oracle);
488 
489  assert(0 <= fromidx);
490  assert(0 <= toidx);
491  assert(fromidx < oracle->nvars);
492  assert(toidx < oracle->nvars);
493 
494  assert(oracle->varlbs[toidx] <= -oracle->infinity);
495  assert(oracle->varubs[toidx] >= oracle->infinity);
496  assert(oracle->varnames == NULL || oracle->varnames[toidx] == NULL);
497  assert(!oracle->vardegreesuptodate || oracle->vardegrees[toidx] == -1);
498 
499  oracle->varlbs[toidx] = oracle->varlbs[fromidx];
500  oracle->varubs[toidx] = oracle->varubs[fromidx];
501 
502  oracle->varlbs[fromidx] = -oracle->infinity;
503  oracle->varubs[fromidx] = oracle->infinity;
504 
505  oracle->vardegrees[toidx] = oracle->vardegrees[fromidx];
506  oracle->vardegrees[fromidx] = -1;
507 
508  if( oracle->varnames != NULL )
509  {
510  oracle->varnames[toidx] = oracle->varnames[fromidx];
511  oracle->varnames[fromidx] = NULL;
512  }
513 
514  return SCIP_OKAY;
515 }
516 
517 /** frees all variables */
518 static
520  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
521  )
522 {
523  int i;
524 
525  assert(oracle != NULL);
526 
527  SCIPdebugMessage("%p free variables\n", (void*)oracle);
528 
529  if( oracle->varnames != NULL )
530  {
531  for( i = 0; i < oracle->nvars; ++i )
532  {
533  if( oracle->varnames[i] != NULL )
534  {
535  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->varnames[i], strlen(oracle->varnames[i])+1); /*lint !e866*/
536  }
537  }
538  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &oracle->varnames, oracle->varssize);
539  }
540  oracle->nvars = 0;
541  oracle->vardegreesuptodate = TRUE;
542 
543  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &oracle->varlbs, oracle->varssize);
544  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &oracle->varubs, oracle->varssize);
545  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &oracle->vardegrees, oracle->varssize);
546 
547  oracle->varssize = 0;
548 }
549 
550 /** increases variable degrees in oracle w.r.t. variables occuring in a single constraint */
551 static
553  SCIP_NLPIORACLE* oracle, /**< oracle data structure */
554  SCIP_NLPIORACLECONS* cons /**< oracle constraint */
555  )
556 {
557  int j;
558 
559  assert(oracle != NULL);
560  assert(oracle->vardegrees != NULL);
561  assert(cons != NULL);
562 
563  for( j = 0; j < cons->nlinidxs; ++j )
564  if( oracle->vardegrees[cons->linidxs[j]] < 1 )
565  oracle->vardegrees[cons->linidxs[j]] = 1;
566 
567  for( j = 0; j < cons->nquadelems; ++j )
568  {
569  if( oracle->vardegrees[cons->quadelems[j].idx1] < 2 )
570  oracle->vardegrees[cons->quadelems[j].idx1] = 2;
571 
572  if( oracle->vardegrees[cons->quadelems[j].idx2] < 2 )
573  oracle->vardegrees[cons->quadelems[j].idx2] = 2;
574  }
575 
576  /* we could use exprtreeGetDegree to get actual degree of a variable in tree,
577  * but so far no solver could make use of this information */
578  if( cons->exprtree != NULL )
579  for( j = SCIPexprtreeGetNVars(cons->exprtree)-1; j >= 0; --j )
580  oracle->vardegrees[cons->exprvaridxs[j]] = INT_MAX;
581 }
582 
583 /** Updates the degrees of all variables. */
584 static
586  SCIP_NLPIORACLE* oracle /**< pointer to store NLPIORACLE data structure */
587  )
588 {
589  int c;
590 
591  assert(oracle != NULL);
592  assert(oracle->nvars == 0 || oracle->vardegrees != NULL);
593  assert(oracle->objective != NULL);
594 
595  SCIPdebugMessage("%p update variable degrees\n", (void*)oracle);
596 
597  if( oracle->vardegreesuptodate || oracle->nvars == 0 )
598  return;
599 
600  /* assume all variables do not appear in NLP */
601  BMSclearMemoryArray(oracle->vardegrees, oracle->nvars);
602 
603  updateVariableDegreesCons(oracle, oracle->objective);
604  for( c = 0; c < oracle->nconss; ++c )
605  updateVariableDegreesCons(oracle, oracle->conss[c]);
606 
607  oracle->vardegreesuptodate = TRUE;
608 }
609 
610 /** applies a mapping of indices to one array of indices */
611 static
613  int* indexmap, /**< mapping from old variable indices to new indices */
614  int nindices, /**< number of indices in indices1 and indices2 */
615  int* indices /**< array of indices to adjust */
616  )
617 {
618  assert(indexmap != NULL);
619  assert(nindices == 0 || indices != NULL);
620 
621  for( ; nindices ; --nindices, ++indices )
622  {
623  assert(indexmap[*indices] >= 0);
624  *indices = indexmap[*indices];
625  }
626 }
627 
628 /** removes entries with index -1 (marked as deleted) from array of linear elements
629  * assumes that array is sorted by index, i.e., all -1 are at the beginning
630  */
631 static
633  BMS_BLKMEM* blkmem, /**< block memory */
634  int** linidxs, /**< variable indices */
635  SCIP_Real** coefs, /**< variable coefficients */
636  int* nidxs /**< number of indices */
637  )
638 {
639  int i;
640  int offset;
641 
642  SCIPdebugMessage("clear deleted linear elements\n");
643 
644  assert(blkmem != NULL);
645  assert(linidxs != NULL);
646  assert(*linidxs != NULL);
647  assert(coefs != NULL);
648  assert(*coefs != NULL);
649  assert(nidxs != NULL);
650  assert(*nidxs > 0);
651 
652  /* search for beginning of non-delete entries @todo binary search? */
653  for( offset = 0; offset < *nidxs; ++offset )
654  if( (*linidxs)[offset] >= 0 )
655  break;
656 
657  /* nothing was deleted */
658  if( offset == 0 )
659  return;
660 
661  /* some or all elements were deleted -> move remaining ones front */
662  for( i = 0; i < *nidxs - offset; ++i )
663  {
664  (*linidxs)[i] = (*linidxs)[i+offset];
665  (*coefs)[i] = (*coefs) [i+offset];
666  }
667  *nidxs -= offset;
668 }
669 
670 /** removes entries with index pair (-1,-1) (marked as deleted) from array of quadratic elements
671  * assumes that array is sorted, i.e., all deleted elements are at the beginning
672  */
673 static
675  BMS_BLKMEM* blkmem, /**< block memory */
676  SCIP_QUADELEM** quadelems, /**< quadratic elements */
677  int* nquadelems /**< number of quadratic elements */
678  )
679 {
680  int i;
681  int offset;
682 
683  SCIPdebugMessage("clear deleted quad elements\n");
684 
685  assert(blkmem != NULL);
686  assert(quadelems != NULL);
687  assert(*quadelems != NULL);
688  assert(nquadelems != NULL);
689  assert(*nquadelems > 0);
690 
691  /* search for beginning of non-delete entries @todo binary search? */
692  for( offset = 0; offset < *nquadelems; ++offset )
693  {
694  /* either both variables are marked as deleted or none of them */
695  assert(((*quadelems)[offset].idx1 >= 0) == ((*quadelems)[offset].idx2 >= 0));
696  if( (*quadelems)[offset].idx1 >= 0 )
697  break;
698  }
699 
700  /* nothing was deleted */
701  if( offset == 0 )
702  return;
703 
704  /* some or all elements were deleted -> move remaining ones front */
705  for( i = 0; i < *nquadelems - offset; ++i )
706  (*quadelems)[i] = (*quadelems)[i+offset];
707  *nquadelems -= offset;
708 }
709 
710 /** applies a mapping of indices to an array of quadratic elements */
711 static
713  int* indexmap, /**< mapping from old variable indices to new indices */
714  int nelems, /**< number of quadratic elements */
715  SCIP_QUADELEM* elems /**< array of quadratic elements to adjust */
716  )
717 {
718  assert(indexmap != NULL);
719  assert(nelems == 0 || elems != NULL);
720 
721  for( ; nelems ; --nelems, ++elems )
722  {
723  assert(indexmap[elems->idx1] >= 0);
724  assert(indexmap[elems->idx2] >= 0);
725  elems->idx1 = indexmap[elems->idx1];
726  elems->idx2 = indexmap[elems->idx2];
727  /* swap indices if not idx1 <= idx2 */
728  if( elems->idx1 > elems->idx2 )
729  {
730  int tmp = elems->idx2;
731  elems->idx2 = elems->idx1;
732  elems->idx1 = tmp;
733  }
734  }
735 }
736 
737 /** computes the value of a function */
738 static
740  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
741  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
742  const SCIP_Real* x, /**< the point where to evaluate */
743  SCIP_Real* val /**< pointer to store function value */
744  )
745 { /*lint --e{715}*/
746  assert(oracle != NULL);
747  assert(cons != NULL);
748  assert(x != NULL || oracle->nvars == 0);
749  assert(val != NULL);
750 
751  SCIPdebugMessage("%p eval function value\n", (void*)oracle);
752 
753  *val = 0.0;
754 
755  if( cons->nlinidxs > 0 )
756  {
757  int* linidxs;
759  int nlin;
760 
761  nlin = cons->nlinidxs;
762  linidxs = cons->linidxs;
763  lincoefs = cons->lincoefs;
764  assert(linidxs != NULL);
765  assert(lincoefs != NULL);
766  assert(x != NULL);
767 
768  for( ; nlin > 0; --nlin, ++linidxs, ++lincoefs )
769  *val += *lincoefs * x[*linidxs];
770  }
771 
772  if( cons->nquadelems > 0 )
773  {
775  int nquadelems;
776 
777  quadelems = cons->quadelems;
778  nquadelems = cons->nquadelems;
779  assert(quadelems != NULL);
780  assert(x != NULL);
781 
782  for( ; nquadelems > 0; --nquadelems, ++quadelems )
783  *val += quadelems->coef * x[quadelems->idx1] * x[quadelems->idx2];
784  }
785 
786  if( cons->exprtree != NULL )
787  {
788  SCIP_Real* xx;
789  int i;
790  SCIP_Real nlval;
791  int nvars;
792 
793  nvars = SCIPexprtreeGetNVars(cons->exprtree);
794 
795  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &xx, nvars) );
796  for( i = 0; i < nvars; ++i )
797  {
798  assert(cons->exprvaridxs[i] >= 0);
799  xx[i] = x[cons->exprvaridxs[i]]; /*lint !e613 !e644*/
800  }
801 
802  SCIP_CALL( SCIPexprintEval(oracle->exprinterpreter, cons->exprtree, xx, &nlval) );
803  if( nlval != nlval || ABS(nlval) >= oracle->infinity ) /*lint !e777*/
804  *val = nlval;
805  else
806  *val += nlval;
807 
808  BMSfreeBlockMemoryArray(oracle->blkmem, &xx, nvars);
809  }
810 
811  return SCIP_OKAY;
812 }
813 
814 /** computes the value and gradient of a function */
815 static
817  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
818  SCIP_NLPIORACLECONS* cons, /**< oracle constraint */
819  const SCIP_Real* x, /**< the point where to evaluate */
820  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
821  SCIP_Real* val, /**< pointer to store function value */
822  SCIP_Real* grad /**< pointer to store function gradient */
823  )
824 { /*lint --e{715}*/
825  assert(oracle != NULL);
826  assert(x != NULL || oracle->nvars == 0);
827  assert(val != NULL);
828  assert(grad != NULL);
829 
830  SCIPdebugMessage("%p eval function gradient\n", (void*)oracle);
831 
832  *val = 0.0;
833  BMSclearMemoryArray(grad, oracle->nvars);
834 
835  if( cons->nlinidxs > 0 )
836  {
837  int* linidxs;
839  int nlin;
840 
841  nlin = cons->nlinidxs;
842  linidxs = cons->linidxs;
843  lincoefs = cons->lincoefs;
844  assert(linidxs != NULL);
845  assert(lincoefs != NULL);
846  assert(x != NULL);
847 
848  for( ; nlin > 0; --nlin, ++linidxs, ++lincoefs )
849  {
850  *val += *lincoefs * x[*linidxs];
851  assert(grad[*linidxs] == 0.0); /* we do not like duplicate indices */
852  grad[*linidxs] = *lincoefs;
853  }
854  }
855 
856  if( cons->nquadelems > 0 )
857  {
858  SCIP_Real tmp;
860  int nquadelems;
861 
862  quadelems = cons->quadelems;
863  nquadelems = cons->nquadelems;
864  assert(quadelems != NULL);
865  assert(x != NULL);
866 
867  for( ; nquadelems > 0; --nquadelems, ++quadelems )
868  {
869  tmp = quadelems->coef * x[quadelems->idx1];
870  *val += tmp * x[quadelems->idx2];
871  grad[quadelems->idx2] += tmp;
872  grad[quadelems->idx1] += quadelems->coef * x[quadelems->idx2];
873  }
874  }
875 
876  if( cons->exprtree != NULL )
877  {
878  SCIP_Real* xx;
879  SCIP_Real* g;
880  int i;
881  SCIP_Real nlval;
882  int nvars;
883 
884  xx = NULL;
885  nvars = SCIPexprtreeGetNVars(cons->exprtree);
886 
887  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &g, nvars) );
888 
889  if( isnewx )
890  {
891  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &xx, nvars) );
892  for( i = 0; i < nvars; ++i )
893  {
894  assert(cons->exprvaridxs[i] >= 0);
895  xx[i] = x[cons->exprvaridxs[i]]; /*lint !e613*/
896  }
897  }
898 
899  SCIPdebugMessage("eval gradient of ");
900  SCIPdebug( if( isnewx ) {printf("\nx ="); for( i = 0; i < nvars; ++i) printf(" %g", xx[i]); /*lint !e613*/ printf("\n");} )
901 
902  SCIP_CALL( SCIPexprintGrad(oracle->exprinterpreter, cons->exprtree, xx, isnewx, &nlval, g) ); /*lint !e644*/
903 
904  SCIPdebug( printf("g ="); for( i = 0; i < nvars; ++i) printf(" %g", g[i]); printf("\n"); )
905 
906  if( nlval != nlval || ABS(nlval) >= oracle->infinity ) /*lint !e777*/
907  {
908  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
909  BMSfreeBlockMemoryArray(oracle->blkmem, &g, nvars);
910  SCIPdebugMessage("gradient evaluation yield invalid function value %g\n", nlval);
911  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
912  }
913  else
914  {
915  *val += nlval;
916  for( i = 0; i < nvars; ++i )
917  if( g[i] != g[i] ) /*lint !e777*/
918  {
919  SCIPdebugMessage("gradient evaluation yield invalid gradient value %g\n", g[i]);
920  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
921  BMSfreeBlockMemoryArray(oracle->blkmem, &g, nvars);
922  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
923  }
924  else
925  {
926  grad[cons->exprvaridxs[i]] += g[i];
927  }
928  }
929 
930  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
931  BMSfreeBlockMemoryArray(oracle->blkmem, &g, nvars);
932  }
933 
934  return SCIP_OKAY;
935 }
936 
937 /** collects nonzeros entries in colnz and increases the nzcount given indices of quadratic terms */
938 static
940  SCIP_NLPIORACLE* oracle, /**< NLPI oracle */
941  int** colnz, /**< indices of nonzero entries for each column */
942  int* collen, /**< space allocated to store indices of nonzeros for each column */
943  int* colnnz, /**< number of nonzero entries for each column */
944  int* nzcount, /**< counter for total number of nonzeros; should be increased whenever some colnnz is increased */
945  int length, /**< length of quadratic part */
946  SCIP_QUADELEM* quadelems /**< quadratic elements */
947  )
948 {
949  int pos;
950 
951  SCIPdebugMessage("%p hess lag sparsity set nzflag for quad\n", (void*)oracle);
952 
953  assert(oracle != NULL);
954  assert(colnz != NULL);
955  assert(collen != NULL);
956  assert(colnnz != NULL);
957  assert(nzcount != NULL);
958  assert(quadelems != NULL);
959  assert(length >= 0);
960 
961  for( ; length > 0; --length, ++quadelems )
962  {
963  assert(quadelems->idx1 <= quadelems->idx2);
964 
965  if( colnz[quadelems->idx2] == NULL || !SCIPsortedvecFindInt(colnz[quadelems->idx2], quadelems->idx1, colnnz[quadelems->idx2], &pos) )
966  {
967  SCIP_CALL( ensureIntArraySize(oracle->blkmem, &colnz[quadelems->idx2], &collen[quadelems->idx2], colnnz[quadelems->idx2]+1) );
968  SCIPsortedvecInsertInt(colnz[quadelems->idx2], quadelems->idx1, &colnnz[quadelems->idx2], NULL);
969  ++(*nzcount);
970  }
971  }
972 
973  return SCIP_OKAY;
974 }
975 
976 /** collects indices of nonzero entries in the lower-left part of the hessian matrix of an expression
977  * adds the indices to a given set of indices, avoiding duplicates */
978 static
980  SCIP_NLPIORACLE* oracle, /**< NLPI oracle */
981  int** colnz, /**< indices of nonzero entries for each column */
982  int* collen, /**< space allocated to store indices of nonzeros for each column */
983  int* colnnz, /**< number of nonzero entries for each column */
984  int* nzcount, /**< counter for total number of nonzeros; should be increased when nzflag is set to 1 the first time */
985  int* exprvaridx, /**< indices of variables from expression tree in NLP */
986  SCIP_EXPRTREE* exprtree, /**< expression tree */
987  int dim /**< dimension of matrix */
988  )
989 {
990  SCIP_Real* x;
991  SCIP_Bool* hesnz;
992  int i;
993  int j;
994  int nvars;
995  int nn;
996  int row;
997  int col;
998  int pos;
999 
1000  assert(oracle != NULL);
1001  assert(colnz != NULL);
1002  assert(collen != NULL);
1003  assert(colnnz != NULL);
1004  assert(nzcount != NULL);
1005  assert(exprvaridx != NULL);
1006  assert(exprtree != NULL);
1007  assert(dim >= 0);
1008 
1009  SCIPdebugMessage("%p hess lag sparsity set nzflag for exprtree\n", (void*)oracle);
1010 
1011  nvars = SCIPexprtreeGetNVars(exprtree);
1012  nn = nvars * nvars;
1013 
1014  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &x, nvars) );
1015  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &hesnz, nn) );
1016 
1017  for( i = 0; i < nvars; ++i )
1018  x[i] = 2.0; /* hope that this value does not make much trouble for the evaluation routines */ /*lint !e644*/
1019 
1020  SCIP_CALL( SCIPexprintHessianSparsityDense(oracle->exprinterpreter, exprtree, x, hesnz) ); /*lint !e644*/
1021 
1022  for( i = 0; i < nvars; ++i ) /* rows */
1023  for( j = 0; j <= i; ++j ) /* cols */
1024  {
1025  if( !hesnz[i*nvars + j] )
1026  continue;
1027 
1028  row = MAX(exprvaridx[i], exprvaridx[j]);
1029  col = MIN(exprvaridx[i], exprvaridx[j]);
1030 
1031  assert(row < dim);
1032  assert(col <= row);
1033 
1034  if( colnz[row] == NULL || !SCIPsortedvecFindInt(colnz[row], col, colnnz[row], &pos) )
1035  {
1036  SCIP_CALL( ensureIntArraySize(oracle->blkmem, &colnz[row], &collen[row], colnnz[row]+1) );
1037  SCIPsortedvecInsertInt(colnz[row], col, &colnnz[row], NULL);
1038  ++(*nzcount);
1039  }
1040  }
1041 
1042  BMSfreeBlockMemoryArray(oracle->blkmem, &x, nvars);
1043  BMSfreeBlockMemoryArray(oracle->blkmem, &hesnz, nn);
1044 
1045  return SCIP_OKAY;
1046 }
1047 
1048 /** adds quadratic part into hessian structure */
1049 static
1051  SCIP_Real weight, /**< weight of quadratic part */
1052  int length, /**< number of elements in matrix of quadratic part */
1053  SCIP_QUADELEM* quadelems, /**< elements in matrix of quadratic part */
1054  int* hesoffset, /**< row offsets in sparse matrix that is to be filled */
1055  int* hescol, /**< column indices in sparse matrix that is to be filled */
1056  SCIP_Real* values /**< buffer for values of sparse matrix that is to be filled */
1057  )
1058 {
1059  int idx;
1060 
1061  SCIPdebugMessage("hess lag add quad\n");
1062 
1063  assert(length >= 0);
1064  assert(quadelems != NULL || length == 0);
1065  assert(hesoffset != NULL);
1066  assert(hescol != NULL);
1067  assert(values != NULL);
1068 
1069  for( ; length > 0; --length, ++quadelems ) /*lint !e613*/
1070  {
1071  assert(quadelems->idx1 <= quadelems->idx2); /*lint !e613*/
1072  if( !SCIPsortedvecFindInt(&hescol[hesoffset[quadelems->idx2]], quadelems->idx1, hesoffset[quadelems->idx2 + 1] - hesoffset[quadelems->idx2], &idx) ) /*lint !e613*/
1073  {
1074  SCIPerrorMessage("Could not find entry in hessian sparsity\n");
1075  return SCIP_ERROR;
1076  }
1077  values[hesoffset[quadelems->idx2] + idx] += weight * ((quadelems->idx1 == quadelems->idx2) ? 2 * quadelems->coef : quadelems->coef); /*lint !e613*/
1078  }
1079 
1080  return SCIP_OKAY;
1081 }
1082 
1083 /** adds hessian of an expression into hessian structure */
1084 static
1086  SCIP_NLPIORACLE* oracle, /**< oracle */
1087  SCIP_Real weight, /**< weight of quadratic part */
1088  const SCIP_Real* x, /**< point for which hessian should be returned */
1089  SCIP_Bool new_x, /**< whether point has been evaluated before */
1090  int* exprvaridx, /**< NLP indices for variables in expression tree */
1091  SCIP_EXPRTREE* exprtree, /**< expression tree */
1092  int* hesoffset, /**< row offsets in sparse matrix that is to be filled */
1093  int* hescol, /**< column indices in sparse matrix that is to be filled */
1094  SCIP_Real* values /**< buffer for values of sparse matrix that is to be filled */
1095  )
1096 {
1097  SCIP_Real* xx;
1098  SCIP_Real* h;
1099  SCIP_Real* hh;
1100  int i;
1101  int j;
1102  int nvars;
1103  int nn;
1104  int row;
1105  int col;
1106  int idx;
1107  SCIP_Real val;
1108 
1109  SCIPdebugMessage("%p hess lag add exprtree\n", (void*)oracle);
1110 
1111  assert(oracle != NULL);
1112  assert(x != NULL || new_x == FALSE);
1113 
1114  nvars = exprtree != NULL ? SCIPexprtreeGetNVars(exprtree) : 0;
1115  if( nvars == 0 )
1116  return SCIP_OKAY;
1117 
1118  assert(exprtree != NULL);
1119  assert(exprvaridx != NULL);
1120  assert(hesoffset != NULL);
1121  assert(hescol != NULL);
1122  assert(values != NULL);
1123 
1124  nn = nvars * nvars;
1125 
1126  xx = NULL;
1127  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &h, nn) );
1128 
1129  if( new_x )
1130  {
1131  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &xx, nvars) );
1132  for( i = 0; i < nvars; ++i )
1133  {
1134  assert(exprvaridx[i] >= 0);
1135  xx[i] = x[exprvaridx[i]]; /*lint !e613*/
1136  }
1137  }
1138 
1139  SCIP_CALL( SCIPexprintHessianDense(oracle->exprinterpreter, exprtree, xx, new_x, &val, h) ); /*lint !e644*/
1140  if( val != val ) /*lint !e777*/
1141  {
1142  SCIPdebugMessage("hessian evaluation yield invalid function value %g\n", val);
1143  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
1144  BMSfreeBlockMemoryArray(oracle->blkmem, &h, nn);
1145  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
1146  }
1147 
1148  hh = h;
1149  for( i = 0; i < nvars; ++i ) /* rows */
1150  {
1151  for( j = 0; j <= i; ++j, ++hh ) /* cols */
1152  {
1153  if( !*hh )
1154  continue;
1155 
1156  if( *hh != *hh ) /*lint !e777*/
1157  {
1158  SCIPdebugMessage("hessian evaluation yield invalid hessian value %g\n", *hh);
1159  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
1160  BMSfreeBlockMemoryArray(oracle->blkmem, &h, nn);
1161  return SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
1162  }
1163 
1164  row = MAX(exprvaridx[i], exprvaridx[j]);
1165  col = MIN(exprvaridx[i], exprvaridx[j]);
1166 
1167  if( !SCIPsortedvecFindInt(&hescol[hesoffset[row]], col, hesoffset[row+1] - hesoffset[row], &idx) )
1168  {
1169  SCIPerrorMessage("Could not find entry (%d, %d) in hessian sparsity\n", row, col);
1170  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
1171  BMSfreeBlockMemoryArray(oracle->blkmem, &h, nn);
1172  return SCIP_ERROR;
1173  }
1174 
1175  values[hesoffset[row] + idx] += weight * *hh;
1176  }
1177  hh += nvars - j;
1178  }
1179 
1180  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, nvars);
1181  BMSfreeBlockMemoryArray(oracle->blkmem, &h, nn);
1182 
1183  return SCIP_OKAY;
1184 }
1185 
1186 /** 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 */
1187 static
1189  char* buffer, /**< buffer to print to, has to be not NULL */
1190  char* name, /**< name, or NULL */
1191  int idx, /**< index of var or cons which the name corresponds to */
1192  char prefix, /**< a letter (typically 'x' or 'e') to distinguish variable and equation names, if names[idx] is not available */
1193  const char* suffix, /**< a suffer to add to the name, or NULL */
1194  SCIP_Bool longnames /**< whether prefixes for long names should be added */
1195  )
1196 {
1197  if( longnames )
1198  {
1199  if( name != NULL )
1200  sprintf(buffer, "%c%05d%.*s%s", prefix, idx, suffix ? (int)(57-strlen(suffix)) : 57, name, suffix ? suffix : "");
1201  else
1202  sprintf(buffer, "%c%05d", prefix, idx);
1203  }
1204  else
1205  {
1206  if( name != NULL )
1207  {
1208  assert(strlen(name) + (suffix ? strlen(suffix) : 0) <= 64);
1209  sprintf(buffer, "%s%s", name, suffix ? suffix : "");
1210  }
1211  else
1212  sprintf(buffer, "%c%d%s", prefix, idx, suffix ? suffix : "");
1213  }
1214 }
1215 
1216 /** prints a function */
1217 static
1219  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1220  SCIP_MESSAGEHDLR* messagehdlr, /**< message handler */
1221  FILE* file, /**< file to print to, has to be not NULL */
1222  SCIP_NLPIORACLECONS* cons, /**< constraint which function to print */
1223  SCIP_Bool longvarnames, /**< whether variable names need to be shorten to 64 characters */
1224  SCIP_Bool longequnames /**< whether equation names need to be shorten to 64 characters */
1225  )
1226 { /*lint --e{715}*/
1227  int i;
1228  int j;
1229  char namebuf[70];
1230 
1231  SCIPdebugMessage("%p print function\n", (void*)oracle);
1232 
1233  assert(oracle != NULL);
1234  assert(file != NULL);
1235  assert(cons != NULL);
1236 
1237  for( i = 0; i < cons->nlinidxs; ++i )
1238  {
1239  printName(namebuf, oracle->varnames != NULL ? oracle->varnames[cons->linidxs[i]] : NULL, cons->linidxs[i], 'x', NULL, longvarnames);
1240  SCIPmessageFPrintInfo(messagehdlr, file, "%+.20g*%s", cons->lincoefs[i], namebuf);
1241  if( i % 10 == 9 )
1242  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
1243  }
1244 
1245  j = 0;
1246  for( i = 0; i < cons->nquadelems; ++i )
1247  {
1248  printName(namebuf, oracle->varnames != NULL ? oracle->varnames[cons->quadelems[i].idx1] : NULL, cons->quadelems[i].idx1, 'x', NULL, longvarnames);
1249  SCIPmessageFPrintInfo(messagehdlr, file, "%+.20g*%s", cons->quadelems[j].coef, namebuf);
1250  printName(namebuf, oracle->varnames != NULL ? oracle->varnames[cons->quadelems[i].idx2] : NULL, cons->quadelems[i].idx2, 'x', NULL, longvarnames);
1251  SCIPmessageFPrintInfo(messagehdlr, file, "*%s", namebuf);
1252  if( i % 10 == 9 )
1253  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
1254  }
1255 
1256  if( cons->exprtree != NULL )
1257  {
1258  char** varnames;
1259  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &varnames, SCIPexprtreeGetNVars(cons->exprtree)) ); /*lint !e666*/
1260 
1261  /* setup variable names */
1262  for( i = 0; i < SCIPexprtreeGetNVars(cons->exprtree); ++i )
1263  {
1264  assert(cons->exprvaridxs[i] < 1e+20);
1265  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &varnames[i], 70) ); /*lint !e866 !e506 !e644*/
1266  printName(varnames[i], oracle->varnames != NULL ? oracle->varnames[cons->exprvaridxs[i]] : NULL, cons->exprvaridxs[i], 'x', NULL, longvarnames);
1267  }
1268 
1269  SCIPmessageFPrintInfo(messagehdlr, file, " +");
1270  SCIPexprtreePrint(cons->exprtree, messagehdlr, file, (const char**)varnames, NULL);
1271 
1272  for( i = 0; i < SCIPexprtreeGetNVars(cons->exprtree); ++i )
1273  {
1274  BMSfreeBlockMemoryArray(oracle->blkmem, &varnames[i], 70); /*lint !e866*/
1275  }
1276  BMSfreeBlockMemoryArray(oracle->blkmem, &varnames, SCIPexprtreeGetNVars(cons->exprtree));
1277  }
1278 
1279  return SCIP_OKAY;
1280 }
1281 
1282 /** returns whether an expression is contains nonsmooth operands (min, max, abs, ...) */
1283 static
1285  SCIP_EXPR* expr /**< expression */
1286  )
1287 {
1288  int i;
1289 
1290  assert(expr != NULL);
1291  assert(SCIPexprGetChildren(expr) != NULL || SCIPexprGetNChildren(expr) == 0);
1292 
1293  for( i = 0; i < SCIPexprGetNChildren(expr); ++i )
1294  {
1295  if( exprIsNonSmooth(SCIPexprGetChildren(expr)[i]) )
1296  return TRUE;
1297  }
1298 
1299  switch( SCIPexprGetOperator(expr) )
1300  {
1301  case SCIP_EXPR_MIN:
1302  case SCIP_EXPR_MAX:
1303  case SCIP_EXPR_ABS:
1304  case SCIP_EXPR_SIGN:
1305  case SCIP_EXPR_SIGNPOWER:
1306  return TRUE;
1307 
1308  default: ;
1309  } /*lint !e788*/
1310 
1311  return FALSE;
1312 }
1313 
1314 /**@} */
1315 
1316 /**@name public function */
1317 /**@{ */
1318 
1319 /** creates an NLPIORACLE data structure */
1321  BMS_BLKMEM* blkmem, /**< block memory */
1322  SCIP_NLPIORACLE** oracle /**< pointer to store NLPIORACLE data structure */
1323  )
1324 {
1325  assert(blkmem != NULL);
1326  assert(oracle != NULL);
1327 
1328  SCIPdebugMessage("%p oracle create\n", (void*)oracle);
1329 
1330  SCIP_ALLOC( BMSallocMemory(oracle) );
1331  BMSclearMemory(*oracle);
1332 
1333  (*oracle)->blkmem = blkmem;
1334  (*oracle)->infinity = SCIP_DEFAULT_INFINITY;
1335  (*oracle)->vardegreesuptodate = TRUE;
1336 
1337  SCIPdebugMessage("Oracle initializes expression interpreter %s\n", SCIPexprintGetName());
1338  SCIP_CALL( SCIPexprintCreate(blkmem, &(*oracle)->exprinterpreter) );
1339 
1340  /* create zero objective function */
1341  SCIP_CALL( createConstraint((*oracle)->blkmem, &(*oracle)->objective, 0, NULL, NULL, 0, NULL, NULL, NULL, 0.0, 0.0, NULL) );
1342 
1343  return SCIP_OKAY;
1344 }
1345 
1346 /** frees an NLPIORACLE data structure */
1348  SCIP_NLPIORACLE** oracle /**< pointer to NLPIORACLE data structure */
1349  )
1350 {
1351  assert(oracle != NULL);
1352  assert(*oracle != NULL);
1353 
1354  SCIPdebugMessage("%p oracle free\n", (void*)oracle);
1355 
1356  invalidateJacobiSparsity(*oracle);
1358 
1359  freeConstraint((*oracle)->blkmem, &(*oracle)->objective);
1360  freeConstraints(*oracle);
1361  freeVariables(*oracle);
1362 
1363  SCIP_CALL( SCIPexprintFree(&(*oracle)->exprinterpreter) );
1364 
1365  if( (*oracle)->name != NULL )
1366  {
1368  }
1369 
1370  BMSfreeMemory(oracle);
1371 
1372  return SCIP_OKAY;
1373 }
1374 
1375 /** sets the value for infinity */
1377  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1378  SCIP_Real infinity /**< value to use for infinity */
1379  )
1380 {
1381  assert(oracle != NULL);
1382  assert(infinity > 0.0);
1383 
1384  SCIPdebugMessage("%p set infinity\n", (void*)oracle);
1385 
1386  oracle->infinity = infinity;
1387 
1388  return SCIP_OKAY;
1389 }
1390 
1391 /** gets the value for infinity */
1393  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1394  )
1395 {
1396  assert(oracle != NULL);
1397 
1398  SCIPdebugMessage("%p get infinity\n", (void*)oracle);
1399 
1400  return oracle->infinity;
1401 }
1402 
1403 /** sets the problem name (used for printing) */
1405  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1406  const char* name /**< name of problem */
1407  )
1408 {
1409  assert(oracle != NULL);
1410 
1411  SCIPdebugMessage("%p set problem name\n", (void*)oracle);
1412 
1413  if( oracle->name != NULL )
1414  {
1415  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->name, strlen(oracle->name)+1);
1416  }
1417 
1418  if( name != NULL )
1419  {
1420  SCIP_ALLOC( BMSduplicateBlockMemoryArray(oracle->blkmem, &oracle->name, name, strlen(name)+1) );
1421  }
1422 
1423  return SCIP_OKAY;
1424 }
1425 
1426 /** gets the problem name, or NULL if none set */
1428  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
1429  )
1430 {
1431  assert(oracle != NULL);
1432 
1433  SCIPdebugMessage("%p get problem name\n", (void*)oracle);
1434 
1435  return oracle->name;
1436 }
1437 
1438 /** adds variables */
1440  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1441  int nvars, /**< number of variables to add */
1442  const SCIP_Real* lbs, /**< array with lower bounds of new variables, or NULL if all -infinity */
1443  const SCIP_Real* ubs, /**< array with upper bounds of new variables, or NULL if all +infinity */
1444  const char** varnames /**< array with names of new variables, or NULL if no names should be stored */
1445  )
1446 {
1447  int i;
1448 
1449  assert(oracle != NULL);
1450 
1451  SCIPdebugMessage("%p add vars\n", (void*)oracle);
1452 
1453  if( nvars == 0 )
1454  return SCIP_OKAY;
1455 
1456  assert(nvars > 0);
1457 
1458  SCIP_CALL( ensureVarsSize(oracle, oracle->nvars + nvars) );
1459 
1460  if( lbs != NULL )
1461  {
1462  BMScopyMemoryArray(&oracle->varlbs[oracle->nvars], lbs, nvars); /*lint !e866*/
1463  }
1464  else
1465  for( i = 0; i < nvars; ++i )
1466  oracle->varlbs[oracle->nvars+i] = -oracle->infinity;
1467 
1468  if( ubs != NULL )
1469  {
1470  BMScopyMemoryArray(&oracle->varubs[oracle->nvars], ubs, nvars); /*lint !e866*/
1471 
1472  /* ensure variable bounds are consistent */
1473  for( i = oracle->nvars; i < oracle->nvars + nvars; ++i )
1474  {
1475  if( oracle->varlbs[i] > oracle->varubs[i] )
1476  {
1477  assert(EPSEQ(oracle->varlbs[i], oracle->varubs[i], SCIP_DEFAULT_EPSILON));
1478  oracle->varlbs[i] = oracle->varubs[i];
1479  }
1480  }
1481  }
1482  else
1483  for( i = 0; i < nvars; ++i )
1484  oracle->varubs[oracle->nvars+i] = oracle->infinity;
1485 
1486  if( varnames != NULL )
1487  {
1488  if( oracle->varnames == NULL )
1489  {
1490  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &oracle->varnames, oracle->varssize) );
1491  BMSclearMemoryArray(oracle->varnames, oracle->nvars);
1492  }
1493 
1494  for( i = 0; i < nvars; ++i )
1495  {
1496  if( varnames[i] != NULL )
1497  {
1498  SCIP_ALLOC( BMSduplicateBlockMemoryArray(oracle->blkmem, &oracle->varnames[oracle->nvars+i], varnames[i], strlen(varnames[i])+1) ); /*lint !e866*/
1499  }
1500  else
1501  oracle->varnames[oracle->nvars+i] = NULL;
1502  }
1503  }
1504  else if( oracle->varnames != NULL )
1505  {
1506  BMSclearMemoryArray(&oracle->varnames[oracle->nvars], nvars); /*lint !e866*/
1507  }
1508 
1509  BMSclearMemoryArray(&oracle->vardegrees[oracle->nvars], nvars); /*lint !e866*/
1510 
1511  /* @TODO update sparsity pattern by extending heslagoffsets */
1513 
1514  oracle->nvars += nvars;
1515 
1516  return SCIP_OKAY;
1517 }
1518 
1519 /** adds constraints
1520  *
1521  * linear coefficients: row(=constraint) oriented matrix;
1522  * quadratic coefficients: row oriented matrix for each constraint
1523  */
1525  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1526  int nconss, /**< number of constraints to add */
1527  const SCIP_Real* lhss, /**< array with left-hand sides of constraints, or NULL if all -infinity */
1528  const SCIP_Real* rhss, /**< array with right-hand sides of constraints, or NULL if all +infinity */
1529  const int* nlininds, /**< number of linear coefficients for each constraint, may be NULL in case of no linear part */
1530  int* const* lininds, /**< indices of variables for linear coefficients for each constraint, may be NULL in case of no linear part */
1531  SCIP_Real* const* linvals, /**< values of linear coefficient for each constraint, may be NULL in case of no linear part */
1532  const int* nquadelems, /**< number of elements in matrix of quadratic part for each constraint,
1533  * may be NULL in case of no quadratic part in any constraint */
1534  SCIP_QUADELEM* const* quadelems, /**< quadratic elements specifying quadratic part for each constraint, entry of array may be NULL in case of no quadratic part,
1535  * may be NULL in case of no quadratic part in any constraint */
1536  int* const* exprvaridxs, /**< NULL if no nonquadratic parts, otherwise epxrvaridxs[.] maps variable indices in expression tree to indices in nlp */
1537  SCIP_EXPRTREE* const* exprtrees, /**< NULL if no nonquadratic parts, otherwise exprtrees[.] gives nonquadratic part,
1538  * or NULL if no nonquadratic part in this constraint */
1539  const char** consnames /**< names of new constraints, or NULL if no names should be stored */
1540  )
1541 { /*lint --e{715}*/
1542  SCIP_NLPIORACLECONS* cons;
1543  SCIP_Bool addednlcon; /* whether a nonlinear constraint was added */
1544  int c;
1545 
1546  assert(oracle != NULL);
1547 
1548  SCIPdebugMessage("%p add constraints\n", (void*)oracle);
1549 
1550  if( nconss == 0 )
1551  return SCIP_OKAY;
1552 
1553  assert(nconss > 0);
1554 
1555  addednlcon = FALSE;
1556 
1557  invalidateJacobiSparsity(oracle); /* @TODO we could also update (extend) the sparsity pattern */
1558 
1559  SCIP_CALL( ensureConssSize(oracle, oracle->nconss + nconss) );
1560  for( c = 0; c < nconss; ++c )
1561  {
1562  SCIP_CALL( createConstraint(oracle->blkmem, &cons,
1563  nlininds != NULL ? nlininds[c] : 0,
1564  lininds != NULL ? lininds[c] : NULL,
1565  linvals != NULL ? linvals[c] : NULL,
1566  nquadelems != NULL ? nquadelems[c] : 0,
1567  quadelems != NULL ? quadelems[c] : NULL,
1568  exprvaridxs != NULL ? exprvaridxs[c] : NULL,
1569  exprtrees != NULL ? exprtrees[c] : NULL,
1570  lhss != NULL ? lhss[c] : -oracle->infinity,
1571  rhss != NULL ? rhss[c] : oracle->infinity,
1572  consnames != NULL ? consnames[c] : NULL
1573  ) );
1574 
1575  if( cons->nquadelems > 0 )
1576  addednlcon = TRUE;
1577 
1578  if( cons->exprtree != NULL )
1579  {
1580  addednlcon = TRUE;
1582  }
1583 
1584  /* keep variable degrees updated */
1585  if( oracle->vardegreesuptodate )
1586  updateVariableDegreesCons(oracle, cons);
1587 
1588  oracle->conss[oracle->nconss+c] = cons;
1589  }
1590  oracle->nconss += nconss;
1591 
1592  if( addednlcon == TRUE )
1594 
1595  return SCIP_OKAY;
1596 }
1597 
1598 /** sets or overwrites objective, a minimization problem is expected
1599  *
1600  * May change sparsity pattern.
1601  */
1603  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1604  const SCIP_Real constant, /**< constant part of objective */
1605  int nlin, /**< number of linear variable coefficients */
1606  const int* lininds, /**< indices of linear variables, or NULL if no linear part */
1607  const SCIP_Real* linvals, /**< coefficients of linear variables, or NULL if no linear part */
1608  int nquadelems, /**< number of entries in matrix of quadratic part */
1609  const SCIP_QUADELEM* quadelems, /**< entries in matrix of quadratic part, may be NULL in case of no quadratic part */
1610  const int* exprvaridxs, /**< maps variable indices in expression tree to indices in nlp, or NULL if no nonquadratic part */
1611  const SCIP_EXPRTREE* exprtree /**< expression tree of nonquadratic part, or NULL if no nonquadratic part */
1612  )
1613 { /*lint --e{715}*/
1614  assert(oracle != NULL);
1615  assert(REALABS(constant) < oracle->infinity);
1616 
1617  SCIPdebugMessage("%p set objective\n", (void*)oracle);
1618 
1619  if( nquadelems > 0 || oracle->objective->quadsize > 0 || exprtree != NULL || oracle->objective->exprtree != NULL )
1621 
1622  /* clear previous objective */
1623  freeConstraint(oracle->blkmem, &oracle->objective);
1624 
1625  SCIP_CALL( createConstraint(oracle->blkmem, &oracle->objective,
1626  nlin, lininds, linvals, nquadelems, quadelems, exprvaridxs, exprtree, constant, constant, NULL) );
1627 
1628  if( oracle->objective->exprtree != NULL )
1629  {
1631  }
1632 
1633  oracle->vardegreesuptodate = FALSE;
1634 
1635  return SCIP_OKAY;
1636 }
1637 
1638 /** change variable bounds */
1640  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1641  int nvars, /**< number of variables to change bounds */
1642  const int* indices, /**< indices of variables to change bounds */
1643  const SCIP_Real* lbs, /**< new lower bounds, or NULL if all should be -infty */
1644  const SCIP_Real* ubs /**< new upper bounds, or NULL if all should be +infty */
1645  )
1646 {
1647  int i;
1648 
1649  assert(oracle != NULL);
1650  assert(indices != NULL || nvars == 0);
1651 
1652  SCIPdebugMessage("%p chg var bounds\n", (void*)oracle);
1653 
1654  for( i = 0; i < nvars; ++i )
1655  {
1656  assert(indices != NULL);
1657  assert(indices[i] >= 0);
1658  assert(indices[i] < oracle->nvars);
1659 
1660  oracle->varlbs[indices[i]] = (lbs != NULL ? lbs[i] : -oracle->infinity);
1661  oracle->varubs[indices[i]] = (ubs != NULL ? ubs[i] : oracle->infinity);
1662 
1663  if( oracle->varlbs[indices[i]] > oracle->varubs[indices[i]] )
1664  {
1665  /* inconsistent bounds; let's assume it's due to rounding and make them equal */
1666  assert(EPSEQ(oracle->varlbs[indices[i]], oracle->varubs[indices[i]], SCIP_DEFAULT_EPSILON));
1667  oracle->varlbs[indices[i]] = oracle->varubs[indices[i]];
1668  }
1669  }
1670 
1671  return SCIP_OKAY;
1672 }
1673 
1674 /** change constraint bounds */
1676  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1677  int nconss, /**< number of constraints to change bounds */
1678  const int* indices, /**< indices of constraints to change bounds */
1679  const SCIP_Real* lhss, /**< new left-hand sides, or NULL if all should be -infty */
1680  const SCIP_Real* rhss /**< new right-hand sides, or NULL if all should be +infty */
1681  )
1682 {
1683  int i;
1684 
1685  assert(oracle != NULL);
1686  assert(indices != NULL || nconss == 0);
1687 
1688  SCIPdebugMessage("%p chg cons sides\n", (void*)oracle);
1689 
1690  for( i = 0; i < nconss; ++i )
1691  {
1692  assert(indices != NULL);
1693  assert(indices[i] >= 0);
1694  assert(indices[i] < oracle->nconss);
1695 
1696  oracle->conss[indices[i]]->lhs = (lhss != NULL ? lhss[i] : -oracle->infinity);
1697  oracle->conss[indices[i]]->rhs = (rhss != NULL ? rhss[i] : oracle->infinity);
1698  if( oracle->conss[indices[i]]->lhs > oracle->conss[indices[i]]->rhs )
1699  {
1700  assert(EPSEQ(oracle->conss[indices[i]]->lhs, oracle->conss[indices[i]]->rhs, SCIP_DEFAULT_EPSILON));
1701  oracle->conss[indices[i]]->lhs = oracle->conss[indices[i]]->rhs;
1702  }
1703  }
1704 
1705  return SCIP_OKAY;
1706 }
1707 
1708 /** deletes a set of variables */
1710  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1711  int* delstats /**< deletion status of vars in input (1 if var should be deleted, 0 if not);
1712  * new position of var in output (-1 if var was deleted) */
1713  )
1714 { /*lint --e{715}*/
1715  int c;
1716  int lastgood; /* index of the last variable that should be kept */
1717  SCIP_NLPIORACLECONS* cons;
1718 
1719  assert(oracle != NULL);
1720 
1721  SCIPdebugMessage("%p del var set\n", (void*)oracle);
1722 
1723  invalidateJacobiSparsity(oracle);
1725 
1726  lastgood = oracle->nvars - 1;
1727  while( lastgood >= 0 && delstats[lastgood] == 1 )
1728  --lastgood;
1729  if( lastgood < 0 )
1730  {
1731  /* all variables should be deleted */
1732  assert(oracle->nconss == 0); /* we could relax this by checking that all constraints are constant */
1733  assert(oracle->objective->exprtree == NULL || SCIPexprtreeGetNVars(oracle->objective->exprtree) == 0);
1734  oracle->objective->nquadelems = 0;
1735  oracle->objective->nlinidxs = 0;
1736  for( c = 0; c < oracle->nvars; ++c )
1737  delstats[c] = -1;
1738  freeVariables(oracle);
1739  return SCIP_OKAY;
1740  }
1741 
1742  /* delete variables at the end */
1743  for( c = oracle->nvars - 1; c > lastgood; --c )
1744  {
1745  if( oracle->varnames && oracle->varnames[c] != NULL )
1746  {
1747  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->varnames[c], strlen(oracle->varnames[c])+1); /*lint !e866*/
1748  }
1749  delstats[c] = -1;
1750  }
1751 
1752  /* go through variables from the beginning on
1753  * if variable should be deleted, free it and move lastgood variable to this position
1754  * then update lastgood */
1755  for( c = 0; c <= lastgood; ++c )
1756  {
1757  if( delstats[c] == 0 )
1758  { /* variable should not be deleted and is kept on position c */
1759  delstats[c] = c;
1760  continue;
1761  }
1762  assert(delstats[c] == 1); /* variable should be deleted */
1763 
1764  if( oracle->varnames && oracle->varnames[c] != NULL )
1765  {
1766  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->varnames[c], strlen(oracle->varnames[c])+1); /*lint !e866*/
1767  }
1768  delstats[c] = -1;
1769 
1770  /* move constraint at position lastgood to position c */
1771  SCIP_CALL( moveVariable(oracle, lastgood, c) );
1772  delstats[lastgood] = c; /* mark that lastgood variable is now at position c */
1773 
1774  /* move lastgood forward, delete variables on the way */
1775  --lastgood;
1776  while( lastgood > c && delstats[lastgood] == 1)
1777  {
1778  if( oracle->varnames && oracle->varnames[lastgood] != NULL )
1779  {
1780  BMSfreeBlockMemoryArray(oracle->blkmem, &oracle->varnames[lastgood], strlen(oracle->varnames[lastgood])+1); /*lint !e866*/
1781  }
1782  delstats[lastgood] = -1;
1783  --lastgood;
1784  }
1785  }
1786  assert(c == lastgood);
1787 
1788  for( c = -1; c < oracle->nconss; ++c )
1789  {
1790  cons = c < 0 ? oracle->objective : oracle->conss[c];
1791  assert(cons != NULL);
1792 
1793  /* update indices in linear part, sort indices, and then clear elements that are marked as deleted */
1794  mapIndices(delstats, cons->nlinidxs, cons->linidxs);
1795  SCIPsortIntReal(cons->linidxs, cons->lincoefs, cons->nlinidxs);
1796  clearDeletedLinearElements(oracle->blkmem, &cons->linidxs, &cons->lincoefs, &cons->nlinidxs);
1797 
1798  /* update indices in quadratic part, sort elements, and then clear elements that are marked as deleted */
1799  mapIndicesQuad(delstats, cons->quadsize, cons->quadelems);
1800  SCIPquadelemSort(cons->quadelems, cons->quadsize);
1801  clearDeletedQuadElements(oracle->blkmem, &cons->quadelems, &cons->quadsize);
1802 
1803  if( cons->exprtree != NULL )
1804  {
1805  mapIndices(delstats, SCIPexprtreeGetNVars(cons->exprtree), cons->exprvaridxs);
1806  /* assert that all variables from this expression have been deleted */
1807  assert(SCIPexprtreeGetNVars(cons->exprtree) == 0 || cons->exprvaridxs[SCIPexprtreeGetNVars(cons->exprtree)-1] == -1);
1809  SCIP_CALL( SCIPexprtreeFree(&cons->exprtree) );
1810  }
1811  }
1812 
1813  oracle->nvars = lastgood+1;
1814 
1815  return SCIP_OKAY;
1816 }
1817 
1818 /** deletes a set of constraints */
1820  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1821  int* delstats /**< array with deletion status of rows in input (1 if row should be deleted, 0 if not);
1822  * new position of row in output (-1 if row was deleted) */
1823  )
1824 { /*lint --e{715}*/
1825  int c;
1826  int lastgood; /* index of the last constraint that should be kept */
1827 
1828  assert(oracle != NULL);
1829 
1830  SCIPdebugMessage("%p del cons set\n", (void*)oracle);
1831 
1832  invalidateJacobiSparsity(oracle);
1834  oracle->vardegreesuptodate = FALSE;
1835 
1836  lastgood = oracle->nconss - 1;
1837  while( lastgood >= 0 && delstats[lastgood] == 1)
1838  --lastgood;
1839  if( lastgood < 0 )
1840  {
1841  /* all constraints should be deleted */
1842  for( c = 0; c < oracle->nconss; ++c )
1843  delstats[c] = -1;
1844  freeConstraints(oracle);
1845  return SCIP_OKAY;
1846  }
1847 
1848  /* delete constraints at the end */
1849  for( c = oracle->nconss - 1; c > lastgood; --c )
1850  {
1851  freeConstraint(oracle->blkmem, &oracle->conss[c]);
1852  assert(oracle->conss[c] == NULL);
1853  delstats[c] = -1;
1854  }
1855 
1856  /* go through constraint from the beginning on
1857  * if constraint should be deleted, free it and move lastgood constraint to this position
1858  * then update lastgood */
1859  for( c = 0; c <= lastgood; ++c )
1860  {
1861  if( delstats[c] == 0 )
1862  {
1863  /* constraint should not be deleted and is kept on position c */
1864  delstats[c] = c;
1865  continue;
1866  }
1867  assert(delstats[c] == 1); /* constraint should be deleted */
1868 
1869  freeConstraint(oracle->blkmem, &oracle->conss[c]);
1870  assert(oracle->conss[c] == NULL);
1871  delstats[c] = -1;
1872 
1873  /* move constraint at position lastgood to position c */
1874  oracle->conss[c] = oracle->conss[lastgood];
1875  assert(oracle->conss[c] != NULL);
1876  delstats[lastgood] = c; /* mark that lastgood constraint is now at position c */
1877  oracle->conss[lastgood] = NULL;
1878 
1879  /* move lastgood forward, delete constraints on the way */
1880  --lastgood;
1881  while( lastgood > c && delstats[lastgood] == 1)
1882  {
1883  freeConstraint(oracle->blkmem, &oracle->conss[lastgood]);
1884  assert(oracle->conss[lastgood] == NULL);
1885  delstats[lastgood] = -1;
1886  --lastgood;
1887  }
1888  }
1889  assert(c == lastgood);
1890 
1891  oracle->nconss = lastgood+1;
1892 
1893  return SCIP_OKAY;
1894 }
1895 
1896 /** changes linear coefficients in one constraint or objective */
1898  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1899  int considx, /**< index of constraint where linear coefficients should be changed, or -1 for objective */
1900  int nentries, /**< number of coefficients to change */
1901  const int* varidxs, /**< array with indices of variables which coefficients should be changed */
1902  const SCIP_Real* newcoefs /**< array with new coefficients of variables */
1903  )
1904 { /*lint --e{715}*/
1905  SCIP_NLPIORACLECONS* cons;
1906  SCIP_Bool needsort;
1907  int i;
1908 
1909  SCIPdebugMessage("%p chg linear coefs\n", (void*)oracle);
1910 
1911  assert(oracle != NULL);
1912  assert(varidxs != NULL || nentries == 0);
1913  assert(newcoefs != NULL || nentries == 0);
1914  assert(considx >= -1);
1915  assert(considx < oracle->nconss);
1916 
1917  if( nentries == 0 )
1918  return SCIP_OKAY;
1919 
1920  SCIPdebugMessage("change %d linear coefficients in cons %d\n", nentries, considx);
1921 
1922  needsort = FALSE;
1923 
1924  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
1925 
1926  if( cons->linsize == 0 )
1927  {
1928  /* first time we have linear coefficients in this constraint (or objective) */
1929  assert(cons->linidxs == NULL);
1930  assert(cons->lincoefs == NULL);
1931 
1932  SCIP_ALLOC( BMSduplicateBlockMemoryArray(oracle->blkmem, &cons->linidxs, varidxs, nentries) );
1933  SCIP_ALLOC( BMSduplicateBlockMemoryArray(oracle->blkmem, &cons->lincoefs, newcoefs, nentries) );
1934  cons->linsize = nentries;
1935  cons->nlinidxs = nentries;
1936 
1937  needsort = TRUE;
1938  }
1939  else
1940  {
1941  int pos;
1942 
1943  for( i = 0; i < nentries; ++i )
1944  {
1945  assert(varidxs[i] >= 0); /*lint !e613*/
1946  assert(varidxs[i] < oracle->nvars); /*lint !e613*/
1947 
1948  if( SCIPsortedvecFindInt(cons->linidxs, varidxs[i], cons->nlinidxs, &pos) ) /*lint !e613*/
1949  {
1950  SCIPdebugMessage("replace coefficient of var %d at pos %d by %g\n", varidxs[i], pos, newcoefs[i]); /*lint !e613*/
1951 
1952  cons->lincoefs[pos] = newcoefs[i]; /*lint !e613*/
1953 
1954  /* remember that we need to sort/merge/squeeze array if coefficient became zero here */
1955  needsort |= (newcoefs[i] == 0.0); /*lint !e613 !e514*/
1956  }
1957  else if( newcoefs[i] != 0.0 ) /*lint !e613*/
1958  {
1959  /* append new entry */
1960  SCIPdebugMessage("add coefficient of var %d at pos %d, value %g\n", varidxs[i], cons->nlinidxs, newcoefs[i]); /*lint !e613*/
1961 
1962  SCIP_CALL( ensureConsLinSize(oracle->blkmem, cons, cons->nlinidxs + (nentries-i)) );
1963  cons->linidxs[cons->nlinidxs] = varidxs[i]; /*lint !e613*/
1964  cons->lincoefs[cons->nlinidxs] = newcoefs[i]; /*lint !e613*/
1965  ++cons->nlinidxs;
1966 
1967  needsort = TRUE;
1968  }
1969  }
1970  }
1971 
1972  if( needsort )
1973  {
1974  int oldlen;
1975 
1976  invalidateJacobiSparsity(oracle);
1977 
1978  oldlen = cons->nlinidxs;
1979  sortLinearCoefficients(&cons->nlinidxs, cons->linidxs, cons->lincoefs);
1980 
1981  /* if sorting removed an entry, then the var degrees are not uptodate anymore */
1982  oracle->vardegreesuptodate &= (cons->nlinidxs == oldlen); /*lint !e514*/
1983 
1984  /* increase variable degrees of variables to 1 */
1985  if( oracle->vardegreesuptodate )
1986  for( i = 0; i < cons->nlinidxs; ++i )
1987  oracle->vardegrees[varidxs[i]] = MAX(1, oracle->vardegrees[varidxs[i]]); /*lint !e613*/
1988  }
1989 
1990  return SCIP_OKAY;
1991 }
1992 
1993 /** changes (or adds) coefficients in the quadratic part of one constraint or objective */
1995  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
1996  int considx, /**< index of constraint where quadratic coefficients should be changed, or -1 for objective */
1997  int nquadelems, /**< number of entries in quadratic constraint to change */
1998  const SCIP_QUADELEM* quadelems /**< new elements in quadratic matrix (replacing already existing ones or adding new ones) */
1999  )
2000 { /*lint --e{715}*/
2001  SCIP_NLPIORACLECONS* cons;
2002  SCIP_Bool needsort;
2003  int i;
2004 
2005  SCIPdebugMessage("%p chg quad coefs\n", (void*)oracle);
2006 
2007  assert(oracle != NULL);
2008  assert(quadelems != NULL || nquadelems == 0);
2009  assert(considx >= -1);
2010  assert(considx < oracle->nconss);
2011 
2012  if( nquadelems == 0 )
2013  return SCIP_OKAY;
2014 
2015  needsort = FALSE;
2016 
2017  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
2018 
2019  if( cons->quadsize == 0 )
2020  {
2021  /* first time we have quadratic coefficients in this constraint (or objective) */
2022  assert(cons->quadelems == NULL);
2023 
2024  SCIP_ALLOC( BMSduplicateBlockMemoryArray(oracle->blkmem, &cons->quadelems, quadelems, nquadelems) );
2025  cons->quadsize = nquadelems;
2026  cons->nquadelems = nquadelems;
2027 
2028  needsort = TRUE;
2029  }
2030  else
2031  {
2032  int pos;
2033 
2034  for( i = 0; i < nquadelems; ++i )
2035  {
2036  assert(quadelems[i].idx1 >= 0); /*lint !e613*/
2037  assert(quadelems[i].idx2 >= 0); /*lint !e613*/
2038  assert(quadelems[i].idx1 < oracle->nvars); /*lint !e613*/
2039  assert(quadelems[i].idx2 < oracle->nvars); /*lint !e613*/
2040 
2041  /* if we already have an entry for quadelems[i], then just replace the coefficient, otherwise append new entry */
2042  if( SCIPquadelemSortedFind(cons->quadelems, quadelems[i].idx1, quadelems[i].idx2, cons->nquadelems, &pos) ) /*lint !e613*/
2043  {
2044  SCIPdebugMessage("replace coefficient of var%d*var%d at pos %d by %g\n", quadelems[i].idx1, quadelems[i].idx2, pos, quadelems[i].coef); /*lint !e613*/
2045 
2046  cons->quadelems[pos].coef = quadelems[i].coef; /*lint !e613*/
2047 
2048  /* remember that we need to sort/merge/squeeze array if coefficient became zero here */
2049  needsort |= (quadelems[i].coef == 0.0); /*lint !e613 !e514*/
2050  }
2051  else
2052  {
2053  /* append new entry */
2054  SCIPdebugMessage("add coefficient of var%d*var%d at pos %d, value %g\n", quadelems[i].idx1, quadelems[i].idx2, cons->nquadelems, quadelems[i].coef); /*lint !e613*/
2055 
2056  SCIP_CALL( ensureConsQuadSize(oracle->blkmem, cons, cons->nquadelems + (nquadelems-i)) );
2057  cons->quadelems[cons->nquadelems] = quadelems[i]; /*lint !e613*/
2058  ++cons->nquadelems;
2059 
2060  needsort = TRUE;
2061  }
2062  }
2063  }
2064 
2065  if( needsort )
2066  {
2067  int oldsize;
2068 
2069  invalidateJacobiSparsity(oracle);
2071 
2072  oldsize = cons->nquadelems;
2073  SCIPquadelemSort(cons->quadelems, cons->nquadelems);
2074  SCIPquadelemSqueeze(cons->quadelems, cons->nquadelems, &cons->nquadelems);
2075 
2076  /* if sorting removed an entry, then the var degrees are not uptodate anymore */
2077  oracle->vardegreesuptodate &= (cons->nquadelems == oldsize); /*lint !e514*/
2078 
2079  /* increase variable degrees of variables to 2 */
2080  if( oracle->vardegreesuptodate )
2081  for( i = 0; i < cons->nquadelems; ++i )
2082  {
2083  oracle->vardegrees[cons->quadelems[i].idx1] = MAX(2, oracle->vardegrees[cons->quadelems[i].idx1]);
2084  oracle->vardegrees[cons->quadelems[i].idx2] = MAX(2, oracle->vardegrees[cons->quadelems[i].idx2]);
2085  }
2086  }
2087 
2088  return SCIP_OKAY;
2089 }
2090 
2091 /** replaces expression tree of one constraint or objective */
2093  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2094  int considx, /**< index of constraint where expression tree should be changed, or -1 for objective */
2095  const int* exprvaridxs, /**< problem indices of variables in expression tree */
2096  const SCIP_EXPRTREE* exprtree /**< new expression tree, or NULL */
2097  )
2098 {
2099  SCIP_NLPIORACLECONS* cons;
2100  int j;
2101 
2102  SCIPdebugMessage("%p chg exprtree\n", (void*)oracle);
2103 
2104  assert(oracle != NULL);
2105  assert(considx >= -1);
2106  assert(considx < oracle->nconss);
2107  assert((exprvaridxs != NULL) == (exprtree != NULL));
2108 
2110  invalidateJacobiSparsity(oracle);
2111 
2112  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
2113 
2114  /* free previous expression tree */
2115  if( cons->exprtree != NULL )
2116  {
2119  oracle->vardegreesuptodate = FALSE;
2120  }
2121 
2122  /* if user did not want to set new tree, then we are done */
2123  if( exprtree == NULL )
2124  return SCIP_OKAY;
2125 
2126  assert(oracle->exprinterpreter != NULL);
2127 
2128  /* install new expression tree */
2129  SCIP_CALL( SCIPexprtreeCopy(oracle->blkmem, &cons->exprtree, (SCIP_EXPRTREE*)exprtree) );
2132 
2133  /* increase variable degree to keep them up to date
2134  * could get more accurate degree via getMaxDegree function in exprtree, but no solver would use this information so far
2135  */
2136  if( oracle->vardegreesuptodate )
2137  for( j = 0; j < SCIPexprtreeGetNVars(cons->exprtree); ++j )
2138  {
2139  assert(cons->exprvaridxs[j] >= 0);
2140  assert(cons->exprvaridxs[j] < oracle->nvars);
2141  oracle->vardegrees[cons->exprvaridxs[j]] = INT_MAX;
2142  }
2143 
2144  return SCIP_OKAY;
2145 }
2146 
2147 /** changes one parameter of expression tree of one constraint or objective
2148  */
2150  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2151  int considx, /**< index of constraint where parameter should be changed in expression tree, or -1 for objective */
2152  int paramidx, /**< index of parameter */
2153  SCIP_Real paramval /**< new value of parameter */
2154  )
2155 {
2156  SCIPdebugMessage("%p chg expr param\n", (void*)oracle);
2157 
2158  assert(oracle != NULL);
2159  assert(considx >= -1);
2160  assert(considx < oracle->nconss);
2161  assert(paramidx >= 0);
2162  assert(considx >= 0 || oracle->objective->exprtree != NULL);
2163  assert(considx >= 0 || paramidx < SCIPexprtreeGetNParams(oracle->objective->exprtree));
2164  assert(considx == -1 || oracle->conss[considx]->exprtree != NULL);
2165  assert(considx == -1 || paramidx < SCIPexprtreeGetNParams(oracle->conss[considx]->exprtree));
2166 
2167  SCIPexprtreeSetParamVal(considx >= 0 ? oracle->conss[considx]->exprtree : oracle->objective->exprtree, paramidx, paramval);
2168 
2169  return SCIP_OKAY;
2170 }
2171 
2172 /** changes the constant value in the objective function
2173  */
2175  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2176  SCIP_Real objconstant /**< new value for objective constant */
2177  )
2178 {
2179  assert(oracle != NULL);
2180 
2181  SCIPdebugMessage("%p chg obj constant\n", (void*)oracle);
2182 
2183  oracle->objective->lhs = objconstant;
2184  oracle->objective->rhs = objconstant;
2185 
2186  return SCIP_OKAY;
2187 }
2188 
2189 /** gives the current number of variables */
2191  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2192  )
2193 {
2194  assert(oracle != NULL);
2195 
2196  return oracle->nvars;
2197 }
2198 
2199 /** gives the current number of constraints */
2201  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2202  )
2203 {
2204  assert(oracle != NULL);
2205 
2206  return oracle->nconss;
2207 }
2208 
2209 /** gives the variables lower bounds */
2211  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2212  )
2213 {
2214  assert(oracle != NULL);
2215 
2216  return oracle->varlbs;
2217 }
2218 
2219 /** gives the variables upper bounds */
2221  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2222  )
2223 {
2224  assert(oracle != NULL);
2225 
2226  return oracle->varubs;
2227 }
2228 
2229 /** gives the variables names, or NULL if not set */
2231  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2232  )
2233 {
2234  assert(oracle != NULL);
2235 
2236  return oracle->varnames;
2237 }
2238 
2239 /** Gives maximum degree of a variable w.r.t. objective and all constraints.
2240  * The degree of a variable is the degree of the summand where it appears in, and is infinity for nonpolynomial terms.
2241  */
2243  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2244  int varidx
2245  )
2246 {
2247  assert(oracle != NULL);
2248  assert(varidx >= 0);
2249  assert(varidx < oracle->nvars);
2250 
2251  updateVariableDegrees(oracle);
2252 
2253  return oracle->vardegrees[varidx];
2254 }
2255 
2256 /** Gives maximum degree of all variables w.r.t. objective and all constraints.
2257  * The degree of a variable is the degree of the summand where it appears in, and is infinity for nonpolynomial terms.
2258  */
2260  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2261  )
2262 {
2263  assert(oracle != NULL);
2264 
2265  updateVariableDegrees(oracle);
2266 
2267  return oracle->vardegrees;
2268 }
2269 
2270 /** gives left-hand side of a constraint */
2272  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2273  int considx /**< constraint index */
2274  )
2275 {
2276  assert(oracle != NULL);
2277  assert(considx >= 0);
2278  assert(considx < oracle->nconss);
2279 
2280  return oracle->conss[considx]->lhs;
2281 }
2282 
2283 /** gives right-hand side of a constraint */
2285  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2286  int considx /**< constraint index */
2287  )
2288 {
2289  assert(oracle != NULL);
2290  assert(considx >= 0);
2291  assert(considx < oracle->nconss);
2292 
2293  return oracle->conss[considx]->rhs;
2294 }
2295 
2296 /** gives name of a constraint, may be NULL */
2298  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2299  int considx /**< constraint index */
2300  )
2301 {
2302  assert(oracle != NULL);
2303  assert(considx >= 0);
2304  assert(considx < oracle->nconss);
2305 
2306  return oracle->conss[considx]->name;
2307 }
2308 
2309 /** gives maximum degree of a constraint or objective
2310  * The degree is the maximal degree of all summands,, and is infinity for nonpolynomial terms.
2311  */
2313  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2314  int considx /**< index of constraint for which the degree is requested, or -1 for objective */
2315  )
2316 {
2317  SCIP_NLPIORACLECONS* cons;
2318 
2319  assert(oracle != NULL);
2320  assert(considx >= -1);
2321  assert(considx < oracle->nconss);
2322 
2323  cons = considx < 0 ? oracle->objective : oracle->conss[considx];
2324 
2325  /* could do something more clever like exprtreeGetMaxDegree, but no solver uses this so far */
2326  if( cons->exprtree != NULL )
2327  return INT_MAX;
2328 
2329  if( cons->nquadelems > 0 )
2330  return 2;
2331 
2332  if( cons->nlinidxs > 0 )
2333  return 1;
2334 
2335  return 0;
2336 }
2337 
2338 /** Gives maximum degree over all constraints and the objective (or over all variables, resp.).
2339  * Thus, if this function returns 0, then the objective and all constraints are constant.
2340  * If it returns 1, then the problem in linear.
2341  * If it returns 2, then its a QP, QCP, or QCQP.
2342  * And if it returns > 2, then it is an NLP.
2343  */
2345  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2346  )
2347 {
2348  int i;
2349  int maxdegree;
2350 
2351  assert(oracle != NULL);
2352 
2353  SCIPdebugMessage("%p get max degree\n", (void*)oracle);
2354 
2355  updateVariableDegrees(oracle);
2356 
2357  maxdegree = 0;
2358  for( i = 0; i < oracle->nvars; ++i )
2359  if( oracle->vardegrees[i] > maxdegree )
2360  {
2361  maxdegree = oracle->vardegrees[i];
2362  if( maxdegree == INT_MAX )
2363  break;
2364  }
2365 
2366  return maxdegree;
2367 }
2368 
2369 /** Gives the evaluation capabilities that are shared among all expression trees in the problem. */
2371  SCIP_NLPIORACLE* oracle /**< pointer to NLPIORACLE data structure */
2372  )
2373 {
2374  int c;
2375  SCIP_EXPRINTCAPABILITY evalcapability;
2376 
2377  assert(oracle != NULL);
2378 
2379  if( oracle->objective->exprtree != NULL )
2380  evalcapability = SCIPexprintGetExprtreeCapability(oracle->exprinterpreter, oracle->objective->exprtree);
2381  else
2382  evalcapability = SCIP_EXPRINTCAPABILITY_ALL;
2383 
2384  for( c = 0; c < oracle->nconss; ++c )
2385  if( oracle->conss[c]->exprtree != NULL )
2386  evalcapability &= SCIPexprintGetExprtreeCapability(oracle->exprinterpreter, oracle->conss[c]->exprtree);
2387 
2388  return evalcapability;
2389 }
2390 
2391 /** evaluates the objective function in a given point */
2393  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2394  const SCIP_Real* x, /**< point where to evaluate */
2395  SCIP_Real* objval /**< pointer to store objective value */
2396  )
2397 {
2398  assert(oracle != NULL);
2399 
2400  SCIPdebugMessage("%p eval obj value\n", (void*)oracle);
2401 
2402  SCIP_CALL_QUIET( evalFunctionValue(oracle, oracle->objective, x, objval) );
2403 
2404  assert(oracle->objective->lhs == oracle->objective->rhs); /*lint !e777*/
2405  *objval += oracle->objective->lhs;
2406 
2407  return SCIP_OKAY;
2408 }
2409 
2410 /** evaluates one constraint function in a given point */
2412  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2413  int considx, /**< index of constraint to evaluate */
2414  const SCIP_Real* x, /**< point where to evaluate */
2415  SCIP_Real* conval /**< pointer to store constraint value */
2416  )
2417 {
2418  assert(oracle != NULL);
2419  assert(x != NULL || oracle->nvars == 0);
2420  assert(conval != NULL);
2421 
2422  SCIPdebugMessage("%p eval cons value\n", (void*)oracle);
2423 
2424  SCIP_CALL_QUIET( evalFunctionValue(oracle, oracle->conss[considx], x, conval) );
2425 
2426  return SCIP_OKAY;
2427 }
2428 
2429 /** evaluates all constraint functions in a given point */
2431  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2432  const SCIP_Real* x, /**< point where to evaluate */
2433  SCIP_Real* convals /**< buffer to store constraint values */
2434  )
2435 {
2436  int i;
2437 
2438  SCIPdebugMessage("%p eval cons values\n", (void*)oracle);
2439 
2440  assert(oracle != NULL);
2441  assert(x != NULL || oracle->nvars == 0);
2442  assert(convals != NULL);
2443 
2444  for( i = 0; i < oracle->nconss; ++i )
2445  {
2446  SCIP_CALL_QUIET( evalFunctionValue(oracle, oracle->conss[i], x, &convals[i]) );
2447  }
2448 
2449  return SCIP_OKAY;
2450 }
2451 
2452 /** computes the objective gradient in a given point */
2454  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2455  const SCIP_Real* x, /**< point where to evaluate */
2456  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
2457  SCIP_Real* objval, /**< pointer to store objective value */
2458  SCIP_Real* objgrad /**< pointer to store (dense) objective gradient */
2459  )
2460 {
2461  assert(oracle != NULL);
2462 
2463  SCIPdebugMessage("%p eval obj grad\n", (void*)oracle);
2464 
2465  SCIP_CALL_QUIET( evalFunctionGradient(oracle, oracle->objective, x, isnewx, objval, objgrad) );
2466 
2467  assert(oracle->objective->lhs == oracle->objective->rhs); /*lint !e777*/
2468  *objval += oracle->objective->lhs;
2469 
2470  return SCIP_OKAY;
2471 }
2472 
2473 /** computes a constraints gradient in a given point */
2475  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2476  const int considx, /**< index of constraint to compute gradient for */
2477  const SCIP_Real* x, /**< point where to evaluate */
2478  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
2479  SCIP_Real* conval, /**< pointer to store constraint value */
2480  SCIP_Real* congrad /**< pointer to store (dense) constraint gradient */
2481  )
2482 {
2483  assert(oracle != NULL);
2484  assert(x != NULL || oracle->nvars == 0);
2485  assert(conval != NULL);
2486 
2487  SCIPdebugMessage("%p eval cons grad\n", (void*)oracle);
2488 
2489  SCIP_CALL_QUIET( evalFunctionGradient(oracle, oracle->conss[considx], x, isnewx, conval, congrad) );
2490 
2491  return SCIP_OKAY;
2492 }
2493 
2494 /** gets sparsity pattern (rowwise) of Jacobian matrix
2495  *
2496  * Note that internal data is returned in *offset and *col, thus the user does not need to allocate memory there.
2497  * Adding or deleting constraints destroys the sparsity structure and make another call to this function necessary.
2498  */
2500  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2501  const int** offset, /**< pointer to store pointer that stores the offsets to each rows sparsity pattern in col, can be NULL */
2502  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 */
2503  )
2504 {
2505  SCIP_NLPIORACLECONS* cons;
2506  SCIP_Bool* nzflag;
2507  int nnz;
2508  int maxnnz;
2509  int i;
2510  int j;
2511 
2512  assert(oracle != NULL);
2513 
2514  SCIPdebugMessage("%p get jacobian sparsity\n", (void*)oracle);
2515 
2516  if( oracle->jacoffsets != NULL )
2517  {
2518  assert(oracle->jaccols != NULL);
2519  if( offset != NULL )
2520  *offset = oracle->jacoffsets;
2521  if( col != NULL )
2522  *col = oracle->jaccols;
2523  return SCIP_OKAY;
2524  }
2525 
2526  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &oracle->jacoffsets, oracle->nconss + 1) );
2527 
2528  maxnnz = MIN(oracle->nvars, 10) * oracle->nconss; /* initial guess */
2529  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &oracle->jaccols, maxnnz) );
2530 
2531  if( maxnnz == 0 )
2532  {
2533  /* no variables */
2534  BMSclearMemoryArray(oracle->jacoffsets, oracle->nconss + 1);
2535  if( offset != NULL )
2536  *offset = oracle->jacoffsets;
2537  if( col != NULL )
2538  *col = oracle->jaccols;
2539  return SCIP_OKAY;
2540  }
2541  nnz = 0;
2542 
2543  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &nzflag, oracle->nvars) );
2544 
2545  for( i = 0; i < oracle->nconss; ++i )
2546  {
2547  oracle->jacoffsets[i] = nnz;
2548 
2549  cons = oracle->conss[i];
2550  assert(cons != NULL);
2551 
2552  if( cons->nquadelems == 0 && cons->exprtree == NULL )
2553  {
2554  /* for a linear constraint, we can just copy the linear indices from the constraint into the sparsity pattern */
2555  if( cons->nlinidxs > 0 )
2556  {
2557  SCIP_CALL( ensureIntArraySize(oracle->blkmem, &oracle->jaccols, &maxnnz, nnz + cons->nlinidxs) );
2558  BMScopyMemoryArray(&oracle->jaccols[nnz], cons->linidxs, cons->nlinidxs); /*lint !e866*/
2559  nnz += cons->nlinidxs;
2560  }
2561  continue;
2562  }
2563  else if( cons->nlinidxs == 0 && cons->nquadelems == 0 )
2564  {
2565  /* for a constraint with exprtree only, we can just copy the exprvaridxs from the constraint into the sparsity pattern */
2566  int nvars;
2567 
2568  assert(cons->exprtree != NULL); /* this had been the first case */
2569 
2570  nvars = SCIPexprtreeGetNVars(cons->exprtree);
2571  assert(cons->exprvaridxs != NULL || nvars == 0);
2572 
2573  if( nvars > 0 )
2574  {
2575  SCIP_CALL( ensureIntArraySize(oracle->blkmem, &oracle->jaccols, &maxnnz, nnz + nvars) );
2576  BMScopyMemoryArray(&oracle->jaccols[nnz], cons->exprvaridxs, nvars); /*lint !e866*/
2577  nnz += nvars;
2578  }
2579  continue;
2580  }
2581 
2582  /* check which variables appear in constraint i
2583  * @todo this could be done faster for very sparse constraint by assembling all appearing variables, sorting, and removing duplicates
2584  */
2585  BMSclearMemoryArray(nzflag, oracle->nvars); /*lint !e644*/
2586 
2587  for( j = 0; j < cons->nlinidxs; ++j )
2588  nzflag[cons->linidxs[j]] = TRUE;
2589 
2590  for( j = 0; j < cons->nquadelems; ++j )
2591  {
2592  nzflag[cons->quadelems[j].idx1] = TRUE;
2593  nzflag[cons->quadelems[j].idx2] = TRUE;
2594  }
2595 
2596  if( cons->exprvaridxs != NULL )
2597  {
2598  assert(cons->exprtree != NULL);
2599  for( j = SCIPexprtreeGetNVars(cons->exprtree)-1; j >= 0; --j )
2600  nzflag[cons->exprvaridxs[j]] = TRUE;
2601  }
2602 
2603  /* store variables indices in jaccols */
2604  for( j = 0; j < oracle->nvars; ++j )
2605  {
2606  if( nzflag[j] == FALSE )
2607  continue;
2608 
2609  SCIP_CALL( ensureIntArraySize(oracle->blkmem, &oracle->jaccols, &maxnnz, nnz + 1) );
2610  oracle->jaccols[nnz] = j;
2611  ++nnz;
2612  }
2613  }
2614 
2615  oracle->jacoffsets[oracle->nconss] = nnz;
2616 
2617  /* shrink jaccols array to nnz */
2618  if( nnz < maxnnz )
2619  {
2620  SCIP_ALLOC( BMSreallocBlockMemoryArray(oracle->blkmem, &oracle->jaccols, maxnnz, nnz) );
2621  }
2622 
2623  BMSfreeBlockMemoryArray(oracle->blkmem, &nzflag, oracle->nvars);
2624 
2625  if( offset != NULL )
2626  *offset = oracle->jacoffsets;
2627  if( col != NULL )
2628  *col = oracle->jaccols;
2629 
2630  return SCIP_OKAY;
2631 }
2632 
2633 /** evaluates the Jacobi matrix in a given point
2634  *
2635  * The values in the Jacobi matrix are returned in the same order as specified by the offset and col arrays obtained by SCIPnlpiOracleGetJacobianSparsity.
2636  * The user need to call SCIPnlpiOracleGetJacobianSparsity at least ones before using this function.
2637  */
2639  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2640  const SCIP_Real* x, /**< point where to evaluate */
2641  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
2642  SCIP_Real* convals, /**< pointer to store constraint values, can be NULL */
2643  SCIP_Real* jacobi /**< pointer to store sparse jacobian values */
2644  )
2645 {
2646  SCIP_NLPIORACLECONS* cons;
2647  SCIP_RETCODE retcode;
2648  SCIP_Real* grad;
2649  SCIP_Real* xx;
2650  SCIP_Real nlval;
2651  int i;
2652  int j;
2653  int k;
2654  int l;
2655 
2656  SCIPdebugMessage("%p eval jacobian\n", (void*)oracle);
2657 
2658  assert(oracle != NULL);
2659  assert(jacobi != NULL);
2660 
2661  assert(oracle->jacoffsets != NULL);
2662  assert(oracle->jaccols != NULL);
2663 
2664  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &grad, oracle->nvars) );
2665  xx = NULL;
2666 
2667  retcode = SCIP_OKAY;
2668 
2669  j = oracle->jacoffsets[0];
2670  k = 0;
2671  for( i = 0; i < oracle->nconss; ++i )
2672  {
2673  cons = oracle->conss[i];
2674  assert(cons != NULL);
2675 
2676  if( cons->nquadelems == 0 && cons->exprtree == NULL )
2677  {
2678  /* for a linear constraint, we can just copy the linear coefs from the constraint into the jacobian */
2679  if( cons->nlinidxs > 0 )
2680  {
2681  BMScopyMemoryArray(&jacobi[k], cons->lincoefs, cons->nlinidxs); /*lint !e866*/
2682  j += cons->nlinidxs;
2683  k += cons->nlinidxs;
2684  }
2685  assert(j == oracle->jacoffsets[i+1]);
2686  continue;
2687  }
2688 
2689  if( cons->nlinidxs == 0 && cons->nquadelems == 0 )
2690  {
2691  /* for a constraint with exprtree only, we can just copy gradient of the exprtree from the constraint into jacobian */
2692  int nvars;
2693 
2694  assert(cons->exprtree != NULL); /* this had been the first case */
2695 
2696  nvars = SCIPexprtreeGetNVars(cons->exprtree);
2697  assert(nvars <= oracle->nvars);
2698  assert(cons->exprvaridxs != NULL || nvars == 0);
2699 
2700  if( nvars > 0 )
2701  {
2702  if( isnewx )
2703  {
2704  if( xx == NULL )
2705  {
2706  /* if no xx yet, alloc it; make it large enough in case we need it again */
2707  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &xx, oracle->nvars) );
2708  }
2709  for( l = 0; l < nvars; ++l )
2710  xx[l] = x[cons->exprvaridxs[l]]; /*lint !e613*/
2711  }
2712 
2713  SCIPdebugMessage("eval gradient of ");
2714  SCIPdebug( if( isnewx ) {printf("\nx ="); for( l = 0; l < nvars; ++l) printf(" %g", xx[l]); /*lint !e613*/ printf("\n");} )
2715 
2716  SCIP_CALL( SCIPexprintGrad(oracle->exprinterpreter, cons->exprtree, xx, isnewx, &nlval, grad) ); /*lint !e644*/
2717 
2718  SCIPdebug( printf("g ="); for( l = 0; l < nvars; ++l) printf(" %g", grad[l]); printf("\n"); )
2719 
2720  if( nlval != nlval || ABS(nlval) >= oracle->infinity ) /*lint !e777*/
2721  {
2722  SCIPdebugMessage("gradient evaluation yield invalid function value %g\n", nlval);
2723  retcode = SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
2724  break;
2725  }
2726  else
2727  {
2728  if( convals != NULL )
2729  convals[i] = nlval;
2730  for( l = 0; l < nvars; ++l )
2731  {
2732  assert(oracle->jaccols[j+l] == cons->exprvaridxs[l]);
2733  if( grad[l] != grad[l] ) /*lint !e777*/
2734  {
2735  SCIPdebugMessage("gradient evaluation yield invalid gradient value %g\n", grad[l]);
2736  retcode = SCIP_INVALIDDATA; /* indicate that the function could not be evaluated at given point */
2737  break;
2738  }
2739  else
2740  jacobi[k++] = grad[l];
2741  }
2742  if( l < nvars )
2743  break;
2744  j += nvars;
2745  }
2746  }
2747  else if( convals != NULL )
2748  {
2749  SCIPdebugMessage("eval value of constant ");
2750 
2751  SCIP_CALL( SCIPexprintEval(oracle->exprinterpreter, cons->exprtree, NULL, &convals[i]) );
2752  }
2753  continue;
2754  }
2755 
2756  /* do dense eval @todo could do it sparse */
2757  retcode = SCIPnlpiOracleEvalConstraintGradient(oracle, i, x, isnewx, (convals ? &convals[i] : &nlval), grad);
2758  if( retcode != SCIP_OKAY )
2759  break;
2760 
2761  while( j < oracle->jacoffsets[i+1] )
2762  jacobi[k++] = grad[oracle->jaccols[j++]];
2763  }
2764 
2765  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &xx, oracle->nvars);
2766  BMSfreeBlockMemoryArray(oracle->blkmem, &grad, oracle->nvars);
2767 
2768  return retcode;
2769 }
2770 
2771 /** gets sparsity pattern of the Hessian matrix of the Lagrangian
2772  *
2773  * Note that internal data is returned in *offset and *col, thus the user must not allocate memory there.
2774  * Adding or deleting variables, objective, or constraints may destroy the sparsity structure and make another call to this function necessary.
2775  * Only elements of the lower left triangle and the diagonal are counted.
2776  */
2778  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2779  const int** offset, /**< pointer to store pointer that stores the offsets to each rows sparsity pattern in col, can be NULL */
2780  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 */
2781  )
2782 {
2783  int** colnz; /* nonzeros in Hessian corresponding to one column */
2784  int* collen; /* collen[i] is length of array colnz[i] */
2785  int* colnnz; /* colnnz[i] is number of entries in colnz[i] (<= collen[i]) */
2786  int nnz;
2787  int i;
2788  int j;
2789  int cnt;
2790 
2791  assert(oracle != NULL);
2792 
2793  SCIPdebugMessage("%p get hessian lag sparsity\n", (void*)oracle);
2794 
2795  if( oracle->heslagoffsets != NULL )
2796  {
2797  assert(oracle->heslagcols != NULL);
2798  if( offset != NULL )
2799  *offset = oracle->heslagoffsets;
2800  if( col != NULL )
2801  *col = oracle->heslagcols;
2802  return SCIP_OKAY;
2803  }
2804 
2805  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &oracle->heslagoffsets, oracle->nvars + 1) );
2806 
2807  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &colnz, oracle->nvars) );
2808  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &collen, oracle->nvars) );
2809  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &colnnz, oracle->nvars) );
2810  BMSclearMemoryArray(colnz, oracle->nvars); /*lint !e644*/
2811  BMSclearMemoryArray(collen, oracle->nvars); /*lint !e644*/
2812  BMSclearMemoryArray(colnnz, oracle->nvars); /*lint !e644*/
2813  nnz = 0;
2814 
2815  if( oracle->objective->nquadelems != 0 )
2816  {
2817  SCIP_CALL( hessLagSparsitySetNzFlagForQuad(oracle, colnz, collen, colnnz, &nnz, oracle->objective->nquadelems, oracle->objective->quadelems) );
2818  }
2819 
2820  if( oracle->objective->exprtree != NULL )
2821  {
2822  SCIP_CALL( hessLagSparsitySetNzFlagForExprtree(oracle, colnz, collen, colnnz, &nnz, oracle->objective->exprvaridxs, oracle->objective->exprtree, oracle->nvars) );
2823  }
2824 
2825  for( i = 0; i < oracle->nconss; ++i )
2826  {
2827  if( oracle->conss[i]->nquadelems != 0 )
2828  {
2829  SCIP_CALL( hessLagSparsitySetNzFlagForQuad(oracle, colnz, collen, colnnz, &nnz, oracle->conss[i]->nquadelems, oracle->conss[i]->quadelems) );
2830  }
2831 
2832  if( oracle->conss[i]->exprtree != NULL )
2833  {
2834  SCIP_CALL( hessLagSparsitySetNzFlagForExprtree(oracle, colnz, collen, colnnz, &nnz, oracle->conss[i]->exprvaridxs, oracle->conss[i]->exprtree, oracle->nvars) );
2835  }
2836  }
2837 
2838  SCIP_ALLOC( BMSallocBlockMemoryArray(oracle->blkmem, &oracle->heslagcols, nnz) );
2839 
2840  /* set hessian sparsity from colnz, colnnz */
2841  cnt = 0;
2842  for( i = 0; i < oracle->nvars; ++i )
2843  {
2844  oracle->heslagoffsets[i] = cnt;
2845  for( j = 0; j < colnnz[i]; ++j )
2846  {
2847  assert(cnt < nnz);
2848  oracle->heslagcols[cnt++] = colnz[i][j];
2849  }
2850  BMSfreeBlockMemoryArrayNull(oracle->blkmem, &colnz[i], collen[i]); /*lint !e866*/
2851  collen[i] = 0;
2852  }
2853  oracle->heslagoffsets[oracle->nvars] = cnt;
2854  assert(cnt == nnz);
2855 
2856  BMSfreeBlockMemoryArray(oracle->blkmem, &colnz, oracle->nvars);
2857  BMSfreeBlockMemoryArray(oracle->blkmem, &colnnz, oracle->nvars);
2858  BMSfreeBlockMemoryArray(oracle->blkmem, &collen, oracle->nvars);
2859 
2860  if( offset != NULL )
2861  *offset = oracle->heslagoffsets;
2862  if( col != NULL )
2863  *col = oracle->heslagcols;
2864 
2865  return SCIP_OKAY;
2866 }
2867 
2868 /** evaluates the Hessian matrix of the Lagrangian in a given point
2869  *
2870  * The values in the Hessian matrix are returned in the same order as specified by the offset and col arrays obtained by SCIPnlpiOracleGetHessianLagSparsity.
2871  * The user must call SCIPnlpiOracleGetHessianLagSparsity at least ones before using this function.
2872  * Only elements of the lower left triangle and the diagonal are computed.
2873  */
2875  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2876  const SCIP_Real* x, /**< point where to evaluate */
2877  SCIP_Bool isnewx, /**< has the point x changed since the last call to some evaluation function? */
2878  SCIP_Real objfactor, /**< weight for objective function */
2879  const SCIP_Real* lambda, /**< weights (Lagrangian multipliers) for the constraints */
2880  SCIP_Real* hessian /**< pointer to store sparse hessian values */
2881  )
2882 { /*lint --e{715}*/
2883  int i;
2884 
2885  assert(oracle != NULL);
2886  assert(x != NULL);
2887  assert(lambda != NULL || oracle->nconss == 0);
2888  assert(hessian != NULL);
2889 
2890  assert(oracle->heslagoffsets != NULL);
2891  assert(oracle->heslagcols != NULL);
2892 
2893  SCIPdebugMessage("%p eval hessian lag\n", (void*)oracle);
2894 
2895  for( i = oracle->heslagoffsets[oracle->nvars] - 1; i >= 0; --i )
2896  hessian[i] = 0.0;
2897 
2898  if( objfactor != 0.0 )
2899  {
2900  SCIP_CALL( hessLagAddQuad(objfactor, oracle->objective->nquadelems, oracle->objective->quadelems, oracle->heslagoffsets, oracle->heslagcols, hessian) );
2901  SCIP_CALL_QUIET( hessLagAddExprtree(oracle, objfactor, x, isnewx, oracle->objective->exprvaridxs, oracle->objective->exprtree, oracle->heslagoffsets, oracle->heslagcols, hessian) );
2902  }
2903 
2904  for( i = 0; i < oracle->nconss; ++i )
2905  {
2906  assert( lambda != NULL ); /* for lint */
2907  if( lambda[i] == 0.0 )
2908  continue;
2909  SCIP_CALL( hessLagAddQuad(lambda[i], oracle->conss[i]->nquadelems, oracle->conss[i]->quadelems, oracle->heslagoffsets, oracle->heslagcols, hessian) );
2910  SCIP_CALL_QUIET( hessLagAddExprtree(oracle, lambda[i], x, isnewx, oracle->conss[i]->exprvaridxs, oracle->conss[i]->exprtree, oracle->heslagoffsets, oracle->heslagcols, hessian) );
2911  }
2912 
2913  return SCIP_OKAY;
2914 }
2915 
2916 /** prints the problem to a file. */
2918  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2919  SCIP_MESSAGEHDLR* messagehdlr, /**< message handler */
2920  FILE* file /**< file to print to, or NULL for standard output */
2921  )
2922 { /*lint --e{777} */
2923  int i;
2924  SCIP_Real lhs;
2925  SCIP_Real rhs;
2926 
2927  assert(oracle != NULL);
2928 
2929  SCIPdebugMessage("%p print problem\n", (void*)oracle);
2930 
2931  if( file == NULL )
2932  file = stdout;
2933 
2934  SCIPmessageFPrintInfo(messagehdlr, file, "NLPI Oracle %s: %d variables and %d constraints\n", oracle->name ? oracle->name : "", oracle->nvars, oracle->nconss);
2935  for( i = 0; i < oracle->nvars; ++i )
2936  {
2937  if( oracle->varnames != NULL && oracle->varnames[i] != NULL )
2938  SCIPmessageFPrintInfo(messagehdlr, file, "%10s", oracle->varnames[i]);
2939  else
2940  SCIPmessageFPrintInfo(messagehdlr, file, "x%09d", i);
2941  SCIPmessageFPrintInfo(messagehdlr, file, ": [%8g, %8g]", oracle->varlbs[i], oracle->varubs[i]);
2942  if( oracle->vardegreesuptodate )
2943  SCIPmessageFPrintInfo(messagehdlr, file, "\t degree: %d", oracle->vardegrees[i]);
2944  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
2945  }
2946 
2947  SCIPmessageFPrintInfo(messagehdlr, file, "objective: ");
2948  SCIP_CALL( printFunction(oracle, messagehdlr, file, oracle->objective, FALSE, FALSE) );
2949  if( oracle->objective->lhs != 0.0 )
2950  SCIPmessageFPrintInfo(messagehdlr, file, "%+.20g", oracle->objective->lhs);
2951  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
2952 
2953  for( i = 0; i < oracle->nconss; ++i )
2954  {
2955  if( oracle->conss[i]->name != NULL )
2956  SCIPmessageFPrintInfo(messagehdlr, file, "%10s", oracle->conss[i]->name);
2957  else
2958  SCIPmessageFPrintInfo(messagehdlr, file, "con%07d", i);
2959 
2960  lhs = oracle->conss[i]->lhs;
2961  rhs = oracle->conss[i]->rhs;
2962  SCIPmessageFPrintInfo(messagehdlr, file, ": ");
2963  if( lhs > -oracle->infinity && rhs < oracle->infinity && lhs != rhs )
2964  SCIPmessageFPrintInfo(messagehdlr, file, "%.20g <= ", lhs);
2965 
2966  SCIP_CALL( printFunction(oracle, messagehdlr, file, oracle->conss[i], FALSE, FALSE) );
2967 
2968  if( lhs == rhs )
2969  SCIPmessageFPrintInfo(messagehdlr, file, " = %.20g", rhs);
2970  else if( rhs < oracle->infinity )
2971  SCIPmessageFPrintInfo(messagehdlr, file, " <= %.20g", rhs);
2972  else if( lhs > -oracle->infinity )
2973  SCIPmessageFPrintInfo(messagehdlr, file, " >= %.20g", lhs);
2974 
2975  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
2976  }
2977 
2978  return SCIP_OKAY;
2979 }
2980 
2981 /** prints the problem to a file in GAMS format
2982  * If there are variable (equation, resp.) names with more than 9 characters, then variable (equation, resp.) names are prefixed with an unique identifier.
2983  * This is to make it easier to identify variables solution output in the listing file.
2984  * Names with more than 64 characters are shorten to 64 letters due to GAMS limits.
2985  */
2987  SCIP_NLPIORACLE* oracle, /**< pointer to NLPIORACLE data structure */
2988  SCIP_Real* initval, /**< starting point values for variables or NULL */
2989  SCIP_MESSAGEHDLR* messagehdlr, /**< message handler */
2990  FILE* file /**< file to print to, or NULL for standard output */
2991  )
2992 { /*lint --e{777} */
2993  int i;
2994  int nllevel; /* level of nonlinearity of problem: linear = 0, quadratic, smooth nonlinear, nonsmooth */
2995  static const char* nllevelname[4] = { "LP", "QCP", "NLP", "DNLP" };
2996  const char* problemname;
2997  char namebuf[70];
2998  SCIP_Bool havelongvarnames;
2999  SCIP_Bool havelongequnames;
3000 
3001  SCIPdebugMessage("%p print problem gams\n", (void*)oracle);
3002 
3003  assert(oracle != NULL);
3004 
3005  if( file == NULL )
3006  file = stdout;
3007 
3008  nllevel = 0;
3009 
3010  havelongvarnames = FALSE;
3011  for( i = 0; i < oracle->nvars; ++i )
3012  if( oracle->varnames != NULL && oracle->varnames[i] != NULL && strlen(oracle->varnames[i]) > 9 )
3013  {
3014  havelongvarnames = TRUE;
3015  break;
3016  }
3017 
3018  havelongequnames = FALSE;
3019  for( i = 0; i < oracle->nconss; ++i )
3020  if( oracle->conss[i]->name && strlen(oracle->conss[i]->name) > 9 )
3021  {
3022  havelongequnames = TRUE;
3023  break;
3024  }
3025 
3026  SCIPmessageFPrintInfo(messagehdlr, file, "$offlisting\n");
3027  SCIPmessageFPrintInfo(messagehdlr, file, "$offdigit\n");
3028  SCIPmessageFPrintInfo(messagehdlr, file, "* NLPI Oracle Problem %s\n", oracle->name ? oracle->name : "");
3029  SCIPmessageFPrintInfo(messagehdlr, file, "Variables ");
3030  for( i = 0; i < oracle->nvars; ++i )
3031  {
3032  printName(namebuf, oracle->varnames != NULL ? oracle->varnames[i] : NULL, i, 'x', NULL, havelongvarnames);
3033  SCIPmessageFPrintInfo(messagehdlr, file, "%s, ", namebuf);
3034  if( i % 10 == 9 )
3035  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
3036  }
3037  SCIPmessageFPrintInfo(messagehdlr, file, "NLPIORACLEOBJVAR;\n\n");
3038  for( i = 0; i < oracle->nvars; ++i )
3039  {
3040  char* name;
3041  name = oracle->varnames != NULL ? oracle->varnames[i] : NULL;
3042  if( oracle->varlbs[i] == oracle->varubs[i] )
3043  {
3044  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
3045  SCIPmessageFPrintInfo(messagehdlr, file, "%s.fx = %.20g;\t", namebuf, oracle->varlbs[i]);
3046  }
3047  else
3048  {
3049  if( oracle->varlbs[i] > -oracle->infinity )
3050  {
3051  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
3052  SCIPmessageFPrintInfo(messagehdlr, file, "%s.lo = %.20g;\t", namebuf, oracle->varlbs[i]);
3053  }
3054  if( oracle->varubs[i] < oracle->infinity )
3055  {
3056  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
3057  SCIPmessageFPrintInfo(messagehdlr, file, "%s.up = %.20g;\t", namebuf, oracle->varubs[i]);
3058  }
3059  }
3060  if( initval != NULL )
3061  {
3062  printName(namebuf, name, i, 'x', NULL, havelongvarnames);
3063  SCIPmessageFPrintInfo(messagehdlr, file, "%s.l = %.20g;\t", namebuf, initval[i]);
3064  }
3065  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
3066  }
3067  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
3068 
3069  SCIPmessageFPrintInfo(messagehdlr, file, "Equations ");
3070  for( i = 0; i < oracle->nconss; ++i )
3071  {
3072  printName(namebuf, oracle->conss[i]->name, i, 'e', NULL, havelongequnames);
3073  SCIPmessageFPrintInfo(messagehdlr, file, "%s, ", namebuf);
3074 
3075  if( oracle->conss[i]->lhs > -oracle->infinity && oracle->conss[i]->rhs < oracle->infinity && oracle->conss[i]->lhs != oracle->conss[i]->rhs )
3076  {
3077  /* ranged row: add second constraint */
3078  printName(namebuf, oracle->conss[i]->name, i, 'e', "_RNG", havelongequnames);
3079  SCIPmessageFPrintInfo(messagehdlr, file, "%s, ", namebuf);
3080  }
3081  if( i % 10 == 9 )
3082  SCIPmessageFPrintInfo(messagehdlr, file, "\n");
3083  }
3084  SCIPmessageFPrintInfo(messagehdlr, file, "NLPIORACLEOBJ;\n\n");
3085 
3086  SCIPmessageFPrintInfo(messagehdlr, file, "NLPIORACLEOBJ.. NLPIORACLEOBJVAR =E= ");
3087  SCIP_CALL( printFunction(oracle, messagehdlr, file, oracle->objective, havelongvarnames, havelongequnames) );
3088  if( oracle->objective->lhs != 0.0 )
3089  SCIPmessageFPrintInfo(messagehdlr, file, "%+.20g", oracle->objective->lhs);
3090  SCIPmessageFPrintInfo(messagehdlr, file, ";\n");
3091 
3092  for( i = 0; i < oracle->nconss; ++i )
3093  {
3094  SCIP_Real lhs;
3095  SCIP_Real rhs;
3096 
3097  printName(namebuf, oracle->conss[i]->name, i, 'e', NULL, havelongequnames);
3098  SCIPmessageFPrintInfo(messagehdlr, file, "%s.. ", namebuf);
3099 
3100  SCIP_CALL( printFunction(oracle, messagehdlr, file, oracle->conss[i], havelongvarnames, havelongequnames) );
3101 
3102  lhs = oracle->conss[i]->lhs;
3103  rhs = oracle->conss[i]->rhs;
3104 
3105  if( lhs == rhs )
3106  SCIPmessageFPrintInfo(messagehdlr, file, " =E= %.20g", rhs);
3107  else if( rhs < oracle->infinity )
3108  SCIPmessageFPrintInfo(messagehdlr, file, " =L= %.20g", rhs);
3109  else if( lhs > -oracle->infinity )
3110  SCIPmessageFPrintInfo(messagehdlr, file, " =G= %.20g", lhs);
3111  else
3112  SCIPmessageFPrintInfo(messagehdlr, file, " =N= 0");
3113  SCIPmessageFPrintInfo(messagehdlr, file, ";\n");
3114 
3115  if( lhs > -oracle->infinity && rhs < oracle->infinity && lhs != rhs )
3116  {
3117  printName(namebuf, oracle->conss[i]->name, i, 'e', "_RNG", havelongequnames);
3118  SCIPmessageFPrintInfo(messagehdlr, file, "%s.. ", namebuf);
3119 
3120  SCIP_CALL( printFunction(oracle, messagehdlr, file, oracle->conss[i], havelongvarnames, havelongequnames) );
3121 
3122  SCIPmessageFPrintInfo(messagehdlr, file, " =G= %.20g;\n", lhs);
3123  }
3124 
3125  if( nllevel <= 0 && oracle->conss[i]->nquadelems > 0 )
3126  nllevel = 1;
3127  if( nllevel <= 1 && oracle->conss[i]->exprtree != NULL )
3128  nllevel = 2;
3129  if( nllevel <= 2 && oracle->conss[i]->exprtree != NULL && exprIsNonSmooth(SCIPexprtreeGetRoot(oracle->conss[i]->exprtree)) )
3130  nllevel = 3;
3131  }
3132 
3133  problemname = oracle->name ? oracle->name : "m";
3134 
3135  SCIPmessageFPrintInfo(messagehdlr, file, "Model %s / all /;\n", problemname);
3136  SCIPmessageFPrintInfo(messagehdlr, file, "option limrow = 0;\n");
3137  SCIPmessageFPrintInfo(messagehdlr, file, "option limcol = 0;\n");
3138  SCIPmessageFPrintInfo(messagehdlr, file, "Solve %s minimizing NLPIORACLEOBJVAR using %s;\n", problemname, nllevelname[nllevel]);
3139 
3140  return SCIP_OKAY;
3141 }
3142 
3143 /**@} */
SCIP_RETCODE SCIPexprintCompile(SCIP_EXPRINT *exprint, SCIP_EXPRTREE *tree)
void SCIPquadelemSort(SCIP_QUADELEM *quadelems, int nquadelems)
Definition: expr.c:9030
SCIP_RETCODE SCIPnlpiOraclePrintProblemGams(SCIP_NLPIORACLE *oracle, SCIP_Real *initval, SCIP_MESSAGEHDLR *messagehdlr, FILE *file)
Definition: nlpioracle.c:2986
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValues(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *convals)
Definition: nlpioracle.c:2430
#define BMSallocBlockMemory(mem, ptr)
Definition: memory.h:406
#define BMSfreeBlockMemoryArrayNull(mem, ptr, num)
Definition: memory.h:422
const SCIP_Real * SCIPnlpiOracleGetVarUbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2220
static void updateVariableDegrees(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:585
methods to interpret (evaluate) an expression tree "fast"
SCIP_RETCODE SCIPnlpiOracleGetJacobianSparsity(SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2499
void SCIPexprtreePrint(SCIP_EXPRTREE *tree, SCIP_MESSAGEHDLR *messagehdlr, FILE *file, const char **varnames, const char **paramnames)
Definition: expr.c:8596
SCIP_EXPRINTCAPABILITY SCIPnlpiOracleGetEvalCapability(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2370
SCIP_RETCODE SCIPnlpiOracleGetHessianLagSparsity(SCIP_NLPIORACLE *oracle, const int **offset, const int **col)
Definition: nlpioracle.c:2777
int SCIPnlpiOracleGetNVars(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2190
SCIP_EXPROP SCIPexprGetOperator(SCIP_EXPR *expr)
Definition: expr.c:5573
static SCIP_Bool exprIsNonSmooth(SCIP_EXPR *expr)
Definition: nlpioracle.c:1284
#define NULL
Definition: lpi_spx.cpp:130
static SCIP_RETCODE ensureConssSize(SCIP_NLPIORACLE *oracle, int minsize)
Definition: nlpioracle.c:141
static SCIP_RETCODE ensureIntArraySize(BMS_BLKMEM *blkmem, int **intarray, int *len, int minsize)
Definition: nlpioracle.c:218
SCIP_RETCODE SCIPnlpiOracleSetObjective(SCIP_NLPIORACLE *oracle, const SCIP_Real constant, int nlin, const int *lininds, const SCIP_Real *linvals, int nquadelems, const SCIP_QUADELEM *quadelems, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree)
Definition: nlpioracle.c:1602
static SCIP_RETCODE hessLagAddExprtree(SCIP_NLPIORACLE *oracle, SCIP_Real weight, const SCIP_Real *x, SCIP_Bool new_x, int *exprvaridx, SCIP_EXPRTREE *exprtree, int *hesoffset, int *hescol, SCIP_Real *values)
Definition: nlpioracle.c:1085
static void freeVariables(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:519
SCIP_Real * lincoefs
Definition: nlpioracle.c:44
methods to store an NLP and request function, gradient, and hessian values
#define BMSclearMemory(ptr)
Definition: memory.h:84
static SCIP_RETCODE createConstraint(BMS_BLKMEM *blkmem, SCIP_NLPIORACLECONS **cons, int nlinidxs, const int *linidxs, const SCIP_Real *lincoefs, int nquadelems, const SCIP_QUADELEM *quadelems, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree, SCIP_Real lhs, SCIP_Real rhs, const char *name)
Definition: nlpioracle.c:341
#define FALSE
Definition: def.h:56
static SCIP_RETCODE ensureVarsSize(SCIP_NLPIORACLE *oracle, int minsize)
Definition: nlpioracle.c:110
SCIP_Real SCIPnlpiOracleGetConstraintLhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2271
static SCIP_RETCODE evalFunctionGradient(SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS *cons, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *val, SCIP_Real *grad)
Definition: nlpioracle.c:816
static void printName(char *buffer, char *name, int idx, char prefix, const char *suffix, SCIP_Bool longnames)
Definition: nlpioracle.c:1188
static void freeConstraint(BMS_BLKMEM *blkmem, SCIP_NLPIORACLECONS **cons)
Definition: nlpioracle.c:420
#define TRUE
Definition: def.h:55
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:53
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveValue(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Real *objval)
Definition: nlpioracle.c:2392
#define SCIP_CALL(x)
Definition: def.h:266
SCIP_RETCODE SCIPexprtreeCopy(BMS_BLKMEM *blkmem, SCIP_EXPRTREE **targettree, SCIP_EXPRTREE *sourcetree)
Definition: expr.c:8652
static SCIP_RETCODE evalFunctionValue(SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS *cons, const SCIP_Real *x, SCIP_Real *val)
Definition: nlpioracle.c:739
static SCIP_RETCODE hessLagAddQuad(SCIP_Real weight, int length, SCIP_QUADELEM *quadelems, int *hesoffset, int *hescol, SCIP_Real *values)
Definition: nlpioracle.c:1050
static void invalidateJacobiSparsity(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:248
static SCIP_RETCODE printFunction(SCIP_NLPIORACLE *oracle, SCIP_MESSAGEHDLR *messagehdlr, FILE *file, SCIP_NLPIORACLECONS *cons, SCIP_Bool longvarnames, SCIP_Bool longequnames)
Definition: nlpioracle.c:1218
#define SCIPdebugMessage
Definition: pub_message.h:77
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:85
static void mapIndices(int *indexmap, int nindices, int *indices)
Definition: nlpioracle.c:612
SCIP_Bool SCIPquadelemSortedFind(SCIP_QUADELEM *quadelems, int idx1, int idx2, int nquadelems, int *pos)
Definition: expr.c:9055
unsigned int SCIP_EXPRINTCAPABILITY
void SCIPexprtreeSetParamVal(SCIP_EXPRTREE *tree, int paramidx, SCIP_Real paramval)
Definition: expr.c:8482
SCIP_RETCODE SCIPnlpiOracleEvalObjectiveGradient(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *objval, SCIP_Real *objgrad)
Definition: nlpioracle.c:2453
const char * SCIPnlpiOracleGetProblemName(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1427
SCIP_RETCODE SCIPnlpiOracleDelConsSet(SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1819
int * heslagoffsets
Definition: nlpioracle.c:80
SCIP_RETCODE SCIPnlpiOracleChgObjConstant(SCIP_NLPIORACLE *oracle, SCIP_Real objconstant)
Definition: nlpioracle.c:2174
SCIP_RETCODE SCIPexprintCreate(BMS_BLKMEM *blkmem, SCIP_EXPRINT **exprint)
int SCIPnlpiOracleGetNConstraints(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2200
static void mapIndicesQuad(int *indexmap, int nelems, SCIP_QUADELEM *elems)
Definition: nlpioracle.c:712
static void clearDeletedQuadElements(BMS_BLKMEM *blkmem, SCIP_QUADELEM **quadelems, int *nquadelems)
Definition: nlpioracle.c:674
public methods for expressions, expression trees, expression graphs, and related stuff ...
#define SCIP_DEFAULT_EPSILON
Definition: def.h:133
static int calcGrowSize(int num)
Definition: nlpioracle.c:94
SCIP_Real coef
Definition: type_expr.h:102
int * jacoffsets
Definition: nlpioracle.c:77
#define BMSduplicateBlockMemoryArray(mem, ptr, source, num)
Definition: memory.h:416
int SCIPnlpiOracleGetConstraintDegree(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2312
SCIP_RETCODE SCIPnlpiOracleChgExprtree(SCIP_NLPIORACLE *oracle, int considx, const int *exprvaridxs, const SCIP_EXPRTREE *exprtree)
Definition: nlpioracle.c:2092
static void sortLinearCoefficients(int *nidxs, int *idxs, SCIP_Real *coefs)
Definition: nlpioracle.c:292
SCIP_RETCODE SCIPnlpiOracleCreate(BMS_BLKMEM *blkmem, SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1320
#define SCIPerrorMessage
Definition: pub_message.h:45
const SCIP_Real * SCIPnlpiOracleGetVarLbs(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2210
SCIP_RETCODE SCIPnlpiOracleChgVarBounds(SCIP_NLPIORACLE *oracle, int nvars, const int *indices, const SCIP_Real *lbs, const SCIP_Real *ubs)
Definition: nlpioracle.c:1639
static void invalidateHessianLagSparsity(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:271
SCIP_RETCODE SCIPnlpiOracleAddVars(SCIP_NLPIORACLE *oracle, int nvars, const SCIP_Real *lbs, const SCIP_Real *ubs, const char **varnames)
Definition: nlpioracle.c:1439
#define BMSallocMemory(ptr)
Definition: memory.h:74
int SCIPexprtreeGetNVars(SCIP_EXPRTREE *tree)
Definition: expr.c:8452
static SCIP_RETCODE ensureConsLinSize(BMS_BLKMEM *blkmem, SCIP_NLPIORACLECONS *cons, int minsize)
Definition: nlpioracle.c:165
#define BMSfreeMemory(ptr)
Definition: memory.h:100
SCIP_RETCODE SCIPexprintHessianDense(SCIP_EXPRINT *exprint, SCIP_EXPRTREE *tree, SCIP_Real *varvals, SCIP_Bool new_varvals, SCIP_Real *val, SCIP_Real *hessian)
SCIP_Real SCIPnlpiOracleGetConstraintRhs(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2284
SCIP_EXPR * SCIPexprtreeGetRoot(SCIP_EXPRTREE *tree)
Definition: expr.c:8442
SCIP_Bool SCIPsortedvecFindInt(int *intarray, int val, int len, int *pos)
SCIP_Real SCIPnlpiOracleGetInfinity(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:1392
SCIP_RETCODE SCIPnlpiOracleFree(SCIP_NLPIORACLE **oracle)
Definition: nlpioracle.c:1347
static SCIP_RETCODE moveVariable(SCIP_NLPIORACLE *oracle, int fromidx, int toidx)
Definition: nlpioracle.c:479
static void clearDeletedLinearElements(BMS_BLKMEM *blkmem, int **linidxs, SCIP_Real **coefs, int *nidxs)
Definition: nlpioracle.c:632
SCIP_RETCODE SCIPnlpiOracleEvalConstraintGradient(SCIP_NLPIORACLE *oracle, const int considx, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *conval, SCIP_Real *congrad)
Definition: nlpioracle.c:2474
#define SCIP_EXPRINTCAPABILITY_ALL
public data structures and miscellaneous methods
static SCIP_RETCODE hessLagSparsitySetNzFlagForQuad(SCIP_NLPIORACLE *oracle, int **colnz, int *collen, int *colnnz, int *nzcount, int length, SCIP_QUADELEM *quadelems)
Definition: nlpioracle.c:939
SCIP_RETCODE SCIPnlpiOraclePrintProblem(SCIP_NLPIORACLE *oracle, SCIP_MESSAGEHDLR *messagehdlr, FILE *file)
Definition: nlpioracle.c:2917
SCIP_EXPR ** SCIPexprGetChildren(SCIP_EXPR *expr)
Definition: expr.c:5593
char ** varnames
Definition: nlpioracle.c:67
#define SCIP_Bool
Definition: def.h:53
SCIP_NLPIORACLECONS * objective
Definition: nlpioracle.c:75
#define BMSallocBlockMemoryArray(mem, ptr, num)
Definition: memory.h:408
SCIP_RETCODE SCIPnlpiOracleSetProblemName(SCIP_NLPIORACLE *oracle, const char *name)
Definition: nlpioracle.c:1404
static void updateVariableDegreesCons(SCIP_NLPIORACLE *oracle, SCIP_NLPIORACLECONS *cons)
Definition: nlpioracle.c:552
SCIP_RETCODE SCIPnlpiOracleAddConstraints(SCIP_NLPIORACLE *oracle, int nconss, const SCIP_Real *lhss, const SCIP_Real *rhss, const int *nlininds, int *const *lininds, SCIP_Real *const *linvals, const int *nquadelems, SCIP_QUADELEM *const *quadelems, int *const *exprvaridxs, SCIP_EXPRTREE *const *exprtrees, const char **consnames)
Definition: nlpioracle.c:1524
#define BMSfreeBlockMemoryArray(mem, ptr, num)
Definition: memory.h:421
void SCIPquadelemSqueeze(SCIP_QUADELEM *quadelems, int nquadelems, int *nquadelemsnew)
Definition: expr.c:9107
#define MAX(x, y)
Definition: tclique_def.h:75
#define SCIP_CALL_QUIET(x)
Definition: def.h:241
int SCIPexprGetNChildren(SCIP_EXPR *expr)
Definition: expr.c:5583
SCIP_Real * varubs
Definition: nlpioracle.c:66
SCIP_RETCODE SCIPnlpiOracleEvalConstraintValue(SCIP_NLPIORACLE *oracle, int considx, const SCIP_Real *x, SCIP_Real *conval)
Definition: nlpioracle.c:2411
SCIP_Real * varlbs
Definition: nlpioracle.c:65
#define BMSfreeBlockMemory(mem, ptr)
Definition: memory.h:419
SCIP_RETCODE SCIPexprtreeFree(SCIP_EXPRTREE **tree)
Definition: expr.c:8692
SCIP_RETCODE SCIPnlpiOracleChgLinearCoefs(SCIP_NLPIORACLE *oracle, int considx, int nentries, const int *varidxs, const SCIP_Real *newcoefs)
Definition: nlpioracle.c:1897
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:89
SCIP_RETCODE SCIPnlpiOracleChgQuadCoefs(SCIP_NLPIORACLE *oracle, int considx, int nquadelems, const SCIP_QUADELEM *quadelems)
Definition: nlpioracle.c:1994
static SCIP_RETCODE hessLagSparsitySetNzFlagForExprtree(SCIP_NLPIORACLE *oracle, int **colnz, int *collen, int *colnnz, int *nzcount, int *exprvaridx, SCIP_EXPRTREE *exprtree, int dim)
Definition: nlpioracle.c:979
SCIP_RETCODE SCIPnlpiOracleChgConsSides(SCIP_NLPIORACLE *oracle, int nconss, const int *indices, const SCIP_Real *lhss, const SCIP_Real *rhss)
Definition: nlpioracle.c:1675
SCIP_Real infinity
Definition: nlpioracle.c:60
static void freeConstraints(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:453
static SCIP_RETCODE ensureConsQuadSize(BMS_BLKMEM *blkmem, SCIP_NLPIORACLECONS *cons, int minsize)
Definition: nlpioracle.c:192
SCIP_RETCODE SCIPnlpiOracleDelVarSet(SCIP_NLPIORACLE *oracle, int *delstats)
Definition: nlpioracle.c:1709
char * SCIPnlpiOracleGetConstraintName(SCIP_NLPIORACLE *oracle, int considx)
Definition: nlpioracle.c:2297
#define SCIP_DEFAULT_INFINITY
Definition: def.h:132
SCIP_QUADELEM * quadelems
Definition: nlpioracle.c:48
#define REALABS(x)
Definition: def.h:151
int * vardegrees
Definition: nlpioracle.c:68
SCIP_RETCODE SCIPnlpiOracleChgExprParam(SCIP_NLPIORACLE *oracle, int considx, int paramidx, SCIP_Real paramval)
Definition: nlpioracle.c:2149
int SCIPnlpiOracleGetVarDegree(SCIP_NLPIORACLE *oracle, int varidx)
Definition: nlpioracle.c:2242
SCIP_RETCODE SCIPnlpiOracleEvalHessianLag(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real objfactor, const SCIP_Real *lambda, SCIP_Real *hessian)
Definition: nlpioracle.c:2874
public methods for message output
SCIP_Bool vardegreesuptodate
Definition: nlpioracle.c:69
void SCIPmessageFPrintInfo(SCIP_MESSAGEHDLR *messagehdlr, FILE *file, const char *formatstr,...)
Definition: message.c:602
#define SCIP_Real
Definition: def.h:127
#define MIN(x, y)
Definition: memory.c:67
SCIP_EXPRINT * exprinterpreter
Definition: nlpioracle.c:84
const char * SCIPexprintGetName(void)
SCIP_RETCODE SCIPnlpiOracleSetInfinity(SCIP_NLPIORACLE *oracle, SCIP_Real infinity)
Definition: nlpioracle.c:1376
int SCIPexprtreeGetNParams(SCIP_EXPRTREE *tree)
Definition: expr.c:8462
BMS_BLKMEM * blkmem
Definition: nlpioracle.c:59
SCIP_RETCODE SCIPnlpiOracleEvalJacobian(SCIP_NLPIORACLE *oracle, const SCIP_Real *x, SCIP_Bool isnewx, SCIP_Real *convals, SCIP_Real *jacobi)
Definition: nlpioracle.c:2638
#define EPSLE(x, y, eps)
Definition: def.h:154
#define SCIPdebug(x)
Definition: pub_message.h:74
SCIP_RETCODE SCIPexprintHessianSparsityDense(SCIP_EXPRINT *exprint, SCIP_EXPRTREE *tree, SCIP_Real *varvals, SCIP_Bool *sparsity)
#define EPSEQ(x, y, eps)
Definition: def.h:152
SCIP_EXPRINTCAPABILITY SCIPexprintGetExprtreeCapability(SCIP_EXPRINT *exprint, SCIP_EXPRTREE *tree)
SCIP_EXPRTREE * exprtree
Definition: nlpioracle.c:51
struct BMS_BlkMem BMS_BLKMEM
Definition: memory.h:392
char ** SCIPnlpiOracleGetVarNames(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2230
#define SCIP_CALL_ABORT(x)
Definition: def.h:245
#define SCIP_ALLOC(x)
Definition: def.h:277
int * SCIPnlpiOracleGetVarDegrees(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2259
SCIP_RETCODE SCIPexprintGrad(SCIP_EXPRINT *exprint, SCIP_EXPRTREE *tree, SCIP_Real *varvals, SCIP_Bool new_varvals, SCIP_Real *val, SCIP_Real *gradient)
void SCIPsortedvecInsertInt(int *intarray, int keyval, int *len, int *pos)
#define BMSreallocBlockMemoryArray(mem, ptr, oldnum, newnum)
Definition: memory.h:412
SCIP_RETCODE SCIPexprintEval(SCIP_EXPRINT *exprint, SCIP_EXPRTREE *tree, SCIP_Real *varvals, SCIP_Real *val)
SCIP_RETCODE SCIPexprintFree(SCIP_EXPRINT **exprint)
int * heslagcols
Definition: nlpioracle.c:81
int SCIPnlpiOracleGetMaxDegree(SCIP_NLPIORACLE *oracle)
Definition: nlpioracle.c:2344
SCIP_NLPIORACLECONS ** conss
Definition: nlpioracle.c:73
void SCIPsortIntReal(int *intarray, SCIP_Real *realarray, int len)