Skip to content

Commit

Permalink
Merge 55db29c into d6ba203
Browse files Browse the repository at this point in the history
  • Loading branch information
gbarter committed Jun 11, 2022
2 parents d6ba203 + 55db29c commit 9496bfe
Show file tree
Hide file tree
Showing 79 changed files with 891 additions and 152 deletions.
12 changes: 9 additions & 3 deletions examples/17_jacket/analysis_options_jacket.yaml
Original file line number Diff line number Diff line change
Expand Up @@ -4,10 +4,10 @@ general:
design_variables:

jacket:
r_foot:
foot_head_ratio:
flag: True
lower_bound: 6.
upper_bound: 16.
lower_bound: 1.1
upper_bound: 5.0
r_head:
flag: True
lower_bound: 2.
Expand Down Expand Up @@ -39,6 +39,12 @@ constraints:
flag: False
shell_buckling:
flag: True
frequency_1:
flag: False #True
lower_bound: 0.13
upper_bound: 0.40
tower_diameter_coupling:
flag: True

driver:
optimization:
Expand Down
30 changes: 27 additions & 3 deletions wisdem/fixed_bottomse/jacket.py
Original file line number Diff line number Diff line change
Expand Up @@ -29,20 +29,24 @@ def setup(self):
n_bays = mod_opt["WISDEM"]["FixedBottomSE"]["n_bays"]

self.add_input("height", val=70.0, units="m")
self.add_input("r_foot", val=10.0, units="m")
self.add_input("r_head", val=6.0, units="m")
self.add_input("foot_head_ratio", val=1.5)
self.add_input("bay_spacing", val=np.linspace(0.1, 0.9, n_bays + 1))
self.add_input("tower_base_diameter", val=0.0, units="m")

self.add_output("r_foot", val=10.0, units="m")
self.add_output("leg_nodes", val=np.zeros((n_legs, n_bays + 2, 3)), units="m")
self.add_output("bay_nodes", val=np.zeros((n_legs, n_bays + 1, 3)), units="m")
self.add_output("constr_diam_consistency", val=0.0)

def compute(self, inputs, outputs):
mod_opt = self.options["modeling_options"]
n_legs = mod_opt["WISDEM"]["FixedBottomSE"]["n_legs"]
n_bays = mod_opt["WISDEM"]["FixedBottomSE"]["n_bays"]

r_foot = inputs["r_foot"]
ratio = inputs["foot_head_ratio"]
r_head = inputs["r_head"]
outputs["r_foot"] = r_foot = r_head * ratio
height = inputs["height"]

leg_spacing = np.linspace(0, 1.0, n_bays + 2)
Expand All @@ -60,6 +64,8 @@ def compute(self, inputs, outputs):
outputs["bay_nodes"][i, :, 1] = x_bay * np.sin(i / n_legs * 2 * np.pi)
outputs["bay_nodes"][i, :, 2] = z_bay

outputs["constr_diam_consistency"] = inputs["tower_base_diameter"] / r_head


class ComputeFrame3DD(om.ExplicitComponent):
"""
Expand Down Expand Up @@ -146,6 +152,14 @@ def setup(self):
self.add_output("jacket_Myy", np.zeros((n_elem, n_dlc)), units="N*m")
self.add_output("jacket_Mzz", np.zeros((n_elem, n_dlc)), units="N*m")

# Frequencies
self.add_output("f1", val=0.0, units="Hz")
self.add_output("f2", val=0.0, units="Hz")
self.add_output("structural_frequencies", np.zeros(NFREQ), units="Hz")

self.add_output("jacket_deflection", np.zeros((n_node, n_dlc)), units="m")
self.add_output("top_deflection", np.zeros(n_dlc), units="m")

self.idx = 0

def compute(self, inputs, outputs, discrete_inputs, discrete_outputs):
Expand Down Expand Up @@ -442,6 +456,15 @@ def add_element(n1_nodes, n1_indices, n2_nodes, n2_indices, itube, idx1, idx2, i
# self.frame.draw()
displacements, forces, reactions, internalForces, mass, modal = self.frame.run()

# natural frequncies
outputs["f1"] = modal.freq[0]
outputs["f2"] = modal.freq[1]
outputs["structural_frequencies"] = modal.freq[:NFREQ]

# deflections due to loading (from cylinder top and wind/wave loads)
outputs["jacket_deflection"] = np.sqrt(displacements.dx**2 + displacements.dy**2).T
outputs["top_deflection"] = outputs["jacket_deflection"][-1, :]

# Determine reaction forces
outputs["jacket_base_F"] = -np.c_[
reactions.Fx.sum(axis=1), reactions.Fy.sum(axis=1), reactions.Fz.sum(axis=1)
Expand Down Expand Up @@ -633,7 +656,8 @@ def compute(self, inputs, outputs, discrete_inputs, discrete_outputs):
leg_L = np.linalg.norm(leg_vec, axis=1)

# Get the angle of intersection between all edges (as vectors) and all legs (as vectors)
leg_alpha = np.arccos(np.dot(edge_vec, leg_vec.T) / np.outer(L, leg_L))
vec_vals = np.dot(edge_vec, leg_vec.T) / np.outer(L, leg_L)
leg_alpha = np.arccos(np.minimum(np.maximum(vec_vals, -1.0), 1.0))
# If angle of intersection is close to 0 or 180, the edge is part of a leg
tol = np.deg2rad(5)
idx_leg = np.any((np.abs(leg_alpha) < tol) | (np.abs(leg_alpha - np.pi) < tol), axis=1)
Expand Down
24 changes: 20 additions & 4 deletions wisdem/glue_code/gc_LoadInputs.py
Original file line number Diff line number Diff line change
@@ -1,4 +1,5 @@
import numpy as np

import wisdem.inputs as sch


Expand Down Expand Up @@ -217,7 +218,7 @@ def set_openmdao_vectors(self):
)

joint_pos = self.wt_init["components"]["blade"]["internal_structure_2d_fem"]["joint"]["position"]
if joint_pos>0.:
if joint_pos > 0.0:
self.modeling_options["WISDEM"]["RotorSE"]["bjs"] = True
# Adjust grid to have grid point at join location
closest_grid_pt = np.argmin(abs(self.modeling_options["WISDEM"]["RotorSE"]["nd_span"] - joint_pos))
Expand Down Expand Up @@ -886,7 +887,7 @@ def write_ontology(self, wt_opt, fname_output):
"grid"
] = wt_opt["blade.internal_structure_2d_fem.s"].tolist()
self.wt_init["components"]["blade"]["internal_structure_2d_fem"]["layers"][i]["width"][
"values"
"values"
] = wt_opt["blade.internal_structure_2d_fem.layer_width"][i, :].tolist()
if (
wt_opt["blade.internal_structure_2d_fem.definition_layer"][i] == 2
Expand All @@ -902,7 +903,7 @@ def write_ontology(self, wt_opt, fname_output):
"grid"
] = wt_opt["blade.internal_structure_2d_fem.s"].tolist()
self.wt_init["components"]["blade"]["internal_structure_2d_fem"]["layers"][i]["offset_y_pa"][
"values"
"values"
] = wt_opt["blade.internal_structure_2d_fem.layer_offset_y_pa"][i, :].tolist()
if (
wt_opt["blade.internal_structure_2d_fem.definition_layer"][i] == 4
Expand All @@ -927,7 +928,9 @@ def write_ontology(self, wt_opt, fname_output):
# TODO assign joint mass to wt_init from rs.bjs
# Elastic properties of the blade
if self.modeling_options["WISDEM"]["RotorSE"]["bjs"]:
self.wt_init["components"]["blade"]["internal_structure_2d_fem"]["joint"]["mass"] = float(wt_opt["rotorse.rs.bjs.joint_mass"])
self.wt_init["components"]["blade"]["internal_structure_2d_fem"]["joint"]["mass"] = float(
wt_opt["rotorse.rs.bjs.joint_mass"]
)

self.wt_init["components"]["blade"]["elastic_properties_mb"] = {}
self.wt_init["components"]["blade"]["elastic_properties_mb"]["six_x_six"] = {}
Expand Down Expand Up @@ -1232,6 +1235,19 @@ def write_ontology(self, wt_opt, fname_output):
"values"
] = wt_opt["monopile.layer_thickness"][i, :].tolist()

# Update jacket
if self.modeling_options["flags"]["jacket"]:
self.wt_init["components"]["jacket"]["r_head"] = float(wt_opt["jacket.r_head"][0])
self.wt_init["components"]["jacket"]["r_foot"] = float(
wt_opt["jacket.r_head"] * wt_opt["jacket.foot_head_ratio"]
)
self.wt_init["components"]["jacket"]["height"] = float(wt_opt["jacket.height"][0])
self.wt_init["components"]["jacket"]["leg_diameter"] = float(wt_opt["jacket.leg_diameter"][0])
self.wt_init["components"]["jacket"]["leg_thickness"] = float(wt_opt["jacket.leg_thickness"][0])
self.wt_init["components"]["jacket"]["brace_diameters"] = wt_opt["jacket.brace_diameters"].tolist()
self.wt_init["components"]["jacket"]["brace_thicknesses"] = wt_opt["jacket.brace_thicknesses"].tolist()
self.wt_init["components"]["jacket"]["bay_spacing"] = wt_opt["jacket.bay_spacing"].tolist()

# Update floating platform and mooring
if self.modeling_options["flags"]["floating"]:
yaml_out = self.wt_init["components"]["floating_platform"]
Expand Down
12 changes: 7 additions & 5 deletions wisdem/glue_code/gc_PoseOptimization.py
Original file line number Diff line number Diff line change
Expand Up @@ -835,12 +835,11 @@ def set_design_variables(self, wt_opt, wt_init):
)

# -- Jacket --
if jacket_opt["r_foot"]["flag"]:
if jacket_opt["foot_head_ratio"]["flag"]:
wt_opt.model.add_design_var(
"jacket.r_foot",
lower=jacket_opt["r_foot"]["lower_bound"],
upper=jacket_opt["r_foot"]["upper_bound"],
ref=5.0,
"jacket.foot_head_ratio",
lower=jacket_opt["foot_head_ratio"]["lower_bound"],
upper=jacket_opt["foot_head_ratio"]["upper_bound"],
)

if jacket_opt["r_head"]["flag"]:
Expand Down Expand Up @@ -1370,6 +1369,9 @@ def set_constraints(self, wt_opt):
if jacket_constr["shell_buckling"]["flag"]:
wt_opt.model.add_constraint("fixedse.constr_shell_buckling", upper=1.0)

if jacket_constr["tower_diameter_coupling"]["flag"]:
wt_opt.model.add_constraint("fixedse.constr_diam_consistency", upper=1.0)

if jacket_constr["frequency_1"]["flag"]:
wt_opt.model.add_constraint(
"fixedse.structural_frequencies",
Expand Down
9 changes: 4 additions & 5 deletions wisdem/glue_code/gc_WT_DataStruc.py
Original file line number Diff line number Diff line change
Expand Up @@ -2171,7 +2171,7 @@ def setup(self):

class Compute_Grid(om.ExplicitComponent):
"""
Compute the non-dimensional grid for a tower or monopile/jacket.
Compute the non-dimensional grid for a tower or monopile.
Using the dimensional `ref_axis` array, this component computes the
non-dimensional grid, height (vertical distance) and length (curve distance)
Expand Down Expand Up @@ -2298,10 +2298,9 @@ def setup(self):

ivc = self.add_subsystem("jacket_indep_vars", om.IndepVarComp(), promotes=["*"])
ivc.add_output(
"r_foot",
val=0.0,
units="m",
desc="Radius of foot (bottom) of jacket, in meters.",
"foot_head_ratio",
val=1.5,
desc="Ratio of radius of foot (bottom) of jacket to head.",
)
ivc.add_output(
"r_head",
Expand Down
2 changes: 1 addition & 1 deletion wisdem/glue_code/gc_WT_InitModel.py
Original file line number Diff line number Diff line change
Expand Up @@ -977,8 +977,8 @@ def assign_jacket_values(wt_opt, modeling_options, jacket):
wt_opt["jacket.transition_piece_mass"] = jacket["transition_piece_mass"]
wt_opt["jacket.transition_piece_cost"] = jacket["transition_piece_cost"]
wt_opt["jacket.gravity_foundation_mass"] = jacket["gravity_foundation_mass"]
wt_opt["jacket.r_foot"] = jacket["r_foot"]
wt_opt["jacket.r_head"] = jacket["r_head"]
wt_opt["jacket.foot_head_ratio"] = jacket["r_foot"] / jacket["r_head"]
wt_opt["jacket.height"] = jacket["height"]
wt_opt["jacket.leg_diameter"] = jacket["leg_diameter"]
wt_opt["jacket.leg_thickness"] = jacket["leg_thickness"]
Expand Down
22 changes: 18 additions & 4 deletions wisdem/glue_code/glue_code.py
Original file line number Diff line number Diff line change
Expand Up @@ -516,8 +516,6 @@ def setup(self):
self.connect("env.Tsig_wave", "fixedse.Tsig_wave")
self.connect("monopile.diameter", "fixedse.monopile_outer_diameter_in")
self.connect("monopile.diameter", "fixedse.monopile_top_diameter", src_indices=[-1])
if modeling_options["flags"]["tower"]:
self.connect("tower.diameter", "fixedse.tower_base_diameter", src_indices=[0])
self.connect("monopile.foundation_height", "fixedse.monopile_foundation_height")
self.connect("monopile.outfitting_factor", "fixedse.outfitting_factor_in")
self.connect("monopile.height", "fixedse.monopile_height")
Expand All @@ -539,6 +537,7 @@ def setup(self):
self.connect("monopile.transition_piece_mass", "fixedse.transition_piece_mass")
self.connect("monopile.gravity_foundation_mass", "fixedse.gravity_foundation_mass")
if modeling_options["flags"]["tower"]:
self.connect("tower.diameter", "fixedse.tower_base_diameter", src_indices=[0])
self.connect("towerse.tower_mass", "fixedse.tower_mass")
self.connect("towerse.tower_cost", "fixedse.tower_cost")
self.connect("towerse.turbine_mass", "fixedse.turbine_mass")
Expand All @@ -564,7 +563,7 @@ def setup(self):
self.connect("materials.G", "fixedse.G_mat")
self.connect("materials.rho", "fixedse.rho_mat")
self.connect("materials.name", "fixedse.material_names")
self.connect("jacket.r_foot", "fixedse.r_foot")
self.connect("jacket.foot_head_ratio", "fixedse.foot_head_ratio")
self.connect("jacket.r_head", "fixedse.r_head")
self.connect("jacket.height", "fixedse.height")
self.connect("jacket.leg_diameter", "fixedse.leg_diameter")
Expand All @@ -580,6 +579,7 @@ def setup(self):
self.connect("towerse.turbine_I_base", "fixedse.turbine_I")
self.connect("towerse.tower.turbine_F", "fixedse.turbine_F")
self.connect("towerse.tower.turbine_M", "fixedse.turbine_M")
self.connect("tower.diameter", "fixedse.tower_base_diameter", src_indices=[0])

if modeling_options["flags"]["floating"]:
self.connect("env.rho_water", "floatingse.rho_water")
Expand Down Expand Up @@ -844,7 +844,14 @@ def setup(self):
self.add_subsystem("wt", WT_RNTA(modeling_options=modeling_options, opt_options=opt_options), promotes=["*"])
if modeling_options["WISDEM"]["BOS"]["flag"]:
if modeling_options["flags"]["offshore"]:
self.add_subsystem("orbit", Orbit(floating=modeling_options["flags"]["floating_platform"]))
self.add_subsystem(
"orbit",
Orbit(
floating=modeling_options["flags"]["floating"],
jacket=modeling_options["flags"]["jacket"],
jacket_legs=modeling_options["WISDEM"]["FixedBottomSE"]["n_legs"],
),
)
else:
self.add_subsystem("landbosse", LandBOSSE())

Expand Down Expand Up @@ -874,6 +881,13 @@ def setup(self):
self.connect("monopile.transition_piece_mass", "orbit.transition_piece_mass")
self.connect("monopile.transition_piece_cost", "orbit.transition_piece_cost")
self.connect("monopile.diameter", "orbit.monopile_diameter", src_indices=[0])
elif modeling_options["flags"]["jacket"]:
self.connect("fixedse.r_foot", "orbit.jacket_r_foot")
self.connect("jacket.height", "orbit.jacket_length")
self.connect("fixedse.jacket_mass", "orbit.jacket_mass")
self.connect("fixedse.jacket_cost", "orbit.jacket_cost")
self.connect("jacket.transition_piece_mass", "orbit.transition_piece_mass")
self.connect("jacket.transition_piece_cost", "orbit.transition_piece_cost")
elif modeling_options["flags"]["floating"]:
self.connect("mooring.n_lines", "orbit.num_mooring_lines")
self.connect("floatingse.line_mass", "orbit.mooring_line_mass", src_indices=[0])
Expand Down
9 changes: 4 additions & 5 deletions wisdem/inputs/analysis_schema.yaml
Original file line number Diff line number Diff line change
Expand Up @@ -589,18 +589,17 @@ properties:
description: Design variables associated with the jacket
default: {}
properties:
r_foot:
foot_head_ratio:
type: object
description: Adjust the radius of the jacket foot.
description: Adjust the ratio of the jacket foot (bottom) radius to that of the head (top)
default: {}
properties:
flag: *flag
lower_bound: &rfootbound
type: number
minimum: 0.1
minimum: 1.0
maximum: 100.0
default: 5.0
unit: m
default: 1.5
description: Design variable bound
upper_bound: *rfootbound
r_head:
Expand Down
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
storage_specs:
max_cargo: 6000 # t
max_cargo: 10000 # t
max_deck_load: 8 # t/m^2
max_deck_space: 600 # m^2
transport_specs:
Expand Down

0 comments on commit 9496bfe

Please sign in to comment.