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 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
53 /**************************************************************/
54 /* Includes */
55 /**************************************************************/
56 
60 
61 // Teuchos includes
62 #include "Teuchos_RCP.hpp"
63 #include "Teuchos_XMLParameterListHelpers.hpp"
64 
65 // SCOREC includes
66 #ifdef HAVE_ZOLTAN2_PARMA
67 #include <parma.h>
68 #include <apf.h>
69 #include <apfMesh.h>
70 #include <apfMDS.h>
71 #include <apfMesh2.h>
72 #include <PCU.h>
73 #include <gmi_mesh.h>
74 #endif
75 
76 using Teuchos::ParameterList;
77 using Teuchos::RCP;
78 
79 #ifdef HAVE_ZOLTAN2_PARMA
80 
81 void runTest(RCP<const Teuchos::Comm<int> >& CommT, apf::Mesh2* m,std::string action,
82  std::string parma_method,int nParts, double imbalance, std::string output_title );
83 
84 #endif
85 /*****************************************************************************/
86 /******************************** MAIN ***************************************/
87 /*****************************************************************************/
88 
89 int main(int narg, char *arg[]) {
90 
91  Tpetra::ScopeGuard tscope(&narg, &arg);
92  Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
93 
94  int me = CommT->getRank();
95  //int numProcs = CommT->getSize();
96 
97  if (me == 0){
98  std::cout
99  << "====================================================================\n"
100  << "| |\n"
101  << "| Example: Partition APF Mesh |\n"
102  << "| |\n"
103  << "| Questions? Contact Karen Devine (kddevin@sandia.gov), |\n"
104  << "| Erik Boman (egboman@sandia.gov), |\n"
105  << "| Siva Rajamanickam (srajama@sandia.gov). |\n"
106  << "| |\n"
107  << "| Pamgen's website: http://trilinos.sandia.gov/packages/pamgen |\n"
108  << "| Zoltan2's website: http://trilinos.sandia.gov/packages/zoltan2 |\n"
109  << "| Trilinos website: http://trilinos.sandia.gov |\n"
110  << "| |\n"
111  << "====================================================================\n";
112  }
113 
114 
115 #ifdef HAVE_MPI
116  if (me == 0) {
117  std::cout << "PARALLEL executable \n";
118  }
119 #else
120  if (me == 0) {
121  std::cout << "SERIAL executable \n";
122  }
123 #endif
124 
125  /***************************************************************************/
126  /******************************* GET INPUTS ********************************/
127  /***************************************************************************/
128 
129  // default values for command-line arguments
130  std::string meshFileName("4/");
131  std::string modelFileName("torus.dmg");
132  std::string action("zoltan_hg");
133  std::string parma_method("VtxElm");
134  std::string output_loc("");
135  int nParts = CommT->getSize();
136  double imbalance=1.1;
137 
138  // Read run-time options.
139  Teuchos::CommandLineProcessor cmdp (false, false);
140  cmdp.setOption("meshfile", &meshFileName,
141  "Mesh file with APF specifications (.smb file(s))");
142  cmdp.setOption("modelfile", &modelFileName,
143  "Model file with APF specifications (.dmg file)");
144  cmdp.setOption("action", &action,
145  "Method to use: mj, scotch, zoltan_rcb, parma or color");
146  cmdp.setOption("parma_method", &parma_method,
147  "Method to use: Vertex, Edge, Element, VtxElm, VtxEdgeElm, ElmLtVtx, Ghost, or Shape ");
148  cmdp.setOption("nparts", &nParts,
149  "Number of parts to create");
150  cmdp.setOption("imbalance", &imbalance,
151  "Target Imbalance for first partitioner");
152  cmdp.setOption("output", &output_loc,
153  "Location of new partitioned apf mesh. Ex: 4/torus.smb");
154  cmdp.parse(narg, arg);
155 
156 
157  /***************************************************************************/
158  /********************** GET CELL TOPOLOGY **********************************/
159  /***************************************************************************/
160 
161  // Get dimensions
162  //int dim = 3;
163 
164  /***************************************************************************/
165  /***************************** GENERATE MESH *******************************/
166  /***************************************************************************/
167 
168 #ifdef HAVE_ZOLTAN2_PARMA
169 
170  if (me == 0) std::cout << "Generating mesh ... \n\n";
171 
172  //Setup for SCOREC
173  PCU_Comm_Init();
174 
175  // Generate mesh with MDS
176  gmi_register_mesh();
177  apf::Mesh2* m = apf::loadMdsMesh(modelFileName.c_str(),meshFileName.c_str());
178 
179  runTest(CommT,m,action,parma_method,nParts,imbalance,"partition");
180 
181  runTest(CommT,m,"parma",parma_method,nParts,imbalance,"parma");
182 
183 
184 
185 
186  if (output_loc!="") {
187  m->writeNative(output_loc.c_str());
188  }
189 
190  // delete mesh
191  if (me == 0) std::cout << "Deleting the mesh ... \n\n";
192 
193  //Delete APF Mesh;
194  m->destroyNative();
195  apf::destroyMesh(m);
196  //End communications
197  PCU_Comm_Free();
198 
199 #endif
200  if (me == 0)
201  std::cout << "PASS" << std::endl;
202 
203  return 0;
204 
205 }
206 /*****************************************************************************/
207 /********************************* END MAIN **********************************/
208 /*****************************************************************************/
209 
210 #ifdef HAVE_ZOLTAN2_PARMA
211 
212 void runTest(RCP<const Teuchos::Comm<int> >& CommT, apf::Mesh2* m,std::string action,
213  std::string parma_method,int nParts, double imbalance,std::string output_title) {
214  //Get rank
215  int me = CommT->getRank();
216 
217  //Data for APF MeshAdapter
218  std::string primary="region";
219  std::string adjacency="face";
220  if (m->getDimension()==2) {
221  primary="face";
222  adjacency="edge";
223  }
224  bool needSecondAdj=false;
225 
226  // Set parameters for partitioning
227  if (me == 0) std::cout << "Creating parameter list ... \n\n";
228 
229  Teuchos::ParameterList params("test params");
230  params.set("timer_output_stream" , "std::cout");
231 
232  if (action == "mj") {
233  params.set("debug_level", "basic_status");
234  params.set("imbalance_tolerance", imbalance);
235  params.set("num_global_parts", nParts);
236  params.set("algorithm", "multijagged");
237  params.set("rectilinear", "yes");
238  }
239  else if (action == "scotch") {
240  params.set("debug_level", "no_status");
241  params.set("imbalance_tolerance", imbalance);
242  params.set("num_global_parts", nParts);
243  params.set("partitioning_approach", "partition");
244  params.set("algorithm", "scotch");
245  needSecondAdj=true;
246  }
247  else if (action == "zoltan_rcb") {
248  params.set("debug_level", "verbose_detailed_status");
249  params.set("imbalance_tolerance", imbalance);
250  params.set("num_global_parts", nParts);
251  params.set("partitioning_approach", "partition");
252  params.set("algorithm", "zoltan");
253  }
254  else if (action == "parma") {
255  params.set("debug_level", "no_status");
256  params.set("imbalance_tolerance", imbalance);
257  params.set("algorithm", "parma");
258  Teuchos::ParameterList &pparams = params.sublist("parma_parameters",false);
259  pparams.set("parma_method",parma_method);
260  pparams.set("step_size",1.1);
261  if (parma_method=="Ghost") {
262  pparams.set("ghost_layers",3);
263  pparams.set("ghost_bridge",m->getDimension()-1);
264  }
265  adjacency="vertex";
266  }
267  else if (action=="zoltan_hg") {
268  params.set("debug_level", "no_status");
269  params.set("imbalance_tolerance", imbalance);
270  params.set("algorithm", "zoltan");
271  params.set("num_global_parts", nParts);
272  Teuchos::ParameterList &zparams = params.sublist("zoltan_parameters",false);
273  zparams.set("LB_METHOD","HYPERGRAPH");
274  //params.set("compute_metrics","yes");
275  adjacency="vertex";
276  }
277 
278  //Print the stats of original mesh
279  Parma_PrintPtnStats(m,output_title+"_before");
280 
281  // Creating mesh adapter
282  if (me == 0) std::cout << "Creating mesh adapter ... \n\n";
283  typedef Zoltan2::APFMeshAdapter<apf::Mesh2*> inputAdapter_t;
285  typedef Zoltan2::MeshAdapter<apf::Mesh2*> baseMeshAdapter_t;
286 
287  double time_1 = PCU_Time();
288  inputAdapter_t *ia =
289  new inputAdapter_t(*CommT, m,primary,adjacency,needSecondAdj);
290  double time_2 = PCU_Time();
291 
292 
293  // create Partitioning problem
294  if (me == 0) std::cout << "Creating partitioning problem ... \n\n";
295  double time_3=PCU_Time();
296  Zoltan2::PartitioningProblem<inputAdapter_t> problem(ia, &params, CommT);
297 
298  // call the partitioner
299  if (me == 0) std::cout << "Calling the partitioner ... \n\n";
300 
301  problem.solve();
302 
303 
304  //apply the partitioning solution to the mesh
305  if (me==0) std::cout << "Applying Solution to Mesh\n\n";
306  apf::Mesh2** new_mesh = &m;
307  ia->applyPartitioningSolution(m,new_mesh,problem.getSolution());
308  new_mesh=NULL;
309  double time_4=PCU_Time();
310 
311  //create metric object
312  RCP<quality_t> metricObject =
313  rcp(new quality_t(ia, &params, CommT, &problem.getSolution()));
314 
315  if (!me) {
316  metricObject->printMetrics(std::cout);
317  }
318 
319  //Print the stats after partitioning
320  Parma_PrintPtnStats(m,output_title+"_after");
321  ia->destroy();
322 
323  time_4-=time_3;
324  time_2-=time_1;
325  PCU_Max_Doubles(&time_2,1);
326  PCU_Max_Doubles(&time_4,1);
327  if (!me) {
328  std::cout<<"\n"<<output_title<<"Construction time: "<<time_2<<"\n"
329  <<output_title<<"Problem time: " << time_4<<"\n\n";
330  }
331 
332 }
333 
334 #endif
Defines the ColoringProblem class.
MeshAdapter defines the interface for mesh input.
int main(int narg, char *arg[])
bool runTest(Adapter &ia, const RCP< const Teuchos::Comm< int > > &comm, std::string hi)
Definition: Mapping.cpp:341
Defines the APFMeshAdapter class.
PartitioningProblem sets up partitioning problems for the user.
#define nParts
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.