15#ifndef dealii_solver_bicgstab_h
16#define dealii_solver_bicgstab_h
78template <
typename VectorType = Vector<
double>>
103 const bool exact_residual =
true,
104 const double breakdown =
105 std::numeric_limits<typename VectorType::value_type>::min())
106 : exact_residual(exact_residual)
107 , breakdown(breakdown)
141 template <
typename MatrixType,
typename PreconditionerType>
145 void solve(const MatrixType &A,
148 const PreconditionerType &preconditioner);
154 template <typename MatrixType>
156 criterion(const MatrixType &A,
167 print_vectors(const
unsigned int step,
170 const VectorType &d) const;
191 const unsigned int last_step,
192 const double last_residual);
199 template <
typename MatrixType,
typename PreconditionerType>
204 const PreconditionerType &preconditioner,
205 const unsigned int step);
215template <
typename VectorType>
218 const bool breakdown,
220 const unsigned int last_step,
221 const double last_residual)
222 : breakdown(breakdown)
224 , last_step(last_step)
225 , last_residual(last_residual)
230template <
typename VectorType>
234 const AdditionalData &data)
236 , additional_data(data)
241template <
typename VectorType>
244 const AdditionalData &data)
246 , additional_data(data)
251template <
typename VectorType>
253template <
typename MatrixType>
260 return std::sqrt(t.add_and_dot(-1.0, b, t));
265template <
typename VectorType>
270 const VectorType &)
const
275template <
typename VectorType>
277template <
typename MatrixType,
typename PreconditionerType>
282 const PreconditionerType &preconditioner,
283 const unsigned int last_step)
296 VectorType &rbar = *Vrbar;
304 rbar.reinit(x,
true);
311 using value_type =
typename VectorType::value_type;
316 value_type res = r.l2_norm();
318 unsigned int step = last_step;
322 return IterationResult(
false, state, step, res);
326 value_type alpha = 1.;
328 value_type omega = 1.;
334 const value_type rhobar = (step == 1 + last_step) ? res * res : r * rbar;
336 if (std::fabs(rhobar) < additional_data.breakdown)
338 return IterationResult(
true, state, step, res);
341 const value_type beta = rhobar * alpha / (rho * omega);
343 if (step == last_step + 1)
350 p.add(-beta * omega, v);
353 preconditioner.vmult(y, p);
355 const value_type rbar_dot_v = rbar * v;
356 if (std::fabs(rbar_dot_v) < additional_data.breakdown)
358 return IterationResult(
true, state, step, res);
361 alpha = rho / rbar_dot_v;
363 res =
std::sqrt(real_type(r.add_and_dot(-alpha, v, r)));
374 print_vectors(step, x, r, y);
378 preconditioner.vmult(z, r);
380 const value_type t_dot_r = t * r;
381 const real_type t_squared = t * t;
382 if (t_squared < additional_data.breakdown)
384 return IterationResult(
true, state, step, res);
386 omega = t_dot_r / t_squared;
387 x.add(alpha, y, omega, z);
389 if (additional_data.exact_residual)
392 res = criterion(A, x, b, t);
395 res =
std::sqrt(real_type(r.add_and_dot(-omega, t, r)));
397 state = this->iteration_status(step, res, x);
398 print_vectors(step, x, r, y);
402 return IterationResult(
false, state, step, res);
407template <
typename VectorType>
409template <
typename MatrixType,
typename PreconditionerType>
416 const PreconditionerType &preconditioner)
423 state = iterate(A, x, b, preconditioner, state.last_step);
430 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
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_CXX20_REQUIRES(condition)
#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)