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> 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_(
"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);
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;
173 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
178 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
183 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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;
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.");
224 Teuchos::RCP<CrsMatrix> A00 = bA->getMatrix(0, 0);
225 Teuchos::RCP<CrsMatrix> A01 = bA->getMatrix(0, 1);
226 Teuchos::RCP<CrsMatrix> A10 = bA->getMatrix(1, 0);
227 Teuchos::RCP<CrsMatrix> A11 = bA->getMatrix(1, 1);
229 Teuchos::RCP<CrsMatrixWrap> Op00 = Teuchos::rcp(
new CrsMatrixWrap(A00));
230 Teuchos::RCP<CrsMatrixWrap> Op01 = Teuchos::rcp(
new CrsMatrixWrap(A01));
231 Teuchos::RCP<CrsMatrixWrap> Op10 = Teuchos::rcp(
new CrsMatrixWrap(A10));
232 Teuchos::RCP<CrsMatrixWrap> Op11 = Teuchos::rcp(
new CrsMatrixWrap(A11));
234 F_ = Teuchos::rcp_dynamic_cast<Matrix>(Op00);
235 G_ = Teuchos::rcp_dynamic_cast<Matrix>(Op01);
236 D_ = Teuchos::rcp_dynamic_cast<Matrix>(Op10);
237 Z_ = Teuchos::rcp_dynamic_cast<Matrix>(Op11);
240 F_->CreateView(
"stridedMaps", bA->getRangeMap(0), bA->getDomainMap(0));
241 G_->CreateView(
"stridedMaps", bA->getRangeMap(0), bA->getDomainMap(1));
242 D_->CreateView(
"stridedMaps", bA->getRangeMap(1), bA->getDomainMap(0));
243 Z_->CreateView(
"stridedMaps", bA->getRangeMap(1), bA->getDomainMap(1));
246 bool bSIMPLEC = pL.get<
bool>(
"UseSIMPLEC");
249 RCP<Vector> diagFVector = VectorFactory::Build(
F_->getRowMap());
251 F_->getLocalDiagCopy(*diagFVector);
252 diagFVector->reciprocal(*diagFVector);
254 const RCP<const Map> rowmap =
F_->getRowMap();
255 size_t locSize = rowmap->getNodeNumElements();
256 Teuchos::ArrayRCP<SC> diag = diagFVector->getDataNonConst(0);
257 Teuchos::ArrayView<const LO> cols;
258 Teuchos::ArrayView<const SC> vals;
259 for (
size_t i=0; i<locSize; ++i) {
260 F_->getLocalRowView(i,cols,vals);
261 Scalar absRowSum = Teuchos::ScalarTraits<Scalar>::zero();
262 for (LO j=0; j<cols.size(); ++j) {
263 absRowSum += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
267 diagFVector->reciprocal(*diagFVector);
274 RCP<const FactoryManagerBase> velpredictFactManager =
FactManager_.at(0);
276 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
279 RCP<const FactoryManagerBase> schurFactManager =
FactManager_.at(1);
281 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
287 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
292 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
294 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
298 LocalOrdinal nSweeps = pL.get<LocalOrdinal>(
"Sweeps");
299 Scalar omega = pL.get<Scalar>(
"Damping factor");
302 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
306 RCP<MultiVector> residual = MultiVectorFactory::Build(B.getMap(), B.getNumVectors());
309 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
311 residual->update(one,B,zero);
312 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
320 RCP<MultiVector> xtilde1 = MultiVectorFactory::Build(
F_->getRowMap(),1);
321 xtilde1->putScalar(zero);
326 RCP<MultiVector> schurCompRHS = MultiVectorFactory::Build(
Z_->getRowMap(),1);
327 D_->apply(*xtilde1,*schurCompRHS);
328 schurCompRHS->update(one,*r2,-one);
332 RCP<MultiVector> xtilde2 = MultiVectorFactory::Build(
Z_->getRowMap(),1);
333 xtilde2->putScalar(zero);
338 RCP<MultiVector> xhat2 = MultiVectorFactory::Build(
Z_->getRowMap(),1);
339 xhat2->update(omega,*xtilde2,zero);
342 RCP<MultiVector> xhat1 = MultiVectorFactory::Build(
F_->getRowMap(),1);
343 RCP<MultiVector> xhat1_temp = MultiVectorFactory::Build(
F_->getRowMap(),1);
344 G_->apply(*xhat2,*xhat1_temp);
345 xhat1->elementWiseMultiply(one,*
diagFinv_,*xhat1_temp,zero);
346 xhat1->update(one,*xtilde1,-one);
354 x1->update(one,*xhat1,one);
355 x2->update( one,*xhat2,one);
364 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
365 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
370 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
372 std::ostringstream out;
374 out <<
"{type = " <<
type_ <<
"}";
378 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
383 out0 <<
"Prec. type: " <<
type_ << std::endl;
386 if (verbLevel &
Debug) {
Teuchos::RCP< SmootherBase > velPredictSmoo_
smoother for velocity prediction
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 type_
smoother type
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::RCP< SmootherBase > schurCompSmoo_
smoother for SchurComplement equation
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
std::vector< Teuchos::RCP< const FactoryManagerBase > > FactManager_
vector of factory managers
Timer to be used in factories. Similar to Monitor but with additional timers.
void Setup(Level ¤tLevel)
Setup routine.
Teuchos::RCP< Matrix > G_
pressure gradient operator
Print additional debugging information.
Namespace for MueLu classes and methods.
virtual const Teuchos::ParameterList & GetParameterList() const
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
RCP< Matrix > A_
block operator
Teuchos::RCP< Matrix > F_
fluid operator
RCP< const ParameterList > GetValidParameterList() const
Input.
std::string description() const
Return a simple one-line description of this object.
Class that holds all level-specific information.
Teuchos::RCP< Vector > diagFinv_
matrices
Teuchos::RCP< Matrix > Z_
pressure stabilization term or null block
Factory for building a thresholded operator.
Teuchos::RCP< Matrix > D_
divergence 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.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Factory for building the Schur Complement for a 2x2 block matrix.
RCP< const MapExtractorClass > domainMapExtractor_
domain map extractor (from A_ generated by AFact)
Exception throws to report errors in the internal logical of the program.
An exception safe way to call the method 'Level::SetFactoryManager()'.
RCP< const MapExtractorClass > rangeMapExtractor_
range map extractor (from A_ generated by AFact)
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.
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory()