JSBSim Flight Dynamics Model 1.0 (23 February 2013)
An Open Source Flight Dynamics and Control Software Library in C++

FGTrimAnalysisControl.cpp

00001 /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00002 
00003  Header:       FGTrimAnalysisControl.cpp
00004  Author:       Agostino De Marco
00005  Date started: Dec/14/2006
00006 
00007  ------------- Copyright (C) 2006  Agostino De Marco (agodemar@unina.it) -------
00008 
00009  This program is free software; you can redistribute it and/or modify it under
00010  the terms of the GNU Lesser General Public License as published by the Free Software
00011  Foundation; either version 2 of the License, or (at your option) any later
00012  version.
00013 
00014  This program is distributed in the hope that it will be useful, but WITHOUT
00015  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00016  FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public License for more
00017  details.
00018 
00019  You should have received a copy of the GNU Lesser General Public License along with
00020  this program; if not, write to the Free Software Foundation, Inc., 59 Temple
00021  Place - Suite 330, Boston, MA  02111-1307, USA.
00022 
00023  Further information about the GNU Lesser General Public License can also be found on
00024  the world wide web at http://www.gnu.org.
00025 
00026 
00027  HISTORY
00028 --------------------------------------------------------------------------------
00029 12/14/06   ADM   Created
00030 
00031 
00032 FUNCTIONAL DESCRIPTION
00033 --------------------------------------------------------------------------------
00034 
00035 
00036 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00037 SENTRY
00038 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
00039 
00040 #ifdef _MSC_VER
00041 #  pragma warning (disable : 4786)
00042 #endif
00043 
00044 #include <string>
00045 #include <stdlib.h>
00046 
00047 #include "FGFDMExec.h"
00048 #include "models/FGAtmosphere.h"
00049 #include "initialization/FGInitialCondition.h"
00050 #include "FGTrimAnalysisControl.h"
00051 #include "models/FGAircraft.h"
00052 #include "models/FGPropulsion.h"
00053 #include "models/FGAerodynamics.h"
00054 
00055 namespace JSBSim {
00056 
00057 static const char *IdSrc = "$Id: FGTrimAnalysisControl.cpp,v 1.5 2012/09/05 21:49:19 bcoconni Exp $";
00058 static const char *IdHdr = ID_TRIMANALYSISCONTROL;
00059 
00060 /*****************************************************************************/
00061 
00062 FGTrimAnalysisControl::FGTrimAnalysisControl(FGFDMExec* fdex, FGInitialCondition* ic,
00063                                              //State st,
00064                                              TaControl ctrl) {
00065   fdmex=fdex;
00066   fgic=ic;
00067   state=taAll;
00068   control=ctrl;
00069   control_initial_value = 0.;
00070   control_value=0;
00071   state_convert=1.0;
00072   control_convert=1.0;
00073   state_value=0;
00074   state_target=0;
00075   control_tolerance = DEFAULT_TOLERANCE;
00076 
00077   switch(control) {
00078   case taThrottle:
00079     control_min=0;
00080     control_max=1;
00081     control_step=0.2;
00082     control_initial_value = 0.5;
00083     control_value=control_initial_value;
00084     control_name = "Throttle (cmd,norm)";
00085     break;
00086   case taBeta:
00087     control_min=-30*degtorad;
00088     control_max=30*degtorad;
00089     control_step=1*degtorad;
00090     control_convert=radtodeg;
00091     break;
00092   case taAlpha:
00093     control_min=fdmex->GetAerodynamics()->GetAlphaCLMin();
00094     control_max=fdmex->GetAerodynamics()->GetAlphaCLMax();
00095     if(control_max <= control_min) {
00096       control_max=20*degtorad;
00097       control_min=-5*degtorad;
00098     }
00099     control_step=1*degtorad;
00100     control_initial_value = (control_min+control_max)/2;
00101     control_value= control_initial_value;
00102     control_convert=radtodeg;
00103     break;
00104   case taPitchTrim:
00105     control_name = "Pitch Trim (cmd,norm)";
00106     control_min=-1;
00107     control_max=1;
00108     control_step=0.1;
00109     state_convert=radtodeg;
00110     break;
00111   case taElevator:
00112     control_name = "Elevator (cmd,norm)";
00113     control_min=-1;
00114     control_max=1;
00115     control_step=0.1;
00116     state_convert=radtodeg;
00117     break;
00118   case taRollTrim:
00119     control_name = "Roll Trim (cmd,norm)";
00120     control_min=-1;
00121     control_max=1;
00122     control_step=0.1;
00123     state_convert=radtodeg;
00124     break;
00125   case taAileron:
00126     control_name = "Ailerons (cmd,norm)";
00127     control_min=-1;
00128     control_max=1;
00129     control_step=0.1;
00130     state_convert=radtodeg;
00131     break;
00132   case taYawTrim:
00133     control_name = "Yaw Trim (cmd,norm)";
00134     control_min=-1;
00135     control_max=1;
00136     control_step=0.1;
00137     state_convert=radtodeg;
00138     break;
00139   case taRudder:
00140     control_name = "Rudder (cmd,norm)";
00141     control_min=-1;
00142     control_max=1;
00143     control_step=0.1;
00144     state_convert=radtodeg;
00145     break;
00146   case taAltAGL:
00147     control_name = "Altitude (ft)";
00148     control_min=0;
00149     control_max=30;
00150     control_step=2;
00151     control_initial_value = fdmex->GetPropagate()->GetDistanceAGL();
00152     control_value = control_initial_value;
00153     break;
00154   case taPhi:
00155     control_name = "Phi (rad)";
00156     control_min=fdmex->GetPropagate()->GetEuler(ePhi) - 30*degtorad;
00157     control_max=fdmex->GetPropagate()->GetEuler(ePhi) + 30*degtorad;
00158     control_step=1*degtorad;
00159     state_convert=radtodeg;
00160     control_convert=radtodeg;
00161     break;
00162   case taTheta:
00163     control_name = "Theta (rad)";
00164     control_min=fdmex->GetPropagate()->GetEuler(eTht) - 5*degtorad;
00165     control_max=fdmex->GetPropagate()->GetEuler(eTht) + 5*degtorad;
00166     control_step=1*degtorad;
00167     state_convert=radtodeg;
00168     break;
00169   case taHeading:
00170     control_name = "Heading (rad)";
00171     control_min=fdmex->GetPropagate()->GetEuler(ePsi) - 30*degtorad;
00172     control_max=fdmex->GetPropagate()->GetEuler(ePsi) + 30*degtorad;
00173     control_step=1*degtorad;
00174     state_convert=radtodeg;
00175     break;
00176   case taGamma:
00177     control_name = "Gamma (rad)";
00178     control_min=-80*degtorad;
00179     control_max=80*degtorad;
00180     control_step=1*degtorad;
00181     control_convert=radtodeg;
00182     break;
00183   }
00184 
00185 //  if (debug_lvl > 0)
00186 //    cout << "FGTrimAnalysisControl created: "<< control_name << endl;
00187 
00188   Debug(0);
00189 }
00190 
00191 /*****************************************************************************/
00192 
00193 FGTrimAnalysisControl::~FGTrimAnalysisControl(void)
00194 {
00195   Debug(1);
00196 }
00197 
00198 /*****************************************************************************/
00199 
00200 void FGTrimAnalysisControl::getState(void) {
00201   switch(state) {
00202   case taUdot: state_value=fdmex->GetPropagate()->GetUVWdot(1)-state_target; break;
00203   case taVdot: state_value=fdmex->GetPropagate()->GetUVWdot(2)-state_target; break;
00204   case taWdot: state_value=fdmex->GetPropagate()->GetUVWdot(3)-state_target; break;
00205   case taPdot: state_value=fdmex->GetPropagate()->GetPQRdot(1)-state_target; break;
00206   case taQdot: state_value=fdmex->GetPropagate()->GetPQRdot(2)-state_target;break;
00207   case taRdot: state_value=fdmex->GetPropagate()->GetPQRdot(3)-state_target; break;
00208   case taHmgt: state_value=computeHmgt()-state_target; break;
00209   case taNlf:  state_value=fdmex->GetAircraft()->GetNlf()-state_target; break;
00210   case taAll: break;
00211   }
00212 }
00213 
00214 /*****************************************************************************/
00215 
00216 //States are not settable
00217 
00218 void FGTrimAnalysisControl::getControl(void) {
00219   switch(control) {
00220   case taThrottle:  control_value=fdmex->GetFCS()->GetThrottleCmd(0); break;
00221   case taBeta:      control_value=fdmex->GetAuxiliary()->Getbeta(); break;
00222   case taAlpha:     control_value=fdmex->GetAuxiliary()->Getalpha();  break;
00223   case taPitchTrim: control_value=fdmex->GetFCS()->GetPitchTrimCmd(); break;
00224   case taElevator:  control_value=fdmex->GetFCS()->GetDeCmd(); break;
00225   case taRollTrim:
00226   case taAileron:   control_value=fdmex->GetFCS()->GetDaCmd(); break;
00227   case taYawTrim:
00228   case taRudder:    control_value=fdmex->GetFCS()->GetDrCmd(); break;
00229   case taAltAGL:    control_value=fdmex->GetPropagate()->GetDistanceAGL();break;
00230   case taTheta:     control_value=fdmex->GetPropagate()->GetEuler(eTht); break;
00231   case taPhi:       control_value=fdmex->GetPropagate()->GetEuler(ePhi); break;
00232   case taGamma:     control_value=fdmex->GetAuxiliary()->GetGamma();break;
00233   case taHeading:   control_value=fdmex->GetPropagate()->GetEuler(ePsi); break;
00234   }
00235 }
00236 
00237 /*****************************************************************************/
00238 
00239 double FGTrimAnalysisControl::computeHmgt(void) {
00240   double diff;
00241 
00242   diff   = fdmex->GetPropagate()->GetEuler(ePsi) -
00243              fdmex->GetAuxiliary()->GetGroundTrack();
00244 
00245   if( diff < -M_PI ) {
00246      return (diff + 2*M_PI);
00247   } else if( diff > M_PI ) {
00248      return (diff - 2*M_PI);
00249   } else {
00250      return diff;
00251   }
00252 
00253 }
00254 
00255 /*****************************************************************************/
00256 
00257 
00258 void FGTrimAnalysisControl::setControl(void) {
00259   switch(control) {
00260   case taThrottle:  setThrottlesPct(); break;
00261   case taBeta:      fgic->SetBetaRadIC(control_value); break;
00262   case taAlpha:     fgic->SetAlphaRadIC(control_value);  break;
00263   case taPitchTrim: fdmex->GetFCS()->SetPitchTrimCmd(control_value); break;
00264   case taElevator:  fdmex->GetFCS()->SetDeCmd(control_value); break;
00265   case taRollTrim:
00266   case taAileron:   fdmex->GetFCS()->SetDaCmd(control_value); break;
00267   case taYawTrim:
00268   case taRudder:    fdmex->GetFCS()->SetDrCmd(control_value); break;
00269   case taAltAGL:    fgic->SetAltitudeAGLFtIC(control_value); break;
00270   case taTheta:     fgic->SetThetaRadIC(control_value); break;
00271   case taPhi:       fgic->SetPhiRadIC(control_value); break;
00272   case taGamma:     fgic->SetFlightPathAngleRadIC(control_value); break;
00273   case taHeading:   fgic->SetPsiRadIC(control_value); break;
00274   }
00275 }
00276 
00277 
00278 
00279 
00280 
00281 /*****************************************************************************/
00282 
00283 // the aircraft center of rotation is no longer the cg once the gear
00284 // contact the ground so the altitude needs to be changed when pitch
00285 // and roll angle are adjusted.  Instead of attempting to calculate the
00286 // new center of rotation, pick a gear unit as a reference and use its
00287 // location vector to calculate the new height change. i.e. new altitude =
00288 // earth z component of that vector (which is in body axes )
00289 void FGTrimAnalysisControl::SetThetaOnGround(double ff) {
00290   int center,i,ref;
00291 
00292   // favor an off-center unit so that the same one can be used for both
00293   // pitch and roll.  An on-center unit is used (for pitch)if that's all
00294   // that's in contact with the ground.
00295   i=0; ref=-1; center=-1;
00296   while( (ref < 0) && (i < fdmex->GetGroundReactions()->GetNumGearUnits()) ) {
00297     if(fdmex->GetGroundReactions()->GetGearUnit(i)->GetWOW()) {
00298       if(fabs(fdmex->GetGroundReactions()->GetGearUnit(i)->GetBodyLocation(2)) > 0.01)
00299         ref=i;
00300       else
00301         center=i;
00302     }
00303     i++;
00304   }
00305   if((ref < 0) && (center >= 0)) {
00306     ref=center;
00307   }
00308   cout << "SetThetaOnGround ref gear: " << ref << endl;
00309   if(ref >= 0) {
00310     double sp = fdmex->GetPropagate()->GetSinEuler(ePhi);
00311     double cp = fdmex->GetPropagate()->GetCosEuler(ePhi);
00312     double lx = fdmex->GetGroundReactions()->GetGearUnit(ref)->GetBodyLocation(1);
00313     double ly = fdmex->GetGroundReactions()->GetGearUnit(ref)->GetBodyLocation(2);
00314     double lz = fdmex->GetGroundReactions()->GetGearUnit(ref)->GetBodyLocation(3);
00315     double hagl = -1*lx*sin(ff) +
00316                     ly*sp*cos(ff) +
00317                     lz*cp*cos(ff);
00318 
00319     fgic->SetAltitudeAGLFtIC(hagl);
00320     cout << "SetThetaOnGround new alt: " << hagl << endl;
00321   }
00322   fgic->SetThetaRadIC(ff);
00323   cout << "SetThetaOnGround new theta: " << ff << endl;
00324 }
00325 
00326 /*****************************************************************************/
00327 
00328 bool FGTrimAnalysisControl::initTheta(void) {
00329   int i,N;
00330   int iForward = 0;
00331   int iAft = 1;
00332   double zAft,zForward,zDiff,theta;
00333   double xAft,xForward,xDiff;
00334   bool level;
00335   double saveAlt;
00336 
00337   saveAlt=fgic->GetAltitudeAGLFtIC();
00338   fgic->SetAltitudeAGLFtIC(100);
00339 
00340 
00341   N=fdmex->GetGroundReactions()->GetNumGearUnits();
00342 
00343   //find the first wheel unit forward of the cg
00344   //the list is short so a simple linear search is fine
00345   for( i=0; i<N; i++ ) {
00346     if(fdmex->GetGroundReactions()->GetGearUnit(i)->GetBodyLocation(1) > 0 ) {
00347         iForward=i;
00348         break;
00349     }
00350   }
00351   //now find the first wheel unit aft of the cg
00352   for( i=0; i<N; i++ ) {
00353     if(fdmex->GetGroundReactions()->GetGearUnit(i)->GetBodyLocation(1) < 0 ) {
00354         iAft=i;
00355         break;
00356     }
00357   }
00358 
00359   // now adjust theta till the wheels are the same distance from the ground
00360   xAft=fdmex->GetGroundReactions()->GetGearUnit(iAft)->GetBodyLocation(1);
00361   xForward=fdmex->GetGroundReactions()->GetGearUnit(iForward)->GetBodyLocation(1);
00362   xDiff = xForward - xAft;
00363   zAft=fdmex->GetGroundReactions()->GetGearUnit(iAft)->GetLocalGear(3);
00364   zForward=fdmex->GetGroundReactions()->GetGearUnit(iForward)->GetLocalGear(3);
00365   zDiff = zForward - zAft;
00366   level=false;
00367   theta=fgic->GetThetaDegIC();
00368   while(!level && (i < 100)) {
00369     theta+=radtodeg*atan(zDiff/xDiff);
00370     fgic->SetThetaDegIC(theta);
00371     fdmex->SuspendIntegration();
00372     fdmex->Initialize(fgic);
00373     fdmex->Run();
00374     fdmex->ResumeIntegration();
00375     zAft=fdmex->GetGroundReactions()->GetGearUnit(iAft)->GetLocalGear(3);
00376     zForward=fdmex->GetGroundReactions()->GetGearUnit(iForward)->GetLocalGear(3);
00377     zDiff = zForward - zAft;
00378     //cout << endl << theta << "  " << zDiff << endl;
00379     //cout << "0: " << fdmex->GetGroundReactions()->GetGearUnit(0)->GetLocalGear() << endl;
00380     //cout << "1: " << fdmex->GetGroundReactions()->GetGearUnit(1)->GetLocalGear() << endl;
00381     if(fabs(zDiff ) < 0.1)
00382         level=true;
00383     i++;
00384   }
00385   //cout << i << endl;
00386   if (debug_lvl > 0) {
00387       cout << "    Initial Theta: " << fdmex->GetPropagate()->GetEuler(eTht)*radtodeg << endl;
00388       cout << "    Used gear unit " << iAft << " as aft and " << iForward << " as forward" << endl;
00389   }
00390   control_min=(theta+5)*degtorad;
00391   control_max=(theta-5)*degtorad;
00392   fgic->SetAltitudeAGLFtIC(saveAlt);
00393   if(i < 100)
00394     return true;
00395   else
00396     return false;
00397 }
00398 
00399 /*****************************************************************************/
00400 
00401 void FGTrimAnalysisControl::SetPhiOnGround(double ff) {
00402   int i,ref;
00403 
00404   i=0; ref=-1;
00405   //must have an off-center unit here
00406   while ( (ref < 0) && (i < fdmex->GetGroundReactions()->GetNumGearUnits()) ) {
00407     if ( (fdmex->GetGroundReactions()->GetGearUnit(i)->GetWOW()) &&
00408       (fabs(fdmex->GetGroundReactions()->GetGearUnit(i)->GetBodyLocation(2)) > 0.01))
00409         ref=i;
00410     i++;
00411   }
00412   if (ref >= 0) {
00413     double st = fdmex->GetPropagate()->GetSinEuler(eTht);
00414     double ct = fdmex->GetPropagate()->GetCosEuler(eTht);
00415     double lx = fdmex->GetGroundReactions()->GetGearUnit(ref)->GetBodyLocation(1);
00416     double ly = fdmex->GetGroundReactions()->GetGearUnit(ref)->GetBodyLocation(2);
00417     double lz = fdmex->GetGroundReactions()->GetGearUnit(ref)->GetBodyLocation(3);
00418     double hagl = -1*lx*st +
00419                     ly*sin(ff)*ct +
00420                     lz*cos(ff)*ct;
00421 
00422     fgic->SetAltitudeAGLFtIC(hagl);
00423   }
00424   fgic->SetPhiRadIC(ff);
00425 
00426 }
00427 
00428 /*****************************************************************************/
00429 
00430 void FGTrimAnalysisControl::Run(void) {
00431 
00432     // ... what's going on here ??
00433 }
00434 
00435 /*****************************************************************************/
00436 
00437 void FGTrimAnalysisControl::setThrottlesPct(void) {
00438   double tMin,tMax;
00439   for(unsigned i=0;i<fdmex->GetPropulsion()->GetNumEngines();i++) {
00440       tMin=fdmex->GetPropulsion()->GetEngine(i)->GetThrottleMin();
00441       tMax=fdmex->GetPropulsion()->GetEngine(i)->GetThrottleMax();
00442       //cout << "setThrottlespct: " << i << ", " << control_min << ", " << control_max << ", " << control_value;
00443       fdmex->GetFCS()->SetThrottleCmd(i,tMin+control_value*(tMax-tMin));
00444       //cout << "setThrottlespct: " << fdmex->GetFCS()->GetThrottleCmd(i) << endl;
00445       fdmex->SuspendIntegration();
00446       fdmex->Initialize(fgic);
00447       fdmex->Run();
00448       fdmex->ResumeIntegration();
00449       fdmex->GetPropulsion()->GetSteadyState();
00450   }
00451 }
00452 
00453 /*****************************************************************************/
00454 
00455 
00456 /*****************************************************************************/
00457 //    The bitmasked value choices are as follows:
00458 //    unset: In this case (the default) JSBSim would only print
00459 //       out the normally expected messages, essentially echoing
00460 //       the config files as they are read. If the environment
00461 //       variable is not set, debug_lvl is set to 1 internally
00462 //    0: This requests JSBSim not to output any messages
00463 //       whatsoever.
00464 //    1: This value explicity requests the normal JSBSim
00465 //       startup messages
00466 //    2: This value asks for a message to be printed out when
00467 //       a class is instantiated
00468 //    4: When this value is set, a message is displayed when a
00469 //       FGModel object executes its Run() method
00470 //    8: When this value is set, various runtime state variables
00471 //       are printed out periodically
00472 //    16: When set various parameters are sanity checked and
00473 //       a message is printed out when they go out of bounds
00474 
00475 void FGTrimAnalysisControl::Debug(int from)
00476 {
00477 
00478   if (debug_lvl <= 0) return;
00479   if (debug_lvl & 1 ) { // Standard console startup message output
00480     if (from == 0) { // Constructor
00481 
00482     }
00483   }
00484   if (debug_lvl & 2 ) { // Instantiation/Destruction notification
00485     if (from == 0) cout << "Instantiated: FGTrimAnalysisControl" << endl;
00486     if (from == 1) cout << "Destroyed:    FGTrimAnalysisControl" << endl;
00487   }
00488   if (debug_lvl & 4 ) { // Run() method entry print for FGModel-derived objects
00489   }
00490   if (debug_lvl & 8 ) { // Runtime state variables
00491   }
00492   if (debug_lvl & 16) { // Sanity checking
00493   }
00494   if (debug_lvl & 64) {
00495     if (from == 0) { // Constructor
00496       cout << IdSrc << endl;
00497       cout << IdHdr << endl;
00498     }
00499   }
00500 }
00501 }