alpaqa no-casadi-dep
Nonconvex constrained optimization
Loading...
Searching...
No Matches
CasADiControlProblem.tpp
Go to the documentation of this file.
1#pragma once
2
10#include "CasADiLoader-util.hpp"
11
12#include <Eigen/Sparse>
13
14#if ALPAQA_WITH_EXTERNAL_CASADI
15#include <casadi/core/external.hpp>
16#endif
17
18#include <fstream>
19#include <memory>
20#include <optional>
21#include <stdexcept>
22#include <type_traits>
23
25
26namespace fs = std::filesystem;
27
28namespace casadi_loader {
29
30using namespace alpaqa::casadi_loader;
31
32template <Config Conf>
35
36 static constexpr bool WithParam = true;
37 length_t nx, nu, nh, nh_N, nc, nc_N, p;
45
52
56
60
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>;
65 }
66 static std::unique_ptr<CasADiControlFunctionsWithParam>
67 load(Loader &&loader) {
68 length_t nx, nu, nh, nh_N, nc, nc_N, p;
69 auto load_f = [&]() -> CasADiFunctionEvaluator<Conf, 3, 1> {
70 casadi::Function ffun = loader("f");
71 using namespace std::literals::string_literals;
72 if (ffun.n_in() != 3)
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));
83 CasADiFunctionEvaluator<Conf, 3, 1> f{std::move(ffun)};
84 f.validate_dimensions({dim(nx, 1), dim(nu, 1), dim(p, 1)},
85 {dim(nx, 1)});
86 return f;
87 };
88 auto load_h = [&]() -> CasADiFunctionEvaluator<Conf, 3, 1> {
89 casadi::Function hfun = loader("h");
90 using namespace std::literals::string_literals;
91 if (hfun.n_in() != 3)
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));
100 CasADiFunctionEvaluator<Conf, 3, 1> h{std::move(hfun)};
101 h.validate_dimensions({dim(nx, 1), dim(nu, 1), dim(p, 1)},
102 {dim(nh, 1)});
103 return h;
104 };
105 auto load_h_N = [&]() -> CasADiFunctionEvaluator<Conf, 2, 1> {
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));
117 CasADiFunctionEvaluator<Conf, 2, 1> h{std::move(hfun)};
118 h.validate_dimensions({dim(nx, 1), dim(p, 1)}, {dim(nh_N, 1)});
119 return h;
120 };
121 auto load_c = [&]() -> CasADiFunctionEvaluator<Conf, 2, 1> {
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));
133 CasADiFunctionEvaluator<Conf, 2, 1> c{std::move(cfun)};
134 c.validate_dimensions({dim(nx, 1), dim(p, 1)}, {dim(nc, 1)});
135 return c;
136 };
137 auto load_c_N = [&]() -> CasADiFunctionEvaluator<Conf, 2, 1> {
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));
149 CasADiFunctionEvaluator<Conf, 2, 1> c{std::move(cfun)};
150 c.validate_dimensions({dim(nx, 1), dim(p, 1)}, {dim(nc_N, 1)});
151 return c;
152 };
153 // Load the functions "f", "h", and "c" to determine the unknown dimensions.
154 auto f = wrap_load(loader, "f", load_f);
155 auto h = wrap_load(loader, "h", load_h);
156 auto h_N = wrap_load(loader, "h_N", load_h_N);
157 auto c = wrap_load(loader, "c", load_c);
158 auto c_N = wrap_load(loader, "c_N", load_c_N);
159
160 auto self = std::make_unique<CasADiControlFunctionsWithParam<Conf>>(
162 .nx = nx,
163 .nu = nu,
164 .nh = nh,
165 .nh_N = nh_N,
166 .nc = nc,
167 .nc_N = nc_N,
168 .p = p,
169 .f = std::move(f),
170 .jac_f = wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
171 loader, "jacobian_f", dims(nx, nu, p),
172 dims(dim(nx, nx + nu))),
173 .grad_f_prod =
175 loader, "grad_f_prod", dims(nx, nu, p, nx),
176 dims(nx + nu)),
177 .h = std::move(h),
178 .h_N = std::move(h_N),
179 .l = wrapped_load<CasADiFunctionEvaluator<Conf, 2, 1>>(
180 loader, "l", dims(nh, p), dims(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})),
195 .c = std::move(c),
196 .grad_c_prod =
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),
202 .grad_c_prod_N =
203 wrapped_load<CasADiFunctionEvaluator<Conf, 3, 1>>(
204 loader, "grad_c_prod_N", dims(nx, p, nc_N), dims(nx)),
205 .gn_hess_c_N =
207 loader, "gn_hess_c_N", dims(nx, p, nc_N),
208 dims(dim{nx, nx})),
209 });
210 return self;
211 }
212};
213
214} // namespace casadi_loader
215
216template <Config Conf>
218 length_t N)
219 : N{N} {
220
221 struct {
222 const std::string &filename;
223 auto operator()(const std::string &name) const {
224 return casadi::external(name, filename);
225 }
226 auto format_name(const std::string &name) const {
227 return filename + ':' + name;
228 }
229 } loader{filename};
231
232 this->nx = impl->nx;
233 this->nu = impl->nu;
234 this->nh = impl->nh;
235 this->nh_N = impl->nh_N;
236 this->nc = impl->nc;
237 this->nc_N = impl->nc_N;
238 this->x_init = vec::Constant(nx, alpaqa::NaN<Conf>);
239 this->param = vec::Constant(impl->p, alpaqa::NaN<Conf>);
240 this->U = Box{nu};
241 this->D = Box{nc};
242 this->D_N = Box{nc_N};
243
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(),
249 });
250 this->work = vec::Constant(static_cast<length_t>(n_work), NaN<Conf>);
251
252 auto bounds_filepath = fs::path{filename}.replace_extension("csv");
253 if (fs::exists(bounds_filepath))
254 load_numerical_data(bounds_filepath);
255}
256
257template <Config Conf>
259 const std::filesystem::path &filepath, char sep) {
260 // Open data file
261 std::ifstream data_file{filepath};
262 if (!data_file)
263 throw std::runtime_error("Unable to open data file \"" +
264 filepath.string() + '"');
265
266 // Helper function for reading single line of (float) data
267 index_t line = 0;
268 auto wrap_data_load = [&](std::string_view name, auto &v,
269 bool fixed_size = true) {
270 try {
271 ++line;
272 if (data_file.peek() == '\n') // Ignore empty lines
273 return static_cast<void>(data_file.get());
274 if (fixed_size) {
275 csv::read_row(data_file, v, sep);
276 } else { // Dynamic size
277 auto s = csv::read_row_std_vector<real_t>(data_file, sep);
278 v = cmvec{s.data(), static_cast<index_t>(s.size())};
279 }
280 } catch (csv::read_error &e) {
281 // Transform any errors in something more readable
282 throw std::runtime_error("Unable to read " + std::string(name) +
283 " from data file \"" + filepath.string() +
284 ':' + std::to_string(line) +
285 "\": " + e.what());
286 }
287 };
288 // Helper function for reading a single value
289 auto read_single = [&](std::string_view name, auto &v) {
290 data_file >> v;
291 if (!data_file)
292 throw std::runtime_error("Unable to read " + std::string(name) +
293 " from data file \"" + filepath.string() +
294 ':' + std::to_string(line) + '"');
295 };
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);
304 // Penalty/ALM split is a single integer
305 read_single("penalty_alm_split", this->penalty_alm_split);
306 read_single("penalty_alm_split_N", this->penalty_alm_split_N);
307}
308
309template <Config Conf>
311 default;
312template <Config Conf>
315
316template <Config Conf>
318 CasADiControlProblem &&) noexcept = default;
319template <Config Conf>
320CasADiControlProblem<Conf> &CasADiControlProblem<Conf>::operator=(
321 CasADiControlProblem &&) noexcept = default;
322
323template <Config Conf>
325
326template <Config Conf>
327void CasADiControlProblem<Conf>::eval_f(index_t, crvec x, crvec u,
328 rvec fxu) const {
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()});
333}
334template <Config Conf>
335void CasADiControlProblem<Conf>::eval_jac_f(index_t, crvec x, crvec u,
336 rmat J_fxu) const {
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()});
342}
343template <Config Conf>
344void CasADiControlProblem<Conf>::eval_grad_f_prod(index_t, crvec x, crvec u,
345 crvec p,
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()});
353}
354template <Config Conf>
355void CasADiControlProblem<Conf>::eval_h(index_t, crvec x, crvec u,
356 rvec h) const {
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()});
361}
362template <Config Conf>
363void CasADiControlProblem<Conf>::eval_h_N(crvec x, rvec h) const {
364 assert(x.size() == nx);
365 assert(h.size() == nh_N);
366 impl->h_N({x.data(), param.data()}, {h.data()});
367}
368template <Config Conf>
369auto CasADiControlProblem<Conf>::eval_l(index_t, crvec h) const -> real_t {
370 assert(h.size() == nh);
371 real_t l;
372 impl->l({h.data(), param.data()}, {&l});
373 return l;
374}
375template <Config Conf>
376auto CasADiControlProblem<Conf>::eval_l_N(crvec h) const -> real_t {
377 assert(h.size() == nh_N);
378 real_t l;
379 impl->l_N({h.data(), param.data()}, {&l});
380 return l;
381}
382template <Config Conf>
383void CasADiControlProblem<Conf>::eval_qr(index_t, crvec xu, crvec h,
384 rvec qr) const {
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()});
389}
390template <Config Conf>
391void CasADiControlProblem<Conf>::eval_q_N(crvec x, crvec h, rvec q) const {
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()});
396}
397template <Config Conf>
398void CasADiControlProblem<Conf>::eval_add_Q(index_t, crvec xu, crvec h,
399 rmat Q) const {
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};
410 else
411 Q += cmspmat{
412 nx,
413 nx,
414 static_cast<length_t>(sparse.nnz()),
415 sparse.colind(),
416 sparse.row(),
417 work.data(),
418 };
419}
420template <Config Conf>
421void CasADiControlProblem<Conf>::eval_add_Q_N(crvec x, crvec h, rmat Q) const {
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};
432 else
433 Q += cmspmat{
434 nx,
435 nx,
436 static_cast<length_t>(sparse.nnz()),
437 sparse.colind(),
438 sparse.row(),
439 work.data(),
440 };
441}
442
443template <Config Conf>
444void CasADiControlProblem<Conf>::eval_add_R_masked(index_t, crvec xu, crvec h,
445 crindexvec mask, rmat R,
446 rvec work) const {
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);
461 } else {
462 cmspmat R_full{
463 nu,
464 nu,
465 static_cast<length_t>(sparse.nnz()),
466 sparse.colind(),
467 sparse.row(),
468 work.data(),
469 };
470 util::sparse_add_masked(R_full, R, mask);
471 }
472}
473
474template <Config Conf>
475void CasADiControlProblem<Conf>::eval_add_S_masked(index_t, crvec xu, crvec h,
476 crindexvec mask, rmat S,
477 rvec work) const {
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);
492 } else {
493 cmspmat S_full{
494 nu,
495 nx,
496 static_cast<length_t>(sparse.nnz()),
497 sparse.colind(),
498 sparse.row(),
499 work.data(),
500 };
501 util::sparse_add_masked_rows(S_full, S, mask);
502 }
503}
504
505template <Config Conf>
507 crindexvec mask_J,
508 crindexvec mask_K,
509 crvec v, rvec out,
510 rvec work) const {
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);
520 } else {
521 cmspmat R{
522 nu,
523 nu,
524 static_cast<length_t>(sparse.nnz()),
525 sparse.colind(),
526 sparse.row(),
527 work.data(),
528 };
529 // out += R_full(mask_J,mask_K) * v(mask_K);
530 util::sparse_matvec_add_masked_rows_cols(R, v, out, mask_J, mask_K);
531 }
532}
533
534template <Config Conf>
536 crindexvec mask_K,
537 crvec v, rvec out,
538 rvec work) const {
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);
549 } else {
550 cmspmat S{
551 nu,
552 nx,
553 static_cast<length_t>(sparse.nnz()),
554 sparse.colind(),
555 sparse.row(),
556 work.data(),
557 };
558 // out += S(mask_K,:)ᵀ * v(mask_K);
559 util::sparse_matvec_add_transpose_masked_rows(S, v, out, mask_K);
560 }
561}
562
563template <Config Conf>
565 auto &&sparse = impl->R.fun.sparsity_out(0);
566 return static_cast<length_t>(sparse.nnz());
567}
568
569template <Config Conf>
571 auto &&sparse = impl->S.fun.sparsity_out(0);
572 return static_cast<length_t>(sparse.nnz());
573}
574
575template <Config Conf>
576void CasADiControlProblem<Conf>::eval_constr(index_t, crvec x, rvec c) const {
577 if (nc == 0)
578 return;
579 assert(x.size() == nx);
580 assert(c.size() == nc);
581 impl->c({x.data(), param.data()}, {c.data()});
582}
583
584template <Config Conf>
586 crvec p,
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()});
592}
593
594template <Config Conf>
596 crvec M,
597 rmat out) const {
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};
609 else
610 out += cmspmat{
611 nx,
612 nx,
613 static_cast<length_t>(sparse.nnz()),
614 sparse.colind(),
615 sparse.row(),
616 work.data(),
617 };
618}
619
620template <Config Conf>
621void CasADiControlProblem<Conf>::eval_constr_N(crvec x, rvec c) const {
622 if (nc_N == 0)
623 return;
624 assert(x.size() == nx);
625 assert(c.size() == nc_N);
626 impl->c_N({x.data(), param.data()}, {c.data()});
627}
628
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()});
636}
637
638template <Config Conf>
640 rmat out) const {
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};
652 else
653 out += cmspmat{
654 nx,
655 nx,
656 static_cast<length_t>(sparse.nnz()),
657 sparse.colind(),
658 sparse.row(),
659 work.data(),
660 };
661}
662
663} // namespace alpaqa::inline ALPAQA_CASADI_LOADER_NAMESPACE
#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
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(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
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)
Definition config.hpp:77
typename Conf::length_t length_t
Definition config.hpp:103
constexpr const auto inf
Definition config.hpp:112
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