299 lines
9.3 KiB
Python
299 lines
9.3 KiB
Python
import contextlib
|
|
import io
|
|
import numpy
|
|
from numpy.linalg import norm
|
|
import pytest
|
|
import warnings
|
|
|
|
from ..fdmath import vec, unvec
|
|
from ..fdfd import waveguide_2d, waveguide_3d, waveguide_cyl
|
|
|
|
|
|
OMEGA = 1 / 1500
|
|
|
|
|
|
def build_waveguide_3d_mode(
|
|
*,
|
|
slice_start: int,
|
|
polarity: int,
|
|
) -> tuple[numpy.ndarray, list[list[numpy.ndarray]], tuple[slice, slice, slice], dict[str, complex | numpy.ndarray]]:
|
|
epsilon = numpy.ones((3, 5, 5, 1), dtype=float)
|
|
dxes = [[numpy.ones(5), numpy.ones(5), numpy.ones(1)] for _ in range(2)]
|
|
slices = (slice(slice_start, slice_start + 1), slice(None), slice(None))
|
|
result = waveguide_3d.solve_mode(
|
|
0,
|
|
omega=OMEGA,
|
|
dxes=dxes,
|
|
axis=0,
|
|
polarity=polarity,
|
|
slices=slices,
|
|
epsilon=epsilon,
|
|
)
|
|
return epsilon, dxes, slices, result
|
|
|
|
|
|
def build_waveguide_cyl_fixture(
|
|
*,
|
|
nonuniform: bool = False,
|
|
) -> tuple[list[list[numpy.ndarray]], numpy.ndarray, float]:
|
|
if nonuniform:
|
|
dxes = [
|
|
[numpy.array([1.0, 1.5, 1.2, 0.8, 1.1]), numpy.ones(5)],
|
|
[numpy.array([0.9, 1.4, 1.0, 0.7, 1.2]), numpy.ones(5)],
|
|
]
|
|
else:
|
|
dxes = [[numpy.ones(5), numpy.ones(5)] for _ in range(2)]
|
|
epsilon = vec(numpy.ones((3, 5, 5), dtype=float))
|
|
return dxes, epsilon, 10.0
|
|
|
|
|
|
def test_waveguide_3d_solve_mode_and_expand_e_are_phase_consistent() -> None:
|
|
epsilon, dxes, slices, result = build_waveguide_3d_mode(slice_start=0, polarity=1)
|
|
axis = 0
|
|
polarity = 1
|
|
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)
|
|
|
|
|
|
@pytest.mark.parametrize(
|
|
('polarity', 'expected_range'),
|
|
[(1, (0, 1)), (-1, (3, 4))],
|
|
)
|
|
def test_waveguide_3d_compute_overlap_e_uses_adjacent_window(
|
|
polarity: int,
|
|
expected_range: tuple[int, int],
|
|
) -> None:
|
|
_epsilon, dxes, slices, result = build_waveguide_3d_mode(slice_start=2, polarity=polarity)
|
|
|
|
with warnings.catch_warnings(record=True) as caught:
|
|
overlap = waveguide_3d.compute_overlap_e(
|
|
E=result['E'],
|
|
wavenumber=result['wavenumber'],
|
|
dxes=dxes,
|
|
axis=0,
|
|
polarity=polarity,
|
|
slices=slices,
|
|
omega=OMEGA,
|
|
)
|
|
|
|
nonzero = numpy.argwhere(numpy.abs(overlap) > 0)
|
|
|
|
assert not caught
|
|
assert numpy.isfinite(overlap).all()
|
|
assert nonzero[:, 1].min() == expected_range[0]
|
|
assert nonzero[:, 1].max() == expected_range[1]
|
|
|
|
|
|
@pytest.mark.parametrize(
|
|
('polarity', 'slice_start', 'expected_index'),
|
|
[(1, 1, 0), (-1, 3, 4)],
|
|
)
|
|
def test_waveguide_3d_compute_overlap_e_warns_when_window_is_clipped(
|
|
polarity: int,
|
|
slice_start: int,
|
|
expected_index: int,
|
|
) -> None:
|
|
_epsilon, dxes, slices, result = build_waveguide_3d_mode(slice_start=slice_start, polarity=polarity)
|
|
|
|
with pytest.warns(RuntimeWarning, match='clipped'):
|
|
overlap = waveguide_3d.compute_overlap_e(
|
|
E=result['E'],
|
|
wavenumber=result['wavenumber'],
|
|
dxes=dxes,
|
|
axis=0,
|
|
polarity=polarity,
|
|
slices=slices,
|
|
omega=OMEGA,
|
|
)
|
|
|
|
nonzero = numpy.argwhere(numpy.abs(overlap) > 0)
|
|
|
|
assert numpy.isfinite(overlap).all()
|
|
assert nonzero[:, 1].min() == expected_index
|
|
assert nonzero[:, 1].max() == expected_index
|
|
|
|
|
|
@pytest.mark.parametrize(
|
|
('polarity', 'slice_start'),
|
|
[(1, 0), (-1, 4)],
|
|
)
|
|
def test_waveguide_3d_compute_overlap_e_rejects_empty_overlap_window(
|
|
polarity: int,
|
|
slice_start: int,
|
|
) -> None:
|
|
_epsilon, dxes, slices, result = build_waveguide_3d_mode(slice_start=slice_start, polarity=polarity)
|
|
|
|
with pytest.raises(ValueError, match='outside the domain'):
|
|
waveguide_3d.compute_overlap_e(
|
|
E=result['E'],
|
|
wavenumber=result['wavenumber'],
|
|
dxes=dxes,
|
|
axis=0,
|
|
polarity=polarity,
|
|
slices=slices,
|
|
omega=OMEGA,
|
|
)
|
|
|
|
|
|
def test_waveguide_3d_compute_overlap_e_rejects_zero_support_window() -> None:
|
|
_epsilon, dxes, slices, result = build_waveguide_3d_mode(slice_start=2, polarity=1)
|
|
|
|
with pytest.raises(ValueError, match='no overlap field support'):
|
|
waveguide_3d.compute_overlap_e(
|
|
E=numpy.zeros_like(result['E']),
|
|
wavenumber=result['wavenumber'],
|
|
dxes=dxes,
|
|
axis=0,
|
|
polarity=1,
|
|
slices=slices,
|
|
omega=OMEGA,
|
|
)
|
|
|
|
|
|
def test_waveguide_cyl_solved_modes_are_ordered_and_low_residual() -> None:
|
|
dxes, epsilon, rmin = build_waveguide_cyl_fixture()
|
|
|
|
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:
|
|
dxes, epsilon, rmin = build_waveguide_cyl_fixture()
|
|
|
|
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=rmin,
|
|
)
|
|
|
|
assert numpy.isfinite(linear_wavenumbers).all()
|
|
assert numpy.all(numpy.real(linear_wavenumbers) > 0)
|
|
assert numpy.all(numpy.diff(numpy.real(linear_wavenumbers)) <= 0)
|
|
|
|
|
|
def test_waveguide_cyl_dxes2t_matches_expected_radius_scaling() -> None:
|
|
dxes, _epsilon, rmin = build_waveguide_cyl_fixture(nonuniform=True)
|
|
Ta, Tb = waveguide_cyl.dxes2T(dxes, rmin)
|
|
|
|
ta = (rmin + numpy.cumsum(dxes[0][0])) / rmin
|
|
tb = (rmin + dxes[0][0] / 2 + numpy.cumsum(dxes[1][0])) / rmin
|
|
|
|
numpy.testing.assert_allclose(Ta.diagonal(), numpy.repeat(ta, dxes[0][1].size))
|
|
numpy.testing.assert_allclose(Tb.diagonal(), numpy.repeat(tb, dxes[1][1].size))
|
|
|
|
|
|
def test_waveguide_cyl_exy2e_and_exy2h_return_finite_full_fields() -> None:
|
|
dxes, epsilon, rmin = build_waveguide_cyl_fixture()
|
|
mu = vec(2 * numpy.ones((3, 5, 5), dtype=float))
|
|
e_xy, angular_wavenumber = waveguide_cyl.solve_mode(
|
|
0,
|
|
omega=OMEGA,
|
|
dxes=dxes,
|
|
epsilon=epsilon,
|
|
rmin=rmin,
|
|
)
|
|
|
|
e_field = waveguide_cyl.exy2e(
|
|
angular_wavenumber=angular_wavenumber,
|
|
omega=OMEGA,
|
|
dxes=dxes,
|
|
rmin=rmin,
|
|
epsilon=epsilon,
|
|
) @ e_xy
|
|
h_field = waveguide_cyl.exy2h(
|
|
angular_wavenumber=angular_wavenumber,
|
|
omega=OMEGA,
|
|
dxes=dxes,
|
|
rmin=rmin,
|
|
epsilon=epsilon,
|
|
mu=mu,
|
|
) @ e_xy
|
|
|
|
assert e_field.shape == (3 * 25,)
|
|
assert h_field.shape == (3 * 25,)
|
|
assert numpy.isfinite(e_field).all()
|
|
assert numpy.isfinite(h_field).all()
|
|
assert unvec(e_field, (5, 5)).shape == (3, 5, 5)
|
|
assert unvec(h_field, (5, 5)).shape == (3, 5, 5)
|
|
|
|
|
|
@pytest.mark.parametrize('use_mu', [False, True])
|
|
def test_waveguide_cyl_normalized_fields_are_unit_norm_and_silent(use_mu: bool) -> None:
|
|
dxes, epsilon, rmin = build_waveguide_cyl_fixture()
|
|
mu = vec(2 * numpy.ones((3, 5, 5), dtype=float)) if use_mu else None
|
|
e_xy, angular_wavenumber = waveguide_cyl.solve_mode(
|
|
0,
|
|
omega=OMEGA,
|
|
dxes=dxes,
|
|
epsilon=epsilon,
|
|
rmin=rmin,
|
|
)
|
|
|
|
output = io.StringIO()
|
|
with contextlib.redirect_stdout(output):
|
|
e_field, h_field = waveguide_cyl.normalized_fields_e(
|
|
e_xy,
|
|
angular_wavenumber=angular_wavenumber,
|
|
omega=OMEGA,
|
|
dxes=dxes,
|
|
rmin=rmin,
|
|
epsilon=epsilon,
|
|
mu=mu,
|
|
)
|
|
|
|
overlap = waveguide_2d.inner_product(e_field, h_field, dxes, conj_h=True)
|
|
|
|
assert output.getvalue() == ''
|
|
assert numpy.isfinite(e_field).all()
|
|
assert numpy.isfinite(h_field).all()
|
|
assert abs(overlap.real - 1.0) < 1e-10
|
|
assert abs(overlap.imag) < 1e-10
|