fdfd¶
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:
operators,functional: General FDFD problem setup.solvers: Solver interface and reference implementation.scpml: Stretched-coordinate perfectly matched layer (SCPML) boundary conditions.waveguide_2d: Operators and mode-solver for waveguides with constant cross-section.waveguide_3d: Functions for transformingwaveguide_2dresults into 3D, including mode-source and overlap-window construction.farfield,bloch,eme: specialized helper modules for near/far transforms, Bloch-periodic problems, and eigenmode expansion.
================================================================
From the "Frequency domain" section of meanas.fdmath, we have
resulting in
Maxwell's equations are then
With \(\Delta_t \to 0\), this simplifies to
and then
Core operator layers¶
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)
e_full
¶
e_full(
omega: complex,
dxes: dx_lists_t,
epsilon: fdfield,
mu: fdfield | None = None,
) -> cfdfield_updater_t
Wave operator for use with E-field. See operators.e_full for details.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
fdfield
|
Dielectric constant |
required |
mu
|
fdfield | None
|
Magnetic permeability (default 1 everywhere) |
None
|
Returns:
| Type | Description |
|---|---|
cfdfield_updater_t
|
Function |
cfdfield_updater_t
|
|
eh_full
¶
eh_full(
omega: complex,
dxes: dx_lists_t,
epsilon: fdfield,
mu: fdfield | None = None,
) -> Callable[
[cfdfield, cfdfield], tuple[cfdfield_t, cfdfield_t]
]
Wave operator for full (both E and H) field representation.
See operators.eh_full.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
fdfield
|
Dielectric constant |
required |
mu
|
fdfield | None
|
Magnetic permeability (default 1 everywhere) |
None
|
Returns:
| Type | Description |
|---|---|
Callable[[cfdfield, cfdfield], tuple[cfdfield_t, cfdfield_t]]
|
Function |
Callable[[cfdfield, cfdfield], tuple[cfdfield_t, cfdfield_t]]
|
|
e2h
¶
Utility operator for converting the E field into the H field.
For use with e_full -- assumes that there is no magnetic current M.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
mu
|
fdfield | None
|
Magnetic permeability (default 1 everywhere) |
None
|
Returns:
| Type | Description |
|---|---|
cfdfield_updater_t
|
Function |
cfdfield_updater_t
|
|
m2j
¶
Utility operator for converting magnetic current M distribution
into equivalent electric current distribution J.
For use with e.g. e_full.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
mu
|
fdfield | None
|
Magnetic permeability (default 1 everywhere) |
None
|
Returns:
| Type | Description |
|---|---|
cfdfield_updater_t
|
Function |
cfdfield_updater_t
|
|
e_tfsf_source
¶
e_tfsf_source(
TF_region: fdfield,
omega: complex,
dxes: dx_lists_t,
epsilon: fdfield,
mu: fdfield | None = None,
) -> cfdfield_updater_t
Operator that turns an E-field distribution into a total-field/scattered-field (TFSF) source.
If A is the full wave operator from e_full(...) and Q is the diagonal
mask selecting the total-field region, then the TFSF source is the commutator
This vanishes in the interior of the total-field and scattered-field regions
and is supported only at their shared boundary, where the mask discontinuity
makes A and Q fail to commute. The returned current is therefore the
distributed source needed to inject the desired total field without also
forcing the scattered-field region.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
TF_region
|
fdfield
|
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. |
required |
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
fdfield
|
Dielectric constant distribution |
required |
mu
|
fdfield | None
|
Magnetic permeability (default 1 everywhere) |
None
|
Returns:
| Type | Description |
|---|---|
cfdfield_updater_t
|
Function |
cfdfield_updater_t
|
|
poynting_e_cross_h
¶
Generates a function that takes the single-frequency E and H fields
and calculates the cross product E x H = \(E \times H\) as required
for the Poynting vector, \(S = E \times H\).
On the Yee grid, the electric and magnetic components are not stored at the
same locations. This helper therefore applies the same one-cell electric-field
shifts used by the sparse operators.poynting_e_cross(...) construction so
that the discrete cross product matches the face-centered energy flux used in
meanas.fdtd.energy.poynting(...).
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.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
dxes
|
dx_lists_t
|
Grid parameters |
required |
Returns:
| Type | Description |
|---|---|
Callable[[cfdfield, cfdfield], cfdfield_t]
|
Function |
Callable[[cfdfield, cfdfield], cfdfield_t]
|
For time-average power, call it as |
meanas.fdfd.operators
¶
Sparse matrix operators for use with electromagnetic wave equations.
These functions return sparse-matrix (scipy.sparse.sparray) representations of
a variety of operators, intended for use with E and H fields vectorized using the
meanas.fdmath.vectorization.vec() and meanas.fdmath.vectorization.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:
- E-only wave operator
- H-only wave operator
- EH wave operator
- Curl for use with E, H fields
- E to H conversion
- M to J conversion
- Poynting cross products
- Circular shifts
- Discrete derivatives
- Averaging operators
- Cross product matrices
e_full
¶
e_full(
omega: complex,
dxes: dx_lists_t,
epsilon: vfdfield | vcfdfield,
mu: vfdfield | None = None,
pec: vfdfield | None = None,
pmc: vfdfield | None = None,
) -> sparse.sparray
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().
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
vfdfield | vcfdfield
|
Vectorized dielectric constant |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere). |
None
|
pec
|
vfdfield | None
|
Vectorized mask specifying PEC cells. Any cells where |
None
|
pmc
|
vfdfield | None
|
Vectorized mask specifying PMC cells. Any cells where |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix containing the wave operator. |
e_full_preconditioners
¶
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
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
dxes
|
dx_lists_t
|
Grid parameters |
required |
Returns:
| Type | Description |
|---|---|
tuple[sparray, sparray]
|
Preconditioner matrices |
h_full
¶
h_full(
omega: complex,
dxes: dx_lists_t,
epsilon: vfdfield,
mu: vfdfield | None = None,
pec: vfdfield | None = None,
pmc: vfdfield | None = None,
) -> sparse.sparray
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
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
vfdfield
|
Vectorized dielectric constant |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere) |
None
|
pec
|
vfdfield | None
|
Vectorized mask specifying PEC cells. Any cells where |
None
|
pmc
|
vfdfield | None
|
Vectorized mask specifying PMC cells. Any cells where |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix containing the wave operator. |
eh_full
¶
eh_full(
omega: complex,
dxes: dx_lists_t,
epsilon: vfdfield,
mu: vfdfield | None = None,
pec: vfdfield | None = None,
pmc: vfdfield | None = None,
) -> sparse.sparray
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)):
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
vfdfield
|
Vectorized dielectric constant |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere) |
None
|
pec
|
vfdfield | None
|
Vectorized mask specifying PEC cells. Any cells where |
None
|
pmc
|
vfdfield | None
|
Vectorized mask specifying PMC cells. Any cells where |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix containing the wave operator. |
e2h
¶
e2h(
omega: complex,
dxes: dx_lists_t,
mu: vfdfield | None = None,
pmc: vfdfield | None = None,
) -> sparse.sparray
Utility operator for converting the E field into the H field.
For use with e_full() -- assumes that there is no magnetic current M.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere) |
None
|
pmc
|
vfdfield | None
|
Vectorized mask specifying PMC cells. Any cells where |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix for converting E to H. |
m2j
¶
Operator for converting a magnetic current M into an electric current J.
For use with eg. e_full().
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere) |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix for converting M to J. |
poynting_e_cross
¶
Operator for computing the staggered-grid (E \times) part of the Poynting vector.
On the Yee grid the E and H components live on different edges, so the
electric field must be shifted by one cell in the appropriate direction
before forming the discrete cross product. This sparse operator encodes that
shifted cross product directly and is the matrix equivalent of
functional.poynting_e_cross_h(...).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
e
|
vcfdfield
|
Vectorized E-field for the ExH cross product |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix containing the |
sparray
|
cross product. |
poynting_h_cross
¶
Operator for computing the staggered-grid (H \times) part of the Poynting vector.
Together with poynting_e_cross(...), this provides the matrix form of the
Yee-grid cross product used in the flux helpers. The two are related by the
usual antisymmetry of the cross product,
once the same staggered field placement is used on both sides.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
h
|
vcfdfield
|
Vectorized H-field for the HxE cross product |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix containing the |
sparray
|
cross product. |
e_tfsf_source
¶
e_tfsf_source(
TF_region: vfdfield,
omega: complex,
dxes: dx_lists_t,
epsilon: vfdfield,
mu: vfdfield | None = None,
) -> sparse.sparray
Operator that turns a desired E-field distribution into a total-field/scattered-field (TFSF) source.
Let A be the full wave operator from e_full(...), and let
Q = \mathrm{diag}(TF_region) be the projector onto the total-field region.
Then the TFSF current operator is the commutator
Inside regions where Q is locally constant, A and Q commute and the
source vanishes. Only cells at the TF/SF boundary contribute nonzero current,
which is exactly the desired distributed source for injecting the chosen
field into the total-field region without directly forcing the
scattered-field region.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
TF_region
|
vfdfield
|
Mask, which is set to 1 inside the total-field region and 0 in the scattered-field region |
required |
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
vfdfield
|
Vectorized dielectric constant |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere). |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix that turns an E-field into a current (J) distribution. |
e_boundary_source
¶
e_boundary_source(
mask: vfdfield,
omega: complex,
dxes: dx_lists_t,
epsilon: vfdfield,
mu: vfdfield | None = None,
periodic_mask_edges: bool = False,
) -> sparse.sparray
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.
Equivalently, this helper first constructs the TFSF commutator source for the
full mask and then zeroes out all cells except the mask boundary. The
boundary is defined as the set of cells whose mask value changes under a
one-cell shift in any Cartesian direction. With periodic_mask_edges=False
the shifts mirror at the domain boundary; with True they wrap periodically.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
mask
|
vfdfield
|
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. |
required |
omega
|
complex
|
Angular frequency of the simulation |
required |
dxes
|
dx_lists_t
|
Grid parameters |
required |
epsilon
|
vfdfield
|
Vectorized dielectric constant |
required |
mu
|
vfdfield | None
|
Vectorized magnetic permeability (default 1 everywhere). |
None
|
Returns:
| Type | Description |
|---|---|
sparray
|
Sparse matrix that turns an E-field into a current (J) distribution. |
meanas.fdfd.solvers
¶
Solvers and solver interface for FDFD problems.
generic
¶
generic(
omega: complex,
dxes: dx_lists_t,
J: vcfdfield,
epsilon: vfdfield,
mu: vfdfield | None = None,
*,
pec: vfdfield | None = None,
pmc: vfdfield | None = None,
adjoint: bool = False,
matrix_solver: Callable[..., ArrayLike] = _scipy_qmr,
matrix_solver_opts: dict[str, Any] | None = None,
E_guess: vcfdfield | None = None,
) -> vcfdfield_t
Conjugate gradient FDFD solver using CSR sparse matrices.
All ndarray arguments should be 1D arrays, as returned by meanas.fdmath.vectorization.vec().
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
complex
|
Complex frequency to solve at. |
required |
dxes
|
dx_lists_t
|
|
required |
J
|
vcfdfield
|
Electric current distribution (at E-field locations) |
required |
epsilon
|
vfdfield
|
Dielectric constant distribution (at E-field locations) |
required |
mu
|
vfdfield | None
|
Magnetic permeability distribution (at H-field locations) |
None
|
pec
|
vfdfield | None
|
Perfect electric conductor distribution (at E-field locations; non-zero value indicates PEC is present) |
None
|
pmc
|
vfdfield | None
|
Perfect magnetic conductor distribution (at H-field locations; non-zero value indicates PMC is present) |
None
|
adjoint
|
bool
|
If true, solves the adjoint problem. |
False
|
matrix_solver
|
Callable[..., ArrayLike]
|
Called as |
_scipy_qmr
|
matrix_solver_opts
|
dict[str, Any] | None
|
Passed as kwargs to |
None
|
E_guess
|
vcfdfield | None
|
Guess at the solution E-field. |
None
|
Returns:
| Type | Description |
|---|---|
vcfdfield_t
|
E-field which solves the system. |
meanas.fdfd.scpml
¶
Functions for creating stretched coordinate perfectly matched layer (PML) absorbers.
s_function_t
module-attribute
¶
Typedef for s-functions, see prepare_s_function()
prepare_s_function
¶
Create an s_function to pass to the SCPML functions. This is used when you would like to customize the PML parameters.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
ln_R
|
float
|
Natural logarithm of the desired reflectance |
-16
|
m
|
float
|
Polynomial order for the PML (imaginary part increases as distance ** m) |
4
|
Returns:
| Type | Description |
|---|---|
s_function_t
|
An s_function, which takes an ndarray (distances) and returns an ndarray (complex part |
s_function_t
|
of the cell width; needs to be divided by |
s_function_t
|
before use. |
uniform_grid_scpml
¶
uniform_grid_scpml(
shape: Sequence[int],
thicknesses: Sequence[int],
omega: float,
epsilon_effective: float = 1.0,
s_function: s_function_t | None = None,
) -> list[list[NDArray[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(...).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
shape
|
Sequence[int]
|
Shape of the grid, including the PMLs (which are 2*thicknesses thick) |
required |
thicknesses
|
Sequence[int]
|
|
required |
omega
|
float
|
Angular frequency for the simulation |
required |
epsilon_effective
|
float
|
Effective epsilon of the PML. Match this to the material at the edge of your grid. Default 1. |
1.0
|
s_function
|
s_function_t | None
|
created by |
None
|
Returns:
| Type | Description |
|---|---|
list[list[NDArray[float64]]]
|
Complex cell widths (dx_lists_mut) as discussed in |
stretch_with_scpml
¶
stretch_with_scpml(
dxes: list[list[NDArray[float64]]],
axis: int,
polarity: int,
omega: float,
epsilon_effective: float = 1.0,
thickness: int = 10,
s_function: s_function_t | None = None,
) -> list[list[NDArray[numpy.float64]]]
Stretch dxes to contain a stretched-coordinate PML (SCPML) in one direction along one axis.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
dxes
|
list[list[NDArray[float64]]]
|
Grid parameters |
required |
axis
|
int
|
axis to stretch (0=x, 1=y, 2=z) |
required |
polarity
|
int
|
direction to stretch (-1 for -ve, +1 for +ve) |
required |
omega
|
float
|
Angular frequency for the simulation |
required |
epsilon_effective
|
float
|
Effective epsilon of the PML. Match this to the material at the edge of your grid. Default 1. |
1.0
|
thickness
|
int
|
number of cells to use for pml (default 10) |
10
|
s_function
|
s_function_t | None
|
Created by |
None
|
Returns:
| Type | Description |
|---|---|
list[list[NDArray[float64]]]
|
Complex cell widths (dx_lists_mut) as discussed in |
list[list[NDArray[float64]]]
|
Multiple calls to this function may be necessary if multiple absorpbing boundaries are needed. |
meanas.fdfd.farfield
¶
Functions for performing near-to-farfield transformation (and the reverse).
near_to_farfield
¶
near_to_farfield(
E_near: cfdfield_t,
H_near: cfdfield_t,
dx: float,
dy: float,
padded_size: list[int] | int | None = None,
) -> dict[str, 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.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
E_near
|
cfdfield_t
|
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). |
required |
H_near
|
cfdfield_t
|
List of 2 ndarrays containing the 2D phasor field slices for the transverse H fields (e.g. [Hx, hy] for calculating the farfield towrad the z-direction). |
required |
dx
|
float
|
Cell size along x-dimension, in units of wavelength. |
required |
dy
|
float
|
Cell size along y-dimension, in units of wavelength. |
required |
padded_size
|
list[int] | int | None
|
Shape of the output. A single integer |
None
|
Returns:
| Type | Description |
|---|---|
dict[str, Any]
|
Dict with keys |
dict[str, Any]
|
|
dict[str, Any]
|
|
dict[str, Any]
|
|
dict[str, Any]
|
|
dict[str, Any]
|
|
dict[str, Any]
|
|
far_to_nearfield
¶
far_to_nearfield(
E_far: cfdfield_t,
H_far: cfdfield_t,
dkx: float,
dky: float,
padded_size: list[int] | int | None = None,
) -> dict[str, 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.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
E_far
|
cfdfield_t
|
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)) |
required |
H_far
|
cfdfield_t
|
List of 2 ndarrays containing the 2D phasor field slices for the transverse H fields (e.g. [Hx, hy] for calculating the nearfield toward the z-direction). Fields should be normalized so that H_far = H_far_actual / (i k exp(-i k r) / (4 pi r)) |
required |
dkx
|
float
|
kx discretization, in units of wavelength. |
required |
dky
|
float
|
ky discretization, in units of wavelength. |
required |
padded_size
|
list[int] | int | None
|
Shape of the output. A single integer |
None
|
Returns:
| Type | Description |
|---|---|
dict[str, Any]
|
Dict with keys |
dict[str, Any]
|
|
dict[str, Any]
|
|
dict[str, Any]
|
|