10 #ifndef _ZOLTAN2_SPHYNXALGORITHM_HPP_
11 #define _ZOLTAN2_SPHYNXALGORITHM_HPP_
37 #include "Zoltan2Sphynx_config.h"
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"
53 #include "BelosLinearProblem.hpp"
54 #include "BelosTpetraOperator.hpp"
56 #include "Ifpack2_Factory.hpp"
58 #include "Teuchos_TimeMonitor.hpp"
60 #ifdef HAVE_ZOLTAN2SPHYNX_MUELU
61 #include "MueLu_CreateTpetraPreconditioner.hpp"
66 template <
typename Adapter>
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>;
86 typedef Anasazi::MultiVecTraits<scalar_t,mvector_t>
MVT;
93 Sphynx(
const RCP<const Environment> &env,
94 const RCP<Teuchos::ParameterList> ¶ms,
95 const RCP<Teuchos::ParameterList> &sphynxParams,
96 const RCP<
const Comm<int> > &comm,
98 env_(env), params_(params), sphynxParams_(sphynxParams), comm_(comm), adapter_(adapter)
102 const Teuchos::ParameterEntry *pe = params_->getEntryPtr(
"num_global_parts");
103 numGlobalParts_ = pe->getValue<
int>(&numGlobalParts_);
104 if(numGlobalParts_ > 1){
106 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer(
"Sphynx::Laplacian"));
109 pe = sphynxParams_->getEntryPtr(
"sphynx_verbosity");
111 verbosity_ = pe->getValue<
int>(&verbosity_);
114 pe = sphynxParams_->getEntryPtr(
"sphynx_skip_preprocessing");
116 skipPrep_ = pe->getValue<
bool>(&skipPrep_);
121 graph_ = adapter_->getUserGraph();
123 throw std::runtime_error(
"\nSphynx Error: Preprocessing has not been implemented yet.\n");
149 template<
typename problem_t>
152 template<
typename problem_t>
155 template<
typename problem_t>
158 template<
typename problem_t>
166 std::vector<int> strides);
169 std::vector<const weight_t *>
weights,
170 std::vector<int> strides,
186 auto rowOffsets = graph_->getLocalGraphDevice().row_map;
187 auto rowOffsets_h = Kokkos::create_mirror_view(rowOffsets);
188 Kokkos::deep_copy(rowOffsets_h, rowOffsets);
191 const size_t numGlobalEntries = graph_->getGlobalNumEntries();
192 const size_t numLocalRows = graph_->getLocalNumRows();
193 const size_t numGlobalRows = graph_->getGlobalNumRows();
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;
203 size_t globalMax = 0;
204 Teuchos::reduceAll<int, size_t> (*comm_, Teuchos::REDUCE_MAX, 1, &localMax, &globalMax);
206 double avg =
static_cast<double>(numGlobalEntries-numGlobalRows)/numGlobalRows;
207 double maxOverAvg =
static_cast<double>(globalMax)/avg;
211 if(maxOverAvg > 10) {
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;
242 precType_ =
"jacobi";
243 #ifdef HAVE_ZOLTAN2SPHYNX_MUELU
248 const Teuchos::ParameterEntry *pe = sphynxParams_->getEntryPtr(
"sphynx_preconditioner_type");
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");
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.");
263 if(precType_ ==
"polynomial")
270 pe = sphynxParams_->getEntryPtr(
"sphynx_problem_type");
272 std::string probType =
"";
273 probType = pe->getValue<std::string>(&probType);
275 if(probType ==
"combinatorial")
277 else if(probType ==
"generalized")
279 else if(probType ==
"normalized")
282 throw std::runtime_error(
"\nSphynx Error: " + probType +
" is not recognized as a problem type.\n"
283 +
" Possible values: combinatorial, generalized, and normalized.\n");
289 if(!irregular_ && (precType_ ==
"jacobi" || precType_ ==
"polynomial"))
294 pe = sphynxParams_->getEntryPtr(
"sphynx_tolerance");
296 tolerance_ = pe->getValue<
scalar_t>(&tolerance_);
305 pe = sphynxParams_->getEntryPtr(
"sphynx_initial_guess");
307 std::string initialGuess =
"";
308 initialGuess = pe->getValue<std::string>(&initialGuess);
310 if(initialGuess ==
"random")
312 else if(initialGuess ==
"constants")
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");
328 if(solverType_ ==
"randomized")
342 auto rowOffsets = graph_->getLocalGraphHost().row_map;
345 Teuchos::RCP<matrix_t> degMat(
new matrix_t (graph_->getRowMap(),
351 lno_t numRows =
static_cast<lno_t>(graph_->getLocalNumRows());
354 for (
lno_t i = 0; i < numRows; ++i) {
355 val[0] = rowOffsets(i+1) - rowOffsets(i) - 1;
357 degMat->insertLocalValues(i, 1, val, ind);
362 degMat->fillComplete(graph_->getDomainMap(), graph_->getRangeMap());
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>;
379 const size_t numEnt = graph_->getLocalNumEntries();
380 const size_t numRows = graph_->getLocalNumRows();
383 values_view_t newVal (Kokkos::view_alloc(
"CombLapl::val", Kokkos::WithoutInitializing), numEnt);
384 Kokkos::deep_copy(newVal, -1);
387 offset_view_t diagOffsets(Kokkos::view_alloc(
"Diag Offsets", Kokkos::WithoutInitializing), numRows);
388 graph_->getLocalDiagOffsets(diagOffsets);
391 auto rowOffsets = graph_->getLocalGraphDevice().row_map;
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;
402 Teuchos::RCP<matrix_t> laplacian (
new matrix_t(graph_, newVal));
403 laplacian->fillComplete (graph_->getDomainMap(), graph_->getRangeMap());
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 = Kokkos::ArithTraits<scalar_t>;
430 const size_t numEnt = graph_->getLocalNumEntries();
431 const size_t numRows = graph_->getLocalNumRows();
434 values_view_t newVal (Kokkos::view_alloc(
"NormLapl::val", Kokkos::WithoutInitializing), numEnt);
436 Kokkos::deep_copy(newVal, 1);
439 Kokkos::deep_copy(newVal, -1);
443 dual_view_t dv (
"MV::DualView", numRows, 1);
444 auto deginvsqrt = dv.d_view;
447 offset_view_t diagOffsets(Kokkos::view_alloc(
"Diag Offsets", Kokkos::WithoutInitializing), numRows);
448 graph_->getLocalDiagOffsets(diagOffsets);
451 auto rowOffsets = graph_->getLocalGraphDevice().row_map;
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);
464 Teuchos::RCP<matrix_t> laplacian (
new matrix_t(graph_, newVal));
465 laplacian->fillComplete (graph_->getDomainMap(), graph_->getRangeMap());
468 vector_t degInvSqrt(graph_->getRowMap(), dv);
471 laplacian->leftScale(degInvSqrt);
472 laplacian->rightScale(degInvSqrt);
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_;
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_;
498 std::string precType_;
499 std::string solverType_;
516 template <
typename Adapter>
519 eigenVectors_ = userEvects;
528 template <
typename Adapter>
532 if(numGlobalParts_ == 1) {
534 size_t numRows =adapter_->getUserGraph()->getLocalNumRows();
535 Teuchos::ArrayRCP<part_t> parts(numRows,0);
536 solution->setParts(parts);
543 int numEigenVectors = (int) log2(numGlobalParts_)+1;
546 if(eigenVectors_ == Teuchos::null){
551 if(computedNumEv <= 1 && (solverType_ ==
"LOBPCG" || solverType_ ==
"GeneralizedDavidson" || solverType_ ==
"BlockDavidson" || solverType_ ==
"BlockKrylovSchur")) {
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");
560 computedNumEv = (int) eigenVectors_->getNumVectors();
561 if(computedNumEv <= numEigenVectors) {
563 std::runtime_error(
"\nSphynx Error: Number of eigenvectors given by user\n"
564 " is less than number of Eigenvectors needed for partition." );
574 std::vector<const weight_t *>
weights;
575 std::vector<int> wstrides;
587 template <
typename Adapter>
591 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer(
"Sphynx::Anasazi"));
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;
605 bool useFullOrtho = sphynxParams_->get(
"sphynx_use_full_ortho",
true);
612 int anasaziVerbosity = Anasazi::Errors + Anasazi::Warnings;
614 anasaziVerbosity += Anasazi::FinalSummary + Anasazi::TimingDetails;
616 anasaziVerbosity += Anasazi::IterationDetails;
618 anasaziVerbosity += Anasazi::StatusTestDetails + Anasazi::OrthoDetails
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);
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 );
643 auto map = laplacian_->getRangeMap();
644 Teuchos::RCP<mvector_t> ivec(
new mvector_t(map, numEigenVectors));
648 Anasazi::MultiVecTraits<scalar_t, mvector_t>::MvRandom(*ivec);
649 ivec->getVectorNonConst(0)->putScalar(1.);
655 ivec->getVectorNonConst(0)->putScalar(1.);
656 for (
int j = 1; j < numEigenVectors; j++)
657 ivec->getVectorNonConst(j)->putScalar(0.);
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\".");
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.);
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);
677 if(solverType_ !=
"randomized"){
684 bool boolret = problem->setProblem();
685 if (boolret !=
true) {
686 throw std::runtime_error(
"\nAnasazi::BasicEigenproblem::setProblem() returned with error.\n");
689 Teuchos::RCP<Anasazi::SolverManager<scalar_t, mvector_t, op_t>> solver;
691 if(solverType_ ==
"LOBPCG"){
692 solver = Teuchos::rcp(
new Anasazi::LOBPCGSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
694 else if(solverType_ ==
"BlockDavidson"){
695 solver = Teuchos::rcp(
new Anasazi::BlockDavidsonSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
697 else if(solverType_ ==
"GeneralizedDavidson"){
698 solver = Teuchos::rcp(
new Anasazi::GeneralizedDavidsonSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
700 else if(solverType_ ==
"BlockKrylovSchur"){
701 solver = Teuchos::rcp(
new Anasazi::BlockKrylovSchurSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
704 solver = Teuchos::rcp(
new Anasazi::Experimental::RandomizedSolMgr<scalar_t, mvector_t, op_t>(problem, anasaziParams));
708 if (verbosity_ > 0 && comm_->getRank() == 0)
709 anasaziParams.print(std::cout);
712 Anasazi::ReturnType returnCode = solver->solve();
715 iter = solver->getNumIters();
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;
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;
732 std::cout <<
"No of comp. vecs. : " << numev << std::endl;
735 std::cout <<
"Solver type: " << solverType_ << std::endl;
737 if(solverType_ ==
"randomized") {
738 std::vector<double> normR(numev);
742 Teuchos::SerialDenseMatrix<int,double> T (numev, numev);
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);
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"
757 <<
"------------------------------------------------------" << std::endl;
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
764 std::cout <<
"------------------------------------------------------" << std::endl;
773 template <
typename Adapter>
774 template <
typename problem_t>
777 if(comm_->getRank() == 0) std::cout <<
"precType_ is: " << precType_ << std::endl;
779 if(precType_ ==
"muelu") {
782 else if(precType_ ==
"polynomial") {
785 else if(precType_ ==
"jacobi") {
792 template <
typename Adapter>
793 template <
typename problem_t>
796 #ifdef HAVE_ZOLTAN2SPHYNX_MUELU
797 Teuchos::ParameterList paramList;
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");
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);
820 paramList.set(
"multigrid algorithm",
"unsmoothed");
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);
828 paramList.set(
"max levels", 5);
829 paramList.set(
"aggregation: drop tol", 0.40);
832 using prec_t = MueLu::TpetraOperator<scalar_t, lno_t, gno_t, node_t>;
833 Teuchos::RCP<prec_t> prec = MueLu::CreateTpetraPreconditioner<
836 problem->setPrec(prec);
839 throw std::runtime_error(
"\nSphynx Error: MueLu requested but not compiled into Trilinos.\n");
844 template <
typename Adapter>
845 template <
typename problem_t>
848 int verbosity2 = Belos::Errors;
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;
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" );
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_);
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),
878 problem->setPrec(polySolver);
882 template <
typename Adapter>
883 template <
typename problem_t>
887 Teuchos::RCP<Ifpack2::Preconditioner<scalar_t, lno_t, gno_t, node_t>> prec;
888 std::string precType =
"RELAXATION";
890 prec = Ifpack2::Factory::create<matrix_t> (precType, laplacian_);
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);
897 prec->setParameters(precParams);
901 problem->setPrec(prec);
907 template <
typename Adapter>
913 Teuchos::Array<size_t> columns (computedNumEv-1);
914 for (
int j = 0; j < computedNumEv-1; ++j) {
917 coordinates = eigenVectors->subCopy (columns());
918 coordinates->setCopyOrView (Teuchos::View);
926 template <
typename Adapter>
928 std::vector<int> strides)
931 int numWeights = adapter_->getNumWeightsPerVertex();
932 int numConstraints = numWeights > 0 ? numWeights : 1;
934 size_t myNumVertices = adapter_->getLocalNumVertices();
936 for(
int j = 0; j < numConstraints; j++)
937 weights[j] =
new weight_t[myNumVertices];
941 const gno_t *columnId;
944 if(numWeights == 0) {
947 adapter_->getEdgesView(offset, columnId);
948 for (
size_t i = 0; i < myNumVertices; i++)
949 weights[0][i] = offset[i+1] - offset[i] - 1;
951 vecweights.push_back(weights[0]);
952 strides.push_back(1);
957 for(
int j = 0; j < numConstraints; j++) {
959 if(adapter_->useDegreeAsVertexWeight(j)) {
961 adapter_->getEdgesView(offset, columnId);
962 for (
size_t i = 0; i < myNumVertices; i++)
963 weights[j][i] = offset[i+1] - offset[i];
968 adapter_->getVertexWeightsView(wgt, stride, j);
970 for (
size_t i = 0; i < myNumVertices; i++)
971 weights[j][i] = wgt[i];
974 vecweights.push_back(weights[j]);
975 strides.push_back(1);
985 template <
typename Adapter>
987 std::vector<const weight_t *>
weights,
988 std::vector<int> strides,
992 Teuchos::TimeMonitor t(*Teuchos::TimeMonitor::getNewTimer(
"Sphynx::MJ"));
999 size_t myNumVertices = coordinates->getLocalLength();
1002 Teuchos::RCP<mvector_adapter_t> adapcoordinates(
new mvector_adapter_t(coordinates,
1005 Teuchos::RCP<const base_adapter_t> baseAdapter =
1006 Teuchos::rcp(dynamic_cast<const base_adapter_t *>(adapcoordinates.get()),
false);
1012 Teuchos::RCP<const Comm<int>> comm2 = comm_;
1013 Teuchos::RCP<mj_t> mj(
new mj_t(env_, comm2, baseAdapter));
1016 Teuchos::RCP<solution_t> vectorsolution(
new solution_t(env_, comm2, 1, mj));
1017 mj->partition(vectorsolution);
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);
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
void partition(const Teuchos::RCP< PartitioningSolution< Adapter > > &solution)
Sphynx(const RCP< const Environment > &env, const RCP< Teuchos::ParameterList > ¶ms, const RCP< Teuchos::ParameterList > &sphynxParams, const RCP< const Comm< int > > &comm, const RCP< const XpetraCrsGraphAdapter< graph_t > > &adapter)
Anasazi::MultiVecTraits< scalar_t, mvector_t > MVT
void computeWeights(std::vector< const weight_t * > vecweights, std::vector< int > strides)
void setMueLuPreconditioner(Teuchos::RCP< problem_t > &problem)
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
void eigenvecsToCoords(Teuchos::RCP< mvector_t > &eigenVectors, int computedNumEv, Teuchos::RCP< mvector_t > &coordinates)
map_t::global_ordinal_type gno_t
void computeDegreeMatrix()
typename Adapter::node_t node_t
void setPolynomialPreconditioner(Teuchos::RCP< problem_t > &problem)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Tpetra::MultiVector< scalar_t, lno_t, gno_t, node_t > mvector_t
void setJacobiPreconditioner(Teuchos::RCP< problem_t > &problem)
Teuchos::RCP< matrix_t > computeCombinatorialLaplacian()
void setUserEigenvectors(const Teuchos::RCP< mvector_t > &userEvects)
Defines the XpetraMultiVectorAdapter.
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
Contains the Multi-jagged algorthm.
Adapter::scalar_t scalar_t
A PartitioningSolution is a solution to a partitioning problem.
Tpetra::CrsGraph< lno_t, gno_t, node_t > graph_t
void determineRegularity()
Tpetra::CrsMatrix< scalar_t, lno_t, gno_t, node_t > matrix_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)
Algorithm defines the base class for all algorithms.
map_t::local_ordinal_type lno_t
An adapter for Xpetra::MultiVector.
int AnasaziWrapper(const int numEigenVectors)
Tpetra::Operator< scalar_t, lno_t, gno_t, node_t > op_t
Defines the CoordinateModel classes.
typename Adapter::offset_t offset_t
void setPreconditioner(Teuchos::RCP< problem_t > &problem)
typename Adapter::scalar_t weight_t
Teuchos::RCP< matrix_t > computeNormalizedLaplacian(bool AHat=false)
Multi Jagged coordinate partitioning algorithm.