16#ifndef dealii_solver_qmrs_h
17#define dealii_solver_qmrs_h
93template <
typename VectorType = Vector<
double>>
177 template <
typename MatrixType,
typename PreconditionerType>
181 const VectorType & b,
182 const PreconditionerType &preconditioner);
191 const VectorType & x,
192 const VectorType & r,
193 const VectorType & d)
const;
219 template <
typename MatrixType,
typename PreconditionerType>
223 const VectorType & b,
224 const PreconditionerType &preconditioner,
243template <
class VectorType>
246 const double last_residual)
248 , last_residual(last_residual)
253template <
class VectorType>
256 const AdditionalData & data)
258 , additional_data(data)
262template <
class VectorType>
264 const AdditionalData &data)
266 , additional_data(data)
270template <
class VectorType>
275 const VectorType &)
const
278template <
class VectorType>
279template <
typename MatrixType,
typename PreconditionerType>
283 const VectorType & b,
284 const PreconditionerType &preconditioner)
314 deallog <<
"Restart step " << step << std::endl;
315 state = iterate(A, x, b, preconditioner, *Vr, *Vu, *Vq, *Vt, *Vd);
326template <
class VectorType>
327template <
typename MatrixType,
typename PreconditionerType>
331 const VectorType & b,
332 const PreconditionerType &preconditioner,
343 double tau, rho, theta = 0;
351 if (additional_data.left_preconditioning)
354 preconditioner.vmult(t, r);
361 preconditioner.vmult(q, t);
380 const double sigma = q * t;
383 if (additional_data.breakdown_testing ==
true &&
384 std::fabs(sigma) < additional_data.breakdown_threshold)
387 const double alpha = rho / sigma;
393 const double theta_old = theta;
396 if (additional_data.left_preconditioning)
399 preconditioner.vmult(t, r);
409 const double psi = 1. / (1. + theta);
413 d.sadd(psi * theta_old, psi * alpha, q);
416 print_vectors(step, x, r, d);
424 if (res < additional_data.solver_tolerance)
430 state = this->iteration_status(step, res, x);
433 return IterationResult(state, res);
438 if (additional_data.breakdown_testing ==
true &&
439 std::fabs(sigma) < additional_data.breakdown_threshold)
442 const double rho_old = rho;
445 if (additional_data.left_preconditioning)
453 preconditioner.vmult(u, t);
458 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)
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)