46#include "models/FGInertial.h"
47#include "models/FGAccelerations.h"
48#include "models/FGMassBalance.h"
49#include "models/FGFCS.h"
50#include "input_output/FGLog.h"
53#pragma warning (disable : 4786 4788)
68 max_sub_iterations=100;
70 A_Tolerance = Tolerance / 10;
74 fgic = *fdmex->
GetIC();
84 log <<
"Instantiated: FGTrim\n";
90FGTrim::~FGTrim(
void) {
93 log <<
"Destroyed: FGTrim\n";
102 out <<
"\n Trim Statistics:\n";
103 out <<
" Total Iterations: " << total_its <<
"\n";
105 out <<
" Sub-iterations:\n";
106 for (
unsigned int current_axis=0; current_axis<TrimAxes.size(); current_axis++) {
107 run_sum += TrimAxes[current_axis].GetRunCount();
108 out <<
" " << setw(5) << TrimAxes[current_axis].GetStateName().c_str()
109 <<
": " << setprecision(3) << sub_iterations[current_axis]
110 <<
" average: " << setprecision(5) << sub_iterations[current_axis]/double(total_its)
111 <<
" successful: " << setprecision(3) << successful[current_axis]
112 <<
" stability: " << setprecision(5) << TrimAxes[current_axis].GetAvgStability()
115 out <<
" Run Count: " << run_sum <<
"\n";
123 out <<
" Trim Results:\n";
124 for(
unsigned int current_axis=0; current_axis<TrimAxes.size(); current_axis++)
125 TrimAxes[current_axis].AxisReport();
141 vector <FGTrimAxis>::iterator iAxes = TrimAxes.begin();
142 for (; iAxes != TrimAxes.end(); ++iAxes) {
143 if (iAxes->GetStateType() == state)
147 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,state,control));
148 sub_iterations.resize(TrimAxes.size());
149 successful.resize(TrimAxes.size());
150 solution.resize(TrimAxes.size());
161 vector <FGTrimAxis>::iterator iAxes = TrimAxes.begin();
162 while (iAxes != TrimAxes.end()) {
163 if( iAxes->GetStateType() == state ) {
164 iAxes = TrimAxes.erase(iAxes);
171 sub_iterations.resize(TrimAxes.size());
172 successful.resize(TrimAxes.size());
173 solution.resize(TrimAxes.size());
182 vector <FGTrimAxis>::iterator iAxes = TrimAxes.begin();
183 while (iAxes != TrimAxes.end()) {
184 if( iAxes->GetStateType() == state ) {
185 *iAxes =
FGTrimAxis(fdmex,&fgic,state,new_control);
196 bool trim_failed=
false;
198 unsigned int axis_count = 0;
199 auto FCS = fdmex->
GetFCS();
201 vector<double> throttle0 = FCS->GetThrottleCmd();
202 double elevator0 = FCS->GetDeCmd();
203 double aileron0 = FCS->GetDaCmd();
204 double rudder0 = FCS->GetDrCmd();
205 double PitchTrim0 = FCS->GetPitchTrimCmd();
207 for(
int i=0;i < GroundReactions->GetNumGearUnits();i++)
208 GroundReactions->GetGearUnit(i)->SetReport(
false);
210 fdmex->SetTrimStatus(
true);
217 if (mode == tGround) {
226 TrimAxes[1].SetControlLimits(theta - 5.0 * degtorad, theta + 5.0 * degtorad);
227 TrimAxes[2].SetControlLimits(phi - 30.0 * degtorad, phi + 30.0 * degtorad);
231 for(
unsigned int current_axis=0;current_axis<TrimAxes.size();current_axis++) {
235 xlo=TrimAxes[current_axis].GetControlMin();
236 xhi=TrimAxes[current_axis].GetControlMax();
237 TrimAxes[current_axis].SetControl((xlo+xhi)/2);
238 TrimAxes[current_axis].Run();
240 sub_iterations[current_axis]=0;
241 successful[current_axis]=0;
242 solution[current_axis]=
false;
245 if(mode == tPullup ) {
247 log <<
"Setting pitch rate and nlf...\n";
249 log <<
"pitch rate done ...\n";
250 TrimAxes[0].SetStateTarget(targetNlf);
252 }
else if (mode == tTurn) {
259 for(
unsigned int current_axis=0;current_axis<TrimAxes.size();current_axis++) {
260 setDebug(TrimAxes[current_axis]);
263 if(!solution[current_axis]) {
264 if(checkLimits(TrimAxes[current_axis])) {
265 solution[current_axis]=
true;
266 solve(TrimAxes[current_axis]);
268 }
else if(findInterval(TrimAxes[current_axis])) {
269 solve(TrimAxes[current_axis]);
271 solution[current_axis]=
false;
273 sub_iterations[current_axis]+=Nsub;
275 for(
unsigned int current_axis=0;current_axis<TrimAxes.size();current_axis++) {
277 if(Debug > 0) TrimAxes[current_axis].AxisReport();
278 if(TrimAxes[current_axis].InTolerance()) {
280 successful[current_axis]++;
284 if((axis_count == TrimAxes.size()-1) && (TrimAxes.size() > 1)) {
292 for(
unsigned int current_axis=0;current_axis<TrimAxes.size();current_axis++) {
294 if(!TrimAxes[current_axis].InTolerance()) {
295 if(!checkLimits(TrimAxes[current_axis])) {
298 if( (gamma_fallback) &&
299 (TrimAxes[current_axis].GetStateType() == tUdot) &&
300 (TrimAxes[current_axis].GetControlType() == tThrottle)) {
302 log <<
" Can't trim udot with throttle, trying flight"
303 <<
" path angle. (" << N <<
")\n";
304 if(TrimAxes[current_axis].GetState() > 0)
305 TrimAxes[current_axis].SetControlToMin();
307 TrimAxes[current_axis].SetControlToMax();
308 TrimAxes[current_axis].Run();
309 TrimAxes[current_axis]=
FGTrimAxis(fdmex,&fgic,tUdot,tGamma);
312 log <<
" Sorry, " << TrimAxes[current_axis].GetStateName()
313 <<
" doesn't appear to be trimmable\n";
322 if(N > max_iterations)
324 }
while((axis_count < TrimAxes.size()) && (!trim_failed));
326 if((!trim_failed) && (axis_count >= TrimAxes.size())) {
330 log <<
"\n Trim successful\n";
336 fgic = *fdmex->
GetIC();
337 FCS->SetDeCmd(elevator0);
338 FCS->SetDaCmd(aileron0);
339 FCS->SetDrCmd(rudder0);
340 FCS->SetPitchTrimCmd(PitchTrim0);
341 for (
unsigned int i=0; i < throttle0.size(); i++)
342 FCS->SetThrottleCmd(i, throttle0[i]);
348 if (GroundReactions->GetWOW())
353 log <<
"\n Trim failed\n";
359 fdmex->SetTrimStatus(
false);
361 for(
int i=0;i < GroundReactions->GetNumGearUnits();i++)
362 GroundReactions->GetGearUnit(i)->SetReport(
true);
391void FGTrim::trimOnGround(
void)
397 vector<ContactPoints> contacts;
398 FGLocation CGLocation = Propagate->GetLocation();
406 for (
int i = 0; i < GroundReactions->GetNumGearUnits(); ++i) {
408 auto gear = GroundReactions->GetGearUnit(i);
411 if (!gear->GetGearUnitDown())
414 c.location = gear->GetBodyLocation();
419 double height = fdmex->
GetInertial()->GetContactPoint(gearLoc, lDummy,
423 if (gear->IsBogey() && !GroundReactions->GetSolid())
426 c.normal = Tec2b * normal;
427 contacts.push_back(c);
431 contactRef = contacts.size() - 1;
435 if (contacts.size() < 3)
442 FGColumnVector3 contact0 = contacts[contactRef].location;
443 contacts.erase(contacts.begin() + contactRef);
453 FGColumnVector3 force = MassBalance->GetMass() * Accelerations->GetUVWdot();
454 FGColumnVector3 moment = MassBalance->GetJ() * Accelerations->GetPQRdot()
456 FGColumnVector3 rotationAxis = moment.Normalize();
460 RotationParameters rParam = calcRotation(contacts, rotationAxis, contact0);
461 FGQuaternion q0(rParam.angleMin, rotationAxis);
464 FGMatrix33 rot = q0.GetTInv();
465 vector<ContactPoints>::iterator iter;
466 for (iter = contacts.begin(); iter != contacts.end(); ++iter)
467 iter->location = contact0 + rot * (iter->location - contact0);
472 FGColumnVector3 contact1 = rParam.contactRef->location;
473 contacts.erase(rParam.contactRef);
478 rotationAxis = contact1 - contact0;
481 rotationAxis = contact0 - contact1;
483 rotationAxis.Normalize();
486 rParam = calcRotation(contacts, rotationAxis, contact0);
487 FGQuaternion q1(rParam.angleMin, rotationAxis);
490 FGColumnVector3 euler = (fgic.
GetOrientation() * q0 * q1).GetEuler();
506FGTrim::RotationParameters FGTrim::calcRotation(vector<ContactPoints>& contacts,
507 const FGColumnVector3& u,
508 const FGColumnVector3& GM0)
510 RotationParameters rParam;
511 vector<ContactPoints>::iterator iter;
513 rParam.angleMin = 3.0 * M_PI;
515 for (iter = contacts.begin(); iter != contacts.end(); ++iter) {
519 FGColumnVector3 t = u * iter->normal;
520 double length = t.Magnitude();
522 FGColumnVector3 v = t * u;
523 FGColumnVector3 MM0 = GM0 - iter->location;
528 double sqrRadius =
DotProduct(MM0, MM0) - d0 * d0;
531 double DistPlane = d0 *
DotProduct(u, iter->normal) / length;
534 double mag = sqrRadius - DistPlane * DistPlane;
536 FGLogging log(LogLevel::WARN);
537 log <<
"FGTrim::calcRotation DistPlane^2 larger than sqrRadius\n";
539 double alpha = sqrt(max(mag, 0.0));
540 FGColumnVector3 CP = alpha * t + DistPlane * v;
544 double cosine = -
DotProduct(MM0, CP) / sqrRadius;
545 double sine =
DotProduct(MM0 * u, CP) / sqrRadius;
546 double angle = atan2(sine, cosine);
547 if (angle < 0.0) angle += 2.0 * M_PI;
548 if (angle < rParam.angleMin) {
549 rParam.angleMin = angle;
550 rParam.contactRef = iter;
559bool FGTrim::solve(FGTrimAxis& axis) {
561 double x1,x2,x3,f1,f2,f3,d,d0;
562 const double relax =0.9;
563 double eps=axis.GetSolverEps();
569 if( solutionDomain != 0) {
580 while ( (axis.InTolerance() ==
false )
581 && (fabs(d) > eps) && (Nsub < max_sub_iterations)) {
584 x2=x1-d*d0*f1/(f3-f1);
589 FGLogging log(LogLevel::DEBUG);
590 log <<
"FGTrim::solve Nsub,x1,x2,x3: " << Nsub <<
", " << x1
591 <<
", " << x2 <<
", " << x3 <<
"\n";
592 log <<
" " << f1 <<
", " << f2 <<
", " << f3 <<
"\n";
601 else if(f2*f3 <= 0.0) {
614 if(Nsub < max_sub_iterations) success=
true;
644bool FGTrim::findInterval(FGTrimAxis& axis) {
647 double current_control=axis.GetControl();
648 double current_accel=axis.GetState();;
649 double xmin=axis.GetControlMin();
650 double xmax=axis.GetControlMax();
651 double lastxlo,lastxhi,lastalo,lastahi;
653 step=0.025*fabs(xmax);
654 xlo=xhi=current_control;
655 alo=ahi=current_accel;
656 lastxlo=xlo;lastxhi=xhi;
657 lastalo=alo;lastahi=ahi;
663 if(xlo < xmin) xlo=xmin;
665 if(xhi > xmax) xhi=xmax;
666 axis.SetControl(xlo);
669 axis.SetControl(xhi);
672 if(fabs(ahi-alo) <= axis.GetTolerance())
continue;
675 if(alo*current_accel <= 0) {
689 lastxlo=xlo;lastxhi=xhi;
690 lastalo=alo;lastahi=ahi;
691 if( !found && xlo==xmin && xhi==xmax )
continue;
693 FGLogging log(LogLevel::DEBUG);
694 log <<
"FGTrim::findInterval: Nsub=" << Nsub <<
" Lo= " << xlo
695 <<
" Hi= " << xhi <<
" alo*ahi: " << alo*ahi <<
"\n";
697 }
while(!found && (Nsub <= max_sub_iterations) );
716bool FGTrim::checkLimits(FGTrimAxis& axis)
719 double current_control=axis.GetControl();
720 double current_accel=axis.GetState();
721 xlo=axis.GetControlMin();
722 xhi=axis.GetControlMax();
724 axis.SetControl(xlo);
727 axis.SetControl(xhi);
731 FGLogging log(LogLevel::DEBUG);
732 log <<
"checkLimits() xlo,xhi,alo,ahi: " << xlo <<
", " << xhi <<
", "
733 << alo <<
", " << ahi <<
"\n";
736 solutionExists=
false;
737 if(fabs(ahi-alo) > axis.GetTolerance()) {
738 if(alo*current_accel <= 0) {
743 }
else if(current_accel*ahi < 0){
750 axis.SetControl(current_control);
752 return solutionExists;
757void FGTrim::setupPullup() {
761 FGLogging log(LogLevel::INFO);
762 log <<
"setPitchRateInPullup(): " << g <<
", " << cgamma <<
", "
765 log << targetNlf <<
", " << q <<
"\n";
767 log <<
"setPitchRateInPullup() complete\n";
773void FGTrim::setupTurn(
void){
776 if( fabs(phi) > 0.001 && fabs(phi) < 1.56 ) {
777 targetNlf = 1 / cos(phi);
778 g = fdmex->
GetInertial()->GetGravity().Magnitude();
780 FGLogging log(LogLevel::INFO);
781 log << targetNlf <<
", " << psidot <<
"\n";
788void FGTrim::updateRates(
void){
794 const double g = fdmex->
GetInertial()->GetGravity().Magnitude();
799 if (fabs(phi) > 0.001 && fabs(phi) < 1.56) {
801 psidot = g * tan(phi) / vTrue;
802 p = -psidot * sin(theta);
803 q = psidot * cos(theta) * sin(phi);
804 r = psidot * cos(theta) * cos(phi);
810 }
else if (mode == tPullup && fabs(targetNlf - 1.0) > 0.01) {
811 const double g = fdmex->
GetInertial()->GetGravity().Magnitude();
813 const double q = g * (targetNlf - cgamma) / vTrue;
820void FGTrim::setDebug(FGTrimAxis& axis) {
821 if(debug_axis == tAll ||
822 axis.GetStateType() == debug_axis ) {
840 log <<
" Full Trim\n";
842 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tWdot,tAlpha));
843 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tUdot,tThrottle ));
844 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tQdot,tPitchTrim ));
846 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tVdot,tPhi ));
847 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tPdot,tAileron ));
848 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tRdot,tRudder ));
853 log <<
" Longitudinal Trim\n";
855 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tWdot,tAlpha ));
856 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tUdot,tThrottle ));
857 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tQdot,tPitchTrim ));
862 log <<
" Ground Trim\n";
864 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tWdot,tAltAGL ));
865 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tQdot,tTheta ));
866 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tPdot,tPhi ));
869 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tNlf,tAlpha ));
870 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tUdot,tThrottle ));
871 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tQdot,tPitchTrim ));
872 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tHmgt,tBeta ));
873 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tVdot,tPhi ));
874 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tPdot,tAileron ));
875 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tRdot,tRudder ));
878 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tWdot,tAlpha ));
879 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tUdot,tThrottle ));
880 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tQdot,tPitchTrim ));
881 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tVdot,tBeta ));
882 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tPdot,tAileron ));
883 TrimAxes.push_back(
FGTrimAxis(fdmex,&fgic,tRdot,tRudder ));
891 sub_iterations.resize(TrimAxes.size());
892 successful.resize(TrimAxes.size());
893 solution.resize(TrimAxes.size());
This class implements a 3 element column vector.
Encapsulates the JSBSim simulation executive.
std::shared_ptr< FGInitialCondition > GetIC(void) const
Returns a pointer to the FGInitialCondition object.
std::shared_ptr< FGFCS > GetFCS(void) const
Returns the FGFCS pointer.
std::shared_ptr< FGPropagate > GetPropagate(void) const
Returns the FGPropagate pointer.
bool Run(void)
This function executes each scheduled model in succession.
std::shared_ptr< FGGroundReactions > GetGroundReactions(void) const
Returns the FGGroundReactions pointer.
std::shared_ptr< FGMassBalance > GetMassBalance(void) const
Returns the FGAircraft pointer.
void SuspendIntegration(void)
Suspends the simulation and sets the delta T to zero.
std::shared_ptr< FGInertial > GetInertial(void) const
Returns the FGInertial pointer.
void Initialize(const FGInitialCondition *FGIC)
Initializes the simulation with initial conditions.
std::shared_ptr< FGAccelerations > GetAccelerations(void) const
Returns the FGAccelerations pointer.
void ResumeIntegration(void)
Resumes the simulation by resetting delta T to the correct value.
void SetPsiRadIC(double psi)
Sets the initial heading angle.
double GetPhiRadIC(void) const
Gets the initial roll angle.
double GetAltitudeASLFtIC(void) const
Gets the initial altitude above sea level.
double GetTargetNlfIC(void) const
Gets the target normal load factor set from IC.
const FGQuaternion & GetOrientation(void) const
Gets the initial orientation.
void SetPRadpsIC(double P)
Sets the initial body axis roll rate.
double GetAltitudeAGLFtIC(void) const
Gets the initial altitude above ground level.
void SetThetaRadIC(double theta)
Sets the initial pitch angle.
double GetVtrueFpsIC(void) const
Gets the initial true velocity.
void SetRRadpsIC(double R)
Sets the initial body axis yaw rate.
void SetQRadpsIC(double Q)
Sets the initial body axis pitch rate.
double GetUBodyFpsIC(void) const
Gets the initial body axis X velocity.
double GetFlightPathAngleRadIC(void) const
Gets the initial flight path angle.
double GetThetaRadIC(void) const
Gets the initial pitch angle.
void SetPhiRadIC(double phi)
Sets the initial roll angle.
void SetAltitudeASLFtIC(double altitudeASL)
Sets the altitude above sea level initial condition in feet.
FGLocation holds an arbitrary location in the Earth centered Earth fixed reference frame (ECEF).
FGLocation LocalToLocation(const FGColumnVector3 &lvec) const
Conversion from Local frame coordinates to a location in the earth centered and fixed frame.
Handles matrix math operations.
bool AddState(State state, Control control)
Add a state-control pair to the current configuration.
void ClearStates(void)
Clear all state-control pairs from the current configuration.
bool RemoveState(State state)
Remove a specific state-control pair from the current configuration.
void Report(void)
Print the results of the trim.
void SetMode(TrimMode tm)
Clear all state-control pairs and set a predefined trim mode.
bool DoTrim(void)
Execute the trim.
bool EditState(State state, Control new_control)
Change the control used to zero a state previously configured.
FGTrim(FGFDMExec *FDMExec, TrimMode tm=tGround)
Initializes the trimming class.
void TrimStats()
Iteration statistics.
Main namespace for the JSBSim Flight Dynamics Model.
State
Models an aircraft axis for purposes of trimming.
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...