Vector Optimized Library of Kernels 3.2.0
Architecture-tuned implementations of math kernels
Loading...
Searching...
No Matches
volk_32f_invsqrt_32f.h
Go to the documentation of this file.
1/* -*- c++ -*- */
2/*
3 * Copyright 2013, 2014 Free Software Foundation, Inc.
4 *
5 * This file is part of VOLK
6 *
7 * SPDX-License-Identifier: LGPL-3.0-or-later
8 */
9
51
52#ifndef INCLUDED_volk_32f_invsqrt_32f_a_H
53#define INCLUDED_volk_32f_invsqrt_32f_a_H
54
55#include <inttypes.h>
56#include <math.h>
57#include <stdio.h>
58#include <string.h>
59
60static inline float Q_rsqrt(float number)
61{
62 float x2;
63 const float threehalfs = 1.5F;
64 union f32_to_i32 {
65 int32_t i;
66 float f;
67 } u;
68
69 x2 = number * 0.5F;
70 u.f = number;
71 u.i = 0x5f3759df - (u.i >> 1); // what the fuck?
72 u.f = u.f * (threehalfs - (x2 * u.f * u.f)); // 1st iteration
73 // u.f = u.f * ( threehalfs - ( x2 * u.f * u.f ) ); // 2nd iteration, this can be
74 // removed
75
76 return u.f;
77}
78
79#ifdef LV_HAVE_AVX
80#include <immintrin.h>
81
82static inline void
83volk_32f_invsqrt_32f_a_avx(float* cVector, const float* aVector, unsigned int num_points)
84{
85 unsigned int number = 0;
86 const unsigned int eighthPoints = num_points / 8;
87
88 float* cPtr = cVector;
89 const float* aPtr = aVector;
90 __m256 aVal, cVal;
91 for (; number < eighthPoints; number++) {
92 aVal = _mm256_load_ps(aPtr);
93 cVal = _mm256_rsqrt_ps(aVal);
94 _mm256_store_ps(cPtr, cVal);
95 aPtr += 8;
96 cPtr += 8;
97 }
98
99 number = eighthPoints * 8;
100 for (; number < num_points; number++) {
101 *cPtr++ = Q_rsqrt(*aPtr++);
102 }
103}
104#endif /* LV_HAVE_AVX */
105
106
107#ifdef LV_HAVE_SSE
108#include <xmmintrin.h>
109
110static inline void
111volk_32f_invsqrt_32f_a_sse(float* cVector, const float* aVector, unsigned int num_points)
112{
113 unsigned int number = 0;
114 const unsigned int quarterPoints = num_points / 4;
115
116 float* cPtr = cVector;
117 const float* aPtr = aVector;
118
119 __m128 aVal, cVal;
120 for (; number < quarterPoints; number++) {
121
122 aVal = _mm_load_ps(aPtr);
123
124 cVal = _mm_rsqrt_ps(aVal);
125
126 _mm_store_ps(cPtr, cVal); // Store the results back into the C container
127
128 aPtr += 4;
129 cPtr += 4;
130 }
131
132 number = quarterPoints * 4;
133 for (; number < num_points; number++) {
134 *cPtr++ = Q_rsqrt(*aPtr++);
135 }
136}
137#endif /* LV_HAVE_SSE */
138
139
140#ifdef LV_HAVE_NEON
141#include <arm_neon.h>
142
143static inline void
144volk_32f_invsqrt_32f_neon(float* cVector, const float* aVector, unsigned int num_points)
145{
146 unsigned int number;
147 const unsigned int quarter_points = num_points / 4;
148
149 float* cPtr = cVector;
150 const float* aPtr = aVector;
151 float32x4_t a_val, c_val;
152 for (number = 0; number < quarter_points; ++number) {
153 a_val = vld1q_f32(aPtr);
154 c_val = vrsqrteq_f32(a_val);
155 vst1q_f32(cPtr, c_val);
156 aPtr += 4;
157 cPtr += 4;
158 }
159
160 for (number = quarter_points * 4; number < num_points; number++) {
161 *cPtr++ = Q_rsqrt(*aPtr++);
162 }
163}
164#endif /* LV_HAVE_NEON */
165
166
167#ifdef LV_HAVE_GENERIC
168
169static inline void volk_32f_invsqrt_32f_generic(float* cVector,
170 const float* aVector,
171 unsigned int num_points)
172{
173 float* cPtr = cVector;
174 const float* aPtr = aVector;
175 unsigned int number = 0;
176 for (number = 0; number < num_points; number++) {
177 *cPtr++ = Q_rsqrt(*aPtr++);
178 }
179}
180#endif /* LV_HAVE_GENERIC */
181
182#ifdef LV_HAVE_AVX
183#include <immintrin.h>
184
185static inline void
186volk_32f_invsqrt_32f_u_avx(float* cVector, const float* aVector, unsigned int num_points)
187{
188 unsigned int number = 0;
189 const unsigned int eighthPoints = num_points / 8;
190
191 float* cPtr = cVector;
192 const float* aPtr = aVector;
193 __m256 aVal, cVal;
194 for (; number < eighthPoints; number++) {
195 aVal = _mm256_loadu_ps(aPtr);
196 cVal = _mm256_rsqrt_ps(aVal);
197 _mm256_storeu_ps(cPtr, cVal);
198 aPtr += 8;
199 cPtr += 8;
200 }
201
202 number = eighthPoints * 8;
203 for (; number < num_points; number++) {
204 *cPtr++ = Q_rsqrt(*aPtr++);
205 }
206}
207#endif /* LV_HAVE_AVX */
208
209#ifdef LV_HAVE_RVV
210#include <riscv_vector.h>
211
212static inline void
213volk_32f_invsqrt_32f_rvv(float* cVector, const float* aVector, unsigned int num_points)
214{
215 size_t n = num_points;
216 for (size_t vl; n > 0; n -= vl, aVector += vl, cVector += vl) {
217 vl = __riscv_vsetvl_e32m8(n);
218 vfloat32m8_t v = __riscv_vle32_v_f32m8(aVector, vl);
219 __riscv_vse32(cVector, __riscv_vfrsqrt7(v, vl), vl);
220 }
221}
222#endif /*LV_HAVE_RVV*/
223
224#endif /* INCLUDED_volk_32f_invsqrt_32f_a_H */
static void volk_32f_invsqrt_32f_neon(float *cVector, const float *aVector, unsigned int num_points)
Definition volk_32f_invsqrt_32f.h:144
static void volk_32f_invsqrt_32f_a_avx(float *cVector, const float *aVector, unsigned int num_points)
Definition volk_32f_invsqrt_32f.h:83
static void volk_32f_invsqrt_32f_generic(float *cVector, const float *aVector, unsigned int num_points)
Definition volk_32f_invsqrt_32f.h:169
static void volk_32f_invsqrt_32f_a_sse(float *cVector, const float *aVector, unsigned int num_points)
Definition volk_32f_invsqrt_32f.h:111
static void volk_32f_invsqrt_32f_u_avx(float *cVector, const float *aVector, unsigned int num_points)
Definition volk_32f_invsqrt_32f.h:186
static float Q_rsqrt(float number)
Definition volk_32f_invsqrt_32f.h:60
for i
Definition volk_config_fixed.tmpl.h:13