solve#

emg3d.solver.solve(model, sfield, sslsolver=True, semicoarsening=True, linerelaxation=True, verb=0, **kwargs)[source]#

Solver for three-dimensional electromagnetic diffusion.

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, semicoarsening, and line relaxation is generally the most robust solution, albeit not necessarily the fastest. It is the default setting for its robustness. Just using traditional multigrid without BiCGSTAB nor semicoarsening nor line relaxation uses the least memory and is often faster but may fail on stretched grids or for models with strong anisotropy. However, only testing the parameters for a given model can give certainty to which parameters are best.

Parameters
modelModel

The model; a emg3d.models.Model instance.

sfieldField

The source field. See emg3d.fields.get_source_field.

sslsolver{str, bool}, default: True

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

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, bool}, default: True

Semicoarsening.

  • 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, bool}, default: True

Line relaxation.

  • 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).

verbint, default: 0

Level of verbosity (the higher the more verbose).

  • -1: Nothing.

  • 0: Warnings.

  • 1: One-liner at the end.

  • 2: One-liner, dynamically updated.

  • 3: Information about the method plus dynamic one-liner.

  • 4: Additional information for each MG-cycle.

  • 5: Everything (slower due to additional error computations).

cycle{str, None}, default: ‘F’

Type of multigrid cycle.

  • 'V': V-cycle, simplest version.

  • 'W': W-cycle, most expensive version.

  • 'F': F-cycle, sort of a compromise between V- and W-cycle.

  • None: Does not use multigrid, only the chosen sslsolver.

A sslsolver has to be provided if cycle=None, and the sslsolver will then 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_     \/       \/\/  \/       \/\/  \/\/
efieldField, default: None

Initial electric field. If is initiated with zeroes if not provided.

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 multigrid cycle helps to stabilize it.

Note that the tangential field at the boundary of a provided efield is set to zero to ensure a PEC boundary (perfect electric conductor).

tolfloat, default: 1e-6

Convergence tolerance.

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.

maxitint, default: 50

Maximum number of multigrid iterations.

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_initint, default: 0

Number of initial smoothing steps, before multigrid cycle.

nu_preint, default: 2

Number of pre-smoothing steps.

nu_coarseint, default: 1

Number of smoothing steps on coarsest grid.

nu_postint, default: 2

Number of post-smoothing steps.

clevelint, default: -1

The maximum coarsening level can be different for each dimension and is, by default, automatically determined (clevel=-1). The parameter clevel restricts the maximum coarsening level by its value.

return_infobool, default: False

If True, a dictionary is returned with runtime info (final norm, number of iterations of multigrid and the sslsolver, log, exit message, etc).

logint, default: 1

Only relevant if return_info=True.

  • -1: LOG ONLY: Only store info in log, do not print on screen.

  • 0: SCREEN only: Only print info to screen, do not store in log.

  • 1: BOTH: Store info in log and print on screen.

plainbool, default: False

Plain multigrid method. This is a shortcut for sslsolver=False, semicoarsening=False, linerelaxation=False. The three parameters remain unchanged if they are set to anything else than True.

Returns
efieldField, returned if efield=None

Resulting electric field. It is not returned but stored in-place if an initial efield was provided.

info_dictdict, returned if return_info=True

Dictionary with solver info. 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.

  • log: Stored log.

Examples

In [1]: import emg3d
   ...: import numpy as np
   ...: 

In [2]: # Create a simple grid, 8 cells of length 1 in each direction,
   ...: # starting at the origin.
   ...: hx = np.ones(8)
   ...: grid = emg3d.TensorMesh([hx, hx, hx], origin=(0, 0, 0))
   ...: 

In [3]: # Create a fullspace model with triaxial anisotropy.
   ...: model = emg3d.Model(grid, property_x=1.5, property_y=1.8,
   ...:                     property_z=3.3, mapping='Resistivity')
   ...: 

In [4]: # The source is a x-directed, horizontal dipole at (4, 4, 4)
   ...: # with a frequency of 10 Hz.
   ...: coo = (4, 4, 4, 0, 0)  # (x, y, z, azimuth, elevation)
   ...: sfield = emg3d.fields.get_source_field(
   ...:             grid, source=coo, frequency=10)
   ...: 

In [5]: # Solve for the electric field.
   ...: efield = emg3d.solve(model, sfield, plain=True, verb=4)
   ...: 

:: emg3d START :: 18:14:27 :: v1.5.0

   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      : 4
   Original grid  :   8 x   8 x   8     => 512 cells
   Coarsest grid  :   2 x   2 x   2     => 8 cells
   Coarsest level :   2 ;   2 ;   2     :: Grid not optimal for MG solver ::

   [hh:mm:ss]  rel. error                  [abs. error, last/prev]   l s

       h_
      2h_ \    /
      4h_  \/\/ 

   [18:14:39]   2.284e-02  after   1 F-cycles   [1.275e-06, 0.023]   0 0
   [18:14:39]   1.565e-03  after   2 F-cycles   [8.739e-08, 0.069]   0 0
   [18:14:39]   1.295e-04  after   3 F-cycles   [7.232e-09, 0.083]   0 0
   [18:14:39]   1.197e-05  after   4 F-cycles   [6.685e-10, 0.092]   0 0
   [18:14:39]   1.233e-06  after   5 F-cycles   [6.886e-11, 0.103]   0 0
   [18:14:39]   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   :: 18:14:39 :: runtime = 0:00:12