Grid to Rod Fretting Wear

One leading cause of failure in the rod assembly is fuel leakage into the coolant. The primary conducive mechanism to that failure mode is wear on the cladding material, whose gradual material loss weakens the ability of cladding material to safely contain nuclear fuel.

This BISON example utilizes a two-dimensional (2D) plane strain LWR section to demonstrate the use of mortar mechanical contact to predict small wear depths in a dynamic simulation, using Newmark-. Wear is only accounted for in the contact equations, so the approach used here is only valid to predict the onset of wear, in which the wear depth is significantly smaller than the cladding thickness.

The contact formulation incorporates dynamic variables (velocity) in the contact equations to stabilize the mortar contact constraints (Recuero and Lindsay, 2023) and a nodal wear depth computed in the weak sense is introduced in the contact equations to modify the contact interface geometry. Wear is based on Archard's law. Note that, due to the difference in time scales between the dynamic excitation caused by flow-induced vibrations (20-30 Hz) and the natural evolution of the rest of the physics, a step-wise procedure to alternate time spans with dynamic excitation and time spans without dynamic excitation needs to be devised (Recuero et al., 2023).

In this problem, we impose relative displacement beween the spacer grid and the cladding material through Dirichlet boundary conditions in and directions at approximately the frequencies and amplitudes at which flow-induced vibrations occur. More accurate multiphysics coupling (with fluid dynamics) and material removal (arbitrary Lagrangian Eulerian, e.g.) formulations can be used.

Archard's Law

Wear is considered to be proportional to the contact pressure, the relative displacement between the contacting surfaces, and a wear coefficient that captures the properties of the materials into contact. Note that, due to alternating "dynamic" and "quasi-static" steps, the wear coefficient in the dynamic step needs to account for the fact that quasi-static steps do not contribute to the generation of wear.

[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
  [worn_depth]
    type = MortarArchardsLawAux<<<{"description": "Returns the weighted gap velocity at a node. This quantity is useful for mortar contact, particularly when dual basis functions are used in contact mechanics", "href": "../source/auxkernels/MortarArchardsLawAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = worn_depth
    primary_boundary<<<{"description": "The name of the primary boundary sideset."}>>> = 101
    secondary_boundary<<<{"description": "The name of the secondary boundary sideset."}>>> = 102
    primary_subdomain<<<{"description": "The name of the primary subdomain."}>>> = 'spacer_clad_mechanical_primary_subdomain'
    secondary_subdomain<<<{"description": "The name of the secondary subdomain."}>>> = 'spacer_clad_mechanical_secondary_subdomain'
    displacements<<<{"description": "The displacement variables. This mortar nodal auxiliary kernel can take two or three displacements"}>>> = 'disp_x disp_y'
    friction_coefficient<<<{"description": "Friction coefficient used to compute wear (to match that of frictional contact)"}>>> = 0.5
    energy_wear_coefficient<<<{"description": "Energy wear coefficient is a surface-dependent parameter used in Archard's wear law"}>>> = 0.1e-9
    normal_pressure<<<{"description": "The name of the Lagrange multiplier that holds the normal contact pressure in mortar formulations"}>>> = spacer_clad_mechanical_normal_lm
    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'
  []
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)

The auxiliary kernel is executed at time step end. This setup lags the local wear depth one time step, which facilitates convergence and is strongly recommended. Said lagging is typically not of significance since dynamic time steps are very small compared to the thermomechanical and wear evolution of the system.

Contact Block

In this example we consider fuel-cladding contact and cladding-spacer grid contact. Fuel cladding contact is considered to be frictionless and does not take place in the setup of this example (see Recuero et al. (2023)) for a 3-year simulation in which fuel-cladding contact is of relevance. Cladding-spacer grid contact is considered to be frictional and coupled to wear defined by Archard's law.

For fuel-cladding contact, no wear is considered:

[Contact<<<{"href": "../syntax/Contact/index.html"}>>>]
  [pellet_clad_mechanical_real]
    formulation<<<{"description": "The contact formulation"}>>> = mortar
    model<<<{"description": "The contact model to use"}>>> = frictionless
    primary<<<{"description": "The list of boundary IDs referring to primary sidesets"}>>> = 7
    secondary<<<{"description": "The list of boundary IDs referring to secondary sidesets"}>>> = 8
    c_normal<<<{"description": "Parameter for balancing the size of the gap and contact pressure for a mortar formulation. This purely numerical parameter affects convergence behavior and, in general, should be larger for stiffer materials. It is recommended that the user tries out various orders of magnitude for this parameter if the default value generates poor contact convergence."}>>> = 1e+16 #
    c_tangential<<<{"description": "Numerical parameter for nonlinear mortar frictional constraints"}>>> = 1e+16
    friction_coefficient<<<{"description": "The friction coefficient"}>>> = 0.4
    # Do not apply dynamic stabilization
    newmark_beta<<<{"description": "Newmark-beta beta parameter for its inclusion in the weighted gap update formula"}>>> = 0.0001
    newmark_gamma<<<{"description": "Newmark-beta gamma parameter for its inclusion in the weighted gap update formula"}>>> = 0.5
    capture_tolerance<<<{"description": "Normal distance from surface within which nodes are captured. This parameter is used for node-face and mortar formulations."}>>> = 0.0
    mortar_dynamics<<<{"description": "Whether to use constraints that account for the persistency condition, giving rise to smoother normal contact pressure evolution. This flag should only be set to yes for dynamic simulations using the Newmark-beta numerical integrator"}>>> = true
    interpolate_normals = false
    generate_mortar_mesh<<<{"description": "Whether to generate the mortar mesh from the action. Typically this will be the case, but one may also want to reuse an existing lower-dimensional mesh prior to a restart."}>>> = false
  []
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)

For spacer grid-cladding contact, we provide wear_depth as an auxiliary variable, which is considered in the normal contact equations. The input parameter generate_mortar_mesh is selected to be false since mortar lower-dimensional domains will be read from the restart files.

Boundary Conditions

Relative displacement between cladding and spacer grid conducive to grid-to-rod fretting wear are captured by time-dependent Dirichlet boundary conditions:

[BCs<<<{"href": "../syntax/NuclearMaterials/BCs/index.html"}>>>]
  [vibration_x]
    # pin pellets and clad along axis of symmetry (y)
    type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.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"}>>> = '112'
    expression = 'if(t < ${end_dynamic_excitation}, 10.0*1.0e-6*sin(2*3.1415926535*20* (t - ${user_start_time})) + 2.0*1.0e-6*sin(2*3.1415926535*35*(t - ${user_start_time})), 0)'
    #expression = '0'
  []
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)

and

[BCs<<<{"href": "../syntax/NuclearMaterials/BCs/index.html"}>>>]
  [vibration_y]
    # pin pellets and clad along axis of symmetry (y)
    type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.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"}>>> = '112'
    expression = 'if(t < ${end_dynamic_excitation}, 10.0*1.0e-6*sin(2*3.1415926535*20*(t-${user_start_time})) + 2.0*1.0e-6*sin(2*3.1415926535*35*(t-${user_start_time})) + 0.9e-4, 0.9e-4)'
    #expression = '5.9e-4'
  []
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)

These conditions approximate flow-induced vibrations.

Restart

In order to switch between simulation time spans with dynamic excitation, which require "dynamic" time steps, and time spans without dynamic excitation, which require typical Bison time steps, we employ restart files. That way, we can modify executioner or algorithmic parameters as needed. Here, only one restart step is shown, but, in practice, to achieve convergence between "dynamic" and "quasi-static" steps, a larger number of restart runs need to be employed.

Restart files must be provided to the mesh and the problem (i.e. restart_file_base), as follows:

[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
  [file]
    type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../source/meshgenerators/FileMeshGenerator.html"}>>>
    file<<<{"description": "The filename to read."}>>> = fretting-wear-initial_out_cp/LATEST
    skip_partitioning<<<{"description": "True to skip partitioning, only after this mesh generator, because the mesh was pre-split for example."}>>> = true
    allow_renumbering<<<{"description": "Whether to allow the mesh to renumber nodes and elements, if not overridden by a global parameter or by a requirement (e.g. an exodus restart or a constraint matrix) that disables renumbering."}>>> = false
  []
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)
[Problem<<<{"href": "../syntax/Problem/index.html"}>>>]
  type = ReferenceResidualProblem
  reference_vector = 'ref'
  extra_tag_vectors = 'ref'
  group_variables = 'disp_x disp_y '
  converge_on = 'disp_x disp_y  temp'
  restart_file_base = ./fretting-wear-initial_out_cp/LATEST
  material_coverage_check = false
  kernel_coverage_check = false
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)

Conclusion

This example partially shows a methodology for the prediction of grid-to-rod fretting wear leveraging Bison models. Assumptions include: 1) Wear depth is small and 2) The action of flow-induced vibrations is known and can be approximated by boundary conditions.

initial_fuel_density = 10431.0

[GlobalParams<<<{"href": "../syntax/GlobalParams/index.html"}>>>]
  temperature = temp
  displacements = 'disp_x disp_y'
  order = FIRST
  family = LAGRANGE
  energy_per_fission = 3.2e-11 # J/fission
  volumetric_locking_correction = true
[]

[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
  [file]
    type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../source/meshgenerators/FileMeshGenerator.html"}>>>
    file<<<{"description": "The filename to read."}>>> = refined_excitation_better_mesh.e
  []
  construct_node_list_from_side_list = true
  patch_size = 100 # For contact algorithm
[]

[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
  [temp]
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 580.0 # set initial temp to ambient
  []
[]

[Problem<<<{"href": "../syntax/Problem/index.html"}>>>]
  type = ReferenceResidualProblem
  reference_vector = 'ref'
  extra_tag_vectors = 'ref'
  group_variables = 'disp_x disp_y '
  converge_on = 'disp_x disp_y  temp'
  material_coverage_check = false
  kernel_coverage_check = false
  #  restart_file_base = planestrain_grid_aux_vars_out_cp/LATEST
[]

[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
  [fission_rate]
    block = pellet_type_1
  []
  [burnup]
    block = pellet_type_1
  []
  [fast_neutron_flux]
    block = 'clad grid'
  []
  [fast_neutron_fluence]
    block = 'clad grid'
  []
  [relocation_strain]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [worn_depth]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
    block = 'spacer_clad_mechanical_secondary_subdomain'
  []
[]

[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
  [power_history]
    type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>> # reads and interpolates an input file containing rod average linear power vs time
    data_file<<<{"description": "File holding CSV data"}>>> = powerhistory.csv
    scale_factor<<<{"description": "Scale factor to be applied to the ordinate values"}>>> = 1
  []
  [axial_peaking_factors]
    type = ConstantFunction<<<{"description": "A function that returns a constant value as defined by an input parameter.", "href": "../source/functions/ConstantFunction.html"}>>>
    value<<<{"description": "The constant value"}>>> = 1
  []
  [pressure_var] # reads and interpolates input data defining amplitude curve for fill gas pressure
    type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
    x<<<{"description": "The abscissa values"}>>> = '0 1e4'
    y<<<{"description": "The ordinate values"}>>> = '0 1'
  []
  [pressure_var_variable] # reads and interpolates input data defining amplitude curve for fill gas pressure
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = 'if(t < 1e4, 1, 1 + sin((t-1e4)*pi/10.0) * (t-1e4))'
  []
[]

[Physics<<<{"href": "../syntax/Physics/index.html"}>>>/SolidMechanics<<<{"href": "../syntax/Physics/SolidMechanics/index.html"}>>>/Dynamic<<<{"href": "../syntax/Physics/SolidMechanics/Dynamic/index.html"}>>>]
  [pellets]
    add_variables<<<{"description": "Add the displacement variables"}>>> = true
    newmark_beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
    newmark_gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
    block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = pellet_type_1
    strain<<<{"description": "Strain formulation"}>>> = FINITE
    planar_formulation<<<{"description": "Out-of-plane stress/strain formulation"}>>> = PLANE_STRAIN
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'fuel_relocation_eigenstrain fuel_thermal_eigenstrain
      fuel_volumetric_eigenstrain'
    generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'vonmises_stress stress_xx stress_yy stress_zz'
    decomposition_method<<<{"description": "Methods to calculate the finite strain and rotation increments"}>>> = EigenSolution
    temperature<<<{"description": "The temperature"}>>> = temp
    extra_vector_tags<<<{"description": "The tag names for extra vectors that residual data should be saved into"}>>> = 'ref'
  []
  [clad]
    add_variables<<<{"description": "Add the displacement variables"}>>> = true
    newmark_beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
    newmark_gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
    block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = clad
    strain<<<{"description": "Strain formulation"}>>> = FINITE
    planar_formulation<<<{"description": "Out-of-plane stress/strain formulation"}>>> = PLANE_STRAIN
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'clad_thermal_eigenstrain clad_irradiation_eigenstrain'
    generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'vonmises_stress stress_xx stress_yy stress_zz'
    decomposition_method<<<{"description": "Methods to calculate the finite strain and rotation increments"}>>> = EigenSolution
    temperature<<<{"description": "The temperature"}>>> = temp
    extra_vector_tags<<<{"description": "The tag names for extra vectors that residual data should be saved into"}>>> = 'ref'
  []
  [grid]
    add_variables<<<{"description": "Add the displacement variables"}>>> = true
    newmark_beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
    newmark_gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
    block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = grid
    strain<<<{"description": "Strain formulation"}>>> = FINITE
    planar_formulation<<<{"description": "Out-of-plane stress/strain formulation"}>>> = PLANE_STRAIN
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'grid_thermal_eigenstrain grid_irradiation_eigenstrain'
    generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'vonmises_stress stress_xx stress_yy stress_zz'
    decomposition_method<<<{"description": "Methods to calculate the finite strain and rotation increments"}>>> = EigenSolution
    temperature<<<{"description": "The temperature"}>>> = temp
    extra_vector_tags<<<{"description": "The tag names for extra vectors that residual data should be saved into"}>>> = 'ref'
  []
[]

[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
  [heat] # gradient term in heat conduction equation
    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"}>>> = temp
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'pellet_type_1 clad grid'
    extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = 'ref'
  []
  [heat_ie] # time term in heat conduction equation
    type = HeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the thermal energy conservation equation.", "href": "../source/kernels/HeatConductionTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temp
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'pellet_type_1 clad'
    extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = 'ref'
  []
  [heat_source] # source term in heat conduction equation
    type = NeutronHeatSource<<<{"description": "Compute heat generation due to fission.", "href": "../source/kernels/NeutronHeatSource.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temp
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    fission_rate<<<{"description": "Coupled Fission Rate"}>>> = fission_rate
    extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = 'ref'
  []
[]

[Contact<<<{"href": "../syntax/Contact/index.html"}>>>]
  # Define mechanical contact between the fuel (sideset=10) and the clad (sideset=5)
  [spacer_clad_mechanical]
    formulation<<<{"description": "The contact formulation"}>>> = mortar
    model<<<{"description": "The contact model to use"}>>> = coulomb
    primary<<<{"description": "The list of boundary IDs referring to primary sidesets"}>>> = 101
    secondary<<<{"description": "The list of boundary IDs referring to secondary sidesets"}>>> = 102
    c_normal<<<{"description": "Parameter for balancing the size of the gap and contact pressure for a mortar formulation. This purely numerical parameter affects convergence behavior and, in general, should be larger for stiffer materials. It is recommended that the user tries out various orders of magnitude for this parameter if the default value generates poor contact convergence."}>>> = 1e+16 # 1e+7
    c_tangential<<<{"description": "Numerical parameter for nonlinear mortar frictional constraints"}>>> = 1e+20
    friction_coefficient<<<{"description": "The friction coefficient"}>>> = 0.4
    # Do not apply dynamic stabilization
    newmark_beta<<<{"description": "Newmark-beta beta parameter for its inclusion in the weighted gap update formula"}>>> = 0.0001
    newmark_gamma<<<{"description": "Newmark-beta gamma parameter for its inclusion in the weighted gap update formula"}>>> = 0.5
    capture_tolerance<<<{"description": "Normal distance from surface within which nodes are captured. This parameter is used for node-face and mortar formulations."}>>> = 0.0
    mortar_dynamics<<<{"description": "Whether to use constraints that account for the persistency condition, giving rise to smoother normal contact pressure evolution. This flag should only be set to yes for dynamic simulations using the Newmark-beta numerical integrator"}>>> = true
    interpolate_normals = false
    generate_mortar_mesh<<<{"description": "Whether to generate the mortar mesh from the action. Typically this will be the case, but one may also want to reuse an existing lower-dimensional mesh prior to a restart."}>>> = true
    wear_depth<<<{"description": "The name of the mortar auxiliary variable that is used to modify the weighted gap definition"}>>> = worn_depth
  []
[]

[Contact]
  # Define mechanical contact between the fuel (sideset=10) and the clad (sideset=5)
  [pellet_clad_mechanical_real]
    formulation<<<{"description": "The contact formulation"}>>> = mortar
    model<<<{"description": "The contact model to use"}>>> = frictionless
    primary<<<{"description": "The list of boundary IDs referring to primary sidesets"}>>> = 7
    secondary<<<{"description": "The list of boundary IDs referring to secondary sidesets"}>>> = 8
    c_normal<<<{"description": "Parameter for balancing the size of the gap and contact pressure for a mortar formulation. This purely numerical parameter affects convergence behavior and, in general, should be larger for stiffer materials. It is recommended that the user tries out various orders of magnitude for this parameter if the default value generates poor contact convergence."}>>> = 1e+16 #
    c_tangential<<<{"description": "Numerical parameter for nonlinear mortar frictional constraints"}>>> = 1e+16
    friction_coefficient<<<{"description": "The friction coefficient"}>>> = 0.4
    # Do not apply dynamic stabilization
    newmark_beta<<<{"description": "Newmark-beta beta parameter for its inclusion in the weighted gap update formula"}>>> = 0.0001
    newmark_gamma<<<{"description": "Newmark-beta gamma parameter for its inclusion in the weighted gap update formula"}>>> = 0.5
    capture_tolerance<<<{"description": "Normal distance from surface within which nodes are captured. This parameter is used for node-face and mortar formulations."}>>> = 0.0
    mortar_dynamics<<<{"description": "Whether to use constraints that account for the persistency condition, giving rise to smoother normal contact pressure evolution. This flag should only be set to yes for dynamic simulations using the Newmark-beta numerical integrator"}>>> = true
    interpolate_normals = false
    generate_mortar_mesh<<<{"description": "Whether to generate the mortar mesh from the action. Typically this will be the case, but one may also want to reuse an existing lower-dimensional mesh prior to a restart."}>>> = true
  []
[]

[ThermalContactMortar<<<{"href": "../syntax/ThermalContactMortar/index.html"}>>>]
  [thermal_contact]
    secondary_variable<<<{"description": "The secondary variable"}>>> = temp
    primary_boundary<<<{"description": "The primary surface"}>>> = 7
    secondary_boundary<<<{"description": "The secondary surface"}>>> = 8
    initial_moles<<<{"description": "The Postprocessor that will give the initial moles of gas."}>>> = initial_moles # coupling to a postprocessor which supplies the initial plenum/gap gas mass
    gas_released<<<{"description": "The postprocessor(s) that will give the gas released."}>>> = fission_gas_released # coupling to a postprocessor which supplies the fission gas addition
  []
[]

[Burnup<<<{"href": "../syntax/Burnup/index.html"}>>>]
  [burnup]
    block<<<{"description": "The blocks where radial power factor should be computed."}>>> = pellet_type_1
    rod_ave_lin_pow<<<{"description": "Rod average linear power function."}>>> = power_history # using the power function defined above
    axial_power_profile<<<{"description": "Axial power peaking function."}>>> = axial_peaking_factors # using the axial power profile function defined above
    num_radial<<<{"description": "Number of radial divisions in secondary grid used to compute radial power profile."}>>> = 80
    num_axial<<<{"description": "Number of axial divisions in secondary grid used to compute radial power profile."}>>> = 21
    axial_axis<<<{"description": "Coordinate axis of the axial direction of the fuel stack (0, 1, or 2 for x, y, or z)"}>>> = 2
    density<<<{"description": "The initial fuel density."}>>> = ${initial_fuel_density}
    a_lower<<<{"description": "The lower axial coordinate of the fuel stack. Required if fuel_pin_geometry is not specified."}>>> = -1e-3 # mesh dependent!
    a_upper<<<{"description": "The upper axial coordinate of the fuel stack. Required if fuel_pin_geometry is not specified."}>>> = 1e-3 # mesh dependent!
    fuel_inner_radius<<<{"description": "The inner radius of the fuel."}>>> = 0
    fuel_outer_radius<<<{"description": "The outer radius of the fuel."}>>> = .0041
    fuel_volume_ratio<<<{"description": "Reduction factor for deviation from right circular cylinder fuel.  The ratio of actual volume to right circular cylinder volume."}>>> = 0.987775 # for use with dished pellets (ratio of actual volume to cylinder volume)

    #N235 = N235 # Activate to write N235 concentration to output file
    #N238 = N238 # Activate to write N238 concentration to output file
    #N239 = N239 # Activate to write N239 concentration to output file
    #N240 = N240 # Activate to write N240 concentration to output file
    #N241 = N241 # Activate to write N241 concentration to output file
    #N242 = N242 # Activate to write N242 concentration to output file
    RPF<<<{"description": "Specifies that the radial power factor is required."}>>> = RPF
  []
[]

[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
  [worn_depth]
    type = MortarArchardsLawAux<<<{"description": "Returns the weighted gap velocity at a node. This quantity is useful for mortar contact, particularly when dual basis functions are used in contact mechanics", "href": "../source/auxkernels/MortarArchardsLawAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = worn_depth
    primary_boundary<<<{"description": "The name of the primary boundary sideset."}>>> = 101
    secondary_boundary<<<{"description": "The name of the secondary boundary sideset."}>>> = 102
    primary_subdomain<<<{"description": "The name of the primary subdomain."}>>> = 'spacer_clad_mechanical_primary_subdomain'
    secondary_subdomain<<<{"description": "The name of the secondary subdomain."}>>> = 'spacer_clad_mechanical_secondary_subdomain'
    displacements<<<{"description": "The displacement variables. This mortar nodal auxiliary kernel can take two or three displacements"}>>> = 'disp_x disp_y'
    friction_coefficient<<<{"description": "Friction coefficient used to compute wear (to match that of frictional contact)"}>>> = 0.5
    energy_wear_coefficient<<<{"description": "Energy wear coefficient is a surface-dependent parameter used in Archard's wear law"}>>> = 0.1e-9
    normal_pressure<<<{"description": "The name of the Lagrange multiplier that holds the normal contact pressure in mortar formulations"}>>> = spacer_clad_mechanical_normal_lm
    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'
  []
  [fast_neutron_flux]
    type = FastNeutronFluxAux<<<{"description": "Computes fast neutron flux.", "href": "../source/auxkernels/FastNeutronFluxAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = fast_neutron_flux
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    rod_ave_lin_pow<<<{"description": "The power history function."}>>> = power_history
    axial_power_profile<<<{"description": "The axial power profile function."}>>> = axial_peaking_factors
    factor<<<{"description": "The fast neutron flux if the function, RALP, and q_variable are not given. A scaling factor if the function, RALP, or q_variable is given.  If RALP or q_variable is given, it is recommended to use a value of 3e13 (n/(m^2s)/(W/m))."}>>> = 3e13
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_begin
  []
  [fast_neutron_fluence]
    type = FastNeutronFluenceAux<<<{"description": "Computes fast neutron fluence from fast_neutron_flux.", "href": "../source/auxkernels/FastNeutronFluenceAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = fast_neutron_fluence
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    fast_neutron_flux<<<{"description": "Coupled Fast Neutron Flux."}>>> = fast_neutron_flux
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_begin
  []
  [relocation_strain]
    type = MaterialRealAux<<<{"description": "Outputs element volume-averaged material properties", "href": "../source/auxkernels/MaterialRealAux.html"}>>>
    property<<<{"description": "The material property name."}>>> = relocation_strain
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = relocation_strain
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
[]

[BCs<<<{"href": "../syntax/NuclearMaterials/BCs/index.html"}>>>]
  # Define boundary conditions
  [no_y_all] # pin pellets and clad along axis of symmetry (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"}>>> = 15
    value<<<{"description": "Value of the BC"}>>> = 0.0
  []
  [no_x_all] # pin pellets and clad along axis of symmetry (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"}>>> = 16
    value<<<{"description": "Value of the BC"}>>> = 0.0
  []
  [no_y_all_grid] # pin pellets and clad along axis of symmetry (y)
    type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.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"}>>> = '112'
    function<<<{"description": "The forcing function."}>>> = 'if(t < 1.0e4,1.0e-4 * t/1.0e4 - 1.0e-5,0.9e-4)'
  []
  [no_x_all_grid] # pin pellets and clad along axis of symmetry (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"}>>> = '112'
    value<<<{"description": "Value of the BC"}>>> = 0.0
  []
  [Pressure<<<{"href": "../syntax/BCs/Pressure/index.html"}>>>] #  apply coolant pressure on clad outer walls
    [coolantPressure]
      boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '2'
      factor<<<{"description": "The magnitude to use in computing the pressure"}>>> = 15.5e6
      function<<<{"description": "The function that describes the pressure"}>>> = pressure_var # use the pressure_ramp function defined above
    []
  []
  [PlenumPressure<<<{"href": "../syntax/BCs/PlenumPressure/index.html"}>>>] #  apply plenum pressure on clad inner walls and pellet surfaces
    [plenumPressure]
      boundary<<<{"description": "The list of boundary IDs from the mesh where the pressure will be applied"}>>> = 9
      initial_pressure<<<{"description": "The initial pressure in the cavity.  If not given, a zero initial pressure will be used."}>>> = 2.0e6
      R<<<{"description": "The universal gas constant for the units used."}>>> = 8.3143
      output_initial_moles<<<{"description": "The name to use when reporting the initial moles of gas"}>>> = initial_moles # coupling to post processor to get initial fill gas mass
      temperature<<<{"description": "The name of the average temperature postprocessor value."}>>> = plenum_temperature # coupling to post processor to get gas temperature approximation
      volume<<<{"description": "The name of the postprocessor(s) that holds the value of the internal volume in the cavity"}>>> = plenum_volume # coupling to post processor to get gas volume
      material_input<<<{"description": "The name of the postprocessor(s) that holds the amount of material injected into the plenum."}>>> = fission_gas_released # coupling to post processor to get fission gas added
      output<<<{"description": "The name to use for the cavity pressure value"}>>> = plenum_pressure # coupling to post processor to output plenum/gap pressure
      displacements<<<{"description": "The nonlinear displacement variables"}>>> = 'disp_x disp_y'
    []
  []
  [convective_clad_surface] # apply convective boundary to clad outer surface
    type = ConvectiveFluxBC<<<{"description": "Determines boundary values via the initial and final values, flux, and exposure duration", "href": "../source/bcs/ConvectiveFluxBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temp
    rate = 38200.0 #convection coefficient (h)
    initial = 580.0
    final = 580.0
    duration = 1.0e4 #duration of initial power ramp
  []
[]

[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
  # Define material behavior models and input material property data
  [fuel_thermal] # temperature and burnup dependent thermal properties of UO2 (BISON kernel)
    type = UO2Thermal<<<{"description": "Computes thermal conductivity and specific heat capacity of uranium dioxide fuel.", "href": "../source/materials/UO2Thermal.html"}>>>
    thermal_conductivity_model<<<{"description": "The thermal conductivity model."}>>> = FINK_LUCUTA
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    burnup<<<{"description": "Coupled burnup"}>>> = burnup
    initial_porosity<<<{"description": "Initial porosity.  Must be between 0.0 and 1.0."}>>> = 0.0
  []
  [fuel_solid_mechanics_swelling] # free expansion strains (swelling and densification) for UO2 (BISON kernel)
    type = UO2VolumetricSwellingEigenstrain<<<{"description": "Computes and sums the change in fuel pellet volume due to densification and fission product release. This class applies a volumetric strain correction before adding the strain from this class to the diagonal entries of the eigenstrain tensor.", "href": "../source/materials/solid_mechanics/UO2VolumetricSwellingEigenstrain.html"}>>>
    gas_swelling_model_type<<<{"description": "Which type of model to use to calculate the gaseous swelling. Choices are SIFGRS MATPRO SIFGRS_IG. If you select SIFGRS or SIFGRS_IG, the SIFGRS model must be included in the input file."}>>> = MATPRO
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    burnup<<<{"description": "Coupled Burnup"}>>> = burnup
    initial_fuel_density<<<{"description": "Initial fuel density in kg-UO2/m^3"}>>> = 10431.0
    temperature<<<{"description": "Coupled Temperature in Kelvin"}>>> = temp
    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."}>>> = 'fuel_volumetric_eigenstrain'
  []
  [fuel_creep]
    type = UO2CreepUpdate<<<{"description": "Computes the secondary thermal and irradiation creep for UO2 LWR fuel. This material must be run in conjunction with ComputeMultipleInelasticStress.", "href": "../source/materials/solid_mechanics/UO2CreepUpdate.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    fission_rate<<<{"description": "Coupled fission rate"}>>> = fission_rate
    density<<<{"description": "Initial fuel density"}>>> = 10431.0

    initial_grain_radius<<<{"description": "Fuel grain radius (m)"}>>> = 10.0e-6
    oxygen_to_metal_ratio<<<{"description": "Oxygen to metal ratio"}>>> = 2.0
  []
  [fuel_elasticity_tensor]
    type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'pellet_type_1'
    youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 906e6
    poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.345
  []
  [fuel_stress]
    type = ComputeMultipleInelasticStress<<<{"description": "Compute state (stress and internal parameters such as plastic strains and internal parameters) using an iterative process.  Combinations of creep models and plastic models may be used.", "href": "../source/materials/ComputeMultipleInelasticStress.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    inelastic_models<<<{"description": "The material objects to use to calculate stress and inelastic strains. Note: specify creep models first and plasticity models second."}>>> = 'fuel_creep'
  []
  [fuel_thermal_expansion]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 10.0e-6
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 580.0
    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."}>>> = 'fuel_thermal_eigenstrain'
  []
  [fuel_relocation]
    type = UO2RelocationEigenstrain<<<{"description": "Accounts for cracking and relocation of fuel pellet fragments in the radial direction and is necessary for accurate modeling of LWR fuel. The linear heat rate must either be a functionor a variable.", "href": "../source/materials/solid_mechanics/UO2RelocationEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    burnup<<<{"description": "Coupled Burnup variable in FIMA"}>>> = burnup
    diameter<<<{"description": "As fabricated cold diameter of pellet in meters"}>>> = 0.0082
    rod_ave_lin_pow<<<{"description": "linear power function."}>>> = power_history
    axial_power_profile<<<{"description": "axial power peaking function."}>>> = axial_peaking_factors
    diametral_gap<<<{"description": "As fabricated cold diametral pellet-cladding gap in meters."}>>> =160e-6
    burnup_relocation_stop<<<{"description": "Burnup at which relocation strain stops (in FIMA)"}>>> = 1.e20
    relocation_activation1<<<{"description": "First activation linear power in W/m term in relocation model"}>>> = 5000
    axial_axis = 2
    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."}>>> = 'fuel_relocation_eigenstrain'
  []
  [clad_thermal]
    type = HeatConductionMaterial<<<{"description": "General-purpose material model for heat conduction", "href": "../source/materials/HeatConductionMaterial.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    thermal_conductivity<<<{"description": "The thermal conductivity value"}>>> = 16.0
    specific_heat<<<{"description": "The specific heat value"}>>> = 330.0
  []
  [clad_elasticity_tensor]
    type = ZryElasticityTensor<<<{"description": "Either provides constant elasticity constants for Zircaloy cladding or calculates the Young's modulus and Poisson's ratio for Zircaloy cladding using MATPRO relations as a function of temperature and fast neutron fluence.", "href": "../source/materials/solid_mechanics/ZryElasticityTensor.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
  []
  [clad_creep_model]
    type = ZryCreepHayesHoppeUpdate<<<{"description": "Computes the secondary thermal Hayes and Kassner secondary creep and the Hoppe irradiation creep for Zircaloy cladding.  This material must be run in conjunction with ComputeMultipleInelasticStress.", "href": "../source/materials/solid_mechanics/ZryCreepHayesHoppeUpdate.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    fast_neutron_flux<<<{"description": "The fast neutron flux"}>>> = fast_neutron_flux
    temperature<<<{"description": "The coupled temperature (K)"}>>> = temp
    zircaloy_material_type<<<{"description": "Type of zircaloy material properties to use in calculating creep. Note: ESCORE_IRRADIATIONGROWTHZR4 is not valid."}>>> = stress_relief_annealed
    model_irradiation_creep<<<{"description": "Set true to activate irradiation induced creep"}>>> = true
    model_thermal_creep<<<{"description": "Set true to activate steady state thermal creep"}>>> = true
  []
  [clad_stress]
    type = ComputeMultipleInelasticStress<<<{"description": "Compute state (stress and internal parameters such as plastic strains and internal parameters) using an iterative process.  Combinations of creep models and plastic models may be used.", "href": "../source/materials/ComputeMultipleInelasticStress.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    tangent_operator<<<{"description": "Type of tangent operator to return.  'elastic': return the elasticity tensor.  'nonlinear': return the full, general consistent tangent operator."}>>> = elastic
    inelastic_models<<<{"description": "The material objects to use to calculate stress and inelastic strains. Note: specify creep models first and plasticity models second."}>>> = 'clad_creep_model'
  []
  [clad_thermal_expansion]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 5.0e-6
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 580.0
    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."}>>> = 'clad_thermal_eigenstrain'
  []
  [clad_irrgrowth]
    type = ZryIrradiationGrowthEigenstrain<<<{"description": "Computes eigenstrain from irradiation growth in Zircaloy cladding using either the Franklin or ESCORE models.", "href": "../source/materials/solid_mechanics/ZryIrradiationGrowthEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    fast_neutron_fluence<<<{"description": "The fast neutron fluence"}>>> = fast_neutron_fluence
    axial_direction<<<{"description": "The direction you want to apply irradiation growth"}>>> = 2
    zircaloy_material_type<<<{"description": "Type of irradiation growth volumetric swelling formulation."}>>> = ESCORE_IrradiationGrowthZr4
    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."}>>> = 'clad_irradiation_eigenstrain'
  []

  [grid_thermal]
    type = HeatConductionMaterial<<<{"description": "General-purpose material model for heat conduction", "href": "../source/materials/HeatConductionMaterial.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    thermal_conductivity<<<{"description": "The thermal conductivity value"}>>> = 16.0
    specific_heat<<<{"description": "The specific heat value"}>>> = 330.0
  []
  [grid_elasticity_tensor]
    type = ZryElasticityTensor<<<{"description": "Either provides constant elasticity constants for Zircaloy cladding or calculates the Young's modulus and Poisson's ratio for Zircaloy cladding using MATPRO relations as a function of temperature and fast neutron fluence.", "href": "../source/materials/solid_mechanics/ZryElasticityTensor.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
  []
  [grid_creep_model]
    type = ZryCreepHayesHoppeUpdate<<<{"description": "Computes the secondary thermal Hayes and Kassner secondary creep and the Hoppe irradiation creep for Zircaloy cladding.  This material must be run in conjunction with ComputeMultipleInelasticStress.", "href": "../source/materials/solid_mechanics/ZryCreepHayesHoppeUpdate.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    fast_neutron_flux<<<{"description": "The fast neutron flux"}>>> = fast_neutron_flux
    temperature<<<{"description": "The coupled temperature (K)"}>>> = temp
    zircaloy_material_type<<<{"description": "Type of zircaloy material properties to use in calculating creep. Note: ESCORE_IRRADIATIONGROWTHZR4 is not valid."}>>> = stress_relief_annealed
    model_irradiation_creep<<<{"description": "Set true to activate irradiation induced creep"}>>> = true
    model_thermal_creep<<<{"description": "Set true to activate steady state thermal creep"}>>> = true
  []
  [grid_stress]
    type = ComputeMultipleInelasticStress<<<{"description": "Compute state (stress and internal parameters such as plastic strains and internal parameters) using an iterative process.  Combinations of creep models and plastic models may be used.", "href": "../source/materials/ComputeMultipleInelasticStress.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    tangent_operator<<<{"description": "Type of tangent operator to return.  'elastic': return the elasticity tensor.  'nonlinear': return the full, general consistent tangent operator."}>>> = elastic
    inelastic_models<<<{"description": "The material objects to use to calculate stress and inelastic strains. Note: specify creep models first and plasticity models second."}>>> = 'grid_creep_model'
  []
  [grid_thermal_expansion]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 5.0e-6
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 580.0
    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."}>>> = 'grid_thermal_eigenstrain'
  []
  [grid_irrgrowth]
    type = ZryIrradiationGrowthEigenstrain<<<{"description": "Computes eigenstrain from irradiation growth in Zircaloy cladding using either the Franklin or ESCORE models.", "href": "../source/materials/solid_mechanics/ZryIrradiationGrowthEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = grid
    fast_neutron_fluence<<<{"description": "The fast neutron fluence"}>>> = fast_neutron_fluence
    axial_direction<<<{"description": "The direction you want to apply irradiation growth"}>>> = 2
    zircaloy_material_type<<<{"description": "Type of irradiation growth volumetric swelling formulation."}>>> = ESCORE_IrradiationGrowthZr4
    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."}>>> = 'grid_irradiation_eigenstrain'
  []

  [fission_gas_release] # Forsberg-Massih fission gas release mode
    type = UO2Sifgrs<<<{"description": "Recommended fission gas model to account for generation of fission gasses in nuclear fuel", "href": "../source/materials/UO2Sifgrs.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    fission_rate<<<{"description": "Coupled fission rate variable (fiss/m^3/s)"}>>> = fission_rate # coupling to fission_rate aux variable
    grain_radius<<<{"description": "Coupled grain Radius"}>>> = 10.0e-6
    #external_pressure = 40e6
  []
  [clad_density]
    type = StrainAdjustedDensity<<<{"description": "Creates density material property", "href": "../source/materials/StrainAdjustedDensity.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    density = 6551.0
  []
  [fuel_density]
    type = StrainAdjustedDensity<<<{"description": "Creates density material property", "href": "../source/materials/StrainAdjustedDensity.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet
    strain_free_density<<<{"description": "Material property for strain-free density"}>>> = 10431.0
  []
  [grid]
    type = StrainAdjustedDensity<<<{"description": "Creates density material property", "href": "../source/materials/StrainAdjustedDensity.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = grid
    density = 6560
  []
[]

[Debug<<<{"href": "../syntax/Debug/index.html"}>>>]
  show_var_residual_norms<<<{"description": "Print the residual norms of the individual solution variables at each nonlinear iteration"}>>> = true
[]

[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
  type = Transient
  solve_type = 'NEWTON'
  petsc_options = '-snes_converged_reason -ksp_converged_reason'
  petsc_options_iname = '-pc_type -pc_factor_mat_solver_package -mat_mffd_err -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = 'lu    superlu_dist   1e-6          NONZERO               1e-10'
  snesmf_reuse_base = true

  line_search = 'none'

  l_max_its = 100
  l_tol = 8e-3

  nl_max_its = 45
  nl_rel_tol = 1e-10 # was -7 and nl 25. Tightening tangential contact forces.
  nl_abs_tol = 1e-12

  [TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
    type = NewmarkBeta
    beta = 0.25
    gamma = 0.5
  []

  start_time = 0.0
  end_time = 1.0e5

  [TimeStepper<<<{"href": "../syntax/Executioner/TimeStepper/index.html"}>>>]
    type = IterationAdaptiveDT
    dt = 2.0e2
    time_t = '1e4 1e5 1e6'
    time_dt = '2e2 1e4 1e5'
    growth_factor = 1.4
    iteration_window = 5.0
    optimal_iterations = 35
  []

  dtmax = 2e5 # Larger causes instabilities 2e6
  dtmin = 1
[]

[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
  # Define postprocessors (some are required as specified above; others are optional; many others are available)
  [average_interior_clad_temperature] # average temperature of cladding interior
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../source/postprocessors/SideAverageValue.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 7
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    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 timestep_end'
  []
  [average_centerline_fuel_temperature] # average temperature of the cladding interior and all pellet exteriors
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../source/postprocessors/SideAverageValue.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 9
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    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 linear'
  []
  [plenum_temperature]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../source/postprocessors/SideAverageValue.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 9
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    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 timestep_end'
  []
  [plenum_volume] # gas volume
    type = InternalVolume<<<{"description": "Computes the volume of an enclosed area by performing an integral over a user-supplied boundary.", "href": "../source/postprocessors/InternalVolume.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 9
    addition<<<{"description": "An additional volume to be included in the internal volume calculation. A time-dependent function is expected."}>>> = 1.3e-5 #rough guess of plenum volume/unit length of fuel
    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 linear'
  []
  [pellet_volume] # fuel pellet total volume
    type = InternalVolume<<<{"description": "Computes the volume of an enclosed area by performing an integral over a user-supplied boundary.", "href": "../source/postprocessors/InternalVolume.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 8
    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 timestep_end'
  []
  [clad_inner_vol] # volume inside of cladding
    type = InternalVolume<<<{"description": "Computes the volume of an enclosed area by performing an integral over a user-supplied boundary.", "href": "../source/postprocessors/InternalVolume.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 7
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = exodus
    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 timestep_end'
  []
  [fission_gas_generated] # fission gas produced (moles)
    type = ElementIntegralFisGasGeneratedSifgrs<<<{"description": "Reports the fission gas that is produced in moles.  To be used in combination with the Sifgrs model.", "href": "../source/postprocessors/ElementIntegralFisGasGeneratedSifgrs.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
  [fission_gas_released] # fission gas released to plenum (moles)
    type = ElementIntegralFisGasReleasedSifgrs<<<{"description": "Reports the fission gas that is released to the plenum in moles.  To be used in combination with the Sifgrs model.", "href": "../source/postprocessors/ElementIntegralFisGasReleasedSifgrs.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
  [flux_from_clad] # area integrated heat flux from the cladding
    type = SideDiffusiveFluxIntegral<<<{"description": "Computes the integral of the diffusive flux over the specified boundary", "href": "../source/postprocessors/SideDiffusiveFluxIntegral.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 5
    diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation. This must be provided if the variable is of finite element type"}>>> = thermal_conductivity
    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
  []
  [flux_from_fuel] # area integrated heat flux from the fuel
    type = SideDiffusiveFluxIntegral<<<{"description": "Computes the integral of the diffusive flux over the specified boundary", "href": "../source/postprocessors/SideDiffusiveFluxIntegral.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 10
    diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation. This must be provided if the variable is of finite element type"}>>> = thermal_conductivity
    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
  []
  [_dt] # time step
    type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../source/postprocessors/TimestepSize.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
  []
  [num_lin_it]
    type = NumLinearIterations<<<{"description": "Compute the number of linear iterations.", "href": "../source/postprocessors/NumLinearIterations.html"}>>>
  []
  [num_nonlin_it]
    type = NumNonlinearIterations<<<{"description": "Outputs the number of nonlinear iterations", "href": "../source/postprocessors/NumNonlinearIterations.html"}>>>
  []
  [tot_lin_it]
    type = CumulativeValuePostprocessor<<<{"description": "Creates a cumulative sum of a Postprocessor value with time.", "href": "../source/postprocessors/CumulativeValuePostprocessor.html"}>>>
    postprocessor<<<{"description": "The name of the postprocessor"}>>> = num_lin_it
  []
  [tot_nonlin_it]
    type = CumulativeValuePostprocessor<<<{"description": "Creates a cumulative sum of a Postprocessor value with time.", "href": "../source/postprocessors/CumulativeValuePostprocessor.html"}>>>
    postprocessor<<<{"description": "The name of the postprocessor"}>>> = num_nonlin_it
  []
  [alive_time]
    type = PerfGraphData<<<{"description": "Retrieves performance information about a section from the PerfGraph.", "href": "../source/postprocessors/PerfGraphData.html"}>>>
    section_name<<<{"description": "The name of the section to get data for"}>>> = Root
    data_type<<<{"description": "The type of data to retrieve for the section_name"}>>> = TOTAL
  []
  [rod_total_power]
    type = ElementIntegralPower<<<{"description": "Computes the power given the fission rate and energy per fission.", "href": "../source/postprocessors/ElementIntegralPower.html"}>>>
    variable<<<{"description": "The name of the variable that this object operates on"}>>> = temp
    fission_rate<<<{"description": "Coupled fission rate"}>>> = fission_rate
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
  [rod_input_power]
    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."}>>> = power_history
    scale_factor<<<{"description": "A scale factor to be applied to the function"}>>> = 0.1186 # rod height
    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
  []
  [fission_gas_released_percentage]
    type = FGRPercent<<<{"description": "Computes the ratio of fission gas released to the fission gas generated in percent.", "href": "../source/postprocessors/FGRPercent.html"}>>>
    fission_gas_released<<<{"description": "The fission gas released postprocessor"}>>> = fission_gas_released
    fission_gas_generated<<<{"description": "The fission gas generated postprocessor"}>>> = fission_gas_generated
  []
[]

[VectorPostprocessors<<<{"href": "../syntax/VectorPostprocessors/index.html"}>>>]
  [contact_pressure]
    type = NodalValueSampler<<<{"description": "Samples values of nodal variable(s).", "href": "../source/vectorpostprocessors/NodalValueSampler.html"}>>>
    sort_by<<<{"description": "What to sort the samples by"}>>> = x
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = spacer_clad_mechanical_normal_lm
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 102
  []
  [frictional_pressure]
    type = NodalValueSampler<<<{"description": "Samples values of nodal variable(s).", "href": "../source/vectorpostprocessors/NodalValueSampler.html"}>>>
    sort_by<<<{"description": "What to sort the samples by"}>>> = x
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = spacer_clad_mechanical_tangential_lm
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 102
  []
  [worn_depth]
    type = NodalValueSampler<<<{"description": "Samples values of nodal variable(s).", "href": "../source/vectorpostprocessors/NodalValueSampler.html"}>>>
    sort_by<<<{"description": "What to sort the samples by"}>>> = x
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = worn_depth
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 102
    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
  []
[]

[Outputs<<<{"href": "../syntax/Outputs/index.html"}>>>]
  perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = true
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true

  [console]
    type = Console<<<{"description": "Object for screen output.", "href": "../source/outputs/Console.html"}>>>
    max_rows<<<{"description": "The maximum number of postprocessor/scalar values displayed on screen during a timestep (set to 0 for unlimited)"}>>> = 25
  []
  checkpoint<<<{"description": "Create checkpoint files using the default options."}>>> = true
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial.i)
user_start_time = 1.0e5
user_end_time = 1.000002e5

end_dynamic_excitation = 1.000002e5
time_step_dynamics = 2.0e-3
step_number = 1

initial_fuel_density = 10431.0

[GlobalParams<<<{"href": "../syntax/GlobalParams/index.html"}>>>]
  temperature = temp
  displacements = 'disp_x disp_y'
  order = FIRST
  family = LAGRANGE
  energy_per_fission = 3.2e-11 # J/fission
  volumetric_locking_correction = true
[]

[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
  [file]
    type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../source/meshgenerators/FileMeshGenerator.html"}>>>
    file<<<{"description": "The filename to read."}>>> = fretting-wear-initial_out_cp/LATEST
    skip_partitioning<<<{"description": "True to skip partitioning, only after this mesh generator, because the mesh was pre-split for example."}>>> = true
    allow_renumbering<<<{"description": "Whether to allow the mesh to renumber nodes and elements, if not overridden by a global parameter or by a requirement (e.g. an exodus restart or a constraint matrix) that disables renumbering."}>>> = false
  []
  patch_size = 100 # For contact algorithm
[]

[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
  [temp]
  []
[]

[Problem<<<{"href": "../syntax/Problem/index.html"}>>>]
  type = ReferenceResidualProblem
  reference_vector = 'ref'
  extra_tag_vectors = 'ref'
  group_variables = 'disp_x disp_y '
  converge_on = 'disp_x disp_y  temp'
  restart_file_base = ./fretting-wear-initial_out_cp/LATEST
  material_coverage_check = false
  kernel_coverage_check = false
[]

[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
  [fission_rate]
    block = pellet_type_1
  []
  [burnup]
    block = pellet_type_1
  []
  [fast_neutron_flux]
    block = 'clad grid'
  []
  [fast_neutron_fluence]
    block = 'clad grid'
  []
  [relocation_strain]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [worn_depth]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
    block = 'spacer_clad_mechanical_secondary_subdomain'
  []
[]

[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
  [power_history]
    type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>> # reads and interpolates an input file containing rod average linear power vs time
    data_file<<<{"description": "File holding CSV data"}>>> = powerhistory.csv
    scale_factor<<<{"description": "Scale factor to be applied to the ordinate values"}>>> = 1
  []
  [axial_peaking_factors]
    type = ConstantFunction<<<{"description": "A function that returns a constant value as defined by an input parameter.", "href": "../source/functions/ConstantFunction.html"}>>>
    value<<<{"description": "The constant value"}>>> = 1
  []
  [pressure_var] # reads and interpolates input data defining amplitude curve for fill gas pressure
    type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
    x<<<{"description": "The abscissa values"}>>> = '0 1e4'
    y<<<{"description": "The ordinate values"}>>> = '0 1'
  []
  [pressure_var_variable] # reads and interpolates input data defining amplitude curve for fill gas pressure
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = 'if(t < 1e4, 1, 1 + sin((t-1e4)*pi/10.0) * (t-1e4))'
  []
[]

[Physics<<<{"href": "../syntax/Physics/index.html"}>>>/SolidMechanics<<<{"href": "../syntax/Physics/SolidMechanics/index.html"}>>>/Dynamic<<<{"href": "../syntax/Physics/SolidMechanics/Dynamic/index.html"}>>>]
  [pellets]
    add_variables<<<{"description": "Add the displacement variables"}>>> = true
    newmark_beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
    newmark_gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
    block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = pellet_type_1
    strain<<<{"description": "Strain formulation"}>>> = FINITE
    planar_formulation<<<{"description": "Out-of-plane stress/strain formulation"}>>> = PLANE_STRAIN
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'fuel_relocation_eigenstrain fuel_thermal_eigenstrain
      fuel_volumetric_eigenstrain'
    generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'vonmises_stress stress_xx stress_yy stress_zz'
    decomposition_method<<<{"description": "Methods to calculate the finite strain and rotation increments"}>>> = EigenSolution
    temperature<<<{"description": "The temperature"}>>> = temp
    extra_vector_tags<<<{"description": "The tag names for extra vectors that residual data should be saved into"}>>> = 'ref'
  []
  [clad]
    add_variables<<<{"description": "Add the displacement variables"}>>> = true
    newmark_beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
    newmark_gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
    block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = clad
    strain<<<{"description": "Strain formulation"}>>> = FINITE
    planar_formulation<<<{"description": "Out-of-plane stress/strain formulation"}>>> = PLANE_STRAIN
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'clad_thermal_eigenstrain clad_irradiation_eigenstrain'
    generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'vonmises_stress stress_xx stress_yy stress_zz'
    decomposition_method<<<{"description": "Methods to calculate the finite strain and rotation increments"}>>> = EigenSolution
    temperature<<<{"description": "The temperature"}>>> = temp
    extra_vector_tags<<<{"description": "The tag names for extra vectors that residual data should be saved into"}>>> = 'ref'
  []
  [grid]
    add_variables<<<{"description": "Add the displacement variables"}>>> = true
    newmark_beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
    newmark_gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
    block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = grid
    strain<<<{"description": "Strain formulation"}>>> = FINITE
    planar_formulation<<<{"description": "Out-of-plane stress/strain formulation"}>>> = PLANE_STRAIN
    eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'grid_thermal_eigenstrain grid_irradiation_eigenstrain'
    generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'vonmises_stress stress_xx stress_yy stress_zz'
    decomposition_method<<<{"description": "Methods to calculate the finite strain and rotation increments"}>>> = EigenSolution
    temperature<<<{"description": "The temperature"}>>> = temp
    extra_vector_tags<<<{"description": "The tag names for extra vectors that residual data should be saved into"}>>> = 'ref'
  []
[]

[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
  [heat] # gradient term in heat conduction equation
    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"}>>> = temp
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'pellet_type_1 clad grid'
    extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = 'ref'
  []
  [heat_ie] # time term in heat conduction equation
    type = HeatConductionTimeDerivative<<<{"description": "Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the thermal energy conservation equation.", "href": "../source/kernels/HeatConductionTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temp
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'pellet_type_1 clad'
    extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = 'ref'
  []
  [heat_source] # source term in heat conduction equation
    type = NeutronHeatSource<<<{"description": "Compute heat generation due to fission.", "href": "../source/kernels/NeutronHeatSource.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temp
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    fission_rate<<<{"description": "Coupled Fission Rate"}>>> = fission_rate
    extra_vector_tags<<<{"description": "The extra tags for the vectors this Kernel should fill"}>>> = 'ref'
  []
[]

[Contact<<<{"href": "../syntax/Contact/index.html"}>>>]
  # Define mechanical contact between the fuel (sideset=10) and the clad (sideset=5)
  [spacer_clad_mechanical]
    formulation<<<{"description": "The contact formulation"}>>> = mortar
    model<<<{"description": "The contact model to use"}>>> = coulomb
    primary<<<{"description": "The list of boundary IDs referring to primary sidesets"}>>> = 101
    secondary<<<{"description": "The list of boundary IDs referring to secondary sidesets"}>>> = 102
    c_normal<<<{"description": "Parameter for balancing the size of the gap and contact pressure for a mortar formulation. This purely numerical parameter affects convergence behavior and, in general, should be larger for stiffer materials. It is recommended that the user tries out various orders of magnitude for this parameter if the default value generates poor contact convergence."}>>> = 1e+12 # 5e13
    c_tangential<<<{"description": "Numerical parameter for nonlinear mortar frictional constraints"}>>> = 1e+18
    friction_coefficient<<<{"description": "The friction coefficient"}>>> = 0.4
    # Do not apply dynamic stabilization
    newmark_beta<<<{"description": "Newmark-beta beta parameter for its inclusion in the weighted gap update formula"}>>> = 0.0001
    newmark_gamma<<<{"description": "Newmark-beta gamma parameter for its inclusion in the weighted gap update formula"}>>> = 0.5
    capture_tolerance<<<{"description": "Normal distance from surface within which nodes are captured. This parameter is used for node-face and mortar formulations."}>>> = 0.0
    mortar_dynamics<<<{"description": "Whether to use constraints that account for the persistency condition, giving rise to smoother normal contact pressure evolution. This flag should only be set to yes for dynamic simulations using the Newmark-beta numerical integrator"}>>> = true
    interpolate_normals = false
    normal_lm_scaling<<<{"description": "Scaling factor to apply to the normal LM variable for a mortar formulation"}>>> = 1.0e-6
    tangential_lm_scaling<<<{"description": "Scaling factor to apply to the tangential LM variable for a mortar formulation"}>>> = 1.0e-6
    generate_mortar_mesh<<<{"description": "Whether to generate the mortar mesh from the action. Typically this will be the case, but one may also want to reuse an existing lower-dimensional mesh prior to a restart."}>>> = false
    wear_depth<<<{"description": "The name of the mortar auxiliary variable that is used to modify the weighted gap definition"}>>> = worn_depth
  []
[]

[Contact]
  # Define mechanical contact between the fuel (sideset=10) and the clad (sideset=5)
  [pellet_clad_mechanical_real]
    formulation<<<{"description": "The contact formulation"}>>> = mortar
    model<<<{"description": "The contact model to use"}>>> = frictionless
    primary<<<{"description": "The list of boundary IDs referring to primary sidesets"}>>> = 7
    secondary<<<{"description": "The list of boundary IDs referring to secondary sidesets"}>>> = 8
    c_normal<<<{"description": "Parameter for balancing the size of the gap and contact pressure for a mortar formulation. This purely numerical parameter affects convergence behavior and, in general, should be larger for stiffer materials. It is recommended that the user tries out various orders of magnitude for this parameter if the default value generates poor contact convergence."}>>> = 1e+16 #
    c_tangential<<<{"description": "Numerical parameter for nonlinear mortar frictional constraints"}>>> = 1e+16
    friction_coefficient<<<{"description": "The friction coefficient"}>>> = 0.4
    # Do not apply dynamic stabilization
    newmark_beta<<<{"description": "Newmark-beta beta parameter for its inclusion in the weighted gap update formula"}>>> = 0.0001
    newmark_gamma<<<{"description": "Newmark-beta gamma parameter for its inclusion in the weighted gap update formula"}>>> = 0.5
    capture_tolerance<<<{"description": "Normal distance from surface within which nodes are captured. This parameter is used for node-face and mortar formulations."}>>> = 0.0
    mortar_dynamics<<<{"description": "Whether to use constraints that account for the persistency condition, giving rise to smoother normal contact pressure evolution. This flag should only be set to yes for dynamic simulations using the Newmark-beta numerical integrator"}>>> = true
    interpolate_normals = false
    generate_mortar_mesh<<<{"description": "Whether to generate the mortar mesh from the action. Typically this will be the case, but one may also want to reuse an existing lower-dimensional mesh prior to a restart."}>>> = false
  []
[]

[ThermalContactMortar<<<{"href": "../syntax/ThermalContactMortar/index.html"}>>>]
  [thermal_contact]
    secondary_variable<<<{"description": "The secondary variable"}>>> = temp
    primary_boundary<<<{"description": "The primary surface"}>>> = 7
    secondary_boundary<<<{"description": "The secondary surface"}>>> = 8
    initial_moles<<<{"description": "The Postprocessor that will give the initial moles of gas."}>>> = initial_moles # coupling to a postprocessor which supplies the initial plenum/gap gas mass
    gas_released<<<{"description": "The postprocessor(s) that will give the gas released."}>>> = fission_gas_released # coupling to a postprocessor which supplies the fission gas addition
    primary_subdomain<<<{"description": "The primary subdomain"}>>> = 'pellet_clad_mechanical_real_primary_subdomain'
    secondary_subdomain<<<{"description": "The secondary subdomain"}>>> = 'pellet_clad_mechanical_real_secondary_subdomain'
  []
[]

[Burnup<<<{"href": "../syntax/Burnup/index.html"}>>>]
  [burnup]
    block<<<{"description": "The blocks where radial power factor should be computed."}>>> = pellet_type_1
    rod_ave_lin_pow<<<{"description": "Rod average linear power function."}>>> = power_history # using the power function defined above
    axial_power_profile<<<{"description": "Axial power peaking function."}>>> = axial_peaking_factors # using the axial power profile function defined above
    num_radial<<<{"description": "Number of radial divisions in secondary grid used to compute radial power profile."}>>> = 80
    num_axial<<<{"description": "Number of axial divisions in secondary grid used to compute radial power profile."}>>> = 21
    axial_axis<<<{"description": "Coordinate axis of the axial direction of the fuel stack (0, 1, or 2 for x, y, or z)"}>>> = 2
    density<<<{"description": "The initial fuel density."}>>> = ${initial_fuel_density}
    a_lower<<<{"description": "The lower axial coordinate of the fuel stack. Required if fuel_pin_geometry is not specified."}>>> = -1e-3 # mesh dependent!
    a_upper<<<{"description": "The upper axial coordinate of the fuel stack. Required if fuel_pin_geometry is not specified."}>>> = 1e-3 # mesh dependent!
    fuel_inner_radius<<<{"description": "The inner radius of the fuel."}>>> = 0
    fuel_outer_radius<<<{"description": "The outer radius of the fuel."}>>> = .0041
    fuel_volume_ratio<<<{"description": "Reduction factor for deviation from right circular cylinder fuel.  The ratio of actual volume to right circular cylinder volume."}>>> = 0.987775 # for use with dished pellets (ratio of actual volume to cylinder volume)

    #N235 = N235 # Activate to write N235 concentration to output file
    #N238 = N238 # Activate to write N238 concentration to output file
    #N239 = N239 # Activate to write N239 concentration to output file
    #N240 = N240 # Activate to write N240 concentration to output file
    #N241 = N241 # Activate to write N241 concentration to output file
    #N242 = N242 # Activate to write N242 concentration to output file
    RPF<<<{"description": "Specifies that the radial power factor is required."}>>> = RPF
  []
[]

[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
  # Define auxilliary kernels for each of the aux variables
  [worn_depth]
    type = MortarArchardsLawAux<<<{"description": "Returns the weighted gap velocity at a node. This quantity is useful for mortar contact, particularly when dual basis functions are used in contact mechanics", "href": "../source/auxkernels/MortarArchardsLawAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = worn_depth
    primary_boundary<<<{"description": "The name of the primary boundary sideset."}>>> = 101
    secondary_boundary<<<{"description": "The name of the secondary boundary sideset."}>>> = 102
    primary_subdomain<<<{"description": "The name of the primary subdomain."}>>> = 'spacer_clad_mechanical_primary_subdomain'
    secondary_subdomain<<<{"description": "The name of the secondary subdomain."}>>> = 'spacer_clad_mechanical_secondary_subdomain'
    displacements<<<{"description": "The displacement variables. This mortar nodal auxiliary kernel can take two or three displacements"}>>> = 'disp_x disp_y'
    friction_coefficient<<<{"description": "Friction coefficient used to compute wear (to match that of frictional contact)"}>>> = 0.5
    energy_wear_coefficient<<<{"description": "Energy wear coefficient is a surface-dependent parameter used in Archard's wear law"}>>> = 0.1e-9
    normal_pressure<<<{"description": "The name of the Lagrange multiplier that holds the normal contact pressure in mortar formulations"}>>> = spacer_clad_mechanical_normal_lm
    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'
  []
  [fast_neutron_flux]
    type = FastNeutronFluxAux<<<{"description": "Computes fast neutron flux.", "href": "../source/auxkernels/FastNeutronFluxAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = fast_neutron_flux
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    rod_ave_lin_pow<<<{"description": "The power history function."}>>> = power_history
    axial_power_profile<<<{"description": "The axial power profile function."}>>> = axial_peaking_factors
    factor<<<{"description": "The fast neutron flux if the function, RALP, and q_variable are not given. A scaling factor if the function, RALP, or q_variable is given.  If RALP or q_variable is given, it is recommended to use a value of 3e13 (n/(m^2s)/(W/m))."}>>> = 3e13
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_begin
  []
  [fast_neutron_fluence]
    type = FastNeutronFluenceAux<<<{"description": "Computes fast neutron fluence from fast_neutron_flux.", "href": "../source/auxkernels/FastNeutronFluenceAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = fast_neutron_fluence
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    fast_neutron_flux<<<{"description": "Coupled Fast Neutron Flux."}>>> = fast_neutron_flux
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_begin
  []
  [relocation_strain]
    type = MaterialRealAux<<<{"description": "Outputs element volume-averaged material properties", "href": "../source/auxkernels/MaterialRealAux.html"}>>>
    property<<<{"description": "The material property name."}>>> = relocation_strain
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = relocation_strain
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
[]

[BCs<<<{"href": "../syntax/NuclearMaterials/BCs/index.html"}>>>]
  # Define boundary conditions
  [no_y_all] # pin pellets and clad along axis of symmetry (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"}>>> = 15
    value<<<{"description": "Value of the BC"}>>> = 0.0
  []
  [no_x_all] # pin pellets and clad along axis of symmetry (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"}>>> = 16
    value<<<{"description": "Value of the BC"}>>> = 0.0
  []

  # Flow-induced vibrations refined_excitation
  [vibration_x] # pin pellets and clad along axis of symmetry (y)
    type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.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"}>>> = '112'
    expression = 'if(t < ${end_dynamic_excitation}, 10.0*1.0e-6*sin(2*3.1415926535*20* (t - ${user_start_time})) + 2.0*1.0e-6*sin(2*3.1415926535*35*(t - ${user_start_time})), 0)'
    #expression = '0'
  []
  [vibration_y] # pin pellets and clad along axis of symmetry (y)
    type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.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"}>>> = '112'
    expression = 'if(t < ${end_dynamic_excitation}, 10.0*1.0e-6*sin(2*3.1415926535*20*(t-${user_start_time})) + 2.0*1.0e-6*sin(2*3.1415926535*35*(t-${user_start_time})) + 0.9e-4, 0.9e-4)'
    #expression = '5.9e-4'
  []

  [Pressure<<<{"href": "../syntax/BCs/Pressure/index.html"}>>>] #  apply coolant pressure on clad outer walls
    [coolantPressure]
      boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = '2'
      factor<<<{"description": "The magnitude to use in computing the pressure"}>>> = 15.5e6
      function<<<{"description": "The function that describes the pressure"}>>> = pressure_var # use the pressure_ramp function defined above
    []
  []
  [PlenumPressure<<<{"href": "../syntax/BCs/PlenumPressure/index.html"}>>>] #  apply plenum pressure on clad inner walls and pellet surfaces
    [plenumPressure]
      boundary<<<{"description": "The list of boundary IDs from the mesh where the pressure will be applied"}>>> = 9
      initial_pressure<<<{"description": "The initial pressure in the cavity.  If not given, a zero initial pressure will be used."}>>> = 2.0e6
      R<<<{"description": "The universal gas constant for the units used."}>>> = 8.3143
      output_initial_moles<<<{"description": "The name to use when reporting the initial moles of gas"}>>> = initial_moles # coupling to post processor to get initial fill gas mass
      temperature<<<{"description": "The name of the average temperature postprocessor value."}>>> = plenum_temperature # coupling to post processor to get gas temperature approximation
      volume<<<{"description": "The name of the postprocessor(s) that holds the value of the internal volume in the cavity"}>>> = plenum_volume # coupling to post processor to get gas volume
      material_input<<<{"description": "The name of the postprocessor(s) that holds the amount of material injected into the plenum."}>>> = fission_gas_released # coupling to post processor to get fission gas added
      output<<<{"description": "The name to use for the cavity pressure value"}>>> = plenum_pressure # coupling to post processor to output plenum/gap pressure
      displacements<<<{"description": "The nonlinear displacement variables"}>>> = 'disp_x disp_y'
    []
  []
  [convective_clad_surface] # apply convective boundary to clad outer surface
    type = ConvectiveFluxBC<<<{"description": "Determines boundary values via the initial and final values, flux, and exposure duration", "href": "../source/bcs/ConvectiveFluxBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '2'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temp
    rate = 38200.0 #convection coefficient (h)
    initial = 580.0
    final = 580.0
    duration = 1.0e4 #duration of initial power ramp
  []
[]

[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
  # Define material behavior models and input material property data
  [fuel_thermal] # temperature and burnup dependent thermal properties of UO2 (BISON kernel)
    type = UO2Thermal<<<{"description": "Computes thermal conductivity and specific heat capacity of uranium dioxide fuel.", "href": "../source/materials/UO2Thermal.html"}>>>
    thermal_conductivity_model<<<{"description": "The thermal conductivity model."}>>> = FINK_LUCUTA
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    burnup<<<{"description": "Coupled burnup"}>>> = burnup
    initial_porosity<<<{"description": "Initial porosity.  Must be between 0.0 and 1.0."}>>> = 0.0
  []
  [fuel_solid_mechanics_swelling] # free expansion strains (swelling and densification) for UO2 (BISON kernel)
    type = UO2VolumetricSwellingEigenstrain<<<{"description": "Computes and sums the change in fuel pellet volume due to densification and fission product release. This class applies a volumetric strain correction before adding the strain from this class to the diagonal entries of the eigenstrain tensor.", "href": "../source/materials/solid_mechanics/UO2VolumetricSwellingEigenstrain.html"}>>>
    gas_swelling_model_type<<<{"description": "Which type of model to use to calculate the gaseous swelling. Choices are SIFGRS MATPRO SIFGRS_IG. If you select SIFGRS or SIFGRS_IG, the SIFGRS model must be included in the input file."}>>> = MATPRO
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    burnup<<<{"description": "Coupled Burnup"}>>> = burnup
    initial_fuel_density<<<{"description": "Initial fuel density in kg-UO2/m^3"}>>> = 10431.0
    temperature<<<{"description": "Coupled Temperature in Kelvin"}>>> = temp
    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."}>>> = 'fuel_volumetric_eigenstrain'
  []
  [fuel_creep]
    type = UO2CreepUpdate<<<{"description": "Computes the secondary thermal and irradiation creep for UO2 LWR fuel. This material must be run in conjunction with ComputeMultipleInelasticStress.", "href": "../source/materials/solid_mechanics/UO2CreepUpdate.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    fission_rate<<<{"description": "Coupled fission rate"}>>> = fission_rate
    density<<<{"description": "Initial fuel density"}>>> = 10431.0

    initial_grain_radius<<<{"description": "Fuel grain radius (m)"}>>> = 10.0e-6
    oxygen_to_metal_ratio<<<{"description": "Oxygen to metal ratio"}>>> = 2.0
  []
  [fuel_elasticity_tensor]
    type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'pellet_type_1'
    youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 906e6
    poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.345
  []
  [fuel_stress]
    type = ComputeMultipleInelasticStress<<<{"description": "Compute state (stress and internal parameters such as plastic strains and internal parameters) using an iterative process.  Combinations of creep models and plastic models may be used.", "href": "../source/materials/ComputeMultipleInelasticStress.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    inelastic_models<<<{"description": "The material objects to use to calculate stress and inelastic strains. Note: specify creep models first and plasticity models second."}>>> = 'fuel_creep'
  []
  [fuel_thermal_expansion]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 10.0e-6
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 580.0
    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."}>>> = 'fuel_thermal_eigenstrain'
  []
  [fuel_relocation]
    type = UO2RelocationEigenstrain<<<{"description": "Accounts for cracking and relocation of fuel pellet fragments in the radial direction and is necessary for accurate modeling of LWR fuel. The linear heat rate must either be a functionor a variable.", "href": "../source/materials/solid_mechanics/UO2RelocationEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    burnup<<<{"description": "Coupled Burnup variable in FIMA"}>>> = burnup
    diameter<<<{"description": "As fabricated cold diameter of pellet in meters"}>>> = 0.0082
    rod_ave_lin_pow<<<{"description": "linear power function."}>>> = power_history
    axial_power_profile<<<{"description": "axial power peaking function."}>>> = axial_peaking_factors
    diametral_gap<<<{"description": "As fabricated cold diametral pellet-cladding gap in meters."}>>> =160e-6
    burnup_relocation_stop<<<{"description": "Burnup at which relocation strain stops (in FIMA)"}>>> = 1.e20
    relocation_activation1<<<{"description": "First activation linear power in W/m term in relocation model"}>>> = 5000
    axial_axis = 2
    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."}>>> = 'fuel_relocation_eigenstrain'
  []
  [clad_thermal]
    type = HeatConductionMaterial<<<{"description": "General-purpose material model for heat conduction", "href": "../source/materials/HeatConductionMaterial.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    thermal_conductivity<<<{"description": "The thermal conductivity value"}>>> = 16.0
    specific_heat<<<{"description": "The specific heat value"}>>> = 330.0
  []
  [clad_elasticity_tensor]
    type = ZryElasticityTensor<<<{"description": "Either provides constant elasticity constants for Zircaloy cladding or calculates the Young's modulus and Poisson's ratio for Zircaloy cladding using MATPRO relations as a function of temperature and fast neutron fluence.", "href": "../source/materials/solid_mechanics/ZryElasticityTensor.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
  []
  [clad_creep_model]
    type = ZryCreepHayesHoppeUpdate<<<{"description": "Computes the secondary thermal Hayes and Kassner secondary creep and the Hoppe irradiation creep for Zircaloy cladding.  This material must be run in conjunction with ComputeMultipleInelasticStress.", "href": "../source/materials/solid_mechanics/ZryCreepHayesHoppeUpdate.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    fast_neutron_flux<<<{"description": "The fast neutron flux"}>>> = fast_neutron_flux
    temperature<<<{"description": "The coupled temperature (K)"}>>> = temp
    zircaloy_material_type<<<{"description": "Type of zircaloy material properties to use in calculating creep. Note: ESCORE_IRRADIATIONGROWTHZR4 is not valid."}>>> = stress_relief_annealed
    model_irradiation_creep<<<{"description": "Set true to activate irradiation induced creep"}>>> = true
    model_thermal_creep<<<{"description": "Set true to activate steady state thermal creep"}>>> = true

  []
  [clad_stress]
    type = ComputeMultipleInelasticStress<<<{"description": "Compute state (stress and internal parameters such as plastic strains and internal parameters) using an iterative process.  Combinations of creep models and plastic models may be used.", "href": "../source/materials/ComputeMultipleInelasticStress.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    tangent_operator<<<{"description": "Type of tangent operator to return.  'elastic': return the elasticity tensor.  'nonlinear': return the full, general consistent tangent operator."}>>> = elastic
    inelastic_models<<<{"description": "The material objects to use to calculate stress and inelastic strains. Note: specify creep models first and plasticity models second."}>>> = 'clad_creep_model'
  []
  [clad_thermal_expansion]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 5.0e-6
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 580.0
    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."}>>> = 'clad_thermal_eigenstrain'
  []
  [clad_irrgrowth]
    type = ZryIrradiationGrowthEigenstrain<<<{"description": "Computes eigenstrain from irradiation growth in Zircaloy cladding using either the Franklin or ESCORE models.", "href": "../source/materials/solid_mechanics/ZryIrradiationGrowthEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = clad
    fast_neutron_fluence<<<{"description": "The fast neutron fluence"}>>> = fast_neutron_fluence
    axial_direction<<<{"description": "The direction you want to apply irradiation growth"}>>> = 2
    zircaloy_material_type<<<{"description": "Type of irradiation growth volumetric swelling formulation."}>>> = ESCORE_IrradiationGrowthZr4
    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."}>>> = 'clad_irradiation_eigenstrain'
  []
  [grid_thermal]
    type = HeatConductionMaterial<<<{"description": "General-purpose material model for heat conduction", "href": "../source/materials/HeatConductionMaterial.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    thermal_conductivity<<<{"description": "The thermal conductivity value"}>>> = 16.0
    specific_heat<<<{"description": "The specific heat value"}>>> = 330.0
  []
  [grid_elasticity_tensor]
    type = ZryElasticityTensor<<<{"description": "Either provides constant elasticity constants for Zircaloy cladding or calculates the Young's modulus and Poisson's ratio for Zircaloy cladding using MATPRO relations as a function of temperature and fast neutron fluence.", "href": "../source/materials/solid_mechanics/ZryElasticityTensor.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
  []
  [grid_creep_model]
    type = ZryCreepHayesHoppeUpdate<<<{"description": "Computes the secondary thermal Hayes and Kassner secondary creep and the Hoppe irradiation creep for Zircaloy cladding.  This material must be run in conjunction with ComputeMultipleInelasticStress.", "href": "../source/materials/solid_mechanics/ZryCreepHayesHoppeUpdate.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    fast_neutron_flux<<<{"description": "The fast neutron flux"}>>> = fast_neutron_flux
    temperature<<<{"description": "The coupled temperature (K)"}>>> = temp
    zircaloy_material_type<<<{"description": "Type of zircaloy material properties to use in calculating creep. Note: ESCORE_IRRADIATIONGROWTHZR4 is not valid."}>>> = stress_relief_annealed
    model_irradiation_creep<<<{"description": "Set true to activate irradiation induced creep"}>>> = true
    model_thermal_creep<<<{"description": "Set true to activate steady state thermal creep"}>>> = true
  []
  [grid_stress]
    type = ComputeMultipleInelasticStress<<<{"description": "Compute state (stress and internal parameters such as plastic strains and internal parameters) using an iterative process.  Combinations of creep models and plastic models may be used.", "href": "../source/materials/ComputeMultipleInelasticStress.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    tangent_operator<<<{"description": "Type of tangent operator to return.  'elastic': return the elasticity tensor.  'nonlinear': return the full, general consistent tangent operator."}>>> = elastic
    inelastic_models<<<{"description": "The material objects to use to calculate stress and inelastic strains. Note: specify creep models first and plasticity models second."}>>> = 'grid_creep_model'
  []
  [grid_thermal_expansion]
    type = ComputeThermalExpansionEigenstrain<<<{"description": "Computes eigenstrain due to thermal expansion with a constant coefficient", "href": "../source/materials/ComputeThermalExpansionEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'grid'
    thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 5.0e-6
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 580.0
    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."}>>> = 'grid_thermal_eigenstrain'
  []
  [grid_irrgrowth]
    type = ZryIrradiationGrowthEigenstrain<<<{"description": "Computes eigenstrain from irradiation growth in Zircaloy cladding using either the Franklin or ESCORE models.", "href": "../source/materials/solid_mechanics/ZryIrradiationGrowthEigenstrain.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = grid
    fast_neutron_fluence<<<{"description": "The fast neutron fluence"}>>> = fast_neutron_fluence
    axial_direction<<<{"description": "The direction you want to apply irradiation growth"}>>> = 2
    zircaloy_material_type<<<{"description": "Type of irradiation growth volumetric swelling formulation."}>>> = ESCORE_IrradiationGrowthZr4
    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."}>>> = 'grid_irradiation_eigenstrain'
  []

  [fission_gas_release] # Forsberg-Massih fission gas release mode
    type = UO2Sifgrs<<<{"description": "Recommended fission gas model to account for generation of fission gasses in nuclear fuel", "href": "../source/materials/UO2Sifgrs.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    temperature<<<{"description": "Coupled temperature"}>>> = temp
    fission_rate<<<{"description": "Coupled fission rate variable (fiss/m^3/s)"}>>> = fission_rate # coupling to fission_rate aux variable
    grain_radius<<<{"description": "Coupled grain Radius"}>>> = 10.0e-6
    #external_pressure = 40e6
  []
  [clad_density]
    type = StrainAdjustedDensity<<<{"description": "Creates density material property", "href": "../source/materials/StrainAdjustedDensity.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'clad'
    density = 6551.0
  []
  [fuel_density]
    type = StrainAdjustedDensity<<<{"description": "Creates density material property", "href": "../source/materials/StrainAdjustedDensity.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet
    strain_free_density<<<{"description": "Material property for strain-free density"}>>> = 10431.0
  []
  [grid]
    type = StrainAdjustedDensity<<<{"description": "Creates density material property", "href": "../source/materials/StrainAdjustedDensity.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = grid
    density = 6560
  []
[]

[Debug<<<{"href": "../syntax/Debug/index.html"}>>>]
  show_var_residual_norms<<<{"description": "Print the residual norms of the individual solution variables at each nonlinear iteration"}>>> = true
[]

[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
  type = Transient
  solve_type = 'NEWTON'
  petsc_options = '-snes_converged_reason -ksp_converged_reason'
  petsc_options_iname = '-pc_type -pc_factor_mat_solver_package -mat_mffd_err -pc_factor_shift_type -pc_factor_shift_amount'
  petsc_options_value = 'lu    superlu_dist   1e-6          NONZERO               1e-14'
  snesmf_reuse_base = true

  line_search = 'basic'

  l_max_its = 100
  l_tol = 8e-3

  nl_max_its = 35
  nl_rel_tol = 1e-7
  nl_abs_tol = 1e-11

  [TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
    type = NewmarkBeta
    beta = 0.25
    gamma = 0.5
  []

  start_time = '${user_start_time}'
  end_time = '${user_end_time}'
  timestep_tolerance = 1e-8

  [TimeStepper<<<{"href": "../syntax/Executioner/TimeStepper/index.html"}>>>]
    type = IterationAdaptiveDT
    dt = '${time_step_dynamics}'
    time_t = '${end_dynamic_excitation}'
    time_dt = '${time_step_dynamics}'
    growth_factor = 1.2
    cutback_factor = 0.75
  []

  dtmax = 3e-2
[]

[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
  # Define postprocessors (some are required as specified above; others are optional; many others are available)
  [average_interior_clad_temperature] # average temperature of cladding interior
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../source/postprocessors/SideAverageValue.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 7
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    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 timestep_end'
  []
  [average_centerline_fuel_temperature] # average temperature of the cladding interior and all pellet exteriors
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../source/postprocessors/SideAverageValue.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 9
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    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 linear'
  []
  [plenum_temperature]
    type = SideAverageValue<<<{"description": "Computes the average value of a variable on a sideset. Note that this cannot be used on the centerline of an axisymmetric model.", "href": "../source/postprocessors/SideAverageValue.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 9
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    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 timestep_end'
  []
  [plenum_volume] # gas volume
    type = InternalVolume<<<{"description": "Computes the volume of an enclosed area by performing an integral over a user-supplied boundary.", "href": "../source/postprocessors/InternalVolume.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 9
    addition<<<{"description": "An additional volume to be included in the internal volume calculation. A time-dependent function is expected."}>>> = 1.3e-5 #rough guess of plenum volume/unit length of fuel
    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 linear'
  []
  [pellet_volume] # fuel pellet total volume
    type = InternalVolume<<<{"description": "Computes the volume of an enclosed area by performing an integral over a user-supplied boundary.", "href": "../source/postprocessors/InternalVolume.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 8
    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 timestep_end'
  []
  [clad_inner_vol] # volume inside of cladding
    type = InternalVolume<<<{"description": "Computes the volume of an enclosed area by performing an integral over a user-supplied boundary.", "href": "../source/postprocessors/InternalVolume.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 7
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = exodus
    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 timestep_end'
  []
  [fission_gas_generated] # fission gas produced (moles)
    type = ElementIntegralFisGasGeneratedSifgrs<<<{"description": "Reports the fission gas that is produced in moles.  To be used in combination with the Sifgrs model.", "href": "../source/postprocessors/ElementIntegralFisGasGeneratedSifgrs.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
  [fission_gas_released] # fission gas released to plenum (moles)
    type = ElementIntegralFisGasReleasedSifgrs<<<{"description": "Reports the fission gas that is released to the plenum in moles.  To be used in combination with the Sifgrs model.", "href": "../source/postprocessors/ElementIntegralFisGasReleasedSifgrs.html"}>>>
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
  [flux_from_clad] # area integrated heat flux from the cladding
    type = SideDiffusiveFluxIntegral<<<{"description": "Computes the integral of the diffusive flux over the specified boundary", "href": "../source/postprocessors/SideDiffusiveFluxIntegral.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 5
    diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation. This must be provided if the variable is of finite element type"}>>> = thermal_conductivity
    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
  []
  [flux_from_fuel] # area integrated heat flux from the fuel
    type = SideDiffusiveFluxIntegral<<<{"description": "Computes the integral of the diffusive flux over the specified boundary", "href": "../source/postprocessors/SideDiffusiveFluxIntegral.html"}>>>
    variable<<<{"description": "The name of the variable which this postprocessor integrates"}>>> = temp
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 10
    diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation. This must be provided if the variable is of finite element type"}>>> = thermal_conductivity
    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
  []
  [_dt] # time step
    type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../source/postprocessors/TimestepSize.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
  []
  [num_lin_it]
    type = NumLinearIterations<<<{"description": "Compute the number of linear iterations.", "href": "../source/postprocessors/NumLinearIterations.html"}>>>
  []
  [num_nonlin_it]
    type = NumNonlinearIterations<<<{"description": "Outputs the number of nonlinear iterations", "href": "../source/postprocessors/NumNonlinearIterations.html"}>>>
  []
  [tot_lin_it]
    type = CumulativeValuePostprocessor<<<{"description": "Creates a cumulative sum of a Postprocessor value with time.", "href": "../source/postprocessors/CumulativeValuePostprocessor.html"}>>>
    postprocessor<<<{"description": "The name of the postprocessor"}>>> = num_lin_it
  []
  [tot_nonlin_it]
    type = CumulativeValuePostprocessor<<<{"description": "Creates a cumulative sum of a Postprocessor value with time.", "href": "../source/postprocessors/CumulativeValuePostprocessor.html"}>>>
    postprocessor<<<{"description": "The name of the postprocessor"}>>> = num_nonlin_it
  []
  [alive_time]
    type = PerfGraphData<<<{"description": "Retrieves performance information about a section from the PerfGraph.", "href": "../source/postprocessors/PerfGraphData.html"}>>>
    section_name<<<{"description": "The name of the section to get data for"}>>> = Root
    data_type<<<{"description": "The type of data to retrieve for the section_name"}>>> = TOTAL
  []
  [rod_total_power]
    type = ElementIntegralPower<<<{"description": "Computes the power given the fission rate and energy per fission.", "href": "../source/postprocessors/ElementIntegralPower.html"}>>>
    variable<<<{"description": "The name of the variable that this object operates on"}>>> = temp
    fission_rate<<<{"description": "Coupled fission rate"}>>> = fission_rate
    block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = pellet_type_1
    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
  []
  [rod_input_power]
    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."}>>> = power_history
    scale_factor<<<{"description": "A scale factor to be applied to the function"}>>> = 0.1186 # rod height
    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
  []
  [fission_gas_released_percentage]
    type = FGRPercent<<<{"description": "Computes the ratio of fission gas released to the fission gas generated in percent.", "href": "../source/postprocessors/FGRPercent.html"}>>>
    fission_gas_released<<<{"description": "The fission gas released postprocessor"}>>> = fission_gas_released
    fission_gas_generated<<<{"description": "The fission gas generated postprocessor"}>>> = fission_gas_generated
  []
[]

[VectorPostprocessors<<<{"href": "../syntax/VectorPostprocessors/index.html"}>>>]
  [contact_pressure]
    type = NodalValueSampler<<<{"description": "Samples values of nodal variable(s).", "href": "../source/vectorpostprocessors/NodalValueSampler.html"}>>>
    sort_by<<<{"description": "What to sort the samples by"}>>> = x
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = spacer_clad_mechanical_normal_lm
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 102
  []
  [frictional_pressure]
    type = NodalValueSampler<<<{"description": "Samples values of nodal variable(s).", "href": "../source/vectorpostprocessors/NodalValueSampler.html"}>>>
    sort_by<<<{"description": "What to sort the samples by"}>>> = x
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = spacer_clad_mechanical_tangential_lm
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 102
  []
  [worn_depth]
    type = NodalValueSampler<<<{"description": "Samples values of nodal variable(s).", "href": "../source/vectorpostprocessors/NodalValueSampler.html"}>>>
    sort_by<<<{"description": "What to sort the samples by"}>>> = x
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    variable<<<{"description": "The names of the variables that this VectorPostprocessor operates on"}>>> = worn_depth
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 102
    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
  []
[]

[Outputs<<<{"href": "../syntax/Outputs/index.html"}>>>]
  perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = true
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'FINAL'
  [console]
    type = Console<<<{"description": "Object for screen output.", "href": "../source/outputs/Console.html"}>>>
    max_rows<<<{"description": "The maximum number of postprocessor/scalar values displayed on screen during a timestep (set to 0 for unlimited)"}>>> = 25
  []
  checkpoint<<<{"description": "Create checkpoint files using the default options."}>>> = true
  file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = 'step_${step_number}'
[]
(examples/2D_plane_strain_fretting_wear/fretting-wear-initial-dyn-exc.i)

References

  1. A. Recuero and A. Lindsay. On practical aspects of variational consistency in contact dynamics. Submitted for publication, -:1 – 10, 2023.[BibTeX]
  2. Antonio Recuero, Alexander Lindsay, and Dewen Yushu. An approach to grid-to-rod fretting wear modeling using dynamic mortar contact. Progress in Nuclear Energy, 163:104793, 2023. doi:10.1016/j.pnucene.2023.104793.[BibTeX]