proxsuite 0.6.7
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
model.hpp
Go to the documentation of this file.
1//
2// Copyright (c) 2022 INRIA
3//
5#ifndef PROXSUITE_PROXQP_SPARSE_MODEL_HPP
6#define PROXSUITE_PROXQP_SPARSE_MODEL_HPP
7
8#include <Eigen/Sparse>
11
12namespace proxsuite {
13namespace proxqp {
14namespace sparse {
18
21template<typename T, typename I>
22struct Model
23{
24 typedef T Scalar;
25 using VectorType = Eigen::Matrix<T, Eigen::Dynamic, 1>;
26
27 isize dim;
28 isize n_eq;
29 isize n_in;
30
31 isize H_nnz;
32 isize A_nnz;
33 isize C_nnz;
34
38
42
47
54 Model(isize dim, isize n_eq, isize n_in)
55 : dim(dim)
56 , n_eq(n_eq)
57 , n_in(n_in)
58 {
60 std::invalid_argument,
61 "wrong argument size: the dimension wrt primal "
62 "variable x should be strictly positive.");
63
64 const T infinite_bound_value = helpers::infinite_bound<T>::value();
65
66 g.setZero();
67 b.setZero();
68 u.setZero();
69 u.fill(+infinite_bound_value); // in case it appears u is nullopt (i.e., the
70 // problem is only lower bounded)
71 l.fill(-infinite_bound_value); // in case it appears l is nullopt (i.e., the
72 // problem is only upper bounded)
73 }
77 auto kkt() const -> proxsuite::linalg::sparse::MatRef<T, I>
78 {
79 auto n_tot = kkt_col_ptrs.len() - 1;
80 auto nnz =
81 isize(proxsuite::linalg::sparse::util::zero_extend(kkt_col_ptrs[n_tot]));
82 return {
83 proxsuite::linalg::sparse::from_raw_parts,
84 n_tot,
85 n_tot,
86 nnz,
88 nullptr,
90 kkt_values.ptr(),
91 };
92 }
97 {
98 auto n_tot = kkt_col_ptrs.len() - 1;
99 auto nnz =
100 isize(proxsuite::linalg::sparse::util::zero_extend(kkt_col_ptrs[n_tot]));
101 return {
102 proxsuite::linalg::sparse::from_raw_parts,
103 n_tot,
104 n_tot,
105 nnz,
107 nullptr,
109 kkt_values.ptr_mut(),
110 };
111 }
115 auto kkt_unscaled() const -> proxsuite::linalg::sparse::MatRef<T, I>
116 {
117 auto n_tot = kkt_col_ptrs_unscaled.len() - 1;
118 auto nnz = isize(proxsuite::linalg::sparse::util::zero_extend(
119 kkt_col_ptrs_unscaled[n_tot]));
120 return {
121 proxsuite::linalg::sparse::from_raw_parts,
122 n_tot,
123 n_tot,
124 nnz,
126 nullptr,
129 };
130 }
135 {
136 auto n_tot = kkt_col_ptrs_unscaled.len() - 1;
137 auto nnz = isize(proxsuite::linalg::sparse::util::zero_extend(
138 kkt_col_ptrs_unscaled[n_tot]));
139 return {
140 proxsuite::linalg::sparse::from_raw_parts,
141 n_tot,
142 n_tot,
143 nnz,
145 nullptr,
147 kkt_values_unscaled.ptr_mut(),
148 };
149 }
150};
151
152template<typename _Scalar>
154{
155 typedef _Scalar Scalar;
156 enum
157 {
158 layout = Eigen::RowMajor
159 };
160
161 using VectorType = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
162
163 Eigen::SparseMatrix<Scalar, 1> H;
165 Eigen::SparseMatrix<Scalar, 1> A;
167 Eigen::SparseMatrix<Scalar, 1> C;
170
171 template<typename Vector_g,
172 typename Vector_b,
173 typename Vector_u,
174 typename Vector_l>
175 SparseModel(const Eigen::SparseMatrix<Scalar, 1>& H,
176 const Eigen::MatrixBase<Vector_g>& g,
177 const Eigen::SparseMatrix<Scalar, 1>& A,
178 const Eigen::MatrixBase<Vector_b>& b,
179 const Eigen::SparseMatrix<Scalar, 1>& C,
180 const Eigen::MatrixBase<Vector_u>& u,
181 const Eigen::MatrixBase<Vector_l>& l) noexcept
182 : H(H)
183 , g(g)
184 , A(A)
185 , b(b)
186 , C(C)
187 , u(u)
188 , l(l)
189 {
190 }
191
193 {
194 return {
195 { proxqp::from_eigen, H },
196 { proxqp::from_eigen, g },
197 { proxqp::from_eigen, A },
198 { proxqp::from_eigen, b },
199 { proxqp::from_ptr_rows_cols_stride,
200 nullptr,
201 0,
202 proxqp::isize(H.rows()),
203 0 },
204 { proxqp::from_ptr_size, nullptr, 0 },
205 };
206 }
208 {
209 return {
210 { proxqp::from_eigen, H },
211 { proxqp::from_eigen, g },
212 { proxqp::from_eigen, A },
213 { proxqp::from_eigen, b },
214 { proxqp::from_ptr_rows_cols_stride,
215 nullptr,
216 0,
217 proxqp::isize(H.rows()),
218 0 },
219 { proxqp::from_ptr_size, nullptr, 0 },
220 };
221 }
222};
223
224} // namespace sparse
225} // namespace proxqp
226} // namespace proxsuite
227
228#endif /* end of include guard PROXSUITE_PROXQP_SPARSE_MODEL_HPP */
#define PROXSUITE_THROW_PRETTY(condition, exception, message)
Definition macros.hpp:18
Eigen::Ref< Eigen::Matrix< T, DYN, DYN > const > MatRef
Definition fwd.hpp:36
VEG_NODISCARD VEG_INLINE auto ptr_mut() VEG_NOEXCEPT -> T *
Definition vec.hpp:966
VEG_NODISCARD VEG_INLINE auto ptr() const VEG_NOEXCEPT -> T const *
Definition vec.hpp:962
VEG_NODISCARD VEG_INLINE auto len() const VEG_NOEXCEPT -> isize
Definition vec.hpp:970
This class stores the model of the QP problem.
Definition model.hpp:23
auto kkt_unscaled() const -> proxsuite::linalg::sparse::MatRef< T, I >
Definition model.hpp:115
proxsuite::linalg::veg::Vec< I > kkt_col_ptrs
Definition model.hpp:35
proxsuite::linalg::veg::Vec< T > kkt_values
Definition model.hpp:37
Eigen::Matrix< T, Eigen::Dynamic, 1 > VectorType
Definition model.hpp:25
Model(isize dim, isize n_eq, isize n_in)
Definition model.hpp:54
proxsuite::linalg::veg::Vec< I > kkt_row_indices_unscaled
Definition model.hpp:40
proxsuite::linalg::veg::Vec< T > kkt_values_unscaled
Definition model.hpp:41
auto kkt_mut() -> proxsuite::linalg::sparse::MatMut< T, I >
Definition model.hpp:96
auto kkt_mut_unscaled() -> proxsuite::linalg::sparse::MatMut< T, I >
Definition model.hpp:134
auto kkt() const -> proxsuite::linalg::sparse::MatRef< T, I >
Definition model.hpp:77
proxsuite::linalg::veg::Vec< I > kkt_row_indices
Definition model.hpp:36
proxsuite::linalg::veg::Vec< I > kkt_col_ptrs_unscaled
Definition model.hpp:39
auto as_view() -> proxqp::dense::QpView< Scalar >
Definition model.hpp:192
SparseModel(const Eigen::SparseMatrix< Scalar, 1 > &H, const Eigen::MatrixBase< Vector_g > &g, const Eigen::SparseMatrix< Scalar, 1 > &A, const Eigen::MatrixBase< Vector_b > &b, const Eigen::SparseMatrix< Scalar, 1 > &C, const Eigen::MatrixBase< Vector_u > &u, const Eigen::MatrixBase< Vector_l > &l) noexcept
Definition model.hpp:175
Eigen::SparseMatrix< Scalar, 1 > C
Definition model.hpp:167
Eigen::SparseMatrix< Scalar, 1 > H
Definition model.hpp:163
Eigen::SparseMatrix< Scalar, 1 > A
Definition model.hpp:165
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > VectorType
Definition model.hpp:161
auto as_mut() -> proxqp::dense::QpViewMut< Scalar >
Definition model.hpp:207