53 #ifndef MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_ 54 #define MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_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> 69 #include "MueLu_Utilities.hpp" 71 #include "MueLu_HierarchyHelpers.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_(
"IndefiniteBlockDiagonalSmoother"), 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 (must be a 2x2 block matrix)");
97 validParamList->set< Scalar > (
"Damping factor", 1.0,
"Damping/Scaling factor");
98 validParamList->set< LocalOrdinal > (
"Sweeps", 1,
"Number of SIMPLE sweeps (default = 1)");
101 return validParamList;
105 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
107 TEUCHOS_TEST_FOR_EXCEPTION(pos < 0,
Exceptions::RuntimeError,
"MueLu::IndefBlockedDiagonalSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
109 size_t myPos = Teuchos::as<size_t>(pos);
118 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
119 *out <<
"Warning: cannot add new FactoryManager at proper position " << pos <<
". The FactoryManager is just appended to the end. Check this!" << std::endl;
126 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
130 TEUCHOS_TEST_FOR_EXCEPTION(
FactManager_.size() != 2,
Exceptions::RuntimeError,
"MueLu::IndefBlockedDiagonalSmoother::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\"!");
133 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
138 currentLevel.
DeclareInput(
"PreSmoother",(*it)->GetFactory(
"Smoother").get());
143 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
145 FactoryMonitor m(*
this,
"Setup for indefinite blocked diagonal smoother", currentLevel);
148 this->
GetOStream(
Warnings0) <<
"MueLu::IndefBlockedDiagonalSmoother::Setup(): Setup() has already been called";
151 A_ = Factory::Get<RCP<Matrix> > (currentLevel,
"A");
153 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
154 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null,
Exceptions::BadCast,
"MueLu::IndefBlockedDiagonalSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
161 Teuchos::RCP<CrsMatrix> A00 = bA->getMatrix(0, 0);
162 Teuchos::RCP<CrsMatrix> A01 = bA->getMatrix(0, 1);
163 Teuchos::RCP<CrsMatrix> A10 = bA->getMatrix(1, 0);
164 Teuchos::RCP<CrsMatrix> A11 = bA->getMatrix(1, 1);
166 Teuchos::RCP<CrsMatrixWrap> Op00 = Teuchos::rcp(
new CrsMatrixWrap(A00));
167 Teuchos::RCP<CrsMatrixWrap> Op01 = Teuchos::rcp(
new CrsMatrixWrap(A01));
168 Teuchos::RCP<CrsMatrixWrap> Op10 = Teuchos::rcp(
new CrsMatrixWrap(A10));
169 Teuchos::RCP<CrsMatrixWrap> Op11 = Teuchos::rcp(
new CrsMatrixWrap(A11));
171 F_ = Teuchos::rcp_dynamic_cast<Matrix>(Op00);
174 Z_ = Teuchos::rcp_dynamic_cast<Matrix>(Op11);
177 F_->CreateView(
"stridedMaps", bA->getRangeMap(0), bA->getDomainMap(0));
180 Z_->CreateView(
"stridedMaps", bA->getRangeMap(1), bA->getDomainMap(1));
211 RCP<const FactoryManagerBase> velpredictFactManager =
FactManager_.at(0);
213 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
216 RCP<const FactoryManagerBase> schurFactManager =
FactManager_.at(1);
218 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
223 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
228 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
230 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
234 LocalOrdinal nSweeps = pL.get<LocalOrdinal>(
"Sweeps");
235 Scalar omega = pL.get<Scalar>(
"Damping factor");
238 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
242 RCP<MultiVector> residual = MultiVectorFactory::Build(B.getMap(), B.getNumVectors());
245 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
247 residual->update(one,B,zero);
248 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
256 RCP<MultiVector> xtilde1 = MultiVectorFactory::Build(
F_->getRowMap(),1);
257 xtilde1->putScalar(zero);
262 RCP<MultiVector> xtilde2 = MultiVectorFactory::Build(
Z_->getRowMap(),1);
263 xtilde2->putScalar(zero);
272 x1->update(omega,*xtilde1,one);
273 x2->update(omega,*xtilde2,one);
283 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
284 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
289 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
291 std::ostringstream out;
293 out <<
"{type = " <<
type_ <<
"}";
297 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
302 out0 <<
"Prec. type: " <<
type_ << std::endl;
305 if (verbLevel &
Debug) {
Important warning messages (one line)
Exception indicating invalid cast attempted.
std::string description() const
Return a simple one-line description of this object.
std::vector< Teuchos::RCP< const FactoryManagerBase > > FactManager_
vector of factory managers
std::string type_
smoother type
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
RCP< Matrix > A_
block operator
IndefBlockedDiagonalSmoother()
Constructor.
void Setup(Level ¤tLevel)
Setup routine.
RCP< Matrix > F_
fluid operator
Timer to be used in factories. Similar to Monitor but with additional timers.
Print additional debugging information.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
Namespace for MueLu classes and methods.
virtual const Teuchos::ParameterList & GetParameterList() const
void DeclareInput(Level ¤tLevel) const
Input.
RCP< const MapExtractorClass > rangeMapExtractor_
range map extractor (from A_ generated by AFact)
Class that holds all level-specific information.
RCP< const ParameterList > GetValidParameterList() const
Input.
bool IsSetup() const
Get the state of a smoother prototype.
RCP< Matrix > Z_
pressure stabilization term or null block
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Teuchos::RCP< SmootherBase > velPredictSmoo_
Block smoothers.
Teuchos::RCP< SmootherBase > schurCompSmoo_
smoother for SchurComplement equation
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos=0)
Add a factory manager for BraessSarazin internal SchurComplement handling.
Exception throws to report errors in the internal logical of the program.
An exception safe way to call the method 'Level::SetFactoryManager()'.
virtual ~IndefBlockedDiagonalSmoother()
Destructor.
RCP< SmootherPrototype > Copy() const
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
RCP< const MapExtractorClass > domainMapExtractor_
domain map extractor (from A_ generated by AFact)
virtual std::string description() const
Return a simple one-line description of this object.
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory()