53 #ifndef MUELU_UZAWASMOOTHER_DEF_HPP_ 54 #define MUELU_UZAWASMOOTHER_DEF_HPP_ 56 #include "Teuchos_ArrayViewDecl.hpp" 57 #include "Teuchos_ScalarTraits.hpp" 61 #include <Xpetra_Matrix.hpp> 62 #include <Xpetra_BlockedCrsMatrix.hpp> 63 #include <Xpetra_MultiVectorFactory.hpp> 64 #include <Xpetra_VectorFactory.hpp> 68 #include "MueLu_Utilities.hpp" 70 #include "MueLu_HierarchyUtils.hpp" 72 #include "MueLu_SubBlockAFactory.hpp" 75 #include "MueLu_SchurComplementFactory.hpp" 76 #include "MueLu_DirectSolver.hpp" 77 #include "MueLu_SmootherFactory.hpp" 78 #include "MueLu_FactoryManager.hpp" 82 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
84 : type_(
"Uzawa"), A_(
Teuchos::null)
88 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
91 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
93 RCP<ParameterList> validParamList = rcp(
new ParameterList());
95 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
96 validParamList->set< Scalar > (
"Damping factor", 1.0,
"Damping/Scaling factor in SIMPLE");
97 validParamList->set< LocalOrdinal > (
"Sweeps", 1,
"Number of SIMPLE sweeps (default = 1)");
99 return validParamList;
103 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
105 TEUCHOS_TEST_FOR_EXCEPTION(pos < 0,
Exceptions::RuntimeError,
"MueLu::UzawaSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
107 size_t myPos = Teuchos::as<size_t>(pos);
116 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
117 *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>
131 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
136 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 TEUCHOS_TEST_FOR_EXCEPTION(
FactManager_.size() != 2,
Exceptions::RuntimeError,
"MueLu::UzawaSmoother::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!");
143 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
148 currentLevel.
DeclareInput(
"PreSmoother",(*it)->GetFactory(
"Smoother").get());
152 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
155 FactoryMonitor m(*
this,
"Setup blocked Uzawa Smoother", currentLevel);
158 this->
GetOStream(
Warnings0) <<
"MueLu::UzawaSmoother::Setup(): Setup() has already been called";
161 A_ = Factory::Get<RCP<Matrix> > (currentLevel,
"A");
163 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
164 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null,
Exceptions::BadCast,
"MueLu::UzawaSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
171 F_ = bA->getMatrix(0, 0);
172 G_ = bA->getMatrix(0, 1);
173 D_ = bA->getMatrix(1, 0);
174 Z_ = bA->getMatrix(1, 1);
179 RCP<const FactoryManagerBase> velpredictFactManager =
FactManager_.at(0);
181 velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother",velpredictFactManager->GetFactory(
"Smoother").get());
184 RCP<const FactoryManagerBase> schurFactManager =
FactManager_.at(1);
186 schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >(
"PreSmoother", schurFactManager->GetFactory(
"Smoother").get());
192 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
197 Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
199 SC zero = Teuchos::ScalarTraits<SC>::zero(), one = Teuchos::ScalarTraits<SC>::one();
203 RCP<BlockedCrsMatrix> bA00 = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(
F_);
204 RCP<BlockedCrsMatrix> bA11 = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(
Z_);
205 bool bFThyraSpecialTreatment =
false;
206 bool bZThyraSpecialTreatment =
false;
207 if (bA00 != Teuchos::null) {
208 if(bA00->Rows() == 1 && bA00->Cols() == 1 &&
rangeMapExtractor_->getThyraMode() ==
true) bFThyraSpecialTreatment =
true;
210 if (bA11 != Teuchos::null) {
211 if(bA11->Rows() == 1 && bA11->Cols() == 1 &&
rangeMapExtractor_->getThyraMode() ==
true) bZThyraSpecialTreatment =
true;
216 LocalOrdinal nSweeps = pL.get<LocalOrdinal>(
"Sweeps");
217 Scalar omega = pL.get<Scalar>(
"Damping factor");
220 RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
224 RCP<MultiVector> residual = MultiVectorFactory::Build(B.getMap(), B.getNumVectors());
227 for (LocalOrdinal run = 0; run < nSweeps; ++run) {
229 residual->update(one,B,zero);
230 A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
238 RCP<MultiVector> xtilde1 = MultiVectorFactory::Build(
F_->getRowMap(),X.getNumVectors(),
true);
242 if(bFThyraSpecialTreatment ==
true) {
243 RCP<MultiVector> xtilde1_thyra =
domainMapExtractor_->getVector(0, X.getNumVectors(),
true);
244 RCP<MultiVector> r1_thyra =
rangeMapExtractor_->getVector(0, B.getNumVectors(),
true);
246 for(
size_t k=0; k < r1->getNumVectors(); k++) {
247 Teuchos::ArrayRCP<const Scalar> xpetraVecData = r1->getData(k);
248 Teuchos::ArrayRCP<Scalar> thyraVecData = r1_thyra->getDataNonConst(k);
249 for(
size_t i=0; i < r1->getLocalLength(); i++) {
250 thyraVecData[i] = xpetraVecData[i];
256 for(
size_t k=0; k < xtilde1_thyra->getNumVectors(); k++) {
257 Teuchos::ArrayRCP<Scalar> xpetraVecData = xtilde1->getDataNonConst(k);
258 Teuchos::ArrayRCP<const Scalar> thyraVecData = xtilde1_thyra->getData(k);
259 for(
size_t i=0; i < xtilde1_thyra->getLocalLength(); i++) {
260 xpetraVecData[i] = thyraVecData[i];
269 RCP<MultiVector> schurCompRHS = MultiVectorFactory::Build(
Z_->getRowMap(),B.getNumVectors());
270 D_->apply(*xtilde1,*schurCompRHS);
271 schurCompRHS->update(one,*r2,-one);
275 RCP<MultiVector> xtilde2 = MultiVectorFactory::Build(
Z_->getRowMap(),X.getNumVectors(),
true);
279 if(bZThyraSpecialTreatment ==
true) {
280 RCP<MultiVector> xtilde2_thyra =
domainMapExtractor_->getVector(1, X.getNumVectors(),
true);
281 RCP<MultiVector> schurCompRHS_thyra =
rangeMapExtractor_->getVector(1, B.getNumVectors(),
true);
283 for(
size_t k=0; k < schurCompRHS->getNumVectors(); k++) {
284 Teuchos::ArrayRCP<const Scalar> xpetraVecData = schurCompRHS->getData(k);
285 Teuchos::ArrayRCP<Scalar> thyraVecData = schurCompRHS_thyra->getDataNonConst(k);
286 for(
size_t i=0; i < schurCompRHS->getLocalLength(); i++) {
287 thyraVecData[i] = xpetraVecData[i];
293 for(
size_t k=0; k < xtilde2_thyra->getNumVectors(); k++) {
294 Teuchos::ArrayRCP<Scalar> xpetraVecData = xtilde2->getDataNonConst(k);
295 Teuchos::ArrayRCP<const Scalar> thyraVecData = xtilde2_thyra->getData(k);
296 for(
size_t i=0; i < xtilde2_thyra->getLocalLength(); i++) {
297 xpetraVecData[i] = thyraVecData[i];
310 x1->update(omega,*xtilde1,one);
311 x2->update(omega,*xtilde2,one);
319 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
320 RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
325 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
327 std::ostringstream out;
329 out <<
"{type = " <<
type_ <<
"}";
333 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
338 out0 <<
"Prec. type: " <<
type_ << std::endl;
341 if (verbLevel &
Debug) {
Important warning messages (one line)
std::string description() const
Return a simple one-line description of this object.
Exception indicating invalid cast attempted.
void Setup(Level ¤tLevel)
Setup routine.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsSetup() const
Get the state of a smoother prototype.
Timer to be used in factories. Similar to Monitor but with additional timers.
Print additional debugging information.
Teuchos::RCP< Matrix > D_
divergence operator
void DeclareInput(Level ¤tLevel) const
Input.
RCP< const MapExtractorClass > domainMapExtractor_
domain map extractor (from A_ generated by AFact)
Namespace for MueLu classes and methods.
RCP< const MapExtractorClass > rangeMapExtractor_
range map extractor (from A_ generated by AFact)
UzawaSmoother()
Constructor.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
RCP< const ParameterList > GetValidParameterList() const
Input.
Teuchos::RCP< Matrix > F_
fluid operator
Class that holds all level-specific information.
RCP< Matrix > A_
block operator
Teuchos::RCP< Matrix > Z_
pressure stabilization term or null block
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
std::vector< Teuchos::RCP< const FactoryManagerBase > > FactManager_
vector of factory managers
virtual const Teuchos::ParameterList & GetParameterList() const
Teuchos::RCP< SmootherBase > schurCompSmoo_
smoother for SchurComplement equation
RCP< SmootherPrototype > Copy() const
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
Teuchos::RCP< Matrix > G_
pressure gradient operator
Teuchos::RCP< SmootherBase > velPredictSmoo_
smoother for velocity prediction
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
std::string type_
smoother type
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory()
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos)
Add a factory manager at a specific position.
Exception throws to report errors in the internal logical of the program.
An exception safe way to call the method 'Level::SetFactoryManager()'.
void Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual ~UzawaSmoother()
Destructor.
virtual std::string description() const
Return a simple one-line description of this object.