Add-ons

fdesign – Digital Linear Filter (DLF) design

The add-on fdesign can be used to design digital linear filters for the Hankel or Fourier transform, or for any linear transform ([Ghos70]). For this included or provided theoretical transform pairs can be used. Alternatively, one can use the EM modeller empymod to use the responses to an arbitrary 1D model as numerical transform pair.

More information can be found in the following places:

This filter designing tool uses the direct matrix inversion method as described in [Kong07] and is based on scripts by [Key12]. The whole project of fdesign started with the Matlab scripts from Kerry Key, which he used to design his filters for [Key09], [Key12]. Fruitful discussions with Evert Slob and Kerry Key improved the add-on substantially.

Note that the use of empymod to create numerical transform pairs is, as of now, only implemented for the Hankel transform.

Implemented analytical transform pairs

The following tables list the transform pairs which are implemented by default. Any other transform pair can be provided as input. A transform pair is defined in the following way:

from empymod.scripts import fdesign

def my_tp_pair(var):
    '''My transform pair.'''

    def lhs(l):
        return func(l, var)

    def rhs(r):
        return func(r, var)

    return fdesign.Ghosh(name, lhs, rhs)

Here, name must be one of j0, j1, sin, or cos, depending what type of transform pair it is. Additional variables are provided with var. The evaluation points of the lhs are denoted by l, and the evaluation points of the rhs are denoted as r. As an example here the implemented transform pair j0_1

def j0_1(a=1):
    '''Hankel transform pair J0_1 ([Ande75]_).'''

    def lhs(l):
        return l*np.exp(-a*l**2)

    def rhs(r):
        return np.exp(-r**2/(4*a))/(2*a)

    return Ghosh('j0', lhs, rhs)

Implemented Hankel transforms

  • j0_1 [Ande75]

    \[\int^\infty_0 l \exp\left(-al^2\right) J_0(lr) dl = \frac{\exp\left(\frac{-r^2}{4a}\right)}{2a}\]
  • j0_2 [Ande75]

    \[\int^\infty_0 \exp\left(-al\right) J_0(lr) dl = \frac{1}{\sqrt{a^2+r^2}}\]
  • j0_3 [GuSi97]

    \[\int^\infty_0 l\exp\left(-al\right) J_0(lr) dl = \frac{a}{(a^2 + r^2)^{3/2}}\]
  • j0_4 [ChCo82]

    \[\int^\infty_0 \frac{l}{\beta} \exp\left(-\beta z_v \right) J_0(lr) dl = \frac{\exp\left(-\gamma R\right)}{R}\]
  • j0_5 [ChCo82]

    \[\int^\infty_0 l \exp\left(-\beta z_v \right) J_0(lr) dl = \frac{ z_v (\gamma R + 1)}{R^3}\exp\left(-\gamma R\right)\]
  • j1_1 [Ande75]

    \[\int^\infty_0 l^2 \exp\left(-al^2\right) J_1(lr) dl = \frac{r}{4a^2} \exp\left(-\frac{r^2}{4a}\right)\]
  • j1_2 [Ande75]

    \[\int^\infty_0 \exp\left(-al\right) J_1(lr) dl = \frac{\sqrt{a^2+r^2}-a}{r\sqrt{a^2 + r^2}}\]
  • j1_3 [Ande75]

    \[\int^\infty_0 l \exp\left(-al\right) J_1(lr) dl = \frac{r}{(a^2 + r^2)^{3/2}}\]
  • j1_4 [ChCo82]

    \[\int^\infty_0 \frac{l^2}{\beta} \exp\left(-\beta z_v \right) J_1(lr) dl = \frac{r(\gamma R+1)}{R^3}\exp\left(-\gamma R\right)\]
  • j1_5 [ChCo82]

    \[\int^\infty_0 l^2 \exp\left(-\beta z_v \right) J_1(lr) dl = \frac{r z_v(\gamma^2R^2+3\gamma R+3)}{R^5}\exp\left(-\gamma R\right)\]

Where

\[a >0, r>0\]
\[z_v = |z_{rec} - z_{src}|\]
\[R = \sqrt{r^2 + z_v^2}\]
\[\gamma = \sqrt{2j\pi\mu_0f/\rho}\]
\[\beta = \sqrt{l^2 + \gamma^2}\]

Implemented Fourier transforms

  • sin_1 [Ande75]

    \[\int^\infty_0 l\exp\left(-a^2l^2\right) \sin(lr) dl = \frac{\sqrt{\pi}r}{4a^3} \exp\left(-\frac{r^2}{4a^2}\right)\]
  • sin_2 [Ande75]

    \[\int^\infty_0 \exp\left(-al\right) \sin(lr) dl = \frac{r}{a^2 + r^2}\]
  • sin_3 [Ande75]

    \[\int^\infty_0 \frac{l}{a^2+l^2} \sin(lr) dl = \frac{\pi}{2} \exp\left(-ar\right)\]
  • cos_1 [Ande75]

    \[\int^\infty_0 \exp\left(-a^2l^2\right) \cos(lr) dl = \frac{\sqrt{\pi}}{2a} \exp\left(-\frac{r^2}{4a^2}\right)\]
  • cos_2 [Ande75]

    \[\int^\infty_0 \exp\left(-al\right) \cos(lr) dl = \frac{a}{a^2 + r^2}\]
  • cos_3 [Ande75]

    \[\int^\infty_0 \frac{1}{a^2+l^2} \cos(lr) dl = \frac{\pi}{2a} \exp\left(-ar\right)\]
empymod.scripts.fdesign.design(n, spacing, shift, fI, fC=False, r=None, r_def=(1, 1, 2), reim=None, cvar='amp', error=0.01, name=None, full_output=False, finish=False, save=True, path='filters', verb=2, plot=1)[source]

Digital linear filter (DLF) design

This routine can be used to design digital linear filters for the Hankel or Fourier transform, or for any linear transform ([Ghos70]). For this included or provided theoretical transform pairs can be used. Alternatively, one can use the EM modeller empymod to use the responses to an arbitrary 1D model as numerical transform pair.

This filter designing tool uses the direct matrix inversion method as described in [Kong07] and is based on scripts by [Key12]. The tool is an add-on to the electromagnetic modeller empymod [Wert17]. Fruitful discussions with Evert Slob and Kerry Key improved the add-on substantially.

Example notebooks of its usage can be found in the repo github.com/empymod/empymod-examples.

Parameters:
n : int

Filter length.

spacing: float or tuple (start, stop, num)

Spacing between filter points. If tuple, it corresponds to the input for np.linspace with endpoint=True.

shift: float or tuple (start, stop, num)

Shift of base from zero. If tuple, it corresponds to the input for np.linspace with endpoint=True.

fI, fC : transform pairs

Theoretical or numerical transform pair(s) for the inversion (I) and for the check of goodness (fC). fC is optional. If not provided, fI is used for both fI and fC.

r : array, optional

Right-hand side evaluation points for the check of goodness (fC). Defaults to r = np.logspace(0, 5, 1000), which are a lot of evaluation points, and depending on the transform pair way too long r’s.

r_def : tuple (add_left, add_right, factor), optional

Definition of the right-hand side evaluation points r of the inversion. r is derived from the base values, default is (1, 1, 2).

  • rmin = log10(1/max(base)) - add_left
  • rmax = log10(1/min(base)) + add_right
  • r = logspace(rmin, rmax, factor*n)
reim : np.real or np.imag, optional

Which part of complex transform pairs is used for the inversion. Defaults to np.real.

cvar : string {‘amp’, ‘r’}, optional

If ‘amp’, the inversion minimizes the amplitude. If ‘r’, the inversion maximizes the right-hand side evaluation point r. Defaults is ‘amp’.

error : float, optional

Up to which relative error the transformation is considered good in the evaluation of the goodness. Default is 0.01 (1 %).

name : str, optional

Name of the filter. Defaults to dlf_+str(n).

full_output : bool, optional

If True, returns best filter and output from scipy.optimize.brute; else only filter. Default is False.

finish : None, True, or callable, optional

If callable, it is passed through to scipy.optimize.brute: minimization function to find minimize best result from brute-force approach. Default is None. You can simply provide True in order to use scipy.optimize.fmin_powell(). Set this to None if you are only interested in the actually provided spacing/shift-values.

save : bool, optional

If True, best filter is saved to plain text files in ./filters/. Can be loaded with fdesign.load_filter(name). If full, the inversion output is stored too. You can add ‘.gz’ to name, which will then save the full inversion output in a compressed file instead of plain text.

path : string, optional

Absolute or relative path where output will be saved if save=True. Default is ‘filters’.

verb : {0, 1, 2}, optional
Level of verbosity, default is 2:
  • 0: Print nothing.
  • 1: Print warnings.
  • 2: Print additional time, progress, and result
plot : {0, 1, 2, 3}, optional
Level of plot-verbosity, default is 1:
  • 0: Plot nothing.
  • 1: Plot brute-force result
  • 2: Plot additional theoretical transform pairs, and best inv.
  • 3: Plot additional inversion result
    (can result in lots of plots depending on spacing and shift) If you are using a notebook, use %matplotlib notebook to have all inversion results appear in the same plot.
Returns:
filter : empymod.filter.DigitalFilter instance

Best filter for the input parameters.

full : tuple

Output from scipy.optimize.brute with full_output=True. (Returned when full_output is True.)

empymod.scripts.fdesign.save_filter(name, filt, full=None, path='filters')[source]

Save DLF-filter and inversion output to plain text files.

empymod.scripts.fdesign.load_filter(name, full=False, path='filters')[source]

Load saved DLF-filter and inversion output from text files.

empymod.scripts.fdesign.plot_result(filt, full, prntres=True)[source]

QC the inversion result.

Parameters:
- filt, full as returned from fdesign.design with full_output=True
- If prntres is True, it calls fdesign.print_result as well.
r
empymod.scripts.fdesign.print_result(filt, full=None)[source]

Print best filter information.

Parameters:
- filt, full as returned from fdesign.design with full_output=True
class empymod.scripts.fdesign.Ghosh(name, lhs, rhs)[source]

Simple Class for Theoretical Transform Pairs.

Named after D. P. Ghosh, honouring his 1970 Ph.D. thesis with which he introduced the digital filter method to geophysics ([Ghos70]).

empymod.scripts.fdesign.j0_1(a=1)[source]

Hankel transform pair J0_1 ([Ande75]).

empymod.scripts.fdesign.j0_2(a=1)[source]

Hankel transform pair J0_2 ([Ande75]).

empymod.scripts.fdesign.j0_3(a=1)[source]

Hankel transform pair J0_3 ([GuSi97]).

empymod.scripts.fdesign.j0_4(f=1, rho=0.3, z=50)[source]

Hankel transform pair J0_4 ([ChCo82]).

Parameters:
f : float

Frequency (Hz)

rho : float

Resistivity (Ohm.m)

z : float

Vertical distance between source and receiver (m)

empymod.scripts.fdesign.j0_5(f=1, rho=0.3, z=50)[source]

Hankel transform pair J0_5 ([ChCo82]).

Parameters:
f : float

Frequency (Hz)

rho : float

Resistivity (Ohm.m)

z : float

Vertical distance between source and receiver (m)

empymod.scripts.fdesign.j1_1(a=1)[source]

Hankel transform pair J1_1 ([Ande75]).

empymod.scripts.fdesign.j1_2(a=1)[source]

Hankel transform pair J1_2 ([Ande75]).

empymod.scripts.fdesign.j1_3(a=1)[source]

Hankel transform pair J1_3 ([Ande75]).

empymod.scripts.fdesign.j1_4(f=1, rho=0.3, z=50)[source]

Hankel transform pair J1_4 ([ChCo82]).

Parameters:
f : float

Frequency (Hz)

rho : float

Resistivity (Ohm.m)

z : float

Vertical distance between source and receiver (m)

empymod.scripts.fdesign.j1_5(f=1, rho=0.3, z=50)[source]

Hankel transform pair J1_5 ([ChCo82]).

Parameters:
f : float

Frequency (Hz)

rho : float

Resistivity (Ohm.m)

z : float

Vertical distance between source and receiver (m)

empymod.scripts.fdesign.sin_1(a=1)[source]

Fourier sine transform pair sin_1 ([Ande75]).

empymod.scripts.fdesign.sin_2(a=1)[source]

Fourier sine transform pair sin_2 ([Ande75]).

empymod.scripts.fdesign.sin_3(a=1)[source]

Fourier sine transform pair sin_3 ([Ande75]).

empymod.scripts.fdesign.cos_1(a=1)[source]

Fourier cosine transform pair cos_1 ([Ande75]).

empymod.scripts.fdesign.cos_2(a=1)[source]

Fourier cosine transform pair cos_2 ([Ande75]).

empymod.scripts.fdesign.cos_3(a=1)[source]

Fourier cosine transform pair cos_3 ([Ande75]).

empymod.scripts.fdesign.empy_hankel(ftype, zsrc, zrec, res, freqtime, depth=None, aniso=None, epermH=None, epermV=None, mpermH=None, mpermV=None, htarg=None, verblhs=0, verbrhs=0)[source]

Numerical transform pair with empymod.

All parameters except ftype, verblhs, and verbrhs correspond to the input parameters to empymod.dipole. See there for more information.

Note that if depth=None or [], the analytical full-space solutions will be used (much faster).

Parameters:
ftype : str or list of strings

Either of: {‘j0’, ‘j1’, ‘j2’, [‘j0’, ‘j1’]}

  • ‘j0’: Analyze J0-term with ab=11, angle=45°
  • ‘j1’: Analyze J1-term with ab=31, angle=0°
  • ‘j2’: Analyze J0- and J1-terms jointly with ab=12, angle=45°
  • [‘j0’, ‘j1’]: Same as calling empy_hankel twice, once with ‘j0’ and
    one with ‘j1’; can be provided like this to fdesign.design.
verblhs, verbrhs: int

verb-values provided to empymod for lhs and rhs.

Note that ftype=’j2’ only works for fC, not for fI.

tmtemod – Calculate up- and down-going TM and TE modes

This add-on for empymod adjusts [HuTS15] for TM/TE-split. The development was initiated by the development of https://github.com/empymod/csem-ziolkowski-and-slob ([ZiSl19]).

This is a stripped-down version of empymod with a lot of simplifications but an important addition. The modeller empymod returns the total field, hence not distinguishing between TM and TE mode, and even less between up- and down-going fields. The reason behind this is simple: The derivation of [HuTS15], on which empymod is based, returns the total field. In this derivation each mode (TM and TE) contains non-physical contributions. The non-physical contributions have opposite signs in TM and TE, so they cancel each other out in the total field. However, in order to obtain the correct TM and TE contributions one has to remove these non-physical parts.

This is what this routine does, but only for an x-directed electric source with an x-directed electric receiver, and in the frequency domain (src and rec in same layer). This version of dipole returns the signal separated into TM++, TM+-, TM-+, TM–, TE++, TE+-, TE-+, and TE– as well as the direct field TM and TE contributions. The first superscript denotes the direction in which the field diffuses towards the receiver and the second superscript denotes the direction in which the field diffuses away from the source. For both the plus-sign indicates the field diffuses in the downward direction and the minus-sign indicates the field diffuses in the upward direction. It uses empymod wherever possible. See the corresponding functions in empymod for more explanation and documentation regarding input parameters. There are important limitations:

  • ab == 11 [=> x-directed el. source & el. receivers]
  • signal == None [=> only frequency domain]
  • xdirect == False [=> direct field calc. in wavenr-domain]
  • ht == ‘fht’
  • htarg == ‘key_201_2012’
  • Options ft, ftarg, opt, and loop are not available.
  • lsrc == lrec [=> src & rec are assumed in same layer!]
  • Model must have more than 1 layer
  • Electric permittivity and magnetic permeability are isotropic.
  • Only one frequency at once.

Theory

The derivation of [HuTS15], on which empymod is based, returns the total field. Internally it also calculates TM and TE modes, and sums these up. However, the separation into TM and TE mode introduces a singularity at \(\kappa = 0\). It has no contribution in the space-frequency domain to the total fields, but it introduces non-physical events in each mode with opposite sign (so they cancel each other out in the total field). In order to obtain the correct TM and TE contributions one has to remove these non-physical parts.

To remove the non-physical part we use the file tmtemod.py in this directory. This routine is basically a heavily simplified version of empymod with the following limitations outlined above.

So tmtemod.py returns the signal separated into TM++, TM+-, TM-+, TM–, TE++, TE+-, TE-+, and TE– as well as the direct field TM and TE contributions. The first superscript denotes the direction in which the field diffuses towards the receiver and the second superscript denotes the direction in which the field diffuses away from the source. For both the plus-sign indicates the field diffuses in the downward direction and the minus-sign indicates the field diffuses in the upward direction. The routine uses empymod wherever possible, see the corresponding functions in empymod for more explanation and documentation regarding input parameters.

Please note that the notation in [HuTS15] differs from the notation in [ZiSl19]. I specify therefore always, which notification applies, either Hun15 or Zio19.

We start with equation (105) in Hun15:

\[\hat{G}^{ee}_{xx}(\boldsymbol{x}, \boldsymbol{x'}, \omega) = \hat{G}^{ee;i}_{xx;s}(\boldsymbol{x}-\boldsymbol{x'}, \omega) + \frac{1}{8\pi}\int^\infty_{\kappa=0} \left(\frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s} - \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s}\right) J_0(\kappa r)\kappa d \kappa\]
\[- \frac{\cos(2\phi)}{8\pi}\int^\infty_{\kappa=0} \left(\frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s} + \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s}\right) J_2(\kappa r)\kappa d \kappa .\]

Ignoring the incident field, and using \(J_2 = \frac{2}{\kappa r}J_1 - J_0\) to avoid \(J_2\)-integrals, we get

\[\hat{G}^{ee}_{xx}(\boldsymbol{x}, \boldsymbol{x'}, \omega) = \frac{1}{8\pi}\int^\infty_{\kappa=0} \left(\frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s}- \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s}\right) J_0(\kappa r)\kappa\,{\mathrm{d}}\kappa\]
\[+ \frac{\cos(2\phi)}{8\pi}\int^\infty_{\kappa=0} \left(\frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s} + \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s}\right) J_0(\kappa r)\kappa\,{\mathrm{d}}\kappa\]
\[- \frac{\cos(2\phi)}{4\pi r}\int^\infty_{\kappa=0} \left(\frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s} + \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s}\right) J_1(\kappa r)\,{\mathrm{d}}\kappa .\]

From this the TM- and TE-parts follow as

\[{\mathrm{TE}} = \frac{\cos(2\phi)-1}{8\pi}\int^\infty_{\kappa=0} \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s} J_0(\kappa r)\kappa\,{\mathrm{d}}\kappa - \frac{\cos(2\phi)}{4\pi r}\int^\infty_{\kappa=0} \frac{\zeta_s \tilde{g}^{te}_{zz;s}}{\bar{\Gamma}_s} J_1(\kappa r)\,{\mathrm{d}}\kappa ,\]
\[ {\mathrm{TM}} = \frac{\cos(2\phi)+1}{8\pi}\int^\infty_{\kappa=0} \frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s} J_0(\kappa r)\kappa\,{\mathrm{d}}\kappa - \frac{\cos(2\phi)}{4\pi r}\int^\infty_{\kappa=0} \frac{\Gamma_s \tilde{g}^{tm}_{hh;s}}{\eta_s} J_1(\kappa r)\,{\mathrm{d}}\kappa .\]

Equations (108) and (109) in Hun15 yield the required parameters \(\tilde{g}^{tm}_{hh;s}\) and \(\tilde{g}^{te}_{zz;s}\),

\[\tilde{g}^{tm}_{hh;s} = P^{u-}_s W^u_s + P^{d-}_s W^d_s ,\]
\[\tilde{g}^{te}_{zz;s} = \bar{P}^{u+}_s \bar{W}^u_s + \bar{P}^{d+}_s \bar{W}^d_s \ .\]

The parameters \(P^{u\pm}_s\) and \(P^{d\pm}_s\) are given in equations (81) and (82), \(\bar{P}^{u\pm}_s\) and \(\bar{P}^{d\pm}_s\) in equations (A-8) and (A-9); \(W^u_s\) and \(W^d_s\) in equation (74) in Hun15. This yields

\[\tilde{g}^{te}_{zz;s} = \frac{\bar{R}_s^+}{\bar{M}_s}\left\{\exp[-\bar{\Gamma}_s(z_s-z+d^+)] + \bar{R}_s^-\exp[-\bar{\Gamma}_s(z_s-z+d_s+d^-)]\right\}\]
\[+ \frac{\bar{R}_s^-}{\bar{M}_s} \left\{\exp[-\bar{\Gamma}_s(z-z_{s-1}+d^-)]+ \bar{R}_s^+\exp[-\bar{\Gamma}_s(z-z_{s-1}+d_s+d^+)]\right\} ,\]
\[=\frac{\bar{R}_s^+}{\bar{M}_s}\left\{\exp[-\bar{\Gamma}_s(2z_s-z-z')] + \bar{R}_s^-\exp[-\bar{\Gamma}_s(z'-z+2d_s)]\right\}\]
\[+ \frac{\bar{R}_s^-}{\bar{M}_s} \left\{\exp[-\bar{\Gamma}_s(z+z'-2z_{s-1})]+ \bar{R}_s^+\exp[-\bar{\Gamma}_s(z-z'+2d_s)]\right\} ,\]

where \(d^\pm\) is taken from the text below equation (67). There are four terms in the right-hand side, two in the first line and two in the second line. The first term in the first line is the integrand of TE+-, the second term in the first line corresponds to TE++, the first term in the second line is TE-+, and the second term in the second line is TE–.

If we look at TE+-, we have

\[\tilde{g}^{te+-}_{zz;s} = \frac{\bar{R}_s^+}{\bar{M}_s}\exp[-\bar{\Gamma}_s(2z_s-z-z')] \ ,\]

and therefore

\[{\mathrm{TE}}^{+-} = \frac{\cos(2\phi)-1}{8\pi}\int^\infty_{\kappa=0} \frac{\zeta_s \bar{R}_s^+}{\bar{\Gamma}_s\bar{M}_s} \exp[-\bar{\Gamma}_s(2z_s-z-z')] J_0(\kappa r)\kappa\,{\mathrm{d}}\kappa\]
\[- \frac{\cos(2\phi)}{4\pi r}\int^\infty_{\kappa=0} \frac{\zeta_s \bar{R}_s^+}{\bar{\Gamma}_s\bar{M}_s} \exp[-\bar{\Gamma}_s(2z_s-z-z')] J_1(\kappa r)\,{\mathrm{d}}\kappa .\]

We can compare this to equation (4.165) in Zio19, with \(\hat{I}^e_x=1\) and slightly re-arranging it to look more alike, we get

\[\hat{E}^{+-}_{xx;H} = \frac{y^2}{4\pi r^2} \int^\infty_{\kappa=0} \frac{\zeta_1}{\Gamma_1} \frac{R^-_{H;1}}{M_{H;1}} \exp(-\Gamma_1 h^{+-})J_0(\kappa r)\kappa d\kappa\]
\[+ \frac{x^2-y^2}{4\pi r^3} \int^\infty_{\kappa=0} \frac{\zeta_1}{\Gamma_1} \left(\frac{R^-_{H;1}}{M_{H;1}} - \frac{R^-_{H;1}(\kappa=0)}{M_{H;1}(\kappa=0)}\right) \exp(-\Gamma_1 h^{+-})J_1(\kappa r) d\kappa\]
\[- \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^-_{H;1}(\kappa=0)}{M_{H;1}(\kappa=0)} \exp(-\gamma_1 R^{+-}) .\]

The notation in this equation follows Zio19.

The difference between the two previous equations is that the first one contains non-physical contributions. These have opposite signs in TM+- and TE+-, and therefore cancel each other out. But if we want to know the specific contributions from TM and TE we have to remove them. The non-physical contributions only affect the \(J_1\)-integrals, and only for \(\kappa = 0\).

The following lists for all 8 cases the term that has to be removed, in the notation of Zio19 (for the notation as in Hun15 see the implementation in tmtemod.py):

\[TE^{++} = + \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{\exp(-\gamma_1 |h^-|) }{M_{H;1}(\kappa=0)} ,\]
\[TE^{-+} = - \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^+_{H;1}(\kappa=0)\exp(-\gamma_1 h^{-+}) }{M_{H;1}(\kappa=0)} ,\]
\[TE^{+-} = - \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^-_{H;1}(\kappa=0)\exp(-\gamma_1 h^{+-}) }{M_{H;1}(\kappa=0)} ,\]
\[TE^{--} = + \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^+_{H;1}(\kappa=0)R^-_{H;1}(\kappa=0)\exp(-\gamma_1 h^{--}) } {M_{H;1}(\kappa=0)} ,\]
\[TM^{++} = - \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{\exp(-\gamma_1 |h^-|) }{M_{V;1}(\kappa=0)} ,\]
\[TM^{-+} = - \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^+_{V;1}(\kappa=0)\exp(-\gamma_1 h^{-+}) }{M_{V;1}(\kappa=0)} ,\]
\[TM^{+-} = - \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^-_{V;1}(\kappa=0)\exp(-\gamma_1 h^{+-}) }{M_{V;1}(\kappa=0)} ,\]
\[TM^{--} = - \frac{\zeta_1 (x^2-y^2)}{4\pi\gamma_1 r^4} \frac{R^+_{V;1}(\kappa=0)R^-_{V;1}(\kappa=0)\exp(-\gamma_1 h^{--}) } {M_{V;1}(\kappa=0)} .\]

Note that in the first and fourth equations the correction terms have opposite sign as those in the fifth and eighth equations because at \(\kappa=0\) the TM and TE mode correction terms are equal. Also note that in the second and third equations the correction terms have the same sign as those in the sixth and seventh equations because at \(\kappa=0\) the TM and TE mode reflection responses in those terms are equal but with opposite sign: \(R^\pm_{V;1}(\kappa=0) = -R^\pm_{V;1}(\kappa=0)\).

Hun15 uses \(\phi\), whereas Zio19 uses \(x\), \(y\), for which we can use

\[\cos(2\phi) = -\frac{x^2-y^2}{r^2} \ .\]
empymod.scripts.tmtemod.dipole(src, rec, depth, res, freqtime, aniso=None, eperm=None, mperm=None, verb=2)[source]

Return the electromagnetic field due to a dipole source.

This is a modified version of empymod.model.dipole(). It returns the separated contributions of TM–, TM-+, TM+-, TM++, TMdirect, TE–, TE-+, TE+-, TE++, and TEdirect.

Parameters:
src, rec : list of floats or arrays

Source and receiver coordinates (m): [x, y, z]. The x- and y-coordinates can be arrays, z is a single value. The x- and y-coordinates must have the same dimension.

Sources or receivers placed on a layer interface are considered in the upper layer.

Sources and receivers must be in the same layer.

depth : list

Absolute layer interfaces z (m); #depth = #res - 1 (excluding +/- infinity).

res : array_like

Horizontal resistivities rho_h (Ohm.m); #res = #depth + 1.

freqtime : float

Frequency f (Hz). (The name freqtime is kept for consistency with empymod.model.dipole(). Only one frequency at once.

aniso : array_like, optional

Anisotropies lambda = sqrt(rho_v/rho_h) (-); #aniso = #res. Defaults to ones.

eperm : array_like, optional

Relative electric permittivities epsilon (-); #eperm = #res. Default is ones.

mperm : array_like, optional

Relative magnetic permeabilities mu (-); #mperm = #res. Default is ones.

verb : {0, 1, 2, 3, 4}, optional
Level of verbosity, default is 2:
  • 0: Print nothing.
  • 1: Print warnings.
  • 2: Print additional runtime and kernel calls
  • 3: Print additional start/stop, condensed parameter information.
  • 4: Print additional full parameter information
Returns:
TM, TE : list of ndarrays, (nfreq, nrec, nsrc)

Frequency-domain EM field [V/m], separated into TM = [TM–, TM-+, TM+-, TM++, TMdirect] and TE = [TE–, TE-+, TE+-, TE++, TEdirect].

However, source and receiver are normalised. So the source strength is 1 A and its length is 1 m. Therefore the electric field could also be written as [V/(A.m2)].

The shape of EM is (nfreq, nrec, nsrc). However, single dimensions are removed.