12#include <Eigen/Sparse>
14#if ALPAQA_WITH_EXTERNAL_CASADI
15#include <casadi/core/external.hpp>
27namespace fs = std::filesystem;
29namespace casadi_loader {
62 template <
class Loader>
64 {
loader(name) } -> std::same_as<casadi::Function>;
65 {
loader.format_name(name) } -> std::same_as<std::string>;
67 static std::unique_ptr<CasADiControlFunctionsWithParam>
72 using namespace std::literals::string_literals;
75 "Invalid number of input arguments: got "s +
76 std::to_string(
ffun.n_in()) +
", should be 3.");
77 if (
ffun.n_out() != 1)
79 "Invalid number of output arguments: got "s +
80 std::to_string(
ffun.n_in()) +
", should be 1.");
91 using namespace std::literals::string_literals;
94 "Invalid number of input arguments: got "s +
95 std::to_string(
hfun.n_in()) +
", should be 3.");
96 if (
hfun.n_out() != 1)
98 "Invalid number of output arguments: got "s +
99 std::to_string(
hfun.n_in()) +
", should be 1.");
108 using namespace std::literals::string_literals;
109 if (
hfun.n_in() != 2)
111 "Invalid number of input arguments: got "s +
112 std::to_string(
hfun.n_in()) +
", should be 2.");
113 if (
hfun.n_out() != 1)
115 "Invalid number of output arguments: got "s +
116 std::to_string(
hfun.n_in()) +
", should be 1.");
124 using namespace std::literals::string_literals;
125 if (
cfun.n_in() != 2)
127 "Invalid number of input arguments: got "s +
128 std::to_string(
cfun.n_in()) +
", should be 2.");
129 if (
cfun.n_out() != 1)
131 "Invalid number of output arguments: got "s +
132 std::to_string(
cfun.n_in()) +
", should be 1.");
140 using namespace std::literals::string_literals;
141 if (
cfun.n_in() != 2)
143 "Invalid number of input arguments: got "s +
144 std::to_string(
cfun.n_in()) +
", should be 2.");
145 if (
cfun.n_out() != 1)
147 "Invalid number of output arguments: got "s +
148 std::to_string(
cfun.n_in()) +
", should be 1.");
161 auto self = std::make_unique<CasADiControlFunctionsWithParam<Conf>>(
179 .h_N = std::move(
h_N),
202 .c_N = std::move(
c_N),
217template <Config Conf>
224 auto operator()(
const std::string &name)
const {
246 impl->Q.fun.sparsity_out(0).nnz(),
247 impl->Q_N.fun.sparsity_out(0).nnz(),
248 impl->gn_hess_c.fun.sparsity_out(0).nnz(),
249 impl->gn_hess_c_N.fun.sparsity_out(0).nnz(),
258template <Config Conf>
260 const std::filesystem::path &
filepath,
char sep) {
264 throw std::runtime_error(
"Unable to open data file \"" +
274 return static_cast<void>(
data_file.get());
278 auto s = csv::read_row_std_vector<real_t>(
data_file, sep);
283 throw std::runtime_error(
"Unable to read " + std::string(name) +
284 " from data file \"" +
filepath.string() +
285 ':' + std::to_string(
line) +
290 auto read_single = [&](std::string_view name,
auto &
v) {
293 throw std::runtime_error(
"Unable to read " + std::string(name) +
294 " from data file \"" +
filepath.string() +
295 ':' + std::to_string(
line) +
'"');
306 read_single(
"penalty_alm_split", this->penalty_alm_split);
307 read_single(
"penalty_alm_split_N", this->penalty_alm_split_N);
310template <Config Conf>
313template <Config Conf>
317template <Config Conf>
333 impl->f({x.data(), u.data(), param.data()}, {
fxu.data()});
335template <Config Conf>
342 impl->jac_f({x.data(), u.data(), param.data()}, {
J_fxu.data()});
344template <Config Conf>
352 impl->grad_f_prod({x.data(), u.data(), param.data(), p.data()},
355template <Config Conf>
361 impl->h({x.data(), u.data(), param.data()}, {h.data()});
363template <Config Conf>
367 impl->h_N({x.data(), param.data()}, {h.data()});
369template <Config Conf>
373 impl->l({h.data(), param.data()}, {&l});
376template <Config Conf>
380 impl->l_N({h.data(), param.data()}, {&l});
383template <Config Conf>
386 assert(xu.size() == nx + nu);
388 assert(qr.size() == nx + nu);
389 impl->qr({xu.data(), h.data(), param.data()}, {qr.data()});
391template <Config Conf>
396 impl->q_N({x.data(), h.data(), param.data()}, {q.data()});
398template <Config Conf>
401 assert(xu.size() == nx + nu);
405 impl->Q({xu.data(), h.data(), param.data()}, {work.data()});
406 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
407 using cmspmat = Eigen::Map<const spmat>;
408 auto &&sparse = impl->Q.fun.sparsity_out(0);
409 if (sparse.is_dense())
410 Q +=
cmmat{work.data(), nx, nx};
415 static_cast<length_t>(sparse.nnz()),
421template <Config Conf>
427 impl->Q_N({x.data(), h.data(), param.data()}, {work.data()});
428 auto &&sparse = impl->Q_N.fun.sparsity_out(0);
429 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
430 using cmspmat = Eigen::Map<const spmat>;
431 if (sparse.is_dense())
432 Q +=
cmmat{work.data(), nx, nx};
437 static_cast<length_t>(sparse.nnz()),
444template <Config Conf>
448 auto &&sparse = impl->R.fun.sparsity_out(0);
449 assert(xu.size() == nx + nu);
456 impl->R({xu.data(), h.data(), param.data()}, {work.data()});
457 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
458 using cmspmat = Eigen::Map<const spmat>;
459 if (sparse.is_dense()) {
466 static_cast<length_t>(sparse.nnz()),
475template <Config Conf>
479 auto &&sparse = impl->S.fun.sparsity_out(0);
480 assert(xu.size() == nx + nu);
486 impl->S({xu.data(), h.data(), param.data()}, {work.data()});
487 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
488 using cmspmat = Eigen::Map<const spmat>;
489 using Eigen::indexing::all;
490 if (sparse.is_dense()) {
497 static_cast<length_t>(sparse.nnz()),
506template <Config Conf>
512 auto &&sparse = impl->R.fun.sparsity_out(0);
516 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
517 using cmspmat = Eigen::Map<const spmat>;
518 if (sparse.is_dense()) {
519 auto R =
cmmat{work.data(), nu, nu};
525 static_cast<length_t>(sparse.nnz()),
535template <Config Conf>
540 auto &&sparse = impl->S.fun.sparsity_out(0);
544 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
545 using cmspmat = Eigen::Map<const spmat>;
546 using Eigen::indexing::all;
547 if (sparse.is_dense()) {
548 auto Sᵀ =
cmmat{work.data(), nu, nx}.transpose();
554 static_cast<length_t>(sparse.nnz()),
564template <Config Conf>
566 auto &&sparse = impl->R.fun.sparsity_out(0);
567 return static_cast<length_t>(sparse.nnz());
570template <Config Conf>
572 auto &&sparse = impl->S.fun.sparsity_out(0);
573 return static_cast<length_t>(sparse.nnz());
576template <Config Conf>
582 impl->c({x.data(), param.data()}, {c.data()});
585template <Config Conf>
592 impl->grad_c_prod({x.data(), param.data(), p.data()}, {
grad_cx_p.data()});
595template <Config Conf>
599 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
605 impl->gn_hess_c({x.data(), param.data(),
M.data()}, {work.data()});
606 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
607 using cmspmat = Eigen::Map<const spmat>;
608 if (sparse.is_dense())
609 out +=
cmmat{work.data(), nx, nx};
614 static_cast<length_t>(sparse.nnz()),
621template <Config Conf>
627 impl->c_N({x.data(), param.data()}, {c.data()});
630template <Config Conf>
636 impl->grad_c_prod_N({x.data(), param.data(), p.data()}, {
grad_cx_p.data()});
639template <Config Conf>
642 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
648 impl->gn_hess_c_N({x.data(), param.data(),
M.data()}, {work.data()});
649 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
650 using cmspmat = Eigen::Map<const spmat>;
651 if (sparse.is_dense())
652 out +=
cmmat{work.data(), nx, nx};
657 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
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
CasADiControlProblem(const std::string &filename, length_t N)
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)
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
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