void_ratio_stress¶
Class summary¶
AvSoilModel (av, siga, ea) |
Linear void ratio-effective stress relationship |
CcCrSoilModel (Cc, Cr, siga, ea) |
Semi-log void ratio-effective stress realationship |
FunctionSoilModel (fn_e_from_stress, ...) |
Functional definition of void-ratio stress relationship |
OneDimensionalVoidRatioEffectiveStress |
Base class for defining 1D void ratio-effective stress relationships |
PwiseLinearSoilModel (siga, ea, Cr[, xlog, ylog]) |
Pwise linear void ratio-effective stress realationship |
YinAndGrahamSoilModel (lam, kap, psi, siga, ...) |
Yin and Graham creep void ratio-effective stress realationship |
Module listing¶
Relationships between void-ratio and effective stress.
-
class
geotecha.constitutive_models.void_ratio_stress.
AvSoilModel
(av, siga, ea)[source]¶ Bases:
geotecha.constitutive_models.void_ratio_stress.OneDimensionalVoidRatioEffectiveStress
Linear void ratio-effective stress relationship
Parameters: av : float
Slope of compression line.
siga, ea : float
Effective stress and void ratio specifying point on av line.
Methods
av_from_stress
(*args, **kwargs)Slope of void ratio from effective stress e_and_stress_for_plotting
(**kwargs)Void ratio and stress values that plot the model e_from_stress
(estress, **kwargs)Void ratio from effective stress plot_model
(**kwargs)plot the void ratio-stress points stress_from_e
(e, **kwargs)Effective stress from void ratio -
e_and_stress_for_plotting
(**kwargs)[source]¶ Void ratio and stress values that plot the model
Parameters: npts : int, optional
Number of points to return. Default npts=100.
xmin, xmax : float, optional
range of x (i.e. effective stress) values from which to return points. Default xmin=1, xmax=100
Returns: x, y : 1d ndarray
npts stress, and void ratio values between xmin and xmax.
-
e_from_stress
(estress, **kwargs)[source]¶ Void ratio from effective stress
Parameters: estress : float
Current effective stress.
Returns: e : float
Void ratio corresponding to current stress state.
Examples
>>> a = AvSoilModel(av=1.5, siga=19, ea=4) >>> a.e_from_stress(20) 2.5
Array inputs:
>>> a = AvSoilModel(av=1.5, siga=19, ea=4) >>> a.e_from_stress(np.array([20, 21])) array([ 2.5, 1. ])
-
stress_from_e
(e, **kwargs)[source]¶ Effective stress from void ratio
Parameters: e : float
Current void ratio.
Returns: estress : float
Effective stress corresponding to current void ratio
Examples
>>> a = AvSoilModel(av=1.5, siga=19, ea=4) >>> a.stress_from_e(1) 21.0
Array inputs:
>>> a = AvSoilModel(av=1.5, siga=19, ea=4) >>> a.stress_from_e(np.array([1, 2.5])) array([ 21., 20.])
-
-
class
geotecha.constitutive_models.void_ratio_stress.
CcCrSoilModel
(Cc, Cr, siga, ea)[source]¶ Bases:
geotecha.constitutive_models.void_ratio_stress.OneDimensionalVoidRatioEffectiveStress
Semi-log void ratio-effective stress realationship
Parameters: Cc : float
Compressibility index, slope of e-log(sig) line.
Cr : float
Recompression index, slope of e-log(sig) line.
siga, ea : float
Point on compression line fixing it in effective stress-void ratio space.
Methods
av_from_stress
(estress, **kwargs)Slope of void ratio from effective stress e_and_stress_for_plotting
(**kwargs)Void ratio and stress values that plot the model e_from_stress
(estress, **kwargs)Void ratio from effective stress plot_model
(**kwargs)plot the void ratio-stress points stress_from_e
(e, **kwargs)Effective stress from void ratio -
av_from_stress
(estress, **kwargs)[source]¶ Slope of void ratio from effective stress
Parameters: estress : float
Current effective stress.
pstress : float, optional
Preconsolidation stress. Default pstress=estress i.e. normally consolidated.
Returns: av : float
Slope of void-ratio vs effective stress plot at current stress state.
Examples
On recompression line:
>>> a = CcCrSoilModel(Cc=3.0, Cr=0.5, siga=10, ea=5) >>> a.av_from_stress(estress=40, pstress=50) 0.00542868...
On compression line:
>>> a = CcCrSoilModel(Cc=3.0, Cr=0.5, siga=10, ea=5) >>> a.av_from_stress(estress=60, pstress=50) 0.02171472...
Array inputs:
>>> a = CcCrSoilModel(Cc=3.0, Cr=0.5, siga=10, ea=5) >>> a.av_from_stress(estress=np.array([40, 60.0]), ... pstress=np.array([50, 55.0])) array([ 0.00542868, 0.02171472])
-
e_and_stress_for_plotting
(**kwargs)[source]¶ Void ratio and stress values that plot the model
Parameters: pstress : float, optional
Preconsolidation stress. Default behaviour is normally consolidated.
npts : int, optional
Number of points to return. Default npts=100.
xmin, xmax : float, optional
range of x (i.e. effective stress) values from which to return points. Default xmin=1, xmax=100
Returns: x, y : 1d ndarray
npts stress, and void ratio values between xmin and xmax.
-
e_from_stress
(estress, **kwargs)[source]¶ Void ratio from effective stress
Parameters: estress : float
Current effective stress.
pstress : float, optional
Preconsolidation stress. Default pstress=estress i.e. normally consolidated.
Returns: e : float
Void ratio corresponding to current stress state.
Examples
On recompression line:
>>> a = CcCrSoilModel(Cc=3, Cr=0.5, siga=10, ea=5) >>> a.e_from_stress(estress=40, pstress=50) 2.95154...
On compression line:
>>> a = CcCrSoilModel(Cc=3, Cr=0.5, siga=10, ea=5) >>> a.e_from_stress(estress=60, pstress=50) 2.66554...
Array inputs:
>>> a = CcCrSoilModel(Cc=3, Cr=0.5, siga=10, ea=5) >>> a.e_from_stress(estress=np.array([40, 60]), ... pstress=np.array([50, 55])) array([ 2.95154499, 2.66554625])
Normally consolidated (pstress not specified):
>>> a = CcCrSoilModel(Cc=3.0, Cr=0.5, siga=10, ea=5) >>> a.e_from_stress(estress=10) 5.0
-
stress_from_e
(e, **kwargs)[source]¶ Effective stress from void ratio
Parameters: e : float
Current void ratio.
pstress : float, optional
Preconsolidation stress. Default pstress=estress i.e. normally consolidated.
Returns: estress : float
Effective stress corresponding to current void ratio.
Examples
On recompression line:
>>> a = CcCrSoilModel(Cc=3.0, Cr=0.5, siga=10, ea=5) >>> a.stress_from_e(e=2.95154499, pstress=50) 40...
On compression line:
>>> a = CcCrSoilModel(Cc=3, Cr=0.5, siga=10, ea=5) >>> a.stress_from_e(e=2.66554625, pstress=50) 59.999...
Array inputs:
>>> a = CcCrSoilModel(Cc=3, Cr=0.5, siga=10, ea=5) >>> a.stress_from_e(e=np.array([ 2.95154499, 2.66554625]), pstress=50) array([ 40.0..., 59.99...])
Normally consolidated:
>>> a = CcCrSoilModel(Cc=3.0, Cr=0.5, siga=10, ea=5) >>> a.stress_from_e(e=5) 10.0
-
-
class
geotecha.constitutive_models.void_ratio_stress.
FunctionSoilModel
(fn_e_from_stress, fn_stress_from_e, fn_av_from_stress, *args, **kwargs)[source]¶ Bases:
geotecha.constitutive_models.void_ratio_stress.OneDimensionalVoidRatioEffectiveStress
Functional definition of void-ratio stress relationship
User provides python functions.
Parameters: fn_e_from_stress: callable object
Function to obtain void ratio from stress. fn_e_from_stress should be the inverse function of fn_stress_from_e.
fn_stress_from_e : callable object
Function to obtain stress from void ratio. fn_stress_from_e should be the inverse function of fn_e_from_stress.
fn_av_from_stress: callable object
Function to obtain slope of void ratio-stress relationship. fn_av_from_stress should be negative the derivative of fn_e_from_stress w.r.t. stress be the inverse function of fn_k_from_e.
*args, **kwargs : anything
Positional and keyword arguments to be passed to the fn_e_from_stress, fn_stress_from_e, fn_av_from_stress functions. Note that any additional args and kwargs passed to the functions will be appended to the args, and kwargs. You may get into a mess for example with e_from_stress because normally the first postional argument passed to such fucntions is estress. If you add your own positional arguments, then k will be after you own arguments. Just be aware that the usual way to call methods of the base object OneDimensionalVoidRatioEffectiveStress a single positonal arguement, e.g. estress, e, followed by any required keyword arguments.
Notes
Any function should be able to accept additonal keywords.
Examples
>>> def efs(s, b=2): ... return s * b >>> def sfe(e, b=2): ... return e / b >>> def afs(s, b=2): ... return -b >>> a = FunctionSoilModel(efs, sfe, afs, b=5) >>> a.e_from_stress(3) 15 >>> a.stress_from_e(15) 3.0 >>> a.av_from_stress(3) -5
Prconsolidation stress:
>>> def efs2(s, pstress=None, b=2): ... if pstress is None: ... pstress=8 ... return s * b + pstress >>> a = FunctionSoilModel(efs2, sfe, afs, b=5) >>> a.e_from_stress(3) 23
Methods
av_from_stress
(estress, **kwargs)Slope of void ratio from effective stress e_and_stress_for_plotting
(**kwargs)Void ratio and stress values that plot the method e_from_stress
(estress, **kwargs)Void ratio from effective stress plot_model
(**kwargs)plot the void ratio-stress points stress_from_e
(e, **kwargs)Effective stress from void ratio
-
class
geotecha.constitutive_models.void_ratio_stress.
OneDimensionalVoidRatioEffectiveStress
[source]¶ Bases:
object
Base class for defining 1D void ratio-effective stress relationships
Methods
av_from_stress
(estress, **kwargs)Slope of void ratio from effective stress e_and_stress_for_plotting
(**kwargs)Void ratio and stress values that plot the method e_from_stress
(estress, **kwargs)Void ratio from effective stress plot_model
(**kwargs)plot the void ratio-stress points stress_from_e
(e, **kwargs)Effective stress from void ratio
-
class
geotecha.constitutive_models.void_ratio_stress.
PwiseLinearSoilModel
(siga, ea, Cr, xlog=False, ylog=False)[source]¶ Bases:
geotecha.constitutive_models.void_ratio_stress.OneDimensionalVoidRatioEffectiveStress
Pwise linear void ratio-effective stress realationship
x and y data can be interpolated natural-natural, natural-log10, log10-natural, or log10, log10.
Parameters: siga, ea : 1d array
Effective stress values and void ratio values defining a one-to-one relationship. Slope of void ratio-effectinve stress plot should never fall below Cr.
Cr : float
Precompression index, slope of void rato-effective stress line. Note that Cr is the slope in whatever scales of slog and ylog that have been chosen.
xlog, ylog : True/False, Optional
If True then interpolation on each axis is assumed to be logarithmic with base 10. Default xlog=ylog=False
Methods
av_from_stress
(estress, **kwargs)Slope of void ratio from effective stress e_and_stress_for_plotting
(**kwargs)Void ratio and stress values that plot the model e_from_stress
(estress, **kwargs)Void ratio from effective stress plot_model
(**kwargs)plot the void ratio-stress points stress_from_e
(e, **kwargs)Effective stress from void ratio -
av_from_stress
(estress, **kwargs)[source]¶ Slope of void ratio from effective stress
Parameters: estress : float
Current effective stress.
pstress : float, optional
Preconsolidation stress. Default pstress=estress i.e. normally consolidated.
Returns: av : float
Slope of void-ratio vs effective stress plot at current stress state.
Examples
On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1) >>> a.av_from_stress(estress=1.25, pstress=2.25) 0.10...
On compression line:
>>> a.av_from_stress(estress=2.25, pstress=2) 3.0...
Normally consolidated (pstress not specified):
>>> a.av_from_stress(estress=2.25) 3.0...
Array inputs:
>>> a.av_from_stress(estress=np.array([1.25, 2.25]), ... pstress=np.array([2.25, 2.])) array([ 0.1, 3. ])
Logarithmic effective stress scale On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, xlog=True) >>> a.av_from_stress(estress=1.25, pstress=2.25) 0.034743...
On compression line:
>>> a.av_from_stress(estress=2.25, pstress=2) 2.987...
Normally consolidated (pstress not specified):
>>> a.av_from_stress(estress=2.25) 2.987...
Array inputs:
>>> a.av_from_stress(estress=np.array([1.25, 2.25]), ... pstress=np.array([2.25, 2.])) array([ 0.034743..., 2.987...])
Logarithmic void ratio scale On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, ylog=True) >>> a.av_from_stress(estress=1.25, pstress=2.25) 0.76694...
On compression line:
>>> a.av_from_stress(estress=2.25, pstress=2) 2.9612...
Normally consolidated (pstress not specified):
>>> a.av_from_stress(estress=2.25) 2.9612...
Array inputs:
>>> a.av_from_stress(estress=np.array([1.25, 2.25]), pstress=2.25) array([ 0.76694..., 2.9612...])
Logarithmic effective stress and void ratio scales On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, xlog=True, ylog=True) >>> a.av_from_stress(estress=1.25, pstress=2.25) 0.2210045...
On compression line:
>>> a.av_from_stress(estress=2.25, pstress=2) 2.9034...
Normally consolidated (pstress not specified):
>>> a.av_from_stress(estress=2.25) 2.9034...
Array inputs:
>>> a.av_from_stress(estress=np.array([1.25, 2.25]), pstress=2.25) array([ 0.2210045..., 2.9034...])
Increasing vs decreasing inputs
>>> ea = np.arange(1,10) >>> siga = 3 * ea >>> np.isclose(PwiseLinearSoilModel(siga, ea, Cr=0.1).av_from_stress(7.2), ... PwiseLinearSoilModel(siga[::-1], ea[::-1], Cr=0.1).av_from_stress(7.2)) True
-
e_and_stress_for_plotting
(**kwargs)[source]¶ Void ratio and stress values that plot the model
Parameters: pstress : float, optional
Preconsolidation stress. Default behaviour is normally consolidated.
npts : int, optional
Number of points to return. Default npts=100
xmin, xmax : float, optional
Range of x (i.e. effective stress) values from which to return points. Default minumum of model siga.
Returns: x, y : 1d ndarray
npts stress, and void ratio values between xmin and xmax.
-
e_from_stress
(estress, **kwargs)[source]¶ Void ratio from effective stress
Parameters: estress : float
Current effective stress. etress must be within the range of the soil model points.
pstress : float, optional
Preconsolidation stress. Default pstress=estress i.e. normally consolidated. pstress must be in the range of the soil model points.
Returns: e : float
Void ratio corresponding to current stress state.
Examples
On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.2) >>> a.e_from_stress(estress=1.25, pstress=2.25) 2.95...
On compression line:
>>> a.e_from_stress(estress=2.25, pstress=2) 2.75...
Normally consolidated (pstress not specified):
>>> a.e_from_stress(estress=2.25) 2.75...
Array inputs:
>>> a.e_from_stress(estress=np.array([1.25, 2.25]), pstress=2.25) array([ 2.95, 2.75])
Logarithmic effective stress scale: On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.2, xlog=True) >>> a.e_from_stress(estress=1.25, pstress=2.25) 2.75930...
On compression line:
>>> a.e_from_stress(estress=2.25, pstress=2) 2.70824...
Normally consolidated (pstress not specified):
>>> a.e_from_stress(estress=2.25) 2.70824...
Array inputs:
>>> a.e_from_stress(estress=np.array([1.25, 2.25]), pstress=2.25) array([ 2.75930..., 2.70824...])
Logarithmic void ratio scale On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, ylog=True) >>> a.e_from_stress(estress=1.25, pstress=2.25) 3.330803...
On compression line:
>>> a.e_from_stress(estress=2.25, pstress=2) 2.64575...
Normally consolidated (pstress not specified):
>>> a.e_from_stress(estress=2.25) 2.64575...
Array inputs:
>>> a.e_from_stress(estress=np.array([1.25, 2.25]), pstress=2.25) array([ 3.330803..., 2.64575...])
Logarithmic effective stress and void ratio scales On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, xlog=True, ylog=True) >>> a.e_from_stress(estress=1.25, pstress=2.25) 2.76255...
On compression line:
>>> a.e_from_stress(estress=2.25, pstress=2) 2.604857...
Normally consolidated (pstress not specified):
>>> a.e_from_stress(estress=2.25) 2.604857...
Array inputs:
>>> a.e_from_stress(estress=np.array([1.25, 2.25]), pstress=2.25) array([ 2.76255..., 2.604857...])
Increasing vs decreasing inputs
>>> ea = np.arange(1,10) >>> siga = 3 * ea >>> np.isclose(PwiseLinearSoilModel(siga, ea, ... Cr=0.1).e_from_stress(7.2), ... PwiseLinearSoilModel(siga[::-1], ea[::-1], ... Cr=0.1).e_from_stress(7.2)) True
-
stress_from_e
(e, **kwargs)[source]¶ Effective stress from void ratio
Parameters: e : float
Current void ratio.
pstress : float, optional
Preconsolidation stress. Default pstress=estress i.e. normally consolidated.
Returns: estress : float
Effective stress corresponding to current void ratio.
Examples
On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1) >>> a.stress_from_e(e=2.8, pstress=2.25) 1.75...
On compression line:
>>> a.stress_from_e(e=2.65, pstress=2) 2.28333...
Normally consolidated (pstress not specified): >>> a.stress_from_e(e=2.65) 2.28333...
Array inputs:
>>> a.stress_from_e(e=np.array([2.8, 2.65]), ... pstress=np.array([2.25, 2.0])) array([ 1.75..., 2.28333...])
Logarithmic effective stress scale On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, xlog=True) >>> a.stress_from_e(e=2.73, pstress=2.25) 1.363494...
On compression line:
>>> a.stress_from_e(e=2.73, pstress=2) 2.2427...
Normally consolidated (pstress not specified):
>>> a.stress_from_e(e=2.73) 2.2427...
Array inputs:
>>> a.stress_from_e(e=2.73, pstress=np.array([2.25, 2.])) array([ 1.363494..., 2.2427...])
Logarithmic void ratio scale On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, ylog=True) >>> a.stress_from_e(e=2.75, pstress=2.25) 2.082163...
On compression line:
>>> a.stress_from_e(e=2.65, pstress=2) 2.24856...
Normally consolidated (pstress not specified):
>>> a.stress_from_e(e=2.65) 2.24856...
Array inputs:
>>> a.stress_from_e(e=np.array([2.75, 2.65]), ... pstress=np.array([2.25, 2])) array([ 2.082163..., 2.24856...])
Logarithmic effective stress and void ratio scales On recompression line:
>>> x, y = np.array([1,2,2.5]), np.array([4, 3.5, 2]) >>> a = PwiseLinearSoilModel(siga=x, ea=y, Cr=0.1, xlog=True, ylog=True) >>> a.stress_from_e(e=2.65, pstress=2.25) 1.8948013...
On compression line:
>>> a.stress_from_e(e=2.65, pstress=2) 2.23463...
Normally consolidated (pstress not specified):
>>> a.stress_from_e(e=2.65) 2.23463...
Array inputs:
>>> a.stress_from_e(e=2.65, pstress=np.array([2.25, 2])) array([ 1.8948013..., 2.23463...])
Increasing vs decreasing inputs
>>> ea = np.arange(1,10) >>> siga = 3 * ea >>> np.isclose(PwiseLinearSoilModel(siga, ea, Cr=0.1).stress_from_e(3.0, pstress=4.0), ... PwiseLinearSoilModel(siga[::-1], ea[::-1], Cr=0.1).stress_from_e(3.0, pstress=4.0)) True
-
-
class
geotecha.constitutive_models.void_ratio_stress.
YinAndGrahamSoilModel
(lam, kap, psi, siga, ea, ta, e0=None, estress0=None, pstress0=None, t0=None)[source]¶ Bases:
geotecha.constitutive_models.void_ratio_stress.OneDimensionalVoidRatioEffectiveStress
Yin and Graham creep void ratio-effective stress realationship
This model uses the integral void ratio form of the Yin and Graham (1996) [R95] soil model as described by Hu et al. (2014) [R96]. Basically you step through time keeping track of the stress state and a cumulatlive integral representing creep strain.
Parameters: lam : float
Compressibility index, slope of the referece time line e-ln(sig) line.
kap : float
Recompression index, slope of e-ln(sig) line.
psi : float
Creep parameter. Slope of the void ratio-ln(time) line.
siga, ea : float
Point on reference time line fixing it in effective stress-void ratio space.
ta : float
Reference time corresponding to the reference time line.
e0 : float, optional
Initial void ratio. Default e0=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
estress0 : float, optional
Initial effective stress. Default estress=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
pstress0 : float, optional
Initial ‘preconsolidation` stress on the reference time line. Default pstress0=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
t0 : float, optional
Initial equivalent time at the initial conditon. Default t0=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
Notes
The void ratio stress relationship is:
\[e=e0-\kappa\ln\left({\frac{\sigma^{\prime}_z}{\sigma^{\prime}_{z}}}\right) -\psi\ln\left[{ \int_{0}^{t}{ \frac{1}{t_0} \left({ \frac{\sigma^{\prime}_z} {\sigma^{\prime}_{0}}}\right)^{\alpha}\,dt }+1 }\right]\]To calculate the integral we keep track of the effective stress and the integration at time t and then add the integration increment between \(t\) and \(t+\bigtriangleup t\):
\[e=e0-\kappa\ln\left({\frac{\sigma^{\prime}_{z,t+\bigtriangleup t}}{\sigma^{\prime}_{z}}}\right) -\psi\ln\left[{ \int_{0}^{t}{ \frac{1}{t_0} \left({ \frac{\sigma^{\prime}_z} {\sigma^{\prime}_{0}}}\right)^{\alpha}\,dt } + \int_{t}^{t+\bigtriangleup t}{ \frac{1}{t_0} \left({ \frac{\sigma^{\prime}_z} {\sigma^{\prime}_{0}}}\right)^{\alpha}\,dt } +1}\right]\]Consider a time increment of length \(\bigtriangleup t\) in which the effective stress changes from \(\sigma^{\prime}_{t}\) to \(\sigma^{\prime}_{t + \bigtriangleup t}\). Assume that the stress changes in a linear fashion over \(\bigtriangleup t\):
\[\sigma^{\prime}\left({t}\right) = \sigma^{\prime}_{t} + t \frac{\sigma^{\prime}_{t+ \bigtriangleup t} -\sigma^{\prime}_{t}} {\bigtriangleup t}\]The time integral becomes:
\[\frac{1}{t_0}\int_{0}^{t}{ \left({ \frac{\sigma^{\prime}_z} {\sigma^{\prime}_{z0}}}\right)^{\alpha}\,dt } + \frac{1}{t_0}\int_{0}^{\bigtriangleup t} {\left({\frac{\sigma^{\prime}_{t} + t \frac{\sigma^{\prime}_{t+ \bigtriangleup t} -\sigma^{\prime}_{t}} {\bigtriangleup t}}{\sigma^{\prime}_0}}\right)^\alpha\,dt}\]The increment that is added to the main integral is:
\[\frac{1}{t_0}\int_{0}^{\bigtriangleup t} {\left({\frac{\sigma^{\prime}_{t} + t \frac{\sigma^{\prime}_{t+ \bigtriangleup t} -\sigma^{\prime}_{t}} {\bigtriangleup t}}{\sigma^{\prime}_0}}\right)^\alpha\,dt} = \frac{1}{t_0} \frac{\sigma^{\prime}_{0} \bigtriangleup t} {\left({\alpha+1}\right) \left({\sigma^{\prime}_{t+ \bigtriangleup t} -\sigma^{\prime}_{t}}\right)} \left({ \left({\frac{\sigma^{\prime}_{t+ \bigtriangleup t}}{\sigma^{\prime}_0}}\right)^{\alpha+1}- \left({\frac{\sigma^{\prime}_{t}}{\sigma^{\prime}_0}}\right)^{\alpha+1} }\right)\]When \(\sigma^{\prime}_{t+ \bigtriangleup t}=\sigma^{\prime}_{t}\) then the increment is:
\[\frac{1}{t_0}\int_{0}^{\bigtriangleup t} {\left({\frac{\sigma^{\prime}_{t} + t \frac{\sigma^{\prime}_{t+ \bigtriangleup t} -\sigma^{\prime}_{t}} {\bigtriangleup t}}{\sigma^{\prime}_0}}\right)^\alpha\,dt} = \frac{1}{t_0}\bigtriangleup t \left({\frac{\sigma^{\prime}_{t}}{\sigma^{\prime}_0}}\right)^{\alpha+1}\]The value of the integral and the effective stress is stored within the YinAndGrahamSoilModel object and updated after each call to e_from_stress or stress_from_e.
e_from_stress
e_from_stress uses the above equations directly (knowing all material parameters, the value of the integral at the start of the increment, and the efffective stresses at the beginning and end of time incremetn dt, then e can be calculated directly).
stress_from_e
stress_from_e recasts the equations as estress=f(estress) and uses fixed point iteration to determine the effective stress at the end of the increment.
initial_conditions
If two of the following four parameters describing the intial state of the soil are known, then the equations below can be combined to solve for the remaining two unknowns: initial void ratio \(e_0\), initial stress \(\sigma^{\prime}_0\), pseudo initial preconsolidation stress \(\sigma^{\prime}_p\) and equivalent time corresponding to initial state \(t_0\).
On the reference time line running through point \((\sigma^{\prime}_a, e_a)\), the psuedo preconsolidation pressure \(\sigma^{\prime}_p\) corresponds to void ratio \(e_p\):
\[e_p = e_a - \lambda \ln\left({ \frac{\sigma^{\prime}_p} {\sigma^{\prime}_a} }\right)\]The instant time line running through point \((\sigma^{\prime}_0, e_0)\) meets the reference time line at \(\sigma^{\prime}_p\) and void ratio \(e_p\):
\[e_p = e_0 - \kappa \ln\left({ \frac{\sigma^{\prime}_p} {\sigma^{\prime}_0} }\right)\]The initial state \((\sigma^{\prime}_0, e_0)\) can be reached from the reference time line by either a) loading to \(\sigma^{\prime}_p\) and then unloading to \(\sigma^{\prime}_0\) or b) creeping from time \(t_a\) to time \(t_0\) at constant stress \(\sigma^{\prime}_0\). equating the two void ratio changes gives:
\[\left({\lambda-\kappa}\right) \ln\left({ \frac{\sigma^{\prime}_p}{\sigma^{\prime}_0} }\right) = \psi \ln\left({ \frac{t_0}{t_a} }\right)\]YinAndGrahamSoilModel()
References
[R95] (1, 2) Yin, J.-H., and Graham, J. (1996). ‘Elastic visco-plastic modelling of one-dimensional consolidation. Geotechnique, 46(3), 515-527. [R96] (1, 2) Hu, Y.-Y., Zhou, W.-H., and Cai, Y.-Q. (2014). ‘Large-strain elastic viscoplastic consolidation analysis of very soft clay layers with vertical drains under preloading. Canadian Geotechnical Journal, 51(2), 144-157. Examples
Combos of initial conditions give the same results
>>> fixed_param = dict(lam=0.2, kap=0.04, psi=0.01, siga=20, ea=1, ta=1) >>> oparam = dict(e0=0.95, estress0=18, pstress0=28.06637954, ... t0=1220.73731) >>> a = YinAndGrahamSoilModel(e0=oparam['e0'], ... estress0=oparam['estress0'], **fixed_param) >>> [a.e0, a.estress0, a.pstress0, a.t0] [0.95..., 18..., 28.0..., 1220.737...] >>> a = YinAndGrahamSoilModel(e0=oparam['e0'], ... pstress0=oparam['pstress0'], **fixed_param) >>> [a.e0, a.estress0, a.pstress0, a.t0] [0.95..., 18..., 28.0..., 1220.737...] >>> a = YinAndGrahamSoilModel(e0=oparam['e0'], ... t0=oparam['t0'], **fixed_param) >>> [a.e0, a.estress0, a.pstress0, a.t0] [0.95..., 18..., 28.0..., 1220.737...] >>> a = YinAndGrahamSoilModel(estress0=oparam['estress0'], ... pstress0=oparam['pstress0'], **fixed_param) >>> [a.e0, a.estress0, a.pstress0, a.t0] [0.95..., 18..., 28.0..., 1220.737...] >>> a = YinAndGrahamSoilModel(estress0=oparam['estress0'], ... t0=oparam['t0'], **fixed_param) >>> [a.e0, a.estress0, a.pstress0, a.t0] [0.95..., 18..., 28.0..., 1220.737...] >>> a = YinAndGrahamSoilModel(pstress0=oparam['pstress0'], ... t0=oparam['t0'], **fixed_param) >>> [a.e0, a.estress0, a.pstress0, a.t0] [0.95..., 18..., 28.0..., 1220.737...]
Attributes
_alpha (float) Composite material parameter (lam - kap / psi) _psi_zero_model (CcCrSoilModel) soil model when psi=0. i.e. CcCrSoilModel following instant time line and reference time line. _igral (float) Current value of the material model integration w.r.t. time. _igral = integrate[1/t0 * (estress/etress0)**alpha, t,0,t] _inc (float) Increment to add to _igral. Helper array to speed up stress_from_e. Only present after initial conditions have been initialized. _inc = integrate[1/t0 * (estress/etress0)**alpha, t,(t),(t + dt)] _es (float) Temp effectve stress. Helper array to speed up stress_from_e calc. Only present after initial conditions have been initialized. _estress (float) Current value of effective stress. Only present after initial conditions have been initialized. Basically when you call method e_from_stress or method stress_from_e with parameter dt, then _estress and _igral correspond to effective stress and the time integral at the start of the increment; the integration increment _inc will be calculated and then _igral and _estress will be updated to reflect values at the end of dt time increment. Methods
CRSN
(tvals, edot, **kwargs)Constant rate of strain simulation CRSS
(tvals, estressdot, **kwargs)Constant rate of stress simulation av_from_stress
(estress, **kwargs)Slope of void ratio from effective stress e_and_stress_for_plotting
(**kwargs)Void ratio and stress values that plot the method e_from_stress
(estress, **kwargs)Void ratio from effective stress initial_conditions
([e0, estress0, pstress0, t0])Calculate the initial conditions plot_model
(**kwargs)plot the void ratio-stress points stress_from_e
(e, **kwargs)Effective stress from void ratio -
CRSN
(tvals, edot, **kwargs)[source]¶ Constant rate of strain simulation
Note that the rate of chage of void ratio is specified. Strain can be back caluclated using strain = (e0-e) / (1 + e0).
Parameters: tvals : 1d array of float
Time values for piecewise linear strain-rate vs time.
edot : 1d array of float
Time rate of change of void ratio corresponding to tvals. Note that you will need a negative edot for void ratio to decrease. Be careful because strain in geotech is usually positive for a decrease in void ratio.
method : [‘step’, ‘ode’] optional
Describes which method to use. ‘step’ will step through time using the stress_from_e function. ‘odeint’ will reframe the constitutive model in terms of stress-rate and use odeint to solve for effective stress. Both methods should produce the same results; the ‘odeint’ method is really just a validity check of the the ‘step’ method. Default method=’step’.
**kwargs : various
Arguments (other than method) will be passed to the subdivide_x_into_segments function. If not given then the following values will be used, dx = (tvals[-1]-tvals[0])/200, min_segments=50. You may have to fiddle around with kwargs to get meaningful results; a small dx value will give good results but could take a long time to compute.
Returns: tvals_ : 1d array of float
Time values
estress : 1d array of float
Effective stress at tvals_.
e : 1d array of float
Void ratio at time tvals_.
edot_ : 1d array of float.
Time rate of change of void ratio at tvals_
See also
geotecha.piecewise.piecewsie_linear_1d.subdivide_x_into_segments
- used to determine time intervals.
-
CRSS
(tvals, estressdot, **kwargs)[source]¶ Constant rate of stress simulation
Note that the output is void ratio not strain. Strain can be back caluclated using strain = (e0-e) / (1 + e0).
Parameters: tvals : 1d array of float
Time values for piecewise linear strain-rate vs time.
estressdot : 1d array of float
Time rate of change of effective stress corresponding to tvals.
method : [‘step’, ‘ode’] optional
Describes which method to use. ‘step’ will step through time using the e_from_stress function. ‘odeint’ will reframe the constitutive model in terms of strain-rate and use odeint to solve for effective stress. Both methods should produce the same results; the ‘odeint’ method is really just a validity check of the the ‘step’ method. Default method=’step’.
**kwargs : various
Arguments (other than method) will be passed to the subdivide_x_into_segments function. If not given then the following values will be used, dx = (tvals[-1]-tvals[0])/200, min_segments=50. You may have to fiddle around with kwargs to get meaningful results; a small dx value will give good results but could take a long time to compute.
Returns: tvals_ : 1d array of float
Time values
estress : 1d array of float
Effective stress at tvals_.
e : 1d array of float
Void ratio at time tvals_.
estressdot_ : 1d array of float.
Time rate of change of stress at tvals_
See also
geotecha.piecewise.piecewsie_linear_1d.subdivide_x_into_segments
- used to determine time intervals.
-
av_from_stress
(estress, **kwargs)[source]¶ Slope of void ratio from effective stress
When estressdot is provided as a kwarg then av is a true representation of de/destress at the current time (i.e. using the current value of self._igral). When estressdot is not provided then a conservative lower value of av based on self.kap and the current estress is calculated.
Parameters: estress : float
Current effective stress.
estressdot: float, optional
Rate of change of effective stress w.r.t. time. Default estressdot=None i.e. av is calculated using self.kap and current effective stress.
Returns: av : float
Slope of void-ratio vs effective stress plot at current stress state, assuming psi=0.
Examples
>>> a = YinAndGrahamSoilModel(lam=0.2, kap=0.04, psi=0.01, siga=20, ... ea=1, ta=1, e0=0.90, estress0=18) >>> a._igral=3.1 #artificially set igral. Usually calculted advancing through time using e_from_stress or stress_from_e >>> a.av_from_stress(estress=10) 0.004 >>> a.av_from_stress(estress=40, estressdot=1.5) 0.00417...
#todo: write up av equations
-
e_from_stress
(estress, **kwargs)[source]¶ Void ratio from effective stress
When dt is None then then self.e0, self.estress0, self.t0 and self.pstress0 will be initialized from estress and pstress; self.e0 will be returned.
Parameters: estress : float
Current effective stress. estress is assumed to be the effective stress at the end of an increment of length dt. i.e. at time t+dt.
dt : float, optional
Time increment. Default dt=None. When dt is None, the void ratio is calculted assuming psi=0. i.e. along the instant time line and then down the reference time line.
pstress : float, optional
Preconsolidation stress. This is only relevant if dt is not None. Default pstress=estress i.e. on the reference time line.
Returns: e : float
Void ratio corresponding to current stress state.
Notes
self._igral and self._estress and self._inc will be updated!
-
initial_conditions
(e0=None, estress0=None, pstress0=None, t0=None)[source]¶ Calculate the initial conditions
Only use two of the paramters, the other two will be calculated.
Parameters: e0 : float, optional
Initial void ratio. Default e0=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
estress0 : float, optional
Initial effective stress. Default estress=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
pstress0 : float, optional
Initial ‘preconsolidation` stress on the reference time line. Default pstress0=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
t0 : float, optional
Initial equivalent time at the initial conditon. Default t0=None. If specifying initial conditions then only use two of the four available (e0, estress0, pstress0, t0).
Returns: e0 : float
Initial void ratio
estress0 : float
Initial effective stress
pstress0 : float
Inital preconsolidation stress
t0 : float
Initial equivalent time.
-
stress_from_e
(e, **kwargs)[source]¶ Effective stress from void ratio
If stress_from_e is called when self.estress0 is not initialized then the initial conditions will be initialized assumming the stress state is on the reference time line.
Parameters: e : float
Current void ratio. e is assumed to be the effective stress at the end of an increment of length dt. i.e. at time t+dt.
dt : float, optional
Time increment. Default dt=0.
Returns: estress : float
Effective stress corresponding to current void ratio.
Notes
self._igral and self._estress and self._inc will be updated!
-