Example 3b: Nonlinear site-response analysis (explicit time integration)

Example 3a 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. Also, displacement input has to be applied to the soil column instead of an acceleration input.

The input file is presented below:

[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  nx = 1
  ny = 1
  nz = 20
  xmin = -0.5
  ymin = -0.5
  zmin = 0.0
  xmax = 0.5
  ymax = 0.5
  zmax = 20.0
  dim = 3
[]

[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"}>>>]
  [./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
  [../]
  [./vel_x]
  [../]
  [./accel_x]
  [../]
  [./vel_y]
  [../]
  [./accel_y]
  [../]
  [./vel_z]
  [../]
  [./accel_z]
  [../]
  [./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.000781
  [../]
  [./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 = 0.64026
  [../]
  [./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 = 0.64026
  [../]
  [./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 = 0.64026
  [../]
  [./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"}>>>]
  [./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"}>>> = strain_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
  [../]
  [./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
  [../]
  [./layer_id]
     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."}>>> = '1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20.2'
     direction<<<{"description": "The direction to apply layering."}>>> = '0.0 0.0 1.0'
     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_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"}>>> = 'back'
    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"}>>> = 'back'
    value<<<{"description": "Value of the BC"}>>> = 0.0
  [../]
  [./bottom_shake]
    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"}>>> = 'back'
    function<<<{"description": "The forcing function."}>>> = disp_bottom
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'disp_x'
  [../]
  [./Periodic<<<{"href": "../syntax/BCs/Periodic/index.html"}>>>]
    [./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."}>>> = 'bottom'
      secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = 'top'
      translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '0.0 1.0 0.0'
    [../]
    [./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."}>>> = 'left'
      secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = 'right'
      translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '1.0 0.0 0.0'
    [../]
  [../]
[]

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

  [./initial_zz]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
    value<<<{"description": "The user defined function."}>>> = '-2000.0 * 9.81 * (20.0 - z)'
  [../]
  [./initial_xx]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../source/functions/MooseParsedFunction.html"}>>>
    value<<<{"description": "The user defined function."}>>> = '-2000.0 * 9.81 * (20.0 - z) * 0.3/0.7'
  [../]
[]

# [Materials]
#   [./sample_isoil]
#     type = ComputeISoilStress
#     soil_type = 'gqh'
#     layer_variable = layer_id
#     layer_ids = '0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19'
#     initial_soil_stress = 'initial_xx 0 0 0 initial_xx 0 0 0 initial_zz'
#     poissons_ratio = '0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3'
#     number_of_points = 100
#     ### GQ/H ####
#     initial_shear_modulus = '125000000 118098000 111392000 103968000 96800000 89888000 83232000 76832000 70688000 64800000 59168000 53792000 48672000 43808000 39200000 34848000 30752000 26912000 23328000 20000000'
#     theta_1 = '-6.66 -6.86 -7.06 -7.35 -7.65 -7.95 -8.3 -8.61 -8.95 -9.3 -9.61 -9.92 -10.0 -10.0 -10.0 -10.0 -10.0 -9.31 -7.17 -5.54'
#     theta_2 = '5.5 5.7 5.9 6.2 6.6 6.9 7.3 7.6 8.0 8.4 8.6 8.8 8.82 8.71 8.55 8.3 7.88 6.4 2.4 -2.28'
#     theta_3 = '1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0'
#     theta_4 = '1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0'
#     theta_5 = '0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99'
#     taumax = '292500 277500 262500 247500 232500 217500 202500 187500 172500 157500 142500 127500 112500 97500 82500 67500 52500 37500 22500 7500'
#     p_ref = '236841 224695 212550 200404 188258 176112 163967 151821 139675 127530 115384 103238 91092 78947 66801 54655 42510 30364 18218 6072'
#     ######
#     a0 = 1.0
#     a1 = 0.0
#     a2 = 0.0
#     b_exp = 0.0
#     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 = '2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0'
#     wave_speed_calculation = false
#     layer_ids = '0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19'
#     poissons_ratio = '0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3'
#     shear_modulus = '125000000 118098000 111392000 103968000 96800000 89888000 83232000 76832000 70688000 64800000 59168000 53792000 48672000 43808000 39200000 34848000 30752000 26912000 23328000 20000000'
#     layer_variable = layer_id
#   [../]
# []

[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
  [./I_Soil<<<{"href": "../syntax/Materials/I_Soil/index.html"}>>>]
    [./soil_all]
       block<<<{"description": "The blocks where this material model is applied."}>>> = 0
       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 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19'
       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_xx 0 0 0 initial_xx 0 0 0 initial_zz'
       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 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3 0.3'
       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."}>>> = 'gqh'
       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)."}>>> = 100
       ### GQ/H ####
       initial_shear_modulus<<<{"description": "The initial shear modulus of the soil layers."}>>> = '125000000 118098000 111392000 103968000 96800000 89888000 83232000 76832000 70688000 64800000 59168000 53792000 48672000 43808000 39200000 34848000 30752000 26912000 23328000 20000000'
       theta_1<<<{"description": "The curve fit coefficient for GQ/H modelfor each soil layer."}>>> = '-6.66 -6.86 -7.06 -7.35 -7.65 -7.95 -8.3 -8.61 -8.95 -9.3 -9.61 -9.92 -10.0 -10.0 -10.0 -10.0 -10.0 -9.31 -7.17 -5.54'
       theta_2<<<{"description": "The curve fit coefficient for GQ/H modelfor each soil layer."}>>> = '5.5 5.7 5.9 6.2 6.6 6.9 7.3 7.6 8.0 8.4 8.6 8.8 8.82 8.71 8.55 8.3 7.88 6.4 2.4 -2.28'
       theta_3<<<{"description": "The curve fit coefficient for GQ/H modelfor each soil layer."}>>> = '1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0'
       theta_4<<<{"description": "The curve fit coefficient for GQ/H modelfor each soil layer."}>>> = '1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0'
       theta_5<<<{"description": "The curve fit coefficient for GQ/H modelfor each soil layer."}>>> = '0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99 0.99'
       taumax<<<{"description": "The ultimate shear strength of the soil layers. Required for GQ/H model"}>>> = '292500 277500 262500 247500 232500 217500 202500 187500 172500 157500 142500 127500 112500 97500 82500 67500 52500 37500 22500 7500'
       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."}>>> = '236841 224695 212550 200404 188258 176112 163967 151821 139675 127530 115384 103238 91092 78947 66801 54655 42510 30364 18218 6072'
       ######
       density<<<{"description": "Vector of density values that map one-to-one with the number 'layer_ids' parameter."}>>> = '2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0 2000.0'
       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."}>>> = 1.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.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."}>>> = 0.0
       b_exp<<<{"description": "The exponential factors for pressure dependent stiffness for all the soil layers."}>>> = 0.0
    [../]
  [../]
[]

[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
  type = Transient
  start_time = 0.0
  end_time = 85.4
  timestep_tolerance = 1e-6
  l_abs_tol = 1e-12
  # dt = 0.001
  [./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"}>>>]
  [./accelx_top]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 20.0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_x
  [../]
  [./accelx_top1]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 19.0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_x
  [../]
  [./accelx_mid]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 10.0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_x
  [../]
  [./accelx_bot]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 0.0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_x
  [../]
  [./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
  [../]
[]

[VectorPostprocessors<<<{"href": "../syntax/VectorPostprocessors/index.html"}>>>]
  [./accel_hist]
    type = ResponseHistoryBuilder<<<{"description": "Calculates response histories for a given node and variable(s).", "href": "../source/vectorpostprocessors/ResponseHistoryBuilder.html"}>>>
    variables<<<{"description": "Variable name for which the response history is requested."}>>> = 'accel_x'
    nodes<<<{"description": "Node number(s) at which the response history is needed."}>>> = '80'
  [../]
  [./accel_spec]
    type = ResponseSpectraCalculator<<<{"description": "Calculate the response spectrum at the requested nodes or points.", "href": "../source/vectorpostprocessors/ResponseSpectraCalculator.html"}>>>
    vectorpostprocessor<<<{"description": "Name of the ResponseHistoryBuilder vectorpostprocessor, for which response spectra are calculated."}>>> = accel_hist
    regularize_dt<<<{"description": "dt for response spectra calculation. The acceleration response will be regularized to this dt prior to the response spectrum calculation."}>>> = 0.001
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out
  [../]
[]

[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)"}>>> = true
  [./out]
   type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../source/outputs/CSV.html"}>>>
   execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_begin'
   file_base<<<{"description": "The desired solution output name without an extension. If not provided, MOOSE sets it with Outputs/file_base when available. Otherwise, MOOSE uses input file name and this object name for a master input or uses master file_base, the subapp name and this object name for a subapp input to set it."}>>> = topsoil_response_explicit
  [../]
  [./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
   interval = 1000
  [../]
[]
(examples/ex03b/shear_beam_Isoil_free_field_explicit.i)

The plot below compares the accelerations recored at the top of the soil column between explicit and implicit integration schemes and also their 5 percent damped response spectra. These results match closely.

Figure 1: Comparison of accelerations at the soil surface generated using explicit and implicit schemes.

Figure 2: Comparison of response spectra at the soil surface 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.