seisflows.tools.model

A Model class for interfacing with SPECFEM velocity models, used to read models into memory, manipulate and write to disk. Also contains plotting functions for SPECFEM2D models

Classes

Model

A container for reading, storing and manipulating model/gradient/kernel

Module Contents

class seisflows.tools.model.Model(path=None, fmt='', parameters=None, regions='123', flavor=None)

A container for reading, storing and manipulating model/gradient/kernel parameters from SPECFEM2D/3D/3D_GLOBE.

acceptable_parameters = ['vp', 'vs', 'vpv', 'vph', 'vsv', 'vsh', 'eta', 'rho', 'rhop', 'kappa', 'mu', 'c11', 'c12',...
path = None
fmt = ''
flavor = None
model = None
coordinates = None
_parameters = None
_ngll = None
_nproc = None
setup(path=None)

Set up the internal Model definition by guessing or ingesting the file format, solver flavor, reading the model values and metadata and assigning internal variables

Parameters:

path (str) – optional path to Model files. If not given, defaults to the internal path attribute.

Raises:

AssertionError – if the path does not exist or the path does not contain expected Model files

fnfmt(i='*', val='*', ext='*')

Expected SPECFEM filename format with some checks to ensure that wildcards and numbers are accepted. An example filename is: ‘proc000001_vs.bin’

Parameters:
  • i (int or str) – processor number or wildcard. If given as an integer, will be converted to a 6 digit zero-leading value to match SPECFEM format

  • val (str) – parameter value (e.g., ‘vs’) or wildcard

  • ext (str) – the file format (e.g., ‘.bin’). If NOT preceded by a ‘.’ will have one prepended

Return type:

str

Returns:

filename formatter for use in model manipulation

property parameters

Returns a list of parameters which defines the model.

property ngll

Provide the number of GLL (Gauss Lobatto Legendre) points per processor chunk. Access hidden attribute _ngll in the case that the model is very large, we only want to count the GLL points once.

Return type:

list of float

Returns:

each float represents the number of GLL points for the chunk number that corresponds to its index

_update_ngll_from_model()

Convenience function to count NGLL points as length of data arrays for each parameter processor chunk

property nproc

Returns the number of processors that define the model/gradient/kernel.

Return type:

int

Returns:

number of processors that define the model

property vector

Conveience property to access the merge() function which creates a linear vector defining all model parameters

Return type:

np.array

Returns:

a linear vector of all model parameters

copy()

Returns a deep copy of self so that models can be transferred

read(parameters=None)

Utility function to load in SPECFEM models/kernels/gradients saved in various formats. Will try to guess format of model. Assumes that models are saved using the following filename format:

proc{num}_{val}.{format} where num is usually a 6 digit number representing the processor number (e.g., 000000), val is the parameter value of the model/kernel/gradient and format is the format of the file (e.g., bin)

Parameters:

parameters (list of str) – unique parameters to load model for, if None will load all available parameters found in path

Return type:

Dict of np arrays

Returns:

Dictionary where keys correspond to model parameters and values are vectors (np.arrays) representing the model

read_coordinates_specfem2d()

Attempt to read coordinate files from the given model definition. This is only really useful for SPECFEM2D, where we can plot the model, kernel and gradient using matplotlib.

Return type:

Dict

Returns:

a dictioanary with the X and Z coordinates read in from a SPECFEM2D model, if applicable

merge(parameter=None)

Convert dictionary representation model to vector representation m where all parameters and processors are stored as a single 1D vector. This vector representation is used by the optimization library during model perturbation.

Parameters:

parameter (str) – single parameter to retrieve model vector from, otherwise returns all parameters merged into single vector

Return type:

np.array

Returns:

vector representation of the model

write(path=None, fmt=None)

Save a SPECFEM model/gradient/kernel vector loaded into memory back to disk in the appropriate format expected by SPECFEM

split(vector=None)

Converts internal vector representation m to dictionary representation model. Does this by separating the vector based on how it was constructed, parameter-wise and processor-wise

Parameters:

vector (np.array) – allow Model to split an input vector. If none given, will split the internal vector representation

Return type:

Dict of np.array

Returns:

dictionary of model parameters split up by number of processors

print_stats()

Print out all model parameters (min, mean, max) in the logger. Useful for checking if models are different, if a model is being satisfactorily udpated, or for debugging purposes.

check(min_pr=-1.0, max_pr=0.5)

Checks parameters in the model. If Vs and Vp present, checks poissons ratio. Checks for negative velocity values.

Para min_pr:

minimum allowable Poisson’s ratio, if applicable

Parameters:

max_pr (float) – maximum allowable Poisson’s ratio, if applicable

Raises:

AssertionError

  • if the input model has no values for any of its parameters

  • if the model contains any NaN values

_check_2d3d_parameters(min_pr=-1.0, max_pr=0.5)

Checks parameters for SPECFEM2D and SPECFEM3D derived models

_check_3dglobe_parameters(min_pr=-1.0, max_pr=0.5)

Checks parameters for SPECFEM3D_GLOBE derived models

update(model=None, vector=None)

Update internal model/vector defitions. Because these two quantities are tied to one another, updating one will update the other. This function simply makes that easier.

plot2d(parameter, vmin=None, vmax=None, cmap=None, levels=None, center_cmap=False, title='', show=True, save=None)

Plot internal model parameters as a 2D image plot.

Warning

This is only available for SPECFEM2D models. SPECFEM3D model coordinates do not match the model vectors (because the grids are irregular) and cannot be visualized like this.

Parameters:
  • parameter (str) – chosen internal parameter value to plot.

  • cmap (str) – colormap which match available matplotlib colormap. If None, will choose default colormap based on parameter choice.

  • show (bool) – show the figure after plotting

  • title (str) – optional title to prepend to some values that are provided by default. Useful for adding information about iteration, step count etc.

  • save (str) – if not None, full path to figure to save the output image

_get_nproc_parameters()

Get the number of processors and the available parameters from a list of output SPECFEM model files.

Return type:

tuple (int, list)

Returns:

(number of processors, list of available parameters in dir)

_guess_file_format()

Guess the file format of model/kernel/gradient files if none provided by the user. Does so by checking file formats against formats expected from SPECFEM2D/3D/3D_GLOBE

Return type:

str

Returns:

file format suffix with a leading ‘.’ e.g., ‘.bin’

_guess_specfem_flavor()

Guess if SPECFEM2D/3D/3D_GLOBE was used to generate the model based on the format of the file. Check based on unique available files or properties

Return type:

str

Returns:

SPECFEM flavor, one of [‘2D’, ‘3D’, ‘3DGLOBE’]

_read_model_fortran_binary(parameter)

Load Fortran binary models into disk. This is the preferred model format for SeisFlows <-> SPECFEM interaction

Parameters:

parameter (str) – chosen parameter to load model for

Return type:

np.array

Returns:

vector of model values for given parameter

abstractmethod _read_model_adios(parameter)

Load ADIOS models into disk

Parameters:

parameter (str) – chosen parameter to load model for

Return type:

np.array

Returns:

vector of model values for given parameter

_read_model_ascii(parameter)

Load ASCII SPECFEM2D models into disk. ASCII models are generally saved all in a single file with all parameters together as a N column ASCII file where columns 1 and 2 are the coordinates of the mesh, and the remainder columns are data corresponding to the filenames e.g., proc000000_rho_vp_vs.dat, rho is column 3, vp is 4 etc.

Parameters:

parameter (str) – chosen parameter to load model for

Return type:

np.array

Returns:

vector of model values for given parameter

_write_model_fortran_binary(path)

Save a SPECFEM model back to Fortran binary format. Data are written as single precision floating point numbers

Note

FORTRAN unformatted binaries are bounded by an INT*4 byte count. This function mimics that behavior by tacking on the boundary data as ‘int32’ at the top and bottom of the data array. https://docs.oracle.com/cd/E19957-01/805-4939/6j4m0vnc4/index.html