1
0
mirror of https://github.com/gnss-sdr/gnss-sdr synced 2025-10-26 13:07:39 +00:00

Added ability to generate real valued codes

Only done for GPS L1 C/A and Galileo E1 OS for now. Also added a
cpu_multicorrelator_real_codes class that performs code correlation
using real-valued local codes
This commit is contained in:
Cillian O'Driscoll
2017-09-11 15:21:05 +01:00
parent 676c1506da
commit e87522880e
11 changed files with 530 additions and 44 deletions

View File

@@ -73,7 +73,7 @@ void galileo_e1_code_gen_int(int* _dest, char _Signal[3], signed int _prn)
void galileo_e1_sinboc_11_gen(std::complex<float>* _dest, int* _prn, unsigned int _length_out) void galileo_e1_sinboc_11_gen_int(int* _dest, int* _prn, unsigned int _length_out)
{ {
const unsigned int _length_in = Galileo_E1_B_CODE_LENGTH_CHIPS; const unsigned int _length_in = Galileo_E1_B_CODE_LENGTH_CHIPS;
unsigned int _period = static_cast<unsigned int>( _length_out / _length_in ); unsigned int _period = static_cast<unsigned int>( _length_out / _length_in );
@@ -81,18 +81,18 @@ void galileo_e1_sinboc_11_gen(std::complex<float>* _dest, int* _prn, unsigned in
{ {
for (unsigned int j = 0; j < (_period / 2); j++) for (unsigned int j = 0; j < (_period / 2); j++)
{ {
_dest[i * _period + j] = std::complex<float>(static_cast<float>(_prn[i]), 0.0); _dest[i * _period + j] = _prn[i];
} }
for (unsigned int j = (_period / 2); j < _period; j++) for (unsigned int j = (_period / 2); j < _period; j++)
{ {
_dest[i * _period + j] = std::complex<float>(static_cast<float>(- _prn[i]), 0.0); _dest[i * _period + j] = - _prn[i];
} }
} }
} }
void galileo_e1_sinboc_61_gen(std::complex<float>* _dest, int* _prn, unsigned int _length_out) void galileo_e1_sinboc_61_gen_int(int* _dest, int* _prn, unsigned int _length_out)
{ {
const unsigned int _length_in = Galileo_E1_B_CODE_LENGTH_CHIPS; const unsigned int _length_in = Galileo_E1_B_CODE_LENGTH_CHIPS;
unsigned int _period = static_cast<unsigned int>(_length_out / _length_in); unsigned int _period = static_cast<unsigned int>(_length_out / _length_in);
@@ -101,42 +101,44 @@ void galileo_e1_sinboc_61_gen(std::complex<float>* _dest, int* _prn, unsigned in
{ {
for (unsigned int j = 0; j < _period; j += 2) for (unsigned int j = 0; j < _period; j += 2)
{ {
_dest[i * _period + j] = std::complex<float>(static_cast<float>(_prn[i]), 0.0); _dest[i * _period + j] = _prn[i];
} }
for (unsigned int j = 1; j < _period; j += 2) for (unsigned int j = 1; j < _period; j += 2)
{ {
_dest[i * _period + j] = std::complex<float>(static_cast<float>(- _prn[i]), 0.0); _dest[i * _period + j] = - _prn[i];
} }
} }
} }
void galileo_e1_gen(std::complex<float>* _dest, int* _prn, char _Signal[3]) void galileo_e1_gen_float(float* _dest, int* _prn, char _Signal[3])
{ {
std::string _galileo_signal = _Signal; std::string _galileo_signal = _Signal;
const unsigned int _codeLength = 12 * Galileo_E1_B_CODE_LENGTH_CHIPS; const unsigned int _codeLength = 12 * Galileo_E1_B_CODE_LENGTH_CHIPS;
const float alpha = sqrt(10.0 / 11.0); const float alpha = sqrt(10.0 / 11.0);
const float beta = sqrt(1.0 / 11.0); const float beta = sqrt(1.0 / 11.0);
std::complex<float> sinboc_11[12 * 4092]; // _codeLength not accepted by Clang int sinboc_11[12 * 4092]; // _codeLength not accepted by Clang
std::complex<float> sinboc_61[12 * 4092]; int sinboc_61[12 * 4092];
galileo_e1_sinboc_11_gen(sinboc_11, _prn, _codeLength); //generate sinboc(1,1) 12 samples per chip galileo_e1_sinboc_11_gen_int(sinboc_11, _prn, _codeLength); //generate sinboc(1,1) 12 samples per chip
galileo_e1_sinboc_61_gen(sinboc_61, _prn, _codeLength); //generate sinboc(6,1) 12 samples per chip galileo_e1_sinboc_61_gen_int(sinboc_61, _prn, _codeLength); //generate sinboc(6,1) 12 samples per chip
if (_galileo_signal.rfind("1B") != std::string::npos && _galileo_signal.length() >= 2) if (_galileo_signal.rfind("1B") != std::string::npos && _galileo_signal.length() >= 2)
{ {
for (unsigned int i = 0; i < _codeLength; i++) for (unsigned int i = 0; i < _codeLength; i++)
{ {
_dest[i] = alpha * sinboc_11[i] + beta * sinboc_61[i]; _dest[i] = alpha * static_cast<float>(sinboc_11[i]) +
beta * static_cast<float>(sinboc_61[i]);
} }
} }
else if (_galileo_signal.rfind("1C") != std::string::npos && _galileo_signal.length() >= 2) else if (_galileo_signal.rfind("1C") != std::string::npos && _galileo_signal.length() >= 2)
{ {
for (unsigned int i = 0; i < _codeLength; i++) for (unsigned int i = 0; i < _codeLength; i++)
{ {
_dest[i] = alpha * sinboc_11[i] - beta * sinboc_61[i]; _dest[i] = alpha * static_cast<float>(sinboc_11[i]) -
beta * static_cast<float>(sinboc_61[i]);
} }
} }
else else
@@ -144,8 +146,7 @@ void galileo_e1_gen(std::complex<float>* _dest, int* _prn, char _Signal[3])
} }
void galileo_e1_code_gen_float_sampled(float* _dest, char _Signal[3],
void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signal[3],
bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift, bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift,
bool _secondary_flag) bool _secondary_flag)
{ {
@@ -164,23 +165,29 @@ void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signa
galileo_e1_code_gen_int(primary_code_E1_chips, _Signal, _prn); //generate Galileo E1 code, 1 sample per chip galileo_e1_code_gen_int(primary_code_E1_chips, _Signal, _prn); //generate Galileo E1 code, 1 sample per chip
std::complex<float>* _signal_E1; float* _signal_E1;
_codeLength = _samplesPerChip * Galileo_E1_B_CODE_LENGTH_CHIPS; _codeLength = _samplesPerChip * Galileo_E1_B_CODE_LENGTH_CHIPS;
_signal_E1 = new std::complex<float>[_codeLength]; _signal_E1 = new float[_codeLength];
if (_cboc == true) if (_cboc == true)
{ {
galileo_e1_gen(_signal_E1, primary_code_E1_chips, _Signal); //generate cboc 12 samples per chip galileo_e1_gen_float(_signal_E1, primary_code_E1_chips, _Signal); //generate cboc 12 samples per chip
} }
else else
{ {
galileo_e1_sinboc_11_gen(_signal_E1, primary_code_E1_chips, _codeLength); //generate sinboc(1,1) 2 samples per chip int _signal_E1_int[_codeLength];
galileo_e1_sinboc_11_gen_int(_signal_E1_int, primary_code_E1_chips, _codeLength); //generate sinboc(1,1) 2 samples per chip
for( unsigned int ii = 0; ii < _codeLength; ++ii )
{
_signal_E1[ii] = static_cast< float >( _signal_E1_int[ii] );
}
} }
if (_fs != _samplesPerChip * _codeFreqBasis) if (_fs != _samplesPerChip * _codeFreqBasis)
{ {
std::complex<float>* _resampled_signal = new std::complex<float>[_samplesPerCode]; float* _resampled_signal = new float[_samplesPerCode];
resampler(_signal_E1, _resampled_signal, _samplesPerChip * _codeFreqBasis, _fs, resampler(_signal_E1, _resampled_signal, _samplesPerChip * _codeFreqBasis, _fs,
_codeLength, _samplesPerCode); //resamples code to fs _codeLength, _samplesPerCode); //resamples code to fs
@@ -192,8 +199,8 @@ void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signa
if (_galileo_signal.rfind("1C") != std::string::npos && _galileo_signal.length() >= 2 && _secondary_flag) if (_galileo_signal.rfind("1C") != std::string::npos && _galileo_signal.length() >= 2 && _secondary_flag)
{ {
std::complex<float>* _signal_E1C_secondary = new std::complex<float> float* _signal_E1C_secondary = new float[
[static_cast<int>(Galileo_E1_C_SECONDARY_CODE_LENGTH) static_cast<int>(Galileo_E1_C_SECONDARY_CODE_LENGTH)
* _samplesPerCode]; * _samplesPerCode];
for (unsigned int i = 0; i < static_cast<unsigned int>(Galileo_E1_C_SECONDARY_CODE_LENGTH); i++) for (unsigned int i = 0; i < static_cast<unsigned int>(Galileo_E1_C_SECONDARY_CODE_LENGTH); i++)
@@ -202,7 +209,7 @@ void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signa
{ {
_signal_E1C_secondary[i*_samplesPerCode + k] = _signal_E1[k] _signal_E1C_secondary[i*_samplesPerCode + k] = _signal_E1[k]
* (Galileo_E1_C_SECONDARY_CODE.at(i) == '0' * (Galileo_E1_C_SECONDARY_CODE.at(i) == '0'
? std::complex<float>(1,0) : std::complex<float>(-1,0)); ? 1.0f : -1.0f);
} }
} }
@@ -220,6 +227,36 @@ void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signa
delete[] _signal_E1; delete[] _signal_E1;
} }
void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signal[3],
bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift,
bool _secondary_flag)
{
std::string _galileo_signal = _Signal;
const int _codeFreqBasis = Galileo_E1_CODE_CHIP_RATE_HZ; //Hz
unsigned int _samplesPerCode = static_cast<unsigned int>(
static_cast<double>(_fs) / (static_cast<double>(_codeFreqBasis )
/ static_cast<double>(Galileo_E1_B_CODE_LENGTH_CHIPS)));
if (_galileo_signal.rfind("1C") != std::string::npos && _galileo_signal.length() >= 2 && _secondary_flag)
{
_samplesPerCode *= static_cast<int>(Galileo_E1_C_SECONDARY_CODE_LENGTH);
}
float real_code[_samplesPerCode];
galileo_e1_code_gen_float_sampled( real_code, _Signal, _cboc, _prn, _fs, _chip_shift, _secondary_flag );
for( unsigned int ii = 0; ii < _samplesPerCode; ++ii )
{
_dest[ii] = std::complex< float >( real_code[ii], 0.0f );
}
}
void galileo_e1_code_gen_float_sampled(float* _dest, char _Signal[3],
bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift)
{
galileo_e1_code_gen_float_sampled(_dest, _Signal, _cboc, _prn, _fs, _chip_shift, false);
}
void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signal[3], void galileo_e1_code_gen_complex_sampled(std::complex<float>* _dest, char _Signal[3],
bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift) bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift)

View File

@@ -36,30 +36,22 @@
/*! /*!
* \brief This function generates Galileo E1 code (one sample per chip). * \brief This function generates Galileo E1 code (can select E1B or E1C, cboc or sinboc
* and the sample frequency _fs).
* *
*/ */
void galileo_e1_code_gen_int(int* _dest, char _Signal[3], signed int _prn); void galileo_e1_code_gen_float_sampled(float* _dest, char _Signal[3],
bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift,
bool _secondary_flag);
/*! /*!
* \brief This function generates Galileo E1 sinboc(1,1) code (minimum 2 samples per chip), * \brief This function generates Galileo E1 code (can select E1B or E1C, cboc or sinboc
* the _codeLength variable must be a multiple of 2*4092. * and the sample frequency _fs).
* *
*/ */
void galileo_e1_sinboc_11_gen(std::complex<float>* _dest, int* _prn, void galileo_e1_code_gen_float_sampled(float* _dest, char _Signal[3],
unsigned int _codeLength); bool _cboc, unsigned int _prn, signed int _fs, unsigned int _chip_shift);
/*!
* \brief This function generates Galileo E1 sinboc(6,1) code (minimum 12 samples per chip),
* the _codeLength variable must be a multiple of 12*4092.
*
*/
void galileo_e1_sinboc_61_gen(std::complex<float>* _dest, int* _prn,
unsigned int _codeLength);
/*!
* \brief This function generates Galileo E1 cboc code (12 samples per chip).
*
*/
void galileo_e1_cboc_gen(std::complex<float>* _dest, int* _prn, char _Signal[3]);
/*! /*!
* \brief This function generates Galileo E1 code (can select E1B or E1C, cboc or sinboc * \brief This function generates Galileo E1 code (can select E1B or E1C, cboc or sinboc
* and the sample frequency _fs). * and the sample frequency _fs).

View File

@@ -156,6 +156,28 @@ void hex_to_binary_converter(int * _dest, char _from)
} }
} }
void resampler(float* _from, float* _dest, float _fs_in,
float _fs_out, unsigned int _length_in, unsigned int _length_out)
{
unsigned int _codeValueIndex;
float aux;
//--- Find time constants --------------------------------------------------
const float _t_in = 1 / _fs_in; // Incoming sampling period in sec
const float _t_out = 1 / _fs_out; // Out sampling period in sec
for (unsigned int i = 0; i < _length_out - 1; i++)
{
//=== Digitizing =======================================================
//--- compute index array to read sampled values -------------------------
//_codeValueIndex = ceil((_t_out * ((float)i + 1)) / _t_in) - 1;
aux = (_t_out * (i + 1)) / _t_in;
_codeValueIndex = auxCeil2(aux) - 1;
//if repeat the chip -> upsample by nearest neighborhood interpolation
_dest[i] = _from[_codeValueIndex];
}
//--- Correct the last index (due to number rounding issues) -----------
_dest[_length_out - 1] = _from[_length_in - 1];
}
void resampler(std::complex<float>* _from, std::complex<float>* _dest, float _fs_in, void resampler(std::complex<float>* _from, std::complex<float>* _dest, float _fs_in,
float _fs_out, unsigned int _length_in, unsigned int _length_out) float _fs_out, unsigned int _length_in, unsigned int _length_out)

View File

@@ -60,6 +60,13 @@ void complex_exp_gen_conj(std::complex<float>* _dest, double _f, double _fs,
*/ */
void hex_to_binary_converter(int * _dest, char _from); void hex_to_binary_converter(int * _dest, char _from);
/*!
* \brief This function resamples a sequence of float values.
*
*/
void resampler(float* _from, float* _dest,
float _fs_in, float _fs_out, unsigned int _length_in,
unsigned int _length_out);
/*! /*!
* \brief This function resamples a sequence of complex values. * \brief This function resamples a sequence of complex values.
* *

View File

@@ -34,7 +34,7 @@
auto auxCeil = [](float x){ return static_cast<int>(static_cast<long>((x)+1)); }; auto auxCeil = [](float x){ return static_cast<int>(static_cast<long>((x)+1)); };
void gps_l1_ca_code_gen_complex(std::complex<float>* _dest, signed int _prn, unsigned int _chip_shift) void gps_l1_ca_code_gen_int(int* _dest, signed int _prn, unsigned int _chip_shift)
{ {
const unsigned int _code_length = 1023; const unsigned int _code_length = 1023;
bool G1[_code_length]; bool G1[_code_length];
@@ -102,17 +102,44 @@ void gps_l1_ca_code_gen_complex(std::complex<float>* _dest, signed int _prn, uns
aux = G1[(lcv + _chip_shift) % _code_length]^G2[delay]; aux = G1[(lcv + _chip_shift) % _code_length]^G2[delay];
if(aux == true) if(aux == true)
{ {
_dest[lcv] = std::complex<float>(1, 0); _dest[lcv] = 1;
} }
else else
{ {
_dest[lcv] = std::complex<float>(-1, 0); _dest[lcv] = -1;
} }
delay++; delay++;
delay %= _code_length; delay %= _code_length;
} }
} }
void gps_l1_ca_code_gen_float(float* _dest, signed int _prn, unsigned int _chip_shift)
{
unsigned int _code_length = 1023;
int ca_code_int[ _code_length ];
gps_l1_ca_code_gen_int( ca_code_int, _prn, _chip_shift );
for( unsigned int ii = 0; ii < _code_length; ++ii )
{
_dest[ii] = static_cast<float>( ca_code_int[ii] );
}
}
void gps_l1_ca_code_gen_complex(std::complex<float>* _dest, signed int _prn, unsigned int _chip_shift)
{
unsigned int _code_length = 1023;
int ca_code_int[ _code_length ];
gps_l1_ca_code_gen_int( ca_code_int, _prn, _chip_shift );
for( unsigned int ii = 0; ii < _code_length; ++ii )
{
_dest[ii] = std::complex< float >( static_cast<float>(ca_code_int[ii]), 0.0f );
}
}
/* /*

View File

@@ -35,6 +35,12 @@
#include <complex> #include <complex>
//!Generates int GPS L1 C/A code for the desired SV ID and code shift
void gps_l1_ca_code_gen_int(int* _dest, signed int _prn, unsigned int _chip_shift);
//!Generates float GPS L1 C/A code for the desired SV ID and code shift
void gps_l1_ca_code_gen_float(float* _dest, signed int _prn, unsigned int _chip_shift);
//!Generates complex GPS L1 C/A code for the desired SV ID and code shift, and sampled to specific sampling frequency //!Generates complex GPS L1 C/A code for the desired SV ID and code shift, and sampled to specific sampling frequency
void gps_l1_ca_code_gen_complex(std::complex<float>* _dest, signed int _prn, unsigned int _chip_shift); void gps_l1_ca_code_gen_complex(std::complex<float>* _dest, signed int _prn, unsigned int _chip_shift);

View File

@@ -33,6 +33,7 @@ endif(ENABLE_CUDA)
set(TRACKING_LIB_SOURCES set(TRACKING_LIB_SOURCES
cpu_multicorrelator.cc cpu_multicorrelator.cc
cpu_multicorrelator_real_codes.cc
cpu_multicorrelator_16sc.cc cpu_multicorrelator_16sc.cc
lock_detectors.cc lock_detectors.cc
tcp_communication.cc tcp_communication.cc

View File

@@ -0,0 +1,147 @@
/*!
* \file cpu_multicorrelator_real_codes.cc
* \brief High optimized CPU vector multiTAP correlator class with real-valued local codes
* \authors <ul>
* <li> Javier Arribas, 2015. jarribas(at)cttc.es
* <li> Cillian O'Driscoll, 2017. cillian.odriscoll(at)gmail.com
* </ul>
*
* Class that implements a high optimized vector multiTAP correlator class for CPUs
*
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2017 (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 <http://www.gnu.org/licenses/>.
*
* -------------------------------------------------------------------------
*/
#include "cpu_multicorrelator_real_codes.h"
#include <cmath>
#include <iostream>
#include <volk_gnsssdr/volk_gnsssdr.h>
cpu_multicorrelator_real_codes::cpu_multicorrelator_real_codes()
{
d_sig_in = nullptr;
d_local_code_in = nullptr;
d_shifts_chips = nullptr;
d_corr_out = nullptr;
d_local_codes_resampled = nullptr;
d_code_length_chips = 0;
d_n_correlators = 0;
}
cpu_multicorrelator_real_codes::~cpu_multicorrelator_real_codes()
{
if(d_local_codes_resampled != nullptr)
{
cpu_multicorrelator_real_codes::free();
}
}
bool cpu_multicorrelator_real_codes::init(
int max_signal_length_samples,
int n_correlators)
{
// ALLOCATE MEMORY FOR INTERNAL vectors
size_t size = max_signal_length_samples * sizeof(float);
d_local_codes_resampled = static_cast<float**>(volk_gnsssdr_malloc(n_correlators * sizeof(float*), volk_gnsssdr_get_alignment()));
for (int n = 0; n < n_correlators; n++)
{
d_local_codes_resampled[n] = static_cast<float*>(volk_gnsssdr_malloc(size, volk_gnsssdr_get_alignment()));
}
d_n_correlators = n_correlators;
return true;
}
bool cpu_multicorrelator_real_codes::set_local_code_and_taps(
int code_length_chips,
const float* local_code_in,
float *shifts_chips)
{
d_local_code_in = local_code_in;
d_shifts_chips = shifts_chips;
d_code_length_chips = code_length_chips;
return true;
}
bool cpu_multicorrelator_real_codes::set_input_output_vectors(std::complex<float>* corr_out, const std::complex<float>* sig_in)
{
// Save CPU pointers
d_sig_in = sig_in;
d_corr_out = corr_out;
return true;
}
void cpu_multicorrelator_real_codes::update_local_code(int correlator_length_samples, float rem_code_phase_chips, float code_phase_step_chips)
{
volk_gnsssdr_32f_xn_resampler_32f_xn(d_local_codes_resampled,
d_local_code_in,
rem_code_phase_chips,
code_phase_step_chips,
d_shifts_chips,
d_code_length_chips,
d_n_correlators,
correlator_length_samples);
}
bool cpu_multicorrelator_real_codes::Carrier_wipeoff_multicorrelator_resampler(
float rem_carrier_phase_in_rad,
float phase_step_rad,
float rem_code_phase_chips,
float code_phase_step_chips,
int signal_length_samples)
{
update_local_code(signal_length_samples, rem_code_phase_chips, code_phase_step_chips);
// Regenerate phase at each call in order to avoid numerical issues
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));
// call VOLK_GNSSSDR kernel
volk_gnsssdr_32fc_32f_rotator_dot_prod_32fc_xn(d_corr_out, d_sig_in, std::exp(lv_32fc_t(0, - phase_step_rad)), phase_offset_as_complex, (const float**)d_local_codes_resampled, d_n_correlators, signal_length_samples);
return true;
}
bool cpu_multicorrelator_real_codes::free()
{
// Free memory
if (d_local_codes_resampled != nullptr)
{
for (int n = 0; n < d_n_correlators; n++)
{
volk_gnsssdr_free(d_local_codes_resampled[n]);
}
volk_gnsssdr_free(d_local_codes_resampled);
d_local_codes_resampled = nullptr;
}
return true;
}

View File

@@ -0,0 +1,70 @@
/*!
* \file cpu_multicorrelator_real_codes.h
* \brief High optimized CPU vector multiTAP correlator class using real-valued local codes
* \authors <ul>
* <li> Javier Arribas, 2015. jarribas(at)cttc.es
* <li> Cillian O'Driscoll, 2017, cillian.odriscoll(at)gmail.com
* </ul>
*
* Class that implements a high optimized vector multiTAP correlator class for CPUs
*
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2017 (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 <http://www.gnu.org/licenses/>.
*
* -------------------------------------------------------------------------
*/
#ifndef GNSS_SDR_CPU_MULTICORRELATOR_REAL_CODES_H_
#define GNSS_SDR_CPU_MULTICORRELATOR_REAL_CODES_H_
#include <complex>
/*!
* \brief Class that implements carrier wipe-off and correlators.
*/
class cpu_multicorrelator_real_codes
{
public:
cpu_multicorrelator_real_codes();
~cpu_multicorrelator_real_codes();
bool init(int max_signal_length_samples, int n_correlators);
bool set_local_code_and_taps(int code_length_chips, const float* local_code_in, float *shifts_chips);
bool set_input_output_vectors(std::complex<float>* corr_out, const std::complex<float>* sig_in);
void update_local_code(int correlator_length_samples, float rem_code_phase_chips, float code_phase_step_chips);
bool Carrier_wipeoff_multicorrelator_resampler(float rem_carrier_phase_in_rad, float phase_step_rad, float rem_code_phase_chips, float code_phase_step_chips, int signal_length_samples);
bool free();
private:
// Allocate the device input vectors
const std::complex<float> *d_sig_in;
float **d_local_codes_resampled;
const float *d_local_code_in;
std::complex<float> *d_corr_out;
float *d_shifts_chips;
int d_code_length_chips;
int d_n_correlators;
};
#endif /* CPU_MULTICORRELATOR_REAL_CODES_H_ */

View File

@@ -0,0 +1,174 @@
/*!
* \file cpu_multicorrelator_real_codes_test.cc
* \brief This file implements timing tests for the FFT.
* \author Carles Fernandez-Prades, 2016. cfernandez(at)cttc.es
*
*
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2016 (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 <http://www.gnu.org/licenses/>.
*
* -------------------------------------------------------------------------
*/
#include <chrono>
#include <complex>
#include <random>
#include <thread>
#include <gtest/gtest.h>
#include <gflags/gflags.h>
#include <gnuradio/gr_complex.h>
#include <volk_gnsssdr/volk_gnsssdr.h>
#include "cpu_multicorrelator_real_codes.h"
#include "gps_sdr_signal_processing.h"
#include "GPS_L1_CA.h"
DEFINE_int32(cpu_multicorrelator_real_codes_iterations_test, 1000, "Number of averaged iterations in CPU multicorrelator test timing test");
DEFINE_int32(cpu_multicorrelator_real_codes_max_threads_test, 12, "Number of maximum concurrent correlators in CPU multicorrelator test timing test");
void run_correlator_cpu(cpu_multicorrelator_real_codes* correlator,
float d_rem_carrier_phase_rad,
float d_carrier_phase_step_rad,
float d_code_phase_step_chips,
float d_rem_code_phase_chips,
int correlation_size)
{
for(int k = 0; k < FLAGS_cpu_multicorrelator_real_codes_iterations_test; k++)
{
correlator->Carrier_wipeoff_multicorrelator_resampler(d_rem_carrier_phase_rad,
d_carrier_phase_step_rad,
d_code_phase_step_chips,
d_rem_code_phase_chips,
correlation_size);
}
}
TEST(CpuMulticorrelatorRealCodesTest, MeasureExecutionTime)
{
std::chrono::time_point<std::chrono::system_clock> start, end;
std::chrono::duration<double> elapsed_seconds(0);
int max_threads = FLAGS_cpu_multicorrelator_real_codes_max_threads_test;
std::vector<std::thread> thread_pool;
cpu_multicorrelator_real_codes* correlator_pool[max_threads];
unsigned int correlation_sizes [3] = { 2048, 4096, 8192};
double execution_times [3];
float* d_ca_code;
gr_complex* in_cpu;
gr_complex* d_correlator_outs;
int d_n_correlator_taps = 3;
int d_vector_length = correlation_sizes[2]; //max correlation size to allocate all the necessary memory
float* d_local_code_shift_chips;
//allocate host memory
// Get space for a vector with the C/A code replica sampled 1x/chip
d_ca_code = static_cast<float*>(volk_gnsssdr_malloc(static_cast<int>(GPS_L1_CA_CODE_LENGTH_CHIPS) * sizeof(float), volk_gnsssdr_get_alignment()));
in_cpu = static_cast<gr_complex*>(volk_gnsssdr_malloc(2 * d_vector_length * sizeof(gr_complex), volk_gnsssdr_get_alignment()));
// correlator outputs (scalar)
d_n_correlator_taps = 3; // Early, Prompt, and Late
d_correlator_outs = static_cast<gr_complex*>(volk_gnsssdr_malloc(d_n_correlator_taps*sizeof(gr_complex), volk_gnsssdr_get_alignment()));
for (int n = 0; n < d_n_correlator_taps; n++)
{
d_correlator_outs[n] = gr_complex(0,0);
}
d_local_code_shift_chips = static_cast<float*>(volk_gnsssdr_malloc(d_n_correlator_taps*sizeof(float), volk_gnsssdr_get_alignment()));
// Set TAPs delay values [chips]
float d_early_late_spc_chips = 0.5;
d_local_code_shift_chips[0] = - d_early_late_spc_chips;
d_local_code_shift_chips[1] = 0.0;
d_local_code_shift_chips[2] = d_early_late_spc_chips;
//--- Perform initializations ------------------------------
//local code resampler on GPU
// generate local reference (1 sample per chip)
gps_l1_ca_code_gen_float(d_ca_code, 1, 0);
// generate inut signal
std::random_device r;
std::default_random_engine e1(r());
std::uniform_real_distribution<float> uniform_dist(0, 1);
for (int n = 0; n < 2 * d_vector_length; n++)
{
in_cpu[n] = std::complex<float>(uniform_dist(e1), uniform_dist(e1));
}
for (int n = 0; n < max_threads; n++)
{
correlator_pool[n] = new cpu_multicorrelator_real_codes();
correlator_pool[n]->init(d_vector_length, d_n_correlator_taps);
correlator_pool[n]->set_input_output_vectors(d_correlator_outs, in_cpu);
correlator_pool[n]->set_local_code_and_taps(static_cast<int>(GPS_L1_CA_CODE_LENGTH_CHIPS), d_ca_code, d_local_code_shift_chips);
}
float d_rem_carrier_phase_rad = 0.0;
float d_carrier_phase_step_rad = 0.1;
float d_code_phase_step_chips = 0.3;
float d_rem_code_phase_chips = 0.4;
EXPECT_NO_THROW(
for(int correlation_sizes_idx = 0; correlation_sizes_idx < 3; correlation_sizes_idx++)
{
for(int current_max_threads = 1; current_max_threads < (max_threads+1); current_max_threads++)
{
std::cout << "Running " << current_max_threads << " concurrent correlators" << std::endl;
start = std::chrono::system_clock::now();
//create the concurrent correlator threads
for (int current_thread = 0; current_thread < current_max_threads; current_thread++)
{
thread_pool.push_back(std::thread(run_correlator_cpu,
correlator_pool[current_thread],
d_rem_carrier_phase_rad,
d_carrier_phase_step_rad,
d_code_phase_step_chips,
d_rem_code_phase_chips,
correlation_sizes[correlation_sizes_idx]));
}
//wait the threads to finish they work and destroy the thread objects
for(auto &t : thread_pool)
{
t.join();
}
thread_pool.clear();
end = std::chrono::system_clock::now();
elapsed_seconds = end - start;
execution_times[correlation_sizes_idx] = elapsed_seconds.count() / static_cast<double>(FLAGS_cpu_multicorrelator_real_codes_iterations_test);
std::cout << "CPU Multicorrelator (real codes) execution time for length=" << correlation_sizes[correlation_sizes_idx]
<< " : " << execution_times[correlation_sizes_idx] << " [s]" << std::endl;
}
}
);
volk_gnsssdr_free(d_local_code_shift_chips);
volk_gnsssdr_free(d_correlator_outs);
volk_gnsssdr_free(d_ca_code);
volk_gnsssdr_free(in_cpu);
for (int n = 0; n< max_threads; n++)
{
correlator_pool[n]->free();
}
}

View File

@@ -33,6 +33,9 @@
#include <complex> #include <complex>
#include <random> #include <random>
#include <thread> #include <thread>
#include <gtest/gtest.h>
#include <gflags/gflags.h>
#include <gnuradio/gr_complex.h>
#include <volk_gnsssdr/volk_gnsssdr.h> #include <volk_gnsssdr/volk_gnsssdr.h>
#include "cpu_multicorrelator.h" #include "cpu_multicorrelator.h"
#include "gps_sdr_signal_processing.h" #include "gps_sdr_signal_processing.h"