Scippy

SCIP

Solving Constraint Integer Programs

cons_quadratic.h File Reference

Detailed Description

constraint handler for quadratic constraints $\textrm{lhs} \leq \sum_{i,j=1}^n a_{i,j} x_ix_j + \sum_{i=1}^n b_i x_i \leq \textrm{rhs}$

Author
Stefan Vigerske

This constraint handler handles constraints of the form

\[ \textrm{lhs} \leq \sum_{i,j=1}^n a_{i,j} x_ix_j + \sum_{i=1}^n b_i x_i \leq \textrm{rhs} \]

Constraints are enforced by separation, domain propagation, and spatial branching.

For semidefinite matrices $A=(a_{i,j})_{i,j}$, cuts based on linearization of $\langle x, Ax\rangle$ are implemented. For underestimating a non-convex term, McCormick underestimators and secants for univariate concave quadratic terms are implemented. If $\langle x, Ax\rangle$ is factorable (i.e., can be written as product of two linear functions), specialized separation techniques (e.g., lifted tangent inequalities) that take the constraint sides into account are applied.

Branching is performed for variables in nonconvex terms, if the relaxation solution cannot be separated. Further, domain propagation is applied.

During presolve, variable products which contain binary variables may be reformulated into linear constraints, thereby introducing new variables.

See also

Timo Berthold and Stefan Heinz and Stefan Vigerske
Extending a CIP framework to solve MIQCPs
In: Jon Lee and Sven Leyffer (eds.), Mixed-integer nonlinear optimization: Algorithmic advances and applications, IMA volumes in Mathematics and its Applications, volume 154, 427-444, 2012.
Stefan Vigerske
Decomposition of Multistage Stochastic Programs and a Constraint Integer Programming Approach to Mixed-Integer Nonlinear Programming
PhD Thesis, Humboldt-University Berlin, 2012, submitted.
Pietro Belotti and Andrew J. Miller and Mahdi Namazifar
Linear inequalities for bounded products of variables
SIAG/OPT Views-and-News 22:1, 1-8, 2011.

Definition in file cons_quadratic.h.

#include "scip/scip.h"
#include "scip/intervalarith.h"
#include "nlpi/type_nlpi.h"

Go to the source code of this file.

Macros

#define SCIP_DECL_QUADCONSUPGD(x)
 

Typedefs

typedef struct
SCIP_QuadVarEventData 
SCIP_QUADVAREVENTDATA
 
typedef struct SCIP_QuadVarTerm SCIP_QUADVARTERM
 
typedef struct SCIP_BilinTerm SCIP_BILINTERM
 

Functions

SCIP_RETCODE SCIPincludeConshdlrQuadratic (SCIP *scip)
 
SCIP_RETCODE SCIPincludeQuadconsUpgrade (SCIP *scip, SCIP_DECL_QUADCONSUPGD((*quadconsupgd)), int priority, SCIP_Bool active, const char *conshdlrname)
 
SCIP_RETCODE SCIPcreateConsQuadratic (SCIP *scip, SCIP_CONS **cons, const char *name, int nlinvars, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nquadterms, SCIP_VAR **quadvars1, SCIP_VAR **quadvars2, SCIP_Real *quadcoeffs, SCIP_Real lhs, SCIP_Real rhs, SCIP_Bool initial, SCIP_Bool separate, SCIP_Bool enforce, SCIP_Bool check, SCIP_Bool propagate, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool dynamic, SCIP_Bool removable)
 
SCIP_RETCODE SCIPcreateConsBasicQuadratic (SCIP *scip, SCIP_CONS **cons, const char *name, int nlinvars, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nquadterms, SCIP_VAR **quadvars1, SCIP_VAR **quadvars2, SCIP_Real *quadcoefs, SCIP_Real lhs, SCIP_Real rhs)
 
SCIP_RETCODE SCIPcreateConsQuadratic2 (SCIP *scip, SCIP_CONS **cons, const char *name, int nlinvars, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nquadvarterms, SCIP_QUADVARTERM *quadvarterms, int nbilinterms, SCIP_BILINTERM *bilinterms, SCIP_Real lhs, SCIP_Real rhs, SCIP_Bool initial, SCIP_Bool separate, SCIP_Bool enforce, SCIP_Bool check, SCIP_Bool propagate, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool dynamic, SCIP_Bool removable)
 
SCIP_RETCODE SCIPcreateConsBasicQuadratic2 (SCIP *scip, SCIP_CONS **cons, const char *name, int nlinvars, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nquadvarterms, SCIP_QUADVARTERM *quadvarterms, int nbilinterms, SCIP_BILINTERM *bilinterms, SCIP_Real lhs, SCIP_Real rhs)
 
void SCIPaddConstantQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_Real constant)
 
SCIP_RETCODE SCIPaddLinearVarQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real coef)
 
SCIP_RETCODE SCIPaddQuadVarQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real lincoef, SCIP_Real sqrcoef)
 
SCIP_RETCODE SCIPaddQuadVarLinearCoefQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real coef)
 
SCIP_RETCODE SCIPaddSquareCoefQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real coef)
 
SCIP_RETCODE SCIPaddBilinTermQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var1, SCIP_VAR *var2, SCIP_Real coef)
 
SCIP_RETCODE SCIPgetNlRowQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_NLROW **nlrow)
 
int SCIPgetNLinearVarsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_VAR ** SCIPgetLinearVarsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_RealSCIPgetCoefsLinearVarsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
int SCIPgetNQuadVarTermsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_QUADVARTERMSCIPgetQuadVarTermsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_RETCODE SCIPsortQuadVarTermsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_RETCODE SCIPfindQuadVarTermQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, int *pos)
 
int SCIPgetNBilinTermsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_BILINTERMSCIPgetBilinTermsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_Real SCIPgetLhsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_Real SCIPgetRhsQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_RETCODE SCIPcheckCurvatureQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_Bool SCIPisConvexQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_Bool SCIPisConcaveQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_RETCODE SCIPgetViolationQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_SOL *sol, SCIP_Real *violation)
 
SCIP_Bool SCIPisLinearLocalQuadratic (SCIP *scip, SCIP_CONS *cons)
 
SCIP_RETCODE SCIPaddToNlpiProblemQuadratic (SCIP *scip, SCIP_CONS *cons, SCIP_NLPI *nlpi, SCIP_NLPIPROBLEM *nlpiprob, SCIP_HASHMAP *scipvar2nlpivar, SCIP_Bool names)
 

Macro Definition Documentation

#define SCIP_DECL_QUADCONSUPGD (   x)
Value:
SCIP_RETCODE x (SCIP* scip, SCIP_CONS* cons, \
int nbinlin, int nbinquad, int nintlin, int nintquad, int nimpllin, int nimplquad, int ncontlin, int ncontquad, \
SCIP_Bool integral, int* nupgdconss, SCIP_CONS** upgdconss, int upgdconsssize)

upgrading method for quadratic constraints into more specific constraints

the method might upgrade a quadratic constraint into a set of quadratic constraints the caller provided an array upgdconss to store upgrade constraints the length of upgdconss is given by upgdconsssize if an upgrade is not possible, set *nupgdconss to zero if more than upgdconsssize many constraints shall replace cons, the function should return the required number as negated value in *nupgdconss i.e., if cons should be replaced by 3 constraints, the function should set *nupgdconss to -3 and return with SCIP_OKAY

input:

  • scip : SCIP main data structure
  • cons : the quadratic constraint to upgrade
  • nbinlin : number of binary variables in linear part
  • nbinquad : number of binary variables in quadratic part
  • nintlin : number of integer variables in linear part
  • nintquad : number of integer variables in quadratic part
  • nimpllin : number of implicit integer variables in linear part
  • nimplquad : number of implicit integer variables in quadratic part
  • ncontlin : number of continuous variables in linear part
  • ncontquad : number of continuous variables in quadratic part
  • integral : TRUE iff constraints activity value is always integral
  • nupgdconss : pointer to store number of constraints that replace this constraint
  • upgdconss : array to store constraints that replace this constraint
  • upgdconsssize : length of the provided upgdconss array

Definition at line 128 of file cons_quadratic.h.

Typedef Documentation

typedef struct SCIP_QuadVarEventData SCIP_QUADVAREVENTDATA

event data for variable bound changes in quadratic constraints

Definition at line 71 of file cons_quadratic.h.

typedef struct SCIP_QuadVarTerm SCIP_QUADVARTERM

Definition at line 87 of file cons_quadratic.h.

typedef struct SCIP_BilinTerm SCIP_BILINTERM

Definition at line 98 of file cons_quadratic.h.

Function Documentation

SCIP_RETCODE SCIPincludeConshdlrQuadratic ( SCIP scip)

creates the handler for quadratic constraints and includes it in SCIP

Parameters
scipSCIP data structure
SCIP_RETCODE SCIPincludeQuadconsUpgrade ( SCIP scip,
SCIP_DECL_QUADCONSUPGD((*quadconsupgd))  ,
int  priority,
SCIP_Bool  active,
const char *  conshdlrname 
)

includes a quadratic constraint upgrade method into the quadratic constraint handler

Parameters
scipSCIP data structure
prioritypriority of upgrading method
activeshould the upgrading method be active by default?
conshdlrnamename of the constraint handler
SCIP_RETCODE SCIPcreateConsQuadratic ( SCIP scip,
SCIP_CONS **  cons,
const char *  name,
int  nlinvars,
SCIP_VAR **  linvars,
SCIP_Real lincoefs,
int  nquadterms,
SCIP_VAR **  quadvars1,
SCIP_VAR **  quadvars2,
SCIP_Real quadcoeffs,
SCIP_Real  lhs,
SCIP_Real  rhs,
SCIP_Bool  initial,
SCIP_Bool  separate,
SCIP_Bool  enforce,
SCIP_Bool  check,
SCIP_Bool  propagate,
SCIP_Bool  local,
SCIP_Bool  modifiable,
SCIP_Bool  dynamic,
SCIP_Bool  removable 
)

Creates and captures a quadratic constraint.

The constraint should be given in the form

\[ \ell \leq \sum_{i=1}^n b_i x_i + \sum_{j=1}^m a_j y_j z_j \leq u, \]

where $x_i = y_j = z_k$ is possible.

Note
the constraint gets captured, hence at one point you have to release it using the method SCIPreleaseCons()
Parameters
scipSCIP data structure
conspointer to hold the created constraint
namename of constraint
nlinvarsnumber of linear terms (n)
linvarsvariables in linear part (x_i) or NULL if nlinvars == 0
lincoefscoefficients of variables in linear part (b_i) or NULL if nlinvars == 0
nquadtermsnumber of quadratic terms (m)
quadvars1array with first variables in quadratic terms (y_j) or NULL if nquadterms == 0
quadvars2array with second variables in quadratic terms (z_j) or NULL if nquadterms == 0
quadcoeffsarray with coefficients of quadratic terms (a_j) or NULL if nquadterms == 0
lhsleft hand side of quadratic equation (l)
rhsright hand side of quadratic equation (u)
initialshould the LP relaxation of constraint be in the initial LP? Usually set to TRUE. Set to FALSE for 'lazy constraints'.
separateshould the constraint be separated during LP processing? Usually set to TRUE.
enforceshould the constraint be enforced during node processing? TRUE for model constraints, FALSE for additional, redundant constraints.
checkshould the constraint be checked for feasibility? TRUE for model constraints, FALSE for additional, redundant constraints.
propagateshould the constraint be propagated during node processing? Usually set to TRUE.
localis constraint only valid locally? Usually set to FALSE. Has to be set to TRUE, e.g., for branching constraints.
modifiableis constraint modifiable (subject to column generation)? Usually set to FALSE. In column generation applications, set to TRUE if pricing adds coefficients to this constraint.
dynamicis constraint subject to aging? Usually set to FALSE. Set to TRUE for own cuts which are separated as constraints.
removableshould the relaxation be removed from the LP due to aging or cleanup? Usually set to FALSE. Set to TRUE for 'lazy constraints' and 'user cuts'.
SCIP_RETCODE SCIPcreateConsBasicQuadratic ( SCIP scip,
SCIP_CONS **  cons,
const char *  name,
int  nlinvars,
SCIP_VAR **  linvars,
SCIP_Real lincoefs,
int  nquadterms,
SCIP_VAR **  quadvars1,
SCIP_VAR **  quadvars2,
SCIP_Real quadcoefs,
SCIP_Real  lhs,
SCIP_Real  rhs 
)

creates and captures a quadratic constraint in its most basic variant, i. e., with all constraint flags set to their default values, which can be set afterwards using SCIPsetConsFLAGNAME() in scip.h

The constraint should be given in the form

\[ \ell \leq \sum_{i=1}^n b_i x_i + \sum_{j=1}^m a_j y_jz_j \leq u, \]

where $x_i = y_j = z_k$ is possible.

See Also
SCIPcreateConsQuadratic() for the default constraint flag configuration
Note
the constraint gets captured, hence at one point you have to release it using the method SCIPreleaseCons()
Parameters
scipSCIP data structure
conspointer to hold the created constraint
namename of constraint
nlinvarsnumber of linear terms (n)
linvarsarray with variables in linear part (x_i)
lincoefsarray with coefficients of variables in linear part (b_i)
nquadtermsnumber of quadratic terms (m)
quadvars1array with first variables in quadratic terms (y_j)
quadvars2array with second variables in quadratic terms (z_j)
quadcoefsarray with coefficients of quadratic terms (a_j)
lhsleft hand side of quadratic equation (ell)
rhsright hand side of quadratic equation (u)
SCIP_RETCODE SCIPcreateConsQuadratic2 ( SCIP scip,
SCIP_CONS **  cons,
const char *  name,
int  nlinvars,
SCIP_VAR **  linvars,
SCIP_Real lincoefs,
int  nquadvarterms,
SCIP_QUADVARTERM quadvarterms,
int  nbilinterms,
SCIP_BILINTERM bilinterms,
SCIP_Real  lhs,
SCIP_Real  rhs,
SCIP_Bool  initial,
SCIP_Bool  separate,
SCIP_Bool  enforce,
SCIP_Bool  check,
SCIP_Bool  propagate,
SCIP_Bool  local,
SCIP_Bool  modifiable,
SCIP_Bool  dynamic,
SCIP_Bool  removable 
)

creates and captures a quadratic constraint.

The constraint should be given in the form

\[ \ell \leq \sum_{i=1}^n b_i x_i + \sum_{j=1}^m (a_j y_j^2 + b_j y_j) + \sum_{k=1}^p c_k v_k w_k \leq u. \]

Note
the constraint gets captured, hence at one point you have to release it using the method SCIPreleaseCons()
Parameters
scipSCIP data structure
conspointer to hold the created constraint
namename of constraint
nlinvarsnumber of linear terms (n)
linvarsarray with variables in linear part (x_i)
lincoefsarray with coefficients of variables in linear part (b_i)
nquadvartermsnumber of quadratic terms (m)
quadvartermsquadratic variable terms
nbilintermsnumber of bilinear terms (p)
bilintermsbilinear terms
lhsconstraint left hand side (ell)
rhsconstraint right hand side (u)
initialshould the LP relaxation of constraint be in the initial LP?
separateshould the constraint be separated during LP processing?
enforceshould the constraint be enforced during node processing?
checkshould the constraint be checked for feasibility?
propagateshould the constraint be propagated during node processing?
localis constraint only valid locally?
modifiableis constraint modifiable (subject to column generation)?
dynamicis constraint dynamic?
removableshould the constraint be removed from the LP due to aging or cleanup?
SCIP_RETCODE SCIPcreateConsBasicQuadratic2 ( SCIP scip,
SCIP_CONS **  cons,
const char *  name,
int  nlinvars,
SCIP_VAR **  linvars,
SCIP_Real lincoefs,
int  nquadvarterms,
SCIP_QUADVARTERM quadvarterms,
int  nbilinterms,
SCIP_BILINTERM bilinterms,
SCIP_Real  lhs,
SCIP_Real  rhs 
)

creates and captures a quadratic constraint in its most basic variant, i. e., with all constraint flags set to their default values, which can be set afterwards using SCIPsetConsFLAGNAME() in scip.h

The constraint should be given in the form

\[ \ell \leq \sum_{i=1}^n b_i x_i + \sum_{j=1}^m (a_j y_j^2 + b_j y_j) + \sum_{k=1}^p c_kv_kw_k \leq u. \]

See Also
SCIPcreateConsQuadratic2() for the default constraint flag configuration
Note
the constraint gets captured, hence at one point you have to release it using the method SCIPreleaseCons()
Parameters
scipSCIP data structure
conspointer to hold the created constraint
namename of constraint
nlinvarsnumber of linear terms (n)
linvarsarray with variables in linear part (x_i)
lincoefsarray with coefficients of variables in linear part (b_i)
nquadvartermsnumber of quadratic terms (m)
quadvartermsquadratic variable terms
nbilintermsnumber of bilinear terms (p)
bilintermsbilinear terms
lhsconstraint left hand side (ell)
rhsconstraint right hand side (u)
void SCIPaddConstantQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_Real  constant 
)

Adds a constant to the constraint function, that is, subtracts a constant from both sides

Parameters
scipSCIP data structure
consconstraint
constantconstant to subtract from both sides
SCIP_RETCODE SCIPaddLinearVarQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_VAR var,
SCIP_Real  coef 
)

Adds a linear variable with coefficient to a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
varvariable
coefcoefficient of variable
SCIP_RETCODE SCIPaddQuadVarQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_VAR var,
SCIP_Real  lincoef,
SCIP_Real  sqrcoef 
)

Adds a quadratic variable with linear and square coefficient to a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
varvariable
lincoeflinear coefficient of variable
sqrcoefsquare coefficient of variable
SCIP_RETCODE SCIPaddQuadVarLinearCoefQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_VAR var,
SCIP_Real  coef 
)

Adds a linear coefficient for a quadratic variable.

Variable will be added with square coefficient 0.0 if not existing yet.

Parameters
scipSCIP data structure
consconstraint
varvariable
coefvalue to add to linear coefficient of variable
SCIP_RETCODE SCIPaddSquareCoefQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_VAR var,
SCIP_Real  coef 
)

Adds a square coefficient for a quadratic variable.

Variable will be added with linear coefficient 0.0 if not existing yet.

Parameters
scipSCIP data structure
consconstraint
varvariable
coefvalue to add to square coefficient of variable
SCIP_RETCODE SCIPaddBilinTermQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_VAR var1,
SCIP_VAR var2,
SCIP_Real  coef 
)

Adds a bilinear term to a quadratic constraint.

Variables will be added with linear and square coefficient 0.0 if not existing yet. If variables are equal, only the square coefficient of the variable is updated.

Parameters
scipSCIP data structure
consconstraint
var1first variable
var2second variable
coefcoefficient of bilinear term
SCIP_RETCODE SCIPgetNlRowQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_NLROW **  nlrow 
)

Gets the quadratic constraint as a nonlinear row representation.

Parameters
scipSCIP data structure
consconstraint
nlrowpointer to store nonlinear row
int SCIPgetNLinearVarsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the number of variables in the linear part of a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
SCIP_VAR** SCIPgetLinearVarsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the variables in the linear part of a quadratic constraint. Length is given by SCIPgetNLinearVarsQuadratic.

Parameters
scipSCIP data structure
consconstraint
SCIP_Real* SCIPgetCoefsLinearVarsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the coefficients in the linear part of a quadratic constraint. Length is given by SCIPgetNQuadVarsQuadratic.

Parameters
scipSCIP data structure
consconstraint
int SCIPgetNQuadVarTermsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the number of quadratic variable terms of a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
SCIP_QUADVARTERM* SCIPgetQuadVarTermsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the quadratic variable terms of a quadratic constraint. Length is given by SCIPgetNQuadVarTermsQuadratic.

Parameters
scipSCIP data structure
consconstraint
SCIP_RETCODE SCIPsortQuadVarTermsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Ensures that quadratic variable terms are sorted.

Parameters
scipSCIP data structure
consconstraint
SCIP_RETCODE SCIPfindQuadVarTermQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_VAR var,
int *  pos 
)

Finds the position of a quadratic variable term for a given variable.

Note
If the quadratic variable terms have not been sorted before, then a search may reorder the current order of the terms.
Parameters
scipSCIP data structure
consconstraint
varvariable to search for
posbuffer to store position of quadvarterm for var, or -1 if not found
int SCIPgetNBilinTermsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the number of bilinear terms of a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
SCIP_BILINTERM* SCIPgetBilinTermsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the bilinear terms of a quadratic constraint. Length is given by SCIPgetNBilinTermQuadratic.

Parameters
scipSCIP data structure
consconstraint
SCIP_Real SCIPgetLhsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the left hand side of a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
SCIP_Real SCIPgetRhsQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Gets the right hand side of a quadratic constraint.

Parameters
scipSCIP data structure
consconstraint
SCIP_RETCODE SCIPcheckCurvatureQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Check the quadratic function of a quadratic constraint for its semi-definiteness, if not done yet.

Parameters
scipSCIP data structure
consconstraint
SCIP_Bool SCIPisConvexQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Indicates whether the quadratic function of a quadratic constraint is (known to be) convex.

Parameters
scipSCIP data structure
consconstraint
SCIP_Bool SCIPisConcaveQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Indicates whether the quadratic function of a quadratic constraint is (known to be) concave.

Parameters
scipSCIP data structure
consconstraint
SCIP_RETCODE SCIPgetViolationQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_SOL sol,
SCIP_Real violation 
)

Gets the violation of a constraint by a solution.

Parameters
scipSCIP data structure
consconstraint
solsolution which violation to calculate, or NULL for LP solution
violationpointer to store violation of constraint
SCIP_Bool SCIPisLinearLocalQuadratic ( SCIP scip,
SCIP_CONS cons 
)

Indicates whether the quadratic constraint is local w.r.t. the current local bounds.

That is, checks whether each variable with a square term is fixed and for each bilinear term at least one variable is fixed.

Parameters
scipSCIP data structure
consconstraint
SCIP_RETCODE SCIPaddToNlpiProblemQuadratic ( SCIP scip,
SCIP_CONS cons,
SCIP_NLPI nlpi,
SCIP_NLPIPROBLEM nlpiprob,
SCIP_HASHMAP scipvar2nlpivar,
SCIP_Bool  names 
)

Adds the constraint to an NLPI problem.

Parameters
scipSCIP data structure
consconstraint
nlpiinterface to NLP solver
nlpiprobNLPI problem where to add constraint
scipvar2nlpivarmapping from SCIP variables to variable indices in NLPI
nameswhether to pass constraint names to NLPI