Example 4: Frictional contact
This example demonstrates a frictional contact using the user-defined backbone I-soil material model in MASTODON. Two brick elements are stacked vertically and a thin layer is modeled to describe the frictional contact (see Figure 1). Top and bottom layers use linear elastic material model whereas thin layer uses I-soil to approximate a contact behavior with elastic - nearly perfectly plastic shear behavior.
This type of backbone curve can now be automated using the 'thin_layer' option in I-soil and ComputeISoilStress.
[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
type = FileMesh
file = '2blockfriction_isoilunit.e'
[]
[GlobalParams<<<{"href": "../syntax/GlobalParams/index.html"}>>>]
displacements = 'disp_x disp_y disp_z'
[]
[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
[./disp_x]
[../]
[./disp_y]
[../]
[./disp_z]
[../]
[]
[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
[./stress_xy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./stress_xx]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./stress_yy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./stress_zz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./stress_yz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./stress_zx]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./strain_yz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./strain_zx]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./layer_id]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[../]
[./vel_x]
[../]
[./accel_x]
[../]
[./vel_y]
[../]
[./accel_y]
[../]
[./vel_z]
[../]
[./accel_z]
[../]
[]
[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
[./TensorMechanics<<<{"href": "../syntax/Kernels/TensorMechanics/index.html"}>>>]
displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
[../]
[./gravity]
type = Gravity<<<{"description": "Apply gravity. Value is in units of acceleration.", "href": "../source/kernels/Gravity.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
value<<<{"description": "Value multiplied against the residual, e.g. gravitational acceleration"}>>> = -9.81
[../]
[./inertia_x]
type = InertialForce<<<{"description": "Calculates the residual for the inertial force ($M \\cdot acceleration$) and the contribution of mass dependent Rayleigh damping and HHT time integration scheme ($\\eta \\cdot M \\cdot ((1+\\alpha)velq2-\\alpha \\cdot vel-old) $)", "href": "../source/kernels/InertialForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
velocity<<<{"description": "velocity variable"}>>> = vel_x
acceleration<<<{"description": "acceleration variable"}>>> = accel_x
beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
eta<<<{"description": "Name of material property or a constant real number defining the eta parameter for the Rayleigh damping."}>>>=0.0
[../]
[./inertia_y]
type = InertialForce<<<{"description": "Calculates the residual for the inertial force ($M \\cdot acceleration$) and the contribution of mass dependent Rayleigh damping and HHT time integration scheme ($\\eta \\cdot M \\cdot ((1+\\alpha)velq2-\\alpha \\cdot vel-old) $)", "href": "../source/kernels/InertialForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
velocity<<<{"description": "velocity variable"}>>> = vel_y
acceleration<<<{"description": "acceleration variable"}>>> = accel_y
beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
eta<<<{"description": "Name of material property or a constant real number defining the eta parameter for the Rayleigh damping."}>>>=0.0
[../]
[./inertia_z]
type = InertialForce<<<{"description": "Calculates the residual for the inertial force ($M \\cdot acceleration$) and the contribution of mass dependent Rayleigh damping and HHT time integration scheme ($\\eta \\cdot M \\cdot ((1+\\alpha)velq2-\\alpha \\cdot vel-old) $)", "href": "../source/kernels/InertialForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
velocity<<<{"description": "velocity variable"}>>> = vel_z
acceleration<<<{"description": "acceleration variable"}>>> = accel_z
beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.25
gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.5
eta<<<{"description": "Name of material property or a constant real number defining the eta parameter for the Rayleigh damping."}>>> = 0.0
[../]
[]
[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
[./stress_zx]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_zx
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 2
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
[../]
[./strain_zx]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = total_strain
variable<<<{"description": "The name of the variable that this object applies to"}>>> = strain_zx
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 2
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
[../]
[./stress_xx]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xx
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
[../]
[./stress_yy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
[../]
[./stress_zz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_zz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 2
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[../]
[./stress_xy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
[../]
[./stress_yz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
[../]
[./layers]
type = UniformLayerAuxKernel<<<{"description": "Computes an AuxVariable for representing a layered structure in an arbitrary direction.", "href": "../source/auxkernels/UniformLayerAuxKernel.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = layer_id
interfaces<<<{"description": "A list of layer interface locations to apply across the domain in the specified direction."}>>> = '11.0'
direction<<<{"description": "The direction to apply layering."}>>> = '0.0 0.0 1.0'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1002
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = initial
[../]
[]
[BCs<<<{"href": "../syntax/BCs/index.html"}>>>]
[./bottom_x]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 109
value<<<{"description": "Value of the BC"}>>> = 0.0
[../]
[./bottom_y]
type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../source/bcs/DirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 109
value<<<{"description": "Value of the BC"}>>> = 0.0
[../]
[./bottom_z]
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"}>>> = 109
value<<<{"description": "Value of the BC"}>>> = 0.0
[../]
[./Periodic<<<{"href": "../syntax/BCs/Periodic/index.html"}>>>]
[./x_dir1]
variable<<<{"description": "Variable for the periodic boundary condition"}>>> = 'disp_x disp_y disp_z'
primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = '101'
secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = '103'
translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '10.0 0.0 0.0'
[../]
[./y_dir1]
variable<<<{"description": "Variable for the periodic boundary condition"}>>> = 'disp_x disp_y disp_z'
primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = '102'
secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = '104'
translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '0.0 10.0 0.0'
[../]
[./x_dir2]
variable<<<{"description": "Variable for the periodic boundary condition"}>>> = 'disp_x disp_y disp_z'
primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = '105'
secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = '107'
translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '10.0 0.0 0.0'
[../]
[./y_dir2]
variable<<<{"description": "Variable for the periodic boundary condition"}>>> = 'disp_x disp_y disp_z'
primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = '106'
secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = '108'
translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '0.0 10.0 0.0'
[../]
[../]
[./X-dir]
type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 110
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
function<<<{"description": "The forcing function."}>>> = loading_bc
[../]
[]
[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
[./loading_bc]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = 'rampload8_unit.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
[../]
[./initial_zztop]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
value<<<{"description": "The user defined function."}>>> = '8792.0 * -9.81 * (21.0 - z)'
[../]
[./initial_ytop]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
value<<<{"description": "The user defined function."}>>> = '8792.0 * -9.81 * (21.0 - z) * 0.2/0.8'
[../]
[./initial_zzmid]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
value<<<{"description": "The user defined function."}>>> = '-862495.2 - 2500.0 * -9.81 * (11.0 - z)'
[../]
[./initial_ymid]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
value<<<{"description": "The user defined function."}>>> = '(-862495.2 - 2500.0 * -9.81 * (11.0 - z)) * 0.45/0.55'
[../]
[./initial_zzbot]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
value<<<{"description": "The user defined function."}>>> = '-887020.2 - 2500.0 * -9.81 * (10.0 - z)'
[../]
[./initial_ybot]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
value<<<{"description": "The user defined function."}>>> = '(-887020.2 - 2500.0 * -9.81 * (10.0 - z)) * 0.3/0.7'
[../]
[]
[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
[./elasticity_tensor_top]
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 6.5e10 #Pa
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.2
type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1003
[../]
[./strain_top]
#Computes the strain, assuming small strains
type = ComputeIncrementalSmallStrain<<<{"description": "Compute a strain increment and rotation increment for small strains.", "href": "../source/materials/ComputeIncrementalStrain.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1003
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = ini_stress
[../]
[./strain_from_initial_stress]
type = ComputeEigenstrainFromInitialStress<<<{"description": "Computes an eigenstrain from an initial stress", "href": "../source/materials/ComputeEigenstrainFromInitialStress.html"}>>>
initial_stress<<<{"description": "A list of functions describing the initial stress. There must be 9 of these, corresponding to the xx, yx, zx, xy, yy, zy, xz, yz, zz components respectively. To compute the eigenstrain correctly, your elasticity tensor should not be time-varying in the first timestep"}>>> = '0 0 0 0 0 0 0 0 initial_zztop'
eigenstrain_name<<<{"description": "Material property name for the eigenstrain tensor computed by this model. IMPORTANT: The name of this property must also be provided to the strain calculator."}>>> = ini_stress
[../]
[./stress_top]
#Computes the stress, using linear elasticity
type = ComputeFiniteStrainElasticStress<<<{"description": "Compute stress using elasticity for finite strains", "href": "../source/materials/ComputeFiniteStrainElasticStress.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1003
[../]
[./den_top]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../source/materials/GenericConstantMaterial.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1003
prop_names<<<{"description": "The names of the properties this material will have"}>>> = density
prop_values<<<{"description": "The values associated with the named properties"}>>> = 8792 #kg/m3
[../]
[./elasticity_tensor_bot]
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 7e9#Pa
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.3
type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1001
[../]
[./strain_bot]
#Computes the strain, assuming small strains
type = ComputeIncrementalSmallStrain<<<{"description": "Compute a strain increment and rotation increment for small strains.", "href": "../source/materials/ComputeIncrementalStrain.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1001
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
eigenstrain_names<<<{"description": "List of eigenstrains to be applied in this strain calculation"}>>> = 'ini_stress_bot'
[../]
[./stress_bot]
#Computes the stress, using linear elasticity
type = ComputeFiniteStrainElasticStress<<<{"description": "Compute stress using elasticity for finite strains", "href": "../source/materials/ComputeFiniteStrainElasticStress.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1001
[../]
[./strain_from_initial_stress_bot]
type = ComputeEigenstrainFromInitialStress<<<{"description": "Computes an eigenstrain from an initial stress", "href": "../source/materials/ComputeEigenstrainFromInitialStress.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1001
initial_stress<<<{"description": "A list of functions describing the initial stress. There must be 9 of these, corresponding to the xx, yx, zx, xy, yy, zy, xz, yz, zz components respectively. To compute the eigenstrain correctly, your elasticity tensor should not be time-varying in the first timestep"}>>> = 'initial_ybot 0 0 0 initial_ybot 0 0 0 initial_zzbot'
eigenstrain_name<<<{"description": "Material property name for the eigenstrain tensor computed by this model. IMPORTANT: The name of this property must also be provided to the strain calculator."}>>> = ini_stress_bot
[../]
[./den_bot]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../source/materials/GenericConstantMaterial.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1001
prop_names<<<{"description": "The names of the properties this material will have"}>>> = density
prop_values<<<{"description": "The values associated with the named properties"}>>> = 2500 #kg/m3
[../]
[./I_Soil<<<{"href": "../syntax/Materials/I_Soil/index.html"}>>>]
[./soil1]
layer_variable<<<{"description": "The auxvariable providing the soil layer identification."}>>> = layer_id
layer_ids<<<{"description": "Vector of layer ids that map one-to-one to the rest of the soil layer parameters provided as input."}>>> = '0'
soil_type<<<{"description": "This parameter determines the type of backbone curve used. Use 'user_defined' for a user defined backbone curve provided in a data file, 'darendeli' if the backbone curve is to be determined using Darandeli equations, 'gqh' if the backbone curve is determined using the GQ/H approach and 'thin_layer' if the soil is being used to simulate a thin-layer friction interface."}>>> = 'user_defined'
backbone_curve_files<<<{"description": "The vector of file names of the files containing stress-strain backbone curves for the different soil layers. The size of the vector should be same as the size of layer_ids. All files should contain the same number of stress-strain points. Headers are not expected and it is assumed that the first column corresponds to strain values and the second column corresponds to the stress values. Additionally, two segments of a backbone curve cannot have the same slope."}>>> = 'backbone_curveunit.csv'
poissons_ratio<<<{"description": "Poissons's ratio for the soil layers. The size of the vector should be same as the size of layer_ids."}>>> = '0.45'
block<<<{"description": "The blocks where this material model is applied."}>>> = 1002
initial_soil_stress<<<{"description": "The function values for the initial stress distribution. 9 function names have to be provided corresponding to stress_xx, stress_xy, stress_xz, stress_yx, stress_yy, stress_yz, stress_zx, stress_zy, stress_zz. Each function can be a function of space."}>>> = 'initial_ymid 0 0 0 initial_ymid 0 0 0 initial_zzmid'
density<<<{"description": "Vector of density values that map one-to-one with the number 'layer_ids' parameter."}>>> = '2500'
#initial_bulk_modulus = '7.83e10'
initial_shear_modulus<<<{"description": "The initial shear modulus of the soil layers."}>>> = '2.70588e11'
#taumax = 275998
# theta_1 = '0'
# theta_2 = '0'
# theta_3 = '0'
# theta_4 = '0'
# theta_5 = '0'
pressure_dependency<<<{"description": "Set to true to turn on pressure dependent stiffness and yield strength calculation."}>>> = true
#b_exp = 0.0
p_ref<<<{"description": "The reference pressure at which the parameters are defined for each soil layer. If 'soil_type = darendeli', then the reference pressure must be input in kilopascals."}>>> = 880500
#tension_pressure_cut_off = -1
a0<<<{"description": "The first coefficient for pressure dependent yield strength calculation for all the soil layers. If a0 = 1, a1 = 0 and a2=0 for one soil layer, then the yield strength of that layer is independent of pressure."}>>> = 0
a1<<<{"description": "The second coefficient for pressure dependent yield strength calculation for all the soil layers. If a0 = 1, a1 = 0, a2 = 0 for one soil layer, then the yield strength of that layer is independent of pressure."}>>> = 0
a2<<<{"description": "The third coefficient for pressure dependent yield strength calculation for all the soil layers. If a0 = 1, a1=0 and a2=0 for one soil layer, then the yield strength of that layer is independent of pressure."}>>> = 1
[../]
[../]
[]
[Preconditioning<<<{"href": "../syntax/Preconditioning/index.html"}>>>]
[./andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
[../]
[]
[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = PJFNK
nl_abs_tol = 1e-7
nl_rel_tol = 1e-7
l_tol = 1e-6
l_max_its = 20
start_time = 0
end_time = 1.2
dt = 0.01
timestep_tolerance = 1e-6
petsc_options = '-snes_ksp_ew'
petsc_options_iname = '-ksp_gmres_restart -pc_type -pc_hypre_type -pc_hypre_boomeramg_max_iter'
petsc_options_value = '201 hypre boomeramg 4'
line_search = 'none'
[]
[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
[./_dt]
type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../source/postprocessors/TimestepSize.html"}>>>
[../]
[./stres_xx_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = stress_xx
[../]
[./stres_xy_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = stress_xy
[../]
[./stres_yy_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = stress_yy
[../]
[./stres_yz_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = stress_yz
[../]
[./stres_zz_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = stress_zz
[../]
[./stress_zx_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = stress_zx
[../]
[./strain_zx_interface]
type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = '1002'
variable<<<{"description": "The name of the variable that this object operates on"}>>> = strain_zx
[../]
[]
[Outputs<<<{"href": "../syntax/Mastodon/Outputs/index.html"}>>>]
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = true
print_linear_residuals<<<{"description": "Enable printing of linear residuals to the screen (Console)"}>>> = false
[./screen]
type = Console<<<{"description": "Object for screen output.", "href": "../source/outputs/Console.html"}>>>
max_rows<<<{"description": "The maximum number of postprocessor/scalar values displayed on screen during a timestep (set to 0 for unlimited)"}>>> = 1
[../]
[]
(examples/ex04/2BlockFriction_Isoilunit.i)
Figure 1: Undeformed shape of the numerical model.
The domain is loaded with gravity and top element is sheared on the top surface. At intermediate steps whereby the contact strength is yet not mobilized, no separation occurs between two elastic blocks as it is shown in Figure 2.

Figure 2: Intermediate step whereby the contact layer strength is not mobilized.
Once the contact layer is yielded, the upper block starts to slide on top of the lower block and thin layer effectively simulates the interface between two elastic layers (see Figure 3).

Figure 3: Deformed shape of the numerical model after the contact is yielded.