13#include <Eigen/Sparse>
15#if ALPAQA_WITH_EXTERNAL_CASADI
16#include <casadi/core/external.hpp>
28namespace fs = std::filesystem;
30namespace casadi_loader {
63 template <
class Loader>
65 {
loader(name) } -> std::same_as<casadi::Function>;
66 {
loader.format_name(name) } -> std::same_as<std::string>;
68 static std::unique_ptr<CasADiControlFunctionsWithParam>
73 using namespace std::literals::string_literals;
76 "Invalid number of input arguments: got "s +
77 std::to_string(
ffun.n_in()) +
", should be 3.");
78 if (
ffun.n_out() != 1)
80 "Invalid number of output arguments: got "s +
81 std::to_string(
ffun.n_in()) +
", should be 1.");
92 using namespace std::literals::string_literals;
95 "Invalid number of input arguments: got "s +
96 std::to_string(
hfun.n_in()) +
", should be 3.");
97 if (
hfun.n_out() != 1)
99 "Invalid number of output arguments: got "s +
100 std::to_string(
hfun.n_in()) +
", should be 1.");
109 using namespace std::literals::string_literals;
110 if (
hfun.n_in() != 2)
112 "Invalid number of input arguments: got "s +
113 std::to_string(
hfun.n_in()) +
", should be 2.");
114 if (
hfun.n_out() != 1)
116 "Invalid number of output arguments: got "s +
117 std::to_string(
hfun.n_in()) +
", should be 1.");
125 using namespace std::literals::string_literals;
126 if (
cfun.n_in() != 2)
128 "Invalid number of input arguments: got "s +
129 std::to_string(
cfun.n_in()) +
", should be 2.");
130 if (
cfun.n_out() != 1)
132 "Invalid number of output arguments: got "s +
133 std::to_string(
cfun.n_in()) +
", should be 1.");
141 using namespace std::literals::string_literals;
142 if (
cfun.n_in() != 2)
144 "Invalid number of input arguments: got "s +
145 std::to_string(
cfun.n_in()) +
", should be 2.");
146 if (
cfun.n_out() != 1)
148 "Invalid number of output arguments: got "s +
149 std::to_string(
cfun.n_in()) +
", should be 1.");
162 auto self = std::make_unique<CasADiControlFunctionsWithParam<Conf>>(
180 .h_N = std::move(
h_N),
203 .c_N = std::move(
c_N),
218template <Config Conf>
227 auto operator()(
const std::string &name)
const {
228#if ALPAQA_WITH_EXTERNAL_CASADI
253 impl->Q.fun.sparsity_out(0).nnz(),
254 impl->Q_N.fun.sparsity_out(0).nnz(),
255 impl->gn_hess_c.fun.sparsity_out(0).nnz(),
256 impl->gn_hess_c_N.fun.sparsity_out(0).nnz(),
265template <Config Conf>
267 const std::filesystem::path &
filepath,
char sep) {
271 throw std::runtime_error(
"Unable to open data file \"" +
281 return static_cast<void>(
data_file.get());
285 auto s = csv::read_row_std_vector<real_t>(
data_file, sep);
290 throw std::runtime_error(
"Unable to read " + std::string(name) +
291 " from data file \"" +
filepath.string() +
292 ':' + std::to_string(
line) +
297 auto read_single = [&](std::string_view name,
auto &
v) {
300 throw std::runtime_error(
"Unable to read " + std::string(name) +
301 " from data file \"" +
filepath.string() +
302 ':' + std::to_string(
line) +
'"');
313 read_single(
"penalty_alm_split", this->penalty_alm_split);
314 read_single(
"penalty_alm_split_N", this->penalty_alm_split_N);
317template <Config Conf>
320template <Config Conf>
324template <Config Conf>
340 impl->f({x.data(), u.data(), param.data()}, {
fxu.data()});
342template <Config Conf>
349 impl->jac_f({x.data(), u.data(), param.data()}, {
J_fxu.data()});
351template <Config Conf>
359 impl->grad_f_prod({x.data(), u.data(), param.data(), p.data()},
362template <Config Conf>
368 impl->h({x.data(), u.data(), param.data()}, {h.data()});
370template <Config Conf>
374 impl->h_N({x.data(), param.data()}, {h.data()});
376template <Config Conf>
380 impl->l({h.data(), param.data()}, {&l});
383template <Config Conf>
387 impl->l_N({h.data(), param.data()}, {&l});
390template <Config Conf>
393 assert(xu.size() == nx + nu);
395 assert(qr.size() == nx + nu);
396 impl->qr({xu.data(), h.data(), param.data()}, {qr.data()});
398template <Config Conf>
403 impl->q_N({x.data(), h.data(), param.data()}, {q.data()});
405template <Config Conf>
408 assert(xu.size() == nx + nu);
412 impl->Q({xu.data(), h.data(), param.data()}, {work.data()});
413 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
414 using cmspmat = Eigen::Map<const spmat>;
415 auto &&sparse = impl->Q.fun.sparsity_out(0);
416 if (sparse.is_dense())
417 Q +=
cmmat{work.data(), nx, nx};
422 static_cast<length_t>(sparse.nnz()),
428template <Config Conf>
434 impl->Q_N({x.data(), h.data(), param.data()}, {work.data()});
435 auto &&sparse = impl->Q_N.fun.sparsity_out(0);
436 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
437 using cmspmat = Eigen::Map<const spmat>;
438 if (sparse.is_dense())
439 Q +=
cmmat{work.data(), nx, nx};
444 static_cast<length_t>(sparse.nnz()),
451template <Config Conf>
455 auto &&sparse = impl->R.fun.sparsity_out(0);
456 assert(xu.size() == nx + nu);
463 impl->R({xu.data(), h.data(), param.data()}, {work.data()});
464 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
465 using cmspmat = Eigen::Map<const spmat>;
466 if (sparse.is_dense()) {
473 static_cast<length_t>(sparse.nnz()),
482template <Config Conf>
486 auto &&sparse = impl->S.fun.sparsity_out(0);
487 assert(xu.size() == nx + nu);
493 impl->S({xu.data(), h.data(), param.data()}, {work.data()});
494 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
495 using cmspmat = Eigen::Map<const spmat>;
496 using Eigen::indexing::all;
497 if (sparse.is_dense()) {
504 static_cast<length_t>(sparse.nnz()),
513template <Config Conf>
519 auto &&sparse = impl->R.fun.sparsity_out(0);
523 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
524 using cmspmat = Eigen::Map<const spmat>;
525 if (sparse.is_dense()) {
526 auto R =
cmmat{work.data(), nu, nu};
532 static_cast<length_t>(sparse.nnz()),
542template <Config Conf>
547 auto &&sparse = impl->S.fun.sparsity_out(0);
551 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
552 using cmspmat = Eigen::Map<const spmat>;
553 using Eigen::indexing::all;
554 if (sparse.is_dense()) {
555 auto Sᵀ =
cmmat{work.data(), nu, nx}.transpose();
561 static_cast<length_t>(sparse.nnz()),
571template <Config Conf>
573 auto &&sparse = impl->R.fun.sparsity_out(0);
574 return static_cast<length_t>(sparse.nnz());
577template <Config Conf>
579 auto &&sparse = impl->S.fun.sparsity_out(0);
580 return static_cast<length_t>(sparse.nnz());
583template <Config Conf>
589 impl->c({x.data(), param.data()}, {c.data()});
592template <Config Conf>
599 impl->grad_c_prod({x.data(), param.data(), p.data()}, {
grad_cx_p.data()});
602template <Config Conf>
606 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
612 impl->gn_hess_c({x.data(), param.data(),
M.data()}, {work.data()});
613 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
614 using cmspmat = Eigen::Map<const spmat>;
615 if (sparse.is_dense())
616 out +=
cmmat{work.data(), nx, nx};
621 static_cast<length_t>(sparse.nnz()),
628template <Config Conf>
634 impl->c_N({x.data(), param.data()}, {c.data()});
637template <Config Conf>
643 impl->grad_c_prod_N({x.data(), param.data(), p.data()}, {
grad_cx_p.data()});
646template <Config Conf>
649 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
655 impl->gn_hess_c_N({x.data(), param.data(),
M.data()}, {work.data()});
656 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
657 using cmspmat = Eigen::Map<const spmat>;
658 if (sparse.is_dense())
659 out +=
cmmat{work.data(), nx, nx};
664 static_cast<length_t>(sparse.nnz()),
#define BEGIN_ALPAQA_CASADI_LOADER_NAMESPACE
#define END_ALPAQA_CASADI_LOADER_NAMESPACE
void eval_jac_f(index_t timestep, crvec x, crvec u, rmat J_fxu) const
void eval_add_R_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_J, crindexvec mask_K, crvec v, rvec out, rvec work) const
CasADiControlProblem & operator=(const CasADiControlProblem &)
void eval_add_gn_hess_constr_N(crvec x, crvec M, rmat out) const
length_t get_S_work_size() const
void eval_qr(index_t timestep, crvec xu, crvec h, rvec qr) const
void eval_add_S_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_K, crvec v, rvec out, rvec work) const
void eval_constr_N(crvec x, rvec c) const
void load_numerical_data(const std::filesystem::path &filepath, char sep=',')
Load the numerical problem data (bounds and parameters) from a CSV file.
void eval_grad_constr_prod_N(crvec x, crvec p, rvec grad_cx_p) const
void eval_add_R_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat R, rvec work) const
real_t eval_l_N(crvec h) const
real_t eval_l(index_t timestep, crvec h) const
void eval_grad_f_prod(index_t timestep, crvec x, crvec u, crvec p, rvec grad_fxu_p) const
void eval_add_gn_hess_constr(index_t timestep, crvec x, crvec M, rmat out) const
void eval_constr(index_t timestep, crvec x, rvec c) const
CasADiControlProblem(const std::string &filename, length_t N, DynamicLoadFlags dl_flags={})
void eval_h(index_t timestep, crvec x, crvec u, rvec h) const
util::copyable_unique_ptr< Functions > impl
length_t get_R_work_size() const
void eval_add_S_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat S, rvec work) const
void eval_q_N(crvec x, crvec h, rvec q) const
void eval_add_Q_N(crvec x, crvec h, rmat Q) const
void eval_grad_constr_prod(index_t timestep, crvec x, crvec p, rvec grad_cx_p) const
void eval_add_Q(index_t timestep, crvec xu, crvec h, rmat Q) const
void eval_h_N(crvec x, rvec h) const
Class that loads and calls pre-compiled CasADi functions in a DLL/SO file.
Class for evaluating CasADi functions, allocating the necessary workspace storage in advance for allo...
static void validate_dimensions(const casadi::Function &fun, const std::array< casadi_dim, N_in > &dim_in={}, const std::array< casadi_dim, N_out > &dim_out={})
#define USING_ALPAQA_CONFIG(Conf)
auto wrapped_load(Loader &&loader, const char *name, Args &&...args)
constexpr auto dims(auto... a)
std::pair< casadi_int, casadi_int > dim
auto wrap_load(Loader &&loader, const char *name, F f)
Function external(const std::string &name, const std::string &bin_name, DynamicLoadFlags dl_flags)
Load the given CasADi function from the given DLL/SO file.
void read_row(std::istream &is, Eigen::Ref< Eigen::VectorX< Eigen::Index > > v, char sep)
void sparse_matvec_add_masked_rows_cols(const SpMat &R, const CVec &v, Vec &&out, const MaskVec &mask_J, const MaskVec &mask_K)
out += R(mask_J,mask_K) * v(mask_K);
void sparse_add_masked_rows(const SpMat &S_full, Mat &&S, const MaskVec &mask)
S += S_full(mask,:)
void sparse_add_masked(const SpMat &R_full, Mat &&R, const MaskVec &mask)
R += R_full(mask,mask)
void sparse_matvec_add_transpose_masked_rows(const SpMat &S, const CVec &v, Vec &&out, const MaskVec &mask)
out += S(mask,:)ᵀ * v(mask);
typename Conf::cmmat cmmat
typename Conf::real_t real_t
typename Conf::index_t index_t
typename Conf::length_t length_t
typename Conf::cmvec cmvec
typename Conf::crvec crvec
typename Conf::crindexvec crindexvec
Flags to be passed to dlopen.
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > qr
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > jac_f
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > f
CasADiFunctionEvaluator< Conf, 1+WithParam, 1 > c
static std::unique_ptr< CasADiControlFunctionsWithParam > load(Loader &&loader)
CasADiFunctionEvaluator< Conf, 3+WithParam, 1 > grad_f_prod
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > R
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > grad_c_prod
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > q_N
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > grad_c_prod_N
CasADiFunctionEvaluator< Conf, 1+WithParam, 1 > h_N
CasADiFunctionEvaluator< Conf, 1+WithParam, 1 > l
static constexpr bool WithParam
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > Q
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > S
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > gn_hess_c_N
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > gn_hess_c
CasADiFunctionEvaluator< Conf, 1+WithParam, 1 > l_N
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > Q_N
CasADiFunctionEvaluator< Conf, 1+WithParam, 1 > c_N
CasADiFunctionEvaluator< Conf, 2+WithParam, 1 > h