fdModel class

The fdModel class

class fdModel

Finite difference wave modelling class.

This class contains everything needed to do finite difference wave forward and adjoint modelling. It contains the entire experimental parameters as fields, which are loaded at runtime from the supplied .ini file. The class contains all necessary functions to perform FWI, but lacks optimization schemes.

Public Functions

explicit fdModel(const char *configuration_file_relative_path)

Constructor for modelling class.

This constructor creates the modelling class from a configuration file supplied at runtime. As such, all fields are created dynamically.

Parameters

configuration_file_relative_path – Relative path to the configuration .ini file. This file should contain all the fields needed for simulation. Arbitrary defaults are hardcoded into the binary as backup within the parse_configuration() method.

fdModel(const int nt, const int nx_inner, const int nz_inner, const int nx_inner_boundary, const int nz_inner_boundary, const real_simulation dx, const real_simulation dz, const real_simulation dt, const int np_boundary, const real_simulation np_factor, const real_simulation scalar_rho, const real_simulation scalar_vp, const real_simulation scalar_vs, const int npx, const int npz, const real_simulation peak_frequency, const real_simulation source_timeshift, const real_simulation delay_cycles_per_shot, const int n_sources, const int n_shots, const std::vector<int> ix_sources_vector, const std::vector<int> iz_sources_vector, const std::vector<real_simulation> moment_angles_vector, const std::vector<std::vector<int>> which_source_to_fire_in_which_shot, const int nr, const std::vector<int> ix_receivers_vector, const std::vector<int> iz_receivers_vector, const int snapshot_interval, const std::string observed_data_folder, const std::string stf_folder)
fdModel(const fdModel &model)
~fdModel()

Destructor for the class.

The destructor properly addresses every used new keyword in the constructor, freeing all memory.

void allocate_memory()
void initialize_arrays()
void copy_arrays(const fdModel &model)
void parse_parameters(const std::vector<int> ix_sources_vector, const std::vector<int> iz_sources_vector, const std::vector<real_simulation> moment_angles_vector, const std::vector<int> ix_receivers_vector, const std::vector<int> iz_receivers_vector)
void parse_configuration_file(const char *configuration_file_relative_path)

Method that parses .ini configuration file. Only used in fdModel().

Parameters

configuration_file_relative_path – Relative path to the configuration .ini file.

void forward_simulate(int i_shot, bool store_fields, bool verbose, bool output_wavefields = false)

Method to forward simulate wavefields for a specific shot.

Forward simulate wavefields of shot i_shot based on currently loaded models. Storage of wavefields and verbosity of run can be toggled.

Parameters
  • i_shot – Integer controlling which shot to simulate.

  • store_fields – Boolean to control storage of wavefields. If storage is not required (i.e. no adjoint modeling), forward simulation should be faster without storage.

  • verbose – Boolean controlling if modelling should be verbose.

void adjoint_simulate(int i_shot, bool verbose)

Method to adjoint simulate wavefields for a specific shot.

Adjoint simulate wavefields of shot i_shot based on currently loaded models and calculate adjoint sources. Verbosity of run can be toggled.

Parameters
  • i_shot – Integer controlling which shot to simulate.

  • verbose – Boolean controlling if modelling should be verbose.

void write_receivers()

Method to write out synthetic seismograms to plaintext.

This method writes out the synthetic seismograms to a plaintext file. Allows one to e.g. subsequently import these files as observed seismograms later using load_receivers(). Every shot generates a separate ux and uz receiver file (rtf_ux/rtf_uz), with every receiver being a single line in these files.

void write_receivers(std::string prefix)
void write_sources()

Method to write out source signals to plaintext.

This method writes out the source time function (without moment tensor) to plaintext file. Useful for e.g. visualizing the source staggering.

void load_receivers(bool verbose)

Method to load receiver files.

This method loads receiver data from observed_data_folder folder into the object. The data has to exactly match the set-up, and be named according to component and shot (as generated by write_receivers() ).

Parameters

verbose – Controls the verbosity of the method during loading.

void map_kernels_to_velocity()

Method to map kernels to velocity parameter set.

This method takes the kernels (lambda, mu, rho) as originally calculated on the Lamé’s parameter set and maps them to the velocity parameter set (vp, vs, rho).

void update_from_velocity()

Method to map velocities into Lamé’s parameter set.

This method updates Lamé’s parameters of the current model based on the velocity parameters of the current model. Typically has to be done every time after updating velocity.

void calculate_l2_misfit()

Method to calculate L2 misfit.

This method calculates L2 misfit between observed seismograms and synthetic seismograms and stores it in the misfit field.

void calculate_l2_adjoint_sources()

Method to calculate L2 adjoint sources This method calculates L2 misfit between observed seismograms and synthetic seismograms and stores it in the misfit field.

void load_model(const std::string &de_path, const std::string &vp_path, const std::string &vs_path, bool verbose)

Method to load models from plaintext into the model.

This methods loads any appropriate model (expressed in density, P-wave velocity, and S-wave velocity) into the class and updates the Lamé fields accordingly.

Parameters
  • de_path – Relative path to plaintext de file.

  • vp_path – Relative path to plaintext vp file.

  • vs_path – Relative path to plaintext vs file.

  • verbose – Boolean controlling the verbosity of the method.

void run_model(bool verbose, bool simulate_adjoint)

Method to perform all steps necessary for FWI, with additional control over adjoint simulation.

This method performs all necessary steps in FWI; forward modelling, misfit calculation and optionally adjoint source calculation, adjoint modelling and kernel projection.

Parameters
  • verbose – Boolean controlling the verbosity of the method.

  • simulate_adjoint – Boolean controlling the execution of the adjoint simulation and kernel computation.

void reset_kernels()

Method to reset all Lamé sensitivity kernels to zero.

This method resets all sensitivity kernels to zero. Essential before performing new adjoint simulations, as otherwise the kernels of subsequent simulations would stack.

void write_kernels()
dynamic_vector get_model_vector()
void set_model_vector(dynamic_vector m)
dynamic_vector get_gradient_vector()
dynamic_vector load_vector(const std::string &vector_path, bool verbose)

Public Members

real_simulation c1 = real_simulation(9.0 / 8.0)
real_simulation c2 = real_simulation(1.0 / 24.0)
bool add_np_to_source_location = true
bool add_np_to_receiver_location = true
real_simulation **vx

Dynamic horizontal velocity field used in the simulations.

real_simulation **vz

Dynamic vertical velocity field used in the simulations.

real_simulation **txx

Dynamic horizontal stress field used in the simulations.

real_simulation **tzz

Dynamic vertical stress field used in the simulations.

real_simulation **txz

Dynamic shear stress field used in the simulations.

real_simulation **lm
real_simulation **la
real_simulation **mu
real_simulation **b_vx
real_simulation **b_vz
real_simulation **rho
real_simulation **vp
real_simulation **vs
real_simulation **lambda_kernel
real_simulation **mu_kernel
real_simulation **density_l_kernel
real_simulation **vp_kernel
real_simulation **vs_kernel
real_simulation **density_v_kernel
real_simulation **starting_rho
real_simulation **starting_vp
real_simulation **starting_vs
real_simulation **taper
real_simulation *t
real_simulation **stf
real_simulation ***moment
real_simulation ***rtf_ux
real_simulation ***rtf_uz
real_simulation ***rtf_ux_true
real_simulation ***rtf_uz_true
real_simulation ***a_stf_ux
real_simulation ***a_stf_uz
real_simulation ****accu_vx
real_simulation ****accu_vz
real_simulation ****accu_txx
real_simulation ****accu_tzz
real_simulation ****accu_txz
int nt
int nx_inner
int nz_inner
int nx_inner_boundary
int nz_inner_boundary
real_simulation dx
real_simulation dz
real_simulation dt
int np_boundary
real_simulation np_factor
real_simulation scalar_rho
real_simulation scalar_vp
real_simulation scalar_vs
int n_sources
int n_shots
std::vector<std::vector<int>> which_source_to_fire_in_which_shot
real_simulation delay_cycles_per_shot
int *ix_sources
int *iz_sources
real_simulation *moment_angles
real_simulation peak_frequency
real_simulation alpha
real_simulation t0
int nr
int *ix_receivers
int *iz_receivers
int snapshot_interval
int snapshots
int nx
int nz
int nx_free_parameters
int nz_free_parameters
int basis_gridpoints_x = 1
int basis_gridpoints_z = 1
int free_parameters
real_simulation misfit
std::string observed_data_folder
std::string stf_folder