Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/solver_cg.h>
Classes | |
struct | AdditionalData |
Public Types | |
using | size_type = types::global_dof_index |
Public Types inherited from SolverBase< VectorType > | |
using | vector_type = VectorType |
Public Member Functions | |
SolverCG (SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData()) | |
SolverCG (SolverControl &cn, const AdditionalData &data=AdditionalData()) | |
virtual | ~SolverCG () override=default |
template<typename MatrixType , typename PreconditionerType > | |
void | solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner) |
boost::signals2::connection | connect_coefficients_slot (const std::function< void(typename VectorType::value_type, typename VectorType::value_type)> &slot) |
boost::signals2::connection | connect_condition_number_slot (const std::function< void(double)> &slot, const bool every_iteration=false) |
boost::signals2::connection | connect_eigenvalues_slot (const std::function< void(const std::vector< double > &)> &slot, const bool every_iteration=false) |
Public Member Functions inherited from SolverBase< VectorType > | |
SolverBase (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory) | |
SolverBase (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 (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Member Functions | |
virtual void | print_vectors (const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const |
Static Protected Member Functions | |
static void | compute_eigs_and_cond (const std::vector< typename VectorType::value_type > &diagonal, const std::vector< typename VectorType::value_type > &offdiagonal, const boost::signals2::signal< void(const std::vector< double > &)> &eigenvalues_signal, const boost::signals2::signal< void(double)> &cond_signal) |
Protected Attributes | |
AdditionalData | additional_data |
boost::signals2::signal< void(typename VectorType::value_type, typename VectorType::value_type)> | coefficients_signal |
boost::signals2::signal< void(double)> | condition_number_signal |
boost::signals2::signal< void(double)> | all_condition_numbers_signal |
boost::signals2::signal< void(const std::vector< double > &)> | eigenvalues_signal |
boost::signals2::signal< void(const std::vector< double > &)> | all_eigenvalues_signal |
Protected Attributes inherited from SolverBase< 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 |
Additional Inherited Members | |
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) |
This class implements the preconditioned Conjugate Gradients (CG) method that can be used to solve linear systems with a symmetric positive definite matrix. This class is used first in step-3 and step-4, but is used in many other tutorial programs as well. Like all other solver classes, it can work on any kind of vector and matrix as long as they satisfy certain requirements (for the requirements on matrices and vectors in order to work with this class, see the documentation of the Solver base class). The type of the solution vector must be passed as template argument, and defaults to Vector<double>.
The cg-method performs an orthogonal projection of the original preconditioned linear system to another system of smaller dimension. Furthermore, the projected matrix T
is tri-diagonal. Since the projection is orthogonal, the eigenvalues of T
approximate those of the original preconditioned matrix PA
. In fact, after n
steps, where n
is the dimension of the original system, the eigenvalues of both matrices are equal. But, even for small numbers of iteration steps, the condition number of T
is a good estimate for the one of PA
.
After m
steps the matrix T_m can be written in terms of the coefficients alpha
and beta
as the tri-diagonal matrix with diagonal elements 1/alpha_0
, 1/alpha_1 + beta_0/alpha_0
, ..., 1/alpha_{m-1
+beta_{m-2}/alpha_{m-2}} and off-diagonal elements sqrt(beta_0)/alpha_0
, ..., sqrt(beta_{m-2
)/alpha_{m-2}}. The eigenvalues of this matrix can be computed by postprocessing.
The coefficients, eigenvalues and condition number (computed as the ratio of the largest over smallest eigenvalue) can be obtained by connecting a function as a slot to the solver using one of the functions connect_coefficients_slot
, connect_eigenvalues_slot
and connect_condition_number_slot
. These slots will then be called from the solver with the estimates as argument.
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 96 of file solver_cg.h.
using SolverCG< VectorType >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 102 of file solver_cg.h.
SolverCG< VectorType >::SolverCG | ( | SolverControl & | cn, |
VectorMemory< VectorType > & | mem, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverCG< VectorType >::SolverCG | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
|
overridevirtualdefault |
Virtual destructor.
void SolverCG< VectorType >::solve | ( | const MatrixType & | A, |
VectorType & | x, | ||
const VectorType & | b, | ||
const PreconditionerType & | preconditioner | ||
) |
Solve the linear system \(Ax=b\) for x.
boost::signals2::connection SolverCG< VectorType >::connect_coefficients_slot | ( | const std::function< void(typename VectorType::value_type, typename VectorType::value_type)> & | slot | ) |
Connect a slot to retrieve the CG coefficients. The slot will be called with alpha as the first argument and with beta as the second argument, where alpha and beta follow the notation in Y. Saad: "Iterative methods for Sparse Linear Systems", section 6.7. Called once per iteration
boost::signals2::connection SolverCG< VectorType >::connect_condition_number_slot | ( | const std::function< void(double)> & | slot, |
const bool | every_iteration = false |
||
) |
Connect a slot to retrieve the estimated condition number. Called on each iteration if every_iteration=true, otherwise called once when iterations are ended (i.e., either because convergence has been achieved, or because divergence has been detected).
boost::signals2::connection SolverCG< VectorType >::connect_eigenvalues_slot | ( | const std::function< void(const std::vector< double > &)> & | slot, |
const bool | every_iteration = false |
||
) |
Connect a slot to retrieve the estimated eigenvalues. Called on each iteration if every_iteration=true, otherwise called once when iterations are ended (i.e., either because convergence has been achieved, or because divergence has been detected).
|
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.
|
staticprotected |
Estimates the eigenvalues from diagonal and offdiagonal. Uses these estimate to compute the condition number. Calls the signals eigenvalues_signal and cond_signal with these estimates as arguments.
|
protected |
Additional parameters.
Definition at line 200 of file solver_cg.h.
|
protected |
Signal used to retrieve the CG coefficients. Called on each iteration.
Definition at line 207 of file solver_cg.h.
|
protected |
Signal used to retrieve the estimated condition number. Called once when all iterations are ended.
Definition at line 213 of file solver_cg.h.
|
protected |
Signal used to retrieve the estimated condition numbers. Called on each iteration.
Definition at line 219 of file solver_cg.h.
|
protected |
Signal used to retrieve the estimated eigenvalues. Called once when all iterations are ended.
Definition at line 225 of file solver_cg.h.
|
protected |
Signal used to retrieve the estimated eigenvalues. Called on each iteration.
Definition at line 232 of file solver_cg.h.