Annulus Shape Optimization — SciPy + StochasticControl
This example demonstrates how to couple SciPy optimization routines with StochasticControl to minimize the maximum temperature in a 2D annular geometry while enforcing a constant volume constraint. It is a practical pattern for black‑box optimization where each function evaluation runs a MOOSE simulation and returns quantities of interest (QoIs).
Problem Description
In this problem, we focus on an annular (ring-shaped) domain subjected to thermal conditions. The domain is defined by two concentric circles with an inner radius and thickness. The physical problem is heat conduction governed by the heat equation with a convective flux type boundary condition at the inner radius and an insulated outer radius. The convective heat transfer coefficient is dependent on the inner radius. A constant source term is distributed throughout the domain.
The optimization goal is to find the inner radius and thickness that minimize the maximum temperature within the domain while maintaining a fixed volume of the annulus. In other words, the design variables are the inner radius and thickness, and the objective function is the maximum temperature within the annulus. This is a constrained optimization problem because the volume of the annulus needs to remain constant during the optimization process.
This example problem is identical to the shape optimization example in the Optimization module. The optimization module version utilizes gradients for its minimization solve, which makes it a faster method, but is more involved to formulate.
Physics Inputs
There are two MOOSE inputs representing the physics. The first utilizes input variables for the inner radius and thickness to directly change the mesh that is generated. The other uses mesh displacement to change the domain of the ring.
Mesh-Based Perturbation Input
The mesh-based input starts with definition of the input variables to-be optimized:
inner_radius = 6
thickness = 4
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)The mesh is then generated based on these values:
[Mesh<<<{"href": "../../../syntax/Mesh/index.html"}>>>]
[disk]
type = ConcentricCircleMeshGenerator<<<{"description": "This ConcentricCircleMeshGenerator source code is to generate concentric circle meshes.", "href": "../../../source/meshgenerators/ConcentricCircleMeshGenerator.html"}>>>
has_outer_square<<<{"description": "It determines if meshes for a outer square are added to concentric circle meshes."}>>> = false
preserve_volumes<<<{"description": "Volume of concentric circles can be preserved using this function."}>>> = false
radii<<<{"description": "Radii of major concentric circles"}>>> = '${inner_radius} ${fparse inner_radius + thickness}'
rings<<<{"description": "Number of rings in each circle or in the enclosing square"}>>> = '16 16'
num_sectors<<<{"description": "num_sectors % 2 = 0, num_sectors > 0Number of azimuthal sectors in each quadrant'num_sectors' must be an even number."}>>> = 16
[]
[ring]
type = BlockDeletionGenerator<<<{"description": "Mesh generator which removes elements from the specified subdomains", "href": "../../../source/meshgenerators/BlockDeletionGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = disk
block<<<{"description": "The list of blocks to be deleted"}>>> = 1
new_boundary<<<{"description": "optional boundary name to assign to the cut surface"}>>> = inner
[]
[]
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)The heat conduction problem is then defined:
[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
[T]
[]
[]
[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
[diffusion]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = k
[]
[src]
type = ADBodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../../../source/kernels/BodyForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 1
[]
[]
[BCs<<<{"href": "../../../syntax/BCs/index.html"}>>>]
[convection]
type = ADMatNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{C \\partial u}{\\partial n}=M*h$, where $h$ is a constant, $M$ is a material property, and $C$ is a coefficient defined by the kernel for $u$.", "href": "../../../source/bcs/MatNeumannBC.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = inner
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
boundary_material<<<{"description": "Material property multiplying the constant that will be enforced by the BC"}>>> = convection
value<<<{"description": "For a Laplacian problem, the value of the gradient dotted with the normals on the boundary."}>>> = 1
[]
[]
[Materials<<<{"href": "../../../syntax/Materials/index.html"}>>>]
[conductivity]
type = ADGenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'k'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '1'
[]
[convection]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../../source/materials/ParsedMaterial.html"}>>>
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = 'h * (100 - T)'
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'T'
constant_names<<<{"description": "Vector of constants used in the parsed function (use this for kB etc.)"}>>> = 'h'
constant_expressions<<<{"description": "Vector of values for the constants in constant_names (can be an FParser expression)"}>>> = '${fparse 10 / (pi * inner_radius^3)}'
property_name<<<{"description": "Name of the parsed material property"}>>> = convection
[]
[]
[Executioner<<<{"href": "../../../syntax/Executioner/index.html"}>>>]
type = Steady
solve_type = NEWTON
line_search = none
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)Finally, postprocessors are defined that represent the objective (Tmax
) and constraint (volume
):
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[Tmax]
type = NodalExtremeValue<<<{"description": "Finds either the min or max elemental value of a variable over the domain.", "href": "../../../source/postprocessors/NodalExtremeValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on"}>>> = T
[]
[volume]
type = VolumePostprocessor<<<{"description": "Computes the volume of a specified block", "href": "../../../source/postprocessors/VolumePostprocessor.html"}>>>
[]
[]
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)Mesh Displacement-Based Input
This input utilizes mesh-displacement to define the ring domain. This strategy is arguably more complex; however, it enables the use of STM's "batch-restore" functionality to prevent re-initialization of the sub-application and to re-use previous results as initial guesses.
Instead of input variables, the inner radius and thickness are defined as ConstantPostprocessors, whose "value" parameter is controllable:
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[inner_radius]
type = ConstantPostprocessor<<<{"description": "Postprocessor that holds a constant value", "href": "../../../source/postprocessors/ConstantPostprocessor.html"}>>>
value<<<{"description": "The value"}>>> = 6
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_begin'
force_preaux<<<{"description": "Forces the UserObject to be executed in PREAUX"}>>> = true
[]
[thickness]
type = ConstantPostprocessor<<<{"description": "Postprocessor that holds a constant value", "href": "../../../source/postprocessors/ConstantPostprocessor.html"}>>>
value<<<{"description": "The value"}>>> = 4
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_begin'
force_preaux<<<{"description": "Forces the UserObject to be executed in PREAUX"}>>> = true
[]
[]
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus_displaced_mesh.i)The base mesh is fixed with an inner radius and thickness of 1:
@@ -1,16 +1,16 @@
[Mesh<<<{"href": "../../../syntax/Mesh/index.html"}>>>]
[disk]
type = ConcentricCircleMeshGenerator<<<{"description": "This ConcentricCircleMeshGenerator source code is to generate concentric circle meshes.", "href": "../../../source/meshgenerators/ConcentricCircleMeshGenerator.html"}>>>
has_outer_square<<<{"description": "It determines if meshes for a outer square are added to concentric circle meshes."}>>> = false
preserve_volumes<<<{"description": "Volume of concentric circles can be preserved using this function."}>>> = false
- radii<<<{"description": "Radii of major concentric circles"}>>> = '${inner_radius} ${fparse inner_radius + thickness}'
+ radii<<<{"description": "Radii of major concentric circles"}>>> = '1 2'
rings<<<{"description": "Number of rings in each circle or in the enclosing square"}>>> = '16 16'
num_sectors<<<{"description": "num_sectors % 2 = 0, num_sectors > 0Number of azimuthal sectors in each quadrant'num_sectors' must be an even number."}>>> = 16
[]
[ring]
type = BlockDeletionGenerator<<<{"description": "Mesh generator which removes elements from the specified subdomains", "href": "../../../source/meshgenerators/BlockDeletionGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = disk
block<<<{"description": "The list of blocks to be deleted"}>>> = 1
new_boundary<<<{"description": "optional boundary name to assign to the cut surface"}>>> = inner
[]
[]
(- modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)(+ modules/stochastic_tools/examples/optimization/annulus_shape/annulus_displaced_mesh.i)
Displacement variables are defined, whose values depend on the postprocessors:
[GlobalParams<<<{"href": "../../../syntax/GlobalParams/index.html"}>>>]
displacements = 'disp_x disp_y'
[]
[AuxVariables<<<{"href": "../../../syntax/AuxVariables/index.html"}>>>]
[disp_x]
[]
[disp_y]
[]
[]
[AuxKernels<<<{"href": "../../../syntax/AuxKernels/index.html"}>>>]
[disp_x_aux]
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../../source/auxkernels/ParsedAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = disp_x
expression<<<{"description": "Parsed function expression to compute"}>>> = 'r:=sqrt(x * x + y * y);theta:=atan(y / x);
dr:=r * (thickness - 1) + inner_radius - thickness;
dr * cos(theta)'
functor_names<<<{"description": "Functors to use in the parsed expression. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'inner_radius thickness'
use_xyzt<<<{"description": "Make coordinate (x,y,z) and time (t) variables available in the function expression."}>>> = true
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_begin'
use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
[]
[disp_y_aux]
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../../source/auxkernels/ParsedAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = disp_y
expression<<<{"description": "Parsed function expression to compute"}>>> = 'r:=sqrt(x * x + y * y);theta:=atan(y / x);
dr:=r * (thickness - 1) + inner_radius - thickness;
dr * sin(theta)'
functor_names<<<{"description": "Functors to use in the parsed expression. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'inner_radius thickness'
use_xyzt<<<{"description": "Make coordinate (x,y,z) and time (t) variables available in the function expression."}>>> = true
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_begin'
use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
[]
[]
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus_displaced_mesh.i)Note that the heat transfer coefficient also depends on the inner radius, so the convection boundary condition material needs adjustment:
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[h]
type = ParsedPostprocessor<<<{"description": "Computes a parsed expression with post-processors", "href": "../../../source/postprocessors/ParsedPostprocessor.html"}>>>
expression<<<{"description": "function expression"}>>> = '${fparse 10 / pi} / inner_radius^3'
pp_names<<<{"description": "Post-processors arguments"}>>> = 'inner_radius'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_begin'
[]
[]
(modules/stochastic_tools/examples/optimization/annulus_shape/annulus_displaced_mesh.i) @@ -1,10 +1,10 @@
[Materials<<<{"href": "../../../syntax/Materials/index.html"}>>>]
[convection]
type = ADParsedMaterial<<<{"description": "Parsed expression Material.", "href": "../../../source/materials/ParsedMaterial.html"}>>>
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = 'h * (100 - T)'
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'T'
- constant_names<<<{"description": "Vector of constants used in the parsed function (use this for kB etc.)"}>>> = 'h'
- constant_expressions<<<{"description": "Vector of values for the constants in constant_names (can be an FParser expression)"}>>> = '${fparse 10 / (pi * inner_radius^3)}'
+ postprocessor_names<<<{"description": "Vector of postprocessor names used in the parsed function"}>>> = 'h'
property_name<<<{"description": "Name of the parsed material property"}>>> = convection
+ use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
[]
[]
(- modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)(+ modules/stochastic_tools/examples/optimization/annulus_shape/annulus_displaced_mesh.i)
Finally, we have to ensure the rest of the heat conduction objects and postprocessors must use the displaced mesh:
@@ -1,32 +1,37 @@
[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
[diffusion]
type = ADMatDiffusion<<<{"description": "Diffusion equation kernel that takes an isotropic diffusivity from a material property", "href": "../../../source/kernels/ADMatDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
diffusivity<<<{"description": "The diffusivity value or material property"}>>> = k
+ use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
[]
[src]
type = ADBodyForce<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../../../source/kernels/BodyForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
value<<<{"description": "Coefficient to multiply by the body force term"}>>> = 1
+ use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
[]
[]
[BCs<<<{"href": "../../../syntax/BCs/index.html"}>>>]
[convection]
type = ADMatNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{C \\partial u}{\\partial n}=M*h$, where $h$ is a constant, $M$ is a material property, and $C$ is a coefficient defined by the kernel for $u$.", "href": "../../../source/bcs/MatNeumannBC.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = inner
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
boundary_material<<<{"description": "Material property multiplying the constant that will be enforced by the BC"}>>> = convection
value<<<{"description": "For a Laplacian problem, the value of the gradient dotted with the normals on the boundary."}>>> = 1
+ use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
[]
[]
+
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[Tmax]
type = NodalExtremeValue<<<{"description": "Finds either the min or max elemental value of a variable over the domain.", "href": "../../../source/postprocessors/NodalExtremeValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on"}>>> = T
[]
[volume]
type = VolumePostprocessor<<<{"description": "Computes the volume of a specified block", "href": "../../../source/postprocessors/VolumePostprocessor.html"}>>>
+ use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
[]
[]
(- modules/stochastic_tools/examples/optimization/annulus_shape/annulus.i)(+ modules/stochastic_tools/examples/optimization/annulus_shape/annulus_displaced_mesh.i)
Optimization Using SciPy
The following presents an example python script that leverages various SciPy optimizers to solve this optimization problem using StochasticControl to drive the simulation. A description of the main function of this script is below:
def optimize_annulus(
input_file: str,
executable: str,
num_procs: int,
volume: float,
mode: int,
optimizer: str,
cli_args: list[str] = [],
) -> OptimizeResult:
"""Perform shape optimization of annular pipe.
Parameters:
input_file (str): Physics input file (options: 'annulus.i' and 'annulus_displaced_mesh.i')
executable (str): Executable with stochastic-tools module included.
num_procs (int): Number of processors to execute the optimization.
volume (float): Volume of annulus to keep constant.
mode (int): MultiApp mode of execution (See moose_stochastic_tools.StochasticRunOptions.MultiAppMode)
optimizer (str): SciPy optimizer to use, e.g. 'shgo', 'dual_annealing', 'slsqp', etc.
cli_args (list[str]): Extra command-line arguments for stochastic input.
Returns
OptimizeResults: Results of the optimization.
"""
(modules/stochastic_tools/examples/optimization/annulus_shape/optimize_annulus.py)Configuring the StochasticControl
The parameters to control are dependent on which input file is being utilized: Next, the postprocessors defining the objective and constraints of the optimization are specified as QoIs of the simulation: A StochasticRunOptions
is then built, to be given as options for the StochasticControl
. The StochasticControl
is created as a context manager which returns a StochasticRunner
object. Since the optimizers often run objective and constraint evaluation separately, but with the same sample, it is recommended to cache the input-output pairs using configCache
. Finally, somewhat arbitrary bounds for the parameters are defined.
# Parameters to optimize depend on which input we're using
if input_file.endswith("annulus.i"):
params = ["inner_radius", "thickness"]
elif input_file.endswith("annulus_displaced_mesh.i"):
params = ["Postprocessors/inner_radius/value", "Postprocessors/thickness/value"]
else:
raise ValueError(f"Unknown input file {input_file}")
# QoIs include the minimization function (max temperature) and constraint (volume)
qois = ["Tmax/value", "volume/value"]
# Build StochasticControl with options
opts = StochasticRunOptions(
num_procs=num_procs,
multiapp_mode=StochasticRunOptions.MultiAppMode(mode),
cli_args=cli_args,
)
with StochasticControl(executable, input_file, params, qois, opts) as runner:
# This will ensure that if the same input is seen twice, it will only run the simulation once.
runner.configCache()
# Somewhat arbitrary bounds to the inner radius and thickness
bounds = [(0.1, 10), (0.01, 10)]
(modules/stochastic_tools/examples/optimization/annulus_shape/optimize_annulus.py)Helper functions (QoI adapters)
def computeTmax(runner: StochasticRunner, x: np.ndarray) -> np.ndarray | float:
"""Use the inputted runner to evaluate the inputted sample and
retrieve the appropriate index for max temperature of the result.
"""
y = runner(x)
return y[:, 0] if y.ndim > 1 else y[0]
def computeVolume(runner: StochasticRunner, x: np.ndarray) -> np.ndarray | float:
"""Use the inputted runner to evaluate the inputted sample and
retrieve the appropriate index for volume of the result.
"""
y = runner(x)
return y[:, 1] if y.ndim > 1 else y[1]
(modules/stochastic_tools/examples/optimization/annulus_shape/optimize_annulus.py)These wrap the runner and select the correct QoI column.
Optimizer behaviors
1) SHGO (global)
# Create lambdas so that 'runner' is not a required argument for the optimizer
objective = lambda x: computeTmax(runner, x)
eq_constraint = lambda x: computeVolume(runner, x) - volume
result = shgo(
objective,
bounds=bounds,
constraints=[{"type": "eq", "fun": eq_constraint}],
workers=1 if num_procs == 1 else runner.parallelWorker,
)
(modules/stochastic_tools/examples/optimization/annulus_shape/optimize_annulus.py)Uses equality constraint via
constraints=[...]
.Sets
workers
torunner.parallelWorker
whennum_procs > 1
to distribute evaluations.
2) Differential Evolution (global, vectorized)
# Differential evolution allows for "vectorization", sampling multiple rows at a time
# The input it provides (and output expected) is actually transposed to what
# the StochasticRunner provides.
# Create lambdas so that 'runner' is not a required argument for the optimizer
objective = lambda x: (
computeTmax(runner, x.T) if x.ndim == 2 else computeTmax(runner, x)
)
eq_constraint = lambda x: (
computeVolume(runner, x.T).reshape((1, -1))
if x.ndim == 2
else computeVolume(runner, x)
)
result = differential_evolution(
objective,
bounds=bounds,
constraints=NonlinearConstraint(eq_constraint, volume, volume),
vectorized=True,
maxiter=10,
)
(modules/stochastic_tools/examples/optimization/annulus_shape/optimize_annulus.py)Vectorized mode samples many candidates at once; SciPy passes a matrix shaped
(n_params, n_points)
. The code transposes to matchrunner
’s(n_points, n_params)
expectation.Equality constraint via
NonlinearConstraint
.maxiter=10
for demonstration; increase for higher‑quality solutions.
3) Local methods via minimize
# Create lambdas so that 'runner' is not a required argument for the optimizer
objective = lambda x: computeTmax(runner, x)
eq_constraint = lambda x: computeVolume(runner, x)
result = minimize(
objective,
x0=np.array([6, 4]),
method=optimizer.upper(),
bounds=bounds,
constraints=NonlinearConstraint(eq_constraint, volume, volume),
options={"workers": 1 if num_procs == 1 else runner.parallelWorker},
)
(modules/stochastic_tools/examples/optimization/annulus_shape/optimize_annulus.py)Starts from
x0=[6, 4]
.Adds equality constraint with
NonlinearConstraint
.Parallel evaluation available through
runner.parallelWorker
in methods that respect theworkers
option.
Running the Example
This example script can be run using the following command-line arguments:
python optimize_annulus.py <args>
Flag | Default | Description |
---|---|---|
-i, --input-file | annulus.i | Use annulus.i or annulus_displaced_mesh.i |
-e, --executable | stochastic_tools-opt | MOOSE executable with stochastic_tools |
-n, --num-procs | 1 | Number of processors for the stochastic run |
-v, --volume | 200.0 | Equality constraint value passed to the optimizer |
-m, --mode | 1 | Integer mapped to StochasticRunOptions.MultiAppMode |
-o, --optimizer | shgo | One of shgo , differential_evolution , or any method supported by scipy.optimize.minimize (e.g., SLSQP , COBYLA , etc.) |
The resulting inner radius is approximately 1.198, the thickness is 6.877, and the max temperature is 160.3.