53 #ifndef MUELU_SIMPLESMOOTHER_DEF_HPP_ 54 #define MUELU_SIMPLESMOOTHER_DEF_HPP_ 56 #include "Teuchos_ArrayViewDecl.hpp" 57 #include "Teuchos_ScalarTraits.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_(
"SIMPLE"), A_(
Teuchos::null)
96 SchurFact->SetParameter(
"omega",Teuchos::ParameterEntry(omega));
97 SchurFact->SetParameter(
"lumping",Teuchos::ParameterEntry(SIMPLEC));
98 SchurFact->SetFactory(
"A", this->
GetFactory(
"A"));
101 Teuchos::ParameterList SCparams;
103 RCP<SmootherPrototype> smoProtoSC = rcp(
new DirectSolver(SCtype,SCparams) );
104 smoProtoSC->SetFactory(
"A", SchurFact);
106 RCP<SmootherFactory> SmooSCFact = rcp(
new SmootherFactory(smoProtoSC) );
109 schurFactManager->SetFactory(
"A", SchurFact);
110 schurFactManager->SetFactory(
"Smoother", SmooSCFact);
111 schurFactManager->SetIgnoreUserData(
true);
115 A00Fact->SetFactory(
"A",this->
GetFactory(
"A"));
116 A00Fact->SetParameter(
"block row",ParameterEntry(0));
117 A00Fact->SetParameter(
"block col",ParameterEntry(0));
118 Teuchos::ParameterList velpredictParams;
119 std::string velpredictType;
120 RCP<SmootherPrototype> smoProtoPredict = rcp(
new DirectSolver(velpredictType,velpredictParams) );
121 smoProtoPredict->SetFactory(
"A", A00Fact);
122 RCP<SmootherFactory> SmooPredictFact = rcp(
new SmootherFactory(smoProtoPredict) );
124 RCP<FactoryManager> velpredictFactManager = rcp(
new FactoryManager());
125 velpredictFactManager->SetFactory(
"A", A00Fact);
126 velpredictFactManager->SetFactory(
"Smoother", SmooPredictFact);
127 velpredictFactManager->SetIgnoreUserData(
true);
134 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
137 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
139 RCP<ParameterList> validParamList = rcp(
new ParameterList());
141 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
142 validParamList->set< Scalar > (
"Damping factor", 1.0,
"Damping/Scaling factor in SIMPLE");
143 validParamList->set< LocalOrdinal > (
"Sweeps", 1,
"Number of SIMPLE sweeps (default = 1)");
144 validParamList->set<
bool > (
"UseSIMPLEC",
false,
"Use SIMPLEC instead of SIMPLE (default = false)");
146 return validParamList;
150 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
152 TEUCHOS_TEST_FOR_EXCEPTION(pos < 0,
Exceptions::RuntimeError,
"MueLu::SimpleSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
154 size_t myPos = Teuchos::as<size_t>(pos);
156 if (myPos < FactManager_.size()) {
158 FactManager_.at(myPos) = FactManager;
159 }
else if( myPos == FactManager_.size()) {
161 FactManager_.push_back(FactManager);
163 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
164 *out <<
"Warning: cannot add new FactoryManager at proper position " << pos <<
". The FactoryManager is just appended to the end. Check this!" << std::endl;
167 FactManager_.push_back(FactManager);
173 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
175 AddFactoryManager(FactManager, 0);
178 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
180 AddFactoryManager(FactManager, 1);
183 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
185 currentLevel.
DeclareInput(
"A",this->GetFactory(
"A").
get());
187 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!");
190 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
191 for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
195 currentLevel.
DeclareInput(
"PreSmoother",(*it)->GetFactory(
"Smoother").get());
199 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
208 FactoryMonitor m(*
this,
"Setup blocked SIMPLE Smoother", currentLevel);
211 this->GetOStream(
Warnings0) <<
"MueLu::SimpleSmoother::Setup(): Setup() has already been called";
214 A_ = Factory::Get<RCP<Matrix> > (currentLevel,
"A");
216 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
217 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null,
Exceptions::BadCast,
"MueLu::SimpleSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
220 rangeMapExtractor_ = bA->getRangeMapExtractor();
221 domainMapExtractor_ = bA->getDomainMapExtractor();
224 F_ = bA->getMatrix(0, 0);
225 G_ = bA->getMatrix(0, 1);
226 D_ = bA->getMatrix(1, 0);
227 Z_ = bA->getMatrix(1, 1);
230 bool bSIMPLEC = pL.get<
bool>(
"UseSIMPLEC");
234 RCP<Vector> diagFVector = VectorFactory::Build(F_->getRowMap());
236 F_->getLocalDiagCopy(*diagFVector);
258 RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
260 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
263 RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
265 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
271 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
275 #ifdef HAVE_MUELU_DEBUG 276 TEUCHOS_TEST_FOR_EXCEPTION(A_->getRangeMap()->isSameAs(*(B.getMap())) ==
false,
Exceptions::RuntimeError,
"MueLu::SimpleSmoother::Apply(): The map of RHS vector B is not the same as range map of the blocked operator A. Please check the map of B and A.");
277 TEUCHOS_TEST_FOR_EXCEPTION(A_->getDomainMap()->isSameAs(*(X.getMap())) ==
false,
Exceptions::RuntimeError,
"MueLu::SimpleSmoother::Apply(): The map of the solution vector X is not the same as domain map of the blocked operator A. Please check the map of X and A.");
280 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
282 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
286 LocalOrdinal nSweeps = pL.get<LocalOrdinal>(
"Sweeps");
287 Scalar omega = pL.get<Scalar>(
"Damping factor");
296 bool bRangeThyraModePredict = rangeMapExtractor_->getThyraMode() && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_) == Teuchos::null);
297 bool bDomainThyraModePredict = domainMapExtractor_->getThyraMode() && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_) == Teuchos::null);
298 bool bRangeThyraModeSchur = rangeMapExtractor_->getThyraMode() && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Z_) == Teuchos::null);
299 bool bDomainThyraModeSchur = domainMapExtractor_->getThyraMode() && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Z_) == Teuchos::null);
303 RCP<BlockedCrsMatrix> bF = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_);
304 RCP<BlockedCrsMatrix> bZ = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Z_);
305 bool bFThyraSpecialTreatment =
false;
306 bool bZThyraSpecialTreatment =
false;
307 if (bF != Teuchos::null) {
308 if(bF->Rows() == 1 && bF->Cols() == 1 && rangeMapExtractor_->getThyraMode() ==
true) bFThyraSpecialTreatment =
true;
310 if (bZ != Teuchos::null) {
311 if(bZ->Rows() == 1 && bZ->Cols() == 1 && rangeMapExtractor_->getThyraMode() ==
true) bZThyraSpecialTreatment =
true;
314 #if 1// new implementation 317 RCP<MultiVector> res = MultiVectorFactory::Build(B.getMap(), B.getNumVectors(),
true);
318 RCP<BlockedMultiVector> residual = Teuchos::rcp(
new BlockedMultiVector(rangeMapExtractor_,res));
321 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
322 RCP<BlockedMultiVector> bX = Teuchos::rcp(
new BlockedMultiVector(domainMapExtractor_,rcpX));
325 RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
326 RCP<const BlockedMultiVector> bB = Teuchos::rcp(
new const BlockedMultiVector(rangeMapExtractor_,rcpB));
330 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
332 residual->update(one,*bB,zero);
333 A_->apply(*bX, *residual, Teuchos::NO_TRANS, -one, one);
336 Teuchos::RCP<MultiVector> r1 = rangeMapExtractor_->ExtractVector(residual, 0, bRangeThyraModePredict);
337 Teuchos::RCP<MultiVector> r2 = rangeMapExtractor_->ExtractVector(residual, 1, bRangeThyraModeSchur);
341 RCP<MultiVector> xtilde1 = domainMapExtractor_->getVector(0, X.getNumVectors(), bDomainThyraModePredict,
true);
344 if(bFThyraSpecialTreatment ==
true) {
345 xtilde1->replaceMap(domainMapExtractor_->getMap(0,
true));
346 r1->replaceMap(rangeMapExtractor_->getMap(0,
true));
347 velPredictSmoo_->Apply(*xtilde1,*r1);
348 xtilde1->replaceMap(domainMapExtractor_->getMap(0,
false));
350 velPredictSmoo_->Apply(*xtilde1,*r1);
355 RCP<MultiVector> schurCompRHS = rangeMapExtractor_->getVector(1, B.getNumVectors(), bRangeThyraModeSchur,
false);
356 if(D_.is_null() ==
false)
357 D_->apply(*xtilde1,*schurCompRHS);
359 schurCompRHS->putScalar(zero);
360 schurCompRHS->update(one,*r2,-one);
364 RCP<MultiVector> xtilde2 = domainMapExtractor_->getVector(1, X.getNumVectors(), bDomainThyraModeSchur,
true);
369 if(bZThyraSpecialTreatment ==
true) {
370 xtilde2->replaceMap(domainMapExtractor_->getMap(1,
true));
371 schurCompRHS->replaceMap(rangeMapExtractor_->getMap(1,
true));
372 schurCompSmoo_->Apply(*xtilde2,*schurCompRHS);
373 xtilde2->replaceMap(domainMapExtractor_->getMap(1,
false));
375 schurCompSmoo_->Apply(*xtilde2,*schurCompRHS);
380 RCP<MultiVector> xhat2 = domainMapExtractor_->getVector(1, X.getNumVectors(), bDomainThyraModeSchur,
false);
381 xhat2->update(omega,*xtilde2,zero);
384 RCP<MultiVector> xhat1 = domainMapExtractor_->getVector(0, X.getNumVectors(), bDomainThyraModePredict,
false);
385 RCP<MultiVector> xhat1_temp = domainMapExtractor_->getVector(0, X.getNumVectors(), bDomainThyraModePredict,
false);
386 if(G_.is_null() ==
false)
387 G_->apply(*xhat2,*xhat1_temp);
389 xhat1_temp->putScalar(zero);
390 xhat1->elementWiseMultiply(one,*diagFinv_,*xhat1_temp,zero);
391 xhat1->update(one,*xtilde1,-one);
394 Teuchos::RCP<MultiVector> x1 = domainMapExtractor_->ExtractVector(bX, 0, bDomainThyraModePredict);
395 Teuchos::RCP<MultiVector> x2 = domainMapExtractor_->ExtractVector(bX, 1, bDomainThyraModeSchur);
399 x1->update(one,*xhat1,one);
400 x2->update( one,*xhat2,one);
402 domainMapExtractor_->InsertVector(x1, 0, bX, bDomainThyraModePredict);
403 domainMapExtractor_->InsertVector(x2, 1, bX, bDomainThyraModeSchur);
407 domainMapExtractor_->InsertVector(bX->getMultiVector(0,bDomainThyraModePredict), 0, rcpX, bDomainThyraModePredict);
408 domainMapExtractor_->InsertVector(bX->getMultiVector(1,bDomainThyraModeSchur), 1, rcpX, bDomainThyraModeSchur);
412 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
416 RCP<MultiVector> residual = MultiVectorFactory::Build(B.getMap(), B.getNumVectors());
419 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
421 residual->update(one,B,zero);
422 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
424 Teuchos::RCP<MultiVector> r1 = rangeMapExtractor_->ExtractVector(residual, 0, bRangeThyraModePredict);
425 Teuchos::RCP<MultiVector> r2 = rangeMapExtractor_->ExtractVector(residual, 1, bRangeThyraModeSchur);
429 RCP<MultiVector> xtilde1 = domainMapExtractor_->getVector(0, X.getNumVectors(), bDomainThyraModePredict);
430 xtilde1->putScalar(zero);
434 if(bFThyraSpecialTreatment ==
true) {
436 RCP<MultiVector> xtilde1_thyra = domainMapExtractor_->getVector(0, X.getNumVectors(),
true);
438 Teuchos::RCP<MultiVector> r1_thyra = rangeMapExtractor_->ExtractVector(residual, 0,
true);
439 velPredictSmoo_->Apply(*xtilde1_thyra,*r1_thyra);
440 for(
size_t k=0; k < xtilde1_thyra->getNumVectors(); k++) {
441 Teuchos::ArrayRCP<Scalar> xpetraVecData = xtilde1->getDataNonConst(k);
442 Teuchos::ArrayRCP<const Scalar> thyraVecData = xtilde1_thyra->getData(k);
443 for(
size_t i=0; i < xtilde1_thyra->getLocalLength(); i++) {
444 xpetraVecData[i] = thyraVecData[i];
448 velPredictSmoo_->Apply(*xtilde1,*r1);
453 RCP<MultiVector> schurCompRHS = rangeMapExtractor_->getVector(1, B.getNumVectors(), bRangeThyraModeSchur);
454 D_->apply(*xtilde1,*schurCompRHS);
455 schurCompRHS->update(one,*r2,-one);
459 RCP<MultiVector> xtilde2 = domainMapExtractor_->getVector(1, X.getNumVectors(), bDomainThyraModeSchur);
460 xtilde2->putScalar(zero);
464 if(bZThyraSpecialTreatment ==
true) {
466 RCP<MultiVector> xtilde2_thyra = domainMapExtractor_->getVector(1, X.getNumVectors(),
true);
468 RCP<MultiVector> schurCompRHS_thyra = rangeMapExtractor_->getVector(1, B.getNumVectors(),
true);
470 for(
size_t k=0; k < schurCompRHS->getNumVectors(); k++) {
471 Teuchos::ArrayRCP<const Scalar> xpetraVecData = schurCompRHS->getData(k);
472 Teuchos::ArrayRCP<Scalar> thyraVecData = schurCompRHS_thyra->getDataNonConst(k);
473 for(
size_t i=0; i < schurCompRHS->getLocalLength(); i++) {
474 thyraVecData[i] = xpetraVecData[i];
478 schurCompSmoo_->Apply(*xtilde2_thyra,*schurCompRHS_thyra);
480 for(
size_t k=0; k < xtilde2_thyra->getNumVectors(); k++) {
481 Teuchos::ArrayRCP<Scalar> xpetraVecData = xtilde2->getDataNonConst(k);
482 Teuchos::ArrayRCP<const Scalar> thyraVecData = xtilde2_thyra->getData(k);
483 for(
size_t i=0; i < xtilde2_thyra->getLocalLength(); i++) {
484 xpetraVecData[i] = thyraVecData[i];
488 schurCompSmoo_->Apply(*xtilde2,*schurCompRHS);
493 RCP<MultiVector> xhat2 = domainMapExtractor_->getVector(1, X.getNumVectors(), bDomainThyraModeSchur);
494 xhat2->update(omega,*xtilde2,zero);
497 RCP<MultiVector> xhat1 = domainMapExtractor_->getVector(0, X.getNumVectors(), bDomainThyraModePredict);
498 RCP<MultiVector> xhat1_temp = domainMapExtractor_->getVector(0, X.getNumVectors(), bDomainThyraModePredict);
499 G_->apply(*xhat2,*xhat1_temp);
500 xhat1->elementWiseMultiply(one,*diagFinv_,*xhat1_temp,zero);
501 xhat1->update(one,*xtilde1,-one);
504 Teuchos::RCP<MultiVector> x1 = domainMapExtractor_->ExtractVector(rcpX, 0, bDomainThyraModePredict);
505 Teuchos::RCP<MultiVector> x2 = domainMapExtractor_->ExtractVector(rcpX, 1, bDomainThyraModeSchur);
509 x1->update(one,*xhat1,one);
510 x2->update( one,*xhat2,one);
512 domainMapExtractor_->InsertVector(x1, 0, rcpX, bDomainThyraModePredict);
513 domainMapExtractor_->InsertVector(x2, 1, rcpX, bDomainThyraModeSchur);
518 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
519 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
524 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
526 std::ostringstream out;
528 out <<
"{type = " << type_ <<
"}";
532 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
537 out0 <<
"Prec. type: " << type_ << std::endl;
540 if (verbLevel &
Debug) {
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...
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 ...
void DeclareInput(Level ¤tLevel) const
Input.
Class that holds all level-specific information.
Factory for building a thresholded operator.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
virtual const Teuchos::ParameterList & GetParameterList() const
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.
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.