API

Input parameters

SeismicWaves.InputParametersAcousticType
struct InputParametersAcoustic{T, N} <: InputParameters{T, N}

Parameters for acoustic wave simulations.

  • ntimesteps::Int64: Number of time steps

  • dt::Any: Time step

  • gridsize::NTuple{N, Int64} where N: Grid size for each dimension

  • gridspacing::NTuple{N, T} where {T, N}: Grid spacing in each direction

  • boundcond::InputBoundaryConditionParameters: Kind of boundary conditions

source
SeismicWaves.InputParametersElasticType
struct InputParametersElastic{T, N} <: InputParameters{T, N}

Parameters for elastic wave simulations

  • ntimesteps::Int64: Number of time steps

  • dt::Any: Time step

  • gridsize::NTuple{N, Int64} where N: Grid size for each dimension

  • gridspacing::NTuple{N, T} where {T, N}: Grid spacing in each direction

  • boundcond::InputBoundaryConditionParameters: Kind of boundary conditions

source

Input BDCs parameters

SeismicWaves.CPMLBoundaryConditionParametersType
struct CPMLBoundaryConditionParameters{T} <: InputBoundaryConditionParameters{T}

CPML boundary conditions parameters for wave simulations.

  • halo::Int64: Number of CPML grid points

  • rcoef::Any: Target reflection coefficient

  • freeboundtop::Bool: Free surface boundary condition at the top

source

Material properties

SeismicWaves.ElasticIsoMaterialPropertiesType
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: Density

source
SeismicWaves.VpRhoAcousticVDMaterialPropertiesType
mutable struct VpRhoAcousticVDMaterialProperties{T, N} <: MaterialProperties{T, N}

Material properties for acoustic variable-density simulation.

  • vp::Array: P-wave velocity

  • rho::Array: Density

  • interp_method::SeismicWaves.InterpolationMethod: Interpolation method

source
SeismicWaves.VpRhoAcousticVDMaterialPropertiesMethod
VpRhoAcousticVDMaterialProperties(
  vp::Array{T, N},
  rho::Array{T, N};
  interp_method::InterpolationMethod=ArithmeticAverageInterpolation(2)
) where {T, N}

Constructor for material properties for acoustic variable-density simulation.

source

Shots, sources and receivers

SeismicWaves.MomentTensorShotType
struct 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.

source
SeismicWaves.MomentTensorSourcesMethod
MomentTensorSources( 
    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.

source
SeismicWaves.ScalarReceiversType
struct ScalarReceivers{T} <: Receivers{T}

Type representing a multi-receiver configuration for a wave propagation shot.

  • positions::Matrix: Receiver positions

  • seismograms::Matrix: Array holding seismograms (as columns)

  • observed::Matrix: Array holding observed seismograms (as columns)

  • invcov::AbstractMatrix: Inverse of the covariance matrix

  • windows::Vector{Pair{Int64, Int64}}: Windows of data used for misfit calculations

source
SeismicWaves.ScalarReceiversMethod
 ScalarReceivers(
   positions::Matrix{T},
   nt::Int;
   observed::Union{Matrix{T}, Nothing}=nothing,
   invcov::Union{AbstractMatrix{T}, Nothing}=nothing,
   windows::Union{Vector{Pair{Int,Int}}, Nothing}=nothing
 ) where {T}
ScalarReceivers(
    positions::Array{T, 2},
    nt::Int64;
    observed,
    invcov,
    windows
) -> ScalarReceivers

Create a single shot wave propagation receivers configuration from receivers positions.

source
SeismicWaves.ScalarShotType
struct 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.

source
SeismicWaves.ScalarSourcesType
struct ScalarSources{T} <: Sources{T}

Type representing a multi-source configuration for a wave propagation shot.

  • positions::Matrix: Source positions

  • tf::Matrix: Source time function

  • domfreq::Any: Dominant frequency

source
SeismicWaves.ScalarSourcesMethod
ScalarSources(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.

source
SeismicWaves.ShotType
abstract type Shot{T}

Shot is the abstract supertype describing a shot composed of sources and receivers.

Currently implemented concrete shot types are:

source
SeismicWaves.VectorReceiversType
struct VectorReceivers{T, N} <: Receivers{T}

Type representing a multi-receiver configuration for a wave propagation shot.

  • positions::Matrix: Receiver positions

  • seismograms::Array{T, 3} where T: Array holding seismograms (as columns)

  • observed::Array{T, 3} where T: Array holding observed seismograms (as columns)

  • invcov::AbstractMatrix: Inverse of the covariance matrix

source

Solver functions

SeismicWaves.swforward!Function
swforward!(
    params::InputParameters{T, N},
    matprop::MaterialProperties{T, N},
    shots::Array{<:Shot{T}, 1};
    parall,
    snapevery,
    infoevery,
    logger
) -> Any

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 are Shot structures. 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.jl GPU backend performing automatic domain decomposition if set to :CUDA
    • the AMDGPU.jl GPU backend performing automatic domain decomposition if set to :AMDGPU
    • Base.Threads CPU threads performing automatic domain decomposition if set to :threads
    • Base.Threads CPU threads sending a group of sources to each thread if set to :threadpersrc
    • otherwise the serial version if set to :serial
  • snapevery::Union{Int, Nothing} = nothing: if specified, saves itermediate snapshots at the specified frequency (one every snapevery time step iteration) and return them as a vector of arrays.
  • infoevery::Union{Int, Nothing} = nothing: if specified, logs info about the current state of simulation every infoevery time steps.
  • logger::Union{Nothing, AbstractLogger} = nothing: if specified, uses the given logger object to print logs, otherwise it uses the logger returned from current_logger().

See also InputParameters, MaterialProperties and Shot.

source
swforward!(
    wavesim::Union{WaveSimulation{T, N}, Array{<:WaveSimulation{T, N}, 1}},
    matprop::MaterialProperties{T, N},
    shots::Array{<:Shot{T}, 1};
    logger,
    kwargs...
) -> 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 are Shot structures. Each shot contains information about both source(s) and receiver(s).

Keyword arguments

  • logger::Union{Nothing, AbstractLogger} = nothing: if specified, uses the given logger object to print logs, otherwise it uses the logger returned from current_logger().

See also InputParameters, MaterialProperties and Shot.

source
SeismicWaves.swmisfit!Function
swmisfit!(
    params::InputParameters{T, N},
    matprop::MaterialProperties{T, N},
    shots::Array{<:Shot{T}, 1};
    parall,
    misfit,
    logger
) -> Any

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 are Shot structures. 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.jl GPU backend performing automatic domain decomposition if set to :CUDA
    • the AMDGPU.jl GPU backend performing automatic domain decomposition if set to :AMDGPU
    • Base.Threads CPU threads performing automatic domain decomposition if set to :threads
    • Base.Threads CPU threads sending a group of sources to each thread if set to :threadpersrc
    • otherwise the serial version if set to :serial
  • logger::Union{Nothing,AbstractLogger}: specifies the logger to be used.

See also InputParameters, MaterialProperties and Shot. See also swforward! and swgradient! and Shot.

source
swmisfit!(
    wavesim::Union{WaveSimulation{T, N}, Array{<:WaveSimulation{T, N}, 1}},
    matprop::MaterialProperties{T, N},
    shots::Array{<:Shot{T}, 1};
    logger,
    kwargs...
) -> 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}}}: input WaveSimulation 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 kind (e.g., Vp only is required for an acoustic constant-density simulation).
  • shots::Vector{<:Shot{T}}: a vector whose elements are Shot structures. 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.jl GPU backend performing automatic domain decomposition if set to :CUDA
    • the AMDGPU.jl GPU backend performing automatic domain decomposition if set to :AMDGPU
    • Base.Threads CPU threads performing automatic domain decomposition if set to :threads
    • Base.Threads CPU threads sending a group of sources to each thread if set to :threadpersrc
    • otherwise the serial version if set to :serial

See also InputParameters, MaterialProperties and Shot. See also swforward! and swgradient! and Shot.

source
SeismicWaves.swgradient!Function
swgradient!(
    params::InputParameters{T, N},
    matprop::MaterialProperties{T, N},
    shots::Array{<:Shot{T}, 1};
    parall,
    check_freq,
    infoevery,
    compute_misfit,
    misfit,
    smooth_radius,
    logger
) -> 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 are Shot structures. 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.jl GPU backend performing automatic domain decomposition if set to :CUDA
    • the AMDGPU.jl GPU backend performing automatic domain decomposition if set to :AMDGPU
    • Base.Threads CPU threads performing automatic domain decomposition if set to :threads
    • Base.Threads CPU threads sending a group of sources to each thread if set to :threadpersrc
    • otherwise the serial version if set to :serial
  • check_freq::Union{Int, Nothing} = nothing: if specified, enables checkpointing and specifies the checkpointing frequency.
  • infoevery::Union{Int, Nothing} = nothing: if specified, logs info about the current state of simulation every infoevery time steps.
  • compute_misfit::Bool = false: if true, also computes and return misfit value.
  • smooth_radius::Int = 5: 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 sources positions.
  • logger::Union{Nothing,AbstractLogger}: specifies the logger to be used.

See also InputParameters, MaterialProperties and Shot. See also swforward! and swmisfit! and Shot.

source
swgradient!(
    wavesim::Union{WaveSimulation{T, N}, Array{<:WaveSimulation{T, N}, 1}},
    matprop::MaterialProperties{T, N},
    shots::Array{<:Shot{T}, 1};
    logger,
    kwargs...
) -> 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}}}: input WaveSimulation 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 kind (e.g., Vp only is required for an acoustic constant-density simulation).
  • shots::Vector{<:Shot{T}}: a vector whose elements are Shot structures. 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.jl GPU backend performing automatic domain decomposition if set to :CUDA
    • the AMDGPU.jl GPU backend performing automatic domain decomposition if set to :AMDGPU
    • Base.Threads CPU threads performing automatic domain decomposition if set to :threads
    • Base.Threads CPU threads sending a group of sources to each thread if set to :threadpersrc
    • otherwise the serial version if set to :serial
  • check_freq::Union{Int, Nothing} = nothing: if specified, enables checkpointing and specifies the checkpointing frequency.
  • infoevery::Union{Int, Nothing} = nothing: if specified, logs info about the current state of simulation every infoevery time steps.
  • compute_misfit::Bool = false: if true, also computes and return misfit value.
  • smooth_radius::Int = 5: 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 sources positions.
  • logger::Union{Nothing,AbstractLogger}: specifies the logger to be used.

See also InputParameters, MaterialProperties and Shot. See also swforward! and swmisfit! and Shot.

source

Utilities

SeismicWaves.build_wavesimFunction
build_wavesim(
    params::InputParameters{T, N},
    matprop::MaterialProperties{T, N};
    parall,
    kwargs...
)

Builds a wave similation 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).

Keyword arguments

  • parall::Symbol = :threads: controls which backend is used for computation:
    • the CUDA.jl GPU backend performing automatic domain decomposition if set to :CUDA
    • the AMDGPU.jl GPU backend performing automatic domain decomposition if set to :AMDGPU
    • Base.Threads CPU threads performing automatic domain decomposition if set to :threads
    • Base.Threads CPU threads sending a group of sources to each thread if set to :threadpersrc
    • otherwise the serial version if set to :serial
  • gradient::Bool = false: whether the wave simulation is used for gradients computations.
  • check_freq::Union{<:Int, Nothing} = nothing: if gradient = true and if specified, enables checkpointing and specifies the checkpointing frequency.
  • snapevery::Union{<:Int, Nothing} = nothing: if specified, saves itermediate snapshots at the specified frequency (one every snapevery time 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 every infoevery time steps.
source
SeismicWaves.gaussstfFunction
gaussstf(t::Real, t0::Real, f0::Real) -> Real

Gaussian source time function for current time t, activation time t0 and dominating frequency f0.

source
SeismicWaves.gaussderivstfFunction
gaussderivstf(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.

source
SeismicWaves.rickerstfFunction
rickerstf(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.

source