MueLu  Version of the Day
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>
67 #include <MueLu_Hierarchy_fwd.hpp>
69 #include <MueLu_PFactory_fwd.hpp>
70 #include <MueLu_PgPFactory_fwd.hpp>
71 #include <MueLu_RAPFactory_fwd.hpp>
73 #include <MueLu_SaPFactory_fwd.hpp>
75 #include <MueLu_ShiftedLaplacianOperator.hpp>
81 #include <MueLu_Utilities_fwd.hpp>
82 
83 // Belos
84 #include <BelosConfigDefs.hpp>
85 #include <BelosLinearProblem.hpp>
86 #include <BelosSolverFactory.hpp>
87 #include <BelosTpetraAdapter.hpp>
88 
89 namespace MueLu {
90 
100  template <class Scalar = Xpetra::Matrix<>::scalar_type,
101  class LocalOrdinal = typename Xpetra::Matrix<Scalar>::local_ordinal_type,
102  class GlobalOrdinal = typename Xpetra::Matrix<Scalar, LocalOrdinal>::global_ordinal_type,
103  class Node = typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
104  class ShiftedLaplacian : public BaseClass {
105 #undef MUELU_SHIFTEDLAPLACIAN_SHORT
106 #include "MueLu_UseShortNames.hpp"
107 
108  typedef Tpetra::Vector<SC,LO,GO,NO> TVEC;
109  typedef Tpetra::MultiVector<SC,LO,GO,NO> TMV;
110  typedef Tpetra::Operator<SC,LO,GO,NO> OP;
111  typedef Belos::LinearProblem<SC,TMV,OP> LinearProblem;
112  typedef Belos::SolverManager<SC,TMV,OP> SolverManager;
113  typedef Belos::SolverFactory<SC,TMV,OP> SolverFactory;
114 
115  public:
116 
117  /*
118  FIXME 26-June-2015 JJH: This contructor is setting numerous defaults. However, they don't match the defaults
119  FIXME int the method setParameters(). There also isn't any parameter validation that I can see.
120  */
121 
124  numPDEs_(1),
125  Smoother_("schwarz"),
126  Aggregation_("uncoupled"),
127  Nullspace_("constant"),
128  numLevels_(5),
129  coarseGridSize_(100),
130  omega_(2.0*M_PI),
131  iters_(500),
132  blksize_(1),
133  tol_(1.0e-4),
134  nsweeps_(5),
135  ncycles_(1),
136  cycles_(8),
137  subiters_(10),
138  option_(1),
139  nproblems_(0),
140  solverType_(1),
141  restart_size_(100),
142  recycle_size_(25),
143  smoother_sweeps_(4),
144  smoother_damping_((SC)1.0),
145  krylov_type_(1),
148  ilu_leveloffill_(5.0),
149  ilu_abs_thresh_(0.0),
150  ilu_rel_thresh_(1.0),
152  ilu_drop_tol_(0.01),
153  ilu_fill_tol_(0.01),
154  ilu_relax_val_(1.0),
155  ilu_rowperm_("LargeDiag"),
156  ilu_colperm_("COLAMD"),
157  ilu_drop_rule_("DROP_BASIC"),
158  ilu_normtype_("INF_NORM"),
159  ilu_milutype_("SILU"),
160  schwarz_overlap_(0),
161  schwarz_usereorder_(true),
162  schwarz_combinemode_(Tpetra::ADD),
163  schwarz_ordermethod_("rcm"),
164  GridTransfersExist_(false),
165  isSymmetric_(true)
166  { }
167 
168  // Destructor
169  virtual ~ShiftedLaplacian();
170 
171  // Parameters
172  void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList);
173 
174  // Set matrices
175  void setProblemMatrix(RCP<Matrix>& A);
176  void setProblemMatrix(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraA);
177  void setPreconditioningMatrix(RCP<Matrix>& P);
178  void setPreconditioningMatrix(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraP);
179  void setstiff(RCP<Matrix>& K);
180  void setstiff(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraK);
181  void setmass(RCP<Matrix>& M);
182  void setmass(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraM);
183  void setcoords(RCP<MultiVector>& Coords);
184  void setNullSpace(RCP<MultiVector> NullSpace);
185  void setLevelShifts(std::vector<Scalar> levelshifts);
186 
187  // initialize: set parameters and factories, construct
188  // prolongation and restriction matrices
189  void initialize();
190  // setupFastRAP: setup hierarchy with
191  // prolongators of the stiffness matrix
192  // constant complex shifts
193  void setupFastRAP();
194  // setupSlowRAP: setup hierarchy with
195  // prolongators of the stiffness matrix
196  // variable complex shifts
197  void setupSlowRAP();
198  // setupNormalRAP: setup hierarchy with
199  // prolongators of the preconditioning matrix
200  void setupNormalRAP();
201  // setupSolver: initialize Belos solver
202  void setupSolver();
203  // resetLinearProblem: for multiple frequencies;
204  // reset the Belos operator if the frequency changes
205  void resetLinearProblem();
206 
207 
208  // Solve phase
209  int solve(const RCP<TMV> B, RCP<TMV>& X);
210  void multigrid_apply(const RCP<MultiVector> B,
211  RCP<MultiVector>& X);
212  void multigrid_apply(const RCP<Tpetra::MultiVector<SC,LO,GO,NO> > B,
213  RCP<Tpetra::MultiVector<SC,LO,GO,NO> >& X);
214  int GetIterations();
215  double GetResidual();
216 
217  RCP<FactoryManager> Manager_;
218 
219  private:
220 
221  // Problem options
222  // numPDEs_ -> number of DOFs at each node
223  int numPDEs_;
224 
225  // Multigrid options
226  // numLevels_ -> number of Multigrid levels
227  // coarseGridSize_ -> size of coarsest grid (if current level has less DOFs, stop coarsening)
230 
231  // Shifted Laplacian/Helmholtz parameters
232  double omega_;
233  std::vector<SC> levelshifts_;
234 
235  // Krylov solver inputs
236  // iters -> max number of iterations
237  // tol -> residual tolerance
239  double tol_;
243 
244  // Smoother parameters
255  Tpetra::CombineMode schwarz_combinemode_;
256  std::string schwarz_ordermethod_;
257 
258  // flags for setup
261 
262  // Xpetra matrices
263  // K_ -> stiffness matrix
264  // M_ -> mass matrix
265  // A_ -> Problem matrix
266  // P_ -> Preconditioning matrix
267  RCP<Matrix> K_, M_, A_, P_;
268  RCP<MultiVector> Coords_, NullSpace_;
269 
270  // Multigrid Hierarchy
271  RCP<Hierarchy> Hierarchy_;
272 
273  // Factories and prototypes
274  RCP<TentativePFactory> TentPfact_;
275  RCP<PFactory> Pfact_;
276  RCP<PgPFactory> PgPfact_;
277  RCP<TransPFactory> TransPfact_;
278  RCP<GenericRFactory> Rfact_;
279  RCP<RAPFactory> Acfact_;
280  RCP<RAPShiftFactory> Acshift_;
281  RCP<CoalesceDropFactory> Dropfact_;
282  RCP<CoupledAggregationFactory> Aggfact_;
283  RCP<UncoupledAggregationFactory> UCaggfact_;
284  RCP<SmootherPrototype> smooProto_, coarsestSmooProto_;
285  RCP<SmootherFactory> smooFact_, coarsestSmooFact_;
286  Teuchos::ParameterList coarsestSmooList_;
287  std::string precType_;
288  Teuchos::ParameterList precList_;
289 
290  // Operator and Preconditioner
291  RCP< MueLu::ShiftedLaplacianOperator<SC,LO,GO,NO> > MueLuOp_;
292  RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> > TpetraA_;
293 
294  // Belos Linear Problem and Solver
295  RCP<LinearProblem> LinearProblem_;
296  RCP<SolverManager> SolverManager_;
297  RCP<SolverFactory> SolverFactory_;
298  RCP<Teuchos::ParameterList> BelosList_;
299 
300  };
301 
302 }
303 
304 #define MUELU_SHIFTEDLAPLACIAN_SHORT
305 
306 #endif //if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
307 
308 #endif // MUELU_SHIFTEDLAPLACIAN_DECL_HPP
RCP< SmootherFactory > coarsestSmooFact_
Belos::SolverManager< SC, TMV, OP > SolverManager
void multigrid_apply(const RCP< MultiVector > B, RCP< MultiVector > &X)
Namespace for MueLu classes and methods.
RCP< CoalesceDropFactory > Dropfact_
void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList)
Tpetra::MultiVector< SC, LO, GO, NO > TMV
RCP< SmootherPrototype > smooProto_
void setLevelShifts(std::vector< Scalar > levelshifts)
RCP< Teuchos::ParameterList > BelosList_
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)
RCP< SmootherPrototype > coarsestSmooProto_
RCP< CoupledAggregationFactory > Aggfact_
Tpetra::Vector< SC, LO, GO, NO > TVEC
void setProblemMatrix(RCP< Matrix > &A)
RCP< TentativePFactory > TentPfact_
Belos::LinearProblem< SC, TMV, OP > LinearProblem
Belos::SolverFactory< SC, TMV, OP > SolverFactory
void setNullSpace(RCP< MultiVector > NullSpace)
void setPreconditioningMatrix(RCP< Matrix > &P)