File: basic.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 (117 lines) | stat: -rw-r--r-- 3,611 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
#include "../Serialize.h"
#include <sstream>
#include <mpi.h>
#include <cassert>
#include <chrono>
#include "boost/serialization/vector.hpp"
#include "boost/serialization/set.hpp"
#include "boost/archive/binary_iarchive.hpp"
#include "boost/archive/binary_oarchive.hpp"

#define NUM_SEND 1000000

using namespace std;

template<class T>
void test_send(const T& t) {
  stringbuf ss;
  chrono::time_point<chrono::high_resolution_clock> start = chrono::high_resolution_clock::now();
  encode(ss, t);
  chrono::time_point<chrono::high_resolution_clock> end = chrono::high_resolution_clock::now();
  chrono::duration<double> dur = end - start;
  double rate = ss.str().size() / dur.count() / 1024 / 1024;
  cout << ss.str().size() << " bytes at " << rate << " MB/s encode" << endl;

  MPI_Send(ss.str().data(), ss.str().size(), MPI_BYTE, 1, 0, MPI_COMM_WORLD);
}

template<class T>
void test_send_boost(const T& t) {
  stringbuf ss;
  boost::archive::binary_oarchive oa(ss);
  chrono::time_point<chrono::high_resolution_clock> start = chrono::high_resolution_clock::now();
  oa << t;
  chrono::time_point<chrono::high_resolution_clock> end = chrono::high_resolution_clock::now();
  chrono::duration<double> dur = end - start;
  double rate = ss.str().size() / dur.count() / 1024 / 1024;
  cout << ss.str().size() << " bytes at " << rate << " MB/s encode (boost)" << endl;

  MPI_Send(ss.str().data(), ss.str().size(), MPI_BYTE, 1, 1, MPI_COMM_WORLD);
}

template<class T>
void test_recv(const T& t) {
  // find out how much data we are recieving
  MPI_Status status;
  int count;
  MPI_Probe(0, 0, MPI_COMM_WORLD, &status);
  MPI_Get_count(&status, MPI_BYTE, &count);

  // recieve data
  string s(count, 0);
  MPI_Recv((void *)s.data(), count, MPI_BYTE, 0, 0, MPI_COMM_WORLD, &status);

  // decode data
  stringbuf iss(s);
  T tt;

  chrono::time_point<chrono::high_resolution_clock> start = chrono::high_resolution_clock::now();
  decode(iss, tt);
  chrono::time_point<chrono::high_resolution_clock> end = chrono::high_resolution_clock::now();
  chrono::duration<double> dur = end - start;
  auto rate = count / dur.count() / 1024 / 1024;
  cout << rate << " MB/s decode" << endl;

  if (t != tt) throw;
}

template<class T>
void test_recv_boost(const T& t) {
  // find out how much data we are recieving
  MPI_Status status;
  int count;
  MPI_Probe(0, 1, MPI_COMM_WORLD, &status);
  MPI_Get_count(&status, MPI_BYTE, &count);

  // recieve data
  string s(count, 0);
  MPI_Recv((void *)s.data(), count, MPI_BYTE, 0, 1, MPI_COMM_WORLD, &status);

  // decode data
  stringbuf iss(s);
  boost::archive::binary_iarchive oa(iss);
  T tt;

  chrono::time_point<chrono::high_resolution_clock> start = chrono::high_resolution_clock::now();
  oa >> tt;
  chrono::time_point<chrono::high_resolution_clock> end = chrono::high_resolution_clock::now();
  chrono::duration<double> dur = end - start;
  auto rate = count / dur.count() / 1024 / 1024;
  cout << rate << " MB/s decode (boost)" << endl;

  if (t != tt) throw;
}

int main(int argc, char* argv[])
{
	int nprocs, myrank;
	MPI_Init(&argc, &argv);
	MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
	MPI_Comm_rank(MPI_COMM_WORLD,&myrank);

  vector<set<int>> myVector;
  for(int i = 0; i < NUM_SEND; i++) {
    myVector.emplace_back(set<int>{i,i,i,i,i,i});
  }

  if (myrank == 0) {
    cout << "Testing " << "std::vector<std::set<int>>" << " size: " << myVector.size() << endl;
    test_send(myVector);
    test_send_boost(myVector);
  } else if (myrank == 1) {
    test_recv(myVector);
    test_recv_boost(myVector);
    cout << "Test passed" << endl;
  }
  MPI_Finalize();
}