16 #ifndef dealii_solver_bicgstab_h 17 #define dealii_solver_bicgstab_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/logstream.h> 22 #include <deal.II/lac/solver.h> 23 #include <deal.II/lac/solver_control.h> 25 #include <deal.II/base/subscriptor.h> 26 #include <deal.II/base/signaling_nan.h> 28 DEAL_II_NAMESPACE_OPEN
120 template <
typename VectorType = Vector<
double> >
181 template <
typename MatrixType,
typename PreconditionerType>
183 solve (
const MatrixType &A,
186 const PreconditionerType &preconditioner);
232 const VectorType *
Vb;
237 template <
typename MatrixType>
238 double criterion (
const MatrixType &A,
const VectorType &x,
const VectorType &b);
248 const VectorType &d)
const;
259 template <
typename MatrixType>
270 unsigned int last_step;
271 double last_residual;
275 const unsigned int last_step,
276 const double last_residual);
283 template <
typename MatrixType,
typename PreconditionerType>
286 const PreconditionerType &preconditioner);
295 template <
typename VectorType>
297 (
const bool breakdown,
299 const unsigned int last_step,
300 const double last_residual)
302 breakdown (breakdown),
304 last_step (last_step),
305 last_residual (last_residual)
310 template <
typename VectorType>
313 const AdditionalData &data)
315 Solver<VectorType>(cn,mem),
323 template <
typename VectorType>
325 const AdditionalData &data)
335 template <
typename VectorType>
341 template <
typename VectorType>
342 template <
typename MatrixType>
357 template <
typename VectorType >
358 template <
typename MatrixType>
363 Vr->sadd(-1.,1.,*
Vb);
371 template <
typename VectorType>
376 const VectorType &)
const 381 template <
typename VectorType>
382 template <
typename MatrixType,
typename PreconditionerType>
385 const PreconditionerType &preconditioner)
392 VectorType &rbar = *
Vrbar;
420 preconditioner.vmult(y,p);
428 if (std::fabs(
alpha) > 1.e10)
429 return IterationResult(
true, state,
step,
res);
431 res = std::sqrt(r.add_and_dot(-
alpha, v, r));
446 preconditioner.vmult(z,r);
458 res = std::sqrt(r.add_and_dot(-
omega, t, r));
464 return IterationResult(
false, state,
step,
res);
468 template <
typename VectorType>
469 template <
typename MatrixType,
typename PreconditionerType>
474 const PreconditionerType &preconditioner)
486 Vrbar->reinit(x,
true);
504 deallog <<
"Restart step " <<
step << std::endl;
510 state =
iterate(A, preconditioner);
513 while (state.breakdown ==
true);
518 state.last_residual));
524 DEAL_II_NAMESPACE_CLOSE
Stop iteration, goal not reached.
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
VectorMemory< VectorType >::Pointer Vrbar
#define AssertThrow(cond, exc)
virtual ~SolverBicgstab()
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)
VectorMemory< VectorType > & memory
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate), StateCombiner > iteration_status
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)
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)