JSBSim Flight Dynamics Model 1.3.0 (09 Apr 2026)
An Open Source Flight Dynamics and Control Software Library in C++
Loading...
Searching...
No Matches
FGPropagate.cpp
1/*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2
3 Module: FGPropagate.cpp
4 Author: Jon S. Berndt
5 Date started: 01/05/99
6 Purpose: Integrate the EOM to determine instantaneous position
7 Called by: FGFDMExec
8
9 ------------- Copyright (C) 1999 Jon S. Berndt (jon@jsbsim.org) -------------
10
11 This program is free software; you can redistribute it and/or modify it under
12 the terms of the GNU Lesser General Public License as published by the Free
13 Software Foundation; either version 2 of the License, or (at your option) any
14 later version.
15
16 This program is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18 FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
19 details.
20
21 You should have received a copy of the GNU Lesser General Public License along
22 with this program; if not, write to the Free Software Foundation, Inc., 59
23 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
24
25 Further information about the GNU Lesser General Public License can also be
26 found on the world wide web at http://www.gnu.org.
27
28FUNCTIONAL DESCRIPTION
29--------------------------------------------------------------------------------
30This class encapsulates the integration of rates and accelerations to get the
31current position of the aircraft.
32
33HISTORY
34--------------------------------------------------------------------------------
3501/05/99 JSB Created
36
37%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
38COMMENTS, REFERENCES, and NOTES
39%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
40[1] Cooke, Zyda, Pratt, and McGhee, "NPSNET: Flight Simulation Dynamic Modeling
41 Using Quaternions", Presence, Vol. 1, No. 4, pp. 404-420 Naval Postgraduate
42 School, January 1994
43[2] D. M. Henderson, "Euler Angles, Quaternions, and Transformation Matrices",
44 JSC 12960, July 1977
45[3] Richard E. McFarland, "A Standard Kinematic Model for Flight Simulation at
46 NASA-Ames", NASA CR-2497, January 1975
47[4] Barnes W. McCormick, "Aerodynamics, Aeronautics, and Flight Mechanics",
48 Wiley & Sons, 1979 ISBN 0-471-03032-5
49[5] Bernard Etkin, "Dynamics of Flight, Stability and Control", Wiley & Sons,
50 1982 ISBN 0-471-08936-2
51[6] S. Buss, "Accurate and Efficient Simulation of Rigid Body Rotations",
52 Technical Report, Department of Mathematics, University of California,
53 San Diego, 1999
54[7] Barker L.E., Bowles R.L. and Williams L.H., "Development and Application of
55 a Local Linearization Algorithm for the Integration of Quaternion Rate
56 Equations in Real-Time Flight Simulation Problems", NASA TN D-7347,
57 December 1973
58[8] Phillips W.F, Hailey C.E and Gebert G.A, "Review of Attitude Representations
59 Used for Aircraft Kinematics", Journal Of Aircraft Vol. 38, No. 4,
60 July-August 2001
61
62%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
63INCLUDES
64%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
65
66#include <iomanip>
67
68#include "FGPropagate.h"
69#include "initialization/FGInitialCondition.h"
70#include "FGFDMExec.h"
71#include "simgear/io/iostreams/sgstream.hxx"
72#include "FGInertial.h"
73#include "input_output/FGLog.h"
74
75using namespace std;
76
77namespace JSBSim {
78
79/*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
80CLASS IMPLEMENTATION
81%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
82
84 : FGModel(fdmex)
85{
86 Name = "FGPropagate";
87
88 Inertial = FDMExec->GetInertial();
89
91 // eNone = 0, eRectEuler, eTrapezoidal, eAdamsBashforth2, eAdamsBashforth3, eAdamsBashforth4};
92
93 integrator_rotational_rate = eRectEuler;
94 integrator_translational_rate = eAdamsBashforth2;
95 integrator_rotational_position = eRectEuler;
96 integrator_translational_position = eAdamsBashforth3;
97
98 VState.dqPQRidot.resize(5, FGColumnVector3(0.0,0.0,0.0));
99 VState.dqUVWidot.resize(5, FGColumnVector3(0.0,0.0,0.0));
100 VState.dqInertialVelocity.resize(5, FGColumnVector3(0.0,0.0,0.0));
101 VState.dqQtrndot.resize(5, FGQuaternion(0.0,0.0,0.0));
102
103 epa = 0.0;
104
105 bind();
106 Debug(0);
107}
108
109//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
110
112{
113 Debug(1);
114}
115
116//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
117
119{
120 if (!FGModel::InitModel()) return false;
121
122 // For initialization ONLY:
123 VState.vLocation.SetEllipse(in.SemiMajor, in.SemiMinor);
124 Inertial->SetAltitudeAGL(VState.vLocation, 4.0);
125
126 VState.dqPQRidot.resize(5, FGColumnVector3(0.0,0.0,0.0));
127 VState.dqUVWidot.resize(5, FGColumnVector3(0.0,0.0,0.0));
128 VState.dqInertialVelocity.resize(5, FGColumnVector3(0.0,0.0,0.0));
129 VState.dqQtrndot.resize(5, FGColumnVector3(0.0,0.0,0.0));
130
131 integrator_rotational_rate = eRectEuler;
132 integrator_translational_rate = eAdamsBashforth2;
133 integrator_rotational_position = eRectEuler;
134 integrator_translational_position = eAdamsBashforth3;
135
136 epa = 0.0;
137
138 return true;
139}
140
141//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
142
143void FGPropagate::SetInitialState(const FGInitialCondition* FGIC)
144{
145 // Initialize the State Vector elements and the transformation matrices
146
147 // Set the position lat/lon/radius
148 VState.vLocation = FGIC->GetPosition();
149
150 epa = FGIC->GetEarthPositionAngleIC();
151 Ti2ec = { cos(epa), sin(epa), 0.0,
152 -sin(epa), cos(epa), 0.0,
153 0.0, 0.0, 1.0 };
154 Tec2i = Ti2ec.Transposed(); // ECEF to ECI frame transform
155
156 VState.vInertialPosition = Tec2i * VState.vLocation;
157
158 UpdateLocationMatrices();
159
160 // Set the orientation from the euler angles (is normalized within the
161 // constructor). The Euler angles represent the orientation of the body
162 // frame relative to the local frame.
163 VState.qAttitudeLocal = FGIC->GetOrientation();
164
165 VState.qAttitudeECI = Ti2l.GetQuaternion()*VState.qAttitudeLocal;
166 UpdateBodyMatrices();
167
168 // Set the velocities in the instantaneus body frame
169 VState.vUVW = FGIC->GetUVWFpsIC();
170
171 // Compute the local frame ECEF velocity
172 vVel = Tb2l * VState.vUVW;
173
174 // Compute local terrain velocity
175 RecomputeLocalTerrainVelocity();
176
177 // Set the angular velocities of the body frame relative to the ECEF frame,
178 // expressed in the body frame.
179 VState.vPQR = FGIC->GetPQRRadpsIC();
180
181 VState.vPQRi = VState.vPQR + Ti2b * in.vOmegaPlanet;
182
183 CalculateInertialVelocity(); // Translational position derivative
184 CalculateQuatdot(); // Angular orientation derivative
185}
186
187//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
188// Initialize the past value deques
189
190void FGPropagate::InitializeDerivatives()
191{
192 VState.dqPQRidot.assign(5, in.vPQRidot);
193 VState.dqUVWidot.assign(5, in.vUVWidot);
194 VState.dqInertialVelocity.assign(5, VState.vInertialVelocity);
195 VState.dqQtrndot.assign(5, VState.vQtrndot);
196}
197
198//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
199/*
200Purpose: Called on a schedule to perform EOM integration
201Notes: [JB] Run in standalone mode, SeaLevelRadius will be reference radius.
202 In FGFS, SeaLevelRadius is stuffed from FGJSBSim in JSBSim.cxx each pass.
203
204At the top of this Run() function, see several "shortcuts" (or, aliases) being
205set up for use later, rather than using the longer class->function() notation.
206
207This propagation is done using the current state values
208and current derivatives. Based on these values we compute an approximation to the
209state values for (now + dt).
210
211In the code below, variables named beginning with a small "v" refer to a
212a column vector, variables named beginning with a "T" refer to a transformation
213matrix. ECEF refers to Earth Centered Earth Fixed. ECI refers to Earth Centered
214Inertial.
215
216*/
217
218bool FGPropagate::Run(bool Holding)
219{
220 if (FGModel::Run(Holding)) return true; // Fast return if we have nothing to do ...
221 if (Holding) return false;
222
223 double dt = in.DeltaT * rate; // The 'stepsize'
224
225 // Propagate rotational / translational velocity, angular /translational position, respectively.
226
227 if (!FDMExec->IntegrationSuspended()) {
228 Integrate(VState.qAttitudeECI, VState.vQtrndot, VState.dqQtrndot, dt, integrator_rotational_position);
229 Integrate(VState.vPQRi, in.vPQRidot, VState.dqPQRidot, dt, integrator_rotational_rate);
230 Integrate(VState.vInertialPosition, VState.vInertialVelocity, VState.dqInertialVelocity, dt, integrator_translational_position);
231 Integrate(VState.vInertialVelocity, in.vUVWidot, VState.dqUVWidot, dt, integrator_translational_rate);
232 }
233
234 // CAUTION : the order of the operations below is very important to get
235 // transformation matrices that are consistent with the new state of the
236 // vehicle
237
238 // 1. Update the Earth position angle (EPA)
239 epa += in.vOmegaPlanet(eZ)*dt;
240
241 // 2. Update the Ti2ec and Tec2i transforms from the updated EPA
242 double cos_epa = cos(epa);
243 double sin_epa = sin(epa);
244 Ti2ec = { cos_epa, sin_epa, 0.0,
245 -sin_epa, cos_epa, 0.0,
246 0.0, 0.0, 1.0 };
247 Tec2i = Ti2ec.Transposed(); // ECEF to ECI frame transform
248
249 // 3. Update the location from the updated Ti2ec and inertial position
250 VState.vLocation = Ti2ec*VState.vInertialPosition;
251
252 // 4. Update the other "Location-based" transformation matrices from the
253 // updated vLocation vector.
254 UpdateLocationMatrices();
255
256 // 5. Update the "Orientation-based" transformation matrices from the updated
257 // orientation quaternion and vLocation vector.
258 UpdateBodyMatrices();
259
260 // Translational position derivative (velocities are integrated in the
261 // inertial frame)
262 CalculateUVW();
263
264 // Set auxilliary state variables
265 RecomputeLocalTerrainVelocity();
266
267 VState.vPQR = VState.vPQRi - Ti2b * in.vOmegaPlanet;
268
269 // Angular orientation derivative
270 CalculateQuatdot();
271
272 VState.qAttitudeLocal = Tl2b.GetQuaternion();
273
274 // Compute vehicle velocity wrt ECEF frame, expressed in Local horizontal
275 // frame.
276 vVel = Tb2l * VState.vUVW;
277
278 // Compute orbital parameters in the inertial frame
279 ComputeOrbitalParameters();
280
281 Debug(2);
282 return false;
283}
284
285//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
286
288{
289 if (hd) {
290 VState.vUVW.InitMatrix();
291 CalculateInertialVelocity();
292 VState.vPQR.InitMatrix();
293 VState.vPQRi = Ti2b * in.vOmegaPlanet;
294 CalculateQuatdot();
295 InitializeDerivatives();
296 }
297}
298
299//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
300// Compute the quaternion orientation derivative
301//
302// vQtrndot is the quaternion derivative.
303// Reference: See Stevens and Lewis, "Aircraft Control and Simulation",
304// Second edition (2004), eqn 1.5-16b (page 50)
305
306void FGPropagate::CalculateQuatdot(void)
307{
308 // Compute quaternion orientation derivative on current body rates
309 VState.vQtrndot = VState.qAttitudeECI.GetQDot(VState.vPQRi);
310}
311
312//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
313 // Transform the velocity vector of the body relative to the origin (Earth
314 // center) to be expressed in the inertial frame, and add the vehicle velocity
315 // contribution due to the rotation of the planet.
316 // Reference: See Stevens and Lewis, "Aircraft Control and Simulation",
317 // Second edition (2004), eqn 1.5-16c (page 50)
318
319void FGPropagate::CalculateInertialVelocity(void)
320{
321 VState.vInertialVelocity = Tb2i * VState.vUVW + (in.vOmegaPlanet * VState.vInertialPosition);
322}
323
324//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
325 // Transform the velocity vector of the inertial frame to be expressed in the
326 // body frame relative to the origin (Earth center), and substract the vehicle
327 // velocity contribution due to the rotation of the planet.
328
329void FGPropagate::CalculateUVW(void)
330{
331 VState.vUVW = Ti2b * (VState.vInertialVelocity - (in.vOmegaPlanet * VState.vInertialPosition));
332}
333
334//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
335
336void FGPropagate::Integrate( FGColumnVector3& Integrand,
337 FGColumnVector3& Val,
338 deque <FGColumnVector3>& ValDot,
339 double dt,
340 eIntegrateType integration_type)
341{
342 ValDot.push_front(Val);
343 ValDot.pop_back();
344
345 switch(integration_type) {
346 case eRectEuler: Integrand += dt*ValDot[0];
347 break;
348 case eTrapezoidal: Integrand += 0.5*dt*(ValDot[0] + ValDot[1]);
349 break;
350 case eAdamsBashforth2: Integrand += dt*(1.5*ValDot[0] - 0.5*ValDot[1]);
351 break;
352 case eAdamsBashforth3: Integrand += (1/12.0)*dt*(23.0*ValDot[0] - 16.0*ValDot[1] + 5.0*ValDot[2]);
353 break;
354 case eAdamsBashforth4: Integrand += (1/24.0)*dt*(55.0*ValDot[0] - 59.0*ValDot[1] + 37.0*ValDot[2] - 9.0*ValDot[3]);
355 break;
356 case eAdamsBashforth5: Integrand += dt*((1901./720.)*ValDot[0] - (1387./360.)*ValDot[1] + (109./30.)*ValDot[2] - (637./360.)*ValDot[3] + (251./720.)*ValDot[4]);
357 break;
358 case eNone: // do nothing, freeze translational rate
359 break;
360 case eBuss1:
361 case eBuss2:
362 case eLocalLinearization:
363 {
364 LogException err;
365 err << "Can only use Buss (1 & 2) or local linearization integration methods in for rotational position!";
366 throw err;
367 }
368 default:
369 break;
370 }
371}
372
373//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
374
375void FGPropagate::Integrate( FGQuaternion& Integrand,
376 FGQuaternion& Val,
377 deque <FGQuaternion>& ValDot,
378 double dt,
379 eIntegrateType integration_type)
380{
381 ValDot.push_front(Val);
382 ValDot.pop_back();
383
384 switch(integration_type) {
385 case eRectEuler: Integrand += dt*ValDot[0];
386 break;
387 case eTrapezoidal: Integrand += 0.5*dt*(ValDot[0] + ValDot[1]);
388 break;
389 case eAdamsBashforth2: Integrand += dt*(1.5*ValDot[0] - 0.5*ValDot[1]);
390 break;
391 case eAdamsBashforth3: Integrand += (1/12.0)*dt*(23.0*ValDot[0] - 16.0*ValDot[1] + 5.0*ValDot[2]);
392 break;
393 case eAdamsBashforth4: Integrand += (1/24.0)*dt*(55.0*ValDot[0] - 59.0*ValDot[1] + 37.0*ValDot[2] - 9.0*ValDot[3]);
394 break;
395 case eAdamsBashforth5: Integrand += dt*((1901./720.)*ValDot[0] - (1387./360.)*ValDot[1] + (109./30.)*ValDot[2] - (637./360.)*ValDot[3] + (251./720.)*ValDot[4]);
396 break;
397 case eBuss1:
398 {
399 // This is the first order method as described in Samuel R. Buss paper[6].
400 // The formula from Buss' paper is transposed below to quaternions and is
401 // actually the exact solution of the quaternion differential equation
402 // qdot = 1/2*w*q when w is constant.
403 Integrand = Integrand * QExp(0.5 * dt * VState.vPQRi);
404 }
405 return; // No need to normalize since the quaternion exponential is always normal
406 case eBuss2:
407 {
408 // This is the 'augmented second-order method' from S.R. Buss paper [6].
409 // Unlike Runge-Kutta or Adams-Bashforth, it is a one-pass second-order
410 // method (see reference [6]).
411 FGColumnVector3 wi = VState.vPQRi;
412 FGColumnVector3 wdoti = in.vPQRidot;
413 FGColumnVector3 omega = wi + 0.5*dt*wdoti + dt*dt/12.*wdoti*wi;
414 Integrand = Integrand * QExp(0.5 * dt * omega);
415 }
416 return; // No need to normalize since the quaternion exponential is always normal
417 case eLocalLinearization:
418 {
419 // This is the local linearization algorithm of Barker et al. (see ref. [7])
420 // It is also a one-pass second-order method. The code below is based on the
421 // more compact formulation issued from equation (107) of ref. [8]. The
422 // constants C1, C2, C3 and C4 have the same value than those in ref. [7] pp. 11
423 FGColumnVector3 wi = 0.5 * VState.vPQRi;
424 FGColumnVector3 wdoti = 0.5 * in.vPQRidot;
425 double omegak2 = DotProduct(VState.vPQRi, VState.vPQRi);
426 double omegak = omegak2 > 1E-6 ? sqrt(omegak2) : 1E-6;
427 double rhok = 0.5 * dt * omegak;
428 double C1 = cos(rhok);
429 double C2 = 2.0 * sin(rhok) / omegak;
430 double C3 = 4.0 * (1.0 - C1) / (omegak*omegak);
431 double C4 = 4.0 * (dt - C2) / (omegak*omegak);
432 FGColumnVector3 Omega = C2*wi + C3*wdoti + C4*wi*wdoti;
433 FGQuaternion q;
434
435 q(1) = C1 - C4*DotProduct(wi, wdoti);
436 q(2) = Omega(eP);
437 q(3) = Omega(eQ);
438 q(4) = Omega(eR);
439
440 Integrand = Integrand * q;
441
442 /* Cross check with ref. [7] pp.11-12 formulas and code pp. 20
443 double pk = VState.vPQRi(eP);
444 double qk = VState.vPQRi(eQ);
445 double rk = VState.vPQRi(eR);
446 double pdotk = in.vPQRidot(eP);
447 double qdotk = in.vPQRidot(eQ);
448 double rdotk = in.vPQRidot(eR);
449 double Ap = -0.25 * (pk*pdotk + qk*qdotk + rk*rdotk);
450 double Bp = 0.25 * (pk*qdotk - qk*pdotk);
451 double Cp = 0.25 * (pdotk*rk - pk*rdotk);
452 double Dp = 0.25 * (qk*rdotk - qdotk*rk);
453 double C2p = sin(rhok) / omegak;
454 double C3p = 2.0 * (1.0 - cos(rhok)) / (omegak*omegak);
455 double H = C1 + C4 * Ap;
456 double G = -C2p*rk - C3p*rdotk + C4*Bp;
457 double J = C2p*qk + C3p*qdotk - C4*Cp;
458 double K = C2p*pk + C3p*pdotk - C4*Dp;
459
460 FGLogging log(LogLevel::INFO);
461 log << "q: " << q << "\n";
462
463 // Warning! In the paper of Barker et al. the quaternion components are not
464 // ordered the same way as in JSBSim (see equations (2) and (3) of ref. [7]
465 // as well as the comment just below equation (3))
466 log << "FORTRAN: " << H << " , " << K << " , " << J << " , " << -G << "\n";*/
467 }
468 break; // The quaternion q is not normal so the normalization needs to be done.
469 case eNone: // do nothing, freeze rotational rate
470 break;
471 default:
472 break;
473 }
474
475 Integrand.Normalize();
476}
477
478//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
479
480void FGPropagate::UpdateLocationMatrices(void)
481{
482 Tl2ec = VState.vLocation.GetTl2ec(); // local to ECEF transform
483 Tec2l = Tl2ec.Transposed(); // ECEF to local frame transform
484 Ti2l = Tec2l * Ti2ec; // ECI to local frame transform
485 Tl2i = Ti2l.Transposed(); // local to ECI transform
486}
487
488//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
489
490void FGPropagate::UpdateBodyMatrices(void)
491{
492 Ti2b = VState.qAttitudeECI.GetT(); // ECI to body frame transform
493 Tb2i = Ti2b.Transposed(); // body to ECI frame transform
494 Tl2b = Ti2b * Tl2i; // local to body frame transform
495 Tb2l = Tl2b.Transposed(); // body to local frame transform
496 Tec2b = Ti2b * Tec2i; // ECEF to body frame transform
497 Tb2ec = Tec2b.Transposed(); // body to ECEF frame tranform
498
499 Qec2b = Tec2b.GetQuaternion();
500}
501
502//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
503
504void FGPropagate::ComputeOrbitalParameters(void)
505{
506 const FGColumnVector3 Z{0., 0., 1.};
507 FGColumnVector3 R = VState.vInertialPosition;
508 FGColumnVector3 angularMomentum = R * VState.vInertialVelocity;
509 h = angularMomentum.Magnitude();
510 Inclination = acos(angularMomentum(eZ)/h)*radtodeg;
511 FGColumnVector3 N;
512 if (abs(Inclination) > 1E-8) {
513 N = Z * angularMomentum;
514 RightAscension = atan2(N(eY), N(eX))*radtodeg;
515 N.Normalize();
516 }
517 else {
518 RightAscension = 0.0;
519 N = {1., 0., 0.};
520 PerigeeArgument = 0.0;
521 }
522 R.Normalize();
523 double vr = DotProduct(R, VState.vInertialVelocity);
524 FGColumnVector3 eVector = (VState.vInertialVelocity*angularMomentum/in.GM - R);
525 Eccentricity = eVector.Magnitude();
526 if (Eccentricity > 1E-8) {
527 eVector /= Eccentricity;
528 if (abs(Inclination) > 1E-8) {
529 PerigeeArgument = acos(DotProduct(N, eVector)) * radtodeg;
530 if (eVector(eZ) < 0) PerigeeArgument = 360. - PerigeeArgument;
531 }
532 }
533 else
534 {
535 eVector = {1., 0., 0.};
536 PerigeeArgument = 0.0;
537 }
538
539 TrueAnomaly = acos(Constrain(-1.0, DotProduct(eVector, R), 1.0)) * radtodeg;
540 if (vr < 0.0) TrueAnomaly = 360. - TrueAnomaly;
541 ApoapsisRadius = h*h / (in.GM * (1-Eccentricity));
542 PeriapsisRadius = h*h / (in.GM * (1+Eccentricity));
543
544 if (Eccentricity < 1.0) {
545 double semimajor = 0.5*(ApoapsisRadius + PeriapsisRadius);
546 OrbitalPeriod = 2.*M_PI*pow(semimajor, 1.5)/sqrt(in.GM);
547 }
548 else
549 OrbitalPeriod = 0.0;
550}
551
552//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
553
554void FGPropagate::SetInertialOrientation(const FGQuaternion& Qi)
555{
556 VState.qAttitudeECI = Qi;
557 VState.qAttitudeECI.Normalize();
558 UpdateBodyMatrices();
559 VState.qAttitudeLocal = Tl2b.GetQuaternion();
560 CalculateQuatdot();
561}
562
563//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
564
565void FGPropagate::SetInertialVelocity(const FGColumnVector3& Vi) {
566 VState.vInertialVelocity = Vi;
567 CalculateUVW();
568 vVel = Tb2l * VState.vUVW;
569}
570
571//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
572
573void FGPropagate::SetInertialRates(const FGColumnVector3& vRates) {
574 VState.vPQRi = Ti2b * vRates;
575 VState.vPQR = VState.vPQRi - Ti2b * in.vOmegaPlanet;
576 CalculateQuatdot();
577}
578
579//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
580
582{
583 return VState.vLocation.GetRadius() - VState.vLocation.GetSeaLevelRadius();
584}
585
586//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
587
588void FGPropagate::SetAltitudeASL(double altASL)
589{
590 double slr = VState.vLocation.GetSeaLevelRadius();
591 VState.vLocation.SetRadius(slr + altASL);
592 UpdateVehicleState();
593}
594
595//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
596
597void FGPropagate::RecomputeLocalTerrainVelocity()
598{
599 FGLocation contact;
600 FGColumnVector3 normal;
601 Inertial->GetContactPoint(VState.vLocation, contact, normal,
602 LocalTerrainVelocity, LocalTerrainAngularVelocity);
603}
604
605//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
606
607double FGPropagate::GetTerrainElevation(void) const
608{
609 FGColumnVector3 vDummy;
610 FGLocation contact;
611 contact.SetEllipse(in.SemiMajor, in.SemiMinor);
612 Inertial->GetContactPoint(VState.vLocation, contact, vDummy, vDummy, vDummy);
613 return contact.GetGeodAltitude();
614}
615
616//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
617
618void FGPropagate::SetTerrainElevation(double terrainElev)
619{
620 Inertial->SetTerrainElevation(terrainElev);
621}
622
623//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
624
626{
627 FGLocation contact;
628 FGColumnVector3 vDummy;
629 Inertial->GetContactPoint(VState.vLocation, contact, vDummy, vDummy, vDummy);
630 return contact.GetRadius();
631}
632
633//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
634
635double FGPropagate::GetDistanceAGL(void) const
636{
637 return Inertial->GetAltitudeAGL(VState.vLocation);
638}
639
640//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
641
642double FGPropagate::GetDistanceAGLKm(void) const
643{
644 return GetDistanceAGL()*0.0003048;
645}
646
647//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
648
649void FGPropagate::SetDistanceAGL(double tt)
650{
651 Inertial->SetAltitudeAGL(VState.vLocation, tt);
652 UpdateVehicleState();
653}
654
655//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
656
657void FGPropagate::SetDistanceAGLKm(double tt)
658{
659 SetDistanceAGL(tt*3280.8399);
660}
661
662//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
663
664void FGPropagate::SetVState(const VehicleState& vstate)
665{
666 //ToDo: Shouldn't all of these be set from the vstate vector passed in?
667 VState.vLocation = vstate.vLocation;
668 UpdateLocationMatrices();
669 SetInertialOrientation(vstate.qAttitudeECI);
670 RecomputeLocalTerrainVelocity();
671 VState.vUVW = vstate.vUVW;
672 vVel = Tb2l * VState.vUVW;
673 VState.vPQR = vstate.vPQR;
674 VState.vPQRi = VState.vPQR + Ti2b * in.vOmegaPlanet;
675 VState.vInertialPosition = vstate.vInertialPosition;
676 CalculateQuatdot();
677}
678
679//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
680
681void FGPropagate::UpdateVehicleState(void)
682{
683 RecomputeLocalTerrainVelocity();
684 VState.vInertialPosition = Tec2i * VState.vLocation;
685 UpdateLocationMatrices();
686 UpdateBodyMatrices();
687 vVel = Tb2l * VState.vUVW;
688 VState.qAttitudeLocal = Tl2b.GetQuaternion();
689}
690
691//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
692
693void FGPropagate::SetLocation(const FGLocation& l)
694{
695 VState.vLocation = l;
696 UpdateVehicleState();
697}
698
699//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
700
702{
703 return VState.qAttitudeLocal.GetEuler() * radtodeg;
704}
705
706//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
707
708void FGPropagate::DumpState(void)
709{
710 FGLogging out(LogLevel::STDOUT);
711 out << "\n";
712 out << LogFormat::BLUE
713 << "------------------------------------------------------------------" << LogFormat::RESET << "\n";
714 out << LogFormat::BOLD << fixed
715 << "State Report at sim time: " << FDMExec->GetSimTime() << " seconds" << LogFormat::RESET << "\n";
716 out << " " << LogFormat::UNDERLINE_ON
717 << "Position" << LogFormat::UNDERLINE_OFF << "\n";
718 out << " ECI: " << VState.vInertialPosition.Dump(", ") << " (x,y,z, in ft)\n";
719 out << " ECEF: " << VState.vLocation << " (x,y,z, in ft)\n";
720 out << " Local: " << VState.vLocation.GetGeodLatitudeDeg()
721 << ", " << VState.vLocation.GetLongitudeDeg()
722 << ", " << GetAltitudeASL() << " (geodetic lat, lon, alt ASL in deg and ft)\n";
723
724 out << "\n " << LogFormat::UNDERLINE_ON
725 << "Orientation" << LogFormat::UNDERLINE_OFF << "\n";
726 out << " ECI: " << VState.qAttitudeECI.GetEulerDeg().Dump(", ") << " (phi, theta, psi in deg)\n";
727 out << " Local: " << VState.qAttitudeLocal.GetEulerDeg().Dump(", ") << " (phi, theta, psi in deg)\n";
728
729 out << "\n " << LogFormat::UNDERLINE_ON
730 << "Velocity" << LogFormat::UNDERLINE_OFF << "\n";
731 out << " ECI: " << VState.vInertialVelocity.Dump(", ") << " (x,y,z in ft/s)\n";
732 out << " ECEF: " << GetECEFVelocity().Dump(", ") << " (x,y,z in ft/s)\n";
733 out << " Local: " << GetVel() << " (n,e,d in ft/sec)\n";
734 out << " Body: " << GetUVW() << " (u,v,w in ft/sec)\n";
735
736 out << "\n" << " " << LogFormat::UNDERLINE_ON
737 << "Body Rates (relative to given frame, expressed in body frame)" << LogFormat::UNDERLINE_OFF << "\n";
738 out << " ECI: " << (VState.vPQRi * radtodeg).Dump(", ") << " (p,q,r in deg/s)" << "\n";
739 out << " ECEF: " << (VState.vPQR * radtodeg).Dump(", ") << " (p,q,r in deg/s)" << "\n";
740}
741
742//******************************************************************************
743
744void FGPropagate::WriteStateFile(int num)
745{
746 sg_ofstream outfile;
747
748 if (num == 0) return;
749
750 SGPath path = FDMExec->GetOutputPath();
751
752 if (path.isNull()) path = SGPath("initfile.");
753 else path.append("initfile.");
754
755 // Append sim time to the filename since there may be more than one created during a simulation run
756 path.concat(to_string((double)FDMExec->GetSimTime())+".xml");
757
758 switch(num) {
759 case 1:
760 outfile.open(path);
761 if (outfile.is_open()) {
762 outfile << "<?xml version=\"1.0\"?>\n";
763 outfile << "<initialize name=\"reset00\">\n";
764 outfile << " <ubody unit=\"FT/SEC\"> " << VState.vUVW(eU) << " </ubody> \n";
765 outfile << " <vbody unit=\"FT/SEC\"> " << VState.vUVW(eV) << " </vbody> \n";
766 outfile << " <wbody unit=\"FT/SEC\"> " << VState.vUVW(eW) << " </wbody> \n";
767 outfile << " <phi unit=\"DEG\"> " << VState.qAttitudeLocal.GetEuler(ePhi)*radtodeg << " </phi>\n";
768 outfile << " <theta unit=\"DEG\"> " << VState.qAttitudeLocal.GetEuler(eTht)*radtodeg << " </theta>\n";
769 outfile << " <psi unit=\"DEG\"> " << VState.qAttitudeLocal.GetEuler(ePsi)*radtodeg << " </psi>\n";
770 outfile << " <longitude unit=\"DEG\"> " << VState.vLocation.GetLongitudeDeg() << " </longitude>\n";
771 outfile << " <latitude unit=\"DEG\"> " << VState.vLocation.GetLatitudeDeg() << " </latitude>\n";
772 outfile << " <altitude unit=\"FT\"> " << GetDistanceAGL() << " </altitude>\n";
773 outfile << "</initialize>\n";
774 outfile.close();
775 } else {
776 FGLogging log(LogLevel::ERROR);
777 log << "Could not open and/or write the state to the initial conditions file: "
778 << path << "\n";
779 }
780 break;
781 case 2:
782 outfile.open(path);
783 if (outfile.is_open()) {
784 outfile << "<?xml version=\"1.0\"?>\n";
785 outfile << "<initialize name=\"IC File\" version=\"2.0\">\n";
786 outfile << "\n";
787 outfile << " <position frame=\"ECEF\">\n";
788 outfile << " <latitude unit=\"DEG\" type=\"geodetic\"> " << VState.vLocation.GetGeodLatitudeDeg() << " </latitude>\n";
789 outfile << " <longitude unit=\"DEG\"> " << VState.vLocation.GetLongitudeDeg() << " </longitude>\n";
790 outfile << " <altitudeMSL unit=\"FT\"> " << GetAltitudeASL() << " </altitudeMSL>\n";
791 outfile << " </position>\n";
792 outfile << "\n";
793 outfile << " <orientation unit=\"DEG\" frame=\"LOCAL\">\n";
794 outfile << " <yaw> " << VState.qAttitudeLocal.GetEulerDeg(eYaw) << " </yaw>\n";
795 outfile << " <pitch> " << VState.qAttitudeLocal.GetEulerDeg(ePitch) << " </pitch>\n";
796 outfile << " <roll> " << VState.qAttitudeLocal.GetEulerDeg(eRoll) << " </roll>\n";
797 outfile << " </orientation>\n";
798 outfile << "\n";
799 outfile << " <velocity unit=\"FT/SEC\" frame=\"LOCAL\">\n";
800 outfile << " <x> " << GetVel(eNorth) << " </x>\n";
801 outfile << " <y> " << GetVel(eEast) << " </y>\n";
802 outfile << " <z> " << GetVel(eDown) << " </z>\n";
803 outfile << " </velocity>\n";
804 outfile << "\n";
805 outfile << " <attitude_rate unit=\"DEG/SEC\" frame=\"BODY\">\n";
806 outfile << " <roll> " << (VState.vPQR*radtodeg)(eRoll) << " </roll>\n";
807 outfile << " <pitch> " << (VState.vPQR*radtodeg)(ePitch) << " </pitch>\n";
808 outfile << " <yaw> " << (VState.vPQR*radtodeg)(eYaw) << " </yaw>\n";
809 outfile << " </attitude_rate>\n";
810 outfile << "\n";
811 outfile << "</initialize>\n";
812 outfile.close();
813 } else {
814 FGLogging log(LogLevel::ERROR);
815 log << "Could not open and/or write the state to the initial conditions file: "
816 << path << "\n";
817 }
818 break;
819 default:
820 FGLogging log(LogLevel::ERROR);
821 log << "When writing a state file, the supplied value must be 1 or 2 for the version number of the resulting IC file\n";
822 }
823}
824
825//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
826
827void FGPropagate::bind(void)
828{
829 PropertyManager->Tie("velocities/h-dot-fps", this, &FGPropagate::Gethdot);
830
831 PropertyManager->Tie("velocities/v-north-fps", this, eNorth, &FGPropagate::GetVel);
832 PropertyManager->Tie("velocities/v-east-fps", this, eEast, &FGPropagate::GetVel);
833 PropertyManager->Tie("velocities/v-down-fps", this, eDown, &FGPropagate::GetVel);
834
835 PropertyManager->Tie("velocities/u-fps", this, eU, &FGPropagate::GetUVW);
836 PropertyManager->Tie("velocities/v-fps", this, eV, &FGPropagate::GetUVW);
837 PropertyManager->Tie("velocities/w-fps", this, eW, &FGPropagate::GetUVW);
838
839 PropertyManager->Tie("velocities/p-rad_sec", this, eP, &FGPropagate::GetPQR);
840 PropertyManager->Tie("velocities/q-rad_sec", this, eQ, &FGPropagate::GetPQR);
841 PropertyManager->Tie("velocities/r-rad_sec", this, eR, &FGPropagate::GetPQR);
842
843 PropertyManager->Tie("velocities/pi-rad_sec", this, eP, &FGPropagate::GetPQRi);
844 PropertyManager->Tie("velocities/qi-rad_sec", this, eQ, &FGPropagate::GetPQRi);
845 PropertyManager->Tie("velocities/ri-rad_sec", this, eR, &FGPropagate::GetPQRi);
846
847 PropertyManager->Tie("velocities/eci-x-fps", this, eX, &FGPropagate::GetInertialVelocity);
848 PropertyManager->Tie("velocities/eci-y-fps", this, eY, &FGPropagate::GetInertialVelocity);
849 PropertyManager->Tie("velocities/eci-z-fps", this, eZ, &FGPropagate::GetInertialVelocity);
850
851 PropertyManager->Tie("velocities/ecef-x-fps", this, eX, &FGPropagate::GetECEFVelocity);
852 PropertyManager->Tie("velocities/ecef-y-fps", this, eY, &FGPropagate::GetECEFVelocity);
853 PropertyManager->Tie("velocities/ecef-z-fps", this, eZ, &FGPropagate::GetECEFVelocity);
854
855 PropertyManager->Tie("velocities/eci-velocity-mag-fps", this, &FGPropagate::GetInertialVelocityMagnitude);
856 PropertyManager->Tie("velocities/ned-velocity-mag-fps", this, &FGPropagate::GetNEDVelocityMagnitude);
857
858 PropertyManager->Tie("position/h-sl-ft", this, &FGPropagate::GetAltitudeASL,
859 &FGPropagate::SetAltitudeASL);
860 PropertyManager->Tie("position/h-sl-meters", this, &FGPropagate::GetAltitudeASLmeters,
861 &FGPropagate::SetAltitudeASLmeters);
862 PropertyManager->Tie("position/lat-gc-rad", this, &FGPropagate::GetLatitude,
863 &FGPropagate::SetLatitude);
864 PropertyManager->Tie("position/long-gc-rad", this, &FGPropagate::GetLongitude,
865 &FGPropagate::SetLongitude);
866 PropertyManager->Tie("position/lat-gc-deg", this, &FGPropagate::GetLatitudeDeg,
867 &FGPropagate::SetLatitudeDeg);
868 PropertyManager->Tie("position/long-gc-deg", this, &FGPropagate::GetLongitudeDeg,
869 &FGPropagate::SetLongitudeDeg);
870 PropertyManager->Tie("position/lat-geod-rad", this, &FGPropagate::GetGeodLatitudeRad);
871 PropertyManager->Tie("position/lat-geod-deg", this, &FGPropagate::GetGeodLatitudeDeg);
872 PropertyManager->Tie("position/geod-alt-ft", this, &FGPropagate::GetGeodeticAltitude);
873 PropertyManager->Tie("position/h-agl-ft", this, &FGPropagate::GetDistanceAGL,
874 &FGPropagate::SetDistanceAGL);
875 PropertyManager->Tie("position/geod-alt-km", this, &FGPropagate::GetGeodeticAltitudeKm);
876 PropertyManager->Tie("position/h-agl-km", this, &FGPropagate::GetDistanceAGLKm,
877 &FGPropagate::SetDistanceAGLKm);
878 PropertyManager->Tie("position/radius-to-vehicle-ft", this, &FGPropagate::GetRadius);
879 PropertyManager->Tie("position/terrain-elevation-asl-ft", this, &FGPropagate::GetTerrainElevation,
880 &FGPropagate::SetTerrainElevation);
881
882 PropertyManager->Tie("position/eci-x-ft", this, eX, &FGPropagate::GetInertialPosition);
883 PropertyManager->Tie("position/eci-y-ft", this, eY, &FGPropagate::GetInertialPosition);
884 PropertyManager->Tie("position/eci-z-ft", this, eZ, &FGPropagate::GetInertialPosition);
885
886 PropertyManager->Tie("position/ecef-x-ft", this, eX, &FGPropagate::GetLocation);
887 PropertyManager->Tie("position/ecef-y-ft", this, eY, &FGPropagate::GetLocation);
888 PropertyManager->Tie("position/ecef-z-ft", this, eZ, &FGPropagate::GetLocation);
889
890 PropertyManager->Tie("position/epa-rad", this, &FGPropagate::GetEarthPositionAngle);
891 PropertyManager->Tie("metrics/terrain-radius", this, &FGPropagate::GetLocalTerrainRadius);
892
893 PropertyManager->Tie("attitude/phi-rad", this, ePhi, &FGPropagate::GetEuler);
894 PropertyManager->Tie("attitude/theta-rad", this, eTht, &FGPropagate::GetEuler);
895 PropertyManager->Tie("attitude/psi-rad", this, ePsi, &FGPropagate::GetEuler);
896
897 PropertyManager->Tie("attitude/phi-deg", this, ePhi, &FGPropagate::GetEulerDeg);
898 PropertyManager->Tie("attitude/theta-deg", this, eTht, &FGPropagate::GetEulerDeg);
899 PropertyManager->Tie("attitude/psi-deg", this, ePsi, &FGPropagate::GetEulerDeg);
900
901 PropertyManager->Tie("attitude/roll-rad", this, ePhi, &FGPropagate::GetEuler);
902 PropertyManager->Tie("attitude/pitch-rad", this, eTht, &FGPropagate::GetEuler);
903 PropertyManager->Tie("attitude/heading-true-rad", this, ePsi, &FGPropagate::GetEuler);
904
905 PropertyManager->Tie("orbital/specific-angular-momentum-ft2_sec", &h);
906 PropertyManager->Tie("orbital/inclination-deg", &Inclination);
907 PropertyManager->Tie("orbital/right-ascension-deg", &RightAscension);
908 PropertyManager->Tie("orbital/eccentricity", &Eccentricity);
909 PropertyManager->Tie("orbital/argument-of-perigee-deg", &PerigeeArgument);
910 PropertyManager->Tie("orbital/true-anomaly-deg", &TrueAnomaly);
911 PropertyManager->Tie("orbital/apoapsis-radius-ft", &ApoapsisRadius);
912 PropertyManager->Tie("orbital/periapsis-radius-ft", &PeriapsisRadius);
913 PropertyManager->Tie("orbital/period-sec", &OrbitalPeriod);
914
915 PropertyManager->Tie("simulation/integrator/rate/rotational", (int*)&integrator_rotational_rate);
916 PropertyManager->Tie("simulation/integrator/rate/translational", (int*)&integrator_translational_rate);
917 PropertyManager->Tie("simulation/integrator/position/rotational", (int*)&integrator_rotational_position);
918 PropertyManager->Tie("simulation/integrator/position/translational", (int*)&integrator_translational_position);
919
920 PropertyManager->Tie<FGPropagate, int>("simulation/write-state-file", this,
921 nullptr, &FGPropagate::WriteStateFile);
922}
923
924//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
925// The bitmasked value choices are as follows:
926// unset: In this case (the default) JSBSim would only print
927// out the normally expected messages, essentially echoing
928// the config files as they are read. If the environment
929// variable is not set, debug_lvl is set to 1 internally
930// 0: This requests JSBSim not to output any messages
931// whatsoever.
932// 1: This value explicity requests the normal JSBSim
933// startup messages
934// 2: This value asks for a message to be printed out when
935// a class is instantiated
936// 4: When this value is set, a message is displayed when a
937// FGModel object executes its Run() method
938// 8: When this value is set, various runtime state variables
939// are printed out periodically
940// 16: When set various parameters are sanity checked and
941// a message is printed out when they go out of bounds
942
943void FGPropagate::Debug(int from)
944{
945 if (debug_lvl <= 0) return;
946
947 if (debug_lvl & 1) { // Standard console startup message output
948 if (from == 0) { // Constructor
949
950 }
951 }
952 if (debug_lvl & 2 ) { // Instantiation/Destruction notification
953 FGLogging log(LogLevel::DEBUG);
954 if (from == 0) log << "Instantiated: FGPropagate\n";
955 if (from == 1) log << "Destroyed: FGPropagate\n";
956 }
957 if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
958 }
959 if (debug_lvl & 8 && from == 2) { // Runtime state variables
960 FGLogging log(LogLevel::DEBUG);
961 log << "\n" << LogFormat::BLUE << LogFormat::BOLD << left << fixed
962 << " Propagation Report (English units: ft, degrees) at simulation time " << FDMExec->GetSimTime() << " seconds"
963 << LogFormat::RESET << "\n\n";
964 log << LogFormat::BOLD << " Earth Position Angle (deg): " << setw(8) << setprecision(3) << LogFormat::RESET
965 << GetEarthPositionAngleDeg() << "\n\n";
966 log << LogFormat::BOLD << " Body velocity (ft/sec): " << setw(8) << setprecision(3) << LogFormat::RESET << VState.vUVW << "\n";
967 log << LogFormat::BOLD << " Local velocity (ft/sec): " << setw(8) << setprecision(3) << LogFormat::RESET << vVel << "\n";
968 log << LogFormat::BOLD << " Inertial velocity (ft/sec): " << setw(8) << setprecision(3) << LogFormat::RESET << VState.vInertialVelocity << "\n";
969 log << LogFormat::BOLD << " Inertial Position (ft): " << setw(10) << setprecision(3) << LogFormat::RESET << VState.vInertialPosition << "\n";
970 log << LogFormat::BOLD << " Latitude (deg): " << setw(8) << setprecision(3) << LogFormat::RESET << VState.vLocation.GetLatitudeDeg() << "\n";
971 log << LogFormat::BOLD << " Longitude (deg): " << setw(8) << setprecision(3) << LogFormat::RESET << VState.vLocation.GetLongitudeDeg() << "\n";
972 log << LogFormat::BOLD << " Altitude ASL (ft): " << setw(8) << setprecision(3) << LogFormat::RESET << GetAltitudeASL() << "\n";
973 // log << LogFormat::BOLD << " Acceleration (NED, ft/sec^2): " << setw(8) << setprecision(3) << LogFormat::RESET << Tb2l*GetUVWdot() << "\n";
974 log << "\n";
975 log << LogFormat::BOLD << " Matrix ECEF to Body (Orientation of Body with respect to ECEF): \n"
976 << LogFormat::RESET << Tec2b.Dump("\t", " ");
977 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
978 << setprecision(3) << LogFormat::RESET << (Tec2b.GetQuaternion().GetEuler()*radtodeg)
979 << "\n\n";
980
981 log << LogFormat::BOLD << " Matrix Body to ECEF (Orientation of ECEF with respect to Body):\n"
982 << LogFormat::RESET << Tb2ec.Dump("\t", " ");
983 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
984 << setprecision(3) << LogFormat::RESET << (Tb2ec.GetQuaternion().GetEuler()*radtodeg)
985 << "\n\n";
986
987 log << LogFormat::BOLD << " Matrix Local to Body (Orientation of Body with respect to Local):\n"
988 << LogFormat::RESET << Tl2b.Dump("\t", " ");
989 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
990 << setprecision(3) << LogFormat::RESET << (Tl2b.GetQuaternion().GetEuler()*radtodeg)
991 << "\n\n";
992
993 log << LogFormat::BOLD << " Matrix Body to Local (Orientation of Local with respect to Body):\n"
994 << LogFormat::RESET << Tb2l.Dump("\t", " ");
995 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
996 << setprecision(3) << LogFormat::RESET << (Tb2l.GetQuaternion().GetEuler()*radtodeg)
997 << "\n\n";
998
999 log << LogFormat::BOLD << " Matrix Local to ECEF (Orientation of ECEF with respect to Local):\n"
1000 << LogFormat::RESET << Tl2ec.Dump("\t", " ");
1001 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1002 << setprecision(3) << LogFormat::RESET << (Tl2ec.GetQuaternion().GetEuler()*radtodeg)
1003 << "\n\n";
1004
1005 log << LogFormat::BOLD << " Matrix ECEF to Local (Orientation of Local with respect to ECEF):\n"
1006 << LogFormat::RESET << Tec2l.Dump("\t", " ");
1007 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1008 << setprecision(3) << LogFormat::RESET << (Tec2l.GetQuaternion().GetEuler()*radtodeg)
1009 << "\n\n";
1010
1011 log << LogFormat::BOLD << " Matrix ECEF to Inertial (Orientation of Inertial with respect to ECEF):\n"
1012 << LogFormat::RESET << Tec2i.Dump("\t", " ");
1013 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1014 << setprecision(3) << LogFormat::RESET << (Tec2i.GetQuaternion().GetEuler()*radtodeg)
1015 << "\n\n";
1016
1017 log << LogFormat::BOLD << " Matrix Inertial to ECEF (Orientation of ECEF with respect to Inertial):\n"
1018 << LogFormat::RESET << Ti2ec.Dump("\t", " ");
1019 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1020 << setprecision(3) << LogFormat::RESET << (Ti2ec.GetQuaternion().GetEuler()*radtodeg)
1021 << "\n\n";
1022
1023 log << LogFormat::BOLD << " Matrix Inertial to Body (Orientation of Body with respect to Inertial):\n"
1024 << LogFormat::RESET << Ti2b.Dump("\t", " ");
1025 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1026 << setprecision(3) << LogFormat::RESET << (Ti2b.GetQuaternion().GetEuler()*radtodeg)
1027 << "\n\n";
1028
1029 log << LogFormat::BOLD << " Matrix Body to Inertial (Orientation of Inertial with respect to Body):\n"
1030 << LogFormat::RESET << Tb2i.Dump("\t", " ");
1031 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1032 << setprecision(3) << LogFormat::RESET << (Tb2i.GetQuaternion().GetEuler()*radtodeg)
1033 << "\n\n";
1034
1035 log << LogFormat::BOLD << " Matrix Inertial to Local (Orientation of Local with respect to Inertial):\n"
1036 << LogFormat::RESET << Ti2l.Dump("\t", " ");
1037 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1038 << setprecision(3) << LogFormat::RESET << (Ti2l.GetQuaternion().GetEuler()*radtodeg)
1039 << "\n\n";
1040
1041 log << LogFormat::BOLD << " Matrix Local to Inertial (Orientation of Inertial with respect to Local):\n"
1042 << LogFormat::RESET << Tl2i.Dump("\t", " ");
1043 log << LogFormat::BOLD << "\n Associated Euler angles (deg): " << setw(8)
1044 << setprecision(3) << LogFormat::RESET << (Tl2i.GetQuaternion().GetEuler()*radtodeg)
1045 << "\n\n";
1046 }
1047 if (debug_lvl & 16) { // Sanity checking
1048 if (from == 2) { // State sanity checking
1049 if (fabs(VState.vPQR.Magnitude()) > 1000.0) {
1050 LogException err;
1051 err << "Vehicle rotation rate is excessive (>1000 rad/sec): " << VState.vPQR.Magnitude() << "\n";
1052 throw err;
1053 }
1054 if (fabs(VState.vUVW.Magnitude()) > 1.0e10) {
1055 LogException err;
1056 err << "Vehicle velocity is excessive (>1e10 ft/sec): " << VState.vUVW.Magnitude() << "\n";
1057 throw err;
1058 }
1059 if (fabs(GetDistanceAGL()) > 1e10) {
1060 LogException err;
1061 err << "Vehicle altitude is excessive (>1e10 ft): " << GetDistanceAGL() << "\n";
1062 throw err;
1063 }
1064 }
1065 }
1066 if (debug_lvl & 64) {
1067 if (from == 0) { // Constructor
1068 }
1069 }
1070}
1071}
This class implements a 3 element column vector.
double Magnitude(void) const
Length of the vector.
std::string Dump(const std::string &delimeter) const
Prints the contents of the vector.
Encapsulates the JSBSim simulation executive.
Definition FGFDMExec.h:185
const SGPath & GetOutputPath(void)
Retrieves the path to the output files.
Definition FGFDMExec.h:404
bool IntegrationSuspended(void) const
Returns the simulation suspension state.
Definition FGFDMExec.h:563
double GetSimTime(void) const
Returns the cumulative simulation time in seconds.
Definition FGFDMExec.h:550
std::shared_ptr< FGInertial > GetInertial(void) const
Returns the FGInertial pointer.
Initializes the simulation run.
const FGColumnVector3 GetUVWFpsIC(void) const
Gets the initial body velocity.
const FGQuaternion & GetOrientation(void) const
Gets the initial orientation.
double GetEarthPositionAngleIC(void) const
Gets the initial Earth position angle.
const FGColumnVector3 GetPQRRadpsIC(void) const
Gets the initial body rotation rate.
const FGLocation & GetPosition(void) const
Gets the initial position.
static constexpr double Constrain(double min, double value, double max)
Constrain a value between a minimum and a maximum value.
Definition FGJSBBase.h:289
FGLocation holds an arbitrary location in the Earth centered Earth fixed reference frame (ECEF).
Definition FGLocation.h:152
void SetEllipse(double semimajor, double semiminor)
Sets the semimajor and semiminor axis lengths for this planet.
const FGMatrix33 & GetTl2ec(void) const
Transform matrix from local horizontal to earth centered frame.
Definition FGLocation.h:296
double GetGeodLatitudeDeg(void) const
Get the GEODETIC latitude in degrees.
Definition FGLocation.h:273
double GetRadius() const
Get the distance from the center of the earth in feet.
Definition FGLocation.h:291
double GetSeaLevelRadius(void) const
Get the sea level radius in feet below the current location.
double GetLatitudeDeg() const
Get the GEOCENTRIC latitude in degrees.
Definition FGLocation.h:267
void SetRadius(double radius)
Set the distance from the center of the earth.
double GetLongitudeDeg() const
Get the longitude.
Definition FGLocation.h:240
FGMatrix33 Transposed(void) const
Transposed matrix.
Definition FGMatrix33.h:221
std::string Dump(const std::string &delimeter) const
Prints the contents of the matrix.
FGQuaternion GetQuaternion(void) const
Returns the quaternion associated with this direction cosine (rotation) matrix.
Base class for all scheduled JSBSim models.
Definition FGModel.h:70
virtual bool Run(bool Holding)
Runs the model; called by the Executive.
Definition FGModel.cpp:90
double GetNEDVelocityMagnitude(void) const
Retrieves the total local NED velocity in ft/sec.
double GetInertialVelocityMagnitude(void) const
Retrieves the total inertial velocity in ft/sec.
double Gethdot(void) const
Returns the current altitude rate.
FGColumnVector3 GetEulerDeg(void) const
Retrieves the Euler angles (in degrees) that define the vehicle orientation.
const FGColumnVector3 & GetPQR(void) const
Retrieves the body angular rates vector, relative to the ECEF frame.
FGColumnVector3 GetECEFVelocity(void) const
Calculates and retrieves the velocity vector relative to the earth centered earth fixed (ECEF) frame.
double GetLocalTerrainRadius(void) const
Returns the "constant" LocalTerrainRadius.
bool Run(bool Holding)
Runs the state propagation model; called by the Executive Can pass in a value indicating if the execu...
double GetEarthPositionAngleDeg(void) const
Returns the Earth position angle in degrees.
double GetAltitudeASL(void) const
Returns the current altitude above sea level.
const FGColumnVector3 & GetPQRi(void) const
Retrieves the body angular rates vector, relative to the ECI (inertial) frame.
const FGColumnVector3 & GetEuler(void) const
Retrieves the Euler angles that define the vehicle orientation.
const FGColumnVector3 & GetUVW(void) const
Retrieves the body frame vehicle velocity vector.
~FGPropagate()
Destructor.
bool InitModel(void)
Initializes the FGPropagate class after instantiation and prior to first execution.
FGPropagate(FGFDMExec *Executive)
Constructor.
const FGColumnVector3 & GetInertialVelocity(void) const
Retrieves the inertial velocity vector in ft/sec.
const FGColumnVector3 & GetVel(void) const
Retrieves the velocity vector.
double GetEarthPositionAngle(void) const
Returns the Earth position angle.
void SetHoldDown(bool hd)
Sets the property forces/hold-down.
double GetAltitudeASLmeters(void) const
Returns the current altitude above sea level.
const FGColumnVector3 & GetInertialPosition(void) const
Retrieves the inertial position vector.
Models the Quaternion representation of rotations.
void Normalize(void)
Normalize.
double GetEulerDeg(int i) const
Retrieves the Euler angles.
FGQuaternion GetQDot(const FGColumnVector3 &PQR) const
Quaternion derivative for given angular rates.
const FGColumnVector3 & GetEuler(void) const
Retrieves the Euler angles.
const FGMatrix33 & GetT(void) const
Transformation matrix.
Main namespace for the JSBSim Flight Dynamics Model.
Definition FGFDMExec.cpp:71
double DotProduct(const FGColumnVector3 &v1, const FGColumnVector3 &v2)
Dot product of two vectors Compute and return the euclidean dot (or scalar) product of two vectors v1...
FGQuaternion QExp(const FGColumnVector3 &omega)
Quaternion exponential.
FGColumnVector3 vUVW
The velocity vector of the vehicle with respect to the ECEF frame, expressed in the body system.
FGQuaternion qAttitudeECI
The current orientation of the vehicle, that is, the orientation of the body frame relative to the in...
FGColumnVector3 vPQRi
The angular velocity vector for the vehicle body frame relative to the ECI frame, expressed in the bo...
FGColumnVector3 vPQR
The angular velocity vector for the vehicle relative to the ECEF frame, expressed in the body frame.
FGLocation vLocation
Represents the current location of the vehicle in Earth centered Earth fixed (ECEF) frame.
FGQuaternion qAttitudeLocal
The current orientation of the vehicle, that is, the orientation of the body frame relative to the lo...