COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
test_mpipspgemm.cpp
Go to the documentation of this file.
1#include <mpi.h>
2#include <sys/time.h>
3#include <iostream>
4#include <functional>
5#include <algorithm>
6#include <vector>
7#include <string>
8#include <sstream>
9#include <stdint.h>
10#include <cmath>
11#include "CombBLAS/CombBLAS.h"
12#include "Glue.h"
13#include "CCGrid.h"
14#include "Reductions.h"
15#include "Multiplier.h"
16#include "SplitMatDist.h"
17
18using namespace std;
19using namespace combblas;
20
30
31#define ITERS 1
32
33int main(int argc, char *argv[])
34{
35 int provided;
36 //MPI_Init_thread(&argc, &argv, MPI_THREAD_SINGLE, &provided);
37
38
41 {
42 printf("ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
44 }
45
46
47 int nprocs, myrank;
50
51 if(argc != 8)
52 {
53 cout << argc << endl;
54 if(myrank == 0)
55 {
56 printf("Usage (input): ./mpipspgemm <GridRows> <GridCols> <Layers> <matA> <matB> <matC> <algo>\n");
57 printf("Example: ./mpipspgemm 4 4 2 matA.mtx matB.mtx matB.mtx threaded\n");
58 printf("algo: outer | column | threaded \n");
59 }
60 return -1;
61 }
62
63 unsigned GRROWS = (unsigned) atoi(argv[1]);
64 unsigned GRCOLS = (unsigned) atoi(argv[2]);
65 unsigned C_FACTOR = (unsigned) atoi(argv[3]);
67
68 if(GRROWS != GRCOLS)
69 {
70 SpParHelper::Print("This version of the Combinatorial BLAS only works on a square logical processor grid\n");
73 }
74
77 {
78 SpParHelper::Print("The product of <GridRows> <GridCols> <Replicas> does not match the number of processes\n");
81 }
82
83
84
87 string type;
88
89
90 string fileA(argv[4]);
91 string fileB(argv[5]);
92 string fileC(argv[6]);
93
94 {
96 layerGrid.reset( new CommGrid(CMG.layerWorld, 0, 0) );
97 FullyDistVec<int32_t, int32_t> p(layerGrid); // permutation vector defined on layers
98
99 double t01 = MPI_Wtime();
100
104
105 SplitMat(CMG, A, splitA, false);
106 SplitMat(CMG, B, splitB, true);
107 SplitMat(CMG, C, controlC, false);
108
109 if(myrank == 0) cout << "Matrices read and replicated along layers : time " << MPI_Wtime() - t01 << endl;
110
111 type = string(argv[7]);
112
113 if(type == string("outer"))
114 {
115 for(int k=0; k<ITERS; k++)
116 {
117 splitB.Transpose(); // locally "transpose" [ABAB: check correctness]
118 splitC = multiply(splitA, splitB, CMG, true, false); // outer product
119 if (controlC == *splitC)
120 SpParHelper::Print("Outer product multiplication working correctly\n");
121 else
122 SpParHelper::Print("ERROR in Outer product multiplication, go fix it!\n");
123 delete splitC;
124 }
125
126 }
127 else if(type == string("column"))
128 {
129
130 for(int k=0; k<ITERS; k++)
131 {
132 splitC = multiply(splitA, splitB, CMG, false, false);
133 if (controlC == *splitC)
134 SpParHelper::Print("Col-heap multiplication working correctly\n");
135 else
136 SpParHelper::Print("ERROR in Col-heap multiplication, go fix it!\n");
137
138 delete splitC;
139 }
140
141 }
142 else // default threaded
143 {
144 for(int k=0; k<ITERS; k++)
145 {
146 splitC = multiply(splitA, splitB, CMG, false, true);
147 if (controlC == *splitC)
148 SpParHelper::Print("Col-heap-threaded multiplication working correctly\n");
149 else
150 SpParHelper::Print("ERROR in Col-heap-threaded multiplication, go fix it!\n");
151 delete splitC;
152 }
153 }
154 }
155
156 MPI_Finalize();
157 return 0;
158}
159
160
int main()
Definition Driver.cpp:12
Definition test.cpp:53
static void Print(const std::string &s)
int nprocs
Definition comms.cpp:55
void SplitMat(CCGrid &CMG, SpDCCols< IT, NT > *localmat, SpDCCols< IT, NT > &splitmat, bool rowsplit=false)
SpDCCols< IT, NT > * multiply(SpDCCols< IT, NT > &splitA, SpDCCols< IT, NT > &splitB, CCGrid &CMG, bool isBT, bool threaded)
Definition Multiplier.h:11
double A
double C
Definition options.h:15
double comm_reduce
double comp_split
#define ITERS
double comm_bcast
double comm_split
double comp_reduce
double comp_result
double comp_trans
double comp_summa
double comp_reduce_layer