16#ifndef dealii_solver_bicgstab_h
17#define dealii_solver_bicgstab_h
122template <
typename VectorType = Vector<
double>>
185 template <
typename MatrixType,
typename PreconditionerType>
189 const VectorType &
b,
190 const PreconditionerType &preconditioner);
236 const VectorType *
Vb;
241 template <
typename MatrixType>
243 criterion(
const MatrixType &
A,
const VectorType &x,
const VectorType &
b);
252 const VectorType & x,
253 const VectorType & r,
254 const VectorType &
d)
const;
283 template <
typename MatrixType,
typename PreconditionerType>
285 iterate(
const MatrixType &
A,
const PreconditionerType &preconditioner);
294template <
typename VectorType>
296 const bool breakdown,
298 const unsigned int last_step,
299 const double last_residual)
300 : breakdown(breakdown)
302 , last_step(last_step)
303 , last_residual(last_residual)
308template <
typename VectorType>
311 const AdditionalData & data)
315 , additional_data(data)
320template <
typename VectorType>
322 const AdditionalData &data)
326 , additional_data(data)
331template <
typename VectorType>
332template <
typename MatrixType>
347template <
typename VectorType>
352 const VectorType &)
const
357template <
typename VectorType>
358template <
typename MatrixType,
typename PreconditionerType>
361 const PreconditionerType &preconditioner)
364 Vr->sadd(-1., 1., *Vb);
368 if (state == SolverControl::State::success)
369 return IterationResult(
false, state, step, res);
371 alpha = omega = rho = 1.;
374 VectorType &rbar = *Vrbar;
389 if (
std::fabs(rhobar) < additional_data.breakdown)
391 return IterationResult(
true, state, step, res);
393 beta = rhobar * alpha / (rho * omega);
403 p.add(-beta * omega, v);
406 preconditioner.vmult(y, p);
409 if (
std::fabs(rhobar) < additional_data.breakdown)
411 return IterationResult(
true, state, step, res);
414 alpha = rho / rhobar;
416 res =
std::sqrt(r.add_and_dot(-alpha, v, r));
427 print_vectors(step, *Vx, r, y);
431 preconditioner.vmult(z, r);
434 auto t_squared = t * t;
435 if (t_squared < additional_data.breakdown)
437 return IterationResult(
true, state, step, res);
439 omega = rhobar / (t * t);
440 Vx->add(alpha, y, omega, z);
442 if (additional_data.exact_residual)
445 res = criterion(
A, *Vx, *Vb);
448 res =
std::sqrt(r.add_and_dot(-omega, t, r));
450 state = this->iteration_status(step, res, *Vx);
451 print_vectors(step, *Vx, r, y);
455 return IterationResult(
false, state, step, res);
460template <
typename VectorType>
461template <
typename MatrixType,
typename PreconditionerType>
465 const VectorType &
b,
466 const PreconditionerType &preconditioner)
480 Vrbar->reinit(x,
true);
495 state = iterate(
A, preconditioner);
512 state.last_residual));
VectorMemory< VectorType >::Pointer Vr
VectorMemory< VectorType >::Pointer Vp
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
VectorMemory< VectorType >::Pointer Vv
VectorMemory< VectorType >::Pointer Vy
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)
VectorMemory< VectorType >::Pointer Vrbar
VectorMemory< VectorType >::Pointer Vt
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
IterationResult iterate(const MatrixType &A, const PreconditionerType &preconditioner)
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
VectorMemory< VectorType >::Pointer Vz
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)
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 > &)
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)