File linearization.h#

namespace sym

Typedefs

template<typename Scalar>
using SparseLinearization = Linearization<Eigen::SparseMatrix<Scalar>>#
using SparseLinearizationd = SparseLinearization<double>#
using SparseLinearizationf = SparseLinearization<float>#
template<typename Scalar>
using DenseLinearization = Linearization<MatrixX<Scalar>>#
using DenseLinearizationd = DenseLinearization<double>#
using DenseLinearizationf = DenseLinearization<float>#

Functions

template<typename Scalar>
sparse_matrix_structure_t GetSparseStructure(const Eigen::SparseMatrix<Scalar> &matrix)#

Returns the sparse matrix structure of matrix.

template<typename Derived>
sparse_matrix_structure_t GetSparseStructure(const Eigen::EigenBase<Derived> &matrix)#

Return a default initialized sparse structure because arg is dense.

template<typename MatrixType, typename Scalar, typename std::enable_if_t<kIsSparseEigenType<MatrixType>, bool> = true>
Eigen::Map<const Eigen::SparseMatrix<Scalar>> MatrixViewFromSparseStructure(const sparse_matrix_structure_t &structure, const MatrixX<Scalar> &values)#

Returns the reconstructed matrix from the stored values and sparse_matrix_structure_t

Useful for reconstructing the jacobian for an iteration in the OptimizationStats

The result contains references to structure and values, so its lifetime must be shorter than those.

Parameters:
  • rows – The number of rows in the matrix (must match structure.shape for sparse matrices)

  • cols – The number of columns in the matrix (must match structure.shape for sparse matrices)

  • structure – The sparsity pattern of the matrix, as obtained from GetSparseStructure or stored in OptimizationStats::iterations. For dense matrices, this should be empty.

  • values – The coefficients of the matrix; flat for sparse matrices, 2d for dense

Template Parameters:
  • MatrixType – The type of the matrix used to decide if the result is sparse or dense. This is not otherwise used

  • Scalar – The scalar type of the result (must match the scalar type of values)

template<typename MatrixType, typename Scalar, typename std::enable_if_t<kIsEigenType<MatrixType>, bool> = true>
Eigen::Map<const MatrixX<Scalar>> MatrixViewFromSparseStructure(const sparse_matrix_structure_t &structure, const MatrixX<Scalar> &values)#

Returns the reconstructed matrix from the stored values and sparse_matrix_structure_t

Useful for reconstructing the jacobian for an iteration in the OptimizationStats

The result contains references to structure and values, so its lifetime must be shorter than those.

Parameters:
  • rows – The number of rows in the matrix (must match structure.shape for sparse matrices)

  • cols – The number of columns in the matrix (must match structure.shape for sparse matrices)

  • structure – The sparsity pattern of the matrix, as obtained from GetSparseStructure or stored in OptimizationStats::iterations. For dense matrices, this should be empty.

  • values – The coefficients of the matrix; flat for sparse matrices, 2d for dense

Template Parameters:
  • MatrixType – The type of the matrix used to decide if the result is sparse or dense. This is not otherwise used

  • Scalar – The scalar type of the result (must match the scalar type of values)

template<typename Scalar>
Eigen::Map<const VectorX<Scalar>> JacobianValues(const Eigen::SparseMatrix<Scalar> &matrix)#

Returns coefficients of matrix. Overloads exist for both dense and sparse matrices to make writing generic code easier.

This version returns the non-zero values of an Eigen::SparseMatrix

Note: it returns a map, so be careful about mutating or disposing of matrix before you are finished with the output.

template<typename Scalar>
const MatrixX<Scalar> &JacobianValues(const MatrixX<Scalar> &matrix)#

Returns coefficients of matrix. Overloads exist for both dense and sparse matrices to make writing generic code easier.

Returns a const-ref to the argument.

template<typename MatrixType>
struct Linearization#
#include <linearization.h>

Class for storing a problem linearization evaluated at a Values (i.e. a residual, jacobian, hessian, and rhs)

MatrixType is expected to be an Eigen MatrixX or SparseMatrix.

Public Types

using Scalar = typename MatrixType::Scalar#
using Vector = VectorX<Scalar>#
using Matrix = MatrixType#

Public Functions

inline void Reset()#

Set to invalid

inline bool IsInitialized() const#

Returns whether the linearization is currently valid for the corresponding values. Accessing any of the members when this is false could result in unexpected behavior

inline void SetInitialized(const bool initialized = true)#
inline Scalar Error() const#
inline Scalar LinearDeltaError(const Vector &x_update, const Vector &damping_vector) const#

Returns the change in error predicted by the Linearization at the given update

Parameters:
  • x_update – The update to the values

  • damping_vector – The vector added to the diagonal of the hessian during the linear solve

Public Members

Vector residual#
Matrix hessian_lower#
Matrix jacobian#
Vector rhs#

Private Members

bool initialized_ = {false}#