16 #ifndef dealii_solver_minres_h 17 #define dealii_solver_minres_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/logstream.h> 23 #include <deal.II/base/signaling_nan.h> 24 #include <deal.II/base/subscriptor.h> 26 #include <deal.II/lac/solver.h> 27 #include <deal.II/lac/solver_control.h> 31 DEAL_II_NAMESPACE_OPEN
71 template <
class VectorType = Vector<
double>>
104 template <
typename MatrixType,
typename PreconditionerType>
106 solve(
const MatrixType & A,
108 const VectorType & b,
109 const PreconditionerType &preconditioner);
136 const VectorType & x,
137 const VectorType & r,
138 const VectorType & d)
const;
154 template <
class VectorType>
157 const AdditionalData &)
159 , res2(
numbers::signaling_nan<double>())
164 template <
class VectorType>
166 const AdditionalData &)
173 template <
class VectorType>
181 template <
class VectorType>
186 const VectorType &)
const 191 template <
class VectorType>
192 template <
typename MatrixType,
typename PreconditionerType>
196 const VectorType & b,
197 const PreconditionerType &preconditioner)
213 using vecptr = VectorType *;
214 vecptr u[3] = {Vu0.get(), Vu1.get(), Vu2.get()};
215 vecptr m[3] = {Vm0.get(), Vm1.get(), Vm2.get()};
220 u[0]->reinit(b,
true);
221 u[1]->reinit(b,
true);
222 u[2]->reinit(b,
true);
223 m[0]->reinit(b,
true);
224 m[1]->reinit(b,
true);
225 m[2]->reinit(b,
true);
229 double delta[3] = {0, 0, 0};
230 double f[2] = {0, 0};
231 double e[2] = {0, 0};
253 preconditioner.vmult(v, *u[1]);
255 delta[1] = v * (*u[1]);
257 Assert(delta[1] >= 0, ExcPreconditionerNotDefinite());
259 r0 = std::sqrt(delta[1]);
273 v *= 1. / std::sqrt(delta[1]);
278 u[2]->add(-std::sqrt(delta[1] / delta[0]), *u[0]);
280 const double gamma = *u[2] * v;
281 u[2]->add(-gamma / std::sqrt(delta[1]), *u[1]);
287 preconditioner.vmult(v, *u[2]);
289 delta[2] = v * (*u[2]);
291 Assert(delta[2] >= 0, ExcPreconditionerNotDefinite());
296 e[1] = std::sqrt(delta[2]);
300 d_ = s *
e[0] - c * gamma;
301 e[0] = c *
e[0] + s * gamma;
302 f[1] = s * std::sqrt(delta[2]);
303 e[1] = -c * std::sqrt(delta[2]);
306 const double d = std::sqrt(d_ * d_ + delta[2]);
313 s = std::sqrt(delta[2]) /
d;
318 m[0]->add(-e[0], *m[1]);
320 m[0]->add(-f[0], *m[2]);
323 r_l2 *= std::fabs(s);
325 conv = this->iteration_status(j, r_l2, x);
358 DEAL_II_NAMESPACE_CLOSE
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverMinRes(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
#define AssertThrow(cond, exc)
void swap(BlockIndices &u, BlockIndices &v)
Stop iteration, goal reached.
#define Assert(cond, exc)
#define DeclException0(Exception0)
virtual ~SolverMinRes() override=default
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
virtual double criterion()
static ::ExceptionBase & ExcPreconditionerNotDefinite()
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const