meanas
meanas is a python package for electromagnetic simulations
** UNSTABLE / WORK IN PROGRESS **
Formerly known as fdfd_tools.
This package is intended for building simulation inputs, analyzing simulation outputs, and running short simulations on unspecialized hardware. It is designed to provide tooling and a baseline for other, high-performance purpose- and hardware-specific solvers.
Contents
This package does not provide a fast matrix solver, though
by default generic()(…)
will call
scipy.sparse.linalg.qmr(…)
to perform a solve. For 2D FDFD
problems this should be fine; likewise, the waveguide mode solver uses
scipy’s eigenvalue solver, with reasonable results.
For solving large (or 3D) FDFD problems, I recommend a GPU-based iterative solver, such as opencl_fdfd or those included in MAGMA. Your solver will need the ability to solve complex symmetric (non-Hermitian) linear systems, ideally with double precision.
Requirements:
Install from PyPI with pip:
pip3 install 'meanas[dev]'
Install python3 and git:
# This is for Debian/Ubuntu/other-apt-based systems; you may need an alternative command
sudo apt install python3 build-essential python3-dev git
In-place development install:
# Download using git
git clone https://mpxd.net/code/jan/meanas.git
# If you'd like to create a virtualenv, do so:
python3 -m venv my_venv
# If you are using a virtualenv, activate it
source my_venv/bin/activate
# Install in-place (-e, editable) from ./meanas, including development dependencies ([dev])
pip3 install --user -e './meanas[dev]'
# Run tests
cd meanas
python3 -m pytest -rsxX | tee test_results.txt
See examples/
for some simple examples; you may need
additional packages such as gridlock to run the
examples.
meanas.eigensolvers
Solvers for eigenvalue / eigenvector problems
power_iteration
def power_iteration(operator: scipy.sparse._matrix.spmatrix, guess_vector: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]] | None = None, iterations: int = 20) -> tuple[complex, numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Use power iteration to estimate the dominant eigenvector of a matrix.
Args —–= operator
: Matrix to
analyze.
guess_vector
iterations
Returns —–= (Largest-magnitude eigenvalue, Corresponding eigenvector estimate)
rayleigh_quotient_iteration
def rayleigh_quotient_iteration(operator: scipy.sparse._matrix.spmatrix | scipy.sparse.linalg._interface.LinearOperator, guess_vector: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], iterations: int = 40, tolerance: float = 1e-13, solver: collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]] | None = None) -> tuple[complex, numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Use Rayleigh quotient iteration to refine an eigenvector guess.
Args —–= operator
: Matrix to
analyze.
guess_vector
iterations
tolerance
(A - I*eigenvalue) @ v < num_vectors * tolerance
,
Default 1e-13.
solver
x = solver(A, b)
. By default,
use scipy.sparse.spsolve for sparse matrices and scipy.sparse.bicgstab
for general LinearOperator instances.
Returns —–= (eigenvalues, eigenvectors)
signed_eigensolve
def signed_eigensolve(operator: scipy.sparse._matrix.spmatrix | scipy.sparse.linalg._interface.LinearOperator, how_many: int, negative: bool = False) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Find the largest-magnitude positive-only (or negative-only) eigenvalues and eigenvectors of the provided matrix.
Args —–= operator
: Matrix to
analyze.
how_many
negative
Returns —–= (sorted list of eigenvalues, 2D ndarray of corresponding
eigenvectors) eigenvectors[:, k]
corresponds to the k-th
eigenvalue
meanas.fdfd
Tools for finite difference frequency-domain (FDFD) simulations and calculations.
These mostly involve picking a single frequency, then setting up and solving a matrix equation (Ax=b) or eigenvalue problem.
Submodules:
meanas.fdfd.operators
, meanas.fdfd.functional
:
General FDFD problem setup.meanas.fdfd.solvers
:
Solver interface and reference implementation.meanas.fdfd.scpml
:
Stretched-coordinate perfectly matched layer (scpml) boundary
conditionsmeanas.fdfd.waveguide_2d
:
Operators and mode-solver for waveguides with constant
cross-section.meanas.fdfd.waveguide_3d
:
Functions for transforming meanas.fdfd.waveguide_2d
results into 3D.================================================================
From the “Frequency domain” section of meanas.fdmath
, we have
resulting in
Maxwell’s equations are then
With
and then
meanas.fdfd.bloch
Bloch eigenmode solver/operators
This module contains functions for generating and solving the 3D Bloch eigenproblem. The approach is to transform the problem into the (spatial) fourier domain, transforming the equation
1/mu * curl(1/eps * curl(H_eigenmode)) = (w/c)^2 H_eigenmode
into
conv(1/mu_k, ik x conv(1/eps_k, ik x H_k)) = (w/c)^2 H_k
where:
_k
subscript denotes a 3D fourier transformed
fieldH_k
corresponds to a plane wave with
wavevector k
x
is the cross productconv()
denotes convolutionSince k
and H
are orthogonal for each plane
wave, we can use each k
to create an orthogonal basis (k,
m, n), with k x m = n
, and |m| = |n| = 1
. The
cross products are then simplified as follows:
h
is shorthand for H_k
(…)_xyz
denotes the (x, y, z)
basis(…)_kmn
denotes the (k, m, n)
basishm
is the component of h
in the
m
direction, etc.We know
k @ h = kx hx + ky hy + kz hz = 0 = hk
h = hk + hm + hn = hm + hn
k = kk + km + kn = kk = |k|
We can write
k x h = (ky hz - kz hy,
kz hx - kx hz,
kx hy - ky hx)_xyz
= ((k x h) @ k, (k x h) @ m, (k x h) @ n)_kmn
= (0, (m x k) @ h, (n x k) @ h)_kmn # triple product ordering
= (0, kk (-n @ h), kk (m @ h))_kmn # (m x k) = -|k| n, etc.
= |k| (0, -h @ n, h @ m)_kmn
which gives us a straightforward way to perform the cross product
while simultaneously transforming into the _kmn
basis. We
can also write
k x h = (km hn - kn hm,
kn hk - kk hn,
kk hm - km hk)_kmn
= (0, -kk hn, kk hm)_kmn
= (-kk hn)(mx, my, mz)_xyz + (kk hm)(nx, ny, nz)_xyz
= |k| (hm * (nx, ny, nz)_xyz
- hn * (mx, my, mz)_xyz)
which gives us a way to perform the cross product while
simultaneously trasnforming back into the _xyz
basis.
We can also simplify conv(X_k, Y_k)
as
fftn(X * ifftn(Y_k))
.
Using these results and storing H_k
as
h = (hm, hn)
, we have
e_xyz = fftn(1/eps * ifftn(|k| (hm * n - hn * m)))
b_mn = |k| (-e_xyz @ n, e_xyz @ m)
h_mn = fftn(1/mu * ifftn(b_m * m + b_n * n))
which forms the operator from the left side of the equation.
We can then use a preconditioned block Rayleigh iteration algorithm, as in SG Johnson and JD Joannopoulos, Block-iterative frequency-domain methods for Maxwell’s equations in a planewave basis, Optics Express 8, 3, 173-190 (2001) (similar to that used in MPB) to find the eigenvectors for this operator.
===
Typically you will want to do something like
recip_lattice = numpy.diag(1/numpy.array(epsilon[0].shape * dx))
n, v = bloch.eigsolve(5, k0, recip_lattice, epsilon)
f = numpy.sqrt(-numpy.real(n[0]))
n_eff = norm(recip_lattice @ k0) / f
v2e = bloch.hmn_2_exyz(k0, recip_lattice, epsilon)
e_field = v2e(v[0])
k, f = find_k(frequency=1/1550,
tolerance=(1/1550 - 1/1551),
direction=[1, 0, 0],
G_matrix=recip_lattice,
epsilon=epsilon,
band=0)
eigsolve
def eigsolve(num_modes: int, k0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, tolerance: float = 1e-07, max_iters: int = 10000, reset_iters: int = 100, y0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]], ForwardRef(None)] = None, callback: collections.abc.Callable[..., None] | None = None) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Find the first (lowest-frequency) num_modes eigenmodes with Bloch wavevector k0 of the specified structure.
Args —–= k0
: Bloch wavevector,
[k0x, k0y, k0z]
.
G_matrix
epsilon
mu
None
(1 everywhere).
tolerance
trace(Z.H @ A @ Z @ inv(Z Z.H))
is smaller than the
tolerance
max_iters
reset_iters
callback
y0
Returns —–= (eigenvalues, eigenvectors)
where
eigenvalues[i]
corresponds to the vector
eigenvectors[i, :]
fftn
def fftn(*args: Any, **kwargs: Any) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]
find_k
def find_k(frequency: float, tolerance: float, direction: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, band: int = 0, k_bounds: tuple[float, float] = (0, 0.5), k_guess: float | None = None, solve_callback: collections.abc.Callable[..., None] | None = None, iter_callback: collections.abc.Callable[..., None] | None = None, v0: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]] | None = None) -> tuple[float, float, numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Search for a bloch vector that has a given frequency.
Args —–= frequency
: Target
frequency.
tolerance
direction
G_matrix
epsilon
mu
band
k_bounds
k_guess
solve_callback
iter_callback
Returns —–= (k, actual_frequency, eigenvalues,
eigenvectors)
The found k-vector and its frequency, along with
all eigenvalues and eigenvectors.
generate_kmn
def generate_kmn(k0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], shape: collections.abc.Sequence[int]) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]
Generate a (k, m, n) orthogonal basis for each k-vector in the simulation grid.
Args —–= k0
: [k0x, k0y, k0z], Bloch
wavevector, in G basis.
G_matrix
shape
Returns —–= (|k|, m, n)
where |k|
has shape
tuple(shape) + (1,)
and m
, n
have
shape tuple(shape) + (3,)
. All are given in the xyz basis
(e.g. |k|[0,0,0] = norm(G_matrix @ k0)
).
hmn_2_exyz
def hmn_2_exyz(k0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Generate an operator which converts a vectorized
spatial-frequency-space h_mn
into an E-field distribution,
i.e.
ifft(conv(1/eps_k, ik x h_mn))
The operator is a function that acts on a vector h_mn
of
size 2 * epsilon[0].size
.
See the meanas.fdfd.bloch
docstring for
more information.
Args —–= k0
: Bloch wavevector,
[k0x, k0y, k0z]
.
G_matrix
epsilon
Returns —–= Function for converting h_mn
into
E_xyz
hmn_2_hxyz
def hmn_2_hxyz(k0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Generate an operator which converts a vectorized
spatial-frequency-space h_mn
into an H-field distribution,
i.e.
ifft(h_mn)
The operator is a function that acts on a vector h_mn
of
size 2 * epsilon[0].size
.
See the meanas.fdfd.bloch
docstring for
more information.
Args —–= k0
: Bloch wavevector,
[k0x, k0y, k0z]
.
G_matrix
epsilon
epsilon[0].shape
is used.
Returns —–= Function for converting h_mn
into
H_xyz
ifftn
def ifftn(*args: Any, **kwargs: Any) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]
inner_product
def inner_product(eL, hL, eR, hR) -> complex
inverse_maxwell_operator_approx
def inverse_maxwell_operator_approx(k0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Generate an approximate inverse of the Maxwell operator,
ik x conv(eps_k, ik x conv(mu_k, ___))
which can be used to improve the speed of ARPACK in shift-invert mode.
See the meanas.fdfd.bloch
docstring for
more information.
Args —–= k0
: Bloch wavevector,
[k0x, k0y, k0z]
.
G_matrix
epsilon
mu
Returns —–= Function which applies the approximate inverse of the
maxwell operator to h_mn
.
maxwell_operator
def maxwell_operator(k0: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], G_matrix: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Generate the Maxwell operator
conv(1/mu_k, ik x conv(1/eps_k, ik x ___))
which is the spatial-frequency-space representation of
1/mu * curl(1/eps * curl(___))
The operator is a function that acts on a vector h_mn of size
2 * epsilon[0].size
See the meanas.fdfd.bloch
docstring for
more information.
Args —–= k0
: Bloch wavevector,
[k0x, k0y, k0z]
.
G_matrix
epsilon
mu
Returns —–= Function which applies the maxwell operator to h_mn.
trq
def trq(eI, hI, eO, hO) -> tuple[complex, complex]
meanas.fdfd.farfield
Functions for performing near-to-farfield transformation (and the reverse).
far_to_nearfield
def far_to_nearfield(E_far: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], H_far: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], dkx: float, dky: float, padded_size: list[int] | int | None = None) -> dict[str, typing.Any]
Compute the farfield, i.e. the distribution of the fields after propagation through several wavelengths of uniform medium.
The input fields should be complex phasors.
Args —–= E_far
: List of 2 ndarrays
containing the 2D phasor field slices for the transverse E fields
(e.g. [Ex, Ey] for calculating the nearfield toward the z-direction).
Fields should be normalized so that E_far = E_far_actual / (i k exp(-i k
r) / (4 pi r))
H_far
dkx
dky
padded_size
n
will be expanded to
(n, n)
. Powers of 2 are most efficient for FFT computation.
Default is the smallest power of 2 larger than the input, for each axis.
Returns —–= Dict with keys
E
: E-field nearfieldH
: H-field nearfielddx
, dy
: spatial discretization, normalized
to wavelength (dimensionless)near_to_farfield
def near_to_farfield(E_near: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], H_near: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], dx: float, dy: float, padded_size: list[int] | int | None = None) -> dict[str, typing.Any]
Compute the farfield, i.e. the distribution of the fields after propagation through several wavelengths of uniform medium.
The input fields should be complex phasors.
Args —–= E_near
: List of 2 ndarrays
containing the 2D phasor field slices for the transverse E fields
(e.g. [Ex, Ey] for calculating the farfield toward the z-direction).
H_near
dx
dy
padded_size
n
will be expanded to
(n, n)
. Powers of 2 are most efficient for FFT computation.
Default is the smallest power of 2 larger than the input, for each axis.
Returns —–= Dict with keys
E_far
: Normalized E-field farfield; multiply by (i k
exp(-i k r) / (4 pi r)) to get the actual field value.H_far
: Normalized H-field farfield; multiply by (i k
exp(-i k r) / (4 pi r)) to get the actual field value.kx
, ky
: Wavevector values corresponding to
the x- and y- axes in E_far and H_far, normalized to wavelength
(dimensionless).dkx
, dky
: step size for kx and ky,
normalized to wavelength.theta
: arctan2(ky, kx) corresponding to each (kx, ky).
This is the angle in the x-y plane, counterclockwise from above,
starting from +x.phi
: arccos(kz / k) corresponding to each (kx, ky).
This is the angle away from +z.meanas.fdfd.functional
Functional versions of many FDFD operators. These can be useful for performing FDFD calculations without needing to construct large matrices in memory.
The functions generated here expect cfdfield_t
inputs
with shape (3, X, Y, Z), e.g. E = [E_x, E_y, E_z] where each (complex)
component has shape (X, Y, Z)
e2h
def e2h(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Utility operator for converting the E
field into the
H
field. For use with e_full()
– assumes that
there is no magnetic current M
.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
mu
Returns —–= Function f
for converting E
to
H
, f(E)
-> H
e_full
def e_full(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Wave operator for use with E-field. See operators.e_full
for details.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
Returns —–= Function f
implementing the wave operator
f(E)
-> -i * omega * J
e_tfsf_source
def e_tfsf_source(TF_region: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Operator that turns an E-field distribution into a total-field/scattered-field (TFSF) source.
Args —–= TF_region
: mask which is set
to 1 in the total-field region, and 0 elsewhere (i.e. in the
scattered-field region). Should have the same shape as the simulation
grid, e.g. epsilon[0].shape
.
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
Returns —–= Function f
which takes an E field and
returns a current distribution, f(E)
->
J
eh_full
def eh_full(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]], tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]]
Wave operator for full (both E and H) field representation. See
operators.eh_full
.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
Returns —–= Function f
implementing the wave operator
f(E, H)
-> (J, -M)
m2j
def m2j(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Utility operator for converting magnetic current M
distribution into equivalent electric current distribution
J
. For use with e.g. e_full()
.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
mu
Returns —–= Function f
for converting M
to
J
, f(M)
-> J
poynting_e_cross_h
def poynting_e_cross_h(dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]]) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Generates a function that takes the single-frequency E
and H
fields and calculates the cross product
E
x H
=
Note —–= This function also shifts the input E
field by
one cell as required for computing the Poynting cross product (see
meanas.fdfd
module docs).
Note —–= If E
and H
are peak amplitudes as
assumed elsewhere in this code, the time-average of the poynting vector
is <S> = Re(S)/2 = Re(E x H*) / 2
. The factor of
1/2
can be omitted if root-mean-square quantities are used
instead.
Args —–= dxes
: Grid parameters
[dx_e, dx_h]
as described in meanas.fdmath.types
Returns —–= Function f
that returns E x H as required
for the poynting vector.
meanas.fdfd.operators
Sparse matrix operators for use with electromagnetic wave equations.
These functions return sparse-matrix
(scipy.sparse.spmatrix
) representations of a variety of
operators, intended for use with E and H fields vectorized using the
vec()
and
unvec()
functions.
E- and H-field values are defined on a Yee cell; epsilon
values should be calculated for cells centered at each E component
(mu
at each H component).
Many of these functions require a dxes
parameter, of
type dx_lists_t
; see the meanas.fdmath.types
submodule for
details.
The following operators are included:
e2h
def e2h(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Utility operator for converting the E field into the H field. For use
with e_full()
–
assumes that there is no magnetic current M.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
mu
pmc
pmc != 0
are interpreted as containing a perfect magnetic
conductor (PMC). The PMC is applied per-field-component
(i.e. pmc.size == epsilon.size
)
Returns —–= Sparse matrix for converting E to H.
e_boundary_source
def e_boundary_source(mask: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, periodic_mask_edges: bool = False) -> scipy.sparse._matrix.spmatrix
Operator that turns an E-field distrubtion into a current (J)
distribution along the edges (external and internal) of the provided
mask. This is just an e_tfsf_source()
with an additional masking step.
Args —–= mask
: The current
distribution is generated at the edges of the mask, i.e. any points
where shifting the mask by one cell in any direction would change its
value.
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
Returns —–= Sparse matrix that turns an E-field into a current (J) distribution.
e_full
def e_full(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Wave operator
del x (1/mu * del x) - omega**2 * epsilon
for use with the E-field, with wave equation
(del x (1/mu * del x) - omega**2 * epsilon) E = -i * omega * J
To make this matrix symmetric, use the preconditioners from e_full_preconditioners()
.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
pec
pec != 0
are interpreted as containing a perfect electrical
conductor (PEC). The PEC is applied per-field-component
(i.e. pec.size == epsilon.size
)
pmc
pmc != 0
are interpreted as containing a perfect magnetic
conductor (PMC). The PMC is applied per-field-component
(i.e. pmc.size == epsilon.size
)
Returns —–= Sparse matrix containing the wave operator.
e_full_preconditioners
def e_full_preconditioners(dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]]) -> tuple[scipy.sparse._matrix.spmatrix, scipy.sparse._matrix.spmatrix]
Left and right preconditioners (Pl, Pr)
for symmetrizing
the e_full()
wave operator.
The preconditioned matrix A_symm = (Pl @ A @ Pr)
is
complex-symmetric (non-Hermitian unless there is no loss or PMLs).
The preconditioner matrices are diagonal and complex, with
Pr = 1 / Pl
Args —–= dxes
: Grid parameters
[dx_e, dx_h]
as described in meanas.fdmath.types
Returns —–= Preconditioner matrices (Pl, Pr)
.
e_tfsf_source
def e_tfsf_source(TF_region: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Operator that turns a desired E-field distribution into a total-field/scattered-field (TFSF) source.
TODO: Reference Rumpf paper
Args —–= TF_region
: Mask, which is set
to 1 inside the total-field region and 0 in the scattered-field
region
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
Returns —–= Sparse matrix that turns an E-field into a current (J) distribution.
eh_full
def eh_full(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Wave operator for [E, H]
field representation. This
operator implements Maxwell’s equations without cancelling out either E
or H. The operator is
[[-i * omega * epsilon, del x ],
[del x, i * omega * mu]]
for use with a field vector of the form cat(vec(E),
vec(H))
:
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
pec
pec != 0
are interpreted as containing a perfect electrical
conductor (PEC). The PEC is applied per-field-component
(i.e. pec.size == epsilon.size
)
pmc
pmc != 0
are interpreted as containing a perfect magnetic
conductor (PMC). The PMC is applied per-field-component
(i.e. pmc.size == epsilon.size
)
Returns —–= Sparse matrix containing the wave operator.
h_full
def h_full(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Wave operator
del x (1/epsilon * del x) - omega**2 * mu
for use with the H-field, with wave equation
(del x (1/epsilon * del x) - omega**2 * mu) E = i * omega * M
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
pec
pec != 0
are interpreted as containing a perfect electrical
conductor (PEC). The PEC is applied per-field-component
(i.e. pec.size == epsilon.size
)
pmc
pmc != 0
are interpreted as containing a perfect magnetic
conductor (PMC). The PMC is applied per-field-component
(i.e. pmc.size == epsilon.size
)
Returns —–= Sparse matrix containing the wave operator.
m2j
def m2j(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Operator for converting a magnetic current M into an electric current
J. For use with eg. e_full()
.
Args —–= omega
: Angular frequency of
the simulation
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
mu
Returns —–= Sparse matrix for converting M to J.
poynting_e_cross
def poynting_e_cross(e: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]]) -> scipy.sparse._matrix.spmatrix
Operator for computing the Poynting vector, containing the (E x) portion of the Poynting vector.
Args —–= e
: Vectorized E-field for the
ExH cross product
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
Returns —–= Sparse matrix containing (E x) portion of Poynting cross product.
poynting_h_cross
def poynting_h_cross(h: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]]) -> scipy.sparse._matrix.spmatrix
Operator for computing the Poynting vector, containing the (H x) portion of the Poynting vector.
Args —–= h
: Vectorized H-field for the
HxE cross product
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
Returns —–= Sparse matrix containing (H x) portion of Poynting cross product.
meanas.fdfd.scpml
Functions for creating stretched coordinate perfectly matched layer (PML) absorbers.
s_function_t
Typedef for s-functions, see prepare_s_function()
prepare_s_function
def prepare_s_function(ln_R: float = -16, m: float = 4) -> collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]
Create an s_function to pass to the SCPML functions. This is used when you would like to customize the PML parameters.
Args —–= ln_R
: Natural logarithm of
the desired reflectance
m
Returns —–= An s_function, which takes an ndarray (distances) and
returns an ndarray (complex part of the cell width; needs to be divided
by sqrt(epilon_effective) * real(omega))
before use.
stretch_with_scpml
def stretch_with_scpml(dxes: list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]], axis: int, polarity: int, omega: float, epsilon_effective: float = 1.0, thickness: int = 10, s_function: collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]] | None = None) -> list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]]
Stretch dxes to contain a stretched-coordinate PML (SCPML) in one direction along one axis.
Args —–= dxes
: Grid parameters
[dx_e, dx_h]
as described in meanas.fdmath.types
axis
polarity
omega
epsilon_effective
thickness
s_function
prepare_s_function()(…)
,
allowing customization of pml parameters. Default uses prepare_s_function()
with no parameters.
Returns —–= Complex cell widths (dx_lists_mut) as discussed in
meanas.fdmath.types
.
Multiple calls to this function may be necessary if multiple absorpbing
boundaries are needed.
uniform_grid_scpml
def uniform_grid_scpml(shape: collections.abc.Sequence[int], thicknesses: collections.abc.Sequence[int], omega: float, epsilon_effective: float = 1.0, s_function: collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]], numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]] | None = None) -> list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]]
Create dx arrays for a uniform grid with a cell width of 1 and a pml.
If you want something more fine-grained, check out stretch_with_scpml()(…)
.
Args —–= shape
: Shape of the grid,
including the PMLs (which are 2*thicknesses thick)
thicknesses
[th_x, th_y, th_z]
Thickness of the PML in each direction.
Both polarities are added. Each th_ of pml is applied twice, once on
each edge of the grid along the given axis. th_*
may be
zero, in which case no pml is added.
omega
epsilon_effective
s_function
prepare_s_function()(…)
,
allowing customization of pml parameters. Default uses prepare_s_function()
with no parameters.
Returns —–= Complex cell widths (dx_lists_mut) as discussed in
meanas.fdmath.types
.
meanas.fdfd.solvers
Solvers and solver interface for FDFD problems.
generic
def generic(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], J: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, adjoint: bool = False, matrix_solver: collections.abc.Callable[..., typing.Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[typing.Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[typing.Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[typing.Union[bool, int, float, complex, str, bytes]]]] = <function _scipy_qmr>, matrix_solver_opts: dict[str, typing.Any] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]
Conjugate gradient FDFD solver using CSR sparse matrices.
All ndarray arguments should be 1D arrays, as returned by vec()
.
Args —–= omega
: Complex frequency to
solve at.
dxes
[[dx_e, dy_e, dz_e], [dx_h, dy_h, dz_h]]
(complex cell
sizes) as discussed in meanas.fdmath.types
J
epsilon
mu
pec
pmc
adjoint
matrix_solver
matrix_solver(A, b, **matrix_solver_opts) -> x
, where
A
: scipy.sparse.csr_matrix
; b
:
ArrayLike
; x
: ArrayLike
; Default
is a wrapped version of scipy.sparse.linalg.qmr()
which
doesn’t return convergence info and logs the residual every 100
iterations.
matrix_solver_opts
matrix_solver(…)
Returns —–= E-field which solves the system.
meanas.fdfd.waveguide_2d
Operators and helper functions for waveguides with unchanging cross-section.
The propagation direction is chosen to be along the z axis, and all
fields are given an implicit z-dependence of the form
exp(-1 * wavenumber * z)
.
As the z-dependence is known, all the functions in this file assume a
2D grid
(i.e. dxes = [[[dx_e[0], dx_e[1], ...], [dy_e[0], ...]], [[dx_h[0], ...], [dy_h[0], ...]]]
).
===============
Consider Maxwell’s equations in continuous space, in the frequency
domain. Assuming a structure with some (x, y) cross-section extending
uniformly into the z dimension, with a diagonal
Expanding the first two equations into vector components, we get
Substituting in our expressions for
Rewrite the last three equations as
Now apply
With a similar approach (but using
We can combine this equation for
and
However, based on our rewritten equation for
and, similarly,
By combining both pairs of expressions, we get
Using these, we can construct the eigenvalue problem
In the literature, meanas
it is allowed to be
complex.
An equivalent eigenvalue problem can be formed using the
Note that
curl_e
def curl_e(wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]]) -> scipy.sparse._matrix.spmatrix
Discretized curl operator for use with the waveguide E field.
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
Returns —–= Sparse matrix representation of the operator.
curl_h
def curl_h(wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]]) -> scipy.sparse._matrix.spmatrix
Discretized curl operator for use with the waveguide H field.
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
Returns —–= Sparse matrix representation of the operator.
e2h
def e2h(wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Returns an operator which, when applied to a vectorized E eigenfield, produces the vectorized H eigenfield.
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
mu
Returns —–= Sparse matrix representation of the operator.
e_err
def e_err(e: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> float
Calculates the relative error in the E field
Args —–= e
: Vectorized E field
wavenumber
exp(-i * wavenumber * z)
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Relative error norm(A_e @ e) / norm(e)
.
exy2e
def exy2e(wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]) -> scipy.sparse._matrix.spmatrix
Operator which transforms the vector e_xy
containing the
vectorized E_x and E_y fields, into a vectorized E containing all three
E components
From the operator derivation (see module docs), we have
as well as the intermediate equations
Combining these, we get
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
It should satisfy
operator_e() @ e_xy == wavenumber**2 * e_xy
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
Returns —–= Sparse matrix representing the operator.
exy2h
def exy2h(wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Operator which transforms the vector e_xy
containing the
vectorized E_x and E_y fields, into a vectorized H containing all three
H components
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
. It should satisfy
operator_e() @ e_xy == wavenumber**2 * e_xy
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Sparse matrix representing the operator.
get_abcd
def get_abcd(eL_xys, wavenumbers_L, eR_xys, wavenumbers_R, **kwargs)
get_s
def get_s(eL_xys, wavenumbers_L, eR_xys, wavenumbers_R, force_nogain: bool = False, force_reciprocal: bool = False, **kwargs)
get_tr
def get_tr(ehL, wavenumbers_L, ehR, wavenumbers_R, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]])
h2e
def h2e(wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]) -> scipy.sparse._matrix.spmatrix
Returns an operator which, when applied to a vectorized H eigenfield, produces the vectorized E eigenfield.
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
Returns —–= Sparse matrix representation of the operator.
h_err
def h_err(h: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> float
Calculates the relative error in the H field
Args —–= h
: Vectorized H field
wavenumber
exp(-i * wavenumber * z)
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Relative error norm(A_h @ h) / norm(h)
.
hxy2e
def hxy2e(wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Operator which transforms the vector h_xy
containing the
vectorized H_x and H_y fields, into a vectorized E containing all three
E components
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
. It should satisfy
operator_h() @ h_xy == wavenumber**2 * h_xy
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Sparse matrix representing the operator.
hxy2h
def hxy2h(wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Operator which transforms the vector h_xy
containing the
vectorized H_x and H_y fields, into a vectorized H containing all three
H components
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
. It should satisfy
operator_h() @ h_xy == wavenumber**2 * h_xy
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
mu
Returns —–= Sparse matrix representing the operator.
inner_product
def inner_product(e1: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], h2: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], prop_phase: float = 0, conj_h: bool = False, trapezoid: bool = False) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
normalized_fields_e
def normalized_fields_e(e_xy: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, prop_phase: float = 0) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Given a vector e_xy
containing the vectorized E_x and
E_y fields, returns normalized, vectorized E and H fields for the
system.
Args —–= e_xy
: Vector containing E_x
and E_y fields
wavenumber
exp(-i * wavenumber * z)
. It should satisfy
operator_e() @ e_xy == wavenumber**2 * e_xy
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
prop_phase
(dz * corrected_wavenumber)
over 1 cell in
propagation direction. Default 0 (continuous propagation direction,
i.e. dz->0).
Returns —–= (e, h)
, where each field is vectorized,
normalized, and contains all three vector components.
normalized_fields_h
def normalized_fields_h(h_xy: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, prop_phase: float = 0) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Given a vector h_xy
containing the vectorized H_x and
H_y fields, returns normalized, vectorized E and H fields for the
system.
Args —–= h_xy
: Vector containing H_x
and H_y fields
wavenumber
exp(-i * wavenumber * z)
. It should satisfy
operator_h() @ h_xy == wavenumber**2 * h_xy
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
prop_phase
(dz * corrected_wavenumber)
over 1 cell in
propagation direction. Default 0 (continuous propagation direction,
i.e. dz->0).
Returns —–= (e, h)
, where each field is vectorized,
normalized, and contains all three vector components.
operator_e
def operator_e(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Waveguide operator of the form
omega**2 * mu * epsilon +
mu * [[-Dy], [Dx]] / mu * [-Dy, Dx] +
[[Dx], [Dy]] / epsilon * [Dx, Dy] * epsilon
for use with a field vector of the form cat([E_x,
E_y])
.
More precisely, the operator is
This operator can be used to form an eigenvalue problem of the form
operator_e(...) @ [E_x, E_y] = wavenumber**2 * [E_x, E_y]
which can then be solved for the eigenmodes of the system (an
exp(-i * wavenumber * z)
z-dependence is assumed for the
fields).
Args —–= omega
: The angular frequency
of the system.
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Sparse matrix representation of the operator.
operator_h
def operator_h(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Waveguide operator of the form
omega**2 * epsilon * mu +
epsilon * [[-Dy], [Dx]] / epsilon * [-Dy, Dx] +
[[Dx], [Dy]] / mu * [Dx, Dy] * mu
for use with a field vector of the form cat([H_x,
H_y])
.
More precisely, the operator is
This operator can be used to form an eigenvalue problem of the form
operator_h(...) @ [H_x, H_y] = wavenumber**2 * [H_x, H_y]
which can then be solved for the eigenmodes of the system (an
exp(-i * wavenumber * z)
z-dependence is assumed for the
fields).
Args —–= omega
: The angular frequency
of the system.
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Sparse matrix representation of the operator.
sensitivity
def sensitivity(e_norm: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], h_norm: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]
Given a waveguide structure (dxes
, epsilon
,
mu
) and mode fields (e_norm
,
h_norm
, wavenumber
, omega
),
calculates the sensitivity of the wavenumber
The output is a vector of the same size as vec(epsilon)
,
with each element specifying the sensitivity of wavenumber
to changes in the corresponding element in vec(epsilon)
,
i.e.
An adjoint approach is used to calculate the sensitivity; the derivation is provided here:
Starting with the eigenvalue equation
where operator_e()
, and
We then multiply by
However, operator_h()
(operator_e()
(
and we can simplify to
This expression can be quickly calculated for all
Args —–= e_norm
: Normalized,
vectorized E_xyz field for the mode. E.g. as returned by normalized_fields_e()
.
h_norm
normalized_fields_e()
.
wavenumber
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Sparse matrix representation of the operator.
solve_mode
def solve_mode(mode_number: int, *args: Any, **kwargs: Any) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], complex]
Wrapper around solve_modes()
that solves for a single mode.
Args —–= mode_number
: 0-indexed mode
number to solve for
*args
solve_modes()
**kwargs
solve_modes()
Returns —–= (e_xy, wavenumber)
solve_modes
def solve_modes(mode_numbers: collections.abc.Sequence[int], omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, mode_margin: int = 2) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
Given a 2D region, attempts to solve for the eigenmode with the specified mode numbers.
Args —–= mode_numbers
: List of
0-indexed mode numbers to solve for
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
epsilon
mu
mode_margin
(max(mode_number) + mode_margin)
modes, but only return the
target mode. Increasing this value can improve the solver’s ability to
find the correct mode. Default 2.
Returns —–= e_xys
: NDArray of vfdfield_t specifying
fields. First dimension is mode number.
wavenumbers
meanas.fdfd.waveguide_3d
Tools for working with waveguide modes in 3D domains.
This module relies heavily on waveguide_2d
and mostly
just transforms its parameters into 2D equivalents and expands the
results back into 3D.
compute_overlap_e
def compute_overlap_e(E: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], axis: int, polarity: int, slices: collections.abc.Sequence[slice]) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]
Given an eigenmode obtained by solve_mode()
,
calculates an overlap_e for the mode orthogonality relation
Integrate(((E x H_mode) + (E_mode x H)) dot dn) [assumes reflection
symmetry].
TODO: add reference
Args —–= E
: E-field of the mode
H
wavenumber
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
axis
polarity
slices
epsilon[tuple(slices)]
is used to select the portion of the
grid to use as the waveguide cross-section. slices[axis] should select
only one item.
mu
Returns —–= overlap_e such that
numpy.sum(overlap_e * other_e.conj())
computes the overlap
integral
compute_source
def compute_source(E: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], axis: int, polarity: int, slices: collections.abc.Sequence[slice], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]
Given an eigenmode obtained by solve_mode()
,
returns the current source distribution necessary to position a
unidirectional source at the slice location.
Args —–= E
: E-field of the mode
wavenumber
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
axis
polarity
slices
epsilon[tuple(slices)]
is used to select the portion of the
grid to use as the waveguide cross-section. slices[axis]
should select only one item.
mu
Returns —–= J distribution for the unidirectional source
expand_e
def expand_e(E: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], axis: int, polarity: int, slices: collections.abc.Sequence[slice]) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]
Given an eigenmode obtained by solve_mode()
,
expands the E-field from the 2D slice where the mode was calculated to
the entire domain (along the propagation axis). This assumes the epsilon
cross-section remains constant throughout the entire domain; it is up to
the caller to truncate the expansion to any regions where it is
valid.
Args —–= E
: E-field of the mode
wavenumber
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
axis
polarity
slices
epsilon[tuple(slices)]
is used to select the portion of the
grid to use as the waveguide cross-section. slices[axis] should select
only one item.
Returns —–= E
, with the original field expanded along
the specified axis
.
solve_mode
def solve_mode(mode_number: int, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], axis: int, polarity: int, slices: collections.abc.Sequence[slice], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> dict[str, complex | numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]]]
Given a 3D grid, selects a slice from the grid and attempts to solve for an eigenmode propagating through that slice.
Args —–= mode_number
: Number of the
mode, 0-indexed
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
axis
polarity
slices
epsilon[tuple(slices)]
is used to select the portion of the
grid to use as the waveguide cross-section. slices[axis]
should select only one item.
epsilon
mu
Returns —–=
{
'E': NDArray[complexfloating],
'H': NDArray[complexfloating],
'wavenumber': complex,
}
meanas.fdfd.waveguide_cyl
Operators and helper functions for cylindrical waveguides with unchanging cross-section.
WORK IN PROGRESS, CURRENTLY BROKEN
As the z-dependence is known, all the functions in this file assume a
2D grid
(i.e. dxes = [[[dr_e_0, dx_e_1, ...], [dy_e_0, ...]], [[dr_h_0, ...], [dy_h_0, ...]]]
).
cylindrical_operator
def cylindrical_operator(omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], rmin: float) -> scipy.sparse._matrix.spmatrix
Cylindrical coordinate waveguide operator of the form
(NOTE: See 10.1364/OL.33.001848) TODO: consider 10.1364/OE.20.021583
TODO
for use with a field vector of the form [E_r, E_y]
.
This operator can be used to form an eigenvalue problem of the form A @ [E_r, E_y] = wavenumber**2 * [E_r, E_y]
which can then be solved for the eigenmodes of the system (an
exp(-i * wavenumber * theta)
theta-dependence is assumed
for the fields).
Args —–= omega
: The angular frequency
of the system
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
rmin
Returns —–= Sparse matrix representation of the operator
dxes2T
def dxes2T(dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], rmin=builtins.float) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]
e2h
def e2h(wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Returns an operator which, when applied to a vectorized E eigenfield, produces the vectorized H eigenfield.
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
mu
Returns —–= Sparse matrix representation of the operator.
exy2e
def exy2e(wavenumber: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]) -> scipy.sparse._matrix.spmatrix
Operator which transforms the vector e_xy
containing the
vectorized E_x and E_y fields, into a vectorized E containing all three
E components
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
It should satisfy
operator_e() @ e_xy == wavenumber**2 * e_xy
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
Returns —–= Sparse matrix representing the operator.
exy2h
def exy2h(wavenumber: complex, omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None) -> scipy.sparse._matrix.spmatrix
Operator which transforms the vector e_xy
containing the
vectorized E_x and E_y fields, into a vectorized H containing all three
H components
Args —–= wavenumber
: Wavenumber
assuming fields have z-dependence of
exp(-i * wavenumber * z)
. It should satisfy
operator_e() @ e_xy == wavenumber**2 * e_xy
omega
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
epsilon
mu
Returns —–= Sparse matrix representing the operator.
linear_wavenumbers
def linear_wavenumbers(e_xys: numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], angular_wavenumbers: Union[collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], rmin: float) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]
Calculate linear wavenumbers (1/distance) based on angular wavenumbers (1/rad) and the mode’s energy distribution.
Args —–= e_xys
: Vectorized mode fields
with shape [num_modes, 2 * x *y)
angular_wavenumbers
e_xys
epsilon
dxes
[dx_e, dx_h]
as described in meanas.fdmath.types
(2D)
rmin
Returns —–= NDArray containing the calculated linear (1/distance) wavenumbers
solve_mode
def solve_mode(mode_number: int, *args: Any, **kwargs: Any) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], complex]
Wrapper around solve_modes()
that solves for a single mode.
Args —–= mode_number
: 0-indexed mode
number to solve for
*args
solve_modes()
**kwargs
solve_modes()
Returns —–= (e_xy, angular_wavenumber)
solve_modes
def solve_modes(mode_numbers: collections.abc.Sequence[int], omega: complex, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], rmin: float, mode_margin: int = 2) -> tuple[numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]]
TODO: fixup Given a 2d (r, y) slice of epsilon, attempts to solve for the eigenmode of the bent waveguide with the specified mode number.
Args —–= mode_number
: Number of the
mode, 0-indexed
omega
dxes
epsilon
rmin
Returns —–= e_xys
: NDArray of vfdfield_t specifying
fields. First dimension is mode number.
angular_wavenumbers
meanas.fdmath
Basic discrete calculus for finite difference (fd) simulations.
Discrete fields are stored in one of two forms:
fdfield_t
form is a multidimensional
numpy.NDArray
U[m, n, p]
, where
m
, n
, and p
are discrete indices
referring to positions on the x, y, and z axes respectively.E[:, m, n, p] = [Ex[m, n, p], Ey[m, n, p], Ez[m, n, p]]
.vfdfield_t
form is simply a vectorzied (i.e. 1D)
version of the fdfield_t
, as obtained by vec()
(effectively
just numpy.ravel
)Operators which act on fields also come in two forms: + Python
functions, created by the functions in meanas.fdmath.functional
.
The generated functions act on fields in the fdfield_t
form. + Linear operators, usually 2D sparse matrices using
scipy.sparse
, created by meanas.fdmath.operators
.
These operators act on vectorized fields in the vfdfield_t
form.
The operations performed should be equivalent:
functional.op(*args)(E)
should be equivalent to
unvec(operators.op(*args) @ vec(E), E.shape[1:])
.
Generally speaking the field_t
form is easier to work
with, but can be harder or less efficient to compose (e.g. it is easy to
generate a single matrix by multiplying a series of other matrices).
This documentation and approach is roughly based on W.C. Chew’s excellent “Electromagnetic Theory on a Lattice” (doi:10.1063/1.355770), which covers a superset of this material with similar notation and more detail.
Define the discrete forward derivative as
If we treat f
as a 1D array of values, with the
i
-th value f[i]
taking up a length
dx[i]
along the x-axis, the forward derivative is
deriv_forward(f)[i] = (f[i + 1] - f[i]) / dx[i]
Likewise, discrete reverse derivative is
deriv_back(f)[i] = (f[i] - f[i - 1]) / dx[i]
The derivatives’ values are shifted by a half-cell relative to the
original function, and will have different cell widths if all the
dx[i]
(
[figure: derivatives and cell sizes]
dx0 dx1 dx2 dx3 cell sizes for function
----- ----- ----------- -----
______________________________
| | | |
f0 | f1 | f2 | f3 | function
_____|_____|___________|_____|
| | | |
| Df0 | Df1 | Df2 | Df3 forward derivative (periodic boundary)
__|_____|________|________|___
dx'3] dx'0 dx'1 dx'2 [dx'3 cell sizes for forward derivative
-- ----- -------- -------- ---
dx'0] dx'1 dx'2 dx'3 [dx'0 cell sizes for reverse derivative
______________________________
| | | |
| df1 | df2 | df3 | df0 reverse derivative (periodic boundary)
__|_____|________|________|___
Periodic boundaries are used here and elsewhere unless otherwise noted.
In the above figure, f0 =
f1 =
Df0 =
Df1 =
df0 =
The fractional subscript
Just as
For the remainder of the Discrete calculus
section, all
figures will show constant-length cells in order to focus on the vector
derivatives themselves. See the Grid description
section
below for additional information on this topic and generalization to
three dimensions.
Expanding to three dimensions, we can define two gradients
or
[code: gradients]
grad_forward(f)[i,j,k] = [Dx_forward(f)[i, j, k],
Dy_forward(f)[i, j, k],
Dz_forward(f)[i, j, k]]
= [(f[i + 1, j, k] - f[i, j, k]) / dx[i],
(f[i, j + 1, k] - f[i, j, k]) / dy[i],
(f[i, j, k + 1] - f[i, j, k]) / dz[i]]
grad_back(f)[i,j,k] = [Dx_back(f)[i, j, k],
Dy_back(f)[i, j, k],
Dz_back(f)[i, j, k]]
= [(f[i, j, k] - f[i - 1, j, k]) / dx[i],
(f[i, j, k] - f[i, j - 1, k]) / dy[i],
(f[i, j, k] - f[i, j, k - 1]) / dz[i]]
The three derivatives in the gradient cause shifts in different directions, so the x/y/z components of the resulting “vector” are defined at different points: the x-component is shifted in the x-direction, y in y, and z in z.
We call the resulting object a “fore-vector” or “back-vector”,
depending on the direction of the shift. We write it as
[figure: gradient / fore-vector]
(m, n+1, p+1) ______________ (m+1, n+1, p+1)
/: /|
/ : / |
/ : / |
(m, n, p+1)/_____________/ | The forward derivatives are defined
| : | | at the Dx, Dy, Dz points,
| :.........|...| but the forward-gradient fore-vector
z y Dz / | / is the set of all three
|/_x | Dy | / and is said to be "located" at (m,n,p)
|/ |/
(m, n, p)|_____Dx______| (m+1, n, p)
There are also two divergences,
or
[code: divergences]
div_forward(g)[i,j,k] = Dx_forward(gx)[i, j, k] +
Dy_forward(gy)[i, j, k] +
Dz_forward(gz)[i, j, k]
= (gx[i + 1, j, k] - gx[i, j, k]) / dx[i] +
(gy[i, j + 1, k] - gy[i, j, k]) / dy[i] +
(gz[i, j, k + 1] - gz[i, j, k]) / dz[i]
div_back(g)[i,j,k] = Dx_back(gx)[i, j, k] +
Dy_back(gy)[i, j, k] +
Dz_back(gz)[i, j, k]
= (gx[i, j, k] - gx[i - 1, j, k]) / dx[i] +
(gy[i, j, k] - gy[i, j - 1, k]) / dy[i] +
(gz[i, j, k] - gz[i, j, k - 1]) / dz[i]
where g = [gx, gy, gz]
is a fore- or back-vector
field.
Since we applied the forward divergence to the back-vector (and
vice-versa), the resulting scalar value is defined at the back-vector’s
(fore-vector’s) location
[figure: divergence]
^^
(m-1/2, n+1/2, p+1/2) _____||_______ (m+1/2, n+1/2, p+1/2)
/: || ,, /|
/ : || // / | The divergence at (m, n, p) (the center
/ : // / | of this cube) of a fore-vector field
(m-1/2, n-1/2, p+1/2)/_____________/ | is the sum of the outward-pointing
| : | | fore-vector components, which are
z y <==|== :.........|.====> located at the face centers.
|/_x | / | /
| / // | / Note that in a nonuniform grid, each
|/ // || |/ dimension is normalized by the cell width.
(m-1/2, n-1/2, p-1/2)|____//_______| (m+1/2, n-1/2, p-1/2)
'' ||
VV
The two curls are then
and
where
[code: curls]
curl_forward(g)[i,j,k] = [Dy_forward(gz)[i, j, k] - Dz_forward(gy)[i, j, k],
Dz_forward(gx)[i, j, k] - Dx_forward(gz)[i, j, k],
Dx_forward(gy)[i, j, k] - Dy_forward(gx)[i, j, k]]
curl_back(g)[i,j,k] = [Dy_back(gz)[i, j, k] - Dz_back(gy)[i, j, k],
Dz_back(gx)[i, j, k] - Dx_back(gz)[i, j, k],
Dx_back(gy)[i, j, k] - Dy_back(gx)[i, j, k]]
For example, consider the forward curl, at (m, n, p), of a
back-vector field g
, defined on a grid containing (m + 1/2,
n + 1/2, p + 1/2). The curl will be a fore-vector, so its z-component
will be defined at (m, n, p + 1/2). Take the nearest x- and y-components
of g
in the xy plane where the curl’s z-component is
located; these are
[curl components]
(m, n + 1/2, p + 1/2) : x-component of back-vector at (m + 1/2, n + 1/2, p + 1/2)
(m + 1, n + 1/2, p + 1/2) : x-component of back-vector at (m + 3/2, n + 1/2, p + 1/2)
(m + 1/2, n , p + 1/2) : y-component of back-vector at (m + 1/2, n + 1/2, p + 1/2)
(m + 1/2, n + 1 , p + 1/2) : y-component of back-vector at (m + 1/2, n + 3/2, p + 1/2)
These four xy-components can be used to form a loop around the curl’s z-component; its magnitude and sign is set by their loop-oriented sum (i.e. two have their signs flipped to complete the loop).
[figure: z-component of curl]
: |
z y : ^^ |
|/_x :....||.<.....| (m+1, n+1, p+1/2)
/ || /
| v || | ^
|/ |/
(m, n, p+1/2) |_____>______| (m+1, n, p+1/2)
If we discretize both space (m,n,p) and time (l), Maxwell’s equations become
with
where the spatial subscripts are abbreviated as
The above is Yee’s algorithm, written in a form analogous to Maxwell’s equations. The time derivatives can be expanded to form the update equations:
[code: Maxwell's equations updates]
H[i, j, k] -= dt * (curl_forward(E)[i, j, k] + M[t, i, j, k]) / mu[i, j, k]
E[i, j, k] += dt * (curl_back( H)[i, j, k] + J[t, i, j, k]) / epsilon[i, j, k]
Note that the E-field fore-vector and H-field back-vector are offset by a half-cell, resulting in distinct locations for all six E- and H-field components:
[figure: Field components]
(m - 1/2,=> ____________Hx__________[H] <= r + 1/2 = (m + 1/2,
n + 1/2, /: /: /| n + 1/2,
z y p + 1/2) / : / : / | p + 1/2)
|/_x / : / : / |
/ : Ez__________Hy | Locations of the E- and
/ : : : /| | H-field components for the
(m - 1/2, / : : Ey...../.|..Hz [E] fore-vector at r = (m,n,p)
n - 1/2, =>/________________________/ | /| (the large cube's center)
p + 1/2) | : : / | | / | and [H] back-vector at r + 1/2
| : :/ | |/ | (the top right corner)
| : [E].......|.Ex |
| :.................|......| <= (m + 1/2, n + 1/2, p + 1/2)
| / | /
| / | /
| / | / This is the Yee discretization
| / | / scheme ("Yee cell").
r - 1/2 = | / | /
(m - 1/2, |/ |/
n - 1/2,=> |________________________| <= (m + 1/2, n - 1/2, p - 1/2)
p - 1/2)
Each component forms its own grid, offset from the others:
[figure: E-fields for adjacent cells]
H1__________Hx0_________H0
z y /: /|
|/_x / : / | This figure shows H back-vector locations
/ : / | H0, H1, etc. and their associated components
Hy1 : Hy0 | H0 = (Hx0, Hy0, Hz0) etc.
/ : / |
/ Hz1 / Hz0
H2___________Hx3_________H3 | The equivalent drawing for E would have
| : | | fore-vectors located at the cube's
| : | | center (and the centers of adjacent cubes),
| : | | with components on the cube's faces.
| H5..........Hx4...|......H4
| / | /
Hz2 / Hz2 /
| / | /
| Hy6 | Hy4
| / | /
|/ |/
H6__________Hx7__________H7
The divergence equations can be derived by taking the divergence of
the curl equations and combining them with charge continuity,
Taking the backward curl of the
We can substitute in a time-harmonic fields
resulting in
This gives the frequency-domain wave equation,
With uniform material distribution and no sources
the frequency domain wave equation simplifies to
Since
and we get
We can convert this to three scalar-wave equations of the form
with
resulting in
with similar expressions for the y and z dimnsions (and
This implies
where
Assuming real
If
As described in the section on scalar discrete derivatives above,
cell widths (dx[i]
, dy[j]
, dz[k]
)
along each axis can be arbitrary and independently defined. Moreover,
all field components are actually defined at “derived” or “dual”
positions, in-between the “base” grid points on one or more axes.
To get a better sense of how this works, let’s start by drawing a
grid with uniform dy
and dz
and nonuniform
dx
. We will only draw one cell in the y and z dimensions to
make the illustration simpler; we need at least two cells in the x
dimension to demonstrate how nonuniform dx
affects the
various components.
Place the E fore-vectors at integer indices
Draw lines to denote the planes on which the H components and back-vectors are defined. For simplicity, don’t draw the equivalent planes for the E components and fore-vectors, except as necessary to show their locations – it’s easiest to just connect them to their associated H-equivalents.
The result looks something like this:
[figure: Component centers]
p=
[H]__________Hx___________[H]_____Hx______[H] __ +1/2
z y /: /: /: /: /| | |
|/_x / : / : / : / : / | | |
/ : / : / : / : / | | |
Hy : Ez...........Hy : Ez......Hy | | |
/: : : : /: : : : /| | | |
/ : Hz : Ey....../.:..Hz : Ey./.|..Hz __ 0 | dz[0]
/ : /: : / / : /: : / / | /| | |
/_________________________/_______________/ | / | | |
| :/ : :/ | :/ : :/ | |/ | | |
| Ex : [E].......|..Ex : [E]..|..Ex | | |
| : | : | | | |
| [H]..........Hx....|......[H].....H|x.....[H] __ --------- (n=+1/2, p=-1/2)
| / | / | / / /
Hz / Hz / Hz / / /
| / | / | / / /
| Hy | Hy | Hy __ 0 / dy[0]
| / | / | / / /
| / | / | / / /
|/ |/ |/ / /
[H]__________Hx___________[H]_____Hx______[H] __ -1/2 /
=n
|------------|------------|-------|-------|
-1/2 0 +1/2 +1 +3/2 = m
------------------------- ----------------
dx[0] dx[1]
Part of a nonuniform "base grid", with labels specifying
positions of the various field components. [E] fore-vectors
are at the cell centers, and [H] back-vectors are at the
vertices. H components along the near (-y) top (+z) edge
have been omitted to make the insides of the cubes easier
to visualize.
The above figure shows where all the components are located; however,
it is also useful to show what volumes those components correspond to.
Consider the Ex component at m = +1/2
: it is shifted in the
x-direction by a half-cell from the E fore-vector at m = 0
(labeled [E]
in the figure). It corresponds to a volume
between m = 0
and m = +1
(the other dimensions
are not shifted, i.e. they are still bounded by
n, p = +-1/2
). (See figure below). Since m
is
an index and not an x-coordinate, the Ex component is not necessarily at
the center of the volume it represents, and the x-length of its volume
is the derived quantity dx'[0] = (dx[0] + dx[1]) / 2
rather
than the base dx
. (See also Scalar derivatives and
cell shifts
).
[figure: Ex volumes]
p=
<_________________________________________> __ +1/2
z y << /: / /: >> | |
|/_x < < / : / / : > > | |
< < / : / / : > > | |
< < / : / / : > > | |
<: < / : : / : >: > | |
< : < / : : / : > : > __ 0 | dz[0]
< : < / : : / :> : > | |
<____________/____________________/_______> : > | |
< : < | : : | > : > | |
< Ex < | : Ex | > Ex > | |
< : < | : : | > : > | |
< : <....|.......:........:...|.......>...:...> __ --------- (n=+1/2, p=-1/2)
< : < | / : /| /> : > / /
< : < | / : / | / > : > / /
< :< | / :/ | / > :> / /
< < | / : | / > > _ 0 / dy[0]
< < | / | / > > / /
< < | / | / > > / /
<< |/ |/ >> / /
<____________|____________________|_______> __ -1/2 /
=n
|------------|------------|-------|-------|
-1/2 0 +1/2 +1 +3/2 = m
~------------ -------------------- -------~
dx'[-1] dx'[0] dx'[1]
The Ex values are positioned on the x-faces of the base
grid. They represent the Ex field in volumes shifted by
a half-cell in the x-dimension, as shown here. Only the
center cell (with width dx'[0]) is fully shown; the
other two are truncated (shown using >< markers).
Note that the Ex positions are the in the same positions
as the previous figure; only the cell boundaries have moved.
Also note that the points at which Ex is defined are not
necessarily centered in the volumes they represent; non-
uniform cell sizes result in off-center volumes like the
center cell here.
The next figure shows the volumes corresponding to the Hy components, which are shifted in two dimensions (x and z) compared to the base grid.
[figure: Hy volumes]
p=
z y mmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmm __ +1/2 s
|/_x << m: m: >> | |
< < m : m : > > | | dz'[1]
< < m : m : > > | |
Hy........... m........Hy...........m......Hy > | |
< < m : m : > > | |
< ______ m_____:_______________m_____:_>______ __ 0
< < m /: m / > > | |
mmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmmm > | |
< < | / : | / > > | | dz'[0]
< < | / : | / > > | |
< < | / : | / > > | |
< wwwww|w/wwwwwwwwwwwwwwwwwww|w/wwwww>wwwwwwww __ s
< < |/ w |/ w> > / /
_____________|_____________________|________ > / /
< < | w | w > > / /
< Hy........|...w........Hy.......|...w...>..Hy _ 0 / dy[0]
< < | w | w > > / /
<< | w | w > > / /
< |w |w >> / /
wwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwwww __ -1/2 /
|------------|------------|--------|-------|
-1/2 0 +1/2 +1 +3/2 = m
~------------ --------------------- -------~
dx'[-1] dx'[0] dx'[1]
The Hy values are positioned on the y-edges of the base
grid. Again here, the 'Hy' labels represent the same points
as in the basic grid figure above; the edges have shifted
by a half-cell along the x- and z-axes.
The grid lines _|:/ are edges of the area represented by
each Hy value, and the lines drawn using <m>.w represent
edges where a cell's faces extend beyond the drawn area
(i.e. where the drawing is truncated in the x- or z-
directions).
In this documentation, the E fore-vectors are placed on the base grid. An equivalent formulation could place the H back-vectors on the base grid instead. However, in the case of a non-uniform grid, the operation to get from the “base” cell widths to “derived” ones is not its own inverse.
The base grid’s cell sizes could be fully described by a list of three 1D arrays, specifying the cell widths along all three axes:
[dx, dy, dz] = [[dx[0], dx[1], ...], [dy[0], ...], [dz[0], ...]]
Note that this is a list-of-arrays rather than a 2D array, as the simulation domain may have a different number of cells along each axis.
Knowing the base grid’s cell widths and the boundary conditions
(periodic unless otherwise noted) is enough information to calculate the
cell widths dx'
, dy'
, and dz'
for
the derived grids.
However, since most operations are trivially generalized to allow either E or H to be defined on the base grid, they are written to take the a full set of base and derived cell widths, distinguished by which field they apply to rather than their “base” or “derived” status. This removes the need for each function to generate the derived widths, and makes the “base” vs “derived” distinction unnecessary in the code.
The resulting data structure containing all the cell widths takes the form of a list-of-lists-of-arrays. The first list-of-arrays provides the cell widths for the E-field fore-vectors, while the second list-of-arrays does the same for the H-field back-vectors:
[[[dx_e[0], dx_e[1], ...], [dy_e[0], ...], [dz_e[0], ...]],
[[dx_h[0], dx_h[1], ...], [dy_h[0], ...], [dz_h[0], ...]]]
where dx_e[0]
is the x-width of the m=0
cells, as used when calculating dE/dx, and dy_h[0]
is the
y-width of the n=0
cells, as used when calculating dH/dy,
etc.
Since each vector component of E and H is defined in a different
location and represents a different volume, the value of the
spatially-discrete epsilon
and mu
can also be
different for all three field components, even when representing a
simple planar interface between two isotropic materials.
As a result, epsilon
and mu
are taken to
have the same dimensions as the field, and composed of the three
diagonal tensor components:
[equations: epsilon_and_mu]
epsilon = [epsilon_xx, epsilon_yy, epsilon_zz]
mu = [mu_xx, mu_yy, mu_zz]
or
where the off-diagonal terms (e.g. epsilon_xy
) are
assumed to be zero.
High-accuracy volumetric integration of shapes on multiple grids can be performed by the gridlock module.
The values of the vacuum permittivity and permability effectively become scaling factors that appear in several locations (e.g. between the E and H fields). In order to limit floating-point inaccuracy and simplify calculations, they are often set to 1 and relative permittivities and permeabilities are used in their places; the true values can be multiplied back in after the simulation is complete if non- normalized results are needed.
meanas.fdmath.functional
Math functions for finite difference simulations
Basic discrete calculus etc.
curl_back
def curl_back(dx_h: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]] | None = None) -> collections.abc.Callable[[~TT], ~TT]
Create a function which takes the backward curl of a field.
Args —–= dx_h
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= Function f
for taking the discrete backward
curl of a field, f(H)
-> curlH
curl_back_parts
def curl_back_parts(dx_h: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]] | None = None) -> collections.abc.Callable
curl_forward
def curl_forward(dx_e: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]] | None = None) -> collections.abc.Callable[[~TT], ~TT]
Curl operator for use with the E field.
Args —–= dx_e
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= Function f
for taking the discrete forward
curl of a field, f(E)
-> curlE
curl_forward_parts
def curl_forward_parts(dx_e: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]] | None = None) -> collections.abc.Callable
deriv_back
def deriv_back(dx_h: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]] | None = None) -> tuple[collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]]
Utility operators for taking discretized derivatives (forward variant).
Args —–= dx_h
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= List of functions for taking forward derivatives along each axis.
deriv_forward
def deriv_forward(dx_e: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]] | None = None) -> tuple[collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]]
Utility operators for taking discretized derivatives (backward variant).
Args —–= dx_e
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= List of functions for taking forward derivatives along each axis.
meanas.fdmath.operators
Matrix operators for finite difference simulations
Basic discrete calculus etc.
avg_back
def avg_back(axis: int, shape: collections.abc.Sequence[int]) -> scipy.sparse._matrix.spmatrix
Backward average operator (x4 = (x4 + x3) / 2)
Args —–= axis
: Axis to average along
(x=0, y=1, z=2)
shape
Returns —–= Sparse matrix for backward average operation.
avg_forward
def avg_forward(axis: int, shape: collections.abc.Sequence[int]) -> scipy.sparse._matrix.spmatrix
Forward average operator (x4 = (x4 + x5) / 2)
Args —–= axis
: Axis to average along
(x=0, y=1, z=2)
shape
Returns —–= Sparse matrix for forward average operation.
cross
def cross(B: collections.abc.Sequence[scipy.sparse._matrix.spmatrix]) -> scipy.sparse._matrix.spmatrix
Cross product operator
Args —–= B
: List [Bx, By,
Bz]
of sparse matrices corresponding to the x, y, z portions of
the operator on the left side of the cross product.
Returns —–= Sparse matrix corresponding to (B x), where x is the cross product.
curl_back
def curl_back(dx_h: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]) -> scipy.sparse._matrix.spmatrix
Curl operator for use with the H field.
Args —–= dx_h
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= Sparse matrix for taking the discretized curl of the H-field
curl_forward
def curl_forward(dx_e: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]) -> scipy.sparse._matrix.spmatrix
Curl operator for use with the E field.
Args —–= dx_e
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= Sparse matrix for taking the discretized curl of the E-field
deriv_back
def deriv_back(dx_h: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]) -> list[scipy.sparse._matrix.spmatrix]
Utility operators for taking discretized derivatives (backward variant).
Args —–= dx_h
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= List of operators for taking forward derivatives along each axis.
deriv_forward
def deriv_forward(dx_e: collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]) -> list[scipy.sparse._matrix.spmatrix]
Utility operators for taking discretized derivatives (forward variant).
Args —–= dx_e
: Lists of cell sizes for
all axes [[dx_0, dx_1, …], [dy_0, dy_1, …], …]
.
Returns —–= List of operators for taking forward derivatives along each axis.
shift_circ
def shift_circ(axis: int, shape: collections.abc.Sequence[int], shift_distance: int = 1) -> scipy.sparse._matrix.spmatrix
Utility operator for performing a circular shift along a specified axis by a specified number of elements.
Args —–= axis
: Axis to shift along.
x=0, y=1, z=2
shape
shift_distance
Returns —–= Sparse matrix for performing the circular shift.
shift_with_mirror
def shift_with_mirror(axis: int, shape: collections.abc.Sequence[int], shift_distance: int = 1) -> scipy.sparse._matrix.spmatrix
Utility operator for performing an n-element shift along a specified axis, with mirror boundary conditions applied to the cells beyond the receding edge.
Args —–= axis
: Axis to shift along.
x=0, y=1, z=2
shape
shift_distance
Returns —–= Sparse matrix for performing the shift-with-mirror.
vec_cross
def vec_cross(b: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]) -> scipy.sparse._matrix.spmatrix
Vector cross product operator
Args —–= b
: Vector on the left side of
the cross product.
Returns —–= Sparse matrix corresponding to (b x), where x is the cross product.
meanas.fdmath.types
Types shared across multiple submodules
cfdfield_t
Complex vector field with shape (3, X, Y, Z) (e.g. [E_x, E_y,
E_z]
)
cfdfield_updater_t
Convenience type for functions which take and return an cfdfield_t
dx_lists_mut
Mutable version of dx_lists_t
dx_lists_t
‘dxes’ datastructure which contains grid cell width information in the following format:
[[[dx_e[0], dx_e[1], ...], [dy_e[0], ...], [dz_e[0], ...]],
[[dx_h[0], dx_h[1], ...], [dy_h[0], ...], [dz_h[0], ...]]]
where dx_e[0]
is the x-width of the x=0
cells, as used when calculating dE/dx, and dy_h[0]
is the
y-width of the y=0
cells, as used when calculating dH/dy,
etc.
fdfield_t
Vector field with shape (3, X, Y, Z) (e.g. [E_x, E_y,
E_z]
)
fdfield_updater_t
Convenience type for functions which take and return an fdfield_t
vcfdfield_t
Linearized complex vector field (single vector of length 3XY*Z)
vfdfield_t
Linearized vector field (single vector of length 3XY*Z)
meanas.fdmath.vectorization
Functions for moving between a vector field (list of 3 ndarrays,
[f_x, f_y, f_z]
) and a 1D array representation of that
field [f_x0, f_x1, f_x2,… f_y0,… f_z0,…]
. Vectorized
versions of the field use row-major (ie., C-style) ordering.
unvec
def unvec(v: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]] | None, shape: collections.abc.Sequence[int], nvdim: int = 3) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]] | None
Perform the inverse of vec(): take a 1D ndarray and output an
nvdim
-component field of form e.g. [f_x, f_y,
f_z]
(nvdim=3
) where each of f_*
is a
len(shape)-dimensional ndarray.
Returns None
if called with v=None
.
Args —–= v
: 1D ndarray representing a
vector field of shape shape (or None)
shape
nvdim
Returns —–= [f_x, f_y, f_z]
where each f_
is a len(shape)
dimensional ndarray (or
None
)
vec
def vec(f: Union[numpy.ndarray[Any, numpy.dtype[numpy.floating]], numpy.ndarray[Any, numpy.dtype[numpy.complexfloating]], collections.abc.Buffer, numpy._typing._array_like._SupportsArray[numpy.dtype[Any]], numpy._typing._nested_sequence._NestedSequence[numpy._typing._array_like._SupportsArray[numpy.dtype[Any]]], bool, int, float, complex, str, bytes, numpy._typing._nested_sequence._NestedSequence[Union[bool, int, float, complex, str, bytes]], ForwardRef(None)]) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | numpy.ndarray[typing.Any, numpy.dtype[numpy.complexfloating]] | None
Create a 1D ndarray from a vector field which spans a 1-3D region.
Returns None
if called with f=None
.
Args —–= f
: A vector field,
e.g. [f_x, f_y, f_z]
where each f_
component
is a 1- to 3-D ndarray (f_*
should all be the same size).
Doesn’t fail with f=None
.
Returns —–= 1D ndarray containing the linearized field (or
None
)
meanas.fdtd
Utilities for running finite-difference time-domain (FDTD) simulations
See the discussion of Maxwell's Equations
in meanas.fdmath
for basic mathematical
background.
From the discussion of “Plane waves and the Dispersion relation” in
meanas.fdmath
, we have
or, if
Based on this, we can set
dt = sqrt(mu.min() * epsilon.min()) / sqrt(1/dx_min**2 + 1/dy_min**2 + 1/dz_min**2)
The dx_min
, dy_min
, dz_min
should be the minimum value across both the base and derived grids.
Let
where
By taking the divergence and rearranging terms, we can show that
where in the last line the spatial subscripts have been dropped to
emphasize the time subscripts
etc. For
and for
These two results form the discrete time-domain analogue to Poynting’s theorem. They hint at the expressions for the energy, which can be calculated at the same time-index as either the E or H field:
Rewriting the Poynting theorem in terms of the energy expressions,
This result is exact and should practically hold to within numerical precision. No time- or spatial-averaging is necessary.
Note that each value of
It is often useful to excite the simulation with an arbitrary broadband pulse and then extract the frequency-domain response by performing an on-the-fly Fourier transform of the time-domain fields.
The Ricker wavelet (normalized second derivative of a Gaussian) is commonly used for the pulse shape. It can be written
with
meanas.fdtd.base
Basic FDTD field updates
maxwell_e
def maxwell_e(dt: float, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]
Build a function which performs a portion the time-domain E-field update,
E += curl_back(H[t]) / epsilon
The full update should be
E += (curl_back(H[t]) + J) / epsilon
which requires an additional step of E += J / epsilon
which is not performed by the generated function.
See meanas.fdmath
for
descriptions of
dxes
: “Datastructure: dx_lists_t” sectionepsilon
: “Permittivity and Permeability” sectionAlso see the “Timestep” section of meanas.fdtd
for a discussion of the
dt
parameter.
Args —–= dt
: Timestep. See meanas.fdtd
for details.
dxes
meanas.fdmath
.
Returns —–= Function
f(E_old, H_old, epsilon) -> E_new
.
maxwell_h
def maxwell_h(dt: float, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]
Build a function which performs part of the time-domain H-field update,
H -= curl_forward(E[t]) / mu
The full update should be
H -= (curl_forward(E[t]) + M) / mu
which requires an additional step of H -= M / mu
which
is not performed by the generated function; this step can be omitted if
there is no magnetic current M
.
See meanas.fdmath
for
descriptions of
dxes
: “Datastructure: dx_lists_t” sectionmu
: “Permittivity and Permeability” sectionAlso see the “Timestep” section of meanas.fdtd
for a discussion of the
dt
parameter.
Args —–= dt
: Timestep. See meanas.fdtd
for details.
dxes
meanas.fdmath
.
Returns —–= Function
f(E_old, H_old, epsilon) -> E_new
.
meanas.fdtd.boundaries
Boundary conditions
#TODO conducting boundary documentation
conducting_boundary
def conducting_boundary(direction: int, polarity: int) -> tuple[collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], collections.abc.Callable[..., numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]]]
meanas.fdtd.energy
delta_energy_e2h
def delta_energy_e2h(dt: float, h0: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], e1: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], h2: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], e3: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Change in energy during the half-step from e1
to
h2
.
This is just from (h2 * h2 + e3 * e1) - (e1 * e1 + h0 * h2)
Args —–= h0
: E-field one half-timestep
before the start of the energy delta.
e1
h2
e1
).
e3
epsilon
mu
dxes
meanas.fdmath
.
Returns —–= Change in energy from the time of e1
to the
time of h2
.
delta_energy_h2e
def delta_energy_h2e(dt: float, e0: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], h1: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], e2: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], h3: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Change in energy during the half-step from h1
to
e2
.
This is just from (e2 * e2 + h3 * h1) - (h1 * h1 + e0 * e2)
Args —–= e0
: E-field one half-timestep
before the start of the energy delta.
h1
e2
h1
).
h3
epsilon
mu
dxes
meanas.fdmath
.
Returns —–= Change in energy from the time of h1
to the
time of e2
.
delta_energy_j
def delta_energy_j(j0: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], e1: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Calculate
Note that each value of
dxmul
def dxmul(ee: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], hh: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | float | None = None, mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | float | None = None, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
energy_estep
def energy_estep(h0: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], e1: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], h2: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Calculate energy U
at the time of the provided E-field
e1
.
TODO: Figure out what this means spatially.
Args —–= h0
: H-field one half-timestep
before the energy.
e1
h2
epsilon
mu
dxes
meanas.fdmath
.
Returns —–= Energy, at the time of the E-field e1
.
energy_hstep
def energy_hstep(e0: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], h1: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], e2: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, mu: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Calculate energy U
at the time of the provided H-field
h1
.
TODO: Figure out what this means spatially.
Args —–= e0
: E-field one half-timestep
before the energy.
h1
e2
epsilon
mu
dxes
meanas.fdmath
.
Returns —–= Energy, at the time of the H-field h1
.
poynting
def poynting(e: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], h: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Calculate the poynting vector S
(
This is the energy transfer rate (amount of energy U
per
dt
transferred between adjacent cells) in each direction
that happens during the half-step bounded by the two provided
fields.
The returned vector field S
is the energy flow across
+x, +y, and +z boundaries of the corresponding U
cell. For
example,
mx = numpy.roll(mask, -1, axis=0)
my = numpy.roll(mask, -1, axis=1)
mz = numpy.roll(mask, -1, axis=2)
u_hstep = fdtd.energy_hstep(e0=es[ii - 1], h1=hs[ii], e2=es[ii], **args)
u_estep = fdtd.energy_estep(h0=hs[ii], e1=es[ii], h2=hs[ii + 1], **args)
delta_j_B = fdtd.delta_energy_j(j0=js[ii], e1=es[ii], dxes=dxes)
du_half_h2e = u_estep - u_hstep - delta_j_B
s_h2e = -fdtd.poynting(e=es[ii], h=hs[ii], dxes=dxes) * dt
planes = [s_h2e[0, mask].sum(), -s_h2e[0, mx].sum(),
s_h2e[1, mask].sum(), -s_h2e[1, my].sum(),
s_h2e[2, mask].sum(), -s_h2e[2, mz].sum()]
assert_close(sum(planes), du_half_h2e[mask])
(see meanas.tests.test_fdtd.test_poynting_planes
)
The full relationship is
These equalities are exact and should practically hold to within
numerical precision. No time- or spatial-averaging is necessary. (See
meanas.fdtd
docs for
derivation.)
Args —–= e
: E-field
h
e
)
dxes
meanas.fdmath
.
Returns —–= s
: Vector field. Components indicate the
energy transfer rate from the corresponding energy cell into its +x, +y,
and +z neighbors during the half-step from the time of the earlier input
field until the time of later input field.
poynting_divergence
def poynting_divergence(s: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, *, e: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, h: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]] | None = None, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]] | None = None) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]
Calculate the divergence of the poynting vector.
This is the net energy flow for each cell, i.e. the change in energy
U
per dt
caused by transfer of energy to
nearby cells (rather than absorption/emission by currents J
or M
).
See poynting()
and meanas.fdtd
for more details.
Args —–= s
: Poynting vector, as
calculated with poynting()
. Optional;
caller can provide e
and h
instead.
e
s
or both e
and
h
)
h
s
or both e
and
h
)
dxes
meanas.fdmath
.
Returns —–= ds
: Divergence of the poynting vector.
Entries indicate the net energy flow out of the corresponding energy
cell.
meanas.fdtd.pml
PML implementations
#TODO discussion of PMLs #TODO cpml documentation
cpml_params
def cpml_params(axis: int, polarity: int, dt: float, thickness: int = 8, ln_R_per_layer: float = -1.6, epsilon_eff: float = 1, mu_eff: float = 1, m: float = 3.5, ma: float = 1, cfs_alpha: float = 0) -> dict[str, typing.Any]
updates_with_cpml
def updates_with_cpml(cpml_params: collections.abc.Sequence[collections.abc.Sequence[dict[str, typing.Any] | None]], dt: float, dxes: collections.abc.Sequence[collections.abc.Sequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], *, dtype: Union[numpy.dtype[Any], ForwardRef(None), type[Any], numpy._typing._dtype_like._SupportsDType[numpy.dtype[Any]], str, tuple[Any, int], tuple[Any, Union[SupportsIndex, collections.abc.Sequence[SupportsIndex]]], list[Any], numpy._typing._dtype_like._DTypeDict, tuple[Any, Any]] = numpy.float32) -> tuple[collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], None], collections.abc.Callable[[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]], numpy.ndarray[typing.Any, numpy.dtype[numpy.floating]]], None]]
meanas.test
Tests (run with
python3 -m pytest -rxPXs | tee results.txt
)
meanas.test.conftest
Test fixtures
dx
def dx(request: Any) -> float
dxes
def dxes(request: Any, shape: tuple[int, ...], dx: float) -> list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]]
epsilon
def epsilon(request: Any, shape: tuple[int, ...], epsilon_bg: float, epsilon_fg: float) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]
epsilon_bg
def epsilon_bg(request: Any) -> float
epsilon_fg
def epsilon_fg(request: Any) -> float
j_mag
def j_mag(request: Any) -> float
shape
def shape(request: Any) -> tuple[int, ...]
meanas.test.test_fdfd
j_distribution
def j_distribution(request: Any, shape: tuple[int, ...], j_mag: float) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]
omega
def omega(request: Any) -> float
pec
def pec(request: Any) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None
pmc
def pmc(request: Any) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None
sim
def sim(request: Any, shape: tuple[int, ...], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], dxes: list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]], j_distribution: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], omega: float, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None) -> meanas.test.test_fdfd.FDResult
Build simulation from parts
test_poynting_planes
def test_poynting_planes(sim: FDResult) -> None
test_residual
def test_residual(sim: FDResult) -> None
FDResult
class FDResult(shape: tuple[int, ...], dxes: list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], omega: complex, j: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], e: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None)
FDResult(shape: tuple[int, …], dxes: list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], omega: complex, j: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], e: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None)
dxes
e
epsilon
j
omega
pec
pmc
shape
meanas.test.test_fdfd_pml
dxes
def dxes(request: Any, shape: tuple[int, ...], dx: float, omega: float, epsilon_fg: float) -> list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]]
epsilon
def epsilon(request: Any, shape: tuple[int, ...], epsilon_bg: float, epsilon_fg: float) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]
j_distribution
def j_distribution(request: Any, shape: tuple[int, ...], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], dxes: collections.abc.MutableSequence[collections.abc.MutableSequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], omega: float, src_polarity: int) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]]
omega
def omega(request: Any) -> float
pec
def pec(request: Any) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None
pmc
def pmc(request: Any) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None
shape
def shape(request: Any) -> tuple[int, int, int]
sim
def sim(request: Any, shape: tuple[int, ...], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], dxes: collections.abc.MutableSequence[collections.abc.MutableSequence[numpy.ndarray[typing.Any, numpy.dtype[numpy.floating | numpy.complexfloating]]]], j_distribution: numpy.ndarray[typing.Any, numpy.dtype[numpy.complex128]], omega: float, pec: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None, pmc: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]] | None) -> meanas.test.test_fdfd.FDResult
src_polarity
def src_polarity(request: Any) -> int
test_pml
def test_pml(sim: meanas.test.test_fdfd.FDResult, src_polarity: int) -> None
meanas.test.test_fdtd
dt
def dt(request: Any) -> float
j_distribution
def j_distribution(request: Any, shape: tuple[int, ...], j_mag: float) -> numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]
j_steps
def j_steps(request: Any) -> tuple[int, ...]
sim
def sim(request: Any, shape: tuple[int, ...], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], dxes: list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]], dt: float, j_distribution: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], j_steps: tuple[int, ...]) -> meanas.test.test_fdtd.TDResult
test_energy_conservation
def test_energy_conservation(sim: TDResult) -> None
Assumes fields start at 0 before J0 is added
test_initial_energy
def test_initial_energy(sim: TDResult) -> None
Assumes fields start at 0 before J0 is added
test_initial_fields
def test_initial_fields(sim: TDResult) -> None
test_poynting_divergence
def test_poynting_divergence(sim: TDResult) -> None
test_poynting_planes
def test_poynting_planes(sim: TDResult) -> None
TDResult
class TDResult(shape: tuple[int, ...], dt: float, dxes: list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]], epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], j_distribution: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]], j_steps: tuple[int, ...], es: list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]] = <factory>, hs: list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]] = <factory>, js: list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]] = <factory>)
TDResult(shape: tuple[int, …], dt: float, dxes:
list[list[numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]]]],
epsilon: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]],
j_distribution: numpy.ndarray[typing.Any, numpy.dtype[numpy.float64]],
j_steps: tuple[int, …], es: list[numpy.ndarray[typing.Any,
numpy.dtype[numpy.float64]]] =
dt
dxes
epsilon
es
hs
j_distribution
j_steps
js
shape
meanas.test.utils
assert_close
def assert_close(x: numpy.ndarray[typing.Any, numpy.dtype[+_ScalarType_co]], y: numpy.ndarray[typing.Any, numpy.dtype[+_ScalarType_co]], *args, **kwargs) -> None
assert_fields_close
def assert_fields_close(x: numpy.ndarray[typing.Any, numpy.dtype[+_ScalarType_co]], y: numpy.ndarray[typing.Any, numpy.dtype[+_ScalarType_co]], *args, **kwargs) -> None
Generated by pdoc 0.11.1 (https://pdoc3.github.io).