rockphypy.EM
¶
Module Contents¶
Classes¶
classical bounds and inclusion models. |
- class rockphypy.EM.EM[source]¶
classical bounds and inclusion models.
- static cripor(K0, G0, phi, phic)[source]¶
Critical porosity model according to Nur’s modified Voigt average.
- static cripor_reuss(M0, Mf, phic, den=False)[source]¶
In the suspension domain, the effective bulk and shear moduli of the rock can be estimated by using the Reuss (isostress) average.
- Parameters:
M0 (float) – The solid phase modulus or density
Mf (float) – The pore filled phase modulus or density
phic (float) – critical porosity
den (bool, optional) – If False: compute the reuss average for effective modulus of two mixing phases. If true, compute avearge density using mass balance, which corresponds to voigt average. Defaults to False.
- Returns:
float or array-like – M (GPa/g.cc): average modulus or average density
References
Section 7.1 Rock physics handbook 2nd edition
- static HS(f, K1, K2, G1, G2, bound='upper')[source]¶
Compute effective moduli of two-phase composite using hashin-strikmann bounds.
- Parameters:
f (float) – 0-1, volume fraction of stiff material
K1 (float or array-like) – bulk modulus of stiff phase
K2 (float or array-like) – bulk modulus of soft phase
G1 (float or array-like) – shear modulus of stiff phase
G2 (float or array-like) – shear modulus of soft phase
bound (str, optional) – upper bound or lower bound. Defaults to ‘upper’.
- Returns:
float or array-like – K, G (GPa): effective moduli of two-phase composite
- static Eshelby_Cheng(K, G, phi, alpha, Kf, mat=False)[source]¶
Compute the effective anisotropic moduli of a cracked isotropic rock with single set fracture using Eshelby–Cheng Model.
- Parameters:
K (float) – bulk modulus of the isotropic matrix GPa
G (float) – shear modulus of the isotropic matrix GPa
phi (float) – (crack) porosity
alpha (float) – aspect ratio of crack
Kf (float) – bulk modulus of the fluid. For dry cracks use fluid bulk modulus 0
mat (bool, optional) – If true: the output is in matrix form, otherwise is numpy array. Defaults to False.
- Returns:
1d or 2d array – C_eff: effective moduli of cracked, transversely isotropic rocks
References
section 4.14 in The Rock Physics Handbook
- static hudson(K, G, Ki, Gi, alpha, crd, order=1, axis=3)[source]¶
Hudson’s effective crack model assuming weak inclusion for media with single crack set with all normals aligned along 1 or 3-axis. First and Second order corrections are both implemented. Notice that the second order correction has limitation. See Cheng (1993).
- Parameters:
K (float) – bulk modulus of isotropic background
G (float) – shear modulus of isotropic background
Ki (float) – bulk modulus of the inclusion material. For dry cracks: Ki=0
Gi (float) – shear modulus of the inclusion material
alpha (float) – crack aspect ratio
crd (float) – crack density
order (int, optional) –
- approximation order.
1: Hudson’s model with first order correction. 2: Hudson’s model with first order correction. Defaults to 1.
axis (int, optional) –
- axis of symmetry.
1: crack normals aligned along 1-axis, output HTI 3: crack normals aligned along 3-axis, output VTI Defaults to 3
- Returns:
2d array – C_eff: effective moduli in 6x6 matrix form.
- static hudson_rand(K, G, Ki, Gi, alpha, crd)[source]¶
Hudson’s crack model of a material containing randomly oriented inclusions. The model results agree with the consistent results of Budiansky and O’Connell (1976).
- Parameters:
K (float or array-like) – bulk modulus of isotropic background
G (float or array-like) – shear modulus of isotropic background
Ki (float) – bulk modulus of the inclusion material. For dry cracks: Ki=0
Gi (float) – shear modulus of the inclusion material, for fluid, Gi=0
alpha (float) – crack aspect ratio
crd (float) – crack density
- Returns:
float or array-like – K_eff, G_eff (GPa): effective moduli of the medium with randomly oriented inclusions
- static hudson_ortho(K, G, Ki, Gi, alpha, crd)[source]¶
Hudson’s first order effective crack model assuming weak inclusion for media with three crack sets with normals aligned along 1 2, and 3-axis respectively. Model is valid for small crack density and aspect ratios.
- Parameters:
K (float) – bulk modulus of isotropic background
G (float) – shear modulus of isotropic background
Ki (float) – bulk modulus of the inclusion material. For dry cracks: Ki=0
Gi (float) – shear modulus of the inclusion material, for fluid, Gi=0
alpha (nd array with size 3) – [alpha1, alpha2,alpha3] aspect ratios of three crack sets
crd (nd array with size 3) – [crd1, crd2, crd3] crack densities of three crack sets
- Returns:
2d array – C_eff: effective moduli in 6x6 matrix form.
- static hudson_cone(K, G, Ki, Gi, alpha, crd, theta)[source]¶
Hudson’s first order effective crack model assuming weak inclusion for media with crack normals randomly distributed at a fixed angle from the TI symmetry axis 3 forming a cone;
- Parameters:
K (float) – bulk modulus of isotropic background
G (float) – shear modulus of isotropic background
Ki (float) – bulk modulus of the inclusion material. For dry cracks: Ki=0
Gi (float) – shear modulus of the inclusion material, for fluid, Gi=0
alpha (float) – aspect ratios of crack sets
crd (float) – total crack density
theta (float) – the fixed angle between the crack normam and the symmetry axis x3. degree unit.
- Returns:
2d array – C_eff: effective moduli of TI medium in 6x6 matrix form.
- static Berryman_sc(K, G, X, Alpha)[source]¶
Effective elastic moduli for multi-component composite using Berryman’s Consistent (Coherent Potential Approximation) method.See also: PQ_vectorize, Berryman_func
- Parameters:
K (array-like) – 1d array of bulk moduli of N constituent phases, [K1,K2,…Kn]
G (array-like) – 1d array of shear moduli of N constituent phases, [G1,G2,…Gn]
X (array-like) – 1d array of volume fractions of N constituent phases, [x1,…xn], Sum(X) = 1.
Alpha (array-like) – aspect ratios of N constituent phases. Note that α <1 for oblate spheroids and α > 1 for prolate spheroids, α = 1 for spherical pores,[α1,α2…αn]
- Returns:
array-like – K_sc,G_sc: Effective bulk and shear moduli of the composite
- static PQ_vectorize(Km, Gm, Ki, Gi, alpha)[source]¶
compute geometric strain concentration factors P and Q for prolate and oblate spheroids according to Berymann (1980).See also: Berryman_sc, Berryman_func
- Parameters:
Km (float) – bulk modulus of matrix phase. For Berryman SC approach, this corresponds to the effective moduli of the composite.
Gm (float) – shear modulus of matrix phase. For Berryman SC approach, this corresponds to the effective moduli of the composite.
Ki (array-like) – 1d array of bulk moduli of N constituent phases, [K1,K2,…Kn]
Gi (array-like) – 1d array of shear moduli of N constituent phases, [G1,G2,…Gn]
alpha (array-like) – aspect ratios of N constituent phases. Note that α <1 for oblate spheroids and α > 1 for prolate spheroids, α = 1 for spherical pores,[α1,α2…αn]
- Returns:
array-like – P,Q (array): geometric strain concentration factors, [P1,,,Pn],[Q1,,,Qn]
- static Berryman_func(params, K, G, X, Alpha)[source]¶
Form the system of equastions to solve. See 4.11.14 and 4.11.15 in Rock physics handbook 2020. See also: Berryman_sc
- Parameters:
params – Parameters to solve, K_sc, G_sc
K (array) – 1d array of bulk moduli of N constituent phases, [K1,K2,…Kn]
G (array) – 1d array of shear moduli of N constituent phases, [G1,G2,…Gn]
X (array) – 1d array of volume fractions of N constituent phases, [x1,…xn]
Alpha (array) – aspect ratios of N constituent phases. Note that α <1 for oblate spheroids and α > 1 for prolate spheroids, α = 1 for spherical pores,[α1,α2…αn]
- Returns:
equation – Eqs to be solved
- static Swiss_cheese(Ks, Gs, phi)[source]¶
Compute effective elastic moduli via “Swiss cheese” model with spherical pores. “Swiss cheese” model assumes a dilute distribution of spherical inclusions embedded in an * unbounded * homogenous solid. It takes the “noninteracting assumption” in which all cavities (pores) are independent so that their contributions can be added.
- static SC(phi, Ks, Gs, iter_n)[source]¶
Self-Consistent(SC) model with spherical pores considering the critical porosity and the interaction effect between inclusions.
- Parameters:
- Returns:
float or array-like – K_eff,G_eff (GPa): effective elastic moduli
- static OConnell_Budiansky(K0, G0, crd)[source]¶
O’Connell and Budiansky (1974) presented equations for effective bulk and shear moduli of a cracked medium with randomly oriented dry penny-shaped cracks (in the limiting case when the aspect ratio α goes to 0)
- static OConnell_Budiansky_fl(K0, G0, Kfl, crd, alpha)[source]¶
Saturated effective elastic moduli using the O’Connell and Budiansky Consistent (SC) formulations under the constraints of small aspect ratio cracks with soft-fluid saturation.
- Parameters:
- Returns:
float – K_sat,G_sat: elastic moduli of cracked background fully saturated by soft fluid.
References – ———-
- O’Connell and Budiansky, (1974)
- static OC_R_funcs(params, crd, nu_0, w)[source]¶
Form the system of equastions to solve. Given crack density and w, solve for the D and nu_eff simulaneously using equations 23 and 25 in O’Connell and Budiansky, (1974)
- static PQ(Km, Gm, Ki, Gi, alpha)[source]¶
compute geometric strain concentration factors P and Q for prolate and oblate spheroids according to Berymann (1980). See also PQ_vectorize
- Parameters:
- Returns:
float – P,Q (unitless): geometric strain concentration factors
- static DEM(y, t, params)[source]¶
ODE solver tutorial: https://physics.nyu.edu/pine/pymanual/html/chap9/chap9_scipy.html.
- static Berryman_DEM(Km, Gm, Ki, Gi, alpha, phi)[source]¶
Compute elastic moduli of two-phase composites by incrementally adding inclusions of one phase (phase 2) to the matrix phase using Berryman DEM theory
- static SC_dilute(Km, Gm, Ki, Gi, f, mode)[source]¶
Elastic solids with elastic micro-inclusions. Random distribution of dilute spherical micro-inclusions in a two phase composite.
- Parameters:
Km (float) – bulk modulus of matrix
Gm (float) – shear modulus of matrix
Ki (float) – bulk modulus of inclusion
Gi (float) – shear modulus of inclusion
f (float or array) – volume fraction of inclusion phases
mode (string) – ‘stress’ if macro stress is prescribed. ‘strain’ if macro strain is prescribed.
References
Nemat-Nasser and M. Hori (book) : Micromechanics: Overall Properties of Heterogeneous Materials. Sec 8
- Returns:
float or array – K, G: effective moduli of the composite
- static SC_flex(f, iter_n, Km, Ki, Gm, Gi)[source]¶
iteratively solving self consistent model for a two phase compposite consisting random distribution of spherical inclusion, not limited to pore.
- Parameters:
f (float or array) – volumetric fraction, f.shape== Km.shape
iter_n (int) – iterations, necessary iterations increases as f increases.
Km (float) – bulk modulus of matrix phase
Ki (float) – bulk modulus of inclusion phase
Gm (float) – shear modulus of matrix phase
Gi (float) – shear modulus of inclusion phase
Reference –
--------- –
(book) (S. Nemat-Nasser and M. Hori) –
- Returns:
float or array – K_eff, G_eff (GPa): Effective elastic moduli
- static MT_average(f, Kmat, Gmat, K1, G1, K2, G2)[source]¶
Compute Two-phase composite without matrix using modified Mori-Takana Scheme according to Iwakuma 2003, one of the inhomogeneities must be considered as a matrix in the limiting model.
- Parameters:
f (float or array) – Volume fraction of matrix/inhomogeneity 1. f1=1-f2, (1-f) can be regarded as pseudo crack density.
Kmat (float) – Bulk modulus of matrix/inhomogeneity 1
Gmat (float) – shear modulus of matrix/ inhomogeneity 1
K1 (float) – Bulk modulus of inhomogeneity 1
G1 (float) – shear modulus of inhomogeneity 1
K2 (float) – Bulk modulus of inhomogeneity 2
G2 (float) – shear modulus of inhomogeneity 2
- Returns:
float or array – K_ave, G_ave [GPa]: MT average bulk and shear modulus