10#ifndef MUELU_REFMAXWELLSMOOTHER_DEF_HPP
11#define MUELU_REFMAXWELLSMOOTHER_DEF_HPP
15#include <Teuchos_ParameterList.hpp>
17#include <Xpetra_CrsMatrix.hpp>
18#include <Xpetra_Matrix.hpp>
19#include <Xpetra_MultiVectorFactory.hpp>
27template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
30 const bool solverSupported = (
type_ ==
"RefMaxwell");
36template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
41template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
43 const ParameterList& pL = this->GetParameterList();
46 if (pL.isType<
int>(
"refmaxwell: space number"))
47 spaceNumber = pL.get<
int>(
"refmaxwell: space number");
48 this->Input(currentLevel,
"A");
49 this->Input(currentLevel,
"Dk_1");
51 this->Input(currentLevel,
"Dk_2");
52 this->Input(currentLevel,
"D0");
53 this->Input(currentLevel,
"Mk_one");
55 this->Input(currentLevel,
"Mk_1_one");
56 this->Input(currentLevel,
"M1_beta");
58 this->Input(currentLevel,
"M1_alpha");
59 this->Input(currentLevel,
"invMk_1_invBeta");
61 this->Input(currentLevel,
"invMk_2_invAlpha");
62 this->Input(currentLevel,
"Coordinates");
65template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 SetupRefMaxwell(currentLevel);
71 this->GetOStream(
Statistics1) << description() << std::endl;
74template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 using coordinateType =
typename Teuchos::ScalarTraits<Scalar>::coordinateType;
77 using RealValuedMultiVector = Xpetra::MultiVector<coordinateType, LO, GO, NO>;
79 ParameterList params = this->GetParameterList();
81 if (params.isType<
int>(
"refmaxwell: space number"))
82 spaceNumber = params.get<
int>(
"refmaxwell: space number");
84 RCP<Matrix> SM, Dk_1, Dk_2, D0, Mk_one, Mk_1_one, M1_beta, M1_alpha, invMk_1_invBeta, invMk_2_invAlpha;
86 SM = currentLevel.
Get<RCP<Matrix> >(
"A");
87 Dk_1 = currentLevel.
Get<RCP<Matrix> >(
"Dk_1");
89 Dk_2 = currentLevel.
Get<RCP<Matrix> >(
"Dk_2");
90 D0 = currentLevel.
Get<RCP<Matrix> >(
"D0");
91 Mk_one = currentLevel.
Get<RCP<Matrix> >(
"Mk_one");
93 Mk_1_one = currentLevel.
Get<RCP<Matrix> >(
"Mk_1_one");
94 M1_beta = currentLevel.
Get<RCP<Matrix> >(
"M1_beta");
96 M1_alpha = currentLevel.
Get<RCP<Matrix> >(
"M1_alpha");
97 invMk_1_invBeta = currentLevel.
Get<RCP<Matrix> >(
"invMk_1_invBeta");
99 invMk_2_invAlpha = currentLevel.
Get<RCP<Matrix> >(
"invMk_2_invAlpha");
100 auto coords = currentLevel.
Get<RCP<RealValuedMultiVector> >(
"Coordinates");
106 invMk_1_invBeta, invMk_2_invAlpha,
107 Teuchos::null, Teuchos::null, coords,
109 std::ostringstream oss;
110 op_->describe(*Teuchos::fancyOStream(Teuchos::rcpFromRef(oss)));
111 cachedDescription_ = oss.str();
114template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 if (InitialGuessIsZero) {
124template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 smoother->SetParameterList(this->GetParameterList());
131template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
133 std::ostringstream out;
135 out << op_->description();
145template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
150 out0 <<
"Parameter list: " << std::endl;
151 Teuchos::OSTab tab2(out);
152 out << this->GetParameterList();
156 if (op_ != Teuchos::null) {
157 Teuchos::OSTab tab2(out);
158 out << *op_ << std::endl
170template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
172 return Teuchos::OrdinalTraits<size_t>::invalid();
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual void SetParameterList(const Teuchos::ParameterList ¶mList)=0
Set parameters from a parameter list and return with default values.
Class that encapsulates Operator smoothers.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void SetupRefMaxwell(Level ¤tLevel)
void Setup(Level ¤tLevel)
Set up the smoother.
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void DeclareInput(Level ¤tLevel) const
Input.
RefMaxwellSmoother(const std::string type, const Teuchos::ParameterList ¶mList)
Constructor.
std::string description() const
Return a simple one-line description of this object.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
RCP< SmootherPrototype > Copy() const
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form.
void declareConstructionOutcome(bool fail, std::string msg)
bool IsSetup() const
Get the state of a smoother prototype.
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.
@ External
Print external lib objects.
@ Parameters1
Print class parameters (more parameters, more verbose)