diff --git a/CMakeLists.txt b/CMakeLists.txt index 2bcb539ea..aaabf7aa8 100644 --- a/CMakeLists.txt +++ b/CMakeLists.txt @@ -345,7 +345,7 @@ endif(NOT Boost_FOUND) ################################################################################ # GNU Radio - http://gnuradio.org/redmine/projects/gnuradio/wiki ################################################################################ -set(GR_REQUIRED_COMPONENTS RUNTIME ANALOG BLOCKS FFT FILTER PMT UHD) +set(GR_REQUIRED_COMPONENTS RUNTIME ANALOG BLOCKS FFT FILTER PMT) find_package(Gnuradio) if(PC_GNURADIO_RUNTIME_VERSION) if(PC_GNURADIO_RUNTIME_VERSION VERSION_LESS 3.7.3) @@ -822,18 +822,17 @@ endif(NOT GNUTLS_OPENSSL_LIBRARY) ################################################################################ -# Universal Hardware Driver (UHD) +# USRP Hardware Driver (UHD) - OPTIONAL ################################################################################ find_package(UHD) if(NOT UHD_FOUND) set(ENABLE_UHD OFF) - message(STATUS "The Universal Hardware Driver (UHD) based signal source will not be built,") - message(STATUS "so all USRP-based front-ends will not be usable.") - message(STATUS "Please check http://code.ettus.com/redmine/ettus/projects/uhd/wiki") + message(STATUS " The USRP Hardware Driver (UHD) signal source will not be built,") + message(STATUS " so all USRP-based front-ends will not be usable.") + message(STATUS " Please check http://files.ettus.com/manual/") else(NOT UHD_FOUND) - if(NOT GNURADIO_UHD_FOUND) - message(FATAL_ERROR "*** gnuradio-uhd 3.7 or later is required to build gnss-sdr") - endif(NOT GNURADIO_UHD_FOUND) + set(GR_REQUIRED_COMPONENTS UHD) + find_package(Gnuradio) set(ENABLE_UHD ON) endif(NOT UHD_FOUND) diff --git a/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/kernels/volk_gnsssdr/volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn.h b/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/kernels/volk_gnsssdr/volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn.h new file mode 100644 index 000000000..a958743a7 --- /dev/null +++ b/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/kernels/volk_gnsssdr/volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn.h @@ -0,0 +1,398 @@ +/*! + * \file volk_gnsssdr_16ic_x2_dot_prod_16ic_xn.h + * \brief Volk protokernel: multiplies N 16 bits vectors by a common vector phase rotated and accumulates the results in N 16 bits short complex outputs. + * \authors + * + * Volk protokernel that multiplies N 16 bits vectors by a common vector, which is phase-rotated by phase offset and phase increment, and accumulates the results in N 16 bits short complex outputs. + * It is optimized to perform the N tap correlation process in GNSS receivers. + * + * ------------------------------------------------------------------------- + * + * Copyright (C) 2010-2015 (see AUTHORS file for a list of contributors) + * + * GNSS-SDR is a software defined Global Navigation + * Satellite Systems receiver + * + * This file is part of GNSS-SDR. + * + * GNSS-SDR is free software: you can redistribute it and/or modify + * it under the terms of the GNU General Public License as published by + * the Free Software Foundation, either version 3 of the License, or + * (at your option) any later version. + * + * GNSS-SDR is distributed in the hope that it will be useful, + * but WITHOUT ANY WARRANTY; without even the implied warranty of + * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the + * GNU General Public License for more details. + * + * You should have received a copy of the GNU General Public License + * along with GNSS-SDR. If not, see . + * + * ------------------------------------------------------------------------- + */ + +#ifndef INCLUDED_volk_gnsssdr_16ic_xn_rotator_dot_prod_16ic_xn_H +#define INCLUDED_volk_gnsssdr_16ic_xn_rotator_dot_prod_16ic_xn_H + + +#include +#include + +#ifdef LV_HAVE_GENERIC +/*! + \brief Multiplies the reference complex vector with multiple versions of another complex vector, accumulates the results and stores them in the output vector + \param[out] result Array of num_a_vectors components with the multiple versions of in_a multiplied and accumulated The vector where the accumulated result will be stored + \param[in] in_common Pointer to one of the vectors to be multiplied and accumulated (reference vector) + \param[in] in_a Pointer to an array of pointers to multiple versions of the other vector to be multiplied and accumulated + \param[in] num_a_vectors Number of vectors to be multiplied by the reference vector and accumulated + \param[in] num_points The Number of complex values to be multiplied together, accumulated and stored into result + */ +static inline void volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn_generic(lv_16sc_t* result, const lv_16sc_t* in_common, const lv_32fc_t phase_inc, lv_32fc_t* phase, const lv_16sc_t** in_a, int num_a_vectors, unsigned int num_points) +{ + lv_16sc_t tmp16; + lv_32fc_t tmp32; + for (int n_vec = 0; n_vec < num_a_vectors; n_vec++) + { + result[n_vec] = lv_cmake(0,0); + } + for (unsigned int n = 0; n < num_points; n++) + { + tmp16 = *in_common++; + tmp32 = lv_cmake((float)lv_creal(tmp16), (float)lv_cimag(tmp16)) * (*phase); + tmp16 = lv_cmake((int16_t)rintf(lv_creal(tmp32)), (int16_t)rintf(lv_cimag(tmp32))); + (*phase) *= phase_inc; + for (int n_vec = 0; n_vec < num_a_vectors; n_vec++) + { + lv_16sc_t tmp = tmp16 * in_a[n_vec][n]; + result[n_vec] = lv_cmake(sat_adds16i(lv_creal(result[n_vec]), lv_creal(tmp)), sat_adds16i(lv_cimag(result[n_vec]), lv_cimag(tmp))); + } + } +} + +#endif /*LV_HAVE_GENERIC*/ + + +#ifdef LV_HAVE_SSE3 +#include + +/*! + \brief Multiplies the reference complex vector with multiple versions of another complex vector, accumulates the results and stores them in the output vector + \param[out] result Array of num_a_vectors components with the multiple versions of in_a multiplied and accumulated The vector where the accumulated result will be stored + \param[in] in_common Pointer to one of the vectors to be multiplied and accumulated (reference vector) + \param[in] in_a Pointer to an array of pointers to multiple versions of the other vector to be multiplied and accumulated + \param[in] num_a_vectors Number of vectors to be multiplied by the reference vector and accumulated + \param[in] num_points The Number of complex values to be multiplied together, accumulated and stored into result + */ +static inline void volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn_a_sse3(lv_16sc_t* out, const lv_16sc_t* in_common, const lv_32fc_t phase_inc, lv_32fc_t* phase, const lv_16sc_t** in_a, int num_a_vectors, unsigned int num_points) +{ + lv_16sc_t dotProduct = lv_cmake(0,0); + + const unsigned int sse_iters = num_points / 4; + + const lv_16sc_t** _in_a = in_a; + const lv_16sc_t* _in_common = in_common; + lv_16sc_t* _out = out; + + __VOLK_ATTR_ALIGNED(16) lv_16sc_t dotProductVector[4]; + + //todo dyn mem reg + + __m128i* realcacc; + __m128i* imagcacc; + + realcacc = (__m128i*)calloc(num_a_vectors, sizeof(__m128i)); //calloc also sets memory to 0 + imagcacc = (__m128i*)calloc(num_a_vectors, sizeof(__m128i)); //calloc also sets memory to 0 + + __m128i a,b,c, c_sr, mask_imag, mask_real, real, imag, imag1,imag2, b_sl, a_sl, result; + + mask_imag = _mm_set_epi8(255, 255, 0, 0, 255, 255, 0, 0, 255, 255, 0, 0, 255, 255, 0, 0); + mask_real = _mm_set_epi8(0, 0, 255, 255, 0, 0, 255, 255, 0, 0, 255, 255, 0, 0, 255, 255); + + // phase rotation registers + __m128 pa, pb, two_phase_acc_reg, two_phase_inc_reg; + __m128i pc1, pc2; + __attribute__((aligned(16))) lv_32fc_t two_phase_inc[2]; + two_phase_inc[0] = phase_inc * phase_inc; + two_phase_inc[1] = phase_inc * phase_inc; + two_phase_inc_reg = _mm_load_ps((float*) two_phase_inc); + __attribute__((aligned(16))) lv_32fc_t two_phase_acc[2]; + two_phase_acc[0] = (*phase); + two_phase_acc[1] = (*phase) * phase_inc; + two_phase_acc_reg = _mm_load_ps((float*)two_phase_acc); + __m128 yl, yh, tmp1, tmp2, tmp3; + lv_16sc_t tmp16; + lv_32fc_t tmp32; + + for(unsigned int number = 0; number < sse_iters; number++) + { + // Phase rotation on operand in_common starts here: + + pa = _mm_set_ps((float)(lv_cimag(_in_common[1])), (float)(lv_creal(_in_common[1])), (float)(lv_cimag(_in_common[0])), (float)(lv_creal(_in_common[0]))); // //load (2 byte imag, 2 byte real) x 2 into 128 bits reg + //complex 32fc multiplication b=a*two_phase_acc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(pa, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + pa = _mm_shuffle_ps(pa, pa, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(pa, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + pb = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + pc1 = _mm_cvtps_epi32(pb); // convert from 32fc to 32ic + + //complex 32fc multiplication two_phase_acc_reg=two_phase_acc_reg*two_phase_inc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(two_phase_inc_reg, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + tmp3 = _mm_shuffle_ps(two_phase_inc_reg, two_phase_inc_reg, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(tmp3, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + two_phase_acc_reg = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + + //next two samples + _in_common += 2; + pa = _mm_set_ps((float)(lv_cimag(_in_common[1])), (float)(lv_creal(_in_common[1])), (float)(lv_cimag(_in_common[0])), (float)(lv_creal(_in_common[0]))); // //load (2 byte imag, 2 byte real) x 2 into 128 bits reg + //complex 32fc multiplication b=a*two_phase_acc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(pa, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + pa = _mm_shuffle_ps(pa, pa, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(pa, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + pb = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + pc2 = _mm_cvtps_epi32(pb); // convert from 32fc to 32ic + + //complex 32fc multiplication two_phase_acc_reg=two_phase_acc_reg*two_phase_inc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(two_phase_inc_reg, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + tmp3 = _mm_shuffle_ps(two_phase_inc_reg, two_phase_inc_reg, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(tmp3, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + two_phase_acc_reg = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + + // store four rotated in_common samples in the register b + b = _mm_packs_epi32(pc1, pc2);// convert from 32ic to 16ic + + //next two samples + _in_common += 2; + + for (int n_vec = 0; n_vec < num_a_vectors; n_vec++) + { + a = _mm_load_si128((__m128i*)&(_in_a[n_vec][number*4])); //load (2 byte imag, 2 byte real) x 4 into 128 bits reg + + c = _mm_mullo_epi16 (a, b); // a3.i*b3.i, a3.r*b3.r, .... + + c_sr = _mm_srli_si128 (c, 2); // Shift a right by imm8 bytes while shifting in zeros, and store the results in dst. + real = _mm_subs_epi16 (c, c_sr); + + b_sl = _mm_slli_si128(b, 2); // b3.r, b2.i .... + a_sl = _mm_slli_si128(a, 2); // a3.r, a2.i .... + + imag1 = _mm_mullo_epi16(a, b_sl); // a3.i*b3.r, .... + imag2 = _mm_mullo_epi16(b, a_sl); // b3.i*a3.r, .... + + imag = _mm_adds_epi16(imag1, imag2); + + realcacc[n_vec] = _mm_adds_epi16 (realcacc[n_vec], real); + imagcacc[n_vec] = _mm_adds_epi16 (imagcacc[n_vec], imag); + + } + } + + for (int n_vec=0;n_vec + +/*! + \brief Multiplies the reference complex vector with multiple versions of another complex vector, accumulates the results and stores them in the output vector + \param[out] result Array of num_a_vectors components with the multiple versions of in_a multiplied and accumulated The vector where the accumulated result will be stored + \param[in] in_common Pointer to one of the vectors to be multiplied and accumulated (reference vector) + \param[in] in_a Pointer to an array of pointers to multiple versions of the other vector to be multiplied and accumulated + \param[in] num_a_vectors Number of vectors to be multiplied by the reference vector and accumulated + \param[in] num_points The Number of complex values to be multiplied together, accumulated and stored into result + */ +static inline void volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn_u_sse3(lv_16sc_t* out, const lv_16sc_t* in_common, const lv_32fc_t phase_inc, lv_32fc_t* phase, const lv_16sc_t** in_a, int num_a_vectors, unsigned int num_points) +{ + lv_16sc_t dotProduct = lv_cmake(0,0); + + const unsigned int sse_iters = num_points / 4; + + const lv_16sc_t** _in_a = in_a; + const lv_16sc_t* _in_common = in_common; + lv_16sc_t* _out = out; + + __VOLK_ATTR_ALIGNED(16) lv_16sc_t dotProductVector[4]; + + //todo dyn mem reg + + __m128i* realcacc; + __m128i* imagcacc; + + realcacc = (__m128i*)calloc(num_a_vectors, sizeof(__m128i)); //calloc also sets memory to 0 + imagcacc = (__m128i*)calloc(num_a_vectors, sizeof(__m128i)); //calloc also sets memory to 0 + + __m128i a,b,c, c_sr, mask_imag, mask_real, real, imag, imag1,imag2, b_sl, a_sl, result; + + mask_imag = _mm_set_epi8(255, 255, 0, 0, 255, 255, 0, 0, 255, 255, 0, 0, 255, 255, 0, 0); + mask_real = _mm_set_epi8(0, 0, 255, 255, 0, 0, 255, 255, 0, 0, 255, 255, 0, 0, 255, 255); + + // phase rotation registers + __m128 pa, pb, two_phase_acc_reg, two_phase_inc_reg; + __m128i pc1, pc2; + __attribute__((aligned(16))) lv_32fc_t two_phase_inc[2]; + two_phase_inc[0] = phase_inc * phase_inc; + two_phase_inc[1] = phase_inc * phase_inc; + two_phase_inc_reg = _mm_load_ps((float*) two_phase_inc); + __attribute__((aligned(16))) lv_32fc_t two_phase_acc[2]; + two_phase_acc[0] = (*phase); + two_phase_acc[1] = (*phase) * phase_inc; + two_phase_acc_reg = _mm_load_ps((float*)two_phase_acc); + __m128 yl, yh, tmp1, tmp2, tmp3; + lv_16sc_t tmp16; + lv_32fc_t tmp32; + + for(unsigned int number = 0; number < sse_iters; number++) + { + // Phase rotation on operand in_common starts here: + + pa = _mm_set_ps((float)(lv_cimag(_in_common[1])), (float)(lv_creal(_in_common[1])), (float)(lv_cimag(_in_common[0])), (float)(lv_creal(_in_common[0]))); // //load (2 byte imag, 2 byte real) x 2 into 128 bits reg + //complex 32fc multiplication b=a*two_phase_acc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(pa, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + pa = _mm_shuffle_ps(pa, pa, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(pa, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + pb = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + pc1 = _mm_cvtps_epi32(pb); // convert from 32fc to 32ic + + //complex 32fc multiplication two_phase_acc_reg=two_phase_acc_reg*two_phase_inc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(two_phase_inc_reg, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + tmp3 = _mm_shuffle_ps(two_phase_inc_reg, two_phase_inc_reg, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(tmp3, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + two_phase_acc_reg = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + + //next two samples + _in_common += 2; + pa = _mm_set_ps((float)(lv_cimag(_in_common[1])), (float)(lv_creal(_in_common[1])), (float)(lv_cimag(_in_common[0])), (float)(lv_creal(_in_common[0]))); // //load (2 byte imag, 2 byte real) x 2 into 128 bits reg + //complex 32fc multiplication b=a*two_phase_acc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(pa, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + pa = _mm_shuffle_ps(pa, pa, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(pa, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + pb = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + pc2 = _mm_cvtps_epi32(pb); // convert from 32fc to 32ic + + //complex 32fc multiplication two_phase_acc_reg=two_phase_acc_reg*two_phase_inc_reg + yl = _mm_moveldup_ps(two_phase_acc_reg); // Load yl with cr,cr,dr,dr + yh = _mm_movehdup_ps(two_phase_acc_reg); // Load yh with ci,ci,di,di + tmp1 = _mm_mul_ps(two_phase_inc_reg, yl); // tmp1 = ar*cr,ai*cr,br*dr,bi*dr + tmp3 = _mm_shuffle_ps(two_phase_inc_reg, two_phase_inc_reg, 0xB1); // Re-arrange x to be ai,ar,bi,br + tmp2 = _mm_mul_ps(tmp3, yh); // tmp2 = ai*ci,ar*ci,bi*di,br*di + two_phase_acc_reg = _mm_addsub_ps(tmp1, tmp2); // ar*cr-ai*ci, ai*cr+ar*ci, br*dr-bi*di, bi*dr+br*di + + // store four rotated in_common samples in the register b + b = _mm_packs_epi32(pc1, pc2);// convert from 32ic to 16ic + + //next two samples + _in_common += 2; + + for (int n_vec = 0; n_vec < num_a_vectors; n_vec++) + { + a = _mm_loadu_si128((__m128i*)&(_in_a[n_vec][number*4])); //load (2 byte imag, 2 byte real) x 4 into 128 bits reg + + c = _mm_mullo_epi16 (a, b); // a3.i*b3.i, a3.r*b3.r, .... + + c_sr = _mm_srli_si128 (c, 2); // Shift a right by imm8 bytes while shifting in zeros, and store the results in dst. + real = _mm_subs_epi16 (c, c_sr); + + b_sl = _mm_slli_si128(b, 2); // b3.r, b2.i .... + a_sl = _mm_slli_si128(a, 2); // a3.r, a2.i .... + + imag1 = _mm_mullo_epi16(a, b_sl); // a3.i*b3.r, .... + imag2 = _mm_mullo_epi16(b, a_sl); // b3.i*a3.r, .... + + imag = _mm_adds_epi16(imag1, imag2); + + realcacc[n_vec] = _mm_adds_epi16 (realcacc[n_vec], real); + imagcacc[n_vec] = _mm_adds_epi16 (imagcacc[n_vec], imag); + + } + } + + for (int n_vec=0;n_vec + *
  • Carles Fernandez Prades 2016 cfernandez at cttc dot cat + * + * + * Volk puppet for integrating the resampler into volk's test system + * + * ------------------------------------------------------------------------- + * + * Copyright (C) 2010-2015 (see AUTHORS file for a list of contributors) + * + * GNSS-SDR is a software defined Global Navigation + * Satellite Systems receiver + * + * This file is part of GNSS-SDR. + * + * GNSS-SDR is free software: you can redistribute it and/or modify + * it under the terms of the GNU General Public License as published by + * the Free Software Foundation, either version 3 of the License, or + * (at your option) any later version. + * + * GNSS-SDR is distributed in the hope that it will be useful, + * but WITHOUT ANY WARRANTY; without even the implied warranty of + * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the + * GNU General Public License for more details. + * + * You should have received a copy of the GNU General Public License + * along with GNSS-SDR. If not, see . + * + * ------------------------------------------------------------------------- + */ + +#ifndef INCLUDED_volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic_H +#define INCLUDED_volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic_H + +#include "volk_gnsssdr/volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn.h" +#include +#include +#include +#include + +#ifdef LV_HAVE_GENERIC +static inline void volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic_generic(lv_16sc_t* result, const lv_16sc_t* local_code, const lv_16sc_t* in, unsigned int num_points) +{ + // phases must be normalized. Phase rotator expects a complex exponential input! + float rem_carrier_phase_in_rad = 0.345; + float phase_step_rad = 0.123; + lv_32fc_t phase[1]; + phase[0] = lv_cmake(cos(rem_carrier_phase_in_rad), -sin(rem_carrier_phase_in_rad)); + lv_32fc_t phase_inc[1]; + phase_inc[0] = lv_cmake(cos(phase_step_rad), -sin(phase_step_rad)); + + int num_a_vectors = 3; + lv_16sc_t** in_a = (lv_16sc_t**)volk_gnsssdr_malloc(sizeof(lv_16sc_t*) * num_a_vectors, volk_gnsssdr_get_alignment()); + for(unsigned int n = 0; n < num_a_vectors; n++) + { + in_a[n] = (lv_16sc_t*)volk_gnsssdr_malloc(sizeof(lv_16sc_t) * num_points, volk_gnsssdr_get_alignment()); + memcpy(in_a[n], in, sizeof(lv_16sc_t) * num_points); + } + volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn_generic(result, local_code, phase_inc[0], phase,(const lv_16sc_t**) in_a, num_a_vectors, num_points); + + for(unsigned int n = 0; n < num_a_vectors; n++) + { + volk_gnsssdr_free(in_a[n]); + } + volk_gnsssdr_free(in_a); +} + +#endif // Generic + +#ifdef LV_HAVE_SSE3 +static inline void volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic_a_sse3(lv_16sc_t* result, const lv_16sc_t* local_code, const lv_16sc_t* in, unsigned int num_points) +{ + // phases must be normalized. Phase rotator expects a complex exponential input! + float rem_carrier_phase_in_rad = 0.345; + float phase_step_rad = 0.123; + lv_32fc_t phase[1]; + phase[0] = lv_cmake(cos(rem_carrier_phase_in_rad), -sin(rem_carrier_phase_in_rad)); + lv_32fc_t phase_inc[1]; + phase_inc[0] = lv_cmake(cos(phase_step_rad), -sin(phase_step_rad)); + + int num_a_vectors = 3; + lv_16sc_t** in_a = (lv_16sc_t**)volk_gnsssdr_malloc(sizeof(lv_16sc_t*) * num_a_vectors, volk_gnsssdr_get_alignment()); + for(unsigned int n = 0; n < num_a_vectors; n++) + { + in_a[n] = (lv_16sc_t*)volk_gnsssdr_malloc(sizeof(lv_16sc_t) * num_points, volk_gnsssdr_get_alignment()); + memcpy((lv_16sc_t*)in_a[n], (lv_16sc_t*)in, sizeof(lv_16sc_t) * num_points); + } + volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn_a_sse3(result, local_code, phase_inc[0], phase, (const lv_16sc_t**) in_a, num_a_vectors, num_points); + + for(unsigned int n = 0; n < num_a_vectors; n++) + { + volk_gnsssdr_free(in_a[n]); + } + volk_gnsssdr_free(in_a); +} + +#endif // SSE3 + +#ifdef LV_HAVE_SSE3 + +static inline void volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic_u_sse3(lv_16sc_t* result, const lv_16sc_t* local_code, const lv_16sc_t* in, unsigned int num_points) +{ + // phases must be normalized. Phase rotator expects a complex exponential input! + float rem_carrier_phase_in_rad = 0.345; + float phase_step_rad = 0.123; + lv_32fc_t phase[1]; + phase[0] = lv_cmake(cos(rem_carrier_phase_in_rad), -sin(rem_carrier_phase_in_rad)); + lv_32fc_t phase_inc[1]; + phase_inc[0] = lv_cmake(cos(phase_step_rad), -sin(phase_step_rad)); + + int num_a_vectors = 3; + lv_16sc_t** in_a = (lv_16sc_t**)volk_gnsssdr_malloc(sizeof(lv_16sc_t*) * num_a_vectors, volk_gnsssdr_get_alignment()); + for(unsigned int n = 0; n < num_a_vectors; n++) + { + in_a[n] = (lv_16sc_t*)volk_gnsssdr_malloc(sizeof(lv_16sc_t)*num_points, volk_gnsssdr_get_alignment()); + memcpy(in_a[n], in, sizeof(lv_16sc_t)*num_points); + } + volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn_u_sse3(result, local_code, phase_inc[0], phase, (const lv_16sc_t**) in_a, num_a_vectors, num_points); + + for(unsigned int n = 0; n < num_a_vectors; n++) + { + volk_gnsssdr_free(in_a[n]); + } + volk_gnsssdr_free(in_a); +} + +#endif // SSE3 + +#endif // INCLUDED_volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic_H + + diff --git a/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/lib/kernel_tests.h b/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/lib/kernel_tests.h index 0235218a1..97bf1ee97 100644 --- a/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/lib/kernel_tests.h +++ b/src/algorithms/libs/volk_gnsssdr_module/volk_gnsssdr/lib/kernel_tests.h @@ -81,6 +81,7 @@ std::vector init_test_list(volk_gnsssdr_test_params_t (VOLK_INIT_PUPP(volk_gnsssdr_16ic_resamplerpuppet_16ic, volk_gnsssdr_16ic_resampler_16ic, test_params)) (VOLK_INIT_PUPP(volk_gnsssdr_16ic_resamplerxnpuppet_16ic, volk_gnsssdr_16ic_xn_resampler_16ic_xn, test_params)) (VOLK_INIT_PUPP(volk_gnsssdr_16ic_x2_dotprodxnpuppet_16ic, volk_gnsssdr_16ic_x2_dot_prod_16ic_xn, test_params)) + (VOLK_INIT_PUPP(volk_gnsssdr_16ic_x2_rotator_dotprodxnpuppet_16ic, volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn, test_params)) ; return test_cases; diff --git a/src/algorithms/signal_source/adapters/osmosdr_signal_source.h b/src/algorithms/signal_source/adapters/osmosdr_signal_source.h index 0fffe7b18..0d1bbbe2e 100644 --- a/src/algorithms/signal_source/adapters/osmosdr_signal_source.h +++ b/src/algorithms/signal_source/adapters/osmosdr_signal_source.h @@ -81,7 +81,7 @@ public: private: std::string role_; - // UHD SETTINGS + // Front-end settings bool AGC_enabled_; double sample_rate_; diff --git a/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.cc b/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.cc index 05a4f20a1..8cbc574dd 100644 --- a/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.cc +++ b/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.cc @@ -44,18 +44,18 @@ bool cpu_multicorrelator_16sc::init( // ALLOCATE MEMORY FOR INTERNAL vectors size_t size = max_signal_length_samples * sizeof(lv_16sc_t); - // NCO signal - d_nco_in = static_cast(volk_gnsssdr_malloc(size, volk_gnsssdr_get_alignment())); - - // Doppler-free signal - d_sig_doppler_wiped = static_cast(volk_gnsssdr_malloc(size, volk_gnsssdr_get_alignment())); + // NCO signal (not needed if the rotator+dot_product kernel is used) + //d_nco_in = static_cast(volk_gnsssdr_malloc(size, volk_gnsssdr_get_alignment())); + // Doppler-free signal (not needed if the rotator+dot_product kernel is used) + //d_sig_doppler_wiped = static_cast(volk_gnsssdr_malloc(size, volk_gnsssdr_get_alignment())); + d_n_correlators = n_correlators; + d_tmp_code_phases_chips = static_cast(volk_gnsssdr_malloc(n_correlators*sizeof(float), volk_gnsssdr_get_alignment())); d_local_codes_resampled = new lv_16sc_t*[n_correlators]; for (int n = 0; n < n_correlators; n++) { d_local_codes_resampled[n] = static_cast(volk_gnsssdr_malloc(size, volk_gnsssdr_get_alignment())); } - d_n_correlators = n_correlators; return true; } @@ -81,26 +81,21 @@ bool cpu_multicorrelator_16sc::set_input_output_vectors(lv_16sc_t* corr_out, con return true; } - - void cpu_multicorrelator_16sc::update_local_code(int correlator_length_samples,float rem_code_phase_chips, float code_phase_step_chips) { - float *tmp_code_phases_chips; - tmp_code_phases_chips = static_cast(volk_gnsssdr_malloc(d_n_correlators*sizeof(float), volk_gnsssdr_get_alignment())); + for (int n = 0; n < d_n_correlators; n++) { - tmp_code_phases_chips[n] = d_shifts_chips[n] - rem_code_phase_chips; + d_tmp_code_phases_chips[n] = d_shifts_chips[n] - rem_code_phase_chips; } volk_gnsssdr_16ic_xn_resampler_16ic_xn(d_local_codes_resampled, d_local_code_in, - tmp_code_phases_chips, + d_tmp_code_phases_chips, code_phase_step_chips, correlator_length_samples, d_n_correlators, d_code_length_chips); - - volk_gnsssdr_free(tmp_code_phases_chips); } @@ -111,11 +106,14 @@ bool cpu_multicorrelator_16sc::Carrier_wipeoff_multicorrelator_resampler( float code_phase_step_chips, int signal_length_samples) { + update_local_code(signal_length_samples, rem_code_phase_chips, code_phase_step_chips); lv_32fc_t phase_offset_as_complex[1]; phase_offset_as_complex[0] = lv_cmake(std::cos(rem_carrier_phase_in_rad), -std::sin(rem_carrier_phase_in_rad)); - volk_gnsssdr_16ic_s32fc_x2_rotator_16ic(d_sig_doppler_wiped, d_sig_in, std::exp(lv_32fc_t(0, -phase_step_rad)), phase_offset_as_complex, signal_length_samples); - update_local_code(signal_length_samples, rem_code_phase_chips, code_phase_step_chips); - volk_gnsssdr_16ic_x2_dot_prod_16ic_xn(d_corr_out, d_sig_doppler_wiped, (const lv_16sc_t**)d_local_codes_resampled, d_n_correlators, signal_length_samples); + //replaced by integrated rotator + dot_product kernel + //volk_gnsssdr_16ic_s32fc_x2_rotator_16ic(d_sig_doppler_wiped, d_sig_in, std::exp(lv_32fc_t(0, -phase_step_rad)), phase_offset_as_complex, signal_length_samples); + //volk_gnsssdr_16ic_x2_dot_prod_16ic_xn(d_corr_out, d_sig_doppler_wiped, (const lv_16sc_t**)d_local_codes_resampled, d_n_correlators, signal_length_samples); + volk_gnsssdr_16ic_x2_rotator_dot_prod_16ic_xn(d_corr_out, d_sig_in, std::exp(lv_32fc_t(0, -phase_step_rad)), phase_offset_as_complex, (const lv_16sc_t**)d_local_codes_resampled, d_n_correlators, signal_length_samples); + return true; } @@ -123,8 +121,8 @@ bool cpu_multicorrelator_16sc::Carrier_wipeoff_multicorrelator_resampler( cpu_multicorrelator_16sc::cpu_multicorrelator_16sc() { d_sig_in = NULL; - d_nco_in = NULL; - d_sig_doppler_wiped = NULL; + //d_nco_in = NULL; + //d_sig_doppler_wiped = NULL; d_local_code_in = NULL; d_shifts_chips = NULL; d_corr_out = NULL; @@ -136,8 +134,9 @@ cpu_multicorrelator_16sc::cpu_multicorrelator_16sc() bool cpu_multicorrelator_16sc::free() { // Free memory - if (d_sig_doppler_wiped != NULL) volk_gnsssdr_free(d_sig_doppler_wiped); - if (d_nco_in != NULL) volk_gnsssdr_free(d_nco_in); + //if (d_sig_doppler_wiped != NULL) volk_gnsssdr_free(d_sig_doppler_wiped); + //if (d_nco_in != NULL) volk_gnsssdr_free(d_nco_in); + if (d_tmp_code_phases_chips != NULL) volk_gnsssdr_free(d_tmp_code_phases_chips); for (int n = 0; n < d_n_correlators; n++) { volk_gnsssdr_free(d_local_codes_resampled[n]); diff --git a/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.h b/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.h index 1c5163cc2..6a6082168 100644 --- a/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.h +++ b/src/algorithms/tracking/libs/cpu_multicorrelator_16sc.h @@ -55,9 +55,10 @@ public: private: // Allocate the device input vectors const lv_16sc_t *d_sig_in; - lv_16sc_t *d_nco_in; + //lv_16sc_t *d_nco_in; + float *d_tmp_code_phases_chips; lv_16sc_t **d_local_codes_resampled; - lv_16sc_t *d_sig_doppler_wiped; + //lv_16sc_t *d_sig_doppler_wiped; const lv_16sc_t *d_local_code_in; lv_16sc_t *d_corr_out; float *d_shifts_chips; diff --git a/src/core/interfaces/gnss_block_interface.h b/src/core/interfaces/gnss_block_interface.h index bb45902c1..dd7f1558e 100644 --- a/src/core/interfaces/gnss_block_interface.h +++ b/src/core/interfaces/gnss_block_interface.h @@ -52,29 +52,29 @@ class GNSSBlockInterface { public: - virtual ~GNSSBlockInterface() - {} - virtual std::string role() = 0; - virtual std::string implementation() = 0; - virtual size_t item_size() = 0; - virtual void connect(gr::top_block_sptr top_block) = 0; - virtual void disconnect(gr::top_block_sptr top_block) = 0; + virtual ~GNSSBlockInterface() + {} + virtual std::string role() = 0; + virtual std::string implementation() = 0; + virtual size_t item_size() = 0; + virtual void connect(gr::top_block_sptr top_block) = 0; + virtual void disconnect(gr::top_block_sptr top_block) = 0; - virtual gr::basic_block_sptr get_left_block() = 0; - virtual gr::basic_block_sptr get_right_block() = 0; + virtual gr::basic_block_sptr get_left_block() = 0; + virtual gr::basic_block_sptr get_right_block() = 0; - virtual gr::basic_block_sptr get_left_block(int RF_channel) - { - assert(RF_channel >= 0); - if (RF_channel == 0){}; // avoid unused param warning - return NULL; // added to support raw array access (non pure virtual to allow left unimplemented)= 0; - } - virtual gr::basic_block_sptr get_right_block(int RF_channel) - { - assert(RF_channel >= 0); - if (RF_channel == 0){}; // avoid unused param warning - return NULL; // added to support raw array access (non pure virtual to allow left unimplemented)= 0; - } + virtual gr::basic_block_sptr get_left_block(int RF_channel) + { + assert(RF_channel >= 0); + if (RF_channel == 0){}; // avoid unused param warning + return nullptr; // added to support raw array access (non pure virtual to allow left unimplemented)= 0; + } + virtual gr::basic_block_sptr get_right_block(int RF_channel) + { + assert(RF_channel >= 0); + if (RF_channel == 0){}; // avoid unused param warning + return nullptr; // added to support raw array access (non pure virtual to allow left unimplemented)= 0; + } }; #endif /*GNSS_SDR_GNSS_BLOCK_INTERFACE_H_*/ diff --git a/src/tests/gnss_block/gnss_block_factory_test.cc b/src/tests/gnss_block/gnss_block_factory_test.cc index 0e425870a..812647e8d 100644 --- a/src/tests/gnss_block/gnss_block_factory_test.cc +++ b/src/tests/gnss_block/gnss_block_factory_test.cc @@ -62,23 +62,6 @@ TEST(GNSS_Block_Factory_Test, InstantiateFileSignalSource) EXPECT_STREQ("File_Signal_Source", signal_source->implementation().c_str()); } -/* -TEST(GNSS_Block_Factory_Test, InstantiateUHDSignalSource) -{ - std::shared_ptr configuration = std::make_shared(); - configuration->set_property("SignalSource.implementation", "UHD_Signal_Source"); - configuration->set_property("SignalSource.item_type", "gr_complex"); - configuration->set_property("SignalSource.device_address", "192.168.40.2"); - gr::msg_queue::sptr queue = gr::msg_queue::make(0); - // Example of a factory created with auto - auto factory = new GNSSBlockFactory(); - // Example of a block created with auto - auto signal_source = factory->GetSignalSource(configuration, queue); - - EXPECT_STREQ("SignalSource", signal_source->role().c_str()); - EXPECT_STREQ("UHD_Signal_Source", signal_source->implementation().c_str()); -} -*/ TEST(GNSS_Block_Factory_Test, InstantiateWrongSignalSource) {