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 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jeremie Gaidamour (jngaida@sandia.gov)
40 // Jonathan Hu (jhu@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_SHIFTEDLAPLACIAN_DECL_HPP
47 #define MUELU_SHIFTEDLAPLACIAN_DECL_HPP
48 
49 // Xpetra
50 #include <Xpetra_Matrix_fwd.hpp>
51 #include <Xpetra_VectorFactory_fwd.hpp>
52 #include <Xpetra_MultiVectorFactory_fwd.hpp>
53 #include <Xpetra_TpetraMultiVector.hpp>
54 
55 // MueLu
56 #include "MueLu.hpp"
57 #include "MueLu_ConfigDefs.hpp"
58 
59 #if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
60 
61 #include <MueLu_BaseClass.hpp>
69 #include <MueLu_Hierarchy_fwd.hpp>
71 #include <MueLu_PFactory_fwd.hpp>
72 #include <MueLu_PgPFactory_fwd.hpp>
73 #include <MueLu_RAPFactory_fwd.hpp>
75 #include <MueLu_SaPFactory_fwd.hpp>
77 #include <MueLu_ShiftedLaplacianOperator.hpp>
83 #include <MueLu_Utilities_fwd.hpp>
84 
85 // Belos
86 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
87 #include <BelosConfigDefs.hpp>
88 #include <BelosLinearProblem.hpp>
89 #include <BelosSolverFactory.hpp>
90 #include <BelosTpetraAdapter.hpp>
91 #endif
92 
93 namespace MueLu {
94 
104  template <class Scalar = DefaultScalar,
107  class Node = DefaultNode>
108  class ShiftedLaplacian : public BaseClass {
109 #undef MUELU_SHIFTEDLAPLACIAN_SHORT
110 #include "MueLu_UseShortNames.hpp"
111 
112  typedef Tpetra::Vector<SC,LO,GO,NO> TVEC;
113  typedef Tpetra::MultiVector<SC,LO,GO,NO> TMV;
114  typedef Tpetra::Operator<SC,LO,GO,NO> OP;
115 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
116  typedef Belos::LinearProblem<SC,TMV,OP> LinearProblem;
117  typedef Belos::SolverManager<SC,TMV,OP> SolverManager;
118  typedef Belos::SolverFactory<SC,TMV,OP> SolverFactory;
119 #endif
120 
121  public:
122 
123  /*
124  FIXME 26-June-2015 JJH: This contructor is setting numerous defaults. However, they don't match the defaults
125  FIXME int the method setParameters(). There also isn't any parameter validation that I can see.
126  */
127 
130  numPDEs_(1),
131  Smoother_("schwarz"),
132  Aggregation_("uncoupled"),
133  Nullspace_("constant"),
134  numLevels_(5),
135  coarseGridSize_(100),
136  omega_(2.0*M_PI),
137  iters_(500),
138  blksize_(1),
139  tol_(1.0e-4),
140  nsweeps_(5),
141  ncycles_(1),
142  cycles_(8),
143  subiters_(10),
144  option_(1),
145  nproblems_(0),
146  solverType_(1),
147  restart_size_(100),
148  recycle_size_(25),
149  smoother_sweeps_(4),
150  smoother_damping_((SC)1.0),
151  krylov_type_(1),
154  ilu_leveloffill_(5.0),
155  ilu_abs_thresh_(0.0),
156  ilu_rel_thresh_(1.0),
158  ilu_drop_tol_(0.01),
159  ilu_fill_tol_(0.01),
160  ilu_relax_val_(1.0),
161  ilu_rowperm_("LargeDiag"),
162  ilu_colperm_("COLAMD"),
163  ilu_drop_rule_("DROP_BASIC"),
164  ilu_normtype_("INF_NORM"),
165  ilu_milutype_("SILU"),
166  schwarz_overlap_(0),
167  schwarz_usereorder_(true),
168  schwarz_combinemode_(Tpetra::ADD),
169  schwarz_ordermethod_("rcm"),
170  GridTransfersExist_(false),
171  isSymmetric_(true)
172  { }
173 
174  // Destructor
175  virtual ~ShiftedLaplacian();
176 
177  // Parameters
179 
180  // Set matrices
181  void setProblemMatrix(RCP<Matrix>& A);
182  void setProblemMatrix(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraA);
184  void setPreconditioningMatrix(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraP);
185  void setstiff(RCP<Matrix>& K);
186  void setstiff(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraK);
187  void setmass(RCP<Matrix>& M);
188  void setmass(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraM);
189  void setcoords(RCP<MultiVector>& Coords);
190  void setNullSpace(RCP<MultiVector> NullSpace);
191  void setLevelShifts(std::vector<Scalar> levelshifts);
192 
193  // initialize: set parameters and factories, construct
194  // prolongation and restriction matrices
195  void initialize();
196  // setupFastRAP: setup hierarchy with
197  // prolongators of the stiffness matrix
198  // constant complex shifts
199  void setupFastRAP();
200  // setupSlowRAP: setup hierarchy with
201  // prolongators of the stiffness matrix
202  // variable complex shifts
203  void setupSlowRAP();
204  // setupNormalRAP: setup hierarchy with
205  // prolongators of the preconditioning matrix
206  void setupNormalRAP();
207  // setupSolver: initialize Belos solver
208  void setupSolver();
209  // resetLinearProblem: for multiple frequencies;
210  // reset the Belos operator if the frequency changes
211  void resetLinearProblem();
212 
213 
214  // Solve phase
215  int solve(const RCP<TMV> B, RCP<TMV>& X);
216  void multigrid_apply(const RCP<MultiVector> B,
217  RCP<MultiVector>& X);
218  void multigrid_apply(const RCP<Tpetra::MultiVector<SC,LO,GO,NO> > B,
219  RCP<Tpetra::MultiVector<SC,LO,GO,NO> >& X);
220  int GetIterations();
222 
224 
225  private:
226 
227  // Problem options
228  // numPDEs_ -> number of DOFs at each node
229  int numPDEs_;
230 
231  // Multigrid options
232  // numLevels_ -> number of Multigrid levels
233  // coarseGridSize_ -> size of coarsest grid (if current level has less DOFs, stop coarsening)
236 
237  // Shifted Laplacian/Helmholtz parameters
238  double omega_;
239  std::vector<SC> levelshifts_;
240 
241  // Krylov solver inputs
242  // iters -> max number of iterations
243  // tol -> residual tolerance
245  double tol_;
249 
250  // Smoother parameters
261  Tpetra::CombineMode schwarz_combinemode_;
262  std::string schwarz_ordermethod_;
263 
264  // flags for setup
267 
268  // Xpetra matrices
269  // K_ -> stiffness matrix
270  // M_ -> mass matrix
271  // A_ -> Problem matrix
272  // P_ -> Preconditioning matrix
275 
276  // Multigrid Hierarchy
278 
279  // Factories and prototypes
295  std::string precType_;
297 
298  // Operator and Preconditioner
301 
302 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
303  // Belos Linear Problem and Solver
304  RCP<LinearProblem> LinearProblem_;
305  RCP<SolverManager> SolverManager_;
306  RCP<SolverFactory> SolverFactory_;
307  RCP<Teuchos::ParameterList> BelosList_;
308 #endif
309 
310  };
311 
312 }
313 
314 #define MUELU_SHIFTEDLAPLACIAN_SHORT
315 
316 #endif //if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
317 
318 #endif // MUELU_SHIFTEDLAPLACIAN_DECL_HPP
RCP< SmootherFactory > coarsestSmooFact_
MueLu::DefaultLocalOrdinal LocalOrdinal
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)
Teuchos::ScalarTraits< Scalar >::magnitudeType GetResidual()
Tpetra::MultiVector< SC, LO, GO, NO > TMV
RCP< SmootherPrototype > smooProto_
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void setLevelShifts(std::vector< Scalar > levelshifts)
RCP< MueLu::ShiftedLaplacianOperator< SC, LO, GO, NO > > MueLuOp_
Teuchos::ParameterList coarsestSmooList_
void setcoords(RCP< MultiVector > &Coords)
RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > TpetraA_
RCP< UncoupledAggregationFactory > UCaggfact_
Tpetra::Operator< SC, LO, GO, NO > OP
int solve(const RCP< TMV > B, RCP< TMV > &X)
Base class for MueLu classes.
RCP< SmootherPrototype > coarsestSmooProto_
RCP< CoupledAggregationFactory > Aggfact_
Tpetra::Vector< SC, LO, GO, NO > TVEC
void setProblemMatrix(RCP< Matrix > &A)
RCP< TentativePFactory > TentPfact_
Shifted Laplacian Helmholtz solver.
RCP< AmalgamationFactory > Amalgfact_
void setNullSpace(RCP< MultiVector > NullSpace)
void setPreconditioningMatrix(RCP< Matrix > &P)