20    std::string s(16, 
'0');
 
   21    auto begin = s.data(), end = begin + s.size();
 
   22    auto [ptr, 
ec] = std::to_chars(begin, end, 
version, 16);
 
   23    if (
ec != std::errc())
 
   24        throw std::logic_error(std::make_error_code(
ec).message());
 
   25    std::rotate(begin, ptr, end);
 
 
   33        throw std::runtime_error(
 
   34            "alpaqa::dl::DLProblem::DLProblem: " 
   35            "Incompatible problem definition (problem ABI version 0x" +
 
 
   41std::shared_ptr<void> load_lib(
const std::string &
so_filename) {
 
   46        throw std::runtime_error(
err);
 
   47    return std::shared_ptr<void>{h, &
::dlclose};
 
 
   54    auto *h = 
::dlsym(handle, name.c_str());
 
   56        throw std::runtime_error(
"Unable to load function '" + name +
 
   59    return reinterpret_cast<F *
>(h);
 
 
   79                .cols     = 
sp.dense.cols,
 
   80                .symmetry = 
static_cast<Symmetry
>(
sp.dense.symmetry),
 
   85                .
rows = 
sp.sparse_csc.rows,
 
   86                .cols = 
sp.sparse_csc.cols,
 
   87                .symmetry = 
static_cast<Symmetry
>(
sp.sparse_csc.symmetry),
 
   88                .inner_idx = 
typename SparseCSC::index_vector_map_t{
sp.sparse_csc.inner_idx, 
sp.sparse_csc.nnz},
 
   89                .outer_ptr = 
typename SparseCSC::index_vector_map_t{
sp.sparse_csc.outer_ptr, 
sp.sparse_csc.cols + 1},
 
   90                .order = 
static_cast<SparseCSC::Order
>(
sp.sparse_csc.order),
 
   95                .
rows = 
sp.sparse_csc_l.rows,
 
   96                .cols = 
sp.sparse_csc_l.cols,
 
   97                .symmetry = 
static_cast<Symmetry
>(
sp.sparse_csc_l.symmetry),
 
   98                .inner_idx = 
typename SparseCSCl::index_vector_map_t{
sp.sparse_csc_l.inner_idx, 
sp.sparse_csc_l.nnz},
 
   99                .outer_ptr = 
typename SparseCSCl::index_vector_map_t{
sp.sparse_csc_l.outer_ptr, 
sp.sparse_csc_l.cols + 1},
 
  100                .order = 
static_cast<SparseCSCl::Order
>(
sp.sparse_csc_l.order),
 
  105                .
rows = 
sp.sparse_csc_ll.rows,
 
  106                .cols = 
sp.sparse_csc_ll.cols,
 
  107                .symmetry = 
static_cast<Symmetry
>(
sp.sparse_csc_ll.symmetry),
 
  108                .inner_idx = 
typename SparseCSCll::index_vector_map_t{
sp.sparse_csc_ll.inner_idx, 
sp.sparse_csc_ll.nnz},
 
  109                .outer_ptr = 
typename SparseCSCll::index_vector_map_t{
sp.sparse_csc_ll.outer_ptr, 
sp.sparse_csc_ll.cols + 1},
 
  110                .order = 
static_cast<SparseCSCll::Order
>(
sp.sparse_csc_ll.order),
 
  115                .
rows = 
sp.sparse_coo.rows,
 
  116                .cols = 
sp.sparse_coo.cols,
 
  117                .symmetry = 
static_cast<Symmetry
>(
sp.sparse_coo.symmetry),
 
  118                .row_indices = 
typename SparseCOO::index_vector_map_t{
sp.sparse_coo.row_indices, 
sp.sparse_coo.nnz},
 
  119                .col_indices = 
typename SparseCOO::index_vector_map_t{
sp.sparse_coo.col_indices, 
sp.sparse_coo.nnz},
 
  120                .order = 
static_cast<SparseCOO::Order
>(
sp.sparse_coo.order),
 
  121                .first_index = 
sp.sparse_coo.first_index,
 
  126                .
rows = 
sp.sparse_coo_l.rows,
 
  127                .cols = 
sp.sparse_coo_l.cols,
 
  128                .symmetry = 
static_cast<Symmetry
>(
sp.sparse_coo_l.symmetry),
 
  129                .row_indices = 
typename SparseCOOl::index_vector_map_t{
sp.sparse_coo_l.row_indices, 
sp.sparse_coo_l.nnz},
 
  130                .col_indices = 
typename SparseCOOl::index_vector_map_t{
sp.sparse_coo_l.col_indices, 
sp.sparse_coo_l.nnz},
 
  131                .order = 
static_cast<SparseCOOl::Order
>(
sp.sparse_coo_l.order),
 
  132                .first_index = 
sp.sparse_coo_l.first_index,
 
  137                .
rows = 
sp.sparse_coo_ll.rows,
 
  138                .cols = 
sp.sparse_coo_ll.cols,
 
  139                .symmetry = 
static_cast<Symmetry
>(
sp.sparse_coo_ll.symmetry),
 
  140                .row_indices = 
typename SparseCOOll::index_vector_map_t{
sp.sparse_coo_ll.row_indices, 
sp.sparse_coo_ll.nnz},
 
  141                .col_indices = 
typename SparseCOOll::index_vector_map_t{
sp.sparse_coo_ll.col_indices, 
sp.sparse_coo_ll.nnz},
 
  142                .order = 
static_cast<SparseCOOll::Order
>(
sp.sparse_coo_ll.order),
 
  143                .first_index = 
sp.sparse_coo_ll.first_index,
 
  145        default: 
throw std::invalid_argument(
"Invalid sparsity kind");
 
 
 
  161    std::unique_ptr<alpaqa_function_dict_t> 
unique_extra{r.extra_functions};
 
  163    check_abi_version(r.abi_version);
 
  178        throw std::logic_error(
"alpaqa::dl::DLProblem::DLProblem: plugin did " 
  179                               "not return any functions");
 
  187    this->
C = 
Box{this->
n};
 
  188    this->
D = 
Box{this->
m};
 
  191                                    this->C.upperbound.data());
 
  194                                    this->D.upperbound.data());
 
 
  208    if (functions->eval_prox_grad_step)
 
  209        return functions->eval_prox_grad_step(
 
  210            instance.get(), γ, x.data(), grad_ψ.data(), x̂.data(), p.data());
 
 
  225auto DLProblem::eval_hess_ψ_prod(
crvec x, 
crvec y, 
crvec Σ, 
real_t scale, 
crvec v, 
rvec Hv) 
const -> 
void { 
return functions->eval_hess_ψ_prod(instance.get(), x.data(), y.data(), Σ.data(), 
scale, D.lowerbound.data(), D.upperbound.data(), 
v.data(), 
Hv.data()); }
 
  226auto DLProblem::eval_hess_ψ(
crvec x, 
crvec y, 
crvec Σ, 
real_t scale, 
rvec H_values) 
const -> 
void { 
return functions->eval_hess_ψ(instance.get(), x.data(), y.data(), Σ.data(), 
scale, D.lowerbound.data(), D.upperbound.data(), 
H_values.size() == 0 ? 
nullptr : 
H_values.data()); }
 
  233auto DLProblem::eval_grad_ψ(
crvec x, 
crvec y, 
crvec Σ, 
rvec grad_ψ, 
rvec work_n, 
rvec work_m) 
const -> 
void { 
return functions->eval_grad_ψ(instance.get(), x.data(), y.data(), Σ.data(), D.lowerbound.data(), D.upperbound.data(), grad_ψ.data(), work_n.data(), work_m.data()); }
 
  234auto DLProblem::eval_ψ_grad_ψ(
crvec x, 
crvec y, 
crvec Σ, 
rvec grad_ψ, 
rvec work_n, 
rvec work_m) 
const -> 
real_t { 
return functions->eval_ψ_grad_ψ(instance.get(), x.data(), y.data(), Σ.data(), D.lowerbound.data(), D.upperbound.data(), grad_ψ.data(), work_n.data(), work_m.data()); }
 
  270    std::unique_ptr<alpaqa_function_dict_t> 
unique_extra{r.extra_functions};
 
  272    check_abi_version(r.abi_version);
 
  287        throw std::logic_error(
"alpaqa::dl::DLControlProblem::DLControlProblem:" 
  288                               " plugin did not return any functions");
 
 
  311auto DLControlProblem::eval_add_R_masked(
index_t timestep, 
crvec xu, 
crvec h, 
crindexvec mask, 
rmat R, 
rvec work) 
const -> 
void { 
return functions->eval_add_R_masked(instance.get(), 
timestep, xu.data(), h.data(), 
mask.data(), R.data(), work.data()); }
 
  312auto DLControlProblem::eval_add_S_masked(
index_t timestep, 
crvec xu, 
crvec h, 
crindexvec mask, 
rmat S, 
rvec work) 
const -> 
void { 
return functions->eval_add_S_masked(instance.get(), 
timestep, xu.data(), h.data(), 
mask.data(), S.data(), work.data()); }
 
  313auto DLControlProblem::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 -> 
void { 
return functions->eval_add_R_prod_masked(instance.get(), 
timestep, xu.data(), h.data(), 
mask_J.data(), 
mask_K.data(), 
v.data(), 
out.data(), work.data()); }
 
  314auto DLControlProblem::eval_add_S_prod_masked(
index_t timestep, 
crvec xu, 
crvec h, 
crindexvec mask_K, 
crvec v, 
rvec out, 
rvec work) 
const -> 
void { 
return functions->eval_add_S_prod_masked(instance.get(), 
timestep, xu.data(), h.data(), 
mask_K.data(), 
v.data(), 
out.data(), work.data()); }
 
Implements common problem functions for minimization problems with box constraints.
Box C
Constraints of the decision variables, .
length_t m
Number of constraints, dimension of g(x) and z.
length_t n
Number of decision variables, dimension of x.
Box D
Other constraints, .
real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x̂, rvec p) const
vec l1_reg
(1-norm) regularization parameter.
void eval_add_S_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat S, rvec work) const
bool provides_eval_grad_constr_prod() const
void eval_add_Q_N(crvec x, crvec h, rmat Q) const
bool provides_eval_add_gn_hess_constr() const
DLControlProblem(const std::string &so_filename, const std::string &function_name="register_alpaqa_control_problem", void *user_param=nullptr)
Load a problem from a shared library.
void eval_q_N(crvec x, crvec h, rvec q) const
void eval_add_gn_hess_constr(index_t timestep, crvec x, crvec M, rmat out) const
void eval_grad_constr_prod(index_t timestep, crvec x, crvec p, rvec grad_cx_p) 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
bool provides_eval_add_R_prod_masked() const
length_t get_S_work_size() const
void eval_constr_N(crvec x, rvec c) const
void get_D_N(Box &D) const
bool provides_get_D() const
bool provides_get_S_work_size() const
void eval_grad_f_prod(index_t timestep, crvec x, crvec u, crvec p, rvec grad_fxu_p) const
bool provides_eval_add_Q_N() const
std::shared_ptr< void > instance
Problem instance created by the registration function, including the deleter to destroy it.
void eval_jac_f(index_t timestep, crvec x, crvec u, rmat J_fxu) const
real_t eval_l_N(crvec h) const
void eval_h_N(crvec x, rvec h) const
real_t eval_l(index_t timestep, crvec h) const
void eval_grad_constr_prod_N(crvec x, crvec p, rvec grad_cx_p) const
bool provides_eval_constr() const
control_problem_functions_t * functions
Pointer to the struct of function pointers for evaluating the objective, constraints,...
bool provides_get_D_N() const
void eval_add_gn_hess_constr_N(crvec x, crvec M, rmat out) const
void eval_f(index_t timestep, crvec x, crvec u, rvec fxu) const
bool provides_eval_grad_constr_prod_N() const
length_t get_R_work_size() const
void eval_add_S_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_K, crvec v, rvec out, rvec work) const
bool provides_eval_add_S_prod_masked() const
void get_x_init(rvec x_init) const
void eval_h(index_t timestep, crvec x, crvec u, rvec h) const
void eval_add_R_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat R, rvec work) const
std::shared_ptr< void > handle
Handle to the shared module defining the problem.
void eval_qr(index_t timestep, crvec xu, crvec h, rvec qr) const
ExtraFuncs extra_funcs
Dictionary of extra functions that were registered by the problem.
bool provides_eval_constr_N() const
void eval_add_Q(index_t timestep, crvec xu, crvec h, rmat Q) const
bool provides_get_R_work_size() const
bool provides_eval_add_gn_hess_constr_N() const
void eval_constr(index_t timestep, crvec x, rvec c) const
bool provides_eval_hess_L() const
real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x̂, rvec p) const
real_t eval_ψ_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const
void eval_grad_gi(crvec x, index_t i, rvec grad_gi) const
bool provides_get_hess_L_sparsity() const
void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const
bool provides_eval_hess_ψ_prod() const
void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const
bool provides_eval_ψ_grad_ψ() const
bool provides_get_box_C() const
Sparsity get_jac_g_sparsity() const
real_t eval_f_g(crvec x, rvec g) const
void eval_grad_f(crvec x, rvec grad_fx) const
Sparsity get_hess_ψ_sparsity() const
bool provides_eval_jac_g() const
Sparsity get_hess_L_sparsity() const
void eval_grad_f_grad_g_prod(crvec x, crvec y, rvec grad_f, rvec grad_gxy) const
real_t eval_ψ(crvec x, crvec y, crvec Σ, rvec ŷ) const
std::shared_ptr< void > instance
Problem instance created by the registration function, including the deleter to destroy it.
bool provides_eval_g() const
void eval_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const
void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const
bool provides_eval_grad_f_grad_g_prod() const
DLProblem(const std::string &so_filename, const std::string &function_name="register_alpaqa_problem", void *user_param=nullptr)
Load a problem from a shared library.
bool provides_eval_f() const
bool provides_get_hess_ψ_sparsity() const
bool provides_eval_hess_L_prod() const
bool provides_eval_grad_g_prod() const
bool provides_get_jac_g_sparsity() const
real_t eval_f_grad_f(crvec x, rvec grad_fx) const
bool provides_eval_f_grad_f() const
bool provides_eval_grad_f() const
bool provides_eval_grad_gi() const
bool provides_eval_f_g() const
problem_functions_t * functions
Pointer to the struct of function pointers for evaluating the objective, constraints,...
void eval_jac_g(crvec x, rvec J_values) const
real_t eval_f(crvec x) const
void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const
bool provides_eval_grad_L() const
void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const
bool provides_eval_grad_ψ() const
std::shared_ptr< void > handle
Handle to the shared module defining the problem.
ExtraFuncs extra_funcs
Dictionary of extra functions that were registered by the problem.
bool provides_eval_hess_ψ() const
void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const
void eval_g(crvec x, rvec gx) const
bool provides_eval_ψ() const
#define USING_ALPAQA_CONFIG(Conf)
#define ALPAQA_DL_ABI_VERSION
void check_abi_version(uint64_t abi_version)
std::string format_abi_version(uint64_t version)
std::list< std::shared_ptr< void > > leaked_modules
F * load_func(void *handle, const std::string &name)
Sparsity< Conf > convert_sparsity(alpaqa_sparsity_t sp)
void leak_lib(std::shared_ptr< void > handle)
std::mutex leaked_modules_mutex
Symmetry
Describes the symmetry of matrices.
typename Conf::real_t real_t
typename Conf::index_t index_t
typename Conf::length_t length_t
typename Conf::crvec crvec
typename Conf::crindexvec crindexvec
Sparse coordinate list structure (COO).
Sparse compressed-column structure (CCS or CSC).
Stores any of the supported sparsity patterns.
void(* eval_add_gn_hess_constr)(void *instance, alpaqa_index_t timestep, const alpaqa_real_t *x, const alpaqa_real_t *M, alpaqa_real_t *out)
alpaqa_length_t(* get_S_work_size)(void *instance)
void(* eval_grad_constr_prod_N)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *p, alpaqa_real_t *grad_cx_p)
void(* get_D)(void *instance, alpaqa_real_t *lb, alpaqa_real_t *ub)
void(* eval_grad_constr_prod)(void *instance, alpaqa_index_t timestep, const alpaqa_real_t *x, const alpaqa_real_t *p, alpaqa_real_t *grad_cx_p)
void(* eval_add_R_prod_masked)(void *instance, alpaqa_index_t timestep, const alpaqa_real_t *xu, const alpaqa_real_t *h, const alpaqa_index_t *mask_J, const alpaqa_index_t *mask_K, const alpaqa_real_t *v, alpaqa_real_t *out, alpaqa_real_t *work)
void(* eval_add_gn_hess_constr_N)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *M, alpaqa_real_t *out)
alpaqa_length_t(* get_R_work_size)(void *instance)
void(* get_D_N)(void *instance, alpaqa_real_t *lb, alpaqa_real_t *ub)
void(* eval_constr)(void *instance, alpaqa_index_t timestep, const alpaqa_real_t *x, alpaqa_real_t *c)
void(* eval_add_Q_N)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *h, alpaqa_real_t *Q)
void(* eval_add_S_prod_masked)(void *instance, alpaqa_index_t timestep, const alpaqa_real_t *xu, const alpaqa_real_t *h, const alpaqa_index_t *mask_K, const alpaqa_real_t *v, alpaqa_real_t *out, alpaqa_real_t *work)
void(* eval_constr_N)(void *instance, const alpaqa_real_t *x, alpaqa_real_t *c)
void(* eval_grad_f)(void *instance, const alpaqa_real_t *x, alpaqa_real_t *grad_fx)
Gradient of the cost function.
alpaqa_real_t(* eval_ψ)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, const alpaqa_real_t *Σ, const alpaqa_real_t *zl, const alpaqa_real_t *zu, alpaqa_real_t *ŷ)
Augmented Lagrangian.
void(* eval_hess_ψ)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, const alpaqa_real_t *Σ, alpaqa_real_t scale, const alpaqa_real_t *zl, const alpaqa_real_t *zu, alpaqa_real_t *H_values)
Hessian of the augmented Lagrangian.
void(* eval_hess_ψ_prod)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, const alpaqa_real_t *Σ, alpaqa_real_t scale, const alpaqa_real_t *zl, const alpaqa_real_t *zu, const alpaqa_real_t *v, alpaqa_real_t *Hv)
Hessian-vector product of the augmented Lagrangian.
void(* initialize_box_C)(void *instance, alpaqa_real_t *lb, alpaqa_real_t *ub)
Provide the initial values for the bounds of alpaqa::BoxConstrProblem::C, i.e.
alpaqa_real_t(* eval_prox_grad_step)(void *instance, alpaqa_real_t γ, const alpaqa_real_t *x, const alpaqa_real_t *grad_ψ, alpaqa_real_t *x̂, alpaqa_real_t *p)
Proximal gradient step.
alpaqa_real_t(* eval_f_g)(void *instance, const alpaqa_real_t *x, alpaqa_real_t *g)
Cost and constraints.
void(* eval_grad_f_grad_g_prod)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, alpaqa_real_t *grad_f, alpaqa_real_t *grad_gxy)
Gradient of the cost and gradient-vector product of the constraints.
void(* eval_grad_L)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, alpaqa_real_t *grad_L, alpaqa_real_t *work_n)
Gradient of the Lagrangian.
alpaqa_real_t(* eval_f_grad_f)(void *instance, const alpaqa_real_t *x, alpaqa_real_t *grad_fx)
Cost and its gradient.
alpaqa_real_t(* eval_f)(void *instance, const alpaqa_real_t *x)
Cost function.
alpaqa_sparsity_t(* get_hess_L_sparsity)(void *instance)
Get the sparsity pattern of the Hessian of the Lagrangian.
alpaqa_sparsity_t(* get_hess_ψ_sparsity)(void *instance)
Get the sparsity pattern of the Hessian of the augmented Lagrangian.
alpaqa_length_t m
Number of constraints.
void(* eval_grad_ψ)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, const alpaqa_real_t *Σ, const alpaqa_real_t *zl, const alpaqa_real_t *zu, alpaqa_real_t *grad_ψ, alpaqa_real_t *work_n, alpaqa_real_t *work_m)
Gradient of the augmented Lagrangian.
alpaqa_length_t n
Number of decision variables.
void(* eval_grad_gi)(void *instance, const alpaqa_real_t *x, alpaqa_index_t i, alpaqa_real_t *grad_gi)
Gradient of specific constraint function.
void(* initialize_l1_reg)(void *instance, alpaqa_real_t *lambda, alpaqa_length_t *size)
Provide the initial values for alpaqa::BoxConstrProblem::l1_reg, the ℓ₁-regularization factor.
void(* eval_hess_L_prod)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, alpaqa_real_t scale, const alpaqa_real_t *v, alpaqa_real_t *Hv)
Hessian-vector product of the Lagrangian.
void(* eval_g)(void *instance, const alpaqa_real_t *x, alpaqa_real_t *gx)
Constraints function.
alpaqa_sparsity_t(* get_jac_g_sparsity)(void *instance)
Get the sparsity pattern of the Jacobian of the constraints function.
void(* initialize_box_D)(void *instance, alpaqa_real_t *lb, alpaqa_real_t *ub)
Provide the initial values for the bounds of alpaqa::BoxConstrProblem::D, i.e.
void(* eval_jac_g)(void *instance, const alpaqa_real_t *x, alpaqa_real_t *J_values)
Jacobian of the constraints function.
alpaqa_real_t(* eval_ψ_grad_ψ)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, const alpaqa_real_t *Σ, const alpaqa_real_t *zl, const alpaqa_real_t *zu, alpaqa_real_t *grad_ψ, alpaqa_real_t *work_n, alpaqa_real_t *work_m)
Augmented Lagrangian and its gradient.
void(* eval_hess_L)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, alpaqa_real_t scale, alpaqa_real_t *H_values)
Hessian of the Lagrangian.
void(* eval_grad_g_prod)(void *instance, const alpaqa_real_t *x, const alpaqa_real_t *y, alpaqa_real_t *grad_gxy)
Gradient-vector product of the constraints function.
@ alpaqa_sparsity_sparse_csc
@ alpaqa_sparsity_sparse_coo_l
@ alpaqa_sparsity_sparse_coo_ll
@ alpaqa_sparsity_sparse_csc_ll
@ alpaqa_sparsity_sparse_csc_l
@ alpaqa_sparsity_sparse_coo