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