16 #ifndef dealii_solver_bicgstab_h 17 #define dealii_solver_bicgstab_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/logstream.h> 23 #include <deal.II/base/signaling_nan.h> 24 #include <deal.II/base/subscriptor.h> 26 #include <deal.II/lac/solver.h> 27 #include <deal.II/lac/solver_control.h> 31 DEAL_II_NAMESPACE_OPEN
123 template <
typename VectorType = Vector<
double>>
183 template <
typename MatrixType,
typename PreconditionerType>
185 solve(
const MatrixType & A,
187 const VectorType & b,
188 const PreconditionerType &preconditioner);
234 const VectorType *
Vb;
239 template <
typename MatrixType>
241 criterion(
const MatrixType &A,
const VectorType &x,
const VectorType &b);
250 const VectorType & x,
251 const VectorType & r,
252 const VectorType & d)
const;
263 template <
typename MatrixType>
265 start(
const MatrixType &A);
275 unsigned int last_step;
276 double last_residual;
280 const unsigned int last_step,
281 const double last_residual);
288 template <
typename MatrixType,
typename PreconditionerType>
290 iterate(
const MatrixType &A,
const PreconditionerType &preconditioner);
299 template <
typename VectorType>
301 const bool breakdown,
303 const unsigned int last_step,
304 const double last_residual)
305 : breakdown(breakdown)
307 , last_step(last_step)
308 , last_residual(last_residual)
313 template <
typename VectorType>
316 const AdditionalData & data)
325 template <
typename VectorType>
327 const AdditionalData &data)
336 template <
typename VectorType>
337 template <
typename MatrixType>
352 template <
typename VectorType>
353 template <
typename MatrixType>
358 Vr->sadd(-1., 1., *
Vb);
366 template <
typename VectorType>
371 const VectorType &)
const 376 template <
typename VectorType>
377 template <
typename MatrixType,
typename PreconditionerType>
380 const PreconditionerType &preconditioner)
388 VectorType &rbar = *
Vrbar;
416 preconditioner.vmult(y, p);
424 if (std::fabs(
alpha) > 1.e10)
425 return IterationResult(
true, state,
step,
res);
427 res = std::sqrt(r.add_and_dot(-
alpha, v, r));
442 preconditioner.vmult(z, r);
454 res = std::sqrt(r.add_and_dot(-
omega, t, r));
460 return IterationResult(
false, state,
step,
res);
464 template <
typename VectorType>
465 template <
typename MatrixType,
typename PreconditionerType>
469 const VectorType & b,
470 const PreconditionerType &preconditioner)
482 Vrbar->reinit(x,
true);
500 deallog <<
"Restart step " <<
step << std::endl;
506 state =
iterate(A, preconditioner);
509 while (state.breakdown ==
true);
514 state.last_residual));
520 DEAL_II_NAMESPACE_CLOSE
Stop iteration, goal not reached.
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
virtual ~SolverBicgstab() override=default
VectorMemory< VectorType >::Pointer Vrbar
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate), StateCombiner > iteration_status
#define AssertThrow(cond, exc)
VectorMemory< VectorType >::Pointer Vv
VectorMemory< VectorType >::Pointer Vp
VectorMemory< VectorType >::Pointer Vr
Stop iteration, goal reached.
AdditionalData additional_data
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
VectorMemory< VectorType >::Pointer Vt
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
AdditionalData(const bool exact_residual=true, const double breakdown=1.e-10)
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
IterationResult iterate(const MatrixType &A, const PreconditionerType &preconditioner)
VectorMemory< VectorType >::Pointer Vy
VectorMemory< VectorType >::Pointer Vz
SolverControl::State start(const MatrixType &A)
VectorMemory< VectorType > & memory
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)