uzawa_solver.hpp
Go to the documentation of this file.
1#ifndef UZAWA_SOLVER_HPP_
2#define UZAWA_SOLVER_HPP_
3//==============================================================================
13//==============================================================================
14
15#include <dune/istl/solvers.hh>
16
19
20namespace Opm {
21 namespace Elasticity {
22
26 template<class X, class Y>
27class UzawaSolver : public Dune::InverseOperator<X,Y>
28{
29 public:
30 typedef std::shared_ptr<Dune::InverseOperator<X,Y> > OperatorPtr;
35 UzawaSolver(OperatorPtr& innersolver_,
36 OperatorPtr& outersolver_,
37 const Matrix& B_) :
38 innersolver(innersolver_), outersolver(outersolver_), B(B_)
39 {
40 }
41
47 void apply(X& x, Y& b, double /* reduction */,
48 Dune::InverseOperatorResult& res) override
49 {
50 apply(x, b, res);
51 }
52
57 void apply(X& x, Y& b, Dune::InverseOperatorResult& res) override
58 {
59 Vector lambda, lambda2, u, u2;
61 Dune::InverseOperatorResult res2;
62 u2.resize(u.size());
63 u2 = 0;
64 innersolver->apply(u2, u, res2);
65 lambda.resize(B.M());
66 B.mtv(u2, lambda);
67 lambda2.resize(lambda.size());
68 lambda2 = 0;
69 outersolver->apply(lambda2, lambda, res2);
71 B.usmv(-1.0, lambda2, u);
72 u2 = 0;
73 innersolver->apply(u2, u, res);
74 MortarUtils::injectBlock(x, u2, B.N());
75 MortarUtils::injectBlock(x, lambda2, B.M(), B.N());
76 }
77
78 Dune::SolverCategory::Category category() const override
79 {
80 return Dune::SolverCategory::sequential;
81 }
82
83 protected:
86 const Matrix& B;
87};
88
89}
90}
91
92#endif
static void extractBlock(Vector &x, const Vector &y, int len, int start=0)
Extract a range of indices from a vector.
Definition: mortar_utils.hpp:27
static void injectBlock(Vector &x, const Vector &y, int len, int start=0)
Inject a range of indices into a vector.
Definition: mortar_utils.hpp:38
Definition: uzawa_solver.hpp:28
OperatorPtr outersolver
The outer solver.
Definition: uzawa_solver.hpp:85
std::shared_ptr< Dune::InverseOperator< X, Y > > OperatorPtr
Definition: uzawa_solver.hpp:30
Dune::SolverCategory::Category category() const override
Definition: uzawa_solver.hpp:78
const Matrix & B
The coupling matrix.
Definition: uzawa_solver.hpp:86
void apply(X &x, Y &b, double, Dune::InverseOperatorResult &res) override
Apply the scheme to a vector.
Definition: uzawa_solver.hpp:47
void apply(X &x, Y &b, Dune::InverseOperatorResult &res) override
Apply the scheme to a vector.
Definition: uzawa_solver.hpp:57
UzawaSolver(OperatorPtr &innersolver_, OperatorPtr &outersolver_, const Matrix &B_)
Default constructor.
Definition: uzawa_solver.hpp:35
OperatorPtr innersolver
The inner solver.
Definition: uzawa_solver.hpp:84
Helper class with some matrix operations.
Mortar helper class.
@ Y
Definition: mpc.hh:24
@ X
Definition: mpc.hh:24
Dune::BCRSMatrix< Dune::FieldMatrix< double, 1, 1 > > Matrix
A sparse matrix holding our operator.
Definition: matrixops.hpp:27
Dune::BlockVector< Dune::FieldVector< double, 1 > > Vector
A vector holding our RHS.
Definition: matrixops.hpp:33
Definition: ImplicitAssembly.hpp:43