Scippy

SCIP

Solving Constraint Integer Programs

nlpi_all.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 2002-2022 Zuse Institute Berlin */
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 nlpi_all.c
26  * @ingroup DEFPLUGINS_NLPI
27  * @brief NLP interface that uses all available NLP interfaces
28  * @author Benjamin Mueller
29  */
30 
31 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
32 
33 #include "scip/nlpi_all.h"
34 #include "scip/scip_mem.h"
35 #include "scip/scip_numerics.h"
36 #include "scip/scip_nlpi.h"
37 
38 #include <string.h>
39 
40 #define NLPI_NAME "all" /**< short concise name of solver */
41 #define NLPI_DESC "NLP interface that uses all available NLP interfaces" /**< description of solver */
42 #define NLPI_PRIORITY -3000 /**< priority of NLP solver */
43 
44 /*
45  * Data structures
46  */
47 
48 struct SCIP_NlpiData
49 {
50  SCIP_NLPI** nlpis; /**< array containing all nlpis */
51  int nnlpis; /**< total number of nlpis */
52 };
53 
55 {
56  SCIP_NLPIPROBLEM** nlpiproblems; /**< array containing all nlpi problems */
57  int nnlpiproblems; /**< total number of nlpi problems */
58  int bestidx; /**< index of NLP solver with the best solution */
59 };
60 
61 #ifdef SCIP_STATISTIC
62 static int _nnlps = 0; /**< number of NLPs that have been solved */
63 #endif
64 
65 /*
66  * Local methods
67  */
68 
69 /*
70  * Callback methods of NLP solver interface
71  */
72 
73 /** copy method of NLP interface (called when SCIP copies plugins) */
74 static
75 SCIP_DECL_NLPICOPY(nlpiCopyAll)
76 {
77  /* include NLPI */
79 
80  return SCIP_OKAY; /*lint !e527*/
81 } /*lint !e715*/
82 
83 /** destructor of NLP interface to free nlpi data */
84 static
85 SCIP_DECL_NLPIFREE(nlpiFreeAll)
86 {
87  assert(nlpi != NULL);
88  assert(nlpidata != NULL);
89  assert(*nlpidata != NULL);
90 
91  SCIPfreeBlockMemoryArrayNull(scip, &(*nlpidata)->nlpis, (*nlpidata)->nnlpis);
92  SCIPfreeBlockMemory(scip, nlpidata);
93  assert(*nlpidata == NULL);
94 
95  return SCIP_OKAY; /*lint !e527*/
96 } /*lint !e715*/
97 
98 /** creates a problem instance */
99 static
100 SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemAll)
101 {
102  SCIP_NLPIDATA* data;
103  int i;
104 
105  assert(nlpi != NULL);
106  assert(problem != NULL);
107 
108  data = SCIPnlpiGetData(nlpi);
109  assert(data != NULL);
110 
112 
113  /* initialize problem */
114  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(*problem)->nlpiproblems, data->nnlpis) );
115  (*problem)->nnlpiproblems = data->nnlpis;
116 
117  for( i = 0; i < data->nnlpis; ++i )
118  {
119  assert(data->nlpis[i] != NULL);
120  SCIP_CALL( SCIPcreateNlpiProblem(scip, data->nlpis[i], &((*problem)->nlpiproblems[i]), name) );
121  }
122 
123  return SCIP_OKAY;
124 } /*lint !e715*/
125 
126 /** free a problem instance */
127 static
128 SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemAll)
129 {
130  SCIP_NLPIDATA* data;
131  int i;
132 
133  assert(nlpi != NULL);
134  assert(problem != NULL);
135  assert(*problem != NULL);
136 
137  data = SCIPnlpiGetData(nlpi);
138  assert(data != NULL);
139 
140  for( i = 0; i < data->nnlpis; ++i )
141  {
142  assert(data->nlpis[i] != NULL);
143  SCIP_CALL( SCIPfreeNlpiProblem(scip, data->nlpis[i], &(*problem)->nlpiproblems[i]) );
144  }
145 
146  SCIPfreeBlockMemoryArrayNull(scip, &(*problem)->nlpiproblems, data->nnlpis);
147  SCIPfreeBlockMemory(scip, problem);
148 
149  return SCIP_OKAY;
150 } /*lint !e715*/
151 
152 /** add variables */
153 static
154 SCIP_DECL_NLPIADDVARS(nlpiAddVarsAll)
155 {
156  SCIP_NLPIDATA* nlpidata;
157  int i;
158 
159  nlpidata = SCIPnlpiGetData(nlpi);
160  assert(nlpidata != NULL);
161 
162  for( i = 0; i < nlpidata->nnlpis; ++i )
163  {
164  assert(nlpidata->nlpis[i] != NULL);
165  assert(problem->nlpiproblems[i] != NULL);
166 
167  SCIP_CALL( SCIPaddNlpiVars(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], nvars, lbs, ubs, varnames) );
168  }
169 
170  return SCIP_OKAY;
171 } /*lint !e715*/
172 
173 
174 /** add constraints */
175 static
176 SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsAll)
177 {
178  SCIP_NLPIDATA* nlpidata;
179  int i;
180 
181  nlpidata = SCIPnlpiGetData(nlpi);
182  assert(nlpidata != NULL);
183 
184  for( i = 0; i < nlpidata->nnlpis; ++i )
185  {
186  assert(nlpidata->nlpis[i] != NULL);
187  assert(problem->nlpiproblems[i] != NULL);
188 
189  SCIP_CALL( SCIPaddNlpiConstraints(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], nconss, lhss, rhss,
190  nlininds, lininds, linvals, exprs, names) );
191  }
192 
193  return SCIP_OKAY;
194 } /*lint !e715*/
195 
196 /** sets or overwrites objective, a minimization problem is expected */
197 static
198 SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveAll)
199 {
200  SCIP_NLPIDATA* nlpidata;
201  int i;
202 
203  nlpidata = SCIPnlpiGetData(nlpi);
204  assert(nlpidata != NULL);
205 
206  for( i = 0; i < nlpidata->nnlpis; ++i )
207  {
208  assert(nlpidata->nlpis[i] != NULL);
209  assert(problem->nlpiproblems[i] != NULL);
210 
211  SCIP_CALL( SCIPsetNlpiObjective(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], nlins, lininds, linvals, expr, constant) );
212  }
213 
214  return SCIP_OKAY;
215 } /*lint !e715*/
216 
217 /** change variable bounds */
218 static
219 SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsAll)
220 {
221  SCIP_NLPIDATA* nlpidata;
222  int i;
223 
224  nlpidata = SCIPnlpiGetData(nlpi);
225  assert(nlpidata != NULL);
226 
227  for( i = 0; i < nlpidata->nnlpis; ++i )
228  {
229  assert(nlpidata->nlpis[i] != NULL);
230  assert(problem->nlpiproblems[i] != NULL);
231 
232  SCIP_CALL( SCIPchgNlpiVarBounds(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], nvars, indices, lbs, ubs) );
233  }
234 
235  return SCIP_OKAY;
236 } /*lint !e715*/
237 
238 /** change constraint bounds */
239 static
240 SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesAll)
241 {
242  SCIP_NLPIDATA* nlpidata;
243  int i;
244 
245  nlpidata = SCIPnlpiGetData(nlpi);
246  assert(nlpidata != NULL);
247 
248  for( i = 0; i < nlpidata->nnlpis; ++i )
249  {
250  assert(nlpidata->nlpis[i] != NULL);
251  assert(problem->nlpiproblems[i] != NULL);
252 
253  SCIP_CALL( SCIPchgNlpiConsSides(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], nconss, indices, lhss, rhss) );
254  }
255 
256  return SCIP_OKAY;
257 } /*lint !e715*/
258 
259 /** delete a set of variables */
260 static
261 SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetAll)
262 {
263  SCIP_NLPIDATA* nlpidata;
264  int* tmpdstats;
265  int i;
266 
267  nlpidata = SCIPnlpiGetData(nlpi);
268  assert(nlpidata != NULL);
269 
270  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &tmpdstats, dstatssize) );
271 
272  for( i = 0; i < nlpidata->nnlpis; ++i )
273  {
274  assert(nlpidata->nlpis[i] != NULL);
275  assert(problem->nlpiproblems[i] != NULL);
276 
277  if( i < nlpidata->nnlpis -1 )
278  {
279  /* restore dstats entries */
280  BMScopyMemoryArray(tmpdstats, dstats, dstatssize);
281 
282  SCIP_CALL( SCIPdelNlpiVarSet(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], tmpdstats, dstatssize) );
283  }
284  else
285  {
286  /* NOTE this works only when all dstats array are the same after calling the nlpidelvarset callback
287  * As long as all solvers use the SCIP NLPI oracle to store the NLP problem data, this is the case.
288  * @TODO Assert that the returned dstats are all the same?
289  */
290  SCIP_CALL( SCIPdelNlpiVarSet(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], dstats, dstatssize) );
291  }
292  }
293 
294  SCIPfreeBlockMemoryArray(scip, &tmpdstats, dstatssize);
295 
296  return SCIP_OKAY;
297 } /*lint !e715*/
298 
299 /** delete a set of constraints */
300 static
301 SCIP_DECL_NLPIDELCONSSET( nlpiDelConstraintSetAll )
302 {
303  SCIP_NLPIDATA* nlpidata;
304  int* tmpdstats;
305  int i;
306 
307  nlpidata = SCIPnlpiGetData(nlpi);
308  assert(nlpidata != NULL);
309 
310  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &tmpdstats, dstatssize) );
311 
312  for( i = 0; i < nlpidata->nnlpis; ++i )
313  {
314  assert(nlpidata->nlpis[i] != NULL);
315  assert(problem->nlpiproblems[i] != NULL);
316 
317  if( i < nlpidata->nnlpis - 1 )
318  {
319  /* restore dstats entries */
320  BMScopyMemoryArray(tmpdstats, dstats, dstatssize);
321 
322  SCIP_CALL( SCIPdelNlpiConsSet(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], tmpdstats, dstatssize) );
323  }
324  else
325  {
326  /* NOTE this works only when all dstats array are the same after calling the nlpidelconsset callback
327  * As long as all solvers use the SCIP NLPI oracle to store the NLP problem data, this is the case.
328  * @TODO Assert that the returned dstats are all the same?
329  */
330  SCIP_CALL( SCIPdelNlpiConsSet(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], dstats, dstatssize) );
331  }
332  }
333 
334  SCIPfreeBlockMemoryArray(scip, &tmpdstats, dstatssize);
335 
336  return SCIP_OKAY;
337 } /*lint !e715*/
338 
339 /** changes (or adds) linear coefficients in a constraint or objective */
340 static
341 SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsAll)
342 {
343  SCIP_NLPIDATA* nlpidata;
344  int i;
345 
346  nlpidata = SCIPnlpiGetData(nlpi);
347  assert(nlpidata != NULL);
348 
349  for( i = 0; i < nlpidata->nnlpis; ++i )
350  {
351  assert(nlpidata->nlpis[i] != NULL);
352  assert(problem->nlpiproblems[i] != NULL);
353 
354  SCIP_CALL( SCIPchgNlpiLinearCoefs(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], idx, nvals, varidxs, vals) );
355  }
356 
357  return SCIP_OKAY;
358 } /*lint !e715*/
359 
360 /** replaces the expression of a constraint or objective */
361 static
362 SCIP_DECL_NLPICHGEXPR(nlpiChgExprAll)
363 {
364  SCIP_NLPIDATA* nlpidata;
365  int i;
366 
367  nlpidata = SCIPnlpiGetData(nlpi);
368  assert(nlpidata != NULL);
369 
370  for( i = 0; i < nlpidata->nnlpis; ++i )
371  {
372  assert(nlpidata->nlpis[i] != NULL);
373  assert(problem->nlpiproblems[i] != NULL);
374 
375  SCIP_CALL( SCIPchgNlpiExpr(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], idxcons, expr) );
376  }
377 
378  return SCIP_OKAY;
379 } /*lint !e715*/
380 
381 /** change the constant offset in the objective */
382 static
383 SCIP_DECL_NLPICHGOBJCONSTANT(nlpiChgObjConstantAll)
384 {
385  SCIP_NLPIDATA* nlpidata;
386  int i;
387 
388  nlpidata = SCIPnlpiGetData(nlpi);
389  assert(nlpidata != NULL);
390 
391  for( i = 0; i < nlpidata->nnlpis; ++i )
392  {
393  assert(nlpidata->nlpis[i] != NULL);
394  assert(problem->nlpiproblems[i] != NULL);
395 
396  SCIP_CALL( SCIPchgNlpiObjConstant(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], objconstant) );
397  }
398 
399  return SCIP_OKAY;
400 } /*lint !e715*/
401 
402 /** sets initial guess for primal variables */
403 static
404 SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessAll)
405 {
406  SCIP_NLPIDATA* nlpidata;
407  int i;
408 
409  nlpidata = SCIPnlpiGetData(nlpi);
410  assert(nlpidata != NULL);
411 
412  for( i = 0; i < nlpidata->nnlpis; ++i )
413  {
414  assert(nlpidata->nlpis[i] != NULL);
415  assert(problem->nlpiproblems[i] != NULL);
416 
417  SCIP_CALL( SCIPsetNlpiInitialGuess(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], primalvalues, consdualvalues,
418  varlbdualvalues, varubdualvalues) );
419  }
420 
421  return SCIP_OKAY;
422 } /*lint !e715*/
423 
424 /** tries to solve NLP */
425 static
426 SCIP_DECL_NLPISOLVE(nlpiSolveAll)
427 {
428  SCIP_NLPIDATA* nlpidata;
429  SCIP_NLPTERMSTAT besttermstat;
430  SCIP_NLPSOLSTAT bestsolstat;
431  SCIP_Real bestsolval;
432  int i;
433 
434  nlpidata = SCIPnlpiGetData(nlpi);
435  assert(nlpidata != NULL);
436 
437  /* use first solver per default */
438  problem->bestidx = 0;
439 
440  /* initialize best solution values */
441  besttermstat = SCIP_NLPTERMSTAT_OTHER;
442  bestsolstat = SCIP_NLPSOLSTAT_UNKNOWN;
443  bestsolval = SCIPinfinity(scip);
444 
445  for( i = 0; i < nlpidata->nnlpis; ++i )
446  {
447  SCIP_NLPTERMSTAT termstat;
448  SCIP_NLPSOLSTAT solstat;
449  SCIP_Real solval;
450  SCIP_Bool update;
451 
452  assert(nlpidata->nlpis[i] != NULL);
453  assert(problem->nlpiproblems[i] != NULL);
454 
455  /* solve NLP */
456  SCIP_CALL( SCIPsolveNlpiParam(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], param) );
457 
458  termstat = SCIPgetNlpiTermstat(scip, nlpidata->nlpis[i], problem->nlpiproblems[i]);
459  solstat = SCIPgetNlpiSolstat(scip, nlpidata->nlpis[i], problem->nlpiproblems[i]);
460  solval = SCIPinfinity(scip);
461  update = FALSE;
462 
463  /* collect solution value */
464  if( solstat <= SCIP_NLPSOLSTAT_FEASIBLE )
465  {
466  SCIP_CALL( SCIPgetNlpiSolution(scip, nlpidata->nlpis[i], problem->nlpiproblems[i],
467  NULL, NULL, NULL, NULL, &solval) );
468  assert(!SCIPisInfinity(scip, solval));
469  }
470 
471  /* better termination status -> update best solver */
472  if( termstat < besttermstat )
473  update = TRUE;
474 
475  /* no feasible solutions have been found so far -> update best solver */
476  else if( bestsolstat >= SCIP_NLPSOLSTAT_LOCINFEASIBLE && solstat <= SCIP_NLPSOLSTAT_LOCINFEASIBLE )
477  update = TRUE;
478 
479  /* use solver with the better solution value */
480  else if( solval < bestsolval )
481  update = TRUE;
482 
483  /* update best solver */
484  if( update )
485  {
486  besttermstat = termstat;
487  bestsolstat = solstat;
488  bestsolval = solval;
489  problem->bestidx = i;
490  }
491 
492 #ifdef SCIP_STATISTIC
493  {
494  SCIP_NLPSTATISTICS stats;
495 
496  SCIP_CALL( SCIPgetNlpiStatistics(scip, nlpidata->nlpis[i], problem->nlpiproblems[i], &stats) );
497 
498  SCIPstatisticMessage("%d solver %s termstat %d solstat %d solval %e iters %d time %g\n",
499  _nnlps, SCIPnlpiGetName(nlpidata->nlpis[i]), termstat, solstat, solval,
500  stats.niterations, stats.totaltime);
501  }
502 #endif
503 
504  /* don't try more NLP solvers if allowed time is exceeded or SCIP is asked to interrupt */
505  if( termstat == SCIP_NLPTERMSTAT_TIMELIMIT || termstat == SCIP_NLPTERMSTAT_INTERRUPT )
506  break;
507  }
508 
509 #ifdef SCIP_STATISTIC
510  ++_nnlps;
511 #endif
512 
513  return SCIP_OKAY;
514 } /*lint !e715*/
515 
516 /** gives solution status */
517 static
518 SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatAll)
519 {
520  SCIP_NLPIDATA* nlpidata;
521 
522  nlpidata = SCIPnlpiGetData(nlpi);
523  assert(nlpidata != NULL);
524  assert(nlpidata->nlpis != NULL);
525  assert(nlpidata->nlpis[problem->bestidx] != NULL);
526  assert(problem->nlpiproblems != NULL);
527  assert(problem->nlpiproblems[problem->bestidx] != NULL);
528 
529  /* return the solution status of the first nlpi */
530  return SCIPgetNlpiSolstat(scip, nlpidata->nlpis[problem->bestidx], problem->nlpiproblems[problem->bestidx]);
531 }
532 
533 /** gives termination reason */
534 static
535 SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatAll)
536 {
537  SCIP_NLPIDATA* nlpidata;
538 
539  nlpidata = SCIPnlpiGetData(nlpi);
540  assert(nlpidata != NULL);
541  assert(nlpidata->nlpis != NULL);
542  assert(nlpidata->nlpis[problem->bestidx] != NULL);
543  assert(problem->nlpiproblems != NULL);
544  assert(problem->nlpiproblems[problem->bestidx] != NULL);
545 
546  /* return the solution status of the first nlpi */
547  return SCIPgetNlpiTermstat(scip, nlpidata->nlpis[problem->bestidx], problem->nlpiproblems[problem->bestidx]);
548 }
549 
550 /** gives primal and dual solution values */
551 static
552 SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionAll)
553 {
554  SCIP_NLPIDATA* nlpidata;
555 
556  nlpidata = SCIPnlpiGetData(nlpi);
557  assert(nlpidata != NULL);
558  assert(nlpidata->nlpis != NULL);
559  assert(nlpidata->nlpis[problem->bestidx] != NULL);
560  assert(problem->nlpiproblems != NULL);
561  assert(problem->nlpiproblems[problem->bestidx] != NULL);
562 
563  /* return the solution status of the first nlpi */
564  SCIP_CALL( SCIPgetNlpiSolution(scip, nlpidata->nlpis[problem->bestidx], problem->nlpiproblems[problem->bestidx],
565  primalvalues, consdualvalues, varlbdualvalues, varubdualvalues, objval) );
566 
567  return SCIP_OKAY;
568 }
569 
570 /** gives solve statistics */
571 static
572 SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsAll)
573 {
574  SCIP_NLPIDATA* nlpidata;
575 
576  nlpidata = SCIPnlpiGetData(nlpi);
577  assert(nlpidata != NULL);
578  assert(nlpidata->nlpis != NULL);
579  assert(nlpidata->nlpis[problem->bestidx] != NULL);
580  assert(problem->nlpiproblems != NULL);
581  assert(problem->nlpiproblems[problem->bestidx] != NULL);
582 
583  /* collect statistics of the first solver */
584  SCIP_CALL( SCIPgetNlpiStatistics(scip, nlpidata->nlpis[problem->bestidx], problem->nlpiproblems[problem->bestidx],
585  statistics) );
586 
587  return SCIP_OKAY;
588 } /*lint !e715*/
589 
590 /*
591  * NLP solver interface specific interface methods
592  */
593 
594 /** create solver interface for the solver "All" and includes it into SCIP, if at least 2 NLPIs have already been included
595  *
596  * This method should be called after all other NLP solver interfaces have been included.
597  */
599  SCIP* scip /**< SCIP data structure */
600  )
601 {
602  SCIP_NLPIDATA* nlpidata;
603  int i;
604 
605  assert(scip != NULL);
606 
607  /* the number of NLPIs so far must be >= 2 */
608  if( SCIPgetNNlpis(scip) < 2 )
609  return SCIP_OKAY;
610 
611  /* create all solver interface data */
612  SCIP_CALL( SCIPallocClearBlockMemory(scip, &nlpidata) );
613 
614  nlpidata->nnlpis = SCIPgetNNlpis(scip);
615  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &nlpidata->nlpis, nlpidata->nnlpis) );
616 
617  /* copy nlpi pointers TODO should not need that */
618  for( i = 0; i < nlpidata->nnlpis; ++i )
619  nlpidata->nlpis[i] = SCIPgetNlpis(scip)[i];
620 
621  /* create solver interface */
624  nlpiCopyAll, nlpiFreeAll, NULL,
625  nlpiCreateProblemAll, nlpiFreeProblemAll, NULL,
626  nlpiAddVarsAll, nlpiAddConstraintsAll, nlpiSetObjectiveAll,
627  nlpiChgVarBoundsAll, nlpiChgConsSidesAll, nlpiDelVarSetAll, nlpiDelConstraintSetAll,
628  nlpiChgLinearCoefsAll, nlpiChgExprAll,
629  nlpiChgObjConstantAll, nlpiSetInitialGuessAll, nlpiSolveAll, nlpiGetSolstatAll, nlpiGetTermstatAll,
630  nlpiGetSolutionAll, nlpiGetStatisticsAll,
631  nlpidata) );
632 
633  return SCIP_OKAY;
634 }
#define NLPI_DESC
Definition: nlpi_all.c:41
#define SCIPfreeBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:110
enum SCIP_NlpTermStat SCIP_NLPTERMSTAT
Definition: type_nlpi.h:194
SCIP_NLPIDATA * SCIPnlpiGetData(SCIP_NLPI *nlpi)
Definition: nlpi.c:712
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:93
static SCIP_DECL_NLPICHGEXPR(nlpiChgExprAll)
Definition: nlpi_all.c:362
public methods for memory management
#define FALSE
Definition: def.h:96
SCIP_Real SCIPinfinity(SCIP *scip)
#define TRUE
Definition: def.h:95
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63
#define SCIPstatisticMessage
Definition: pub_message.h:123
SCIP_NLPIPROBLEM ** nlpiproblems
Definition: nlpi_all.c:56
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:108
const char * SCIPnlpiGetName(SCIP_NLPI *nlpi)
Definition: nlpi.c:722
SCIP_NLPI ** SCIPgetNlpis(SCIP *scip)
Definition: scip_nlpi.c:186
static SCIP_DECL_NLPICOPY(nlpiCopyAll)
Definition: nlpi_all.c:75
public methods for numerical tolerances
int SCIPgetNNlpis(SCIP *scip)
Definition: scip_nlpi.c:199
public methods for NLPI solver interfaces
static SCIP_DECL_NLPIGETSOLSTAT(nlpiGetSolstatAll)
Definition: nlpi_all.c:518
#define NLPI_PRIORITY
Definition: nlpi_all.c:42
struct SCIP_NlpiData SCIP_NLPIDATA
Definition: type_nlpi.h:52
enum SCIP_NlpSolStat SCIP_NLPSOLSTAT
Definition: type_nlpi.h:168
static SCIP_DECL_NLPISOLVE(nlpiSolveAll)
Definition: nlpi_all.c:426
SCIP_Real totaltime
Definition: type_nlpi.h:200
#define NULL
Definition: lpi_spx1.cpp:164
#define NLPI_NAME
Definition: nlpi_all.c:40
static SCIP_DECL_NLPISETOBJECTIVE(nlpiSetObjectiveAll)
Definition: nlpi_all.c:198
#define SCIP_CALL(x)
Definition: def.h:393
static SCIP_DECL_NLPIGETSOLUTION(nlpiGetSolutionAll)
Definition: nlpi_all.c:552
static SCIP_DECL_NLPICHGLINEARCOEFS(nlpiChgLinearCoefsAll)
Definition: nlpi_all.c:341
static SCIP_DECL_NLPIGETSTATISTICS(nlpiGetStatisticsAll)
Definition: nlpi_all.c:572
#define SCIP_Bool
Definition: def.h:93
static SCIP_DECL_NLPIDELCONSSET(nlpiDelConstraintSetAll)
Definition: nlpi_all.c:301
static SCIP_DECL_NLPICHGOBJCONSTANT(nlpiChgObjConstantAll)
Definition: nlpi_all.c:383
static SCIP_DECL_NLPIADDCONSTRAINTS(nlpiAddConstraintsAll)
Definition: nlpi_all.c:176
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:136
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
NLP interface that uses all available NLP interfaces.
static SCIP_DECL_NLPIGETTERMSTAT(nlpiGetTermstatAll)
Definition: nlpi_all.c:535
static SCIP_DECL_NLPICREATEPROBLEM(nlpiCreateProblemAll)
Definition: nlpi_all.c:100
static SCIP_DECL_NLPISETINITIALGUESS(nlpiSetInitialGuessAll)
Definition: nlpi_all.c:404
static SCIP_DECL_NLPIFREE(nlpiFreeAll)
Definition: nlpi_all.c:85
#define SCIP_Real
Definition: def.h:186
SCIP_RETCODE SCIPincludeNlpSolverAll(SCIP *scip)
Definition: nlpi_all.c:598
static SCIP_DECL_NLPICHGVARBOUNDS(nlpiChgVarBoundsAll)
Definition: nlpi_all.c:219
static SCIP_DECL_NLPIFREEPROBLEM(nlpiFreeProblemAll)
Definition: nlpi_all.c:128
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:111
SCIP_RETCODE SCIPincludeNlpi(SCIP *scip, const char *name, const char *description, int priority, SCIP_DECL_NLPICOPY((*nlpicopy)), SCIP_DECL_NLPIFREE((*nlpifree)), SCIP_DECL_NLPIGETSOLVERPOINTER((*nlpigetsolverpointer)), SCIP_DECL_NLPICREATEPROBLEM((*nlpicreateproblem)), SCIP_DECL_NLPIFREEPROBLEM((*nlpifreeproblem)), SCIP_DECL_NLPIGETPROBLEMPOINTER((*nlpigetproblempointer)), SCIP_DECL_NLPIADDVARS((*nlpiaddvars)), SCIP_DECL_NLPIADDCONSTRAINTS((*nlpiaddconstraints)), SCIP_DECL_NLPISETOBJECTIVE((*nlpisetobjective)), SCIP_DECL_NLPICHGVARBOUNDS((*nlpichgvarbounds)), SCIP_DECL_NLPICHGCONSSIDES((*nlpichgconssides)), SCIP_DECL_NLPIDELVARSET((*nlpidelvarset)), SCIP_DECL_NLPIDELCONSSET((*nlpidelconsset)), SCIP_DECL_NLPICHGLINEARCOEFS((*nlpichglinearcoefs)), SCIP_DECL_NLPICHGEXPR((*nlpichgexpr)), SCIP_DECL_NLPICHGOBJCONSTANT((*nlpichgobjconstant)), SCIP_DECL_NLPISETINITIALGUESS((*nlpisetinitialguess)), SCIP_DECL_NLPISOLVE((*nlpisolve)), SCIP_DECL_NLPIGETSOLSTAT((*nlpigetsolstat)), SCIP_DECL_NLPIGETTERMSTAT((*nlpigettermstat)), SCIP_DECL_NLPIGETSOLUTION((*nlpigetsolution)), SCIP_DECL_NLPIGETSTATISTICS((*nlpigetstatistics)), SCIP_NLPIDATA *nlpidata)
Definition: scip_nlpi.c:107
#define SCIPallocClearBlockMemory(scip, ptr)
Definition: scip_mem.h:91
static SCIP_DECL_NLPICHGCONSSIDES(nlpiChgConsSidesAll)
Definition: nlpi_all.c:240
static SCIP_DECL_NLPIDELVARSET(nlpiDelVarSetAll)
Definition: nlpi_all.c:261
static SCIP_DECL_NLPIADDVARS(nlpiAddVarsAll)
Definition: nlpi_all.c:154