Commercial Aircraft Range Maximization by Differential Inclusion#
In this example we seek to maximize the range of a commercial aircraft flying at a constant Mach number with a fixed initial fuel load. The initial and final altitudes are constrained to 10000 ft (that is, in this simplified example we’re ignoring the aircraft behavior near the ground).
This example is simplified in several ways. The Mach number is fixed at 0.8. When only considering fuel burn, the aircraft will tend to fly slower to conserve fuel, ignoring the running costs of being airborne.
For a more detailed optimization of commercial aircraft performance, see Betts [CBC95].
Differential Inclusion#
In this case we also use the steady flight assumption to allow the problem to be parameterized with a different control set. Rather than providing throttle/thrust and angle of attack as control variables, we use a differential inclusion approach. Whereas altitude and velocity are often treated as state variables to be integrated, here we take a slightly different approach.
Rather than using the differential equations to compute rates for the flight path angle and true airspeed, we use a nonlinear solver to determine the alpha and thrust time-history which is needed to make the given altitude and airspeed time-history possible. This technique is known as differential inclusion. It was demonstrated by Seywald [CSey94]. Since the collocation techniques in Dymos are based on polynomial approximations, rates of the control variables can be easily calculated, which can then be applied to the reordered dynamics equations to solve for angle of attack and thrust.
The use of collocation/psuedospectral techniques in solving problems via differential inclusion has been examined by Fahroo and Ross [CFR01] and Kumar and Seywald [CKS96]. Since the differential inclusion approach demonstrated here relies on nonlinear solvers running within the ODE model, obtaining accurate derivatives is paramount. Computing derivatives with finite differences in the presence of nonlinear solvers can lead to inaccurate derivatives based on the properties of the solver. With OpenMDAO’s unified derivatives framework, the sensitivy to solver tolerances is greatly reduced.
Pitfalls of Differential Inclusion#
Differential inclusion works well in this case but it is highly dependent on the parameterization of the problem. For instance, it is possible for the optimizer to suggest a flight profile that requires thrust well beyond the capabilities of the vehicle. In our case, the throttle parameter \(\tau\) is path constrained within the optimizer, so that the optimizer can make the solution physically attainable by our model.
We have two options for how to tackle this problem. We could specify altitude as a control variable, and use it’s approximate derivative (\(alt\_rate\)) to provide the altitude rate for the flight equilibrium conditions. In practice, however, this noisy derivative of altitude can cause numerical difficulties. Other times, the optimizer will shape the control to satisfy our equilibrium conditions at the collocation nodes but be wildly wrong in between these nodes. One simple way to get around this is to still integrate altitude, but to provide \(climb\_rate\) as the control variable. This makes it easy to constraint \(climb\_rate\) with bounds constraints, and the time-history of altitude will be smoother since it is the time integral of the rate of climb.
Solver Options#
In general, openmdao.api.AnalysisError should be raised when the design variables provided by the optimizer result in nonconvergence of the solvers. This informs pyoptsparse to backtrack along its line search for a solution. This error can be raised by the compute method of a component. The nonlinear solver can also raise this error if it fails to converge by setting the following options:
self.nonlinear_solver.options['maxiter'] = 150
self.nonlinear_solver.options['err_on_maxiter'] = True
In the course of solving, an outer NewtonSolver won’t necessarily force subsystems to be solved, instead relying on the derivative at the at the previous solution. In the context of OpenMDAO, even explicit components constitute subsystems. This can lead to poorly converged solutions.
Setting the solve_subystems option to True will force subsystems to be updated, making sure that all derivative information is up-to-date.
self.nonlinear_solver.options['solve_subsystems'] = True
self.nonlinear_solver.options['max_sub_solves'] = 10
The Impact of Choice of Transcription#
The high-order Gauss-Lobatto method used by Dymos evaluates the ordinary differential equations in a two-step process. First, the ODEs are evaluated at the state discretization nodes. Using the state values and derivatives, values and rates at the collocation nodes are then interpolated. If a large rate is found at one of the state discretization nodes, the interpolation at the collocation node can be outside of the expected range, leading to nonconvergence of solvers within the ODEs.
There is also an overhead cost to invoking solvers. While differential inclusion reduces the number of design variables for the optimizer, it also slows down each evaluation of the ODE somewhat. Since evaluating the collocation constraints in the Gauss-Lobatto method requires two evaluations of the ODE, it typically suffers greater slowdown that the Radau Pseudospectral method, which evaluates the defects of the ODE with a single evaluation that encompasses all the nodes.
For these reasons the Radau Pseudospectral transcription is typically preferable when using a differential inclusion approach to problems.
Problem Formulation#
State Variables#
Name |
Description |
Fixed Initial Value |
Fixed Final Value |
---|---|---|---|
range |
distance flown along ground |
True (0 NM) |
False |
mass_fuel |
mass of fuel aboard aircraft |
True (30000 lbm) |
True (0 lbm) |
alt |
aircraft altitude |
True (10000 ft) |
True (10000 ft) |
Dynamic Controls#
Name |
Description |
Optimized |
Fixed Initial Value |
Fixed Final Value |
---|---|---|---|---|
climb_rate |
aircraft rate of climb |
True |
False |
False |
Parameters#
Name |
Description |
Optimized |
---|---|---|
mach |
mach number |
False (0.8) |
S |
aerodynamic reference area |
False (427.8 m**2) |
mass_empty |
aircraft empty mass |
False (330693.393 lbm) |
mass_payload |
aircraft payload mass |
False (74100 lbm) |
Objective#
Name |
Description |
Location |
Minimized or Maximized |
---|---|---|---|
r |
range |
final |
Maximized |
Nonlinear Path Constraints#
Name |
Description |
Lower |
Upper |
---|---|---|---|
tau |
engine throttle parameter |
0.01 |
1.0 |
Nonlinear Boundary Constraints#
None
Models#
Atmosphere#
This problem uses an analytic fit to the 1976 standard atmosphere.
Name |
Description |
Input or Output |
---|---|---|
alt |
altitude (m) |
input |
pres |
static pressure (Pa) |
output |
temp |
temperature (K) |
output |
sos |
speed of sound (m/s) |
output |
rho |
density (kg/m**3) |
output |
True Airspeed#
TrueAirspeedComp uses the Mach number, provided as a control, and the speed of sound from the atmosphere model to compute the true airspeed of the aircraft.
Name |
Description |
Input or Output |
---|---|---|
mach |
Mach number |
input |
sos |
speed of sound (m/s) |
input |
TAS |
true airspeed (m/s) |
output |
Flight Path Angle#
SteadyFlightPathAngleComp uses the true airspeed and the climb rate, obtained by differentiating the altitude time history at the nodes, to compute the flight path angle.
Name |
Description |
Input or Output |
---|---|---|
TAS |
true airspeed (m/s) |
input |
alt_rate |
climb rate (m/s) |
input |
gamma |
flight path angle (rad) |
output |
Range Rate#
RangeRateComp uses the true airspeed and the flight path angle to determine the velocity projected along the ground. This is the derivative of the state variable range.
Name |
Description |
Input or Output |
---|---|---|
TAS |
true airspeed (m/s) |
input |
gamma |
flight path angle (rad) |
input |
dXdt:range |
range rate (m/s) |
output |
Mass#
The component MassComp defined in
mass_comp.py
computes the aircraft total mass based on
its empty mass, payload mass, and current fuel mass. It also computes
total weight which simplifies some equations later on.
Name |
Description |
Input or Output |
---|---|---|
mass_empty |
aircraft empty mass (kg) |
input |
mass_payload |
payload mass (kg) |
input |
mass_fuel |
fuel mass (kg) |
input |
mass_total |
total aircraft mass (kg) |
output |
W_total |
total aircraft weight (N) |
output |
Dynamic Pressure#
The DynamicPressureComp computes the dynamic pressure from true airspeed and atmospheric density.
Name |
Description |
Input or Output |
---|---|---|
TAS |
true airspeed (m/s) |
input |
rho |
atmospheric density (kg/m**3) |
input |
q |
dynamic pressure (Pa) |
output |
Aerodynamics#
The aerodynamics group computes the aerodynamic coefficients and forces on the vehicle. It consists of an interpolation component which outputs lift, drag, and moment coefficients as a function of Mach number, angle of attack, altitude, and tail rotation angle. A second component then uses these coefficients, along with dynamic pressure and aerodynamic reference area, to compute the lift and drag forces on the vehicle.
The aerodynamics group resides within the flight equilibrium group. As that group iterates to find the combination of thrust coefficient, angle of attack, and tail rotation angle, aerodynamics needs to update the values of the interpolated coefficients and resulting forces.
Organizationally speaking, we logically could have put the dynamic pressure component within the aerodynamics group. However, since that group doesn’t need to be updated with changes in alpha and tail angle, it’s more efficient to leave it outside of flight equilibrium group.
Name |
Description |
Input or Output |
---|---|---|
mach |
mach number |
input |
alt |
altitude (m) |
input |
alpha |
angle of attack (deg) |
input |
eta |
tail rotation angle (deg) |
input |
CL |
lift coefficient |
output |
CD |
drag coefficient |
output |
CM |
moment coefficient |
output |
L |
lift force (N) |
output |
D |
drag force (N) |
output |
Note
This example uses openmdao.api.MetaModelStructuredComp to interpolate aerodynamic properties of the vehicle. This component is somewhat easier to use since it is distributed as part of OpenMDAO, but it can be significantly slower than alternatives such as SMT.
Flight Equilibrium#
The steady flight equilibrium group uses balances to solve for the angle of attack and tail plane rotation angle such that the aircraft is in steady flight (the rates of change in flight path angle and true airspeed are zero) and the aerodynamic moment in the pitch axis is zero.
Of course, in reality the vehicle will accelerate, but the flight profile being modeled is so benign that assuming steady flight at discrete points (nodes) in the trajectory is not terribly inaccurate.
The thrust necessary for steady flight is computed by balancing the drag equation
The lift coefficient required for steady flight is found by balancing lift and weight:
Using coefficients in the balance equations is better scaled from a numerical standpoint.
Propulsion#
Having determined the thrust, the propulsion group then computes the rate of fuel burn. In addition, by normalizing thrust at any point by the maximum possible thrust, we obtain the throttle parameter \(\tau\). The propulsion group uses a number of components to perform these calculations.
Maximum thrust is computed by multiplying sea-level thrust by the ratio of pressure to sea-level atmospheric pressure.
The throttle parameter is then the ratio current thrust to maximum possible thrust.
The thrust specific fuel consumption is computed as follows:
Finally, fuel burn rate is:
The ODE System: aircraft_ode.py#
class AircraftODE(om.Group):
def initialize(self):
self.options.declare('num_nodes', types=int,
desc='Number of nodes to be evaluated in the RHS')
def setup(self):
nn = self.options['num_nodes']
self.add_subsystem(name='mass_comp',
subsys=MassComp(num_nodes=nn),
promotes_inputs=['mass_fuel'])
self.connect('mass_comp.W_total', 'flight_equilibrium.W_total')
self.add_subsystem(name='atmos',
subsys=USatm1976Comp(num_nodes=nn),
promotes_inputs=[('h', 'alt')])
self.connect('atmos.pres', 'propulsion.pres')
self.connect('atmos.sos', 'tas_comp.sos')
self.connect('atmos.rho', 'q_comp.rho')
self.add_subsystem(name='tas_comp',
subsys=TrueAirspeedComp(num_nodes=nn))
self.connect('tas_comp.TAS',
('gam_comp.TAS', 'q_comp.TAS', 'range_rate_comp.TAS'))
self.add_subsystem(name='gam_comp',
subsys=SteadyFlightPathAngleComp(num_nodes=nn))
self.connect('gam_comp.gam', ('flight_equilibrium.gam', 'range_rate_comp.gam'))
self.add_subsystem(name='q_comp',
subsys=DynamicPressureComp(num_nodes=nn))
self.connect('q_comp.q', ('aero.q', 'flight_equilibrium.q', 'propulsion.q'))
self.add_subsystem(name='flight_equilibrium',
subsys=SteadyFlightEquilibriumGroup(num_nodes=nn),
promotes_inputs=['aero.*', 'alt'],
promotes_outputs=['aero.*'])
self.connect('flight_equilibrium.CT', 'propulsion.CT')
self.add_subsystem(name='propulsion', subsys=PropulsionGroup(num_nodes=nn),
promotes_inputs=['alt'])
self.add_subsystem(name='range_rate_comp', subsys=RangeRateComp(num_nodes=nn))
# We promoted multiple inputs to the group named 'alt'.
# In order to automatically create a source variable for 'alt', we must specify their units
# since they have different units in different components.
self.set_input_defaults('alt', val=np.ones(nn,), units='m')
In this case the system has only two integrated states: range and mass_fuel. There are six parameters. Two of them (alt and climb_rate) will be varied dynamically in the phase The other four (mach, S, mass_empty, and mass_payload), will be set to fixed values as non-optimized parameters. More details on the various models involved can be found in the examples code.
Building the problem#
In the following code we define and solve the optimal control problem. Note that we demonstrate the use of externally-connected design parameters in this case. The four parameters are connected to a source provided by the assumptions IndepVarComp.
Show code cell outputs
import numpy as np
import openmdao.api as om
from dymos.models.atmosphere import USatm1976Comp
from .steady_flight_path_angle_comp import SteadyFlightPathAngleComp
from .dynamic_pressure_comp import DynamicPressureComp
from .flight_equlibrium.steady_flight_equilibrium_group import SteadyFlightEquilibriumGroup
from .propulsion.propulsion_group import PropulsionGroup
from .range_rate_comp import RangeRateComp
from .true_airspeed_comp import TrueAirspeedComp
from .mass_comp import MassComp
class AircraftODE(om.Group):
def initialize(self):
self.options.declare('num_nodes', types=int,
desc='Number of nodes to be evaluated in the RHS')
def setup(self):
nn = self.options['num_nodes']
self.add_subsystem(name='mass_comp',
subsys=MassComp(num_nodes=nn),
promotes_inputs=['mass_fuel'])
self.connect('mass_comp.W_total', 'flight_equilibrium.W_total')
self.add_subsystem(name='atmos',
subsys=USatm1976Comp(num_nodes=nn),
promotes_inputs=[('h', 'alt')])
self.connect('atmos.pres', 'propulsion.pres')
self.connect('atmos.sos', 'tas_comp.sos')
self.connect('atmos.rho', 'q_comp.rho')
self.add_subsystem(name='tas_comp',
subsys=TrueAirspeedComp(num_nodes=nn))
self.connect('tas_comp.TAS',
('gam_comp.TAS', 'q_comp.TAS', 'range_rate_comp.TAS'))
self.add_subsystem(name='gam_comp',
subsys=SteadyFlightPathAngleComp(num_nodes=nn))
self.connect('gam_comp.gam', ('flight_equilibrium.gam', 'range_rate_comp.gam'))
self.add_subsystem(name='q_comp',
subsys=DynamicPressureComp(num_nodes=nn))
self.connect('q_comp.q', ('aero.q', 'flight_equilibrium.q', 'propulsion.q'))
self.add_subsystem(name='flight_equilibrium',
subsys=SteadyFlightEquilibriumGroup(num_nodes=nn),
promotes_inputs=['aero.*', 'alt'],
promotes_outputs=['aero.*'])
self.connect('flight_equilibrium.CT', 'propulsion.CT')
self.add_subsystem(name='propulsion', subsys=PropulsionGroup(num_nodes=nn),
promotes_inputs=['alt'])
self.add_subsystem(name='range_rate_comp', subsys=RangeRateComp(num_nodes=nn))
# We promoted multiple inputs to the group named 'alt'.
# In order to automatically create a source variable for 'alt', we must specify their units
# since they have different units in different components.
self.set_input_defaults('alt', val=np.ones(nn,), units='m')
import matplotlib.pyplot as plt
import openmdao.api as om
import dymos as dm
from dymos.examples.aircraft_steady_flight.aircraft_ode import AircraftODE
from dymos.examples.plotting import plot_results
from dymos.utils.lgl import lgl
p = om.Problem(model=om.Group())
p.driver = om.pyOptSparseDriver()
p.driver.options['optimizer'] = 'IPOPT'
p.driver.opt_settings['print_level'] = 0
p.driver.declare_coloring()
num_seg = 15
seg_ends, _ = lgl(num_seg + 1)
traj = p.model.add_subsystem('traj', dm.Trajectory())
phase = traj.add_phase('phase0',
dm.Phase(ode_class=AircraftODE,
transcription=dm.Radau(num_segments=num_seg,
segment_ends=seg_ends,
order=3, compressed=False)))
# Pass Reference Area from an external source
assumptions = p.model.add_subsystem('assumptions', om.IndepVarComp())
assumptions.add_output('S', val=427.8, units='m**2')
assumptions.add_output('mass_empty', val=1.0, units='kg')
assumptions.add_output('mass_payload', val=1.0, units='kg')
phase.set_time_options(fix_initial=True,
duration_bounds=(300, 10000),
duration_ref=1000)
phase.add_state('range', units='NM',
rate_source='range_rate_comp.dXdt:range',
fix_initial=True, fix_final=False, ref=1e-3,
defect_ref=1e-3, lower=0, upper=2000)
phase.add_state('mass_fuel', units='lbm',
rate_source='propulsion.dXdt:mass_fuel',
fix_initial=True, fix_final=True,
upper=1.5E5, lower=0.0, ref=1e2, defect_ref=1e2)
phase.add_state('alt', units='kft',
rate_source='climb_rate',
fix_initial=True, fix_final=True,
lower=0.0, upper=60, ref=1e-3, defect_ref=1e-3)
phase.add_control('climb_rate', units='ft/min', opt=True, lower=-3000, upper=3000,
targets=['gam_comp.climb_rate'],
rate_continuity=True, rate2_continuity=False)
phase.add_control('mach', targets=['tas_comp.mach', 'aero.mach'], units=None, opt=False)
phase.add_parameter('S',
targets=['aero.S', 'flight_equilibrium.S', 'propulsion.S'],
units='m**2')
phase.add_parameter('mass_empty', targets=['mass_comp.mass_empty'], units='kg')
phase.add_parameter('mass_payload', targets=['mass_comp.mass_payload'], units='kg')
phase.add_path_constraint('propulsion.tau', lower=0.01, upper=2.0)
p.model.connect('assumptions.S', 'traj.phase0.parameters:S')
p.model.connect('assumptions.mass_empty', 'traj.phase0.parameters:mass_empty')
p.model.connect('assumptions.mass_payload', 'traj.phase0.parameters:mass_payload')
phase.add_objective('range', loc='final', ref=-1.0)
phase.add_timeseries_output('aero.CL')
phase.add_timeseries_output('aero.CD')
p.setup()
--- Constraint Report [traj] ---
--- phase0 ---
[path] 1.0000e-02 <= propulsion.tau <= 2.0000e+00 [None]
<openmdao.core.problem.Problem at 0x7fb3f85fc670>
Set the initial guesses and assumptions#
p['traj.phase0.t_initial'] = 0.0
p['traj.phase0.t_duration'] = 3600.0
p['traj.phase0.states:range'] = phase.interp('range', ys=(0, 724.0))
p['traj.phase0.states:mass_fuel'] = phase.interp('mass_fuel', ys=(30000, 1e-3))
p['traj.phase0.states:alt'][:] = 10.0
p['traj.phase0.controls:mach'][:] = 0.8
p['assumptions.S'] = 427.8
p['assumptions.mass_empty'] = 0.15E6
p['assumptions.mass_payload'] = 84.02869 * 400
Run the problem#
The dymos run_problem
method is used to automatically run the driver and follow the driver execution with a simulation run to help verify the results.
Plots of the resulting solution and simulation will be generated.
dm.run_problem(p, simulate=True, make_plots=True)
Model viewer data has already been recorded for Driver.
Full total jacobian was computed 3 times, taking 1.365613 seconds.
Total jacobian shape: (224, 221)
Jacobian shape: (224, 221) ( 2.90% nonzero)
FWD solves: 12 REV solves: 0
Total colors vs. total size: 12 vs 221 (94.6% improvement)
Sparsity computed using tolerance: 1e-25
Time to compute sparsity: 1.365613 sec.
Time to compute coloring: 0.232396 sec.
Memory to compute coloring: 0.250000 MB.
Optimization Problem -- Optimization using pyOpt_sparse
================================================================================
Objective Function: _objfunc
Solution:
--------------------------------------------------------------------------------
Total Time: 28.0042
User Objective Time : 25.0847
User Sensitivity Time : 2.4292
Interface Time : 0.2170
Opt Solver Time: 0.2732
Calls to Objective Function : 29
Calls to Sens Function : 28
Objectives
Index Name Value
0 traj.phases.phase0.indep_states.states:range -7.260052E+02
Variables (c - continuous, i - integer, d - discrete)
Index Name Type Lower Bound Value Upper Bound Status
0 traj.phase0.t_duration_0 c 3.000000E-01 5.637564E+00 1.000000E+01
1 traj.phases.phase0.indep_states.states:range_0 c 0.000000E+00 4.300778E+03 2.000000E+06
2 traj.phases.phase0.indep_states.states:range_1 c 0.000000E+00 1.019478E+04 2.000000E+06
3 traj.phases.phase0.indep_states.states:range_2 c 0.000000E+00 1.205044E+04 2.000000E+06
4 traj.phases.phase0.indep_states.states:range_3 c 0.000000E+00 1.205044E+04 2.000000E+06
5 traj.phases.phase0.indep_states.states:range_4 c 0.000000E+00 2.179716E+04 2.000000E+06
6 traj.phases.phase0.indep_states.states:range_5 c 0.000000E+00 3.502521E+04 2.000000E+06
7 traj.phases.phase0.indep_states.states:range_6 c 0.000000E+00 3.915789E+04 2.000000E+06
8 traj.phases.phase0.indep_states.states:range_7 c 0.000000E+00 3.915789E+04 2.000000E+06
9 traj.phases.phase0.indep_states.states:range_8 c 0.000000E+00 5.336791E+04 2.000000E+06
10 traj.phases.phase0.indep_states.states:range_9 c 0.000000E+00 7.250027E+04 2.000000E+06
11 traj.phases.phase0.indep_states.states:range_10 c 0.000000E+00 7.845561E+04 2.000000E+06
12 traj.phases.phase0.indep_states.states:range_11 c 0.000000E+00 7.845561E+04 2.000000E+06
13 traj.phases.phase0.indep_states.states:range_12 c 0.000000E+00 9.625935E+04 2.000000E+06
14 traj.phases.phase0.indep_states.states:range_13 c 0.000000E+00 1.208337E+05 2.000000E+06
15 traj.phases.phase0.indep_states.states:range_14 c 0.000000E+00 1.286109E+05 2.000000E+06
16 traj.phases.phase0.indep_states.states:range_15 c 0.000000E+00 1.286109E+05 2.000000E+06
17 traj.phases.phase0.indep_states.states:range_16 c 0.000000E+00 1.498279E+05 2.000000E+06
18 traj.phases.phase0.indep_states.states:range_17 c 0.000000E+00 1.791031E+05 2.000000E+06
19 traj.phases.phase0.indep_states.states:range_18 c 0.000000E+00 1.883686E+05 2.000000E+06
20 traj.phases.phase0.indep_states.states:range_19 c 0.000000E+00 1.883686E+05 2.000000E+06
21 traj.phases.phase0.indep_states.states:range_20 c 0.000000E+00 2.121241E+05 2.000000E+06
22 traj.phases.phase0.indep_states.states:range_21 c 0.000000E+00 2.449018E+05 2.000000E+06
23 traj.phases.phase0.indep_states.states:range_22 c 0.000000E+00 2.552758E+05 2.000000E+06
24 traj.phases.phase0.indep_states.states:range_23 c 0.000000E+00 2.552758E+05 2.000000E+06
25 traj.phases.phase0.indep_states.states:range_24 c 0.000000E+00 2.805971E+05 2.000000E+06
26 traj.phases.phase0.indep_states.states:range_25 c 0.000000E+00 3.155354E+05 2.000000E+06
27 traj.phases.phase0.indep_states.states:range_26 c 0.000000E+00 3.265932E+05 2.000000E+06
28 traj.phases.phase0.indep_states.states:range_27 c 0.000000E+00 3.265932E+05 2.000000E+06
29 traj.phases.phase0.indep_states.states:range_28 c 0.000000E+00 3.524437E+05 2.000000E+06
30 traj.phases.phase0.indep_states.states:range_29 c 0.000000E+00 3.881121E+05 2.000000E+06
31 traj.phases.phase0.indep_states.states:range_30 c 0.000000E+00 3.994010E+05 2.000000E+06
32 traj.phases.phase0.indep_states.states:range_31 c 0.000000E+00 3.994010E+05 2.000000E+06
33 traj.phases.phase0.indep_states.states:range_32 c 0.000000E+00 4.247223E+05 2.000000E+06
34 traj.phases.phase0.indep_states.states:range_33 c 0.000000E+00 4.596605E+05 2.000000E+06
35 traj.phases.phase0.indep_states.states:range_34 c 0.000000E+00 4.707183E+05 2.000000E+06
36 traj.phases.phase0.indep_states.states:range_35 c 0.000000E+00 4.707183E+05 2.000000E+06
37 traj.phases.phase0.indep_states.states:range_36 c 0.000000E+00 4.944738E+05 2.000000E+06
38 traj.phases.phase0.indep_states.states:range_37 c 0.000000E+00 5.272515E+05 2.000000E+06
39 traj.phases.phase0.indep_states.states:range_38 c 0.000000E+00 5.376255E+05 2.000000E+06
40 traj.phases.phase0.indep_states.states:range_39 c 0.000000E+00 5.376255E+05 2.000000E+06
41 traj.phases.phase0.indep_states.states:range_40 c 0.000000E+00 5.588425E+05 2.000000E+06
42 traj.phases.phase0.indep_states.states:range_41 c 0.000000E+00 5.881176E+05 2.000000E+06
43 traj.phases.phase0.indep_states.states:range_42 c 0.000000E+00 5.973831E+05 2.000000E+06
44 traj.phases.phase0.indep_states.states:range_43 c 0.000000E+00 5.973831E+05 2.000000E+06
45 traj.phases.phase0.indep_states.states:range_44 c 0.000000E+00 6.151931E+05 2.000000E+06
46 traj.phases.phase0.indep_states.states:range_45 c 0.000000E+00 6.397625E+05 2.000000E+06
47 traj.phases.phase0.indep_states.states:range_46 c 0.000000E+00 6.475361E+05 2.000000E+06
48 traj.phases.phase0.indep_states.states:range_47 c 0.000000E+00 6.475361E+05 2.000000E+06
49 traj.phases.phase0.indep_states.states:range_48 c 0.000000E+00 6.612680E+05 2.000000E+06
50 traj.phases.phase0.indep_states.states:range_49 c 0.000000E+00 6.805947E+05 2.000000E+06
51 traj.phases.phase0.indep_states.states:range_50 c 0.000000E+00 6.868473E+05 2.000000E+06
52 traj.phases.phase0.indep_states.states:range_51 c 0.000000E+00 6.868473E+05 2.000000E+06
53 traj.phases.phase0.indep_states.states:range_52 c 0.000000E+00 6.963491E+05 2.000000E+06
54 traj.phases.phase0.indep_states.states:range_53 c 0.000000E+00 7.096820E+05 2.000000E+06
55 traj.phases.phase0.indep_states.states:range_54 c 0.000000E+00 7.139547E+05 2.000000E+06
56 traj.phases.phase0.indep_states.states:range_55 c 0.000000E+00 7.139547E+05 2.000000E+06
57 traj.phases.phase0.indep_states.states:range_56 c 0.000000E+00 7.182109E+05 2.000000E+06
58 traj.phases.phase0.indep_states.states:range_57 c 0.000000E+00 7.241240E+05 2.000000E+06
59 traj.phases.phase0.indep_states.states:range_58 c 0.000000E+00 7.260052E+05 2.000000E+06
60 traj.phases.phase0.indep_states.states:mass_fuel_0 c 0.000000E+00 2.957208E+02 1.500000E+03
61 traj.phases.phase0.indep_states.states:mass_fuel_1 c 0.000000E+00 2.900371E+02 1.500000E+03
62 traj.phases.phase0.indep_states.states:mass_fuel_2 c 0.000000E+00 2.882890E+02 1.500000E+03
63 traj.phases.phase0.indep_states.states:mass_fuel_3 c 0.000000E+00 2.882890E+02 1.500000E+03
64 traj.phases.phase0.indep_states.states:mass_fuel_4 c 0.000000E+00 2.794117E+02 1.500000E+03
65 traj.phases.phase0.indep_states.states:mass_fuel_5 c 0.000000E+00 2.680940E+02 1.500000E+03
66 traj.phases.phase0.indep_states.states:mass_fuel_6 c 0.000000E+00 2.647067E+02 1.500000E+03
67 traj.phases.phase0.indep_states.states:mass_fuel_7 c 0.000000E+00 2.647067E+02 1.500000E+03
68 traj.phases.phase0.indep_states.states:mass_fuel_8 c 0.000000E+00 2.535688E+02 1.500000E+03
69 traj.phases.phase0.indep_states.states:mass_fuel_9 c 0.000000E+00 2.401846E+02 1.500000E+03
70 traj.phases.phase0.indep_states.states:mass_fuel_10 c 0.000000E+00 2.364961E+02 1.500000E+03
71 traj.phases.phase0.indep_states.states:mass_fuel_11 c 0.000000E+00 2.364961E+02 1.500000E+03
72 traj.phases.phase0.indep_states.states:mass_fuel_12 c 0.000000E+00 2.271299E+02 1.500000E+03
73 traj.phases.phase0.indep_states.states:mass_fuel_13 c 0.000000E+00 2.167574E+02 1.500000E+03
74 traj.phases.phase0.indep_states.states:mass_fuel_14 c 0.000000E+00 2.137292E+02 1.500000E+03
75 traj.phases.phase0.indep_states.states:mass_fuel_15 c 0.000000E+00 2.137292E+02 1.500000E+03
76 traj.phases.phase0.indep_states.states:mass_fuel_16 c 0.000000E+00 2.055422E+02 1.500000E+03
77 traj.phases.phase0.indep_states.states:mass_fuel_17 c 0.000000E+00 1.942336E+02 1.500000E+03
78 traj.phases.phase0.indep_states.states:mass_fuel_18 c 0.000000E+00 1.906363E+02 1.500000E+03
79 traj.phases.phase0.indep_states.states:mass_fuel_19 c 0.000000E+00 1.906363E+02 1.500000E+03
80 traj.phases.phase0.indep_states.states:mass_fuel_20 c 0.000000E+00 1.813581E+02 1.500000E+03
81 traj.phases.phase0.indep_states.states:mass_fuel_21 c 0.000000E+00 1.684992E+02 1.500000E+03
82 traj.phases.phase0.indep_states.states:mass_fuel_22 c 0.000000E+00 1.644355E+02 1.500000E+03
83 traj.phases.phase0.indep_states.states:mass_fuel_23 c 0.000000E+00 1.644355E+02 1.500000E+03
84 traj.phases.phase0.indep_states.states:mass_fuel_24 c 0.000000E+00 1.545567E+02 1.500000E+03
85 traj.phases.phase0.indep_states.states:mass_fuel_25 c 0.000000E+00 1.410094E+02 1.500000E+03
86 traj.phases.phase0.indep_states.states:mass_fuel_26 c 0.000000E+00 1.367383E+02 1.500000E+03
87 traj.phases.phase0.indep_states.states:mass_fuel_27 c 0.000000E+00 1.367383E+02 1.500000E+03
88 traj.phases.phase0.indep_states.states:mass_fuel_28 c 0.000000E+00 1.267738E+02 1.500000E+03
89 traj.phases.phase0.indep_states.states:mass_fuel_29 c 0.000000E+00 1.130185E+02 1.500000E+03
90 traj.phases.phase0.indep_states.states:mass_fuel_30 c 0.000000E+00 1.086503E+02 1.500000E+03
91 traj.phases.phase0.indep_states.states:mass_fuel_31 c 0.000000E+00 1.086503E+02 1.500000E+03
92 traj.phases.phase0.indep_states.states:mass_fuel_32 c 0.000000E+00 9.882576E+01 1.500000E+03
93 traj.phases.phase0.indep_states.states:mass_fuel_33 c 0.000000E+00 8.534597E+01 1.500000E+03
94 traj.phases.phase0.indep_states.states:mass_fuel_34 c 0.000000E+00 8.113325E+01 1.500000E+03
95 traj.phases.phase0.indep_states.states:mass_fuel_35 c 0.000000E+00 8.113325E+01 1.500000E+03
96 traj.phases.phase0.indep_states.states:mass_fuel_36 c 0.000000E+00 7.218114E+01 1.500000E+03
97 traj.phases.phase0.indep_states.states:mass_fuel_37 c 0.000000E+00 5.977772E+01 1.500000E+03
98 traj.phases.phase0.indep_states.states:mass_fuel_38 c 0.000000E+00 5.576977E+01 1.500000E+03
99 traj.phases.phase0.indep_states.states:mass_fuel_39 c 0.000000E+00 5.576977E+01 1.500000E+03
100 traj.phases.phase0.indep_states.states:mass_fuel_40 c 0.000000E+00 4.741802E+01 1.500000E+03
101 traj.phases.phase0.indep_states.states:mass_fuel_41 c 0.000000E+00 3.592792E+01 1.500000E+03
102 traj.phases.phase0.indep_states.states:mass_fuel_42 c 0.000000E+00 3.239650E+01 1.500000E+03
103 traj.phases.phase0.indep_states.states:mass_fuel_43 c 0.000000E+00 3.239650E+01 1.500000E+03
104 traj.phases.phase0.indep_states.states:mass_fuel_44 c 0.000000E+00 2.595940E+01 1.500000E+03
105 traj.phases.phase0.indep_states.states:mass_fuel_45 c 0.000000E+00 1.877799E+01 1.500000E+03
106 traj.phases.phase0.indep_states.states:mass_fuel_46 c 0.000000E+00 1.714586E+01 1.500000E+03
107 traj.phases.phase0.indep_states.states:mass_fuel_47 c 0.000000E+00 1.714586E+01 1.500000E+03
108 traj.phases.phase0.indep_states.states:mass_fuel_48 c 0.000000E+00 1.500522E+01 1.500000E+03
109 traj.phases.phase0.indep_states.states:mass_fuel_49 c 0.000000E+00 1.243555E+01 1.500000E+03
110 traj.phases.phase0.indep_states.states:mass_fuel_50 c 0.000000E+00 1.141787E+01 1.500000E+03
111 traj.phases.phase0.indep_states.states:mass_fuel_51 c 0.000000E+00 1.141787E+01 1.500000E+03
112 traj.phases.phase0.indep_states.states:mass_fuel_52 c 0.000000E+00 9.500576E+00 1.500000E+03
113 traj.phases.phase0.indep_states.states:mass_fuel_53 c 0.000000E+00 5.938689E+00 1.500000E+03
114 traj.phases.phase0.indep_states.states:mass_fuel_54 c 0.000000E+00 4.560755E+00 1.500000E+03
115 traj.phases.phase0.indep_states.states:mass_fuel_55 c 0.000000E+00 4.560755E+00 1.500000E+03
116 traj.phases.phase0.indep_states.states:mass_fuel_56 c 0.000000E+00 3.066294E+00 1.500000E+03
117 traj.phases.phase0.indep_states.states:mass_fuel_57 c 0.000000E+00 7.801232E-01 1.500000E+03
118 traj.phases.phase0.indep_states.states:alt_0 c 0.000000E+00 1.152283E+04 6.000000E+04
119 traj.phases.phase0.indep_states.states:alt_1 c 0.000000E+00 1.362403E+04 6.000000E+04
120 traj.phases.phase0.indep_states.states:alt_2 c 0.000000E+00 1.428905E+04 6.000000E+04
121 traj.phases.phase0.indep_states.states:alt_3 c 0.000000E+00 1.428905E+04 6.000000E+04
122 traj.phases.phase0.indep_states.states:alt_4 c 0.000000E+00 1.781028E+04 6.000000E+04
123 traj.phases.phase0.indep_states.states:alt_5 c 0.000000E+00 2.266886E+04 6.000000E+04
124 traj.phases.phase0.indep_states.states:alt_6 c 0.000000E+00 2.420659E+04 6.000000E+04
125 traj.phases.phase0.indep_states.states:alt_7 c 0.000000E+00 2.420659E+04 6.000000E+04
126 traj.phases.phase0.indep_states.states:alt_8 c 0.000000E+00 2.947162E+04 6.000000E+04
127 traj.phases.phase0.indep_states.states:alt_9 c 0.000000E+00 3.563844E+04 6.000000E+04
128 traj.phases.phase0.indep_states.states:alt_10 c 0.000000E+00 3.710397E+04 6.000000E+04
129 traj.phases.phase0.indep_states.states:alt_11 c 0.000000E+00 3.710397E+04 6.000000E+04
130 traj.phases.phase0.indep_states.states:alt_12 c 0.000000E+00 3.977281E+04 6.000000E+04
131 traj.phases.phase0.indep_states.states:alt_13 c 0.000000E+00 4.064815E+04 6.000000E+04
132 traj.phases.phase0.indep_states.states:alt_14 c 0.000000E+00 4.061936E+04 6.000000E+04
133 traj.phases.phase0.indep_states.states:alt_15 c 0.000000E+00 4.061936E+04 6.000000E+04
134 traj.phases.phase0.indep_states.states:alt_16 c 0.000000E+00 4.049114E+04 6.000000E+04
135 traj.phases.phase0.indep_states.states:alt_17 c 0.000000E+00 4.035335E+04 6.000000E+04
136 traj.phases.phase0.indep_states.states:alt_18 c 0.000000E+00 4.033653E+04 6.000000E+04
137 traj.phases.phase0.indep_states.states:alt_19 c 0.000000E+00 4.033653E+04 6.000000E+04
138 traj.phases.phase0.indep_states.states:alt_20 c 0.000000E+00 4.037083E+04 6.000000E+04
139 traj.phases.phase0.indep_states.states:alt_21 c 0.000000E+00 4.051554E+04 6.000000E+04
140 traj.phases.phase0.indep_states.states:alt_22 c 0.000000E+00 4.056201E+04 6.000000E+04
141 traj.phases.phase0.indep_states.states:alt_23 c 0.000000E+00 4.056201E+04 6.000000E+04
142 traj.phases.phase0.indep_states.states:alt_24 c 0.000000E+00 4.064374E+04 6.000000E+04
143 traj.phases.phase0.indep_states.states:alt_25 c 0.000000E+00 4.069400E+04 6.000000E+04
144 traj.phases.phase0.indep_states.states:alt_26 c 0.000000E+00 4.069891E+04 6.000000E+04
145 traj.phases.phase0.indep_states.states:alt_27 c 0.000000E+00 4.069891E+04 6.000000E+04
146 traj.phases.phase0.indep_states.states:alt_28 c 0.000000E+00 4.070253E+04 6.000000E+04
147 traj.phases.phase0.indep_states.states:alt_29 c 0.000000E+00 4.075182E+04 6.000000E+04
148 traj.phases.phase0.indep_states.states:alt_30 c 0.000000E+00 4.079373E+04 6.000000E+04
149 traj.phases.phase0.indep_states.states:alt_31 c 0.000000E+00 4.079373E+04 6.000000E+04
150 traj.phases.phase0.indep_states.states:alt_32 c 0.000000E+00 4.093549E+04 6.000000E+04
151 traj.phases.phase0.indep_states.states:alt_33 c 0.000000E+00 4.107641E+04 6.000000E+04
152 traj.phases.phase0.indep_states.states:alt_34 c 0.000000E+00 4.106576E+04 6.000000E+04
153 traj.phases.phase0.indep_states.states:alt_35 c 0.000000E+00 4.106576E+04 6.000000E+04
154 traj.phases.phase0.indep_states.states:alt_36 c 0.000000E+00 4.093919E+04 6.000000E+04
155 traj.phases.phase0.indep_states.states:alt_37 c 0.000000E+00 4.085650E+04 6.000000E+04
156 traj.phases.phase0.indep_states.states:alt_38 c 0.000000E+00 4.093624E+04 6.000000E+04
157 traj.phases.phase0.indep_states.states:alt_39 c 0.000000E+00 4.093624E+04 6.000000E+04
158 traj.phases.phase0.indep_states.states:alt_40 c 0.000000E+00 4.129894E+04 6.000000E+04
159 traj.phases.phase0.indep_states.states:alt_41 c 0.000000E+00 4.178378E+04 6.000000E+04
160 traj.phases.phase0.indep_states.states:alt_42 c 0.000000E+00 4.181540E+04 6.000000E+04
161 traj.phases.phase0.indep_states.states:alt_43 c 0.000000E+00 4.181540E+04 6.000000E+04
162 traj.phases.phase0.indep_states.states:alt_44 c 0.000000E+00 4.146282E+04 6.000000E+04
163 traj.phases.phase0.indep_states.states:alt_45 c 0.000000E+00 3.890556E+04 6.000000E+04
164 traj.phases.phase0.indep_states.states:alt_46 c 0.000000E+00 3.730605E+04 6.000000E+04
165 traj.phases.phase0.indep_states.states:alt_47 c 0.000000E+00 3.730605E+04 6.000000E+04
166 traj.phases.phase0.indep_states.states:alt_48 c 0.000000E+00 3.341148E+04 6.000000E+04
167 traj.phases.phase0.indep_states.states:alt_49 c 0.000000E+00 2.654155E+04 6.000000E+04
168 traj.phases.phase0.indep_states.states:alt_50 c 0.000000E+00 2.420659E+04 6.000000E+04
169 traj.phases.phase0.indep_states.states:alt_51 c 0.000000E+00 2.420659E+04 6.000000E+04
170 traj.phases.phase0.indep_states.states:alt_52 c 0.000000E+00 2.068536E+04 6.000000E+04
171 traj.phases.phase0.indep_states.states:alt_53 c 0.000000E+00 1.582678E+04 6.000000E+04
172 traj.phases.phase0.indep_states.states:alt_54 c 0.000000E+00 1.428905E+04 6.000000E+04
173 traj.phases.phase0.indep_states.states:alt_55 c 0.000000E+00 1.428905E+04 6.000000E+04
174 traj.phases.phase0.indep_states.states:alt_56 c 0.000000E+00 1.276622E+04 6.000000E+04
175 traj.phases.phase0.indep_states.states:alt_57 c 0.000000E+00 1.066502E+04 6.000000E+04
176 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_0 c -3.000000E+03 3.000000E+03 3.000000E+03 U
177 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_1 c -3.000000E+03 3.000000E+03 3.000000E+03 U
178 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_2 c -3.000000E+03 3.000000E+03 3.000000E+03 U
179 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_3 c -3.000000E+03 3.000000E+03 3.000000E+03 U
180 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_4 c -3.000000E+03 3.000000E+03 3.000000E+03 U
181 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_5 c -3.000000E+03 3.000000E+03 3.000000E+03 U
182 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_6 c -3.000000E+03 3.000000E+03 3.000000E+03 U
183 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_7 c -3.000000E+03 2.834019E+03 3.000000E+03
184 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_8 c -3.000000E+03 2.059977E+03 3.000000E+03
185 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_9 c -3.000000E+03 1.683328E+03 3.000000E+03
186 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_10 c -3.000000E+03 6.801696E+02 3.000000E+03
187 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_11 c -3.000000E+03 8.660699E-01 3.000000E+03
188 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_12 c -3.000000E+03 -4.385213E+01 3.000000E+03
189 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_13 c -3.000000E+03 -4.577801E+01 3.000000E+03
190 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_14 c -3.000000E+03 -2.088393E+01 3.000000E+03
191 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_15 c -3.000000E+03 -6.349127E+00 3.000000E+03
192 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_16 c -3.000000E+03 2.513465E+01 3.000000E+03
193 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_17 c -3.000000E+03 3.612047E+01 3.000000E+03
194 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_18 c -3.000000E+03 3.175682E+01 3.000000E+03
195 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_19 c -3.000000E+03 1.816953E+01 3.000000E+03
196 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_20 c -3.000000E+03 4.887492E+00 3.000000E+03
197 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_21 c -3.000000E+03 2.004171E+00 3.000000E+03
198 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_22 c -3.000000E+03 2.209239E+00 3.000000E+03
199 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_23 c -3.000000E+03 2.286212E+01 3.000000E+03
200 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_24 c -3.000000E+03 3.431969E+01 3.000000E+03
201 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_25 c -3.000000E+03 4.562497E+01 3.000000E+03
202 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_26 c -3.000000E+03 5.244817E+00 3.000000E+03
203 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_27 c -3.000000E+03 -2.105889E+01 3.000000E+03
204 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_28 c -3.000000E+03 -4.837126E+01 3.000000E+03
205 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_29 c -3.000000E+03 3.274564E+01 3.000000E+03
206 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_30 c -3.000000E+03 8.711934E+01 3.000000E+03
207 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_31 c -3.000000E+03 1.552976E+02 3.000000E+03
208 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_32 c -3.000000E+03 6.174642E+01 3.000000E+03
209 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_33 c -3.000000E+03 -1.318852E+01 3.000000E+03
210 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_34 c -3.000000E+03 -3.453103E+02 3.000000E+03
211 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_35 c -3.000000E+03 -1.352390E+03 3.000000E+03
212 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_36 c -3.000000E+03 -1.803712E+03 3.000000E+03
213 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_37 c -3.000000E+03 -2.502393E+03 3.000000E+03
214 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_38 c -3.000000E+03 -2.971240E+03 3.000000E+03
215 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_39 c -3.000000E+03 -3.000000E+03 3.000000E+03 L
216 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_40 c -3.000000E+03 -3.000000E+03 3.000000E+03 L
217 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_41 c -3.000000E+03 -3.000000E+03 3.000000E+03 L
218 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_42 c -3.000000E+03 -3.000000E+03 3.000000E+03 L
219 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_43 c -3.000000E+03 -3.000000E+03 3.000000E+03 L
220 traj.phases.phase0.control_group.indep_controls.controls:climb_rate_44 c -3.000000E+03 -3.000000E+03 3.000000E+03 L
Constraints (i - inequality, e - equality)
Index Name Type Lower Value Upper Status Lagrange Multiplier (N/A)
0 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.380902E-12 0.000000E+00 9.00000E+100
1 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -1.190451E-12 0.000000E+00 9.00000E+100
2 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -3.571353E-12 0.000000E+00 9.00000E+100
3 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.651601E-11 0.000000E+00 9.00000E+100
4 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.752669E-12 0.000000E+00 9.00000E+100
5 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
6 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -4.612515E-11 0.000000E+00 9.00000E+100
7 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 8.386390E-12 0.000000E+00 9.00000E+100
8 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 8.805710E-11 0.000000E+00 9.00000E+100
9 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.966213E-10 0.000000E+00 9.00000E+100
10 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 4.533214E-10 0.000000E+00 9.00000E+100
11 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.638511E-11 0.000000E+00 9.00000E+100
12 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -3.253272E-11 0.000000E+00 9.00000E+100
13 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 4.489516E-10 0.000000E+00 9.00000E+100
14 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 5.530563E-10 0.000000E+00 9.00000E+100
15 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 3.715345E-10 0.000000E+00 9.00000E+100
16 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.913996E-11 0.000000E+00 9.00000E+100
17 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -2.185497E-11 0.000000E+00 9.00000E+100
18 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.251902E-10 0.000000E+00 9.00000E+100
19 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 4.736759E-10 0.000000E+00 9.00000E+100
20 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.009473E-10 0.000000E+00 9.00000E+100
21 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 3.805178E-10 0.000000E+00 9.00000E+100
22 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 9.750770E-10 0.000000E+00 9.00000E+100
23 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.140413E-10 0.000000E+00 9.00000E+100
24 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 3.882589E-11 0.000000E+00 9.00000E+100
25 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.048299E-09 0.000000E+00 9.00000E+100
26 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.630687E-10 0.000000E+00 9.00000E+100
27 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 8.013488E-11 0.000000E+00 9.00000E+100
28 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.216593E-09 0.000000E+00 9.00000E+100
29 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 3.351095E-10 0.000000E+00 9.00000E+100
30 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 4.554581E-10 0.000000E+00 9.00000E+100
31 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 3.773796E-10 0.000000E+00 9.00000E+100
32 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 7.092134E-10 0.000000E+00 9.00000E+100
33 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -4.260128E-10 0.000000E+00 9.00000E+100
34 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.403149E-10 0.000000E+00 9.00000E+100
35 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.293915E-10 0.000000E+00 9.00000E+100
36 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 9.225029E-11 0.000000E+00 9.00000E+100
37 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.515917E-10 0.000000E+00 9.00000E+100
38 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 2.809441E-10 0.000000E+00 9.00000E+100
39 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 5.505337E-12 0.000000E+00 9.00000E+100
40 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 5.505337E-11 0.000000E+00 9.00000E+100
41 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 1.101067E-11 0.000000E+00 9.00000E+100
42 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 -7.142707E-11 0.000000E+00 9.00000E+100
43 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 8.214113E-11 0.000000E+00 9.00000E+100
44 traj.phases.phase0.collocation_constraint.defects:range e 0.000000E+00 3.809444E-11 0.000000E+00 9.00000E+100
45 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -2.491376E-13 0.000000E+00 9.00000E+100
46 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -4.952277E-14 0.000000E+00 9.00000E+100
47 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -4.876088E-14 0.000000E+00 9.00000E+100
48 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.321281E-13 0.000000E+00 9.00000E+100
49 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -3.875757E-14 0.000000E+00 9.00000E+100
50 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -8.808540E-14 0.000000E+00 9.00000E+100
51 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 5.904019E-14 0.000000E+00 9.00000E+100
52 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 3.220374E-14 0.000000E+00 9.00000E+100
53 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 2.683645E-15 0.000000E+00 9.00000E+100
54 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.922519E-14 0.000000E+00 9.00000E+100
55 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 7.454132E-12 0.000000E+00 9.00000E+100
56 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 2.995984E-11 0.000000E+00 9.00000E+100
57 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 3.969929E-11 0.000000E+00 9.00000E+100
58 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 2.360887E-11 0.000000E+00 9.00000E+100
59 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -2.061273E-13 0.000000E+00 9.00000E+100
60 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 6.084423E-13 0.000000E+00 9.00000E+100
61 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 3.193739E-12 0.000000E+00 9.00000E+100
62 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 6.434103E-13 0.000000E+00 9.00000E+100
63 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 6.410931E-13 0.000000E+00 9.00000E+100
64 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 8.413726E-12 0.000000E+00 9.00000E+100
65 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 3.986953E-11 0.000000E+00 9.00000E+100
66 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 4.058096E-11 0.000000E+00 9.00000E+100
67 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 8.660586E-12 0.000000E+00 9.00000E+100
68 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 9.720963E-12 0.000000E+00 9.00000E+100
69 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.362944E-11 0.000000E+00 9.00000E+100
70 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 7.206085E-14 0.000000E+00 9.00000E+100
71 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 3.107562E-11 0.000000E+00 9.00000E+100
72 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 3.704738E-11 0.000000E+00 9.00000E+100
73 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.145550E-11 0.000000E+00 9.00000E+100
74 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.045309E-11 0.000000E+00 9.00000E+100
75 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.674004E-11 0.000000E+00 9.00000E+100
76 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 7.560085E-12 0.000000E+00 9.00000E+100
77 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 4.207913E-12 0.000000E+00 9.00000E+100
78 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 8.840094E-12 0.000000E+00 9.00000E+100
79 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.067173E-11 0.000000E+00 9.00000E+100
80 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 2.171573E-12 0.000000E+00 9.00000E+100
81 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 4.783597E-13 0.000000E+00 9.00000E+100
82 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 8.386390E-15 0.000000E+00 9.00000E+100
83 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -1.341822E-14 0.000000E+00 9.00000E+100
84 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 4.404270E-15 0.000000E+00 9.00000E+100
85 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -4.404270E-15 0.000000E+00 9.00000E+100
86 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -3.523416E-15 0.000000E+00 9.00000E+100
87 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 -1.904722E-15 0.000000E+00 9.00000E+100
88 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 4.190388E-15 0.000000E+00 9.00000E+100
89 traj.phases.phase0.collocation_constraint.defects:mass_fuel e 0.000000E+00 1.904722E-15 0.000000E+00 9.00000E+100
90 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -5.952256E-13 0.000000E+00 9.00000E+100
91 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -2.678515E-12 0.000000E+00 9.00000E+100
92 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -5.952256E-13 0.000000E+00 9.00000E+100
93 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 6.193504E-12 0.000000E+00 9.00000E+100
94 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 6.881672E-12 0.000000E+00 9.00000E+100
95 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.376334E-12 0.000000E+00 9.00000E+100
96 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -3.669046E-11 0.000000E+00 9.00000E+100
97 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.677278E-11 0.000000E+00 9.00000E+100
98 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.467618E-11 0.000000E+00 9.00000E+100
99 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.297154E-11 0.000000E+00 9.00000E+100
100 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.433697E-11 0.000000E+00 9.00000E+100
101 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.707816E-11 0.000000E+00 9.00000E+100
102 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 2.816114E-11 0.000000E+00 9.00000E+100
103 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.743551E-11 0.000000E+00 9.00000E+100
104 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 2.484433E-11 0.000000E+00 9.00000E+100
105 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -2.258774E-11 0.000000E+00 9.00000E+100
106 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.209422E-12 0.000000E+00 9.00000E+100
107 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -2.014755E-11 0.000000E+00 9.00000E+100
108 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 2.238555E-11 0.000000E+00 9.00000E+100
109 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 9.494143E-12 0.000000E+00 9.00000E+100
110 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.027521E-11 0.000000E+00 9.00000E+100
111 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -3.836184E-11 0.000000E+00 9.00000E+100
112 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.425432E-11 0.000000E+00 9.00000E+100
113 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -2.136697E-11 0.000000E+00 9.00000E+100
114 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.665267E-11 0.000000E+00 9.00000E+100
115 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.826030E-11 0.000000E+00 9.00000E+100
116 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.320232E-11 0.000000E+00 9.00000E+100
117 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 5.975968E-13 0.000000E+00 9.00000E+100
118 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -6.915049E-12 0.000000E+00 9.00000E+100
119 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.141125E-11 0.000000E+00 9.00000E+100
120 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 4.015758E-11 0.000000E+00 9.00000E+100
121 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 1.016648E-12 0.000000E+00 9.00000E+100
122 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 2.236625E-12 0.000000E+00 9.00000E+100
123 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 4.272289E-12 0.000000E+00 9.00000E+100
124 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -6.827129E-12 0.000000E+00 9.00000E+100
125 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.160612E-11 0.000000E+00 9.00000E+100
126 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -6.813942E-12 0.000000E+00 9.00000E+100
127 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
128 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 2.096598E-12 0.000000E+00 9.00000E+100
129 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.582784E-11 0.000000E+00 9.00000E+100
130 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -1.376334E-12 0.000000E+00 9.00000E+100
131 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -4.129003E-12 0.000000E+00 9.00000E+100
132 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 7.737932E-12 0.000000E+00 9.00000E+100
133 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -5.952256E-13 0.000000E+00 9.00000E+100
134 traj.phases.phase0.collocation_constraint.defects:alt e 0.000000E+00 -2.678515E-12 0.000000E+00 9.00000E+100
135 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
136 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
137 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
138 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
139 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
140 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
141 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
142 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
143 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
144 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
145 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
146 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
147 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
148 traj.phases.phase0.continuity_comp.defect_states:range e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
149 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
150 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
151 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
152 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
153 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
154 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
155 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
156 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
157 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
158 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
159 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
160 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
161 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
162 traj.phases.phase0.continuity_comp.defect_states:mass_fuel e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
163 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
164 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
165 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
166 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
167 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
168 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
169 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
170 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
171 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
172 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
173 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
174 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
175 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
176 traj.phases.phase0.continuity_comp.defect_states:alt e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
177 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 5.052444E-11 0.000000E+00 9.00000E+100
178 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 2.929924E-11 0.000000E+00 9.00000E+100
179 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -1.001433E-11 0.000000E+00 9.00000E+100
180 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 8.253995E-12 0.000000E+00 9.00000E+100
181 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 1.564738E-13 0.000000E+00 9.00000E+100
182 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -1.564738E-13 0.000000E+00 9.00000E+100
183 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -3.911846E-14 0.000000E+00 9.00000E+100
184 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 7.823692E-14 0.000000E+00 9.00000E+100
185 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -3.129477E-13 0.000000E+00 9.00000E+100
186 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
187 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -2.503581E-12 0.000000E+00 9.00000E+100
188 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -3.004298E-11 0.000000E+00 9.00000E+100
189 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 -1.042787E-11 0.000000E+00 9.00000E+100
190 traj.phases.phase0.continuity_comp.defect_control_rates:climb_rate_rate e 0.000000E+00 4.756028E-11 0.000000E+00 9.00000E+100
191 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
192 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 -4.547474E-13 0.000000E+00 9.00000E+100
193 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 -2.273737E-13 0.000000E+00 9.00000E+100
194 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 3.552714E-14 0.000000E+00 9.00000E+100
195 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 8.881784E-16 0.000000E+00 9.00000E+100
196 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
197 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 -8.881784E-16 0.000000E+00 9.00000E+100
198 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
199 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
200 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
201 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 1.065814E-14 0.000000E+00 9.00000E+100
202 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 -2.273737E-13 0.000000E+00 9.00000E+100
203 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 0.000000E+00 0.000000E+00 9.00000E+100
204 traj.phases.phase0.continuity_comp.defect_controls:climb_rate e 0.000000E+00 -4.547474E-13 0.000000E+00 9.00000E+100
205 traj.phases.phase0->path_constraint->tau i 1.000000E-02 5.293064E-01 2.000000E+00 9.00000E+100
206 traj.phases.phase0->path_constraint->tau i 1.000000E-02 5.459977E-01 2.000000E+00 9.00000E+100
207 traj.phases.phase0->path_constraint->tau i 1.000000E-02 5.710745E-01 2.000000E+00 9.00000E+100
208 traj.phases.phase0->path_constraint->tau i 1.000000E-02 5.795573E-01 2.000000E+00 9.00000E+100
209 traj.phases.phase0->path_constraint->tau i 1.000000E-02 5.795573E-01 2.000000E+00 9.00000E+100
210 traj.phases.phase0->path_constraint->tau i 1.000000E-02 6.296025E-01 2.000000E+00 9.00000E+100
211 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.166094E-01 2.000000E+00 9.00000E+100
212 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.497818E-01 2.000000E+00 9.00000E+100
213 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.497818E-01 2.000000E+00 9.00000E+100
214 traj.phases.phase0->path_constraint->tau i 1.000000E-02 8.706197E-01 2.000000E+00 9.00000E+100
215 traj.phases.phase0->path_constraint->tau i 1.000000E-02 9.783659E-01 2.000000E+00 9.00000E+100
216 traj.phases.phase0->path_constraint->tau i 1.000000E-02 9.719173E-01 2.000000E+00 9.00000E+100
217 traj.phases.phase0->path_constraint->tau i 1.000000E-02 9.719173E-01 2.000000E+00 9.00000E+100
218 traj.phases.phase0->path_constraint->tau i 1.000000E-02 8.846864E-01 2.000000E+00 9.00000E+100
219 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.714098E-01 2.000000E+00 9.00000E+100
220 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.598734E-01 2.000000E+00 9.00000E+100
221 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.598734E-01 2.000000E+00 9.00000E+100
222 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.532675E-01 2.000000E+00 9.00000E+100
223 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.518367E-01 2.000000E+00 9.00000E+100
224 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.538609E-01 2.000000E+00 9.00000E+100
225 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.538609E-01 2.000000E+00 9.00000E+100
226 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.607274E-01 2.000000E+00 9.00000E+100
227 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.668862E-01 2.000000E+00 9.00000E+100
228 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.671535E-01 2.000000E+00 9.00000E+100
229 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.671535E-01 2.000000E+00 9.00000E+100
230 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.659281E-01 2.000000E+00 9.00000E+100
231 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.630010E-01 2.000000E+00 9.00000E+100
232 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.619392E-01 2.000000E+00 9.00000E+100
233 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.619392E-01 2.000000E+00 9.00000E+100
234 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.606804E-01 2.000000E+00 9.00000E+100
235 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.651014E-01 2.000000E+00 9.00000E+100
236 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.685947E-01 2.000000E+00 9.00000E+100
237 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.685947E-01 2.000000E+00 9.00000E+100
238 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.751753E-01 2.000000E+00 9.00000E+100
239 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.698307E-01 2.000000E+00 9.00000E+100
240 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.629945E-01 2.000000E+00 9.00000E+100
241 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.629945E-01 2.000000E+00 9.00000E+100
242 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.508298E-01 2.000000E+00 9.00000E+100
243 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.635160E-01 2.000000E+00 9.00000E+100
244 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.778113E-01 2.000000E+00 9.00000E+100
245 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.778113E-01 2.000000E+00 9.00000E+100
246 traj.phases.phase0->path_constraint->tau i 1.000000E-02 8.060096E-01 2.000000E+00 9.00000E+100
247 traj.phases.phase0->path_constraint->tau i 1.000000E-02 8.034679E-01 2.000000E+00 9.00000E+100
248 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.873301E-01 2.000000E+00 9.00000E+100
249 traj.phases.phase0->path_constraint->tau i 1.000000E-02 7.873301E-01 2.000000E+00 9.00000E+100
250 traj.phases.phase0->path_constraint->tau i 1.000000E-02 6.985485E-01 2.000000E+00 9.00000E+100
251 traj.phases.phase0->path_constraint->tau i 1.000000E-02 4.201175E-01 2.000000E+00 9.00000E+100
252 traj.phases.phase0->path_constraint->tau i 1.000000E-02 3.127165E-01 2.000000E+00 9.00000E+100
253 traj.phases.phase0->path_constraint->tau i 1.000000E-02 3.127165E-01 2.000000E+00 9.00000E+100
254 traj.phases.phase0->path_constraint->tau i 1.000000E-02 1.836616E-01 2.000000E+00 9.00000E+100
255 traj.phases.phase0->path_constraint->tau i 1.000000E-02 1.531111E-01 2.000000E+00 9.00000E+100
256 traj.phases.phase0->path_constraint->tau i 1.000000E-02 1.633905E-01 2.000000E+00 9.00000E+100
257 traj.phases.phase0->path_constraint->tau i 1.000000E-02 1.633905E-01 2.000000E+00 9.00000E+100
258 traj.phases.phase0->path_constraint->tau i 1.000000E-02 1.816829E-01 2.000000E+00 9.00000E+100
259 traj.phases.phase0->path_constraint->tau i 1.000000E-02 2.026853E-01 2.000000E+00 9.00000E+100
260 traj.phases.phase0->path_constraint->tau i 1.000000E-02 2.082759E-01 2.000000E+00 9.00000E+100
261 traj.phases.phase0->path_constraint->tau i 1.000000E-02 2.082759E-01 2.000000E+00 9.00000E+100
262 traj.phases.phase0->path_constraint->tau i 1.000000E-02 2.133263E-01 2.000000E+00 9.00000E+100
263 traj.phases.phase0->path_constraint->tau i 1.000000E-02 2.195263E-01 2.000000E+00 9.00000E+100
264 traj.phases.phase0->path_constraint->tau i 1.000000E-02 2.213099E-01 2.000000E+00 9.00000E+100
Exit Status
Inform Description
0 Solve Succeeded
--------------------------------------------------------------------------------
/usr/share/miniconda/envs/test/lib/python3.10/site-packages/openmdao/visualization/opt_report/opt_report.py:634: UserWarning: Attempting to set identical low and high ylims makes transformation singular; automatically expanding.
ax.set_ylim([ymin_plot, ymax_plot])
Simulating trajectory traj
Model viewer data has already been recorded for Driver.
Done simulating trajectory traj
False
from IPython.display import Image
Image(p.get_reports_dir() / 'plots/states_alt.png', width=600)
---------------------------------------------------------------------------
TypeError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1045, in Image._data_and_metadata(self, always_both)
1044 try:
-> 1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
TypeError: a bytes-like object is required, not 'str'
The above exception was the direct cause of the following exception:
FileNotFoundError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/formatters.py:974, in MimeBundleFormatter.__call__(self, obj, include, exclude)
971 method = get_real_method(obj, self.print_method)
973 if method is not None:
--> 974 return method(include=include, exclude=exclude)
975 return None
976 else:
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1035, in Image._repr_mimebundle_(self, include, exclude)
1033 if self.embed:
1034 mimetype = self._mimetype
-> 1035 data, metadata = self._data_and_metadata(always_both=True)
1036 if metadata:
1037 metadata = {mimetype: metadata}
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1047, in Image._data_and_metadata(self, always_both)
1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
-> 1047 raise FileNotFoundError(
1048 "No such file or directory: '%s'" % (self.data)) from e
1049 md = {}
1050 if self.metadata:
FileNotFoundError: No such file or directory: 'reports/problem/plots/states_alt.png'
---------------------------------------------------------------------------
TypeError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1045, in Image._data_and_metadata(self, always_both)
1044 try:
-> 1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
TypeError: a bytes-like object is required, not 'str'
The above exception was the direct cause of the following exception:
FileNotFoundError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/formatters.py:344, in BaseFormatter.__call__(self, obj)
342 method = get_real_method(obj, self.print_method)
343 if method is not None:
--> 344 return method()
345 return None
346 else:
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1067, in Image._repr_png_(self)
1065 def _repr_png_(self):
1066 if self.embed and self.format == self._FMT_PNG:
-> 1067 return self._data_and_metadata()
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1047, in Image._data_and_metadata(self, always_both)
1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
-> 1047 raise FileNotFoundError(
1048 "No such file or directory: '%s'" % (self.data)) from e
1049 md = {}
1050 if self.metadata:
FileNotFoundError: No such file or directory: 'reports/problem/plots/states_alt.png'
<IPython.core.display.Image object>
Image(p.get_reports_dir() / 'plots/states_mass_fuel.png', width=600)
---------------------------------------------------------------------------
TypeError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1045, in Image._data_and_metadata(self, always_both)
1044 try:
-> 1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
TypeError: a bytes-like object is required, not 'str'
The above exception was the direct cause of the following exception:
FileNotFoundError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/formatters.py:974, in MimeBundleFormatter.__call__(self, obj, include, exclude)
971 method = get_real_method(obj, self.print_method)
973 if method is not None:
--> 974 return method(include=include, exclude=exclude)
975 return None
976 else:
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1035, in Image._repr_mimebundle_(self, include, exclude)
1033 if self.embed:
1034 mimetype = self._mimetype
-> 1035 data, metadata = self._data_and_metadata(always_both=True)
1036 if metadata:
1037 metadata = {mimetype: metadata}
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1047, in Image._data_and_metadata(self, always_both)
1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
-> 1047 raise FileNotFoundError(
1048 "No such file or directory: '%s'" % (self.data)) from e
1049 md = {}
1050 if self.metadata:
FileNotFoundError: No such file or directory: 'reports/problem/plots/states_mass_fuel.png'
---------------------------------------------------------------------------
TypeError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1045, in Image._data_and_metadata(self, always_both)
1044 try:
-> 1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
TypeError: a bytes-like object is required, not 'str'
The above exception was the direct cause of the following exception:
FileNotFoundError Traceback (most recent call last)
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/formatters.py:344, in BaseFormatter.__call__(self, obj)
342 method = get_real_method(obj, self.print_method)
343 if method is not None:
--> 344 return method()
345 return None
346 else:
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1067, in Image._repr_png_(self)
1065 def _repr_png_(self):
1066 if self.embed and self.format == self._FMT_PNG:
-> 1067 return self._data_and_metadata()
File /usr/share/miniconda/envs/test/lib/python3.10/site-packages/IPython/core/display.py:1047, in Image._data_and_metadata(self, always_both)
1045 b64_data = b2a_base64(self.data, newline=False).decode("ascii")
1046 except TypeError as e:
-> 1047 raise FileNotFoundError(
1048 "No such file or directory: '%s'" % (self.data)) from e
1049 md = {}
1050 if self.metadata:
FileNotFoundError: No such file or directory: 'reports/problem/plots/states_mass_fuel.png'
<IPython.core.display.Image object>
References#
John T Betts and Evin J Cramer. Application of direct transcription to commercial aircraft trajectory optimization. Journal of Guidance, Control, and Dynamics, 18(1):151–159, 1995.
Fariba Fahroo and I Michael Ross. Second look at approximating differential inclusions. Journal of Guidance, Control, and Dynamics, 24(1):131–133, 2001.
Renjith R Kumar and Hans Seywald. Should controls be eliminated while solving optimal control problems via direct methods? Journal of Guidance, Control, and Dynamics, 19(2):418–423, mar 1996. URL: https://doi.org/10.2514/3.21634, doi:10.2514/3.21634.
Hans Seywald. Trajectory optimization based on differential inclusion (Revised). Journal of Guidance, Control, and Dynamics, 17(3):480–487, may 1994. URL: https://doi.org/10.2514/3.21224, doi:10.2514/3.21224.