16 #include <Teuchos_RCP.hpp>
17 #include <Teuchos_Array.hpp>
18 #include <Teuchos_ParameterList.hpp>
19 #include "Teuchos_XMLParameterListHelpers.hpp"
21 #include "Tpetra_MultiVector.hpp"
22 #include <Tpetra_CrsGraph.hpp>
23 #include <Tpetra_Map.hpp>
34 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>
54 int nx,
int ny,
int nz,
55 int numProcs, Teuchos::RCP<
const Teuchos::Comm<int> > tcomm,
56 RCP<Zoltan2::Environment> env,
zscalar_t ** &partCenters,
zgno_t & myTasks) {
58 int rank = tcomm->getRank();
59 using namespace Teuchos;
62 zgno_t numGlobalTasks = nx*ny*nz;
65 myTasks = numGlobalTasks / numProcs;
66 zgno_t taskLeftOver = numGlobalTasks % numProcs;
67 if (
zgno_t(rank) < taskLeftOver ) ++myTasks;
69 zgno_t myTaskBegin = (numGlobalTasks / numProcs) * rank;
70 myTaskBegin += (taskLeftOver <
zgno_t(rank) ? taskLeftOver : rank);
71 zgno_t myTaskEnd = myTaskBegin + myTasks;
75 for(
int i = 0; i < coordDim; ++i) {
79 zgno_t *task_communication_xadj_ =
new zgno_t [myTasks + 1];
80 zgno_t *task_communication_adj_ =
new zgno_t [myTasks * 6];
84 task_communication_xadj_[0] = 0;
85 for (
zgno_t i = myTaskBegin; i < myTaskEnd; ++i) {
87 int y = (i / (nx)) % ny;
88 int z = (i / (nx)) / ny;
89 partCenters[0][i - myTaskBegin] = x;
90 partCenters[1][i - myTaskBegin] = y;
91 partCenters[2][i - myTaskBegin] = z;
94 task_communication_adj_[prevNCount++] = i - 1;
97 task_communication_adj_[prevNCount++] = i + 1;
100 task_communication_adj_[prevNCount++] = i - nx;
103 task_communication_adj_[prevNCount++] = i + nx;
106 task_communication_adj_[prevNCount++] = i - nx * ny;
109 task_communication_adj_[prevNCount++] = i + nx * ny;
112 task_communication_xadj_[i + 1 - myTaskBegin] = prevNCount;
117 using namespace Teuchos;
118 RCP<const mytest_map_t> map = rcp (
new mytest_map_t (numGlobalTasks,
121 Teuchos::Array<size_t> adjPerTask(myTasks);
122 for (
zgno_t lclRow = 0; lclRow < myTasks; ++lclRow)
123 adjPerTask[lclRow] = task_communication_xadj_[lclRow+1]
124 - task_communication_xadj_[lclRow];
130 for (
zgno_t lclRow = 0; lclRow < myTasks; ++lclRow) {
131 const zgno_t gblRow = map->getGlobalElement (lclRow);
132 zgno_t begin = task_communication_xadj_[lclRow];
133 zgno_t end = task_communication_xadj_[lclRow + 1];
134 const ArrayView< const zgno_t > indices(task_communication_adj_ + begin,
136 TpetraCrsGraph->insertGlobalIndices(gblRow, indices);
138 TpetraCrsGraph->fillComplete ();
140 delete [] task_communication_xadj_;
141 delete [] task_communication_adj_;
144 return TpetraCrsGraph;
148 RCP<Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> >
153 const int coord_dim = 3;
154 Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
157 Teuchos::ArrayView<const zscalar_t> a(partCenters[0], myTasks);
159 Teuchos::ArrayView<const zscalar_t> b(partCenters[1], myTasks);
161 Teuchos::ArrayView<const zscalar_t> c(partCenters[2], myTasks);
165 Teuchos::ArrayView<const zscalar_t> a;
172 RCP<mytest_tMVector_t> coords(
174 RCP<const mytest_tMVector_t> const_coords =
176 RCP <Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> > adapter(
183 int nx,
int ny,
int nz,
184 Teuchos::RCP<
const Teuchos::Comm<int> > global_tcomm)
186 Teuchos::RCP<const Teuchos::Comm<int> > tcomm = global_tcomm;
190 Teuchos::ParameterList distributed_problemParams;
193 distributed_problemParams.set(
"Machine_Optimization_Level", 10);
194 distributed_problemParams.set(
"mapping_algorithm",
"geometric");
195 distributed_problemParams.set(
"distributed_input_adapter",
true);
196 distributed_problemParams.set(
"mj_enable_rcb",
true);
197 distributed_problemParams.set(
"algorithm",
"multijagged");
198 distributed_problemParams.set(
"num_global_parts", numProcs);
200 RCP<Zoltan2::Environment> env(
202 RCP<Zoltan2::TimerManager> timer(
205 env->setTimer(timer);
211 RCP<mytest_tcrsGraph_t> distributed_tpetra_graph =
213 tcomm, env, partCenters,
215 RCP<const mytest_map_t> distributed_map =
216 distributed_tpetra_graph->getMap();
217 global_tcomm->barrier();
221 RCP<const mytest_tcrsGraph_t> const_tpetra_graph =
226 RCP<Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> >
228 partCenters, myTasks);
229 ia->setCoordinateInput(distributed_adapter.getRawPtr());
231 global_tcomm->barrier();
252 ia.getRawPtr(), &distributed_problemParams, global_tcomm);
254 partitioning_problem.
solve();
257 partitioning_problem.getSolution();
265 single_phase_mapping_solution(env, global_tcomm, 0);
266 Teuchos::ArrayView< const zgno_t> gids =
267 distributed_map->getLocalElementList();
269 ArrayRCP<int> initial_part_ids(myTasks);
270 for (
zgno_t i = 0; i < myTasks; ++i) {
271 initial_part_ids[i] = gids[i];
273 single_phase_mapping_solution.
setParts(initial_part_ids);
279 ia.getRawPtr(), &distributed_problemParams,
280 global_tcomm, &partition_solution);
287 ia.getRawPtr(), &distributed_problemParams, global_tcomm);
292 ia.getRawPtr(), &distributed_problemParams,
293 global_tcomm, &single_phase_mapping_solution);
299 distributed_map_problem_1.
solve(
true);
301 distributed_map_problem_2.solve(
true);
303 distributed_map_problem_3.solve(
true);
310 distributed_map_problem_1.getSolution();
313 distributed_map_problem_2.getSolution();
316 distributed_map_problem_3.getSolution();
318 timer->printAndResetToZero();
324 RCP<quality_t> metricObject_1 =
325 rcp(
new quality_t(ia.getRawPtr(), &distributed_problemParams,
326 global_tcomm, msoln1,
327 distributed_map_problem_1.getMachine().getRawPtr()));
330 RCP<quality_t> metricObject_2 =
331 rcp(
new quality_t(ia.getRawPtr(), &distributed_problemParams,
332 global_tcomm, msoln2,
333 distributed_map_problem_2.getMachine().getRawPtr()));
336 RCP<quality_t> metricObject_3 =
337 rcp(
new quality_t(ia.getRawPtr(), &distributed_problemParams,
338 global_tcomm, msoln3,
339 distributed_map_problem_3.getMachine().getRawPtr()));
342 if (global_tcomm->getRank() == 0) {
343 std::cout <<
"METRICS FOR THE FIRST CASE - TWO PHASE MAPPING"
345 metricObject_1->printMetrics(std::cout);
346 std::cout <<
"METRICS FOR THE SECOND CASE - TWO PHASE MAPPING"
347 <<
" - INITIAL ASSIGNMENT ARE ASSUMED TO BE A PART" << std::endl;
348 metricObject_2->printMetrics(std::cout);
349 std::cout <<
"METRICS FOR THE THIRD CASE - ONE PHASE MAPPING"
350 <<
" - EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING"
352 metricObject_3->printMetrics(std::cout);
355 for (
int i = 0; i < 3; i++)
356 delete [] partCenters[i];
357 delete [] partCenters;
363 int nx,
int ny,
int nz,
364 Teuchos::RCP<
const Teuchos::Comm<int> > global_tcomm)
367 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm =
368 Teuchos::createSerialComm<int>();
373 Teuchos::ParameterList serial_problemParams;
375 serial_problemParams.set(
"Machine_Optimization_Level", 10);
376 serial_problemParams.set(
"mapping_algorithm",
"geometric");
377 serial_problemParams.set(
"distributed_input_adapter",
false);
378 serial_problemParams.set(
"algorithm",
"multijagged");
379 serial_problemParams.set(
"num_global_parts", numProcs);
381 RCP<Zoltan2::Environment> env(
383 RCP<Zoltan2::TimerManager> timer(
386 env->setTimer(timer);
392 RCP<mytest_tcrsGraph_t> serial_tpetra_graph =
394 env, partCenters, myTasks);
395 RCP<const mytest_map_t> serial_map = serial_tpetra_graph->getMap();
396 global_tcomm->barrier();
400 RCP<const mytest_tcrsGraph_t> const_tpetra_graph =
405 RCP <Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t>> serial_adapter =
407 ia->setCoordinateInput(serial_adapter.getRawPtr());
409 global_tcomm->barrier();
437 single_phase_mapping_solution(env, global_tcomm, 0);
438 Teuchos::ArrayView< const zgno_t> gids = serial_map->getLocalElementList();
440 ArrayRCP<int> initial_part_ids(myTasks);
441 for (
zgno_t i = 0; i < myTasks; ++i) {
442 initial_part_ids[i] = gids[i];
444 single_phase_mapping_solution.
setParts(initial_part_ids);
454 ia.getRawPtr(), &serial_problemParams,
455 global_tcomm, &single_phase_mapping_solution);
461 serial_map_problem.
solve(
true);
466 serial_map_problem.getSolution();
468 timer->printAndResetToZero();
476 RCP<quality_t> metricObject_3 =
478 &serial_problemParams, serial_comm,msoln3,
479 serial_map_problem.getMachine().getRawPtr()));
481 if (global_tcomm->getRank() == 0) {
482 std::cout <<
"METRICS FOR THE SERIAL CASE - ONE PHASE MAPPING "
483 <<
"- EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING"
485 metricObject_3->printMetrics(std::cout);
488 for (
int i = 0; i < 3; i++)
489 delete [] partCenters[i];
490 delete [] partCenters;
493 int main(
int narg,
char *arg[]) {
495 Tpetra::ScopeGuard tscope(&narg, &arg);
496 Teuchos::RCP<const Teuchos::Comm<int> > global_tcomm =
497 Tpetra::getDefaultComm();
499 int nx = 16, ny = 16, nz = 16;
500 for (
int i = 1 ; i < narg ; ++i) {
501 if (0 == strcasecmp(arg[i] ,
"NX")) {
502 nx = atoi( arg[++i] );
504 else if (0 == strcasecmp( arg[i] ,
"NY")) {
505 ny = atoi( arg[++i] );
507 else if (0 == strcasecmp( arg[i] ,
"NZ")) {
508 nz = atoi( arg[++i] );
511 std::cerr <<
"Unrecognized command line argument #"
512 << i <<
": " << arg[i] << std::endl ;
519 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm =
520 Teuchos::createSerialComm<int>();
526 part_t my_parts = 0, *my_result_parts;
529 std::cout <<
"me:" << global_tcomm->getRank()
530 <<
" my_parts:" << my_parts
531 <<
" myTasks:" << myTasks << std::endl;
532 if (global_tcomm->getRank() == 0) {
536 FILE *f2 = fopen(
"plot.gnuplot",
"w");
537 for (
int i = 0; i< global_tcomm->getSize(); ++i) {
539 sprintf(str,
"coords%d.txt", i);
541 fprintf(f2,
"splot \"%s\"\n", str);
544 fprintf(f2,
"replot \"%s\"\n", str);
547 fprintf(f2,
"pause-1\n");
551 int myrank = global_tcomm->getRank();
552 sprintf(str,
"coords%d.txt", myrank);
553 FILE *coord_files = fopen(str,
"w");
556 for (
int j = 0; j < my_parts; ++j) {
557 int findex = my_result_parts[j];
558 std::cout <<
"findex " << findex << std::endl;
559 fprintf(coord_files,
"%lf %lf %lf\n",
560 partCenters[0][findex],
561 partCenters[1][findex],
562 partCenters[2][findex]);
568 if (global_tcomm->getRank() == 0) {
569 std::cout <<
"PASS" << std::endl;
572 catch(std::string &s) {
573 std::cerr << s << std::endl;
577 std::cerr << s << std::endl;
RCP< Zoltan2::XpetraMultiVectorAdapter< mytest_tMVector_t > > create_multi_vector_adapter(RCP< const mytest_map_t > map, zscalar_t **partCenters, zgno_t myTasks)
Time an algorithm (or other entity) as a whole.
Zoltan2::XpetraCrsGraphAdapter< mytest_tcrsGraph_t, mytest_tMVector_t > mytest_adapter_t
MappingInputDistributution
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > mytest_tcrsGraph_t
PartitionMapping maps a solution or an input distribution to ranks.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
mytest_adapter_t::part_t mytest_part_t
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
common code used by tests
typename InputTraits< User >::part_t part_t
Tpetra::Map::node_type mytest_znode_t
void test_distributed_input_adapter(int nx, int ny, int nz, Teuchos::RCP< const Teuchos::Comm< int > > global_tcomm)
Tpetra::Map< zlno_t, zgno_t, mytest_znode_t > mytest_map_t
Defines the XpetraMultiVectorAdapter.
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
void test_serial_input_adapter(int nx, int ny, int nz, Teuchos::RCP< const Teuchos::Comm< int > > global_tcomm)
A PartitioningSolution is a solution to a partitioning problem.
Defines the EvaluatePartition class.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
An adapter for Xpetra::MultiVector.
Tpetra::Map::local_ordinal_type zlno_t
void setParts(ArrayRCP< part_t > &partList)
The algorithm uses setParts to set the solution.
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > mytest_tMVector_t
Defines the MappingSolution class.
PartitioningProblem sets up partitioning problems for the user.
RCP< mytest_tcrsGraph_t > create_tpetra_input_matrix(int nx, int ny, int nz, int numProcs, Teuchos::RCP< const Teuchos::Comm< int > > tcomm, RCP< Zoltan2::Environment > env, zscalar_t **&partCenters, zgno_t &myTasks)
Defines the PartitioningProblem class.
MappingProblem enables mapping of a partition (either computed or input) to MPI ranks.
A class that computes and returns quality metrics.
Tpetra::Map::global_ordinal_type zgno_t
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Declarations for TimerManager.
Defines the MappingProblem class.