53 #ifndef MUELU_SIMPLESMOOTHER_DEF_HPP_ 54 #define MUELU_SIMPLESMOOTHER_DEF_HPP_ 56 #include "Teuchos_ArrayViewDecl.hpp" 57 #include "Teuchos_ScalarTraits.hpp" 61 #include <Xpetra_Matrix.hpp> 62 #include <Xpetra_CrsMatrixWrap.hpp> 63 #include <Xpetra_BlockedCrsMatrix.hpp> 64 #include <Xpetra_MultiVectorFactory.hpp> 65 #include <Xpetra_VectorFactory.hpp> 66 #include <Xpetra_ReorderedBlockedCrsMatrix.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)
97 SchurFact->SetParameter(
"omega",Teuchos::ParameterEntry(omega));
98 SchurFact->SetParameter(
"lumping",Teuchos::ParameterEntry(SIMPLEC));
99 SchurFact->SetFactory(
"A", this->
GetFactory(
"A"));
102 Teuchos::ParameterList SCparams;
104 RCP<SmootherPrototype> smoProtoSC = rcp(
new DirectSolver(SCtype,SCparams) );
105 smoProtoSC->SetFactory(
"A", SchurFact);
107 RCP<SmootherFactory> SmooSCFact = rcp(
new SmootherFactory(smoProtoSC) );
110 schurFactManager->SetFactory(
"A", SchurFact);
111 schurFactManager->SetFactory(
"Smoother", SmooSCFact);
112 schurFactManager->SetIgnoreUserData(
true);
116 A00Fact->SetFactory(
"A",this->
GetFactory(
"A"));
117 A00Fact->SetParameter(
"block row",ParameterEntry(0));
118 A00Fact->SetParameter(
"block col",ParameterEntry(0));
119 Teuchos::ParameterList velpredictParams;
120 std::string velpredictType;
121 RCP<SmootherPrototype> smoProtoPredict = rcp(
new DirectSolver(velpredictType,velpredictParams) );
122 smoProtoPredict->SetFactory(
"A", A00Fact);
123 RCP<SmootherFactory> SmooPredictFact = rcp(
new SmootherFactory(smoProtoPredict) );
125 RCP<FactoryManager> velpredictFactManager = rcp(
new FactoryManager());
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>
140 RCP<ParameterList> validParamList = rcp(
new ParameterList());
142 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
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>
153 TEUCHOS_TEST_FOR_EXCEPTION(pos < 0,
Exceptions::RuntimeError,
"MueLu::SimpleSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
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);
164 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
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");
217 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
218 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null,
Exceptions::BadCast,
"MueLu::SimpleSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
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);
258 RCP<BlockedVector> bdiagFinv = Teuchos::rcp_dynamic_cast<BlockedVector>(diagFinv_);
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);
267 RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
269 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
272 RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
274 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
280 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
286 RCP<MultiVector> rcpDebugX = Teuchos::rcpFromRef(X);
287 RCP<const MultiVector> rcpDebugB = Teuchos::rcpFromRef(B);
288 RCP<BlockedMultiVector> rcpBDebugX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpDebugX);
289 RCP<const BlockedMultiVector> rcpBDebugB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpDebugB);
290 if(rcpBDebugB.is_null() ==
false) {
295 if(rcpBDebugX.is_null() ==
false) {
302 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
304 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
312 bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
313 bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
318 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
319 RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
322 bool bCopyResultX =
false;
323 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
325 RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
326 RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
328 if(bX.is_null() ==
true) {
329 RCP<MultiVector> test = Teuchos::rcp(
new BlockedMultiVector(bA->getBlockedDomainMap(),rcpX));
334 if(bB.is_null() ==
true) {
335 RCP<const MultiVector> test = Teuchos::rcp(
new BlockedMultiVector(bA->getBlockedRangeMap(),rcpB));
340 bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
341 bB = Teuchos::rcp_dynamic_cast<
const BlockedMultiVector>(rcpB);
344 RCP<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > rbA = Teuchos::rcp_dynamic_cast<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA);
345 if(rbA.is_null() ==
false) {
347 Teuchos::RCP<const Xpetra::BlockReorderManager > brm = rbA->getBlockReorderManager();
350 if(bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps()) {
352 Teuchos::RCP<MultiVector> test =
353 buildReorderedBlockedMultiVector(brm, bX);
356 if(bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps()) {
358 Teuchos::RCP<const MultiVector> test =
359 buildReorderedBlockedMultiVector(brm, bB);
368 RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
369 RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
370 Teuchos::RCP<MultiVector> r1 = bresidual->getMultiVector(0,bRangeThyraMode);
371 Teuchos::RCP<MultiVector> r2 = bresidual->getMultiVector(1,bRangeThyraMode);
374 RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
375 RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
376 RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0,bDomainThyraMode);
377 RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1,bDomainThyraMode);
380 RCP<MultiVector> xhat = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
381 RCP<BlockedMultiVector> bxhat = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xhat);
382 RCP<MultiVector> xhat1 = bxhat->getMultiVector(0,bDomainThyraMode);
383 RCP<MultiVector> xhat2 = bxhat->getMultiVector(1,bDomainThyraMode);
389 residual->update(one,*rcpB,zero);
390 if(InitialGuessIsZero ==
false || run > 0)
391 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
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) {
425 RCP<MultiVector> Xmerged = bX->Merge();
426 X.update(one, *Xmerged, zero);
431 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
432 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, 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>
461 return Teuchos::OrdinalTraits<size_t>::invalid();
Important warning messages (one line)
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
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...
MueLu::DefaultLocalOrdinal LocalOrdinal
SimpleSmoother()
Constructor.
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.
void Setup(Level ¤tLevel)
Setup routine.
Print additional debugging information.
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Namespace for MueLu classes and methods.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
virtual ~SimpleSmoother()
Destructor.
RCP< SmootherPrototype > Copy() const
Class that encapsulates direct solvers. Autoselection of AmesosSmoother or Amesos2Smoother according ...
MueLu::DefaultScalar Scalar
void DeclareInput(Level ¤tLevel) const
Input.
Class that holds all level-specific information.
virtual const Teuchos::ParameterList & GetParameterList() const
Factory for building a thresholded operator.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Factory for building the Schur Complement for a 2x2 block matrix.
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())
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
RCP< const ParameterList > GetValidParameterList() const
Input.
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory()
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()
std::string description() const
Return a simple one-line description of this object.
virtual std::string description() const
Return a simple one-line description of this object.
SIMPLE smoother for 2x2 block matrices.
void Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.