sepa_eccuts.c
Go to the documentation of this file.
32 /**@todo only add nonlinear row aggregations where at least ...% of the variables (bilinear terms?) are in edge concave
34 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
51 #define SEPA_DELAY FALSE /**< should separation method be delayed, if other separators found cuts? */
53 #define CLIQUE_MAXFIRSTNODEWEIGHT 1000 /**< maximum weight of branching nodes in level 0; 0 if not used for cliques
57 #define CLIQUE_BACKTRACKFREQ 10000 /**< frequency to backtrack to first level of tree (0: no premature backtracking) */
59 #define DEFAULT_DYNAMICCUTS TRUE /**< should generated cuts be removed from the LP if they are no longer tight? */
60 #define DEFAULT_MAXROUNDS 10 /**< maximal number of separation rounds per node (-1: unlimited) */
61 #define DEFAULT_MAXROUNDSROOT 250 /**< maximal number of separation rounds in the root node (-1: unlimited) */
63 #define DEFAULT_MAXSEPACUTS 10 /**< maximal number of e.c. cuts separated per separation round */
64 #define DEFAULT_MAXSEPACUTSROOT 50 /**< maximal number of e.c. cuts separated per separation round in root node */
65 #define DEFAULT_CUTMAXRANGE 1e+7 /**< maximal coefficient range of a cut (maximal coefficient divided by minimal
68 #define DEFAULT_MINAGGRSIZE 3 /**< search for e.c. aggregation of at least this size (has to be >= 3) */
69 #define DEFAULT_MAXAGGRSIZE 4 /**< search for e.c. aggregation of at most this size (has to be >= minaggrsize) */
70 #define DEFAULT_MAXBILINTERMS 500 /**< maximum number of bilinear terms allowed to be in a quadratic constraint */
71 #define DEFAULT_MAXSTALLROUNDS 5 /**< maximum number of unsuccessful rounds in the e.c. aggregation search */
75 #define ADJUSTFACETTOL 1e-6 /**< adjust resulting facets in checkRikun() up to a violation of this value */
79 static const int poweroftwo[] = { 1, 2, 4, 8, 16, 32, 64, 128, 256, 512, 1024, 2048, 4096, 8192 };
85 /** data to store a single edge-concave aggregations; an edge-concave aggregation of a quadratic constraint is a subset
99 };
102 /** data to store all edge-concave aggregations and the remaining part of a nonlinear row of the form g(x) <= rhs */
131 };
141 int minaggrsize; /**< only search for e.c. aggregations of at least this size (has to be >= 3) */
142 int maxaggrsize; /**< only search for e.c. aggregations of at most this size (has to be >= minaggrsize) */
144 int maxbilinterms; /**< maximum number of bilinear terms allowed to be in a quadratic constraint */
145 int maxstallrounds; /**< maximum number of unsuccessful rounds in the e.c. aggregation search */
147 SCIP_LPI* lpi; /**< LP interface to solve the LPs to compute the facets of the convex envelopes */
152 SCIP_Bool dynamiccuts; /**< should generated cuts be removed from the LP if they are no longer tight? */
159 int maxsepacutsroot; /**< maximal number of e.c. cuts separated per separation round in root node */
163 int nlhsnlrowaggrs; /**< number of found nonlinear row aggregations for SCIP_NLROWs of the form g(x) <= rhs */
164 int nrhsnlrowaggrs; /**< number of found nonlinear row aggregations for SCIP_NLROWs of the form g(x) >= lhs */
303 SCIPdebugMsgPrint(scip, "%e %s * %s ", ecaggr->termcoefs[i], SCIPvarGetName(x), SCIPvarGetName(y) );
366 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &nlrowaggr->linvars, nlrowaggr->linvarssize, newsize) );
367 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &nlrowaggr->lincoefs, nlrowaggr->linvarssize, newsize) );
395 SCIP_CALL( SCIPensureBlockMemoryArray(scip, &nlrowaggr->quadvars, &nlrowaggr->quadvarssize, nlrowaggr->nquadvars+1) );
479 SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(*nlrowaggr)->quadvar2aggr, quadvar2aggr, nquadvars) );
482 SCIP_CALL( nlrowaggrStoreLinearTerms(scip, *nlrowaggr, SCIPnlrowGetLinearVars(nlrow), SCIPnlrowGetLinearCoefs(nlrow),
529 /* only handle qterm1 == qterm here; the other case will be handled when its turn for qterm2 to be qterm */
620 /* only handle qterm1 == qterm here; the other case will be handled when its turn for qterm2 to be qterm */
711 SCIPdebugMsgPrint(scip, "%e %s * %s + ", nlrowaggr->remtermcoefs[i], SCIPvarGetName(nlrowaggr->remtermvars1[i]),
714 SCIPdebugMsgPrint(scip, "%e %s + ", nlrowaggr->lincoefs[i], SCIPvarGetName(nlrowaggr->linvars[i]) );
819 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &sepadata->nlrowaggrs, sepadata->nlrowaggrssize, 2 * sepadata->nlrowaggrssize) ); /*lint !e506 !e647*/
861 /** creates an MIP to search for cycles with an odd number of positive edges in the graph representation of a nonlinear row
865 * If the term is quadratic (e.g., loop in the graph) we fix the corresponding variables to zero.
880 int* narcs /**< pointer to store the number of created arc variables (number of square and bilinear terms) */
938 SCIP_CALL( SCIPcreateVarBasic(subscip, &forwardarcs[*narcs], name, 0.0, 0.0, 0.01, SCIP_VARTYPE_BINARY) );
941 SCIP_CALL( SCIPcreateVarBasic(subscip, &backwardarcs[*narcs], name, 0.0, 0.0, 0.01, SCIP_VARTYPE_BINARY) );
967 SCIP_CALL( SCIPcreateVarBasic(subscip, &forwardarcs[*narcs], name, 0.0, 1.0, 0.01 + edgeweight, SCIP_VARTYPE_BINARY) );
971 SCIP_CALL( SCIPcreateVarBasic(subscip, &backwardarcs[*narcs], name, 0.0, 1.0, 0.01 + edgeweight, SCIP_VARTYPE_BINARY) );
993 SCIP_CALL( SCIPcreateConsBasicLinear(subscip, &noparallelcons, name, 0, NULL, NULL, 0.0, 1.0) );
1006 SCIP_CALL( SCIPcreateConsBasicXor(subscip, &oddcyclecons, name, TRUE, noddcyclearcs, oddcyclearcs) );
1081 /** fixed all arc variables (u,v) for which u or v is already in an edge-concave aggregation */
1308 * SCIP's clique code can only handle integer node weights; all node weights are scaled by a factor of 100; since the
1342 /* SCIPdebugMsg(scip, "%d (%s): nodeweight %d \n", i, SCIPvarGetName(SCIPnlrowGetQuadVars(nlrow)[i]), nodeweight); */
1397 /** searches for edge-concave aggregations by computing cliques in the graph representation of a given nonlinear row
1399 * update graph, compute clique, store clique; after computing a clique we heuristically check if the clique contains
1455 tcliqueMaxClique(tcliqueGetNNodes, tcliqueGetWeights, tcliqueIsEdge, tcliqueSelectAdjnodes, graph, NULL, NULL,
1456 maxcliquenodes, &nmaxcliquenodes, &maxcliqueweight, CLIQUE_MAXFIRSTNODEWEIGHT, CLIQUE_MINWEIGHT,
1468 SCIP_CALL( SCIPhashmapInsertInt(cliquemap, (void*) (size_t) maxcliquenodes[i], 0) ); /*lint !e571*/
1489 isoddcycleedge = (rhsaggr && SCIPisPositive(scip, coef)) || (!rhsaggr && SCIPisNegative(scip, coef));
1509 isoddcycleedge = (rhsaggr && SCIPisPositive(scip, coef)) || (!rhsaggr && SCIPisNegative(scip, coef));
1528 if( noddcycleedges == 0 || (aggrsize == 3 && noddcycleedges == 2) || (aggrsize == 4 && ntwodegrees == 4) )
1533 /* add the found clique as an edge-concave aggregation or exclude the nodes from the remaining search */
1537 SCIPdebugMsg(scip, "%s %d\n", *foundaggr ? "aggregate node: " : "exclude node: ", maxcliquenodes[i]);
1558 SCIP_Bool rhsaggr, /**< consider nonlinear row aggregation for g(x) <= rhs (TRUE) or g(x) >= lhs (FALSE) */
1589 /* arrays to store all arc variables of the MIP model; note that we introduce variables even for loops in the graph
1597 /* initialize mapping from quadvars to e.c. aggregation index (-1: quadvar is in no aggregation); compute node
1610 nodeweights[i] = phi(scip, SCIPgetSolVal(scip, sol, var), SCIPvarGetLbLocal(var), SCIPvarGetUbLocal(var));
1611 SCIPdebugMsg(scip, "%s = %e (%e in [%e, %e])\n", SCIPvarGetName(var), nodeweights[i], SCIPgetSolVal(scip, sol, var),
1615 SCIP_CALL( createMIP(scip, subscip, sepadata, nlrow, rhsaggr, forwardarcs, backwardarcs, nodeweights, &nedges, &narcs) );
1658 SCIP_CALL( storeAggrFromMIP(subscip, nlrow, forwardarcs, backwardarcs, quadvar2aggr, *nfound) );
1676 /* 2.a - search and store a single edge-concave aggregation by computing a clique with a good cycle */
1677 SCIP_CALL_FINALLY( searchEcAggrWithCliques(scip, graph, sepadata, nlrow, quadvar2aggr, *nfound, rhsaggr,
1706 SCIP_CALL_FINALLY( updateMIP(subscip, nlrow, forwardarcs, backwardarcs, quadvar2aggr, &nedges), tcliqueFree(&graph) );
1737 /** computes a partitioning into edge-concave aggregations for a given (quadratic) nonlinear row
1739 * Each aggregation has to contain a cycle with an odd number of positive weighted edges (good cycles) in the corresponding graph representation.
1741 * -# use a MIP model based on binary flow variables to compute good cycles and store the implied subgraphs as an e.c. aggr.
1742 * -# if we find a good cycle, store the implied subgraph, delete it from the graph representation and go to 1)
1745 * -# if we find a good cycle in C, store the implied subgraph of C, delete it from the graph representation and go to 1)
1754 SCIP_Bool rhsaggr, /**< consider nonlinear row aggregation for g(x) <= rhs (TRUE) or g(x) >= lhs (FALSE) */
1775 /** returns whether a given nonlinear row can be used to compute edge-concave aggregations for which their convex
1779 * an odd number of positive edges in the corresponding graph representation of the nonlinear row.
1786 SCIP_Bool* rhscandidate, /**< pointer to store if we should compute edge-concave aggregations for
1788 SCIP_Bool* lhscandidate /**< pointer to store if we should compute edge-concave aggregations for
1807 /* check whether nlrow is in the NLP, is quadratic in variables, and there are enough quadratic variables */
1894 if( nposedges + nnegedges > sepadata->maxbilinterms || nposedges + nnegedges < sepadata->minaggrsize
1902 /* check if there are enough positive/negative edges; for a 3-clique there has to be an odd number of those edges */
1934 /* check obvious conditions for existing cycles with an odd number of positive/negative edges */
1941 SCIPexprGetQuadraticData(SCIPnlrowGetExpr(nlrow), NULL, NULL, NULL, NULL, &nquadvars, NULL, NULL, NULL);
2015 * To make the facet valid we subtract the maximum violation from the constant part of the facet.
2064 val += facet[pos] * (SCIPvarGetUbLocal(ecaggr->vars[pos]) - SCIPvarGetLbLocal(ecaggr->vars[pos]));
2066 val -= facet[pos] * (SCIPvarGetUbLocal(ecaggr->vars[pos]) - SCIPvarGetLbLocal(ecaggr->vars[pos]));
2077 /* there seem to be numerical problems if the violation is too large; in this case we reject the facet */
2122 SCIP_CALL( SCIPlpiCreate(&(sepadata->lpi), SCIPgetMessagehdlr(scip), "e.c. LP", SCIP_OBJSEN_MINIMIZE) );
2228 bound1 = ((poweroftwo[idx1]) & k) == 0 ? SCIPvarGetLbLocal(ecaggr->vars[idx1]) : SCIPvarGetUbLocal(ecaggr->vars[idx1]); /*lint !e661*/
2229 bound2 = ((poweroftwo[idx2]) & k) == 0 ? SCIPvarGetLbLocal(ecaggr->vars[idx2]) : SCIPvarGetUbLocal(ecaggr->vars[idx2]); /*lint !e661*/
2250 assert(!SCIPisFeasEQ(scip, ub - lb, 0.0)); /* this would mean that a variable has been fixed */
2271 * where \f$f\f$ is an edge concave function, \f$x\in [l,u]\f$ is a solution of the current relaxation, and \f$v_i\f$ are the vertices of \f$[l,u]\f$.
2272 * The method transforms the problem to the domain \f$[0,1]^n\f$, computes a facet, and transforms this facet to the
2288 SCIP_Real* facet, /**< array to store the coefficients of the resulting facet; size has to be at least (ecaggr->nvars + 1) */
2289 SCIP_Real* facetval, /**< pointer to store the value of the facet evaluated at the current solution */
2363 side[i] = transformValue(scip, SCIPvarGetLbLocal(x), SCIPvarGetUbLocal(x), SCIPgetSolVal(scip, sol, x));
2400 /* the dual solution corresponds to the coefficients of the facet in the transformed problem; note that it might be
2436 * we now have the linear underestimator L(x) = beta^T x + beta_0, which needs to be transform to the original space
2437 * the underestimator in the original space, G(x) = alpha^T x + alpha_0, is given by G(x) = L(T(x)), where T(.) is
2489 /** method to add a facet of the convex envelope of an edge-concave aggregation to a given cut */
2652 assert(SCIPisLE(scip, refpoint, SCIPvarGetUbLocal(x)) && SCIPisGE(scip, refpoint, SCIPvarGetLbLocal(x)));
2655 SCIPaddSquareLinearization(scip, coeff, refpoint, SCIPvarIsIntegral(x), &lincoef, &linconst, success);
2657 SCIPaddSquareSecant(scip, coeff, SCIPvarGetLbLocal(x), SCIPvarGetUbLocal(x), &lincoef, &linconst, success);
2683 refpointx = SCIPisLT(scip, refpointx, SCIPvarGetLbLocal(x)) ? SCIPvarGetLbLocal(x) : refpointx;
2684 refpointx = SCIPisGT(scip, refpointx, SCIPvarGetUbLocal(x)) ? SCIPvarGetUbLocal(x) : refpointx;
2685 refpointy = SCIPisLT(scip, refpointy, SCIPvarGetLbLocal(y)) ? SCIPvarGetLbLocal(y) : refpointy;
2686 refpointy = SCIPisGT(scip, refpointy, SCIPvarGetUbLocal(y)) ? SCIPvarGetUbLocal(y) : refpointy;
2687 assert(SCIPisLE(scip, refpointx, SCIPvarGetUbLocal(x)) && SCIPisGE(scip, refpointx, SCIPvarGetLbLocal(x)));
2688 assert(SCIPisLE(scip, refpointy, SCIPvarGetUbLocal(y)) && SCIPisGE(scip, refpointy, SCIPvarGetLbLocal(y)));
2690 SCIPaddBilinMcCormick(scip, coeff, SCIPvarGetLbLocal(x), SCIPvarGetUbLocal(x), refpointx, SCIPvarGetLbLocal(y),
2700 SCIPdebugMsg(scip, "add to cut: %e * %s + %e * %s + %e\n", lincoefx, SCIPvarGetName(x), lincoefy,
2710 * the remaining bilinear terms will be underestimated with McCormick, secants or linearizations;
2743 /* we use SCIPgetDepth because we add the cut to the global cut pool if cut is globally valid */
2748 SCIP_CALL( SCIPcreateEmptyRowSepa(scip, &cut, sepa, cutname, -SCIPinfinity(scip), SCIPinfinity(scip), islocalcut, FALSE,
2767 SCIPdebugMsg(scip, " %s = [%e, %e] solval = %e\n", SCIPvarGetName(var), SCIPvarGetLbLocal(var),
2782 SCIP_CALL( computeConvexEnvelopeFacet(scip, sepadata, sol, ecaggr, bestfacet, &bestfacetval, &found) );
2784 SCIPdebugMsg(scip, "found facet for edge-concave aggregation %d/%d ? %s\n", i, nlrowaggr->necaggr,
2792 /* do not add any cut because we did not found a facet for at least one edge-concave aggregation */
2797 SCIP_CALL( addFacetToCut(scip, sol, cut, bestfacet, ecaggr->vars, ecaggr->nvars, &cutconstant, &cutactivity,
2804 /* compute an underestimate for each bilinear term which is not in any edge-concave aggregation */
2816 SCIP_CALL( addBilinearTermToCut(scip, sol, cut, x, y, nlrowaggr->remtermcoefs[i], &cutconstant, &cutactivity,
2835 SCIP_CALL( addLinearTermToCut(scip, sol, cut, x, coef, &cutconstant, &cutactivity, &success) );
2848 /* assert(SCIPisFeasEQ(scip, cutactivity - cutconstant, SCIPgetRowSolActivity(scip, cut, sol)) ); */
2855 SCIProwGetName(cut), SCIPgetRowSolActivity(scip, cut, sol), SCIPgetCutEfficacy(scip, sol, cut),
2858 /* try to add the cut has a finite rhs, is efficacious, and does not exceed the maximum cut range */
2859 if( !SCIPisInfinity(scip, nlrowaggr->rhs - cutconstant) && SCIPisCutEfficacious(scip, sol, cut)
2912 if( SCIPisInfinity(scip, REALABS(SCIPvarGetLbLocal(var))) || SCIPisInfinity(scip, REALABS(SCIPvarGetUbLocal(var)))
2920 if( nlrowaggr->quadvar2aggr[i] >= 0 && SCIPisFeasEQ(scip, SCIPvarGetLbLocal(var), SCIPvarGetUbLocal(var)) )
2980 SCIPdebugMsg(scip, "found a cut: %s cutoff: %s\n", separated ? "yes" : "no", cutoff ? "yes" : "no");
3035 /** solving process deinitialization method of separator (called before branch and bound process data is freed) */
3136 )
3145 SCIP_CALL( SCIPincludeSepaBasic(scip, &sepa, SEPA_NAME, SEPA_DESC, SEPA_PRIORITY, SEPA_FREQ, SEPA_MAXBOUNDDIST,
3187 "maximal coef. range of a cut (max coef. divided by min coef.) in order to be added to LP relaxation",
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:4067
#define SCIPreallocBlockMemoryArray(scip, ptr, oldnum, newnum)
Definition: scip_mem.h:99
Definition: type_result.h:42
SCIP_RETCODE SCIPcacheRowExtensions(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:1635
static SCIP_RETCODE addFacetToCut(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut, SCIP_Real *facet, SCIP_VAR **vars, int nvars, SCIP_Real *cutconstant, SCIP_Real *cutactivity, SCIP_Bool *success)
Definition: sepa_eccuts.c:2493
SCIP_Bool SCIPisFeasEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:780
SCIP_RETCODE SCIPcreateConsBasicLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs)
Definition: cons_linear.c:18053
Definition: struct_scip.h:68
TCLIQUE_Bool tcliqueCreate(TCLIQUE_GRAPH **tcliquegraph)
Definition: tclique_graph.c:185
Definition: type_prob.h:47
SCIP_RETCODE SCIPflushRowExtensions(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:1658
#define SCIPallocClearBufferArray(scip, ptr, num)
Definition: scip_mem.h:126
SCIP_RETCODE SCIPlpiGetSol(SCIP_LPI *lpi, SCIP_Real *objval, SCIP_Real *primsol, SCIP_Real *dualsol, SCIP_Real *activity, SCIP_Real *redcost)
Definition: lpi_clp.cpp:2788
TCLIQUE_Bool tcliqueAddNode(TCLIQUE_GRAPH *tcliquegraph, int node, TCLIQUE_WEIGHT weight)
Definition: tclique_graph.c:340
SCIP_RETCODE SCIPgetRealParam(SCIP *scip, const char *name, SCIP_Real *value)
Definition: scip_param.c:307
SCIP_RETCODE SCIPaddVarToRow(SCIP *scip, SCIP_ROW *row, SCIP_VAR *var, SCIP_Real val)
Definition: scip_lp.c:1698
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:497
static SCIP_RETCODE sepadataFreeNlrows(SCIP *scip, SCIP_SEPADATA *sepadata)
Definition: sepa_eccuts.c:746
SCIP_RETCODE SCIPlpiChgSides(SCIP_LPI *lpi, int nrows, const int *ind, const SCIP_Real *lhs, const SCIP_Real *rhs)
Definition: lpi_clp.cpp:1167
Definition: struct_var.h:207
static SCIP_RETCODE searchEcAggrWithCliques(SCIP *scip, TCLIQUE_GRAPH *graph, SCIP_SEPADATA *sepadata, SCIP_NLROW *nlrow, int *quadvar2aggr, int nfoundsofar, SCIP_Bool rhsaggr, SCIP_Bool *foundaggr, SCIP_Bool *foundclique)
Definition: sepa_eccuts.c:1405
static SCIP_RETCODE searchEcAggrWithMIP(SCIP *subscip, SCIP_Real timelimit, int nedges, SCIP_Bool *aggrleft, SCIP_Bool *found)
Definition: sepa_eccuts.c:1248
SCIP_RETCODE SCIPsetHeuristics(SCIP *scip, SCIP_PARAMSETTING paramsetting, SCIP_Bool quiet)
Definition: scip_param.c:906
SCIP_RETCODE SCIPcheckExprQuadratic(SCIP *scip, SCIP_EXPR *expr, SCIP_Bool *isquadratic)
Definition: scip_expr.c:2342
SCIP_RETCODE SCIPhashmapCreate(SCIP_HASHMAP **hashmap, BMS_BLKMEM *blkmem, int mapsize)
Definition: misc.c:3024
static SCIP_RETCODE createTcliqueGraph(SCIP_NLROW *nlrow, TCLIQUE_GRAPH **graph, SCIP_Real *nodeweights)
Definition: sepa_eccuts.c:1314
static SCIP_RETCODE computeConvexEnvelopeFacet(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_SOL *sol, SCIP_ECAGGR *ecaggr, SCIP_Real *facet, SCIP_Real *facetval, SCIP_Bool *success)
Definition: sepa_eccuts.c:2285
Definition: type_var.h:62
SCIP_RETCODE SCIPhashmapInsertInt(SCIP_HASHMAP *hashmap, void *origin, int image)
Definition: misc.c:3142
SCIP_RETCODE SCIPcreateVarBasic(SCIP *scip, SCIP_VAR **var, const char *name, SCIP_Real lb, SCIP_Real ub, SCIP_Real obj, SCIP_VARTYPE vartype)
Definition: scip_var.c:194
static SCIP_RETCODE sepadataFree(SCIP *scip, SCIP_SEPADATA **sepadata)
Definition: sepa_eccuts.c:776
edge concave cut separator
Definition: type_result.h:49
Definition: tclique.h:66
tclique user interface
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:445
Definition: struct_sepa.h:46
SCIP_RETCODE SCIPsetSepaCopy(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPACOPY((*sepacopy)))
Definition: scip_sepa.c:151
SCIP_RETCODE SCIPsetRealParam(SCIP *scip, const char *name, SCIP_Real value)
Definition: scip_param.c:603
SCIP_Bool SCIPexprAreQuadraticExprsVariables(SCIP_EXPR *expr)
Definition: expr.c:4195
SCIP_RETCODE SCIPaddIntParam(SCIP *scip, const char *name, const char *desc, int *valueptr, SCIP_Bool isadvanced, int defaultvalue, int minvalue, int maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:83
static SCIP_RETCODE nlrowaggrAddRemBilinTerm(SCIP_NLROWAGGR *nlrowaggr, SCIP_VAR *x, SCIP_VAR *y, SCIP_Real coef)
Definition: sepa_eccuts.c:407
SCIP_RETCODE SCIPlpiCreate(SCIP_LPI **lpi, SCIP_MESSAGEHDLR *messagehdlr, const char *name, SCIP_OBJSEN objsen)
Definition: lpi_clp.cpp:531
SCIP_RETCODE SCIPaddCoefLinear(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real val)
Definition: cons_linear.c:18205
SCIP_RETCODE SCIPcreateProbBasic(SCIP *scip, const char *name)
Definition: scip_prob.c:180
SCIP_RETCODE SCIPlpiAddCols(SCIP_LPI *lpi, int ncols, const SCIP_Real *obj, const SCIP_Real *lb, const SCIP_Real *ub, char **colnames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
Definition: lpi_clp.cpp:758
static SCIP_RETCODE ecaggrCreateEmpty(SCIP *scip, SCIP_ECAGGR **ecaggr, int nquadvars, int nquadterms)
Definition: sepa_eccuts.c:177
Definition: struct_sol.h:73
SCIP_Bool SCIPhashmapExists(SCIP_HASHMAP *hashmap, void *origin)
Definition: misc.c:3373
static SCIP_RETCODE ecaggrFree(SCIP *scip, SCIP_ECAGGR **ecaggr)
Definition: sepa_eccuts.c:207
static SCIP_Bool isPossibleToComputeCut(SCIP *scip, SCIP_SOL *sol, SCIP_NLROWAGGR *nlrowaggr)
Definition: sepa_eccuts.c:2891
SCIP_RETCODE SCIPcreateConsBasicXor(SCIP *scip, SCIP_CONS **cons, const char *name, SCIP_Bool rhs, int nvars, SCIP_VAR **vars)
Definition: cons_xor.c:5929
#define SCIPduplicateBlockMemoryArray(scip, ptr, source, num)
Definition: scip_mem.h:105
Definition: struct_misc.h:137
SCIP_Bool SCIPisCutEfficacious(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut)
Definition: scip_cut.c:117
SCIP_RETCODE SCIPsetObjsense(SCIP *scip, SCIP_OBJSENSE objsense)
Definition: scip_prob.c:1250
SCIP_RETCODE SCIPlpiAddRows(SCIP_LPI *lpi, int nrows, const SCIP_Real *lhs, const SCIP_Real *rhs, char **rownames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
Definition: lpi_clp.cpp:914
TCLIQUE_Bool tcliqueAddEdge(TCLIQUE_GRAPH *tcliquegraph, int node1, int node2)
Definition: tclique_graph.c:380
static SCIP_RETCODE doSeachEcAggr(SCIP *scip, SCIP *subscip, SCIP_SEPADATA *sepadata, SCIP_NLROW *nlrow, SCIP_SOL *sol, SCIP_Bool rhsaggr, int *quadvar2aggr, int *nfound)
Definition: sepa_eccuts.c:1554
Definition: type_result.h:44
Definition: struct_cons.h:46
SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:458
void tcliqueChangeWeight(TCLIQUE_GRAPH *tcliquegraph, int node, TCLIQUE_WEIGHT weight)
Definition: tclique_graph.c:362
Definition: type_retcode.h:57
Definition: type_lpi.h:43
SCIP_RETCODE SCIPchgVarUb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound)
Definition: scip_var.c:4766
static SCIP_Real phi(SCIP *scip, SCIP_Real val, SCIP_Real lb, SCIP_Real ub)
Definition: sepa_eccuts.c:846
static SCIP_RETCODE sepadataAddNlrowaggr(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_NLROWAGGR *nlrowaggr)
Definition: sepa_eccuts.c:802
static SCIP_RETCODE storeAggrFromMIP(SCIP *subscip, SCIP_NLROW *nlrow, SCIP_VAR **forwardarcs, SCIP_VAR **backwardarcs, int *quadvar2aggr, int nfoundsofar)
Definition: sepa_eccuts.c:1164
internal methods for NLP management
Definition: sepa_eccuts.c:90
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:4114
static SCIP_RETCODE addBilinearTermToCut(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut, SCIP_VAR *x, SCIP_VAR *y, SCIP_Real coeff, SCIP_Real *cutconstant, SCIP_Real *cutactivity, SCIP_Bool *success)
Definition: sepa_eccuts.c:2602
Definition: type_retcode.h:42
void tcliqueMaxClique(TCLIQUE_GETNNODES((*getnnodes)), TCLIQUE_GETWEIGHTS((*getweights)), TCLIQUE_ISEDGE((*isedge)), TCLIQUE_SELECTADJNODES((*selectadjnodes)), TCLIQUE_GRAPH *tcliquegraph, TCLIQUE_NEWSOL((*newsol)), TCLIQUE_DATA *tcliquedata, int *maxcliquenodes, int *nmaxcliquenodes, TCLIQUE_WEIGHT *maxcliqueweight, TCLIQUE_WEIGHT maxfirstnodeweight, TCLIQUE_WEIGHT minweight, int maxntreenodes, int backtrackfreq, int maxnzeroextensions, int fixednode, int *ntreenodes, TCLIQUE_STATUS *status)
Definition: tclique_branch.c:1010
#define SCIPensureBlockMemoryArray(scip, ptr, arraysizeptr, minsize)
Definition: scip_mem.h:107
SCIP_RETCODE SCIPaddRow(SCIP *scip, SCIP_ROW *row, SCIP_Bool forcecut, SCIP_Bool *infeasible)
Definition: scip_cut.c:250
Definition: type_retcode.h:43
SCIP_RETCODE SCIPincludeSepaBasic(SCIP *scip, SCIP_SEPA **sepa, const char *name, const char *desc, int priority, int freq, SCIP_Real maxbounddist, SCIP_Bool usessubscip, SCIP_Bool delay, SCIP_DECL_SEPAEXECLP((*sepaexeclp)), SCIP_DECL_SEPAEXECSOL((*sepaexecsol)), SCIP_SEPADATA *sepadata)
Definition: scip_sepa.c:109
Definition: struct_expr.h:104
SCIP_RETCODE SCIPsetSepaExitsol(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAEXITSOL((*sepaexitsol)))
Definition: scip_sepa.c:231
SCIP_RETCODE SCIPincludeDefaultPlugins(SCIP *scip)
Definition: scipdefplugins.c:37
SCIP_RETCODE SCIPchgRowRhs(SCIP *scip, SCIP_ROW *row, SCIP_Real rhs)
Definition: scip_lp.c:1607
static SCIP_RETCODE searchEcAggr(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_NLROW *nlrow, SCIP_SOL *sol, SCIP_Bool rhsaggr, int *quadvar2aggr, int *nfound)
Definition: sepa_eccuts.c:1751
Definition: struct_lp.h:201
static SCIP_Real transformValue(SCIP *scip, SCIP_Real lb, SCIP_Real ub, SCIP_Real val)
Definition: sepa_eccuts.c:2241
SCIP_RETCODE SCIPsetIntParam(SCIP *scip, const char *name, int value)
Definition: scip_param.c:487
SCIP_RETCODE SCIPcreateEmptyRowSepa(SCIP *scip, SCIP_ROW **row, SCIP_SEPA *sepa, const char *name, SCIP_Real lhs, SCIP_Real rhs, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool removable)
Definition: scip_lp.c:1453
SCIP_Real SCIPgetCutEfficacy(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut)
Definition: scip_cut.c:94
void SCIPaddSquareSecant(SCIP *scip, SCIP_Real sqrcoef, SCIP_Real lb, SCIP_Real ub, SCIP_Real *lincoef, SCIP_Real *linconstant, SCIP_Bool *success)
Definition: expr_pow.c:3303
static SCIP_RETCODE nlrowaggrCreate(SCIP *scip, SCIP_NLROW *nlrow, SCIP_NLROWAGGR **nlrowaggr, int *quadvar2aggr, int nfound, SCIP_Bool rhsaggr)
Definition: sepa_eccuts.c:432
void SCIPaddBilinMcCormick(SCIP *scip, SCIP_Real bilincoef, SCIP_Real lbx, SCIP_Real ubx, SCIP_Real refpointx, SCIP_Real lby, SCIP_Real uby, SCIP_Real refpointy, SCIP_Bool overestimate, SCIP_Real *lincoefx, SCIP_Real *lincoefy, SCIP_Real *linconstant, SCIP_Bool *success)
Definition: nlhdlr_bilinear.c:1762
static SCIP_RETCODE isCandidate(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_NLROW *nlrow, SCIP_Bool *rhscandidate, SCIP_Bool *lhscandidate)
Definition: sepa_eccuts.c:1784
static SCIP_RETCODE nlrowaggrAddQuadraticVar(SCIP *scip, SCIP_NLROWAGGR *nlrowaggr, SCIP_VAR *quadvar)
Definition: sepa_eccuts.c:387
SCIP_Real SCIPgetRowSolActivity(SCIP *scip, SCIP_ROW *row, SCIP_SOL *sol)
Definition: scip_lp.c:2141
static SCIP_RETCODE nlrowaggrStoreLinearTerms(SCIP *scip, SCIP_NLROWAGGR *nlrowaggr, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nlinvars)
Definition: sepa_eccuts.c:313
SCIP_RETCODE SCIPlpiChgBounds(SCIP_LPI *lpi, int ncols, const int *ind, const SCIP_Real *lb, const SCIP_Real *ub)
Definition: lpi_clp.cpp:1084
static SCIP_RETCODE nlrowaggrFree(SCIP *scip, SCIP_NLROWAGGR **nlrowaggr)
Definition: sepa_eccuts.c:654
static SCIP_RETCODE separateCuts(SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, int depth, SCIP_SOL *sol, SCIP_RESULT *result)
Definition: sepa_eccuts.c:2934
Constraint handler for XOR constraints, .
static SCIP_RETCODE addLinearTermToCut(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut, SCIP_VAR *x, SCIP_Real coeff, SCIP_Real *cutconstant, SCIP_Real *cutactivity, SCIP_Bool *success)
Definition: sepa_eccuts.c:2551
SCIP_RETCODE SCIPsetSepaFree(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAFREE((*sepafree)))
Definition: scip_sepa.c:167
SCIP_Bool SCIPisGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:484
static SCIP_RETCODE createLP(SCIP *scip, SCIP_SEPADATA *sepadata)
Definition: sepa_eccuts.c:2091
Definition: sepa_eccuts.c:105
static SCIP_RETCODE findAndStoreEcAggregations(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_NLROW *nlrow, SCIP_SOL *sol)
Definition: sepa_eccuts.c:1915
static SCIP_RETCODE nlrowaggrAddLinearTerm(SCIP *scip, SCIP_NLROWAGGR *nlrowaggr, SCIP_VAR *linvar, SCIP_Real lincoef)
Definition: sepa_eccuts.c:354
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip_cons.c:1119
void SCIPexprGetQuadraticBilinTerm(SCIP_EXPR *expr, int termidx, SCIP_EXPR **expr1, SCIP_EXPR **expr2, SCIP_Real *coef, int *pos2, SCIP_EXPR **prodexpr)
Definition: expr.c:4157
Definition: lpi_clp.cpp:104
Definition: struct_nlp.h:64
static SCIP_RETCODE createMIP(SCIP *scip, SCIP *subscip, SCIP_SEPADATA *sepadata, SCIP_NLROW *nlrow, SCIP_Bool rhsaggr, SCIP_VAR **forwardarcs, SCIP_VAR **backwardarcs, SCIP_Real *nodeweights, int *nedges, int *narcs)
Definition: sepa_eccuts.c:871
SCIP_RETCODE SCIPprintRow(SCIP *scip, SCIP_ROW *row, FILE *file)
Definition: scip_lp.c:2209
SCIP_RETCODE SCIPprintNlRow(SCIP *scip, SCIP_NLROW *nlrow, FILE *file)
Definition: scip_nlp.c:1556
static SCIP_RETCODE sepadataCreate(SCIP *scip, SCIP_SEPADATA **sepadata)
Definition: sepa_eccuts.c:730
SCIP_RETCODE SCIPlpiChgObj(SCIP_LPI *lpi, int ncols, const int *ind, const SCIP_Real *obj)
Definition: lpi_clp.cpp:1240
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:471
static SCIP_RETCODE computeCut(SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_NLROWAGGR *nlrowaggr, SCIP_SOL *sol, SCIP_Bool *separated, SCIP_Bool *cutoff)
Definition: sepa_eccuts.c:2716
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:111
void SCIPaddSquareLinearization(SCIP *scip, SCIP_Real sqrcoef, SCIP_Real refpoint, SCIP_Bool isint, SCIP_Real *lincoef, SCIP_Real *linconstant, SCIP_Bool *success)
Definition: expr_pow.c:3235
static SCIP_RETCODE ecaggrAddBilinTerm(SCIP *scip, SCIP_ECAGGR *ecaggr, SCIP_VAR *x, SCIP_VAR *y, SCIP_Real coef)
Definition: sepa_eccuts.c:239
Definition: objbenders.h:43
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1361
default SCIP plugins
Definition: type_stat.h:62
SCIP_RETCODE SCIPaddRealParam(SCIP *scip, const char *name, const char *desc, SCIP_Real *valueptr, SCIP_Bool isadvanced, SCIP_Real defaultvalue, SCIP_Real minvalue, SCIP_Real maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:139
static SCIP_RETCODE ecaggrAddQuadvar(SCIP_ECAGGR *ecaggr, SCIP_VAR *x)
Definition: sepa_eccuts.c:228
static SCIP_RETCODE updateMIP(SCIP *subscip, SCIP_NLROW *nlrow, SCIP_VAR **forwardarcs, SCIP_VAR **backwardarcs, int *quadvar2aggr, int *nedges)
Definition: sepa_eccuts.c:1085
SCIP_RETCODE SCIPsetLongintParam(SCIP *scip, const char *name, SCIP_Longint value)
Definition: scip_param.c:545
Definition: type_result.h:48
Definition: type_paramset.h:61
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:57
static SCIP_Bool checkRikun(SCIP *scip, SCIP_ECAGGR *ecaggr, SCIP_Real *fvals, SCIP_Real *facet)
Definition: sepa_eccuts.c:2020
SCIP_RETCODE SCIPprintSol(SCIP *scip, SCIP_SOL *sol, FILE *file, SCIP_Bool printzeros)
Definition: scip_sol.c:1775