Main solver routine

emg3d.solver.solve(grid, model, sfield, efield=None, cycle='F', sslsolver=False, semicoarsening=False, linerelaxation=False, verb=2, **kwargs)[source]

Solver for 3D CSEM data with tri-axial electrical anisotropy.

The principal solver of emg3d is using the multigrid method as presented in [Muld06]. Multigrid can be used as a standalone solver, or as a preconditioner for an iterative solver from the scipy.sparse.linalg-library, e.g., scipy.sparse.linalg.bicgstab(). Alternatively, these Krylov subspace solvers can also be used without multigrid at all. See the cycle and sslsolver parameters.

Implemented are the F-, V-, and W-cycle schemes for multigrid (cycle parameter), and the amount of smoothing steps (initial smoothing, pre-smoothing, coarsest-grid smoothing, and post-smoothing) can be set individually (nu_init, nu_pre, nu_coarse, and nu_post, respectively). The maximum level of coarsening can be restricted with the clevel parameter.

Semicoarsening and line relaxation, as presented in [Muld07], are implemented, see the semicoarsening and linerelaxation parameters. Using the BiCGSTAB solver together with multigrid preconditioning with semicoarsening and line relaxation is slow but generally the most robust. Not using BiCGSTAB nor semicoarsening nor line relaxation is fast but may fail on stretched grids.

Parameters:
grid : emg3d.meshes.TensorMesh

The grid. See emg3d.meshes.TensorMesh.

model : emg3d.models.Model

The model. See emg3d.models.Model.

sfield : emg3d.fields.SourceField

The source field. See emg3d.fields.get_source_field().

efield : emg3d.fields.Field, optional

Initial electric field. It is initiated with zeroes if not provided. A provided efield MUST have frequency information (initiated with emg3d.fields.Field(..., freq)).

If an initial efield is provided nothing is returned, but the final efield is directly put into the provided efield.

If an initial field is provided and a sslsolver is used, then it first carries out one multigrid cycle without semicoarsening nor line relaxation. The sslsolver is at times unstable with an initial guess, carrying out one MG cycle helps to stabilize it.

cycle : str; optional.

Type of multigrid cycle. Default is ‘F’.

  • ‘V’: V-cycle, simplest version;
  • ‘W’: W-cycle, most expensive version;
  • ‘F’: F-cycle, sort of a compromise between ‘V’ and ‘W’;
  • None: Does not use multigrid, only sslsolver.

If None, sslsolver must be provided, and the sslsolver will be used without multigrid pre-conditioning.

Comparison of V (left), F (middle), and W (right) cycles for the case of four grids (three relaxation and prolongation steps):

 h_
2h_   \    /   \          /   \            /
4h_    \  /     \    /\  /     \    /\    /
8h_     \/       \/\/  \/       \/\/  \/\/
sslsolver : str, optional

A scipy.sparse.linalg-solver, to use with MG as pre-conditioner or on its own (if cycle=None). Default is False.

Current possibilities:

It does currently not work with ‘cg’, ‘bicg’, ‘qmr’, and ‘minres’ for various reasons (e.g., some require rmatvec in addition to matvec).

semicoarsening : int; optional

Semicoarsening. Default is False.

  • True: Cycling over 1, 2, 3.
  • 0 or False: No semicoarsening.
  • 1: Semicoarsening in x direction.
  • 2: Semicoarsening in y direction.
  • 3: Semicoarsening in z direction.
  • Multi-digit number containing digits from 0 to 3. Multigrid will cycle over these values, e.g., semicoarsening=1213 will cycle over [1, 2, 1, 3].
linerelaxation : int; optional

Line relaxation. Default is False.

This parameter is not respected on the coarsest grid, except if it is set to 0. If it is bigger than zero line relaxation on the coarsest grid is carried out along all dimensions which have more than 2 cells.

  • True: Cycling over [4, 5, 6].
  • 0 or False: No line relaxation.
  • 1: line relaxation in x direction.
  • 2: line relaxation in y direction.
  • 3: line relaxation in z direction.
  • 4: line relaxation in y and z directions.
  • 5: line relaxation in x and z directions.
  • 6: line relaxation in x and y directions.
  • 7: line relaxation in x, y, and z directions.
  • Multi-digit number containing digits from 0 to 7. Multigrid will cycle over these values, e.g., linerelaxation=1213 will cycle over [1, 2, 1, 3].

Note: Smoothing is generally done in lexicographical order, except for line relaxation in y direction; the reason is speed (memory access).

verb : int; optional

Level of verbosity (the higher the more verbose). Default is 2.

  • 0: Print nothing.
  • 1: Print warnings.
  • 2: Print runtime and information about the method.
  • 3: Print additional information for each MG-cycle.
  • 4: Print everything (slower due to additional error calculations).
  • -1: Print one-liner (dynamically updated).
**kwargs : Optional solver options:
  • tol : float

    Convergence tolerance. Default is 1e-6.

    Iterations stop as soon as the norm of the residual has decreased by this factor, relative to the residual norm obtained for a zero electric field.

  • maxit : int

    Maximum number of multigrid iterations. Default is 50.

    If sslsolver is used, this applies to the sslsolver.

    In the case that multigrid is used as a pre-conditioner for the sslsolver, the maximum iteration for multigrid is defined by the maximum length of the linerelaxation and semicoarsening-cycles.

  • nu_init : int

    Number of initial smoothing steps, before MG cycle. Default is 0.

  • nu_pre : int

    Number of pre-smoothing steps. Default is 2.

  • nu_coarse : int

    Number of smoothing steps on coarsest grid. Default is 1.

  • nu_post : int

    Number of post-smoothing steps. Default is 2.

  • clevel : int

    The maximum coarsening level can be different for each dimension and is, by default, automatically determined (clevel=-1). The parameter clevel can be used to restrict the maximum coarsening level in any direction by its value. Default is -1.

  • return_info : bool

    If True, a dictionary is returned with runtime info (final norm and number of iterations of MG and the sslsolver).

Returns:
efield : emg3d.fields.Field

Resulting electric field. Is not returned but replaced in-place if an initial efield was provided.

info_dict : dict

Dictionary with runtime info; only if return_info=True.

Keys:

  • exit: Exit status, 0=Success, 1=Failure;
  • exit_message: Exit message, check this if exit=1;
  • abs_error: Absolute error;
  • rel_error: Relative error;
  • ref_error: Reference error [norm(sfield)];
  • tol: Tolerance (abs_error<ref_error*tol);
  • it_mg: Number of multigrid iterations;
  • it_ssl: Number of SSL iterations;
  • time: Runtime (s).
  • runtime_at_cycle: Runtime after each cycle (s).
  • error_at_cycle: Absolute error after each cycle.

Examples

>>> import emg3d
>>> import numpy as np
>>> # Create a simple grid, 8 cells of length 1 in each direction,
>>> # starting at the origin.
>>> grid = emg3d.meshes.TensorMesh(
>>>         [np.ones(8), np.ones(8), np.ones(8)],
>>>         x0=np.array([0, 0, 0]))
>>> # The model is a fullspace with tri-axial anisotropy.
>>> model = emg3d.models.Model(grid, res_x=1.5, res_y=1.8, res_z=3.3)
>>> # The source is a x-directed, horizontal dipole at (4, 4, 4)
>>> # with a frequency of 10 Hz.
>>> sfield = emg3d.fields.get_source_field(
>>>         grid, src=[4, 4, 4, 0, 0], freq=10)
>>> # Calculate the electric signal.
>>> efield = emg3d.solve(grid, model, sfield, verb=3)
>>> # Get the corresponding magnetic signal.
>>> hfield = emg3d.fields.get_h_field(grid, model, efield)
.
:: emg3d START :: 10:27:25 :: v0.9.1
.
   MG-cycle       : 'F'                 sslsolver : False
   semicoarsening : False [0]           tol       : 1e-06
   linerelaxation : False [0]           maxit     : 50
   nu_{i,1,c,2}   : 0, 2, 1, 2          verb      : 3
   Original grid  :   8 x   8 x   8     => 512 cells
   Coarsest grid  :   2 x   2 x   2     => 8 cells
   Coarsest level :   2 ;   2 ;   2
.
   [hh:mm:ss]  rel. error                  [abs. error, last/prev]   l s
.
       h_
      2h_ \    /
      4h_  \/\/
.
   [10:27:25]   2.284e-02  after   1 F-cycles   [1.275e-06, 0.023]   0 0
   [10:27:25]   1.565e-03  after   2 F-cycles   [8.739e-08, 0.069]   0 0
   [10:27:25]   1.295e-04  after   3 F-cycles   [7.232e-09, 0.083]   0 0
   [10:27:25]   1.197e-05  after   4 F-cycles   [6.685e-10, 0.092]   0 0
   [10:27:25]   1.233e-06  after   5 F-cycles   [6.886e-11, 0.103]   0 0
   [10:27:25]   1.415e-07  after   6 F-cycles   [7.899e-12, 0.115]   0 0
.
   > CONVERGED
   > MG cycles        : 6
   > Final rel. error : 1.415e-07
.
:: emg3d END   :: 10:27:25 :: runtime = 0:00:00