Conjugate Heat Transfer
This tutorial uses Modulus to study the conjugate heat transfer between the heat sink and the surrounding fluid. The temperature variations inside solid and fluid would be solved in a coupled manner with appropriate interface boundary conditions. In this tutorial, you will learn:
How to generate a 3D geometry using the geometry module in Modulus.
How to set up a conjugate heat transfer problem using the interface boundary conditions in Modulus.
How to use the Multi-Phase training approach in Modulus for one way coupled problems.
This tutorial assumes that you have completed tutorial Introductory Example on and have familiarized yourself with the basics of the Modulus APIs. Also, you should review the Scalar Transport: 2D Advection Diffusion tutorial for additional details on writing some of the thermal boundary conditions.
The scripts used in this problem can be found at examples/three_fin_3d/
. The scripts are made configurable such that the same three fin problem can be solved for higher reynolds numbers,
parameterized geometry, etc. very easily by only configuring the Hydra configs. This example focuses on the laminar variant of this as the purpose of this tutorial is to cover
the ideas related to multi-phase training and conjugate heat transfer. The tutorial Parameterized 3D Heat Sink covers the parameterization aspect of this problem.
Therefore in this problem, the custom
configs for turbulent
and parameterized
are both set to false
.
The geometry for a 3-fin heat sink placed inside a channel is shown in Fig. 142. The inlet to the channel is at 1 \(m/s\). The pressure at the outlet is specified as 0 \(Pa\). All the other surfaces of the geometry are treated as no-slip walls.
Fig. 142 Three fin heat sink geometry (All dimensions in \(m\))
The inlet is at 273.15 \(K\). The channel walls are adiabatic. The heat sink has a heat source of \(0.2 \times 0.4\) \(m\) at the bottom of the heat sink situated centrally on the bottom surface. The heat source generates heat such that the temperature gradient on the source surface is 360 \(K/m\) in the normal direction. Conjugate heat transfer takes place between the fluid-solid contact surface.
The properties fluid and thermal properties of the fluid and the solid are as follows:
Property |
Fluid |
Solid |
Kinematic Viscosity \((m^2/s)\) |
0.02 |
NA |
Thermal Diffusivity \((m^2/s)\) |
0.02 |
0.0625 |
Thermal Conductivity \((W/m.K)\) |
1.0 |
5.0 |
In this tutorial, since you are dealing with only incompressible flow, there is a one way coupling between the heat and flow equations. This means that it is possible to train the temperature field after the flow field is trained and converged. Such an approach is useful while training the multiphysics problems which are one way coupled as it is possible to achieve significant speed-up, as well as simulate cases with same flow boundary conditions but different thermal boundary conditions. One can easily use the same flow field as in input to train for different thermal boundary conditions.
Therefore, for this problem you have three separate files for the geometry,
flow solver, and heat solver. The three_fin_geometry.py
will contain
all the definitions of geometry. three_fin_flow.py
and three_fin_thermal.py
would then use this geometry to setup the relevant flow and heat constraints and
solve them individually. The basic idea would be to train the flow model to convergence
and then start the heat training after by initializing from the trained flow model to solve
for the temperature distributions in fluid and solid simultaneously.
In this problem you will nondimensionalize the temperature according to the following equation:
(207)\[\theta= T/273.15 - 1.0\]
Creating Geometry
The three_fin_geometry.py
script contains all the details relevant to the geometry generation.
We will use the Box
primitive to create the
3-fin geometry and Channel
primitive to generate the channel.
Similar to 2D, Channel
and Box
are defined by using the two
corner points. Like 2D, the Channel
geometry has no bounding planes
in the x-direction. You will also make use of the repeat
method to
create the fins. This speeds up the generation of repetitive
structures in comparison to constructing the same geometry separately and doing boolean operations to assemble them.
Use the Plane
geometry to create the planes at the inlet and
outlet. The code for generating the required geometries is shown below.
Please note the normal directions for the inlet and outlet planes.
Additionally, the parameters required for solving the heat part as also defined upfront, ex. dimensions and locations of source etc.
The script contains a few extra definitions that are only relevant for parameterized geometry. These are not relevant for this tutorial.
# geometry params for domain
channel_origin = (-2.5, -0.5, -0.5)
channel_dim = (5.0, 1.0, 1.0)
heat_sink_base_origin = (-1.0, -0.5, -0.3)
heat_sink_base_dim = (1.0, 0.2, 0.6)
fin_origin = (heat_sink_base_origin[0] + 0.5 - fin_length_s / 2, -0.3, -0.3)
fin_dim = (fin_length_s, fin_height_s, fin_thickness_s) # two side fins
total_fins = 2 # two side fins
flow_box_origin = (-1.1, -0.5, -0.5)
flow_box_dim = (1.6, 1.0, 1.0)
source_origin = (-0.7, -0.5, -0.1)
source_dim = (0.4, 0.0, 0.2)
source_area = 0.08
# define geometry
class ThreeFin(object):
def __init__(self, parameterized: bool = False):
# set param ranges
if parameterized:
pr = Parameterization(param_ranges)
self.pr = param_ranges
else:
pr = Parameterization(fixed_param_ranges)
self.pr = fixed_param_ranges
# channel
self.channel = Channel(
channel_origin,
(
channel_origin[0] + channel_dim[0],
channel_origin[1] + channel_dim[1],
channel_origin[2] + channel_dim[2],
),
parameterization=pr,
)
# three fin heat sink
heat_sink_base = Box(
heat_sink_base_origin,
(
heat_sink_base_origin[0] + heat_sink_base_dim[0], # base of heat sink
heat_sink_base_origin[1] + heat_sink_base_dim[1],
heat_sink_base_origin[2] + heat_sink_base_dim[2],
),
parameterization=pr,
)
fin_center = (
fin_origin[0] + fin_dim[0] / 2,
fin_origin[1] + fin_dim[1] / 2,
fin_origin[2] + fin_dim[2] / 2,
)
fin = Box(
fin_origin,
(
fin_origin[0] + fin_dim[0],
fin_origin[1] + fin_dim[1],
fin_origin[2] + fin_dim[2],
),
parameterization=pr,
)
gap = (heat_sink_base_dim[2] - fin_dim[2]) / (
total_fins - 1
) # gap between fins
fin_2 = fin.translate([0, 0, gap])
fin = fin + fin_2
three_fin = heat_sink_base + fin
# parameterized center fin
center_fin_origin = (
heat_sink_base_origin[0] + 0.5 - fin_length_m / 2,
fin_origin[1],
-fin_thickness_m / 2,
)
center_fin_dim = (fin_length_m, fin_height_m, fin_thickness_m)
center_fin = Box(
center_fin_origin,
(
center_fin_origin[0] + center_fin_dim[0],
center_fin_origin[1] + center_fin_dim[1],
center_fin_origin[2] + center_fin_dim[2],
),
parameterization=pr,
)
self.three_fin = three_fin + center_fin
# entire geometry
self.geo = self.channel - self.three_fin
# low and high resultion geo away and near the heat sink
flow_box = Box(
flow_box_origin,
(
flow_box_origin[0] + flow_box_dim[0], # base of heat sink
flow_box_origin[1] + flow_box_dim[1],
flow_box_origin[2] + flow_box_dim[2],
),
)
self.lr_geo = self.geo - flow_box
self.hr_geo = self.geo & flow_box
lr_bounds_x = (channel_origin[0], channel_origin[0] + channel_dim[0])
lr_bounds_y = (channel_origin[1], channel_origin[1] + channel_dim[1])
lr_bounds_z = (channel_origin[2], channel_origin[2] + channel_dim[2])
self.lr_bounds = {x: lr_bounds_x, y: lr_bounds_y, z: lr_bounds_z}
hr_bounds_x = (flow_box_origin[0], flow_box_origin[0] + flow_box_dim[0])
hr_bounds_y = (flow_box_origin[1], flow_box_origin[1] + flow_box_dim[1])
hr_bounds_z = (flow_box_origin[2], flow_box_origin[2] + flow_box_dim[2])
self.hr_bounds = {x: hr_bounds_x, y: hr_bounds_y, z: hr_bounds_z}
# inlet and outlet
self.inlet = Plane(
channel_origin,
(
channel_origin[0],
channel_origin[1] + channel_dim[1],
channel_origin[2] + channel_dim[2],
),
-1,
parameterization=pr,
)
self.outlet = Plane(
(channel_origin[0] + channel_dim[0], channel_origin[1], channel_origin[2]),
(
channel_origin[0] + channel_dim[0],
channel_origin[1] + channel_dim[1],
channel_origin[2] + channel_dim[2],
),
1,
parameterization=pr,
)
# planes for integral continuity
self.integral_plane = Plane(
(x_pos, channel_origin[1], channel_origin[2]),
(
x_pos,
channel_origin[1] + channel_dim[1],
channel_origin[2] + channel_dim[2],
),
1,
)
Neural network, Nodes and Multi-Phase training
Let’s have a look at the networks and nodes required to solve the flow and heat for this problem. The architectures and nodes for flow problem are very similar to previous tutorials. You will add the nodes for NavierStokes
and NormalDotVec
and create a single flow network that has the coordinates as inputs and the velocity components and the pressure as output. The code for the flow nodes can be found here:
# make navier stokes equations
if cfg.custom.turbulent:
ze = ZeroEquation(nu=0.002, dim=3, time=False, max_distance=0.5)
ns = NavierStokes(nu=ze.equations["nu"], rho=1.0, dim=3, time=False)
navier_stokes_nodes = ns.make_nodes() + ze.make_nodes()
else:
ns = NavierStokes(nu=0.01, rho=1.0, dim=3, time=False)
navier_stokes_nodes = ns.make_nodes()
normal_dot_vel = NormalDotVec()
# make network arch
if cfg.custom.parameterized:
input_keys = [
Key("x"),
Key("y"),
Key("z"),
Key("fin_height_m"),
Key("fin_height_s"),
Key("fin_length_m"),
Key("fin_length_s"),
Key("fin_thickness_m"),
Key("fin_thickness_s"),
]
else:
input_keys = [Key("x"), Key("y"), Key("z")]
flow_net = FullyConnectedArch(
input_keys=input_keys, output_keys=[Key("u"), Key("v"), Key("w"), Key("p")]
)
# make list of nodes to unroll graph on
flow_nodes = (
navier_stokes_nodes
+ normal_dot_vel.make_nodes()
+ [flow_net.make_node(name="flow_network")]
)
For the thermal nodes, start by adding nodes for relevant equations like AdvectionDiffusion
, Diffusion
, DiffusionInterface
and GradNormal
that will be used to define the various thermal boundary conditions relevant to this problem. Also, create 3 separate neural networks flow_net
, thermal_f_net
and thermal_s_net
. The first one is the same flow network defined in the flow scripts. This network architecture definition in heat script must exactly match to that of the flow script for successful initialization of the flow model during heat training. Set the optimize
argument as False
while making the nodes of flow network to avoid optimizing the flow network during the heat training. Finally, separate networks to predict the temperatures in fluid and solid are created.
# make thermal equations
ad = AdvectionDiffusion(T="theta_f", rho=1.0, D=0.02, dim=3, time=False)
dif = Diffusion(T="theta_s", D=0.0625, dim=3, time=False)
dif_inteface = DiffusionInterface("theta_f", "theta_s", 1.0, 5.0, dim=3, time=False)
f_grad = GradNormal("theta_f", dim=3, time=False)
s_grad = GradNormal("theta_s", dim=3, time=False)
# make network arch
if cfg.custom.parameterized:
input_keys = [
Key("x"),
Key("y"),
Key("z"),
Key("fin_height_m"),
Key("fin_height_s"),
Key("fin_length_m"),
Key("fin_length_s"),
Key("fin_thickness_m"),
Key("fin_thickness_s"),
]
else:
input_keys = [Key("x"), Key("y"), Key("z")]
flow_net = FullyConnectedArch(
input_keys=input_keys,
output_keys=[Key("u"), Key("v"), Key("w"), Key("p")],
)
thermal_f_net = FullyConnectedArch(
input_keys=input_keys, output_keys=[Key("theta_f")]
)
thermal_s_net = FullyConnectedArch(
input_keys=input_keys, output_keys=[Key("theta_s")]
)
# make list of nodes to unroll graph on
thermal_nodes = (
ad.make_nodes()
+ dif.make_nodes()
+ dif_inteface.make_nodes()
+ f_grad.make_nodes()
+ s_grad.make_nodes()
+ [flow_net.make_node(name="flow_network", optimize=False)]
+ [thermal_f_net.make_node(name="thermal_f_network")]
+ [thermal_s_net.make_node(name="thermal_s_network")]
)
Setting up Flow Domain and Constraints
The contents of three_fin_flow.py
script are described below.
Inlet, Outlet and Channel and Heat Sink walls
For inlet boundary conditions, specify the velocity to be a constant velocity of 1.0 \(m/s\) in x-direction. Like in tutorial Introductory Example, weight the velocity by the SDF of the channel to avoid sharp discontinuity at the boundaries. For outlet, specify the pressure to be 0. All the channel walls and heat sink walls are treated as no slip boundaries.
Interior
The flow equations can be specified in the low resolution and high
resolution domains of the problem by using PointwiseInteriorConstraint
. This allows
independent point densities in these two areas to be controlled easily.
Integral Continuity
The inlet volumetric flow is 1 \(m^3/s\) so,
specify 1.0 as the value for integral_continuity
key.
The code for flow domain is shown below.
# make flow domain
flow_domain = Domain()
# inlet
u_profile = inlet_vel * tanh((0.5 - Abs(y)) / 0.02) * tanh((0.5 - Abs(z)) / 0.02)
constraint_inlet = PointwiseBoundaryConstraint(
nodes=flow_nodes,
geometry=geo.inlet,
outvar={"u": u_profile, "v": 0, "w": 0},
batch_size=cfg.batch_size.Inlet,
criteria=Eq(x, channel_origin[0]),
lambda_weighting={
"u": 1.0,
"v": 1.0,
"w": 1.0,
}, # weight zero on edges
parameterization=geo.pr,
batch_per_epoch=5000,
)
flow_domain.add_constraint(constraint_inlet, "inlet")
# outlet
constraint_outlet = PointwiseBoundaryConstraint(
nodes=flow_nodes,
geometry=geo.outlet,
outvar={"p": 0},
batch_size=cfg.batch_size.Outlet,
criteria=Eq(x, channel_origin[0] + channel_dim[0]),
lambda_weighting={"p": 1.0},
parameterization=geo.pr,
batch_per_epoch=5000,
)
flow_domain.add_constraint(constraint_outlet, "outlet")
# no slip for channel walls
no_slip = PointwiseBoundaryConstraint(
nodes=flow_nodes,
geometry=geo.geo,
outvar={"u": 0, "v": 0, "w": 0},
batch_size=cfg.batch_size.NoSlip,
lambda_weighting={
"u": 1.0,
"v": 1.0,
"w": 1.0,
}, # weight zero on edges
parameterization=geo.pr,
batch_per_epoch=5000,
)
flow_domain.add_constraint(no_slip, "no_slip")
# flow interior low res away from three fin
lr_interior = PointwiseInteriorConstraint(
nodes=flow_nodes,
geometry=geo.geo,
outvar={"continuity": 0, "momentum_x": 0, "momentum_y": 0, "momentum_z": 0},
batch_size=cfg.batch_size.InteriorLR,
lambda_weighting={
"continuity": Symbol("sdf"),
"momentum_x": Symbol("sdf"),
"momentum_y": Symbol("sdf"),
"momentum_z": Symbol("sdf"),
},
compute_sdf_derivatives=True,
parameterization=geo.pr,
batch_per_epoch=5000,
criteria=Or(x < -1.1, x > 0.5),
)
flow_domain.add_constraint(lr_interior, "lr_interior")
# flow interiror high res near three fin
hr_interior = PointwiseInteriorConstraint(
nodes=flow_nodes,
geometry=geo.geo,
outvar={"continuity": 0, "momentum_x": 0, "momentum_z": 0, "momentum_y": 0},
batch_size=cfg.batch_size.InteriorLR,
lambda_weighting={
"continuity": Symbol("sdf"),
"momentum_x": Symbol("sdf"),
"momentum_y": Symbol("sdf"),
"momentum_z": Symbol("sdf"),
},
compute_sdf_derivatives=True,
parameterization=geo.pr,
batch_per_epoch=5000,
criteria=And(x > -1.1, x < 0.5),
)
flow_domain.add_constraint(hr_interior, "hr_interior")
# integral continuity
def integral_criteria(invar, params):
sdf = geo.geo.sdf(invar, params)
return np.greater(sdf["sdf"], 0)
integral_continuity = IntegralBoundaryConstraint(
nodes=flow_nodes,
geometry=geo.integral_plane,
outvar={"normal_dot_vel": volumetric_flow},
batch_size=5,
integral_batch_size=cfg.batch_size.IntegralContinuity,
criteria=integral_criteria,
lambda_weighting={"normal_dot_vel": 1.0},
parameterization={**geo.pr, **{x_pos: (-1.1, 0.1)}},
fixed_dataset=False,
num_workers=4,
)
flow_domain.add_constraint(integral_continuity, "integral_continuity")
The addition of integral continuity planes and separate flow box for dense sampling are examples of adding more training data and user knowledge/bias of the problem to the training. This addition helps to improve the accuracy and convergence to a great extent and it is recommended wherever possible.
Setting up Thermal Multi-Phase Domain and Constraints
The contents of three_fin_thermal.py
are described below.
Inlet, Outlet and Channel walls:
For the heat part, specify temperature at the inlet. All the
outlet and the channel walls will have a zero gradient boundary
condition which will be enforced by setting
'normal_gradient_theta_f'
equal to 0. We will use 'theta_f'
for
defining the temperatures in fluid and 'theta_s'
for defining the
temperatures in solid.
Fluid and Solid Interior:
Just like the Scalar Transport: 2D Advection Diffusion tutorial,
set 'advection_diffusion'
equal to 0 in both, low and high
resolution fluid domains. For solid interior, set
'diffusion'
equal to 0.
Fluid-Solid Interface:
For the interface between fluid and solid, enforce both Neumann
and Dirichlet boundary condition by setting
'diffusion_interface_dirichlet_theta_f_theta_s'
and
'diffusion_interface_neumann_theta_f_theta_s'
both equal to 0.
The order in which you define 'theta_f'
and 'theta_s'
in the interface boundary condition must match the one in PDE definition DiffusionInterface
to avoid any graph unroll errors.
The corresponding conductivities must be also specified in the same order in the PDE definition.
Heat Source:
Apply a \(tanh\) smoothing while defining the heat source on
the bottom wall of the heat sink. Smoothing out the
sharp boundaries helps in training the Neural Network converge faster.
The 'normal_gradient_theta_s'
is set equal to grad_t
in the area
of source and 0 everywhere else on the bottom surface of heat sink.
The code for heat domain is shown below.
# make flow domain
thermal_domain = Domain()
# inlet
constraint_inlet = PointwiseBoundaryConstraint(
nodes=thermal_nodes,
geometry=geo.inlet,
outvar={"theta_f": inlet_t},
batch_size=cfg.batch_size.Inlet,
criteria=Eq(x, channel_origin[0]),
lambda_weighting={"theta_f": 1.0}, # weight zero on edges
parameterization=geo.pr,
)
thermal_domain.add_constraint(constraint_inlet, "inlet")
# outlet
constraint_outlet = PointwiseBoundaryConstraint(
nodes=thermal_nodes,
geometry=geo.outlet,
outvar={"normal_gradient_theta_f": 0},
batch_size=cfg.batch_size.Outlet,
criteria=Eq(x, channel_origin[0] + channel_dim[0]),
lambda_weighting={"normal_gradient_theta_f": 1.0}, # weight zero on edges
parameterization=geo.pr,
)
thermal_domain.add_constraint(constraint_outlet, "outlet")
# channel walls insulating
def wall_criteria(invar, params):
sdf = geo.three_fin.sdf(invar, params)
return np.less(sdf["sdf"], -1e-5)
channel_walls = PointwiseBoundaryConstraint(
nodes=thermal_nodes,
geometry=geo.channel,
outvar={"normal_gradient_theta_f": 0},
batch_size=cfg.batch_size.ChannelWalls,
criteria=wall_criteria,
lambda_weighting={"normal_gradient_theta_f": 1.0},
parameterization=geo.pr,
)
thermal_domain.add_constraint(channel_walls, "channel_walls")
# fluid solid interface
def interface_criteria(invar, params):
sdf = geo.channel.sdf(invar, params)
return np.greater(sdf["sdf"], 0)
fluid_solid_interface = PointwiseBoundaryConstraint(
nodes=thermal_nodes,
geometry=geo.three_fin,
outvar={
"diffusion_interface_dirichlet_theta_f_theta_s": 0,
"diffusion_interface_neumann_theta_f_theta_s": 0,
},
batch_size=cfg.batch_size.SolidInterface,
criteria=interface_criteria,
parameterization=geo.pr,
)
thermal_domain.add_constraint(fluid_solid_interface, "fluid_solid_interface")
# heat source
sharpen_tanh = 60.0
source_func_xl = (tanh(sharpen_tanh * (x - source_origin[0])) + 1.0) / 2.0
source_func_xh = (
tanh(sharpen_tanh * ((source_origin[0] + source_dim[0]) - x)) + 1.0
) / 2.0
source_func_zl = (tanh(sharpen_tanh * (z - source_origin[2])) + 1.0) / 2.0
source_func_zh = (
tanh(sharpen_tanh * ((source_origin[2] + source_dim[2]) - z)) + 1.0
) / 2.0
gradient_normal = (
grad_t * source_func_xl * source_func_xh * source_func_zl * source_func_zh
)
heat_source = PointwiseBoundaryConstraint(
nodes=thermal_nodes,
geometry=geo.three_fin,
outvar={"normal_gradient_theta_s": gradient_normal},
batch_size=cfg.batch_size.HeatSource,
criteria=Eq(y, source_origin[1]),
)
thermal_domain.add_constraint(heat_source, "heat_source")
# flow interior low res away from three fin
lr_flow_interior = PointwiseInteriorConstraint(
nodes=thermal_nodes,
geometry=geo.geo,
outvar={"advection_diffusion_theta_f": 0},
batch_size=cfg.batch_size.InteriorLR,
criteria=Or(x < -1.1, x > 0.5),
)
thermal_domain.add_constraint(lr_flow_interior, "lr_flow_interior")
# flow interiror high res near three fin
hr_flow_interior = PointwiseInteriorConstraint(
nodes=thermal_nodes,
geometry=geo.geo,
outvar={"advection_diffusion_theta_f": 0},
batch_size=cfg.batch_size.InteriorHR,
criteria=And(x > -1.1, x < 0.5),
)
thermal_domain.add_constraint(hr_flow_interior, "hr_flow_interior")
# solid interior
solid_interior = PointwiseInteriorConstraint(
nodes=thermal_nodes,
geometry=geo.three_fin,
outvar={"diffusion_theta_s": 0},
batch_size=cfg.batch_size.SolidInterior,
lambda_weighting={"diffusion_theta_s": 100.0},
)
thermal_domain.add_constraint(solid_interior, "solid_interior")
Adding Validators and Monitors
In this tutorial you will monitors for pressure drops during the flow field simulation and monitors for peak temperature reached in the source chip during the heat simulation.
Similarly, respective validation data are added to the flow and heat scripts. Only flow monitors and validators are shown for brevity.
Monitors and validators in flow script:
# flow data
mapping = {
"Points:0": "x",
"Points:1": "y",
"Points:2": "z",
"U:0": "u",
"U:1": "v",
"U:2": "w",
"p_rgh": "p",
}
if cfg.custom.turbulent:
openfoam_var = csv_to_dict(
to_absolute_path("openfoam/threeFin_extend_zeroEq_re500_fluid.csv"), mapping
)
else:
openfoam_var = csv_to_dict(
to_absolute_path("openfoam/threeFin_extend_fluid0.csv"), mapping
)
openfoam_var = {key: value[0::4] for key, value in openfoam_var.items()}
openfoam_var["x"] = openfoam_var["x"] + channel_origin[0]
openfoam_var["y"] = openfoam_var["y"] + channel_origin[1]
openfoam_var["z"] = openfoam_var["z"] + channel_origin[2]
openfoam_var.update({"fin_height_m": np.full_like(openfoam_var["x"], 0.4)})
openfoam_var.update({"fin_height_s": np.full_like(openfoam_var["x"], 0.4)})
openfoam_var.update({"fin_thickness_m": np.full_like(openfoam_var["x"], 0.1)})
openfoam_var.update({"fin_thickness_s": np.full_like(openfoam_var["x"], 0.1)})
openfoam_var.update({"fin_length_m": np.full_like(openfoam_var["x"], 1.0)})
openfoam_var.update({"fin_length_s": np.full_like(openfoam_var["x"], 1.0)})
openfoam_invar_numpy = {
key: value
for key, value in openfoam_var.items()
if key
in [
"x",
"y",
"z",
"fin_height_m",
"fin_height_s",
"fin_thickness_m",
"fin_thickness_s",
"fin_length_m",
"fin_length_s",
]
}
openfoam_outvar_numpy = {
key: value for key, value in openfoam_var.items() if key in ["u", "v", "w", "p"]
}
openfoam_validator = PointwiseValidator(
nodes=flow_nodes,
invar=openfoam_invar_numpy,
true_outvar=openfoam_outvar_numpy,
)
flow_domain.add_validator(openfoam_validator)
# add pressure monitor
invar_inlet_pressure = geo.integral_plane.sample_boundary(
1024, parameterization={**fixed_param_ranges, **{x_pos: -2}}
)
pressure_monitor = PointwiseMonitor(
invar_inlet_pressure,
output_names=["p"],
metrics={"inlet_pressure": lambda var: torch.mean(var["p"])},
nodes=flow_nodes,
)
flow_domain.add_monitor(pressure_monitor)
Once both the flow and heat scripts are defined, run the
three_fin_flow.py
first to solve for the flow field. Once a
desired level of convergence is achieved you can run three_fin_thermal.py
to solve for heat.
The table and figures below show the results of Pressure drop and Peak temperatures obtained from the Modulus and compare it with the results from OpenFOAM solver.
Modulus |
OpenFOAM |
|
Pressure Drop \((Pa)\) |
7.51 |
7.49 |
Peak Temperature \((^{\circ} C)\) |
78.35 |
78.05 |
Fig. 143 Left: Modulus. Center: OpenFOAM. Right: Difference. Top: Temperature distribution in Fluid. Bottom: Temperature distribution in Solid (Temperature scales in C)
Fig. 144 Left: Modulus. Center: OpenFOAM. Right: Difference. (Temperature scales in C)
Plotting gradient quantities: Wall Velocity Gradients
In a variety of applications, it is desirable to plot the gradients of
some quantities inside the domain. One such example relevant to fluid
flows is the wall velocity gradients and wall shear stresses. These
quantities are often plotted to compute frictional forces, etc. You can
visualize such quantities in Modulus by outputting the \(x\),
\(y\) and \(z\) derivatives of the desired variables using an
PointwiseInferencer
.
...
# add inferencer
inferencer = PointwiseInferencer(
geo.sample_boundary(4000, parameterization=pr),
["u__x", "u__y", "u__z",
"v__x", "v__y", "v__z",
"w__x", "w__y", "w__z"],
nodes=flow_nodes,
)
flow_domain.add_inferencer(inferencer, "inf_data")
You can then post-process these quantities based on your choice to visualize the desired variables. Paraview’s Calculator Filter was used for the plot shown below. The wall velocity gradients comparison between OpenFOAM and Modulus is shown in Fig. 145.
Fig. 145 Comparison of magnitude of wall velocity gradients. Left: Modulus. Right: OpenFOAM