seis
Reflection Interface Methods
sim2x.seis_interface
Functions for interface reflectivity modelling.
This module contains functions for.
Zoeppritz full scattering matrix.
akirichards(thetai, vp1, vs1, rho1, vp2, vs2, rho2, method='avseth')
Aki-Richards forumlation of reflectivity functions.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
P-wave angle of incidence for wavefront in radians |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
method |
Optional |
Defaults to 'avseth' - avseth formulation 'ar' - original aki-richards |
'avseth' |
Returns:
Type | Description |
---|---|
(numpy.ndarray) |
Rp(theta) |
bortfeld(thetai, vp1, vs1, rho1, vp2, vs2, rho2)
Calculates the solution to the full Bortfeld equations (1961)
These are approximations which work well when the interval velocity is well defined.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
P-wave angle of incidence for wavefront in radians |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
Returns:
Type | Description |
---|---|
ndarray |
Rp(thetai) P-wave reflectivity of angle theta |
calcreflp(vp1, vs1, rho1, vp2, vs2, rho2)
Calculates the reflectivity parameters of an interface.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
Returns:
Type | Description |
---|---|
Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray] |
[rVp, rVs, rrho, rVsVp, dVp, dVs, drho] |
shuey(thetai, vp1, vs1, rho1, vp2, vs2, rho2, mode='rtheta')
Shuey approximation to the Aki-Richards equations.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
P-wave angle of incidence for wavefront in radians |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
mode |
Literal['rtheta', 'R0_G'] |
what to return 'rtheta' returns Rp(theta) 'R0_G' returns [R0,G] aka [A,B] |
'rtheta' |
Returns:
Type | Description |
---|---|
ndarray |
if |
snellrr(thetai, vp1, vs1, vp2, vs2)
Snell's reflection and refraction angles for a two layered half space.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
downgoing p-wave angle of incidence for wavefront (radians) |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
Returns:
Type | Description |
---|---|
Tuple[float, numpy.ndarray, numpy.ndarray, numpy.ndarray] |
returns a list of calculated angles using Snell's Law thetai = input angle of P-wave incidence and output angle P-wave reflection thetat = output angle of P-wave transmission phir = output angle of S-wave reflection phit = output angle of S-wave transmission |
zoeppritz_pdpu_only(thetai, vp1, vs1, rho1, vp2, vs2, rho2)
Calculates the solution to an incident down-going P-wave ray and upgoing P-wave only.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
P-wave angle of incidence for wavefront in radians |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
Returns:
Type | Description |
---|---|
Rp |
zoeppritz_ponly(thetai, vp1, vs1, rho1, vp2, vs2, rho2)
Calculates the solution to an incident down-going P-wave ray.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
P-wave angle of incidence for wavefront in radians |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
Returns:
Type | Description |
---|---|
Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray] |
(Rp, Rs, Tp, Ts) of amplitudes for reflected and transmitted rays |
zoeppritzfull(thetai, vp1, vs1, rho1, vp2, vs2, rho2)
Full Zoeppritz scattering matrix solution
Parameters:
Name | Type | Description | Default |
---|---|---|---|
thetai |
float |
P-wave angle of incidence for wavefront in radians |
required |
vp1 |
ndarray |
p-velocity layer 1 |
required |
vs1 |
ndarray |
s-velocity layer 1 |
required |
rho1 |
ndarray |
density layer 1 |
required |
vp2 |
ndarray |
p-velocity layer 2 |
required |
vs2 |
ndarray |
s-veloicty layer 2 |
required |
rho2 |
ndarray |
density layer 2 |
required |
Returns:
Type | Description |
---|---|
ndarray |
[Rp, Rs, Tp, Ts] of amplitudes for reflected and transmitted rays |
Domain Conversion Methods
sim2x.twt_conv
calc_domain_range(domain_vol, srate=0.001)
Calculate the domain range so we convieniently land on whole numbers starting from near the limit of our data and incrementing at srate.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
vol |
Volume mapping dims of |
required | |
srate |
float |
The sample rate in the target domain |
0.001 |
Returns:
Type | Description |
---|---|
appropriate domain_stick for domain_vol to convert to |
domain_convert_vol(domain_vol, ds, mapping_dims=('xline', 'iline'), interpolation_kind='linear', to_dim='twt', from_dim='depth', domain_vol_direction='forward', output_samples=None)
Convert volumes between different domains (e.g. TWT to DEPTH)
Parameters:
Name | Type | Description | Default |
---|---|---|---|
domain_vol |
DataArray |
A DataArray with the |
required |
ds |
Dataset |
A Dataset with the properties to depth convert. Has same dims s |
required |
mapping_dims |
Sequence[str] |
The dimensions over which to map the funciton |
('xline', 'iline') |
interpolation_kind |
As for scipy.interpolate.interp1d(kind=) or |
'linear' |
|
to_dim |
str |
The domain dimension for the current input |
'twt' |
from_dim |
str |
The domain dimension to map to |
'depth' |
domain_vol_direction |
Literal['forward', 'reverse'] |
If domain volume maps |
'forward' |
output_samples |
Union[NoneType, Sequence[Sequence[Sequence[Sequence[Sequence[Any]]]]], numpy._array_like._SupportsArray[numpy.dtype], Sequence[numpy._array_like._SupportsArray[numpy.dtype]], Sequence[Sequence[numpy._array_like._SupportsArray[numpy.dtype]]], Sequence[Sequence[Sequence[numpy._array_like._SupportsArray[numpy.dtype]]]], Sequence[Sequence[Sequence[Sequence[numpy._array_like._SupportsArray[numpy.dtype]]]]], bool, int, float, complex, str, bytes, Sequence[Union[bool, int, float, complex, str, bytes]], Sequence[Sequence[Union[bool, int, float, complex, str, bytes]]], Sequence[Sequence[Sequence[Union[bool, int, float, complex, str, bytes]]]], Sequence[Sequence[Sequence[Sequence[Union[bool, int, float, complex, str, bytes]]]]]] |
Specify the samples at which ds will be interpolated to using |
None |
Returns:
Type | Description |
---|---|
Dataset |
The properties of |
peg_time_interval(twt_vol, depth_surf, twt_surf, mapping_dims=('xline', 'iline'))
Shifts twt_vol to create such that the depth conversion of depth_surf
with twt_vol
matches
twt_surf
Parameters:
Name | Type | Description | Default |
---|---|---|---|
twt_vol |
DataArray |
TWT volume in depth |
required |
depth_surf |
DataArray |
The depth surface at |
required |
twt_surf |
DataArray |
The TWT surface at |
required |
time_integral(vp, depth_dim='depth')
Get the time integral from a vp Dataset with axis depth
Synthetic Seismic
sim2x.seis_synth
Functions necessary for synthetic seismic modelling (e.g. sim2seis).
convolution1d_vol(ds, reflectivity_key, wavelet_amp, mapping_dims=('xline', 'iline'))
Convolved a reflectivity array from a Dataset with a wavelet.
Wavelet must have same sample rate as seismic twt
Parameters:
Name | Type | Description | Default |
---|---|---|---|
ds |
Dataset |
A Dataset with the properties to depth convert. Has same dims s |
required |
mapping_dims |
Sequence[str] |
The dimensions over which to map the funciton |
('xline', 'iline') |
Returns:
Type | Description |
---|---|
Dataset |
The properties of |
reflectivity(theta, velp, vels, rho, method='full')
Reflectivity for theta from acoustic vectors
Parameters:
Name | Type | Description | Default |
---|---|---|---|
theta |
float |
P-wave incidence angle |
required |
velp |
ndarray |
P-wave velcoities |
required |
vels |
ndarray |
S-wave velocities |
required |
rho |
ndarray |
density values |
required |
method |
Literal['full', 'zoep_pud', 'ar'] |
The modelling method to use |
'full' |
Returns:
Type | Description |
---|---|
ndarray |
interface reflectivity arrays (padded 0 at front to keep same length as input) |
reflectivity_vol(ds, theta, mapping_dims=('xline', 'iline'))
Convert elastic properties "vp", "vs" and "density" to reflectivity
Parameters:
Name | Type | Description | Default |
---|---|---|---|
ds |
Dataset |
A Dataset with the properties to depth convert. Has same dims s |
required |
mapping_dims |
Sequence[str] |
The dimensions over which to map the funciton |
('xline', 'iline') |
Returns:
Type | Description |
---|---|
Dataset |
The properties of |