Example 12b: Fluid-structure interaction analysis of a tank subjected to uni-directional earthquake acceleration
Giga kilograms, Giga Newtons, meters, and seconds were used for the modeling. However, the post-processed results maybe presented in different units.
Modeling
This example uses the same model and kernels as discussed in example 12a. The major difference is from example 12a is, instead of using a harmonic acceleration, an earthquake acceleration corresponding to the ElCentro earthquake is used. The input file is presented below.
[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
[file]
type = FileMeshGenerator<<<{"description": "Read a mesh from a file.", "href": "../source/meshgenerators/FileMeshGenerator.html"}>>>
file<<<{"description": "The filename to read."}>>> = '../ex12a/Model_1.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"}>>> = file
primary_block<<<{"description": "The primary set of blocks for which to draw a sideset between"}>>> = '2'
paired_block<<<{"description": "The paired set of blocks for which to draw a sideset between"}>>> = '1'
new_boundary<<<{"description": "The list of boundary names to create on the supplied subdomain"}>>> = 'Interface'
[]
[]
[GlobalParams<<<{"href": "../syntax/GlobalParams/index.html"}>>>]
[]
[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
[p]
block = 2
[]
[disp_x]
block = 1
[]
[disp_y]
block = 1
[]
[disp_z]
block = 1
[]
[]
[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
[Wave1]
block = 2
[]
[vel_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
block = 1
[]
[accel_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
block = 1
[]
[vel_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
block = 1
[]
[accel_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
block = 1
[]
[vel_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
block = 1
[]
[accel_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
block = 1
[]
[stress_xx]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
block = 1
[]
[stress_yy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
block = 1
[]
[stress_xy]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
block = 1
[]
[stress_zz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
block = 1
[]
[stress_yz]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
block = 1
[]
[stress_xz]
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
block = 1
[]
[]
[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"}>>> = 2
[]
[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"}>>> = 2
[]
[DynamicTensorMechanics<<<{"href": "../syntax/Kernels/DynamicTensorMechanics/index.html"}>>>]
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"}>>> = 1
[]
[inertia_x1]
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
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[inertia_y1]
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
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[inertia_z1]
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
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[../]
[]
[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
[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"}>>> = 2
[]
[accel_x]
type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_x
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_x
first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[vel_x]
type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_x
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_x
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[accel_y]
type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_y
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_y
first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[vel_y]
type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_y
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[accel_z]
type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_z
variable<<<{"description": "The name of the variable that this object applies to"}>>> = accel_z
first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[vel_z]
type = TestNewmarkTI<<<{"description": "Assigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.", "href": "../source/auxkernels/TestNewmarkTI.html"}>>>
displacement<<<{"description": "The variable whose first/second derivative needs to be copied to the provided auxvariable."}>>> = disp_z
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_z
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[stress_xx]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xx
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 0
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[stress_yy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[stress_xy]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xy
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 0
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[stress_zz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_zz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 2
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[stress_yz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yz
index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 2
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[stress_xz]
type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_xz
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
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[]
[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
[]
[]
[BCs<<<{"href": "../syntax/BCs/index.html"}>>>]
[bottom_accel1]
type = PresetAcceleration<<<{"description": "Prescribe acceleration on a given boundary in a given direction", "href": "../source/bcs/PresetAcceleration.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
velocity<<<{"description": "The velocity variable."}>>> = vel_x
acceleration<<<{"description": "The acceleration variable."}>>> = accel_x
beta<<<{"description": "beta parameter for Newmark time integration."}>>> = 0.25
function<<<{"description": "Function describing the velocity."}>>> = accel_bottom
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'Bottom'
[]
[disp_x2]
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"}>>> = 'Bottom'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[disp_x3]
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"}>>> = 'Bottom'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[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'
alpha<<<{"description": "Inverse of the acceleration due to gravity."}>>> = '0.1'
[]
[]
[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
[accel_bottom]
type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
data_file<<<{"description": "File holding CSV data"}>>> = SC_ElCentro_input.csv
scale_factor<<<{"description": "Scale factor to be applied to the output values"}>>> = 9.81
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = 'columns'
[]
[]
[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
[density]
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"}>>> = 2
[]
[density0]
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"}>>> = 1
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
[]
[elasticity_base]
type = ComputeIsotropicElasticityTensor<<<{"description": "Compute a constant isotropic elasticity tensor.", "href": "../source/materials/ComputeIsotropicElasticityTensor.html"}>>>
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 2e2
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.27
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[strain]
type = ComputeFiniteStrain<<<{"description": "Compute a strain increment and rotation increment for finite strains.", "href": "../source/materials/ComputeFiniteStrain.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
[]
[stress]
type = ComputeFiniteStrainElasticStress<<<{"description": "Compute stress using elasticity for finite strains", "href": "../source/materials/ComputeFiniteStrainElasticStress.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 1
[]
[]
[Preconditioning<<<{"href": "../syntax/Preconditioning/index.html"}>>>]
[andy]
type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
[]
[]
[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
start_time = 0.0
end_time = 7.0
dt = 0.005
dtmin = 0.00001
nl_abs_tol = 1e-14
nl_rel_tol = 1e-14
l_tol = 1e-14
l_max_its = 25
timestep_tolerance = 1e-8
automatic_scaling = true
[TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
type = NewmarkBeta
[]
[]
[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
[PBottom]
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.79 0.0 0.0079'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = p
[]
[Wave]
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.71 0.0 1.8079'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Wave1
[]
[]
[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
file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = Ex_SC_ElCentro
[out]
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'
type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../source/outputs/CSV.html"}>>>
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."}>>> = SC_ElCentro
[]
[]
(examples/ex12b/FSI_ex12b.i)Input motion
The bottom of the tank is subjected to an ElCentro earthquake acceleration with peak amplitude of 0.2 . The length of the input motion is about 7 . Figure 1 presents a plot of the input motion.

Figure 1: Input harmonic motion to the tank.
Results: Pressures on the tank and wave heights over the fluid
Pressures acting on the tank wall as a function of time are recorded at the circular red dot as illustrated in example 12a. Wave heights as the function of time are also recorded at the triangular red dot as illustrated in example 12a. Figure 3 and Figure 4 present the pressures and wave heights, respectively.

Figure 3: Pressures recorded at the circular dot as illustrated in example 12a.

Figure 4: Wave heights recorded at the triangular dot as illustrated in example 12a.
Results: Comparison to a reference solution
Figure 2 presents the peak pressure and peak wave height in comparison to a reference solution. It is observed that the modeled results match quite well with the reference solution from Yu and Whittaker (2020).

Figure 2: Comparison of the modeled results to a reference solution from Yu and Whittaker (2020).
References
- C. C. Yu and A. S. Whittaker.
Analytical solutions for seismic fluid-structure interaction of head-supported cylindrical tanks.
Journal of Engineering Mechanics, 146(10):04020112, 2020.[BibTeX]