Austin Schuh | 9a24b37 | 2018-01-28 16:12:29 -0800 | [diff] [blame^] | 1 | /************************************************************************************************** |
| 2 | * * |
| 3 | * This file is part of BLASFEO. * |
| 4 | * * |
| 5 | * BLASFEO -- BLAS For Embedded Optimization. * |
| 6 | * Copyright (C) 2016-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, giaf (at) dtu.dk * |
| 25 | * gianluca.frison (at) imtek.uni-freiburg.de * |
| 26 | * * |
| 27 | **************************************************************************************************/ |
| 28 | |
| 29 | #include <stdlib.h> |
| 30 | #include <stdio.h> |
| 31 | #include <sys/time.h> |
| 32 | |
| 33 | #include "tools.h" |
| 34 | |
| 35 | #include "../include/blasfeo_common.h" |
| 36 | #include "../include/blasfeo_i_aux_ext_dep.h" |
| 37 | #include "../include/blasfeo_d_aux_ext_dep.h" |
| 38 | #include "../include/blasfeo_d_aux.h" |
| 39 | #include "../include/blasfeo_d_kernel.h" |
| 40 | #include "../include/blasfeo_d_blas.h" |
| 41 | |
| 42 | |
| 43 | |
| 44 | static void d_back_ric_sv_libstr(int N, int *nx, int *nu, struct d_strmat *hsBAbt, struct d_strmat *hsRSQrq, struct d_strmat *hsL, struct d_strvec *hsux, struct d_strvec *hspi, struct d_strmat *hswork_mat, struct d_strvec *hswork_vec) |
| 45 | { |
| 46 | |
| 47 | int nn; |
| 48 | |
| 49 | // factorization and backward substitution |
| 50 | |
| 51 | // last stage |
| 52 | dpotrf_l_libstr(nx[N]+1, nx[N], &hsRSQrq[N], 0, 0, &hsL[N], 0, 0); |
| 53 | |
| 54 | // middle stages |
| 55 | for(nn=0; nn<N; nn++) |
| 56 | { |
| 57 | dtrmm_rlnn_libstr(nu[N-nn-1]+nx[N-nn-1]+1, nx[N-nn], 1.0, &hsL[N-nn], nu[N-nn], nu[N-nn], &hsBAbt[N-nn-1], 0, 0, &hswork_mat[0], 0, 0); |
| 58 | dgead_libstr(1, nx[N-nn], 1.0, &hsL[N-nn], nu[N-nn]+nx[N-nn], nu[N-nn], &hswork_mat[0], nu[N-nn-1]+nx[N-nn-1], 0); |
| 59 | #if 1 |
| 60 | dsyrk_dpotrf_ln_libstr(nu[N-nn-1]+nx[N-nn-1]+1, nu[N-nn-1]+nx[N-nn-1], nx[N-nn], &hswork_mat[0], 0, 0, &hswork_mat[0], 0, 0, &hsRSQrq[N-nn-1], 0, 0, &hsL[N-nn-1], 0, 0); |
| 61 | #else |
| 62 | dsyrk_ln_libstr(nu[N-nn-1]+nx[N-nn-1]+1, nu[N-nn-1]+nx[N-nn-1], nx[N-nn], 1.0, &hswork_mat[0], 0, 0, &hswork_mat[0], 0, 0, 1.0, &hsRSQrq[N-nn-1], 0, 0, &hsL[N-nn-1], 0, 0); |
| 63 | dpotrf_l_libstr(nu[N-nn-1]+nx[N-nn-1]+1, nu[N-nn-1]+nx[N-nn-1], &hsL[N-nn-1], 0, 0, &hsL[N-nn-1], 0, 0); |
| 64 | #endif |
| 65 | } |
| 66 | |
| 67 | // forward substitution |
| 68 | |
| 69 | // first stage |
| 70 | nn = 0; |
| 71 | drowex_libstr(nu[nn]+nx[nn], -1.0, &hsL[nn], nu[nn]+nx[nn], 0, &hsux[nn], 0); |
| 72 | dtrsv_ltn_libstr(nu[nn]+nx[nn], nu[nn]+nx[nn], &hsL[nn], 0, 0, &hsux[nn], 0, &hsux[nn], 0); |
| 73 | drowex_libstr(nx[nn+1], 1.0, &hsBAbt[nn], nu[nn]+nx[nn], 0, &hsux[nn+1], nu[nn+1]); |
| 74 | dgemv_t_libstr(nu[nn]+nx[nn], nx[nn+1], 1.0, &hsBAbt[nn], 0, 0, &hsux[nn], 0, 1.0, &hsux[nn+1], nu[nn+1], &hsux[nn+1], nu[nn+1]); |
| 75 | dveccp_libstr(nx[nn+1], 1.0, &hsux[nn+1], nu[nn+1], &hspi[nn], 0); |
| 76 | drowex_libstr(nx[nn+1], 1.0, &hsL[nn+1], nu[nn+1]+nx[nn+1], nu[nn+1], &hswork_vec[0], 0); |
| 77 | dtrmv_ltn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hspi[nn], 0, &hspi[nn], 0); |
| 78 | daxpy_libstr(nx[nn+1], 1.0, &hswork_vec[0], 0, &hspi[nn], 0); |
| 79 | dtrmv_lnn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hspi[nn], 0, &hspi[nn], 0); |
| 80 | |
| 81 | // middle stages |
| 82 | for(nn=1; nn<N; nn++) |
| 83 | { |
| 84 | drowex_libstr(nu[nn], -1.0, &hsL[nn], nu[nn]+nx[nn], 0, &hsux[nn], 0); |
| 85 | dtrsv_ltn_libstr(nu[nn]+nx[nn], nu[nn], &hsL[nn], 0, 0, &hsux[nn], 0, &hsux[nn], 0); |
| 86 | drowex_libstr(nx[nn+1], 1.0, &hsBAbt[nn], nu[nn]+nx[nn], 0, &hsux[nn+1], nu[nn+1]); |
| 87 | dgemv_t_libstr(nu[nn]+nx[nn], nx[nn+1], 1.0, &hsBAbt[nn], 0, 0, &hsux[nn], 0, 1.0, &hsux[nn+1], nu[nn+1], &hsux[nn+1], nu[nn+1]); |
| 88 | dveccp_libstr(nx[nn+1], 1.0, &hsux[nn+1], nu[nn+1], &hspi[nn], 0); |
| 89 | drowex_libstr(nx[nn+1], 1.0, &hsL[nn+1], nu[nn+1]+nx[nn+1], nu[nn+1], &hswork_vec[0], 0); |
| 90 | dtrmv_ltn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hspi[nn], 0, &hspi[nn], 0); |
| 91 | daxpy_libstr(nx[nn+1], 1.0, &hswork_vec[0], 0, &hspi[nn], 0); |
| 92 | dtrmv_lnn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hspi[nn], 0, &hspi[nn], 0); |
| 93 | } |
| 94 | |
| 95 | return; |
| 96 | |
| 97 | } |
| 98 | |
| 99 | |
| 100 | |
| 101 | static void d_back_ric_trf_libstr(int N, int *nx, int *nu, struct d_strmat *hsBAbt, struct d_strmat *hsRSQrq, struct d_strmat *hsL, struct d_strmat *hswork_mat) |
| 102 | { |
| 103 | |
| 104 | int nn; |
| 105 | |
| 106 | // factorization |
| 107 | |
| 108 | // last stage |
| 109 | dpotrf_l_libstr(nx[N], nx[N], &hsRSQrq[N], 0, 0, &hsL[N], 0, 0); |
| 110 | |
| 111 | // middle stages |
| 112 | for(nn=0; nn<N; nn++) |
| 113 | { |
| 114 | dtrmm_rlnn_libstr(nu[N-nn-1]+nx[N-nn-1], nx[N-nn], 1.0, &hsL[N-nn], nu[N-nn], nu[N-nn], &hsBAbt[N-nn-1], 0, 0, &hswork_mat[0], 0, 0); |
| 115 | #if 1 |
| 116 | dsyrk_dpotrf_ln_libstr(nu[N-nn-1]+nx[N-nn-1], nu[N-nn-1]+nx[N-nn-1], nx[N-nn], &hswork_mat[0], 0, 0, &hswork_mat[0], 0, 0, &hsRSQrq[N-nn-1], 0, 0, &hsL[N-nn-1], 0, 0); |
| 117 | #else |
| 118 | dsyrk_ln_libstr(nu[N-nn-1]+nx[N-nn-1], nu[N-nn-1]+nx[N-nn-1], nx[N-nn], 1.0, &hswork_mat[0], 0, 0, &hswork_mat[0], 0, 0, 1.0, &hsRSQrq[N-nn-1], 0, 0, &hsL[N-nn-1], 0, 0); |
| 119 | dpotrf_l_libstr(nu[N-nn-1]+nx[N-nn-1], nu[N-nn-1]+nx[N-nn-1], &hsL[N-nn-1], 0, 0, &hsL[N-nn-1], 0, 0); |
| 120 | #endif |
| 121 | } |
| 122 | |
| 123 | return; |
| 124 | |
| 125 | } |
| 126 | |
| 127 | |
| 128 | |
| 129 | static void d_back_ric_trs_libstr(int N, int *nx, int *nu, struct d_strmat *hsBAbt, struct d_strvec *hsb, struct d_strvec *hsrq, struct d_strmat *hsL, struct d_strvec *hsPb, struct d_strvec *hsux, struct d_strvec *hspi, struct d_strvec *hswork_vec) |
| 130 | { |
| 131 | |
| 132 | int nn; |
| 133 | |
| 134 | // backward substitution |
| 135 | |
| 136 | // last stage |
| 137 | dveccp_libstr(nu[N]+nx[N], 1.0, &hsrq[N], 0, &hsux[N], 0); |
| 138 | |
| 139 | // middle stages |
| 140 | for(nn=0; nn<N-1; nn++) |
| 141 | { |
| 142 | // compute Pb |
| 143 | dtrmv_ltn_libstr(nx[N-nn], nx[N-nn], &hsL[N-nn], nu[N-nn], nu[N-nn], &hsb[N-nn-1], 0, &hsPb[N-nn-1], 0); |
| 144 | dtrmv_lnn_libstr(nx[N-nn], nx[N-nn], &hsL[N-nn], nu[N-nn], nu[N-nn], &hsPb[N-nn-1], 0, &hsPb[N-nn-1], 0); |
| 145 | dveccp_libstr(nu[N-nn-1]+nx[N-nn-1], 1.0, &hsrq[N-nn-1], 0, &hsux[N-nn-1], 0); |
| 146 | dveccp_libstr(nx[N-nn], 1.0, &hsPb[N-nn-1], 0, &hswork_vec[0], 0); |
| 147 | daxpy_libstr(nx[N-nn], 1.0, &hsux[N-nn], nu[N-nn], &hswork_vec[0], 0); |
| 148 | dgemv_n_libstr(nu[N-nn-1]+nx[N-nn-1], nx[N-nn], 1.0, &hsBAbt[N-nn-1], 0, 0, &hswork_vec[0], 0, 1.0, &hsux[N-nn-1], 0, &hsux[N-nn-1], 0); |
| 149 | dtrsv_lnn_libstr(nu[N-nn-1]+nx[N-nn-1], nu[N-nn-1], &hsL[N-nn-1], 0, 0, &hsux[N-nn-1], 0, &hsux[N-nn-1], 0); |
| 150 | } |
| 151 | |
| 152 | // first stage |
| 153 | nn = N-1; |
| 154 | dtrmv_ltn_libstr(nx[N-nn], nx[N-nn], &hsL[N-nn], nu[N-nn], nu[N-nn], &hsb[N-nn-1], 0, &hsPb[N-nn-1], 0); |
| 155 | dtrmv_lnn_libstr(nx[N-nn], nx[N-nn], &hsL[N-nn], nu[N-nn], nu[N-nn], &hsPb[N-nn-1], 0, &hsPb[N-nn-1], 0); |
| 156 | dveccp_libstr(nu[N-nn-1]+nx[N-nn-1], 1.0, &hsrq[N-nn-1], 0, &hsux[N-nn-1], 0); |
| 157 | dveccp_libstr(nx[N-nn], 1.0, &hsPb[N-nn-1], 0, &hswork_vec[0], 0); |
| 158 | daxpy_libstr(nx[N-nn], 1.0, &hsux[N-nn], nu[N-nn], &hswork_vec[0], 0); |
| 159 | dgemv_n_libstr(nu[N-nn-1]+nx[N-nn-1], nx[N-nn], 1.0, &hsBAbt[N-nn-1], 0, 0, &hswork_vec[0], 0, 1.0, &hsux[N-nn-1], 0, &hsux[N-nn-1], 0); |
| 160 | dtrsv_lnn_libstr(nu[N-nn-1]+nx[N-nn-1], nu[N-nn-1]+nx[N-nn-1], &hsL[N-nn-1], 0, 0, &hsux[N-nn-1], 0, &hsux[N-nn-1], 0); |
| 161 | |
| 162 | // forward substitution |
| 163 | |
| 164 | // first stage |
| 165 | nn = 0; |
| 166 | dveccp_libstr(nx[nn+1], 1.0, &hsux[nn+1], nu[nn+1], &hspi[nn], 0); |
| 167 | dveccp_libstr(nu[nn]+nx[nn], -1.0, &hsux[nn], 0, &hsux[nn], 0); |
| 168 | dtrsv_ltn_libstr(nu[nn]+nx[nn], nu[nn]+nx[nn], &hsL[nn], 0, 0, &hsux[nn], 0, &hsux[nn], 0); |
| 169 | dgemv_t_libstr(nu[nn]+nx[nn], nx[nn+1], 1.0, &hsBAbt[nn], 0, 0, &hsux[nn], 0, 1.0, &hsb[nn], 0, &hsux[nn+1], nu[nn+1]); |
| 170 | dveccp_libstr(nx[nn+1], 1.0, &hsux[nn+1], nu[nn+1], &hswork_vec[0], 0); |
| 171 | dtrmv_ltn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hswork_vec[0], 0, &hswork_vec[0], 0); |
| 172 | dtrmv_lnn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hswork_vec[0], 0, &hswork_vec[0], 0); |
| 173 | daxpy_libstr(nx[nn+1], 1.0, &hswork_vec[0], 0, &hspi[nn], 0); |
| 174 | |
| 175 | // middle stages |
| 176 | for(nn=1; nn<N; nn++) |
| 177 | { |
| 178 | dveccp_libstr(nx[nn+1], 1.0, &hsux[nn+1], nu[nn+1], &hspi[nn], 0); |
| 179 | dveccp_libstr(nu[nn], -1.0, &hsux[nn], 0, &hsux[nn], 0); |
| 180 | dtrsv_ltn_libstr(nu[nn]+nx[nn], nu[nn], &hsL[nn], 0, 0, &hsux[nn], 0, &hsux[nn], 0); |
| 181 | dgemv_t_libstr(nu[nn]+nx[nn], nx[nn+1], 1.0, &hsBAbt[nn], 0, 0, &hsux[nn], 0, 1.0, &hsb[nn], 0, &hsux[nn+1], nu[nn+1]); |
| 182 | dveccp_libstr(nx[nn+1], 1.0, &hsux[nn+1], nu[nn+1], &hswork_vec[0], 0); |
| 183 | dtrmv_ltn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hswork_vec[0], 0, &hswork_vec[0], 0); |
| 184 | dtrmv_lnn_libstr(nx[nn+1], nx[nn+1], &hsL[nn+1], nu[nn+1], nu[nn+1], &hswork_vec[0], 0, &hswork_vec[0], 0); |
| 185 | daxpy_libstr(nx[nn+1], 1.0, &hswork_vec[0], 0, &hspi[nn], 0); |
| 186 | } |
| 187 | |
| 188 | return; |
| 189 | |
| 190 | } |
| 191 | |
| 192 | |
| 193 | |
| 194 | /************************************************ |
| 195 | Mass-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. |
| 196 | ************************************************/ |
| 197 | static void d_mass_spring_system(double Ts, int nx, int nu, int N, double *A, double *B, double *b, double *x0) |
| 198 | { |
| 199 | |
| 200 | int nx2 = nx*nx; |
| 201 | |
| 202 | int info = 0; |
| 203 | |
| 204 | int pp = nx/2; // number of masses |
| 205 | |
| 206 | /************************************************ |
| 207 | * build the continuous time system |
| 208 | ************************************************/ |
| 209 | |
| 210 | double *T; d_zeros(&T, pp, pp); |
| 211 | int ii; |
| 212 | for(ii=0; ii<pp; ii++) T[ii*(pp+1)] = -2; |
| 213 | for(ii=0; ii<pp-1; ii++) T[ii*(pp+1)+1] = 1; |
| 214 | for(ii=1; ii<pp; ii++) T[ii*(pp+1)-1] = 1; |
| 215 | |
| 216 | double *Z; d_zeros(&Z, pp, pp); |
| 217 | double *I; d_zeros(&I, pp, pp); for(ii=0; ii<pp; ii++) I[ii*(pp+1)]=1.0; // = eye(pp); |
| 218 | double *Ac; d_zeros(&Ac, nx, nx); |
| 219 | dmcopy(pp, pp, Z, pp, Ac, nx); |
| 220 | dmcopy(pp, pp, T, pp, Ac+pp, nx); |
| 221 | dmcopy(pp, pp, I, pp, Ac+pp*nx, nx); |
| 222 | dmcopy(pp, pp, Z, pp, Ac+pp*(nx+1), nx); |
| 223 | free(T); |
| 224 | free(Z); |
| 225 | free(I); |
| 226 | |
| 227 | d_zeros(&I, nu, nu); for(ii=0; ii<nu; ii++) I[ii*(nu+1)]=1.0; //I = eye(nu); |
| 228 | double *Bc; d_zeros(&Bc, nx, nu); |
| 229 | dmcopy(nu, nu, I, nu, Bc+pp, nx); |
| 230 | free(I); |
| 231 | |
| 232 | /************************************************ |
| 233 | * compute the discrete time system |
| 234 | ************************************************/ |
| 235 | |
| 236 | double *bb; d_zeros(&bb, nx, 1); |
| 237 | dmcopy(nx, 1, bb, nx, b, nx); |
| 238 | |
| 239 | dmcopy(nx, nx, Ac, nx, A, nx); |
| 240 | dscal_3l(nx2, Ts, A); |
| 241 | expm(nx, A); |
| 242 | |
| 243 | d_zeros(&T, nx, nx); |
| 244 | d_zeros(&I, nx, nx); for(ii=0; ii<nx; ii++) I[ii*(nx+1)]=1.0; //I = eye(nx); |
| 245 | dmcopy(nx, nx, A, nx, T, nx); |
| 246 | daxpy_3l(nx2, -1.0, I, T); |
| 247 | dgemm_nn_3l(nx, nu, nx, T, nx, Bc, nx, B, nx); |
| 248 | free(T); |
| 249 | free(I); |
| 250 | |
| 251 | int *ipiv = (int *) malloc(nx*sizeof(int)); |
| 252 | dgesv_3l(nx, nu, Ac, nx, ipiv, B, nx, &info); |
| 253 | free(ipiv); |
| 254 | |
| 255 | free(Ac); |
| 256 | free(Bc); |
| 257 | free(bb); |
| 258 | |
| 259 | |
| 260 | /************************************************ |
| 261 | * initial state |
| 262 | ************************************************/ |
| 263 | |
| 264 | if(nx==4) |
| 265 | { |
| 266 | x0[0] = 5; |
| 267 | x0[1] = 10; |
| 268 | x0[2] = 15; |
| 269 | x0[3] = 20; |
| 270 | } |
| 271 | else |
| 272 | { |
| 273 | int jj; |
| 274 | for(jj=0; jj<nx; jj++) |
| 275 | x0[jj] = 1; |
| 276 | } |
| 277 | |
| 278 | } |
| 279 | |
| 280 | |
| 281 | |
| 282 | int main() |
| 283 | { |
| 284 | |
| 285 | printf("\nExample of LU factorization and backsolve\n\n"); |
| 286 | |
| 287 | #if defined(LA_HIGH_PERFORMANCE) |
| 288 | |
| 289 | printf("\nLA provided by BLASFEO\n\n"); |
| 290 | |
| 291 | #elif defined(LA_BLAS) |
| 292 | |
| 293 | printf("\nLA provided by BLAS\n\n"); |
| 294 | |
| 295 | #elif defined(LA_REFERENCE) |
| 296 | |
| 297 | printf("\nLA provided by REFERENCE\n\n"); |
| 298 | |
| 299 | #else |
| 300 | |
| 301 | printf("\nLA provided by ???\n\n"); |
| 302 | exit(2); |
| 303 | |
| 304 | #endif |
| 305 | |
| 306 | // loop index |
| 307 | int ii; |
| 308 | |
| 309 | /************************************************ |
| 310 | * problem size |
| 311 | ************************************************/ |
| 312 | |
| 313 | // problem size |
| 314 | int N = 4; |
| 315 | int nx_ = 4; |
| 316 | int nu_ = 1; |
| 317 | |
| 318 | // stage-wise variant size |
| 319 | int nx[N+1]; |
| 320 | nx[0] = 0; |
| 321 | for(ii=1; ii<=N; ii++) |
| 322 | nx[ii] = nx_; |
| 323 | nx[N] = nx_; |
| 324 | |
| 325 | int nu[N+1]; |
| 326 | for(ii=0; ii<N; ii++) |
| 327 | nu[ii] = nu_; |
| 328 | nu[N] = 0; |
| 329 | |
| 330 | /************************************************ |
| 331 | * dynamical system |
| 332 | ************************************************/ |
| 333 | |
| 334 | double *A; d_zeros(&A, nx_, nx_); // states update matrix |
| 335 | |
| 336 | double *B; d_zeros(&B, nx_, nu_); // inputs matrix |
| 337 | |
| 338 | double *b; d_zeros(&b, nx_, 1); // states offset |
| 339 | double *x0; d_zeros(&x0, nx_, 1); // initial state |
| 340 | |
| 341 | double Ts = 0.5; // sampling time |
| 342 | d_mass_spring_system(Ts, nx_, nu_, N, A, B, b, x0); |
| 343 | |
| 344 | for(ii=0; ii<nx_; ii++) |
| 345 | b[ii] = 0.1; |
| 346 | |
| 347 | for(ii=0; ii<nx_; ii++) |
| 348 | x0[ii] = 0; |
| 349 | x0[0] = 2.5; |
| 350 | x0[1] = 2.5; |
| 351 | |
| 352 | d_print_mat(nx_, nx_, A, nx_); |
| 353 | d_print_mat(nx_, nu_, B, nx_); |
| 354 | d_print_mat(1, nx_, b, 1); |
| 355 | d_print_mat(1, nx_, x0, 1); |
| 356 | |
| 357 | /************************************************ |
| 358 | * cost function |
| 359 | ************************************************/ |
| 360 | |
| 361 | double *R; d_zeros(&R, nu_, nu_); |
| 362 | for(ii=0; ii<nu_; ii++) R[ii*(nu_+1)] = 2.0; |
| 363 | |
| 364 | double *S; d_zeros(&S, nu_, nx_); |
| 365 | |
| 366 | double *Q; d_zeros(&Q, nx_, nx_); |
| 367 | for(ii=0; ii<nx_; ii++) Q[ii*(nx_+1)] = 1.0; |
| 368 | |
| 369 | double *r; d_zeros(&r, nu_, 1); |
| 370 | for(ii=0; ii<nu_; ii++) r[ii] = 0.2; |
| 371 | |
| 372 | double *q; d_zeros(&q, nx_, 1); |
| 373 | for(ii=0; ii<nx_; ii++) q[ii] = 0.1; |
| 374 | |
| 375 | d_print_mat(nu_, nu_, R, nu_); |
| 376 | d_print_mat(nu_, nx_, S, nu_); |
| 377 | d_print_mat(nx_, nx_, Q, nx_); |
| 378 | d_print_mat(1, nu_, r, 1); |
| 379 | d_print_mat(1, nx_, q, 1); |
| 380 | |
| 381 | /************************************************ |
| 382 | * matrices as strmat |
| 383 | ************************************************/ |
| 384 | |
| 385 | struct d_strmat sA; |
| 386 | d_allocate_strmat(nx_, nx_, &sA); |
| 387 | d_cvt_mat2strmat(nx_, nx_, A, nx_, &sA, 0, 0); |
| 388 | struct d_strvec sb; |
| 389 | d_allocate_strvec(nx_, &sb); |
| 390 | d_cvt_vec2strvec(nx_, b, &sb, 0); |
| 391 | struct d_strvec sx0; |
| 392 | d_allocate_strvec(nx_, &sx0); |
| 393 | d_cvt_vec2strvec(nx_, x0, &sx0, 0); |
| 394 | struct d_strvec sb0; |
| 395 | d_allocate_strvec(nx_, &sb0); |
| 396 | double *b0; d_zeros(&b0, nx_, 1); // states offset |
| 397 | dgemv_n_libstr(nx_, nx_, 1.0, &sA, 0, 0, &sx0, 0, 1.0, &sb, 0, &sb0, 0); |
| 398 | d_print_tran_strvec(nx_, &sb0, 0); |
| 399 | |
| 400 | struct d_strmat sBbt0; |
| 401 | d_allocate_strmat(nu_+nx_+1, nx_, &sBbt0); |
| 402 | d_cvt_tran_mat2strmat(nx_, nx_, B, nx_, &sBbt0, 0, 0); |
| 403 | drowin_libstr(nx_, 1.0, &sb0, 0, &sBbt0, nu_, 0); |
| 404 | d_print_strmat(nu_+1, nx_, &sBbt0, 0, 0); |
| 405 | |
| 406 | struct d_strmat sBAbt1; |
| 407 | d_allocate_strmat(nu_+nx_+1, nx_, &sBAbt1); |
| 408 | d_cvt_tran_mat2strmat(nx_, nu_, B, nx_, &sBAbt1, 0, 0); |
| 409 | d_cvt_tran_mat2strmat(nx_, nx_, A, nx_, &sBAbt1, nu_, 0); |
| 410 | d_cvt_tran_mat2strmat(nx_, 1, b, nx_, &sBAbt1, nu_+nx_, 0); |
| 411 | d_print_strmat(nu_+nx_+1, nx_, &sBAbt1, 0, 0); |
| 412 | |
| 413 | struct d_strvec sr0; // XXX no need to update r0 since S=0 |
| 414 | d_allocate_strvec(nu_, &sr0); |
| 415 | d_cvt_vec2strvec(nu_, r, &sr0, 0); |
| 416 | |
| 417 | struct d_strmat sRr0; |
| 418 | d_allocate_strmat(nu_+1, nu_, &sRr0); |
| 419 | d_cvt_mat2strmat(nu_, nu_, R, nu_, &sRr0, 0, 0); |
| 420 | drowin_libstr(nu_, 1.0, &sr0, 0, &sRr0, nu_, 0); |
| 421 | d_print_strmat(nu_+1, nu_, &sRr0, 0, 0); |
| 422 | |
| 423 | struct d_strvec srq1; |
| 424 | d_allocate_strvec(nu_+nx_, &srq1); |
| 425 | d_cvt_vec2strvec(nu_, r, &srq1, 0); |
| 426 | d_cvt_vec2strvec(nx_, q, &srq1, nu_); |
| 427 | |
| 428 | struct d_strmat sRSQrq1; |
| 429 | d_allocate_strmat(nu_+nx_+1, nu_+nx_, &sRSQrq1); |
| 430 | d_cvt_mat2strmat(nu_, nu_, R, nu_, &sRSQrq1, 0, 0); |
| 431 | d_cvt_tran_mat2strmat(nu_, nx_, S, nu_, &sRSQrq1, nu_, 0); |
| 432 | d_cvt_mat2strmat(nx_, nx_, Q, nx_, &sRSQrq1, nu_, nu_); |
| 433 | drowin_libstr(nu_+nx_, 1.0, &srq1, 0, &sRSQrq1, nu_+nx_, 0); |
| 434 | d_print_strmat(nu_+nx_+1, nu_+nx_, &sRSQrq1, 0, 0); |
| 435 | |
| 436 | struct d_strvec sqN; |
| 437 | d_allocate_strvec(nx_, &sqN); |
| 438 | d_cvt_vec2strvec(nx_, q, &sqN, 0); |
| 439 | |
| 440 | struct d_strmat sQqN; |
| 441 | d_allocate_strmat(nx_+1, nx_, &sQqN); |
| 442 | d_cvt_mat2strmat(nx_, nx_, Q, nx_, &sQqN, 0, 0); |
| 443 | drowin_libstr(nx_, 1.0, &sqN, 0, &sQqN, nx_, 0); |
| 444 | d_print_strmat(nx_+1, nx_, &sQqN, 0, 0); |
| 445 | |
| 446 | /************************************************ |
| 447 | * array of matrices |
| 448 | ************************************************/ |
| 449 | |
| 450 | struct d_strmat hsBAbt[N]; |
| 451 | struct d_strvec hsb[N]; |
| 452 | struct d_strmat hsRSQrq[N+1]; |
| 453 | struct d_strvec hsrq[N+1]; |
| 454 | struct d_strmat hsL[N+1]; |
| 455 | struct d_strvec hsPb[N]; |
| 456 | struct d_strvec hsux[N+1]; |
| 457 | struct d_strvec hspi[N]; |
| 458 | struct d_strmat hswork_mat[1]; |
| 459 | struct d_strvec hswork_vec[1]; |
| 460 | |
| 461 | hsBAbt[0] = sBbt0; |
| 462 | hsb[0] = sb0; |
| 463 | hsRSQrq[0] = sRr0; |
| 464 | hsrq[0] = sr0; |
| 465 | d_allocate_strmat(nu_+1, nu_, &hsL[0]); |
| 466 | d_allocate_strvec(nx_, &hsPb[0]); |
| 467 | d_allocate_strvec(nx_+nu_+1, &hsux[0]); |
| 468 | d_allocate_strvec(nx_, &hspi[0]); |
| 469 | for(ii=1; ii<N; ii++) |
| 470 | { |
| 471 | hsBAbt[ii] = sBAbt1; |
| 472 | hsb[ii] = sb; |
| 473 | hsRSQrq[ii] = sRSQrq1; |
| 474 | hsrq[ii] = srq1; |
| 475 | d_allocate_strmat(nu_+nx_+1, nu_+nx_, &hsL[ii]); |
| 476 | d_allocate_strvec(nx_, &hsPb[ii]); |
| 477 | d_allocate_strvec(nx_+nu_+1, &hsux[ii]); |
| 478 | d_allocate_strvec(nx_, &hspi[ii]); |
| 479 | } |
| 480 | hsRSQrq[N] = sQqN; |
| 481 | hsrq[N] = sqN; |
| 482 | d_allocate_strmat(nx_+1, nx_, &hsL[N]); |
| 483 | d_allocate_strvec(nx_+nu_+1, &hsux[N]); |
| 484 | d_allocate_strmat(nu_+nx_+1, nx_, &hswork_mat[0]); |
| 485 | d_allocate_strvec(nx_, &hswork_vec[0]); |
| 486 | |
| 487 | // for(ii=0; ii<N; ii++) |
| 488 | // d_print_strmat(nu[ii]+nx[ii]+1, nx[ii+1], &hsBAbt[ii], 0, 0); |
| 489 | // return 0; |
| 490 | |
| 491 | /************************************************ |
| 492 | * call Riccati solver |
| 493 | ************************************************/ |
| 494 | |
| 495 | // timing |
| 496 | struct timeval tv0, tv1, tv2, tv3; |
| 497 | int nrep = 1000; |
| 498 | int rep; |
| 499 | |
| 500 | gettimeofday(&tv0, NULL); // time |
| 501 | |
| 502 | for(rep=0; rep<nrep; rep++) |
| 503 | { |
| 504 | d_back_ric_sv_libstr(N, nx, nu, hsBAbt, hsRSQrq, hsL, hsux, hspi, hswork_mat, hswork_vec); |
| 505 | } |
| 506 | |
| 507 | gettimeofday(&tv1, NULL); // time |
| 508 | |
| 509 | for(rep=0; rep<nrep; rep++) |
| 510 | { |
| 511 | d_back_ric_trf_libstr(N, nx, nu, hsBAbt, hsRSQrq, hsL, hswork_mat); |
| 512 | } |
| 513 | |
| 514 | gettimeofday(&tv2, NULL); // time |
| 515 | |
| 516 | for(rep=0; rep<nrep; rep++) |
| 517 | { |
| 518 | d_back_ric_trs_libstr(N, nx, nu, hsBAbt, hsb, hsrq, hsL, hsPb, hsux, hspi, hswork_vec); |
| 519 | } |
| 520 | |
| 521 | gettimeofday(&tv3, NULL); // time |
| 522 | |
| 523 | float time_sv = (float) (tv1.tv_sec-tv0.tv_sec)/(nrep+0.0)+(tv1.tv_usec-tv0.tv_usec)/(nrep*1e6); |
| 524 | float time_trf = (float) (tv2.tv_sec-tv1.tv_sec)/(nrep+0.0)+(tv2.tv_usec-tv1.tv_usec)/(nrep*1e6); |
| 525 | float time_trs = (float) (tv3.tv_sec-tv2.tv_sec)/(nrep+0.0)+(tv3.tv_usec-tv2.tv_usec)/(nrep*1e6); |
| 526 | |
| 527 | // print sol |
| 528 | printf("\nux = \n\n"); |
| 529 | for(ii=0; ii<=N; ii++) |
| 530 | d_print_tran_strvec(nu[ii]+nx[ii], &hsux[ii], 0); |
| 531 | |
| 532 | printf("\npi = \n\n"); |
| 533 | for(ii=0; ii<N; ii++) |
| 534 | d_print_tran_strvec(nx[ii+1], &hspi[ii], 0); |
| 535 | |
| 536 | // printf("\nL = \n\n"); |
| 537 | // for(ii=0; ii<=N; ii++) |
| 538 | // d_print_strmat(nu[ii]+nx[ii]+1, nu[ii]+nx[ii], &hsL[ii], 0, 0); |
| 539 | |
| 540 | printf("\ntime sv\t\ttime trf\t\ttime trs\n"); |
| 541 | printf("\n%e\t%e\t%e\n", time_sv, time_trf, time_trs); |
| 542 | printf("\n"); |
| 543 | |
| 544 | /************************************************ |
| 545 | * free memory |
| 546 | ************************************************/ |
| 547 | |
| 548 | d_free(A); |
| 549 | d_free(B); |
| 550 | d_free(b); |
| 551 | d_free(x0); |
| 552 | d_free(R); |
| 553 | d_free(S); |
| 554 | d_free(Q); |
| 555 | d_free(r); |
| 556 | d_free(q); |
| 557 | d_free(b0); |
| 558 | d_free_strmat(&sA); |
| 559 | d_free_strvec(&sb); |
| 560 | d_free_strmat(&sBbt0); |
| 561 | d_free_strvec(&sb0); |
| 562 | d_free_strmat(&sBAbt1); |
| 563 | d_free_strmat(&sRr0); |
| 564 | d_free_strvec(&sr0); |
| 565 | d_free_strmat(&sRSQrq1); |
| 566 | d_free_strvec(&srq1); |
| 567 | d_free_strmat(&sQqN); |
| 568 | d_free_strvec(&sqN); |
| 569 | d_free_strmat(&hsL[0]); |
| 570 | d_free_strvec(&hsPb[0]); |
| 571 | d_free_strvec(&hsux[0]); |
| 572 | d_free_strvec(&hspi[0]); |
| 573 | for(ii=1; ii<N; ii++) |
| 574 | { |
| 575 | d_free_strmat(&hsL[ii]); |
| 576 | d_free_strvec(&hsPb[ii]); |
| 577 | d_free_strvec(&hsux[ii]); |
| 578 | d_free_strvec(&hspi[ii]); |
| 579 | } |
| 580 | d_free_strmat(&hsL[N]); |
| 581 | d_free_strvec(&hsux[N]); |
| 582 | d_free_strmat(&hswork_mat[0]); |
| 583 | d_free_strvec(&hswork_vec[0]); |
| 584 | |
| 585 | |
| 586 | /************************************************ |
| 587 | * return |
| 588 | ************************************************/ |
| 589 | |
| 590 | return 0; |
| 591 | |
| 592 | } |
| 593 | |
| 594 | |
| 595 | |