COMBINATORIAL_BLAS
1.6
Loading...
Searching...
No Matches
SpAsgnTest.cpp
Go to the documentation of this file.
1
/****************************************************************/
2
/* Parallel Combinatorial BLAS Library (for Graph Computations) */
3
/* version 1.5 -------------------------------------------------*/
4
/* date: 10/09/2015 ---------------------------------------------*/
5
/* authors: Ariful Azad, Aydin Buluc, Adam Lugowski ------------*/
6
/****************************************************************/
7
/*
8
Copyright (c) 2010-2015, 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
#include <mpi.h>
30
#include <
sys/time.h
>
31
#include <iostream>
32
#include <functional>
33
#include <algorithm>
34
#include <vector>
35
#include <sstream>
36
#include "CombBLAS/CombBLAS.h"
37
38
using namespace
std;
39
using namespace
combblas
;
40
41
template
<
typename
IT,
typename
NT>
42
pair< FullyDistVec<IT,IT>
,
FullyDistVec<IT,NT>
>
TopK
(
FullyDistSpVec<IT,NT>
& v,
IT
k)
43
{
44
// FullyDistVec::FullyDistVec(IT glen, NT initval)
45
FullyDistVec<IT,IT>
sel
(v.getcommgrid(), k, 0);
46
47
//void FullyDistVec::iota(IT globalsize, NT first)
48
sel
.iota(k, v.TotalLength() - k);
49
50
FullyDistSpVec<IT,NT>
sorted
(v);
51
FullyDistSpVec<IT,IT>
perm
=
sorted
.sort();
52
53
// FullyDistVec FullyDistSpVec::operator(FullyDistVec & v)
54
FullyDistVec<IT,IT>
topkind
=
perm
(
sel
);
55
FullyDistVec<IT,NT>
topkele
= v(
topkind
);
56
return
make_pair
(
topkind
,
topkele
);
57
}
58
59
struct
mypair
60
{
61
mypair
(
double
rhs
)
62
{
63
val
=
make_pair
(
rhs
, -
rhs
);
64
}
65
mypair
()
66
{
67
val
=
make_pair
(1,-1);
68
}
69
pair<double, double>
val
;
70
};
71
72
73
int
main
(
int
argc
,
char
*
argv
[])
74
{
75
int
nprocs
, myrank;
76
MPI_Init
(&
argc
, &
argv
);
77
MPI_Comm_size
(
MPI_COMM_WORLD
,&
nprocs
);
78
MPI_Comm_rank
(
MPI_COMM_WORLD
,&myrank);
79
80
if
(
argc
< 8)
81
{
82
if
(myrank == 0)
83
{
84
cout
<<
"Usage: ./SpAsgnTest <BASEADDRESS> <Matrix> <PrunedMatrix> <RHSMatrix> <AssignedMatrix> <VectorRowIndices> <VectorColIndices>"
<<
endl
;
85
cout
<<
"Example: ./SpAsgnTest TESTDATA/ A_100x100.txt A_with20x30hole.txt dense_20x30matrix.txt A_wdenseblocks.txt 20outta100.txt 30outta100.txt"
<<
endl
;
86
cout
<<
"Input files should be under <BASEADDRESS> in triples format"
<<
endl
;
87
}
88
MPI_Finalize
();
89
return
-1;
90
}
91
{
92
string
directory
(
argv
[1]);
93
string
normalname
(
argv
[2]);
94
string
prunedname
(
argv
[3]);
95
string
rhsmatname
(
argv
[4]);
96
string
assignname
(
argv
[5]);
97
string
vec1name
(
argv
[6]);
98
string
vec2name
(
argv
[7]);
99
normalname
=
directory
+
"/"
+
normalname
;
100
prunedname
=
directory
+
"/"
+
prunedname
;
101
rhsmatname
=
directory
+
"/"
+
rhsmatname
;
102
assignname
=
directory
+
"/"
+
assignname
;
103
vec1name
=
directory
+
"/"
+
vec1name
;
104
vec2name
=
directory
+
"/"
+
vec2name
;
105
106
ifstream
inputvec1
(
vec1name
.c_str());
107
ifstream
inputvec2
(
vec2name
.c_str());
108
109
if
(myrank == 0)
110
{
111
if
(
inputvec1
.fail() ||
inputvec2
.fail())
112
{
113
cout
<<
"One of the input vector files do not exist, aborting"
<<
endl
;
114
MPI_Abort
(
MPI_COMM_WORLD
,
NOFILE
);
115
return
-1;
116
}
117
}
118
MPI_Barrier
(
MPI_COMM_WORLD
);
119
typedef
SpParMat <int64_t, double , SpDCCols<int64_t,double>
>
PARDBMAT
;
120
typedef
SpParMat <int64_t, mypair , SpDCCols<int64_t, mypair >
>
PARPAIRMAT
;
121
shared_ptr<CommGrid>
fullWorld
;
122
fullWorld
.reset(
new
CommGrid
(
MPI_COMM_WORLD
, 0, 0) );
123
124
PARDBMAT
A
(
fullWorld
);
125
PARDBMAT
Apr
(
fullWorld
);
126
PARDBMAT
B
(
fullWorld
);
127
PARDBMAT
C
(
fullWorld
);
128
FullyDistVec<int64_t,int64_t>
vec1
(
fullWorld
);
129
FullyDistVec<int64_t,int64_t>
vec2
(
fullWorld
);
130
131
A
.ReadDistribute(
normalname
, 0);
132
Apr
.ReadDistribute(
prunedname
, 0);
133
B
.ReadDistribute(
rhsmatname
, 0);
134
C
.ReadDistribute(
assignname
, 0);
135
vec1
.ReadDistribute(
inputvec1
, 0);
136
vec2
.ReadDistribute(
inputvec2
, 0);
137
138
vec1
.Apply(
bind2nd
(
minus<int64_t>
(), 1));
// For 0-based indexing
139
vec2
.Apply(
bind2nd
(
minus<int64_t>
(), 1));
140
141
PARDBMAT
Atemp
=
A
;
142
Atemp
.Prune(
vec1
,
vec2
);
143
144
PARPAIRMAT
Apair
=
A
;
145
Apair
.Prune(
vec1
,
vec2
);
146
147
PARDBMAT
Apruned
=
A
;
148
Apruned
.PruneFull(
vec1
,
vec2
);
149
Apruned
.ParallelWriteMM(
"ArowscolsPruned.mtx"
,
true
);
150
151
// We should get the original A back.
152
if
(
Atemp
==
Apr
)
153
{
154
SpParHelper::Print
(
"Pruning is working\n"
);
155
}
156
else
157
{
158
SpParHelper::Print
(
"Error in pruning, go fix it\n"
);
159
}
160
161
A
.SpAsgn(
vec1
,
vec2
,
B
);
162
if
(
A
==
C
)
163
{
164
SpParHelper::Print
(
"SpAsgn working correctly\n"
);
165
}
166
else
167
{
168
SpParHelper::Print
(
"ERROR in SpAsgn, go fix it!\n"
);
169
A
.SaveGathered(
"Erroneous_SpAsgnd.txt"
);
170
}
171
172
FullyDistVec<int64_t,int64_t>
crow
(
fullWorld
);
173
FullyDistVec<int64_t,int64_t>
ccol
(
fullWorld
);
174
FullyDistVec<int64_t,double>
cval
(
fullWorld
);
175
A
.Find(
crow
,
ccol
,
cval
);
176
FullyDistSpVec<int64_t, double>
sval
=
cval
;
177
// sval.DebugPrint();
178
179
pair< FullyDistVec<int64_t,int64_t>
,
FullyDistVec<int64_t,double>
>
ptopk
;
180
ptopk
=
TopK
(
sval
, (
int64_t
) 3);
181
//ptopk.first.DebugPrint();
182
//ptopk.second.DebugPrint();
183
184
185
inputvec1
.clear();
186
inputvec1
.close();
187
inputvec2
.clear();
188
inputvec2
.close();
189
}
190
MPI_Finalize
();
191
return
0;
192
}
IT
int64_t IT
Definition
BlockedSpGEMM.cpp:9
main
int main()
Definition
Driver.cpp:12
TopK
pair< FullyDistVec< IT, IT >, FullyDistVec< IT, NT > > TopK(FullyDistSpVec< IT, NT > &v, IT k)
Definition
SpAsgnTest.cpp:42
B
Definition
test.cpp:53
combblas::CommGrid
Definition
CommGrid.h:45
combblas::DistEdgeList
Definition
DistEdgeList.h:82
combblas::SpParHelper::Print
static void Print(const std::string &s)
Definition
SpParHelper.cpp:836
nprocs
int nprocs
Definition
comms.cpp:55
int64_t
long int64_t
Definition
compat.h:21
NOFILE
#define NOFILE
Definition
SpDefs.h:75
combblas
Definition
CCGrid.h:4
A
double A
C
double C
Definition
options.h:15
mypair
Definition
SpAsgnTest.cpp:60
mypair::mypair
mypair(double rhs)
Definition
SpAsgnTest.cpp:61
mypair::val
pair< double, double > val
Definition
SpAsgnTest.cpp:69
mypair::mypair
mypair()
Definition
SpAsgnTest.cpp:65
time.h
ReleaseTests
SpAsgnTest.cpp
Generated by
1.9.8