COMBINATORIAL_BLAS  1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GalerkinNew.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 < 5)
61  {
62  if(myrank == 0)
63  {
64  cout << "Usage: ./GalerkinNew <Matrix> <OffDiagonal> <Diagonal> <T(right hand side restriction matrix)>" << endl;
65  cout << "<Matrix> <OffDiagonal> <Diagonal> <T> are absolute addresses, and files should be in triples format" << endl;
66  cout << "Example: ./GalerkinNew TESTDATA/grid3d_k5.txt TESTDATA/offdiag_grid3d_k5.txt TESTDATA/diag_grid3d_k5.txt TESTDATA/restrict_T_grid3d_k5.txt" << endl;
67  }
68  MPI_Finalize();
69  return -1;
70  }
71  {
72  string Aname(argv[1]);
73  string Aoffd(argv[2]);
74  string Adiag(argv[3]);
75  string Tname(argv[4]);
76 
77  // A = L+D
78  // A*T = L*T + D*T;
79  // S*(A*T) = S*L*T + S*D*T;
80  ifstream inputD(Adiag.c_str());
81 
82  MPI_Barrier(MPI_COMM_WORLD);
83  typedef PlusTimesSRing<double, double> PTDD;
84 
85  PSpMat<double>::MPI_DCCols A, L, T; // construct objects
87 
88  // For matrices, passing the file names as opposed to fstream objects
89  A.ReadDistribute(Aname, 0);
90  L.ReadDistribute(Aoffd, 0);
91  T.ReadDistribute(Tname, 0);
92  dvec.ReadDistribute(inputD,0);
93  SpParHelper::Print("Data read\n");
94 
96  S.Transpose();
97 
98  // force the calling of C's destructor; warm up instruction cache - also check correctness
99  {
100  PSpMat<double>::MPI_DCCols AT = PSpGEMM<PTDD>(A, T);
101  PSpMat<double>::MPI_DCCols SAT = PSpGEMM<PTDD>(S, AT);
102 
103  PSpMat<double>::MPI_DCCols LT = PSpGEMM<PTDD>(L, T);
104  PSpMat<double>::MPI_DCCols SLT = PSpGEMM<PTDD>(S, LT);
106  SD.DimApply(Column, dvec, multiplies<double>()); // scale columns of S to get SD
107  PSpMat<double>::MPI_DCCols SDT = PSpGEMM<PTDD>(SD, T);
108  SLT += SDT; // now this is SAT
109 
110  if(SLT == SAT)
111  {
112  SpParHelper::Print("Splitting approach is correct\n");
113  }
114  else
115  {
116  SpParHelper::Print("Error in splitting, go fix it\n");
117  SLT.PrintInfo();
118  SAT.PrintInfo();
119  //SLT.SaveGathered("SLT.txt");
120  //SAT.SaveGathered("SAT.txt");
121  }
122  }
123  MPI_Barrier(MPI_COMM_WORLD);
124  double t1 = MPI_Wtime(); // initilize (wall-clock) timer
125  for(int i=0; i<ITERATIONS; i++)
126  {
127  PSpMat<double>::MPI_DCCols AT = PSpGEMM<PTDD>(A, T);
128  PSpMat<double>::MPI_DCCols SAT = PSpGEMM<PTDD>(S, AT);
129  }
130  MPI_Barrier(MPI_COMM_WORLD);
131  double t2 = MPI_Wtime();
132  if(myrank == 0)
133  {
134  cout<<"Full restriction (without splitting) finished"<<endl;
135  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
136  }
137 
138  MPI_Barrier(MPI_COMM_WORLD);
139  t1 = MPI_Wtime(); // initilize (wall-clock) timer
140  for(int i=0; i<ITERATIONS; i++)
141  {
142 
143  PSpMat<double>::MPI_DCCols LT = PSpGEMM<PTDD>(L, T);
144  PSpMat<double>::MPI_DCCols SLT = PSpGEMM<PTDD>(S, LT);
146  SD.DimApply(Column, dvec, multiplies<double>()); // scale columns of S to get SD
147  PSpMat<double>::MPI_DCCols SDT = PSpGEMM<PTDD>(SD, T);
148  SLT += SDT;
149  }
150  MPI_Barrier(MPI_COMM_WORLD);
151  t2 = MPI_Wtime();
152  if(myrank == 0)
153  {
154  cout<<"Full restriction (with splitting) finished"<<endl;
155  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
156  }
157  inputD.clear();inputD.close();
158  }
159  MPI_Finalize();
160  return 0;
161 }
162