3D laser welding

commentnote

This input requires tests objects. It can only be run with the Navier Stokes module executable or an application that also compiled the Navier Stokes test objects. The --allow-test-objects argument must be passed on the command line as well.

The input file below can be used to model a full rotation of a laser spot around the surface of a cubic representation of a welding material. This input, whose results are published in Lindsay et al. (2021), reproduces a model outlined in Noble et al. (2007). A simple Laplacian equation is used to model the displacement field:

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [disp_x]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
  []

  [disp_y]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
  []

  [disp_z]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

The incompressible Navier-Stokes equations are solved for mass, momentum, and energy. These equations are run on the displaced mesh such that a mesh convection term (see INSADMomentumMeshAdvection) must be added to correct the material velocity, as shown in equation 2 of Kong and Cai (2017). Both SUPG and PSPG stabilizations are used in this input. The kernels used to model the Navier-Stokes equations are shown below:

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [mass]
    type = INSADMass<<<{"description": "This class computes the mass equation residual and Jacobian contributions (the latter using automatic differentiation) for the incompressible Navier-Stokes equations.", "href": "../../source/kernels/INSADMass.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
    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
  []
  [mass_pspg]
    type = INSADMassPSPG<<<{"description": "This class adds PSPG stabilization to the mass equation, enabling use of equal order shape functions for pressure and velocity variables", "href": "../../source/kernels/INSADMassPSPG.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
    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
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative<<<{"description": "This class computes the time derivative for the incompressible Navier-Stokes momentum equation.", "href": "../../source/kernels/INSADMomentumTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    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
  []
  [momentum_advection]
    type = INSADMomentumAdvection<<<{"description": "Adds the advective term to the INS momentum equation", "href": "../../source/kernels/INSADMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    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
  []
  [momentum_mesh_advection]
    type = INSADMomentumMeshAdvection<<<{"description": "Corrects the convective derivative for situations in which the fluid mesh is dynamic.", "href": "../../source/kernels/INSADMomentumMeshAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    disp_x<<<{"description": "The x displacement"}>>> = disp_x
    disp_y<<<{"description": "The y displacement"}>>> = disp_y
    disp_z<<<{"description": "The z displacement"}>>> = disp_z
    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
  []
  [momentum_viscous]
    type = INSADMomentumViscous<<<{"description": "Adds the viscous term to the INS momentum equation", "href": "../../source/kernels/INSADMomentumViscous.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    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
  []
  [momentum_pressure]
    type = INSADMomentumPressure<<<{"description": "Adds the pressure term to the INS momentum equation", "href": "../../source/kernels/INSADMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    pressure<<<{"description": "The pressure"}>>> = p
    integrate_p_by_parts<<<{"description": "Whether to integrate the pressure term by parts"}>>> = true
    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
  []
  [momentum_supg]
    type = INSADMomentumSUPG<<<{"description": "Adds the supg stabilization to the INS momentum equation", "href": "../../source/kernels/INSADMomentumSUPG.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    material_velocity<<<{"description": "A material property describing the velocity"}>>> = relative_velocity
    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
  []
  [temperature_time]
    type = INSADHeatConductionTimeDerivative<<<{"description": "AD Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation for quasi-constant specific heat $c_p$ and the density $\\rho$.", "href": "../../source/kernels/INSADHeatConductionTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    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
  []
  [temperature_advection]
    type = INSADEnergyAdvection<<<{"description": "This class computes the residual and Jacobian contributions for temperature advection for a divergence free velocity field.", "href": "../../source/kernels/INSADEnergyAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    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
  []
  [temperature_mesh_advection]
    type = INSADEnergyMeshAdvection<<<{"description": "This class computes the residual and Jacobian contributions for temperature advection from mesh velocity in an ALE simulation.", "href": "../../source/kernels/INSADEnergyMeshAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    disp_x<<<{"description": "The x displacement"}>>> = disp_x
    disp_y<<<{"description": "The y displacement"}>>> = disp_y
    disp_z<<<{"description": "The z displacement"}>>> = disp_z
    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
  []
  [temperature_conduction]
    type = ADHeatConduction<<<{"description": "Same as `Diffusion` in terms of physics/residual, but the Jacobian is computed using forward automatic differentiation", "href": "../../source/kernels/ADHeatConduction.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    thermal_conductivity<<<{"description": "the name of the thermal conductivity material property"}>>> = 'k'
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
  [temperature_supg]
    type = INSADEnergySUPG<<<{"description": "Adds the supg stabilization to the INS temperature/energy equation", "href": "../../source/kernels/INSADEnergySUPG.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    velocity<<<{"description": "The velocity variable."}>>> = vel
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

This simulation is driven by the rotating laser spot, whose effect is introduced via the GaussianEnergyFluxBC object

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [weld_flux]
    type = GaussianEnergyFluxBC<<<{"description": "Describes an incoming heat flux beam with a Gaussian profile", "href": "../../source/bcs/GaussianEnergyFluxBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    P0<<<{"description": "The total power of the beam."}>>> = 159.96989792079225
    R<<<{"description": "The radius at which the beam intensity falls to $1/e^2$ of its axis value."}>>> = 1.8257418583505537e-4
    x_beam_coord<<<{"description": "The x coordinate of the center of the beam as a function of time. Note that we will pass the origin as the spatial argument to the function; any spatial dependence in the passed-in function will be ignored"}>>> = '2e-4 * cos(t * 2 * pi / ${period})'
    y_beam_coord<<<{"description": "The y coordinate of the center of the beam as a function of time. Note that we will pass the origin as the spatial argument to the function; any spatial dependence in the passed-in function will be ignored"}>>> = '2e-4 * sin(t * 2 * pi / ${period})'
    z_beam_coord<<<{"description": "The z coordinate of the center of the beam as a function of time. Note that we will pass the origin as the spatial argument to the function; any spatial dependence in the passed-in function will be ignored"}>>> = 0
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

In addition to this incoming heat flux, an outgoing radiative heat flux is modeled using FunctionRadiativeBC

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [radiation_flux]
    type = FunctionRadiativeBC<<<{"description": "Boundary condition for radiative heat exchange where the emissivity function is supplied by a Function.", "href": "../../source/bcs/FunctionRadiativeBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    emissivity_function<<<{"description": "Function describing emissivity for radiative boundary condition"}>>> = '1'
    Tinfinity<<<{"description": "Temperature of the body in radiative heat transfer."}>>> = 300
    stefan_boltzmann_constant<<<{"description": "The Stefan-Boltzmann constant."}>>> = ${sb}
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

Evaporation of material from the liquefied material surface helps drive momentum changes at the surface of the condensed phase; this effect is incorporated via the INSADVaporRecoilPressureMomentumFluxBC object:

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [vapor_recoil]
    type = INSADVaporRecoilPressureMomentumFluxBC<<<{"description": "Vapor recoil pressure momentum flux", "href": "../../source/bcs/INSADVaporRecoilPressureMomentumFluxBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

The surface tension of the liquid also contributes to the forces that determine the deformation of the surface. This effect is added using the INSADSurfaceTensionBC object:

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [surface_tension]
    type = INSADSurfaceTensionBC<<<{"description": "Surface tension stresses.", "href": "../../source/bcs/INSADSurfaceTensionBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation.  Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

These surface momentum and velocity changes are then translated into mesh displacement through INSADDisplaceBoundaryBC objects.

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [displace_x_top]
    type = INSADDisplaceBoundaryBC<<<{"description": "Boundary condition for displacing a boundary", "href": "../../source/bcs/INSADDisplaceBoundaryBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_x'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 0
    associated_subdomain<<<{"description": "The subdomain that the boundary nodeset is associated with. This will be passed to the coupled functor for unambigious evaluation (e.g. at the edge of the node-patch where we might run into the intersection of subdomains"}>>> = 0
  []

  [displace_y_top]
    type = INSADDisplaceBoundaryBC<<<{"description": "Boundary condition for displacing a boundary", "href": "../../source/bcs/INSADDisplaceBoundaryBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_y'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 1
    associated_subdomain<<<{"description": "The subdomain that the boundary nodeset is associated with. This will be passed to the coupled functor for unambigious evaluation (e.g. at the edge of the node-patch where we might run into the intersection of subdomains"}>>> = 0
  []

  [displace_z_top]
    type = INSADDisplaceBoundaryBC<<<{"description": "Boundary condition for displacing a boundary", "href": "../../source/bcs/INSADDisplaceBoundaryBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_z'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 2
    associated_subdomain<<<{"description": "The subdomain that the boundary nodeset is associated with. This will be passed to the coupled functor for unambigious evaluation (e.g. at the edge of the node-patch where we might run into the intersection of subdomains"}>>> = 0
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

We also introduce INSADDummyDisplaceBoundaryIntegratedBC objects in order to fill the sparsity dependence of the surface displacement degrees of freedom on the surface velocity degrees of freedom before the Jacobian matrix is assembled prior to executing nodal boundary conditions. This sparsity filling is necessary in order to prevent new nonzero allocations from occurring when the INSADDisplaceBoundaryBC nodal boundary conditions are executed.

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [displace_x_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC<<<{"description": "This object adds Jacobian entries for the boundary displacement dependence on the velocity", "href": "../../source/bcs/INSADDummyDisplaceBoundaryIntegratedBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_x'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 0
  []

  [displace_y_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC<<<{"description": "This object adds Jacobian entries for the boundary displacement dependence on the velocity", "href": "../../source/bcs/INSADDummyDisplaceBoundaryIntegratedBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_y'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 1
  []

  [displace_z_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC<<<{"description": "This object adds Jacobian entries for the boundary displacement dependence on the velocity", "href": "../../source/bcs/INSADDummyDisplaceBoundaryIntegratedBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_z'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 2
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

No-slip boundary conditions are applied at all surfaces other than at the front surface where the laser spot is applied

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [no_slip]
    type = ADVectorFunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $\\vec{u}=\\vec{g}$, where $\\vec{g}$ components are calculated with functions.", "href": "../../source/bcs/ADVectorFunctionDirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom right left top back'
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

The back surface is held at a constant temperature of 300 Kelvin

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [T_cold]
    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"}>>> = T
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 300
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

Zero displacements are applied at the back surface

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [x_no_disp]
    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"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0
  []

  [y_no_disp]
    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"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0
  []

  [z_no_disp]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

Material properties are based on 304L stainless steel

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [steel]
    type = AriaLaserWeld304LStainlessSteel
    temperature = T
    beta = 1e7
    use_displaced_mesh = true
  []

  [steel_boundary]
    type = AriaLaserWeld304LStainlessSteelBoundary
    boundary = 'front'
    temperature = T
    use_displaced_mesh = true
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

The full input is listed below

period=1.25e-3
endtime=${period}
timestep=1.25e-5
surfacetemp=300
sb=5.67e-8

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  temperature = T
[]

[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 3
  xmin = -.35e-3
  xmax = 0.35e-3
  ymin = -.35e-3
  ymax = .35e-3
  zmin = -.7e-3
  zmax = 0
  nx = 2
  ny = 2
  nz = 2
  displacements<<<{"description": "The variables corresponding to the x y z displacements of the mesh.  If this is provided then the displacements will be taken into account during the computation. Creation of the displaced mesh can be suppressed even if this is set by setting 'use_displaced_mesh = false'."}>>> = 'disp_x disp_y disp_z'
  uniform_refine<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 2
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [vel]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE_VEC
  []
  [T]
  []
  [p]
  []
  [disp_x]
  []
  [disp_y]
  []
  [disp_z]
  []
[]

[ICs<<<{"href": "../../syntax/ICs/index.html"}>>>]
  [T]
    type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../source/ics/FunctionIC.html"}>>>
    variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = T
    function<<<{"description": "The initial condition function."}>>> = '(${surfacetemp} - 300) / .7e-3 * z + ${surfacetemp}'
  []
[]

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [disp_x]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
  []
  [disp_y]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
  []
  [disp_z]
    type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../../source/kernels/Diffusion.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
  []
  [mass]
    type = INSADMass<<<{"description": "This class computes the mass equation residual and Jacobian contributions (the latter using automatic differentiation) for the incompressible Navier-Stokes equations.", "href": "../../source/kernels/INSADMass.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
    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
  []
  [mass_pspg]
    type = INSADMassPSPG<<<{"description": "This class adds PSPG stabilization to the mass equation, enabling use of equal order shape functions for pressure and velocity variables", "href": "../../source/kernels/INSADMassPSPG.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
    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
  []
  [momentum_time]
    type = INSADMomentumTimeDerivative<<<{"description": "This class computes the time derivative for the incompressible Navier-Stokes momentum equation.", "href": "../../source/kernels/INSADMomentumTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    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
  []
  [momentum_advection]
    type = INSADMomentumAdvection<<<{"description": "Adds the advective term to the INS momentum equation", "href": "../../source/kernels/INSADMomentumAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    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
  []
  [momentum_mesh_advection]
    type = INSADMomentumMeshAdvection<<<{"description": "Corrects the convective derivative for situations in which the fluid mesh is dynamic.", "href": "../../source/kernels/INSADMomentumMeshAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    disp_x<<<{"description": "The x displacement"}>>> = disp_x
    disp_y<<<{"description": "The y displacement"}>>> = disp_y
    disp_z<<<{"description": "The z displacement"}>>> = disp_z
    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
  []
  [momentum_viscous]
    type = INSADMomentumViscous<<<{"description": "Adds the viscous term to the INS momentum equation", "href": "../../source/kernels/INSADMomentumViscous.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    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
  []
  [momentum_pressure]
    type = INSADMomentumPressure<<<{"description": "Adds the pressure term to the INS momentum equation", "href": "../../source/kernels/INSADMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    pressure<<<{"description": "The pressure"}>>> = p
    integrate_p_by_parts<<<{"description": "Whether to integrate the pressure term by parts"}>>> = true
    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
  []
  [momentum_supg]
    type = INSADMomentumSUPG<<<{"description": "Adds the supg stabilization to the INS momentum equation", "href": "../../source/kernels/INSADMomentumSUPG.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    material_velocity<<<{"description": "A material property describing the velocity"}>>> = relative_velocity
    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
  []
  [temperature_time]
    type = INSADHeatConductionTimeDerivative<<<{"description": "AD Time derivative term $\\rho c_p \\frac{\\partial T}{\\partial t}$ of the heat equation for quasi-constant specific heat $c_p$ and the density $\\rho$.", "href": "../../source/kernels/INSADHeatConductionTimeDerivative.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    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
  []
  [temperature_advection]
    type = INSADEnergyAdvection<<<{"description": "This class computes the residual and Jacobian contributions for temperature advection for a divergence free velocity field.", "href": "../../source/kernels/INSADEnergyAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    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
  []
  [temperature_mesh_advection]
    type = INSADEnergyMeshAdvection<<<{"description": "This class computes the residual and Jacobian contributions for temperature advection from mesh velocity in an ALE simulation.", "href": "../../source/kernels/INSADEnergyMeshAdvection.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    disp_x<<<{"description": "The x displacement"}>>> = disp_x
    disp_y<<<{"description": "The y displacement"}>>> = disp_y
    disp_z<<<{"description": "The z displacement"}>>> = disp_z
    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
  []
  [temperature_conduction]
    type = ADHeatConduction<<<{"description": "Same as `Diffusion` in terms of physics/residual, but the Jacobian is computed using forward automatic differentiation", "href": "../../source/kernels/ADHeatConduction.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    thermal_conductivity<<<{"description": "the name of the thermal conductivity material property"}>>> = 'k'
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
  [temperature_supg]
    type = INSADEnergySUPG<<<{"description": "Adds the supg stabilization to the INS temperature/energy equation", "href": "../../source/kernels/INSADEnergySUPG.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    velocity<<<{"description": "The velocity variable."}>>> = vel
    use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = true
  []
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  [x_no_disp]
    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"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0
  []
  [y_no_disp]
    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"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0
  []
  [z_no_disp]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0
  []
  [no_slip]
    type = ADVectorFunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $\\vec{u}=\\vec{g}$, where $\\vec{g}$ components are calculated with functions.", "href": "../../source/bcs/ADVectorFunctionDirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom right left top back'
  []
  [T_cold]
    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"}>>> = T
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 300
  []
  [radiation_flux]
    type = FunctionRadiativeBC<<<{"description": "Boundary condition for radiative heat exchange where the emissivity function is supplied by a Function.", "href": "../../source/bcs/FunctionRadiativeBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    emissivity_function<<<{"description": "Function describing emissivity for radiative boundary condition"}>>> = '1'
    Tinfinity<<<{"description": "Temperature of the body in radiative heat transfer."}>>> = 300
    stefan_boltzmann_constant<<<{"description": "The Stefan-Boltzmann constant."}>>> = ${sb}
    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
  []
  [weld_flux]
    type = GaussianEnergyFluxBC<<<{"description": "Describes an incoming heat flux beam with a Gaussian profile", "href": "../../source/bcs/GaussianEnergyFluxBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    P0<<<{"description": "The total power of the beam."}>>> = 159.96989792079225
    R<<<{"description": "The radius at which the beam intensity falls to $1/e^2$ of its axis value."}>>> = 1.8257418583505537e-4
    x_beam_coord<<<{"description": "The x coordinate of the center of the beam as a function of time. Note that we will pass the origin as the spatial argument to the function; any spatial dependence in the passed-in function will be ignored"}>>> = '2e-4 * cos(t * 2 * pi / ${period})'
    y_beam_coord<<<{"description": "The y coordinate of the center of the beam as a function of time. Note that we will pass the origin as the spatial argument to the function; any spatial dependence in the passed-in function will be ignored"}>>> = '2e-4 * sin(t * 2 * pi / ${period})'
    z_beam_coord<<<{"description": "The z coordinate of the center of the beam as a function of time. Note that we will pass the origin as the spatial argument to the function; any spatial dependence in the passed-in function will be ignored"}>>> = 0
    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
  []
  [vapor_recoil]
    type = INSADVaporRecoilPressureMomentumFluxBC<<<{"description": "Vapor recoil pressure momentum flux", "href": "../../source/bcs/INSADVaporRecoilPressureMomentumFluxBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    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
  []
  [surface_tension]
    type = INSADSurfaceTensionBC<<<{"description": "Surface tension stresses.", "href": "../../source/bcs/INSADSurfaceTensionBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    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
  []
  [displace_x_top]
    type = INSADDisplaceBoundaryBC<<<{"description": "Boundary condition for displacing a boundary", "href": "../../source/bcs/INSADDisplaceBoundaryBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_x'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 0
    associated_subdomain<<<{"description": "The subdomain that the boundary nodeset is associated with. This will be passed to the coupled functor for unambigious evaluation (e.g. at the edge of the node-patch where we might run into the intersection of subdomains"}>>> = 0
  []
  [displace_y_top]
    type = INSADDisplaceBoundaryBC<<<{"description": "Boundary condition for displacing a boundary", "href": "../../source/bcs/INSADDisplaceBoundaryBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_y'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 1
    associated_subdomain<<<{"description": "The subdomain that the boundary nodeset is associated with. This will be passed to the coupled functor for unambigious evaluation (e.g. at the edge of the node-patch where we might run into the intersection of subdomains"}>>> = 0
  []
  [displace_z_top]
    type = INSADDisplaceBoundaryBC<<<{"description": "Boundary condition for displacing a boundary", "href": "../../source/bcs/INSADDisplaceBoundaryBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_z'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 2
    associated_subdomain<<<{"description": "The subdomain that the boundary nodeset is associated with. This will be passed to the coupled functor for unambigious evaluation (e.g. at the edge of the node-patch where we might run into the intersection of subdomains"}>>> = 0
  []
  [displace_x_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC<<<{"description": "This object adds Jacobian entries for the boundary displacement dependence on the velocity", "href": "../../source/bcs/INSADDummyDisplaceBoundaryIntegratedBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_x'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 0
  []
  [displace_y_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC<<<{"description": "This object adds Jacobian entries for the boundary displacement dependence on the velocity", "href": "../../source/bcs/INSADDummyDisplaceBoundaryIntegratedBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_y'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 1
  []
  [displace_z_top_dummy]
    type = INSADDummyDisplaceBoundaryIntegratedBC<<<{"description": "This object adds Jacobian entries for the boundary displacement dependence on the velocity", "href": "../../source/bcs/INSADDummyDisplaceBoundaryIntegratedBC.html"}>>>
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'front'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_z'
    velocity<<<{"description": "The velocity at which to displace. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'vel'
    component<<<{"description": "What component of velocity/displacement this object is acting on."}>>> = 2
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [ins_mat]
    type = INSADStabilized3Eqn<<<{"description": "This is the material class used to compute the stabilization parameter tau for momentum and tau_energy for the energy equation.", "href": "../../source/materials/INSADStabilized3Eqn.html"}>>>
    velocity<<<{"description": "The velocity"}>>> = vel
    pressure<<<{"description": "The pressure"}>>> = p
    temperature<<<{"description": "The temperature"}>>> = T
    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
  []
  [steel]
    type = AriaLaserWeld304LStainlessSteel
    temperature = T
    beta = 1e7
    use_displaced_mesh = true
  []
  [steel_boundary]
    type = AriaLaserWeld304LStainlessSteelBoundary
    boundary = 'front'
    temperature = T
    use_displaced_mesh = true
  []
  [const]
    type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../source/materials/GenericConstantMaterial.html"}>>>
    prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'abs sb_constant'
    prop_values<<<{"description": "The values associated with the named properties"}>>> = '1 ${sb}'
    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
  []
[]

[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
    petsc_options_iname<<<{"description": "Names of PETSc name/value pairs"}>>> = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type'
    petsc_options_value<<<{"description": "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\""}>>> = 'lu       NONZERO               strumpack'
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = Transient
  end_time = ${endtime}
  dtmin = 1e-8
  dtmax = ${timestep}
  petsc_options = '-snes_converged_reason -ksp_converged_reason -options_left'
  solve_type = 'NEWTON'
  line_search = 'none'
  nl_max_its = 12
  l_max_its = 100
  [TimeStepper<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
    type = IterationAdaptiveDT
    optimal_iterations = 7
    dt = ${timestep}
    linear_iteration_ratio = 1e6
    growth_factor = 1.5
  []
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  [exodus]
    type = Exodus<<<{"description": "Object for output data in the Exodus format", "href": "../../source/outputs/Exodus.html"}>>>
    output_material_properties<<<{"description": "Flag indicating if material properties should be output"}>>> = true
    show_material_properties<<<{"description": "List of material properties that should be written to the output"}>>> = 'mu'
  []
  checkpoint<<<{"description": "Create checkpoint files using the default options."}>>> = true
  perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = true
[]

[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
[]

[Adaptivity<<<{"href": "../../syntax/Adaptivity/index.html"}>>>]
  marker<<<{"description": "The name of the Marker to use to actually adapt the mesh."}>>> = combo
  max_h_level<<<{"description": "Maximum number of times a single element can be refined. If 0 then infinite."}>>> = 4

  [Indicators<<<{"href": "../../syntax/Adaptivity/Indicators/index.html"}>>>]
    [error_T]
      type = GradientJumpIndicator<<<{"description": "Compute the jump of the solution gradient across element boundaries.", "href": "../../source/indicators/GradientJumpIndicator.html"}>>>
      variable<<<{"description": "The name of the variable that this side indicator applies to"}>>> = T
    []
    [error_dispz]
      type = GradientJumpIndicator<<<{"description": "Compute the jump of the solution gradient across element boundaries.", "href": "../../source/indicators/GradientJumpIndicator.html"}>>>
      variable<<<{"description": "The name of the variable that this side indicator applies to"}>>> = disp_z
    []
  []

  [Markers<<<{"href": "../../syntax/Adaptivity/Markers/index.html"}>>>]
    [errorfrac_T]
      type = ErrorFractionMarker<<<{"description": "Marks elements for refinement or coarsening based on the fraction of the min/max error from the supplied indicator.", "href": "../../source/markers/ErrorFractionMarker.html"}>>>
      refine<<<{"description": "Elements within this percentage of the max error will be refined.  Must be between 0 and 1!"}>>> = 0.4
      coarsen<<<{"description": "Elements within this percentage of the min error will be coarsened.  Must be between 0 and 1!"}>>> = 0.2
      indicator<<<{"description": "The name of the Indicator that this Marker uses."}>>> = error_T
    []
    [errorfrac_dispz]
      type = ErrorFractionMarker<<<{"description": "Marks elements for refinement or coarsening based on the fraction of the min/max error from the supplied indicator.", "href": "../../source/markers/ErrorFractionMarker.html"}>>>
      refine<<<{"description": "Elements within this percentage of the max error will be refined.  Must be between 0 and 1!"}>>> = 0.4
      coarsen<<<{"description": "Elements within this percentage of the min error will be coarsened.  Must be between 0 and 1!"}>>> = 0.2
      indicator<<<{"description": "The name of the Indicator that this Marker uses."}>>> = error_dispz
    []
    [combo]
      type = ComboMarker<<<{"description": "A marker that converts many markers into a single marker by considering the maximum value of the listed markers (i.e., refinement takes precedent).", "href": "../../source/markers/ComboMarker.html"}>>>
      markers<<<{"description": "A list of marker names to combine into a single marker."}>>> = 'errorfrac_T errorfrac_dispz'
    []
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [num_dofs]
    type = NumDOFs<<<{"description": "Return the number of Degrees of freedom from either the NL, Aux or both systems.", "href": "../../source/postprocessors/NumDOFs.html"}>>>
    system<<<{"description": "The system(s) for which you want to retrieve the number of DOFs (NL, AUX, ALL). Default == ALL"}>>> = 'NL'
  []
  [nl]
    type = NumNonlinearIterations<<<{"description": "Outputs the number of nonlinear iterations", "href": "../../source/postprocessors/NumNonlinearIterations.html"}>>>
  []
  [tot_nl]
    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"}>>> = 'nl'
  []
[]
(modules/navier_stokes/examples/laser-welding/3d.i)

References

  1. Fande Kong and Xiao-Chuan Cai. A scalable nonlinear fluid–structure interaction solver based on a schwarz preconditioner with isogeometric unstructured coarse spaces in 3d. Journal of Computational Physics, 340:498–518, 2017.[BibTeX]
  2. Alexander Lindsay, Roy Stogner, Derek Gaston, Daniel Schwen, Christopher Matthews, Wen Jiang, Larry K Aagesen, Robert Carlsen, Fande Kong, Andrew Slaughter, and others. Automatic differentiation in metaphysicl and its applications in moose. Nuclear Technology, 207(7):905–922, 2021.[BibTeX]
  3. David R Noble, Patrick K Notz, Mario J Martinez, and Andrew Michael Kraynik. Use of aria to simulate laser weld pool dynamics for neutron generator production. Technical Report, Sandia National Laboratories (SNL), Albuquerque, NM, and Livermore, CA …, 2007.[BibTeX]