blob: cc0d96f29f5c941a8a099c888aec638be3f87b38 [file] [log] [blame]
Austin Schuh16ce3c72018-01-28 16:17:08 -08001/**************************************************************************************************
2* *
3* This file is part of HPIPM. *
4* *
5* HPIPM -- High Performance Interior Point Method. *
6* Copyright (C) 2017 by Gianluca Frison. *
7* Developed at IMTEK (University of Freiburg) under the supervision of Moritz Diehl. *
8* All rights reserved. *
9* *
10* HPMPC is free software; you can redistribute it and/or *
11* modify it under the terms of the GNU Lesser General Public *
12* License as published by the Free Software Foundation; either *
13* version 2.1 of the License, or (at your option) any later version. *
14* *
15* HPMPC is distributed in the hope that it will be useful, *
16* but WITHOUT ANY WARRANTY; without even the implied warranty of *
17* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. *
18* See the GNU Lesser General Public License for more details. *
19* *
20* You should have received a copy of the GNU Lesser General Public *
21* License along with HPMPC; if not, write to the Free Software *
22* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA *
23* *
24* Author: Gianluca Frison, gianluca.frison (at) imtek.uni-freiburg.de *
25* *
26**************************************************************************************************/
27
28#include <stdlib.h>
29#include <stdio.h>
30#include <math.h>
31#include <sys/time.h>
32
33#include <blasfeo_target.h>
34#include <blasfeo_common.h>
35#include <blasfeo_v_aux_ext_dep.h>
36#include <blasfeo_d_aux_ext_dep.h>
37#include <blasfeo_i_aux_ext_dep.h>
38#include <blasfeo_d_aux.h>
39#include <blasfeo_d_blas.h>
40
41#include "../include/hpipm_d_ocp_qp.h"
42#include "../include/hpipm_d_ocp_qp_sol.h"
43#include "../include/hpipm_d_ocp_qp_ipm_hard.h"
44
45#include "d_tools.h"
46
47
48
49#if ! defined(EXT_DEP)
50/* creates a zero matrix */
51void d_zeros(double **pA, int row, int col)
52 {
53 *pA = malloc((row*col)*sizeof(double));
54 double *A = *pA;
55 int i;
56 for(i=0; i<row*col; i++) A[i] = 0.0;
57 }
58/* frees matrix */
59void d_free(double *pA)
60 {
61 free( pA );
62 }
63/* prints a matrix in column-major format */
64void d_print_mat(int m, int n, double *A, int lda)
65 {
66 int i, j;
67 for(i=0; i<m; i++)
68 {
69 for(j=0; j<n; j++)
70 {
71 printf("%9.5f ", A[i+lda*j]);
72 }
73 printf("\n");
74 }
75 printf("\n");
76 }
77/* prints the transposed of a matrix in column-major format */
78void d_print_tran_mat(int row, int col, double *A, int lda)
79 {
80 int i, j;
81 for(j=0; j<col; j++)
82 {
83 for(i=0; i<row; i++)
84 {
85 printf("%9.5f ", A[i+lda*j]);
86 }
87 printf("\n");
88 }
89 printf("\n");
90 }
91/* prints a matrix in column-major format (exponential notation) */
92void d_print_e_mat(int m, int n, double *A, int lda)
93 {
94 int i, j;
95 for(i=0; i<m; i++)
96 {
97 for(j=0; j<n; j++)
98 {
99 printf("%e\t", A[i+lda*j]);
100 }
101 printf("\n");
102 }
103 printf("\n");
104 }
105/* prints the transposed of a matrix in column-major format (exponential notation) */
106void d_print_e_tran_mat(int row, int col, double *A, int lda)
107 {
108 int i, j;
109 for(j=0; j<col; j++)
110 {
111 for(i=0; i<row; i++)
112 {
113 printf("%e\t", A[i+lda*j]);
114 }
115 printf("\n");
116 }
117 printf("\n");
118 }
119/* creates a zero matrix aligned */
120void int_zeros(int **pA, int row, int col)
121 {
122 void *temp = malloc((row*col)*sizeof(int));
123 *pA = temp;
124 int *A = *pA;
125 int i;
126 for(i=0; i<row*col; i++) A[i] = 0;
127 }
128/* frees matrix */
129void int_free(int *pA)
130 {
131 free( pA );
132 }
133/* prints a matrix in column-major format */
134void int_print_mat(int row, int col, int *A, int lda)
135 {
136 int i, j;
137 for(i=0; i<row; i++)
138 {
139 for(j=0; j<col; j++)
140 {
141 printf("%d ", A[i+lda*j]);
142 }
143 printf("\n");
144 }
145 printf("\n");
146 }
147#endif
148
149
150
151#define KEEP_X0 0
152
153// printing
154#define PRINT 1
155
156/************************************************
157Mass-spring system: nx/2 masses connected each other with springs (in a row), and the first and the last one to walls. nu (<=nx) controls act on the first nu masses. The system is sampled with sampling time Ts.
158************************************************/
159void mass_spring_system(double Ts, int nx, int nu, int N, double *A, double *B, double *b, double *x0)
160 {
161
162 int nx2 = nx*nx;
163
164 int info = 0;
165
166 int pp = nx/2; // number of masses
167
168/************************************************
169* build the continuous time system
170************************************************/
171
172 double *T; d_zeros(&T, pp, pp);
173 int ii;
174 for(ii=0; ii<pp; ii++) T[ii*(pp+1)] = -2;
175 for(ii=0; ii<pp-1; ii++) T[ii*(pp+1)+1] = 1;
176 for(ii=1; ii<pp; ii++) T[ii*(pp+1)-1] = 1;
177
178 double *Z; d_zeros(&Z, pp, pp);
179 double *I; d_zeros(&I, pp, pp); for(ii=0; ii<pp; ii++) I[ii*(pp+1)]=1.0; // = eye(pp);
180 double *Ac; d_zeros(&Ac, nx, nx);
181 dmcopy(pp, pp, Z, pp, Ac, nx);
182 dmcopy(pp, pp, T, pp, Ac+pp, nx);
183 dmcopy(pp, pp, I, pp, Ac+pp*nx, nx);
184 dmcopy(pp, pp, Z, pp, Ac+pp*(nx+1), nx);
185 free(T);
186 free(Z);
187 free(I);
188
189 d_zeros(&I, nu, nu); for(ii=0; ii<nu; ii++) I[ii*(nu+1)]=1.0; //I = eye(nu);
190 double *Bc; d_zeros(&Bc, nx, nu);
191 dmcopy(nu, nu, I, nu, Bc+pp, nx);
192 free(I);
193
194/************************************************
195* compute the discrete time system
196************************************************/
197
198 double *bb; d_zeros(&bb, nx, 1);
199 dmcopy(nx, 1, bb, nx, b, nx);
200
201 dmcopy(nx, nx, Ac, nx, A, nx);
202 dscal_3l(nx2, Ts, A);
203 expm(nx, A);
204
205 d_zeros(&T, nx, nx);
206 d_zeros(&I, nx, nx); for(ii=0; ii<nx; ii++) I[ii*(nx+1)]=1.0; //I = eye(nx);
207 dmcopy(nx, nx, A, nx, T, nx);
208 daxpy_3l(nx2, -1.0, I, T);
209 dgemm_nn_3l(nx, nu, nx, T, nx, Bc, nx, B, nx);
210 free(T);
211 free(I);
212
213 int *ipiv = (int *) malloc(nx*sizeof(int));
214 dgesv_3l(nx, nu, Ac, nx, ipiv, B, nx, &info);
215 free(ipiv);
216
217 free(Ac);
218 free(Bc);
219 free(bb);
220
221
222/************************************************
223* initial state
224************************************************/
225
226 if(nx==4)
227 {
228 x0[0] = 5;
229 x0[1] = 10;
230 x0[2] = 15;
231 x0[3] = 20;
232 }
233 else
234 {
235 int jj;
236 for(jj=0; jj<nx; jj++)
237 x0[jj] = 1;
238 }
239
240 }
241
242
243
244int main()
245 {
246
247
248 // local variables
249
250 int ii, jj;
251
252 int rep, nrep=1000;
253
254 struct timeval tv0, tv1;
255
256
257
258 // problem size
259
260 int nx_ = 8; // number of states (it has to be even for the mass-spring system test problem)
261 int nu_ = 3; // number of inputs (controllers) (it has to be at least 1 and at most nx/2 for the mass-spring system test problem)
262 int N = 5; // horizon lenght
263
264
265
266 // stage-wise variant size
267
268 int nx[N+1];
269#if KEEP_X0
270 nx[0] = nx_;
271#else
272 nx[0] = 0;
273#endif
274 for(ii=1; ii<=N; ii++)
275 nx[ii] = nx_;
276// nx[N] = 0;
277
278 int nu[N+1];
279 for(ii=0; ii<N; ii++)
280 nu[ii] = nu_;
281 nu[N] = 0;
282
283#if 1
284 int nb[N+1];
285#if KEEP_X0
286 nb[0] = nu[0]+nx[0]/2;
287#else
288 nb[0] = nu[0];
289#endif
290 for(ii=1; ii<N; ii++)
291 nb[ii] = nu[1]+nx[1]/2;
292 nb[N] = nx[N]/2;
293
294 int ng[N+1];
295 ng[0] = 0;
296 for(ii=1; ii<N; ii++)
297 ng[ii] = 0;
298 ng[N] = 0;
299#elif 0
300 int nb[N+1];
301 nb[0] = 0;
302 for(ii=1; ii<N; ii++)
303 nb[ii] = 0;
304 nb[N] = 0;
305
306 int ng[N+1];
307#if KEEP_X0
308 ng[0] = nu[0]+nx[0]/2;
309#else
310 ng[0] = nu[0];
311#endif
312 for(ii=1; ii<N; ii++)
313 ng[ii] = nu[1]+nx[1]/2;
314 ng[N] = nx[N]/2;
315#else
316 int nb[N+1];
317 nb[0] = nu[0] + nx[0]/2;
318 for(ii=1; ii<N; ii++)
319 nb[ii] = nu[ii] + nx[ii]/2;
320 nb[N] = nu[N] + nx[N]/2;
321
322 int ng[N+1];
323#if KEEP_X0
324 ng[0] = nx[0]/2;
325#else
326 ng[0] = 0;
327#endif
328 for(ii=1; ii<N; ii++)
329 ng[ii] = nx[1]/2;
330 ng[N] = nx[N]/2;
331#endif
332
333/************************************************
334* dynamical system
335************************************************/
336
337 double *A; d_zeros(&A, nx_, nx_); // states update matrix
338
339 double *B; d_zeros(&B, nx_, nu_); // inputs matrix
340
341 double *b; d_zeros(&b, nx_, 1); // states offset
342 double *x0; d_zeros(&x0, nx_, 1); // initial state
343
344 double Ts = 0.5; // sampling time
345 mass_spring_system(Ts, nx_, nu_, N, A, B, b, x0);
346
347 for(jj=0; jj<nx_; jj++)
348 b[jj] = 0.1;
349
350 for(jj=0; jj<nx_; jj++)
351 x0[jj] = 0;
352 x0[0] = 2.5;
353 x0[1] = 2.5;
354
355 double *b0; d_zeros(&b0, nx_, 1);
356 dgemv_n_3l(nx_, nx_, A, nx_, x0, b0);
357 daxpy_3l(nx_, 1.0, b, b0);
358
359#if PRINT
360 d_print_mat(nx_, nx_, A, nx_);
361 d_print_mat(nx_, nu_, B, nu_);
362 d_print_mat(1, nx_, b, 1);
363 d_print_mat(1, nx_, x0, 1);
364 d_print_mat(1, nx_, b0, 1);
365#endif
366
367/************************************************
368* cost function
369************************************************/
370
371 double *Q; d_zeros(&Q, nx_, nx_);
372 for(ii=0; ii<nx_; ii++) Q[ii*(nx_+1)] = 1.0;
373
374 double *R; d_zeros(&R, nu_, nu_);
375 for(ii=0; ii<nu_; ii++) R[ii*(nu_+1)] = 2.0;
376
377 double *S; d_zeros(&S, nu_, nx_);
378
379 double *q; d_zeros(&q, nx_, 1);
380 for(ii=0; ii<nx_; ii++) q[ii] = 0.1;
381
382 double *r; d_zeros(&r, nu_, 1);
383 for(ii=0; ii<nu_; ii++) r[ii] = 0.2;
384
385 double *r0; d_zeros(&r0, nu_, 1);
386 dgemv_n_3l(nu_, nx_, S, nu_, x0, r0);
387 daxpy_3l(nu_, 1.0, r, r0);
388
389#if PRINT
390 d_print_mat(nx_, nx_, Q, nx_);
391 d_print_mat(nu_, nu_, R, nu_);
392 d_print_mat(nu_, nx_, S, nu_);
393 d_print_mat(1, nx_, q, 1);
394 d_print_mat(1, nu_, r, 1);
395 d_print_mat(1, nu_, r0, 1);
396#endif
397
398 // maximum element in cost functions
399 double mu0 = 2.0;
400
401/************************************************
402* box & general constraints
403************************************************/
404
405 int *idxb0; int_zeros(&idxb0, nb[0], 1);
406 double *d_lb0; d_zeros(&d_lb0, nb[0], 1);
407 double *d_ub0; d_zeros(&d_ub0, nb[0], 1);
408 double *d_lg0; d_zeros(&d_lg0, ng[0], 1);
409 double *d_ug0; d_zeros(&d_ug0, ng[0], 1);
410 for(ii=0; ii<nb[0]; ii++)
411 {
412 if(ii<nu[0]) // input
413 {
414 d_lb0[ii] = - 0.5; // umin
415 d_ub0[ii] = 0.5; // umax
416 }
417 else // state
418 {
419 d_lb0[ii] = - 4.0; // xmin
420 d_ub0[ii] = 4.0; // xmax
421 }
422 idxb0[ii] = ii;
423 }
424 for(ii=0; ii<ng[0]; ii++)
425 {
426 if(ii<nu[0]-nb[0]) // input
427 {
428 d_lg0[ii] = - 0.5; // umin
429 d_ug0[ii] = 0.5; // umax
430 }
431 else // state
432 {
433 d_lg0[ii] = - 4.0; // xmin
434 d_ug0[ii] = 4.0; // xmax
435 }
436 }
437
438 int *idxb1; int_zeros(&idxb1, nb[1], 1);
439 double *d_lb1; d_zeros(&d_lb1, nb[1], 1);
440 double *d_ub1; d_zeros(&d_ub1, nb[1], 1);
441 double *d_lg1; d_zeros(&d_lg1, ng[1], 1);
442 double *d_ug1; d_zeros(&d_ug1, ng[1], 1);
443 for(ii=0; ii<nb[1]; ii++)
444 {
445 if(ii<nu[1]) // input
446 {
447 d_lb1[ii] = - 0.5; // umin
448 d_ub1[ii] = 0.5; // umax
449 }
450 else // state
451 {
452 d_lb1[ii] = - 4.0; // xmin
453 d_ub1[ii] = 4.0; // xmax
454 }
455 idxb1[ii] = ii;
456 }
457 for(ii=0; ii<ng[1]; ii++)
458 {
459 if(ii<nu[1]-nb[1]) // input
460 {
461 d_lg1[ii] = - 0.5; // umin
462 d_ug1[ii] = 0.5; // umax
463 }
464 else // state
465 {
466 d_lg1[ii] = - 4.0; // xmin
467 d_ug1[ii] = 4.0; // xmax
468 }
469 }
470
471
472 int *idxbN; int_zeros(&idxbN, nb[N], 1);
473 double *d_lbN; d_zeros(&d_lbN, nb[N], 1);
474 double *d_ubN; d_zeros(&d_ubN, nb[N], 1);
475 double *d_lgN; d_zeros(&d_lgN, ng[N], 1);
476 double *d_ugN; d_zeros(&d_ugN, ng[N], 1);
477 for(ii=0; ii<nb[N]; ii++)
478 {
479 d_lbN[ii] = - 4.0; // xmin
480 d_ubN[ii] = 4.0; // xmax
481 idxbN[ii] = ii;
482 }
483 for(ii=0; ii<ng[N]; ii++)
484 {
485 d_lgN[ii] = - 4.0; // dmin
486 d_ugN[ii] = 4.0; // dmax
487 }
488
489 double *C0; d_zeros(&C0, ng[0], nx[0]);
490 double *D0; d_zeros(&D0, ng[0], nu[0]);
491 for(ii=0; ii<nu[0]-nb[0] & ii<ng[0]; ii++)
492 D0[ii+(nb[0]+ii)*ng[0]] = 1.0;
493 for(; ii<ng[0]; ii++)
494 C0[ii+(nb[0]+ii-nu[0])*ng[0]] = 1.0;
495
496 double *C1; d_zeros(&C1, ng[1], nx[1]);
497 double *D1; d_zeros(&D1, ng[1], nu[1]);
498 for(ii=0; ii<nu[1]-nb[1] & ii<ng[1]; ii++)
499 D1[ii+(nb[1]+ii)*ng[1]] = 1.0;
500 for(; ii<ng[1]; ii++)
501 C1[ii+(nb[1]+ii-nu[1])*ng[1]] = 1.0;
502
503 double *CN; d_zeros(&CN, ng[N], nx[N]);
504 double *DN; d_zeros(&DN, ng[N], nu[N]);
505 for(ii=0; ii<nu[N]-nb[N] & ii<ng[N]; ii++)
506 DN[ii+(nb[N]+ii)*ng[N]] = 1.0;
507 for(; ii<ng[N]; ii++)
508 CN[ii+(nb[N]+ii-nu[N])*ng[N]] = 1.0;
509
510#if PRINT
511 // box constraints
512 int_print_mat(1, nb[0], idxb0, 1);
513 d_print_mat(1, nb[0], d_lb0, 1);
514 d_print_mat(1, nb[0], d_ub0, 1);
515 int_print_mat(1, nb[1], idxb1, 1);
516 d_print_mat(1, nb[1], d_lb1, 1);
517 d_print_mat(1, nb[1], d_ub1, 1);
518 int_print_mat(1, nb[N], idxbN, 1);
519 d_print_mat(1, nb[N], d_lbN, 1);
520 d_print_mat(1, nb[N], d_ubN, 1);
521 // general constraints
522 d_print_mat(1, ng[0], d_lg0, 1);
523 d_print_mat(1, ng[0], d_ug0, 1);
524 d_print_mat(ng[0], nu[0], D0, ng[0]);
525 d_print_mat(ng[0], nx[0], C0, ng[0]);
526 d_print_mat(1, ng[1], d_lg1, 1);
527 d_print_mat(1, ng[1], d_ug1, 1);
528 d_print_mat(ng[1], nu[1], D1, ng[1]);
529 d_print_mat(ng[1], nx[1], C1, ng[1]);
530 d_print_mat(1, ng[N], d_lgN, 1);
531 d_print_mat(1, ng[N], d_ugN, 1);
532 d_print_mat(ng[N], nu[N], DN, ng[N]);
533 d_print_mat(ng[N], nx[N], CN, ng[N]);
534#endif
535
536/************************************************
537* array of matrices
538************************************************/
539
540 double *hA[N];
541 double *hB[N];
542 double *hb[N];
543 double *hQ[N+1];
544 double *hS[N+1];
545 double *hR[N+1];
546 double *hq[N+1];
547 double *hr[N+1];
548 double *hd_lb[N+1];
549 double *hd_ub[N+1];
550 double *hd_lg[N+1];
551 double *hd_ug[N+1];
552 double *hC[N+1];
553 double *hD[N+1];
554 int *hidxb[N+1];
555
556 hA[0] = A;
557 hB[0] = B;
558 hb[0] = b0;
559 hQ[0] = Q;
560 hS[0] = S;
561 hR[0] = R;
562 hq[0] = q;
563 hr[0] = r0;
564 hidxb[0] = idxb0;
565 hd_lb[0] = d_lb0;
566 hd_ub[0] = d_ub0;
567 hd_lg[0] = d_lg0;
568 hd_ug[0] = d_ug0;
569 hC[0] = C0;
570 hD[0] = D0;
571 for(ii=1; ii<N; ii++)
572 {
573 hA[ii] = A;
574 hB[ii] = B;
575 hb[ii] = b;
576 hQ[ii] = Q;
577 hS[ii] = S;
578 hR[ii] = R;
579 hq[ii] = q;
580 hr[ii] = r;
581 hidxb[ii] = idxb1;
582 hd_lb[ii] = d_lb1;
583 hd_ub[ii] = d_ub1;
584 hd_lg[ii] = d_lg1;
585 hd_ug[ii] = d_ug1;
586 hC[ii] = C1;
587 hD[ii] = D1;
588 }
589 hQ[N] = Q;
590 hS[N] = S;
591 hR[N] = R;
592 hq[N] = q;
593 hr[N] = r;
594 hidxb[N] = idxbN;
595 hd_lb[N] = d_lbN;
596 hd_ub[N] = d_ubN;
597 hd_lg[N] = d_lgN;
598 hd_ug[N] = d_ugN;
599 hC[N] = CN;
600 hD[N] = DN;
601
602/************************************************
603* ocp qp
604************************************************/
605
606 int qp_size = d_memsize_ocp_qp(N, nx, nu, nb, ng);
607 printf("\nqp size = %d\n", qp_size);
608 void *qp_mem = malloc(qp_size);
609
610 struct d_ocp_qp qp;
611 d_create_ocp_qp(N, nx, nu, nb, ng, &qp, qp_mem);
612 d_cvt_colmaj_to_ocp_qp(hA, hB, hb, hQ, hS, hR, hq, hr, hidxb, hd_lb, hd_ub, hC, hD, hd_lg, hd_ug, &qp);
613#if 0
614 printf("\nN = %d\n", qp.N);
615 for(ii=0; ii<N; ii++)
616 d_print_strmat(qp.nu[ii]+qp.nx[ii]+1, qp.nx[ii+1], qp.BAbt+ii, 0, 0);
617 for(ii=0; ii<N; ii++)
618 d_print_tran_strvec(qp.nx[ii+1], qp.b+ii, 0);
619 for(ii=0; ii<=N; ii++)
620 d_print_strmat(qp.nu[ii]+qp.nx[ii]+1, qp.nu[ii]+qp.nx[ii], qp.RSQrq+ii, 0, 0);
621 for(ii=0; ii<=N; ii++)
622 d_print_tran_strvec(qp.nu[ii]+qp.nx[ii], qp.rq+ii, 0);
623 for(ii=0; ii<=N; ii++)
624 int_print_mat(1, nb[ii], qp.idxb[ii], 1);
625 for(ii=0; ii<=N; ii++)
626 d_print_tran_strvec(qp.nb[ii], qp.d_lb+ii, 0);
627 for(ii=0; ii<=N; ii++)
628 d_print_tran_strvec(qp.nb[ii], qp.d_ub+ii, 0);
629 for(ii=0; ii<=N; ii++)
630 d_print_strmat(qp.nu[ii]+qp.nx[ii], qp.ng[ii], qp.DCt+ii, 0, 0);
631 for(ii=0; ii<=N; ii++)
632 d_print_tran_strvec(qp.ng[ii], qp.d_lg+ii, 0);
633 for(ii=0; ii<=N; ii++)
634 d_print_tran_strvec(qp.ng[ii], qp.d_ug+ii, 0);
635 return;
636#endif
637
638/************************************************
639* ocp qp sol
640************************************************/
641
642 int qp_sol_size = d_memsize_ocp_qp_sol(N, nx, nu, nb, ng);
643 printf("\nqp sol size = %d\n", qp_sol_size);
644 void *qp_sol_mem = malloc(qp_sol_size);
645
646 struct d_ocp_qp_sol qp_sol;
647 d_create_ocp_qp_sol(N, nx, nu, nb, ng, &qp_sol, qp_sol_mem);
648
649/************************************************
650* ipm
651************************************************/
652
653 struct d_ipm_hard_ocp_qp_arg arg;
654 arg.alpha_min = 1e-8;
655 arg.mu_max = 1e-12;
656 arg.iter_max = 20;
657 arg.mu0 = 2.0;
658
659 int ipm_size = d_memsize_ipm_hard_ocp_qp(&qp, &arg);
660 printf("\nipm size = %d\n", ipm_size);
661 void *ipm_mem = malloc(ipm_size);
662
663 struct d_ipm_hard_ocp_qp_workspace workspace;
664 d_create_ipm_hard_ocp_qp(&qp, &arg, &workspace, ipm_mem);
665
666 gettimeofday(&tv0, NULL); // start
667
668 for(rep=0; rep<nrep; rep++)
669 {
670// d_solve_ipm_hard_ocp_qp(&qp, &qp_sol, &workspace);
671 d_solve_ipm2_hard_ocp_qp(&qp, &qp_sol, &workspace);
672 }
673
674 gettimeofday(&tv1, NULL); // stop
675
676 double time_ocp_ipm = (tv1.tv_sec-tv0.tv_sec)/(nrep+0.0)+(tv1.tv_usec-tv0.tv_usec)/(nrep*1e6);
677
678/************************************************
679* extract and print solution
680************************************************/
681
682 double *u[N+1]; for(ii=0; ii<=N; ii++) d_zeros(u+ii, nu[ii], 1);
683 double *x[N+1]; for(ii=0; ii<=N; ii++) d_zeros(x+ii, nx[ii], 1);
684 double *pi[N]; for(ii=0; ii<N; ii++) d_zeros(pi+ii, nx[ii+1], 1);
685 double *lam_lb[N+1]; for(ii=0; ii<=N; ii++) d_zeros(lam_lb+ii, nb[ii], 1);
686 double *lam_ub[N+1]; for(ii=0; ii<=N; ii++) d_zeros(lam_ub+ii, nb[ii], 1);
687 double *lam_lg[N+1]; for(ii=0; ii<=N; ii++) d_zeros(lam_lg+ii, ng[ii], 1);
688 double *lam_ug[N+1]; for(ii=0; ii<=N; ii++) d_zeros(lam_ug+ii, ng[ii], 1);
689
690 d_cvt_ocp_qp_sol_to_colmaj(&qp, &qp_sol, u, x, pi, lam_lb, lam_ub, lam_lg, lam_ug);
691
692#if 1
693 printf("\nsolution\n\n");
694 printf("\nu\n");
695 for(ii=0; ii<=N; ii++)
696 d_print_mat(1, nu[ii], u[ii], 1);
697 printf("\nx\n");
698 for(ii=0; ii<=N; ii++)
699 d_print_mat(1, nx[ii], x[ii], 1);
700 printf("\npi\n");
701 for(ii=0; ii<N; ii++)
702 d_print_mat(1, nx[ii+1], pi[ii], 1);
703 printf("\nlam_lb\n");
704 for(ii=0; ii<=N; ii++)
705 d_print_mat(1, nb[ii], lam_lb[ii], 1);
706 printf("\nlam_ub\n");
707 for(ii=0; ii<=N; ii++)
708 d_print_mat(1, nb[ii], lam_ub[ii], 1);
709 printf("\nlam_lg\n");
710 for(ii=0; ii<=N; ii++)
711 d_print_mat(1, ng[ii], lam_lg[ii], 1);
712 printf("\nlam_ug\n");
713 for(ii=0; ii<=N; ii++)
714 d_print_mat(1, ng[ii], lam_ug[ii], 1);
715
716 printf("\nt_lb\n");
717 for(ii=0; ii<=N; ii++)
718 d_print_mat(1, nb[ii], (qp_sol.t_lb+ii)->pa, 1);
719 printf("\nt_ub\n");
720 for(ii=0; ii<=N; ii++)
721 d_print_mat(1, nb[ii], (qp_sol.t_ub+ii)->pa, 1);
722 printf("\nt_lg\n");
723 for(ii=0; ii<=N; ii++)
724 d_print_mat(1, ng[ii], (qp_sol.t_lg+ii)->pa, 1);
725 printf("\nt_ug\n");
726 for(ii=0; ii<=N; ii++)
727 d_print_mat(1, ng[ii], (qp_sol.t_ug+ii)->pa, 1);
728
729 printf("\nresiduals\n\n");
730 printf("\nres_g\n");
731 for(ii=0; ii<=N; ii++)
732 d_print_e_mat(1, nu[ii]+nx[ii], (workspace.res_g+ii)->pa, 1);
733 printf("\nres_b\n");
734 for(ii=0; ii<N; ii++)
735 d_print_e_mat(1, nx[ii+1], (workspace.res_b+ii)->pa, 1);
736 printf("\nres_m_lb\n");
737 for(ii=0; ii<=N; ii++)
738 d_print_e_mat(1, nb[ii], (workspace.res_m_lb+ii)->pa, 1);
739 printf("\nres_m_ub\n");
740 for(ii=0; ii<=N; ii++)
741 d_print_e_mat(1, nb[ii], (workspace.res_m_ub+ii)->pa, 1);
742 printf("\nres_m_lg\n");
743 for(ii=0; ii<=N; ii++)
744 d_print_e_mat(1, ng[ii], (workspace.res_m_lg+ii)->pa, 1);
745 printf("\nres_m_ug\n");
746 for(ii=0; ii<=N; ii++)
747 d_print_e_mat(1, ng[ii], (workspace.res_m_ug+ii)->pa, 1);
748 printf("\nres_d_lb\n");
749 for(ii=0; ii<=N; ii++)
750 d_print_e_mat(1, nb[ii], (workspace.res_d_lb+ii)->pa, 1);
751 printf("\nres_d_ub\n");
752 for(ii=0; ii<=N; ii++)
753 d_print_e_mat(1, nb[ii], (workspace.res_d_ub+ii)->pa, 1);
754 printf("\nres_d_lg\n");
755 for(ii=0; ii<=N; ii++)
756 d_print_e_mat(1, ng[ii], (workspace.res_d_lg+ii)->pa, 1);
757 printf("\nres_d_ug\n");
758 for(ii=0; ii<=N; ii++)
759 d_print_e_mat(1, ng[ii], (workspace.res_d_ug+ii)->pa, 1);
760 printf("\nres_mu\n");
761 printf("\n%e\n\n", workspace.res_mu);
762#endif
763
764 printf("\nipm iter = %d\n", workspace.iter);
765 printf("\nalpha_aff\tmu_aff\t\tsigma\t\talpha\t\tmu\n");
766 d_print_e_tran_mat(5, workspace.iter, workspace.stat, 5);
767
768 printf("\nocp ipm time = %e [s]\n\n", time_ocp_ipm);
769
770/************************************************
771* free memory
772************************************************/
773
774 d_free(A);
775 d_free(B);
776 d_free(b);
777 d_free(x0);
778 d_free(Q);
779 d_free(R);
780 d_free(S);
781 d_free(q);
782 d_free(r);
783 d_free(r0);
784 int_free(idxb0);
785 d_free(d_lb0);
786 d_free(d_ub0);
787 int_free(idxb1);
788 d_free(d_lb1);
789 d_free(d_ub1);
790 int_free(idxbN);
791 d_free(d_lbN);
792 d_free(d_ubN);
793 d_free(C0);
794 d_free(D0);
795 d_free(d_lg0);
796 d_free(d_ug0);
797 d_free(C1);
798 d_free(D1);
799 d_free(d_lg1);
800 d_free(d_ug1);
801 d_free(CN);
802 d_free(DN);
803 d_free(d_lgN);
804 d_free(d_ugN);
805
806 for(ii=0; ii<N; ii++)
807 {
808 d_free(u[ii]);
809 d_free(x[ii]);
810 d_free(pi[ii]);
811 d_free(lam_lb[ii]);
812 d_free(lam_ub[ii]);
813 d_free(lam_lg[ii]);
814 d_free(lam_ug[ii]);
815 }
816 d_free(u[ii]);
817 d_free(x[ii]);
818 d_free(lam_lb[ii]);
819 d_free(lam_ub[ii]);
820 d_free(lam_lg[ii]);
821 d_free(lam_ug[ii]);
822
823 free(qp_mem);
824 free(qp_sol_mem);
825 free(ipm_mem);
826
827/************************************************
828* return
829************************************************/
830
831 return 0;
832
833 }