WPILibC++ 2027.0.0-alpha-4
Loading...
Searching...
No Matches
sqp_matrix_callbacks.hpp
Go to the documentation of this file.
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <functional>
6
7#include <Eigen/Core>
8#include <Eigen/SparseCore>
9
10namespace slp {
11
12/// Matrix callbacks for the Sequential Quadratic Programming (SQP) solver.
13///
14/// @tparam Scalar Scalar type.
15template <typename Scalar>
17 /// Type alias for dense vector.
18 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
19 /// Type alias for sparse matrix.
20 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
21 /// Type alias for sparse vector.
22 using SparseVector = Eigen::SparseVector<Scalar>;
23
24 /// Cost function value f(x) getter.
25 ///
26 /// <table>
27 /// <tr>
28 /// <th>Variable</th>
29 /// <th>Rows</th>
30 /// <th>Columns</th>
31 /// </tr>
32 /// <tr>
33 /// <td>x</td>
34 /// <td>num_decision_variables</td>
35 /// <td>1</td>
36 /// </tr>
37 /// <tr>
38 /// <td>f(x)</td>
39 /// <td>1</td>
40 /// <td>1</td>
41 /// </tr>
42 /// </table>
43 std::function<Scalar(const DenseVector& x)> f;
44
45 /// Cost function gradient ∇f(x) getter.
46 ///
47 /// <table>
48 /// <tr>
49 /// <th>Variable</th>
50 /// <th>Rows</th>
51 /// <th>Columns</th>
52 /// </tr>
53 /// <tr>
54 /// <td>x</td>
55 /// <td>num_decision_variables</td>
56 /// <td>1</td>
57 /// </tr>
58 /// <tr>
59 /// <td>∇f(x)</td>
60 /// <td>num_decision_variables</td>
61 /// <td>1</td>
62 /// </tr>
63 /// </table>
64 std::function<SparseVector(const DenseVector& x)> g;
65
66 /// Lagrangian Hessian ∇ₓₓ²L(x, y) getter.
67 ///
68 /// L(xₖ, yₖ) = f(xₖ) − yₖᵀcₑ(xₖ)
69 ///
70 /// <table>
71 /// <tr>
72 /// <th>Variable</th>
73 /// <th>Rows</th>
74 /// <th>Columns</th>
75 /// </tr>
76 /// <tr>
77 /// <td>x</td>
78 /// <td>num_decision_variables</td>
79 /// <td>1</td>
80 /// </tr>
81 /// <tr>
82 /// <td>y</td>
83 /// <td>num_equality_constraints</td>
84 /// <td>1</td>
85 /// </tr>
86 /// <tr>
87 /// <td>∇ₓₓ²L(x, y)</td>
88 /// <td>num_decision_variables</td>
89 /// <td>num_decision_variables</td>
90 /// </tr>
91 /// </table>
92 std::function<SparseMatrix(const DenseVector& x, const DenseVector& y)> H;
93
94 /// Equality constraint value cₑ(x) getter.
95 ///
96 /// <table>
97 /// <tr>
98 /// <th>Variable</th>
99 /// <th>Rows</th>
100 /// <th>Columns</th>
101 /// </tr>
102 /// <tr>
103 /// <td>x</td>
104 /// <td>num_decision_variables</td>
105 /// <td>1</td>
106 /// </tr>
107 /// <tr>
108 /// <td>cₑ(x)</td>
109 /// <td>num_equality_constraints</td>
110 /// <td>1</td>
111 /// </tr>
112 /// </table>
113 std::function<DenseVector(const DenseVector& x)> c_e;
114
115 /// Equality constraint Jacobian ∂cₑ/∂x getter.
116 ///
117 /// ```
118 /// [∇ᵀcₑ₁(xₖ)]
119 /// Aₑ(x) = [∇ᵀcₑ₂(xₖ)]
120 /// [ ⋮ ]
121 /// [∇ᵀcₑₘ(xₖ)]
122 /// ```
123 ///
124 /// <table>
125 /// <tr>
126 /// <th>Variable</th>
127 /// <th>Rows</th>
128 /// <th>Columns</th>
129 /// </tr>
130 /// <tr>
131 /// <td>x</td>
132 /// <td>num_decision_variables</td>
133 /// <td>1</td>
134 /// </tr>
135 /// <tr>
136 /// <td>Aₑ(x)</td>
137 /// <td>num_equality_constraints</td>
138 /// <td>num_decision_variables</td>
139 /// </tr>
140 /// </table>
141 std::function<SparseMatrix(const DenseVector& x)> A_e;
142};
143
144} // namespace slp
Definition expression_graph.hpp:11
Matrix callbacks for the Sequential Quadratic Programming (SQP) solver.
Definition sqp_matrix_callbacks.hpp:16
std::function< SparseVector(const DenseVector &x)> g
Cost function gradient ∇f(x) getter.
Definition sqp_matrix_callbacks.hpp:64
std::function< Scalar(const DenseVector &x)> f
Cost function value f(x) getter.
Definition sqp_matrix_callbacks.hpp:43
Eigen::SparseVector< Scalar > SparseVector
Type alias for sparse vector.
Definition sqp_matrix_callbacks.hpp:22
std::function< DenseVector(const DenseVector &x)> c_e
Equality constraint value cₑ(x) getter.
Definition sqp_matrix_callbacks.hpp:113
Eigen::SparseMatrix< Scalar > SparseMatrix
Type alias for sparse matrix.
Definition sqp_matrix_callbacks.hpp:20
Eigen::Vector< Scalar, Eigen::Dynamic > DenseVector
Type alias for dense vector.
Definition sqp_matrix_callbacks.hpp:18
std::function< SparseMatrix(const DenseVector &x)> A_e
Equality constraint Jacobian ∂cₑ/∂x getter.
Definition sqp_matrix_callbacks.hpp:141
std::function< SparseMatrix(const DenseVector &x, const DenseVector &y)> H
Lagrangian Hessian ∇ₓₓ²L(x, y) getter.
Definition sqp_matrix_callbacks.hpp:92