pyavo.seismodel package¶
Submodules¶
pyavo.seismodel.angle_stack module¶
Functions to import segy data and compute AVO attributes from Near and Far Stacks.
- pyavo.seismodel.angle_stack.avo_attributes(horizon_data: str, near_stack: xarray.core.dataarray.DataArray, far_stack: xarray.core.dataarray.DataArray, theta_near: int, theta_far: int, TWT_slice: tuple, XL_slice: tuple, inline: int, well_XL=None, well_display=True, robust=True, interpolation='spline16', add_colorbar=True, cmap='jet_r') → dict¶
Computes intercept, gradient values and display images of AVO attributes & crossplot.
- Parameters
horizon_data – file path, any valid string path to a txt file is acceptable.
near_stack – DataArray of Near angle stack
far_stack – DataArray of Far angle stack
theta_near – Incidence angle of near stack
theta_far – Incidence angle of far stack
TWT_slice – Two-way time slice
XL_slice – Crossline time slice
inline – Inline number
well_XL – Well crossline number to define well position
- Returns
c: Intercept, m: Gradient
- pyavo.seismodel.angle_stack.nf_attributes(horizon_data: str, near_stack: xarray.core.dataarray.DataArray, far_stack: xarray.core.dataarray.DataArray, TWT_slice: tuple, XL_slice: tuple, inline: int, well_XL=None, well_display=True, robust=True, interpolation='spline16', add_colorbar=True, cmap='jet_r')¶
Display images of near and far angle stack attributes. Allowed file extension for horizon_data is .txt. Both near and far stack data are read as DataArray to enhance computation speed.
- Parameters
horizon_data – file path, any valid string path to a txt file is acceptable.
near_stack – DataArray of Near angle stack
far_stack – DataArray of Far angle stack
TWT_slice – Two-way time slice
XL_slice – Crossline time slice
inline – Inline number
well_XL – Well crossline number to define well position
well_display – Show well trajectory on plot, vertical
- Reference:
Avseth, P., Mukerji, T., & Mavko, G., 2010. Quantitative seismic interpretation: Applying rock physics tools to reduce interpretation risk.Cambridge university press.
- pyavo.seismodel.angle_stack.nfstack(horizon_data: str, near_stack: xarray.core.dataarray.DataArray, far_stack: xarray.core.dataarray.DataArray, inline: int, robust=True, interpolation='spline16', cmap='RdBu', well_XL=None, well_display=True)¶
Display images of Near and Far angle stacks from DataArray. Allowed file extension for horizon_data is .txt. Both near and far stack data are read as DataArray to enhance computation speed.
- Parameters
horizon_data – file path, any valid string path to a txt file is acceptable.
near_stack – DataArray of Near angle stack
far_stack – DataArray of Far angle stack
inline – Inline number
well_XL – Well crossline number to define well position
well_display – Show well trajectory on plot, vertical
- pyavo.seismodel.angle_stack.read_segy(filepath: str, byte_il=14, byte_xl=21, ignore_geometry=False) → xarray.core.dataarray.DataArray¶
Read a segyfile into a Data array.
- Parameters
filepath – filepath, any valid string path is acceptable
byte_il – inline byte
byte_xl – crossline byte
ignore_geometry – Opt out on building geometry information, useful for e.g. shot organised files. Defaults to False.
- Returns
data_array: xarray
pyavo.seismodel.rpm module¶
Functions to create rock physics models and apply them to well log analysis to quantify seismic responses.
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
- pyavo.seismodel.rpm.cp_model(k_min: float, mu_min: float, phi: float, cp=0.4) → tuple¶
Build rock physics models using critical porosity concept.
- Parameters
k_min – Bulk modulus of mineral (Gpa)
mu_min – Shear modulus of mineral (Gpa)
phi – Porosity (fraction)
cp – critical porosity, default is 40%
- Returns
Bulk and Shear modulus of rock frame
- Reference
Critical porosity, Nur et al. (1991, 1995) Mavko, G, T. Mukerji and J. Dvorkin (2009), The Rock Physics Handbook: Cambridge University Press. p.353
- pyavo.seismodel.rpm.hz_mindlin(k_min, mu_min, phi: numpy.ndarray, cd_num=8.6, cp=0.4, P=10, f=1) → tuple¶
Computes the bulk modulus and shear modulus of the Hertz-Mindlin rock physics model.
- Parameters
k_min – bulk modulus of mineral (Gpa)
mu_min – shear modulus of mineral (Gpa)
phi – porosity (fraction) expressed as vector values
cd_num – coordination number
cp – critical porosity (random close pack of well-sorted rounded quartz grains)
P – Confining pressure/rock stress (Mpa)
f – shear modulus correction factor 0 = dry frictionless packing 1 = dry pack with perfect adhesion
- Returns
Bulk modulus and Shear modulus
- Reference
The Rock Physics Handbook: Tools for Seismic Analysis of Porous Media - Gary M. Mavko, Jack Dvorkin, and Tapan Mukerji
- pyavo.seismodel.rpm.plot_log(data: pandas.core.frame.DataFrame, vsh: pandas.core.series.Series, vp: pandas.core.series.Series, vs: pandas.core.series.Series, imp: pandas.core.series.Series, vp_vs: pandas.core.series.Series, phi: pandas.core.series.Series, rho: pandas.core.series.Series, z_init: float, z_final: float, shale_cutoff: float, sand_cutoff: float)¶
Explore key input logs for Rock physics modeling. Generates a log plot of the shale volume, acoustic impedance, Vp/Vs and a cross plot of Porosity against P-wave velocity and acoustic impedance against Vp/VS.
- Parameters
data – pandas DataFrame
vsh – Shale volume as input from vshale log
vp – Compressional wave velocity log
vs – Shear wave velocity log
imp – Acoustic impedance from acoustic impedance log
vp_vs – Vp/Vs ratios
phi – porosity values form PHI log
rho – density values from RHOB log
z_init – minimum depth from depth log
z_final – maximum depth from depth log
shale_cutoff – Shale cutoff
sand_cutoff – Sand cutoff
- Returns
log plots(Vp/Vs, AI, Vsh) and cross plots(AI,Vp/Vs & PHI,Vp)
- pyavo.seismodel.rpm.plot_rpm(df: pandas.core.frame.DataFrame, k_qtz: float, mu_qtz: float, k_sh: float, mu_sh: float, rho_sh: float, rho_qtz: float, k_br: float, rho_br: float, phi: numpy.ndarray, NG: numpy.ndarray, vsh: pandas.core.series.Series, z_init: float, z_final: float, sand_cut: float, shale_cut: float, P: Union[int, float], cp_ssm: float, cn_ssm: float, cp_stm: float, cn_stm: float, eff_phi: pandas.core.series.Series, vp: pandas.core.series.Series)¶
Plot soft sand and stiff sand rock physics models for different mineralogies.
- Parameters
df – Pandas DataFrame
k_qtz – Bulk modulus of quartz (GPa)
mu_qtz – Shear modulus of quartz (GPa)
k_sh – Bulk modulus of shale (Gpa)
mu_sh – Shear modulus of shale (Gpa)
rho_sh – Density of shale (g/cc)
rho_qtz – Density of quartz (g/cc)
k_br – Bulk modulus of brine (GPa)
rho_br – Density of brine (g/cc)
phi – Porosity (fraction)
NG – Net-to-Gross
vsh – Volume of shale from VSH log
z_init – minimum depth withing range of input depth log (ft)
z_final – maximum depth withing range of input depth log (ft)
sand_cut – sand cutoff
shale_cut – shale cutoff
P – Confining pressure (Mpa)
cp_ssm – Critical porosity for soft sand model
cn_ssm – Coordination number for soft sand model
cp_stm – Critical porosity for stiff sand model
cn_stm – Coordination number for stiff sand model
eff_phi – effective porosity values from porosity log (fraction)
vp – P-wave velocity from VP log (m/s)
- Returns
Plot of soft sand and stiff sand models with curves showing varying net-to-gross.
- pyavo.seismodel.rpm.plot_rpt(model='soft', fluid='oil', display=True, vsh=0.0, k_sh=15, k_qtz=37, mu_sh=5, mu_qtz=44, rho_sh=2.8, rho_qtz=2.6, k_gas=0.06, k_oil=0.9, k_br=2.8, rho_gas=0.2, rho_oil=0.8, rho_br=1.1, cp=0.4, cd_num=8.6, P=10, f=1)¶
Plot RPTs showing variations in porosity, mineralogy and fluid content defined by water saturation.
- Parameters
model – type of rock physics model, ‘soft’ or ‘stiff’. default = ‘soft’
fluid – desired fluid after Gassmann’s fluid substitution, default = ‘oil’
display – show plot with calculated parameters, default = True
vsh – volume of shale, default = 0
k_sh – bulk modulus of shale, default = 15GPa
k_qtz – bulk modulus of quartz, default = 37GPa
mu_sh – shear modulus of shale, default = 5GPa
mu_qtz – shear modulus of quartz, default = 44GPa
rho_sh – density of shale, default = 2.8g/cc
rho_qtz – density of quartz, default = 2.6g/cc
k_gas – bulk modulus of gas (GPa), default = 0.06GPa
k_oil – bulk modulus of oil (GPa), default = 0.9GPa
k_br – bulk modulus of brine (Gpa), default = 2.8GPa
rho_gas – density of gas (g/cc), default = 0.2g/cc
rho_oil – density of oil (g/cc), default = 0.8g/cc
rho_br – density of brine (g/cc) default = 1.1g/cc
cp – critical porosity, default = 0.4
cd_num – coordination number, default = 8.6
P – confining pressure(MPa), default = 10MPa
f – shear modulus correction factor, default=1
- Returns
plot of RPT with labels, Array of Acoustic impedance and Vp/Vs
- Reference:
Quantitative Seismic Interpretation: Applying Rock Physics Tools to Reduce Interpretation Risk (Per Avseth, T. Mukerji and G.mavko, 2015)
- pyavo.seismodel.rpm.soft_sand(k_min: float, mu_min: float, phi: numpy.ndarray, cd_num=8.6, cp=0.4, P=10, f=1)¶
Computes the bulk modulus and shear modulus of soft-sand (uncemented) rock physics model.
- Parameters
k_min – bulk modulus of mineral (Gpa)
mu_min – shear modulus of mineral (Gpa)
phi – porosity (fraction) expressed as vectors
cd_num – coordination number
cp – critical porosity (random close pack of well-sorted rounded quartz grains)
P – Confining pressure/rock stress (Gpa)
f – shear modulus correction factor 0 = dry frictionless packing 1 = dry pack with perfect adhesion
- Returns
Bulk modulus and Shear modulus
- Reference:
The Rock Physics Handbook: Tools for Seismic Analysis of Porous Media - Gary M. Mavko, Jack Dvorkin, and Tapan Mukerji
- pyavo.seismodel.rpm.stiff_sand(k_min: float, mu_min: float, phi: numpy.ndarray, cd_num=8.6, cp=0.4, P=10, f=1)¶
Computes the bulk modulus and shear modulus of stiff-sand(uncemented) rock physics model.
- Parameters
k_min – bulk modulus of mineral (Gpa)
mu_min – shear modulus of mineral (Gpa)
phi – porosity (fraction) expressed as vectors
cd_num – coordination number
cp – critical porosity (random close pack of well-sorted rounded quartz grains)
P – Confining pressure/rock stress (Gpa)
f – shear modulus correction factor 0 = dry frictionless packing 1 = dry pack with perfect adhesion
- Returns
Bulk modulus and Shear modulus
- Reference:
The Rock Physics Handbook: Tools for Seismic Analysis of Porous Media - Gary M. Mavko, Jack Dvorkin, and Tapan Mukerji
- pyavo.seismodel.rpm.twolayer(n_samples: int, vp1: float, vs1: float, rho1: float, vp2: float, vs2: float, rho2: float)¶
Display seismic signatures at Near and Far offset traces and AVO curve for a two layer model.
- Parameters
n_samples – Number of samples
vp1 – P-wave velocity in top layer
vs1 – S-Wave velocity in top layer
rho1 – Density of top layer
vp2 – P-wave velocity in bottom layer
vs2 – S-wave velocity in bottom layer
rho2 – Density of second layer
- Reference
Alessandro aadm (2015)
- pyavo.seismodel.rpm.vel_sat(k_min: float, rho_min: float, k_fl: float, rho_fl: float, k_frame: float, mu_frame: float, phi: numpy.ndarray) → tuple¶
Computes velocities and densities of saturated rock using Gassmann’s fluid substitution equations.
- Parameters
k_min – bulk modulus of mineral (Gpa)
rho_min – density of mineral (g/cc)
k_fl – bulk modulus of fluid (Gpa)
rho_fl – density of fluid (g/cc)
k_frame – bulk modulus of rock frame/dry rock (Gpa)
mu_frame – shear modulus of rock frame/dry rock (Gpa)
phi – porosity expressed as vector values
- Returns
k_sat: Bulk modulus rho_sat: Density vp_sat: P-wave velocity vs_sat: S-wave velocity
- pyavo.seismodel.rpm.voigt_reuss(mod1: float, mod2: float, vfrac: float) → tuple¶
Computes the Voigt-Reuss-Hill bounds and average for a 2-component mixture. The Voigt Bound is obtained by assuming the strain field remains constant throughout the material when subjected to an arbitrary average stress field. The Reuss Bound is obtained by assuming the stress field remains constant throughout the material in an arbitrary average strain field.
- Parameters
mod1 – elastic modulus of first mineral
mod2 – elastic modulus of second mineral
vfrac – volumetric fraction of first mineral (net-to-gross)
- Returns
upper bound, lower bound, Voigt-Ruess-Hill average
- Reference
Mavko, G, T. Mukerji and J. Dvorkin (2009), The Rock Physics Handbook: Cambridge University Press.
- pyavo.seismodel.rpm.well_rpt(df, z_init: float, z_final: float, sand_cut: float, vsh: pandas.core.series.Series, vp: pandas.core.series.Series, vs: pandas.core.series.Series, rho: pandas.core.series.Series, ip_rpt0: numpy.ndarray, vpvs_rpt0: numpy.ndarray, ip_rpt1: numpy.ndarray, vpvs_rpt1: numpy.ndarray)¶
Plot Soft and Stiff Sand RPT models and superimpose on well data.
- Parameters
df – Pandas Dataframe
z_init – minimum depth values withing depth log range
z_final – maximum depth values within depth log range
sand_cut – cut off sand
vsh – volume of shale from VSH log
vp – P-wave velocity values from VP log
vs – S-wave velocity values from VS log
rho – density values from RHOB log
ip_rpt0 – Acoustic impedance values of first RPT model
vpvs_rpt0 – Vp/Vs values of first RPT model
ip_rpt1 – Acoustic impedance values of second RPT model
vpvs_rpt1 – Vp/Vs values of second RPT model
pyavo.seismodel.tuning_prestack module¶
Functions to generate a synthetic angle gather from a 3-layer property model to examine pre-stack tuning effects.
- pyavo.seismodel.tuning_prestack.avo_inv(rc_zoep, ntrc: int, top: list) → dict¶
AVO inversion for INTERCEPT and GRADIENT from analytic and convolved reflectivity values. Linear least squares method is used for estimating INTERCEPT and GRADIENT coefficients.
- Parameters
rc_zoep – Zoeppritz reflection coefficient values
ntrc – Number of traces
top – Convolved top reflectivity values
- Returns
Zoeppritz and Convolved reflectivities
- pyavo.seismodel.tuning_prestack.calc_theta_rc(theta1_min: float, theta1_step: float, vp: list, vs: list, rho: list, ang) → tuple¶
Computes the reflection coefficients for given incidence angles in a three layer model.
- Parameters
theta1_min – minimum incidence angle
theta1_step – maximum incidence angle
vp – Layer P-wave velocities
vs – Layer S-wave velocities
rho – Density of layers
ang – Angle of incidence
- Returns
theta1_samp : n-samples of incidence angles, rc_1: Reflection coefficient at layer1 / layer2, rc_2: Reflection coefficient at layer2 / layer3.
- pyavo.seismodel.tuning_prestack.int_depth(h_int: list, thickness: float)¶
Computes depth to the interface for a three layer model.
- Parameters
h_int – depth to to first interface
thickness – thickness
- Returns
d_interface: interface depths
- pyavo.seismodel.tuning_prestack.layer_index(lyr_times: list, dt=0.0001) → tuple¶
Calculate array indices corresponding to top/base interfaces.
- Parameters
lyr_times – Interface times
dt – change in time, default is 0.0001. Changing this from 0.00005 can affect the display quality.
- Returns
- lyr1_indx: array
Layer 1
- lyr2_indx: array
Layer 2
- pyavo.seismodel.tuning_prestack.n_angles(theta1_min: float, theta1_max: float, theta1_step=1) → int¶
Computes number of traces for given angles on incidence.
- Parameters
theta1_min – minimum incidence angle
theta1_max – maximum incidence angle
theta1_step – angle of incidence steps, default is 1
- Returns
n_trace: number of traces
- pyavo.seismodel.tuning_prestack.plot_vawig(axhdl, data, t, excursion)¶
Format the display of synthetic angle gather.
- Parameters
axhdl – Plot axes
data – Synthetic traces generated from convolved zoeppritz.
t – Regularly spaced ime samples
excursion – Adjust plot width
- pyavo.seismodel.tuning_prestack.ray_param(v_int: float, theta: float) → float¶
Calculates the ray parameter P.
- Parameters
v_int – Interval velocity
theta – Angle of incidence (deg)
- Returns
p: Ray parameter
- pyavo.seismodel.tuning_prestack.rc_zoep(vp1: float, vs1: float, vp2: float, vs2: float, rho1: float, rho2: float, theta1: float) → numpy.ndarray¶
Calculate the Reflection & Transmission coefficients using Zoeppritz equations.
- Parameters
vp1 – P-wave velocity in upper layer
vs1 – S-wave velocity in upper layer
vp2 – P-wave velocity in lower layer
vs2 – S-wave velocity in lower layer
rho1 – Density of upper layer
rho2 – Density of lower layer
theta1 – Angle of incidence of ray
- Returns
R: Reflection coefficients
- Reference:
The Rock Physics Handbook, Dvorkin et al.
- pyavo.seismodel.tuning_prestack.syn_angle_gather(min_time: float, max_time: float, lyr_times: numpy.ndarray, thickness: float, top_ref: list, base_ref: list, vp_dig: numpy.ndarray, vs_dig: numpy.ndarray, rho_dig: numpy.ndarray, syn_zoep: numpy.ndarray, rc_zoep: numpy.ndarray, t: numpy.ndarray, excursion: int)¶
Plot synthetic angle gather for three layer model displayed in normal polarity and amplitudes extracted along the upper and lower interfaces.
- Parameters
min_time – Minimum plot time (s)
max_time – Maximum plot time (s)
lyr_times – interface times
thickness – apparent thickness of second layer
top_ref – convolved top reflectivity values
base_ref – convolved base reflectivity values
vp_dig – P-wave velocity in digital time domain
vs_dig – S-wave velocity in digital time domain
rho_dig – Density of layers in digital time domain
syn_zoep – Synthetic seismogram
rc_zoep – Zoeppritz reflectivities
t – regularly spaced time samples
excursion – Adjust plot width
- pyavo.seismodel.tuning_prestack.t_domain(t, vp: list, vs: list, rho: list, lyr1_index: list, lyr2_index: list) → tuple¶
Create a “digital” time domain version of an input property model.
- Parameters
t – array of regularly spaced time samples
vp – P-wave velocity
vs – S-wave velocity
rho – Density of layers
lyr1_index – Array indices of top interface
lyr2_index – Array indices of base interface
- Returns
vp_dig: P-wave velocity in digital t-domain, vs_dig: S-wave velocity in digital t-domain, rho_dig: Density in digital t-domain.
pyavo.seismodel.tuning_wedge module¶
Functions to generate a three layer wedge model using the reflection coefficients and travel times
Reference:
Mavko, G., T. Mukerji, and J. Dvorkin, 2009, The rock physics handbook: Tools for seismic analysis of porous media, 2nd ed.: Cambridge University Press. CrossRef
Shuey, R. T., 1985, A simplification of the Zoeppritz equations: Geophysics, 50, no. 4, 609–614, http://dx.doi.org/10.1190/1.1441936.
Widess, M., 1973, How thin is a thin bed?: Geophysics, 38, no. 6, 1176–1180, http://dx.doi.org/10.1190/1.1440403.
- pyavo.seismodel.tuning_wedge.calc_times(z_int: list, vp: list) → list¶
Calculate the travel time to a reflector.
- Parameters
z_int – Depth to a reflector
vp – P-Wave Velocity
- Returns
t_int: Interface times
- Usage
t_int = calc_times(z_int, vp)
- pyavo.seismodel.tuning_wedge.get_rc(Vp: list, rho: list) → list¶
Calculates the Reflection Coefficient at the interface separating two layers. The ratio of amplitude of the reflected wave to the incident wave, or how much energy is reflected. Typical values of R are approximately −1 from water to air, meaning that nearly 100% of the energy is reflected and none is transmitted; ~0.5 from water to rock; and ~0.2 for shale to sand.
- Parameters
Vp – P-wave velocity (m/s)
rho – Layer density (g/cc)
- Returns
rc_int: Reflection Coefficient
- pyavo.seismodel.tuning_wedge.int_depth(h_int: list, dh_min: float, dh_step: float, mod: int) → list¶
Computes the depth to an interface.
- Parameters
h_int – depth to first interface
dh_min – minimum thickness of layer 2
dh_step – Thickness step from trace-to-trace (usually 1.0m)
mod – model traces generated within a thickness or interval
- Returns
d_inteface: depth to interface
- pyavo.seismodel.tuning_wedge.mod_digitize(rc_int: list, t_int: list, t: list) → list¶
Digitize a 3 layer model using reflection coefficients and interface times.
- Parameters
rc_int – reflection coefficients corresponding to interface times
t_int – interface times
t – regularly sampled time series defining model sampling
- Returns
rc: Reflection Coefficient
- Usage
rc = digitize_model(rc_int, t_int, t)
- pyavo.seismodel.tuning_wedge.n_model(h_min: float, h_max: float, h_step=1) → int¶
Computes number of traces within an interval or thickness.
- Parameters
h_min – minimum thickness
h_max – maximum thickness
h_step – thickness steps, default is 1
- Returns
n_trace: number of traces
- pyavo.seismodel.tuning_wedge.ricker(sample_rate: float, length: float, c_freq: float)¶
Generate time and amplitude values for a zero-phase wavelet. The second derivative of the Gaussian function or the third derivative of the normal-probability density function.
- Parameters
sample_rate – sample rate in seconds
length – length of time (dt) in seconds
c_freq – central frequency of wavelet (cycles/seconds or Hz).
- Returns
wv_time: zero-phase wavelet time, wv_amp: zero-phase wavelet amplitude.
- Usage:
time, wavelet = (sample_rate,duration,c_freq)
- Reference:
A Ricker wavelet is often used as a zero-phase embedded wavelet in modeling and synthetic seismogram manufacture. Norman H. Ricker (1896–1980), American geophysicist.
- pyavo.seismodel.tuning_wedge.syn_seis(ref_coef: list, wav_amp: numpy.ndarray) → list¶
Generate synthetic seismogram from convolved reflectivities and wavelet.
- Parameters
ref_coef – Reflection coefficient
wav_amp – wavelet amplitude
- Returns
- smg: array
Synthetic seismogram
- pyavo.seismodel.tuning_wedge.time_samples(t_min: float, t_max: float, dt=0.0001) → list¶
Create regularly sampled time series defining model sampling.
- Parameters
t_min – Minimum time duration
t_max – Maximum time duration
dt – Change in time, default = 0.0001
- Returns
time
- pyavo.seismodel.tuning_wedge.wedge_model(syn_zo: numpy.ndarray, layer_times: numpy.ndarray, t: numpy.ndarray, t_min: float, t_max: float, h_max: float, h_step: float, excursion: int, dt=0.0001)¶
Plot a three layer wedge model amplitudes and zero offset seismogram.
- Parameters
syn_zo – Synthetic Seismogram
layer_times – travel time to reflectors
t – regularly sampled time series defining model sampling
t_min – minimum time duration
t_max – maximum time duration
h_max – maximum depth
h_step – depth steps, default is 1
excursion – adjust plot width
dt – trace parameter, changing this from 0.0001 can affect the display quality.
pyavo.seismodel.wavelet module¶
Functions to generate a ricker and trapezoidal bandpass wavelet.
- pyavo.seismodel.wavelet.bandpass(f1: float, f2: float, f3: float, f4: float, phase: float, dt: float, wvlt_length: float) → tuple¶
Calculate a trapezoidal bandpass wavelet.
f1: Low truncation frequency of wavelet in Hz f2: Low cut frequency of wavelet in Hz f3: High cut frequency of wavelet in Hz f4: High truncation frequency of wavelet in Hz phase: wavelet phase in degrees dt: sample rate in seconds wvlt_length: length of wavelet in seconds
- Usage:
t, wvlt = wvlt_ricker(f1, f2, f3, f4, phase, dt, wvlt_length)
- Reference
Wes Hamlyn, 2014.
- pyavo.seismodel.wavelet.plot_ricker(sample_rate=0.001, length=0.128, c_freq=25)¶
Computes and display a zero-phase wavelet.
- Parameters
sample_rate – sample rate in seconds.
length – length of time (dt) in seconds.
c_freq – central frequency of wavelet (cycles/seconds or Hz)
- Usage:
plot_ricker(sample_rate,duration,c_freq)
- pyavo.seismodel.wavelet.ricker(sample_rate=0.001, length=0.128, c_freq=25) → tuple¶
Generate time and amplitude values for a zero-phase wavelet. The second derivative of the Gaussian function or the third derivative of the normal-probability density function. A Ricker wavelet is often used as a zero-phase embedded wavelet in modeling and synthetic seismogram manufacture.
- Reference:
Norman H. Ricker (1896–1980), American geophysicist.
- Parameters
sample_rate – sample rate in seconds (float, int)
length – length of time (dt) in seconds (float, int)
c_freq – central frequency of wavelet (cycles/seconds or Hz). (float, int)
- Returns
ndarray
- Usage:
time, wavelet = (sample_rate,duration,c_freq)