Example 14: Seismic analysis of a base-isolated nuclear power plant building
GN, GPa, m, and sec
This example demonstrates the seismic analysis of a nuclear power plant (NPP) building that houses a reactor vessel and four steam generators. This NPP building is identical to that of Example 13 and is a hypothetical molten-chloride fast reactor (NPP). Details of the NPP building are provided in Example 13 as well as Bolisetti et al. (2020). This example is presented in two parts: (1) seismic analysis of just the foundation basemat on seismic isolators, and (2) seismic analysis of the full building on seismic isolators. The foundation basemat and the isolation system in both cases is identical. In each case, the modeling is briefly described and the results are presented. The mesh for both the parts of this example is generated in Cubit and the Cubit journal file (fixed_base_with_isolators_new.jou
) for these meshes is also included in the folder of this example in the repository.
Part 1: Foundation basemat analysis
Modeling
Before performing the seismic analysis of the full NPP building, the base-isolated basemat is first analyzed. The isolation system of the building comprises 70 Friction PendulumTM (FP) isolators arranged in a grid shown in Figure 3 below. Note that the figure shows a view of the foundation from the bottom, and therefore the +Z axis is pointing inwards into the basemat. The basemat is 98.75m x 68.75m in plan and is 1m thick.

Figure 3: Isolators under the basemat of the NPP building analyzed in Part 1 of this example

Figure 4: Isometric view of the NPP building mesh analyzed in Part 2 of this example.
The basemat is modeled to be almost rigid, with an elastic modulus of almost 99.2 GPa, which is four times stiffer than concrete. Each of the isolators is independently attached to the basemat with rigid beam elements, thereby simulating a rigid connection between them. The isolators themselves are 0.3m long and have a friction coefficient, mu_ref=0.06
, radius of curvature, and r_eff=1.0
. These properties result in a sliding period of 2 sec, which is reasonable for the NPP building. Two methods of modeling the isolation system are demonstrated here. In method 1, the friction coefficients of the isolators are assumed to be independent of pressure, temperature, and velocity and the mu_ref
value is used throughout the simulations. (See the theory manual for a detailed description of pressure, temperature, velocity dependency of the friction coefficient of FP isolators). Additionally, a unidirectional ground motion is applied in the X direction. In method 2, the pressure, temperature, and velocity dependencies are switched on and ground motions are applied in all three directions. The Materials blocks for the isolators for the two mwethods are listed below, and the acceleration histories of the input motions are presented in Figure 1 below.
[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
[elasticity]
type = ComputeFPIsolatorElasticity<<<{"description": "Compute the forces and the stiffness matrix for a single concave Friction Pendulum isolator element.", "href": "../source/materials/ComputeFPIsolatorElasticity.html"}>>>
mu_ref<<<{"description": "Reference co-efficient of friction."}>>> = 0.06
p_ref<<<{"description": "Reference axial pressure."}>>> = 0.05 # GPa
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
diffusivity<<<{"description": "Thermal diffusivity of steel."}>>> = 4.4e-6
conductivity<<<{"description": "Thermal conductivity of steel."}>>> = 18
a<<<{"description": "Rate parameter."}>>> = 100
r_eff<<<{"description": "Effective radius of curvature of sliding surface."}>>> = 1.0 # meters. 2sec sliding period
r_contact<<<{"description": "Radius of contact area at sliding surface."}>>> = 0.2
uy<<<{"description": "Yield displacement of the bearing in shear."}>>> = 0.001
unit<<<{"description": "Tag for conversion in the pressure factor computation when different unit systems are used. Enter 1 for N m s C; 2 for kN m s C; 3 for N mm s C; 4 for kN mm s C; 5 for lb in s C; 6 for kip in s C; 7 for lb ft s C; 8 for kip ft s C. "}>>> = 4
beta<<<{"description": "Beta parameter of Newmark algorithm."}>>> = 0.275625
gamma<<<{"description": "Gamma parameter of Newmark algorithm."}>>> = 0.55
pressure_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous pressure."}>>> = false
temperature_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous temperature."}>>> = false
velocity_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous sliding velocity."}>>> = false
k_x<<<{"description": "Stiffness of the bearing in translation along x axis. Defaults to 10e13 (SI units)"}>>> = 78.53 # 7.853e10 N
k_xx<<<{"description": "Stiffness of the bearing in rotation about x axis. Defaults to 10e13 (SI units) "}>>> = 0.62282 # 622820743.6 N
k_yy<<<{"description": "Stiffness of the bearing in rotation about y axis. Defaults to 10e13 (SI units) "}>>> = 0.3114 # 311410371.8 N
k_zz<<<{"description": "Stiffness of the bearing in rotation about z axis. Defaults to 10e13 (SI units) "}>>> = 0.3114 # 311410371.8 N
[]
[]
(examples/ex14/basemat_with_isolators_new.i)[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
[elasticity]
type = ComputeFPIsolatorElasticity<<<{"description": "Compute the forces and the stiffness matrix for a single concave Friction Pendulum isolator element.", "href": "../source/materials/ComputeFPIsolatorElasticity.html"}>>>
mu_ref<<<{"description": "Reference co-efficient of friction."}>>> = 0.06
p_ref<<<{"description": "Reference axial pressure."}>>> = 0.1 # GPa
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
diffusivity<<<{"description": "Thermal diffusivity of steel."}>>> = 4.4e-6
conductivity<<<{"description": "Thermal conductivity of steel."}>>> = 18
a<<<{"description": "Rate parameter."}>>> = 100
r_eff<<<{"description": "Effective radius of curvature of sliding surface."}>>> = 1.0 # meters. 2sec sliding period
r_contact<<<{"description": "Radius of contact area at sliding surface."}>>> = 0.2
uy<<<{"description": "Yield displacement of the bearing in shear."}>>> = 0.001
unit<<<{"description": "Tag for conversion in the pressure factor computation when different unit systems are used. Enter 1 for N m s C; 2 for kN m s C; 3 for N mm s C; 4 for kN mm s C; 5 for lb in s C; 6 for kip in s C; 7 for lb ft s C; 8 for kip ft s C. "}>>> = 4
beta<<<{"description": "Beta parameter of Newmark algorithm."}>>> = 0.275625
gamma<<<{"description": "Gamma parameter of Newmark algorithm."}>>> = 0.55
pressure_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous pressure."}>>> = true
temperature_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous temperature."}>>> = true
velocity_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous sliding velocity."}>>> = true
k_x<<<{"description": "Stiffness of the bearing in translation along x axis. Defaults to 10e13 (SI units)"}>>> = 78.53 # 7.853e10 N
k_xx<<<{"description": "Stiffness of the bearing in rotation about x axis. Defaults to 10e13 (SI units) "}>>> = 0.62282 # 622820743.6 N
k_yy<<<{"description": "Stiffness of the bearing in rotation about y axis. Defaults to 10e13 (SI units) "}>>> = 0.3114 # 311410371.8 N
k_zz<<<{"description": "Stiffness of the bearing in rotation about z axis. Defaults to 10e13 (SI units) "}>>> = 0.3114 # 311410371.8 N
[]
[]
(examples/ex14/basemat_with_isolators_new_3D.i)
Figure 1: Acceleration histories of the input ground motions used in this example.
Continuum elements such as HEX8 in 3D and QUAD4 in 2D, do not have rotational degrees of freedom. When connecting other elements such as beams or link elements (such as isolators) or shell elements, which do have rotational dofs, to continuum elements, these connection nodes are at a risk of having zero stiffness unless they involve another element with a rotational dof connected to them. In the case of this example, where each of the rigid beam elements (that connect the isolators to the basemat) are directly connected to the solid elements the rotational stiffness of this node (i.e., the Jacobian for the rotational dof of the node) is zero. Having zeros in the Jacobian can lead to significant convergence issues, and with enough zeros, non convergence. In this example, this is avoided by constraining the rotational dofs of these specific nodes as can be seen in the input file. Users are recommended to constrain such rotational dofs when their stiffness is zero to avoid convergence issues.
Gravity analysis is an integral part of seismic analysis of base-isolated structures, especially with FP bearings. The shear strength of the FP bearings depends on the normal force on the surface, which comes from gravity. In this example, gravity loads are initiated through a static initialization performed in the first 3 time steps. Static analysis is achieved by (a) switching off the velocity and acceleration calculations in the time integrator for the first couple of timesteps in the NewmarkBeta
TimeIntegrator block, and (b) switching off the inertia kernels using the controls block. Static initialization also requires that the stiffness damping is turned off in the stress divergence kernels, typically by setting static_initialization=true
when available. This option is currently only available for continuum elements through the DynamicTensorMechanics
block and the isolators through the StressDivergenceIsolator
block. Static initialization is turned on for both of these kernels in this example. Note from the input file, that time derivative calculations are turned off for the first two timesteps and the inertia kernels are switched on in the 3rd time step. It is highly recommended that users go through these parts of the input file that enable gravity simulation.
Outputs
An important output in seismic analysis of isolated systems is the force-displacement relationship from the isolators. The forces and deformations in the isolators can be recorded directly using MaterialRealCMMAux
, which can retrieve the isolator deformations, forces, stiffnesses (all of which, are of MOOSE ColumnMajorMatrix
type in their implementation; see source code) and store them in an AuxVariable
. In this example, the isolator deformations and forces in the basic co-ordinate system are recorded and presented. A sample MaterialRealCMMAux
definition for the calculation and storage of the axial forces in the isolators is shown below. On examination of the source code for ComputeFPIsolatorElasticity
it can be seen that the forces in the basic co-ordinate system are stored in the material property named basic_forces
. This is 6 x 1 matrix that stores the forces and moments of the isolators. The first element (row=0
and column=0
) corresponds to the axial forces in the isolators and is calculated and stored here. Similarly, the deformations are also evaluated in other blocks.
[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
[Fb_x]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = basic_forces
row<<<{"description": "The row component to consider for this kernel"}>>> = 0
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Fb_x
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[]
(examples/ex14/basemat_with_isolators_new.i)Both the acceleration responses of the basemat and the force-displacement relationship of one of the isolators is presented in this section. In Figure 5 and Figure 6, the orange and blue lines correspond to methods 1 and 2, respectively. Note that when the pressure, temperature, and velocity dependencies are switched on (method 2, blue) the hysteresis is a lot more complicated and the stiffness varies significantly during the earthquake. Note that the results for method 1 do not appear in Figure 6 because they are almost zero, given that the ground motion for method 1 is unidirectional. The results presented here are for an isolator at the center of the basemat plan.

Figure 5: Shear force-displacement history in the XZ direction. Orange is method 1 and blue is method 2.

Figure 6: Shear force-displacement history in the YZ direction. Orange is method 1 and blue is method 2.
The acceleration response spectra calculated at the center of the basemat are presented in the figures below. The top row presents the input acceleration spectra and the bottom row presents the output spectra at the center of the basemat in X, Y, and Z directions. Again, in these figures, blue plots present method 1 and orange plots present method 2. The figures show that the accelerations at the top of the basemat reduce drastically due to seismic isolation, except in the Z direction, in which, the response almost stays the same, demonstrating a typical isolation response with FP isolation systems.

Figure 7: 5% damped acceleration response spectra in X direction (input)

Figure 8: 5% damped acceleration response spectra in X direction (basemat center)

Figure 9: 5% damped acceleration response spectra in Y direction (input)

Figure 10: 5% damped acceleration response spectra in Y direction (basemat center)

Figure 11: 5% damped acceleration response spectra in Z direction (input)

Figure 12: 5% damped acceleration response spectra in Z direction (basemat center)
Part 2: NPP building analysis
Now that the seismic response with just the basemat is shown to be reasonable, the modeling and response of the seismically isolated building is presented in this section.
Modeling
The finite-element mesh of the building, developed in Cubit is presented in Figure 4 above. The building is founded on the basemat and the isolation system presented in Figure 3. The figures below show the plan view and an isometric view of the internal components of the building. The building is a one-story shear wall structure with buttresses on all four sides and the roof. It houses four steam generators (shown in purple and modeled as solid cylinders for simplicity) supported by individual concrete bases (in yellow), and head-supported reactor vessel is suspended by an internal concrete slab. The reactor vessel is contained in a concrete cylindrical housing that can be seen in red the figures below. The reactor vessel contains molten salt with properties that are assumed to be the same as water. The reactor vessel, as seen from the bottom without its housing or the internal walls, is shown in Figure 2 below.

Figure 13: Isometric view of the internal components of the NPP building.

Figure 14: Plan view of the internal components of the NPP building.

Figure 2: Head-supported reactor vessel as seen from the bottom.
All the materials in the building are modeled using a linear elastic material and 3D 8-noded HEX elements, except for the FP isolators, which are modeled using two noded link elements with the FP isolator material model (see theory and user manuals for more information). The reactor vessel is modeled as a thin cylindrical vessel using continuum elements. The fluid inside the reactor vessel is modeled using a linear elastic material with a very small shear modulus to reproduce fluid-like behavior. A more comprehensive, fluid-structure interaction (FSI) approach is also possible in MASTODON to model the fluid. Examples describing the FSI approach are presented in Example 13. More information on FSI modeling in MOOSE and MASTODON is provided here. Further information on the building model is also provided in Bolisetti et al. (2020).
In this demonstration, the building is subjected to ground motions in X, Y, and Z directions at the base of the isolation system. These ground motions are presented in Figure 1 above, and their spectral accelerations are shown in Figure 7, Figure 9, and Figure 11 above. Pressure, temperature, and velocity dependencies of the isolators are switched on, as described in method 2 of Part 1 of this example. The first three timesteps of the analysis involve a static gravity analysis as described in the note above.
The input file for the simulation of Part 2 is listed below.
[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."}>>> = full_structure_with_isolators_new.e
[]
[]
[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
[disp_x]
[]
[disp_y]
[]
[disp_z]
[]
[rot_x]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_y]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_z]
block = 'isolator_elems upper_rigid_elems'
[]
[]
[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
[vel_x]
[]
[accel_x]
[]
[vel_y]
[]
[accel_y]
[]
[vel_z]
[]
[accel_z]
[]
[rot_vel_x]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_vel_y]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_vel_z]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_accel_x]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_accel_y]
block = 'isolator_elems upper_rigid_elems'
[]
[rot_accel_z]
block = 'isolator_elems upper_rigid_elems'
[]
[Fb_x]
block = 'isolator_elems'
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
[]
[Fb_y]
block = 'isolator_elems'
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
[]
[Fb_z]
block = 'isolator_elems'
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
[]
[Defb_x]
block = 'isolator_elems'
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
[]
[Velb_x]
block = 'isolator_elems'
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
[]
[Defb_y]
block = 'isolator_elems'
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
[]
[Defb_z]
block = 'isolator_elems'
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
[]
[]
[Physics<<<{"href": "../syntax/Physics/index.html"}>>>/SolidMechanics<<<{"href": "../syntax/Physics/SolidMechanics/index.html"}>>>/LineElement<<<{"href": "../syntax/Physics/SolidMechanics/LineElement/index.html"}>>>/QuasiStatic<<<{"href": "../syntax/Physics/SolidMechanics/LineElement/QuasiStatic/index.html"}>>>]
displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotations appropriate for the simulation geometry and coordinate system"}>>> = 'rot_x rot_y rot_z'
velocities<<<{"description": "Translational velocity variables"}>>> = 'vel_x vel_y vel_z'
accelerations<<<{"description": "Translational acceleration variables"}>>> = 'accel_x accel_y accel_z'
rotational_velocities<<<{"description": "Rotational velocity variables"}>>> = 'rot_vel_x rot_vel_y rot_vel_z'
rotational_accelerations<<<{"description": "Rotational acceleration variables"}>>> = 'rot_accel_x rot_accel_y rot_accel_z'
beta<<<{"description": "beta parameter for Newmark Time integration"}>>> = 0.275625
gamma<<<{"description": "gamma parameter for Newmark Time integration"}>>> = 0.55
alpha<<<{"description": "alpha parameter for mass dependent numerical damping induced by HHT time integration scheme"}>>> = -0.05
[rigid_beams]
block<<<{"description": "The list of ids of the blocks (subdomain) that the stress divergence, inertia kernels and materials will be applied to"}>>> = 'upper_rigid_elems'
area<<<{"description": "Cross-section area of the beam. Can be supplied as either a number or a variable name."}>>> = 130.06
Iy<<<{"description": "Second moment of area of the beam about y axis. Can be supplied as either a number or a variable name."}>>> = 24166.729
Iz<<<{"description": "Second moment of area of the beam about z axis. Can be supplied as either a number or a variable name."}>>> = 24166.729
y_orientation<<<{"description": "Orientation of the y direction along which Iyy is provided. This should be perpendicular to the axis of the beam."}>>> = '0.0 0.0 1.0'
[]
[]
[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
[DynamicTensorMechanics<<<{"href": "../syntax/Kernels/DynamicTensorMechanics/index.html"}>>>]
displacements<<<{"description": "The nonlinear displacement variables for the problem"}>>> = 'disp_x disp_y disp_z'
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 SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
hht_alpha<<<{"description": "alpha parameter for mass dependent numerical damping induced by HHT time integration scheme"}>>> = -0.05
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
stiffness_damping_coefficient<<<{"description": "Name of material property or a constant real number defining stiffness Rayleigh parameter (zeta)."}>>> = 0.0019
[]
[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"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
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."}>>> = 0.038
alpha<<<{"description": "alpha parameter for mass dependent numerical damping induced by HHT time integration scheme"}>>> = -0.05
[]
[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"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
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."}>>> = 0.038
alpha<<<{"description": "alpha parameter for mass dependent numerical damping induced by HHT time integration scheme"}>>> = -0.05
[]
[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"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
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."}>>> = 0.038
alpha<<<{"description": "alpha parameter for mass dependent numerical damping induced by HHT time integration scheme"}>>> = -0.05
[]
[lr_disp_x]
type = StressDivergenceIsolator<<<{"description": "Kernel for isolator element", "href": "../source/kernels/StressDivergenceIsolator.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
displacements<<<{"description": "The displacement variables for isolator."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables for the isolator."}>>> = 'rot_x rot_y rot_z'
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z, 3 for rot_x, 4 for rot_y and 5 for rot_z)."}>>> = 0
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
static_initialization<<<{"description": "Set to true to get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step"}>>> = true
zeta<<<{"description": "Name of material property or a constant real number defining the zeta parameter for the Rayleigh damping."}>>> = 0.0019
alpha<<<{"description": "alpha parameter for HHT time integration"}>>> = -0.05
[]
[lr_disp_y]
type = StressDivergenceIsolator<<<{"description": "Kernel for isolator element", "href": "../source/kernels/StressDivergenceIsolator.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
displacements<<<{"description": "The displacement variables for isolator."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables for the isolator."}>>> = 'rot_x rot_y rot_z'
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z, 3 for rot_x, 4 for rot_y and 5 for rot_z)."}>>> = 1
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
static_initialization<<<{"description": "Set to true to get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step"}>>> = true
zeta<<<{"description": "Name of material property or a constant real number defining the zeta parameter for the Rayleigh damping."}>>> = 0.0019
alpha<<<{"description": "alpha parameter for HHT time integration"}>>> = -0.05
[]
[lr_disp_z]
type = StressDivergenceIsolator<<<{"description": "Kernel for isolator element", "href": "../source/kernels/StressDivergenceIsolator.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
displacements<<<{"description": "The displacement variables for isolator."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables for the isolator."}>>> = 'rot_x rot_y rot_z'
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z, 3 for rot_x, 4 for rot_y and 5 for rot_z)."}>>> = 2
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
static_initialization<<<{"description": "Set to true to get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step"}>>> = true
zeta<<<{"description": "Name of material property or a constant real number defining the zeta parameter for the Rayleigh damping."}>>> = 0.0019
alpha<<<{"description": "alpha parameter for HHT time integration"}>>> = -0.05
[]
[lr_rot_x]
type = StressDivergenceIsolator<<<{"description": "Kernel for isolator element", "href": "../source/kernels/StressDivergenceIsolator.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
displacements<<<{"description": "The displacement variables for isolator."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables for the isolator."}>>> = 'rot_x rot_y rot_z'
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z, 3 for rot_x, 4 for rot_y and 5 for rot_z)."}>>> = 3
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = rot_x
static_initialization<<<{"description": "Set to true to get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step"}>>> = true
zeta<<<{"description": "Name of material property or a constant real number defining the zeta parameter for the Rayleigh damping."}>>> = 0.0019
alpha<<<{"description": "alpha parameter for HHT time integration"}>>> = -0.05
[]
[lr_rot_y]
type = StressDivergenceIsolator<<<{"description": "Kernel for isolator element", "href": "../source/kernels/StressDivergenceIsolator.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
displacements<<<{"description": "The displacement variables for isolator."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables for the isolator."}>>> = 'rot_x rot_y rot_z'
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z, 3 for rot_x, 4 for rot_y and 5 for rot_z)."}>>> = 4
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = rot_y
static_initialization<<<{"description": "Set to true to get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step"}>>> = true
zeta<<<{"description": "Name of material property or a constant real number defining the zeta parameter for the Rayleigh damping."}>>> = 0.0019
alpha<<<{"description": "alpha parameter for HHT time integration"}>>> = -0.05
[]
[lr_rot_z]
type = StressDivergenceIsolator<<<{"description": "Kernel for isolator element", "href": "../source/kernels/StressDivergenceIsolator.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
displacements<<<{"description": "The displacement variables for isolator."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables for the isolator."}>>> = 'rot_x rot_y rot_z'
component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z, 3 for rot_x, 4 for rot_y and 5 for rot_z)."}>>> = 5
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = rot_z
static_initialization<<<{"description": "Set to true to get the system to equilibrium under gravity by running a quasi-static analysis (by solving Ku = F) in the first time step"}>>> = true
zeta<<<{"description": "Name of material property or a constant real number defining the zeta parameter for the Rayleigh damping."}>>> = 0.0019
alpha<<<{"description": "alpha parameter for HHT time integration"}>>> = -0.05
[]
[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
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
alpha<<<{"description": "alpha parameter required for HHT time integration scheme"}>>> = -0.05
[]
[]
[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
[Fb_x]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = basic_forces
row<<<{"description": "The row component to consider for this kernel"}>>> = 0
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Fb_x
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[Fb_y]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = basic_forces
row<<<{"description": "The row component to consider for this kernel"}>>> = 1
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Fb_y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[Fb_z]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = basic_forces
row<<<{"description": "The row component to consider for this kernel"}>>> = 2
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Fb_z
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[Defb_x]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = deformations
row<<<{"description": "The row component to consider for this kernel"}>>> = 0
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Defb_x
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[Velb_x]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = deformation_rates
row<<<{"description": "The row component to consider for this kernel"}>>> = 0
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Velb_x
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[Defb_y]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = deformations
row<<<{"description": "The row component to consider for this kernel"}>>> = 1
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Defb_y
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[Defb_z]
type = MaterialRealCMMAux<<<{"description": "Populate an auxiliary variable with an entry from a ColumnMajorMatrix material property.", "href": "../source/auxkernels/MaterialRealCMMAux.html"}>>>
property<<<{"description": "The material property name."}>>> = deformations
row<<<{"description": "The row component to consider for this kernel"}>>> = 2
column<<<{"description": "The column component to consider for this kernel"}>>> = 0
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Defb_z
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[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
[]
[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
[]
[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
[]
[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
[]
[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
[]
[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
[]
[rot_accel_x]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems upper_rigid_elems'
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."}>>> = rot_x
variable<<<{"description": "The name of the variable that this object applies to"}>>> = rot_accel_x
first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
[]
[rot_vel_x]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems upper_rigid_elems'
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."}>>> = rot_x
variable<<<{"description": "The name of the variable that this object applies to"}>>> = rot_vel_x
[]
[rot_accel_y]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems upper_rigid_elems'
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."}>>> = rot_y
variable<<<{"description": "The name of the variable that this object applies to"}>>> = rot_accel_y
first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
[]
[rot_vel_y]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems upper_rigid_elems'
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."}>>> = rot_y
variable<<<{"description": "The name of the variable that this object applies to"}>>> = rot_vel_y
[]
[rot_accel_z]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems upper_rigid_elems'
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."}>>> = rot_z
variable<<<{"description": "The name of the variable that this object applies to"}>>> = rot_accel_z
first<<<{"description": "Set to true to copy the first derivative to the auxvariable. If false, the second derivative is copied."}>>> = false
[]
[rot_vel_z]
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems upper_rigid_elems'
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."}>>> = rot_z
variable<<<{"description": "The name of the variable that this object applies to"}>>> = rot_vel_z
[]
[]
[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
[elasticity_concrete]
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 RV_slab'
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 24.8 #GPa
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.2
[]
[elasticity_rigid_concrete]
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"}>>> = 'upper_basemat'
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 99.2 #GPa # 4 x concrete for rigid basemat
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.2
[]
[elasticity_steel_316]
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"}>>> = 'SGs RV'
youngs_modulus<<<{"description": "Young's modulus of the material."}>>> = 170 #GPa
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.3
[]
[elasticity_fluid]
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"}>>> = 'fluid_material'
bulk_modulus<<<{"description": "The bulk modulus for the material."}>>> = 2 #GPa #water
poissons_ratio<<<{"description": "Poisson's ratio for the material."}>>> = 0.45 #water
[]
[strain_1]
type = ComputeFiniteStrain<<<{"description": "Compute a strain increment and rotation increment for finite strains.", "href": "../source/materials/ComputeFiniteStrain.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 SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
[]
[stress_1]
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"}>>> = 'roof ext_buttresses ext_walls int_buttresses SG_bases SGs int_wall int_slab RV_housing RV small_walls upper_basemat fluid_material RV_slab'
[]
[concrete_density]
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 upper_basemat 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"}>>> = 2.4e-6 #e9kg/m3
[]
[steel_density]
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 RV'
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 #e9kg/m3
[]
[fluid_density]
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"}>>> = 1.0e-6 #e9kg/m3 #water
[]
[isolator_deformation]
type = ComputeIsolatorDeformation<<<{"description": "Compute the deformations and rotations in a two-noded isolator element.", "href": "../source/materials/ComputeIsolatorDeformation.html"}>>>
sd_ratio<<<{"description": "Shear distance ratio."}>>> = 0.5
y_orientation<<<{"description": "Orientation of the local Y direction along which, Ky is provided. This should be perpendicular to the axis of the isolator."}>>> = '1.0 0.0 0.0'
displacements<<<{"description": "The displacement variables appropriate for the simulation geometry and coordinate system."}>>> = 'disp_x disp_y disp_z'
rotations<<<{"description": "The rotation variables appropriate for the simulation geometry and coordinate system."}>>> = 'rot_x rot_y rot_z'
velocities<<<{"description": "Translational velocity variables."}>>> = 'vel_x vel_y vel_z'
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
[]
[elasticity]
type = ComputeFPIsolatorElasticity<<<{"description": "Compute the forces and the stiffness matrix for a single concave Friction Pendulum isolator element.", "href": "../source/materials/ComputeFPIsolatorElasticity.html"}>>>
mu_ref<<<{"description": "Reference co-efficient of friction."}>>> = 0.06
p_ref<<<{"description": "Reference axial pressure."}>>> = 0.006 # GPa
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'isolator_elems'
diffusivity<<<{"description": "Thermal diffusivity of steel."}>>> = 4.4e-6
conductivity<<<{"description": "Thermal conductivity of steel."}>>> = 18
a<<<{"description": "Rate parameter."}>>> = 100
r_eff<<<{"description": "Effective radius of curvature of sliding surface."}>>> = 1.0 # meters. 2sec sliding period
r_contact<<<{"description": "Radius of contact area at sliding surface."}>>> = 0.2
uy<<<{"description": "Yield displacement of the bearing in shear."}>>> = 0.001
unit<<<{"description": "Tag for conversion in the pressure factor computation when different unit systems are used. Enter 1 for N m s C; 2 for kN m s C; 3 for N mm s C; 4 for kN mm s C; 5 for lb in s C; 6 for kip in s C; 7 for lb ft s C; 8 for kip ft s C. "}>>> = 4
beta<<<{"description": "Beta parameter of Newmark algorithm."}>>> = 0.275625
gamma<<<{"description": "Gamma parameter of Newmark algorithm."}>>> = 0.55
pressure_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous pressure."}>>> = false
temperature_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous temperature."}>>> = false
velocity_dependent<<<{"description": "Switch for modeling friction dependence on the instantaneous sliding velocity."}>>> = false
k_x<<<{"description": "Stiffness of the bearing in translation along x axis. Defaults to 10e13 (SI units)"}>>> = 78.53 # 7.853e10 N
k_xx<<<{"description": "Stiffness of the bearing in rotation about x axis. Defaults to 10e13 (SI units) "}>>> = 0.62282 # 622820743.6 N
k_yy<<<{"description": "Stiffness of the bearing in rotation about y axis. Defaults to 10e13 (SI units) "}>>> = 0.3114 # 311410371.8 N
k_zz<<<{"description": "Stiffness of the bearing in rotation about z axis. Defaults to 10e13 (SI units) "}>>> = 0.3114 # 311410371.8 N
[]
[elasticity_beam_rigid]
type = ComputeElasticityBeam<<<{"description": "Computes the equivalent of the elasticity tensor for the beam element, which are vectors of material translational and flexural stiffness.", "href": "../source/materials/ComputeElasticityBeam.html"}>>>
youngs_modulus<<<{"description": "Young's modulus of the material. Can be supplied as either a number or a variable name."}>>> = 2e4
poissons_ratio<<<{"description": "Poisson's ratio of the material. Can be supplied as either a number or a variable name."}>>> = 0.27
shear_coefficient<<<{"description": "Scale factor for the shear modulus. Can be supplied as either a number or a variable name."}>>> = 0.85
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'upper_rigid_elems'
[]
[stress_beam_rigid]
type = ComputeBeamResultants<<<{"description": "Compute forces and moments using elasticity", "href": "../source/materials/ComputeBeamResultants.html"}>>>
block<<<{"description": "The list of blocks (ids or names) that this object will be applied"}>>> = 'upper_rigid_elems'
[]
[]
[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"}>>> = 'case2_scaled.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
scale_factor<<<{"description": "Scale factor to be applied to the output values"}>>> = 9.81
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 1
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
[]
[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"}>>> = 'case2_scaled.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
scale_factor<<<{"description": "Scale factor to be applied to the output values"}>>> = 9.81
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 2
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
[]
[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"}>>> = 'case2_scaled.csv'
format<<<{"description": "Format of csv data file that is in either in columns or rows"}>>> = columns
scale_factor<<<{"description": "Scale factor to be applied to the output values"}>>> = 9.81
y_index_in_file<<<{"description": "The ordinate index in the data file"}>>> = 3
xy_in_file_only<<<{"description": "If the data file only contains abscissa and ordinate data"}>>> = false
[]
[]
[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.2725625
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_isolators'
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.2725625
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_isolators'
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.2725625
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_isolators'
function<<<{"description": "Function describing the velocity."}>>> = 'input_motion_z'
[]
[fixrxbot]
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"}>>> = rot_x
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_isolators'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[fixrybot]
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"}>>> = rot_y
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_isolators'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[fixrzbot]
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"}>>> = rot_z
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom_isolators'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[fixrxcon]
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"}>>> = rot_x
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'connections'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[fixrycon]
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"}>>> = rot_y
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'connections'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[fixrzcon]
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"}>>> = rot_z
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'connections'
value<<<{"description": "Value of the BC"}>>> = 0.0
[]
[]
[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-15
dt = 0.01
end_time = 28
timestep_tolerance = 1e-6
automatic_scaling = true
[TimeIntegrator<<<{"href": "../syntax/Executioner/TimeIntegrator/index.html"}>>>]
type = NewmarkBeta
beta = 0.275625
gamma = 0.55
inactive_tsteps = 2
[]
[]
[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.03
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
*/rot_vel_x */rot_vel_y */rot_vel_z
*/rot_accel_x */rot_accel_y */rot_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'
[]
[]
[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
[inp_accel_x]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_x
[]
[inp_accel_y]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_y
[]
[inp_accel_z]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_z
[]
[basemat_accel_x]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 0.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_x
[]
[basemat_accel_y]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 0.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_y
[]
[basemat_accel_z]
type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0.0 0.0'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = accel_z
[]
[iso1_fb_axial]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Fb_x
[]
[iso1_defb_axial]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Defb_x
[]
[iso1_fb_shear1]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Fb_y
[]
[iso1_defb_shear1]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Defb_y
[]
[iso1_fb_shear2]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Fb_z
[]
[iso1_defb_shear2]
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."}>>> = '5.0 0.0 -1.3'
variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = Defb_z
[]
[]
[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'
nodes<<<{"description": "Node number(s) at which the response history is needed."}>>> = '5252 2767 59044 24207 44503 41781 59152 38767 59100'
# locations:
# 5252-roof_edge
# 2767-roof_center
# 59044-RV_slab_center
# 24207-SG_base
# 44503-basemat_center-(0.35,-0.75,-1)(approx)
# 41781-center_isolator_top-(5,0,-1)
# 59152-center_isolator_bottom-(5,0,-1.3)
# 38767-edge_isolator_top-(-45,30,-1)
# 59100-edge_isolator_bottom(-45,30,-1.3)
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.01
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."}>>> = 50
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'
nodes<<<{"description": "Node number(s) at which the response history is needed."}>>> = '5252 2767 59044 24207 44503 41781 59152 38767 59100'
# locations:
# 5252-roof_edge
# 2767-roof_center
# 59044-RV_slab_center
# 24207-SG_base
# 44503-basemat_center-(0.35,-0.75,-1)(approx)
# 41781-center_isolator_top-(5,0,-1)
# 59152-center_isolator_bottom-(5,0,-1.3)
# 38767-edge_isolator_top-(-45,30,-1)
# 59100-edge_isolator_bottom(-45,30,-1.3)
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.01
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."}>>> = 50
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'
nodes<<<{"description": "Node number(s) at which the response history is needed."}>>> = '5252 2767 59044 24207 44503 41781 59152 38767 59100'
# locations:
# 5252-roof_edge
# 2767-roof_center
# 59044-RV_slab_center
# 24207-SG_base
# 44503-basemat_center-(0.35,-0.75,-1)(approx)
# 41781-center_isolator_top-(5,0,-1)
# 59152-center_isolator_bottom-(5,0,-1.3)
# 38767-edge_isolator_top-(-45,30,-1)
# 59100-edge_isolator_bottom(-45,30,-1.3)
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.01
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."}>>> = 50
outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = out1
[]
[]
[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/ex14/building_basemat_with_isolators_new.i)Outputs
The output locations and responses for Part 2 include the same as those of Part 1 (an isolator at the center of the isolation system and the center of the basemat). In addition to these responses, the acceleration response at the roof, at the center of the reactor head, and the base of one of the steam generators is shown below.
The figures below present the isolator shear responses in XZ and YZ directions calculated in the same manner described in Part 1.

Figure 15: Shear force-displacement history in the XZ direction for an isolator at the center of the isolation system.

Figure 16: Shear force-displacement history in the YZ direction for an isolator at the center of the isolation system.
The figure below present the building responses for the same input ground motions shown in Part 1. The first row presents the spectral accelerations in all three directions at a node located at the center of the basemat of the building (same node as Part 1). The second, third, and fourth rows present these responses calculated at the center of the roof of the building, center of the reactor head, and at the base of one of the steam generators, respectively.

Figure 17: 5% damped acceleration response spectra in X direction (basemat center)

Figure 18: 5% damped acceleration response spectra in Y direction (basemat center)

Figure 19: 5% damped acceleration response spectra in Z direction (basemat center)

Figure 20: 5% damped acceleration response spectra in X direction (center of the building roof)

Figure 21: 5% damped acceleration response spectra in Y direction (center of the building roof)

Figure 22: 5% damped acceleration response spectra in Z direction (center of the building roof)

Figure 23: 5% damped acceleration response spectra in X direction (center of reactor vessel head)

Figure 24: 5% damped acceleration response spectra in Y direction (center of reactor vessel head)

Figure 25: 5% damped acceleration response spectra in Z direction (center of reactor vessel head)

Figure 26: 5% damped acceleration response spectra in X direction (base of one of the steam generators)

Figure 27: 5% damped acceleration response spectra in Y direction (base of one of the steam generators)

Figure 28: 5% damped acceleration response spectra in Z direction (base of one of the steam generators)
References
- C. Bolisetti, W. Hoffman, J. L. Coleman, S. S. Parsi, K. Lal, A. S. Whittaker, M. Cohen, K. Kramer, P. Kirchman, H. Bowers, and J. Redd.
Seismic isolation of major advanced reactor systems for economic improvement and safety assurance.
Technical Report INL/EXT-20-59608, Idaho National Laboratory, 2020.[BibTeX]