Example 13: Fixed base analysis of nuclear power plant structure with fluid structure interaction
GigaPascals, meters, and seconds.
Model description
This example demonstrates a fixed base analysis of a representative nuclear power-plant structure complete with an internal reactor vessel and generic fluid material where fluid structure interaction (FSI) effects are modeled. Details regarding the model can be seen in Figure 1 below.

Figure 1: Top Left: View of inner structure. Bottom Left: View of reactor vessel and supporting components. Right: Full fixed base structural model.
The structure shown in the figure above sits on top of a 1m thick basemat. The building is surrounded by external buttresses that also span the top of the structure beneath the roof slab. The cylindrical reactor vessel has a spherical bottom and is suspended from a mid level slab, which is supported by internal walls as well as a cylindrical citadel. A more detailed look at the reactor vessel is shown in Figure 2 below.

Figure 2: Section view of the support citadel, reactor vessel and internal fluid.
The reactor vessel shown in Figure 2 is 2cm thick, the inner diameter of the vessel is 5 m and the vessel is approximately 6m tall from the bottom of the spherical end to the top of the vessel.
Finite element mesh generation
The three-dimensional finite element mesh is generated using Cubit 15.6. The procedure for generating this mesh is not trivial, so the journal script used to create the model is provided in the example 13 directory along with other necessary input files. In order to adequately capture bending of the structure walls, a minimum of 2 linear elements through the thickness of the walls is required. Alternatively, 1 or 2 second order elements could be used, but this will result in an increase in computational costs. The model mesh can be seen in Figure 1.
Modeling in MASTODON
This example consists of two individual analysis efforts for: (i) site response analysis (no structure); (ii) fixed-base structure (no soil)
Site response analysis
A single stack of 3D volume elements are generated within the mastodon input file in order to define the site response model domain. Nonlinear soil is modeled using the Isoil material block, where the backbone curve is supplied as a csv file. Periodic boundary conditions are applied to the two opposite faces of the soil domain, which are perpendicular to the Z-axis to ensure that they move together. This is achieved using the Periodic
option in the boundary conditions. The translational motion is fully defined in X Y and Z directions as shown in Figure 3 and is applied at the bottom of the soil block. The acceleration response is measured at the top using ResponseHistoryBuilder
. This represents the free field soil response and is used as the input motion for the fixed base model. For more information on how to set up a site response model in MASTODON, see Example 3a

Figure 3: Ormsby wavelet input ground motion used for the analysis.
The input file for simulating the free field soil response is presented below.
[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
nx = 1
ny = 1
nz = 17
xmin = -0.05
ymin = -0.05
zmin = 0.0
xmax = 0.05
ymax = 0.05
zmax = 32.0
dim = 3
[]
[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
[disp_x]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
[]
[disp_y]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
[]
[disp_z]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
[]
[]
[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
[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"}>>>]
stiffness_damping_coefficient<<<{"description": "Name of material property or a constant real number defining stiffness Rayleigh parameter (zeta)."}>>> = 0.0017490 #1.5% 1-20hz
displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
static_initialization<<<{"description": "Set to true get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step."}>>> = true
[]
[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.0607069 #1.5% 1-20hz
[]
[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.0607069 #1.5% 1-20hz
[]
[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.0607069 #1.5% 1-20hz
[]
[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 = NewmarkAccelAux<<<{"description": "Computes the current acceleration using the Newmark method.", "href": "../source/auxkernels/NewmarkAccelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_x
displacement<<<{"description": "displacement variable"}>>> = disp_x
velocity<<<{"description": "velocity variable"}>>> = vel_x
beta<<<{"description": "beta parameter for Newmark method"}>>> = 0.25
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_end
[]
[vel_x]
type = NewmarkVelAux<<<{"description": "Calculates the current velocity using Newmark method.", "href": "../source/auxkernels/NewmarkVelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_x
acceleration<<<{"description": "acceleration variable"}>>> = accel_x
gamma<<<{"description": "gamma parameter for Newmark method"}>>> = 0.5
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_end
[]
[accel_y]
type = NewmarkAccelAux<<<{"description": "Computes the current acceleration using the Newmark method.", "href": "../source/auxkernels/NewmarkAccelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_y
displacement<<<{"description": "displacement variable"}>>> = disp_y
velocity<<<{"description": "velocity variable"}>>> = vel_y
beta<<<{"description": "beta parameter for Newmark method"}>>> = 0.25
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_end
[]
[vel_y]
type = NewmarkVelAux<<<{"description": "Calculates the current velocity using Newmark method.", "href": "../source/auxkernels/NewmarkVelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_y
acceleration<<<{"description": "acceleration variable"}>>> = accel_y
gamma<<<{"description": "gamma parameter for Newmark method"}>>> = 0.5
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_end
[]
[accel_z]
type = NewmarkAccelAux<<<{"description": "Computes the current acceleration using the Newmark method.", "href": "../source/auxkernels/NewmarkAccelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_z
displacement<<<{"description": "displacement variable"}>>> = disp_z
velocity<<<{"description": "velocity variable"}>>> = vel_z
beta<<<{"description": "beta parameter for Newmark method"}>>> = 0.25
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_end
[]
[vel_z]
type = NewmarkVelAux<<<{"description": "Calculates the current velocity using Newmark method.", "href": "../source/auxkernels/NewmarkVelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_z
acceleration<<<{"description": "acceleration variable"}>>> = accel_z
gamma<<<{"description": "gamma parameter for Newmark method"}>>> = 0.5
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_end
[]
[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."}>>> = '32'
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
[]
[]
[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'
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.csv'
displacements<<<{"description": "The vector of displacement variables in the problem."}>>> = 'disp_x disp_y disp_z'
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'
initial_shear_modulus<<<{"description": "The initial shear modulus of the soil layers."}>>> = '0.523325'
density<<<{"description": "Vector of density values that map one-to-one with the number 'layer_ids' parameter."}>>> = '1.730e-6'
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
[]
[]
[]
[Controls<<<{"href": "../syntax/Controls/index.html"}>>>]
[inertia_switch]
type = TimePeriod<<<{"description": "Control the enabled/disabled state of objects with time.", "href": "../source/controls/TimePeriod.html"}>>>
start_time<<<{"description": "The time at which the objects are to be enabled/disabled."}>>> = 0.0
end_time<<<{"description": "The time at which the objects are to be enable/disabled."}>>> = 0.01
disable_objects<<<{"description": "A list of object tags to disable."}>>> = '*/inertia_x */inertia_y */inertia_z */vel_x */vel_y */vel_z */accel_x */accel_y */accel_z'
set_sync_times<<<{"description": "Set the start and end time as execute sync times."}>>> = true
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 timestep_end'
[]
[]
[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
[input_motion_x]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = '../ground_motion.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 1
[]
[input_motion_y]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = '../ground_motion.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 2
[]
[input_motion_z]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = '../ground_motion.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 3
[]
[]
[BCs<<<{"href": "../syntax/BCs/index.html"}>>>]
[x_motion]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
acceleration<<<{"description": "The acceleration variable."}>>> = accel_x
velocity<<<{"description": "The velocity variable."}>>> = vel_x
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_x'
[]
[y_motion]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
acceleration<<<{"description": "The acceleration variable."}>>> = accel_y
velocity<<<{"description": "The velocity variable."}>>> = vel_y
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_y'
[]
[z_motion]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
acceleration<<<{"description": "The acceleration variable."}>>> = accel_z
velocity<<<{"description": "The velocity variable."}>>> = vel_z
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_z'
[]
[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 0.1 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."}>>> = '0.1 0.0 0.0'
[]
[]
[]
[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 accel_y accel_z'
nodes<<<{"description": "Node number(s) at which the response history is needed."}>>> = '0 71'
[]
[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.002
damping_ratio<<<{"description": "Damping ratio for response spectra calculation."}>>> = 0.05
start_frequency<<<{"description": "Start frequency for the response spectra calculation."}>>> = 0.1
end_frequency<<<{"description": "End frequency for the response spectra calculation."}>>> = 1000
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[]
[Preconditioning<<<{"href": "../syntax/Preconditioning/index.html"}>>>]
[smp]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
[]
[]
[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
type = Transient
petsc_options = '-ksp_snes_ew'
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
solve_type = 'NEWTON'
nl_rel_tol = 1e-7
nl_abs_tol = 1e-12
# l_tol = 1e-2
dt = 0.01
end_time = 28.0
timestep_tolerance = 1e-6
[TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
type = NewmarkBeta
beta = 0.25
gamma = 0.5
[]
[]
[Outputs<<<{"href": "../syntax/Mastodon/Outputs/index.html"}>>>]
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
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[out1]
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."}>>> = 'final'
[]
[]
(examples/ex13/site_response/sr_input.i)Fixed-base structure analysis
Once the free field soil response has been collected from the site response model, that acceleration history can be used as an acceleration input for the fixed base structure model. The bottom surface is fixed to the accelerations collected from the site response model using the PresetAcceleration
boundary condition.

Figure 4: Meshed fixed base model used in this analysis.
[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
[mesh_gen]
type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../source/meshgenerators/FileMeshGenerator.html"}>>>
file<<<{"description": "The filename to read."}>>> = fixed_base_mesh.e
[]
[interface1]
type = SideSetsBetweenSubdomainsGenerator<<<{"description": "MeshGenerator that creates a sideset composed of the nodes located between two or more subdomains.", "href": "../source/meshgenerators/SideSetsBetweenSubdomainsGenerator.html"}>>>
input<<<{"description": "The mesh we want to modify"}>>> = mesh_gen
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = 'fluid_material'
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = 'RV'
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = 'Interface'
[]
displacements<<<{"description": "The variables corresponding to the x y z displacements of the mesh. If this is provided then the displacements will be taken into account during the computation. Creation of the displaced mesh can be suppressed even if this is set by setting 'use_displaced_mesh = false'."}>>> = 'disp_x disp_y disp_z'
[]
[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
[disp_x]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' #' # 99 98 97 96 95 94'
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
[]
[disp_y]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
[]
[disp_z]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
[]
[p]
block = 'fluid_material'
[]
[]
[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
[Wave1]
block = 'fluid_material'
[]
[vel_x]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[accel_x]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[vel_y]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[accel_y]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[vel_z]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[accel_z]
block = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[]
[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
[diffusion]
type = Diffusion<<<{"description": "The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak form of $(\\nabla \\phi_i, \\nabla u_h)$.", "href": "../source/kernels/Diffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = 'p'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'fluid_material'
use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
[]
[inertia]
type = AcousticInertia<<<{"description": "Calculates the residual for the inertial force which is the double time derivative of pressure.", "href": "../source/kernels/AcousticInertia.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'fluid_material'
use_displaced_mesh<<<{"description": "Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used."}>>> = false
[]
[DynamicTensorMechanics<<<{"href": "../syntax/Kernels/DynamicTensorMechanics/index.html"}>>>]
stiffness_damping_coefficient<<<{"description": "Name of material property or a constant real number defining stiffness Rayleigh parameter (zeta)."}>>> = 0.00006366
displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence kernels will be applied to"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
[]
[inertia_x]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
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]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
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]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
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
[]
[]
[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
[accel_x]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
type = NewmarkAccelAux<<<{"description": "Computes the current acceleration using the Newmark method.", "href": "../source/auxkernels/NewmarkAccelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_x
displacement<<<{"description": "displacement variable"}>>> = disp_x
velocity<<<{"description": "velocity variable"}>>> = vel_x
beta<<<{"description": "beta parameter for Newmark method"}>>> = 0.25
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_end
[]
[vel_x]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
type = NewmarkVelAux<<<{"description": "Calculates the current velocity using Newmark method.", "href": "../source/auxkernels/NewmarkVelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_x
acceleration<<<{"description": "acceleration variable"}>>> = accel_x
gamma<<<{"description": "gamma parameter for Newmark method"}>>> = 0.5
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_end
[]
[accel_y]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
type = NewmarkAccelAux<<<{"description": "Computes the current acceleration using the Newmark method.", "href": "../source/auxkernels/NewmarkAccelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_y
displacement<<<{"description": "displacement variable"}>>> = disp_y
velocity<<<{"description": "velocity variable"}>>> = vel_y
beta<<<{"description": "beta parameter for Newmark method"}>>> = 0.25
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_end
[]
[vel_y]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
type = NewmarkVelAux<<<{"description": "Calculates the current velocity using Newmark method.", "href": "../source/auxkernels/NewmarkVelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_y
acceleration<<<{"description": "acceleration variable"}>>> = accel_y
gamma<<<{"description": "gamma parameter for Newmark method"}>>> = 0.5
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_end
[]
[accel_z]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
type = NewmarkAccelAux<<<{"description": "Computes the current acceleration using the Newmark method.", "href": "../source/auxkernels/NewmarkAccelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_z
displacement<<<{"description": "displacement variable"}>>> = disp_z
velocity<<<{"description": "velocity variable"}>>> = vel_z
beta<<<{"description": "beta parameter for Newmark method"}>>> = 0.25
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_end
[]
[vel_z]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat RV SGs RV_slab' # 99 98 97 96 95 94'
type = NewmarkVelAux<<<{"description": "Calculates the current velocity using Newmark method.", "href": "../source/auxkernels/NewmarkVelAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_z
acceleration<<<{"description": "acceleration variable"}>>> = accel_z
gamma<<<{"description": "gamma parameter for Newmark method"}>>> = 0.5
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_end
[]
[waves]
type = WaveHeightAuxKernel<<<{"description": "Calculates the wave heights given pressures.", "href": "../source/auxkernels/WaveHeightAuxKernel.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = 'Wave1'
pressure<<<{"description": "pressure variable"}>>> = p
density<<<{"description": "fluid density"}>>> = 1e-6
gravity<<<{"description": "acceleration due to gravity"}>>> = 9.81
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_end
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'fluid_material'
[]
[]
[InterfaceKernels<<<{"href": "../syntax/InterfaceKernels/index.html"}>>>]
[interface1]
type = StructureAcousticInterface<<<{"description": "Enforces displacement and stress/pressure continuity between the fluid and structural domains. Element is always the structure and neighbor is always the fluid.", "href": "../source/interfacekernels/StructureAcousticInterface.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = disp_x
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'Interface'
D<<<{"description": "Fluid density."}>>> = 1e-6
component<<<{"description": "The desired component of displacement."}>>> = 0
[]
[interface2]
type = StructureAcousticInterface<<<{"description": "Enforces displacement and stress/pressure continuity between the fluid and structural domains. Element is always the structure and neighbor is always the fluid.", "href": "../source/interfacekernels/StructureAcousticInterface.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = disp_y
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'Interface'
D<<<{"description": "Fluid density."}>>> = 1e-6
component<<<{"description": "The desired component of displacement."}>>> = 1
[]
[interface3]
type = StructureAcousticInterface<<<{"description": "Enforces displacement and stress/pressure continuity between the fluid and structural domains. Element is always the structure and neighbor is always the fluid.", "href": "../source/interfacekernels/StructureAcousticInterface.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
neighbor_var<<<{"description": "The variable on the other side of the interface."}>>> = disp_z
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'Interface'
D<<<{"description": "Fluid density."}>>> = 1e-6
component<<<{"description": "The desired component of displacement."}>>> = 2
[]
[]
[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
[elasticity_1]
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"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat'
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 24.8
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.2
[]
[strain_1]
type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../source/materials/ComputeSmallStrain.html"}>>>
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat'
[]
[stress_1]
type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../source/materials/ComputeLinearElasticStress.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat'
[]
[den_1]
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"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases int_wall int_slab RV_housing small_walls basemat'
prop_names<<<{"description": "The names of the properties this material will have"}>>> = density
prop_values<<<{"description": "The values associated with the named properties"}>>> = 2.4e-6 #kg/m3
[]
[elasticity_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"}>>> = 'RV SGs RV_slab'
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 170
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.3
[]
[strain_2]
type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../source/materials/ComputeSmallStrain.html"}>>>
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'RV SGs RV_slab'
[]
[stress_2]
type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../source/materials/ComputeLinearElasticStress.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'RV SGs RV_slab'
[]
[den_rv]
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"}>>> = 'RV RV_slab'
prop_names<<<{"description": "The names of the properties this material will have"}>>> = density
prop_values<<<{"description": "The values associated with the named properties"}>>> = 7.85e-6 #kg/m3
[]
[den_sg]
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"}>>> = 'SGs'
prop_names<<<{"description": "The names of the properties this material will have"}>>> = density
prop_values<<<{"description": "The values associated with the named properties"}>>> = 7.85e-6 #kg/m3
[]
[density_fluid]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = inv_co_sq
prop_values<<<{"description": "The values associated with the named properties"}>>> = 4.65e-7
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'fluid_material'
[]
[dens2]
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"}>>> = 'fluid_material'
prop_names<<<{"description": "The names of the properties this material will have"}>>> = density
prop_values<<<{"description": "The values associated with the named properties"}>>> = 1e-6
[]
[]
[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
[input_motion_x]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = 'fixed_base_input_motion.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 1
[]
[input_motion_y]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = 'fixed_base_input_motion.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 2
[]
[input_motion_z]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = 'fixed_base_input_motion.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 3
[]
[]
[BCs<<<{"href": "../syntax/BCs/index.html"}>>>]
[x_motion]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
acceleration<<<{"description": "The acceleration variable."}>>> = accel_x
velocity<<<{"description": "The velocity variable."}>>> = vel_x
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_surface'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_x'
[]
[y_motion]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
acceleration<<<{"description": "The acceleration variable."}>>> = accel_y
velocity<<<{"description": "The velocity variable."}>>> = vel_y
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_surface'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_y'
[]
[z_motion]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
acceleration<<<{"description": "The acceleration variable."}>>> = accel_z
velocity<<<{"description": "The velocity variable."}>>> = vel_z
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_surface'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_z'
[]
[free]
type = FluidFreeSurfaceBC<<<{"description": "Applies a mixed Dirichlet-Neumann BC on the fluid surface.", "href": "../source/bcs/FluidFreeSurfaceBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = p
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'fluid_top'
[]
[]
[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
[Wave1]
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."}>>> = '-17.235 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave2]
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."}>>> = '-17.0 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave3]
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."}>>> = '-16.75 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave4]
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."}>>> = '-16.5 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave5]
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."}>>> = '-16.25 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave6]
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."}>>> = '-16.0 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave7]
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."}>>> = '-15.75 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave8]
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."}>>> = '-15.5 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave9]
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."}>>> = '-15.25 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave10]
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."}>>> = '-15.0 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave11]
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."}>>> = '-14.75 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave12]
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."}>>> = '-14.5 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave13]
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."}>>> = '-14.25 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave14]
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."}>>> = '-14.0 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave15]
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."}>>> = '-13.75 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave16]
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."}>>> = '-13.5 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave17]
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."}>>> = '-13.25 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave18]
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."}>>> = '-13.0 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[Wave19]
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."}>>> = '-12.765 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[p1_n]
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."}>>> = '-17.48 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p2_n]
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."}>>> = '-17.48 0.0 6.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p3_n]
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."}>>> = '-17.48 0.0 6.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p4_n]
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."}>>> = '-17.48 0.0 5.7'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p5_n]
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."}>>> = '-17.48 0.0 5.4'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p6_n]
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."}>>> = '-17.48 0.0 5.1'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p7_n]
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."}>>> = '-17.48 0.0 4.8'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p8_n]
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."}>>> = '-17.48 0.0 4.5'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p9_n]
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."}>>> = '-17.48 0.0 4.2'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p10_n]
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."}>>> = '-17.48 0.0 3.9'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p11_n]
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."}>>> = '-17.48 0.0 3.6'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p12_n]
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."}>>> = '-17.48 0.0 3.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p13_n]
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."}>>> = '-17.48 0.0 3.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p14_n]
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."}>>> = '-17.48 0.0 2.855'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p1_p]
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."}>>> = '-12.52 0.0 6.655'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p2_p]
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."}>>> = '-12.52 0.0 6.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p3_p]
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."}>>> = '-12.52 0.0 6.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p4_p]
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."}>>> = '-12.52 0.0 5.7'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p5_p]
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."}>>> = '-12.52 0.0 5.4'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p6_p]
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."}>>> = '-12.52 0.0 5.1'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p7_p]
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."}>>> = '-12.52 0.0 4.8'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p8_p]
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."}>>> = '-12.52 0.0 4.5'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p9_p]
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."}>>> = '-12.52 0.0 4.2'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p10_p]
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."}>>> = '-12.52 0.0 3.9'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p11_p]
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."}>>> = '-12.52 0.0 3.6'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p12_p]
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."}>>> = '-12.52 0.0 3.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p13_p]
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."}>>> = '-12.52 0.0 3.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[p14_p]
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."}>>> = '-12.52 0.0 2.855'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[moment_x]
type = SidesetMoment<<<{"description": "Computes the product of integrated reaction force and lever arm in a user-specified direction on a sideset from the surface traction", "href": "../source/postprocessors/SidesetMoment.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'fluid_interface'
reference_point<<<{"description": "Reference point on the sideset about which the moment is computed"}>>> = '-15.0 0.0 7.625'
moment_direction<<<{"description": "Moment direction"}>>> = '1 0 0'
pressure<<<{"description": "The scalar pressure"}>>> = p
[]
[moment_y]
type = SidesetMoment<<<{"description": "Computes the product of integrated reaction force and lever arm in a user-specified direction on a sideset from the surface traction", "href": "../source/postprocessors/SidesetMoment.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'fluid_interface'
reference_point<<<{"description": "Reference point on the sideset about which the moment is computed"}>>> = '-15.0 0.0 7.625'
moment_direction<<<{"description": "Moment direction"}>>> = '0 1 0'
pressure<<<{"description": "The scalar pressure"}>>> = p
[]
[moment_z]
type = SidesetMoment<<<{"description": "Computes the product of integrated reaction force and lever arm in a user-specified direction on a sideset from the surface traction", "href": "../source/postprocessors/SidesetMoment.html"}>>>
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'fluid_interface'
reference_point<<<{"description": "Reference point on the sideset about which the moment is computed"}>>> = '-15.0 0.0 7.625'
moment_direction<<<{"description": "Moment direction"}>>> = '0 0 1'
pressure<<<{"description": "The scalar pressure"}>>> = p
[]
[]
[VectorPostprocessors<<<{"href": "../syntax/VectorPostprocessors/index.html"}>>>]
[accel_hist_x]
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'
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'rv_head_nodes'
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[accel_spec_x]
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_x
regularize_dt<<<{"description": "dt for response spectra calculation. The acceleration response will be regularized to this dt prior to the response spectrum calculation."}>>> = 0.002
damping_ratio<<<{"description": "Damping ratio for response spectra calculation."}>>> = 0.05
start_frequency<<<{"description": "Start frequency for the response spectra calculation."}>>> = 0.1
end_frequency<<<{"description": "End frequency for the response spectra calculation."}>>> = 1000
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[accel_hist_y]
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_y'
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'rv_head_nodes'
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[accel_spec_y]
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_y
regularize_dt<<<{"description": "dt for response spectra calculation. The acceleration response will be regularized to this dt prior to the response spectrum calculation."}>>> = 0.002
damping_ratio<<<{"description": "Damping ratio for response spectra calculation."}>>> = 0.05
start_frequency<<<{"description": "Start frequency for the response spectra calculation."}>>> = 0.1
end_frequency<<<{"description": "End frequency for the response spectra calculation."}>>> = 1000
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[accel_hist_z]
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_z'
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies"}>>> = 'rv_head_nodes'
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[accel_spec_z]
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_z
regularize_dt<<<{"description": "dt for response spectra calculation. The acceleration response will be regularized to this dt prior to the response spectrum calculation."}>>> = 0.002
damping_ratio<<<{"description": "Damping ratio for response spectra calculation."}>>> = 0.05
start_frequency<<<{"description": "Start frequency for the response spectra calculation."}>>> = 0.1
end_frequency<<<{"description": "End frequency for the response spectra calculation."}>>> = 1000
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[]
[Preconditioning<<<{"href": "../syntax/Preconditioning/index.html"}>>>]
[smp]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
[]
[]
[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
type = Transient
petsc_options = '-ksp_snes_ew'
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
solve_type = 'NEWTON'
nl_rel_tol = 1e-7
nl_abs_tol = 1e-12
dt = 0.01
end_time = 28.0
timestep_tolerance = 1e-6
[TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
type = NewmarkBeta
beta = 0.25
gamma = 0.5
[]
[]
[Outputs<<<{"href": "../syntax/Mastodon/Outputs/index.html"}>>>]
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
csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
[out1]
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."}>>> = 'final'
[]
[]
(examples/ex13/fixed_base/fixed_base.i)The response history builder vectorpostprocessor was used to capture the response at various of points of interest in this model, specifcally around the reactor vessel. Figure 5 shows the acceleration response spectrum at a point on the slab directly above the center of the reactor vessel. A nodeset was created during the meshing proccess to gather struture response information like this at various points around the model, additionally node ids can be supplied to the vectorpostprocessor.

Figure 5: plot of the spectral acceleration at a point on the slab directly above the center of the reactor vessel.
The fluid pressures were measured at various locations along the wall of the reactor vessel. The locations of each measurement point are shown in figure Figure 6 below. The pressure was measured through the depth of the fluid at three instances in time, as shown in figure Figure 7, additionally it was measured over time at three points on the vessel wall as shown in figure Figure 8.

Figure 6: diagram of the reactor vessel showing the various fluid pressure measurement points.

Figure 7: Pressure profiles with depth at times 7.5, 20, and 24.66 seconds on the (a) west and (b) east faces of the vessel. (c) Wave height profiles with the radial distance from center.

Figure 8: Pressure time histories at the locations of interest on the (a) west and (b) east faces of the reactor vessel. Wave height time histories at the locations of interest on the (c) west and (d) east sides of the fluid surface.