daqp_runtime.hpp
1 //
2 // MIT No Attribution
3 //
4 // Copyright (C) 2010-2023 Joel Andersson, Joris Gillis, Moritz Diehl, KU Leuven.
5 //
6 // Permission is hereby granted, free of charge, to any person obtaining a copy of this
7 // software and associated documentation files (the "Software"), to deal in the Software
8 // without restriction, including without limitation the rights to use, copy, modify,
9 // merge, publish, distribute, sublicense, and/or sell copies of the Software, and to
10 // permit persons to whom the Software is furnished to do so.
11 //
12 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED,
13 // INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A
14 // PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT
15 // HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION
16 // OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
17 // SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
18 //
19 
20 // C-REPLACE "casadi_qp_prob<T1>" "struct casadi_qp_prob"
21 // C-REPLACE "casadi_qp_data<T1>" "struct casadi_qp_data"
22 
23 // C-REPLACE "reinterpret_cast<int**>" "(int**) "
24 // C-REPLACE "reinterpret_cast<int*>" "(int*) "
25 // C-REPLACE "const_cast<DAQPSettings*>" "(DAQPSettings*) "
26 // C-REPLACE "const_cast<T1*>" "(casadi_real*) "
27 
28 template<typename T1>
31 
32  DAQPSettings settings;
33 };
34 // C-REPLACE "casadi_daqp_prob<T1>" "struct casadi_daqp_prob"
35 
36 // SYMBOL "daqp_setup"
37 template<typename T1>
38 void casadi_daqp_setup(casadi_daqp_prob<T1>* p) {
39 
40 }
41 
42 
43 
44 // SYMBOL "daqp_data"
45 template<typename T1>
47  // Problem structure
49  // Problem structure
51 
52  DAQPWorkspace work;
53  DAQPProblem daqp;
54  DAQPResult res;
55 
57 
58 };
59 // C-REPLACE "casadi_daqp_data<T1>" "struct casadi_daqp_data"
60 
61 // SYMBOL "daqp_init_mem"
62 template<typename T1>
63 int daqp_init_mem(casadi_daqp_data<T1>* d) {
64  return 0;
65 }
66 
67 // SYMBOL "daqp_free_mem"
68 template<typename T1>
69 void daqp_free_mem(casadi_daqp_data<T1>* d) {
70 
71 }
72 
73 // SYMBOL "daqp_work"
74 template<typename T1>
75 void casadi_daqp_work(const casadi_daqp_prob<T1>* p, casadi_int* sz_arg, casadi_int* sz_res, casadi_int* sz_iw, casadi_int* sz_w) {
76  casadi_qp_work(p->qp, sz_arg, sz_res, sz_iw, sz_w);
77  // Local variables
78  casadi_int nx, na;
79  const casadi_qp_prob<T1>* p_qp = p->qp;
80  // Get matrix number of nonzeros
81  na = p_qp->na;
82  nx = p_qp->nx;
83 
84  // Reset sz_w, sz_iw
85  *sz_w = *sz_iw = 0;
86  // Temporary work vectors
87  // Persistent work vectors
88  *sz_w += nx*nx; // H
89  *sz_w += na*nx; // A
90  *sz_w += p->qp->nz; // lbz
91  *sz_w += p->qp->nz; // ubz
92  *sz_w += p->qp->nz; // lam
93  *sz_iw += p->qp->nz; // sense
94 }
95 
96 // SYMBOL "daqp_init"
97 template<typename T1>
98 void casadi_daqp_init(casadi_daqp_data<T1>* d, const T1*** arg, T1*** res, casadi_int** iw, T1** w) {
99  // Local variables
100  casadi_int nx, na;
101  const casadi_daqp_prob<T1>* p = d->prob;
102  const casadi_qp_prob<T1>* p_qp = p->qp;
103  // Get matrix number of nonzeros
104  na = p_qp->na;
105  nx = p_qp->nx;
106 
107  d->daqp.H = *w; *w += nx*nx;
108  d->daqp.A = *w; *w += na*nx;
109  d->daqp.blower = *w; *w += p->qp->nz;
110  d->daqp.bupper = *w; *w += p->qp->nz;
111  d->res.lam = *w; *w += p->qp->nz;
112  d->daqp.sense = (int*) *iw; *iw += p->qp->nz;
113 }
114 
115 
116 // C-REPLACE "SOLVER_RET_SUCCESS" "0"
117 // C-REPLACE "SOLVER_RET_LIMITED" "2"
118 
119 // SYMBOL "daqp_solve"
120 template<typename T1>
121 int casadi_daqp_solve(casadi_daqp_data<T1>* d, const double** arg, double** res, casadi_int* iw, double* w) {
122  // Local variables
123  casadi_int i;
124  int flag;
125  const casadi_daqp_prob<T1>* p = d->prob;
126  const casadi_qp_prob<T1>* p_qp = p->qp;
127  casadi_qp_data<T1>* d_qp = d->qp;
128 
129  for (i=0;i<p_qp->nx;++i) {
130  d->daqp.sense[i] = d_qp->lbx[i]==d_qp->ubx[i] ? 5 : 0;
131  }
132  for (i=0;i<p_qp->na;++i) {
133  d->daqp.sense[p_qp->nx+i] = d_qp->lba[i]==d_qp->uba[i] ? 5 : 0;
134  }
135 
136  d->daqp.n = p_qp->nx;
137  d->daqp.m = p_qp->nx + p_qp->na;
138  d->daqp.ms = p_qp->nx;
139  d->work.settings = const_cast<DAQPSettings*>(&p->settings);
140  d->res.x = d_qp->x;
141  casadi_densify(d->qp->h, d->prob->qp->sp_h, d->daqp.H, 0);
142  casadi_densify(d->qp->a, d->prob->qp->sp_a, d->daqp.A, 1);
143  casadi_copy(d_qp->lbx, p_qp->nx, d->daqp.blower);
144  casadi_copy(d_qp->lba, p_qp->na, d->daqp.blower+p_qp->nx);
145  casadi_copy(d_qp->ubx, p_qp->nx, d->daqp.bupper);
146  casadi_copy(d_qp->uba, p_qp->na, d->daqp.bupper+p_qp->nx);
147  d->daqp.f = const_cast<T1*>(d_qp->g);
148 
149  flag = setup_daqp(&d->daqp,&d->work,&(d->res.setup_time));
150  if (flag<0) return 1;
151  daqp_solve(&d->res,&d->work);
152  casadi_copy(d->res.lam, p_qp->nx, d_qp->lam_x);
153  casadi_copy(d->res.lam+p_qp->nx, p_qp->na, d_qp->lam_a);
154  if (d_qp->f) *d_qp->f = d->res.fval;
155  d->work.settings = 0;
156  free_daqp_workspace(&d->work);
157  free_daqp_ldp(&d->work);
158 
159  d->return_status = d->res.exitflag;
160  d_qp->success = d->res.exitflag==EXIT_OPTIMAL;
161 
173  return 0;
174 }
DAQPProblem daqp
DAQPWorkspace work
const casadi_daqp_prob< T1 > * prob
casadi_qp_data< T1 > * qp
DAQPSettings settings
const casadi_qp_prob< T1 > * qp
const T1 * lba
Definition: casadi_qp.hpp:64
const T1 * lbx
Definition: casadi_qp.hpp:64
const T1 * uba
Definition: casadi_qp.hpp:64
const T1 * ubx
Definition: casadi_qp.hpp:64
const T1 * g
Definition: casadi_qp.hpp:64
casadi_int nx
Definition: casadi_qp.hpp:33
casadi_int na
Definition: casadi_qp.hpp:33