MueLu  Version of the Day
MueLu_ShiftedLaplacian_def.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_DEF_HPP
47 #define MUELU_SHIFTEDLAPLACIAN_DEF_HPP
48 
50 
51 #if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
52 
53 #include <MueLu_CoalesceDropFactory.hpp>
54 #include <MueLu_CoupledAggregationFactory.hpp>
55 #include <MueLu_CoupledRBMFactory.hpp>
56 #include <MueLu_DirectSolver.hpp>
57 #include <MueLu_GenericRFactory.hpp>
58 #include <MueLu_Hierarchy.hpp>
59 #include <MueLu_Ifpack2Smoother.hpp>
60 #include <MueLu_PFactory.hpp>
61 #include <MueLu_PgPFactory.hpp>
62 #include <MueLu_RAPFactory.hpp>
63 #include <MueLu_RAPShiftFactory.hpp>
64 #include <MueLu_SaPFactory.hpp>
65 #include <MueLu_ShiftedLaplacian.hpp>
66 #include <MueLu_ShiftedLaplacianOperator.hpp>
67 #include <MueLu_SmootherFactory.hpp>
68 #include <MueLu_SmootherPrototype.hpp>
69 #include <MueLu_TentativePFactory.hpp>
70 #include <MueLu_TransPFactory.hpp>
71 #include <MueLu_UncoupledAggregationFactory.hpp>
72 #include <MueLu_Utilities.hpp>
73 
74 namespace MueLu {
75 
76 // Destructor
77 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79 
80 // Input
81 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82 void ShiftedLaplacian<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList) {
83 
84  // Parameters
85  coarseGridSize_ = paramList->get("MueLu: coarse size", 1000);
86  numLevels_ = paramList->get("MueLu: levels", 3);
87  int stype = paramList->get("MueLu: smoother", 8);
88  if(stype==1) { Smoother_="jacobi"; }
89  else if(stype==2) { Smoother_="gauss-seidel"; }
90  else if(stype==3) { Smoother_="symmetric gauss-seidel"; }
91  else if(stype==4) { Smoother_="chebyshev"; }
92  else if(stype==5) { Smoother_="krylov"; }
93  else if(stype==6) { Smoother_="ilut"; }
94  else if(stype==7) { Smoother_="riluk"; }
95  else if(stype==8) { Smoother_="schwarz"; }
96  else if(stype==9) { Smoother_="superilu"; }
97  else if(stype==10) { Smoother_="superlu"; }
98  else { Smoother_="schwarz"; }
99  smoother_sweeps_ = paramList->get("MueLu: sweeps", 5);
100  smoother_damping_ = paramList->get("MueLu: relax val", 1.0);
101  ncycles_ = paramList->get("MueLu: cycles", 1);
102  iters_ = paramList->get("MueLu: iterations", 500);
103  solverType_ = paramList->get("MueLu: solver type", 1);
104  restart_size_ = paramList->get("MueLu: restart size", 100);
105  recycle_size_ = paramList->get("MueLu: recycle size", 25);
106  isSymmetric_ = paramList->get("MueLu: symmetric", true);
107  ilu_leveloffill_ = paramList->get("MueLu: level-of-fill", 5);
108  ilu_abs_thresh_ = paramList->get("MueLu: abs thresh", 0.0);
109  ilu_rel_thresh_ = paramList->get("MueLu: rel thresh", 1.0);
110  ilu_diagpivotthresh_ = paramList->get("MueLu: piv thresh", 0.1);
111  ilu_drop_tol_ = paramList->get("MueLu: drop tol", 0.01);
112  ilu_fill_tol_ = paramList->get("MueLu: fill tol", 0.01);
113  schwarz_overlap_ = paramList->get("MueLu: overlap", 0);
114  schwarz_usereorder_ = paramList->get("MueLu: use reorder", true);
115  int combinemode = paramList->get("MueLu: combine mode", 1);
116  if(combinemode==0) { schwarz_combinemode_ = Tpetra::ZERO; }
117  else { schwarz_combinemode_ = Tpetra::ADD; }
118  tol_ = paramList->get("MueLu: tolerance", 0.001);
119 
120 }
121 
122 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
124 
125  A_=A;
126  if(A_!=Teuchos::null)
127  TpetraA_ = Utils::Op2NonConstTpetraCrs(A_);
128  if(LinearProblem_!=Teuchos::null)
129  LinearProblem_ -> setOperator ( TpetraA_ );
130 
131 }
132 
133 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134 void ShiftedLaplacian<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setProblemMatrix(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraA) {
135 
136  TpetraA_=TpetraA;
137  if(LinearProblem_!=Teuchos::null)
138  LinearProblem_ -> setOperator ( TpetraA_ );
139 
140 }
141 
142 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
144 
145  P_=P;
146  GridTransfersExist_=false;
147 
148 }
149 
150 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152 
153  RCP< Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Atmp
154  = rcp( new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(TpetraP) );
155  P_= rcp( new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Atmp) );
156  GridTransfersExist_=false;
157 
158 }
159 
160 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
162 
163  K_=K;
164 
165 }
166 
167 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
168 void ShiftedLaplacian<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setstiff(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraK) {
169 
170  RCP< Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Atmp
171  = rcp( new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(TpetraK) );
172  K_= rcp( new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Atmp) );
173 
174 }
175 
176 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
178 
179  M_=M;
180 
181 }
182 
183 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
184 void ShiftedLaplacian<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setmass(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraM) {
185 
186  RCP< Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Atmp
187  = rcp( new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(TpetraM) );
188  M_= rcp( new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Atmp) );
189 
190 }
191 
192 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
194 
195  Coords_=Coords;
196 
197 }
198 
199 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201 
202  NullSpace_=NullSpace;
203 
204 }
205 
206 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
208 
209  levelshifts_=levelshifts;
210  numLevels_=levelshifts_.size();
211 
212 }
213 
214 // initialize
215 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
217 
218  TentPfact_ = rcp( new TentativePFactory );
219  Pfact_ = rcp( new SaPFactory );
220  PgPfact_ = rcp( new PgPFactory );
221  TransPfact_= rcp( new TransPFactory );
222  Rfact_ = rcp( new GenericRFactory );
223  Acfact_ = rcp( new RAPFactory );
224  Acshift_ = rcp( new RAPShiftFactory );
225  Dropfact_ = rcp( new CoalesceDropFactory );
226  Aggfact_ = rcp( new CoupledAggregationFactory );
227  UCaggfact_ = rcp( new UncoupledAggregationFactory );
228  Manager_ = rcp( new FactoryManager );
229  if(isSymmetric_==true) {
230  Manager_ -> SetFactory("P", Pfact_);
231  Manager_ -> SetFactory("R", TransPfact_);
232  }
233  else {
234  Manager_ -> SetFactory("P", PgPfact_);
235  Manager_ -> SetFactory("R", Rfact_);
236  solverType_ = 10;
237  }
238  Manager_ -> SetFactory("Ptent", TentPfact_);
239  Teuchos::ParameterList params;
240  params.set("lightweight wrap",true);
241  params.set("aggregation: drop scheme","classical");
242  Dropfact_ -> SetParameterList(params);
243  Manager_ -> SetFactory("Graph", Dropfact_);
244  if(Aggregation_=="coupled") {
245  Manager_ -> SetFactory("Aggregates", Aggfact_ );
246  }
247  else {
248  Manager_ -> SetFactory("Aggregates", UCaggfact_ );
249  }
250 
251  // choose smoother
252  if(Smoother_=="jacobi") {
253  precType_ = "RELAXATION";
254  precList_.set("relaxation: type", "Jacobi");
255  precList_.set("relaxation: sweeps", smoother_sweeps_);
256  precList_.set("relaxation: damping factor", smoother_damping_);
257  }
258  else if(Smoother_=="gauss-seidel") {
259  precType_ = "RELAXATION";
260  precList_.set("relaxation: type", "Gauss-Seidel");
261  precList_.set("relaxation: sweeps", smoother_sweeps_);
262  precList_.set("relaxation: damping factor", smoother_damping_);
263  }
264  else if(Smoother_=="symmetric gauss-seidel") {
265  precType_ = "RELAXATION";
266  precList_.set("relaxation: type", "Symmetric Gauss-Seidel");
267  precList_.set("relaxation: sweeps", smoother_sweeps_);
268  precList_.set("relaxation: damping factor", smoother_damping_);
269  }
270  else if(Smoother_=="chebyshev") {
271  precType_ = "CHEBYSHEV";
272  }
273  else if(Smoother_=="krylov") {
274  precType_ = "KRYLOV";
275  precList_.set("krylov: iteration type", krylov_type_);
276  precList_.set("krylov: number of iterations", krylov_iterations_);
277  precList_.set("krylov: residual tolerance",1.0e-8);
278  precList_.set("krylov: block size",1);
279  precList_.set("krylov: preconditioner type", krylov_preconditioner_);
280  precList_.set("relaxation: sweeps",1);
281  solverType_=10;
282  }
283  else if(Smoother_=="ilut") {
284  precType_ = "ILUT";
285  precList_.set("fact: ilut level-of-fill", ilu_leveloffill_);
286  precList_.set("fact: absolute threshold", ilu_abs_thresh_);
287  precList_.set("fact: relative threshold", ilu_rel_thresh_);
288  precList_.set("fact: drop tolerance", ilu_drop_tol_);
289  precList_.set("fact: relax value", ilu_relax_val_);
290  }
291  else if(Smoother_=="riluk") {
292  precType_ = "RILUK";
293  precList_.set("fact: iluk level-of-fill", ilu_leveloffill_);
294  precList_.set("fact: absolute threshold", ilu_abs_thresh_);
295  precList_.set("fact: relative threshold", ilu_rel_thresh_);
296  precList_.set("fact: drop tolerance", ilu_drop_tol_);
297  precList_.set("fact: relax value", ilu_relax_val_);
298  }
299  else if(Smoother_=="schwarz") {
300  precType_ = "SCHWARZ";
301  precList_.set("schwarz: overlap level", schwarz_overlap_);
302  precList_.set("schwarz: combine mode", schwarz_combinemode_);
303  precList_.set("schwarz: use reordering", schwarz_usereorder_);
304  precList_.set("schwarz: filter singletons", true);
305  precList_.set("order_method",schwarz_ordermethod_);
306  precList_.sublist("schwarz: reordering list").set("order_method",schwarz_ordermethod_);
307  precList_.sublist("schwarz: subdomain solver parameters").set("fact: ilut level-of-fill", ilu_leveloffill_);
308  precList_.sublist("schwarz: subdomain solver parameters").set("fact: absolute threshold", ilu_abs_thresh_);
309  precList_.sublist("schwarz: subdomain solver parameters").set("fact: relative threshold", ilu_rel_thresh_);
310  precList_.sublist("schwarz: subdomain solver parameters").set("fact: drop tolerance", ilu_drop_tol_);
311  precList_.sublist("schwarz: subdomain solver parameters").set("fact: relax value", ilu_relax_val_);
312  }
313  else if(Smoother_=="superilu") {
314  precType_ = "superlu";
315  precList_.set("RowPerm", ilu_rowperm_);
316  precList_.set("ColPerm", ilu_colperm_);
317  precList_.set("DiagPivotThresh", ilu_diagpivotthresh_);
318  precList_.set("ILU_DropRule",ilu_drop_rule_);
319  precList_.set("ILU_DropTol",ilu_drop_tol_);
320  precList_.set("ILU_FillFactor",ilu_leveloffill_);
321  precList_.set("ILU_Norm",ilu_normtype_);
322  precList_.set("ILU_MILU",ilu_milutype_);
323  precList_.set("ILU_FillTol",ilu_fill_tol_);
324  precList_.set("ILU_Flag",true);
325  }
326  else if(Smoother_=="superlu") {
327  precType_ = "superlu";
328  precList_.set("ColPerm", ilu_colperm_);
329  precList_.set("DiagPivotThresh", ilu_diagpivotthresh_);
330  }
331  // construct smoother
332  smooProto_ = rcp( new Ifpack2Smoother(precType_,precList_) );
333  smooFact_ = rcp( new SmootherFactory(smooProto_) );
334 #if defined(HAVE_MUELU_AMESOS2) and defined(HAVE_AMESOS2_SUPERLU)
335  coarsestSmooProto_ = rcp( new DirectSolver("Superlu",coarsestSmooList_) );
336 #elif defined(HAVE_MUELU_AMESOS2) and defined(HAVE_AMESOS2_KLU2)
337  coarsestSmooProto_ = rcp( new DirectSolver("Klu",coarsestSmooList_) );
338 #elif defined(HAVE_MUELU_AMESOS2) and defined(HAVE_AMESOS2_SUPERLUDIST)
339  coarsestSmooProto_ = rcp( new DirectSolver("Superludist",coarsestSmooList_) );
340 #else
341  coarsestSmooProto_ = rcp( new Ifpack2Smoother(precType_,precList_) );
342 #endif
343  coarsestSmooFact_ = rcp( new SmootherFactory(coarsestSmooProto_, Teuchos::null) );
344 
345  // For setupSlowRAP and setupFastRAP, the prolongation/restriction matrices
346  // are constructed with the stiffness matrix. These matrices are kept for future
347  // setup calls; this is achieved by calling Hierarchy->Keep(). It is particularly
348  // useful for multiple frequency problems - when the frequency/preconditioner
349  // changes, you only compute coarse grids (RAPs) and setup level smoothers when
350  // you call Hierarchy->Setup().
351  if(K_!=Teuchos::null) {
352  Manager_ -> SetFactory("Smoother", Teuchos::null);
353  Manager_ -> SetFactory("CoarseSolver", Teuchos::null);
354  Hierarchy_ = rcp( new Hierarchy(K_) );
355  if(NullSpace_!=Teuchos::null)
356  Hierarchy_ -> GetLevel(0) -> Set("Nullspace", NullSpace_);
357  if(isSymmetric_==true) {
358  Hierarchy_ -> Keep("P", Pfact_.get());
359  Hierarchy_ -> Keep("R", TransPfact_.get());
360  Hierarchy_ -> SetImplicitTranspose(true);
361  }
362  else {
363  Hierarchy_ -> Keep("P", PgPfact_.get());
364  Hierarchy_ -> Keep("R", Rfact_.get());
365  }
366  Hierarchy_ -> Keep("Ptent", TentPfact_.get());
367  Hierarchy_ -> SetMaxCoarseSize( coarseGridSize_ );
368  Hierarchy_ -> Setup(*Manager_, 0, numLevels_);
369  GridTransfersExist_=true;
370  }
371  // Use preconditioning matrix to setup prolongation/restriction operators
372  else {
373  Manager_ -> SetFactory("Smoother", smooFact_);
374  Manager_ -> SetFactory("CoarseSolver", coarsestSmooFact_);
375  Hierarchy_ = rcp( new Hierarchy(P_) );
376  if(NullSpace_!=Teuchos::null)
377  Hierarchy_ -> GetLevel(0) -> Set("Nullspace", NullSpace_);
378  if(isSymmetric_==true)
379  Hierarchy_ -> SetImplicitTranspose(true);
380  Hierarchy_ -> SetMaxCoarseSize( coarseGridSize_ );
381  Hierarchy_ -> Setup(*Manager_, 0, numLevels_);
382  GridTransfersExist_=true;
383  }
384 
385  // Belos Linear Problem and Solver Manager
386  BelosList_ = rcp( new Teuchos::ParameterList("GMRES") );
387  BelosList_ -> set("Maximum Iterations",iters_ );
388  BelosList_ -> set("Convergence Tolerance",tol_ );
389  BelosList_ -> set("Verbosity", Belos::Errors + Belos::Warnings + Belos::StatusTestDetails);
390  BelosList_ -> set("Output Frequency",1);
391  BelosList_ -> set("Output Style",Belos::Brief);
392  BelosList_ -> set("Num Blocks",restart_size_);
393  BelosList_ -> set("Num Recycled Blocks",recycle_size_);
394 
395 }
396 
397 // setup coarse grids for new frequency
398 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
400 
401  int numLevels = Hierarchy_ -> GetNumLevels();
402  Manager_ -> SetFactory("Smoother", smooFact_);
403  Manager_ -> SetFactory("CoarseSolver", coarsestSmooFact_);
404  Hierarchy_ -> GetLevel(0) -> Set("A", P_);
405  Hierarchy_ -> Setup(*Manager_, 0, numLevels);
406  setupSolver();
407 
408 }
409 
410 // setup coarse grids for new frequency
411 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
413 
414  int numLevels = Hierarchy_ -> GetNumLevels();
415  Acshift_->SetShifts(levelshifts_);
416  Manager_ -> SetFactory("Smoother", smooFact_);
417  Manager_ -> SetFactory("CoarseSolver", coarsestSmooFact_);
418  Manager_ -> SetFactory("A", Acshift_);
419  Manager_ -> SetFactory("K", Acshift_);
420  Manager_ -> SetFactory("M", Acshift_);
421  Hierarchy_ -> GetLevel(0) -> Set("A", P_);
422  Hierarchy_ -> GetLevel(0) -> Set("K", K_);
423  Hierarchy_ -> GetLevel(0) -> Set("M", M_);
424  Hierarchy_ -> Setup(*Manager_, 0, numLevels);
425  setupSolver();
426 
427 }
428 
429 // setup coarse grids for new frequency
430 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
432 
433  // Only setup hierarchy again if preconditioning matrix has changed
434  if( GridTransfersExist_ == false ) {
435  Hierarchy_ = rcp( new Hierarchy(P_) );
436  if(NullSpace_!=Teuchos::null)
437  Hierarchy_ -> GetLevel(0) -> Set("Nullspace", NullSpace_);
438  if(isSymmetric_==true)
439  Hierarchy_ -> SetImplicitTranspose(true);
440  Hierarchy_ -> SetMaxCoarseSize( coarseGridSize_ );
441  Hierarchy_ -> Setup(*Manager_, 0, numLevels_);
442  GridTransfersExist_=true;
443  }
444  setupSolver();
445 
446 }
447 
448 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
450 
451  // Define Preconditioner and Operator
453  (Hierarchy_, A_, ncycles_, subiters_, option_, tol_) );
454  // Belos Linear Problem
455  if(LinearProblem_==Teuchos::null)
456  LinearProblem_ = rcp( new LinearProblem );
457  LinearProblem_ -> setOperator ( TpetraA_ );
458  LinearProblem_ -> setRightPrec( MueLuOp_ );
459  if(SolverManager_==Teuchos::null) {
460  std::string solverName;
461  SolverFactory_= rcp( new SolverFactory() );
462  if(solverType_==1) { solverName="Block GMRES"; }
463  else if(solverType_==2) { solverName="Recycling GMRES"; }
464  else { solverName="Flexible GMRES"; }
465  SolverManager_ = SolverFactory_->create( solverName, BelosList_ );
466  SolverManager_ -> setProblem( LinearProblem_ );
467  }
468 
469 }
470 
471 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
473 {
474  LinearProblem_ -> setOperator ( TpetraA_ );
475 }
476 
477 // Solve phase
478 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
480 {
481  // Set left and right hand sides for Belos
482  LinearProblem_ -> setProblem(X, B);
483  // iterative solve
484  SolverManager_ -> solve();
485  return 0;
486 }
487 
488 // Solve phase
489 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
491  RCP<MultiVector>& X)
492 {
493  // Set left and right hand sides for Belos
494  Hierarchy_ -> Iterate(*B, *X, 1, true, 0);
495 }
496 
497 // Solve phase
498 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
499 void ShiftedLaplacian<Scalar,LocalOrdinal,GlobalOrdinal,Node>::multigrid_apply(const RCP<Tpetra::MultiVector<SC,LO,GO,NO> > B,
500  RCP<Tpetra::MultiVector<SC,LO,GO,NO> >& X)
501 {
502  Teuchos::RCP< Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > XpetraX
503  = Teuchos::rcp( new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(X) );
504  Teuchos::RCP< Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > XpetraB
505  = Teuchos::rcp( new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(B) );
506  // Set left and right hand sides for Belos
507  Hierarchy_ -> Iterate(*XpetraB, *XpetraX, 1, true, 0);
508 }
509 
510 // Get most recent iteration count
511 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
513 {
514  int numiters = SolverManager_ -> getNumIters();
515  return numiters;
516 }
517 
518 // Get most recent solver tolerance achieved
519 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
521 {
522  double residual = SolverManager_ -> achievedTol();
523  return residual;
524 }
525 
526 }
527 
528 #define MUELU_SHIFTEDLAPLACIAN_SHORT
529 
530 #endif //if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA)
531 #endif // MUELU_SHIFTEDLAPLACIAN_DEF_HPP
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
This class specifies the default factory that should generate some data on a Level if the data does n...
Factory for building coarse grid matrices, when the matrix is of the form K+a*M. Useful when you want...
void multigrid_apply(const RCP< MultiVector > B, RCP< MultiVector > &X)
Namespace for MueLu classes and methods.
Factory for building tentative prolongator.
void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList)
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
Factory for coarsening a graph with uncoupled aggregation.
Factory for building restriction operators using a prolongator factory.
void setLevelShifts(std::vector< Scalar > levelshifts)
Always keep data, even accross run. This flag is set by Level::Keep(). This flag is propagated to coa...
static RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2NonConstTpetraCrs(RCP< Matrix > Op)
void setcoords(RCP< MultiVector > &Coords)
Wraps an existing MueLu::Hierarchy as a Tpetra::Operator, with an optional two-level correction...
int solve(const RCP< TMV > B, RCP< TMV > &X)
Factory for creating a graph base on a given matrix.
void setProblemMatrix(RCP< Matrix > &A)
Factory for building restriction operators.
Print all warning messages.
Factory for building coarse matrices.
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
Belos::LinearProblem< SC, TMV, OP > LinearProblem
Class that encapsulates Ifpack2 smoothers.
Factory for building Smoothed Aggregation prolongators.Input/output of SaPFactory
Factory for building uncoupled aggregates.
Belos::SolverFactory< SC, TMV, OP > SolverFactory
void setNullSpace(RCP< MultiVector > NullSpace)
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
void setPreconditioningMatrix(RCP< Matrix > &P)