Skip to content

Commit

Permalink
Fix bug introduced in recent refactor
Browse files Browse the repository at this point in the history
* x, u are in the col direction, y, v in the row direction python is row, col order
  • Loading branch information
samaloney committed Jul 12, 2024
1 parent 0d737e5 commit 3cfd052
Show file tree
Hide file tree
Showing 3 changed files with 20 additions and 18 deletions.
14 changes: 7 additions & 7 deletions xrayvision/tests/test_imaging.py
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,7 @@
import numpy as np
import pytest
from astropy.convolution.kernels import Gaussian2DKernel
from numpy.testing import assert_allclose
from numpy.testing import assert_allclose, assert_array_equal
from sunpy.map import Map

from xrayvision.imaging import image_to_vis, map_to_vis, vis_psf_image, vis_to_image, vis_to_map
Expand Down Expand Up @@ -161,8 +161,8 @@ def test_map_to_vis(pos, pixel):
pixel = pixel * apu.arcsec / apu.pix

# Calculate full u, v coverage so will be equivalent to a discrete Fourier transform (DFT)
u = generate_uv(m * apu.pix, phase_center=pos[0], pixel_size=pixel[0])
v = generate_uv(n * apu.pix, phase_center=pos[1], pixel_size=pixel[1])
v = generate_uv(m * apu.pix, phase_center=pos[0], pixel_size=pixel[0])
u = generate_uv(n * apu.pix, phase_center=pos[1], pixel_size=pixel[1])
u, v = np.meshgrid(u, v, indexing="ij")
u, v = np.array([u, v]).reshape(2, size) / apu.arcsec

Expand All @@ -181,10 +181,10 @@ def test_map_to_vis(pos, pixel):
mp = Map((data, header))
vis = map_to_vis(mp, u=u, v=v)

assert np.array_equal(vis.phase_center, pos)
assert_array_equal(vis.phase_center, pos)

res = vis_to_image(vis, shape=(m, n) * apu.pixel, pixel_size=pixel)
assert np.allclose(res, data)
assert_allclose(res, data)


def test_vis_to_image():
Expand Down Expand Up @@ -255,8 +255,8 @@ def test_vis_to_image_invalid_pixel_size():
def test_vis_to_map(m, n, pos, pixel):
pos = pos * apu.arcsec
pixel = pixel * apu.arcsec / apu.pix
u = generate_uv(m * apu.pix, phase_center=pos[0], pixel_size=pixel[0])
v = generate_uv(n * apu.pix, phase_center=pos[1], pixel_size=pixel[1])
v = generate_uv(m * apu.pix, phase_center=pos[0], pixel_size=pixel[0])
u = generate_uv(n * apu.pix, phase_center=pos[1], pixel_size=pixel[1])
u, v = np.meshgrid(u, v, indexing="ij")
uv = np.array([u, v]).reshape(2, m * n) / apu.arcsec
u, v = uv
Expand Down
10 changes: 5 additions & 5 deletions xrayvision/tests/test_transform.py
Original file line number Diff line number Diff line change
Expand Up @@ -113,10 +113,10 @@ def test_generate_uv_offset_size(phase_center, pixel_size):
@pytest.mark.parametrize("center", [(0, 0), (-1, 1), (1, -1), (1, 1)])
def test_dft_idft(shape, pixel_size, center):
data = np.arange(np.prod(shape)).reshape(shape)
uu = generate_uv(
vv = generate_uv(
shape[0] * apu.pix, phase_center=center[0] * apu.arcsec, pixel_size=pixel_size[0] * apu.arcsec / apu.pix
)
vv = generate_uv(
uu = generate_uv(
shape[1] * apu.pix, phase_center=center[1] * apu.arcsec, pixel_size=pixel_size[1] * apu.arcsec / apu.pix
)
u, v = np.meshgrid(uu, vv, indexing="ij")
Expand Down Expand Up @@ -386,13 +386,13 @@ def test_fft_equivalence():
center = (1, 1)

data = np.arange(np.prod(shape)).reshape(shape)
uu = generate_uv(
vv = generate_uv(
shape[0] * apu.pix, phase_center=center[0] * apu.arcsec, pixel_size=pixel[0] * apu.arcsec / apu.pix
)
vv = generate_uv(
uu = generate_uv(
shape[1] * apu.pix, phase_center=center[1] * apu.arcsec, pixel_size=pixel[1] * apu.arcsec / apu.pix
)
u, v = np.meshgrid(uu, vv, indexing="ij")
u, v = np.meshgrid(uu, vv)
u = u.flatten()
v = v.flatten()

Expand Down
14 changes: 8 additions & 6 deletions xrayvision/transform.py
Original file line number Diff line number Diff line change
Expand Up @@ -154,10 +154,11 @@ def dft_map(
"""
m, n = input_array.shape * apu.pix
x = generate_xy(m, phase_center=phase_center[0], pixel_size=pixel_size[0]) # type: ignore
y = generate_xy(n, phase_center=phase_center[1], pixel_size=pixel_size[1]) # type: ignore
# python array index in row, column hence y, x
y = generate_xy(m, phase_center=phase_center[0], pixel_size=pixel_size[0]) # type: ignore
x = generate_xy(n, phase_center=phase_center[1], pixel_size=pixel_size[1]) # type: ignore

x, y = np.meshgrid(x, y, indexing="ij")
x, y = np.meshgrid(x, y)
uv = np.vstack([u, v])
# Check units are correct for exp need to be dimensionless and then remove units for speed
if (uv[0, :] * x[0, 0]).unit == apu.dimensionless_unscaled and (
Expand Down Expand Up @@ -220,10 +221,11 @@ def idft_map(
"""
m, n = shape
x = generate_xy(m, phase_center=phase_center[0], pixel_size=pixel_size[0]) # type: ignore
y = generate_xy(n, phase_center=phase_center[1], pixel_size=pixel_size[1]) # type: ignore
# python array index in row, column hence y, x
y = generate_xy(m, phase_center=phase_center[0], pixel_size=pixel_size[0]) # type: ignore
x = generate_xy(n, phase_center=phase_center[1], pixel_size=pixel_size[1]) # type: ignore

x, y = np.meshgrid(x, y, indexing="ij")
x, y = np.meshgrid(x, y)

if weights is None:
weights = np.ones(input_vis.shape)
Expand Down

0 comments on commit 3cfd052

Please sign in to comment.