16#ifndef dealii_solver_qmrs_h
17#define dealii_solver_qmrs_h
91template <
typename VectorType = Vector<
double>>
175 template <
typename MatrixType,
typename PreconditionerType>
179 const VectorType &
b,
180 const PreconditionerType &preconditioner);
189 const VectorType & x,
190 const VectorType & r,
191 const VectorType &
d)
const;
217 template <
typename MatrixType,
typename PreconditionerType>
221 const VectorType &
b,
222 const PreconditionerType &preconditioner,
241template <
class VectorType>
244 const double last_residual)
246 , last_residual(last_residual)
251template <
class VectorType>
254 const AdditionalData & data)
256 , additional_data(data)
260template <
class VectorType>
262 const AdditionalData &data)
264 , additional_data(data)
268template <
class VectorType>
273 const VectorType &)
const
276template <
class VectorType>
277template <
typename MatrixType,
typename PreconditionerType>
281 const VectorType &
b,
282 const PreconditionerType &preconditioner)
312 deallog <<
"Restart step " << step << std::endl;
313 state = iterate(
A, x,
b, preconditioner, *Vr, *Vu, *Vq, *Vt, *Vd);
324template <
class VectorType>
325template <
typename MatrixType,
typename PreconditionerType>
329 const VectorType &
b,
330 const PreconditionerType &preconditioner,
341 double tau, rho, theta = 0;
349 if (additional_data.left_preconditioning)
352 preconditioner.vmult(t, r);
359 preconditioner.vmult(q, t);
378 const double sigma = q * t;
381 if (additional_data.breakdown_testing ==
true &&
382 std::fabs(sigma) < additional_data.breakdown_threshold)
385 const double alpha = rho / sigma;
391 const double theta_old = theta;
394 if (additional_data.left_preconditioning)
397 preconditioner.vmult(t, r);
407 const double psi = 1. / (1. + theta);
411 d.sadd(psi * theta_old, psi * alpha, q);
414 print_vectors(step, x, r,
d);
422 if (res < additional_data.solver_tolerance)
428 state = this->iteration_status(step, res, x);
431 return IterationResult(state, res);
436 if (additional_data.breakdown_testing ==
true &&
437 std::fabs(sigma) < additional_data.breakdown_threshold)
440 const double rho_old = rho;
443 if (additional_data.left_preconditioning)
451 preconditioner.vmult(u, t);
456 const double beta = rho / rho_old;
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
@ failure
Stop iteration, goal not reached.
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
IterationResult iterate(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner, VectorType &r, VectorType &u, VectorType &q, VectorType &t, VectorType &d)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverQMRS(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
AdditionalData additional_data
SolverQMRS(SolverControl &cn, const AdditionalData &data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertThrow(cond, exc)
Expression fabs(const Expression &x)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
double breakdown_threshold
AdditionalData(const bool left_preconditioning=false, const double solver_tolerance=1.e-9, const bool breakdown_testing=true, const double breakdown_threshold=1.e-16)
bool left_preconditioning
SolverControl::State state
IterationResult(const SolverControl::State state, const double last_residual)