COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
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/CombBLAS.h"
37
38using namespace std;
39using 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"
45template <class NT>
46class PSpMat
47{
48public:
51};
52
53
54int main(int argc, char* argv[])
55{
56 int nprocs, myrank;
57 MPI_Init(&argc, &argv);
60
61 if(argc < 4)
62 {
63 if(myrank == 0)
64 {
65 cout << "Usage: ./Galerkin <Matrix> <S> <STranspose>" << endl;
66 cout << "<Matrix>,<S>,<STranspose> are absolute addresses, and files should be in triples format" << endl;
67 }
69 return -1;
70 }
71 {
72 string Aname(argv[1]);
73 string Sname(argv[2]);
74 string STname(argv[3]);
75
78
79 PSpMat<double>::MPI_DCCols A, S, ST; // construct objects
80
81 A.ReadDistribute(Aname, 0);
82 S.ReadDistribute(Sname, 0);
83 ST.ReadDistribute(STname, 0);
84 SpParHelper::Print("Data read\n");
85
86 // force the calling of C's destructor
87 {
90 SpParHelper::Print("Warmed up for DoubleBuff (right evaluate)\n");
91 }
93 MPI_Pcontrol(1,"SpGEMM_DoubleBuff_right");
94 double t1 = MPI_Wtime(); // initilize (wall-clock) timer
95 for(int i=0; i<ITERATIONS; i++)
96 {
99 }
101 double t2 = MPI_Wtime();
102 MPI_Pcontrol(-1,"SpGEMM_DoubleBuff_right");
103 if(myrank == 0)
104 {
105 cout<<"Double buffered multiplications (right evaluate) finished"<<endl;
106 printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
107 }
108
109 // force the calling of C's destructor
110 {
113 SpParHelper::Print("Warmed up for DoubleBuff (left evaluate)\n");
114 }
116 MPI_Pcontrol(1,"SpGEMM_DoubleBuff_left");
117 t1 = MPI_Wtime(); // initilize (wall-clock) timer
118 for(int i=0; i<ITERATIONS; i++)
119 {
122 }
124 t2 = MPI_Wtime();
125 MPI_Pcontrol(-1,"SpGEMM_DoubleBuff_left");
126 if(myrank == 0)
127 {
128 cout<<"Double buffered multiplications (left evaluate) finished"<<endl;
129 printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
130 }
131
132 // force the calling of C's destructor
133 {
136 }
137 SpParHelper::Print("Warmed up for Synch (right evaluate)\n");
139 MPI_Pcontrol(1,"SpGEMM_Synch_right");
140 t1 = MPI_Wtime(); // initilize (wall-clock) timer
141 for(int i=0; i<ITERATIONS; i++)
142 {
145 }
147 MPI_Pcontrol(-1,"SpGEMM_Synch_right");
148 t2 = MPI_Wtime();
149 if(myrank == 0)
150 {
151 cout<<"Synchronous multiplications (right evaluate) finished"<<endl;
152 printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
153 }
154
155 // force the calling of C's destructor
156 {
159 }
160 SpParHelper::Print("Warmed up for Synch (left evaluate)\n");
162 MPI_Pcontrol(1,"SpGEMM_Synch_left");
163 t1 = MPI_Wtime(); // initilize (wall-clock) timer
164 for(int i=0; i<ITERATIONS; i++)
165 {
168 }
170 MPI_Pcontrol(-1,"SpGEMM_Synch_left");
171 t2 = MPI_Wtime();
172 if(myrank == 0)
173 {
174 cout<<"Synchronous multiplications (left evaluate) finished"<<endl;
175 printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
176 }
177 }
178 MPI_Finalize();
179 return 0;
180}
181
int main()
Definition Driver.cpp:12
#define ITERATIONS
Definition Galerkin.cpp:40
SpParMat< int, NT, DCCols > MPI_DCCols
Definition Galerkin.cpp:50
SpDCCols< int, NT > DCCols
Definition Galerkin.cpp:49
static void Print(const std::string &s)
int nprocs
Definition comms.cpp:55
double A
double C
Definition options.h:15
double D
Definition options.h:15