COMBINATORIAL_BLAS  1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Galerkin.cpp
Go to the documentation of this file.
1 /****************************************************************/
2 /* Parallel Combinatorial BLAS Library (for Graph Computations) */
3 /* version 1.3 -------------------------------------------------*/
4 /* date: 2/1/2013 ----------------------------------------------*/
5 /* authors: Aydin Buluc (abuluc@lbl.gov), Adam Lugowski --------*/
6 /****************************************************************/
7 /*
8  Copyright (c) 2010-, Aydin Buluc
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.h"
37 
38 using namespace std;
39 #define ITERATIONS 10
40 
41 // Simple helper class for declarations: Just the numerical type is templated
42 // The index type and the sequential matrix type stays the same for the whole code
43 // In this case, they are "int" and "SpDCCols"
44 template <class NT>
45 class PSpMat
46 {
47 public:
50 };
51 
52 
53 int main(int argc, char* argv[])
54 {
55  int nprocs, myrank;
56  MPI_Init(&argc, &argv);
57  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
58  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
59 
60  if(argc < 4)
61  {
62  if(myrank == 0)
63  {
64  cout << "Usage: ./Galerkin <Matrix> <S> <STranspose>" << endl;
65  cout << "<Matrix>,<S>,<STranspose> are absolute addresses, and files should be in triples format" << endl;
66  }
67  MPI_Finalize();
68  return -1;
69  }
70  {
71  string Aname(argv[1]);
72  string Sname(argv[2]);
73  string STname(argv[3]);
74 
75  ifstream inputA(Aname.c_str());
76  ifstream inputS(Sname.c_str());
77  ifstream inputST(STname.c_str());
78 
79  MPI_Barrier(MPI_COMM_WORLD);
80  typedef PlusTimesSRing<double, double> PTDOUBLEDOUBLE;
81 
82  PSpMat<double>::MPI_DCCols A, S, ST; // construct objects
83 
84  A.ReadDistribute(inputA, 0);
85  S.ReadDistribute(inputS, 0);
86  ST.ReadDistribute(inputST, 0);
87  SpParHelper::Print("Data read\n");
88 
89  // force the calling of C's destructor
90  {
91  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(A, ST);
92  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(S, C);
93  SpParHelper::Print("Warmed up for DoubleBuff (right evaluate)\n");
94  }
95  MPI_Barrier(MPI_COMM_WORLD);
96  MPI_Pcontrol(1,"SpGEMM_DoubleBuff_right");
97  double t1 = MPI_Wtime(); // initilize (wall-clock) timer
98  for(int i=0; i<ITERATIONS; i++)
99  {
100  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(A, ST);
101  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(S, C);
102  }
103  MPI_Barrier(MPI_COMM_WORLD);
104  double t2 = MPI_Wtime();
105  MPI_Pcontrol(-1,"SpGEMM_DoubleBuff_right");
106  if(myrank == 0)
107  {
108  cout<<"Double buffered multiplications (right evaluate) finished"<<endl;
109  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
110  }
111 
112  // force the calling of C's destructor
113  {
114  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(S, A);
115  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(C, ST);
116  SpParHelper::Print("Warmed up for DoubleBuff (left evaluate)\n");
117  }
118  MPI_Barrier(MPI_COMM_WORLD);
119  MPI_Pcontrol(1,"SpGEMM_DoubleBuff_left");
120  t1 = MPI_Wtime(); // initilize (wall-clock) timer
121  for(int i=0; i<ITERATIONS; i++)
122  {
123  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(S, A);
124  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_DoubleBuff<PTDOUBLEDOUBLE>(C, ST);
125  }
126  MPI_Barrier(MPI_COMM_WORLD);
127  t2 = MPI_Wtime();
128  MPI_Pcontrol(-1,"SpGEMM_DoubleBuff_left");
129  if(myrank == 0)
130  {
131  cout<<"Double buffered multiplications (left evaluate) finished"<<endl;
132  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
133  }
134 
135  // force the calling of C's destructor
136  {
137  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(A, ST);
138  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(S, C);
139  }
140  SpParHelper::Print("Warmed up for Synch (right evaluate)\n");
141  MPI_Barrier(MPI_COMM_WORLD);
142  MPI_Pcontrol(1,"SpGEMM_Synch_right");
143  t1 = MPI_Wtime(); // initilize (wall-clock) timer
144  for(int i=0; i<ITERATIONS; i++)
145  {
146  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(A, ST);
147  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(S, C);
148  }
149  MPI_Barrier(MPI_COMM_WORLD);
150  MPI_Pcontrol(-1,"SpGEMM_Synch_right");
151  t2 = MPI_Wtime();
152  if(myrank == 0)
153  {
154  cout<<"Synchronous multiplications (right evaluate) finished"<<endl;
155  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
156  }
157 
158  // force the calling of C's destructor
159  {
160  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(S, A);
161  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(C, ST);
162  }
163  SpParHelper::Print("Warmed up for Synch (left evaluate)\n");
164  MPI_Barrier(MPI_COMM_WORLD);
165  MPI_Pcontrol(1,"SpGEMM_Synch_left");
166  t1 = MPI_Wtime(); // initilize (wall-clock) timer
167  for(int i=0; i<ITERATIONS; i++)
168  {
169  PSpMat<double>::MPI_DCCols C = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(S, A);
170  PSpMat<double>::MPI_DCCols D = Mult_AnXBn_Synch<PTDOUBLEDOUBLE>(C, ST);
171  }
172  MPI_Barrier(MPI_COMM_WORLD);
173  MPI_Pcontrol(-1,"SpGEMM_Synch_left");
174  t2 = MPI_Wtime();
175  if(myrank == 0)
176  {
177  cout<<"Synchronous multiplications (left evaluate) finished"<<endl;
178  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
179  }
180  inputA.clear();
181  inputA.close();
182  inputS.clear();
183  inputS.close();
184  inputST.clear();
185  inputST.clear();
186  }
187  MPI_Finalize();
188  return 0;
189 }
190