11#include <Eigen/Sparse>
12#include <casadi/core/external.hpp>
21namespace fs = std::filesystem;
23namespace casadi_loader {
56 template <
class Loader>
58 {
loader(name) } -> std::same_as<casadi::Function>;
59 {
loader.format_name(name) } -> std::same_as<std::string>;
61 static std::unique_ptr<CasADiControlFunctionsWithParam>
66 using namespace std::literals::string_literals;
69 "Invalid number of input arguments: got "s +
70 std::to_string(
ffun.n_in()) +
", should be 3.");
71 if (
ffun.n_out() != 1)
73 "Invalid number of output arguments: got "s +
74 std::to_string(
ffun.n_in()) +
", should be 1.");
85 using namespace std::literals::string_literals;
88 "Invalid number of input arguments: got "s +
89 std::to_string(
hfun.n_in()) +
", should be 3.");
90 if (
hfun.n_out() != 1)
92 "Invalid number of output arguments: got "s +
93 std::to_string(
hfun.n_in()) +
", should be 1.");
102 using namespace std::literals::string_literals;
103 if (
hfun.n_in() != 2)
105 "Invalid number of input arguments: got "s +
106 std::to_string(
hfun.n_in()) +
", should be 2.");
107 if (
hfun.n_out() != 1)
109 "Invalid number of output arguments: got "s +
110 std::to_string(
hfun.n_in()) +
", should be 1.");
118 using namespace std::literals::string_literals;
119 if (
cfun.n_in() != 2)
121 "Invalid number of input arguments: got "s +
122 std::to_string(
cfun.n_in()) +
", should be 2.");
123 if (
cfun.n_out() != 1)
125 "Invalid number of output arguments: got "s +
126 std::to_string(
cfun.n_in()) +
", should be 1.");
134 using namespace std::literals::string_literals;
135 if (
cfun.n_in() != 2)
137 "Invalid number of input arguments: got "s +
138 std::to_string(
cfun.n_in()) +
", should be 2.");
139 if (
cfun.n_out() != 1)
141 "Invalid number of output arguments: got "s +
142 std::to_string(
cfun.n_in()) +
", should be 1.");
155 auto self = std::make_unique<CasADiControlFunctionsWithParam<Conf>>(
173 .h_N = std::move(
h_N),
196 .c_N = std::move(
c_N),
211template <Config Conf>
218 auto operator()(
const std::string &name)
const {
219 return casadi::external(name,
filename);
240 impl->Q.fun.sparsity_out(0).nnz(),
241 impl->Q_N.fun.sparsity_out(0).nnz(),
242 impl->gn_hess_c.fun.sparsity_out(0).nnz(),
243 impl->gn_hess_c_N.fun.sparsity_out(0).nnz(),
252template <Config Conf>
254 const std::filesystem::path &
filepath,
char sep) {
258 throw std::runtime_error(
"Unable to open data file \"" +
268 return static_cast<void>(
data_file.get());
272 auto s = csv::read_row_std_vector<real_t>(
data_file, sep);
277 throw std::runtime_error(
"Unable to read " + std::string(name) +
278 " from data file \"" +
filepath.string() +
279 ':' + std::to_string(
line) +
284 auto read_single = [&](std::string_view name,
auto &
v) {
287 throw std::runtime_error(
"Unable to read " + std::string(name) +
288 " from data file \"" +
filepath.string() +
289 ':' + std::to_string(
line) +
'"');
300 read_single(
"penalty_alm_split", this->penalty_alm_split);
301 read_single(
"penalty_alm_split_N", this->penalty_alm_split_N);
304template <Config Conf>
307template <Config Conf>
311template <Config Conf>
327 impl->f({x.data(), u.data(), param.data()}, {
fxu.data()});
329template <Config Conf>
336 impl->jac_f({x.data(), u.data(), param.data()}, {
J_fxu.data()});
338template <Config Conf>
346 impl->grad_f_prod({x.data(), u.data(), param.data(), p.data()},
349template <Config Conf>
355 impl->h({x.data(), u.data(), param.data()}, {h.data()});
357template <Config Conf>
361 impl->h_N({x.data(), param.data()}, {h.data()});
363template <Config Conf>
367 impl->l({h.data(), param.data()}, {&l});
370template <Config Conf>
374 impl->l_N({h.data(), param.data()}, {&l});
377template <Config Conf>
380 assert(xu.size() == nx + nu);
382 assert(qr.size() == nx + nu);
383 impl->qr({xu.data(), h.data(), param.data()}, {qr.data()});
385template <Config Conf>
390 impl->q_N({x.data(), h.data(), param.data()}, {q.data()});
392template <Config Conf>
395 assert(xu.size() == nx + nu);
399 impl->Q({xu.data(), h.data(), param.data()}, {work.data()});
400 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
401 using cmspmat = Eigen::Map<const spmat>;
402 auto &&sparse = impl->Q.fun.sparsity_out(0);
403 if (sparse.is_dense())
404 Q +=
cmmat{work.data(), nx, nx};
409 static_cast<length_t>(sparse.nnz()),
415template <Config Conf>
421 impl->Q_N({x.data(), h.data(), param.data()}, {work.data()});
422 auto &&sparse = impl->Q_N.fun.sparsity_out(0);
423 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
424 using cmspmat = Eigen::Map<const spmat>;
425 if (sparse.is_dense())
426 Q +=
cmmat{work.data(), nx, nx};
431 static_cast<length_t>(sparse.nnz()),
438template <Config Conf>
442 auto &&sparse = impl->R.fun.sparsity_out(0);
443 assert(xu.size() == nx + nu);
450 impl->R({xu.data(), h.data(), param.data()}, {work.data()});
451 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
452 using cmspmat = Eigen::Map<const spmat>;
453 if (sparse.is_dense()) {
460 static_cast<length_t>(sparse.nnz()),
469template <Config Conf>
473 auto &&sparse = impl->S.fun.sparsity_out(0);
474 assert(xu.size() == nx + nu);
480 impl->S({xu.data(), h.data(), param.data()}, {work.data()});
481 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
482 using cmspmat = Eigen::Map<const spmat>;
483 using Eigen::indexing::all;
484 if (sparse.is_dense()) {
491 static_cast<length_t>(sparse.nnz()),
500template <Config Conf>
506 auto &&sparse = impl->R.fun.sparsity_out(0);
510 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
511 using cmspmat = Eigen::Map<const spmat>;
512 if (sparse.is_dense()) {
513 auto R =
cmmat{work.data(), nu, nu};
519 static_cast<length_t>(sparse.nnz()),
529template <Config Conf>
534 auto &&sparse = impl->S.fun.sparsity_out(0);
538 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
539 using cmspmat = Eigen::Map<const spmat>;
540 using Eigen::indexing::all;
541 if (sparse.is_dense()) {
542 auto Sᵀ =
cmmat{work.data(), nu, nx}.transpose();
548 static_cast<length_t>(sparse.nnz()),
558template <Config Conf>
560 auto &&sparse = impl->R.fun.sparsity_out(0);
561 return static_cast<length_t>(sparse.nnz());
564template <Config Conf>
566 auto &&sparse = impl->S.fun.sparsity_out(0);
567 return static_cast<length_t>(sparse.nnz());
570template <Config Conf>
576 impl->c({x.data(), param.data()}, {c.data()});
579template <Config Conf>
586 impl->grad_c_prod({x.data(), param.data(), p.data()}, {
grad_cx_p.data()});
589template <Config Conf>
593 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
599 impl->gn_hess_c({x.data(), param.data(),
M.data()}, {work.data()});
600 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
601 using cmspmat = Eigen::Map<const spmat>;
602 if (sparse.is_dense())
608 static_cast<length_t>(sparse.nnz()),
615template <Config Conf>
621 impl->c_N({x.data(), param.data()}, {c.data()});
624template <Config Conf>
630 impl->grad_c_prod_N({x.data(), param.data(), p.data()}, {
grad_cx_p.data()});
633template <Config Conf>
636 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
642 impl->gn_hess_c_N({x.data(), param.data(),
M.data()}, {work.data()});
643 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
644 using cmspmat = Eigen::Map<const spmat>;
645 if (sparse.is_dense())
651 static_cast<length_t>(sparse.nnz()),
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 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)
void read_row(std::istream &is, Eigen::Ref< Eigen::VectorX< float > > 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