1#ifndef PROXSUITE_PROXQP_UTILS_RANDOM_QP_PROBLEMS_HPP
2#define PROXSUITE_PROXQP_UTILS_RANDOM_QP_PROBLEMS_HPP
5#include <Eigen/SparseCore>
6#include <Eigen/Cholesky>
7#include <Eigen/Eigenvalues>
25template<
typename T, proxqp::Layout L>
32using Vec = Eigen::Matrix<T, Eigen::Dynamic, 1>;
34template<
typename Scalar>
35using SparseMat = Eigen::SparseMatrix<Scalar, Eigen::ColMajor, c_int>;
75std::mt19937 gen(1234);
76std::uniform_real_distribution<> uniform_dist(0.0, 1.0);
77std::normal_distribution<double> normal_dist;
82 double output = double(uniform_dist(gen));
88 static u64 output = gen();
101 return normal_dist(gen);
111 return g_lehmer64_state;
133 return double(a) / double(u64(1) << 53U);
138 static const double pi2 = std::atan(
static_cast<double>(1)) * 8;
143 double ln = std::log(u1);
144 double sqrt = std::sqrt(-2 * ln);
146 return sqrt * std::cos(pi2 * u2);
150template<
typename Scalar>
156 for (isize i = 0; i < nrows; ++i) {
162template<
typename Scalar>
168 for (isize i = 0; i < nrows; ++i) {
169 for (isize j = 0; j < ncols; ++j) {
178template<
typename Scalar>
182 auto mat = rand::matrix_rand<Scalar>(n, n);
183 auto qr = mat.householderQr();
187using Input = std::pair<u128, isize>;
190template<
typename Scalar>
195 static auto cache = std::map<detail::Input, Mat<Scalar, colmajor>>{};
197 auto it = cache.find(input);
198 if (it == cache.end()) {
199 auto res = cache.insert({
201 detail::orthonormal_rand_impl<Scalar>(n),
208template<
typename Scalar>
212 auto const& q = rand::orthonormal_rand<Scalar>(n);
218 Scalar diff = log(cond);
219 for (isize i = 0; i < n; ++i) {
220 d(i) = exp(Scalar(i) / Scalar(n) * diff);
224 return q * d.asDiagonal() * q.transpose();
227template<
typename Scalar>
235 for (isize i = 0; i < n; ++i) {
237 H.insert(i, i) = random;
240 for (isize i = 0; i < n; ++i) {
241 for (isize j = i + 1; j < n; ++j) {
243 if (urandom < p / 2) {
245 H.insert(i, j) = random;
252 H_dense.template selfadjointView<Eigen::Upper>().eigenvalues();
254 Scalar min = eigh.minCoeff();
255 Scalar max = eigh.maxCoeff();
265 Scalar rho = (max - min) / (cond - 1) - min;
270 for (isize i = 0; i < n; ++i) {
271 H.coeffRef(i, i) += rho;
278template<
typename Scalar>
288 for (isize i = 0; i < n; ++i) {
289 for (isize j = i + 1; j < n; ++j) {
293 H.insert(i, j) = random;
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));
309template<
typename Scalar>
318 for (isize i = 0; i < n; ++i) {
319 for (isize j = 0; j < n; ++j) {
321 if (urandom < p / 2) {
328 H = ((H + H.transpose()) * 0.5)
332 Vec<Scalar> eigh = H.template selfadjointView<Eigen::Upper>().eigenvalues();
333 Scalar min = eigh.minCoeff();
334 H.diagonal().array() += (rho + abs(min));
339template<
typename Scalar>
345 for (isize i = 0; i < nrows; ++i) {
346 for (isize j = 0; j < ncols; ++j) {
356template<
typename Scalar>
364 for (isize i = 0; i < nrows; ++i) {
365 for (isize j = 0; j < ncols; ++j) {
390 return lhs.operator*(rhs);
392template<
typename To,
typename From>
396 return from.template cast<To>();
402template<
typename MatLhs,
typename MatRhs,
typename T =
typename MatLhs::Scalar>
406 using Upscaled =
typename std::
407 conditional<std::is_floating_point<T>::value,
long double, T>::type;
409 return mat_cast<T, Upscaled>(matmul_impl<Upscaled>(
414template<
typename MatLhs,
417 typename T =
typename MatLhs::Scalar>
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); }
443template<
typename Scalar>
446 Scalar sparsity_factor,
447 Scalar strong_convexity_factor = Scalar(1e-2))
451 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
452 dim, strong_convexity_factor, sparsity_factor);
455 rand::sparse_matrix_rand_not_compressed<Scalar>(0, dim, sparsity_factor);
458 rand::sparse_matrix_rand_not_compressed<Scalar>(0, dim, sparsity_factor);
467template<
typename Scalar>
472 Scalar sparsity_factor,
473 Scalar strong_convexity_factor = Scalar(1e-2))
477 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
478 dim, strong_convexity_factor, sparsity_factor);
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);
485 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
488 for (proxqp::isize i = 0; i < n_in; ++i) {
509template<
typename Scalar>
514 Scalar sparsity_factor)
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);
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);
530 for (proxqp::isize i = 0; i < n_in; ++i) {
537 Vec<Scalar> g = -(H * x_sol + C.transpose() * z_sol + A.transpose() * y_sol);
550template<
typename Scalar>
555 Scalar sparsity_factor,
556 Scalar strong_convexity_factor = Scalar(1e-2))
560 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
561 dim, strong_convexity_factor, sparsity_factor);
564 rand::sparse_matrix_rand_not_compressed<Scalar>(n_eq, dim, sparsity_factor);
567 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
570 for (proxqp::isize i = 0; i < 2 * n_in; ++i) {
576 rand::sparse_matrix_rand_not_compressed<Scalar>(n_in, dim, sparsity_factor);
578 C.block(0, 0, n_in, dim) = C_;
579 C.block(n_in, 0, n_in, dim) = C_;
596template<
typename Scalar>
601 Scalar sparsity_factor,
602 Scalar strong_convexity_factor = Scalar(1e-2))
606 rand::sparse_positive_definite_rand_not_compressed<Scalar>(
607 dim, strong_convexity_factor, sparsity_factor);
610 rand::sparse_matrix_rand_not_compressed<Scalar>(n_eq, dim, sparsity_factor);
613 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
616 for (proxqp::isize i = 0; i < n_in; ++i) {
621 C.diagonal().array() += 1;
635template<
typename Scalar>
640 Scalar sparsity_factor,
641 Scalar strong_convexity_factor = Scalar(1e-2))
645 dim, strong_convexity_factor, sparsity_factor);
648 rand::sparse_matrix_rand<Scalar>(n_eq, dim, sparsity_factor);
650 rand::sparse_matrix_rand<Scalar>(n_in, dim, sparsity_factor);
652 Vec<Scalar> x_sol = rand::vector_rand<Scalar>(dim);
655 for (proxqp::isize i = 0; i < n_in; ++i) {
#define LDLT_EXPLICIT_TPL_DECL(NParams,...)
#define VEG_TAG(Name, Type)
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 >
std::pair< u128, isize > Input
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 >
auto normal_rand() -> double
auto lehmer_global() -> u128 &
auto uniform_rand() -> double
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
constexpr Layout colmajor
This class stores the model of the QP problem.
EigenNoAlloc(EigenNoAlloc const &)=delete
auto operator=(EigenNoAlloc const &) -> EigenNoAlloc &=delete
EigenNoAlloc(EigenNoAlloc &&)=delete
auto operator=(EigenNoAlloc &&) -> EigenNoAlloc &=delete