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_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_(
"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);
111 if (myPos < FactManager_.size()) {
113 FactManager_.at(myPos) = FactManager;
114 }
else if( myPos == FactManager_.size()) {
116 FactManager_.push_back(FactManager);
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;
122 FactManager_.push_back(FactManager);
126 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 currentLevel.
DeclareInput(
"A",this->GetFactory(
"A").
get());
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;
134 for(it = FactManager_.begin(); it!=FactManager_.end(); ++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.");
157 rangeMapExtractor_ = bA->getRangeMapExtractor();
158 domainMapExtractor_ = bA->getDomainMapExtractor();
161 Teuchos::RCP<Matrix> A00 = bA->getMatrix(0, 0);
162 Teuchos::RCP<Matrix> A01 = bA->getMatrix(0, 1);
163 Teuchos::RCP<Matrix> A10 = bA->getMatrix(1, 0);
164 Teuchos::RCP<Matrix> A11 = bA->getMatrix(1, 1);
198 RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
200 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
203 RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
205 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
210 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
215 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
217 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
221 RCP<BlockedCrsMatrix> bA00 = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_);
222 RCP<BlockedCrsMatrix> bA11 = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Z_);
223 bool bA00ThyraSpecialTreatment =
false;
224 bool bA11ThyraSpecialTreatment =
false;
225 if (bA00 != Teuchos::null) {
226 if(bA00->Rows() == 1 && bA00->Cols() == 1 && rangeMapExtractor_->getThyraMode() ==
true) bA00ThyraSpecialTreatment =
true;
228 if (bA11 != Teuchos::null) {
229 if(bA11->Rows() == 1 && bA11->Cols() == 1 && rangeMapExtractor_->getThyraMode() ==
true) bA11ThyraSpecialTreatment =
true;
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);
251 Teuchos::RCP<MultiVector> r1 = rangeMapExtractor_->ExtractVector(residual, 0);
252 Teuchos::RCP<MultiVector> r2 = rangeMapExtractor_->ExtractVector(residual, 1);
256 RCP<MultiVector> xtilde1 = MultiVectorFactory::Build(F_->getRowMap(),X.getNumVectors(),
true);
260 if(bA00ThyraSpecialTreatment ==
true) {
261 RCP<MultiVector> xtilde1_thyra = domainMapExtractor_->getVector(0, X.getNumVectors(),
true);
262 RCP<MultiVector> r1_thyra = rangeMapExtractor_->getVector(0, B.getNumVectors(),
true);
264 for(
size_t k=0; k < r1->getNumVectors(); k++) {
265 Teuchos::ArrayRCP<const Scalar> xpetraVecData = r1->getData(k);
266 Teuchos::ArrayRCP<Scalar> thyraVecData = r1_thyra->getDataNonConst(k);
267 for(
size_t i=0; i < r1->getLocalLength(); i++) {
268 thyraVecData[i] = xpetraVecData[i];
272 velPredictSmoo_->Apply(*xtilde1_thyra,*r1_thyra);
274 for(
size_t k=0; k < xtilde1_thyra->getNumVectors(); k++) {
275 Teuchos::ArrayRCP<Scalar> xpetraVecData = xtilde1->getDataNonConst(k);
276 Teuchos::ArrayRCP<const Scalar> thyraVecData = xtilde1_thyra->getData(k);
277 for(
size_t i=0; i < xtilde1_thyra->getLocalLength(); i++) {
278 xpetraVecData[i] = thyraVecData[i];
282 velPredictSmoo_->Apply(*xtilde1,*r1);
287 RCP<MultiVector> xtilde2 = MultiVectorFactory::Build(Z_->getRowMap(),X.getNumVectors(),
true);
291 if(bA11ThyraSpecialTreatment ==
true) {
292 RCP<MultiVector> xtilde2_thyra = domainMapExtractor_->getVector(1, X.getNumVectors(),
true);
293 RCP<MultiVector> r2_thyra = rangeMapExtractor_->getVector(1, B.getNumVectors(),
true);
295 for(
size_t k=0; k < r2->getNumVectors(); k++) {
296 Teuchos::ArrayRCP<const Scalar> xpetraVecData = r2->getData(k);
297 Teuchos::ArrayRCP<Scalar> thyraVecData = r2_thyra->getDataNonConst(k);
298 for(
size_t i=0; i < r2->getLocalLength(); i++) {
299 thyraVecData[i] = xpetraVecData[i];
303 schurCompSmoo_->Apply(*xtilde2_thyra,*r2_thyra);
305 for(
size_t k=0; k < xtilde2_thyra->getNumVectors(); k++) {
306 Teuchos::ArrayRCP<Scalar> xpetraVecData = xtilde2->getDataNonConst(k);
307 Teuchos::ArrayRCP<const Scalar> thyraVecData = xtilde2_thyra->getData(k);
308 for(
size_t i=0; i < xtilde2_thyra->getLocalLength(); i++) {
309 xpetraVecData[i] = thyraVecData[i];
313 schurCompSmoo_->Apply(*xtilde2,*r2);
317 Teuchos::RCP<MultiVector> x1 = domainMapExtractor_->ExtractVector(rcpX, 0);
318 Teuchos::RCP<MultiVector> x2 = domainMapExtractor_->ExtractVector(rcpX, 1);
322 x1->update(omega,*xtilde1,one);
323 x2->update(omega,*xtilde2,one);
326 domainMapExtractor_->InsertVector(x1, 0, rcpX);
327 domainMapExtractor_->InsertVector(x2, 1, rcpX);
331 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
332 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
337 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
339 std::ostringstream out;
341 out <<
"{type = " << type_ <<
"}";
345 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
350 out0 <<
"Prec. type: " << type_ << std::endl;
353 if (verbLevel &
Debug) {
Important warning messages (one line)
Exception indicating invalid cast attempted.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
IndefBlockedDiagonalSmoother()
Constructor.
bool IsSetup() const
Get the state of a smoother prototype.
void Setup(Level ¤tLevel)
Setup routine.
Timer to be used in factories. Similar to Monitor but with additional timers.
RCP< const ParameterList > GetValidParameterList() const
Input.
Print additional debugging information.
Namespace for MueLu classes and methods.
void DeclareInput(Level ¤tLevel) const
Input.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
RCP< SmootherPrototype > Copy() const
Class that holds all level-specific information.
Cheap Blocked diagonal smoother for indefinite 2x2 block matrices.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
virtual const Teuchos::ParameterList & GetParameterList() const
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=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.
std::string description() const
Return a simple one-line description of this object.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual std::string description() const
Return a simple one-line description of this object.