0.3.0
Loading...
Searching...
No Matches
ObservationEstimator.hpp
Go to the documentation of this file.
1// This file is part of INSTINCT, the INS Toolkit for Integrated
2// Navigation Concepts and Training by the Institute of Navigation of
3// the University of Stuttgart, Germany.
4//
5// This Source Code Form is subject to the terms of the Mozilla Public
6// License, v. 2.0. If a copy of the MPL was not distributed with this
7// file, You can obtain one at https://mozilla.org/MPL/2.0/.
8
13
14#pragma once
15
16#include <algorithm>
17#include <cmath>
18#include <cstddef>
19#include <cstdint>
20#include <memory>
21#include <unordered_map>
22#include <utility>
23#include <vector>
24
25#include <imgui.h>
26#include <application.h>
27#include <imgui_stdlib.h>
36
43
45#include "util/Json.hpp"
46#include <Eigen/src/Core/MatrixBase.h>
47#include <fmt/core.h>
48#include <fmt/format.h>
49
50namespace NAV
51{
52
55{
56 public:
59 explicit ObservationEstimator(size_t receiverCount)
60 : _receiverCount(receiverCount)
61 {
62 for (size_t i = 0; i < receiverCount; i++)
63 {
64 _antenna.emplace(i, Antenna{});
65 }
66 }
67
75
82 template<typename ReceiverType>
84 const Receiver<ReceiverType>& receiver,
85 const std::shared_ptr<const IonosphericCorrections>& ionosphericCorrections,
86 [[maybe_unused]] const std::string& nameId,
87 ObservationDifference obsDiff) const
88 {
89 LOG_DATA("{}: Calculating all observation estimates:", nameId);
90 for (auto& [satSigId, observation] : observations.signals)
91 {
92 calcObservationEstimate(satSigId, observation, receiver, ionosphericCorrections, nameId, obsDiff);
93 }
94 }
95
103 template<typename ReceiverType>
104 void calcObservationEstimate(const SatSigId& satSigId,
106 const Receiver<ReceiverType>& receiver,
107 const std::shared_ptr<const IonosphericCorrections>& ionosphericCorrections,
108 [[maybe_unused]] const std::string& nameId,
109 ObservationDifference obsDiff) const
110 {
111 LOG_DATA("{}: Calculating observation estimates for [{}]", nameId, satSigId);
112 const Frequency freq = satSigId.freq();
113 const SatelliteSystem satSys = freq.getSatSys();
114
115 auto& recvObs = observation.recvObs.at(receiver.type);
116
117 Eigen::Vector3d e_recvPosAPC = e_calcRecvPosAPC(freq,
118 receiver.type,
119 receiver.e_posMarker,
120 receiver.gnssObs, nameId);
121 Eigen::Vector3d lla_recvPosAPC = trafo::ecef2lla_WGS84(e_recvPosAPC);
122
123 double elevation = recvObs->satElevation(e_recvPosAPC, lla_recvPosAPC);
124 double azimuth = recvObs->satAzimuth(e_recvPosAPC, lla_recvPosAPC);
125 Eigen::Vector3d e_pLOS = recvObs->e_pLOS(e_recvPosAPC);
126
127 // Receiver-Satellite Range [m]
128 double rho_r_s = (recvObs->e_satPos() - e_recvPosAPC).norm();
129 recvObs->terms.rho_r_s = rho_r_s;
130 // Troposphere
131 auto tropo_r_s = calcTroposphericDelayAndMapping(receiver.gnssObs->insTime, lla_recvPosAPC,
132 elevation, azimuth, _troposphereModels, nameId);
133 recvObs->terms.tropoZenithDelay = tropo_r_s;
134 // Estimated troposphere propagation error [m]
135 double dpsr_T_r_s = tropo_r_s.ZHD * tropo_r_s.zhdMappingFactor + tropo_r_s.ZWD * tropo_r_s.zwdMappingFactor;
136 recvObs->terms.dpsr_T_r_s = dpsr_T_r_s;
137 // Estimated ionosphere propagation error [m]
138 double dpsr_I_r_s = calcIonosphericDelay(static_cast<double>(receiver.gnssObs->insTime.toGPSweekTow().tow),
139 freq, observation.freqNum(), lla_recvPosAPC, elevation, azimuth,
140 _ionosphereModel, ionosphericCorrections.get());
141 recvObs->terms.dpsr_I_r_s = dpsr_I_r_s;
142 // Sagnac correction [m]
143 double dpsr_ie_r_s = calcSagnacCorrection(e_recvPosAPC, recvObs->e_satPos());
144 recvObs->terms.dpsr_ie_r_s = dpsr_ie_r_s;
145
146 // Earth's gravitational field causes relativistic signal delay due to space-time curvature (Shapiro effect) [s]
147 // double posNorm = recvObs->e_satPos().norm() + e_recvPosAPC.norm();
148 // double dt_rel_stc = e_recvPosAPC.norm() > InsConst::WGS84::a / 2.0 ? 2.0 * InsConst::WGS84::MU / std::pow(InsConst::C, 3) * std::log((posNorm + rho_r_s) / (posNorm - rho_r_s))
149 // : 0.0;
150
151 double cn0 = recvObs->gnssObsData().CN0.value_or(1.0);
152
153 double pcv = calcPCV(freq, receiver.gnssObs->insTime, elevation, std::nullopt, receiver.type, receiver.gnssObs, nameId); // TODO: use 'azimuth', but need antenna azimuth
154
155 for (auto& [obsType, obsData] : recvObs->obs)
156 {
157 LOG_DATA("{}: [{}][{:11}][{:5}] Observation estimate", nameId, satSigId, obsType, receiver.type);
158 switch (obsType)
159 {
161 obsData.estimate = rho_r_s
162 + dpsr_ie_r_s
163 + dpsr_T_r_s
164 + dpsr_I_r_s
166 * (*receiver.recvClk.biasFor(satSys) - recvObs->satClock().bias
167 // + dt_rel_stc
168 + (receiver.interFrequencyBias.contains(freq) ? receiver.interFrequencyBias.at(freq).value : 0.0))
169 + pcv;
170 obsData.measVar = _gnssMeasurementErrorModel.psrMeasErrorVar(satSys, elevation, cn0);
171 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m] Geometrical range", nameId, satSigId, obsType, receiver.type, rho_r_s);
172 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Sagnac correction", nameId, satSigId, obsType, receiver.type, dpsr_ie_r_s);
173 if (dpsr_T_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Tropospheric delay", nameId, satSigId, obsType, receiver.type, dpsr_T_r_s); }
174 if (dpsr_I_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Ionospheric delay", nameId, satSigId, obsType, receiver.type, dpsr_I_r_s); }
175 if (*receiver.recvClk.biasFor(satSys) != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Receiver clock bias", nameId, satSigId, obsType, receiver.type, InsConst::C * *receiver.recvClk.biasFor(satSys)); }
176 LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m] Satellite clock bias", nameId, satSigId, obsType, receiver.type, InsConst::C * recvObs->satClock().bias);
177 // LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Shapiro effect", nameId, satSigId, obsType, receiver.type, InsConst::C * dt_rel_stc);
178 if (receiver.interFrequencyBias.contains(freq)) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Inter-frequency bias", nameId, satSigId, obsType, receiver.type, InsConst::C * receiver.interFrequencyBias.at(freq).value); }
179 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4f} [m] Pseudorange estimate", nameId, satSigId, obsType, receiver.type, obsData.estimate);
180 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m] Pseudorange measurement", nameId, satSigId, obsType, receiver.type, obsData.measurement);
181 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4e} [m] Difference to measurement", nameId, satSigId, obsType, receiver.type, obsData.measurement - obsData.estimate);
182 break;
183 case GnssObs::Carrier:
184 obsData.estimate = rho_r_s
185 + dpsr_ie_r_s
186 + dpsr_T_r_s
187 - dpsr_I_r_s
189 * (*receiver.recvClk.biasFor(satSys) - recvObs->satClock().bias
190 // + dt_rel_stc
191 )
192 + pcv;
193 obsData.measVar = _gnssMeasurementErrorModel.carrierMeasErrorVar(satSys, elevation, cn0);
194 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m] Geometrical range", nameId, satSigId, obsType, receiver.type, rho_r_s);
195 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Sagnac correction", nameId, satSigId, obsType, receiver.type, dpsr_ie_r_s);
196 if (dpsr_T_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Tropospheric delay", nameId, satSigId, obsType, receiver.type, dpsr_T_r_s); }
197 if (dpsr_I_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m] Ionospheric delay", nameId, satSigId, obsType, receiver.type, dpsr_I_r_s); }
198 if (*receiver.recvClk.biasFor(satSys) != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Receiver clock bias", nameId, satSigId, obsType, receiver.type, InsConst::C * *receiver.recvClk.biasFor(satSys)); }
199 LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m] Satellite clock bias", nameId, satSigId, obsType, receiver.type, InsConst::C * recvObs->satClock().bias);
200 // LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Shapiro effect", nameId, satSigId, obsType, receiver.type, InsConst::C * dt_rel_stc);
201 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4f} [m] Carrier-phase estimate", nameId, satSigId, obsType, receiver.type, obsData.estimate);
202 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m] Carrier-phase measurement", nameId, satSigId, obsType, receiver.type, obsData.measurement);
203 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4e} [m] Difference to measurement", nameId, satSigId, obsType, receiver.type, obsData.measurement - obsData.estimate);
204 break;
205 case GnssObs::Doppler:
206 obsData.estimate = e_pLOS.dot(recvObs->e_satVel() - receiver.e_vel)
207 - calcSagnacRateCorrection(e_recvPosAPC, recvObs->e_satPos(), receiver.e_vel, recvObs->e_satVel())
209 * (*receiver.recvClk.driftFor(satSys)
210 - recvObs->satClock().drift);
211 obsData.measVar = _gnssMeasurementErrorModel.psrRateMeasErrorVar(freq, observation.freqNum(), elevation, cn0);
212 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m/s] ", nameId, satSigId, obsType, receiver.type, e_pLOS.dot(recvObs->e_satVel() - receiver.e_vel));
213 LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m/s] Sagnac rate correction", nameId, satSigId, obsType, receiver.type, calcSagnacRateCorrection(e_recvPosAPC, recvObs->e_satPos(), receiver.e_vel, recvObs->e_satVel()));
214 if (*receiver.recvClk.driftFor(satSys) != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m/s] Receiver clock drift", nameId, satSigId, obsType, receiver.type, InsConst::C * *receiver.recvClk.driftFor(satSys)); }
215 LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m/s] Satellite clock drift", nameId, satSigId, obsType, receiver.type, InsConst::C * recvObs->satClock().drift);
216 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4f} [m/s] Doppler estimate", nameId, satSigId, obsType, receiver.type, obsData.estimate);
217 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m/s] Doppler measurement", nameId, satSigId, obsType, receiver.type, obsData.measurement);
218 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4e} [m/s] Difference to measurement", nameId, satSigId, obsType, receiver.type, obsData.measurement - obsData.estimate);
219 break;
221 break;
222 }
223 LOG_DATA("{}: [{}][{:11}][{:5}] Observation error variance", nameId, satSigId, obsType, receiver.type);
224 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4g} [{}] Measurement error variance", nameId, satSigId, obsType, receiver.type, obsData.measVar,
225 obsType == GnssObs::Doppler ? "m^2/s^2" : "m^2");
226
227 if (obsDiff == NoDifference)
228 {
229 if (obsType == GnssObs::Pseudorange || obsType == GnssObs::Carrier)
230 {
231 obsData.measVar += observation.navData()->calcSatellitePositionVariance()
232 + ionoErrorVar(dpsr_I_r_s, freq, observation.freqNum())
233 + tropoErrorVar(dpsr_T_r_s, elevation);
234 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Satellite position variance", nameId, satSigId, obsType, receiver.type, observation.navData()->calcSatellitePositionVariance());
235 if (dpsr_I_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Ionosphere variance", nameId, satSigId, obsType, receiver.type, ionoErrorVar(dpsr_I_r_s, freq, observation.freqNum())); }
236 if (dpsr_T_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Troposphere variance", nameId, satSigId, obsType, receiver.type, tropoErrorVar(dpsr_T_r_s, elevation)); }
237 }
238 if (obsType == GnssObs::Pseudorange)
239 {
240 obsData.measVar += _gnssMeasurementErrorModel.codeBiasErrorVar();
241 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Code bias variance", nameId, satSigId, obsType, receiver.type, _gnssMeasurementErrorModel.codeBiasErrorVar());
242 }
243 }
244 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4g} [{}] Observation error variance", nameId, satSigId, obsType, receiver.type, obsData.measVar,
245 obsType == GnssObs::Doppler ? "m^2/s^2" : "m^2");
246 }
247 }
248
251 {
259
260 Eigen::Vector3d e_satPos;
261 Eigen::Vector3d e_satVel;
262 double satClockBias{};
263 double satClockDrift{};
264 };
265
278 template<typename Derived>
280 size_t receiverType,
281 const Eigen::MatrixBase<Derived>& e_recvPosMarker,
282 auto recvClockBias,
283 auto interFreqBias,
284 auto dpsr_T_r_s,
285 auto dpsr_I_r_s,
286 const std::shared_ptr<const SatelliteInfo>& satInfo,
287 const std::shared_ptr<const GnssObs>& gnssObs,
288 [[maybe_unused]] const std::string& nameId) const
289 {
290 Eigen::Vector3d e_recvPosAPC = e_calcRecvPosAPC(freq, receiverType, e_recvPosMarker, gnssObs, nameId);
291 Eigen::Vector3d lla_recvPosAPC = trafo::ecef2lla_WGS84(e_recvPosAPC);
292
293 Eigen::Vector3d e_pLOS = e_calcLineOfSightUnitVector(e_recvPosAPC, satInfo->e_satPos);
294 Eigen::Vector3d n_lineOfSightUnitVector = trafo::n_Quat_e(lla_recvPosAPC(0), lla_recvPosAPC(1)) * e_pLOS;
295 double elevation = calcSatElevation(n_lineOfSightUnitVector);
296 // double azimuth = calcSatAzimuth(n_lineOfSightUnitVector);
297
298 // Receiver-Satellite Range [m]
299 auto rho_r_s = (satInfo->e_satPos - e_recvPosAPC).norm();
300
301 // Sagnac correction [m]
302 auto dpsr_ie_r_s = calcSagnacCorrection(e_recvPosAPC, satInfo->e_satPos);
303
304 // Earth's gravitational field causes relativistic signal delay due to space-time curvature (Shapiro effect) [s]
305 // double posNorm = e_satPos.norm() + e_recvPosAPC.norm();
306 // double dt_rel_stc = e_recvPosAPC.norm() > InsConst::WGS84::a / 2.0 ? 2.0 * InsConst::WGS84::MU / std::pow(InsConst::C, 3) * std::log((posNorm + rho_r_s) / (posNorm - rho_r_s))
307 // : 0.0;
308
309 double pcv = calcPCV(freq, gnssObs->insTime, elevation, std::nullopt, receiverType, gnssObs, nameId); // TODO: use 'azimuth', but need antenna azimuth
310
311 return rho_r_s
312 + dpsr_ie_r_s
313 + dpsr_T_r_s
314 + dpsr_I_r_s
315 + InsConst::C * (recvClockBias - satInfo->satClockBias /* + dt_rel_stc */ + interFreqBias)
316 + pcv;
317 }
318
330 template<typename Derived>
331 std::pair<double, Eigen::Vector3d>
333 int8_t freqNum,
334 size_t receiverType,
335 const Eigen::MatrixBase<Derived>& e_recvPosMarker,
336 double recvClockBias,
337 const std::shared_ptr<const SatelliteInfo>& satInfo,
338 const std::shared_ptr<const GnssObs>& gnssObs,
339 const std::shared_ptr<const IonosphericCorrections>& ionosphericCorrections,
340 [[maybe_unused]] const std::string& nameId) const
341 {
342 Eigen::Vector3d e_recvPosAPC = e_calcRecvPosAPC(freq, receiverType, e_recvPosMarker, gnssObs, nameId);
343 Eigen::Vector3d lla_recvPosAPC = trafo::ecef2lla_WGS84(e_recvPosAPC);
344
345 Eigen::Vector3d e_pLOS = e_calcLineOfSightUnitVector(e_recvPosAPC, satInfo->e_satPos);
346 Eigen::Vector3d n_lineOfSightUnitVector = trafo::n_Quat_e(lla_recvPosAPC(0), lla_recvPosAPC(1)) * e_pLOS;
347 double elevation = calcSatElevation(n_lineOfSightUnitVector);
348 double azimuth = calcSatAzimuth(n_lineOfSightUnitVector);
349
350 // Receiver-Satellite Range [m]
351 double rho_r_s = (satInfo->e_satPos - e_recvPosAPC).norm();
352 // Estimated troposphere propagation error [m]
353 double dpsr_T_r_s = calculateTroposphericDelay(lla_recvPosAPC, gnssObs, elevation, azimuth, nameId);
354 // Estimated ionosphere propagation error [m]
355 double dpsr_I_r_s = calculateIonosphericDelay(freq, freqNum, lla_recvPosAPC, gnssObs, ionosphericCorrections, elevation, azimuth);
356 // Sagnac correction [m]
357 double dpsr_ie_r_s = calcSagnacCorrection(e_recvPosAPC, satInfo->e_satPos);
358
359 // Earth's gravitational field causes relativistic signal delay due to space-time curvature (Shapiro effect) [s]
360 // double posNorm = e_satPos.norm() + e_recvPosAPC.norm();
361 // double dt_rel_stc = e_recvPosAPC.norm() > InsConst::WGS84::a / 2.0 ? 2.0 * InsConst::WGS84::MU / std::pow(InsConst::C, 3) * std::log((posNorm + rho_r_s) / (posNorm - rho_r_s))
362 // : 0.0;
363
364 double pcv = calcPCV(freq, gnssObs->insTime, elevation, std::nullopt, receiverType, gnssObs, nameId); // TODO: use 'azimuth', but need antenna azimuth
365
366 double estimate = rho_r_s
367 + dpsr_ie_r_s
368 + dpsr_T_r_s
369 - dpsr_I_r_s
370 + InsConst::C * (recvClockBias - satInfo->satClockBias /* + dt_rel_stc */)
371 + pcv;
372
373 return { estimate, e_pLOS };
374 }
375
386 template<typename DerivedPos, typename DerivedVel>
388 size_t receiverType,
389 const Eigen::MatrixBase<DerivedPos>& e_recvPosMarker,
390 const Eigen::MatrixBase<DerivedVel>& e_recvVel,
391 auto recvClockDrift,
392 const std::shared_ptr<const SatelliteInfo>& satInfo,
393 const std::shared_ptr<const GnssObs>& gnssObs,
394 [[maybe_unused]] const std::string& nameId) const
395 {
396 auto e_recvPosAPC = e_calcRecvPosAPC(freq, receiverType, e_recvPosMarker, gnssObs, nameId);
397
398 auto e_pLOS = e_calcLineOfSightUnitVector(e_recvPosAPC, satInfo->e_satPos);
399
400 return e_pLOS.dot(satInfo->e_satVel - e_recvVel)
401 - calcSagnacRateCorrection(e_recvPosAPC, satInfo->e_satPos, e_recvVel, satInfo->e_satVel)
402 + InsConst::C * (recvClockDrift - satInfo->satClockDrift);
403 }
404
420 int8_t freqNum,
421 size_t receiverType,
423 const Eigen::Vector3d& e_recvPosMarker,
424 const Eigen::Vector3d& e_satPos,
425 double cn0,
426 const std::shared_ptr<const GnssObs>& gnssObs,
427 const std::shared_ptr<const SatNavData>& navData,
428 const std::shared_ptr<const IonosphericCorrections>& ionosphericCorrections,
429 [[maybe_unused]] const std::string& nameId,
430 ObservationDifference obsDiff) const
431 {
432 const SatelliteSystem satSys = freq.getSatSys();
433
434 Eigen::Vector3d e_recvPosAPC = e_calcRecvPosAPC(freq, receiverType, e_recvPosMarker, gnssObs, nameId);
435 Eigen::Vector3d lla_recvPosAPC = trafo::ecef2lla_WGS84(e_recvPosAPC);
436
437 Eigen::Vector3d e_pLOS = e_calcLineOfSightUnitVector(e_recvPosAPC, e_satPos);
438 Eigen::Vector3d n_lineOfSightUnitVector = trafo::n_Quat_e(lla_recvPosAPC(0), lla_recvPosAPC(1)) * e_pLOS;
439 double elevation = calcSatElevation(n_lineOfSightUnitVector);
440 double azimuth = calcSatAzimuth(n_lineOfSightUnitVector);
441
442 // Estimated troposphere propagation error [m]
443 double dpsr_T_r_s = calculateTroposphericDelay(lla_recvPosAPC, gnssObs, elevation, azimuth, nameId);
444 // Estimated ionosphere propagation error [m]
445 double dpsr_I_r_s = calculateIonosphericDelay(freq, freqNum, lla_recvPosAPC, gnssObs, ionosphericCorrections, elevation, azimuth);
446
447 double variance = 0.0;
448
449 switch (obsType)
450 {
452 variance = _gnssMeasurementErrorModel.psrMeasErrorVar(satSys, elevation, cn0);
453 break;
454 case GnssObs::Carrier:
455 variance = _gnssMeasurementErrorModel.carrierMeasErrorVar(satSys, elevation, cn0);
456 break;
457 case GnssObs::Doppler:
458 variance = _gnssMeasurementErrorModel.psrRateMeasErrorVar(freq, freqNum, elevation, cn0);
459 break;
461 break;
462 }
463
464 if (obsDiff == NoDifference)
465 {
466 if (obsType == GnssObs::Pseudorange || obsType == GnssObs::Carrier)
467 {
468 variance += navData->calcSatellitePositionVariance()
469 + ionoErrorVar(dpsr_I_r_s, freq, freqNum)
470 + tropoErrorVar(dpsr_T_r_s, elevation);
471 }
472 if (obsType == GnssObs::Pseudorange)
473 {
475 }
476 }
477
478 return variance;
479 }
480
487 template<typename Derived>
488 [[nodiscard]] Eigen::Vector3<typename Derived::Scalar> e_calcRecvPosAPC(const Frequency& freq,
489 size_t receiverType,
490 const Eigen::MatrixBase<Derived>& e_recvPosMarker,
491 const std::shared_ptr<const GnssObs>& gnssObs,
492 [[maybe_unused]] const std::string& nameId) const
493 {
494 const Antenna& antenna = _antenna.at(receiverType);
495
496 auto e_recvPosAPC = trafo::e_posMarker2ARP(e_recvPosMarker,
497 gnssObs,
498 antenna.hen_delta);
499 if (antenna.correctPCO)
500 {
501 std::string antennaType;
502 if (antenna.autoDetermine && gnssObs->receiverInfo)
503 {
504 antennaType = gnssObs->receiverInfo->get().antennaType;
505 }
506 else if (!antenna.autoDetermine)
507 {
508 antennaType = antenna.name;
509 }
510 e_recvPosAPC = trafo::e_posARP2APC(e_recvPosAPC, gnssObs, freq, antennaType, nameId);
511 }
512
513 return e_recvPosAPC;
514 }
515
525 [[nodiscard]] double calcPCV(const Frequency& freq,
526 const InsTime& insTime,
527 double elevation,
528 std::optional<double> azimuth,
529 size_t receiverType,
530 const std::shared_ptr<const GnssObs>& gnssObs,
531 [[maybe_unused]] const std::string& nameId) const
532 {
533 if (std::isnan(elevation)
534 || (azimuth.has_value() && std::isnan(*azimuth))) { return 0.0; }
535 const Antenna& antenna = _antenna.at(receiverType);
536 if (antenna.correctPCV)
537 {
538 std::string antennaType;
539 if (antenna.autoDetermine && gnssObs->receiverInfo)
540 {
541 antennaType = gnssObs->receiverInfo->get().antennaType;
542 }
543 else if (!antenna.autoDetermine)
544 {
545 antennaType = antenna.name;
546 }
548 Frequency_(freq),
549 insTime,
550 elevation,
551 azimuth,
552 nameId)
553 .value_or(0.0);
554 }
555
556 return 0.0;
557 }
558
569 int8_t freqNum,
570 const Eigen::Vector3d& lla_recvPosAPC,
571 const std::shared_ptr<const GnssObs>& gnssObs,
572 const std::shared_ptr<const IonosphericCorrections>& ionosphericCorrections,
573 double elevation,
574 double azimuth) const
575 {
576 return calcIonosphericDelay(static_cast<double>(gnssObs->insTime.toGPSweekTow().tow),
577 freq, freqNum, lla_recvPosAPC,
578 elevation, azimuth,
579 _ionosphereModel, ionosphericCorrections.get());
580 }
581
589 double calculateTroposphericDelay(const Eigen::Vector3d& lla_recvPosAPC,
590 const std::shared_ptr<const GnssObs>& gnssObs,
591 double elevation,
592 double azimuth,
593 [[maybe_unused]] const std::string& nameId) const
594 {
595 auto tropo_r_s = calcTroposphericDelayAndMapping(gnssObs->insTime, lla_recvPosAPC,
596 elevation, azimuth, _troposphereModels, nameId);
597
598 return tropo_r_s.ZHD * tropo_r_s.zhdMappingFactor + tropo_r_s.ZWD * tropo_r_s.zwdMappingFactor;
599 }
600
604 template<typename ReceiverType>
605 bool ShowGuiWidgets(const char* id, float itemWidth)
606 {
607 bool changed = false;
608
609 ImGui::SetNextItemOpen(true, ImGuiCond_FirstUseEver);
610 if (ImGui::TreeNode(fmt::format("Compensation models##{}", id).c_str()))
611 {
612 ImGui::SetNextItemWidth(itemWidth - ImGui::GetStyle().IndentSpacing);
613 if (ComboIonosphereModel(fmt::format("Ionosphere Model##{}", id).c_str(), _ionosphereModel))
614 {
615 LOG_DEBUG("{}: Ionosphere Model changed to {}", id, NAV::to_string(_ionosphereModel));
616 changed = true;
617 }
618 if (ComboTroposphereModel(fmt::format("Troposphere Model##{}", id).c_str(), _troposphereModels, itemWidth - ImGui::GetStyle().IndentSpacing))
619 {
620 changed = true;
621 }
622 ImGui::TreePop();
623 }
624
625 ImGui::SetNextItemOpen(true, ImGuiCond_FirstUseEver);
626 if (ImGui::TreeNode(fmt::format("GNSS Measurement Error##{}", id).c_str()))
627 {
628 if (_gnssMeasurementErrorModel.ShowGuiWidgets(id, itemWidth - ImGui::GetStyle().IndentSpacing))
629 {
630 LOG_DEBUG("{}: GNSS Measurement Error Model changed.", id);
631 changed = true;
632 }
633 ImGui::TreePop();
634 }
635
636 if (ImGui::TreeNode(fmt::format("Antenna settings##{}", id).c_str()))
637 {
638 for (size_t i = 0; i < _receiverCount; ++i)
639 {
640 if (_receiverCount > 1)
641 {
642 ImGui::SetNextItemOpen(true, ImGuiCond_FirstUseEver);
643 if (!ImGui::TreeNode(fmt::format("{}", static_cast<ReceiverType>(i)).c_str())) { continue; }
644 }
645
646 if (ImGui::Checkbox(fmt::format("Correct Phase-center offset##antenna {} {}", id, i).c_str(), &_antenna.at(i).correctPCO))
647 {
648 LOG_DEBUG("{}: Antenna {} correctPCO: {}", id, i, _antenna.at(i).correctPCO);
649 changed = true;
650 }
651 ImGui::SameLine();
652 if (ImGui::Checkbox(fmt::format("Correct Phase-center variation##antenna {} {}", id, i).c_str(), &_antenna.at(i).correctPCV))
653 {
654 LOG_DEBUG("{}: Antenna {} correctPCV: {}", id, i, _antenna.at(i).correctPCV);
655 changed = true;
656 }
657
658 if (!_antenna.at(i).correctPCO && !_antenna.at(i).correctPCV) { ImGui::BeginDisabled(); }
659 if (ImGui::Checkbox(fmt::format("Auto determine##antenna {} {}", id, i).c_str(), &_antenna.at(i).autoDetermine))
660 {
661 LOG_DEBUG("{}: Antenna {} auto: {}", id, i, _antenna.at(i).autoDetermine);
662 changed = true;
663 }
664 ImGui::SameLine();
665 if (_antenna.at(i).autoDetermine) { ImGui::BeginDisabled(); }
666 ImGui::SetNextItemWidth(itemWidth - (1.0F + static_cast<float>(_receiverCount > 1)) * ImGui::GetStyle().IndentSpacing);
667
668 if (ImGui::BeginCombo(fmt::format("Type##Combo {} {}", id, i).c_str(), _antenna.at(i).name.c_str()))
669 {
670 ImGui::PushFont(Application::MonoFont());
671 if (ImGui::Selectable(fmt::format("##Selectable {} {}", id, i).c_str(), _antenna.at(i).name.empty()))
672 {
673 _antenna.at(i).name = "";
674 }
675
676 for (const auto& iter : AntexReader::Get().antennas())
677 {
678 const bool is_selected = (_antenna.at(i).name == iter);
679 if (ImGui::Selectable(iter.c_str(), is_selected))
680 {
681 _antenna.at(i).name = iter;
682 }
683 if (is_selected) { ImGui::SetItemDefaultFocus(); }
684 }
685 ImGui::PopFont();
686 ImGui::EndCombo();
687 }
688 if (_antenna.at(i).autoDetermine) { ImGui::EndDisabled(); }
689 if (!_antenna.at(i).correctPCO && !_antenna.at(i).correctPCV) { ImGui::EndDisabled(); }
690 ImGui::SameLine();
691 gui::widgets::HelpMarker("Used for lookup in ANTEX file.");
692
693 ImGui::SetNextItemWidth(itemWidth - (1.0F + static_cast<float>(_receiverCount > 1)) * ImGui::GetStyle().IndentSpacing);
694 if (ImGui::InputDouble3(fmt::format("Delta H/E/N [m]##{} {}", id, i).c_str(), _antenna.at(i).hen_delta.data(), "%.3f"))
695 {
696 LOG_DEBUG("{}: Antenna {} delta: {}", id, i, _antenna.at(i).hen_delta.transpose());
697 changed = true;
698 }
699 ImGui::SameLine();
700 gui::widgets::HelpMarker("- Antenna height: Height of the antenna reference point (ARP) above the marker\n"
701 "- Horizontal eccentricity of ARP relative to the marker (east/north)");
702
703 if (_receiverCount > 1) { ImGui::TreePop(); }
704 }
705 ImGui::TreePop();
706 }
707
708 return changed;
709 }
710
712 [[nodiscard]] const IonosphereModel& getIonosphereModel() const { return _ionosphereModel; }
717
718 private:
719 IonosphereModel _ionosphereModel = IonosphereModel::Klobuchar;
722
724
726 struct Antenna
727 {
728 bool correctPCO = true;
729 bool correctPCV = true;
730 bool autoDetermine = true;
731 std::string name;
732 Eigen::Vector3d hen_delta = Eigen::Vector3d::Zero();
733 };
734
735 std::unordered_map<size_t, Antenna> _antenna;
736
740 friend void to_json(json& j, const ObservationEstimator::Antenna& obj)
741 {
742 j = json{
743 { "correctPCO", obj.correctPCO },
744 { "correctPCV", obj.correctPCV },
745 { "autoDetermine", obj.autoDetermine },
746 { "name", obj.name },
747 { "hen_delta", obj.hen_delta },
748 };
749 }
753 friend void from_json(const json& j, ObservationEstimator::Antenna& obj)
754 {
755 if (j.contains("correctPCO")) { j.at("correctPCO").get_to(obj.correctPCO); }
756 if (j.contains("correctPCV")) { j.at("correctPCV").get_to(obj.correctPCV); }
757 if (j.contains("autoDetermine")) { j.at("autoDetermine").get_to(obj.autoDetermine); }
758 if (j.contains("name")) { j.at("name").get_to(obj.name); }
759 if (j.contains("hen_delta")) { j.at("hen_delta").get_to(obj.hen_delta); }
760 }
761
765 friend void to_json(json& j, const ObservationEstimator& obj)
766 {
767 j = json{
768 { "ionosphereModel", Frequency_(obj._ionosphereModel) },
769 { "troposphereModels", obj._troposphereModels },
770 { "gnssMeasurementError", obj._gnssMeasurementErrorModel },
771 { "antenna", obj._antenna },
772 };
773 }
777 friend void from_json(const json& j, ObservationEstimator& obj)
778 {
779 if (j.contains("ionosphereModel")) { j.at("ionosphereModel").get_to(obj._ionosphereModel); }
780 if (j.contains("troposphereModels")) { j.at("troposphereModels").get_to(obj._troposphereModels); }
781 if (j.contains("gnssMeasurementError")) { j.at("gnssMeasurementError").get_to(obj._gnssMeasurementErrorModel); }
782 if (j.contains("antenna")) { j.at("antenna").get_to(obj._antenna); }
783 }
784};
785
786} // namespace NAV
GNSS Antenna related transformations.
ANTEX file reader.
Transformation collection.
nlohmann::json json
json namespace
Definition FlowManager.hpp:21
Frequency definition for different satellite systems.
Frequency_
Enumerate for GNSS frequencies.
Definition Frequency.hpp:26
Galileo Ephemeris information.
Derived::Scalar calcSatElevation(const Eigen::MatrixBase< Derived > &n_lineOfSightUnitVector)
Calculates the elevation of the satellite from the antenna.
Definition Functions.hpp:43
DerivedA::Scalar calcSagnacCorrection(const Eigen::MatrixBase< DerivedA > &e_posAnt, const Eigen::MatrixBase< DerivedB > &e_satPos)
Calculates the Earth rotation/Sagnac correction.
Definition Functions.hpp:66
DerivedA::Scalar calcSagnacRateCorrection(const Eigen::MatrixBase< DerivedA > &e_recvPos, const Eigen::MatrixBase< DerivedB > &e_satPos, const Eigen::MatrixBase< DerivedC > &e_recvVel, const Eigen::MatrixBase< DerivedD > &e_satVel)
Calculates the Range-rate Earth rotation/Sagnac correction.
Definition Functions.hpp:81
Derived::Scalar calcSatAzimuth(const Eigen::MatrixBase< Derived > &n_lineOfSightUnitVector)
Calculates the azimuth of the satellite from the antenna.
Definition Functions.hpp:54
Eigen::Vector3< typename DerivedA::Scalar > e_calcLineOfSightUnitVector(const Eigen::MatrixBase< DerivedA > &e_posAnt, const Eigen::MatrixBase< DerivedB > &e_posSat)
Calculates the line-of-sight unit vector from the user antenna to the satellite.
Definition Functions.hpp:31
GNSS Observation messages.
Text Help Marker (?) with Tooltip.
Ionosphere Models.
IonosphereModel
Available Ionosphere Models.
Definition Ionosphere.hpp:28
double ionoErrorVar(double dpsr_I, Frequency freq, int8_t num)
Calculates the ionospheric error variance.
bool ComboIonosphereModel(const char *label, IonosphereModel &ionosphereModel)
Shows a ComboBox to select the ionosphere model.
double calcIonosphericDelay(double tow, Frequency freq, int8_t freqNum, const Eigen::Vector3d &lla_pos, double elevation, double azimuth, IonosphereModel ionosphereModel=IonosphereModel::None, const IonosphericCorrections *corrections=nullptr)
Calculates the ionospheric delay.
Defines how to save certain datatypes to json.
#define LOG_DEBUG
Debug information. Should not be called on functions which receive observations (spamming)
Definition Logger.hpp:67
#define LOG_DATA
All output which occurs repeatedly every time observations are received.
Definition Logger.hpp:29
Errors concerning GNSS observations.
Observation data used for calculations.
const char * to_string(gui::widgets::PositionWithFrame::ReferenceFrame refFrame)
Converts the enum to a string.
Receiver information.
Structs identifying a unique satellite.
Troposphere Models.
bool ComboTroposphereModel(const char *label, TroposphereModelSelection &troposphereModelSelection, float width=0.0F)
Shows a ComboBox and button for advanced configuration to select the troposphere models.
double tropoErrorVar(double dpsr_T, double elevation)
Calculates the tropospheric error variance.
ZenithDelay calcTroposphericDelayAndMapping(const InsTime &insTime, const Eigen::Vector3d &lla_pos, double elevation, double azimuth, const TroposphereModelSelection &troposphereModels, const std::string &nameId)
Calculates the tropospheric zenith hydrostatic and wet delays and corresponding mapping factors.
Values with an uncertainty (Standard Deviation)
void move(std::vector< T > &v, size_t sourceIdx, size_t targetIdx)
Moves an element within a vector to a new position.
Definition Vector.hpp:26
static AntexReader & Get()
Get the static Instance of the reader.
Definition AntexReader.hpp:98
std::optional< double > getAntennaPhaseCenterVariation(const std::string &antennaType, Frequency_ freq, const InsTime &insTime, double elevation, std::optional< double > azimuth, const std::string &nameId) const
Gets the phase center variation for given elevation and azimuth.
Definition AntexReader.hpp:394
Frequency definition for different satellite systems.
Definition Frequency.hpp:59
SatelliteSystem getSatSys() const
Get the satellite system for which this frequency is defined.
Definition Frequency.hpp:152
Errors concerning GNSS observations.
Definition MeasurementErrors.hpp:31
double carrierMeasErrorVar(const SatelliteSystem &satSys, double elevation, double cn0) const
Calculates the measurement Error Variance for carrier-phase observations.
double codeBiasErrorVar() const
Returns the Code Bias Error Variance.
bool ShowGuiWidgets(const char *id, float width)
Shows a GUI widgets.
double psrRateMeasErrorVar(const Frequency &freq, int8_t num, double elevation, double cn0) const
Returns the Pseudo-range rate Error Variance.
double psrMeasErrorVar(const SatelliteSystem &satSys, double elevation, double cn0) const
Calculates the measurement Error Variance for pseudorange observations.
ObservationType
Observation types.
Definition GnssObs.hpp:37
@ Doppler
Doppler (Pseudorange rate)
Definition GnssObs.hpp:40
@ ObservationType_COUNT
Count.
Definition GnssObs.hpp:41
@ Carrier
Carrier-Phase.
Definition GnssObs.hpp:39
@ Pseudorange
Pseudorange.
Definition GnssObs.hpp:38
static constexpr double C
Speed of light [m/s].
Definition Constants.hpp:34
The class is responsible for all time-related tasks.
Definition InsTime.hpp:668
Calculates Observation estimates.
Definition ObservationEstimator.hpp:55
friend void to_json(json &j, const ObservationEstimator::Antenna &obj)
Converts the provided object into json.
Definition ObservationEstimator.hpp:740
const IonosphereModel & getIonosphereModel() const
Get the Ionosphere Model selected.
Definition ObservationEstimator.hpp:712
TroposphereModelSelection _troposphereModels
Troposphere Models used for the calculation.
Definition ObservationEstimator.hpp:720
friend void from_json(const json &j, ObservationEstimator &obj)
Converts the provided json object into a node object.
Definition ObservationEstimator.hpp:777
ObservationDifference
How the observation gets used. Influenced the measurement variance.
Definition ObservationEstimator.hpp:70
@ NoDifference
Estimation is not differenced.
Definition ObservationEstimator.hpp:71
@ DoubleDifference
Double Difference.
Definition ObservationEstimator.hpp:73
@ SingleDifference
Single Difference.
Definition ObservationEstimator.hpp:72
ObservationEstimator(size_t receiverCount)
Constructor.
Definition ObservationEstimator.hpp:59
GnssMeasurementErrorModel _gnssMeasurementErrorModel
GNSS measurement error model to use.
Definition ObservationEstimator.hpp:721
auto calcPseudorangeEstimate(Frequency freq, size_t receiverType, const Eigen::MatrixBase< Derived > &e_recvPosMarker, auto recvClockBias, auto interFreqBias, auto dpsr_T_r_s, auto dpsr_I_r_s, const std::shared_ptr< const SatelliteInfo > &satInfo, const std::shared_ptr< const GnssObs > &gnssObs, const std::string &nameId) const
Calculates the pseudorange estimate.
Definition ObservationEstimator.hpp:279
double calcPCV(const Frequency &freq, const InsTime &insTime, double elevation, std::optional< double > azimuth, size_t receiverType, const std::shared_ptr< const GnssObs > &gnssObs, const std::string &nameId) const
Calculates the phase center variation.
Definition ObservationEstimator.hpp:525
Eigen::Vector3< typename Derived::Scalar > e_calcRecvPosAPC(const Frequency &freq, size_t receiverType, const Eigen::MatrixBase< Derived > &e_recvPosMarker, const std::shared_ptr< const GnssObs > &gnssObs, const std::string &nameId) const
Calculates the antenna phase center position for the marker.
Definition ObservationEstimator.hpp:488
double calculateIonosphericDelay(const Frequency &freq, int8_t freqNum, const Eigen::Vector3d &lla_recvPosAPC, const std::shared_ptr< const GnssObs > &gnssObs, const std::shared_ptr< const IonosphericCorrections > &ionosphericCorrections, double elevation, double azimuth) const
Calculate the ionospheric delay with the model selected in the GUI.
Definition ObservationEstimator.hpp:568
const GnssMeasurementErrorModel & getGnssMeasurementErrorModel() const
Get the GNSS measurement error model.
Definition ObservationEstimator.hpp:716
size_t _receiverCount
Number of receivers (used for GUI)
Definition ObservationEstimator.hpp:723
auto calcDopplerEstimate(Frequency freq, size_t receiverType, const Eigen::MatrixBase< DerivedPos > &e_recvPosMarker, const Eigen::MatrixBase< DerivedVel > &e_recvVel, auto recvClockDrift, const std::shared_ptr< const SatelliteInfo > &satInfo, const std::shared_ptr< const GnssObs > &gnssObs, const std::string &nameId) const
Calculates the doppler estimate.
Definition ObservationEstimator.hpp:387
IonosphereModel _ionosphereModel
Ionosphere Model used for the calculation.
Definition ObservationEstimator.hpp:719
double calcObservationVariance(Frequency freq, int8_t freqNum, size_t receiverType, GnssObs::ObservationType obsType, const Eigen::Vector3d &e_recvPosMarker, const Eigen::Vector3d &e_satPos, double cn0, const std::shared_ptr< const GnssObs > &gnssObs, const std::shared_ptr< const SatNavData > &navData, const std::shared_ptr< const IonosphericCorrections > &ionosphericCorrections, const std::string &nameId, ObservationDifference obsDiff) const
Calculates the observation variance.
Definition ObservationEstimator.hpp:419
void calcObservationEstimate(const SatSigId &satSigId, Observations::SignalObservation &observation, const Receiver< ReceiverType > &receiver, const std::shared_ptr< const IonosphericCorrections > &ionosphericCorrections, const std::string &nameId, ObservationDifference obsDiff) const
Calculates the observation estimate for the given signal.
Definition ObservationEstimator.hpp:104
double calculateTroposphericDelay(const Eigen::Vector3d &lla_recvPosAPC, const std::shared_ptr< const GnssObs > &gnssObs, double elevation, double azimuth, const std::string &nameId) const
Calculate the tropospheric delay with the model selected in the GUI.
Definition ObservationEstimator.hpp:589
std::pair< double, Eigen::Vector3d > calcCarrierEstimate(Frequency freq, int8_t freqNum, size_t receiverType, const Eigen::MatrixBase< Derived > &e_recvPosMarker, double recvClockBias, const std::shared_ptr< const SatelliteInfo > &satInfo, const std::shared_ptr< const GnssObs > &gnssObs, const std::shared_ptr< const IonosphericCorrections > &ionosphericCorrections, const std::string &nameId) const
Calculates the carrier-phase estimate.
Definition ObservationEstimator.hpp:332
bool ShowGuiWidgets(const char *id, float itemWidth)
Shows the GUI input to select the options.
Definition ObservationEstimator.hpp:605
const TroposphereModelSelection & getTroposphereModels() const
Get the Troposphere Model selection.
Definition ObservationEstimator.hpp:714
friend void from_json(const json &j, ObservationEstimator::Antenna &obj)
Converts the provided json object into a node object.
Definition ObservationEstimator.hpp:753
std::unordered_map< size_t, Antenna > _antenna
User antenna selection. Key is receiver type.
Definition ObservationEstimator.hpp:735
friend void to_json(json &j, const ObservationEstimator &obj)
Converts the provided object into json.
Definition ObservationEstimator.hpp:765
void calcObservationEstimates(Observations &observations, const Receiver< ReceiverType > &receiver, const std::shared_ptr< const IonosphericCorrections > &ionosphericCorrections, const std::string &nameId, ObservationDifference obsDiff) const
Calculates the observation estimates.
Definition ObservationEstimator.hpp:83
ImGui extensions.
Antenna information.
Definition ObservationEstimator.hpp:727
bool autoDetermine
Try determine antenna type from e.g. RINEX file.
Definition ObservationEstimator.hpp:730
Eigen::Vector3d hen_delta
Delta values of marker to antenna base.
Definition ObservationEstimator.hpp:732
bool correctPCV
Correct phase center variation.
Definition ObservationEstimator.hpp:729
std::string name
Type of the antenna for the Antex lookup.
Definition ObservationEstimator.hpp:731
bool correctPCO
Correct phase center offset.
Definition ObservationEstimator.hpp:728
Satellite Info used for estimate calculation.
Definition ObservationEstimator.hpp:251
Eigen::Vector3d e_satVel
Satellite velocity in ECEF frame [m/s].
Definition ObservationEstimator.hpp:261
SatelliteInfo(Eigen::Vector3d e_satPos, Eigen::Vector3d e_satVel, double satClockBias, double satClockDrift)
Constructor.
Definition ObservationEstimator.hpp:257
Eigen::Vector3d e_satPos
Satellite position in ECEF coordinates [m].
Definition ObservationEstimator.hpp:260
double satClockBias
Satellite clock bias [s].
Definition ObservationEstimator.hpp:262
double satClockDrift
Satellite clock drift [s/s].
Definition ObservationEstimator.hpp:263
Receiver specific observation of the signal.
Definition Observation.hpp:44
int8_t freqNum() const
Frequency number. Only used for GLONASS G1 and G2.
Definition Observation.hpp:198
std::unordered_map< size_t, std::shared_ptr< ReceiverSpecificData > > recvObs
Receiver specific data.
Definition Observation.hpp:193
const std::shared_ptr< const SatNavData > & navData() const
Satellite Navigation data.
Definition Observation.hpp:196
Observation storage type.
Definition Observation.hpp:41
unordered_map< SatSigId, SignalObservation > signals
Observations and calculated data for each signal.
Definition Observation.hpp:205
const double * biasFor(SatelliteSystem satSys) const
Get the bias for the given satellite system.
Definition ReceiverClock.hpp:81
const double * driftFor(SatelliteSystem satSys) const
Get the drift for the given satellite system.
Definition ReceiverClock.hpp:115
Receiver information.
Definition Receiver.hpp:34
std::shared_ptr< const GnssObs > gnssObs
Latest GNSS observation.
Definition Receiver.hpp:54
std::unordered_map< Frequency, UncertainValue< double > > interFrequencyBias
Inter frequency biases [s].
Definition Receiver.hpp:52
ReceiverType type
Receiver Type.
Definition Receiver.hpp:42
Eigen::Vector3d e_posMarker
Marker Position in ECEF frame [m].
Definition Receiver.hpp:44
Eigen::Vector3d e_vel
Velocity in ECEF frame [m/s].
Definition Receiver.hpp:48
ReceiverClock recvClk
Estimated receiver clock parameters.
Definition Receiver.hpp:50
Identifies a satellite signal (satellite frequency and number)
Definition SatelliteIdentifier.hpp:67
Frequency freq() const
Returns the frequency of the satellite signal.
Definition SatelliteIdentifier.hpp:108
Satellite System type.
Definition SatelliteSystem.hpp:44
Collection of troposphere model selections.
Definition Troposphere.hpp:63