Zoltan2
Loading...
Searching...
No Matches
Sphynx_Research_Driver.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
10#include "Teuchos_CommandLineProcessor.hpp"
11#include "Tpetra_CrsMatrix.hpp"
12#include "Tpetra_Core.hpp"
13#include "Tpetra_KokkosCompat_DefaultNode.hpp"
15
18
19#include "Teuchos_TimeMonitor.hpp"
20#include "Teuchos_StackedTimer.hpp"
22
23#include <Galeri_MultiVectorTraits.hpp>
24#include <Galeri_XpetraProblemFactory.hpp>
25#include <Galeri_XpetraParameters.hpp>
26#include <MatrixMarket_Tpetra.hpp>
27
29// This is a driver with many available options that can be used to test
30// a variety of features of Sphynx.
31// Note: This is research code. We do not guarantee it is without bugs.
33
34/* -------------------------------------------------------------------------
35 * Function: buildCrsMatrix
36 * Purpose: When the user does not input a matrix file, buildCrsMatrix will
37 * construct a Brick3D matrix in either 1, 2, or 3 dimensions. The
38 * Brick3D matrix is a 27-point difference stencil for the Laplace
39 * operator on a hex mesh.
40 * ------------------------------------------------------------------------- */
41
42template <typename lno_t, typename gno_t, typename scalar_t, typename nod_t>
43int buildCrsMatrix(int xdim, int ydim, int zdim, std::string problemType,
44 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
45 Teuchos::RCP<Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, nod_t>> &M_)
46{
47 /* Print size of the mesh being constructed */
48 if (comm->getRank() == 0){
49 std::cout << "Create matrix with " << problemType;
50 std::cout << " (a " << xdim;
51 if (zdim > 0)
52 std::cout << " x " << ydim << " x " << zdim << " ";
53 else if (ydim > 0)
54 std::cout << " x" << ydim << " x 1 ";
55 else
56 std::cout << "x 1 x 1 ";
57 std::cout << " mesh)" << std::endl;
58 }
59
60 Teuchos::CommandLineProcessor tclp;
61 Galeri::Xpetra::Parameters<gno_t> params(tclp, xdim, ydim, zdim, problemType);
62
63 Teuchos::RCP<const Tpetra::Map<lno_t, gno_t> > map =
64 Teuchos::rcp(new Tpetra::Map<lno_t, gno_t>(params.GetNumGlobalElements(), 0, comm));
65
66 /* Build the Brick3D matrix */
67 try{
68 Teuchos::RCP<Galeri::Xpetra::Problem<Tpetra::Map<lno_t, gno_t>,
69 Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, nod_t>,
70 Tpetra::MultiVector<scalar_t, lno_t, gno_t> > > Pr=
71 Galeri::Xpetra::BuildProblem<scalar_t, lno_t, gno_t,
72 Tpetra::Map<lno_t, gno_t>,
73 Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, nod_t>,
74 Tpetra::MultiVector<scalar_t, lno_t, gno_t, nod_t> >
75 (params.GetMatrixType(), map, params.GetParameterList());
76
77 M_ = Pr->BuildMatrix();
78 }
79 catch (std::exception &e) { // Probably not enough memory
80 if(comm->getRank() == 0) std::cout << "Error returned from Galeri " << e.what() << std::endl;
81 exit(-1);
82 }
83 if (M_.is_null())
84 return 1;
85 else
86 return 0;
87}
88
89template <typename adapter_type>
90 void
91compute_edgecut(Teuchos::RCP<adapter_type> &adapter,
93{
94 typedef typename adapter_type::user_t graph_type;
95 typedef typename graph_type::global_ordinal_type GO;
96 typedef typename graph_type::local_ordinal_type LO;
97 typedef typename graph_type::node_type NO;
98 typedef typename adapter_type::part_t PT;
99
100 using ordinal_view_t = Kokkos::View<GO*, typename NO::device_type>;
101 using part_view_t = Kokkos::View<PT*, typename NO::device_type>;
102
103 auto graph = adapter->getUserGraph();
104 auto rowMap = graph->getRowMap();
105 auto colMap = graph->getColMap();
106
107 size_t numLclRows = rowMap->getLocalNumElements();
108 size_t numGblRows = rowMap->getGlobalNumElements();
109 size_t numLclCols = colMap->getLocalNumElements();
110
111 ordinal_view_t colLocalToGlobal(Kokkos::view_alloc("colLocalToGlobal", Kokkos::WithoutInitializing), numLclCols);
112 auto colMapHost = Kokkos::create_mirror_view (Kokkos::HostSpace (), colLocalToGlobal);
113 for(size_t i = 0; i < numLclCols; ++i)
114 colMapHost[i] = colMap->getGlobalElement(i);
115 Kokkos::deep_copy (colLocalToGlobal, colMapHost);
116
117 ordinal_view_t rowLocalToGlobal(Kokkos::view_alloc("rowLocalToGlobal", Kokkos::WithoutInitializing), numLclRows);
118 auto rowMapHost = Kokkos::create_mirror_view (Kokkos::HostSpace (), rowLocalToGlobal);
119 for(size_t i = 0; i < numLclRows; ++i)
120 rowMapHost[i] = rowMap->getGlobalElement(i);
121 Kokkos::deep_copy (rowLocalToGlobal, rowMapHost);
122
123 part_view_t localParts("localParts", numGblRows);
124 part_view_t globalParts("globalParts", numGblRows);
125 auto localPartsHost = Kokkos::create_mirror_view(Kokkos::HostSpace(), localParts);
126
127 auto parts = solution.getPartListView();
128 for(size_t i = 0; i < numLclRows; i++){
129 GO gi = rowMap->getGlobalElement(i);
130 localPartsHost(gi) = parts[i];
131 }
132 Kokkos::deep_copy(localParts, localPartsHost);
133
134 auto comm = graph->getComm();
135 Teuchos::reduceAll<int, PT> (*comm, Teuchos::REDUCE_SUM, numGblRows, localParts.data(), globalParts.data());
136
137 auto rowPtr = graph->getLocalGraphHost().row_map;
138 auto colInd = graph->getLocalGraphHost().entries;
139
140 size_t localtotalcut = 0, totalcut = 0;
141
142 using execution_space = typename NO::device_type::execution_space;
143 using range_policy = Kokkos::RangePolicy<execution_space, Kokkos::IndexType<LO>>;
144 Kokkos::parallel_reduce("Compute cut", range_policy(0, numLclRows),
145 KOKKOS_LAMBDA(const LO i, size_t &cut){
146
147 const GO gRid = rowLocalToGlobal(i);
148 const PT gi = globalParts(gRid);
149
150 const size_t start = rowPtr(i);
151 const size_t end = rowPtr(i+1);
152 for(size_t j = start; j < end; ++j) {
153
154 const GO gCid = colLocalToGlobal(colInd(j));
155 PT gj = globalParts(gCid);
156 if(gi != gj)
157 cut += 1;
158 }
159 }, localtotalcut);
160
161 Teuchos::reduceAll (*comm, Teuchos::REDUCE_SUM, 1, &localtotalcut, &totalcut);
162
163 // compute imbalance
164 auto rowPtr_h = Kokkos::create_mirror_view(rowPtr);
165 Kokkos::deep_copy(rowPtr_h, rowPtr);
166 int nparts = (int)solution.getTargetGlobalNumberOfParts();
167
168 size_t *partw = new size_t[nparts];
169 size_t *partc = new size_t[nparts];
170
171 size_t *gpartw = new size_t[nparts];
172 size_t *gpartc = new size_t[nparts];
173
174 for(int i = 0; i < nparts; i++){
175 partw[i] = 0; partc[i] = 0;
176 gpartw[i] = 0; gpartc[i] = 0;
177 }
178
179 for(size_t i = 0; i < numLclRows; i++){
180 partw[parts[i]] += rowPtr_h(i+1) - rowPtr_h(i) - 1;
181 partc[parts[i]] ++;
182 }
183
184 Teuchos::reduceAll (*comm, Teuchos::REDUCE_SUM, nparts, partw, gpartw);
185 Teuchos::reduceAll (*comm, Teuchos::REDUCE_SUM, nparts, partc, gpartc);
186
187 size_t maxc = 0, totc = 0;
188 size_t maxw = 0, totw = 0;
189
190 for(int i = 0; i < nparts; i++){
191 if(gpartw[i] > maxw)
192 maxw = gpartw[i];
193 if(gpartc[i] > maxc)
194 maxc = gpartc[i];
195 totw += gpartw[i];
196 totc += gpartc[i];
197 }
198
199 double imbw = (double)maxw/((double)totw/nparts);
200 double imbc = (double)maxc/((double)totc/nparts);
201
202 if(comm->getRank() == 0) {
203
204 std::cout << "\n\n************************************************" << std::endl;
205 std::cout << " EDGECUT: " << totalcut << std::endl;
206 std::cout << " MAX/AVG WEIGHT: " << imbw << std::endl;
207 std::cout << " MAX/AVG COUNT: " << imbc << std::endl;
208 std::cout << "************************************************\n\n" << std::endl;
209
210 }
211
212}
213
214int main(int narg, char *arg[])
215{
216
217 Tpetra::ScopeGuard tpetraScope (&narg, &arg);
218 {
219
220 const Teuchos::RCP<const Teuchos::Comm<int>> pComm= Tpetra::getDefaultComm();
221
222 int me = pComm->getRank();
223
224 // Parameters
225 int nparts = 64;
226 int max_iters = 1000;
227 int block_size = -1;
228 int rand_seed =1;
229 int resFreq = 0;
230 int orthoFreq = 0;
231 std::string matrix_file = "";
232 std::string eigensolve = "LOBPCG";
233 bool parmetis = false;
234 bool pulp = false;
235
236 int verbosity = 1;
237
238 std::string ptype = "";
239 std::string prec = "jacobi";
240 std::string init = "random";
241 double tol = 1e-1;
242
243 // Echo the command line
244 if (me == 0) {
245 for (int i = 0; i < narg; i++)
246 std::cout << arg[i] << " ";
247 std::cout << std::endl;
248 }
249
250 Teuchos::CommandLineProcessor cmdp(false,true);
251 cmdp.setOption("matrix_file",&matrix_file,
252 "Path and filename of the matrix to be read.");
253 cmdp.setOption("nparts",&nparts,
254 "Number of global parts desired in the resulting partition.");
255 cmdp.setOption("rand_seed",&rand_seed,
256 "Seed for the random multivector.");
257 cmdp.setOption("max_iters",&max_iters,
258 "Maximum iters for the eigensolver.");
259 cmdp.setOption("block_size",&block_size,
260 "Block size.");
261 cmdp.setOption("verbosity", &verbosity,
262 "Verbosity level");
263 cmdp.setOption("parmetis", "sphynx", &parmetis,
264 "Whether to use parmetis.");
265 cmdp.setOption("pulp", "sphynx", &pulp,
266 "Whether to use pulp.");
267 cmdp.setOption("prec", &prec,
268 "Prec type to use.");
269 cmdp.setOption("eigensolve", &eigensolve,
270 "Eigensolver to use: LOBPCG, BlockDavidson, GeneralizedDavidson, BlockKrylovSchur or randomized.");
271 cmdp.setOption("prob", &ptype,
272 "Problem type to use. Options are combinatorial, normalized or generalized.");
273 cmdp.setOption("tol", &tol,
274 "Tolerance to use.");
275 cmdp.setOption("init", &init,
276 "Sphynx Initial guess. Options: random or constants. Default: random if randomized solver is used.");
277 cmdp.setOption("resFreq", &resFreq,
278 "(For randomized) Specify how often to check the residuals. Orthogonalization of the basis is also done.");
279 cmdp.setOption("orthoFreq", &orthoFreq,
280 "(For randomized) Specify how often to orthogonalize the basis.");
281
282 if (cmdp.parse(narg,arg)!=Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) return -1;
283
284 // Print the most essential options (not in the MyPL parameters later)
285 if (me == 0 ){
286 std::cout << std::endl << "--------------------------------------------------" << matrix_file << std::endl;
287 std::cout << "| Sphynx Parameters |" << matrix_file << std::endl;
288 std::cout << "--------------------------------------------------" << matrix_file << std::endl;
289 std::cout << " matrix file = " << matrix_file << std::endl;
290 std::cout << " nparts = " << nparts << std::endl;
291 std::cout << " verbosity = " << verbosity << std::endl;
292 std::cout << " parmetis = " << parmetis << std::endl;
293 std::cout << " pulp = " << pulp << std::endl;
294 std::cout << " prec = " << prec << std::endl;
295 std::cout << " eigensolver = " << eigensolve << std::endl;
296 std::cout << " prob = " << ptype << std::endl;
297 std::cout << " tol = " << tol << std::endl;
298 std::cout << " init = " << init << std::endl;
299 std::cout << " resFreq = " << resFreq << std::endl;
300 std::cout << " orthoFreq = " << orthoFreq << std::endl;
301 std::cout << "--------------------------------------------------" << matrix_file << std::endl << std::endl;
302 }
303
304 using scalar_type = Tpetra::Details::DefaultTypes::scalar_type;
305 using local_ordinal_type = Tpetra::Details::DefaultTypes::local_ordinal_type;
306 using global_ordinal_type = Tpetra::Details::DefaultTypes::global_ordinal_type;
307 using node_type = Tpetra::Details::DefaultTypes::node_type;
308
309 using crs_matrix_type = Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type>;
311 using solution_type = Zoltan2::PartitioningSolution<adapter_type>;
312
313 // Read the input matrix
314 Teuchos::RCP<adapter_type> adapter;
315 Teuchos::RCP<crs_matrix_type> tmatrix;
316
317 // Set the random seed and hope it goes through to Tpetra.
318 std::srand(rand_seed);
319
320 // Read in user-specified matrix or create default Brick3D matrix (100^3)
321 std::string mtx = ".mtx", mm = ".mm", lc = ".largestComp", lc2 = ".bin";
322 if(std::equal(lc.rbegin(), lc.rend(), matrix_file.rbegin()) || std::equal(lc2.rbegin(), lc2.rend(), matrix_file.rbegin())) {
323 tmatrix = readMatrixFromBinaryFile<crs_matrix_type>(matrix_file, pComm, true, verbosity>0);
324 if(me == 0 && verbosity > 1) std::cout << "Used Seher's reader for Largest Comp." << std::endl;
325 }
326 else if(std::equal(mtx.rbegin(), mtx.rend(), matrix_file.rbegin()) || std::equal(mm.rbegin(), mm.rend(), matrix_file.rbegin())) {
327 typedef Tpetra::MatrixMarket::Reader<crs_matrix_type> reader_type;
328 reader_type r;
329 tmatrix = r.readSparseFile(matrix_file, pComm);
330 if(me == 0 && verbosity > 1) std::cout << "Used standard Matrix Market reader." << std::endl;
331 }
332 else { // Build Brick3D matrix
333 /* If the user entered something that was not an .mtx or .largestComp file,
334 * check if it is a mesh size. If invalid input, default to 100^3 Brick3D.
335 */
336 int meshdim = 100;
337 if(matrix_file != "")
338 {
339 std::string::const_iterator it = matrix_file.begin();
340 while(it != matrix_file.end() && std::isdigit(*it)) ++it;
341 if(it == matrix_file.end())
342 {
343 meshdim = std::stoi(matrix_file);
344 }
345 else {
346 if(me == 0) std::cout << "Invalid matrix file entered. Reverting to default matrix." << std::endl;
347 }
348 }
349 int ierr = buildCrsMatrix<local_ordinal_type, global_ordinal_type, scalar_type>
350 (meshdim, meshdim, meshdim, "Brick3D", pComm, tmatrix);
351 if (ierr != 0) std::cout << "Error! Brick3D failed to build!" << std::endl;
352 }
353
354 if(me == 0 && verbosity > 0) std::cout << "Done with reading/creating the matrix." << std::endl;
355
356 Teuchos::RCP<const Tpetra::Map<> > map = tmatrix->getMap();
357
358 adapter = Teuchos::rcp(new adapter_type(tmatrix->getCrsGraph(), 1));
359 adapter->setVertexWeightIsDegree(0);
360
361 // Set the Sphynx parameters
362 Teuchos::RCP<Teuchos::ParameterList> params = Teuchos::rcp(new Teuchos::ParameterList());
363 Teuchos::RCP<Teuchos::ParameterList> sphynxParams(new Teuchos::ParameterList);
364 params->set("num_global_parts", nparts);
365
366 Teuchos::RCP<Teuchos::StackedTimer> stacked_timer;
367 stacked_timer = Teuchos::rcp(new Teuchos::StackedTimer("SphynxDriver"));
368 Teuchos::TimeMonitor::setStackedTimer(stacked_timer);
369 if(parmetis || pulp) {
370
371 params->set("partitioning_approach", "partition");
372 params->set("imbalance_tolerance", 1.01);
373 if(parmetis) {
374 params->set("algorithm", "parmetis");
375 params->set("imbalance_tolerance", 1.01);
376 }
377 else {
378 params->set("algorithm", "pulp");
379 params->set("pulp_vert_imbalance", 1.01);
380 }
381
383 Teuchos::RCP<problem_type> problem;
384 pComm->barrier();
385 {
386 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer("Partitioning::All"));
387 {
388 Teuchos::TimeMonitor t2(*Teuchos::TimeMonitor::getNewTimer("Partitioning::Problem"));
389 problem = Teuchos::rcp(new problem_type(adapter.getRawPtr(), params.getRawPtr(), Tpetra::getDefaultComm()));
390 }
391 {
392 Teuchos::TimeMonitor t3(*Teuchos::TimeMonitor::getNewTimer("Partitioning::Solve"));
393 problem->solve();
394 }
395 }
396 pComm->barrier();
397
398 solution_type solution = problem->getSolution();
399 compute_edgecut<adapter_type>(adapter, solution);
400
401 } else {
402
403 sphynxParams->set("sphynx_verbosity", verbosity);
404 sphynxParams->set("sphynx_max_iterations", max_iters);
405 sphynxParams->set("sphynx_ortho_freq", orthoFreq);
406 sphynxParams->set("sphynx_res_freq", resFreq);
407 sphynxParams->set("sphynx_skip_preprocessing", true);
408 sphynxParams->set("sphynx_eigensolver", eigensolve);
409 if (block_size > 0) sphynxParams->set("sphynx_block_size", block_size);
410 if (ptype != "") sphynxParams->set("sphynx_problem_type", ptype);
411 if (init != "") sphynxParams->set("sphynx_initial_guess", init);
412 if (prec != "") sphynxParams->set("sphynx_preconditioner_type", prec);
413 if (tol != -1) sphynxParams->set("sphynx_tolerance", tol);
414
415 using problem_type = Zoltan2::SphynxProblem<adapter_type>; //We found sphynx
416 Teuchos::RCP<problem_type> problem;
417 pComm->barrier();
418 {
419 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer("Partitioning::All"));
420 {
421 Teuchos::TimeMonitor t2(*Teuchos::TimeMonitor::getNewTimer("Partitioning::Problem"));
422 problem = Teuchos::rcp(new problem_type(adapter.get(), params.get(), sphynxParams, Tpetra::getDefaultComm()));
423 }
424 {
425 if(me == 0) std::cout << eigensolve << " will be used to solve the partitioning problem." << std::endl;
426 problem->solve();
427 }
428 pComm->barrier();
429 }
430 solution_type solution = problem->getSolution();
431 compute_edgecut<adapter_type>(adapter, solution);
432 }
433 stacked_timer->stopBaseTimer();
434
435 Teuchos::RCP<Teuchos::FancyOStream> fancy2 = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
436 Teuchos::FancyOStream& out2 = *fancy2;
437 Teuchos::StackedTimer::OutputOptions options;
438 options.output_fraction = options.output_histogram = options.output_minmax = true;
439 stacked_timer->report(out2, pComm, options);
440
441 Teuchos::TimeMonitor::summarize();
442
443 } //End Tpetra scope guard
444 return 0;
445}
int buildCrsMatrix(int xdim, int ydim, int zdim, std::string problemType, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, Teuchos::RCP< Tpetra::CrsMatrix< scalar_t, lno_t, gno_t, nod_t > > &M_)
void compute_edgecut(Teuchos::RCP< adapter_type > &adapter, Zoltan2::PartitioningSolution< adapter_type > &solution)
Defines the PartitioningProblem class.
Defines XpetraCrsGraphAdapter class.
int main()
PartitioningProblem sets up partitioning problems for the user.
A PartitioningSolution is a solution to a partitioning problem.
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
size_t getTargetGlobalNumberOfParts() const
Returns the global number of parts desired in the solution.
Provides access for Zoltan2 to Xpetra::CrsGraph data.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t