62 #include "Teuchos_RCP.hpp"
63 #include "Teuchos_XMLParameterListHelpers.hpp"
66 #include "create_inline_mesh.h"
68 using Teuchos::ParameterList;
69 using Teuchos::ArrayRCP;
77 int main(
int narg,
char *arg[]) {
79 Tpetra::ScopeGuard tscope(&narg, &arg);
80 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
82 int me = CommT->getRank();
83 int numProcs = CommT->getSize();
87 <<
"====================================================================\n"
89 <<
"| Example: Partition Pamgen Hexahedral Mesh |\n"
91 <<
"| Questions? Contact Karen Devine (kddevin@sandia.gov), |\n"
92 <<
"| Erik Boman (egboman@sandia.gov), |\n"
93 <<
"| Siva Rajamanickam (srajama@sandia.gov). |\n"
95 <<
"| Pamgen's website: http://trilinos.sandia.gov/packages/pamgen |\n"
96 <<
"| Zoltan2's website: http://trilinos.sandia.gov/packages/zoltan2 |\n"
97 <<
"| Trilinos website: http://trilinos.sandia.gov |\n"
99 <<
"====================================================================\n";
105 std::cout <<
"PARALLEL executable \n";
109 std::cout <<
"SERIAL executable \n";
118 std::string xmlMeshInFileName(
"Poisson.xml");
119 std::string action(
"mj");
120 int nParts = CommT->getSize();
123 Teuchos::CommandLineProcessor cmdp (
false,
false);
124 cmdp.setOption(
"xmlfile", &xmlMeshInFileName,
125 "XML file with PamGen specifications");
126 cmdp.setOption(
"action", &action,
127 "Method to use: mj, scotch, zoltan_rcb, zoltan_hg, hg_ghost, "
129 cmdp.setOption(
"nparts", &nParts,
130 "Number of parts to create");
131 cmdp.parse(narg, arg);
134 ParameterList inputMeshList;
136 if(xmlMeshInFileName.length()) {
138 std::cout <<
"\nReading parameter list from the XML file \""
139 <<xmlMeshInFileName<<
"\" ...\n\n";
141 Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
142 Teuchos::inoutArg(inputMeshList));
144 inputMeshList.print(std::cout,2,
true,
true);
149 std::cout <<
"Cannot read input file: " << xmlMeshInFileName <<
"\n";
154 std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
168 if (me == 0) std::cout <<
"Generating mesh ... \n\n";
171 long long maxInt = 9223372036854775807LL;
172 Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
175 if (me == 0) std::cout <<
"Creating mesh adapter ... \n\n";
180 inputAdapter_t *ia =
new inputAdapter_t(*CommT,
"region");
184 if (me == 0) std::cout <<
"Creating parameter list ... \n\n";
186 Teuchos::ParameterList params(
"test params");
187 params.set(
"timer_output_stream" ,
"std::cout");
189 bool do_partitioning =
false;
190 if (action ==
"mj") {
191 do_partitioning =
true;
192 params.set(
"debug_level",
"basic_status");
193 params.set(
"imbalance_tolerance", 1.1);
194 params.set(
"num_global_parts", nParts);
195 params.set(
"algorithm",
"multijagged");
196 params.set(
"rectilinear",
true);
198 else if (action ==
"scotch") {
199 do_partitioning =
true;
200 params.set(
"debug_level",
"verbose_detailed_status");
201 params.set(
"imbalance_tolerance", 1.1);
202 params.set(
"num_global_parts", nParts);
203 params.set(
"partitioning_approach",
"partition");
204 params.set(
"algorithm",
"scotch");
206 else if (action ==
"zoltan_rcb") {
207 do_partitioning =
true;
208 params.set(
"debug_level",
"verbose_detailed_status");
209 params.set(
"imbalance_tolerance", 1.1);
210 params.set(
"num_global_parts", nParts);
211 params.set(
"partitioning_approach",
"partition");
212 params.set(
"algorithm",
"zoltan");
214 else if (action ==
"zoltan_hg") {
215 do_partitioning =
true;
216 params.set(
"debug_level",
"verbose_detailed_status");
217 params.set(
"imbalance_tolerance", 1.1);
218 params.set(
"num_global_parts", nParts);
219 params.set(
"partitioning_approach",
"partition");
220 params.set(
"algorithm",
"zoltan");
221 Teuchos::ParameterList &zparams = params.sublist(
"zoltan_parameters",
false);
222 zparams.set(
"LB_METHOD",
"phg");
223 zparams.set(
"FINAL_OUTPUT",
"1");
225 else if (action==
"hg_ghost") {
226 do_partitioning =
true;
227 params.set(
"debug_level",
"no_status");
228 params.set(
"imbalance_tolerance", 1.1);
229 params.set(
"algorithm",
"zoltan");
230 params.set(
"num_global_parts", nParts);
231 params.set(
"hypergraph_model_type",
"ghosting");
232 params.set(
"ghost_layers",2);
233 Teuchos::ParameterList &zparams = params.sublist(
"zoltan_parameters",
false);
234 zparams.set(
"LB_METHOD",
"HYPERGRAPH");
235 zparams.set(
"LB_APPROACH",
"PARTITION");
236 zparams.set(
"PHG_EDGE_SIZE_THRESHOLD",
"1.0");
239 else if (action ==
"parma") {
240 do_partitioning =
true;
241 params.set(
"debug_level",
"basic_status");
242 params.set(
"imbalance_tolerance", 1.05);
243 params.set(
"algorithm",
"parma");
244 Teuchos::ParameterList &pparams = params.sublist(
"parma_parameters",
false);
245 pparams.set(
"parma_method",
"VtxElm");
247 else if (action==
"zoltan_hg") {
248 do_partitioning =
true;
249 params.set(
"debug_level",
"no_status");
250 params.set(
"imbalance_tolerance", 1.1);
251 params.set(
"algorithm",
"zoltan");
252 params.set(
"num_global_parts", nParts);
253 Teuchos::ParameterList &zparams = params.sublist(
"zoltan_parameters",
false);
254 zparams.set(
"LB_METHOD",
"HYPERGRAPH");
258 else if (action ==
"color") {
259 params.set(
"debug_level",
"verbose_detailed_status");
260 params.set(
"debug_output_file",
"kdd");
261 params.set(
"debug_procs",
"all");
264 if(me == 0) std::cout <<
"Action: " << action << std::endl;
266 if (do_partitioning) {
267 if (me == 0) std::cout <<
"Creating partitioning problem ... \n\n";
272 if (me == 0) std::cout <<
"Calling the partitioner ... \n\n";
278 Teuchos::RCP<quality_t> metricObject =
279 rcp(
new quality_t(ia, ¶ms, CommT, &problem.
getSolution()));
282 metricObject->printMetrics(std::cout);
286 if (me == 0) std::cout <<
"Creating coloring problem ... \n\n";
291 if (me == 0) std::cout <<
"Calling the coloring algorithm ... \n\n";
299 if (me == 0) std::cout <<
"Deleting the mesh ... \n\n";
301 Delete_Pamgen_Mesh();
304 std::cout <<
"PASS" << std::endl;
ColoringProblem sets up coloring problems for the user.
Defines the ColoringProblem class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
int main(int narg, char *arg[])
A simple class that can be the User template argument for an InputAdapter.
Defines the PamgenMeshAdapter class.
void printTimers() const
Return the communicator passed to the problem.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.
This class represents a mesh.
void solve(bool updateInputData=true)
Direct the problem to create a solution.