8 #include <Teuchos_RCP.hpp> 9 #include <Teuchos_Array.hpp> 10 #include <Teuchos_GlobalMPISession.hpp> 11 #include <Teuchos_ParameterList.hpp> 12 #include "Teuchos_XMLParameterListHelpers.hpp" 14 #include "Tpetra_MultiVector.hpp" 15 #include <Tpetra_CrsGraph.hpp> 16 #include <Tpetra_Map.hpp> 36 typedef Tpetra::MultiVector<test_scalar_t, test_lno_t, test_gno_t, znode_t>
mytest_tMVector_t;
39 typedef Tpetra::Map<test_lno_t, test_gno_t, mytest_znode_t>
mytest_map_t;
54 char *input_binary_graph,
55 Teuchos::RCP<
const Teuchos::Comm<int> > tcomm,
57 std::vector <int> &task_communication_xadj_,
58 std::vector <int> &task_communication_adj_,
59 std::vector <double> &task_communication_adjw_){
61 int rank = tcomm->getRank();
69 FILE *f2 = fopen(input_binary_graph,
"rb");
72 fread(&num_vertices,
sizeof(
int),1,f2);
73 fread(&num_edges,
sizeof(
int),1,f2);
75 myTasks = num_vertices;
77 std::cout <<
"numParts:" << num_vertices <<
" ne:" << num_edges << std::endl;
79 task_communication_xadj_.resize(num_vertices+1);
80 task_communication_adj_.resize(num_edges);
81 task_communication_adjw_.resize(num_edges);
83 fread((
void *)&(task_communication_xadj_[0]),
sizeof(
int),num_vertices + 1,f2);
84 fread((
void *)&(task_communication_adj_[0]),
sizeof(
int),num_edges ,f2);
85 fread((
void *)&(task_communication_adjw_[0]),
sizeof(
double),num_edges,f2);
91 tcomm->broadcast(0,
sizeof(
test_lno_t), (
char *) &myTasks);
92 tcomm->broadcast(0,
sizeof(
test_lno_t), (
char *) &myEdges);
95 task_communication_xadj_.resize(myTasks+1);
96 task_communication_adj_.resize(myEdges);
97 task_communication_adjw_.resize(myEdges);
100 tcomm->broadcast(0,
sizeof(
test_lno_t) * (myTasks+1), (
char *) &(task_communication_xadj_[0]));
101 tcomm->broadcast(0,
sizeof(
test_lno_t)* (myEdges), (
char *) &(task_communication_adj_[0]));
102 tcomm->broadcast(0,
sizeof(
test_scalar_t)* (myEdges), (
char *) &(task_communication_adjw_[0]));
106 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm = Teuchos::createSerialComm<int>();
107 RCP<const mytest_map_t> map = rcp (
new mytest_map_t (myTasks, myTasks, 0, serial_comm));
112 std::vector<test_gno_t> tmp(myEdges);
113 for (
test_lno_t lclRow = 0; lclRow < myTasks; ++lclRow) {
114 const test_gno_t gblRow = map->getGlobalElement (lclRow);
115 test_lno_t begin = task_communication_xadj_[gblRow];
116 test_lno_t end = task_communication_xadj_[gblRow + 1];
118 tmp[m - begin] = task_communication_adj_[m];
120 const ArrayView< const test_gno_t > indices(&(tmp[0]), end-begin);
121 TpetraCrsGraph->insertGlobalIndices(gblRow, indices);
123 TpetraCrsGraph->fillComplete ();
126 return TpetraCrsGraph;
131 RCP<const mytest_map_t> map,
int coord_dim,
135 Teuchos::Array<Teuchos::ArrayView<const test_scalar_t> > coordView(coord_dim);
138 for (
int i = 0; i < coord_dim; ++i){
139 Teuchos::ArrayView<const test_scalar_t> a(partCenters[i], myTasks);
144 for (
int i = 0; i < coord_dim; ++i){
145 Teuchos::ArrayView<const test_scalar_t> a;
149 RCP<mytest_tMVector_t> coords(
new mytest_tMVector_t(map, coordView.view(0, coord_dim), coord_dim));
150 RCP<const mytest_tMVector_t> const_coords = rcp_const_cast<
const mytest_tMVector_t>(coords);
157 char *input_binary_graph,
char *input_binary_coordinate,
char *input_machine_file,
158 int machine_optimization_level,
bool divide_prime_first,
int rank_per_node,
bool visualize_mapping,
int reduce_best_mapping){
160 if (input_binary_graph == NULL || input_binary_coordinate == NULL || input_machine_file == NULL){
161 throw "Binary files is mandatory";
164 Teuchos::RCP<const Teuchos::Comm<int> > serial_comm = Teuchos::createSerialComm<int>();
168 Teuchos::ParameterList serial_problemParams;
170 serial_problemParams.set(
"mapping_algorithm",
"geometric");
171 serial_problemParams.set(
"distributed_input_adapter",
false);
172 serial_problemParams.set(
"algorithm",
"multijagged");
173 serial_problemParams.set(
"Machine_Optimization_Level", machine_optimization_level);
174 serial_problemParams.set(
"Input_RCA_Machine_Coords", input_machine_file);
175 serial_problemParams.set(
"divide_prime_first", divide_prime_first);
176 serial_problemParams.set(
"ranks_per_node", rank_per_node);
177 if (reduce_best_mapping)
178 serial_problemParams.set(
"reduce_best_mapping",
true);
180 serial_problemParams.set(
"reduce_best_mapping",
false);
185 serial_problemParams.set(
"num_global_parts", numProcs);
188 env->setTimer(timer);
191 std::vector <double> task_communication_adjw_;
193 std::vector <int> task_communication_xadj_;
194 std::vector <int> task_communication_adj_;
203 task_communication_xadj_, task_communication_adj_,
204 task_communication_adjw_);
205 RCP<const mytest_map_t> serial_map = serial_tpetra_graph->getMap();
206 global_tcomm->barrier();
210 RCP<const mytest_tcrsGraph_t> const_tpetra_graph = rcp_const_cast<
const mytest_tcrsGraph_t>(serial_tpetra_graph);
213 int rank = global_tcomm->getRank();
215 int numParts = 0;
int coordDim = 0;
219 FILE *f2 = fopen(input_binary_coordinate,
"rb");
220 fread((
void *)&numParts,
sizeof(
int),1,f2);
221 fread((
void *)&coordDim,
sizeof(
int),1,f2);
225 for(
int i = 0; i < coordDim; ++i){
227 fread((
void *) partCenters[i],
sizeof(
double),numParts, f2);
232 global_tcomm->broadcast(0,
sizeof(
test_lno_t), (
char *) &numParts);
233 global_tcomm->broadcast(0,
sizeof(
test_lno_t), (
char *) &coordDim);
235 if (numParts!= myTasks){
236 throw "myTasks is different than numParts";
240 for(
int i = 0; i < coordDim; ++i){
245 for(
int i = 0; i < coordDim; ++i){
246 global_tcomm->broadcast(0,
sizeof(
test_scalar_t)* (numParts), (
char *) partCenters[i]);
250 RCP <Zoltan2::XpetraMultiVectorAdapter<mytest_tMVector_t> > serial_adapter =
create_multi_vector_adapter(serial_map, coordDim, partCenters, myTasks);
251 ia->setCoordinateInput(serial_adapter.getRawPtr());
253 ia->setEdgeWeights(&(task_communication_adjw_[0]), 1, 0);
265 global_tcomm->barrier();
287 Teuchos::ArrayView< const test_gno_t> gids = serial_map->getNodeElementList();
289 ArrayRCP<int> initial_part_ids(myTasks);
291 initial_part_ids[i] = gids[i];
293 single_phase_mapping_solution.
setParts(initial_part_ids);
301 ia.getRawPtr(), &serial_problemParams, global_tcomm, &single_phase_mapping_solution, &transformed_machine);
306 serial_map_problem.
solve(
true);
312 timer->printAndResetToZero();
320 RCP<quality_t> metricObject_3 = rcp(
321 new quality_t(ia.getRawPtr(),&serial_problemParams,serial_comm,msoln3, serial_map_problem.getMachine().getRawPtr()));
323 if (global_tcomm->getRank() == 0){
324 std::cout <<
"METRICS FOR THE SERIAL CASE - ONE PHASE MAPPING - EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING" << std::endl;
325 metricObject_3->printMetrics(cout);
327 if (machine_optimization_level > 0){
329 Teuchos::ParameterList serial_problemParams_2;
330 serial_problemParams_2.set(
"Input_RCA_Machine_Coords", input_machine_file);
334 RCP<quality_t> metricObject_4 = rcp(
335 new quality_t(ia.getRawPtr(),&serial_problemParams_2,serial_comm,msoln3, &actual_machine));
337 if (global_tcomm->getRank() == 0){
338 std::cout <<
"METRICS FOR THE SERIAL CASE - ONE PHASE MAPPING - EACH ELEMENT IS ASSUMED TO BE IN UNIQUE PART AT THE BEGINNING" << std::endl;
339 metricObject_4->printMetrics(cout);
343 if (visualize_mapping && global_tcomm->getRank() == 0){
345 Teuchos::ParameterList serial_problemParams_2;
346 serial_problemParams_2.set(
"Input_RCA_Machine_Coords", input_machine_file);
350 Zoltan2::visualize_mapping<zscalar_t, int> (0, actual_machine.
getMachineDim(), actual_machine.
getNumRanks(), coords,
351 int (task_communication_xadj_.size())-1, &(task_communication_xadj_[0]), &(task_communication_adj_[0]), msoln3->getPartListView());
357 int main(
int argc,
char *argv[]){
360 Teuchos::GlobalMPISession session(&argc, &argv);
361 Teuchos::RCP<const Teuchos::Comm<int> > global_tcomm = Teuchos::DefaultComm<int>::getComm();
363 char *input_binary_graph = NULL;
364 char *input_binary_coordinate = NULL;
365 char *input_machine_file = NULL;
366 int machine_optimization_level = 10;
367 bool divide_prime_first =
false;
368 int rank_per_node = 1;
369 int reduce_best_mapping = 1;
371 for (
int i = 1 ; i < argc ; ++i ) {
372 if ( 0 == strcasecmp( argv[i] ,
"BG" ) ) {
374 input_binary_graph = argv[++i];
376 else if ( 0 == strcasecmp( argv[i] ,
"BC" ) ) {
377 input_binary_coordinate = argv[++i];
379 else if ( 0 == strcasecmp( argv[i] ,
"MF" ) ) {
381 input_machine_file = argv[++i];
383 else if ( 0 == strcasecmp( argv[i] ,
"OL" ) ) {
384 machine_optimization_level = atoi( argv[++i] );
386 else if ( 0 == strcasecmp( argv[i] ,
"DP" ) ) {
387 if (atoi( argv[++i] )){
388 divide_prime_first =
true;
391 else if ( 0 == strcasecmp( argv[i] ,
"RPN" ) ) {
392 rank_per_node = atoi( argv[++i] );
394 else if ( 0 == strcasecmp( argv[i] ,
"VM" ) ) {
395 visualize_mapping =
true;
397 else if ( 0 == strcasecmp( argv[i] ,
"RBM" ) ) {
398 reduce_best_mapping = atoi( argv[++i] );
401 std::cerr <<
"Unrecognized command line argument #" << i <<
": " << argv[i] << std::endl ;
410 machine_optimization_level, divide_prime_first, rank_per_node, visualize_mapping, reduce_best_mapping);
414 part_t my_parts = 0, *my_result_parts;
417 std::cout <<
"me:" << global_tcomm->getRank() <<
" my_parts:" << my_parts <<
" myTasks:" << myTasks << std::endl;
418 if (global_tcomm->getRank() == 0) {
422 FILE *f2 = fopen(
"plot.gnuplot",
"w");
423 for (
int i = 0; i< global_tcomm->getSize(); ++i){
425 sprintf(str,
"coords%d.txt", i);
427 fprintf(f2,
"splot \"%s\"\n", str);
430 fprintf(f2,
"replot \"%s\"\n", str);
433 fprintf(f2,
"pause-1\n");
437 int myrank = global_tcomm->getRank();
438 sprintf(str,
"coords%d.txt", myrank);
439 FILE *coord_files = fopen(str,
"w");
442 for (
int j = 0; j < my_parts; ++j){
443 int findex = my_result_parts[j];
444 std::cout <<
"findex " << findex << std::endl;
445 fprintf(coord_files,
"%lf %lf %lf\n", partCenters[0][findex], partCenters[1][findex], partCenters[2][findex]);
451 if (global_tcomm->getRank() == 0){
452 cout <<
"PASS" << endl;
455 catch(std::string &s){
Time an algorithm (or other entity) as a whole.
MappingInputDistributution
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > mytest_tcrsGraph_t
void visualize_mapping(int myRank, const int machine_coord_dim, const int num_ranks, proc_coord_t **machine_coords, const v_lno_t num_tasks, const v_lno_t *task_communication_xadj, const v_lno_t *task_communication_adj, const int *task_to_rank)
PartitionMapping maps a solution or an input distribution to ranks.
mytest_adapter_t::part_t mytest_part_t
int main(int argc, char *argv[])
Tpetra::MultiVector< test_scalar_t, test_lno_t, test_gno_t, znode_t > mytest_tMVector_t
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Provides access for Zoltan2 to Xpetra::CrsGraph data.
void test_serial_input_adapter(Teuchos::RCP< const Teuchos::Comm< int > > global_tcomm, char *input_binary_graph, char *input_binary_coordinate, char *input_machine_file, int machine_optimization_level, bool divide_prime_first, int rank_per_node, bool visualize_mapping, int reduce_best_mapping)
common code used by tests
RCP< mytest_tcrsGraph_t > create_tpetra_input_matrix(char *input_binary_graph, Teuchos::RCP< const Teuchos::Comm< int > > tcomm, test_gno_t &myTasks, std::vector< int > &task_communication_xadj_, std::vector< int > &task_communication_adj_, std::vector< double > &task_communication_adjw_)
Tpetra::Map ::node_type mytest_znode_t
Tpetra::Map< zlno_t, zgno_t, mytest_znode_t > mytest_map_t
Defines the XpetraMultiVectorAdapter.
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
Defines the EvaluatePartition class.
Zoltan2::XpetraCrsGraphAdapter< mytest_tcrsGraph_t, mytest_tMVector_t > mytest_adapter_t
int getMachineDim() const
returns the dimension (number of coords per node) in the machine
InputTraits< User >::part_t part_t
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
bool getAllMachineCoordinatesView(pcoord_t **&allCoords) const
getProcDim function set the coordinates of all ranks allCoords[i][j], i=0,...,getMachineDim(), j=0,...,getNumRanks(), is the i-th dimensional coordinate for rank j. return true if coordinates are available for all ranks
An adapter for Xpetra::MultiVector.
MachineRepresentation Class Base class for representing machine coordinates, networks, etc.
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.
RCP< Zoltan2::XpetraMultiVectorAdapter< mytest_tMVector_t > > create_multi_vector_adapter(RCP< const mytest_map_t > map, int coord_dim, test_scalar_t **partCenters, test_gno_t myTasks)
int getNumRanks() const
return the number of ranks.
Defines the PartitioningProblem class.
Tpetra::Map< test_lno_t, test_gno_t, mytest_znode_t > mytest_map_t
MappingProblem enables mapping of a partition (either computed or input) to MPI ranks.
A class that computes and returns quality metrics.
Declarations for TimerManager.
Defines the MappingProblem class.
Tpetra::CrsGraph< test_lno_t, test_gno_t, znode_t > mytest_tcrsGraph_t