- variableThe name of the variable that this object applies to
C++ Type:AuxVariableName
Controllable:No
Description:The name of the variable that this object applies to
CellMaterialIDAux
Display the OpenMC fluid material ID mapped to each MOOSE element
Description
Displays the OpenMC fluid material ID mapped to the MOOSE elements which are providing density feedback.
If a MOOSE element is not fluid or did not map at all to an OpenMC cell, then this auxiliary kernel returns .
Example Input Syntax
As an example, the material_id
auxiliary kernel will extract the OpenMC material ID and map it to the MOOSE elements corresponding to each OpenMC cell.
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
(test/tests/neutronics/feedback/lattice/openmc.i)Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- check_boundary_restrictedTrueWhether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
Default:True
C++ Type:bool
Controllable:No
Description:Whether to check for multiple element sides on the boundary in the case of a boundary restricted, element aux variable. Setting this to false will allow contribution to a single element's elemental value(s) from multiple boundary sides on the same element (example: when the restricted boundary exists on two or more sides of an element, such as at a corner of a mesh
- execute_onLINEAR TIMESTEP_ENDThe list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, POSTCHECK, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, PRE_DISPLACE.
Default:LINEAR TIMESTEP_END
C++ Type:ExecFlagEnum
Options:NONE, INITIAL, LINEAR, NONLINEAR, POSTCHECK, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, PRE_DISPLACE
Controllable:No
Description:The list of flag(s) indicating when this object should be executed, the available options include NONE, INITIAL, LINEAR, NONLINEAR, POSTCHECK, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, FINAL, CUSTOM, PRE_DISPLACE.
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether 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.
Default:False
C++ Type:bool
Controllable:No
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.
Advanced Parameters
Input Files
- (test/tests/neutronics/feedback/unmapped_moose/openmc.i)
- (test/tests/neutronics/feedback/interchangeable/no_data_flow.i)
- (test/tests/neutronics/dagmc/cell_tallies/csg_step_1/openmc.i)
- (test/tests/neutronics/feedback/different_units/openmc_cm.i)
- (tutorials/gas_assembly/openmc.i)
- (test/tests/neutronics/feedback/interchangeable/no_tally.i)
- (tutorials/pincell_multiphysics/openmc.i)
- (test/tests/neutronics/feedback/materials/openmc.i)
- (test/tests/neutronics/feedback/single_level/openmc.i)
- (test/tests/neutronics/dagmc/density_skin/openmc.i)
- (test/tests/neutronics/density/openmc.i)
- (test/tests/neutronics/dagmc/mesh_tallies/csg_step_1/openmc.i)
- (test/tests/neutronics/feedback/different_units/openmc.i)
- (test/tests/neutronics/dagmc/density_skin/csg_step_2/openmc.i)
- (test/tests/neutronics/feedback/lattice/openmc.i)
- (test/tests/neutronics/feedback/interchangeable/openmc.i)
(test/tests/neutronics/feedback/lattice/openmc.i)
[Mesh]
type = FileMesh
file = ../../meshes/pincell.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
power = 500.0
temperature_blocks = '1 2 3'
density_blocks = '2'
cell_level = 1
[Tallies]
[Cell]
type = CellTally
blocks = '1'
name = heat_source
[]
[]
[]
[Executioner]
type = Transient
[]
[Postprocessors]
[heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
[]
[fluid_heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
block = '2'
[]
[solid_heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
block = '1 3'
[]
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/feedback/unmapped_moose/openmc.i)
[Mesh]
[pincell]
type = FileMeshGenerator
file = ../../meshes/pincell.e
[]
[combine]
type = CombinerGenerator
inputs = 'pincell'
positions = '0 0 0
5 0 0'
[]
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
power = 500.0
temperature_blocks = '1 2 3'
density_blocks = '2'
verbose = true
cell_level = 0
[Tallies]
[Cell]
type = CellTally
blocks = '1'
[]
[]
[]
[Executioner]
type = Transient
[]
[Postprocessors]
[heat_source]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
[]
[fluid_heat_source]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '2'
[]
[solid_heat_source]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '1 3'
[]
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/feedback/interchangeable/no_data_flow.i)
[Mesh]
[3d]
type = GeneratedMeshGenerator
dim = 3
nx = 4
ny = 4
nz = 2
xmin = -12.5
xmax = 87.5
ymin = -12.5
ymax = 37.5
zmin = -12.5
zmax = 12.5
[]
[upper_block]
type = ParsedSubdomainMeshGenerator
input = 3d
combinatorial_geometry = 'y > 12.5'
block_id = 1
[]
[]
[AuxVariables]
[cell_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_instance]
family = MONOMIAL
order = CONSTANT
[]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[cell_id]
type = CellIDAux
variable = cell_id
[]
[cell_instance]
type = CellInstanceAux
variable = cell_instance
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[]
[ICs]
active = ''
[temp]
type = ConstantIC
variable = temp
value = 300.0
[]
[density]
type = ConstantIC
variable = density
value = 15.0e3
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[id_0]
type = ElementAverageValue
variable = cell_id
block = '0'
[]
[id_1]
type = ElementAverageValue
variable = cell_id
block = '1'
[]
[instance_0]
type = ElementAverageValue
variable = cell_instance
block = '0'
[]
[instance_1]
type = ElementAverageValue
variable = cell_instance
block = '1'
[]
[temp_0]
type = ElementAverageValue
variable = cell_temperature
block = '0'
[]
[temp_1]
type = ElementAverageValue
variable = cell_temperature
block = '1'
[]
[rho_0]
type = ElementAverageValue
variable = cell_density
block = '0'
[]
[rho_1]
type = ElementAverageValue
variable = cell_density
block = '1'
[]
[mat_0]
type = ElementAverageValue
variable = material_id
block = '0'
[]
[mat_1]
type = ElementAverageValue
variable = material_id
block = '1'
[]
[]
[Executioner]
type = Transient
num_steps = 1
[]
[Outputs]
csv = true
[]
(test/tests/neutronics/dagmc/cell_tallies/csg_step_1/openmc.i)
[Mesh]
type = FileMesh
file = ../../mesh_tallies/slab.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[]
[AuxKernels]
[temp]
type = ConstantAux
variable = temp
value = 500.0
execute_on = timestep_begin
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
temperature_blocks = '1 2'
cell_level = 0
power = 100.0
[Tallies]
[Cell]
type = CellTally
blocks = '1 2'
[]
[]
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[]
[Executioner]
type = Transient
num_steps = 2
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/feedback/different_units/openmc_cm.i)
[Mesh]
type = FileMesh
file = sphere_in_cm.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
power = 1e6
temperature_blocks = '1'
density_blocks = '1'
cell_level = 0
[Tallies]
[Cell]
type = CellTally
blocks = '1'
name = 'heat_source'
[]
[]
[]
[Executioner]
type = Transient
[]
[Postprocessors]
[heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
[]
[]
[Outputs]
exodus = true
[]
(tutorials/gas_assembly/openmc.i)
num_layers_for_THM = 50 # number of elements in the THM model; for the converged
# case, we set this to 150
[Mesh]
# mesh mirror for the solid regions
[solid]
type = FileMeshGenerator
file = solid_mesh_in.e
[]
# create a mesh for a single coolant channel; because we will receive uniform
# temperatures and densities from THM on each x-y plane, we can use a very coarse
# mesh in the radial direction
[coolant_face]
type = AnnularMeshGenerator
nr = 1
nt = 8
rmin = 0.0
rmax = ${fparse channel_diameter / 2.0}
[]
[extrude]
type = AdvancedExtruderGenerator
input = coolant_face
num_layers = ${num_layers_for_THM}
direction = '0 0 1'
heights = '${height}'
top_boundary = '300' # inlet
bottom_boundary = '400' # outlet
[]
[rename]
type = RenameBlockGenerator
input = extrude
old_block = '1'
new_block = '101'
[]
# repeat the coolant channels and then combine together to get a combined mesh mirror
[repeat]
type = CombinerGenerator
inputs = rename
positions_file = coolant_channel_positions.txt
[]
[add]
type = CombinerGenerator
inputs = 'solid repeat'
[]
[]
[AuxVariables]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[thm_temp_wall]
block = '101'
[]
[flux]
[]
# just for postprocessing purposes
[thm_pressure]
block = '101'
[]
[thm_velocity]
block = '101'
[]
[z]
family = MONOMIAL
order = CONSTANT
block = 'compacts'
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[density]
type = FluidDensityAux
variable = density
p = ${outlet_P}
T = thm_temp
fp = helium
execute_on = 'timestep_begin linear'
[]
[z]
type = ParsedAux
variable = z
use_xyzt = true
expression = 'z'
[]
[]
[ICs]
[fluid_temp_wall]
type = FunctionIC
variable = thm_temp_wall
function = temp_ic
[]
[fluid_temp]
type = FunctionIC
variable = thm_temp
function = temp_ic
[]
[heat_source]
type = ConstantIC
variable = heat_source
block = 'compacts'
value = ${fparse power / (n_bundles * n_fuel_compacts_per_block) / (pi * compact_diameter * compact_diameter / 4.0 * height)}
[]
[]
[Functions]
[temp_ic]
type = ParsedFunction
expression = '${inlet_T} + (${height} - z) / ${height} * ${power} / ${mdot} / ${fluid_Cp}'
[]
[]
[FluidProperties]
[helium]
type = IdealGasFluidProperties
molar_mass = 4e-3
gamma = 1.668282 # should correspond to Cp = 5189 J/kg/K
k = 0.2556
mu = 3.22639e-5
[]
[]
[Problem]
type = OpenMCCellAverageProblem
identical_cell_fills = '2'
power = ${fparse power / n_bundles}
scaling = 100.0
cell_level = 1
relaxation = constant
relaxation_factor = 0.5
# to get a faster-running tutorial, we use only 1000 particles per batch; converged
# results are instead obtained by increasing this parameter to 10000. We also use fewer
# batches to speed things up; the converged results were obtained with 500 inactive batches
# and 2000 active batches
particles = 1000
inactive_batches = 200
batches = 1000
# we will read temperature from THM (for the fluid) and MOOSE (for the solid)
# into variables we name as 'solid_temp' and 'thm_temp'. This syntax will automatically
# create those variabes for us
temperature_variables = 'solid_temp; thm_temp'
temperature_blocks = '1 2 4; 101'
density_blocks = '101'
[Tallies]
[heat_source]
type = CellTally
blocks = '2'
name = heat_source
check_equal_mapped_tally_volumes = true
output = 'unrelaxed_tally_std_dev'
[]
[]
[]
[MultiApps]
[bison]
type = TransientMultiApp
input_files = 'solid.i'
execute_on = timestep_begin
[]
[thm]
type = FullSolveMultiApp
input_files = 'thm.i'
execute_on = timestep_end
max_procs_per_app = 1
bounding_box_padding = '0.1 0.1 0'
positions_file = coolant_channel_positions.txt
output_in_position = true
[]
[]
[Transfers]
[solid_temp_to_openmc]
type = MultiAppGeometricInterpolationTransfer
source_variable = T
variable = solid_temp
from_multi_app = bison
[]
[heat_flux_to_openmc]
type = MultiAppGeneralFieldNearestLocationTransfer
source_variable = flux
variable = flux
from_multi_app = bison
from_boundaries = 'fluid_solid_interface'
to_boundaries = 'fluid_solid_interface'
from_postprocessors_to_be_preserved = flux_integral
to_postprocessors_to_be_preserved = flux_integral
[]
[source_to_bison]
type = MultiAppGeneralFieldShapeEvaluationTransfer
source_variable = heat_source
variable = power
to_multi_app = bison
from_postprocessors_to_be_preserved = heat_source
to_postprocessors_to_be_preserved = power
[]
[thm_temp_to_bison]
type = MultiAppGeometricInterpolationTransfer
source_variable = thm_temp_wall
variable = thm_temp
to_multi_app = bison
[]
[q_wall_to_thm]
type = MultiAppGeneralFieldUserObjectTransfer
variable = q_wall
to_multi_app = thm
source_user_object = q_wall_avg
[]
[T_wall_from_thm]
type = MultiAppGeneralFieldNearestLocationTransfer
source_variable = T_wall
from_multi_app = thm
variable = thm_temp_wall
to_boundaries = 'fluid_solid_interface'
[]
[T_bulk_from_thm]
type = MultiAppGeneralFieldNearestLocationTransfer
source_variable = T
from_multi_app = thm
variable = thm_temp
[]
# just for postprocessing purposes
[pressure_from_thm]
type = MultiAppGeneralFieldNearestLocationTransfer
source_variable = p
from_multi_app = thm
variable = thm_pressure
[]
[velocity_from_thm]
type = MultiAppGeneralFieldNearestLocationTransfer
source_variable = vel_z
from_multi_app = thm
variable = thm_velocity
[]
[]
[UserObjects]
[q_wall_avg]
type = NearestPointLayeredSideAverage
boundary = 'fluid_solid_interface'
variable = flux
# Note: make this to match the num_elems in the channel
direction = z
num_layers = ${num_layers_for_THM}
points_file = coolant_channel_positions.txt
direction_min = 0.0
direction_max = ${height}
[]
[average_power_axial]
type = LayeredAverage
variable = heat_source
direction = z
num_layers = ${num_layers_for_plots}
block = 'compacts'
[]
[average_fluid_axial]
type = LayeredAverage
variable = thm_temp
direction = z
num_layers = ${num_layers_for_plots}
block = '101'
[]
[average_pressure]
type = LayeredAverage
variable = thm_pressure
direction = z
num_layers = ${num_layers_for_plots}
block = '101'
[]
[average_axial_velocity]
type = LayeredAverage
variable = thm_velocity
direction = z
num_layers = ${num_layers_for_plots}
block = '101'
[]
[]
[VectorPostprocessors]
[power_avg]
type = SpatialUserObjectVectorPostprocessor
userobject = average_power_axial
[]
[fluid_avg]
type = SpatialUserObjectVectorPostprocessor
userobject = average_fluid_axial
[]
[pressure_avg]
type = SpatialUserObjectVectorPostprocessor
userobject = average_pressure
[]
[velocity_avg]
type = SpatialUserObjectVectorPostprocessor
userobject = average_axial_velocity
[]
[]
[Postprocessors]
[flux_integral]
type = SideIntegralVariablePostprocessor
variable = flux
boundary = 'fluid_solid_interface'
execute_on = 'transfer linear'
[]
[heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
execute_on = 'transfer initial timestep_end'
[]
[max_tally_rel_err]
type = TallyRelativeError
value_type = max
[]
[k]
type = KEigenvalue
[]
[k_std_dev]
type = KStandardDeviation
[]
[min_power]
type = ElementExtremeValue
variable = heat_source
value_type = min
block = 'compacts'
[]
[max_power]
type = ElementExtremeValue
variable = heat_source
value_type = max
block = 'compacts'
[]
[z_max_power]
type = ElementExtremeValue
proxy_variable = heat_source
variable = z
block = 'compacts'
[]
[max_Tf]
type = ElementExtremeValue
variable = thm_temp
block = '101'
[]
[P_in]
type = SideAverageValue
variable = thm_pressure
boundary = '300'
[]
[pressure_drop]
type = LinearCombinationPostprocessor
pp_names = 'P_in'
pp_coefs = '1.0'
b = '${fparse -outlet_P}'
[]
[]
[Executioner]
type = Transient
num_steps = 10
[]
[Outputs]
exodus = true
csv = true
hide = 'P_in flux_integral z'
[]
(test/tests/neutronics/feedback/interchangeable/no_tally.i)
[Mesh]
[3d]
type = GeneratedMeshGenerator
dim = 3
nx = 4
ny = 4
nz = 2
xmin = -12.5
xmax = 87.5
ymin = -12.5
ymax = 37.5
zmin = -12.5
zmax = 12.5
[]
[upper_block]
type = ParsedSubdomainMeshGenerator
input = 3d
combinatorial_geometry = 'y > 12.5'
block_id = 1
[]
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[]
[ICs]
active = ''
[temp]
type = ConstantIC
variable = temp
value = 300.0
[]
[density]
type = ConstantIC
variable = density
value = 15.0e3
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
cell_level = 0
power = 1.0
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[power_0]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '0'
[]
[power_1]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '1'
[]
[id_0]
type = ElementAverageValue
variable = cell_id
block = '0'
[]
[id_1]
type = ElementAverageValue
variable = cell_id
block = '1'
[]
[instance_0]
type = ElementAverageValue
variable = cell_instance
block = '0'
[]
[instance_1]
type = ElementAverageValue
variable = cell_instance
block = '1'
[]
[temp_0]
type = ElementAverageValue
variable = cell_temperature
block = '0'
[]
[temp_1]
type = ElementAverageValue
variable = cell_temperature
block = '1'
[]
[rho_0]
type = ElementAverageValue
variable = cell_density
block = '0'
[]
[rho_1]
type = ElementAverageValue
variable = cell_density
block = '1'
[]
[mat_0]
type = ElementAverageValue
variable = material_id
block = '0'
[]
[mat_1]
type = ElementAverageValue
variable = material_id
block = '1'
[]
[]
[Executioner]
type = Transient
num_steps = 1
[]
[Outputs]
csv = true
[]
(tutorials/pincell_multiphysics/openmc.i)
inlet_T = 573.0 # inlet temperature
power = 1e3 # total power (W)
Re = 500.0 # Reynolds number
outlet_P = 1e6
height = 0.5 # total height of the domain
Df = 0.825e-2 # fuel diameter
pin_diameter = 0.97e-2 # pin outer diameter
pin_pitch = 1.28e-2 # pin pitch
mu = 8.8e-5 # fluid dynamic viscosity
rho = 723.6 # fluid density
Cp = 5512.0 # fluid isobaric specific heat capacity
Rf = ${fparse Df / 2.0}
flow_area = ${fparse pin_pitch * pin_pitch - pi * pin_diameter * pin_diameter / 4.0}
wetted_perimeter = ${fparse pi * pin_diameter}
hydraulic_diameter = ${fparse 4.0 * flow_area / wetted_perimeter}
U_ref = ${fparse Re * mu / rho / hydraulic_diameter}
mdot = ${fparse rho * U_ref * flow_area}
dT = ${fparse power / mdot / Cp}
[Mesh]
[solid]
type = FileMeshGenerator
file = solid_in.e
[]
[]
[AuxVariables]
# These auxiliary variables are all just for visualizing the solution and
# the mapping - none of these are part of the calculation sequence
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[density]
type = FluidDensityAux
variable = density
p = ${outlet_P}
T = nek_temp
fp = sodium
execute_on = 'INITIAL TIMESTEP_END'
[]
[]
[FluidProperties]
[sodium]
type = SodiumSaturationFluidProperties
[]
[]
[ICs]
[nek_temp]
type = FunctionIC
variable = nek_temp
function = temp_ic
[]
[solid_temp]
type = FunctionIC
variable = solid_temp
function = temp_ic
[]
[heat_source]
type = ConstantIC
variable = heat_source
block = '2'
value = ${fparse power / (pi * Rf * Rf * height)}
[]
[]
[Functions]
[temp_ic]
type = ParsedFunction
expression = '${inlet_T} + z / ${height} * ${dT}'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
power = ${power}
scaling = 100.0
density_blocks = '1'
cell_level = 0
# This automatically creates these variables and will read from the non-default choice of 'temp'
temperature_variables = 'solid_temp; nek_temp'
temperature_blocks = '2 3; 1'
relaxation = robbins_monro
# Set some parameters for when we terminate the OpenMC solve in each iteration;
# this will run a minimum of 30 batches, and after that, terminate once reaching
# the specified std. dev. of k and rel. err. of the fission tally
inactive_batches = 20
batches = 30
k_trigger = std_dev
k_trigger_threshold = 7.5e-4
batch_interval = 50
max_batches = 1000
[Tallies]
[heat_source]
type = CellTally
blocks = '2'
name = heat_source
check_equal_mapped_tally_volumes = true
trigger = rel_err
trigger_threshold = 2e-2
output = unrelaxed_tally_std_dev
[]
[]
[]
[MultiApps]
[bison]
type = TransientMultiApp
input_files = 'bison.i'
execute_on = timestep_begin
sub_cycling = true
[]
[]
[Transfers]
[solid_temp_to_openmc]
type = MultiAppGeneralFieldShapeEvaluationTransfer
source_variable = T
variable = solid_temp
from_multi_app = bison
[]
[source_to_bison]
type = MultiAppCopyTransfer
source_variable = heat_source
variable = power
to_multi_app = bison
[]
[temp_from_nek]
type = MultiAppGeneralFieldShapeEvaluationTransfer
source_variable = nek_temp
from_multi_app = bison
variable = nek_temp
[]
[]
[Outputs]
exodus = true
csv = true
[]
[Executioner]
type = Transient
dt = 0.5
num_steps = 10
[]
[Postprocessors]
[max_tally_rel_err]
type = TallyRelativeError
value_type = max
[]
[k]
type = KEigenvalue
[]
[k_std_dev]
type = KStandardDeviation
[]
[]
(test/tests/neutronics/feedback/materials/openmc.i)
[Mesh]
type = FileMesh
file = ../../meshes/pincell.e
[]
[ICs]
[temp]
type = FunctionIC
variable = temp
function = temp
[]
[density]
type = FunctionIC
variable = density
function = density
[]
[]
[Functions]
[temp]
type = ParsedFunction
expression = '500'
[]
[density]
type = ParsedFunction
expression = '800+z*10'
[]
[]
[AuxVariables]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
power = 500.0
temperature_blocks = '1 2 3'
density_blocks = '2'
cell_level = 0
map_density_by_cell = false
[Tallies]
[Cell]
type = CellTally
blocks = '1'
[]
[]
[]
[Postprocessors]
[avg_density]
type = ElementAverageValue
variable = density
block = '2'
[]
[rho1]
type = PointValue
variable = cell_density
point = '0.5 0 1.5'
[]
[rho2]
type = PointValue
variable = cell_density
point = '0.5 0 8.5'
[]
[mat1]
type = PointValue
variable = material_id
point = '0.5 0 1.5'
[]
[mat2]
type = PointValue
variable = material_id
point = '0.5 0 8.5'
[]
[k]
type = KEigenvalue
value_type = collision
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
csv = true
[]
(test/tests/neutronics/feedback/single_level/openmc.i)
[Mesh]
type = FileMesh
file = ../../meshes/pincell.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
power = 500.0
temperature_blocks = '1 2 3'
density_blocks = '2'
verbose = true
cell_level = 0
[Tallies]
[Cell]
type = CellTally
blocks = '1'
[]
[]
[]
[Executioner]
type = Transient
[]
[Postprocessors]
[kappa_fission]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
[]
[fluid_kappa_fission]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '2'
[]
[solid_kappa_fission]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '1 3'
[]
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/dagmc/density_skin/openmc.i)
[Mesh]
type = FileMesh
file = ../mesh_tallies/slab.e
allow_renumbering = false
[]
x0 = 12.5
x1 = 37.5
x2 = 62.5
[AuxVariables]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[]
rho0 = 600.0
drho = 50.0
[Functions]
[density]
type = ParsedFunction
expression = 'if (x <= ${x0}, ${fparse rho0 - drho}, if (x <= ${x1}, ${rho0}, if (x <= ${x2}, ${fparse rho0 + drho}, ${fparse rho0 + 2 * drho})))'
[]
[]
[ICs]
[temp]
type = ConstantIC
variable = temp
value = 500.0
[]
[]
[AuxKernels]
[density]
type = FunctionAux
variable = density
function = density
execute_on = timestep_begin
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
temperature_blocks = '1 2'
density_blocks = '1 2'
power = 100.0
skinner = moab
cell_level = 0
[Tallies]
[Mesh]
type = CellTally
blocks = '1 2'
[]
[]
[]
[UserObjects]
[moab]
type = MoabSkinner
# just one temperature bin
temperature_min = 0.0
temperature_max = 1000.0
n_temperature_bins = 1
temperature = temp
density_min = ${fparse rho0 - drho}
density_max = ${fparse rho0 + 2 * drho}
n_density_bins = 4
density = density
build_graveyard = true
[]
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[]
[Executioner]
type = Transient
num_steps = 2
[]
[Outputs]
exodus = true
hide = 'temp density cell_instance cell_id'
[]
(test/tests/neutronics/density/openmc.i)
[Mesh]
[file]
type = FileMeshGenerator
file = ../dagmc/mesh_tallies/slab.e
[]
[split_right]
type = ParsedSubdomainMeshGenerator
input = file
combinatorial_geometry = 'x > 37.5'
block_id = 100
excluded_subdomains = '1'
[]
[split_left]
type = ParsedSubdomainMeshGenerator
input = split_right
combinatorial_geometry = 'x < 37.5'
block_id = 200
excluded_subdomains = '2'
[]
[scale]
type = TransformGenerator
input = split_left
transform = translate
vector_value = '0.0 0.0 0.0'
[]
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[cell_material]
family = MONOMIAL
order = CONSTANT
[]
[]
x0 = 12.5
x1 = 37.5
x2 = 62.5
T0 = 600.0
dT = 50.0
drho = 100.0
rho0 = 1000.0
[Functions]
[temp]
type = ParsedFunction
expression = 'if (x <= ${x0}, ${fparse T0 - dT}, if (x <= ${x1}, ${T0}, if (x <= ${x2}, ${fparse T0 + dT}, ${fparse T0 + 2 * dT})))'
[]
[rho]
type = ParsedFunction
expression = 'if (x <= ${x0}, ${fparse rho0 - drho}, if (x <= ${x1}, ${rho0}, if (x <= ${x2}, ${fparse rho0 + drho}, ${fparse rho0 + 2 * drho})))'
[]
[]
[AuxKernels]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[cell_material]
type = CellMaterialIDAux
variable = cell_material
[]
[]
[AuxKernels]
[temp]
type = FunctionAux
variable = temp
function = temp
execute_on = timestep_begin
[]
[rho]
type = FunctionAux
variable = density
function = rho
execute_on = timestep_begin
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
density_blocks = '2 1'
temperature_blocks = '200 1'
cell_level = 0
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[rho_2]
type = ElementAverageValue
variable = density
block = '2'
[]
[rho_1]
type = ElementAverageValue
variable = density
block = '1'
[]
[T_200]
type = ElementAverageValue
variable = temp
block = '200'
[]
[T_1]
type = ElementAverageValue
variable = temp
block = '1'
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/dagmc/mesh_tallies/csg_step_1/openmc.i)
[Mesh]
type = FileMesh
file = ../slab.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[]
[AuxKernels]
[temp]
type = ConstantAux
variable = temp
value = 500.0
execute_on = timestep_begin
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
temperature_blocks = '1 2'
cell_level = 0
power = 100.0
[Tallies]
[Mesh]
type = MeshTally
mesh_template = ../slab.e
[]
[]
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[]
[Executioner]
type = Transient
num_steps = 2
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/feedback/different_units/openmc.i)
[Mesh]
type = FileMesh
file = sphere_in_m.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
power = 1e6
temperature_blocks = '1'
density_blocks = '1'
cell_level = 0
scaling = 100.0
[Tallies]
[Cell]
type = CellTally
blocks = '1'
[]
[]
[]
[Executioner]
type = Transient
[]
[Postprocessors]
[heat_source]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
[]
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/dagmc/density_skin/csg_step_2/openmc.i)
[Mesh]
type = FileMesh
file = ../../mesh_tallies/slab.e
[]
x0 = 12.5
x1 = 37.5
x2 = 62.5
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[]
T0 = 600.0
dT = 50.0
[Functions]
[density]
type = ParsedFunction
expression = 'if (x <= ${x0}, ${fparse T0 - dT}, if (x <= ${x1}, ${T0}, if (x <= ${x2}, ${fparse T0 + dT}, ${fparse T0 + 2 * dT})))'
[]
[]
[ICs]
[temp]
type = ConstantIC
variable = temp
value = 500.0
[]
[]
[AuxKernels]
[density]
type = FunctionAux
variable = density
function = density
execute_on = timestep_begin
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
temperature_blocks = '1 2'
density_blocks = '1 2'
cell_level = 0
power = 100.0
[Tallies]
[Mesh]
type = CellTally
blocks = '1 2'
[]
[]
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[]
[Executioner]
type = Transient
num_steps = 2
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/feedback/lattice/openmc.i)
[Mesh]
type = FileMesh
file = ../../meshes/pincell.e
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
execute_on = 'timestep_end'
[]
[cell_density]
type = CellDensityAux
variable = cell_density
execute_on = 'timestep_end'
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
power = 500.0
temperature_blocks = '1 2 3'
density_blocks = '2'
cell_level = 1
[Tallies]
[Cell]
type = CellTally
blocks = '1'
name = heat_source
[]
[]
[]
[Executioner]
type = Transient
[]
[Postprocessors]
[heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
[]
[fluid_heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
block = '2'
[]
[solid_heat_source]
type = ElementIntegralVariablePostprocessor
variable = heat_source
block = '1 3'
[]
[]
[Outputs]
exodus = true
[]
(test/tests/neutronics/feedback/interchangeable/openmc.i)
[Mesh]
[3d]
type = GeneratedMeshGenerator
dim = 3
nx = 4
ny = 4
nz = 2
xmin = -12.5
xmax = 87.5
ymin = -12.5
ymax = 37.5
zmin = -12.5
zmax = 12.5
[]
[upper_block]
type = ParsedSubdomainMeshGenerator
input = 3d
combinatorial_geometry = 'y > 12.5'
block_id = 1
[]
[]
[AuxVariables]
[cell_temperature]
family = MONOMIAL
order = CONSTANT
[]
[material_id]
family = MONOMIAL
order = CONSTANT
[]
[cell_density]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[cell_temperature]
type = CellTemperatureAux
variable = cell_temperature
[]
[material_id]
type = CellMaterialIDAux
variable = material_id
[]
[cell_density]
type = CellDensityAux
variable = cell_density
[]
[]
[ICs]
active = ''
[temp]
type = ConstantIC
variable = temp
value = 300.0
[]
[density]
type = ConstantIC
variable = density
value = 15.0e3
[]
[]
[Problem]
type = OpenMCCellAverageProblem
verbose = true
cell_level = 0
power = 1.0
[Tallies]
[Cell]
type = CellTally
blocks = '0 1'
[]
[]
[]
[Postprocessors]
[k]
type = KEigenvalue
[]
[power_0]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '0'
[]
[power_1]
type = ElementIntegralVariablePostprocessor
variable = kappa_fission
block = '1'
[]
[id_0]
type = ElementAverageValue
variable = cell_id
block = '0'
[]
[id_1]
type = ElementAverageValue
variable = cell_id
block = '1'
[]
[instance_0]
type = ElementAverageValue
variable = cell_instance
block = '0'
[]
[instance_1]
type = ElementAverageValue
variable = cell_instance
block = '1'
[]
[temp_0]
type = ElementAverageValue
variable = cell_temperature
block = '0'
[]
[temp_1]
type = ElementAverageValue
variable = cell_temperature
block = '1'
[]
[rho_0]
type = ElementAverageValue
variable = cell_density
block = '0'
[]
[rho_1]
type = ElementAverageValue
variable = cell_density
block = '1'
[]
[mat_0]
type = ElementAverageValue
variable = material_id
block = '0'
[]
[mat_1]
type = ElementAverageValue
variable = material_id
block = '1'
[]
[]
[Executioner]
type = Transient
num_steps = 1
[]
[Outputs]
csv = true
[]