10 #ifndef MUELU_UZAWASMOOTHER_DEF_HPP_
11 #define MUELU_UZAWASMOOTHER_DEF_HPP_
13 #include <Teuchos_ArrayViewDecl.hpp>
20 #include <Xpetra_MultiVectorFactory.hpp>
27 #include "MueLu_HierarchyUtils.hpp"
30 #include "MueLu_FactoryManager.hpp"
34 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
39 validParamList->
set<
Scalar>(
"Damping factor", 1.0,
"Damping/Scaling factor in SIMPLE");
40 validParamList->
set<
LocalOrdinal>(
"Sweeps", 1,
"Number of SIMPLE sweeps (default = 1)");
42 return validParamList;
45 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
49 size_t myPos = Teuchos::as<size_t>(pos);
51 if (myPos < FactManager_.size()) {
53 FactManager_.at(myPos) = FactManager;
54 }
else if (myPos == FactManager_.size()) {
56 FactManager_.push_back(FactManager);
59 *out <<
"Warning: cannot add new FactoryManager at proper position " << pos <<
". The FactoryManager is just appended to the end. Check this!" << std::endl;
62 FactManager_.push_back(FactManager);
66 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 AddFactoryManager(FactManager, 0);
71 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 AddFactoryManager(FactManager, 1);
76 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
78 currentLevel.
DeclareInput(
"A", this->GetFactory(
"A").
get());
80 TEUCHOS_TEST_FOR_EXCEPTION(FactManager_.size() != 2,
Exceptions::RuntimeError,
"MueLu::UzawaSmoother::DeclareInput: You have to declare two FactoryManagers with a \"Smoother\" object: One for predicting the primary variable and one for the SchurComplement system. The smoother for the SchurComplement system needs a SchurComplementFactory as input for variable \"A\". make sure that you use the same proper damping factors for omega both in the SchurComplementFactory and in the SIMPLE smoother!");
83 std::vector<Teuchos::RCP<const FactoryManagerBase>>::const_iterator it;
84 for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
88 currentLevel.
DeclareInput(
"PreSmoother", (*it)->GetFactory(
"Smoother").get());
92 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
94 FactoryMonitor m(*
this,
"Setup blocked Uzawa Smoother", currentLevel);
97 this->GetOStream(
Warnings0) <<
"MueLu::UzawaSmoother::Setup(): Setup() has already been called";
100 A_ = Factory::Get<RCP<Matrix>>(currentLevel,
"A");
106 rangeMapExtractor_ = bA->getRangeMapExtractor();
107 domainMapExtractor_ = bA->getDomainMapExtractor();
110 F_ = bA->getMatrix(0, 0);
111 G_ = bA->getMatrix(0, 1);
112 D_ = bA->getMatrix(1, 0);
113 Z_ = bA->getMatrix(1, 1);
131 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
139 bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
140 bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
152 bool bCopyResultX =
false;
153 bool bReorderX =
false;
154 bool bReorderB =
false;
164 if (rbA.is_null() ==
false) {
169 if (bX.is_null() ==
true) {
171 rcpX.swap(vectorWithBlockedMap);
177 if (bB.is_null() ==
true) {
179 rcpB.
swap(vectorWithBlockedMap);
184 if (bX.is_null() ==
true) {
186 rcpX.swap(vectorWithBlockedMap);
190 if (bB.is_null() ==
true) {
192 rcpB.
swap(vectorWithBlockedMap);
197 bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
198 bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
201 if (rbA.is_null() ==
false) {
205 if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
208 rcpX.swap(reorderedBlockedVector);
211 if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
214 rcpB.
swap(reorderedBlockedVector);
228 RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
236 residual->update(one, *rcpB, zero);
241 bxtilde->putScalar(zero);
242 velPredictSmoo_->Apply(*xtilde1, *r1);
246 RCP<MultiVector> schurCompRHS = MultiVectorFactory::Build(r2->getMap(), rcpB->getNumVectors());
247 D_->apply(*xtilde1, *schurCompRHS);
248 schurCompRHS->update(one, *r2, -one);
252 schurCompSmoo_->Apply(*xtilde2, *schurCompRHS);
254 rcpX->update(omega, *bxtilde, one);
257 if (bCopyResultX ==
true) {
259 X.update(one, *Xmerged, zero);
263 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
269 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
271 std::ostringstream out;
273 out <<
"{type = " << type_ <<
"}";
277 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
282 out0 <<
"Prec. type: " << type_ << std::endl;
285 if (verbLevel &
Debug) {
290 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
Important warning messages (one line)
virtual const Teuchos::ParameterList & GetParameterList() const
Exception indicating invalid cast attempted.
void Setup(Level ¤tLevel)
Setup routine.
MueLu::DefaultLocalOrdinal LocalOrdinal
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void swap(RCP< T > &r_ptr)
Print additional debugging information.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void DeclareInput(Level ¤tLevel) const
Input.
void Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
RCP< const ParameterList > GetValidParameterList() const
Input.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
Set factory manager for internal velocity prediction.
virtual const RCP< const FactoryBase > GetFactory(const std::string &varName) const =0
Get.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
Class that holds all level-specific information.
Block triangle Uzawa smoother for 2x2 block matrices.
std::string description() const
Return a simple one-line description of this object.
RCP< SmootherPrototype > Copy() const
bool IsSetup() const
Get the state of a smoother prototype.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
Set factory manager for internal SchurComplement handling.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos)
Add a factory manager at a specific position.
void residual(const Operator< SC, LO, GO, NO > &Aop, const MultiVector< SC, LO, GO, NO > &X_in, const MultiVector< SC, LO, GO, NO > &B_in, MultiVector< SC, LO, GO, NO > &R_in)
Exception throws to report errors in the internal logical of the program.
#define MUELU_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
An exception safe way to call the method 'Level::SetFactoryManager()'.
Teuchos::RCP< const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > buildReorderedBlockedMultiVector(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bvec)
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual std::string description() const
Return a simple one-line description of this object.
std::string toString(const T &t)