MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_ShiftedLaplacian_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_SHIFTEDLAPLACIAN_DECL_HPP
11#define MUELU_SHIFTEDLAPLACIAN_DECL_HPP
12
13// Xpetra
14#include <Xpetra_Matrix_fwd.hpp>
15#include <Xpetra_VectorFactory_fwd.hpp>
16#include <Xpetra_MultiVectorFactory_fwd.hpp>
17#include <Xpetra_TpetraMultiVector.hpp>
18
19// MueLu
20#include "MueLu.hpp"
21#include "MueLu_ConfigDefs.hpp"
22
23#if defined(HAVE_MUELU_IFPACK2)
24
25#include <MueLu_BaseClass.hpp>
40#include <MueLu_ShiftedLaplacianOperator.hpp>
47
48// Belos
49#ifdef HAVE_MUELU_TPETRA_INST_INT_INT
50#include <BelosConfigDefs.hpp>
51#include <BelosLinearProblem.hpp>
52#include <BelosSolverFactory.hpp>
53#include <BelosTpetraAdapter.hpp>
54#endif
55
56#include "Kokkos_Core.hpp"
57
58namespace MueLu {
59
69template <class Scalar = DefaultScalar,
72 class Node = DefaultNode>
74#undef MUELU_SHIFTEDLAPLACIAN_SHORT
76
77 typedef Tpetra::Vector<SC, LO, GO, NO> TVEC;
78 typedef Tpetra::MultiVector<SC, LO, GO, NO> TMV;
79 typedef Tpetra::Operator<SC, LO, GO, NO> OP;
80#ifdef HAVE_MUELU_TPETRA_INST_INT_INT
81 typedef Belos::LinearProblem<SC, TMV, OP> LinearProblem;
82 typedef Belos::SolverManager<SC, TMV, OP> SolverManager;
83 typedef Belos::SolverFactory<SC, TMV, OP> SolverFactory;
84#endif
85
86 public:
87 /*
88 FIXME 26-June-2015 JJH: This contructor is setting numerous defaults. However, they don't match the defaults
89 FIXME int the method setParameters(). There also isn't any parameter validation that I can see.
90 */
91
94 : numPDEs_(1)
95 , Smoother_("schwarz")
96 , Aggregation_("uncoupled")
97 , Nullspace_("constant")
98 , numLevels_(5)
99 , coarseGridSize_(100)
100 , omega_(2.0 * Kokkos::numbers::pi_v<double>)
101 , iters_(500)
102 , blksize_(1)
103 , tol_(1.0e-4)
104 , nsweeps_(5)
105 , ncycles_(1)
106 , cycles_(8)
107 , subiters_(10)
108 , option_(1)
109 , nproblems_(0)
110 , solverType_(1)
111 , restart_size_(100)
112 , recycle_size_(25)
114 , smoother_damping_((SC)1.0)
115 , krylov_type_(1)
118 , ilu_leveloffill_(5.0)
119 , ilu_abs_thresh_(0.0)
120 , ilu_rel_thresh_(1.0)
122 , ilu_drop_tol_(0.01)
123 , ilu_fill_tol_(0.01)
124 , ilu_relax_val_(1.0)
125 , ilu_rowperm_("LargeDiag")
126 , ilu_colperm_("COLAMD")
127 , ilu_drop_rule_("DROP_BASIC")
128 , ilu_normtype_("INF_NORM")
129 , ilu_milutype_("SILU")
131 , schwarz_usereorder_(true)
132 , schwarz_combinemode_(Tpetra::ADD)
133 , schwarz_ordermethod_("rcm")
134 , GridTransfersExist_(false)
135 , isSymmetric_(true) {}
136
137 // Destructor
138 virtual ~ShiftedLaplacian();
139
140 // Parameters
141 void setParameters(Teuchos::RCP<Teuchos::ParameterList> paramList);
142
143 // Set matrices
144 void setProblemMatrix(RCP<Matrix>& A);
145 void setProblemMatrix(RCP<Tpetra::CrsMatrix<SC, LO, GO, NO> >& TpetraA);
146 void setPreconditioningMatrix(RCP<Matrix>& P);
147 void setPreconditioningMatrix(RCP<Tpetra::CrsMatrix<SC, LO, GO, NO> >& TpetraP);
148 void setstiff(RCP<Matrix>& K);
149 void setstiff(RCP<Tpetra::CrsMatrix<SC, LO, GO, NO> >& TpetraK);
150 void setmass(RCP<Matrix>& M);
151 void setmass(RCP<Tpetra::CrsMatrix<SC, LO, GO, NO> >& TpetraM);
152 void setcoords(RCP<MultiVector>& Coords);
153 void setNullSpace(RCP<MultiVector> NullSpace);
154 void setLevelShifts(std::vector<Scalar> levelshifts);
155
156 // initialize: set parameters and factories, construct
157 // prolongation and restriction matrices
158 void initialize();
159 // setupFastRAP: setup hierarchy with
160 // prolongators of the stiffness matrix
161 // constant complex shifts
162 void setupFastRAP();
163 // setupSlowRAP: setup hierarchy with
164 // prolongators of the stiffness matrix
165 // variable complex shifts
166 void setupSlowRAP();
167 // setupNormalRAP: setup hierarchy with
168 // prolongators of the preconditioning matrix
169 void setupNormalRAP();
170 // setupSolver: initialize Belos solver
171 void setupSolver();
172 // resetLinearProblem: for multiple frequencies;
173 // reset the Belos operator if the frequency changes
174 void resetLinearProblem();
175
176 // Solve phase
177 int solve(const RCP<TMV> B, RCP<TMV>& X);
178 void multigrid_apply(const RCP<MultiVector> B,
179 RCP<MultiVector>& X);
180 void multigrid_apply(const RCP<Tpetra::MultiVector<SC, LO, GO, NO> > B,
181 RCP<Tpetra::MultiVector<SC, LO, GO, NO> >& X);
182 int GetIterations();
183 typename Teuchos::ScalarTraits<Scalar>::magnitudeType GetResidual();
184
185 RCP<FactoryManager> Manager_;
186
187 private:
188 // Problem options
189 // numPDEs_ -> number of DOFs at each node
191
192 // Multigrid options
193 // numLevels_ -> number of Multigrid levels
194 // coarseGridSize_ -> size of coarsest grid (if current level has less DOFs, stop coarsening)
197
198 // Shifted Laplacian/Helmholtz parameters
199 double omega_;
200 std::vector<SC> levelshifts_;
201
202 // Krylov solver inputs
203 // iters -> max number of iterations
204 // tol -> residual tolerance
206 double tol_;
210
211 // Smoother parameters
222 Tpetra::CombineMode schwarz_combinemode_;
224
225 // flags for setup
228
229 // Xpetra matrices
230 // K_ -> stiffness matrix
231 // M_ -> mass matrix
232 // A_ -> Problem matrix
233 // P_ -> Preconditioning matrix
234 RCP<Matrix> K_, M_, A_, P_;
235 RCP<MultiVector> Coords_, NullSpace_;
236
237 // Multigrid Hierarchy
238 RCP<Hierarchy> Hierarchy_;
239
240 // Factories and prototypes
241 RCP<TentativePFactory> TentPfact_;
242 RCP<PFactory> Pfact_;
243 RCP<PgPFactory> PgPfact_;
244 RCP<TransPFactory> TransPfact_;
245 RCP<GenericRFactory> Rfact_;
246 RCP<RAPFactory> Acfact_;
247 RCP<RAPShiftFactory> Acshift_;
248 RCP<AmalgamationFactory> Amalgfact_;
249 RCP<CoalesceDropFactory> Dropfact_;
250 RCP<UncoupledAggregationFactory> UCaggfact_;
251 RCP<CoarseMapFactory> CoarseMapfact_;
252 RCP<SmootherPrototype> smooProto_, coarsestSmooProto_;
253 RCP<SmootherFactory> smooFact_, coarsestSmooFact_;
254 Teuchos::ParameterList coarsestSmooList_;
255 std::string precType_;
256 Teuchos::ParameterList precList_;
257
258 // Operator and Preconditioner
259 RCP<MueLu::ShiftedLaplacianOperator<SC, LO, GO, NO> > MueLuOp_;
260 RCP<Tpetra::CrsMatrix<SC, LO, GO, NO> > TpetraA_;
261
262#ifdef HAVE_MUELU_TPETRA_INST_INT_INT
263 // Belos Linear Problem and Solver
264 RCP<LinearProblem> LinearProblem_;
265 RCP<SolverManager> SolverManager_;
266 RCP<SolverFactory> SolverFactory_;
267 RCP<Teuchos::ParameterList> BelosList_;
268#endif
269};
270
271} // namespace MueLu
272
273#define MUELU_SHIFTEDLAPLACIAN_SHORT
274
275#endif // if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
276
277#endif // MUELU_SHIFTEDLAPLACIAN_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Base class for MueLu classes.
Tpetra::MultiVector< SC, LO, GO, NO > TMV
void setPreconditioningMatrix(RCP< Matrix > &P)
RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > TpetraA_
Tpetra::Operator< SC, LO, GO, NO > OP
void setNullSpace(RCP< MultiVector > NullSpace)
RCP< MueLu::ShiftedLaplacianOperator< SC, LO, GO, NO > > MueLuOp_
RCP< SmootherPrototype > coarsestSmooProto_
void setcoords(RCP< MultiVector > &Coords)
Teuchos::ScalarTraits< Scalar >::magnitudeType GetResidual()
RCP< UncoupledAggregationFactory > UCaggfact_
int solve(const RCP< TMV > B, RCP< TMV > &X)
void setLevelShifts(std::vector< Scalar > levelshifts)
void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList)
void setProblemMatrix(RCP< Matrix > &A)
void multigrid_apply(const RCP< MultiVector > B, RCP< MultiVector > &X)
RCP< AmalgamationFactory > Amalgfact_
Tpetra::Vector< SC, LO, GO, NO > TVEC
RCP< CoalesceDropFactory > Dropfact_
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar