autogenu-jupyter
An automatic code generator and the continuation/GMRES (C/GMRES) based numerical solvers for nonlinear MPC
Loading...
Searching...
No Matches
cgmres::SingleShootingCGMRESSolver< OCP, N, kmax > Class Template Reference

Single-shooting C/GMRES solver for nonlinear MPC. More...

#include <single_shooting_cgmres_solver.hpp>

Public Types

using SingleShootingNLP_ = detail::SingleShootingNLP< OCP, N >
 
using ContinuationGMRES_ = detail::ContinuationGMRES< SingleShootingNLP_ >
 
using MatrixFreeGMRES_ = detail::MatrixFreeGMRES< ContinuationGMRES_, kmax >
 

Public Member Functions

 SingleShootingCGMRESSolver (const OCP &ocp, const Horizon &horizon, const SolverSettings &settings)
 Constructs the single-shooting C/GMRES solver. More...
 
 SingleShootingCGMRESSolver ()=default
 Default constructor. More...
 
 ~SingleShootingCGMRESSolver ()=default
 Default destructor. More...
 
template<typename VectorType >
void set_u (const MatrixBase< VectorType > &u)
 Sets the control input vector. More...
 
template<typename VectorType >
void set_uc (const MatrixBase< VectorType > &uc)
 Sets the control input vector and Lagrange multiplier with respect to the equality constraints. More...
 
template<typename VectorType >
void set_dummy (const MatrixBase< VectorType > &dummy)
 Sets the dummy input vector with respect to the control input bounds constraint. More...
 
template<typename VectorType >
void set_mu (const MatrixBase< VectorType > &mu)
 Sets the Lagrange multiplier with respect to the control input bounds constraint. More...
 
template<typename T >
void set_u_array (const T &u_array)
 Sets the array of control input vector. More...
 
template<typename T >
void set_uc_array (const T &uc_array)
 Sets the array of control input vector and Lagrange multiplier with respect to the equality constraints. More...
 
template<typename T >
void set_dummy_array (const T &dummy_array)
 Sets the array of the dummy input vector with respect to the control input bounds constraint. More...
 
template<typename T >
void set_mu_array (const T &mu_array)
 Sets the array of the Lagrange multiplier with respect to the control input bounds constraint. More...
 
void init_dummy_mu ()
 Initializes the dummy input vectors and Lagrange multipliers with respect to the control input bounds constraint. More...
 
const std::array< Vector< nu >, N > & uopt () const
 Getter of the optimal solution. More...
 
const std::array< Vector< nuc >, N > & ucopt () const
 Getter of the optimal solution. More...
 
const std::array< Vector< nx >, N+1 > & xopt () const
 Getter of the optimal solution. More...
 
const std::array< Vector< nx >, N+1 > & lmdopt () const
 Getter of the optimal solution. More...
 
const std::array< Vector< nub >, N > & dummyopt () const
 Getter of the optimal solution. More...
 
const std::array< Vector< nub >, N > & muopt () const
 Getter of the optimal solution. More...
 
Scalar optError () const
 Gets the l2-norm of the current optimality errors. More...
 
template<typename VectorType >
Scalar optError (const Scalar t, const MatrixBase< VectorType > &x)
 Computes and gets the l2-norm of the current optimality errors. More...
 
template<typename VectorType >
void update (const Scalar t, const MatrixBase< VectorType > &x)
 Updates the solution by performing C/GMRES method. More...
 
TimingProfile getProfile () const
 Get timing result as TimingProfile. More...
 
void disp (std::ostream &os) const
 

Static Public Attributes

static constexpr int nx = OCP::nx
 Dimension of the state. More...
 
static constexpr int nu = OCP::nu
 Dimension of the control input. More...
 
static constexpr int nc = OCP::nc
 Dimension of the equality constraints. More...
 
static constexpr int nuc = nu + nc
 Dimension of the concatenation of the control input and equality constraints. More...
 
static constexpr int nub = OCP::nub
 Dimension of the bound constraints on the control input. More...
 
static constexpr int dim = nuc * N + 2 * N * nub
 Dimension of the linear problem solved by the GMRES solver. More...
 

Friends

std::ostream & operator<< (std::ostream &os, const SingleShootingCGMRESSolver &solver)
 

Detailed Description

template<class OCP, int N, int kmax>
class cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >

Single-shooting C/GMRES solver for nonlinear MPC.

Template Parameters
OCPA definition of the optimal control problem (OCP).
NNumber of discretizationn grids of the horizon. Must be positive.
kmaxMaximum number of the GMRES iterations. Must be positive.

Member Typedef Documentation

◆ ContinuationGMRES_

template<class OCP , int N, int kmax>
using cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::ContinuationGMRES_ = detail::ContinuationGMRES<SingleShootingNLP_>

◆ MatrixFreeGMRES_

template<class OCP , int N, int kmax>
using cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::MatrixFreeGMRES_ = detail::MatrixFreeGMRES<ContinuationGMRES_, kmax>

◆ SingleShootingNLP_

template<class OCP , int N, int kmax>
using cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::SingleShootingNLP_ = detail::SingleShootingNLP<OCP, N>

Constructor & Destructor Documentation

◆ SingleShootingCGMRESSolver() [1/2]

template<class OCP , int N, int kmax>
cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::SingleShootingCGMRESSolver ( const OCP &  ocp,
const Horizon horizon,
const SolverSettings settings 
)
inline

Constructs the single-shooting C/GMRES solver.

Parameters
[in]ocpA definition of the optimal control problem (OCP).
[in]horizonPrediction horizon of MPC.
[in]settingsSolver settings.

◆ SingleShootingCGMRESSolver() [2/2]

template<class OCP , int N, int kmax>
cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::SingleShootingCGMRESSolver ( )
default

Default constructor.

◆ ~SingleShootingCGMRESSolver()

template<class OCP , int N, int kmax>
cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::~SingleShootingCGMRESSolver ( )
default

Default destructor.

Member Function Documentation

◆ disp()

template<class OCP , int N, int kmax>
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::disp ( std::ostream &  os) const
inline

◆ dummyopt()

template<class OCP , int N, int kmax>
const std::array< Vector< nub >, N > & cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::dummyopt ( ) const
inline

Getter of the optimal solution.

Returns
const reference to the optimal dummy input vectors with respect to the control input bounds constraint over the horizon.

◆ getProfile()

template<class OCP , int N, int kmax>
TimingProfile cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::getProfile ( ) const
inline

Get timing result as TimingProfile.

Returns
Timing profile.

◆ init_dummy_mu()

template<class OCP , int N, int kmax>
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::init_dummy_mu ( )
inline

Initializes the dummy input vectors and Lagrange multipliers with respect to the control input bounds constraint.

◆ lmdopt()

template<class OCP , int N, int kmax>
const std::array< Vector< nx >, N+1 > & cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::lmdopt ( ) const
inline

Getter of the optimal solution.

Returns
const reference to the optimal costate vectors over the horizon.

◆ muopt()

template<class OCP , int N, int kmax>
const std::array< Vector< nub >, N > & cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::muopt ( ) const
inline

Getter of the optimal solution.

Returns
const reference to the Lagrange multipliers with respect to the control input bounds constraint over the horizon.

◆ optError() [1/2]

template<class OCP , int N, int kmax>
Scalar cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::optError ( ) const
inline

Gets the l2-norm of the current optimality errors.

Returns
The l2-norm of the current optimality errors.

◆ optError() [2/2]

template<class OCP , int N, int kmax>
template<typename VectorType >
Scalar cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::optError ( const Scalar  t,
const MatrixBase< VectorType > &  x 
)
inline

Computes and gets the l2-norm of the current optimality errors.

Parameters
[in]tInitial time of the horizon.
[in]xInitial state of the horizon. Size must be SingleShootingCGMRESSolver::nx.
Returns
The l2-norm of the current optimality errors.

◆ set_dummy()

template<class OCP , int N, int kmax>
template<typename VectorType >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_dummy ( const MatrixBase< VectorType > &  dummy)
inline

Sets the dummy input vector with respect to the control input bounds constraint.

Parameters
[in]dummyThe dummy input vector. Size must be SingleShootingCGMRESSolver::nub.

◆ set_dummy_array()

template<class OCP , int N, int kmax>
template<typename T >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_dummy_array ( const T &  dummy_array)
inline

Sets the array of the dummy input vector with respect to the control input bounds constraint.

Parameters
[in]dummy_arrayA container (std::vector, std::array, etc.) of the dummy input vector. Size must be N and size of each element must be SingleShootingCGMRESSolver::nub.

◆ set_mu()

template<class OCP , int N, int kmax>
template<typename VectorType >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_mu ( const MatrixBase< VectorType > &  mu)
inline

Sets the Lagrange multiplier with respect to the control input bounds constraint.

Parameters
[in]muThe Lagrange multiplier. Size must be SingleShootingCGMRESSolver::nub.

◆ set_mu_array()

template<class OCP , int N, int kmax>
template<typename T >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_mu_array ( const T &  mu_array)
inline

Sets the array of the Lagrange multiplier with respect to the control input bounds constraint.

Parameters
[in]mu_arrayA container (std::vector, std::array, etc.) of the Lagrange multiplier. Size must be N and size of each element must be SingleShootingCGMRESSolver::nub.

◆ set_u()

template<class OCP , int N, int kmax>
template<typename VectorType >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_u ( const MatrixBase< VectorType > &  u)
inline

Sets the control input vector.

Parameters
[in]uThe control input vector. Size must be SingleShootingCGMRESSolver::nu.

◆ set_u_array()

template<class OCP , int N, int kmax>
template<typename T >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_u_array ( const T &  u_array)
inline

Sets the array of control input vector.

Parameters
[in]u_arrayA container (std::vector, std::array, etc.) of the control input vector. Size must be N and size of each element must be SingleShootingCGMRESSolver::nu.

◆ set_uc()

template<class OCP , int N, int kmax>
template<typename VectorType >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_uc ( const MatrixBase< VectorType > &  uc)
inline

Sets the control input vector and Lagrange multiplier with respect to the equality constraints.

Parameters
[in]ucConcatenatin of the control input vector and Lagrange multiplier with respect to the equality constraints. Size must be SingleShootingCGMRESSolver::nuc.

◆ set_uc_array()

template<class OCP , int N, int kmax>
template<typename T >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::set_uc_array ( const T &  uc_array)
inline

Sets the array of control input vector and Lagrange multiplier with respect to the equality constraints.

Parameters
[in]uc_arrayA container (std::vector, std::array, etc.) of the concatenatin of the control input vector and Lagrange multiplier with respect to the equality constraints. Size must be N and size of each element must be SingleShootingCGMRESSolver::nuc.

◆ ucopt()

template<class OCP , int N, int kmax>
const std::array< Vector< nuc >, N > & cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::ucopt ( ) const
inline

Getter of the optimal solution.

Returns
const reference to the optimal concatenatins of the control input vector and Lagrange multiplier with respect to the equality constraints.

◆ uopt()

template<class OCP , int N, int kmax>
const std::array< Vector< nu >, N > & cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::uopt ( ) const
inline

Getter of the optimal solution.

Returns
const reference to the optimal control input vectors over the horizon.

◆ update()

template<class OCP , int N, int kmax>
template<typename VectorType >
void cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::update ( const Scalar  t,
const MatrixBase< VectorType > &  x 
)
inline

Updates the solution by performing C/GMRES method.

Parameters
[in]tInitial time of the horizon.
[in]xInitial state of the horizon. Size must be SingleShootingCGMRESSolver::nx.

◆ xopt()

template<class OCP , int N, int kmax>
const std::array< Vector< nx >, N+1 > & cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::xopt ( ) const
inline

Getter of the optimal solution.

Returns
const reference to the optimal state vectors over the horizon.

Friends And Related Function Documentation

◆ operator<<

template<class OCP , int N, int kmax>
std::ostream & operator<< ( std::ostream &  os,
const SingleShootingCGMRESSolver< OCP, N, kmax > &  solver 
)
friend

Member Data Documentation

◆ dim

template<class OCP , int N, int kmax>
constexpr int cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::dim = nuc * N + 2 * N * nub
staticconstexpr

Dimension of the linear problem solved by the GMRES solver.

◆ nc

template<class OCP , int N, int kmax>
constexpr int cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::nc = OCP::nc
staticconstexpr

Dimension of the equality constraints.

◆ nu

template<class OCP , int N, int kmax>
constexpr int cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::nu = OCP::nu
staticconstexpr

Dimension of the control input.

◆ nub

template<class OCP , int N, int kmax>
constexpr int cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::nub = OCP::nub
staticconstexpr

Dimension of the bound constraints on the control input.

◆ nuc

template<class OCP , int N, int kmax>
constexpr int cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::nuc = nu + nc
staticconstexpr

Dimension of the concatenation of the control input and equality constraints.

◆ nx

template<class OCP , int N, int kmax>
constexpr int cgmres::SingleShootingCGMRESSolver< OCP, N, kmax >::nx = OCP::nx
staticconstexpr

Dimension of the state.


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