1#ifndef UZAWA_SOLVER_HPP_
2#define UZAWA_SOLVER_HPP_
15#include <dune/istl/solvers.hh>
21 namespace Elasticity {
26 template<
class X,
class Y>
30 typedef std::shared_ptr<Dune::InverseOperator<X,Y> > OperatorPtr;
36 OperatorPtr& outersolver_,
48 Dune::InverseOperatorResult& res)
override
57 void apply(X& x, Y& b, Dune::InverseOperatorResult& res)
override
59 Vector lambda, lambda2, u, u2;
61 Dune::InverseOperatorResult res2;
67 lambda2.resize(lambda.size());
71 B.usmv(-1.0, lambda2, u);
78 Dune::SolverCategory::Category category()
const override
80 return Dune::SolverCategory::sequential;
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
OperatorPtr outersolver
The outer solver.
Definition uzawa_solver.hpp:85
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.
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
Inverting small matrices.
Definition ImplicitAssembly.hpp:43