-
Notifications
You must be signed in to change notification settings - Fork 85
/
runTests_MphysAeroThermal.py
executable file
·334 lines (296 loc) · 11.1 KB
/
runTests_MphysAeroThermal.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
#!/usr/bin/env python
"""
Run Python tests for optimization integration
"""
from mpi4py import MPI
import os
import numpy as np
from testFuncs import *
import openmdao.api as om
from mphys.multipoint import Multipoint
from dafoam.mphys import DAFoamBuilder, OptFuncs
from funtofem.mphys import MeldThermalBuilder
from pygeo import geo_utils
from mphys.scenario_aerothermal import ScenarioAeroThermal
from pygeo.mphys import OM_DVGEOCOMP
gcomm = MPI.COMM_WORLD
os.chdir("./reg_test_files-main/ChannelConjugateHeat")
if gcomm.rank == 0:
os.system("rm -rf */processor*")
# aero setup
U0 = 10.0
daOptionsAero = {
"designSurfaces": [
"hot_air_inner",
"hot_air_outer",
"hot_air_sides",
"cold_air_outer",
"cold_air_inner",
"cold_air_sides",
],
"solverName": "DARhoSimpleFoam",
"primalMinResTol": 1.0e-17,
"primalMinResTolDiff": 1.0e17,
"discipline": "aero",
"primalBC": {
"UHot": {"variable": "U", "patches": ["hot_air_in"], "value": [U0, 0.0, 0.0]},
"UCold": {"variable": "U", "patches": ["cold_air_in"], "value": [-U0, 0.0, 0.0]},
"useWallFunction": False,
},
"objFunc": {
"PL": {
"part1": {
"type": "totalPressure",
"source": "patchToFace",
"patches": ["hot_air_in"],
"scale": 1.0,
"addToAdjoint": True,
},
"part2": {
"type": "totalPressure",
"source": "patchToFace",
"patches": ["hot_air_out"],
"scale": -1.0,
"addToAdjoint": True,
},
"part3": {
"type": "totalPressure",
"source": "patchToFace",
"patches": ["cold_air_in"],
"scale": 1.0,
"addToAdjoint": True,
},
"part4": {
"type": "totalPressure",
"source": "patchToFace",
"patches": ["cold_air_out"],
"scale": -1.0,
"addToAdjoint": True,
},
},
"TOUT": {
"part1": {
"type": "patchMean",
"source": "patchToFace",
"patches": ["hot_air_out"],
"varName": "T",
"varType": "scalar",
"component": 0,
"scale": 1.0,
"addToAdjoint": True,
}
},
"HFH": {
"part1": {
"type": "wallHeatFlux",
"source": "patchToFace",
"patches": ["hot_air_inner"],
"scale": 1,
"addToAdjoint": True,
}
},
"HFC": {
"part1": {
"type": "wallHeatFlux",
"source": "patchToFace",
"patches": ["cold_air_outer"],
"scale": 1,
"addToAdjoint": True,
}
},
},
"couplingInfo": {
"aerothermal": {
"active": True,
"couplingSurfaceGroups": {
"wallGroup": ["hot_air_inner", "cold_air_outer"],
},
}
},
"adjStateOrdering": "cell",
"adjEqnOption": {
"gmresRelTol": 1.0e-3,
"pcFillLevel": 1,
"jacMatReOrdering": "natural",
"useNonZeroInitGuess": True,
},
"normalizeStates": {
"U": U0,
"p": 101325,
"nuTilda": 1e-3,
"T": 300,
"phi": 1.0,
},
"designVar": {
"shape": {"designVarType": "FFD"},
},
}
daOptionsThermal = {
"designSurfaces": ["channel_outer", "channel_inner", "channel_sides"],
"solverName": "DAHeatTransferFoam",
"primalMinResTol": 1.0e-17,
"primalMinResTolDiff": 1.0e17,
"discipline": "thermal",
"objFunc": {
"HF_INNER": {
"part1": {
"type": "wallHeatFlux",
"source": "patchToFace",
"patches": ["channel_inner"],
"scale": 1,
"addToAdjoint": True,
}
},
"HF_OUTER": {
"part1": {
"type": "wallHeatFlux",
"source": "patchToFace",
"patches": ["channel_outer"],
"scale": 1,
"addToAdjoint": True,
}
},
},
"couplingInfo": {
"aerothermal": {
"active": True,
"couplingSurfaceGroups": {
"wallGroup": ["channel_outer", "channel_inner"],
},
}
},
"adjStateOrdering": "cell",
"adjEqnOption": {
"gmresRelTol": 1.0e-3,
"pcFillLevel": 1,
"jacMatReOrdering": "natural",
"useNonZeroInitGuess": True,
},
"normalizeStates": {
"T": 300.0,
},
"designVar": {
"shape": {"designVarType": "FFD"},
"power": {
"designVarType": "HSC",
"heatSourceName": "source1",
"comps": [6, 7, 8],
},
},
"fvSource": {
"source1": {
"type": "heatSource",
"source": "cylinderSmooth",
"center": [0.2, 0.055, 0.025],
"axis": [1.0, 0.0, 0.0],
"radius": 0.005,
"length": 0.4,
"power": 10.0,
"eps": 0.001,
},
},
}
# Mesh deformation setup
meshOptions = {
"gridFile": os.getcwd(),
"fileType": "OpenFOAM",
# point and normal for the symmetry plane
"symmetryPlanes": [],
}
def power(val, DASolver):
radius = float(val[0])
length = float(val[1])
power = float(val[2])
# only change design variables
DASolver.setOption("fvSource", {"source1": {"power": power, "radius": radius, "length": length}})
DASolver.updateDAOption()
class Top(Multipoint):
def setup(self):
# create the builder to initialize the DASolvers for both cases (they share the same mesh option)
dafoam_builder_aero = DAFoamBuilder(daOptionsAero, meshOptions, scenario="aerothermal", run_directory="aero")
dafoam_builder_aero.initialize(self.comm)
dafoam_builder_thermal = DAFoamBuilder(
daOptionsThermal, meshOptions, scenario="aerothermal", run_directory="thermal"
)
dafoam_builder_thermal.initialize(self.comm)
thermalxfer_builder = MeldThermalBuilder(dafoam_builder_aero, dafoam_builder_thermal, n=1, beta=0.5)
thermalxfer_builder.initialize(self.comm)
# add the design variable component to keep the top level design variables
self.add_subsystem("dvs", om.IndepVarComp(), promotes=["*"])
# add the mesh component
self.add_subsystem("mesh_aero", dafoam_builder_aero.get_mesh_coordinate_subsystem())
self.add_subsystem("mesh_thermal", dafoam_builder_thermal.get_mesh_coordinate_subsystem())
# add the geometry component (FFD). Note that the aero and thermal use the exact same FFD file
self.add_subsystem("geometry_aero", OM_DVGEOCOMP(file="aero/FFD/channelFFD.xyz", type="ffd"))
self.add_subsystem("geometry_thermal", OM_DVGEOCOMP(file="aero/FFD/channelFFD.xyz", type="ffd"))
# add a scenario (flow condition) for optimization, we pass the builder
# to the scenario to actually run the flow and adjoint
self.mphys_add_scenario(
"scenario",
ScenarioAeroThermal(
aero_builder=dafoam_builder_aero,
thermal_builder=dafoam_builder_thermal,
thermalxfer_builder=thermalxfer_builder,
),
om.NonlinearBlockGS(maxiter=10, iprint=2, use_aitken=True, rtol=1e-10, atol=1e-6),
om.LinearBlockGS(maxiter=10, iprint=2, use_aitken=True, rtol=1e-8, atol=1e-4),
)
# need to manually connect the x_aero0 between the mesh and geometry components
self.connect("mesh_aero.x_aero0", "geometry_aero.x_aero_in")
self.connect("geometry_aero.x_aero0", "scenario.x_aero")
self.connect("mesh_thermal.x_thermal0", "geometry_thermal.x_thermal_in")
self.connect("geometry_thermal.x_thermal0", "scenario.x_thermal")
def configure(self):
super().configure()
# add the objective function to the cruise scenario
self.scenario.aero_post.mphys_add_funcs()
self.scenario.thermal_post.mphys_add_funcs()
# get the surface coordinates from the mesh component
points_aero = self.mesh_aero.mphys_get_surface_mesh()
points_thermal = self.mesh_thermal.mphys_get_surface_mesh()
# add pointset to the geometry component
self.geometry_aero.nom_add_discipline_coords("aero", points_aero)
self.geometry_thermal.nom_add_discipline_coords("thermal", points_thermal)
# self.scenario.coupling._mphys_promote_coupling_variables()
# select the FFD points to move
pts = self.geometry_aero.DVGeo.getLocalIndex(0)
indexList = pts[3:6, :, :].flatten()
PS = geo_utils.PointSelect("list", indexList)
nShapes = self.geometry_aero.nom_addLocalDV(dvName="shape", axis="y", pointSelect=PS)
nShapes = self.geometry_thermal.nom_addLocalDV(dvName="shape", axis="y", pointSelect=PS)
self.scenario.coupling.thermal.solver.add_dv_func("power", power)
self.scenario.thermal_post.add_dv_func("power", power)
# add the design variables to the dvs component's output
self.dvs.add_output("shape", val=np.array([0] * nShapes))
self.dvs.add_output("power", val=np.array([0.005, 0.4, 0.0]))
# manually connect the dvs output to the geometry and cruise
# sa and sst cases share the same shape
self.connect("shape", "geometry_aero.shape")
self.connect("shape", "geometry_thermal.shape")
self.connect("power", "scenario.power")
# define the design variables to the top level
self.add_design_var("shape", lower=-1, upper=1, scaler=1.0)
self.add_design_var("power", lower=-100.0, upper=100.0, scaler=1.0)
# add objective and constraints to the top level
self.add_objective("scenario.aero_post.HFH", scaler=1.0)
self.add_constraint("scenario.aero_post.PL", upper=1000, scaler=1.0)
self.add_constraint("scenario.thermal_post.HF_INNER", upper=1000, scaler=1.0)
prob = om.Problem(reports=None)
prob.model = Top()
prob.setup(mode="rev")
om.n2(prob, show_browser=False, outfile="mphys_aerothermal.html")
optFuncs = OptFuncs([daOptionsAero, daOptionsThermal], prob)
prob.run_model()
totals = prob.compute_totals()
if gcomm.rank == 0:
derivDict = {}
derivDict["PL"] = {}
derivDict["PL"]["shape"] = totals[("scenario.aero_post.functionals.PL", "dvs.shape")][0]
derivDict["PL"]["power"] = totals[("scenario.aero_post.functionals.PL", "dvs.power")][0]
derivDict["HFH"] = {}
derivDict["HFH"]["shape"] = totals[("scenario.aero_post.functionals.HFH", "dvs.shape")][0]
derivDict["HFH"]["power"] = totals[("scenario.aero_post.functionals.HFH", "dvs.power")][0]
derivDict["HF_INNER"] = {}
derivDict["HF_INNER"]["shape"] = totals[("scenario.thermal_post.functionals.HF_INNER", "dvs.shape")][0]
derivDict["HF_INNER"]["power"] = totals[("scenario.thermal_post.functionals.HF_INNER", "dvs.power")][0]
reg_write_dict(derivDict, 1e-4, 1e-6)