1 #include <mpi.h>
2 #include <sys/time.h>
3 #include <iostream>
4 #include <functional>
5 #include <algorithm>
6 #include <vector>
7 #include <sstream>
8 #include "CombBLAS/CombBLAS.h"
9 
10 using namespace std;
11 using namespace combblas;
12 #define ITERATIONS 1
13 
14 // Simple helper class for declarations: Just the numerical type is templated
15 // The index type and the sequential matrix type stays the same for the whole code
16 // In this case, they are "int" and "SpDCCols"
17 template <class NT>
18 class PSpMat
19 {
20 public:
21 	typedef SpDCCols < int, NT > DCCols;
22 	typedef SpParMat < int, NT, DCCols > MPI_DCCols;
23 };
24 
25 #define ElementType double
26 
27 
main(int argc,char * argv[])28 int main(int argc, char* argv[])
29 {
30 	int nprocs, myrank;
31 	MPI_Init(&argc, &argv);
32 	MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
33 	MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
34 
35 	if(argc < 3)
36 	{
37 		if(myrank == 0)
38 		{
39 			cout << "Usage: ./MultTest <MatrixA> <MatrixB>" << endl;
40 			cout << "<MatrixA>,<MatrixB> are absolute addresses, and files should be in triples format" << endl;
41 		}
42 		MPI_Finalize();
43 		return -1;
44 	}
45 	{
46 		string Aname(argv[1]);
47 		string Bname(argv[2]);
48 		typedef PlusTimesSRing<ElementType, ElementType> PTDOUBLEDOUBLE;
49 		PSpMat<ElementType>::MPI_DCCols A, B;	// construct objects
50 
51 		A.ReadDistribute(Aname, 0);
52 		A.PrintInfo();
53 		B.ReadDistribute(Bname, 0);
54 		B.PrintInfo();
55 		SpParHelper::Print("Data read\n");
56 
57 		{ // force the calling of C's destructor
58 			PSpMat<ElementType>::MPI_DCCols C = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
59 			int64_t cnnz = C.getnnz();
60 			ostringstream tinfo;
61 			tinfo << "C has a total of " << cnnz << " nonzeros" << endl;
62 			SpParHelper::Print(tinfo.str());
63 			SpParHelper::Print("Warmed up for DoubleBuff\n");
64 			C.PrintInfo();
65 		}
66 		MPI_Barrier(MPI_COMM_WORLD);
67 		MPI_Pcontrol(1,"SpGEMM_DoubleBuff");
68 		double t1 = MPI_Wtime(); 	// initilize (wall-clock) timer
69 		for(int i=0; i<ITERATIONS; i++)
70 		{
71 			PSpMat<ElementType>::MPI_DCCols C = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
72 		}
73 		MPI_Barrier(MPI_COMM_WORLD);
74 		double t2 = MPI_Wtime();
75 		MPI_Pcontrol(-1,"SpGEMM_DoubleBuff");
76 		if(myrank == 0)
77 		{
78 			cout<<"Double buffered multiplications finished"<<endl;
79 			printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
80 		}
81 
82 		{// force the calling of C's destructor
83 			PSpMat<ElementType>::MPI_DCCols C = Mult_AnXBn_Synch<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
84 			C.PrintInfo();
85 		}
86 		SpParHelper::Print("Warmed up for Synch\n");
87 		MPI_Barrier(MPI_COMM_WORLD);
88 		MPI_Pcontrol(1,"SpGEMM_Synch");
89 		t1 = MPI_Wtime(); 	// initilize (wall-clock) timer
90 		for(int i=0; i<ITERATIONS; i++)
91 		{
92 			PSpMat<ElementType>::MPI_DCCols C = Mult_AnXBn_Synch<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
93 		}
94 		MPI_Barrier(MPI_COMM_WORLD);
95 		MPI_Pcontrol(-1,"SpGEMM_Synch");
96 		t2 = MPI_Wtime();
97 		if(myrank == 0)
98 		{
99 			cout<<"Synchronous multiplications finished"<<endl;
100 			printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
101 		}
102 
103 
104 		/*
105 		C = Mult_AnXBn_ActiveTarget<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
106 		SpParHelper::Print("Warmed up for ActiveTarget\n");
107 		MPI_Barrier(MPI_COMM_WORLD);
108 		MPI_Pcontrol(1,"SpGEMM_ActiveTarget");
109 		t1 = MPI_Wtime(); 	// initilize (wall-clock) timer
110 		for(int i=0; i<ITERATIONS; i++)
111 		{
112 			C = Mult_AnXBn_ActiveTarget<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
113 		}
114 		MPI_Barrier(MPI_COMM_WORLD);
115 		MPI_Pcontrol(-1,"SpGEMM_ActiveTarget");
116 		t2 = MPI_Wtime();
117 		if(myrank == 0)
118 		{
119 			cout<<"Active target multiplications finished"<<endl;
120 			printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
121 		}
122 
123 		C = Mult_AnXBn_PassiveTarget<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
124 		SpParHelper::Print("Warmed up for PassiveTarget\n");
125 		MPI_Barrier(MPI_COMM_WORLD);
126 		MPI_Pcontrol(1,"SpGEMM_PassiveTarget");
127 		t1 = MPI_Wtime(); 	// initilize (wall-clock) timer
128 		for(int i=0; i<ITERATIONS; i++)
129 		{
130 			C = Mult_AnXBn_PassiveTarget<PTDOUBLEDOUBLE, ElementType, PSpMat<ElementType>::DCCols >(A, B);
131 		}
132 		MPI_Barrier(MPI_COMM_WORLD);
133 		MPI_Pcontrol(-1,"SpGEMM_PassiveTarget");
134 		t2 = MPI_Wtime();
135 		if(myrank == 0)
136 		{
137 			cout<<"Passive target multiplications finished"<<endl;
138 			printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
139 		}
140 		*/
141 	}
142 	MPI_Finalize();
143 	return 0;
144 }
145 
146