53 #ifndef MUELU_UZAWASMOOTHER_DEF_HPP_
54 #define MUELU_UZAWASMOOTHER_DEF_HPP_
56 #include "Teuchos_ArrayViewDecl.hpp"
61 #include <Xpetra_Matrix.hpp>
62 #include <Xpetra_BlockedCrsMatrix.hpp>
63 #include <Xpetra_MultiVectorFactory.hpp>
64 #include <Xpetra_VectorFactory.hpp>
65 #include <Xpetra_ReorderedBlockedCrsMatrix.hpp>
69 #include "MueLu_Utilities.hpp"
71 #include "MueLu_HierarchyUtils.hpp"
73 #include "MueLu_SubBlockAFactory.hpp"
76 #include "MueLu_SchurComplementFactory.hpp"
77 #include "MueLu_DirectSolver.hpp"
78 #include "MueLu_SmootherFactory.hpp"
79 #include "MueLu_FactoryManager.hpp"
83 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
85 : type_(
"Uzawa"), A_(Teuchos::null)
89 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
92 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
97 validParamList->
set<
Scalar > (
"Damping factor", 1.0,
"Damping/Scaling factor in SIMPLE");
98 validParamList->
set<
LocalOrdinal > (
"Sweeps", 1,
"Number of SIMPLE sweeps (default = 1)");
100 return validParamList;
104 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
108 size_t myPos = Teuchos::as<size_t>(pos);
110 if (myPos < FactManager_.size()) {
112 FactManager_.at(myPos) = FactManager;
113 }
else if( myPos == FactManager_.size()) {
115 FactManager_.push_back(FactManager);
118 *out <<
"Warning: cannot add new FactoryManager at proper position " << pos <<
". The FactoryManager is just appended to the end. Check this!" << std::endl;
121 FactManager_.push_back(FactManager);
127 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
129 AddFactoryManager(FactManager, 0);
132 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 AddFactoryManager(FactManager, 1);
137 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
139 currentLevel.
DeclareInput(
"A",this->GetFactory(
"A").
get());
141 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!");
144 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
145 for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
149 currentLevel.
DeclareInput(
"PreSmoother",(*it)->GetFactory(
"Smoother").get());
153 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
156 FactoryMonitor m(*
this,
"Setup blocked Uzawa Smoother", currentLevel);
159 this->GetOStream(
Warnings0) <<
"MueLu::UzawaSmoother::Setup(): Setup() has already been called";
162 A_ = Factory::Get<RCP<Matrix> > (currentLevel,
"A");
168 rangeMapExtractor_ = bA->getRangeMapExtractor();
169 domainMapExtractor_ = bA->getDomainMapExtractor();
172 F_ = bA->getMatrix(0, 0);
173 G_ = bA->getMatrix(0, 1);
174 D_ = bA->getMatrix(1, 0);
175 Z_ = bA->getMatrix(1, 1);
193 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
202 bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
203 bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
215 bool bCopyResultX =
false;
221 if(bX.is_null() ==
true) {
227 if(bB.is_null() ==
true) {
233 bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
234 bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
238 if(rbA.is_null() ==
false) {
243 if(bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps()) {
246 buildReorderedBlockedMultiVector(brm, bX);
249 if(bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps()) {
252 buildReorderedBlockedMultiVector(brm, bB);
261 RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
267 RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
275 residual->update(one,*rcpB,zero);
280 bxtilde->putScalar(zero);
281 velPredictSmoo_->Apply(*xtilde1,*r1);
285 RCP<MultiVector> schurCompRHS = MultiVectorFactory::Build(r2->getMap(),rcpB->getNumVectors());
286 D_->apply(*xtilde1,*schurCompRHS);
287 schurCompRHS->update(one,*r2,-one);
291 schurCompSmoo_->Apply(*xtilde2,*schurCompRHS);
293 rcpX->update(omega,*bxtilde,one);
296 if (bCopyResultX ==
true) {
298 X.update(one, *Xmerged, zero);
302 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
308 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
310 std::ostringstream out;
312 out <<
"{type = " << type_ <<
"}";
316 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
321 out0 <<
"Prec. type: " << type_ << std::endl;
324 if (verbLevel &
Debug) {
329 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)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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.
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.
UzawaSmoother()
Constructor.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
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)
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.
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()'.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual ~UzawaSmoother()
Destructor.
virtual std::string description() const
Return a simple one-line description of this object.
std::string toString(const T &t)