23 namespace seqan3::detail
28 template <simd::simd_concept simd_t>
29 constexpr simd_t load_avx512(
void const * mem_addr);
34 template <simd::simd_concept simd_t>
35 inline void transpose_matrix_avx512(
std::array<simd_t, simd_traits<simd_t>::length> & matrix);
40 template <simd::simd_concept target_simd_t, simd::simd_concept source_simd_t>
41 constexpr target_simd_t upcast_signed_avx512(source_simd_t
const & src);
46 template <simd::simd_concept target_simd_t, simd::simd_concept source_simd_t>
47 constexpr target_simd_t upcast_unsigned_avx512(source_simd_t
const & src);
52 template <u
int8_t index, simd::simd_concept simd_t>
53 constexpr simd_t extract_halve_avx512(simd_t
const & src);
58 template <u
int8_t index, simd::simd_concept simd_t>
59 constexpr simd_t extract_quarter_avx512(simd_t
const & src);
64 template <u
int8_t index, simd::simd_concept simd_t>
65 constexpr simd_t extract_eighth_avx512(simd_t
const & src);
75 namespace seqan3::detail
78 template <simd::simd_concept simd_t>
79 constexpr simd_t load_avx512(
void const * mem_addr)
81 return reinterpret_cast<simd_t>(_mm512_loadu_si512(mem_addr));
85 template <simd::simd_concept simd_t>
86 inline void transpose_matrix_avx512(
std::array<simd_t, simd_traits<simd_t>::length> & matrix);
88 template <simd::simd_concept target_simd_t, simd::simd_concept source_simd_t>
89 constexpr target_simd_t upcast_signed_avx512(source_simd_t
const & src)
91 __m512i
const & tmp = reinterpret_cast<__m512i const &>(src);
92 if constexpr (simd_traits<source_simd_t>::length == 64)
94 if constexpr (simd_traits<target_simd_t>::length == 32)
95 return reinterpret_cast<target_simd_t>(_mm512_cvtepi8_epi16(_mm512_castsi512_si256(tmp)));
96 if constexpr (simd_traits<target_simd_t>::length == 16)
97 return reinterpret_cast<target_simd_t>(_mm512_cvtepi8_epi32(_mm512_castsi512_si128(tmp)));
98 if constexpr (simd_traits<target_simd_t>::length == 8)
99 return reinterpret_cast<target_simd_t>(_mm512_cvtepi8_epi64(_mm512_castsi512_si128(tmp)));
101 else if constexpr (simd_traits<source_simd_t>::length == 32)
103 if constexpr (simd_traits<target_simd_t>::length == 16)
104 return reinterpret_cast<target_simd_t>(_mm512_cvtepi16_epi32(_mm512_castsi512_si256(tmp)));
105 if constexpr (simd_traits<target_simd_t>::length == 8)
106 return reinterpret_cast<target_simd_t>(_mm512_cvtepi16_epi64(_mm512_castsi512_si128(tmp)));
110 static_assert(simd_traits<source_simd_t>::length == 16,
"Expected 32 bit scalar type.");
111 return reinterpret_cast<target_simd_t>(_mm512_cvtepi32_epi64(_mm512_castsi512_si256(tmp)));
115 template <simd::simd_concept target_simd_t, simd::simd_concept source_simd_t>
116 constexpr target_simd_t upcast_unsigned_avx512(source_simd_t
const & src)
118 __m512i
const & tmp = reinterpret_cast<__m512i const &>(src);
119 if constexpr (simd_traits<source_simd_t>::length == 64)
121 if constexpr (simd_traits<target_simd_t>::length == 32)
122 return reinterpret_cast<target_simd_t>(_mm512_cvtepu8_epi16(_mm512_castsi512_si256(tmp)));
123 if constexpr (simd_traits<target_simd_t>::length == 16)
124 return reinterpret_cast<target_simd_t>(_mm512_cvtepu8_epi32(_mm512_castsi512_si128(tmp)));
125 if constexpr (simd_traits<target_simd_t>::length == 8)
126 return reinterpret_cast<target_simd_t>(_mm512_cvtepu8_epi64(_mm512_castsi512_si128(tmp)));
128 else if constexpr (simd_traits<source_simd_t>::length == 32)
130 if constexpr (simd_traits<target_simd_t>::length == 16)
131 return reinterpret_cast<target_simd_t>(_mm512_cvtepu16_epi32(_mm512_castsi512_si256(tmp)));
132 if constexpr (simd_traits<target_simd_t>::length == 8)
133 return reinterpret_cast<target_simd_t>(_mm512_cvtepu16_epi64(_mm512_castsi512_si128(tmp)));
137 static_assert(simd_traits<source_simd_t>::length == 16,
"Expected 32 bit scalar type.");
138 return reinterpret_cast<target_simd_t>(_mm512_cvtepu32_epi64(_mm512_castsi512_si256(tmp)));
143 template <u
int8_t index, simd::simd_concept simd_t>
144 constexpr simd_t extract_halve_avx512(simd_t
const & src);
147 template <u
int8_t index, simd::simd_concept simd_t>
148 constexpr simd_t extract_quarter_avx512(simd_t
const & src);
151 template <u
int8_t index, simd::simd_concept simd_t>
152 constexpr simd_t extract_eighth_avx512(simd_t
const & src);
156 #endif // __AVX512F__