API
Input parameters
SeismicWaves.InputParameters — Typeabstract type InputParameters{T, N}InputParameters is the abstract supertype describing input parameters for wave simulations.
Currently implemented concrete parameters are InputParametersAcoustic and InputParametersElastic, for acoustic and elasic wave simulation respectively.
SeismicWaves.InputParametersAcoustic — Typestruct InputParametersAcoustic{T, N} <: InputParameters{T, N}Parameters for acoustic wave simulations.
ntimesteps::Int64: Number of time stepsdt::Any: Time stepgridsize::NTuple{N, Int64} where N: Grid size for each dimensiongridspacing::NTuple{N, T} where {T, N}: Grid spacing in each directionboundcond::InputBoundaryConditionParameters: Kind of boundary conditions
SeismicWaves.InputParametersElastic — Typestruct InputParametersElastic{T, N} <: InputParameters{T, N}Parameters for elastic wave simulations
ntimesteps::Int64: Number of time stepsdt::Any: Time stepgridsize::NTuple{N, Int64} where N: Grid size for each dimensiongridspacing::NTuple{N, T} where {T, N}: Grid spacing in each directionboundcond::InputBoundaryConditionParameters: Kind of boundary conditions
Input BDCs parameters
SeismicWaves.CPMLBoundaryConditionParameters — Typestruct CPMLBoundaryConditionParameters{T} <: InputBoundaryConditionParameters{T}CPML boundary conditions parameters for wave simulations.
halo::Int64: Number of CPML grid pointsrcoef::Any: Target reflection coefficientfreeboundtop::Bool: Free surface boundary condition at the topvel_max::Union{Nothing, T} where T: Maximum velocity for C-PML coefficients
SeismicWaves.InputBoundaryConditionParameters — Typeabstract type InputBoundaryConditionParameters{T}InputBoundaryConditionParameters is the abstract supertype describing boundary conditions input parameters for wave simulations.
Currently implemented concrete parameters are CPMLBoundaryConditionParameters and ReflectiveBoundaryConditionParameters.
SeismicWaves.ReflectiveBoundaryConditionParameters — Typestruct ReflectiveBoundaryConditionParameters{T, N} <: InputBoundaryConditionParameters{T}Reflective boundary conditions parameters for wave simulations.
Material properties
SeismicWaves.ElasticIsoMaterialProperties — Typemutable struct ElasticIsoMaterialProperties{T, N} <: SeismicWaves.AbstrElasticIsoMaterialProperties{T, N}Material properties for elastic isotropic simulation.
λ::Array: First Lamé parameterμ::Array: Second Lamé parameter (shear modulus)ρ::Array: Densityinterp_method_ρ::SeismicWaves.AbstractInterpolationMethod: Interpolation method for densityinterp_method_λ::SeismicWaves.AbstractInterpolationMethod: Interpolation method for densityinterp_method_μ::SeismicWaves.AbstractInterpolationMethod: Interpolation method for density
SeismicWaves.MaterialProperties — Typeabstract type MaterialProperties{T, N}MaterialProperties is the abstract supertype describing material properties for wave simulations. It defines which type of wave equation is solved.
Currently implemented concrete properties are:
VpAcousticCDMaterialPropertiesfor acoustic constant density wave equationVpRhoAcousticVDMaterialPropertiesfor acoustic variable density wave equationElasticIsoMaterialPropertiesfor elasitic isotropic wave equation
SeismicWaves.VpAcousticCDMaterialProperties — Typestruct VpAcousticCDMaterialProperties{T, N} <: MaterialProperties{T, N}Material properties for acoustic constant-density simulation.
vp::Array: P-wave velocity
SeismicWaves.VpRhoAcousticVDMaterialProperties — Typemutable struct VpRhoAcousticVDMaterialProperties{T, N} <: MaterialProperties{T, N}Material properties for acoustic variable-density simulation.
vp::Array: P-wave velocityrho::Array: Densityinterp_method::SeismicWaves.AbstractInterpolationMethod: Interpolation method
SeismicWaves.VpRhoAcousticVDMaterialProperties — MethodVpRhoAcousticVDMaterialProperties(
vp::Array{T, N},
rho::Array{T, N};
interp_method::InterpolationMethod=ArithmeticAverageInterpolation()
) where {T, N}Constructor for material properties for acoustic variable-density simulation.
Shots, sources and receivers
SeismicWaves.ExternalForceShot — Typestruct ExternalForceShot{T, N} <: Shot{T}Type representing a shot with external force sources and multi-component receivers.
srcs::ExternalForceSources: Structure containing the ExternalForcesSources for a given simulation.recs::VectorReceivers: Structure containing the VectorReceivers for a given simulation.
SeismicWaves.MomentTensorShot — Typestruct MomentTensorShot{T, N, M<:SeismicWaves.MomentTensor{T, N}} <: Shot{T}Type representing a shot with moment tensor sources and multi-component receivers.
srcs::MomentTensorSources: Structure containing the MomentTensorSources for a given simulation.recs::VectorReceivers: Structure containing the VectorReceivers for a given simulation.
SeismicWaves.MomentTensorSources — TypeType representing a multi-source configuration for a wave propagation shot.
SeismicWaves.MomentTensorSources — MethodMomentTensorSources(
positions::Matrix{T},
tf::Matrix{T},
momtens::Vector{M}
domfreq::T
) where {T, N, M <: MomentTensor{T}}Create a single shot wave propagation source configuration from source positions, time-functions and a dominant frequency.
SeismicWaves.Receivers — Typeabstract type Receivers{T}Sources is the abstract supertype describing seismic receivers.
Currently implemented concrete receivers types are:
ScalarReceiversfor scalar field receiversVectorReceiversfor vector field receivers
SeismicWaves.ScalarReceivers — Typestruct ScalarReceivers{T} <: Receivers{T}Type representing a multi-receiver configuration for a wave propagation shot.
positions::Matrix: Receiver positionsseismograms::Matrix: Array holding seismograms (as columns)
SeismicWaves.ScalarReceivers — MethodScalarReceivers(
positions::Array{T, 2},
nt::Int64
) -> ScalarReceivers
Create a single shot wave propagation receivers configuration from receiver positions.
SeismicWaves.ScalarShot — Typestruct ScalarShot{T} <: Shot{T}Type representing a shot with scalar sources and receivers.
srcs::ScalarSources: Structure containing the ScalarSources for a given simulation.recs::ScalarReceivers: Structure containing the ScalarReceivers for a given simulation.
SeismicWaves.ScalarSources — Typestruct ScalarSources{T} <: Sources{T}Type representing a multi-source configuration for a wave propagation shot.
positions::Matrix: Source positionstf::Matrix: Source time functiondomfreq::Any: Dominant frequency
SeismicWaves.ScalarSources — MethodScalarSources(positions::Matrix{T}, tf::Matrix{T}, domfreq::T) where {T}Create a single shot wave propagation source configuration from source positions, time-functions and a dominant frequency.
SeismicWaves.Shot — Typeabstract type Shot{T}Shot is the abstract supertype describing a shot composed of sources and receivers.
Currently implemented concrete shot types are:
ScalarShotfor scalar field sources and receiversMomentTensorShotfor moment tensor sources and vector field receivers
SeismicWaves.Sources — Typeabstract type Sources{T}Sources is the abstract supertype describing seismic sources.
Currently implemented concrete sources types are:
ScalarSourcesfor scalar field sourcesMomentTensorSourcesfor moment tensor sources
SeismicWaves.VectorReceivers — Typestruct VectorReceivers{T, N} <: Receivers{T}Type representing a multi-receiver configuration for a wave propagation shot.
positions::Matrix: Receiver positionsseismograms::Array{T, 3} where T: Array holding seismograms (as columns)
SeismicWaves.VectorReceivers — MethodVectorReceivers(
positions::Array{T, 2},
nt::Int64
) -> VectorReceivers
VectorReceivers(
positions::Array{T, 2},
nt::Int64,
ndim::Int64
) -> VectorReceivers
Create a single shot wave propagation receivers configuration from receiver positions.
Solver functions
SeismicWaves.swforward! — Functionswforward!(
params::InputParameters{T, N},
matprop::MaterialProperties{T, N},
shots::Array{<:Shot{T}, 1};
runparams
)
Compute forward simulation using the given input parameters params and material properties matprop on multiple shots. Receivers traces are stored in the receivers for each shot.
Return a vector of Dict containing for each shot the snapshots of the fields computed in the simulation for each timestep.
Positional arguments
params::InputParameters{T, N}: input parameters for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation type.matprop::MaterialProperties{T, N}: material properties for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation type.shots::Vector{<:Shot{T}}: a vector whose elements areShotstructures. Each shot contains information about both source(s) and receiver(s).
Keyword arguments
runparams::RunParameters: a struct containing parameters related to forward calculations. SeeRunParametersfor details. In case of a forward simulation,gradparamsis set tonothing.
See also InputParameters, MaterialProperties and Shot.
swforward!(
wavesim::Union{WaveSimulation{T, N}, Array{<:WaveSimulation{T, N}, 1}},
matprop::MaterialProperties{T, N},
shots::Array{<:Shot{T}, 1}
) -> Any
Compute forward simulation using a previously constructed WaveSimulation object. See also build_wavesim on how to build the WaveSimulation. Receivers traces are stored in the receivers for each shot.
Return a vector of Dict containing for each shot the snapshots of the fields computed in the simulation for each timestep.
Positional arguments
wavesim: wave simulation object containing all required information to run the simulation.matprop::MaterialProperties{T, N}: material properties for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation type.shots::Vector{<:Shot{T}}: a vector whose elements areShotstructures. Each shot contains information about both source(s) and receiver(s).
See also InputParameters, MaterialProperties and Shot.
SeismicWaves.swmisfit! — Functionswmisfit!(
params::InputParameters{T, N},
matprop::MaterialProperties{T, N},
shots::Array{<:Shot{T}, 1},
misfit::Array{<:AbstractMisfit{T}, 1};
runparams
)
Return the misfit w.r.t. observed data by running a forward simulation using the given input parameters params and material properties matprop on multiple shots. Receivers traces are stored in the receivers for each shot.
Positional arguments
params::InputParameters{T, N}: input parameters for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., acoustic variable-density).matprop::MaterialProperties{T, N}: material properties for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., Vp only is required for an acoustic constant-density simulation).shots::Vector{<:Shot{T}}: a vector whose elements areShotstructures. Each shot contains information about both source(s) and receiver(s).
Keyword arguments
parall::Symbol = :threads: controls which backend is used for computation:- the
CUDA.jlGPU backend performing automatic domain decomposition if set to:CUDA - the
AMDGPU.jlGPU backend performing automatic domain decomposition if set to:AMDGPU - the
Metal.jlGPU backend performing automatic domain decomposition if set to:Metal Base.ThreadsCPU threads performing automatic domain decomposition if set to:threadsBase.ThreadsCPU threads sending a group of sources to each thread if set to:threadpersrc- otherwise the serial version if set to
:serial
- the
logger::Union{Nothing,AbstractLogger}: specifies the logger to be used.
See also InputParameters, MaterialProperties and Shot. See also swforward! and swgradient! and Shot.
swmisfit!(
wavesim::Union{WaveSimulation{T, N}, Array{<:WaveSimulation{T, N}, 1}},
matprop::MaterialProperties{T, N},
shots::Array{<:Shot{T}, 1},
misfit::Array{<:AbstractMisfit{T}, 1}
) -> Any
Return the misfit w.r.t. observed data by running a forward simulation using the given WaveSimulation object as an input. Receivers traces are stored in the receivers for each shot. See also build_wavesim on how to build the WaveSimulation.
Positional arguments
wavesim::Union{WaveSimulation{T,N},Vector{<:WaveSimulation{T,N}}}: inputWaveSimulationobject containing all required information to run the simulation.matprop::MaterialProperties{T, N}: material properties for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., Vp only is required for an acoustic constant-density simulation).shots::Vector{<:Shot{T}}: a vector whose elements areShotstructures. Each shot contains information about both source(s) and receiver(s).
Keyword arguments
parall::Symbol = :threads: controls which backend is used for computation:- the
CUDA.jlGPU backend performing automatic domain decomposition if set to:CUDA - the
AMDGPU.jlGPU backend performing automatic domain decomposition if set to:AMDGPU - the
Metal.jlGPU backend performing automatic domain decomposition if set to:Metal Base.ThreadsCPU threads performing automatic domain decomposition if set to:threadsBase.ThreadsCPU threads sending a group of sources to each thread if set to:threadpersrc- otherwise the serial version if set to
:serial
- the
See also InputParameters, MaterialProperties and Shot. See also swforward! and swgradient! and Shot.
SeismicWaves.swgradient! — Functionswgradient!(
params::InputParameters{T, N},
matprop::MaterialProperties{T, N},
shots::Array{<:Shot{T}, 1},
misfit::Array{<:AbstractMisfit{T}, 1};
runparams,
gradparams
) -> Any
Compute gradients w.r.t. model parameters using the given input parameters params and material parameters matprop on multiple shots.
The check_freq parameter controls the checkpoiting frequency for adjoint computation. If nothing, no checkpointing is performed. If greater than 2, a checkpoint is saved every check_freq time step. The optimal tradeoff value is check_freq = sqrt(nt) where nt is the number of time steps of the forward simulation. Bigger values speed up computation at the cost of using more memory.
Positional arguments
params::InputParameters{T, N}: input parameters for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., acoustic variable-density).matprop::MaterialProperties{T, N}: material properties for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., Vp only is required for an acoustic constant-density simulation).shots::Vector{<:Shot{T}}: a vector whose elements areShotstructures. Each shot contains information about both source(s) and receiver(s).
Keyword arguments
runparams::RunParameters: a struct containing parameters related to forward calculations. SeeRunParametersfor details. In case of a forward simulation,gradparamsis set tonothing.gradparams::Union{GradParameters,Nothing}: a struct containing parameters related to gradient calculations. SeeGradParametersfor details.
Return
The gradient as a dictionary, one key per each material property. If compute_misfit in gradparams is set to true, additionally returns also the misfit value.
See also InputParameters, MaterialProperties and Shot. See also swforward! and swmisfit! and Shot.
swgradient!(
wavesim::Union{WaveSimulation{T, N}, Array{<:WaveSimulation{T, N}, 1}},
matprop::MaterialProperties{T, N},
shots::Array{<:Shot{T}, 1},
misfit::Array{<:AbstractMisfit{T}, 1}
) -> Any
Compute gradients w.r.t. model parameters using the previously built WaveSimulation. This avoids re-initializing and re-allocating several arrays in case of multiple gradient calculations. See also build_wavesim on how to build the WaveSimulation.
The check_freq parameter controls the checkpoiting frequency for adjoint computation. If nothing, no checkpointing is performed. If greater than 2, a checkpoint is saved every check_freq time step. The optimal tradeoff value is check_freq = sqrt(nt) where nt is the number of time steps of the forward simulation. Bigger values speed up computation at the cost of using more memory.
Positional arguments
wavesim::Union{WaveSimulation{T,N},Vector{<:WaveSimulation{T,N}}}: inputWaveSimulationobject containing all required information to run the simulation.matprop::MaterialProperties{T, N}: material properties for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., Vp only is required for an acoustic constant-density simulation).shots::Vector{<:Shot{T}}: a vector whose elements areShotstructures. Each shot contains information about both source(s) and receiver(s).
Keyword arguments
runparams::RunParameters: a struct containing parameters related to forward calculations. SeeRunParametersfor details. In case of a forward simulation,gradparamsis set tonothing.gradparams::Union{GradParameters,Nothing}: a struct containing parameters related to gradient calculations. SeeGradParametersfor details.
Return
The gradient as a dictionary, one key per each material property. If compute_misfit in gradparams is set to true, additionally returns also the misfit value.
See also InputParameters, MaterialProperties and Shot. See also swforward! and swmisfit! and Shot.
SeismicWaves.RunParameters — Typestruct RunParametersRunParameters is a struct containing various parameters related to forward (and gradient) simulations.
Keyword arguments
parall::Symbol = :threads: controls which backend is used for computation:- the
CUDA.jlGPU backend performing automatic domain decomposition if set to:CUDA - the
AMDGPU.jlGPU backend performing automatic domain decomposition if set to:AMDGPU - the
Metal.jlGPU backend performing automatic domain decomposition if set to:Metal Base.ThreadsCPU threads performing automatic domain decomposition if set to:threadsBase.ThreadsCPU threads sending a group of sources to each thread if set to:threadpersrc- otherwise the serial version if set to
:serial
- the
snapevery::Union{<:Int, Nothing} = nothing: if specified, saves itermediate snapshots at the specified frequency (one everysnapeverytime step iteration) and return them as a vector of arrays (only for forward simulations).infoevery::Union{<:Int, Nothing} = nothing: if specified, logs info about the current state of simulation everyinfoeverytime steps.logger::AbstractLogger=current_logger(): specifies the logger to be used.erroronCFL::Bool=true: throw an error if the CFL condition is not met.minPPW::Int=10: minimum number of points per wavelength (PPW).erroronPPW::Bool=true: throw an error if the minimum number of points per wavelength (PPW) is not achieved (otherwise warn about it).
SeismicWaves.GradParameters — Typestruct GradParametersGradParameters is a struct containing various parameters related specifically to gradient simulations.
Keyword arguments
mute_radius_src::Int=0: grid points inside a ball with radius specified by the parameter (in grid points) will have their gradient smoothed by a factor inversely proportional to their distance from source positions.mute_radius_rec::Int=0: grid points inside a ball with radius specified by the parameter (in grid points) will have their gradient smoothed by a factor inversely proportional to their distance from receiver positions.compute_misfit::Bool=false: default false. If true, also computes and returns the misfit value.check_freq::Int=1: if specified, enables checkpointing and defines the checkpointing frequency. Deafults to no checkpointing (=1).
Utilities
SeismicWaves.build_wavesim — Functionbuild_wavesim(
params::InputParameters{T, N},
matprop::MaterialProperties{T, N};
runparams,
gradparams,
kwargs...
)
Builds a wave simulation object based on the input paramters params and keyword arguments kwargs.
Positional arguments
params::InputParameters{T,N}: input parameters for the simulation, where T represents the data type and N represents the number of dimensions. They vary depending on the simulation kind (e.g., acoustic variable-density).matprop::MaterialProperties{T, N}: material properties.
Keyword arguments
gradient::Bool = false: whether the wave simulation is used for gradients computations or not.runparams::RunParameters: a struct containing parameters related to forward calculations. SeeRunParametersfor details.gradparams::Union{GradParameters,Nothing}: a struct containing parameters related to gradient calculations. SeeGradParametersfor details. In case of a forward simulation,gradparamsis set tonothing.
SeismicWaves.gaussstf — Functiongaussstf(t::Real, t0::Real, f0::Real) -> Real
Gaussian source time function for current time t, activation time t0 and dominating frequency f0.
SeismicWaves.gaussderivstf — Functiongaussderivstf(t::Real, t0::Real, f0::Real) -> Real
First derivative of gaussian source time function for current time t, activation time t0 and dominating frequency f0.
SeismicWaves.rickerstf — Functionrickerstf(t::Real, t0::Real, f0::Real) -> Real
Ricker source (second derivative of gaussian) source time function for current time t, activation time t0 and dominating frequency f0.