Zoltan2
taskMappingTest3.cpp
Go to the documentation of this file.
3 #include "Tpetra_MultiVector_decl.hpp"
4 
5 #include <GeometricGenerator.hpp>
6 #include <string>
7 #include "Teuchos_XMLParameterListHelpers.hpp"
8 //#include "Teuchos_MPIComm.hpp"
9 #define partDIM 3
10 #define procDIM 3
11 #define nProcs 200;
12 #define nParts 200;
13 
14 
16  const string& s,
17  const string& delimiters = " \f\n\r\t\v" )
18 {
19  return s.substr( 0, s.find_last_not_of( delimiters ) + 1 );
20 }
21 
23  const string& s,
24  const string& delimiters = " \f\n\r\t\v" )
25 {
26  return s.substr( s.find_first_not_of( delimiters ) );
27 }
28 
29 string trim_copy(
30  const string& s,
31  const string& delimiters = " \f\n\r\t\v" )
32 {
33  return trim_left_copy( trim_right_copy( s, delimiters ), delimiters );
34 }
35 
36 const char param_comment = '#';
37 void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP<const Teuchos::Comm<int> > & comm){
38  std::string input = "";
39  char inp[25000];
40  for(int i = 0; i < 25000; ++i){
41  inp[i] = 0;
42  }
43 
44  bool fail = false;
45  if(comm->getRank() == 0){
46 
47  std::fstream inParam(paramFileName.c_str());
48  if (inParam.fail())
49  {
50  fail = true;
51  }
52  if(!fail)
53  {
54  std::string tmp = "";
55  getline (inParam,tmp);
56  while (!inParam.eof()){
57  if(tmp != ""){
58  tmp = trim_copy(tmp);
59  if(tmp != ""){
60  input += tmp + "\n";
61  }
62  }
63  getline (inParam,tmp);
64  }
65  inParam.close();
66  for (size_t i = 0; i < input.size(); ++i){
67  inp[i] = input[i];
68  }
69  }
70  }
71 
72 
73 
74  int size = input.size();
75  if(fail){
76  size = -1;
77  }
78  comm->broadcast(0, sizeof(int), (char*) &size);
79  if(size == -1){
80  throw "File " + paramFileName + " cannot be opened.";
81  }
82  comm->broadcast(0, size, inp);
83  std::istringstream inParam(inp);
84  std::string str;
85  getline (inParam,str);
86  while (!inParam.eof()){
87  if(str[0] != param_comment){
88  size_t pos = str.find('=');
89  if(pos == std::string::npos){
90  throw "Invalid Line:" + str + " in parameter file";
91  }
92  string paramname = trim_copy(str.substr(0,pos));
93  string paramvalue = trim_copy(str.substr(pos + 1));
94  geoparams.set(paramname, paramvalue);
95  }
96  getline (inParam,str);
97  }
98 }
99 string convert_to_string(char *args){
100  string tmp = "";
101  for(int i = 0; args[i] != 0; i++)
102  tmp += args[i];
103  return tmp;
104 }
105 bool getArgumentValue(string &argumentid, double &argumentValue, string argumentline){
106  std::stringstream stream(std::stringstream::in | std::stringstream::out);
107  stream << argumentline;
108  getline(stream, argumentid, '=');
109  if (stream.eof()){
110  return false;
111  }
112  stream >> argumentValue;
113  return true;
114 }
115 
116 template <typename part_t>
118  int argc,
119  char **argv,
120  std::string &procF,
121  part_t &nx,
122  part_t &ny,
123  part_t &nz, bool &divide_prime, int &ranks_per_node, string &taskGraphFile, string &taskCoordFile){
124 
125  bool isprocset = false;
126  int ispartset = 0;
127 
128  for(int i = 0; i < argc; ++i){
129  string tmp = convert_to_string(argv[i]);
130  string tmp2 = "";
131  string identifier = "";
132  double fval = -1;
133  if(!getArgumentValue(identifier, fval, tmp)) continue;
134 
135  if(identifier == "PROC"){
136  std::stringstream stream(std::stringstream::in | std::stringstream::out);
137  stream << tmp;
138  getline(stream, procF, '=');
139 
140  stream >> procF;
141  isprocset = true;
142  }
143  else if(identifier == "NX"){
144  std::stringstream stream(std::stringstream::in | std::stringstream::out);
145  stream << tmp;
146  getline(stream, tmp2, '=');
147 
148  stream >> nx;
149  ispartset++;
150  }
151  else if(identifier == "NY"){
152  std::stringstream stream(std::stringstream::in | std::stringstream::out);
153  stream << tmp;
154  getline(stream, tmp2, '=');
155 
156  stream >> ny;
157  ispartset++;
158  }
159  else if(identifier == "NZ"){
160  std::stringstream stream(std::stringstream::in | std::stringstream::out);
161  stream << tmp;
162  getline(stream, tmp2, '=');
163 
164  stream >> nz;
165  ispartset++;
166  }
167  else if(identifier == "DP"){
168  std::stringstream stream(std::stringstream::in | std::stringstream::out);
169  stream << tmp;
170  getline(stream, tmp2, '=');
171  int val;
172  stream >> val;
173  if (val) divide_prime = true;
174  ispartset++;
175  }
176  else if(identifier == "RPN"){
177  std::stringstream stream(std::stringstream::in | std::stringstream::out);
178  stream << tmp;
179  getline(stream, tmp2, '=');
180  stream >> ranks_per_node;
181  ispartset++;
182  } else if(identifier == "TG"){
183  std::stringstream stream(std::stringstream::in | std::stringstream::out);
184  stream << tmp;
185  getline(stream, taskGraphFile, '=');
186 
187  stream >> taskGraphFile;
188  }
189  else if(identifier == "TC"){
190  std::stringstream stream(std::stringstream::in | std::stringstream::out);
191  stream << tmp;
192  getline(stream, taskCoordFile, '=');
193  stream >> taskCoordFile;
194  }
195 
196  else {
197  throw "Invalid argument at " + tmp;
198  }
199 
200  }
201  if(!(ispartset >= 3&& isprocset)){
202  throw "(PROC && PART) are mandatory arguments.";
203  }
204 
205 }
206 int main(int argc, char *argv[]){
207 
208  typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tMVector_t;
209  typedef Zoltan2::XpetraMultiVectorAdapter<tMVector_t> inputAdapter_t;
211 
212  Teuchos::GlobalMPISession session(&argc, &argv);
213  //if (argc != 3){
214  // cout << "Usage: " << argv[0] << " PART=partGeoParams.txt PROC=procGeoParams.txt" << endl;
215  // exit(1);
216  //}
217  part_t numParts = 0;
218  zscalar_t **partCenters = NULL;
219  int coordDim = 0;
220 
221  part_t numProcs = 0;
222  zscalar_t **procCoordinates = NULL;
223  int procDim = 0;
224 
225 
226  string taskGraphFile = "";
227  string taskCoordFile = "";
228  bool divide_prime = false;
229  part_t jobX = 1, jobY = 1, jobZ = 1;
230  string procfile = "";
231  int rank_per_node = 16;
232 
233  const RCP<Comm<int> > commN;
234  RCP<Comm<int> >comm = Teuchos::rcp_const_cast<Comm<int> >
235  (Teuchos::DefaultComm<int>::getDefaultSerialComm(commN));
236 
237  part_t *task_communication_xadj_ = NULL;
238  part_t *task_communication_adj_ = NULL;
239  zscalar_t *task_communication_adjw_ = NULL;
240  try {
241 
242  getArgVals<part_t>(
243  argc,
244  argv,
245  procfile ,
246  jobX, jobY, jobZ, divide_prime, rank_per_node, taskGraphFile, taskCoordFile);
247 
248  coordDim = 3;
249  procDim = 3;
250  numParts = jobZ*jobY*jobX;
251 
252  //numProcs = numParts;
253  //cout << "part:" << numParts << " proc:" << procfile << endl;
254  if (taskGraphFile == "" || taskCoordFile == "")
255  {
256 
257  partCenters = new zscalar_t * [coordDim];
258  for(int i = 0; i < coordDim; ++i){
259  partCenters[i] = new zscalar_t[numParts];
260  }
261 
262 
263  task_communication_xadj_ = new part_t [numParts+1];
264  task_communication_adj_ = new part_t [numParts * 6];
265 
266  int prevNCount = 0;
267  task_communication_xadj_[0] = 0;
268  for (part_t i = 0; i < numParts; ++i) {
269  int x = i % jobX;
270  int y = (i / (jobX)) % jobY;
271  int z = (i / (jobX)) / jobY;
272  partCenters[0][i] = x;
273  partCenters[1][i] = y;
274  partCenters[2][i] = z;
275 
276  if (x > 0){
277  task_communication_adj_[prevNCount++] = i - 1;
278  }
279  if (x < jobX - 1){
280  task_communication_adj_[prevNCount++] = i + 1;
281  }
282  if (y > 0){
283  task_communication_adj_[prevNCount++] = i - jobX;
284  }
285  if (y < jobY - 1){
286  task_communication_adj_[prevNCount++] = i + jobX;
287  }
288  if (z > 0){
289  task_communication_adj_[prevNCount++] = i - jobX * jobY;
290  }
291  if (z < jobZ - 1){
292  task_communication_adj_[prevNCount++] = i + jobX * jobY;
293  }
294  task_communication_xadj_[i+1] = prevNCount;
295  }
296  }
297  else {
298  int ne = 0;
299  FILE *f2 = fopen(taskGraphFile.c_str(), "rb");
300 
301  fread(&numParts,sizeof(int),1,f2); // write 10 bytes to our buffer
302  fread(&ne,sizeof(int),1,f2); // write 10 bytes to our buffer
303 
304 
305  std::cout << "numParts:" << numParts << " ne:" << ne << std::endl;
306 
307  task_communication_xadj_ = new part_t [numParts+1];
308  task_communication_adj_ = new part_t [ne];
309  task_communication_adjw_ = new zscalar_t [ne];
310 
311  fread((void *)task_communication_xadj_,sizeof(int),numParts + 1,f2); // write 10 bytes to our buffer
312  fread((void *)task_communication_adj_,sizeof(int),ne ,f2); // write 10 bytes to our buffer
313  fread((void *)task_communication_adjw_,sizeof(double),ne,f2); // write 10 bytes to our buffer
314  fclose(f2);
315 
316  f2 = fopen(taskCoordFile.c_str(), "rb");
317  fread((void *)&numParts,sizeof(int),1,f2); // write 10 bytes to our buffer
318  fread((void *)&coordDim,sizeof(int),1,f2); // write 10 bytes to our buffer
319 
320  std::cout << "numParts:" << numParts << " coordDim:" << coordDim << std::endl;
321 
322  partCenters = new zscalar_t * [coordDim];
323  for(int i = 0; i < coordDim; ++i){
324  partCenters[i] = new zscalar_t[numParts];
325  fread((void *) partCenters[i],sizeof(double),numParts, f2); // write 10 bytes to our buffer
326  }
327  fclose(f2);
328  }
329 
330 
331 
332 
333  {
334  std::vector < std::vector <zscalar_t> > proc_coords(procDim);
335  std::fstream m(procfile.c_str());
336  procCoordinates = new zscalar_t * [procDim];
337 
338  part_t i = 0;
339  zscalar_t a,b, c;
340  m >> a >> b >> c;
341  while(!m.eof()){
342  proc_coords[0].push_back(a);
343  proc_coords[1].push_back(b);
344  proc_coords[2].push_back(c);
345  ++i;
346  m >> a >> b >> c;
347  }
348 
349  m.close();
350  numProcs = i;
351  for(int ii = 0; ii < procDim; ++ii){
352  procCoordinates[ii] = new zscalar_t[numProcs];
353  for (int j = 0; j < numProcs; ++j){
354  procCoordinates[ii][j] = proc_coords[ii][j];
355  }
356  }
357  }
358  std::cout << "numProcs:" << numProcs << std::endl;
359 
360  /*
361  Zoltan2::CoordinateCommunicationModel<zscalar_t,zscalar_t,int> *cm =
362  new Zoltan2::CoordinateCommunicationModel<zscalar_t,zscalar_t,int>(
363  procDim, procCoordinates,
364  coordDim, partCenters,
365  numProcs, numParts);
366 
367  Zoltan2::Environment *env = new Zoltan2::Environment();
368  Zoltan2::CoordinateTaskMapper <inputAdapter_t, int> *ctm=
369  new Zoltan2::CoordinateTaskMapper<inputAdapter_t,int>(env, cm);
370 
371  */
372  RCP<const Teuchos::Comm<int> > tcomm = Teuchos::DefaultComm<int>::getComm();
373  part_t *proc_to_task_xadj_ = new part_t[numProcs+1];
374  part_t *proc_to_task_adj_ = new part_t[numParts];
375 /*
376  cout << "procDim:" << procDim <<
377  " numProcs:" << numProcs <<
378  " coordDim:" << coordDim <<
379  " numParts" << numParts << endl;
380 
381  for(part_t j = 0; j < numProcs; ++j){
382  cout << "proc - coord:" << j << " " << procCoordinates[0][j]<< " " << procCoordinates[1][j]<< " " << procCoordinates[2][j] << endl;
383  }
384 
385  for(part_t j = 0; j < numParts; ++j){
386  cout << "part - coord:" << j << " " << partCenters[0][j]<< " " << partCenters[1][j]<< " " << partCenters[2][j] << endl;
387  }
388 */
389  /*
390  int partArray[3];
391  partArray[0] = 8;
392  partArray[1] = 4;
393  partArray[2] = 16;
394  */
395  part_t *partArray = NULL;
396  int partArraysize = -1;
397  //part_t hopper[3];
398  //hopper[0] = 17;
399  //hopper[1] = 8;
400  //hopper[2] = 24;
401  part_t *machineDimensions = NULL;
402  //machineDimensions = hopper;
403  Zoltan2::coordinateTaskMapperInterface<part_t, zscalar_t, zscalar_t>(
404  tcomm,
405  procDim,
406  numProcs,
407  procCoordinates,
408 
409  coordDim,
410  numParts,
411  partCenters,
412 
413  task_communication_xadj_,
414  task_communication_adj_,
415  task_communication_adjw_,
416 
417  proc_to_task_xadj_, /*output*/
418  proc_to_task_adj_, /*output*/
419 
420  partArraysize,
421  partArray,
422  machineDimensions, rank_per_node, divide_prime
423  );
424 
425  if (tcomm->getRank() == 0){
426  cout << "PASS" << endl;
427  }
428  /*
429  delete ctm;
430  delete cm;
431  delete env;
432  */
433  delete [] proc_to_task_xadj_;
434  delete [] proc_to_task_adj_;
435  delete [] task_communication_xadj_;
436  delete [] task_communication_adj_;
437  delete [] task_communication_adjw_;
438 
439  for (int i = 0; i < coordDim; i++) delete [] partCenters[i];
440  delete [] partCenters;
441  for (int i = 0; i < procDim; i++) delete [] procCoordinates[i];
442  delete [] procCoordinates;
443  }
444  catch(std::string &s){
445  cerr << s << endl;
446  }
447 
448  catch(char * s){
449  cerr << s << endl;
450  }
451 }
452 
int main(int argc, char *argv[])
double zscalar_t
string convert_to_string(char *args)
common code used by tests
void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP< const Teuchos::Comm< int > > &comm)
SparseMatrixAdapter_t::part_t part_t
const char param_comment
string trim_copy(const string &s, const string &delimiters=" \\\)
An adapter for Xpetra::MultiVector.
static const std::string fail
string trim_left_copy(const string &s, const string &delimiters=" \\\)
Tpetra::MultiVector< double, int, int > tMVector_t
bool getArgumentValue(string &argumentid, double &argumentValue, string argumentline)
string trim_right_copy(const string &s, const string &delimiters=" \\\)
void getArgVals(int argc, char **argv, std::string &procF, part_t &nx, part_t &ny, part_t &nz, bool &divide_prime, int &ranks_per_node, string &taskGraphFile, string &taskCoordFile)