53 #ifndef MUELU_UZAWASMOOTHER_DEF_HPP_ 54 #define MUELU_UZAWASMOOTHER_DEF_HPP_ 56 #include "Teuchos_ArrayViewDecl.hpp" 57 #include "Teuchos_ScalarTraits.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>
94 RCP<ParameterList> validParamList = rcp(
new ParameterList());
96 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
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>
106 TEUCHOS_TEST_FOR_EXCEPTION(pos < 0,
Exceptions::RuntimeError,
"MueLu::UzawaSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
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);
117 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
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");
164 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
165 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null,
Exceptions::BadCast,
"MueLu::UzawaSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
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);
180 RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
182 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
185 RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
187 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
193 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
198 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
200 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
202 bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
203 bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
211 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
212 RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
215 bool bCopyResultX =
false;
216 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
218 RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
219 RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
221 if(bX.is_null() ==
true) {
222 RCP<MultiVector> test = Teuchos::rcp(
new BlockedMultiVector(bA->getBlockedDomainMap(),rcpX));
227 if(bB.is_null() ==
true) {
228 RCP<const MultiVector> test = Teuchos::rcp(
new BlockedMultiVector(bA->getBlockedRangeMap(),rcpB));
233 bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
234 bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
237 RCP<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > rbA = Teuchos::rcp_dynamic_cast<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA);
238 if(rbA.is_null() ==
false) {
240 Teuchos::RCP<const Xpetra::BlockReorderManager > brm = rbA->getBlockReorderManager();
243 if(bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps()) {
245 Teuchos::RCP<MultiVector> test =
246 buildReorderedBlockedMultiVector(brm, bX);
249 if(bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps()) {
251 Teuchos::RCP<const MultiVector> test =
252 buildReorderedBlockedMultiVector(brm, bB);
261 RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
262 RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
263 Teuchos::RCP<MultiVector> r1 = bresidual->getMultiVector(0,bRangeThyraMode);
264 Teuchos::RCP<MultiVector> r2 = bresidual->getMultiVector(1,bRangeThyraMode);
267 RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
268 RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
269 RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0,bDomainThyraMode);
270 RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1,bDomainThyraMode);
275 residual->update(one,*rcpB,zero);
276 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
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) {
297 RCP<MultiVector> Xmerged = bX->Merge();
298 X.update(one, *Xmerged, zero);
302 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
303 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, 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>
332 return Teuchos::OrdinalTraits<size_t>::invalid();
Important warning messages (one line)
std::string description() const
Return a simple one-line description of this object.
Exception indicating invalid cast attempted.
void Setup(Level ¤tLevel)
Setup routine.
MueLu::DefaultLocalOrdinal LocalOrdinal
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
bool IsSetup() const
Get the state of a smoother prototype.
Timer to be used in factories. Similar to Monitor but with additional timers.
Print additional debugging information.
void DeclareInput(Level ¤tLevel) const
Input.
Namespace for MueLu classes and methods.
UzawaSmoother()
Constructor.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
RCP< const ParameterList > GetValidParameterList() const
Input.
MueLu::DefaultScalar Scalar
Class that holds all level-specific information.
Block triangle Uzawa smoother for 2x2 block matrices.
virtual const Teuchos::ParameterList & GetParameterList() const
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
RCP< SmootherPrototype > Copy() const
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
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 Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual ~UzawaSmoother()
Destructor.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
virtual std::string description() const
Return a simple one-line description of this object.