Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
PartitionAndParMATest.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 
17 /**************************************************************/
18 /* Includes */
19 /**************************************************************/
20 
24 
25 // Teuchos includes
26 #include "Teuchos_RCP.hpp"
27 #include "Teuchos_XMLParameterListHelpers.hpp"
28 
29 // SCOREC includes
30 #ifdef HAVE_ZOLTAN2_PARMA
31 #include <parma.h>
32 #include <apf.h>
33 #include <apfMesh.h>
34 #include <apfMDS.h>
35 #include <apfMesh2.h>
36 #include <PCU.h>
37 #include <gmi_mesh.h>
38 #endif
39 
40 using Teuchos::ParameterList;
41 using Teuchos::RCP;
42 
43 #ifdef HAVE_ZOLTAN2_PARMA
44 
45 void runTest(RCP<const Teuchos::Comm<int> >& CommT, apf::Mesh2* m,std::string action,
46  std::string parma_method,int nParts, double imbalance, std::string output_title );
47 
48 #endif
49 /*****************************************************************************/
50 /******************************** MAIN ***************************************/
51 /*****************************************************************************/
52 
53 int main(int narg, char *arg[]) {
54 
55  Tpetra::ScopeGuard tscope(&narg, &arg);
56  Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
57 
58  int me = CommT->getRank();
59  //int numProcs = CommT->getSize();
60 
61  if (me == 0){
62  std::cout
63  << "====================================================================\n"
64  << "| |\n"
65  << "| Example: Partition APF Mesh |\n"
66  << "| |\n"
67  << "| Questions? Contact Karen Devine (kddevin@sandia.gov), |\n"
68  << "| Erik Boman (egboman@sandia.gov), |\n"
69  << "| Siva Rajamanickam (srajama@sandia.gov). |\n"
70  << "| |\n"
71  << "| Pamgen's website: http://trilinos.sandia.gov/packages/pamgen |\n"
72  << "| Zoltan2's website: http://trilinos.sandia.gov/packages/zoltan2 |\n"
73  << "| Trilinos website: http://trilinos.sandia.gov |\n"
74  << "| |\n"
75  << "====================================================================\n";
76  }
77 
78 
79 #ifdef HAVE_MPI
80  if (me == 0) {
81  std::cout << "PARALLEL executable \n";
82  }
83 #else
84  if (me == 0) {
85  std::cout << "SERIAL executable \n";
86  }
87 #endif
88 
89  /***************************************************************************/
90  /******************************* GET INPUTS ********************************/
91  /***************************************************************************/
92 
93  // default values for command-line arguments
94  std::string meshFileName("4/");
95  std::string modelFileName("torus.dmg");
96  std::string action("zoltan_hg");
97  std::string parma_method("VtxElm");
98  std::string output_loc("");
99  int nParts = CommT->getSize();
100  double imbalance=1.1;
101 
102  // Read run-time options.
103  Teuchos::CommandLineProcessor cmdp (false, false);
104  cmdp.setOption("meshfile", &meshFileName,
105  "Mesh file with APF specifications (.smb file(s))");
106  cmdp.setOption("modelfile", &modelFileName,
107  "Model file with APF specifications (.dmg file)");
108  cmdp.setOption("action", &action,
109  "Method to use: mj, scotch, zoltan_rcb, parma or color");
110  cmdp.setOption("parma_method", &parma_method,
111  "Method to use: Vertex, Edge, Element, VtxElm, VtxEdgeElm, ElmLtVtx, Ghost, or Shape ");
112  cmdp.setOption("nparts", &nParts,
113  "Number of parts to create");
114  cmdp.setOption("imbalance", &imbalance,
115  "Target Imbalance for first partitioner");
116  cmdp.setOption("output", &output_loc,
117  "Location of new partitioned apf mesh. Ex: 4/torus.smb");
118  cmdp.parse(narg, arg);
119 
120 
121  /***************************************************************************/
122  /********************** GET CELL TOPOLOGY **********************************/
123  /***************************************************************************/
124 
125  // Get dimensions
126  //int dim = 3;
127 
128  /***************************************************************************/
129  /***************************** GENERATE MESH *******************************/
130  /***************************************************************************/
131 
132 #ifdef HAVE_ZOLTAN2_PARMA
133 
134  if (me == 0) std::cout << "Generating mesh ... \n\n";
135 
136  //Setup for SCOREC
137  PCU_Comm_Init();
138 
139  // Generate mesh with MDS
140  gmi_register_mesh();
141  apf::Mesh2* m = apf::loadMdsMesh(modelFileName.c_str(),meshFileName.c_str());
142 
143  runTest(CommT,m,action,parma_method,nParts,imbalance,"partition");
144 
145  runTest(CommT,m,"parma",parma_method,nParts,imbalance,"parma");
146 
147 
148 
149 
150  if (output_loc!="") {
151  m->writeNative(output_loc.c_str());
152  }
153 
154  // delete mesh
155  if (me == 0) std::cout << "Deleting the mesh ... \n\n";
156 
157  //Delete APF Mesh;
158  m->destroyNative();
159  apf::destroyMesh(m);
160  //End communications
161  PCU_Comm_Free();
162 
163 #endif
164  if (me == 0)
165  std::cout << "PASS" << std::endl;
166 
167  return 0;
168 
169 }
170 /*****************************************************************************/
171 /********************************* END MAIN **********************************/
172 /*****************************************************************************/
173 
174 #ifdef HAVE_ZOLTAN2_PARMA
175 
176 void runTest(RCP<const Teuchos::Comm<int> >& CommT, apf::Mesh2* m,std::string action,
177  std::string parma_method,int nParts, double imbalance,std::string output_title) {
178  //Get rank
179  int me = CommT->getRank();
180 
181  //Data for APF MeshAdapter
182  std::string primary="region";
183  std::string adjacency="face";
184  if (m->getDimension()==2) {
185  primary="face";
186  adjacency="edge";
187  }
188  bool needSecondAdj=false;
189 
190  // Set parameters for partitioning
191  if (me == 0) std::cout << "Creating parameter list ... \n\n";
192 
193  Teuchos::ParameterList params("test params");
194  params.set("timer_output_stream" , "std::cout");
195 
196  if (action == "mj") {
197  params.set("debug_level", "basic_status");
198  params.set("imbalance_tolerance", imbalance);
199  params.set("num_global_parts", nParts);
200  params.set("algorithm", "multijagged");
201  params.set("rectilinear", "yes");
202  }
203  else if (action == "scotch") {
204  params.set("debug_level", "no_status");
205  params.set("imbalance_tolerance", imbalance);
206  params.set("num_global_parts", nParts);
207  params.set("partitioning_approach", "partition");
208  params.set("algorithm", "scotch");
209  needSecondAdj=true;
210  }
211  else if (action == "zoltan_rcb") {
212  params.set("debug_level", "verbose_detailed_status");
213  params.set("imbalance_tolerance", imbalance);
214  params.set("num_global_parts", nParts);
215  params.set("partitioning_approach", "partition");
216  params.set("algorithm", "zoltan");
217  }
218  else if (action == "parma") {
219  params.set("debug_level", "no_status");
220  params.set("imbalance_tolerance", imbalance);
221  params.set("algorithm", "parma");
222  Teuchos::ParameterList &pparams = params.sublist("parma_parameters",false);
223  pparams.set("parma_method",parma_method);
224  pparams.set("step_size",1.1);
225  if (parma_method=="Ghost") {
226  pparams.set("ghost_layers",3);
227  pparams.set("ghost_bridge",m->getDimension()-1);
228  }
229  adjacency="vertex";
230  }
231  else if (action=="zoltan_hg") {
232  params.set("debug_level", "no_status");
233  params.set("imbalance_tolerance", imbalance);
234  params.set("algorithm", "zoltan");
235  params.set("num_global_parts", nParts);
236  Teuchos::ParameterList &zparams = params.sublist("zoltan_parameters",false);
237  zparams.set("LB_METHOD","HYPERGRAPH");
238  //params.set("compute_metrics","yes");
239  adjacency="vertex";
240  }
241 
242  //Print the stats of original mesh
243  Parma_PrintPtnStats(m,output_title+"_before");
244 
245  // Creating mesh adapter
246  if (me == 0) std::cout << "Creating mesh adapter ... \n\n";
247  typedef Zoltan2::APFMeshAdapter<apf::Mesh2*> inputAdapter_t;
249  typedef Zoltan2::MeshAdapter<apf::Mesh2*> baseMeshAdapter_t;
250 
251  double time_1 = PCU_Time();
252  inputAdapter_t *ia =
253  new inputAdapter_t(*CommT, m,primary,adjacency,needSecondAdj);
254  double time_2 = PCU_Time();
255 
256 
257  // create Partitioning problem
258  if (me == 0) std::cout << "Creating partitioning problem ... \n\n";
259  double time_3=PCU_Time();
260  Zoltan2::PartitioningProblem<inputAdapter_t> problem(ia, &params, CommT);
261 
262  // call the partitioner
263  if (me == 0) std::cout << "Calling the partitioner ... \n\n";
264 
265  problem.solve();
266 
267 
268  //apply the partitioning solution to the mesh
269  if (me==0) std::cout << "Applying Solution to Mesh\n\n";
270  apf::Mesh2** new_mesh = &m;
271  ia->applyPartitioningSolution(m,new_mesh,problem.getSolution());
272  new_mesh=NULL;
273  double time_4=PCU_Time();
274 
275  //create metric object
276  RCP<quality_t> metricObject =
277  rcp(new quality_t(ia, &params, CommT, &problem.getSolution()));
278 
279  if (!me) {
280  metricObject->printMetrics(std::cout);
281  }
282 
283  //Print the stats after partitioning
284  Parma_PrintPtnStats(m,output_title+"_after");
285  ia->destroy();
286 
287  time_4-=time_3;
288  time_2-=time_1;
289  PCU_Max_Doubles(&time_2,1);
290  PCU_Max_Doubles(&time_4,1);
291  if (!me) {
292  std::cout<<"\n"<<output_title<<"Construction time: "<<time_2<<"\n"
293  <<output_title<<"Problem time: " << time_4<<"\n\n";
294  }
295 
296 }
297 
298 #endif
Defines the ColoringProblem class.
#define nParts
MeshAdapter defines the interface for mesh input.
bool runTest(Adapter &ia, const RCP< const Teuchos::Comm< int > > &comm, std::string hi)
Definition: Mapping.cpp:306
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
Definition: coloring1.cpp:164
Defines the APFMeshAdapter class.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.