46 #ifndef MUELU_SHIFTEDLAPLACIAN_DECL_HPP 47 #define MUELU_SHIFTEDLAPLACIAN_DECL_HPP 50 #include <Xpetra_Matrix_fwd.hpp> 51 #include <Xpetra_VectorFactory_fwd.hpp> 52 #include <Xpetra_MultiVectorFactory_fwd.hpp> 53 #include <Xpetra_TpetraMultiVector.hpp> 59 #if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA) 77 #include <MueLu_ShiftedLaplacianOperator.hpp> 86 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT 87 #include <BelosConfigDefs.hpp> 88 #include <BelosLinearProblem.hpp> 89 #include <BelosSolverFactory.hpp> 90 #include <BelosTpetraAdapter.hpp> 109 #undef MUELU_SHIFTEDLAPLACIAN_SHORT 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;
178 void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList);
186 void setstiff(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraK);
188 void setmass(RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> >& TpetraM);
189 void setcoords(RCP<MultiVector>& Coords);
215 int solve(
const RCP<TMV> B, RCP<TMV>& X);
217 RCP<MultiVector>& X);
219 RCP<Tpetra::MultiVector<SC,LO,GO,NO> >& X);
221 typename Teuchos::ScalarTraits<Scalar>::magnitudeType
GetResidual();
299 RCP< MueLu::ShiftedLaplacianOperator<SC,LO,GO,NO> >
MueLuOp_;
302 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT 304 RCP<LinearProblem> LinearProblem_;
305 RCP<SolverManager> SolverManager_;
306 RCP<SolverFactory> SolverFactory_;
307 RCP<Teuchos::ParameterList> BelosList_;
314 #define MUELU_SHIFTEDLAPLACIAN_SHORT 316 #endif //if defined(HAVE_MUELU_IFPACK2) and defined(HAVE_MUELU_TPETRA) 318 #endif // MUELU_SHIFTEDLAPLACIAN_DECL_HPP
ShiftedLaplacian()
Constructors.
Tpetra::CombineMode schwarz_combinemode_
RCP< SmootherFactory > coarsestSmooFact_
MueLu::DefaultLocalOrdinal LocalOrdinal
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
std::string ilu_normtype_
RCP< TransPFactory > TransPfact_
void multigrid_apply(const RCP< MultiVector > B, RCP< MultiVector > &X)
RCP< RAPShiftFactory > Acshift_
void setmass(RCP< Matrix > &M)
std::string ilu_drop_rule_
RCP< CoarseMapFactory > CoarseMapfact_
RCP< SmootherFactory > smooFact_
Namespace for MueLu classes and methods.
double ilu_diagpivotthresh_
RCP< CoalesceDropFactory > Dropfact_
void resetLinearProblem()
std::string schwarz_ordermethod_
void setParameters(Teuchos::RCP< Teuchos::ParameterList > paramList)
Teuchos::ParameterList precList_
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< PgPFactory > PgPfact_
RCP< MultiVector > Coords_
virtual ~ShiftedLaplacian()
RCP< MueLu::ShiftedLaplacianOperator< SC, LO, GO, NO > > MueLuOp_
Teuchos::ParameterList coarsestSmooList_
void setcoords(RCP< MultiVector > &Coords)
RCP< GenericRFactory > Rfact_
RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > TpetraA_
RCP< MultiVector > NullSpace_
RCP< UncoupledAggregationFactory > UCaggfact_
Tpetra::Operator< SC, LO, GO, NO > OP
RCP< RAPFactory > Acfact_
int solve(const RCP< TMV > B, RCP< TMV > &X)
Base class for MueLu classes.
RCP< SmootherPrototype > coarsestSmooProto_
RCP< Hierarchy > Hierarchy_
RCP< CoupledAggregationFactory > Aggfact_
Tpetra::Vector< SC, LO, GO, NO > TVEC
void setProblemMatrix(RCP< Matrix > &A)
void setstiff(RCP< Matrix > &K)
int krylov_preconditioner_
RCP< TentativePFactory > TentPfact_
Shifted Laplacian Helmholtz solver.
RCP< AmalgamationFactory > Amalgfact_
std::vector< SC > levelshifts_
RCP< FactoryManager > Manager_
void setNullSpace(RCP< MultiVector > NullSpace)
std::string ilu_milutype_
void setPreconditioningMatrix(RCP< Matrix > &P)