Meep contains a feature to decompose arbitrary fields into a superposition of the harmonic modes of a given structure via its integration with the eigenmode solver MPB. This section provides an overview of the theory and implementation of this feature. For examples, see Tutorial/Mode Decomposition.
[TOC]
The theory underlying mode decomposition is described in Chapter 31 ("Modal methods for Maxwell's equations") of Optical Waveguide Theory by Snyder and Love.
Consider a waveguide with propagation axis along the
Any arbitrary fields
-
In Meep, compute the Fourier-transformed transverse fields
$\psi$ on a surface that is transverse to the waveguide, stored in adft_flux
object. -
In MPB, compute the eigenmodes
$\psi^\pm_n$ as well as the propagation constants$\beta_n$ for the same cross-sectional structure. -
Compute the coefficients
$\alpha_n^\pm$ for any number of eigenmodes $n=1,2,...$.
This is all done automatically in Meep using the get_eigenmode_coefficients
routine.
Meep normalizes the modes
The mode-decomposition feature is available via the meep::fields::get_eigenmode_coefficients
function callable from Python or C++. This function makes use of several lower-level functions which are described in more detail below. The C++ header for this function is:
void fields::get_eigenmode_coefficients(dft_flux flux,
const volume &eig_vol,
int *bands,
int num_bands,
int parity,
double eig_resolution,
double eigensolver_tol,
std::complex<double> *coeffs,
double *vgrp,
kpoint_func user_kpoint_func,
void *user_kpoint_data,
vec *kpoint_list,
vec *kdom_list,
double *cscale,
direction d,
diffractedplanewave *dp)
The following are the parameters:
-
flux
is adft_flux
object containing the frequency-domain fields on a cross-sectional slice perpendicular to the waveguide axis. -
eig_vol
is thevolume
passed to MPB for the eigenmode calculation. In most cases this will simply be the volume over which the frequency-domain fields are tabulated (i.e.flux.where
). -
bands
is an array of integers corresponding to the mode indices (equivalent to$n$ in the two formulas above). -
num_bands
is the length of thebands
array. -
parity
is the parity of the mode to calculate, assuming the structure has$z$ and/or$y$ mirror symmetry in the source region. If the structure has both$y$ and$z$ mirror symmetry, you can combine more than one of these, e.g.ODD_Z+EVEN_Y
. This is especially useful in 2d simulations to restrict yourself to a desired polarization. -
eig_resolution
is the spatial resolution to use in MPB for the eigenmode calculations. -
eigensolver_tol
is the tolerance to use in the MPB eigensolver. MPB terminates when the eigenvalues stop changing by less than this fractional tolerance. -
coeffs
is a user-allocated array of typestd::complex<double>
(shortened hereafter tocdouble
) of length2*num_freqs*num_bands
wherenum_freqs
is the number of frequencies stored in theflux
object (equivalent toflux->Nfreq
) andnum_bands
is the length of thebands
input array. The expansion coefficients for the mode with frequencynf
and band indexnb
are stored sequentially as$\alpha^+$ ,$\alpha^-$ starting at slot2*nb*num_freqs+nf
of this array -
vgrp
is an optional user-allocateddouble
array of lengthnum_freqs*num_bands
. On return,vgrp[nb*num_freqs + nf]
is the group velocity of the mode with frequencynf
and band indexnb.
If you do not need this information, simply passNULL
for this parameter. -
user_kpoint_func
is an optional function you supply to provide an initial guess of the wavevector of a mode with given frequency and band index having the following prototype:
vec (*kpoint_func)(double freq, int mode, void *user_data);
-
user_kpoint_data
is the user data passed to theuser_kpoint_func
. -
kpoint_list
is a user allocated array ofmeep::vec
objects of length (num_bands*num_freqs
). If non-null, this array is filled in with the wavevectors of the eigenmode for each band from 1 tonum_bands*num_freqs
. -
kdom_list
is a user allocated array ofmeep::vec
objects of length (num_bands*num_freqs
). If non-null, this array is filled in with the wavevectors of the dominant planewave in the Fourier series expansion for each band from 1 tonum_bands*num_freqs
.kdom_list[nb*num_freqs + nf]
is the dominant planewave of the mode with frequencynf
and band indexnb
which defaults toNULL
. This is especially useful for interpreting the modes computed in a uniform medium, because those modes are exactly planewaves proportional to $\exp(2\pi i \vec{k}{dom}\cdot \vec{x})$ where $\vec{k}{dom}$ (kdom
) is the wavevector. -
cscale
is a user allocated array ofdouble
objects of length (num_bands*num_freqs
). If non-null, this array is filled in with scalar coefficients used for adjoint calculations. -
d
is the direction normal to the monitor plane. -
dp
is a user allocateddiffractedplanewave
used to specify a diffracted planewave in homogeneous media.
The following is a demonstration of typical usage:
int num_bands = bands.size();
int num_freqs = Flux->Nfreq;
std::vector<cdouble> coeffs(2*num_bands*num_freqs);
f.get_eigenmode_coefficients(...);
for(int nb=0; nb<num_bands; nb++)
for(int nf=0; nf<num_freqs++; nf++)
{
// get coefficients of forward- and backward-traveling
// waves in eigenmode bands[nb] at frequency nf
cdouble AlphaPlus = coeffs[2*nb*num_freqs+nf+0];
cdouble AlphaMinus = coeffs[2*nb*num_freqs+nf+1];
...
The get_eigenmode_coefficients
are normalized such that their squared magnitude equals the power carried by the corresponding eigenmode:
where
Besides get_eigenmode_coefficients
, there are a few computational routines in libmeep
that you may find useful for problems like those considered above.
void *fields::get_eigenmode(double frequency,
direction d,
const volume where,
const volume eig_vol,
int band_num,
const vec kpoint,
bool match_frequency,
int parity,
double resolution,
double eigensolver_tol,
double *kdom,
void **user_mdata,
diffractedplanewave *dp);
Calls MPB to compute the band_num
th eigenmode at frequency frequency
for the portion of your geometry lying in where
which is typically a cross-sectional slice of a waveguide. kpoint
is an initial starting guess for what the propagation vector of the waveguide mode will be. kdom
, if non-NULL and length 3, is filled in with the dominant planewave for the current band (see above). This is implemented in src/mpb.cpp.
The return value of get_eigenmode
is an opaque pointer to a data structure storing information about the computed eigenmode, which may be passed to the following routines:
// get a single component of the eigenmode field at a given point in space
std::complex<double> eigenmode_amplitude(const vec &p, void *vedata, component c);
// get the group velocity of the eigenmode
double get_group_velocity(void *vedata);
// free all memory associated with the eigenmode
void destroy_eigenmode_data(void *vedata);
These functions are implemented in src/mpb.cpp.
void output_dft(dft_flux flux, const char *HDF5FileName);
output_dft
exports the components of the frequency-domain fields stored in flux
to an HDF5 file with the given filename In general, flux
will store data for fields at multiple frequencies.
This function is implemented in src/dft.cpp.
std::complex<double> get_mode_flux_overlap(void *mode_data,
dft_flux *flux,
int num_freq,
std::complex<double>overlap[2]);
std::complex<double> get_mode_mode_overlap(void *mode1_data,
void *mode2_data,
dft_flux *flux,
std::complex<double>overlap[2]);
get_mode_flux_overlap
computes the overlap integral between the eigenmode described by mode_data
and the fields stored in flux
for the num_freq
th stored frequency, where num_freq
ranges from 0 to flux->Nfreq-1
. mode_data
should be the return value of a previous call to get_eigenmode.
get_mode_mode_overlap
is similar, but computes the overlap integral between two eigenmodes. mode1_data
and mode2_data
may be identical, in which case you get the inner product of the mode with itself. This should equal the group velocity of the mode based on the MPB's normalization convention.
These functions are implemented in src/dft.cpp.
The theoretical basis of the mode-decomposition algorithm is an orthogonality relation satisfied by the eigenmodes:
where the (indefinite) inner product between two vectors of the transverse field components, e.g.
$$ \left\langle \psi , \psi' \right\rangle \equiv \int_{S} \Big[ \mathbf{E}^(\vec \rho) \times \mathbf{H}'(\vec \rho) + \mathbf{E}'(\vec \rho) \times \mathbf{H}^(\vec \rho)\Big] \cdot \hat{\mathbf{n}} , dA \tag{5} $$
where
(There is some subtlety about the use of complex conjugation for orthogonality in the inner product above that arises for evanescent modes, which we don't consider here because Meep only computes coefficients of propagating modes. The above definition has the nice property that
Now consider a Meep calculation in which we have accumulated tranverse frequency-domain fields $\psi^{\text{meep}} = (\mathbf E^{\text{meep}}\parallel, \mathbf H^{\text{meep}}\parallel)$ on a dft_flux
object located on a cross-sectional surface
Taking the inner product of this equation with the
Some simplifications arise from the fact that, in a constant cross-section waveguide, the tangential components of the forward- and backward-traveling propagating modes are related by a mirror reflection in fields::get_mode_flux_overlap
. Although simple in principle, the implementation is complicated by the fact that, in multi-processor calculations, the Meep fields needed to evaluate the integrals are generally not all present on any one processor, but are instead distributed over multiple processors, requiring some interprocess communication to evaluate the full integral.
As mentioned above, the Poynting flux
where the right-hand-side is obtained by substituting the modal expansion and
using the mode orthogonality and the normalization
The get_eigenmode_coefficients
.