Library
Documentation for ReferenceFrameRotations.jl.
ReferenceFrameRotations.ReferenceFrameRotation — TypeReferenceFramerotationA Union of all supported rotation types.
ReferenceFrameRotations.DCM — TypeDCM{T}Direction Cosine Matrix (DCM) of type T, which is a 3x3 static matrix of type T.
Examples
julia> DCM(1.0I)
DCM{Float64}:
1.0 0.0 0.0
0.0 1.0 0.0
0.0 0.0 1.0
julia> DCM([1 0 0; 0 -1 0; 0 0 -1])
DCM{Int64}:
1 0 0
0 -1 0
0 0 -1ReferenceFrameRotations.EulerAngleAxis — TypeEulerAngleAxis{T}The definition of Euler Angle and Axis to represent a 3D rotation.
Fields
a::T: The Euler angle [rad].v::SVector{3, T}: The unitary vector aligned with the Euler axis.
Constructor
EulerAngleAxis(a::T1, v::AbstractVector{T2}) where {T1,T2}Create an Euler Angle and Axis representation structure with angle a [rad] and vector v.
The vector v will not be normalized.
The returned structure type will be selected according to the input types.
Examples
julia> EulerAngleAxis(pi / 3, [sqrt(2), sqrt(2), 0])
EulerAngleAxis{Float64}:
Euler angle: 1.0472 rad ( 60.0000 deg)
Euler axis: [ 1.4142, 1.4142, 0.0000]ReferenceFrameRotations.EulerAngles — TypeEulerAngles{T}The definition of Euler Angles, which is composed of three angles a1, a2, and a3 together with a rotation sequence rot_seq.
Fields
a1::T: First rotation [rad].a2::T: Second rotation [rad].a3::T: Third rotation [rad].rot_seq::Symbol: Rotation sequence.
rot_seq is provided by a symbol with three characters, each one indicating the rotation axis of the corresponding angle, e.g. :ZYX. The valid values for rot_seq are:
:XYX,:XYZ,:XZX,:XZY,:YXY,:YXZ,:YZX,:YZY,:ZXY,:ZXZ,:ZYX, andZYZ.
Constructor
EulerAngles(a1::T1, a2::T2, a3::T3, rot_seq::Symbol = :ZYX) where {T1, T2, T3}Create a new instance of EulerAngles with the angles a1, a2, and a3 and the rotation sequence rot_seq.
The type will be inferred from T1, T2, and T3.
If rot_seq is not provided, then it defaults to :ZYX.
Examples
julia> EulerAngles(pi / 2, pi / 4, -pi, :XYZ)
EulerAngles{Float64}:
R(X) : 1.5707963267948966 rad ( 90.0°)
R(Y) : 0.7853981633974483 rad ( 45.0°)
R(Z) : -3.141592653589793 rad (-180.0°)ReferenceFrameRotations.Quaternion — TypeQuaternion{T}The definition of the quaternion.
Fields
q0::T: Quaternion real part.q1::T: X component of the quaternion imaginary part.q2::T: Y component of the quaternion imaginary part.q3::T: Z component of the quaternion imaginary part.
Example
julia> Quaternion(cosd(45), sind(45), 0, 0)
Quaternion{Float64}:
+ 0.7071067811865476 + 0.7071067811865476.i + 0.0.j + 0.0.kReferenceFrameRotations.Quaternion — MethodQuaternion(q0::T0, q1::T1, q2::T2, q3::T3) where {T0, T1, T2, T3}Create the following quaternion:
q0 + q1.i + q2.j + q3.kin which:
q0is the real part of the quaternion.q1is the X component of the quaternion vectorial part.q2is the Y component of the quaternion vectorial part.q3is the Z component of the quaternion vectorial part.
Examples
julia> Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> Quaternion(1, 0, 0, 0.0)
Quaternion{Float64}:
+ 1.0 + 0.0⋅i + 0.0⋅j + 0.0⋅kQuaternion(v::AbstractVector)If the vector v has 3 components, then create a quaternion in which the real part is 0 and the vectorial or imaginary part has the same components of the vector v. In other words:
q = 0 + v[1].i + v[2].j + v[3].kOtherwise, if the vector v has 4 components, then create a quaternion in which the elements match those of the input vector:
q = v[1] + v[2].i + v[3].j + v[4].kExamples
julia> Quaternion([0, cosd(45), sind(45)])
Quaternion{Float64}:
+ 0.0 + 0.0⋅i + 0.707107⋅j + 0.707107⋅k
julia> Quaternion([cosd(45), 0, sind(45), 0])
Quaternion{Float64}:
+ 0.707107 + 0.0⋅i + 0.707107⋅j + 0.0⋅kQuaternion(r::Number, v::AbstractVector)Create a quaternion with real part r and vectorial or imaginary part v:
r + v[1].i + v[2].j + v[3].kExamples
julia> Quaternion(cosd(45), [0, sind(45), 0])
Quaternion{Float64}:
+ 0.707107 + 0.0⋅i + 0.707107⋅j + 0.0⋅kQuaternion(u::UniformScaling{T}) where T
Quaternion{T}(u::UniformScaling) where T
Quaternion(u::UniformScaling, Q::Quaternion{T}) where TCreate the quaternion u.λ + 0.i + 0.j + 0.k.
If a quaternion is passed as in the third signature, then the new quaternion will have the same type.
Examples
julia> Quaternion(I)
Quaternion{Bool}:
+ true + false⋅i + false⋅j + false⋅k
julia> Quaternion(1.0I)
Quaternion{Float64}:
+ 1.0 + 0.0⋅i + 0.0⋅j + 0.0⋅k
julia> q = Quaternion{Float32}(I)
Quaternion{Float32}:
+ 1.0 + 0.0⋅i + 0.0⋅j + 0.0⋅k
julia> Quaternion(I, q)
Quaternion{Float32}:
+ 1.0 + 0.0⋅i + 0.0⋅j + 0.0⋅kReferenceFrameRotations._EulerAngleConversion — Typestruct _EulerAngleConversion{R}This private structure is used only to enable the rotation conversion to Euler angles using the Julia API.
Base.:* — Method*(v::AbstractVector, q::Quaternion) -> Quaternion
*(q::Quaternion, v::AbstractVector) -> QuaternionCompute the multiplication qv * q or q * qv in which qv is a quaternion with real part 0 and vectorial/imaginary part v (Hamilton product).
Examples
julia> q = Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> v = [0, cosd(60), sind(60)]
3-element Vector{Float64}:
0.0
0.5
0.8660254037844386
julia> q * v
Quaternion{Float64}:
+ 0.0 + 0.0⋅i + 0.5⋅j + 0.866025⋅kBase.:* — Method*(Θ₂::EulerAngles, Θ₁::EulerAngles) -> EulerAnglesCompute the composed rotation of Θ₁ -> Θ₂.
The rotation will be represented by Euler angles (see EulerAngles) with the same rotation sequence as Θ₂.
Examples
julia> ea1 = EulerAngles(deg2rad(35), 0, 0, :XYZ)
EulerAngles{Float64}:
R(X) : 0.610865 rad ( 35.0°)
R(Y) : 0.0 rad ( 0.0°)
R(Z) : 0.0 rad ( 0.0°)
julia> ea2 = EulerAngles(0, 0, deg2rad(25), :ZYX)
EulerAngles{Float64}:
R(Z) : 0.0 rad ( 0.0°)
R(Y) : 0.0 rad ( 0.0°)
R(X) : 0.436332 rad ( 25.0°)
julia> ea2 * ea1
EulerAngles{Float64}:
R(Z) : 0.0 rad ( 0.0°)
R(Y) : -0.0 rad (-0.0°)
R(X) : 1.0472 rad ( 60.0°)Base.:* — Method*(λ::Number, q::Quaternion) -> Quaternion
*(q::Quaternion, λ::Number) -> QuaternionCompute λ * q or q * λ, in which λ is a scalar.
Examples
julia> q = Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> 2 * q
Quaternion{Int64}:
+ 2 + 0⋅i + 0⋅j + 0⋅kBase.:* — Method*(q1::Quaternion, q2::Quaternion) -> QuaternionCompute the quaternion multiplication q1 * q2 (Hamilton product).
If one of the operands is a UniformScaling:
*(u::UniformScaling, q::Quaternion)
*(q::Quaternion, u::UniformScaling)then it is considered as the quaternion u.λ + 0 ⋅ i + 0 ⋅ j + 0 ⋅ k.
Examples
julia> q1 = Quaternion(cosd(30), 0, sind(30), 0)
Quaternion{Float64}:
+ 0.866025 + 0.0⋅i + 0.5⋅j + 0.0⋅k
julia> q2 = Quaternion(cosd(60), 0, sind(60), 0)
Quaternion{Float64}:
+ 0.5 + 0.0⋅i + 0.866025⋅j + 0.0⋅k
julia> q1 * q2
Quaternion{Float64}:
+ 0.0 + 0.0⋅i + 1.0⋅j + 0.0⋅k
julia> I * q1
Quaternion{Float64}:
+ 0.866025 + 0.0⋅i + 0.5⋅j + 0.0⋅kBase.:* — Method*(av₂::EulerAngleAxis{T1}, av₁::EulerAngleAxis{T2}) where {T1, T2} -> EulerAngleAxis{T3}Compute the composed rotation of av₁ -> av₂.
The rotation will be represented by a Euler angle and axis (see EulerAngleAxis). By convention, the output angle will always be in the range [0, π] rad.
Notice that the vector representing the axis in av₁ and av₂ must be unitary. This function neither verifies this nor normalizes the vector.
Examples
julia> av1 = EulerAngleAxis(deg2rad(45), [sqrt(2)/2, sqrt(2)/2, 0])
EulerAngleAxis{Float64}:
Euler angle : 0.785398 rad (45.0°)
Euler axis : [0.707107, 0.707107, 0.0]
julia> av2 = EulerAngleAxis(deg2rad(22.5), [sqrt(2)/2, sqrt(2)/2, 0])
EulerAngleAxis{Float64}:
Euler angle : 0.392699 rad (22.5°)
Euler axis : [0.707107, 0.707107, 0.0]
julia> av1 * av2
EulerAngleAxis{Float64}:
Euler angle : 1.1781 rad (67.5°)
Euler axis : [0.707107, 0.707107, 0.0]Base.:+ — Method+(qa::Quaternion, qb::Quaternion) -> QuaternionCompute qa + qb.
If one of the operands is a UniformScaling:
+(u::UniformScaling, q::Quaternion)
+(q::Quaternion, u::UniformScaling)then it is considered as the quaternion u.λ + 0 ⋅ i + 0 ⋅ j + 0 ⋅ k.
Examples
julia> q1 = Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> q2 = Quaternion(0, cosd(45), 0, sind(45))
Quaternion{Float64}:
+ 0.0 + 0.707107⋅i + 0.0⋅j + 0.707107⋅k
julia> q1 + q2
Quaternion{Float64}:
+ 1.0 + 0.707107⋅i + 0.0⋅j + 0.707107⋅k
julia> q1 + I
Quaternion{Int64}:
+ 2 + 0⋅i + 0⋅j + 0⋅kBase.:- — Method-(qa::Quaternion, qb::Quaternion) -> QuaternionCompute qa - qb.
If one of the operands is a UniformScaling:
-(u::UniformScaling, q::Quaternion)
-(q::Quaternion, u::UniformScaling)then it is considered as the quaternion u.λ + 0 ⋅ i + 0 ⋅ j + 0 ⋅ k.
Examples
julia> q1 = Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> q2 = Quaternion(0, cosd(45), 0, sind(45))
Quaternion{Float64}:
+ 0.0 + 0.707107⋅i + 0.0⋅j + 0.707107⋅k
julia> q1 - q2
Quaternion{Float64}:
+ 1.0 - 0.707107⋅i + 0.0⋅j - 0.707107⋅k
julia> q1 - I
Quaternion{Int64}:
+ 0 + 0⋅i + 0⋅j + 0⋅kBase.:- — Method-(q::Quaternion) -> QuaternionReturn the quaterion -q.
Examples
julia> q = Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> -q
Quaternion{Int64}:
- 1 + 0⋅i + 0⋅j + 0⋅kBase.:/ — Method/(λ::Number, q::Quaternion) -> Quaternion
/(q::Quaternion, λ::Number) -> QuaternionCompute the division λ / q or q / λ, in which λ is a scalar.
Examples
julia> q = Quaternion(2, 0, 0, 0)
Quaternion{Int64}:
+ 2 + 0⋅i + 0⋅j + 0⋅k
julia> q / 2
Quaternion{Float64}:
+ 1.0 + 0.0⋅i + 0.0⋅j + 0.0⋅kBase.:/ — Method/(q1::Quaternion, q2::Quaternion) -> QuaternionCompute q1 * inv(q2) (Hamilton product).
If one of the operands is a UniformScaling:
/(u::UniformScaling, q::Quaternion)
/(q::Quaternion, u::UniformScaling)then it is considered as the quaternion u.λ + 0 ⋅ i + 0 ⋅ j + 0 ⋅ k.
Examples
julia> q1 = Quaternion(cosd(75), 0, sind(75), 0)
Quaternion{Float64}:
+ 0.258819 + 0.0⋅i + 0.965926⋅j + 0.0⋅k
julia> q2 = Quaternion(cosd(30), 0, sind(30), 0)
Quaternion{Float64}:
+ 0.866025 + 0.0⋅i + 0.5⋅j + 0.0⋅k
julia> q1 / q2
Quaternion{Float64}:
+ 0.707107 + 0.0⋅i + 0.707107⋅j + 0.0⋅k
julia> q1 / (2 * I)
Quaternion{Float64}:
+ 0.12941 + 0.0⋅i + 0.482963⋅j + 0.0⋅kBase.:\ — Method\(v::AbstractVector, q::Quaternion) -> Quaternion
\(q::Quaternion, v::AbstractVector) -> QuaternionCompute the division qv \ q or q \ qv in which qv is a quaternion with real part 0 and vectorial/imaginary part v (Hamilton product).
Examples
julia> q = Quaternion(1, 0, 0, 0)
Quaternion{Int64}:
+ 1 + 0⋅i + 0⋅j + 0⋅k
julia> v = [0, cosd(60), sind(60)]
3-element Vector{Float64}:
0.0
0.5
0.8660254037844386
julia> v \ q
Quaternion{Float64}:
+ 0.0 + 0.0⋅i - 0.5⋅j - 0.866025⋅kBase.:\ — Method\(q1::Quaternion, q2::Quaternion) -> QuaternionCompute inv(q1) * q2.
If one of the operands is a UniformScaling:
\(u::UniformScaling, q::Quaternion)
\(q::Quaternion, u::UniformScaling)then it is considered as the quaternion u.λ + 0 ⋅ i + 0 ⋅ j + 0 ⋅ k.
Examples
julia> q1 = Quaternion(cosd(75), 0, sind(75), 0)
Quaternion{Float64}:
+ 0.258819 + 0.0⋅i + 0.965926⋅j + 0.0⋅k
julia> q2 = Quaternion(cosd(30), 0, sind(30), 0)
Quaternion{Float64}:
+ 0.866025 + 0.0⋅i + 0.5⋅j + 0.0⋅k
julia> q2 \ q1
Quaternion{Float64}:
+ 0.707107 + 0.0⋅i + 0.707107⋅j + 0.0⋅kBase.conj — Methodconj(q::Quaternion) -> QuaternionCompute the complex conjugate of the quaternion q:
q0 - q1.i - q2.j - q3.kSee also: inv
Examples
julia> q = Quaternion(1, cosd(75), 0, sind(75))
Quaternion{Float64}:
+ 1.0 + 0.258819⋅i + 0.0⋅j + 0.965926⋅k
julia> conj(q)
Quaternion{Float64}:
+ 1.0 - 0.258819⋅i - 0.0⋅j - 0.965926⋅kBase.copy — Methodcopy(q::Quaternion{T}) where T -> QuaternionCreate a copy of the quaternion q.
Base.imag — Methodimag(q::Quaternion{T}) -> SVector{3, T}Return the vectorial or imaginary part of the quaternion q represented by a 3 × 1 vector of type SVector{3}.
Examples
julia> q = Quaternion(cosd(75), 0, sind(75), 0)
Quaternion{Float64}:
+ 0.258819 + 0.0⋅i + 0.965926⋅j + 0.0⋅k
julia> imag(q)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
0.0
0.9659258262890683
0.0Base.inv — Methodinv(Θ::EulerAngles) -> EulerAnglesReturn the Euler angles that represent the inverse rotation of Θ.
The rotation sequence of the result will be the inverse of the input. Hence, if the input rotation sequence is, for example, :XYZ, then the result will be represented using :ZYX.
Examples
julia> ea = EulerAngles(π / 3, π / 6, 2 / 3 * π, :ZYX)
EulerAngles{Float64}:
R(Z) : 1.0472 rad ( 60.0°)
R(Y) : 0.523599 rad ( 30.0°)
R(X) : 2.0944 rad ( 120.0°)
julia> inv(ea)
EulerAngles{Float64}:
R(X) : -2.0944 rad (-120.0°)
R(Y) : -0.523599 rad (-30.0°)
R(Z) : -1.0472 rad (-60.0°)Base.inv — Methodinv(q::Quaternion) -> QuaternionCompute the inverse of the quaternion q:
conj(q)
───────
|q|²See also: conj
Examples
julia> q = Quaternion(1, 0, cosd(75), sind(75))
Quaternion{Float64}:
+ 1.0 + 0.0⋅i + 0.258819⋅j + 0.965926⋅k
julia> inv(q)
Quaternion{Float64}:
+ 0.5 - 0.0⋅i - 0.12941⋅j - 0.482963⋅kBase.inv — Methodinv(av::EulerAngleAxis) -> EulerAngleAxisCompute the inverse rotation of the Euler angle and axis av.
The Euler angle returned by this function will always be in the interval [0, π] rad.
Examples
julia> av = EulerAngleAxis(deg2rad(20), [sqrt(2) / 2, 0, sqrt(2) / 2])
EulerAngleAxis{Float64}:
Euler angle : 0.349066 rad (20.0°)
Euler axis : [0.707107, 0.0, 0.707107]
julia> inv(av)
EulerAngleAxis{Float64}:
Euler angle : 0.349066 rad (20.0°)
Euler axis : [-0.707107, -0.0, -0.707107]
julia> av = EulerAngleAxis(deg2rad(-20), [sqrt(2) / 2, 0, sqrt(2) / 2])
EulerAngleAxis{Float64}:
Euler angle : -0.349066 rad (-20.0°)
Euler axis : [0.707107, 0.0, 0.707107]
julia> inv(av)
EulerAngleAxis{Float64}:
Euler angle : 0.349066 rad (20.0°)
Euler axis : [0.707107, 0.0, 0.707107]Base.real — Methodreal(q::Quaternion{T}) -> TReturn the real part of the quaternion q: q0.
Examples
julia> q = Quaternion(cosd(75), 0, sind(75), 0)
Quaternion{Float64}:
+ 0.258819 + 0.0⋅i + 0.965926⋅j + 0.0⋅k
julia> real(q)
0.25881904510252074LinearAlgebra.norm — Methodnorm(q::Quaternion{T}) -> float(T)Compute the Euclidean norm of the quaternion q:
√(q0² + q1² + q2² + q3²)Examples
julia> q = Quaternion(cosd(75), 0, sind(75), 0)
Quaternion{Float64}:
+ 0.258819 + 0.0⋅i + 0.965926⋅j + 0.0⋅k
julia> norm(q)
1.0ReferenceFrameRotations.angle_to_angle — Methodangle_to_angle(θ₁::Number, θ₂::Number, θ₃::Number, rot_seq_orig::Symbol, rot_seq_dest::Symbol) -> EulerAngles
angle_to_angle(Θ::EulerAngles, rot_seq_dest::Symbol) -> EulerAnglesConvert the Euler angles θ₁, θ₂, and θ₃ [rad] with the rotation sequence rot_seq_orig to a new set of Euler angles with rotation sequence rot_seq_dest.
The input values of the origin Euler angles can also be passed inside the structure Θ (see EulerAngles).
The rotation sequence is defined by a :Symbol. The possible values are: :XYX, XYZ, :XZX, :XZY, :YXY, :YXZ, :YZX, :YZY, :ZXY, :ZXZ, :ZYX, and :ZYZ.
Example
julia> angle_to_angle(-pi / 2, -pi / 3, -pi / 4, :ZYX, :XYZ)
EulerAngles{Float64}:
R(X) : -1.0472 rad (-60.0°)
R(Y) : 0.785398 rad ( 45.0°)
R(Z) : -1.5708 rad (-90.0°)
julia> angle_to_angle(-pi / 2, 0, 0, :ZYX, :XYZ)
EulerAngles{Float64}:
R(X) : 0.0 rad ( 0.0°)
R(Y) : 0.0 rad ( 0.0°)
R(Z) : -1.5708 rad (-90.0°)
julia> Θ = EulerAngles(1, 2, 3, :XYX)
EulerAngles{Int64}:
R(X) : 1 rad ( 57.2958°)
R(Y) : 2 rad ( 114.592°)
R(X) : 3 rad ( 171.887°)
julia> angle_to_angle(Θ, :ZYZ)
EulerAngles{Float64}:
R(Z) : -2.70239 rad (-154.836°)
R(Y) : 1.46676 rad ( 84.0393°)
R(Z) : -1.05415 rad (-60.3984°)ReferenceFrameRotations.angle_to_angleaxis — Functionangle_to_angleaxis(θ₁::Number, θ₂::Number, θ₃::Number, rot_seq::Symbol = :ZYX) -> EulerAngleAxis
angle_to_angleaxis(Θ::EulerAngles) -> EulerAngleAxisConvert the Euler angles θ₁, θ₂, and θ₃ [rad] with the rotation sequence rot_seq to an Euler angle and axis representation.
Those values can also be passed inside the structure Θ (see EulerAngles).
The rotation sequence is defined by a :Symbol. The possible values are: :XYX, XYZ, :XZX, :XZY, :YXY, :YXZ, :YZX, :YZY, :ZXY, :ZXZ, :ZYX, and :ZYZ. If no value is specified, it defaults to :ZYX.
Example
julia> angle_to_angleaxis(1, 0, 0, :XYZ)
EulerAngleAxis{Float64}:
Euler angle : 1.0 rad (57.2958°)
Euler axis : [1.0, 0.0, 0.0]
julia> Θ = EulerAngles(1, 1, 1, :XYZ);
julia> angle_to_angleaxis(Θ)
EulerAngleAxis{Float64}:
Euler angle : 1.93909 rad (111.102°)
Euler axis : [0.692363, 0.203145, 0.692363]ReferenceFrameRotations.angle_to_dcm — Methodangle_to_dcm(θ₁::Number[, θ₂::Number[, θ₃::Number]], rot_seq::Symbol = :ZYX) -> DCM
angle_to_dcm(Θ::EulerAngles) -> DCMCreate a direction cosine matrix that perform a set of rotations (θ₁, θ₂, θ₃) about the coordinate axes specified in rot_seq.
The input values of the origin Euler angles can also be passed inside the structure Θ (see EulerAngles).
The rotation sequence is defined by a Symbol specifing the rotation axes. The possible values depends on the number of rotations as follows:
- 1 rotation (
θ₁)::X,:Y, or:Z. - 2 rotations (
θ₁,θ₂)::XY,:XZ,:YX,:YZ,:ZX, or:ZY. - 3 rotations (
θ₁,θ₂,θ₃)::XYX,XYZ,:XZX,:XZY,:YXY,:YXZ,:YZX,:YZY,:ZXY,:ZXZ,:ZYX, or:ZYZ
Remarks
This function assigns dcm = A3 * A2 * A1 in which Ai is the DCM related with the i-th rotation, i ∈ [1,2,3]. If the i-th rotation is not specified, then Ai = I.
Example
julia> angle_to_dcm(pi / 2, :X)
DCM{Float64}:
1.0 0.0 0.0
0.0 6.12323e-17 1.0
0.0 -1.0 6.12323e-17
julia> angle_to_dcm(pi / 5, pi / 7, :YZ)
DCM{Float64}:
0.728899 0.433884 -0.529576
-0.351019 0.900969 0.25503
0.587785 0.0 0.809017
julia> angle_to_dcm(pi / 5, pi / 7, 0, :YZY)
DCM{Float64}:
0.728899 0.433884 -0.529576
-0.351019 0.900969 0.25503
0.587785 0.0 0.809017
julia> dcm = angle_to_dcm(pi / 2, pi / 3, pi / 4, :ZYX)
DCM{Float64}:
3.06162e-17 0.5 -0.866025
-0.707107 0.612372 0.353553
0.707107 0.612372 0.353553ReferenceFrameRotations.angle_to_quat — Methodangle_to_quat(θ₁::T1[, θ₂::T2[, θ₃::T3]], rot_seq::Symbol = :ZYX) where {T1<:Number, T2<:Number, T3<:Number} -> Quaternion
angle_to_quat(eulerang::EulerAngles) -> QuaternionCreate a quaternion that perform a set of rotations (θ₁, θ₂, θ₃) about the coordinate axes specified in rot_seq.
The input values of the origin Euler angles can also be passed inside the structure Θ (see EulerAngles).
The rotation sequence is defined by a Symbol specifing the rotation axes. The possible values depends on the number of rotations as follows:
- 1 rotation (
θ₁)::X,:Y, or:Z. - 2 rotations (
θ₁,θ₂)::XY,:XZ,:YX,:YZ,:ZX, or:ZY. - 3 rotations (
θ₁,θ₂,θ₃)::XYX,XYZ,:XZX,:XZY,:YXY,:YXZ,:YZX,:YZY,:ZXY,:ZXZ,:ZYX, or:ZYZ
Remarks
This function assigns q = q1 * q2 * q3 in which qi is the quaternion related with the i-th rotation, i ∈ [1,2,3]. If the i-th rotation is not specified, then qi = Quaternion(I).
Example
julia> angle_to_quat(pi / 2, :X)
Quaternion{Float64}:
+ 0.707107 + 0.707107⋅i + 0.0⋅j + 0.0⋅k
julia> angle_to_quat(pi / 5, pi / 7, :YZ)
Quaternion{Float64}:
+ 0.927212 + 0.0687628⋅i + 0.301269⋅j + 0.21163⋅k
julia> angle_to_quat(pi / 5, pi / 7, 0, :YZX)
Quaternion{Float64}:
+ 0.927212 + 0.0687628⋅i + 0.301269⋅j + 0.21163⋅k
julia> angle_to_quat(pi / 2, pi / 3, pi / 4, :ZYX)
Quaternion{Float64}:
+ 0.701057 - 0.092296⋅i + 0.560986⋅j + 0.430459⋅kReferenceFrameRotations.angle_to_rot — Methodangle_to_rot([T,] θ₁::Number[, θ₂::Number[, θ₃::Number]], rot_seq::Symbol) -> T
angle_to_rot([T,] Θ::EulerAngles) -> TCreate a rotation description of type T that perform a set of rotations (θ₁, θ₂, θ₃) about the coordinate axes specified in rot_seq.
The input values of the origin Euler angles can also be passed inside the structure Θ (see EulerAngles).
The rotation sequence is defined by a Symbol specifing the rotation axes. The possible values depends on the number of rotations as follows:
- 1 rotation (
θ₁)::X,:Y, or:Z. - 2 rotations (
θ₁,θ₂)::XY,:XZ,:YX,:YZ,:ZX, or:ZY. - 3 rotations (
θ₁,θ₂,θ₃)::XYX,XYZ,:XZX,:XZY,:YXY,:YXZ,:YZX,:YZY,:ZXY,:ZXZ,:ZYX, or:ZYZ
Example
julia> dcm = angle_to_rot(pi / 5, :Z)
3×3 StaticArrays.SMatrix{3, 3, Float64, 9} with indices SOneTo(3)×SOneTo(3):
0.809017 0.587785 0.0
-0.587785 0.809017 0.0
0.0 0.0 1.0
julia> quat = angle_to_rot(Quaternion, pi / 5, :Z)
Quaternion{Float64}:
+ 0.951057 + 0.0⋅i + 0.0⋅j + 0.309017⋅k
julia> dcm = angle_to_rot(pi / 5, pi / 7, :YZ)
3×3 StaticArrays.SMatrix{3, 3, Float64, 9} with indices SOneTo(3)×SOneTo(3):
0.728899 0.433884 -0.529576
-0.351019 0.900969 0.25503
0.587785 0.0 0.809017
julia> quat = angle_to_rot(Quaternion, pi / 5, pi / 7, :YZ)
Quaternion{Float64}:
+ 0.927212 + 0.0687628⋅i + 0.301269⋅j + 0.21163⋅k
julia> dcm = angle_to_rot(pi / 2, pi / 3, pi / 4, :ZYX)
3×3 StaticArrays.SMatrix{3, 3, Float64, 9} with indices SOneTo(3)×SOneTo(3):
3.06162e-17 0.5 -0.866025
-0.707107 0.612372 0.353553
0.707107 0.612372 0.353553
julia> q = angle_to_rot(Quaternion, pi / 2, pi / 3, pi / 4, :ZYX)
Quaternion{Float64}:
+ 0.701057 - 0.092296⋅i + 0.560986⋅j + 0.430459⋅kReferenceFrameRotations.angleaxis_to_angle — Methodangleaxis_to_angle(θ::Number, v::AbstractVector, rot_seq::Symbol) -> EulerAngles
angleaxis_to_angle(av::EulerAngleAxis, rot_seq::Symbol) -> EulerAnglesConvert the Euler angle θ [rad] and Euler axis v to Euler angles with rotation sequence rot_seq.
Those values can also be passed inside the structure av (see EulerAngleAxis).
The rotation sequence is defined by a :Symbol. The possible values are: :XYX, XYZ, :XZX, :XZY, :YXY, :YXZ, :YZX, :YZY, :ZXY, :ZXZ, :ZYX, and :ZYZ. If no value is specified, it defaults to :ZYX.
It is expected that the vector v is unitary. However, no verification is performed inside the function. The user must handle this situation.
Example
julia> av = EulerAngleAxis(deg2rad(45), [1, 0, 0]);
julia> angleaxis_to_angle(av, :ZXY)
EulerAngles{Float64}:
R(Z) : 0.0 rad ( 0.0°)
R(X) : 0.785398 rad ( 45.0°)
R(Y) : 0.0 rad ( 0.0°)ReferenceFrameRotations.angleaxis_to_dcm — Methodangleaxis_to_dcm(a::Number, v::AbstractVector) -> DCM
angleaxis_to_dcm(av::EulerAngleAxis) -> DCMConvert the Euler angle a [rad] and Euler axis v to a DCM.
Those values can also be passed inside the structure ea (see EulerAngleAxis).
It is expected that the vector v is unitary. However, no verification is performed inside the function. The user must handle this situation.
Example
julia> v = [1, 1, 1];
julia> v /= norm(v);
julia> angleaxis_to_dcm(pi / 2, v)
DCM{Float64}:
0.333333 0.910684 -0.244017
-0.244017 0.333333 0.910684
0.910684 -0.244017 0.333333
julia> ea = EulerAngleAxis(pi / 2, v);
julia> angleaxis_to_dcm(ea)
DCM{Float64}:
0.333333 0.910684 -0.244017
-0.244017 0.333333 0.910684
0.910684 -0.244017 0.333333ReferenceFrameRotations.angleaxis_to_quat — Methodangleaxis_to_quat(θ::Number, v::AbstractVector) -> Quaternion
angleaxis_to_quat(angleaxis::EulerAngleAxis) -> QuaternionConvert the Euler angle θ [rad] and Euler axis v to a quaternion.
Those values can also be passed inside the structure ea (see EulerAngleAxis).
It is expected that the vector v is unitary. However, no verification is performed inside the function. The user must handle this situation.
Example
julia> v = [1, 1, 1];
julia> v /= norm(v);
julia> angleaxis_to_quat(pi / 2, v)
Quaternion{Float64}:
+ 0.707107 + 0.408248⋅i + 0.408248⋅j + 0.408248⋅kReferenceFrameRotations.compose_rotation — Methodcompose_rotation(R1::T, [, R2::T, R3::T, R4::T, R5::T, ...]) -> TCompute a composed rotation using the rotations R1, R2, R3, R4, ..., in the following order:
First rotation
|
|
R1 => R2 => R3 => R4 => ...
|
|
Second rotationThe rotations can be described by:
- A direction cosine matrix (
DCM); - An Euler angle and axis (
EulerAngleAxis); - A set of Euler angles (
EulerAngles); or - A quaternion (
Quaternion).
Notice, however, that all rotations must be of the same type (DCM or quaternion).
The output will have the same type as the inputs.
Example
julia> D1 = angle_to_dcm(pi / 3, pi / 4, pi / 5, :ZYX);
julia> D2 = angle_to_dcm(-pi / 5, -pi / 4, -pi / 3, :XYZ);
julia> compose_rotation(D1, D2)
DCM{Float64}:
1.0 1.08801e-17 3.54837e-17
1.08801e-17 1.0 2.88714e-17
3.54837e-17 2.88714e-17 1.0
julia> ea1 = EulerAngleAxis(30 * pi / 180, [0, 1, 0]);
julia> ea2 = EulerAngleAxis(45 * pi / 180, [0, 1, 0]);
julia> compose_rotation(ea1, ea2)
EulerAngleAxis{Float64}:
Euler angle : 1.309 rad (75.0°)
Euler axis : [0.0, 1.0, 0.0]
julia> Θ1 = EulerAngles(1, 2, 3, :ZYX);
julia> Θ2 = EulerAngles(-3, -2, -1, :XYZ);
julia> compose_rotation(Θ1, Θ2)
EulerAngles{Float64}:
R(X) : -1.66533e-16 rad (-9.54166e-15°)
R(Y) : 9.24446e-33 rad ( 5.29669e-31°)
R(Z) : -1.11022e-16 rad (-6.36111e-15°)
julia> q1 = angle_to_quat(pi / 3, pi / 4, pi / 5, :ZYX);
julia> q2 = angle_to_quat(-pi / 5, -pi / 4, -pi / 3, :XYZ);
julia> compose_rotation(q1, q2)
Quaternion{Float64}:
+ 1.0 + 0.0⋅i + 2.08167e-17⋅j + 5.55112e-17⋅kReferenceFrameRotations.dcm_to_angle — Methoddcm_to_angle(dcm::DCM, rot_seq::Symbol=:ZYX) -> EulerAnglesConvert the dcm to Euler Angles (see EulerAngles) given a rotation sequence rot_seq.
The rotation sequence is defined by a :Symbol. The possible values are: :XYX, XYZ, :XZX, :XZY, :YXY, :YXZ, :YZX, :YZY, :ZXY, :ZXZ, :ZYX, and :ZYZ. If no value is specified, it defaults to :ZYX.
Gimbal-lock and special cases
If the rotations are about three different axes, e.g. :XYZ, :ZYX, etc., then a second rotation of ±90˚ yields a gimbal-lock. This means that the rotations between the first and third axes have the same effect. In this case, the net rotation angle is assigned to the first rotation, and the angle of the third rotation is set to 0.
If the rotations are about two different axes, e.g. :XYX, :YXY, etc., then a rotation about the duplicated axis yields multiple representations. In this case, the entire angle is assigned to the first rotation and the third rotation is set to 0.
Example
julia> D = DCM([1. 0. 0.; 0. 0. -1; 0. -1 0.]);
julia> dcm_to_angle(D,:XYZ)
EulerAngles{Float64}:
R(X) : 1.5708 rad ( 90.0°)
R(Y) : 0.0 rad ( 0.0°)
R(Z) : 0.0 rad ( 0.0°)
julia> D = angle_to_dcm(1, -pi / 2, 2, :ZYX);
julia> dcm_to_angle(D, :ZYX)
EulerAngles{Float64}:
R(Z) : 3.0 rad ( 171.887°)
R(Y) : -1.5708 rad (-90.0°)
R(X) : 0.0 rad ( 0.0°)
julia> D = angle_to_dcm(1, :X) * angle_to_dcm(2, :X);
julia> dcm_to_angle(D, :XYX)
EulerAngles{Float64}:
R(X) : 3.0 rad ( 171.887°)
R(Y) : 0.0 rad ( 0.0°)
R(X) : 0.0 rad ( 0.0°)ReferenceFrameRotations.dcm_to_angleaxis — Methoddcm_to_angleaxis(dcm::DCM{T}) where T<:Number -> EulerAngleAxisConvert the dcm to an Euler angle and axis representation.
By convention, the returned Euler angle will always be in the interval [0, π].
ReferenceFrameRotations.dcm_to_quat — Methoddcm_to_quat(dcm::DCM) -> QuaternionConvert the dcm to a quaternion.
The type of the quaternion will be automatically selected by the constructor Quaternion to avoid InexactError.
Remarks
By convention, the real part of the quaternion will always be positive. Moreover, the function does not check if dcm is a valid direction cosine matrix. This must be handle by the user.
This algorithm was obtained from [1].
Example
julia> dcm = angle_to_dcm(pi / 2, 0.0, 0.0, :XYZ);
julia> q = dcm_to_quat(dcm)
Quaternion{Float64}:
+ 0.707107 + 0.707107⋅i + 0.0⋅j + 0.0⋅kReferences
- [1]: http://www.euclideanspace.com/maths/geometry/rotations/conversions/matrixToQuaternion/
ReferenceFrameRotations.ddcm — Methodddcm(Dba::DCM, wba_b::AbstractArray) -> SMatrix{3, 3}Compute the time-derivative of the dcm that rotates a reference frame a into alignment with the reference frame b in which the angular velocity of b with respect to a, and represented in b, is wba_b.
Example
julia> D = DCM(1.0I);
julia> ddcm(D, [1, 0, 0])
3×3 StaticArrays.SArray{Tuple{3,3},Float64,2,9}:
0.0 0.0 0.0
0.0 0.0 1.0
0.0 -1.0 0.0ReferenceFrameRotations.dquat — Methoddquat(qba::Quaternion, wba_b::AbstractVector) -> QuaternionCompute the time-derivative of the quaternion qba that rotates a reference frame a into alignment to the reference frame b in which the angular velocity of b with respect to a, and represented in b, is wba_b.
Examples
julia> q = Quaternion(1.0I);
julia> dquat(q,[1;0;0])
Quaternion{Float64}:
- 0.0 + 0.5⋅i + 0.0⋅j + 0.0⋅kReferenceFrameRotations.inv_rotation — Methodinv_rotation(R::T) -> TCompute the inverse rotation of R, which can be:
- A direction cosine matrix (
DCM); - An Euler angle and axis (
EulerAngleAxis); - A set of Euler anlges (
EulerAngles); or - A quaternion (
Quaternion).
The output will have the same type as R.
If R is a DCM, than its transpose is computed instead of its inverse to reduce the computational burden. The both are equal if the DCM has unit norm. This must be verified by the user.
If R is a quaternion, than its conjugate is computed instead of its inverse to reduce the computational burden. The both are equal if the quaternion has unit norm. This must be verified by the used.
Example
julia> D = angle_to_dcm(pi / 3, pi / 4, pi / 5, :ZYX);
julia> inv_rotation(D)
DCM{Float64}:
0.353553 -0.492816 0.795068
0.612372 0.764452 0.201527
-0.707107 0.415627 0.572061
julia> ea = EulerAngleAxis(30 * pi / 180, [1, 0, 0]);
julia> inv_rotation(ea)
EulerAngleAxis{Float64}:
Euler angle : 0.523599 rad (30.0°)
Euler axis : [-1.0, -0.0, -0.0]
julia> Θ = EulerAngles(-pi / 3, -pi / 2, -pi, :YXZ);
julia> inv_rotation(Θ)
EulerAngles{Float64}:
R(Z) : 3.14159 rad ( 180.0°)
R(X) : 1.5708 rad ( 90.0°)
R(Y) : 1.0472 rad ( 60.0°)
julia> q = angle_to_quat(pi / 3, pi / 4, pi / 5, :ZYX);
julia> inv_rotation(q)
Quaternion{Float64}:
+ 0.820071 - 0.0652687⋅i - 0.45794⋅j - 0.336918⋅kReferenceFrameRotations.orthonormalize — Methodorthonormalize(dcm::DCM) -> DCMPerform the Gram-Schmidt orthonormalization process in the dcm and return the new matrix.
This function does not check if the columns of the input matrix span a three-dimensional space. If not, then the returned matrix should have NaN. Notice, however, that such input matrix is not a valid direction cosine matrix.
Example
julia> D = DCM(3I)
julia> orthonormalize(D)
3×3 StaticArrays.SArray{Tuple{3,3},Float64,2,9}:
1.0 0.0 0.0
0.0 1.0 0.0
0.0 0.0 1.0ReferenceFrameRotations.quat_to_angle — Functionquat_to_angle(q::Quaternion, rot_seq::Symbol = :ZYX) -> EulerAnglesConvert the quaternion q to Euler Angles (see EulerAngles) given a rotation sequence rot_seq.
The rotation sequence is defined by a :Symbol. The possible values are: :XYX, XYZ, :XZX, :XZY, :YXY, :YXZ, :YZX, :YZY, :ZXY, :ZXZ, :ZYX, and :ZYZ. If no value is specified, it defaults to :ZYX.
Examples
julia> q = Quaternion(cosd(45/2), sind(45/2), 0, 0);
julia> quat_to_angle(q, :XYZ)
EulerAngles{Float64}:
R(X) : 0.785398 rad ( 45.0°)
R(Y) : 0.0 rad ( 0.0°)
R(Z) : 0.0 rad ( 0.0°)ReferenceFrameRotations.quat_to_angleaxis — Methodquat_to_angleaxis(q::Quaternion{T}) where T -> EulerAngleAxisConvert the quaternion q to a Euler angle and axis representation (see EulerAngleAxis). By convention, the Euler angle will be kept between [0, π] rad.
Remarks
This function will not fail if the quaternion norm is not 1. However, the meaning of the results will not be defined, because the input quaternion does not represent a 3D rotation. The user must handle such situations.
Examples
julia> q = Quaternion(cosd(45/2), sind(45/2), 0, 0);
julia> quat_to_angleaxis(q)
EulerAngleAxis{Float64}:
Euler angle : 0.785398 rad (45.0°)
Euler axis : [1.0, 0.0, 0.0]ReferenceFrameRotations.quat_to_dcm — Methodquat_to_dcm(q::Quaternion) -> DCMConvert the quaternion q to a Direction Cosine Matrix (DCM).
Examples
julia> q = Quaternion(cosd(45/2), sind(45/2), 0, 0);
julia> quat_to_dcm(q)
DCM{Float64}:
1.0 0.0 0.0
0.0 0.707107 0.707107
0.0 -0.707107 0.707107ReferenceFrameRotations.smallangle_to_dcm — Methodsmallangle_to_dcm(θx::Number, θy::Number, θz::Number; normalize = true) -> DCMCreate a direction cosine matrix from three small rotations of angles θx, θy, and θz [rad] about the axes X, Y, and Z, respectively.
If the keyword normalize is true, the matrix will be normalized using the function orthonormalize.
Example
julia> smallangle_to_dcm(+0.01, -0.01, +0.01)
DCM{Float64}:
0.9999 0.00989903 0.010098
-0.009999 0.999901 0.00989802
-0.009999 -0.009998 0.9999
julia> smallangle_to_dcm(+0.01, -0.01, +0.01; normalize = false)
DCM{Float64}:
1.0 0.01 0.01
-0.01 1.0 0.01
-0.01 -0.01 1.0ReferenceFrameRotations.smallangle_to_quat — Methodsmallangle_to_quat(θx::Number, θy::Number, θz::Number) -> QuaternionCreate a quaternion from three small rotations of angles θx, θy, and θz [rad] about the axes X, Y, and Z, respectively.
Example
julia> smallangle_to_quat(+0.01, -0.01, +0.01)
Quaternion{Float64}:
+ 0.999963 + 0.00499981⋅i - 0.00499981⋅j + 0.00499981⋅kReferenceFrameRotations.smallangle_to_rot — Methodsmallangle_to_rot([T,] θx::Number, θy::Number, θz::Number[; normalize = true]) -> TCreate a rotation description of type T from three small rotations of angles θx, θy, and θz [rad] about the axes X, Y, and Z, respectively.
The type T of the rotation description can be DCM or Quaternion. If the type T is not specified, if defaults to DCM.
If T is DCM, the resulting matrix will be orthonormalized using the orthonormalize function if the keyword normalize is true.
Example
julia> dcm = smallangle_to_rot(+0.01, -0.01, +0.01)
DCM{Float64}:
0.9999 0.00989903 0.010098
-0.009999 0.999901 0.00989802
-0.009999 -0.009998 0.9999
julia> dcm = smallangle_to_rot(+0.01, -0.01, +0.01; normalize = false)
DCM{Float64}:
1.0 0.01 0.01
-0.01 1.0 0.01
-0.01 -0.01 1.0
julia> q = smallangle_to_rot(Quaternion, +0.01, -0.01, +0.01)
Quaternion{Float64}:
+ 0.999963 + 0.00499981⋅i - 0.00499981⋅j + 0.00499981⋅kReferenceFrameRotations.vect — Methodvect(q::Quaternion{T}) -> SVector{3, T}Return the vectorial or imaginary part of the quaternion q represented by a 3 × 1 vector of type SVector{3, T}.
Examples
julia> q = Quaternion(cosd(75), 0, sind(75), 0)
Quaternion{Float64}:
+ 0.258819 + 0.0⋅i + 0.965926⋅j + 0.0⋅k
julia> vect(q)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
0.0
0.9659258262890683
0.0