C++ API Reference#
hyhound#
-
template<class T = real_t>
using hyhound::MatrixView = guanaqo::MatrixView<T, index_t># Non-owning view of a dense matrix (column-major storage order with unit inner stride).
- Template Parameters:
T – The matrix element type.
-
template<class T, Config<T> Conf = {}, class UpDown>
void hyhound::update_cholesky(MatrixView<T> L, MatrixView<T> A, UpDown signs, MatrixView<T> Ws = {{}})# Update the Cholesky factorization \( LL^\top = H \) such that the updated factorization satisfies \( \tilde L \tilde L^\top = H + A S A^\top \).
\( S \) is a diagonal matrix with the diagonal elements determined by the
signsparameter:If
signsis of type Update, then \( S = \mathrm{I} \);If
signsis of type Downdate, then \( S = -\mathrm{I} \);If
signsis of type UpDowndate, then \( S = \mathrm{diag}(s) \), where \( s = \mathrm{copysign}(1, \texttt{signs.signs}) \), and wheresigns.signscontains only the values+0.0and-0.0;If
signsis of type DownUpdate, then \( S = -\mathrm{diag}(s) \), where \( s = \mathrm{copysign}(1, \texttt{signs.signs}) \), and wheresigns.signscontains only the values+0.0and-0.0;If
signsis of type DiagonalUpDowndate, then \( S = \mathrm{diag}(\texttt{signs.diag}) \).
If \( L \) is a tall matrix, the function also returns a matrix \( \tilde A \) such that \( \tilde L \tilde L^\top + \begin{pmatrix} 0 \\ \tilde A \end{pmatrix} S \begin{pmatrix} 0 & \tilde A^\top \end{pmatrix} = L L^\top + A S A^\top \).
- Parameters:
L – [inout] (k × n, lower trapezoidal) On input, the lower-triangular Cholesky factor \( L \) of \( H \). On output, the lower-triangular Cholesky factor \( \tilde L \) of \( H + A S A^\top \).
A – [inout] (k × m) On input, the matrix \( A \). On output, the top n rows are overwritten, and the bottom k-n rows contain the matrix \( \tilde A \).
signs – [in] The update/downdate specification. See above for details.
Ws – [out] (r × n) If specified, this matrix will contain part of the block Householder representation that was applied to perform the update/downdate, and the top n rows of A will contain the corresponding Householder reflectors. The number of rows r depends on the block size specified by Config::block_size_r.
-
template<class T, Config<T> Conf = {}, class UpDown>
void hyhound::apply_householder(MatrixView<T> L, MatrixView<T> A, UpDown signs, std::type_identity_t<MatrixView<const T>> B, std::type_identity_t<MatrixView<const T>> Ws)# Apply a block Householder transformation \( \breve Q \), returning \( \begin{pmatrix} \tilde L & \tilde A \end{pmatrix} = \begin{pmatrix} L & A \end{pmatrix} \breve Q \).
- Parameters:
L – [inout] (l × n) On input, the matrix \( L \). On output, the matrix \( \tilde L \).
A – [inout] (l × m) On input, the matrix \( A \). On output, the matrix \( \tilde A \).
signs – [in] The update/downdate specification. See update_cholesky for details.
B – [in] (n × m) The Householder reflectors generated by update_cholesky.
Ws – [in] (r × n) The upper triangular block Householder representation generated by update_cholesky.
-
template<class T>
struct Config# Blocking and packing parameters for the update_cholesky and apply_householder functions.
Public Members
-
index_t block_size_r = detail::DefaultMicroKernelSizes<T>::DefaultSizeR#
Block size of the block column of L to process in the micro-kernels.
-
index_t block_size_s = detail::DefaultMicroKernelSizes<T>::DefaultSizeS#
Block size of the block row of L to process in the micro-kernels.
-
index_t num_blocks_r = 1#
Number of block columns per cache block.
-
index_t prefetch_dist_col_a = 4#
Column prefetch distance for the matrix A.
-
bool enable_packing = true#
Enable cache blocking by copying the current block row of A to a temporary buffer.
-
index_t block_size_r = detail::DefaultMicroKernelSizes<T>::DefaultSizeR#
-
struct Update#
Perform a factorization update, i.e. given the factorization of \( H \), compute the factorization of \( H + A A^\top \).
-
struct Downdate#
Perform a factorization downdate, i.e. given the factorization of \( H \), compute the factorization of \( H - A A^\top \).
-
template<class T>
struct UpDowndate# Perform a factorization update or downdate, depending on the given signs, i.e. given the factorization of \( H \), compute the factorization of \( H + A S A^\top \), where \( S \) is a diagonal matrix with \( S_{jj} = 1 \) if
signs[j] is+0.0, and \( S_{jj} = -1 \) ifsigns[j] is-0.0. Other values forsignsare not allowed.Public Functions
-
inline index_t size() const#
-
inline index_t size() const#
-
template<class T>
struct DownUpdate# Perform a factorization downdate or update, depending on the given signs, i.e. given the factorization of \( H \), compute the factorization of \( H - A S A^\top \), where \( S \) is a diagonal matrix with \( S_{jj} = 1 \) if
signs[j] is+0.0, and \( S_{jj} = -1 \) ifsigns[j] is-0.0. Other values forsignsare not allowed.Public Functions
-
inline index_t size() const#
-
inline index_t size() const#
OCP#
-
namespace ocp#
-
Functions
-
void factor(RiccatiFactor &factor, Eigen::Ref<const mat> Sigma)#
-
void update(RiccatiFactor &factor, Eigen::Ref<const mat> DeltaSigma)#
-
void solve(RiccatiFactor &factor)#
-
void factor(SchurFactor &factor, Eigen::Ref<const mat> Sigma)#
-
void update(SchurFactor &factor, Eigen::Ref<const mat> DeltaSigma)#
-
void solve(SchurFactor &factor)#
-
static auto with_signature_matrix(auto &&S)#
-
static void update_schur_rank_one(SchurFactor &factor, index_t j, index_t i, real_t Sigma_ji)#
-
struct OCPDataRiccati#
- #include <riccati.hpp>
-
struct RiccatiFactor#
- #include <riccati.hpp>
-
struct OCPDataSchur#
- #include <schur.hpp>
-
struct SchurFactor#
- #include <schur.hpp>
-
void factor(RiccatiFactor &factor, Eigen::Ref<const mat> Sigma)#