proxsuite 0.6.7
The Advanced Proximal Optimization Toolbox
Loading...
Searching...
No Matches
random_qp_problems.hpp
Go to the documentation of this file.
1#ifndef PROXSUITE_PROXQP_UTILS_RANDOM_QP_PROBLEMS_HPP
2#define PROXSUITE_PROXQP_UTILS_RANDOM_QP_PROBLEMS_HPP
3
4#include <Eigen/Core>
5#include <Eigen/SparseCore>
6#include <Eigen/Cholesky>
7#include <Eigen/Eigenvalues>
8#include <Eigen/QR>
9#include <utility>
13#include <map>
14#include <random>
15
16namespace proxsuite {
17namespace proxqp {
18namespace utils {
19
20using c_int = long long;
21using c_float = double;
22
23namespace proxqp = proxsuite::proxqp;
24
25template<typename T, proxqp::Layout L>
26using Mat =
27 Eigen::Matrix<T,
28 Eigen::Dynamic,
29 Eigen::Dynamic,
30 (L == proxqp::colmajor) ? Eigen::ColMajor : Eigen::RowMajor>;
31template<typename T>
32using Vec = Eigen::Matrix<T, Eigen::Dynamic, 1>;
33
34template<typename Scalar>
35using SparseMat = Eigen::SparseMatrix<Scalar, Eigen::ColMajor, c_int>;
36
37using namespace proxqp;
38namespace eigen {
39template<typename T>
40void
42 Eigen::LLT<T>& out,
43 T const& mat)
44{
45 out.compute(mat);
46}
47template<typename T>
48void
50 Eigen::LDLT<T>& out,
51 T const& mat)
52{
53 out.compute(mat);
54}
59
64} // namespace eigen
65namespace rand {
66
67using proxqp::u32;
68using proxqp::u64;
69
70#ifdef _MSC_VER
71/* Using the MSCV compiler on Windows causes problems because the type uint128
72is not available. Therefore, we use a random number generator from the stdlib
73instead of our custom Lehmer random number generator. The necessary lehmer
74functions used in in our code are remplaced with calls to the stdlib.*/
75std::mt19937 gen(1234);
76std::uniform_real_distribution<> uniform_dist(0.0, 1.0);
77std::normal_distribution<double> normal_dist;
78using u128 = u64;
79inline auto
80uniform_rand() -> double
81{
82 double output = double(uniform_dist(gen));
83 return output;
84}
85inline auto
87{
88 static u64 output = gen();
89 return output;
90}
91
92inline void
93set_seed(u64 seed)
94{
95 gen.seed(seed);
96}
97
98inline auto
99normal_rand() -> double
100{
101 return normal_dist(gen);
102}
103#else
104using u128 = __uint128_t;
105
106constexpr u128 lehmer64_constant(0xda942042e4dd58b5);
107inline auto
109{
110 static u128 g_lehmer64_state = lehmer64_constant * lehmer64_constant;
111 return g_lehmer64_state;
112}
113
114inline auto
115lehmer64() -> u64
116{ // [0, 2^64)
118 return u64(lehmer_global() >> u128(64U));
119}
120
121inline void
122set_seed(u64 seed)
123{
124 lehmer_global() = u128(seed) + 1;
125 lehmer64();
126 lehmer64();
127}
128
129inline auto
130uniform_rand() -> double
131{ // [0, 2^53]
132 u64 a = lehmer64() / (1U << 11U);
133 return double(a) / double(u64(1) << 53U);
134}
135inline auto
136normal_rand() -> double
137{
138 static const double pi2 = std::atan(static_cast<double>(1)) * 8;
139
140 double u1 = uniform_rand();
141 double u2 = uniform_rand();
142
143 double ln = std::log(u1);
144 double sqrt = std::sqrt(-2 * ln);
145
146 return sqrt * std::cos(pi2 * u2);
147}
148#endif
149
150template<typename Scalar>
151auto
153{
154 auto v = Vec<Scalar>(nrows);
155
156 for (isize i = 0; i < nrows; ++i) {
157 v(i) = Scalar(rand::normal_rand());
158 }
159
160 return v;
161}
162template<typename Scalar>
163auto
164matrix_rand(isize nrows, isize ncols) -> Mat<Scalar, colmajor>
165{
166 auto v = Mat<Scalar, colmajor>(nrows, ncols);
167
168 for (isize i = 0; i < nrows; ++i) {
169 for (isize j = 0; j < ncols; ++j) {
170 v(i, j) = Scalar(rand::normal_rand());
171 }
172 }
173
174 return v;
175}
176
177namespace detail {
178template<typename Scalar>
179auto
181{
182 auto mat = rand::matrix_rand<Scalar>(n, n);
183 auto qr = mat.householderQr();
184 Mat<Scalar, colmajor> q = qr.householderQ();
185 return q;
186}
187using Input = std::pair<u128, isize>;
188} // namespace detail
189
190template<typename Scalar>
191auto
193{
194
195 static auto cache = std::map<detail::Input, Mat<Scalar, colmajor>>{};
196 auto input = detail::Input{ lehmer_global(), n };
197 auto it = cache.find(input);
198 if (it == cache.end()) {
199 auto res = cache.insert({
200 input,
201 detail::orthonormal_rand_impl<Scalar>(n),
202 });
203 it = res.first;
204 }
205 return (*it).second;
206}
207
208template<typename Scalar>
209auto
211{
212 auto const& q = rand::orthonormal_rand<Scalar>(n);
213 auto d = Vec<Scalar>(n);
214
215 {
216 using std::exp;
217 using std::log;
218 Scalar diff = log(cond);
219 for (isize i = 0; i < n; ++i) {
220 d(i) = exp(Scalar(i) / Scalar(n) * diff);
221 }
222 }
223
224 return q * d.asDiagonal() * q.transpose();
225}
226
227template<typename Scalar>
228auto
230 Scalar cond,
231 Scalar p) -> SparseMat<Scalar>
232{
233 auto H = SparseMat<Scalar>(n, n);
234
235 for (isize i = 0; i < n; ++i) {
236 auto random = Scalar(rand::normal_rand());
237 H.insert(i, i) = random;
238 }
239
240 for (isize i = 0; i < n; ++i) {
241 for (isize j = i + 1; j < n; ++j) {
242 auto urandom = rand::uniform_rand();
243 if (urandom < p / 2) {
244 auto random = Scalar(rand::normal_rand());
245 H.insert(i, j) = random;
246 }
247 }
248 }
249
250 Mat<Scalar, colmajor> H_dense = H.toDense();
251 Vec<Scalar> eigh =
252 H_dense.template selfadjointView<Eigen::Upper>().eigenvalues();
253
254 Scalar min = eigh.minCoeff();
255 Scalar max = eigh.maxCoeff();
256
257 // new_min = min + rho
258 // new_max = max + rho
259 //
260 // (max + rho)/(min + rho) = cond
261 // 1 + (max - min) / (min + rho) = cond
262 // (max - min) / (min + rho) = cond - 1
263 // min + rho = (max - min) / (cond - 1)
264 // rho = (max - min)/(cond - 1) - min
265 Scalar rho = (max - min) / (cond - 1) - min;
266 if (max == min) {
267 rho += 1;
268 }
269
270 for (isize i = 0; i < n; ++i) {
271 H.coeffRef(i, i) += rho;
272 }
273
274 H.makeCompressed();
275 return H;
276}
277
278template<typename Scalar>
279auto
281 Scalar rho,
282 Scalar p) -> SparseMat<Scalar>
283{
284 auto H = SparseMat<Scalar>(n, n);
285
286 H.setZero();
287
288 for (isize i = 0; i < n; ++i) {
289 for (isize j = i + 1; j < n; ++j) {
290 auto urandom = rand::uniform_rand();
291 if (urandom < p) {
292 auto random = Scalar(rand::normal_rand());
293 H.insert(i, j) = random;
294 }
295 }
296 }
297 Mat<Scalar, colmajor> H_dense = H.toDense();
298 Vec<Scalar> eigh =
299 H_dense.template selfadjointView<Eigen::Upper>().eigenvalues();
300 Scalar min = eigh.minCoeff();
301 for (isize i = 0; i < n; ++i) {
302 H.coeffRef(i, i) += (rho + abs(min));
303 }
304
305 H.makeCompressed();
306 return H;
307}
308
309template<typename Scalar>
310auto
312 Scalar rho,
313 Scalar p) -> Mat<Scalar, colmajor>
314{
315 auto H = Mat<Scalar, colmajor>(n, n);
316 H.setZero();
317
318 for (isize i = 0; i < n; ++i) {
319 for (isize j = 0; j < n; ++j) {
320 auto urandom = rand::uniform_rand();
321 if (urandom < p / 2) {
322 auto random = Scalar(rand::normal_rand());
323 H(i, j) = random;
324 }
325 }
326 }
327
328 H = ((H + H.transpose()) * 0.5)
329 .eval(); // safe no aliasing :
330 // https://eigen.tuxfamily.org/dox/group__TopicAliasing.html
331 // H.array() /= 2.;
332 Vec<Scalar> eigh = H.template selfadjointView<Eigen::Upper>().eigenvalues();
333 Scalar min = eigh.minCoeff();
334 H.diagonal().array() += (rho + abs(min));
335
336 return H;
337}
338
339template<typename Scalar>
340auto
341sparse_matrix_rand(isize nrows, isize ncols, Scalar p) -> SparseMat<Scalar>
342{
343 auto A = SparseMat<Scalar>(nrows, ncols);
344
345 for (isize i = 0; i < nrows; ++i) {
346 for (isize j = 0; j < ncols; ++j) {
347 if (rand::uniform_rand() < p) {
348 A.insert(i, j) = Scalar(rand::normal_rand());
349 }
350 }
351 }
352 A.makeCompressed();
353 return A;
354}
355
356template<typename Scalar>
357auto
359 isize ncols,
360 Scalar p) -> Mat<Scalar, colmajor>
361{
362 auto A = Mat<Scalar, colmajor>(nrows, ncols);
363 A.setZero();
364 for (isize i = 0; i < nrows; ++i) {
365 for (isize j = 0; j < ncols; ++j) {
366 if (rand::uniform_rand() < p) {
367 A(i, j) = Scalar(rand::normal_rand());
368 }
369 }
370 }
371 return A;
372}
373
374} // namespace rand
375using proxqp::usize;
376
377namespace osqp {
378auto
380auto
382} // namespace osqp
383
384template<typename T>
385auto
387 Mat<T, proxqp::colmajor> const& lhs,
389{
390 return lhs.operator*(rhs);
391}
392template<typename To, typename From>
393auto
395{
396 return from.template cast<To>();
397}
398LDLT_EXPLICIT_TPL_DECL(2, matmul_impl<long double>);
399LDLT_EXPLICIT_TPL_DECL(1, mat_cast<proxqp::f64, long double>);
400LDLT_EXPLICIT_TPL_DECL(1, mat_cast<proxqp::f32, long double>);
401
402template<typename MatLhs, typename MatRhs, typename T = typename MatLhs::Scalar>
403auto
404matmul(MatLhs const& a, MatRhs const& b) -> Mat<T, proxqp::colmajor>
405{
406 using Upscaled = typename std::
407 conditional<std::is_floating_point<T>::value, long double, T>::type;
408
409 return mat_cast<T, Upscaled>(matmul_impl<Upscaled>(
410 Mat<T, proxqp::colmajor>(a).template cast<Upscaled>(),
411 Mat<T, proxqp::colmajor>(b).template cast<Upscaled>()));
412}
413
414template<typename MatLhs,
415 typename MatMid,
416 typename MatRhs,
417 typename T = typename MatLhs::Scalar>
418auto
419matmul3(MatLhs const& a,
420 MatMid const& b,
421 MatRhs const& c) -> Mat<T, proxqp::colmajor>
422{
423 return matmul(matmul(a, b), c);
424}
425
426VEG_TAG(from_data, FromData);
427
429{
431 EigenNoAlloc(EigenNoAlloc const&) = delete;
432 auto operator=(EigenNoAlloc&&) -> EigenNoAlloc& = delete;
433 auto operator=(EigenNoAlloc const&) -> EigenNoAlloc& = delete;
434
435#if defined(EIGEN_RUNTIME_NO_MALLOC)
436 EigenNoAlloc() noexcept { Eigen::internal::set_is_malloc_allowed(false); }
437 ~EigenNoAlloc() noexcept { Eigen::internal::set_is_malloc_allowed(true); }
438#else
439 EigenNoAlloc() = default;
440#endif
441};
442
443template<typename Scalar>
445dense_unconstrained_qp(proxqp::isize dim,
446 Scalar sparsity_factor,
447 Scalar strong_convexity_factor = Scalar(1e-2))
448{
449
451 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
452 dim, strong_convexity_factor, sparsity_factor);
453 Vec<Scalar> g = rand::vector_rand<Scalar>(dim);
455 rand::sparse_matrix_rand_not_compressed<Scalar>(0, dim, sparsity_factor);
456 Vec<Scalar> b = rand::vector_rand<Scalar>(0);
458 rand::sparse_matrix_rand_not_compressed<Scalar>(0, dim, sparsity_factor);
459 Vec<Scalar> u = rand::vector_rand<Scalar>(0);
460 Vec<Scalar> l = rand::vector_rand<Scalar>(0);
462 model.H = H;
463 model.g = g;
464 return model;
465}
466
467template<typename Scalar>
469dense_strongly_convex_qp(proxqp::isize dim,
470 proxqp::isize n_eq,
471 proxqp::isize n_in,
472 Scalar sparsity_factor,
473 Scalar strong_convexity_factor = Scalar(1e-2))
474{
475
477 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
478 dim, strong_convexity_factor, sparsity_factor);
479 Vec<Scalar> g = rand::vector_rand<Scalar>(dim);
481 rand::sparse_matrix_rand_not_compressed<Scalar>(n_eq, dim, sparsity_factor);
483 rand::sparse_matrix_rand_not_compressed<Scalar>(n_in, dim, sparsity_factor);
484
485 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
486 auto delta = Vec<Scalar>(n_in);
487
488 for (proxqp::isize i = 0; i < n_in; ++i) {
489 delta(i) = rand::uniform_rand();
490 }
491
492 Vec<Scalar> u = C * x_sol + delta;
493 Vec<Scalar> b = A * x_sol;
494 Vec<Scalar> l(n_in);
495 l.setZero();
496 l.array() -= 1.e20;
497
498 proxsuite::proxqp::dense::Model<Scalar> model(dim, n_eq, n_in);
499 model.H = H;
500 model.g = g;
501 model.A = A;
502 model.b = b;
503 model.C = C;
504 model.u = u;
505 model.l = l;
506 return model;
507}
508
509template<typename Scalar>
512 proxqp::isize n_eq,
513 proxqp::isize n_in,
514 Scalar sparsity_factor)
515{
516
518 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
519 dim, Scalar(0), sparsity_factor);
521 rand::sparse_matrix_rand_not_compressed<Scalar>(n_eq, dim, sparsity_factor);
523 rand::sparse_matrix_rand_not_compressed<Scalar>(n_in, dim, sparsity_factor);
524
525 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
526 Vec<Scalar> y_sol = rand::vector_rand<Scalar>(n_eq);
527 Vec<Scalar> z_sol = rand::vector_rand<Scalar>(n_in);
528 auto delta = Vec<Scalar>(n_in);
529
530 for (proxqp::isize i = 0; i < n_in; ++i) {
531 delta(i) = rand::uniform_rand();
532 }
533 auto Cx = C * x_sol;
534 Vec<Scalar> u = Cx + delta;
535 Vec<Scalar> b = A * x_sol;
536 Vec<Scalar> l = Cx - delta;
537 Vec<Scalar> g = -(H * x_sol + C.transpose() * z_sol + A.transpose() * y_sol);
538
539 proxsuite::proxqp::dense::Model<Scalar> model(dim, n_eq, n_in);
540 model.H = H;
541 model.g = g;
542 model.A = A;
543 model.b = b;
544 model.C = C;
545 model.u = u;
546 model.l = l;
547 return model;
548}
549
550template<typename Scalar>
552dense_degenerate_qp(proxqp::isize dim,
553 proxqp::isize n_eq,
554 proxqp::isize n_in,
555 Scalar sparsity_factor,
556 Scalar strong_convexity_factor = Scalar(1e-2))
557{
558
560 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
561 dim, strong_convexity_factor, sparsity_factor);
562 Vec<Scalar> g = rand::vector_rand<Scalar>(dim);
564 rand::sparse_matrix_rand_not_compressed<Scalar>(n_eq, dim, sparsity_factor);
566
567 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
568 auto delta = Vec<Scalar>(2 * n_in);
569
570 for (proxqp::isize i = 0; i < 2 * n_in; ++i) {
571 delta(i) = rand::uniform_rand();
572 }
573 Vec<Scalar> b = A * x_sol;
574
575 auto C_ =
576 rand::sparse_matrix_rand_not_compressed<Scalar>(n_in, dim, sparsity_factor);
577 C.setZero();
578 C.block(0, 0, n_in, dim) = C_;
579 C.block(n_in, 0, n_in, dim) = C_;
580 Vec<Scalar> u = C * x_sol + delta;
581 Vec<Scalar> l(2 * n_in);
582 l.setZero();
583 l.array() -= 1.e20;
584
585 proxsuite::proxqp::dense::Model<Scalar> model(dim, n_eq, n_in);
586 model.H = H;
587 model.g = g;
588 model.A = A;
589 model.b = b;
590 model.C = C;
591 model.u = u;
592 model.l = l;
593 return model;
594}
595
596template<typename Scalar>
598dense_box_constrained_qp(proxqp::isize dim,
599 proxqp::isize n_eq,
600 proxqp::isize n_in,
601 Scalar sparsity_factor,
602 Scalar strong_convexity_factor = Scalar(1e-2))
603{
604
606 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
607 dim, strong_convexity_factor, sparsity_factor);
608 Vec<Scalar> g = rand::vector_rand<Scalar>(dim);
610 rand::sparse_matrix_rand_not_compressed<Scalar>(n_eq, dim, sparsity_factor);
612
613 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
614 auto delta = Vec<Scalar>(n_in);
615
616 for (proxqp::isize i = 0; i < n_in; ++i) {
617 delta(i) = rand::uniform_rand();
618 }
619 Vec<Scalar> b = A * x_sol;
620 C.setZero();
621 C.diagonal().array() += 1;
622 Vec<Scalar> u = x_sol + delta;
623 Vec<Scalar> l = x_sol - delta;
624 proxsuite::proxqp::dense::Model<Scalar> model(dim, n_eq, n_in);
625 model.H = H;
626 model.g = g;
627 model.A = A;
628 model.b = b;
629 model.C = C;
630 model.u = u;
631 model.l = l;
632 return model;
633}
634
635template<typename Scalar>
638 proxqp::isize n_eq,
639 proxqp::isize n_in,
640 Scalar sparsity_factor,
641 Scalar strong_convexity_factor = Scalar(1e-2))
642{
643
644 SparseMat<Scalar> H = rand::sparse_positive_definite_rand_compressed<Scalar>(
645 dim, strong_convexity_factor, sparsity_factor);
646 Vec<Scalar> g = rand::vector_rand<Scalar>(dim);
648 rand::sparse_matrix_rand<Scalar>(n_eq, dim, sparsity_factor);
650 rand::sparse_matrix_rand<Scalar>(n_in, dim, sparsity_factor);
651
652 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
653 auto delta = Vec<Scalar>(n_in);
654
655 for (proxqp::isize i = 0; i < n_in; ++i) {
656 delta(i) = rand::uniform_rand();
657 }
658
659 Vec<Scalar> u = C * x_sol + delta;
660 Vec<Scalar> b = A * x_sol;
661 Vec<Scalar> l(n_in);
662 l.setZero();
663 l.array() -= 1.e20;
664
665 proxsuite::proxqp::sparse::SparseModel<Scalar> res{ H, g, A, b, C, u, l };
666 return res;
667}
668
669} // namespace utils
670} // namespace proxqp
671} // namespace proxsuite
672
673#endif /* end of include guard PROXSUITE_PROXQP_UTILS_RANDOM_QP_PROBLEMS_HPP \
674 */
#define LDLT_EXPLICIT_TPL_DECL(NParams,...)
Definition views.hpp:58
#define VEG_TAG(Name, Type)
Definition macros.hpp:550
void ldlt_compute(Eigen::LDLT< T > &out, T const &mat)
void llt_compute(Eigen::LLT< T > &out, T const &mat)
auto to_sparse_sym(Mat< c_float, colmajor > const &mat) -> SparseMat< c_float >
auto to_sparse(Mat< c_float, colmajor > const &mat) -> SparseMat< c_float >
auto orthonormal_rand_impl(isize n) -> Mat< Scalar, colmajor >
auto orthonormal_rand(isize n) -> Mat< Scalar, colmajor > const &
auto sparse_matrix_rand_not_compressed(isize nrows, isize ncols, Scalar p) -> Mat< Scalar, colmajor >
auto sparse_positive_definite_rand_compressed(isize n, Scalar rho, Scalar p) -> SparseMat< Scalar >
auto matrix_rand(isize nrows, isize ncols) -> Mat< Scalar, colmajor >
auto sparse_positive_definite_rand(isize n, Scalar cond, Scalar p) -> SparseMat< Scalar >
auto positive_definite_rand(isize n, Scalar cond) -> Mat< Scalar, colmajor >
auto vector_rand(isize nrows) -> Vec< Scalar >
auto sparse_positive_definite_rand_not_compressed(isize n, Scalar rho, Scalar p) -> Mat< Scalar, colmajor >
constexpr u128 lehmer64_constant(0xda942042e4dd58b5)
auto sparse_matrix_rand(isize nrows, isize ncols, Scalar p) -> SparseMat< Scalar >
Eigen::Matrix< T, Eigen::Dynamic, 1 > Vec
proxsuite::proxqp::dense::Model< Scalar > dense_strongly_convex_qp(proxqp::isize dim, proxqp::isize n_eq, proxqp::isize n_in, Scalar sparsity_factor, Scalar strong_convexity_factor=Scalar(1e-2))
proxsuite::proxqp::dense::Model< Scalar > dense_box_constrained_qp(proxqp::isize dim, proxqp::isize n_eq, proxqp::isize n_in, Scalar sparsity_factor, Scalar strong_convexity_factor=Scalar(1e-2))
proxsuite::proxqp::dense::Model< Scalar > dense_degenerate_qp(proxqp::isize dim, proxqp::isize n_eq, proxqp::isize n_in, Scalar sparsity_factor, Scalar strong_convexity_factor=Scalar(1e-2))
Eigen::SparseMatrix< Scalar, Eigen::ColMajor, c_int > SparseMat
proxsuite::proxqp::dense::Model< Scalar > dense_unconstrained_qp(proxqp::isize dim, Scalar sparsity_factor, Scalar strong_convexity_factor=Scalar(1e-2))
auto matmul3(MatLhs const &a, MatMid const &b, MatRhs const &c) -> Mat< T, proxqp::colmajor >
auto mat_cast(Mat< From, proxqp::colmajor > const &from) -> Mat< To, proxqp::colmajor >
auto matmul_impl(Mat< T, proxqp::colmajor > const &lhs, Mat< T, proxqp::colmajor > const &rhs) -> Mat< T, proxqp::colmajor >
proxsuite::proxqp::sparse::SparseModel< Scalar > sparse_strongly_convex_qp(proxqp::isize dim, proxqp::isize n_eq, proxqp::isize n_in, Scalar sparsity_factor, Scalar strong_convexity_factor=Scalar(1e-2))
proxsuite::proxqp::dense::Model< Scalar > dense_not_strongly_convex_qp(proxqp::isize dim, proxqp::isize n_eq, proxqp::isize n_in, Scalar sparsity_factor)
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic,(L==proxqp::colmajor) ? Eigen::ColMajor :Eigen::RowMajor > Mat
auto matmul(MatLhs const &a, MatRhs const &b) -> Mat< T, proxqp::colmajor >
decltype(sizeof(0)) usize
Definition views.hpp:24
constexpr Layout colmajor
Definition views.hpp:244
This class stores the model of the QP problem.
Definition model.hpp:24
EigenNoAlloc(EigenNoAlloc const &)=delete
auto operator=(EigenNoAlloc const &) -> EigenNoAlloc &=delete
EigenNoAlloc(EigenNoAlloc &&)=delete
auto operator=(EigenNoAlloc &&) -> EigenNoAlloc &=delete