Caribou
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SofaCaribou::solver::ConjugateGradientSolver Class Reference

Detailed Description

Implements a Conjugate Gradient method for solving a linear system of the type Ax = b.

When using no preconditioner, no actual system is built by this component, meaning that the matrix A, and the vectors x and b are never accumulated into memory. Instead, this solver relies on the vectors stored in the mechanical objects of the current scene context graph. The multiplication of the matrix A and a given vector x is computed by accumulating the result b=Ax from ff.addDForce(x, b) of every forcefields acting on a given mechanical object.

When using a preconditioner, the matrix has to be assembled since the preconditioner needs to factorize it. In this case, the complete system matrix A and dense vector b are first accumulated from the mechanical objects of the current scene context graph. Once the dense vector x is found, it is propagated back to the mechanical object's vectors.

#include <ConjugateGradientSolver.h>

Inheritance diagram for SofaCaribou::solver::ConjugateGradientSolver:
SofaCaribou::solver::LinearSolver

Public Types

enum  PreconditioningMethod : unsigned int { PreconditioningMethod::None = 0, PreconditioningMethod::Identity = 1, PreconditioningMethod::Diagonal = 2, PreconditioningMethod::IncompleteLU = 5 }
 Preconditioning methods. More...
 
using SparseMatrix = Eigen::SparseMatrix< FLOATING_POINT_TYPE, Eigen::ColMajor >
 
using Vector = Eigen::Matrix< FLOATING_POINT_TYPE, Eigen::Dynamic, 1 >
 

Public Member Functions

 SOFA_CLASS (ConjugateGradientSolver, LinearSolver)
 
void resetSystem () final
 Reset the complete system (A, x and b are cleared). More...
 
void setSystemMBKMatrix (const sofa::core::MechanicalParams *mparams) final
 Set the linear system matrix A = (mM + bB + kK), storing the coefficients m, b and k of the mechanical M,B,K matrices. More...
 
void setSystemRHVector (sofa::core::MultiVecDerivId b_id) final
 Gives the identifier of the right-hand side vector b. More...
 
void setSystemLHVector (sofa::core::MultiVecDerivId x_id) final
 Gives the identifier of the left-hand side vector x. More...
 
void solveSystem () final
 Solves the system by the conjugate gradient method using the coefficients m, b and k; and the vectors x and b.
 
auto squared_residuals () const -> const std::vector< FLOATING_POINT_TYPE > &
 List of squared residual norms (||r||^2) of every CG iterations of the last solve call.
 
auto squared_initial_residual () const -> const FLOATING_POINT_TYPE &
 Squared residual norm (||r||^2) of the last right-hand side term (b in Ax=b) of the last solve call.
 
auto system_matrix () const -> const SparseMatrix &
 Get the system matrix.
 
void assemble (const sofa::core::MechanicalParams *mparams)
 Assemble the system matrix A = (mM + bB + kK). More...
 
- Public Member Functions inherited from SofaCaribou::solver::LinearSolver
virtual sofa::defaulttype::BaseMatrix * create_new_matrix (unsigned int rows, unsigned int cols) const =0
 Creates a new BaseMatrix of size rows x cols. More...
 
virtual sofa::defaulttype::BaseVector * create_new_vector (unsigned int n) const =0
 Creates a new BaseVector of size n. More...
 
virtual bool solve (const sofa::defaulttype::BaseMatrix *A, const sofa::defaulttype::BaseVector *F, sofa::defaulttype::BaseVector *X) const =0
 Solve the linear system A [X] = F. More...
 

Protected Member Functions

 ConjugateGradientSolver ()
 Constructor.
 
void solve (sofa::core::behavior::MultiVecDeriv &b, sofa::core::behavior::MultiVecDeriv &x)
 Methods. More...
 
template<typename Matrix , typename Preconditioner >
void solve (const Preconditioner &precond, const Matrix &A, const Vector &b, Vector &x)
 Solve the linear system Ax = b using a preconditioner. More...
 

Protected Attributes

Data< bool > d_verbose
 INPUTS.
 
Data< unsigned int > d_maximum_number_of_iterations
 
Data< FLOATING_POINT_TYPE > d_residual_tolerance_threshold
 
Data< sofa::helper::OptionsGroup > d_preconditioning_method
 

Member Enumeration Documentation

◆ PreconditioningMethod

Preconditioning methods.

Enumerator
None 

No preconditioning, hence the complete matrix won't be built. (default)

Identity 

A naive preconditioner which approximates any matrix as the identity matrix.

Diagonal 

Preconditioning using an approximation of A.x = b by ignoring all off-diagonal entries of A.

IncompleteLU 

Preconditioning based on the incomplete LU factorization.

Member Function Documentation

◆ assemble()

void SofaCaribou::solver::ConjugateGradientSolver::assemble ( const sofa::core::MechanicalParams *  mparams)

Assemble the system matrix A = (mM + bB + kK).

Parameters
mparamsMechanical parameters containing the m, b and k factors.

◆ resetSystem()

void SofaCaribou::solver::ConjugateGradientSolver::resetSystem ( )
final

Reset the complete system (A, x and b are cleared).

When using no preconditioner (None), this does absolutely nothing here since the complete system is never built.

This method is called by the MultiMatrix::clear() and MultiMatrix::reset() methods.

◆ setSystemLHVector()

void SofaCaribou::solver::ConjugateGradientSolver::setSystemLHVector ( sofa::core::MultiVecDerivId  x_id)
final

Gives the identifier of the left-hand side vector x.

This identifier will be used to find the actual vector in the mechanical objects of the system. When using a preconditioner (other than None), the complete dense vector is accumulated from the mechanical objects found in the graph subtree of the current context.

◆ setSystemMBKMatrix()

void SofaCaribou::solver::ConjugateGradientSolver::setSystemMBKMatrix ( const sofa::core::MechanicalParams *  mparams)
final

Set the linear system matrix A = (mM + bB + kK), storing the coefficients m, b and k of the mechanical M,B,K matrices.

When using no preconditioner (None), no actual matrix is built by this method since it is not required by the Conjugate Gradient algorithm. Only the coefficients m, b and k are stored.

When using a preconditioner, the complete system A is accumulated into a sparse matrix, and the preconditioner factorize this resulting matrix for later use during the solve step.

Parameters
mparamsContains the coefficients m, b and k of the matrices M, B and K

◆ setSystemRHVector()

void SofaCaribou::solver::ConjugateGradientSolver::setSystemRHVector ( sofa::core::MultiVecDerivId  b_id)
final

Gives the identifier of the right-hand side vector b.

This identifier will be used to find the actual vector in the mechanical objects of the system.When using a preconditioner (other than None), the complete dense vector is accumulated from the mechanical objects found in the graph subtree of the current context.

◆ solve() [1/2]

template<typename Matrix , typename Preconditioner >
void SofaCaribou::solver::ConjugateGradientSolver::solve ( const Preconditioner &  precond,
const Matrix &  A,
const Vector &  b,
Vector &  x 
)
protected

Solve the linear system Ax = b using a preconditioner.

Template Parameters
DerivedThe type of the matrix A, can be a dense or a sparse matrix.
PreconditionerThe type of the preconditioner.
Parameters
precondThe preconditioner (must be a
AThe system matrix as an Eigen matrix
bThe right-hand side vector of the system
xThe solution vector of the system. It should be filled with an initial guess or the previous solution.

◆ solve() [2/2]

void SofaCaribou::solver::ConjugateGradientSolver::solve ( sofa::core::behavior::MultiVecDeriv &  b,
sofa::core::behavior::MultiVecDeriv &  x 
)
protected

Methods.

Solve the linear system Ax = b using Sofa's graph scene.

Here, the matrix A is never built. The matrix-vector product Ax is done by going down in the subgraph of the current context and doing the multiplication directly in the vectors of mechanical objects found.

Warning
Since the matrix A is never built, no preconditioning method can be used with this method.
Parameters
bThe right-hand side vector of the system
xThe solution vector of the system. It should be filled with an initial guess or the previous solution.

The documentation for this class was generated from the following files: