Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

eigenmodes of diffracted planewaves #1316

Merged
merged 15 commits into from
Aug 26, 2020
Merged
Show file tree
Hide file tree
Changes from 3 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
34 changes: 34 additions & 0 deletions python/meep.i
Original file line number Diff line number Diff line change
Expand Up @@ -512,6 +512,23 @@ kpoint_list get_eigenmode_coefficients_and_kpoints(meep::fields *f, meep::dft_fl
return res;
}

kpoint_list get_eigenmode_coefficients_and_kpoints(meep::fields *f, meep::dft_flux flux, const meep::volume &eig_vol,
meep::diffractedplanewave dp, int parity, double eig_resolution,
double eigensolver_tol, std::complex<double> *coeffs,
double *vgrp, meep::kpoint_func user_kpoint_func,
void *user_kpoint_data, double *cscale, meep::direction d) {

size_t num_kpoints = flux.freq.size();
meep::vec *kpoints = new meep::vec[num_kpoints];
meep::vec *kdom = new meep::vec[num_kpoints];
f->get_eigenmode_coefficients(flux, eig_vol, NULL, 1, parity, eig_resolution, eigensolver_tol,
coeffs, vgrp, user_kpoint_func, user_kpoint_data, kpoints, kdom, cscale, d, &dp);

kpoint_list res = {kpoints, num_kpoints, kdom, num_kpoints};

return res;
}

PyObject *_get_array_slice_dimensions(meep::fields *f, const meep::volume &where, size_t dims[3],
bool collapse_empty_dimensions, bool snap_empty_dimensions,
meep::component cgrid = Centered, PyObject *min_max_loc = NULL) {
Expand Down Expand Up @@ -1034,6 +1051,16 @@ void _get_gradient(PyObject *grad, PyObject *fields_a, PyObject *fields_f, PyObj
$1 = (std::complex<double> *)array_data($input);
}

// typemaps for diffractedplanewave

%typecheck(SWIG_TYPECHECK_POINTER, fragment="NumPy_Fragments") double axis[3] {
$1 = is_array($input);
}

%typemap(in) double axis[3] {
$1 = (double *)array_data($input);
}

// typemaps for gaussianbeam

%typecheck(SWIG_TYPECHECK_POINTER, fragment="NumPy_Fragments") std::complex<double> E0[3] {
Expand Down Expand Up @@ -1476,6 +1503,12 @@ kpoint_list get_eigenmode_coefficients_and_kpoints(meep::fields *f, meep::dft_fl
std::complex<double> *coeffs, double *vgrp,
meep::kpoint_func user_kpoint_func, void *user_kpoint_data,
double *cscale, meep::direction d);
kpoint_list get_eigenmode_coefficients_and_kpoints(meep::fields *f, meep::dft_flux flux,
const meep::volume &eig_vol, meep::diffractedplanewave dp,
int parity, double eig_resolution, double eigensolver_tol,
std::complex<double> *coeffs, double *vgrp,
meep::kpoint_func user_kpoint_func, void *user_kpoint_data,
double *cscale, meep::direction d);
PyObject *_get_array_slice_dimensions(meep::fields *f, const meep::volume &where, size_t dims[3],
bool collapse_empty_dimensions, bool snap_empty_dimensions,
meep::component cgrid = Centered, PyObject *min_max_loc = NULL);
Expand Down Expand Up @@ -1562,6 +1595,7 @@ PyObject *_get_array_slice_dimensions(meep::fields *f, const meep::volume &where
DftEnergy,
DftFields,
Volume,
DiffractedPlanewave,
after_sources,
after_sources_and_time,
after_time,
Expand Down
94 changes: 75 additions & 19 deletions python/simulation.py
Original file line number Diff line number Diff line change
Expand Up @@ -99,6 +99,32 @@ def py_v3_to_vec(dims, iterable, is_cylindrical=False):
else:
raise ValueError("Invalid dimensions in Volume: {}".format(dims))

class DiffractedPlanewave(object):
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Note: the current plan is to optionally pass this for the band_num parameter to eigenmode sources etcetera.

def __init__(self,
g=None,
axis=None,
Copy link
Collaborator

@stevengj stevengj Aug 14, 2020

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Eventually, we'll have the EigenmodeSource etcetera choose a reasonable default if axis=None when the source is added.

(You'll generally want axis to lie in the plane of the source, I think. In 2d you'll also want it to be in the xy plane, i.e. it should be along the source line.)

s=None,
p=None):
self._g = g
self._axis = Vector3(*axis)
self._s = complex(s)
self._p = complex(p)

@property
def g(self):
return self._g

@property
def axis(self):
return self._axis

@property
def s(self):
return self._s

@property
def p(self):
return self._p

class PML(object):
"""
Expand Down Expand Up @@ -3257,25 +3283,55 @@ def get_eigenmode_coefficients(self, flux, bands, eig_parity=mp.NO_PARITY, eig_v
if direction is None or direction == mp.AUTOMATIC:
direction = flux.normal_direction

num_bands = len(bands)
coeffs = np.zeros(2 * num_bands * flux.freq.size(), dtype=np.complex128)
vgrp = np.zeros(num_bands * flux.freq.size())
cscale = np.zeros(num_bands * flux.freq.size())

kpoints, kdom = mp.get_eigenmode_coefficients_and_kpoints(
self.fields,
flux.swigobj,
eig_vol,
np.array(bands, dtype=np.intc),
eig_parity,
eig_resolution,
eig_tolerance,
coeffs,
vgrp,
kpoint_func,
cscale,
direction
)
if isinstance(bands,(list,range)):
num_bands = len(bands)
coeffs = np.zeros(2 * num_bands * flux.freq.size(), dtype=np.complex128)
vgrp = np.zeros(num_bands * flux.freq.size())
cscale = np.zeros(num_bands * flux.freq.size())

kpoints, kdom = mp.get_eigenmode_coefficients_and_kpoints(
self.fields,
flux.swigobj,
eig_vol,
np.array(bands, dtype=np.intc),
eig_parity,
eig_resolution,
eig_tolerance,
coeffs,
vgrp,
kpoint_func,
cscale,
direction
)
elif isinstance(bands,DiffractedPlanewave):
num_bands = 1
coeffs = np.zeros(2 * num_bands * flux.freq.size(), dtype=np.complex128)
vgrp = np.zeros(num_bands * flux.freq.size())
cscale = np.zeros(num_bands * flux.freq.size())
diffractedplanewave_args = [
np.array(bands.g, dtype=np.intc),
np.array([bands.axis.x, bands.axis.y, bands.axis.z], dtype=np.float64),
bands.s * 1.0,
bands.p * 1.0
]
diffractedplanewave = mp.diffractedplanewave(*diffractedplanewave_args)

kpoints, kdom = mp.get_eigenmode_coefficients_and_kpoints(
self.fields,
flux.swigobj,
eig_vol,
diffractedplanewave,
eig_parity,
eig_resolution,
eig_tolerance,
coeffs,
vgrp,
kpoint_func,
cscale,
direction
)
else:
raise TypeError("get_eigenmode_coefficients: bands must be either a list or DiffractedPlanewave object")

return EigCoeffsResult(np.reshape(coeffs, (num_bands, flux.freq.size(), 2)), vgrp, kpoints, kdom, cscale)

Expand Down
24 changes: 21 additions & 3 deletions src/meep.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1484,6 +1484,22 @@ class chunkloop_field_components {
std::vector<std::complex<double> > values; // updated by update_values(idx)
};

class diffractedplanewave {

public:
diffractedplanewave(int g[3], double axis[3], std::complex<double> s, std::complex<double> p);
int *get_g() { return g; };
double *get_axis() { return axis; }
std::complex<double> get_s() const { return s; };
std::complex<double> get_p() const { return p; };

private:
int g[3]; // diffraction order
double axis[3]; // axis vector
std::complex<double> s; // s polarization amplitude
std::complex<double> p; // p polarization ampiltude
};

class gaussianbeam {

public:
Expand Down Expand Up @@ -1710,7 +1726,7 @@ class fields {
void *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 = 0,
void **user_mdata = 0);
void **user_mdata = 0, diffractedplanewave *dp = 0);

void add_eigenmode_source(component c, const src_time &src, direction d, const volume &where,
const volume &eig_vol, int band_num, const vec &kpoint,
Expand All @@ -1722,12 +1738,14 @@ class fields {
int parity, double eig_resolution, double eigensolver_tol,
std::complex<double> *coeffs, double *vgrp,
kpoint_func user_kpoint_func, void *user_kpoint_data,
vec *kpoints, vec *kdom, double *cscale, direction d);
vec *kpoints, vec *kdom, double *cscale, direction d,
diffractedplanewave *dp = 0);
void 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 = 0, void *user_kpoint_data = 0,
vec *kpoints = 0, vec *kdom = 0, double *cscale = 0);
vec *kpoints = 0, vec *kdom = 0, double *cscale = 0,
diffractedplanewave *dp = 0);

// initialize.cpp:
void initialize_field(component, std::complex<double> f(const vec &));
Expand Down
46 changes: 38 additions & 8 deletions src/mpb.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -292,6 +292,15 @@ static int nextpow2357(int n) {
}
}

diffractedplanewave::diffractedplanewave(int g_[3], double axis_[3], std::complex<double> s_, std::complex<double> p_) {
for (int j = 0; j < 3; ++j) {
g[j] = g_[j];
axis[j] = axis_[j];
}
s = s_;
p = p_;
};

/****************************************************************/
/* call MPB to get the band_numth eigenmode at freq frequency. */
/* */
Expand All @@ -313,7 +322,7 @@ static int nextpow2357(int n) {
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) {
void **user_mdata, diffractedplanewave *dp) {
/*--------------------------------------------------------------*/
/*- part 1: preliminary setup for calling MPB -----------------*/
/*--------------------------------------------------------------*/
Expand Down Expand Up @@ -538,7 +547,7 @@ void *fields::get_eigenmode(double frequency, direction d, const volume where, c
/*- part 2: newton iteration loop with call to MPB on each step */
/*- until eigenmode converged to requested tolerance */
/*--------------------------------------------------------------*/
if (am_master()) do {
if (am_master() && !dp) do {
eigensolver(H, eigvals, maxwell_operator, (void *)mdata,
#if MPB_VERSION_MAJOR > 1 || (MPB_VERSION_MAJOR == 1 && MPB_VERSION_MINOR >= 6)
NULL, NULL, /* eventually, we can support mu here */
Expand Down Expand Up @@ -591,6 +600,28 @@ void *fields::get_eigenmode(double frequency, direction d, const volume where, c
} while (match_frequency &&
fabs(sqrt(eigvals[band_num - 1]) - frequency) > frequency * match_tol);

if (am_master() && dp) {
scalar_complex s = {real(dp->get_s()),imag(dp->get_s())};
scalar_complex p = {real(dp->get_p()),imag(dp->get_p())};

master_printf("diffracted planewave: k = %g, %g, %g for omega = %g\n", k[0],k[1],k[2], frequency);
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It needs to first compute the correct k vector from the frequency and get_g().

update_maxwell_data_k(mdata, k, G[0], G[1], G[2]);
maxwell_set_planewave(mdata, H, band_num, dp->get_g(), s, p, dp->get_axis());

eigvals[band_num - 1] = frequency;
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

or if match_frequency is false then it needs to compute the frequency from k and get_g().


evectmatrix_resize(&W[0], 1, 0);
evectmatrix_resize(&W[1], 1, 0);
for (int i = 0; i < H.n; ++i)
W[0].data[i] = H.data[H.p - 1 + i * H.p];
maxwell_ucross_op(W[0], W[1], mdata, kdir); // W[1] = (dTheta/dk) W[0]
mpb_real vscratch;
evectmatrix_XtY_diag_real(W[0], W[1], &vgrp, &vscratch);
vgrp /= sqrt(eigvals[band_num - 1]);

master_printf("diffracted planewave: group velocity v=%g\n", vgrp);
}

double eigval = eigvals[band_num - 1];

// cleanup temporary storage
Expand Down Expand Up @@ -828,10 +859,9 @@ void fields::get_eigenmode_coefficients(dft_flux flux, const volume &eig_vol, in
double eigensolver_tol, std::complex<double> *coeffs,
double *vgrp, kpoint_func user_kpoint_func,
void *user_kpoint_data, vec *kpoints, vec *kdom_list,
double *cscale, direction d) {
double *cscale, direction d, diffractedplanewave *dp) {
int num_freqs = flux.freq.size();
bool match_frequency = true;

if (flux.use_symmetry && S.multiplicity() > 1 && parity == 0)
abort("flux regions for eigenmode projection with symmetry should be created by "
"add_mode_monitor()");
Expand All @@ -847,13 +877,13 @@ void fields::get_eigenmode_coefficients(dft_flux flux, const volume &eig_vol, in
/*--------------------------------------------------------------*/
/*- call mpb to compute the eigenmode --------------------------*/
/*--------------------------------------------------------------*/
int band_num = bands[nb];
int band_num = bands ? bands[nb] : 1;
double kdom[3];
if (user_kpoint_func) kpoint = user_kpoint_func(flux.freq[nf], band_num, user_kpoint_data);
am_now_working_on(MPBTime);
void *mode_data =
get_eigenmode(flux.freq[nf], d, flux.where, eig_vol, band_num, kpoint, match_frequency,
parity, eig_resolution, eigensolver_tol, kdom, (void **)&mdata);
parity, eig_resolution, eigensolver_tol, kdom, (void **)&mdata, dp);
finished_working();
if (!mode_data) { // mode not found, assume evanescent
coeffs[2 * nb * num_freqs + 2 * nf] = coeffs[2 * nb * num_freqs + 2 * nf + 1] = 0;
Expand Down Expand Up @@ -1001,10 +1031,10 @@ void fields::get_eigenmode_coefficients(dft_flux flux, const volume &eig_vol, in
double eigensolver_tol, std::complex<double> *coeffs,
double *vgrp, kpoint_func user_kpoint_func,
void *user_kpoint_data, vec *kpoints, vec *kdom,
double *cscale) {
double *cscale, diffractedplanewave *dp) {
get_eigenmode_coefficients(flux, eig_vol, bands, num_bands, parity, eig_resolution,
eigensolver_tol, coeffs, vgrp, user_kpoint_func, user_kpoint_data,
kpoints, kdom, cscale, flux.normal_direction);
kpoints, kdom, cscale, flux.normal_direction, dp);
}

} // namespace meep