ModelSpace
All Classes Namespaces Functions Variables Enumerations Pages
Euler321.hpp
1/******************************************************************************
2* Copyright (c) ATTX LLC 2024. All Rights Reserved.
3*
4* This software and associated documentation (the "Software") are the
5* proprietary and confidential information of ATTX, LLC. The Software is
6* furnished under a license agreement between ATTX and the user organization
7* and may be used or copied only in accordance with the terms of the agreement.
8* Refer to 'license/attx_license.adoc' for standard license terms.
9*
10* EXPORT CONTROL NOTICE: THIS SOFTWARE MAY INCLUDE CONTENT CONTROLLED UNDER THE
11* INTERNATIONAL TRAFFIC IN ARMS REGULATIONS (ITAR) OR THE EXPORT ADMINISTRATION
12* REGULATIONS (EAR99). No part of the Software may be used, reproduced, or
13* transmitted in any form or by any means, for any purpose, without the express
14* written permission of ATTX, LLC.
15******************************************************************************/
16/*
17Euler Angle 321 set header file
18
19Author: Alex Reynolds
20*/
21#ifndef EULER321_HPP
22#define EULER321_HPP
23
24#include <cmath>
25
26#include "core/Matrix.hpp"
27#include "core/CartesianVector.hpp"
28#include "core/safemath.hpp"
29#include "core/matrixmath.hpp"
30
31#include "DCM.hpp"
32
33namespace clockwerk {
34
35 /// @brief Class defining a 3-2-1 Euler angle sequence
36 ///
37 /// This file defines a simple 321 Euler angle attitude representation for cartesian
38 /// coordinate systems. It is largely the same as the base vector,
39 /// with the following exceptions:
40 /// - Set to be 3 elements
41 /// - Functions defined to convert to other attitude representations
42 /// - Function defined to calculate rate of change as a function of omega
43 ///
44 /// Naming conventions
45 /// ------------------
46 /// - All naming conventions and equations are per Analytical Mechanics of
47 /// Space Systems by Schaub and Junkins
48 /// - The omega vector is sometimes denoted by w and assumed frame consistent with
49 /// the attitude representation. For example, if a DCM represents the
50 /// relative orientation between two frames [BN] the omega vector is
51 /// assumed to be w_(B/N)
52 /// - The naming convention for all vectors is:
53 /// variablename_obj1_obj2__frameN (note the double underscore preceding frame)
54 /// and reads back as:
55 /// "variable variablename representing the relationship between obj1 and obj2
56 /// represented in frameN"
57 /// - Unless otherwise noted angles are in RADIANS
58 ///
59 /// NOTE: All attitude representations, including DCMs, are assumed to represent a
60 /// three dimensional, cartesian coordinate system because that is what they
61 /// are used, and in many cases defined, for.
62 template<typename T>
63 class Euler321 : public CartesianVector<T, 3> {
64 public:
65 /// @brief Default constructor generates Euler 321 sequence as all zeroes
66 Euler321<T>();
67
68 /// Constructor for Matrix class with initialization.
69 /// Initializes matrix to values passed in via array
70 Euler321<T>(const T(&initial)[3]);
71
72 /// Copy constructor for Matrix class. Copies data from matrix
73 /// object to the current instance
74 Euler321<T>(const Euler321<T> &initial);
75
76 /// Constructor for Matrix class with initialization.
77 /// Initializes matrix to values passed in via array
78 Euler321<T>(const std::array<T, 3> &initial);
79
80 /// Destructor -- doesn't do anything because we don't dynamically
81 /// allocate
82 ~Euler321<T>() {}
83
84 /// ---------------------------------------------------------------------------
85 /// Dynamics and kinematics functions
86 /// ---------------------------------------------------------------------------
87 /// @brief Function to calculate the rate of change in the current
88 /// representation based on the omega vector
89 /// @param omega_f1_f2__f1 Angular velocity vector for the attitude
90 /// representation - omega of frame 1 wrt frame 2
91 /// expressed in frame 1
92 /// @param eulerdot_f1_f2 Rate of change in the current Euler angles
93 /// @return Integer error code corresponding to errors in clockwerkerrors.h
94 int rate(const CartesianVector<T, 3> &omega_f1_f2__f1, Matrix<T, 3, 1> &eulerdot_f1_f2) const;
95
96 /// ---------------------------------------------------------------------------
97 /// Conversions to other attitude representations
98 /// ---------------------------------------------------------------------------
99 /// @brief Function to convert current attitude to DCM
100 /// @param PBR return of DCM in PBR case
101 /// @return Integer error code corresponding to errors in clockwerkerrors.h
102 void toDCM(DCM<T> &dcm_f1_f2) const;
103 DCM<T> toDCM() const;
104 };
105
106 template <typename T>
107 Euler321<T>::Euler321() :
108 CartesianVector<T, 3>() {}
109
110 template <typename T>
111 Euler321<T>::Euler321(const T(&initial)[3]) :
112 CartesianVector<T, 3>(initial) {}
113
114 template <typename T>
115 Euler321<T>::Euler321(const Euler321<T> &initial) :
116 CartesianVector<T, 3>(initial) {}
117
118 template <typename T>
119 Euler321<T>::Euler321(const std::array<T, 3> &initial) :
120 CartesianVector<T, 3>(initial) {}
121
122 template <typename T>
123 int Euler321<T>::rate(const CartesianVector<T, 3> &omega_f1_f2__f1, Matrix<T, 3, 1> &eulerdot_f1_f2) const {
124 T multiplier;
125 int err = safeDivide(1, cos(CartesianVector<T, 3>::values[1][0]), multiplier);
126 if(err) {
127 return err;
128 }
129 eulerdot_f1_f2.values[0][0] = sin(CartesianVector<T, 3>::values[2][0])*omega_f1_f2__f1.values[1][0] + cos(CartesianVector<T, 3>::values[2][0])*omega_f1_f2__f1.values[2][0];
130 eulerdot_f1_f2.values[0][1] = cos(CartesianVector<T, 3>::values[2][0])*cos(CartesianVector<T, 3>::values[1][0])*omega_f1_f2__f1.values[1][0] - sin(CartesianVector<T, 3>::values[2][0])*cos(CartesianVector<T, 3>::values[1][0])*omega_f1_f2__f1.values[2][0];
131 eulerdot_f1_f2.values[0][2] = cos(CartesianVector<T, 3>::values[1][0])*omega_f1_f2__f1.values[0][0] + sin(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[1][0])*omega_f1_f2__f1.values[1][0] + cos(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[1][0])*omega_f1_f2__f1.values[2][0];
132 multiply(multiplier, eulerdot_f1_f2, eulerdot_f1_f2);
133
134 return NO_ERROR;
135 }
136
137 template <typename T>
138 void Euler321<T>::toDCM(DCM<T> &dcm_f1_f2) const {
139 dcm_f1_f2.values[0][0] = cos(CartesianVector<T, 3>::values[1][0])*cos(CartesianVector<T, 3>::values[0][0]);
140 dcm_f1_f2.values[0][1] = cos(CartesianVector<T, 3>::values[1][0])*sin(CartesianVector<T, 3>::values[0][0]);
141 dcm_f1_f2.values[0][2] = -sin(CartesianVector<T, 3>::values[1][0]);
142 dcm_f1_f2.values[1][0] = sin(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[1][0])*cos(CartesianVector<T, 3>::values[0][0]) - cos(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[0][0]);
143 dcm_f1_f2.values[1][1] = sin(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[1][0])*sin(CartesianVector<T, 3>::values[0][0]) + cos(CartesianVector<T, 3>::values[2][0])*cos(CartesianVector<T, 3>::values[0][0]);
144 dcm_f1_f2.values[1][2] = sin(CartesianVector<T, 3>::values[2][0])*cos(CartesianVector<T, 3>::values[1][0]);
145 dcm_f1_f2.values[2][0] = cos(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[1][0])*cos(CartesianVector<T, 3>::values[0][0]) + sin(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[0][0]);
146 dcm_f1_f2.values[2][1] = cos(CartesianVector<T, 3>::values[2][0])*sin(CartesianVector<T, 3>::values[1][0])*sin(CartesianVector<T, 3>::values[0][0]) - sin(CartesianVector<T, 3>::values[2][0])*cos(CartesianVector<T, 3>::values[0][0]);
147 dcm_f1_f2.values[2][2] = cos(CartesianVector<T, 3>::values[2][0])*cos(CartesianVector<T, 3>::values[1][0]);
148 }
149 template <typename T>
150 DCM<T> Euler321<T>::toDCM() const {
151 DCM<T> tmp;
152 toDCM(tmp);
153 return tmp;
154 }
155
156}
157
158#endif
Standard vector class derived from Matrix.
Definition CartesianVector.hpp:39
Class defining a direction cosine matrix inherited from Matrix.
Definition DCM.hpp:71
void toDCM(DCM< T > &dcm_f1_f2) const
Function to convert current attitude to DCM.
Definition Euler321.hpp:138
int rate(const CartesianVector< T, 3 > &omega_f1_f2__f1, Matrix< T, 3, 1 > &eulerdot_f1_f2) const
Function to calculate the rate of change in the current representation based on the omega vector.
Definition Euler321.hpp:123
Matrix math implementation.
Definition Matrix.hpp:54
#define NO_ERROR
Definition clockwerkerrors.h:31