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
|
/****************************************************************/
/* Parallel Combinatorial BLAS Library (for Graph Computations) */
/* version 1.6 -------------------------------------------------*/
/* date: 6/15/2017 ---------------------------------------------*/
/* authors: Ariful Azad, Aydin Buluc --------------------------*/
/****************************************************************/
/*
Copyright (c) 2010-2017, The Regents of the University of California
Permission is hereby granted, free of charge, to any person obtaining a copy
of this software and associated documentation files (the "Software"), to deal
in the Software without restriction, including without limitation the rights
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
copies of the Software, and to permit persons to whom the Software is
furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
THE SOFTWARE.
*/
#include "SpCCols.h"
#include "Deleter.h"
#include <algorithm>
#include <functional>
#include <vector>
#include <climits>
#include <iomanip>
#include <cassert>
namespace combblas {
/****************************************************************************/
/********************* PUBLIC CONSTRUCTORS/DESTRUCTORS **********************/
/****************************************************************************/
template <class IT, class NT>
const IT SpCCols<IT,NT>::esscount = static_cast<IT>(3);
template <class IT, class NT>
SpCCols<IT,NT>::SpCCols():csc(NULL), m(0), n(0), nnz(0), splits(0){
}
// Allocate all the space necessary
template <class IT, class NT>
SpCCols<IT,NT>::SpCCols(IT size, IT nRow, IT nCol)
:m(nRow), n(nCol), nnz(size), splits(0)
{
if(nnz > 0)
csc = new Csc<IT,NT>(nnz, n);
else
csc = NULL;
}
template <class IT, class NT>
SpCCols<IT,NT>::~SpCCols()
{
if(nnz > 0)
{
if(csc != NULL)
{
if(splits > 0)
{
for(int i=0; i<splits; ++i)
delete cscarr[i];
delete [] cscarr;
}
else
{
delete csc;
}
}
}
}
// Copy constructor (constructs a new object. i.e. this is NEVER called on an existing object)
// Derived's copy constructor can safely call Base's default constructor as base has no data members
template <class IT, class NT>
SpCCols<IT,NT>::SpCCols(const SpCCols<IT,NT> & rhs)
: m(rhs.m), n(rhs.n), nnz(rhs.nnz), splits(rhs.splits)
{
if(splits > 0)
{
for(int i=0; i<splits; ++i)
CopyCsc(rhs.cscarr[i]);
}
else
{
CopyCsc(rhs.csc);
}
}
/**
* Constructor for converting SpTuples matrix -> SpCCols
* @param[in] rhs if transpose=true,
* \n then rhs is assumed to be a row sorted SpTuples object
* \n else rhs is assumed to be a column sorted SpTuples object
**/
template <class IT, class NT>
SpCCols<IT,NT>::SpCCols(const SpTuples<IT, NT> & rhs, bool transpose)
: m(rhs.m), n(rhs.n), nnz(rhs.nnz), splits(0)
{
if(nnz == 0) // m by n matrix of complete zeros
{
if(transpose) std::swap(m,n);
csc = NULL;
}
else
{
if(transpose)
{
std::swap(m,n);
csc = new Csc<IT,NT>(nnz,n); // the swap is already done here
std::vector< std::pair<IT,NT> > tosort (nnz);
std::vector<IT> work(n+1, (IT) 0 ); // workspace, zero initialized, first entry stays zero
for (IT k = 0 ; k < nnz ; ++k)
{
IT tmp = rhs.rowindex(k);
work [ tmp+1 ]++ ; // column counts (i.e, w holds the "col difference array")
}
if(nnz > 0)
{
std::partial_sum(work.begin(), work.end(), work.begin());
std::copy(work.begin(), work.end(), csc->jc);
IT last;
for (IT k = 0 ; k < nnz ; ++k)
{
tosort[ work[ rhs.rowindex(k) ]++] = std::make_pair( rhs.colindex(k), rhs.numvalue(k));
}
#ifdef _OPENMP
#pragma omp parallel for
#endif
for(int i=0; i< n; ++i)
{
sort(tosort.begin() + csc->jc[i], tosort.begin() + csc->jc[i+1]);
IT ind;
typename std::vector<std::pair<IT,NT> >::iterator itr; // iterator is a dependent name
for(itr = tosort.begin() + csc->jc[i], ind = csc->jc[i]; itr != tosort.begin() + csc->jc[i+1]; ++itr, ++ind)
{
csc->ir[ind] = itr->first;
csc->num[ind] = itr->second;
}
}
}
}
else
{
csc = new Csc<IT,NT>(nnz,n); // the swap is already done here
std::vector< std::pair<IT,NT> > tosort (nnz);
std::vector<IT> work(n+1, (IT) 0 ); // workspace, zero initialized, first entry stays zero
for (IT k = 0 ; k < nnz ; ++k)
{
IT tmp = rhs.colindex(k);
work [ tmp+1 ]++ ; // column counts (i.e, w holds the "col difference array")
}
if(nnz > 0)
{
std::partial_sum(work.begin(), work.end(), work.begin());
std::copy(work.begin(), work.end(), csc->jc);
IT last;
for (IT k = 0 ; k < nnz ; ++k)
{
tosort[ work[ rhs.colindex(k) ]++] = std::make_pair( rhs.rowindex(k), rhs.numvalue(k));
}
#ifdef _OPENMP
#pragma omp parallel for
#endif
for(int i=0; i< n; ++i)
{
sort(tosort.begin() + csc->jc[i], tosort.begin() + csc->jc[i+1]);
IT ind;
typename std::vector<std::pair<IT,NT> >::iterator itr; // iterator is a dependent name
for(itr = tosort.begin() + csc->jc[i], ind = csc->jc[i]; itr != tosort.begin() + csc->jc[i+1]; ++itr, ++ind)
{
csc->ir[ind] = itr->first;
csc->num[ind] = itr->second;
}
}
}
}
}
}
/****************************************************************************/
/************************** PUBLIC OPERATORS ********************************/
/****************************************************************************/
/**
* The assignment operator operates on an existing object
* The assignment operator is the only operator that is not inherited.
* But there is no need to call base's assigment operator as it has no data members
*/
template <class IT, class NT>
SpCCols<IT,NT> & SpCCols<IT,NT>::operator=(const SpCCols<IT,NT> & rhs)
{
// this pointer stores the address of the class instance
// check for self assignment using address comparison
if(this != &rhs)
{
if(csc != NULL && nnz > 0)
{
delete csc;
}
if(rhs.csc != NULL)
{
csc = new Csc<IT,NT>(*(rhs.csc));
nnz = rhs.nnz;
}
else
{
csc = NULL;
nnz = 0;
}
m = rhs.m;
n = rhs.n;
splits = rhs.splits;
}
return *this;
}
template <class IT, class NT>
void SpCCols<IT,NT>::RowSplit(int numsplits)
{
splits = numsplits;
IT perpiece = m / splits;
std::vector<IT> nnzs(splits, 0);
std::vector < std::vector < std::tuple<IT,IT,NT> > > colrowpairs(splits);
std::vector< std::vector<IT> > colcnts(splits);
for(int i=0; i< splits; ++i)
colcnts[i].resize(n, 0);
if(nnz > 0 && csc != NULL)
{
for(IT i=0; i< csc->n; ++i)
{
for(IT j = csc->jc[i]; j< csc->jc[i+1]; ++j)
{
IT rowid = csc->ir[j]; // colid=i
IT owner = std::min(rowid / perpiece, static_cast<IT>(splits-1));
colrowpairs[owner].push_back(std::make_tuple(i, rowid - owner*perpiece, csc->num[i]));
++(colcnts[owner][i]);
++(nnzs[owner]);
}
}
}
delete csc; // claim memory
cscarr = new Csc<IT,NT>*[splits];
#ifdef _OPENMP
#pragma omp parallel for
#endif
for(int i=0; i< splits; ++i) // i iterates over splits
{
cscarr[i] = new Csc<IT,NT>(nnzs[i],n);
sort(colrowpairs[i].begin(), colrowpairs[i].end()); // sort w.r.t. columns first and rows second
cscarr[i]->jc[0] = 0;
std::partial_sum(colcnts[i].begin(), colcnts[i].end(), cscarr[i]->jc+1);
std::copy(cscarr[i]->jc, cscarr[i]->jc+n, colcnts[i].begin()); // reuse the colcnts as "current column pointers"
for(IT k=0; k<nnzs[i]; ++k) // k iterates over all nonzeros
{
IT cindex = std::get<0>(colrowpairs[i][k]);
IT rindex = std::get<1>(colrowpairs[i][k]);
NT value = std::get<2>(colrowpairs[i][k]);
IT curcptr = (colcnts[i][cindex])++; // fetch the pointer and post increment
cscarr[i]->ir[curcptr] = rindex;
cscarr[i]->num[curcptr] = value;
}
}
}
template<class IT, class NT>
void SpCCols<IT,NT>::PrintInfo() const
{
std::cout << "m: " << m ;
std::cout << ", n: " << n ;
std::cout << ", nnz: "<< nnz ;
if(splits > 0)
{
std::cout << ", local splits: " << splits << std::endl;
#ifdef _OPENMP
if(omp_get_thread_num() == 0)
{
SubPrintInfo(cscarr[0]);
}
#endif
}
else
{
std::cout << std::endl;
SubPrintInfo(csc);
}
}
template <class IT, class NT>
template <typename UnaryOperation, typename GlobalIT>
SpCCols<IT, NT> *
SpCCols<IT, NT>::PruneI (UnaryOperation unary_op,
bool inPlace,
GlobalIT rowOffset,
GlobalIT colOffset
)
{
SpCCols<IT, NT> *retcols = NULL;
if (nnz > 0)
{
Csc<IT, NT> *ret = csc->PruneI(unary_op, inPlace, rowOffset, colOffset);
if (inPlace)
{
nnz = csc->nz;
if (nnz == 0)
{
delete csc;
csc = NULL;
}
retcols = NULL;
}
else
{
// wrap the pruned csc into a new SpCCols
retcols = new SpCCols<IT, NT>();
retcols->csc = ret;
retcols->nnz = retcols->csc->nz;
retcols->n = n;
retcols->m = m;
}
}
else
{
if (inPlace)
retcols = NULL;
else
{
retcols = new SpCCols<IT, NT>();
retcols->csc = NULL;
retcols->nnz = 0;
retcols->n = 0;
retcols->m = 0;
}
}
return retcols;
}
template <class IT, class NT>
std::vector<IT>
SpCCols<IT, NT>::GetEssentials (void) const
{
std::vector<IT> essentials(esscount);
essentials[0] = nnz;
essentials[1] = m;
essentials[2] = n;
return essentials;
}
template <class IT, class NT>
void
SpCCols<IT, NT>::CreateImpl (const std::vector<IT> &essentials)
{
assert(essentials.size() == esscount);
nnz = essentials[0];
m = essentials[1];
n = essentials[2];
if (nnz > 0)
csc = new Csc<IT, NT>(nnz, n);
else
csc = NULL;
}
template <class IT, class NT>
void
SpCCols<IT, NT>::CreateImpl (IT size,
IT nRow,
IT nCol,
std::tuple<IT, IT, NT> *mytuples)
{
SpTuples<IT,NT> tuples(size, nRow, nCol, mytuples);
tuples.SortColBased();
SpCCols<IT, NT> tmp(tuples, false);
*this = tmp;
}
template <class IT, class NT>
Arr<IT, NT>
SpCCols<IT, NT>::GetArrays (void) const
{
Arr<IT, NT> arr(2, 1);
if (nnz > 0)
{
arr.indarrs[0] = LocArr<IT, IT>(csc->jc, csc->n + 1);
arr.indarrs[1] = LocArr<IT, IT>(csc->ir, csc->nz);
arr.numarrs[0] = LocArr<NT, IT>(csc->num, csc->nz);
}
else
{
arr.indarrs[0] = LocArr<IT, IT>(NULL, 0);
arr.indarrs[1] = LocArr<IT, IT>(NULL, 0);
arr.numarrs[0] = LocArr<NT, IT>(NULL, 0);
}
return arr;
}
template <class IT, class NT>
void
SpCCols<IT, NT>::Transpose (void)
{
if (nnz > 0)
{
SpTuples<IT, NT> tuples(*this);
tuples.SortRowBased();
*this = SpCCols<IT, NT>(tuples, true);
}
else
*this = SpCCols<IT, NT>(0, n, m);
}
template <class IT, class NT>
SpCCols<IT, NT>
SpCCols<IT, NT>::TransposeConst (void) const
{
SpTuples<IT, NT> tuples(*this);
tuples.SortRowBased();
return SpCCols<IT, NT>(tuples, true);
}
template <class IT, class NT>
SpCCols<IT, NT> *
SpCCols<IT, NT>::TransposeConstPtr (void) const
{
SpTuples<IT, NT> tuples(*this);
tuples.SortRowBased();
return new SpCCols<IT, NT>(tuples, true);
}
template <class IT, class NT>
void
SpCCols<IT, NT>::Split (SpCCols<IT, NT> &partA,
SpCCols<IT, NT> &partB
)
{
IT cut = n/2;
if (cut == 0)
{
std::cout<< "Matrix is too small to be splitted" << std::endl;
return;
}
Csc<IT, NT> *Acsc = NULL;
Csc<IT, NT> *Bcsc = NULL;
if (nnz != 0)
csc->Split(Acsc, Bcsc, cut);
partA = SpCCols<IT, NT>(m, cut, Acsc);
partB = SpCCols<IT, NT>(m, n - cut, Bcsc);
*this = SpCCols<IT, NT>();
}
template <class IT, class NT>
void
SpCCols<IT, NT>::Merge (SpCCols<IT, NT> &partA,
SpCCols<IT, NT> &partB
)
{
assert (partA.m == partB.m);
Csc<IT, NT> *Ccsc = new Csc<IT, NT>();
if (partA.nnz == 0 && partB.nnz == 0)
Ccsc = NULL;
else if (partA.nnz == 0)
{
Ccsc = new Csc<IT, NT>(partB.nnz, partA.n + partB.n);
std::fill(Ccsc->jc, Ccsc->jc + partA.n, 0);
std::copy(partB.csc->jc, partB.csc->jc + partB.n + 1,
Ccsc->jc + partA.n);
std::copy(partB.csc->ir, partB.csc->ir + partB.nnz, Ccsc->ir);
std::copy(partB.csc->num, partB.csc->num + partB.nnz, Ccsc->num);
// Ccsc = new Csc<IT,NT>(*(partB.csc));
}
else if (partB.nnz == 0)
{
Ccsc = new Csc<IT, NT>(partA.nnz, partA.n + partB.n);
std::copy(partA.csc->jc, partA.csc->jc + partA.n + 1, Ccsc->jc);
std::fill(Ccsc->jc + partA.n + 1, Ccsc->jc + partA.n + partB.n + 1,
partA.csc->jc[partA.n]);
std::copy(partA.csc->ir, partA.csc->ir + partA.nnz, Ccsc->ir);
std::copy(partA.csc->num, partA.csc->num + partA.nnz, Ccsc->num);
// Ccsc = new Csc<IT,NT>(*(partA.csc));
}
else
Ccsc->Merge(partA.csc, partB.csc, partA.n); // 3rd param not used
*this = SpCCols<IT, NT>(partA.m, partA.n + partB.n, Ccsc);
partA = SpCCols<IT, NT>();
partB = SpCCols<IT, NT>();
}
template <class IT, class NT>
std::ofstream &
SpCCols<IT, NT>::put (std::ofstream &outfile) const
{
if (nnz == 0)
{
outfile << "Matrix doesn't have any nonzeros" << std::endl;
return outfile;
}
SpTuples<IT, NT> tuples(*this);
outfile << tuples << std::endl;
return outfile;
}
/****************************************************************************/
/************************* PRIVATE MEMBER FUNCTIONS *************************/
/****************************************************************************/
template <class IT, class NT>
SpCCols<IT, NT>::SpCCols (IT nRow,
IT nCol,
Csc<IT, NT> *mycsc
) :
csc(mycsc), m(nRow), n(nCol), splits(0)
{
if (mycsc == NULL)
nnz = 0;
else
nnz = mycsc->nz;
}
template <class IT, class NT>
void SpCCols<IT,NT>::SubPrintInfo(Csc<IT,NT> * mycsc) const
{
#ifdef _OPENMP
std::cout << "Printing for thread " << omp_get_thread_num() << std::endl;
#endif
if(m < PRINT_LIMIT && n < PRINT_LIMIT) // small enough to print
{
NT ** A = SpHelper::allocate2D<NT>(m,n);
for(IT i=0; i< m; ++i)
for(IT j=0; j<n; ++j)
A[i][j] = NT();
if(mycsc != NULL)
{
for(IT i=0; i< n; ++i)
{
for(IT j = mycsc->jc[i]; j< mycsc->jc[i+1]; ++j)
{
IT rowid = mycsc->ir[j];
A[rowid][i] = mycsc->num[j];
}
}
}
for(IT i=0; i< m; ++i)
{
for(IT j=0; j<n; ++j)
{
std::cout << std::setiosflags(std::ios::fixed) << std::setprecision(2) << A[i][j];
std::cout << " ";
}
std::cout << std::endl;
}
SpHelper::deallocate2D(A,m);
}
}
template <class IT, class NT>
inline void SpCCols<IT,NT>::CopyCsc(Csc<IT,NT> * source)
{
// source csc will be NULL if number of nonzeros is zero
if(source != NULL)
csc = new Csc<IT,NT>(*source);
else
csc = NULL;
}
}
|