10#ifndef MUELU_SHIFTEDLAPLACIAN_DEF_HPP
11#define MUELU_SHIFTEDLAPLACIAN_DEF_HPP
15#include <MueLu_AmalgamationFactory.hpp>
16#include <MueLu_CoalesceDropFactory.hpp>
17#include <MueLu_CoarseMapFactory.hpp>
18#include <MueLu_CoupledRBMFactory.hpp>
19#include <MueLu_DirectSolver.hpp>
20#include <MueLu_GenericRFactory.hpp>
21#include <MueLu_Hierarchy.hpp>
22#include <MueLu_Ifpack2Smoother.hpp>
23#include <MueLu_PFactory.hpp>
24#include <MueLu_PgPFactory.hpp>
25#include <MueLu_RAPFactory.hpp>
26#include <MueLu_RAPShiftFactory.hpp>
27#include <MueLu_SaPFactory.hpp>
28#include <MueLu_ShiftedLaplacian.hpp>
29#include <MueLu_ShiftedLaplacianOperator.hpp>
30#include <MueLu_SmootherFactory.hpp>
31#include <MueLu_SmootherPrototype.hpp>
32#include <MueLu_TentativePFactory.hpp>
33#include <MueLu_TransPFactory.hpp>
34#include <MueLu_UncoupledAggregationFactory.hpp>
35#include <MueLu_Utilities.hpp>
40template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
44template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
47 coarseGridSize_ = paramList->get(
"MueLu: coarse size", 1000);
48 numLevels_ = paramList->get(
"MueLu: levels", 3);
49 int stype = paramList->get(
"MueLu: smoother", 8);
52 }
else if (stype == 2) {
53 Smoother_ =
"gauss-seidel";
54 }
else if (stype == 3) {
55 Smoother_ =
"symmetric gauss-seidel";
56 }
else if (stype == 4) {
57 Smoother_ =
"chebyshev";
58 }
else if (stype == 5) {
60 }
else if (stype == 6) {
62 }
else if (stype == 7) {
64 }
else if (stype == 8) {
65 Smoother_ =
"schwarz";
66 }
else if (stype == 9) {
67 Smoother_ =
"superilu";
68 }
else if (stype == 10) {
69 Smoother_ =
"superlu";
71 Smoother_ =
"schwarz";
73 smoother_sweeps_ = paramList->get(
"MueLu: sweeps", 5);
74 smoother_damping_ = paramList->get(
"MueLu: relax val", 1.0);
75 ncycles_ = paramList->get(
"MueLu: cycles", 1);
76 iters_ = paramList->get(
"MueLu: iterations", 500);
77 solverType_ = paramList->get(
"MueLu: solver type", 1);
78 restart_size_ = paramList->get(
"MueLu: restart size", 100);
79 recycle_size_ = paramList->get(
"MueLu: recycle size", 25);
80 isSymmetric_ = paramList->get(
"MueLu: symmetric",
true);
81 ilu_leveloffill_ = paramList->get(
"MueLu: level-of-fill", 5);
82 ilu_abs_thresh_ = paramList->get(
"MueLu: abs thresh", 0.0);
83 ilu_rel_thresh_ = paramList->get(
"MueLu: rel thresh", 1.0);
84 ilu_diagpivotthresh_ = paramList->get(
"MueLu: piv thresh", 0.1);
85 ilu_drop_tol_ = paramList->get(
"MueLu: drop tol", 0.01);
86 ilu_fill_tol_ = paramList->get(
"MueLu: fill tol", 0.01);
87 schwarz_overlap_ = paramList->get(
"MueLu: overlap", 0);
88 schwarz_usereorder_ = paramList->get(
"MueLu: use reorder",
true);
89 int combinemode = paramList->get(
"MueLu: combine mode", 1);
90 if (combinemode == 0) {
91 schwarz_combinemode_ = Tpetra::ZERO;
93 schwarz_combinemode_ = Tpetra::ADD;
95 tol_ = paramList->get(
"MueLu: tolerance", 0.001);
98template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
101 if (A_ != Teuchos::null)
102 TpetraA_ = toTpetra(A_);
103 if (LinearProblem_ != Teuchos::null)
104 LinearProblem_->setOperator(TpetraA_);
107template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 if (LinearProblem_ != Teuchos::null)
111 LinearProblem_->setOperator(TpetraA_);
114template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
117 GridTransfersExist_ =
false;
120template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
122 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Atmp = rcp(
new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(TpetraP));
123 P_ = rcp(
new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Atmp));
124 GridTransfersExist_ =
false;
127template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
132template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Atmp = rcp(
new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(TpetraK));
135 K_ = rcp(
new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Atmp));
138template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
143template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
145 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Atmp = rcp(
new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(TpetraM));
146 M_ = rcp(
new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Atmp));
149template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
154template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
156 NullSpace_ = NullSpace;
159template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
161 levelshifts_ = levelshifts;
162 numLevels_ = levelshifts_.size();
166template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
180 Manager_->SetFactory(
"UnAmalgamationInfo", Amalgfact_);
181 Teuchos::ParameterList params;
182 params.set(
"lightweight wrap",
true);
183 params.set(
"aggregation: drop scheme",
"classical");
184 Dropfact_->SetParameterList(params);
185 Manager_->SetFactory(
"Graph", Dropfact_);
186 Manager_->SetFactory(
"Aggregates", UCaggfact_);
187 Manager_->SetFactory(
"CoarseMap", CoarseMapfact_);
188 Manager_->SetFactory(
"Ptent", TentPfact_);
189 if (isSymmetric_ ==
true) {
190 Manager_->SetFactory(
"P", Pfact_);
191 Manager_->SetFactory(
"R", TransPfact_);
193 Manager_->SetFactory(
"P", PgPfact_);
194 Manager_->SetFactory(
"R", Rfact_);
199 if (Smoother_ ==
"jacobi") {
200 precType_ =
"RELAXATION";
201 precList_.set(
"relaxation: type",
"Jacobi");
202 precList_.set(
"relaxation: sweeps", smoother_sweeps_);
203 precList_.set(
"relaxation: damping factor", smoother_damping_);
204 }
else if (Smoother_ ==
"gauss-seidel") {
205 precType_ =
"RELAXATION";
206 precList_.set(
"relaxation: type",
"Gauss-Seidel");
207 precList_.set(
"relaxation: sweeps", smoother_sweeps_);
208 precList_.set(
"relaxation: damping factor", smoother_damping_);
209 }
else if (Smoother_ ==
"symmetric gauss-seidel") {
210 precType_ =
"RELAXATION";
211 precList_.set(
"relaxation: type",
"Symmetric Gauss-Seidel");
212 precList_.set(
"relaxation: sweeps", smoother_sweeps_);
213 precList_.set(
"relaxation: damping factor", smoother_damping_);
214 }
else if (Smoother_ ==
"chebyshev") {
215 precType_ =
"CHEBYSHEV";
216 }
else if (Smoother_ ==
"krylov") {
217 precType_ =
"KRYLOV";
218 precList_.set(
"krylov: iteration type", krylov_type_);
219 precList_.set(
"krylov: number of iterations", krylov_iterations_);
220 precList_.set(
"krylov: residual tolerance", 1.0e-8);
221 precList_.set(
"krylov: block size", 1);
222 precList_.set(
"krylov: preconditioner type", krylov_preconditioner_);
223 precList_.set(
"relaxation: sweeps", 1);
225 }
else if (Smoother_ ==
"ilut") {
227 precList_.set(
"fact: ilut level-of-fill", ilu_leveloffill_);
228 precList_.set(
"fact: absolute threshold", ilu_abs_thresh_);
229 precList_.set(
"fact: relative threshold", ilu_rel_thresh_);
230 precList_.set(
"fact: drop tolerance", ilu_drop_tol_);
231 precList_.set(
"fact: relax value", ilu_relax_val_);
232 }
else if (Smoother_ ==
"riluk") {
234 precList_.set(
"fact: iluk level-of-fill", ilu_leveloffill_);
235 precList_.set(
"fact: absolute threshold", ilu_abs_thresh_);
236 precList_.set(
"fact: relative threshold", ilu_rel_thresh_);
237 precList_.set(
"fact: drop tolerance", ilu_drop_tol_);
238 precList_.set(
"fact: relax value", ilu_relax_val_);
239 }
else if (Smoother_ ==
"schwarz") {
240 precType_ =
"SCHWARZ";
241 precList_.set(
"schwarz: overlap level", schwarz_overlap_);
242 precList_.set(
"schwarz: combine mode", schwarz_combinemode_);
243 precList_.set(
"schwarz: use reordering", schwarz_usereorder_);
245 precList_.set(
"order_method", schwarz_ordermethod_);
246 precList_.sublist(
"schwarz: reordering list").set(
"order_method", schwarz_ordermethod_);
247 precList_.sublist(
"schwarz: subdomain solver parameters").set(
"fact: ilut level-of-fill", ilu_leveloffill_);
248 precList_.sublist(
"schwarz: subdomain solver parameters").set(
"fact: absolute threshold", ilu_abs_thresh_);
249 precList_.sublist(
"schwarz: subdomain solver parameters").set(
"fact: relative threshold", ilu_rel_thresh_);
250 precList_.sublist(
"schwarz: subdomain solver parameters").set(
"fact: drop tolerance", ilu_drop_tol_);
251 precList_.sublist(
"schwarz: subdomain solver parameters").set(
"fact: relax value", ilu_relax_val_);
252 }
else if (Smoother_ ==
"superilu") {
253 precType_ =
"superlu";
254 precList_.set(
"RowPerm", ilu_rowperm_);
255 precList_.set(
"ColPerm", ilu_colperm_);
256 precList_.set(
"DiagPivotThresh", ilu_diagpivotthresh_);
257 precList_.set(
"ILU_DropRule", ilu_drop_rule_);
258 precList_.set(
"ILU_DropTol", ilu_drop_tol_);
259 precList_.set(
"ILU_FillFactor", ilu_leveloffill_);
260 precList_.set(
"ILU_Norm", ilu_normtype_);
261 precList_.set(
"ILU_MILU", ilu_milutype_);
262 precList_.set(
"ILU_FillTol", ilu_fill_tol_);
263 precList_.set(
"ILU_Flag",
true);
264 }
else if (Smoother_ ==
"superlu") {
265 precType_ =
"superlu";
266 precList_.set(
"ColPerm", ilu_colperm_);
267 precList_.set(
"DiagPivotThresh", ilu_diagpivotthresh_);
272#if defined(HAVE_AMESOS2_SUPERLU)
273 coarsestSmooProto_ = rcp(
new DirectSolver(
"Superlu", coarsestSmooList_));
274#elif defined(HAVE_AMESOS2_KLU2)
275 coarsestSmooProto_ = rcp(
new DirectSolver(
"Klu", coarsestSmooList_));
276#elif defined(HAVE_AMESOS2_SUPERLUDIST)
277 coarsestSmooProto_ = rcp(
new DirectSolver(
"Superludist", coarsestSmooList_));
281 coarsestSmooFact_ = rcp(
new SmootherFactory(coarsestSmooProto_, Teuchos::null));
289 if (K_ != Teuchos::null) {
290 Manager_->SetFactory(
"Smoother", Teuchos::null);
291 Manager_->SetFactory(
"CoarseSolver", Teuchos::null);
293 if (NullSpace_ != Teuchos::null)
294 Hierarchy_->GetLevel(0)->Set(
"Nullspace", NullSpace_);
295 if (isSymmetric_ ==
true) {
296 Hierarchy_->Keep(
"P", Pfact_.get());
297 Hierarchy_->Keep(
"R", TransPfact_.get());
298 Hierarchy_->SetImplicitTranspose(
true);
300 Hierarchy_->Keep(
"P", PgPfact_.get());
301 Hierarchy_->Keep(
"R", Rfact_.get());
303 Hierarchy_->Keep(
"Ptent", TentPfact_.get());
304 Hierarchy_->SetMaxCoarseSize(coarseGridSize_);
305 Hierarchy_->Setup(*Manager_, 0, numLevels_);
306 GridTransfersExist_ =
true;
310 Manager_->SetFactory(
"Smoother", smooFact_);
311 Manager_->SetFactory(
"CoarseSolver", coarsestSmooFact_);
313 if (NullSpace_ != Teuchos::null)
314 Hierarchy_->GetLevel(0)->Set(
"Nullspace", NullSpace_);
315 if (isSymmetric_ ==
true)
316 Hierarchy_->SetImplicitTranspose(
true);
317 Hierarchy_->SetMaxCoarseSize(coarseGridSize_);
318 Hierarchy_->Setup(*Manager_, 0, numLevels_);
319 GridTransfersExist_ =
true;
323 BelosList_ = rcp(
new Teuchos::ParameterList(
"GMRES"));
324 BelosList_->set(
"Maximum Iterations", iters_);
325 BelosList_->set(
"Convergence Tolerance", tol_);
326 BelosList_->set(
"Verbosity", Belos::Errors + Belos::Warnings + Belos::StatusTestDetails);
327 BelosList_->set(
"Output Frequency", 1);
328 BelosList_->set(
"Output Style", Belos::Brief);
329 BelosList_->set(
"Num Blocks", restart_size_);
330 BelosList_->set(
"Num Recycled Blocks", recycle_size_);
334template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
336 int numLevels = Hierarchy_->GetNumLevels();
337 Manager_->SetFactory(
"Smoother", smooFact_);
338 Manager_->SetFactory(
"CoarseSolver", coarsestSmooFact_);
339 Hierarchy_->GetLevel(0)->Set(
"A", P_);
340 Hierarchy_->Setup(*Manager_, 0, numLevels);
345template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
347 int numLevels = Hierarchy_->GetNumLevels();
348 Acshift_->SetShifts(levelshifts_);
349 Manager_->SetFactory(
"Smoother", smooFact_);
350 Manager_->SetFactory(
"CoarseSolver", coarsestSmooFact_);
351 Manager_->SetFactory(
"A", Acshift_);
352 Manager_->SetFactory(
"K", Acshift_);
353 Manager_->SetFactory(
"M", Acshift_);
354 Hierarchy_->GetLevel(0)->Set(
"A", P_);
355 Hierarchy_->GetLevel(0)->Set(
"K", K_);
356 Hierarchy_->GetLevel(0)->Set(
"M", M_);
357 Hierarchy_->Setup(*Manager_, 0, numLevels);
362template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
365 if (GridTransfersExist_ ==
false) {
367 if (NullSpace_ != Teuchos::null)
368 Hierarchy_->GetLevel(0)->Set(
"Nullspace", NullSpace_);
369 if (isSymmetric_ ==
true)
370 Hierarchy_->SetImplicitTranspose(
true);
371 Hierarchy_->SetMaxCoarseSize(coarseGridSize_);
372 Hierarchy_->Setup(*Manager_, 0, numLevels_);
373 GridTransfersExist_ =
true;
378template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
383 if (LinearProblem_ == Teuchos::null)
385 LinearProblem_->setOperator(TpetraA_);
386 LinearProblem_->setRightPrec(MueLuOp_);
387 if (SolverManager_ == Teuchos::null) {
388 std::string solverName;
390 if (solverType_ == 1) {
391 solverName =
"Block GMRES";
392 }
else if (solverType_ == 2) {
393 solverName =
"Recycling GMRES";
395 solverName =
"Flexible GMRES";
397 SolverManager_ = SolverFactory_->create(solverName, BelosList_);
398 SolverManager_->setProblem(LinearProblem_);
402template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
404 LinearProblem_->setOperator(TpetraA_);
408template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
411 LinearProblem_->setProblem(X, B);
413 return SolverManager_->solve();
417template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
419 RCP<MultiVector>& X) {
421 Hierarchy_->Iterate(*B, *X, 1,
true, 0);
425template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
427 RCP<Tpetra::MultiVector<SC, LO, GO, NO> >& X) {
428 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > XpetraX = Teuchos::rcp(
new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(X));
429 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > XpetraB = Teuchos::rcp(
new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(B));
431 Hierarchy_->Iterate(*XpetraB, *XpetraX, 1,
true, 0);
435template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
437 int numiters = SolverManager_->getNumIters();
442template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
443typename Teuchos::ScalarTraits<Scalar>::magnitudeType
445 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType MT;
446 MT residual = SolverManager_->achievedTol();
452#define MUELU_SHIFTEDLAPLACIAN_SHORT
AmalgamationFactory for subblocks of strided map based amalgamation data.
Factory for creating a graph based on a given matrix.
Factory for generating coarse level map. Used by TentativePFactory.
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
This class specifies the default factory that should generate some data on a Level if the data does n...
Factory for building restriction operators using a prolongator factory.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Class that encapsulates Ifpack2 smoothers.
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
Factory for building coarse matrices.
Factory for building coarse grid matrices, when the matrix is of the form K+a*M. Useful when you want...
Factory for building Smoothed Aggregation prolongators.
void setPreconditioningMatrix(RCP< Matrix > &P)
void setmass(RCP< Matrix > &M)
void resetLinearProblem()
void setNullSpace(RCP< MultiVector > NullSpace)
void setcoords(RCP< MultiVector > &Coords)
Belos::LinearProblem< SC, TMV, OP > LinearProblem
Teuchos::ScalarTraits< Scalar >::magnitudeType GetResidual()
void setLevelShifts(std::vector< Scalar > levelshifts)
void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList)
virtual ~ShiftedLaplacian()
void setstiff(RCP< Matrix > &K)
void setProblemMatrix(RCP< Matrix > &A)
void multigrid_apply(const RCP< MultiVector > B, RCP< MultiVector > &X)
Belos::SolverFactory< SC, TMV, OP > SolverFactory
Belos::ReturnType solve(const RCP< TMV > B, RCP< TMV > &X)
Wraps an existing MueLu::Hierarchy as a Tpetra::Operator, with an optional two-level correction....
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
Factory for building tentative prolongator.
Factory for building restriction operators.
Factory for building uncoupled aggregates.
Namespace for MueLu classes and methods.