- disp_xThe x displacement variable
C++ Type:std::vector<VariableName>
Description:The x displacement variable
- disp_yThe y displacement variable
C++ Type:std::vector<VariableName>
Description:The y displacement variable
- primary_boundaryThe name of the primary boundary sideset.
C++ Type:BoundaryName
Description:The name of the primary boundary sideset.
- primary_subdomainThe name of the primary subdomain.
C++ Type:SubdomainName
Description:The name of the primary subdomain.
- secondary_boundaryThe name of the secondary boundary sideset.
C++ Type:BoundaryName
Description:The name of the secondary boundary sideset.
- secondary_subdomainThe name of the secondary subdomain.
C++ Type:SubdomainName
Description:The name of the secondary subdomain.
ComputeWeightedGapLMMechanicalContact
The Karush-Kuhn-Tucker conditions of mechanical contact are:
where is the gap and is the contact pressure, a Lagrange multipler variable living on the secondary face. Per (Wohlmuth, 2011) and (Popp and Wall, 2014), the variationally consistent, discretized version of the KKT conditions are:
where indicates the normal direction, denotes the j'th secondary contact interface node, and is the discrete weighted gap, computed by:
where denotes the secondary contact interface, is the j'th lagrange multiplier test function, and is the discretized version of the gap function.
The ComputeWeightedGapLMMechanicalContact
object computes the weighted gap. It does not apply** the KKT conditions. That is done with the ApplyPenetrationConstraintLMMechanicalContact object. Consequently, the two objects must always be used in conjunction.
Computes the weighted gap that will later be used to enforce the zero-penetration mechanical contact conditions
Input Parameters
- compute_lm_residualsTrueWhether to compute Lagrange Multiplier residuals
Default:True
C++ Type:bool
Options:
Description:Whether to compute Lagrange Multiplier residuals
- compute_primal_residualsTrueWhether to compute residuals for the primal variable.
Default:True
C++ Type:bool
Options:
Description:Whether to compute residuals for the primal variable.
- periodicFalseWhether this constraint is going to be used to enforce a periodic condition. This has the effect of changing the normals vector for projection from outward to inward facing
Default:False
C++ Type:bool
Options:
Description:Whether this constraint is going to be used to enforce a periodic condition. This has the effect of changing the normals vector for projection from outward to inward facing
- variableThe name of the lagrange multiplier variable that this constraint is applied to. This parameter may not be supplied in the case of using penalty methods for example
C++ Type:NonlinearVariableName
Options:
Description:The name of the lagrange multiplier variable that this constraint is applied to. This parameter may not be supplied in the case of using penalty methods for example
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Options:
Description:Determines whether this object is calculated using an implicit or explicit form
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Options:
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Options:
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Options:
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Options:
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- (modules/combined/test/tests/gap_heat_transfer_mortar/small-2d/small.i)
- (modules/contact/test/tests/bouncing-block-contact/frictionless-weighted-gap.i)
- (modules/combined/test/tests/gap_heat_transfer_mortar/finite-2d-rz/finite.i)
- (modules/combined/test/tests/gap_heat_transfer_mortar/finite-2d/finite.i)
- (modules/combined/test/tests/gap_heat_transfer_mortar/small-2d-rz/small.i)
- (modules/contact/test/tests/bouncing-block-contact/frictionless-weighted-gap-mixed-basis.i)
References
- Alexander Popp and WA Wall.
Dual mortar methods for computational contact mechanics–overview and recent developments.
GAMM-Mitteilungen, 37(1):66–84, 2014.[BibTeX]
@article{popp2014dual, author = "Popp, Alexander and Wall, WA", title = "Dual mortar methods for computational contact mechanics--overview and recent developments", journal = "GAMM-Mitteilungen", volume = "37", number = "1", pages = "66--84", year = "2014", publisher = "Wiley Online Library" }
- Barbara Wohlmuth.
Variationally consistent discretization schemes and numerical algorithms for contact problems.
Acta Numerica, 20:569–734, 2011.[BibTeX]
@article{wohlmuth2011variationally, author = "Wohlmuth, Barbara", title = "Variationally consistent discretization schemes and numerical algorithms for contact problems", journal = "Acta Numerica", volume = "20", pages = "569--734", year = "2011" }
(modules/combined/test/tests/gap_heat_transfer_mortar/small-2d/small.i)
E_block = 1e7
E_plank = 1e7
elem = QUAD4
order = FIRST
name = 'small'
[Mesh]
patch_size = 80
patch_update_strategy = auto
[plank]
type = GeneratedMeshGenerator
dim = 2
xmin = -0.3
xmax = 0.3
ymin = -10
ymax = 10
nx = 2
ny = 67
elem_type = ${elem}
boundary_name_prefix = plank
[]
[plank_id]
type = SubdomainIDGenerator
input = plank
subdomain_id = 1
[]
[block]
type = GeneratedMeshGenerator
dim = 2
xmin = 0.31
xmax = 0.91
ymin = 7.7
ymax = 8.5
nx = 3
ny = 4
elem_type = ${elem}
boundary_name_prefix = block
boundary_id_offset = 10
[]
[block_id]
type = SubdomainIDGenerator
input = block
subdomain_id = 2
[]
[combined]
type = MeshCollectionGenerator
inputs = 'plank_id block_id'
[]
[block_rename]
type = RenameBlockGenerator
input = combined
old_block_id = '1 2'
new_block_name = 'plank block'
[]
[secondary]
input = block_rename
type = LowerDBlockFromSidesetGenerator
sidesets = 'block_left'
new_block_id = '30'
new_block_name = 'frictionless_secondary_subdomain'
[]
[primary]
input = secondary
type = LowerDBlockFromSidesetGenerator
sidesets = 'plank_right'
new_block_id = '20'
new_block_name = 'frictionless_primary_subdomain'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Variables]
[disp_x]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[disp_y]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[temp]
order = ${order}
block = 'plank block'
scaling = 1e-1
[]
[thermal_lm]
order = ${order}
block = 'frictionless_secondary_subdomain'
scaling = 1e-7
[]
[frictionless_normal_lm]
order = FIRST
block = 'frictionless_secondary_subdomain'
scaling = 1e3
use_dual = true
[]
[]
[Modules/TensorMechanics/Master]
[action]
generate_output = 'stress_xx stress_yy stress_zz vonmises_stress hydrostatic_stress strain_xx strain_yy strain_zz'
block = 'plank block'
use_automatic_differentiation = true
[]
[]
[Kernels]
[hc]
type = ADHeatConduction
variable = temp
use_displaced_mesh = true
block = 'plank block'
[]
[]
[Constraints]
[weighted_gap_lm]
type = ComputeWeightedGapLMMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[]
[ncp_lm]
type = ApplyPenetrationConstraintLMMechanicalContact
secondary = block_left
primary = plank_right
variable = frictionless_normal_lm
primary_variable = disp_x
[]
[normal_x]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_x
component = x
use_displaced_mesh = true
compute_lm_residuals = false
[]
[normal_y]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_y
component = y
use_displaced_mesh = true
compute_lm_residuals = false
[]
[thermal_contact]
type = GapConductanceConstraint
variable = thermal_lm
secondary_variable = temp
k = 1
use_displaced_mesh = true
primary_boundary = plank_right
primary_subdomain = frictionless_primary_subdomain
secondary_boundary = block_left
secondary_subdomain = frictionless_secondary_subdomain
displacements = 'disp_x disp_y'
[]
[]
[BCs]
[left_temp]
type = DirichletBC
variable = temp
boundary = 'plank_left'
value = 400
[]
[right_temp]
type = DirichletBC
variable = temp
boundary = 'block_right'
value = 300
[]
[left_x]
type = DirichletBC
variable = disp_x
boundary = plank_left
value = 0.0
[]
[left_y]
type = DirichletBC
variable = disp_y
boundary = plank_bottom
value = 0.0
[]
[right_x]
type = ADFunctionDirichletBC
variable = disp_x
boundary = block_right
function = '-0.04*sin(4*(t+1.5))+0.02'
[]
[right_y]
type = ADFunctionDirichletBC
variable = disp_y
boundary = block_right
function = '-t'
[]
[]
[Materials]
[plank]
type = ADComputeIsotropicElasticityTensor
block = 'plank'
poissons_ratio = 0.3
youngs_modulus = ${E_plank}
[]
[block]
type = ADComputeIsotropicElasticityTensor
block = 'block'
poissons_ratio = 0.3
youngs_modulus = ${E_block}
[]
[stress]
type = ADComputeLinearElasticStress
block = 'plank block'
[]
[heat_plank]
type = ADHeatConductionMaterial
block = plank
thermal_conductivity = 2
specific_heat = 1
[]
[heat_block]
type = ADHeatConductionMaterial
block = block
thermal_conductivity = 1
specific_heat = 1
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options = '-snes_converged_reason -ksp_converged_reason'
petsc_options_iname = '-pc_type -mat_mffd_err -pc_factor_shift_type -pc_factor_shift_amount -snes_max_it'
petsc_options_value = 'lu 1e-5 NONZERO 1e-15 20'
end_time = 13.5
dt = 0.1
dtmin = 0.1
timestep_tolerance = 1e-6
line_search = 'none'
[]
[Postprocessors]
[nl_its]
type = NumNonlinearIterations
[]
[total_nl_its]
type = CumulativeValuePostprocessor
postprocessor = nl_its
[]
[l_its]
type = NumLinearIterations
[]
[total_l_its]
type = CumulativeValuePostprocessor
postprocessor = l_its
[]
[contact]
type = ContactDOFSetSize
variable = frictionless_normal_lm
subdomain = frictionless_secondary_subdomain
[]
[avg_hydro]
type = ElementAverageValue
variable = hydrostatic_stress
block = 'block'
[]
[avg_temp]
type = ElementAverageValue
variable = temp
block = 'block'
[]
[max_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
[]
[min_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
value_type = min
[]
[avg_vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 'block'
[]
[max_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
[]
[min_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
value_type = min
[]
[]
[Outputs]
exodus = true
file_base = ${name}
checkpoint = true
[comp]
type = CSV
show = 'contact avg_temp'
[]
[out]
type = CSV
file_base = '${name}_out'
[]
[]
[Debug]
show_var_residual_norms = true
[]
(modules/contact/test/tests/bouncing-block-contact/frictionless-weighted-gap.i)
starting_point = 2e-1
offset = 1e-2
[GlobalParams]
displacements = 'disp_x disp_y'
diffusivity = 1e0
scaling = 1e0
[]
[Mesh]
file = long-bottom-block-1elem-blocks.e
[]
[Variables]
[./disp_x]
block = '1 2'
[../]
[./disp_y]
block = '1 2'
[../]
[./normal_lm]
block = 3
[../]
[]
[ICs]
[./disp_y]
block = 2
variable = disp_y
value = ${fparse starting_point + offset}
type = ConstantIC
[../]
[]
[Kernels]
[./disp_x]
type = MatDiffusion
variable = disp_x
[../]
[./disp_y]
type = MatDiffusion
variable = disp_y
[../]
[]
[Constraints]
[./weighted_gap_lm]
type = ComputeWeightedGapLMMechanicalContact
primary_boundary = 20
secondary_boundary = 10
primary_subdomain = 4
secondary_subdomain = 3
variable = normal_lm
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[../]
[./ncp_lm]
type = ApplyPenetrationConstraintLMMechanicalContact
secondary = 10
primary = 20
variable = normal_lm
primary_variable = disp_x
c = 1
[../]
[normal_x]
type = NormalMortarMechanicalContact
primary_boundary = 20
secondary_boundary = 10
primary_subdomain = 4
secondary_subdomain = 3
variable = normal_lm
secondary_variable = disp_x
component = x
use_displaced_mesh = true
compute_lm_residuals = false
[]
[normal_y]
type = NormalMortarMechanicalContact
primary_boundary = 20
secondary_boundary = 10
primary_subdomain = 4
secondary_subdomain = 3
variable = normal_lm
secondary_variable = disp_y
component = y
use_displaced_mesh = true
compute_lm_residuals = false
[]
[]
[BCs]
[./botx]
type = DirichletBC
variable = disp_x
boundary = 40
value = 0.0
[../]
[./boty]
type = DirichletBC
variable = disp_y
boundary = 40
value = 0.0
[../]
[./topy]
type = FunctionDirichletBC
variable = disp_y
boundary = 30
function = '${starting_point} * cos(2 * pi / 40 * t) + ${offset}'
[../]
[./leftx]
type = FunctionDirichletBC
variable = disp_x
boundary = 50
function = '1e-2 * t'
[../]
[]
[Executioner]
type = Transient
end_time = 200
dt = 5
dtmin = .1
solve_type = 'PJFNK'
petsc_options = '-snes_converged_reason -ksp_converged_reason -pc_svd_monitor -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount -mat_mffd_err'
petsc_options_value = 'lu NONZERO 1e-15 1e-5'
l_max_its = 30
nl_max_its = 20
line_search = 'none'
snesmf_reuse_base = false
[]
[Debug]
show_var_residual_norms = true
[]
[Outputs]
exodus = true
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
active = 'num_nl cumulative contact'
[./num_nl]
type = NumNonlinearIterations
[../]
[./cumulative]
type = CumulativeValuePostprocessor
postprocessor = num_nl
[../]
[contact]
type = ContactDOFSetSize
variable = normal_lm
subdomain = '3'
execute_on = 'nonlinear timestep_end'
[]
[]
(modules/combined/test/tests/gap_heat_transfer_mortar/finite-2d-rz/finite.i)
E_block = 1e7
E_plank = 1e7
elem = QUAD4
order = FIRST
name = 'finite'
[Mesh]
patch_size = 80
patch_update_strategy = auto
[plank]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 0.6
ymin = -10
ymax = 10
nx = 2
ny = 67
elem_type = ${elem}
boundary_name_prefix = plank
[]
[plank_id]
type = SubdomainIDGenerator
input = plank
subdomain_id = 1
[]
[block]
type = GeneratedMeshGenerator
dim = 2
xmin = 0.61
xmax = 1.21
ymin = 7.7
ymax = 8.5
nx = 3
ny = 4
elem_type = ${elem}
boundary_name_prefix = block
boundary_id_offset = 10
[]
[block_id]
type = SubdomainIDGenerator
input = block
subdomain_id = 2
[]
[combined]
type = MeshCollectionGenerator
inputs = 'plank_id block_id'
[]
[block_rename]
type = RenameBlockGenerator
input = combined
old_block_id = '1 2'
new_block_name = 'plank block'
[]
[secondary]
input = block_rename
type = LowerDBlockFromSidesetGenerator
sidesets = 'block_left'
new_block_id = '30'
new_block_name = 'frictionless_secondary_subdomain'
[]
[primary]
input = secondary
type = LowerDBlockFromSidesetGenerator
sidesets = 'plank_right'
new_block_id = '20'
new_block_name = 'frictionless_primary_subdomain'
[]
[]
[Problem]
coord_type = RZ
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Variables]
[disp_x]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[disp_y]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[temp]
order = ${order}
block = 'plank block'
scaling = 1e-1
[]
[thermal_lm]
order = ${order}
block = 'frictionless_secondary_subdomain'
scaling = 1e-7
[]
[frictionless_normal_lm]
order = FIRST
block = 'frictionless_secondary_subdomain'
scaling = 1e3
use_dual = true
[]
[]
[Modules/TensorMechanics/Master]
[action]
generate_output = 'stress_xx stress_yy stress_zz vonmises_stress hydrostatic_stress strain_xx strain_yy strain_zz'
block = 'plank block'
use_automatic_differentiation = true
strain = FINITE
[]
[]
[Kernels]
[hc]
type = ADHeatConduction
variable = temp
use_displaced_mesh = true
block = 'plank block'
[]
[]
[Constraints]
[weighted_gap_lm]
type = ComputeWeightedGapLMMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[]
[ncp_lm]
type = ApplyPenetrationConstraintLMMechanicalContact
secondary = block_left
primary = plank_right
variable = frictionless_normal_lm
primary_variable = disp_x
[]
[normal_x]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_x
component = x
use_displaced_mesh = true
compute_lm_residuals = false
[]
[normal_y]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_y
component = y
use_displaced_mesh = true
compute_lm_residuals = false
[]
[thermal_contact]
type = GapConductanceConstraint
variable = thermal_lm
secondary_variable = temp
k = 1
use_displaced_mesh = true
primary_boundary = plank_right
primary_subdomain = frictionless_primary_subdomain
secondary_boundary = block_left
secondary_subdomain = frictionless_secondary_subdomain
displacements = 'disp_x disp_y'
[]
[]
[BCs]
[left_temp]
type = DirichletBC
variable = temp
boundary = 'plank_left'
value = 400
[]
[right_temp]
type = DirichletBC
variable = temp
boundary = 'block_right'
value = 300
[]
[left_x]
type = DirichletBC
variable = disp_x
boundary = plank_left
value = 0.0
[]
[left_y]
type = DirichletBC
variable = disp_y
boundary = plank_bottom
value = 0.0
[]
[right_x]
type = ADFunctionDirichletBC
variable = disp_x
boundary = block_right
function = '-0.04*sin(4*(t+1.5))+0.02'
preset = false
[]
[right_y]
type = ADFunctionDirichletBC
variable = disp_y
boundary = block_right
function = '-t'
preset = false
[]
[]
[Materials]
[plank]
type = ADComputeIsotropicElasticityTensor
block = 'plank'
poissons_ratio = 0.3
youngs_modulus = ${E_plank}
[]
[block]
type = ADComputeIsotropicElasticityTensor
block = 'block'
poissons_ratio = 0.3
youngs_modulus = ${E_block}
[]
[stress]
type = ADComputeFiniteStrainElasticStress
block = 'plank block'
[]
[heat_plank]
type = ADHeatConductionMaterial
block = plank
thermal_conductivity = 2
specific_heat = 1
[]
[heat_block]
type = ADHeatConductionMaterial
block = block
thermal_conductivity = 1
specific_heat = 1
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options = '-snes_converged_reason -ksp_converged_reason'
petsc_options_iname = '-pc_type -mat_mffd_err -pc_factor_shift_type -pc_factor_shift_amount -snes_max_it'
petsc_options_value = 'lu 1e-5 NONZERO 1e-15 20'
end_time = 13.5
dt = 0.1
dtmin = 0.1
timestep_tolerance = 1e-6
line_search = 'none'
[]
[Postprocessors]
[nl_its]
type = NumNonlinearIterations
[]
[total_nl_its]
type = CumulativeValuePostprocessor
postprocessor = nl_its
[]
[l_its]
type = NumLinearIterations
[]
[total_l_its]
type = CumulativeValuePostprocessor
postprocessor = l_its
[]
[contact]
type = ContactDOFSetSize
variable = frictionless_normal_lm
subdomain = frictionless_secondary_subdomain
[]
[avg_hydro]
type = ElementAverageValue
variable = hydrostatic_stress
block = 'block'
[]
[avg_temp]
type = ElementAverageValue
variable = temp
block = 'block'
[]
[max_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
[]
[min_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
value_type = min
[]
[avg_vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 'block'
[]
[max_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
[]
[min_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
value_type = min
[]
[]
[Outputs]
exodus = true
file_base = ${name}
checkpoint = true
[comp]
type = CSV
show = 'contact avg_temp'
[]
[out]
type = CSV
file_base = '${name}_out'
[]
[]
[Debug]
show_var_residual_norms = true
[]
(modules/combined/test/tests/gap_heat_transfer_mortar/finite-2d/finite.i)
E_block = 1e7
E_plank = 1e7
elem = QUAD4
order = FIRST
name = 'finite'
[Mesh]
patch_size = 80
patch_update_strategy = auto
[plank]
type = GeneratedMeshGenerator
dim = 2
xmin = -0.3
xmax = 0.3
ymin = -10
ymax = 10
nx = 2
ny = 67
elem_type = ${elem}
boundary_name_prefix = plank
[]
[plank_id]
type = SubdomainIDGenerator
input = plank
subdomain_id = 1
[]
[block]
type = GeneratedMeshGenerator
dim = 2
xmin = 0.31
xmax = 0.91
ymin = 7.7
ymax = 8.5
nx = 3
ny = 4
elem_type = ${elem}
boundary_name_prefix = block
boundary_id_offset = 10
[]
[block_id]
type = SubdomainIDGenerator
input = block
subdomain_id = 2
[]
[combined]
type = MeshCollectionGenerator
inputs = 'plank_id block_id'
[]
[block_rename]
type = RenameBlockGenerator
input = combined
old_block_id = '1 2'
new_block_name = 'plank block'
[]
[secondary]
input = block_rename
type = LowerDBlockFromSidesetGenerator
sidesets = 'block_left'
new_block_id = '30'
new_block_name = 'frictionless_secondary_subdomain'
[]
[primary]
input = secondary
type = LowerDBlockFromSidesetGenerator
sidesets = 'plank_right'
new_block_id = '20'
new_block_name = 'frictionless_primary_subdomain'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Variables]
[disp_x]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[disp_y]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[temp]
order = ${order}
block = 'plank block'
scaling = 1e-1
[]
[thermal_lm]
order = ${order}
block = 'frictionless_secondary_subdomain'
scaling = 1e-7
[]
[frictionless_normal_lm]
order = FIRST
block = 'frictionless_secondary_subdomain'
scaling = 1e3
use_dual = true
[]
[]
[Modules/TensorMechanics/Master]
[action]
generate_output = 'stress_xx stress_yy stress_zz vonmises_stress hydrostatic_stress strain_xx strain_yy strain_zz'
block = 'plank block'
use_automatic_differentiation = true
strain = FINITE
[]
[]
[Kernels]
[hc]
type = ADHeatConduction
variable = temp
use_displaced_mesh = true
block = 'plank block'
[]
[]
[Constraints]
[weighted_gap_lm]
type = ComputeWeightedGapLMMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[]
[ncp_lm]
type = ApplyPenetrationConstraintLMMechanicalContact
secondary = block_left
primary = plank_right
variable = frictionless_normal_lm
primary_variable = disp_x
[]
[normal_x]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_x
component = x
use_displaced_mesh = true
compute_lm_residuals = false
[]
[normal_y]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_y
component = y
use_displaced_mesh = true
compute_lm_residuals = false
[]
[thermal_contact]
type = GapConductanceConstraint
variable = thermal_lm
secondary_variable = temp
k = 1
use_displaced_mesh = true
primary_boundary = plank_right
primary_subdomain = frictionless_primary_subdomain
secondary_boundary = block_left
secondary_subdomain = frictionless_secondary_subdomain
displacements = 'disp_x disp_y'
[]
[]
[BCs]
[left_temp]
type = DirichletBC
variable = temp
boundary = 'plank_left'
value = 400
[]
[right_temp]
type = DirichletBC
variable = temp
boundary = 'block_right'
value = 300
[]
[left_x]
type = DirichletBC
variable = disp_x
boundary = plank_left
value = 0.0
[]
[left_y]
type = DirichletBC
variable = disp_y
boundary = plank_bottom
value = 0.0
[]
[right_x]
type = ADFunctionDirichletBC
variable = disp_x
boundary = block_right
function = '-0.04*sin(4*(t+1.5))+0.02'
preset = false
[]
[right_y]
type = ADFunctionDirichletBC
variable = disp_y
boundary = block_right
function = '-t'
preset = false
[]
[]
[Materials]
[plank]
type = ADComputeIsotropicElasticityTensor
block = 'plank'
poissons_ratio = 0.3
youngs_modulus = ${E_plank}
[]
[block]
type = ADComputeIsotropicElasticityTensor
block = 'block'
poissons_ratio = 0.3
youngs_modulus = ${E_block}
[]
[stress]
type = ADComputeFiniteStrainElasticStress
block = 'plank block'
[]
[heat_plank]
type = ADHeatConductionMaterial
block = plank
thermal_conductivity = 2
specific_heat = 1
[]
[heat_block]
type = ADHeatConductionMaterial
block = block
thermal_conductivity = 1
specific_heat = 1
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options = '-snes_converged_reason -ksp_converged_reason'
petsc_options_iname = '-pc_type -mat_mffd_err -pc_factor_shift_type -pc_factor_shift_amount -snes_max_it'
petsc_options_value = 'lu 1e-5 NONZERO 1e-15 20'
end_time = 13.5
dt = 0.1
dtmin = 0.1
timestep_tolerance = 1e-6
line_search = 'none'
[]
[Postprocessors]
[nl_its]
type = NumNonlinearIterations
[]
[total_nl_its]
type = CumulativeValuePostprocessor
postprocessor = nl_its
[]
[l_its]
type = NumLinearIterations
[]
[total_l_its]
type = CumulativeValuePostprocessor
postprocessor = l_its
[]
[contact]
type = ContactDOFSetSize
variable = frictionless_normal_lm
subdomain = frictionless_secondary_subdomain
[]
[avg_hydro]
type = ElementAverageValue
variable = hydrostatic_stress
block = 'block'
[]
[avg_temp]
type = ElementAverageValue
variable = temp
block = 'block'
[]
[max_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
[]
[min_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
value_type = min
[]
[avg_vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 'block'
[]
[max_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
[]
[min_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
value_type = min
[]
[]
[Outputs]
exodus = true
file_base = ${name}
checkpoint = true
[comp]
type = CSV
show = 'contact avg_temp'
[]
[out]
type = CSV
file_base = '${name}_out'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Debug]
show_var_residual_norms = true
[]
(modules/combined/test/tests/gap_heat_transfer_mortar/small-2d-rz/small.i)
E_block = 1e7
E_plank = 1e7
elem = QUAD4
order = FIRST
name = 'small'
[Mesh]
patch_size = 80
patch_update_strategy = auto
[plank]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 0.6
ymin = -10
ymax = 10
nx = 2
ny = 67
elem_type = ${elem}
boundary_name_prefix = plank
[]
[plank_id]
type = SubdomainIDGenerator
input = plank
subdomain_id = 1
[]
[block]
type = GeneratedMeshGenerator
dim = 2
xmin = 0.61
xmax = 1.21
ymin = 7.7
ymax = 8.5
nx = 3
ny = 4
elem_type = ${elem}
boundary_name_prefix = block
boundary_id_offset = 10
[]
[block_id]
type = SubdomainIDGenerator
input = block
subdomain_id = 2
[]
[combined]
type = MeshCollectionGenerator
inputs = 'plank_id block_id'
[]
[block_rename]
type = RenameBlockGenerator
input = combined
old_block_id = '1 2'
new_block_name = 'plank block'
[]
[secondary]
input = block_rename
type = LowerDBlockFromSidesetGenerator
sidesets = 'block_left'
new_block_id = '30'
new_block_name = 'frictionless_secondary_subdomain'
[]
[primary]
input = secondary
type = LowerDBlockFromSidesetGenerator
sidesets = 'plank_right'
new_block_id = '20'
new_block_name = 'frictionless_primary_subdomain'
[]
[]
[Problem]
coord_type = RZ
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Variables]
[disp_x]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[disp_y]
order = ${order}
block = 'plank block'
scaling = ${fparse 2.0 / (E_plank + E_block)}
[]
[temp]
order = ${order}
block = 'plank block'
scaling = 1e-1
[]
[thermal_lm]
order = ${order}
block = 'frictionless_secondary_subdomain'
scaling = 1e-7
[]
[frictionless_normal_lm]
order = FIRST
block = 'frictionless_secondary_subdomain'
scaling = 1e3
use_dual = true
[]
[]
[Modules/TensorMechanics/Master]
[action]
generate_output = 'stress_xx stress_yy stress_zz vonmises_stress hydrostatic_stress strain_xx strain_yy strain_zz'
block = 'plank block'
use_automatic_differentiation = true
[]
[]
[Kernels]
[hc]
type = ADHeatConduction
variable = temp
use_displaced_mesh = true
block = 'plank block'
[]
[]
[Constraints]
[weighted_gap_lm]
type = ComputeWeightedGapLMMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[]
[ncp_lm]
type = ApplyPenetrationConstraintLMMechanicalContact
secondary = block_left
primary = plank_right
variable = frictionless_normal_lm
primary_variable = disp_x
[]
[normal_x]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_x
component = x
use_displaced_mesh = true
compute_lm_residuals = false
[]
[normal_y]
type = NormalMortarMechanicalContact
primary_boundary = plank_right
secondary_boundary = block_left
primary_subdomain = frictionless_primary_subdomain
secondary_subdomain = frictionless_secondary_subdomain
variable = frictionless_normal_lm
secondary_variable = disp_y
component = y
use_displaced_mesh = true
compute_lm_residuals = false
[]
[thermal_contact]
type = GapConductanceConstraint
variable = thermal_lm
secondary_variable = temp
k = 1
use_displaced_mesh = true
primary_boundary = plank_right
primary_subdomain = frictionless_primary_subdomain
secondary_boundary = block_left
secondary_subdomain = frictionless_secondary_subdomain
displacements = 'disp_x disp_y'
[]
[]
[BCs]
[left_temp]
type = DirichletBC
variable = temp
boundary = 'plank_left'
value = 400
[]
[right_temp]
type = DirichletBC
variable = temp
boundary = 'block_right'
value = 300
[]
[left_x]
type = DirichletBC
variable = disp_x
boundary = plank_left
value = 0.0
[]
[left_y]
type = DirichletBC
variable = disp_y
boundary = plank_bottom
value = 0.0
[]
[right_x]
type = ADFunctionDirichletBC
variable = disp_x
boundary = block_right
function = '-0.04*sin(4*(t+1.5))+0.02'
preset = false
[]
[right_y]
type = ADFunctionDirichletBC
variable = disp_y
boundary = block_right
function = '-t'
preset = false
[]
[]
[Materials]
[plank]
type = ADComputeIsotropicElasticityTensor
block = 'plank'
poissons_ratio = 0.3
youngs_modulus = ${E_plank}
[]
[block]
type = ADComputeIsotropicElasticityTensor
block = 'block'
poissons_ratio = 0.3
youngs_modulus = ${E_block}
[]
[stress]
type = ADComputeLinearElasticStress
block = 'plank block'
[]
[heat_plank]
type = ADHeatConductionMaterial
block = plank
thermal_conductivity = 2
specific_heat = 1
[]
[heat_block]
type = ADHeatConductionMaterial
block = block
thermal_conductivity = 1
specific_heat = 1
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options = '-snes_converged_reason -ksp_converged_reason'
petsc_options_iname = '-pc_type -mat_mffd_err -pc_factor_shift_type -pc_factor_shift_amount -snes_max_it'
petsc_options_value = 'lu 1e-5 NONZERO 1e-15 20'
end_time = 13.5
dt = 0.1
dtmin = 0.1
timestep_tolerance = 1e-6
line_search = 'none'
[]
[Postprocessors]
[nl_its]
type = NumNonlinearIterations
[]
[total_nl_its]
type = CumulativeValuePostprocessor
postprocessor = nl_its
[]
[l_its]
type = NumLinearIterations
[]
[total_l_its]
type = CumulativeValuePostprocessor
postprocessor = l_its
[]
[contact]
type = ContactDOFSetSize
variable = frictionless_normal_lm
subdomain = frictionless_secondary_subdomain
[]
[avg_hydro]
type = ElementAverageValue
variable = hydrostatic_stress
block = 'block'
[]
[avg_temp]
type = ElementAverageValue
variable = temp
block = 'block'
[]
[max_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
[]
[min_hydro]
type = ElementExtremeValue
variable = hydrostatic_stress
block = 'block'
value_type = min
[]
[avg_vonmises]
type = ElementAverageValue
variable = vonmises_stress
block = 'block'
[]
[max_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
[]
[min_vonmises]
type = ElementExtremeValue
variable = vonmises_stress
block = 'block'
value_type = min
[]
[]
[Outputs]
exodus = true
file_base = ${name}
checkpoint = true
[comp]
type = CSV
show = 'contact avg_temp'
[]
[out]
type = CSV
file_base = '${name}_out'
[]
[]
[Debug]
show_var_residual_norms = true
[]
(modules/contact/test/tests/bouncing-block-contact/frictionless-weighted-gap-mixed-basis.i)
starting_point = 2e-1
offset = 1e-2
[GlobalParams]
displacements = 'disp_x disp_y'
diffusivity = 1e0
scaling = 1e0
[]
[Mesh]
file = long-bottom-block-1elem-blocks.e
second_order = true
patch_update_strategy = always
[]
[Variables]
[./disp_x]
block = '1 2'
order = SECOND
[../]
[./disp_y]
block = '1 2'
order = SECOND
[../]
[./normal_lm]
block = 3
[../]
[]
[ICs]
[./disp_y]
block = 2
variable = disp_y
value = ${fparse starting_point + offset}
type = ConstantIC
[../]
[]
[Kernels]
[./disp_x]
type = MatDiffusion
variable = disp_x
[../]
[./disp_y]
type = MatDiffusion
variable = disp_y
[../]
[]
[Constraints]
[./weighted_gap_lm]
type = ComputeWeightedGapLMMechanicalContact
primary_boundary = 20
secondary_boundary = 10
primary_subdomain = 4
secondary_subdomain = 3
variable = normal_lm
disp_x = disp_x
disp_y = disp_y
use_displaced_mesh = true
[../]
[./ncp_lm]
type = ApplyPenetrationConstraintLMMechanicalContact
secondary = 10
primary = 20
variable = normal_lm
primary_variable = disp_x
c = 1
[../]
[normal_x]
type = NormalMortarMechanicalContact
primary_boundary = 20
secondary_boundary = 10
primary_subdomain = 4
secondary_subdomain = 3
variable = normal_lm
secondary_variable = disp_x
component = x
use_displaced_mesh = true
compute_lm_residuals = false
[]
[normal_y]
type = NormalMortarMechanicalContact
primary_boundary = 20
secondary_boundary = 10
primary_subdomain = 4
secondary_subdomain = 3
variable = normal_lm
secondary_variable = disp_y
component = y
use_displaced_mesh = true
compute_lm_residuals = false
[]
[]
[BCs]
[./botx]
type = DirichletBC
variable = disp_x
boundary = 40
value = 0.0
preset = false
[../]
[./boty]
type = DirichletBC
variable = disp_y
boundary = 40
value = 0.0
preset = false
[../]
[./topy]
type = FunctionDirichletBC
variable = disp_y
boundary = 30
function = '${starting_point} * cos(2 * pi / 40 * t) + ${offset}'
preset = false
[../]
[./leftx]
type = FunctionDirichletBC
variable = disp_x
boundary = 50
function = '1e-2 * t'
preset = false
[../]
[]
[Executioner]
type = Transient
end_time = 200
dt = 5
dtmin = .1
solve_type = 'PJFNK'
petsc_options = '-snes_converged_reason -ksp_converged_reason -pc_svd_monitor -snes_linesearch_monitor'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_shift_amount -mat_mffd_err'
petsc_options_value = 'lu NONZERO 1e-15 1e-5'
l_max_its = 30
nl_max_its = 20
line_search = 'none'
snesmf_reuse_base = false
abort_on_solve_fail = true
[]
[Debug]
show_var_residual_norms = true
[]
[Outputs]
exodus = true
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Postprocessors]
active = 'num_nl cumulative contact'
[./num_nl]
type = NumNonlinearIterations
[../]
[./cumulative]
type = CumulativeValuePostprocessor
postprocessor = num_nl
[../]
[contact]
type = ContactDOFSetSize
variable = normal_lm
subdomain = '3'
execute_on = 'nonlinear timestep_end'
[]
[]