You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

omatcopy_rt_bench.c 18 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444
  1. /* +r: %0 = src, %1 = dst, %2 = src_ld, %3 = dst_ld, %4 = dst_tmp */
  2. /* m: %5 = num_rows, %6 = alpha */
  3. /* xmm15 = alpha */
  4. #ifdef HAVE_AVX
  5. #define TRANS_4x4(a1_no,a2_no,a3_no,a4_no,t1_no,t2_no,t3_no,t4_no)\
  6. "vunpcklps %%xmm"#a2_no",%%xmm"#a1_no",%%xmm"#t1_no"; vunpckhps %%xmm"#a2_no",%%xmm"#a1_no",%%xmm"#t2_no";"\
  7. "vunpcklps %%xmm"#a4_no",%%xmm"#a3_no",%%xmm"#t3_no"; vunpckhps %%xmm"#a4_no",%%xmm"#a3_no",%%xmm"#t4_no";"\
  8. "vunpcklpd %%xmm"#t3_no",%%xmm"#t1_no",%%xmm"#a1_no"; vunpckhpd %%xmm"#t3_no",%%xmm"#t1_no",%%xmm"#a2_no";"\
  9. "vunpcklpd %%xmm"#t4_no",%%xmm"#t2_no",%%xmm"#a3_no"; vunpckhpd %%xmm"#t4_no",%%xmm"#t2_no",%%xmm"#a4_no";"
  10. #define TRANS_4x8(a1_no,a2_no,a3_no,a4_no,t1_no,t2_no,t3_no,t4_no)\
  11. "vunpcklps %%ymm"#a2_no",%%ymm"#a1_no",%%ymm"#t1_no"; vunpckhps %%ymm"#a2_no",%%ymm"#a1_no",%%ymm"#t2_no";"\
  12. "vunpcklps %%ymm"#a4_no",%%ymm"#a3_no",%%ymm"#t3_no"; vunpckhps %%ymm"#a4_no",%%ymm"#a3_no",%%ymm"#t4_no";"\
  13. "vunpcklpd %%ymm"#t3_no",%%ymm"#t1_no",%%ymm"#a1_no"; vunpckhpd %%ymm"#t3_no",%%ymm"#t1_no",%%ymm"#a2_no";"\
  14. "vunpcklpd %%ymm"#t4_no",%%ymm"#t2_no",%%ymm"#a3_no"; vunpckhpd %%ymm"#t4_no",%%ymm"#t2_no",%%ymm"#a4_no";"
  15. #define SAVE_4x4(b1_no,b2_no,b3_no,b4_no)\
  16. "vmovups %%xmm"#b1_no",(%4); vmovups %%xmm"#b2_no",(%4,%3,1); leaq (%4,%3,2),%4;"\
  17. "vmovups %%xmm"#b3_no",(%4); vmovups %%xmm"#b4_no",(%4,%3,1); leaq (%4,%3,2),%4;"
  18. #define SAVE_4x8(b1_no,b2_no,b3_no,b4_no) SAVE_4x4(b1_no,b2_no,b3_no,b4_no)\
  19. "vextractf128 $1,%%ymm"#b1_no",(%4); vextractf128 $1,%%ymm"#b2_no",(%4,%3,1); leaq (%4,%3,2),%4;"\
  20. "vextractf128 $1,%%ymm"#b3_no",(%4); vextractf128 $1,%%ymm"#b4_no",(%4,%3,1); leaq (%4,%3,2),%4;"
  21. #define COPY_4x16 "movq %1,%4; addq $16,%1;"\
  22. "vmulps (%0),%%ymm15,%%ymm0; vmulps 32(%0),%%ymm15,%%ymm4; vmulps (%0,%2,1),%%ymm15,%%ymm1; vmulps 32(%0,%2,1),%%ymm15,%%ymm5; leaq (%0,%2,2),%0;"\
  23. "vmulps (%0),%%ymm15,%%ymm2; vmulps 32(%0),%%ymm15,%%ymm6; vmulps (%0,%2,1),%%ymm15,%%ymm3; vmulps 32(%0,%2,1),%%ymm15,%%ymm7; leaq (%0,%2,2),%0;"\
  24. TRANS_4x8(0,1,2,3,8,9,10,11) SAVE_4x8(0,1,2,3)\
  25. TRANS_4x8(4,5,6,7,8,9,10,11) SAVE_4x8(4,5,6,7)
  26. #define COPY_4x8 "movq %1,%4; addq $16,%1;"\
  27. "vmulps (%0),%%ymm15,%%ymm0; vmulps (%0,%2,1),%%ymm15,%%ymm1; leaq (%0,%2,2),%0;"\
  28. "vmulps (%0),%%ymm15,%%ymm2; vmulps (%0,%2,1),%%ymm15,%%ymm3; leaq (%0,%2,2),%0;"\
  29. TRANS_4x8(0,1,2,3,8,9,10,11) SAVE_4x8(0,1,2,3)
  30. #define COPY_4x4 "movq %1,%4; addq $16,%1;"\
  31. "vmulps (%0),%%xmm15,%%xmm0; vmulps (%0,%2,1),%%xmm15,%%xmm1; leaq (%0,%2,2),%0;"\
  32. "vmulps (%0),%%xmm15,%%xmm2; vmulps (%0,%2,1),%%xmm15,%%xmm3; leaq (%0,%2,2),%0;"\
  33. TRANS_4x4(0,1,2,3,8,9,10,11) SAVE_4x4(0,1,2,3)
  34. #define COPY_4x2 \
  35. "vmovsd (%0),%%xmm0; vmovhpd (%0,%2,1),%%xmm0,%%xmm0; vmulps %%xmm15,%%xmm0,%%xmm0; leaq (%0,%2,2),%0;"\
  36. "vmovsd (%0),%%xmm1; vmovhpd (%0,%2,1),%%xmm1,%%xmm1; vmulps %%xmm15,%%xmm1,%%xmm1; leaq (%0,%2,2),%0;"\
  37. "vpermilps $216,%%xmm0,%%xmm0; vpermilps $216,%%xmm1,%%xmm1; vunpcklpd %%xmm1,%%xmm0,%%xmm2; vunpckhpd %%xmm1,%%xmm0,%%xmm3;"\
  38. "vmovups %%xmm2,(%1); vmovups %%xmm3,(%1,%3,1); addq $16,%1;"
  39. #define COPY_4x1 \
  40. "vmovss (%0),%%xmm0; vinsertps $16,(%0,%2,1),%%xmm0,%%xmm0; leaq (%0,%2,2),%0;"\
  41. "vinsertps $32,(%0),%%xmm0,%%xmm0; vinsertps $48,(%0,%2,1),%%xmm0,%%xmm0; leaq (%0,%2,2),%0;"\
  42. "vmulps %%xmm15,%%xmm0,%%xmm0; vmovups %%xmm0,(%1); addq $16,%1;"
  43. #define SAVE_2x4(c1_no,c2_no,t1_no,t2_no) \
  44. "vunpcklps %%xmm"#c2_no",%%xmm"#c1_no",%%xmm"#t1_no"; vmulps %%xmm15,%%xmm"#t1_no",%%xmm"#t1_no";"\
  45. "vmovsd %%xmm"#t1_no",(%4); vmovhpd %%xmm"#t1_no",(%4,%3,1); leaq (%4,%3,2),%4;"\
  46. "vunpckhps %%xmm"#c2_no",%%xmm"#c1_no",%%xmm"#t2_no"; vmulps %%xmm15,%%xmm"#t2_no",%%xmm"#t2_no";"\
  47. "vmovsd %%xmm"#t2_no",(%4); vmovhpd %%xmm"#t2_no",(%4,%3,1); leaq (%4,%3,2),%4;"
  48. #define COPY_2x16 "movq %1,%4; addq $8,%1;"\
  49. "vmovups (%0),%%ymm0; vmovups 32(%0),%%ymm2; vmovups (%0,%2,1),%%ymm1; vmovups 32(%0,%2,1),%%ymm3; leaq (%0,%2,2),%0;"\
  50. "vextractf128 $1,%%ymm0,%%xmm4; vextractf128 $1,%%ymm2,%%xmm6; vextractf128 $1,%%ymm1,%%xmm5; vextractf128 $1,%%ymm3,%%xmm7;"\
  51. SAVE_2x4(0,1,8,9) SAVE_2x4(4,5,8,9) SAVE_2x4(2,3,8,9) SAVE_2x4(6,7,8,9)
  52. #define COPY_2x8 "movq %1,%4; addq $8,%1;"\
  53. "vmovups (%0),%%ymm0; vmovups (%0,%2,1),%%ymm1; leaq (%0,%2,2),%0;"\
  54. "vextractf128 $1,%%ymm0,%%xmm2; vextractf128 $1,%%ymm1,%%xmm3;"\
  55. SAVE_2x4(0,1,4,5) SAVE_2x4(2,3,4,5)
  56. #define COPY_2x4 "movq %1,%4; addq $8,%1;"\
  57. "vmovups (%0),%%xmm0; vmovups (%0,%2,1),%%xmm1; leaq (%0,%2,2),%0;"\
  58. SAVE_2x4(0,1,4,5)
  59. #define COPY_2x2 \
  60. "vmovsd (%0),%%xmm0; vmovhpd (%0,%2,1),%%xmm0,%%xmm0; vmulps %%xmm15,%%xmm0,%%xmm0; leaq (%0,%2,2),%0; vpermilps $216,%%xmm0,%%xmm0;"\
  61. "vmovsd %%xmm0,(%1); vmovhpd %%xmm0,(%1,%3,1); addq $8,%1;"
  62. #define COPY_2x1 \
  63. "vmovss (%0),%%xmm0; vinsertps $16,(%0,%2,1),%%xmm0,%%xmm0; vmulps %%xmm15,%%xmm0,%%xmm0; leaq (%0,%2,2),%0; vmovsd %%xmm0,(%1); addq $8,%1;"
  64. #define SAVE_1x4(c1_no)\
  65. "vmulps %%xmm15,%%xmm"#c1_no",%%xmm"#c1_no"; vmovss %%xmm"#c1_no",(%4); vextractps $1,%%xmm"#c1_no",(%4,%3,1); leaq (%4,%3,2),%4;"\
  66. "vextractps $2,%%xmm"#c1_no",(%4); vextractps $3,%%xmm"#c1_no",(%4,%3,1); leaq (%4,%3,2),%4;"
  67. #define COPY_1x16 "movq %1,%4; addq $4,%1;"\
  68. "vmovups (%0),%%xmm1;" SAVE_1x4(1) "vmovups 16(%0),%%xmm2;" SAVE_1x4(2)\
  69. "vmovups 32(%0),%%xmm1;" SAVE_1x4(1) "vmovups 48(%0),%%xmm2;" SAVE_1x4(2) "addq %2,%0;"
  70. #define COPY_1x8 "movq %1,%4; addq $4,%1;"\
  71. "vmovups (%0),%%xmm1;" SAVE_1x4(1) "vmovups 16(%0),%%xmm2;" SAVE_1x4(2) "addq %2,%0;"
  72. #define COPY_1x4 "movq %1,%4; addq $4,%1; vmovups (%0),%%xmm1;" SAVE_1x4(1) "addq %2,%0;"
  73. #define COPY_1x2 "vmovsd (%0),%%xmm1; addq %2,%0; vmulps %%xmm15,%%xmm1,%%xmm1; vmovss %%xmm1,(%1); vextractps $1,%%xmm1,(%1,%3,1); addq $4,%1;"
  74. #define COPY_1x1 "vmulss (%0),%%xmm15,%%xmm1; vmovss %%xmm1,(%1); addq %2,%0; addq $4,%1;"
  75. #define COMPUTE(ndim){\
  76. src = src_base; dst = dst_base;\
  77. __asm__ __volatile__(\
  78. "vbroadcastss %6,%%ymm15; movq %5,%%r11; cmpq $4,%%r11; jb "#ndim"32f;"\
  79. #ndim"31:\n\t"\
  80. COPY_4x##ndim "subq $4,%%r11; cmpq $4,%%r11; jnb "#ndim"31b;"\
  81. #ndim"32:\n\t"\
  82. "cmpq $2,%%r11; jb "#ndim"33f;"\
  83. COPY_2x##ndim "subq $2,%%r11;"\
  84. #ndim"33:\n\t"\
  85. "testq %%r11,%%r11; jz "#ndim"34f;"\
  86. COPY_1x##ndim "subq $1,%%r11;"\
  87. #ndim"34:\n\t"\
  88. :"+r"(src),"+r"(dst),"+r"(src_ld_bytes),"+r"(dst_ld_bytes),"+r"(dst_tmp):"m"(num_rows),"m"(ALPHA):"r11","cc","memory"\
  89. ,"xmm0","xmm1","xmm2","xmm3","xmm4","xmm5","xmm6","xmm7","xmm8","xmm9","xmm10","xmm11","xmm12","xmm13","xmm14","xmm15");\
  90. }
  91. #endif
  92. #ifndef FLOAT
  93. #define FLOAT float
  94. #endif
  95. #ifdef CNAME
  96. #include "common.h"
  97. #endif
  98. #define ROWS_OF_BLOCK 384
  99. #ifndef BLASLONG
  100. #define BLASLONG long
  101. #endif
  102. #include <stdint.h>
  103. #include <string.h>
  104. #ifdef HAVE_AVX
  105. int CNAME(BLASLONG rows, BLASLONG cols, FLOAT alpha, FLOAT *a, BLASLONG lda, FLOAT *b, BLASLONG ldb){
  106. float *src, *dst, *dst_tmp, *src_base, *dst_base;
  107. uint64_t src_ld_bytes = (uint64_t)lda * sizeof(float), dst_ld_bytes = (uint64_t)ldb * sizeof(float), num_rows = 0;
  108. BLASLONG cols_left, rows_done; float ALPHA = alpha;
  109. if(ALPHA==0.0){
  110. dst_base = b;
  111. for(cols_left=cols;cols_left>0;cols_left--) {memset(dst_base,0,rows*sizeof(float)); dst_base += ldb;}
  112. return 0;
  113. }
  114. for(rows_done=0;rows_done<rows;rows_done+=num_rows){
  115. num_rows = rows-rows_done;
  116. if(num_rows > ROWS_OF_BLOCK) num_rows = ROWS_OF_BLOCK;
  117. cols_left = cols; src_base = a + (int64_t)lda * (int64_t)rows_done; dst_base = b + rows_done;
  118. if(ldb%1024>3 && ldb%1024<1021) for(;cols_left>15;cols_left-=16){COMPUTE(16) src_base += 16; dst_base += 16 * ldb;}
  119. for(;cols_left>7;cols_left-=8){COMPUTE(8) src_base += 8; dst_base += 8 * ldb;}
  120. for(;cols_left>3;cols_left-=4){COMPUTE(4) src_base += 4; dst_base += 4 * ldb;}
  121. for(;cols_left>1;cols_left-=2){COMPUTE(2) src_base += 2; dst_base += 2 * ldb;}
  122. if(cols_left>0){COMPUTE(1) src_base ++; dst_base += ldb;}
  123. }
  124. return 0;
  125. }
  126. #endif
  127. int CNAME2(BLASLONG rows, BLASLONG cols, FLOAT alpha, FLOAT *a, BLASLONG lda, FLOAT *b, BLASLONG ldb)
  128. {
  129. BLASLONG i, j;
  130. FLOAT *a_offset, *a_offset1, *a_offset2, *a_offset3, *a_offset4;
  131. FLOAT *b_offset, *b_offset1, *b_offset2, *b_offset3, *b_offset4;
  132. if (rows <= 0) return 0;
  133. if (cols <= 0) return 0;
  134. a_offset = a;
  135. b_offset = b;
  136. i = (rows >> 2);
  137. if (i > 0) {
  138. do {
  139. a_offset1 = a_offset;
  140. a_offset2 = a_offset1 + lda;
  141. a_offset3 = a_offset2 + lda;
  142. a_offset4 = a_offset3 + lda;
  143. a_offset += 4 * lda;
  144. b_offset1 = b_offset;
  145. b_offset2 = b_offset1 + ldb;
  146. b_offset3 = b_offset2 + ldb;
  147. b_offset4 = b_offset3 + ldb;
  148. b_offset += 4;
  149. j = (cols >> 2);
  150. if (j > 0) {
  151. do {
  152. /* Column 1 of MAT_B */
  153. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha; // Row 1 of MAT_A
  154. *(b_offset2 + 0) = *(a_offset1 + 1)*alpha;
  155. *(b_offset3 + 0) = *(a_offset1 + 2)*alpha;
  156. *(b_offset4 + 0) = *(a_offset1 + 3)*alpha;
  157. /* Column 2 of MAT_B */
  158. *(b_offset1 + 1) = *(a_offset2 + 0)*alpha; // Row 2 of MAT_A
  159. *(b_offset2 + 1) = *(a_offset2 + 1)*alpha;
  160. *(b_offset3 + 1) = *(a_offset2 + 2)*alpha;
  161. *(b_offset4 + 1) = *(a_offset2 + 3)*alpha;
  162. /* Column 3 of MAT_B */
  163. *(b_offset1 + 2) = *(a_offset3 + 0)*alpha; // Row 3 of MAT_A
  164. *(b_offset2 + 2) = *(a_offset3 + 1)*alpha;
  165. *(b_offset3 + 2) = *(a_offset3 + 2)*alpha;
  166. *(b_offset4 + 2) = *(a_offset3 + 3)*alpha;
  167. /* Column 4 of MAT_B */
  168. *(b_offset1 + 3) = *(a_offset4 + 0)*alpha; // Row 4 of MAT_A
  169. *(b_offset2 + 3) = *(a_offset4 + 1)*alpha;
  170. *(b_offset3 + 3) = *(a_offset4 + 2)*alpha;
  171. *(b_offset4 + 3) = *(a_offset4 + 3)*alpha;
  172. a_offset1 += 4;
  173. a_offset2 += 4;
  174. a_offset3 += 4;
  175. a_offset4 += 4;
  176. b_offset1 += ldb * 4;
  177. b_offset2 += ldb * 4;
  178. b_offset3 += ldb * 4;
  179. b_offset4 += ldb * 4;
  180. j--;
  181. } while (j > 0);
  182. } // if(j > 0)
  183. if (cols & 2) {
  184. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  185. *(b_offset2 + 0) = *(a_offset1 + 1)*alpha;
  186. *(b_offset1 + 1) = *(a_offset2 + 0)*alpha;
  187. *(b_offset2 + 1) = *(a_offset2 + 1)*alpha;
  188. *(b_offset1 + 2) = *(a_offset3 + 0)*alpha;
  189. *(b_offset2 + 2) = *(a_offset3 + 1)*alpha;
  190. *(b_offset1 + 3) = *(a_offset4 + 0)*alpha;
  191. *(b_offset2 + 3) = *(a_offset4 + 1)*alpha;
  192. a_offset1 += 2;
  193. a_offset2 += 2;
  194. a_offset3 += 2;
  195. a_offset4 += 2;
  196. b_offset1 += ldb*2;
  197. }
  198. if (cols & 1) {
  199. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  200. *(b_offset1 + 1) = *(a_offset2 + 0)*alpha;
  201. *(b_offset1 + 2) = *(a_offset3 + 0)*alpha;
  202. *(b_offset1 + 3) = *(a_offset4 + 0)*alpha;
  203. }
  204. i--;
  205. } while (i > 0);
  206. }
  207. if (rows & 2) {
  208. a_offset1 = a_offset;
  209. a_offset2 = a_offset1 + lda;
  210. a_offset += 2 * lda;
  211. b_offset1 = b_offset;
  212. b_offset2 = b_offset1 + ldb;
  213. b_offset3 = b_offset2 + ldb;
  214. b_offset4 = b_offset3 + ldb;
  215. b_offset += 2;
  216. j = (cols >> 2);
  217. if (j > 0){
  218. do {
  219. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  220. *(b_offset2 + 0) = *(a_offset1 + 1)*alpha;
  221. *(b_offset3 + 0) = *(a_offset1 + 2)*alpha;
  222. *(b_offset4 + 0) = *(a_offset1 + 3)*alpha;
  223. *(b_offset1 + 1) = *(a_offset2 + 0)*alpha;
  224. *(b_offset2 + 1) = *(a_offset2 + 1)*alpha;
  225. *(b_offset3 + 1) = *(a_offset2 + 2)*alpha;
  226. *(b_offset4 + 1) = *(a_offset2 + 3)*alpha;
  227. a_offset1 += 4;
  228. a_offset2 += 4;
  229. b_offset1 += ldb * 4;
  230. b_offset2 += ldb * 4;
  231. b_offset3 += ldb * 4;
  232. b_offset4 += ldb * 4;
  233. j--;
  234. } while (j > 0);
  235. }
  236. if (cols & 2){
  237. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  238. *(b_offset2 + 0) = *(a_offset1 + 1)*alpha;
  239. *(b_offset1 + 1) = *(a_offset2 + 0)*alpha;
  240. *(b_offset2 + 1) = *(a_offset2 + 1)*alpha;
  241. a_offset1 += 2;
  242. a_offset2 += 2;
  243. b_offset1 += ldb*2;
  244. }
  245. if (cols & 1){
  246. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  247. *(b_offset1 + 1) = *(a_offset2 + 0)*alpha;
  248. }
  249. } // if (rows & 2)
  250. if (rows & 1) {
  251. a_offset1 = a_offset;
  252. a_offset += lda;
  253. b_offset1 = b_offset;
  254. b_offset2 = b_offset1 + ldb;
  255. b_offset3 = b_offset2 + ldb;
  256. b_offset4 = b_offset3 + ldb;
  257. j = (cols >> 2);
  258. if (j > 0){
  259. do {
  260. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  261. *(b_offset2 + 0) = *(a_offset1 + 1)*alpha;
  262. *(b_offset3 + 0) = *(a_offset1 + 2)*alpha;
  263. *(b_offset4 + 0) = *(a_offset1 + 3)*alpha;
  264. a_offset1 += 4;
  265. b_offset1 += ldb * 4;
  266. b_offset2 += ldb * 4;
  267. b_offset3 += ldb * 4;
  268. b_offset4 += ldb * 4;
  269. j--;
  270. } while (j > 0);
  271. }
  272. if (cols & 2){
  273. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  274. *(b_offset2 + 0) = *(a_offset1 + 1)*alpha;
  275. a_offset1 += 2;
  276. b_offset1 += ldb * 2;
  277. }
  278. if (cols & 1){
  279. *(b_offset1 + 0) = *(a_offset1 + 0)*alpha;
  280. }
  281. }
  282. return 0;
  283. }
  284. #ifndef CNAME
  285. int REF(BLASLONG rows, BLASLONG cols, FLOAT alpha, FLOAT *a, BLASLONG lda, FLOAT *b, BLASLONG ldb){
  286. BLASLONG i,j;
  287. FLOAT *aptr,*bptr;
  288. if (rows <= 0) return 0;
  289. if (cols <= 0) return 0;
  290. aptr = a;
  291. for(i=0;i<rows;i++){
  292. bptr = &b[i];
  293. for(j=0;j<cols;j++) bptr[j*ldb] = alpha * aptr[j];
  294. aptr += lda;
  295. }
  296. return 0;
  297. }
  298. #include <stdio.h>
  299. #include <stdlib.h>
  300. #include <sys/time.h>
  301. int main(){
  302. struct timeval time1,time2,time3,time4;
  303. int testdim[5]={2000,4000,6000,8000,10000}; long length = testdim[4]*testdim[4];
  304. const float alpha = 1.5; const long lda = testdim[4] - testdim[4]%60, ldb = testdim[4] - testdim[4]%70;
  305. float *A = (float *)malloc(length*sizeof(float));
  306. float *B1 = (float *)malloc(length*sizeof(float));
  307. float *B2 = (float *)malloc(length*sizeof(float));
  308. float *B3 = (float *)malloc(length*sizeof(float));
  309. long count, rows, cols, iter, usec1, usec2, usec3; float tmpdif,maxdif,tmpdif2,maxdif2;
  310. for(count=0;count<length;count++) A[count]=B1[count]=B2[count]=B3[count]=(float)rand()/(float)RAND_MAX;
  311. for(iter=0;iter<5;iter++){
  312. rows = testdim[iter] - testdim[iter] % 23;
  313. cols = testdim[iter] - testdim[iter] % 53;
  314. #ifdef HAVE_AVX
  315. CNAME(rows,cols,alpha,A,lda,B2,ldb);
  316. #endif
  317. REF(rows,cols,alpha,A,lda,B1,ldb);
  318. CNAME2(rows,cols,alpha,A,lda,B3,ldb);
  319. printf("\nTest omatcopy_rt with rows = %ld and cols = %ld.\n",rows,cols);
  320. gettimeofday(&time1,0);
  321. #ifdef HAVE_AVX
  322. CNAME(rows,cols,alpha,A,lda,B2,ldb);
  323. #endif
  324. gettimeofday(&time2,0);
  325. REF(rows,cols,alpha,A,lda,B1,ldb);
  326. gettimeofday(&time3,0);
  327. CNAME2(rows,cols,alpha,A,lda,B3,ldb);
  328. gettimeofday(&time4,0);
  329. usec1 = (time2.tv_sec - time1.tv_sec) * 1000000 + (time2.tv_usec - time1.tv_usec);
  330. usec2 = (time3.tv_sec - time2.tv_sec) * 1000000 + (time3.tv_usec - time2.tv_usec);
  331. usec3 = (time4.tv_sec - time3.tv_sec) * 1000000 + (time4.tv_usec - time3.tv_usec);
  332. printf("Reference C implementation: %.2e M elements per second.\n",(double)rows*(double)cols/(double)usec2);
  333. #ifdef HAVE_AVX
  334. printf("AVX-enhanced implementation: %.2e M elements per second.\n",(double)rows*(double)cols/(double)usec1);
  335. #endif
  336. printf("4x4 blocked C implementation: %.2e M elements per second.\n",(double)rows*(double)cols/(double)usec3);
  337. maxdif=0.0;
  338. for(count=0;count<length;count++){
  339. tmpdif = B2[count] - B1[count];
  340. tmpdif2 = B3[count] - B1[count];
  341. if(tmpdif<0) tmpdif*=-1.0;
  342. if(tmpdif>maxdif) maxdif = tmpdif;
  343. if(tmpdif2<0) tmpdif2 *=-1.0;
  344. if(tmpdif2>maxdif2) maxdif2 = tmpdif2;
  345. }
  346. #ifdef HAVE_AVX
  347. printf("Maxdiff between output matrices: (AVX vx Ref %.2e\n",maxdif);
  348. #endif
  349. printf("Maxdiff between output matrices: (4x4 vx Ref %.2e\n",maxdif2);
  350. #ifdef HAVE_AVX
  351. CNAME(rows/2,cols/2,0.0,A,lda,B2,ldb);
  352. #endif
  353. REF(rows/2,cols/2,0.0,A,lda,B1,ldb);
  354. CNAME2(rows/2,cols/2,0.0,A,lda,B3,ldb);
  355. maxdif=0.0;
  356. maxdif2=0.0;
  357. for(count=0;count<length;count++){
  358. tmpdif = B2[count] - B1[count];
  359. tmpdif2 = B3[count] - B1[count];
  360. if(tmpdif<0) tmpdif*=-1.0;
  361. if(tmpdif>maxdif) maxdif = tmpdif;
  362. if(tmpdif2<0) tmpdif2 *=-1.0;
  363. if(tmpdif2>maxdif2) maxdif2 = tmpdif2;
  364. }
  365. #ifdef HAVE_AVX
  366. printf("Maxdiff between output matrices after subsequent alpha=zero runs: AVX vs Ref %.2e\n",maxdif);
  367. #endif
  368. printf("Maxdiff between output matrices after subsequent alpha=zero runs: 4x4 vs Ref %.2e\n",maxdif2);
  369. }
  370. free(A); free(B1); free(B2); free(B3); A=B1=B2=B3=NULL;
  371. }
  372. #endif