def getkernel(i, inplambd, inpoff, inpfang): r"""Return wavenumber-domain-kernel as a fct of interval i.""" # Indices and factor for this interval iB = i*nquad + np.arange(nquad) # PJ0 and PJ1 for this interval PJ0, PJ1, PJ0b = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH[None, :], etaV[None, :], zetaH[None, :], zetaV[None, :], np.atleast_2d(inplambd)[:, iB], ab, xdirect, msrc, mrec) # Carry out and return the Hankel transform for this interval gEM = np.zeros_like(inpoff, dtype=np.complex128) if k_used[1]: gEM += inpfang*np.dot(PJ1[0, :], BJ1[iB]) if ab in [11, 12, 21, 22, 14, 24, 15, 25]: # Because of J2 # J2(kr) = 2/(kr)*J1(kr) - J0(kr) gEM /= np.atleast_1d(inpoff) if k_used[2]: gEM += inpfang*np.dot(PJ0b[0, :], BJ0[iB]) if k_used[0]: gEM += np.dot(PJ0[0, :], BJ0[iB]) return gEM
def test_wavenumber(): # 1. wavenumber dat = DATAKERNEL['wave'][()] for key, val in dat.items(): out = kernel.wavenumber(ab=val[0], msrc=val[1], mrec=val[2], **val[3]) assert_allclose(out[0], val[4][0]) assert_allclose(out[1], val[4][1]) assert_allclose(out[2], val[4][2])
def getkernel(i, inplambd, inpoff, inpfang): iB = i * dat['nquad'] + np.arange(dat['nquad']) PJ = kernel.wavenumber(lambd=np.atleast_2d(inplambd)[:, iB], **dat['nsinp']) fEM = inpfang * np.dot(PJ[1][0, :], dat['BJ1'][iB]) if dat['ab'] in [11, 12, 21, 22, 14, 24, 15, 25]: fEM /= np.atleast_1d(inpoff) fEM += inpfang * np.dot(PJ[2][0, :], dat['BJ0'][iB]) fEM += np.dot(PJ[0][0, :], dat['BJ0'][iB]) return fEM
def test_wavenumber(): # 1. wavenumber dat = DATAKERNEL['wave'][()] for _, val in dat.items(): out = kernel.wavenumber(ab=val[0], msrc=val[1], mrec=val[2], **val[3]) if val[0] in [11, 22, 24, 15, 33]: assert_allclose(out[0], val[4][0]) else: assert out[0] is None if val[0] == 33: assert out[1] is None else: assert_allclose(out[1], val[4][1]) if val[0] in [11, 22, 24, 15, 12, 21, 14, 25]: assert_allclose(out[2], val[4][2]) else: assert out[2] is None
break xint = np.concatenate((np.array([1e-20]), b_zero)) dx = np.repeat(np.diff(xint) / 2, nquad) Bx = dx * (np.tile(g_x, maxint) + 1) + np.repeat(xint[:-1], nquad) BJ0 = special.j0(Bx) * np.tile(g_w, maxint) BJ1 = special.j1(Bx) * np.tile(g_w, maxint) intervals = xint / off[:, None] lambd = Bx / off[:, None] factAng = kernel.angle_factor(angle, ab, msrc, mrec) # 1 Spline version start = np.log(lambd.min()) stop = np.log(lambd.max()) ilambd = np.logspace(start, stop, (stop - start) * pts_per_dec + 1, 10) PJ0, PJ1, PJ0b = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH[None, :], etaV[None, :], zetaH[None, :], zetaV[None, :], np.atleast_2d(ilambd), ab, False, msrc, mrec, False) si_PJ0r = iuSpline(np.log(ilambd), PJ0.real) si_PJ0i = iuSpline(np.log(ilambd), PJ0.imag) si_PJ1r = iuSpline(np.log(ilambd), PJ1.real) si_PJ1i = iuSpline(np.log(ilambd), PJ1.imag) si_PJ0br = iuSpline(np.log(ilambd), PJ0b.real) si_PJ0bi = iuSpline(np.log(ilambd), PJ0b.imag) sPJ0 = si_PJ0r(np.log(lambd)) + 1j * si_PJ0i(np.log(lambd)) sPJ1 = si_PJ1r(np.log(lambd)) + 1j * si_PJ1i(np.log(lambd)) sPJ0b = si_PJ0br(np.log(lambd)) + 1j * si_PJ0bi(np.log(lambd)) sEM = np.sum(np.reshape(sPJ1 * BJ1, (off.size, nquad, -1), order='F'), 1) if ab in [11, 12, 21, 22, 14, 24, 15, 25]: # Because of J2 # J2(kr) = 2/(kr)*J1(kr) - J0(kr) sEM /= np.atleast_1d(off[:, np.newaxis])
def test_dlf(): # 10. dlf # DLF is integral of hankel_dlf and fourier_dlf, and therefore tested a lot # through those. Here we just ensure status quo. And if a problem arises in # hankel_dlf or fourier_dlf, it would make it obvious if the problem arises # from dlf or not. # Check DLF for Fourier t = DATA['t'][()] for i in [0, 1, 2]: dat = DATA['fourier_dlf' + str(i)][()] tres = DATA['tEM' + str(i)][()] finp = dat['fEM'] ftarg = dat['ftarg'] if i > 0: finp /= 2j * np.pi * dat['f'] if i > 1: finp *= -1 if ftarg['pts_per_dec'] == 0: finp = finp.reshape(t.size, -1) tEM = transform.dlf(finp, 2 * np.pi * dat['f'], t, ftarg['dlf'], ftarg['pts_per_dec'], kind=ftarg['kind']) assert_allclose(tEM * 2 / np.pi, tres, rtol=1e-3) # Check DLF for Hankel for ab in [12, 22, 13, 33]: model = utils.check_model([], 10, 2, 2, 5, 1, 10, True, 0) depth, res, aniso, epermH, epermV, mpermH, mpermV, _ = model frequency = utils.check_frequency(1, res, aniso, epermH, epermV, mpermH, mpermV, 0) _, etaH, etaV, zetaH, zetaV = frequency src = [0, 0, 0] src, nsrc = utils.check_dipole(src, 'src', 0) ab, msrc, mrec = utils.check_ab(ab, 0) ht, htarg = utils.check_hankel('dlf', {}, 0) xdirect = False # Important, as we want to comp. wavenumber-frequency! rec = [np.arange(1, 11) * 500, np.zeros(10), 300] rec, nrec = utils.check_dipole(rec, 'rec', 0) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, 0) lsrc, zsrc = utils.get_layer_nr(src, depth) lrec, zrec = utils.get_layer_nr(rec, depth) dlf = htarg['dlf'] pts_per_dec = htarg['pts_per_dec'] # # # 0. No Spline # # # # dlf calculation lambd = dlf.base / off[:, None] PJ = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) # Angle factor, one example with None instead of 1's. if ab != 13: ang_fact = kernel.angle_factor(angle, ab, msrc, mrec) else: ang_fact = None # dlf calculation fEM0 = transform.dlf(PJ, lambd, off, dlf, 0, ang_fact=ang_fact, ab=ab) # Analytical frequency-domain solution freq1 = kernel.fullspace(off, angle, zsrc, zrec, etaH, etaV, zetaH, zetaV, ab, msrc, mrec) # Compare assert_allclose(np.squeeze(fEM0), np.squeeze(freq1)) # # # 1. Spline; One angle # # # # dlf calculation lambd, _ = transform.get_dlf_points(dlf, off, pts_per_dec) PJ1 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) # dlf calculation fEM1 = transform.dlf(PJ1, lambd, off, dlf, pts_per_dec, ang_fact=ang_fact, ab=ab) # Compare assert_allclose(np.squeeze(fEM1), np.squeeze(freq1), rtol=1e-4) # # # 2.a Lagged; One angle # # # rec = [np.arange(1, 11) * 500, np.arange(-5, 5) * 0, 300] rec, nrec = utils.check_dipole(rec, 'rec', 0) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, 0) # dlf calculation lambd, _ = transform.get_dlf_points(dlf, off, -1) PJ2 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) ang_fact = kernel.angle_factor(angle, ab, msrc, mrec) # dlf calculation fEM2 = transform.dlf(PJ2, lambd, off, dlf, -1, ang_fact=ang_fact, ab=ab) # Analytical frequency-domain solution freq2 = kernel.fullspace(off, angle, zsrc, zrec, etaH, etaV, zetaH, zetaV, ab, msrc, mrec) # Compare assert_allclose(np.squeeze(fEM2), np.squeeze(freq2), rtol=1e-4) # # # 2.b Lagged; Multi angle # # # rec = [np.arange(1, 11) * 500, np.arange(-5, 5) * 200, 300] rec, nrec = utils.check_dipole(rec, 'rec', 0) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, 0) # dlf calculation lambd, _ = transform.get_dlf_points(dlf, off, -1) PJ2 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) ang_fact = kernel.angle_factor(angle, ab, msrc, mrec) # dlf calculation fEM2 = transform.dlf(PJ2, lambd, off, dlf, -1, ang_fact=ang_fact, ab=ab) # Analytical frequency-domain solution freq2 = kernel.fullspace(off, angle, zsrc, zrec, etaH, etaV, zetaH, zetaV, ab, msrc, mrec) # Compare assert_allclose(np.squeeze(fEM2), np.squeeze(freq2), rtol=1e-4) # # # 3. Spline; Multi angle # # # lambd, _ = transform.get_dlf_points(dlf, off, 30) # dlf calculation PJ3 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) # dlf calculation fEM3 = transform.dlf(PJ3, lambd, off, dlf, 30, ang_fact=ang_fact, ab=ab) # Compare assert_allclose(np.squeeze(fEM3), np.squeeze(freq2), rtol=1e-3)
def test_dlf(): # 10. dlf # DLF is integral of fht and ffht, and therefore tested a lot through # those. Here we just ensure status quo. And if a problem arises in fht or # ffht, it would make it obvious if the problem arises from dlf or not. # Check DLF for Fourier t = DATA['t'][()] for i in [0, 1, 2]: dat = DATA['ffht' + str(i)][()] tres = DATA['tEM' + str(i)][()] finp = dat['fEM'] ftarg = dat['ftarg'] if i > 0: finp /= 2j * np.pi * dat['f'] if i > 1: finp *= -1 if ftarg[1] == 0: finp = finp.reshape(t.size, -1) tEM = transform.dlf(finp, 2 * np.pi * dat['f'], t, ftarg[0], ftarg[1], kind=ftarg[2]) assert_allclose(tEM * 2 / np.pi, tres, rtol=1e-3) # Check DLF for Hankel for ab in [12, 22, 13, 33]: model = utils.check_model([], 10, 2, 2, 5, 1, 10, True, 0) depth, res, aniso, epermH, epermV, mpermH, mpermV, isfullspace = model frequency = utils.check_frequency(1, res, aniso, epermH, epermV, mpermH, mpermV, 0) freq, etaH, etaV, zetaH, zetaV = frequency src = [0, 0, 0] src, nsrc = utils.check_dipole(src, 'src', 0) ab, msrc, mrec = utils.check_ab(ab, 0) ht, htarg = utils.check_hankel('fht', None, 0) options = utils.check_opt(None, None, ht, htarg, 0) use_ne_eval, loop_freq, loop_off = options xdirect = False # Important, as we want to comp. wavenumber-frequency! rec = [np.arange(1, 11) * 500, np.zeros(10), 300] rec, nrec = utils.check_dipole(rec, 'rec', 0) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, 0) lsrc, zsrc = utils.get_layer_nr(src, depth) lrec, zrec = utils.get_layer_nr(rec, depth) fhtfilt = htarg[0] pts_per_dec = htarg[1] # # # 0. No Spline # # # # fht calculation lambd = fhtfilt.base / off[:, None] PJ = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec, use_ne_eval) factAng = kernel.angle_factor(angle, ab, msrc, mrec) # dlf calculation fEM0 = transform.dlf(PJ, lambd, off, fhtfilt, 0, factAng=factAng, ab=ab) # Analytical frequency-domain solution freq1 = kernel.fullspace(off, angle, zsrc, zrec, etaH, etaV, zetaH, zetaV, ab, msrc, mrec) # Compare assert_allclose(np.squeeze(fEM0), np.squeeze(freq1)) # # # 1. Spline; One angle # # # options = utils.check_opt('spline', None, ht, htarg, 0) use_ne_eval, loop_freq, loop_off = options # fht calculation lambd, _ = transform.get_spline_values(fhtfilt, off, pts_per_dec) PJ1 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec, use_ne_eval) # dlf calculation fEM1 = transform.dlf(PJ1, lambd, off, fhtfilt, pts_per_dec, factAng=factAng, ab=ab) # Compare assert_allclose(np.squeeze(fEM1), np.squeeze(freq1), rtol=1e-4) # # # 2.a Lagged; One angle # # # rec = [np.arange(1, 11) * 500, np.arange(-5, 5) * 0, 300] rec, nrec = utils.check_dipole(rec, 'rec', 0) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, 0) # fht calculation lambd, _ = transform.get_spline_values(fhtfilt, off, -1) PJ2 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec, use_ne_eval) factAng = kernel.angle_factor(angle, ab, msrc, mrec) # dlf calculation fEM2 = transform.dlf(PJ2, lambd, off, fhtfilt, -1, factAng=factAng, ab=ab) # Analytical frequency-domain solution freq2 = kernel.fullspace(off, angle, zsrc, zrec, etaH, etaV, zetaH, zetaV, ab, msrc, mrec) # Compare assert_allclose(np.squeeze(fEM2), np.squeeze(freq2), rtol=1e-4) # # # 2.b Lagged; Multi angle # # # rec = [np.arange(1, 11) * 500, np.arange(-5, 5) * 200, 300] rec, nrec = utils.check_dipole(rec, 'rec', 0) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, 0) # fht calculation lambd, _ = transform.get_spline_values(fhtfilt, off, -1) PJ2 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec, use_ne_eval) factAng = kernel.angle_factor(angle, ab, msrc, mrec) # dlf calculation fEM2 = transform.dlf(PJ2, lambd, off, fhtfilt, -1, factAng=factAng, ab=ab) # Analytical frequency-domain solution freq2 = kernel.fullspace(off, angle, zsrc, zrec, etaH, etaV, zetaH, zetaV, ab, msrc, mrec) # Compare assert_allclose(np.squeeze(fEM2), np.squeeze(freq2), rtol=1e-4) # # # 3. Spline; Multi angle # # # lambd, _ = transform.get_spline_values(fhtfilt, off, 10) # fht calculation PJ3 = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec, use_ne_eval) # dlf calculation fEM3 = transform.dlf(PJ3, lambd, off, fhtfilt, 10, factAng=factAng, ab=ab) # Compare assert_allclose(np.squeeze(fEM3), np.squeeze(freq2), rtol=1e-3)
lambd = filters.key_51_2012().base/np.array([0.001, 1, 100, 10000])[:, None] depth = np.array([-np.infty, 0, 150, 300, 500, 600]) inp1 = {'zsrc': 100., 'zrec': 650., 'lsrc': np.array(1), 'lrec': np.array(5), 'depth': depth, 'etaH': etaH, 'etaV': etaV, 'zetaH': zetaH, 'zetaV': zetaV, 'lambd': lambd, 'xdirect': False} wave = {} for key, val in iab.items(): res = kernel.wavenumber(ab=val[2], msrc=val[0], mrec=val[1], **inp1) wave[key] = (val[2], val[0], val[1], inp1, res) # # C -- GREENFCT # # # Standard example inp2 = deepcopy(inp1) # Source and receiver in same layer (last) inp3 = deepcopy(inp1) inp3['zsrc'] = 610. inp3['lsrc'] = np.array(5) # Receiver in first layer inp4 = deepcopy(inp1) inp4['zrec'] = -30. inp4['lrec'] = np.array(0) green = {}
def hankel_dlf(zsrc, zrec, lsrc, lrec, off, ang_fact, depth, ab, etaH, etaV, zetaH, zetaV, xdirect, htarg, msrc, mrec): r"""Hankel Transform using the Digital Linear Filter method. The *Digital Linear Filter* method was introduced to geophysics by [Ghos70]_, and made popular and wide-spread by [Ande75]_, [Ande79]_, [Ande82]_. The DLF is sometimes referred to as the *Fast Hankel Transform* FHT, from which this routine has its name. This implementation of the DLF follows [Key12]_, equation 6. Without going into the mathematical details (which can be found in any of the above papers) and following [Key12]_, the DLF method rewrites the Hankel transform of the form .. math:: :label: dlf1 F(r) = \int^\infty_0 f(\lambda)J_v(\lambda r)\ \mathrm{d}\lambda as .. math:: :label: dlf2 F(r) = \sum^n_{i=1} f(b_i/r)h_i/r \ , where :math:`h` is the digital filter.The Filter abscissae b is given by .. math:: :label: dlf3 b_i = \lambda_ir = e^{ai}, \qquad i = -l, -l+1, \cdots, l \ , with :math:`l=(n-1)/2`, and :math:`a` is the spacing coefficient. This function is loosely based on `get_CSEM1D_FD_FHT.m` from the source code distributed with [Key12]_. The function is called from one of the modelling routines in :mod:`empymod.model`. Consult these modelling routines for a description of the input and output parameters. Returns ------- fEM : array Returns frequency-domain EM response. kcount : int Kernel count. For DLF, this is 1. conv : bool Only relevant for QWE/QUAD. """ # Compute required lambdas for given Hankel-filter-base lambd, int_pts = get_dlf_points(htarg['dlf'], off, htarg['pts_per_dec']) # Call the kernel PJ = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) # Carry out the dlf fEM = dlf(PJ, lambd, off, htarg['dlf'], htarg['pts_per_dec'], ang_fact=ang_fact, ab=ab, int_pts=int_pts) return fEM, 1, True
def hankel_quad(zsrc, zrec, lsrc, lrec, off, ang_fact, depth, ab, etaH, etaV, zetaH, zetaV, xdirect, htarg, msrc, mrec): r"""Hankel Transform using the `QUADPACK` library. This routine uses the :func:`scipy.integrate.quad` module, which in turn makes use of the Fortran library `QUADPACK` (`qagse`). It is massively (orders of magnitudes) slower than either :func:`hankel_dlf` or :func:`hankel_qwe`, and is mainly here for completeness and comparison purposes. It always uses interpolation in the wavenumber domain, hence it generally will not be as precise as the other methods. However, it might work in some areas where the others fail. The function is called from one of the modelling routines in :mod:`empymod.model`. Consult these modelling routines for a description of the input and output parameters. Returns ------- fEM : array Returns frequency-domain EM response. kcount : int Kernel count. For HQUAD, this is 1. conv : bool If true, QUAD converged. If not, `htarg` might have to be adjusted. """ # Get required lambdas la = np.log10(htarg['a']) lb = np.log10(htarg['b']) ilambd = np.logspace(la, lb, int((lb-la)*htarg['pts_per_dec'] + 1)) # Call the kernel PJ0, PJ1, PJ0b = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, np.atleast_2d(ilambd), ab, xdirect, msrc, mrec) # Interpolation in wavenumber domain: Has to be done separately on each PJ, # in order to work with multiple offsets which have different angles. # We check if the kernels are zero, to avoid unnecessary calculations. if PJ0 is not None: sPJ0r = iuSpline(np.log(ilambd), PJ0.real) sPJ0i = iuSpline(np.log(ilambd), PJ0.imag) else: sPJ0r = None sPJ0i = None if PJ1 is not None: sPJ1r = iuSpline(np.log(ilambd), PJ1.real) sPJ1i = iuSpline(np.log(ilambd), PJ1.imag) else: sPJ1r = None sPJ1i = None if PJ0b is not None: sPJ0br = iuSpline(np.log(ilambd), PJ0b.real) sPJ0bi = iuSpline(np.log(ilambd), PJ0b.imag) else: sPJ0br = None sPJ0bi = None # Pre-allocate output array fEM = np.zeros(off.size, dtype=np.complex128) conv = True # Input-dictionary for quad iinp = {'a': htarg['a'], 'b': htarg['b'], 'epsabs': htarg['atol'], 'epsrel': htarg['rtol'], 'limit': htarg['limit']} # Loop over offsets for i in range(off.size): fEM[i], tc = quad(sPJ0r, sPJ0i, sPJ1r, sPJ1i, sPJ0br, sPJ0bi, ab, off[i], ang_fact[i], iinp) conv *= tc # Return the electromagnetic field # Second argument (1) is the kernel count, last argument is only for QWE. return fEM, 1, conv
def hankel_qwe(zsrc, zrec, lsrc, lrec, off, ang_fact, depth, ab, etaH, etaV, zetaH, zetaV, xdirect, htarg, msrc, mrec): r"""Hankel Transform using Quadrature-With-Extrapolation. *Quadrature-With-Extrapolation* was introduced to geophysics by [Key12]_. It is one of many so-called *ISE* methods to solve Hankel Transforms, where *ISE* stands for Integration, Summation, and Extrapolation. Following [Key12]_, but without going into the mathematical details here, the QWE method rewrites the Hankel transform of the form .. math:: :label: qwe1 F(r) = \int^\infty_0 f(\lambda)J_v(\lambda r)\ \mathrm{d}\lambda as a quadrature sum which form is similar to the DLF (equation 15), .. math:: :label: qwe2 F_i \approx \sum^m_{j=1} f(x_j/r)w_j g(x_j) = \sum^m_{j=1} f(x_j/r)\hat{g}(x_j) \ , but with various bells and whistles applied (using the so-called Shanks transformation in the form of a routine called :math:`\epsilon`-algorithm ([Shan55]_, [Wynn56]_; implemented with algorithms from [Tref00]_ and [Weni89]_). This function is based on `get_CSEM1D_FD_QWE.m`, `qwe.m`, and `getBesselWeights.m` from the source code distributed with [Key12]_. In the spline-version, :func:`hankel_qwe` checks how steep the decay of the wavenumber-domain result is, and calls QUAD for the very steep interval, for which QWE is not suited. The function is called from one of the modelling routines in :mod:`empymod.model`. Consult these modelling routines for a description of the input and output parameters. Returns ------- fEM : array Returns frequency-domain EM response. kcount : int Kernel count. conv : bool If true, QWE/QUAD converged. If not, `htarg` might have to be adjusted. """ # Input params have an additional dimension for frequency, reduce here etaH = etaH[0, :] etaV = etaV[0, :] zetaH = zetaH[0, :] zetaV = zetaV[0, :] # Get rtol, atol, nquad, maxint, and pts_per_dec rtol = htarg['rtol'] atol = htarg['atol'] nquad = htarg['nquad'] maxint = htarg['maxint'] pts_per_dec = htarg['pts_per_dec'] # 1. PRE-COMPUTE THE BESSEL FUNCTIONS # at fixed quadrature points for each interval and multiply by the # corresponding Gauss quadrature weights # Get Gauss quadrature weights g_x, g_w = special.roots_legendre(nquad) # Compute n zeros of the Bessel function of the first kind of order 1 using # the Newton-Raphson method, which is fast enough for our purposes. Could # be done with a loop for (but it is slower): # b_zero[i] = optimize.newton(special.j1, b_zero[i]) # Initial guess using asymptotic zeros b_zero = np.pi*np.arange(1.25, maxint+1) # Newton-Raphson iterations for i in range(10): # 10 is more than enough, usually stops in 5 # Evaluate b_x0 = special.j1(b_zero) # j0 and j1 have faster versions b_x1 = special.jv(2, b_zero) # j2 does not have a faster version # The step length b_h = -b_x0/(b_x0/b_zero - b_x1) # Take the step b_zero += b_h # Check for convergence if all(np.abs(b_h) < 8*np.finfo(float).eps*b_zero): break # 2. COMPUTE THE QUADRATURE INTERVALS AND BESSEL FUNCTION WEIGHTS # Lower limit of integrand, a small but non-zero value xint = np.concatenate((np.array([1e-20]), b_zero)) # Assemble the output arrays dx = np.repeat(np.diff(xint)/2, nquad) Bx = dx*(np.tile(g_x, maxint) + 1) + np.repeat(xint[:-1], nquad) BJ0 = special.j0(Bx)*np.tile(g_w, maxint) BJ1 = special.j1(Bx)*np.tile(g_w, maxint) # 3. START QWE # Intervals and lambdas for all offset intervals = xint/off[:, None] lambd = Bx/off[:, None] # The following lines until # "Call and return QWE, depending if spline or not" # are part of the splined routine. However, we calculate it here to get # the non-zero kernels, `k_used`. # New lambda, from min to max required lambda with pts_per_dec start = np.log10(lambd.min()) stop = np.log10(lambd.max()) # If not spline, we just calculate three lambdas to check if pts_per_dec == 0: ilambd = np.logspace(start, stop, 3) else: ilambd = np.logspace(start, stop, int((stop-start)*pts_per_dec + 1)) # Call the kernel PJ0, PJ1, PJ0b = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH[None, :], etaV[None, :], zetaH[None, :], zetaV[None, :], np.atleast_2d(ilambd), ab, xdirect, msrc, mrec) # Check which kernels have information k_used = [True, True, True] for i, val in enumerate((PJ0, PJ1, PJ0b)): if val is None: k_used[i] = False # Call and return QWE, depending if spline or not if pts_per_dec != 0: # If spline, we calculate all kernels here # Interpolation : Has to be done separately on each PJ, # in order to work with multiple offsets which have different angles. if k_used[0]: sPJ0r = iuSpline(np.log(ilambd), PJ0.real) sPJ0i = iuSpline(np.log(ilambd), PJ0.imag) else: sPJ0r = None sPJ0i = None if k_used[1]: sPJ1r = iuSpline(np.log(ilambd), PJ1.real) sPJ1i = iuSpline(np.log(ilambd), PJ1.imag) else: sPJ1r = None sPJ1i = None if k_used[2]: sPJ0br = iuSpline(np.log(ilambd), PJ0b.real) sPJ0bi = iuSpline(np.log(ilambd), PJ0b.imag) else: sPJ0br = None sPJ0bi = None # Get htarg: diff_quad, a, b, limit diff_quad = htarg['diff_quad'] a = htarg['a'] b = htarg['b'] limit = htarg['limit'] # Set htarg if not given: if not limit: limit = maxint if not a: a = intervals[:, 0] else: a = a*np.ones(off.shape) if not b: b = intervals[:, -1] else: b = b*np.ones(off.shape) # Check if we use QWE or SciPy's QUAD # If there are any steep decays within an interval we have to use QUAD, # as QWE is not designed for these intervals. check0 = np.log(intervals[:, :-1]) check1 = np.log(intervals[:, 1:]) numerator = np.zeros((off.size, maxint), dtype=np.complex128) denominator = np.zeros((off.size, maxint), dtype=np.complex128) if k_used[0]: numerator += sPJ0r(check0) + 1j*sPJ0i(check0) denominator += sPJ0r(check1) + 1j*sPJ0i(check1) if k_used[1]: numerator += sPJ1r(check0) + 1j*sPJ1i(check0) denominator += sPJ1r(check1) + 1j*sPJ1i(check1) if k_used[2]: numerator += sPJ0br(check0) + 1j*sPJ0bi(check0) denominator += sPJ0br(check1) + 1j*sPJ0bi(check1) doqwe = np.all((np.abs(numerator)/np.abs(denominator) < diff_quad), 1) # Pre-allocate output array fEM = np.zeros(off.size, dtype=np.complex128) conv = True # Carry out SciPy's Quad if required if np.any(~doqwe): # Loop over offsets that require Quad for i in np.where(~doqwe)[0]: # Input-dictionary for quad iinp = {'a': a[i], 'b': b[i], 'epsabs': atol, 'epsrel': rtol, 'limit': limit} fEM[i], tc = quad(sPJ0r, sPJ0i, sPJ1r, sPJ1i, sPJ0br, sPJ0bi, ab, off[i], ang_fact[i], iinp) # Update conv conv *= tc # Return kcount=1 in case no QWE is calculated kcount = 1 if np.any(doqwe): # Get EM-field at required offsets if k_used[0]: sPJ0 = sPJ0r(np.log(lambd)) + 1j*sPJ0i(np.log(lambd)) if k_used[1]: sPJ1 = sPJ1r(np.log(lambd)) + 1j*sPJ1i(np.log(lambd)) if k_used[2]: sPJ0b = sPJ0br(np.log(lambd)) + 1j*sPJ0bi(np.log(lambd)) # Carry out and return the Hankel transform for this interval sEM = np.zeros_like(numerator, dtype=np.complex128) if k_used[1]: sEM += np.sum(np.reshape(sPJ1*BJ1, (off.size, nquad, -1), order='F'), 1) if ab in [11, 12, 21, 22, 14, 24, 15, 25]: # Because of J2 # J2(kr) = 2/(kr)*J1(kr) - J0(kr) sEM /= np.atleast_1d(off[:, np.newaxis]) if k_used[2]: sEM += np.sum(np.reshape(sPJ0b*BJ0, (off.size, nquad, -1), order='F'), 1) if k_used[1] or k_used[2]: sEM *= ang_fact[:, np.newaxis] if k_used[0]: sEM += np.sum(np.reshape(sPJ0*BJ0, (off.size, nquad, -1), order='F'), 1) getkernel = sEM[doqwe, :] # Get QWE fEM[doqwe], kcount, tc = qwe(rtol, atol, maxint, getkernel, intervals[doqwe, :], None, None, None) conv *= tc else: # If not spline, we define the wavenumber-kernel here def getkernel(i, inplambd, inpoff, inpfang): r"""Return wavenumber-domain-kernel as a fct of interval i.""" # Indices and factor for this interval iB = i*nquad + np.arange(nquad) # PJ0 and PJ1 for this interval PJ0, PJ1, PJ0b = kernel.wavenumber(zsrc, zrec, lsrc, lrec, depth, etaH[None, :], etaV[None, :], zetaH[None, :], zetaV[None, :], np.atleast_2d(inplambd)[:, iB], ab, xdirect, msrc, mrec) # Carry out and return the Hankel transform for this interval gEM = np.zeros_like(inpoff, dtype=np.complex128) if k_used[1]: gEM += inpfang*np.dot(PJ1[0, :], BJ1[iB]) if ab in [11, 12, 21, 22, 14, 24, 15, 25]: # Because of J2 # J2(kr) = 2/(kr)*J1(kr) - J0(kr) gEM /= np.atleast_1d(inpoff) if k_used[2]: gEM += inpfang*np.dot(PJ0b[0, :], BJ0[iB]) if k_used[0]: gEM += np.dot(PJ0[0, :], BJ0[iB]) return gEM # Get QWE fEM, kcount, conv = qwe(rtol, atol, maxint, getkernel, intervals, lambd, off, ang_fact) return fEM, kcount, conv
def setup_cache(self): """setup_cache is not parametrized, so we do it manually. """ data = {} for size in self.params[0]: # size data[size] = {} # One big, one small model if size == 'Small': # Small; Total size: 5*1*1*1 = 5 x = np.array([500., 1000.]) else: # Big; Total size: 5*100*100*201 = 10'050'000 x = np.arange(1, 101)*200. # Define model parameters freq = np.array([1]) src = [0, 0, 250] rec = [x, np.zeros(x.shape), 300] depth = np.array([-np.infty, 0, 300, 2000, 2100]) res = np.array([2e14, .3, 1, 50, 1]) ab = 11 xdirect = False verb = 0 if not VERSION2: use_ne_eval = False # Checks (since DLF exists the `utils`-checks haven't changed, so # we just use them here. model = utils.check_model(depth, res, None, None, None, None, None, xdirect, verb) depth, res, aniso, epermH, epermV, mpermH, mpermV, _ = model frequency = utils.check_frequency(freq, res, aniso, epermH, epermV, mpermH, mpermV, verb) freq, etaH, etaV, zetaH, zetaV = frequency ab, msrc, mrec = utils.check_ab(ab, verb) src, nsrc = utils.check_dipole(src, 'src', verb) rec, nrec = utils.check_dipole(rec, 'rec', verb) off, angle = utils.get_off_ang(src, rec, nsrc, nrec, verb) lsrc, zsrc = utils.get_layer_nr(src, depth) lrec, zrec = utils.get_layer_nr(rec, depth) for htype in self.params[1]: # htype # pts_per_dec depending on htype if htype == 'Standard': pts_per_dec = 0 elif htype == 'Lagged': pts_per_dec = -1 else: pts_per_dec = 10 # Compute kernels for dlf if VERSION2: # HT arguments _, fhtarg = utils.check_hankel( 'dlf', {'dlf': 'key_201_2009', 'pts_per_dec': pts_per_dec}, 0) inp = (fhtarg['dlf'], off, fhtarg['pts_per_dec']) lambd, _ = transform.get_dlf_points(*inp) else: # HT arguments _, fhtarg = utils.check_hankel( 'fht', ['key_201_2009', pts_per_dec], 0) inp = (fhtarg[0], off, fhtarg[1]) lambd, _ = transform.get_spline_values(*inp) if VERSION2: inp = (zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec) else: inp = (zsrc, zrec, lsrc, lrec, depth, etaH, etaV, zetaH, zetaV, lambd, ab, xdirect, msrc, mrec, use_ne_eval) PJ = kernel.wavenumber(*inp) factAng = kernel.angle_factor(angle, ab, msrc, mrec) # Signature changed at commit a15af07 (20/05/2018; before # v1.6.2) try: dlf = {'signal': PJ, 'points': lambd, 'out_pts': off, 'ab': ab} if VERSION2: dlf['ang_fact'] = factAng dlf['filt'] = fhtarg['dlf'] dlf['pts_per_dec'] = fhtarg['pts_per_dec'] else: dlf['factAng'] = factAng dlf['filt'] = fhtarg[0] dlf['pts_per_dec'] = fhtarg[1] transform.dlf(**dlf) except VariableCatch: dlf = {'signal': PJ, 'points': lambd, 'out_pts': off, 'targ': fhtarg, 'factAng': factAng} data[size][htype] = dlf return data