mars_lib 0.1.0.2abe2576fe7f
Modular and Robust Sensor-Fusion
Loading...
Searching...
No Matches
gps_sensor_class.h
Go to the documentation of this file.
1// Copyright (C) 2021 Christian Brommer, Control of Networked Systems, University of Klagenfurt, Austria.
2//
3// All rights reserved.
4//
5// This software is licensed under the terms of the BSD-2-Clause-License with
6// no commercial use allowed, the full terms of which are made available
7// in the LICENSE file. No license in patents is granted.
8//
9// You can contact the author at <christian.brommer@ieee.org>
10
11#ifndef GPSSENSORCLASS_H
12#define GPSSENSORCLASS_H
13
14#include <mars/core_state.h>
15#include <mars/ekf.h>
21#include <mars/time.h>
23#include <iostream>
24#include <memory>
25#include <string>
26#include <utility>
27
28namespace mars
29{
31
33{
34public:
35 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
36
40
41 GpsSensorClass(const std::string& name, std::shared_ptr<CoreState> core_states)
42 {
43 name_ = name;
44 core_states_ = std::move(core_states);
45 const_ref_to_nav_ = false;
49
50 // chi2
51 chi2_.set_dof(3);
52
53 std::cout << "Created: [" << this->name_ << "] Sensor" << std::endl;
54 }
55
56 virtual ~GpsSensorClass() = default;
57
58 GpsSensorStateType get_state(const std::shared_ptr<void>& sensor_data)
59 {
60 GpsSensorData data = *static_cast<GpsSensorData*>(sensor_data.get());
61 return data.state_;
62 }
63
64 Eigen::MatrixXd get_covariance(const std::shared_ptr<void>& sensor_data)
65 {
66 GpsSensorData data = *static_cast<GpsSensorData*>(sensor_data.get());
67 return data.get_full_cov();
68 }
69
70 void set_initial_calib(std::shared_ptr<void> calibration)
71 {
72 initial_calib_ = calibration;
74 }
75
76 void set_gps_reference_coordinates(const double& latitude, const double& longitude, const double& altitude)
77 {
78 set_gps_reference_coordinates(mars::GpsCoordinates(latitude, longitude, altitude));
79 }
80
82 {
84 {
88 std::cout << "Info: [" << name_ << "] Set External GPS Reference: \n" << gps_reference << std::endl;
89 }
90 else
91 {
92 std::cout << "Warning: [" << name_ << "] "
93 << "Trying to set GPS reference but reference was already set. Action has no effect." << std::endl;
94 }
95 }
96
97 BufferDataType Initialize(const Time& timestamp, std::shared_ptr<void> sensor_data,
98 std::shared_ptr<CoreType> latest_core_data)
99 {
100 GpsMeasurementType measurement = *static_cast<GpsMeasurementType*>(sensor_data.get());
101
103 {
104 GpsCoordinates gps_reference(measurement.coordinates_.latitude_, measurement.coordinates_.longitude_,
105 measurement.coordinates_.altitude_);
106
107 gps_conversion_.set_gps_reference(gps_reference);
109
110 std::cout << "Info: [" << name_ << "] Set Internal GPS Reference: \n" << gps_reference << std::endl;
111 }
112
113 GpsSensorData sensor_state;
114 std::string calibration_type;
115
116 if (this->initial_calib_provided_)
117 {
118 calibration_type = "Given";
119
120 GpsSensorData calib = *static_cast<GpsSensorData*>(initial_calib_.get());
121
122 sensor_state.state_ = calib.state_;
123 sensor_state.sensor_cov_ = calib.sensor_cov_;
124 }
125 else
126 {
127 calibration_type = "Auto";
128 std::cout << "GPS calibration AUTO init not implemented yet" << std::endl;
129 exit(EXIT_FAILURE);
130 }
131
132 // Bypass core state for the returned object
133 BufferDataType result(std::make_shared<CoreType>(*latest_core_data.get()),
134 std::make_shared<GpsSensorData>(sensor_state));
135
136 is_initialized_ = true;
137
138 std::cout << "Info: Initialized [" << name_ << "] with [" << calibration_type << "] Calibration at t=" << timestamp
139 << std::endl;
140
141 std::cout << "Info: [" << name_ << "] Calibration(rounded):" << std::endl;
142 std::cout << "\tPosition[m]: [" << sensor_state.state_.p_ig_.transpose() << " ]" << std::endl;
143 std::cout << "\tReference: \n" << gps_conversion_.get_gps_reference() << std::endl;
144
145 return result;
146 }
147
148 bool CalcUpdate(const Time& /*timestamp*/, std::shared_ptr<void> measurement, const CoreStateType& prior_core_state,
149 std::shared_ptr<void> latest_sensor_data, const Eigen::MatrixXd& prior_cov,
150 BufferDataType* new_state_data)
151 {
152 // Cast the sensor measurement and prior state information
153 GpsMeasurementType* meas = static_cast<GpsMeasurementType*>(measurement.get());
154 GpsSensorData* prior_sensor_data = static_cast<GpsSensorData*>(latest_sensor_data.get());
155
156 // Decompose sensor measurement
157 Eigen::Vector3d p_meas = gps_conversion_.get_enu(meas->coordinates_);
158
159 // Extract sensor state
160 GpsSensorStateType prior_sensor_state(prior_sensor_data->state_);
161
162 // Generate measurement noise matrix and check
163 // if noisevalues from the measurement object should be used
164 Eigen::MatrixXd R_meas_dyn;
166 {
167 meas->get_meas_noise(&R_meas_dyn);
168 }
169 else
170 {
171 R_meas_dyn = this->R_.asDiagonal();
172 }
173 const Eigen::Matrix<double, 3, 3> R_meas = R_meas_dyn;
174
175 const int size_of_core_state = CoreStateType::size_error_;
176 const int size_of_sensor_state = prior_sensor_state.cov_size_;
177 const int size_of_full_error_state = size_of_core_state + size_of_sensor_state;
178 const Eigen::MatrixXd P = prior_cov;
179 assert(P.size() == size_of_full_error_state * size_of_full_error_state);
180
181 // Calculate the measurement jacobian H
182 // const Eigen::Matrix3d I_3 = Eigen::Matrix3d::Identity();
183 const Eigen::Vector3d P_wi = prior_core_state.p_wi_;
184 const Eigen::Matrix3d R_wi = prior_core_state.q_wi_.toRotationMatrix();
185 const Eigen::Vector3d P_ig = prior_sensor_state.p_ig_;
186
187 const Eigen::Vector3d P_gw_w = prior_sensor_state.p_gw_w_;
188 const Eigen::Matrix3d R_gw_w = prior_sensor_state.q_gw_w_.toRotationMatrix();
189
190 // Position
191 const Eigen::Matrix3d Hp_pwi = R_gw_w;
192 const Eigen::Matrix3d Hp_vwi = Eigen::Matrix3d::Zero();
193 const Eigen::Matrix3d Hp_rwi = -R_gw_w * R_wi * Utils::Skew(P_ig);
194 const Eigen::Matrix3d Hp_bw = Eigen::Matrix3d::Zero();
195 const Eigen::Matrix3d Hp_ba = Eigen::Matrix3d::Zero();
196
197 const Eigen::Matrix3d Hp_pig = R_gw_w * R_wi;
198 const Eigen::Matrix3d Hp_pgw_w = Eigen::Matrix3d::Zero(); // I_3;
199 const Eigen::Matrix3d Hp_rgw_w = Eigen::Matrix3d::Zero(); // R_gw_w * Utils::Skew(P_wi + R_wi * P_ig);
200
201 // Assemble the jacobian for the position (horizontal)
202 // H_p = [Hp_pwi Hp_vwi Hp_rwi Hp_bw Hp_ba Hp_pig Hp_pgw_w Hp_rgw_w ];
203 Eigen::MatrixXd H(3, Hp_pwi.cols() + Hp_vwi.cols() + Hp_rwi.cols() + Hp_bw.cols() + Hp_ba.cols() + Hp_pig.cols() +
204 Hp_pgw_w.cols() + Hp_rgw_w.cols());
205
206 H << Hp_pwi, Hp_vwi, Hp_rwi, Hp_bw, Hp_ba, Hp_pig, Hp_pgw_w, Hp_rgw_w;
207
208 // Calculate the residual z = z~ - (estimate)
209 // Position
210 const Eigen::Vector3d p_est = P_gw_w + R_gw_w * (P_wi + R_wi * P_ig);
211 residual_ = Eigen::MatrixXd(p_est.rows(), 1);
212 residual_ = p_meas - p_est;
213
214 // Perform EKF calculations
215 mars::Ekf ekf(H, R_meas, residual_, P);
216 const Eigen::MatrixXd correction = ekf.CalculateCorrection(&chi2_);
217 assert(correction.size() == size_of_full_error_state * 1);
218
219 // Perform Chi2 test
220 if (!chi2_.passed_ && chi2_.do_test_)
221 {
223 return false;
224 }
225
226 Eigen::MatrixXd P_updated = ekf.CalculateCovUpdate();
227 assert(P_updated.size() == size_of_full_error_state * size_of_full_error_state);
228 P_updated = Utils::EnforceMatrixSymmetry(P_updated);
229
230 // Apply Core Correction
231 CoreStateVector core_correction = correction.block(0, 0, CoreStateType::size_error_, 1);
232 CoreStateType corrected_core_state = CoreStateType::ApplyCorrection(prior_core_state, core_correction);
233
234 // Apply Sensor Correction
235 const Eigen::MatrixXd sensor_correction = correction.block(size_of_core_state, 0, size_of_sensor_state, 1);
236 const GpsSensorStateType corrected_sensor_state = ApplyCorrection(prior_sensor_state, sensor_correction);
237
238 // Return Results
239 // CoreState data
240 CoreType core_data;
241 core_data.cov_ = P_updated.block(0, 0, CoreStateType::size_error_, CoreStateType::size_error_);
242 core_data.state_ = corrected_core_state;
243
244 // SensorState data
245 std::shared_ptr<GpsSensorData> sensor_data(std::make_shared<GpsSensorData>());
246 sensor_data->set_cov(P_updated);
247 sensor_data->state_ = corrected_sensor_state;
248
249 BufferDataType state_entry(std::make_shared<CoreType>(core_data), sensor_data);
250
252 {
253 // corrected_sensor_data.ref_to_nav = prior_ref_to_nav;
254 }
255 else
256 {
257 // TODO also estimate ref to nav
258 }
259
260 *new_state_data = state_entry;
261
262 return true;
263 }
264
265 GpsSensorStateType ApplyCorrection(const GpsSensorStateType& prior_sensor_state, const Eigen::MatrixXd& correction)
266 {
267 // state + error state correction
268 // with quaternion from small angle approx -> new state
269
270 GpsSensorStateType corrected_sensor_state;
271 corrected_sensor_state.p_ig_ = prior_sensor_state.p_ig_ + correction.block(0, 0, 3, 1);
272 corrected_sensor_state.p_gw_w_ = prior_sensor_state.p_gw_w_ + correction.block(3, 0, 3, 1);
273 corrected_sensor_state.q_gw_w_ =
274 Utils::ApplySmallAngleQuatCorr(prior_sensor_state.q_gw_w_, correction.block(6, 0, 3, 1));
275 return corrected_sensor_state;
276 }
277};
278} // namespace mars
279
280#endif // GPSSENSORCLASS_H
bool has_meas_noise
Definition measurement_base_class.h:23
bool get_meas_noise(Eigen::MatrixXd *meas_noise)
get the measurement noise associated with the current sensor measurement
Definition measurement_base_class.h:25
int cov_size_
Definition base_states.h:25
The BaseSensorData class binds the sensor state and covariance matrix.
Definition bind_sensor_data.h:29
EIGEN_MAKE_ALIGNED_OPERATOR_NEW T state_
Definition bind_sensor_data.h:35
Eigen::MatrixXd get_full_cov() const
get_full_cov builds the full covariance matrix
Definition bind_sensor_data.h:63
Eigen::MatrixXd sensor_cov_
covariance of the sensor states
Definition bind_sensor_data.h:37
The BufferDataType binds the core and sensor state in form of a shared void pointer.
Definition buffer_data_type.h:36
bool passed_
Determine if the test is performed or not.
Definition ekf.h:84
bool do_test_
Upper critival value.
Definition ekf.h:83
void set_dof(const int &value)
set_dof Set degree of freedom for the X2 distribution
void PrintReport(const std::string &name)
PrintReport Print a formated report e.g. if the test did not pass.
Definition core_state_type.h:21
static constexpr int size_error_
Definition core_state_type.h:38
Eigen::Vector3d p_wi_
Definition core_state_type.h:27
Eigen::Quaternion< double > q_wi_
Definition core_state_type.h:29
static CoreStateType ApplyCorrection(CoreStateType state_prior, Eigen::Matrix< double, CoreStateType::size_error_, 1 > correction)
ApplyCorrection.
Definition core_state_type.h:46
Definition core_type.h:19
CoreStateMatrix cov_
Definition core_type.h:22
CoreStateType state_
Definition core_type.h:21
Definition ekf.h:92
Eigen::MatrixXd CalculateCovUpdate()
CalculateCovUpdate Updating the state covariance after the state update.
Eigen::MatrixXd CalculateCorrection()
Kalman gain.
The GpsConversion class.
Definition gps_conversion.h:66
void set_gps_reference(mars::GpsCoordinates coordinates)
set_gps_reference
mars::GpsCoordinates get_gps_reference()
get_gps_reference
Eigen::Matrix< double, 3, 1 > get_enu(mars::GpsCoordinates coordinates)
get_enu get current GPS reference coordinates
Definition gps_measurement_type.h:21
GpsCoordinates coordinates_
Definition gps_measurement_type.h:23
Definition gps_sensor_class.h:33
GpsSensorStateType ApplyCorrection(const GpsSensorStateType &prior_sensor_state, const Eigen::MatrixXd &correction)
Definition gps_sensor_class.h:265
EIGEN_MAKE_ALIGNED_OPERATOR_NEW GpsConversion gps_conversion_
Definition gps_sensor_class.h:37
GpsSensorStateType get_state(const std::shared_ptr< void > &sensor_data)
Definition gps_sensor_class.h:58
bool CalcUpdate(const Time &, std::shared_ptr< void > measurement, const CoreStateType &prior_core_state, std::shared_ptr< void > latest_sensor_data, const Eigen::MatrixXd &prior_cov, BufferDataType *new_state_data)
CalcUpdate Calculates the update for an individual sensor definition.
Definition gps_sensor_class.h:148
void set_gps_reference_coordinates(const double &latitude, const double &longitude, const double &altitude)
Definition gps_sensor_class.h:76
void set_initial_calib(std::shared_ptr< void > calibration)
set_initial_calib Sets the calibration of an individual sensor
Definition gps_sensor_class.h:70
bool gps_reference_is_set_
Definition gps_sensor_class.h:39
GpsSensorClass(const std::string &name, std::shared_ptr< CoreState > core_states)
Definition gps_sensor_class.h:41
void set_gps_reference_coordinates(const mars::GpsCoordinates &gps_reference)
Definition gps_sensor_class.h:81
Eigen::MatrixXd get_covariance(const std::shared_ptr< void > &sensor_data)
get_covariance Resolves a void pointer to the covariance matrix of the corresponding sensor type Each...
Definition gps_sensor_class.h:64
virtual ~GpsSensorClass()=default
bool using_external_gps_reference_
Definition gps_sensor_class.h:38
BufferDataType Initialize(const Time &timestamp, std::shared_ptr< void > sensor_data, std::shared_ptr< CoreType > latest_core_data)
Initialize the state of an individual sensor.
Definition gps_sensor_class.h:97
Definition gps_sensor_state_type.h:20
Eigen::Quaterniond q_gw_w_
Definition gps_sensor_state_type.h:26
Eigen::Vector3d p_gw_w_
Definition gps_sensor_state_type.h:25
EIGEN_MAKE_ALIGNED_OPERATOR_NEW Eigen::Vector3d p_ig_
Definition gps_sensor_state_type.h:24
std::string name_
Name of the individual sensor instance.
Definition sensor_abs_class.h:23
bool is_initialized_
True if the sensor has been initialized.
Definition sensor_abs_class.h:24
bool const_ref_to_nav_
True if the reference should not be estimated.
Definition sensor_abs_class.h:27
bool use_dynamic_meas_noise_
True if dynamic noise values from measurements should be used.
Definition sensor_abs_class.h:29
Definition time.h:20
Definition update_sensor_abs_class.h:24
bool initial_calib_provided_
True if an initial calibration was provided.
Definition update_sensor_abs_class.h:38
std::shared_ptr< CoreState > core_states_
Definition update_sensor_abs_class.h:42
Eigen::VectorXd R_
Measurement noise "squared".
Definition update_sensor_abs_class.h:32
std::shared_ptr< void > initial_calib_
Definition update_sensor_abs_class.h:37
Chi2 chi2_
Definition update_sensor_abs_class.h:40
Eigen::MatrixXd residual_
Definition update_sensor_abs_class.h:31
static Eigen::MatrixXd EnforceMatrixSymmetry(const Eigen::Ref< const Eigen::MatrixXd > &mat_in)
EnforceMatrixSymmetry.
static Eigen::Matrix3d Skew(const Eigen::Vector3d &v)
skew generate the skew symmetric matrix of v
static Eigen::Quaterniond ApplySmallAngleQuatCorr(const Eigen::Quaterniond &q_prior, const Eigen::Vector3d &correction)
ApplySmallAngleQuatCorr.
Definition buffer.h:27
Eigen::Matrix< double, CoreStateType::size_error_, 1 > CoreStateVector
Definition core_state_type.h:135
The GpsCoordinates struct.
Definition gps_conversion.h:22
double longitude_
Definition gps_conversion.h:29
double altitude_
Definition gps_conversion.h:30
double latitude_
Definition gps_conversion.h:28