WPILibC++ 2027.0.0-alpha-5
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 /// Number of decision variables.
26
27 /// Number of equality constraints.
29
30 /// Cost function value f(x) getter.
31 ///
32 /// <table>
33 /// <tr>
34 /// <th>Variable</th>
35 /// <th>Rows</th>
36 /// <th>Columns</th>
37 /// </tr>
38 /// <tr>
39 /// <td>x</td>
40 /// <td>num_decision_variables</td>
41 /// <td>1</td>
42 /// </tr>
43 /// <tr>
44 /// <td>f(x)</td>
45 /// <td>1</td>
46 /// <td>1</td>
47 /// </tr>
48 /// </table>
49 std::function<Scalar(const DenseVector& x)> f;
50
51 /// Cost function gradient ∇f(x) getter.
52 ///
53 /// <table>
54 /// <tr>
55 /// <th>Variable</th>
56 /// <th>Rows</th>
57 /// <th>Columns</th>
58 /// </tr>
59 /// <tr>
60 /// <td>x</td>
61 /// <td>num_decision_variables</td>
62 /// <td>1</td>
63 /// </tr>
64 /// <tr>
65 /// <td>∇f(x)</td>
66 /// <td>num_decision_variables</td>
67 /// <td>1</td>
68 /// </tr>
69 /// </table>
70 std::function<SparseVector(const DenseVector& x)> g;
71
72 /// Lagrangian Hessian ∇ₓₓ²L(x, y) getter.
73 ///
74 /// L(x, y) = f(x) − yᵀcₑ(x)
75 ///
76 /// <table>
77 /// <tr>
78 /// <th>Variable</th>
79 /// <th>Rows</th>
80 /// <th>Columns</th>
81 /// </tr>
82 /// <tr>
83 /// <td>x</td>
84 /// <td>num_decision_variables</td>
85 /// <td>1</td>
86 /// </tr>
87 /// <tr>
88 /// <td>y</td>
89 /// <td>num_equality_constraints</td>
90 /// <td>1</td>
91 /// </tr>
92 /// <tr>
93 /// <td>∇ₓₓ²L(x, y)</td>
94 /// <td>num_decision_variables</td>
95 /// <td>num_decision_variables</td>
96 /// </tr>
97 /// </table>
98 std::function<SparseMatrix(const DenseVector& x, const DenseVector& y)> H;
99
100 /// Constraint part of Lagrangian Hessian ∇ₓₓ²(−yᵀcₑ(x)) getter.
101 ///
102 /// <table>
103 /// <tr>
104 /// <th>Variable</th>
105 /// <th>Rows</th>
106 /// <th>Columns</th>
107 /// </tr>
108 /// <tr>
109 /// <td>x</td>
110 /// <td>num_decision_variables</td>
111 /// <td>1</td>
112 /// </tr>
113 /// <tr>
114 /// <td>y</td>
115 /// <td>num_equality_constraints</td>
116 /// <td>1</td>
117 /// </tr>
118 /// <tr>
119 /// <td>∇ₓₓ²(−yᵀcₑ(x))</td>
120 /// <td>num_decision_variables</td>
121 /// <td>num_decision_variables</td>
122 /// </tr>
123 /// </table>
124 std::function<SparseMatrix(const DenseVector& x, const DenseVector& y)> H_c;
125
126 /// Equality constraint value cₑ(x) getter.
127 ///
128 /// <table>
129 /// <tr>
130 /// <th>Variable</th>
131 /// <th>Rows</th>
132 /// <th>Columns</th>
133 /// </tr>
134 /// <tr>
135 /// <td>x</td>
136 /// <td>num_decision_variables</td>
137 /// <td>1</td>
138 /// </tr>
139 /// <tr>
140 /// <td>cₑ(x)</td>
141 /// <td>num_equality_constraints</td>
142 /// <td>1</td>
143 /// </tr>
144 /// </table>
145 std::function<DenseVector(const DenseVector& x)> c_e;
146
147 /// Equality constraint Jacobian ∂cₑ/∂x getter.
148 ///
149 /// ```
150 /// [∇ᵀcₑ₁(x)]
151 /// Aₑ(x) = [∇ᵀcₑ₂(x)]
152 /// [ ⋮ ]
153 /// [∇ᵀcₑₘ(x)]
154 /// ```
155 ///
156 /// <table>
157 /// <tr>
158 /// <th>Variable</th>
159 /// <th>Rows</th>
160 /// <th>Columns</th>
161 /// </tr>
162 /// <tr>
163 /// <td>x</td>
164 /// <td>num_decision_variables</td>
165 /// <td>1</td>
166 /// </tr>
167 /// <tr>
168 /// <td>Aₑ(x)</td>
169 /// <td>num_equality_constraints</td>
170 /// <td>num_decision_variables</td>
171 /// </tr>
172 /// </table>
173 std::function<SparseMatrix(const DenseVector& x)> A_e;
174};
175
176} // namespace slp
Definition to_underlying.hpp:7
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:70
int num_decision_variables
Number of decision variables.
Definition sqp_matrix_callbacks.hpp:25
std::function< Scalar(const DenseVector &x)> f
Cost function value f(x) getter.
Definition sqp_matrix_callbacks.hpp:49
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:145
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:173
int num_equality_constraints
Number of equality constraints.
Definition sqp_matrix_callbacks.hpp:28
std::function< SparseMatrix(const DenseVector &x, const DenseVector &y)> H_c
Constraint part of Lagrangian Hessian ∇ₓₓ²(−yᵀcₑ(x)) getter.
Definition sqp_matrix_callbacks.hpp:124
std::function< SparseMatrix(const DenseVector &x, const DenseVector &y)> H
Lagrangian Hessian ∇ₓₓ²L(x, y) getter.
Definition sqp_matrix_callbacks.hpp:98