12#include <Eigen/Sparse>
14#if ALPAQA_WITH_EXTERNAL_CASADI
15#include <casadi/core/external.hpp>
26namespace fs = std::filesystem;
28namespace casadi_loader {
30using namespace alpaqa::casadi_loader;
36 static constexpr bool WithParam =
true;
37 length_t nx, nu, nh, nh_N,
nc, nc_N, p;
61 template <
class Loader>
62 requires requires(Loader &&loader,
const char *name) {
63 { loader(name) } -> std::same_as<casadi::Function>;
64 { loader.format_name(name) } -> std::same_as<std::string>;
66 static std::unique_ptr<CasADiControlFunctionsWithParam>
68 length_t nx, nu, nh, nh_N, nc, nc_N, p;
70 casadi::Function ffun = loader(
"f");
71 using namespace std::literals::string_literals;
74 "Invalid number of input arguments: got "s +
75 std::to_string(ffun.n_in()) +
", should be 3.");
76 if (ffun.n_out() != 1)
78 "Invalid number of output arguments: got "s +
79 std::to_string(ffun.n_in()) +
", should be 1.");
80 nx =
static_cast<length_t
>(ffun.size1_in(0));
81 nu =
static_cast<length_t
>(ffun.size1_in(1));
82 p =
static_cast<length_t
>(ffun.size1_in(2));
89 casadi::Function hfun = loader(
"h");
90 using namespace std::literals::string_literals;
93 "Invalid number of input arguments: got "s +
94 std::to_string(hfun.n_in()) +
", should be 3.");
95 if (hfun.n_out() != 1)
97 "Invalid number of output arguments: got "s +
98 std::to_string(hfun.n_in()) +
", should be 1.");
99 nh =
static_cast<length_t
>(hfun.size1_out(0));
106 casadi::Function hfun = loader(
"h_N");
107 using namespace std::literals::string_literals;
108 if (hfun.n_in() != 2)
110 "Invalid number of input arguments: got "s +
111 std::to_string(hfun.n_in()) +
", should be 2.");
112 if (hfun.n_out() != 1)
114 "Invalid number of output arguments: got "s +
115 std::to_string(hfun.n_in()) +
", should be 1.");
116 nh_N =
static_cast<length_t
>(hfun.size1_out(0));
122 casadi::Function cfun = loader(
"c");
123 using namespace std::literals::string_literals;
124 if (cfun.n_in() != 2)
126 "Invalid number of input arguments: got "s +
127 std::to_string(cfun.n_in()) +
", should be 2.");
128 if (cfun.n_out() != 1)
130 "Invalid number of output arguments: got "s +
131 std::to_string(cfun.n_in()) +
", should be 1.");
132 nc =
static_cast<length_t
>(cfun.size1_out(0));
138 casadi::Function cfun = loader(
"c_N");
139 using namespace std::literals::string_literals;
140 if (cfun.n_in() != 2)
142 "Invalid number of input arguments: got "s +
143 std::to_string(cfun.n_in()) +
", should be 2.");
144 if (cfun.n_out() != 1)
146 "Invalid number of output arguments: got "s +
147 std::to_string(cfun.n_in()) +
", should be 1.");
148 nc_N =
static_cast<length_t
>(cfun.size1_out(0));
156 auto h_N =
wrap_load(loader,
"h_N", load_h_N);
158 auto c_N =
wrap_load(loader,
"c_N", load_c_N);
160 auto self = std::make_unique<CasADiControlFunctionsWithParam<Conf>>(
170 .jac_f = wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
171 loader,
"jacobian_f",
dims(nx, nu, p),
175 loader,
"grad_f_prod",
dims(nx, nu, p, nx),
178 .h_N = std::move(h_N),
179 .l = wrapped_load<CasADiFunctionEvaluator<Conf, 2, 1>>(
182 loader,
"l_N",
dims(nh_N, p),
dims(1)),
184 loader,
"qr",
dims(nx + nu, nh, p),
dims(nx + nu)),
186 loader,
"q_N",
dims(nx, nh_N, p),
dims(nx)),
188 loader,
"Q",
dims(nx + nu, nh, p),
dims(
dim{nx, nx})),
189 .Q_N = wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
190 loader,
"Q_N",
dims(nx, nh_N, p),
dims(
dim{nx, nx})),
191 .R = wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
192 loader,
"R",
dims(nx + nu, nh, p),
dims(
dim{nu, nu})),
193 .S = wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
194 loader,
"S",
dims(nx + nu, nh, p),
dims(
dim{nu, nx})),
197 wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
198 loader,
"grad_c_prod",
dims(nx, p, nc),
dims(nx)),
200 loader,
"gn_hess_c",
dims(nx, p, nc),
dims(
dim{nx, nx})),
201 .c_N = std::move(c_N),
203 wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
204 loader,
"grad_c_prod_N",
dims(nx, p, nc_N),
dims(nx)),
207 loader,
"gn_hess_c_N",
dims(nx, p, nc_N),
216template <Config Conf>
222 const std::string &filename;
223 auto operator()(
const std::string &name)
const {
224 return casadi::external(name, filename);
226 auto format_name(
const std::string &name)
const {
227 return filename +
':' + name;
244 auto n_work = std::max({
245 impl->Q.fun.sparsity_out(0).nnz(),
246 impl->Q_N.fun.sparsity_out(0).nnz(),
247 impl->gn_hess_c.fun.sparsity_out(0).nnz(),
248 impl->gn_hess_c_N.fun.sparsity_out(0).nnz(),
250 this->
work = vec::Constant(
static_cast<length_t
>(n_work), NaN<Conf>);
252 auto bounds_filepath = fs::path{filename}.replace_extension(
"csv");
253 if (fs::exists(bounds_filepath))
257template <Config Conf>
259 const std::filesystem::path &filepath,
char sep) {
261 std::ifstream data_file{filepath};
263 throw std::runtime_error(
"Unable to open data file \"" +
264 filepath.string() +
'"');
268 auto wrap_data_load = [&](std::string_view name,
auto &v,
269 bool fixed_size =
true) {
272 if (data_file.peek() ==
'\n')
273 return static_cast<void>(data_file.get());
275 csv::read_row(data_file, v, sep);
277 auto s = csv::read_row_std_vector<real_t>(data_file, sep);
278 v = cmvec{s.data(),
static_cast<index_t
>(s.size())};
280 }
catch (csv::read_error &e) {
282 throw std::runtime_error(
"Unable to read " + std::string(name) +
283 " from data file \"" + filepath.string() +
284 ':' + std::to_string(line) +
289 auto read_single = [&](std::string_view name,
auto &v) {
292 throw std::runtime_error(
"Unable to read " + std::string(name) +
293 " from data file \"" + filepath.string() +
294 ':' + std::to_string(line) +
'"');
296 wrap_data_load(
"U.lowerbound", this->U.lowerbound);
297 wrap_data_load(
"U.upperbound", this->U.upperbound);
298 wrap_data_load(
"D.lowerbound", this->D.lowerbound);
299 wrap_data_load(
"D.upperbound", this->D.upperbound);
300 wrap_data_load(
"D_N.lowerbound", this->D_N.lowerbound);
301 wrap_data_load(
"D_N.upperbound", this->D_N.upperbound);
302 wrap_data_load(
"x_init", this->x_init);
303 wrap_data_load(
"param", this->param);
305 read_single(
"penalty_alm_split", this->penalty_alm_split);
306 read_single(
"penalty_alm_split_N", this->penalty_alm_split_N);
309template <Config Conf>
312template <Config Conf>
316template <Config Conf>
319template <Config Conf>
323template <Config Conf>
326template <Config Conf>
329 assert(x.size() == nx);
330 assert(u.size() == nu);
331 assert(fxu.size() == nx);
332 impl->f({x.data(), u.data(), param.data()}, {fxu.data()});
334template <Config Conf>
337 assert(x.size() == nx);
338 assert(u.size() == nu);
339 assert(J_fxu.rows() == nx);
340 assert(J_fxu.cols() == nx + nu);
341 impl->jac_f({x.data(), u.data(), param.data()}, {J_fxu.data()});
343template <Config Conf>
346 rvec grad_fxu_p)
const {
347 assert(x.size() == nx);
348 assert(u.size() == nu);
349 assert(p.size() == nx);
350 assert(grad_fxu_p.size() == nx + nu);
351 impl->grad_f_prod({x.data(), u.data(), param.data(), p.data()},
352 {grad_fxu_p.data()});
354template <Config Conf>
357 assert(x.size() == nx);
358 assert(u.size() == nu);
359 assert(h.size() == nh);
360 impl->h({x.data(), u.data(), param.data()}, {h.data()});
362template <Config Conf>
364 assert(x.size() == nx);
365 assert(h.size() == nh_N);
366 impl->h_N({x.data(), param.data()}, {h.data()});
368template <Config Conf>
370 assert(h.size() == nh);
372 impl->l({h.data(), param.data()}, {&l});
375template <Config Conf>
377 assert(h.size() == nh_N);
379 impl->l_N({h.data(), param.data()}, {&l});
382template <Config Conf>
385 assert(xu.size() == nx + nu);
386 assert(h.size() == nh);
387 assert(qr.size() == nx + nu);
388 impl->qr({xu.data(), h.data(), param.data()}, {qr.data()});
390template <Config Conf>
392 assert(x.size() == nx);
393 assert(h.size() == nh_N);
394 assert(q.size() == nx);
395 impl->q_N({x.data(), h.data(), param.data()}, {q.data()});
397template <Config Conf>
400 assert(xu.size() == nx + nu);
401 assert(h.size() == nh);
402 assert(Q.rows() == nx);
403 assert(Q.cols() == nx);
404 impl->Q({xu.data(), h.data(), param.data()}, {work.data()});
405 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
406 using cmspmat = Eigen::Map<const spmat>;
407 auto &&sparse = impl->Q.fun.sparsity_out(0);
408 if (sparse.is_dense())
409 Q += cmmat{work.data(), nx, nx};
414 static_cast<length_t
>(sparse.nnz()),
420template <Config Conf>
422 assert(x.size() == nx);
423 assert(h.size() == nh_N);
424 assert(Q.rows() == nx);
425 assert(Q.cols() == nx);
426 impl->Q_N({x.data(), h.data(), param.data()}, {work.data()});
427 auto &&sparse = impl->Q_N.fun.sparsity_out(0);
428 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
429 using cmspmat = Eigen::Map<const spmat>;
430 if (sparse.is_dense())
431 Q += cmmat{work.data(), nx, nx};
436 static_cast<length_t
>(sparse.nnz()),
443template <Config Conf>
445 crindexvec mask, rmat R,
447 auto &&sparse = impl->R.fun.sparsity_out(0);
448 assert(xu.size() == nx + nu);
449 assert(h.size() == nh);
450 assert(R.rows() <= nu);
451 assert(R.cols() <= nu);
452 assert(R.rows() == mask.size());
453 assert(R.cols() == mask.size());
454 assert(work.size() >=
static_cast<length_t
>(sparse.nnz()));
455 impl->R({xu.data(), h.data(), param.data()}, {work.data()});
456 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
457 using cmspmat = Eigen::Map<const spmat>;
458 if (sparse.is_dense()) {
459 cmmat R_full{work.data(), nu, nu};
460 R += R_full(mask, mask);
465 static_cast<length_t
>(sparse.nnz()),
470 util::sparse_add_masked(R_full, R, mask);
474template <Config Conf>
476 crindexvec mask, rmat S,
478 auto &&sparse = impl->S.fun.sparsity_out(0);
479 assert(xu.size() == nx + nu);
480 assert(h.size() == nh);
481 assert(S.rows() <= nu);
482 assert(S.rows() == mask.size());
483 assert(S.cols() == nx);
484 assert(work.size() >=
static_cast<length_t
>(sparse.nnz()));
485 impl->S({xu.data(), h.data(), param.data()}, {work.data()});
486 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
487 using cmspmat = Eigen::Map<const spmat>;
488 using Eigen::indexing::all;
489 if (sparse.is_dense()) {
490 cmmat S_full{work.data(), nu, nx};
491 S += S_full(mask, all);
496 static_cast<length_t
>(sparse.nnz()),
501 util::sparse_add_masked_rows(S_full, S, mask);
505template <Config Conf>
511 auto &&sparse = impl->R.fun.sparsity_out(0);
512 assert(v.size() == nu);
513 assert(out.size() == mask_J.size());
514 assert(work.size() >=
static_cast<length_t
>(sparse.nnz()));
515 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
516 using cmspmat = Eigen::Map<const spmat>;
517 if (sparse.is_dense()) {
518 auto R = cmmat{work.data(), nu, nu};
519 out.noalias() += R(mask_J, mask_K) * v(mask_K);
524 static_cast<length_t
>(sparse.nnz()),
530 util::sparse_matvec_add_masked_rows_cols(R, v, out, mask_J, mask_K);
534template <Config Conf>
539 auto &&sparse = impl->S.fun.sparsity_out(0);
540 assert(v.size() == nu);
541 assert(out.size() == nx);
542 assert(work.size() >=
static_cast<length_t
>(sparse.nnz()));
543 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
544 using cmspmat = Eigen::Map<const spmat>;
545 using Eigen::indexing::all;
546 if (sparse.is_dense()) {
547 auto Sᵀ = cmmat{work.data(), nu, nx}.transpose();
548 out.noalias() += Sᵀ(all, mask_K) * v(mask_K);
553 static_cast<length_t
>(sparse.nnz()),
559 util::sparse_matvec_add_transpose_masked_rows(S, v, out, mask_K);
563template <Config Conf>
565 auto &&sparse = impl->R.fun.sparsity_out(0);
566 return static_cast<length_t
>(sparse.nnz());
569template <Config Conf>
571 auto &&sparse = impl->S.fun.sparsity_out(0);
572 return static_cast<length_t
>(sparse.nnz());
575template <Config Conf>
579 assert(x.size() == nx);
580 assert(c.size() == nc);
581 impl->c({x.data(), param.data()}, {c.data()});
584template <Config Conf>
587 rvec grad_cx_p)
const {
588 assert(x.size() == nx);
589 assert(p.size() == nc);
590 assert(grad_cx_p.size() == nx);
591 impl->grad_c_prod({x.data(), param.data(), p.data()}, {grad_cx_p.data()});
594template <Config Conf>
598 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
599 assert(x.size() == nx);
600 assert(M.size() == nc);
601 assert(out.rows() == nx);
602 assert(out.cols() == nx);
603 assert(work.size() >=
static_cast<length_t
>(sparse.nnz()));
604 impl->gn_hess_c({x.data(), param.data(), M.data()}, {work.data()});
605 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
606 using cmspmat = Eigen::Map<const spmat>;
607 if (sparse.is_dense())
608 out += cmmat{work.data(), nx, nx};
613 static_cast<length_t
>(sparse.nnz()),
620template <Config Conf>
624 assert(x.size() == nx);
625 assert(c.size() == nc_N);
626 impl->c_N({x.data(), param.data()}, {c.data()});
629template <Config Conf>
631 rvec grad_cx_p)
const {
632 assert(x.size() == nx);
633 assert(p.size() == nc_N);
634 assert(grad_cx_p.size() == nx);
635 impl->grad_c_prod_N({x.data(), param.data(), p.data()}, {grad_cx_p.data()});
638template <Config Conf>
641 auto &&sparse = impl->gn_hess_c.fun.sparsity_out(0);
642 assert(x.size() == nx);
643 assert(M.size() == nc_N);
644 assert(out.rows() == nx);
645 assert(out.cols() == nx);
646 assert(work.size() >=
static_cast<length_t
>(sparse.nnz()));
647 impl->gn_hess_c_N({x.data(), param.data(), M.data()}, {work.data()});
648 using spmat = Eigen::SparseMatrix<real_t, Eigen::ColMajor, casadi_int>;
649 using cmspmat = Eigen::Map<const spmat>;
650 if (sparse.is_dense())
651 out += cmmat{work.data(), nx, nx};
656 static_cast<length_t
>(sparse.nnz()),
#define 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 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)
typename Conf::length_t length_t
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)
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
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