Description
Base class for iterative solvers aimed at solving complementarity problems arising from QP optimization problems.
The problem is described by a variational inequality VI(Z*x-d,K):
| M -Cq'|*|q|- | f|= |0| | Cq -E | |l| |-b| |c|
with l \(\in\) Y, C \(\in\) Ny, normal cone to Y
Also Z symmetric by flipping sign of l_i:
|M Cq'|*| q|-| f|=|0| |Cq E | |-l| |-b| |c|
- case linear problem: all Y_i = R, Ny=0, ex. all bilaterals
- case LCP: all Y_i = R+: c>=0, l>=0, l*c=0
- case CCP: Y_i are friction cones
The default maximum number of iterations is 50.
The 'tolerance' threshold value used in the stopping criteria may have different meaning for different solvers (threshold on maximum constraint violation, threshold on projected residual, etc). See the GetError method for an individual iterative VI solver for details.
Diagonal preconditioning is enabled by default, but may not supported by all iterative VI solvers.
#include <ChIterativeSolverVI.h>
Public Member Functions | |
void | SetOmega (double mval) |
Set the overrelaxation factor (default: 1.0). More... | |
void | SetSharpnessLambda (double mval) |
Set the sharpness factor (default: 1.0). More... | |
void | SetRecordViolation (bool mval) |
Enable/disable recording of the constraint violation history. More... | |
double | GetOmega () const |
Return the current value of the overrelaxation factor. | |
double | GetSharpnessLambda () const |
Return the current value of the sharpness factor. | |
virtual int | GetIterations () const override |
Return the number of iterations performed during the last solve. | |
const std::vector< double > & | GetViolationHistory () const |
Access the vector of constraint violation history. More... | |
const std::vector< double > & | GetDeltalambdaHistory () const |
Access the vector with history of maximum change in Lagrange multipliers Note that collection of constraint violations must be enabled through SetRecordViolation. | |
Public Member Functions inherited from chrono::ChIterativeSolver | |
void | SetMaxIterations (int max_iterations) |
Set the maximum number of iterations. | |
void | SetTolerance (double tolerance) |
Set the tolerance threshold used by the stopping criteria. | |
void | EnableDiagonalPreconditioner (bool val) |
Enable/disable use of a simple diagonal preconditioner (default: true). More... | |
void | EnableWarmStart (bool val) |
Enable/disable warm starting by providing an initial guess (default: false). More... | |
int | GetMaxIterations () const |
Get the current maximum number of iterations. | |
double | GetTolerance () const |
Get the current tolerance value. | |
virtual double | GetError () const =0 |
Return the tolerance error reached during the last solve. | |
Public Member Functions inherited from chrono::ChSolver | |
virtual Type | GetType () const |
Return type of the solver. | |
virtual double | Solve (ChSystemDescriptor &sysd)=0 |
Performs the solution of the problem. More... | |
virtual bool | Setup (ChSystemDescriptor &sysd) |
This function does the setup operations for the solver. More... | |
void | SetVerbose (bool mv) |
Set verbose output from solver. | |
Protected Member Functions | |
void | AtIterationEnd (double mmaxviolation, double mdeltalambda, unsigned int iternum) |
This method MUST be called by all iterative methods INSIDE their iteration loops (at the end). More... | |
virtual void | ArchiveOUT (ChArchiveOut &marchive) override |
Method to allow serialization of transient data to archives. | |
virtual void | ArchiveIN (ChArchiveIn &marchive) override |
Method to allow de-serialization of transient data from archives. | |
virtual bool | SolveRequiresMatrix () const override |
Indicate whether ot not the Solve() phase requires an up-to-date problem matrix. More... | |
Protected Member Functions inherited from chrono::ChIterativeSolver | |
ChIterativeSolver (int max_iterations, double tolerance, bool use_precond, bool warm_start) | |
void | SaveMatrix (ChSystemDescriptor &sysd) |
double | CheckSolution (ChSystemDescriptor &sysd, const ChVectorDynamic<> &x) |
Protected Attributes | |
int | m_iterations |
total number of iterations performed by the solver | |
double | m_omega |
over-relaxation factor | |
double | m_shlambda |
sharpness factor | |
bool | record_violation_history |
std::vector< double > | violation_history |
std::vector< double > | dlambda_history |
Protected Attributes inherited from chrono::ChIterativeSolver | |
bool | m_use_precond |
use diagonal preconditioning? | |
bool | m_warm_start |
use initial guesss? | |
int | m_max_iterations |
maximum number of iterations | |
double | m_tolerance |
tolerance threshold in stopping criteria | |
Protected Attributes inherited from chrono::ChSolver | |
bool | verbose |
Additional Inherited Members | |
Public Types inherited from chrono::ChSolver | |
enum | Type { Type::PSOR = 0, Type::PSSOR, Type::PJACOBI, Type::PMINRES, Type::BARZILAIBORWEIN, Type::APGD, Type::ADDM, Type::SPARSE_LU, Type::SPARSE_QR, Type::PARDISO_MKL, Type::PARDISO_PROJECT, Type::MUMPS, Type::GMRES, Type::MINRES, Type::BICGSTAB, CUSTOM } |
Available types of solvers. More... | |
Member Function Documentation
◆ AtIterationEnd()
|
protected |
This method MUST be called by all iterative methods INSIDE their iteration loops (at the end).
If history recording is enabled, this function will store the current values as passed as arguments. Note: 'iternum' starts at 0 for the first iteration.
◆ GetViolationHistory()
|
inline |
Access the vector of constraint violation history.
Note that collection of constraint violations must be enabled through SetRecordViolation.
◆ SetOmega()
void chrono::ChIterativeSolverVI::SetOmega | ( | double | mval | ) |
Set the overrelaxation factor (default: 1.0).
This factor may be used by PSOR-like methods. A good value for Jacobi solver is 0.2; for other iterative solvers it can be up to 1.0
◆ SetRecordViolation()
|
inline |
Enable/disable recording of the constraint violation history.
If enabled, the maximum constraint violation at the end of each iteration is stored in a vector (see GetViolationHistory).
◆ SetSharpnessLambda()
void chrono::ChIterativeSolverVI::SetSharpnessLambda | ( | double | mval | ) |
Set the sharpness factor (default: 1.0).
This factor may be used by PSOR-like methods with projection (see Mangasarian LCP method). A good sharpness value is in the 0.8 ... 1.0 range (lower values improve accuracy but at the cost of slower convergence)
◆ SolveRequiresMatrix()
|
inlineoverrideprotectedvirtual |
Indicate whether ot not the Solve() phase requires an up-to-date problem matrix.
Typically, this is the case for iterative solvers (as the matrix is needed for the matrix-vector operations).
Implements chrono::ChSolver.
The documentation for this class was generated from the following files:
- /builds/uwsbel/chrono/src/chrono/solver/ChIterativeSolverVI.h
- /builds/uwsbel/chrono/src/chrono/solver/ChIterativeSolverVI.cpp