53 #ifndef MUELU_SIMPLESMOOTHER_DEF_HPP_
54 #define MUELU_SIMPLESMOOTHER_DEF_HPP_
56 #include "Teuchos_ArrayViewDecl.hpp"
70 #include "MueLu_Utilities.hpp"
72 #include "MueLu_HierarchyUtils.hpp"
74 #include "MueLu_SubBlockAFactory.hpp"
77 #include "MueLu_SchurComplementFactory.hpp"
78 #include "MueLu_DirectSolver.hpp"
79 #include "MueLu_SmootherFactory.hpp"
80 #include "MueLu_FactoryManager.hpp"
84 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
86 : type_(
"SIMPLE"), A_(Teuchos::null)
99 SchurFact->SetFactory(
"A", this->
GetFactory(
"A"));
105 smoProtoSC->SetFactory(
"A", SchurFact);
110 schurFactManager->SetFactory(
"A", SchurFact);
111 schurFactManager->SetFactory(
"Smoother", SmooSCFact);
112 schurFactManager->SetIgnoreUserData(
true);
116 A00Fact->SetFactory(
"A",this->
GetFactory(
"A"));
120 std::string velpredictType;
122 smoProtoPredict->SetFactory(
"A", A00Fact);
126 velpredictFactManager->SetFactory(
"A", A00Fact);
127 velpredictFactManager->SetFactory(
"Smoother", SmooPredictFact);
128 velpredictFactManager->SetIgnoreUserData(
true);
135 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
138 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
143 validParamList->
set< Scalar > (
"Damping factor", 1.0,
"Damping/Scaling factor in SIMPLE");
144 validParamList->
set< LocalOrdinal > (
"Sweeps", 1,
"Number of SIMPLE sweeps (default = 1)");
145 validParamList->
set<
bool > (
"UseSIMPLEC",
false,
"Use SIMPLEC instead of SIMPLE (default = false)");
147 return validParamList;
151 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
155 size_t myPos = Teuchos::as<size_t>(pos);
157 if (myPos < FactManager_.size()) {
159 FactManager_.at(myPos) = FactManager;
160 }
else if( myPos == FactManager_.size()) {
162 FactManager_.push_back(FactManager);
165 *out <<
"Warning: cannot add new FactoryManager at proper position " << pos <<
". The FactoryManager is just appended to the end. Check this!" << std::endl;
168 FactManager_.push_back(FactManager);
174 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
176 AddFactoryManager(FactManager, 0);
179 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
181 AddFactoryManager(FactManager, 1);
184 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
186 currentLevel.
DeclareInput(
"A",this->GetFactory(
"A").
get());
188 TEUCHOS_TEST_FOR_EXCEPTION(FactManager_.size() != 2,
Exceptions::RuntimeError,
"MueLu::SimpleSmoother::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!");
191 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
192 for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
196 currentLevel.
DeclareInput(
"PreSmoother",(*it)->GetFactory(
"Smoother").get());
200 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
209 FactoryMonitor m(*
this,
"Setup blocked SIMPLE Smoother", currentLevel);
212 this->GetOStream(
Warnings0) <<
"MueLu::SimpleSmoother::Setup(): Setup() has already been called";
215 A_ = Factory::Get<RCP<Matrix> > (currentLevel,
"A");
221 rangeMapExtractor_ = bA->getRangeMapExtractor();
222 domainMapExtractor_ = bA->getDomainMapExtractor();
225 F_ = bA->getMatrix(0, 0);
226 G_ = bA->getMatrix(0, 1);
227 D_ = bA->getMatrix(1, 0);
228 Z_ = bA->getMatrix(1, 1);
231 bool bSIMPLEC = pL.
get<
bool>(
"UseSIMPLEC");
235 RCP<Vector> diagFVector = VectorFactory::Build(F_->getRowMap());
237 F_->getLocalDiagCopy(*diagFVector);
259 if(bdiagFinv.
is_null() ==
false && bdiagFinv->getBlockedMap()->getNumMaps() == 1) {
260 RCP<Vector> nestedVec = bdiagFinv->getMultiVector(0,bdiagFinv->getBlockedMap()->getThyraMode())->getVectorNonConst(0);
261 diagFinv_.
swap(nestedVec);
280 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
290 if(rcpBDebugB.is_null() ==
false) {
295 if(rcpBDebugX.
is_null() ==
false) {
308 LocalOrdinal nSweeps = pL.
get<LocalOrdinal>(
"Sweeps");
309 Scalar omega = pL.
get<Scalar>(
"Damping factor");
312 bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
313 bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
322 bool bCopyResultX =
false;
328 if(bX.is_null() ==
true) {
334 if(bB.is_null() ==
true) {
340 bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
341 bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
345 if(rbA.is_null() ==
false) {
350 if(bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps()) {
356 if(bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps()) {
368 RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
374 RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
380 RCP<MultiVector> xhat = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
387 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
389 residual->update(one,*rcpB,zero);
390 if(InitialGuessIsZero ==
false || run > 0)
395 xtilde1->putScalar(zero);
396 xtilde2->putScalar(zero);
397 velPredictSmoo_->Apply(*xtilde1,*r1);
401 RCP<MultiVector> schurCompRHS = rangeMapExtractor_->getVector(1, rcpB->getNumVectors(), bRangeThyraMode);
402 D_->apply(*xtilde1,*schurCompRHS);
404 schurCompRHS->update(one,*r2,-one);
408 schurCompSmoo_->Apply(*xtilde2,*schurCompRHS);
412 xhat2->update(omega,*xtilde2,zero);
415 RCP<MultiVector> xhat1_temp = domainMapExtractor_->getVector(0, rcpX->getNumVectors(), bDomainThyraMode);
416 G_->apply(*xhat2,*xhat1_temp);
418 xhat1->elementWiseMultiply(one,*diagFinv_,*xhat1_temp,zero);
419 xhat1->update(one,*xtilde1,-one);
421 rcpX->update(one,*bxhat,one);
424 if (bCopyResultX ==
true) {
426 X.update(one, *Xmerged, zero);
431 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
437 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
439 std::ostringstream out;
441 out <<
"{type = " << type_ <<
"}";
445 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
450 out0 <<
"Prec. type: " << type_ << std::endl;
453 if (verbLevel &
Debug) {
458 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
Important warning messages (one line)
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
virtual const Teuchos::ParameterList & GetParameterList() const
Exception indicating invalid cast attempted.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos)
Add a factory manager at a specific position.
This class specifies the default factory that should generate some data on a Level if the data does n...
SimpleSmoother()
Constructor.
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)
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void Setup(Level ¤tLevel)
Setup routine.
void swap(RCP< T > &r_ptr)
Print additional debugging information.
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
void DeclareInput(Level ¤tLevel) const
Input.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
virtual ~SimpleSmoother()
Destructor.
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
RCP< SmootherPrototype > Copy() const
virtual const RCP< const FactoryBase > GetFactory(const std::string &varName) const =0
Get.
RCP< const ParameterList > GetValidParameterList() const
Input.
std::string description() const
Return a simple one-line description of this object.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Factory for building a thresholded operator.
bool IsSetup() const
Get the state of a smoother prototype.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
void Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
Factory for building the Schur Complement for a 2x2 block matrix.
virtual const RCP< const FactoryBase > GetFactory(const std::string &varName) const =0
static Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetInverse(Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero())
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()
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.
SIMPLE smoother for 2x2 block matrices.
std::string toString(const T &t)