Zoltan2
Loading...
Searching...
No Matches
Zoltan2_Sphynx.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Sphynx
4//
5// Copyright 2020 NTESS and the Sphynx contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef _ZOLTAN2_SPHYNXALGORITHM_HPP_
11#define _ZOLTAN2_SPHYNXALGORITHM_HPP_
12
13
15// This file contains the Sphynx algorithm.
16//
17// Sphynx is a graph partitioning algorithm that is based on a spectral method.
18// It has three major steps:
19// 1) compute the Laplacian matrix of the input graph,
20// 2) compute logK+1 eigenvectors on the Laplacian matrix,
21// 3) use eigenvector entries as coordinates and compute a K-way partition on
22// them using a geometric method.
23//
24// Step1: Sphynx provides three eigenvalue problems and hence Laplacian matrix:
25// i) combinatorial (Lx = \lambdax, where L = A - D)
26// ii) generalized (Lx = \lambdaDx, where L = A - D)
27// iii) normalized (L_nx, \lambdax, where Ln = D^{-1/2}LD^{-1/2}
28// and L = A - D)
29//
30// Step2: Sphynx calls the LOBPCG algorithm provided in Anasazi to obtain
31// logK+1 eigenvectors.
32// (Alternately, uses the experimental Randomized eigensolver.)
33// Step3: Sphynx calls the MJ algorithm provided in Zoltan2Core to compute the
34// partition.
36
37#include "Zoltan2Sphynx_config.h"
38
41
44
45#include "AnasaziLOBPCGSolMgr.hpp"
46#include "AnasaziBlockDavidsonSolMgr.hpp"
47#include "AnasaziGeneralizedDavidsonSolMgr.hpp"
48#include "AnasaziBlockKrylovSchurSolMgr.hpp"
49#include "AnasaziRandomizedSolMgr.hpp"
50#include "AnasaziBasicEigenproblem.hpp"
51#include "AnasaziTpetraAdapter.hpp"
52
53#include "BelosLinearProblem.hpp"
54#include "BelosTpetraOperator.hpp"
55
56#include "Ifpack2_Factory.hpp"
57
58#include "Teuchos_TimeMonitor.hpp"
59
60#ifdef HAVE_ZOLTAN2SPHYNX_MUELU
61#include "MueLu_CreateTpetraPreconditioner.hpp"
62#endif
63
64namespace Zoltan2 {
65
66 template <typename Adapter>
67 class Sphynx : public Algorithm<Adapter>
68 {
69
70 public:
71
72 using scalar_t = double; // Sphynx with scalar_t=double obtains better cutsize
73 using lno_t = typename Adapter::lno_t;
74 using gno_t = typename Adapter::gno_t;
75 using node_t = typename Adapter::node_t;
76 using offset_t = typename Adapter::offset_t;
77 using part_t = typename Adapter::part_t;
78 using weight_t = typename Adapter::scalar_t;
79
80 using graph_t = Tpetra::CrsGraph<lno_t, gno_t, node_t>;
81 using matrix_t = Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t>;
82 using mvector_t = Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t>;
83 using op_t = Tpetra::Operator<scalar_t, lno_t, gno_t, node_t>;
84
86 typedef Anasazi::MultiVecTraits<scalar_t,mvector_t> MVT;
87
91
92 // Takes the graph from the input adapter and computes the Laplacian matrix
93 Sphynx(const RCP<const Environment> &env,
94 const RCP<Teuchos::ParameterList> &params,
95 const RCP<Teuchos::ParameterList> &sphynxParams,
96 const RCP<const Comm<int> > &comm,
97 const RCP<const XpetraCrsGraphAdapter<graph_t> > &adapter) :
98 env_(env), params_(params), sphynxParams_(sphynxParams), comm_(comm), adapter_(adapter)
99 {
100
101 // Don't compute the Laplacian if the number of requested parts is 1
102 const Teuchos::ParameterEntry *pe = params_->getEntryPtr("num_global_parts");
103 numGlobalParts_ = pe->getValue<int>(&numGlobalParts_);
104 if(numGlobalParts_ > 1){
105
106 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer("Sphynx::Laplacian"));
107
108 // The verbosity is common throughout the algorithm, better to check and set now
109 pe = sphynxParams_->getEntryPtr("sphynx_verbosity");
110 if (pe)
111 verbosity_ = pe->getValue<int>(&verbosity_);
112
113 // Do we need to pre-process the input?
114 pe = sphynxParams_->getEntryPtr("sphynx_skip_preprocessing");
115 if (pe)
116 skipPrep_ = pe->getValue<bool>(&skipPrep_);
117
118 // Get the graph from XpetraCrsGraphAdapter if skipPrep_ is true
119 // We assume the graph has all of the symmetric and diagonal entries
120 if(skipPrep_)
121 graph_ = adapter_->getUserGraph();
122 else {
123 throw std::runtime_error("\nSphynx Error: Preprocessing has not been implemented yet.\n");
124 }
125
126 // Check if the graph is regular
128
129 // Set Sphynx defaults: preconditioner, problem type, tolerance, initial vectors.
130 setDefaults();
131
132 // Compute the Laplacian matrix
134
135 if(problemType_ == GENERALIZED)
137
138 }
139 }
140
144
145 void partition(const Teuchos::RCP<PartitioningSolution<Adapter> > &solution);
146
147 int AnasaziWrapper(const int numEigenVectors);
148
149 template<typename problem_t>
150 void setPreconditioner(Teuchos::RCP<problem_t> &problem);
151
152 template<typename problem_t>
153 void setMueLuPreconditioner(Teuchos::RCP<problem_t> &problem);
154
155 template<typename problem_t>
156 void setJacobiPreconditioner(Teuchos::RCP<problem_t> &problem);
157
158 template<typename problem_t>
159 void setPolynomialPreconditioner(Teuchos::RCP<problem_t> &problem);
160
161 void eigenvecsToCoords(Teuchos::RCP<mvector_t> &eigenVectors,
162 int computedNumEv,
163 Teuchos::RCP<mvector_t> &coordinates);
164
165 void computeWeights(std::vector<const weight_t *> vecweights,
166 std::vector<int> strides);
167
168 void MJwrapper(const Teuchos::RCP<const mvector_t> &coordinates,
169 std::vector<const weight_t *> weights,
170 std::vector<int> strides,
171 const Teuchos::RCP<PartitioningSolution<Adapter>> &solution);
172
173 void setUserEigenvectors(const Teuchos::RCP<mvector_t> &userEvects);
177
179 // Determine input graph's regularity = maxDegree/AvgDegree < 10.
180 // If Laplacian type is not specified, then use combinatorial for regular
181 // and generalized for irregular.
182 // MueLu settings will differ depending on the regularity, too.
184 {
185 // Get the row pointers in the host
186 auto rowOffsets = graph_->getLocalGraphDevice().row_map;
187 auto rowOffsets_h = Kokkos::create_mirror_view(rowOffsets);
188 Kokkos::deep_copy(rowOffsets_h, rowOffsets);
189
190 // Get size information
191 const size_t numGlobalEntries = graph_->getGlobalNumEntries();
192 const size_t numLocalRows = graph_->getLocalNumRows();
193 const size_t numGlobalRows = graph_->getGlobalNumRows();
194
195 // Compute local maximum degree
196 size_t localMax = 0;
197 for(size_t i = 0; i < numLocalRows; i++){
198 if(rowOffsets_h(i+1) - rowOffsets_h(i) - 1 > localMax)
199 localMax = rowOffsets_h(i+1) - rowOffsets_h(i) - 1;
200 }
201
202 // Compute global maximum degree
203 size_t globalMax = 0;
204 Teuchos::reduceAll<int, size_t> (*comm_, Teuchos::REDUCE_MAX, 1, &localMax, &globalMax);
205
206 double avg = static_cast<double>(numGlobalEntries-numGlobalRows)/numGlobalRows;
207 double maxOverAvg = static_cast<double>(globalMax)/avg;
208
209 // Use generalized Laplacian on irregular graphs
210 irregular_ = false;
211 if(maxOverAvg > 10) {
212 irregular_ = true;
213 }
214
215 // Let the user know about what we determined if verbose
216 if(verbosity_ > 0) {
217 if(comm_->getRank() == 0) {
218 std::cout << "Regularity of Graph ----------------" << std::endl;
219 std::cout << " Maximum degree: " << globalMax << std::endl;
220 std::cout << " Average degree: " << avg << std::endl;
221 std::cout << " Max/Avg: " << globalMax/avg << std::endl;
222 std::cout << " Regular graph?: " << !irregular_ << std::endl;
223 }
224 }
225 }
226
228 // If preconditioner type is not specified:
229 // use muelu if it is enabled, and jacobi otherwise.
230 // If eigenvalue problem type is not specified:
231 // use combinatorial for regular and
232 // normalized for irregular with polynomial preconditioner,
233 // generalized for irregular with other preconditioners.
234 // If convergence tolerance is not specified:
235 // use 1e-3 for regular with jacobi and polynomial, and 1e-2 otherwise.
236 // If how to decide the initial vectors is not specified:
237 // use random for regular and constant for irregular
239 {
240
241 // Set the default preconditioner to muelu if it is enabled, jacobi otherwise.
242 precType_ = "jacobi";
243#ifdef HAVE_ZOLTAN2SPHYNX_MUELU
244 precType_ = "muelu";
245#endif
246
247 // Override the preconditioner type with the user's preference
248 const Teuchos::ParameterEntry *pe = sphynxParams_->getEntryPtr("sphynx_preconditioner_type");
249 if (pe) {
250 precType_ = pe->getValue<std::string>(&precType_);
251 if(precType_ != "muelu" && precType_ != "jacobi" && precType_ != "polynomial")
252 throw std::runtime_error("\nSphynx Error: " + precType_ + " is not recognized as a preconditioner.\n"
253 + " Possible values: muelu (if enabled), jacobi, and polynomial\n");
254 }
255
256 solverType_ = sphynxParams_->get("sphynx_eigensolver","LOBPCG");
257 TEUCHOS_TEST_FOR_EXCEPTION(!(solverType_ == "LOBPCG" || solverType_ == "randomized" || solverType_ == "BlockDavidson" || solverType_ == "GeneralizedDavidson" || solverType_ == "BlockKrylovSchur" ),
258 std::invalid_argument, "Sphynx: sphynx_eigensolver must be set to LOBPCG, randomized, BlockDavidson, GeneralizedDavidson, or BlockKrylovSchur.");
259
260 // Set the default problem type
261 problemType_ = COMBINATORIAL;
262 if(irregular_) {
263 if(precType_ == "polynomial")
264 problemType_ = NORMALIZED;
265 else
266 problemType_ = GENERALIZED;
267 }
268
269 // Override the problem type with the user's preference
270 pe = sphynxParams_->getEntryPtr("sphynx_problem_type");
271 if (pe) {
272 std::string probType = "";
273 probType = pe->getValue<std::string>(&probType);
274
275 if(probType == "combinatorial")
276 problemType_ = COMBINATORIAL;
277 else if(probType == "generalized")
278 problemType_ = GENERALIZED;
279 else if(probType == "normalized")
280 problemType_ = NORMALIZED;
281 else
282 throw std::runtime_error("\nSphynx Error: " + probType + " is not recognized as a problem type.\n"
283 + " Possible values: combinatorial, generalized, and normalized.\n");
284 }
285
286
287 // Set the default for tolerance
288 tolerance_ = 1.0e-2;
289 if(!irregular_ && (precType_ == "jacobi" || precType_ == "polynomial"))
290 tolerance_ = 1.0e-3;
291
292
293 // Override the tolerance with the user's preference
294 pe = sphynxParams_->getEntryPtr("sphynx_tolerance");
295 if (pe)
296 tolerance_ = pe->getValue<scalar_t>(&tolerance_);
297
298
299 // Set the default for initial vectors
300 randomInit_ = true;
301 if(irregular_)
302 randomInit_ = false;
303
304 // Override the initialization method with the user's preference
305 pe = sphynxParams_->getEntryPtr("sphynx_initial_guess");
306 if (pe) {
307 std::string initialGuess = "";
308 initialGuess = pe->getValue<std::string>(&initialGuess);
309
310 if(initialGuess == "random")
311 randomInit_ = true;
312 else if(initialGuess == "constants")
313 randomInit_ = false;
314 else
315 throw std::runtime_error("\nSphynx Error: " + initialGuess + " is not recognized as an option for initial guess.\n"
316 + " Possible values: random and constants.\n");
317 }
318
319 }
320
322 // Compute the Laplacian matrix depending on the eigenvalue problem type.
323 // There are 3 options for the type: combinatorial, generalized, and normalized.
324 // Combinatorial and generalized share the same Laplacian but generalized
325 // also needs a degree matrix.
327 {
328 if(solverType_ == "randomized")
329 laplacian_ = computeNormalizedLaplacian(true);
330 else if(problemType_ == NORMALIZED)
331 laplacian_ = computeNormalizedLaplacian();
332 else
333 laplacian_ = computeCombinatorialLaplacian();
334 }
335
337 // Compute a diagonal matrix with the vertex degrees in the input graph
339 {
340
341 // Get the row pointers in the host
342 auto rowOffsets = graph_->getLocalGraphHost().row_map;
343
344 // Create the degree matrix with max row size set to 1
345 Teuchos::RCP<matrix_t> degMat(new matrix_t (graph_->getRowMap(),
346 graph_->getRowMap(),
347 1));
348
349 scalar_t *val = new scalar_t[1];
350 lno_t *ind = new lno_t[1];
351 lno_t numRows = static_cast<lno_t>(graph_->getLocalNumRows());
352
353 // Insert the diagonal entries as the degrees
354 for (lno_t i = 0; i < numRows; ++i) {
355 val[0] = rowOffsets(i+1) - rowOffsets(i) - 1;
356 ind[0] = i;
357 degMat->insertLocalValues(i, 1, val, ind);
358 }
359 delete [] val;
360 delete [] ind;
361
362 degMat->fillComplete(graph_->getDomainMap(), graph_->getRangeMap());
363
364 degMatrix_ = degMat;
365 }
366
368 // Compute the combinatorial Laplacian: L = D - A.
369 // l_ij = degree of vertex i if i = j
370 // l_ij = -1 if i != j and a_ij != 0
371 // l_ij = 0 if i != j and a_ij = 0
372 Teuchos::RCP<matrix_t> computeCombinatorialLaplacian()
373 {
374 using range_policy = Kokkos::RangePolicy<
375 typename node_t::device_type::execution_space, Kokkos::IndexType<lno_t>>;
376 using values_view_t = Kokkos::View<scalar_t*, typename node_t::device_type>;
377 using offset_view_t = Kokkos::View<size_t*, typename node_t::device_type>;
378
379 const size_t numEnt = graph_->getLocalNumEntries();
380 const size_t numRows = graph_->getLocalNumRows();
381
382 // Create new values for the Laplacian, initialize to -1
383 values_view_t newVal (Kokkos::view_alloc("CombLapl::val", Kokkos::WithoutInitializing), numEnt);
384 Kokkos::deep_copy(newVal, -1);
385
386 // Get the diagonal offsets
387 offset_view_t diagOffsets(Kokkos::view_alloc("Diag Offsets", Kokkos::WithoutInitializing), numRows);
388 graph_->getLocalDiagOffsets(diagOffsets);
389
390 // Get the row pointers in the host
391 auto rowOffsets = graph_->getLocalGraphDevice().row_map;
392
393 // Compute the diagonal entries as the vertex degrees
394 Kokkos::parallel_for("Combinatorial Laplacian", range_policy(0, numRows),
395 KOKKOS_LAMBDA(const lno_t i){
396 newVal(rowOffsets(i) + diagOffsets(i)) = rowOffsets(i+1) - rowOffsets(i) - 1;
397 }
398 );
399 Kokkos::fence ();
400
401 // Create the Laplacian matrix using the input graph and with the new values
402 Teuchos::RCP<matrix_t> laplacian (new matrix_t(graph_, newVal));
403 laplacian->fillComplete (graph_->getDomainMap(), graph_->getRangeMap());
404
405 // Create the Laplacian maatrix using the input graph and with the new values
406 return laplacian;
407 }
408
410 // For AHat = false:
411 // Compute the normalized Laplacian: L_N = D^{-1/2} L D^{-1/2}, where L = D - A.
412 // l_ij = 1
413 // l_ij = -1/(sqrt(deg(v_i))sqrt(deg(v_j)) if i != j and a_ij != 0
414 // l_ij = 0 if i != j and a_ij = 0
415 //
416 // For AHat = true:
417 // AHat is turned to true if (and only if) using the randomized Eigensolver.
418 // For the randomized Eigensolver, we find eigenvalues of the matrix
419 // AHat = 2*I - L_N, and this is the matrix computed by this function.
420 Teuchos::RCP<matrix_t> computeNormalizedLaplacian(bool AHat = false)
421 {
422 using range_policy = Kokkos::RangePolicy<
423 typename node_t::device_type::execution_space, Kokkos::IndexType<lno_t>>;
424 using values_view_t = Kokkos::View<scalar_t*, typename node_t::device_type>;
425 using offset_view_t = Kokkos::View<size_t*, typename node_t::device_type>;
426 using vector_t = Tpetra::Vector<scalar_t, lno_t, gno_t, node_t>;
427 using dual_view_t = typename vector_t::dual_view_type;
428 using KAT = KokkosKernels::ArithTraits<scalar_t>;
429
430 const size_t numEnt = graph_->getLocalNumEntries();
431 const size_t numRows = graph_->getLocalNumRows();
432
433 // Create new values for the Laplacian, initialize to -1
434 values_view_t newVal (Kokkos::view_alloc("NormLapl::val", Kokkos::WithoutInitializing), numEnt);
435 if(AHat){
436 Kokkos::deep_copy(newVal, 1);
437 }
438 else{
439 Kokkos::deep_copy(newVal, -1);
440 }
441
442 // D^{-1/2}
443 dual_view_t dv ("MV::DualView", numRows, 1);
444 auto deginvsqrt = dv.view_device();
445
446 // Get the diagonal offsets
447 offset_view_t diagOffsets(Kokkos::view_alloc("Diag Offsets", Kokkos::WithoutInitializing), numRows);
448 graph_->getLocalDiagOffsets(diagOffsets);
449
450 // Get the row pointers
451 auto rowOffsets = graph_->getLocalGraphDevice().row_map;
452
453 // if(!AHat){
454 // Compute the diagonal entries as the vertex degrees
455 Kokkos::parallel_for("Combinatorial Laplacian", range_policy(0, numRows),
456 KOKKOS_LAMBDA(const lno_t i){
457 newVal(rowOffsets(i) + diagOffsets(i)) = rowOffsets(i+1) - rowOffsets(i) - 1;
458 deginvsqrt(i,0) = 1.0/KAT::sqrt(rowOffsets(i+1) - rowOffsets(i) - 1);
459 }
460 );
461 Kokkos::fence ();
462
463 // Create the Laplacian graph using the same graph structure with the new values
464 Teuchos::RCP<matrix_t> laplacian (new matrix_t(graph_, newVal));
465 laplacian->fillComplete (graph_->getDomainMap(), graph_->getRangeMap());
466
467 // Create the vector for D^{-1/2}
468 vector_t degInvSqrt(graph_->getRowMap(), dv);
469
470 // Normalize the laplacian matrix as D^{-1/2} L D^{-1/2}
471 laplacian->leftScale(degInvSqrt);
472 laplacian->rightScale(degInvSqrt);
473
474 return laplacian;
475 }
476
480
481 private:
482
483 // User-provided members
484 const Teuchos::RCP<const Environment> env_;
485 Teuchos::RCP<Teuchos::ParameterList> params_;
486 Teuchos::RCP<Teuchos::ParameterList> sphynxParams_;
487 const Teuchos::RCP<const Teuchos::Comm<int> > comm_;
488 const Teuchos::RCP<const Adapter> adapter_;
489
490 // Internal members
491 int numGlobalParts_;
492 Teuchos::RCP<const graph_t> graph_;
493 Teuchos::RCP<matrix_t> laplacian_;
494 Teuchos::RCP<const matrix_t> degMatrix_;
495 Teuchos::RCP<mvector_t> eigenVectors_;
496
497 bool irregular_; // decided internally
498 std::string precType_; // obtained from user params or decided internally
499 std::string solverType_; // obtained from user params or decided internally
500 problemType problemType_; // obtained from user params or decided internally
501 double tolerance_; // obtained from user params or decided internally
502 bool randomInit_; // obtained from user params or decided internally
503 int verbosity_; // obtained from user params
504 bool skipPrep_; // obtained from user params
505 };
506
510
512 // Allows the user to manually set eigenvectors for the Sphynx partitioner
513 // to use rather than solving for them with Anasazi. Mainly intended
514 // for debugging purposes.
516 template <typename Adapter>
517 void Sphynx<Adapter>::setUserEigenvectors(const Teuchos::RCP<mvector_t> &userEvects)
518 {
519 eigenVectors_ = userEvects;
520 }
521
523 // Compute a partition using the Laplacian matrix (and possibly the degree
524 // matrix as well). First call LOBPCG (or Randomized solver)
525 // to compute logK+1 eigenvectors, then
526 // transform the eigenvectors to coordinates, and finally call MJ to compute
527 // a partition on the coordinates.
528 template <typename Adapter>
530 {
531 // Return a trivial solution if only one part is requested
532 if(numGlobalParts_ == 1) {
533
534 size_t numRows =adapter_->getUserGraph()->getLocalNumRows();
535 Teuchos::ArrayRCP<part_t> parts(numRows,0);
536 solution->setParts(parts);
537
538 return;
539
540 }
541
542 // The number of eigenvectors to be computed
543 int numEigenVectors = (int) log2(numGlobalParts_)+1;
544 int computedNumEv;
545
546 if(eigenVectors_ == Teuchos::null){
547 // Compute the eigenvectors using LOBPCG
548 // or Randomized eigensolver
549 computedNumEv = Sphynx::AnasaziWrapper(numEigenVectors);
550
551 if(computedNumEv <= 1 && (solverType_ == "LOBPCG" || solverType_ == "GeneralizedDavidson" || solverType_ == "BlockDavidson" || solverType_ == "BlockKrylovSchur")) {
552 throw
553 std::runtime_error("\nAnasazi Error: LOBPCGSolMgr::solve() returned unconverged.\n"
554 "Sphynx Error: LOBPCG could not compute any eigenvectors.\n"
555 " Increase either max iters or tolerance.\n");
556 }
557 }
558 else{
559 // Use eigenvectors provided by user.
560 computedNumEv = (int) eigenVectors_->getNumVectors();
561 if(computedNumEv <= numEigenVectors) {
562 throw
563 std::runtime_error("\nSphynx Error: Number of eigenvectors given by user\n"
564 " is less than number of Eigenvectors needed for partition." );
565 }
566 }
567
568 // Transform the eigenvectors into coordinates
569 Teuchos::RCP<mvector_t> coordinates;
570
571 Sphynx::eigenvecsToCoords(eigenVectors_, computedNumEv, coordinates);
572
573 // Get the weights from the adapter
574 std::vector<const weight_t *> weights;
575 std::vector<int> wstrides;
577
578
579 // Compute the partition using MJ on coordinates
580 Sphynx::MJwrapper(coordinates, weights, wstrides, solution);
581
582 }
583
585 // Call LOBPCG on the Laplacian matrix.
586 // Or use the randomized eigensolver.
587 template <typename Adapter>
588 int Sphynx<Adapter>::AnasaziWrapper(const int numEigenVectors)
589 {
590
591 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer("Sphynx::Anasazi"));
592
593 // Set defaults for the parameters
594 // and get user-set values.
595 std::string which = (solverType_ == "randomized" ? "LM" : "SR");
596 std::string ortho = "SVQB";
597 bool relConvTol = false;
598 int maxIterations = sphynxParams_->get("sphynx_max_iterations",1000);
599 int blockSize = sphynxParams_->get("sphynx_block_size",numEigenVectors);
600 int orthoFreq = sphynxParams_->get("sphynx_ortho_freq", 0);
601 int resFreq = sphynxParams_->get("sphynx_res_freq", 0);
602 bool isHermitian = true;
603 bool relLockTol = false;
604 bool lock = false;
605 bool useFullOrtho = sphynxParams_->get("sphynx_use_full_ortho",true);
606
607 // Information to output in a verbose run
608 int numfailed = 0;
609 int iter = 0;
610
611 // Set Anasazi verbosity level
612 int anasaziVerbosity = Anasazi::Errors + Anasazi::Warnings;
613 if (verbosity_ >= 1) // low
614 anasaziVerbosity += Anasazi::FinalSummary + Anasazi::TimingDetails;
615 if (verbosity_ >= 2) // medium
616 anasaziVerbosity += Anasazi::IterationDetails;
617 if (verbosity_ >= 3) // high
618 anasaziVerbosity += Anasazi::StatusTestDetails + Anasazi::OrthoDetails
619 + Anasazi::Debug;
620
621 // Create the parameter list to pass into solver
622 Teuchos::ParameterList anasaziParams;
623 anasaziParams.set("Verbosity", anasaziVerbosity);
624 anasaziParams.set("Which", which);
625 anasaziParams.set("Convergence Tolerance", tolerance_);
626 anasaziParams.set("Maximum Iterations", maxIterations);
627 anasaziParams.set("Block Size", blockSize);
628 anasaziParams.set("Relative Convergence Tolerance", relConvTol);
629 anasaziParams.set("Orthogonalization", ortho);
630 anasaziParams.set("Use Locking", lock);
631 anasaziParams.set("Relative Locking Tolerance", relLockTol);
632 anasaziParams.set("Full Ortho", useFullOrtho);
633 anasaziParams.set("Orthogonalization Frequency", orthoFreq);
634 anasaziParams.set("Residual Frequency", resFreq);
635
636 if(solverType_ == "GeneralizedDavidson" || solverType_ == "BlockKrylovSchur" || solverType_ == "BlockDavidson"){
637 anasaziParams.set( "Num Blocks", maxIterations+1 );
638 anasaziParams.set( "Maximum Restarts", 0 );
639 anasaziParams.set( "Maximum Subspace Dimension", (maxIterations+1)*blockSize );
640 }
641
642 // Create and set initial vectors
643 auto map = laplacian_->getRangeMap();
644 Teuchos::RCP<mvector_t> ivec( new mvector_t(map, numEigenVectors));
645
646 if (randomInit_) {
647 // 0-th vector constant 1, rest random
648 Anasazi::MultiVecTraits<scalar_t, mvector_t>::MvRandom(*ivec);
649 ivec->getVectorNonConst(0)->putScalar(1.);
650 }
651 else { // This implies we will use constant initial vectors.
652 // 0-th vector constant 1, other vectors constant per block
653 // Create numEigenVectors blocks, but only use numEigenVectors-1 of them.
654 // This assures orthogonality.
655 ivec->getVectorNonConst(0)->putScalar(1.);
656 for (int j = 1; j < numEigenVectors; j++)
657 ivec->getVectorNonConst(j)->putScalar(0.);
658
659 gno_t blkSize = map->getGlobalNumElements() / numEigenVectors;
660 TEUCHOS_TEST_FOR_EXCEPTION(blkSize <= 0, std::runtime_error, "Blocksize too small for \"constants\" initial guess. Try \"random\".");
661
662 for (size_t lid = 0; lid < ivec->getLocalLength(); lid++) {
663 gno_t gid = map->getGlobalElement(lid);
664 for (int j = 1; j < numEigenVectors; j++){
665 if (((j-1)*blkSize <= gid) && (j*blkSize > gid))
666 ivec->replaceLocalValue(lid,j,1.);
667 }
668 }
669 }
670
671 // Create the eigenproblem to be solved
672 using problem_t = Anasazi::BasicEigenproblem<scalar_t, mvector_t, op_t>;
673 Teuchos::RCP<problem_t> problem (new problem_t(laplacian_, ivec));
674 problem->setHermitian(isHermitian);
675 problem->setNEV(numEigenVectors);
676
677 if(solverType_ != "randomized"){
678 // Set preconditioner
680 if(problemType_ == Sphynx::GENERALIZED) problem->setM(degMatrix_);
681 }
682
683 // Inform the eigenproblem that you are finished passing it information
684 bool boolret = problem->setProblem();
685 if (boolret != true) {
686 throw std::runtime_error("\nAnasazi::BasicEigenproblem::setProblem() returned with error.\n");
687 }
688 // Set Eigensolver
689 Teuchos::RCP<Anasazi::SolverManager<scalar_t, mvector_t, op_t>> solver;
690
691 if(solverType_ == "LOBPCG"){
692 solver = Teuchos::rcp(new Anasazi::LOBPCGSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
693 }
694 else if(solverType_ == "BlockDavidson"){
695 solver = Teuchos::rcp(new Anasazi::BlockDavidsonSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
696 }
697 else if(solverType_ == "GeneralizedDavidson"){
698 solver = Teuchos::rcp(new Anasazi::GeneralizedDavidsonSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
699 }
700 else if(solverType_ == "BlockKrylovSchur"){
701 solver = Teuchos::rcp(new Anasazi::BlockKrylovSchurSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
702 }
703 else{
704 solver = Teuchos::rcp(new Anasazi::Experimental::RandomizedSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
705 //numFailed = solver->getNumFailed();
706 }
707
708 if (verbosity_ > 0 && comm_->getRank() == 0)
709 anasaziParams.print(std::cout);
710
711 // Solve the problem
712 Anasazi::ReturnType returnCode = solver->solve();
713
714 // Check convergence, niters, and solvetime
715 iter = solver->getNumIters();
716 //solvetime = (solver->getTimers()[0])->totalElapsedTime();
717
718 // Retrieve the solution
719 using solution_t = Anasazi::Eigensolution<scalar_t, mvector_t>;
720 solution_t sol = problem->getSolution();
721 eigenVectors_ = sol.Evecs;
722 int numev = sol.numVecs;
723 std::vector<Anasazi::Value<double>> eigenValues_ = sol.Evals;
724
725 // Summarize iteration counts and solve time
726 if (verbosity_ > 0 && comm_->getRank() == 0) {
727 std::cout << std::endl;
728 std::cout << "ANASAZI SUMMARY" << std::endl;
729 std::cout << "Failed to converge: " << numfailed << std::endl;
730 std::cout << "No of iterations : " << iter << std::endl;
731 std::cout << "Solve time: " << std::endl; // solvetime << std::endl;
732 std::cout << "No of comp. vecs. : " << numev << std::endl;
733 }
734
735 std::cout << "Solver type: " << solverType_ << std::endl;
736 // Compute residuals (LOBPCG does this internally)
737 if(solverType_ == "randomized") {
738 std::vector<double> normR(numev);
739 mvector_t Aevec (map, numev);
740
741 if (numev > 0) {
742 Teuchos::SerialDenseMatrix<int,double> T (numev, numev);
743 T.putScalar(0.0);
744 for (int i = 0; i < numev; ++i) T(i,i) = eigenValues_[i].realpart;
745 laplacian_->apply(*eigenVectors_, Aevec);
746 MVT::MvTimesMatAddMv(-1.0, *eigenVectors_, T, 1.0, Aevec);
747 MVT::MvNorm(Aevec, normR);
748 } // End residual computation
749
750 if(comm_->getRank() == 0 && verbosity_ > 0) {
751 std::cout << std::endl << "Solver manager returned "
752 << (returnCode == Anasazi::Converged ? "converged." : "unconverged.")
753 << std::endl << std::endl
754 << std::setw(16) << "Eigenvalue"
755 << std::setw(18) << "Direct Residual"
756 << std::endl
757 << "------------------------------------------------------" << std::endl;
758
759 for (int i = 0; i < numev; ++i) {
760 std::cout << std::setw(16) << 2.0-eigenValues_[i].realpart
761 << std::setw(18) << normR[i] / eigenValues_[i].realpart
762 << std::endl;
763 }
764 std::cout << "------------------------------------------------------" << std::endl;
765 }
766 }
767
768 return numev;
769
770 }
771
773 template <typename Adapter>
774 template <typename problem_t>
775 void Sphynx<Adapter>::setPreconditioner(Teuchos::RCP<problem_t> &problem)
776 {
777 if(comm_->getRank() == 0) std::cout << "precType_ is: " << precType_ << std::endl;
778 // Set the preconditioner
779 if(precType_ == "muelu") {
781 }
782 else if(precType_ == "polynomial") {
784 }
785 else if(precType_ == "jacobi") {
787 }
788 // else: do we want a case where no preconditioning is applied?
789 }
790
792 template <typename Adapter>
793 template <typename problem_t>
794 void Sphynx<Adapter>::setMueLuPreconditioner(Teuchos::RCP<problem_t> &problem)
795 {
796#ifdef HAVE_ZOLTAN2SPHYNX_MUELU
797 Teuchos::ParameterList paramList;
798 if(verbosity_ == 0)
799 paramList.set("verbosity", "none");
800 else if(verbosity_ == 1)
801 paramList.set("verbosity", "low");
802 else if(verbosity_ == 2)
803 paramList.set("verbosity", "medium");
804 else if(verbosity_ == 3)
805 paramList.set("verbosity", "high");
806 else if(verbosity_ >= 4)
807 paramList.set("verbosity", "extreme");
808
809 paramList.set("smoother: type", "CHEBYSHEV");
810 Teuchos::ParameterList smootherParamList;
811 smootherParamList.set("chebyshev: degree", 3);
812 smootherParamList.set("chebyshev: ratio eigenvalue", 7.0);
813 smootherParamList.set("chebyshev: eigenvalue max iterations", irregular_ ? 100 : 10);
814 paramList.set("smoother: params", smootherParamList);
815 paramList.set("use kokkos refactor", true);
816 paramList.set("transpose: use implicit", true);
817
818 if(irregular_) {
819
820 paramList.set("multigrid algorithm", "unsmoothed");
821
822 paramList.set("coarse: type", "CHEBYSHEV");
823 Teuchos::ParameterList coarseParamList;
824 coarseParamList.set("chebyshev: degree", 3);
825 coarseParamList.set("chebyshev: ratio eigenvalue", 7.0);
826 paramList.set("coarse: params", coarseParamList);
827
828 paramList.set("max levels", 5);
829 paramList.set("aggregation: drop tol", 0.40);
830
831 }
832 using prec_t = MueLu::TpetraOperator<scalar_t, lno_t, gno_t, node_t>;
833 Teuchos::RCP<prec_t> prec = MueLu::CreateTpetraPreconditioner<
834 scalar_t, lno_t, gno_t, node_t>(laplacian_, paramList);
835
836 problem->setPrec(prec);
837
838#else
839 throw std::runtime_error("\nSphynx Error: MueLu requested but not compiled into Trilinos.\n");
840#endif
841 }
842
844 template <typename Adapter>
845 template <typename problem_t>
846 void Sphynx<Adapter>::setPolynomialPreconditioner(Teuchos::RCP<problem_t> &problem)
847 {
848 int verbosity2 = Belos::Errors;
849 if(verbosity_ == 1)
850 verbosity2 += Belos::Warnings;
851 else if(verbosity_ == 2)
852 verbosity2 += Belos::Warnings + Belos::FinalSummary;
853 else if(verbosity_ == 3)
854 verbosity2 += Belos::Warnings + Belos::FinalSummary + Belos::TimingDetails;
855 else if(verbosity_ >= 4)
856 verbosity2 += Belos::Warnings + Belos::FinalSummary + Belos::TimingDetails
857 + Belos::StatusTestDetails;
858
859 Teuchos::ParameterList paramList;
860 paramList.set("Polynomial Type", "Roots");
861 paramList.set("Orthogonalization","ICGS");
862 paramList.set("Maximum Degree", laplacian_->getGlobalNumRows() > 100 ? 25 : 5);
863 paramList.set("Polynomial Tolerance", 1.0e-6 );
864 paramList.set("Verbosity", verbosity2 );
865 paramList.set("Random RHS", false );
866 paramList.set("Outer Solver", "");
867 paramList.set("Timer Label", "Belos Polynomial Solve" );
868
869 // Construct a linear problem for the polynomial solver manager
870 using lproblem_t = Belos::LinearProblem<scalar_t, mvector_t, op_t>;
871 Teuchos::RCP<lproblem_t> innerPolyProblem(new lproblem_t());
872 innerPolyProblem->setOperator(laplacian_);
873
874 using btop_t = Belos::TpetraOperator<scalar_t, lno_t, gno_t, node_t>;
875 Teuchos::RCP<btop_t> polySolver(new btop_t(innerPolyProblem,
876 Teuchos::rcpFromRef(paramList),
877 "GmresPoly", true));
878 problem->setPrec(polySolver);
879 }
880
882 template <typename Adapter>
883 template <typename problem_t>
884 void Sphynx<Adapter>::setJacobiPreconditioner(Teuchos::RCP<problem_t> &problem)
885 {
886
887 Teuchos::RCP<Ifpack2::Preconditioner<scalar_t, lno_t, gno_t, node_t>> prec;
888 std::string precType = "RELAXATION";
889
890 prec = Ifpack2::Factory::create<matrix_t> (precType, laplacian_);
891
892 Teuchos::ParameterList precParams;
893 precParams.set("relaxation: type", "Jacobi");
894 precParams.set("relaxation: fix tiny diagonal entries", true);
895 precParams.set("relaxation: min diagonal value", 1.0e-8);
896
897 prec->setParameters(precParams);
898 prec->initialize();
899 prec->compute();
900
901 problem->setPrec(prec);
902
903 }
904
906 // Transform the computed eigenvectors into coordinates
907 template <typename Adapter>
908 void Sphynx<Adapter>::eigenvecsToCoords(Teuchos::RCP<mvector_t> &eigenVectors,
909 int computedNumEv,
910 Teuchos::RCP<mvector_t> &coordinates)
911 {
912 // Extract the meaningful eigenvectors by getting rid of the first one
913 Teuchos::Array<size_t> columns (computedNumEv-1);
914 for (int j = 0; j < computedNumEv-1; ++j) {
915 columns[j] = j+1;
916 }
917 coordinates = eigenVectors->subCopy (columns());
918 coordinates->setCopyOrView (Teuchos::View);
919 }
920
921
923 // If user provided some weights, use them by getting them from the adapter.
924 // If user didn't provide weights but told us to use degree as weight, do so.
925 // If user neither provided weights nor told us what to do, use degree as weight.
926 template <typename Adapter>
927 void Sphynx<Adapter>::computeWeights(std::vector<const weight_t *> vecweights,
928 std::vector<int> strides)
929 {
930
931 int numWeights = adapter_->getNumWeightsPerVertex();
932 int numConstraints = numWeights > 0 ? numWeights : 1;
933
934 size_t myNumVertices = adapter_->getLocalNumVertices();
935 weight_t ** weights = new weight_t*[numConstraints];
936 for(int j = 0; j < numConstraints; j++)
937 weights[j] = new weight_t[myNumVertices];
938
939 // Will be needed if we use degree as weight
940 const offset_t *offset;
941 const gno_t *columnId;
942
943 // If user hasn't set any weights, use vertex degrees as weights
944 if(numWeights == 0) {
945
946 // Compute the weight of vertex i as the number of nonzeros in row i.
947 adapter_->getEdgesView(offset, columnId);
948 for (size_t i = 0; i < myNumVertices; i++)
949 weights[0][i] = offset[i+1] - offset[i] - 1;
950
951 vecweights.push_back(weights[0]);
952 strides.push_back(1);
953 }
954 else {
955
956 // Use the weights if there are any already set in the adapter
957 for(int j = 0; j < numConstraints; j++) {
958
959 if(adapter_->useDegreeAsVertexWeight(j)) {
960 // Compute the weight of vertex i as the number of nonzeros in row i.
961 adapter_->getEdgesView(offset, columnId);
962 for (size_t i = 0; i < myNumVertices; i++)
963 weights[j][i] = offset[i+1] - offset[i];
964 }
965 else{
966 int stride;
967 const weight_t *wgt = NULL;
968 adapter_->getVertexWeightsView(wgt, stride, j);
969
970 for (size_t i = 0; i < myNumVertices; i++)
971 weights[j][i] = wgt[i];
972 }
973
974 vecweights.push_back(weights[j]);
975 strides.push_back(1);
976
977 }
978 }
979
980 }
981
982
984 // Compute a partition by calling MJ on given coordinates with given weights
985 template <typename Adapter>
986 void Sphynx<Adapter>::MJwrapper(const Teuchos::RCP<const mvector_t> &coordinates,
987 std::vector<const weight_t *> weights,
988 std::vector<int> strides,
989 const Teuchos::RCP<PartitioningSolution<Adapter>> &solution)
990 {
991
992 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer("Sphynx::MJ"));
993
994 using mvector_adapter_t = Zoltan2::XpetraMultiVectorAdapter<mvector_t>;
995 using base_adapter_t = typename mvector_adapter_t::base_adapter_t;
998
999 size_t myNumVertices = coordinates->getLocalLength();
1000
1001 // Create the base adapter for the multivector adapter
1002 Teuchos::RCP<mvector_adapter_t> adapcoordinates(new mvector_adapter_t(coordinates,
1003 weights,
1004 strides));
1005 Teuchos::RCP<const base_adapter_t> baseAdapter =
1006 Teuchos::rcp(dynamic_cast<const base_adapter_t *>(adapcoordinates.get()), false);
1007
1008 // Create the coordinate model using the base multivector adapter
1010
1011 // Create the MJ object
1012 Teuchos::RCP<const Comm<int>> comm2 = comm_;
1013 Teuchos::RCP<mj_t> mj(new mj_t(env_, comm2, baseAdapter));
1014
1015 // Partition with MJ
1016 Teuchos::RCP<solution_t> vectorsolution( new solution_t(env_, comm2, 1, mj));
1017 mj->partition(vectorsolution);
1018
1019 // Transform the solution
1020 Teuchos::ArrayRCP<part_t> parts(myNumVertices);
1021 for(size_t i = 0; i < myNumVertices; i++) parts[i] = vectorsolution->getPartListView()[i];
1022 solution->setParts(parts);
1023 }
1024
1025} // namespace Zoltan2
1026
1027#endif
Contains the Multi-jagged algorthm.
Defines the CoordinateModel classes.
Defines XpetraCrsGraphAdapter class.
Defines the XpetraMultiVectorAdapter.
Algorithm defines the base class for all algorithms.
A PartitioningSolution is a solution to a partitioning problem.
Sphynx(const RCP< const Environment > &env, const RCP< Teuchos::ParameterList > &params, const RCP< Teuchos::ParameterList > &sphynxParams, const RCP< const Comm< int > > &comm, const RCP< const XpetraCrsGraphAdapter< graph_t > > &adapter)
void setMueLuPreconditioner(Teuchos::RCP< problem_t > &problem)
typename Adapter::offset_t offset_t
void MJwrapper(const Teuchos::RCP< const mvector_t > &coordinates, std::vector< const weight_t * > weights, std::vector< int > strides, const Teuchos::RCP< PartitioningSolution< Adapter > > &solution)
Tpetra::CrsGraph< lno_t, gno_t, node_t > graph_t
void eigenvecsToCoords(Teuchos::RCP< mvector_t > &eigenVectors, int computedNumEv, Teuchos::RCP< mvector_t > &coordinates)
Tpetra::MultiVector< scalar_t, lno_t, gno_t, node_t > mvector_t
void setUserEigenvectors(const Teuchos::RCP< mvector_t > &userEvects)
typename Adapter::node_t node_t
Tpetra::Operator< scalar_t, lno_t, gno_t, node_t > op_t
typename Adapter::lno_t lno_t
int AnasaziWrapper(const int numEigenVectors)
Anasazi::MultiVecTraits< scalar_t, mvector_t > MVT
typename Adapter::scalar_t weight_t
void partition(const Teuchos::RCP< PartitioningSolution< Adapter > > &solution)
Teuchos::RCP< matrix_t > computeCombinatorialLaplacian()
Teuchos::RCP< matrix_t > computeNormalizedLaplacian(bool AHat=false)
typename Adapter::gno_t gno_t
void setPreconditioner(Teuchos::RCP< problem_t > &problem)
void setPolynomialPreconditioner(Teuchos::RCP< problem_t > &problem)
void setJacobiPreconditioner(Teuchos::RCP< problem_t > &problem)
Tpetra::CrsMatrix< scalar_t, lno_t, gno_t, node_t > matrix_t
typename Adapter::part_t part_t
void computeWeights(std::vector< const weight_t * > vecweights, std::vector< int > strides)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Multi Jagged coordinate partitioning algorithm.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Created by mbenlioglu on Aug 31, 2020.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
static RCP< tMVector_t > coordinates
static ArrayRCP< ArrayRCP< zscalar_t > > weights