MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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>
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>
32 #include <MueLu_Hierarchy_fwd.hpp>
34 #include <MueLu_PFactory_fwd.hpp>
35 #include <MueLu_PgPFactory_fwd.hpp>
36 #include <MueLu_RAPFactory_fwd.hpp>
38 #include <MueLu_SaPFactory_fwd.hpp>
40 #include <MueLu_ShiftedLaplacianOperator.hpp>
46 #include <MueLu_Utilities_fwd.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 namespace MueLu {
57 
67 template <class Scalar = DefaultScalar,
70  class Node = DefaultNode>
71 class ShiftedLaplacian : public BaseClass {
72 #undef MUELU_SHIFTEDLAPLACIAN_SHORT
73 #include "MueLu_UseShortNames.hpp"
74 
78 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
79  typedef Belos::LinearProblem<SC, TMV, OP> LinearProblem;
80  typedef Belos::SolverManager<SC, TMV, OP> SolverManager;
81  typedef Belos::SolverFactory<SC, TMV, OP> SolverFactory;
82 #endif
83 
84  public:
85  /*
86  FIXME 26-June-2015 JJH: This contructor is setting numerous defaults. However, they don't match the defaults
87  FIXME int the method setParameters(). There also isn't any parameter validation that I can see.
88  */
89 
92  : numPDEs_(1)
93  , Smoother_("schwarz")
94  , Aggregation_("uncoupled")
95  , Nullspace_("constant")
96  , numLevels_(5)
97  , coarseGridSize_(100)
98  , omega_(2.0 * M_PI)
99  , iters_(500)
100  , blksize_(1)
101  , tol_(1.0e-4)
102  , nsweeps_(5)
103  , ncycles_(1)
104  , cycles_(8)
105  , subiters_(10)
106  , option_(1)
107  , nproblems_(0)
108  , solverType_(1)
109  , restart_size_(100)
110  , recycle_size_(25)
111  , smoother_sweeps_(4)
112  , smoother_damping_((SC)1.0)
113  , krylov_type_(1)
114  , krylov_iterations_(5)
116  , ilu_leveloffill_(5.0)
117  , ilu_abs_thresh_(0.0)
118  , ilu_rel_thresh_(1.0)
119  , ilu_diagpivotthresh_(0.1)
120  , ilu_drop_tol_(0.01)
121  , ilu_fill_tol_(0.01)
122  , ilu_relax_val_(1.0)
123  , ilu_rowperm_("LargeDiag")
124  , ilu_colperm_("COLAMD")
125  , ilu_drop_rule_("DROP_BASIC")
126  , ilu_normtype_("INF_NORM")
127  , ilu_milutype_("SILU")
128  , schwarz_overlap_(0)
129  , schwarz_usereorder_(true)
130  , schwarz_combinemode_(Tpetra::ADD)
131  , schwarz_ordermethod_("rcm")
132  , GridTransfersExist_(false)
133  , isSymmetric_(true) {}
134 
135  // Destructor
136  virtual ~ShiftedLaplacian();
137 
138  // Parameters
140 
141  // Set matrices
142  void setProblemMatrix(RCP<Matrix>& A);
146  void setstiff(RCP<Matrix>& K);
148  void setmass(RCP<Matrix>& M);
150  void setcoords(RCP<MultiVector>& Coords);
151  void setNullSpace(RCP<MultiVector> NullSpace);
152  void setLevelShifts(std::vector<Scalar> levelshifts);
153 
154  // initialize: set parameters and factories, construct
155  // prolongation and restriction matrices
156  void initialize();
157  // setupFastRAP: setup hierarchy with
158  // prolongators of the stiffness matrix
159  // constant complex shifts
160  void setupFastRAP();
161  // setupSlowRAP: setup hierarchy with
162  // prolongators of the stiffness matrix
163  // variable complex shifts
164  void setupSlowRAP();
165  // setupNormalRAP: setup hierarchy with
166  // prolongators of the preconditioning matrix
167  void setupNormalRAP();
168  // setupSolver: initialize Belos solver
169  void setupSolver();
170  // resetLinearProblem: for multiple frequencies;
171  // reset the Belos operator if the frequency changes
172  void resetLinearProblem();
173 
174  // Solve phase
175  int solve(const RCP<TMV> B, RCP<TMV>& X);
176  void multigrid_apply(const RCP<MultiVector> B,
177  RCP<MultiVector>& X);
180  int GetIterations();
182 
184 
185  private:
186  // Problem options
187  // numPDEs_ -> number of DOFs at each node
188  int numPDEs_;
189 
190  // Multigrid options
191  // numLevels_ -> number of Multigrid levels
192  // coarseGridSize_ -> size of coarsest grid (if current level has less DOFs, stop coarsening)
195 
196  // Shifted Laplacian/Helmholtz parameters
197  double omega_;
198  std::vector<SC> levelshifts_;
199 
200  // Krylov solver inputs
201  // iters -> max number of iterations
202  // tol -> residual tolerance
204  double tol_;
208 
209  // Smoother parameters
221  std::string schwarz_ordermethod_;
222 
223  // flags for setup
226 
227  // Xpetra matrices
228  // K_ -> stiffness matrix
229  // M_ -> mass matrix
230  // A_ -> Problem matrix
231  // P_ -> Preconditioning matrix
234 
235  // Multigrid Hierarchy
237 
238  // Factories and prototypes
253  std::string precType_;
255 
256  // Operator and Preconditioner
259 
260 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
261  // Belos Linear Problem and Solver
262  RCP<LinearProblem> LinearProblem_;
263  RCP<SolverManager> SolverManager_;
264  RCP<SolverFactory> SolverFactory_;
265  RCP<Teuchos::ParameterList> BelosList_;
266 #endif
267 };
268 
269 } // namespace MueLu
270 
271 #define MUELU_SHIFTEDLAPLACIAN_SHORT
272 
273 #endif // if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
274 
275 #endif // MUELU_SHIFTEDLAPLACIAN_DECL_HPP
RCP< SmootherFactory > coarsestSmooFact_
MueLu::DefaultLocalOrdinal LocalOrdinal
Tpetra::MultiVector< SC, LO, GO, NO > TMV
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
void multigrid_apply(const RCP< MultiVector > B, RCP< MultiVector > &X)
RCP< CoarseMapFactory > CoarseMapfact_
RCP< CoalesceDropFactory > Dropfact_
MueLu::DefaultNode Node
void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList)
Tpetra::Operator< SC, LO, GO, NO > OP
Teuchos::ScalarTraits< Scalar >::magnitudeType GetResidual()
RCP< SmootherPrototype > smooProto_
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void setLevelShifts(std::vector< Scalar > levelshifts)
Tpetra::Vector< SC, LO, GO, NO > TVEC
Teuchos::ParameterList coarsestSmooList_
void setcoords(RCP< MultiVector > &Coords)
RCP< UncoupledAggregationFactory > UCaggfact_
int solve(const RCP< TMV > B, RCP< TMV > &X)
Base class for MueLu classes.
RCP< SmootherPrototype > coarsestSmooProto_
void setProblemMatrix(RCP< Matrix > &A)
Scalar SC
RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > TpetraA_
RCP< TentativePFactory > TentPfact_
Shifted Laplacian Helmholtz solver.
RCP< MueLu::ShiftedLaplacianOperator< SC, LO, GO, NO > > MueLuOp_
RCP< AmalgamationFactory > Amalgfact_
void setNullSpace(RCP< MultiVector > NullSpace)
void setPreconditioningMatrix(RCP< Matrix > &P)