specbeam

Class summary

MovingPointLoads(x, p[, offset, v, L, t0, x0]) Multi-axle vehicle represented by point loads in a line
SpecBeam([BC, nterms, E, I, rho, A, L, k1, …]) Finite elastic Euler-Bernoulli beam resting on viscoelastic foundation subjected to a moving load, piecewise-linear properties.

Function summary

DAFinfinite(alp, bet[, kratio]) Dynamic amplification factor for infinite beam on winkler foundation
DAFrelative(alp, bet, k2_k1) Relative Dynamic Amplification Factor, based on ratio of beam stiffness to some reference beam stiffness.
FIGURE_DAF_constant_prop([saveas, nterms, …]) Dynamic amplification plot, esveld fig 6.18.
FIGURE_DAF_transition_length([saveas, …]) Dynamic amplification factor length of transition
FIGURE_DAF_transition_length_esveld([…]) Dynamic amplification factor length of transition, reproduce esveld figure
FIGURE_DAFwrtinf_transition_length([saveas, …]) Dynamic amplification factor with rest to length of transition
FIGURE_DAFwrtinf_transition_length_betaconst([…]) Dynamic amplification factor with rest to length of transition with beta term constant
FIGURE_DAFwrtinf_transition_lengthrev([…]) Dynamic amplification factor with rest to length of transition
FIGURE_convergence_step_change([v_nterms, …]) Check convergence with step change in stiffness at midpoint
FIGURE_defl_envelope_vary_velocity([…]) Deflection envelopes with different velocities
FIGURE_verify_against_DingEtAl2012_fig8([…]) Ding et al Figure 8, displacement vs time at beam midpoint (using the runme method) (50 and 200 terms)
FIG_DAF_envelope_vary_beta([ax, betas, …]) DAF envelope for whole beam with different alpha & beta
SpecBeam_const_mat_midpoint_defl_runme_analytical_play() Test SpecBeam for constant mat: close to Ding et al Figure 8, displacement vs time at beam midpoint (using the runme method) but with k3=0
SpecBeam_stationary_load_const_mat_midpoint_defl_runme_analytical_play() Test SpecBeam for constant mat: close to Ding et al Figure 8, displacement vs time at beam midpoint (using the runme method) but with k3=0
align_yaxis(ax1, ax2) Adjust y-axis limits so zeros of the two axes align, zooming them out by same ratio.
article_figure_01() Figure 1 Mid point deflection comparison with Ding et al 2012, figure from Walker and Indraratna (in press).
article_figure_02() Figure 2 Deflection amplification factor compared with analytical infinite beam solution, figure from Walker and Indraratna (in press).
article_figure_03() Figure 3 Low damping DAF for velocity ratio alp=0p5, figure from Walker and Indraratna (in press).
article_figure_04() Figure 4 Relative DAF for various stiffness ratios, figure from Walker and Indraratna (in press).
article_figure_05a_05b_05c_07() Figure 5 a),b),&c) Max deflection envelope for bet=0p1 k2_k1=10, and Figure 7 Deflection vs time at point which experiences maximum deflection alp=0p8, bet=0p1, k2_k1=10, soft to stiff=forward, figure from Walker and Indraratna (in press).
article_figure_06() Additional deflection for transition effects, bet=0p1, alp=10, (dashed=stiff to soft, solid=soft to stiff), figure from Walker and Indraratna (in press).
article_figure_08() Figure 8 Deflection vs time under four moving loads alp=0p8, bet=0p1, k2/k1=10, soft to stiff, figure from Walker and Indraratna (in press).
article_figure_09() Figure 9 is drawing in a word file, this routine is just a marker to indicated that figure 9 has not been forgotten, figure from Walker and Indraratna (in press).
article_figure_10_11() Figure 10 Measured vs modelled deflections at various locations in the transition zone, Figure 11 Peak displacement and stiffness distribution, figure from Walker and Indraratna (in press).
article_figure_12() Figure 12 python code for relative dynamic amplification factor calculation, figure from Walker and Indraratna (in press).
case_study_PaixaoEtAl2014([saveas, nterms, …]) Try and replicate Paixao et al 2014 for 6 rail cars going over transition using specbeam.
dingetal_figure_8([v_nterms]) Reproduce Ding Et Al 2012 Figure 8 (might take a while).
esveld_fig6p17() Displacement shapes before and after a moving point load on infinite beam on elastic foundation for a sub crtical , critical, super critcal velocity ratio (alpha) and damping ratios (beta) , Esveld Figure 6p17.
esveld_fig6p18() Dynamic Amplification Factor vs velocity ratio (alpha) for variuos damping ratios (beta) for moving point load on infinite beam on elastic foundation, Esveld Figure 6p18.
moving_load_on_beam_on_elastic_foundation(…) Displacement vs distance from moving point load on infinite beam on elastic foundation.
multi_static(xc, frel[, wsingle, xlim, n, L]) Deflection shape for multiple static point loads
point_load_static_defl_shape_function(x, L) Deflection shape function for static point load on beam on foundation.
reverse_curve_transition(x, x1, y1, x2, y2) Transition curve between (x1, y1) and (x2, y2)
specbeam_vs_analytical_single_static_load() Compare specbeam vs analytical solution for single static load at mid point of ‘long’ beam.
static_point_load_shape(x, L) Deflection shape parameter for point load at x=0 on infinite beam on elastic foundation with characteristic length L.
tanhtran(x, a, b) Transition between two functions using tanh
tanhtrani(y, a, b) Inverse of Transition between two functions using tanh
tanhtranib(x, y, a) Smoothing factor of tansition curve to get known (x,y)
transition_DAF1([saveas, ax, kratios, …]) Plot of Dynamic Amplification Factor (DAF) vs velocity ratio for different stiffness and damping values.
transition_DAF1_at_specific_alp([saveas, …]) Plot of Dynamic Amplification Factor (DAF) vs stiffness ratio for different velocity ratios and damping values.
transition_zones1([Lt_Lcs, alphas, …]) Exploration of deflection envelopes caused by point loads travelling through a stiffness transition zone on beam on visco-elastic foundation using Specbeam.
transition_zones2([animateme, saveas, xlim, …]) Essentially a duplicate of transition_zones1 to get a DAF* vs Lt/Lc plot for forward and reverse plots

Module listing

Finite elastic Euler-Bernoulli beam resting on viscoelastic foundation with piecewise-linear properties, subjected to multiple moving point loads.

geotecha.beam_on_foundation.specbeam.DAFinfinite(alp, bet, kratio=1.0)[source]

Dynamic amplification factor for infinite beam on winkler foundation

Parameters:
alp : float

Velocity ratio. alp = v/vcr where critical velocity vcr=sqrt(2*sqrt(k*E*I)/(rho*A))**0.5

bet : float

Damping ratio. bet=c/ccr where critical damping ccr=2*sqrt(rho*A*k). The damping ratio is with respect to the reference foundation stiffness.

kratio : float, optional

Ratio of foundation stiffness to reference value i.e. k2/k1 where k1 is the reference value. Default kratio=1

Returns:
DAF : float

Dynamic/deflection amplification factor with respect to the reference foundation stiffness, k1. i.e. deflection = DAF*[static deflection of point load on beam on k1 ]

Notes
—–
Only works for scalar input
geotecha.beam_on_foundation.specbeam.DAFrelative(alp, bet, k2_k1)[source]

Relative Dynamic Amplification Factor, based on ratio of beam stiffness to some reference beam stiffness.

Parameters:
alp : float

velocity ratio. alp=v/vcr; vcr= SQRT(2*SQRT(k*EI)/(rho*A))

bet : float

damping ratio. bet=c/ccr; ccr=2*SQRT(rho*A*k)

k2_k1 : float

stiffness ratio of beam to current beam

Returns:
DAF : float

Relative Dynamic Amplification Factor

Notes

DAFrelative is a terse implementation for an appendix of code in Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.

geotecha.beam_on_foundation.specbeam.FIGURE_DAF_constant_prop(saveas=None, nterms=100, force_calc=True, prefix='', alphas=None, betas=None, nt=20, nx=55, end_damp=0, xwindow=(0.3, 0.7), xeval=(0.3, 0.7), numerical_DAF=True, article_formatting=False)[source]

Dynamic amplification plot, esveld fig 6.18. p120

geotecha.beam_on_foundation.specbeam.FIGURE_DAF_transition_length(saveas=None, k_ratio=4, beta=0.1, nterms=100, force_calc=True)[source]

Dynamic amplification factor length of transition

geotecha.beam_on_foundation.specbeam.FIGURE_DAF_transition_length_esveld(saveas=None, k_ratio=4, beta=0.1, nterms=100, force_calc=True)[source]

Dynamic amplification factor length of transition, reproduce esveld figure

geotecha.beam_on_foundation.specbeam.FIGURE_DAFwrtinf_transition_length(saveas=None, k_ratio=4, beta=0.1, nterms=100, force_calc=True)[source]

Dynamic amplification factor with rest to length of transition

geotecha.beam_on_foundation.specbeam.FIGURE_DAFwrtinf_transition_length_betaconst(saveas=None, k_ratio=4, beta=0.1, nterms=100, force_calc=True)[source]

Dynamic amplification factor with rest to length of transition with beta term constant

geotecha.beam_on_foundation.specbeam.FIGURE_DAFwrtinf_transition_lengthrev(saveas=None, k_ratio=4, beta=0.1, nterms=100, force_calc=True)[source]

Dynamic amplification factor with rest to length of transition

geotecha.beam_on_foundation.specbeam.FIGURE_convergence_step_change(v_nterms=(50, 75, 150, 200), saveas=None)[source]

Check convergence with step change in stiffness at midpoint

geotecha.beam_on_foundation.specbeam.FIGURE_defl_envelope_vary_velocity(velocities_norm=(0.01165, ), saveas=None, nterms=100, force_calc=True, nx=400, nt=400)[source]

Deflection envelopes with different velocities

geotecha.beam_on_foundation.specbeam.FIGURE_verify_against_DingEtAl2012_fig8(v_nterms=(50, 75, 150, 200), saveas=None)[source]

Ding et al Figure 8, displacement vs time at beam midpoint (using the runme method) (50 and 200 terms)

geotecha.beam_on_foundation.specbeam.FIG_DAF_envelope_vary_beta(ax=None, betas=None, alphas=None, saveas=None, nx=200, nt=200, nterms=100, force_calc=True, t_extend=1, ylim=None, numerical_DAF=False, article_formatting=False)[source]

DAF envelope for whole beam with different alpha & beta

A plot will be produced of DAF vs normalised distance along beam, where the dynamic deflection in the DAF calculation is the maximum experienced each x/L point.

Parameters:
ax : matplotlib.Axes object, optional

Axes object to plot on. Default ax=None i.e. fig and axis will be created.

betas : list of float, optional

list of damping ratios, Default betas=None which will use [0, 0.05, 0.1].

alphas : list of float, optional

list of velocity ratios, Default alphas=None which will use [0.7].

saveas : string, optional

Filename stem (not including extension) to save figure to. Default saveas=None, i.e. figure won’t be saved to disk.

nx : int, optional

Number of points between xlim[0] and xlim[1] at which to calculate deflection. Default nx = 200.

nt : int, optional

Number of evaluation times from when load enters xlim[0] and leaves xlim[1]. Default nt=200

nterms : int, optional

Number of series terms in approximation. More terms leads to higher accuracy but also longer computation times. Start of with a small enough value to observe behaviour and then increase to confirm. Default nterms=100 .

force_calc : [True, False], optional

When True all calcuations will be done from scratch. When False Current working directory will be searched for a results file named according to “saveas” and if found it will be loaded. Be careful as name of results file may not be unique to your analysis parameters. Safest to use force_calc = True, but if you are simply regenerating graphs then consider force_calc=False. Default force_calc=True.

t_extend : float, optional

factor to extend maximum evalutated time. As per nt’ evaluation times are normally between first point entering `xlim[0]’ and leaving `xlim[1]. However, with mulitple loads, loads can still be on the beam at that point. max eval time can be extended fby a factor. t_extend = 1 .

ylim : float, optional

Y axis limit. default ylim=None i.e. auto axis scale

numerical_DAF : [False, True], optional

If True then the static deflection for use in Dynamic amplification factor calcuations (DAF = w/wstatic) willbe determined from a very slow moving load travelling over the beam/foundation. If False then the analytical expression for wstatic will be used. Default numerical_DAF=False.

article_formatting : [False, True], optional

You should never need this. It was a quick hack to do some custom formatting for a journal article. Dont even try messing with it. will ony work with 3 Lt_Lcs values values or less. Default=False, i.e. no special formatting.

class geotecha.beam_on_foundation.specbeam.MovingPointLoads(x, p, offset=0, v=None, L=None, t0=0, x0=0.0)[source]

Bases: object

Multi-axle vehicle represented by point loads in a line

Parameters:
x : 1d array or list

x-coord of point loads. Vehicle assumed to start at 0.

p : 1d array or list

Magnitude of point loads at corresponding x value.

offset : float, optional

x-coords will be offset by this distance. Positive values will move x to the right (i.e. increase in x).

v : float, optional

Velocity. Default v = None, i.e. you’ll need to enter in various functions.

L : float, optional

Length of Beam. Default L=None i.e. you’ll need to explicitly enter in various functions. You can use a psuedo L with x0 to have the loads begin and end midspan.

t0 : float, optional

Time the point loads begin at x0 on the beam. default t0=0 ie. you’ll have to explicitly enter in various functions.

x0 : float, optional

Start position of loads on beam. Defualt x0=0.0 .

Attributes:
xt : 1d array

transformed x coords after mirror and offset are applied.

start_time, end_time : 1d array

start and end time for each point load traversing length L.

Methods

animateme([x0, t0, tf, v0, xlim, ylim]) Animate moving vehicle
axle_positions([t, v0, t0, x0]) Positions of axles at given time
convert_to_specbeam([v, L, t0, x0]) Convert vehicle at velocity v into surcharge_vs_time and surcharge_omega_phase PolyLines for use with spectral method.
plot_specbeam_PolyLines([ax, v, L, t0, x0]) Plot load_vs_time for each vehicle axis.
plotme([ax, t, v0, t0, x0, xlim, ylim]) Plot the load at time t
point_loads_on_beam(t) Position and magnitude of point loads between [0, L] at time t
setv(v) Set the vehicle’s velocity, then calc start_time and end_time when vehicle’s axles enter and leave the beam
time_axles_pass_beam_point(x) Time values when each axle passes a particular point on a beam
time_axles_pass_point(x[, v0, t0, x0]) time values when each moving load passes a certain point
delv  
getv  
animateme(x0=0, t0=0, tf=1, v0=1, xlim=None, ylim=None)[source]

Animate moving vehicle

Parameters:
x0 : float, optional

Start position of vehicle, default=0.

t0 : float, optional

Start time of animation. Default t0=0.

tf : float, optional

End time of animation. Default tf=1.

v0 : float, optional

Initial velocity of vehicle. Default v0=1. If v0<0 then positions will be multiplied by -1 before adjusting position.

xlim, ylim : tuple of size 2, default xlim=ylim=None.

axes limits.

Returns:
ani : matplotlib animation

Animation of moving vehicle

axle_positions(t=0, v0=0, t0=0, x0=0)[source]

Positions of axles at given time

Parameters:
t : float

Time of interest

v0 : float, optional

Initial speed of vehicle. Default v0=0. If v0<0 then axle coords will be multiplied by -1 before moving.

t0 : float, optional

Time when vehicle starts moving at velocity v0. Default t0=0

x0 : float, optional

Start postion of vehicle

Returns:
xnew : 1d array

xcoords of axles at time t.

convert_to_specbeam(v=None, L=None, t0=None, x0=None)[source]

Convert vehicle at velocity v into surcharge_vs_time and surcharge_omega_phase PolyLines for use with spectral method.

Parameters:
v : float, optional

Velocity of vehicle. Default v=None in which case the v used during object initialization will be uses

L : float, optional

Length of beam. Default L=None in which case the L used during object initialization will be uses

t0 : float, optional

time when vehicle reference point begins entering the beam. Default t0=0. If v>0 vehicle will enter from the right. If v<0 vehicle will enter from the right. v=0 will result in an error.

x0 : float, optional

Position where loads start. Default x0=none, i.e. value used during initialization will be used.

Returns:
plines : list of PolyLine
omega_phase : list of two element tuples

frequency and phase angle for use in spec beam. as per the notes above will need to multiply both omega and phase by M for use in specbeam.

Notes

Vehicle is made up of axles at a relative distance dxi from a reference point. dxi should be negative. If velocity v is positive then the reference point enters the beam (at x=x0) at time t0. The position of the ith axle is:

xi = x0 + dxi + v * (t - t0)
The ith axle enters the beam at position x0 at time
tstarti = t0 - dxi / v
and leaves the beam (x = x0 + L) at time
tendi = t0+(L-dxi)/v

If v is negative then the reference point enters the beam from the right, then we first flip the vehicle by making dx_i=-dxi. The positiion and times entering and leaving the beam are

xi = x0 + L + dx_i + v * (t - t0)
and
tstart = t0 - dx_i / v tend = t0 - (L + dx_i)/v

For use in spec beam we will need sin(M*xi) in a form sin(omega*t+phase).

for v>0 xi = v * t + (x0 + dxi - v * t0) so omega/M = v and phase/M = x0 + dxi-v*t0

for v<0 xi = v * t + (x0 + L + dx_i - v * t0) so omega/M = v and phase/M = x0 + L + dx_i - v*t0

delv()[source]
getv()[source]
plot_specbeam_PolyLines(ax=None, v=1, L=1, t0=0, x0=0)[source]

Plot load_vs_time for each vehicle axis.

Parameters:
ax : Matplotlib.Axes, optional

Axes to plot on. Default ax=None i.e. a new fig and ax will be created.

v : float, optional

velocity, default v=1.

L : float, optional

Lenght of beam default L=1

x0, t0 : float, optional

Vehicle starts at positon x0 at time t0. Default x0=0, t0=0

plotme(ax=None, t=0, v0=0, t0=0, x0=0, xlim=None, ylim=None)[source]

Plot the load at time t

Each axle load is represented by an arrow, with head at 0 and tail at the load magnitude.

Parameters:
ax : pyplot.Axes, optional

Axes object to plot on. Default ax=None i.e. new figure and axes will be created.

t : float, optional

Time at which to plot the vehicle positoin. Default t=0.

v0 : float, optional

Initial velocity of vehicle. default v0=0. If v0<0 then positions will be multiplied by -1 before adjusting position.

t0 : float, optional

Time the vehicle starts moving. Default t0=0.

x0 : float, optional

Start position of vehicle

xlim, ylim : tuple of size 2, default xlim=ylim=None.

axes limits.

Returns:
ax : pyplot.Axes opbect

Axes object with plotted vehicle position.

point_loads_on_beam(t)[source]

Position and magnitude of point loads between [0, L] at time t

Parameters:
t : float

Current time.

Returns:
x : 1d array of float

position of each load at time t.

p : 1d array of float

Load at time t. Will be set to zero if load is not beteen 0 and L.

Notes

Position can be offset by self.x0

setv(v)[source]

Set the vehicle’s velocity, then calc start_time and end_time when vehicle’s axles enter and leave the beam

time_axles_pass_beam_point(x)[source]

Time values when each axle passes a particular point on a beam

Parameters:
x : float

x coord of interest

Returns:
t_pass : 1d aray of float

time values when each axle passes the particualr point

Notes

If self.v < 0 then loads enter beam from the left, i.e. at time t0 axles will be at L+x0+xt

Not really sure if I’ve accountered for self.x0!=0 for the negative velocity. Just make self.x0 =0 for simplicity

time_axles_pass_point(x, v0=0, t0=0, x0=0)[source]

time values when each moving load passes a certain point

Parameters:
x : float

x-coord of interest

v0 : float, optional

Initial speed of vehicle. Default v0=0. If v0<0 then axle coords will be multiplied by -1 before moving.

t0 : float, optional

Time when vehicle starts moving at velocity v0. Default t0=0

x0 : float, optional

Start postion of vehicle

L : float, optional

Length of beam. Default L=None in which case the L used during object initialization will be uses.

Returns:
t_pass : 1d array of float

list of times corresponding to when each axle passes the given x value.

v

Velocity

class geotecha.beam_on_foundation.specbeam.SpecBeam(BC='SS', nterms=50, E=None, I=None, rho=None, A=None, L=None, k1=None, k3=None, mu=None, moving_loads_x=None, moving_loads_Fz=None, moving_loads_v=None, moving_loads_offset=None, moving_loads_t0=None, moving_loads_x0=None, moving_loads_L=None, tvals=None, xvals=None, v_norm=None, moving_loads_x_norm=None, moving_loads_Fz_norm=None, moving_loads_v_norm=None, moving_loads_offset_norm=None, moving_loads_t0_norm=None, moving_loads_x0_norm=None, moving_loads_L_norm=None, tvals_norm=None, xvals_norm=None, stationary_loads_x=None, stationary_loads_vs_t=None, stationary_loads_omega_phase=None, stationary_loads_x_norm=None, stationary_loads_vs_t_norm=None, stationary_loads_omega_phase_norm=None, kf=None, mu_norm=None, k1_norm=None, k3_norm=None, nquad=30, Ebar=PolyLine(array([[0., 1.], [1., 1.]])), rhobar=PolyLine(array([[0., 1.], [1., 1.]])), Ibar=PolyLine(array([[0., 1.], [1., 1.]])), Abar=PolyLine(array([[0., 1.], [1., 1.]])), k1bar=PolyLine(array([[0., 1.], [1., 1.]])), k3bar=PolyLine(array([[0., 1.], [1., 1.]])), mubar=PolyLine(array([[0., 1.], [1., 1.]])), implementation='vectorized', use_analytical=False, file_stem='specbeam_', force_calc=False)[source]

Bases: object

Finite elastic Euler-Bernoulli beam resting on

viscoelastic foundation subjected to a moving load, piecewise-linear properties.

An extension of Ding et al. (2012) [Rb51e5289a983-1] with piecewise linear material properties and non-linear foundation stiffness k3=0 in k3*w**3.

You don’t need all the parameters. Basically if normalised values are ‘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.

Note: some of the code is a relic of geotecha.beam_on_foundation.dingetal2012 which is a straight up implementation of [Rb51e5289a983-1] . Sorry for any confusion.

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

Boundary condition. Can only have Simply Supported (SS). Clamped Clamped (CC) and Free Free (FF) are a relic from a differnet program. 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 [Rb51e5289a983-1]).

rho : float, optional

Mass density of beam.

I : float, optional

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

A : float, optional

Cross sectional area of beam

L : float, optional

Length of beam. (L in [Rb51e5289a983-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. This is a relic; loads now specified with moving_loads_…
v : list of float, optional

Velocity of each group of moving point loads.

moving_loads_x : list of lists of floats, optional

Each sublist contains x-coords of each point load in the moving load group. The coords are relative to a reference point which will enter the beam at the corresponding time value in moving_loads_t0. Coords are usually negative for vehicle facing left to right.

moving_loads_Fz : list of lists of floats, optional

Each sublist contains magnitude of of each point load in moving load group.

moving_loads_v : list of float, optional

Velocity of each moving point load group. When velocity is negative the x-coords in the corresponding member of moving_loads_x will be multiplied by -1 before the beam enters the beam from the right.

moving_loads_offset : list of float, optional

An initial distance by which the corresponding moving_loads_x is offset. Will default to zero if not inputted.

moving_loads_t0 : list of lists, optional

Time that moving load group enters the beam (strictyly speaking it is when the 0 coord in moving_loads_x enters the beam). If v>0 moving load group starts at x=0, t=0 and moves to the right. If v<0 then moving load group starts at x=L,t=0 and moves to the left. Defaults to 0.

moving_loads_x0 : list of float, optional

Position where moving load group appears on the beam. Defaults to zero if ommited.

moving_loads_L : list of float, optional

Distance the load group travels on beam. Defaults to L if omitted. Used in conjunction with moving loads_x0.

tvals : 1d array of float, optional

Output times

xvals : 1d array of float, optional

Output x-coord

v_norm : list of float, optional

normalised velocity = V * sqrt(rho / E) for each group of moving point loads. 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.

moving_loads_x_norm : list of lists of floats, optional

As per moving_loads_x but nomalised. x_norm = x / L

moving_loads_Fz_norm : list of lists of floats, optional

As per moving_loads_Fz but nomalised. Fz_norm = Fz / (E * A).

moving_loads_v_norm : list of float, optional

As per moving_loads_v but nomalised. v_norm= v * sqrt(rho / E) for each group of moving point loads. 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.

moving_loads_offset_norm : list of float, optional

As per moving_loads_offset but nomalised. offset_norm = offet / L

moving_loads_t0_norm : list of lists, optional

As per moving_loads_t0 but nomalised. t0_norm = t0 / L * np.sqrt(E / rho)

moving_loads_x0_norm : list of float, optional

As per moving_loads_x0. x0_norm=x0/L.

moving_loads_L_norm : list of float, optional

As per moving_loads_L_norm but normalised.L_norm = L’ / L

tvals_norm : 1d array of float, optional

Normalised output times. t_norm = t / L * np.sqrt(E / rho)

xvals_norm : 1d array of float, optional

Normalised output coords. x_norm = x / L

stationary_loads_x : list of float, optional

Location of stationary point loads.

stationary_loads_vs_t : list of Polyline, optional

Stationary load magnitude variation with time. PolyLine(time, magnitude).

stationary_loads_omega_phase : list of 2 element tuples, optional

(omega, phase) to define cyclic variation of stationary point loads. i.e. mag_vs_time * cos(omega*t + phase). If stationary_loads_omega_phase is None then cyclic component will be ignored. If stationary_loads_omega_phase is a list then if any member is None then cyclic component will not be applied for that load combo.

stationary_loads_x_norm : list of float, optional

As per stationary_loads_x but normalised

stationary_loads_vs_t_norm : list of Polyline, optional

As per stationary_loads_vs_t but normalised

stationary_loads_omega_phase_norm : list of 2 element tuples, optional

As per stationary_loads_omega_phase but normalised.

kf : float, optional

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

# Fz_norm : float, optional
# Normalised load = Fz / (E * A). This is a relic; loads now specified
# with moving_loads_…
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.

Ebar : PolyLine, optional

PolyLine object representing piecewise linear E/Eref vs x. Default Ebar = PolyLine([0],[1],[1],[1]) i.e. uniform.

rhobar : PolyLine, optional

PolyLine object representing piecewise linear rho/rhoref vs x. Default rhobar = PolyLine([0],[1],[1],[1]) i.e. uniform.

Ibar : PolyLine, optional

PolyLine object representing piecewise linear I/Iref vs x. Default Ibar = PolyLine([0],[1],[1],[1]) i.e. uniform.

Abar : PolyLine, optional

PolyLine object representing piecewise linear A/Aref vs x. Default Abar = PolyLine([0],[1],[1],[1]) i.e. uniform.

k1bar : PolyLine, optional

PolyLine object representing piecewise linear k1/k1ref vs x. Default k1bar = PolyLine([0],[1],[1],[1]) i.e. uniform.

k3bar : PolyLine, optional

PolyLine object representing piecewise linear k3/k3ref vs x. Default k3bar = PolyLine([0],[1],[1],[1]) i.e. uniform.

mubar : PolyLine, optional

PolyLine object representing piecewise linear mu/muref vs x. Default mubar = PolyLine([0],[1],[1],[1]) i.e. uniform.

use_analytical : [False,true]

If true then system of second order ode resulting from spectral galerkin method will be solved analytically. If False then the ode system will be solved numerically (potentially very slow). analytical approach is only avaialble when BC=’SS’ and k3=0 (or k3_norm=0).

file_stem : string, optional

File stem for saving deflection data. Default file_stem=”specbeam”.

force_calc : [False, True]

If True then calcualtion will happen regardless of presence of output file (file_stem)

Notes

defl and def_norm will be 1d arrays if len(xvals)=1.

Attributes:
phi : function

Relevant Galerkin trial function. Depends on BC. See [Rb51e5289a983-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

I don’t think this is relevant anymore. 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 [Rb51e5289a983-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 [Rb51e5289a983-1]. w(x) = sum(qk * phi_k(x)). Only valid after running calulate_qk.

lam_qdotdot : 2d array of shape (nterms, nterms)

spectral coefficeient matrix for qdotdot

psi_q : 2d array of shape (nterms, nterms)

spectral coefficeient matrix for q

psi_qdot : 2d array of shape (nterms, nterms)

spectral coefficeient matrix for qdot.

k3barj : 1d array of float with len(nquad)

value of k3 interpolated from k3bar at the quadrature points. Used in the k3*w**3 term.

L_norm : float

Normalised length. Always =1.

defl_norm : array of float, size (len(xvals_norm), len(tvals_norm))

normalised defelction at xval,tval

defl : array of float, size (len(xvals_norm), len(tvals_norm))

deflection at xval,tval. ONly if self.L is not None.

beta_block : 1d array float with size (2*nterms)

Block column vector of repeated beta values.

has_moving_loads : boolean

If True then there are moving point loads.

has_stationary_loads : boolean

If True then there are stationary loads.

moving_load_w_norm : list of ndarray

Normalised deflections under each moving load at each time value. Each element of the list corresponds to each moving load. The array is of size (len(mvpl.x), len(tvals)) i.e. for each axle in the moving load group.

moving_load_w : list of ndarray

Non normalised deflections under each moving load.

Methods

animateme([xlim, ylim, norm, saveme, interval]) Animate the deflection vs distance over time plot
calulate_qk([t, t_norm]) Calculate the nterm Galerkin coefficients qk at each time value
defl_vs_time_under_loads([xlim_norm]) Deflection vs time under each moving load
load_defl() Loads defl_norm from txt file if it exists
make_E_Igamv_the() sum contributions from all loads
make_output() make all output
make_time_dependent_arrays() make all time dependent arrays
plot_w_vs_x_overtime([ax, norm]) Plot the normalised deflection vs distance for all times
saveme() Save deflection vs time to file
time_independant_matrices() Make all time independent matrices
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]
onClick  
phiCC  
phiFF  
phiSS  
runme  
animateme(xlim=None, ylim=None, norm=True, saveme=False, interval=50)[source]

Animate the deflection vs distance over time plot

Will display beam osciallations and evolving max and min deflection envelopes. Position and magnitude of load will also be shown.

Parameters:
xlim : tuple of 2 floats, optional

Limits of x-axis , default xlim=None, i.e. use Beam ends

ylim : tuple of 2 floats, optional

Limits of y-axis , default ylim=None, i.e. limits will be calculated ylim=(1.5*defl_min,1.1*defl_max)

norm : True/False, optional

If norm=True normalised values will be plotted. Default norm=True.

saveme : False/True

Whether to save the animation to disk. Default saveme=False.

interval : float

Number of miliseconds for each frame. Default interval=50.

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.

defl_vs_time_under_loads(xlim_norm=(0, 1))[source]

Deflection vs time under each moving load

i.e. What is deflection under each axle

Creates and populates self.moving_load_w_norm . Normalised deflections under each moving load at each time value. Each element of the list corresponds to each moving load. The array is of size (len(mvpl.x), len(tvals)) i.e. for each axle in the moving load group.

Parameters:
xlim_norm : tuple, optional

won’t show defelction untill within these normalised bounds Default xlim_norm=(0,1)

load_defl()[source]

Loads defl_norm from txt file if it exists

Returns:
Loaded : boolean

If True then self.defl either already exists or has been successfully loaded.

make_E_Igamv_the()[source]

sum contributions from all loads

Calculates all contributions to E*inverse(gam*v)*theta part of solution defl=phi*vE*inverse(gam*v)*theta. i.e. surcharge, vacuum, top and bottom pore pressure boundary conditions. make_E_Igamv_the will create self.E_Igamv_the. self.E_Igamv_the is an array of size (nterms, len(tvals)). So the columns are the column array E*inverse(gam*v)*theta calculated at each output time. This will allow us later to do defl = phi*v*self.E_Igamv_the

See also

_make_E_Igamv_the_mvl
moiving point load contribution
make_output()[source]

make all output

make_time_dependent_arrays()[source]

make all time dependent arrays

See also

self.make_E_Igamv_the

onClick(event)[source]
phiCC(x, beta)[source]
phiFF(x, beta)[source]
phiSS(x, beta)[source]
plot_w_vs_x_overtime(ax=None, norm=True)[source]

Plot the normalised deflection vs distance for all times

runme()[source]
saveme()[source]

Save deflection vs time to file

Deflection will be saved to : ‘self.file_stem + “_defl.csv’ Normalised deflection will be saved to: ‘self.file_stem + “_defl_norm.csv”

Might not handle singel x value and single t vlue.

time_independant_matrices()[source]

Make all time independent matrices

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, shape(len(x),len(t))

Deflections at x and self.t_norm[tslice]

geotecha.beam_on_foundation.specbeam.SpecBeam_const_mat_midpoint_defl_runme_analytical_play()[source]

Test SpecBeam for constant mat: close to Ding et al Figure 8, displacement vs time at beam midpoint (using the runme method) but with k3=0

geotecha.beam_on_foundation.specbeam.SpecBeam_stationary_load_const_mat_midpoint_defl_runme_analytical_play()[source]

Test SpecBeam for constant mat: close to Ding et al Figure 8, displacement vs time at beam midpoint (using the runme method) but with k3=0

geotecha.beam_on_foundation.specbeam.align_yaxis(ax1, ax2)[source]

Adjust y-axis limits so zeros of the two axes align, zooming them out by same ratio.

Parameters:
ax1, ax2 : matpolotlib.Axes

Axes to align.

geotecha.beam_on_foundation.specbeam.article_figure_01()[source]

Figure 1 Mid point deflection comparison with Ding et al 2012, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_02()[source]

Figure 2 Deflection amplification factor compared with analytical infinite beam solution, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_03()[source]

Figure 3 Low damping DAF for velocity ratio alp=0p5, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_04()[source]

Figure 4 Relative DAF for various stiffness ratios, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_05a_05b_05c_07()[source]

Figure 5 a),b),&c) Max deflection envelope for bet=0p1 k2_k1=10, and Figure 7 Deflection vs time at point which experiences maximum deflection alp=0p8, bet=0p1, k2_k1=10, soft to stiff=forward, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_06()[source]

Additional deflection for transition effects, bet=0p1, alp=10, (dashed=stiff to soft, solid=soft to stiff), figure from Walker and Indraratna (in press).

This takes forever to run.

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_08()[source]

Figure 8 Deflection vs time under four moving loads alp=0p8, bet=0p1, k2/k1=10, soft to stiff, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_09()[source]

Figure 9 is drawing in a word file, this routine is just a marker to indicated that figure 9 has not been forgotten, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_10_11()[source]

Figure 10 Measured vs modelled deflections at various locations in the transition zone, Figure 11 Peak displacement and stiffness distribution, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.article_figure_12()[source]

Figure 12 python code for relative dynamic amplification factor calculation, figure from Walker and Indraratna (in press).

References

[1]Walker, R.T.R. and Indraratna, B, (in press) “Moving loads on a viscoelastic foundation with special reference to railway transition zones”. International Journal of Geomechanics.
geotecha.beam_on_foundation.specbeam.case_study_PaixaoEtAl2014(saveas=None, nterms=100, force_calc=True, nx=100, nt=100)[source]

Try and replicate Paixao et al 2014 for 6 rail cars going over transition using specbeam.

Some eplots will be produced.

You’ll have to read the source code for this one.

Parameters:
saveas : string, optional

Filename stem (not including extension) to save figure to. Default saveas=None, i.e. figure won’t be saved to disk.

nx : int, optional

Number of points between xlim[0] and xlim[1] at which to calculate deflection. Default nx = 100.

nt : int, optional

Number of evaluation times from when load enters xlim[0] and leaves xlim[1]. Default nt=100

nterms : int, optional

Number of series terms in approximation. More terms leads to higher accuracy but also longer computation times. Start of with a small enough value to observe behaviour and then increase to confirm. Default nterms=100 .

force_calc : [True, False], optional

When True all calcuations will be done from scratch. When False Current working directory will be searched for a results file named according to “saveas” and if found it will be loaded. Be careful as name of results file may not be unique to your analysis parameters. Safest to use force_calc = True, but if you are simply regenerating graphs then consider force_calc=False. Default force_calc=True.

geotecha.beam_on_foundation.specbeam.dingetal_figure_8(v_nterms=(50, 75, 150, 200))[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.

geotecha.beam_on_foundation.specbeam.esveld_fig6p17()[source]

Displacement shapes before and after a moving point load on infinite beam on elastic foundation for a sub crtical , critical, super critcal velocity ratio (alpha) and damping ratios (beta) , Esveld Figure 6p17.

Esveld, C. (2001). Modern railway track. MRT-productions Zaltbommel, The Netherlands, Netherlands.

Returns:
fig : matplotlib figure

the deflection shape plot.

geotecha.beam_on_foundation.specbeam.esveld_fig6p18()[source]

Dynamic Amplification Factor vs velocity ratio (alpha) for variuos damping ratios (beta) for moving point load on infinite beam on elastic foundation, Esveld Figure 6p18.

Esveld, C. (2001). Modern railway track. MRT-productions Zaltbommel, The Netherlands, Netherlands.

Returns:
fig : matplotlib figure

the Dynamic amplification plotplot of

geotecha.beam_on_foundation.specbeam.moving_load_on_beam_on_elastic_foundation(alp, bet, slim=(-6, 6), npts=1001)[source]

Displacement vs distance from moving point load on infinite beam on elastic foundation.

This is to generate data for a a single subplot of esveld Figure 6.17

Esveld, C. (2001). Modern railway track. MRT-productions Zaltbommel, The Netherlands, Netherlands.0

s is in multiples of lambda.

Parameters:
alp : float

velocity ratio.

bet : float

damping ratio

slim : two element tuple of float, optional

Distance, in multiples of characteristic length, on either side of moving load to assess deflection at Default slim = (-6,6).

npts : int, optional

number of points witin slim to evaluate deflection. Use an odd value to capture the mid point directly under the mmoving load. Default npts=1001.

Returns:
s : 1D numpy array of float

Distances, in multiples of characteristic length, from moving point load. w.r.t. normalised by characteristic lengths

disp : 1D numpy array of float

Displacements at s.

geotecha.beam_on_foundation.specbeam.multi_static(xc, frel, wsingle=1, xlim=(-10, 10), n=500, L=1)[source]

Deflection shape for multiple static point loads

Parameters:
xc : float

x-coord of each point load

frel : float

relative magnitude of each point load

wsingle : float, optional

maximum deflection under unit load

xlim : tuple with two values

xmin, xmax of evaluation

L : float, optional

Characterisic length

Returns:
x,w : 1d ndarray

x coords, and static deflection

See also

point_load_static_defl_shape_function
single load function
geotecha.beam_on_foundation.specbeam.point_load_static_defl_shape_function(x, L, xc=0.0)[source]

Deflection shape function for static point load on beam on foundation.

Parameters:
x : float

x value

L : float

Characteristc length (1/lam)

x0 : float, optional

position of load.

Returns:
nu : float

shape function at x

geotecha.beam_on_foundation.specbeam.reverse_curve_transition(x, x1, y1, x2, y2)[source]

Transition curve between (x1, y1) and (x2, y2)

Connects two parallel horizontal lines by two circular arcs.

Parameters:
x : array of float

value(s) at which to caluclate transition curve

x1, y1 : float

start point of transition

x2, y2 : float

end point of transition

Returns:
y : float

value of transition curve at x.

Notes

Two circular curves transition in a rectangle h high by L wide. Each curve has radius R and traces out a curve of delta radians. h=2*R*(1-cos(delta)) & L=2*R*(sin(delta))

Solve these in wolfram alpha for delta :Solve[(1-Cos[x])/Sin[x]-h/L,x] to get : delta = 2*arctan(h/L) Then use equation of circle in particualr quadrant to get points on curve. adjust for x1,y1 and x2,y2 starting points.

geotecha.beam_on_foundation.specbeam.specbeam_vs_analytical_single_static_load()[source]

Compare specbeam vs analytical solution for single static load at mid point of ‘long’ beam.

Will produce a plot and an animation.

Note that you can’t really have a static load in specbeam. You can have a slowly applied load that then remains constant. Thus depending on the beam properties e.g. damping, you can get a maximum displacement envelope that is higher than the static case in that the momentum of the beam carries the deflection past the ‘analytical” maximum and then oscillations are created untill they are damped out.

geotecha.beam_on_foundation.specbeam.static_point_load_shape(x, L)[source]

Deflection shape parameter for point load at x=0 on infinite beam on elastic foundation with characteristic length L.

Deflection is relative to the maximum value

Parameters:
x : float

x coordinate(s)

L : float

Beam characteristic Length.

Returns:
out : float

Relative deflection at x.

Notes

FYI the defelctions are relative to w_max which is :
L = Characteristic length = k1 / ((4 * E * I)**0.25) Q = Point load w_max = static = Q * lam / (2 * k1)
geotecha.beam_on_foundation.specbeam.tanhtran(x, a, b)[source]

Transition between two functions using tanh

t(x) Smoothly transitions between 0 and 1;

Used with:
[transition f(x) to g(x)](x) = g(x)*t(x) + (1-t(x))*f(x)
Parameters:
x : float

x values to evaluate at

a : float

Position of tranistion. Often a=0.

b : float,

b close to zero gives a sharp change.

Returns:
tfunc: float

0.5 + 0.5 * tanh((x-a)/b)

See also

tanhtrani
inverse of tanhtran
tanhtranib
get smoothing factor from speccific x,y
geotecha.beam_on_foundation.specbeam.tanhtrani(y, a, b)[source]

Inverse of Transition between two functions using tanh

x value that matches y

Parameters:
y : float

y value to evaluate at

a : float

Position of tranistion. Often a=0.

b : float,

b close to zero gives a sharp change.

Returns:
tfunc: float

b * arctanh((y - 0.5) / 0.5) + a

See also

tanhtran
inverse of tanhtrani
tanhtranib
get smoothing factor from specific point.
geotecha.beam_on_foundation.specbeam.tanhtranib(x, y, a)[source]

Smoothing factor of tansition curve to get known (x,y)

b value such that t(x) goes through (x,y)

Parameters:
x,y : float

x values to evaluate at

a : float

Position of tranistion. Often a=0.

Returns:
tfunc: float

(x-a) / arctanh((y - 0.5) / 0.5)

See also

tanhtran
transition function
tanhtrani
inverse of tanhtran
geotecha.beam_on_foundation.specbeam.transition_DAF1(saveas=None, ax=None, kratios=None, alpmax=3, npts=500, betas=None, article_formatting=False, DAFmax=4)[source]

Plot of Dynamic Amplification Factor (DAF) vs velocity ratio for different stiffness and damping values.

Parameters:
saveas : string, optional

Filename (not including extension) to save figure to. Default saveas=None, i.e. figure won’t be saved to disk.

ax : matplotlib.Axes object, optional

Axes object to plot on. Default ax=None, i.e. figure and axis will be created. Nb this might return an error if used as the function will try and return then figure.

kratios : list of float, optional

striffness ratio values. Default kratios=None which will mean kratios=[1,2].

alpmax : float, optional

maximum velocit ratio considered (i.e x-axis max). Default alpmax=3

npts : int, optional

velocity ratios (i.e x-axis values) will be npts between 0 and alpmax. Default npts=500.

betas : list of float, optional

list of damping ratios to consider. Default betas=None which will mean betas=[0,0.05, 0.1, 0.3, 1.1, 2.0]

article_formatting : [False, True], optional

You should never need this. It was a quick hack to do some custom formatting for a journal article. Dont even try messing with it. will ony work with 6 beta values or less. Default=False, i.e. no special formatting.

DAFmax : float, optional

Y-axis max. Default DAFmax=4

Returns:
fig : matplotlib figure

the figure

geotecha.beam_on_foundation.specbeam.transition_DAF1_at_specific_alp(saveas=None, ax=None, alphas=None, kratios=None, betas=None, DAFmax=4)[source]

Plot of Dynamic Amplification Factor (DAF) vs stiffness ratio for different velocity ratios and damping values.

Parameters:
saveas : string, optional

Filename (not including extension) to save figure to. Default saveas=None, i.e. figure won’t be saved to disk.

ax : matplotlib.Axes object, optional

Axes object to plot on. Default ax=None, i.e. figure and axis will be created. Nb this might return an error if used as the function will try and return then figure.

alphas : list of float,optional

Specific velocity ratios to calc at. Default alphas=None which will use [0.1, 0.5, 0.8].

kratios : list of float, optional

striffness ratio values. Default kratios=None which will mean kratios=np.logspace(np.log10(0.1), np.log10(100), 100), i.e 100 points logspaced between 0.1 and 100.

betas : list of float, optional

list of damping ratios to consider. Default betas=None which will

DAFmax : float, optional

Y-axis max. Default DAFmax=4

Returns:
fig : matplotlib figure

the figure

geotecha.beam_on_foundation.specbeam.transition_zones1(Lt_Lcs=[0.1, 4], alphas=[0.8], reverses=[False, True], xi_Lc=None, kratio=50, beta=0.2, alpha_design=0.8, ntrans=50, DAF_distribution='linear', tanhend=0.01, xlim=(-8, 8), nx=100, nt=100, t_extend=1.2, nterms=150, ninterp=15, nloads=1, load_spacing=0.03, animateme=False, saveas=None, force_calc=True, article_formatting=False)[source]

Exploration of deflection envelopes caused by point loads travelling through a stiffness transition zone on beam on visco-elastic foundation using Specbeam.

  • Transition between reference beam/foundation with stiffness k to beam on

foundation with stiffness k*`kratio`. - Transition zone length Lt is a multiple of the characteristic length Lc of the reference beam/foundation. Lt = Lt_Lc`*Lc - Transition zone divided into ntrans points. - Change in stiffness over transition zone is a calculated value based on a goal displacement profile. The shape of the goal displacement profile is controlled be `DAF_distribution which is usually linear. - At each evaluation point in transition zone the stiffness needed to achieve the desired displacement value is back calcualted based on infinite beam theory with a point load travelling at velocity ratio of alpha_design with damping ratio beta (alpha_design and beta are w.r.t. the reference beam/foundation). - for input into specbeam the design stiffness profile piecewise linear interpolated at ninterp points to produce piecewise linear analysis profile. ninterp should be less thatn or equal to ntrans - Moving load will then traverse the design stiffness profile and depending on dynamic effects the moving load displacement may or may not match the design displacements. You can only design for one speed so deflection response might be differnt when actual speed is different to design speed. - Note only the stiffness changes. We assume that damping determined by beta of the reference beam/foundation does not change in the transition zone. Perhaps not realistic but simply an assumption to enable exploration.

3 plots will be produced: - Line will be plotted based on permututions of Lt_Lcs, alphas, reverses. So number of lines can quickly get out of hand. Usually interested in varyin one of the parameters at a time. - Plot 1a) shows the design transtion stiffness profile Plot 1b) shows max deflection experienced at each point - Plot 2) is deflection vs time at the point which experiences the

maximum deflection. also plots a marker at the time when load passes the point.
  • Plot 3) is the deflection vs time curve under the moving point load.
Parameters:
Lt_Lcs : list of float, optional

Transtion zone lengths as a mulitple of characteristic length. Default Lt_Lcs = [0.1,4].

alphas : list of float, optional

Velocity ratios to analyze. Velocity rations are with respect to the reference beam/foundation. Default alphas = [0.8]

reverses : list of Bool, optional

Directions of travel. False means left to right from reference beam/foundation to other beam/foundation. True means the reverse direction from right to left from the other beam/foundation to the reference beam foundation. Default reverses=[False,True].

xi_Lc : [False, True], optional

If True then deflection vs time will be plotted at the point where the maximum deflection was experienced. Default xi_Lc=None i.e. False and no plot will be created. Also plots a marker at the time when load passes the point.

kratio : float optional

stiffness ratio k2/k1 where k1 is the stiffness of the reference beam/foundation. Default kratio = 50.

beta : float, optional

Reference value of damping ratio. Default beta = 0.2

alpha_design : float, optional

Design velocity ratio used to back calculate the design stiffness ratio. Default alpha_design = 0.8.

ntrans : int, optional

Number of points to approximate the desired displacement transition profile into, to then back calculate the stiffness profile used for analysis. Default ntrans = 50 .

DAF_distribution : [“linear”, “tanh”, “reverse”], optional

Shape of goal transition deflection distrribution. “linear” isself explanatory and simplest. “tanh” is a hyperbolic tangent, with steepness of transtion controlled by tanend. Not well tested. “reverse” is a compund reverse curve of two curclar arcs tangent to the end points and tangent at the arc join. Not well tested Default DAF_distribution = “linear”.

tanend : float, optional

value controling the steepness of a “tanh” DAF_distribution. Say transitioning between points (-0.5, -1) and (0.5, 1). Tanh curve will never actually hit valus of -1 and +1, it will approach those values as x goes to += infinity. Using tanend=0.01 will control transtion such that the tanh curve evalutes to y=1+0.01 and y=1-0.01 at the start/end points. This evaluated value is vever actually used because we have given start and end points but it does control the position of points just inside the start/end points and hence the steepness of transition.

xlim : two element tuple of float, optional

Distances on either side of the transtion zone midpoint to evaluate deflection at. The mid point of the transition would be at 0. For example default xlim=(-8,8) which would be 8 characteristic lenghts on either side of the transition zone midpoint. Note: Characteristic lenth Lc is with respect to the reference beam/foundation. Default xlim=(-8, 8).

nx : int, optional

Number of points between xlim[0] and xlim[1] at which to calculate deflection. Default nx = 100.

nt : int, optional

Number of evaluation times from when load enters xlim[0] and leaves xlim[1]. Default nt=100

t_extend : float, optional

factor to extend maximum evalutated time. As per nt’ evaluation times are normally between first point entering `xlim[0]’ and leaving `xlim[1]. However, with mulitple loads, loads can still be on the beam at that point. max eval time can be extended fby a factor. t_extend = 1.2 .

ninterp : int, optional

For use in specbeam the back calculated stiffness and damping transitions profiles are interpoalted at ninterp points. A coarser discretization will improve computation times but might no capture the changing properties approraitely. Note: ninterp should be less thatn or equal to ntrans. Default ninterp = 15 .

nterms : int, optional

Number of series terms in approximation. More terms leads to higher accuracy but also longer computation times. Start of with a small enough value to observe behaviour and then increase to confirm. Default nterms=150 .

nloads : int, optional

Number of eqully spaced moving point loads. Default nloads = 1 .

load_spacing : float, optional

Spacing between the multiple points loads. Units are multiples of Lc characteristic length. Default load_spacing = 0.03 .

animateme : [False, True], optional

If True the deflectons and moving loads will be animated. Default anmateme=False.

saveas : string, optional

Filename stem (not including extension) to save figure to. Default saveas=None, i.e. figure won’t be saved to disk.

force_calc : [True, False], optional

When True all calcuations will be done from scratch. When False Current working directory will be searched for a results file named according to “saveas” and if found it will be loaded. Be careful as name of results file may not be unique to your analysis parameters. Safest to use force_calc = True, but if you are simply regenerating graphs then consider force_calc=False. Default force_calc=True.

article_formatting : [False, True], optional

You should never need this. It was a quick hack to do some custom formatting for a journal article. Dont even try messing with it. will ony work with 3 Lt_Lcs values values or less. Default=False, i.e. no special formatting.

Notes

While all parameters are normalised the raw properties that all others are based on are: - E = 6.998*1e9, #Pa - rho = 2373, #kg/m3 - L = 160, #m - A = 0.3*0.1, # m^2 - I = 0.00022477900799999998, - k1 = 800002.6125,

geotecha.beam_on_foundation.specbeam.transition_zones2(animateme=False, saveas=None, xlim=(-8, 8), ntrans=50, kratio=50, beta=0.2, alphas=[0.8], nx=100, nt=100, nterms=150, force_calc=True, ninterp=15, nloads=1, load_spacing=0.03, alpha_design=0.8, Lt_Lcs=[0.1, 4], t_extend=1.1, DAF_distribution='linear', tanhend=0.01, reverses=[False, True], xi_Lc=None, prefix='', article_formatting=True)[source]

Essentially a duplicate of transition_zones1 to get a DAF* vs Lt/Lc plot for forward and reverse plots

Be warned that this can produce hundreds of plots and take hours to run when all you are interested is the single resulting plot.

Basically a transtion_zones1 analysis will be run. The maximum defelction (w.r.t. to the infinite beam deflection at same speed and damping) for that run will be recorded. That will be one point for one alpha value for one Lt/Lc value. Points will be joined to form a line for for each alphas and reverses combination.

I used this once for a journal article figure.