26 #include "Teuchos_RCP.hpp"
27 #include "Teuchos_XMLParameterListHelpers.hpp"
30 #include "create_inline_mesh.h"
32 using Teuchos::ParameterList;
33 using Teuchos::ArrayRCP;
41 int main(
int narg,
char *arg[]) {
43 Tpetra::ScopeGuard tscope(&narg, &arg);
44 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
46 int me = CommT->getRank();
47 int numProcs = CommT->getSize();
51 <<
"====================================================================\n"
53 <<
"| Example: Partition Pamgen Hexahedral Mesh |\n"
55 <<
"| Questions? Contact Karen Devine (kddevin@sandia.gov), |\n"
56 <<
"| Erik Boman (egboman@sandia.gov), |\n"
57 <<
"| Siva Rajamanickam (srajama@sandia.gov). |\n"
59 <<
"| Pamgen's website: http://trilinos.sandia.gov/packages/pamgen |\n"
60 <<
"| Zoltan2's website: http://trilinos.sandia.gov/packages/zoltan2 |\n"
61 <<
"| Trilinos website: http://trilinos.sandia.gov |\n"
63 <<
"====================================================================\n";
69 std::cout <<
"PARALLEL executable \n";
73 std::cout <<
"SERIAL executable \n";
82 std::string xmlMeshInFileName(
"Poisson.xml");
83 std::string action(
"mj");
84 int nParts = CommT->getSize();
87 Teuchos::CommandLineProcessor cmdp (
false,
false);
88 cmdp.setOption(
"xmlfile", &xmlMeshInFileName,
89 "XML file with PamGen specifications");
90 cmdp.setOption(
"action", &action,
91 "Method to use: mj, scotch, zoltan_rcb, zoltan_hg, hg_ghost, "
93 cmdp.setOption(
"nparts", &nParts,
94 "Number of parts to create");
95 cmdp.parse(narg, arg);
98 ParameterList inputMeshList;
100 if(xmlMeshInFileName.length()) {
102 std::cout <<
"\nReading parameter list from the XML file \""
103 <<xmlMeshInFileName<<
"\" ...\n\n";
105 Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
106 Teuchos::inoutArg(inputMeshList));
108 inputMeshList.print(std::cout,2,
true,
true);
113 std::cout <<
"Cannot read input file: " << xmlMeshInFileName <<
"\n";
118 std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
132 if (me == 0) std::cout <<
"Generating mesh ... \n\n";
135 long long maxInt = 9223372036854775807LL;
136 Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
139 if (me == 0) std::cout <<
"Creating mesh adapter ... \n\n";
144 inputAdapter_t *ia =
new inputAdapter_t(*CommT,
"region");
148 if (me == 0) std::cout <<
"Creating parameter list ... \n\n";
150 Teuchos::ParameterList params(
"test params");
151 params.set(
"timer_output_stream" ,
"std::cout");
153 bool do_partitioning =
false;
154 if (action ==
"mj") {
155 do_partitioning =
true;
156 params.set(
"debug_level",
"basic_status");
157 params.set(
"imbalance_tolerance", 1.1);
158 params.set(
"num_global_parts", nParts);
159 params.set(
"algorithm",
"multijagged");
160 params.set(
"rectilinear",
true);
162 else if (action ==
"scotch") {
163 do_partitioning =
true;
164 params.set(
"debug_level",
"verbose_detailed_status");
165 params.set(
"imbalance_tolerance", 1.1);
166 params.set(
"num_global_parts", nParts);
167 params.set(
"partitioning_approach",
"partition");
168 params.set(
"algorithm",
"scotch");
170 else if (action ==
"zoltan_rcb") {
171 do_partitioning =
true;
172 params.set(
"debug_level",
"verbose_detailed_status");
173 params.set(
"imbalance_tolerance", 1.1);
174 params.set(
"num_global_parts", nParts);
175 params.set(
"partitioning_approach",
"partition");
176 params.set(
"algorithm",
"zoltan");
178 else if (action ==
"zoltan_hg") {
179 do_partitioning =
true;
180 params.set(
"debug_level",
"verbose_detailed_status");
181 params.set(
"imbalance_tolerance", 1.1);
182 params.set(
"num_global_parts", nParts);
183 params.set(
"partitioning_approach",
"partition");
184 params.set(
"algorithm",
"zoltan");
185 Teuchos::ParameterList &zparams = params.sublist(
"zoltan_parameters",
false);
186 zparams.set(
"LB_METHOD",
"phg");
187 zparams.set(
"FINAL_OUTPUT",
"1");
189 else if (action==
"hg_ghost") {
190 do_partitioning =
true;
191 params.set(
"debug_level",
"no_status");
192 params.set(
"imbalance_tolerance", 1.1);
193 params.set(
"algorithm",
"zoltan");
194 params.set(
"num_global_parts", nParts);
195 params.set(
"hypergraph_model_type",
"ghosting");
196 params.set(
"ghost_layers",2);
197 Teuchos::ParameterList &zparams = params.sublist(
"zoltan_parameters",
false);
198 zparams.set(
"LB_METHOD",
"HYPERGRAPH");
199 zparams.set(
"LB_APPROACH",
"PARTITION");
200 zparams.set(
"PHG_EDGE_SIZE_THRESHOLD",
"1.0");
203 else if (action ==
"parma") {
204 do_partitioning =
true;
205 params.set(
"debug_level",
"basic_status");
206 params.set(
"imbalance_tolerance", 1.05);
207 params.set(
"algorithm",
"parma");
208 Teuchos::ParameterList &pparams = params.sublist(
"parma_parameters",
false);
209 pparams.set(
"parma_method",
"VtxElm");
211 else if (action==
"zoltan_hg") {
212 do_partitioning =
true;
213 params.set(
"debug_level",
"no_status");
214 params.set(
"imbalance_tolerance", 1.1);
215 params.set(
"algorithm",
"zoltan");
216 params.set(
"num_global_parts", nParts);
217 Teuchos::ParameterList &zparams = params.sublist(
"zoltan_parameters",
false);
218 zparams.set(
"LB_METHOD",
"HYPERGRAPH");
222 else if (action ==
"color") {
223 params.set(
"debug_level",
"verbose_detailed_status");
224 params.set(
"debug_output_file",
"kdd");
225 params.set(
"debug_procs",
"all");
228 if(me == 0) std::cout <<
"Action: " << action << std::endl;
230 if (do_partitioning) {
231 if (me == 0) std::cout <<
"Creating partitioning problem ... \n\n";
236 if (me == 0) std::cout <<
"Calling the partitioner ... \n\n";
242 Teuchos::RCP<quality_t> metricObject =
246 metricObject->printMetrics(std::cout);
250 if (me == 0) std::cout <<
"Creating coloring problem ... \n\n";
255 if (me == 0) std::cout <<
"Calling the coloring algorithm ... \n\n";
263 if (me == 0) std::cout <<
"Deleting the mesh ... \n\n";
265 Delete_Pamgen_Mesh();
268 std::cout <<
"PASS" << std::endl;
ColoringProblem sets up coloring problems for the user.
Defines the ColoringProblem class.
Zoltan2::BasicUserTypes basic_user_t
void solve(bool updateInputData=true)
Direct the problem to create a solution.
A simple class that can be the User template argument for an InputAdapter.
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
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.