Using Stochastic Tools with Multiphysics Models
The purpose of this document is to present a multiphysics example using the Stochastic Tools Module. The intention is to showcase the capabilities of the module to produce statistically relevant results including uncertainty propagation and sensitivity, as well as the module's surrogate modeling infrastructure. The problem of interest is a thermomechanics model using a combination of the Heat Transfer and Solid Mechanics modules. The problem involves multiple uncertain material properties and multiple quantities of interest (QoI). Using both Monte Carlo sampling and polynomial chaos surrogate modeling, the effect of these properties' uncertainties are quantified with uncertainty propagation and global sensitivity analysis.
Problem Description
The problem of interest involves a steady-state thermomechanics model. The geometry is a 3-D finite hollow cylinder with two concentric layers of different material properties, seen in Figure 1. Due to symmetry, only 1/8 of the cylinder is represented by the mesh. The inner surface of the ring is exposed to a surface heat source that has a cosine shape along the axis of the cylinder with its peak being at the center. The end and outside of the cylinder have convective boundary conditions. The cylinder is free to displace in all directions due to the thermal expansion. The relevant material and geometric properties are listed in Table 1. This table also lists the "uncertain" parameters that will be described in the following section. For reference, the temperature and displacement profiles are shown in Figure 4 to Figure 8 where the uncertain properties are set to some arbitrary values. The full input file for this model is shown in Listing 1.

Figure 1: Model Problem Geometry
Table 1: Material Properties for Thermomechanics Cylinder
Property | Symbol | Value | Units |
---|---|---|---|
Half Cylinder Height | 3 | m | |
Inner Radius | 1.0 | m | |
Inner Ring Width | 0.1 | m | |
Outer Ring Width | 0.1 | m | |
Outer Heat Transfer Coef. | 10 | W/mK | |
Outer Free Temperature | 300 | K | |
End Heat Transfer Coef. | 10 | W/mK | |
End Free Temperature | 300 | K | |
Heat Source Magnitude | Uncertain | W | |
Inner Thermal Conductivity | Uncertain | W/mK | |
Outer Thermal Conductivity | Uncertain | W/mK | |
Inner Young's Modulus | Uncertain | Pa | |
Outer Young's Modulus | Uncertain | Pa | |
Inner Poisson's Ratio | Uncertain | ||
Outer Poisson's Ratio | Uncertain | ||
Inner Thermal Expansion Coef. | Uncertain | 1/K | |
Outer Thermal Expansion Coef. | Uncertain | 1/K | |
At Rest Temperature | 300 | K |

Figure 4: Temperature (K)

Figure 5: Displacement (m)

Figure 6: X-Displacement (m)

Figure 7: Y-Displacement (m)

Figure 8: Z-Displacement (m)
Listing 1: Thermomechanics model input file
# Generate 1/4 of a 2-ring disk and extrude it by half to obtain
# 1/8 of a 3D tube. Mirror boundary conditions will exist on the
# cut portions.
[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"}>>>
num_sectors<<<{"description": "num_sectors % 2 = 0, num_sectors > 0Number of azimuthal sectors in each quadrant'num_sectors' must be an even number."}>>> = 10
radii<<<{"description": "Radii of major concentric circles"}>>> = '1.0 1.1 1.2'
rings<<<{"description": "Number of rings in each circle or in the enclosing square"}>>> = '1 1 1'
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
portion<<<{"description": "Control of which part of mesh is created"}>>> = top_right
[]
[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'
[]
[cylinder]
type = MeshExtruderGenerator<<<{"description": "Takes a 1D or 2D mesh and extrudes the entire structure along the specified axis increasing the dimensionality of the mesh.", "href": "../../../source/meshgenerators/MeshExtruderGenerator.html"}>>>
input<<<{"description": "the mesh we want to extrude"}>>> = ring
extrusion_vector<<<{"description": "The direction and length of the extrusion"}>>> = '0 0 1.5'
num_layers<<<{"description": "The number of layers in the extruded mesh"}>>> = 15
bottom_sideset<<<{"description": "The boundary that will be applied to the bottom of the extruded mesh"}>>> = 'back'
top_sideset<<<{"description": "The boundary that will be to the top of the extruded mesh"}>>> = 'front'
[]
[]
[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
[T]
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 300
[]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[]
[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
[hc]
type = HeatConduction<<<{"description": "Diffusive heat conduction term $-\\nabla\\cdot(k\\nabla T)$ of the thermal energy conservation equation", "href": "../../../source/kernels/HeatConduction.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
[]
[TensorMechanics<<<{"href": "../../../syntax/Kernels/TensorMechanics/index.html"}>>>]
displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
[]
[]
[BCs<<<{"href": "../../../syntax/BCs/index.html"}>>>]
[temp_inner]
type = FunctionNeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, where $h$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../../../source/bcs/FunctionNeumannBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'inner'
function<<<{"description": "The function."}>>> = surface_source
[]
[temp_front]
type = ConvectiveHeatFluxBC<<<{"description": "Convective heat transfer boundary condition with temperature and heat transfer coefficent given by material properties.", "href": "../../../source/bcs/ConvectiveHeatFluxBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
T_infinity<<<{"description": "Material property for far-field temperature"}>>> = 300
heat_transfer_coefficient<<<{"description": "Material property for heat transfer coefficient"}>>> = 10
[]
[temp_outer]
type = ConvectiveHeatFluxBC<<<{"description": "Convective heat transfer boundary condition with temperature and heat transfer coefficent given by material properties.", "href": "../../../source/bcs/ConvectiveHeatFluxBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'outer'
T_infinity<<<{"description": "Material property for far-field temperature"}>>> = 300
heat_transfer_coefficient<<<{"description": "Material property for heat transfer coefficient"}>>> = 10
[]
# mirror boundary conditions.
[disp_x]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[disp_y]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[disp_z]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[]
[Materials<<<{"href": "../../../syntax/Materials/index.html"}>>>]
[cond_inner]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../../source/materials/GenericConstantMaterial.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
prop_names<<<{"description": "The names of the properties this material will have"}>>> = thermal_conductivity
prop_values<<<{"description": "The values associated with the named properties"}>>> = 25
[]
[cond_outer]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../../source/materials/GenericConstantMaterial.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
prop_names<<<{"description": "The names of the properties this material will have"}>>> = thermal_conductivity
prop_values<<<{"description": "The values associated with the named properties"}>>> = 100
[]
[elasticity_tensor_inner]
type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../../../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 2.1e5
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.3
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[elasticity_tensor_outer]
type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../../../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 3.1e5
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.2
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[thermal_strain_inner]
type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../../../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 2e-6
temperature<<<{"description": "Coupled temperature"}>>> = T
stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 300
eigenstrain_name<<<{"description": "Material property name for the eigenstrain tensor computed by this model. IMPORTANT: The name of this property must also be provided to the strain calculator."}>>> = eigenstrain_inner
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[thermal_strain_outer]
type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../../../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 1e-6
temperature<<<{"description": "Coupled temperature"}>>> = T
stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 300
eigenstrain_name<<<{"description": "Material property name for the eigenstrain tensor computed by this model. IMPORTANT: The name of this property must also be provided to the strain calculator."}>>> = eigenstrain_outer
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[strain_inner] #We use small deformation mechanics
type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../../../source/materials/ComputeSmallStrain.html"}>>>
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'eigenstrain_inner'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 2
[]
[strain_outer] #We use small deformation mechanics
type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../../../source/materials/ComputeSmallStrain.html"}>>>
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'eigenstrain_outer'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 3
[]
[stress] #We use linear elasticity
type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../../../source/materials/ComputeLinearElasticStress.html"}>>>
[]
[]
[Functions<<<{"href": "../../../syntax/Functions/index.html"}>>>]
[surface_source]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = 'Q_t*pi/2.0/3.0 * cos(pi/3.0*z)'
symbol_names<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = 'Q_t'
symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = heat_source
[]
[]
[Executioner<<<{"href": "../../../syntax/Executioner/index.html"}>>>]
type = Steady
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart'
petsc_options_value = 'hypre boomeramg 101'
l_max_its = 30
nl_max_its = 100
nl_abs_tol = 1e-9
l_tol = 1e-04
[]
[Preconditioning<<<{"href": "../../../syntax/Preconditioning/index.html"}>>>]
[SMP]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
[]
[]
[Controls<<<{"href": "../../../syntax/Controls/index.html"}>>>]
[stochastic]
type = SamplerReceiver<<<{"description": "Control for receiving data from a Sampler via SamplerParameterTransfer.", "href": "../../../source/controls/SamplerReceiver.html"}>>>
[]
[]
[VectorPostprocessors<<<{"href": "../../../syntax/VectorPostprocessors/index.html"}>>>]
[temp_center]
type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = T
start_point<<<{"description": "The beginning of the line"}>>> = '1 0 0'
end_point<<<{"description": "The ending of the line"}>>> = '1.2 0 0'
num_points<<<{"description": "The number of points to sample along the line"}>>> = 11
sort_by<<<{"description": "What to sort the samples by"}>>> = 'x'
[]
[temp_end]
type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = T
start_point<<<{"description": "The beginning of the line"}>>> = '1 0 1.5'
end_point<<<{"description": "The ending of the line"}>>> = '1.2 0 1.5'
num_points<<<{"description": "The number of points to sample along the line"}>>> = 11
sort_by<<<{"description": "What to sort the samples by"}>>> = 'x'
[]
[dispx_center]
type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = disp_x
start_point<<<{"description": "The beginning of the line"}>>> = '1 0 0'
end_point<<<{"description": "The ending of the line"}>>> = '1.2 0 0'
num_points<<<{"description": "The number of points to sample along the line"}>>> = 11
sort_by<<<{"description": "What to sort the samples by"}>>> = 'x'
[]
[dispx_end]
type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = disp_x
start_point<<<{"description": "The beginning of the line"}>>> = '1 0 1.5'
end_point<<<{"description": "The ending of the line"}>>> = '1.2 0 1.5'
num_points<<<{"description": "The number of points to sample along the line"}>>> = 11
sort_by<<<{"description": "What to sort the samples by"}>>> = 'x'
[]
[dispz_end]
type = LineValueSampler<<<{"description": "Samples variable(s) along a specified line", "href": "../../../source/vectorpostprocessors/LineValueSampler.html"}>>>
variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = disp_z
start_point<<<{"description": "The beginning of the line"}>>> = '1 0 1.5'
end_point<<<{"description": "The ending of the line"}>>> = '1.2 0 1.5'
num_points<<<{"description": "The number of points to sample along the line"}>>> = 11
sort_by<<<{"description": "What to sort the samples by"}>>> = 'x'
[]
[]
[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
[heat_source]
type = FunctionValuePostprocessor<<<{"description": "Computes the value of a supplied function at a single point (scalable)", "href": "../../../source/postprocessors/FunctionValuePostprocessor.html"}>>>
function<<<{"description": "The function which supplies the postprocessor value."}>>> = 1
scale_factor<<<{"description": "A scale factor to be applied to the function"}>>> = 10000
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."}>>> = linear
[]
[temp_center_inner]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = T
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0 0'
[]
[temp_center_outer]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = T
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1.2 0 0'
[]
[temp_end_inner]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = T
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0 1.5'
[]
[temp_end_outer]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = T
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1.2 0 1.5'
[]
[dispx_center_inner]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_x
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0 0'
[]
[dispx_center_outer]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_x
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1.2 0 0'
[]
[dispx_end_inner]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_x
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0 1.5'
[]
[dispx_end_outer]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_x
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1.2 0 1.5'
[]
[dispz_inner]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_z
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0 1.5'
[]
[dispz_outer]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../../../source/postprocessors/PointValue.html"}>>>
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_z
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1.2 0 1.5'
[]
[]
[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = false
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = false
[]
(modules/combined/examples/stochastic/thermomech/graphite_ring_thermomechanics.i)Uncertain Parameters
A total of nine properties of the model are not known exactly, but with some known probability of values occurring. The probabilities are represented by each parameter's probability density function. All parameters have a independent uniform distribution, , where and are the lower and upper limit values for the property, respectively. Table 2 lists the details of each of the property's distribution.
Table 2: Uniform distribution parameters for uncertain properties
Index | Property | ||
---|---|---|---|
1 | 20 | 30 | |
2 | 90 | 110 | |
3 | 9000 | 11000 | |
4 | 1.0 | 3.0 | |
5 | 0.5 | 1.5 | |
6 | 2.0 | 2.2 | |
7 | 3.0 | 3.2 | |
8 | 0.29 | 0.31 | |
9 | 0.19 | 0.21 |
Quantities of Interest
There are a total of ten QoIs for the model, which involve temperature and displacement:
Temperature at center of inner surface —
Temperature at center of outer surface —
Temperature at end of inner surface —
Temperature at end of outer surface —
x-displacement at center of inner surface —
x-displacement at center of outer surface —
x-displacement at end of inner surface —
x-displacement at end of outer surface —
z-displacement at end of inner surface —
z-displacement at end of outer surface —
Figure 2 shows geometrically where these QoIs are located in the model.

Figure 2: Geometric description of quantities of interest
Results
In this exercise, we will use the statistics and Sobol sensitivity capabilities available in the stochastic tools module. The goal of this exercise is to understand how the uncertainty in the parameters affects the resulting QoIs. This is done through sampling the model at different perturbations of the parameters and performing statistical calculations on resulting QoI values. Two methods are used to perform this analysis. First is using the sampler system to perturb the uncertain properties and retrieve the QoIs which will undergo the analysis. The second is training a polynomial chaos surrogate and using that reduced order model to sample and perform the analysis. The idea is that many evaluations of the model are necessary to compute accurate statistical quantities and surrogate modeling speeds up this computation by requiring much fewer full model evaluations for training and is significantly faster to evaluate once trained.
Using latin hypercube sampling, the thermomechanics model was run with a total of 100,000 samples, the input file is shown by Listing 2. A order four polynomial chaos surrogate was training using a Smolyak sparse quadrature for a total of 7,344 runs of the full model. The training input is shown by Listing 3 and the evaluation input is shown by Listing 4. Table 3 shows the run-time for sampling the full order model and training and evaluating the surrogate. We see here that cumulative time for training and evaluating the surrogate is much smaller than just sampling the full order model, this is because building the surrogate required far fewer evaluations of the full model and evaluating the surrogate is much faster than evaluating the full model.
Table 3: Stochastic run-time results
Simulation | Samples | CPU Time |
---|---|---|
Full-Order Sampling | 100,000 | 176 hr |
Polynomial Chaos — Training | 7,344 | 13.7 hr |
Polynomial Chaos — Evaluation | 100,000 | 6.8 s |
Listing 2: Latin hypercube sampling and statistics input file
[StochasticTools<<<{"href": "../../../syntax/StochasticTools/index.html"}>>>]
[]
[Distributions<<<{"href": "../../../syntax/Distributions/index.html"}>>>]
[cond_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 20
upper_bound<<<{"description": "Distribution upper bound"}>>> = 30
[]
[cond_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 90
upper_bound<<<{"description": "Distribution upper bound"}>>> = 110
[]
[heat_source]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 9000
upper_bound<<<{"description": "Distribution upper bound"}>>> = 11000
[]
[alpha_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 1e-6
upper_bound<<<{"description": "Distribution upper bound"}>>> = 3e-6
[]
[alpha_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 5e-7
upper_bound<<<{"description": "Distribution upper bound"}>>> = 1.5e-6
[]
[ymod_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 2e5
upper_bound<<<{"description": "Distribution upper bound"}>>> = 2.2e5
[]
[ymod_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 3e5
upper_bound<<<{"description": "Distribution upper bound"}>>> = 3.2e5
[]
[prat_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 0.29
upper_bound<<<{"description": "Distribution upper bound"}>>> = 0.31
[]
[prat_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 0.19
upper_bound<<<{"description": "Distribution upper bound"}>>> = 0.21
[]
[]
[Samplers<<<{"href": "../../../syntax/Samplers/index.html"}>>>]
[sample]
type = LatinHypercube<<<{"description": "Latin Hypercube Sampler.", "href": "../../../source/samplers/LatinHypercubeSampler.html"}>>>
num_rows<<<{"description": "The size of the square matrix to generate."}>>> = 100000
distributions<<<{"description": "The distribution names to be sampled, the number of distributions provided defines the number of columns per matrix."}>>> = 'cond_inner cond_outer heat_source alpha_inner alpha_outer ymod_inner ymod_outer prat_inner prat_outer'
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."}>>> = INITIAL
[]
[]
[MultiApps<<<{"href": "../../../syntax/MultiApps/index.html"}>>>]
[sub]
type = SamplerFullSolveMultiApp<<<{"description": "Creates a full-solve type sub-application for each row of each Sampler matrix.", "href": "../../../source/multiapps/SamplerFullSolveMultiApp.html"}>>>
input_files<<<{"description": "The input file for each App. If this parameter only contains one input file it will be used for all of the Apps. When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = graphite_ring_thermomechanics.i
sampler<<<{"description": "The Sampler object to utilize for creating MultiApps."}>>> = sample
mode<<<{"description": "The operation mode, 'normal' creates one sub-application for each row in the Sampler and 'batch-reset' and 'batch-restore' creates N sub-applications, where N is the minimum of 'num_rows' in the Sampler and floor(number of processes / min_procs_per_app). To run the rows in the Sampler, 'batch-reset' will destroy and re-create sub-apps as needed, whereas the 'batch-restore' will backup and restore sub-apps to the initial state prior to execution, without destruction."}>>> = batch-reset
[]
[]
[Transfers<<<{"href": "../../../syntax/Transfers/index.html"}>>>]
[sub]
type = SamplerParameterTransfer<<<{"description": "Copies Sampler data to a SamplerReceiver object.", "href": "../../../source/transfers/SamplerParameterTransfer.html"}>>>
to_multi_app<<<{"description": "The name of the MultiApp to transfer the data to"}>>> = sub
sampler<<<{"description": "A the Sampler object that Transfer is associated.."}>>> = sample
parameters<<<{"description": "A list of parameters (on the sub application) to control with the Sampler data. The order of the parameters listed here should match the order of the items in the Sampler."}>>> = 'Materials/cond_inner/prop_values Materials/cond_outer/prop_values
Postprocessors/heat_source/scale_factor
Materials/thermal_strain_inner/thermal_expansion_coeff Materials/thermal_strain_outer/thermal_expansion_coeff
Materials/elasticity_tensor_inner/youngs_modulus Materials/elasticity_tensor_outer/youngs_modulus
Materials/elasticity_tensor_inner/poissons_ratio Materials/elasticity_tensor_outer/poissons_ratio'
check_multiapp_execute_on<<<{"description": "When false the check between the multiapp and transfer execute on flags is not performed."}>>> = false
[]
[data]
type = SamplerReporterTransfer<<<{"description": "Transfers data from Reporters on the sub-application to a StochasticReporter on the main application.", "href": "../../../source/transfers/SamplerReporterTransfer.html"}>>>
from_multi_app<<<{"description": "The name of the MultiApp to receive data from"}>>> = sub
sampler<<<{"description": "A the Sampler object that Transfer is associated.."}>>> = sample
stochastic_reporter<<<{"description": "The name of the StochasticReporter object to transfer values to."}>>> = storage
from_reporter<<<{"description": "The name(s) of the Reporter(s) on the sub-app to transfer from."}>>> = 'temp_center_inner/value temp_center_outer/value temp_end_inner/value temp_end_outer/value
dispx_center_inner/value dispx_center_outer/value dispx_end_inner/value dispx_end_outer/value
dispz_inner/value dispz_outer/value'
[]
[]
[Reporters<<<{"href": "../../../syntax/Reporters/index.html"}>>>]
[storage]
type = StochasticReporter<<<{"description": "Storage container for stochastic simulation results coming from Reporters.", "href": "../../../source/reporters/StochasticReporter.html"}>>>
parallel_type<<<{"description": "This parameter will determine how the stochastic data is gathered. It is common for outputting purposes that this parameter be set to ROOT, otherwise, many files will be produced showing the values on each processor. However, if there are lot of samples, gathering on root may be memory restrictive."}>>> = ROOT
[]
[stats]
type = StatisticsReporter<<<{"description": "Compute statistical values of a given VectorPostprocessor objects and vectors.", "href": "../../../source/reporters/StatisticsReporter.html"}>>>
reporters<<<{"description": "List of Reporter values to utilized for statistic computations."}>>> = 'storage/data:temp_center_inner:value storage/data:temp_center_outer:value storage/data:temp_end_inner:value storage/data:temp_end_outer:value
storage/data:dispx_center_inner:value storage/data:dispx_center_outer:value storage/data:dispx_end_inner:value storage/data:dispx_end_outer:value
storage/data:dispz_inner:value storage/data:dispz_outer:value'
compute<<<{"description": "The statistic(s) to compute for each of the supplied vector postprocessors."}>>> = 'mean stddev'
ci_method<<<{"description": "The method to use for computing confidence level intervals."}>>> = 'percentile'
ci_levels<<<{"description": "A vector of confidence levels to consider, values must be in (0, 1)."}>>> = '0.05 0.95'
[]
[]
[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
[out]
type = JSON<<<{"description": "Output for Reporter values using JSON format.", "href": "../../../source/outputs/JSONOutput.html"}>>>
[]
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_END
[]
(modules/combined/examples/stochastic/thermomech/lhs_uniform.i)Listing 3: Polynomial chaos training input file
[StochasticTools<<<{"href": "../../../syntax/StochasticTools/index.html"}>>>]
[]
[Distributions<<<{"href": "../../../syntax/Distributions/index.html"}>>>]
[cond_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 20
upper_bound<<<{"description": "Distribution upper bound"}>>> = 30
[]
[cond_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 90
upper_bound<<<{"description": "Distribution upper bound"}>>> = 110
[]
[heat_source]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 9000
upper_bound<<<{"description": "Distribution upper bound"}>>> = 11000
[]
[alpha_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 1e-6
upper_bound<<<{"description": "Distribution upper bound"}>>> = 3e-6
[]
[alpha_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 5e-7
upper_bound<<<{"description": "Distribution upper bound"}>>> = 1.5e-6
[]
[ymod_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 2e5
upper_bound<<<{"description": "Distribution upper bound"}>>> = 2.2e5
[]
[ymod_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 3e5
upper_bound<<<{"description": "Distribution upper bound"}>>> = 3.2e5
[]
[prat_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 0.29
upper_bound<<<{"description": "Distribution upper bound"}>>> = 0.31
[]
[prat_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 0.19
upper_bound<<<{"description": "Distribution upper bound"}>>> = 0.21
[]
[]
[GlobalParams<<<{"href": "../../../syntax/GlobalParams/index.html"}>>>]
distributions = 'cond_inner cond_outer heat_source alpha_inner alpha_outer ymod_inner ymod_outer prat_inner prat_outer'
[]
[Samplers<<<{"href": "../../../syntax/Samplers/index.html"}>>>]
[sample]
type = Quadrature<<<{"description": "Quadrature sampler for Polynomial Chaos.", "href": "../../../source/samplers/QuadratureSampler.html"}>>>
sparse_grid<<<{"description": "Type of sparse grid to use, if none, full tensor product is used."}>>> = smolyak
order<<<{"description": "Specify the maximum order of the polynomials in the expansion."}>>> = 5
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."}>>> = INITIAL
[]
[]
[MultiApps<<<{"href": "../../../syntax/MultiApps/index.html"}>>>]
[sub]
type = SamplerFullSolveMultiApp<<<{"description": "Creates a full-solve type sub-application for each row of each Sampler matrix.", "href": "../../../source/multiapps/SamplerFullSolveMultiApp.html"}>>>
input_files<<<{"description": "The input file for each App. If this parameter only contains one input file it will be used for all of the Apps. When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = graphite_ring_thermomechanics.i
sampler<<<{"description": "The Sampler object to utilize for creating MultiApps."}>>> = sample
mode<<<{"description": "The operation mode, 'normal' creates one sub-application for each row in the Sampler and 'batch-reset' and 'batch-restore' creates N sub-applications, where N is the minimum of 'num_rows' in the Sampler and floor(number of processes / min_procs_per_app). To run the rows in the Sampler, 'batch-reset' will destroy and re-create sub-apps as needed, whereas the 'batch-restore' will backup and restore sub-apps to the initial state prior to execution, without destruction."}>>> = batch-reset
[]
[]
[Transfers<<<{"href": "../../../syntax/Transfers/index.html"}>>>]
[sub]
type = SamplerParameterTransfer<<<{"description": "Copies Sampler data to a SamplerReceiver object.", "href": "../../../source/transfers/SamplerParameterTransfer.html"}>>>
to_multi_app<<<{"description": "The name of the MultiApp to transfer the data to"}>>> = sub
sampler<<<{"description": "A the Sampler object that Transfer is associated.."}>>> = sample
parameters<<<{"description": "A list of parameters (on the sub application) to control with the Sampler data. The order of the parameters listed here should match the order of the items in the Sampler."}>>> = 'Materials/cond_inner/prop_values Materials/cond_outer/prop_values
Postprocessors/heat_source/scale_factor
Materials/thermal_strain_inner/thermal_expansion_coeff Materials/thermal_strain_outer/thermal_expansion_coeff
Materials/elasticity_tensor_inner/youngs_modulus Materials/elasticity_tensor_outer/youngs_modulus
Materials/elasticity_tensor_inner/poissons_ratio Materials/elasticity_tensor_outer/poissons_ratio'
check_multiapp_execute_on<<<{"description": "When false the check between the multiapp and transfer execute on flags is not performed."}>>> = false
[]
[data]
type = SamplerReporterTransfer<<<{"description": "Transfers data from Reporters on the sub-application to a StochasticReporter on the main application.", "href": "../../../source/transfers/SamplerReporterTransfer.html"}>>>
from_multi_app<<<{"description": "The name of the MultiApp to receive data from"}>>> = sub
sampler<<<{"description": "A the Sampler object that Transfer is associated.."}>>> = sample
stochastic_reporter<<<{"description": "The name of the StochasticReporter object to transfer values to."}>>> = storage
from_reporter<<<{"description": "The name(s) of the Reporter(s) on the sub-app to transfer from."}>>> = 'temp_center_inner/value temp_center_outer/value temp_end_inner/value temp_end_outer/value
dispx_center_inner/value dispx_center_outer/value dispx_end_inner/value dispx_end_outer/value
dispz_inner/value dispz_outer/value'
[]
[]
[Reporters<<<{"href": "../../../syntax/Reporters/index.html"}>>>]
[storage]
type = StochasticReporter<<<{"description": "Storage container for stochastic simulation results coming from Reporters.", "href": "../../../source/reporters/StochasticReporter.html"}>>>
[]
[]
[Trainers<<<{"href": "../../../syntax/Trainers/index.html"}>>>]
[temp_center_inner]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:temp_center_inner:value
[]
[temp_center_outer]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:temp_center_outer:value
[]
[temp_end_inner]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:temp_end_inner:value
[]
[temp_end_outer]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:temp_end_outer:value
[]
[dispx_center_inner]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:dispx_center_inner:value
[]
[dispx_center_outer]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:dispx_center_outer:value
[]
[dispx_end_inner]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:dispx_end_inner:value
[]
[dispx_end_outer]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:dispx_end_outer:value
[]
[dispz_inner]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:dispz_inner:value
[]
[dispz_outer]
type = PolynomialChaosTrainer<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/trainers/PolynomialChaosTrainer.html"}>>>
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_end
order<<<{"description": "Maximum polynomial order."}>>> = 4
sampler<<<{"description": "Sampler used to create predictor and response data."}>>> = sample
response<<<{"description": "Reporter value of response results, can be vpp with <vpp_name>/<vector_name> or sampler column with 'sampler/col_<index>'."}>>> = storage/data:dispz_outer:value
[]
[]
[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
[out]
type = SurrogateTrainerOutput<<<{"description": "Output for trained surrogate model data.", "href": "../../../source/outputs/SurrogateTrainerOutput.html"}>>>
trainers<<<{"description": "A list of SurrogateTrainer objects to output."}>>> = 'temp_center_inner temp_center_outer temp_end_inner temp_end_outer
dispx_center_inner dispx_center_outer dispx_end_inner dispx_end_outer
dispz_inner dispz_outer'
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."}>>> = FINAL
[]
[]
(modules/combined/examples/stochastic/thermomech/poly_chaos_train_uniform.i)Listing 4: Polynomial chaos evaluation input file
[StochasticTools<<<{"href": "../../../syntax/StochasticTools/index.html"}>>>]
[]
[Distributions<<<{"href": "../../../syntax/Distributions/index.html"}>>>]
[cond_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 20
upper_bound<<<{"description": "Distribution upper bound"}>>> = 30
[]
[cond_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 90
upper_bound<<<{"description": "Distribution upper bound"}>>> = 110
[]
[heat_source]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 9000
upper_bound<<<{"description": "Distribution upper bound"}>>> = 11000
[]
[alpha_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 1e-6
upper_bound<<<{"description": "Distribution upper bound"}>>> = 3e-6
[]
[alpha_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 5e-7
upper_bound<<<{"description": "Distribution upper bound"}>>> = 1.5e-6
[]
[ymod_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 2e5
upper_bound<<<{"description": "Distribution upper bound"}>>> = 2.2e5
[]
[ymod_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 3e5
upper_bound<<<{"description": "Distribution upper bound"}>>> = 3.2e5
[]
[prat_inner]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 0.29
upper_bound<<<{"description": "Distribution upper bound"}>>> = 0.31
[]
[prat_outer]
type = Uniform<<<{"description": "Continuous uniform distribution.", "href": "../../../source/distributions/Uniform.html"}>>>
lower_bound<<<{"description": "Distribution lower bound"}>>> = 0.19
upper_bound<<<{"description": "Distribution upper bound"}>>> = 0.21
[]
[]
[Samplers<<<{"href": "../../../syntax/Samplers/index.html"}>>>]
[sample]
type = MonteCarlo<<<{"description": "Monte Carlo Sampler.", "href": "../../../source/samplers/MonteCarloSampler.html"}>>>
num_rows<<<{"description": "The number of rows per matrix to generate."}>>> = 100000
distributions<<<{"description": "The distribution names to be sampled, the number of distributions provided defines the number of columns per matrix."}>>> = 'cond_inner cond_outer heat_source alpha_inner alpha_outer ymod_inner ymod_outer prat_inner prat_outer'
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."}>>> = INITIAL
[]
[]
[Surrogates<<<{"href": "../../../syntax/Surrogates/index.html"}>>>]
[temp_center_inner]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_temp_center_inner.rd'
[]
[temp_center_outer]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_temp_center_outer.rd'
[]
[temp_end_inner]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_temp_end_inner.rd'
[]
[temp_end_outer]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_temp_end_outer.rd'
[]
[dispx_center_inner]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_dispx_center_inner.rd'
[]
[dispx_center_outer]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_dispx_center_outer.rd'
[]
[dispx_end_inner]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_dispx_end_inner.rd'
[]
[dispx_end_outer]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_dispx_end_outer.rd'
[]
[dispz_inner]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_dispz_inner.rd'
[]
[dispz_outer]
type = PolynomialChaos<<<{"description": "Computes and evaluates polynomial chaos surrogate model.", "href": "../../../source/surrogates/PolynomialChaos.html"}>>>
filename<<<{"description": "The name of the file which will be associated with the saved/loaded data."}>>> = 'poly_chaos_train_uniform_out_dispz_outer.rd'
[]
[]
[Reporters<<<{"href": "../../../syntax/Reporters/index.html"}>>>]
[storage]
type = EvaluateSurrogate<<<{"description": "Tool for sampling surrogate models.", "href": "../../../source/reporters/EvaluateSurrogate.html"}>>>
sampler<<<{"description": "Sampler to use for evaluating surrogate models."}>>> = sample
model<<<{"description": "Name of surrogate models."}>>> = 'temp_center_inner temp_center_outer temp_end_inner temp_end_outer
dispx_center_inner dispx_center_outer dispx_end_inner dispx_end_outer
dispz_inner dispz_outer'
parallel_type<<<{"description": "This parameter will determine how the stochastic data is gathered. It is common for outputting purposes that this parameter be set to ROOT, otherwise, many files will be produced showing the values on each processor. However, if there are lot of samples, gathering on root may be memory restrictive."}>>> = ROOT
[]
[stats]
type = PolynomialChaosReporter<<<{"description": "Tool for extracting data from PolynomialChaos surrogates and computing statistics.", "href": "../../../source/reporters/PolynomialChaosReporter.html"}>>>
pc_name<<<{"description": "Name(s) of PolynomialChaos surrogate object(s)."}>>> = 'temp_center_inner temp_center_outer temp_end_inner temp_end_outer
dispx_center_inner dispx_center_outer dispx_end_inner dispx_end_outer
dispz_inner dispz_outer'
statistics<<<{"description": "Statistics to compute."}>>> = 'mean stddev'
include_sobol<<<{"description": "True to compute Sobol indices."}>>> = true
[]
[]
[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
[out]
type = JSON<<<{"description": "Output for Reporter values using JSON format.", "href": "../../../source/outputs/JSONOutput.html"}>>>
[]
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_END
[]
(modules/combined/examples/stochastic/thermomech/poly_chaos_uniform.i)Statistics
Table 4 shows the statistical results of sampling the thermomechanics model and the polynomial chaos surrogate. and represent the mean and standard deviation of the QoI, and CI is the confidence interval. Note that the confidence interval for the polynomial chaos statistics is not relevant since these values were found analytically using integration techniques. Figure 9 to Figure 11 compares several of the probability distributions of the QoIs between sampling the full-order model and the polynomial chaos surrogate.
Table 4: Statistics Results
QoI | 95% CI | 95% CI | PC – | PC – | ||
---|---|---|---|---|---|---|
609.97 | (609.87, 610.06) | 18.22 | (18.18, 18.27) | 609.97 | 18.23 | |
586.85 | (586.77, 586.94) | 16.64 | (16.60, 16.68) | 586.85 | 16.64 | |
506.31 | (506.25, 506.37) | 12.05 | (12.03, 12.08) | 506.31 | 12.05 | |
507.92 | (507.85, 507.98) | 12.13 | (12.10, 12.15) | 507.92 | 12.12 | |
4.032E-04 | (4.028E-04, 4.037E-04) | 8.608E-05 | (8.581E-05, 8.636E-05) | 4.032E-04 | 8.625E-05 | |
4.996E-04 | (4.990E-04, 5.001E-04) | 1.078E-04 | (1.075E-04, 1.082E-04) | 4.996E-04 | 1.081E-04 | |
4.220E-04 | (4.213E-04, 4.226E-04) | 1.255E-04 | (1.252E-04, 1.258E-04) | 4.220E-04 | 1.256E-04 | |
4.793E-04 | (4.786E-04, 4.800E-04) | 1.352E-04 | (1.349E-04, 1.356E-04) | 4.794E-04 | 1.354E-04 | |
6.139E-04 | (6.131E-04, 6.146E-04) | 1.466E-04 | (1.462E-04, 1.470E-04) | 6.139E-04 | 1.469E-04 | |
4.995E-04 | (4.989E-04, 5.000E-04) | 1.067E-04 | (1.065E-04, 1.072E-04) | 4.995E-04 | 1.071E-04 |

Figure 9: Histogram

Figure 10: Histogram

Figure 11: Histogram
Sobol Sensitivities
Sobol sensitivities, or Sobol indices, are a metric to compare the global sensitivity a parameter has on a QoI. This examples demonstrates several different types of the sensitivities. The first is total sensitivity, which measure the total sensitivity from a parameter, Figure 3 shows these values for each QoI and parameter. The second is a correlative sensitivity, which measures the sensitivity due to a combination of parameters, Figure 12 to Figure 14 show heat maps of these values for several QoIs.

Figure 3: Total Sobol sensitivities

Figure 12: Sobol sensitivities

Figure 13: Sobol sensitivities

Figure 14: Sobol sensitivities
Supplementary Figures
Probability Distributions for All QoIs










Sobol Sensitivity Heatmaps for All QoIs









