0.2.0
Loading...
Searching...
No Matches
Mechanization.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 <Eigen/Core>
17#include <Eigen/Dense>
18
27#include "Navigation/Transformations/Units.hpp"
28#include "util/Logger.hpp"
29
30#include <cmath>
31
32namespace NAV
33{
51template<typename DerivedA, typename DerivedB>
52Eigen::Vector4<typename DerivedA::Scalar> calcTimeDerivativeFor_n_Quat_b(const Eigen::MatrixBase<DerivedA>& b_omega_nb,
53 const Eigen::MatrixBase<DerivedB>& n_Quat_b_coeffs)
54{
55 // Angular rates in matrix form (Titterton (2005), eq. (11.35))
56 Eigen::Matrix4<typename DerivedA::Scalar> A;
57
58 // clang-format off
59 A << 0.0 , -b_omega_nb(0), -b_omega_nb(1), -b_omega_nb(2),
60 b_omega_nb(0), 0.0 , b_omega_nb(2), -b_omega_nb(1),
61 b_omega_nb(1), -b_omega_nb(2), 0.0 , b_omega_nb(0),
62 b_omega_nb(2), b_omega_nb(1), -b_omega_nb(0), 0.0 ;
63 // clang-format on
64
65 // Propagation of an attitude Quaternion with time (Titterton, ch. 11.2.5, eq. 11.33-11.35, p. 319)
66 return 0.5 * A * n_Quat_b_coeffs; // (w, x, y, z)
67}
68
86template<typename DerivedA, typename DerivedB, typename DerivedC, typename DerivedD>
87Eigen::Vector3<typename DerivedA::Scalar> n_calcTimeDerivativeForVelocity(const Eigen::MatrixBase<DerivedA>& n_measuredForce,
88 const Eigen::MatrixBase<DerivedB>& n_coriolisAcceleration,
89 const Eigen::MatrixBase<DerivedC>& n_gravitation,
90 const Eigen::MatrixBase<DerivedD>& n_centrifugalAcceleration)
91{
92 return n_measuredForce
93 - n_coriolisAcceleration
94 + n_gravitation
95 - n_centrifugalAcceleration;
96}
97
116template<typename Derived>
117Eigen::Vector3<typename Derived::Scalar> lla_calcTimeDerivativeForPosition(const Eigen::MatrixBase<Derived>& n_velocity,
118 const typename Derived::Scalar& phi, const typename Derived::Scalar& h,
119 const typename Derived::Scalar& R_N, const typename Derived::Scalar& R_E)
120{
121 // Velocity North in [m/s]
122 const auto& v_N = n_velocity(0);
123 // Velocity East in [m/s]
124 const auto& v_E = n_velocity(1);
125 // Velocity Down in [m/s]
126 const auto& v_D = n_velocity(2);
127
128 return { v_N / (R_N + h),
129 v_E / ((R_E + h) * std::cos(phi)),
130 -v_D };
131}
132
138template<typename Scalar, typename = std::enable_if_t<std::is_floating_point_v<Scalar>>>
139Eigen::Vector<Scalar, 10> n_calcPosVelAttDerivative(const Eigen::Vector<Scalar, 10>& y, const Eigen::Vector<Scalar, 6>& z, const PosVelAttDerivativeConstants<Scalar>& c, Scalar /* t */ = 0.0)
140{
141 // 0 1 2 3 4 5 6 7 8 9
142 // [w, x, y, z, v_N, v_E, v_D, 𝜙, λ, h]^T
143 Eigen::Vector<Scalar, 10> y_dot = Eigen::Vector<Scalar, 10>::Zero();
144
145 Eigen::Quaternion<Scalar> n_Quat_b{ y(0), y(1), y(2), y(3) };
146 n_Quat_b.normalize();
147 Eigen::Quaternion<Scalar> n_Quat_e = trafo::n_Quat_e(y(7), y(8));
148
149 LOG_DATA("rollPitchYaw = {} [°]", rad2deg(trafo::quat2eulerZYX(n_Quat_b)).transpose());
150 LOG_DATA("n_velocity = {} [m/s]", y.template segment<3>(4).transpose());
151 LOG_DATA("lla_position = {}°, {}°, {} m", rad2deg(y(7)), rad2deg(y(8)), y(9));
152
153 auto R_N = calcEarthRadius_N(y(7));
154 LOG_DATA("R_N = {} [m]", R_N);
155 auto R_E = calcEarthRadius_E(y(7));
156 LOG_DATA("R_E = {} [m]", R_E);
157
158 // ω_ie_n Turn rate of the Earth expressed in local-navigation frame coordinates
159 Eigen::Vector3<Scalar> n_omega_ie = n_Quat_e * InsConst<>::e_omega_ie;
160 LOG_DATA("n_omega_ie = {} [rad/s]", n_omega_ie.transpose());
161 // ω_en_n Turn rate of the local frame with respect to the Earth-fixed frame, called the transport rate, expressed in local-navigation frame coordinates
162 Eigen::Vector3<Scalar> n_omega_en = n_calcTransportRate(y.template segment<3>(7), y.template segment<3>(4), R_N, R_E);
163 LOG_DATA("n_omega_en = {} [rad/s]", n_omega_en.transpose());
164 // ω_nb_b = ω_ib_b - C_bn * (ω_ie_n + ω_en_n)
165 Eigen::Vector3<Scalar> b_omega_nb = z.template segment<3>(3)
166 - n_Quat_b.conjugate()
167 * ((c.angularRateEarthRotationCompensationEnabled ? n_omega_ie : Eigen::Vector3<Scalar>::Zero())
168 + (c.angularRateTransportRateCompensationEnabled ? n_omega_en : Eigen::Vector3<Scalar>::Zero()));
169 LOG_DATA("b_omega_nb = {} [rad/s]", b_omega_nb.transpose());
170
171 // Coriolis acceleration in [m/s^2] (acceleration due to motion in rotating reference frame)
172 Eigen::Vector3<Scalar> n_coriolisAcceleration = c.coriolisAccelerationCompensationEnabled
173 ? n_calcCoriolisAcceleration(n_omega_ie, n_omega_en, y.template segment<3>(4))
174 : Eigen::Vector3<Scalar>::Zero();
175 LOG_DATA("n_coriolisAcceleration = {} [m/s^2]", n_coriolisAcceleration.transpose());
176 // Centrifugal acceleration in [m/s^2] (acceleration that makes a body follow a curved path)
177 Eigen::Vector3<Scalar> n_centrifugalAcceleration = c.centrifgalAccelerationCompensationEnabled
178 ? n_Quat_e * e_calcCentrifugalAcceleration(trafo::lla2ecef_WGS84(y.template segment<3>(7)), InsConst<>::e_omega_ie)
179 : Eigen::Vector3<Scalar>::Zero();
180 LOG_DATA("n_centrifugalAcceleration = {} [m/s^2]", n_centrifugalAcceleration.transpose());
181
182 Eigen::Vector3<Scalar> n_gravitation = n_calcGravitation(y.template segment<3>(7), c.gravitationModel);
183 LOG_DATA("n_gravitation = {} [m/s^2] ({})", n_gravitation.transpose(), to_string(c.gravitationModel));
184
185 y_dot.template segment<4>(0) = calcTimeDerivativeFor_n_Quat_b(b_omega_nb, // ω_nb_b Body rate with respect to the navigation frame, expressed in the body frame
186 y.template segment<4>(0)); // n_Quat_b_coeffs Coefficients of the quaternion n_Quat_b in order w, x, y, z (q = w + ix + jy + kz)
187
188 y_dot.template segment<3>(4) = n_calcTimeDerivativeForVelocity(n_Quat_b * z.template segment<3>(0), // f_n Specific force vector as measured by a triad of accelerometers and resolved into local-navigation frame coordinates
189 n_coriolisAcceleration, // Coriolis acceleration in local-navigation coordinates in [m/s^2]
190 n_gravitation, // Local gravitation vector (caused by effects of mass attraction) in local-navigation frame coordinates [m/s^2]
191 n_centrifugalAcceleration); // Centrifugal acceleration in local-navigation coordinates in [m/s^2]
192
193 y_dot.template segment<3>(7) = lla_calcTimeDerivativeForPosition(y.template segment<3>(4), // Velocity with respect to the Earth in local-navigation frame coordinates [m/s]
194 y(7), // 𝜙 Latitude in [rad]
195 y(9), // h Altitude in [m]
196 R_N, // North/South (meridian) earth radius [m]
197 R_E); // East/West (prime vertical) earth radius [m]
198
199 LOG_DATA("n_Quat_b_dot = {} ", y_dot.template segment<4>(0).transpose());
200 LOG_DATA("n_velocity_dot = {} [m/s^2]", y_dot.template segment<3>(4).transpose());
201 LOG_DATA("lla_position_dot = {} [rad/s, rad/s, m/s]", y_dot.template segment<3>(7).transpose());
202
203 return y_dot;
204}
205
206} // namespace NAV
Holds all Constants.
Transformation collection.
Functions concerning the ellipsoid model.
Scalar calcEarthRadius_E(const Scalar &latitude, const Scalar &a=InsConst< Scalar >::WGS84::a, const Scalar &e_squared=InsConst<>::WGS84::e_squared)
Calculates the East/West (prime vertical) earth radius.
Definition Ellipsoid.hpp:57
Scalar calcEarthRadius_N(const Scalar &latitude, const Scalar &a=InsConst<>::WGS84::a, const Scalar &e_squared=InsConst<>::WGS84::e_squared)
Calculates the North/South (meridian) earth radius.
Definition Ellipsoid.hpp:41
Different Gravity Models.
Eigen::Vector3< typename Derived::Scalar > n_calcGravitation(const Eigen::MatrixBase< Derived > &lla_position, GravitationModel gravitationModel=GravitationModel::EGM96)
Calculates the gravitation (acceleration due to mass attraction of the Earth)
Definition Gravity.hpp:193
Inertial Navigation Helper Functions.
Eigen::Vector3< typename DerivedA::Scalar > n_calcCoriolisAcceleration(const Eigen::MatrixBase< DerivedA > &n_omega_ie, const Eigen::MatrixBase< DerivedB > &n_omega_en, const Eigen::MatrixBase< DerivedC > &n_velocity)
Calculates the coriolis acceleration in [m/s^2] (acceleration due to motion in rotating reference fra...
Definition Functions.hpp:90
Eigen::Vector3< typename DerivedA::Scalar > e_calcCentrifugalAcceleration(const Eigen::MatrixBase< DerivedA > &e_position, const Eigen::MatrixBase< DerivedB > &e_omega_ie=InsConst< typename DerivedB::Scalar >::e_omega_ie)
Calculates the centrifugal acceleration in [m/s^2] (acceleration that makes a body follow a curved pa...
Definition Functions.hpp:72
Eigen::Vector3< typename DerivedA::Scalar > n_calcTransportRate(const Eigen::MatrixBase< DerivedA > &lla_position, const Eigen::MatrixBase< DerivedB > &n_velocity, const typename DerivedA::Scalar &R_N, const typename DerivedA::Scalar &R_E)
Calculates the transport rate ω_en_n (rotation rate of the Earth frame relative to the navigation fra...
Definition Functions.hpp:39
Eigen::Vector3< typename DerivedA::Scalar > n_calcTimeDerivativeForVelocity(const Eigen::MatrixBase< DerivedA > &n_measuredForce, const Eigen::MatrixBase< DerivedB > &n_coriolisAcceleration, const Eigen::MatrixBase< DerivedC > &n_gravitation, const Eigen::MatrixBase< DerivedD > &n_centrifugalAcceleration)
Calculates the time derivative of the velocity in local-navigation frame coordinates.
Definition Mechanization.hpp:87
Eigen::Vector3< typename Derived::Scalar > lla_calcTimeDerivativeForPosition(const Eigen::MatrixBase< Derived > &n_velocity, const typename Derived::Scalar &phi, const typename Derived::Scalar &h, const typename Derived::Scalar &R_N, const typename Derived::Scalar &R_E)
Calculates the time derivative of the curvilinear position.
Definition Mechanization.hpp:117
Eigen::Vector< Scalar, 10 > n_calcPosVelAttDerivative(const Eigen::Vector< Scalar, 10 > &y, const Eigen::Vector< Scalar, 6 > &z, const PosVelAttDerivativeConstants< Scalar > &c, Scalar=0.0)
Calculates the derivative of the quaternion, velocity and curvilinear position.
Definition Mechanization.hpp:139
Eigen::Vector4< typename DerivedA::Scalar > calcTimeDerivativeFor_n_Quat_b(const Eigen::MatrixBase< DerivedA > &b_omega_nb, const Eigen::MatrixBase< DerivedB > &n_Quat_b_coeffs)
Calculates the time derivative of the quaternion n_Quat_b.
Definition Mechanization.hpp:52
Utility class for logging to console and file.
#define LOG_DATA
All output which occurs repeatedly every time observations are received.
Definition Logger.hpp:29
Simple Math functions.
Inertial Navigation Mechanization Functions.
Provides Numerical integration methods.
const char * to_string(gui::widgets::PositionWithFrame::ReferenceFrame refFrame)
Converts the enum to a string.
Constants.
Definition Constants.hpp:26
Values needed to calculate the PosVelAttDerivative for the local-navigation frame.
Definition Mechanization.hpp:27
bool centrifgalAccelerationCompensationEnabled
Apply the centrifugal acceleration compensation to the measured accelerations.
Definition Mechanization.hpp:30
bool angularRateTransportRateCompensationEnabled
Apply the transport rate compensation to the measured angular rates.
Definition Mechanization.hpp:32
bool angularRateEarthRotationCompensationEnabled
Apply the Earth rotation rate compensation to the measured angular rates.
Definition Mechanization.hpp:31
GravitationModel gravitationModel
Gravity Model to use.
Definition Mechanization.hpp:28
bool coriolisAccelerationCompensationEnabled
Apply the Coriolis acceleration compensation to the measured accelerations.
Definition Mechanization.hpp:29