26 #include "Teuchos_RCP.hpp"
27 #include "Teuchos_XMLParameterListHelpers.hpp"
30 #ifdef HAVE_ZOLTAN2_PARMA
40 using Teuchos::ParameterList;
43 #ifdef HAVE_ZOLTAN2_PARMA
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 );
53 int main(
int narg,
char *arg[]) {
55 Tpetra::ScopeGuard tscope(&narg, &arg);
56 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
58 int me = CommT->getRank();
63 <<
"====================================================================\n"
65 <<
"| Example: Partition APF Mesh |\n"
67 <<
"| Questions? Contact Karen Devine (kddevin@sandia.gov), |\n"
68 <<
"| Erik Boman (egboman@sandia.gov), |\n"
69 <<
"| Siva Rajamanickam (srajama@sandia.gov). |\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"
75 <<
"====================================================================\n";
81 std::cout <<
"PARALLEL executable \n";
85 std::cout <<
"SERIAL executable \n";
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;
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);
132 #ifdef HAVE_ZOLTAN2_PARMA
134 if (me == 0) std::cout <<
"Generating mesh ... \n\n";
141 apf::Mesh2* m = apf::loadMdsMesh(modelFileName.c_str(),meshFileName.c_str());
143 runTest(CommT,m,action,parma_method,nParts,imbalance,
"partition");
145 runTest(CommT,m,
"parma",parma_method,nParts,imbalance,
"parma");
150 if (output_loc!=
"") {
151 m->writeNative(output_loc.c_str());
155 if (me == 0) std::cout <<
"Deleting the mesh ... \n\n";
165 std::cout <<
"PASS" << std::endl;
174 #ifdef HAVE_ZOLTAN2_PARMA
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) {
179 int me = CommT->getRank();
182 std::string primary=
"region";
183 std::string adjacency=
"face";
184 if (m->getDimension()==2) {
188 bool needSecondAdj=
false;
191 if (me == 0) std::cout <<
"Creating parameter list ... \n\n";
193 Teuchos::ParameterList params(
"test params");
194 params.set(
"timer_output_stream" ,
"std::cout");
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");
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");
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");
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);
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");
243 Parma_PrintPtnStats(m,output_title+
"_before");
246 if (me == 0) std::cout <<
"Creating mesh adapter ... \n\n";
251 double time_1 = PCU_Time();
253 new inputAdapter_t(*CommT, m,primary,adjacency,needSecondAdj);
254 double time_2 = PCU_Time();
258 if (me == 0) std::cout <<
"Creating partitioning problem ... \n\n";
259 double time_3=PCU_Time();
263 if (me == 0) std::cout <<
"Calling the partitioner ... \n\n";
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());
273 double time_4=PCU_Time();
276 RCP<quality_t> metricObject =
277 rcp(
new quality_t(ia, ¶ms, CommT, &problem.getSolution()));
280 metricObject->printMetrics(std::cout);
284 Parma_PrintPtnStats(m,output_title+
"_after");
289 PCU_Max_Doubles(&time_2,1);
290 PCU_Max_Doubles(&time_4,1);
292 std::cout<<
"\n"<<output_title<<
"Construction time: "<<time_2<<
"\n"
293 <<output_title<<
"Problem time: " << time_4<<
"\n\n";
Defines the ColoringProblem class.
MeshAdapter defines the interface for mesh input.
bool runTest(Adapter &ia, const RCP< const Teuchos::Comm< int > > &comm, std::string hi)
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
Defines the APFMeshAdapter class.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.