pyavo.avotools package¶
Submodules¶
pyavo.avotools.approx module¶
Functions to compute Aki-Richards and Shuey 2-term & 3-term approximations.
- pyavo.avotools.approx.aki_richards(vp1: numpy.ndarray, vs1: numpy.ndarray, rho1: numpy.ndarray, vp2: numpy.ndarray, vs2: numpy.ndarray, rho2: numpy.ndarray, theta1: numpy.ndarray) → numpy.ndarray¶
Computes the Reflection Coefficient with Aki and Richard’s (1980) equation for a two-layered model.
- Parameters
vp1 – P-wave velocity in upper layer
vs1 – S-wave velocity in upper layer
rho1 – Density of upper layer
vp2 – P-wave velocity in lower layer
vs2 – S-wave velocity in lower layer
rho2 – Density in lower layer
theta1 – Angle of incidence
- Returns
rc: Reflection coefficient
Reference: AVO - Chopra and Castagna, 2014, Page 62.
- pyavo.avotools.approx.ref_coeff(imp: Union[numpy.ndarray, pandas.core.series.Series])¶
Computes the reflection coefficient for a plane incident P-wave.
- Parameters
imp – Acoustic Impedance
- Returns
rc: reflection coefficient
- pyavo.avotools.approx.shuey(vp1: numpy.ndarray, vs1: numpy.ndarray, rho1: numpy.ndarray, vp2: numpy.ndarray, vs2: numpy.ndarray, rho2: numpy.ndarray, theta1: numpy.ndarray) → tuple¶
Computes the Reflectiviy parameters with Shuey (1985) 2 and 3 terms for a two-layered model.
- Parameters
vp1 – P-wave velocity in upper layer
vs1 – S-wave velocity in lower layer
rho1 – Density in upper later
vp2 – P-wave velocity in lower layer
vs2 – S-wave velocity in lower layer
rho2 – Density in lower layer
theta1 – Angles of incidence
- Returns
c : Intercept. m : Gradient. rc2 : Reflection coefficient for the 2-term approximation. rc3 : Reflection coefficient for the 3-term approximation.
Reference: Avseth et al., Quantitative seismic interpretation, 2006, Page 182.
- pyavo.avotools.approx.shueyrc(vp0: Union[numpy.ndarray, pandas.core.series.Series, float], vs0: Union[numpy.ndarray, pandas.core.series.Series, float], rho0: Union[numpy.ndarray, pandas.core.series.Series, float], theta1: Union[numpy.ndarray, pandas.core.series.Series, float]) → tuple¶
Computes the P-wave reflectivity with Shuey (1985) 2 terms for a given well log.
- Parameters
vp0 – P-wave velocities from Vp log
vs0 – S-Wave velocities from Vs log
rho0 – Density from RHOB log
theta1 – Angles of incidence
- Returns
RC: Reflection coefficient for the 2-term approximation. c: Intercept. m: Gradient.
Reference: Avseth et al., Quantitative seismic interpretation, 2006, Page 182.
- pyavo.avotools.approx.snell(vp1: numpy.ndarray, vp2: numpy.ndarray, theta1: float) → tuple¶
Computes the angles of refraction for an incident P-wave in a two-layered model.
- Parameters
vp1 – P-wave velocity in upper layer
vp2 – P-wave velocity in lower layer
theta1 – Angle of incidence
- Returns
theta2: Angle of Refraction p: Ray Parameter
Reference: AVO - Chopra and Castagna, 2014, Page 6.
pyavo.avotools.gassmann module¶
Class for modelling fluid properties for gas, oil and brine saturated sands.
References: Wang(2001), Batzle and Wang(1992), Geophysics.
- class pyavo.avotools.gassmann.GassmannSub(vp: Union[int, float], vs: Union[int, float], rho: Union[int, float], rho_o: Union[int, float], rho_g: Union[int, float], vsh: float, phi: float, swi: float, swt: float, S: float, T: Union[int, float], P: Union[int, float], init_fluid: str, final_fluid: str, GOR: float)¶
Bases:
object
Class to model Gassmann fluid substitution for brine sands, oil sands, and gas sands. It generates the P and S wave velocities and density after fluid substitution according to input parameters.
- Arguments:
vp = P-wave velocity from log (ft/s)
vs = S-wave velocity from log (ft/s)
rho = Bulk density form log (g/cc)
rho_o = Oil gravity (deg API)
rho_g = Gas gravity (API)
vsh = Shale volume from log
phi = Porosity
swi = Initial water saturation from log
swt = Target water saturation
S = Salinity (ppm)
T = Temperature (deg)
P = Pressure (psi)
init_fluid = Fluid type of initial hydrocarbon (gas or oil)
final_fluid = Fluid type of desired output where (gas or oil)
GOR = Gas-Oil ratio
- final_hc() → tuple¶
Computes the bulk modulus and density of the desired hydrocarbon (oil or gas).
- Returns
k_hyc: bulk modulus, rho_hyc: density
- init_hyc() → tuple¶
Computes Bulk modulus and density of initial hydrocarbon.
- Returns
k_hyc: Bulk modulus, rho_hyc: Density
- insitu_moduli(rho_fluid: float, rho_matrix: float, d_phi=True) → tuple¶
Computes the initial original moduli for saturated insitu rock. Density of the insitu saturated rock is calculated from the porosity log using the mass balance equation.
- Parameters
rho_fluid – density of insitu pore fluid phase (g/cc)
rho_matrix – density of rock matrix (g/cc)
d_phi – density derived from input RHOB log.
- Returns
k_sat_init: Bulk modulus, mu_sat_init: Shear modulus
- k_frame(k_mat: float, k_fld: float, k_sat: float) → float¶
Computes the bulk modulus of porous rock frame.
- Parameters
k_mat – Bulk modulus of rock matrix
k_fld – Bulk modulus of pore fluid phase
k_sat – Bulk modulus of saturated insitu rock
- Returns
k_frame: Bulk modulus
- k_rho_brine() → tuple¶
Computes the bulk modulus(Gpa) and density(g/cc) of brine.
- Returns
Bulk Modulus, Density
- k_rho_fluid() → tuple¶
Computes the Bulk modulus and density of the mixed pore fluid phase (Initial insitu model)
- Returns
k_fld: Bulk modulus of pore fluid, rho_fld: Density of pore fluid
- k_rho_sat(k_mat: float, rho_mat: float, k_frame: float) → tuple¶
Computes the Gassmann bulk modulus and density of saturated rock after fluid substitution. The desired fluid is defined by the target water saturation and type of hydrocarbon.
- Parameters
k_mat – Bulk modulus of rock matrix (GPa)
rho_mat – Density of rock matrix (g/cc)
k_frame – Bulk modulus of rock frame (Gpa)
- Returns
k_sat_new: Bulk modulus, rho_sat_new: Density
- pyavo.avotools.gassmann.k_rho_matrix(v_cly: float, k_cly: float, k_qtz: float, rho_cly: float, rho_qtz: float) → tuple¶
Calculate the Bulk modulus and Density of rock matrix.
- Parameters
v_cly – Volume of clay assumed to be 70% of Shale volume
k_cly – Bulk modulus of clay (Gpa)
k_qtz – Bulk modulus of quartz (Gpa)
rho_cly – Density of clay (g/cc)
rho_qtz – Density of quartz (g/cc)
- Returns
k_mat : Bulk modulus of rock matrix, rho_mat : Density of rock matrix
- pyavo.avotools.gassmann.vel_sat(k_sat: float, rho_sat: float, mu: float) → tuple¶
Estimate the seismic velocities after Gassmann fluid substitution using density and elastic moduli of saturated rock.
- Returns
vp_new : P-wave velocity, vs_new : S-wave velocity
pyavo.avotools.impedance module¶
Functions to calculate acoustic and elastic impedance from well logs
- pyavo.avotools.impedance.ai(vp: float, rho: float)¶
Calculates the acoustic impedance of a layer given velocities and densities.
- Parameters
vp – P-wave velocity (m/s)
rho – Density (g/cc)
- Returns
z: Acoustic Impedance
- pyavo.avotools.impedance.ei(vp: float, vs: float, rho: float, ang: int)¶
Computes the elastic impedance of a layer given velocities, densities and incidence angle.
- Parameters
vp – P-wave velocity (m/s)
vs – S-wave velocity (m/s)
rho – Density (g/cc)
ang – Angle of incidence (deg)
- Returns
ei: Elastic impedance
Reference: Connolly, P., 1999, Elastic impedance: The Leading Edge, 18, 438–452.
- pyavo.avotools.impedance.lame(vp: float, vs: float, rho: float)¶
Computes Lamé parameters - lambda_rho and mu_rho.
- Parameters
vp – P-wave velocity (m/s)
vs – S-wave velocity (m/s)
rho – Density (g/cc)
- Returns
lambda_rho: Lamé first parameter, mu_rho: Lamé second parameter
Reference: Goodway, B., T. Chen, and J. Downton, 1997, Improved AVO fluid detection and lithology discrimination using Lamé petrophysical parameters; “λρ”, “μρ”, & “λ/μ fluid stack” from P and S inversions: 67th Annual International Meeting, SEG, Expanded Abstracts, 183-186.
- pyavo.avotools.impedance.norm_ei(vp: float, vs: float, rho: float, vp_sh: float, vs_sh: float, rho_sh: float, ang: int)¶
Computes the normalized elastic impedance.
- Parameters
vp – P-wave velocity. (m/s)
vs – S-wave velocity. (m/s)
rho – Density. (g/cc)
vp_sh – P-wave velocity in shale (m/s)
vs_sh – S-wave velocity in shale (m/s)
rho_sh – Shale reference density.
ang – Angle of incidence
- Returns
nei: Normalized elastic impedance.
- Reference:
Whitcombe, D, 2002, Elastic impedance normalization, Geophysics, 67 (1), 60–62.
pyavo.avotools.log_crossplot module¶
Functions to create well log plots & crossplots of P-wave reflectivity vs angles, and gradient vs intercept.
- pyavo.avotools.log_crossplot.avo_plot(angle: numpy.ndarray, shuey: float, intercept: float, gradient: float, lim=1.5)¶
Generates crossplots of P-wave reflectivity vs angles, and gradient vs intercept.
- Parameters
angle – Angle of incidence
shuey – Shuey 2-term approximation of reflectivity
gradient – Gradient from shuey approximation
intercept – Intercept from shuey approximation
- pyavo.avotools.log_crossplot.crossplot(vp, vs, vpvs, rho, phi, GR, AI, NEI, lambda_rho, mu_rho)¶
Generates cross plots of elastic impedance, porosity, density and velocities.
- Parameters
vp – P-wave velocity from Vp log (m/s)
vs – S-wave velocity from Vs log (m/s)
vpvs – VpVS ratio from VP-VS log (m/s)
rho – Density from RHOB log (g/cc)
phi – Porosity from PHI log
GR – Gamma ray from GR log
AI – Acoustic impedance
NEI – Normalized acoustic impedance
lambda_rho – Lame’s first parameter
mu_rho – Lame’s second parameter
- pyavo.avotools.log_crossplot.plot_imp(vpvs: Union[numpy.ndarray, pandas.core.series.Series, float], vp: Union[numpy.ndarray, pandas.core.series.Series, float], vs: Union[numpy.ndarray, pandas.core.series.Series, float], rho: Union[numpy.ndarray, pandas.core.series.Series, float], angle: Union[float, int], h_well: Union[numpy.ndarray, pandas.core.series.Series], h_ref: Union[int, float])¶
Creates log plots of Poisson ratio, Acoustic impedance-Normalized Elastic Impedance(AI-NEI) and Lame’s Parameters (lambda-rho & mu-rho).
- Parameters
vpvs – Vp/Vs values form VP/VS Log (m/s)
vp – P-wave velocities form Vp log (m/s)
vs – S-wave velocities from Vs Log (m/s)
rho – Density form RHOB log (g/cc)
angle – Angle of incidence (deg)
h_well – Depth of well from Depth log (m)
h_ref – closest location of the top on the depth array as reference to normal elastic impedance
- Returns
Poisson ratio, Acoustic impedance, lambda_rho, mu_rho, NEI
pyavo.avotools.moduli module¶
Functions to convert between various acoustic/elastic parameters, and provides a way to calculate all the elastic moduli from Vp, Vs, and rho.
Using equations http://www.subsurfwiki.org/wiki/Elastic_modulus from Mavko, G, T Mukerji and J Dvorkin (2003), The Rock Physics Handbook, Cambridge University Press
Created by Matt Hall, 2014. Modified by Tola Abiodun, 2021.
- pyavo.avotools.moduli.bulk(vp=None, vs=None, rho=None, mu=None, lam=None, youngs=None, pr=None, pmod=None)¶
Calculate the Bulk modulus given either P-wave Velocity, S-wave velocity, and Density, or any two elastic moduli (e.g. lambda and mu, or bulk and P moduli). SI units only.
- Args:
vp, vs, and rho or any 2 from lam, mu, youngs, pr, and pmod
- Returns
Bulk modulus in pascals, Pa
- pyavo.avotools.moduli.comp_vel(youngs=None, vs=None, rho=None, mu=None, lam=None, bulk=None, pr=None, pmod=None)¶
Calculates the P-wave velocity given bulk density and any two elastic moduli (e.g. lambda and mu, or Young’s and P moduli). SI units only.
- Args:
Any 2 from lam, mu, youngs, pr, pmod, bulk and Rho
- Returns
Vp in m/s
- pyavo.avotools.moduli.elastic_mod(vp, vs, rho) → dict¶
Computes elastic moduli given P-wave velocity, S-wave velocity, and Density. SI units only.
- Args:
Vp, Vs, and rho
- Returns
A dict of elastic moduli, plus P-wave impedance.
- pyavo.avotools.moduli.lame_param(vp=None, vs=None, rho=None, pr=None, mu=None, youngs=None, bulk=None, pmod=None)¶
Computes Lame’s first parameter given either P-wave velocity, S-wave velocity, and Density, or any two elastic moduli (e.g. bulk and mu, or Young’s and P moduli). SI units only.
- Args:
vp, vs, and rho or any 2 from bulk, mu, youngs, pr, and pmod
- Returns
Lambda in pascals, Pa
- pyavo.avotools.moduli.mu(vp=None, vs=None, rho=None, pr=None, lam=None, youngs=None, bulk=None, pmod=None)¶
Computes shear modulus given either Vp, Vs, and rho, or any two elastic moduli (e.g. lambda and bulk modulus, or Young’s and P moduli). SI units only.
- Args:
vp, vs, and rho or any 2 from lam, bulk, youngs, pr, and pmod
- Returns
Shear modulus in pascals, Pa
- pyavo.avotools.moduli.p_mod(vp=None, vs=None, rho=None, pr=None, mu=None, lam=None, youngs=None, bulk=None)¶
Computes P-wave modulus given either P-wave velocity, S-wave velocity, and Density, or any two elastic moduli (e.g. lambda and mu, or Young’s and bulk moduli). SI units only.
- Args:
vp, vs, and rho or any 2 from lam, mu, youngs, pr, and bulk
- Returns
P-wave modulus in pascals, Pa
- pyavo.avotools.moduli.poisson(vp=None, vs=None, rho=None, mu=None, lam=None, youngs=None, bulk=None, pmod=None)¶
Calculate the Poisson ratio given either P-wave Velocity, S-wave velocity, and Density, or any two elastic moduli (e.g. lambda and mu, or bulk and P moduli). SI units only.
- Args:
vp, vs, and rho or any 2 from lam, mu, youngs, bulk, and pmod
- Returns
Poisson’s ratio, dimensionless
- pyavo.avotools.moduli.shear_vel(rho=None, mu=None)¶
Computes the Shear wave velocity given bulk density and shear modulus. SI units only.
- Args:
Mu, Rho
- Returns
Vs in m/s
- pyavo.avotools.moduli.youngs(vp=None, vs=None, rho=None, mu=None, lam=None, bulk=None, pr=None, pmod=None)¶
Computes Young’s modulus given either P-wave Velocity, S-wave velocity, and Density, or any two elastic moduli (e.g. lambda and mu, or bulk and P moduli). SI units only.
- Args:
vp, vs, and rho or any 2 from lam, mu, bulk, pr, and pmod
- Returns
Young’s modulus in pascals, Pa