blob: 1e71494328f30e311ebb4e5768e4ed09643f7b0d [file] [log] [blame]
Austin Schuh9a24b372018-01-28 16:12:29 -08001/**************************************************************************************************
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//#if defined(TARGET_X64_INTEL_SANDY_BRIDGE)
34//#include <xmmintrin.h> // needed to flush to zero sub-normals with _MM_SET_FLUSH_ZERO_MODE (_MM_FLUSH_ZERO_ON); in the main()
35//#endif
36
37#include "../include/blasfeo_common.h"
38#include "../include/blasfeo_d_aux_ext_dep.h"
39#include "../include/blasfeo_d_aux.h"
40#include "../include/blasfeo_i_aux_ext_dep.h"
41#include "../include/blasfeo_v_aux_ext_dep.h"
42#include "../include/blasfeo_d_kernel.h"
43#include "../include/blasfeo_d_blas.h"
44
45#ifndef D_PS
46#define D_PS 1
47#endif
48#ifndef D_NC
49#define D_NC 1
50#endif
51
52
53
54#if defined(REF_BLAS_OPENBLAS)
55void openblas_set_num_threads(int num_threads);
56#endif
57#if defined(REF_BLAS_BLIS)
58void omp_set_num_threads(int num_threads);
59#endif
60#if defined(REF_BLAS_MKL)
61#include "mkl.h"
62#endif
63
64
65
66#include "cpu_freq.h"
67
68
69
70int main()
71 {
72
73#if defined(REF_BLAS_OPENBLAS)
74 openblas_set_num_threads(1);
75#endif
76#if defined(REF_BLAS_BLIS)
77 omp_set_num_threads(1);
78#endif
79#if defined(REF_BLAS_MKL)
80 mkl_set_num_threads(1);
81#endif
82
83//#if defined(TARGET_X64_INTEL_SANDY_BRIDGE)
84// _MM_SET_FLUSH_ZERO_MODE(_MM_FLUSH_ZERO_ON); // flush to zero subnormals !!! works only with one thread !!!
85//#endif
86
87 printf("\n");
88 printf("\n");
89 printf("\n");
90
91 printf("BLAS performance test - double precision\n");
92 printf("\n");
93
94 // maximum frequency of the processor
95 const float GHz_max = GHZ_MAX;
96 printf("Frequency used to compute theoretical peak: %5.1f GHz (edit test_param.h to modify this value).\n", GHz_max);
97 printf("\n");
98
99 // maximum flops per cycle, double precision
100#if defined(TARGET_X64_INTEL_HASWELL)
101 const float flops_max = 16;
102 printf("Testing BLAS version for AVX2 and FMA instruction sets, 64 bit (optimized for Intel Haswell): theoretical peak %5.1f Gflops\n", flops_max*GHz_max);
103#elif defined(TARGET_X64_INTEL_SANDY_BRIDGE)
104 const float flops_max = 8;
105 printf("Testing BLAS version for AVX instruction set, 64 bit (optimized for Intel Sandy Bridge): theoretical peak %5.1f Gflops\n", flops_max*GHz_max);
106#elif defined(TARGET_X64_INTEL_CORE)
107 const float flops_max = 4;
108 printf("Testing BLAS version for SSE3 instruction set, 64 bit (optimized for Intel Core): theoretical peak %5.1f Gflops\n", flops_max*GHz_max);
109#elif defined(TARGET_X64_AMD_BULLDOZER)
110 const float flops_max = 8;
111 printf("Testing BLAS version for SSE3 and FMA instruction set, 64 bit (optimized for AMD Bulldozer): theoretical peak %5.1f Gflops\n", flops_max*GHz_max);
112#elif defined(TARGET_ARMV8A_ARM_CORTEX_A57)
113 const float flops_max = 4;
114 printf("Testing BLAS version for NEONv2 instruction set, 64 bit (optimized for ARM Cortex A57): theoretical peak %5.1f Gflops\n", flops_max*GHz_max);
115#elif defined(TARGET_ARMV7A_ARM_CORTEX_A15)
116 const float flops_max = 2;
117 printf("Testing BLAS version for VFPv4 instruction set, 32 bit (optimized for ARM Cortex A15): theoretical peak %5.1f Gflops\n", flops_max*GHz_max);
118#elif defined(TARGET_GENERIC)
119 const float flops_max = 2;
120 printf("Testing BLAS version for generic scalar instruction set: theoretical peak %5.1f Gflops ???\n", flops_max*GHz_max);
121#endif
122
123// FILE *f;
124// f = fopen("./test_problems/results/test_blas.m", "w"); // a
125
126#if defined(TARGET_X64_INTEL_HASWELL)
127// fprintf(f, "C = 'd_x64_intel_haswell';\n");
128// fprintf(f, "\n");
129#elif defined(TARGET_X64_INTEL_SANDY_BRIDGE)
130// fprintf(f, "C = 'd_x64_intel_sandybridge';\n");
131// fprintf(f, "\n");
132#elif defined(TARGET_X64_INTEL_CORE)
133// fprintf(f, "C = 'd_x64_intel_core';\n");
134// fprintf(f, "\n");
135#elif defined(TARGET_X64_AMD_BULLDOZER)
136// fprintf(f, "C = 'd_x64_amd_bulldozer';\n");
137// fprintf(f, "\n");
138#elif defined(TARGET_ARMV8A_ARM_CORTEX_A57)
139// fprintf(f, "C = 'd_armv8a_arm_cortex_a57';\n");
140// fprintf(f, "\n");
141#elif defined(TARGET_ARMV7A_ARM_CORTEX_A15)
142// fprintf(f, "C = 'd_armv7a_arm_cortex_a15';\n");
143// fprintf(f, "\n");
144#elif defined(TARGET_GENERIC)
145// fprintf(f, "C = 'd_generic';\n");
146// fprintf(f, "\n");
147#endif
148
149// fprintf(f, "A = [%f %f];\n", GHz_max, flops_max);
150// fprintf(f, "\n");
151
152// fprintf(f, "B = [\n");
153
154
155
156 int i, j, rep, ll;
157
158 const int bsd = D_PS;
159 const int ncd = D_NC;
160
161/* int info = 0;*/
162
163 printf("\nn\t dgemm_blasfeo\t dgemm_blas\n");
164 printf("\nn\t Gflops\t %%\t Gflops\n\n");
165
166#if 1
167 int nn[] = {4, 8, 12, 16, 20, 24, 28, 32, 36, 40, 44, 48, 52, 56, 60, 64, 68, 72, 76, 80, 84, 88, 92, 96, 100, 104, 108, 112, 116, 120, 124, 128, 132, 136, 140, 144, 148, 152, 156, 160, 164, 168, 172, 176, 180, 184, 188, 192, 196, 200, 204, 208, 212, 216, 220, 224, 228, 232, 236, 240, 244, 248, 252, 256, 260, 264, 268, 272, 276, 280, 284, 288, 292, 296, 300, 304, 308, 312, 316, 320, 324, 328, 332, 336, 340, 344, 348, 352, 356, 360, 364, 368, 372, 376, 380, 384, 388, 392, 396, 400, 404, 408, 412, 416, 420, 424, 428, 432, 436, 440, 444, 448, 452, 456, 460, 500, 550, 600, 650, 700};
168 int nnrep[] = {10000, 10000, 10000, 10000, 10000, 10000, 10000, 10000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 1000, 400, 400, 400, 400, 400, 200, 200, 200, 200, 200, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 40, 40, 40, 40, 40, 40, 40, 40, 40, 40, 20, 20, 20, 20, 20, 20, 20, 10, 10, 10, 10, 10, 10, 10, 10, 10, 10, 10, 10, 10, 10, 10, 4, 4, 4, 4, 4};
169
170// for(ll=0; ll<24; ll++)
171 for(ll=0; ll<75; ll++)
172// for(ll=0; ll<115; ll++)
173// for(ll=0; ll<120; ll++)
174
175 {
176
177 int n = nn[ll];
178 int nrep = nnrep[ll];
179// int n = ll+1;
180// int nrep = nnrep[0];
181// n = n<12 ? 12 : n;
182// n = n<8 ? 8 : n;
183
184#else
185 int nn[] = {1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24};
186
187 for(ll=0; ll<24; ll++)
188
189 {
190
191 int n = nn[ll];
192 int nrep = 40000; //nnrep[ll];
193#endif
194
195 double *A; d_zeros(&A, n, n);
196 double *B; d_zeros(&B, n, n);
197 double *C; d_zeros(&C, n, n);
198 double *M; d_zeros(&M, n, n);
199
200 char c_n = 'n';
201 char c_l = 'l';
202 char c_r = 'r';
203 char c_t = 't';
204 char c_u = 'u';
205 int i_1 = 1;
206 int i_t;
207 double d_1 = 1;
208 double d_0 = 0;
209
210 for(i=0; i<n*n; i++)
211 A[i] = i;
212
213 for(i=0; i<n; i++)
214 B[i*(n+1)] = 1;
215
216 for(i=0; i<n*n; i++)
217 M[i] = 1;
218
219 int n2 = n*n;
220 double *B2; d_zeros(&B2, n, n);
221 for(i=0; i<n*n; i++)
222 B2[i] = 1e-15;
223 for(i=0; i<n; i++)
224 B2[i*(n+1)] = 1;
225
226 int pnd = ((n+bsd-1)/bsd)*bsd;
227 int cnd = ((n+ncd-1)/ncd)*ncd;
228 int cnd2 = 2*((n+ncd-1)/ncd)*ncd;
229
230 double *x; d_zeros_align(&x, pnd, 1);
231 double *y; d_zeros_align(&y, pnd, 1);
232 double *x2; d_zeros_align(&x2, pnd, 1);
233 double *y2; d_zeros_align(&y2, pnd, 1);
234 double *diag; d_zeros_align(&diag, pnd, 1);
235 int *ipiv; int_zeros(&ipiv, n, 1);
236
237 for(i=0; i<pnd; i++) x[i] = 1;
238 for(i=0; i<pnd; i++) x2[i] = 1;
239
240 // matrix struct
241#if 0
242 struct d_strmat sA; d_allocate_strmat(n+4, n+4, &sA);
243 struct d_strmat sB; d_allocate_strmat(n+4, n+4, &sB);
244 struct d_strmat sC; d_allocate_strmat(n+4, n+4, &sC);
245 struct d_strmat sD; d_allocate_strmat(n+4, n+4, &sD);
246 struct d_strmat sE; d_allocate_strmat(n+4, n+4, &sE);
247#else
248 struct d_strmat sA; d_allocate_strmat(n, n, &sA);
249 struct d_strmat sB; d_allocate_strmat(n, n, &sB);
250 struct d_strmat sB2; d_allocate_strmat(n, n, &sB2);
251 struct d_strmat sB3; d_allocate_strmat(n, n, &sB3);
252 struct d_strmat sC; d_allocate_strmat(n, n, &sC);
253 struct d_strmat sD; d_allocate_strmat(n, n, &sD);
254 struct d_strmat sE; d_allocate_strmat(n, n, &sE);
255#endif
256 struct d_strvec sx; d_allocate_strvec(n, &sx);
257 struct d_strvec sy; d_allocate_strvec(n, &sy);
258 struct d_strvec sz; d_allocate_strvec(n, &sz);
259
260 d_cvt_mat2strmat(n, n, A, n, &sA, 0, 0);
261 d_cvt_mat2strmat(n, n, B, n, &sB, 0, 0);
262 d_cvt_mat2strmat(n, n, B2, n, &sB2, 0, 0);
263 d_cvt_vec2strvec(n, x, &sx, 0);
264 int ii;
265 for(ii=0; ii<n; ii++)
266 {
267 DMATEL_LIBSTR(&sB3, ii, ii) = 1.0;
268// DMATEL_LIBSTR(&sB3, n-1, ii) = 1.0;
269 DMATEL_LIBSTR(&sB3, ii, n-1) = 1.0;
270 DVECEL_LIBSTR(&sx, ii) = 1.0;
271 }
272// d_print_strmat(n, n, &sB3, 0, 0);
273// if(n==20) return;
274
275 int qr_work_size = 0;//dgeqrf_work_size_libstr(n, n);
276 void *qr_work;
277 v_zeros_align(&qr_work, qr_work_size);
278
279 int lq_work_size = 0;//dgelqf_work_size_libstr(n, n);
280 void *lq_work;
281 v_zeros_align(&lq_work, lq_work_size);
282
283 // create matrix to pivot all the time
284// dgemm_nt_libstr(n, n, n, 1.0, &sA, 0, 0, &sA, 0, 0, 1.0, &sB, 0, 0, &sD, 0, 0);
285
286 double *dummy;
287
288 int info;
289
290 /* timing */
291 struct timeval tvm1, tv0, tv1, tv2, tv3, tv4, tv5, tv6, tv7, tv8, tv9, tv10, tv11, tv12, tv13, tv14, tv15, tv16;
292
293 /* warm up */
294 for(rep=0; rep<nrep; rep++)
295 {
296 dgemm_nt_libstr(n, n, n, 1.0, &sA, 0, 0, &sA, 0, 0, 1.0, &sB, 0, 0, &sC, 0, 0);
297 }
298
299 double alpha = 1.0;
300 double beta = 0.0;
301
302 gettimeofday(&tv0, NULL); // stop
303
304 for(rep=0; rep<nrep; rep++)
305 {
306
307// dgemm_nt_lib(n, n, n, 1.0, pA, cnd, pB, cnd, 0.0, pC, cnd, pC, cnd);
308// dgemm_nn_lib(n, n, n, 1.0, pA, cnd, pB, cnd, 0.0, pC, cnd, pC, cnd);
309// dsyrk_nt_l_lib(n, n, n, 1.0, pA, cnd, pB, cnd, 1.0, pC, cnd, pD, cnd);
310// dtrmm_nt_ru_lib(n, n, pA, cnd, pB, cnd, 0, pC, cnd, pD, cnd);
311// dpotrf_nt_l_lib(n, n, pB, cnd, pD, cnd, diag);
312// dsyrk_dpotrf_nt_l_lib(n, n, n, pA, cnd, pA, cnd, 1, pB, cnd, pD, cnd, diag);
313// dsyrk_nt_l_lib(n, n, n, pA, cnd, pA, cnd, 1, pB, cnd, pD, cnd);
314// dpotrf_nt_l_lib(n, n, pD, cnd, pD, cnd, diag);
315// dgetrf_nn_nopivot_lib(n, n, pB, cnd, pB, cnd, diag);
316// dgetrf_nn_lib(n, n, pB, cnd, pB, cnd, diag, ipiv);
317// dtrsm_nn_ll_one_lib(n, n, pD, cnd, pB, cnd, pB, cnd);
318// dtrsm_nn_lu_inv_lib(n, n, pD, cnd, diag, pB, cnd, pB, cnd);
319 }
320
321 gettimeofday(&tv1, NULL); // stop
322
323 for(rep=0; rep<nrep; rep++)
324 {
325// kernel_dgemm_nt_12x4_lib4(n, &alpha, sA.pA, sA.cn, sB.pA, &beta, sD.pA, sD.cn, sD.pA, sD.cn);
326// kernel_dgemm_nt_8x8_lib4(n, &alpha, sA.pA, sA.cn, sB.pA, sB.cn, &beta, sD.pA, sD.cn, sD.pA, sD.cn);
327// kernel_dsyrk_nt_l_8x8_lib4(n, &alpha, sA.pA, sA.cn, sB.pA, sB.cn, &beta, sD.pA, sD.cn, sD.pA, sD.cn);
328// kernel_dgemm_nt_8x4_lib4(n, &alpha, sA.pA, sA.cn, sB.pA, &beta, sD.pA, sD.cn, sD.pA, sD.cn);
329// kernel_dgemm_nt_4x8_lib4(n, &alpha, sA.pA, sB.pA, sB.cn, &beta, sD.pA, sD.pA);
330// kernel_dgemm_nt_4x4_lib4(n, &alpha, sA.pA, sB.pA, &beta, sD.pA, sD.pA);
331// kernel_dger4_12_sub_lib4(n, sA.pA, sA.cn, sB.pA, sD.pA, sD.cn);
332// kernel_dger4_sub_12r_lib4(n, sA.pA, sA.cn, sB.pA, sD.pA, sD.cn);
333// kernel_dger4_sub_8r_lib4(n, sA.pA, sA.cn, sB.pA, sD.pA, sD.cn);
334// kernel_dger12_add_4r_lib4(n, sA.pA, sB.pA, sB.cn, sD.pA);
335// kernel_dger8_add_4r_lib4(n, sA.pA, sB.pA, sB.cn, sD.pA);
336// kernel_dger4_sub_4r_lib4(n, sA.pA, sB.pA, sD.pA);
337// kernel_dger2_sub_4r_lib4(n, sA.pA, sB.pA, sD.pA);
338// kernel_dger4_sub_8c_lib4(n, sA.pA, sA.cn, sB.pA, sD.pA, sD.cn);
339// kernel_dger4_sub_4c_lib4(n, sA.pA, sA.cn, sB.pA, sD.pA, sD.cn);
340// kernel_dgemm_nn_4x12_lib4(n, &alpha, sA.pA, 0, sB.pA, sB.cn, &beta, sD.pA, sD.pA);
341// kernel_dgemm_nn_4x8_lib4(n, &alpha, sA.pA, 0, sB.pA, sB.cn, &beta, sD.pA, sD.pA);
342// kernel_dgemm_nn_4x4_lib4(n, &alpha, sA.pA, 0, sB.pA, sB.cn, &beta, sD.pA, sD.pA);
343
344 dgemm_nt_libstr(n, n, n, 1.0, &sA, 0, 0, &sB, 0, 0, 0.0, &sC, 0, 0, &sD, 0, 0);
345// dgemm_nn_libstr(n, n, n, 1.0, &sA, 0, 0, &sB, 0, 0, 0.0, &sC, 0, 0, &sD, 0, 0);
346// dsyrk_ln_libstr(n, n, 1.0, &sA, 0, 0, &sA, 0, 0, 0.0, &sC, 0, 0, &sD, 0, 0);
347// dsyrk_ln_mn_libstr(n, n, n, 1.0, &sA, 0, 0, &sA, 0, 0, 0.0, &sC, 0, 0, &sD, 0, 0);
348// dpotrf_l_mn_libstr(n, n, &sB, 0, 0, &sB, 0, 0);
349// dpotrf_l_libstr(n, &sB, 0, 0, &sB, 0, 0);
350// dgetrf_nopivot_libstr(n, n, &sB, 0, 0, &sB, 0, 0);
351// dgetrf_libstr(n, n, &sB, 0, 0, &sB, 0, 0, ipiv);
352// dgeqrf_libstr(n, n, &sC, 0, 0, &sD, 0, 0, qr_work);
353// dcolin_libstr(n, &sx, 0, &sB3, 0, n-1);
354// dgelqf_libstr(n, n, &sB3, 0, 0, &sB3, 0, 0, lq_work);
355// dtrmm_rlnn_libstr(n, n, 1.0, &sA, 0, 0, &sD, 0, 0, &sD, 0, 0); //
356// dtrmm_rutn_libstr(n, n, 1.0, &sA, 0, 0, &sB, 0, 0, &sD, 0, 0);
357// dtrsm_llnu_libstr(n, n, 1.0, &sD, 0, 0, &sB, 0, 0, &sB, 0, 0);
358// dtrsm_lunn_libstr(n, n, 1.0, &sD, 0, 0, &sB, 0, 0, &sB, 0, 0);
359// dtrsm_rltn_libstr(n, n, 1.0, &sB2, 0, 0, &sD, 0, 0, &sD, 0, 0); //
360// dtrsm_rltu_libstr(n, n, 1.0, &sD, 0, 0, &sB, 0, 0, &sB, 0, 0);
361// dtrsm_rutn_libstr(n, n, 1.0, &sD, 0, 0, &sB, 0, 0, &sB, 0, 0);
362// dgemv_n_libstr(n, n, 1.0, &sA, 0, 0, &sx, 0, 0.0, &sy, 0, &sz, 0);
363// dgemv_t_libstr(n, n, 1.0, &sA, 0, 0, &sx, 0, 0.0, &sy, 0, &sz, 0);
364// dsymv_l_libstr(n, n, 1.0, &sA, 0, 0, &sx, 0, 0.0, &sy, 0, &sz, 0);
365// dgemv_nt_libstr(n, n, 1.0, 1.0, &sA, 0, 0, &sx, 0, &sx, 0, 0.0, 0.0, &sy, 0, &sy, 0, &sz, 0, &sz, 0);
366 }
367
368// d_print_strmat(n, n, &sD, 0, 0);
369
370 gettimeofday(&tv2, NULL); // stop
371
372 for(rep=0; rep<nrep; rep++)
373 {
374#if defined(REF_BLAS_OPENBLAS) || defined(REF_BLAS_NETLIB) || defined(REF_BLAS_MKL)
375// dgemm_(&c_n, &c_t, &n, &n, &n, &d_1, A, &n, M, &n, &d_0, C, &n);
376// dpotrf_(&c_l, &n, B2, &n, &info);
377// dgemm_(&c_n, &c_n, &n, &n, &n, &d_1, A, &n, M, &n, &d_0, C, &n);
378// dsyrk_(&c_l, &c_n, &n, &n, &d_1, A, &n, &d_0, C, &n);
379// dtrmm_(&c_r, &c_u, &c_t, &c_n, &n, &n, &d_1, A, &n, C, &n);
380// dgetrf_(&n, &n, B2, &n, ipiv, &info);
381// dtrsm_(&c_l, &c_l, &c_n, &c_u, &n, &n, &d_1, B2, &n, B, &n);
382// dtrsm_(&c_l, &c_u, &c_n, &c_n, &n, &n, &d_1, B2, &n, B, &n);
383// dtrtri_(&c_l, &c_n, &n, B2, &n, &info);
384// dlauum_(&c_l, &n, B, &n, &info);
385// dgemv_(&c_n, &n, &n, &d_1, A, &n, x, &i_1, &d_0, y, &i_1);
386// dgemv_(&c_t, &n, &n, &d_1, A, &n, x2, &i_1, &d_0, y2, &i_1);
387// dtrmv_(&c_l, &c_n, &c_n, &n, B, &n, x, &i_1);
388// dtrsv_(&c_l, &c_n, &c_n, &n, B, &n, x, &i_1);
389// dsymv_(&c_l, &n, &d_1, A, &n, x, &i_1, &d_0, y, &i_1);
390
391// for(i=0; i<n; i++)
392// {
393// i_t = n-i;
394// dcopy_(&i_t, &B[i*(n+1)], &i_1, &C[i*(n+1)], &i_1);
395// }
396// dsyrk_(&c_l, &c_n, &n, &n, &d_1, A, &n, &d_1, C, &n);
397// dpotrf_(&c_l, &n, C, &n, &info);
398
399#endif
400
401#if defined(REF_BLAS_BLIS)
402// dgemm_(&c_n, &c_t, &n77, &n77, &n77, &d_1, A, &n77, B, &n77, &d_0, C, &n77);
403// dgemm_(&c_n, &c_n, &n77, &n77, &n77, &d_1, A, &n77, B, &n77, &d_0, C, &n77);
404// dsyrk_(&c_l, &c_n, &n77, &n77, &d_1, A, &n77, &d_0, C, &n77);
405// dtrmm_(&c_r, &c_u, &c_t, &c_n, &n77, &n77, &d_1, A, &n77, C, &n77);
406// dpotrf_(&c_l, &n77, B, &n77, &info);
407// dtrtri_(&c_l, &c_n, &n77, B, &n77, &info);
408// dlauum_(&c_l, &n77, B, &n77, &info);
409#endif
410 }
411
412 gettimeofday(&tv3, NULL); // stop
413
414 float Gflops_max = flops_max * GHz_max;
415
416// float flop_operation = 4*16.0*2*n; // kernel 12x4
417// float flop_operation = 3*16.0*2*n; // kernel 12x4
418// float flop_operation = 2*16.0*2*n; // kernel 8x4
419// float flop_operation = 1*16.0*2*n; // kernel 4x4
420// float flop_operation = 0.5*16.0*2*n; // kernel 2x4
421
422 float flop_operation = 2.0*n*n*n; // dgemm
423// float flop_operation = 1.0*n*n*n; // dsyrk dtrmm dtrsm
424// float flop_operation = 1.0/3.0*n*n*n; // dpotrf dtrtri
425// float flop_operation = 2.0/3.0*n*n*n; // dgetrf
426// float flop_operation = 4.0/3.0*n*n*n; // dgeqrf
427// float flop_operation = 2.0*n*n; // dgemv dsymv
428// float flop_operation = 1.0*n*n; // dtrmv dtrsv
429// float flop_operation = 4.0*n*n; // dgemv_nt
430
431// float flop_operation = 4.0/3.0*n*n*n; // dsyrk+dpotrf
432
433 float time_hpmpc = (float) (tv1.tv_sec-tv0.tv_sec)/(nrep+0.0)+(tv1.tv_usec-tv0.tv_usec)/(nrep*1e6);
434 float time_blasfeo = (float) (tv2.tv_sec-tv1.tv_sec)/(nrep+0.0)+(tv2.tv_usec-tv1.tv_usec)/(nrep*1e6);
435 float time_blas = (float) (tv3.tv_sec-tv2.tv_sec)/(nrep+0.0)+(tv3.tv_usec-tv2.tv_usec)/(nrep*1e6);
436
437 float Gflops_hpmpc = 1e-9*flop_operation/time_hpmpc;
438 float Gflops_blasfeo = 1e-9*flop_operation/time_blasfeo;
439 float Gflops_blas = 1e-9*flop_operation/time_blas;
440
441
442// printf("%d\t%7.2f\t%7.2f\t%7.2f\t%7.2f\t%7.2f\t%7.2f\n", n, Gflops_hpmpc, 100.0*Gflops_hpmpc/Gflops_max, Gflops_blasfeo, 100.0*Gflops_blasfeo/Gflops_max, Gflops_blas, 100.0*Gflops_blas/Gflops_max);
443// fprintf(f, "%d\t%7.2f\t%7.2f\t%7.2f\t%7.2f\t%7.2f\t%7.2f\n", n, Gflops_hpmpc, 100.0*Gflops_hpmpc/Gflops_max, Gflops_blasfeo, 100.0*Gflops_blasfeo/Gflops_max, Gflops_blas, 100.0*Gflops_blas/Gflops_max);
444 printf("%d\t%7.2f\t%7.2f\t%7.2f\t%7.2f\n", n, Gflops_blasfeo, 100.0*Gflops_blasfeo/Gflops_max, Gflops_blas, 100.0*Gflops_blas/Gflops_max);
445// fprintf(f, "%d\t%7.2f\t%7.2f\t%7.2f\t%7.2f\n", n, Gflops_blasfeo, 100.0*Gflops_blasfeo/Gflops_max, Gflops_blas, 100.0*Gflops_blas/Gflops_max);
446
447
448 d_free(A);
449 d_free(B);
450 d_free(B2);
451 d_free(M);
452 d_free_align(x);
453 d_free_align(y);
454 d_free_align(x2);
455 d_free_align(y2);
456 int_free(ipiv);
457 free(qr_work);
458 free(lq_work);
459
460 d_free_strmat(&sA);
461 d_free_strmat(&sB);
462 d_free_strmat(&sB2);
463 d_free_strmat(&sB3);
464 d_free_strmat(&sC);
465 d_free_strmat(&sD);
466 d_free_strmat(&sE);
467 d_free_strvec(&sx);
468 d_free_strvec(&sy);
469 d_free_strvec(&sz);
470
471 }
472
473 printf("\n");
474
475// fprintf(f, "];\n");
476// fclose(f);
477
478 return 0;
479
480 }