mirror of
https://github.com/gnss-sdr/gnss-sdr
synced 2025-04-07 03:07:02 +00:00
Adding GPS L2C M code generators and PCPS Acquisition algorithms
skeletons
This commit is contained in:
parent
e72a743e22
commit
75983a0471
372
src/algorithms/acquisition/adapters/gps_l2_m_pcps_acquisition.cc
Normal file
372
src/algorithms/acquisition/adapters/gps_l2_m_pcps_acquisition.cc
Normal file
@ -0,0 +1,372 @@
|
||||
/*!
|
||||
* \file gps_l1_ca_pcps_acquisition.cc
|
||||
* \brief Adapts a PCPS acquisition block to an AcquisitionInterface for
|
||||
* GPS L1 C/A signals
|
||||
* \authors <ul>
|
||||
* <li> Javier Arribas, 2011. jarribas(at)cttc.es
|
||||
* </ul>
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*
|
||||
* 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 <http://www.gnu.org/licenses/>.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*/
|
||||
|
||||
#include "gps_l2_m_pcps_acquisition.h"
|
||||
#include <iostream>
|
||||
#include <stdexcept>
|
||||
#include <boost/math/distributions/exponential.hpp>
|
||||
#include <glog/logging.h>
|
||||
#include <gnuradio/msg_queue.h>
|
||||
#include "gps_sdr_signal_processing.h"
|
||||
#include "GPS_L1_CA.h"
|
||||
#include "configuration_interface.h"
|
||||
|
||||
|
||||
using google::LogMessage;
|
||||
|
||||
GpsL2MPcpsAcquisition::GpsL2MPcpsAcquisition(
|
||||
ConfigurationInterface* configuration, std::string role,
|
||||
unsigned int in_streams, unsigned int out_streams,
|
||||
gr::msg_queue::sptr queue) :
|
||||
role_(role), in_streams_(in_streams), out_streams_(out_streams), queue_(queue)
|
||||
{
|
||||
configuration_ = configuration;
|
||||
std::string default_item_type = "gr_complex";
|
||||
std::string default_dump_filename = "./data/acquisition.dat";
|
||||
|
||||
DLOG(INFO) << "role " << role;
|
||||
|
||||
item_type_ = configuration_->property(role + ".item_type",
|
||||
default_item_type);
|
||||
//float pfa = configuration_->property(role + ".pfa", 0.0);
|
||||
|
||||
fs_in_ = configuration_->property("GNSS-SDR.internal_fs_hz", 2048000);
|
||||
if_ = configuration_->property(role + ".ifreq", 0);
|
||||
dump_ = configuration_->property(role + ".dump", false);
|
||||
shift_resolution_ = configuration_->property(role + ".doppler_max", 15);
|
||||
sampled_ms_ = configuration_->property(role + ".coherent_integration_time_ms", 1);
|
||||
|
||||
bit_transition_flag_ = configuration_->property(role + ".bit_transition_flag", false);
|
||||
|
||||
if (!bit_transition_flag_)
|
||||
{
|
||||
max_dwells_ = configuration_->property(role + ".max_dwells", 1);
|
||||
}
|
||||
else
|
||||
{
|
||||
max_dwells_ = 2;
|
||||
}
|
||||
|
||||
dump_filename_ = configuration_->property(role + ".dump_filename",
|
||||
default_dump_filename);
|
||||
|
||||
//--- Find number of samples per spreading code -------------------------
|
||||
code_length_ = round(fs_in_
|
||||
/ (GPS_L1_CA_CODE_RATE_HZ / GPS_L1_CA_CODE_LENGTH_CHIPS));
|
||||
|
||||
vector_length_ = code_length_ * sampled_ms_;
|
||||
|
||||
code_= new gr_complex[vector_length_];
|
||||
|
||||
// if (item_type_.compare("gr_complex") == 0 )
|
||||
// {
|
||||
item_size_ = sizeof(gr_complex);
|
||||
acquisition_cc_ = pcps_make_acquisition_cc(sampled_ms_, max_dwells_,
|
||||
shift_resolution_, if_, fs_in_, code_length_, code_length_,
|
||||
bit_transition_flag_, queue_, dump_, dump_filename_);
|
||||
|
||||
stream_to_vector_ = gr::blocks::stream_to_vector::make(item_size_, vector_length_);
|
||||
|
||||
DLOG(INFO) << "stream_to_vector(" << stream_to_vector_->unique_id() << ")";
|
||||
DLOG(INFO) << "acquisition(" << acquisition_cc_->unique_id() << ")";
|
||||
// }
|
||||
|
||||
if (item_type_.compare("cshort") == 0)
|
||||
{
|
||||
cshort_to_float_x2_ = make_cshort_to_float_x2();
|
||||
float_to_complex_ = gr::blocks::float_to_complex::make();
|
||||
}
|
||||
|
||||
if (item_type_.compare("cbyte") == 0)
|
||||
{
|
||||
cbyte_to_float_x2_ = make_complex_byte_to_float_x2();
|
||||
float_to_complex_ = gr::blocks::float_to_complex::make();
|
||||
}
|
||||
//}
|
||||
//else
|
||||
// {
|
||||
// LOG(WARNING) << item_type_
|
||||
// << " unknown acquisition item type";
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
GpsL2MPcpsAcquisition::~GpsL2MPcpsAcquisition()
|
||||
{
|
||||
delete[] code_;
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_channel(unsigned int channel)
|
||||
{
|
||||
channel_ = channel;
|
||||
//if (item_type_.compare("gr_complex") == 0)
|
||||
//{
|
||||
acquisition_cc_->set_channel(channel_);
|
||||
//}
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_threshold(float threshold)
|
||||
{
|
||||
float pfa = configuration_->property(role_ + boost::lexical_cast<std::string>(channel_) + ".pfa", 0.0);
|
||||
|
||||
if(pfa == 0.0)
|
||||
{
|
||||
pfa = configuration_->property(role_+".pfa", 0.0);
|
||||
}
|
||||
if(pfa == 0.0)
|
||||
{
|
||||
threshold_ = threshold;
|
||||
}
|
||||
else
|
||||
{
|
||||
threshold_ = calculate_threshold(pfa);
|
||||
}
|
||||
|
||||
DLOG(INFO) <<"Channel "<<channel_<<" Threshold = " << threshold_;
|
||||
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_threshold(threshold_);
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_doppler_max(unsigned int doppler_max)
|
||||
{
|
||||
doppler_max_ = doppler_max;
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_doppler_max(doppler_max_);
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_doppler_step(unsigned int doppler_step)
|
||||
{
|
||||
doppler_step_ = doppler_step;
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_doppler_step(doppler_step_);
|
||||
// }
|
||||
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_channel_queue(
|
||||
concurrent_queue<int> *channel_internal_queue)
|
||||
{
|
||||
channel_internal_queue_ = channel_internal_queue;
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_channel_queue(channel_internal_queue_);
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_gnss_synchro(Gnss_Synchro* gnss_synchro)
|
||||
{
|
||||
gnss_synchro_ = gnss_synchro;
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_gnss_synchro(gnss_synchro_);
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
signed int GpsL2MPcpsAcquisition::mag()
|
||||
{
|
||||
// // if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
return acquisition_cc_->mag();
|
||||
// }
|
||||
// else
|
||||
// {
|
||||
// return 0;
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::init()
|
||||
{
|
||||
acquisition_cc_->init();
|
||||
set_local_code();
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_local_code()
|
||||
{
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
std::complex<float>* code = new std::complex<float>[code_length_];
|
||||
|
||||
gps_l1_ca_code_gen_complex_sampled(code, gnss_synchro_->PRN, fs_in_, 0);
|
||||
|
||||
for (unsigned int i = 0; i < sampled_ms_; i++)
|
||||
{
|
||||
memcpy(&(code_[i*code_length_]), code,
|
||||
sizeof(gr_complex)*code_length_);
|
||||
}
|
||||
|
||||
acquisition_cc_->set_local_code(code_);
|
||||
|
||||
delete[] code;
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::reset()
|
||||
{
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_active(true);
|
||||
// }
|
||||
}
|
||||
|
||||
void GpsL2MPcpsAcquisition::set_state(int state)
|
||||
{
|
||||
// if (item_type_.compare("gr_complex") == 0)
|
||||
// {
|
||||
acquisition_cc_->set_state(state);
|
||||
// }
|
||||
}
|
||||
|
||||
|
||||
|
||||
float GpsL2MPcpsAcquisition::calculate_threshold(float pfa)
|
||||
{
|
||||
//Calculate the threshold
|
||||
unsigned int frequency_bins = 0;
|
||||
for (int doppler = (int)(-doppler_max_); doppler <= (int)doppler_max_; doppler += doppler_step_)
|
||||
{
|
||||
frequency_bins++;
|
||||
}
|
||||
DLOG(INFO) << "Channel " << channel_<< " Pfa = " << pfa;
|
||||
unsigned int ncells = vector_length_ * frequency_bins;
|
||||
double exponent = 1 / static_cast<double>(ncells);
|
||||
double val = pow(1.0 - pfa, exponent);
|
||||
double lambda = double(vector_length_);
|
||||
boost::math::exponential_distribution<double> mydist (lambda);
|
||||
float threshold = (float)quantile(mydist,val);
|
||||
|
||||
return threshold;
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::connect(gr::top_block_sptr top_block)
|
||||
{
|
||||
if (item_type_.compare("gr_complex") == 0)
|
||||
{
|
||||
top_block->connect(stream_to_vector_, 0, acquisition_cc_, 0);
|
||||
}
|
||||
else if (item_type_.compare("cshort") == 0)
|
||||
{
|
||||
top_block->connect(cshort_to_float_x2_, 0, float_to_complex_, 0);
|
||||
top_block->connect(cshort_to_float_x2_, 1, float_to_complex_, 1);
|
||||
top_block->connect(float_to_complex_, 0, stream_to_vector_, 0);
|
||||
top_block->connect(stream_to_vector_, 0, acquisition_cc_, 0);
|
||||
}
|
||||
else if (item_type_.compare("cbyte") == 0)
|
||||
{
|
||||
top_block->connect(cbyte_to_float_x2_, 0, float_to_complex_, 0);
|
||||
top_block->connect(cbyte_to_float_x2_, 1, float_to_complex_, 1);
|
||||
top_block->connect(float_to_complex_, 0, stream_to_vector_, 0);
|
||||
top_block->connect(stream_to_vector_, 0, acquisition_cc_, 0);
|
||||
}
|
||||
else
|
||||
{
|
||||
LOG(WARNING) << item_type_ << " unknown acquisition item type";
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
|
||||
void GpsL2MPcpsAcquisition::disconnect(gr::top_block_sptr top_block)
|
||||
{
|
||||
if (item_type_.compare("gr_complex") == 0)
|
||||
{
|
||||
top_block->disconnect(stream_to_vector_, 0, acquisition_cc_, 0);
|
||||
}
|
||||
else if (item_type_.compare("cshort") == 0)
|
||||
{
|
||||
// Since a short-based acq implementation is not available,
|
||||
// we just convert cshorts to gr_complex
|
||||
top_block->disconnect(cshort_to_float_x2_, 0, float_to_complex_, 0);
|
||||
top_block->disconnect(cshort_to_float_x2_, 1, float_to_complex_, 1);
|
||||
top_block->disconnect(float_to_complex_, 0, stream_to_vector_, 0);
|
||||
top_block->disconnect(stream_to_vector_, 0, acquisition_cc_, 0);
|
||||
}
|
||||
else if (item_type_.compare("cbyte") == 0)
|
||||
{
|
||||
// Since a byte-based acq implementation is not available,
|
||||
// we just convert cshorts to gr_complex
|
||||
top_block->disconnect(cbyte_to_float_x2_, 0, float_to_complex_, 0);
|
||||
top_block->disconnect(cbyte_to_float_x2_, 1, float_to_complex_, 1);
|
||||
top_block->disconnect(float_to_complex_, 0, stream_to_vector_, 0);
|
||||
top_block->disconnect(stream_to_vector_, 0, acquisition_cc_, 0);
|
||||
}
|
||||
else
|
||||
{
|
||||
LOG(WARNING) << item_type_ << " unknown acquisition item type";
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
gr::basic_block_sptr GpsL2MPcpsAcquisition::get_left_block()
|
||||
{
|
||||
if (item_type_.compare("gr_complex") == 0)
|
||||
{
|
||||
return stream_to_vector_;
|
||||
}
|
||||
else if (item_type_.compare("cshort") == 0)
|
||||
{
|
||||
return cshort_to_float_x2_;
|
||||
}
|
||||
else if (item_type_.compare("cbyte") == 0)
|
||||
{
|
||||
return cbyte_to_float_x2_;
|
||||
}
|
||||
else
|
||||
{
|
||||
LOG(WARNING) << item_type_ << " unknown acquisition item type";
|
||||
return nullptr;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
gr::basic_block_sptr GpsL2MPcpsAcquisition::get_right_block()
|
||||
{
|
||||
return acquisition_cc_;
|
||||
}
|
||||
|
177
src/algorithms/acquisition/adapters/gps_l2_m_pcps_acquisition.h
Normal file
177
src/algorithms/acquisition/adapters/gps_l2_m_pcps_acquisition.h
Normal file
@ -0,0 +1,177 @@
|
||||
/*!
|
||||
* \file gps_l2_m_pcps_acquisition.h
|
||||
* \brief Adapts a PCPS acquisition block to an AcquisitionInterface for
|
||||
* GPS L2 M signals
|
||||
* \authors <ul>
|
||||
* <li> Javier Arribas, 2015. jarribas(at)cttc.es
|
||||
* </ul>
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*
|
||||
* 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 <http://www.gnu.org/licenses/>.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*/
|
||||
|
||||
#ifndef GNSS_SDR_GPS_L2_M_PCPS_ACQUISITION_H_
|
||||
#define GNSS_SDR_GPS_L2_M_PCPS_ACQUISITION_H_
|
||||
|
||||
#include <string>
|
||||
#include <gnuradio/msg_queue.h>
|
||||
#include <gnuradio/blocks/stream_to_vector.h>
|
||||
#include <gnuradio/blocks/float_to_complex.h>
|
||||
#include "gnss_synchro.h"
|
||||
#include "acquisition_interface.h"
|
||||
#include "pcps_acquisition_cc.h"
|
||||
#include "cshort_to_float_x2.h"
|
||||
#include "complex_byte_to_float_x2.h"
|
||||
|
||||
|
||||
|
||||
class ConfigurationInterface;
|
||||
|
||||
/*!
|
||||
* \brief This class adapts a PCPS acquisition block to an AcquisitionInterface
|
||||
* for GPS L2 M signals
|
||||
*/
|
||||
class GpsL2MPcpsAcquisition: public AcquisitionInterface
|
||||
{
|
||||
public:
|
||||
GpsL2MPcpsAcquisition(ConfigurationInterface* configuration,
|
||||
std::string role, unsigned int in_streams,
|
||||
unsigned int out_streams, boost::shared_ptr<gr::msg_queue> queue);
|
||||
|
||||
virtual ~GpsL2MPcpsAcquisition();
|
||||
|
||||
std::string role()
|
||||
{
|
||||
return role_;
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Returns "GPS_L2_M_PCPS_Acquisition"
|
||||
*/
|
||||
std::string implementation()
|
||||
{
|
||||
return "GPS_L2_M_PCPS_Acquisition";
|
||||
}
|
||||
size_t item_size()
|
||||
{
|
||||
return item_size_;
|
||||
}
|
||||
|
||||
void connect(gr::top_block_sptr top_block);
|
||||
void disconnect(gr::top_block_sptr top_block);
|
||||
gr::basic_block_sptr get_left_block();
|
||||
gr::basic_block_sptr get_right_block();
|
||||
|
||||
/*!
|
||||
* \brief Set acquisition/tracking common Gnss_Synchro object pointer
|
||||
* to efficiently exchange synchronization data between acquisition and
|
||||
* tracking blocks
|
||||
*/
|
||||
void set_gnss_synchro(Gnss_Synchro* p_gnss_synchro);
|
||||
|
||||
/*!
|
||||
* \brief Set acquisition channel unique ID
|
||||
*/
|
||||
void set_channel(unsigned int channel);
|
||||
|
||||
/*!
|
||||
* \brief Set statistics threshold of PCPS algorithm
|
||||
*/
|
||||
void set_threshold(float threshold);
|
||||
|
||||
/*!
|
||||
* \brief Set maximum Doppler off grid search
|
||||
*/
|
||||
void set_doppler_max(unsigned int doppler_max);
|
||||
|
||||
/*!
|
||||
* \brief Set Doppler steps for the grid search
|
||||
*/
|
||||
void set_doppler_step(unsigned int doppler_step);
|
||||
|
||||
/*!
|
||||
* \brief Set tracking channel internal queue
|
||||
*/
|
||||
void set_channel_queue(concurrent_queue<int> *channel_internal_queue);
|
||||
|
||||
/*!
|
||||
* \brief Initializes acquisition algorithm.
|
||||
*/
|
||||
void init();
|
||||
|
||||
/*!
|
||||
* \brief Sets local code for GPS L2/M PCPS acquisition algorithm.
|
||||
*/
|
||||
void set_local_code();
|
||||
|
||||
/*!
|
||||
* \brief Returns the maximum peak of grid search
|
||||
*/
|
||||
signed int mag();
|
||||
|
||||
/*!
|
||||
* \brief Restart acquisition algorithm
|
||||
*/
|
||||
void reset();
|
||||
|
||||
/*!
|
||||
* \brief If state = 1, it forces the block to start acquiring from the first sample
|
||||
*/
|
||||
void set_state(int state);
|
||||
|
||||
private:
|
||||
ConfigurationInterface* configuration_;
|
||||
pcps_acquisition_cc_sptr acquisition_cc_;
|
||||
gr::blocks::stream_to_vector::sptr stream_to_vector_;
|
||||
gr::blocks::float_to_complex::sptr float_to_complex_;
|
||||
cshort_to_float_x2_sptr cshort_to_float_x2_;
|
||||
complex_byte_to_float_x2_sptr cbyte_to_float_x2_;
|
||||
size_t item_size_;
|
||||
std::string item_type_;
|
||||
unsigned int vector_length_;
|
||||
unsigned int code_length_;
|
||||
bool bit_transition_flag_;
|
||||
unsigned int channel_;
|
||||
float threshold_;
|
||||
unsigned int doppler_max_;
|
||||
unsigned int doppler_step_;
|
||||
unsigned int shift_resolution_;
|
||||
unsigned int sampled_ms_;
|
||||
unsigned int max_dwells_;
|
||||
long fs_in_;
|
||||
long if_;
|
||||
bool dump_;
|
||||
std::string dump_filename_;
|
||||
std::complex<float> * code_;
|
||||
Gnss_Synchro * gnss_synchro_;
|
||||
std::string role_;
|
||||
unsigned int in_streams_;
|
||||
unsigned int out_streams_;
|
||||
boost::shared_ptr<gr::msg_queue> queue_;
|
||||
concurrent_queue<int> *channel_internal_queue_;
|
||||
|
||||
float calculate_threshold(float pfa);
|
||||
};
|
||||
|
||||
#endif /* GNSS_SDR_GPS_L2_M_PCPS_ACQUISITION_H_ */
|
107
src/algorithms/libs/gps_l2c_signal.cc
Normal file
107
src/algorithms/libs/gps_l2c_signal.cc
Normal file
@ -0,0 +1,107 @@
|
||||
/*!
|
||||
* \file gps_l2c_signal.cc
|
||||
* \brief This class implements signal generators for the GPS L2C signals
|
||||
* \author Javier Arribas, 2015. jarribas(at)cttc.es
|
||||
*
|
||||
* Detailed description of the file here if needed.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*
|
||||
* 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 <http://www.gnu.org/licenses/>.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*/
|
||||
|
||||
#include "GPS_L2C.h"
|
||||
#include <stdlib.h>
|
||||
#include <stdint.h>
|
||||
#include <cmath>
|
||||
|
||||
|
||||
int8_t gps_l2c_m_shift(int x)
|
||||
{
|
||||
int8_t a;
|
||||
return (int8_t)((x>>1)^(x&1)*0445112474);
|
||||
}
|
||||
|
||||
void gps_l2c_m_code(int8_t * _dest, int _prn)
|
||||
{
|
||||
int x;
|
||||
x= GPS_L2C_M_INIT_REG[_prn];
|
||||
for (int n=0; n<GPS_L2_M_CODE_LENGTH_CHIPS; n++)
|
||||
{
|
||||
x= gps_l2c_m_shift(x);
|
||||
_dest[n]=(int8_t)(x&1);
|
||||
}
|
||||
}
|
||||
|
||||
/*
|
||||
* Generates complex GPS L1 C/A code for the desired SV ID and sampled to specific sampling frequency
|
||||
*/
|
||||
void gps_l1_ca_code_gen_complex_sampled(std::complex<float>* _dest, unsigned int _prn, signed int _fs)
|
||||
{
|
||||
int8_t _code[GPS_L2_M_CODE_LENGTH_CHIPS];
|
||||
|
||||
|
||||
signed int _samplesPerCode, _codeValueIndex;
|
||||
float _ts;
|
||||
float _tc;
|
||||
const signed int _codeFreqBasis = GPS_L2_M_CODE_RATE_HZ; //Hz
|
||||
const signed int _codeLength = GPS_L2_M_CODE_LENGTH_CHIPS;
|
||||
|
||||
//--- Find number of samples per spreading code ----------------------------
|
||||
_samplesPerCode = round(_fs / (_codeFreqBasis / _codeLength));
|
||||
|
||||
//--- Find time constants --------------------------------------------------
|
||||
_ts = 1/(float)_fs; // Sampling period in sec
|
||||
_tc = 1/(float)_codeFreqBasis; // C/A chip period in sec
|
||||
//gps_l1_ca_code_gen_complex(_code,_prn); //generate C/A code 1 sample per chip
|
||||
|
||||
for (signed int i=0; i<_samplesPerCode; i++)
|
||||
{
|
||||
//=== Digitizing =======================================================
|
||||
|
||||
//--- Make index array to read C/A code values -------------------------
|
||||
// The length of the index array depends on the sampling frequency -
|
||||
// number of samples per millisecond (because one C/A code period is one
|
||||
// millisecond).
|
||||
|
||||
_codeValueIndex = ceil((_ts * ((float)i + 1)) / _tc) - 1;
|
||||
|
||||
//--- Make the digitized version of the C/A code -----------------------
|
||||
// The "upsampled" code is made by selecting values form the CA code
|
||||
// chip array (caCode) for the time instances of each sample.
|
||||
if (i == _samplesPerCode - 1)
|
||||
{
|
||||
//--- Correct the last index (due to number rounding issues) -----------
|
||||
_dest[i] = std::complex(1.0-2.0*_code[_codeLength - 1],0);
|
||||
|
||||
}
|
||||
else
|
||||
{
|
||||
_dest[i] = std::complex(1.0-2.0*_code[_codeValueIndex],0);; //repeat the chip -> upsample
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
|
49
src/algorithms/libs/gps_l2c_signal.h
Normal file
49
src/algorithms/libs/gps_l2c_signal.h
Normal file
@ -0,0 +1,49 @@
|
||||
/*!
|
||||
* \file gps_l2c_signal.h
|
||||
* \brief This class implements signal generators for the GPS L2C signals
|
||||
* \author Javier Arribas, 2015. jarribas(at)cttc.es
|
||||
*
|
||||
* Detailed description of the file here if needed.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*
|
||||
* 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 <http://www.gnu.org/licenses/>.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*/
|
||||
|
||||
#ifndef GNSS_GPS_L2C_SIGNAL_H_
|
||||
#define GNSS_GPS_L2C_SIGNAL_H_
|
||||
|
||||
#include <complex>
|
||||
#include <iostream>
|
||||
#include "GPS_L2C.h"
|
||||
|
||||
//!Generates complex GPS L2C M code for the desired SV ID and code shift, and sampled to specific sampling frequency
|
||||
void gps_l2_m_code_gen_complex(std::complex<float>* _dest, signed int _prn, unsigned int _chip_shift);
|
||||
|
||||
//! Generates N complex GPS L2C M codes for the desired SV ID and code shift
|
||||
void gps_l2_m_code_gen_complex_sampled(std::complex<float>* _dest, unsigned int _prn, signed int _fs, unsigned int _chip_shift, unsigned int _ncodes);
|
||||
|
||||
//! Generates complex GPS L2C M code for the desired SV ID and code shift
|
||||
void gps_l2_m_code_gen_complex_sampled(std::complex<float>* _dest, unsigned int _prn, signed int _fs, unsigned int _chip_shift);
|
||||
|
||||
#endif /* GNSS_GPS_L2C_SIGNAL_H_ */
|
85
src/core/system_parameters/GPS_L2C.h
Normal file
85
src/core/system_parameters/GPS_L2C.h
Normal file
@ -0,0 +1,85 @@
|
||||
/*!
|
||||
* \file GPS_L2C.h
|
||||
* \brief Defines system parameters for GPS L1 C/A signal and NAV data
|
||||
* \author Javier Arribas, 2015. jarribas(at)cttc.es
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*
|
||||
* 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 <http://www.gnu.org/licenses/>.
|
||||
*
|
||||
* -------------------------------------------------------------------------
|
||||
*/
|
||||
|
||||
|
||||
#ifndef GNSS_SDR_GPS_L2C_H_
|
||||
#define GNSS_SDR_GPS_L2C_H_
|
||||
|
||||
#include <complex>
|
||||
#include <vector>
|
||||
#include <utility> // std::pair
|
||||
#include <gnss_satellite.h>
|
||||
#include "MATH_CONSTANTS.h"
|
||||
|
||||
|
||||
// carrier and code frequencies
|
||||
const double GPS_L2_FREQ_HZ = 1.2276e9; //!< L2 [Hz]
|
||||
|
||||
const double GPS_L2_M_CODE_RATE_HZ = 0.5115e6; //!< GPS L2 M code rate [chips/s]
|
||||
const int GPS_L2_M_CODE_LENGTH_CHIPS = 10230; //!< GPS L2 M code length [chips]
|
||||
const double GPS_L2_M_PERIOD = 0.02; //!< GPS L2 M code period [seconds]
|
||||
|
||||
const double GPS_L2_L_CODE_RATE_HZ = 0.5115e6; //!< GPS L2 L code rate [chips/s]
|
||||
const int GPS_L2_L_CODE_LENGTH_CHIPS = 767250; //!< GPS L2 L code length [chips]
|
||||
const double GPS_L2_L_PERIOD = 1.5; //!< GPS L2 L code period [seconds]
|
||||
|
||||
|
||||
const int GPS_L2C_M_INIT_REG[115] =
|
||||
{0742417664, 0756014035,0002747144,0066265724, // 1:4
|
||||
0601403471, 0703232733, 0124510070,0617316361, // 5:8
|
||||
0047541621, 0733031046, 0713512145, 0024437606,
|
||||
0021264003, 0230655351, 0001314400, 0222021506,
|
||||
0540264026, 0205521705, 0064022144, 0120161274,
|
||||
0044023533, 0724744327, 0045743577, 0741201660,
|
||||
0700274134, 0010247261, 0713433445, 0737324162,
|
||||
0311627434, 0710452007, 0722462133, 0050172213,
|
||||
0500653703, 0755077436, 0136717361, 0756675453,
|
||||
0435506112, 0771353753, 0226107701, 0022025110,
|
||||
0402466344, 0752566114, 0702011164, 0041216771,
|
||||
0047457275, 0266333164, 0713167356, 0060546335,
|
||||
0355173035, 0617201036, 0157465571, 0767360553,
|
||||
0023127030, 0431343777, 0747317317, 0045706125,
|
||||
0002744276, 0060036467, 0217744147, 0603340174,//57:60
|
||||
0326616775, 0063240065, 0111460621, //61:63
|
||||
0604055104, 0157065232, 0013305707, 0603552017,//159:162
|
||||
0230461355, 0603653437, 0652346475, 0743107103,
|
||||
0401521277, 0167335110, 0014013575, 0362051132,
|
||||
0617753265, 0216363634, 0755561123, 0365304033,
|
||||
0625025543, 0054420334, 0415473671, 0662364360,
|
||||
0373446602, 0417564100, 0000526452, 0226631300,
|
||||
0113752074, 0706134401, 0041352546, 0664630154,
|
||||
0276524255, 0714720530, 0714051771, 0044526647,
|
||||
0207164322, 0262120161, 0204244652, 0202133131,
|
||||
0714351204, 0657127260, 0130567507, 0670517677,
|
||||
0607275514, 0045413633, 0212645405, 0613700455,
|
||||
0706202440, 0705056276, 0020373522, 0746013617,
|
||||
0132720621, 0434015513, 0566721727, 0140633660};
|
||||
|
||||
#endif /* GNSS_SDR_GPS_L2C_H_ */
|
Loading…
x
Reference in New Issue
Block a user