16#ifndef dealii_solver_bicgstab_h
17#define dealii_solver_bicgstab_h
78template <
typename VectorType = Vector<
double>>
104 std::numeric_limits<typename VectorType::value_type>::min())
140 template <
typename MatrixType,
typename PreconditionerType>
144 const VectorType & b,
145 const PreconditionerType &preconditioner);
151 template <
typename MatrixType>
165 const VectorType & x,
166 const VectorType & r,
167 const VectorType & d)
const;
196 template <
typename MatrixType,
typename PreconditionerType>
200 const VectorType & b,
201 const PreconditionerType &preconditioner,
202 const unsigned int step);
212template <
typename VectorType>
214 const bool breakdown,
216 const unsigned int last_step,
217 const double last_residual)
218 : breakdown(breakdown)
220 , last_step(last_step)
221 , last_residual(last_residual)
226template <
typename VectorType>
229 const AdditionalData & data)
231 , additional_data(data)
236template <
typename VectorType>
238 const AdditionalData &data)
240 , additional_data(data)
245template <
typename VectorType>
246template <
typename MatrixType>
254 return std::sqrt(t.add_and_dot(-1.0, b, t));
259template <
typename VectorType>
264 const VectorType &)
const
269template <
typename VectorType>
270template <
typename MatrixType,
typename PreconditionerType>
274 const VectorType & b,
275 const PreconditionerType &preconditioner,
276 const unsigned int last_step)
289 VectorType &rbar = *Vrbar;
297 rbar.reinit(x,
true);
304 using value_type =
typename VectorType::value_type;
309 value_type res = r.l2_norm();
311 unsigned int step = last_step;
315 return IterationResult(
false, state, step, res);
319 value_type alpha = 1.;
321 value_type omega = 1.;
327 const value_type rhobar = (step == 1 + last_step) ? res * res : r * rbar;
329 if (std::fabs(rhobar) < additional_data.breakdown)
331 return IterationResult(
true, state, step, res);
334 const value_type beta = rhobar * alpha / (rho * omega);
336 if (step == last_step + 1)
343 p.add(-beta * omega, v);
346 preconditioner.vmult(y, p);
348 const value_type rbar_dot_v = rbar * v;
349 if (std::fabs(rbar_dot_v) < additional_data.breakdown)
351 return IterationResult(
true, state, step, res);
354 alpha = rho / rbar_dot_v;
356 res =
std::sqrt(real_type(r.add_and_dot(-alpha, v, r)));
367 print_vectors(step, x, r, y);
371 preconditioner.vmult(z, r);
373 const value_type t_dot_r = t * r;
374 const real_type t_squared = t * t;
375 if (t_squared < additional_data.breakdown)
377 return IterationResult(
true, state, step, res);
379 omega = t_dot_r / t_squared;
380 x.add(alpha, y, omega, z);
382 if (additional_data.exact_residual)
385 res = criterion(A, x, b, t);
388 res =
std::sqrt(real_type(r.add_and_dot(-omega, t, r)));
390 state = this->iteration_status(step, res, x);
391 print_vectors(step, x, r, y);
395 return IterationResult(
false, state, step, res);
400template <
typename VectorType>
401template <
typename MatrixType,
typename PreconditionerType>
405 const VectorType & b,
406 const PreconditionerType &preconditioner)
413 state = iterate(A, x, b, preconditioner, state.last_step);
420 state.last_residual));
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b, VectorType &t)
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
virtual ~SolverBicgstab() override=default
AdditionalData additional_data
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
IterationResult iterate(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner, const unsigned int step)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
@ failure
Stop iteration, goal not reached.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertThrow(cond, exc)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
AdditionalData(const bool exact_residual=true, const double breakdown=std::numeric_limits< typename VectorType::value_type >::min())
SolverControl::State state
IterationResult(const bool breakdown, const SolverControl::State state, const unsigned int last_step, const double last_residual)