Scippy

SCIP

Solving Constraint Integer Programs

nlhdlr_convex.c
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2 /* */
3 /* This file is part of the program and library */
4 /* SCIP --- Solving Constraint Integer Programs */
5 /* */
6 /* Copyright (C) 2002-2022 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not visit scip.zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file nlhdlr_convex.c
17  * @ingroup DEFPLUGINS_NLHDLR
18  * @brief nonlinear handlers for convex and concave expressions
19  * @author Benjamin Mueller
20  * @author Stefan Vigerske
21  *
22  * TODO convex: perturb reference point if separation fails due to too large numbers
23  */
24 
25 #include <string.h>
26 
27 #include "scip/nlhdlr_convex.h"
28 #include "scip/pub_nlhdlr.h"
29 #include "scip/scip_expr.h"
30 #include "scip/cons_nonlinear.h"
31 #include "scip/expr_var.h"
32 #include "scip/pub_misc_rowprep.h"
33 #include "scip/dbldblarith.h"
34 
35 /* fundamental nonlinear handler properties */
36 #define CONVEX_NLHDLR_NAME "convex"
37 #define CONVEX_NLHDLR_DESC "handler that identifies and estimates convex expressions"
38 #define CONVEX_NLHDLR_DETECTPRIORITY 50
39 #define CONVEX_NLHDLR_ENFOPRIORITY 50
40 
41 #define CONCAVE_NLHDLR_NAME "concave"
42 #define CONCAVE_NLHDLR_DESC "handler that identifies and estimates concave expressions"
43 #define CONCAVE_NLHDLR_DETECTPRIORITY 40
44 #define CONCAVE_NLHDLR_ENFOPRIORITY 40
45 
46 #define DEFAULT_DETECTSUM FALSE
47 #define DEFAULT_EXTENDEDFORM TRUE
48 #define DEFAULT_CVXQUADRATIC_CONVEX TRUE
49 #define DEFAULT_CVXQUADRATIC_CONCAVE FALSE
50 #define DEFAULT_CVXSIGNOMIAL TRUE
51 #define DEFAULT_CVXPRODCOMP TRUE
52 #define DEFAULT_HANDLETRIVIAL FALSE
53 
54 #define INITLPMAXVARVAL 1000.0 /**< maximal absolute value of variable for still generating a linearization cut at that point in initlp */
55 
56 /*lint -e440*/
57 /*lint -e441*/
58 /*lint -e666*/
59 /*lint -e777*/
60 
61 /*
62  * Data structures
63  */
64 
65 /** nonlinear handler expression data */
66 struct SCIP_NlhdlrExprData
67 {
68  SCIP_EXPR* nlexpr; /**< expression (copy) for which this nlhdlr estimates */
69  SCIP_HASHMAP* nlexpr2origexpr; /**< mapping of our copied expression to original expression */
70 
71  int nleafs; /**< number of distinct leafs of nlexpr, i.e., number of distinct (auxiliary) variables handled */
72  SCIP_EXPR** leafexprs; /**< distinct leaf expressions (excluding value-expressions), thus variables */
73 };
74 
75 /** nonlinear handler data */
76 struct SCIP_NlhdlrData
77 {
78  SCIP_Bool isnlhdlrconvex; /**< whether this data is used for the convex nlhdlr (TRUE) or the concave one (FALSE) */
79  SCIP_SOL* evalsol; /**< solution used for evaluating expression in a different point,
80  e.g., for facet computation of vertex-polyhedral function */
81 
82  /* parameters */
83  SCIP_Bool detectsum; /**< whether to run detection when the root of an expression is a non-quadratic sum */
84  SCIP_Bool extendedform; /**< whether to create extended formulations instead of looking for maximal possible subexpression */
85 
86  /* advanced parameters (maybe remove some day) */
87  SCIP_Bool cvxquadratic; /**< whether to use convexity check on quadratics */
88  SCIP_Bool cvxsignomial; /**< whether to use convexity check on signomials */
89  SCIP_Bool cvxprodcomp; /**< whether to use convexity check on product composition f(h)*h */
90  SCIP_Bool handletrivial; /**< whether to handle trivial expressions, i.e., those where all children are variables */
91 };
92 
93 /** data struct to be be passed on to vertexpoly-evalfunction (see SCIPcomputeFacetVertexPolyhedralNonlinear) */
94 typedef struct
95 {
100 
101 /** stack used in constructExpr to store expressions that need to be investigated ("to do list") */
102 typedef struct
103 {
104  SCIP_EXPR** stack; /**< stack elements */
105  int stacksize; /**< allocated space (in number of pointers) */
106  int stackpos; /**< position of top element of stack */
107 } EXPRSTACK;
108 
109 #define DECL_CURVCHECK(x) SCIP_RETCODE x( \
110  SCIP* scip, /**< SCIP data structure */ \
111  SCIP_EXPR* nlexpr, /**< nlhdlr-expr to check */ \
112  SCIP_Bool isrootexpr, /**< whether nlexpr is the root from where detection has been started */ \
113  EXPRSTACK* stack, /**< stack where to add generated leafs */ \
114  SCIP_HASHMAP* nlexpr2origexpr, /**< mapping from our expression copy to original expression */ \
115  SCIP_NLHDLRDATA* nlhdlrdata, /**< data of nlhdlr */ \
116  SCIP_HASHMAP* assumevarfixed, /**< hashmap containing variables that should be assumed to be fixed, or NULL */ \
117  SCIP_Bool* success /**< whether we found something */ \
118  )
119 
120 /*
121  * static methods
122  */
123 
124 /** create nlhdlr-expression
125  *
126  * does not create children, i.e., assumes that this will be a leaf
127  */
128 static
130  SCIP* scip, /**< SCIP data structure */
131  SCIP_HASHMAP* nlexpr2origexpr, /**< mapping from copied to original expression */
132  SCIP_EXPR** nlhdlrexpr, /**< buffer to store created expr */
133  SCIP_EXPR* origexpr, /**< original expression to be copied */
134  SCIP_EXPRCURV curv /**< curvature to achieve */
135  )
136 {
137  assert(scip != NULL);
138  assert(nlexpr2origexpr != NULL);
139  assert(nlhdlrexpr != NULL);
140  assert(origexpr != NULL);
141 
142  if( SCIPexprGetNChildren(origexpr) == 0 )
143  {
144  /* for leaves, do not copy */
145  *nlhdlrexpr = origexpr;
146  SCIPcaptureExpr(*nlhdlrexpr);
147  if( !SCIPhashmapExists(nlexpr2origexpr, (void*)*nlhdlrexpr) )
148  {
149  SCIP_CALL( SCIPhashmapInsert(nlexpr2origexpr, (void*)*nlhdlrexpr, (void*)origexpr) );
150  }
151  return SCIP_OKAY;
152  }
153 
154  /* create copy of expression, but without children */
155  SCIP_CALL( SCIPduplicateExprShallow(scip, origexpr, nlhdlrexpr, NULL, NULL) );
156  assert(*nlhdlrexpr != NULL); /* copies within the same SCIP must always work */
157 
158  /* store the curvature we want to get in the curvature flag of the copied expression
159  * it's a bit of a misuse, but once we are done with everything, this is actually correct
160  */
161  SCIPexprSetCurvature(*nlhdlrexpr, curv);
162 
163  /* remember which the original expression was */
164  SCIP_CALL( SCIPhashmapInsert(nlexpr2origexpr, (void*)*nlhdlrexpr, (void*)origexpr) );
165 
166  return SCIP_OKAY;
167 }
168 
169 /** expand nlhdlr-expression by adding children according to original expression */
170 static
172  SCIP* scip, /**< SCIP data structure */
173  SCIP_HASHMAP* nlexpr2origexpr, /**< mapping from copied to original expression */
174  SCIP_EXPR* nlhdlrexpr, /**< expression for which to create children */
175  SCIP_EXPRCURV* childrencurv /**< curvature required for children, or NULL if to set to UNKNOWN */
176  )
177 {
178  SCIP_EXPR* origexpr;
179  SCIP_EXPR* child;
180  int nchildren;
181  int i;
182 
183  assert(scip != NULL);
184  assert(nlhdlrexpr != NULL);
185  assert(SCIPexprGetNChildren(nlhdlrexpr) == 0);
186 
187  origexpr = (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlhdlrexpr);
188 
189  nchildren = SCIPexprGetNChildren(origexpr);
190  if( nchildren == 0 )
191  return SCIP_OKAY;
192 
193  for( i = 0; i < nchildren; ++i )
194  {
195  SCIP_CALL( nlhdlrExprCreate(scip, nlexpr2origexpr, &child, SCIPexprGetChildren(origexpr)[i],
196  childrencurv != NULL ? childrencurv[i] : SCIP_EXPRCURV_UNKNOWN) );
197  SCIP_CALL( SCIPappendExprChild(scip, nlhdlrexpr, child) );
198  /* append captures child, so we can release the capture from nlhdlrExprCreate */
199  SCIP_CALL( SCIPreleaseExpr(scip, &child) );
200  }
201 
202  assert(SCIPexprGetNChildren(nlhdlrexpr) == SCIPexprGetNChildren(origexpr));
203 
204  return SCIP_OKAY;
205 }
206 
207 /** evaluate expression at solution w.r.t. auxiliary variables */
208 static
209 SCIP_DECL_VERTEXPOLYFUN(nlhdlrExprEvalConcave)
210 {
211  VERTEXPOLYFUN_EVALDATA* evaldata = (VERTEXPOLYFUN_EVALDATA*)funcdata;
212  int i;
213 
214  assert(args != NULL);
215  assert(nargs == evaldata->nlhdlrexprdata->nleafs);
216  assert(evaldata != NULL);
217 
218 #ifdef SCIP_MORE_DEBUG
219  SCIPdebugMsg(evaldata->scip, "eval vertexpolyfun at\n");
220 #endif
221  for( i = 0; i < nargs; ++i )
222  {
223 #ifdef SCIP_MORE_DEBUG
224  SCIPdebugMsg(evaldata->scip, " <%s> = %g\n",
225  SCIPvarGetName(SCIPgetVarExprVar(evaldata->nlhdlrexprdata->leafexprs[i])), args[i]);
226 #endif
227  SCIP_CALL_ABORT( SCIPsetSolVal(evaldata->scip, evaldata->evalsol,
228  SCIPgetVarExprVar(evaldata->nlhdlrexprdata->leafexprs[i]), args[i]) );
229  }
230 
231  SCIP_CALL_ABORT( SCIPevalExpr(evaldata->scip, evaldata->nlhdlrexprdata->nlexpr, evaldata->evalsol, 0L) );
232 
233  return SCIPexprGetEvalValue(evaldata->nlhdlrexprdata->nlexpr);
234 }
235 
236 /** initialize expression stack */
237 static
239  SCIP* scip, /**< SCIP data structure */
240  EXPRSTACK* exprstack, /**< stack to initialize */
241  int initsize /**< initial size */
242  )
243 {
244  assert(scip != NULL);
245  assert(exprstack != NULL);
246  assert(initsize > 0);
247 
248  SCIP_CALL( SCIPallocBufferArray(scip, &exprstack->stack, initsize) );
249  exprstack->stacksize = initsize;
250  exprstack->stackpos = -1;
251 
252  return SCIP_OKAY;
253 }
254 
255 /** free expression stack */
256 static
258  SCIP* scip, /**< SCIP data structure */
259  EXPRSTACK* exprstack /**< free expression stack */
260  )
261 {
262  assert(scip != NULL);
263  assert(exprstack != NULL);
264 
265  SCIPfreeBufferArray(scip, &exprstack->stack);
266 }
267 
268 /** add expressions to expression stack */
269 static
271  SCIP* scip, /**< SCIP data structure */
272  EXPRSTACK* exprstack, /**< expression stack */
273  int nexprs, /**< number of expressions to push */
274  SCIP_EXPR** exprs /**< expressions to push */
275  )
276 {
277  assert(scip != NULL);
278  assert(exprstack != NULL);
279 
280  if( nexprs == 0 )
281  return SCIP_OKAY;
282 
283  assert(exprs != NULL);
284 
285  if( exprstack->stackpos+1 + nexprs > exprstack->stacksize ) /*lint !e644*/
286  {
287  exprstack->stacksize = SCIPcalcMemGrowSize(scip, exprstack->stackpos+1 + nexprs); /*lint !e644*/
288  SCIP_CALL( SCIPreallocBufferArray(scip, &exprstack->stack, exprstack->stacksize) );
289  }
290 
291  memcpy(exprstack->stack + (exprstack->stackpos+1), exprs, nexprs * sizeof(SCIP_EXPR*)); /*lint !e679*/ /*lint !e737*/
292  exprstack->stackpos += nexprs;
293 
294  return SCIP_OKAY;
295 }
296 
297 /** gives expression from top of expression stack and removes it from stack */
298 static
300  EXPRSTACK* exprstack /**< expression stack */
301  )
302 {
303  assert(exprstack != NULL);
304  assert(exprstack->stackpos >= 0);
305 
306  return exprstack->stack[exprstack->stackpos--];
307 }
308 
309 /** indicate whether expression stack is empty */
310 static
312  EXPRSTACK* exprstack /**< expression stack */
313  )
314 {
315  assert(exprstack != NULL);
316 
317  return exprstack->stackpos < 0;
318 }
319 
320 /** looks whether given expression is (proper) quadratic and has a given curvature
321  *
322  * If having a given curvature, currently require all arguments of quadratic to be linear.
323  * Hence, not using this for a simple square term, as curvCheckExprhdlr may provide a better condition on argument curvature then.
324  * Also we wouldn't do anything useful for a single bilinear term.
325  * Thus, run on sum's only.
326  */
327 static
328 DECL_CURVCHECK(curvCheckQuadratic)
329 { /*lint --e{715}*/
330  SCIP_EXPR* expr;
331  SCIP_EXPRCURV presentcurv;
332  SCIP_EXPRCURV wantedcurv;
333  SCIP_HASHSET* lonelysquares = NULL;
334  SCIP_Bool isquadratic;
335  int nbilinexprs;
336  int nquadexprs;
337  int i;
338 
339  assert(nlexpr != NULL);
340  assert(stack != NULL);
341  assert(nlexpr2origexpr != NULL);
342  assert(success != NULL);
343 
344  *success = FALSE;
345 
346  if( !nlhdlrdata->cvxquadratic )
347  return SCIP_OKAY;
348 
349  if( !SCIPisExprSum(scip, nlexpr) )
350  return SCIP_OKAY;
351 
352  wantedcurv = SCIPexprGetCurvature(nlexpr);
353  if( wantedcurv == SCIP_EXPRCURV_LINEAR )
354  return SCIP_OKAY;
355  assert(wantedcurv == SCIP_EXPRCURV_CONVEX || wantedcurv == SCIP_EXPRCURV_CONCAVE);
356 
357  expr = (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlexpr);
358  assert(expr != NULL);
359 
360  /* check whether quadratic */
361  SCIP_CALL( SCIPcheckExprQuadratic(scip, expr, &isquadratic) );
362 
363  /* if not quadratic, then give up here */
364  if( !isquadratic )
365  return SCIP_OKAY;
366 
367  SCIPexprGetQuadraticData(expr, NULL, NULL, NULL, NULL, &nquadexprs, &nbilinexprs, NULL, NULL);
368 
369  /* if only single square term (+linear), then give up here (let curvCheckExprhdlr handle this) */
370  if( nquadexprs <= 1 )
371  return SCIP_OKAY;
372 
373  /* if root expression is only sum of squares (+linear) and detectsum is disabled, then give up here, too */
374  if( isrootexpr && !nlhdlrdata->detectsum && nbilinexprs == 0 )
375  return SCIP_OKAY;
376 
377  /* get curvature of quadratic
378  * TODO as we know what curvature we want, we could first do some simple checks like computing xQx for a random x
379  */
380  SCIP_CALL( SCIPcomputeExprQuadraticCurvature(scip, expr, &presentcurv, assumevarfixed, FALSE) );
381 
382  /* if not having desired curvature, return */
383  if( presentcurv != wantedcurv )
384  return SCIP_OKAY;
385 
386  *success = TRUE;
387 
388  if( !nlhdlrdata->detectsum )
389  {
390  /* first step towards block-decomposition of quadratic term:
391  * collect all square-expressions (in original expr) which have no adjacent bilinear term
392  * we will treat these x^2 as linear, i.e., add an auxvar for them, so x^2 maybe linearized
393  * more efficiently (in particular if x is discrete)
394  */
395  SCIP_CALL( SCIPhashsetCreate(&lonelysquares, SCIPblkmem(scip), nquadexprs) );
396  for( i = 0; i < nquadexprs; ++i )
397  {
398  int nadjbilin;
399  SCIP_EXPR* sqrexpr;
400 
401  SCIPexprGetQuadraticQuadTerm(expr, i, NULL, NULL, NULL, &nadjbilin, NULL, &sqrexpr);
402  if( nadjbilin == 0 )
403  {
404  assert(sqrexpr != NULL);
405  SCIP_CALL( SCIPhashsetInsert(lonelysquares, SCIPblkmem(scip), (void*)sqrexpr) );
406  }
407  }
408  }
409 
410  /* add immediate children to nlexpr */
411  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, nlexpr, NULL) );
412  assert(SCIPexprGetNChildren(nlexpr) == SCIPexprGetNChildren(expr));
413 
414  /* put children that are not square or product on stack
415  * grow child for children that are square or product and put this child on stack
416  * require all children to be linear
417  */
418  for( i = 0; i < SCIPexprGetNChildren(nlexpr); ++i )
419  {
420  SCIP_EXPR* child;
422 
423  child = SCIPexprGetChildren(nlexpr)[i];
424  assert(child != NULL);
425 
426  assert(SCIPhashmapGetImage(nlexpr2origexpr, (void*)child) == SCIPexprGetChildren(expr)[i]);
427 
428  if( SCIPisExprPower(scip, child) && SCIPgetExponentExprPow(child) == 2.0 &&
429  (lonelysquares == NULL || !SCIPhashsetExists(lonelysquares, SCIPexprGetChildren(expr)[i])) )
430  {
431  /* square term that isn't lonely, i.e., orig-version of child is a square-expr and nadjbilin>0 */
432  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, child, curvlinear) );
433  assert(SCIPexprGetNChildren(child) == 1);
434  SCIP_CALL( exprstackPush(scip, stack, 1, SCIPexprGetChildren(child)) );
435  }
436  else if( SCIPisExprProduct(scip, child) && SCIPexprGetNChildren(SCIPexprGetChildren(expr)[i]) == 2 )
437  /* using original version of child here as NChildren(child)==0 atm */
438  {
439  /* bilinear term */
440  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, child, curvlinear) );
441  assert(SCIPexprGetNChildren(child) == 2);
442  SCIP_CALL( exprstackPush(scip, stack, 2, SCIPexprGetChildren(child)) );
443  }
444  else
445  {
446  /* linear term (or term to be considered as linear) */
447  SCIPexprSetCurvature(child, SCIP_EXPRCURV_LINEAR);
448  SCIP_CALL( exprstackPush(scip, stack, 1, &child) );
449  }
450  }
451 
452  if( lonelysquares != NULL )
453  SCIPhashsetFree(&lonelysquares, SCIPblkmem(scip));
454 
455  return SCIP_OKAY;
456 }
457 
458 /** looks whether top of given expression looks like a signomial that can have a given curvature
459  *
460  * e.g., sqrt(x)*sqrt(y) is convex if x,y >= 0 and x and y are convex
461  *
462  * unfortunately, doesn't work for tls, because i) it's originally sqrt(x*y), and ii) it is expanded into some sqrt(z*y+y);
463  * but works for cvxnonsep_nsig
464  */
465 static
466 DECL_CURVCHECK(curvCheckSignomial)
467 { /*lint --e{715}*/
468  SCIP_EXPR* expr;
469  SCIP_EXPR* child;
470  SCIP_Real* exponents;
471  SCIP_INTERVAL* bounds;
472  SCIP_EXPRCURV* curv;
473  int nfactors;
474  int i;
475 
476  assert(nlexpr != NULL);
477  assert(stack != NULL);
478  assert(nlexpr2origexpr != NULL);
479  assert(success != NULL);
480 
481  *success = FALSE;
482 
483  if( !nlhdlrdata->cvxsignomial )
484  return SCIP_OKAY;
485 
486  if( !SCIPisExprProduct(scip, nlexpr) )
487  return SCIP_OKAY;
488 
489  expr = (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlexpr);
490  assert(expr != NULL);
491 
492  nfactors = SCIPexprGetNChildren(expr);
493  if( nfactors <= 1 ) /* boooring */
494  return SCIP_OKAY;
495 
496  SCIP_CALL( SCIPallocBufferArray(scip, &exponents, nfactors) );
497  SCIP_CALL( SCIPallocBufferArray(scip, &bounds, nfactors) );
498  SCIP_CALL( SCIPallocBufferArray(scip, &curv, nfactors) );
499 
500  for( i = 0; i < nfactors; ++i )
501  {
502  child = SCIPexprGetChildren(expr)[i];
503  assert(child != NULL);
504 
505  if( !SCIPisExprPower(scip, child) )
506  {
507  exponents[i] = 1.0;
509  bounds[i] = SCIPexprGetActivity(child);
510  }
511  else
512  {
513  exponents[i] = SCIPgetExponentExprPow(child);
515  bounds[i] = SCIPexprGetActivity(SCIPexprGetChildren(child)[0]);
516  }
517  }
518 
520  nfactors, exponents, bounds, curv) )
521  goto TERMINATE;
522 
523  /* add immediate children to nlexpr
524  * some entries in curv actually apply to arguments of pow's, will correct this next
525  */
526  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, nlexpr, curv) );
527  assert(SCIPexprGetNChildren(nlexpr) == nfactors);
528 
529  /* put children that are not power on stack
530  * grow child for children that are power and put this child on stack
531  * if extendedform, then require children to be linear
532  * unless they are linear, an auxvar will be introduced for them and thus they will be handled as var here
533  */
534  for( i = 0; i < nfactors; ++i )
535  {
536  child = SCIPexprGetChildren(nlexpr)[i];
537  assert(child != NULL);
538 
539  if( SCIPisExprPower(scip, child) )
540  {
541  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, child, &curv[i]) );
542  assert(SCIPexprGetNChildren(child) == 1);
543  child = SCIPexprGetChildren(child)[0];
544  }
545  assert(SCIPexprGetNChildren(child) == 0);
546 
547  if( nlhdlrdata->extendedform )
548  {
550 #ifdef SCIP_DEBUG
551  SCIPinfoMessage(scip, NULL, "Extendedform: Require linearity for ");
552  SCIPprintExpr(scip, child, NULL);
553  SCIPinfoMessage(scip, NULL, "\n");
554 #endif
555  }
556 
557  SCIP_CALL( exprstackPush(scip, stack, 1, &child) );
558  }
559 
560  *success = TRUE;
561 
562 TERMINATE:
563  SCIPfreeBufferArray(scip, &curv);
564  SCIPfreeBufferArray(scip, &bounds);
565  SCIPfreeBufferArray(scip, &exponents);
566 
567  return SCIP_OKAY;
568 }
569 
570 /** looks for \f$f(c h(x)+d) h(x) \cdot \text{constant}\f$ and tries to conclude conditions on curvature
571  *
572  * Assume \f$h\f$ is univariate:
573  * - First derivative is \f$f'(c h + d) c h' h + f(c h + d) h'\f$.
574  * - Second derivative is \f{align}{&f''(c h + d) c h' c h' h + f'(c h + d) (c h'' h + c h' h') + f'(c h + d) c h' h' + f(c h + d) h'' \\
575  * =& f''(c h + d) c^2 h'^2 h + f'(c h + d) c h'' h + 2 f'(c h + d) c h'^2 + f(c h + d) h''.\f}
576  * Remove always positive factors leaves \f[f''(c h + d) h,\quad f'(c h + d) c h'' h,\quad f'(c h + d) c,\quad f(c h + d) h''.\f]
577  * For convexity we want all these terms to be nonnegative. For concavity we want all of them to be nonpositive.
578  * Note, that in each term either both \f$f'(c h + d)\f$ and \f$c\f$ occur, or none of them.
579  * - Thus, \f$f(c h(x) + d)h(x)\f$ is convex if \f$cf\f$ is monotonically increasing \f$(c f' \geq 0)\f$ and either
580  * - \f$f\f$ is convex \f$(f'' \geq 0)\f$ and \f$h\f$ is nonnegative \f$(h \geq 0)\f$ and \f$h\f$ is convex \f$(h'' \geq 0)\f$ and [\f$f\f$ is nonnegative \f$(f \geq 0)\f$ or \f$h\f$ is linear \f$(h''=0)\f$], or
581  * - \f$f\f$ is concave \f$(f'' \leq 0)\f$ and \f$h\f$ is nonpositive \f$(h \leq 0)\f$ and \f$h\f$ is concave \f$(h'' \leq 0)\f$ and [\f$f\f$ is nonpositive \f$(f \leq 0)\f$ or \f$h\f$ is linear \f$(h''=0)\f$].
582  * - Further, \f$f(c h(x) + d)h(x)\f$ is concave if \f$cf\f$ is monotonically decreasing \f$(c f' \leq 0)\f$ and either
583  * - f is convex \f$(f'' \geq 0)\f$ and \f$h\f$ is nonpositive \f$(h \leq 0)\f$ and \f$h\f$ is concave \f$(h'' \leq 0)\f$ and [\f$f\f$ is nonnegative \f$(f \geq 0)\f$ or \f$h\f$ is linear \f$(h''=0)\f$], or
584  * - f is concave \f$(f'' \leq 0)\f$ and \f$h\f$ is nonnegative \f$(h >= 0)\f$ and \f$h\f$ is convex \f$(h'' \geq 0)\f$ and [\f$f\f$ is nonpositive \f$(f \leq 0)\f$ or \f$h\f$ is linear \f$(h''=0)\f$].
585  *
586  * This should hold also for multivariate and linear \f$h\f$, as things are invariant under linear transformations.
587  * Similar to signomial, I'll assume that this will also hold for other multivariate \f$h\f$ (someone has a formal proof?).
588  */
589 static
590 DECL_CURVCHECK(curvCheckProductComposite)
591 { /*lint --e{715}*/
592  SCIP_EXPR* expr;
593  SCIP_EXPR* f;
594  SCIP_EXPR* h = NULL;
595  SCIP_Real c = 0.0;
596  SCIP_EXPR* ch = NULL; /* c * h */
597  SCIP_INTERVAL fbounds;
598  SCIP_INTERVAL hbounds;
599  SCIP_MONOTONE fmonotonicity;
600  SCIP_EXPRCURV desiredcurv;
601  SCIP_EXPRCURV hcurv;
602  SCIP_EXPRCURV dummy;
603  int fidx;
604 
605  assert(nlexpr != NULL);
606  assert(stack != NULL);
607  assert(nlexpr2origexpr != NULL);
608  assert(success != NULL);
609 
610  *success = FALSE;
611 
612  if( !nlhdlrdata->cvxprodcomp )
613  return SCIP_OKAY;
614 
615  if( !SCIPisExprProduct(scip, nlexpr) )
616  return SCIP_OKAY;
617 
618  expr = (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlexpr);
619  assert(expr != NULL);
620 
621  if( SCIPexprGetNChildren(expr) != 2 )
622  return SCIP_OKAY;
623 
624  /* check whether we have f(c * h(x)) * h(x) or h(x) * f(c * h(x)) */
625  for( fidx = 0; fidx <= 1; ++fidx )
626  {
627  f = SCIPexprGetChildren(expr)[fidx];
628 
629  if( SCIPexprGetNChildren(f) != 1 )
630  continue;
631 
632  ch = SCIPexprGetChildren(f)[0];
633  c = 1.0;
634  h = ch;
635 
636  /* check whether ch is of the form c*h(x), then switch h to child ch */
637  if( SCIPisExprSum(scip, ch) && SCIPexprGetNChildren(ch) == 1 )
638  {
639  c = SCIPgetCoefsExprSum(ch)[0];
640  h = SCIPexprGetChildren(ch)[0];
641  assert(c != 1.0 || SCIPgetConstantExprSum(ch) != 0.0); /* we could handle this, but it should have been simplified away */
642  }
643 
644 #ifndef NLHDLR_CONVEX_UNITTEST
645  /* can assume that duplicate subexpressions have been identified and comparing pointer is sufficient */
646  if( SCIPexprGetChildren(expr)[1-fidx] == h )
647 #else
648  /* called from unittest -> duplicate subexpressions were not identified -> compare more expensively */
649  if( SCIPcompareExpr(scip, SCIPexprGetChildren(expr)[1-fidx], h) == 0 )
650 #endif
651  break;
652  }
653  if( fidx == 2 )
654  return SCIP_OKAY;
655 
656 #ifdef SCIP_MORE_DEBUG
657  SCIPinfoMessage(scip, NULL, "f(c*h+d)*h with f = %s, c = %g, d = %g, h = ", SCIPexprhdlrGetName(SCIPexprGetHdlr(f)),
658  c, h != ch ? SCIPgetConstantExprSum(ch) : 0.0);
659  SCIPprintExpr(scip, h, NULL);
660  SCIPinfoMessage(scip, NULL, "\n");
661 #endif
662 
663  assert(c != 0.0);
664 
667  fbounds = SCIPexprGetActivity(f);
668  hbounds = SCIPexprGetActivity(h);
669 
670  /* if h has mixed sign, then cannot conclude anything */
671  if( hbounds.inf < 0.0 && hbounds.sup > 0.0 )
672  return SCIP_OKAY;
673 
674  SCIP_CALL( SCIPcallExprMonotonicity(scip, f, 0, &fmonotonicity) );
675 
676  /* if f is not monotone, then cannot conclude anything */
677  if( fmonotonicity == SCIP_MONOTONE_UNKNOWN )
678  return SCIP_OKAY;
679 
680  /* curvature we want to achieve (negate if product has negative coef) */
681  desiredcurv = SCIPexprcurvMultiply(SCIPgetCoefExprProduct(nlexpr), SCIPexprGetCurvature(nlexpr));
682 
683  /* now check the conditions as stated above */
684  if( desiredcurv == SCIP_EXPRCURV_CONVEX )
685  {
686  /* f(c h(x)+d)h(x) is convex if c*f is monotonically increasing (c f' >= 0) and either
687  * - f is convex (f'' >= 0) and h is nonnegative (h >= 0) and h is convex (h'' >= 0) and [f is nonnegative (f >= 0) or h is linear (h''=0)], or
688  * - f is concave (f'' <= 0) and h is nonpositive (h <= 0) and h is concave (h'' <= 0) and [f is nonpositive (f <= 0) or h is linear (h''=0)]
689  * as the curvature requirements on f are on f only and not the composition f(h), we can ignore the requirements returned by SCIPcallExprCurvature (last arg)
690  */
691  if( (c > 0.0 && fmonotonicity != SCIP_MONOTONE_INC) || (c < 0.0 && fmonotonicity != SCIP_MONOTONE_DEC) )
692  return SCIP_OKAY;
693 
694  /* check whether f can be convex (h>=0) or concave (h<=0), resp., and derive requirements for h */
695  if( hbounds.inf >= 0 )
696  {
697  SCIP_CALL( SCIPcallExprCurvature(scip, f, SCIP_EXPRCURV_CONVEX, success, &dummy) );
698 
699  /* now h also needs to be convex; and if f < 0, then h actually needs to be linear */
700  if( fbounds.inf < 0.0 )
701  hcurv = SCIP_EXPRCURV_LINEAR;
702  else
703  hcurv = SCIP_EXPRCURV_CONVEX;
704  }
705  else
706  {
707  SCIP_CALL( SCIPcallExprCurvature(scip, f, SCIP_EXPRCURV_CONCAVE, success, &dummy) );
708 
709  /* now h also needs to be concave; and if f > 0, then h actually needs to be linear */
710  if( fbounds.sup > 0.0 )
711  hcurv = SCIP_EXPRCURV_LINEAR;
712  else
713  hcurv = SCIP_EXPRCURV_CONCAVE;
714  }
715  }
716  else
717  {
718  /* f(c h(x)+d)*h(x) is concave if c*f is monotonically decreasing (c f' <= 0) and either
719  * - f is convex (f'' >= 0) and h is nonpositive (h <= 0) and h is concave (h'' <= 0) and [f is nonnegative (f >= 0) or h is linear (h''=0)], or
720  * - f is concave (f'' <= 0) and h is nonnegative (h >= 0) and h is convex (h'' >= 0) and [f is nonpositive (f <= 0) or h is linear (h''=0)]
721  * as the curvature requirements on f are on f only and not the composition f(h), we can ignore the requirements returned by SCIPcallExprCurvature (last arg)
722  */
723  if( (c > 0.0 && fmonotonicity != SCIP_MONOTONE_DEC) || (c < 0.0 && fmonotonicity != SCIP_MONOTONE_INC) )
724  return SCIP_OKAY;
725 
726  /* check whether f can be convex (h<=0) or concave (h>=0), resp., and derive requirements for h */
727  if( hbounds.sup <= 0 )
728  {
729  SCIP_CALL( SCIPcallExprCurvature(scip, f, SCIP_EXPRCURV_CONVEX, success, &dummy) );
730 
731  /* now h also needs to be concave; and if f < 0, then h actually needs to be linear */
732  if( fbounds.inf < 0.0 )
733  hcurv = SCIP_EXPRCURV_LINEAR;
734  else
735  hcurv = SCIP_EXPRCURV_CONCAVE;
736  }
737  else
738  {
739  SCIP_CALL( SCIPcallExprCurvature(scip, f, SCIP_EXPRCURV_CONCAVE, success, &dummy) );
740 
741  /* now h also needs to be convex; and if f > 0, then h actually needs to be linear */
742  if( fbounds.sup > 0.0 )
743  hcurv = SCIP_EXPRCURV_LINEAR;
744  else
745  hcurv = SCIP_EXPRCURV_CONVEX;
746  }
747  }
748 
749  if( !*success )
750  return SCIP_OKAY;
751 
752  /* add immediate children (f and ch) to nlexpr; we set required curvature for h further below */
753  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, nlexpr, NULL) );
754  assert(SCIPexprGetNChildren(nlexpr) == 2);
755 
756  /* copy of f (and h) should have same child position in nlexpr as f (and h) has on expr (resp) */
757  assert(SCIPhashmapGetImage(nlexpr2origexpr, (void*)SCIPexprGetChildren(nlexpr)[fidx]) == (void*)f);
758 #ifndef NLHDLR_CONVEX_UNITTEST
759  assert(SCIPhashmapGetImage(nlexpr2origexpr, (void*)SCIPexprGetChildren(nlexpr)[1-fidx]) == (void*)h);
760 #endif
761  /* push this h onto stack for further checking */
762  SCIP_CALL( exprstackPush(scip, stack, 1, &(SCIPexprGetChildren(nlexpr)[1-fidx])) );
763 
764  /* if we prefer extended formulations, then we always want h() to be linear */
765  if( nlhdlrdata->extendedform )
766  hcurv = SCIP_EXPRCURV_LINEAR;
767 
768  /* h-child of product should have curvature hcurv */
769  SCIPexprSetCurvature(SCIPexprGetChildren(nlexpr)[1-fidx], hcurv);
770 
771  if( h != ch )
772  {
773  /* add copy of ch as child to copy of f */
774  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, SCIPexprGetChildren(nlexpr)[fidx], NULL) );
775  assert(SCIPexprGetNChildren(SCIPexprGetChildren(nlexpr)[fidx]) == 1);
776  assert(SCIPhashmapGetImage(nlexpr2origexpr, (void*)SCIPexprGetChildren(SCIPexprGetChildren(nlexpr)[fidx])[0]) == (void*)ch);
777 
778  /* add copy of h (created above as child of product) as child in copy of ch */
780  SCIPexprGetChildren(SCIPexprGetChildren(nlexpr)[fidx])[0] /* copy of ch */,
781  SCIPexprGetChildren(nlexpr)[1-fidx] /* copy of h */) );
782  }
783  else
784  {
785  /* add copy of h (created above as child of product) as child in copy of f */
787  SCIPexprGetChildren(nlexpr)[fidx] /* copy of f */,
788  SCIPexprGetChildren(nlexpr)[1-fidx] /* copy of h */) );
789  }
790 
791  return SCIP_OKAY;
792 }
793 
794 /** use expression handlers curvature callback to check whether given curvature can be achieved */
795 static
796 DECL_CURVCHECK(curvCheckExprhdlr)
797 { /*lint --e{715}*/
798  SCIP_EXPR* origexpr;
799  int nchildren;
800  SCIP_EXPRCURV* childcurv;
801 
802  assert(nlexpr != NULL);
803  assert(stack != NULL);
804  assert(nlexpr2origexpr != NULL);
805  assert(success != NULL);
806 
807  origexpr = (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, nlexpr);
808  assert(origexpr != NULL);
809  nchildren = SCIPexprGetNChildren(origexpr);
810 
811  if( nchildren == 0 )
812  {
813  /* if originally no children, then should be var or value, which should have every curvature,
814  * so should always be success
815  */
816  SCIP_CALL( SCIPcallExprCurvature(scip, origexpr, SCIPexprGetCurvature(nlexpr), success, NULL) );
817  assert(*success);
818 
819  return SCIP_OKAY;
820  }
821 
822  /* ignore sums if > 1 children
823  * NOTE: this means that for something like 1+f(x), even if f is a trivial convex expression, we would handle 1+f(x)
824  * with this nlhdlr, instead of formulating this as 1+z and handling z=f(x) with the default nlhdlr, i.e., the exprhdlr
825  * today, I prefer handling this here, as it avoids introducing an extra auxiliary variable
826  */
827  if( isrootexpr && !nlhdlrdata->detectsum && SCIPisExprSum(scip, nlexpr) && nchildren > 1 )
828  return SCIP_OKAY;
829 
830  SCIP_CALL( SCIPallocBufferArray(scip, &childcurv, nchildren) );
831 
832  /* check whether and under which conditions origexpr can have desired curvature */
833  SCIP_CALL( SCIPcallExprCurvature(scip, origexpr, SCIPexprGetCurvature(nlexpr), success, childcurv) );
834 #ifdef SCIP_MORE_DEBUG
835  SCIPprintExpr(scip, origexpr, NULL);
836  SCIPinfoMessage(scip, NULL, " is %s? %d\n", SCIPexprcurvGetName(SCIPexprGetCurvature(nlexpr)), *success);
837 #endif
838  if( !*success )
839  goto TERMINATE;
840 
841  /* if origexpr can have curvature curv, then don't treat it as leaf, but include its children */
842  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, nlexpr, childcurv) );
843  assert(SCIPexprGetChildren(nlexpr) != NULL);
844  assert(SCIPexprGetNChildren(nlexpr) == nchildren);
845 
846  /* If we prefer extended formulations, then require all children to be linear.
847  * Unless they are, auxvars will be introduced and they will be handles as variables, which can be an
848  * advantage in the context of extended formulations.
849  */
850  if( nlhdlrdata->extendedform )
851  {
852  int i;
853  for( i = 0; i < nchildren; ++i )
855 #ifdef SCIP_DEBUG
856  SCIPinfoMessage(scip, NULL, "require linearity for children of ");
857  SCIPprintExpr(scip, origexpr, NULL);
858  SCIPinfoMessage(scip, NULL, "\n");
859 #endif
860  }
861 
862  /* add children expressions to to-do list (stack) */
863  SCIP_CALL( exprstackPush(scip, stack, nchildren, SCIPexprGetChildren(nlexpr)) );
864 
865 TERMINATE:
866  SCIPfreeBufferArray(scip, &childcurv);
867 
868  return SCIP_OKAY;
869 }
870 
871 /** curvature check and expression-growing methods
872  *
873  * some day this could be plugins added by users at runtime, but for now we have a fixed list here
874  * @note curvCheckExprhdlr should be last
875  */
876 static DECL_CURVCHECK((*CURVCHECKS[])) = { curvCheckProductComposite, curvCheckSignomial, curvCheckQuadratic, curvCheckExprhdlr };
877 /** number of curvcheck methods */
878 static const int NCURVCHECKS = sizeof(CURVCHECKS) / sizeof(void*);
879 
880 /** checks whether expression is a sum with more than one child and each child being a variable or going to be a variable if `expr` is a nlhdlr-specific copy
881  *
882  * Within constructExpr(), we can have an expression of any type which is a copy of an original expression,
883  * but without children. At the end of constructExpr() (after the loop with the stack), these expressions
884  * will remain as leafs and will eventually be turned into variables in collectLeafs(). Thus, we treat
885  * every child that has no children as if it were a variable. Theoretically, there is still the possibility
886  * that it could be a constant (value-expression), but simplify should have removed these.
887  */
888 static
890  SCIP* scip, /**< SCIP data structure */
891  SCIP_EXPR* expr /**< expression to check */
892  )
893 {
894  int nchildren;
895  int c;
896 
897  assert(expr != NULL);
898 
899  if( !SCIPisExprSum(scip, expr) )
900  return FALSE;
901 
902  nchildren = SCIPexprGetNChildren(expr);
903  if( nchildren <= 1 )
904  return FALSE;
905 
906  for( c = 0; c < nchildren; ++c )
907  /*if( !SCIPisExprVar(scip, SCIPexprGetChildren(expr)[c]) ) */
908  if( SCIPexprGetNChildren(SCIPexprGetChildren(expr)[c]) > 0 )
909  return FALSE;
910 
911  return TRUE;
912 }
913 
914 /** constructs a subexpression (as nlhdlr-expression) of maximal size that has a given curvature
915  *
916  * If the curvature cannot be achieved for an expression in the original expression graph,
917  * then this expression becomes a leaf in the nlhdlr-expression.
918  *
919  * Sets `*rootnlexpr` to NULL if failed.
920  */
921 static
923  SCIP* scip, /**< SCIP data structure */
924  SCIP_NLHDLRDATA* nlhdlrdata, /**< nonlinear handler data */
925  SCIP_EXPR** rootnlexpr, /**< buffer to store created expression */
926  SCIP_HASHMAP* nlexpr2origexpr, /**< mapping from our expression copy to original expression */
927  int* nleafs, /**< number of leafs in constructed expression */
928  SCIP_EXPR* rootexpr, /**< expression */
929  SCIP_EXPRCURV curv, /**< curvature to achieve */
930  SCIP_HASHMAP* assumevarfixed, /**< hashmap containing variables that should be assumed to be fixed, or NULL */
931  SCIP_Bool assumecurvature, /**< whether to assume that desired curvature is given (skips curvature checks) */
932  SCIP_Bool* curvsuccess /**< pointer to store whether the curvature could be achieved
933  w.r.t. the original variables (might be NULL) */
934  )
935 {
936  SCIP_EXPR* nlexpr;
937  EXPRSTACK stack; /* to do list: expressions where to check whether they can have the desired curvature when taking their children into account */
938  int oldstackpos;
939  SCIP_Bool isrootexpr = TRUE;
940 
941  assert(scip != NULL);
942  assert(nlhdlrdata != NULL);
943  assert(rootnlexpr != NULL);
944  assert(nlexpr2origexpr != NULL);
945  assert(nleafs != NULL);
946  assert(rootexpr != NULL);
947  assert(curv == SCIP_EXPRCURV_CONVEX || curv == SCIP_EXPRCURV_CONCAVE);
948 
949  /* create root expression */
950  SCIP_CALL( nlhdlrExprCreate(scip, nlexpr2origexpr, rootnlexpr, rootexpr, curv) );
951 
952  *nleafs = 0;
953  if( curvsuccess != NULL )
954  *curvsuccess = TRUE;
955 
956  SCIP_CALL( exprstackInit(scip, &stack, 20) );
957  SCIP_CALL( exprstackPush(scip, &stack, 1, rootnlexpr) );
958  while( !exprstackIsEmpty(&stack) )
959  {
960  /* take expression from stack */
961  nlexpr = exprstackPop(&stack);
962  assert(nlexpr != NULL);
963  assert(SCIPexprGetNChildren(nlexpr) == 0);
964 
965  oldstackpos = stack.stackpos;
966  if( nlhdlrdata->isnlhdlrconvex && !SCIPexprhdlrHasBwdiff(SCIPexprGetHdlr(nlexpr)) )
967  {
968  /* if bwdiff is not implemented, then we could not generate cuts in the convex nlhdlr, so "stop" (treat nlexpr as variable) */
969  }
970  else if( !nlhdlrdata->isnlhdlrconvex && exprIsMultivarLinear(scip, (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlexpr)) )
971  {
972  /* if we are in the concave handler, we would like to treat linear multivariate subexpressions by a new auxvar always,
973  * e.g., handle log(x+y) as log(z), z=x+y, because the estimation problem will be smaller then without making the estimator worse
974  * (cons_nonlinear does this, too)
975  * this check takes care of this when x and y are original variables
976  * however, it isn't unlikely that we will have sums that become linear after we add auxvars for some children
977  * this will be handled in a postprocessing below
978  * for now, the check is performed on the original expression since there is not enough information in nlexpr yet
979  */
980 #ifdef SCIP_MORE_DEBUG
981  SCIPprintExpr(scip, SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlexpr), NULL);
982  SCIPinfoMessage(scip, NULL, "... is a multivariate linear sum that we'll treat as auxvar\n");
983 #endif
984  }
985  else if( SCIPexprGetCurvature(nlexpr) != SCIP_EXPRCURV_UNKNOWN && !assumecurvature )
986  {
987  /* if we are here, either convexity or concavity is required; try to check for this curvature */
988  SCIP_Bool success;
989  int method;
990 
991  /* try through curvature check methods until one succeeds */
992  for( method = 0; method < NCURVCHECKS; ++method )
993  {
994  SCIP_CALL( CURVCHECKS[method](scip, nlexpr, isrootexpr, &stack, nlexpr2origexpr, nlhdlrdata, assumevarfixed, &success) );
995  if( success )
996  break;
997  }
998  }
999  else
1000  {
1001  /* if we don't care about curvature in this subtree anymore (very unlikely),
1002  * or we are told to assume that the desired curvature is present (assumecurvature==TRUE),
1003  * then only continue iterating this subtree to assemble leaf expressions
1004  */
1005  SCIP_CALL( nlhdlrExprGrowChildren(scip, nlexpr2origexpr, nlexpr, NULL) );
1006 
1007  /* add children expressions, if any, to to-do list (stack) */
1008  SCIP_CALL( exprstackPush(scip, &stack, SCIPexprGetNChildren(nlexpr), SCIPexprGetChildren(nlexpr)) );
1009  }
1010  assert(stack.stackpos >= oldstackpos); /* none of the methods above should have removed something from the stack */
1011 
1012  isrootexpr = FALSE;
1013 
1014  /* if nothing was added, then none of the successors of nlexpr were added to the stack
1015  * this is either because nlexpr was already a variable or value expressions, thus a leaf,
1016  * or because the desired curvature could not be achieved, so it will be handled as variables, thus a leaf
1017  */
1018  if( stack.stackpos == oldstackpos )
1019  {
1020  ++*nleafs;
1021 
1022  /* check whether the new leaf is not an original variable (or constant) */
1023  if( curvsuccess != NULL && !SCIPisExprVar(scip, nlexpr) && !SCIPisExprValue(scip, nlexpr) )
1024  *curvsuccess = FALSE;
1025  }
1026  }
1027 
1028  exprstackFree(scip, &stack);
1029 
1030  if( !nlhdlrdata->isnlhdlrconvex && *rootnlexpr != NULL )
1031  {
1032  /* remove multivariate linear subexpressions, that is, change some f(z1+z2) into f(z3) (z3=z1+z2 will be done by nlhdlr_default)
1033  * this handles the case that was not covered by the above check, which could recognize f(x+y) for x, y original variables
1034  */
1035  SCIP_EXPRITER* it;
1036 
1037  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
1038  SCIP_CALL( SCIPexpriterInit(it, *rootnlexpr, SCIP_EXPRITER_DFS, FALSE) );
1040 
1041  while( !SCIPexpriterIsEnd(it) )
1042  {
1043  SCIP_EXPR* child;
1044 
1045  child = SCIPexpriterGetChildExprDFS(it);
1046  assert(child != NULL);
1047 
1048  /* We want to change some f(x+y+z) into just f(), where f is the expression the iterator points to
1049  * and x+y+z is child. A child of a child, e.g., z, may not be a variable yet (these are added in collectLeafs later),
1050  * but an expression of some nonlinear type without children.
1051  */
1052  if( exprIsMultivarLinear(scip, child) )
1053  {
1054  /* turn child (x+y+z) into a sum without children
1055  * collectLeafs() should then replace this by an auxvar
1056  */
1057 #ifdef SCIP_MORE_DEBUG
1058  SCIPprintExpr(scip, child, NULL);
1059  SCIPinfoMessage(scip, NULL, "... is a multivariate linear sum that we'll treat as auxvar instead (postprocess)\n");
1060 #endif
1061 
1062  /* TODO remove children from nlexpr2origexpr ?
1063  * should also do this if they are not used somewhere else; we could check nuses for this
1064  * however, it shouldn't matter to have some stray entries in the hashmap either
1065  */
1066  SCIP_CALL( SCIPremoveExprChildren(scip, child) );
1067  assert(SCIPexprGetNChildren(child) == 0);
1068 
1069  (void) SCIPexpriterSkipDFS(it);
1070  }
1071  else
1072  {
1073  (void) SCIPexpriterGetNext(it);
1074  }
1075  }
1076 
1077  SCIPfreeExpriter(&it);
1078  }
1079 
1080  if( *rootnlexpr != NULL )
1081  {
1082  SCIP_Bool istrivial = TRUE;
1083 
1084  /* if handletrivial is enabled, then only require that rootnlexpr itself has required curvature (so has children; see below) and
1085  * that we are not a trivial sum (because the previous implementation of this nlhdlr didn't allow this, either)
1086  */
1087  if( !nlhdlrdata->handletrivial || SCIPisExprSum(scip, *rootnlexpr) )
1088  {
1089  /* if all children do not have children, i.e., are variables, or will be replaced by auxvars, then free
1090  * also if rootnlexpr has no children, then free
1091  */
1092  int i;
1093  for( i = 0; i < SCIPexprGetNChildren(*rootnlexpr); ++i )
1094  {
1095  if( SCIPexprGetNChildren(SCIPexprGetChildren(*rootnlexpr)[i]) > 0 )
1096  {
1097  istrivial = FALSE;
1098  break;
1099  }
1100  }
1101  }
1102  else if( SCIPexprGetNChildren(*rootnlexpr) > 0 ) /* if handletrivial, then just require children */
1103  istrivial = FALSE;
1104 
1105  if( istrivial )
1106  {
1107  SCIP_CALL( SCIPreleaseExpr(scip, rootnlexpr) );
1108  }
1109  }
1110 
1111  return SCIP_OKAY;
1112 }
1113 
1114 /** collects (non-value) leaf expressions and ensure that they correspond to a variable (original or auxiliary)
1115  *
1116  * For children where we could not achieve the desired curvature, get the auxvar and replace the child by a
1117  * var-expression that points to this auxvar.
1118  * Collect all leaf expressions (if not a value-expression) and index them.
1119  */
1120 static
1122  SCIP* scip, /**< SCIP data structure */
1123  SCIP_NLHDLREXPRDATA* nlhdlrexprdata /**< nlhdlr expression data */
1124  )
1125 {
1126  SCIP_EXPRITER* it;
1127  SCIP_EXPR* nlexpr;
1128  SCIP_HASHMAP* leaf2index;
1129  int i;
1130 
1131  assert(nlhdlrexprdata != NULL);
1132  assert(nlhdlrexprdata->nlexpr != NULL);
1133  assert(nlhdlrexprdata->nlexpr2origexpr != NULL);
1134  /* nleafs should be the upper bound on the number of variables given by constructExpr
1135  * leafexprs should be NULL, as this is what we want to setup here
1136  */
1137  assert(nlhdlrexprdata->nleafs > 0);
1138  assert(nlhdlrexprdata->leafexprs == NULL);
1139 
1140  /* collect all auxvars and collect all variables */
1141  SCIP_CALL( SCIPhashmapCreate(&leaf2index, SCIPblkmem(scip), nlhdlrexprdata->nleafs) );
1142  nlhdlrexprdata->nleafs = 0; /* we start a new count, this time skipping value-expressions */
1143 
1144  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
1145  SCIP_CALL( SCIPexpriterInit(it, nlhdlrexprdata->nlexpr, SCIP_EXPRITER_DFS, FALSE) );
1147 
1148  for( nlexpr = SCIPexpriterGetCurrent(it); !SCIPexpriterIsEnd(it); nlexpr = SCIPexpriterGetNext(it) )
1149  {
1150  SCIP_EXPR* child;
1151  SCIP_EXPR* origexpr;
1152 
1153  assert(nlexpr != NULL);
1154 
1155  child = SCIPexpriterGetChildExprDFS(it);
1156 
1157  /* if the to-be-visited child has children, then it doesn't need to be replaced by a new expression (representing the auxvar) */
1158  if( SCIPexprGetNChildren(child) > 0 )
1159  continue;
1160 
1161  origexpr = (SCIP_EXPR*)SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, (void*)child);
1162  assert(origexpr != NULL);
1163 
1164  if( SCIPexprGetNChildren(origexpr) > 0 )
1165  {
1166  SCIP_EXPR* newchild;
1167  int childidx;
1168  SCIP_VAR* var;
1169 
1170  /* having a child that had children in original but not in copy means that we could not achieve the desired curvature
1171  * thus, replace by a new child that points to the auxvar of the original expression
1172  * we registered in createNlhdlrExprData that we need an auxvar, so it should exist now
1173  */
1174  var = SCIPgetExprAuxVarNonlinear(origexpr);
1175  assert(var != NULL);
1176 
1177  SCIP_CALL( SCIPcreateExprVar(scip, &newchild, var, NULL, NULL) ); /* this captures newchild once */
1178 
1179  childidx = SCIPexpriterGetChildIdxDFS(it);
1180  SCIP_CALL( SCIPreplaceExprChild(scip, nlexpr, childidx, newchild) ); /* this captures newchild again */
1181 
1182  /* do not remove child->origexpr from hashmap, as child may appear again due to common subexprs
1183  * (created by curvCheckProductComposite, for example)
1184  * if it doesn't reappear, though, but the memory address is reused, we need to make sure it
1185  * points to the right origexpr
1186  */
1187  /* SCIP_CALL( SCIPhashmapRemove(nlexpr2origexpr, (void*)child) ); */
1188  SCIP_CALL( SCIPhashmapSetImage(nlhdlrexprdata->nlexpr2origexpr, (void*)newchild, (void*)origexpr) );
1189 
1190  if( !SCIPhashmapExists(leaf2index, (void*)newchild) )
1191  {
1192  /* new leaf -> new index and remember in hashmap */
1193  SCIP_CALL( SCIPhashmapInsertInt(leaf2index, (void*)newchild, nlhdlrexprdata->nleafs++) );
1194  }
1195 
1196  child = newchild;
1197  SCIP_CALL( SCIPreleaseExpr(scip, &newchild) ); /* because it was captured by both create and replace */
1198  }
1199  else if( SCIPisExprVar(scip, child) )
1200  {
1201  /* if variable, then add to hashmap, if not already there */
1202  if( !SCIPhashmapExists(leaf2index, (void*)child) )
1203  {
1204  SCIP_CALL( SCIPhashmapInsertInt(leaf2index, (void*)child, nlhdlrexprdata->nleafs++) );
1205  }
1206  }
1207  /* else: it's probably a value-expression, nothing to do */
1208 
1209  /* update integrality flag for future leaf expressions: convex nlhdlr may use this information */
1210  SCIP_CALL( SCIPcomputeExprIntegrality(scip, child) );
1211  }
1212  assert(nlhdlrexprdata->nleafs > 0);
1213 
1214  SCIPfreeExpriter(&it);
1215 
1216  /* assemble auxvars array */
1217  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(nlhdlrexprdata->leafexprs), nlhdlrexprdata->nleafs) );
1218  for( i = 0; i < SCIPhashmapGetNEntries(leaf2index); ++i )
1219  {
1220  SCIP_HASHMAPENTRY* entry;
1221  SCIP_EXPR* leaf;
1222  int idx;
1223 
1224  entry = SCIPhashmapGetEntry(leaf2index, i);
1225  if( entry == NULL )
1226  continue;
1227 
1228  leaf = (SCIP_EXPR*) SCIPhashmapEntryGetOrigin(entry);
1229  assert(leaf != NULL);
1230  assert(SCIPisExprVar(scip, leaf));
1231 
1232  idx = SCIPhashmapEntryGetImageInt(entry);
1233  assert(idx >= 0);
1234  assert(idx < nlhdlrexprdata->nleafs);
1235 
1236  nlhdlrexprdata->leafexprs[idx] = leaf;
1237 
1238  SCIPdebugMsg(scip, "leaf %d: <%s>\n", idx, SCIPvarGetName(SCIPgetVarExprVar(leaf)));
1239  }
1240 
1241  SCIPhashmapFree(&leaf2index);
1242 
1243  return SCIP_OKAY;
1244 }
1245 
1246 /** creates nonlinear handler expression data structure and registers expr usage */
1247 static
1249  SCIP* scip, /**< SCIP data structure */
1250  SCIP_NLHDLRDATA* nlhdlrdata, /**< nlhdlr data */
1251  SCIP_NLHDLREXPRDATA** nlhdlrexprdata, /**< pointer to store nlhdlr expression data */
1252  SCIP_EXPR* expr, /**< original expression */
1253  SCIP_EXPR* nlexpr, /**< our copy of expression */
1254  SCIP_HASHMAP* nlexpr2origexpr, /**< mapping of expression copy to original */
1255  int nleafs, /**< number of leafs as counted by constructExpr */
1256  SCIP_NLHDLR_METHOD participating /**< the enfo methods in which we plan to participate */
1257  )
1258 {
1259  SCIP_EXPRITER* it;
1260  SCIP_Bool usingaux;
1261 
1262  assert(scip != NULL);
1263  assert(expr != NULL);
1264  assert(nlhdlrexprdata != NULL);
1265  assert(*nlhdlrexprdata == NULL);
1266  assert(nlexpr != NULL);
1267  assert(nlexpr2origexpr != NULL);
1268 
1269  assert(SCIPexprGetNChildren(nlexpr) > 0);
1270  assert(SCIPexprGetChildren(nlexpr) != NULL);
1271 
1272  SCIP_CALL( SCIPallocClearBlockMemory(scip, nlhdlrexprdata) );
1273  (*nlhdlrexprdata)->nlexpr = nlexpr;
1274  (*nlhdlrexprdata)->nlexpr2origexpr = nlexpr2origexpr;
1275  (*nlhdlrexprdata)->nleafs = nleafs;
1276 
1277  usingaux = FALSE;
1278 
1279  SCIP_CALL( SCIPcreateExpriter(scip, &it) );
1282 
1283  for( ; !SCIPexpriterIsEnd(it); (void) SCIPexpriterGetNext(it) )
1284  {
1285  SCIP_EXPR* child;
1286  SCIP_EXPR* origexpr;
1287 
1288  /* check whether to-be-visited child needs to be replaced by a new expression (representing the auxvar)
1289  * if child has children, then that is not the case
1290  * if child has no children, but also corresponding origexpr has no chilren, then this is also not the case
1291  */
1292  child = SCIPexpriterGetChildExprDFS(it);
1293  if( SCIPexprGetNChildren(child) > 0 )
1294  continue;
1295 
1296  origexpr = (SCIP_EXPR*)SCIPhashmapGetImage(nlexpr2origexpr, (void*)child);
1297  assert(origexpr != NULL);
1298 
1299  /* if child had children in original but not in copy means that we could not achieve the desired curvature
1300  * thus, we will later replace by a new child that points to the auxvar of the original expression
1301  * as we do not have the auxvar now, we will only register that we will need the auxvar later (if origexpr isn't a variable or constant)
1302  * if we are working for the concave nlhdlr, then we also indicate interest on the exprs activity for estimate (distinguish below or above)
1303  */
1304  SCIP_CALL( SCIPregisterExprUsageNonlinear(scip, origexpr,
1305  SCIPexprGetNChildren(origexpr) > 0, FALSE,
1306  !nlhdlrdata->isnlhdlrconvex && (participating & SCIP_NLHDLR_METHOD_SEPABELOW),
1307  !nlhdlrdata->isnlhdlrconvex && (participating & SCIP_NLHDLR_METHOD_SEPAABOVE)) );
1308 
1309  /* remember that we use an auxvar */
1310  if( SCIPexprGetNChildren(origexpr) > 0 )
1311  usingaux = TRUE;
1312  }
1313 
1314  SCIPfreeExpriter(&it);
1315 
1316 #ifdef SCIP_DEBUG
1317  SCIPprintExpr(scip, nlexpr, NULL);
1318  SCIPinfoMessage(scip, NULL, " (%p) is handled as %s\n", SCIPhashmapGetImage(nlexpr2origexpr, (void*)nlexpr),
1320 #endif
1321 
1322  /* If we don't work on the extended formulation, then set curvature also in original expression
1323  * (in case someone wants to pick this up; this might be removed again).
1324  * This doesn't ensure that every convex or concave original expression is actually marked here.
1325  * Not only because our tests are incomprehensive, but also because we may not detect on sums,
1326  * prefer extended formulations (in nlhdlr_convex), or introduce auxvars for linear subexpressions
1327  * on purpose (in nlhdlr_concave).
1328  */
1329  if( !usingaux )
1331 
1332  return SCIP_OKAY;
1333 }
1334 
1335 /** adds an estimator for a vertex-polyhedral (e.g., concave) function to a given rowprep
1336  *
1337  * Calls \ref SCIPcomputeFacetVertexPolyhedralNonlinear() for given function and
1338  * box set to local bounds of auxiliary variables.
1339  */
1340 static
1342  SCIP* scip, /**< SCIP data structure */
1343  SCIP_CONSHDLR* conshdlr, /**< nonlinear constraint handler */
1344  SCIP_NLHDLR* nlhdlr, /**< nonlinear handler */
1345  SCIP_NLHDLREXPRDATA* nlhdlrexprdata, /**< nonlinear handler expression data */
1346  SCIP_SOL* sol, /**< solution to use, unless usemidpoint is TRUE */
1347  SCIP_Bool usemidpoint, /**< whether to use the midpoint of the domain instead of sol */
1348  SCIP_Bool overestimate, /**< whether over- or underestimating */
1349  SCIP_Real targetvalue, /**< a target value to achieve; if not reachable, then can give up early */
1350  SCIP_ROWPREP* rowprep, /**< rowprep where to store estimator */
1351  SCIP_Bool* success /**< buffer to store whether successful */
1352  )
1353 {
1354  SCIP_NLHDLRDATA* nlhdlrdata;
1355  VERTEXPOLYFUN_EVALDATA evaldata;
1356  SCIP_Real* xstar;
1357  SCIP_Real* box;
1358  SCIP_Real facetconstant;
1359  SCIP_VAR* var;
1360  int i;
1361  SCIP_Bool allfixed;
1362 
1363  assert(scip != NULL);
1364  assert(nlhdlr != NULL);
1365  assert(nlhdlrexprdata != NULL);
1366  assert(rowprep != NULL);
1367  assert(success != NULL);
1368 
1369  *success = FALSE;
1370 
1371  /* caller is responsible to have checked whether we can estimate, i.e., expression curvature and overestimate flag match */
1372  assert( overestimate || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONCAVE); /* if underestimate, then must be concave */
1373  assert(!overestimate || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONVEX); /* if overestimate, then must be convex */
1374 
1375 #ifdef SCIP_DEBUG
1376  SCIPinfoMessage(scip, NULL, "%sestimate expression ", overestimate ? "over" : "under");
1377  SCIPprintExpr(scip, nlhdlrexprdata->nlexpr, NULL);
1378  SCIPinfoMessage(scip, NULL, " at point\n");
1379  for( i = 0; i < nlhdlrexprdata->nleafs; ++i )
1380  {
1381  var = SCIPgetVarExprVar(nlhdlrexprdata->leafexprs[i]);
1382  assert(var != NULL);
1383 
1384  SCIPinfoMessage(scip, NULL, " <%s> = %g [%g,%g]\n", SCIPvarGetName(var),
1385  usemidpoint ? 0.5 * (SCIPvarGetLbLocal(var) + SCIPvarGetUbLocal(var)) : SCIPgetSolVal(scip, sol, var),
1387  }
1388 #endif
1389 
1390  nlhdlrdata = SCIPnlhdlrGetData(nlhdlr);
1391  assert(nlhdlrdata != NULL);
1392 
1393  if( nlhdlrdata->evalsol == NULL )
1394  {
1395  SCIP_CALL( SCIPcreateSol(scip, &nlhdlrdata->evalsol, NULL) );
1396  }
1397 
1398  evaldata.nlhdlrexprdata = nlhdlrexprdata;
1399  evaldata.evalsol = nlhdlrdata->evalsol;
1400  evaldata.scip = scip;
1401 
1402  SCIP_CALL( SCIPallocBufferArray(scip, &xstar, nlhdlrexprdata->nleafs) );
1403  SCIP_CALL( SCIPallocBufferArray(scip, &box, 2*nlhdlrexprdata->nleafs) );
1404 
1405  allfixed = TRUE;
1406  for( i = 0; i < nlhdlrexprdata->nleafs; ++i )
1407  {
1408  var = SCIPgetVarExprVar(nlhdlrexprdata->leafexprs[i]);
1409  assert(var != NULL);
1410 
1411  box[2*i] = SCIPvarGetLbLocal(var);
1412  if( SCIPisInfinity(scip, -box[2*i]) )
1413  {
1414  SCIPdebugMsg(scip, "lower bound at -infinity, no estimate possible\n");
1415  goto TERMINATE;
1416  }
1417 
1418  box[2*i+1] = SCIPvarGetUbLocal(var);
1419  if( SCIPisInfinity(scip, box[2*i+1]) )
1420  {
1421  SCIPdebugMsg(scip, "upper bound at +infinity, no estimate possible\n");
1422  goto TERMINATE;
1423  }
1424 
1425  if( !SCIPisRelEQ(scip, box[2*i], box[2*i+1]) )
1426  allfixed = FALSE;
1427 
1428  if( usemidpoint )
1429  xstar[i] = 0.5 * (box[2*i] + box[2*i+1]);
1430  else
1431  xstar[i] = SCIPgetSolVal(scip, sol, var);
1432  assert(xstar[i] != SCIP_INVALID);
1433  }
1434 
1435  if( allfixed )
1436  {
1437  /* SCIPcomputeFacetVertexPolyhedralNonlinear prints a warning and does not succeed if all is fixed */
1438  SCIPdebugMsg(scip, "all variables fixed, skip estimate\n");
1439  goto TERMINATE;
1440  }
1441 
1442  SCIP_CALL( SCIPensureRowprepSize(scip, rowprep, nlhdlrexprdata->nleafs + 1) );
1443 
1444  SCIP_CALL( SCIPcomputeFacetVertexPolyhedralNonlinear(scip, conshdlr, overestimate, nlhdlrExprEvalConcave, (void*)&evaldata,
1445  xstar, box, nlhdlrexprdata->nleafs, targetvalue, success, SCIProwprepGetCoefs(rowprep), &facetconstant) );
1446 
1447  if( !*success )
1448  {
1449  SCIPdebugMsg(scip, "failed to compute facet of convex hull\n");
1450  goto TERMINATE;
1451  }
1452 
1453  SCIProwprepSetLocal(rowprep, TRUE);
1454  SCIProwprepAddConstant(rowprep, facetconstant);
1455  for( i = 0; i < nlhdlrexprdata->nleafs; ++i )
1456  {
1457  SCIP_CALL( SCIPaddRowprepTerm(scip, rowprep, SCIPgetVarExprVar(nlhdlrexprdata->leafexprs[i]), SCIProwprepGetCoefs(rowprep)[i]) );
1458  }
1459 
1460 #ifdef SCIP_DEBUG
1461  SCIPinfoMessage(scip, NULL, "computed estimator: ");
1462  SCIPprintRowprep(scip, rowprep, NULL);
1463 #endif
1464 
1465  TERMINATE:
1466  SCIPfreeBufferArray(scip, &box);
1467  SCIPfreeBufferArray(scip, &xstar);
1468 
1469  return SCIP_OKAY;
1470 }
1471 
1472 /** adds an estimator computed via a gradient to a given rowprep */
1473 static
1475  SCIP* scip, /**< SCIP data structure */
1476  SCIP_NLHDLREXPRDATA* nlhdlrexprdata, /**< nonlinear handler expression data */
1477  SCIP_SOL* sol, /**< solution to use */
1478  SCIP_Real auxvalue, /**< value of nlexpr in sol - we may not be able to take this value
1479  from nlexpr if it was evaluated at a different sol recently */
1480  SCIP_ROWPREP* rowprep, /**< rowprep where to store estimator */
1481  SCIP_Bool* success /**< buffer to store whether successful */
1482  )
1483 {
1484  SCIP_EXPR* nlexpr;
1485  SCIP_Real QUAD(constant);
1486  int i;
1487 
1488  assert(nlhdlrexprdata != NULL);
1489  assert(rowprep != NULL);
1490  assert(success != NULL);
1491 
1492  nlexpr = nlhdlrexprdata->nlexpr;
1493  assert(nlexpr != NULL);
1494 
1495 #ifdef SCIP_DEBUG
1496  SCIPinfoMessage(scip, NULL, "estimate expression ");
1497  SCIPprintExpr(scip, nlexpr, NULL);
1498  SCIPinfoMessage(scip, NULL, " by gradient\n");
1499 #endif
1500 
1501  *success = FALSE;
1502 
1503  /* evaluation error -> skip */
1504  if( auxvalue == SCIP_INVALID )
1505  {
1506  SCIPdebugMsg(scip, "evaluation error / too large value (%g) for %p\n", auxvalue, (void*)nlexpr);
1507  return SCIP_OKAY;
1508  }
1509 
1510  /* compute gradient (TODO: this also re-evaluates (soltag=0), which shouldn't be necessary unless we tried ConvexSecant before) */
1511  SCIP_CALL( SCIPevalExprGradient(scip, nlexpr, sol, 0L) );
1512 
1513  /* gradient evaluation error -> skip */
1514  if( SCIPexprGetDerivative(nlexpr) == SCIP_INVALID )
1515  {
1516  SCIPdebugMsg(scip, "gradient evaluation error for %p\n", (void*)nlexpr);
1517  return SCIP_OKAY;
1518  }
1519 
1520  /* add gradient underestimator to rowprep: f(sol) + (x - sol) \nabla f(sol)
1521  * constant will store f(sol) - sol * \nabla f(sol)
1522  * to avoid some cancellation errors when linear variables take huge values (like 1e20),
1523  * we use double-double arithemtic here
1524  */
1525  QUAD_ASSIGN(constant, SCIPexprGetEvalValue(nlexpr)); /* f(sol) */
1526  for( i = 0; i < nlhdlrexprdata->nleafs; ++i )
1527  {
1528  SCIP_VAR* var;
1529  SCIP_Real deriv;
1530  SCIP_Real varval;
1531 
1532  assert(SCIPexprGetDiffTag(nlhdlrexprdata->leafexprs[i]) == SCIPexprGetDiffTag(nlexpr));
1533  deriv = SCIPexprGetDerivative(nlhdlrexprdata->leafexprs[i]);
1534  if( deriv == SCIP_INVALID )
1535  {
1536  SCIPdebugMsg(scip, "gradient evaluation error for component %d of %p\n", i, (void*)nlexpr);
1537  return SCIP_OKAY;
1538  }
1539 
1540  var = SCIPgetVarExprVar(nlhdlrexprdata->leafexprs[i]);
1541  assert(var != NULL);
1542 
1543  varval = SCIPgetSolVal(scip, sol, var);
1544 
1545  SCIPdebugMsg(scip, "add %g * (<%s> - %g) to rowprep\n", deriv, SCIPvarGetName(var), varval);
1546 
1547  /* add deriv * var to rowprep and deriv * (-varval) to constant */
1548  SCIP_CALL( SCIPaddRowprepTerm(scip, rowprep, var, deriv) );
1549  SCIPquadprecSumQD(constant, constant, -deriv * varval);
1550  }
1551 
1552  SCIProwprepAddConstant(rowprep, QUAD_TO_DBL(constant));
1553  SCIProwprepSetLocal(rowprep, FALSE);
1554 
1555  *success = TRUE;
1556 
1557  return SCIP_OKAY;
1558 }
1559 
1560 /** adds an estimator generated by putting a secant through the coordinates given by the two closest integer points */
1561 static
1563  SCIP* scip, /**< SCIP data structure */
1564  SCIP_NLHDLR* nlhdlr, /**< nonlinear handler */
1565  SCIP_NLHDLREXPRDATA* nlhdlrexprdata, /**< nonlinear handler expression data */
1566  SCIP_SOL* sol, /**< solution to use, unless usemidpoint is TRUE */
1567  SCIP_ROWPREP* rowprep, /**< rowprep where to store estimator */
1568  SCIP_Bool* success /**< buffer to store whether successful */
1569  )
1570 {
1571  SCIP_NLHDLRDATA* nlhdlrdata;
1572  SCIP_EXPR* nlexpr;
1573  SCIP_VAR* var;
1574  SCIP_Real x;
1575  SCIP_Real left, right;
1576  SCIP_Real fleft, fright;
1577 
1578  assert(nlhdlrexprdata != NULL);
1579  assert(nlhdlrexprdata->nleafs == 1);
1580  assert(rowprep != NULL);
1581  assert(success != NULL);
1582 
1583  nlexpr = nlhdlrexprdata->nlexpr;
1584  assert(nlexpr != NULL);
1585 
1586  *success = FALSE;
1587 
1588  nlhdlrdata = SCIPnlhdlrGetData(nlhdlr);
1589  assert(nlhdlrdata != NULL);
1590 
1591  var = SCIPgetVarExprVar(nlhdlrexprdata->leafexprs[0]);
1592  assert(var != NULL);
1593 
1594  x = SCIPgetSolVal(scip, sol, var);
1595 
1596 #ifdef SCIP_DEBUG
1597  SCIPinfoMessage(scip, NULL, "estimate expression ");
1598  SCIPprintExpr(scip, nlexpr, NULL);
1599  SCIPinfoMessage(scip, NULL, " by secant\n");
1600  SCIPinfoMessage(scip, NULL, "integral variable <%s> = %g [%g,%g]\n", SCIPvarGetName(var),
1601  x, SCIPvarGetLbGlobal(var), SCIPvarGetUbGlobal(var));
1602 #endif
1603 
1604  /* find out coordinates of var left and right to sol */
1605  if( SCIPisIntegral(scip, x) )
1606  {
1607  x = SCIPround(scip, x);
1608  if( SCIPisEQ(scip, x, SCIPvarGetLbGlobal(var)) )
1609  {
1610  left = x;
1611  right = left + 1.0;
1612  }
1613  else
1614  {
1615  right = x;
1616  left = right - 1.0;
1617  }
1618  }
1619  else
1620  {
1621  left = SCIPfloor(scip, x);
1622  right = SCIPceil(scip, x);
1623  }
1624  assert(left != right);
1625 
1626  /* now evaluate at left and right */
1627  if( nlhdlrdata->evalsol == NULL )
1628  {
1629  SCIP_CALL( SCIPcreateSol(scip, &nlhdlrdata->evalsol, NULL) );
1630  }
1631 
1632  SCIP_CALL( SCIPsetSolVal(scip, nlhdlrdata->evalsol, var, left) );
1633  SCIP_CALL( SCIPevalExpr(scip, nlexpr, nlhdlrdata->evalsol, 0L) );
1634 
1635  /* evaluation error or a too large constant -> skip */
1636  fleft = SCIPexprGetEvalValue(nlexpr);
1637  if( SCIPisInfinity(scip, REALABS(fleft)) )
1638  {
1639  SCIPdebugMsg(scip, "evaluation error / too large value (%g) for %p\n", SCIPexprGetEvalValue(nlexpr), (void*)nlexpr);
1640  return SCIP_OKAY;
1641  }
1642 
1643  SCIP_CALL( SCIPsetSolVal(scip, nlhdlrdata->evalsol, var, right) );
1644  SCIP_CALL( SCIPevalExpr(scip, nlexpr, nlhdlrdata->evalsol, 0L) );
1645 
1646  /* evaluation error or a too large constant -> skip */
1647  fright = SCIPexprGetEvalValue(nlexpr);
1648  if( SCIPisInfinity(scip, REALABS(fright)) )
1649  {
1650  SCIPdebugMsg(scip, "evaluation error / too large value (%g) for %p\n", SCIPexprGetEvalValue(nlexpr), (void*)nlexpr);
1651  return SCIP_OKAY;
1652  }
1653 
1654  SCIPdebugMsg(scip, "f(%g)=%g, f(%g)=%g\n", left, fleft, right, fright);
1655 
1656  /* skip if too steep
1657  * for clay0204h, this resulted in a wrong cut from f(0)=1e12 f(1)=0.99998,
1658  * since due to limited precision, this was handled as if f(1)=1
1659  */
1660  if( (!SCIPisZero(scip, fleft) && REALABS(fright/fleft)*SCIPepsilon(scip) > 1.0) ||
1661  (!SCIPisZero(scip, fright) && REALABS(fleft/fright)*SCIPepsilon(scip) > 1.0) )
1662  {
1663  SCIPdebugMsg(scip, "function is too steep, abandoning\n");
1664  return SCIP_OKAY;
1665  }
1666 
1667  /* now add f(left) + (f(right) - f(left)) * (x - left) as estimator to rowprep */
1668  SCIP_CALL( SCIPaddRowprepTerm(scip, rowprep, var, fright - fleft) );
1669  SCIProwprepAddConstant(rowprep, fleft - (fright - fleft) * left);
1670  SCIProwprepSetLocal(rowprep, FALSE);
1671 
1672  *success = TRUE;
1673 
1674  return SCIP_OKAY;
1675 }
1676 
1677 /*
1678  * Callback methods of convex nonlinear handler
1679  */
1680 
1681 /** free handler data of convex or concave nlhdlr */
1682 static
1683 SCIP_DECL_NLHDLRFREEHDLRDATA(nlhdlrfreeHdlrDataConvexConcave)
1684 { /*lint --e{715}*/
1685  assert(scip != NULL);
1686  assert(nlhdlrdata != NULL);
1687  assert(*nlhdlrdata != NULL);
1688  assert((*nlhdlrdata)->evalsol == NULL);
1689 
1690  SCIPfreeBlockMemory(scip, nlhdlrdata);
1691 
1692  return SCIP_OKAY;
1693 }
1694 
1695 /** callback to free expression specific data */
1696 static
1697 SCIP_DECL_NLHDLRFREEEXPRDATA(nlhdlrfreeExprDataConvexConcave)
1698 { /*lint --e{715}*/
1699  assert(scip != NULL);
1700  assert(nlhdlrexprdata != NULL);
1701  assert(*nlhdlrexprdata != NULL);
1702 
1703  SCIPfreeBlockMemoryArrayNull(scip, &(*nlhdlrexprdata)->leafexprs, (*nlhdlrexprdata)->nleafs);
1704  SCIP_CALL( SCIPreleaseExpr(scip, &(*nlhdlrexprdata)->nlexpr) );
1705  SCIPhashmapFree(&(*nlhdlrexprdata)->nlexpr2origexpr);
1706 
1707  SCIPfreeBlockMemory(scip, nlhdlrexprdata);
1708 
1709  return SCIP_OKAY;
1710 }
1711 
1712 /** deinitialization of problem-specific data */
1713 static
1714 SCIP_DECL_NLHDLREXIT(nlhdlrExitConvex)
1715 {
1716  SCIP_NLHDLRDATA* nlhdlrdata;
1717 
1718  nlhdlrdata = SCIPnlhdlrGetData(nlhdlr);
1719  assert(nlhdlrdata != NULL);
1720 
1721  if( nlhdlrdata->evalsol != NULL )
1722  {
1723  SCIP_CALL( SCIPfreeSol(scip, &nlhdlrdata->evalsol) );
1724  }
1725 
1726  return SCIP_OKAY;
1727 }
1728 
1729 /** checks whether expression (or -expression) is convex, possibly after introducing auxiliary variables */
1730 static
1731 SCIP_DECL_NLHDLRDETECT(nlhdlrDetectConvex)
1732 { /*lint --e{715}*/
1733  SCIP_NLHDLRDATA* nlhdlrdata;
1734  SCIP_EXPR* nlexpr = NULL;
1735  SCIP_HASHMAP* nlexpr2origexpr;
1736  int nleafs = 0;
1737 
1738  assert(scip != NULL);
1739  assert(nlhdlr != NULL);
1740  assert(expr != NULL);
1741  assert(enforcing != NULL);
1742  assert(participating != NULL);
1743  assert(nlhdlrexprdata != NULL);
1744 
1745  /* we currently do not participate if only activity computation is required */
1746  if( (*enforcing & SCIP_NLHDLR_METHOD_SEPABOTH) == SCIP_NLHDLR_METHOD_SEPABOTH )
1747  return SCIP_OKAY;
1748 
1749  /* ignore pure constants and variables */
1750  if( SCIPexprGetNChildren(expr) == 0 )
1751  return SCIP_OKAY;
1752 
1753  nlhdlrdata = SCIPnlhdlrGetData(nlhdlr);
1754  assert(nlhdlrdata != NULL);
1755  assert(nlhdlrdata->isnlhdlrconvex);
1756 
1757  SCIPdebugMsg(scip, "nlhdlr_convex detect for expr %p\n", (void*)expr);
1758 
1759  /* initialize mapping from copied expression to original one
1760  * 20 is not a bad estimate for the size of convex subexpressions that we can usually discover
1761  * when expressions will be allowed to store "user"data, we could get rid of this hashmap (TODO)
1762  */
1763  SCIP_CALL( SCIPhashmapCreate(&nlexpr2origexpr, SCIPblkmem(scip), 20) );
1764 
1765  if( (*enforcing & SCIP_NLHDLR_METHOD_SEPABELOW) == 0 ) /* if no separation below yet */
1766  {
1767  SCIP_CALL( constructExpr(scip, nlhdlrdata, &nlexpr, nlexpr2origexpr, &nleafs, expr,
1769  if( nlexpr != NULL )
1770  {
1771  assert(SCIPexprGetNChildren(nlexpr) > 0); /* should not be trivial */
1772 
1773  *participating |= SCIP_NLHDLR_METHOD_SEPABELOW;
1774 
1775  SCIPdebugMsg(scip, "detected expr %p to be convex -> can enforce expr <= auxvar\n", (void*)expr);
1776  }
1777  else
1778  {
1779  SCIP_CALL( SCIPhashmapRemoveAll(nlexpr2origexpr) );
1780  }
1781  }
1782 
1783  if( (*enforcing & SCIP_NLHDLR_METHOD_SEPAABOVE) == 0 && nlexpr == NULL ) /* if no separation above and not convex */
1784  {
1785  SCIP_CALL( constructExpr(scip, nlhdlrdata, &nlexpr, nlexpr2origexpr, &nleafs, expr,
1787  if( nlexpr != NULL )
1788  {
1789  assert(SCIPexprGetNChildren(nlexpr) > 0); /* should not be trivial */
1790 
1791  *participating |= SCIP_NLHDLR_METHOD_SEPAABOVE;
1792 
1793  SCIPdebugMsg(scip, "detected expr %p to be concave -> can enforce expr >= auxvar\n", (void*)expr);
1794  }
1795  }
1796 
1797  /* everything we participate in we also enforce */
1798  *enforcing |= *participating;
1799 
1800  assert(*participating || nlexpr == NULL);
1801  if( !*participating )
1802  {
1803  SCIPhashmapFree(&nlexpr2origexpr);
1804  return SCIP_OKAY;
1805  }
1806 
1807  /* create the expression data of the nonlinear handler
1808  * notify conshdlr about expr for which we will require auxiliary variables
1809  */
1810  SCIP_CALL( createNlhdlrExprData(scip, nlhdlrdata, nlhdlrexprdata, expr, nlexpr, nlexpr2origexpr, nleafs, *participating) );
1811 
1812  return SCIP_OKAY;
1813 }
1814 
1815 /** auxiliary evaluation callback */
1816 static
1817 SCIP_DECL_NLHDLREVALAUX(nlhdlrEvalAuxConvexConcave)
1818 { /*lint --e{715}*/
1819  assert(nlhdlrexprdata != NULL);
1820  assert(nlhdlrexprdata->nlexpr != NULL);
1821  assert(auxvalue != NULL);
1822 
1823  SCIP_CALL( SCIPevalExpr(scip, nlhdlrexprdata->nlexpr, sol, 0L) );
1824  *auxvalue = SCIPexprGetEvalValue(nlhdlrexprdata->nlexpr);
1825 
1826  return SCIP_OKAY;
1827 }
1828 
1829 /** init sepa callback that initializes LP */
1830 static
1831 SCIP_DECL_NLHDLRINITSEPA(nlhdlrInitSepaConvex)
1832 { /*lint --e{715}*/
1833  SCIP_EXPR* nlexpr;
1834  SCIP_EXPRCURV curvature;
1835  SCIP_Bool success;
1836  SCIP_ROWPREP* rowprep = NULL;
1837  SCIP_ROW* row;
1838  SCIP_Real lb;
1839  SCIP_Real ub;
1840  SCIP_Real lambda;
1841  SCIP_SOL* sol;
1842  int k;
1843 
1844  assert(scip != NULL);
1845  assert(expr != NULL);
1846  assert(nlhdlrexprdata != NULL);
1847 
1848  /* setup nlhdlrexprdata->leafexprs */
1849  SCIP_CALL( collectLeafs(scip, nlhdlrexprdata) );
1850 
1851  nlexpr = nlhdlrexprdata->nlexpr;
1852  assert(nlexpr != NULL);
1853  assert(SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, (void*)nlexpr) == expr);
1854 
1855  curvature = SCIPexprGetCurvature(nlexpr);
1856  assert(curvature == SCIP_EXPRCURV_CONVEX || curvature == SCIP_EXPRCURV_CONCAVE);
1857 
1858  /* we can only be estimating on the convex side */
1859  if( curvature == SCIP_EXPRCURV_CONVEX )
1860  overestimate = FALSE;
1861  else if( curvature == SCIP_EXPRCURV_CONCAVE )
1862  underestimate = FALSE;
1863  if( !overestimate && !underestimate )
1864  return SCIP_OKAY;
1865 
1866  /* linearizes at 5 different points obtained as convex combination of the lower and upper bound of the variables
1867  * present in the convex expression; whether more weight is given to the lower or upper bound of a variable depends
1868  * on whether the fixing of the variable to that value is better for the objective function
1869  */
1870  SCIP_CALL( SCIPcreateSol(scip, &sol, NULL) );
1871 
1872  *infeasible = FALSE;
1873 
1874  for( k = 0; k < 5; ++k )
1875  {
1876  int i;
1877  lambda = 0.1 * (k+1); /* lambda = 0.1, 0.2, 0.3, 0.4, 0.5 */
1878 
1879  for( i = 0; i < nlhdlrexprdata->nleafs; ++i )
1880  {
1881  SCIP_VAR* var;
1882 
1883  var = SCIPgetVarExprVar(nlhdlrexprdata->leafexprs[i]);
1884 
1885  lb = SCIPvarGetLbGlobal(var);
1886  ub = SCIPvarGetUbGlobal(var);
1887 
1888  /* make sure the absolute values of bounds are not too large */
1889  if( ub > -INITLPMAXVARVAL )
1890  lb = MAX(lb, -INITLPMAXVARVAL);
1891  if( lb < INITLPMAXVARVAL )
1892  ub = MIN(ub, INITLPMAXVARVAL);
1893 
1894  /* in the case when ub < -maxabsbnd or lb > maxabsbnd, we still want to at least make bounds finite */
1895  if( SCIPisInfinity(scip, -lb) )
1896  lb = MIN(-10.0, ub - 0.1*REALABS(ub));
1897  if( SCIPisInfinity(scip, ub) )
1898  ub = MAX( 10.0, lb + 0.1*REALABS(lb));
1899 
1901  SCIP_CALL( SCIPsetSolVal(scip, sol, var, lambda * ub + (1.0 - lambda) * lb) );
1902  else
1903  SCIP_CALL( SCIPsetSolVal(scip, sol, var, lambda * lb + (1.0 - lambda) * ub) );
1904  }
1905 
1907  SCIP_CALL( estimateGradient(scip, nlhdlrexprdata, sol, 0.0, rowprep, &success) );
1908  if( !success )
1909  {
1910  SCIPdebugMsg(scip, "failed to linearize for k = %d\n", k);
1911  if( rowprep != NULL )
1912  SCIPfreeRowprep(scip, &rowprep);
1913  continue;
1914  }
1915 
1916  /* add auxiliary variable */
1917  SCIP_CALL( SCIPaddRowprepTerm(scip, rowprep, SCIPgetExprAuxVarNonlinear(expr), -1.0) );
1918 
1919  /* straighten out numerics */
1920  SCIP_CALL( SCIPcleanupRowprep2(scip, rowprep, NULL, SCIPgetHugeValue(scip), &success) );
1921  if( !success )
1922  {
1923  SCIPdebugMsg(scip, "failed to cleanup rowprep numerics for k = %d\n", k);
1924  if( rowprep != NULL )
1925  SCIPfreeRowprep(scip, &rowprep);
1926  continue;
1927  }
1928 
1929  (void) SCIPsnprintf(SCIProwprepGetName(rowprep), SCIP_MAXSTRLEN, "%sestimate_gradient%p_initsepa_%d",
1930  overestimate ? "over" : "under", (void*)expr, k);
1931  SCIP_CALL( SCIPgetRowprepRowCons(scip, &row, rowprep, cons) );
1932  SCIPfreeRowprep(scip, &rowprep);
1933 
1934 #ifdef SCIP_DEBUG
1935  SCIPinfoMessage(scip, NULL, "initsepa computed row: ");
1936  SCIPprintRow(scip, row, NULL);
1937 #endif
1938 
1939  SCIP_CALL( SCIPaddRow(scip, row, FALSE, infeasible) );
1940  SCIP_CALL( SCIPreleaseRow(scip, &row) );
1941 
1942  if( *infeasible )
1943  break;
1944  }
1945 
1946  SCIP_CALL( SCIPfreeSol(scip, &sol) );
1947 
1948  return SCIP_OKAY;
1949 }
1950 
1951 /** estimator callback */
1952 static
1953 SCIP_DECL_NLHDLRESTIMATE(nlhdlrEstimateConvex)
1954 { /*lint --e{715}*/
1955  SCIP_ROWPREP* rowprep;
1956 
1957  assert(scip != NULL);
1958  assert(expr != NULL);
1959  assert(nlhdlrexprdata != NULL);
1960  assert(nlhdlrexprdata->nlexpr != NULL);
1961  assert(rowpreps != NULL);
1962  assert(success != NULL);
1963 
1964  assert(SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, (void*)nlhdlrexprdata->nlexpr) == expr);
1965 
1966  /* we must be called only for the side that we indicated to participate in during DETECT */
1967  assert(SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONVEX
1968  || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONCAVE);
1969  assert(!overestimate || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONCAVE);
1970  assert( overestimate || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONVEX);
1971 
1972  *success = FALSE;
1973  *addedbranchscores = FALSE;
1974 
1975  /* we can skip eval as nlhdlrEvalAux should have been called for same solution before */
1976  /* SCIP_CALL( nlhdlrExprEval(scip, nlexpr, sol) ); */
1977  assert(auxvalue == SCIPexprGetEvalValue(nlhdlrexprdata->nlexpr)); /* given value (originally from
1978  nlhdlrEvalAuxConvexConcave) should coincide with the one stored in nlexpr */
1979 
1981 
1982  if( nlhdlrexprdata->nleafs == 1 && SCIPexprIsIntegral(nlhdlrexprdata->leafexprs[0]) )
1983  {
1984  SCIP_CALL( estimateConvexSecant(scip, nlhdlr, nlhdlrexprdata, sol, rowprep, success) );
1985 
1986  (void) SCIPsnprintf(SCIProwprepGetName(rowprep), SCIP_MAXSTRLEN, "%sestimate_convexsecant%p_%s%d",
1987  overestimate ? "over" : "under",
1988  (void*)expr,
1989  sol != NULL ? "sol" : "lp",
1990  sol != NULL ? SCIPsolGetIndex(sol) : SCIPgetNLPs(scip));
1991  }
1992 
1993  /* if secant method was not used or failed, then try with gradient */
1994  if( !*success )
1995  {
1996  SCIP_CALL( estimateGradient(scip, nlhdlrexprdata, sol, auxvalue, rowprep, success) );
1997 
1998  (void) SCIPsnprintf(SCIProwprepGetName(rowprep), SCIP_MAXSTRLEN, "%sestimate_convexgradient%p_%s%d",
1999  overestimate ? "over" : "under",
2000  (void*)expr,
2001  sol != NULL ? "sol" : "lp",
2002  sol != NULL ? SCIPsolGetIndex(sol) : SCIPgetNLPs(scip));
2003  }
2004 
2005  if( *success )
2006  {
2007  SCIP_CALL( SCIPsetPtrarrayVal(scip, rowpreps, 0, rowprep) );
2008  }
2009  else
2010  {
2011  SCIPfreeRowprep(scip, &rowprep);
2012  }
2013 
2014  return SCIP_OKAY;
2015 }
2016 
2017 /** include nlhdlr in another scip instance */
2018 static
2019 SCIP_DECL_NLHDLRCOPYHDLR(nlhdlrCopyhdlrConvex)
2020 { /*lint --e{715}*/
2021  assert(targetscip != NULL);
2022  assert(sourcenlhdlr != NULL);
2023  assert(strcmp(SCIPnlhdlrGetName(sourcenlhdlr), CONVEX_NLHDLR_NAME) == 0);
2024 
2025  SCIP_CALL( SCIPincludeNlhdlrConvex(targetscip) );
2026 
2027  return SCIP_OKAY;
2028 }
2029 
2030 /** includes convex nonlinear handler in nonlinear constraint handler */
2032  SCIP* scip /**< SCIP data structure */
2033  )
2034 {
2035  SCIP_NLHDLR* nlhdlr;
2036  SCIP_NLHDLRDATA* nlhdlrdata;
2037 
2038  assert(scip != NULL);
2039 
2040  SCIP_CALL( SCIPallocBlockMemory(scip, &nlhdlrdata) );
2041  nlhdlrdata->isnlhdlrconvex = TRUE;
2042  nlhdlrdata->evalsol = NULL;
2043 
2045  CONVEX_NLHDLR_DETECTPRIORITY, CONVEX_NLHDLR_ENFOPRIORITY, nlhdlrDetectConvex, nlhdlrEvalAuxConvexConcave, nlhdlrdata) );
2046  assert(nlhdlr != NULL);
2047 
2048  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONVEX_NLHDLR_NAME "/detectsum",
2049  "whether to run convexity detection when the root of an expression is a non-quadratic sum",
2050  &nlhdlrdata->detectsum, FALSE, DEFAULT_DETECTSUM, NULL, NULL) );
2051 
2052  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONVEX_NLHDLR_NAME "/extendedform",
2053  "whether to create extended formulations instead of looking for maximal convex expressions",
2054  &nlhdlrdata->extendedform, FALSE, DEFAULT_EXTENDEDFORM, NULL, NULL) );
2055 
2056  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONVEX_NLHDLR_NAME "/cvxquadratic",
2057  "whether to use convexity check on quadratics",
2058  &nlhdlrdata->cvxquadratic, TRUE, DEFAULT_CVXQUADRATIC_CONVEX, NULL, NULL) );
2059 
2060  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONVEX_NLHDLR_NAME "/cvxsignomial",
2061  "whether to use convexity check on signomials",
2062  &nlhdlrdata->cvxsignomial, TRUE, DEFAULT_CVXSIGNOMIAL, NULL, NULL) );
2063 
2064  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONVEX_NLHDLR_NAME "/cvxprodcomp",
2065  "whether to use convexity check on product composition f(h)*h",
2066  &nlhdlrdata->cvxprodcomp, TRUE, DEFAULT_CVXPRODCOMP, NULL, NULL) );
2067 
2068  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONVEX_NLHDLR_NAME "/handletrivial",
2069  "whether to also handle trivial convex expressions",
2070  &nlhdlrdata->handletrivial, TRUE, DEFAULT_HANDLETRIVIAL, NULL, NULL) );
2071 
2072  SCIPnlhdlrSetFreeHdlrData(nlhdlr, nlhdlrfreeHdlrDataConvexConcave);
2073  SCIPnlhdlrSetCopyHdlr(nlhdlr, nlhdlrCopyhdlrConvex);
2074  SCIPnlhdlrSetFreeExprData(nlhdlr, nlhdlrfreeExprDataConvexConcave);
2075  SCIPnlhdlrSetSepa(nlhdlr, nlhdlrInitSepaConvex, NULL, nlhdlrEstimateConvex, NULL);
2076  SCIPnlhdlrSetInitExit(nlhdlr, NULL, nlhdlrExitConvex);
2077 
2078  return SCIP_OKAY;
2079 }
2080 
2081 /*
2082  * Callback methods of concave nonlinear handler
2083  */
2084 
2085 /** deinitialization of problem-specific data */
2086 static
2087 SCIP_DECL_NLHDLREXIT(nlhdlrExitConcave)
2088 {
2089  SCIP_NLHDLRDATA* nlhdlrdata;
2090 
2091  nlhdlrdata = SCIPnlhdlrGetData(nlhdlr);
2092  assert(nlhdlrdata != NULL);
2093 
2094  if( nlhdlrdata->evalsol != NULL )
2095  {
2096  SCIP_CALL( SCIPfreeSol(scip, &nlhdlrdata->evalsol) );
2097  }
2098 
2099  return SCIP_OKAY;
2100 }
2101 
2102 /** checks whether expression (or -expression) is concave, possibly after introducing auxiliary variables */
2103 static
2104 SCIP_DECL_NLHDLRDETECT(nlhdlrDetectConcave)
2105 { /*lint --e{715}*/
2106  SCIP_NLHDLRDATA* nlhdlrdata;
2107  SCIP_EXPR* nlexpr = NULL;
2108  SCIP_HASHMAP* nlexpr2origexpr;
2109  int nleafs = 0;
2110 
2111  assert(scip != NULL);
2112  assert(nlhdlr != NULL);
2113  assert(expr != NULL);
2114  assert(enforcing != NULL);
2115  assert(participating != NULL);
2116  assert(nlhdlrexprdata != NULL);
2117 
2118  /* we currently do not participate if only activity computation is required */
2119  if( (*enforcing & SCIP_NLHDLR_METHOD_SEPABOTH) == SCIP_NLHDLR_METHOD_SEPABOTH )
2120  return SCIP_OKAY;
2121 
2122  /* ignore pure constants and variables */
2123  if( SCIPexprGetNChildren(expr) == 0 )
2124  return SCIP_OKAY;
2125 
2126  nlhdlrdata = SCIPnlhdlrGetData(nlhdlr);
2127  assert(nlhdlrdata != NULL);
2128  assert(!nlhdlrdata->isnlhdlrconvex);
2129 
2130  SCIPdebugMsg(scip, "nlhdlr_concave detect for expr %p\n", (void*)expr);
2131 
2132  /* initialize mapping from copied expression to original one
2133  * 20 is not a bad estimate for the size of concave subexpressions that we can usually discover
2134  * when expressions will be allowed to store "user"data, we could get rid of this hashmap (TODO)
2135  */
2136  SCIP_CALL( SCIPhashmapCreate(&nlexpr2origexpr, SCIPblkmem(scip), 20) );
2137 
2138  if( (*enforcing & SCIP_NLHDLR_METHOD_SEPABELOW) == 0 ) /* if no separation below yet */
2139  {
2140  SCIP_CALL( constructExpr(scip, nlhdlrdata, &nlexpr, nlexpr2origexpr, &nleafs, expr,
2142 
2143  if( nlexpr != NULL && nleafs > SCIP_MAXVERTEXPOLYDIM )
2144  {
2145  SCIPdebugMsg(scip, "Too many variables (%d) in constructed expression. Will not be able to estimate. Rejecting.\n", nleafs);
2146  SCIP_CALL( SCIPreleaseExpr(scip, &nlexpr) );
2147  }
2148 
2149  if( nlexpr != NULL )
2150  {
2151  assert(SCIPexprGetNChildren(nlexpr) > 0); /* should not be trivial */
2152 
2153  *participating |= SCIP_NLHDLR_METHOD_SEPABELOW;
2154 
2155  SCIPdebugMsg(scip, "detected expr %p to be concave -> can enforce expr <= auxvar\n", (void*)expr);
2156  }
2157  else
2158  {
2159  SCIP_CALL( SCIPhashmapRemoveAll(nlexpr2origexpr) );
2160  }
2161  }
2162 
2163  if( (*enforcing & SCIP_NLHDLR_METHOD_SEPAABOVE) == 0 && nlexpr == NULL ) /* if no separation above and not concave */
2164  {
2165  SCIP_CALL( constructExpr(scip, nlhdlrdata, &nlexpr, nlexpr2origexpr, &nleafs, expr,
2167 
2168  if( nlexpr != NULL && nleafs > SCIP_MAXVERTEXPOLYDIM )
2169  {
2170  SCIPdebugMsg(scip, "Too many variables (%d) in constructed expression. Will not be able to estimate. Rejecting.\n", nleafs);
2171  SCIP_CALL( SCIPreleaseExpr(scip, &nlexpr) );
2172  }
2173 
2174  if( nlexpr != NULL )
2175  {
2176  assert(SCIPexprGetNChildren(nlexpr) > 0); /* should not be trivial */
2177 
2178  *participating |= SCIP_NLHDLR_METHOD_SEPAABOVE;
2179 
2180  SCIPdebugMsg(scip, "detected expr %p to be convex -> can enforce expr >= auxvar\n", (void*)expr);
2181  }
2182  }
2183 
2184  /* everything we participate in we also enforce (at the moment) */
2185  *enforcing |= *participating;
2186 
2187  assert(*participating || nlexpr == NULL);
2188  if( !*participating )
2189  {
2190  SCIPhashmapFree(&nlexpr2origexpr);
2191  return SCIP_OKAY;
2192  }
2193 
2194  /* create the expression data of the nonlinear handler
2195  * notify conshdlr about expr for which we will require auxiliary variables and use activity
2196  */
2197  SCIP_CALL( createNlhdlrExprData(scip, nlhdlrdata, nlhdlrexprdata, expr, nlexpr, nlexpr2origexpr, nleafs, *participating) );
2198 
2199  return SCIP_OKAY;
2200 }
2201 
2202 /** init sepa callback that initializes LP */
2203 static
2204 SCIP_DECL_NLHDLRINITSEPA(nlhdlrInitSepaConcave)
2205 {
2206  SCIP_EXPR* nlexpr;
2207  SCIP_EXPRCURV curvature;
2208  SCIP_Bool success;
2209  SCIP_ROWPREP* rowprep = NULL;
2210  SCIP_ROW* row;
2211 
2212  assert(scip != NULL);
2213  assert(expr != NULL);
2214  assert(nlhdlrexprdata != NULL);
2215 
2216  nlexpr = nlhdlrexprdata->nlexpr;
2217  assert(nlexpr != NULL);
2218  assert(SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, (void*)nlexpr) == expr);
2219 
2220  /* setup nlhdlrexprdata->leafexprs */
2221  SCIP_CALL( collectLeafs(scip, nlhdlrexprdata) );
2222 
2223  curvature = SCIPexprGetCurvature(nlexpr);
2224  assert(curvature == SCIP_EXPRCURV_CONVEX || curvature == SCIP_EXPRCURV_CONCAVE);
2225  /* we can only be estimating on non-convex side */
2226  if( curvature == SCIP_EXPRCURV_CONCAVE )
2227  overestimate = FALSE;
2228  else if( curvature == SCIP_EXPRCURV_CONVEX )
2229  underestimate = FALSE;
2230  if( !overestimate && !underestimate )
2231  return SCIP_OKAY;
2232 
2233  /* compute estimator and store in rowprep */
2235  SCIP_CALL( estimateVertexPolyhedral(scip, conshdlr, nlhdlr, nlhdlrexprdata, NULL, TRUE, overestimate,
2236  overestimate ? SCIPinfinity(scip) : -SCIPinfinity(scip), rowprep, &success) );
2237  if( !success )
2238  {
2239  SCIPdebugMsg(scip, "failed to compute facet of convex hull\n");
2240  goto TERMINATE;
2241  }
2242 
2243  /* add auxiliary variable */
2244  SCIP_CALL( SCIPaddRowprepTerm(scip, rowprep, SCIPgetExprAuxVarNonlinear(expr), -1.0) );
2245 
2246  /* straighten out numerics */
2247  SCIP_CALL( SCIPcleanupRowprep2(scip, rowprep, NULL, SCIPgetHugeValue(scip), &success) );
2248  if( !success )
2249  {
2250  SCIPdebugMsg(scip, "failed to cleanup rowprep numerics\n");
2251  goto TERMINATE;
2252  }
2253 
2254  (void) SCIPsnprintf(SCIProwprepGetName(rowprep), SCIP_MAXSTRLEN, "%sestimate_concave%p_initsepa",
2255  overestimate ? "over" : "under", (void*)expr);
2256  SCIP_CALL( SCIPgetRowprepRowCons(scip, &row, rowprep, cons) );
2257 
2258 #ifdef SCIP_DEBUG
2259  SCIPinfoMessage(scip, NULL, "initsepa computed row: ");
2260  SCIPprintRow(scip, row, NULL);
2261 #endif
2262 
2263  SCIP_CALL( SCIPaddRow(scip, row, FALSE, infeasible) );
2264  SCIP_CALL( SCIPreleaseRow(scip, &row) );
2265 
2266  TERMINATE:
2267  if( rowprep != NULL )
2268  SCIPfreeRowprep(scip, &rowprep);
2269 
2270  return SCIP_OKAY;
2271 }
2272 
2273 /** estimator callback */
2274 static
2275 SCIP_DECL_NLHDLRESTIMATE(nlhdlrEstimateConcave)
2276 { /*lint --e{715}*/
2277  SCIP_ROWPREP* rowprep;
2278 
2279  assert(scip != NULL);
2280  assert(expr != NULL);
2281  assert(nlhdlrexprdata != NULL);
2282  assert(nlhdlrexprdata->nlexpr != NULL);
2283  assert(rowpreps != NULL);
2284  assert(success != NULL);
2285 
2286  assert(SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, (void*)nlhdlrexprdata->nlexpr) == expr);
2287 
2288  /* we must be called only for the side that we indicated to participate in during DETECT */
2289  assert(SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONVEX
2290  || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONCAVE);
2291  assert(!overestimate || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONVEX);
2292  assert( overestimate || SCIPexprGetCurvature(nlhdlrexprdata->nlexpr) == SCIP_EXPRCURV_CONCAVE);
2293 
2294  *success = FALSE;
2295  *addedbranchscores = FALSE;
2296 
2298 
2299  SCIP_CALL( estimateVertexPolyhedral(scip, conshdlr, nlhdlr, nlhdlrexprdata, sol, FALSE, overestimate, targetvalue, rowprep, success) );
2300 
2301  if( *success )
2302  {
2303  SCIP_CALL( SCIPsetPtrarrayVal(scip, rowpreps, 0, rowprep) );
2304 
2305  (void) SCIPsnprintf(SCIProwprepGetName(rowprep), SCIP_MAXSTRLEN, "%sestimate_concave%p_%s%d",
2306  overestimate ? "over" : "under",
2307  (void*)expr,
2308  sol != NULL ? "sol" : "lp",
2309  sol != NULL ? SCIPsolGetIndex(sol) : SCIPgetNLPs(scip));
2310  }
2311  else
2312  {
2313  SCIPfreeRowprep(scip, &rowprep);
2314  }
2315 
2316  if( addbranchscores )
2317  {
2318  SCIP_Real violation;
2319 
2320  /* check how much is the violation on the side that we estimate */
2321  if( auxvalue == SCIP_INVALID )
2322  {
2323  /* if cannot evaluate, then always branch */
2324  violation = SCIPinfinity(scip);
2325  }
2326  else
2327  {
2328  SCIP_Real auxval;
2329 
2330  /* get value of auxiliary variable of this expression */
2331  assert(SCIPgetExprAuxVarNonlinear(expr) != NULL);
2332  auxval = SCIPgetSolVal(scip, sol, SCIPgetExprAuxVarNonlinear(expr));
2333 
2334  /* compute the violation
2335  * if we underestimate, then we enforce expr <= auxval, so violation is (positive part of) auxvalue - auxval
2336  * if we overestimate, then we enforce expr >= auxval, so violation is (positive part of) auxval - auxvalue
2337  */
2338  if( !overestimate )
2339  violation = MAX(0.0, auxvalue - auxval);
2340  else
2341  violation = MAX(0.0, auxval - auxvalue);
2342  }
2343  assert(violation >= 0.0);
2344 
2345  /* add violation as branching-score to expressions; the core will take care distributing this onto variables */
2346  if( nlhdlrexprdata->nleafs == 1 )
2347  {
2348  SCIP_EXPR* e;
2349  e = (SCIP_EXPR*)SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, nlhdlrexprdata->leafexprs[0]);
2350  SCIP_CALL( SCIPaddExprsViolScoreNonlinear(scip, &e, 1, violation, sol, addedbranchscores) );
2351  }
2352  else
2353  {
2354  SCIP_EXPR** exprs;
2355  int c;
2356 
2357  /* map leaf expressions back to original expressions
2358  * TODO do this once at end of detect and store in nlhdlrexprdata
2359  */
2360  SCIP_CALL( SCIPallocBufferArray(scip, &exprs, nlhdlrexprdata->nleafs) );
2361  for( c = 0; c < nlhdlrexprdata->nleafs; ++c )
2362  exprs[c] = (SCIP_EXPR*)SCIPhashmapGetImage(nlhdlrexprdata->nlexpr2origexpr, nlhdlrexprdata->leafexprs[c]);
2363 
2364  SCIP_CALL( SCIPaddExprsViolScoreNonlinear(scip, exprs, nlhdlrexprdata->nleafs, violation, sol, addedbranchscores) );
2365 
2366  SCIPfreeBufferArray(scip, &exprs);
2367  }
2368  }
2369 
2370  return SCIP_OKAY;
2371 }
2372 
2373 /** includes nonlinear handler in another scip instance */
2374 static
2375 SCIP_DECL_NLHDLRCOPYHDLR(nlhdlrCopyhdlrConcave)
2376 { /*lint --e{715}*/
2377  assert(targetscip != NULL);
2378  assert(sourcenlhdlr != NULL);
2379  assert(strcmp(SCIPnlhdlrGetName(sourcenlhdlr), CONCAVE_NLHDLR_NAME) == 0);
2380 
2381  SCIP_CALL( SCIPincludeNlhdlrConcave(targetscip) );
2382 
2383  return SCIP_OKAY;
2384 }
2385 
2386 /** includes concave nonlinear handler in nonlinear constraint handler */
2387 SCIP_EXPORT
2389  SCIP* scip /**< SCIP data structure */
2390  )
2391 {
2392  SCIP_NLHDLR* nlhdlr;
2393  SCIP_NLHDLRDATA* nlhdlrdata;
2394 
2395  assert(scip != NULL);
2396 
2397  SCIP_CALL( SCIPallocBlockMemory(scip, &nlhdlrdata) );
2398  nlhdlrdata->isnlhdlrconvex = FALSE;
2399  nlhdlrdata->evalsol = NULL;
2400 
2402  CONCAVE_NLHDLR_DETECTPRIORITY, CONCAVE_NLHDLR_ENFOPRIORITY, nlhdlrDetectConcave, nlhdlrEvalAuxConvexConcave, nlhdlrdata) );
2403  assert(nlhdlr != NULL);
2404 
2405  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONCAVE_NLHDLR_NAME "/detectsum",
2406  "whether to run convexity detection when the root of an expression is a sum",
2407  &nlhdlrdata->detectsum, FALSE, DEFAULT_DETECTSUM, NULL, NULL) );
2408 
2409  /* "extended" formulations of a concave expressions can give worse estimators */
2410  nlhdlrdata->extendedform = FALSE;
2411 
2412  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONCAVE_NLHDLR_NAME "/cvxquadratic",
2413  "whether to use convexity check on quadratics",
2414  &nlhdlrdata->cvxquadratic, TRUE, DEFAULT_CVXQUADRATIC_CONCAVE, NULL, NULL) );
2415 
2416  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONCAVE_NLHDLR_NAME "/cvxsignomial",
2417  "whether to use convexity check on signomials",
2418  &nlhdlrdata->cvxsignomial, TRUE, DEFAULT_CVXSIGNOMIAL, NULL, NULL) );
2419 
2420  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONCAVE_NLHDLR_NAME "/cvxprodcomp",
2421  "whether to use convexity check on product composition f(h)*h",
2422  &nlhdlrdata->cvxprodcomp, TRUE, DEFAULT_CVXPRODCOMP, NULL, NULL) );
2423 
2424  SCIP_CALL( SCIPaddBoolParam(scip, "nlhdlr/" CONCAVE_NLHDLR_NAME "/handletrivial",
2425  "whether to also handle trivial convex expressions",
2426  &nlhdlrdata->handletrivial, TRUE, DEFAULT_HANDLETRIVIAL, NULL, NULL) );
2427 
2428  SCIPnlhdlrSetFreeHdlrData(nlhdlr, nlhdlrfreeHdlrDataConvexConcave);
2429  SCIPnlhdlrSetCopyHdlr(nlhdlr, nlhdlrCopyhdlrConcave);
2430  SCIPnlhdlrSetFreeExprData(nlhdlr, nlhdlrfreeExprDataConvexConcave);
2431  SCIPnlhdlrSetSepa(nlhdlr, nlhdlrInitSepaConcave, NULL, nlhdlrEstimateConcave, NULL);
2432  SCIPnlhdlrSetInitExit(nlhdlr, NULL, nlhdlrExitConcave);
2433 
2434  return SCIP_OKAY;
2435 }
2436 
2437 /** checks whether a given expression is convex or concave w.r.t. the original variables
2438  *
2439  * This function uses the methods that are used in the detection algorithm of the convex nonlinear handler.
2440  */
2442  SCIP* scip, /**< SCIP data structure */
2443  SCIP_EXPR* expr, /**< expression */
2444  SCIP_EXPRCURV curv, /**< curvature to check for */
2445  SCIP_Bool* success, /**< buffer to store whether expression has curvature curv (w.r.t. original variables) */
2446  SCIP_HASHMAP* assumevarfixed /**< hashmap containing variables that should be assumed to be fixed, or NULL */
2447  )
2448 {
2449  SCIP_NLHDLRDATA nlhdlrdata;
2450  SCIP_EXPR* rootnlexpr;
2451  SCIP_HASHMAP* nlexpr2origexpr;
2452  int nleafs;
2453 
2454  assert(expr != NULL);
2455  assert(curv != SCIP_EXPRCURV_UNKNOWN);
2456  assert(success != NULL);
2457 
2458  /* create temporary hashmap */
2459  SCIP_CALL( SCIPhashmapCreate(&nlexpr2origexpr, SCIPblkmem(scip), 20) );
2460 
2461  /* prepare nonlinear handler data */
2462  nlhdlrdata.isnlhdlrconvex = TRUE;
2463  nlhdlrdata.evalsol = NULL;
2464  nlhdlrdata.detectsum = TRUE;
2465  nlhdlrdata.extendedform = FALSE;
2466  nlhdlrdata.cvxquadratic = TRUE;
2467  nlhdlrdata.cvxsignomial = TRUE;
2468  nlhdlrdata.cvxprodcomp = TRUE;
2469  nlhdlrdata.handletrivial = TRUE;
2470 
2471  SCIP_CALL( constructExpr(scip, &nlhdlrdata, &rootnlexpr, nlexpr2origexpr, &nleafs, expr, curv, assumevarfixed, FALSE, success) );
2472 
2473  /* free created expression */
2474  if( rootnlexpr != NULL )
2475  {
2476  SCIP_CALL( SCIPreleaseExpr(scip, &rootnlexpr) );
2477  }
2478 
2479  /* free hashmap */
2480  SCIPhashmapFree(&nlexpr2origexpr);
2481 
2482  return SCIP_OKAY;
2483 }
static SCIP_RETCODE estimateConvexSecant(SCIP *scip, SCIP_NLHDLR *nlhdlr, SCIP_NLHDLREXPRDATA *nlhdlrexprdata, SCIP_SOL *sol, SCIP_ROWPREP *rowprep, SCIP_Bool *success)
void SCIPexprGetQuadraticData(SCIP_EXPR *expr, SCIP_Real *constant, int *nlinexprs, SCIP_EXPR ***linexprs, SCIP_Real **lincoefs, int *nquadexprs, int *nbilinexprs, SCIP_Real **eigenvalues, SCIP_Real **eigenvectors)
Definition: expr.c:4055
static SCIP_Bool exprstackIsEmpty(EXPRSTACK *exprstack)
#define DECL_CURVCHECK(x)
const char * SCIPexprcurvGetName(SCIP_EXPRCURV curv)
Definition: exprcurv.c:576
SCIP_RETCODE SCIPexpriterInit(SCIP_EXPRITER *iterator, SCIP_EXPR *expr, SCIP_EXPRITER_TYPE type, SCIP_Bool allowrevisit)
Definition: expriter.c:491
SCIP_RETCODE SCIPduplicateExprShallow(SCIP *scip, SCIP_EXPR *expr, SCIP_EXPR **copyexpr, SCIP_DECL_EXPR_OWNERCREATE((*ownercreate)), void *ownercreatedata)
Definition: scip_expr.c:1291
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:84
SCIP_MONOTONE
Definition: type_expr.h:57
SCIP_RETCODE SCIPprintExpr(SCIP *scip, SCIP_EXPR *expr, FILE *file)
Definition: scip_expr.c:1476
SCIP_Bool SCIPisRelEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
#define SCIPquadprecSumQD(r, a, b)
Definition: dbldblarith.h:53
SCIP_RETCODE SCIPcomputeExprQuadraticCurvature(SCIP *scip, SCIP_EXPR *expr, SCIP_EXPRCURV *curv, SCIP_HASHMAP *assumevarfixed, SCIP_Bool storeeigeninfo)
Definition: scip_expr.c:2549
SCIP_RETCODE SCIPevalExprActivity(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1706
#define CONCAVE_NLHDLR_NAME
Definition: nlhdlr_convex.c:41
#define SCIP_NLHDLR_METHOD_SEPABELOW
Definition: type_nlhdlr.h:42
int SCIPexprGetNChildren(SCIP_EXPR *expr)
Definition: expr.c:3798
#define DEFAULT_EXTENDEDFORM
Definition: nlhdlr_convex.c:47
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:17910
#define SCIP_MAXSTRLEN
Definition: def.h:293
void SCIPnlhdlrSetFreeExprData(SCIP_NLHDLR *nlhdlr, SCIP_DECL_NLHDLRFREEEXPRDATA((*freeexprdata)))
Definition: nlhdlr.c:85
#define DEFAULT_CVXPRODCOMP
Definition: nlhdlr_convex.c:51
SCIP_NLHDLREXPRDATA * nlhdlrexprdata
Definition: nlhdlr_convex.c:96
const char * SCIPexprhdlrGetName(SCIP_EXPRHDLR *exprhdlr)
Definition: expr.c:525
SCIP_RETCODE SCIPcleanupRowprep2(SCIP *scip, SCIP_ROWPREP *rowprep, SCIP_SOL *sol, SCIP_Real maxcoefbound, SCIP_Bool *success)
int SCIPcalcMemGrowSize(SCIP *scip, int num)
Definition: scip_mem.c:130
SCIP_Real SCIPvarGetLbLocal(SCIP_VAR *var)
Definition: var.c:17966
SCIP_Real SCIPgetConstantExprSum(SCIP_EXPR *expr)
Definition: expr_sum.c:1172
#define CONCAVE_NLHDLR_DESC
Definition: nlhdlr_convex.c:42
SCIP_EXPR * SCIPexpriterSkipDFS(SCIP_EXPRITER *iterator)
Definition: expriter.c:920
void SCIPprintRowprep(SCIP *scip, SCIP_ROWPREP *rowprep, FILE *file)
Definition: misc_rowprep.c:769
void * SCIPhashmapEntryGetOrigin(SCIP_HASHMAPENTRY *entry)
Definition: misc.c:3500
SCIP_RETCODE SCIPregisterExprUsageNonlinear(SCIP *scip, SCIP_EXPR *expr, SCIP_Bool useauxvar, SCIP_Bool useactivityforprop, SCIP_Bool useactivityforsepabelow, SCIP_Bool useactivityforsepaabove)
SCIP_Bool SCIPassumeConvexNonlinear(SCIP_CONSHDLR *conshdlr)
SCIP_RETCODE SCIPcheckExprQuadratic(SCIP *scip, SCIP_EXPR *expr, SCIP_Bool *isquadratic)
Definition: scip_expr.c:2340
#define FALSE
Definition: def.h:87
SCIP_RETCODE SCIPhashmapCreate(SCIP_HASHMAP **hashmap, BMS_BLKMEM *blkmem, int mapsize)
Definition: misc.c:3014
void SCIPexprSetCurvature(SCIP_EXPR *expr, SCIP_EXPRCURV curvature)
Definition: expr.c:4006
void SCIPnlhdlrSetFreeHdlrData(SCIP_NLHDLR *nlhdlr, SCIP_DECL_NLHDLRFREEHDLRDATA((*freehdlrdata)))
Definition: nlhdlr.c:74
SCIP_Bool SCIPexprhdlrHasBwdiff(SCIP_EXPRHDLR *exprhdlr)
Definition: expr.c:575
#define CONVEX_NLHDLR_DESC
Definition: nlhdlr_convex.c:37
char * SCIProwprepGetName(SCIP_ROWPREP *rowprep)
Definition: misc_rowprep.c:664
SCIP_Real SCIPinfinity(SCIP *scip)
int SCIPsnprintf(char *t, int len, const char *s,...)
Definition: misc.c:10755
SCIP_Real SCIPgetExponentExprPow(SCIP_EXPR *expr)
Definition: expr_pow.c:3343
#define TRUE
Definition: def.h:86
SCIP_RETCODE SCIPhashsetCreate(SCIP_HASHSET **hashset, BMS_BLKMEM *blkmem, int size)
Definition: misc.c:3699
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:54
SCIP_EXPR ** stack
SCIP_RETCODE SCIPhashmapInsertInt(SCIP_HASHMAP *hashmap, void *origin, int image)
Definition: misc.c:3132
SCIP_Bool SCIPhashsetExists(SCIP_HASHSET *hashset, void *element)
Definition: misc.c:3757
static SCIP_RETCODE estimateVertexPolyhedral(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_NLHDLR *nlhdlr, SCIP_NLHDLREXPRDATA *nlhdlrexprdata, SCIP_SOL *sol, SCIP_Bool usemidpoint, SCIP_Bool overestimate, SCIP_Real targetvalue, SCIP_ROWPREP *rowprep, SCIP_Bool *success)
SCIP_INTERVAL SCIPexprGetActivity(SCIP_EXPR *expr)
Definition: expr.c:3954
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:99
#define CONVEX_NLHDLR_ENFOPRIORITY
Definition: nlhdlr_convex.c:39
void * SCIPhashmapGetImage(SCIP_HASHMAP *hashmap, void *origin)
Definition: misc.c:3201
SCIP_RETCODE SCIPincludeNlhdlrConvex(SCIP *scip)
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_RETCODE SCIPcreateExprVar(SCIP *scip, SCIP_EXPR **expr, SCIP_VAR *var, SCIP_DECL_EXPR_OWNERCREATE((*ownercreate)), void *ownercreatedata)
Definition: expr_var.c:377
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:127
#define QUAD_ASSIGN(a, constant)
Definition: dbldblarith.h:42
defines macros for basic operations in double-double arithmetic giving roughly twice the precision of...
void SCIPcaptureExpr(SCIP_EXPR *expr)
Definition: scip_expr.c:1399
variable expression handler
SCIP_RETCODE SCIPensureRowprepSize(SCIP *scip, SCIP_ROWPREP *rowprep, int size)
Definition: misc_rowprep.c:855
#define SCIPdebugMsg
Definition: scip_message.h:69
SCIP_EXPR * SCIPexpriterGetCurrent(SCIP_EXPRITER *iterator)
Definition: expriter.c:673
SCIP_VAR ** x
Definition: circlepacking.c:54
SCIP_RETCODE SCIPevalExprGradient(SCIP *scip, SCIP_EXPR *expr, SCIP_SOL *sol, SCIP_Longint soltag)
Definition: scip_expr.c:1656
void SCIPinfoMessage(SCIP *scip, FILE *file, const char *formatstr,...)
Definition: scip_message.c:199
static SCIP_DECL_NLHDLREVALAUX(nlhdlrEvalAuxConvexConcave)
SCIP_Real SCIPepsilon(SCIP *scip)
SCIP_Bool SCIPisExprProduct(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1454
#define QUAD_TO_DBL(x)
Definition: dbldblarith.h:40
SCIP_RETCODE SCIPgetRowprepRowCons(SCIP *scip, SCIP_ROW **row, SCIP_ROWPREP *rowprep, SCIP_CONS *cons)
SCIP_RETCODE SCIPevalExpr(SCIP *scip, SCIP_EXPR *expr, SCIP_SOL *sol, SCIP_Longint soltag)
Definition: scip_expr.c:1623
int SCIPcompareExpr(SCIP *scip, SCIP_EXPR *expr1, SCIP_EXPR *expr2)
Definition: scip_expr.c:1723
SCIP_Bool SCIPhashmapExists(SCIP_HASHMAP *hashmap, void *origin)
Definition: misc.c:3363
SCIP_Real * SCIProwprepGetCoefs(SCIP_ROWPREP *rowprep)
Definition: misc_rowprep.c:624
SCIP_EXPR ** SCIPexprGetChildren(SCIP_EXPR *expr)
Definition: expr.c:3808
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:17920
SCIP_Real inf
Definition: intervalarith.h:46
static SCIP_RETCODE exprstackInit(SCIP *scip, EXPRSTACK *exprstack, int initsize)
#define CONVEX_NLHDLR_NAME
Definition: nlhdlr_convex.c:36
static SCIP_DECL_NLHDLRESTIMATE(nlhdlrEstimateConvex)
SCIP_Real * SCIPgetCoefsExprSum(SCIP_EXPR *expr)
Definition: expr_sum.c:1157
SCIP_Real SCIPexprGetEvalValue(SCIP_EXPR *expr)
Definition: expr.c:3872
SCIP_Real SCIPexprGetDerivative(SCIP_EXPR *expr)
Definition: expr.c:3898
int SCIPhashmapGetNEntries(SCIP_HASHMAP *hashmap)
Definition: misc.c:3481
#define CONVEX_NLHDLR_DETECTPRIORITY
Definition: nlhdlr_convex.c:38
SCIP_HASHMAPENTRY * SCIPhashmapGetEntry(SCIP_HASHMAP *hashmap, int entryidx)
Definition: misc.c:3489
SCIP_RETCODE SCIPcomputeFacetVertexPolyhedralNonlinear(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_Bool overestimate, SCIP_DECL_VERTEXPOLYFUN((*function)), void *fundata, SCIP_Real *xstar, SCIP_Real *box, int nallvars, SCIP_Real targetvalue, SCIP_Bool *success, SCIP_Real *facetcoefs, SCIP_Real *facetconstant)
const char * SCIPnlhdlrGetName(SCIP_NLHDLR *nlhdlr)
Definition: nlhdlr.c:141
SCIP_VAR * SCIPgetVarExprVar(SCIP_EXPR *expr)
Definition: expr_var.c:403
static SCIP_RETCODE estimateGradient(SCIP *scip, SCIP_NLHDLREXPRDATA *nlhdlrexprdata, SCIP_SOL *sol, SCIP_Real auxvalue, SCIP_ROWPREP *rowprep, SCIP_Bool *success)
void SCIPhashsetFree(SCIP_HASHSET **hashset, BMS_BLKMEM *blkmem)
Definition: misc.c:3730
SCIP_RETCODE SCIPincludeNlhdlrNonlinear(SCIP *scip, SCIP_NLHDLR **nlhdlr, const char *name, const char *desc, int detectpriority, int enfopriority, SCIP_DECL_NLHDLRDETECT((*detect)), SCIP_DECL_NLHDLREVALAUX((*evalaux)), SCIP_NLHDLRDATA *nlhdlrdata)
BMS_BLKMEM * SCIPblkmem(SCIP *scip)
Definition: scip_mem.c:48
public functions to work with algebraic expressions
#define QUAD(x)
Definition: dbldblarith.h:38
SCIP_Bool SCIPisExprValue(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1432
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:17251
void SCIPhashmapFree(SCIP_HASHMAP **hashmap)
Definition: misc.c:3048
void SCIPexprGetQuadraticQuadTerm(SCIP_EXPR *quadexpr, int termidx, SCIP_EXPR **expr, SCIP_Real *lincoef, SCIP_Real *sqrcoef, int *nadjbilin, int **adjbilin, SCIP_EXPR **sqrexpr)
Definition: expr.c:4102
#define NULL
Definition: lpi_spx1.cpp:155
static SCIP_RETCODE collectLeafs(SCIP *scip, SCIP_NLHDLREXPRDATA *nlhdlrexprdata)
#define REALABS(x)
Definition: def.h:201
SCIP_RETCODE SCIPreplaceExprChild(SCIP *scip, SCIP_EXPR *expr, int childidx, SCIP_EXPR *newchild)
Definition: scip_expr.c:1238
SCIP_Longint SCIPexprGetDiffTag(SCIP_EXPR *expr)
Definition: expr.c:3941
SCIP_Bool SCIPisExprSum(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1443
void SCIPfreeRowprep(SCIP *scip, SCIP_ROWPREP **rowprep)
Definition: misc_rowprep.c:558
#define SCIP_CALL(x)
Definition: def.h:384
int SCIPexpriterGetChildIdxDFS(SCIP_EXPRITER *iterator)
Definition: expriter.c:697
SCIP_VAR * h
Definition: circlepacking.c:59
SCIP_Real sup
Definition: intervalarith.h:47
static SCIP_Bool exprIsMultivarLinear(SCIP *scip, SCIP_EXPR *expr)
static SCIP_DECL_VERTEXPOLYFUN(nlhdlrExprEvalConcave)
SCIP_RETCODE SCIPhasExprCurvature(SCIP *scip, SCIP_EXPR *expr, SCIP_EXPRCURV curv, SCIP_Bool *success, SCIP_HASHMAP *assumevarfixed)
#define CONCAVE_NLHDLR_ENFOPRIORITY
Definition: nlhdlr_convex.c:44
SCIP_RETCODE SCIPaddRow(SCIP *scip, SCIP_ROW *row, SCIP_Bool forcecut, SCIP_Bool *infeasible)
Definition: scip_cut.c:241
SCIP_EXPRCURV SCIPexprGetCurvature(SCIP_EXPR *expr)
Definition: expr.c:3996
static SCIP_RETCODE nlhdlrExprGrowChildren(SCIP *scip, SCIP_HASHMAP *nlexpr2origexpr, SCIP_EXPR *nlhdlrexpr, SCIP_EXPRCURV *childrencurv)
static SCIP_DECL_NLHDLRDETECT(nlhdlrDetectConvex)
SCIP_RETCODE SCIPcreateExpriter(SCIP *scip, SCIP_EXPRITER **iterator)
Definition: scip_expr.c:2300
SCIP_Real SCIPgetCoefExprProduct(SCIP_EXPR *expr)
nonlinear handlers for convex and concave expressions, respectively
void SCIPnlhdlrSetSepa(SCIP_NLHDLR *nlhdlr, SCIP_DECL_NLHDLRINITSEPA((*initsepa)), SCIP_DECL_NLHDLRENFO((*enfo)), SCIP_DECL_NLHDLRESTIMATE((*estimate)), SCIP_DECL_NLHDLREXITSEPA((*exitsepa)))
Definition: nlhdlr.c:123
static SCIP_DECL_NLHDLREXIT(nlhdlrExitConvex)
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:115
SCIP_RETCODE SCIPsetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var, SCIP_Real val)
Definition: scip_sol.c:1212
#define SCIP_Bool
Definition: def.h:84
static SCIP_DECL_NLHDLRCOPYHDLR(nlhdlrCopyhdlrConvex)
void SCIProwprepAddConstant(SCIP_ROWPREP *rowprep, SCIP_Real constant)
Definition: misc_rowprep.c:728
SCIP_RETCODE SCIPhashmapRemoveAll(SCIP_HASHMAP *hashmap)
Definition: misc.c:3573
SCIP_EXPRCURV
Definition: type_expr.h:48
unsigned int SCIP_NLHDLR_METHOD
Definition: type_nlhdlr.h:48
SCIP_RETCODE SCIPreleaseExpr(SCIP *scip, SCIP_EXPR **expr)
Definition: scip_expr.c:1407
#define DEFAULT_CVXQUADRATIC_CONCAVE
Definition: nlhdlr_convex.c:49
constraint handler for nonlinear constraints specified by algebraic expressions
static SCIP_DECL_NLHDLRINITSEPA(nlhdlrInitSepaConvex)
#define MAX(x, y)
Definition: tclique_def.h:83
SCIP_Bool SCIPexprcurvMonomialInv(SCIP_EXPRCURV monomialcurv, int nfactors, SCIP_Real *exponents, SCIP_INTERVAL *factorbounds, SCIP_EXPRCURV *factorcurv)
Definition: exprcurv.c:447
SCIP_RETCODE SCIPappendExprChild(SCIP *scip, SCIP_EXPR *expr, SCIP_EXPR *child)
Definition: scip_expr.c:1220
SCIP_RETCODE SCIPfreeSol(SCIP *scip, SCIP_SOL **sol)
Definition: scip_sol.c:976
SCIP_Bool SCIPexprIsIntegral(SCIP_EXPR *expr)
Definition: expr.c:4017
SCIP_RETCODE SCIPcomputeExprIntegrality(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1988
SCIP_EXPR * SCIPexpriterGetNext(SCIP_EXPRITER *iterator)
Definition: expriter.c:848
SCIP_EXPR * SCIPexpriterGetChildExprDFS(SCIP_EXPRITER *iterator)
Definition: expriter.c:711
static SCIP_EXPR * exprstackPop(EXPRSTACK *exprstack)
SCIP_EXPRHDLR * SCIPexprGetHdlr(SCIP_EXPR *expr)
Definition: expr.c:3821
SCIP_Bool SCIPisExprPower(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1465
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
static void exprstackFree(SCIP *scip, EXPRSTACK *exprstack)
#define DEFAULT_CVXQUADRATIC_CONVEX
Definition: nlhdlr_convex.c:48
#define CONCAVE_NLHDLR_DETECTPRIORITY
Definition: nlhdlr_convex.c:43
void SCIPexpriterSetStagesDFS(SCIP_EXPRITER *iterator, SCIP_EXPRITER_STAGE stopstages)
Definition: expriter.c:654
void SCIPfreeExpriter(SCIP_EXPRITER **iterator)
Definition: scip_expr.c:2314
#define INITLPMAXVARVAL
Definition: nlhdlr_convex.c:54
SCIP_RETCODE SCIPreleaseRow(SCIP *scip, SCIP_ROW **row)
Definition: scip_lp.c:1553
SCIP_Bool SCIPisIntegral(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPgetHugeValue(SCIP *scip)
static const int NCURVCHECKS
#define SCIP_MAXVERTEXPOLYDIM
SCIP_Bool SCIPisExprVar(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1421
SCIP_RETCODE SCIPhashsetInsert(SCIP_HASHSET *hashset, BMS_BLKMEM *blkmem, void *element)
Definition: misc.c:3740
SCIP_RETCODE SCIPhashmapSetImage(SCIP_HASHMAP *hashmap, void *origin, void *image)
Definition: misc.c:3263
void SCIPnlhdlrSetInitExit(SCIP_NLHDLR *nlhdlr, SCIP_DECL_NLHDLRINIT((*init)), SCIP_DECL_NLHDLREXIT((*exit_)))
Definition: nlhdlr.c:97
SCIP_RETCODE SCIPaddExprsViolScoreNonlinear(SCIP *scip, SCIP_EXPR **exprs, int nexprs, SCIP_Real violscore, SCIP_SOL *sol, SCIP_Bool *success)
#define SCIP_Real
Definition: def.h:177
#define SCIP_INVALID
Definition: def.h:197
int SCIPhashmapEntryGetImageInt(SCIP_HASHMAPENTRY *entry)
Definition: misc.c:3520
SCIP_RETCODE SCIPprintRow(SCIP *scip, SCIP_ROW *row, FILE *file)
Definition: scip_lp.c:2197
SCIP_NLHDLRDATA * SCIPnlhdlrGetData(SCIP_NLHDLR *nlhdlr)
Definition: nlhdlr.c:191
#define DEFAULT_HANDLETRIVIAL
Definition: nlhdlr_convex.c:52
SCIP_BOUNDTYPE SCIPvarGetBestBoundType(SCIP_VAR *var)
Definition: var.c:18022
static SCIP_RETCODE nlhdlrExprCreate(SCIP *scip, SCIP_HASHMAP *nlexpr2origexpr, SCIP_EXPR **nlhdlrexpr, SCIP_EXPR *origexpr, SCIP_EXPRCURV curv)
static SCIP_DECL_NLHDLRFREEEXPRDATA(nlhdlrfreeExprDataConvexConcave)
struct SCIP_NlhdlrExprData SCIP_NLHDLREXPRDATA
Definition: type_nlhdlr.h:404
SCIP_RETCODE SCIPcreateRowprep(SCIP *scip, SCIP_ROWPREP **rowprep, SCIP_SIDETYPE sidetype, SCIP_Bool local)
Definition: misc_rowprep.c:538
struct SCIP_NlhdlrData SCIP_NLHDLRDATA
Definition: type_nlhdlr.h:403
SCIP_Bool SCIPisZero(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPvarGetUbLocal(SCIP_VAR *var)
Definition: var.c:17976
SCIP_RETCODE SCIPaddRowprepTerm(SCIP *scip, SCIP_ROWPREP *rowprep, SCIP_VAR *var, SCIP_Real coef)
Definition: misc_rowprep.c:881
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:102
#define DEFAULT_DETECTSUM
Definition: nlhdlr_convex.c:46
SCIP_EXPRCURV SCIPexprcurvMultiply(SCIP_Real factor, SCIP_EXPRCURV curvature)
Definition: exprcurv.c:78
#define DEFAULT_CVXSIGNOMIAL
Definition: nlhdlr_convex.c:50
static SCIP_RETCODE exprstackPush(SCIP *scip, EXPRSTACK *exprstack, int nexprs, SCIP_EXPR **exprs)
#define SCIP_EXPRITER_VISITINGCHILD
Definition: type_expr.h:668
SCIP_VAR * SCIPgetExprAuxVarNonlinear(SCIP_EXPR *expr)
SCIP_RETCODE SCIPhashmapInsert(SCIP_HASHMAP *hashmap, void *origin, void *image)
Definition: misc.c:3096
static SCIP_RETCODE createNlhdlrExprData(SCIP *scip, SCIP_NLHDLRDATA *nlhdlrdata, SCIP_NLHDLREXPRDATA **nlhdlrexprdata, SCIP_EXPR *expr, SCIP_EXPR *nlexpr, SCIP_HASHMAP *nlexpr2origexpr, int nleafs, SCIP_NLHDLR_METHOD participating)
static SCIP_RETCODE constructExpr(SCIP *scip, SCIP_NLHDLRDATA *nlhdlrdata, SCIP_EXPR **rootnlexpr, SCIP_HASHMAP *nlexpr2origexpr, int *nleafs, SCIP_EXPR *rootexpr, SCIP_EXPRCURV curv, SCIP_HASHMAP *assumevarfixed, SCIP_Bool assumecurvature, SCIP_Bool *curvsuccess)
#define SCIPallocClearBlockMemory(scip, ptr)
Definition: scip_mem.h:82
SCIP_Bool SCIPexpriterIsEnd(SCIP_EXPRITER *iterator)
Definition: expriter.c:959
SCIPallocBlockMemory(scip, subsol))
void SCIProwprepSetLocal(SCIP_ROWPREP *rowprep, SCIP_Bool islocal)
Definition: misc_rowprep.c:748
public functions of nonlinear handlers of nonlinear constraints
#define SCIP_CALL_ABORT(x)
Definition: def.h:363
SCIP_Real SCIPceil(SCIP *scip, SCIP_Real val)
SCIP_RETCODE SCIPincludeNlhdlrConcave(SCIP *scip)
SCIP_Longint SCIPgetNLPs(SCIP *scip)
SCIP_Real SCIPround(SCIP *scip, SCIP_Real val)
#define SCIP_NLHDLR_METHOD_SEPABOTH
Definition: type_nlhdlr.h:44
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1352
#define SCIP_NLHDLR_METHOD_SEPAABOVE
Definition: type_nlhdlr.h:43
static SCIP_DECL_NLHDLRFREEHDLRDATA(nlhdlrfreeHdlrDataConvexConcave)
int SCIPsolGetIndex(SCIP_SOL *sol)
Definition: sol.c:2635
SCIP_Real SCIPfloor(SCIP *scip, SCIP_Real val)
preparation of a linear inequality to become a SCIP_ROW
SCIP_RETCODE SCIPremoveExprChildren(SCIP *scip, SCIP_EXPR *expr)
Definition: scip_expr.c:1257
SCIP_RETCODE SCIPsetPtrarrayVal(SCIP *scip, SCIP_PTRARRAY *ptrarray, int idx, void *val)
SCIP_RETCODE SCIPaddBoolParam(SCIP *scip, const char *name, const char *desc, SCIP_Bool *valueptr, SCIP_Bool isadvanced, SCIP_Bool defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:48
SCIP_RETCODE SCIPcreateSol(SCIP *scip, SCIP_SOL **sol, SCIP_HEUR *heur)
Definition: scip_sol.c:319
#define SCIPreallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:119
void SCIPnlhdlrSetCopyHdlr(SCIP_NLHDLR *nlhdlr, SCIP_DECL_NLHDLRCOPYHDLR((*copy)))
Definition: nlhdlr.c:63