FlightGear next
FGPropagate.cpp
Go to the documentation of this file.
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"
70#include "FGFDMExec.h"
71#include "simgear/io/iostreams/sgstream.hxx"
72#include "FGInertial.h"
73
74using namespace std;
75
76namespace JSBSim {
77
78/*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
79CLASS IMPLEMENTATION
80%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
81
83 : FGModel(fdmex)
84{
85 Debug(0);
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
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
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
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
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 Debug(2);
279 return false;
280}
281
282//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
283
285{
286 if (hd) {
287 VState.vUVW.InitMatrix();
288 CalculateInertialVelocity();
289 VState.vPQR.InitMatrix();
290 VState.vPQRi = Ti2b * in.vOmegaPlanet;
291 CalculateQuatdot();
293 }
294}
295
296//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
297// Compute the quaternion orientation derivative
298//
299// vQtrndot is the quaternion derivative.
300// Reference: See Stevens and Lewis, "Aircraft Control and Simulation",
301// Second edition (2004), eqn 1.5-16b (page 50)
302
303void FGPropagate::CalculateQuatdot(void)
304{
305 // Compute quaternion orientation derivative on current body rates
306 VState.vQtrndot = VState.qAttitudeECI.GetQDot(VState.vPQRi);
307}
308
309//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
310 // Transform the velocity vector of the body relative to the origin (Earth
311 // center) to be expressed in the inertial frame, and add the vehicle velocity
312 // contribution due to the rotation of the planet.
313 // Reference: See Stevens and Lewis, "Aircraft Control and Simulation",
314 // Second edition (2004), eqn 1.5-16c (page 50)
315
316void FGPropagate::CalculateInertialVelocity(void)
317{
318 VState.vInertialVelocity = Tb2i * VState.vUVW + (in.vOmegaPlanet * VState.vInertialPosition);
319}
320
321//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
322 // Transform the velocity vector of the inertial frame to be expressed in the
323 // body frame relative to the origin (Earth center), and substract the vehicle
324 // velocity contribution due to the rotation of the planet.
325
326void FGPropagate::CalculateUVW(void)
327{
328 VState.vUVW = Ti2b * (VState.vInertialVelocity - (in.vOmegaPlanet * VState.vInertialPosition));
329}
330
331//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
332
333void FGPropagate::Integrate( FGColumnVector3& Integrand,
334 FGColumnVector3& Val,
335 deque <FGColumnVector3>& ValDot,
336 double dt,
337 eIntegrateType integration_type)
338{
339 ValDot.push_front(Val);
340 ValDot.pop_back();
341
342 switch(integration_type) {
343 case eRectEuler: Integrand += dt*ValDot[0];
344 break;
345 case eTrapezoidal: Integrand += 0.5*dt*(ValDot[0] + ValDot[1]);
346 break;
347 case eAdamsBashforth2: Integrand += dt*(1.5*ValDot[0] - 0.5*ValDot[1]);
348 break;
349 case eAdamsBashforth3: Integrand += (1/12.0)*dt*(23.0*ValDot[0] - 16.0*ValDot[1] + 5.0*ValDot[2]);
350 break;
351 case eAdamsBashforth4: Integrand += (1/24.0)*dt*(55.0*ValDot[0] - 59.0*ValDot[1] + 37.0*ValDot[2] - 9.0*ValDot[3]);
352 break;
353 case eAdamsBashforth5: Integrand += dt*((1901./720.)*ValDot[0] - (1387./360.)*ValDot[1] + (109./30.)*ValDot[2] - (637./360.)*ValDot[3] + (251./720.)*ValDot[4]);
354 break;
355 case eNone: // do nothing, freeze translational rate
356 break;
357 case eBuss1:
358 case eBuss2:
360 throw("Can only use Buss (1 & 2) or local linearization integration methods in for rotational position!");
361 default:
362 break;
363 }
364}
365
366//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
367
368void FGPropagate::Integrate( FGQuaternion& Integrand,
369 FGQuaternion& Val,
370 deque <FGQuaternion>& ValDot,
371 double dt,
372 eIntegrateType integration_type)
373{
374 ValDot.push_front(Val);
375 ValDot.pop_back();
376
377 switch(integration_type) {
378 case eRectEuler: Integrand += dt*ValDot[0];
379 break;
380 case eTrapezoidal: Integrand += 0.5*dt*(ValDot[0] + ValDot[1]);
381 break;
382 case eAdamsBashforth2: Integrand += dt*(1.5*ValDot[0] - 0.5*ValDot[1]);
383 break;
384 case eAdamsBashforth3: Integrand += (1/12.0)*dt*(23.0*ValDot[0] - 16.0*ValDot[1] + 5.0*ValDot[2]);
385 break;
386 case eAdamsBashforth4: Integrand += (1/24.0)*dt*(55.0*ValDot[0] - 59.0*ValDot[1] + 37.0*ValDot[2] - 9.0*ValDot[3]);
387 break;
388 case eAdamsBashforth5: Integrand += dt*((1901./720.)*ValDot[0] - (1387./360.)*ValDot[1] + (109./30.)*ValDot[2] - (637./360.)*ValDot[3] + (251./720.)*ValDot[4]);
389 break;
390 case eBuss1:
391 {
392 // This is the first order method as described in Samuel R. Buss paper[6].
393 // The formula from Buss' paper is transposed below to quaternions and is
394 // actually the exact solution of the quaternion differential equation
395 // qdot = 1/2*w*q when w is constant.
396 Integrand = Integrand * QExp(0.5 * dt * VState.vPQRi);
397 }
398 return; // No need to normalize since the quaternion exponential is always normal
399 case eBuss2:
400 {
401 // This is the 'augmented second-order method' from S.R. Buss paper [6].
402 // Unlike Runge-Kutta or Adams-Bashforth, it is a one-pass second-order
403 // method (see reference [6]).
404 FGColumnVector3 wi = VState.vPQRi;
405 FGColumnVector3 wdoti = in.vPQRidot;
406 FGColumnVector3 omega = wi + 0.5*dt*wdoti + dt*dt/12.*wdoti*wi;
407 Integrand = Integrand * QExp(0.5 * dt * omega);
408 }
409 return; // No need to normalize since the quaternion exponential is always normal
411 {
412 // This is the local linearization algorithm of Barker et al. (see ref. [7])
413 // It is also a one-pass second-order method. The code below is based on the
414 // more compact formulation issued from equation (107) of ref. [8]. The
415 // constants C1, C2, C3 and C4 have the same value than those in ref. [7] pp. 11
416 FGColumnVector3 wi = 0.5 * VState.vPQRi;
417 FGColumnVector3 wdoti = 0.5 * in.vPQRidot;
418 double omegak2 = DotProduct(VState.vPQRi, VState.vPQRi);
419 double omegak = omegak2 > 1E-6 ? sqrt(omegak2) : 1E-6;
420 double rhok = 0.5 * dt * omegak;
421 double C1 = cos(rhok);
422 double C2 = 2.0 * sin(rhok) / omegak;
423 double C3 = 4.0 * (1.0 - C1) / (omegak*omegak);
424 double C4 = 4.0 * (dt - C2) / (omegak*omegak);
425 FGColumnVector3 Omega = C2*wi + C3*wdoti + C4*wi*wdoti;
426 FGQuaternion q;
427
428 q(1) = C1 - C4*DotProduct(wi, wdoti);
429 q(2) = Omega(eP);
430 q(3) = Omega(eQ);
431 q(4) = Omega(eR);
432
433 Integrand = Integrand * q;
434
435 /* Cross check with ref. [7] pp.11-12 formulas and code pp. 20
436 double pk = VState.vPQRi(eP);
437 double qk = VState.vPQRi(eQ);
438 double rk = VState.vPQRi(eR);
439 double pdotk = in.vPQRidot(eP);
440 double qdotk = in.vPQRidot(eQ);
441 double rdotk = in.vPQRidot(eR);
442 double Ap = -0.25 * (pk*pdotk + qk*qdotk + rk*rdotk);
443 double Bp = 0.25 * (pk*qdotk - qk*pdotk);
444 double Cp = 0.25 * (pdotk*rk - pk*rdotk);
445 double Dp = 0.25 * (qk*rdotk - qdotk*rk);
446 double C2p = sin(rhok) / omegak;
447 double C3p = 2.0 * (1.0 - cos(rhok)) / (omegak*omegak);
448 double H = C1 + C4 * Ap;
449 double G = -C2p*rk - C3p*rdotk + C4*Bp;
450 double J = C2p*qk + C3p*qdotk - C4*Cp;
451 double K = C2p*pk + C3p*pdotk - C4*Dp;
452
453 cout << "q: " << q << endl;
454
455 // Warning! In the paper of Barker et al. the quaternion components are not
456 // ordered the same way as in JSBSim (see equations (2) and (3) of ref. [7]
457 // as well as the comment just below equation (3))
458 cout << "FORTRAN: " << H << " , " << K << " , " << J << " , " << -G << endl;*/
459 }
460 break; // The quaternion q is not normal so the normalization needs to be done.
461 case eNone: // do nothing, freeze rotational rate
462 break;
463 default:
464 break;
465 }
466
467 Integrand.Normalize();
468}
469
470//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
471
472void FGPropagate::UpdateLocationMatrices(void)
473{
474 Tl2ec = VState.vLocation.GetTl2ec(); // local to ECEF transform
475 Tec2l = Tl2ec.Transposed(); // ECEF to local frame transform
476 Ti2l = Tec2l * Ti2ec; // ECI to local frame transform
477 Tl2i = Ti2l.Transposed(); // local to ECI transform
478}
479
480//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
481
482void FGPropagate::UpdateBodyMatrices(void)
483{
484 Ti2b = VState.qAttitudeECI.GetT(); // ECI to body frame transform
485 Tb2i = Ti2b.Transposed(); // body to ECI frame transform
486 Tl2b = Ti2b * Tl2i; // local to body frame transform
487 Tb2l = Tl2b.Transposed(); // body to local frame transform
488 Tec2b = Ti2b * Tec2i; // ECEF to body frame transform
489 Tb2ec = Tec2b.Transposed(); // body to ECEF frame tranform
490
491 Qec2b = Tec2b.GetQuaternion();
492}
493
494//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
495
497{
498 VState.qAttitudeECI = Qi;
499 VState.qAttitudeECI.Normalize();
500 UpdateBodyMatrices();
501 VState.qAttitudeLocal = Tl2b.GetQuaternion();
502 CalculateQuatdot();
503}
504
505//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
506
508 VState.vInertialVelocity = Vi;
509 CalculateUVW();
510 vVel = Tb2l * VState.vUVW;
511}
512
513//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
514
516 VState.vPQRi = Ti2b * vRates;
517 VState.vPQR = VState.vPQRi - Ti2b * in.vOmegaPlanet;
518 CalculateQuatdot();
519}
520
521//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
522
524{
525 return VState.vLocation.GetRadius() - VState.vLocation.GetSeaLevelRadius();
526}
527
528//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
529
531{
532 double slr = VState.vLocation.GetSeaLevelRadius();
533 VState.vLocation.SetRadius(slr + altASL);
534 UpdateVehicleState();
535}
536
537//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
538
540{
541 FGLocation contact;
542 FGColumnVector3 normal;
543 Inertial->GetContactPoint(VState.vLocation, contact, normal,
544 LocalTerrainVelocity, LocalTerrainAngularVelocity);
545}
546
547//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
548
550{
551 FGColumnVector3 vDummy;
552 FGLocation contact;
553 contact.SetEllipse(in.SemiMajor, in.SemiMinor);
554 Inertial->GetContactPoint(VState.vLocation, contact, vDummy, vDummy, vDummy);
555 return contact.GetGeodAltitude();
556}
557
558//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
559
560void FGPropagate::SetTerrainElevation(double terrainElev)
561{
562 Inertial->SetTerrainElevation(terrainElev);
563}
564
565//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
566
568{
569 FGLocation contact;
570 FGColumnVector3 vDummy;
571 Inertial->GetContactPoint(VState.vLocation, contact, vDummy, vDummy, vDummy);
572 return contact.GetRadius();
573}
574
575//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
576
578{
579 return Inertial->GetAltitudeAGL(VState.vLocation);
580}
581
582//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
583
585{
586 return GetDistanceAGL()*0.0003048;
587}
588
589//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
590
592{
593 Inertial->SetAltitudeAGL(VState.vLocation, tt);
594 UpdateVehicleState();
595}
596
597//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
598
600{
601 SetDistanceAGL(tt*3280.8399);
602}
603
604//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
605
607{
608 //ToDo: Shouldn't all of these be set from the vstate vector passed in?
609 VState.vLocation = vstate.vLocation;
610 UpdateLocationMatrices();
613 VState.vUVW = vstate.vUVW;
614 vVel = Tb2l * VState.vUVW;
615 VState.vPQR = vstate.vPQR;
616 VState.vPQRi = VState.vPQR + Ti2b * in.vOmegaPlanet;
617 VState.vInertialPosition = vstate.vInertialPosition;
618 CalculateQuatdot();
619}
620
621//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
622
623void FGPropagate::UpdateVehicleState(void)
624{
626 VState.vInertialPosition = Tec2i * VState.vLocation;
627 UpdateLocationMatrices();
628 UpdateBodyMatrices();
629 vVel = Tb2l * VState.vUVW;
630 VState.qAttitudeLocal = Tl2b.GetQuaternion();
631}
632
633//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
634
636{
637 VState.vLocation = l;
638 UpdateVehicleState();
639}
640
641//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
642
644{
645 return VState.qAttitudeLocal.GetEuler() * radtodeg;
646}
647
648//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
649
651{
652 cout << endl;
653 cout << fgblue
654 << "------------------------------------------------------------------" << reset << endl;
655 cout << highint
656 << "State Report at sim time: " << FDMExec->GetSimTime() << " seconds" << reset << endl;
657 cout << " " << underon
658 << "Position" << underoff << endl;
659 cout << " ECI: " << VState.vInertialPosition.Dump(", ") << " (x,y,z, in ft)" << endl;
660 cout << " ECEF: " << VState.vLocation << " (x,y,z, in ft)" << endl;
661 cout << " Local: " << VState.vLocation.GetGeodLatitudeDeg()
662 << ", " << VState.vLocation.GetLongitudeDeg()
663 << ", " << GetAltitudeASL() << " (geodetic lat, lon, alt ASL in deg and ft)" << endl;
664
665 cout << endl << " " << underon
666 << "Orientation" << underoff << endl;
667 cout << " ECI: " << VState.qAttitudeECI.GetEulerDeg().Dump(", ") << " (phi, theta, psi in deg)" << endl;
668 cout << " Local: " << VState.qAttitudeLocal.GetEulerDeg().Dump(", ") << " (phi, theta, psi in deg)" << endl;
669
670 cout << endl << " " << underon
671 << "Velocity" << underoff << endl;
672 cout << " ECI: " << VState.vInertialVelocity.Dump(", ") << " (x,y,z in ft/s)" << endl;
673 cout << " ECEF: " << (Tb2ec * VState.vUVW).Dump(", ") << " (x,y,z in ft/s)" << endl;
674 cout << " Local: " << GetVel() << " (n,e,d in ft/sec)" << endl;
675 cout << " Body: " << GetUVW() << " (u,v,w in ft/sec)" << endl;
676
677 cout << endl << " " << underon
678 << "Body Rates (relative to given frame, expressed in body frame)" << underoff << endl;
679 cout << " ECI: " << (VState.vPQRi*radtodeg).Dump(", ") << " (p,q,r in deg/s)" << endl;
680 cout << " ECEF: " << (VState.vPQR*radtodeg).Dump(", ") << " (p,q,r in deg/s)" << endl;
681}
682
683//******************************************************************************
684
685void FGPropagate::WriteStateFile(int num)
686{
687 sg_ofstream outfile;
688
689 if (num == 0) return;
690
691 SGPath path = FDMExec->GetOutputPath();
692
693 if (path.isNull()) path = SGPath("initfile.");
694 else path.append("initfile.");
695
696 // Append sim time to the filename since there may be more than one created during a simulation run
697 path.concat(to_string((double)FDMExec->GetSimTime())+".xml");
698
699 switch(num) {
700 case 1:
701 outfile.open(path);
702 if (outfile.is_open()) {
703 outfile << "<?xml version=\"1.0\"?>" << endl;
704 outfile << "<initialize name=\"reset00\">" << endl;
705 outfile << " <ubody unit=\"FT/SEC\"> " << VState.vUVW(eU) << " </ubody> " << endl;
706 outfile << " <vbody unit=\"FT/SEC\"> " << VState.vUVW(eV) << " </vbody> " << endl;
707 outfile << " <wbody unit=\"FT/SEC\"> " << VState.vUVW(eW) << " </wbody> " << endl;
708 outfile << " <phi unit=\"DEG\"> " << VState.qAttitudeLocal.GetEuler(ePhi)*radtodeg << " </phi>" << endl;
709 outfile << " <theta unit=\"DEG\"> " << VState.qAttitudeLocal.GetEuler(eTht)*radtodeg << " </theta>" << endl;
710 outfile << " <psi unit=\"DEG\"> " << VState.qAttitudeLocal.GetEuler(ePsi)*radtodeg << " </psi>" << endl;
711 outfile << " <longitude unit=\"DEG\"> " << VState.vLocation.GetLongitudeDeg() << " </longitude>" << endl;
712 outfile << " <latitude unit=\"DEG\"> " << VState.vLocation.GetLatitudeDeg() << " </latitude>" << endl;
713 outfile << " <altitude unit=\"FT\"> " << GetDistanceAGL() << " </altitude>" << endl;
714 outfile << "</initialize>" << endl;
715 outfile.close();
716 } else {
717 cerr << "Could not open and/or write the state to the initial conditions file: "
718 << path << endl;
719 }
720 break;
721 case 2:
722 outfile.open(path);
723 if (outfile.is_open()) {
724 outfile << "<?xml version=\"1.0\"?>" << endl;
725 outfile << "<initialize name=\"IC File\" version=\"2.0\">" << endl;
726 outfile << "" << endl;
727 outfile << " <position frame=\"ECEF\">" << endl;
728 outfile << " <latitude unit=\"DEG\" type=\"geodetic\"> " << VState.vLocation.GetGeodLatitudeDeg() << " </latitude>" << endl;
729 outfile << " <longitude unit=\"DEG\"> " << VState.vLocation.GetLongitudeDeg() << " </longitude>" << endl;
730 outfile << " <altitudeMSL unit=\"FT\"> " << GetAltitudeASL() << " </altitudeMSL>" << endl;
731 outfile << " </position>" << endl;
732 outfile << "" << endl;
733 outfile << " <orientation unit=\"DEG\" frame=\"LOCAL\">" << endl;
734 outfile << " <yaw> " << VState.qAttitudeLocal.GetEulerDeg(eYaw) << " </yaw>" << endl;
735 outfile << " <pitch> " << VState.qAttitudeLocal.GetEulerDeg(ePitch) << " </pitch>" << endl;
736 outfile << " <roll> " << VState.qAttitudeLocal.GetEulerDeg(eRoll) << " </roll>" << endl;
737 outfile << " </orientation>" << endl;
738 outfile << "" << endl;
739 outfile << " <velocity unit=\"FT/SEC\" frame=\"LOCAL\">" << endl;
740 outfile << " <x> " << GetVel(eNorth) << " </x>" << endl;
741 outfile << " <y> " << GetVel(eEast) << " </y>" << endl;
742 outfile << " <z> " << GetVel(eDown) << " </z>" << endl;
743 outfile << " </velocity>" << endl;
744 outfile << "" << endl;
745 outfile << " <attitude_rate unit=\"DEG/SEC\" frame=\"BODY\">" << endl;
746 outfile << " <roll> " << (VState.vPQR*radtodeg)(eRoll) << " </roll>" << endl;
747 outfile << " <pitch> " << (VState.vPQR*radtodeg)(ePitch) << " </pitch>" << endl;
748 outfile << " <yaw> " << (VState.vPQR*radtodeg)(eYaw) << " </yaw>" << endl;
749 outfile << " </attitude_rate>" << endl;
750 outfile << "" << endl;
751 outfile << "</initialize>" << endl;
752 outfile.close();
753 } else {
754 cerr << "Could not open and/or write the state to the initial conditions file: "
755 << path << endl;
756 }
757 break;
758 default:
759 cerr << "When writing a state file, the supplied value must be 1 or 2 for the version number of the resulting IC file" << endl;
760 }
761}
762
763//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
764
765void FGPropagate::bind(void)
766{
767 typedef double (FGPropagate::*PMF)(int) const;
768 typedef int (FGPropagate::*iPMF)(void) const;
769
770 PropertyManager->Tie("velocities/h-dot-fps", this, &FGPropagate::Gethdot);
771
772 PropertyManager->Tie("velocities/v-north-fps", this, eNorth, (PMF)&FGPropagate::GetVel);
773 PropertyManager->Tie("velocities/v-east-fps", this, eEast, (PMF)&FGPropagate::GetVel);
774 PropertyManager->Tie("velocities/v-down-fps", this, eDown, (PMF)&FGPropagate::GetVel);
775
776 PropertyManager->Tie("velocities/u-fps", this, eU, (PMF)&FGPropagate::GetUVW);
777 PropertyManager->Tie("velocities/v-fps", this, eV, (PMF)&FGPropagate::GetUVW);
778 PropertyManager->Tie("velocities/w-fps", this, eW, (PMF)&FGPropagate::GetUVW);
779
780 PropertyManager->Tie("velocities/p-rad_sec", this, eP, (PMF)&FGPropagate::GetPQR);
781 PropertyManager->Tie("velocities/q-rad_sec", this, eQ, (PMF)&FGPropagate::GetPQR);
782 PropertyManager->Tie("velocities/r-rad_sec", this, eR, (PMF)&FGPropagate::GetPQR);
783
784 PropertyManager->Tie("velocities/pi-rad_sec", this, eP, (PMF)&FGPropagate::GetPQRi);
785 PropertyManager->Tie("velocities/qi-rad_sec", this, eQ, (PMF)&FGPropagate::GetPQRi);
786 PropertyManager->Tie("velocities/ri-rad_sec", this, eR, (PMF)&FGPropagate::GetPQRi);
787
788 PropertyManager->Tie("velocities/eci-x-fps", this, eX, (PMF)&FGPropagate::GetInertialVelocity);
789 PropertyManager->Tie("velocities/eci-y-fps", this, eY, (PMF)&FGPropagate::GetInertialVelocity);
790 PropertyManager->Tie("velocities/eci-z-fps", this, eZ, (PMF)&FGPropagate::GetInertialVelocity);
791
792 PropertyManager->Tie("velocities/eci-velocity-mag-fps", this, &FGPropagate::GetInertialVelocityMagnitude);
793 PropertyManager->Tie("velocities/ned-velocity-mag-fps", this, &FGPropagate::GetNEDVelocityMagnitude);
794
797 PropertyManager->Tie("position/lat-gc-rad", this, &FGPropagate::GetLatitude, &FGPropagate::SetLatitude);
798 PropertyManager->Tie("position/long-gc-rad", this, &FGPropagate::GetLongitude, &FGPropagate::SetLongitude);
799 PropertyManager->Tie("position/lat-gc-deg", this, &FGPropagate::GetLatitudeDeg, &FGPropagate::SetLatitudeDeg);
801 PropertyManager->Tie("position/lat-geod-rad", this, &FGPropagate::GetGeodLatitudeRad);
802 PropertyManager->Tie("position/lat-geod-deg", this, &FGPropagate::GetGeodLatitudeDeg);
803 PropertyManager->Tie("position/geod-alt-ft", this, &FGPropagate::GetGeodeticAltitude);
805 PropertyManager->Tie("position/geod-alt-km", this, &FGPropagate::GetGeodeticAltitudeKm);
807 PropertyManager->Tie("position/radius-to-vehicle-ft", this, &FGPropagate::GetRadius);
808 PropertyManager->Tie("position/terrain-elevation-asl-ft", this,
811
812 PropertyManager->Tie("position/eci-x-ft", this, eX, (PMF)&FGPropagate::GetInertialPosition);
813 PropertyManager->Tie("position/eci-y-ft", this, eY, (PMF)&FGPropagate::GetInertialPosition);
814 PropertyManager->Tie("position/eci-z-ft", this, eZ, (PMF)&FGPropagate::GetInertialPosition);
815
816 PropertyManager->Tie("position/ecef-x-ft", this, eX, (PMF)&FGPropagate::GetLocation);
817 PropertyManager->Tie("position/ecef-y-ft", this, eY, (PMF)&FGPropagate::GetLocation);
818 PropertyManager->Tie("position/ecef-z-ft", this, eZ, (PMF)&FGPropagate::GetLocation);
819
820 PropertyManager->Tie("position/epa-rad", this, &FGPropagate::GetEarthPositionAngle);
821 PropertyManager->Tie("metrics/terrain-radius", this, &FGPropagate::GetLocalTerrainRadius);
822
823 PropertyManager->Tie("attitude/phi-rad", this, (int)ePhi, (PMF)&FGPropagate::GetEuler);
824 PropertyManager->Tie("attitude/theta-rad", this, (int)eTht, (PMF)&FGPropagate::GetEuler);
825 PropertyManager->Tie("attitude/psi-rad", this, (int)ePsi, (PMF)&FGPropagate::GetEuler);
826
827 PropertyManager->Tie("attitude/phi-deg", this, (int)ePhi, (PMF)&FGPropagate::GetEulerDeg);
828 PropertyManager->Tie("attitude/theta-deg", this, (int)eTht, (PMF)&FGPropagate::GetEulerDeg);
829 PropertyManager->Tie("attitude/psi-deg", this, (int)ePsi, (PMF)&FGPropagate::GetEulerDeg);
830
831 PropertyManager->Tie("attitude/roll-rad", this, (int)ePhi, (PMF)&FGPropagate::GetEuler);
832 PropertyManager->Tie("attitude/pitch-rad", this, (int)eTht, (PMF)&FGPropagate::GetEuler);
833 PropertyManager->Tie("attitude/heading-true-rad", this, (int)ePsi, (PMF)&FGPropagate::GetEuler);
834
835 PropertyManager->Tie("simulation/integrator/rate/rotational", (int*)&integrator_rotational_rate);
836 PropertyManager->Tie("simulation/integrator/rate/translational", (int*)&integrator_translational_rate);
837 PropertyManager->Tie("simulation/integrator/position/rotational", (int*)&integrator_rotational_position);
838 PropertyManager->Tie("simulation/integrator/position/translational", (int*)&integrator_translational_position);
839
840 PropertyManager->Tie("simulation/write-state-file", this, (iPMF)0, &FGPropagate::WriteStateFile);
841}
842
843//%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
844// The bitmasked value choices are as follows:
845// unset: In this case (the default) JSBSim would only print
846// out the normally expected messages, essentially echoing
847// the config files as they are read. If the environment
848// variable is not set, debug_lvl is set to 1 internally
849// 0: This requests JSBSim not to output any messages
850// whatsoever.
851// 1: This value explicity requests the normal JSBSim
852// startup messages
853// 2: This value asks for a message to be printed out when
854// a class is instantiated
855// 4: When this value is set, a message is displayed when a
856// FGModel object executes its Run() method
857// 8: When this value is set, various runtime state variables
858// are printed out periodically
859// 16: When set various parameters are sanity checked and
860// a message is printed out when they go out of bounds
861
862void FGPropagate::Debug(int from)
863{
864 if (debug_lvl <= 0) return;
865
866 if (debug_lvl & 1) { // Standard console startup message output
867 if (from == 0) { // Constructor
868
869 }
870 }
871 if (debug_lvl & 2 ) { // Instantiation/Destruction notification
872 if (from == 0) cout << "Instantiated: FGPropagate" << endl;
873 if (from == 1) cout << "Destroyed: FGPropagate" << endl;
874 }
875 if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
876 }
877 if (debug_lvl & 8 && from == 2) { // Runtime state variables
878 cout << endl << fgblue << highint << left
879 << " Propagation Report (English units: ft, degrees) at simulation time " << FDMExec->GetSimTime() << " seconds"
880 << reset << endl;
881 cout << endl;
882 cout << highint << " Earth Position Angle (deg): " << setw(8) << setprecision(3) << reset
883 << GetEarthPositionAngleDeg() << endl;
884 cout << endl;
885 cout << highint << " Body velocity (ft/sec): " << setw(8) << setprecision(3) << reset << VState.vUVW << endl;
886 cout << highint << " Local velocity (ft/sec): " << setw(8) << setprecision(3) << reset << vVel << endl;
887 cout << highint << " Inertial velocity (ft/sec): " << setw(8) << setprecision(3) << reset << VState.vInertialVelocity << endl;
888 cout << highint << " Inertial Position (ft): " << setw(10) << setprecision(3) << reset << VState.vInertialPosition << endl;
889 cout << highint << " Latitude (deg): " << setw(8) << setprecision(3) << reset << VState.vLocation.GetLatitudeDeg() << endl;
890 cout << highint << " Longitude (deg): " << setw(8) << setprecision(3) << reset << VState.vLocation.GetLongitudeDeg() << endl;
891 cout << highint << " Altitude ASL (ft): " << setw(8) << setprecision(3) << reset << GetAltitudeASL() << endl;
892// cout << highint << " Acceleration (NED, ft/sec^2): " << setw(8) << setprecision(3) << reset << Tb2l*GetUVWdot() << endl;
893 cout << endl;
894 cout << highint << " Matrix ECEF to Body (Orientation of Body with respect to ECEF): "
895 << reset << endl << Tec2b.Dump("\t", " ") << endl;
896 cout << highint << " Associated Euler angles (deg): " << setw(8)
897 << setprecision(3) << reset << (Tec2b.GetQuaternion().GetEuler()*radtodeg)
898 << endl << endl;
899
900 cout << highint << " Matrix Body to ECEF (Orientation of ECEF with respect to Body):"
901 << reset << endl << Tb2ec.Dump("\t", " ") << endl;
902 cout << highint << " Associated Euler angles (deg): " << setw(8)
903 << setprecision(3) << reset << (Tb2ec.GetQuaternion().GetEuler()*radtodeg)
904 << endl << endl;
905
906 cout << highint << " Matrix Local to Body (Orientation of Body with respect to Local):"
907 << reset << endl << Tl2b.Dump("\t", " ") << endl;
908 cout << highint << " Associated Euler angles (deg): " << setw(8)
909 << setprecision(3) << reset << (Tl2b.GetQuaternion().GetEuler()*radtodeg)
910 << endl << endl;
911
912 cout << highint << " Matrix Body to Local (Orientation of Local with respect to Body):"
913 << reset << endl << Tb2l.Dump("\t", " ") << endl;
914 cout << highint << " Associated Euler angles (deg): " << setw(8)
915 << setprecision(3) << reset << (Tb2l.GetQuaternion().GetEuler()*radtodeg)
916 << endl << endl;
917
918 cout << highint << " Matrix Local to ECEF (Orientation of ECEF with respect to Local):"
919 << reset << endl << Tl2ec.Dump("\t", " ") << endl;
920 cout << highint << " Associated Euler angles (deg): " << setw(8)
921 << setprecision(3) << reset << (Tl2ec.GetQuaternion().GetEuler()*radtodeg)
922 << endl << endl;
923
924 cout << highint << " Matrix ECEF to Local (Orientation of Local with respect to ECEF):"
925 << reset << endl << Tec2l.Dump("\t", " ") << endl;
926 cout << highint << " Associated Euler angles (deg): " << setw(8)
927 << setprecision(3) << reset << (Tec2l.GetQuaternion().GetEuler()*radtodeg)
928 << endl << endl;
929
930 cout << highint << " Matrix ECEF to Inertial (Orientation of Inertial with respect to ECEF):"
931 << reset << endl << Tec2i.Dump("\t", " ") << endl;
932 cout << highint << " Associated Euler angles (deg): " << setw(8)
933 << setprecision(3) << reset << (Tec2i.GetQuaternion().GetEuler()*radtodeg)
934 << endl << endl;
935
936 cout << highint << " Matrix Inertial to ECEF (Orientation of ECEF with respect to Inertial):"
937 << reset << endl << Ti2ec.Dump("\t", " ") << endl;
938 cout << highint << " Associated Euler angles (deg): " << setw(8)
939 << setprecision(3) << reset << (Ti2ec.GetQuaternion().GetEuler()*radtodeg)
940 << endl << endl;
941
942 cout << highint << " Matrix Inertial to Body (Orientation of Body with respect to Inertial):"
943 << reset << endl << Ti2b.Dump("\t", " ") << endl;
944 cout << highint << " Associated Euler angles (deg): " << setw(8)
945 << setprecision(3) << reset << (Ti2b.GetQuaternion().GetEuler()*radtodeg)
946 << endl << endl;
947
948 cout << highint << " Matrix Body to Inertial (Orientation of Inertial with respect to Body):"
949 << reset << endl << Tb2i.Dump("\t", " ") << endl;
950 cout << highint << " Associated Euler angles (deg): " << setw(8)
951 << setprecision(3) << reset << (Tb2i.GetQuaternion().GetEuler()*radtodeg)
952 << endl << endl;
953
954 cout << highint << " Matrix Inertial to Local (Orientation of Local with respect to Inertial):"
955 << reset << endl << Ti2l.Dump("\t", " ") << endl;
956 cout << highint << " Associated Euler angles (deg): " << setw(8)
957 << setprecision(3) << reset << (Ti2l.GetQuaternion().GetEuler()*radtodeg)
958 << endl << endl;
959
960 cout << highint << " Matrix Local to Inertial (Orientation of Inertial with respect to Local):"
961 << reset << endl << Tl2i.Dump("\t", " ") << endl;
962 cout << highint << " Associated Euler angles (deg): " << setw(8)
963 << setprecision(3) << reset << (Tl2i.GetQuaternion().GetEuler()*radtodeg)
964 << endl << endl;
965
966 cout << setprecision(6); // reset the output stream
967 }
968 if (debug_lvl & 16) { // Sanity checking
969 if (from == 2) { // State sanity checking
970 if (fabs(VState.vPQR.Magnitude()) > 1000.0) {
971 stringstream s;
972 s << "Vehicle rotation rate is excessive (>1000 rad/sec): " << VState.vPQR.Magnitude();
973 cerr << endl << s.str() << endl;
974 throw BaseException(s.str());
975 }
976 if (fabs(VState.vUVW.Magnitude()) > 1.0e10) {
977 stringstream s;
978 s << "Vehicle velocity is excessive (>1e10 ft/sec): " << VState.vUVW.Magnitude();
979 cerr << endl << s.str() << endl;
980 throw BaseException(s.str());
981 }
982 if (fabs(GetDistanceAGL()) > 1e10) {
983 stringstream s;
984 s << "Vehicle altitude is excessive (>1e10 ft): " << GetDistanceAGL();
985 cerr << endl << s.str() << endl;
986 throw BaseException(s.str());
987 }
988 }
989 }
990 if (debug_lvl & 64) {
991 if (from == 0) { // Constructor
992 }
993 }
994}
995}
JSBSim::FGFDMExec * FDMExec
Definition JSBSim.cpp:88
This class implements a 3 element column vector.
FGColumnVector3 & Normalize(void)
Normalize.
const FGQuaternion & GetOrientation(void) const
Gets the initial orientation.
const FGColumnVector3 GetPQRRadpsIC(void) const
Gets the initial body rotation rate.
const FGLocation & GetPosition(void) const
Gets the initial position.
double GetEarthPositionAngleIC(void) const
Gets the initial Earth position angle.
const FGColumnVector3 GetUVWFpsIC(void) const
Gets the initial body velocity.
static constexpr double radtodeg
Definition FGJSBBase.h:348
static char highint[5]
highlights text
Definition FGJSBBase.h:123
static char reset[5]
resets text properties
Definition FGJSBBase.h:129
static char fgblue[6]
blue text
Definition FGJSBBase.h:135
static char underoff[6]
underline off
Definition FGJSBBase.h:133
static short debug_lvl
Definition FGJSBBase.h:190
static char underon[5]
underlines text
Definition FGJSBBase.h:131
FGLocation holds an arbitrary location in the Earth centered Earth fixed reference frame (ECEF).
Definition FGLocation.h:152
double GetGeodAltitude(void) const
Gets the geodetic altitude in feet.
Definition FGLocation.h:279
double GetLongitudeDeg() const
Get the longitude.
Definition FGLocation.h:240
double GetRadius() const
Get the distance from the center of the earth in feet.
Definition FGLocation.h:291
void SetEllipse(double semimajor, double semiminor)
Sets the semimajor and semiminor axis lengths for this planet.
double GetLatitudeDeg() const
Get the GEOCENTRIC latitude in degrees.
Definition FGLocation.h:267
FGQuaternion GetQuaternion(void) const
Returns the quaternion associated with this direction cosine (rotation) matrix.
unsigned int rate
Definition FGModel.h:102
FGPropertyManager * PropertyManager
Definition FGModel.h:117
bool InitModel(void) override
Definition FGModel.cpp:81
FGFDMExec * FDMExec
Definition FGModel.h:116
FGModel(FGFDMExec *)
Constructor.
Definition FGModel.cpp:57
std::string Name
Definition FGModel.h:103
virtual bool Run(bool Holding)
Runs the model; called by the Executive.
Definition FGModel.cpp:89
void SetDistanceAGLKm(double tt)
double GetLatitude(void) const
void SetVState(const VehicleState &vstate)
double GetTerrainElevation(void) const
~FGPropagate()
Destructor.
void RecomputeLocalTerrainVelocity()
void SetInertialRates(const FGColumnVector3 &vRates)
double GetEarthPositionAngle(void) const
Returns the Earth position angle.
double GetAltitudeASLmeters(void) const
Returns the current altitude above sea level.
void SetDistanceAGL(double tt)
double GetGeodLatitudeRad(void) const
void SetAltitudeASLmeters(double altASL)
void SetTerrainElevation(double tt)
double GetDistanceAGLKm(void) const
void SetLatitude(double lat)
double GetNEDVelocityMagnitude(void) const
Retrieves the total local NED velocity in ft/sec.
double GetLongitude(void) const
const FGColumnVector3 & GetEuler(void) const
Retrieves the Euler angles that define the vehicle orientation.
void SetLongitudeDeg(double lon)
void SetLatitudeDeg(double lat)
void SetLongitude(double lon)
void SetInitialState(const FGInitialCondition *)
FGColumnVector3 GetEulerDeg(void) const
Retrieves the Euler angles (in degrees) that define the vehicle orientation.
double GetAltitudeASL(void) const
Returns the current altitude above sea level.
double GetLocalTerrainRadius(void) const
Returns the "constant" LocalTerrainRadius.
void SetInertialVelocity(const FGColumnVector3 &Vi)
const FGColumnVector3 & GetInertialVelocity(void) const
Retrieves the inertial velocity vector in ft/sec.
double GetRadius(void) const
double GetDistanceAGL(void) const
double Gethdot(void) const
Returns the current altitude rate.
const FGColumnVector3 & GetPQRi(void) const
Retrieves the body angular rates vector, relative to the ECI (inertial) frame.
const FGColumnVector3 & GetInertialPosition(void) const
Retrieves the inertial position vector.
void SetHoldDown(bool hd)
Sets the property forces/hold-down.
void SetInertialOrientation(const FGQuaternion &Qi)
FGPropagate(FGFDMExec *Executive)
Constructor.
const FGLocation & GetLocation(void) const
void SetLocation(const FGLocation &l)
void SetAltitudeASL(double altASL)
double GetInertialVelocityMagnitude(void) const
Retrieves the total inertial velocity in ft/sec.
const FGColumnVector3 & GetPQR(void) const
Retrieves the body angular rates vector, relative to the ECEF frame.
double GetLatitudeDeg(void) const
const FGColumnVector3 & GetUVW(void) const
Retrieves the body frame vehicle velocity vector.
bool Run(bool Holding)
Runs the state propagation model; called by the Executive Can pass in a value indicating if the execu...
const FGColumnVector3 & GetVel(void) const
Retrieves the velocity vector.
struct JSBSim::FGPropagate::Inputs in
double GetGeodLatitudeDeg(void) const
double GetLongitudeDeg(void) const
double GetGeodeticAltitudeKm(void) const
double GetGeodeticAltitude(void) const
bool InitModel(void)
Initializes the FGPropagate class after instantiation and prior to first execution.
double GetEarthPositionAngleDeg(void) const
Returns the Earth position angle in degrees.
Models the Quaternion representation of rotations.
const FGColumnVector3 & GetEuler(void) const
Retrieves the Euler angles.
FGQuaternion GetQDot(const FGColumnVector3 &PQR) const
Quaternion derivative for given angular rates.
void Normalize(void)
Normalize.
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 vOmegaPlanet
The current vehicle state vector structure contains the translational and angular position,...
Definition FGPropagate.h:98
FGColumnVector3 vUVW
The velocity vector of the vehicle with respect to the ECEF frame, expressed in the body system.
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.
FGQuaternion qAttitudeECI
The current orientation of the vehicle, that is, the orientation of the body frame relative to the in...
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...