speccon1d_vrc

Class summary

Speccon1dVRC([reader, pkg_for_version]) Multilayer consolidation with stone columns including vertical and radial drainage in the column, using the spectral Galerkin method.

Function summary

main() Run speccon1d_vrc as a script

Module listing

Multilayer consolidation with stone columns including vertical and radial drainage using the spectral Galerkin method.

class geotecha.speccon.speccon1d_vrc.Speccon1dVRC(reader=None, pkg_for_version='geotecha')[source]

Bases: geotecha.speccon.speccon1d.Speccon1d

Multilayer consolidation with stone columns including vertical and radial drainage in the column, using the spectral Galerkin method.

Features:

  • Multiple layers.
  • Vertical and radial drainage in a unit cell. (radial drainage uses the eta method).
  • Finite vertical and radial permeability in the column.
  • Material properties that are constant in time but piecewsie linear with depth (including column permeability and stiffness). Soil-column stiffness is modelled using a lumped soil/column volume compressibilty (mv).
  • Surcharge loading (vacuum loading is via specifying non-zero negative pore pressure at top and bottom of drain.
  • Non-zero top and bottom pore pressure boundary conditions (top and bottom pore pressures are the same in the column and soil).
  • Surcharge/Boundary Conditions/vary with time in a piecewise-linear function multiplied by a cosine function of time.
    • Surcharge can also vary piecewise linear with depth. The depth dependence does not vary with time.
    • Mulitple loads will be combined using superposition.
  • Subset of Python syntax available in input files/strings allowing basic calculations within input files.
  • Output:
    • Excess pore pressure at depth in soil and column.
    • Average excess pore pressure between two depths (in soil, column and overall).
    • Settlement between two depths.
    • Charts and csv output available.
  • Program can be run as script or in a python interpreter.
  • Note there is no pumping or fixed pore pressure functionality.

Warning

The ‘Parameters’ and ‘Attributes’ sections below require further explanation. The parameters listed below are not used to explicitly initialize the object. Rather they are defined in either a multi-line string or a file-like object using python syntax; the file/string is then used to initialize the object using the reader parameter. As well as simple assignment statements (H = 1, drn = 0 etc.), the input file/string can contain basic calculations (z = np.linspace(0, H, 20) etc.). Not all of the listed parameters are needed. The user should pick an appropriate combination of attributes for their analysis (minimal explicit checks on input data will be performed). Each ‘parameter’ will be turned into an attribute that can be accessed using conventional python dot notation, after the object has been initialised. The attributes listed below are calculated values (i.e. they could be interpreted as results) which are accessible using dot notation after all calculations are complete.

Parameters:
H : float, optional

Total height of soil profile. Default H=1.0. Note that even though this program deals with normalised depth values it is important to enter the correct H valu, as it is used when plotting, outputing data and in normalising gradient boundary conditions (see bot_vs_time below).

mvref : float, optional

Reference value of volume compressibility mv (used with H in settlement calculations). Default mvref=1.0.

kvref : float, optional

Reference value of vertical permeability kv in soil (only used for pretty output). Default kvref=1.0.

khref : float, optional

Reference value of horizontal permeability kh in soil (only used for pretty output). Default khref=1.0.

kvcref : float, optional

Reference value of vertical permeability kvc in column (only used for pretty output). Default kvcref=1.0.

khcref : float, optional

Reference value of horizontal permeability khc in column (only used for pretty output). Default khcref=1.0.

etref : float, optional

Reference value of lumped drain parameter et (only used for pretty output). Default etref=1.0. et = 2 / (mu * re^2) where mu is smear-zone/geometry parameter and re is radius of influence of vertical drain.

drn : {0, 1}, optional

drainage boundary condition. Default drn=0. 0 = Pervious top pervious bottom (PTPB). 1 = Pervious top impoervious bottom (PTIB).

dT : float, optional

Convienient normaliser for time factor multiplier. Default dT=1.0.

neig : int, optional

Number of series terms to use in solution. Default neig=2. Don’t use neig=1.

dTv : float, optional

Vertical reference time factor multiplier. dTv is calculated with the chosen reference values of kv and mv: dTv = kv /(mv*gamw) / H ^ 2

dTvc : float, optional

Well reference time factor multiplier. dTvc is calculated with the chosen reference values of kvc and mv: dTvc = kvc /(mv*gamw) / H ^ 2

dTh : float, optional

horizontal reference time factor multiplier. dTh is calculated with the reference values of kh, et, and mv: dTh = kh / (mv * gamw) * et

dThc : float, optional

Horizontal reference time factor multiplier in column. dThc is calculated with the reference values of khc, and mv: dThc = 8*khc / (mv * gamw) /rc**2.

mv : PolyLine

Normalised volume compressibility PolyLine(depth, mv). The mv here is the value of ms / (1 + alp * (Y - 1)) normalised by mvref. Y is column to soil stiffness ratio ms/mc or Ec/Es. alp = 1/n**2 where n is the ratio between influence radius and column radius n=re/rw.

kh : PolyLine, optional

Normalised horizontal permeability in soil PolyLine(depth, kh).

kv : PolyLine , optional

Normalised vertical permeability in soil PolyLine(depth, kv).

khc : PolyLine, optional

Normalised horizontal permeability in column PolyLine(depth, kh).

kvc : PolyLine , optional

Normalised vertical permeability in column PolyLine(depth, kv).

et : PolyLine, optional

Normalised vertical drain parameter PolyLine(depth, et). et = 2 / (mu * re^2) where mu is smear-zone/geometry parameter and re is radius of influence of vertical drain.

surcharge_vs_depth : list of Polyline, optional

Surcharge variation with depth. PolyLine(depth, multiplier).

surcharge_vs_time : list of Polyline, optional

Surcharge magnitude variation with time. PolyLine(time, magnitude).

surcharge_omega_phase : list of 2 element tuples, optional

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

top_vs_time : list of Polyline, optional

Top p.press variation with time. Polyline(time, magnitude).

top_omega_phase : list of 2 element tuples, optional

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

bot_vs_time : list of Polyline, optional

Bottom p.press variation with time. Polyline(time, magnitude). When drn=1, i.e. PTIB, bot_vs_time is equivilent to saying D[u(H,t), z] = bot_vs_time. Within the program the actual gradient will be normalised with depth by multiplying H.

bot_omega_phase : list of 2 element tuples, optional

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

ppress_z : list_like of float, optional

Normalised z to calculate pore pressure at.

avg_ppress_z_pairs : list of two element list of float, optional

Nomalised zs to calculate average pore pressure between e.g. average of all profile is [[0,1]].

settlement_z_pairs : list of two element list of float, optional

Normalised depths to calculate normalised settlement between. e.g. surface settlement would be [[0, 1]].

tvals : list of float

Times to calculate output at.

ppress_z_tval_indexes : list/array of int, slice, optional

Indexes of tvals at which to calculate ppress_z. i.e. only calculate ppress_z at a subset of the tvals values. Default ppress_z_tval_indexes=slice(None, None) i.e. use all the tvals.

avg_ppress_z_pairs_tval_indexes : list/array of int, slice, optional

Indexes of tvals at which to calculate avg_ppress_z_pairs. i.e. only calc avg_ppress_z_pairs at a subset of the tvals values. Default avg_ppress_z_pairs_tval_indexes=slice(None, None) i.e. use all the tvals.

settlement_z_pairs_tval_indexes : list/array of int, slice, optional

Indexes of tvals at which to calculate settlement_z_pairs. i.e. only calc settlement_z_pairs at a subset of the tvals values. Default settlement_z_pairs_tval_indexes=slice(None, None) i.e. use all the tvals.

implementation : [‘scalar’, ‘vectorized’,’fortran’], optional

Where possible use the stated implementation type. ‘scalar’= python loops (slowest), ‘vectorized’ = numpy (fast), ‘fortran’ = fortran extension (fastest). Note only some functions have multiple implementations.

RLzero : float, optional

Reduced level of the top of the soil layer. If RLzero is not None then all depths (in plots and results) will be transformed to an RL by RL = RLzero - z*H. If RLzero is None (i.e. the default) then all depths will be reported z*H (i.e. positive numbers).

plot_properties : dict of dict, optional

Dictionary that overrides some of the plot properties. Each member of plot_properties will correspond to one of the plots.

plot_properties description
por dict of prop to pass to overall pore pressure plot.
porc dict of prop to pass to soil pore pressure plot.
pors dict of prop to pass to average pore pressure plot.
avp dict of prop to pass to overall average pore pressure plot.
avps dict of prop to pass to soil average pore pressure plot.
avpc dict of prop to pass to c average pore pressure plot.
set dict of prop to pass to settlement plot.
load dict of prop to pass to loading plot.
material dict of prop to pass to materials plot.

see geotecha.plotting.one_d.plot_vs_depth and geotecha.plotting.one_d.plot_vs_time for options to specify in each plot dict.

save_data_to_file : True/False, optional

If True data will be saved to file. Default save_data_to_file=False

save_figures_to_file : True/False

If True then figures will be saved to file. Default save_figures_to_file=False

show_figures : True/False, optional

If True the after calculation figures will be shown on screen. Default show_figures=False.

directory : string, optional

Path to directory where files should be stored. Default directory=None which will use the current working directory. Note if you keep getting directory does not exist errors then try putting an r before the string definition. i.e. directory = r’C:Users…’

overwrite : True/False, optional

If True then existing files will be overwritten. Default overwrite=False.

prefix : string, optional

Filename prefix for all output files. Default prefix= ‘out’

create_directory : True/Fase, optional

If True a new sub-folder with name based on prefix and an incremented number will contain the output files. Default create_directory=True.

data_ext : string, optional

File extension for data files. Default data_ext=’.csv’

input_ext : string, optional

File extension for original and parsed input files. default = “.py”

figure_ext : string, optional

File extension for figures. Can be any valid matplotlib option for savefig. Default figure_ext=”.eps”. Others include ‘pdf’, ‘png’.

title : str, optional

A title for the input file. This will appear at the top of data files. Default title=None, i.e. no title.

author : str, optional

Author of analysis. Default author=’unknown’.

Notes

Gotchas

All the loading terms e.g. surcharge_vs_time, surcharge_vs_depth, surcharge_omega_phase can be either a single value or a list of values. The corresponding lists that define a load must have the same length e.g. if specifying multiple surcharge loads then surcharge_vs_time and surcharge_vs_depth must be lists of the same length such that surcharge_vs_time[0] can be paired with surcharge_vs_depth[0], surcharge_vs_time[1] can be paired with surcharge_vs_depth[1], etc.

Material and geometric properties

  • \(k_v\) is vertical soil permeability.
  • \(k_h\) is horizontal soilpermeability.
  • \(k_{vc}\) is vertical permeability in column.
  • \(k_{hc}\) is horizontal permeability in column.
  • \(m_s\) is the volume compressibility in the soil.
  • \(m_c\) is the volume compressibility in the column.
  • \(m_v\) is lumped volume compressibility. \(m_v=m_s/(1+\alpha(Y-1))\). Stiffness ratio \(Y=m_s/m_c=E_c/E_s\). Also \(\alpha=1/n^2\).
  • \(\eta\) is the radial drainage parameter \(\eta = \frac{2}{r_e^2 \mu}\).
  • \(r_e\) is influence radius of drain.
  • \(r_w\) is drain radius.
  • \(n=r_e/r_w\) is ratio of influence radius to drain radius.
  • \(\mu\) is any of the smear zone geometry parameters dependent on the distribution of permeabilit in the smear zone (see geotecha.consolidation.smear_zones).
  • \(\gamma_w\) is the unit weight of water.
  • \(Z\) is the nomalised depth (\(Z=z/H\)).
  • \(H\) is the total height of the soil profile.

Governing equation Radially averaged soil and column pore pressures, \(u_s\), \(u_c\) are related to the overall average pore pressure, \(u\) by:

\[u=\left({1-\alpha}\right)u_s + \alpha u_c\]

where \(\alpha=1-n^2=1-\left({r_e/r_w}\right)\). Equal strain in the soil and column requires:

\[\frac{\varepsilon}{m_{vref}}= \overline{m}_v \left[{\sigma -\left({1-\alpha}\right)u_s -\alpha u_c}\right]\]

Continuity of flow in the soil, column, and across the soil-column boundary respectively requires:

\[\frac{\dot{\varepsilon}}{m_{vref}}= dT_h\overline{k}_h\overline{\eta}\left({u_s-\beta}\right) - dT_v\left({\overline{k}_v u_s,_Z}\right),_Z\]
\[\frac{\dot{\varepsilon}}{m_{vref}}= dT_{hc}\overline{k}_{hc}\left({u_c-\beta}\right) - dT_{vc}\left({\overline{k}_{vc} u_c,_Z}\right),_Z\]
\[\frac{\dot{\varepsilon}}{m_{vref}}= \left({1 - \alpha}\right) dT_v\left({\overline{k}_v u,_Z}\right),_Z - \alpha dT_{vc}\left({\overline{k}_{vc} u_c,_Z}\right),_Z\]

Each pore pressure value is a function of normalised depth \(Z\) and time \(t\). \(\beta\) is the pore pressure at the soil/column interface at a particular depth.

where

\[dT_v = \frac{k_{v\textrm{ref}}} {H^2 m_{v\textrm{ref}} \gamma_w}\]
\[dT_{vc} = \frac{k_{vc\textrm{ref}}} {H^2 m_{v\textrm{ref}} \gamma_w}\]
\[dT_{hc} = \frac{8 k_{hc\textrm{ref}}} {m_{v\textrm{ref}} \gamma_w r_c^2}\]
\[dT_h = \frac{k_{h\textrm{ref}} \eta_{\textrm{ref}}} {m_{v\textrm{ref}} \gamma_w}\]
\[\eta = \frac{2}{r_e^2 \mu}\]

\(\mu\) is any of the smear zone geometry parameters dependent on the distribution of permeabilit in the smear zone (see geotecha.consolidation.smear_zones).

The overline notation represents a depth dependent property normalised by the relevant reference property. e.g. \(\overline{k}_v = k_v\left({z}\right) / k_{v\textrm{ref}}\).

A comma followed by a subscript represents differentiation with respect to the subscripted variable e.g. \(u,_Z = u\left({Z,t}\right) / \partial Z\).

Non-zero Boundary conditions

The following two sorts of boundary conditions (identical in in both the soil and column) can be modelled:

\[\left.u\left({Z,t}\right)\right|_{Z=0} = u^{\textrm{top}}\left({t}\right) \textrm{ and } \left.u\left({Z,t}\right)\right|_{Z=1} = u^{\textrm{bot}}\left({t}\right)\]
\[\left.u\left({Z,t}\right)\right|_{Z=0} = u^{\textrm{top}}\left({t}\right) \textrm{ and } \left.u\left({Z,t}\right),_Z\right|_{Z=1} = u^{\textrm{bot}}\left({t}\right)\]

The boundary conditions are incorporated by homogenising the governing equation with the following substitution (same process for the column and beta):

\[u\left({Z,t}\right) = \hat{u}\left({Z,t}\right) + u_b\left({Z,t}\right)\]

where for the two types of non zero boundary boundary conditions:

\[u_b\left({Z,t}\right) = u^{\textrm{top}}\left({t}\right) \left({1-Z}\right) + u^{\textrm{bot}}\left({t}\right) Z\]
\[u_b\left({Z,t}\right) = u^{\textrm{top}}\left({t}\right) + u^{\textrm{bot}}\left({t}\right) Z\]

Time and depth dependence of loads/material properties

Soil properties do not vary with time.

Loads are formulated as the product of separate time and depth dependant functions as well as a cyclic component:

\[\sigma\left({Z,t}\right)= \sigma\left({Z}\right) \sigma\left({t}\right) \cos\left(\omega t + \phi\right)\]

\(\sigma\left(t\right)\) is a piecewise linear function of time that within the kth loading stage is defined by the load magnitude at the start and end of the stage:

\[\sigma\left(t\right) = \sigma_k^{\textrm{start}} + \frac{\sigma_k^{\textrm{end}} - \sigma_k^{\textrm{start}}} {t_k^{\textrm{end}} - t_k^{\textrm{start}}} \left(t - t_k^{\textrm{start}}\right)\]

The depth dependence of loads and material property \(a\left(Z\right)\) is a piecewise linear function with respect to \(Z\), that within a layer are defined by:

\[a\left(z\right) = a_t + \frac{a_b - a_t}{z_b - z_t}\left(z - z_t\right)\]

with \(t\) and \(b\) subscripts representing ‘top’ and ‘bottom’ of each layer respectively.

References

The genesis of this work is from research carried out by Dr. Rohan Walker, Prof. Buddhima Indraratna and others at the University of Wollongong, NSW, Austrlia, [1], [2], [3], [4].

[1](1, 2) Walker, Rohan. 2006. ‘Analytical Solutions for Modeling Soft Soil Consolidation by Vertical Drains’. PhD Thesis, Wollongong, NSW, Australia: University of Wollongong.
[2](1, 2) Walker, R., and B. Indraratna. 2009. ‘Consolidation Analysis of a Stratified Soil with Vertical and Horizontal Drainage Using the Spectral Method’. Geotechnique 59 (5) (January): 439-449. doi:10.1680/geot.2007.00019.
[3](1, 2) Walker, Rohan, Buddhima Indraratna, and Nagaratnam Sivakugan. 2009. ‘Vertical and Radial Consolidation Analysis of Multilayered Soil Using the Spectral Method’. Journal of Geotechnical and Geoenvironmental Engineering 135 (5) (May): 657-663. doi:10.1061/(ASCE)GT.1943-5606.0000075.
[4](1, 2) Walker, Rohan T. 2011. Vertical Drain Consolidation Analysis in One, Two and Three Dimensions’. Computers and Geotechnics 38 (8) (December): 1069-1077. doi:10.1016/j.compgeo.2011.07.006.
Attributes:
por, porc, pors : ndarray, only present if ppress_z is input

Calculated pore pressure at depths corresponding to ppress_z and times corresponding to tvals. This is an output array of size (len(ppress_z), len(tvals[ppress_z_tval_indexes])).

avp, avpc, avps : ndarray, only present if avg_ppress_z_pairs is input

Calculated average pore pressure between depths corresponding to avg_ppress_z_pairs and times corresponding to tvals. This is an output array of size (len(avg_ppress_z_pairs), len(tvals[avg_ppress_z_pairs_tval_indexes])).

set : ndarray, only present if settlement_z_pairs is input

Settlement between depths corresponding to settlement_z_pairs and times corresponding to tvals. This is an output array of size (len(avg_ppress_z_pairs), len(tvals[settlement_z_pairs_tval_indexes]))

Methods

check_input_attributes() Perform checks on attributes
make_E_Igamv_the() sum contributions from all loads
make_all() Run checks, make all arrays, make output
make_output() make all output
make_time_dependent_arrays() make all time dependent arrays
make_time_independent_arrays() make all time independent arrays
produce_plots() produce plots of analysis
make_E_Igamv_the()[source]

sum contributions from all loads

Calculates all contributions to E*inverse(gam*v)*theta part of solution u=phi*vE*inverse(gam*v)*theta. i.e. surcharge, vacuum, top and bottom pore pressure boundary conditions. make_load_matrices will create `self.E_Igamv_the. self.E_Igamv_the is an array of size (neig, 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 u = phi*v*self.E_Igamv_the

See also

_make_E_Igamv_the_surcharge
surchage contribution
_make_E_Igamv_the_BC
top boundary pore pressure contribution
_make_E_Igamv_the_bot
bottom boundary pore pressure contribution
make_output()[source]

make all output

make_time_dependent_arrays()[source]

make all time dependent arrays

See also

self.make_E_Igamv_the

make_time_independent_arrays()[source]

make all time independent arrays

See also

self._make_m
make the basis function eigenvalues
self._make_gam
make the mv dependent gamma matrix
self._make_psi
make the kv, kh, et dependent psi matrix
self._make_eigs_and_v
make eigenvalues, eigenvectors and I_gamv
produce_plots()[source]

produce plots of analysis

geotecha.speccon.speccon1d_vrc.main()[source]

Run speccon1d_vrc as a script