Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
AVX512/PacketMath.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2016 Benoit Steiner ([email protected])
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_PACKET_MATH_AVX512_H
11#define EIGEN_PACKET_MATH_AVX512_H
12
13// IWYU pragma: private
14#include "../../InternalHeaderCheck.h"
15
16namespace Eigen {
17
18namespace internal {
19
20#ifndef EIGEN_CACHEFRIENDLY_PRODUCT_THRESHOLD
21#define EIGEN_CACHEFRIENDLY_PRODUCT_THRESHOLD 8
22#endif
23
24#ifndef EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS
25#define EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS 32
26#endif
27
28#ifdef EIGEN_VECTORIZE_FMA
29#ifndef EIGEN_HAS_SINGLE_INSTRUCTION_MADD
30#define EIGEN_HAS_SINGLE_INSTRUCTION_MADD
31#endif
32#endif
33
34typedef __m512 Packet16f;
35typedef __m512i Packet16i;
36typedef __m512d Packet8d;
37typedef eigen_packet_wrapper<__m512i, 1> Packet8l;
38#ifndef EIGEN_VECTORIZE_AVX512FP16
39typedef eigen_packet_wrapper<__m256i, 1> Packet16h;
40#endif
41typedef eigen_packet_wrapper<__m256i, 2> Packet16bf;
42
43template <>
44struct is_arithmetic<__m512> {
45 enum { value = true };
46};
47template <>
48struct is_arithmetic<__m512i> {
49 enum { value = true };
50};
51template <>
52struct is_arithmetic<__m512d> {
53 enum { value = true };
54};
55template <>
56struct is_arithmetic<Packet8l> {
57 enum { value = true };
58};
59
60#ifndef EIGEN_VECTORIZE_AVX512FP16
61template <>
62struct is_arithmetic<Packet16h> {
63 enum { value = true };
64};
65
66template <>
67struct packet_traits<half> : default_packet_traits {
68 typedef Packet16h type;
69 // There is no half-size packet for Packet16h.
70 typedef Packet16h half;
71 enum {
72 Vectorizable = 1,
73 AlignedOnScalar = 1,
74 size = 16,
75
76 HasCmp = 1,
77 HasAdd = 1,
78 HasSub = 1,
79 HasMul = 1,
80 HasDiv = 1,
81 HasNegate = 1,
82 HasAbs = 1,
83 HasAbs2 = 0,
84 HasMin = 1,
85 HasMax = 1,
86 HasConj = 1,
87 HasSetLinear = 0,
88 HasSqrt = 1,
89 HasRsqrt = 1,
90 HasLog = 1,
91 HasLog1p = 1,
92 HasExp = 1,
93 HasExpm1 = 1,
94 HasBessel = 1,
95 HasNdtri = 1,
96 HasSin = EIGEN_FAST_MATH,
97 HasCos = EIGEN_FAST_MATH,
98 HasTanh = EIGEN_FAST_MATH,
99 HasErf = EIGEN_FAST_MATH,
100 HasBlend = 0
101 };
102};
103#endif
104
105template <>
106struct packet_traits<float> : default_packet_traits {
107 typedef Packet16f type;
108 typedef Packet8f half;
109 enum {
110 Vectorizable = 1,
111 AlignedOnScalar = 1,
112 size = 16,
113
114 HasAbs = 1,
115 HasMin = 1,
116 HasMax = 1,
117 HasConj = 1,
118 HasBlend = 1,
119 HasSin = EIGEN_FAST_MATH,
120 HasCos = EIGEN_FAST_MATH,
121 HasACos = 1,
122 HasASin = 1,
123 HasATan = 1,
124 HasATanh = 1,
125 HasSqrt = 1,
126 HasRsqrt = 1,
127 HasLog = 1,
128 HasLog1p = 1,
129 HasExpm1 = 1,
130 HasNdtri = 1,
131 HasBessel = 1,
132 HasExp = 1,
133 HasReciprocal = EIGEN_FAST_MATH,
134 HasTanh = EIGEN_FAST_MATH,
135 HasErf = EIGEN_FAST_MATH,
136 HasCmp = 1,
137 HasDiv = 1
138 };
139};
140template <>
141struct packet_traits<double> : default_packet_traits {
142 typedef Packet8d type;
143 typedef Packet4d half;
144 enum {
145 Vectorizable = 1,
146 AlignedOnScalar = 1,
147 size = 8,
148 HasBlend = 1,
149 HasSqrt = 1,
150 HasRsqrt = 1,
151 HasSin = EIGEN_FAST_MATH,
152 HasCos = EIGEN_FAST_MATH,
153 HasLog = 1,
154 HasExp = 1,
155 HasATan = 1,
156 HasCmp = 1,
157 HasDiv = 1
158 };
159};
160
161template <>
162struct packet_traits<int> : default_packet_traits {
163 typedef Packet16i type;
164 typedef Packet8i half;
165 enum { Vectorizable = 1, AlignedOnScalar = 1, HasBlend = 0, HasCmp = 1, HasDiv = 1, size = 16 };
166};
167
168template <>
169struct packet_traits<int64_t> : default_packet_traits {
170 typedef Packet8l type;
171 typedef Packet4l half;
172 enum { Vectorizable = 1, AlignedOnScalar = 1, HasCmp = 1, size = 8 };
173};
174
175template <>
176struct unpacket_traits<Packet16f> {
177 typedef float type;
178 typedef Packet8f half;
179 typedef Packet16i integer_packet;
180 typedef uint16_t mask_t;
181 enum {
182 size = 16,
183 alignment = Aligned64,
184 vectorizable = true,
185 masked_load_available = true,
186 masked_store_available = true,
187 masked_fpops_available = true
188 };
189};
190template <>
191struct unpacket_traits<Packet8d> {
192 typedef double type;
193 typedef Packet4d half;
194 typedef Packet8l integer_packet;
195 typedef uint8_t mask_t;
196 enum {
197 size = 8,
198 alignment = Aligned64,
199 vectorizable = true,
200 masked_load_available = true,
201 masked_store_available = true,
202 masked_fpops_available = true
203 };
204};
205template <>
206struct unpacket_traits<Packet16i> {
207 typedef int type;
208 typedef Packet8i half;
209 enum {
210 size = 16,
211 alignment = Aligned64,
212 vectorizable = true,
213 masked_load_available = false,
214 masked_store_available = false
215 };
216};
217
218template <>
219struct unpacket_traits<Packet8l> {
220 typedef int64_t type;
221 typedef Packet4l half;
222 enum {
223 size = 8,
224 alignment = Aligned64,
225 vectorizable = true,
226 masked_load_available = false,
227 masked_store_available = false
228 };
229};
230
231#ifndef EIGEN_VECTORIZE_AVX512FP16
232template <>
233struct unpacket_traits<Packet16h> {
234 typedef Eigen::half type;
235 typedef Packet8h half;
236 enum {
237 size = 16,
238 alignment = Aligned32,
239 vectorizable = true,
240 masked_load_available = false,
241 masked_store_available = false
242 };
243};
244#endif
245
246template <>
247EIGEN_STRONG_INLINE Packet16f pset1<Packet16f>(const float& from) {
248 return _mm512_set1_ps(from);
249}
250template <>
251EIGEN_STRONG_INLINE Packet8d pset1<Packet8d>(const double& from) {
252 return _mm512_set1_pd(from);
253}
254template <>
255EIGEN_STRONG_INLINE Packet16i pset1<Packet16i>(const int& from) {
256 return _mm512_set1_epi32(from);
257}
258template <>
259EIGEN_STRONG_INLINE Packet8l pset1<Packet8l>(const int64_t& from) {
260 return _mm512_set1_epi64(from);
261}
262
263template <>
264EIGEN_STRONG_INLINE Packet16f pset1frombits<Packet16f>(unsigned int from) {
265 return _mm512_castsi512_ps(_mm512_set1_epi32(from));
266}
267
268template <>
269EIGEN_STRONG_INLINE Packet8d pset1frombits<Packet8d>(const numext::uint64_t from) {
270 return _mm512_castsi512_pd(_mm512_set1_epi64(from));
271}
272
273template <>
274EIGEN_STRONG_INLINE Packet16f pzero(const Packet16f& /*a*/) {
275 return _mm512_setzero_ps();
276}
277template <>
278EIGEN_STRONG_INLINE Packet8d pzero(const Packet8d& /*a*/) {
279 return _mm512_setzero_pd();
280}
281template <>
282EIGEN_STRONG_INLINE Packet16i pzero(const Packet16i& /*a*/) {
283 return _mm512_setzero_si512();
284}
285
286template <>
287EIGEN_STRONG_INLINE Packet8l pzero(const Packet8l& /*a*/) {
288 return _mm512_setzero_si512();
289}
290
291template <>
292EIGEN_STRONG_INLINE Packet16f peven_mask(const Packet16f& /*a*/) {
293 return _mm512_castsi512_ps(_mm512_set_epi32(0, -1, 0, -1, 0, -1, 0, -1, 0, -1, 0, -1, 0, -1, 0, -1));
294}
295template <>
296EIGEN_STRONG_INLINE Packet16i peven_mask(const Packet16i& /*a*/) {
297 return _mm512_set_epi32(0, -1, 0, -1, 0, -1, 0, -1, 0, -1, 0, -1, 0, -1, 0, -1);
298}
299template <>
300EIGEN_STRONG_INLINE Packet8d peven_mask(const Packet8d& /*a*/) {
301 return _mm512_castsi512_pd(_mm512_set_epi32(0, 0, -1, -1, 0, 0, -1, -1, 0, 0, -1, -1, 0, 0, -1, -1));
302}
303template <>
304EIGEN_STRONG_INLINE Packet8l peven_mask(const Packet8l& /*a*/) {
305 return _mm512_set_epi32(0, 0, -1, -1, 0, 0, -1, -1, 0, 0, -1, -1, 0, 0, -1, -1);
306}
307
308template <>
309EIGEN_STRONG_INLINE Packet16f pload1<Packet16f>(const float* from) {
310#if (EIGEN_COMP_GNUC != 0) || (EIGEN_COMP_CLANG != 0)
311 // Inline asm here helps reduce some register spilling in TRSM kernels.
312 // See note in unrolls::gemm::microKernel in TrsmKernel.h
313 Packet16f ret;
314 __asm__("vbroadcastss %[mem], %[dst]" : [dst] "=v"(ret) : [mem] "m"(*from));
315 return ret;
316#else
317 return _mm512_broadcastss_ps(_mm_load_ps1(from));
318#endif
319}
320template <>
321EIGEN_STRONG_INLINE Packet8d pload1<Packet8d>(const double* from) {
322#if (EIGEN_COMP_GNUC != 0) || (EIGEN_COMP_CLANG != 0)
323 Packet8d ret;
324 __asm__("vbroadcastsd %[mem], %[dst]" : [dst] "=v"(ret) : [mem] "m"(*from));
325 return ret;
326#else
327 return _mm512_set1_pd(*from);
328#endif
329}
330
331template <>
332EIGEN_STRONG_INLINE Packet16f plset<Packet16f>(const float& a) {
333 return _mm512_add_ps(_mm512_set1_ps(a), _mm512_set_ps(15.0f, 14.0f, 13.0f, 12.0f, 11.0f, 10.0f, 9.0f, 8.0f, 7.0f,
334 6.0f, 5.0f, 4.0f, 3.0f, 2.0f, 1.0f, 0.0f));
335}
336template <>
337EIGEN_STRONG_INLINE Packet8d plset<Packet8d>(const double& a) {
338 return _mm512_add_pd(_mm512_set1_pd(a), _mm512_set_pd(7.0, 6.0, 5.0, 4.0, 3.0, 2.0, 1.0, 0.0));
339}
340template <>
341EIGEN_STRONG_INLINE Packet16i plset<Packet16i>(const int& a) {
342 return _mm512_add_epi32(_mm512_set1_epi32(a), _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0));
343}
344template <>
345EIGEN_STRONG_INLINE Packet8l plset<Packet8l>(const int64_t& a) {
346 return _mm512_add_epi64(_mm512_set1_epi64(a), _mm512_set_epi64(7, 6, 5, 4, 3, 2, 1, 0));
347}
348
349template <>
350EIGEN_STRONG_INLINE Packet16f padd<Packet16f>(const Packet16f& a, const Packet16f& b) {
351 return _mm512_add_ps(a, b);
352}
353template <>
354EIGEN_STRONG_INLINE Packet8d padd<Packet8d>(const Packet8d& a, const Packet8d& b) {
355 return _mm512_add_pd(a, b);
356}
357template <>
358EIGEN_STRONG_INLINE Packet16i padd<Packet16i>(const Packet16i& a, const Packet16i& b) {
359 return _mm512_add_epi32(a, b);
360}
361template <>
362EIGEN_STRONG_INLINE Packet8l padd<Packet8l>(const Packet8l& a, const Packet8l& b) {
363 return _mm512_add_epi64(a, b);
364}
365
366template <>
367EIGEN_STRONG_INLINE Packet16f padd<Packet16f>(const Packet16f& a, const Packet16f& b, uint16_t umask) {
368 __mmask16 mask = static_cast<__mmask16>(umask);
369 return _mm512_maskz_add_ps(mask, a, b);
370}
371template <>
372EIGEN_STRONG_INLINE Packet8d padd<Packet8d>(const Packet8d& a, const Packet8d& b, uint8_t umask) {
373 __mmask8 mask = static_cast<__mmask8>(umask);
374 return _mm512_maskz_add_pd(mask, a, b);
375}
376
377template <>
378EIGEN_STRONG_INLINE Packet16f psub<Packet16f>(const Packet16f& a, const Packet16f& b) {
379 return _mm512_sub_ps(a, b);
380}
381template <>
382EIGEN_STRONG_INLINE Packet8d psub<Packet8d>(const Packet8d& a, const Packet8d& b) {
383 return _mm512_sub_pd(a, b);
384}
385template <>
386EIGEN_STRONG_INLINE Packet16i psub<Packet16i>(const Packet16i& a, const Packet16i& b) {
387 return _mm512_sub_epi32(a, b);
388}
389template <>
390EIGEN_STRONG_INLINE Packet8l psub<Packet8l>(const Packet8l& a, const Packet8l& b) {
391 return _mm512_sub_epi64(a, b);
392}
393
394template <>
395EIGEN_STRONG_INLINE Packet16f pnegate(const Packet16f& a) {
396 // NOTE: MSVC seems to struggle with _mm512_set1_epi32, leading to random results.
397 // The intel docs give it a relatively high latency as well, so we're probably
398 // better off with using _mm512_set_epi32 directly anyways.
399 const __m512i mask =
400 _mm512_set_epi32(0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000,
401 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000, 0x80000000);
402 return _mm512_castsi512_ps(_mm512_xor_epi32(_mm512_castps_si512(a), mask));
403}
404template <>
405EIGEN_STRONG_INLINE Packet8d pnegate(const Packet8d& a) {
406 const __m512i mask =
407 _mm512_set_epi64(0x8000000000000000ULL, 0x8000000000000000ULL, 0x8000000000000000ULL, 0x8000000000000000ULL,
408 0x8000000000000000ULL, 0x8000000000000000ULL, 0x8000000000000000ULL, 0x8000000000000000ULL);
409 return _mm512_castsi512_pd(_mm512_xor_epi64(_mm512_castpd_si512(a), mask));
410}
411template <>
412EIGEN_STRONG_INLINE Packet16i pnegate(const Packet16i& a) {
413 return _mm512_sub_epi32(_mm512_setzero_si512(), a);
414}
415template <>
416EIGEN_STRONG_INLINE Packet8l pnegate(const Packet8l& a) {
417 return _mm512_sub_epi64(_mm512_setzero_si512(), a);
418}
419
420template <>
421EIGEN_STRONG_INLINE Packet16f pconj(const Packet16f& a) {
422 return a;
423}
424template <>
425EIGEN_STRONG_INLINE Packet8d pconj(const Packet8d& a) {
426 return a;
427}
428template <>
429EIGEN_STRONG_INLINE Packet16i pconj(const Packet16i& a) {
430 return a;
431}
432template <>
433EIGEN_STRONG_INLINE Packet8l pconj(const Packet8l& a) {
434 return a;
435}
436
437template <>
438EIGEN_STRONG_INLINE Packet16f pmul<Packet16f>(const Packet16f& a, const Packet16f& b) {
439 return _mm512_mul_ps(a, b);
440}
441template <>
442EIGEN_STRONG_INLINE Packet8d pmul<Packet8d>(const Packet8d& a, const Packet8d& b) {
443 return _mm512_mul_pd(a, b);
444}
445template <>
446EIGEN_STRONG_INLINE Packet16i pmul<Packet16i>(const Packet16i& a, const Packet16i& b) {
447 return _mm512_mullo_epi32(a, b);
448}
449template <>
450EIGEN_STRONG_INLINE Packet8l pmul<Packet8l>(const Packet8l& a, const Packet8l& b) {
451#ifdef EIGEN_VECTORIZE_AVX512DQ
452 return _mm512_mullo_epi64(a, b);
453#else
454 return _mm512_mullox_epi64(a, b);
455#endif
456}
457
458template <>
459EIGEN_STRONG_INLINE Packet16f pdiv<Packet16f>(const Packet16f& a, const Packet16f& b) {
460 return _mm512_div_ps(a, b);
461}
462
463template <>
464EIGEN_STRONG_INLINE Packet8d pdiv<Packet8d>(const Packet8d& a, const Packet8d& b) {
465 return _mm512_div_pd(a, b);
466}
467
468template <>
469EIGEN_STRONG_INLINE Packet16i pdiv<Packet16i>(const Packet16i& a, const Packet16i& b) {
470 Packet8i q_lo = pdiv<Packet8i>(_mm512_extracti64x4_epi64(a, 0), _mm512_extracti64x4_epi64(b, 0));
471 Packet8i q_hi = pdiv<Packet8i>(_mm512_extracti64x4_epi64(a, 1), _mm512_extracti64x4_epi64(b, 1));
472 return _mm512_inserti64x4(_mm512_castsi256_si512(q_lo), q_hi, 1);
473}
474
475#ifdef EIGEN_VECTORIZE_FMA
476template <>
477EIGEN_STRONG_INLINE Packet16f pmadd(const Packet16f& a, const Packet16f& b, const Packet16f& c) {
478 return _mm512_fmadd_ps(a, b, c);
479}
480template <>
481EIGEN_STRONG_INLINE Packet8d pmadd(const Packet8d& a, const Packet8d& b, const Packet8d& c) {
482 return _mm512_fmadd_pd(a, b, c);
483}
484
485template <>
486EIGEN_STRONG_INLINE Packet16f pmsub(const Packet16f& a, const Packet16f& b, const Packet16f& c) {
487 return _mm512_fmsub_ps(a, b, c);
488}
489template <>
490EIGEN_STRONG_INLINE Packet8d pmsub(const Packet8d& a, const Packet8d& b, const Packet8d& c) {
491 return _mm512_fmsub_pd(a, b, c);
492}
493
494template <>
495EIGEN_STRONG_INLINE Packet16f pnmadd(const Packet16f& a, const Packet16f& b, const Packet16f& c) {
496 return _mm512_fnmadd_ps(a, b, c);
497}
498template <>
499EIGEN_STRONG_INLINE Packet8d pnmadd(const Packet8d& a, const Packet8d& b, const Packet8d& c) {
500 return _mm512_fnmadd_pd(a, b, c);
501}
502
503template <>
504EIGEN_STRONG_INLINE Packet16f pnmsub(const Packet16f& a, const Packet16f& b, const Packet16f& c) {
505 return _mm512_fnmsub_ps(a, b, c);
506}
507template <>
508EIGEN_STRONG_INLINE Packet8d pnmsub(const Packet8d& a, const Packet8d& b, const Packet8d& c) {
509 return _mm512_fnmsub_pd(a, b, c);
510}
511#endif
512
513template <>
514EIGEN_DEVICE_FUNC inline Packet16f pselect(const Packet16f& mask, const Packet16f& a, const Packet16f& b) {
515 __mmask16 mask16 = _mm512_cmpeq_epi32_mask(_mm512_castps_si512(mask), _mm512_setzero_epi32());
516 return _mm512_mask_blend_ps(mask16, a, b);
517}
518
519template <>
520EIGEN_DEVICE_FUNC inline Packet16i pselect(const Packet16i& mask, const Packet16i& a, const Packet16i& b) {
521 __mmask16 mask16 = _mm512_cmpeq_epi32_mask(mask, _mm512_setzero_epi32());
522 return _mm512_mask_blend_epi32(mask16, a, b);
523}
524
525template <>
526EIGEN_DEVICE_FUNC inline Packet8l pselect(const Packet8l& mask, const Packet8l& a, const Packet8l& b) {
527 __mmask8 mask8 = _mm512_cmpeq_epi64_mask(mask, _mm512_setzero_si512());
528 return _mm512_mask_blend_epi64(mask8, a, b);
529}
530
531template <>
532EIGEN_DEVICE_FUNC inline Packet8d pselect(const Packet8d& mask, const Packet8d& a, const Packet8d& b) {
533 __mmask8 mask8 = _mm512_cmp_epi64_mask(_mm512_castpd_si512(mask), _mm512_setzero_epi32(), _MM_CMPINT_EQ);
534 return _mm512_mask_blend_pd(mask8, a, b);
535}
536
537template <>
538EIGEN_STRONG_INLINE Packet16f pmin<Packet16f>(const Packet16f& a, const Packet16f& b) {
539 // Arguments are reversed to match NaN propagation behavior of std::min.
540 return _mm512_min_ps(b, a);
541}
542template <>
543EIGEN_STRONG_INLINE Packet8d pmin<Packet8d>(const Packet8d& a, const Packet8d& b) {
544 // Arguments are reversed to match NaN propagation behavior of std::min.
545 return _mm512_min_pd(b, a);
546}
547template <>
548EIGEN_STRONG_INLINE Packet16i pmin<Packet16i>(const Packet16i& a, const Packet16i& b) {
549 return _mm512_min_epi32(b, a);
550}
551template <>
552EIGEN_STRONG_INLINE Packet8l pmin<Packet8l>(const Packet8l& a, const Packet8l& b) {
553 return _mm512_min_epi64(b, a);
554}
555
556template <>
557EIGEN_STRONG_INLINE Packet16f pmax<Packet16f>(const Packet16f& a, const Packet16f& b) {
558 // Arguments are reversed to match NaN propagation behavior of std::max.
559 return _mm512_max_ps(b, a);
560}
561template <>
562EIGEN_STRONG_INLINE Packet8d pmax<Packet8d>(const Packet8d& a, const Packet8d& b) {
563 // Arguments are reversed to match NaN propagation behavior of std::max.
564 return _mm512_max_pd(b, a);
565}
566template <>
567EIGEN_STRONG_INLINE Packet16i pmax<Packet16i>(const Packet16i& a, const Packet16i& b) {
568 return _mm512_max_epi32(b, a);
569}
570template <>
571EIGEN_STRONG_INLINE Packet8l pmax<Packet8l>(const Packet8l& a, const Packet8l& b) {
572 return _mm512_max_epi64(b, a);
573}
574
575// Add specializations for min/max with prescribed NaN progation.
576template <>
577EIGEN_STRONG_INLINE Packet16f pmin<PropagateNumbers, Packet16f>(const Packet16f& a, const Packet16f& b) {
578 return pminmax_propagate_numbers(a, b, pmin<Packet16f>);
579}
580template <>
581EIGEN_STRONG_INLINE Packet8d pmin<PropagateNumbers, Packet8d>(const Packet8d& a, const Packet8d& b) {
582 return pminmax_propagate_numbers(a, b, pmin<Packet8d>);
583}
584template <>
585EIGEN_STRONG_INLINE Packet16f pmax<PropagateNumbers, Packet16f>(const Packet16f& a, const Packet16f& b) {
586 return pminmax_propagate_numbers(a, b, pmax<Packet16f>);
587}
588template <>
589EIGEN_STRONG_INLINE Packet8d pmax<PropagateNumbers, Packet8d>(const Packet8d& a, const Packet8d& b) {
590 return pminmax_propagate_numbers(a, b, pmax<Packet8d>);
591}
592template <>
593EIGEN_STRONG_INLINE Packet16f pmin<PropagateNaN, Packet16f>(const Packet16f& a, const Packet16f& b) {
594 return pminmax_propagate_nan(a, b, pmin<Packet16f>);
595}
596template <>
597EIGEN_STRONG_INLINE Packet8d pmin<PropagateNaN, Packet8d>(const Packet8d& a, const Packet8d& b) {
598 return pminmax_propagate_nan(a, b, pmin<Packet8d>);
599}
600template <>
601EIGEN_STRONG_INLINE Packet16f pmax<PropagateNaN, Packet16f>(const Packet16f& a, const Packet16f& b) {
602 return pminmax_propagate_nan(a, b, pmax<Packet16f>);
603}
604template <>
605EIGEN_STRONG_INLINE Packet8d pmax<PropagateNaN, Packet8d>(const Packet8d& a, const Packet8d& b) {
606 return pminmax_propagate_nan(a, b, pmax<Packet8d>);
607}
608
609#ifdef EIGEN_VECTORIZE_AVX512DQ
610template <int I_>
611EIGEN_STRONG_INLINE Packet8f extract256(Packet16f x) {
612 return _mm512_extractf32x8_ps(x, I_);
613}
614template <int I_>
615EIGEN_STRONG_INLINE Packet2d extract128(Packet8d x) {
616 return _mm512_extractf64x2_pd(x, I_);
617}
618EIGEN_STRONG_INLINE Packet16f cat256(Packet8f a, Packet8f b) {
619 return _mm512_insertf32x8(_mm512_castps256_ps512(a), b, 1);
620}
621EIGEN_STRONG_INLINE Packet16i cat256i(Packet8i a, Packet8i b) {
622 return _mm512_inserti32x8(_mm512_castsi256_si512(a), b, 1);
623}
624#else
625// AVX512F does not define _mm512_extractf32x8_ps to extract _m256 from _m512
626template <int I_>
627EIGEN_STRONG_INLINE Packet8f extract256(Packet16f x) {
628 return _mm256_castsi256_ps(_mm512_extracti64x4_epi64(_mm512_castps_si512(x), I_));
629}
630
631// AVX512F does not define _mm512_extractf64x2_pd to extract _m128 from _m512
632template <int I_>
633EIGEN_STRONG_INLINE Packet2d extract128(Packet8d x) {
634 return _mm_castsi128_pd(_mm512_extracti32x4_epi32(_mm512_castpd_si512(x), I_));
635}
636
637EIGEN_STRONG_INLINE Packet16f cat256(Packet8f a, Packet8f b) {
638 return _mm512_castsi512_ps(
639 _mm512_inserti64x4(_mm512_castsi256_si512(_mm256_castps_si256(a)), _mm256_castps_si256(b), 1));
640}
641EIGEN_STRONG_INLINE Packet16i cat256i(Packet8i a, Packet8i b) {
642 return _mm512_inserti64x4(_mm512_castsi256_si512(a), b, 1);
643}
644#endif
645
646// Helper function for bit packing snippet of low precision comparison.
647// It packs the flags from 32x16 to 16x16.
648EIGEN_STRONG_INLINE __m256i Pack32To16(Packet16f rf) {
649 // Split data into small pieces and handle with AVX instructions
650 // to guarantee internal order of vector.
651 // Operation:
652 // dst[15:0] := Saturate16(rf[31:0])
653 // dst[31:16] := Saturate16(rf[63:32])
654 // ...
655 // dst[255:240] := Saturate16(rf[255:224])
656 __m256i lo = _mm256_castps_si256(extract256<0>(rf));
657 __m256i hi = _mm256_castps_si256(extract256<1>(rf));
658 __m128i result_lo = _mm_packs_epi32(_mm256_extractf128_si256(lo, 0), _mm256_extractf128_si256(lo, 1));
659 __m128i result_hi = _mm_packs_epi32(_mm256_extractf128_si256(hi, 0), _mm256_extractf128_si256(hi, 1));
660 return _mm256_insertf128_si256(_mm256_castsi128_si256(result_lo), result_hi, 1);
661}
662
663template <>
664EIGEN_STRONG_INLINE Packet16f pisnan(const Packet16f& a) {
665 __mmask16 mask = _mm512_cmp_ps_mask(a, a, _CMP_UNORD_Q);
666 return _mm512_castsi512_ps(_mm512_maskz_set1_epi32(mask, int32_t(-1)));
667}
668
669template <>
670EIGEN_STRONG_INLINE Packet16f pcmp_eq(const Packet16f& a, const Packet16f& b) {
671 __mmask16 mask = _mm512_cmp_ps_mask(a, b, _CMP_EQ_OQ);
672 return _mm512_castsi512_ps(_mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1)));
673}
674template <>
675EIGEN_STRONG_INLINE Packet16f pcmp_le(const Packet16f& a, const Packet16f& b) {
676 __mmask16 mask = _mm512_cmp_ps_mask(a, b, _CMP_LE_OQ);
677 return _mm512_castsi512_ps(_mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1)));
678}
679
680template <>
681EIGEN_STRONG_INLINE Packet16f pcmp_lt(const Packet16f& a, const Packet16f& b) {
682 __mmask16 mask = _mm512_cmp_ps_mask(a, b, _CMP_LT_OQ);
683 return _mm512_castsi512_ps(_mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1)));
684}
685
686template <>
687EIGEN_STRONG_INLINE Packet16f pcmp_lt_or_nan(const Packet16f& a, const Packet16f& b) {
688 __mmask16 mask = _mm512_cmp_ps_mask(a, b, _CMP_NGE_UQ);
689 return _mm512_castsi512_ps(_mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1)));
690}
691
692template <>
693EIGEN_STRONG_INLINE Packet16i pcmp_eq(const Packet16i& a, const Packet16i& b) {
694 __mmask16 mask = _mm512_cmp_epi32_mask(a, b, _MM_CMPINT_EQ);
695 return _mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1));
696}
697template <>
698EIGEN_STRONG_INLINE Packet16i pcmp_le(const Packet16i& a, const Packet16i& b) {
699 __mmask16 mask = _mm512_cmp_epi32_mask(a, b, _MM_CMPINT_LE);
700 return _mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1));
701}
702template <>
703EIGEN_STRONG_INLINE Packet16i pcmp_lt(const Packet16i& a, const Packet16i& b) {
704 __mmask16 mask = _mm512_cmp_epi32_mask(a, b, _MM_CMPINT_LT);
705 return _mm512_mask_set1_epi32(_mm512_setzero_epi32(), mask, int32_t(-1));
706}
707
708template <>
709EIGEN_STRONG_INLINE Packet8l pcmp_eq(const Packet8l& a, const Packet8l& b) {
710 __mmask8 mask = _mm512_cmp_epi64_mask(a, b, _MM_CMPINT_EQ);
711 return _mm512_mask_set1_epi64(_mm512_setzero_si512(), mask, int64_t(-1));
712}
713template <>
714EIGEN_STRONG_INLINE Packet8l pcmp_le(const Packet8l& a, const Packet8l& b) {
715 __mmask8 mask = _mm512_cmp_epi64_mask(a, b, _MM_CMPINT_LE);
716 return _mm512_mask_set1_epi64(_mm512_setzero_si512(), mask, int64_t(-1));
717}
718template <>
719EIGEN_STRONG_INLINE Packet8l pcmp_lt(const Packet8l& a, const Packet8l& b) {
720 __mmask8 mask = _mm512_cmp_epi64_mask(a, b, _MM_CMPINT_LT);
721 return _mm512_mask_set1_epi64(_mm512_setzero_si512(), mask, int64_t(-1));
722}
723
724template <>
725EIGEN_STRONG_INLINE Packet8d pcmp_eq(const Packet8d& a, const Packet8d& b) {
726 __mmask8 mask = _mm512_cmp_pd_mask(a, b, _CMP_EQ_OQ);
727 return _mm512_castsi512_pd(_mm512_mask_set1_epi64(_mm512_setzero_epi32(), mask, 0xffffffffffffffffu));
728}
729template <>
730EIGEN_STRONG_INLINE Packet8d pcmp_le(const Packet8d& a, const Packet8d& b) {
731 __mmask8 mask = _mm512_cmp_pd_mask(a, b, _CMP_LE_OQ);
732 return _mm512_castsi512_pd(_mm512_mask_set1_epi64(_mm512_setzero_epi32(), mask, 0xffffffffffffffffu));
733}
734template <>
735EIGEN_STRONG_INLINE Packet8d pcmp_lt(const Packet8d& a, const Packet8d& b) {
736 __mmask8 mask = _mm512_cmp_pd_mask(a, b, _CMP_LT_OQ);
737 return _mm512_castsi512_pd(_mm512_mask_set1_epi64(_mm512_setzero_epi32(), mask, 0xffffffffffffffffu));
738}
739template <>
740EIGEN_STRONG_INLINE Packet8d pcmp_lt_or_nan(const Packet8d& a, const Packet8d& b) {
741 __mmask8 mask = _mm512_cmp_pd_mask(a, b, _CMP_NGE_UQ);
742 return _mm512_castsi512_pd(_mm512_mask_set1_epi64(_mm512_setzero_epi32(), mask, 0xffffffffffffffffu));
743}
744
745template <>
746EIGEN_STRONG_INLINE Packet16f print<Packet16f>(const Packet16f& a) {
747 return _mm512_roundscale_ps(a, _MM_FROUND_CUR_DIRECTION);
748}
749template <>
750EIGEN_STRONG_INLINE Packet8d print<Packet8d>(const Packet8d& a) {
751 return _mm512_roundscale_pd(a, _MM_FROUND_CUR_DIRECTION);
752}
753
754template <>
755EIGEN_STRONG_INLINE Packet16f pceil<Packet16f>(const Packet16f& a) {
756 return _mm512_roundscale_ps(a, _MM_FROUND_TO_POS_INF);
757}
758template <>
759EIGEN_STRONG_INLINE Packet8d pceil<Packet8d>(const Packet8d& a) {
760 return _mm512_roundscale_pd(a, _MM_FROUND_TO_POS_INF);
761}
762
763template <>
764EIGEN_STRONG_INLINE Packet16f pfloor<Packet16f>(const Packet16f& a) {
765 return _mm512_roundscale_ps(a, _MM_FROUND_TO_NEG_INF);
766}
767template <>
768EIGEN_STRONG_INLINE Packet8d pfloor<Packet8d>(const Packet8d& a) {
769 return _mm512_roundscale_pd(a, _MM_FROUND_TO_NEG_INF);
770}
771
772template <>
773EIGEN_STRONG_INLINE Packet16f ptrunc<Packet16f>(const Packet16f& a) {
774 return _mm512_roundscale_ps(a, _MM_FROUND_TO_ZERO);
775}
776template <>
777EIGEN_STRONG_INLINE Packet8d ptrunc<Packet8d>(const Packet8d& a) {
778 return _mm512_roundscale_pd(a, _MM_FROUND_TO_ZERO);
779}
780
781template <>
782EIGEN_STRONG_INLINE Packet16i ptrue<Packet16i>(const Packet16i& /*a*/) {
783 return _mm512_set1_epi32(int32_t(-1));
784}
785
786template <>
787EIGEN_STRONG_INLINE Packet8l ptrue<Packet8l>(const Packet8l& /*a*/) {
788 return _mm512_set1_epi64(int64_t(-1));
789}
790
791template <>
792EIGEN_STRONG_INLINE Packet16f ptrue<Packet16f>(const Packet16f& a) {
793 return _mm512_castsi512_ps(ptrue<Packet16i>(_mm512_castps_si512(a)));
794}
795
796template <>
797EIGEN_STRONG_INLINE Packet8d ptrue<Packet8d>(const Packet8d& a) {
798 return _mm512_castsi512_pd(ptrue<Packet16i>(_mm512_castpd_si512(a)));
799}
800
801template <>
802EIGEN_STRONG_INLINE Packet16i pand<Packet16i>(const Packet16i& a, const Packet16i& b) {
803 return _mm512_and_si512(a, b);
804}
805
806template <>
807EIGEN_STRONG_INLINE Packet8l pand<Packet8l>(const Packet8l& a, const Packet8l& b) {
808 return _mm512_and_si512(a, b);
809}
810
811template <>
812EIGEN_STRONG_INLINE Packet16f pand<Packet16f>(const Packet16f& a, const Packet16f& b) {
813#ifdef EIGEN_VECTORIZE_AVX512DQ
814 return _mm512_and_ps(a, b);
815#else
816 return _mm512_castsi512_ps(pand(_mm512_castps_si512(a), _mm512_castps_si512(b)));
817#endif
818}
819template <>
820EIGEN_STRONG_INLINE Packet8d pand<Packet8d>(const Packet8d& a, const Packet8d& b) {
821#ifdef EIGEN_VECTORIZE_AVX512DQ
822 return _mm512_and_pd(a, b);
823#else
824 Packet8d res = _mm512_undefined_pd();
825 Packet4d lane0_a = _mm512_extractf64x4_pd(a, 0);
826 Packet4d lane0_b = _mm512_extractf64x4_pd(b, 0);
827 res = _mm512_insertf64x4(res, _mm256_and_pd(lane0_a, lane0_b), 0);
828
829 Packet4d lane1_a = _mm512_extractf64x4_pd(a, 1);
830 Packet4d lane1_b = _mm512_extractf64x4_pd(b, 1);
831 return _mm512_insertf64x4(res, _mm256_and_pd(lane1_a, lane1_b), 1);
832#endif
833}
834
835template <>
836EIGEN_STRONG_INLINE Packet16i por<Packet16i>(const Packet16i& a, const Packet16i& b) {
837 return _mm512_or_si512(a, b);
838}
839
840template <>
841EIGEN_STRONG_INLINE Packet8l por<Packet8l>(const Packet8l& a, const Packet8l& b) {
842 return _mm512_or_si512(a, b);
843}
844
845template <>
846EIGEN_STRONG_INLINE Packet16f por<Packet16f>(const Packet16f& a, const Packet16f& b) {
847#ifdef EIGEN_VECTORIZE_AVX512DQ
848 return _mm512_or_ps(a, b);
849#else
850 return _mm512_castsi512_ps(por(_mm512_castps_si512(a), _mm512_castps_si512(b)));
851#endif
852}
853
854template <>
855EIGEN_STRONG_INLINE Packet8d por<Packet8d>(const Packet8d& a, const Packet8d& b) {
856#ifdef EIGEN_VECTORIZE_AVX512DQ
857 return _mm512_or_pd(a, b);
858#else
859 return _mm512_castsi512_pd(por(_mm512_castpd_si512(a), _mm512_castpd_si512(b)));
860#endif
861}
862
863template <>
864EIGEN_STRONG_INLINE Packet16i pxor<Packet16i>(const Packet16i& a, const Packet16i& b) {
865 return _mm512_xor_si512(a, b);
866}
867
868template <>
869EIGEN_STRONG_INLINE Packet8l pxor<Packet8l>(const Packet8l& a, const Packet8l& b) {
870 return _mm512_xor_si512(a, b);
871}
872
873template <>
874EIGEN_STRONG_INLINE Packet16f pxor<Packet16f>(const Packet16f& a, const Packet16f& b) {
875#ifdef EIGEN_VECTORIZE_AVX512DQ
876 return _mm512_xor_ps(a, b);
877#else
878 return _mm512_castsi512_ps(pxor(_mm512_castps_si512(a), _mm512_castps_si512(b)));
879#endif
880}
881
882template <>
883EIGEN_STRONG_INLINE Packet8d pxor<Packet8d>(const Packet8d& a, const Packet8d& b) {
884#ifdef EIGEN_VECTORIZE_AVX512DQ
885 return _mm512_xor_pd(a, b);
886#else
887 return _mm512_castsi512_pd(pxor(_mm512_castpd_si512(a), _mm512_castpd_si512(b)));
888#endif
889}
890
891template <>
892EIGEN_STRONG_INLINE Packet16i pandnot<Packet16i>(const Packet16i& a, const Packet16i& b) {
893 return _mm512_andnot_si512(b, a);
894}
895
896template <>
897EIGEN_STRONG_INLINE Packet8l pandnot<Packet8l>(const Packet8l& a, const Packet8l& b) {
898 return _mm512_andnot_si512(b, a);
899}
900
901template <>
902EIGEN_STRONG_INLINE Packet16f pandnot<Packet16f>(const Packet16f& a, const Packet16f& b) {
903#ifdef EIGEN_VECTORIZE_AVX512DQ
904 return _mm512_andnot_ps(b, a);
905#else
906 return _mm512_castsi512_ps(pandnot(_mm512_castps_si512(a), _mm512_castps_si512(b)));
907#endif
908}
909template <>
910EIGEN_STRONG_INLINE Packet8d pandnot<Packet8d>(const Packet8d& a, const Packet8d& b) {
911#ifdef EIGEN_VECTORIZE_AVX512DQ
912 return _mm512_andnot_pd(b, a);
913#else
914 return _mm512_castsi512_pd(pandnot(_mm512_castpd_si512(a), _mm512_castpd_si512(b)));
915#endif
916}
917
918template <>
919EIGEN_STRONG_INLINE Packet16f pround<Packet16f>(const Packet16f& a) {
920 // Work-around for default std::round rounding mode.
921 const Packet16f mask = pset1frombits<Packet16f>(static_cast<numext::uint32_t>(0x80000000u));
922 const Packet16f prev0dot5 = pset1frombits<Packet16f>(static_cast<numext::uint32_t>(0x3EFFFFFFu));
923 return _mm512_roundscale_ps(padd(por(pand(a, mask), prev0dot5), a), _MM_FROUND_TO_ZERO);
924}
925template <>
926EIGEN_STRONG_INLINE Packet8d pround<Packet8d>(const Packet8d& a) {
927 // Work-around for default std::round rounding mode.
928 const Packet8d mask = pset1frombits<Packet8d>(static_cast<numext::uint64_t>(0x8000000000000000ull));
929 const Packet8d prev0dot5 = pset1frombits<Packet8d>(static_cast<numext::uint64_t>(0x3FDFFFFFFFFFFFFFull));
930 return _mm512_roundscale_pd(padd(por(pand(a, mask), prev0dot5), a), _MM_FROUND_TO_ZERO);
931}
932
933template <int N>
934EIGEN_STRONG_INLINE Packet16i parithmetic_shift_right(Packet16i a) {
935 return _mm512_srai_epi32(a, N);
936}
937
938template <int N>
939EIGEN_STRONG_INLINE Packet16i plogical_shift_right(Packet16i a) {
940 return _mm512_srli_epi32(a, N);
941}
942
943template <int N>
944EIGEN_STRONG_INLINE Packet16i plogical_shift_left(Packet16i a) {
945 return _mm512_slli_epi32(a, N);
946}
947
948template <int N>
949EIGEN_STRONG_INLINE Packet8l parithmetic_shift_right(Packet8l a) {
950 return _mm512_srai_epi64(a, N);
951}
952
953template <int N>
954EIGEN_STRONG_INLINE Packet8l plogical_shift_right(Packet8l a) {
955 return _mm512_srli_epi64(a, N);
956}
957
958template <int N>
959EIGEN_STRONG_INLINE Packet8l plogical_shift_left(Packet8l a) {
960 return _mm512_slli_epi64(a, N);
961}
962
963template <>
964EIGEN_STRONG_INLINE Packet16f pload<Packet16f>(const float* from) {
965 EIGEN_DEBUG_ALIGNED_LOAD return _mm512_load_ps(from);
966}
967template <>
968EIGEN_STRONG_INLINE Packet8d pload<Packet8d>(const double* from) {
969 EIGEN_DEBUG_ALIGNED_LOAD return _mm512_load_pd(from);
970}
971template <>
972EIGEN_STRONG_INLINE Packet16i pload<Packet16i>(const int* from) {
973 EIGEN_DEBUG_ALIGNED_LOAD return _mm512_load_epi64(from);
974}
975template <>
976EIGEN_STRONG_INLINE Packet8l pload<Packet8l>(const int64_t* from) {
977 EIGEN_DEBUG_ALIGNED_LOAD return _mm512_load_epi64(from);
978}
979
980template <>
981EIGEN_STRONG_INLINE Packet16f ploadu<Packet16f>(const float* from) {
982 EIGEN_DEBUG_UNALIGNED_LOAD return _mm512_loadu_ps(from);
983}
984template <>
985EIGEN_STRONG_INLINE Packet8d ploadu<Packet8d>(const double* from) {
986 EIGEN_DEBUG_UNALIGNED_LOAD return _mm512_loadu_pd(from);
987}
988template <>
989EIGEN_STRONG_INLINE Packet16i ploadu<Packet16i>(const int* from) {
990 EIGEN_DEBUG_UNALIGNED_LOAD return _mm512_loadu_epi32(from);
991}
992template <>
993EIGEN_STRONG_INLINE Packet8l ploadu<Packet8l>(const int64_t* from) {
994 EIGEN_DEBUG_UNALIGNED_LOAD return _mm512_loadu_epi64(from);
995}
996
997template <>
998EIGEN_STRONG_INLINE Packet16f ploadu<Packet16f>(const float* from, uint16_t umask) {
999 __mmask16 mask = static_cast<__mmask16>(umask);
1000 EIGEN_DEBUG_UNALIGNED_LOAD return _mm512_maskz_loadu_ps(mask, from);
1001}
1002template <>
1003EIGEN_STRONG_INLINE Packet8d ploadu<Packet8d>(const double* from, uint8_t umask) {
1004 __mmask8 mask = static_cast<__mmask8>(umask);
1005 EIGEN_DEBUG_UNALIGNED_LOAD return _mm512_maskz_loadu_pd(mask, from);
1006}
1007
1008// Loads 8 floats from memory a returns the packet
1009// {a0, a0 a1, a1, a2, a2, a3, a3, a4, a4, a5, a5, a6, a6, a7, a7}
1010template <>
1011EIGEN_STRONG_INLINE Packet16f ploaddup<Packet16f>(const float* from) {
1012 // an unaligned load is required here as there is no requirement
1013 // on the alignment of input pointer 'from'
1014 __m256i low_half = _mm256_castps_si256(_mm256_loadu_ps(from));
1015 __m512 even_elements = _mm512_castsi512_ps(_mm512_cvtepu32_epi64(low_half));
1016 __m512 pairs = _mm512_permute_ps(even_elements, _MM_SHUFFLE(2, 2, 0, 0));
1017 return pairs;
1018}
1019
1020// Loads 4 doubles from memory a returns the packet {a0, a0, a1, a1, a2, a2, a3,
1021// a3}
1022template <>
1023EIGEN_STRONG_INLINE Packet8d ploaddup<Packet8d>(const double* from) {
1024 Packet8d tmp = _mm512_castpd256_pd512(ploadu<Packet4d>(from));
1025 const Packet8l scatter_mask = _mm512_set_epi64(3, 3, 2, 2, 1, 1, 0, 0);
1026 return _mm512_permutexvar_pd(scatter_mask, tmp);
1027}
1028
1029// Loads 4 int64_t from memory a returns the packet {a0, a0, a1, a1, a2, a2, a3,
1030// a3}
1031template <>
1032EIGEN_STRONG_INLINE Packet8l ploaddup<Packet8l>(const int64_t* from) {
1033 Packet8l tmp = _mm512_castsi256_si512(ploadu<Packet4l>(from));
1034 const Packet8l scatter_mask = _mm512_set_epi64(3, 3, 2, 2, 1, 1, 0, 0);
1035 return _mm512_permutexvar_epi64(scatter_mask, tmp);
1036}
1037
1038// Loads 8 integers from memory and returns the packet
1039// {a0, a0 a1, a1, a2, a2, a3, a3, a4, a4, a5, a5, a6, a6, a7, a7}
1040template <>
1041EIGEN_STRONG_INLINE Packet16i ploaddup<Packet16i>(const int* from) {
1042 __m256i low_half = _mm256_load_si256(reinterpret_cast<const __m256i*>(from));
1043 __m512 even_elements = _mm512_castsi512_ps(_mm512_cvtepu32_epi64(low_half));
1044 __m512 pairs = _mm512_permute_ps(even_elements, _MM_SHUFFLE(2, 2, 0, 0));
1045 return _mm512_castps_si512(pairs);
1046}
1047
1048// Loads 4 floats from memory a returns the packet
1049// {a0, a0 a0, a0, a1, a1, a1, a1, a2, a2, a2, a2, a3, a3, a3, a3}
1050template <>
1051EIGEN_STRONG_INLINE Packet16f ploadquad<Packet16f>(const float* from) {
1052 Packet16f tmp = _mm512_castps128_ps512(ploadu<Packet4f>(from));
1053 const Packet16i scatter_mask = _mm512_set_epi32(3, 3, 3, 3, 2, 2, 2, 2, 1, 1, 1, 1, 0, 0, 0, 0);
1054 return _mm512_permutexvar_ps(scatter_mask, tmp);
1055}
1056
1057// Loads 2 doubles from memory a returns the packet
1058// {a0, a0 a0, a0, a1, a1, a1, a1}
1059template <>
1060EIGEN_STRONG_INLINE Packet8d ploadquad<Packet8d>(const double* from) {
1061 __m256d lane0 = _mm256_set1_pd(*from);
1062 __m256d lane1 = _mm256_set1_pd(*(from + 1));
1063 __m512d tmp = _mm512_undefined_pd();
1064 tmp = _mm512_insertf64x4(tmp, lane0, 0);
1065 return _mm512_insertf64x4(tmp, lane1, 1);
1066}
1067
1068// Loads 2 int64_t from memory a returns the packet
1069// {a0, a0 a0, a0, a1, a1, a1, a1}
1070template <>
1071EIGEN_STRONG_INLINE Packet8l ploadquad<Packet8l>(const int64_t* from) {
1072 __m256i lane0 = _mm256_set1_epi64x(*from);
1073 __m256i lane1 = _mm256_set1_epi64x(*(from + 1));
1074 __m512i tmp = _mm512_undefined_epi32();
1075 tmp = _mm512_inserti64x4(tmp, lane0, 0);
1076 return _mm512_inserti64x4(tmp, lane1, 1);
1077}
1078
1079// Loads 4 integers from memory and returns the packet
1080// {a0, a0 a0, a0, a1, a1, a1, a1, a2, a2, a2, a2, a3, a3, a3, a3}
1081template <>
1082EIGEN_STRONG_INLINE Packet16i ploadquad<Packet16i>(const int* from) {
1083 Packet16i tmp = _mm512_castsi128_si512(ploadu<Packet4i>(from));
1084 const Packet16i scatter_mask = _mm512_set_epi32(3, 3, 3, 3, 2, 2, 2, 2, 1, 1, 1, 1, 0, 0, 0, 0);
1085 return _mm512_permutexvar_epi32(scatter_mask, tmp);
1086}
1087
1088template <>
1089EIGEN_STRONG_INLINE void pstore<float>(float* to, const Packet16f& from) {
1090 EIGEN_DEBUG_ALIGNED_STORE _mm512_store_ps(to, from);
1091}
1092template <>
1093EIGEN_STRONG_INLINE void pstore<double>(double* to, const Packet8d& from) {
1094 EIGEN_DEBUG_ALIGNED_STORE _mm512_store_pd(to, from);
1095}
1096template <>
1097EIGEN_STRONG_INLINE void pstore<int>(int* to, const Packet16i& from) {
1098 EIGEN_DEBUG_ALIGNED_STORE _mm512_store_epi32(to, from);
1099}
1100template <>
1101EIGEN_STRONG_INLINE void pstore<int64_t>(int64_t* to, const Packet8l& from) {
1102 EIGEN_DEBUG_ALIGNED_STORE _mm512_store_epi64(to, from);
1103}
1104
1105template <>
1106EIGEN_STRONG_INLINE void pstoreu<float>(float* to, const Packet16f& from) {
1107 EIGEN_DEBUG_UNALIGNED_STORE _mm512_storeu_ps(to, from);
1108}
1109template <>
1110EIGEN_STRONG_INLINE void pstoreu<double>(double* to, const Packet8d& from) {
1111 EIGEN_DEBUG_UNALIGNED_STORE _mm512_storeu_pd(to, from);
1112}
1113template <>
1114EIGEN_STRONG_INLINE void pstoreu<int>(int* to, const Packet16i& from) {
1115 EIGEN_DEBUG_UNALIGNED_STORE _mm512_storeu_epi32(to, from);
1116}
1117template <>
1118EIGEN_STRONG_INLINE void pstoreu<int64_t>(int64_t* to, const Packet8l& from) {
1119 EIGEN_DEBUG_UNALIGNED_STORE _mm512_storeu_epi64(to, from);
1120}
1121template <>
1122EIGEN_STRONG_INLINE void pstoreu<float>(float* to, const Packet16f& from, uint16_t umask) {
1123 __mmask16 mask = static_cast<__mmask16>(umask);
1124 EIGEN_DEBUG_UNALIGNED_STORE return _mm512_mask_storeu_ps(to, mask, from);
1125}
1126template <>
1127EIGEN_STRONG_INLINE void pstoreu<double>(double* to, const Packet8d& from, uint8_t umask) {
1128 __mmask8 mask = static_cast<__mmask8>(umask);
1129 EIGEN_DEBUG_UNALIGNED_STORE return _mm512_mask_storeu_pd(to, mask, from);
1130}
1131
1132template <typename Scalar, typename Packet>
1133EIGEN_DEVICE_FUNC inline Packet pgather(const Packet& src, const Scalar* from, Index stride,
1134 typename unpacket_traits<Packet>::mask_t umask);
1135template <>
1136EIGEN_DEVICE_FUNC inline Packet16f pgather<float, Packet16f>(const Packet16f& src, const float* from, Index stride,
1137 uint16_t umask) {
1138 Packet16i stride_vector = _mm512_set1_epi32(convert_index<int>(stride));
1139 Packet16i stride_multiplier = _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0);
1140 Packet16i indices = _mm512_mullo_epi32(stride_vector, stride_multiplier);
1141 __mmask16 mask = static_cast<__mmask16>(umask);
1142
1143 return _mm512_mask_i32gather_ps(src, mask, indices, from, 4);
1144}
1145template <>
1146EIGEN_DEVICE_FUNC inline Packet8d pgather<double, Packet8d>(const Packet8d& src, const double* from, Index stride,
1147 uint8_t umask) {
1148 Packet8i stride_vector = _mm256_set1_epi32(convert_index<int>(stride));
1149 Packet8i stride_multiplier = _mm256_set_epi32(7, 6, 5, 4, 3, 2, 1, 0);
1150 Packet8i indices = _mm256_mullo_epi32(stride_vector, stride_multiplier);
1151 __mmask8 mask = static_cast<__mmask8>(umask);
1152
1153 return _mm512_mask_i32gather_pd(src, mask, indices, from, 8);
1154}
1155
1156template <>
1157EIGEN_DEVICE_FUNC inline Packet16f pgather<float, Packet16f>(const float* from, Index stride) {
1158 Packet16i stride_vector = _mm512_set1_epi32(convert_index<int>(stride));
1159 Packet16i stride_multiplier = _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0);
1160 Packet16i indices = _mm512_mullo_epi32(stride_vector, stride_multiplier);
1161
1162 return _mm512_i32gather_ps(indices, from, 4);
1163}
1164template <>
1165EIGEN_DEVICE_FUNC inline Packet8d pgather<double, Packet8d>(const double* from, Index stride) {
1166 Packet8i stride_vector = _mm256_set1_epi32(convert_index<int>(stride));
1167 Packet8i stride_multiplier = _mm256_set_epi32(7, 6, 5, 4, 3, 2, 1, 0);
1168 Packet8i indices = _mm256_mullo_epi32(stride_vector, stride_multiplier);
1169
1170 return _mm512_i32gather_pd(indices, from, 8);
1171}
1172template <>
1173EIGEN_DEVICE_FUNC inline Packet8l pgather<int64_t, Packet8l>(const int64_t* from, Index stride) {
1174 Packet8i stride_vector = _mm256_set1_epi32(convert_index<int>(stride));
1175 Packet8i stride_multiplier = _mm256_set_epi32(7, 6, 5, 4, 3, 2, 1, 0);
1176 Packet8i indices = _mm256_mullo_epi32(stride_vector, stride_multiplier);
1177
1178 return _mm512_i32gather_epi64(indices, from, 8);
1179}
1180template <>
1181EIGEN_DEVICE_FUNC inline Packet16i pgather<int, Packet16i>(const int* from, Index stride) {
1182 Packet16i stride_vector = _mm512_set1_epi32(convert_index<int>(stride));
1183 Packet16i stride_multiplier = _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0);
1184 Packet16i indices = _mm512_mullo_epi32(stride_vector, stride_multiplier);
1185 return _mm512_i32gather_epi32(indices, from, 4);
1186}
1187
1188template <typename Scalar, typename Packet>
1189EIGEN_DEVICE_FUNC inline void pscatter(Scalar* to, const Packet& from, Index stride,
1190 typename unpacket_traits<Packet>::mask_t umask);
1191template <>
1192EIGEN_DEVICE_FUNC inline void pscatter<float, Packet16f>(float* to, const Packet16f& from, Index stride,
1193 uint16_t umask) {
1194 Packet16i stride_vector = _mm512_set1_epi32(convert_index<int>(stride));
1195 Packet16i stride_multiplier = _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0);
1196 Packet16i indices = _mm512_mullo_epi32(stride_vector, stride_multiplier);
1197 __mmask16 mask = static_cast<__mmask16>(umask);
1198 _mm512_mask_i32scatter_ps(to, mask, indices, from, 4);
1199}
1200template <>
1201EIGEN_DEVICE_FUNC inline void pscatter<double, Packet8d>(double* to, const Packet8d& from, Index stride,
1202 uint8_t umask) {
1203 Packet8i stride_vector = _mm256_set1_epi32(convert_index<int>(stride));
1204 Packet8i stride_multiplier = _mm256_set_epi32(7, 6, 5, 4, 3, 2, 1, 0);
1205 Packet8i indices = _mm256_mullo_epi32(stride_vector, stride_multiplier);
1206 __mmask8 mask = static_cast<__mmask8>(umask);
1207 _mm512_mask_i32scatter_pd(to, mask, indices, from, 8);
1208}
1209template <>
1210EIGEN_DEVICE_FUNC inline void pscatter<float, Packet16f>(float* to, const Packet16f& from, Index stride) {
1211 Packet16i stride_vector = _mm512_set1_epi32(convert_index<int>(stride));
1212 Packet16i stride_multiplier = _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0);
1213 Packet16i indices = _mm512_mullo_epi32(stride_vector, stride_multiplier);
1214 _mm512_i32scatter_ps(to, indices, from, 4);
1215}
1216template <>
1217EIGEN_DEVICE_FUNC inline void pscatter<double, Packet8d>(double* to, const Packet8d& from, Index stride) {
1218 Packet8i stride_vector = _mm256_set1_epi32(convert_index<int>(stride));
1219 Packet8i stride_multiplier = _mm256_set_epi32(7, 6, 5, 4, 3, 2, 1, 0);
1220 Packet8i indices = _mm256_mullo_epi32(stride_vector, stride_multiplier);
1221 _mm512_i32scatter_pd(to, indices, from, 8);
1222}
1223template <>
1224EIGEN_DEVICE_FUNC inline void pscatter<int64_t, Packet8l>(int64_t* to, const Packet8l& from, Index stride) {
1225 Packet8i stride_vector = _mm256_set1_epi32(convert_index<int>(stride));
1226 Packet8i stride_multiplier = _mm256_set_epi32(7, 6, 5, 4, 3, 2, 1, 0);
1227 Packet8i indices = _mm256_mullo_epi32(stride_vector, stride_multiplier);
1228 _mm512_i32scatter_epi64(to, indices, from, 8);
1229}
1230template <>
1231EIGEN_DEVICE_FUNC inline void pscatter<int, Packet16i>(int* to, const Packet16i& from, Index stride) {
1232 Packet16i stride_vector = _mm512_set1_epi32(convert_index<int>(stride));
1233 Packet16i stride_multiplier = _mm512_set_epi32(15, 14, 13, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0);
1234 Packet16i indices = _mm512_mullo_epi32(stride_vector, stride_multiplier);
1235 _mm512_i32scatter_epi32(to, indices, from, 4);
1236}
1237
1238template <>
1239EIGEN_STRONG_INLINE void pstore1<Packet16f>(float* to, const float& a) {
1240 Packet16f pa = pset1<Packet16f>(a);
1241 pstore(to, pa);
1242}
1243template <>
1244EIGEN_STRONG_INLINE void pstore1<Packet8d>(double* to, const double& a) {
1245 Packet8d pa = pset1<Packet8d>(a);
1246 pstore(to, pa);
1247}
1248template <>
1249EIGEN_STRONG_INLINE void pstore1<Packet16i>(int* to, const int& a) {
1250 Packet16i pa = pset1<Packet16i>(a);
1251 pstore(to, pa);
1252}
1253template <>
1254EIGEN_STRONG_INLINE void pstore1<Packet8l>(int64_t* to, const int64_t& a) {
1255 Packet8l pa = pset1<Packet8l>(a);
1256 pstore(to, pa);
1257}
1258
1259template <>
1260EIGEN_STRONG_INLINE void prefetch<float>(const float* addr) {
1261 _mm_prefetch((SsePrefetchPtrType)(addr), _MM_HINT_T0);
1262}
1263template <>
1264EIGEN_STRONG_INLINE void prefetch<double>(const double* addr) {
1265 _mm_prefetch((SsePrefetchPtrType)(addr), _MM_HINT_T0);
1266}
1267template <>
1268EIGEN_STRONG_INLINE void prefetch<int>(const int* addr) {
1269 _mm_prefetch((SsePrefetchPtrType)(addr), _MM_HINT_T0);
1270}
1271
1272template <>
1273EIGEN_STRONG_INLINE float pfirst<Packet16f>(const Packet16f& a) {
1274 return _mm512_cvtss_f32(a);
1275}
1276template <>
1277EIGEN_STRONG_INLINE double pfirst<Packet8d>(const Packet8d& a) {
1278 return _mm512_cvtsd_f64(a);
1279}
1280template <>
1281EIGEN_STRONG_INLINE int64_t pfirst<Packet8l>(const Packet8l& a) {
1282 int64_t x = _mm_extract_epi64_0(_mm512_extracti32x4_epi32(a, 0));
1283 return x;
1284}
1285template <>
1286EIGEN_STRONG_INLINE int pfirst<Packet16i>(const Packet16i& a) {
1287#if EIGEN_GNUC_STRICT_LESS_THAN(11, 0, 0)
1288 return _mm_cvtsi128_si32(_mm512_castsi512_si128(a));
1289#else
1290 return _mm512_cvtsi512_si32(a);
1291#endif
1292}
1293
1294template <>
1295EIGEN_STRONG_INLINE Packet16f preverse(const Packet16f& a) {
1296 return _mm512_permutexvar_ps(_mm512_set_epi32(0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15), a);
1297}
1298
1299template <>
1300EIGEN_STRONG_INLINE Packet8d preverse(const Packet8d& a) {
1301 return _mm512_permutexvar_pd(_mm512_set_epi32(0, 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7), a);
1302}
1303
1304template <>
1305EIGEN_STRONG_INLINE Packet16i preverse(const Packet16i& a) {
1306 return _mm512_permutexvar_epi32(_mm512_set_epi32(0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15), a);
1307}
1308
1309template <>
1310EIGEN_STRONG_INLINE Packet8l preverse(const Packet8l& a) {
1311 return _mm512_permutexvar_epi64(_mm512_set_epi64(0, 1, 2, 3, 4, 5, 6, 7), a);
1312}
1313
1314template <>
1315EIGEN_STRONG_INLINE Packet16f pabs(const Packet16f& a) {
1316 // _mm512_abs_ps intrinsic not found, so hack around it
1317 return _mm512_castsi512_ps(_mm512_and_si512(_mm512_castps_si512(a), _mm512_set1_epi32(0x7fffffff)));
1318}
1319template <>
1320EIGEN_STRONG_INLINE Packet8d pabs(const Packet8d& a) {
1321 // _mm512_abs_ps intrinsic not found, so hack around it
1322 return _mm512_castsi512_pd(_mm512_and_si512(_mm512_castpd_si512(a), _mm512_set1_epi64(0x7fffffffffffffff)));
1323}
1324template <>
1325EIGEN_STRONG_INLINE Packet16i pabs(const Packet16i& a) {
1326 return _mm512_abs_epi32(a);
1327}
1328template <>
1329EIGEN_STRONG_INLINE Packet8l pabs(const Packet8l& a) {
1330 return _mm512_abs_epi64(a);
1331}
1332
1333template <>
1334EIGEN_STRONG_INLINE Packet16h psignbit(const Packet16h& a) {
1335 return _mm256_srai_epi16(a, 15);
1336}
1337template <>
1338EIGEN_STRONG_INLINE Packet16bf psignbit(const Packet16bf& a) {
1339 return _mm256_srai_epi16(a, 15);
1340}
1341template <>
1342EIGEN_STRONG_INLINE Packet16f psignbit(const Packet16f& a) {
1343 return _mm512_castsi512_ps(_mm512_srai_epi32(_mm512_castps_si512(a), 31));
1344}
1345template <>
1346EIGEN_STRONG_INLINE Packet8d psignbit(const Packet8d& a) {
1347 return _mm512_castsi512_pd(_mm512_srai_epi64(_mm512_castpd_si512(a), 63));
1348}
1349
1350template <>
1351EIGEN_STRONG_INLINE Packet16f pfrexp<Packet16f>(const Packet16f& a, Packet16f& exponent) {
1352 return pfrexp_generic(a, exponent);
1353}
1354
1355// Extract exponent without existence of Packet8l.
1356template <>
1357EIGEN_STRONG_INLINE Packet8d pfrexp_generic_get_biased_exponent(const Packet8d& a) {
1358 const Packet8d cst_exp_mask = pset1frombits<Packet8d>(static_cast<uint64_t>(0x7ff0000000000000ull));
1359#ifdef EIGEN_VECTORIZE_AVX512DQ
1360 return _mm512_cvtepi64_pd(_mm512_srli_epi64(_mm512_castpd_si512(pand(a, cst_exp_mask)), 52));
1361#else
1362 return _mm512_cvtepi32_pd(_mm512_cvtepi64_epi32(_mm512_srli_epi64(_mm512_castpd_si512(pand(a, cst_exp_mask)), 52)));
1363#endif
1364}
1365
1366template <>
1367EIGEN_STRONG_INLINE Packet8d pfrexp<Packet8d>(const Packet8d& a, Packet8d& exponent) {
1368 return pfrexp_generic(a, exponent);
1369}
1370
1371template <>
1372EIGEN_STRONG_INLINE Packet16f pldexp<Packet16f>(const Packet16f& a, const Packet16f& exponent) {
1373 return pldexp_generic(a, exponent);
1374}
1375
1376template <>
1377EIGEN_STRONG_INLINE Packet8d pldexp<Packet8d>(const Packet8d& a, const Packet8d& exponent) {
1378 // Clamp exponent to [-2099, 2099]
1379 const Packet8d max_exponent = pset1<Packet8d>(2099.0);
1380 const Packet8i e = _mm512_cvtpd_epi32(pmin(pmax(exponent, pnegate(max_exponent)), max_exponent));
1381
1382 // Split 2^e into four factors and multiply.
1383 const Packet8i bias = pset1<Packet8i>(1023);
1384 Packet8i b = parithmetic_shift_right<2>(e); // floor(e/4)
1385
1386 // 2^b
1387 const Packet8i permute_idx = _mm256_setr_epi32(0, 4, 1, 5, 2, 6, 3, 7);
1388 Packet8i hi = _mm256_permutevar8x32_epi32(padd(b, bias), permute_idx);
1389 Packet8i lo = _mm256_slli_epi64(hi, 52);
1390 hi = _mm256_slli_epi64(_mm256_srli_epi64(hi, 32), 52);
1391 Packet8d c = _mm512_castsi512_pd(_mm512_inserti64x4(_mm512_castsi256_si512(lo), hi, 1));
1392 Packet8d out = pmul(pmul(pmul(a, c), c), c); // a * 2^(3b)
1393
1394 // 2^(e - 3b)
1395 b = psub(psub(psub(e, b), b), b); // e - 3b
1396 hi = _mm256_permutevar8x32_epi32(padd(b, bias), permute_idx);
1397 lo = _mm256_slli_epi64(hi, 52);
1398 hi = _mm256_slli_epi64(_mm256_srli_epi64(hi, 32), 52);
1399 c = _mm512_castsi512_pd(_mm512_inserti64x4(_mm512_castsi256_si512(lo), hi, 1));
1400 out = pmul(out, c); // a * 2^e
1401 return out;
1402}
1403
1404#ifdef EIGEN_VECTORIZE_AVX512DQ
1405// AVX512F does not define _mm512_extractf32x8_ps to extract _m256 from _m512
1406#define EIGEN_EXTRACT_8f_FROM_16f(INPUT, OUTPUT) \
1407 __m256 OUTPUT##_0 = _mm512_extractf32x8_ps(INPUT, 0); \
1408 __m256 OUTPUT##_1 = _mm512_extractf32x8_ps(INPUT, 1)
1409
1410// AVX512F does not define _mm512_extracti32x8_epi32 to extract _m256i from _m512i
1411#define EIGEN_EXTRACT_8i_FROM_16i(INPUT, OUTPUT) \
1412 __m256i OUTPUT##_0 = _mm512_extracti32x8_epi32(INPUT, 0); \
1413 __m256i OUTPUT##_1 = _mm512_extracti32x8_epi32(INPUT, 1)
1414#else
1415#define EIGEN_EXTRACT_8f_FROM_16f(INPUT, OUTPUT) \
1416 __m256 OUTPUT##_0 = _mm256_insertf128_ps(_mm256_castps128_ps256(_mm512_extractf32x4_ps(INPUT, 0)), \
1417 _mm512_extractf32x4_ps(INPUT, 1), 1); \
1418 __m256 OUTPUT##_1 = _mm256_insertf128_ps(_mm256_castps128_ps256(_mm512_extractf32x4_ps(INPUT, 2)), \
1419 _mm512_extractf32x4_ps(INPUT, 3), 1)
1420
1421#define EIGEN_EXTRACT_8i_FROM_16i(INPUT, OUTPUT) \
1422 __m256i OUTPUT##_0 = _mm256_insertf128_si256(_mm256_castsi128_si256(_mm512_extracti32x4_epi32(INPUT, 0)), \
1423 _mm512_extracti32x4_epi32(INPUT, 1), 1); \
1424 __m256i OUTPUT##_1 = _mm256_insertf128_si256(_mm256_castsi128_si256(_mm512_extracti32x4_epi32(INPUT, 2)), \
1425 _mm512_extracti32x4_epi32(INPUT, 3), 1)
1426#endif
1427
1428#ifdef EIGEN_VECTORIZE_AVX512DQ
1429#define EIGEN_INSERT_8f_INTO_16f(OUTPUT, INPUTA, INPUTB) \
1430 OUTPUT = _mm512_insertf32x8(_mm512_castps256_ps512(INPUTA), INPUTB, 1);
1431
1432#define EIGEN_INSERT_8i_INTO_16i(OUTPUT, INPUTA, INPUTB) \
1433 OUTPUT = _mm512_inserti32x8(_mm512_castsi256_si512(INPUTA), INPUTB, 1);
1434#else
1435#define EIGEN_INSERT_8f_INTO_16f(OUTPUT, INPUTA, INPUTB) \
1436 OUTPUT = _mm512_undefined_ps(); \
1437 OUTPUT = _mm512_insertf32x4(OUTPUT, _mm256_extractf128_ps(INPUTA, 0), 0); \
1438 OUTPUT = _mm512_insertf32x4(OUTPUT, _mm256_extractf128_ps(INPUTA, 1), 1); \
1439 OUTPUT = _mm512_insertf32x4(OUTPUT, _mm256_extractf128_ps(INPUTB, 0), 2); \
1440 OUTPUT = _mm512_insertf32x4(OUTPUT, _mm256_extractf128_ps(INPUTB, 1), 3);
1441
1442#define EIGEN_INSERT_8i_INTO_16i(OUTPUT, INPUTA, INPUTB) \
1443 OUTPUT = _mm512_undefined_epi32(); \
1444 OUTPUT = _mm512_inserti32x4(OUTPUT, _mm256_extractf128_si256(INPUTA, 0), 0); \
1445 OUTPUT = _mm512_inserti32x4(OUTPUT, _mm256_extractf128_si256(INPUTA, 1), 1); \
1446 OUTPUT = _mm512_inserti32x4(OUTPUT, _mm256_extractf128_si256(INPUTB, 0), 2); \
1447 OUTPUT = _mm512_inserti32x4(OUTPUT, _mm256_extractf128_si256(INPUTB, 1), 3);
1448#endif
1449
1450template <>
1451EIGEN_STRONG_INLINE float predux<Packet16f>(const Packet16f& a) {
1452#ifdef EIGEN_VECTORIZE_AVX512DQ
1453 __m256 lane0 = _mm512_extractf32x8_ps(a, 0);
1454 __m256 lane1 = _mm512_extractf32x8_ps(a, 1);
1455 Packet8f x = _mm256_add_ps(lane0, lane1);
1456 return predux<Packet8f>(x);
1457#else
1458 __m128 lane0 = _mm512_extractf32x4_ps(a, 0);
1459 __m128 lane1 = _mm512_extractf32x4_ps(a, 1);
1460 __m128 lane2 = _mm512_extractf32x4_ps(a, 2);
1461 __m128 lane3 = _mm512_extractf32x4_ps(a, 3);
1462 __m128 sum = _mm_add_ps(_mm_add_ps(lane0, lane1), _mm_add_ps(lane2, lane3));
1463 return predux<Packet4f>(sum);
1464#endif
1465}
1466template <>
1467EIGEN_STRONG_INLINE double predux<Packet8d>(const Packet8d& a) {
1468 __m256d lane0 = _mm512_extractf64x4_pd(a, 0);
1469 __m256d lane1 = _mm512_extractf64x4_pd(a, 1);
1470 __m256d sum = _mm256_add_pd(lane0, lane1);
1471 return predux<Packet4d>(sum);
1472}
1473
1474template <>
1475EIGEN_STRONG_INLINE int64_t predux<Packet8l>(const Packet8l& a) {
1476 return _mm512_reduce_add_epi64(a);
1477}
1478
1479template <>
1480EIGEN_STRONG_INLINE int predux<Packet16i>(const Packet16i& a) {
1481 return _mm512_reduce_add_epi32(a);
1482}
1483
1484template <>
1485EIGEN_STRONG_INLINE Packet8f predux_half_dowto4<Packet16f>(const Packet16f& a) {
1486#ifdef EIGEN_VECTORIZE_AVX512DQ
1487 __m256 lane0 = _mm512_extractf32x8_ps(a, 0);
1488 __m256 lane1 = _mm512_extractf32x8_ps(a, 1);
1489 return _mm256_add_ps(lane0, lane1);
1490#else
1491 __m128 lane0 = _mm512_extractf32x4_ps(a, 0);
1492 __m128 lane1 = _mm512_extractf32x4_ps(a, 1);
1493 __m128 lane2 = _mm512_extractf32x4_ps(a, 2);
1494 __m128 lane3 = _mm512_extractf32x4_ps(a, 3);
1495 __m128 sum0 = _mm_add_ps(lane0, lane2);
1496 __m128 sum1 = _mm_add_ps(lane1, lane3);
1497 return _mm256_insertf128_ps(_mm256_castps128_ps256(sum0), sum1, 1);
1498#endif
1499}
1500template <>
1501EIGEN_STRONG_INLINE Packet4d predux_half_dowto4<Packet8d>(const Packet8d& a) {
1502 __m256d lane0 = _mm512_extractf64x4_pd(a, 0);
1503 __m256d lane1 = _mm512_extractf64x4_pd(a, 1);
1504 return _mm256_add_pd(lane0, lane1);
1505}
1506template <>
1507EIGEN_STRONG_INLINE Packet8i predux_half_dowto4<Packet16i>(const Packet16i& a) {
1508#ifdef EIGEN_VECTORIZE_AVX512DQ
1509 __m256i lane0 = _mm512_extracti32x8_epi32(a, 0);
1510 __m256i lane1 = _mm512_extracti32x8_epi32(a, 1);
1511 return _mm256_add_epi32(lane0, lane1);
1512#else
1513 __m128i lane0 = _mm512_extracti32x4_epi32(a, 0);
1514 __m128i lane1 = _mm512_extracti32x4_epi32(a, 1);
1515 __m128i lane2 = _mm512_extracti32x4_epi32(a, 2);
1516 __m128i lane3 = _mm512_extracti32x4_epi32(a, 3);
1517 __m128i sum0 = _mm_add_epi32(lane0, lane2);
1518 __m128i sum1 = _mm_add_epi32(lane1, lane3);
1519 return _mm256_inserti128_si256(_mm256_castsi128_si256(sum0), sum1, 1);
1520#endif
1521}
1522
1523template <>
1524EIGEN_STRONG_INLINE Packet4l predux_half_dowto4<Packet8l>(const Packet8l& a) {
1525 __m256i lane0 = _mm512_extracti64x4_epi64(a, 0);
1526 __m256i lane1 = _mm512_extracti64x4_epi64(a, 1);
1527 return _mm256_add_epi64(lane0, lane1);
1528}
1529
1530template <>
1531EIGEN_STRONG_INLINE float predux_mul<Packet16f>(const Packet16f& a) {
1532// #ifdef EIGEN_VECTORIZE_AVX512DQ
1533#if 0
1534 Packet8f lane0 = _mm512_extractf32x8_ps(a, 0);
1535 Packet8f lane1 = _mm512_extractf32x8_ps(a, 1);
1536 Packet8f res = pmul(lane0, lane1);
1537 res = pmul(res, _mm256_permute2f128_ps(res, res, 1));
1538 res = pmul(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 3, 2)));
1539 return pfirst(pmul(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 0, 1))));
1540#else
1541 __m128 lane0 = _mm512_extractf32x4_ps(a, 0);
1542 __m128 lane1 = _mm512_extractf32x4_ps(a, 1);
1543 __m128 lane2 = _mm512_extractf32x4_ps(a, 2);
1544 __m128 lane3 = _mm512_extractf32x4_ps(a, 3);
1545 __m128 res = pmul(pmul(lane0, lane1), pmul(lane2, lane3));
1546 res = pmul(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 3, 2)));
1547 return pfirst(pmul(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 0, 1))));
1548#endif
1549}
1550template <>
1551EIGEN_STRONG_INLINE double predux_mul<Packet8d>(const Packet8d& a) {
1552 __m256d lane0 = _mm512_extractf64x4_pd(a, 0);
1553 __m256d lane1 = _mm512_extractf64x4_pd(a, 1);
1554 __m256d res = pmul(lane0, lane1);
1555 res = pmul(res, _mm256_permute2f128_pd(res, res, 1));
1556 return pfirst(pmul(res, _mm256_shuffle_pd(res, res, 1)));
1557}
1558template <>
1559EIGEN_STRONG_INLINE int predux_mul<Packet16i>(const Packet16i& a) {
1560 return _mm512_reduce_mul_epi32(a);
1561}
1562
1563#if EIGEN_COMP_MSVC
1564// MSVC's _mm512_reduce_mul_epi64 is borked, at least up to and including 1939.
1565// alignas(64) int64_t data[] = { 1,1,-1,-1,1,-1,-1,-1 };
1566// int64_t out = _mm512_reduce_mul_epi64(_mm512_load_epi64(data));
1567// produces garbage: 4294967295. It seems to happen whenever the output is supposed to be negative.
1568// Fall back to a manual approach:
1569template <>
1570EIGEN_STRONG_INLINE int64_t predux_mul<Packet8l>(const Packet8l& a) {
1571 Packet4l lane0 = _mm512_extracti64x4_epi64(a, 0);
1572 Packet4l lane1 = _mm512_extracti64x4_epi64(a, 1);
1573 Packet4l res = pmul(lane0, lane1);
1574 res = pmul(res, Packet4l(_mm256_permute2x128_si256(res, res, 1)));
1575 res = pmul(res, Packet4l(_mm256_shuffle_epi32(res, 0xE)));
1576 return pfirst(res);
1577}
1578#else
1579template <>
1580EIGEN_STRONG_INLINE int64_t predux_mul<Packet8l>(const Packet8l& a) {
1581 return _mm512_reduce_mul_epi64(a);
1582}
1583#endif
1584
1585template <>
1586EIGEN_STRONG_INLINE float predux_min<Packet16f>(const Packet16f& a) {
1587 __m128 lane0 = _mm512_extractf32x4_ps(a, 0);
1588 __m128 lane1 = _mm512_extractf32x4_ps(a, 1);
1589 __m128 lane2 = _mm512_extractf32x4_ps(a, 2);
1590 __m128 lane3 = _mm512_extractf32x4_ps(a, 3);
1591 __m128 res = _mm_min_ps(_mm_min_ps(lane0, lane1), _mm_min_ps(lane2, lane3));
1592 res = _mm_min_ps(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 3, 2)));
1593 return pfirst(_mm_min_ps(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 0, 1))));
1594}
1595template <>
1596EIGEN_STRONG_INLINE double predux_min<Packet8d>(const Packet8d& a) {
1597 __m256d lane0 = _mm512_extractf64x4_pd(a, 0);
1598 __m256d lane1 = _mm512_extractf64x4_pd(a, 1);
1599 __m256d res = _mm256_min_pd(lane0, lane1);
1600 res = _mm256_min_pd(res, _mm256_permute2f128_pd(res, res, 1));
1601 return pfirst(_mm256_min_pd(res, _mm256_shuffle_pd(res, res, 1)));
1602}
1603template <>
1604EIGEN_STRONG_INLINE int predux_min<Packet16i>(const Packet16i& a) {
1605 return _mm512_reduce_min_epi32(a);
1606}
1607template <>
1608EIGEN_STRONG_INLINE int64_t predux_min<Packet8l>(const Packet8l& a) {
1609 return _mm512_reduce_min_epi64(a);
1610}
1611
1612template <>
1613EIGEN_STRONG_INLINE float predux_max<Packet16f>(const Packet16f& a) {
1614 __m128 lane0 = _mm512_extractf32x4_ps(a, 0);
1615 __m128 lane1 = _mm512_extractf32x4_ps(a, 1);
1616 __m128 lane2 = _mm512_extractf32x4_ps(a, 2);
1617 __m128 lane3 = _mm512_extractf32x4_ps(a, 3);
1618 __m128 res = _mm_max_ps(_mm_max_ps(lane0, lane1), _mm_max_ps(lane2, lane3));
1619 res = _mm_max_ps(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 3, 2)));
1620 return pfirst(_mm_max_ps(res, _mm_permute_ps(res, _MM_SHUFFLE(0, 0, 0, 1))));
1621}
1622
1623template <>
1624EIGEN_STRONG_INLINE double predux_max<Packet8d>(const Packet8d& a) {
1625 __m256d lane0 = _mm512_extractf64x4_pd(a, 0);
1626 __m256d lane1 = _mm512_extractf64x4_pd(a, 1);
1627 __m256d res = _mm256_max_pd(lane0, lane1);
1628 res = _mm256_max_pd(res, _mm256_permute2f128_pd(res, res, 1));
1629 return pfirst(_mm256_max_pd(res, _mm256_shuffle_pd(res, res, 1)));
1630}
1631template <>
1632EIGEN_STRONG_INLINE int predux_max<Packet16i>(const Packet16i& a) {
1633 return _mm512_reduce_max_epi32(a);
1634}
1635template <>
1636EIGEN_STRONG_INLINE int64_t predux_max<Packet8l>(const Packet8l& a) {
1637 return _mm512_reduce_max_epi64(a);
1638}
1639
1640template <>
1641EIGEN_STRONG_INLINE bool predux_any(const Packet16f& x) {
1642 Packet16i xi = _mm512_castps_si512(x);
1643 __mmask16 tmp = _mm512_test_epi32_mask(xi, xi);
1644 return !_mm512_kortestz(tmp, tmp);
1645}
1646
1647template <>
1648EIGEN_STRONG_INLINE bool predux_any(const Packet16i& x) {
1649 __mmask16 tmp = _mm512_test_epi32_mask(x, x);
1650 return !_mm512_kortestz(tmp, tmp);
1651}
1652
1653#define PACK_OUTPUT(OUTPUT, INPUT, INDEX, STRIDE) \
1654 EIGEN_INSERT_8f_INTO_16f(OUTPUT[INDEX], INPUT[INDEX], INPUT[INDEX + STRIDE]);
1655
1656EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet16f, 16>& kernel) {
1657 __m512 T0 = _mm512_unpacklo_ps(kernel.packet[0], kernel.packet[1]);
1658 __m512 T1 = _mm512_unpackhi_ps(kernel.packet[0], kernel.packet[1]);
1659 __m512 T2 = _mm512_unpacklo_ps(kernel.packet[2], kernel.packet[3]);
1660 __m512 T3 = _mm512_unpackhi_ps(kernel.packet[2], kernel.packet[3]);
1661 __m512 T4 = _mm512_unpacklo_ps(kernel.packet[4], kernel.packet[5]);
1662 __m512 T5 = _mm512_unpackhi_ps(kernel.packet[4], kernel.packet[5]);
1663 __m512 T6 = _mm512_unpacklo_ps(kernel.packet[6], kernel.packet[7]);
1664 __m512 T7 = _mm512_unpackhi_ps(kernel.packet[6], kernel.packet[7]);
1665 __m512 T8 = _mm512_unpacklo_ps(kernel.packet[8], kernel.packet[9]);
1666 __m512 T9 = _mm512_unpackhi_ps(kernel.packet[8], kernel.packet[9]);
1667 __m512 T10 = _mm512_unpacklo_ps(kernel.packet[10], kernel.packet[11]);
1668 __m512 T11 = _mm512_unpackhi_ps(kernel.packet[10], kernel.packet[11]);
1669 __m512 T12 = _mm512_unpacklo_ps(kernel.packet[12], kernel.packet[13]);
1670 __m512 T13 = _mm512_unpackhi_ps(kernel.packet[12], kernel.packet[13]);
1671 __m512 T14 = _mm512_unpacklo_ps(kernel.packet[14], kernel.packet[15]);
1672 __m512 T15 = _mm512_unpackhi_ps(kernel.packet[14], kernel.packet[15]);
1673 __m512 S0 = _mm512_shuffle_ps(T0, T2, _MM_SHUFFLE(1, 0, 1, 0));
1674 __m512 S1 = _mm512_shuffle_ps(T0, T2, _MM_SHUFFLE(3, 2, 3, 2));
1675 __m512 S2 = _mm512_shuffle_ps(T1, T3, _MM_SHUFFLE(1, 0, 1, 0));
1676 __m512 S3 = _mm512_shuffle_ps(T1, T3, _MM_SHUFFLE(3, 2, 3, 2));
1677 __m512 S4 = _mm512_shuffle_ps(T4, T6, _MM_SHUFFLE(1, 0, 1, 0));
1678 __m512 S5 = _mm512_shuffle_ps(T4, T6, _MM_SHUFFLE(3, 2, 3, 2));
1679 __m512 S6 = _mm512_shuffle_ps(T5, T7, _MM_SHUFFLE(1, 0, 1, 0));
1680 __m512 S7 = _mm512_shuffle_ps(T5, T7, _MM_SHUFFLE(3, 2, 3, 2));
1681 __m512 S8 = _mm512_shuffle_ps(T8, T10, _MM_SHUFFLE(1, 0, 1, 0));
1682 __m512 S9 = _mm512_shuffle_ps(T8, T10, _MM_SHUFFLE(3, 2, 3, 2));
1683 __m512 S10 = _mm512_shuffle_ps(T9, T11, _MM_SHUFFLE(1, 0, 1, 0));
1684 __m512 S11 = _mm512_shuffle_ps(T9, T11, _MM_SHUFFLE(3, 2, 3, 2));
1685 __m512 S12 = _mm512_shuffle_ps(T12, T14, _MM_SHUFFLE(1, 0, 1, 0));
1686 __m512 S13 = _mm512_shuffle_ps(T12, T14, _MM_SHUFFLE(3, 2, 3, 2));
1687 __m512 S14 = _mm512_shuffle_ps(T13, T15, _MM_SHUFFLE(1, 0, 1, 0));
1688 __m512 S15 = _mm512_shuffle_ps(T13, T15, _MM_SHUFFLE(3, 2, 3, 2));
1689
1690 EIGEN_EXTRACT_8f_FROM_16f(S0, S0);
1691 EIGEN_EXTRACT_8f_FROM_16f(S1, S1);
1692 EIGEN_EXTRACT_8f_FROM_16f(S2, S2);
1693 EIGEN_EXTRACT_8f_FROM_16f(S3, S3);
1694 EIGEN_EXTRACT_8f_FROM_16f(S4, S4);
1695 EIGEN_EXTRACT_8f_FROM_16f(S5, S5);
1696 EIGEN_EXTRACT_8f_FROM_16f(S6, S6);
1697 EIGEN_EXTRACT_8f_FROM_16f(S7, S7);
1698 EIGEN_EXTRACT_8f_FROM_16f(S8, S8);
1699 EIGEN_EXTRACT_8f_FROM_16f(S9, S9);
1700 EIGEN_EXTRACT_8f_FROM_16f(S10, S10);
1701 EIGEN_EXTRACT_8f_FROM_16f(S11, S11);
1702 EIGEN_EXTRACT_8f_FROM_16f(S12, S12);
1703 EIGEN_EXTRACT_8f_FROM_16f(S13, S13);
1704 EIGEN_EXTRACT_8f_FROM_16f(S14, S14);
1705 EIGEN_EXTRACT_8f_FROM_16f(S15, S15);
1706
1707 PacketBlock<Packet8f, 32> tmp;
1708
1709 tmp.packet[0] = _mm256_permute2f128_ps(S0_0, S4_0, 0x20);
1710 tmp.packet[1] = _mm256_permute2f128_ps(S1_0, S5_0, 0x20);
1711 tmp.packet[2] = _mm256_permute2f128_ps(S2_0, S6_0, 0x20);
1712 tmp.packet[3] = _mm256_permute2f128_ps(S3_0, S7_0, 0x20);
1713 tmp.packet[4] = _mm256_permute2f128_ps(S0_0, S4_0, 0x31);
1714 tmp.packet[5] = _mm256_permute2f128_ps(S1_0, S5_0, 0x31);
1715 tmp.packet[6] = _mm256_permute2f128_ps(S2_0, S6_0, 0x31);
1716 tmp.packet[7] = _mm256_permute2f128_ps(S3_0, S7_0, 0x31);
1717
1718 tmp.packet[8] = _mm256_permute2f128_ps(S0_1, S4_1, 0x20);
1719 tmp.packet[9] = _mm256_permute2f128_ps(S1_1, S5_1, 0x20);
1720 tmp.packet[10] = _mm256_permute2f128_ps(S2_1, S6_1, 0x20);
1721 tmp.packet[11] = _mm256_permute2f128_ps(S3_1, S7_1, 0x20);
1722 tmp.packet[12] = _mm256_permute2f128_ps(S0_1, S4_1, 0x31);
1723 tmp.packet[13] = _mm256_permute2f128_ps(S1_1, S5_1, 0x31);
1724 tmp.packet[14] = _mm256_permute2f128_ps(S2_1, S6_1, 0x31);
1725 tmp.packet[15] = _mm256_permute2f128_ps(S3_1, S7_1, 0x31);
1726
1727 // Second set of _m256 outputs
1728 tmp.packet[16] = _mm256_permute2f128_ps(S8_0, S12_0, 0x20);
1729 tmp.packet[17] = _mm256_permute2f128_ps(S9_0, S13_0, 0x20);
1730 tmp.packet[18] = _mm256_permute2f128_ps(S10_0, S14_0, 0x20);
1731 tmp.packet[19] = _mm256_permute2f128_ps(S11_0, S15_0, 0x20);
1732 tmp.packet[20] = _mm256_permute2f128_ps(S8_0, S12_0, 0x31);
1733 tmp.packet[21] = _mm256_permute2f128_ps(S9_0, S13_0, 0x31);
1734 tmp.packet[22] = _mm256_permute2f128_ps(S10_0, S14_0, 0x31);
1735 tmp.packet[23] = _mm256_permute2f128_ps(S11_0, S15_0, 0x31);
1736
1737 tmp.packet[24] = _mm256_permute2f128_ps(S8_1, S12_1, 0x20);
1738 tmp.packet[25] = _mm256_permute2f128_ps(S9_1, S13_1, 0x20);
1739 tmp.packet[26] = _mm256_permute2f128_ps(S10_1, S14_1, 0x20);
1740 tmp.packet[27] = _mm256_permute2f128_ps(S11_1, S15_1, 0x20);
1741 tmp.packet[28] = _mm256_permute2f128_ps(S8_1, S12_1, 0x31);
1742 tmp.packet[29] = _mm256_permute2f128_ps(S9_1, S13_1, 0x31);
1743 tmp.packet[30] = _mm256_permute2f128_ps(S10_1, S14_1, 0x31);
1744 tmp.packet[31] = _mm256_permute2f128_ps(S11_1, S15_1, 0x31);
1745
1746 // Pack them into the output
1747 PACK_OUTPUT(kernel.packet, tmp.packet, 0, 16);
1748 PACK_OUTPUT(kernel.packet, tmp.packet, 1, 16);
1749 PACK_OUTPUT(kernel.packet, tmp.packet, 2, 16);
1750 PACK_OUTPUT(kernel.packet, tmp.packet, 3, 16);
1751
1752 PACK_OUTPUT(kernel.packet, tmp.packet, 4, 16);
1753 PACK_OUTPUT(kernel.packet, tmp.packet, 5, 16);
1754 PACK_OUTPUT(kernel.packet, tmp.packet, 6, 16);
1755 PACK_OUTPUT(kernel.packet, tmp.packet, 7, 16);
1756
1757 PACK_OUTPUT(kernel.packet, tmp.packet, 8, 16);
1758 PACK_OUTPUT(kernel.packet, tmp.packet, 9, 16);
1759 PACK_OUTPUT(kernel.packet, tmp.packet, 10, 16);
1760 PACK_OUTPUT(kernel.packet, tmp.packet, 11, 16);
1761
1762 PACK_OUTPUT(kernel.packet, tmp.packet, 12, 16);
1763 PACK_OUTPUT(kernel.packet, tmp.packet, 13, 16);
1764 PACK_OUTPUT(kernel.packet, tmp.packet, 14, 16);
1765 PACK_OUTPUT(kernel.packet, tmp.packet, 15, 16);
1766}
1767#define PACK_OUTPUT_2(OUTPUT, INPUT, INDEX, STRIDE) \
1768 EIGEN_INSERT_8f_INTO_16f(OUTPUT[INDEX], INPUT[2 * INDEX], INPUT[2 * INDEX + STRIDE]);
1769
1770EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet16f, 8>& kernel) {
1771 __m512 T0 = _mm512_unpacklo_ps(kernel.packet[0], kernel.packet[1]);
1772 __m512 T1 = _mm512_unpackhi_ps(kernel.packet[0], kernel.packet[1]);
1773 __m512 T2 = _mm512_unpacklo_ps(kernel.packet[2], kernel.packet[3]);
1774 __m512 T3 = _mm512_unpackhi_ps(kernel.packet[2], kernel.packet[3]);
1775 __m512 T4 = _mm512_unpacklo_ps(kernel.packet[4], kernel.packet[5]);
1776 __m512 T5 = _mm512_unpackhi_ps(kernel.packet[4], kernel.packet[5]);
1777 __m512 T6 = _mm512_unpacklo_ps(kernel.packet[6], kernel.packet[7]);
1778 __m512 T7 = _mm512_unpackhi_ps(kernel.packet[6], kernel.packet[7]);
1779
1780 kernel.packet[0] = _mm512_castpd_ps(_mm512_unpacklo_pd(_mm512_castps_pd(T0), _mm512_castps_pd(T2)));
1781 kernel.packet[1] = _mm512_castpd_ps(_mm512_unpackhi_pd(_mm512_castps_pd(T0), _mm512_castps_pd(T2)));
1782 kernel.packet[2] = _mm512_castpd_ps(_mm512_unpacklo_pd(_mm512_castps_pd(T1), _mm512_castps_pd(T3)));
1783 kernel.packet[3] = _mm512_castpd_ps(_mm512_unpackhi_pd(_mm512_castps_pd(T1), _mm512_castps_pd(T3)));
1784 kernel.packet[4] = _mm512_castpd_ps(_mm512_unpacklo_pd(_mm512_castps_pd(T4), _mm512_castps_pd(T6)));
1785 kernel.packet[5] = _mm512_castpd_ps(_mm512_unpackhi_pd(_mm512_castps_pd(T4), _mm512_castps_pd(T6)));
1786 kernel.packet[6] = _mm512_castpd_ps(_mm512_unpacklo_pd(_mm512_castps_pd(T5), _mm512_castps_pd(T7)));
1787 kernel.packet[7] = _mm512_castpd_ps(_mm512_unpackhi_pd(_mm512_castps_pd(T5), _mm512_castps_pd(T7)));
1788
1789 T0 = _mm512_shuffle_f32x4(kernel.packet[0], kernel.packet[4], 0x44);
1790 T1 = _mm512_shuffle_f32x4(kernel.packet[0], kernel.packet[4], 0xee);
1791 T2 = _mm512_shuffle_f32x4(kernel.packet[1], kernel.packet[5], 0x44);
1792 T3 = _mm512_shuffle_f32x4(kernel.packet[1], kernel.packet[5], 0xee);
1793 T4 = _mm512_shuffle_f32x4(kernel.packet[2], kernel.packet[6], 0x44);
1794 T5 = _mm512_shuffle_f32x4(kernel.packet[2], kernel.packet[6], 0xee);
1795 T6 = _mm512_shuffle_f32x4(kernel.packet[3], kernel.packet[7], 0x44);
1796 T7 = _mm512_shuffle_f32x4(kernel.packet[3], kernel.packet[7], 0xee);
1797
1798 kernel.packet[0] = _mm512_shuffle_f32x4(T0, T2, 0x88);
1799 kernel.packet[2] = _mm512_shuffle_f32x4(T0, T2, 0xdd);
1800 kernel.packet[1] = _mm512_shuffle_f32x4(T4, T6, 0x88);
1801 kernel.packet[3] = _mm512_shuffle_f32x4(T4, T6, 0xdd);
1802 kernel.packet[4] = _mm512_shuffle_f32x4(T1, T3, 0x88);
1803 kernel.packet[6] = _mm512_shuffle_f32x4(T1, T3, 0xdd);
1804 kernel.packet[5] = _mm512_shuffle_f32x4(T5, T7, 0x88);
1805 kernel.packet[7] = _mm512_shuffle_f32x4(T5, T7, 0xdd);
1806}
1807
1808EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet16f, 4>& kernel) {
1809 __m512 T0 = _mm512_unpacklo_ps(kernel.packet[0], kernel.packet[1]);
1810 __m512 T1 = _mm512_unpackhi_ps(kernel.packet[0], kernel.packet[1]);
1811 __m512 T2 = _mm512_unpacklo_ps(kernel.packet[2], kernel.packet[3]);
1812 __m512 T3 = _mm512_unpackhi_ps(kernel.packet[2], kernel.packet[3]);
1813
1814 __m512 S0 = _mm512_shuffle_ps(T0, T2, _MM_SHUFFLE(1, 0, 1, 0));
1815 __m512 S1 = _mm512_shuffle_ps(T0, T2, _MM_SHUFFLE(3, 2, 3, 2));
1816 __m512 S2 = _mm512_shuffle_ps(T1, T3, _MM_SHUFFLE(1, 0, 1, 0));
1817 __m512 S3 = _mm512_shuffle_ps(T1, T3, _MM_SHUFFLE(3, 2, 3, 2));
1818
1819 EIGEN_EXTRACT_8f_FROM_16f(S0, S0);
1820 EIGEN_EXTRACT_8f_FROM_16f(S1, S1);
1821 EIGEN_EXTRACT_8f_FROM_16f(S2, S2);
1822 EIGEN_EXTRACT_8f_FROM_16f(S3, S3);
1823
1824 PacketBlock<Packet8f, 8> tmp;
1825
1826 tmp.packet[0] = _mm256_permute2f128_ps(S0_0, S1_0, 0x20);
1827 tmp.packet[1] = _mm256_permute2f128_ps(S2_0, S3_0, 0x20);
1828 tmp.packet[2] = _mm256_permute2f128_ps(S0_0, S1_0, 0x31);
1829 tmp.packet[3] = _mm256_permute2f128_ps(S2_0, S3_0, 0x31);
1830
1831 tmp.packet[4] = _mm256_permute2f128_ps(S0_1, S1_1, 0x20);
1832 tmp.packet[5] = _mm256_permute2f128_ps(S2_1, S3_1, 0x20);
1833 tmp.packet[6] = _mm256_permute2f128_ps(S0_1, S1_1, 0x31);
1834 tmp.packet[7] = _mm256_permute2f128_ps(S2_1, S3_1, 0x31);
1835
1836 PACK_OUTPUT_2(kernel.packet, tmp.packet, 0, 1);
1837 PACK_OUTPUT_2(kernel.packet, tmp.packet, 1, 1);
1838 PACK_OUTPUT_2(kernel.packet, tmp.packet, 2, 1);
1839 PACK_OUTPUT_2(kernel.packet, tmp.packet, 3, 1);
1840}
1841
1842#define PACK_OUTPUT_SQ_D(OUTPUT, INPUT, INDEX, STRIDE) \
1843 OUTPUT[INDEX] = _mm512_insertf64x4(OUTPUT[INDEX], INPUT[INDEX], 0); \
1844 OUTPUT[INDEX] = _mm512_insertf64x4(OUTPUT[INDEX], INPUT[INDEX + STRIDE], 1);
1845
1846#define PACK_OUTPUT_D(OUTPUT, INPUT, INDEX, STRIDE) \
1847 OUTPUT[INDEX] = _mm512_insertf64x4(OUTPUT[INDEX], INPUT[(2 * INDEX)], 0); \
1848 OUTPUT[INDEX] = _mm512_insertf64x4(OUTPUT[INDEX], INPUT[(2 * INDEX) + STRIDE], 1);
1849
1850#define PACK_OUTPUT_L(OUTPUT, INPUT, INDEX, STRIDE) \
1851 OUTPUT[INDEX] = _mm512_inserti64x4(OUTPUT[INDEX], INPUT[(2 * INDEX)], 0); \
1852 OUTPUT[INDEX] = _mm512_inserti64x4(OUTPUT[INDEX], INPUT[(2 * INDEX) + STRIDE], 1);
1853
1854EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet8d, 4>& kernel) {
1855 __m512d T0 = _mm512_shuffle_pd(kernel.packet[0], kernel.packet[1], 0);
1856 __m512d T1 = _mm512_shuffle_pd(kernel.packet[0], kernel.packet[1], 0xff);
1857 __m512d T2 = _mm512_shuffle_pd(kernel.packet[2], kernel.packet[3], 0);
1858 __m512d T3 = _mm512_shuffle_pd(kernel.packet[2], kernel.packet[3], 0xff);
1859
1860 PacketBlock<Packet4d, 8> tmp;
1861
1862 tmp.packet[0] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T0, 0), _mm512_extractf64x4_pd(T2, 0), 0x20);
1863 tmp.packet[1] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T1, 0), _mm512_extractf64x4_pd(T3, 0), 0x20);
1864 tmp.packet[2] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T0, 0), _mm512_extractf64x4_pd(T2, 0), 0x31);
1865 tmp.packet[3] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T1, 0), _mm512_extractf64x4_pd(T3, 0), 0x31);
1866
1867 tmp.packet[4] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T0, 1), _mm512_extractf64x4_pd(T2, 1), 0x20);
1868 tmp.packet[5] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T1, 1), _mm512_extractf64x4_pd(T3, 1), 0x20);
1869 tmp.packet[6] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T0, 1), _mm512_extractf64x4_pd(T2, 1), 0x31);
1870 tmp.packet[7] = _mm256_permute2f128_pd(_mm512_extractf64x4_pd(T1, 1), _mm512_extractf64x4_pd(T3, 1), 0x31);
1871
1872 PACK_OUTPUT_D(kernel.packet, tmp.packet, 0, 1);
1873 PACK_OUTPUT_D(kernel.packet, tmp.packet, 1, 1);
1874 PACK_OUTPUT_D(kernel.packet, tmp.packet, 2, 1);
1875 PACK_OUTPUT_D(kernel.packet, tmp.packet, 3, 1);
1876}
1877
1878EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet8d, 8>& kernel) {
1879 __m512d T0 = _mm512_unpacklo_pd(kernel.packet[0], kernel.packet[1]);
1880 __m512d T1 = _mm512_unpackhi_pd(kernel.packet[0], kernel.packet[1]);
1881 __m512d T2 = _mm512_unpacklo_pd(kernel.packet[2], kernel.packet[3]);
1882 __m512d T3 = _mm512_unpackhi_pd(kernel.packet[2], kernel.packet[3]);
1883 __m512d T4 = _mm512_unpacklo_pd(kernel.packet[4], kernel.packet[5]);
1884 __m512d T5 = _mm512_unpackhi_pd(kernel.packet[4], kernel.packet[5]);
1885 __m512d T6 = _mm512_unpacklo_pd(kernel.packet[6], kernel.packet[7]);
1886 __m512d T7 = _mm512_unpackhi_pd(kernel.packet[6], kernel.packet[7]);
1887
1888 kernel.packet[0] = _mm512_permutex_pd(T2, 0x4E);
1889 kernel.packet[0] = _mm512_mask_blend_pd(0xCC, T0, kernel.packet[0]);
1890 kernel.packet[2] = _mm512_permutex_pd(T0, 0x4E);
1891 kernel.packet[2] = _mm512_mask_blend_pd(0xCC, kernel.packet[2], T2);
1892 kernel.packet[1] = _mm512_permutex_pd(T3, 0x4E);
1893 kernel.packet[1] = _mm512_mask_blend_pd(0xCC, T1, kernel.packet[1]);
1894 kernel.packet[3] = _mm512_permutex_pd(T1, 0x4E);
1895 kernel.packet[3] = _mm512_mask_blend_pd(0xCC, kernel.packet[3], T3);
1896 kernel.packet[4] = _mm512_permutex_pd(T6, 0x4E);
1897 kernel.packet[4] = _mm512_mask_blend_pd(0xCC, T4, kernel.packet[4]);
1898 kernel.packet[6] = _mm512_permutex_pd(T4, 0x4E);
1899 kernel.packet[6] = _mm512_mask_blend_pd(0xCC, kernel.packet[6], T6);
1900 kernel.packet[5] = _mm512_permutex_pd(T7, 0x4E);
1901 kernel.packet[5] = _mm512_mask_blend_pd(0xCC, T5, kernel.packet[5]);
1902 kernel.packet[7] = _mm512_permutex_pd(T5, 0x4E);
1903 kernel.packet[7] = _mm512_mask_blend_pd(0xCC, kernel.packet[7], T7);
1904
1905 T0 = _mm512_shuffle_f64x2(kernel.packet[4], kernel.packet[4], 0x4E);
1906 T0 = _mm512_mask_blend_pd(0xF0, kernel.packet[0], T0);
1907 T4 = _mm512_shuffle_f64x2(kernel.packet[0], kernel.packet[0], 0x4E);
1908 T4 = _mm512_mask_blend_pd(0xF0, T4, kernel.packet[4]);
1909 T1 = _mm512_shuffle_f64x2(kernel.packet[5], kernel.packet[5], 0x4E);
1910 T1 = _mm512_mask_blend_pd(0xF0, kernel.packet[1], T1);
1911 T5 = _mm512_shuffle_f64x2(kernel.packet[1], kernel.packet[1], 0x4E);
1912 T5 = _mm512_mask_blend_pd(0xF0, T5, kernel.packet[5]);
1913 T2 = _mm512_shuffle_f64x2(kernel.packet[6], kernel.packet[6], 0x4E);
1914 T2 = _mm512_mask_blend_pd(0xF0, kernel.packet[2], T2);
1915 T6 = _mm512_shuffle_f64x2(kernel.packet[2], kernel.packet[2], 0x4E);
1916 T6 = _mm512_mask_blend_pd(0xF0, T6, kernel.packet[6]);
1917 T3 = _mm512_shuffle_f64x2(kernel.packet[7], kernel.packet[7], 0x4E);
1918 T3 = _mm512_mask_blend_pd(0xF0, kernel.packet[3], T3);
1919 T7 = _mm512_shuffle_f64x2(kernel.packet[3], kernel.packet[3], 0x4E);
1920 T7 = _mm512_mask_blend_pd(0xF0, T7, kernel.packet[7]);
1921
1922 kernel.packet[0] = T0;
1923 kernel.packet[1] = T1;
1924 kernel.packet[2] = T2;
1925 kernel.packet[3] = T3;
1926 kernel.packet[4] = T4;
1927 kernel.packet[5] = T5;
1928 kernel.packet[6] = T6;
1929 kernel.packet[7] = T7;
1930}
1931
1932EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet8l, 4>& kernel) {
1933 __m512i T0 = _mm512_castpd_si512(
1934 _mm512_shuffle_pd(_mm512_castsi512_pd(kernel.packet[0]), _mm512_castsi512_pd(kernel.packet[1]), 0));
1935 __m512i T1 = _mm512_castpd_si512(
1936 _mm512_shuffle_pd(_mm512_castsi512_pd(kernel.packet[0]), _mm512_castsi512_pd(kernel.packet[1]), 0xff));
1937 __m512i T2 = _mm512_castpd_si512(
1938 _mm512_shuffle_pd(_mm512_castsi512_pd(kernel.packet[2]), _mm512_castsi512_pd(kernel.packet[3]), 0));
1939 __m512i T3 = _mm512_castpd_si512(
1940 _mm512_shuffle_pd(_mm512_castsi512_pd(kernel.packet[2]), _mm512_castsi512_pd(kernel.packet[3]), 0xff));
1941
1942 PacketBlock<Packet4l, 8> tmp;
1943
1944 tmp.packet[0] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T0, 0), _mm512_extracti64x4_epi64(T2, 0), 0x20);
1945 tmp.packet[1] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T1, 0), _mm512_extracti64x4_epi64(T3, 0), 0x20);
1946 tmp.packet[2] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T0, 0), _mm512_extracti64x4_epi64(T2, 0), 0x31);
1947 tmp.packet[3] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T1, 0), _mm512_extracti64x4_epi64(T3, 0), 0x31);
1948
1949 tmp.packet[4] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T0, 1), _mm512_extracti64x4_epi64(T2, 1), 0x20);
1950 tmp.packet[5] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T1, 1), _mm512_extracti64x4_epi64(T3, 1), 0x20);
1951 tmp.packet[6] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T0, 1), _mm512_extracti64x4_epi64(T2, 1), 0x31);
1952 tmp.packet[7] = _mm256_permute2x128_si256(_mm512_extracti64x4_epi64(T1, 1), _mm512_extracti64x4_epi64(T3, 1), 0x31);
1953
1954 PACK_OUTPUT_L(kernel.packet, tmp.packet, 0, 1);
1955 PACK_OUTPUT_L(kernel.packet, tmp.packet, 1, 1);
1956 PACK_OUTPUT_L(kernel.packet, tmp.packet, 2, 1);
1957 PACK_OUTPUT_L(kernel.packet, tmp.packet, 3, 1);
1958}
1959
1960EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet8l, 8>& kernel) {
1961 __m512i T0 = _mm512_unpacklo_epi64(kernel.packet[0], kernel.packet[1]);
1962 __m512i T1 = _mm512_unpackhi_epi64(kernel.packet[0], kernel.packet[1]);
1963 __m512i T2 = _mm512_unpacklo_epi64(kernel.packet[2], kernel.packet[3]);
1964 __m512i T3 = _mm512_unpackhi_epi64(kernel.packet[2], kernel.packet[3]);
1965 __m512i T4 = _mm512_unpacklo_epi64(kernel.packet[4], kernel.packet[5]);
1966 __m512i T5 = _mm512_unpackhi_epi64(kernel.packet[4], kernel.packet[5]);
1967 __m512i T6 = _mm512_unpacklo_epi64(kernel.packet[6], kernel.packet[7]);
1968 __m512i T7 = _mm512_unpackhi_epi64(kernel.packet[6], kernel.packet[7]);
1969
1970 kernel.packet[0] = _mm512_permutex_epi64(T2, 0x4E);
1971 kernel.packet[0] = _mm512_mask_blend_epi64(0xCC, T0, kernel.packet[0]);
1972 kernel.packet[2] = _mm512_permutex_epi64(T0, 0x4E);
1973 kernel.packet[2] = _mm512_mask_blend_epi64(0xCC, kernel.packet[2], T2);
1974 kernel.packet[1] = _mm512_permutex_epi64(T3, 0x4E);
1975 kernel.packet[1] = _mm512_mask_blend_epi64(0xCC, T1, kernel.packet[1]);
1976 kernel.packet[3] = _mm512_permutex_epi64(T1, 0x4E);
1977 kernel.packet[3] = _mm512_mask_blend_epi64(0xCC, kernel.packet[3], T3);
1978 kernel.packet[4] = _mm512_permutex_epi64(T6, 0x4E);
1979 kernel.packet[4] = _mm512_mask_blend_epi64(0xCC, T4, kernel.packet[4]);
1980 kernel.packet[6] = _mm512_permutex_epi64(T4, 0x4E);
1981 kernel.packet[6] = _mm512_mask_blend_epi64(0xCC, kernel.packet[6], T6);
1982 kernel.packet[5] = _mm512_permutex_epi64(T7, 0x4E);
1983 kernel.packet[5] = _mm512_mask_blend_epi64(0xCC, T5, kernel.packet[5]);
1984 kernel.packet[7] = _mm512_permutex_epi64(T5, 0x4E);
1985 kernel.packet[7] = _mm512_mask_blend_epi64(0xCC, kernel.packet[7], T7);
1986
1987 T0 = _mm512_shuffle_i64x2(kernel.packet[4], kernel.packet[4], 0x4E);
1988 T0 = _mm512_mask_blend_epi64(0xF0, kernel.packet[0], T0);
1989 T4 = _mm512_shuffle_i64x2(kernel.packet[0], kernel.packet[0], 0x4E);
1990 T4 = _mm512_mask_blend_epi64(0xF0, T4, kernel.packet[4]);
1991 T1 = _mm512_shuffle_i64x2(kernel.packet[5], kernel.packet[5], 0x4E);
1992 T1 = _mm512_mask_blend_epi64(0xF0, kernel.packet[1], T1);
1993 T5 = _mm512_shuffle_i64x2(kernel.packet[1], kernel.packet[1], 0x4E);
1994 T5 = _mm512_mask_blend_epi64(0xF0, T5, kernel.packet[5]);
1995 T2 = _mm512_shuffle_i64x2(kernel.packet[6], kernel.packet[6], 0x4E);
1996 T2 = _mm512_mask_blend_epi64(0xF0, kernel.packet[2], T2);
1997 T6 = _mm512_shuffle_i64x2(kernel.packet[2], kernel.packet[2], 0x4E);
1998 T6 = _mm512_mask_blend_epi64(0xF0, T6, kernel.packet[6]);
1999 T3 = _mm512_shuffle_i64x2(kernel.packet[7], kernel.packet[7], 0x4E);
2000 T3 = _mm512_mask_blend_epi64(0xF0, kernel.packet[3], T3);
2001 T7 = _mm512_shuffle_i64x2(kernel.packet[3], kernel.packet[3], 0x4E);
2002 T7 = _mm512_mask_blend_epi64(0xF0, T7, kernel.packet[7]);
2003
2004 kernel.packet[0] = T0;
2005 kernel.packet[1] = T1;
2006 kernel.packet[2] = T2;
2007 kernel.packet[3] = T3;
2008 kernel.packet[4] = T4;
2009 kernel.packet[5] = T5;
2010 kernel.packet[6] = T6;
2011 kernel.packet[7] = T7;
2012}
2013
2014#define PACK_OUTPUT_I32(OUTPUT, INPUT, INDEX, STRIDE) \
2015 EIGEN_INSERT_8i_INTO_16i(OUTPUT[INDEX], INPUT[INDEX], INPUT[INDEX + STRIDE]);
2016
2017#define PACK_OUTPUT_I32_2(OUTPUT, INPUT, INDEX, STRIDE) \
2018 EIGEN_INSERT_8i_INTO_16i(OUTPUT[INDEX], INPUT[2 * INDEX], INPUT[2 * INDEX + STRIDE]);
2019
2020#define SHUFFLE_EPI32(A, B, M) _mm512_castps_si512(_mm512_shuffle_ps(_mm512_castsi512_ps(A), _mm512_castsi512_ps(B), M))
2021
2022EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet16i, 16>& kernel) {
2023 __m512i T0 = _mm512_unpacklo_epi32(kernel.packet[0], kernel.packet[1]);
2024 __m512i T1 = _mm512_unpackhi_epi32(kernel.packet[0], kernel.packet[1]);
2025 __m512i T2 = _mm512_unpacklo_epi32(kernel.packet[2], kernel.packet[3]);
2026 __m512i T3 = _mm512_unpackhi_epi32(kernel.packet[2], kernel.packet[3]);
2027 __m512i T4 = _mm512_unpacklo_epi32(kernel.packet[4], kernel.packet[5]);
2028 __m512i T5 = _mm512_unpackhi_epi32(kernel.packet[4], kernel.packet[5]);
2029 __m512i T6 = _mm512_unpacklo_epi32(kernel.packet[6], kernel.packet[7]);
2030 __m512i T7 = _mm512_unpackhi_epi32(kernel.packet[6], kernel.packet[7]);
2031 __m512i T8 = _mm512_unpacklo_epi32(kernel.packet[8], kernel.packet[9]);
2032 __m512i T9 = _mm512_unpackhi_epi32(kernel.packet[8], kernel.packet[9]);
2033 __m512i T10 = _mm512_unpacklo_epi32(kernel.packet[10], kernel.packet[11]);
2034 __m512i T11 = _mm512_unpackhi_epi32(kernel.packet[10], kernel.packet[11]);
2035 __m512i T12 = _mm512_unpacklo_epi32(kernel.packet[12], kernel.packet[13]);
2036 __m512i T13 = _mm512_unpackhi_epi32(kernel.packet[12], kernel.packet[13]);
2037 __m512i T14 = _mm512_unpacklo_epi32(kernel.packet[14], kernel.packet[15]);
2038 __m512i T15 = _mm512_unpackhi_epi32(kernel.packet[14], kernel.packet[15]);
2039 __m512i S0 = SHUFFLE_EPI32(T0, T2, _MM_SHUFFLE(1, 0, 1, 0));
2040 __m512i S1 = SHUFFLE_EPI32(T0, T2, _MM_SHUFFLE(3, 2, 3, 2));
2041 __m512i S2 = SHUFFLE_EPI32(T1, T3, _MM_SHUFFLE(1, 0, 1, 0));
2042 __m512i S3 = SHUFFLE_EPI32(T1, T3, _MM_SHUFFLE(3, 2, 3, 2));
2043 __m512i S4 = SHUFFLE_EPI32(T4, T6, _MM_SHUFFLE(1, 0, 1, 0));
2044 __m512i S5 = SHUFFLE_EPI32(T4, T6, _MM_SHUFFLE(3, 2, 3, 2));
2045 __m512i S6 = SHUFFLE_EPI32(T5, T7, _MM_SHUFFLE(1, 0, 1, 0));
2046 __m512i S7 = SHUFFLE_EPI32(T5, T7, _MM_SHUFFLE(3, 2, 3, 2));
2047 __m512i S8 = SHUFFLE_EPI32(T8, T10, _MM_SHUFFLE(1, 0, 1, 0));
2048 __m512i S9 = SHUFFLE_EPI32(T8, T10, _MM_SHUFFLE(3, 2, 3, 2));
2049 __m512i S10 = SHUFFLE_EPI32(T9, T11, _MM_SHUFFLE(1, 0, 1, 0));
2050 __m512i S11 = SHUFFLE_EPI32(T9, T11, _MM_SHUFFLE(3, 2, 3, 2));
2051 __m512i S12 = SHUFFLE_EPI32(T12, T14, _MM_SHUFFLE(1, 0, 1, 0));
2052 __m512i S13 = SHUFFLE_EPI32(T12, T14, _MM_SHUFFLE(3, 2, 3, 2));
2053 __m512i S14 = SHUFFLE_EPI32(T13, T15, _MM_SHUFFLE(1, 0, 1, 0));
2054 __m512i S15 = SHUFFLE_EPI32(T13, T15, _MM_SHUFFLE(3, 2, 3, 2));
2055
2056 EIGEN_EXTRACT_8i_FROM_16i(S0, S0);
2057 EIGEN_EXTRACT_8i_FROM_16i(S1, S1);
2058 EIGEN_EXTRACT_8i_FROM_16i(S2, S2);
2059 EIGEN_EXTRACT_8i_FROM_16i(S3, S3);
2060 EIGEN_EXTRACT_8i_FROM_16i(S4, S4);
2061 EIGEN_EXTRACT_8i_FROM_16i(S5, S5);
2062 EIGEN_EXTRACT_8i_FROM_16i(S6, S6);
2063 EIGEN_EXTRACT_8i_FROM_16i(S7, S7);
2064 EIGEN_EXTRACT_8i_FROM_16i(S8, S8);
2065 EIGEN_EXTRACT_8i_FROM_16i(S9, S9);
2066 EIGEN_EXTRACT_8i_FROM_16i(S10, S10);
2067 EIGEN_EXTRACT_8i_FROM_16i(S11, S11);
2068 EIGEN_EXTRACT_8i_FROM_16i(S12, S12);
2069 EIGEN_EXTRACT_8i_FROM_16i(S13, S13);
2070 EIGEN_EXTRACT_8i_FROM_16i(S14, S14);
2071 EIGEN_EXTRACT_8i_FROM_16i(S15, S15);
2072
2073 PacketBlock<Packet8i, 32> tmp;
2074
2075 tmp.packet[0] = _mm256_permute2f128_si256(S0_0, S4_0, 0x20);
2076 tmp.packet[1] = _mm256_permute2f128_si256(S1_0, S5_0, 0x20);
2077 tmp.packet[2] = _mm256_permute2f128_si256(S2_0, S6_0, 0x20);
2078 tmp.packet[3] = _mm256_permute2f128_si256(S3_0, S7_0, 0x20);
2079 tmp.packet[4] = _mm256_permute2f128_si256(S0_0, S4_0, 0x31);
2080 tmp.packet[5] = _mm256_permute2f128_si256(S1_0, S5_0, 0x31);
2081 tmp.packet[6] = _mm256_permute2f128_si256(S2_0, S6_0, 0x31);
2082 tmp.packet[7] = _mm256_permute2f128_si256(S3_0, S7_0, 0x31);
2083
2084 tmp.packet[8] = _mm256_permute2f128_si256(S0_1, S4_1, 0x20);
2085 tmp.packet[9] = _mm256_permute2f128_si256(S1_1, S5_1, 0x20);
2086 tmp.packet[10] = _mm256_permute2f128_si256(S2_1, S6_1, 0x20);
2087 tmp.packet[11] = _mm256_permute2f128_si256(S3_1, S7_1, 0x20);
2088 tmp.packet[12] = _mm256_permute2f128_si256(S0_1, S4_1, 0x31);
2089 tmp.packet[13] = _mm256_permute2f128_si256(S1_1, S5_1, 0x31);
2090 tmp.packet[14] = _mm256_permute2f128_si256(S2_1, S6_1, 0x31);
2091 tmp.packet[15] = _mm256_permute2f128_si256(S3_1, S7_1, 0x31);
2092
2093 // Second set of _m256 outputs
2094 tmp.packet[16] = _mm256_permute2f128_si256(S8_0, S12_0, 0x20);
2095 tmp.packet[17] = _mm256_permute2f128_si256(S9_0, S13_0, 0x20);
2096 tmp.packet[18] = _mm256_permute2f128_si256(S10_0, S14_0, 0x20);
2097 tmp.packet[19] = _mm256_permute2f128_si256(S11_0, S15_0, 0x20);
2098 tmp.packet[20] = _mm256_permute2f128_si256(S8_0, S12_0, 0x31);
2099 tmp.packet[21] = _mm256_permute2f128_si256(S9_0, S13_0, 0x31);
2100 tmp.packet[22] = _mm256_permute2f128_si256(S10_0, S14_0, 0x31);
2101 tmp.packet[23] = _mm256_permute2f128_si256(S11_0, S15_0, 0x31);
2102
2103 tmp.packet[24] = _mm256_permute2f128_si256(S8_1, S12_1, 0x20);
2104 tmp.packet[25] = _mm256_permute2f128_si256(S9_1, S13_1, 0x20);
2105 tmp.packet[26] = _mm256_permute2f128_si256(S10_1, S14_1, 0x20);
2106 tmp.packet[27] = _mm256_permute2f128_si256(S11_1, S15_1, 0x20);
2107 tmp.packet[28] = _mm256_permute2f128_si256(S8_1, S12_1, 0x31);
2108 tmp.packet[29] = _mm256_permute2f128_si256(S9_1, S13_1, 0x31);
2109 tmp.packet[30] = _mm256_permute2f128_si256(S10_1, S14_1, 0x31);
2110 tmp.packet[31] = _mm256_permute2f128_si256(S11_1, S15_1, 0x31);
2111
2112 // Pack them into the output
2113 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 0, 16);
2114 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 1, 16);
2115 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 2, 16);
2116 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 3, 16);
2117
2118 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 4, 16);
2119 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 5, 16);
2120 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 6, 16);
2121 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 7, 16);
2122
2123 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 8, 16);
2124 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 9, 16);
2125 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 10, 16);
2126 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 11, 16);
2127
2128 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 12, 16);
2129 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 13, 16);
2130 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 14, 16);
2131 PACK_OUTPUT_I32(kernel.packet, tmp.packet, 15, 16);
2132}
2133
2134EIGEN_DEVICE_FUNC inline void ptranspose(PacketBlock<Packet16i, 4>& kernel) {
2135 __m512i T0 = _mm512_unpacklo_epi32(kernel.packet[0], kernel.packet[1]);
2136 __m512i T1 = _mm512_unpackhi_epi32(kernel.packet[0], kernel.packet[1]);
2137 __m512i T2 = _mm512_unpacklo_epi32(kernel.packet[2], kernel.packet[3]);
2138 __m512i T3 = _mm512_unpackhi_epi32(kernel.packet[2], kernel.packet[3]);
2139
2140 __m512i S0 = SHUFFLE_EPI32(T0, T2, _MM_SHUFFLE(1, 0, 1, 0));
2141 __m512i S1 = SHUFFLE_EPI32(T0, T2, _MM_SHUFFLE(3, 2, 3, 2));
2142 __m512i S2 = SHUFFLE_EPI32(T1, T3, _MM_SHUFFLE(1, 0, 1, 0));
2143 __m512i S3 = SHUFFLE_EPI32(T1, T3, _MM_SHUFFLE(3, 2, 3, 2));
2144
2145 EIGEN_EXTRACT_8i_FROM_16i(S0, S0);
2146 EIGEN_EXTRACT_8i_FROM_16i(S1, S1);
2147 EIGEN_EXTRACT_8i_FROM_16i(S2, S2);
2148 EIGEN_EXTRACT_8i_FROM_16i(S3, S3);
2149
2150 PacketBlock<Packet8i, 8> tmp;
2151
2152 tmp.packet[0] = _mm256_permute2f128_si256(S0_0, S1_0, 0x20);
2153 tmp.packet[1] = _mm256_permute2f128_si256(S2_0, S3_0, 0x20);
2154 tmp.packet[2] = _mm256_permute2f128_si256(S0_0, S1_0, 0x31);
2155 tmp.packet[3] = _mm256_permute2f128_si256(S2_0, S3_0, 0x31);
2156
2157 tmp.packet[4] = _mm256_permute2f128_si256(S0_1, S1_1, 0x20);
2158 tmp.packet[5] = _mm256_permute2f128_si256(S2_1, S3_1, 0x20);
2159 tmp.packet[6] = _mm256_permute2f128_si256(S0_1, S1_1, 0x31);
2160 tmp.packet[7] = _mm256_permute2f128_si256(S2_1, S3_1, 0x31);
2161
2162 PACK_OUTPUT_I32_2(kernel.packet, tmp.packet, 0, 1);
2163 PACK_OUTPUT_I32_2(kernel.packet, tmp.packet, 1, 1);
2164 PACK_OUTPUT_I32_2(kernel.packet, tmp.packet, 2, 1);
2165 PACK_OUTPUT_I32_2(kernel.packet, tmp.packet, 3, 1);
2166}
2167
2168template <size_t N>
2169EIGEN_STRONG_INLINE int avx512_blend_mask(const Selector<N>& ifPacket) {
2170 alignas(__m128i) uint8_t aux[sizeof(__m128i)];
2171 for (size_t i = 0; i < N; i++) aux[i] = static_cast<uint8_t>(ifPacket.select[i]);
2172 __m128i paux = _mm_sub_epi8(_mm_setzero_si128(), _mm_load_si128(reinterpret_cast<const __m128i*>(aux)));
2173 return _mm_movemask_epi8(paux);
2174}
2175
2176template <>
2177EIGEN_STRONG_INLINE Packet16f pblend(const Selector<16>& ifPacket, const Packet16f& thenPacket,
2178 const Packet16f& elsePacket) {
2179 __mmask16 m = avx512_blend_mask(ifPacket);
2180 return _mm512_mask_blend_ps(m, elsePacket, thenPacket);
2181}
2182template <>
2183EIGEN_STRONG_INLINE Packet8d pblend(const Selector<8>& ifPacket, const Packet8d& thenPacket,
2184 const Packet8d& elsePacket) {
2185 __mmask8 m = avx512_blend_mask(ifPacket);
2186 return _mm512_mask_blend_pd(m, elsePacket, thenPacket);
2187}
2188
2189// Packet math for Eigen::half
2190template <>
2191EIGEN_STRONG_INLINE Packet16h pset1<Packet16h>(const Eigen::half& from) {
2192 return _mm256_set1_epi16(from.x);
2193}
2194
2195template <>
2196EIGEN_STRONG_INLINE Eigen::half pfirst<Packet16h>(const Packet16h& from) {
2197 return half_impl::raw_uint16_to_half(static_cast<unsigned short>(_mm256_extract_epi16(from, 0)));
2198}
2199
2200template <>
2201EIGEN_STRONG_INLINE Packet16h pload<Packet16h>(const Eigen::half* from) {
2202 return _mm256_load_si256(reinterpret_cast<const __m256i*>(from));
2203}
2204
2205template <>
2206EIGEN_STRONG_INLINE Packet16h ploadu<Packet16h>(const Eigen::half* from) {
2207 return _mm256_loadu_si256(reinterpret_cast<const __m256i*>(from));
2208}
2209
2210template <>
2211EIGEN_STRONG_INLINE void pstore<half>(Eigen::half* to, const Packet16h& from) {
2212 // (void*) -> workaround clang warning:
2213 // cast from 'Eigen::half *' to '__m256i *' increases required alignment from 2 to 32
2214 _mm256_store_si256((__m256i*)(void*)to, from);
2215}
2216
2217template <>
2218EIGEN_STRONG_INLINE void pstoreu<half>(Eigen::half* to, const Packet16h& from) {
2219 // (void*) -> workaround clang warning:
2220 // cast from 'Eigen::half *' to '__m256i *' increases required alignment from 2 to 32
2221 _mm256_storeu_si256((__m256i*)(void*)to, from);
2222}
2223
2224template <>
2225EIGEN_STRONG_INLINE Packet16h ploaddup<Packet16h>(const Eigen::half* from) {
2226 unsigned short a = from[0].x;
2227 unsigned short b = from[1].x;
2228 unsigned short c = from[2].x;
2229 unsigned short d = from[3].x;
2230 unsigned short e = from[4].x;
2231 unsigned short f = from[5].x;
2232 unsigned short g = from[6].x;
2233 unsigned short h = from[7].x;
2234 return _mm256_set_epi16(h, h, g, g, f, f, e, e, d, d, c, c, b, b, a, a);
2235}
2236
2237template <>
2238EIGEN_STRONG_INLINE Packet16h ploadquad(const Eigen::half* from) {
2239 unsigned short a = from[0].x;
2240 unsigned short b = from[1].x;
2241 unsigned short c = from[2].x;
2242 unsigned short d = from[3].x;
2243 return _mm256_set_epi16(d, d, d, d, c, c, c, c, b, b, b, b, a, a, a, a);
2244}
2245
2246EIGEN_STRONG_INLINE Packet16f half2float(const Packet16h& a) { return _mm512_cvtph_ps(a); }
2247
2248EIGEN_STRONG_INLINE Packet16h float2half(const Packet16f& a) {
2249 return _mm512_cvtps_ph(a, _MM_FROUND_TO_NEAREST_INT | _MM_FROUND_NO_EXC);
2250}
2251
2252template <>
2253EIGEN_STRONG_INLINE Packet16h ptrue(const Packet16h& a) {
2254 return Packet16h(ptrue(Packet8i(a)));
2255}
2256
2257template <>
2258EIGEN_STRONG_INLINE Packet16h pabs(const Packet16h& a) {
2259 const __m256i sign_mask = _mm256_set1_epi16(static_cast<numext::uint16_t>(0x8000));
2260 return _mm256_andnot_si256(sign_mask, a);
2261}
2262
2263template <>
2264EIGEN_STRONG_INLINE Packet16h pmin<Packet16h>(const Packet16h& a, const Packet16h& b) {
2265 return float2half(pmin<Packet16f>(half2float(a), half2float(b)));
2266}
2267
2268template <>
2269EIGEN_STRONG_INLINE Packet16h pmax<Packet16h>(const Packet16h& a, const Packet16h& b) {
2270 return float2half(pmax<Packet16f>(half2float(a), half2float(b)));
2271}
2272
2273template <>
2274EIGEN_STRONG_INLINE Packet16h plset<Packet16h>(const half& a) {
2275 return float2half(plset<Packet16f>(static_cast<float>(a)));
2276}
2277
2278template <>
2279EIGEN_STRONG_INLINE Packet16h por(const Packet16h& a, const Packet16h& b) {
2280 // in some cases Packet8i is a wrapper around __m256i, so we need to
2281 // cast to Packet8i to call the correct overload.
2282 return Packet16h(por(Packet8i(a), Packet8i(b)));
2283}
2284template <>
2285EIGEN_STRONG_INLINE Packet16h pxor(const Packet16h& a, const Packet16h& b) {
2286 return Packet16h(pxor(Packet8i(a), Packet8i(b)));
2287}
2288template <>
2289EIGEN_STRONG_INLINE Packet16h pand(const Packet16h& a, const Packet16h& b) {
2290 return Packet16h(pand(Packet8i(a), Packet8i(b)));
2291}
2292template <>
2293EIGEN_STRONG_INLINE Packet16h pandnot(const Packet16h& a, const Packet16h& b) {
2294 return Packet16h(pandnot(Packet8i(a), Packet8i(b)));
2295}
2296
2297template <>
2298EIGEN_STRONG_INLINE Packet16h pselect(const Packet16h& mask, const Packet16h& a, const Packet16h& b) {
2299 return _mm256_blendv_epi8(b, a, mask);
2300}
2301
2302template <>
2303EIGEN_STRONG_INLINE Packet16h pround<Packet16h>(const Packet16h& a) {
2304 return float2half(pround<Packet16f>(half2float(a)));
2305}
2306
2307template <>
2308EIGEN_STRONG_INLINE Packet16h print<Packet16h>(const Packet16h& a) {
2309 return float2half(print<Packet16f>(half2float(a)));
2310}
2311
2312template <>
2313EIGEN_STRONG_INLINE Packet16h pceil<Packet16h>(const Packet16h& a) {
2314 return float2half(pceil<Packet16f>(half2float(a)));
2315}
2316
2317template <>
2318EIGEN_STRONG_INLINE Packet16h pfloor<Packet16h>(const Packet16h& a) {
2319 return float2half(pfloor<Packet16f>(half2float(a)));
2320}
2321
2322template <>
2323EIGEN_STRONG_INLINE Packet16h ptrunc<Packet16h>(const Packet16h& a) {
2324 return float2half(ptrunc<Packet16f>(half2float(a)));
2325}
2326
2327template <>
2328EIGEN_STRONG_INLINE Packet16h pcmp_eq(const Packet16h& a, const Packet16h& b) {
2329 Packet16f af = half2float(a);
2330 Packet16f bf = half2float(b);
2331 return Pack32To16(pcmp_eq(af, bf));
2332}
2333
2334template <>
2335EIGEN_STRONG_INLINE Packet16h pcmp_le(const Packet16h& a, const Packet16h& b) {
2336 return Pack32To16(pcmp_le(half2float(a), half2float(b)));
2337}
2338
2339template <>
2340EIGEN_STRONG_INLINE Packet16h pcmp_lt(const Packet16h& a, const Packet16h& b) {
2341 return Pack32To16(pcmp_lt(half2float(a), half2float(b)));
2342}
2343
2344template <>
2345EIGEN_STRONG_INLINE Packet16h pcmp_lt_or_nan(const Packet16h& a, const Packet16h& b) {
2346 return Pack32To16(pcmp_lt_or_nan(half2float(a), half2float(b)));
2347}
2348
2349template <>
2350EIGEN_STRONG_INLINE Packet16h pconj(const Packet16h& a) {
2351 return a;
2352}
2353
2354template <>
2355EIGEN_STRONG_INLINE Packet16h pnegate(const Packet16h& a) {
2356 Packet16h sign_mask = _mm256_set1_epi16(static_cast<unsigned short>(0x8000));
2357 return _mm256_xor_si256(a, sign_mask);
2358}
2359
2360#ifndef EIGEN_VECTORIZE_AVX512FP16
2361template <>
2362EIGEN_STRONG_INLINE Packet16h padd<Packet16h>(const Packet16h& a, const Packet16h& b) {
2363 Packet16f af = half2float(a);
2364 Packet16f bf = half2float(b);
2365 Packet16f rf = padd(af, bf);
2366 return float2half(rf);
2367}
2368
2369template <>
2370EIGEN_STRONG_INLINE Packet16h psub<Packet16h>(const Packet16h& a, const Packet16h& b) {
2371 Packet16f af = half2float(a);
2372 Packet16f bf = half2float(b);
2373 Packet16f rf = psub(af, bf);
2374 return float2half(rf);
2375}
2376
2377template <>
2378EIGEN_STRONG_INLINE Packet16h pmul<Packet16h>(const Packet16h& a, const Packet16h& b) {
2379 Packet16f af = half2float(a);
2380 Packet16f bf = half2float(b);
2381 Packet16f rf = pmul(af, bf);
2382 return float2half(rf);
2383}
2384
2385template <>
2386EIGEN_STRONG_INLINE Packet16h pdiv<Packet16h>(const Packet16h& a, const Packet16h& b) {
2387 Packet16f af = half2float(a);
2388 Packet16f bf = half2float(b);
2389 Packet16f rf = pdiv(af, bf);
2390 return float2half(rf);
2391}
2392
2393template <>
2394EIGEN_STRONG_INLINE half predux<Packet16h>(const Packet16h& from) {
2395 Packet16f from_float = half2float(from);
2396 return half(predux(from_float));
2397}
2398
2399#endif
2400
2401template <>
2402EIGEN_STRONG_INLINE Packet8h predux_half_dowto4<Packet16h>(const Packet16h& a) {
2403 Packet8h lane0 = _mm256_extractf128_si256(a, 0);
2404 Packet8h lane1 = _mm256_extractf128_si256(a, 1);
2405 return padd<Packet8h>(lane0, lane1);
2406}
2407
2408template <>
2409EIGEN_STRONG_INLINE Eigen::half predux_max<Packet16h>(const Packet16h& a) {
2410 Packet16f af = half2float(a);
2411 float reduced = predux_max<Packet16f>(af);
2412 return Eigen::half(reduced);
2413}
2414
2415template <>
2416EIGEN_STRONG_INLINE Eigen::half predux_min<Packet16h>(const Packet16h& a) {
2417 Packet16f af = half2float(a);
2418 float reduced = predux_min<Packet16f>(af);
2419 return Eigen::half(reduced);
2420}
2421
2422template <>
2423EIGEN_STRONG_INLINE half predux_mul<Packet16h>(const Packet16h& from) {
2424 Packet16f from_float = half2float(from);
2425 return half(predux_mul(from_float));
2426}
2427
2428template <>
2429EIGEN_STRONG_INLINE Packet16h preverse(const Packet16h& a) {
2430 __m128i m = _mm_setr_epi8(14, 15, 12, 13, 10, 11, 8, 9, 6, 7, 4, 5, 2, 3, 0, 1);
2431 return _mm256_insertf128_si256(_mm256_castsi128_si256(_mm_shuffle_epi8(_mm256_extractf128_si256(a, 1), m)),
2432 _mm_shuffle_epi8(_mm256_extractf128_si256(a, 0), m), 1);
2433}
2434
2435template <>
2436EIGEN_STRONG_INLINE Packet16h pgather<Eigen::half, Packet16h>(const Eigen::half* from, Index stride) {
2437 return _mm256_set_epi16(from[15 * stride].x, from[14 * stride].x, from[13 * stride].x, from[12 * stride].x,
2438 from[11 * stride].x, from[10 * stride].x, from[9 * stride].x, from[8 * stride].x,
2439 from[7 * stride].x, from[6 * stride].x, from[5 * stride].x, from[4 * stride].x,
2440 from[3 * stride].x, from[2 * stride].x, from[1 * stride].x, from[0 * stride].x);
2441}
2442
2443template <>
2444EIGEN_STRONG_INLINE void pscatter<half, Packet16h>(half* to, const Packet16h& from, Index stride) {
2445 EIGEN_ALIGN64 half aux[16];
2446 pstore(aux, from);
2447 to[stride * 0] = aux[0];
2448 to[stride * 1] = aux[1];
2449 to[stride * 2] = aux[2];
2450 to[stride * 3] = aux[3];
2451 to[stride * 4] = aux[4];
2452 to[stride * 5] = aux[5];
2453 to[stride * 6] = aux[6];
2454 to[stride * 7] = aux[7];
2455 to[stride * 8] = aux[8];
2456 to[stride * 9] = aux[9];
2457 to[stride * 10] = aux[10];
2458 to[stride * 11] = aux[11];
2459 to[stride * 12] = aux[12];
2460 to[stride * 13] = aux[13];
2461 to[stride * 14] = aux[14];
2462 to[stride * 15] = aux[15];
2463}
2464
2465EIGEN_STRONG_INLINE void ptranspose(PacketBlock<Packet16h, 16>& kernel) {
2466 __m256i a = kernel.packet[0];
2467 __m256i b = kernel.packet[1];
2468 __m256i c = kernel.packet[2];
2469 __m256i d = kernel.packet[3];
2470 __m256i e = kernel.packet[4];
2471 __m256i f = kernel.packet[5];
2472 __m256i g = kernel.packet[6];
2473 __m256i h = kernel.packet[7];
2474 __m256i i = kernel.packet[8];
2475 __m256i j = kernel.packet[9];
2476 __m256i k = kernel.packet[10];
2477 __m256i l = kernel.packet[11];
2478 __m256i m = kernel.packet[12];
2479 __m256i n = kernel.packet[13];
2480 __m256i o = kernel.packet[14];
2481 __m256i p = kernel.packet[15];
2482
2483 __m256i ab_07 = _mm256_unpacklo_epi16(a, b);
2484 __m256i cd_07 = _mm256_unpacklo_epi16(c, d);
2485 __m256i ef_07 = _mm256_unpacklo_epi16(e, f);
2486 __m256i gh_07 = _mm256_unpacklo_epi16(g, h);
2487 __m256i ij_07 = _mm256_unpacklo_epi16(i, j);
2488 __m256i kl_07 = _mm256_unpacklo_epi16(k, l);
2489 __m256i mn_07 = _mm256_unpacklo_epi16(m, n);
2490 __m256i op_07 = _mm256_unpacklo_epi16(o, p);
2491
2492 __m256i ab_8f = _mm256_unpackhi_epi16(a, b);
2493 __m256i cd_8f = _mm256_unpackhi_epi16(c, d);
2494 __m256i ef_8f = _mm256_unpackhi_epi16(e, f);
2495 __m256i gh_8f = _mm256_unpackhi_epi16(g, h);
2496 __m256i ij_8f = _mm256_unpackhi_epi16(i, j);
2497 __m256i kl_8f = _mm256_unpackhi_epi16(k, l);
2498 __m256i mn_8f = _mm256_unpackhi_epi16(m, n);
2499 __m256i op_8f = _mm256_unpackhi_epi16(o, p);
2500
2501 __m256i abcd_03 = _mm256_unpacklo_epi32(ab_07, cd_07);
2502 __m256i abcd_47 = _mm256_unpackhi_epi32(ab_07, cd_07);
2503 __m256i efgh_03 = _mm256_unpacklo_epi32(ef_07, gh_07);
2504 __m256i efgh_47 = _mm256_unpackhi_epi32(ef_07, gh_07);
2505 __m256i ijkl_03 = _mm256_unpacklo_epi32(ij_07, kl_07);
2506 __m256i ijkl_47 = _mm256_unpackhi_epi32(ij_07, kl_07);
2507 __m256i mnop_03 = _mm256_unpacklo_epi32(mn_07, op_07);
2508 __m256i mnop_47 = _mm256_unpackhi_epi32(mn_07, op_07);
2509
2510 __m256i abcd_8b = _mm256_unpacklo_epi32(ab_8f, cd_8f);
2511 __m256i abcd_cf = _mm256_unpackhi_epi32(ab_8f, cd_8f);
2512 __m256i efgh_8b = _mm256_unpacklo_epi32(ef_8f, gh_8f);
2513 __m256i efgh_cf = _mm256_unpackhi_epi32(ef_8f, gh_8f);
2514 __m256i ijkl_8b = _mm256_unpacklo_epi32(ij_8f, kl_8f);
2515 __m256i ijkl_cf = _mm256_unpackhi_epi32(ij_8f, kl_8f);
2516 __m256i mnop_8b = _mm256_unpacklo_epi32(mn_8f, op_8f);
2517 __m256i mnop_cf = _mm256_unpackhi_epi32(mn_8f, op_8f);
2518
2519 __m256i abcdefgh_01 = _mm256_unpacklo_epi64(abcd_03, efgh_03);
2520 __m256i abcdefgh_23 = _mm256_unpackhi_epi64(abcd_03, efgh_03);
2521 __m256i ijklmnop_01 = _mm256_unpacklo_epi64(ijkl_03, mnop_03);
2522 __m256i ijklmnop_23 = _mm256_unpackhi_epi64(ijkl_03, mnop_03);
2523 __m256i abcdefgh_45 = _mm256_unpacklo_epi64(abcd_47, efgh_47);
2524 __m256i abcdefgh_67 = _mm256_unpackhi_epi64(abcd_47, efgh_47);
2525 __m256i ijklmnop_45 = _mm256_unpacklo_epi64(ijkl_47, mnop_47);
2526 __m256i ijklmnop_67 = _mm256_unpackhi_epi64(ijkl_47, mnop_47);
2527 __m256i abcdefgh_89 = _mm256_unpacklo_epi64(abcd_8b, efgh_8b);
2528 __m256i abcdefgh_ab = _mm256_unpackhi_epi64(abcd_8b, efgh_8b);
2529 __m256i ijklmnop_89 = _mm256_unpacklo_epi64(ijkl_8b, mnop_8b);
2530 __m256i ijklmnop_ab = _mm256_unpackhi_epi64(ijkl_8b, mnop_8b);
2531 __m256i abcdefgh_cd = _mm256_unpacklo_epi64(abcd_cf, efgh_cf);
2532 __m256i abcdefgh_ef = _mm256_unpackhi_epi64(abcd_cf, efgh_cf);
2533 __m256i ijklmnop_cd = _mm256_unpacklo_epi64(ijkl_cf, mnop_cf);
2534 __m256i ijklmnop_ef = _mm256_unpackhi_epi64(ijkl_cf, mnop_cf);
2535
2536 // NOTE: no unpacklo/hi instr in this case, so using permute instr.
2537 __m256i a_p_0 = _mm256_permute2x128_si256(abcdefgh_01, ijklmnop_01, 0x20);
2538 __m256i a_p_1 = _mm256_permute2x128_si256(abcdefgh_23, ijklmnop_23, 0x20);
2539 __m256i a_p_2 = _mm256_permute2x128_si256(abcdefgh_45, ijklmnop_45, 0x20);
2540 __m256i a_p_3 = _mm256_permute2x128_si256(abcdefgh_67, ijklmnop_67, 0x20);
2541 __m256i a_p_4 = _mm256_permute2x128_si256(abcdefgh_89, ijklmnop_89, 0x20);
2542 __m256i a_p_5 = _mm256_permute2x128_si256(abcdefgh_ab, ijklmnop_ab, 0x20);
2543 __m256i a_p_6 = _mm256_permute2x128_si256(abcdefgh_cd, ijklmnop_cd, 0x20);
2544 __m256i a_p_7 = _mm256_permute2x128_si256(abcdefgh_ef, ijklmnop_ef, 0x20);
2545 __m256i a_p_8 = _mm256_permute2x128_si256(abcdefgh_01, ijklmnop_01, 0x31);
2546 __m256i a_p_9 = _mm256_permute2x128_si256(abcdefgh_23, ijklmnop_23, 0x31);
2547 __m256i a_p_a = _mm256_permute2x128_si256(abcdefgh_45, ijklmnop_45, 0x31);
2548 __m256i a_p_b = _mm256_permute2x128_si256(abcdefgh_67, ijklmnop_67, 0x31);
2549 __m256i a_p_c = _mm256_permute2x128_si256(abcdefgh_89, ijklmnop_89, 0x31);
2550 __m256i a_p_d = _mm256_permute2x128_si256(abcdefgh_ab, ijklmnop_ab, 0x31);
2551 __m256i a_p_e = _mm256_permute2x128_si256(abcdefgh_cd, ijklmnop_cd, 0x31);
2552 __m256i a_p_f = _mm256_permute2x128_si256(abcdefgh_ef, ijklmnop_ef, 0x31);
2553
2554 kernel.packet[0] = a_p_0;
2555 kernel.packet[1] = a_p_1;
2556 kernel.packet[2] = a_p_2;
2557 kernel.packet[3] = a_p_3;
2558 kernel.packet[4] = a_p_4;
2559 kernel.packet[5] = a_p_5;
2560 kernel.packet[6] = a_p_6;
2561 kernel.packet[7] = a_p_7;
2562 kernel.packet[8] = a_p_8;
2563 kernel.packet[9] = a_p_9;
2564 kernel.packet[10] = a_p_a;
2565 kernel.packet[11] = a_p_b;
2566 kernel.packet[12] = a_p_c;
2567 kernel.packet[13] = a_p_d;
2568 kernel.packet[14] = a_p_e;
2569 kernel.packet[15] = a_p_f;
2570}
2571
2572EIGEN_STRONG_INLINE void ptranspose(PacketBlock<Packet16h, 8>& kernel) {
2573 EIGEN_ALIGN64 half in[8][16];
2574 pstore<half>(in[0], kernel.packet[0]);
2575 pstore<half>(in[1], kernel.packet[1]);
2576 pstore<half>(in[2], kernel.packet[2]);
2577 pstore<half>(in[3], kernel.packet[3]);
2578 pstore<half>(in[4], kernel.packet[4]);
2579 pstore<half>(in[5], kernel.packet[5]);
2580 pstore<half>(in[6], kernel.packet[6]);
2581 pstore<half>(in[7], kernel.packet[7]);
2582
2583 EIGEN_ALIGN64 half out[8][16];
2584
2585 for (int i = 0; i < 8; ++i) {
2586 for (int j = 0; j < 8; ++j) {
2587 out[i][j] = in[j][2 * i];
2588 }
2589 for (int j = 0; j < 8; ++j) {
2590 out[i][j + 8] = in[j][2 * i + 1];
2591 }
2592 }
2593
2594 kernel.packet[0] = pload<Packet16h>(out[0]);
2595 kernel.packet[1] = pload<Packet16h>(out[1]);
2596 kernel.packet[2] = pload<Packet16h>(out[2]);
2597 kernel.packet[3] = pload<Packet16h>(out[3]);
2598 kernel.packet[4] = pload<Packet16h>(out[4]);
2599 kernel.packet[5] = pload<Packet16h>(out[5]);
2600 kernel.packet[6] = pload<Packet16h>(out[6]);
2601 kernel.packet[7] = pload<Packet16h>(out[7]);
2602}
2603
2604EIGEN_STRONG_INLINE void ptranspose(PacketBlock<Packet16h, 4>& kernel) {
2605 EIGEN_ALIGN64 half in[4][16];
2606 pstore<half>(in[0], kernel.packet[0]);
2607 pstore<half>(in[1], kernel.packet[1]);
2608 pstore<half>(in[2], kernel.packet[2]);
2609 pstore<half>(in[3], kernel.packet[3]);
2610
2611 EIGEN_ALIGN64 half out[4][16];
2612
2613 for (int i = 0; i < 4; ++i) {
2614 for (int j = 0; j < 4; ++j) {
2615 out[i][j] = in[j][4 * i];
2616 }
2617 for (int j = 0; j < 4; ++j) {
2618 out[i][j + 4] = in[j][4 * i + 1];
2619 }
2620 for (int j = 0; j < 4; ++j) {
2621 out[i][j + 8] = in[j][4 * i + 2];
2622 }
2623 for (int j = 0; j < 4; ++j) {
2624 out[i][j + 12] = in[j][4 * i + 3];
2625 }
2626 }
2627
2628 kernel.packet[0] = pload<Packet16h>(out[0]);
2629 kernel.packet[1] = pload<Packet16h>(out[1]);
2630 kernel.packet[2] = pload<Packet16h>(out[2]);
2631 kernel.packet[3] = pload<Packet16h>(out[3]);
2632}
2633
2634template <>
2635struct is_arithmetic<Packet16bf> {
2636 enum { value = true };
2637};
2638
2639template <>
2640struct packet_traits<bfloat16> : default_packet_traits {
2641 typedef Packet16bf type;
2642 typedef Packet8bf half;
2643 enum {
2644 Vectorizable = 1,
2645 AlignedOnScalar = 1,
2646 size = 16,
2647 HasBlend = 0,
2648 HasInsert = 1,
2649 HasSin = EIGEN_FAST_MATH,
2650 HasCos = EIGEN_FAST_MATH,
2651 HasSqrt = 1,
2652 HasRsqrt = 1,
2653#ifdef EIGEN_VECTORIZE_AVX512DQ
2654 HasLog = 1, // Currently fails test with bad accuracy.
2655 HasLog1p = 1,
2656 HasExpm1 = 1,
2657 HasNdtri = 1,
2658 HasBessel = 1,
2659#endif
2660 HasExp = 1,
2661 HasTanh = EIGEN_FAST_MATH,
2662 HasErf = EIGEN_FAST_MATH,
2663 HasCmp = 1,
2664 HasDiv = 1
2665 };
2666};
2667
2668template <>
2669struct unpacket_traits<Packet16bf> {
2670 typedef bfloat16 type;
2671 enum {
2672 size = 16,
2673 alignment = Aligned32,
2674 vectorizable = true,
2675 masked_load_available = false,
2676 masked_store_available = false
2677 };
2678 typedef Packet8bf half;
2679};
2680
2681template <>
2682EIGEN_STRONG_INLINE Packet16bf pset1<Packet16bf>(const bfloat16& from) {
2683 return _mm256_set1_epi16(from.value);
2684}
2685
2686template <>
2687EIGEN_STRONG_INLINE bfloat16 pfirst<Packet16bf>(const Packet16bf& from) {
2688 bfloat16 t;
2689 t.value = static_cast<unsigned short>(_mm256_extract_epi16(from, 0));
2690 return t;
2691}
2692
2693template <>
2694EIGEN_STRONG_INLINE Packet16bf pload<Packet16bf>(const bfloat16* from) {
2695 return _mm256_load_si256(reinterpret_cast<const __m256i*>(from));
2696}
2697
2698template <>
2699EIGEN_STRONG_INLINE Packet16bf ploadu<Packet16bf>(const bfloat16* from) {
2700 return _mm256_loadu_si256(reinterpret_cast<const __m256i*>(from));
2701}
2702
2703template <>
2704EIGEN_STRONG_INLINE void pstore<bfloat16>(bfloat16* to, const Packet16bf& from) {
2705 _mm256_store_si256(reinterpret_cast<__m256i*>(to), from);
2706}
2707
2708template <>
2709EIGEN_STRONG_INLINE void pstoreu<bfloat16>(bfloat16* to, const Packet16bf& from) {
2710 _mm256_storeu_si256(reinterpret_cast<__m256i*>(to), from);
2711}
2712
2713template <>
2714EIGEN_STRONG_INLINE Packet16bf ploaddup<Packet16bf>(const bfloat16* from) {
2715 unsigned short a = from[0].value;
2716 unsigned short b = from[1].value;
2717 unsigned short c = from[2].value;
2718 unsigned short d = from[3].value;
2719 unsigned short e = from[4].value;
2720 unsigned short f = from[5].value;
2721 unsigned short g = from[6].value;
2722 unsigned short h = from[7].value;
2723 return _mm256_set_epi16(h, h, g, g, f, f, e, e, d, d, c, c, b, b, a, a);
2724}
2725
2726template <>
2727EIGEN_STRONG_INLINE Packet16bf ploadquad(const bfloat16* from) {
2728 unsigned short a = from[0].value;
2729 unsigned short b = from[1].value;
2730 unsigned short c = from[2].value;
2731 unsigned short d = from[3].value;
2732 return _mm256_set_epi16(d, d, d, d, c, c, c, c, b, b, b, b, a, a, a, a);
2733}
2734
2735EIGEN_STRONG_INLINE Packet16f Bf16ToF32(const Packet16bf& a) {
2736 return _mm512_castsi512_ps(_mm512_slli_epi32(_mm512_cvtepu16_epi32(a), 16));
2737}
2738
2739// Convert float to bfloat16 according to round-to-nearest-even/denormals algorithm.
2740EIGEN_STRONG_INLINE Packet16bf F32ToBf16(const Packet16f& a) {
2741 Packet16bf r;
2742
2743#if defined(EIGEN_VECTORIZE_AVX512BF16) && EIGEN_GNUC_STRICT_AT_LEAST(10, 1, 0)
2744 // Since GCC 10.1 supports avx512bf16 and C style explicit cast
2745 // (C++ static_cast is not supported yet), do conversion via intrinsic
2746 // and register path for performance.
2747 r = (__m256i)(_mm512_cvtneps_pbh(a));
2748
2749#else
2750 __m512i t;
2751 __m512i input = _mm512_castps_si512(a);
2752 __m512i nan = _mm512_set1_epi32(0x7fc0);
2753
2754 // uint32_t lsb = (input >> 16) & 1;
2755 t = _mm512_and_si512(_mm512_srli_epi32(input, 16), _mm512_set1_epi32(1));
2756 // uint32_t rounding_bias = 0x7fff + lsb;
2757 t = _mm512_add_epi32(t, _mm512_set1_epi32(0x7fff));
2758 // input += rounding_bias;
2759 t = _mm512_add_epi32(t, input);
2760 // input = input >> 16;
2761 t = _mm512_srli_epi32(t, 16);
2762
2763 // Check NaN before converting back to bf16
2764 __mmask16 mask = _mm512_cmp_ps_mask(a, a, _CMP_ORD_Q);
2765
2766 t = _mm512_mask_blend_epi32(mask, nan, t);
2767 // output.value = static_cast<uint16_t>(input);
2768 r = _mm512_cvtepi32_epi16(t);
2769#endif // EIGEN_VECTORIZE_AVX512BF16
2770
2771 return r;
2772}
2773
2774template <>
2775EIGEN_STRONG_INLINE Packet16bf ptrue(const Packet16bf& a) {
2776 return Packet16bf(ptrue<Packet8i>(Packet8i(a)));
2777}
2778
2779template <>
2780EIGEN_STRONG_INLINE Packet16bf por(const Packet16bf& a, const Packet16bf& b) {
2781 return Packet16bf(por<Packet8i>(Packet8i(a), Packet8i(b)));
2782}
2783
2784template <>
2785EIGEN_STRONG_INLINE Packet16bf pxor(const Packet16bf& a, const Packet16bf& b) {
2786 return Packet16bf(pxor<Packet8i>(Packet8i(a), Packet8i(b)));
2787}
2788
2789template <>
2790EIGEN_STRONG_INLINE Packet16bf pand(const Packet16bf& a, const Packet16bf& b) {
2791 return Packet16bf(pand<Packet8i>(Packet8i(a), Packet8i(b)));
2792}
2793
2794template <>
2795EIGEN_STRONG_INLINE Packet16bf pandnot(const Packet16bf& a, const Packet16bf& b) {
2796 return Packet16bf(pandnot<Packet8i>(Packet8i(a), Packet8i(b)));
2797}
2798
2799template <>
2800EIGEN_STRONG_INLINE Packet16bf pselect(const Packet16bf& mask, const Packet16bf& a, const Packet16bf& b) {
2801 // Input mask is expected to be all 0/1, handle it with 8-bit
2802 // intrinsic for performance.
2803 return _mm256_blendv_epi8(b, a, mask);
2804}
2805
2806template <>
2807EIGEN_STRONG_INLINE Packet16bf pround<Packet16bf>(const Packet16bf& a) {
2808 return F32ToBf16(pround<Packet16f>(Bf16ToF32(a)));
2809}
2810
2811template <>
2812EIGEN_STRONG_INLINE Packet16bf print<Packet16bf>(const Packet16bf& a) {
2813 return F32ToBf16(print<Packet16f>(Bf16ToF32(a)));
2814}
2815
2816template <>
2817EIGEN_STRONG_INLINE Packet16bf pceil<Packet16bf>(const Packet16bf& a) {
2818 return F32ToBf16(pceil<Packet16f>(Bf16ToF32(a)));
2819}
2820
2821template <>
2822EIGEN_STRONG_INLINE Packet16bf pfloor<Packet16bf>(const Packet16bf& a) {
2823 return F32ToBf16(pfloor<Packet16f>(Bf16ToF32(a)));
2824}
2825
2826template <>
2827EIGEN_STRONG_INLINE Packet16bf ptrunc<Packet16bf>(const Packet16bf& a) {
2828 return F32ToBf16(ptrunc<Packet16f>(Bf16ToF32(a)));
2829}
2830
2831template <>
2832EIGEN_STRONG_INLINE Packet16bf pcmp_eq(const Packet16bf& a, const Packet16bf& b) {
2833 return Pack32To16(pcmp_eq(Bf16ToF32(a), Bf16ToF32(b)));
2834}
2835
2836template <>
2837EIGEN_STRONG_INLINE Packet16bf pcmp_le(const Packet16bf& a, const Packet16bf& b) {
2838 return Pack32To16(pcmp_le(Bf16ToF32(a), Bf16ToF32(b)));
2839}
2840
2841template <>
2842EIGEN_STRONG_INLINE Packet16bf pcmp_lt(const Packet16bf& a, const Packet16bf& b) {
2843 return Pack32To16(pcmp_lt(Bf16ToF32(a), Bf16ToF32(b)));
2844}
2845
2846template <>
2847EIGEN_STRONG_INLINE Packet16bf pcmp_lt_or_nan(const Packet16bf& a, const Packet16bf& b) {
2848 return Pack32To16(pcmp_lt_or_nan(Bf16ToF32(a), Bf16ToF32(b)));
2849}
2850
2851template <>
2852EIGEN_STRONG_INLINE Packet16bf pnegate(const Packet16bf& a) {
2853 Packet16bf sign_mask = _mm256_set1_epi16(static_cast<unsigned short>(0x8000));
2854 return _mm256_xor_si256(a, sign_mask);
2855}
2856
2857template <>
2858EIGEN_STRONG_INLINE Packet16bf pconj(const Packet16bf& a) {
2859 return a;
2860}
2861
2862template <>
2863EIGEN_STRONG_INLINE Packet16bf pabs(const Packet16bf& a) {
2864 const __m256i sign_mask = _mm256_set1_epi16(static_cast<numext::uint16_t>(0x8000));
2865 return _mm256_andnot_si256(sign_mask, a);
2866}
2867
2868template <>
2869EIGEN_STRONG_INLINE Packet16bf padd<Packet16bf>(const Packet16bf& a, const Packet16bf& b) {
2870 return F32ToBf16(padd<Packet16f>(Bf16ToF32(a), Bf16ToF32(b)));
2871}
2872
2873template <>
2874EIGEN_STRONG_INLINE Packet16bf psub<Packet16bf>(const Packet16bf& a, const Packet16bf& b) {
2875 return F32ToBf16(psub<Packet16f>(Bf16ToF32(a), Bf16ToF32(b)));
2876}
2877
2878template <>
2879EIGEN_STRONG_INLINE Packet16bf pmul<Packet16bf>(const Packet16bf& a, const Packet16bf& b) {
2880 return F32ToBf16(pmul<Packet16f>(Bf16ToF32(a), Bf16ToF32(b)));
2881}
2882
2883template <>
2884EIGEN_STRONG_INLINE Packet16bf pdiv<Packet16bf>(const Packet16bf& a, const Packet16bf& b) {
2885 return F32ToBf16(pdiv<Packet16f>(Bf16ToF32(a), Bf16ToF32(b)));
2886}
2887
2888template <>
2889EIGEN_STRONG_INLINE Packet16bf pmin<Packet16bf>(const Packet16bf& a, const Packet16bf& b) {
2890 return F32ToBf16(pmin<Packet16f>(Bf16ToF32(a), Bf16ToF32(b)));
2891}
2892
2893template <>
2894EIGEN_STRONG_INLINE Packet16bf pmax<Packet16bf>(const Packet16bf& a, const Packet16bf& b) {
2895 return F32ToBf16(pmax<Packet16f>(Bf16ToF32(a), Bf16ToF32(b)));
2896}
2897
2898template <>
2899EIGEN_STRONG_INLINE Packet16bf plset<Packet16bf>(const bfloat16& a) {
2900 return F32ToBf16(plset<Packet16f>(static_cast<float>(a)));
2901}
2902
2903template <>
2904EIGEN_STRONG_INLINE Packet8bf predux_half_dowto4<Packet16bf>(const Packet16bf& a) {
2905 Packet8bf lane0 = _mm256_extractf128_si256(a, 0);
2906 Packet8bf lane1 = _mm256_extractf128_si256(a, 1);
2907 return padd<Packet8bf>(lane0, lane1);
2908}
2909
2910template <>
2911EIGEN_STRONG_INLINE bfloat16 predux<Packet16bf>(const Packet16bf& p) {
2912 return static_cast<bfloat16>(predux<Packet16f>(Bf16ToF32(p)));
2913}
2914
2915template <>
2916EIGEN_STRONG_INLINE bfloat16 predux_mul<Packet16bf>(const Packet16bf& from) {
2917 return static_cast<bfloat16>(predux_mul<Packet16f>(Bf16ToF32(from)));
2918}
2919
2920template <>
2921EIGEN_STRONG_INLINE bfloat16 predux_min<Packet16bf>(const Packet16bf& from) {
2922 return static_cast<bfloat16>(predux_min<Packet16f>(Bf16ToF32(from)));
2923}
2924
2925template <>
2926EIGEN_STRONG_INLINE bfloat16 predux_max<Packet16bf>(const Packet16bf& from) {
2927 return static_cast<bfloat16>(predux_max<Packet16f>(Bf16ToF32(from)));
2928}
2929
2930template <>
2931EIGEN_STRONG_INLINE Packet16bf preverse(const Packet16bf& a) {
2932 __m256i m = _mm256_setr_epi8(14, 15, 12, 13, 10, 11, 8, 9, 6, 7, 4, 5, 2, 3, 0, 1, 14, 15, 12, 13, 10, 11, 8, 9, 6, 7,
2933 4, 5, 2, 3, 0, 1);
2934
2935 Packet16bf res;
2936 // Swap hi and lo first because shuffle is in 128-bit lanes.
2937 res = _mm256_permute2x128_si256(a, a, 1);
2938 // Shuffle 8-bit values in src within 2*128-bit lanes.
2939 return _mm256_shuffle_epi8(res, m);
2940}
2941
2942template <>
2943EIGEN_STRONG_INLINE Packet16bf pgather<bfloat16, Packet16bf>(const bfloat16* from, Index stride) {
2944 return _mm256_set_epi16(
2945 from[15 * stride].value, from[14 * stride].value, from[13 * stride].value, from[12 * stride].value,
2946 from[11 * stride].value, from[10 * stride].value, from[9 * stride].value, from[8 * stride].value,
2947 from[7 * stride].value, from[6 * stride].value, from[5 * stride].value, from[4 * stride].value,
2948 from[3 * stride].value, from[2 * stride].value, from[1 * stride].value, from[0 * stride].value);
2949}
2950
2951template <>
2952EIGEN_STRONG_INLINE void pscatter<bfloat16, Packet16bf>(bfloat16* to, const Packet16bf& from, Index stride) {
2953 EIGEN_ALIGN64 bfloat16 aux[16];
2954 pstore(aux, from);
2955 to[stride * 0] = aux[0];
2956 to[stride * 1] = aux[1];
2957 to[stride * 2] = aux[2];
2958 to[stride * 3] = aux[3];
2959 to[stride * 4] = aux[4];
2960 to[stride * 5] = aux[5];
2961 to[stride * 6] = aux[6];
2962 to[stride * 7] = aux[7];
2963 to[stride * 8] = aux[8];
2964 to[stride * 9] = aux[9];
2965 to[stride * 10] = aux[10];
2966 to[stride * 11] = aux[11];
2967 to[stride * 12] = aux[12];
2968 to[stride * 13] = aux[13];
2969 to[stride * 14] = aux[14];
2970 to[stride * 15] = aux[15];
2971}
2972
2973EIGEN_STRONG_INLINE void ptranspose(PacketBlock<Packet16bf, 16>& kernel) {
2974 __m256i a = kernel.packet[0];
2975 __m256i b = kernel.packet[1];
2976 __m256i c = kernel.packet[2];
2977 __m256i d = kernel.packet[3];
2978 __m256i e = kernel.packet[4];
2979 __m256i f = kernel.packet[5];
2980 __m256i g = kernel.packet[6];
2981 __m256i h = kernel.packet[7];
2982 __m256i i = kernel.packet[8];
2983 __m256i j = kernel.packet[9];
2984 __m256i k = kernel.packet[10];
2985 __m256i l = kernel.packet[11];
2986 __m256i m = kernel.packet[12];
2987 __m256i n = kernel.packet[13];
2988 __m256i o = kernel.packet[14];
2989 __m256i p = kernel.packet[15];
2990
2991 __m256i ab_07 = _mm256_unpacklo_epi16(a, b);
2992 __m256i cd_07 = _mm256_unpacklo_epi16(c, d);
2993 __m256i ef_07 = _mm256_unpacklo_epi16(e, f);
2994 __m256i gh_07 = _mm256_unpacklo_epi16(g, h);
2995 __m256i ij_07 = _mm256_unpacklo_epi16(i, j);
2996 __m256i kl_07 = _mm256_unpacklo_epi16(k, l);
2997 __m256i mn_07 = _mm256_unpacklo_epi16(m, n);
2998 __m256i op_07 = _mm256_unpacklo_epi16(o, p);
2999
3000 __m256i ab_8f = _mm256_unpackhi_epi16(a, b);
3001 __m256i cd_8f = _mm256_unpackhi_epi16(c, d);
3002 __m256i ef_8f = _mm256_unpackhi_epi16(e, f);
3003 __m256i gh_8f = _mm256_unpackhi_epi16(g, h);
3004 __m256i ij_8f = _mm256_unpackhi_epi16(i, j);
3005 __m256i kl_8f = _mm256_unpackhi_epi16(k, l);
3006 __m256i mn_8f = _mm256_unpackhi_epi16(m, n);
3007 __m256i op_8f = _mm256_unpackhi_epi16(o, p);
3008
3009 __m256i abcd_03 = _mm256_unpacklo_epi32(ab_07, cd_07);
3010 __m256i abcd_47 = _mm256_unpackhi_epi32(ab_07, cd_07);
3011 __m256i efgh_03 = _mm256_unpacklo_epi32(ef_07, gh_07);
3012 __m256i efgh_47 = _mm256_unpackhi_epi32(ef_07, gh_07);
3013 __m256i ijkl_03 = _mm256_unpacklo_epi32(ij_07, kl_07);
3014 __m256i ijkl_47 = _mm256_unpackhi_epi32(ij_07, kl_07);
3015 __m256i mnop_03 = _mm256_unpacklo_epi32(mn_07, op_07);
3016 __m256i mnop_47 = _mm256_unpackhi_epi32(mn_07, op_07);
3017
3018 __m256i abcd_8b = _mm256_unpacklo_epi32(ab_8f, cd_8f);
3019 __m256i abcd_cf = _mm256_unpackhi_epi32(ab_8f, cd_8f);
3020 __m256i efgh_8b = _mm256_unpacklo_epi32(ef_8f, gh_8f);
3021 __m256i efgh_cf = _mm256_unpackhi_epi32(ef_8f, gh_8f);
3022 __m256i ijkl_8b = _mm256_unpacklo_epi32(ij_8f, kl_8f);
3023 __m256i ijkl_cf = _mm256_unpackhi_epi32(ij_8f, kl_8f);
3024 __m256i mnop_8b = _mm256_unpacklo_epi32(mn_8f, op_8f);
3025 __m256i mnop_cf = _mm256_unpackhi_epi32(mn_8f, op_8f);
3026
3027 __m256i abcdefgh_01 = _mm256_unpacklo_epi64(abcd_03, efgh_03);
3028 __m256i abcdefgh_23 = _mm256_unpackhi_epi64(abcd_03, efgh_03);
3029 __m256i ijklmnop_01 = _mm256_unpacklo_epi64(ijkl_03, mnop_03);
3030 __m256i ijklmnop_23 = _mm256_unpackhi_epi64(ijkl_03, mnop_03);
3031 __m256i abcdefgh_45 = _mm256_unpacklo_epi64(abcd_47, efgh_47);
3032 __m256i abcdefgh_67 = _mm256_unpackhi_epi64(abcd_47, efgh_47);
3033 __m256i ijklmnop_45 = _mm256_unpacklo_epi64(ijkl_47, mnop_47);
3034 __m256i ijklmnop_67 = _mm256_unpackhi_epi64(ijkl_47, mnop_47);
3035 __m256i abcdefgh_89 = _mm256_unpacklo_epi64(abcd_8b, efgh_8b);
3036 __m256i abcdefgh_ab = _mm256_unpackhi_epi64(abcd_8b, efgh_8b);
3037 __m256i ijklmnop_89 = _mm256_unpacklo_epi64(ijkl_8b, mnop_8b);
3038 __m256i ijklmnop_ab = _mm256_unpackhi_epi64(ijkl_8b, mnop_8b);
3039 __m256i abcdefgh_cd = _mm256_unpacklo_epi64(abcd_cf, efgh_cf);
3040 __m256i abcdefgh_ef = _mm256_unpackhi_epi64(abcd_cf, efgh_cf);
3041 __m256i ijklmnop_cd = _mm256_unpacklo_epi64(ijkl_cf, mnop_cf);
3042 __m256i ijklmnop_ef = _mm256_unpackhi_epi64(ijkl_cf, mnop_cf);
3043
3044 // NOTE: no unpacklo/hi instr in this case, so using permute instr.
3045 kernel.packet[0] = _mm256_permute2x128_si256(abcdefgh_01, ijklmnop_01, 0x20);
3046 kernel.packet[1] = _mm256_permute2x128_si256(abcdefgh_23, ijklmnop_23, 0x20);
3047 kernel.packet[2] = _mm256_permute2x128_si256(abcdefgh_45, ijklmnop_45, 0x20);
3048 kernel.packet[3] = _mm256_permute2x128_si256(abcdefgh_67, ijklmnop_67, 0x20);
3049 kernel.packet[4] = _mm256_permute2x128_si256(abcdefgh_89, ijklmnop_89, 0x20);
3050 kernel.packet[5] = _mm256_permute2x128_si256(abcdefgh_ab, ijklmnop_ab, 0x20);
3051 kernel.packet[6] = _mm256_permute2x128_si256(abcdefgh_cd, ijklmnop_cd, 0x20);
3052 kernel.packet[7] = _mm256_permute2x128_si256(abcdefgh_ef, ijklmnop_ef, 0x20);
3053 kernel.packet[8] = _mm256_permute2x128_si256(abcdefgh_01, ijklmnop_01, 0x31);
3054 kernel.packet[9] = _mm256_permute2x128_si256(abcdefgh_23, ijklmnop_23, 0x31);
3055 kernel.packet[10] = _mm256_permute2x128_si256(abcdefgh_45, ijklmnop_45, 0x31);
3056 kernel.packet[11] = _mm256_permute2x128_si256(abcdefgh_67, ijklmnop_67, 0x31);
3057 kernel.packet[12] = _mm256_permute2x128_si256(abcdefgh_89, ijklmnop_89, 0x31);
3058 kernel.packet[13] = _mm256_permute2x128_si256(abcdefgh_ab, ijklmnop_ab, 0x31);
3059 kernel.packet[14] = _mm256_permute2x128_si256(abcdefgh_cd, ijklmnop_cd, 0x31);
3060 kernel.packet[15] = _mm256_permute2x128_si256(abcdefgh_ef, ijklmnop_ef, 0x31);
3061}
3062
3063EIGEN_STRONG_INLINE void ptranspose(PacketBlock<Packet16bf, 4>& kernel) {
3064 __m256i a = kernel.packet[0];
3065 __m256i b = kernel.packet[1];
3066 __m256i c = kernel.packet[2];
3067 __m256i d = kernel.packet[3];
3068
3069 __m256i ab_07 = _mm256_unpacklo_epi16(a, b);
3070 __m256i cd_07 = _mm256_unpacklo_epi16(c, d);
3071 __m256i ab_8f = _mm256_unpackhi_epi16(a, b);
3072 __m256i cd_8f = _mm256_unpackhi_epi16(c, d);
3073
3074 __m256i abcd_03 = _mm256_unpacklo_epi32(ab_07, cd_07);
3075 __m256i abcd_47 = _mm256_unpackhi_epi32(ab_07, cd_07);
3076 __m256i abcd_8b = _mm256_unpacklo_epi32(ab_8f, cd_8f);
3077 __m256i abcd_cf = _mm256_unpackhi_epi32(ab_8f, cd_8f);
3078
3079 // NOTE: no unpacklo/hi instr in this case, so using permute instr.
3080 kernel.packet[0] = _mm256_permute2x128_si256(abcd_03, abcd_47, 0x20);
3081 kernel.packet[1] = _mm256_permute2x128_si256(abcd_8b, abcd_cf, 0x20);
3082 kernel.packet[2] = _mm256_permute2x128_si256(abcd_03, abcd_47, 0x31);
3083 kernel.packet[3] = _mm256_permute2x128_si256(abcd_8b, abcd_cf, 0x31);
3084}
3085
3086} // end namespace internal
3087
3088} // end namespace Eigen
3089
3090#endif // EIGEN_PACKET_MATH_AVX512_H
@ Aligned64
Definition Constants.h:239
@ Aligned32
Definition Constants.h:238
Namespace containing all symbols from the Eigen library.
Definition Core:137