Scippy

UG

Ubiquity Generator framework

scipParaSolver.cpp
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the program and software framework */
4/* UG --- Ubquity Generator Framework */
5/* */
6/* Copyright Written by Yuji Shinano <shinano@zib.de>, */
7/* Copyright (C) 2021-2024 by Zuse Institute Berlin, */
8/* licensed under LGPL version 3 or later. */
9/* Commercial licenses are available through <licenses@zib.de> */
10/* */
11/* This code is free software; you can redistribute it and/or */
12/* modify it under the terms of the GNU Lesser General Public License */
13/* as published by the Free Software Foundation; either version 3 */
14/* of the License, or (at your option) any later version. */
15/* */
16/* This program is distributed in the hope that it will be useful, */
17/* but WITHOUT ANY WARRANTY; without even the implied warranty of */
18/* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
19/* GNU Lesser General Public License for more details. */
20/* */
21/* You should have received a copy of the GNU Lesser General Public License */
22/* along with this program. If not, see <http://www.gnu.org/licenses/>. */
23/* */
24/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
25
26/**@file scipParaSolver.cpp
27 * @brief ParaSolver extension for SCIP: Parallelized solver implementation for SCIP.
28 * @author Yuji Shinano
29 *
30 *
31 *
32 */
33
34/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
35
36
37#include <cfloat>
38#include <cstring>
39#include <cstdlib>
40#ifndef _MSC_VER
41#include <unistd.h>
42#endif
43#include <typeinfo>
44#include <string>
45#include <sstream>
46#include "ug/paraInitialStat.h"
47#include "ug_bb/bbParaComm.h"
48#include "ug_bb/bbParaNode.h"
50#include "ug_bb/bbParaSolver.h"
53#include "objscip/objscip.h"
54#include "scipParaTagDef.h"
55#include "scipParaParamSet.h"
59#include "scipParaObjProp.h"
61#include "scipParaInitialStat.h"
63#include "scipParaObjNodesel.h"
65#include "scip/scip.h"
66#ifdef UG_DEBUG_SOLUTION
67#ifndef WITH_DEBUG_SOLUTION
68#define WITH_DEBUG_SOLUTION
69#endif
70#include "scip/debug.h"
71#include "scip/struct_scip.h"
72#include "scip/struct_set.h"
73#endif
74// #include "scip/scipdefplugins.h"
75
76using namespace ParaSCIP;
77
78#if ( defined(_COMM_PTH) || defined(_COMM_CPP11) )
79extern long long virtualMemUsedAtLc;
80extern double memoryLimitOfSolverSCIP;
81#endif
82
83extern void
85extern void
87
88/*
89 * Callback methods of conflict handler
90 */
91#define CONFLICTHDLR_NAME "conflictCollector"
92#define CONFLICTHDLR_DESC "conflict handler to collect conflicts"
93#define CONFLICTHDLR_PRIORITY +100000000
94static
95SCIP_DECL_CONFLICTEXEC(conflictExecCollector)
96{ /*lint --e{715}*/
97 SCIP_VAR** vars;
98 SCIP_Real* vals;
99 SCIP_Real lhs;
100 int i;
101
102 assert(conflicthdlr != NULL);
103 assert(strcmp(SCIPconflicthdlrGetName(conflicthdlr), CONFLICTHDLR_NAME) == 0);
104 assert(bdchginfos != NULL || nbdchginfos == 0);
105 assert(result != NULL);
106
107 /* don't process already resolved conflicts */
108 if( resolved )
109 {
110 *result = SCIP_DIDNOTRUN;
111 return SCIP_OKAY;
112 }
113
114 *result = SCIP_DIDNOTFIND;
115
116 /* create array of variables and coefficients: sum_{i \in P} x_i - sum_{i \in N} x_i >= 1 - |N| */
117 SCIP_CALL( SCIPallocBufferArray(scip, &vars, nbdchginfos) );
118 SCIP_CALL( SCIPallocBufferArray(scip, &vals, nbdchginfos) );
119 lhs = 1.0;
120 for( i = 0; i < nbdchginfos; ++i )
121 {
122 assert(bdchginfos != NULL);
123
124 vars[i] = SCIPbdchginfoGetVar(bdchginfos[i]);
125
126 /* we can only treat binary variables */
127 /**@todo extend linear conflict constraints to some non-binary cases */
128 if( !SCIPvarIsBinary(vars[i]) )
129 break;
130
131 /* check whether the variable is fixed to zero (P) or one (N) in the conflict set */
132 if( SCIPbdchginfoGetNewbound(bdchginfos[i]) < 0.5 )
133 vals[i] = 1.0;
134 else
135 {
136 vals[i] = -1.0;
137 lhs -= 1.0;
138 }
139 }
140
141 if( i == nbdchginfos )
142 {
143 ScipParaSolver *scipParaSolver = reinterpret_cast<ScipParaSolver *>(SCIPconflicthdlrGetData(conflicthdlr));
144 std::list<LocalNodeInfoPtr> *conflictConsList = scipParaSolver->getConflictConsList();
145 LocalNodeInfo *localNodeInfo = new LocalNodeInfo;
146 localNodeInfo->linearRhs = SCIPinfinity(scip);
147 localNodeInfo->nLinearCoefs = nbdchginfos;
148 localNodeInfo->idxLinearCoefsVars = new int[nbdchginfos];
149 localNodeInfo->linearCoefs = new double[nbdchginfos];
150 for( i = 0; i < nbdchginfos; ++i )
151 {
152 SCIP_VAR *transformVar = vars[i];
153 SCIP_Real scalar = vals[i];
154 SCIP_Real constant = 0.0;
155 if( SCIPvarGetOrigvarSum(&transformVar, &scalar, &constant ) == SCIP_INVALIDDATA )
156 break;
157 // assert(transformVar != NULL);
158 if( transformVar )
159 {
160 lhs -= constant;
161 if( scipParaSolver->isOriginalIndeciesMap() )
162 {
163 localNodeInfo->idxLinearCoefsVars[i] = scipParaSolver->getOriginalIndex(SCIPvarGetIndex(transformVar));
164 }
165 else
166 {
167 localNodeInfo->idxLinearCoefsVars[i] = SCIPvarGetIndex(transformVar);
168 }
169 localNodeInfo->linearCoefs[i] = scalar;
170 }
171 else
172 {
173 break;
174 }
175 }
176 if( i == nbdchginfos )
177 {
178 localNodeInfo->linearLhs = lhs;
179 conflictConsList->push_back(localNodeInfo);
180 }
181 else
182 {
183 delete [] localNodeInfo->idxLinearCoefsVars;
184 delete [] localNodeInfo->linearCoefs;
185 delete localNodeInfo;
186 }
187 }
188
189 /* free temporary memory */
190 SCIPfreeBufferArray(scip, &vals);
191 SCIPfreeBufferArray(scip, &vars);
192
193 return SCIP_OKAY;
194}
195
196void
197ScipParaSolver::setWinnerRacingParams(
198 UG::ParaRacingRampUpParamSet *inRacingParams /**< winner solver pramset */
199 )
200{
201 if( !userPlugins )
202 {
203 SCIP_CALL_ABORT( SCIPresetParams(scip) );
204 }
205
207 {
208 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
209 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
210 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
211 if( inRacingParams )
212 {
213#if (SCIP_VERSION < 321 || ( SCIP_VERSION == 321 && SCIP_SUBVERSION < 2) )
214 ScipParaRacingRampUpParamSet *scipRacingParams = dynamic_cast< ScipParaRacingRampUpParamSet * >(inRacingParams);
215 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "misc/permutationseed", scipRacingParams->getPermuteProbSeed()) );
216#endif
217 }
218 }
219 else
220 {
221 ScipParaRacingRampUpParamSet *scipRacingParams = dynamic_cast< ScipParaRacingRampUpParamSet * >(inRacingParams);
222 setRacingParams(scipRacingParams, true);
223 }
224
225#if SCIP_VERSION >= 320
227#endif
228
229}
230
231void
233 UG::ParaRacingRampUpParamSet *inRacingParams,
234 bool winnerParam
235 )
236{
237 ScipParaRacingRampUpParamSet *scipRacingParams = dynamic_cast< ScipParaRacingRampUpParamSet * >(inRacingParams);
238
239 if( !winnerParam && !userPlugins )
240 {
241 SCIP_CALL_ABORT( SCIPresetParams(scip) );
242 }
243
244 if ( std::string(paraParams->getStringParamValue(UG::RacingParamsDirPath)) != std::string("") )
245 {
246 assert( scipRacingParams->getScipDiffParamSet() );
247 if( !winnerParam )
248 {
249 scipRacingParams->getScipDiffParamSet()->setParametersInScip(scip);
250 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "randomization/randomseedshift", scipRacingParams->getScipRacingParamSeed()) );
251 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "randomization/permutationseed", scipRacingParams->getPermuteProbSeed()) );
252 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "randomization/lpseed", scipRacingParams->getPermuteProbSeed()) );
253 // int tempInt = 0;
254 // SCIP_CALL_ABORT( SCIPgetIntParam(scip, "lp/solvefreq", &tempInt) );
255 // std::cout << "R." << paraComm->getRank() << " lp/solvefreq = " << tempInt << std::endl;
256 }
257 }
258 else
259 {
261 {
262 if( !winnerParam && scipRacingParams->getScipDiffParamSet() )
263 {
264 scipRacingParams->getScipDiffParamSet()->setParametersInScip(scip);
265 }
266 }
267 else
268 {
269 int nHeuristics;
271 {
272 nHeuristics = scipRacingParams->getScipRacingParamSeed() % 2;
273 }
274 else
275 {
276 nHeuristics = scipRacingParams->getScipRacingParamSeed() % 4;
277 }
278 int nPresolving = (scipRacingParams->getScipRacingParamSeed()/4) % 4;
279 int nSeparating = (scipRacingParams->getScipRacingParamSeed()/(4*4)) % 4;
280
281 switch( nHeuristics )
282 {
283 case 0:
284 {
285 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
286 break;
287 }
288 case 1:
289 {
290 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
291 break;
292 }
293 case 2:
294 {
295 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_FAST, TRUE) );
296 break;
297 }
298 case 3:
299 {
300 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_OFF, TRUE) );
301 break;
302 }
303 default:
304 THROW_LOGICAL_ERROR1("invalid nHeuristics");
305 }
306 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "heuristics/ObjLim/freq", 1) );
307
308 switch( nPresolving )
309 {
310 case 0:
311 {
312 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
313 break;
314 }
315 case 1:
316 {
317#ifdef _COMM_PTH
319 {
320 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
321 }
322 else
323 {
324 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
325 }
326#else
327 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
328#endif
329 break;
330 }
331 case 2:
332 {
333 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_FAST, TRUE) );
334 break;
335 }
336 case 3:
337 {
338 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_OFF, TRUE) );
339 break;
340 }
341 default:
342 THROW_LOGICAL_ERROR1("invalid nPresolving");
343 }
344
345 switch( nSeparating )
346 {
347 case 0:
348 {
349 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
350 break;
351 }
352 case 1:
353 {
354#ifdef _COMM_PTH
356 {
357 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
358 }
359 else
360 {
362 {
363 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
364 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
365 }
366 else
367 {
368 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
369 }
370 }
371#else
373 {
374 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
375 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
376 }
377 else
378 {
379 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_AGGRESSIVE, TRUE) );
380 }
381#endif
382 break;
383 }
384 case 2:
385 {
386 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_FAST, TRUE) );
387 break;
388 }
389 case 3:
390 {
391 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_OFF, TRUE) );
392 break;
393 }
394 default:
395 THROW_LOGICAL_ERROR1("invalid nSeparating");
396 }
397 }
398
399 assert(SCIPgetStage(scip) <= SCIP_STAGE_TRANSFORMED);
400 // make sure that the permutation works on transformed problem
401#if (SCIP_VERSION < 321 || ( SCIP_VERSION == 321 && SCIP_SUBVERSION < 2) )
402 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "misc/permutationseed", scipRacingParams->getPermuteProbSeed()) );
403#endif
404
408 {
410 }
411
412 if( !winnerParam && scipRacingParams->getPermuteProbSeed() >= 64 ) // after all parameters tested, random branchig variable selection
413 {
414 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "branching/random/maxdepth", 2) );
415 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "branching/random/priority", 100000) );
416 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "branching/random/seed", scipRacingParams->getGenerateBranchOrderSeed()) );
417 }
418 }
419
421 {
422 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "randomization/randomseedshift", scipRacingParams->getScipRacingParamSeed()) );
423 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "randomization/permutationseed", scipRacingParams->getPermuteProbSeed()) );
424 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "randomization/lpseed", scipRacingParams->getPermuteProbSeed()) );
425 if( scipRacingParams->getScipDiffParamSet() )
426 {
427 scipRacingParams->getScipDiffParamSet()->setParametersInScip(scip);
428 }
429 }
430
431#if SCIP_VERSION >= 320
433#endif
434
435 // writeSubproblem();
436
437}
438
439void
441 )
442{
443 assert(currentTask);
444
445 UG::BbParaNode *bbCurrentNode = dynamic_cast<UG::BbParaNode *>(currentTask);
446
447#ifdef UG_DEBUG_SOLUTION
448 if( scip->set->debugsoldata == NULL )
449 {
450 SCIP_CALL_ABORT( SCIPdebugSolDataCreate(&((scip->set)->debugsoldata)));
451 // SCIPdebugSetMainscipset(scip->set);
452 }
453#endif
454
455 /** set instance specific parameters */
456 // if( currentNode->isRootNode() && !(paraParams->getBoolParamValue(UseRootNodeCuts)) )
457 // Probably, in order to avoid root node settings twice. Once root nodes is solved in LC
458 // when UseRootNodeCuts is specified. However, for racing ramp-up, this is too bad.
459 // So, I changed the specification. The root node parameter settings is applied twice now
460
461 // Do not reset here, because racing parameters might be set already.
462 // SCIP_CALL_ABORT( SCIPresetParams(scip) );
463
464 /** set original node selection strategy */
465 // setOriginalNodeSelectionStrategy();
467 nodesel->reset();
468 if( bbCurrentNode->isRootTask() )
469 {
470 if ( std::string(paraParams->getStringParamValue(UG::RacingParamsDirPath)) == std::string("") )
471 {
473 }
476 isRacingStage() &&
477 paraComm->getRank() == 1 ) // rank 1 should be all default
478 )
479 {
480 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_OFF, TRUE) );
481 }
482 }
483 else
484 {
486 {
487 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_DEFAULT, TRUE) );
488 }
490 {
491 if( scipDiffParamSet ) // this may not be necessary , check setWinnerRacingParams(0)
492 {
494 }
495 }
496
497 /*
498 if( paraParams->getBoolParamValue(UG::ControlCollectingModeOnSolverSide) )
499 {
500 int maxrestarts;
501 SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/maxrestarts", &maxrestarts) );
502 if( maxrestarts < 0 )
503 {
504 std::cerr << "presolving/maxrestarts >= 0 when you specify ControlCollectingModeOnSolverSide = TRUE."
505 << std::endl;
506 exit(1);
507 }
508 }
509 else
510 {
511 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/maxrestarts", 0 ) );
512 }
513 */
514 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/maxrestarts", getOriginalMaxRestart()) );
515 }
516
517 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "misc/usesymmetry", 0 ) ); // Symmetry handling technique is explicitly turn off in Solever for this version (SCIP 5.0)
518
519#if SCIP_VERSION >= 320
521#endif
522
523 double dualBoundValue = bbCurrentNode->getDualBoundValue();
524
525 ScipParaDiffSubproblem *scipParaDiffSubproblem = dynamic_cast< ScipParaDiffSubproblem* >(currentTask->getDiffSubproblem());
526
527 SCIP_VAR **orgVars = SCIPgetOrigVars(scip); // variables are indexed by index
528 int nOrg = SCIPgetNOrigVars(scip); // the number of original variables
529 if( scipParaDiffSubproblem )
530 {
532 {
533 assert( mapToSolverLocalIndecies );
534 for(int v = 0; v < scipParaDiffSubproblem->getNBoundChanges(); v++)
535 {
536 assert(mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]] >= 0);
537 if( mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]] < nOrg )
538 {
539 if( scipParaDiffSubproblem->getBoundType(v) == SCIP_BOUNDTYPE_LOWER )
540 {
541 SCIP_CALL_ABORT(
542 SCIPchgVarLbGlobal(
543 scip,
544 orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]],
545 scipParaDiffSubproblem->getBranchBound(v) )
546 );
547 if( scipParaDiffSubproblem->getIndex(v) < nOrgVars )
548 {
549 assert(SCIPisEQ(scip,SCIPvarGetLbGlobal(orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]]),scipParaDiffSubproblem->getBranchBound(v)));
550 assert(SCIPisLE(scip,SCIPvarGetLbGlobal(orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]]),SCIPvarGetUbGlobal(orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]])));
551 }
552 }
553 else if (scipParaDiffSubproblem->getBoundType(v) == SCIP_BOUNDTYPE_UPPER)
554 {
555 SCIP_CALL_ABORT(SCIPchgVarUbGlobal(
556 scip,
557 orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]],
558 scipParaDiffSubproblem->getBranchBound(v) )
559 );
560 if( scipParaDiffSubproblem->getIndex(v) < nOrgVars )
561 {
562 assert(SCIPisEQ(scip,SCIPvarGetUbGlobal(orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]]),scipParaDiffSubproblem->getBranchBound(v)));
563 assert(SCIPisLE(scip,SCIPvarGetLbGlobal(orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]]),SCIPvarGetUbGlobal(orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]]])));
564 }
565 }
566 else
567 {
568 THROW_LOGICAL_ERROR2("Invalid bound type: type = ", static_cast<int>(scipParaDiffSubproblem->getBoundType(v))) ;
569 }
570 }
571 else
572 {
573 std::cout << "fixing branching variable index = " << mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIndex(v)]] << " is omitted!" << std::endl;
574 }
575 }
576 }
577 else
578 {
579 for(int v = 0; v < scipParaDiffSubproblem->getNBoundChanges(); v++)
580 {
581 if( scipParaDiffSubproblem->getIndex(v) < nOrg )
582 {
583 if( scipParaDiffSubproblem->getBoundType(v) == SCIP_BOUNDTYPE_LOWER )
584 {
585 SCIP_CALL_ABORT(
586 SCIPchgVarLbGlobal(
587 scip,
588 orgVars[scipParaDiffSubproblem->getIndex(v)],
589 scipParaDiffSubproblem->getBranchBound(v) )
590 );
591 if( scipParaDiffSubproblem->getIndex(v) < nOrgVars )
592 {
593 assert(SCIPisEQ(scip,SCIPvarGetLbGlobal(orgVars[scipParaDiffSubproblem->getIndex(v)]),scipParaDiffSubproblem->getBranchBound(v)));
594 assert(SCIPisLE(scip,SCIPvarGetLbGlobal(orgVars[scipParaDiffSubproblem->getIndex(v)]),SCIPvarGetUbGlobal(orgVars[scipParaDiffSubproblem->getIndex(v)])));
595 }
596 }
597 else if (scipParaDiffSubproblem->getBoundType(v) == SCIP_BOUNDTYPE_UPPER)
598 {
599 SCIP_CALL_ABORT(SCIPchgVarUbGlobal(
600 scip,
601 orgVars[scipParaDiffSubproblem->getIndex(v)],
602 scipParaDiffSubproblem->getBranchBound(v) )
603 );
604 if( scipParaDiffSubproblem->getIndex(v) < nOrgVars )
605 {
606 assert(SCIPisEQ(scip,SCIPvarGetUbGlobal(orgVars[scipParaDiffSubproblem->getIndex(v)]),scipParaDiffSubproblem->getBranchBound(v)));
607 assert(SCIPisLE(scip,SCIPvarGetLbGlobal(orgVars[scipParaDiffSubproblem->getIndex(v)]),SCIPvarGetUbGlobal(orgVars[scipParaDiffSubproblem->getIndex(v)])));
608 }
609 }
610 else
611 {
612 THROW_LOGICAL_ERROR2("Invalid bound type: type = ", static_cast<int>(scipParaDiffSubproblem->getBoundType(v))) ;
613 }
614 }
615 else
616 {
617 std::cout << "fixing branching variable index = " << scipParaDiffSubproblem->getIndex(v) << " is omitted!" << std::endl;
618 }
619 }
620 }
621
622 if( scipParaDiffSubproblem->getNBranchConsLinearConss() > 0 ||
623 scipParaDiffSubproblem->getNBranchConsSetppcConss() > 0 ||
624 scipParaDiffSubproblem->getNLinearConss() > 0 ||
625 scipParaDiffSubproblem->getNBendersLinearConss() > 0 ||
626 scipParaDiffSubproblem->getNBoundDisjunctions() )
627 {
628 assert(addedConss == 0);
629 addedConss = new SCIP_CONS*[scipParaDiffSubproblem->getNBranchConsLinearConss()
630 + scipParaDiffSubproblem->getNBranchConsSetppcConss()
631 + scipParaDiffSubproblem->getNLinearConss()
632 + scipParaDiffSubproblem->getNBendersLinearConss()
633 + scipParaDiffSubproblem->getNBoundDisjunctions()];
634 }
635
636 SCIP_CONS* cons;
637 char consname[SCIP_MAXSTRLEN];
638
639 int c = 0;
640 for(; c < scipParaDiffSubproblem->getNBranchConsLinearConss() ; c++ )
641 {
642 SCIP_VAR** vars;
643 SCIP_Real* vals;
644 int nVars = scipParaDiffSubproblem->getBranchConsNLinearCoefs(c);
645
646 /* create array of variables and coefficients */
647 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vars, nVars) );
648 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vals, nVars) );
649
651 {
652 for( int v = 0; v < nVars; ++v )
653 {
654 vars[v] = orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getBranchConsLinearIdxCoefsVars(c,v)]]];
655 vals[v] = scipParaDiffSubproblem->getBranchConsLinearCoefs(c,v);
656 }
657 }
658 else
659 {
660 for( int v = 0; v < nVars; ++v )
661 {
662 vars[v] = orgVars[scipParaDiffSubproblem->getBranchConsLinearIdxCoefsVars(c,v)];
663 vals[v] = scipParaDiffSubproblem->getBranchConsLinearCoefs(c,v);
664 }
665 }
666
667 /* create a constraint */
668 (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "%s", scipParaDiffSubproblem->getBranchConsLinearConsNames(c));
669 SCIP_CALL_ABORT( SCIPcreateConsLinear(scip, &cons, consname, nVars, vars, vals,
670 scipParaDiffSubproblem->getBranchConsLinearLhs(c), scipParaDiffSubproblem->getBranchConsLinearRhs(c),
671 TRUE, TRUE, FALSE, FALSE, TRUE, FALSE, FALSE, FALSE, TRUE, FALSE) );
672 /** only a constraint whose "enforce is TRUE can be written in transformed problem */
673
674 /* add constraint to SCIP */
675 SCIP_CALL_ABORT( SCIPaddCons(scip, cons) );
676 assert(cons);
677 addedConss[c] = cons;
678 SCIP_CALL_ABORT( SCIPreleaseCons(scip, &cons) );
679 /* free temporary memory */
680 SCIPfreeBufferArray(scip, &vals);
681 SCIPfreeBufferArray(scip, &vars);
682 }
683
684 int i = 0;
685 for(; c < (scipParaDiffSubproblem->getNBranchConsLinearConss()
686 + scipParaDiffSubproblem->getNBranchConsSetppcConss()) ; c++ )
687 {
688 SCIP_VAR** vars;
689
690 int nVars = scipParaDiffSubproblem->getBranchConsSetppcNVars(i);
691 /* create array of variables, types and bounds */
692 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vars, nVars) );
693
695 {
696 for( int v = 0; v < nVars; ++v )
697 {
698 vars[v] = orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getBranchConsSetppcVars(i,v)]]];
699 }
700 }
701 else
702 {
703 for( int v = 0; v < nVars; ++v )
704 {
705 vars[v] = orgVars[scipParaDiffSubproblem->getBranchConsSetppcVars(i,v)];
706 }
707 }
708
709 /* create a constraint */
710 assert( scipParaDiffSubproblem->getBranchConsSetppcType(i) == SCIP_SETPPCTYPE_PARTITIONING ); // currently, only this should be used
711 (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "%s", scipParaDiffSubproblem->getBranchConsSetppcConsNames(i));
712 if( scipParaDiffSubproblem->getBranchConsSetppcType(i) == SCIP_SETPPCTYPE_PARTITIONING )
713 {
714 SCIP_CALL_ABORT( SCIPcreateConsSetpart(scip, &cons, consname, nVars, vars,
715 TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, FALSE, FALSE, FALSE ) );
716 } else if ( scipParaDiffSubproblem->getBranchConsSetppcType(i) == SCIP_SETPPCTYPE_PACKING )
717 {
718 SCIP_CALL_ABORT( SCIPcreateConsSetpack(scip, &cons, consname, nVars, vars,
719 TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, FALSE, FALSE, FALSE ) );
720 } else if ( scipParaDiffSubproblem->getBranchConsSetppcType(i) == SCIP_SETPPCTYPE_COVERING )
721 {
722 SCIP_CALL_ABORT( SCIPcreateConsSetcover(scip, &cons, consname, nVars, vars,
723 TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, FALSE, FALSE, FALSE ) );
724 } else {
725 THROW_LOGICAL_ERROR2("Unknown setppc constraint is received: type = ", scipParaDiffSubproblem->getBranchConsSetppcType(i));
726 }
727
728 /* add constraint to SCIP */
729 SCIP_CALL_ABORT( SCIPaddCons(scip, cons) );
730 assert(cons);
731 addedConss[c] = cons;
732 SCIP_CALL_ABORT( SCIPreleaseCons(scip, &cons) );
733 /* free temporary memory */
734 SCIPfreeBufferArray(scip, &vars);
735 i++;
736 }
737
738 i = 0;
739 for(; c < (scipParaDiffSubproblem->getNBranchConsLinearConss()
740 + scipParaDiffSubproblem->getNBranchConsSetppcConss()
741 + scipParaDiffSubproblem->getNLinearConss()) ; c++ )
742 {
743 SCIP_VAR** vars;
744 SCIP_Real* vals;
745 int nVars = scipParaDiffSubproblem->getNLinearCoefs(i);
746
747 /* create array of variables and coefficients */
748 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vars, nVars) );
749 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vals, nVars) );
750
752 {
753 for( int v = 0; v < nVars; ++v )
754 {
755 vars[v] = orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIdxLinearCoefsVars(i,v)]]];
756 vals[v] = scipParaDiffSubproblem->getLinearCoefs(i,v);
757 }
758 }
759 else
760 {
761 for( int v = 0; v < nVars; ++v )
762 {
763 vars[v] = orgVars[scipParaDiffSubproblem->getIdxLinearCoefsVars(i,v)];
764 vals[v] = scipParaDiffSubproblem->getLinearCoefs(i,v);
765 }
766 }
767
768 /* create a constraint */
769 (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "cli%d", i);
770 SCIP_CALL_ABORT( SCIPcreateConsLinear(scip, &cons, consname, nVars, vars, vals,
771 scipParaDiffSubproblem->getLinearLhs(i), scipParaDiffSubproblem->getLinearRhs(i),
772 TRUE, TRUE, FALSE, FALSE, TRUE, FALSE, FALSE, FALSE, TRUE, FALSE) );
773 /** only a constraint whose "enforce is TRUE can be written in transformed problem */
774
775 /* add constraint to SCIP */
776 SCIP_CALL_ABORT( SCIPaddCons(scip, cons) );
777 assert(cons);
778 addedConss[c] = cons;
779 SCIP_CALL_ABORT( SCIPreleaseCons(scip, &cons) );
780 /* free temporary memory */
781 SCIPfreeBufferArray(scip, &vals);
782 SCIPfreeBufferArray(scip, &vars);
783 i++;
784 }
785
786 i = 0;
787 for(; c < (scipParaDiffSubproblem->getNBranchConsLinearConss()
788 + scipParaDiffSubproblem->getNBranchConsSetppcConss()
789 + scipParaDiffSubproblem->getNLinearConss()
790 + scipParaDiffSubproblem->getNBendersLinearConss()) ; c++ )
791 {
792 SCIP_VAR** vars;
793 SCIP_Real* vals;
794 int nVars = scipParaDiffSubproblem->getNBendersLinearCoefs(i);
795
796 /* create array of variables and coefficients */
797 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vars, nVars) );
798 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vals, nVars) );
799
801 {
802 for( int v = 0; v < nVars; ++v )
803 {
804 vars[v] = orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIdxBendersLinearCoefsVars(i,v)]]];
805 vals[v] = scipParaDiffSubproblem->getBendersLinearCoefs(i,v);
806 }
807 }
808 else
809 {
810 for( int v = 0; v < nVars; ++v )
811 {
812 vars[v] = orgVars[scipParaDiffSubproblem->getIdxBendersLinearCoefsVars(i,v)];
813 vals[v] = scipParaDiffSubproblem->getBendersLinearCoefs(i,v);
814 }
815 }
816
817 /* create a constraint */
818 (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "cli%d", i);
819 SCIP_CALL_ABORT( SCIPcreateConsLinear(scip, &cons, consname, nVars, vars, vals,
820 scipParaDiffSubproblem->getBendersLinearLhs(i), scipParaDiffSubproblem->getBendersLinearRhs(i),
821 TRUE, TRUE, FALSE, FALSE, TRUE, FALSE, FALSE, FALSE, TRUE, FALSE) );
822 /** only a constraint whose "enforce is TRUE can be written in transformed problem */
823
824 /* add constraint to SCIP */
825 SCIP_CALL_ABORT( SCIPaddCons(scip, cons) );
826 assert(cons);
827 addedConss[c] = cons;
828 SCIP_CALL_ABORT( SCIPreleaseCons(scip, &cons) );
829 /* free temporary memory */
830 SCIPfreeBufferArray(scip, &vals);
831 SCIPfreeBufferArray(scip, &vars);
832 i++;
833 }
834
835 i = 0;
836 for(; c < (scipParaDiffSubproblem->getNBranchConsLinearConss()
837 + scipParaDiffSubproblem->getNBranchConsSetppcConss()
838 + scipParaDiffSubproblem->getNLinearConss()
839 + scipParaDiffSubproblem->getNBendersLinearConss()
840 + scipParaDiffSubproblem->getNBoundDisjunctions()) ; c++ )
841 {
842 SCIP_VAR** vars;
843 SCIP_BOUNDTYPE *types;
844 SCIP_Real* bounds;
845 int nVars = scipParaDiffSubproblem->getNVarsBoundDisjunction(i);
846 /* create array of variables, types and bounds */
847 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &vars, nVars) );
848 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &types, nVars) );
849 SCIP_CALL_ABORT( SCIPallocBufferArray(scip, &bounds, nVars) );
850
852 {
853 for( int v = 0; v < nVars; ++v )
854 {
855 vars[v] = orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIdxBoundDisjunctionVars(i,v)]]];
856 types[v] = scipParaDiffSubproblem->getBoundTypesBoundDisjunction(i,v);
857 bounds[v] = scipParaDiffSubproblem->getBoundsBoundDisjunction(i,v);
858 }
859 }
860 else
861 {
862 for( int v = 0; v < nVars; ++v )
863 {
864 vars[v] = orgVars[scipParaDiffSubproblem->getIdxBoundDisjunctionVars(i,v)];
865 types[v] = scipParaDiffSubproblem->getBoundTypesBoundDisjunction(i,v);
866 bounds[v] = scipParaDiffSubproblem->getBoundsBoundDisjunction(i,v);
867 }
868 }
869
870 /* create a constraint */
871 (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bdj%d", i);
872 SCIP_CALL_ABORT( SCIPcreateConsBounddisjunction(scip, &cons, consname, nVars, vars, types, bounds,
873 scipParaDiffSubproblem->getFlagBoundDisjunctionInitial(i),
874 scipParaDiffSubproblem->getFlagBoundDisjunctionSeparate(i),
875 scipParaDiffSubproblem->getFlagBoundDisjunctionEnforce(i),
876 scipParaDiffSubproblem->getFlagBoundDisjunctionCheck(i),
877 scipParaDiffSubproblem->getFlagBoundDisjunctionPropagate(i),
878 scipParaDiffSubproblem->getFlagBoundDisjunctionLocal(i),
879 scipParaDiffSubproblem->getFlagBoundDisjunctionModifiable(i),
880 scipParaDiffSubproblem->getFlagBoundDisjunctionDynamic(i),
881 scipParaDiffSubproblem->getFlagBoundDisjunctionRemovable(i),
882 scipParaDiffSubproblem->getFlagBoundDisjunctionStickingatnode(i) ) );
883 /* add constraint to SCIP */
884 SCIP_CALL_ABORT( SCIPaddCons(scip, cons) );
885 assert(cons);
886 addedConss[c] = cons;
887 SCIP_CALL_ABORT( SCIPreleaseCons(scip, &cons) );
888 /* free temporary memory */
889 SCIPfreeBufferArray(scip, &bounds);
890 SCIPfreeBufferArray(scip, &types);
891 SCIPfreeBufferArray(scip, &vars);
892 i++;
893 }
894 nAddedConss = c;
895 }
896
897
898 int addingConsParam = getParaParamSet()->getIntParamValue(AddDualBoundCons);
899 addedDualCons = 0;
900 if( addingConsParam != 0 )
901 {
902 if( ( addingConsParam == 1 && !SCIPisGT(scip, bbCurrentNode->getDualBoundValue(), bbCurrentNode->getInitialDualBoundValue()) )
903 || addingConsParam == 2 || addingConsParam == 3 )
904 {
905 SCIP_CONS* cons;
906 int nvars = SCIPgetNVars(scip);
907 SCIP_VAR **vars = SCIPgetVars(scip);
908 SCIP_Real* vals = new SCIP_Real[nvars];
909 for(int v = 0; v < nvars; ++v )
910 {
911 vals[v] = SCIPvarGetObj(vars[v]);
912 }
913 SCIP_CALL_ABORT( SCIPcreateConsLinear(scip, &cons, "objective", nvars, vars, vals, dualBoundValue, SCIPinfinity(scip),
914 TRUE, TRUE, FALSE, FALSE, TRUE, FALSE, FALSE, FALSE, TRUE, FALSE) );
915 assert( SCIPisEQ( scip, SCIPgetTransObjscale(scip), 1.0 ) );
916 assert( SCIPisZero( scip, SCIPgetTransObjoffset(scip) ) );
917
918 /** try to automatically convert a linear constraint into a more specific and more specialized constraint */
919
920 /* add constraint to SCIP */
921 SCIP_CALL_ABORT( SCIPaddCons(scip, cons) );
922 SCIP_CALL_ABORT( SCIPreleaseCons(scip, &cons) );
923 addedDualCons = cons;
924 }
925 }
926
927 if( userPlugins && scipParaDiffSubproblem )
928 {
929 userPlugins->newSubproblem(scip, scipParaDiffSubproblem->getBranchLinearConss(), scipParaDiffSubproblem->getBranchSetppcConss());
930 }
931
933 {
934 setWinnerRacingParams(winnerRacingParams); // winner parameters are set, again
935 // std::cout << winnerRacingParams->toString() << std::endl;
936 }
937
938 /// do not save solutions to original problem space
939 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "misc/transsolsorig", FALSE) );
940
941 if( SCIPgetStage(scip) == SCIP_STAGE_PROBLEM)
942 {
943 SCIP_CALL_ABORT( SCIPtransformProb(scip));
944 }
945
946 if( scipParaDiffSubproblem && scipParaDiffSubproblem->getNVarBranchStats() > 0 )
947 {
948 orgVars = SCIPgetOrigVars(scip); /* original problem's variables */
950 {
951 for( int i = 0; i < scipParaDiffSubproblem->getNVarBranchStats(); i++ )
952 {
953 SCIP_CALL_ABORT( SCIPinitVarBranchStats(scip, orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIdxLBranchStatsVars(i)]]],
954 scipParaDiffSubproblem->getDownpscost(i),
955 scipParaDiffSubproblem->getUppscost(i),
956 scipParaDiffSubproblem->getDownvsids(i),
957 scipParaDiffSubproblem->getUpvsids(i),
958 scipParaDiffSubproblem->getDownconflen(i),
959 scipParaDiffSubproblem->getUpconflen(i),
960 scipParaDiffSubproblem->getDowninfer(i),
961 scipParaDiffSubproblem->getUpinfer(i),
962 scipParaDiffSubproblem->getDowncutoff(i),
963 scipParaDiffSubproblem->getUpcutoff(i)
964 )
965 );
966 }
967 }
968 else
969 {
970 for( int i = 0; i < scipParaDiffSubproblem->getNVarBranchStats(); i++ )
971 {
972 SCIP_CALL_ABORT( SCIPinitVarBranchStats(scip, orgVars[scipParaDiffSubproblem->getIdxLBranchStatsVars(i)],
973 scipParaDiffSubproblem->getDownpscost(i),
974 scipParaDiffSubproblem->getUppscost(i),
975 scipParaDiffSubproblem->getDownvsids(i),
976 scipParaDiffSubproblem->getUpvsids(i),
977 scipParaDiffSubproblem->getDownconflen(i),
978 scipParaDiffSubproblem->getUpconflen(i),
979 scipParaDiffSubproblem->getDowninfer(i),
980 scipParaDiffSubproblem->getUpinfer(i),
981 scipParaDiffSubproblem->getDowncutoff(i),
982 scipParaDiffSubproblem->getUpcutoff(i)
983 )
984 );
985 }
986 }
987 }
988
989#if SCIP_VERSION >= 312
990 // std::cout << " VERSION >= 320 " << std::endl;
991
992 if( scipParaDiffSubproblem && scipParaDiffSubproblem->getNVarValueVars() > 0 )
993 {
994 orgVars = SCIPgetOrigVars(scip); /* original problem's variables */
996 {
997 for( int i = 0; i < scipParaDiffSubproblem->getNVarValueVars(); i++ )
998 {
999 for( int j = 0; j < scipParaDiffSubproblem->getNVarValueValues(i); j++ )
1000 {
1001 SCIP_CALL_ABORT( SCIPinitVarValueBranchStats(scip, orgVars[mapToProbIndecies[mapToSolverLocalIndecies[scipParaDiffSubproblem->getIdxLBranchStatsVars(i)]]],
1002 scipParaDiffSubproblem->getVarValue(i,j),
1003 scipParaDiffSubproblem->getVarValueDownvsids(i,j),
1004 scipParaDiffSubproblem->getVarVlaueUpvsids(i,j),
1005 scipParaDiffSubproblem->getVarValueDownconflen(i,j),
1006 scipParaDiffSubproblem->getVarValueUpconflen(i,j),
1007 scipParaDiffSubproblem->getVarValueDowninfer(i,j),
1008 scipParaDiffSubproblem->getVarValueUpinfer(i,j),
1009 scipParaDiffSubproblem->getVarValueDowncutoff(i,j),
1010 scipParaDiffSubproblem->getVarValueUpcutoff(i,j)
1011 )
1012 );
1013 /*
1014 std::cout << mapToOriginalIndecies[scipParaDiffSubproblem->getIdxLBranchStatsVars(i)]
1015 << ", "
1016 << scipParaDiffSubproblem->getVarValue(i,j)
1017 << ", "
1018 << scipParaDiffSubproblem->getVarValueDownvsids(i,j)
1019 << ", "
1020 << scipParaDiffSubproblem->getVarVlaueUpvsids(i,j)
1021 << ", "
1022 << scipParaDiffSubproblem->getVarValueDownconflen(i,j)
1023 << ", "
1024 << scipParaDiffSubproblem->getVarValueUpconflen(i,j)
1025 << ", "
1026 << scipParaDiffSubproblem->getVarValueDowninfer(i,j)
1027 << ", "
1028 << scipParaDiffSubproblem->getVarValueUpinfer(i,j)
1029 << ", "
1030 << scipParaDiffSubproblem->getVarValueDowncutoff(i,j)
1031 << ", "
1032 << scipParaDiffSubproblem->getVarValueUpcutoff(i,j)
1033 << std::endl;
1034 */
1035 }
1036 }
1037 }
1038 else
1039 {
1040 for( int i = 0; i < scipParaDiffSubproblem->getNVarValueVars(); i++ )
1041 {
1042 for( int j = 0; j < scipParaDiffSubproblem->getNVarValueValues(i); j++ )
1043 {
1044 SCIP_CALL_ABORT( SCIPinitVarValueBranchStats(scip, orgVars[scipParaDiffSubproblem->getIdxLBranchStatsVars(i)],
1045 scipParaDiffSubproblem->getVarValue(i,j),
1046 scipParaDiffSubproblem->getVarValueDownvsids(i,j),
1047 scipParaDiffSubproblem->getVarVlaueUpvsids(i,j),
1048 scipParaDiffSubproblem->getVarValueDownconflen(i,j),
1049 scipParaDiffSubproblem->getVarValueUpconflen(i,j),
1050 scipParaDiffSubproblem->getVarValueDowninfer(i,j),
1051 scipParaDiffSubproblem->getVarValueUpinfer(i,j),
1052 scipParaDiffSubproblem->getVarValueDowncutoff(i,j),
1053 scipParaDiffSubproblem->getVarValueUpcutoff(i,j)
1054 )
1055 );
1056 }
1057 }
1058 }
1059 /** for debug *********************
1060 int nvars; ** number of variables
1061 int nbinvars; ** number of binary variables
1062 int nintvars; ** number of integer variables
1063 SCIP_VAR** vars; ** transformed problem's variables
1064 SCIP_CALL_ABORT( SCIPgetVarsData(scip, &vars, &nvars, &nbinvars, &nintvars, NULL, NULL) );
1065 int ngenvars = nbinvars+nintvars;
1066 int nOrgVarst = 0;
1067
1068 if( ngenvars > 0 )
1069 {
1070 std::cout << "R." << paraComm->getRank() << ", ngenvars = " << ngenvars << std::endl;;
1071 for( int i = 0; i < ngenvars; i++ )
1072 {
1073 assert( SCIPvarGetType(vars[i]) == SCIP_VARTYPE_BINARY || SCIPvarGetType(vars[i]) == SCIP_VARTYPE_INTEGER );
1074
1075 SCIP_VAR *transformVar = vars[i];
1076 SCIP_Real scalar = 1.0;
1077 SCIP_Real constant = 0.0;
1078 SCIP_CALL_ABORT( SCIPvarGetOrigvarSum(&transformVar, &scalar, &constant ) );
1079 assert(transformVar != NULL);
1080
1081 if( transformVar ) // The variable in the transformed space
1082 {
1083 SCIP_VALUEHISTORY* valuehistory = SCIPvarGetValuehistory(vars[i]);
1084 if( valuehistory != NULL )
1085 {
1086 nOrgVarst++;
1087 }
1088 else
1089 {
1090 std::cout << "R." << paraComm->getRank() << ", no history for var i = " << i << std::endl;
1091 std::cout << "R." << paraComm->getRank() << ", org = " << SCIPvarGetValuehistory(transformVar) << std::endl;
1092 }
1093 }
1094 else
1095 {
1096 std::cout << "R." << paraComm->getRank() << ", no transfrom var i = " << i << std::endl;;
1097 }
1098 }
1099 }
1100 if( nOrgVarst == 0 )
1101 {
1102 std::cout << "Failed to set Var Value stat. R." << paraComm->getRank() << std::endl;;
1103 abort();
1104 }
1105 else
1106 {
1107 std::cout << "Set " << nOrgVarst << " Var Value stat. R." << paraComm->getRank() << std::endl;;
1108 }
1109 *** end of debug */
1110 }
1111#endif
1112
1113// if( userPlugins && scipParaDiffSubproblem)
1114// {
1115// userPlugins->newSubproblem(scip, scipParaDiffSubproblem->getBranchLinearConss(), scipParaDiffSubproblem->getBranchSetppcConss());
1116// }
1117
1118}
1119
1120void
1122 )
1123{
1125 {
1126 DEF_SCIP_PARA_COMM( scipParaComm, paraComm);
1127 ScipParaInitialStat *initialStat = scipParaComm->createScipParaInitialStat(scip);
1128 initialStat->send(paraComm, 0);
1129 delete initialStat;
1130 }
1131 SCIP_CALL_ABORT( SCIPfreeTransform(scip) );
1132
1133 int c = 0;
1135 {
1136 ScipParaDiffSubproblem *scipParaDiffSubproblem = dynamic_cast< ScipParaDiffSubproblem* >(currentTask->getDiffSubproblem());
1137 if( scipParaDiffSubproblem->getNLinearConss() > 0 )
1138 {
1139 for(; c < scipParaDiffSubproblem->getNLinearConss() ; c++ )
1140 {
1141 if( !SCIPconsIsDeleted(addedConss[c]) )
1142 {
1143 SCIP_CALL_ABORT( SCIPdelCons(scip, addedConss[c]) );
1144 }
1145 }
1146 }
1147 }
1148
1149 if( addedConss )
1150 {
1151 for(; c < nAddedConss; c++ )
1152 {
1153 if( !SCIPconsIsDeleted(addedConss[c]) )
1154 {
1155 SCIP_CALL_ABORT( SCIPdelCons(scip, addedConss[c]) );
1156 }
1157 }
1158 delete [] addedConss;
1159 addedConss = 0;
1160 }
1161
1162 if( addedDualCons )
1163 {
1164 SCIP_CALL_ABORT( SCIPdelCons(scip, addedDualCons) );
1165 addedDualCons = 0;
1166 }
1167
1168 SCIP_VAR **orgVars = SCIPgetOrigVars(scip); // variables are indexed by index
1169 // Taking into account multi-aggregate vars.
1170 // int n = SCIPgetNOrigVars(scip); // the number of original variables
1171 // assert( n == nOrgVars );
1172 assert( nOrgVarsInSolvers == SCIPgetNOrigVars(scip));
1173 // for( int v = 0; v < n; v++ )
1174 for( int v = 0; v < nOrgVars; v++ )
1175 {
1176 SCIP_CALL_ABORT( SCIPchgVarLbGlobal( scip,orgVars[v], orgVarLbs[v] ) );
1177 SCIP_CALL_ABORT( SCIPchgVarUbGlobal( scip,orgVars[v], orgVarUbs[v] ) );
1178 }
1179
1180 if( racingWinner )
1181 {
1184 // std::cout << "Winner: nTightened = " << nTightened << ", nTightenedInt = " << nTightenedInt << std::endl;
1185 }
1186
1190 {
1192 }
1193
1194}
1195
1196void
1198 )
1199{
1200 DEF_SCIP_PARA_COMM( scipParaComm, paraComm);
1201 SCIP_SOL *sol = SCIPgetBestSol(scip);
1202 int nVars = SCIPgetNOrigVars(scip);
1203 SCIP_VAR **vars = SCIPgetOrigVars(scip);
1204 SCIP_Real *vals = new SCIP_Real[nVars];
1205 SCIP_CALL_ABORT( SCIPgetSolVals(scip, sol, nVars, vars, vals) );
1207 {
1208 SCIP_VAR **varsInOrig = new SCIP_VAR*[nVars];
1209 SCIP_Real *valsInOrig = new SCIP_Real[nVars]();
1210 int nVarsInOrig = 0;
1211 for( int i = 0; i < nVars; i++ )
1212 {
1213 if( getOriginalIndex(SCIPvarGetIndex(vars[i])) >= 0 )
1214 {
1215 varsInOrig[nVarsInOrig] = vars[i];
1216 valsInOrig[nVarsInOrig] = vals[i];
1217 nVarsInOrig++;
1218 }
1219 }
1221 scipParaComm->createScipParaSolution(
1222 this,
1223 SCIPgetSolOrigObj(scip, sol),
1224 nVarsInOrig,
1225 varsInOrig,
1226 valsInOrig
1227 )
1228 );
1229 delete [] varsInOrig;
1230 delete [] valsInOrig;
1231 }
1232 else
1233 {
1235 scipParaComm->createScipParaSolution(
1236 this,
1237 SCIPgetSolOrigObj(scip, sol),
1238 nVars,
1239 vars,
1240 vals
1241 )
1242 );
1243 }
1244 delete [] vals;
1245}
1246
1247void
1249 )
1250{
1251
1253 {
1255 }
1256
1257 // if( paraParams->getBoolParamValue(UG::CheckGapInLC) )
1258 // {
1259 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "limits/gap", 0.0 ) );
1260 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "limits/absgap", 0.0 ) );
1261 // }
1262
1263 // if( paraParams->getBoolParamValue(UG::CheckFeasibilityInLC) )
1264 // {
1265 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "numerics/feastol", (orgFeastol/10.0) ) );
1266 if( SCIP_APIVERSION < 61 )
1267 {
1268 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "numerics/lpfeastol", (orgLpfeastol/10.0) ) );
1269 }
1270 /*
1271 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "numerics/feastol", (orgFeastol/100.0) ) );
1272 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "numerics/lpfeastol", (orgLpfeastol/100.0) ) );
1273 */
1274 // }
1275
1276 /** solve */
1278 {
1279 assert(conflictConsList->size() == 0);
1280 }
1281
1282 /* don't catch control+c */
1283 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "misc/catchctrlc", FALSE) );
1284
1285 /** set cutoff value */
1286#ifdef UG_DEBUG_SOLUTION
1287 const double limit = DBL_MAX;
1288 if( limit < globalBestIncumbentValue )
1289 {
1290 SCIP_CALL_ABORT( SCIPsetObjlimit(scip, limit) );
1291 }
1292 else
1293 {
1294 SCIP_CALL_ABORT( SCIPsetObjlimit(scip, globalBestIncumbentValue) );
1295 }
1296 // if( ( !currentNode->getDiffSubproblem() ) ||
1297 // currentNode->getDiffSubproblem()->isOptimalSolIncluded() )
1298 // {
1299 // writeSubproblem();
1300 // }
1302#else
1303 /** set cutoff value */
1304 SCIP_CALL_ABORT( SCIPsetObjlimit(scip, globalBestIncumbentValue) );
1305#endif
1307
1308#ifdef _DEBUG_DET
1310#endif
1311
1312 // if( paraComm->getRank() == 1) checkVarsAndIndex("***before solve***",scip);
1313
1314 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "misc/usesymmetry", 0 ) ); // Symmetry handling technique is explicitly turn off in Solever for this version (SCIP 5.0)
1315
1316 SCIP_CALL_ABORT( SCIPsetIntParam(scip,"timing/clocktype", 2) ); // always wall clock time, racing may change clocktype
1317 // writeSubproblem();
1318
1319#ifdef UG_DEBUG_SOLUTION
1320 // assert( paraNode == currentNode );
1321 if( ( !currentTask->getDiffSubproblem() ) ||
1322 currentTask->getDiffSubproblem()->isOptimalSolIncluded() )
1323 {
1324 SCIPdebugSolEnable(scip);
1325 std::cout << "R." << paraComm->getRank() << ": enable debug" << std::endl;
1326 assert( SCIPdebugSolIsEnabled(scip) == TRUE );
1327 }
1328 else
1329 {
1330 SCIPdebugSolDisable(scip);
1331 std::cout << "R." << paraComm->getRank() << ": disable debug" << std::endl;
1332 assert( SCIPdebugSolIsEnabled(scip) == FALSE );
1333 }
1334#endif
1335
1336#if (SCIP_VERSION >= 700)
1338 {
1339 if( !isRacingStage() )
1340 {
1341 SCIP_CALL_ABORT( SCIPsetCharParam(scip, "estimation/restarts/restartpolicy", 'n' ) );
1342 }
1343 }
1344#endif
1345
1347 {
1348 double timeRemains = std::max(0.0, (paraParams->getRealParamValue(UG::TimeLimit) - paraTimer->getElapsedTime()) );
1349 SCIP_CALL_ABORT( SCIPsetIntParam(scip,"timing/clocktype", 2) ); // to confirm this for the scip environment actually to work
1350 SCIP_CALL_ABORT( SCIPsetRealParam(scip,"limits/time", timeRemains) );
1351 }
1352
1353#if SCIP_APIVERSION >= 101
1354 if( SCIPgetParam(scip, "presolving/milp/threads") != NULL )
1355 {
1356 int nmilpthreads = 0;
1357 SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/milp/threads", &nmilpthreads) );
1358 assert( nmilpthreads == 1 );
1359 }
1360#endif
1361
1362 terminationMode = UG::CompTerminatedNormally; // assume that it terminates normally
1363
1364 // writeSubproblem();
1365 SCIP_RETCODE ret = SCIPsolve(scip);
1366 if( ret != SCIP_OKAY )
1367 {
1368#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
1369 SCIPprintError(ret, NULL);
1370#else
1371 SCIPprintError(ret);
1372#endif
1374 if( userPlugins )
1375 {
1377 }
1378#ifdef _MSC_VER
1379 _sleep(10);
1380#else
1381 sleep(10);
1382#endif
1383 std::cout << "ret = " << (int)ret << std::endl;
1384 THROW_LOGICAL_ERROR1("SCIP terminate with NO SCIP_OKAY");
1385 }
1386
1387 // Notification message has to complete
1389 {
1391 }
1392
1393 // Then, solver status should be checked
1394 SCIP_STATUS status = SCIPgetStatus(scip);
1395
1396#if SCIP_APIVERSION >= 101
1397 if( paraParams->getIntParamValue(UG::RampUpPhaseProcess) == 3 && // self-split ramp-up
1398 currentTask->isRootTask() // the following should do only for the self-split ramp-up procedure
1399 )
1400 {
1401 if( status == SCIP_STATUS_OPTIMAL )
1402 {
1404 }
1405 if( selfSplitNodesel->inSampling() && status == SCIP_STATUS_OPTIMAL ) // it is still sampling. This means the poroblem was solved
1406 {
1407 if( paraComm->getRank() == 1 )
1408 {
1410 }
1411 return;
1412 }
1413 else
1414 {
1415 sendLocalSolution(); // solution was not sent during sampling to generate the same search tree
1416 }
1417 }
1418#endif
1419
1420// if( status == SCIP_STATUS_OPTIMAL ) // when sub-MIP is solved at root node, the solution may not be saved
1421 if( SCIPgetNSols(scip) > 0 )
1422 {
1425 }
1426// else
1427// {
1428#ifdef UG_DEBUG_SOLUTION
1429 if( status != SCIP_STATUS_OPTIMAL )
1430 {
1431 if( SCIPdebugSolIsEnabled(scip) &&
1432 // currentTask->getMergingStatus() != 3 &&
1435 {
1436 std::cout << "R" << paraComm->getRank() << " solver lost optimal solution." << std::endl;
1437 throw "Optimal solution lost!";
1438 }
1439 }
1440#endif
1441// if( status == SCIP_STATUS_MEMLIMIT )
1442// {
1443// std::cout << "Warning: SCIP was interrupted because the memory limit was reached" << std::endl;
1444// abort();
1445// }
1446// }
1447
1448// std::cout << "R" << paraComm->getRank() << " status = " << (int)status << std::endl;
1449
1450 if( status == SCIP_STATUS_OPTIMAL ||
1451 status == SCIP_STATUS_GAPLIMIT )
1452
1453 {
1454 solverDualBound = SCIPgetDualbound(scip);
1455 // if( status == SCIP_STATUS_OPTIMAL )
1456 // {
1457 // current best solution may not be accepted in LC
1458 // solverDualBound = std::max(solverDualBound, getGlobalBestIncumbentValue() );
1459 // should not do above, since globalBestIncumbentValue might not be updated
1460 // }
1461
1462 /*
1463 double dualBound = SCIPgetDualbound(scip);
1464 if( isRacingStage() )
1465 {
1466 if( dualBound > solverDualBound )
1467 {
1468 solverDualBound = dualBound;
1469 }
1470 }
1471 else
1472 {
1473 solverDualBound = dualBound;
1474 // maximalDualBound = std::max(dualBound, currentNode->getDualBoundValue());
1475 if( EPSEQ( maximalDualBound, -DBL_MAX, eps ) )
1476 {
1477 maximalDualBound = std::max(dualBound, currentNode->getDualBoundValue());
1478 }
1479 else
1480 {
1481 if( !SCIPisInfinity(scip, -dualBound) &&
1482 dualBound > currentNode->getDualBoundValue() &&
1483 dualBound < maximalDualBound )
1484 {
1485 maximalDualBound = dualBound;
1486 }
1487 else
1488 {
1489 if( maximalDualBound < currentNode->getDualBoundValue() )
1490 {
1491 maximalDualBound = currentNode->getDualBoundValue();
1492 }
1493 }
1494 }
1495 }
1496 */
1497 }
1498 else if( status == SCIP_STATUS_INFEASIBLE )
1499 {
1500 if( EPSEQ(globalBestIncumbentValue, DBL_MAX, eps ) )
1501 {
1502 solverDualBound = SCIPgetDualbound(scip);
1503 }
1504 else
1505 {
1507 }
1508 }
1509 else
1510 {
1511 if( status == SCIP_STATUS_NODELIMIT )
1512 {
1513 throw "SCIP terminated with SCIP_STATUS_NODELIMIT";
1514 }
1515 else if( status == SCIP_STATUS_TOTALNODELIMIT )
1516 {
1517 throw "SCIP terminated with SCIP_STATUS_TOTALNODELIMIT";
1518 }
1519 else if( status == SCIP_STATUS_STALLNODELIMIT )
1520 {
1521 throw "SCIP terminated with SCIP_STATUS_STALLNODELIMIT";
1522 }
1523 else if( status == SCIP_STATUS_TIMELIMIT )
1524 {
1525 // throw "SCIP terminated with SCIP_STATUS_TIMELIMIT";
1526 solverDualBound = SCIPgetDualbound(scip);
1528 if( isRacingStage() )
1529 {
1530 racingIsInterrupted = true;
1531 }
1532 }
1533 else if( status == SCIP_STATUS_MEMLIMIT )
1534 {
1535 memoryLimitIsReached = true;
1536 }
1537 else if( status == SCIP_STATUS_SOLLIMIT )
1538 {
1539 throw "SCIP terminated with SCIP_STATUS_SOLLIMIT";
1540 }
1541 else if( status == SCIP_STATUS_BESTSOLLIMIT )
1542 {
1543 throw "SCIP terminated with SCIP_STATUS_BESTSOLLIMIT";
1544 }
1545 else if( status == SCIP_STATUS_USERINTERRUPT && SCIPisObjIntegral(scip) )
1546 {
1547 if( SCIPfeasCeil(scip, dynamic_cast<UG::BbParaNode *>(getCurrentNode())->getDualBoundValue()) >= getGlobalBestIncumbentValue() )
1548 {
1549 solverDualBound = SCIPfeasCeil(scip, dynamic_cast<UG::BbParaNode *>(getCurrentNode())->getDualBoundValue());
1550 }
1551 else
1552 {
1553 solverDualBound = std::max(dynamic_cast<UG::BbParaNode *>(getCurrentNode())->getDualBoundValue(), SCIPgetDualbound(scip));
1554 }
1555 }
1556 else if( status == SCIP_STATUS_USERINTERRUPT )
1557 {
1558 solverDualBound = std::max(dynamic_cast<UG::BbParaNode *>(getCurrentNode())->getDualBoundValue(), SCIPgetDualbound(scip));
1559 }
1560 else
1561 {
1562 solverDualBound = -DBL_MAX;
1563 }
1564 /*
1565 if( maximalDualBound < currentNode->getDualBoundValue() )
1566 {
1567 maximalDualBound = currentNode->getDualBoundValue();
1568 }
1569 */
1570 }
1571
1572 /*
1573 double dualBound = SCIPgetDualbound(scip);
1574 if( isRacingStage() )
1575 {
1576 if( dualBound > maximalDualBound )
1577 {
1578 maximalDualBound = dualBound;
1579 }
1580 }
1581 else
1582 {
1583 if( !SCIPisInfinity(scip, -dualBound) &&
1584 dualBound < maximalDualBound )
1585 {
1586 maximalDualBound = dualBound;
1587 }
1588 }
1589 */
1590
1591
1592 if( conflictConsList && conflictConsList->size() > 0 )
1593 {
1594 int nConfilcts = conflictConsList->size();
1595 for(int i = 0; i < nConfilcts; i++ )
1596 {
1597 assert(!conflictConsList->empty());
1598 LocalNodeInfo *info= conflictConsList->front();
1599 conflictConsList->pop_front();
1600 if( info->linearCoefs ) delete[] info->linearCoefs;
1601 if( info->idxLinearCoefsVars ) delete[] info->idxLinearCoefsVars;
1602 delete info;
1603 }
1604 }
1605
1606#if SCIP_APIVERSION >= 101
1607 if( paraParams->getIntParamValue(UG::RampUpPhaseProcess) == 3 && // self-split ramp-up
1609 {
1610 int numnodesels = SCIPgetNNodesels( scip );
1611 SCIP_NODESEL** nodesels = SCIPgetNodesels( scip );
1612 int i;
1613 for( i = 0; i < numnodesels; ++i )
1614 {
1615 std::string nodeselname(SCIPnodeselGetName(nodesels[i]));
1616 if( std::string(nodeselname) == std::string("ScipParaObjSelfSplitNodeSel") )
1617 {
1618 break;
1619 }
1620 }
1621 assert( i != numnodesels );
1622 SCIP_CALL_ABORT( SCIPsetNodeselStdPriority(scip, nodesels[i], -INT_MAX/4 ) );
1623 }
1624#endif
1625
1626}
1627
1628long long
1630 )
1631{
1632 if( SCIPgetStage(scip) == SCIP_STAGE_SOLVING || SCIPgetStage(scip) == SCIP_STAGE_SOLVED )
1633 {
1634 return SCIPgetNTotalNodes(scip);
1635 }
1636 else
1637 {
1638 return 0;
1639 }
1640}
1641
1642int
1644 )
1645{
1646 if( SCIPgetStage(scip) == SCIP_STAGE_SOLVING || SCIPgetStage(scip) == SCIP_STAGE_SOLVED )
1647 {
1648 return SCIPgetNNodesLeft(scip);
1649 }
1650 else
1651 {
1652 if( SCIPgetStage(scip) >= SCIP_STAGE_PRESOLVING && SCIPgetStage(scip) <= SCIP_STAGE_INITSOLVE )
1653 {
1654 return 1;
1655 }
1656 else
1657 {
1658 return 0;
1659 }
1660 }
1661}
1662
1663double
1665 )
1666{
1667 if( SCIPgetStage(scip) == SCIP_STAGE_PRESOLVING || SCIPgetStage(scip) == SCIP_STAGE_INITSOLVE )
1668 {
1669 return dynamic_cast<UG::BbParaNode *>(currentTask)->getDualBoundValue();
1670 }
1671 else
1672 {
1673 return SCIPgetDualbound(scip);
1674 }
1675}
1676
1678 int argc,
1679 char **argv,
1680 int inNhanders, // number of valid message handlers
1682 UG::ParaParamSet *inParaParamSet,
1683 UG::ParaInstance *inParaInstance,
1684 UG::ParaDeterministicTimer *inDetTimer
1685 ) : BbParaSolver(argc, argv, inNhanders, comm, inParaParamSet, inParaInstance, inDetTimer),
1686 messagehdlr(0),
1687 logfile(0),
1688 originalParamSet(0),
1689 conflictConsList(0),
1690 userPlugins(0),
1691 commPointHdlr(0),
1692 nodesel(0),
1693#if SCIP_APIVERSION >= 101
1694 selfSplitNodesel(0),
1695#endif
1696 scipPropagator(0),
1697 interruptMsgMonitor(0),
1698 nPreviousNodesLeft(0),
1699 originalPriority(0),
1700 nOrgVars(0),
1701 nOrgVarsInSolvers(0),
1702 orgVarLbs(0),
1703 orgVarUbs(0),
1704 tightenedVarLbs(0),
1705 tightenedVarUbs(0),
1706 mapToOriginalIndecies(0),
1707 mapToSolverLocalIndecies(0),
1708 mapToProbIndecies(0),
1709 // stuffingMaxrounds(0),
1710 // domcolMaxrounds(0),
1711 // dualcompMaxrounds(0),
1712 // dualinferMaxrounds(0),
1713 // dualaggMaxrounds(0),
1714 // abspowerDualpresolve(0),
1715 // andDualpresolving(0),
1716 // cumulativeDualpresolve(0),
1717 // knapsackDualpresolving(0),
1718 // linearDualpresolving(0),
1719 // setppcDualpresolving(0),
1720 // logicorDualpresolving(0),
1721 miscAllowdualreds(0),
1722 nAddedConss(0),
1723 addedConss(0),
1724 addedDualCons(0),
1725 settingsNameLC(0),
1726 fiberSCIP(false),
1727 quiet(false),
1728 collectingModeIsProhibited(false),
1729 problemFileName(0),
1730 orgFeastol(0.0),
1731 orgLpfeastol(0.0),
1732 copyIncreasedVariables(false)
1733{
1734
1735 // ScipMessageHandlerFunctionPointer *scipMessageHandler = reinterpret_cast<ScipMessageHandlerFunctionPointer *>(messageHandler);
1736 // no additional message handlers
1737
1738 char* logname = NULL;
1739
1740 ScipParaInstance *scipParaInstance = dynamic_cast< ScipParaInstance *>(paraInstance);
1741
1742 /* Initialize the SCIP environment */
1743 /*********
1744 * Setup *
1745 *********/
1746 /* initialize SCIP */
1747 SCIP_CALL_ABORT( SCIPcreate(&scip) );
1748 SCIP_CALL_ABORT( SCIPsetIntParam(scip,"timing/clocktype", 2) ); // always wall clock time
1750 {
1751 double timeRemains = std::max(0.0, (paraParams->getRealParamValue(UG::TimeLimit) - paraTimer->getElapsedTime() + 3.0) ); // 3.0: timming issue
1752 SCIP_CALL_ABORT( SCIPsetRealParam(scip,"limits/time", timeRemains) );
1753 }
1754 /** set user plugins */
1755 ::setUserPlugins(this);
1757 {
1758 /* include default SCIP plugins */
1759 SCIP_CALL_ABORT( SCIPincludeDefaultPlugins(scip) );
1760 /** include user plugins */
1762 }
1763 /* include communication point handler */
1765 commPointHdlr = new ScipParaObjCommPointHdlr(paraComm, this, updator);
1766 SCIP_CALL_ABORT( SCIPincludeObjEventhdlr(scip, commPointHdlr, TRUE) );
1767 SCIP_CALL_ABORT( SCIPincludeObjHeur(scip, updator, TRUE) );
1768
1769 /* include propagator */
1773 {
1775 SCIP_CALL_ABORT( SCIPincludeObjProp(scip, scipPropagator, TRUE) );
1778 }
1779
1780 /* include node selector */
1781 nodesel = new ScipParaObjNodesel(this);
1782 SCIP_CALL_ABORT( SCIPincludeObjNodesel(scip, nodesel, TRUE) );
1783#if SCIP_APIVERSION >= 101
1784 if( paraParams->getIntParamValue(UG::RampUpPhaseProcess) == 3 ) // self-split ramp-up
1785 {
1786 selfSplitNodesel = new ScipParaObjSelfSplitNodesel(
1787 paraComm->getRank() - 1,
1788 paraComm->getSize() - 1,
1790 paraComm,
1791 this,
1792 scip
1793 );
1794 SCIP_CALL_ABORT( SCIPincludeObjNodesel(scip, selfSplitNodesel, TRUE) );
1795 }
1796#endif
1797 /* include branch rule plugins */
1798 SCIP_CALL_ABORT( SCIPincludeObjBranchrule(scip, new ScipParaObjBranchRule(this), TRUE) );
1799
1800 if( inParaParamSet->getBoolParamValue(UG::TransferConflictCuts) )
1801 {
1802 conflictConsList = new std::list<LocalNodeInfoPtr>;
1803 SCIP_CONFLICTHDLRDATA *conflictHdrData = reinterpret_cast< SCIP_CONFLICTHDLRDATA * >(this);
1804 /* create conflict handler to collects conflicts */
1805#if SCIP_VERSION == 211 && SCIP_SUBVERSION == 0
1806 SCIP_CALL_ABORT( SCIPincludeConflicthdlr(scip, CONFLICTHDLR_NAME, CONFLICTHDLR_DESC, CONFLICTHDLR_PRIORITY,
1807 NULL, NULL, NULL, NULL, NULL, NULL, conflictExecCollector, conflictHdrData) );
1808#else
1809 SCIP_CALL_ABORT( SCIPincludeConflicthdlrBasic(scip, NULL, CONFLICTHDLR_NAME, CONFLICTHDLR_DESC, CONFLICTHDLR_PRIORITY,
1810 conflictExecCollector, conflictHdrData) );
1811#endif
1812 }
1813
1814 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/feastol", &orgFeastol ) );
1815 if( SCIP_APIVERSION < 61 )
1816 {
1817 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/lpfeastol", &orgLpfeastol ) );
1818 }
1819
1820 /********************
1821 * Parse parameters *
1822 ********************/
1823 for( int i = 3; i < argc; ++i ) /** the first argument is runtime parameter file for ParaSCIP */
1824 {
1825 if( strcmp(argv[i], "-l") == 0 )
1826 {
1827 i++;
1828 if( i < argc )
1829 logname = argv[i];
1830 else
1831 {
1832 THROW_LOGICAL_ERROR1("missing log filename after parameter '-l'");
1833 }
1834 }
1835 else if( strcmp(argv[i], "-q") == 0 )
1836 quiet = true;
1837 // other arguments are omitted in Solver
1838 }
1839
1840 /***********************************
1841 * create log file message handler *
1842 ***********************************/
1844 {
1845 // SCIP_CALL_ABORT( SCIPsetMessagehdlr(NULL) );
1846#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
1847 SCIP_CALL_ABORT( SCIPcreateObjMessagehdlr(&messagehdlr, new ScipParaObjMessageHdlr(paraComm, NULL, TRUE, FALSE), TRUE) );
1848 SCIP_CALL_ABORT( SCIPsetMessagehdlr(messagehdlr) );
1849#else
1850 SCIPsetMessagehdlrQuiet(scip, TRUE);
1851#endif
1852 }
1853 else
1854 {
1855 if( logname != NULL || quiet )
1856 {
1857 if( logname != NULL )
1858 {
1859 std::ostringstream os;
1860 os << logname << comm->getRank();
1861 logfile = fopen(os.str().c_str(), "a"); // append to log file */
1862 if( logfile == NULL )
1863 {
1864 THROW_LOGICAL_ERROR3("cannot open log file <", logname, "> for writing");
1865 }
1866 }
1867 SCIP_CALL_ABORT( SCIPcreateObjMessagehdlr(&messagehdlr, new ScipParaObjMessageHdlr(paraComm, logfile, quiet, FALSE), TRUE) );
1868#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
1869 SCIP_CALL_ABORT( SCIPsetMessagehdlr(messagehdlr) );
1870#else
1871 SCIP_CALL_ABORT( SCIPsetMessagehdlr(scip, messagehdlr) );
1872 SCIP_CALL_ABORT( SCIPmessagehdlrRelease(&messagehdlr));
1873#endif
1874 }
1875 }
1876
1877 DEF_SCIP_PARA_COMM( scipParaComm, paraComm );
1878 scipDiffParamSetRoot = scipParaComm->createScipDiffParamSet();
1879 scipDiffParamSetRoot->bcast(comm, 0); /** this bcast is sent as SolverInitializationMessage */
1880 scipDiffParamSet = scipParaComm->createScipDiffParamSet();
1881 scipDiffParamSet->bcast(comm, 0); /** this bcast is sent as SolverInitializationMessage */
1882 int tempIsWarmStarted;
1883 comm->bcast(&tempIsWarmStarted, 1, UG::ParaINT, 0);
1884 warmStarted = (tempIsWarmStarted == 1);
1887 {
1888 int solutionExists = 0;
1889 comm->bcast(&solutionExists, 1, UG::ParaINT, 0);
1890 if( solutionExists )
1891 {
1894 }
1895 }
1896
1897 /** set parameters for SCIP: this values are reseted before solving */
1898 /* move to the place after problem has been created because the parameters may be changed.
1899 scipDiffParamSet->setParametersInScip(scip);
1900 SCIP_Real epsilon;
1901 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon));
1902 eps = epsilon;
1903 */
1904
1905 char *isolname = 0;
1906 for( int i = 3; i < argc; ++i ) /** the first argument is runtime parameter file for ParaSCIP */
1907 /** the second argument is problem file name */
1908 {
1909 if( strcmp(argv[i], "-sl") == 0 )
1910 {
1911 i++;
1912 if( i < argc )
1913 {
1914 settingsNameLC = argv[i];
1915 break;
1916 }
1917 else
1918 {
1919 std::cerr << "missing settings filename after parameter '-sl'" << std::endl;
1920 exit(1);
1921 }
1922 }
1923 else if ( strcmp(argv[i], "-isol") == 0 )
1924 {
1925 i++;
1926 if( i < argc )
1927 {
1928 isolname = argv[i];
1929 }
1930 else
1931 {
1932 std::cerr << "missing settings filename after parameter '-isol'" << std::endl;
1933 exit(1);
1934 }
1935 }
1936 }
1937
1938#if( !defined(UG_QUBO) && !defined(UG_SMOOTHIE) )
1939
1940 /** create problem */
1941 scipParaInstance->createProblem(scip,
1948 isolname
1949 );
1950
1951 if( SCIPgetStage(scip) == SCIP_STAGE_INIT && paraParams->getIntParamValue(UG::InstanceTransferMethod) == 2)
1952 {
1953 delete paraComm;
1954 exit(0); // the problem should be solved in LC
1955 }
1956
1957 /** set parameters for SCIP: this values are reseted before solving */
1960 SCIP_Real epsilon;
1961 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon));
1962 eps = epsilon;
1963
1966 {
1967 /* initialize SCIP to check root solvability */
1968 SCIP_CALL_ABORT( SCIPcreate(&scipToCheckEffectOfRootNodeProcesses) );
1969 /* include default SCIP plugins */
1970 SCIP_CALL_ABORT( SCIPincludeDefaultPlugins(scipToCheckEffectOfRootNodeProcesses) );
1971 /* include scipParaConshdlr plugins */
1979 isolname
1980 );
1982 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "presolving/maxrestarts", 0) );
1983 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "presolving/maxrounds", 0) );
1984 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/linear/presolpairwise", FALSE) );
1985 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/and/presolpairwise", FALSE) );
1986 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/logicor/presolpairwise", FALSE) );
1987 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/setppc/presolpairwise", FALSE) );
1988 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "propagating/probing/maxprerounds", 0) );
1989 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "heuristics/feaspump/freq", -1) );
1990 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "heuristics/rens/freq", -1) );
1991 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "separating/maxcutsroot", 100) );
1992 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "separating/maxroundsroot", 5) );
1993 }
1994 delete paraInstance;
1995 paraInstance = 0;
1996
1997 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "limits/memory", dynamic_cast<ScipParaParamSet *>(paraParams)->getRealParamValue(MemoryLimit)) );
1998// if( paraComm->getRank() == 1 )
1999// {
2000// std::cout << "*** set memory limit to " << dynamic_cast<ScipParaParamSet *>(paraParams)->getRealParamValue(MemoryLimit) << " for each SCIP ***" << std::endl;
2001// }
2002
2003 /** save original priority of changing node selector */
2004 assert(scip);
2006
2007#ifndef SCIP_EVENTTYPE_COMM
2008 // std::cout << "#### SCIP_EVENTTYPE_COMM is not defined!" << std::endl;
2010 {
2011 interruptMsgMonitor = new ScipParaInterruptMsgMonitor(scipParaComm, this);
2012 // interruptMsgMonitorThread = std::thread( runInterruptMsgMonitorThread, interruptMsgMonitor );
2014 t.detach();
2015 }
2016#endif
2017
2018#endif
2019}
2020
2022 int argc,
2023 char **argv,
2024 int inNhanders, // number of valid message handlers
2026 UG::ParaParamSet *inParaParamSet,
2027 UG::ParaInstance *inParaInstance,
2028 UG::ParaDeterministicTimer *inDetTimer,
2029 double timeOffset,
2030 bool thread
2031 ) : BbParaSolver(argc, argv, inNhanders, comm, inParaParamSet, inParaInstance, inDetTimer),
2032 messagehdlr(0),
2033 logfile(0),
2034 originalParamSet(0),
2035 conflictConsList(0),
2036 userPlugins(0),
2037 commPointHdlr(0),
2038 nodesel(0),
2039#if SCIP_APIVERSION >= 101
2040 selfSplitNodesel(0),
2041#endif
2042 scipPropagator(0),
2043 interruptMsgMonitor(0),
2044 originalPriority(0),
2045 orgMaxRestart(0),
2046 nOrgVars(0),
2047 nOrgVarsInSolvers(0),
2048 orgVarLbs(0),
2049 orgVarUbs(0),
2050 tightenedVarLbs(0),
2051 tightenedVarUbs(0),
2052 mapToOriginalIndecies(0),
2053 mapToSolverLocalIndecies(0),
2054 mapToProbIndecies(0),
2055 // stuffingMaxrounds(0),
2056 // domcolMaxrounds(0),
2057 // dualcompMaxrounds(0),
2058 // dualinferMaxrounds(0),
2059 // dualaggMaxrounds(0),
2060 // abspowerDualpresolve(0),
2061 // andDualpresolving(0),
2062 // cumulativeDualpresolve(0),
2063 // knapsackDualpresolving(0),
2064 // linearDualpresolving(0),
2065 // setppcDualpresolving(0),
2066 // logicorDualpresolving(0),
2067 miscAllowdualreds(0),
2068 nAddedConss(0),
2069 addedConss(0),
2070 addedDualCons(0),
2071 settingsNameLC(0),
2072 fiberSCIP(true),
2073 quiet(false),
2074 collectingModeIsProhibited(false),
2075 problemFileName(0),
2076 orgFeastol(0.0),
2077 orgLpfeastol(0.0),
2078 copyIncreasedVariables(false)
2079{
2080 assert(thread); // This is a constructor for threads parallel version
2081
2082 // ScipMessageHandlerFunctionPointer *scipMessageHandler = reinterpret_cast<ScipMessageHandlerFunctionPointer *>(messageHandler);
2083 // no additional message handlers
2084
2085 char* logname = NULL;
2086
2087 ScipParaInstance *scipParaInstance = dynamic_cast<ScipParaInstance *>(paraInstance);
2088
2089 paraTimer->setOffset(timeOffset);
2090
2091 /* Initialize the SCIP environment */
2092 /*********
2093 * Setup *
2094 *********/
2096 {
2097 /* copy SCIP environment */
2098 scip = scipParaInstance->getScip();
2099 SCIP_CALL_ABORT( SCIPresetParams(scip) ); // if LC parameter settings are applied,
2100 // it is necessary to reset them
2101 }
2102 else
2103 {
2104 SCIP_CALL_ABORT( SCIPcreate(&scip) );
2105 SCIP_CALL_ABORT( SCIPsetIntParam(scip,"timing/clocktype", 2) ); // always wall clock time
2107 {
2108 double timeRemains = std::max( 0.0, (paraParams->getRealParamValue(UG::TimeLimit) - paraTimer->getElapsedTime() + 3.0) ); // 3.0: timming issue
2109 SCIP_CALL_ABORT( SCIPsetRealParam(scip,"limits/time", timeRemains) );
2110 }
2111 SCIP_CALL_ABORT( SCIPincludeDefaultPlugins(scip) );
2112 /** user include plugins */
2114 }
2115
2116 /* include communication point handler */
2118 commPointHdlr = new ScipParaObjCommPointHdlr(paraComm, this, updator);
2119 SCIP_CALL_ABORT( SCIPincludeObjEventhdlr(scip, commPointHdlr, TRUE) );
2120 SCIP_CALL_ABORT( SCIPincludeObjHeur(scip, updator, TRUE) );
2121
2122 /* include propagator */
2126 {
2128 SCIP_CALL_ABORT( SCIPincludeObjProp(scip, scipPropagator, TRUE) );
2131 }
2132
2133 /* include node selector */
2134 nodesel = new ScipParaObjNodesel(this);
2135 SCIP_CALL_ABORT( SCIPincludeObjNodesel(scip, nodesel, TRUE) );
2136#if SCIP_APIVERSION >= 101
2137 if( paraParams->getIntParamValue(UG::RampUpPhaseProcess) == 3 ) // self-split ramp-up
2138 {
2139 selfSplitNodesel = new ScipParaObjSelfSplitNodesel(
2140 paraComm->getRank() - 1,
2141 paraComm->getSize() - 1,
2143 paraComm,
2144 this,
2145 scip
2146 );
2147 SCIP_CALL_ABORT( SCIPincludeObjNodesel(scip, selfSplitNodesel, TRUE) );
2148 }
2149#endif
2150 /* include branch rule plugins */
2151 SCIP_CALL_ABORT( SCIPincludeObjBranchrule(scip, new ScipParaObjBranchRule(this), TRUE) );
2152
2153 if( inParaParamSet->getBoolParamValue(UG::TransferConflictCuts) )
2154 {
2155 conflictConsList = new std::list<LocalNodeInfoPtr>;
2156 SCIP_CONFLICTHDLRDATA *conflictHdrData = reinterpret_cast< SCIP_CONFLICTHDLRDATA * >(this);
2157 /* create conflict handler to collects conflicts */
2158#if SCIP_VERSION == 211 && SCIP_SUBVERSION == 0
2159 SCIP_CALL_ABORT( SCIPincludeConflicthdlr(scip, CONFLICTHDLR_NAME, CONFLICTHDLR_DESC, CONFLICTHDLR_PRIORITY,
2160 NULL, NULL, NULL, NULL, NULL, NULL, conflictExecCollector, conflictHdrData) );
2161#else
2162 SCIP_CALL_ABORT( SCIPincludeConflicthdlrBasic(scip, NULL, CONFLICTHDLR_NAME, CONFLICTHDLR_DESC, CONFLICTHDLR_PRIORITY,
2163 conflictExecCollector, conflictHdrData) );
2164#endif
2165 }
2166
2167 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/feastol", &orgFeastol ) );
2168 if( SCIP_APIVERSION < 61 )
2169 {
2170 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/lpfeastol", &orgLpfeastol ) );
2171 }
2172
2173 /********************
2174 * Parse parameters *
2175 ********************/
2176 for( int i = 3; i < argc; ++i ) /** the first argument is runtime parameter file for ParaSCIP */
2177 {
2178 if( strcmp(argv[i], "-l") == 0 )
2179 {
2180 i++;
2181 if( i < argc )
2182 logname = argv[i];
2183 else
2184 {
2185 THROW_LOGICAL_ERROR1("missing log filename after parameter '-l'");
2186 }
2187 }
2188 else if( strcmp(argv[i], "-q") == 0 )
2189 quiet = true;
2190 // other arguments are omitted in Solver
2191 }
2192
2193 /***********************************
2194 * create log file message handler *
2195 ***********************************/
2197 {
2198#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
2199 SCIP_CALL_ABORT( SCIPcreateObjMessagehdlr(&messagehdlr, new ScipParaObjMessageHdlr(paraComm, NULL, TRUE, FALSE), TRUE) );
2200 SCIP_CALL_ABORT( SCIPsetMessagehdlr(messagehdlr) );
2201#else
2202 SCIPsetMessagehdlrQuiet(scip, TRUE );
2203#endif
2204 }
2205 else
2206 {
2207 if( logname != NULL || quiet )
2208 {
2209 paraComm->lockApp(); // if solver runs as thread, this lock is necessary
2210 if( logname != NULL )
2211 {
2212 std::ostringstream os;
2213 os << logname << comm->getRank();
2214 logfile = fopen(os.str().c_str(), "a"); // append to log file */
2215 if( logfile == NULL )
2216 {
2217 THROW_LOGICAL_ERROR3("cannot open log file <", logname, "> for writing");
2218 }
2219 }
2220 SCIP_CALL_ABORT( SCIPcreateObjMessagehdlr(&messagehdlr, new ScipParaObjMessageHdlr(paraComm, logfile, quiet, TRUE), TRUE) );
2221#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
2222 SCIP_CALL_ABORT( SCIPsetMessagehdlr(messagehdlr) );
2223#else
2224 SCIP_CALL_ABORT( SCIPsetMessagehdlr(scip, messagehdlr) );
2225 SCIP_CALL_ABORT( SCIPmessagehdlrRelease(&messagehdlr));
2226#endif
2228 }
2229 }
2230
2231 DEF_SCIP_PARA_COMM( scipParaComm, paraComm );
2232 scipDiffParamSetRoot = scipParaComm->createScipDiffParamSet();
2233 scipDiffParamSetRoot->bcast(comm, 0); /** this bcast is sent as SolverInitializationMessage */
2234 scipDiffParamSet = scipParaComm->createScipDiffParamSet();
2235 scipDiffParamSet->bcast(comm, 0); /** this bcast is sent as SolverInitializationMessage */
2236 int tempIsWarmStarted;
2237 comm->bcast(&tempIsWarmStarted, 1, UG::ParaINT, 0);
2238 warmStarted = (tempIsWarmStarted == 1);
2240
2242 {
2243 int solutionExists = 0;
2244 paraComm->bcast(&solutionExists, 1, UG::ParaINT, 0);
2245 }
2246 else
2247 {
2249 {
2250 int solutionExists = 0;
2251 comm->bcast(&solutionExists, 1, UG::ParaINT, 0);
2252 if( solutionExists )
2253 {
2256 }
2257 }
2258 }
2259
2260 /** set parameters for SCIP: this values are reseted before solving */
2261 /* move to the place after problem has been created because the parameters may be changed.
2262 scipDiffParamSet->setParametersInScip(scip);
2263 SCIP_Real epsilon;
2264 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon));
2265 eps = epsilon;
2266 */
2267
2268 char *isolname = 0;
2269 for( int i = 3; i < argc; ++i ) /** the first argument is runtime parameter file for ParaSCIP */
2270 /** the second argument is problem file name */
2271 {
2272 if( strcmp(argv[i], "-sl") == 0 )
2273 {
2274 i++;
2275 if( i < argc )
2276 {
2277 settingsNameLC = argv[i];
2278 break;
2279 }
2280 else
2281 {
2282 std::cerr << "missing settings filename after parameter '-sl'" << std::endl;
2283 exit(1);
2284 }
2285 }
2286 else if ( strcmp(argv[i], "-isol") == 0 )
2287 {
2288 i++;
2289 if( i < argc )
2290 {
2291 isolname = argv[i];
2292 }
2293 else
2294 {
2295 std::cerr << "missing settings filename after parameter '-isol'" << std::endl;
2296 exit(1);
2297 }
2298 }
2299 }
2300
2301#if( !defined(UG_QUBO) && !defined(UG_SMOOTHIE) )
2302
2303 /** create problem */
2304 scipParaInstance->createProblem(scip,
2311 isolname
2312 );
2313
2314 /** set parameters for SCIP: this values are reseted before solving */
2316 SCIP_Real epsilon;
2317 SCIP_CALL_ABORT( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon));
2318 eps = epsilon;
2319
2322 {
2323 /* initialize SCIP to check root solvability */
2324 SCIP_CALL_ABORT( SCIPcreate(&scipToCheckEffectOfRootNodeProcesses) );
2325 /* include default SCIP plugins */
2326 SCIP_CALL_ABORT( SCIPincludeDefaultPlugins(scipToCheckEffectOfRootNodeProcesses) );
2327 /* include scipParaConshdlr plugins */
2335 isolname
2336 );
2338 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "presolving/maxrestarts", 0) );
2339 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "presolving/maxrounds", 0) );
2340 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/linear/presolpairwise", FALSE) );
2341 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/and/presolpairwise", FALSE) );
2342 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/logicor/presolpairwise", FALSE) );
2343 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/setppc/presolpairwise", FALSE) );
2344 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "propagating/probing/maxprerounds", 0) );
2345 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "heuristics/feaspump/freq", -1) );
2346 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "heuristics/rens/freq", -1) );
2347 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "separating/maxcutsroot", 100) );
2348 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "separating/maxroundsroot", 5) );
2349 }
2350 delete paraInstance;
2351 paraInstance = 0;
2352
2353 // std::cout << "##################### paraInstance is deleted ####################" << std::endl;
2354
2355#if ( defined(_COMM_PTH) || defined(_COMM_CPP11) )
2356 assert( memoryLimitOfSolverSCIP > 0.0 );
2357 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "limits/memory", memoryLimitOfSolverSCIP) );
2358#else
2359 SCIP_CALL_ABORT( SCIPsetRealParam(scip, "limits/memory", (dynamic_cast<ScipParaParamSet *>(paraParams)->getRealParamValue(MemoryLimit))/(paraComm->getSize()*SCIP_MEMORY_COPY_FACTOR))); // LC has SCIP env.
2360#endif
2361// if( paraComm->getRank() == 1 )
2362// {
2363// std::cout << "*** set memory limit to " << dynamic_cast<ScipParaParamSet *>(paraParams)->getRealParamValue(MemoryLimit)/(paraComm->getSize()*SCIP_MEMORY_COPY_FACTOR) << " for each SCIP ***" << std::endl;
2364// }
2365
2366 /** save original priority of changing node selector */
2368
2369#ifndef SCIP_EVENTTYPE_COMM
2370 // std::cout << "#### SCIP_EVENTTYPE_COMM is not defined!" << std::endl;
2372 {
2373 interruptMsgMonitor = new ScipParaInterruptMsgMonitor(scipParaComm, this);
2374 // interruptMsgMonitorThread = std::thread( runInterruptMsgMonitorThread, interruptMsgMonitor );
2376 t.detach();
2377 }
2378#endif
2379
2380#endif
2381}
2382
2384 )
2385{
2387 {
2389 }
2390
2391 /** delete scip diff param sets */
2394 if( userPlugins ) delete userPlugins;
2395
2396 /** reset message handler */
2397 // message handler is mangaed within scip. It is freed at SCIPfree
2398#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
2399 if( messagehdlr )
2400 {
2401 SCIP_CALL_ABORT( SCIPsetDefaultMessagehdlr() );
2402 SCIP_CALL_ABORT( SCIPfreeObjMessagehdlr(&messagehdlr) );
2403 }
2404#endif
2405
2406 /* free SCIP environment */
2408 {
2409 SCIP_CALL_ABORT( SCIPfree(&scipToCheckEffectOfRootNodeProcesses) );
2410 }
2411 SCIP_CALL_ABORT( SCIPfree(&scip) );
2412
2413 /** close log file */
2414 if( logfile ) fclose(logfile);
2415
2416 if( conflictConsList && conflictConsList->size() > 0 )
2417 {
2418 int nConfilcts = conflictConsList->size();
2419 for(int i = 0; i < nConfilcts; i++ )
2420 {
2421 assert(!conflictConsList->empty());
2422 LocalNodeInfo *info= conflictConsList->front();
2423 conflictConsList->pop_front();
2424 if( info->linearCoefs ) delete[] info->linearCoefs;
2425 if( info->idxLinearCoefsVars ) delete[] info->idxLinearCoefsVars;
2426 delete info;
2427 }
2428 }
2429
2431
2432 if( orgVarLbs ) delete [] orgVarLbs;
2433 if( orgVarUbs ) delete [] orgVarUbs;
2434 if( tightenedVarLbs ) delete [] tightenedVarLbs;
2435 if( tightenedVarUbs ) delete [] tightenedVarUbs;
2438 if( mapToProbIndecies) delete [] mapToProbIndecies;
2439 if( addedConss ) delete [] addedConss;
2440
2441 // DON'T FREE problemFileName
2442
2443 ///
2444 /// The following code was in Solver base class before
2445 ///
2448 {
2450 }
2451
2453 {
2454 for(;;)
2455 {
2456 // if( !( paraParams->getRealParamValue(UG::TimeLimit) > 0.0 && paraParams->getRealParamValue(UG::TimeLimit) <= paraTimer->getElapsedTime() ) )
2457 // {
2458 while( !waitToken(paraComm->getRank()) )
2459 {
2461 // if( ( paraParams->getRealParamValue(UG::TimeLimit) > 0.0
2462 // && paraParams->getRealParamValue(UG::TimeLimit) <= paraTimer->getElapsedTime() ) )
2463 // {
2464 // break;
2465 // }
2466 }
2467 // paraDetTimer->update(1.0);
2468 // previousCommTime = paraDetTimer->getElapsedTime();
2469 if( paraComm->passTermToken(paraComm->getRank()) ) break;
2470 // }
2471 }
2472 }
2473
2474 double stopTime = paraTimer->getElapsedTime();
2478
2479 double detTime = -1.0;
2480 if( paraDetTimer )
2481 {
2482 detTime = paraDetTimer->getElapsedTime();
2483 }
2484
2485 DEF_BB_PARA_COMM(bbParaComm, paraComm);
2486
2487 UG::BbParaSolverTerminationState *paraSolverTerminationState = dynamic_cast<UG::BbParaSolverTerminationState *>(bbParaComm->createParaSolverTerminationState(
2489 paraComm->getRank(),
2491 minNSolved,
2492 maxNSolved,
2493 totalNSent,
2508 nTightened,
2510 calcTermState,
2511 stopTime,
2521 detTime ));
2522 paraSolverTerminationState->send(paraComm, 0, UG::TagTerminated);
2523
2524 delete paraSolverTerminationState;
2525
2526 // if( interruptMsgMonitor )
2527 // {
2528 // interruptMsgMonitorThread.join();
2529 // }
2530
2531}
2532
2533void
2535 const std::string& filename
2536 )
2537{
2538 FILE *file = fopen(filename.c_str(),"a");
2539 if( !file )
2540 {
2541 std::cout << "file : " << filename << "cannot open." << std::endl;
2542 abort();
2543 }
2544 SCIP_CALL_ABORT( SCIPprintTransProblem(scip, file, "cip", FALSE) );
2545}
2546
2547void
2549 )
2550{
2551 SCIP *backupScip = scip;
2554 /** set cutoff value */
2555 SCIP_CALL_ABORT( SCIPsetObjlimit(scip, globalBestIncumbentValue) );
2556 /** solve */
2557 SCIP_CALL_ABORT( SCIPsolve(scip) );
2558 nSolvedWithNoPreprocesses = SCIPgetNNodes(scip);
2560 scip = backupScip;
2561}
2562
2563void
2565 UG::ParaSolution *sol
2566 )
2567{
2568
2569 if( SCIPgetStage(scip) <= SCIP_STAGE_TRANSFORMING || SCIPgetStage(scip) >= SCIP_STAGE_SOLVED )
2570 {
2571 THROW_LOGICAL_ERROR1("invalid tyrNewSolution");
2572 }
2573
2574 ScipParaSolution *tempSol = dynamic_cast< ScipParaSolution * >(sol);
2575 SCIP_SOL* newsol; /* solution to be created for the original problem */
2576
2578 {
2579 return;
2580 //
2581 // It would be good to install the solution to SCIP as the PartialSolution.
2582 // However, currently SCIP does not support this.
2583 // In the future, this might be better to change.
2584 //
2585 // if( SCIPcreatePartialSol(scip, &newsol, 0) != SCIP_OKAY ) // SCIP_CALL_ABORT ????
2586 // {
2587 // return;
2588 // }
2589 }
2590 else
2591 {
2592 if( SCIPcreateOrigSol(scip, &newsol, 0) != SCIP_OKAY ) // SCIP_CALL_ABORT ????
2593 {
2594 return;
2595 }
2596 }
2597
2598
2599 SCIP_VAR** vars = SCIPgetOrigVars(scip);
2600 SCIP_Real* vals = new SCIP_Real[tempSol->getNVars()]();
2601 int i;
2602 for(i = 0; i < tempSol->getNVars(); i++ )
2603 {
2604 int probindex = tempSol->indexAmongSolvers(i);
2605 if( mapToProbIndecies )
2606 {
2607 probindex = mapToProbIndecies[tempSol->indexAmongSolvers(i)];
2608 }
2609 // assert( probindex >= 0 );
2610 if( probindex >= 0 )
2611 {
2612 vals[probindex] = tempSol->getValues()[i];
2613 }
2614 // /* skip inactive varibales */
2615 // if( probindex < 0 )
2616 // continue;
2617 // SCIP_CALL_ABORT( SCIPsetSolVal(scip, newsol, vars[probindex], tempSol->getValues()[i]) );
2618 }
2619 SCIP_CALL_ABORT( SCIPsetSolVals(scip, newsol, tempSol->getNVars(), vars, vals) );
2620 delete [] vals;
2621 // if( i != tempSol->getNVars() )
2622 // {
2623 // /** the given solution should be generated in original space,
2624 // * therefore the solution values cannot use for ParaSCIP
2625 // */
2626 // SCIP_CALL_ABORT( SCIPfreeSol(scip, &newsol) );
2627 // // delete tempSol; // this case, DO NOT DELETE tempSol.
2628 // return;
2629 // }
2630
2631#if SCIP_VERSION == 211 && SCIP_SUBVERSION == 0
2632 if( SCIPgetStage(scip) == SCIP_STAGE_TRANSFORMED || SCIPgetStage(scip) == SCIP_STAGE_PRESOLVED )
2633#else
2634 if( // SCIPgetStage(scip) == SCIP_STAGE_TRANSFORMED ||
2635 SCIPgetStage(scip) == SCIP_STAGE_PRESOLVED ||
2636 SCIPgetStage(scip) == SCIP_STAGE_INITPRESOLVE )
2637#endif
2638 {
2639 SCIP_Bool success;
2640#if (SCIP_VERSION < 321 || ( SCIP_VERSION == 321 && SCIP_SUBVERSION < 2) )
2641 SCIP_CALL_ABORT( SCIPtrySolFree(scip, &newsol, FALSE, TRUE, TRUE, TRUE, &success) );
2642#else
2643 SCIP_CALL_ABORT( SCIPtrySolFree(scip, &newsol, FALSE, TRUE, TRUE, TRUE, TRUE, &success) );
2644#endif
2645 // std::cout << "Rank " << paraComm->getRank() << ": success = " << success << std::endl;
2646 // if( !success ) abort();
2647 }
2648 else
2649 {
2650 SCIP_CALL_ABORT( SCIPheurPassSolTrySol(scip, SCIPfindHeur(scip, "trysol"), newsol) );
2651 SCIP_CALL_ABORT( SCIPfreeSol(scip, &newsol) );
2652 }
2653}
2654
2655void
2657 )
2658{
2660 if( !originalParamSet )
2661 {
2662 DEF_SCIP_PARA_COMM( scipParaComm, paraComm );
2663 originalParamSet = scipParaComm->createScipDiffParamSet(scip);;
2664 }
2665 SCIP_CALL_ABORT( SCIPsetHeuristics(scip, SCIP_PARAMSETTING_FAST, TRUE) );
2666 SCIP_CALL_ABORT( SCIPsetPresolving(scip, SCIP_PARAMSETTING_FAST, TRUE) );
2667 SCIP_CALL_ABORT( SCIPsetSeparating(scip, SCIP_PARAMSETTING_FAST, TRUE) );
2668}
2669
2670void
2672 )
2673{
2674 assert(originalParamSet);
2677
2678#if SCIP_VERSION >= 320
2680#endif
2681
2682}
2683
2684// static int id = 0;
2685
2686void
2688 )
2689{
2690 std::string subcipprefix("SolverCip");
2691 std::string subcipfilename;
2692 std::ostringstream oss;
2693 oss << subcipprefix << paraComm->getRank();
2694 // oss << subcipprefix << paraComm->getRank() << "." << id++;;
2695 subcipfilename = oss.str();
2696 subcipfilename += ".cip";
2697 SCIP_CALL_ABORT( SCIPwriteOrigProblem(scip, subcipfilename.c_str(), "cip", FALSE) );
2698#ifdef UG_DEBUG_SOLUTION
2699 if( ( !currentTask->getDiffSubproblem() ) ||
2700 currentTask->getDiffSubproblem()->isOptimalSolIncluded() )
2701 {
2702 std::cout << "** " << subcipfilename << " contains optimal solution." << std::endl;
2703 }
2704 else
2705 {
2706 std::cout << "** " << subcipfilename << " does NOT contain optimal solution." << std::endl;
2707 }
2708#endif
2709 subcipfilename = oss.str();
2710 subcipfilename += "-t.cip";
2711 if( SCIPgetStage(scip) >= SCIP_STAGE_TRANSFORMED )
2712 {
2713 SCIP_CALL_ABORT( SCIPwriteTransProblem(scip, subcipfilename.c_str(), "cip", FALSE) );
2714 }
2715 char name[SCIP_MAXSTRLEN];
2716 (void)SCIPsnprintf(name, SCIP_MAXSTRLEN, "SolverCip%d.set", paraComm->getRank());
2717 SCIP_CALL_ABORT( SCIPwriteParams(scip, name, TRUE, FALSE) );
2718}
2719
2720void
2722 )
2723{
2724 assert(paraInstance);
2725 ScipParaInstance *scipParaInstance = dynamic_cast<ScipParaInstance*>(paraInstance);
2726 nOrgVars = scipParaInstance->getNVars(); // number of original variables in LC
2727 nOrgVarsInSolvers = SCIPgetNOrigVars(scip); // maybe increaded
2728 assert( nOrgVarsInSolvers == scipParaInstance->getVarIndexRange() );
2729 // if( nOrgVars == 0 ) nOrgVars = nOrgVarsInSolver;
2730 // nOrgVars = nOrgVarsInSolver;
2731 // assert( nOrgVars <= paraInstance->getVarIndexRange() ); // when copy generated additional variables, this does not hold
2732 orgVarLbs = new SCIP_Real[nOrgVarsInSolvers];
2733 orgVarUbs = new SCIP_Real[nOrgVarsInSolvers];
2734 if( scipParaInstance->isOriginalIndeciesMap() )
2735 {
2736 SCIP_VAR **vars = SCIPgetVars(scip);
2740 {
2741 tightenedVarLbs = new double[nOrgVarsInSolvers];
2742 tightenedVarUbs = new double[nOrgVarsInSolvers];
2743 // for(int v = 0; v < paraInstance->getVarIndexRange(); v++)
2744 // for( int v = 0; v < nOrgVarsInSolvers ; v++ )
2745 for( int v = 0; v < nOrgVars ; v++ )
2746 {
2747 // int orgIndex = scipParaInstance->getOrigProbIndex(SCIPvarGetIndex(vars[v]));
2748 // assert(orgIndex >= 0);
2749 tightenedVarLbs[v] = orgVarLbs[v] = scipParaInstance->getVarLb(v);
2750 tightenedVarUbs[v] = orgVarUbs[v] = scipParaInstance->getVarUb(v);
2751 // std::cout << scipParaInstance->getVarName(v) << " orgVarLbs[" << v << "] = " << orgVarLbs[v] << std::endl;
2752 // std::cout << scipParaInstance->getVarName(v) << " orgVarUbs[" << v << "] = " << orgVarUbs[v] << std::endl;
2753 }
2754 }
2755 else
2756 {
2757 // for(int v = 0; v < paraInstance->getVarIndexRange(); v++)
2758 // for( int v = 0; v < nOrgVarsInSolvers ; v++ )
2759 for( int v = 0; v < nOrgVars; v++ )
2760 {
2761 // int orgIndex = scipParaInstance->getOrigProbIndex(SCIPvarGetIndex(vars[v]));
2762 // assert(orgIndex >= 0);
2763 orgVarLbs[v] = scipParaInstance->getVarLb(v);
2764 orgVarUbs[v] = scipParaInstance->getVarUb(v);
2765 // std::cout << scipParaInstance->getVarName(v) << " orgVarLbs[" << v << "] = " << orgVarLbs[v] << std::endl;
2766 // std::cout << scipParaInstance->getVarName(v) << " orgVarUbs[" << v << "] = " << orgVarUbs[v] << std::endl;
2767 }
2768 }
2769 assert( mapToOriginalIndecies == 0 && mapToProbIndecies == 0 );
2770 mapToOriginalIndecies = scipParaInstance->extractOrigProbIndexMap();
2772 mapToProbIndecies = new int[SCIPgetNTotalVars(scip)]; // need to allocate enough for SCIPvarGetIndex(scip)
2773 for( int i = 0; i < SCIPgetNTotalVars(scip); i++ )
2774 {
2775 mapToProbIndecies[i] = -1;
2776 }
2777 for( int i = 0; i < nOrgVarsInSolvers; i++ )
2778 {
2779 mapToProbIndecies[SCIPvarGetIndex(vars[i])] = i;
2780 }
2781 if ( dynamic_cast<ScipParaInstance *>(paraInstance)->isCopyIncreasedVariables() )
2782 {
2784 }
2785 }
2786 else
2787 {
2791 {
2792 tightenedVarLbs = new double[nOrgVars];
2793 tightenedVarUbs = new double[nOrgVars];
2794 // for(int v = 0; v < paraInstance->getVarIndexRange(); v++)
2795 // for( int v = 0; v < nOrgVarsInSolvers; v++ )
2796 for( int v = 0; v < nOrgVars; v++ )
2797 {
2798 tightenedVarLbs[v] = orgVarLbs[v] = scipParaInstance->getVarLb(v);
2799 tightenedVarUbs[v] = orgVarUbs[v] = scipParaInstance->getVarUb(v);
2800 // std::cout << scipParaInstance->getVarName(v) << " orgVarLbs[" << v << "] = " << orgVarLbs[v] << std::endl;
2801 // std::cout << scipParaInstance->getVarName(v) << " orgVarUbs[" << v << "] = " << orgVarUbs[v] << std::endl;
2802 }
2803 }
2804 else
2805 {
2806 // for(int v = 0; v < paraInstance->getVarIndexRange(); v++)
2807 // for( int v = 0; v < nOrgVarsInSolvers; v++ )
2808 for( int v = 0; v < nOrgVars; v++ )
2809 {
2810 orgVarLbs[v] = scipParaInstance->getVarLb(v);
2811 orgVarUbs[v] = scipParaInstance->getVarUb(v);
2812 // std::cout << scipParaInstance->getVarName(v) << " orgVarLbs[" << v << "] = " << orgVarLbs[v] << std::endl;
2813 // std::cout << scipParaInstance->getVarName(v) << " orgVarUbs[" << v << "] = " << orgVarUbs[v] << std::endl;
2814 }
2815 }
2816 }
2817}
2818
2819void
2821 )
2822{
2823 /*************************************
2824 ** This function does not work well **
2825 **************************************/
2826 /****************************
2827 ** reset original instance *
2828 *****************************/
2829 /* Reinitialize the SCIP environment */
2830 /* free SCIP environment */
2832 {
2833 SCIP_CALL_ABORT( SCIPfree(&scipToCheckEffectOfRootNodeProcesses) );
2834 }
2835 SCIP_CALL_ABORT( SCIPfree(&scip) );
2836 if( orgVarLbs ) delete [] orgVarLbs;
2837 if( orgVarUbs ) delete [] orgVarUbs;
2838 /*********
2839 * Setup *
2840 *********/
2841 if( fiberSCIP )
2842 {
2844 // setUserPlugins(paraInstance); // instance data should not be read from original data file
2846 /* copy SCIP environment */
2847 ScipParaInstance *scipParaInstance = dynamic_cast<ScipParaInstance*>(paraInstance);
2848 scip = scipParaInstance->getScip();
2849 SCIP_CALL_ABORT( SCIPresetParams(scip) ); // if LC parameter settings are applied,
2850 // it is necessary to reset them
2851 }
2852 else
2853 {
2856 ScipParaInstance *scipParaInstance = dynamic_cast<ScipParaInstance*>(paraInstance);
2857 scipParaInstance->setFileName(problemFileName);
2859 /* initialize SCIP */
2860 SCIP_CALL_ABORT( SCIPcreate(&scip) );
2861 SCIP_CALL_ABORT( SCIPsetIntParam(scip,"timing/clocktype", 2) ); // always wall clock time
2863 {
2864 double timeRemains = std::max( 0.0, (paraParams->getRealParamValue(UG::TimeLimit) - paraTimer->getElapsedTime() + 3.0) ); // 3.0: timming issue
2865 SCIP_CALL_ABORT( SCIPsetRealParam(scip,"limits/time", timeRemains) );
2866 }
2867 /* include default SCIP plugins */
2868 SCIP_CALL_ABORT( SCIPincludeDefaultPlugins(scip) );
2869 /** user include plugins */
2871 }
2872
2873 /* include communication point handler */
2875 commPointHdlr = new ScipParaObjCommPointHdlr(paraComm, this, updator);
2876 SCIP_CALL_ABORT( SCIPincludeObjEventhdlr(scip, commPointHdlr, TRUE) );
2877 SCIP_CALL_ABORT( SCIPincludeObjHeur(scip, updator, TRUE) );
2878
2879 /* include propagator */
2883 {
2885 SCIP_CALL_ABORT( SCIPincludeObjProp(scip, scipPropagator, TRUE) );
2888 }
2889
2890 /* include node selector */
2891 nodesel = new ScipParaObjNodesel(this);
2892 SCIP_CALL_ABORT( SCIPincludeObjNodesel(scip, nodesel, TRUE) );
2893 /* include branch rule plugins */
2894 SCIP_CALL_ABORT( SCIPincludeObjBranchrule(scip, new ScipParaObjBranchRule(this), TRUE) );
2895
2897 {
2898 assert(conflictConsList);
2899 delete conflictConsList;
2900 conflictConsList = new std::list<LocalNodeInfoPtr>;
2901 SCIP_CONFLICTHDLRDATA *conflictHdrData = reinterpret_cast< SCIP_CONFLICTHDLRDATA * >(this);
2902 /* create conflict handler to collects conflicts */
2903#if SCIP_VERSION == 211 && SCIP_SUBVERSION == 0
2904 SCIP_CALL_ABORT( SCIPincludeConflicthdlr(scip, CONFLICTHDLR_NAME, CONFLICTHDLR_DESC, CONFLICTHDLR_PRIORITY,
2905 NULL, NULL, NULL, NULL, NULL, NULL, conflictExecCollector, conflictHdrData) );
2906#else
2907 SCIP_CALL_ABORT( SCIPincludeConflicthdlrBasic(scip, NULL, CONFLICTHDLR_NAME, CONFLICTHDLR_DESC, CONFLICTHDLR_PRIORITY,
2908 conflictExecCollector, conflictHdrData) );
2909#endif
2910 }
2911
2913 {
2914 // SCIP_CALL_ABORT( SCIPsetMessagehdlr(NULL) );
2915#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
2916 SCIP_CALL_ABORT( SCIPcreateObjMessagehdlr(&messagehdlr, new ScipParaObjMessageHdlr(paraComm, NULL, TRUE, FALSE), TRUE) );
2917 SCIP_CALL_ABORT( SCIPsetMessagehdlr(messagehdlr) );
2918#else
2919 SCIPsetMessagehdlrQuiet(scip, TRUE);
2920#endif
2921 }
2922 else
2923 {
2924 if( logfile != NULL || quiet )
2925 {
2926 SCIP_CALL_ABORT( SCIPcreateObjMessagehdlr(&messagehdlr, new ScipParaObjMessageHdlr(paraComm, logfile, quiet, FALSE), TRUE) );
2927#ifndef SCIP_THREADSAFE_MESSAGEHDLRS
2928 SCIP_CALL_ABORT( SCIPsetMessagehdlr(messagehdlr) );
2929#else
2930 SCIP_CALL_ABORT( SCIPsetMessagehdlr(scip, messagehdlr) );
2931 SCIP_CALL_ABORT( SCIPmessagehdlrRelease(&messagehdlr));
2932#endif
2933 }
2934 }
2935
2936 ScipParaInstance *scipParaInstance = dynamic_cast<ScipParaInstance*>(paraInstance);
2937 scipParaInstance->createProblem(scip,
2944 NULL
2945 );
2948 {
2956 NULL
2957 );
2959 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "presolving/maxrestarts", 0) );
2960 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "presolving/maxrounds", 0) );
2961 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/linear/presolpairwise", FALSE) );
2962 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/and/presolpairwise", FALSE) );
2963 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/logicor/presolpairwise", FALSE) );
2964 SCIP_CALL_ABORT( SCIPsetBoolParam(scipToCheckEffectOfRootNodeProcesses, "constraints/setppc/presolpairwise", FALSE) );
2965 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "propagating/probing/maxprerounds", 0) );
2966 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "heuristics/feaspump/freq", -1) );
2967 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "heuristics/rens/freq", -1) );
2968 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "separating/maxcutsroot", 100) );
2969 SCIP_CALL_ABORT( SCIPsetIntParam(scipToCheckEffectOfRootNodeProcesses, "separating/maxroundsroot", 5) );
2970 }
2971 delete paraInstance;
2972 paraInstance = 0;
2973}
2974
2975void
2977 )
2978{
2979 if( SCIPgetStage(scip) != SCIP_STAGE_SOLVING )
2980 {
2982 {
2983 // setWinnerRacingParams(winnerRacingParams); // winner parameters are set, again
2984 // std::cout << winnerRacingParams->toString() << std::endl;
2985 }
2986 else
2987 {
2989 }
2990 }
2992
2993#if SCIP_VERSION >= 320
2994 if( currentTask )
2995 {
2997 }
2998#endif
2999}
3000
3001void
3003 )
3004{
3005 char *bakFileName = NULL;
3006 SCIP_CALL_ABORT( SCIPgetStringParam(scip,"visual/bakfilename", &bakFileName) );
3007 if( strcmp(bakFileName,"-") != 0 )
3008 {
3009 std::ostringstream os;
3010 os << bakFileName << "_" << paraComm->getRank();
3011 SCIP_CALL_ABORT( SCIPsetStringParam(scip,"visual/bakfilename", os.str().c_str() ) );
3012 SCIP_CALL_ABORT( SCIPsetStringParam(scip,"visual/baknodeprefix", ((currentTask->getTaskId()).toString()+":").c_str() ) );
3013 SCIP_CALL_ABORT( SCIPsetStringParam(scip,"visual/bakrootinfo", (currentTask->getGeneratorTaskId()).toString().c_str() ) );
3014 SCIP_CALL_ABORT( SCIPsetRealParam(scip,"visual/baktimeoffset", paraTimer->getElapsedTime() ) );
3015 }
3016}
3017
3018int
3020 int source,
3021 int tag
3022 )
3023{
3024 int tightenedIdex;
3025 double tightenedBound;
3027 paraComm->receive( (void *)&tightenedIdex, 1, UG::ParaINT, source, UG::TagLbBoundTightenedIndex )
3028 );
3030 paraComm->receive( (void *)&tightenedBound, 1, UG::ParaDOUBLE, source, UG::TagLbBoundTightenedBound )
3031 );
3032
3034 {
3035 return 0;
3036 }
3037
3038 if( mapToProbIndecies )
3039 {
3040 assert( mapToSolverLocalIndecies );
3041 assert(mapToProbIndecies[mapToSolverLocalIndecies[tightenedIdex]] >= 0);
3042 tightenedIdex = mapToProbIndecies[mapToSolverLocalIndecies[tightenedIdex]];
3043 }
3044 assert( SCIPisLE(scip,orgVarLbs[tightenedIdex], tightenedBound) &&
3045 SCIPisGE(scip,orgVarUbs[tightenedIdex], tightenedBound) );
3046 // std::cout << "Rank " << paraComm->getRank() << ": receive tightened lower bond. idx = " << tightenedIdex << ", bound = " << tightenedBound << std::endl;
3047 SCIP_Var* var = SCIPvarGetTransVar(SCIPgetOrigVars(scip)[tightenedIdex]);
3048 if( var && SCIPisLT(scip,tightenedVarLbs[tightenedIdex], tightenedBound) && SCIPvarGetStatus(var) != SCIP_VARSTATUS_MULTAGGR )
3049 {
3050 // std::cout << "Solver Lb = " << tightenedVarLbs[tightenedIdex]
3051 // << ", Rank " << paraComm->getRank() << ": try to tighten lower bond. idx = " << tightenedIdex << ", bound = " << tightenedBound << std::endl;
3052 scipPropagator->addBoundChange(scip, SCIP_BOUNDTYPE_LOWER, tightenedIdex, tightenedBound );
3053 tightenedVarLbs[tightenedIdex] = tightenedBound;
3054
3055 }
3056
3057 return 0;
3058}
3059
3060int
3062 int source,
3063 int tag
3064 )
3065{
3066 int tightenedIdex;
3067 double tightenedBound;
3069 paraComm->receive( (void *)&tightenedIdex, 1, UG::ParaINT, source, UG::TagUbBoundTightenedIndex )
3070 );
3072 paraComm->receive( (void *)&tightenedBound, 1, UG::ParaDOUBLE, source, UG::TagUbBoundTightenedBound )
3073 );
3074
3076 {
3077 return 0;
3078 }
3079
3080 if( mapToProbIndecies )
3081 {
3082 assert( mapToSolverLocalIndecies );
3083 assert(mapToProbIndecies[mapToSolverLocalIndecies[tightenedIdex]] >= 0);
3084 tightenedIdex = mapToProbIndecies[mapToSolverLocalIndecies[tightenedIdex]];
3085 }
3086 assert( SCIPisLE(scip,orgVarLbs[tightenedIdex], tightenedBound) &&
3087 SCIPisGE(scip,orgVarUbs[tightenedIdex], tightenedBound) );
3088 // std::cout << "Rank " << paraComm->getRank() << ": receive tightened upper bond. idx = " << tightenedIdex << ", bound = " << tightenedBound << std::endl;
3089 SCIP_Var* var = SCIPvarGetTransVar(SCIPgetOrigVars(scip)[tightenedIdex]);
3090 if( var && SCIPisGT(scip,tightenedVarUbs[tightenedIdex], tightenedBound) && SCIPvarGetStatus(var) != SCIP_VARSTATUS_MULTAGGR )
3091 {
3092 // std::cout << "Solver Ub = " << tightenedVarUbs[tightenedIdex]
3093 // << ", Rank " << paraComm->getRank() << ": try to tighten upper bond. idx = " << tightenedIdex << ", bound = " << tightenedBound << std::endl;
3094 scipPropagator->addBoundChange(scip, SCIP_BOUNDTYPE_UPPER, tightenedIdex, tightenedBound );
3095 tightenedVarUbs[tightenedIdex] = tightenedBound;
3096 }
3097
3098 return 0;
3099}
3100
3101/** get number of tightened variables during racing */
3102int
3104 )
3105{
3106 if( scipPropagator )
3107 {
3108 return scipPropagator->getNtightened();
3109 }
3110 else
3111 {
3112 if( nTightened > 0 )
3113 {
3114 return nTightened;
3115 }
3116 else
3117 {
3118 return 0;
3119 }
3120 }
3121}
3122
3123/** get number of tightened integral variables during racing */
3124int
3126 )
3127{
3128 if( scipPropagator )
3129 {
3131 }
3132 else
3133 {
3134 if( nTightened > 0 )
3135 {
3136 return nTightenedInt;
3137 }
3138 else
3139 {
3140 return 0;
3141 }
3142 }
3143}
3144
3145void
3147 )
3148{
3149 /* presolvers */
3150#if 0
3151 SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/stuffing/maxrounds", &stuffingMaxrounds) );
3152 SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/domcol/maxrounds", &domcolMaxrounds) );
3153#if ( SCIP_VERSION >= 322 || (SCIP_VERSION == 321 && SCIP_SUBVERSION >= 2) )
3154 SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/dualcomp/maxrounds", &dualcompMaxrounds) ); /*TODO: ok? */
3155#endif
3156 SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/dualinfer/maxrounds", &dualinferMaxrounds) ); /*TODO: probably fine */
3157 // SCIP_CALL_ABORT( SCIPgetIntParam(scip, "presolving/dualagg/maxrounds", &dualaggMaxrounds ) ); // TODO: seems to have no copy callback */
3158 /* constraint handlers */
3159 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/abspower/dualpresolve", &abspowerDualpresolve) );
3160 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/and/dualpresolving", &andDualpresolving) );
3161 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/cumulative/dualpresolve", &cumulativeDualpresolve) );
3162 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/knapsack/dualpresolving", &knapsackDualpresolving) );
3163 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/linear/dualpresolving", &linearDualpresolving) );
3164 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/setppc/dualpresolving", &setppcDualpresolving) );
3165 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "constraints/logicor/dualpresolving", &logicorDualpresolving) );
3166#endif
3167 if ( isRacingStage() && paraComm->getRank() != 1 )
3168 {
3169#if SCIP_APIVERSION > 34
3170 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "misc/allowstrongdualreds", &miscAllowdualreds) );
3171#else
3172 SCIP_CALL_ABORT( SCIPgetBoolParam(scip, "misc/allowdualreds", &miscAllowdualreds) );
3173#endif
3174 }
3175}
3176
3177void
3179 )
3180{
3181 /* presolvers */
3182#if 0
3183 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/stuffing/maxrounds", 0) );
3184 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/domcol/maxrounds", 0) );
3185#if ( SCIP_VERSION >= 322 || (SCIP_VERSION == 321 && SCIP_SUBVERSION >= 2) )
3186 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/dualcomp/maxrounds", 0) ); /*TODO: ok? */
3187#endif
3188 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/dualinfer/maxrounds", 0) ); /*TODO: probably fine */
3189 // SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/dualagg/maxrounds", 0) ); // TODO: seems to have no copy callback */
3190 /* constraint handlers */
3191 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/abspower/dualpresolve", FALSE) );
3192 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/and/dualpresolving", FALSE) );
3193 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/cumulative/dualpresolve", FALSE) );
3194 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/knapsack/dualpresolving", FALSE) );
3195 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/linear/dualpresolving", FALSE) );
3196 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/setppc/dualpresolving", FALSE) );
3197 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/logicor/dualpresolving", FALSE) );
3198#endif
3199 if ( isRacingStage() && paraComm->getRank() != 1 )
3200 {
3201#if SCIP_APIVERSION > 34
3202 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "misc/allowstrongdualreds", FALSE) );
3203#else
3204 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "misc/allowdualreds", FALSE) );
3205#endif
3206 }
3207}
3208
3209void
3211 )
3212{
3213 /* presolvers */
3214#if 0
3215 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/stuffing/maxrounds", stuffingMaxrounds) );
3216 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/domcol/maxrounds", domcolMaxrounds) );
3217#if ( SCIP_VERSION >= 322 || (SCIP_VERSION == 321 && SCIP_SUBVERSION >= 2) )
3218 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/dualcomp/maxrounds", dualcompMaxrounds) ); /*TODO: ok? */
3219#endif
3220 SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/dualinfer/maxrounds", dualinferMaxrounds) ); /*TODO: probably fine */
3221 // SCIP_CALL_ABORT( SCIPsetIntParam(scip, "presolving/dualagg/maxrounds", dualaggMaxrounds ) ); // TODO: seems to have no copy callback */
3222 /* constraint handlers */
3223 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/setppc/dualpresolving", setppcDualpresolving) );
3224 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "constraints/logicor/dualpresolving", logicorDualpresolving) );
3225#endif
3227 paraParams->getIntParamValue(UG::RampUpPhaseProcess) == 2 ) // may not be in racing statge
3228 && paraComm->getRank() != 1 )
3229 {
3230#if SCIP_APIVERSION > 34
3231 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "misc/allowstrongdualreds", miscAllowdualreds) );
3232#else
3233 SCIP_CALL_ABORT( SCIPsetBoolParam(scip, "misc/allowdualreds", miscAllowdualreds) );
3234#endif
3235 }
3236}
3237
3238void
3240{
3241 SCIP_CALL_ABORT( SCIPinterruptSolve(scip) );
3243}
3244
3245bool
3247{
3248 return commPointHdlr->isInterrupting();
3249}
#define DEF_BB_PARA_COMM(para_comm, comm)
Base class of communicator for UG Framework.
Base class for instance data.
Base class for BbParaNode.
Base class for solution.
This class contains solver termination state which is transferred form Solver to LC.
Base class for solver: Generic parallelized solver.
void setParametersInScip(SCIP *scip)
virtual int bcast(UG::ParaComm *comm, int root)=0
SCIP_Real getBoundsBoundDisjunction(int i, int j)
ScipParaDiffSubproblemBranchLinearCons * getBranchLinearConss()
SCIP_Real getVarValueUpinfer(int i, int j)
SCIP_Real getBendersLinearCoefs(int i, int j)
SCIP_Real getVarValueDownconflen(int i, int j)
SCIP_Real getVarValueDownvsids(int i, int j)
SCIP_BOUNDTYPE getBoundTypesBoundDisjunction(int i, int j)
SCIP_Real getVarValueUpcutoff(int i, int j)
SCIP_Real getBranchConsLinearCoefs(int i, int j)
SCIP_Real getVarValueDowninfer(int i, int j)
SCIP_Real getVarValueUpconflen(int i, int j)
ScipParaDiffSubproblemBranchSetppcCons * getBranchSetppcConss()
SCIP_Real getVarValueDowncutoff(int i, int j)
SCIP_Real getVarVlaueUpvsids(int i, int j)
virtual void setFileName(const char *fileName)=0
void createProblem(SCIP *scip, int method, bool noPreprocessingInLC, bool usetRootNodeCuts, ScipDiffParamSet *scipDiffParamSetRoot, ScipDiffParamSet *scipDiffParamSet, char *settingsNameLC, char *isolname)
void addBoundChange(SCIP *scip, SCIP_BOUNDTYPE boundType, int index, SCIP_Real bound)
class BbParaParamSet
void tryNewSolution(UG::ParaSolution *sol)
int lbBoundTightened(int source, int tag)
virtual void createSubproblem()
int getOriginalIndex(int index)
ScipParaSolver(int argc, char **argv, UG::ParaComm *comm, UG::ParaParamSet *paraParamSet, UG::ParaInstance *inParaInstance, UG::ParaDeterministicTimer *detTimer)
ScipDiffParamSet * originalParamSet
SCIP_MESSAGEHDLR * messagehdlr
std::list< LocalNodeInfoPtr > * conflictConsList
void setRacingParams(UG::ParaRacingRampUpParamSet *inRacingParams, bool winnerParam)
ScipDiffParamSet * scipDiffParamSet
ScipParaObjNodesel * nodesel
ScipParaInterruptMsgMonitor * interruptMsgMonitor
interrupt message monitor
ScipUserPlugins * userPlugins
ScipParaObjCommPointHdlr * commPointHdlr
void setWinnerRacingParams(UG::ParaRacingRampUpParamSet *inRacingParams)
void setUserPlugins(ScipUserPlugins *inUi)
unsigned int miscAllowdualreds
ScipParaObjProp * scipPropagator
void writeCurrentTaskProblem(const std::string &filename)
int ubBoundTightened(int source, int tag)
static void runInterruptMsgMonitorThread(void *threadData)
std::list< LocalNodeInfoPtr > * getConflictConsList()
SCIP * scipToCheckEffectOfRootNodeProcesses
ScipDiffParamSet * scipDiffParamSetRoot
void includeUserPlugins(SCIP *inScip)
virtual void writeSubproblem(SCIP *scip)
virtual void newSubproblem(SCIP *scip, const ScipParaDiffSubproblemBranchLinearCons *linearConss, const ScipParaDiffSubproblemBranchSetppcCons *setppcConss)
class BbParaNode
Definition: bbParaNode.h:62
double getInitialDualBoundValue()
getter of initial dual bound value
Definition: bbParaNode.h:314
double getDualBoundValue()
getter of dual bound value
Definition: bbParaNode.h:303
class BbParaSolverTerminationState (Solver termination state in a ParaSolver)
int minNSolved
minimum number of subtree nodes rooted from ParaNode
Definition: bbParaSolver.h:103
bool restartingRacing
indicate that this solver is restarting racing
Definition: bbParaSolver.h:146
int maxNSolved
maximum number of subtree nodes rooted from ParaNode
Definition: bbParaSolver.h:104
bool isRacingInterruptRequested()
check if racing interrupt was requested or not
Definition: bbParaSolver.h:834
int nTransferredLocalCutsFromSolver
number of local cuts transferred from this Solver
Definition: bbParaSolver.h:105
int nTightened
the number of tightened variable bounds in racing
Definition: bbParaSolver.h:126
bool isCollectingAllNodes()
check if Solver is sending all nodes to LoadCoordinaor or not
int maxRestarts
maximum number of restarts
Definition: bbParaSolver.h:113
ParaTask * getCurrentNode()
get current ParaNode object
Definition: bbParaSolver.h:981
int totalNSolved
Counters related to this BbParaSolver.
Definition: bbParaSolver.h:102
ParaParamSet * getParaParamSet()
get ParaParamSet object
int nParaNodesSolvedAtRoot
number of ParaNodes solved at root node
Definition: bbParaSolver.h:116
bool lightWeightRootNodeComputation
indicate that fast root node computation is required
Definition: bbParaSolver.h:80
virtual void waitNotificationIdMessage()
wait notification id message to synchronized with LoadCoordinator
double maxRootNodeTime
maximum time consumed by root node process
Definition: bbParaSolver.h:90
virtual void sendLocalSolution()
send solution found in this Solver
double getGlobalBestIncumbentValue()
get global best incumbent value
Definition: bbParaSolver.h:971
virtual void iReceiveMessages()
non-blocking receive messages
int nSolvedWithNoPreprocesses
number of nodes solved when it is solved with no preprocesses
Definition: bbParaSolver.h:97
int totalNSent
accumulated number of nodes sent from this BbParaSolver
Definition: bbParaSolver.h:114
int minRestarts
minimum number of restarts
Definition: bbParaSolver.h:112
double solverDualBound
dual bound value achieved for a subproblem
Definition: bbParaSolver.h:137
int nTightenedInt
the number of tightened integral variable bounds in racing
Definition: bbParaSolver.h:127
int nTransferredBendersCutsFromSolver
number of benders cuts transferred from this Solver
Definition: bbParaSolver.h:108
bool isRacingStage()
check if Solver is in racing stage or not
double minRootNodeTime
minimum time consumed by root node process
Definition: bbParaSolver.h:89
int minTransferredLocalCutsFromSolver
minimum number of local cuts transferred from this Solver
Definition: bbParaSolver.h:106
int minTransferredBendersCutsFromSolver
minimum number of benders cuts transferred from this Solver
Definition: bbParaSolver.h:109
int nParaNodesSolvedAtPreCheck
number of ParaNodes solved at pre-checking of root node solvability
Definition: bbParaSolver.h:117
int maxTransferredBendersCutsFromSolver
maximum number of benders cuts transferred from this Solver
Definition: bbParaSolver.h:110
virtual bool saveIfImprovedSolutionWasFound(ParaSolution *sol)
save improved solution if it was found in this Solver
int maxTransferredLocalCutsFromSolver
maximum number of local cuts transferred from this Solver
Definition: bbParaSolver.h:107
double totalRootNodeTime
accumulated root node process time solved by this solver so far
Definition: bbParaSolver.h:88
int totalNImprovedIncumbent
accumulated number of improvements of incumbent value in this BbParaSolver
Definition: bbParaSolver.h:115
int nTotalRestarts
number of total restarts
Definition: bbParaSolver.h:111
bool waitToken(int rank)
wait token for deterministic mode
int getRank()
get rank of caller's thread
int bcast(void *buffer, int count, const int datatypeId, int root)
broadcast function for standard ParaData types
Base class of communicator object.
Definition: paraComm.h:102
virtual void lockApp()=0
lock UG application to synchronize with other threads
virtual ParaSolution * createParaSolution()=0
create ParaSolution object by default constructor
virtual int getSize()=0
get number of UG processes or UG threads depending on run-time environment
virtual bool passTermToken(int rank)
pass termination token from the rank to the next
Definition: paraComm.h:211
virtual void unlockApp()=0
unlock UG application to synchronize with other threads
virtual int receive(void *bufer, int count, const int datatypeId, int source, const int tag)=0
receive function for standard ParaData types
virtual int bcast(void *buffer, int count, const int datatypeId, int root)=0
Some action need to be taken for fault tolerant, when the functions return. So, they rerun status val...
virtual ParaInstance * createParaInstance()=0
create ParaInstance object by default constructor
virtual int getRank()=0
get rank of this process or this thread depending on run-time environment
class for deterministic timer
virtual double getElapsedTime()=0
getter of the deterministic time
virtual void send(ParaComm *comm, int dest)=0
send function for ParaInitialStat object
class for instance data
Definition: paraInstance.h:51
virtual int bcast(ParaComm *comm, int rank, int method)=0
broadcast function to all solvers
class ParaParamSet
Definition: paraParamSet.h:850
bool getBoolParamValue(int param)
get bool parameter value
bool getBoolParamDefaultValue(int param)
get default value of bool parameter
double getRealParamValue(int param)
get real parameter value
int getIntParamValue(int param)
get int parameter value
const char * getStringParamValue(int param)
get string parameter value
class ParaRacingRampUpParamSet (parameter set for racing ramp-up)
class for solution
Definition: paraSolution.h:54
virtual void bcast(ParaComm *comm, int root)=0
broadcast solution data
virtual void send(ParaComm *comm, int destination, int tag)=0
send this object
class ParaSolver
Definition: paraSolver.h:71
double idleTimeToWaitNotificationId
idle time to wait a message within collecting mode
Definition: paraSolver.h:121
ParaParamSet * paraParams
ParaParamSet object.
Definition: paraSolver.h:83
ParaComm * paraComm
ParaCommunicator object.
Definition: paraSolver.h:82
void setTerminationMode(int tm)
set termination mode
Definition: paraSolver.h:837
ParaTask * currentTask
solving task
Definition: paraSolver.h:97
bool notificationProcessed
if true, notification is issued but not receive the corresponding LCB
Definition: paraSolver.h:140
int nParaTasksReceived
Counters related to this ParaSolver.
Definition: paraSolver.h:135
double idleTimeBetweenParaTasks
idle time between ParaTasks processing
Definition: paraSolver.h:119
bool racingIsInterrupted
indicate whether if racing phases is interrupted or not: true - interrupted
Definition: paraSolver.h:106
ParaSolution * globalBestIncumbentSolution
global best solution. However, this is not always feasible for the current sub-MIP
Definition: paraSolver.h:92
double previousIdleTimeToWaitToken
previous idle time to wait token
Definition: paraSolver.h:124
bool memoryLimitIsReached
indicate if memory limit is reached or not, when base solver has memory management feature
Definition: paraSolver.h:109
ParaRacingRampUpParamSet * winnerRacingParams
Winner ParaRacingRampUpParamSet object.
Definition: paraSolver.h:86
bool warmStarted
indicate whether if system is warm started or not
Definition: paraSolver.h:103
double eps
absolute values smaller than this are considered zero esp should be set in the constructor of the der...
Definition: paraSolver.h:142
double idleTimeToWaitToken
idle time to wait token
Definition: paraSolver.h:123
ParaTimer * paraTimer
timer for this ParaSolver
Definition: paraSolver.h:88
double idleTimeToWaitAckCompletion
idle time to wait acknowledgment of completion
Definition: paraSolver.h:122
bool isRacingWinner()
check if this solver is in racing ramp-up or not
Definition: paraSolver.h:530
double previousStopTime
Idle Times.
Definition: paraSolver.h:117
double idleTimeAfterLastParaTask
idle time after the last ParaTask was solved
Definition: paraSolver.h:120
double globalBestIncumbentValue
global best incumbent value
Definition: paraSolver.h:91
bool racingWinner
indicate racing ramp-up winner or not: true - winner
Definition: paraSolver.h:107
ParaDeterministicTimer * paraDetTimer
deterministic timer for this ParaSolver
Definition: paraSolver.h:89
double idleTimeToFirstParaTask
idle time to start solving the first ParaTask
Definition: paraSolver.h:118
int nParaTasksSolved
number of ParaTasks solved ( received ) in this ParaSolver
Definition: paraSolver.h:136
int terminationMode
indicate that termination mode 0: no termination mode 1: normal termination mode 2: interrupted termi...
Definition: paraSolver.h:99
ParaInstance * paraInstance
root problem instance
Definition: paraSolver.h:96
bool isRootTask()
check if root task or not
Definition: paraTask.h:624
TaskId getGeneratorTaskId()
getter of generator task id
Definition: paraTask.h:814
TaskId getTaskId()
getter of task id
Definition: paraTask.h:794
ParaDiffSubproblem * getDiffSubproblem()
getter of diffSubproblem
Definition: paraTask.h:854
void setOffset(double time)
Definition: paraTimer.h:74
virtual double getElapsedTime()=0
get elapsed time
static ScipParaCommTh * comm
Definition: fscip.cpp:73
static bool interrupted
Definition: fscip.cpp:78
static const int CustomizedToSharedMemory
static const int AddDualBoundCons
struct ParaSCIP::LocalNodeInfo_t LocalNodeInfo
static const int MemoryLimit
static const int RacingParamsDirPath
Definition: paraParamSet.h:131
static const int NoSolverPresolvingAtRootDefaultSet
static const int CompTerminatedByInterruptRequest
Definition: paraDef.h:184
static const int TagLbBoundTightenedBound
Definition: bbParaTagDef.h:63
static const int NoUpperBoundTransferInRacing
static const int UseRootNodeCuts
static const int NoTerminationMode
termination mode
Definition: paraSolver.h:60
static const int NoAggressiveSeparatorInRacing
static const int TagLbBoundTightenedIndex
Definition: bbParaTagDef.h:62
static const int ProvingRun
static const int ParaINT
Definition: paraComm.h:66
static const int InstanceTransferMethod
static const int TagTerminated
Definition: paraTagDef.h:58
static const int TagUbBoundTightenedIndex
Definition: bbParaTagDef.h:64
static const int TransferConflictCuts
static const int CheckEffectOfRootNodePreprocesses
static const int NoPreprocessingInLC
static const int TimeLimit
Definition: paraParamSet.h:106
static const int TimeLimitTerminationMode
Definition: paraSolver.h:63
static const int ControlCollectingModeOnSolverSide
static const int Deterministic
Definition: paraParamSet.h:76
static const int RampUpPhaseProcess
static const int CompTerminatedNormally
Definition: paraDef.h:182
static const int InterruptedTerminationMode
Definition: paraSolver.h:62
static const int TagUbBoundTightenedBound
Definition: bbParaTagDef.h:65
static const int RacingStatBranching
static const int AllowTreeSearchRestart
static const int Quiet
Definition: paraParamSet.h:71
static const int NoSolverPresolvingAtRoot
static const int DistributeBestPrimalSolution
static const int ParaDOUBLE
Definition: paraComm.h:76
static const int SelfSplitTreeDepth
static const int SetAllDefaultsAfterRacing
static const int CommunicateTighterBoundsInRacing
#define PARA_COMM_CALL(paracommcall)
Definition: paraComm.h:47
#define THROW_LOGICAL_ERROR1(msg1)
Definition: paraDef.h:52
#define THROW_LOGICAL_ERROR2(msg1, msg2)
Definition: paraDef.h:69
#define EPSEQ(x, y, eps)
Definition: paraDef.h:166
#define THROW_LOGICAL_ERROR3(msg1, msg2, msg3)
Definition: paraDef.h:86
Base class for initial statistics collecting class.
#define DEF_SCIP_PARA_COMM(scip_para_comm, comm)
ParaInitialStat extension for SCIP solver.
Branching rule plug-in for SCIP solver.
Event handlr for communication point.
heuristic to update objlimit
SCIP message handler for ParaSCIP and FiberSCIP.
C++ wrapper for propagators.
node selector for self-split ramp-up
#define SCIP_MEMORY_COPY_FACTOR
ParaRacingRampUpParamSet extension for SCIP solver.
void setUserPlugins(UG::ParaInstance *instance)
static SCIP_DECL_CONFLICTEXEC(conflictExecCollector)
#define CONFLICTHDLR_PRIORITY
#define CONFLICTHDLR_NAME
#define CONFLICTHDLR_DESC
double memoryLimitOfSolverSCIP
Definition: fscip.cpp:71
long long virtualMemUsedAtLc
Definition: fscip.cpp:70