Example 1b: Direct shear test using the Darandeli backbone curve formulation (explicit time integration)

Example 1a can also be run using an explicit time integration through the CentralDifference scheme. Unlike implicit, Newmark-Beta integration Central Difference integration does not require nonlinear iterations (i.e., calculating the Jacobian, etc.). However, explicit schemes are not unconditionally stable. For stability, the timestep should be less than the so-called critical time step. CentralDifference and CriticalTimeStep pages provide more information on explicit integration and its stability, respectively.

When the using the explicit scheme, while the user can specify the time step for the solver, it is recommended to use the CriticalTimeStep postprocessor which automatically calculates and applies the critical time step based on the mesh size and material properties. A factor smaller than unity is recommended to be multiplied to the critical time step to ensure stability.

The input file is presented below:

# One element test to test the auto-generated Darendeli backbone curve.
# The back surface of the element (z=0) is fixed and the front surface (z=1)
# is moved by applying a cyclic preset displacement.

# The resulting shear stress vs shear strain curve was verified against that obtained
# from DEEPSOIL.

[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh # Can generate simple lines, rectangles and rectangular prisms
  dim = 3 # Dimension of the mesh
  nx = 1 # Number of elements in the x direction
  ny = 1 # Number of elements in the y direction
  nz = 1 # Number of elements in the z direction
  xmin = 0.0
  xmax = 1
  ymin = 0.0
  ymax = 1
  zmin = 0.0
  zmax = 1
[]

[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"}>>>]
  [./vel_x]
  [../]
  [./accel_x]
  [../]
  [./vel_y]
  [../]
  [./accel_y]
  [../]
  [./vel_z]
  [../]
  [./accel_z]
  [../]
  [./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_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_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
  [../]
  [./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
  [../]
  [./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
  [../]
  [./strain_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
  [../]
  [./strain_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
  [../]
  [./strain_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
  [../]
  [./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
  [../]
[]

[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
  [./DynamicTensorMechanics<<<{"href": "../syntax/Kernels/DynamicTensorMechanics/index.html"}>>>]
    displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
    # stiffness_damping_coefficient = 0.00006366
  [../]
  [./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
    eta<<<{"description": "Name of material property or a constant real number defining the eta parameter for the Rayleigh damping."}>>> = 7.854
  [../]
  [./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
    eta<<<{"description": "Name of material property or a constant real number defining the eta parameter for the Rayleigh damping."}>>> = 7.854
  [../]
  [./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
    eta<<<{"description": "Name of material property or a constant real number defining the eta parameter for the Rayleigh damping."}>>> = 7.854
  [../]
  [./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
  [../]
[]

[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
  [./accel_x]
    type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_x
    displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_x
    first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
  [../]
  [./vel_x]
    type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_x
    displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_x
  [../]
  [./accel_y]
    type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_y
    displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_y
    first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
  [../]
  [./vel_y]
    type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_y
    displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_y
  [../]
  [./accel_z]
    type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_z
    displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_z
    first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
  [../]
  [./vel_z]
    type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_z
    displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_z
  [../]
  [./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)"}>>> = 1
    index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
  [../]
  [./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)"}>>> = 2
    index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
  [../]
  [./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)"}>>> = 0
    index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
  [../]
  [./strain_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"}>>> = total_strain
    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)"}>>> = 1
    index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
  [../]
  [./strain_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"}>>> = total_strain
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = strain_yz
    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)"}>>> = 1
  [../]
  [./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)"}>>> = 0
    index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
  [../]
  [./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
  [../]
  [./strain_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"}>>> = total_strain
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = strain_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
  [../]
  [./strain_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"}>>> =total_strain
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = strain_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
  [../]
  [./strain_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"}>>> = total_strain
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = strain_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
  [../]
  [./layer]
    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."}>>> = '2.0'
    direction<<<{"description": "The direction to apply layering."}>>> = '0 0 1'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = initial
  [../]
[]

[BCs<<<{"href": "../syntax/BCs/index.html"}>>>]
  [./x_bot]
    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"}>>> = 0
    value<<<{"description": "Value of the BC"}>>> = 0.0
  [../]
  [./y_bot]
    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"}>>> = 0
    value<<<{"description": "Value of the BC"}>>> = 0.0
  [../]
  [./z_bot]
    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"}>>> = 0
    value<<<{"description": "Value of the BC"}>>> = 0.0
  [../]
  [./Periodic<<<{"href": "../syntax/BCs/Periodic/index.html"}>>>]
    [./x_dir]
      variable<<<{"description": "Variable for the periodic boundary condition"}>>> = 'disp_x disp_y disp_z'
      primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = '4'
      secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = '2'
      translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '1.0 0.0 0.0'
    [../]
    [./y_dir]
      variable<<<{"description": "Variable for the periodic boundary condition"}>>> = 'disp_x disp_y disp_z'
      primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = '1'
      secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = '3'
      translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '0.0 1.0 0.0'
    [../]
  [../]
  [./top_x]
    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"}>>> = 5
    function<<<{"description": "The forcing function."}>>> = top_disp
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_x'
  [../]
[]

[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
  [./top_disp]
    type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
    data_file<<<{"description": "File holding CSV data"}>>> = '../ex01a/Displacement2.csv'
    format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
  [../]
[]

# [Materials]
#   [./sample_isoil]
#     type = ComputeISoilStress
#     soil_type = 'darendeli'
#     layer_variable = layer_id
#     layer_ids = '0'
#     initial_soil_stress = '-4.204286 0 0  0 -4.204286 0  0 0 -9.810'
#     poissons_ratio = '0.3'
#     number_of_points = 10
#     over_consolidation_ratio = '1'
#     plasticity_index = '0'
#     initial_shear_modulus = '20000'
#     p_ref = '6.07286'
#     block = '0'
#   [../]
#   [./sample_isoil_strain]
#     type = ComputeIncrementalSmallStrain
#     block = '0'
#     displacements = 'disp_x disp_y disp_z'
#     implicit = false
#   [../]
#   [./sample_isoil_elasticitytensor]
#     type = ComputeIsotropicElasticityTensorSoil
#     block = '0'
#     density = '2'
#     wave_speed_calculation = false
#     layer_ids = '0'
#     poissons_ratio = '0.3'
#     shear_modulus = '20000'
#     layer_variable = layer_id
#   [../]
# []

[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
  [./I_Soil<<<{"href": "../syntax/Materials/I_Soil/index.html"}>>>]
    [./soil_1]
      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."}>>> = 'darendeli'
      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'
      over_consolidation_ratio<<<{"description": "The over consolidation ratio of the soil layers. Required for Darandeli backbone curve."}>>> = '1'
      plasticity_index<<<{"description": "The plasticity index of the soil layers. Required for Darandeli backbone curve."}>>> = '0'
      initial_shear_modulus<<<{"description": "The initial shear modulus of the soil layers."}>>> = '20000'
      number_of_points<<<{"description": "The total number of data points in which the backbone curve needs to be split for all soil layers (required for Darandeli or GQ/H type backbone curves)."}>>> = 10
      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.3'
      block<<<{"description": "The blocks where this material model is applied."}>>> = 0
      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."}>>> = '-4.204286 0 0  0 -4.204286 0  0 0 -9.810'
      density<<<{"description": "Vector of density values that map one-to-one with the number 'layer_ids' parameter."}>>> = '2'
      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."}>>> = '6.07286'
    [../]
  [../]
[]

[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
  start_time = 0.0
  end_time = 8.0
  timestep_tolerance = 1e-6
  l_abs_tol = 1e-12
  [./TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
    type = CentralDifference
  [../]
  [./TimeStepper<<<{"href": "../syntax/Executioner/TimeStepper/index.html"}>>>]
    type = PostprocessorDT
    postprocessor = time_step
    dt = 0.001
  [../]
[]

[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
  [./_dt]
    type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../source/postprocessors/TimestepSize.html"}>>>
  [../]
  [./disp_6x]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = disp_x
  [../]
  [./disp_6y]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = disp_y
  [../]
  [./disp_6z]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = disp_z
  [../]
  [./vel_6x]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = vel_x
  [../]
  [./vel_6y]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = vel_y
  [../]
  [./vel_6z]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = vel_z
  [../]
  [./accel_6x]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = accel_x
  [../]
  [./accel_6y]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = accel_y
  [../]
  [./accel_6z]
    type = NodalVariableValue<<<{"description": "Outputs values of a nodal variable at a particular location", "href": "../source/postprocessors/NodalVariableValue.html"}>>>
    nodeid<<<{"description": "The ID of the node where we monitor"}>>> = 6
    variable<<<{"description": "The variable to be monitored"}>>> = accel_z
  [../]
  [./stress_xy_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = stress_xy
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./stress_yz_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = stress_yz
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./stress_zx_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = stress_zx
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./strain_xy_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = strain_xy
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./strain_yz_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = strain_yz
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./strain_zx_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = strain_zx
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./stress_xx_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = stress_xx
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./stress_yy_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = stress_yy
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./stress_zz_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = stress_zz
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./strain_xx_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = strain_xx
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./strain_yy_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = strain_yy
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./strain_zz_el]
    type = ElementalVariableValue<<<{"description": "Outputs an elemental variable value at a particular location", "href": "../source/postprocessors/ElementalVariableValue.html"}>>>
    variable<<<{"description": "The variable to be monitored"}>>> = strain_zz
    elementid<<<{"description": "The ID of the element where we monitor"}>>> = 0
  [../]
  [./time_step]
    type = CriticalTimeStep<<<{"description": "Computes and reports the critical time step for the explicit solver.", "href": "../source/postprocessors/CriticalTimeStep.html"}>>>
    factor<<<{"description": "Factor to mulitply to the critical time step."}>>> = 0.9
  [../]
[]

[Outputs<<<{"href": "../syntax/Mastodon/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = false
  file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = Test # HYS_Darendeli_explicit
[]
(examples/ex01b/HYS_darendeli_explicit.i)

The plot below compares the hysteretic response between explicit and implicit integration schemes. These results match closely.

Figure 1: Comparison of hysteretic responses generated using explicit and implicit schemes.

commentnote:Current limitations of explicit integration in MASTODON

Explicit integration using the CentralDifference scheme: (1) does not account for element stiffness damping; (2) cannot be used for element types other than solid elements; and (3) does not work with input acceleration (input displacement should be applied). However, the MASTODON team is working to remove these limitations in the near future.