Scippy

SCIP

Solving Constraint Integer Programs

stptest_extutils.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_extutils.c
17  * @brief tests for Steiner tree extended reduction utilites
18  * @author Daniel Rehfeldt
19  *
20  * This file implements tests for Steiner problem reductions.
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 
28 
29 #include <stdio.h>
30 #include <assert.h>
31 #include "stptest.h"
32 #include "portab.h"
33 #include "graph.h"
34 #include "reduce.h"
35 #include "extreduce.h"
36 
37 
38 /** tests building and un-building */
39 static
41  SCIP* scip /**< SCIP data structure */
42 )
43 {
44  const int maxnlevels = 3;
45  const int maxnslots = 12;
46  const int maxntargets = 7;
47 #ifndef NDEBUG
48  const SCIP_Real* dists;
49 #endif
50 
51  MLDISTS* mldists;
52 
53  SCIP_CALL( extreduce_mldistsInit(scip, maxnlevels, maxnslots, maxntargets, 0, TRUE, &mldists) );
54 
55  /* add */
56 
57  extreduce_mldistsLevelAddTop(2, 7, mldists);
58 
59 #ifndef NDEBUG
61  assert(dists);
62 #endif
63 
64  for( int i = 0; i < 2; i++ )
65  {
66  extreduce_mldistsEmptySlotSetBase(i + 10, mldists);
68  }
69 
70  extreduce_mldistsLevelAddTop(12, 5, mldists);
71 
72  for( int i = 0; i < 12; i++ )
73  {
76  }
77 
78  extreduce_mldistsLevelAddTop(2, 7, mldists);
79 
80  for( int i = 0; i < 2; i++ )
81  {
84  }
85 
86  /* remove */
87 
91 
92  STPTEST_ASSERT_MSG(extreduce_mldistsIsEmpty(mldists), "MLDISTS not empty \n");
93 
94  extreduce_mldistsFree(scip, &mldists);
95 
96  return SCIP_OKAY;
97 }
98 
99 
100 
101 /** test close nodes are computed correctly */
102 static
104  SCIP* scip /**< SCIP data structure */
105 )
106 {
107  DISTDATA* distdata;
108  GRAPH* graph;
109  const int nnodes = 4;
110  const int nedges = 8;
111  const int nclosenodes = 3;
112 
113  assert(scip);
114 
115  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
116 
117  /* build tree */
118  graph_knot_add(graph, STP_TERM); /* node 0 */
119  graph_knot_add(graph, STP_TERM); /* node 1 */
120  graph_knot_add(graph, STP_TERM); /* node 2 */
121  graph_knot_add(graph, STP_TERM); /* node 3 */
122 
123  graph->source = 0;
124 
125  graph_edge_addBi(scip, graph, 0, 1, 1.0);
126  graph_edge_addBi(scip, graph, 0, 2, 1.0);
127  graph_edge_addBi(scip, graph, 0, 3, 2.1);
128  graph_edge_addBi(scip, graph, 2, 3, 1.0);
129 
130  SCIP_CALL( graph_pc_initPrizes(scip, graph, nnodes) );
131 
132  graph->prize[0] = 5.0;
133  graph->prize[1] = 0.6;
134  graph->prize[2] = 0.5;
135  graph->prize[3] = 0.5;
136 
137  SCIP_CALL( stptest_graphSetUpPcOrg(scip, graph, NULL, NULL) );
138 
139  SCIP_CALL( graph_init_dcsr(scip, graph) );
140  SCIP_CALL( extreduce_distDataInit(scip, graph, nclosenodes, FALSE, FALSE, &distdata) );
141 
142  /* 2. do the actual test */
143 
144  /* test close nodes */
145  {
146  const RANGE* node_range = distdata->closenodes_range;
147  const int* node_idx = distdata->closenodes_indices;
148  const SCIP_Real* node_dist = distdata->closenodes_distances;
149 
150  STPTEST_ASSERT(node_idx[node_range[0].start] == 1);
151  STPTEST_ASSERT(node_idx[node_range[0].start + 1] == 2);
152  STPTEST_ASSERT(node_idx[node_range[0].start + 2] == 3);
153 
154  STPTEST_ASSERT(EQ(node_dist[node_range[0].start], 1.0));
155  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 1], 1.0));
156  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 2], 1.5));
157  }
158 
159  extreduce_distDataFree(scip, graph, &distdata);
160  graph_free_dcsr(scip, graph);
161  stptest_extreduceTearDown(scip, graph, NULL);
162 
163  return SCIP_OKAY;
164 }
165 
166 
167 /** test close nodes are computed correctly */
168 static
170  SCIP* scip /**< SCIP data structure */
171 )
172 {
173  DISTDATA* distdata;
174  GRAPH* graph;
175  const int nnodes = 5;
176  const int nedges = 10;
177  const int nclosenodes = 4;
178 
179  assert(scip);
180 
181  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
182 
183  /* build tree */
184  graph_knot_add(graph, STP_TERM); /* node 0 */
185  graph_knot_add(graph, STP_TERM); /* node 1 */
186  graph_knot_add(graph, STP_TERM); /* node 2 */
187  graph_knot_add(graph, STP_TERM); /* node 3 */
188  graph_knot_add(graph, STP_TERM_NONE); /* node 4 */
189 
190  graph->source = 0;
191 
192  graph_edge_addBi(scip, graph, 0, 1, 1.5);
193  graph_edge_addBi(scip, graph, 1, 2, 1.0);
194  graph_edge_addBi(scip, graph, 1, 3, 2.0);
195  graph_edge_addBi(scip, graph, 1, 4, 1.0);
196  graph_edge_addBi(scip, graph, 3, 2, 0.5);
197 
198  SCIP_CALL( graph_pc_initPrizes(scip, graph, nnodes) );
199 
200  graph->prize[0] = 5.0;
201  graph->prize[1] = 2.0;
202  graph->prize[2] = 1.0;
203  graph->prize[3] = 1.0;
204  graph->prize[4] = 0.0;
205 
206  SCIP_CALL( stptest_graphSetUpPcOrg(scip, graph, NULL, NULL) );
207 
208  SCIP_CALL( graph_init_dcsr(scip, graph) );
209  SCIP_CALL( extreduce_distDataInit(scip, graph, nclosenodes, FALSE, FALSE, &distdata) );
210 
211  /* 2. do the actual test */
212 
213  /* test close nodes */
214  {
215  const RANGE* node_range = distdata->closenodes_range;
216  const SCIP_Real* node_dist = distdata->closenodes_distances;
217 
218  STPTEST_ASSERT(EQ(node_dist[node_range[0].start], 1.5));
219  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 1], 1.5));
220  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 2], 1.5));
221  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 3], 1.5));
222  }
223 
224  extreduce_distDataFree(scip, graph, &distdata);
225  graph_free_dcsr(scip, graph);
226  stptest_extreduceTearDown(scip, graph, NULL);
227 
228  return SCIP_OKAY;
229 }
230 
231 
232 /** test close nodes are computed correctly */
233 static
235  SCIP* scip /**< SCIP data structure */
236 )
237 {
238  DISTDATA* distdata;
239  GRAPH* graph;
240  const int nnodes = 5;
241  const int nedges = 12;
242  const int nclosenodes = 4;
243 
244  assert(scip);
245 
246  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
247 
248  /* build tree */
249  graph_knot_add(graph, STP_TERM_NONE); /* node 0 */
250  graph_knot_add(graph, STP_TERM); /* node 1 */
251  graph_knot_add(graph, STP_TERM); /* node 2 */
252  graph_knot_add(graph, STP_TERM); /* node 3 */
253  graph_knot_add(graph, STP_TERM); /* node 4 */
254 
255  graph->source = 0;
256 
257  graph_edge_addBi(scip, graph, 0, 1, 1.0);
258  graph_edge_addBi(scip, graph, 0, 2, 1.0);
259  graph_edge_addBi(scip, graph, 1, 3, 1.5);
260  graph_edge_addBi(scip, graph, 2, 3, 1.0);
261  graph_edge_addBi(scip, graph, 2, 4, 0.2);
262  graph_edge_addBi(scip, graph, 3, 4, 0.7);
263 
264  SCIP_CALL( graph_pc_initPrizes(scip, graph, nnodes) );
265 
266  graph->prize[0] = 0.0;
267  graph->prize[1] = 1.2;
268  graph->prize[2] = 0.4;
269  graph->prize[3] = 1.0;
270  graph->prize[4] = 2.0;
271 
272  SCIP_CALL( stptest_graphSetUpPcOrg(scip, graph, NULL, NULL) );
273 
274  SCIP_CALL( graph_init_dcsr(scip, graph) );
275  SCIP_CALL( extreduce_distDataInit(scip, graph, nclosenodes, FALSE, FALSE, &distdata) );
276 
277  /* 2. do the actual test */
278 
279  /* test close nodes */
280  {
281  const RANGE* node_range = distdata->closenodes_range;
282  const int* node_idx = distdata->closenodes_indices;
283  const SCIP_Real* node_dist = distdata->closenodes_distances;
284 
285  STPTEST_ASSERT(EQ(node_dist[node_range[0].start], 1.0));
286  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 1], 1.0));
287  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 2], 1.5));
288  STPTEST_ASSERT(EQ(node_dist[node_range[0].start + 3], 1.0));
289 
290  STPTEST_ASSERT(node_idx[node_range[0].start + 3] == 4);
291 
292  extreduce_edgeRemove(scip, 2, graph, distdata, NULL);
293  extreduce_edgeRemove(scip, 4, graph, distdata, NULL);
294  extreduce_edgeRemove(scip, 8, graph, distdata, NULL);
295 
296  STPTEST_ASSERT(extreduce_distCloseNodesAreValid(scip, graph, distdata));
297  }
298 
299  extreduce_distDataFree(scip, graph, &distdata);
300  graph_free_dcsr(scip, graph);
301  stptest_extreduceTearDown(scip, graph, NULL);
302 
303  return SCIP_OKAY;
304 }
305 
306 
307 /** test close nodes are computed correctly */
308 static
310  SCIP* scip /**< SCIP data structure */
311 )
312 {
313  GRAPH* graph;
314  const int nnodes = 55;
315  const int nedges = 28;
316  const int root = 0;
317  const int nclosenodes = 2; /**< max. number of close nodes to each node */
318 
319  DISTDATA* distdata;
320 
321  assert(scip);
322 
323  /* 1. build a test graph */
324 
325  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
326 
327  graph_knot_add(graph, 0);
328 
329  /* also add dummy nodes to avoid full stack */
330  for( int i = 1; i < nnodes; i++ )
332 
333  graph->source = root;
334 
335  graph_edge_addBi(scip, graph, 0, 1, 1.0);
336  graph_edge_addBi(scip, graph, 1, 2, 0.4);
337  graph_edge_addBi(scip, graph, 1, 3, 1.0);
338  graph_edge_addBi(scip, graph, 1, 4, 1.0);
339  graph_edge_addBi(scip, graph, 2, 5, 0.5);
340  graph_edge_addBi(scip, graph, 3, 6, 1.6);
341  graph_edge_addBi(scip, graph, 3, 7, 1.6);
342  graph_edge_addBi(scip, graph, 4, 8, 1.5);
343  graph_edge_addBi(scip, graph, 4, 9, 1.5);
344  graph_edge_addBi(scip, graph, 5, 10, 1.0);
345  graph_edge_addBi(scip, graph, 10, 11, 1.0);
346  graph_edge_addBi(scip, graph, 11, 12, 1.0);
347 
348  graph_mark(graph);
349  SCIP_CALL( graph_initHistory(scip, graph) );
350  SCIP_CALL( graph_path_init(scip, graph) );
351 
352  SCIP_CALL( graph_init_dcsr(scip, graph) );
353  SCIP_CALL( extreduce_distDataInit(scip, graph, nclosenodes, FALSE, FALSE, &distdata) );
354 
355  /* 2. do the actual test */
356 
357  /* test close nodes */
358  {
359  const RANGE* node_range = distdata->closenodes_range;
360  const int* node_idx = distdata->closenodes_indices;
361 #if 0
362  const SCIP_Real* node_dist = distdata->closenodes_distances;
363 
364  for( int node = 0; node < 6; node++ )
365  {
366  printf("node %d: \n", node);
367 
368  for( int i = node_range[node].start; i < node_range[node].end; i++ )
369  {
370  printf("...closenode=%d dist=%f\n", node_idx[i], node_dist[i]);
371  }
372  }
373 #endif
374 
375  STPTEST_ASSERT(node_idx[node_range[1].start] == 2);
376  STPTEST_ASSERT(node_idx[node_range[1].start + 1] == 5);
377  STPTEST_ASSERT(node_idx[node_range[5].start] == 1);
378  STPTEST_ASSERT(node_idx[node_range[5].start + 1] == 2);
379  }
380 
381  extreduce_distDataFree(scip, graph, &distdata);
382  graph_free_dcsr(scip, graph);
383 
384  graph_path_exit(scip, graph);
385  graph_free(scip, &graph, TRUE);
386  assert(graph == NULL);
387 
388  return SCIP_OKAY;
389 }
390 
391 
392 
393 /** test root paths are computed correctly */
394 static
396  SCIP* scip /**< SCIP data structure */
397 )
398 {
399  GRAPH* graph;
400  const int nnodes = 55;
401  const int nedges = 28;
402  const int root = 0;
403  const int nclosenodes = 2; /**< max. number of close nodes to each node */
404 
405  DISTDATA* distdata;
406 
407  assert(scip);
408 
409  /* 1. build a test graph */
410 
411  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
412 
413  graph_knot_add(graph, 0);
414 
415  /* also add dummy nodes to avoid full stack */
416  for( int i = 1; i < nnodes; i++ )
418 
419  graph->source = root;
420 
421  graph_edge_addBi(scip, graph, 0, 1, 1.0);
422  graph_edge_addBi(scip, graph, 1, 2, 0.4);
423  graph_edge_addBi(scip, graph, 1, 3, 1.0);
424  graph_edge_addBi(scip, graph, 1, 4, 1.0);
425  graph_edge_addBi(scip, graph, 2, 5, 0.5);
426  graph_edge_addBi(scip, graph, 3, 6, 1.6);
427  graph_edge_addBi(scip, graph, 3, 7, 1.6);
428  graph_edge_addBi(scip, graph, 4, 8, 1.5);
429  graph_edge_addBi(scip, graph, 4, 9, 1.5);
430  graph_edge_addBi(scip, graph, 5, 10, 1.0);
431  graph_edge_addBi(scip, graph, 10, 11, 1.0);
432  graph_edge_addBi(scip, graph, 11, 12, 1.0);
433 
434  graph_mark(graph);
435  SCIP_CALL( graph_initHistory(scip, graph) );
436  SCIP_CALL( graph_path_init(scip, graph) );
437 
438  SCIP_CALL( graph_init_dcsr(scip, graph) );
439  SCIP_CALL( extreduce_distDataInit(scip, graph, nclosenodes, FALSE, FALSE, &distdata) );
440 
441  /* 2. do the actual test */
442 
443  /* test edge root paths */
444  {
445  const int edge = 2;
446  PRSTATE** pathroot_blocks = distdata->pathroot_blocks;
447  int* pathroot_blocksizes = distdata->pathroot_blocksizes;
448 
449 
450 #ifdef SCIP_DEBUG
451  printf("edge ");
452  graph_edge_printInfo(graph, edge);
453 #endif
454 
455  STPTEST_ASSERT(graph->head[edge] == 2 && graph->tail[edge] == 1);
456  STPTEST_ASSERT(pathroot_blocksizes[edge / 2] == 6);
457 
458  for( int i = 0; i < pathroot_blocksizes[edge / 2]; i++ )
459  SCIPdebugMessage("...root=%d \n", pathroot_blocks[edge / 2][i].pathroot_id);
460  }
461 
462  extreduce_distDataFree(scip, graph, &distdata);
463  graph_free_dcsr(scip, graph);
464 
465  graph_path_exit(scip, graph);
466  graph_free(scip, &graph, TRUE);
467  assert(graph == NULL);
468 
469  return SCIP_OKAY;
470 }
471 
472 
473 /** test that distances are computed correctly */
474 static
476  SCIP* scip /**< SCIP data structure */
477 )
478 {
479  GRAPH* graph;
480  const int nnodes = 55;
481  const int nedges = 28;
482  const int root = 0;
483  const int nclosenodes = 2; /**< max. number of close nodes to each node */
484 
485  DISTDATA* distdata;
486 
487  assert(scip);
488 
489  /* 1. build a test graph */
490 
491  SCIP_CALL( graph_init(scip, &graph, nnodes, nedges, 1) );
492 
493  graph_knot_add(graph, STP_TERM);
494 
495  for( int i = 1; i < nnodes; i++ )
497 
498  graph->source = root;
499 
500  graph_edge_addBi(scip, graph, 0, 1, 1.0);
501  graph_edge_addBi(scip, graph, 1, 2, 0.4);
502  graph_edge_addBi(scip, graph, 1, 3, 1.0);
503  graph_edge_addBi(scip, graph, 1, 4, 1.0);
504  graph_edge_addBi(scip, graph, 2, 5, 0.5);
505  graph_edge_addBi(scip, graph, 3, 6, 1.6);
506  graph_edge_addBi(scip, graph, 3, 7, 1.6);
507  graph_edge_addBi(scip, graph, 4, 8, 1.5);
508  graph_edge_addBi(scip, graph, 4, 9, 1.5);
509  graph_edge_addBi(scip, graph, 5, 10, 1.0);
510  graph_edge_addBi(scip, graph, 10, 11, 1.0);
511  graph_edge_addBi(scip, graph, 11, 12, 1.0);
512 
513  graph_mark(graph);
514  SCIP_CALL( graph_initHistory(scip, graph) );
515  SCIP_CALL( graph_path_init(scip, graph) );
516 
517  SCIP_CALL( graph_init_dcsr(scip, graph) );
518  SCIP_CALL( extreduce_distDataInit(scip, graph, nclosenodes, FALSE, FALSE, &distdata) );
519 
520  /* test distances */
521  {
522  const int edge = 2;
523  const SCIP_Real dist1_2 = extreduce_distDataGetSd(scip, graph, 1, 2, distdata);
524  const SCIP_Real dist1_3 = extreduce_distDataGetSd(scip, graph, 1, 3, distdata);
525  const SCIP_Real dist1_5 = extreduce_distDataGetSd(scip, graph, 1, 5, distdata);
526  const SCIP_Real dist2_1 = extreduce_distDataGetSd(scip, graph, 2, 1, distdata);
527  const SCIP_Real dist2_5 = extreduce_distDataGetSd(scip, graph, 2, 5, distdata);
528  const SCIP_Real dist2_6 = extreduce_distDataGetSd(scip, graph, 2, 6, distdata);
529 
530  STPTEST_ASSERT(dist1_2 == 0.4);
531  STPTEST_ASSERT(dist1_3 == -1.0);
532  STPTEST_ASSERT(dist1_5 == 0.9);
533  STPTEST_ASSERT(dist2_1 == 0.4);
534  STPTEST_ASSERT(dist2_5 == 0.5);
535  STPTEST_ASSERT(dist2_6 == -1.0);
536 
537  assert(graph->head[edge] == 2 && graph->tail[edge] == 1);
538  graph_edge_delFull(scip, graph, edge, TRUE);
539  extreduce_distDataDeleteEdge(scip, graph, edge, distdata);
540 
541  {
542  const SCIP_Real dist1_2_b = extreduce_distDataGetSd(scip, graph, 1, 2, distdata);
543  const SCIP_Real dist2_6_b = extreduce_distDataGetSd(scip, graph, 2, 6, distdata);
544  const SCIP_Real dist2_5_b = extreduce_distDataGetSd(scip, graph, 2, 5, distdata);
545 
546  STPTEST_ASSERT(dist1_2_b == -1.0);
547  STPTEST_ASSERT(dist2_6_b == -1.0);
548  STPTEST_ASSERT(dist2_5_b == 0.5);
549  }
550  }
551 
552  extreduce_distDataFree(scip, graph, &distdata);
553  graph_free_dcsr(scip, graph);
554 
555  graph_path_exit(scip, graph);
556  graph_free(scip, &graph, TRUE);
557  assert(graph == NULL);
558 
559  return SCIP_OKAY;
560 }
561 
562 
563 /** tests for utilits */
565  SCIP* scip /**< SCIP data structure */
566 )
567 {
575 
576  printf("extended utilities test: all ok \n");
577 
578  return SCIP_OKAY;
579 }
#define STPTEST_ASSERT(cond)
Definition: stptest.h:63
int *RESTRICT head
Definition: graphdefs.h:224
static SCIP_RETCODE testDistCloseNodesPcAreValid1(SCIP *scip)
SCIP_RETCODE stptest_extmldists(SCIP *scip)
int source
Definition: graphdefs.h:195
void graph_edge_addBi(SCIP *, GRAPH *, int, int, double)
void extreduce_distDataFree(SCIP *, const GRAPH *, DISTDATA **)
void graph_mark(GRAPH *)
Definition: graph_base.c:1130
void graph_free_dcsr(SCIP *, GRAPH *)
Definition: graph_util.c:1807
void extreduce_mldistsEmptySlotSetBase(int, MLDISTS *)
void graph_free(SCIP *, GRAPH **, SCIP_Bool)
Definition: graph_base.c:796
SCIP_RETCODE graph_path_init(SCIP *, GRAPH *)
Definition: graph_path.c:480
#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
SCIP_RETCODE stptest_graphSetUpPcOrg(SCIP *, GRAPH *, int *, int *)
void stptest_extreduceTearDown(SCIP *, GRAPH *, REDCOST **)
#define SCIPdebugMessage
Definition: pub_message.h:87
RANGE * closenodes_range
Definition: extreducedefs.h:84
int * closenodes_indices
Definition: extreducedefs.h:85
SCIP_RETCODE graph_init_dcsr(SCIP *, GRAPH *)
Definition: graph_util.c:1721
int * pathroot_blocksizes
Definition: extreducedefs.h:88
#define STP_TERM_NONE
Definition: graphdefs.h:62
SCIP_RETCODE graph_pc_initPrizes(SCIP *, GRAPH *, int)
Definition: graph_pcbase.c:741
This file implements extended reduction techniques for several Steiner problems.
PRSTATE ** pathroot_blocks
Definition: extreducedefs.h:90
#define STPTEST_ASSERT_MSG(cond, msg)
Definition: stptest.h:51
SCIP_Bool extreduce_mldistsIsEmpty(const MLDISTS *)
static SCIP_RETCODE testMldistsBuilding(SCIP *scip)
void extreduce_mldistsEmptySlotSetFilled(MLDISTS *)
static SCIP_RETCODE testDistRootPathsAreValid(SCIP *scip)
SCIP_Real * prize
Definition: graphdefs.h:210
SCIP_RETCODE graph_initHistory(SCIP *, GRAPH *)
Definition: graph_base.c:695
void graph_knot_add(GRAPH *, int)
Definition: graph_node.c:64
SCIP_Real * closenodes_distances
Definition: extreducedefs.h:87
#define NULL
Definition: lpi_spx1.cpp:155
static SCIP_RETCODE testDistCloseNodesAreValid(SCIP *scip)
void extreduce_edgeRemove(SCIP *, int, GRAPH *, DISTDATA *, EXTPERMA *)
void extreduce_mldistsLevelRemoveTop(MLDISTS *)
#define EQ(a, b)
Definition: portab.h:79
#define SCIP_CALL(x)
Definition: def.h:384
void graph_edge_delFull(SCIP *, GRAPH *, int, SCIP_Bool)
Definition: graph_edge.c:418
void extreduce_distDataDeleteEdge(SCIP *, const GRAPH *, int, DISTDATA *)
SCIP_Real extreduce_distDataGetSd(SCIP *, const GRAPH *, int, int, DISTDATA *)
int *RESTRICT tail
Definition: graphdefs.h:223
SCIP_Bool extreduce_distCloseNodesAreValid(SCIP *, const GRAPH *, const DISTDATA *)
#define STP_TERM
Definition: graphdefs.h:61
SCIP_Real * extreduce_mldistsEmptySlotTargetDists(const MLDISTS *)
SCIP_RETCODE extreduce_mldistsInit(SCIP *, int, int, int, int, SCIP_Bool, MLDISTS **)
void graph_edge_printInfo(const GRAPH *, int)
Definition: graph_stats.c:133
void graph_path_exit(SCIP *, GRAPH *)
Definition: graph_path.c:515
Portable definitions.
static SCIP_RETCODE testDistDistancesAreValid(SCIP *scip)
static SCIP_RETCODE testDistCloseNodesPcAreValidAfterDeletion(SCIP *scip)
void extreduce_mldistsLevelAddTop(int, int, MLDISTS *)
SCIP_RETCODE extreduce_distDataInit(SCIP *, GRAPH *, int, SCIP_Bool, SCIP_Bool, DISTDATA **)
includes various testing methods for Steiner tree problems
void extreduce_mldistsFree(SCIP *, MLDISTS **)
#define SCIP_Real
Definition: def.h:177
#define nnodes
Definition: gastrans.c:65
SCIP_RETCODE graph_init(SCIP *, GRAPH **, int, int, int)
Definition: graph_base.c:607
static SCIP_RETCODE testDistCloseNodesPcAreValid2(SCIP *scip)
includes various reduction methods for Steiner tree problems