Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/solver_bicgstab.h>
Classes | |
struct | AdditionalData |
struct | IterationResult |
Public Member Functions | |
SolverBicgstab (SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData()) | |
SolverBicgstab (SolverControl &cn, const AdditionalData &data=AdditionalData()) | |
virtual | ~SolverBicgstab () |
template<typename MatrixType , typename PreconditionerType > | |
void | solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner) |
Public Member Functions inherited from Solver< VectorType > | |
Solver (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory) | |
Solver (SolverControl &solver_control) | |
boost::signals2::connection | connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate)> &slot) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Member Functions | |
template<typename MatrixType > | |
double | criterion (const MatrixType &A, const VectorType &x, const VectorType &b) |
virtual void | print_vectors (const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const |
Protected Member Functions inherited from internal::SolverBicgstabData | |
SolverBicgstabData () | |
Protected Attributes | |
VectorType * | Vx |
VectorMemory< VectorType >::Pointer | Vr |
VectorMemory< VectorType >::Pointer | Vrbar |
VectorMemory< VectorType >::Pointer | Vp |
VectorMemory< VectorType >::Pointer | Vy |
VectorMemory< VectorType >::Pointer | Vz |
VectorMemory< VectorType >::Pointer | Vt |
VectorMemory< VectorType >::Pointer | Vv |
const VectorType * | Vb |
AdditionalData | additional_data |
Protected Attributes inherited from Solver< VectorType > | |
GrowingVectorMemory< VectorType > | static_vector_memory |
VectorMemory< VectorType > & | memory |
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate), StateCombiner > | iteration_status |
Protected Attributes inherited from internal::SolverBicgstabData | |
double | alpha |
double | beta |
double | omega |
double | rho |
double | rhobar |
unsigned int | step |
double | res |
Private Member Functions | |
template<typename MatrixType > | |
SolverControl::State | start (const MatrixType &A) |
template<typename MatrixType , typename PreconditionerType > | |
IterationResult | iterate (const MatrixType &A, const PreconditionerType &preconditioner) |
Additional Inherited Members | |
Public Types inherited from Solver< VectorType > | |
typedef VectorType | vector_type |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Bicgstab algorithm by van der Vorst.
For the requirements on matrices and vectors in order to work with this class, see the documentation of the Solver base class.
Like all other solver classes, this class has a local structure called AdditionalData
which is used to pass additional parameters to the solver, like damping parameters or the number of temporary vectors. We use this additional structure instead of passing these values directly to the constructor because this makes the use of the SolverSelector
and other classes much easier and guarantees that these will continue to work even if number or type of the additional parameters for a certain solver changes.
The Bicgstab-method has two additional parameters: the first is a boolean, deciding whether to compute the actual residual in each step (true
) or to use the length of the computed orthogonal residual (false
). Note that computing the residual causes a third matrix-vector-multiplication, though no additional preconditioning, in each step. The reason for doing this is, that the size of the orthogonalized residual computed during the iteration may be larger by orders of magnitude than the true residual. This is due to numerical instabilities related to badly conditioned matrices. Since this instability results in a bad stopping criterion, the default for this parameter is true
. Whenever the user knows that the estimated residual works reasonably as well, the flag should be set to false
in order to increase the performance of the solver.
The second parameter is the size of a breakdown criterion. It is difficult to find a general good criterion, so if things do not work for you, try to change this value.
The solve() function of this class uses the mechanism described in the Solver base class to determine convergence. This mechanism can also be used to observe the progress of the iteration.
Definition at line 121 of file solver_bicgstab.h.
SolverBicgstab< VectorType >::SolverBicgstab | ( | SolverControl & | cn, |
VectorMemory< VectorType > & | mem, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverBicgstab< VectorType >::SolverBicgstab | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
|
virtual |
Virtual destructor.
void SolverBicgstab< VectorType >::solve | ( | const MatrixType & | A, |
VectorType & | x, | ||
const VectorType & | b, | ||
const PreconditionerType & | preconditioner | ||
) |
Solve primal problem only.
|
protected |
Computation of the stopping criterion.
|
protectedvirtual |
Interface for derived class. This function gets the current iteration vector, the residual and the update vector in each step. It can be used for graphical output of the convergence history.
|
private |
Everything before the iteration loop.
|
private |
The iteration loop itself. The function returns a structure indicating what happened in this function.
|
protected |
Auxiliary vector.
Definition at line 192 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 197 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 202 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 207 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 212 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 217 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 222 of file solver_bicgstab.h.
|
protected |
Auxiliary vector.
Definition at line 227 of file solver_bicgstab.h.
|
protected |
Right hand side vector.
Definition at line 232 of file solver_bicgstab.h.
|
protected |
Additional parameters.
Definition at line 253 of file solver_bicgstab.h.