Scippy

SCIP

Solving Constraint Integer Programs

stptest_graphpath.c
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2 /* */
3 /* This file is part of the program and library */
4 /* SCIP --- Solving Constraint Integer Programs */
5 /* */
6 /* Copyright (C) 2002-2022 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not visit scip.zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file stptest_graphpath.c
17  * @brief graph utility tests for Steiner tree problem methods
18  * @author Daniel Rehfeldt
19  *
20  * This file implements graph utility tests for Steiner tree problems.
21  *
22  * A list of all interface methods can be found in stptest.h.
23  *
24  */
25 
26 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
27 //#define SCIP_DEBUG
28 #include "scip/scip.h"
29 #include "stptest.h"
30 #include "graph.h"
31 #include "portab.h"
32 
33 
34 
35 /** tests that (unbiased) terminal paths are found */
36 static
38  SCIP* scip /**< SCIP data structure */
39 )
40 {
41  SCIP_Real dists[4];
42  int closeterms[4];
43  TPATHS* tpaths;
44  GRAPH* graph;
45  int nnodes = 6;
46  int nedges = 10;
47  int ncloseterms = -1;
48 
49  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
50 
51  for( int i = 0; i < nnodes; i++ )
53 
54  graph->source = 2;
55  graph_knot_chg(graph, 2, STP_TERM);
56  graph_knot_chg(graph, 3, STP_TERM);
57  graph_knot_chg(graph, 4, STP_TERM);
58  graph_knot_chg(graph, 5, STP_TERM);
59 
60  graph_edge_addBi(scip, graph, 0, 1, 1.0); // 0
61  graph_edge_addBi(scip, graph, 0, 3, 2.1); // 2
62  graph_edge_addBi(scip, graph, 0, 4, 1.6);
63  graph_edge_addBi(scip, graph, 1, 2, 1.0);
64  graph_edge_addBi(scip, graph, 4, 5, 1.0);
65 
66  SCIP_CALL( stptest_graphSetUp(scip, graph) );
67  SCIP_CALL( graph_tpathsInit(scip, graph, &tpaths) );
68  graph_tpathsSetAll4(graph, graph->cost, graph->cost, NULL, tpaths);
69 
70  graph_tpathsGet4CloseTerms(graph, tpaths, 0, FARAWAY, closeterms, NULL, dists, &ncloseterms);
71  STPTEST_ASSERT(ncloseterms == 3);
72  STPTEST_ASSERT(closeterms[0] == 4);
73  STPTEST_ASSERT(closeterms[1] == 2);
74  STPTEST_ASSERT(closeterms[2] == 3);
75  STPTEST_ASSERT(EQ(dists[0], 1.6));
76  STPTEST_ASSERT(EQ(dists[1], 2.0));
77  STPTEST_ASSERT(EQ(dists[2], 2.1));
78 
79  graph_tpathsGet4CloseTerms(graph, tpaths, 1, FARAWAY, closeterms, NULL, dists, &ncloseterms);
80  STPTEST_ASSERT(ncloseterms == 3);
81  STPTEST_ASSERT(closeterms[0] == 2);
82  STPTEST_ASSERT(closeterms[1] == 4);
83  STPTEST_ASSERT(closeterms[2] == 3);
84  STPTEST_ASSERT(EQ(dists[0], 1.0));
85  STPTEST_ASSERT(EQ(dists[1], 2.6));
86  STPTEST_ASSERT(EQ(dists[2], 3.1));
87 
88  stptest_graphTearDown(scip, graph);
89  graph_tpathsFree(scip, &tpaths);
90 
91  return SCIP_OKAY;
92 }
93 
94 
95 
96 /** tests that (unbiased) terminal paths are repaired properly (level 1) */
97 static
99  SCIP* scip /**< SCIP data structure */
100 )
101 {
102  SCIP_Real dists[4];
103  int closeterms[4];
104  TPATHS* tpaths;
105  GRAPH* graph;
106  int nnodes = 5;
107  int nedges = 10;
108  int ncloseterms = -1;
109 
110  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
111 
112  for( int i = 0; i < nnodes; i++ )
114 
115  graph->source = 0;
116  graph_knot_chg(graph, 0, STP_TERM);
117  graph_knot_chg(graph, 2, STP_TERM);
118  graph_knot_chg(graph, 4, STP_TERM);
119 
120  graph_edge_addBi(scip, graph, 0, 1, 1.1); // 0
121  graph_edge_addBi(scip, graph, 1, 2, 2.1); // 2
122  graph_edge_addBi(scip, graph, 1, 3, 1.0); //
123  graph_edge_addBi(scip, graph, 3, 4, 2.2); //
124  graph_edge_addBi(scip, graph, 0, 3, 2.15); //
125 
126  SCIP_CALL( stptest_graphSetUp(scip, graph) );
127  SCIP_CALL( graph_tpathsInit(scip, graph, &tpaths) );
128  graph_tpathsSetAll4(graph, graph->cost, graph->cost, NULL, tpaths);
129  SCIP_CALL( graph_tpathsRepairSetUp(graph, tpaths) );
130 
131  graph_tpathsGet4CloseTerms(graph, tpaths, 1, FARAWAY, closeterms, NULL, dists, &ncloseterms);
132  STPTEST_ASSERT(ncloseterms == 3);
133  STPTEST_ASSERT(closeterms[0] == 0);
134  STPTEST_ASSERT(EQ(dists[0], 1.1));
135 
136  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
137  STPTEST_ASSERT(closeterms[0] == 0);
138  STPTEST_ASSERT(EQ(dists[0], 2.1));
139 
140  SCIP_CALL( graph_tpathsRepair(scip, 0, FALSE, graph, tpaths) );
141 
142  graph_tpathsGet4CloseTerms(graph, tpaths, 1, FARAWAY, closeterms, NULL, dists, &ncloseterms);
143 // STPTEST_ASSERT(ncloseterms == 1);
144  STPTEST_ASSERT(closeterms[0] == 2);
145  STPTEST_ASSERT(EQ(dists[0], 2.1));
146 
147  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
148  STPTEST_ASSERT(closeterms[0] == 0);
149  STPTEST_ASSERT(EQ(dists[0], 2.15));
150 
151  stptest_graphTearDown(scip, graph);
152  graph_tpathsFree(scip, &tpaths);
153 
154  return SCIP_OKAY;
155 }
156 
157 
158 /** tests that (unbiased) terminal paths are repaired properly (level 1 + 2) */
159 static
161  SCIP* scip /**< SCIP data structure */
162 )
163 {
164  SCIP_Real dists[4];
165  int closeterms[4];
166  TPATHS* tpaths;
167  GRAPH* graph;
168  int nnodes = 5;
169  int nedges = 10;
170  int ncloseterms = -1;
171 
172  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
173 
174  for( int i = 0; i < nnodes; i++ )
176 
177  graph->source = 0;
178  graph_knot_chg(graph, 0, STP_TERM);
179  graph_knot_chg(graph, 2, STP_TERM);
180  graph_knot_chg(graph, 4, STP_TERM);
181 
182 
183  graph_edge_addBi(scip, graph, 0, 1, 1.1); // 0
184  graph_edge_addBi(scip, graph, 1, 2, 2.0); // 2
185  graph_edge_addBi(scip, graph, 1, 3, 1.0); // 4
186  graph_edge_addBi(scip, graph, 3, 4, 2.05); // 6
187  graph_edge_addBi(scip, graph, 0, 3, 2.15); // 8
188 
189  SCIP_CALL( stptest_graphSetUp(scip, graph) );
190  SCIP_CALL( graph_tpathsInit(scip, graph, &tpaths) );
191  graph_tpathsSetAll4(graph, graph->cost, graph->cost, NULL, tpaths);
192  SCIP_CALL( graph_tpathsRepairSetUp(graph, tpaths) );
193 
194  graph_tpathsGet4CloseTerms(graph, tpaths, 1, FARAWAY, closeterms, NULL, dists, &ncloseterms);
195  STPTEST_ASSERT(ncloseterms == 3);
196  STPTEST_ASSERT(closeterms[0] == 0);
197  STPTEST_ASSERT(EQ(dists[0], 1.1));
198  STPTEST_ASSERT(closeterms[1] == 2);
199  STPTEST_ASSERT(EQ(dists[1], 2.0));
200 
201  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
202  STPTEST_ASSERT(closeterms[0] == 4);
203  STPTEST_ASSERT(EQ(dists[0], 2.05));
204  STPTEST_ASSERT(closeterms[1] == 0);
205  STPTEST_ASSERT(EQ(dists[1], 2.1));
206 
207  /* delete edge 0-1 */
208  SCIP_CALL( graph_tpathsRepair(scip, 0, FALSE, graph, tpaths) );
209 
210  graph_tpathsGet4CloseTerms(graph, tpaths, 1, FARAWAY, closeterms, NULL, dists, &ncloseterms);
211  STPTEST_ASSERT(closeterms[0] == 2);
212  STPTEST_ASSERT(EQ(dists[0], 2.0));
213 
214  STPTEST_ASSERT(closeterms[1] == 4);
215  STPTEST_ASSERT(EQ(dists[1], 3.05));
216 
217  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
218  STPTEST_ASSERT(closeterms[0] == 4);
219  STPTEST_ASSERT(EQ(dists[0], 2.05));
220  STPTEST_ASSERT(closeterms[1] == 0);
221  STPTEST_ASSERT(EQ(dists[1], 2.15));
222 
223  graph_edge_del(scip, graph, 0, TRUE);
224 
225  /* delete edge 0-3 */
226  SCIP_CALL( graph_tpathsRepair(scip, 8, FALSE, graph, tpaths) );
227 
228  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
229  STPTEST_ASSERT(closeterms[0] == 4);
230  STPTEST_ASSERT(EQ(dists[0], 2.05));
231  STPTEST_ASSERT(closeterms[1] == 2);
232  STPTEST_ASSERT(EQ(dists[1], 3.0));
233 
234  stptest_graphTearDown(scip, graph);
235  graph_tpathsFree(scip, &tpaths);
236 
237  return SCIP_OKAY;
238 }
239 
240 
241 /** tests that (unbiased) terminal paths are repaired properly (level 1 + 2 + 3) */
242 static
244  SCIP* scip /**< SCIP data structure */
245 )
246 {
247  SCIP_Real dists[4];
248  int closeterms[4];
249  TPATHS* tpaths;
250  GRAPH* graph;
251  int nnodes = 6;
252  int nedges = 16;
253  int ncloseterms = -1;
254 
255  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
256 
257  for( int i = 0; i < nnodes; i++ )
259 
260  graph->source = 0;
261  graph_knot_chg(graph, 0, STP_TERM);
262  graph_knot_chg(graph, 1, STP_TERM);
263  graph_knot_chg(graph, 2, STP_TERM);
264 
265  graph_edge_addBi(scip, graph, 0, 3, 1.0); // 0
266  graph_edge_addBi(scip, graph, 0, 4, 1.0); // 2
267  graph_edge_addBi(scip, graph, 1, 5, 1.1); // 4
268  graph_edge_addBi(scip, graph, 2, 4, 1.2); // 6
269  graph_edge_addBi(scip, graph, 3, 4, 1.0); // 8
270  graph_edge_addBi(scip, graph, 3, 5, 1.0); // 10
271  graph_edge_addBi(scip, graph, 4, 5, 1.2); // 12
272  graph_edge_addBi(scip, graph, 2, 5, 2.1); // 12
273 
274 
275  SCIP_CALL( stptest_graphSetUp(scip, graph) );
276  SCIP_CALL( graph_tpathsInit(scip, graph, &tpaths) );
277  graph_tpathsSetAll4(graph, graph->cost, graph->cost, NULL, tpaths);
278  SCIP_CALL( graph_tpathsRepairSetUp(graph, tpaths) );
279 
280  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
281  STPTEST_ASSERT(ncloseterms == 3);
282  STPTEST_ASSERT(closeterms[0] == 0);
283  STPTEST_ASSERT(EQ(dists[0], 1.0));
284  STPTEST_ASSERT(closeterms[1] == 1);
285  STPTEST_ASSERT(EQ(dists[1], 2.1));
286 
287  graph_tpathsGet4CloseTerms(graph, tpaths, 5, FARAWAY, closeterms, NULL, dists, &ncloseterms);
288  STPTEST_ASSERT(closeterms[0] == 1);
289  STPTEST_ASSERT(EQ(dists[0], 1.1));
290  STPTEST_ASSERT(closeterms[1] == 0);
291  STPTEST_ASSERT(EQ(dists[1], 2.0));
292 
293  /* delete edge 0-3 */
294  SCIP_CALL( graph_tpathsRepair(scip, 0, FALSE, graph, tpaths) );
295 
296  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
297  STPTEST_ASSERT(closeterms[0] == 0);
298  STPTEST_ASSERT(EQ(dists[0], 2.0));
299  STPTEST_ASSERT(closeterms[1] == 1);
300  STPTEST_ASSERT(EQ(dists[1], 2.1));
301 
302  graph_tpathsGet4CloseTerms(graph, tpaths, 5, FARAWAY, closeterms, NULL, dists, &ncloseterms);
303  STPTEST_ASSERT(closeterms[0] == 1);
304  STPTEST_ASSERT(EQ(dists[0], 1.1));
305  STPTEST_ASSERT(closeterms[1] == 2);
306  STPTEST_ASSERT(EQ(dists[1], 2.1));
307  STPTEST_ASSERT(closeterms[2] == 0);
308  STPTEST_ASSERT(EQ(dists[2], 2.2));
309 
310  graph_tpathsGet4CloseTerms(graph, tpaths, 4, FARAWAY, closeterms, NULL, dists, &ncloseterms);
311  STPTEST_ASSERT(closeterms[0] == 0);
312  STPTEST_ASSERT(EQ(dists[0], 1.0));
313  STPTEST_ASSERT(closeterms[1] == 2);
314  STPTEST_ASSERT(EQ(dists[1], 1.2));
315  STPTEST_ASSERT(closeterms[2] == 1);
316  STPTEST_ASSERT(EQ(dists[2], 2.3));
317 
318  graph_edge_del(scip, graph, 0, TRUE);
319 
320  /* delete edge 0-4 */
321  SCIP_CALL( graph_tpathsRepair(scip, 2, FALSE, graph, tpaths) );
322 
323  graph_tpathsGet4CloseTerms(graph, tpaths, 4, FARAWAY, closeterms, NULL, dists, &ncloseterms);
324  assert(ncloseterms == 2);
325  STPTEST_ASSERT(closeterms[0] == 2);
326  STPTEST_ASSERT(EQ(dists[0], 1.2));
327  STPTEST_ASSERT(closeterms[1] == 1);
328  STPTEST_ASSERT(EQ(dists[1], 2.3));
329 
330 
331  graph_tpathsGet4CloseTerms(graph, tpaths, 3, FARAWAY, closeterms, NULL, dists, &ncloseterms);
332  STPTEST_ASSERT(closeterms[0] == 1);
333  STPTEST_ASSERT(EQ(dists[0], 2.1));
334  STPTEST_ASSERT(closeterms[1] == 2);
335  STPTEST_ASSERT(EQ(dists[1], 2.2));
336 
337  stptest_graphTearDown(scip, graph);
338  graph_tpathsFree(scip, &tpaths);
339 
340  return SCIP_OKAY;
341 }
342 
343 
344 /** tests that SD is repaired properly */
345 static
347  SCIP* scip /**< SCIP data structure */
348 )
349 {
350  SD* sd;
351  GRAPH* graph;
352  int nnodes = 7;
353  int nedges = 22;
354 
355  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
356 
357  for( int i = 0; i < nnodes; i++ )
359 
360  graph->source = 4;
361  graph_knot_chg(graph, 4, STP_TERM);
362  graph_knot_chg(graph, 5, STP_TERM);
363  graph_knot_chg(graph, 6, STP_TERM);
364 
365  graph_edge_addBi(scip, graph, 0, 4, 1.0); // 0
366  graph_edge_addBi(scip, graph, 0, 5, 3.0); // 2
367  graph_edge_addBi(scip, graph, 1, 4, 1.1); // 4
368  graph_edge_addBi(scip, graph, 2, 5, 2.2); // 6
369  graph_edge_addBi(scip, graph, 3, 5, 2.3); // 8
370  graph_edge_addBi(scip, graph, 4, 5, 1.0); // 10
371  graph_edge_addBi(scip, graph, 1, 5, 1.3); // 12
372  graph_edge_addBi(scip, graph, 2, 6, 2.2); // 14
373  graph_edge_addBi(scip, graph, 0, 6, 2.4); // 16
374  graph_edge_addBi(scip, graph, 5, 6, 2.0); // 18
375  graph_edge_addBi(scip, graph, 1, 0, 2.0); // 20
376 
377  SCIP_CALL( stptest_graphSetUp(scip, graph) );
378  SCIP_CALL( reduce_sdInit(scip, graph, &sd) );
379  SCIP_CALL( reduce_sdRepairSetUp(scip, graph, sd) );
380 
381 
382  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 0, 3, FARAWAY, 0.0, sd), 2.3));
383  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 1, 2, FARAWAY, 0.0, sd), 2.2));
384 
385  SCIP_CALL( reduce_sdRepair(scip, 0, FALSE, graph, sd) );
386  graph_edge_del(scip, graph, 0, TRUE);
387 
388  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 0, 3, FARAWAY, 0.0, sd), 2.4));
389  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 1, 2, FARAWAY, 0.0, sd), 2.2));
390  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 0, 3, FARAWAY, 0.0, sd), 2.4));
391 
392  SCIP_CALL( reduce_sdRepair(scip, 4, FALSE, graph, sd) );
393  graph_edge_del(scip, graph, 4, TRUE);
394  SCIP_CALL( reduce_sdRepair(scip, 12, FALSE, graph, sd) );
395  graph_edge_del(scip, graph, 12, TRUE);
396 
397  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 1, 2, FARAWAY, 0.0, sd), 4.4));
398  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 3, 4, FARAWAY, 0.0, sd), 2.3));
399 
400  SCIP_CALL( reduce_sdRepair(scip, 8, FALSE, graph, sd) );
401  graph_edge_del(scip, graph, 8, TRUE);
402  STPTEST_ASSERT(EQ(reduce_sdGetSd(graph, 3, 4, FARAWAY, 0.0, sd), FARAWAY));
403 
404  reduce_sdFree(scip, &sd);
405  stptest_graphTearDown(scip, graph);
406 
407  return SCIP_OKAY;
408 }
409 
410 /** tests that (unbiased) terminal paths are found */
411 static
413  SCIP* scip /**< SCIP data structure */
414 )
415 {
416  SCIP_Real dists[4];
417  int closeterms[4];
418  TPATHS* tpaths;
419  SDPROFIT* sdprofit;
420  GRAPH* graph;
421  int nnodes = 7;
422  int nedges = 16;
423  int ncloseterms = -1;
424 
425  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
426 
427  for( int i = 0; i < nnodes; i++ )
429 
430  graph->source = 2;
431  graph_knot_chg(graph, 2, STP_TERM);
432  graph_knot_chg(graph, 3, STP_TERM);
433  graph_knot_chg(graph, 4, STP_TERM);
434  graph_knot_chg(graph, 5, STP_TERM);
435  graph_knot_chg(graph, 6, STP_TERM);
436 
437  graph_edge_addBi(scip, graph, 0, 1, 1.0); // 0
438  graph_edge_addBi(scip, graph, 0, 3, 2.1); // 2
439  graph_edge_addBi(scip, graph, 0, 4, 1.7);
440  graph_edge_addBi(scip, graph, 1, 2, 1.1);
441  graph_edge_addBi(scip, graph, 1, 6, 1.0);
442  graph_edge_addBi(scip, graph, 4, 5, 1.8);
443  graph_edge_addBi(scip, graph, 5, 6, 1.5);
444 
445  SCIP_CALL( stptest_graphSetUp(scip, graph) );
446  SCIP_CALL( reduce_sdprofitInit(scip, graph, &sdprofit) );
447  SCIP_CALL( graph_tpathsInit(scip, graph, &tpaths) );
448  graph_tpathsSetAll4(graph, graph->cost, graph->cost, sdprofit, tpaths);
449 
450  graph_tpathsGet4CloseTerms(graph, tpaths, 0, FARAWAY, closeterms, NULL, dists, &ncloseterms);
451  STPTEST_ASSERT(ncloseterms == 4);
452  STPTEST_ASSERT(closeterms[0] == 6);
453  STPTEST_ASSERT(closeterms[1] == 2);
454  STPTEST_ASSERT(closeterms[2] == 4);
455  STPTEST_ASSERT(closeterms[3] == 3);
456  STPTEST_ASSERT(EQ(dists[0], 1.0));
457  STPTEST_ASSERT(EQ(dists[1], 1.6));
458  STPTEST_ASSERT(EQ(dists[2], 1.7));
459  STPTEST_ASSERT(EQ(dists[3], 2.1));
460 
461  graph_tpathsGet4CloseTerms(graph, tpaths, 1, FARAWAY, closeterms, NULL, dists, &ncloseterms);
462  STPTEST_ASSERT(ncloseterms == 4);
463  STPTEST_ASSERT(closeterms[2] == 4);
464  STPTEST_ASSERT(closeterms[3] == 3);
465  STPTEST_ASSERT(EQ(dists[0], 1.0));
466  STPTEST_ASSERT(EQ(dists[1], 1.1));
467  STPTEST_ASSERT(EQ(dists[2], 1.7));
468  STPTEST_ASSERT(EQ(dists[3], 3.0));
469 
470  stptest_graphTearDown(scip, graph);
471  graph_tpathsFree(scip, &tpaths);
472  reduce_sdprofitFree(scip, &sdprofit);
473 
474  return SCIP_OKAY;
475 }
476 
477 
478 /** tests terminal paths */
480  SCIP* scip /**< SCIP data structure */
481 )
482 {
483  assert(scip);
484 
485  SCIP_CALL( testSdRepair(scip) );
486 
490 
493 
494  printf("TPATHS test: all ok \n");
495 
496  return SCIP_OKAY;
497 }
SCIP_RETCODE reduce_sdInit(SCIP *, GRAPH *, SD **)
void graph_knot_chg(GRAPH *, int, int)
Definition: graph_node.c:86
static SCIP_RETCODE testTerminalPathsRepair(SCIP *scip)
#define STPTEST_ASSERT(cond)
Definition: stptest.h:63
void graph_edge_del(SCIP *, GRAPH *, int, SCIP_Bool)
Definition: graph_edge.c:368
int source
Definition: graphdefs.h:195
void graph_edge_addBi(SCIP *, GRAPH *, int, int, double)
static SCIP_RETCODE testTerminalPathsTo3NextFound(SCIP *scip)
#define FALSE
Definition: def.h:87
#define TRUE
Definition: def.h:86
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:54
includes various files containing graph methods used for Steiner tree problems
#define FARAWAY
Definition: graphdefs.h:89
static SCIP_RETCODE testTerminalPathsRepair2(SCIP *scip)
#define STP_TERM_NONE
Definition: graphdefs.h:62
void graph_tpathsSetAll4(GRAPH *, const SCIP_Real *, const SCIP_Real *, const SDPROFIT *, TPATHS *)
Definition: graph_tpath.c:2260
static SCIP_RETCODE testTerminalPathsRepair3(SCIP *scip)
SCIP_Real reduce_sdGetSd(const GRAPH *, int, int, SCIP_Real, SCIP_Real, SD *)
Definition: reduce_sd.c:2433
SCIP_RETCODE graph_tpathsRepair(SCIP *, int, SCIP_Bool, const GRAPH *, TPATHS *)
Definition: graph_tpath.c:1744
void graph_knot_add(GRAPH *, int)
Definition: graph_node.c:64
#define NULL
Definition: lpi_spx1.cpp:155
SCIP_RETCODE stptest_tpaths(SCIP *scip)
#define EQ(a, b)
Definition: portab.h:79
void reduce_sdprofitFree(SCIP *, SDPROFIT **)
#define SCIP_CALL(x)
Definition: def.h:384
void reduce_sdFree(SCIP *, SD **)
SCIP_RETCODE graph_tpathsRepairSetUp(const GRAPH *, TPATHS *)
Definition: graph_tpath.c:1719
static SCIP_RETCODE testSdRepair(SCIP *scip)
void stptest_graphTearDown(SCIP *, GRAPH *)
SCIP_RETCODE reduce_sdRepair(SCIP *, int, SCIP_Bool, GRAPH *, SD *)
#define STP_TERM
Definition: graphdefs.h:61
void graph_tpathsGet4CloseTerms(const GRAPH *, const TPATHS *, int, SCIP_Real, int *RESTRICT, int *RESTRICT, SCIP_Real *RESTRICT, int *RESTRICT)
static SCIP_RETCODE testBiasedTerminalPathsTo4NextFound(SCIP *scip)
Portable definitions.
includes various testing methods for Steiner tree problems
SCIP_Real * cost
Definition: graphdefs.h:221
SCIP_RETCODE stptest_graphSetUp(SCIP *, GRAPH *)
SCIP_RETCODE graph_tpathsInit(SCIP *, GRAPH *, TPATHS **)
Definition: graph_tpath.c:1682
#define SCIP_Real
Definition: def.h:177
SCIP_RETCODE reduce_sdprofitInit(SCIP *, const GRAPH *, SDPROFIT **)
#define nnodes
Definition: gastrans.c:65
SCIP_RETCODE graph_init(SCIP *, GRAPH **, int, int, int)
Definition: graph_base.c:607
SCIP_RETCODE reduce_sdRepairSetUp(SCIP *, const GRAPH *, SD *)
void graph_tpathsFree(SCIP *, TPATHS **)
Definition: graph_tpath.c:2300
SCIP callable library.