File: ssp.cpp

package info (click to toggle)
combblas 2.0.0-7
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 190,488 kB
  • sloc: cpp: 55,918; ansic: 25,134; sh: 3,691; makefile: 548; csh: 66; python: 49; perl: 21
file content (147 lines) | stat: -rw-r--r-- 3,567 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
#include <mpi.h>
#include <sys/time.h> 
#include <iostream>
#include <stdio.h>
#include <functional>
#include <algorithm>
#include <vector>
#include <sstream>
#include <ctime>

#include "../../trunk/CombBLAS/CombBLAS.h"

#include "Node.h"

using namespace std;

class BinaryOp {
  public: 
    Node operator() (Node a, Node b) {
      if (a.dist >= b.dist)
        return Node(b);
      else 
        return Node(a);
    }
};

class DoOp {
  public: 
    bool operator() (Node a, Node b) {
      if (a.dist >= b.dist)
        return false;
      else
        return true;
    }
};

int main(int argc, char* argv[])
{
  MPI::Init(argc, argv);
  int nprocs = MPI::COMM_WORLD.Get_size();
  int myrank = MPI::COMM_WORLD.Get_rank();
  extern MPI_Op staticMPIop;

  {
    // int id, double dist, int parent
    MPI::Datatype types[3] = {MPI::INT, MPI::DOUBLE, MPI::INT};
    int lengths[3] = {1, 1, 1};
    Node n;
    MPI::Aint disp[3];
    disp[0] = MPI::Get_address(&n.id) - MPI::Get_address(&n);
    disp[1] = MPI::Get_address(&n.dist) - MPI::Get_address(&n);
    disp[2] = MPI::Get_address(&n.parent) - MPI::Get_address(&n);

    Node_MPI_datatype = MPI::Datatype::Create_struct(3, lengths, disp, types);
    Node_MPI_datatype.Commit();
    MPI_Op_create(apply, true , &staticMPIop);
  }

  {
    if (argc != 4) {
      cout << endl << "Require 3 args..." << endl <<
        "fileName startV testV" << endl;
      MPI::Finalize();
      return -1;
    }

    char* fileName = argv[1];
    stringstream sstr(argv[2]);
    int startVert;
    sstr >> startVert;
    stringstream sstr2(argv[3]);
    int testVert;
    sstr2 >> testVert;

    if (myrank == 0)
      cout << "startV: " << startVert << endl;

    MPI::COMM_WORLD.Barrier();

    // the graph
    SpParMat<int, double, SpDCCols <int, double> > G;
    G.ReadDistribute(fileName, 0);
    int numVerts = G.getncol();
    if (myrank == 0)
      cout << "numVerts: " << numVerts << endl;

    if (startVert > numVerts || startVert <= 0) {
      cout << "Invalid start vertex id." << endl;
      return -1;
    }

    G.Transpose();
    Node zero(double(0), -1);

    time_t startTime, endTime;
    double elapsedTime;
    if (myrank == 0) {
      startTime = time(NULL);
      cout << "start computing" << endl;
    }

    int iteration;
    bool finished = false;

    FullyDistVec<int, Node> result(G.getcommgrid(), G.getncol(), Node());
    FullyDistSpVec<int, Node> frontier(G.getcommgrid(), G.getncol());

    frontier.SetElement(startVert - 1, zero);
    frontier.setNumToInd();

    BinaryOp binaryOp;
    DoOp doOp;

    frontier = EWiseApply<Node>(frontier, result, binaryOp, doOp, false, Node());
    result.EWiseApply(frontier, binaryOp, false, Node());

    for(iteration = 1; iteration < numVerts; iteration++) {
      frontier = SpMV<SPSRing>(G, frontier);
      frontier.setNumToInd();
      frontier = EWiseApply<Node>(frontier, result, binaryOp, doOp, false, Node());
      if (frontier.getnnz() == 0) {
        finished = true;
        break;
      }
      result.EWiseApply(frontier, binaryOp, false, Node());
    }

    Node res = result[testVert - 1];
    if (myrank == 0) {
      endTime = time(NULL);
      elapsedTime = difftime(endTime, startTime);
      if(finished) {
        cout << "finished" << endl;
        cout << res << endl;
      } else {
        cout << "negative loop" << endl;
      }
      cout << "number of iterations: " << iteration << endl;
      cout << "running time: " << elapsedTime << "s" << endl;
    }
    // G.Transpose();
  }

  MPI::Finalize();
  return 0;
}