dingetal2012

Class summary

DingEtAl2012([BC, nterms, E, I, rho, A, L, …]) Finite elastic Euler-Bernoulli beam resting on non-linear viscoelastic foundation subjected to a moving load.

Function summary

dingetal_figure_8() Reproduce Ding Et Al 2012 Figure 8 (might take a while).

Module listing

Ding et al (2012) “Convergence of Galerkin truncation for dynamic response of finite beams on nonlinear foundations under a moving load”.

class geotecha.beam_on_foundation.dingetal2012.DingEtAl2012(BC='SS', nterms=50, E=None, I=None, rho=None, A=None, L=None, k1=None, k3=None, mu=None, Fz=None, v=None, v_norm=None, kf=None, Fz_norm=None, mu_norm=None, k1_norm=None, k3_norm=None, nquad=30)[source]

Bases: object

Finite elastic Euler-Bernoulli beam resting on non-linear viscoelastic foundation subjected to a moving load.

An implementation of Ding et al. (2012) [1].

You don’t need all the parameters. Basically if normalised values are equal None (i.e. default) then those properties will be calculated from the non normalised quantities. All calculations are done in normalised space. You only need the non-normalised variables if you want non-normalised output.

Parameters:
BC : [“SS”, “CC”, “FF”], optional

Boundary condition. Simply Supported (SS), Clamped Clamped (“CC”), Free Free (FF). Default BC=”SS”.

nterms : integer, optional

Number of terms for Galerkin truncation. Default nterms=50

E : float, optional

Young’s modulus of beam.(E in [1]).

rho : float, optional

Mass density of beam.

I : float, optional

Second moment of area of beam (I in [1]).

A : float, optional

Cross sectional area of beam

L : float, optional

Length of beam. (L in [1]).

k1 : float, optional

Mean stiffness of foundation.

k3 : float, optional

Non linear stiffness of foundation.

mu : float, optional

Viscous damping of foundation.

Fz : float, optional

Load.

v : float, optional

Speed of the moving force.

v_norm : float, optional

normalised velocity = V * sqrt(rho / E). An example of a consistent set of units to get the correct v_norm is rho in kg/m^3, L in m, and E in Pa.

kf : float, optional

Normalised modulus of elasticity = 1 / L * sqrt(I / A).

Fz_norm : float, optional

Normalised load = Fz / (E * A).

mu_norm : float, optional

Normalised damping = mu / A * sqrt(L**2 / (rho * E)).

k1_norm : float, optional

Normalised mean stiffness = k1 * L**2 / (E * A).

k3_norm : float, optional

Normalised non-linear stiffness = k3 * L**4 / (E * A)

nquad : integer, optional

Number of quadrature points for numerical integration of non-linear k3*w**3*w_ term. Default nquad=30. Note I’ve had errors when n>35. For the special case of nquad=None then integration will be performed using scipy.integrate.quad; this is slower.

References

[1](1, 2, 3, 4, 5, 6, 7, 8) Ding, H., Chen, L.-Q., and Yang, S.-P. (2012). “Convergence of Galerkin truncation for dynamic response of finite beams on nonlinear foundations under a moving load.” Journal of Sound and Vibration, 331(10), 2426-2442.
Attributes:
phi : function

Relevant Galerkin trial function. Depends on BC. See [1] for details.

beta : 1d ndarray of nterms float

beta terms in Galerkin trial function.

xj : 1d ndarray of nquad float

Quadrature integration points.

Ij : 1d ndarray of nquad float

Weighting coefficients for numerical integration.

BC_coeff : int

Coefficent to multiply the Fz and k3 terms in the ode. For `BC`=”SS” BC_coeff=2, for `BC`=”CC” or “FF” BC_coeff=2. See [1] Equation (31) and (32).

t : 1d ndarray of float

Raw time values. Only valid after running calulate_qk.

t_norm : 1d ndarray of float

Normlised time values = t / L * sqrt(E / rho). Only valid after running calulate_qk.

qsol : 2d ndarray of shape(len(t), 2* nterms) float

Values of Galerkin coefficients at each time i.e. qk(t) in [1]. w(x) = sum(qk * phi_k(x)). Only valid after running calulate_qk.

Methods

calulate_qk([t, t_norm]) Calculate the nterm Galerkin coefficients qk at each time value
vectorfield(vsv, tnow[, p])
Parameters:
w(qk, x) Nomalised vertcal deformation at x
w_cubed_wi(x, q, i) non-linear cube term for numerical integration
wofx([x, x_norm, tslice, normalise_w]) Deflection at distance x, and times t[tslice]
phiCC  
phiFF  
phiSS  
calulate_qk(t=None, t_norm=None, **odeint_kwargs)[source]

Calculate the nterm Galerkin coefficients qk at each time value

Parameters:
t : float or array of float, optional

Raw time values.

t_norm : float or array of float, optional

Normalised time values. If t_norm==None then it will be calculated from raw t values and other params t_norm = t / L * sqrt(E / rho).

Notes

This method determines initializes self.t and self.t_norm and calculates self.qsol.

phiCC(x, beta)[source]
phiFF(x, beta)[source]
phiSS(x, beta)[source]
vectorfield(vsv, tnow, p=())[source]
Parameters:
vsv : float

Vector of the state variables. vsv = [q1, q2, …qk, q1dot, q2dot, …, qkdot] where qk is the kth galerkin coefficient and qkdot is the time derivative of the kth Galerkin coefficient.

tnow : float

Current time.

p : various

Vector of parameters

Returns:
vsvdot : vector of state variables first derivatives

vsvdot = [q1dot, q2dot, …qkdot, q1dotdot, q2dotdot, …, qkdotdot]

w(qk, x)[source]

Nomalised vertcal deformation at x

Parameters:
qk : 1d ndarray of nterm floats

Galerkin coefficients.

x : float

nomalised distances to calculate deflection at.

Returns:
w : float

vertical deformation at x value.

w_cubed_wi(x, q, i)[source]

non-linear cube term for numerical integration

wofx(x=None, x_norm=None, tslice=slice(None, None, None), normalise_w=True)[source]

Deflection at distance x, and times t[tslice]

Parameters:
x : float or ndarray of float, optional

Raw values of x to calculate deflection at. Default x=None.

x_norm : float or array of float, optional

Normalised x values to calc deflection at. If x_norm==None then it will be calculated frm x and other properties : x_norm = x / L. Default x_norm=None.

tslice : slice object, optional

slice to select subset of time values. Default tslice=slice(None, None, None) i.e. all time values. Note the array of time values is already in the object (it was used to calc the qk galerkin coefficients).

normalise_w : True/False, optional

If True then output is normalised deflection. Default nomalise_w=True.

Returns:
w : array of float

Deflections at x and self.t_norm[tslice]

geotecha.beam_on_foundation.dingetal2012.dingetal_figure_8()[source]

Reproduce Ding Et Al 2012 Figure 8 (might take a while).

Note that a plot will be be saved to disk in current working directory as well as a timing file.