Thermo-mechanical coupling of an one pin problem

This is an example case of coupling a subchannel calculation (SCM) with a 2D-RZ thermo_mechanical model of a fuel pin.

Example Description

The example models 4 subchannels around one fuel pin (quadrilateral geometry). SCM calculates coolant temperature distribution and an average pin surface temperature Tpin that gets transfered to the pin model. The pin model calculates heat conduction within the pin and thermal expansion. Then the pin model returns linear heat flux at the surface of the pin q_prime and the updated pin diameter Dpin to SCM. In order for the updated pin diameter to be considered by the SCM model the user must have activated the flag: "deformation" = true in the problem/SubChannel block of the SCM input file.

Input files

The input files are the following:

SCM Model.

######## BC's #################
T_in = 359.15
# [1e+6 kg/m^2-hour] turns into kg/m^2-sec
mass_flux_in = '${fparse 1e+6 * 17.00 / 3600.}'
P_out = 4.923e6 # Pa
heated_length = 1.0
######## Geometry #
[GlobalParams<<<{"href": "../../../syntax/GlobalParams/index.html"}>>>]
    nx = 2
    ny = 2
    n_cells = 25
    pitch = 0.0126
    pin_diameter = 0.00950
    gap = 0.00095
    heated_length = ${heated_length}
    spacer_z = '0.0'
    spacer_k = '0.0'
    power = 100000.0 # W
[]

[QuadSubChannelMesh]
    [subchannel]
        type = SCMQuadSubChannelMeshGenerator
    []

    [fuel_pins]
        type = SCMQuadPinMeshGenerator
        input = subchannel
    []
[]

[Functions<<<{"href": "../../../syntax/Functions/index.html"}>>>]
    [axial_heat_rate]
        type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../source/functions/MooseParsedFunction.html"}>>>
        expression<<<{"description": "The user defined function."}>>> = '(pi/2)*sin(pi*z/L)'
        symbol_names<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = 'L'
        symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = '${heated_length}'
    []
[]

[AuxVariables<<<{"href": "../../../syntax/AuxVariables/index.html"}>>>]
    [mdot]
        block = subchannel
    []
    [SumWij]
        block = subchannel
    []
    [P]
        block = subchannel
    []
    [DP]
        block = subchannel
    []
    [h]
        block = subchannel
    []
    [T]
        block = subchannel
    []
    [rho]
        block = subchannel
    []
    [mu]
        block = subchannel
    []
    [S]
        block = subchannel
    []
    [w_perim]
        block = subchannel
    []
    [q_prime]
        block = fuel_pins
    []
    [Tpin]
        block = fuel_pins
    []
    [Dpin]
        block = fuel_pins
    []
[]

[FluidProperties<<<{"href": "../../../syntax/FluidProperties/index.html"}>>>]
    [water]
        type = Water97FluidProperties<<<{"description": "Fluid properties for water and steam (H2O) using IAPWS-IF97", "href": "../../../source/fluidproperties/Water97FluidProperties.html"}>>>
    []
[]

[SubChannel]
    type = QuadSubChannel1PhaseProblem
    fp = water
    n_blocks = 1
    beta = 0.006
    CT = 2.6
    compute_density = true
    compute_viscosity = true
    compute_power = true
    P_out = ${P_out}
    verbose_subchannel = true
    deformation = true
[]

[ICs<<<{"href": "../../../syntax/ICs/index.html"}>>>]
    [S_IC]
        type = SCMQuadFlowAreaIC<<<{"description": "Computes subchannel flow area in the square lattice subchannel arrangement", "href": "../../../source/ics/SCMQuadFlowAreaIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = S
    []

    [w_perim_IC]
        type = SCMQuadWettedPerimIC<<<{"description": "Computes wetted perimeter of subchannels in a square lattice arrangement", "href": "../../../source/ics/SCMQuadWettedPerimIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = w_perim
    []

    [q_prime_IC]
        type = SCMQuadPowerIC<<<{"description": "Computes axial heat rate (W/m) that goes into the subchannel cells or is assigned to the fuel pins, in a square lattice arrangement", "href": "../../../source/ics/SCMQuadPowerIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = q_prime
        filename<<<{"description": "name of radial power profile .txt file (should be a single column) [UnitLess]."}>>> = "power_profile.txt"
        axial_heat_rate<<<{"description": "user provided normalized function of axial heat rate [Unitless]. The integral over pin heated length should equal the heated length"}>>> = axial_heat_rate
    []

    [T_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = T
        value<<<{"description": "The value to be set in IC"}>>> = ${T_in}
    []

    [Dpin_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = Dpin
        value<<<{"description": "The value to be set in IC"}>>> = 0.00950
    []

    [P_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = P
        value<<<{"description": "The value to be set in IC"}>>> = 0.0
    []

    [DP_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = DP
        value<<<{"description": "The value to be set in IC"}>>> = 0.0
    []

    [Viscosity_ic]
        type = ViscosityIC<<<{"description": "Computes viscosity from specified pressure and temperature", "href": "../../../source/ics/ViscosityIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = mu
        p<<<{"description": "Pressure [Pa]"}>>> = ${P_out}
        T<<<{"description": "Temperature [K]"}>>> = T
        fp<<<{"description": "Fluid properties user object name"}>>> = water
    []

    [rho_ic]
        type = RhoFromPressureTemperatureIC<<<{"description": "Computes the density from pressure and temperature.", "href": "../../../source/ics/RhoFromPressureTemperatureIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = rho
        p<<<{"description": "The pressure [Pa]"}>>> = ${P_out}
        T<<<{"description": "The temperature [K]"}>>> = T
        fp<<<{"description": "The name of fluid properties user object."}>>> = water
    []

    [h_ic]
        type = SpecificEnthalpyFromPressureTemperatureIC<<<{"description": "Computes the specific enthalpy from pressure and temperature.", "href": "../../../source/ics/SpecificEnthalpyFromPressureTemperatureIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = h
        p<<<{"description": "The pressure [Pa]"}>>> = ${P_out}
        T<<<{"description": "The temperature [K]"}>>> = T
        fp<<<{"description": "The name of fluid properties user object."}>>> = water
    []

    [mdot_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = mdot
        value<<<{"description": "The value to be set in IC"}>>> = 0.0
    []
[]

[AuxKernels<<<{"href": "../../../syntax/AuxKernels/index.html"}>>>]
    [T_in_bc]
        type = ConstantAux<<<{"description": "Creates a constant field in the domain.", "href": "../../../source/auxkernels/ConstantAux.html"}>>>
        variable<<<{"description": "The name of the variable that this object applies to"}>>> = T
        boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = inlet
        value<<<{"description": "Some constant value that can be read from the input file"}>>> = ${T_in}
        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'
    []
    [mdot_in_bc]
        type = SCMMassFlowRateAux<<<{"description": "Computes mass flow rate from specified mass flux and subchannel cross-sectional area. Can read either PostprocessorValue or Real", "href": "../../../source/auxkernels/SCMMassFlowRateAux.html"}>>>
        variable<<<{"description": "The name of the variable that this object applies to"}>>> = mdot
        boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = inlet
        area<<<{"description": "Cross sectional area [m^2]"}>>> = S
        mass_flux<<<{"description": "The postprocessor or Real to use for the value of mass_flux"}>>> = ${mass_flux_in}
        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'
    []
[]

[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
    exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
    [Temp_Out_MATRIX]
        type = QuadSubChannelNormalSliceValues<<<{"description": "Prints out a user selected value at a user selected axial height in a matrix format to be used for post-processing", "href": "../../../source/outputs/QuadSubChannelNormalSliceValues.html"}>>>
        variable<<<{"description": "Variable you want the value of"}>>> = T
        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
        file_base<<<{"description": "The desired solution output name without an extension. If not provided, MOOSE sets it with Outputs/file_base when available. Otherwise, MOOSE uses input file name and this object name for a master input or uses master file_base, the subapp name and this object name for a subapp input to set it."}>>> = "Temp_Out.txt"
        height<<<{"description": "Axial location of normal slice [m]"}>>> = 1.0
    []
    [mdot_Out_MATRIX]
        type = QuadSubChannelNormalSliceValues<<<{"description": "Prints out a user selected value at a user selected axial height in a matrix format to be used for post-processing", "href": "../../../source/outputs/QuadSubChannelNormalSliceValues.html"}>>>
        variable<<<{"description": "Variable you want the value of"}>>> = mdot
        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
        file_base<<<{"description": "The desired solution output name without an extension. If not provided, MOOSE sets it with Outputs/file_base when available. Otherwise, MOOSE uses input file name and this object name for a master input or uses master file_base, the subapp name and this object name for a subapp input to set it."}>>> = "mdot_Out.txt"
        height<<<{"description": "Axial location of normal slice [m]"}>>> = 1.0
    []
    [mdot_In_MATRIX]
        type = QuadSubChannelNormalSliceValues<<<{"description": "Prints out a user selected value at a user selected axial height in a matrix format to be used for post-processing", "href": "../../../source/outputs/QuadSubChannelNormalSliceValues.html"}>>>
        variable<<<{"description": "Variable you want the value of"}>>> = mdot
        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
        file_base<<<{"description": "The desired solution output name without an extension. If not provided, MOOSE sets it with Outputs/file_base when available. Otherwise, MOOSE uses input file name and this object name for a master input or uses master file_base, the subapp name and this object name for a subapp input to set it."}>>> = "mdot_In.txt"
        height<<<{"description": "Axial location of normal slice [m]"}>>> = 0.0
    []
[]

[Executioner<<<{"href": "../../../syntax/Executioner/index.html"}>>>]
    type = Steady
    petsc_options_iname = '-pc_type -pc_hypre_type'
    petsc_options_value = 'hypre boomeramg'
    fixed_point_max_its = 8
    fixed_point_min_its = 6
    fixed_point_rel_tol = 1e-6
[]

################################################################################
[MultiApps<<<{"href": "../../../syntax/MultiApps/index.html"}>>>]
    ################################################################################
    # Couple to Thermo-Mechanical Pin model (uses kernels available in MOOSE)
    ################################################################################
    [sub]
        type = FullSolveMultiApp<<<{"description": "Performs a complete simulation during each execution.", "href": "../../../source/multiapps/FullSolveMultiApp.html"}>>>
        input_files<<<{"description": "The input file for each App.  If this parameter only contains one input file it will be used for all of the Apps.  When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = one_pin_problem_sub.i
        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'
        positions<<<{"description": "The positions of the App locations.  Each set of 3 values will represent a Point.  This and 'positions_file' cannot be both supplied. If this and 'positions_file'/'_objects' are not supplied, a single position (0,0,0) will be used"}>>> = '0   0   0 '
        output_in_position<<<{"description": "If true this will cause the output from the MultiApp to be 'moved' by its position vector"}>>> = true
        bounding_box_padding<<<{"description": "Additional padding added to the dimensions of the bounding box. The values are added to the x, y and z dimension respectively."}>>> = '0 0 0.01'
    []

    ################################################################################
    # A multiapp that projects the solution to a detailed mesh for visualization purposes
    ################################################################################
    [viz]
        type = FullSolveMultiApp<<<{"description": "Performs a complete simulation during each execution.", "href": "../../../source/multiapps/FullSolveMultiApp.html"}>>>
        input_files<<<{"description": "The input file for each App.  If this parameter only contains one input file it will be used for all of the Apps.  When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = '3d.i'
        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'
    []
[]

[Transfers<<<{"href": "../../../syntax/Transfers/index.html"}>>>]
    [Tpin] # send pin surface temperature to pin model,
        type = MultiAppGeneralFieldNearestLocationTransfer<<<{"description": "Transfers field data at the MultiApp position by finding the value at the nearest neighbor(s) in the origin application.", "href": "../../../source/transfers/MultiAppGeneralFieldNearestLocationTransfer.html"}>>>
        to_multi_app<<<{"description": "The name of the MultiApp to transfer the data to"}>>> = sub
        variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = Pin_surface_temperature
        source_variable<<<{"description": "The variable to transfer from."}>>> = Tpin
        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'
        greedy_search<<<{"description": "Whether or not to send a point to all the domains. If true, all the processors will be checked for a given point.The code will be slow if this flag is on but it will give a better solution."}>>> = true
    []

    [diameter] # send diameter information from pin model to subchannel
        type = MultiAppGeneralFieldNearestLocationTransfer<<<{"description": "Transfers field data at the MultiApp position by finding the value at the nearest neighbor(s) in the origin application.", "href": "../../../source/transfers/MultiAppGeneralFieldNearestLocationTransfer.html"}>>>
        from_multi_app<<<{"description": "The name of the MultiApp to receive data from"}>>> = sub
        variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = Dpin
        source_variable<<<{"description": "The variable to transfer from."}>>> = pin_diameter_deformed
        from_boundaries<<<{"description": "The boundary we are transferring from (if not specified, whole domain is used)."}>>> = right
        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'
        greedy_search<<<{"description": "Whether or not to send a point to all the domains. If true, all the processors will be checked for a given point.The code will be slow if this flag is on but it will give a better solution."}>>> = true
    []

    [q_prime] # send heat flux from pin model to subchannel
        type = MultiAppGeneralFieldNearestLocationTransfer<<<{"description": "Transfers field data at the MultiApp position by finding the value at the nearest neighbor(s) in the origin application.", "href": "../../../source/transfers/MultiAppGeneralFieldNearestLocationTransfer.html"}>>>
        from_multi_app<<<{"description": "The name of the MultiApp to receive data from"}>>> = sub
        variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = q_prime
        source_variable<<<{"description": "The variable to transfer from."}>>> = q_prime
        from_boundaries<<<{"description": "The boundary we are transferring from (if not specified, whole domain is used)."}>>> = right
        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'
        greedy_search<<<{"description": "Whether or not to send a point to all the domains. If true, all the processors will be checked for a given point.The code will be slow if this flag is on but it will give a better solution."}>>> = true
    []

    [subchannel_transfer]
        type = SCMSolutionTransfer<<<{"description": "Transfers subchannel solution from computational mesh onto visualization mesh", "href": "../../../source/transfers/SCMSolutionTransfer.html"}>>>
        to_multi_app<<<{"description": "The name of the MultiApp to transfer the data to"}>>> = viz
        variable<<<{"description": "The auxiliary variables to transfer."}>>> = 'mdot SumWij P DP h T rho mu S w_perim'
    []

    [pin_transfer]
        type = SCMPinSolutionTransfer<<<{"description": "Transfers pin solution from computational mesh onto visualization mesh", "href": "../../../source/transfers/SCMPinSolutionTransfer.html"}>>>
        to_multi_app<<<{"description": "The name of the MultiApp to transfer the data to"}>>> = viz
        variable<<<{"description": "The auxiliary variables to transfer."}>>> = 'Tpin Dpin q_prime'
    []
[]
(modules/subchannel/examples/coupling/thermo_mech/quad/one_pin_problem.i)

Fuel Pin Model.

pin_diameter = 0.00950
heated_length = 1.0
T_in = 359.15
[GlobalParams<<<{"href": "../../../syntax/GlobalParams/index.html"}>>>]
    displacements = 'disp_x disp_y'
[]

[Mesh<<<{"href": "../../../syntax/Mesh/index.html"}>>>]
    second_order<<<{"description": "Converts a first order mesh to a second order mesh.  Note: This is NOT needed if you are reading an actual first order mesh."}>>> = true
    [PinMesh]
        type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
        dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
        xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = '${fparse pin_diameter / 2.0}'
        bias_x<<<{"description": "The amount by which to grow (or shrink) the cells in the x-direction."}>>> = 1.0
        nx<<<{"description": "Number of elements in the X direction"}>>> = 20
        ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = ${heated_length}
        ny<<<{"description": "Number of elements in the Y direction"}>>> = 100 # number of axial cells
    []
    coord_type = RZ
    rz_coord_axis = Y
    beta_rotation = 90
[]

[Functions<<<{"href": "../../../syntax/Functions/index.html"}>>>]
    [volumetric_heat_rate] # Defined such as to be consistent with the IC in SCM
        type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../source/functions/MooseParsedFunction.html"}>>>
        expression<<<{"description": "The user defined function."}>>> = '(4.0 * 100000.0 / (pi * D * D * L)) * (pi/2)*sin(pi*y/L)'
        symbol_names<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = 'L D'
        symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = '${heated_length} ${pin_diameter}'
    []
[]

[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
    [temperature]
        order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = SECOND
        family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
    []
[]

[AuxVariables<<<{"href": "../../../syntax/AuxVariables/index.html"}>>>]
    [Pin_surface_temperature]
    []
    [pin_diameter_deformed]
        # order = CONSTANT
        # family = MONOMIAL
    []
    [q_prime]
        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
    []
[]

[Physics<<<{"href": "../../../syntax/Physics/index.html"}>>>]
    [SolidMechanics<<<{"href": "../../../syntax/Physics/SolidMechanics/index.html"}>>>]
        [QuasiStatic<<<{"href": "../../../syntax/Physics/SolidMechanics/QuasiStatic/index.html"}>>>]
            add_variables<<<{"description": "Add the displacement variables"}>>> = true
            strain<<<{"description": "Strain formulation"}>>> = SMALL
            incremental<<<{"description": "Use incremental or total strain (if not explicitly specified this defaults to incremental for finite strain and total for small strain)"}>>> = true
            generate_output<<<{"description": "Add scalar quantity output for stress and/or strain"}>>> = 'stress_xx stress_yy stress_xy'
            temperature<<<{"description": "The temperature"}>>> = temperature

            [block0]
                eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = eigenstrain
                block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = 0
            []
        []
    []
[]

[AuxKernels<<<{"href": "../../../syntax/AuxKernels/index.html"}>>>]
    [QPrime]
        type = SCMRZPinQPrimeAux<<<{"description": "Axial heat rate on pin surface for a 2D-RZ axi-symmetric fuel pin model", "href": "../../../source/auxkernels/SCMRZPinQPrimeAux.html"}>>> # This kernel calculates linear heat rate for the 2D-RZ pin model
        diffusivity<<<{"description": "The name of the diffusivity material property that will be used in the flux computation."}>>> = 'thermal_conductivity'
        variable<<<{"description": "The name of the variable that this object applies to"}>>> = q_prime
        diffusion_variable<<<{"description": "The name of the variable"}>>> = temperature
        component<<<{"description": "The desired component of flux."}>>> = normal
        boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'right'
        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'
        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
    []
    [Deformation]
        type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../../source/auxkernels/ParsedAux.html"}>>>
        variable<<<{"description": "The name of the variable that this object applies to"}>>> = pin_diameter_deformed
        coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'disp_x'
        expression<<<{"description": "Parsed function expression to compute"}>>> = '2.0 * disp_x + ${pin_diameter}'
        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'
    []
[]

[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
    [heat_conduction]
        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"}>>> = temperature
        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
    []
    [heat_source]
        type = HeatSource<<<{"description": "Demonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form $(\\psi_i, -f)$.", "href": "../../../source/kernels/HeatSource.html"}>>>
        variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
        function<<<{"description": "Function describing the volumetric heat source"}>>> = volumetric_heat_rate
        use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
    []
[]

[Materials<<<{"href": "../../../syntax/Materials/index.html"}>>>]
    [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"}>>> = 0
        bulk_modulus<<<{"description": "The bulk modulus for the material."}>>> = 0.333333333333e6
        poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.0
    []
    [thermal_strain]
        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"}>>> = 0
        temperature<<<{"description": "Coupled temperature"}>>> = temperature
        stress_free_temperature<<<{"description": "Reference temperature at which there is no thermal expansion for thermal eigenstrain calculation"}>>> = 117.56
        thermal_expansion_coeff<<<{"description": "Thermal expansion coefficient"}>>> = 1e-5
        eigenstrain_name<<<{"description": "Material property name for the eigenstrain tensor computed by this model. IMPORTANT: The name of this property must also be provided to the strain calculator."}>>> = eigenstrain
    []
    [stress]
        type = ComputeStrainIncrementBasedStress<<<{"description": "Compute stress after subtracting inelastic strain increments", "href": "../../../source/materials/ComputeStrainIncrementBasedStress.html"}>>>
        block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 0
    []
    [heat_conductor]
        type = HeatConductionMaterial<<<{"description": "General-purpose material model for heat conduction", "href": "../../../source/materials/HeatConductionMaterial.html"}>>>
        thermal_conductivity<<<{"description": "The thermal conductivity value"}>>> = 1.0
        block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 0
    []
    [density]
        type = Density<<<{"description": "Creates density material property. This class is deprecated, and its functionalityis replaced by StrainAdjustedDensity for cases when the density should be adjustedto account for material deformation. If it is not desired to adjust the density fordeformation, a variety of general-purpose Materials, such as GenericConstantMaterialor ParsedMaterial can be used to define the density.", "href": "../../../source/materials/Density.html"}>>>
        block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 0
        density<<<{"description": "Density"}>>> = 1.0
    []
[]

[BCs<<<{"href": "../../../syntax/BCs/index.html"}>>>]
    [left]
        type = NeumannBC<<<{"description": "Imposes the integrated boundary condition $\\frac{\\partial u}{\\partial n}=h$, where $h$ is a constant, controllable value.", "href": "../../../source/bcs/NeumannBC.html"}>>>
        variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
        boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
    []
    [right]
        type = MatchedValueBC<<<{"description": "Implements a NodalBC which equates two different Variables' values on a specified boundary.", "href": "../../../source/bcs/MatchedValueBC.html"}>>>
        variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = temperature
        boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
        v<<<{"description": "The variable whose value we are to match."}>>> = Pin_surface_temperature
    []
    [no_x]
        type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../source/bcs/DirichletBC.html"}>>>
        variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
        boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
        value<<<{"description": "Value of the BC"}>>> = 0.0
    []
    [no_y]
        type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../../source/bcs/DirichletBC.html"}>>>
        variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
        boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom top'
        value<<<{"description": "Value of the BC"}>>> = 0.0
    []
[]

[ICs<<<{"href": "../../../syntax/ICs/index.html"}>>>]
    [temperature_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = temperature
        value<<<{"description": "The value to be set in IC"}>>> = ${T_in}
    []
    [q_prime_ic]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = q_prime
        value<<<{"description": "The value to be set in IC"}>>> = 0.0
    []
    [RD_IC]
        type = ConstantIC<<<{"description": "Sets a constant field value.", "href": "../../../source/ics/ConstantIC.html"}>>>
        variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = pin_diameter_deformed
        value<<<{"description": "The value to be set in IC"}>>> = ${pin_diameter}
    []
[]

[Preconditioning<<<{"href": "../../../syntax/Preconditioning/index.html"}>>>]
    [smp]
        type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
        full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
    []
[]

[Executioner<<<{"href": "../../../syntax/Executioner/index.html"}>>>]
    type = Transient
    solve_type = 'PJFNK'
    nl_abs_tol = 1e-7
    l_max_its = 20
    start_time = 0.0
    dt = 0.5
    num_steps = 2
    end_time = 1.0
    petsc_options_iname = '-pc_type -mat_fd_coloring_err -mat_fd_type'
    petsc_options_value = 'lu       1e-6                 ds'

    [Quadrature<<<{"href": "../../../syntax/Executioner/Quadrature/index.html"}>>>]
        order<<<{"description": "Order of the quadrature"}>>> = FIFTH
        side_order<<<{"description": "Order of the quadrature for sides"}>>> = SEVENTH
    []
[]

[Postprocessors<<<{"href": "../../../syntax/Postprocessors/index.html"}>>>]
    [total_power]
        type = ElementIntegralFunctorPostprocessor<<<{"description": "Computes a volume integral of the specified functor", "href": "../../../source/postprocessors/ElementIntegralFunctorPostprocessor.html"}>>>
        functor<<<{"description": "The name of the functor that this object operates on. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = volumetric_heat_rate
        use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
    []
    [total_power_disp]
        type = ElementIntegralFunctorPostprocessor<<<{"description": "Computes a volume integral of the specified functor", "href": "../../../source/postprocessors/ElementIntegralFunctorPostprocessor.html"}>>>
        functor<<<{"description": "The name of the functor that this object operates on. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = volumetric_heat_rate
        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
    []
    [volume]
        type = VolumePostprocessor<<<{"description": "Computes the volume of a specified block", "href": "../../../source/postprocessors/VolumePostprocessor.html"}>>>
    []
    [volume_disp]
        type = VolumePostprocessor<<<{"description": "Computes the volume of a specified block", "href": "../../../source/postprocessors/VolumePostprocessor.html"}>>>
        use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
    []
[]

[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
    exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(modules/subchannel/examples/coupling/thermo_mech/quad/one_pin_problem_sub.i)

Caveat

The SCM input file utilizes the MultiApp functionality and runs the fuel model as a SubApp. The executable to be used is the combined-opt that is compiled with the make file in combined folder of MOOSE.