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