COMBINATORIAL_BLAS  1.6
FullyDistSpVec.h
Go to the documentation of this file.
1 /****************************************************************/
2 /* Parallel Combinatorial BLAS Library (for Graph Computations) */
3 /* version 1.6 -------------------------------------------------*/
4 /* date: 6/15/2017 ---------------------------------------------*/
5 /* authors: Ariful Azad, Aydin Buluc --------------------------*/
6 /****************************************************************/
7 /*
8  Copyright (c) 2010-2017, The Regents of the University of California
9 
10  Permission is hereby granted, free of charge, to any person obtaining a copy
11  of this software and associated documentation files (the "Software"), to deal
12  in the Software without restriction, including without limitation the rights
13  to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
14  copies of the Software, and to permit persons to whom the Software is
15  furnished to do so, subject to the following conditions:
16 
17  The above copyright notice and this permission notice shall be included in
18  all copies or substantial portions of the Software.
19 
20  THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
21  IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
22  FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
23  AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
24  LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
25  OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
26  THE SOFTWARE.
27  */
28 
29 
30 #ifndef _FULLY_DIST_SP_VEC_H_
31 #define _FULLY_DIST_SP_VEC_H_
32 
33 #include <iostream>
34 #include <vector>
35 #include <utility>
36 #include "CommGrid.h"
37 #include "promote.h"
38 #include "SpParMat.h"
39 #include "FullyDist.h"
40 #include "Exception.h"
41 #include "OptBuf.h"
42 #include "CombBLAS.h"
43 
44 namespace combblas {
45 
46 template <class IT, class NT, class DER>
47 class SpParMat;
48 
49 template <class IT>
50 class DistEdgeList;
51 
52 template <class IU, class NU>
53 class FullyDistVec;
54 
55 template <class IU, class NU>
56 class SparseVectorLocalIterator;
57 
72 template <class IT, class NT>
73 class FullyDistSpVec: public FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>
74 {
75 public:
76  FullyDistSpVec ( );
77  explicit FullyDistSpVec ( IT glen );
78  FullyDistSpVec ( std::shared_ptr<CommGrid> grid);
79  FullyDistSpVec ( std::shared_ptr<CommGrid> grid, IT glen);
80 
81  template <typename _UnaryOperation>
82  FullyDistSpVec (const FullyDistVec<IT,NT> & rhs, _UnaryOperation unop);
83  FullyDistSpVec (const FullyDistVec<IT,NT> & rhs); // Conversion copy-constructor
84  FullyDistSpVec (IT globalsize, const FullyDistVec<IT,IT> & inds, const FullyDistVec<IT,NT> & vals, bool SumDuplicates = false);
85  FullyDistSpVec (std::shared_ptr<CommGrid> grid, IT globallen, const std::vector<IT>& indvec, const std::vector<NT> & numvec, bool SumDuplicates = false, bool sorted=false);
86 
87  IT NnzUntil() const;
88 
89  FullyDistSpVec<IT,NT> Invert (IT globallen);
90  template <typename _BinaryOperationIdx, typename _BinaryOperationVal, typename _BinaryOperationDuplicate>
91  FullyDistSpVec<IT,NT> Invert (IT globallen, _BinaryOperationIdx __binopIdx, _BinaryOperationVal __binopVal, _BinaryOperationDuplicate __binopDuplicate);
92  template <typename _BinaryOperationIdx, typename _BinaryOperationVal>
93  FullyDistSpVec<IT,NT> InvertRMA (IT globallen, _BinaryOperationIdx __binopIdx, _BinaryOperationVal __binopVal);
94 
95 
96  template <typename NT1, typename _UnaryOperation>
97  void Select (const FullyDistVec<IT,NT1> & denseVec, _UnaryOperation unop);
98  template <typename _UnaryOperation>
99  void FilterByVal (FullyDistSpVec<IT,IT> Selector, _UnaryOperation __unop, bool filterByIndex);
100  template <typename NT1>
101  void Setminus (const FullyDistSpVec<IT,NT1> & other);
102 
103  //template <typename NT1, typename _UnaryOperation>
104  //void Set (FullyDistSpVec<IT,NT1> Selector, _UnaryOperation __unop);
105 
106  template <typename NT1, typename _UnaryOperation, typename _BinaryOperation>
107  void SelectApply (const FullyDistVec<IT,NT1> & denseVec, _UnaryOperation __unop, _BinaryOperation __binop);
108 
109 
110 
113  void stealFrom(FullyDistSpVec<IT,NT> & victim);
114  FullyDistSpVec<IT,NT> & operator=(const FullyDistSpVec< IT,NT > & rhs);
115  FullyDistSpVec<IT,NT> & operator=(const FullyDistVec< IT,NT > & rhs); // convert from dense
116  FullyDistSpVec<IT,NT> & operator=(NT fixedval) // assign fixed value
117  {
118 #ifdef _OPENMP
119 #pragma omp parallel for
120 #endif
121  for(size_t i=0; i < ind.size(); ++i)
122  num[i] = fixedval;
123  return *this;
124  }
127 
128  class ScalarReadSaveHandler
129  {
130  public:
131  NT getNoNum(IT index) { return static_cast<NT>(1); }
132 
133  template <typename c, typename t>
134  NT read(std::basic_istream<c,t>& is, IT index)
135  {
136  NT v;
137  is >> v;
138  return v;
139  }
140 
141  template <typename c, typename t>
142  void save(std::basic_ostream<c,t>& os, const NT& v, IT index)
143  {
144  os << v;
145  }
146  };
147 
148  template <class HANDLER>
149  void ParallelWrite(const std::string & filename, bool onebased, HANDLER handler, bool includeindices = true, bool includeheader = false);
150  void ParallelWrite(const std::string & filename, bool onebased, bool includeindices = true) { ParallelWrite(filename, onebased, ScalarReadSaveHandler(), includeindices); };
151 
152 
153  template <typename _BinaryOperation>
154  void ParallelRead (const std::string & filename, bool onebased, _BinaryOperation BinOp);
155 
156 
158  template <class HANDLER>
159  std::ifstream& ReadDistribute (std::ifstream& infile, int master, HANDLER handler);
160  std::ifstream& ReadDistribute (std::ifstream& infile, int master) { return ReadDistribute(infile, master, ScalarReadSaveHandler()); }
161 
162  template <class HANDLER>
163  void SaveGathered(std::ofstream& outfile, int master, HANDLER handler, bool printProcSplits = false);
164  void SaveGathered(std::ofstream& outfile, int master) { SaveGathered(outfile, master, ScalarReadSaveHandler()); }
165 
166 
167  template <typename NNT> operator FullyDistSpVec< IT,NNT > () const
168  {
169  FullyDistSpVec<IT,NNT> CVT(commGrid);
170  CVT.ind = std::vector<IT>(ind.begin(), ind.end());
171  CVT.num = std::vector<NNT>(num.begin(), num.end());
172  CVT.glen = glen;
173  return CVT;
174  }
175 
176  bool operator==(const FullyDistSpVec<IT,NT> & rhs) const
177  {
178  FullyDistVec<IT,NT> v = *this;
179  FullyDistVec<IT,NT> w = rhs;
180  return (v == w);
181  }
182 
183  void PrintInfo(std::string vecname) const;
184  void iota(IT globalsize, NT first);
185  void nziota(NT first);
187  void SetElement (IT indx, NT numx); // element-wise assignment
188  void DelElement (IT indx); // element-wise deletion
189  NT operator[](IT indx);
190  bool WasFound() const { return wasFound; }
191 
194 
195 #if __cplusplus > 199711L
196  template <typename _BinaryOperation = minimum<NT> >
197  FullyDistSpVec<IT, NT> Uniq(_BinaryOperation __binary_op = _BinaryOperation(), MPI_Op mympiop = MPI_MIN);
198 #else
199  template <typename _BinaryOperation >
200  FullyDistSpVec<IT, NT> Uniq(_BinaryOperation __binary_op, MPI_Op mympiop);
201 #endif
202 
203 
204  IT getlocnnz() const
205  {
206  return ind.size();
207  }
208  IT getnnz() const
209  {
210  IT totnnz = 0;
211  IT locnnz = ind.size();
212  MPI_Allreduce( &locnnz, &totnnz, 1, MPIType<IT>(), MPI_SUM, commGrid->GetWorld());
213  return totnnz;
214  }
221 
222  void setNumToInd()
223  {
224  IT offset = LengthUntil();
225  IT spsize = ind.size();
226  #ifdef _OPENMP
227  #pragma omp parallel for
228  #endif
229  for(IT i=0; i< spsize; ++i)
230  num[i] = ind[i] + offset;
231  }
232 
233  template <typename _Predicate>
234  IT Count(_Predicate pred) const;
235 
236  template <typename _UnaryOperation>
237  void Apply(_UnaryOperation __unary_op)
238  {
239  //transform(num.begin(), num.end(), num.begin(), __unary_op);
240  IT spsize = num.size();
241 #ifdef _OPENMP
242 #pragma omp parallel for
243 #endif
244  for(IT i=0; i < spsize; ++i)
245  num[i] = __unary_op(num[i]);
246  }
247 
248  template <typename _BinaryOperation>
249  void ApplyInd(_BinaryOperation __binary_op)
250  {
251  IT offset = LengthUntil();
252  IT spsize = ind.size();
253  #ifdef _OPENMP
254  #pragma omp parallel for
255  #endif
256  for(IT i=0; i < spsize; ++i)
257  num[i] = __binary_op(num[i], ind[i] + offset);
258  }
259 
260 
261 
262  template <typename _BinaryOperation>
263  NT Reduce(_BinaryOperation __binary_op, NT init) const;
264 
265  template <typename OUT, typename _BinaryOperation, typename _UnaryOperation>
266  OUT Reduce(_BinaryOperation __binary_op, OUT default_val, _UnaryOperation __unary_op) const;
267 
268  void DebugPrint();
269  std::shared_ptr<CommGrid> getcommgrid() const { return commGrid; }
270 
271  void Reset();
272  NT GetLocalElement(IT indx);
273  void BulkSet(IT inds[], int count);
274  std::vector<IT> GetLocalInd (){std::vector<IT> rind = ind; return rind;};
275  std::vector<NT> GetLocalNum (){std::vector<NT> rnum = num; return rnum;};
276 
277  template <typename _Predicate>
278  FullyDistVec<IT,IT> FindInds(_Predicate pred) const;
279  template <typename _Predicate>
280  FullyDistVec<IT,NT> FindVals(_Predicate pred) const;
281 
282 
283 protected:
286 
287 private:
288  std::vector< IT > ind; // ind.size() give the number of nonzeros
289  std::vector< NT > num;
290  bool wasFound; // true if the last GetElement operation returned an actual value
291 
292 
293  template <typename _BinaryOperation>
294  void SparseCommon(std::vector< std::vector < std::pair<IT,NT> > > & data, _BinaryOperation BinOp);
295 
296 
297 #if __cplusplus > 199711L
298  template <typename _BinaryOperation = minimum<NT> >
299  FullyDistSpVec<IT, NT> UniqAll2All(_BinaryOperation __binary_op = _BinaryOperation(), MPI_Op mympiop = MPI_MIN);
300 #else
301  template <typename _BinaryOperation >
302  FullyDistSpVec<IT, NT> UniqAll2All(_BinaryOperation __binary_op, MPI_Op mympiop);
303 #endif
304 
305 
306  template <class IU, class NU>
307  friend class FullyDistSpVec;
308 
309  template <class IU, class NU>
310  friend class FullyDistVec;
311 
312  template <class IU, class NU, class UDER>
313  friend class SpParMat;
314 
315  template <class IU, class NU>
317 
318  template <typename SR, typename IU, typename NUM, typename NUV, typename UDER>
320  SpMV (const SpParMat<IU,NUM,UDER> & A, const FullyDistSpVec<IU,NUV> & x );
321 
322  template <typename SR, typename IU, typename NUM, typename UDER>
324  SpMV (const SpParMat<IU,NUM,UDER> & A, const FullyDistSpVec<IU,IU> & x, bool indexisvalue);
325 
326  template <typename VT, typename IU, typename UDER> // NoSR version (in BFSFriends.h)
328 
329  template <typename SR, typename IVT, typename OVT, typename IU, typename NUM, typename UDER>
330  friend void SpMV (const SpParMat<IU,NUM,UDER> & A, const FullyDistSpVec<IU,IVT> & x, FullyDistSpVec<IU,OVT> & y,bool indexisvalue, OptBuf<int32_t, OVT > & optbuf);
331 
332  template <typename SR, typename IVT, typename OVT, typename IU, typename NUM, typename UDER>
333  friend void SpMV (const SpParMat<IU,NUM,UDER> & A, const FullyDistSpVec<IU,IVT> & x, FullyDistSpVec<IU,OVT> & y,bool indexisvalue, OptBuf<int32_t, OVT > & optbuf, PreAllocatedSPA<OVT> & SPA);
334 
335  template <typename IU, typename NU1, typename NU2>
337  EWiseMult (const FullyDistSpVec<IU,NU1> & V, const FullyDistVec<IU,NU2> & W , bool exclude, NU2 zero);
338 
339  template <typename RET, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
340  friend FullyDistSpVec<IU,RET>
341  EWiseApply (const FullyDistSpVec<IU,NU1> & V, const FullyDistVec<IU,NU2> & W , _BinaryOperation _binary_op, _BinaryPredicate _doOp, bool allowVNulls, NU1 Vzero, const bool useExtendedBinOp);
342 
343  template <typename RET, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
345  EWiseApply_threaded (const FullyDistSpVec<IU,NU1> & V, const FullyDistVec<IU,NU2> & W , _BinaryOperation _binary_op, _BinaryPredicate _doOp, bool allowVNulls, NU1 Vzero, const bool useExtendedBinOp);
346 
347 
348  template <typename RET, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
350  EWiseApply (const FullyDistSpVec<IU,NU1> & V, const FullyDistSpVec<IU,NU2> & W , _BinaryOperation _binary_op, _BinaryPredicate _doOp, bool allowVNulls, bool allowWNulls, NU1 Vzero, NU2 Wzero, const bool allowIntersect, const bool useExtendedBinOp);
351 
352  template <typename IU>
353  friend void RandPerm(FullyDistSpVec<IU,IU> & V); // called on an existing object, randomly permutes it
354 
355  template <typename IU>
356  friend void RenameVertices(DistEdgeList<IU> & DEL);
357 
359  // Ariful: I made this an internal function in ParFriends.h
360  //template <typename SR, typename IU, typename OVT>
361  //friend void MergeContributions(FullyDistSpVec<IU,OVT> & y, int * & recvcnt, int * & rdispls, int32_t * & recvindbuf, OVT * & recvnumbuf, int rowneighs);
362 
363  template <typename IU, typename VT>
364  friend void MergeContributions(FullyDistSpVec<IU,VT> & y, int * & recvcnt, int * & rdispls, int32_t * & recvindbuf, VT * & recvnumbuf, int rowneighs);
365 
366  template<typename IU, typename NV>
367  friend void TransposeVector(MPI_Comm & World, const FullyDistSpVec<IU,NV> & x, int32_t & trxlocnz, IU & lenuntil, int32_t * & trxinds, NV * & trxnums, bool indexisvalue);
368 
369  template <class IU, class NU, class DER, typename _UnaryOperation>
370  friend SpParMat<IU, bool, DER> PermMat1 (const FullyDistSpVec<IU,NU> & ri, const IU ncol, _UnaryOperation __unop);
371 };
372 
373 }
374 
375 #include "FullyDistSpVec.cpp"
376 
377 #endif
combblas::FullyDistSpVec::GetLocalInd
std::vector< IT > GetLocalInd()
Definition: FullyDistSpVec.h:274
combblas::FullyDistVec
Definition: FullyDistSpVec.h:53
combblas::FullyDistSpVec::NnzUntil
IT NnzUntil() const
Returns the number of nonzeros until this processor.
Definition: FullyDistSpVec.cpp:696
combblas::FullyDistSpVec::Select
void Select(const FullyDistVec< IT, NT1 > &denseVec, _UnaryOperation unop)
Definition: FullyDistSpVec.cpp:2358
combblas::FullyDistSpVec::Setminus
void Setminus(const FullyDistSpVec< IT, NT1 > &other)
Definition: FullyDistSpVec.cpp:2400
combblas::FullyDistSpVec::ScalarReadSaveHandler::getNoNum
NT getNoNum(IT index)
Definition: FullyDistSpVec.h:131
combblas::FullyDistSpVec::nziota
void nziota(NT first)
iota over existing nonzero entries
Definition: FullyDistSpVec.cpp:689
combblas
Definition: CCGrid.h:4
combblas::FullyDistSpVec::GetLocalElement
NT GetLocalElement(IT indx)
Definition: FullyDistSpVec.cpp:491
combblas::FullyDistSpVec::PrintInfo
void PrintInfo(std::string vecname) const
Definition: FullyDistSpVec.cpp:1686
combblas::FullyDistSpVec::setNumToInd
void setNumToInd()
Definition: FullyDistSpVec.h:222
combblas::FullyDistSpVec::ScalarReadSaveHandler::save
void save(std::basic_ostream< c, t > &os, const NT &v, IT index)
Definition: FullyDistSpVec.h:142
combblas::FullyDistSpVec::GetLocalNum
std::vector< NT > GetLocalNum()
Definition: FullyDistSpVec.h:275
FullyDist.h
combblas::FullyDistSpVec::FindInds
FullyDistVec< IT, IT > FindInds(_Predicate pred) const
Definition: FullyDistSpVec.cpp:394
combblas::FullyDistSpVec::SetElement
void SetElement(IT indx, NT numx)
Indexing is performed 0-based.
Definition: FullyDistSpVec.cpp:511
combblas::FullyDistSpVec::Reset
void Reset()
Definition: FullyDistSpVec.cpp:1775
combblas::FullyDistSpVec::ApplyInd
void ApplyInd(_BinaryOperation __binary_op)
Definition: FullyDistSpVec.h:249
combblas::FullyDistSpVec::ReadDistribute
std::ifstream & ReadDistribute(std::ifstream &infile, int master, HANDLER handler)
Totally obsolete version that only accepts an ifstream object and ascii files.
Definition: FullyDistSpVec.cpp:1397
Exception.h
combblas::FullyDistSpVec::operator[]
NT operator[](IT indx)
Definition: FullyDistSpVec.cpp:464
combblas::FullyDistSpVec::ReadDistribute
std::ifstream & ReadDistribute(std::ifstream &infile, int master)
Definition: FullyDistSpVec.h:160
combblas::FullyDistSpVec::TransposeVector
friend void TransposeVector(MPI_Comm &World, const FullyDistSpVec< IU, NV > &x, int32_t &trxlocnz, IU &lenuntil, int32_t *&trxinds, NV *&trxnums, bool indexisvalue)
Definition: ParFriends.h:1057
combblas::FullyDistSpVec::EWiseApply_threaded
friend FullyDistSpVec< IU, RET > EWiseApply_threaded(const FullyDistSpVec< IU, NU1 > &V, const FullyDistVec< IU, NU2 > &W, _BinaryOperation _binary_op, _BinaryPredicate _doOp, bool allowVNulls, NU1 Vzero, const bool useExtendedBinOp)
Definition: ParFriends.h:1981
combblas::FullyDistSpVec::WasFound
bool WasFound() const
Definition: FullyDistSpVec.h:190
promote.h
FullyDistSpVec.cpp
combblas::FullyDistSpVec::SelectApply
void SelectApply(const FullyDistVec< IT, NT1 > &denseVec, _UnaryOperation __unop, _BinaryOperation __binop)
Definition: FullyDistSpVec.cpp:2454
combblas::FullyDistSpVec::operator-=
FullyDistSpVec< IT, NT > & operator-=(const FullyDistSpVec< IT, NT > &rhs)
Definition: FullyDistSpVec.cpp:1094
OptBuf.h
combblas::FullyDistSpVec::Uniq
FullyDistSpVec< IT, NT > Uniq(_BinaryOperation __binary_op, MPI_Op mympiop)
Definition: FullyDistSpVec.cpp:1029
combblas::FullyDistSpVec::PermMat1
friend SpParMat< IU, bool, DER > PermMat1(const FullyDistSpVec< IU, NU > &ri, const IU ncol, _UnaryOperation __unop)
combblas::FullyDistSpVec
Definition: FullyDistSpVec.h:73
combblas::FullyDistSpVec::FindVals
FullyDistVec< IT, NT > FindVals(_Predicate pred) const
Definition: FullyDistSpVec.cpp:330
combblas::FullyDistSpVec::operator==
bool operator==(const FullyDistSpVec< IT, NT > &rhs) const
Definition: FullyDistSpVec.h:176
combblas::FullyDistSpVec::Apply
void Apply(_UnaryOperation __unary_op)
Definition: FullyDistSpVec.h:237
combblas::FullyDistSpVec::MergeContributions
friend void MergeContributions(FullyDistSpVec< IU, VT > &y, int *&recvcnt, int *&rdispls, int32_t *&recvindbuf, VT *&recvnumbuf, int rowneighs)
Helper functions for sparse matrix X sparse vector.
Definition: BFSFriends.h:224
combblas::SpParMat
Definition: BFSFriends.h:47
combblas::FullyDistSpVec::RandPerm
friend void RandPerm(FullyDistSpVec< IU, IU > &V)
combblas::FullyDistSpVec::operator+=
FullyDistSpVec< IT, NT > & operator+=(const FullyDistSpVec< IT, NT > &rhs)
Definition: FullyDistSpVec.cpp:1035
combblas::FullyDistSpVec::stealFrom
void stealFrom(FullyDistSpVec< IT, NT > &victim)
Definition: FullyDistSpVec.cpp:456
combblas::FullyDistSpVec::SaveGathered
void SaveGathered(std::ofstream &outfile, int master, HANDLER handler, bool printProcSplits=false)
Definition: FullyDistSpVec.cpp:1529
A
double A
combblas::FullyDistSpVec::operator=
FullyDistSpVec< IT, NT > & operator=(const FullyDistSpVec< IT, NT > &rhs)
Definition: FullyDistSpVec.cpp:70
combblas::FullyDistSpVec::SaveGathered
void SaveGathered(std::ofstream &outfile, int master)
Definition: FullyDistSpVec.h:164
combblas::SparseVectorLocalIterator
Definition: FullyDistSpVec.h:56
combblas::FullyDistSpVec::BulkSet
void BulkSet(IT inds[], int count)
Definition: FullyDistSpVec.cpp:1783
combblas::FullyDistSpVec::sort
FullyDistSpVec< IT, IT > sort()
sort the vector itself, return the permutation vector (0-based)
Definition: FullyDistSpVec.cpp:769
combblas::FullyDistSpVec::ScalarReadSaveHandler::read
NT read(std::basic_istream< c, t > &is, IT index)
Definition: FullyDistSpVec.h:134
CombBLAS.h
combblas::FullyDistSpVec::SpMV
friend FullyDistSpVec< IU, typename promote_trait< NUM, NUV >::T_promote > SpMV(const SpParMat< IU, NUM, UDER > &A, const FullyDistSpVec< IU, NUV > &x)
Definition: ParFriends.h:1675
combblas::FullyDistSpVec::ScalarReadSaveHandler
Definition: FullyDistSpVec.h:128
combblas::FullyDistSpVec::FilterByVal
void FilterByVal(FullyDistSpVec< IT, IT > Selector, _UnaryOperation __unop, bool filterByIndex)
Definition: FullyDistSpVec.cpp:2517
combblas::FullyDist
Definition: FullyDist.h:40
combblas::FullyDistSpVec::InvertRMA
FullyDistSpVec< IT, NT > InvertRMA(IT globallen, _BinaryOperationIdx __binopIdx, _BinaryOperationVal __binopVal)
Definition: FullyDistSpVec.cpp:2172
combblas::FullyDistSpVec::getcommgrid
std::shared_ptr< CommGrid > getcommgrid() const
Definition: FullyDistSpVec.h:269
combblas::FullyDistSpVec::getnnz
IT getnnz() const
Definition: FullyDistSpVec.h:208
combblas::FullyDistSpVec::getlocnnz
IT getlocnnz() const
Definition: FullyDistSpVec.h:204
combblas::FullyDistSpVec::FullyDistSpVec
friend class FullyDistSpVec
Definition: FullyDistSpVec.h:307
combblas::FullyDistSpVec::operator()
FullyDistVec< IT, NT > operator()(const FullyDistVec< IT, IT > &ri) const
SpRef (expects ri to be 0-based)
Definition: FullyDistSpVec.cpp:565
SpParMat.h
combblas::FullyDistSpVec::iota
void iota(IT globalsize, NT first)
Definition: FullyDistSpVec.cpp:676
combblas::FullyDistSpVec::Reduce
NT Reduce(_BinaryOperation __binary_op, NT init) const
Definition: FullyDistSpVec.cpp:1650
combblas::FullyDistSpVec::EWiseApply
friend FullyDistSpVec< IU, RET > EWiseApply(const FullyDistSpVec< IU, NU1 > &V, const FullyDistVec< IU, NU2 > &W, _BinaryOperation _binary_op, _BinaryPredicate _doOp, bool allowVNulls, NU1 Vzero, const bool useExtendedBinOp)
Definition: ParFriends.h:2127
combblas::PreAllocatedSPA
Definition: PreAllocatedSPA.h:42
combblas::DistEdgeList
Definition: DistEdgeList.h:81
combblas::FullyDistSpVec::RenameVertices
friend void RenameVertices(DistEdgeList< IU > &DEL)
Definition: DistEdgeList.cpp:364
combblas::FullyDistSpVec::EWiseMult
friend FullyDistSpVec< IU, typename promote_trait< NU1, NU2 >::T_promote > EWiseMult(const FullyDistSpVec< IU, NU1 > &V, const FullyDistVec< IU, NU2 > &W, bool exclude, NU2 zero)
Definition: ParFriends.h:1892
combblas::OptBuf
Definition: OptBuf.h:43
init
int init
Definition: ApproxWeightPerfectMatching.cpp:32
combblas::FullyDistSpVec::DelElement
void DelElement(IT indx)
Definition: FullyDistSpVec.cpp:541
CommGrid.h
combblas::FullyDistSpVec::operator=
FullyDistSpVec< IT, NT > & operator=(NT fixedval)
Definition: FullyDistSpVec.h:116
combblas::FullyDistSpVec::ParallelWrite
void ParallelWrite(const std::string &filename, bool onebased, bool includeindices=true)
Definition: FullyDistSpVec.h:150
combblas::FullyDistSpVec::Count
IT Count(_Predicate pred) const
Return the number of elements for which pred is true.
Definition: FullyDistSpVec.cpp:1639
combblas::FullyDistSpVec::DebugPrint
void DebugPrint()
Definition: FullyDistSpVec.cpp:1694
combblas::FullyDistSpVec::ParallelWrite
void ParallelWrite(const std::string &filename, bool onebased, HANDLER handler, bool includeindices=true, bool includeheader=false)
Definition: FullyDistSpVec.cpp:1310
combblas::FullyDistSpVec::Invert
FullyDistSpVec< IT, NT > Invert(IT globallen)
Definition: FullyDistSpVec.cpp:1901
combblas::FullyDistSpVec::ParallelRead
void ParallelRead(const std::string &filename, bool onebased, _BinaryOperation BinOp)
Definition: FullyDistSpVec.cpp:1209