213 lines
8.5 KiB
Python
213 lines
8.5 KiB
Python
import numpy
|
|
|
|
from ..fdmath import vec, unvec
|
|
from ..fdmath import functional as fd_functional
|
|
from ..fdfd import operators, scpml
|
|
from ._fdfd_case import (
|
|
BOUNDARY_DXES,
|
|
BOUNDARY_EPSILON,
|
|
BOUNDARY_FIELD,
|
|
BOUNDARY_SHAPE,
|
|
DXES,
|
|
EPSILON,
|
|
E_FIELD,
|
|
MU,
|
|
H_FIELD,
|
|
OMEGA,
|
|
PEC,
|
|
PMC,
|
|
SHAPE,
|
|
apply_fdfd_matrix,
|
|
)
|
|
from .utils import assert_close, assert_fields_close
|
|
|
|
|
|
def _dense_e_full(mu: numpy.ndarray | None) -> numpy.ndarray:
|
|
ce = fd_functional.curl_forward(DXES[0])
|
|
ch = fd_functional.curl_back(DXES[1])
|
|
pe = numpy.where(PEC, 0.0, 1.0)
|
|
pm = numpy.where(PMC, 0.0, 1.0)
|
|
|
|
masked_e = pe * E_FIELD
|
|
curl_term = ce(masked_e)
|
|
if mu is not None:
|
|
curl_term = curl_term / mu
|
|
curl_term = pm * curl_term
|
|
curl_term = ch(curl_term)
|
|
return pe * (curl_term - OMEGA**2 * EPSILON * masked_e)
|
|
|
|
|
|
def _dense_h_full(mu: numpy.ndarray | None) -> numpy.ndarray:
|
|
ce = fd_functional.curl_forward(DXES[0])
|
|
ch = fd_functional.curl_back(DXES[1])
|
|
pe = numpy.where(PEC, 0.0, 1.0)
|
|
pm = numpy.where(PMC, 0.0, 1.0)
|
|
magnetic = numpy.ones_like(EPSILON) if mu is None else mu
|
|
|
|
masked_h = pm * H_FIELD
|
|
curl_term = ch(masked_h)
|
|
curl_term = pe * (curl_term / EPSILON)
|
|
curl_term = ce(curl_term)
|
|
return pm * (curl_term - OMEGA**2 * magnetic * masked_h)
|
|
|
|
|
|
def _normalized_distance(u: numpy.ndarray, size: int, thickness: int) -> numpy.ndarray:
|
|
return ((thickness - u).clip(0) + (u - (size - thickness)).clip(0)) / thickness
|
|
|
|
|
|
def test_h_full_matches_dense_reference_with_and_without_mu() -> None:
|
|
for mu in (None, MU):
|
|
matrix_result = apply_fdfd_matrix(
|
|
operators.h_full(OMEGA, DXES, vec(EPSILON), None if mu is None else vec(mu), vec(PEC), vec(PMC)),
|
|
H_FIELD,
|
|
)
|
|
dense_result = _dense_h_full(mu)
|
|
assert_fields_close(matrix_result, dense_result, atol=1e-10, rtol=1e-10)
|
|
|
|
|
|
def test_e_full_matches_dense_reference_with_masks() -> None:
|
|
for mu in (None, MU):
|
|
matrix_result = apply_fdfd_matrix(
|
|
operators.e_full(OMEGA, DXES, vec(EPSILON), None if mu is None else vec(mu), vec(PEC), vec(PMC)),
|
|
E_FIELD,
|
|
)
|
|
dense_result = _dense_e_full(mu)
|
|
assert_fields_close(matrix_result, dense_result, atol=1e-10, rtol=1e-10)
|
|
|
|
|
|
def test_h_full_without_masks_matches_dense_reference() -> None:
|
|
ce = fd_functional.curl_forward(DXES[0])
|
|
ch = fd_functional.curl_back(DXES[1])
|
|
dense_result = ce(ch(H_FIELD) / EPSILON) - OMEGA**2 * MU * H_FIELD
|
|
matrix_result = apply_fdfd_matrix(
|
|
operators.h_full(OMEGA, DXES, vec(EPSILON), vec(MU)),
|
|
H_FIELD,
|
|
)
|
|
assert_fields_close(matrix_result, dense_result, atol=1e-10, rtol=1e-10)
|
|
|
|
|
|
def test_eh_full_matches_manual_block_operator_with_masks() -> None:
|
|
pe = numpy.where(PEC, 0.0, 1.0)
|
|
pm = numpy.where(PMC, 0.0, 1.0)
|
|
ce = fd_functional.curl_forward(DXES[0])
|
|
ch = fd_functional.curl_back(DXES[1])
|
|
|
|
matrix_result = operators.eh_full(OMEGA, DXES, vec(EPSILON), vec(MU), vec(PEC), vec(PMC)) @ numpy.concatenate(
|
|
[vec(E_FIELD), vec(H_FIELD)],
|
|
)
|
|
matrix_e, matrix_h = (unvec(part, SHAPE) for part in numpy.split(matrix_result, 2))
|
|
|
|
dense_e = pe * ch(pm * H_FIELD) - pe * (1j * OMEGA * EPSILON * (pe * E_FIELD))
|
|
dense_h = pm * ce(pe * E_FIELD) + pm * (1j * OMEGA * MU * (pm * H_FIELD))
|
|
|
|
assert_fields_close(matrix_e, dense_e, atol=1e-10, rtol=1e-10)
|
|
assert_fields_close(matrix_h, dense_h, atol=1e-10, rtol=1e-10)
|
|
|
|
|
|
def test_e2h_pmc_mask_matches_masked_unmasked_result() -> None:
|
|
pmc_complement = numpy.where(PMC, 0.0, 1.0)
|
|
unmasked = apply_fdfd_matrix(operators.e2h(OMEGA, DXES, vec(MU)), E_FIELD)
|
|
masked = apply_fdfd_matrix(operators.e2h(OMEGA, DXES, vec(MU), vec(PMC)), E_FIELD)
|
|
|
|
assert_fields_close(masked, pmc_complement * unmasked, atol=1e-10, rtol=1e-10)
|
|
|
|
|
|
def test_poynting_h_cross_matches_negative_e_cross_relation() -> None:
|
|
h_cross_e = apply_fdfd_matrix(operators.poynting_h_cross(vec(H_FIELD), DXES), E_FIELD)
|
|
e_cross_h = apply_fdfd_matrix(operators.poynting_e_cross(vec(E_FIELD), DXES), H_FIELD)
|
|
|
|
assert_fields_close(h_cross_e, -e_cross_h, atol=1e-10, rtol=1e-10)
|
|
|
|
|
|
def test_e_boundary_source_interior_mask_is_independent_of_periodic_edges() -> None:
|
|
mask = numpy.zeros((3, *BOUNDARY_SHAPE), dtype=float)
|
|
mask[:, 1, 1, 1] = 1.0
|
|
|
|
periodic = operators.e_boundary_source(vec(mask), OMEGA, BOUNDARY_DXES, vec(BOUNDARY_EPSILON), periodic_mask_edges=True)
|
|
mirrored = operators.e_boundary_source(vec(mask), OMEGA, BOUNDARY_DXES, vec(BOUNDARY_EPSILON), periodic_mask_edges=False)
|
|
|
|
assert_close(periodic.toarray(), mirrored.toarray())
|
|
|
|
|
|
def test_e_boundary_source_periodic_edges_add_opposite_face_response() -> None:
|
|
mask = numpy.zeros((3, *BOUNDARY_SHAPE), dtype=float)
|
|
mask[:, 0, 1, 1] = 1.0
|
|
|
|
periodic = operators.e_boundary_source(vec(mask), OMEGA, BOUNDARY_DXES, vec(BOUNDARY_EPSILON), periodic_mask_edges=True)
|
|
mirrored = operators.e_boundary_source(vec(mask), OMEGA, BOUNDARY_DXES, vec(BOUNDARY_EPSILON), periodic_mask_edges=False)
|
|
diff = unvec((periodic - mirrored) @ vec(BOUNDARY_FIELD), BOUNDARY_SHAPE)
|
|
|
|
assert numpy.isfinite(diff).all()
|
|
assert_close(diff[:, 1:-1, :, :], 0.0)
|
|
assert numpy.linalg.norm(diff[:, -1, :, :]) > 0
|
|
|
|
|
|
def test_prepare_s_function_matches_closed_form_polynomial() -> None:
|
|
ln_r = -12.0
|
|
order = 3.0
|
|
distances = numpy.array([0.0, 0.25, 0.5, 1.0])
|
|
s_function = scpml.prepare_s_function(ln_R=ln_r, m=order)
|
|
expected = (order + 1) * ln_r / 2 * distances**order
|
|
|
|
assert_close(s_function(distances), expected)
|
|
|
|
|
|
def test_uniform_grid_scpml_matches_expected_stretch_profile() -> None:
|
|
s_function = scpml.prepare_s_function(ln_R=-12.0, m=3.0)
|
|
dxes = scpml.uniform_grid_scpml((6, 4, 3), (2, 0, 1), omega=2.0, epsilon_effective=4.0, s_function=s_function)
|
|
correction = numpy.sqrt(4.0) * 2.0
|
|
|
|
for axis, size, thickness in ((0, 6, 2), (2, 3, 1)):
|
|
grid = numpy.arange(size, dtype=float)
|
|
expected_a = 1 + 1j * s_function(_normalized_distance(grid, size, thickness)) / correction
|
|
expected_b = 1 + 1j * s_function(_normalized_distance(grid + 0.5, size, thickness)) / correction
|
|
assert_close(dxes[0][axis], expected_a)
|
|
assert_close(dxes[1][axis], expected_b)
|
|
|
|
assert_close(dxes[0][1], 1.0)
|
|
assert_close(dxes[1][1], 1.0)
|
|
assert numpy.isfinite(dxes[0][0]).all()
|
|
assert numpy.isfinite(dxes[1][0]).all()
|
|
|
|
|
|
def test_uniform_grid_scpml_default_s_function_matches_explicit_default() -> None:
|
|
implicit = scpml.uniform_grid_scpml((6, 4, 3), (2, 0, 1), omega=2.0)
|
|
explicit = scpml.uniform_grid_scpml((6, 4, 3), (2, 0, 1), omega=2.0, s_function=scpml.prepare_s_function())
|
|
|
|
for implicit_group, explicit_group in zip(implicit, explicit, strict=True):
|
|
for implicit_axis, explicit_axis in zip(implicit_group, explicit_group, strict=True):
|
|
assert_close(implicit_axis, explicit_axis)
|
|
|
|
|
|
def test_stretch_with_scpml_only_modifies_requested_front_edge() -> None:
|
|
s_function = scpml.prepare_s_function(ln_R=-12.0, m=3.0)
|
|
base = [[numpy.ones(6), numpy.ones(4), numpy.ones(3)] for _ in range(2)]
|
|
stretched = scpml.stretch_with_scpml(base, axis=0, polarity=1, omega=2.0, epsilon_effective=4.0, thickness=2, s_function=s_function)
|
|
|
|
assert_close(stretched[0][0][2:], 1.0)
|
|
assert_close(stretched[1][0][2:], 1.0)
|
|
assert_close(stretched[0][0][-2:], 1.0)
|
|
assert_close(stretched[1][0][-2:], 1.0)
|
|
assert numpy.linalg.norm(stretched[0][0][:2] - 1.0) > 0
|
|
assert numpy.linalg.norm(stretched[1][0][:2] - 1.0) > 0
|
|
|
|
|
|
def test_stretch_with_scpml_only_modifies_requested_back_edge() -> None:
|
|
s_function = scpml.prepare_s_function(ln_R=-12.0, m=3.0)
|
|
base = [[numpy.ones(6), numpy.ones(4), numpy.ones(3)] for _ in range(2)]
|
|
stretched = scpml.stretch_with_scpml(base, axis=0, polarity=-1, omega=2.0, epsilon_effective=4.0, thickness=2, s_function=s_function)
|
|
|
|
assert_close(stretched[0][0][:4], 1.0)
|
|
assert_close(stretched[1][0][:4], 1.0)
|
|
assert numpy.linalg.norm(stretched[0][0][-2:] - 1.0) > 0
|
|
assert numpy.linalg.norm(stretched[1][0][-2:] - 1.0) > 0
|
|
|
|
|
|
def test_stretch_with_scpml_thickness_zero_is_noop() -> None:
|
|
s_function = scpml.prepare_s_function(ln_R=-12.0, m=3.0)
|
|
base = [[numpy.ones(6), numpy.ones(4), numpy.ones(3)] for _ in range(2)]
|
|
stretched = scpml.stretch_with_scpml(base, axis=0, polarity=-1, omega=2.0, epsilon_effective=4.0, thickness=0, s_function=s_function)
|
|
|
|
for grid_group in stretched:
|
|
for axis_grid in grid_group:
|
|
assert_close(axis_grid, 1.0)
|