COMBINATORIAL_BLAS  1.6
GalerkinNew.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 #define ITERATIONS 10
41 
42 // Simple helper class for declarations: Just the numerical type is templated
43 // The index type and the sequential matrix type stays the same for the whole code
44 // In this case, they are "int" and "SpDCCols"
45 template <class NT>
46 class PSpMat
47 {
48 public:
51 };
52 
53 
54 int main(int argc, char* argv[])
55 {
56  int nprocs, myrank;
57  MPI_Init(&argc, &argv);
58  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
59  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
60 
61  if(argc < 5)
62  {
63  if(myrank == 0)
64  {
65  cout << "Usage: ./GalerkinNew <Matrix> <OffDiagonal> <Diagonal> <T(right hand side restriction matrix)>" << endl;
66  cout << "<Matrix> <OffDiagonal> <Diagonal> <T> are absolute addresses, and files should be in triples format" << endl;
67  cout << "Example: ./GalerkinNew TESTDATA/grid3d_k5.txt TESTDATA/offdiag_grid3d_k5.txt TESTDATA/diag_grid3d_k5.txt TESTDATA/restrict_T_grid3d_k5.txt" << endl;
68  }
69  MPI_Finalize();
70  return -1;
71  }
72  {
73  string Aname(argv[1]);
74  string Aoffd(argv[2]);
75  string Adiag(argv[3]);
76  string Tname(argv[4]);
77 
78  // A = L+D
79  // A*T = L*T + D*T;
80  // S*(A*T) = S*L*T + S*D*T;
81  ifstream inputD(Adiag.c_str());
82 
83  MPI_Barrier(MPI_COMM_WORLD);
84  typedef PlusTimesSRing<double, double> PTDD;
85  shared_ptr<CommGrid> fullWorld;
86  fullWorld.reset( new CommGrid(MPI_COMM_WORLD, 0, 0) );
87 
88  PSpMat<double>::MPI_DCCols A(fullWorld); // construct objects
89  PSpMat<double>::MPI_DCCols L(fullWorld);
90  PSpMat<double>::MPI_DCCols T(fullWorld);
91  FullyDistVec<int,double> dvec(fullWorld);
92 
93  // For matrices, passing the file names as opposed to fstream objects
94  A.ReadDistribute(Aname, 0);
95  L.ReadDistribute(Aoffd, 0);
96  T.ReadDistribute(Tname, 0);
97  dvec.ReadDistribute(inputD,0);
98  SpParHelper::Print("Data read\n");
99 
101  S.Transpose();
102 
103  // force the calling of C's destructor; warm up instruction cache - also check correctness
104  {
105  PSpMat<double>::MPI_DCCols AT = PSpGEMM<PTDD>(A, T);
106  PSpMat<double>::MPI_DCCols SAT = PSpGEMM<PTDD>(S, AT);
107 
108  PSpMat<double>::MPI_DCCols LT = PSpGEMM<PTDD>(L, T);
109  PSpMat<double>::MPI_DCCols SLT = PSpGEMM<PTDD>(S, LT);
111  SD.DimApply(Column, dvec, multiplies<double>()); // scale columns of S to get SD
112  PSpMat<double>::MPI_DCCols SDT = PSpGEMM<PTDD>(SD, T);
113  SLT += SDT; // now this is SAT
114 
115  if(SLT == SAT)
116  {
117  SpParHelper::Print("Splitting approach is correct\n");
118  }
119  else
120  {
121  SpParHelper::Print("Error in splitting, go fix it\n");
122  SLT.PrintInfo();
123  SAT.PrintInfo();
124  //SLT.SaveGathered("SLT.txt");
125  //SAT.SaveGathered("SAT.txt");
126  }
127  }
128  MPI_Barrier(MPI_COMM_WORLD);
129  double t1 = MPI_Wtime(); // initilize (wall-clock) timer
130  for(int i=0; i<ITERATIONS; i++)
131  {
132  PSpMat<double>::MPI_DCCols AT = PSpGEMM<PTDD>(A, T);
133  PSpMat<double>::MPI_DCCols SAT = PSpGEMM<PTDD>(S, AT);
134  }
135  MPI_Barrier(MPI_COMM_WORLD);
136  double t2 = MPI_Wtime();
137  if(myrank == 0)
138  {
139  cout<<"Full restriction (without splitting) finished"<<endl;
140  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
141  }
142 
143  MPI_Barrier(MPI_COMM_WORLD);
144  t1 = MPI_Wtime(); // initilize (wall-clock) timer
145  for(int i=0; i<ITERATIONS; i++)
146  {
147 
148  PSpMat<double>::MPI_DCCols LT = PSpGEMM<PTDD>(L, T);
149  PSpMat<double>::MPI_DCCols SLT = PSpGEMM<PTDD>(S, LT);
151  SD.DimApply(Column, dvec, multiplies<double>()); // scale columns of S to get SD
152  PSpMat<double>::MPI_DCCols SDT = PSpGEMM<PTDD>(SD, T);
153  SLT += SDT;
154  }
155  MPI_Barrier(MPI_COMM_WORLD);
156  t2 = MPI_Wtime();
157  if(myrank == 0)
158  {
159  cout<<"Full restriction (with splitting) finished"<<endl;
160  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
161  }
162  inputD.clear();inputD.close();
163  }
164  MPI_Finalize();
165  return 0;
166 }
167 
#define ITERATIONS
Definition: GalerkinNew.cpp:40
void DimApply(Dim dim, const FullyDistVec< IT, NT > &v, _BinaryOperation __binary_op)
Definition: SpParMat.cpp:704
void ReadDistribute(const std::string &filename, int master, bool nonum, HANDLER handler, bool transpose=false, bool pario=false)
Definition: SpParMat.cpp:3560
double A
int main(int argc, char *argv[])
Definition: GalerkinNew.cpp:54
SpDCCols< int, NT > DCCols
Definition: GalerkinNew.cpp:49
Definition: CCGrid.h:4
std::ifstream & ReadDistribute(std::ifstream &infile, int master, HANDLER handler)
SpParMat< int, NT, DCCols > MPI_DCCols
Definition: GalerkinNew.cpp:50