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>
44 int nx,
int ny,
int nz,
45 int numProcs, Teuchos::RCP<
const Teuchos::Comm<int> > tcomm,
46 RCP<Zoltan2::Environment> env,
zscalar_t ** &partCenters,
zgno_t & myTasks)
48 int rank = tcomm->getRank();
49 using namespace Teuchos;
52 zgno_t numGlobalTasks = nx*ny*nz;
55 myTasks = numGlobalTasks / numProcs;
56 zgno_t taskLeftOver = numGlobalTasks % numProcs;
57 if (
zgno_t(rank) < taskLeftOver ) ++myTasks;
59 zgno_t myTaskBegin = (numGlobalTasks / numProcs) * rank;
60 myTaskBegin += (taskLeftOver <
zgno_t(rank) ? taskLeftOver : rank);
61 zgno_t myTaskEnd = myTaskBegin + myTasks;
65 for(
int i = 0; i < coordDim; ++i){
69 zgno_t *task_communication_xadj_ =
new zgno_t [myTasks+1];
70 zgno_t *task_communication_adj_ =
new zgno_t [myTasks * 6];
74 task_communication_xadj_[0] = 0;
75 for (
zgno_t i = myTaskBegin; i < myTaskEnd; ++i) {
77 int y = (i / (nx)) % ny;
78 int z = (i / (nx)) / ny;
79 partCenters[0][i - myTaskBegin] = x;
80 partCenters[1][i - myTaskBegin] = y;
81 partCenters[2][i - myTaskBegin] = z;
84 task_communication_adj_[prevNCount++] = i - 1;
87 task_communication_adj_[prevNCount++] = i + 1;
90 task_communication_adj_[prevNCount++] = i - nx;
93 task_communication_adj_[prevNCount++] = i + nx;
96 task_communication_adj_[prevNCount++] = i - nx * ny;
99 task_communication_adj_[prevNCount++] = i + nx * ny;
101 task_communication_xadj_[i+1 - myTaskBegin] = prevNCount;
106 using namespace Teuchos;
107 RCP<const mytest_map_t> map = rcp (
new mytest_map_t (numGlobalTasks, myTasks, 0, tcomm));
113 for (
zgno_t lclRow = 0; lclRow < myTasks; ++lclRow) {
114 const zgno_t gblRow = map->getGlobalElement (lclRow);
115 zgno_t begin = task_communication_xadj_[lclRow];
116 zgno_t end = task_communication_xadj_[lclRow + 1];
117 const ArrayView< const zgno_t > indices(task_communication_adj_+begin, end-begin);
118 TpetraCrsGraph->insertGlobalIndices(gblRow, indices);
120 TpetraCrsGraph->fillComplete ();
122 delete [] task_communication_xadj_;
123 delete [] task_communication_adj_;
126 return TpetraCrsGraph;
131 RCP<const mytest_map_t> map,
134 const int coord_dim = 3;
135 Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
138 Teuchos::ArrayView<const zscalar_t> a(partCenters[0], myTasks);
140 Teuchos::ArrayView<const zscalar_t> b(partCenters[1], myTasks);
142 Teuchos::ArrayView<const zscalar_t> c(partCenters[2], myTasks);
146 Teuchos::ArrayView<const zscalar_t> a;
151 RCP<mytest_tMVector_t> coords(
new mytest_tMVector_t(map, coordView.view(0, coord_dim), coord_dim));
152 RCP<const mytest_tMVector_t> const_coords = rcp_const_cast<
const mytest_tMVector_t>(coords);
159 int nx,
int ny,
int nz, Teuchos::RCP<
const Teuchos::Comm<int> > global_tcomm)
161 Teuchos::RCP<const Teuchos::Comm<int> > tcomm = global_tcomm;
163 Teuchos::ParameterList distributed_problemParams;
165 distributed_problemParams.set(
"Machine_Optimization_Level", 10);
166 distributed_problemParams.set(
"mapping_algorithm",
"geometric");
167 distributed_problemParams.set(
"distributed_input_adapter",
true);
168 distributed_problemParams.set(
"mj_enable_rcb",
true);
169 distributed_problemParams.set(
"algorithm",
"multijagged");
170 distributed_problemParams.set(
"num_global_parts", numProcs);
172 RCP<Zoltan2::Environment> env(
new Zoltan2::Environment(distributed_problemParams, global_tcomm));
174 env->setTimer(timer);
180 RCP<mytest_tcrsGraph_t> distributed_tpetra_graph =
182 tcomm, env, partCenters,
184 RCP<const mytest_map_t> distributed_map = distributed_tpetra_graph->getMap();
185 global_tcomm->barrier();
189 RCP<const mytest_tcrsGraph_t> const_tpetra_graph =
194 RCP<Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> > distributed_adapter =
create_multi_vector_adapter(distributed_map, partCenters, myTasks);
195 ia->setCoordinateInput(distributed_adapter.getRawPtr());
197 global_tcomm->barrier();
212 partitioning_problem.
solve();
220 Teuchos::ArrayView< const zgno_t> gids = distributed_map->getNodeElementList();
222 ArrayRCP<int> initial_part_ids(myTasks);
223 for (
zgno_t i = 0; i < myTasks; ++i){
224 initial_part_ids[i] = gids[i];
226 single_phase_mapping_solution.
setParts(initial_part_ids);
243 distributed_map_problem_1.
solve(
true);
244 distributed_map_problem_2.solve(
true);
245 distributed_map_problem_3.solve(
true);
253 timer->printAndResetToZero();
259 RCP<quality_t> metricObject_1 = rcp(
new quality_t(ia.getRawPtr(),&distributed_problemParams,global_tcomm,msoln1, distributed_map_problem_1.getMachine().getRawPtr()));
261 RCP<quality_t> metricObject_2 = rcp(
new quality_t(ia.getRawPtr(),&distributed_problemParams,global_tcomm,msoln2, distributed_map_problem_2.getMachine().getRawPtr()));
263 RCP<quality_t> metricObject_3 = rcp(
new quality_t(ia.getRawPtr(),&distributed_problemParams,global_tcomm,msoln3, distributed_map_problem_3.getMachine().getRawPtr()));
266 if (global_tcomm->getRank() == 0){
267 std::cout <<
"METRICS FOR THE FIRST CASE - TWO PHASE MAPPING" << std::endl;
268 metricObject_1->printMetrics(std::cout);
269 std::cout <<
"METRICS FOR THE SECOND CASE - TWO PHASE MAPPING - INITIAL ASSIGNMENT ARE ASSUMED TO BE A PART" << std::endl;
270 metricObject_2->printMetrics(std::cout);
271 std::cout <<
"METRICS FOR THE THIRD CASE - ONE PHASE MAPPING - EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING" << std::endl;
272 metricObject_3->printMetrics(std::cout);
275 for (
int i = 0; i < 3; i++)
delete [] partCenters[i];
276 delete [] partCenters;
282 int nx,
int ny,
int nz, Teuchos::RCP<
const Teuchos::Comm<int> > global_tcomm)
285 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm =
286 Teuchos::createSerialComm<int>();
290 Teuchos::ParameterList serial_problemParams;
292 serial_problemParams.set(
"Machine_Optimization_Level", 10);
293 serial_problemParams.set(
"mapping_algorithm",
"geometric");
294 serial_problemParams.set(
"distributed_input_adapter",
false);
295 serial_problemParams.set(
"algorithm",
"multijagged");
296 serial_problemParams.set(
"num_global_parts", numProcs);
300 env->setTimer(timer);
306 RCP<mytest_tcrsGraph_t> serial_tpetra_graph =
create_tpetra_input_matrix(nx, ny, nz, numProcs, serial_comm, env, partCenters, myTasks);
307 RCP<const mytest_map_t> serial_map = serial_tpetra_graph->getMap();
308 global_tcomm->barrier();
312 RCP<const mytest_tcrsGraph_t> const_tpetra_graph = rcp_const_cast<
const mytest_tcrsGraph_t>(serial_tpetra_graph);
316 RCP <Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> > serial_adapter =
create_multi_vector_adapter(serial_map, partCenters, myTasks);
317 ia->setCoordinateInput(serial_adapter.getRawPtr());
319 global_tcomm->barrier();
341 Teuchos::ArrayView< const zgno_t> gids = serial_map->getNodeElementList();
343 ArrayRCP<int> initial_part_ids(myTasks);
344 for (
zgno_t i = 0; i < myTasks; ++i){
345 initial_part_ids[i] = gids[i];
347 single_phase_mapping_solution.
setParts(initial_part_ids);
359 serial_map_problem.
solve(
true);
365 timer->printAndResetToZero();
373 RCP<quality_t> metricObject_3 = rcp(
new quality_t(ia.getRawPtr(),&serial_problemParams,serial_comm,msoln3, serial_map_problem.getMachine().getRawPtr()));
375 if (global_tcomm->getRank() == 0){
376 std::cout <<
"METRICS FOR THE SERIAL CASE - ONE PHASE MAPPING - EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING" << std::endl;
377 metricObject_3->printMetrics(std::cout);
380 for (
int i = 0; i < 3; i++)
delete [] partCenters[i];
381 delete [] partCenters;
384 int main(
int narg,
char *arg[]){
386 Tpetra::ScopeGuard tscope(&narg, &arg);
387 Teuchos::RCP<const Teuchos::Comm<int> > global_tcomm = Tpetra::getDefaultComm();
389 int nx = 16, ny = 16, nz = 16;
390 for (
int i = 1 ; i < narg ; ++i ) {
391 if ( 0 == strcasecmp( arg[i] ,
"NX" ) ) {
392 nx = atoi( arg[++i] );
394 else if ( 0 == strcasecmp( arg[i] ,
"NY" ) ) {
395 ny = atoi( arg[++i] );
397 else if ( 0 == strcasecmp( arg[i] ,
"NZ" ) ) {
398 nz = atoi( arg[++i] );
401 std::cerr <<
"Unrecognized command line argument #" << i <<
": " << arg[i] << std::endl ;
410 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm = Teuchos::createSerialComm<int>();
416 part_t my_parts = 0, *my_result_parts;
419 std::cout <<
"me:" << global_tcomm->getRank() <<
" my_parts:" << my_parts <<
" myTasks:" << myTasks << std::endl;
420 if (global_tcomm->getRank() == 0) {
424 FILE *f2 = fopen(
"plot.gnuplot",
"w");
425 for (
int i = 0; i< global_tcomm->getSize(); ++i){
427 sprintf(str,
"coords%d.txt", i);
429 fprintf(f2,
"splot \"%s\"\n", str);
432 fprintf(f2,
"replot \"%s\"\n", str);
435 fprintf(f2,
"pause-1\n");
439 int myrank = global_tcomm->getRank();
440 sprintf(str,
"coords%d.txt", myrank);
441 FILE *coord_files = fopen(str,
"w");
444 for (
int j = 0; j < my_parts; ++j){
445 int findex = my_result_parts[j];
446 std::cout <<
"findex " << findex << std::endl;
447 fprintf(coord_files,
"%lf %lf %lf\n", partCenters[0][findex], partCenters[1][findex], partCenters[2][findex]);
453 if (global_tcomm->getRank() == 0){
454 std::cout <<
"PASS" << std::endl;
457 catch(std::string &s){
458 std::cerr << s << std::endl;
462 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.