1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145
|
/************************************************************************************
Copyright (C) 2005-2008 Assefaw H. Gebremedhin, Arijit Tarafdar, Duc Nguyen,
Alex Pothen
This file is part of ColPack.
ColPack is free software: you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published
by the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
ColPack is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with ColPack. If not, see <http://www.gnu.org/licenses/>.
************************************************************************************/
#include "ColPackHeaders.h"
using namespace std;
namespace ColPack
{
//Private Function 2401
int BipartiteGraphPartialColoring::CalculateVertexColorClasses()
{
if(m_s_VertexColoringVariant.empty())
{
return(_FALSE);
}
if(m_i_LeftVertexColorCount != _UNKNOWN)
{
int i_TotalLeftVertexColors = STEP_UP(m_i_LeftVertexColorCount);
m_vi_LeftVertexColorFrequency.clear();
m_vi_LeftVertexColorFrequency.resize((unsigned) i_TotalLeftVertexColors, _FALSE);
int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
for(int i = 0; i < i_LeftVertexCount; i++)
{
m_vi_LeftVertexColorFrequency[m_vi_LeftVertexColors[i]]++;
}
for(int i = 0; i < i_TotalLeftVertexColors; i++)
{
if(m_i_LargestLeftVertexColorClassSize < m_vi_LeftVertexColorFrequency[i])
{
m_i_LargestLeftVertexColorClass = i;
m_i_LargestLeftVertexColorClassSize = m_vi_LeftVertexColorFrequency[i];
}
if(m_i_SmallestLeftVertexColorClassSize == _UNKNOWN)
{
m_i_SmallestLeftVertexColorClass = i;
m_i_SmallestLeftVertexColorClassSize = m_vi_LeftVertexColorFrequency[i];
}
else
if(m_i_SmallestLeftVertexColorClassSize > m_vi_LeftVertexColorFrequency[i])
{
m_i_SmallestLeftVertexColorClass = i;
m_i_SmallestLeftVertexColorClassSize = m_vi_LeftVertexColorFrequency[i];
}
}
m_d_AverageLeftVertexColorClassSize = i_LeftVertexCount / i_TotalLeftVertexColors;
}
if(m_i_RightVertexColorCount != _UNKNOWN)
{
int i_TotalRightVertexColors = STEP_UP(m_i_RightVertexColorCount);
m_vi_RightVertexColorFrequency.clear();
m_vi_RightVertexColorFrequency.resize((unsigned) i_TotalRightVertexColors, _FALSE);
int i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
for(int i = 0; i < i_RightVertexCount; i++)
{
m_vi_RightVertexColorFrequency[m_vi_RightVertexColors[i]]++;
}
for(int i = 0; i < i_TotalRightVertexColors; i++)
{
if(m_i_LargestRightVertexColorClassSize < m_vi_RightVertexColorFrequency[i])
{
m_i_LargestRightVertexColorClass = i;
m_i_LargestRightVertexColorClassSize = m_vi_RightVertexColorFrequency[i];
}
if(m_i_SmallestRightVertexColorClassSize == _UNKNOWN)
{
m_i_SmallestRightVertexColorClass = i;
m_i_SmallestRightVertexColorClassSize = m_vi_RightVertexColorFrequency[i];
}
else
if(m_i_SmallestRightVertexColorClassSize > m_vi_RightVertexColorFrequency[i])
{
m_i_SmallestRightVertexColorClass = i;
m_i_SmallestRightVertexColorClassSize = m_vi_RightVertexColorFrequency[i];
}
}
m_d_AverageRightVertexColorClassSize = i_RightVertexCount / i_TotalRightVertexColors;
}
return(_TRUE);
}
//Private Function 2402
int BipartiteGraphPartialColoring::CheckVertexColoring(string s_VertexColoringVariant)
{
if(m_s_VertexColoringVariant.compare(s_VertexColoringVariant) == 0)
{
return(_TRUE);
}
if(m_s_VertexColoringVariant.compare("ALL") != 0)
{
m_s_VertexColoringVariant = s_VertexColoringVariant;
}
if(m_s_VertexColoringVariant.compare("ROW_PARTIAL_DISTANCE_TWO") == 0)
{
if(m_s_VertexOrderingVariant.empty())
{
RowNaturalOrdering();
}
}
else
if(m_s_VertexColoringVariant.compare("COLUMN_PARTIAL_DISTANCE_TWO") == 0)
{
if(m_s_VertexOrderingVariant.empty())
{
ColumnNaturalOrdering();
}
}
else
{
if(m_s_VertexOrderingVariant.empty())
{
RowNaturalOrdering();
}
}
return(_FALSE);
}
//Public Constructor 2451
BipartiteGraphPartialColoring::BipartiteGraphPartialColoring()
{
Clear();
Seed_init();
}
//Public Destructor 2452
BipartiteGraphPartialColoring::~BipartiteGraphPartialColoring()
{
Clear();
Seed_reset();
}
void BipartiteGraphPartialColoring::Seed_init() {
seed_available = false;
i_seed_rowCount = 0;
dp2_Seed = NULL;
}
void BipartiteGraphPartialColoring::Seed_reset() {
if(seed_available) {
seed_available = false;
free_2DMatrix(dp2_Seed, i_seed_rowCount);
dp2_Seed = NULL;
i_seed_rowCount = 0;
}
}
//Virtual Function 2453
void BipartiteGraphPartialColoring::Clear()
{
BipartiteGraphPartialOrdering::Clear();
m_i_LeftVertexColorCount = _UNKNOWN;
m_i_RightVertexColorCount = _UNKNOWN;
m_i_VertexColorCount = _UNKNOWN;
m_i_ViolationCount =_UNKNOWN;
m_i_ColoringUnits = _UNKNOWN;
m_i_LargestLeftVertexColorClass = _UNKNOWN;
m_i_LargestRightVertexColorClass = _UNKNOWN;
m_i_LargestLeftVertexColorClassSize = _UNKNOWN;
m_i_LargestRightVertexColorClassSize = _UNKNOWN;
m_i_SmallestLeftVertexColorClass = _UNKNOWN;
m_i_SmallestRightVertexColorClass = _UNKNOWN;
m_i_SmallestLeftVertexColorClassSize = _UNKNOWN;
m_i_SmallestRightVertexColorClassSize = _UNKNOWN;
m_d_AverageLeftVertexColorClassSize = _UNKNOWN;
m_d_AverageRightVertexColorClassSize = _UNKNOWN;
m_d_ColoringTime = _UNKNOWN;
m_d_CheckingTime = _UNKNOWN;
m_s_VertexColoringVariant.clear();
m_vi_LeftVertexColors.clear();
m_vi_RightVertexColors.clear();
m_vi_LeftVertexColorFrequency.clear();
m_vi_RightVertexColorFrequency.clear();
return;
}
//Virtual Function 2454
void BipartiteGraphPartialColoring::Reset()
{
BipartiteGraphPartialOrdering::Reset();
m_i_LeftVertexColorCount = _UNKNOWN;
m_i_RightVertexColorCount = _UNKNOWN;
m_i_VertexColorCount = _UNKNOWN;
m_i_ViolationCount = _UNKNOWN;
m_i_ColoringUnits = _UNKNOWN;
m_i_LargestLeftVertexColorClass = _UNKNOWN;
m_i_LargestRightVertexColorClass = _UNKNOWN;
m_i_LargestLeftVertexColorClassSize = _UNKNOWN;
m_i_LargestRightVertexColorClassSize = _UNKNOWN;
m_i_SmallestLeftVertexColorClass = _UNKNOWN;
m_i_SmallestRightVertexColorClass = _UNKNOWN;
m_i_SmallestLeftVertexColorClassSize = _UNKNOWN;
m_i_SmallestRightVertexColorClassSize = _UNKNOWN;
m_d_AverageLeftVertexColorClassSize = _UNKNOWN;
m_d_AverageRightVertexColorClassSize = _UNKNOWN;
m_d_ColoringTime = _UNKNOWN;
m_d_CheckingTime = _UNKNOWN;
m_s_VertexColoringVariant.clear();
m_vi_LeftVertexColors.clear();
m_vi_RightVertexColors.clear();
m_vi_LeftVertexColorFrequency.clear();
m_vi_RightVertexColorFrequency.clear();
return;
}
int f(int x) {return x;}
int BipartiteGraphPartialColoring::PartialDistanceTwoRowColoring_OMP()
{
if(CheckVertexColoring("ROW_PARTIAL_DISTANCE_TWO"))
{
return(_TRUE);
}
int i_LeftVertexCount, i_RightVertexCount, i_CurrentVertex;
bool cont=false;
vector<int> vi_forbiddenColors, vi_VerticesToBeColored, vi_verticesNeedNewColor;
i_LeftVertexCount = (int) m_vi_LeftVertices.size() - 1;
i_RightVertexCount = (int)m_vi_RightVertices.size () - 1;
m_i_LeftVertexColorCount = m_i_RightVertexColorCount = m_i_VertexColorCount = 0;
// !!! do sections for this part ? forbiddenColors may need to be private for each thread
// resize the colors
m_vi_LeftVertexColors.resize ( i_LeftVertexCount, _UNKNOWN );
// resize the forbidden colors. !!! should be resize to max D2 degree instead
vi_forbiddenColors.resize ( i_LeftVertexCount, _UNKNOWN );
//Algo 4 - Line 2: U <- V . U is vi_VerticesToBeColored
vi_VerticesToBeColored.reserve(i_LeftVertexCount);
for(int i=0; i<i_LeftVertexCount; i++) {
vi_VerticesToBeColored.push_back(m_vi_OrderedVertices[i]);
//vi_VerticesToBeColored.push_back(i);
}
// !!! pick reasonable amount to reserve
vi_verticesNeedNewColor.reserve(i_LeftVertexCount);
int i_NumOfVerticesToBeColored = vi_VerticesToBeColored.size();
// cout<<"m_vi_Edges.size() = "<<m_vi_Edges.size()<<endl;
// cout<<"m_vi_LeftVertexColors.size() = "<<m_vi_LeftVertexColors.size()<<endl;
//Algo 4 - Line 3: while U != 0 ; do
while(i_NumOfVerticesToBeColored!=0) {
// cout<<"i_NumOfVerticesToBeColored = "<<i_NumOfVerticesToBeColored<<endl;
//Phase 1: tentative coloring
//Algo 4 - Line 4: for each right vertex v in U (in parallel) do
#pragma omp parallel for default(none) schedule(dynamic) shared(cout, i_NumOfVerticesToBeColored, vi_VerticesToBeColored) firstprivate(vi_forbiddenColors)
for(int i=0; i<i_NumOfVerticesToBeColored; i++) {
int v = vi_VerticesToBeColored[i];
// CoutLock::set(); cout<<"t"<< omp_get_thread_num() <<": i="<<i<<", left vertex v="<<v<<endl;CoutLock::unset();
//Algo 4 - Line 5: for each left vertex w in adj (v) do
for (int w=m_vi_LeftVertices [v]; w<m_vi_LeftVertices [v+1]; w++ ) {
// CoutLock::set();
// cout<<"\t t"<< omp_get_thread_num() <<": w="<<w<<", right vertex m_vi_Edges [w]="<<m_vi_Edges [w]<<endl;
// CoutLock::unset();
//Algo 4 - Line 6: mark color [w] as forbidden to vertex v. NOTE: !!! Not needed
//Algo 4 - Line 7: for each right vertex x in adj (w) and x != v do
for (int x=m_vi_RightVertices [m_vi_Edges [w]]; x<m_vi_RightVertices [m_vi_Edges [w]+1]; x++ ) {
//Algo 4 - Line 8: mark color [x] as forbidden to vertex v
if ( m_vi_LeftVertexColors [m_vi_Edges [x]] != _UNKNOWN ) {
// CoutLock::set();
// cout<<"\t\t t"<< omp_get_thread_num() <<": x="<<x<<endl;
// if(x>=m_vi_Edges.size()) {
// cout<<"ERR: x>=m_vi_Edges.size()"<<endl;
// PrintBipartiteGraph();
// Pause();
// }
// cout<<"\t\t left vertex m_vi_Edges [x] = "<<m_vi_Edges [x]<<endl;
// cout<<"\t\t m_vi_LeftVertexColors [m_vi_Edges [x]] = "<<m_vi_LeftVertexColors [m_vi_Edges [x]]<<endl;
// cout<<"\t\t v="<<v<<endl;
// CoutLock::unset();
// !!! each thread should has their own vi_forbiddenColors[] vector just to make sure the don't override each other
vi_forbiddenColors [m_vi_LeftVertexColors [m_vi_Edges [x]]] = v;
}
}
}
//Algo 4 - Line 9: Pick a permissible color c for vertex v using some strategy
int i_cadidateColor;
// First fit
{
i_cadidateColor = 0;
while(vi_forbiddenColors[i_cadidateColor]==v) i_cadidateColor++;
}
m_vi_LeftVertexColors[v] = i_cadidateColor;
if(m_i_LeftVertexColorCount < i_cadidateColor) {
m_i_LeftVertexColorCount = i_cadidateColor;
}
}
//Algo 4 - Line 10: R.clear() ; R denotes the set of vertices to be recolored
vi_verticesNeedNewColor.clear();
//Phase 2: conflict detection. For each vertex v in U, check and see if v need to be recolored
//Algo 4 - Line 11: for each vertex v in U (in parallel) do
#pragma omp parallel for default(none) schedule(dynamic) shared(i_NumOfVerticesToBeColored, vi_VerticesToBeColored,vi_verticesNeedNewColor, cout) private(cont)
for(int i=0; i<i_NumOfVerticesToBeColored; i++) {
//Algo 4 - Line 12: cont <- true ; cont is used to break from the outer loop below
cont = true;
int v = vi_VerticesToBeColored[i];
//Algo 4 - Line 13: for each vertex w in adj (v) and cont = true do
for (int w=m_vi_LeftVertices [v]; (w<m_vi_LeftVertices [v+1]) && (cont == true); w++ ) {
//Algo 4 - Line 14: if color [v] = color [w] and f (v) > f (w) then . NOTE: !!! Not needed
//Algo 4 - Line 15: add [v] to R ; break . NOTE: !!! Not needed
//Algo 4 - Line 16: for each vertex x in adj (w) and v != x do
for (int x=m_vi_RightVertices [m_vi_Edges [w]]; x<m_vi_RightVertices [m_vi_Edges [w]+1]; x++ ) {
//cout<<" m_vi_LeftVertexColors [m_vi_Edges [x]="<<x<<"]"<< m_vi_LeftVertexColors [m_vi_Edges [x]] <<endl;
// cout<<" m_vi_LeftVertexColors [v="<<v<<"]"<< m_vi_LeftVertexColors [v] <<endl;
//cout<<"f(v="<<v<<")="<<f(v)<<endl;
//cout<<"f(m_vi_Edges [x]="<<x<<")="<<f(m_vi_Edges [x])<<endl;
//Algo 4 - Line 17: if color [v] = color [x] and f (v) > f (x) then
if ( m_vi_LeftVertexColors [m_vi_Edges [x]] == m_vi_LeftVertexColors[v] && f(v) > f(m_vi_Edges [x]) ) {
//Algo 4 - Line 18: add [v] to R ; cont <- false; break
#pragma omp critical
{
vi_verticesNeedNewColor.push_back(v);
}
cont = false;
break;
}
}
}
}
//Algo 4 - Line 19: U <- R , i.e., vi_VerticesToBeColored <- vi_verticesNeedNewColor
vi_VerticesToBeColored.clear();
i_NumOfVerticesToBeColored = vi_verticesNeedNewColor.size();
vi_VerticesToBeColored.reserve(vi_verticesNeedNewColor.size());
for(int i=0; i<vi_verticesNeedNewColor.size(); i++) {
vi_VerticesToBeColored.push_back(vi_verticesNeedNewColor[i]);
}
/*
vi_VerticesToBeColored.resize(vi_verticesNeedNewColor.size());
#pragma omp parallel for default(none) schedule(static) shared(i_NumOfVerticesToBeColored,vi_VerticesToBeColored,vi_verticesNeedNewColor)
for(int i=0; i<i_NumOfVerticesToBeColored; i++) {
vi_VerticesToBeColored[i]=vi_verticesNeedNewColor[i];
}
//*/
}
// Note that m_i_LeftVertexColorCount has not been updated yet
m_i_VertexColorCount = m_i_LeftVertexColorCount;
return _TRUE;
}
int BipartiteGraphPartialColoring::PartialDistanceTwoRowColoring() {
#ifdef _OPENMP
return BipartiteGraphPartialColoring::PartialDistanceTwoRowColoring_OMP();
#else
return BipartiteGraphPartialColoring::PartialDistanceTwoRowColoring_serial();
#endif
}
//Public Function 2455
int BipartiteGraphPartialColoring::PartialDistanceTwoRowColoring_serial()
{
if(CheckVertexColoring("ROW_PARTIAL_DISTANCE_TWO"))
{
return(_TRUE);
}
int i, w, x, c;
int i_LeftVertexCount, i_CurrentVertex;
vector<int> vi_forbiddenColors;
i_LeftVertexCount = (int)m_vi_LeftVertices.size () - 1;
// resize the colors
m_vi_LeftVertexColors.resize ( i_LeftVertexCount, _UNKNOWN );
// resize the forbidden colors
vi_forbiddenColors.resize ( i_LeftVertexCount, _UNKNOWN );
m_i_LeftVertexColorCount = m_i_RightVertexColorCount = m_i_VertexColorCount = 0;
for ( i=0; i<i_LeftVertexCount; ++i )
{
i_CurrentVertex = m_vi_OrderedVertices[i];
for ( w=m_vi_LeftVertices [i_CurrentVertex]; w<m_vi_LeftVertices [i_CurrentVertex+1]; ++w )
{
for ( x=m_vi_RightVertices [m_vi_Edges [w]]; x<m_vi_RightVertices [m_vi_Edges [w]+1]; ++x )
{
if ( m_vi_LeftVertexColors [m_vi_Edges [x]] != _UNKNOWN )
{
vi_forbiddenColors [m_vi_LeftVertexColors [m_vi_Edges [x]]] = i_CurrentVertex;
}
}
}
// do color[vi] <-min {c>0:forbiddenColors[c]=/=vi
for ( c=0; c<i_LeftVertexCount; ++c )
{
if ( vi_forbiddenColors [c] != i_CurrentVertex )
{
m_vi_LeftVertexColors [i_CurrentVertex] = c;
if(m_i_LeftVertexColorCount < c)
{
m_i_LeftVertexColorCount = c;
}
break;
}
}
//
}
m_i_VertexColorCount = m_i_LeftVertexColorCount;
return ( _TRUE );
}
int BipartiteGraphPartialColoring::PartialDistanceTwoColumnColoring_OMP() {
if(CheckVertexColoring("COLUMN_PARTIAL_DISTANCE_TWO"))
{
return(_TRUE);
}
int i_LeftVertexCount, i_RightVertexCount, i_CurrentVertex;
bool cont=false;
vector<int> vi_forbiddenColors, vi_VerticesToBeColored, vi_verticesNeedNewColor;
i_LeftVertexCount = (int) m_vi_LeftVertices.size() - 1;
i_RightVertexCount = (int)m_vi_RightVertices.size () - 1;
m_i_LeftVertexColorCount = m_i_RightVertexColorCount = m_i_VertexColorCount = 0;
// !!! do sections for this part ? forbiddenColors may need to be private for each thread
// resize the colors
m_vi_RightVertexColors.resize ( i_RightVertexCount, _UNKNOWN );
// resize the forbidden colors. !!! should be resize to max D2 degree instead
vi_forbiddenColors.resize ( i_RightVertexCount, _UNKNOWN );
//Algo 4 - Line 2: U <- V . U is vi_VerticesToBeColored
vi_VerticesToBeColored.reserve(i_RightVertexCount);
for(int i=0; i<i_RightVertexCount; i++) {
vi_VerticesToBeColored.push_back(m_vi_OrderedVertices[i] - i_LeftVertexCount);
//vi_VerticesToBeColored.push_back(i);
}
// !!! pick reasonable amount to reserve
vi_verticesNeedNewColor.reserve(i_RightVertexCount);
int i_NumOfVerticesToBeColored = vi_VerticesToBeColored.size();
//Algo 4 - Line 3: while U != 0 ; do
while(i_NumOfVerticesToBeColored!=0) {
//cout<<"i_NumOfVerticesToBeColored = "<<i_NumOfVerticesToBeColored<<endl;
//Phase 1: tentative coloring
//Algo 4 - Line 4: for each right vertex v in U (in parallel) do
#pragma omp parallel for default(none) schedule(dynamic) shared(i_NumOfVerticesToBeColored, vi_VerticesToBeColored) firstprivate(vi_forbiddenColors)
for(int i=0; i<i_NumOfVerticesToBeColored; i++) {
int v = vi_VerticesToBeColored[i];
//Algo 4 - Line 5: for each left vertex w in adj (v) do
for (int w=m_vi_RightVertices [v]; w<m_vi_RightVertices [v+1]; w++ ) {
//Algo 4 - Line 6: mark color [w] as forbidden to vertex v. NOTE: !!! Not needed
//Algo 4 - Line 7: for each right vertex x in adj (w) and x != v do
for (int x=m_vi_LeftVertices [m_vi_Edges [w]]; x<m_vi_LeftVertices [m_vi_Edges [w]+1]; x++ ) {
//Algo 4 - Line 8: mark color [x] as forbidden to vertex v
if ( m_vi_RightVertexColors [m_vi_Edges [x]] != _UNKNOWN ) {
// !!! each thread should has their own vi_forbiddenColors[] vector just to make sure the don't override each other
vi_forbiddenColors [m_vi_RightVertexColors [m_vi_Edges [x]]] = v;
}
}
}
//Algo 4 - Line 9: Pick a permissible color c for vertex v using some strategy
int i_cadidateColor;
// First fit
{
i_cadidateColor = 0;
while(vi_forbiddenColors[i_cadidateColor]==v) i_cadidateColor++;
}
m_vi_RightVertexColors[v] = i_cadidateColor;
if(m_i_RightVertexColorCount < i_cadidateColor) {
m_i_RightVertexColorCount = i_cadidateColor;
}
}
//Algo 4 - Line 10: R.clear() ; R denotes the set of vertices to be recolored
vi_verticesNeedNewColor.clear();
//Phase 2: conflict detection. For each vertex v in U, check and see if v need to be recolored
//Algo 4 - Line 11: for each vertex v in U (in parallel) do
#pragma omp parallel for default(none) schedule(dynamic) shared(i_NumOfVerticesToBeColored, vi_VerticesToBeColored,vi_verticesNeedNewColor, cout) private(cont)
for(int i=0; i<i_NumOfVerticesToBeColored; i++) {
//Algo 4 - Line 12: cont <- true ; cont is used to break from the outer loop below
cont = true;
int v = vi_VerticesToBeColored[i];
//Algo 4 - Line 13: for each vertex w in adj (v) and cont = true do
for (int w=m_vi_RightVertices [v]; (w<m_vi_RightVertices [v+1]) && (cont == true); w++ ) {
//Algo 4 - Line 14: if color [v] = color [w] and f (v) > f (w) then . NOTE: !!! Not needed
//Algo 4 - Line 15: add [v] to R ; break . NOTE: !!! Not needed
//Algo 4 - Line 16: for each vertex x in adj (w) and v != x do
for (int x=m_vi_LeftVertices [m_vi_Edges [w]]; x<m_vi_LeftVertices [m_vi_Edges [w]+1]; x++ ) {
//cout<<" m_vi_RightVertexColors [m_vi_Edges [x]="<<x<<"]"<< m_vi_RightVertexColors [m_vi_Edges [x]] <<endl;
// cout<<" m_vi_RightVertexColors [v="<<v<<"]"<< m_vi_RightVertexColors [v] <<endl;
//cout<<"f(v="<<v<<")="<<f(v)<<endl;
//cout<<"f(m_vi_Edges [x]="<<x<<")="<<f(m_vi_Edges [x])<<endl;
//Algo 4 - Line 17: if color [v] = color [x] and f (v) > f (x) then
if ( m_vi_RightVertexColors [m_vi_Edges [x]] == m_vi_RightVertexColors[v] && f(v) > f(m_vi_Edges [x]) ) {
//Algo 4 - Line 18: add [v] to R ; cont <- false; break
#pragma omp critical
{
vi_verticesNeedNewColor.push_back(v);
}
cont = false;
break;
}
}
}
}
//Algo 4 - Line 19: U <- R , i.e., vi_VerticesToBeColored <- vi_verticesNeedNewColor
vi_VerticesToBeColored.clear();
i_NumOfVerticesToBeColored = vi_verticesNeedNewColor.size();
vi_VerticesToBeColored.reserve(vi_verticesNeedNewColor.size());
for(int i=0; i<vi_verticesNeedNewColor.size(); i++) {
vi_VerticesToBeColored.push_back(vi_verticesNeedNewColor[i]);
}
/*
vi_VerticesToBeColored.resize(vi_verticesNeedNewColor.size());
#pragma omp parallel for default(none) schedule(static) shared(i_NumOfVerticesToBeColored,vi_VerticesToBeColored,vi_verticesNeedNewColor)
for(int i=0; i<i_NumOfVerticesToBeColored; i++) {
vi_VerticesToBeColored[i]=vi_verticesNeedNewColor[i];
}
//*/
}
//note that m_i_RightVertexColorCount has not been updated yet
m_i_VertexColorCount = m_i_RightVertexColorCount;
return _TRUE;
}
int BipartiteGraphPartialColoring::PartialDistanceTwoColumnColoring() {
#ifdef _OPENMP
return BipartiteGraphPartialColoring::PartialDistanceTwoColumnColoring_OMP();
#else
return BipartiteGraphPartialColoring::PartialDistanceTwoColumnColoring_serial();
#endif
}
//Public Function 2456
int BipartiteGraphPartialColoring::PartialDistanceTwoColumnColoring_serial()
{
if(CheckVertexColoring("COLUMN_PARTIAL_DISTANCE_TWO"))
{
return(_TRUE);
}
int i, w, x, c;
int i_LeftVertexCount, i_RightVertexCount, i_CurrentVertex;
vector<int> vi_forbiddenColors;
i_LeftVertexCount = (int) m_vi_LeftVertices.size() - 1;
i_RightVertexCount = (int)m_vi_RightVertices.size () - 1;
// resize the colors
m_vi_RightVertexColors.resize ( i_RightVertexCount, _UNKNOWN );
// resize the forbidden colors
vi_forbiddenColors.resize ( i_RightVertexCount, _UNKNOWN );
m_i_LeftVertexColorCount = m_i_RightVertexColorCount = m_i_VertexColorCount = 0;
//cout<<" i_RightVertexCount = " <<i_RightVertexCount<<endl;
for ( i=0; i<i_RightVertexCount; ++i )
{
i_CurrentVertex = m_vi_OrderedVertices[i] - i_LeftVertexCount;
for ( w=m_vi_RightVertices [i_CurrentVertex]; w<m_vi_RightVertices [i_CurrentVertex+1]; ++w )
{
for ( x=m_vi_LeftVertices [m_vi_Edges [w]]; x<m_vi_LeftVertices [m_vi_Edges [w]+1]; ++x )
{
if ( m_vi_RightVertexColors [m_vi_Edges [x]] != _UNKNOWN )
{
vi_forbiddenColors [m_vi_RightVertexColors [m_vi_Edges [x]]] = i_CurrentVertex;
}
}
}
// do color[vi] <-min {c>0:forbiddenColors[c]=/=vi
for ( c=0; c<i_RightVertexCount; ++c )
{
if ( vi_forbiddenColors [c] != i_CurrentVertex )
{
m_vi_RightVertexColors [i_CurrentVertex] = c;
if(m_i_RightVertexColorCount < c)
{
m_i_RightVertexColorCount = c;
}
break;
}
}
//
}
m_i_VertexColorCount = m_i_RightVertexColorCount;
return ( _TRUE );
}
//Public Function 2457
int BipartiteGraphPartialColoring::CheckPartialDistanceTwoRowColoring()
{
for(int i=0;i<(signed) m_vi_LeftVertices.size()-1;i++)
//for each of left vertices, find its D1 neighbour (right vertices)
{
for(int j=m_vi_LeftVertices[i];j<m_vi_LeftVertices[i+1];j++)
{
for(int k=m_vi_RightVertices[m_vi_Edges[j]]; k<m_vi_RightVertices[m_vi_Edges[j]+1];k++)
//for each of the right vertices, find its D1 neighbour (left vertices exclude the original left)
{
if(m_vi_Edges[k]==i) continue;
if(m_vi_LeftVertexColors[m_vi_Edges[k]]==m_vi_LeftVertexColors[i])
{
cout<<"Left vertices "<<i+1<<" and "<<m_vi_Edges[k]+1<< " (connected by right vectex "<<m_vi_Edges[j]+1<<") have the same color ("<<m_vi_LeftVertexColors[i]<<")"<<endl;
return _FALSE;
}
}
}
}
return _TRUE;
}
//Public Function 2458
int BipartiteGraphPartialColoring::CheckPartialDistanceTwoColumnColoring()
{
for(int i=0;i<(signed) m_vi_RightVertices.size()-1;i++)
//for each of right vertices, find its D1 neighbour (left vertices)
{
for(int j=m_vi_RightVertices[i];j<m_vi_RightVertices[i+1];j++)
{
for(int k=m_vi_LeftVertices[m_vi_Edges[j]]; k<m_vi_LeftVertices[m_vi_Edges[j]+1];k++)
//for each of the left vertices, find its D1 neighbour (right vertices exclude the original right)
{
if(m_vi_Edges[k]==i) continue;
if(m_vi_RightVertexColors[m_vi_Edges[k]]==m_vi_RightVertexColors[i])
{
cout<<"Right vertices "<<i+1<<" and "<<m_vi_Edges[k]+1<< " (connected by left vectex "<<m_vi_Edges[j]+1<<") have the same color ("<<m_vi_RightVertexColors[i]<<")"<<endl;
return _FALSE;
}
}
}
}
return (_TRUE);
}
//Public Function 2459
int BipartiteGraphPartialColoring::GetLeftVertexColorCount()
{
if(m_i_LeftVertexColorCount<0 && GetVertexColoringVariant() == "Row Partial Distance Two" ) {
for(int i=0; i<m_vi_LeftVertexColors.size();i++) {
if(m_i_LeftVertexColorCount<m_vi_LeftVertexColors[i]) m_i_LeftVertexColorCount = m_vi_LeftVertexColors[i];
}
}
return(STEP_UP(m_i_LeftVertexColorCount));
}
//Public Function 2460
int BipartiteGraphPartialColoring::GetRightVertexColorCount()
{
if(m_i_RightVertexColorCount<0 && GetVertexColoringVariant() == "Column Partial Distance Two" ) {
for(int i=0; i<m_vi_RightVertexColors.size();i++) {
if(m_i_RightVertexColorCount<m_vi_RightVertexColors[i]) m_i_RightVertexColorCount = m_vi_RightVertexColors[i];
}
}
return(STEP_UP(m_i_RightVertexColorCount));
}
//Public Function 2461
int BipartiteGraphPartialColoring::GetVertexColorCount()
{
if(m_i_VertexColorCount<0 && GetVertexColoringVariant() != "Unknown" ) {
if(GetVertexColoringVariant() == "Row Partial Distance Two") {
m_i_VertexColorCount = GetLeftVertexColorCount() - 1;
}
else { // GetVertexColoringVariant() == "Column Partial Distance Two"
m_i_VertexColorCount = GetRightVertexColorCount() - 1;
}
}
return(STEP_UP(m_i_VertexColorCount));
}
//Public Function 2462
string BipartiteGraphPartialColoring::GetVertexColoringVariant()
{
if(m_s_VertexColoringVariant.compare("ROW_PARTIAL_DISTANCE_TWO") == 0)
{
return("Row Partial Distance Two");
}
else
if(m_s_VertexColoringVariant.compare("COLUMN_PARTIAL_DISTANCE_TWO") == 0)
{
return("Column Partial Distance Two");
}
else
{
return("Unknown");
}
}
//Public Function 2463
void BipartiteGraphPartialColoring::GetLeftVertexColors(vector<int> &output)
{
output = (m_vi_LeftVertexColors);
}
//Public Function 2464
void BipartiteGraphPartialColoring::GetRightVertexColors(vector<int> &output)
{
output = (m_vi_RightVertexColors);
}
//Public Function 2465
void BipartiteGraphPartialColoring::PrintRowPartialColors()
{
int i;
int i_LeftVertexCount;
string _SLASH("/");
StringTokenizer SlashTokenizer(m_s_InputFile, _SLASH);
m_s_InputFile = SlashTokenizer.GetLastToken();
i_LeftVertexCount = (signed) m_vi_LeftVertexColors.size();
cout<<endl;
cout<<"Bipartite Graph | Row Partial Coloring | Row Vertices | Vertex Colors "<<m_s_InputFile<<endl;
cout<<endl;
for(i=0; i<i_LeftVertexCount; i++)
{
cout<<STEP_UP(i)<<"\t"<<" : "<<STEP_UP(m_vi_LeftVertexColors[i])<<endl;
}
cout<<endl;
cout<<"[Total Row Colors = "<<GetLeftVertexColorCount()<<"]"<<endl;
cout<<endl;
return;
}
//Public Function 2466
void BipartiteGraphPartialColoring::PrintColumnPartialColors()
{
int i;
int i_RightVertexCount;
string _SLASH("/");
StringTokenizer SlashTokenizer(m_s_InputFile, _SLASH);
m_s_InputFile = SlashTokenizer.GetLastToken();
i_RightVertexCount = (signed) m_vi_RightVertexColors.size();
cout<<endl;
cout<<"Bipartite Graph | Column Partial Coloring | Column Vertices | Vertex Colors | "<<m_s_InputFile<<endl;
cout<<endl;
for(i=0; i<i_RightVertexCount; i++)
{
cout<<STEP_UP(i)<<"\t"<<" : "<<STEP_UP(m_vi_RightVertexColors[i])<<endl;
}
cout<<endl;
cout<<"[Total Column Colors = "<<GetRightVertexColorCount()<<"]"<<endl;
cout<<endl;
return;
}
//Public Function 2467
void BipartiteGraphPartialColoring::PrintRowPartialColoringMetrics()
{
string _SLASH("/");
StringTokenizer SlashTokenizer(m_s_InputFile, _SLASH);
string s_InputFile = SlashTokenizer.GetLastToken();
cout<<endl;
cout<<GetVertexColoringVariant()<<" Bicoloring | "<<GetVertexOrderingVariant()<<" Ordering | "<<s_InputFile<<endl;
cout<<endl;
cout<<endl;
cout<<"[Total Row Colors = "<<STEP_UP(m_i_VertexColorCount)<<"; Violation Count = "<<m_i_ViolationCount<<"]"<<endl;
cout<<"[Row Vertex Count = "<<STEP_DOWN(m_vi_LeftVertices.size())<<"; Column Vertex Count = "<<STEP_DOWN(m_vi_RightVertices.size())<<endl;
cout<<"[Ordering Time = "<<m_d_OrderingTime<<"; Coloring Time = "<<m_d_ColoringTime<<"; Checking Time = "<<m_d_CheckingTime<<"]"<<endl;
cout<<endl;
}
//Public Function 2468
void BipartiteGraphPartialColoring::PrintColumnPartialColoringMetrics()
{
string _SLASH("/");
StringTokenizer SlashTokenizer(m_s_InputFile, _SLASH);
string s_InputFile = SlashTokenizer.GetLastToken();
cout<<endl;
cout<<GetVertexColoringVariant()<<" Bicoloring | "<<GetVertexOrderingVariant()<<" Ordering | "<<s_InputFile<<endl;
cout<<endl;
cout<<endl;
cout<<"[Total Column Colors = "<<STEP_UP(m_i_VertexColorCount)<<"; Violation Count = "<<m_i_ViolationCount<<"]"<<endl;
cout<<"[Row Vertex Count = "<<STEP_DOWN(m_vi_LeftVertices.size())<<"; Column Vertex Count = "<<STEP_DOWN(m_vi_RightVertices.size())<<endl;
cout<<"[Ordering Time = "<<m_d_OrderingTime<<"; Coloring Time = "<<m_d_ColoringTime<<"; Checking Time = "<<m_d_CheckingTime<<"]"<<endl;
cout<<endl;
}
//Public Function 2469
void BipartiteGraphPartialColoring::PrintVertexPartialColorClasses()
{
if(CalculateVertexColorClasses() != _TRUE)
{
cout<<endl;
cout<<"Vertex Partial Color Classes | "<<m_s_VertexColoringVariant<<" Coloring | "<<m_s_VertexOrderingVariant<<" Ordering | "<<m_s_InputFile<<" | Vertex Partial Colors Not Set"<<endl;
cout<<endl;
return;
}
if(m_i_LeftVertexColorCount != _UNKNOWN)
{
cout<<endl;
cout<<"Row Color Classes | "<<m_s_VertexColoringVariant<<" Coloring | "<<m_s_VertexOrderingVariant<<" Ordering | "<<m_s_InputFile<<endl;
cout<<endl;
int i_TotalLeftVertexColors = STEP_UP(m_i_LeftVertexColorCount);
for(int i = 0; i < i_TotalLeftVertexColors; i++)
{
if(m_vi_LeftVertexColorFrequency[i] <= 0)
{
continue;
}
cout<<"Color "<<STEP_UP(i)<<" : "<<m_vi_LeftVertexColorFrequency[i]<<endl;
}
cout<<endl;
cout<<"[Largest Row Color Class : "<<STEP_UP(m_i_LargestLeftVertexColorClass)<<"; Largest Row Color Class Size : "<<m_i_LargestLeftVertexColorClassSize<<"]"<<endl;
cout<<"[Smallest Row Color Class : "<<STEP_UP(m_i_SmallestLeftVertexColorClass)<<"; Smallest Row Color Class Size : "<<m_i_SmallestLeftVertexColorClassSize<<"]"<<endl;
cout<<"[Average Row Color Class Size : "<<m_d_AverageLeftVertexColorClassSize<<"]"<<endl;
cout<<endl;
}
if(m_i_RightVertexColorCount != _UNKNOWN)
{
cout<<endl;
cout<<"Column Color Classes | "<<m_s_VertexColoringVariant<<" Coloring | "<<m_s_VertexOrderingVariant<<" Ordering | "<<m_s_InputFile<<endl;
cout<<endl;
int i_TotalRightVertexColors = STEP_UP(m_i_RightVertexColorCount);
for(int i = 0; i < i_TotalRightVertexColors; i++)
{
if(m_vi_RightVertexColorFrequency[i] <= 0)
{
continue;
}
cout<<"Color "<<STEP_UP(i)<<" : "<<m_vi_RightVertexColorFrequency[i]<<endl;
}
cout<<endl;
cout<<"[Largest Column Color Class : "<<STEP_UP(m_i_LargestRightVertexColorClass)<<"; Largest Column Color Class Size : "<<m_i_LargestRightVertexColorClassSize<<"]"<<endl;
cout<<"[Smallest Column Color Class : "<<STEP_UP(m_i_SmallestRightVertexColorClass)<<"; Smallest Column Color Class Size : "<<m_i_SmallestRightVertexColorClassSize<<"]"<<endl;
cout<<"[Average Column Color Class Size : "<<m_d_AverageRightVertexColorClassSize<<"]"<<endl;
cout<<endl;
}
return;
}
double** BipartiteGraphPartialColoring::GetLeftSeedMatrix(int* i_SeedRowCount, int* i_SeedColumnCount) {
if(seed_available) Seed_reset();
dp2_Seed = GetLeftSeedMatrix_unmanaged(i_SeedRowCount, i_SeedColumnCount);
i_seed_rowCount = *i_SeedRowCount;
seed_available = true;
return dp2_Seed;
}
double** BipartiteGraphPartialColoring::GetRightSeedMatrix(int* i_SeedRowCount, int* i_SeedColumnCount) {
if(seed_available) Seed_reset();
dp2_Seed = GetRightSeedMatrix_unmanaged(i_SeedRowCount, i_SeedColumnCount);
i_seed_rowCount = *i_SeedRowCount;
seed_available = true;
return dp2_Seed;
}
double** BipartiteGraphPartialColoring::GetLeftSeedMatrix_unmanaged(int* i_SeedRowCount, int* i_SeedColumnCount) {
int i_size = m_vi_LeftVertexColors.size();
int i_num_of_colors = GetLeftVertexColorCount();
(*i_SeedRowCount) = i_num_of_colors;
(*i_SeedColumnCount) = i_size;
if(i_num_of_colors == 0 || i_size == 0) return NULL;
double** Seed = new double*[i_num_of_colors];
// allocate and initialize Seed matrix
for (int i=0; i<i_num_of_colors; i++) {
Seed[i] = new double[i_size];
for(int j=0; j<i_size; j++) Seed[i][j]=0.;
}
// populate Seed matrix
for (int i=0; i < i_size; i++) {
Seed[m_vi_LeftVertexColors[i]][i] = 1.;
}
return Seed;
}
double** BipartiteGraphPartialColoring::GetRightSeedMatrix_unmanaged(int* i_SeedRowCount, int* i_SeedColumnCount) {
int i_size = m_vi_RightVertexColors.size();
int i_num_of_colors = GetRightVertexColorCount();
(*i_SeedRowCount) = i_size;
(*i_SeedColumnCount) = i_num_of_colors;
if(i_num_of_colors == 0 || i_size == 0) return NULL;
double** Seed = new double*[i_size];
// allocate and initialize Seed matrix
for (int i=0; i<i_size; i++) {
Seed[i] = new double[i_num_of_colors];
for(int j=0; j<i_num_of_colors; j++) Seed[i][j]=0.;
}
// populate Seed matrix
for (int i=0; i < i_size; i++) {
// if(m_vi_RightVertexColors[i]>=i_num_of_colors) {
// cout<<"i="<<i<<endl;
// cout<<"m_vi_RightVertexColors[i]="<<m_vi_RightVertexColors[i]<<endl;
// cout<<"i_num_of_colors="<<i_num_of_colors<<endl;
// }
Seed[i][m_vi_RightVertexColors[i]] = 1.;
}
return Seed;
}
void BipartiteGraphPartialColoring::PrintPartialColoringMetrics() {
if ( m_s_VertexColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") {
PrintColumnPartialColoringMetrics();
}
else if (m_s_VertexColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") {
PrintRowPartialColoringMetrics();
}
else { // Unrecognized Coloring Method
cerr<<" Unknown Partial Distance Two Coloring Method "<<m_s_VertexColoringVariant
<<". Please use a legal Method before calling PrintPartialColors()."<<endl;
}
}
double** BipartiteGraphPartialColoring::GetSeedMatrix(int* i_SeedRowCount, int* i_SeedColumnCount) {
if ( m_s_VertexColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") {
return GetRightSeedMatrix(i_SeedRowCount, i_SeedColumnCount);
}
else if (m_s_VertexColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") {
return GetLeftSeedMatrix(i_SeedRowCount, i_SeedColumnCount);
}
else { // Unrecognized Coloring Method
cerr<<" Unknown Partial Distance Two Coloring Method "<<m_s_VertexColoringVariant
<<". Please use a legal Method before calling PrintPartialColors()."<<endl;
}
return NULL;
}
double** BipartiteGraphPartialColoring::GetSeedMatrix_unmanaged(int* i_SeedRowCount, int* i_SeedColumnCount) {
if ( m_s_VertexColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") {
return GetRightSeedMatrix_unmanaged(i_SeedRowCount, i_SeedColumnCount);
}
else if (m_s_VertexColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") {
return GetLeftSeedMatrix_unmanaged(i_SeedRowCount, i_SeedColumnCount);
}
else { // Unrecognized Coloring Method
cerr<<" Unknown Partial Distance Two Coloring Method "<<m_s_VertexColoringVariant
<<". Please use a legal Method before calling PrintPartialColors()."<<endl;
}
return NULL;
}
void BipartiteGraphPartialColoring::GetVertexPartialColors(vector<int> &output)
{
if ( m_s_VertexColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") {
GetRightVertexColors(output);
}
else if (m_s_VertexColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") {
GetLeftVertexColors(output);
}
else { // Unrecognized Coloring Method
cerr<<" Unknown Partial Distance Two Coloring Method: "<<m_s_VertexColoringVariant
<<". Please use a legal Method before calling GetVertexColors()."<<endl;
}
}
void BipartiteGraphPartialColoring::PrintPartialColors() {
if ( m_s_VertexColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") {
PrintColumnPartialColors();
}
else if (m_s_VertexColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") {
PrintRowPartialColors();
}
else { // Unrecognized Coloring Method
cerr<<" Unknown Partial Distance Two Coloring Method "<<m_s_VertexColoringVariant
<<". Please use a legal Method before calling PrintPartialColors()."<<endl;
}
}
int BipartiteGraphPartialColoring::CheckPartialDistanceTwoColoring() {
if ( m_s_VertexColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") {
return CheckPartialDistanceTwoColumnColoring();
}
else if (m_s_VertexColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") {
return CheckPartialDistanceTwoRowColoring();
}
else { // Unrecognized Coloring Method
cerr<<" Unknown Partial Distance Two Coloring Method: "<<m_s_VertexColoringVariant
<<". Please use a legal Method before calling CheckPartialDistanceTwoColoring()."<<endl;
return _FALSE;
}
}
double BipartiteGraphPartialColoring::GetVertexColoringTime() {
return m_d_ColoringTime;
}
void BipartiteGraphPartialColoring::SetVertexColoringVariant(string s_VertexColoringVariant) {
m_s_VertexColoringVariant = s_VertexColoringVariant;
}
}
|