<ConjugateGradientSolver />

Doxygen: SofaCaribou::solver::ConjugateGradientSolver

Implementation of a Conjugate Gradient linear solver for selfadjoint (symmetric) positive-definite matrices.

Attribute

Format

Default

Description

printLog

bool

false

Output informative messages at the initialization and during the simulation.

verbose

bool

false

Output convergence status at each iterations of the CG.

maximum_number_of_iterations

int

25

Maximum number of iterations before diverging.

residual_tolerance_threshold

float

1e-5

Convergence criterion: The CG iterations will stop when the residual norm of the residual \(\frac{|r_{k}|}{|r_0|} = \frac{|r_{k} - a_k A p_k|}{|b|}\) at iteration k is lower than this threshold (here \(b\) is the right-hand side vector).

preconditioning_method

option

None

Preconditioning method used.

  • None: No preconditioning, hence the complete matrix won’t be built. (default)

  • Identity: A naive preconditioner which approximates any matrix as the identity matrix. This is equivalent as using None, except the system matrix is built.

  • Diagonal: Preconditioning using an approximation of A.x = b by ignoring all off-diagonal entries of A. Also called Jacobi preconditioner, work very well on diagonal dominant matrices. This preconditioner is suitable for both selfadjoint and general problems. The diagonal entries are pre-inverted and stored into a dense vector. See here for more details.

  • IncompleteCholesky: Preconditioning based on an modified incomplete Cholesky with dual threshold. See here for more details.

  • IncompleteLU: Preconditioning based on the incomplete LU factorization. See here for more details.

Quick example

XML

<Node>
    <StaticODESolver newton_iterations="10" correction_tolerance_threshold="1e-8" residual_tolerance_threshold="1e-8" printLog="1" />
    <ConjugateGradientSolver maximum_number_of_iterations="2500" residual_tolerance_threshold="1e-12" preconditioning_method="Diagonal" printLog="0" />
</Node>

Python

node.addObject('StaticODESolver', newton_iterations=10, correction_tolerance_threshold=1e-8, residual_tolerance_threshold=1e-8, printLog=True)
node.addObject('ConjugateGradientSolver', maximum_number_of_iterations=2500, residual_tolerance_threshold=1e-12, preconditioning_method="Diagonal", printLog=False)

Available python bindings

class ConjugateGradientSolver
K()
Returns

Reference to the system matrix

Return type

scipy.sparse.csc_matrix

Note

No copy involved.

Get the system matrix A = (mM + bB + kK) as a compressed sparse column major matrix.

assemble(m, b, k)
Parameters
  • m (float) – Factor for the mass (M) matrix.

  • b (float) – Factor for the damping (b) matrix.

  • k (float) – Factor for the stiffness (K) matrix.

Returns

Reference to the system matrix

Return type

scipy.sparse.csc_matrix

Note

No copy involved.

Get the system matrix A = (mM + bB + kK) as a compressed sparse column major matrix.