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