8 #include <Teuchos_RCP.hpp>
9 #include <Teuchos_Array.hpp>
10 #include <Teuchos_ParameterList.hpp>
11 #include "Teuchos_XMLParameterListHelpers.hpp"
13 #include "Tpetra_MultiVector.hpp"
14 #include <Tpetra_CrsGraph.hpp>
15 #include <Tpetra_Map.hpp>
26 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>
46 int nx,
int ny,
int nz,
47 int numProcs, Teuchos::RCP<
const Teuchos::Comm<int> > tcomm,
48 RCP<Zoltan2::Environment> env,
zscalar_t ** &partCenters,
zgno_t & myTasks) {
50 int rank = tcomm->getRank();
51 using namespace Teuchos;
54 zgno_t numGlobalTasks = nx*ny*nz;
57 myTasks = numGlobalTasks / numProcs;
58 zgno_t taskLeftOver = numGlobalTasks % numProcs;
59 if (
zgno_t(rank) < taskLeftOver ) ++myTasks;
61 zgno_t myTaskBegin = (numGlobalTasks / numProcs) * rank;
62 myTaskBegin += (taskLeftOver <
zgno_t(rank) ? taskLeftOver : rank);
63 zgno_t myTaskEnd = myTaskBegin + myTasks;
67 for(
int i = 0; i < coordDim; ++i) {
71 zgno_t *task_communication_xadj_ =
new zgno_t [myTasks + 1];
72 zgno_t *task_communication_adj_ =
new zgno_t [myTasks * 6];
76 task_communication_xadj_[0] = 0;
77 for (
zgno_t i = myTaskBegin; i < myTaskEnd; ++i) {
79 int y = (i / (nx)) % ny;
80 int z = (i / (nx)) / ny;
81 partCenters[0][i - myTaskBegin] = x;
82 partCenters[1][i - myTaskBegin] = y;
83 partCenters[2][i - myTaskBegin] = z;
86 task_communication_adj_[prevNCount++] = i - 1;
89 task_communication_adj_[prevNCount++] = i + 1;
92 task_communication_adj_[prevNCount++] = i - nx;
95 task_communication_adj_[prevNCount++] = i + nx;
98 task_communication_adj_[prevNCount++] = i - nx * ny;
101 task_communication_adj_[prevNCount++] = i + nx * ny;
104 task_communication_xadj_[i + 1 - myTaskBegin] = prevNCount;
109 using namespace Teuchos;
110 RCP<const mytest_map_t> map = rcp (
new mytest_map_t (numGlobalTasks,
113 Teuchos::Array<size_t> adjPerTask(myTasks);
114 for (
zgno_t lclRow = 0; lclRow < myTasks; ++lclRow)
115 adjPerTask[lclRow] = task_communication_xadj_[lclRow+1]
116 - task_communication_xadj_[lclRow];
122 for (
zgno_t lclRow = 0; lclRow < myTasks; ++lclRow) {
123 const zgno_t gblRow = map->getGlobalElement (lclRow);
124 zgno_t begin = task_communication_xadj_[lclRow];
125 zgno_t end = task_communication_xadj_[lclRow + 1];
126 const ArrayView< const zgno_t > indices(task_communication_adj_ + begin,
128 TpetraCrsGraph->insertGlobalIndices(gblRow, indices);
130 TpetraCrsGraph->fillComplete ();
132 delete [] task_communication_xadj_;
133 delete [] task_communication_adj_;
136 return TpetraCrsGraph;
140 RCP<Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> >
145 const int coord_dim = 3;
146 Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
149 Teuchos::ArrayView<const zscalar_t> a(partCenters[0], myTasks);
151 Teuchos::ArrayView<const zscalar_t> b(partCenters[1], myTasks);
153 Teuchos::ArrayView<const zscalar_t> c(partCenters[2], myTasks);
157 Teuchos::ArrayView<const zscalar_t> a;
164 RCP<mytest_tMVector_t> coords(
166 RCP<const mytest_tMVector_t> const_coords =
168 RCP <Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> > adapter(
175 int nx,
int ny,
int nz,
176 Teuchos::RCP<
const Teuchos::Comm<int> > global_tcomm)
178 Teuchos::RCP<const Teuchos::Comm<int> > tcomm = global_tcomm;
182 Teuchos::ParameterList distributed_problemParams;
185 distributed_problemParams.set(
"Machine_Optimization_Level", 10);
186 distributed_problemParams.set(
"mapping_algorithm",
"geometric");
187 distributed_problemParams.set(
"distributed_input_adapter",
true);
188 distributed_problemParams.set(
"mj_enable_rcb",
true);
189 distributed_problemParams.set(
"algorithm",
"multijagged");
190 distributed_problemParams.set(
"num_global_parts", numProcs);
192 RCP<Zoltan2::Environment> env(
194 RCP<Zoltan2::TimerManager> timer(
197 env->setTimer(timer);
203 RCP<mytest_tcrsGraph_t> distributed_tpetra_graph =
205 tcomm, env, partCenters,
207 RCP<const mytest_map_t> distributed_map =
208 distributed_tpetra_graph->getMap();
209 global_tcomm->barrier();
213 RCP<const mytest_tcrsGraph_t> const_tpetra_graph =
218 RCP<Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> >
220 partCenters, myTasks);
221 ia->setCoordinateInput(distributed_adapter.getRawPtr());
223 global_tcomm->barrier();
244 ia.getRawPtr(), &distributed_problemParams, global_tcomm);
246 partitioning_problem.
solve();
249 partitioning_problem.getSolution();
257 single_phase_mapping_solution(env, global_tcomm, 0);
258 Teuchos::ArrayView< const zgno_t> gids =
259 distributed_map->getNodeElementList();
261 ArrayRCP<int> initial_part_ids(myTasks);
262 for (
zgno_t i = 0; i < myTasks; ++i) {
263 initial_part_ids[i] = gids[i];
265 single_phase_mapping_solution.
setParts(initial_part_ids);
271 ia.getRawPtr(), &distributed_problemParams,
272 global_tcomm, &partition_solution);
279 ia.getRawPtr(), &distributed_problemParams, global_tcomm);
284 ia.getRawPtr(), &distributed_problemParams,
285 global_tcomm, &single_phase_mapping_solution);
291 distributed_map_problem_1.
solve(
true);
293 distributed_map_problem_2.solve(
true);
295 distributed_map_problem_3.solve(
true);
302 distributed_map_problem_1.getSolution();
305 distributed_map_problem_2.getSolution();
308 distributed_map_problem_3.getSolution();
310 timer->printAndResetToZero();
316 RCP<quality_t> metricObject_1 =
317 rcp(
new quality_t(ia.getRawPtr(), &distributed_problemParams,
318 global_tcomm, msoln1,
319 distributed_map_problem_1.getMachine().getRawPtr()));
322 RCP<quality_t> metricObject_2 =
323 rcp(
new quality_t(ia.getRawPtr(), &distributed_problemParams,
324 global_tcomm, msoln2,
325 distributed_map_problem_2.getMachine().getRawPtr()));
328 RCP<quality_t> metricObject_3 =
329 rcp(
new quality_t(ia.getRawPtr(), &distributed_problemParams,
330 global_tcomm, msoln3,
331 distributed_map_problem_3.getMachine().getRawPtr()));
334 if (global_tcomm->getRank() == 0) {
335 std::cout <<
"METRICS FOR THE FIRST CASE - TWO PHASE MAPPING"
337 metricObject_1->printMetrics(std::cout);
338 std::cout <<
"METRICS FOR THE SECOND CASE - TWO PHASE MAPPING"
339 <<
" - INITIAL ASSIGNMENT ARE ASSUMED TO BE A PART" << std::endl;
340 metricObject_2->printMetrics(std::cout);
341 std::cout <<
"METRICS FOR THE THIRD CASE - ONE PHASE MAPPING"
342 <<
" - EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING"
344 metricObject_3->printMetrics(std::cout);
347 for (
int i = 0; i < 3; i++)
348 delete [] partCenters[i];
349 delete [] partCenters;
355 int nx,
int ny,
int nz,
356 Teuchos::RCP<
const Teuchos::Comm<int> > global_tcomm)
359 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm =
360 Teuchos::createSerialComm<int>();
365 Teuchos::ParameterList serial_problemParams;
367 serial_problemParams.set(
"Machine_Optimization_Level", 10);
368 serial_problemParams.set(
"mapping_algorithm",
"geometric");
369 serial_problemParams.set(
"distributed_input_adapter",
false);
370 serial_problemParams.set(
"algorithm",
"multijagged");
371 serial_problemParams.set(
"num_global_parts", numProcs);
373 RCP<Zoltan2::Environment> env(
375 RCP<Zoltan2::TimerManager> timer(
378 env->setTimer(timer);
384 RCP<mytest_tcrsGraph_t> serial_tpetra_graph =
386 env, partCenters, myTasks);
387 RCP<const mytest_map_t> serial_map = serial_tpetra_graph->getMap();
388 global_tcomm->barrier();
392 RCP<const mytest_tcrsGraph_t> const_tpetra_graph =
397 RCP <Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t>> serial_adapter =
399 ia->setCoordinateInput(serial_adapter.getRawPtr());
401 global_tcomm->barrier();
429 single_phase_mapping_solution(env, global_tcomm, 0);
430 Teuchos::ArrayView< const zgno_t> gids = serial_map->getNodeElementList();
432 ArrayRCP<int> initial_part_ids(myTasks);
433 for (
zgno_t i = 0; i < myTasks; ++i) {
434 initial_part_ids[i] = gids[i];
436 single_phase_mapping_solution.
setParts(initial_part_ids);
446 ia.getRawPtr(), &serial_problemParams,
447 global_tcomm, &single_phase_mapping_solution);
453 serial_map_problem.
solve(
true);
458 serial_map_problem.getSolution();
460 timer->printAndResetToZero();
468 RCP<quality_t> metricObject_3 =
469 rcp(
new quality_t(ia.getRawPtr(),
470 &serial_problemParams, serial_comm,msoln3,
471 serial_map_problem.getMachine().getRawPtr()));
473 if (global_tcomm->getRank() == 0) {
474 std::cout <<
"METRICS FOR THE SERIAL CASE - ONE PHASE MAPPING "
475 <<
"- EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING"
477 metricObject_3->printMetrics(std::cout);
480 for (
int i = 0; i < 3; i++)
481 delete [] partCenters[i];
482 delete [] partCenters;
485 int main(
int narg,
char *arg[]) {
487 Tpetra::ScopeGuard tscope(&narg, &arg);
488 Teuchos::RCP<const Teuchos::Comm<int> > global_tcomm =
489 Tpetra::getDefaultComm();
491 int nx = 16, ny = 16, nz = 16;
492 for (
int i = 1 ; i < narg ; ++i) {
493 if (0 == strcasecmp(arg[i] ,
"NX")) {
494 nx = atoi( arg[++i] );
496 else if (0 == strcasecmp( arg[i] ,
"NY")) {
497 ny = atoi( arg[++i] );
499 else if (0 == strcasecmp( arg[i] ,
"NZ")) {
500 nz = atoi( arg[++i] );
503 std::cerr <<
"Unrecognized command line argument #"
504 << i <<
": " << arg[i] << std::endl ;
511 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm =
512 Teuchos::createSerialComm<int>();
518 part_t my_parts = 0, *my_result_parts;
521 std::cout <<
"me:" << global_tcomm->getRank()
522 <<
" my_parts:" << my_parts
523 <<
" myTasks:" << myTasks << std::endl;
524 if (global_tcomm->getRank() == 0) {
528 FILE *f2 = fopen(
"plot.gnuplot",
"w");
529 for (
int i = 0; i< global_tcomm->getSize(); ++i) {
531 sprintf(str,
"coords%d.txt", i);
533 fprintf(f2,
"splot \"%s\"\n", str);
536 fprintf(f2,
"replot \"%s\"\n", str);
539 fprintf(f2,
"pause-1\n");
543 int myrank = global_tcomm->getRank();
544 sprintf(str,
"coords%d.txt", myrank);
545 FILE *coord_files = fopen(str,
"w");
548 for (
int j = 0; j < my_parts; ++j) {
549 int findex = my_result_parts[j];
550 std::cout <<
"findex " << findex << std::endl;
551 fprintf(coord_files,
"%lf %lf %lf\n",
552 partCenters[0][findex],
553 partCenters[1][findex],
554 partCenters[2][findex]);
560 if (global_tcomm->getRank() == 0) {
561 std::cout <<
"PASS" << std::endl;
564 catch(std::string &s) {
565 std::cerr << s << std::endl;
569 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.
int main(int narg, char *arg[])
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.
common code used by tests
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.
InputTraits< User >::part_t part_t
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
An adapter for Xpetra::MultiVector.
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.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Declarations for TimerManager.
Defines the MappingProblem class.