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
|
/****************************************************************/
/* 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.
*/
#ifndef _SP_DCCOLS_H
#define _SP_DCCOLS_H
#include <iostream>
#include <fstream>
#include <cmath>
#include "SpMat.h" // Best to include the base class first
#include "SpHelper.h"
#include "StackEntry.h"
#include "dcsc.h"
#include "Isect.h"
#include "Semirings.h"
#include "MemoryPool.h"
#include "LocArr.h"
#include "Friends.h"
#include "CombBLAS.h"
#include "FullyDist.h"
namespace combblas {
template <class IT, class NT>
class SpDCCols: public SpMat<IT, NT, SpDCCols<IT, NT> >
{
public:
typedef IT LocalIT;
typedef NT LocalNT;
// Constructors :
SpDCCols ();
SpDCCols (IT size, IT nRow, IT nCol, IT nzc);
SpDCCols (const SpTuples<IT,NT> & rhs, bool transpose);
SpDCCols (IT nRow, IT nCol, IT nnz1, const std::tuple<IT, IT, NT> * rhs, bool transpose);
SpDCCols (const SpDCCols<IT,NT> & rhs); // Actual copy constructor
~SpDCCols();
template <typename NNT> operator SpDCCols<IT,NNT> () const; //!< NNT: New numeric type
template <typename NIT, typename NNT> operator SpDCCols<NIT,NNT> () const; //!< NNT: New numeric type, NIT: New index type
// Member Functions and Operators:
SpDCCols<IT,NT> & operator= (const SpDCCols<IT, NT> & rhs);
SpDCCols<IT,NT> & operator+= (const SpDCCols<IT, NT> & rhs);
SpDCCols<IT,NT> operator() (IT ri, IT ci) const;
SpDCCols<IT,NT> operator() (const std::vector<IT> & ri, const std::vector<IT> & ci) const;
bool operator== (const SpDCCols<IT, NT> & rhs) const
{
if(rhs.nnz == 0 && nnz == 0)
return true;
if(nnz != rhs.nnz || m != rhs.m || n != rhs.n)
return false;
return ((*dcsc) == (*(rhs.dcsc)));
}
class SpColIter //! Iterate over (sparse) columns of the sparse matrix
{
public:
class NzIter //! Iterate over the nonzeros of the sparse column
{
public:
NzIter(IT * ir = NULL, NT * num = NULL) : rid(ir), val(num) {}
bool operator==(const NzIter & other)
{
return(rid == other.rid); // compare pointers
}
bool operator!=(const NzIter & other)
{
return(rid != other.rid);
}
bool operator<(const NzIter & other)
{
return(rid < other.rid);
}
NzIter operator+(IT inc)
{
rid+=inc;
val+=inc;
return(*this);
}
NzIter operator-(IT inc)
{
rid-=inc;
val-=inc;
return(*this);
}
NzIter & operator++() // prefix operator
{
++rid;
++val;
return(*this);
}
NzIter operator++(int) // postfix operator
{
NzIter tmp(*this);
++(*this);
return(tmp);
}
IT rowid() const //!< Return the "local" rowid of the current nonzero entry.
{
return (*rid);
}
NT & value() //!< value is returned by reference for possible updates
{
return (*val);
}
private:
IT * rid;
NT * val;
};
SpColIter(IT * cp = NULL, IT * jc = NULL) : cptr(cp), cid(jc) {}
bool operator==(const SpColIter& other)
{
return(cptr == other.cptr); // compare pointers
}
bool operator!=(const SpColIter& other)
{
return(cptr != other.cptr);
}
SpColIter& operator++() // prefix operator
{
++cptr;
++cid;
return(*this);
}
SpColIter operator++(int) // postfix operator (common)
{
SpColIter tmp(*this);
++(*this);
return(tmp);
}
SpColIter operator+(IT inc)
{
cptr+=inc;
cid+=inc;
return(*this);
}
SpColIter operator-(IT inc)
{
cptr-=inc;
cid-=inc;
return(*this);
}
IT colid() const //!< Return the "local" colid of the current column.
{
return (*cid);
}
IT colptr() const
{
return (*cptr);
}
IT colptrnext() const
{
return (*(cptr+1));
}
IT nnz() const
{
return (colptrnext() - colptr());
}
private:
IT * cptr;
IT * cid;
};
SpColIter begcol()
{
if( nnz > 0 )
return SpColIter(dcsc->cp, dcsc->jc);
else
return SpColIter(NULL, NULL);
}
SpColIter begcol(int i) // multithreaded version
{
if( dcscarr[i] )
return SpColIter(dcscarr[i]->cp, dcscarr[i]->jc);
else
return SpColIter(NULL, NULL);
}
SpColIter endcol()
{
if( nnz > 0 )
return SpColIter(dcsc->cp + dcsc->nzc, NULL);
else
return SpColIter(NULL, NULL);
}
SpColIter endcol(int i) //multithreaded version
{
if( dcscarr[i] )
return SpColIter(dcscarr[i]->cp + dcscarr[i]->nzc, NULL);
else
return SpColIter(NULL, NULL);
}
typename SpColIter::NzIter begnz(const SpColIter & ccol) //!< Return the beginning iterator for the nonzeros of the current column
{
return typename SpColIter::NzIter( dcsc->ir + ccol.colptr(), dcsc->numx + ccol.colptr() );
}
typename SpColIter::NzIter endnz(const SpColIter & ccol) //!< Return the ending iterator for the nonzeros of the current column
{
return typename SpColIter::NzIter( dcsc->ir + ccol.colptrnext(), NULL );
}
typename SpColIter::NzIter begnz(const SpColIter & ccol, int i) //!< multithreaded version
{
return typename SpColIter::NzIter( dcscarr[i]->ir + ccol.colptr(), dcscarr[i]->numx + ccol.colptr() );
}
typename SpColIter::NzIter endnz(const SpColIter & ccol, int i) //!< multithreaded version
{
return typename SpColIter::NzIter( dcscarr[i]->ir + ccol.colptrnext(), NULL );
}
template <typename _UnaryOperation>
void Apply(_UnaryOperation __unary_op)
{
if(nnz > 0)
dcsc->Apply(__unary_op);
}
template <typename _UnaryOperation, typename GlobalIT>
SpDCCols<IT,NT>* PruneI(_UnaryOperation __unary_op, bool inPlace, GlobalIT rowOffset, GlobalIT colOffset);
template <typename _UnaryOperation>
SpDCCols<IT,NT>* Prune(_UnaryOperation __unary_op, bool inPlace);
template <typename _BinaryOperation>
SpDCCols<IT,NT>* PruneColumn(NT* pvals, _BinaryOperation __binary_op, bool inPlace);
template <typename _BinaryOperation>
SpDCCols<IT,NT>* PruneColumn(IT* pinds, NT* pvals, _BinaryOperation __binary_op, bool inPlace);
void PruneColumnByIndex(const std::vector<IT>& ci);
template <typename _BinaryOperation>
void UpdateDense(NT ** array, _BinaryOperation __binary_op) const
{
if(nnz > 0 && dcsc != NULL)
dcsc->UpdateDense(array, __binary_op);
}
void EWiseScale(NT ** scaler, IT m_scaler, IT n_scaler);
void EWiseMult (const SpDCCols<IT,NT> & rhs, bool exclude);
void SetDifference (const SpDCCols<IT,NT> & rhs);
void Transpose(); //!< Mutator version, replaces the calling object
SpDCCols<IT,NT> TransposeConst() const; //!< Const version, doesn't touch the existing object
SpDCCols<IT,NT> * TransposeConstPtr() const;
void RowSplit(int numsplits)
{
BooleanRowSplit(*this, numsplits); // only works with boolean arrays
}
void ColSplit(int parts, std::vector< SpDCCols<IT,NT> > & matrices); //!< \attention Destroys calling object (*this)
void ColSplit(int parts, std::vector< SpDCCols<IT,NT>* > & matrices); //!< \attention Destroys calling object (*this)
void ColSplit(std::vector<IT> & cutSizes, std::vector< SpDCCols<IT,NT> > & matrices); //!< \attention Destroys calling object (*this)
void ColSplit(std::vector<IT> & cutSizes, std::vector< SpDCCols<IT,NT>* > & matrices); //!< \attention Destroys calling object (*this)
void ColConcatenate(std::vector< SpDCCols<IT,NT> > & matrices); //!< \attention Destroys its parameters (matrices)
void ColConcatenate(std::vector< SpDCCols<IT,NT>* > & matrices); //!< \attention Destroys its parameters (matrices)
void Split(SpDCCols<IT,NT> & partA, SpDCCols<IT,NT> & partB); //!< \attention Destroys calling object (*this)
void Merge(SpDCCols<IT,NT> & partA, SpDCCols<IT,NT> & partB); //!< \attention Destroys its parameters (partA & partB)
void CreateImpl(const std::vector<IT> & essentials);
void CreateImpl(IT size, IT nRow, IT nCol, std::tuple<IT, IT, NT> * mytuples);
void CreateImpl(IT * _cp, IT * _jc, IT * _ir, NT * _numx, IT _nz, IT _nzc, IT _m, IT _n);
Arr<IT,NT> GetArrays() const;
std::vector<IT> GetEssentials() const;
const static IT esscount;
bool isZero() const { return (nnz == 0); }
IT getnrow() const { return m; }
IT getncol() const { return n; }
IT getnnz() const { return nnz; }
IT getnzc() const { return (nnz == 0) ? 0: dcsc->nzc; }
int getnsplit() const { return splits; }
std::ofstream& put(std::ofstream & outfile) const;
std::ifstream& get(std::ifstream & infile);
void PrintInfo() const;
void PrintInfo(std::ofstream & out) const;
template <typename SR>
int PlusEq_AtXBt(const SpDCCols<IT,NT> & A, const SpDCCols<IT,NT> & B);
template <typename SR>
int PlusEq_AtXBn(const SpDCCols<IT,NT> & A, const SpDCCols<IT,NT> & B);
template <typename SR>
int PlusEq_AnXBt(const SpDCCols<IT,NT> & A, const SpDCCols<IT,NT> & B);
template <typename SR>
int PlusEq_AnXBn(const SpDCCols<IT,NT> & A, const SpDCCols<IT,NT> & B);
Dcsc<IT, NT> * GetDCSC() const // only for single threaded matrices
{
return dcsc;
}
Dcsc<IT, NT> * GetDCSC(int i) const // only for split (multithreaded) matrices
{
return dcscarr[i];
}
auto GetInternal() const { return GetDCSC(); }
auto GetInternal(int i) const { return GetDCSC(i); }
private:
void CopyDcsc(Dcsc<IT,NT> * source);
SpDCCols<IT,NT> ColIndex(const std::vector<IT> & ci) const; //!< col indexing without multiplication
template <typename SR, typename NTR>
SpDCCols< IT, typename promote_trait<NT,NTR>::T_promote > OrdOutProdMult(const SpDCCols<IT,NTR> & rhs) const;
template <typename SR, typename NTR>
SpDCCols< IT, typename promote_trait<NT,NTR>::T_promote > OrdColByCol(const SpDCCols<IT,NTR> & rhs) const;
SpDCCols (IT size, IT nRow, IT nCol, const std::vector<IT> & indices, bool isRow); // Constructor for indexing
SpDCCols (IT nRow, IT nCol, Dcsc<IT,NT> * mydcsc); // Constructor for multiplication
// Anonymous union
union {
Dcsc<IT, NT> * dcsc;
Dcsc<IT, NT> ** dcscarr;
};
IT m;
IT n;
IT nnz;
int splits; // for multithreading
template <class IU, class NU>
friend class SpDCCols; // Let other template instantiations (of the same class) access private members
template <class IU, class NU>
friend class SpTuples;
// AL: removed this because it appears illegal and causes this compiler warning:
// warning: dependent nested name specifier 'SpDCCols<IU, NU>::' for friend class declaration is not supported; turning off access control for 'SpDCCols'
//template <class IU, class NU>
//friend class SpDCCols<IU, NU>::SpColIter;
template<typename IU>
friend void BooleanRowSplit(SpDCCols<IU, bool> & A, int numsplits);
template<typename IU, typename NU1, typename NU2>
friend SpDCCols<IU, typename promote_trait<NU1,NU2>::T_promote > EWiseMult (const SpDCCols<IU,NU1> & A, const SpDCCols<IU,NU2> & B, bool exclude);
template<typename N_promote, typename IU, typename NU1, typename NU2, typename _BinaryOperation>
friend SpDCCols<IU, N_promote > EWiseApply (const SpDCCols<IU,NU1> & A, const SpDCCols<IU,NU2> & B, _BinaryOperation __binary_op, bool notB, const NU2& defaultBVal);
template <typename RETT, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
friend SpDCCols<IU,RETT> EWiseApply (const SpDCCols<IU,NU1> & A, const SpDCCols<IU,NU2> & B, _BinaryOperation __binary_op, _BinaryPredicate do_op, bool allowANulls, bool allowBNulls, const NU1& ANullVal, const NU2& BNullVal, const bool allowIntersect);
template<class SR, class NUO, class IU, class NU1, class NU2>
friend SpTuples<IU, NUO> * Tuples_AnXBn
(const SpDCCols<IU, NU1> & A, const SpDCCols<IU, NU2> & B, bool clearA, bool clearB);
template<class SR, class NUO, class IU, class NU1, class NU2>
friend SpTuples<IU, NUO> * Tuples_AnXBt
(const SpDCCols<IU, NU1> & A, const SpDCCols<IU, NU2> & B, bool clearA, bool clearB);
template<class SR, class NUO, class IU, class NU1, class NU2>
friend SpTuples<IU, NUO> * Tuples_AtXBn
(const SpDCCols<IU, NU1> & A, const SpDCCols<IU, NU2> & B, bool clearA, bool clearB);
template<class SR, class NUO, class IU, class NU1, class NU2>
friend SpTuples<IU, NUO> * Tuples_AtXBt
(const SpDCCols<IU, NU1> & A, const SpDCCols<IU, NU2> & B, bool clearA, bool clearB);
template <typename SR, typename IU, typename NU, typename RHS, typename LHS>
friend void dcsc_gespmv (const SpDCCols<IU, NU> & A, const RHS * x, LHS * y);
template <typename SR, typename IU, typename NU, typename RHS, typename LHS>
friend void dcsc_gespmv_threaded (const SpDCCols<IU, NU> & A, const RHS * x, LHS * y);
template <typename SR, typename IU, typename NU, typename RHS, typename LHS>
friend void dcsc_gespmv_threaded_nosplit (const SpDCCols<IU, NU> & A, const RHS * x, LHS * y);
template <typename SR, typename IU, typename NUM, typename DER, typename IVT, typename OVT>
friend int generic_gespmv_threaded (const SpMat<IU,NUM,DER> & A, const int32_t * indx, const IVT * numx, int32_t nnzx,
int32_t * & sendindbuf, OVT * & sendnumbuf, int * & sdispls, int p_c);
};
// At this point, complete type of of SpDCCols is known, safe to declare these specialization (but macros won't work as they are preprocessed)
// General case #1: When both NT is the same
template <class IT, class NT> struct promote_trait< SpDCCols<IT,NT> , SpDCCols<IT,NT> >
{
typedef SpDCCols<IT,NT> T_promote;
};
// General case #2: First is boolean the second is anything except boolean (to prevent ambiguity)
template <class IT, class NT> struct promote_trait< SpDCCols<IT,bool> , SpDCCols<IT,NT>, typename combblas::disable_if< combblas::is_boolean<NT>::value >::type >
{
typedef SpDCCols<IT,NT> T_promote;
};
// General case #3: Second is boolean the first is anything except boolean (to prevent ambiguity)
template <class IT, class NT> struct promote_trait< SpDCCols<IT,NT> , SpDCCols<IT,bool>, typename combblas::disable_if< combblas::is_boolean<NT>::value >::type >
{
typedef SpDCCols<IT,NT> T_promote;
};
template <class IT> struct promote_trait< SpDCCols<IT,int> , SpDCCols<IT,float> >
{
typedef SpDCCols<IT,float> T_promote;
};
template <class IT> struct promote_trait< SpDCCols<IT,float> , SpDCCols<IT,int> >
{
typedef SpDCCols<IT,float> T_promote;
};
template <class IT> struct promote_trait< SpDCCols<IT,int> , SpDCCols<IT,double> >
{
typedef SpDCCols<IT,double> T_promote;
};
template <class IT> struct promote_trait< SpDCCols<IT,double> , SpDCCols<IT,int> >
{
typedef SpDCCols<IT,double> T_promote;
};
// Below are necessary constructs to be able to define a SpMat<NT,IT> where
// all we know is DER (say SpDCCols<int, double>) and NT,IT
// in other words, we infer the templated SpDCCols<> type
// This is not a type conversion from an existing object,
// but a type inference for the newly created object
// NIT: New IT, NNT: New NT
template <class DER, class NIT, class NNT>
struct create_trait
{
// none
};
// Capture everything of the form SpDCCols<OIT, ONT>
// it may come as a surprise that the partial specializations can
// involve more template parameters than the primary template
template <class NIT, class NNT, class OIT, class ONT>
struct create_trait< SpDCCols<OIT, ONT> , NIT, NNT >
{
typedef SpDCCols<NIT,NNT> T_inferred;
};
}
#include "SpDCCols.cpp"
#endif
|