def test_im_to_vis_phase_centre():
    """
    The simplest test here is to see if a single source at the phase centre
    returns simply the flux everywhere with zero imaginary part
    """
    from africanus.dft.kernels import im_to_vis

    nrow = 100
    uvw = np.random.random(size=(nrow, 3))
    npix = 35  # must be odd for this test to work
    x = np.linspace(-0.1, 0.1, npix)
    ll, mm = np.meshgrid(x, x)
    lm = np.vstack((ll.flatten(), mm.flatten())).T
    nchan = 11
    frequency = np.linspace(1.0, 2.0, nchan, endpoint=True)

    image = np.zeros([npix, npix, nchan], dtype=np.float64)
    I0 = 1.0
    ref_freq = frequency[nchan//2]
    Inu = I0*(frequency/ref_freq)**(-0.7)
    image[npix//2, npix//2, :] = Inu
    image = image.reshape(npix**2, nchan)

    vis = im_to_vis(image, uvw, lm, frequency)

    for i in range(nchan):
        tmp = vis[:, i] - Inu[i]
        assert np.all(tmp.real < 1e-13)
        assert np.all(tmp.imag < 1e-13)
def test_im_to_vis_zero_w():
    """
    This test checks that the result matches the analytic result in the case when 
    w = 0 for multiple channels and sources.
    """
    from africanus.dft.kernels import im_to_vis
    from africanus.constants import minus_two_pi_over_c
    np.random.seed(123)
    nrow = 100
    uvw = np.random.random(size=(nrow, 3))
    uvw[:, 2] = 0.0
    nchan = 3
    frequency = np.linspace(1.e9, 2.e9, nchan, endpoint=True)
    nsource = 5
    I0 = np.random.randn(nsource)
    ref_freq = frequency[nchan//2]
    image = I0[:, None] * (frequency/ref_freq)**(-0.7)
    l = 0.001 + 0.1*np.random.random(nsource)
    m = 0.001 + 0.1*np.random.random(nsource)
    lm = np.vstack((l,m)).T
    vis = im_to_vis(image, uvw, lm, frequency)

    vis_true = np.zeros([nrow, nchan], dtype=np.complex128)


    for ch in range(nchan):
        for source in range(nsource):
            vis_true[:, ch] += image[source, ch]*np.exp(
                                 minus_two_pi_over_c*frequency[ch] * 1.0j *
                                 (uvw[:,0]*lm[source,0]+uvw[:,1]*lm[source,1]))

    assert np.allclose(vis, vis_true)
Beispiel #3
0
def test_im_to_vis_simple():
    """
    This test checks that the result matches the analytic result
    w = 0 for multiple channels and sources but a single correlation
    """
    from africanus.dft.kernels import im_to_vis
    from africanus.constants import minus_two_pi_over_c
    np.random.seed(123)
    nrow = 100
    uvw = np.random.random(size=(nrow, 3))
    nchan = 3
    frequency = np.linspace(1.e9, 2.e9, nchan, endpoint=True)
    nsource = 5
    I0 = np.random.randn(nsource)
    ref_freq = frequency[nchan // 2]
    image = I0[:, None] * (frequency / ref_freq)**(-0.7)
    # add correlation axis
    image = image[:, :, None]
    l = 0.001 + 0.1 * np.random.random(nsource)  # noqa
    m = 0.001 + 0.1 * np.random.random(nsource)
    lm = np.vstack((l, m)).T
    vis = im_to_vis(image, uvw, lm, frequency).squeeze()

    vis_true = np.zeros([nrow, nchan], dtype=np.complex128)

    for ch in range(nchan):
        for source in range(nsource):
            l, m = lm[source]
            n = np.sqrt(1.0 - l**2 - m**2)
            phase = minus_two_pi_over_c*frequency[ch] * 1.0j * \
                (uvw[:, 0]*l + uvw[:, 1]*m + uvw[:, 2]*(n-1))

            vis_true[:, ch] += image[source, ch, 0] * np.exp(phase)
    assert_array_almost_equal(vis, vis_true, decimal=14)
Beispiel #4
0
def test_symmetric_covariance():
    """
    Test that the image plane precision matrix R^H Sigma^-1R is Hermitian
    (symmetric since its real).
    """
    from africanus.dft.kernels import vis_to_im, im_to_vis
    np.random.seed(123)

    nsource = 25
    ncorr = 1
    nchan = 1

    lmmax = 0.05
    ll = -lmmax + 2 * lmmax * np.random.random(nsource)  # noqa
    mm = -lmmax + 2 * lmmax * np.random.random(nsource)
    lm = np.vstack((ll, mm)).T

    nrows = 1000
    uvw = np.random.randn(nrows, 3) * 1000
    uvw[:, 2] = 0.0

    freq = np.array([1.0e9])

    flags = np.zeros((nrows, nchan, ncorr), dtype=np.bool)

    # get the "psf" matrix at source locations
    psf_source = np.zeros((nsource, nsource), dtype=np.float64)
    point_source = np.ones((1, nchan, ncorr), dtype=np.float64)
    for source in range(nsource):
        lm_i = lm[source].reshape(1, 2)
        Ki = im_to_vis(point_source, uvw, lm_i, freq)
        psf_source[:, source] = vis_to_im(Ki, uvw, lm, freq, flags).squeeze()

    assert_array_almost_equal(psf_source, psf_source.T, decimal=14)
Beispiel #5
0
def test_degrid_dft_packed_dask_dft_check():
    da = pytest.importorskip("dask.array")

    # construct kernel
    W = 5
    OS = 3
    kern = kernels.pack_kernel(kernels.kbsinc(W, oversample=OS),
                               W,
                               oversample=OS)
    nrow = 100
    nrow_chunk = nrow // 8
    uvw = np.column_stack(
        (5000.0 * np.cos(np.linspace(0, 2 * np.pi, nrow)),
         5000.0 * np.sin(np.linspace(0, 2 * np.pi, nrow)), np.zeros(nrow)))

    pxacrossbeam = 10
    nchan = 16
    frequency = np.linspace(1.0e9, 1.4e9, nchan)
    wavelength = lightspeed / frequency

    cell = np.rad2deg(
        wavelength[0] /
        (2 * max(np.max(np.abs(uvw[:, 0])), np.max(np.abs(uvw[:, 1]))) *
         pxacrossbeam))
    npix = 512
    mod = np.zeros((1, npix, npix), dtype=np.complex64)
    mod[0, npix // 2 - 5, npix // 2 - 5] = 1.0

    ftmod = np.fft.ifftshift(np.fft.fft2(np.fft.fftshift(
        mod[0, :, :]))).reshape((1, 1, npix, npix))
    chanmap = np.zeros(nchan, dtype=np.int64)
    dec, ra = np.meshgrid(
        np.arange(-npix // 2, npix // 2) * np.deg2rad(cell),
        np.arange(-npix // 2, npix // 2) * np.deg2rad(cell))
    radec = np.column_stack((ra.flatten(), dec.flatten()))
    vis_dft = im_to_vis(mod[0, :, :].reshape(1, 1, npix * npix).T.copy(), uvw,
                        radec, frequency)

    vis_degrid = dwrap.degridder(
        da.from_array(uvw, chunks=(nrow_chunk, 3)),
        da.from_array(ftmod, chunks=(1, 1, npix, npix)),
        da.from_array(wavelength, chunks=(nchan, )),
        da.from_array(chanmap, chunks=(nchan, )),
        cell * 3600.0,
        da.from_array(np.array([[0, np.pi / 4.0]]), chunks=(1, 2)),
        (0, np.pi / 4.0),
        kern,
        W,
        OS,
        "None",  # no faceting
        "None",  # no faceting
        "XXYY_FROM_I",
        "conv_1d_axisymmetric_packed_gather")

    vis_degrid = vis_degrid.compute()

    assert np.percentile(
        np.abs(vis_dft[:, 0, 0].real - vis_degrid[:, 0, 0].real), 99.0) < 0.05
    assert np.percentile(
        np.abs(vis_dft[:, 0, 0].imag - vis_degrid[:, 0, 0].imag), 99.0) < 0.05
Beispiel #6
0
def test_im_to_vis_fft(convention):
    """
    Test against the fft when uv on regular and w is zero.
    """
    from africanus.dft.kernels import im_to_vis
    np.random.seed(123)
    Fs = np.fft.fftshift
    iFs = np.fft.ifftshift

    # set image and take fft
    npix = 29
    ncorr = 1
    image = np.zeros((npix, npix, ncorr), dtype=np.float64)
    fft_image = np.zeros((npix, npix, ncorr), dtype=np.complex128)
    nsource = 25
    for corr in range(ncorr):
        Ix = np.random.randint(5, npix - 5, nsource)
        Iy = np.random.randint(5, npix - 5, nsource)
        image[Ix, Iy, corr] = np.random.randn(nsource)
        fft_image[:, :, corr] = Fs(np.fft.fft2(iFs(image[:, :, corr])))

    # image space coords
    deltal = 0.001
    # this assumes npix is odd
    l_coord = np.arange(-(npix // 2), npix // 2 + 1) * deltal
    ll, mm = np.meshgrid(l_coord, l_coord)
    lm = np.vstack((ll.flatten(), mm.flatten())).T
    # uv-space coords
    u = Fs(np.fft.fftfreq(npix, d=deltal))
    uu, vv = np.meshgrid(u, u)
    uvw = np.zeros((npix**2, 3), dtype=np.float64)
    uvw[:, 0] = uu.flatten()
    uvw[:, 1] = vv.flatten()
    nchan = 1
    image = image.reshape(npix**2, nchan, ncorr)
    frequency = np.ones(nchan, dtype=np.float64)
    from africanus.constants import c as lightspeed
    frequency *= lightspeed  # makes result independent of frequency

    # take DFT and compare
    vis = im_to_vis(image, uvw, lm, frequency, convention=convention)
    fft_image = fft_image.reshape(npix**2, nchan, ncorr)
    fft_image = np.conj(fft_image) if convention == 'casa' else fft_image

    assert_array_almost_equal(vis, fft_image, decimal=13)
def test_im_to_vis_single_baseline_and_chan():
    """
    Here we check that the result is consistent for a single baseline, source and 
    channel. 
    """
    from africanus.dft.kernels import im_to_vis
    from africanus.constants import minus_two_pi_over_c
    nrow = 1
    uvw = np.random.random(size=(nrow, 3))
    frequency = np.array([1.5e9])
    l = 0.015
    m = -0.0123
    lm = np.array([[l, m]])
    n = np.sqrt(1 - l**2 - m**2)
    image = np.array([[1.0]])
    vis = im_to_vis(image, uvw, lm, frequency)

    vis_true = image*np.exp(minus_two_pi_over_c * frequency * 1.0j *
                            (uvw[:,0]*l + uvw[:,1]*m + uvw[:,2]*(n - 1.0)))

    assert np.allclose(vis, vis_true)
Beispiel #8
0
def gen_MODEL_DATA(sel_FREQ, sel_UVW, sel_LSM, 
                    gen_lm, n_chan, n_row, n_dir):

    # Create initial model array
    model = np.zeros((n_dir, n_chan, 4), dtype=np.float64)

    # Cycle coordinates creating a source with flux
    for d, source in enumerate(sel_LSM.sources):
        flux = getattr(source, 'flux')
        for i, stokes in enumerate(['I', 'Q', 'U', 'V']):
            stokes_val = getattr(flux, stokes)
            if stokes_val:
                # Get spectrum (only spi currently supported)
                tmp_spec = source.spectrum
                spi = [tmp_spec.spi if tmp_spec 
                            is not None else 0.0]
                ref_freq = [tmp_spec.freq0 if tmp_spec 
                                    is not None else 1.0]
    
                # Generate model flux
                model[d, :, i] = stokes_val\
                    * (sel_FREQ/ref_freq)**spi

    # Get model visibilities
    model_vis = np.zeros((n_row, n_chan, n_dir, 4), 
                            dtype=np.complex128)
    for s in range(n_dir):
        model_vis[:, :, s] = im_to_vis(
            model[s].reshape((1, n_chan, 4)),
            sel_UVW, 
            gen_lm[s].reshape((1, 2)), 
            sel_FREQ, 
            dtype=np.complex64)

    # Convert Stokes to corr
    in_schema = ['I', 'Q', 'U', 'V']
    out_schema = [['RR', 'RL'], ['LR', 'LL']]
    model_vis = convert(model_vis, in_schema, out_schema)

    return model_vis
Beispiel #9
0
def test_wcorrection_faceting_forward(tmp_path_factory):
    # construct kernel
    W = 5
    OS = 9
    kern = kernels.pack_kernel(kernels.kbsinc(W, oversample=OS), W, OS)
    nrow = 5000
    np.random.seed(0)
    # simulate some ficticious baselines rotated by an hour angle
    uvw = np.zeros((nrow, 3), dtype=np.float64)
    blpos = np.random.uniform(26, 10000, size=(25, 3))
    ntime = int(nrow / 25.0)
    d0 = np.pi / 4.0
    for n in range(25):
        for ih0, h0 in enumerate(
                np.linspace(np.deg2rad(-20), np.deg2rad(20), ntime)):
            s = np.sin
            c = np.cos
            R = np.array([[s(h0), c(h0), 0],
                          [-s(d0) * c(h0),
                           s(d0) * s(h0), c(d0)],
                          [c(d0) * c(h0), -c(d0) * s(h0),
                           s(d0)]])
            uvw[n * ntime + ih0, :] = np.dot(R, blpos[n, :].T)

    pxacrossbeam = 5
    frequency = np.array([1.4e9])
    wavelength = lightspeed / frequency

    cell = np.rad2deg(
        wavelength[0] /
        (max(np.max(np.abs(uvw[:, 0])), np.max(np.abs(uvw[:, 1]))) *
         pxacrossbeam))
    npixfacet = 100
    mod = np.ones((1, 1, 1), dtype=np.complex64)
    deltaradec = np.array([[20 * np.deg2rad(cell), 20 * np.deg2rad(cell)]])
    lm = radec_to_lmn(deltaradec + np.array([[0, d0]]),
                      phase_centre=np.array([0, d0]))

    vis_dft = im_to_vis(mod, uvw, lm[:, 0:2],
                        frequency).repeat(2).reshape(nrow, 1, 2)
    chanmap = np.array([0])
    ftmod = np.ones((1, npixfacet, npixfacet),
                    dtype=np.complex64)  # point source at centre of facet
    vis_degrid = degridder.degridder(
        uvw,
        ftmod,
        wavelength,
        chanmap,
        cell * 3600.0,
        (deltaradec + np.array([[0, d0]]))[0, :],
        (0, d0),
        kern,
        W,
        OS,
        "rotate",  # no faceting
        "phase_rotate",  # no faceting
        "XXYY_FROM_I",
        "conv_1d_axisymmetric_packed_gather")

    try:
        import matplotlib
    except ImportError:
        pass
    else:
        matplotlib.use("agg")
        from matplotlib import pyplot as plt
        plot_dir = tmp_path_factory.mktemp("wcorrection_forward")

        plt.figure()
        plt.plot(vis_degrid[:, 0, 0].real,
                 label=r"$\Re(\mathtt{degrid facet})$")
        plt.plot(vis_dft[:, 0, 0].real, label=r"$\Re(\mathtt{dft})$")
        plt.plot(np.abs(vis_dft[:, 0, 0].real - vis_degrid[:, 0, 0].real),
                 label="Error")
        plt.legend()
        plt.xlabel("sample")
        plt.ylabel("Real of predicted")
        plt.savefig(plot_dir / "facet_degrid_vs_dft_re_packed.png")
        plt.figure()
        plt.plot(vis_degrid[:, 0, 0].imag,
                 label=r"$\Im(\mathtt{degrid facet})$")
        plt.plot(vis_dft[:, 0, 0].imag, label=r"$\Im(\mathtt{dft})$")
        plt.plot(np.abs(vis_dft[:, 0, 0].imag - vis_degrid[:, 0, 0].imag),
                 label="Error")
        plt.legend()
        plt.xlabel("sample")
        plt.ylabel("Imag of predicted")
        plt.savefig(plot_dir / "facet_degrid_vs_dft_im_packed.png")

    assert np.percentile(
        np.abs(vis_dft[:, 0, 0].real - vis_degrid[:, 0, 0].real), 99.0) < 0.05
    assert np.percentile(
        np.abs(vis_dft[:, 0, 0].imag - vis_degrid[:, 0, 0].imag), 99.0) < 0.05
Beispiel #10
0
def test_wcorrection_faceting_backward(tmp_path_factory):
    # construct kernel
    W = 5
    OS = 9
    kern = kernels.pack_kernel(kernels.kbsinc(W, oversample=OS), W, OS)
    nrow = 5000
    np.random.seed(0)
    # simulate some ficticious baselines rotated by an hour angle
    uvw = np.zeros((nrow, 3), dtype=np.float64)
    blpos = np.random.uniform(26, 10000, size=(25, 3))
    ntime = int(nrow / 25.0)
    d0 = np.pi / 4.0
    for n in range(25):
        for ih0, h0 in enumerate(
                np.linspace(np.deg2rad(-20), np.deg2rad(20), ntime)):
            s = np.sin
            c = np.cos
            R = np.array([[s(h0), c(h0), 0],
                          [-s(d0) * c(h0),
                           s(d0) * s(h0), c(d0)],
                          [c(d0) * c(h0), -c(d0) * s(h0),
                           s(d0)]])
            uvw[n * ntime + ih0, :] = np.dot(R, blpos[n, :].T)

    pxacrossbeam = 5
    frequency = np.array([1.4e9])
    wavelength = lightspeed / frequency

    cell = np.rad2deg(
        wavelength[0] /
        (max(np.max(np.abs(uvw[:, 0])), np.max(np.abs(uvw[:, 1]))) *
         pxacrossbeam))
    npix = 2048
    npixfacet = 100
    fftpad = 1.1
    mod = np.ones((1, 1, 1), dtype=np.complex64)
    deltaradec = np.array([[600 * np.deg2rad(cell), 600 * np.deg2rad(cell)]])
    lm = radec_to_lmn(deltaradec + np.array([[0, d0]]),
                      phase_centre=np.array([0, d0]))

    vis_dft = im_to_vis(mod, uvw, lm[:, 0:2],
                        frequency).repeat(2).reshape(nrow, 1, 2)
    chanmap = np.array([0])

    detaper = kernels.compute_detaper_dft_seperable(
        int(npix * fftpad), kernels.unpack_kernel(kern, W, OS), W, OS)
    vis_grid_nofacet = gridder.gridder(
        uvw,
        vis_dft,
        wavelength,
        chanmap,
        int(npix * fftpad),
        cell * 3600.0,
        (0, d0),
        (0, d0),
        kern,
        W,
        OS,
        "None",  # no faceting
        "None",  # no faceting
        "I_FROM_XXYY",
        "conv_1d_axisymmetric_packed_scatter",
        do_normalize=True)
    ftvis = (np.fft.fftshift(
        np.fft.ifft2(np.fft.ifftshift(vis_grid_nofacet[0, :, :]))).reshape(
            (1, int(npix * fftpad), int(npix * fftpad)))).real / detaper * int(
                npix * fftpad)**2
    ftvis = ftvis[:,
                  int(npix * fftpad) // 2 - npix // 2:int(npix * fftpad) // 2 -
                  npix // 2 + npix,
                  int(npix * fftpad) // 2 - npix // 2:int(npix * fftpad) // 2 -
                  npix // 2 + npix]

    detaper_facet = kernels.compute_detaper_dft_seperable(
        int(npixfacet * fftpad), kernels.unpack_kernel(kern, W, OS), W, OS)
    vis_grid_facet = gridder.gridder(uvw,
                                     vis_dft,
                                     wavelength,
                                     chanmap,
                                     int(npixfacet * fftpad),
                                     cell * 3600.0,
                                     (deltaradec + np.array([[0, d0]]))[0, :],
                                     (0, d0),
                                     kern,
                                     W,
                                     OS,
                                     "rotate",
                                     "phase_rotate",
                                     "I_FROM_XXYY",
                                     "conv_1d_axisymmetric_packed_scatter",
                                     do_normalize=True)
    ftvisfacet = (np.fft.fftshift(
        np.fft.ifft2(np.fft.ifftshift(vis_grid_facet[0, :, :]))).reshape(
            (1, int(npixfacet * fftpad), int(
                npixfacet * fftpad)))).real / detaper_facet * int(
                    npixfacet * fftpad)**2
    ftvisfacet = ftvisfacet[:,
                            int(npixfacet * fftpad) // 2 -
                            npixfacet // 2:int(npixfacet * fftpad) // 2 -
                            npixfacet // 2 + npixfacet,
                            int(npixfacet * fftpad) // 2 -
                            npixfacet // 2:int(npixfacet * fftpad) // 2 -
                            npixfacet // 2 + npixfacet]

    try:
        import matplotlib
    except ImportError:
        pass
    else:
        matplotlib.use("agg")
        from matplotlib import pyplot as plt
        plot_dir = tmp_path_factory.mktemp("wcorrection_backward")

        plt.figure()
        plt.subplot(121)
        plt.imshow(ftvis[0, 1624 - 50:1624 + 50, 1447 - 50:1447 + 50])
        plt.colorbar()
        plt.title("Offset FFT (peak={0:.1f})".format(np.max(ftvis)))
        plt.subplot(122)
        plt.imshow(ftvisfacet[0, :, :])
        plt.colorbar()
        plt.title("Faceted FFT (peak={0:.1f})".format(np.max(ftvisfacet)))
        plt.savefig(plot_dir / "facet_imaging.png")

    assert (np.abs(np.max(ftvisfacet[0, :, :]) - 1.0) < 1.0e-6)
Beispiel #11
0
def test_grid_dft_packed(tmp_path_factory):
    # construct kernel
    W = 7
    OS = 1009
    kern = kernels.pack_kernel(kernels.kbsinc(W, oversample=OS), W, OS)
    nrow = 5000
    np.random.seed(0)
    uvw = np.random.normal(scale=6000, size=(nrow, 3))
    uvw[:, 2] = 0.0  # ignore widefield effects for now

    pxacrossbeam = 10
    frequency = np.array([30.0e9])
    wavelength = lightspeed / frequency

    cell = np.rad2deg(
        wavelength[0] /
        (2 * max(np.max(np.abs(uvw[:, 0])), np.max(np.abs(uvw[:, 1]))) *
         pxacrossbeam))
    npix = 256
    fftpad = 1.25
    mod = np.zeros((1, npix, npix), dtype=np.complex64)
    for n in [int(n) for n in np.linspace(npix // 8, 2 * npix // 5, 5)]:
        mod[0, npix // 2 + n, npix // 2 + n] = 1.0
        mod[0, npix // 2 + n, npix // 2 - n] = 1.0
        mod[0, npix // 2 - n, npix // 2 - n] = 1.0
        mod[0, npix // 2 - n, npix // 2 + n] = 1.0
        mod[0, npix // 2, npix // 2 + n] = 1.0
        mod[0, npix // 2, npix // 2 - n] = 1.0
        mod[0, npix // 2 - n, npix // 2] = 1.0
        mod[0, npix // 2 + n, npix // 2] = 1.0

    dec, ra = np.meshgrid(
        np.arange(-npix // 2, npix // 2) * np.deg2rad(cell),
        np.arange(-npix // 2, npix // 2) * np.deg2rad(cell))
    radec = np.column_stack((ra.flatten(), dec.flatten()))

    vis_dft = im_to_vis(mod[0, :, :].reshape(1, 1, npix * npix).T.copy(), uvw,
                        radec, frequency).repeat(2).reshape(nrow, 1, 2)
    chanmap = np.array([0])
    detaper = kernels.compute_detaper_dft_seperable(
        int(npix * fftpad), kernels.unpack_kernel(kern, W, OS), W, OS)
    vis_grid = gridder.gridder(
        uvw,
        vis_dft,
        wavelength,
        chanmap,
        int(npix * fftpad),
        cell * 3600.0,
        (0, np.pi / 4.0),
        (0, np.pi / 4.0),
        kern,
        W,
        OS,
        "None",  # no faceting
        "None",  # no faceting
        "I_FROM_XXYY",
        "conv_1d_axisymmetric_packed_scatter",
        do_normalize=True)

    ftvis = (np.fft.fftshift(np.fft.ifft2(np.fft.ifftshift(
        vis_grid[0, :, :]))).reshape((1, int(npix * fftpad), int(
            npix * fftpad)))).real / detaper * int(npix * fftpad)**2
    ftvis = ftvis[:,
                  int(npix * fftpad) // 2 - npix // 2:int(npix * fftpad) // 2 -
                  npix // 2 + npix,
                  int(npix * fftpad) // 2 - npix // 2:int(npix * fftpad) // 2 -
                  npix // 2 + npix]
    dftvis = vis_to_im(vis_dft, uvw, radec, frequency,
                       np.zeros(vis_dft.shape,
                                dtype=np.bool)).T.copy().reshape(
                                    2, 1, npix, npix) / nrow

    try:
        import matplotlib
    except ImportError:
        pass
    else:
        matplotlib.use("agg")
        from matplotlib import pyplot as plt
        plt.figure()
        plt.subplot(131)
        plt.title("FFT")
        plt.imshow(ftvis[0, :, :])
        plt.colorbar()
        plt.subplot(132)
        plt.title("DFT")
        plt.imshow(dftvis[0, 0, :, :])
        plt.colorbar()
        plt.subplot(133)
        plt.title("ABS diff")
        plt.imshow(np.abs(ftvis[0, :, :] - dftvis[0, 0, :, :]))
        plt.colorbar()
        plt.savefig(
            tmp_path_factory.mktemp("grid_dft_packed") /
            "grid_diff_dft_packed.png")

    assert (np.percentile(np.abs(ftvis[0, :, :] - dftvis[0, 0, :, :]), 95.0) <
            0.15)
Beispiel #12
0
def test_degrid_dft_packed(tmp_path_factory):
    # construct kernel
    W = 5
    OS = 3
    kern = kernels.pack_kernel(kernels.kbsinc(W, oversample=OS),
                               W,
                               oversample=OS)
    uvw = np.column_stack(
        (5000.0 * np.cos(np.linspace(0, 2 * np.pi, 1000)),
         5000.0 * np.sin(np.linspace(0, 2 * np.pi, 1000)), np.zeros(1000)))

    pxacrossbeam = 10
    frequency = np.array([1.4e9])
    wavelength = lightspeed / frequency

    cell = np.rad2deg(
        wavelength[0] /
        (2 * max(np.max(np.abs(uvw[:, 0])), np.max(np.abs(uvw[:, 1]))) *
         pxacrossbeam))
    npix = 512
    mod = np.zeros((1, npix, npix), dtype=np.complex64)
    mod[0, npix // 2 - 5, npix // 2 - 5] = 1.0

    ftmod = np.fft.ifftshift(np.fft.fft2(np.fft.fftshift(
        mod[0, :, :]))).reshape((1, npix, npix))
    chanmap = np.array([0])
    vis_degrid = degridder.degridder(
        uvw,
        ftmod,
        wavelength,
        chanmap,
        cell * 3600.0,
        (0, np.pi / 4.0),
        (0, np.pi / 4.0),
        kern,
        W,
        OS,
        "None",  # no faceting
        "None",  # no faceting
        "XXYY_FROM_I",
        "conv_1d_axisymmetric_packed_gather")

    dec, ra = np.meshgrid(
        np.arange(-npix // 2, npix // 2) * np.deg2rad(cell),
        np.arange(-npix // 2, npix // 2) * np.deg2rad(cell))
    radec = np.column_stack((ra.flatten(), dec.flatten()))

    vis_dft = im_to_vis(mod[0, :, :].reshape(1, 1, npix * npix).T.copy(), uvw,
                        radec, frequency)

    try:
        import matplotlib
    except ImportError:
        pass
    else:
        matplotlib.use("agg")
        from matplotlib import pyplot as plt
        plt.figure()
        plt.plot(vis_degrid[:, 0, 0].real, label=r"$\Re(\mathtt{degrid})$")
        plt.plot(vis_dft[:, 0, 0].real, label=r"$\Re(\mathtt{dft})$")
        plt.plot(np.abs(vis_dft[:, 0, 0].real - vis_degrid[:, 0, 0].real),
                 label="Error")
        plt.legend()
        plt.xlabel("sample")
        plt.ylabel("Real of predicted")
        plt.savefig(
            os.path.join(os.environ.get("TMPDIR", "/tmp"),
                         "degrid_vs_dft_re_packed.png"))
        plt.figure()
        plt.plot(vis_degrid[:, 0, 0].imag, label=r"$\Im(\mathtt{degrid})$")
        plt.plot(vis_dft[:, 0, 0].imag, label=r"$\Im(\mathtt{dft})$")
        plt.plot(np.abs(vis_dft[:, 0, 0].imag - vis_degrid[:, 0, 0].imag),
                 label="Error")
        plt.legend()
        plt.xlabel("sample")
        plt.ylabel("Imag of predicted")
        plt.savefig(
            tmp_path_factory.mktemp("degrid_dft_packed") /
            "degrid_vs_dft_im_packed.png")

    assert np.percentile(
        np.abs(vis_dft[:, 0, 0].real - vis_degrid[:, 0, 0].real), 99.0) < 0.05
    assert np.percentile(
        np.abs(vis_dft[:, 0, 0].imag - vis_degrid[:, 0, 0].imag), 99.0) < 0.05