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
60
int
main
(
int
argc
,
char
*
argv
[])
61
{
62
int
nprocs, myrank;
63
MPI_Init
(&
argc
, &
argv
);
64
MPI_Comm_size
(
MPI_COMM_WORLD
,&nprocs);
65
MPI_Comm_rank
(
MPI_COMM_WORLD
,&myrank);
66
67
if
(
argc
< 8)
68
{
69
if
(myrank == 0)
70
{
71
cout
<<
"Usage: ./SpAsgnTest <BASEADDRESS> <Matrix> <PrunedMatrix> <RHSMatrix> <AssignedMatrix> <VectorRowIndices> <VectorColIndices>"
<<
endl
;
72
cout
<<
"Example: ./SpAsgnTest ../mfiles A_100x100.txt A_with20x30hole.txt dense_20x30matrix.txt A_wdenseblocks.txt 20outta100.txt 30outta100.txt"
<<
endl
;
73
cout
<<
"Input files should be under <BASEADDRESS> in tuples format"
<<
endl
;
74
}
75
MPI_Finalize
();
76
return
-1;
77
}
78
{
79
string
directory
(
argv
[1]);
80
string
normalname
(
argv
[2]);
81
string
prunedname
(
argv
[3]);
82
string
rhsmatname
(
argv
[4]);
83
string
assignname
(
argv
[5]);
84
string
vec1name
(
argv
[6]);
85
string
vec2name
(
argv
[7]);
86
normalname
=
directory
+
"/"
+
normalname
;
87
prunedname
=
directory
+
"/"
+
prunedname
;
88
rhsmatname
=
directory
+
"/"
+
rhsmatname
;
89
assignname
=
directory
+
"/"
+
assignname
;
90
vec1name
=
directory
+
"/"
+
vec1name
;
91
vec2name
=
directory
+
"/"
+
vec2name
;
92
93
ifstream
inputvec1
(
vec1name
.c_str());
94
ifstream
inputvec2
(
vec2name
.c_str());
95
96
if
(myrank == 0)
97
{
98
if
(
inputvec1
.fail() ||
inputvec2
.fail())
99
{
100
cout
<<
"One of the input vector files do not exist, aborting"
<<
endl
;
101
MPI_Abort
(
MPI_COMM_WORLD
,
NOFILE
);
102
return
-1;
103
}
104
}
105
MPI_Barrier
(
MPI_COMM_WORLD
);
106
typedef
SpParMat <int, double , SpDCCols<int,double>
>
PARDBMAT
;
107
shared_ptr<CommGrid>
fullWorld
;
108
fullWorld
.reset(
new
CommGrid
(
MPI_COMM_WORLD
, 0, 0) );
109
110
PARDBMAT
A
(
fullWorld
);
111
PARDBMAT
Apr
(
fullWorld
);
112
PARDBMAT
B
(
fullWorld
);
113
PARDBMAT
C
(
fullWorld
);
114
FullyDistVec<int,int>
vec1
(
fullWorld
);
115
FullyDistVec<int,int>
vec2
(
fullWorld
);
116
117
A
.ReadDistribute(
normalname
, 0);
118
Apr
.ReadDistribute(
prunedname
, 0);
119
B
.ReadDistribute(
rhsmatname
, 0);
120
C
.ReadDistribute(
assignname
, 0);
121
vec1
.ReadDistribute(
inputvec1
, 0);
122
vec2
.ReadDistribute(
inputvec2
, 0);
123
124
vec1
.Apply(
bind2nd
(
minus<int>
(), 1));
// For 0-based indexing
125
vec2
.Apply(
bind2nd
(
minus<int>
(), 1));
126
127
PARDBMAT
Atemp
=
A
;
128
Atemp
.Prune(
vec1
,
vec2
);
129
130
// We should get the original A back.
131
if
(
Atemp
==
Apr
)
132
{
133
SpParHelper::Print
(
"Pruning is working\n"
);
134
}
135
else
136
{
137
SpParHelper::Print
(
"Error in pruning, go fix it\n"
);
138
}
139
140
A
.SpAsgn(
vec1
,
vec2
,
B
);
141
if
(
A
==
C
)
142
{
143
SpParHelper::Print
(
"SpAsgn working correctly\n"
);
144
}
145
else
146
{
147
SpParHelper::Print
(
"ERROR in SpAsgn, go fix it!\n"
);
148
A
.SaveGathered(
"Erroneous_SpAsgnd.txt"
);
149
}
150
151
FullyDistVec<int,int>
crow
(
fullWorld
);
152
FullyDistVec<int,int>
ccol
(
fullWorld
);
153
FullyDistVec<int,double>
cval
(
fullWorld
);
154
A
.Find(
crow
,
ccol
,
cval
);
155
FullyDistSpVec<int, double>
sval
=
cval
;
156
// sval.DebugPrint();
157
158
pair< FullyDistVec<int,int>
,
FullyDistVec<int,double>
>
ptopk
;
159
ptopk
=
TopK
(
sval
, 3);
160
//ptopk.first.DebugPrint();
161
//ptopk.second.DebugPrint();
162
163
164
inputvec1
.clear();
165
inputvec1
.close();
166
inputvec2
.clear();
167
inputvec2
.close();
168
}
169
MPI_Finalize
();
170
return
0;
171
}
TopK
pair< FullyDistVec< IT, IT >, FullyDistVec< IT, NT > > TopK(FullyDistSpVec< IT, NT > &v, IT k)
Definition
SpAsgnTest.cpp:42
main
int main(int argc, char *argv[])
Definition
SpAsgnTest.cpp:60
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:811
NOFILE
#define NOFILE
Definition
SpDefs.h:75
combblas
Definition
CCGrid.h:4
A
double A
C
double C
Definition
options.h:15
B
double B
Definition
options.h:15
time.h
ReleaseTests
SpAsgnTest.cpp
Generated by
1.9.8