- beta0.25beta value
Default:0.25
C++ Type:double
Description:beta value
- gamma0.5gamma value
Default:0.5
C++ Type:double
Description:gamma value
NewmarkBeta
Computes the first and second time derivative of variable using Newmark-Beta method.
Description
Newmark time integration (Newmark, 1959) is one of the commonly used time integration methods in structural dynamics problems. In this method, the second () and first () time derivatives of a variable at are written in terms of the , and at time , and at as shown below:
In the above equations, and are Newmark time integration parameters.
For and , the Newmark time integration method is implicit, unconditionally stable and second order accurate in time. This is the constant average acceleration method with no numerical damping.
and results in the linear acceleration method where the acceleration is linearly varying between and . This method is also implicit, unconditionally stable and second order accurate in time. However, there is a small numerical damping when the linear acceleration method is used.
For , the method is second order accurate and it is unconditionally stable for .
When using the constant average acceleration method that has no numerical damping, high frequency noise can sometimes be observed in the velocity and acceleration time histories for a problem with prescribed displacement (Bathe and Noh, 2012). Using other parameters for and results in non-zero numerical damping that damps out part of the high frequency noise but not all of it. Hilber-Hughes-Taylor (HHT) time integration is a variation of the Newmark method that damps out high frequency noise especially in structural dynamics problems. More details about this Newmark and HHT time integration schemes can be found in these lecture notes. HHT time integration requires modification to the equation of motion and is currently implemented only for structural dynamics problems in tensor mechanics module.
Input Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
Advanced Parameters
Input Files
- modules/tensor_mechanics/test/tests/beam/dynamic/dyn_euler_small_rayleigh_hht_ti.i
- modules/tensor_mechanics/test/tests/dynamics/wave_1D/wave_rayleigh_hht_ti.i
- test/tests/kernels/vector_dot_dot/vector_test.i
- test/tests/variables/time_derivatives_neighbor/test.i
- test/tests/time_integrators/newmark-beta/newmark_beta_prescribed_parameters.i
modules/tensor_mechanics/test/tests/beam/dynamic/dyn_euler_small_rayleigh_hht_ti.i
# Test for damped small strain euler beam vibration in y direction
# An impulse load is applied at the end of a cantilever beam of length 4m.
# The properties of the cantilever beam are as follows:
# Young's modulus (E) = 1e4
# Shear modulus (G) = 4e7
# Shear coefficient (k) = 1.0
# Cross-section area (A) = 0.01
# Iy = 1e-4 = Iz
# Length (L)= 4 m
# density (rho) = 1.0
# mass proportional rayleigh damping(eta) = 0.1
# stiffness proportional rayleigh damping(eta) = 0.1
# HHT time integration parameter (alpha) = -0.3
# Corresponding Newmark beta time integration parameters beta = 0.4225 and gamma = 0.8
# For this beam, the dimensionless parameter alpha = kAGL^2/EI = 6.4e6
# Therefore, the behaves like a Euler-Bernoulli beam.
# The displacement time history from this analysis matches with that obtained from Abaqus.
# Values from the first few time steps are as follows:
# time disp_y vel_y accel_y
# 0.0 0.0 0.0 0.0
# 0.2 0.019898364318588 0.18838688112273 1.1774180070171
# 0.4 0.045577003505278 0.087329917525455 -0.92596052423724
# 0.6 0.063767907208218 0.084330765885995 0.21274543331268
# 0.8 0.073602908614573 0.020029576220975 -0.45506879373455
# 1.0 0.06841704414745 -0.071840076837194 -0.46041813317992
[Mesh]
type = GeneratedMesh
nx = 10
dim = 1
xmin = 0.0
xmax = 4.0
displacements = 'disp_x disp_y disp_z'
[]
[Variables]
[./disp_x]
order = FIRST
family = LAGRANGE
[../]
[./disp_y]
order = FIRST
family = LAGRANGE
[../]
[./disp_z]
order = FIRST
family = LAGRANGE
[../]
[./rot_x]
order = FIRST
family = LAGRANGE
[../]
[./rot_y]
order = FIRST
family = LAGRANGE
[../]
[./rot_z]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./vel_x]
order = FIRST
family = LAGRANGE
[../]
[./vel_y]
order = FIRST
family = LAGRANGE
[../]
[./vel_z]
order = FIRST
family = LAGRANGE
[../]
[./accel_x]
order = FIRST
family = LAGRANGE
[../]
[./accel_y]
order = FIRST
family = LAGRANGE
[../]
[./accel_z]
order = FIRST
family = LAGRANGE
[../]
[./rot_vel_x]
order = FIRST
family = LAGRANGE
[../]
[./rot_vel_y]
order = FIRST
family = LAGRANGE
[../]
[./rot_vel_z]
order = FIRST
family = LAGRANGE
[../]
[./rot_accel_x]
order = FIRST
family = LAGRANGE
[../]
[./rot_accel_y]
order = FIRST
family = LAGRANGE
[../]
[./rot_accel_z]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxKernels]
[./accel_x] # These auxkernels are only to check output
type = TestNewmarkTI
displacement = disp_x
variable = accel_x
first = false
[../]
[./accel_y]
type = TestNewmarkTI
displacement = disp_y
variable = accel_y
first = false
[../]
[./accel_z]
type = TestNewmarkTI
displacement = disp_z
variable = accel_z
first = false
[../]
[./vel_x]
type = TestNewmarkTI
displacement = disp_x
variable = vel_x
[../]
[./vel_y]
type = TestNewmarkTI
displacement = disp_y
variable = vel_y
[../]
[./vel_z]
type = TestNewmarkTI
displacement = disp_z
variable = vel_z
[../]
[./rot_accel_x]
type = TestNewmarkTI
displacement = rot_x
variable = rot_accel_x
first = false
[../]
[./rot_accel_y]
type = TestNewmarkTI
displacement = rot_y
variable = rot_accel_y
first = false
[../]
[./rot_accel_z]
type = TestNewmarkTI
displacement = rot_z
variable = rot_accel_z
first = false
[../]
[./rot_vel_x]
type = TestNewmarkTI
displacement = rot_x
variable = rot_vel_x
[../]
[./rot_vel_y]
type = TestNewmarkTI
displacement = rot_y
variable = rot_vel_y
[../]
[./rot_vel_z]
type = TestNewmarkTI
displacement = rot_z
variable = rot_vel_z
[../]
[]
[BCs]
[./fixx1]
type = DirichletBC
variable = disp_x
boundary = left
value = 0.0
[../]
[./fixy1]
type = DirichletBC
variable = disp_y
boundary = left
value = 0.0
[../]
[./fixz1]
type = DirichletBC
variable = disp_z
boundary = left
value = 0.0
[../]
[./fixr1]
type = DirichletBC
variable = rot_x
boundary = left
value = 0.0
[../]
[./fixr2]
type = DirichletBC
variable = rot_y
boundary = left
value = 0.0
[../]
[./fixr3]
type = DirichletBC
variable = rot_z
boundary = left
value = 0.0
[../]
[]
[NodalKernels]
[./force_y2]
type = UserForcingFunctionNodalKernel
variable = disp_y
boundary = right
function = force
[../]
[]
[Functions]
[./force]
type = PiecewiseLinear
x = '0.0 0.2 0.4 10.0'
y = '0.0 0.01 0.0 0.0'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
line_search = 'none'
l_tol = 1e-11
nl_max_its = 15
nl_rel_tol = 1e-10
nl_abs_tol = 1e-10
start_time = 0.0
dt = 0.2
end_time = 5.0
timestep_tolerance = 1e-6
# Time integrator
[./TimeIntegrator]
type = NewmarkBeta
beta = 0.4225
gamma = 0.8
[../]
[]
[Kernels]
[./solid_disp_x]
type = StressDivergenceBeam
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
component = 0
variable = disp_x
zeta = 0.1
alpha = -0.3
[../]
[./solid_disp_y]
type = StressDivergenceBeam
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
component = 1
variable = disp_y
zeta = 0.1
alpha = -0.3
[../]
[./solid_disp_z]
type = StressDivergenceBeam
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
component = 2
variable = disp_z
zeta = 0.1
alpha = -0.3
[../]
[./solid_rot_x]
type = StressDivergenceBeam
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
component = 3
variable = rot_x
zeta = 0.1
alpha = -0.3
[../]
[./solid_rot_y]
type = StressDivergenceBeam
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
component = 4
variable = rot_y
zeta = 0.1
alpha = -0.3
[../]
[./solid_rot_z]
type = StressDivergenceBeam
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
component = 5
variable = rot_z
zeta = 0.1
alpha = -0.3
[../]
[./inertial_force_x]
type = InertialForceBeam
block = 0
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
eta = 0.1
area = 0.01
Iy = 1e-4
Iz = 1e-4
Ay = 0.0
Az = 0.0
component = 0
variable = disp_x
alpha = -0.3
[../]
[./inertial_force_y]
type = InertialForceBeam
block = 0
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
eta = 0.1
area = 0.01
Iy = 1e-4
Iz = 1e-4
Ay = 0.0
Az = 0.0
component = 1
variable = disp_y
alpha = -0.3
[../]
[./inertial_force_z]
type = InertialForceBeam
block = 0
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
eta = 0.1
area = 0.01
Iy = 1e-4
Iz = 1e-4
Ay = 0.0
Az = 0.0
component = 2
variable = disp_z
alpha = -0.3
[../]
[./inertial_force_rot_x]
type = InertialForceBeam
block = 0
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
eta = 0.1
area = 0.01
Iy = 1e-4
Iz = 1e-4
Ay = 0.0
Az = 0.0
component = 3
variable = rot_x
alpha = -0.3
[../]
[./inertial_force_rot_y]
type = InertialForceBeam
block = 0
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
eta = 0.1
area = 0.01
Iy = 1e-4
Iz = 1e-4
Ay = 0.0
Az = 0.0
component = 4
variable = rot_y
alpha = -0.3
[../]
[./inertial_force_rot_z]
type = InertialForceBeam
block = 0
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
eta = 0.1
area = 0.01
Iy = 1e-4
Iz = 1e-4
Ay = 0.0
Az = 0.0
component = 5
variable = rot_z
alpha = -0.3
[../]
[]
[Materials]
[./elasticity]
type = ComputeElasticityBeam
youngs_modulus = 1.0e4
poissons_ratio = -0.999875
shear_coefficient = 1.0
block = 0
[../]
[./strain]
type = ComputeIncrementalBeamStrain
block = '0'
displacements = 'disp_x disp_y disp_z'
rotations = 'rot_x rot_y rot_z'
area = 0.01
Ay = 0.0
Az = 0.0
Iy = 1.0e-4
Iz = 1.0e-4
y_orientation = '0.0 1.0 0.0'
[../]
[./stress]
type = ComputeBeamResultants
block = 0
[../]
[./density]
type = GenericConstantMaterial
block = 0
prop_names = 'density'
prop_values = '1.0'
[../]
[]
[Postprocessors]
[./disp_x]
type = PointValue
point = '4.0 0.0 0.0'
variable = disp_x
[../]
[./disp_y]
type = PointValue
point = '4.0 0.0 0.0'
variable = disp_y
[../]
[./vel_y]
type = PointValue
point = '4.0 0.0 0.0'
variable = vel_y
[../]
[./accel_y]
type = PointValue
point = '4.0 0.0 0.0'
variable = accel_y
[../]
[]
[Outputs]
file_base = 'dyn_euler_small_rayleigh_hht_out'
exodus = true
csv = true
perf_graph = true
[]
modules/tensor_mechanics/test/tests/dynamics/wave_1D/wave_rayleigh_hht_ti.i
# Wave propogation in 1D using HHT time integration in the presence of Rayleigh damping
#
# The test is for an 1D bar element of length 4m fixed on one end
# with a sinusoidal pulse dirichlet boundary condition applied to the other end.
# alpha, beta and gamma are HHT time integration parameters
# eta and zeta are mass dependent and stiffness dependent Rayleigh damping
# coefficients, respectively.
# The equation of motion in terms of matrices is:
#
# M*accel + (eta*M+zeta*K)*((1+alpha)*vel-alpha*vel_old)
# +(1+alpha)*K*disp-alpha*K*disp_old = 0
#
# Here M is the mass matrix, K is the stiffness matrix
#
# The displacement at the first, second, third and fourth node at t = 0.1 are
# -7.787499960311491942e-02, 1.955566679096475483e-02 and -4.634888180231294501e-03, respectively.
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 4
nz = 1
xmin = 0.0
xmax = 0.1
ymin = 0.0
ymax = 4.0
zmin = 0.0
zmax = 0.1
[]
[Variables]
[./disp_x]
[../]
[./disp_y]
[../]
[./disp_z]
[../]
[]
[AuxVariables]
[./vel_x]
[../]
[./accel_x]
[../]
[./vel_y]
[../]
[./accel_y]
[../]
[./vel_z]
[../]
[./accel_z]
[../]
[]
[Kernels]
[./DynamicTensorMechanics]
displacements = 'disp_x disp_y disp_z'
alpha = -0.3
zeta = 0.1
[../]
[./inertia_x]
type = InertialForce
variable = disp_x
eta=0.1
alpha = -0.3
[../]
[./inertia_y]
type = InertialForce
variable = disp_y
eta=0.1
alpha = -0.3
[../]
[./inertia_z]
type = InertialForce
variable = disp_z
eta = 0.1
alpha = -0.3
[../]
[]
[AuxKernels]
[./accel_x] # These auxkernels are only to check output
type = TestNewmarkTI
displacement = disp_x
variable = accel_x
first = false
[../]
[./accel_y]
type = TestNewmarkTI
displacement = disp_y
variable = accel_y
first = false
[../]
[./accel_z]
type = TestNewmarkTI
displacement = disp_z
variable = accel_z
first = false
[../]
[./vel_x]
type = TestNewmarkTI
displacement = disp_x
variable = vel_x
[../]
[./vel_y]
type = TestNewmarkTI
displacement = disp_y
variable = vel_y
[../]
[./vel_z]
type = TestNewmarkTI
displacement = disp_z
variable = vel_z
[../]
[]
[BCs]
[./top_y]
type = DirichletBC
variable = disp_y
boundary = top
value=0.0
[../]
[./top_x]
type = DirichletBC
variable = disp_x
boundary = top
value=0.0
[../]
[./top_z]
type = DirichletBC
variable = disp_z
boundary = top
value=0.0
[../]
[./right_x]
type = DirichletBC
variable = disp_x
boundary = right
value=0.0
[../]
[./right_z]
type = DirichletBC
variable = disp_z
boundary = right
value=0.0
[../]
[./left_x]
type = DirichletBC
variable = disp_x
boundary = left
value=0.0
[../]
[./left_z]
type = DirichletBC
variable = disp_z
boundary = left
value=0.0
[../]
[./front_x]
type = DirichletBC
variable = disp_x
boundary = front
value=0.0
[../]
[./front_z]
type = DirichletBC
variable = disp_z
boundary = front
value=0.0
[../]
[./back_x]
type = DirichletBC
variable = disp_x
boundary = back
value=0.0
[../]
[./back_z]
type = DirichletBC
variable = disp_z
boundary = back
value=0.0
[../]
[./bottom_x]
type = DirichletBC
variable = disp_x
boundary = bottom
value=0.0
[../]
[./bottom_z]
type = DirichletBC
variable = disp_z
boundary = bottom
value=0.0
[../]
[./bottom_y]
type = FunctionDirichletBC
variable = disp_y
boundary = bottom
function = displacement_bc
[../]
[]
[Materials]
[./Elasticity_tensor]
type = ComputeElasticityTensor
block = 0
fill_method = symmetric_isotropic
C_ijkl = '1 0'
[../]
[./strain]
type = ComputeSmallStrain
block = 0
displacements = 'disp_x disp_y disp_z'
[../]
[./stress]
type = ComputeLinearElasticStress
block = 0
[../]
[./density]
type = GenericConstantMaterial
block = 0
prop_names = 'density'
prop_values = '1'
[../]
[]
[Executioner]
type = Transient
start_time = 0
end_time = 6.0
l_tol = 1e-12
nl_rel_tol = 1e-12
dt = 0.1
[./TimeIntegrator]
type = NewmarkBeta
beta = 0.422
gamma = 0.8
[../]
[]
[Functions]
[./displacement_bc]
type = PiecewiseLinear
data_file = 'sine_wave.csv'
format = columns
[../]
[]
[Postprocessors]
[./_dt]
type = TimestepSize
[../]
[./disp_1]
type = NodalVariableValue
nodeid = 1
variable = disp_y
[../]
[./disp_2]
type = NodalVariableValue
nodeid = 3
variable = disp_y
[../]
[./disp_3]
type = NodalVariableValue
nodeid = 10
variable = disp_y
[../]
[./disp_4]
type = NodalVariableValue
nodeid = 14
variable = disp_y
[../]
[]
[Outputs]
file_base = 'wave_rayleigh_hht_out'
exodus = true
perf_graph = true
[]
test/tests/kernels/vector_dot_dot/vector_test.i
# Tests calculation of first and second time derivative
# of a coupled vector variable in a material
# a_vec(x,y,z,t) = [t*(t*x + y), t*y, 0]
# a_vec_dot(x,y,z,t) = [2*t*x + y, y, 0]
# a_vec_dot_dot(x,y,z,t) = [2*x, 0, 0]
#
# IMPORTANT NOTE:
# Currently, this test produces a_vec_dot and a_vec_dot_dot that contains oscillations over time.
# This is a known by-product of Newmark Beta time integration (see the Newmark Beta documentation),
# but as of Summer 2019, there is no alternative time integrator in MOOSE that can dampen these
# oscillations. This test is used as coverage for the function call coupledVectorDotDot.
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 4
ymin = 0
ymax = 4
nx = 8
ny = 8
[]
[Functions]
[a_fn]
type = ParsedVectorFunction
value_x = 't * (t * x + y)'
value_y = 't * y'
value_z = 0
[]
[]
[AuxVariables]
[a]
family = LAGRANGE_VEC
order = FIRST
[]
[]
[AuxKernels]
[a_ak]
type = VectorFunctionAux
variable = a
function = a_fn
[]
[]
[Materials]
[cm]
type = VectorCoupledValuesMaterial
variable = a
[]
[]
[Variables]
[u] # u is zero
family = LAGRANGE_VEC
order = FIRST
[]
[]
[Kernels]
[td]
type = VectorTimeDerivative
variable = u
[]
[]
[Executioner]
type = Transient
dt = 1
num_steps = 3
[TimeIntegrator]
type = NewmarkBeta
[]
[]
[Outputs]
[./out]
type = Exodus
output_material_properties = true
show_material_properties = 'a_value a_dot a_dot_dot a_dot_du a_dot_dot_du'
[../]
[]
test/tests/variables/time_derivatives_neighbor/test.i
[Mesh]
type = GeneratedMesh
dim = 1
xmin = 0
xmax = 4
nx = 2
[]
[Functions]
[a_fn]
type = ParsedFunction
value = 't*(t+x)'
[]
[]
[AuxVariables]
[a]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[a_ak]
type = FunctionAux
variable = a
function = a_fn
[]
[]
[Materials]
[cm]
type = CoupledValuesMaterial
variable = a
[]
[]
[Variables]
[u]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[td]
type = TimeDerivative
variable = u
[]
[]
[DGKernels]
[dgk]
type = MatDGKernel
variable = u
mat_prop = a_value
[]
[]
[Executioner]
type = Transient
dt = 1
num_steps = 3
[TimeIntegrator]
type = NewmarkBeta
[]
[Quadrature]
type = GAUSS
order = FIRST
[]
[]
[Outputs]
[./out]
type = Exodus
output_material_properties = true
show_material_properties = 'a_value a_dot a_dot_dot a_dot_du a_dot_dot_du'
[../]
[]
test/tests/time_integrators/newmark-beta/newmark_beta_prescribed_parameters.i
###########################################################
# This is a simple test with a time-dependent problem
# demonstrating the use of the TimeIntegrator system.
#
# Testing that the first and second time derivatives
# are calculated correctly using the Newmark-Beta method
#
# @Requirement F1.30
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 1
ny = 1
[]
[Variables]
[u]
[]
[]
[Functions]
[forcing_fn]
type = PiecewiseLinear
x = '0.0 0.1 0.2 0.3 0.4 0.5 0.6'
y = '0.0 0.0 0.0025 0.01 0.0175 0.02 0.02'
[]
[]
[Kernels]
[ie]
type = TimeDerivative
variable = u
[]
[diff]
type = Diffusion
variable = u
[]
[]
[BCs]
[left]
type = FunctionDirichletBC
variable = u
boundary = 'left'
function = forcing_fn
[]
[right]
type = FunctionDirichletBC
variable = u
boundary = 'right'
function = forcing_fn
[]
[]
[Executioner]
type = Transient
start_time = 0.0
num_steps = 6
dt = 0.1
[TimeIntegrator]
type = NewmarkBeta
beta = 0.4225
gamma = 0.8
[]
[]
[Postprocessors]
[udot]
type = ElementAverageTimeDerivative
variable = u
[]
[udotdot]
type = ElementAverageSecondTimeDerivative
variable = u
[]
[u]
type = ElementAverageValue
variable = u
[]
[]
[Outputs]
csv = true
[]
References
- K. J. Bathe and G. Noh.
Insight into an implicit time integration scheme for structural dynamics.
Computers and Structures, 98-99:1–6, 2012.[BibTeX]
@article{bathe2012insight, author = "Bathe, K. J. and Noh, G.", title = "Insight into an implicit time integration scheme for structural dynamics", journal = "Computers and Structures", volume = "98-99", pages = "1-6", year = "2012" }
- N. M. Newmark.
A method of computation for structural dynamics.
Journal of Engineering Mechanics, 85(EM3):67–94, 1959.[BibTeX]
@article{newmark1959amethod, author = "Newmark, N. M.", title = "A method of computation for structural dynamics", journal = "Journal of Engineering Mechanics", publisher = "ASCE", volume = "85", number = "EM3", pages = "67--94", year = "1959" }