chrono::ChSolverPMINRES Class Reference

Description

An iterative solver based on modified Krylov iteration of MINRES type with gradient projections (similar to nonlinear CG with Polyak-Ribiere).

The problem is described by a variational inequality VI(Z*x-d,K):

| M -Cq'|*|q|- | f|= |0| , l Y, C Ny, normal cone to Y | Cq -E | |l| |-b| |c|

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

#include <ChSolverPMINRES.h>

Inheritance diagram for chrono::ChSolverPMINRES:
Collaboration diagram for chrono::ChSolverPMINRES:

Public Member Functions

 ChSolverPMINRES (int mmax_iters=50, bool mwarm_start=false, double mtolerance=0.0)
 
virtual Type GetType () const override
 Return type of the solver. More...
 
virtual double Solve (ChSystemDescriptor &sysd) override
 Performs the solution of the problem. More...
 
virtual double Solve_SupportingStiffness (ChSystemDescriptor &sysd)
 Same as Solve(), but this also supports the presence of ChKblock blocks. More...
 
void SetGradDiffStep (double mf)
 For the case where inequalities are introduced, the gradient is projected. More...
 
double GetGradDiffStep ()
 
void SetRelTolerance (double mrt)
 Set relative tolerance. More...
 
double GetRelTolerance ()
 
void SetDiagonalPreconditioning (bool mp)
 Enable diagonal preconditioning. More...
 
bool GetDiagonalPreconditioning ()
 
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.
 
- Public Member Functions inherited from chrono::ChIterativeSolver
 ChIterativeSolver (int mmax_iters=50, bool mwarm_start=false, double mtolerance=0.0, double momega=1.0, double mshlambda=1.0)
 
virtual bool SolveRequiresMatrix () const override
 Indicate whether ot not the Solve() phase requires an up-to-date problem matrix. More...
 
void SetMaxIterations (int mval)
 Set the maximum number of iterations. More...
 
int GetMaxIterations () const
 Get the current value for the maximum allowed number of iterations.
 
int GetTotalIterations () const
 Get the number of total iterations taken by the solver.
 
void SetOmega (double mval)
 Set the overrelaxation factor. More...
 
double GetOmega () const
 Return the current value of the overrelaxation factor.
 
virtual void SetSharpnessLambda (double mval)
 Set the sharpness factor. More...
 
virtual double GetSharpnessLambda () const
 Return the current value of the sharpness factor.
 
void SetWarmStart (bool mval)
 Enable/disable 'warm start'. More...
 
bool GetWarmStart () const
 Return a flag indicating whether or not warm start is enabled.
 
void SetTolerance (double mval)
 Set the tolerance for stopping criterion. More...
 
double GetTolerance () const
 Return the current value of the solver tolerance.
 
void SetRecordViolation (bool mval)
 Enable/disable recording of the constraint violation history. More...
 
bool GetRecordViolation () const
 Return a flag indicating whether or not constraint violations are recorded.
 
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. More...
 
- Public Member Functions inherited from chrono::ChSolver
virtual bool Setup (ChSystemDescriptor &sysd)
 This function does the setup operations for the solver. More...
 
void SetVerbose (bool mv)
 Set verbose output from solver.
 
bool GetVerbose () const
 

Protected Attributes

double grad_diffstep
 
double rel_tolerance
 
bool diag_preconditioning
 
- Protected Attributes inherited from chrono::ChIterativeSolver
int max_iterations
 maximum allowed iterations
 
int tot_iterations
 total number of iterations performed by the solver
 
bool warm_start
 indicate whether or not to use warm starting
 
double tolerance
 tolerance for termination criteria
 
double omega
 over-relaxation factor
 
double shlambda
 sharpness factor
 
bool record_violation_history
 
std::vector< double > violation_history
 
std::vector< double > dlambda_history
 
- Protected Attributes inherited from chrono::ChSolver
bool verbose
 

Additional Inherited Members

- Public Types inherited from chrono::ChSolver
enum  Type {
  SOR = 0, SYMMSOR, JACOBI, SOR_MULTITHREAD,
  PMINRES, BARZILAIBORWEIN, PCG, APGD,
  MINRES, SOLVER_DEM, CUSTOM
}
 Available types of solvers.
 
- Protected Member Functions inherited from chrono::ChIterativeSolver
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...
 

Constructor & Destructor Documentation

chrono::ChSolverPMINRES::ChSolverPMINRES ( int  mmax_iters = 50,
bool  mwarm_start = false,
double  mtolerance = 0.0 
)
Parameters
mmax_itersmax.number of iterations
mwarm_startuses warm start?
mtolerancetolerance for termination criterion

Member Function Documentation

virtual Type chrono::ChSolverPMINRES::GetType ( ) const
overridevirtual

Return type of the solver.

Default is CUSTOM. Derived classes should override this function.

Reimplemented from chrono::ChSolver.

void chrono::ChSolverPMINRES::SetDiagonalPreconditioning ( bool  mp)

Enable diagonal preconditioning.

It a simple but fast preconditioning technique that is expecially useful to fix slow convergence in case variables have very different orders of magnitude.

void chrono::ChSolverPMINRES::SetGradDiffStep ( double  mf)

For the case where inequalities are introduced, the gradient is projected.

A numerical differentiation is used, this is the delta.

void chrono::ChSolverPMINRES::SetRelTolerance ( double  mrt)

Set relative tolerance.

the iteration stops when the (projected) residual r is smaller than absolute tolerance, that you set via SetTolerance(), OR it is smaller than 'rhs_rel_tol' i.e. the norm of the right hand side multiplied by this relative tolerance. Set to 0 if you do not want relative tolerance to enter into play.

double chrono::ChSolverPMINRES::Solve ( ChSystemDescriptor sysd)
overridevirtual

Performs the solution of the problem.

Returns
the maximum constraint violation after termination.
Parameters
sysdsystem description with constraints and variables

Implements chrono::ChSolver.

double chrono::ChSolverPMINRES::Solve_SupportingStiffness ( ChSystemDescriptor sysd)
virtual

Same as Solve(), but this also supports the presence of ChKblock blocks.

If Solve() is called and stiffness is present, Solve() automatically falls back to this function. It does not solve the Schur complement N*l-r=0 as Solve does, here the entire system KKT matrix with duals l and primals q is used.

Parameters
sysdsystem description with constraints and variables