User-Provided Mesh Example
Here is an example script with a custom mesh provided as an array of coordinates. This should help you understand how meshes are defined in OpenAeroStruct and how to create them for your own custom planform shapes. This is an alternative to the helper-function approach described in Aerodynamic Optimization.
The following shows the portion of the example script in which the user provides the coordinates for the mesh. This example is for a wing with a kink and two distinct trapezoidal segments.
# -----------------------------------------------------------------------------
# CUSTOM MESH: Example mesh for a 2-segment wing with sweep
# -----------------------------------------------------------------------------
# Planform specifications
half_span = 12.0 # wing half-span in m
kink_location = 4.0 # spanwise location of the kink in m
root_chord = 6.0 # root chord in m
kink_chord = 3.0 # kink chord in m
tip_chord = 2.0 # tip chord in m
inboard_LE_sweep = 10.0 # inboard leading-edge sweep angle in deg
outboard_LE_sweep = -10.0 # outboard leading-edge sweep angle in deg
# Mesh specifications
nx = 5 # number of chordwise nodal points (should be odd)
ny_outboard = 9 # number of spanwise nodal points for the outboard segment
ny_inboard = 7 # number of spanwise nodal points for the inboard segment
# Initialize the 3-D mesh object. Indexing: Chordwise, spanwise, then the 3-D coordinates.
# We use ny_inboard+ny_outboard-1 because the 2 segments share the nodes where they connect.
mesh = np.zeros((nx, ny_inboard + ny_outboard - 1, 3))
# The form of this 3-D array can be confusing initially.
# For each node, we are providing the x, y, and z coordinates.
# x is streamwise, y is spanwise, and z is up.
# For example, the node for the leading edge at the tip would be specified as mesh[0, 0, :] = np.array([x, y, z]).
# And the node at the trailing edge at the root would be mesh[nx-1, ny-1, :] = np.array([x, y, z]).
# We only provide the right half of the wing here because we use symmetry.
# Print elements of the mesh to better understand the form.
####### THE Z-COORDINATES ######
# Assume no dihedral, so set the z-coordinate for all the points to 0.
mesh[:, :, 2] = 0.0
####### THE Y-COORDINATES ######
# Using uniform spacing for the spanwise locations of all the nodes within each of the two trapezoidal segments:
# Outboard
mesh[:, :ny_outboard, 1] = np.linspace(half_span, kink_location, ny_outboard)
# Inboard
mesh[:, ny_outboard : ny_outboard + ny_inboard, 1] = np.linspace(kink_location, 0, ny_inboard)[1:]
###### THE X-COORDINATES ######
# Start with the leading edge and create some intermediate arrays that we will use
x_LE = np.zeros(ny_inboard + ny_outboard - 1)
array_for_inboard_leading_edge_x_coord = np.linspace(0, kink_location, ny_inboard) * np.tan(
inboard_LE_sweep / 180.0 * np.pi
)
array_for_outboard_leading_edge_x_coord = (
np.linspace(0, half_span - kink_location, ny_outboard) * np.tan(outboard_LE_sweep / 180.0 * np.pi)
+ np.ones(ny_outboard) * array_for_inboard_leading_edge_x_coord[-1]
)
x_LE[:ny_inboard] = array_for_inboard_leading_edge_x_coord
x_LE[ny_inboard : ny_inboard + ny_outboard] = array_for_outboard_leading_edge_x_coord[1:]
# Then the trailing edge
x_TE = np.zeros(ny_inboard + ny_outboard - 1)
array_for_inboard_trailing_edge_x_coord = np.linspace(
array_for_inboard_leading_edge_x_coord[0] + root_chord,
array_for_inboard_leading_edge_x_coord[-1] + kink_chord,
ny_inboard,
)
array_for_outboard_trailing_edge_x_coord = np.linspace(
array_for_outboard_leading_edge_x_coord[0] + kink_chord,
array_for_outboard_leading_edge_x_coord[-1] + tip_chord,
ny_outboard,
)
x_TE[:ny_inboard] = array_for_inboard_trailing_edge_x_coord
x_TE[ny_inboard : ny_inboard + ny_outboard] = array_for_outboard_trailing_edge_x_coord[1:]
# # Quick plot to check leading and trailing edge x-coords
# plt.plot(x_LE, np.arange(0, ny_inboard+ny_outboard-1), marker='*')
# plt.plot(x_TE, np.arange(0, ny_inboard+ny_outboard-1), marker='*')
# plt.show()
# exit()
for i in range(0, ny_inboard + ny_outboard - 1):
mesh[:, i, 0] = np.linspace(np.flip(x_LE)[i], np.flip(x_TE)[i], nx)
# -----------------------------------------------------------------------------
# END MESH
# -----------------------------------------------------------------------------
The following shows a visualization of the mesh.
The complete script for the optimization is as follows. Make sure you go through the Aerostructural Optimization with Wingbox before trying to understand this setup.
import warnings
import matplotlib
warnings.filterwarnings('ignore')
matplotlib.use('Agg')
"""
This example script can be used to run a multipoint aerostructural (w/ wingbox) optimization for a custom user-provided mesh.
The fuel burn from the cruise flight-point is the objective function and a 2.5g
maneuver flight-point is used for the structural sizing.
After running the optimization, use the 'plot_wingbox.py' script in the utils/
directory (e.g., as 'python ../utils/plot_wingbox.py aerostruct.db' if running
from this directory) to visualize the results.
'plot_wingbox.py' is still under development and will probably not work as it is for other types of cases for now.
"""
import numpy as np
from openaerostruct.integration.aerostruct_groups import AerostructGeometry, AerostructPoint
from openaerostruct.structures.wingbox_fuel_vol_delta import WingboxFuelVolDelta
import openmdao.api as om
# docs checkpoint 0
# -----------------------------------------------------------------------------
# CUSTOM MESH: Example mesh for a 2-segment wing with sweep
# -----------------------------------------------------------------------------
# Planform specifications
half_span = 12.0 # wing half-span in m
kink_location = 4.0 # spanwise location of the kink in m
root_chord = 6.0 # root chord in m
kink_chord = 3.0 # kink chord in m
tip_chord = 2.0 # tip chord in m
inboard_LE_sweep = 10.0 # inboard leading-edge sweep angle in deg
outboard_LE_sweep = -10.0 # outboard leading-edge sweep angle in deg
# Mesh specifications
nx = 5 # number of chordwise nodal points (should be odd)
ny_outboard = 9 # number of spanwise nodal points for the outboard segment
ny_inboard = 7 # number of spanwise nodal points for the inboard segment
# Initialize the 3-D mesh object. Indexing: Chordwise, spanwise, then the 3-D coordinates.
# We use ny_inboard+ny_outboard-1 because the 2 segments share the nodes where they connect.
mesh = np.zeros((nx, ny_inboard + ny_outboard - 1, 3))
# The form of this 3-D array can be confusing initially.
# For each node, we are providing the x, y, and z coordinates.
# x is streamwise, y is spanwise, and z is up.
# For example, the node for the leading edge at the tip would be specified as mesh[0, 0, :] = np.array([x, y, z]).
# And the node at the trailing edge at the root would be mesh[nx-1, ny-1, :] = np.array([x, y, z]).
# We only provide the right half of the wing here because we use symmetry.
# Print elements of the mesh to better understand the form.
####### THE Z-COORDINATES ######
# Assume no dihedral, so set the z-coordinate for all the points to 0.
mesh[:, :, 2] = 0.0
####### THE Y-COORDINATES ######
# Using uniform spacing for the spanwise locations of all the nodes within each of the two trapezoidal segments:
# Outboard
mesh[:, :ny_outboard, 1] = np.linspace(half_span, kink_location, ny_outboard)
# Inboard
mesh[:, ny_outboard : ny_outboard + ny_inboard, 1] = np.linspace(kink_location, 0, ny_inboard)[1:]
###### THE X-COORDINATES ######
# Start with the leading edge and create some intermediate arrays that we will use
x_LE = np.zeros(ny_inboard + ny_outboard - 1)
array_for_inboard_leading_edge_x_coord = np.linspace(0, kink_location, ny_inboard) * np.tan(
inboard_LE_sweep / 180.0 * np.pi
)
array_for_outboard_leading_edge_x_coord = (
np.linspace(0, half_span - kink_location, ny_outboard) * np.tan(outboard_LE_sweep / 180.0 * np.pi)
+ np.ones(ny_outboard) * array_for_inboard_leading_edge_x_coord[-1]
)
x_LE[:ny_inboard] = array_for_inboard_leading_edge_x_coord
x_LE[ny_inboard : ny_inboard + ny_outboard] = array_for_outboard_leading_edge_x_coord[1:]
# Then the trailing edge
x_TE = np.zeros(ny_inboard + ny_outboard - 1)
array_for_inboard_trailing_edge_x_coord = np.linspace(
array_for_inboard_leading_edge_x_coord[0] + root_chord,
array_for_inboard_leading_edge_x_coord[-1] + kink_chord,
ny_inboard,
)
array_for_outboard_trailing_edge_x_coord = np.linspace(
array_for_outboard_leading_edge_x_coord[0] + kink_chord,
array_for_outboard_leading_edge_x_coord[-1] + tip_chord,
ny_outboard,
)
x_TE[:ny_inboard] = array_for_inboard_trailing_edge_x_coord
x_TE[ny_inboard : ny_inboard + ny_outboard] = array_for_outboard_trailing_edge_x_coord[1:]
# # Quick plot to check leading and trailing edge x-coords
# plt.plot(x_LE, np.arange(0, ny_inboard+ny_outboard-1), marker='*')
# plt.plot(x_TE, np.arange(0, ny_inboard+ny_outboard-1), marker='*')
# plt.show()
# exit()
for i in range(0, ny_inboard + ny_outboard - 1):
mesh[:, i, 0] = np.linspace(np.flip(x_LE)[i], np.flip(x_TE)[i], nx)
# -----------------------------------------------------------------------------
# END MESH
# -----------------------------------------------------------------------------
# docs checkpoint 1
# -----------------------------------------------------------------------------
# On to the problem setup (this is the same setup used for the Q400 example)
# -----------------------------------------------------------------------------
# Provide coordinates for a portion of an airfoil for the wingbox cross-section as an nparray with dtype=complex (to work with the complex-step approximation for derivatives).
# These should be for an airfoil with the chord scaled to 1.
# We use the 10% to 60% portion of the NASA SC2-0612 airfoil for this case
# We use the coordinates available from airfoiltools.com. Using such a large number of coordinates is not necessary.
# The first and last x-coordinates of the upper and lower surfaces must be the same
# fmt: off
upper_x = np.array([0.1, 0.11, 0.12, 0.13, 0.14, 0.15, 0.16, 0.17, 0.18, 0.19, 0.2, 0.21, 0.22, 0.23, 0.24, 0.25, 0.26, 0.27, 0.28, 0.29, 0.3, 0.31, 0.32, 0.33, 0.34, 0.35, 0.36, 0.37, 0.38, 0.39, 0.4, 0.41, 0.42, 0.43, 0.44, 0.45, 0.46, 0.47, 0.48, 0.49, 0.5, 0.51, 0.52, 0.53, 0.54, 0.55, 0.56, 0.57, 0.58, 0.59, 0.6], dtype="complex128")
lower_x = np.array([0.1, 0.11, 0.12, 0.13, 0.14, 0.15, 0.16, 0.17, 0.18, 0.19, 0.2, 0.21, 0.22, 0.23, 0.24, 0.25, 0.26, 0.27, 0.28, 0.29, 0.3, 0.31, 0.32, 0.33, 0.34, 0.35, 0.36, 0.37, 0.38, 0.39, 0.4, 0.41, 0.42, 0.43, 0.44, 0.45, 0.46, 0.47, 0.48, 0.49, 0.5, 0.51, 0.52, 0.53, 0.54, 0.55, 0.56, 0.57, 0.58, 0.59, 0.6], dtype="complex128")
upper_y = np.array([ 0.0447, 0.046, 0.0472, 0.0484, 0.0495, 0.0505, 0.0514, 0.0523, 0.0531, 0.0538, 0.0545, 0.0551, 0.0557, 0.0563, 0.0568, 0.0573, 0.0577, 0.0581, 0.0585, 0.0588, 0.0591, 0.0593, 0.0595, 0.0597, 0.0599, 0.06, 0.0601, 0.0602, 0.0602, 0.0602, 0.0602, 0.0602, 0.0601, 0.06, 0.0599, 0.0598, 0.0596, 0.0594, 0.0592, 0.0589, 0.0586, 0.0583, 0.058, 0.0576, 0.0572, 0.0568, 0.0563, 0.0558, 0.0553, 0.0547, 0.0541], dtype="complex128") # noqa: E201, E241
lower_y = np.array([-0.0447, -0.046, -0.0473, -0.0485, -0.0496, -0.0506, -0.0515, -0.0524, -0.0532, -0.054, -0.0547, -0.0554, -0.056, -0.0565, -0.057, -0.0575, -0.0579, -0.0583, -0.0586, -0.0589, -0.0592, -0.0594, -0.0595, -0.0596, -0.0597, -0.0598, -0.0598, -0.0598, -0.0598, -0.0597, -0.0596, -0.0594, -0.0592, -0.0589, -0.0586, -0.0582, -0.0578, -0.0573, -0.0567, -0.0561, -0.0554, -0.0546, -0.0538, -0.0529, -0.0519, -0.0509, -0.0497, -0.0485, -0.0472, -0.0458, -0.0444], dtype="complex128")
# fmt: on
surf_dict = {
# Wing definition
"name": "wing", # name of the surface
"symmetry": True, # if true, model one half of wing
"S_ref_type": "wetted", # how we compute the wing area,
# can be 'wetted' or 'projected'
"mesh": mesh,
"twist_cp": np.array([6.0, 7.0, 7.0, 7.0]),
"fem_model_type": "wingbox",
"data_x_upper": upper_x,
"data_x_lower": lower_x,
"data_y_upper": upper_y,
"data_y_lower": lower_y,
"spar_thickness_cp": np.array([0.004, 0.004, 0.004, 0.004]), # [m]
"skin_thickness_cp": np.array([0.003, 0.006, 0.010, 0.012]), # [m]
"original_wingbox_airfoil_t_over_c": 0.12,
# Aerodynamic deltas.
# These CL0 and CD0 values are added to the CL and CD
# obtained from aerodynamic analysis of the surface to get
# the total CL and CD.
# These CL0 and CD0 values do not vary wrt alpha.
# They can be used to account for things that are not included, such as contributions from the fuselage, nacelles, tail surfaces, etc.
"CL0": 0.0,
"CD0": 0.0142,
"with_viscous": True, # if true, compute viscous drag
"with_wave": True, # if true, compute wave drag
# Airfoil properties for viscous drag calculation
"k_lam": 0.05, # percentage of chord with laminar
# flow, used for viscous drag
"c_max_t": 0.38, # chordwise location of maximum thickness
"t_over_c_cp": np.array([0.1, 0.1, 0.15, 0.15]),
# Structural values are based on aluminum 7075
"E": 73.1e9, # [Pa] Young's modulus
"G": (73.1e9 / 2 / 1.33), # [Pa] shear modulus (calculated using E and the Poisson's ratio here)
"yield": (420.0e6 / 1.5), # [Pa] allowable yield stress
"mrho": 2.78e3, # [kg/m^3] material density
"strength_factor_for_upper_skin": 1.0, # the yield stress is multiplied by this factor for the upper skin
"wing_weight_ratio": 1.25,
"exact_failure_constraint": False, # if false, use KS function
"struct_weight_relief": True,
"distributed_fuel_weight": True,
"fuel_density": 803.0, # [kg/m^3] fuel density (only needed if the fuel-in-wing volume constraint is used)
"Wf_reserve": 500.0, # [kg] reserve fuel mass
}
surfaces = [surf_dict]
# Create the problem and assign the model group
prob = om.Problem()
# Add problem information as an independent variables component
indep_var_comp = om.IndepVarComp()
indep_var_comp.add_output("v", val=np.array([0.5 * 310.95, 0.3 * 340.294]), units="m/s")
indep_var_comp.add_output("alpha", val=0.0, units="deg")
indep_var_comp.add_output("alpha_maneuver", val=0.0, units="deg")
indep_var_comp.add_output("Mach_number", val=np.array([0.5, 0.3]))
indep_var_comp.add_output(
"re",
val=np.array([0.569 * 310.95 * 0.5 * 1.0 / (1.56 * 1e-5), 1.225 * 340.294 * 0.3 * 1.0 / (1.81206 * 1e-5)]),
units="1/m",
)
indep_var_comp.add_output("rho", val=np.array([0.569, 1.225]), units="kg/m**3")
indep_var_comp.add_output("CT", val=0.43 / 3600, units="1/s")
indep_var_comp.add_output("R", val=2e6, units="m")
indep_var_comp.add_output("W0", val=25400 + surf_dict["Wf_reserve"], units="kg")
indep_var_comp.add_output("speed_of_sound", val=np.array([310.95, 340.294]), units="m/s")
indep_var_comp.add_output("load_factor", val=np.array([1.0, 2.5]))
indep_var_comp.add_output("empty_cg", val=np.zeros((3)), units="m")
indep_var_comp.add_output("fuel_mass", val=3000.0, units="kg")
prob.model.add_subsystem("prob_vars", indep_var_comp, promotes=["*"])
# Loop over each surface in the surfaces list
for surface in surfaces:
# Get the surface name and create a group to contain components
# only for this surface
name = surface["name"]
aerostruct_group = AerostructGeometry(surface=surface)
# Add group to the problem with the name of the surface.
prob.model.add_subsystem(name, aerostruct_group)
# Loop through and add a certain number of aerostruct points
for i in range(2):
point_name = "AS_point_{}".format(i)
# Connect the parameters within the model for each aerostruct point
# Create the aero point group and add it to the model
AS_point = AerostructPoint(surfaces=surfaces, internally_connect_fuelburn=False)
prob.model.add_subsystem(point_name, AS_point)
# Connect flow properties to the analysis point
prob.model.connect("v", point_name + ".v", src_indices=[i])
prob.model.connect("Mach_number", point_name + ".Mach_number", src_indices=[i])
prob.model.connect("re", point_name + ".re", src_indices=[i])
prob.model.connect("rho", point_name + ".rho", src_indices=[i])
prob.model.connect("CT", point_name + ".CT")
prob.model.connect("R", point_name + ".R")
prob.model.connect("W0", point_name + ".W0")
prob.model.connect("speed_of_sound", point_name + ".speed_of_sound", src_indices=[i])
prob.model.connect("empty_cg", point_name + ".empty_cg")
prob.model.connect("load_factor", point_name + ".load_factor", src_indices=[i])
prob.model.connect("fuel_mass", point_name + ".total_perf.L_equals_W.fuelburn")
prob.model.connect("fuel_mass", point_name + ".total_perf.CG.fuelburn")
for surface in surfaces:
name = surface["name"]
if surf_dict["distributed_fuel_weight"]:
prob.model.connect("load_factor", point_name + ".coupled.load_factor", src_indices=[i])
com_name = point_name + "." + name + "_perf."
prob.model.connect(
name + ".local_stiff_transformed", point_name + ".coupled." + name + ".local_stiff_transformed"
)
prob.model.connect(name + ".nodes", point_name + ".coupled." + name + ".nodes")
# Connect aerodynamic mesh to coupled group mesh
prob.model.connect(name + ".mesh", point_name + ".coupled." + name + ".mesh")
if surf_dict["struct_weight_relief"]:
prob.model.connect(name + ".element_mass", point_name + ".coupled." + name + ".element_mass")
# Connect performance calculation variables
prob.model.connect(name + ".nodes", com_name + "nodes")
prob.model.connect(name + ".cg_location", point_name + "." + "total_perf." + name + "_cg_location")
prob.model.connect(name + ".structural_mass", point_name + "." + "total_perf." + name + "_structural_mass")
# Connect wingbox properties to von Mises stress calcs
prob.model.connect(name + ".Qz", com_name + "Qz")
prob.model.connect(name + ".J", com_name + "J")
prob.model.connect(name + ".A_enc", com_name + "A_enc")
prob.model.connect(name + ".htop", com_name + "htop")
prob.model.connect(name + ".hbottom", com_name + "hbottom")
prob.model.connect(name + ".hfront", com_name + "hfront")
prob.model.connect(name + ".hrear", com_name + "hrear")
prob.model.connect(name + ".spar_thickness", com_name + "spar_thickness")
prob.model.connect(name + ".t_over_c", com_name + "t_over_c")
prob.model.connect("alpha", "AS_point_0" + ".alpha")
prob.model.connect("alpha_maneuver", "AS_point_1" + ".alpha")
# Here we add the fuel volume constraint componenet to the model
prob.model.add_subsystem("fuel_vol_delta", WingboxFuelVolDelta(surface=surface))
prob.model.connect("wing.struct_setup.fuel_vols", "fuel_vol_delta.fuel_vols")
prob.model.connect("AS_point_0.fuelburn", "fuel_vol_delta.fuelburn")
if surf_dict["distributed_fuel_weight"]:
prob.model.connect("wing.struct_setup.fuel_vols", "AS_point_0.coupled.wing.struct_states.fuel_vols")
prob.model.connect("fuel_mass", "AS_point_0.coupled.wing.struct_states.fuel_mass")
prob.model.connect("wing.struct_setup.fuel_vols", "AS_point_1.coupled.wing.struct_states.fuel_vols")
prob.model.connect("fuel_mass", "AS_point_1.coupled.wing.struct_states.fuel_mass")
comp = om.ExecComp("fuel_diff = (fuel_mass - fuelburn) / fuelburn", units="kg")
prob.model.add_subsystem("fuel_diff", comp, promotes_inputs=["fuel_mass"], promotes_outputs=["fuel_diff"])
prob.model.connect("AS_point_0.fuelburn", "fuel_diff.fuelburn")
## Use these settings if you do not have pyOptSparse or SNOPT
prob.driver = om.ScipyOptimizeDriver()
prob.driver.options["optimizer"] = "SLSQP"
prob.driver.options["tol"] = 1e-4
recorder = om.SqliteRecorder("aerostruct.db")
prob.driver.add_recorder(recorder)
# We could also just use prob.driver.recording_options['includes']=['*'] here, but for large meshes the database file becomes extremely large. So we just select the variables we need.
prob.driver.recording_options["includes"] = [
"alpha",
"rho",
"v",
"cg",
"AS_point_1.cg",
"AS_point_0.cg",
"AS_point_0.coupled.wing_loads.loads",
"AS_point_1.coupled.wing_loads.loads",
"AS_point_0.coupled.wing.normals",
"AS_point_1.coupled.wing.normals",
"AS_point_0.coupled.wing.widths",
"AS_point_1.coupled.wing.widths",
"AS_point_0.coupled.aero_states.wing_sec_forces",
"AS_point_1.coupled.aero_states.wing_sec_forces",
"AS_point_0.wing_perf.CL1",
"AS_point_1.wing_perf.CL1",
"AS_point_0.coupled.wing.S_ref",
"AS_point_1.coupled.wing.S_ref",
"wing.geometry.twist",
"wing.mesh",
"wing.skin_thickness",
"wing.spar_thickness",
"wing.t_over_c",
"wing.structural_mass",
"AS_point_0.wing_perf.vonmises",
"AS_point_1.wing_perf.vonmises",
"AS_point_0.coupled.wing.def_mesh",
"AS_point_1.coupled.wing.def_mesh",
]
prob.driver.recording_options["record_objectives"] = True
prob.driver.recording_options["record_constraints"] = True
prob.driver.recording_options["record_desvars"] = True
prob.driver.recording_options["record_inputs"] = True
prob.model.add_objective("AS_point_0.fuelburn", scaler=1e-5)
prob.model.add_design_var("wing.twist_cp", lower=-15.0, upper=15.0, scaler=0.1)
prob.model.add_design_var("wing.spar_thickness_cp", lower=0.003, upper=0.1, scaler=1e2)
prob.model.add_design_var("wing.skin_thickness_cp", lower=0.003, upper=0.1, scaler=1e2)
prob.model.add_design_var("wing.geometry.t_over_c_cp", lower=0.07, upper=0.2, scaler=10.0)
prob.model.add_design_var("fuel_mass", lower=0.0, upper=2e5, scaler=1e-5)
prob.model.add_design_var("alpha_maneuver", lower=-15.0, upper=15)
prob.model.add_constraint("AS_point_0.CL", equals=0.6)
prob.model.add_constraint("AS_point_1.L_equals_W", equals=0.0)
prob.model.add_constraint("AS_point_1.wing_perf.failure", upper=0.0)
prob.model.add_constraint("fuel_vol_delta.fuel_vol_delta", lower=0.0)
prob.model.add_constraint("fuel_diff", equals=0.0)
# Set up the problem
prob.setup()
/home/docs/checkouts/readthedocs.org/user_builds/mdolab-openaerostruct/envs/latest/lib/python3.11/site-packages/openmdao/core/system.py:2415: PromotionWarning:'AS_point_0.coupled.wing.struct_states' <class SpatialBeamStates>: input variable 'nodes', promoted using 'nodes', was already promoted using 'nodes'. /home/docs/checkouts/readthedocs.org/user_builds/mdolab-openaerostruct/envs/latest/lib/python3.11/site-packages/openmdao/core/system.py:2415: PromotionWarning:'AS_point_0.coupled.wing.struct_states' <class SpatialBeamStates>: input variable 'load_factor', promoted using 'load_factor', was already promoted using 'load_factor'. /home/docs/checkouts/readthedocs.org/user_builds/mdolab-openaerostruct/envs/latest/lib/python3.11/site-packages/openmdao/core/system.py:2415: PromotionWarning:'AS_point_1.coupled.wing.struct_states' <class SpatialBeamStates>: input variable 'nodes', promoted using 'nodes', was already promoted using 'nodes'. /home/docs/checkouts/readthedocs.org/user_builds/mdolab-openaerostruct/envs/latest/lib/python3.11/site-packages/openmdao/core/system.py:2415: PromotionWarning:'AS_point_1.coupled.wing.struct_states' <class SpatialBeamStates>: input variable 'load_factor', promoted using 'load_factor', was already promoted using 'load_factor'. /home/docs/checkouts/readthedocs.org/user_builds/mdolab-openaerostruct/envs/latest/lib/python3.11/site-packages/openmdao/core/driver.py:550: OpenMDAOWarning:ScipyOptimizeDriver: No matches for pattern 'cg' in recording_options['includes'].
# change linear solver for aerostructural coupled adjoint
prob.model.AS_point_0.coupled.linear_solver = om.LinearBlockGS(iprint=0, maxiter=30, use_aitken=True)
prob.model.AS_point_1.coupled.linear_solver = om.LinearBlockGS(iprint=0, maxiter=30, use_aitken=True)
# om.view_model(prob)
# prob.check_partials(form='central', compact_print=True)
prob.run_driver()
================== AS_point_0.coupled ================== NL: NLBGS 1 ; 59474.4786 1 NL: NLBGS 2 ; 58887.5619 0.990131622 NL: NLBGS 3 ; 2197.7072 0.0369521054 NL: NLBGS 4 ; 67.9021348 0.00114170206 NL: NLBGS 5 ; 0.185879282 3.12536212e-06 NL: NLBGS 6 ; 0.00905010794 1.52167924e-07 NL: NLBGS 7 ; 7.21557797e-05 1.21322257e-09 NL: NLBGS 8 ; 1.5466389e-07 2.60050855e-12 NL: NLBGS 9 ; 2.80240759e-09 4.71194984e-14 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 63459.5261 1 NL: NLBGS 2 ; 54460.64 0.858194874 NL: NLBGS 3 ; 1868.49116 0.0294438246 NL: NLBGS 4 ; 53.4050567 0.000841560913 NL: NLBGS 5 ; 0.136184357 2.14600337e-06 NL: NLBGS 6 ; 0.00600916274 9.46928399e-08 NL: NLBGS 7 ; 5.31043825e-05 8.36822867e-10 NL: NLBGS 8 ; 8.3244288e-08 1.31176977e-12 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 1.55051266e-10 1 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 1.85127692e-09 1 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 12750.5196 1 NL: NLBGS 2 ; 11862.3875 0.930345417 NL: NLBGS 3 ; 1908.38371 0.149671054 NL: NLBGS 4 ; 88.3895614 0.0069322321 NL: NLBGS 5 ; 0.962907613 7.5519088e-05 NL: NLBGS 6 ; 0.0841114197 6.59670525e-06 NL: NLBGS 7 ; 0.0023381442 1.83376386e-07 NL: NLBGS 8 ; 2.70106133e-05 2.11839314e-09 NL: NLBGS 9 ; 1.64728483e-06 1.29193545e-10 NL: NLBGS 10 ; 1.5370671e-07 1.2054937e-11 NL: NLBGS 11 ; 5.86333262e-09 4.59850484e-13 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 99441.4903 1 NL: NLBGS 2 ; 97674.0839 0.982226669 NL: NLBGS 3 ; 16423.3169 0.165155579 NL: NLBGS 4 ; 936.878493 0.00942140438 NL: NLBGS 5 ; 35.1251723 0.000353224516 NL: NLBGS 6 ; 1.74058079 1.7503567e-05 NL: NLBGS 7 ; 0.0617140596 6.20606745e-07 NL: NLBGS 8 ; 0.00344670525 3.46606355e-08 NL: NLBGS 9 ; 6.40780037e-05 6.44378955e-10 NL: NLBGS 10 ; 1.18287731e-05 1.1895209e-10 NL: NLBGS 11 ; 2.12569858e-06 2.13763749e-11 NL: NLBGS 12 ; 5.24892486e-08 5.27840527e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 5539.43468 1 NL: NLBGS 2 ; 5477.32024 0.988786864 NL: NLBGS 3 ; 590.748547 0.106644194 NL: NLBGS 4 ; 34.2115619 0.00617600241 NL: NLBGS 5 ; 0.168620199 3.04399653e-05 NL: NLBGS 6 ; 0.0113633451 2.05135465e-06 NL: NLBGS 7 ; 0.000326410906 5.8924949e-08 NL: NLBGS 8 ; 1.5720212e-06 2.83787298e-10 NL: NLBGS 9 ; 7.44351779e-08 1.34373239e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 11777.2144 1 NL: NLBGS 2 ; 10952.5649 0.929979244 NL: NLBGS 3 ; 1788.45454 0.151857178 NL: NLBGS 4 ; 60.9781265 0.00517763578 NL: NLBGS 5 ; 0.437799145 3.71734038e-05 NL: NLBGS 6 ; 0.0508806373 4.32026078e-06 NL: NLBGS 7 ; 0.000849916796 7.21661991e-08 NL: NLBGS 8 ; 1.59619689e-05 1.35532634e-09 NL: NLBGS 9 ; 7.85334203e-07 6.66825091e-11 NL: NLBGS 10 ; 3.41791072e-08 2.9021385e-12 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 8615.56355 1 NL: NLBGS 2 ; 8655.28974 1.00461098 NL: NLBGS 3 ; 663.946766 0.0770636491 NL: NLBGS 4 ; 53.7985308 0.00624434264 NL: NLBGS 5 ; 0.203596967 2.36313e-05 NL: NLBGS 6 ; 0.00470022774 5.45550818e-07 NL: NLBGS 7 ; 0.000316431831 3.67279319e-08 NL: NLBGS 8 ; 9.10786959e-07 1.05714148e-10 NL: NLBGS 9 ; 1.62236906e-08 1.88306784e-12 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 8615.70264 1 NL: NLBGS 2 ; 8051.10149 0.934468358 NL: NLBGS 3 ; 913.855255 0.10606857 NL: NLBGS 4 ; 63.16661 0.00733156803 NL: NLBGS 5 ; 0.444283501 5.15667172e-05 NL: NLBGS 6 ; 0.0142047939 1.64870987e-06 NL: NLBGS 7 ; 0.00076555045 8.88552544e-08 NL: NLBGS 8 ; 4.09016357e-06 4.74733605e-10 NL: NLBGS 9 ; 1.12087446e-07 1.30096697e-11 NL: NLBGS 10 ; 8.74768891e-09 1.01531927e-12 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 7540.08452 1 NL: NLBGS 2 ; 7559.68166 1.00259906 NL: NLBGS 3 ; 531.752285 0.0705233851 NL: NLBGS 4 ; 44.7494365 0.005934872 NL: NLBGS 5 ; 0.177640554 2.3559491e-05 NL: NLBGS 6 ; 0.00265966629 3.52736934e-07 NL: NLBGS 7 ; 0.000195422982 2.59178768e-08 NL: NLBGS 8 ; 4.46395023e-07 5.92029203e-11 NL: NLBGS 9 ; 6.74856433e-09 8.95025023e-13 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 7724.76613 1 NL: NLBGS 2 ; 7175.74992 0.928927789 NL: NLBGS 3 ; 551.431117 0.07138483 NL: NLBGS 4 ; 53.1207267 0.00687667766 NL: NLBGS 5 ; 0.323619426 4.18937506e-05 NL: NLBGS 6 ; 0.00299123921 3.87227155e-07 NL: NLBGS 7 ; 0.000234576535 3.03668139e-08 NL: NLBGS 8 ; 3.75142897e-07 4.85636575e-11 NL: NLBGS 9 ; 2.41490015e-09 3.12617898e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 2014.41317 1 NL: NLBGS 2 ; 1978.12375 0.981985117 NL: NLBGS 3 ; 139.189496 0.0690967964 NL: NLBGS 4 ; 11.4813476 0.00569959919 NL: NLBGS 5 ; 0.0460692172 2.28697955e-05 NL: NLBGS 6 ; 0.000750100698 3.72366857e-07 NL: NLBGS 7 ; 5.38851181e-05 2.67497845e-08 NL: NLBGS 8 ; 1.13957998e-07 5.65713133e-11 NL: NLBGS 9 ; 1.80139007e-09 8.94250546e-13 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 2312.19683 1 NL: NLBGS 2 ; 1944.22814 0.84085754 NL: NLBGS 3 ; 129.459481 0.0559898184 NL: NLBGS 4 ; 12.893806 0.00557643097 NL: NLBGS 5 ; 0.0841375873 3.63885921e-05 NL: NLBGS 6 ; 0.00106515598 4.60668385e-07 NL: NLBGS 7 ; 8.60360144e-05 3.72096412e-08 NL: NLBGS 8 ; 8.42602936e-08 3.64416613e-11 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 183.212481 1 NL: NLBGS 2 ; 181.59268 0.991158896 NL: NLBGS 3 ; 11.7045833 0.0638852942 NL: NLBGS 4 ; 1.03329787 0.00563988797 NL: NLBGS 5 ; 0.00410609139 2.24116358e-05 NL: NLBGS 6 ; 2.76986324e-05 1.51183108e-07 NL: NLBGS 7 ; 2.23586963e-06 1.22036972e-08 NL: NLBGS 8 ; 4.3823043e-09 2.39192454e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 195.409124 1 NL: NLBGS 2 ; 168.768785 0.863668909 NL: NLBGS 3 ; 12.7712212 0.0653563196 NL: NLBGS 4 ; 1.36292825 0.00697474214 NL: NLBGS 5 ; 0.00816984045 4.18088996e-05 NL: NLBGS 6 ; 0.000186725362 9.5556112e-07 NL: NLBGS 7 ; 1.61939077e-05 8.28718092e-08 NL: NLBGS 8 ; 4.3231067e-08 2.21233615e-10 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 1343.47867 1 NL: NLBGS 2 ; 1346.98623 1.00261081 NL: NLBGS 3 ; 91.1437478 0.067841604 NL: NLBGS 4 ; 7.67892275 0.00571570129 NL: NLBGS 5 ; 0.0304483752 2.26638323e-05 NL: NLBGS 6 ; 0.000360670882 2.68460446e-07 NL: NLBGS 7 ; 2.70141243e-05 2.01075945e-08 NL: NLBGS 8 ; 5.7636712e-08 4.29010995e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 1370.3534 1 NL: NLBGS 2 ; 1272.6822 0.928725544 NL: NLBGS 3 ; 92.39475 0.0674240312 NL: NLBGS 4 ; 9.40052849 0.00685993009 NL: NLBGS 5 ; 0.0597843705 4.3626973e-05 NL: NLBGS 6 ; 0.000974536059 7.11156744e-07 NL: NLBGS 7 ; 8.10426379e-05 5.91399548e-08 NL: NLBGS 8 ; 1.65072587e-07 1.20459867e-10 NL: NLBGS 9 ; 8.67706865e-10 6.33199339e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 4822.62203 1 NL: NLBGS 2 ; 4842.50586 1.00412303 NL: NLBGS 3 ; 324.50954 0.0672890263 NL: NLBGS 4 ; 26.4175138 0.00547783212 NL: NLBGS 5 ; 0.100280525 2.07937765e-05 NL: NLBGS 6 ; 0.00118738247 2.46210974e-07 NL: NLBGS 7 ; 8.59497537e-05 1.7822204e-08 NL: NLBGS 8 ; 1.87186772e-07 3.88143153e-11 NL: NLBGS 9 ; 2.22078532e-09 4.60493339e-13 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 4860.41105 1 NL: NLBGS 2 ; 4534.53956 0.932953925 NL: NLBGS 3 ; 335.124981 0.0689499258 NL: NLBGS 4 ; 33.0028568 0.00679013699 NL: NLBGS 5 ; 0.214345699 4.41003234e-05 NL: NLBGS 6 ; 0.0030550078 6.28549268e-07 NL: NLBGS 7 ; 0.000247604507 5.09431208e-08 NL: NLBGS 8 ; 6.09803127e-07 1.25463283e-10 NL: NLBGS 9 ; 2.66849728e-09 5.49027079e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 3485.25768 1 NL: NLBGS 2 ; 3500.46703 1.00436391 NL: NLBGS 3 ; 235.64132 0.0676108747 NL: NLBGS 4 ; 19.6262611 0.0056312224 NL: NLBGS 5 ; 0.0724877317 2.07983852e-05 NL: NLBGS 6 ; 0.000907449966 2.60368113e-07 NL: NLBGS 7 ; 6.70855168e-05 1.92483664e-08 NL: NLBGS 8 ; 1.51352021e-07 4.34263503e-11 NL: NLBGS 9 ; 1.90054938e-09 5.45311007e-13 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 3511.55641 1 NL: NLBGS 2 ; 3277.6173 0.933380222 NL: NLBGS 3 ; 242.93978 0.0691829354 NL: NLBGS 4 ; 24.3704278 0.00694006444 NL: NLBGS 5 ; 0.146858967 4.18216169e-05 NL: NLBGS 6 ; 0.0023281715 6.63002733e-07 NL: NLBGS 7 ; 0.000189994096 5.41053807e-08 NL: NLBGS 8 ; 3.98646693e-07 1.13524217e-10 NL: NLBGS 9 ; 2.15855747e-09 6.14701066e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 664.060481 1 NL: NLBGS 2 ; 665.207743 1.00172765 NL: NLBGS 3 ; 44.893418 0.0676044114 NL: NLBGS 4 ; 3.72008872 0.00560203298 NL: NLBGS 5 ; 0.0144829441 2.18096763e-05 NL: NLBGS 6 ; 0.000183233148 2.75928403e-07 NL: NLBGS 7 ; 1.34410271e-05 2.0240667e-08 NL: NLBGS 8 ; 2.87785254e-08 4.33372053e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 694.408037 1 NL: NLBGS 2 ; 641.696794 0.92409183 NL: NLBGS 3 ; 45.6395911 0.0657244568 NL: NLBGS 4 ; 4.41769808 0.00636181877 NL: NLBGS 5 ; 0.028068028 4.04200793e-05 NL: NLBGS 6 ; 0.000310838863 4.4763143e-07 NL: NLBGS 7 ; 2.4509455e-05 3.52954656e-08 NL: NLBGS 8 ; 2.3369394e-08 3.36536917e-11 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 655.050457 1 NL: NLBGS 2 ; 659.006221 1.00603887 NL: NLBGS 3 ; 43.2608947 0.06604208 NL: NLBGS 4 ; 3.56473056 0.00544191752 NL: NLBGS 5 ; 0.0129012738 1.96950841e-05 NL: NLBGS 6 ; 0.000148694262 2.26996655e-07 NL: NLBGS 7 ; 1.0965564e-05 1.67400296e-08 NL: NLBGS 8 ; 2.32941556e-08 3.55608569e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 634.621514 1 NL: NLBGS 2 ; 598.263059 0.942708442 NL: NLBGS 3 ; 44.82813 0.0706375833 NL: NLBGS 4 ; 4.63589998 0.00730498396 NL: NLBGS 5 ; 0.0275378923 4.33926233e-05 NL: NLBGS 6 ; 0.000570957062 8.99681226e-07 NL: NLBGS 7 ; 4.82788208e-05 7.60749828e-08 NL: NLBGS 8 ; 1.43492117e-07 2.26106607e-10 NL: NLBGS 9 ; 7.82655301e-10 1.23326311e-12 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 2960.222 1 NL: NLBGS 2 ; 2980.21055 1.00675238 NL: NLBGS 3 ; 187.696189 0.0634061194 NL: NLBGS 4 ; 14.5587701 0.00491813455 NL: NLBGS 5 ; 0.0499636812 1.68783562e-05 NL: NLBGS 6 ; 0.000566507807 1.9137342e-07 NL: NLBGS 7 ; 3.94359404e-05 1.33219537e-08 NL: NLBGS 8 ; 8.39264344e-08 2.83513988e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 2848.59411 1 NL: NLBGS 2 ; 2691.45476 0.944836176 NL: NLBGS 3 ; 195.880671 0.0687639811 NL: NLBGS 4 ; 19.0753199 0.00669639799 NL: NLBGS 5 ; 0.112759199 3.9584158e-05 NL: NLBGS 6 ; 0.00202812291 7.11973288e-07 NL: NLBGS 7 ; 0.000163275998 5.73180986e-08 NL: NLBGS 8 ; 5.03378208e-07 1.76711103e-10 NL: NLBGS 9 ; 2.59469079e-09 9.10867148e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 6471.9453 1 NL: NLBGS 2 ; 6529.67288 1.00891966 NL: NLBGS 3 ; 375.035847 0.0579479321 NL: NLBGS 4 ; 25.4106985 0.0039262845 NL: NLBGS 5 ; 0.0768211936 1.18698768e-05 NL: NLBGS 6 ; 0.00085326348 1.31840342e-07 NL: NLBGS 7 ; 5.22200246e-05 8.06867521e-09 NL: NLBGS 8 ; 1.21330463e-07 1.87471397e-11 NL: NLBGS 9 ; 1.16034893e-09 1.7928905e-13 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 6087.07023 1 NL: NLBGS 2 ; 5797.62261 0.952448778 NL: NLBGS 3 ; 406.616453 0.0668000266 NL: NLBGS 4 ; 34.31692 0.0056376744 NL: NLBGS 5 ; 0.196216134 3.22349056e-05 NL: NLBGS 6 ; 0.00243557383 4.00122512e-07 NL: NLBGS 7 ; 0.000173606935 2.85206066e-08 NL: NLBGS 8 ; 5.7236871e-07 9.40302458e-11 NL: NLBGS 9 ; 2.61780747e-09 4.30060337e-13 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 2057.66048 1 NL: NLBGS 2 ; 2078.23617 1.00999955 NL: NLBGS 3 ; 113.732228 0.0552725917 NL: NLBGS 4 ; 7.27351709 0.00353484802 NL: NLBGS 5 ; 0.0187756501 9.12475614e-06 NL: NLBGS 6 ; 0.000222350741 1.08059975e-07 NL: NLBGS 7 ; 1.28941088e-05 6.26639277e-09 NL: NLBGS 8 ; 3.44857088e-08 1.67596691e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 1877.21453 1 NL: NLBGS 2 ; 1778.10123 0.947201931 NL: NLBGS 3 ; 138.538803 0.0738001975 NL: NLBGS 4 ; 10.6368758 0.00566630806 NL: NLBGS 5 ; 0.0541784652 2.88610941e-05 NL: NLBGS 6 ; 0.00042094052 2.24236768e-07 NL: NLBGS 7 ; 2.73666583e-05 1.45783328e-08 NL: NLBGS 8 ; 9.11444311e-08 4.85530181e-11 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 1497.98558 1 NL: NLBGS 2 ; 1515.29925 1.01155797 NL: NLBGS 3 ; 80.902143 0.0540072908 NL: NLBGS 4 ; 5.01742194 0.00334944608 NL: NLBGS 5 ; 0.0127501103 8.51150405e-06 NL: NLBGS 6 ; 0.000146538169 9.78234842e-08 NL: NLBGS 7 ; 8.26042821e-06 5.51435762e-09 NL: NLBGS 8 ; 2.32460288e-08 1.55181926e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 1366.21726 1 NL: NLBGS 2 ; 1309.71385 0.958642445 NL: NLBGS 3 ; 94.3270529 0.0690424984 NL: NLBGS 4 ; 7.07182244 0.00517620634 NL: NLBGS 5 ; 0.0368336883 2.69603448e-05 NL: NLBGS 6 ; 0.000294931586 2.15874587e-07 NL: NLBGS 7 ; 1.87322768e-05 1.37110527e-08 NL: NLBGS 8 ; 6.39670159e-08 4.68205298e-11 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 731.499093 1 NL: NLBGS 2 ; 740.121272 1.011787 NL: NLBGS 3 ; 39.0067248 0.0533243652 NL: NLBGS 4 ; 2.3769829 0.003249468 NL: NLBGS 5 ; 0.00585368333 8.00231112e-06 NL: NLBGS 6 ; 6.83052099e-05 9.3377026e-08 NL: NLBGS 7 ; 3.78495627e-06 5.17424601e-09 NL: NLBGS 8 ; 1.09510881e-08 1.49707473e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 664.841379 1 NL: NLBGS 2 ; 636.909943 0.95798782 NL: NLBGS 3 ; 46.2035734 0.0694956343 NL: NLBGS 4 ; 3.36143401 0.00505599398 NL: NLBGS 5 ; 0.0172914521 2.60083874e-05 NL: NLBGS 6 ; 0.000117098111 1.76129397e-07 NL: NLBGS 7 ; 7.18132547e-06 1.08015621e-08 NL: NLBGS 8 ; 2.42947524e-08 3.65421786e-11 NL: NLBGS Converged ================== AS_point_0.coupled ================== NL: NLBGS 1 ; 26.8045091 1 NL: NLBGS 2 ; 27.1582821 1.01319826 NL: NLBGS 3 ; 1.37361708 0.0512457467 NL: NLBGS 4 ; 0.0835572222 0.00311728232 NL: NLBGS 5 ; 0.000181945209 6.78785828e-06 NL: NLBGS 6 ; 2.20226058e-06 8.21600789e-08 NL: NLBGS 7 ; 1.22292281e-07 4.56237721e-09 NL: NLBGS 8 ; 3.37861909e-10 1.26046669e-11 NL: NLBGS Converged ================== AS_point_1.coupled ================== NL: NLBGS 1 ; 24.4803462 1 NL: NLBGS 2 ; 23.3242428 0.952774224 NL: NLBGS 3 ; 1.73922447 0.0710457465 NL: NLBGS 4 ; 0.113499497 0.00463635178 NL: NLBGS 5 ; 0.000568308863 2.3214903e-05 NL: NLBGS 6 ; 2.82165022e-06 1.1526186e-07 NL: NLBGS 7 ; 1.6185974e-07 6.611824e-09 NL: NLBGS 8 ; 5.12382513e-10 2.09303622e-11 NL: NLBGS Converged Optimization terminated successfully (Exit mode 0) Current function value: 0.026705653242651523 Iterations: 16 Function evaluations: 18 Gradient evaluations: 16 Optimization Complete -----------------------------------
# prob.run_model()
print("The fuel burn value is", prob["AS_point_0.fuelburn"][0], "[kg]")
The fuel burn value is 2670.565324265152 [kg]
print(
"The wingbox mass (excluding the wing_weight_ratio) is",
prob["wing.structural_mass"][0] / surf_dict["wing_weight_ratio"],
"[kg]",
)
The wingbox mass (excluding the wing_weight_ratio) is 1124.0854960417034 [kg]
There is plenty of room for improvement. A finer mesh and a tighter optimization tolerance should be used.