fixing coverity issues

This commit is contained in:
Carles Fernandez 2015-05-13 23:50:21 +02:00
parent cf9945899f
commit 79192a0bbd
14 changed files with 76 additions and 60 deletions

View File

@ -107,7 +107,7 @@ arma::vec galileo_e1_ls_pvt::rotateSatellite(double traveltime, arma::vec X_sat)
R3(1, 2) = 0.0;
R3(2, 0) = 0.0;
R3(2, 1) = 0.0;
R3(2, 2) = 1;
R3(2, 2) = 1.0;
//--- Do the rotation ------------------------------------------------------
arma::vec X_sat_rot;
@ -188,7 +188,7 @@ arma::vec galileo_e1_ls_pvt::leastSquarePos(arma::mat satpos, arma::vec obs, arm
togeod(&dphi, &dlambda, &h, 6378137.0, 298.257223563, pos(0), pos(1), pos(2));
//--- Find delay due to troposphere (in meters)
tropo(&trop, sin(d_visible_satellites_El[i] * GALILEO_PI/180.0), h/1000, 1013.0, 293.0, 50.0, 0.0, 0.0, 0.0);
tropo(&trop, sin(d_visible_satellites_El[i] * GALILEO_PI / 180.0), h / 1000.0, 1013.0, 293.0, 50.0, 0.0, 0.0, 0.0);
if(trop > 50.0 ) trop = 0.0;
}
}
@ -238,9 +238,9 @@ bool galileo_e1_ls_pvt::get_PVT(std::map<int,Gnss_Synchro> gnss_pseudoranges_map
arma::mat satpos = arma::zeros(3, valid_pseudoranges); // satellite positions matrix
int Galileo_week_number = 0;
double utc = 0;
double TX_time_corrected_s;
double SV_clock_bias_s = 0;
double utc = 0.0;
double TX_time_corrected_s = 0.0;
double SV_clock_bias_s = 0.0;
d_flag_averaging = flag_averaging;
@ -260,7 +260,7 @@ bool galileo_e1_ls_pvt::get_PVT(std::map<int,Gnss_Synchro> gnss_pseudoranges_map
/*!
* \todo Place here the satellite CN0 (power level, or weight factor)
*/
W(obs_counter, obs_counter) = 1;
W(obs_counter, obs_counter) = 1.0;
// COMMON RX TIME PVT ALGORITHM MODIFICATION (Like RINEX files)
// first estimate of transmit time
@ -519,14 +519,14 @@ void galileo_e1_ls_pvt::cart2geo(double X, double Y, double Z, int elipsoid_sele
double phi = atan(Z / ((sqrt(X * X + Y * Y) * (1.0 - (2.0 - f[elipsoid_selection])) * f[elipsoid_selection])));
double h = 0.1;
double oldh = 0;
double oldh = 0.0;
double N;
int iterations = 0;
do
{
oldh = h;
N = c / sqrt(1 + ex2 * (cos(phi) * cos(phi)));
phi = atan(Z / ((sqrt(X * X + Y * Y) * (1 - (2 - f[elipsoid_selection]) * f[elipsoid_selection] * N / (N + h) ))));
phi = atan(Z / ((sqrt(X * X + Y * Y) * (1.0 - (2.0 - f[elipsoid_selection]) * f[elipsoid_selection] * N / (N + h) ))));
h = sqrt(X * X + Y * Y) / cos(phi) - N;
iterations = iterations + 1;
if (iterations > 100)
@ -568,17 +568,17 @@ void galileo_e1_ls_pvt::togeod(double *dphi, double *dlambda, double *h, double
*h = 0;
double tolsq = 1.e-10; // tolerance to accept convergence
int maxit = 10; // max number of iterations
double rtd = 180/GPS_PI;
double rtd = 180.0 / GPS_PI;
// compute square of eccentricity
double esq;
if (finv < 1.0E-20)
{
esq = 0;
esq = 0.0;
}
else
{
esq = (2 - 1/finv) / finv;
esq = (2.0 - 1.0 / finv) / finv;
}
// first guess
@ -591,7 +591,7 @@ void galileo_e1_ls_pvt::togeod(double *dphi, double *dlambda, double *h, double
}
else
{
*dlambda = 0;
*dlambda = 0.0;
}
// correct longitude bound
@ -609,7 +609,7 @@ void galileo_e1_ls_pvt::togeod(double *dphi, double *dlambda, double *h, double
}
else
{
sinphi = 0;
sinphi = 0.0;
}
*dphi = asin(sinphi);
@ -628,7 +628,7 @@ void galileo_e1_ls_pvt::togeod(double *dphi, double *dlambda, double *h, double
double N_phi;
double dP;
double dZ;
double oneesq = 1 - esq;
double oneesq = 1.0 - esq;
for (int i = 0; i < maxit; i++)
{
@ -700,7 +700,7 @@ void galileo_e1_ls_pvt::topocent(double *Az, double *El, double *D, arma::vec x,
F(1,1) = -sb * sl;
F(1,2) = cb * sl;
F(2,0) = 0;
F(2,0) = 0.0;
F(2,1) = cb;
F(2,2) = sb;
@ -717,8 +717,8 @@ void galileo_e1_ls_pvt::topocent(double *Az, double *El, double *D, arma::vec x,
if (hor_dis < 1.0E-20)
{
*Az = 0;
*El = 90;
*Az = 0.0;
*El = 90.0;
}
else
{

View File

@ -175,7 +175,8 @@ gr::basic_block_sptr Channel::get_right_block()
void Channel::set_signal(Gnss_Signal gnss_signal)
{
gnss_signal_ = gnss_signal;
const char * str = gnss_signal_.get_signal_str().c_str(); // get a C style null terminated string
std::string str_aux = gnss_signal_.get_signal_str();
const char * str = str_aux.c_str(); // get a C style null terminated string
std::memcpy((void*)gnss_synchro_.Signal, str, 3); // copy string into synchro char array: 2 char + null
gnss_synchro_.Signal[2] = 0; // make sure that string length is only two characters
gnss_synchro_.PRN = gnss_signal_.get_satellite().get_PRN();

View File

@ -115,6 +115,9 @@ public:
GpsL1CaChannelFsm::GpsL1CaChannelFsm()
{
acq_ = nullptr;
trk_ = nullptr;
channel_ = 0;
initiate(); //start the FSM
}
@ -123,6 +126,8 @@ GpsL1CaChannelFsm::GpsL1CaChannelFsm()
GpsL1CaChannelFsm::GpsL1CaChannelFsm(AcquisitionInterface *acquisition) :
acq_(acquisition)
{
trk_ = nullptr;
channel_ = 0;
initiate(); //start the FSM
}

View File

@ -85,7 +85,6 @@ public:
private:
AcquisitionInterface *acq_;
TrackingInterface *trk_;
TelemetryDecoderInterface *nav_;
boost::shared_ptr<gr::msg_queue> queue_;
unsigned int channel_;
};

View File

@ -101,11 +101,10 @@ void galileo_e5_a_code_gen_complex_sampled(std::complex<float>* _dest, char _Sig
// This function is based on the GNU software GPS for MATLAB in the Kay Borre book
unsigned int _samplesPerCode;
unsigned int delay;
unsigned int _codeLength = Galileo_E5a_CODE_LENGTH_CHIPS;
const unsigned int _codeLength = Galileo_E5a_CODE_LENGTH_CHIPS;
const int _codeFreqBasis = Galileo_E5a_CODE_CHIP_RATE_HZ; //Hz
std::complex<float>* _code;
_code = new std::complex<float>[_codeLength];
std::complex<float>* _code = new std::complex<float>[_codeLength];
galileo_e5_a_code_gen_complex_primary(_code , _prn , _Signal);
@ -118,8 +117,8 @@ void galileo_e5_a_code_gen_complex_sampled(std::complex<float>* _dest, char _Sig
std::complex<float>* _resampled_signal;
if (posix_memalign((void**)&_resampled_signal, 16, _samplesPerCode * sizeof(gr_complex)) == 0){};
resampler(_code, _resampled_signal, _codeFreqBasis, _fs, _codeLength, _samplesPerCode); //resamples code to fs
delete[] _code;
_code = _resampled_signal;
//delete[] _code;
memcpy (_code, _resampled_signal, _codeLength);
}
for (unsigned int i = 0; i < _samplesPerCode; i++)
@ -127,5 +126,5 @@ void galileo_e5_a_code_gen_complex_sampled(std::complex<float>* _dest, char _Sig
_dest[(i + delay) % _samplesPerCode] = _code[i];
}
free(_code);
delete[] _code;
}

View File

@ -38,8 +38,6 @@ 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)
{
const unsigned int _code_length = 1023;
bool G1[_code_length];
bool G2[_code_length];
@ -135,14 +133,14 @@ void gps_l1_ca_code_gen_complex_sampled(std::complex<float>* _dest, unsigned int
const signed int _codeLength = 1023;
//--- Find number of samples per spreading code ----------------------------
_samplesPerCode = round(_fs / (_codeFreqBasis / _codeLength));
_samplesPerCode = static_cast<signed int>(static_cast<double>(_fs) / static_cast<double>(_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, _chip_shift); //generate C/A code 1 sample per chip
_ts = 1.0 / static_cast<float>(_fs); // Sampling period in sec
_tc = 1.0 / static_cast<float>(_codeFreqBasis); // C/A chip period in sec
gps_l1_ca_code_gen_complex(_code, _prn, _chip_shift); //generate C/A code 1 sample per chip
for (signed int i=0; i<_samplesPerCode; i++)
for (signed int i = 0; i < _samplesPerCode; i++)
{
//=== Digitizing =======================================================

View File

@ -75,7 +75,6 @@ private:
unsigned int out_stream_;
std::string item_type_;
size_t item_size_;
unsigned long long samples_;
bool dump_;
std::string dump_filename_;
double sample_freq_in_;

View File

@ -76,6 +76,7 @@ GalileoE1BTelemetryDecoder::GalileoE1BTelemetryDecoder(ConfigurationInterface* c
telemetry_decoder_->set_iono_queue(&global_galileo_iono_queue);
telemetry_decoder_->set_almanac_queue(&global_galileo_almanac_queue);
telemetry_decoder_->set_utc_model_queue(&global_galileo_utc_model_queue);
channel_ = 0;
}

View File

@ -82,6 +82,7 @@ GalileoE5aTelemetryDecoder::GalileoE5aTelemetryDecoder(ConfigurationInterface* c
telemetry_decoder_->set_utc_model_queue(&global_galileo_utc_model_queue);
DLOG(INFO) << "global navigation message queue assigned to telemetry_decoder ("<< telemetry_decoder_->unique_id() << ")";
channel_ = 0;
}

View File

@ -80,6 +80,7 @@ GpsL1CaTelemetryDecoder::GpsL1CaTelemetryDecoder(ConfigurationInterface* configu
int decimation_factor = configuration->property(role + ".decimation_factor", 1);
telemetry_decoder_->set_decimation(decimation_factor);
DLOG(INFO) << "global navigation message queue assigned to telemetry_decoder ("<< telemetry_decoder_->unique_id() << ")";
channel_ = 0;
}

View File

@ -49,20 +49,20 @@ class Tracking_2nd_DLL_filter
{
private:
// PLL filter parameters
float d_tau1_code;
float d_tau2_code;
float d_pdi_code;
float d_dllnoisebandwidth;
float d_dlldampingratio;
float d_old_code_error;
float d_old_code_nco;
float d_tau1_code = 0;
float d_tau2_code = 0;
float d_pdi_code = 0;
float d_dllnoisebandwidth = 0;
float d_dlldampingratio = 0;
float d_old_code_error = 0;
float d_old_code_nco = 0;
void calculate_lopp_coef(float* tau1,float* tau2, float lbw, float zeta, float k);
public:
void set_DLL_BW(float dll_bw_hz); //! Set DLL filter bandwidth [Hz]
void set_pdi(float pdi_code); //! Set Summation interval for code [s]
void initialize(); //! Start tracking with acquisition information
float get_code_nco(float DLL_discriminator); //! Numerically controlled oscillator
void set_DLL_BW(float dll_bw_hz); //! Set DLL filter bandwidth [Hz]
void set_pdi(float pdi_code); //! Set Summation interval for code [s]
void initialize(); //! Start tracking with acquisition information
float get_code_nco(float DLL_discriminator); //! Numerically controlled oscillator
Tracking_2nd_DLL_filter(float pdi_code);
Tracking_2nd_DLL_filter();
~Tracking_2nd_DLL_filter();

View File

@ -48,26 +48,26 @@ class Tracking_2nd_PLL_filter
{
private:
// PLL filter parameters
float d_tau1_carr;
float d_tau2_carr;
float d_pdi_carr;
float d_tau1_carr = 0;
float d_tau2_carr = 0;
float d_pdi_carr = 0;
float d_pllnoisebandwidth;
float d_plldampingratio;
float d_pllnoisebandwidth = 0;
float d_plldampingratio = 0;
float d_old_carr_error;
float d_old_carr_nco;
float d_old_carr_error = 0;
float d_old_carr_nco = 0;
void calculate_lopp_coef(float* tau1,float* tau2, float lbw, float zeta, float k);
public:
void set_PLL_BW(float pll_bw_hz); //! Set PLL loop bandwidth [Hz]
void set_pdi(float pdi_carr); //! Set Summation interval for code [s]
void initialize();
float get_carrier_nco(float PLL_discriminator);
Tracking_2nd_PLL_filter(float pdi_carr);
Tracking_2nd_PLL_filter();
~Tracking_2nd_PLL_filter();
void set_PLL_BW(float pll_bw_hz); //! Set PLL loop bandwidth [Hz]
void set_pdi(float pdi_carr); //! Set Summation interval for code [s]
void initialize();
float get_carrier_nco(float PLL_discriminator);
Tracking_2nd_PLL_filter(float pdi_carr);
Tracking_2nd_PLL_filter();
~Tracking_2nd_PLL_filter();
};
#endif

View File

@ -124,7 +124,19 @@ float Tracking_FLL_PLL_filter::get_carrier_error(float FLL_discriminator, float
Tracking_FLL_PLL_filter::Tracking_FLL_PLL_filter ()
{}
{
d_order = 0;
d_pll_w = 0;
d_pll_w0p3 = 0;
d_pll_w0f2 = 0;
d_pll_x = 0;
d_pll_a2 = 0;
d_pll_w0f = 0;
d_pll_a3 = 0;
d_pll_w0p2 = 0;
d_pll_b3 = 0;
d_pll_w0p = 0;
}
Tracking_FLL_PLL_filter::~Tracking_FLL_PLL_filter ()
{}

View File

@ -98,7 +98,7 @@ void GpsL2MPcpsAcquisitionTest::init()
gnss_synchro.Channel_ID = 0;
gnss_synchro.System = 'G';
std::string signal = "2S";
strcpy(gnss_synchro.Signal,signal.c_str());
strncpy(gnss_synchro.Signal, signal.c_str(), 3);
gnss_synchro.PRN = 7;
sampling_freqeuncy_hz = 5000000;