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