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
72 template <
typename VectorType = Vector<
double> >
132 template<
typename MatrixType,
typename PreconditionerType>
134 solve (
const MatrixType &A,
137 const PreconditionerType &precondition);
143 template <
typename MatrixType>
144 double criterion (
const MatrixType &A,
const VectorType &x,
const VectorType &b);
154 const VectorType &d)
const;
191 const VectorType *
Vb;
233 template <
typename MatrixType>
244 unsigned int last_step;
245 double last_residual;
249 const unsigned int last_step,
250 const double last_residual);
257 template<
typename MatrixType,
typename PreconditionerType>
260 const PreconditionerType &precondition);
269 template<
typename VectorType>
271 (
const bool breakdown,
273 const unsigned int last_step,
274 const double last_residual)
276 breakdown (breakdown),
278 last_step (last_step),
279 last_residual (last_residual)
283 template<
typename VectorType>
286 const AdditionalData &data)
288 Solver<VectorType>(cn,mem),
294 template<
typename VectorType>
296 const AdditionalData &data)
320 template<
typename VectorType>
326 template <
typename VectorType>
327 template <
typename MatrixType>
342 template <
typename VectorType >
343 template <
typename MatrixType>
348 Vr->sadd(-1.,1.,*
Vb);
356 template<
typename VectorType>
361 const VectorType &)
const 366 template<
typename VectorType>
367 template<
typename MatrixType,
typename PreconditionerType>
370 const PreconditionerType &precondition)
377 VectorType &rbar = *
Vrbar;
405 precondition.vmult(y,p);
413 if (std::fabs(
alpha) > 1.e10)
414 return IterationResult(
true, state,
step,
res);
416 res = std::sqrt(r.add_and_dot(-
alpha, v, r));
431 precondition.vmult(z,r);
443 res = std::sqrt(r.add_and_dot(-
omega, t, r));
449 return IterationResult(
false, state,
step,
res);
453 template<
typename VectorType>
454 template<
typename MatrixType,
typename PreconditionerType>
459 const PreconditionerType &precondition)
461 deallog.
push(
"Bicgstab");
465 Vrbar->reinit(x,
true);
488 deallog <<
"Restart step " <<
step << std::endl;
494 state =
iterate(A, precondition);
497 while (state.breakdown ==
true);
512 state.last_residual));
518 DEAL_II_NAMESPACE_CLOSE
Stop iteration, goal not reached.
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
IterationResult iterate(const MatrixType &A, const PreconditionerType &precondition)
#define AssertThrow(cond, exc)
virtual ~SolverBicgstab()
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
Stop iteration, goal reached.
AdditionalData additional_data
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)
void push(const std::string &text)
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
SolverControl::State start(const MatrixType &A)
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)