sepa_mcf.c
Go to the documentation of this file.
34 * (1) sum_{a in delta^+(v)} f_a^k - sum_{a in delta^-(v)} f_a^k <= -d_v^k for all v in V and k in K
37 * Constraints (1) are flow conservation constraints, which say that for each commodity k and node v the
38 * outflow (delta^+(v)) minus the inflow (delta^-(v)) of a node v must not exceed the negative of the demand of
39 * node v in commodity k. To say it the other way around, inflow minus outflow must be at least equal to the demand.
40 * Constraints (2) are the arc capacity constraints, which say that the sum of all flow over an arc a must not
42 * c_a x_a does not need to be a single product of a capacity and an integer variable; we also accept general scalar
46/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
90#define SEPA_DELAY FALSE /**< should separation method be delayed, if other separators found cuts? */
94#define DEFAULT_MAXWEIGHTRANGE 1e+06 /**< maximal valid range max(|weights|)/min(|weights|) of row weights for CMIR */
95#define DEFAULT_MAXTESTDELTA 20 /**< maximal number of different deltas to try (-1: unlimited) for CMIR */
96#define DEFAULT_TRYNEGSCALING FALSE /**< should negative values also be tested in scaling? for CMIR */
97#define DEFAULT_FIXINTEGRALRHS TRUE /**< should an additional variable be complemented if f0 = 0? for CMIR */
98#define DEFAULT_DYNAMICCUTS TRUE /**< should generated cuts be removed from the LP if they are no longer tight? */
100#define DEFAULT_MAXSEPACUTS 100 /**< maximal number of cuts separated per separation round (-1: unlimited) */
101#define DEFAULT_MAXSEPACUTSROOT 200 /**< maximal number of cuts separated per separation round in root node (-1: unlimited) */
102#define DEFAULT_MAXINCONSISTENCYRATIO 0.02 /**< maximum inconsistency ratio (inconsistencies/(arcs*commodities)) at all */
103#define DEFAULT_MAXARCINCONSISTENCYRATIO 0.5 /**< maximum inconsistency ratio for arcs not to be deleted */
104#define DEFAULT_CHECKCUTSHORECONNECTIVITY TRUE /**< should we only separate if the cuts shores are connected */
105#define DEFAULT_SEPARATESINGLENODECUTS TRUE /**< should we separate inequalities based on single node cuts */
120#define MAXFLOWVARFLOWROWRATIO 100.0 /**< maximum ratio flowvars/flowrows for the separator to be switched on */
121#define MAXARCNODERATIO 100.0 /**< maximum ratio arcs/nodes for the separator to be switched on */
124#define MINCOMNODESFRACTION 0.5 /**< minimal size of commodity relative to largest commodity to keep it in the network */
127#define MAXCAPACITYSLACK 0.1 /**< maximal slack of weighted capacity constraints to use in aggregation */
128#define UNCAPACITATEDARCSTRESHOLD 0.8 /**< threshold for the percentage of commodities an uncapacitated arc should appear in */
131/* #define OUTPUTGRAPH should a .gml graph of the network be generated for debugging purposes? */
157 SCIP_MCFMODELTYPE_UNDIRECTED = 2 /**< directed network where anti-parallel arcs share the capacity */
158};
167};
175 SCIP_Real** nodeflowscales; /**< scaling factors to convert nodeflowrows[v][k] into a +/-1 <= row */
176 SCIP_Bool** nodeflowinverted; /**< does nodeflowrows[v][k] have to be inverted to fit the network structure? */
197 int nclusters; /**< number of clusters to generate in the shrunken network -- default separation */
198 SCIP_Real maxweightrange; /**< maximal valid range max(|weights|)/min(|weights|) of row weights */
199 int maxtestdelta; /**< maximal number of different deltas to try (-1: unlimited) -- default separation */
202 SCIP_Bool dynamiccuts; /**< should generated cuts be removed from the LP if they are no longer tight? */
205 int maxsepacutsroot; /**< maximal number of cmir cuts separated per separation round in root node -- default separation */
206 SCIP_Real maxinconsistencyratio; /**< maximum inconsistency ratio (inconsistencies/(arcs*commodities)) for separation at all*/
207 SCIP_Real maxarcinconsistencyratio; /**< maximum inconsistency ratio for arcs not to be deleted */
208 SCIP_Bool checkcutshoreconnectivity;/**< should we only separate if the cuts shores are connected */
209 SCIP_Bool separatesinglenodecuts; /**< should we separate inequalities based on single node cuts ? */
210 SCIP_Bool separateflowcutset; /**< should we separate flowcutset inequalities on the network cuts ? */
211 SCIP_Bool separateknapsack; /**< should we separate knapsack cover inequalities on the network cuts ? */
219 unsigned char* flowrowsigns; /**< potential or actual sides of rows to be used as flow conservation constraint */
221 SCIP_Real* flowrowscores; /**< score value indicating how sure we are that this is indeed a flow conservation constraint */
222 unsigned char* capacityrowsigns; /**< potential or actual sides of rows to be used as capacity constraint */
223 SCIP_Real* capacityrowscores; /**< score value indicating how sure we are that this is indeed a capacity constraint */
224 int* flowcands; /**< list of row indices that are candidates for flow conservation constraints */
231 int nemptycommodities; /**< number of commodities that have been discarded but still counted in 'ncommodities' */
248 int* newcols; /**< columns of current commodity that have to be inspected for incident flow conservation rows */
353 SCIPfreeBlockMemoryArrayNull(scip, &(*mcfnetwork)->nodeflowrows[v], (*mcfnetwork)->ncommodities);
354 SCIPfreeBlockMemoryArrayNull(scip, &(*mcfnetwork)->nodeflowscales[v], (*mcfnetwork)->ncommodities);
355 SCIPfreeBlockMemoryArrayNull(scip, &(*mcfnetwork)->nodeflowinverted[v], (*mcfnetwork)->ncommodities);
489 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &mcfnetwork->nodeflowinverted, mcfnetwork->nnodes) );
492 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &mcfnetwork->nodeflowrows[v], mcfnetwork->ncommodities) ); /*lint !e866*/
493 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &mcfnetwork->nodeflowscales[v], mcfnetwork->ncommodities) ); /*lint !e866*/
494 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &mcfnetwork->nodeflowinverted[v], mcfnetwork->ncommodities) ); /*lint !e866*/
504 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &mcfnetwork->arccapacityscales, mcfnetwork->narcs) );
582 /* Scale constraint such that the coefficients for the flow variables are normalized in such a way that coefficients in
584 * This is needed for binary flow variable models in which the demand of each commodity is stored as the coefficient in
585 * the capacity constraints. Since it may happen (e.g., due to presolve) that different capacity constraints are scaled
586 * differently, we need to find scaling factors to make the demand coefficients of each commodity equal.
587 * To do this, we scale the first capacity constraint with +1 and then store the coefficients of the flow variables
588 * as target demands for the commodities. Then, we scale the other constraints in such a way that these demands are hit, if possible.
589 * Note that for continuous flow variable models, the coefficients in the capacity constraints are usually +1.0.
611 /* use negative scaling if we use the left hand side, use positive scaling if we use the right hand side */
637 assert(0 <= compnodeid[mcfdata->arcsources[globala]] && compnodeid[mcfdata->arcsources[globala]] < mcfnetwork->nnodes);
643 assert(0 <= compnodeid[mcfdata->arctargets[globala]] && compnodeid[mcfdata->arctargets[globala]] < mcfnetwork->nnodes);
708 MCFdebugMessage("arc %2d [%2d -> %2d]: ", a, mcfnetwork->arcsources[a], mcfnetwork->arctargets[a]);
711 MCFdebugMessage("<%s> [%+g]\n", SCIProwGetName(mcfnetwork->arccapacityrows[a]), mcfnetwork->arccapacityscales[a]);
759 MCFdebugMessage(" col <%s>: arc %d\n", SCIPvarGetName(SCIPcolGetVar(cols[c])), colarcid != NULL ? colarcid[c] : -1);
764 MCFdebugMessage(" row <%s>: node %d [sign:%+d, inv:%+d]\n", SCIProwGetName(rows[r]), rownodeid != NULL ? rownodeid[r] : -1,
817 MCFdebugMessage(" col <%s> [%g,%g]\n", SCIPvarGetName(var), SCIPvarGetLbGlobal(var), SCIPvarGetUbGlobal(var));
827 if( rowcommodity[r] == -1 && (capacityrowsigns == NULL || (capacityrowsigns[r] & (LHSASSIGNED | RHSASSIGNED)) == 0) )
1021 /** @todo go through list of several model types, depending on the current model type throw away invalid constraints
1059 SCIPdebugMsg(scip, "%4d [score: %2g]: %s\n", mcfdata->flowcands[r], flowrowscores[mcfdata->flowcands[r]],
1111 /* identify model type and set the maximal number of flow variables per capacity constraint and commodity */
1115 maxcolspercommoditylimit = 2; /* will be set to 1 later if we detect that the network is directed */
1282 ABS((nposflowcoefs + nnegflowcoefs)/(SCIP_Real)ncoveredcommodities - maxcolspercommoditylimit);
1291 if( (capacityrowsigns[r] & RHSPOSSIBLE) != 0 && nnegflowcoefs == 0 && nposcapacitycoefs == 0 && nnegcapacitycoefs > 0 )
1293 if( (capacityrowsigns[r] & LHSPOSSIBLE) != 0 && nposflowcoefs == 0 && nposcapacitycoefs > 0 && nnegcapacitycoefs == 0 )
1301 * use slightly increased denominator in order to not increase score too much for very few commodities
1323 capacityrowscores[r] += 20.0 * MAX(nposflowcoefs, nnegflowcoefs)/MAX(1.0,(SCIP_Real)(nposflowcoefs + nnegflowcoefs));
1325 /* capacity coefficients are mostly of the same sign: score +10*max(npos,nneg)/(npos+nneg+1) */
1326 capacityrowscores[r] += 10.0 * MAX(nposcapacitycoefs, nnegcapacitycoefs)/(SCIP_Real)(nposcapacitycoefs+nnegcapacitycoefs+1.0);
1337 SCIPdebugMsg(scip, "row <%s>: maxcolspercommodity=%d capacityrowsign=%d nposflowcoefs=%d nnegflowcoefs=%d nposcapacitycoefs=%d nnegcapacitycoefs=%d nbadcoefs=%d nactivecommodities=%d sameflowcoef=%g -> score=%g\n",
1338 SCIProwGetName(row), maxcolspercommodity[r], capacityrowsigns[r], nposflowcoefs, nnegflowcoefs, nposcapacitycoefs, nnegcapacitycoefs, nbadcoefs, nactivecommodities, sameflowcoef, capacityrowscores[r]);
1343 /* if the model type should be detected automatically, count the number of directed and undirected capacity candidates */
1355 SCIPdebugMsg(scip, "row <%s>: rowsign = %d nposflowcoefs = %d nnegflowcoefs = %d -> discard\n",
1369 modeltype == SCIP_MCFMODELTYPE_DIRECTED ? "directed" : "undirected", directedcandsscore, undirectedcandsscore);
1418 SCIPsortInd(mcfdata->capacitycands, compCands, (void*)capacityrowscores, mcfdata->ncapacitycands);
1425 capacityrowscores[mcfdata->capacitycands[r]], SCIProwGetName(rows[mcfdata->capacitycands[r]]));
1449 SCIP_CALL( SCIPreallocMemoryArray(scip, &mcfdata->commoditysigns, mcfdata->commoditysignssize) );
1499 SCIPdebugMsg(scip, "**** creating new arc %d: %d -> %d ****\n", mcfdata->narcs, source, target);
1751 SCIPdebugMsg(scip, "deleting commodity %d (%d total commodities) with %d flow rows\n", k, ncommodities, nrows);
1809/** checks whether the given row fits into the given commodity and returns the possible flow row signs */
1817 SCIP_Bool* invertcommodity /**< pointer to store whether the commodity has to be inverted to accommodate the row */
1936 SCIPdebugMsg(scip, " -> discard flow row %d <%s>, comoditysign=%d\n", r, SCIProwGetName(row), commoditysigns[k]);
1948 SCIP_Bool* nextinvertcommodity /**< pointer to store whether current commodity has to be inverted to accommodate the next row */
1972 /* check if there are any columns left in the commodity that have not yet been inspected for incident flow rows */
1995 /* check whether column is incident to a valid flow row that fits into the current commodity */
2043 * Note: This is not the overall best row, only the one for the first column that has a valid row.
2128 * (iii) For the newly added columns search for an incident flow conservation constraint. Pick the one of highest ranking.
2154 SCIPdebugMsg(scip, " -------------------start new commodity %d--------------------- \n", mcfdata->ncommodities );
2164 invertCommodity(scip, mcfdata, mcfdata->ncommodities-1, comrows, nnodes, comcolids, ncomcolids);
2181 SCIPdebugMsg(scip, " -> finished commodity %d: identified %d nodes, maxnnodes=%d\n", mcfdata->ncommodities-1, nnodes, maxnnodes);
2183 /* if the commodity has too few nodes, or if it has much fewer nodes than the largest commodity, discard it */
2189 deleteCommodity(scip, mcfdata, mcfdata->ncommodities-1, comrows, nnodes, &ndelflowrows, &ndelflowvars);
2237 MCFdebugMessage("identified %d commodities (%d empty) with a maximum of %d nodes and %d flowrows, %d flowvars \n",
2252/** Arc-Detection -- identifies capacity constraints for the arcs and assigns arc ids to columns and capacity constraints */
2367 SCIPdebugMsg(scip, "discarding capacity candidate row %d <%s> [score:%g]: %d assigned flowvars, %d unassigned flowvars\n",
2368 r, SCIProwGetName(capacityrow), mcfdata->capacityrowscores[r], nassignedflowvars, nunassignedflowvars);
2398 SCIPdebugMsg(scip, "assigning capacity row %d <%s> with sign %+d to arc %d [score:%g]: %d assigned flowvars, %d unassigned flowvars\n",
2399 r, SCIProwGetName(capacityrow), (capacityrowsigns[r] & RHSASSIGNED) != 0 ? +1 : -1, mcfdata->narcs,
2409 /* due to aggregations in preprocessing it may happen that a flow variable appears in multiple capacity constraints;
2423/** collects all flow columns of all commodities (except the one of the base row) that are incident to the node described by the given flow row */
2488 assert(colarcid[caprowc] <= arcid); /* colarcid < arcid if column belongs to multiple arcs, for example, due to an aggregation in presolving */
2515 int basenposuncap, /**< number of uncapacitated vars in base node flow row with positive coeff*/
2516 int basenneguncap, /**< number of uncapacitated vars in base node flow row with negative coeff*/
2519 SCIP_Bool* invertcommodity /**< pointer to store whether the arcs in the commodity of the row have
2651 assert(modeltype == SCIP_MCFMODELTYPE_UNDIRECTED || rowcomsign * basearcpattern[arcid] * arcpattern[arcid] > 0);
2666/* calculate the score: maximize overlap and use minimal number of non-overlapping entries as tie breaker */
2693 *score += 1.0 - (ABS(nneguncap - basenposuncap) + ABS(nposuncap - basenneguncap))/(2.0 * ncols + 1.0);
2695 *score += 1.0 - (ABS(nposuncap - basenposuncap) + ABS(nneguncap - basenneguncap))/(2.0 * ncols + 1.0);
2697 *score = MAX(*score, 1e-6); /* score may get negative due to many columns having the same arcid */
2700 SCIPdebugMsg(scip, " -> node similarity: row <%s>: incompatible=%u overlap=%g rowlen=%d baserowlen=%d score=%g\n",
2850 /* we also count variables that have no arc -- these have no capacity constraint --> uncapacitated */
3203 /* for each node s, count for each other node t the number of flow variables that are not yet assigned
3206 for( v = 0; n < nflowcands; v++ ) /*lint !e440*/ /* for flexelint: n is used as abort criterion for loop */
3218 SCIPdebugMsg(scip, " node %d starts with flowcand %d: <%s>\n", v, n, SCIProwGetName(rows[sortedflowcands[n]]));
3232 /* update sourcecount and targetcount for all flow columns in the row that are not yet assigned to an arc */
3285 * Note: each column can be collected at most once for node v, because each column can appear in at most one
3286 * commodity, and in each commodity a column can have at most one +1 and one -1 entry. One of the two +/-1 entries
3305 /* if other node has been seen the first time, store it in adjlist for sparse access of count arrays */
3353 SCIPdebugMsg(scip, " -> assign arcid:%i to column <%s>\n", arcid, SCIPvarGetName(SCIPcolGetVar(cols[c])));
3390 SCIPdebugMsg(scip, " -> assign arcid:%i to column <%s>\n", arcid, SCIPvarGetName(SCIPcolGetVar(cols[c])));
3410 /* mark the incident columns that could not be assigned to a new arc such that we do not inspect them again */
3435 MCFdebugMessage("network after finding uncapacitated arcs has %d nodes, %d arcs, and %d commodities\n",
3478 MCFdebugMessage("network before cleanup has %d nodes, %d arcs, and %d commodities\n", nnodes, narcs, ncommodities);
3572 SCIPdebugMsg(scip, " -> commodity %d: %d nodes, %d arcs\n", k, nnodespercom[k], narcspercom[k]);
3595 SCIPdebugMsg(scip, " -> discarding %d of %d commodities\n", ncommodities - newncommodities, ncommodities);
3731 assert(rownodeid[r] >= 0); /* otherwise we would have deleted the commodity in the rowcommodity update above */
3754 MCFdebugMessage("network after cleanup has %d nodes, %d arcs, and %d commodities\n", nnodes, narcs, ncommodities);
4008 /* in the undirected model, we use source for the maximum and target for the second largest number of total hits */
4057 SCIPdebugMsg(scip, "arc %d: %d -> %d (len=%d, sourcecnt=%g/%g, targetcnt=%g/%g, %g/%g inconsistencies)\n",
4091 MCFdebugMessage("extracted network has %g inconsistencies (ratio %g) -> separating with effort %d\n",
4274 * 2. Sort flow conservation candidates by a ranking on how sure we are that it is indeed a constraint of the desired type.
4280 * (iii) For the newly added columns search for an incident flow conservation constraint. Pick the one of highest ranking.
4285 * 5. Sort capacity constraint candidates by a ranking on how sure we are that it is indeed a constraint of the desired type.
4286 * 6. Identify capacity constraints for the arcs and assign arc ids to columns and capacity constraints.
4341 * 2. Sort flow conservation candidates by a ranking on how sure we are that it is indeed a constraint of the desired type
4371 * 5. sort capacity constraint candidates by a ranking on how sure we are that it is indeed a constraint of the desired type
4384 /* 6. arc-detection -- identify capacity constraints for the arcs and assign arc ids to columns and capacity constraints */
4455 SCIP_CALL( identifyComponent(scip, &mcfdata, nodevisited, v, compnodes, &ncompnodes, comparcs, &ncomparcs) );
4480 SCIP_CALL( mcfnetworkFill(scip, mcfnetwork, &mcfdata, compnodeid, compnodes, ncompnodes, comparcs, ncomparcs) );
4496 SCIPdebugMsg(scip, " -> discarded network with %d nodes and %d arcs due to maxnetworks (minnodes=%d)\n",
4505 SCIPdebugMsg(scip, " -> discarded component with %d nodes and %d arcs\n", ncompnodes, ncomparcs);
4695 MCFdebugMessage(" nof flowvars: %5d of which [ %d , %d , %d ] are continuous, integer, binary\n",
4697 MCFdebugMessage(" nof capvars: %5d of which [ %d , %d , %d ] are continuous, integer, binary\n",
4884 /* For every nodepair that is used in the network (at least one arc exists having this nodepair as endnodes)
4949 SCIPdebugMsg(scip, "arc %i = (%i %i)\n", a, mcfnetwork->arcsources[a], mcfnetwork->arctargets[a]);
5004 SCIPdebugMsg(scip, " col <%s>: %g\n", SCIPvarGetName(SCIPcolGetVar(rowcols[i])), SCIPcolGetPrimsol(rowcols[i]) );
5008 SCIPdebugMsg(scip, " flow col <%s>: %g\n", SCIPvarGetName(SCIPcolGetVar(rowcols[i])), REALABS(SCIPcolGetPrimsol(rowcols[i])) );
5013 SCIPdebugMsg(scip, " cap col <%s>: %g\n", SCIPvarGetName(SCIPcolGetVar(rowcols[i])), REALABS(SCIPcolGetPrimsol(rowcols[i])) );
5018 SCIPdebugMsg(scip, "cap arc -- slack:%g -- dual:%g -- flow:%g -- cap:%g \n", scale * slack, dualsol/scale, totalflow * scale, totalcap * scale);
5054 SCIPdebugMsg(scip, "nodepair known [%d,%d] -- old weight:%g -- new weight:%g\n", nodepair.node1,nodepair.node2,nodepairptr->weight,
5067 SCIPdebugMsg(scip, "new nodepair [%d,%d]-- weight:%g\n", nodepair.node1, nodepair.node2, nodepair.weight);
5100 SCIP_CALL( SCIPpqueueCreate(&(*nodepairqueue)->pqueue, nnodepairs, 2.0, compNodepairs, NULL) );
5182 SCIP_CALL( SCIPpqueueInsert((*nodepairqueue)->pqueue, (void*)&(*nodepairqueue)->nodepairs[n]) );
5277 SCIP_CALL( SCIPallocBufferArray(scip, &(*nodepartition)->representatives, mcfnetwork->nnodes) );
5319 SCIPdebugMsg(scip, "shrinking nodepair (%d,%d) with weight %g: join representatives %d and %d\n",
5328 /* if there have been too few arcs to shrink the graph to the required number of clusters, join clusters with first cluster
5491 if( nodeInPartition(nodepartition, partition, FALSE, s) == nodeInPartition(nodepartition, partition, FALSE, t) )
5599 (void) SCIPsnprintf(label, SCIP_MAXSTRLEN, "%s", SCIProwGetName(mcfnetwork->nodeflowrows[v][0]));
5724 (void) SCIPsnprintf(label, SCIP_MAXSTRLEN, "%s", SCIProwGetName(mcfnetwork->arccapacityrows[a]));
5756 fprintf(file, " source %d\n", mcfnetwork->arcsources[a] >= 0 ? mcfnetwork->arcsources[a] : mcfnetwork->nnodes);
5757 fprintf(file, " target %d\n", mcfnetwork->arctargets[a] >= 0 ? mcfnetwork->arctargets[a] : mcfnetwork->nnodes);
5833 (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "mcf%" SCIP_LONGINT_FORMAT "_%d", SCIPgetNLPs(scip), *ncuts);
5834 SCIP_CALL( SCIPcreateEmptyRowSepa(scip, &cut, sepa, cutname, -SCIPinfinity(scip), cutrhs, cutislocal, FALSE,
5846 cutname, cutrhs, SCIPgetRowSolActivity(scip, cut, sol), SCIPgetCutEfficacy(scip, sol, cut), SCIProwGetRank(cut));
5865 SCIP_CALL( SCIPseparateRelaxedKnapsack(scip, NULL, sepa, cutnnz, cutvars, cutcoefs, +1.0, cutrhs, sol, cutoff, ncuts) );
5948 * (1) \sum_{a \in \delta^+(v)} f_a^k - \sum_{a \in \delta^-(v)} f_a^k <= -d_v^k for all v \in V and k \in K
5961 * Then, the base constraint of the c-MIR cut is the sum of those flow conservation constraints and the
5964 * sum_{k \in K^+} sum_{v \in S} (sum_{a \in \delta^+(v)} f_a^k - sum_{a \in \delta^-(v)} f_a^k) <= sum_{k \in K^+} sum_{v \in S} -d_v^k
5984 SCIPdebugMsg(scip, "maxtestdelta: %d, maxsepacuts: %d, nnodes: %d \n", maxtestdelta, maxsepacuts, nnodes);
5994 SCIPdebugMsg(scip, "maxtestdelta: %d, maxsepacuts: %d, nclusters: %d \n", maxtestdelta, maxsepacuts, nclusters);
5996 /* We fix the last cluster to belong to partition T. In the undirected case, this is sufficient,
5997 * because S-T is equivalent to T-S. In the directed case, the loop below is conducted two times,
6005 for( partition = startpartition; partition <= allpartitions-1 && !SCIPisStopped(scip) && *ncuts < maxsepacuts && !*cutoff; partition++ )
6021 if( nodepartition != NULL && !nodepartitionIsConnected(scip, mcfnetwork, nodepartition, partition ) )
6036 SCIPdebugMsg(scip, "generating single-node cuts for node %u (inverted: %u)\n", partition, inverted);
6040 SCIPdebugMsg(scip, "generating cluster cuts for partition 0x%x (inverted: %u)\n", partition, inverted);
6094 if( sepadata->separatesinglenodecuts && nodepartition != NULL && (nnodesinS == 1 || nnodesinS == nnodes-1) )
6196 SCIPdebugMsg(scip, " -> arc %d, r=%d, capacity row <%s>: weight=%g slack=%g dual=%g\n", a, r, SCIProwGetName(arccapacityrows[a]), rowweights[r],
6203 baserhs += rowweights[r] * (SCIProwGetRhs(arccapacityrows[a]) - SCIProwGetConstant(arccapacityrows[a]));
6205 baserhs += rowweights[r] * (SCIProwGetLhs(arccapacityrows[a]) - SCIProwGetConstant(arccapacityrows[a]));
6208 /* extract useful deltas for c-MIR scaling and update the demand value for commodities (in binary flow model) */
6229 if( !SCIPisZero(scip, 1.0 / coef) && !SCIPisFeasZero(scip, coef) && SCIPcolIsIntegral(rowcols[j]) )
6286 /* if commodity was not hit by the capacity constraints of the cut in the graph, ignore the commodity */
6350 SCIPdebugMsg(scip, " -> node %d, commodity %d, r=%d, flow row <%s>: scale=%g weight=%g slack=%g dual=%g\n",
6357 baserhs += rowweights[r] * (SCIProwGetRhs(nodeflowrows[v][k]) - SCIProwGetConstant(nodeflowrows[v][k]));
6359 baserhs += rowweights[r] * (SCIProwGetLhs(nodeflowrows[v][k]) - SCIProwGetConstant(nodeflowrows[v][k]));
6401 SCIP_CALL( SCIPcalcMIR(scip, sol, POSTPROCESS, BOUNDSWITCH, USEVBDS, allowlocal, sepadata->fixintegralrhs, NULL, NULL, MINFRAC, MAXFRAC,
6402 1.0/deltas[d], aggrrow, cutcoefs, &cutrhs, cutinds, &cutnnz, &cutefficacy, &cutrank, &cutislocal, &success) );
6422 SCIPdebugMsg(scip, "success -> delta = %g -> rhs: %g, efficacy: %g\n", deltas[d], cutrhs, cutefficacy);
6423 SCIP_CALL( addCut(scip, sepa, sepadata, sol, cutcoefs, cutrhs, cutinds, cutnnz, cutislocal, cutrank, ncuts, cutoff) );
6441 SCIPgetRowFeasibility(scip, arccapacityrows[a]) * arccapacityscales[a], SCIProwGetDualsol(arccapacityrows[a]));
6560 SCIPdebugMsg(scip, " -> discarding capacity row <%s> of weight %g and slack %g: increases MIR violation by %g\n",
6561 SCIProwGetName(arccapacityrows[a]), rowweights[r], SCIPgetRowFeasibility(scip, arccapacityrows[a]),
6582 SCIPdebugMsg(scip, "applying MIR with delta = %g to flowcut inequality (violation improvement: %g)\n", bestdelta, totalviolationdelta);
6589 SCIP_CALL( SCIPcalcMIR(scip, sol, POSTPROCESS, BOUNDSWITCH, USEVBDS, allowlocal, sepadata->fixintegralrhs, NULL, NULL, MINFRAC, MAXFRAC,
6590 1.0/bestdelta, aggrrow, cutcoefs, &cutrhs, cutinds, &cutnnz, &cutefficacy, &cutrank, &cutislocal, &success) );
6596 SCIPdebugMsg(scip, " -> delta = %g -> rhs: %g, efficacy: %g\n", bestdelta, cutrhs, cutefficacy);
6597 SCIP_CALL( addCut(scip, sepa, sepadata, sol, cutcoefs, cutrhs, cutinds, cutnnz, cutislocal, cutrank, ncuts, cutoff) );
6653 /* Notice that in principle we should have most of the columns in the LP to be able to detect a structure */
6675 /* if separation was not delayed before and we had no success in previous round then delay the separation*/
6684 MCFdebugMessage("%d columns, %d rows, ratio %g is not in [%g,%g] -> exit\n", ncols, nrows, colrowratio, MINCOLROWRATIO, MAXCOLROWRATIO);
6695 SCIP_CALL( mcfnetworkExtract(scip, sepadata, &sepadata->mcfnetworks, &sepadata->nmcfnetworks, &sepadata->effortlevel) );
6701 MCFdebugMessage(" -> extracted network %d has %d nodes, %d (%d) arcs (uncapacitated), and %d commodities (modeltype: %s)\n",
6702 i, sepadata->mcfnetworks[i]->nnodes, sepadata->mcfnetworks[i]->narcs, sepadata->mcfnetworks[i]->nuncapacitatedarcs,
6704 sepadata->mcfnetworks[i]->modeltype == SCIP_MCFMODELTYPE_DIRECTED ? "directed" : "undirected");
6739 /* do not allow networks exceeding the arcs/nodes ratio ( = average node degree / 2 (directed)) */
6750 SCIP_CALL( generateClusterCuts(scip, sepa, sepadata, sol, allowlocal, depth, mcfnetwork, NULL, &ncuts, &cutoff) );
6757 sepadata->effortlevel == MCFEFFORTLEVEL_DEFAULT ? sepadata->nclusters : 2 * sepadata->nclusters) );
6763 SCIP_CALL( generateClusterCuts(scip, sepa, sepadata, sol, allowlocal, depth, mcfnetwork, nodepartition, &ncuts, &cutoff) );
6769 MCFdebugMessage("MCF network has %d nodes, %d arcs, %d commodities. Found %d MCF network cuts, cutoff = %u.\n",
6828/** solving process initialization method of separator (called when branch and bound process is about to begin) */
6846/** solving process deinitialization method of separator (called before branch and bound process data is freed) */
6933 SCIP_CALL( SCIPincludeSepaBasic(scip, &sepa, SEPA_NAME, SEPA_DESC, SEPA_PRIORITY, SEPA_FREQ, SEPA_MAXBOUNDDIST,
6951 &sepadata->nclusters, TRUE, DEFAULT_NCLUSTERS, 2, (int) (8*sizeof(unsigned int)), NULL, NULL) );
6982 "maximal number of mcf cuts separated per separation round in the root node -- default separation",
6987 &sepadata->maxinconsistencyratio, TRUE, DEFAULT_MAXINCONSISTENCYRATIO, 0.0, SCIP_REAL_MAX, NULL, NULL) );
6991 &sepadata->maxarcinconsistencyratio, TRUE, DEFAULT_MAXARCINCONSISTENCYRATIO, 0.0, SCIP_REAL_MAX, NULL, NULL) );
Constraint handler for knapsack constraints of the form , x binary and .
methods for the aggregation rows
SCIP_RETCODE SCIPseparateRelaxedKnapsack(SCIP *scip, SCIP_CONS *cons, SCIP_SEPA *sepa, int nknapvars, SCIP_VAR **knapvars, SCIP_Real *knapvals, SCIP_Real valscale, SCIP_Real rhs, SCIP_SOL *sol, SCIP_Bool *cutoff, int *ncuts)
Definition: cons_knapsack.c:5777
SCIP_RETCODE SCIPgetVarsData(SCIP *scip, SCIP_VAR ***vars, int *nvars, int *nbinvars, int *nintvars, int *nimplvars, int *ncontvars)
Definition: scip_prob.c:1866
SCIP_RETCODE SCIPhashtableCreate(SCIP_HASHTABLE **hashtable, BMS_BLKMEM *blkmem, int tablesize, SCIP_DECL_HASHGETKEY((*hashgetkey)), SCIP_DECL_HASHKEYEQ((*hashkeyeq)), SCIP_DECL_HASHKEYVAL((*hashkeyval)), void *userptr)
Definition: misc.c:2296
void * SCIPhashtableRetrieve(SCIP_HASHTABLE *hashtable, void *key)
Definition: misc.c:2608
SCIP_RETCODE SCIPhashtableInsert(SCIP_HASHTABLE *hashtable, void *element)
Definition: misc.c:2547
void SCIPinfoMessage(SCIP *scip, FILE *file, const char *formatstr,...)
Definition: scip_message.c:208
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
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
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
SCIP_RETCODE SCIPpqueueInsert(SCIP_PQUEUE *pqueue, void *elem)
Definition: misc.c:1394
SCIP_RETCODE SCIPpqueueCreate(SCIP_PQUEUE **pqueue, int initsize, SCIP_Real sizefac, SCIP_DECL_SORTPTRCOMP((*ptrcomp)), SCIP_DECL_PQUEUEELEMCHGPOS((*elemchgpos)))
Definition: misc.c:1295
SCIP_Real SCIPgetCutEfficacy(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut)
Definition: scip_cut.c:94
SCIP_RETCODE SCIPaggrRowCreate(SCIP *scip, SCIP_AGGRROW **aggrrow)
Definition: cuts.c:1723
SCIP_Bool SCIPisCutEfficacious(SCIP *scip, SCIP_SOL *sol, SCIP_ROW *cut)
Definition: scip_cut.c:117
SCIP_Bool SCIPisEfficacious(SCIP *scip, SCIP_Real efficacy)
Definition: scip_cut.c:135
SCIP_RETCODE SCIPaddRow(SCIP *scip, SCIP_ROW *row, SCIP_Bool forcecut, SCIP_Bool *infeasible)
Definition: scip_cut.c:250
SCIP_RETCODE SCIPaggrRowSumRows(SCIP *scip, SCIP_AGGRROW *aggrrow, SCIP_Real *weights, int *rowinds, int nrowinds, SCIP_Bool sidetypebasis, SCIP_Bool allowlocal, int negslack, int maxaggrlen, SCIP_Bool *valid)
Definition: cuts.c:2279
SCIP_RETCODE SCIPcalcMIR(SCIP *scip, SCIP_SOL *sol, SCIP_Bool postprocess, SCIP_Real boundswitch, SCIP_Bool usevbds, SCIP_Bool allowlocal, SCIP_Bool fixintegralrhs, int *boundsfortrans, SCIP_BOUNDTYPE *boundtypesfortrans, SCIP_Real minfrac, SCIP_Real maxfrac, SCIP_Real scale, SCIP_AGGRROW *aggrrow, SCIP_Real *cutcoefs, SCIP_Real *cutrhs, int *cutinds, int *cutnnz, SCIP_Real *cutefficacy, int *cutrank, SCIP_Bool *cutislocal, SCIP_Bool *success)
Definition: cuts.c:3873
SCIP_RETCODE SCIPgetLPColsData(SCIP *scip, SCIP_COL ***cols, int *ncols)
Definition: scip_lp.c:471
SCIP_RETCODE SCIPgetLPRowsData(SCIP *scip, SCIP_ROW ***rows, int *nrows)
Definition: scip_lp.c:570
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:111
SCIP_Real SCIPgetRowFeasibility(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:2124
SCIP_Real SCIPgetRowSolFeasibility(SCIP *scip, SCIP_ROW *row, SCIP_SOL *sol)
Definition: scip_lp.c:2167
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 SCIPgetRowLPFeasibility(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:2010
SCIP_RETCODE SCIPaddVarsToRow(SCIP *scip, SCIP_ROW *row, int nvars, SCIP_VAR **vars, SCIP_Real *vals)
Definition: scip_lp.c:1727
SCIP_Real SCIPgetRowSolActivity(SCIP *scip, SCIP_ROW *row, SCIP_SOL *sol)
Definition: scip_lp.c:2144
SCIP_RETCODE SCIPsetSepaInitsol(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAINITSOL((*sepainitsol)))
Definition: scip_sepa.c:215
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
SCIP_RETCODE SCIPsetSepaFree(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAFREE((*sepafree)))
Definition: scip_sepa.c:167
SCIP_RETCODE SCIPsetSepaExitsol(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAEXITSOL((*sepaexitsol)))
Definition: scip_sepa.c:231
SCIP_RETCODE SCIPsetSepaCopy(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPACOPY((*sepacopy)))
Definition: scip_sepa.c:151
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1217
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:497
SCIP_Bool SCIPisFeasIntegral(SCIP *scip, SCIP_Real val)
Definition: scip_numerics.c:881
SCIP_Bool SCIPisFeasGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:819
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:445
SCIP_Bool SCIPisFeasPositive(SCIP *scip, SCIP_Real val)
Definition: scip_numerics.c:857
void SCIPsortInd(int *indarray, SCIP_DECL_SORTINDCOMP((*indcomp)), void *dataptr, int len)
void SCIPsortIntInt(int *intarray1, int *intarray2, int len)
memory allocation routines
Definition: objbenders.h:44
public methods for LP management
public methods for message output
public data structures and miscellaneous methods
methods for sorting joint arrays of various types
public methods for separators
public methods for problem variables
public methods for branching rule plugins and branching
public methods for cuts and aggregation rows
general public methods
public methods for the LP relaxation, rows and columns
public methods for memory management
public methods for message handling
public methods for numerical tolerances
public methods for SCIP parameter handling
public methods for global and local (sub)problems
public methods for separator plugins
public methods for solutions
public methods for querying solving statistics
public methods for the branch-and-bound tree
static void addFlowrowToCommodity(SCIP *scip, MCFDATA *mcfdata, SCIP_ROW *row, unsigned char rowsign, int *comcolids, int *ncomcolids)
Definition: sepa_mcf.c:1516
static SCIP_RETCODE extractFlow(SCIP *scip, MCFDATA *mcfdata, SCIP_Real maxflowvarflowrowratio, SCIP_Bool *failed)
Definition: sepa_mcf.c:2060
static SCIP_RETCODE nodepartitionCreate(SCIP *scip, SCIP_MCFNETWORK *mcfnetwork, NODEPARTITION **nodepartition, int nclusters)
Definition: sepa_mcf.c:5257
static SCIP_RETCODE mcfnetworkFree(SCIP *scip, SCIP_MCFNETWORK **mcfnetwork)
Definition: sepa_mcf.c:330
static SCIP_RETCODE mcfnetworkFill(SCIP *scip, SCIP_MCFNETWORK *mcfnetwork, MCFDATA *mcfdata, int *compnodeid, int *compnodes, int ncompnodes, int *comparcs, int ncomparcs)
Definition: sepa_mcf.c:381
static void getNextFlowrow(SCIP *scip, MCFDATA *mcfdata, SCIP_ROW **nextrow, unsigned char *nextrowsign, SCIP_Bool *nextinvertcommodity)
Definition: sepa_mcf.c:1943
static SCIP_RETCODE findUncapacitatedArcs(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:3103
static SCIP_RETCODE getNodeSimilarityScore(SCIP *scip, MCFDATA *mcfdata, int baserowlen, int *basearcpattern, int basenposuncap, int basenneguncap, SCIP_ROW *row, SCIP_Real *score, SCIP_Bool *invertcommodity)
Definition: sepa_mcf.c:2510
static void nodepartitionJoin(NODEPARTITION *nodepartition, int rep1, int rep2)
Definition: sepa_mcf.c:5246
static SCIP_RETCODE identifySourcesTargets(SCIP *scip, MCFDATA *mcfdata, SCIP_SEPADATA *sepadata, MCFEFFORTLEVEL *effortlevel)
Definition: sepa_mcf.c:3761
static SCIP_RETCODE mcfnetworkCreate(SCIP *scip, SCIP_MCFNETWORK **mcfnetwork)
Definition: sepa_mcf.c:304
static void nodepairqueueFree(SCIP *scip, NODEPAIRQUEUE **nodepairqueue)
Definition: sepa_mcf.c:5191
static SCIP_RETCODE mcfnetworkExtract(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_MCFNETWORK ***mcfnetworks, int *nmcfnetworks, MCFEFFORTLEVEL *effortlevel)
Definition: sepa_mcf.c:4241
static SCIP_RETCODE extractNodes(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:2711
static void unionfindJoinSets(int *representatives, int rep1, int rep2)
Definition: sepa_mcf.c:4748
static void deleteCommodity(SCIP *scip, MCFDATA *mcfdata, int k, SCIP_ROW **comrows, int nrows, int *ndelflowrows, int *ndelflowvars)
Definition: sepa_mcf.c:1727
static SCIP_RETCODE nodepairqueueCreate(SCIP *scip, SCIP_MCFNETWORK *mcfnetwork, NODEPAIRQUEUE **nodepairqueue)
Definition: sepa_mcf.c:4878
static int nodepartitionIsConnected(SCIP *scip, SCIP_MCFNETWORK *mcfnetwork, NODEPARTITION *nodepartition, unsigned int partition)
Definition: sepa_mcf.c:5452
#define DEFAULT_CHECKCUTSHORECONNECTIVITY
Definition: sepa_mcf.c:104
static SCIP_RETCODE extractFlowRows(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:853
static void getFlowrowFit(SCIP *scip, MCFDATA *mcfdata, SCIP_ROW *row, int k, unsigned char *rowsign, SCIP_Bool *invertcommodity)
Definition: sepa_mcf.c:1811
static void getIncidentNodes(SCIP *scip, MCFDATA *mcfdata, SCIP_COL *col, int *sourcenode, int *targetnode)
Definition: sepa_mcf.c:3004
static SCIP_Bool nodeInPartition(NODEPARTITION *nodepartition, unsigned int partition, SCIP_Bool inverted, int v)
Definition: sepa_mcf.c:5422
static SCIP_RETCODE extractCapacities(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:2254
static SCIP_RETCODE createNewArc(SCIP *scip, MCFDATA *mcfdata, int source, int target, int *newarcid)
Definition: sepa_mcf.c:1463
static SCIP_RETCODE identifyComponent(SCIP *scip, MCFDATA *mcfdata, int *nodevisited, int startv, int *compnodes, int *ncompnodes, int *comparcs, int *ncomparcs)
Definition: sepa_mcf.c:4110
static SCIP_RETCODE generateClusterCuts(SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_SOL *sol, SCIP_Bool allowlocal, int depth, SCIP_MCFNETWORK *mcfnetwork, NODEPARTITION *nodepartition, int *ncuts, SCIP_Bool *cutoff)
Definition: sepa_mcf.c:5876
static SCIP_RETCODE cleanupNetwork(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:3443
static void nodepartitionFree(SCIP *scip, NODEPARTITION **nodepartition)
Definition: sepa_mcf.c:5405
static SCIP_RETCODE extractCapacityRows(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:1069
static SCIP_Bool nodepairqueueIsEmpty(NODEPAIRQUEUE *nodepairqueue)
Definition: sepa_mcf.c:5207
static int nodepartitionGetRepresentative(NODEPARTITION *nodepartition, int v)
Definition: sepa_mcf.c:5236
static SCIP_RETCODE addCut(SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_SOL *sol, SCIP_Real *cutcoefs, SCIP_Real cutrhs, int *cutinds, int cutnnz, SCIP_Bool cutislocal, int cutrank, int *ncuts, SCIP_Bool *cutoff)
Definition: sepa_mcf.c:5789
static void unionfindInitSets(int *representatives, int nelems)
Definition: sepa_mcf.c:4716
static void collectIncidentFlowCols(SCIP *scip, MCFDATA *mcfdata, SCIP_ROW *flowrow, int basecommodity)
Definition: sepa_mcf.c:2425
static void invertCommodity(SCIP *scip, MCFDATA *mcfdata, int k, SCIP_ROW **comrows, int ncomrows, int *comcolids, int ncomcolids)
Definition: sepa_mcf.c:1662
static NODEPAIRENTRY * nodepairqueueRemove(NODEPAIRQUEUE *nodepairqueue)
Definition: sepa_mcf.c:5219
static SCIP_RETCODE separateCuts(SCIP *scip, SCIP_SEPA *sepa, SCIP_SOL *sol, SCIP_Bool allowlocal, int depth, SCIP_RESULT *result)
Definition: sepa_mcf.c:6620
static SCIP_RETCODE createNewCommodity(SCIP *scip, MCFDATA *mcfdata)
Definition: sepa_mcf.c:1439
static int unionfindGetRepresentative(int *representatives, int v)
Definition: sepa_mcf.c:4730
multi-commodity-flow network cut separator
Definition: struct_cuts.h:41
Definition: struct_lp.h:136
Definition: struct_misc.h:90
Definition: sepa_mcf.c:172
Definition: struct_misc.h:79
Definition: struct_lp.h:202
Definition: struct_sepa.h:47
Definition: struct_sol.h:74
Definition: struct_var.h:208
Definition: struct_scip.h:70