MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RefMaxwellSmoother_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_REFMAXWELLSMOOTHER_DEF_HPP
11#define MUELU_REFMAXWELLSMOOTHER_DEF_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14
15#include <Teuchos_ParameterList.hpp>
16
17#include <Xpetra_CrsMatrix.hpp>
18#include <Xpetra_Matrix.hpp>
19#include <Xpetra_MultiVectorFactory.hpp>
20
22#include "MueLu_Level.hpp"
23#include "MueLu_Monitor.hpp"
24
25namespace MueLu {
26
27template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28RefMaxwellSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::RefMaxwellSmoother(const std::string type, const Teuchos::ParameterList& paramList)
29 : type_(type) {
30 const bool solverSupported = (type_ == "RefMaxwell");
31 this->declareConstructionOutcome(!solverSupported, "RefMaxwellSmoother does not provide the smoother '" + type_ + "'.");
32 if (solverSupported)
33 SetParameterList(paramList);
34}
35
36template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40
41template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43 const ParameterList& pL = this->GetParameterList();
44
45 int spaceNumber = 1;
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");
50 if (spaceNumber > 1)
51 this->Input(currentLevel, "Dk_2");
52 this->Input(currentLevel, "D0");
53 this->Input(currentLevel, "Mk_one");
54 if (spaceNumber > 1)
55 this->Input(currentLevel, "Mk_1_one");
56 this->Input(currentLevel, "M1_beta");
57 if (spaceNumber > 1)
58 this->Input(currentLevel, "M1_alpha");
59 this->Input(currentLevel, "invMk_1_invBeta");
60 if (spaceNumber > 1)
61 this->Input(currentLevel, "invMk_2_invAlpha");
62 this->Input(currentLevel, "Coordinates");
63}
64
65template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67 FactoryMonitor m(*this, "Setup Smoother", currentLevel);
68
69 SetupRefMaxwell(currentLevel);
71 this->GetOStream(Statistics1) << description() << std::endl;
72}
73
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>;
78
79 ParameterList params = this->GetParameterList();
80 int spaceNumber = 1;
81 if (params.isType<int>("refmaxwell: space number"))
82 spaceNumber = params.get<int>("refmaxwell: space number");
83
84 RCP<Matrix> SM, Dk_1, Dk_2, D0, Mk_one, Mk_1_one, M1_beta, M1_alpha, invMk_1_invBeta, invMk_2_invAlpha;
85
86 SM = currentLevel.Get<RCP<Matrix> >("A");
87 Dk_1 = currentLevel.Get<RCP<Matrix> >("Dk_1");
88 if (spaceNumber > 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");
92 if (spaceNumber > 1)
93 Mk_1_one = currentLevel.Get<RCP<Matrix> >("Mk_1_one");
94 M1_beta = currentLevel.Get<RCP<Matrix> >("M1_beta");
95 if (spaceNumber > 1)
96 M1_alpha = currentLevel.Get<RCP<Matrix> >("M1_alpha");
97 invMk_1_invBeta = currentLevel.Get<RCP<Matrix> >("invMk_1_invBeta");
98 if (spaceNumber > 1)
99 invMk_2_invAlpha = currentLevel.Get<RCP<Matrix> >("invMk_2_invAlpha");
100 auto coords = currentLevel.Get<RCP<RealValuedMultiVector> >("Coordinates");
101
103 Dk_1, Dk_2, D0,
104 M1_beta, M1_alpha,
105 Mk_one, Mk_1_one,
106 invMk_1_invBeta, invMk_2_invAlpha,
107 Teuchos::null, Teuchos::null, coords,
108 params));
109 std::ostringstream oss;
110 op_->describe(*Teuchos::fancyOStream(Teuchos::rcpFromRef(oss)));
111 cachedDescription_ = oss.str();
112}
113
114template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115void RefMaxwellSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
116 TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::RefMaxwellSmoother::Apply(): Setup() has not been called");
117
118 if (InitialGuessIsZero) {
119 X.putScalar(0.0);
120 }
121 op_->apply(B, X);
122}
123
124template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> > RefMaxwellSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
126 RCP<RefMaxwellSmoother> smoother = rcp(new RefMaxwellSmoother(*this));
127 smoother->SetParameterList(this->GetParameterList());
128 return smoother;
129}
130
131template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
133 std::ostringstream out;
135 out << op_->description();
136 out << std::endl
137 << std::endl;
138 // out << cachedDescription_;
139 } else {
140 out << "RefMaxwell";
141 }
142 return out.str();
143}
144
145template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
146void RefMaxwellSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream& out, const VerbLevel verbLevel) const {
148
149 if (verbLevel & Parameters1) {
150 out0 << "Parameter list: " << std::endl;
151 Teuchos::OSTab tab2(out);
152 out << this->GetParameterList();
153 }
154
155 if (verbLevel & External) {
156 if (op_ != Teuchos::null) {
157 Teuchos::OSTab tab2(out);
158 out << *op_ << std::endl
159 << std::endl;
160 }
161 }
162
163 // if (verbLevel & Debug) {
164 // out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
165 // << "-" << std::endl
166 // << "RCP<solver_>: " << tSolver_ << std::endl;
167 // }
168}
169
170template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
172 return Teuchos::OrdinalTraits<size_t>::invalid();
173}
174
175} // namespace MueLu
176
177#endif // MUELU_REFMAXWELLSMOOTHER_DEF_HPP
#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 &paramList)=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 &currentLevel)
void Setup(Level &currentLevel)
Set up the smoother.
void SetParameterList(const Teuchos::ParameterList &paramList)
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 &currentLevel) const
Input.
RefMaxwellSmoother(const std::string type, const Teuchos::ParameterList &paramList)
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)