COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
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
44namespace combblas {
45
46template <class IT, class NT, class DER>
47class SpParMat;
48
49template <class IT>
50class DistEdgeList;
51
52template <class IU, class NU>
53class FullyDistVec;
54
55template <class IU, class NU>
56class SparseVectorLocalIterator;
57
72template <class IT, class NT>
73class FullyDistSpVec: public FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>
74{
75public:
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>
83 FullyDistSpVec (const FullyDistVec<IT,NT> & rhs); // Conversion copy-constructor
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
90 template <typename _BinaryOperationIdx, typename _BinaryOperationVal, typename _BinaryOperationDuplicate>
92 template <typename _BinaryOperationIdx, typename _BinaryOperationVal>
94
95
96 template <typename NT1, typename _UnaryOperation>
98 template <typename _UnaryOperation>
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>
108
109
110
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
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);
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);
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
177 {
178 FullyDistVec<IT,NT> v = *this;
180 return (v == w);
181 }
182
183 void PrintInfo(std::string vecname) const;
187 void SetElement (IT indx, NT numx); // element-wise assignment
188 void DelElement (IT indx); // element-wise deletion
190 bool WasFound() const { return wasFound; }
191
194
195#if __cplusplus > 199711L
196 template <typename _BinaryOperation = minimum<NT> >
198#else
199 template <typename _BinaryOperation >
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 }
215 using FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>::LengthUntil;
216 using FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>::MyLocLength;
217 using FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>::MyRowLength;
218 using FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>::TotalLength;
220 using FullyDist<IT,NT,typename combblas::disable_if< combblas::is_boolean<NT>::value, NT >::type>::RowLenUntil;
221
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>
235
236 template <typename _UnaryOperation>
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>
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>
264
265 template <typename OUT, typename _BinaryOperation, typename _UnaryOperation>
267
269 std::shared_ptr<CommGrid> getcommgrid() const { return commGrid; }
270
271 void Reset();
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>
279 template <typename _Predicate>
281
282
283protected:
286
287private:
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> >
300#else
301 template <typename _BinaryOperation >
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>
321
322 template <typename SR, typename IU, typename NUM, typename UDER>
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>
331
332 template <typename SR, typename IVT, typename OVT, typename IU, typename NUM, typename UDER>
334
335 template <typename IU, typename NU1, typename NU2>
338
339 template <typename RET, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
342
343 template <typename RET, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
346
347
348 template <typename RET, typename IU, typename NU1, typename NU2, typename _BinaryOperation, typename _BinaryPredicate>
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>
368
369 template <class IU, class NU, class DER, typename _UnaryOperation>
371};
372
373}
374
375#include "FullyDistSpVec.cpp"
376
377#endif
NT read(std::basic_istream< c, t > &is, IT index)
void save(std::basic_ostream< c, t > &os, const NT &v, IT index)
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)
friend FullyDistSpVec< IU, VT > SpMV(const SpParMat< IU, bool, UDER > &A, const FullyDistSpVec< IU, VT > &x, OptBuf< int32_t, VT > &optbuf)
FullyDistSpVec< IT, NT > & operator=(const FullyDistSpVec< IT, NT > &rhs)
void Select(const FullyDistVec< IT, NT1 > &denseVec, _UnaryOperation unop)
FullyDistSpVec< IT, IT > sort()
sort the vector itself, return the permutation vector (0-based)
std::ifstream & ReadDistribute(std::ifstream &infile, int master, HANDLER handler)
Totally obsolete version that only accepts an ifstream object and ascii files.
FullyDistSpVec< IT, NT > & operator-=(const FullyDistSpVec< IT, NT > &rhs)
NT Reduce(_BinaryOperation __binary_op, NT init) const
FullyDistSpVec< IT, NT > Uniq(_BinaryOperation __binary_op, MPI_Op mympiop)
FullyDistSpVec(std::shared_ptr< CommGrid > grid, IT globallen, const std::vector< IT > &indvec, const std::vector< NT > &numvec, bool SumDuplicates=false, bool sorted=false)
friend FullyDistSpVec< IU, typename promote_trait< NUM, NUV >::T_promote > SpMV(const SpParMat< IU, NUM, UDER > &A, const FullyDistSpVec< IU, NUV > &x)
std::vector< NT > GetLocalNum()
void Apply(_UnaryOperation __unary_op)
FullyDistVec< IT, NT > FindVals(_Predicate pred) const
friend void RandPerm(FullyDistSpVec< IU, IU > &V)
void ParallelWrite(const std::string &filename, bool onebased, HANDLER handler, bool includeindices=true, bool includeheader=false)
FullyDistSpVec< IT, NT > Invert(IT globallen, _BinaryOperationIdx __binopIdx, _BinaryOperationVal __binopVal, _BinaryOperationDuplicate __binopDuplicate)
void FilterByVal(FullyDistSpVec< IT, IT > Selector, _UnaryOperation __unop, bool filterByIndex)
FullyDistSpVec< IT, NT > Invert(IT globallen)
std::shared_ptr< CommGrid > getcommgrid() const
void DelElement(IT indx)
void SelectApply(const FullyDistVec< IT, NT1 > &denseVec, _UnaryOperation __unop, _BinaryOperation __binop)
void ParallelWrite(const std::string &filename, bool onebased, bool includeindices=true)
IT Count(_Predicate pred) const
Return the number of elements for which pred is true.
FullyDistVec< IT, NT > operator()(const FullyDistVec< IT, IT > &ri) const
SpRef (expects ri to be 0-based)
void SaveGathered(std::ofstream &outfile, int master)
FullyDistSpVec< IT, NT > & operator=(NT fixedval)
void ParallelRead(const std::string &filename, bool onebased, _BinaryOperation BinOp)
FullyDistSpVec(std::shared_ptr< CommGrid > grid)
OUT Reduce(_BinaryOperation __binary_op, OUT default_val, _UnaryOperation __unary_op) const
bool operator==(const FullyDistSpVec< IT, NT > &rhs) const
std::ifstream & ReadDistribute(std::ifstream &infile, int master)
FullyDistSpVec(const FullyDistVec< IT, NT > &rhs)
friend void RenameVertices(DistEdgeList< IU > &DEL)
FullyDistSpVec(IT globalsize, const FullyDistVec< IT, IT > &inds, const FullyDistVec< IT, NT > &vals, bool SumDuplicates=false)
FullyDistVec< IT, IT > FindInds(_Predicate pred) const
void SaveGathered(std::ofstream &outfile, int master, HANDLER handler, bool printProcSplits=false)
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)
void BulkSet(IT inds[], int count)
friend void TransposeVector(MPI_Comm &World, const FullyDistSpVec< IU, NV > &x, int32_t &trxlocnz, IU &lenuntil, int32_t *&trxinds, NV *&trxnums, bool indexisvalue)
void ApplyInd(_BinaryOperation __binary_op)
FullyDistSpVec< IT, NT > & operator=(const FullyDistVec< IT, NT > &rhs)
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)
FullyDistSpVec< IT, NT > & operator+=(const FullyDistSpVec< IT, NT > &rhs)
std::vector< IT > GetLocalInd()
void PrintInfo(std::string vecname) const
FullyDistSpVec(const FullyDistVec< IT, NT > &rhs, _UnaryOperation unop)
void Setminus(const FullyDistSpVec< IT, NT1 > &other)
FullyDistSpVec(std::shared_ptr< CommGrid > grid, IT glen)
friend SpParMat< IU, bool, DER > PermMat1(const FullyDistSpVec< IU, NU > &ri, const IU ncol, _UnaryOperation __unop)
void iota(IT globalsize, NT first)
void SetElement(IT indx, NT numx)
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
void stealFrom(FullyDistSpVec< IT, NT > &victim)
FullyDistSpVec< IT, NT > InvertRMA(IT globallen, _BinaryOperationIdx __binopIdx, _BinaryOperationVal __binopVal)
double A
signed int int32_t
Definition stdint.h:77