[tests] add some more tests around numerical self-consistency

This commit is contained in:
Jan Petykiewicz 2026-04-17 20:16:16 -07:00
commit 38a5c1a9aa
3 changed files with 258 additions and 0 deletions

View file

@ -0,0 +1,103 @@
import numpy
from numpy.linalg import norm
from ..fdmath import vec
from ..fdfd import waveguide_3d, waveguide_cyl
OMEGA = 1 / 1500
def test_waveguide_3d_solve_mode_and_expand_e_are_phase_consistent() -> None:
epsilon = numpy.ones((3, 5, 5, 1), dtype=float)
dxes = [[numpy.ones(5), numpy.ones(5), numpy.ones(1)] for _ in range(2)]
axis = 0
polarity = 1
slices = (slice(0, 1), slice(None), slice(None))
result = waveguide_3d.solve_mode(
0,
omega=OMEGA,
dxes=dxes,
axis=axis,
polarity=polarity,
slices=slices,
epsilon=epsilon,
)
expanded = waveguide_3d.expand_e(
E=result['E'],
wavenumber=result['wavenumber'],
dxes=dxes,
axis=axis,
polarity=polarity,
slices=slices,
)
dx_prop = 0.5 * numpy.array([dx[2][slices[2]] for dx in dxes]).sum()
expected_wavenumber = 2 / dx_prop * numpy.arcsin(result['wavenumber_2d'] * dx_prop / 2)
solved_slice = (slice(None), *slices)
assert result['E'].shape == epsilon.shape
assert result['H'].shape == epsilon.shape
assert numpy.isfinite(result['E']).all()
assert numpy.isfinite(result['H']).all()
assert abs(result['wavenumber'] - expected_wavenumber) < 1e-12
assert numpy.allclose(expanded[solved_slice], result['E'][solved_slice])
component, _x, y_index, z_index = numpy.unravel_index(
numpy.abs(result['E']).argmax(),
result['E'].shape,
)
values = expanded[component, :, y_index, z_index]
ratios = values[1:] / values[:-1]
expected_ratio = numpy.exp(-1j * result['wavenumber'])
numpy.testing.assert_allclose(ratios, expected_ratio, rtol=1e-6, atol=1e-9)
def test_waveguide_cyl_solved_modes_are_ordered_and_low_residual() -> None:
shape = (5, 5)
dxes = [[numpy.ones(shape[0]), numpy.ones(shape[1])] for _ in range(2)]
epsilon = vec(numpy.ones((3, *shape), dtype=float))
rmin = 10.0
e_xys, angular_wavenumbers = waveguide_cyl.solve_modes(
[0, 1],
omega=OMEGA,
dxes=dxes,
epsilon=epsilon,
rmin=rmin,
)
operator = waveguide_cyl.cylindrical_operator(OMEGA, dxes, epsilon, rmin=rmin)
assert numpy.all(numpy.diff(numpy.real(angular_wavenumbers)) <= 0)
for e_xy, angular_wavenumber in zip(e_xys, angular_wavenumbers, strict=True):
eigenvalue = (angular_wavenumber / rmin) ** 2
residual = norm(operator @ e_xy - eigenvalue * e_xy) / norm(e_xy)
assert residual < 1e-6
def test_waveguide_cyl_linear_wavenumbers_are_finite_and_ordered() -> None:
shape = (5, 5)
dxes = [[numpy.ones(shape[0]), numpy.ones(shape[1])] for _ in range(2)]
epsilon = vec(numpy.ones((3, *shape), dtype=float))
e_xys, angular_wavenumbers = waveguide_cyl.solve_modes(
[0, 1],
omega=OMEGA,
dxes=dxes,
epsilon=epsilon,
rmin=10.0,
)
linear_wavenumbers = waveguide_cyl.linear_wavenumbers(
e_xys,
angular_wavenumbers,
epsilon=epsilon,
dxes=dxes,
rmin=10.0,
)
assert numpy.isfinite(linear_wavenumbers).all()
assert numpy.all(numpy.real(linear_wavenumbers) > 0)
assert numpy.all(numpy.diff(numpy.real(linear_wavenumbers)) <= 0)