Zoltan2
Loading...
Searching...
No Matches
pamgenMeshAdapterTest.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// Pamgen includes
30#include "create_inline_mesh.h"
31
32using Teuchos::ParameterList;
33using Teuchos::ArrayRCP;
34
36
37/*****************************************************************************/
38/******************************** MAIN ***************************************/
39/*****************************************************************************/
40
41int main(int narg, char *arg[]) {
42
43 Tpetra::ScopeGuard tscope(&narg, &arg);
44 Teuchos::RCP<const Teuchos::Comm<int> > CommT = Tpetra::getDefaultComm();
45
46 int me = CommT->getRank();
47 int numProcs = CommT->getSize();
48
49 if (me == 0){
50 std::cout
51 << "====================================================================\n"
52 << "| |\n"
53 << "| Example: Partition Pamgen Hexahedral Mesh |\n"
54 << "| |\n"
55 << "| Questions? Contact Karen Devine (kddevin@sandia.gov), |\n"
56 << "| Erik Boman (egboman@sandia.gov), |\n"
57 << "| Siva Rajamanickam (srajama@sandia.gov). |\n"
58 << "| |\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"
62 << "| |\n"
63 << "====================================================================\n";
64 }
65
66
67#ifdef HAVE_MPI
68 if (me == 0) {
69 std::cout << "PARALLEL executable \n";
70 }
71#else
72 if (me == 0) {
73 std::cout << "SERIAL executable \n";
74 }
75#endif
76
77 /***************************************************************************/
78 /*************************** GET XML INPUTS ********************************/
79 /***************************************************************************/
80
81 // default values for command-line arguments
82 std::string xmlMeshInFileName("Poisson.xml");
83 std::string action("mj");
84 int nParts = CommT->getSize();
85
86 // Read run-time options.
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, "
92 "parma or color");
93 cmdp.setOption("nparts", &nParts,
94 "Number of parts to create");
95 cmdp.parse(narg, arg);
96
97 // Read xml file into parameter list
98 ParameterList inputMeshList;
99
100 if(xmlMeshInFileName.length()) {
101 if (me == 0) {
102 std::cout << "\nReading parameter list from the XML file \""
103 <<xmlMeshInFileName<<"\" ...\n\n";
104 }
105 Teuchos::updateParametersFromXmlFile(xmlMeshInFileName,
106 Teuchos::inoutArg(inputMeshList));
107 if (me == 0) {
108 inputMeshList.print(std::cout,2,true,true);
109 std::cout << "\n";
110 }
111 }
112 else {
113 std::cout << "Cannot read input file: " << xmlMeshInFileName << "\n";
114 return 5;
115 }
116
117 // Get pamgen mesh definition
118 std::string meshInput = Teuchos::getParameter<std::string>(inputMeshList,
119 "meshInput");
120
121 /***************************************************************************/
122 /********************** GET CELL TOPOLOGY **********************************/
123 /***************************************************************************/
124
125 // Get dimensions
126 int dim = 3;
127
128 /***************************************************************************/
129 /***************************** GENERATE MESH *******************************/
130 /***************************************************************************/
131
132 if (me == 0) std::cout << "Generating mesh ... \n\n";
133
134 // Generate mesh with Pamgen
135 long long maxInt = 9223372036854775807LL;
136 Create_Pamgen_Mesh(meshInput.c_str(), dim, me, numProcs, maxInt);
137
138 // Creating mesh adapter
139 if (me == 0) std::cout << "Creating mesh adapter ... \n\n";
140
141 typedef Zoltan2::PamgenMeshAdapter<basic_user_t> inputAdapter_t;
143
144 inputAdapter_t *ia = new inputAdapter_t(*CommT, "region");
145 ia->print(me);
146
147 // Set parameters for partitioning
148 if (me == 0) std::cout << "Creating parameter list ... \n\n";
149
150 Teuchos::ParameterList params("test params");
151 params.set("timer_output_stream" , "std::cout");
152
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); // bool parameter
161 }
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");
169 }
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");
177 }
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");
188 }
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");
201 }
202
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");
210 }
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");
219
220 }
221
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");
226 }
227
228 if(me == 0) std::cout << "Action: " << action << std::endl;
229 // create Partitioning problem
230 if (do_partitioning) {
231 if (me == 0) std::cout << "Creating partitioning problem ... \n\n";
232
233 Zoltan2::PartitioningProblem<inputAdapter_t> problem(ia, &params, CommT);
234
235 // call the partitioner
236 if (me == 0) std::cout << "Calling the partitioner ... \n\n";
237
238 problem.solve();
239
240 // create metric object
241
242 Teuchos::RCP<quality_t> metricObject =
243 rcp(new quality_t(ia, &params, CommT, &problem.getSolution()));
244
245 if (!me) {
246 metricObject->printMetrics(std::cout);
247 }
248 }
249 else {
250 if (me == 0) std::cout << "Creating coloring problem ... \n\n";
251
252 Zoltan2::ColoringProblem<inputAdapter_t> problem(ia, &params, CommT);
253
254 // call the partitioner
255 if (me == 0) std::cout << "Calling the coloring algorithm ... \n\n";
256
257 problem.solve();
258
259 problem.printTimers();
260 }
261
262 // delete mesh
263 if (me == 0) std::cout << "Deleting the mesh ... \n\n";
264
265 Delete_Pamgen_Mesh();
266
267 if (me == 0)
268 std::cout << "PASS" << std::endl;
269
270 Kokkos::finalize(); // Test uses directory with kokkos
271
272 return 0;
273}
274/*****************************************************************************/
275/********************************* END MAIN **********************************/
276/*****************************************************************************/
#define nParts
Defines the ColoringProblem class.
Defines the PamgenMeshAdapter class.
Defines the PartitioningProblem class.
int main()
A simple class that can be the User template argument for an InputAdapter.
ColoringProblem sets up coloring problems for the user.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
A class that computes and returns quality metrics.
This class represents a mesh.
PartitioningProblem sets up partitioning problems for the user.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
void printTimers() const
Return the communicator passed to the problem.
Zoltan2::BasicUserTypes basic_user_t
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t