COMBINATORIAL_BLAS  1.6
VectorIndexing.cpp
Go to the documentation of this file.
1 #include <mpi.h>
2 #include <sys/time.h>
3 #include <iostream>
4 #include <functional>
5 #include <algorithm>
6 #include <vector>
7 #include <sstream>
8 #include "CombBLAS/CombBLAS.h"
9 
10 using namespace std;
11 using namespace combblas;
12 
13 int main(int argc, char* argv[])
14 {
15  int nprocs, myrank;
16  MPI_Init(&argc, &argv);
17  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
18  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
19 
20  if(argc < 4)
21  {
22  if(myrank == 0)
23  {
24  cout << "Usage: ./VectorIndexing <BASEADDRESS> <VectorOne> <VectorTwo>" << endl;
25  cout << "Example: ./VectorIndexing ../TESTDATA sp10outta100.txt sp30outta100.txt" << endl;
26  cout << "Input files should be under <BASEADDRESS> in tuples format" << endl;
27  }
28  MPI_Finalize();
29  return -1;
30  }
31  {
32  string directory(argv[1]);
33  string vec1name(argv[2]);
34  string vec2name(argv[3]);
35  vec1name = directory+"/"+vec1name;
36  vec2name = directory+"/"+vec2name;
37 
38  ifstream inputvec1(vec1name.c_str());
39  ifstream inputvec2(vec2name.c_str());
40 
41  if(myrank == 0)
42  {
43  if(inputvec1.fail() || inputvec2.fail())
44  {
45  cout << "One of the input vector files do not exist, aborting" << endl;
46  MPI_Abort(MPI_COMM_WORLD, NOFILE);
47  return -1;
48  }
49  }
50  MPI_Barrier(MPI_COMM_WORLD);
51  FullyDistSpVec<int,int> vec1, vec2;
52 
53  vec1.ReadDistribute(inputvec1, 0);
54  vec2.ReadDistribute(inputvec2, 0);
55 
56  vec1.PrintInfo("vec1");
57  vec1.DebugPrint();
58  vec2.PrintInfo("vec2");
59  vec2.DebugPrint();
60 
62  dvec.iota(100, 1001); // with 100 entries, first being 1001
63  dvec.PrintInfo("dvec");
64  dvec.DebugPrint();
65 
66  auto subvec1 = dvec(vec1);
67  subvec1.DebugPrint();
68 
69  auto subvec2 = dvec(vec2);
70  subvec2.DebugPrint();
71 
72  FullyDistSpVec<int,int> vecA(12);
73  for(int i=0; i<12; i+=3) vecA.SetElement(i,i);
74  MPI_Barrier(MPI_COMM_WORLD);
75  vecA.DebugPrint();
77  dvecA.iota(12, 0); // with 12 entries, first being 0
78 
79  auto subvecA = dvecA(vecA);
80  subvecA.DebugPrint();
81 
82 
83  inputvec1.clear();
84  inputvec1.close();
85  inputvec2.clear();
86  inputvec2.close();
87  }
88  MPI_Finalize();
89  return 0;
90 }
void SetElement(IT indx, NT numx)
Indexing is performed 0-based.
int main(int argc, char *argv[])
#define NOFILE
Definition: SpDefs.h:75
void PrintInfo(std::string vecname) const
void iota(IT globalsize, NT first)
std::ifstream & ReadDistribute(std::ifstream &infile, int master, HANDLER handler)
Totally obsolete version that only accepts an ifstream object and ascii files.
void PrintInfo(std::string vectorname) const
Definition: CCGrid.h:4