0.2.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 <cstddef>
17#include <utility>
18#include <vector>
19
20#include <imgui.h>
21#include <application.h>
22#include <imgui_stdlib.h>
27
34
35#include "util/Json.hpp"
36#include <fmt/core.h>
37#include <fmt/format.h>
38
39namespace NAV
40{
41
44{
45 public:
48 explicit ObservationEstimator(size_t receiverCount)
49 : _antenna(receiverCount) {}
50
58
65 template<typename ReceiverType>
67 const std::array<Receiver<ReceiverType>, ReceiverType::ReceiverType_COUNT>& receivers,
68 const IonosphericCorrections& ionosphericCorrections,
69 [[maybe_unused]] const std::string& nameId,
71 {
72 LOG_DATA("{}: Calculating observation estimates:", nameId);
73
74 for (auto& [satSigId, observation] : observations.signals)
75 {
76 const Frequency freq = satSigId.freq();
77 const SatelliteSystem satSys = freq.getSatSys();
78
79 for (size_t r = 0; r < observation.recvObs.size(); r++)
80 {
81 auto& recvObs = observation.recvObs.at(r);
82 [[maybe_unused]] auto recv = static_cast<ReceiverType>(r);
83 const Receiver<ReceiverType>& receiver = receivers.at(r);
84 const Antenna& antenna = _antenna.at(r);
85
86 Eigen::Vector3d hen_delta(0.0, 0.0, antenna.hen_delta(0));
87 Eigen::Vector3d e_recvPosAPC;
88 Eigen::Vector3d lla_recvPosAPC;
89 if (antenna.enabled)
90 {
91 std::string antennaType;
92 if (antenna.autoDetermine && receiver.gnssObs->receiverInfo)
93 {
94 antennaType = receiver.gnssObs->receiverInfo->get().antennaType;
95 }
96 else if (!antenna.autoDetermine)
97 {
98 antennaType = antenna.name;
99 }
100 e_recvPosAPC = receiver.e_posAntennaPhaseCenter(freq, antennaType, nameId);
101 lla_recvPosAPC = receiver.lla_posAntennaPhaseCenter(freq, antennaType, nameId);
102 }
103 else
104 {
105 e_recvPosAPC = receiver.e_posARP();
106 lla_recvPosAPC = receiver.lla_posARP();
107 }
108 e_recvPosAPC += trafo::e_Quat_n(receiver.lla_posMarker(0), receiver.lla_posMarker(1)) * hen_delta;
109 lla_recvPosAPC.z() += hen_delta.z();
110
111 // Receiver-Satellite Range [m]
112 double rho_r_s = (recvObs.e_satPos() - e_recvPosAPC).norm();
113 recvObs.terms.rho_r_s = rho_r_s;
114 // Troposphere
115 auto tropo_r_s = calcTroposphericDelayAndMapping(receiver.gnssObs->insTime, lla_recvPosAPC,
116 recvObs.satElevation(), recvObs.satAzimuth(), _troposphereModels);
117 recvObs.terms.tropoZenithDelay = tropo_r_s;
118 // Estimated troposphere propagation error [m]
119 double dpsr_T_r_s = tropo_r_s.ZHD * tropo_r_s.zhdMappingFactor + tropo_r_s.ZWD * tropo_r_s.zwdMappingFactor;
120 recvObs.terms.dpsr_T_r_s = dpsr_T_r_s;
121 // Estimated ionosphere propagation error [m]
122 double dpsr_I_r_s = calcIonosphericDelay(static_cast<double>(receiver.gnssObs->insTime.toGPSweekTow().tow),
123 freq, observation.freqNum(), lla_recvPosAPC, recvObs.satElevation(), recvObs.satAzimuth(),
124 _ionosphereModel, &ionosphericCorrections);
125 recvObs.terms.dpsr_I_r_s = dpsr_I_r_s;
126 // Sagnac correction [m]
127 double dpsr_ie_r_s = calcSagnacCorrection(e_recvPosAPC, recvObs.e_satPos());
128 recvObs.terms.dpsr_ie_r_s = dpsr_ie_r_s;
129
130 // Earth's gravitational field causes relativistic signal delay due to space-time curvature (Shapiro effect) [s]
131 // double posNorm = recvObs.e_satPos().norm() + e_recvPosAPC.norm();
132 // 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))
133 // : 0.0;
134
135 double cn0 = recvObs.gnssObsData().CN0.value_or(1.0);
136
137 for (auto& [obsType, obsData] : recvObs.obs)
138 {
139 LOG_DATA("{}: [{}][{:11}][{:5}] Observation estimate", nameId, satSigId, obsType, recv);
140 switch (obsType)
141 {
143 obsData.estimate = rho_r_s
144 + dpsr_ie_r_s
145 + dpsr_T_r_s
146 + dpsr_I_r_s
148 * (receiver.recvClk.bias.value * (obsDiff != DoubleDifference)
149 - recvObs.satClock().bias * (obsDiff == NoDifference)
150 + receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).value
151 // + dt_rel_stc
152 + (receiver.interFrequencyBias.contains(freq) ? receiver.interFrequencyBias.at(freq).value : 0.0));
153 obsData.measVar = _gnssMeasurementErrorModel.psrMeasErrorVar(satSys, recvObs.satElevation(), cn0);
154 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m] Geometrical range", nameId, satSigId, obsType, recv, rho_r_s);
155 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Sagnac correction", nameId, satSigId, obsType, recv, dpsr_ie_r_s);
156 if (dpsr_T_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Tropospheric delay", nameId, satSigId, obsType, recv, dpsr_T_r_s); }
157 if (dpsr_I_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Ionospheric delay", nameId, satSigId, obsType, recv, dpsr_I_r_s); }
158 if (obsDiff != DoubleDifference) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Receiver clock bias", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.recvClk.bias.value); }
159 if (obsDiff == NoDifference) { LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m] Satellite clock bias", nameId, satSigId, obsType, recv, InsConst<>::C * recvObs.satClock().bias); }
160 if (receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).value != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Inter-system clock bias", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).value); }
161 // LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Shapiro effect", nameId, satSigId, obsType, recv, InsConst<>::C * dt_rel_stc);
162 if (receiver.interFrequencyBias.contains(freq)) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Inter-frequency bias", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.interFrequencyBias.at(freq).value); }
163 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4f} [m] Pseudorange estimate", nameId, satSigId, obsType, recv, obsData.estimate);
164 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4e} [m] Difference to measurement", nameId, satSigId, obsType, recv, obsData.measurement - obsData.estimate);
165 break;
166 case GnssObs::Carrier:
167 obsData.estimate = rho_r_s
168 + dpsr_ie_r_s
169 + dpsr_T_r_s
170 - dpsr_I_r_s
172 * (receiver.recvClk.bias.value * (obsDiff != DoubleDifference)
173 - recvObs.satClock().bias * (obsDiff == NoDifference)
174 + receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).value
175 // + dt_rel_stc
176 );
177 obsData.measVar = _gnssMeasurementErrorModel.carrierMeasErrorVar(satSys, recvObs.satElevation(), cn0);
178 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m] Geometrical range", nameId, satSigId, obsType, recv, rho_r_s);
179 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Sagnac correction", nameId, satSigId, obsType, recv, dpsr_ie_r_s);
180 if (dpsr_T_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Tropospheric delay", nameId, satSigId, obsType, recv, dpsr_T_r_s); }
181 if (dpsr_I_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m] Ionospheric delay", nameId, satSigId, obsType, recv, dpsr_I_r_s); }
182 if (obsDiff != DoubleDifference) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Receiver clock bias", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.recvClk.bias.value); }
183 if (obsDiff == NoDifference) { LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m] Satellite clock bias", nameId, satSigId, obsType, recv, InsConst<>::C * recvObs.satClock().bias); }
184 if (receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).value != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Inter-system clock bias", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).value); }
185 // LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m] Shapiro effect", nameId, satSigId, obsType, recv, InsConst<>::C * dt_rel_stc);
186 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4f} [m] Carrier-phase estimate", nameId, satSigId, obsType, recv, obsData.estimate);
187 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4e} [m] Difference to measurement", nameId, satSigId, obsType, recv, obsData.measurement - obsData.estimate);
188 break;
189 case GnssObs::Doppler:
190 obsData.estimate = recvObs.e_pLOS().dot(recvObs.e_satVel() - receiver.e_vel)
191 - calcSagnacRateCorrection(e_recvPosAPC, recvObs.e_satPos(), receiver.e_vel, recvObs.e_satVel())
193 * (receiver.recvClk.drift.value * (obsDiff != DoubleDifference)
194 - recvObs.satClock().drift * (obsDiff == NoDifference)
195 + receiver.recvClk.sysTimeDiffDrift.at(satSys.toEnumeration()).value * (obsDiff == NoDifference));
196 obsData.measVar = _gnssMeasurementErrorModel.psrRateMeasErrorVar(freq, observation.freqNum(), recvObs.satElevation(), cn0);
197 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4f} [m/s] ", nameId, satSigId, obsType, recv, recvObs.e_pLOS().dot(recvObs.e_satVel() - receiver.e_vel));
198 LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m/s] Sagnac rate correction", nameId, satSigId, obsType, recv, calcSagnacRateCorrection(e_recvPosAPC, recvObs.e_satPos(), receiver.e_vel, recvObs.e_satVel()));
199 if (obsDiff != DoubleDifference) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m/s] Receiver clock drift", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.recvClk.drift.value); }
200 if (obsDiff == NoDifference) { LOG_DATA("{}: [{}][{:11}][{:5}] - {:.4f} [m/s] Satellite clock drift", nameId, satSigId, obsType, recv, InsConst<>::C * recvObs.satClock().drift); }
201 if (receiver.recvClk.sysTimeDiffDrift.at(satSys.toEnumeration()).value != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m/s] Inter-system clock drift", nameId, satSigId, obsType, recv, InsConst<>::C * receiver.recvClk.sysTimeDiffDrift.at(satSys.toEnumeration()).value); }
202 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4f} [m/s] Doppler estimate", nameId, satSigId, obsType, recv, obsData.estimate);
203 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4e} [m/s] Difference to measurement", nameId, satSigId, obsType, recv, obsData.measurement - obsData.estimate);
204 break;
206 break;
207 }
208 LOG_DATA("{}: [{}][{:11}][{:5}] Observation error variance", nameId, satSigId, obsType, recv);
209 LOG_DATA("{}: [{}][{:11}][{:5}] {:.4g} [{}] Measurement error variance", nameId, satSigId, obsType, recv, obsData.measVar,
210 obsType == GnssObs::Doppler ? "m^2/s^2" : "m^2");
211
212 if (obsDiff == NoDifference)
213 {
214 if (obsType == GnssObs::Pseudorange || obsType == GnssObs::Carrier)
215 {
216 obsData.measVar += observation.navData()->calcSatellitePositionVariance()
217 + ionoErrorVar(dpsr_I_r_s, freq, observation.freqNum())
218 + tropoErrorVar(dpsr_T_r_s, recvObs.satElevation());
219 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Satellite position variance", nameId, satSigId, obsType, recv, observation.navData()->calcSatellitePositionVariance());
220 if (dpsr_I_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Ionosphere variance", nameId, satSigId, obsType, recv, ionoErrorVar(dpsr_I_r_s, freq, observation.freqNum())); }
221 if (dpsr_T_r_s != 0.0) { LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Troposphere variance", nameId, satSigId, obsType, recv, tropoErrorVar(dpsr_T_r_s, recvObs.satElevation())); }
222 }
223 if (obsType == GnssObs::Pseudorange)
224 {
225 obsData.measVar += _gnssMeasurementErrorModel.codeBiasErrorVar();
226 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Code bias variance", nameId, satSigId, obsType, recv, _gnssMeasurementErrorModel.codeBiasErrorVar());
227
228 if (receiver.interFrequencyBias.contains(freq))
229 {
230 obsData.measVar += std::pow(receiver.interFrequencyBias.at(freq).stdDev, 2);
231 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Inter-frequency bias variance", nameId, satSigId, obsType, recv,
232 std::pow(InsConst<>::C * receiver.interFrequencyBias.at(freq).stdDev, 2));
233 }
234 }
235 }
236 if (obsDiff != DoubleDifference)
237 {
238 if (obsType == GnssObs::Pseudorange || obsType == GnssObs::Carrier)
239 {
240 double recvClockVariance = std::pow(InsConst<>::C, 2)
241 * (std::pow(receiver.recvClk.bias.stdDev, 2)
242 + std::pow(receiver.recvClk.sysTimeDiffBias.at(satSys.toEnumeration()).stdDev, 2));
243 obsData.measVar += recvClockVariance;
244 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2] Receiver clock bias variance", nameId, satSigId, obsType, recv, recvClockVariance);
245 }
246 else if (obsType == GnssObs::Doppler)
247 {
248 double recvClockVariance = std::pow(InsConst<>::C, 2)
249 * (std::pow(receiver.recvClk.drift.stdDev, 2)
250 + std::pow(receiver.recvClk.sysTimeDiffDrift.at(satSys.toEnumeration()).stdDev, 2));
251 obsData.measVar += recvClockVariance;
252 LOG_DATA("{}: [{}][{:11}][{:5}] + {:.4f} [m^2/s^2] Receiver clock drift variance", nameId, satSigId, obsType, recv, recvClockVariance);
253 }
254 }
255 LOG_DATA("{}: [{}][{:11}][{:5}] = {:.4g} [{}] Observation error variance", nameId, satSigId, obsType, recv, obsData.measVar,
256 obsType == GnssObs::Doppler ? "m^2/s^2" : "m^2");
257 }
258 }
259 }
260 }
261
265 template<typename ReceiverType>
266 bool ShowGuiWidgets(const char* id, float itemWidth)
267 {
268 bool changed = false;
269
270 ImGui::SetNextItemOpen(true, ImGuiCond_FirstUseEver);
271 if (ImGui::TreeNode(fmt::format("Compensation models##{}", id).c_str()))
272 {
273 ImGui::SetNextItemWidth(itemWidth - ImGui::GetStyle().IndentSpacing);
274 if (ComboIonosphereModel(fmt::format("Ionosphere Model##{}", id).c_str(), _ionosphereModel))
275 {
276 LOG_DEBUG("{}: Ionosphere Model changed to {}", id, NAV::to_string(_ionosphereModel));
277 changed = true;
278 }
279 if (ComboTroposphereModel(fmt::format("Troposphere Model##{}", id).c_str(), _troposphereModels, itemWidth - ImGui::GetStyle().IndentSpacing))
280 {
281 changed = true;
282 }
283 ImGui::TreePop();
284 }
285
286 ImGui::SetNextItemOpen(true, ImGuiCond_FirstUseEver);
287 if (ImGui::TreeNode(fmt::format("GNSS Measurement Error##{}", id).c_str()))
288 {
289 if (_gnssMeasurementErrorModel.ShowGuiWidgets(id, itemWidth - ImGui::GetStyle().IndentSpacing))
290 {
291 LOG_DEBUG("{}: GNSS Measurement Error Model changed.", id);
292 changed = true;
293 }
294 ImGui::TreePop();
295 }
296
297 if (ImGui::TreeNode(fmt::format("Antenna settings##{}", id).c_str()))
298 {
299 for (size_t i = 0; i < ReceiverType::ReceiverType_COUNT; ++i)
300 {
301 if (ReceiverType::ReceiverType_COUNT > 1)
302 {
303 ImGui::SetNextItemOpen(true, ImGuiCond_FirstUseEver);
304 if (!ImGui::TreeNode(fmt::format("{}", static_cast<ReceiverType>(i)).c_str())) { continue; }
305 }
306
307 if (ImGui::Checkbox(fmt::format("Correct Phase-center offset##antenna {} {}", id, i).c_str(), &_antenna.at(i).enabled))
308 {
309 LOG_DEBUG("{}: Antenna {} enabled: {}", id, i, _antenna.at(i).enabled);
310 changed = true;
311 }
312 ImGui::SameLine();
313 if (!_antenna.at(i).enabled) { ImGui::BeginDisabled(); }
314 if (ImGui::Checkbox(fmt::format("Auto determine##antenna {} {}", id, i).c_str(), &_antenna.at(i).autoDetermine))
315 {
316 LOG_DEBUG("{}: Antenna {} auto: {}", id, i, _antenna.at(i).autoDetermine);
317 changed = true;
318 }
319 if (_antenna.at(i).autoDetermine) { ImGui::BeginDisabled(); }
320 ImGui::SetNextItemWidth(itemWidth - (1.0F + static_cast<float>(ReceiverType::ReceiverType_COUNT > 1)) * ImGui::GetStyle().IndentSpacing);
321
322 if (ImGui::BeginCombo(fmt::format("Type##Combo {} {}", id, i).c_str(), _antenna.at(i).name.c_str()))
323 {
324 ImGui::PushFont(Application::MonoFont());
325 if (ImGui::Selectable(fmt::format("##Selectable {} {}", id, i).c_str(), _antenna.at(i).name.empty()))
326 {
327 _antenna.at(i).name = "";
328 }
329
330 for (const auto& iter : AntexReader::Get().antennas())
331 {
332 const bool is_selected = (_antenna.at(i).name == iter);
333 if (ImGui::Selectable(iter.c_str(), is_selected))
334 {
335 _antenna.at(i).name = iter;
336 }
337 if (is_selected) { ImGui::SetItemDefaultFocus(); }
338 }
339 ImGui::PopFont();
340 ImGui::EndCombo();
341 }
342 if (_antenna.at(i).autoDetermine) { ImGui::EndDisabled(); }
343 if (!_antenna.at(i).enabled) { ImGui::EndDisabled(); }
344 ImGui::SameLine();
345 gui::widgets::HelpMarker("Used for lookup in ANTEX file.");
346
347 ImGui::SetNextItemWidth(itemWidth - (1.0F + static_cast<float>(ReceiverType::ReceiverType_COUNT > 1)) * ImGui::GetStyle().IndentSpacing);
348 if (ImGui::InputDouble3(fmt::format("Delta H/E/N [m]##{} {}", id, i).c_str(), _antenna.at(i).hen_delta.data(), "%.3f"))
349 {
350 LOG_DEBUG("{}: Antenna {} delta: {}", id, i, _antenna.at(i).hen_delta.transpose());
351 changed = true;
352 }
353 ImGui::SameLine();
354 gui::widgets::HelpMarker("- Antenna height: Height of the antenna reference point (ARP) above the marker\n"
355 "- Horizontal eccentricity of ARP relative to the marker (east/north)");
356
357 if (ReceiverType::ReceiverType_COUNT > 1) { ImGui::TreePop(); }
358 }
359 ImGui::TreePop();
360 }
361
362 return changed;
363 }
364
365 private:
366 IonosphereModel _ionosphereModel = IonosphereModel::Klobuchar;
367 TroposphereModelSelection _troposphereModels;
368 GnssMeasurementErrorModel _gnssMeasurementErrorModel;
369
371 struct Antenna
372 {
373 bool enabled = true;
374 bool autoDetermine = true;
375 std::string name;
376 Eigen::Vector3d hen_delta;
377 };
378
379 std::vector<Antenna> _antenna;
380
384 friend void to_json(json& j, const ObservationEstimator::Antenna& obj)
385 {
386 j = json{
387 { "enabled", obj.enabled },
388 { "autoDetermine", obj.autoDetermine },
389 { "name", obj.name },
390 { "hen_delta", obj.hen_delta },
391 };
392 }
396 friend void from_json(const json& j, ObservationEstimator::Antenna& obj)
397 {
398 if (j.contains("enabled")) { j.at("enabled").get_to(obj.enabled); }
399 if (j.contains("autoDetermine")) { j.at("autoDetermine").get_to(obj.autoDetermine); }
400 if (j.contains("name")) { j.at("name").get_to(obj.name); }
401 if (j.contains("hen_delta")) { j.at("hen_delta").get_to(obj.hen_delta); }
402 }
403
407 friend void to_json(json& j, const ObservationEstimator& obj)
408 {
409 j = json{
410 { "ionosphereModel", Frequency_(obj._ionosphereModel) },
411 { "troposphereModels", obj._troposphereModels },
412 { "gnssMeasurementError", obj._gnssMeasurementErrorModel },
413 { "antenna", obj._antenna },
414 };
415 }
419 friend void from_json(const json& j, ObservationEstimator& obj)
420 {
421 if (j.contains("ionosphereModel")) { j.at("ionosphereModel").get_to(obj._ionosphereModel); }
422 if (j.contains("troposphereModels")) { j.at("troposphereModels").get_to(obj._troposphereModels); }
423 if (j.contains("gnssMeasurementError")) { j.at("gnssMeasurementError").get_to(obj._gnssMeasurementErrorModel); }
424 if (j.contains("antenna")) { j.at("antenna").get_to(obj._antenna); }
425 }
426};
427
428} // namespace NAV
ANTEX file reader.
Transformation collection.
nlohmann::json json
json namespace
Definition FlowManager.hpp:21
Frequency_
Enumerate for GNSS frequencies.
Definition Frequency.hpp:26
Galileo Ephemeris information.
double calcSagnacRateCorrection(const Eigen::Vector3d &e_posAnt, const Eigen::Vector3d &e_satPos, const Eigen::Vector3d &e_velAnt, const Eigen::Vector3d &e_satVel)
Calculates the Range-rate Earth rotation/Sagnac correction.
double calcSagnacCorrection(const Eigen::Vector3d &e_posAnt, const Eigen::Vector3d &e_satPos)
Calculates the Earth rotation/Sagnac correction.
Text Help Marker (?) with Tooltip.
Ionosphere Models.
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.
IonosphereModel
Available Ionosphere Models.
Definition Ionosphere.hpp:26
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.
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.
ZenithDelay calcTroposphericDelayAndMapping(const InsTime &insTime, const Eigen::Vector3d &lla_pos, double elevation, double azimuth, const TroposphereModelSelection &troposphereModels)
Calculates the tropospheric zenith hydrostatic and wet delays and corresponding mapping factors.
double tropoErrorVar(double dpsr_T, double elevation)
Calculates the tropospheric error variance.
const std::set< std::string > & antennas() const
Antennas read from the ANTEX files.
Definition AntexReader.hpp:294
static AntexReader & Get()
Get the static Instance of the reader.
Definition AntexReader.hpp:80
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.
@ Doppler
Doppler (Pseudorange rate)
Definition GnssObs.hpp:39
@ ObservationType_COUNT
Count.
Definition GnssObs.hpp:40
@ Carrier
Carrier-Phase.
Definition GnssObs.hpp:38
@ Pseudorange
Pseudorange.
Definition GnssObs.hpp:37
Constants.
Definition Constants.hpp:26
Ionospheric Corrections.
Definition IonosphericCorrections.hpp:30
Calculates Observation estimates.
Definition ObservationEstimator.hpp:44
friend void to_json(json &j, const ObservationEstimator::Antenna &obj)
Converts the provided object into json.
Definition ObservationEstimator.hpp:384
friend void from_json(const json &j, ObservationEstimator &obj)
Converts the provided json object into a node object.
Definition ObservationEstimator.hpp:419
ObservationEstimator(size_t receiverCount)
Constructor.
Definition ObservationEstimator.hpp:48
void calcObservationEstimates(Observations &observations, const std::array< Receiver< ReceiverType >, ReceiverType::ReceiverType_COUNT > &receivers, const IonosphericCorrections &ionosphericCorrections, const std::string &nameId, ObservationDifference obsDiff)
Calculates the observation estimates.
Definition ObservationEstimator.hpp:66
ObservationDifference
How the observation gets used. Influenced the measurement variance.
Definition ObservationEstimator.hpp:53
@ NoDifference
Estimation is not differenced.
Definition ObservationEstimator.hpp:54
@ DoubleDifference
Double Difference.
Definition ObservationEstimator.hpp:56
@ SingleDifference
Single Difference.
Definition ObservationEstimator.hpp:55
bool ShowGuiWidgets(const char *id, float itemWidth)
Shows the GUI input to select the options.
Definition ObservationEstimator.hpp:266
friend void from_json(const json &j, ObservationEstimator::Antenna &obj)
Converts the provided json object into a node object.
Definition ObservationEstimator.hpp:396
friend void to_json(json &j, const ObservationEstimator &obj)
Converts the provided object into json.
Definition ObservationEstimator.hpp:407
ImGui extensions.
Observation storage type.
Definition Observation.hpp:38
unordered_map< SatSigId, SignalObservation > signals
Observations and calculated data for each signal.
Definition Observation.hpp:148
UncertainValue< double > drift
Estimated receiver clock drift [s/s].
Definition ReceiverClock.hpp:30
std::array< UncertainValue< double >, SatelliteSystem::Enum::Enum_COUNT > sysTimeDiffDrift
System time difference drift [s/s].
Definition ReceiverClock.hpp:37
std::array< UncertainValue< double >, SatelliteSystem::Enum::Enum_COUNT > sysTimeDiffBias
System time difference bias [s].
Definition ReceiverClock.hpp:35
UncertainValue< double > bias
Estimated receiver clock bias [s].
Definition ReceiverClock.hpp:28
Receiver information.
Definition Receiver.hpp:35
std::shared_ptr< const GnssObs > gnssObs
Latest GNSS observation.
Definition Receiver.hpp:53
Eigen::Vector3d e_posARP() const
Antenna Reference Point position in ECEF frame [m] (Marker + antennaDeltaNEU)
Definition Receiver.hpp:56
Eigen::Vector3d lla_posMarker
Marker Position in LLA frame [rad, rad, m].
Definition Receiver.hpp:45
Eigen::Vector3d lla_posAntennaPhaseCenter(Frequency freq, const std::string &antennaType, const std::string &nameId) const
Marker position in LLA frame [rad, rad, m] (ARP + antenna phase center)
Definition Receiver.hpp:108
Eigen::Vector3d lla_posARP() const
Antenna Reference Point position in LLA frame [rad, rad, m] (Marker + antennaDeltaNEU)
Definition Receiver.hpp:70
std::unordered_map< Frequency, UncertainValue< double > > interFrequencyBias
Inter frequency biases.
Definition Receiver.hpp:51
Eigen::Vector3d e_vel
Velocity in ECEF frame [m/s].
Definition Receiver.hpp:47
Eigen::Vector3d e_posAntennaPhaseCenter(Frequency freq, const std::string &antennaType, const std::string &nameId) const
Marker position in ECEF frame [m] (ARP + antenna phase center)
Definition Receiver.hpp:84
ReceiverClock recvClk
Estimated receiver clock parameters.
Definition Receiver.hpp:49
Satellite System type.
Definition SatelliteSystem.hpp:43
Enum toEnumeration() const
Returns a continuous enumeration of the object.
Collection of troposphere model selections.
Definition Troposphere.hpp:61
T stdDev
Standard deviation.
Definition UncertainValue.hpp:25
T value
Value.
Definition UncertainValue.hpp:24