Using the mode solver for optical mode analysis#

This tutorial shows how to use the mode solver plugin in Tidy3D.

[1]:
# standard python imports
import numpy as np
import matplotlib.pylab as plt

# tidy3D import
import tidy3d as td
from tidy3d.constants import C_0
import tidy3d.web as web
from tidy3d.plugins.mode import ModeSolver

Setup#

We first set up the mode solver with information about our system. We start by setting parameters

[2]:
# size of simulation domain
Lx, Ly, Lz = 6, 6, 6
dl = 0.05

# waveguide information
wg_width = 1.5
wg_height = 1.0
wg_permittivity = 4.0

# central frequency
wvl_um = 2.0
freq0 = C_0 / wvl_um
fwidth = freq0 / 3

# run_time in ps
run_time = 1e-12

# automatic grid specification
grid_spec = td.GridSpec.auto(min_steps_per_wvl=20, wavelength=wvl_um)

Then we set up a simulation, in this case including a straight waveguide and periodic boundary conditions. Note that Tidy3D warns us that we have not added any sources in our Simulation object, however for purposes of mode solving it is not necessary.

[3]:
waveguide = td.Structure(
    geometry=td.Box(size=(wg_width, td.inf, wg_height)),
    medium=td.Medium(permittivity=wg_permittivity),
)

sim = td.Simulation(
    size=(Lx, Ly, Lz),
    grid_spec=grid_spec,
    structures=[waveguide],
    run_time=run_time,
    boundary_spec=td.BoundarySpec.all_sides(boundary=td.Periodic()),
)

ax = sim.plot(z=0)
plt.show()

../_images/notebooks_ModeSolver_5_0.png

Initialize Mode Solver#

With our system defined, we can now create our mode solver. We first need to specify on what plane we want to solve the modes using a td.Box() object.

[4]:
plane = td.Box(center=(0, 0, 0), size=(4, 0, 3.5))

The mode solver can now compute the modes given a ModeSpec object that specifies everything about the modes we’re looking for, for example:

  • num_modes: how many modes to compute.

  • target_neff: float, default=None, initial guess for the effective index of the mode; if not specified, the modes with the largest real part of the effective index are computed.

The full list of specification parameters can be found here.

[5]:
mode_spec = td.ModeSpec(
    num_modes=3,
    target_neff=2.0,
)

We can also specify a list of frequencies at which to solve for the modes.

[6]:
num_freqs = 11
f0_ind = num_freqs // 2
freqs = np.linspace(freq0 - fwidth / 2, freq0 + fwidth / 2, num_freqs)

Finally, we can initialize the ModeSolver, and call the solve method.

[7]:
mode_solver = ModeSolver(
    simulation=sim,
    plane=plane,
    mode_spec=mode_spec,
    freqs=freqs,
)
mode_data = mode_solver.solve()

[15:52:36] WARNING: Use the remote mode solver with subpixel averaging for      
           better accuracy through 'tidy3d.plugins.mode.web.run(...)'.          

We can also summarize useful mode information using to_dataframe(). Note that the group index was not computed; this can be included by setting group_index_step=True in the ModeSpec.

[8]:
mode_data.to_dataframe()
[8]:
wavelength n eff k eff TE (Ex) fraction wg TE fraction wg TM fraction mode area
f mode_index
1.249135e+14 0 2.400000 1.686130 0.0 0.995963 0.852047 0.891376 1.617417
1 2.400000 1.610272 0.0 0.010245 0.747363 0.926638 2.003674
2 2.400000 1.294575 0.0 0.226875 0.896082 0.622542 2.695679
1.299101e+14 0 2.307692 1.706504 0.0 0.996620 0.862301 0.895654 1.539832
1 2.307692 1.637271 0.0 0.008747 0.759892 0.930604 1.882750
2 2.307692 1.339798 0.0 0.203097 0.888992 0.647088 2.457125
1.349066e+14 0 2.222222 1.724923 0.0 0.997149 0.871672 0.899721 1.473565
1 2.222222 1.661620 0.0 0.007504 0.772128 0.934256 1.778317
2 2.222222 1.381435 0.0 0.180304 0.883673 0.670054 2.257858
1.399031e+14 0 2.142857 1.741628 0.0 0.997578 0.880228 0.903585 1.416463
1 2.142857 1.683633 0.0 0.006470 0.783911 0.937626 1.687731
2 2.142857 1.419700 0.0 0.159137 0.880026 0.691309 2.090148
1.448997e+14 0 2.068966 1.756829 0.0 0.997929 0.888036 0.907257 1.366859
1 2.068966 1.703585 0.0 0.005608 0.795151 0.940743 1.608796
2 2.068966 1.454838 0.0 0.139923 0.877877 0.710836 1.948158
1.498962e+14 0 2.000000 1.770701 0.0 0.998217 0.895166 0.910747 1.323446
1 2.000000 1.721718 0.0 0.004884 0.805802 0.943631 1.539695
2 2.000000 1.487099 0.0 0.122765 0.877010 0.728694 1.827302
1.548928e+14 0 1.935484 1.783397 0.0 0.998457 0.901682 0.914064 1.285192
1 1.935484 1.738240 0.0 0.004275 0.815850 0.946311 1.478919
2 1.935484 1.516732 0.0 0.107622 0.877209 0.744987 1.723899
1.598893e+14 0 1.875000 1.795048 0.0 0.998656 0.907645 0.917216 1.251271
1 1.875000 1.753333 0.0 0.003760 0.825298 0.948805 1.425223
2 1.875000 1.543974 0.0 0.094364 0.878263 0.759842 1.634964
1.648859e+14 0 1.818182 1.805767 0.0 0.998824 0.913108 0.920214 1.221021
1 1.818182 1.767156 0.0 0.003321 0.834165 0.951128 1.377568
2 1.818182 1.569047 0.0 0.082819 0.879989 0.773393 1.558063
1.698824e+14 0 1.764706 1.815652 0.0 0.998965 0.918121 0.923064 1.193899
1 1.764706 1.779845 0.0 0.002946 0.842474 0.953296 1.335090
2 1.764706 1.592156 0.0 0.072798 0.882226 0.785770 1.491208
1.748789e+14 0 1.714286 1.824788 0.0 0.999086 0.922729 0.925776 1.169463
1 1.714286 1.791520 0.0 0.002624 0.850255 0.955322 1.297071
2 1.714286 1.613486 0.0 0.064114 0.884840 0.797097 1.432774

Visualizing Mode Data#

The mode_info object contains information about the effective index of the mode and the field profiles, as well as the mode_spec that was used in the solver. The effective index data and the field profile data is in the form of xarray DataArrays.

We can for example plot the real part of the effective index for all three modes as follows.

[9]:
fig, ax = plt.subplots(1)
n_eff = mode_data.n_eff  # real part of the effective mode index
n_eff.plot.line(x="f")
plt.show()

../_images/notebooks_ModeSolver_17_0.png

The raw data can also be accessed.

[10]:
n_complex = mode_data.n_complex  # complex effective index as a DataArray
n_eff = mode_data.n_eff.values  # real part of the effective index as numpy array
k_eff = mode_data.k_eff.values  # imag part of the effective index as numpy array

print(
    f"first mode effective index at freq0: n_eff = {n_eff[f0_ind, 0]:.2f}, k_eff = {k_eff[f0_ind, 0]:.2e}"
)

first mode effective index at freq0: n_eff = 1.77, k_eff = 0.00e+00

The fields stored in mode_data can be visualized using in-built xarray methods.

[11]:
f, (ax1, ax2) = plt.subplots(1, 2, tight_layout=True, figsize=(10, 3))
abs(mode_data.Ex.isel(mode_index=0, f=f0_ind)).plot(x="x", y="z", ax=ax1, cmap="magma")
abs(mode_data.Ez.isel(mode_index=0, f=f0_ind)).plot(x="x", y="z", ax=ax2, cmap="magma")

ax1.set_title("|Ex(x, y)|")
ax1.set_aspect("equal")
ax2.set_title("|Ez(x, y)|")
ax2.set_aspect("equal")
plt.show()

../_images/notebooks_ModeSolver_21_0.png

Alternatively, we can use the in-built plot_field method of mode_data, which also allows us to overlay the structures in the simulation. The image also looks slightly different because the plot_field method uses robust=True option by default, which scales the colorbar to between the 2nd and 98th percentile of the data.

[12]:
f, (ax1, ax2) = plt.subplots(1, 2, tight_layout=True, figsize=(10, 3))
mode_solver.plot_field("Ex", "abs", mode_index=0, f=freq0, ax=ax1)
mode_solver.plot_field("Ez", "abs", mode_index=0, f=freq0, ax=ax2)
plt.show()

../_images/notebooks_ModeSolver_23_0.png

Choosing the mode of interest#

We can also look at the other modes that were computed.

[13]:
mode_index = 1
f, (ax1, ax2) = plt.subplots(1, 2, tight_layout=True, figsize=(10, 3))
mode_solver.plot_field("Ex", "abs", mode_index=mode_index, f=freq0, ax=ax1)
mode_solver.plot_field("Ez", "abs", mode_index=mode_index, f=freq0, ax=ax2)
plt.show()

../_images/notebooks_ModeSolver_25_0.png

This looks like an Ez-dominant mode. Finally, next-order mode has mixed polarization.

[14]:
mode_index = 2
f, (ax1, ax2) = plt.subplots(1, 2, tight_layout=True, figsize=(10, 3))
mode_solver.plot_field("Ex", "abs", mode_index=mode_index, f=freq0, ax=ax1)
mode_solver.plot_field("Ez", "abs", mode_index=mode_index, f=freq0, ax=ax2)
plt.show()

../_images/notebooks_ModeSolver_27_0.png

Exporting Results#

This looks promising!

Now we can choose the mode specifications to use in our mode source and mode monitors. These can be created separately, can be exported directly from the mode solver, for example:

[15]:
# Makes a modal source with geometry of `plane` with modes specified by `mode_spec` and a selected `mode_index`
source_time = td.GaussianPulse(freq0=freq0, fwidth=fwidth)
mode_src = mode_solver.to_source(mode_index=2, source_time=source_time, direction="-")

# Makes a mode monitor with geometry of `plane`.
mode_mon = mode_solver.to_monitor(name="mode", freqs=freqs)
# Offset the monitor along the propagation direction
mode_mon = mode_mon.copy(update=dict(center=(0, -2, 0)))

[16]:
# In-plane field monitor, slightly offset along x
monitor = td.FieldMonitor(
    center=(0, 0, 0.1), size=(td.inf, td.inf, 0), freqs=[freq0], name="field"
)

sim = td.Simulation(
    size=(Lx, Ly, Lz),
    grid_spec=grid_spec,
    run_time=run_time,
    boundary_spec=td.BoundarySpec.all_sides(boundary=td.PML()),
    structures=[waveguide],
    sources=[mode_src],
    monitors=[monitor, mode_mon],
)

sim.plot(z=0)
plt.show()

[15:54:13] WARNING: Default value for the field monitor 'colocate' setting has  
           changed to 'True' in Tidy3D 2.4.0. All field components will be      
           colocated to the grid boundaries. Set to 'False' to get the raw      
           fields on the Yee grid instead.                                      
../_images/notebooks_ModeSolver_30_1.png
[17]:
job = web.Job(simulation=sim, task_name="mode_simulation", verbose=True)
sim_data = job.run(path="data/simulation_data.hdf5")

[15:54:14] Created task 'mode_simulation' with task_id
           'fdve-bc0cf21a-7395-4ea6-9794-5c3668537f92v1'.
[15:54:20] status = queued
[15:54:24] status = preprocess
[15:54:30] Maximum FlexCredit cost: 0.025. Use 'web.real_cost(task_id)' to get
           the billed FlexCredit cost after a simulation run.
           starting up solver
           running solver
           To cancel the simulation, use 'web.abort(task_id)' or
           'web.delete(task_id)' or abort/delete the task in the web UI.
           Terminating the Python script will not stop the job running on the
           cloud.
[15:54:38] early shutoff detected, exiting.
           status = postprocess
[15:54:42] status = success
[15:54:43] loading SimulationData from data/simulation_data.hdf5

We can now plot the in-plane field and the modal amplitudes. Since we injected mode 2 and we just have a straight waveguide, all the power recorded by the modal monitor is in mode 2, going backwards.

[18]:
fig, ax = plt.subplots(1, 2, figsize=(10, 4))
sim_data.plot_field("field", "Ez", f=freq0, ax=ax[0])
sim_data["mode"].amps.sel(direction="-").abs.plot.line(x="f", ax=ax[1])
plt.show()

../_images/notebooks_ModeSolver_33_0.png

Storing server-side computed modes#

We can also use a ModeSolverMonitor to store the modes as they are computed server-side. This is illustrated below. We will also request in the mode specification that the modes are filtered by their tm polarization. In this particular simulation, TM refers to Ez polarization. The effect of the filtering is that modes with a tm polarization fraction larger than or equal to 0.5 will come first in the list of modes (while still ordered by decreasing effective index). After that, the set of predominantly te-polarized modes (tm fraction < 0.5) follows.

[23]:
mode_spec = mode_spec.copy(update=dict(filter_pol="tm"))

# Update mode source to use the highest-tm-fraction mode
mode_src = mode_src.copy(update=dict(mode_spec=mode_spec))
mode_src = mode_src.copy(update=dict(mode_index=0))

# Update mode monitor to use the tm_fraction ordered mode_spec
mode_mon = mode_mon.copy(update=dict(mode_spec=mode_spec))

# New monitor to record the modes computed at the mode decomposition monitor location
mode_solver_mon = td.ModeSolverMonitor(
    center=mode_mon.center,
    size=mode_mon.size,
    freqs=mode_mon.freqs,
    mode_spec=mode_spec,
    name="mode_solver",
)

sim = td.Simulation(
    size=(Lx, Ly, Lz),
    grid_spec=grid_spec,
    run_time=run_time,
    boundary_spec=td.BoundarySpec.all_sides(boundary=td.PML()),
    structures=[waveguide],
    sources=[mode_src],
    monitors=[monitor, mode_mon, mode_solver_mon],
)

[24]:
job = web.Job(simulation=sim, task_name="mode_simulation", verbose=True)
sim_data = job.run(path="data/simulation_data.hdf5")

[15:14:53] Created task 'mode_simulation' with task_id
           'fdve-58775b76-2c3f-4ede-9c56-ff331ad4eb07v1'.
[15:14:54] status = queued
[15:15:00] status = preprocess
[15:15:06] Maximum FlexCredit cost: 0.025. Use 'web.real_cost(task_id)' to get
           the billed FlexCredit cost after a simulation run.
           starting up solver
           running solver
           To cancel the simulation, use 'web.abort(task_id)' or
           'web.delete(task_id)' or abort/delete the task in the web UI.
           Terminating the Python script will not stop the job running on the
           cloud.
[15:15:14] early shutoff detected, exiting.
           status = postprocess
[15:15:20] status = success
[15:15:22] loading SimulationData from data/simulation_data.hdf5

Note the different ordering of the recorded modes compared to what we saw above, with the fundamental TM mode coming first.

[25]:
fig, ax = plt.subplots(1)
n_eff = sim_data["mode"].n_eff  # real part of the effective mode index
n_eff.plot.line(x="f")
plt.show()

../_images/notebooks_ModeSolver_38_0.png

Now the fundamental Ez-polarized mode is injected, and as before it is the only one that the mode monitor records any intensity in.

[26]:
fig, ax = plt.subplots(1, 2, figsize=(10, 4))
sim_data.plot_field("field", "Ez", f=freq0, ax=ax[0])
sim_data["mode"].amps.sel(direction="-").abs.plot.line(x="f", ax=ax[1])
plt.show()

../_images/notebooks_ModeSolver_40_0.png

We can also have a look at the mode fields stored in the ModeFieldMonitor either directly using xarray methods as above, or using the Tidy3D SimulationData in-built field plotting.

[27]:
fig, ax = plt.subplots(1, 2, figsize=(12, 4))
sim_data.plot_field("mode_solver", "Ex", f=freq0, val="abs", mode_index=0, ax=ax[0])
sim_data.plot_field("mode_solver", "Ez", f=freq0, val="abs", mode_index=0, ax=ax[1])
plt.show()

../_images/notebooks_ModeSolver_42_0.png

Computing modes in the server#

The advantage of the server-side mode calculation is that it can use subpixel averaging, which increases the computation accuracy, something that the local solver cannot do. The web module from the mode solver plugin allows us to skip the local solver entirely and use the server-side computation instead. We start by importing the run function, which will automate the complete process for us. If we want finer control over the server interaction, the class ModeSolverTask can also be used.

(Note that using the mode solver on the server will cost a small amount of FlexCredits.)

[28]:
from tidy3d.plugins.mode.web import run as run_mode_solver

We use the previous mode solver as an example:

[24]:
mode_data = run_mode_solver(mode_solver)
[13:59:35] Mode solver created with                                    web.py:80
           task_id='fdve-ea498b47-738e-4a8c-aea3-21d5f1a83878v1',               
           solver_id='mo-edf74f79-26da-45d6-81f1-8ed812458d87'.                 
[13:59:36] Mode solver status: queued                                  web.py:93
[13:59:39] Mode solver status: running                                 web.py:93
[13:59:57] Mode solver status: success                                web.py:103
[25]:
fig, ax = plt.subplots(1)
n_eff = mode_data.n_eff
n_eff.plot.line(x="f")
plt.show()

../_images/notebooks_ModeSolver_47_0.png

To demonstrate the accuracy improvement, we can show a convergence plot comparing the local and remote solvers on a waveguide that doesn’t align with the mesh, in this case a dielectric rod.

[26]:
resolutions = range(10, 31, 2)
n_eff = ([], [])

cylinder = td.Structure(
    geometry=td.Cylinder(axis=1, center=(0, 0, 0), radius=0.7, length=2 * Lz),
    medium=td.Medium(permittivity=wg_permittivity),
)

dummy_source = td.UniformCurrentSource(
    size=(0,0,0),
    source_time=td.GaussianPulse(freq0=freq0, fwidth=fwidth),
    polarization='Ex',
)

mode_spec = td.ModeSpec(
    num_modes=1,
    target_neff=2.0,
)

def create_mode_solver(resolution):
    grid_spec = td.GridSpec.auto(min_steps_per_wvl=resolution, wavelength=wvl_um)

    sim = td.Simulation(
        size=(Lx, Ly, Lz),
        grid_spec=grid_spec,
        structures=[cylinder],
        sources=[dummy_source],
        run_time=run_time,
        boundary_spec=td.BoundarySpec.all_sides(boundary=td.Periodic()),
    )

    return ModeSolver(
        simulation=sim,
        plane=td.Box(center=(0, 0, 0), size=(4, 0, 4)),
        mode_spec=mode_spec,
        freqs=[freq0],
    )


for resolution in resolutions:
    print("Solving for resolution =", resolution, flush=True)
    mode_solver = create_mode_solver(resolution)
    local_data = mode_solver.solve()
    n_eff[0].append(local_data.n_eff.item())

    mode_solver = create_mode_solver(resolution)
    server_data = run_mode_solver(mode_solver, verbose=False)
    n_eff[1].append(server_data.n_eff.item())
Solving for resolution = 10
Solving for resolution = 12
Solving for resolution = 14
Solving for resolution = 16
Solving for resolution = 18
Solving for resolution = 20
Solving for resolution = 22
Solving for resolution = 24
Solving for resolution = 26
Solving for resolution = 28
Solving for resolution = 30
[27]:
fig, ax = plt.subplots(1)
ax.plot(resolutions, n_eff[0], label="Local")
ax.plot(resolutions, n_eff[1], label="Server")
ax.set(xlabel="Resolution", ylabel="Progagation index")
ax.legend()
ax.grid()
../_images/notebooks_ModeSolver_50_0.png

Mode tracking#

As typical for eigenvalue type solvers the Tidy3D’s mode solver returns calculated modes in the order corresponding to the magnitudes of the found eigenvalues (effective index). Since the effective index of a mode generally depends on frequency and can become larger or smaller in magnitude compared to the effective index values of other modes, the same mode_index value in the returned mode solver data may correspond to physically different modes for different frequencies. To avoid such a mismatch ModeSolver by default performs sorting of modes based on their overlap values (see discussion on mode decomposition, that is, modes at one frequency are matched to the most similar modes at the other frequency. This behavior is controlled by the parameter track_freq of ModeSpec object. It can be set to either None, "lowest", "central" (default), or "highest", where in the first case no sorting is performed while in the other three cases sorting is performed starting from the specified frequency. Additionally, any unsorted mode solver data (objects of type ModeSolverData) returned by ModeSolver.solve() can be sorted afterwards by using function ModeSolverData.overlap_sort(). This function takes in arguments track_freq (default: "central"), which has the same meaning as in ModeSpec, and overlap_thresh (default: 0.9) that defines the threshold overlap value above which a pair of spatial fields is considered to correspond to physically the same mode. While this functionality is not expected to be used in the most cases, we still demonstrate it here to provide a more clear understanding of ModeSolver results. Note that for mode sources and mode monitors sorting is always performed internally.

We start with calculating twelve modes at eleven different frequencies for a layered waveguide with the default sorting turned off.

[28]:
# bottom layer
bottom = td.Structure(
    geometry=td.Box(
        center=(0, 0, -0.8 * wg_height), size=(wg_width, td.inf, wg_height)
    ),
    medium=td.Medium(permittivity=wg_permittivity * 2),
)

# top layer
top = td.Structure(
    geometry=td.Box(
        center=(0, 0, 0.7 * wg_height), size=(wg_width, td.inf, wg_height * 0.8)
    ),
    medium=td.Medium(permittivity=wg_permittivity * 2),
)

# new simulation object
sim = td.Simulation(
    size=(Lx, Ly, Lz),
    grid_spec=grid_spec,
    structures=[waveguide, bottom, top],
    run_time=run_time,
    boundary_spec=td.BoundarySpec.all_sides(boundary=td.Periodic()),
)

# visualize
sim.plot(y=0)
plt.show()

../_images/notebooks_ModeSolver_52_0.png

Again, Tidy3D is warning us that the simulation does not contain a source. However, since this simulation is used to construct the mode solver and will not be run directly, we can ignore this warning.

[29]:
# set mode spec to find 12 modes and turn off automatic mode tracking
mode_spec = td.ModeSpec(num_modes=12, track_freq=None)

# we will use a larger solver plane for this structure
plane = td.Box(center=(0, 0, 0), size=(10, 0, 10))

# 11 frequencies in +/- fwidth range
num_freqs = 11
freqs = np.linspace(freq0 - fwidth, freq0 + fwidth, num_freqs)

# make new mode solver object and solve for modes
mode_solver = ModeSolver(
    simulation=sim,
    plane=plane,
    mode_spec=mode_spec,
    freqs=freqs,
)
mode_data = mode_solver.solve()

By inspecting the values of effective indices we note that some of them come very close to each other, which is a sign that some mode index values might not correspond to the same modes throughout the considered frequency range.

[30]:
fig, ax = plt.subplots(1)
n_eff = mode_data.n_eff  # real part of the effective mode index
n_eff.plot.line(".-", x="f")
plt.show()

../_images/notebooks_ModeSolver_56_0.png

Indeed, by plotting field distributions for the calculated modes we see mixing of modes for mode_index > 4.

[31]:
fig, ax = plt.subplots(
    mode_spec.num_modes,
    len(freqs),
    tight_layout=True,
    figsize=(len(freqs), mode_spec.num_modes),
)
for j in range(mode_spec.num_modes):
    for i, freq in enumerate(freqs):
        ax[j, i].imshow(
            mode_data.Ex.isel(mode_index=j, f=i, drop=True).squeeze().abs.T,
            cmap="magma",
            origin="lower",
        )
        ax[j, i].axis("off")
plt.show()

../_images/notebooks_ModeSolver_58_0.png

This situation can be remedied by using overlap_sort() function as follows.

[32]:
mode_data_sorted = mode_data.overlap_sort(track_freq="central")

Plotting the same information as previously one can confirm that modes are matching now as closely as possible.

[33]:
n_eff_sorted = mode_data_sorted.n_eff  # real part of the effective mode index
n_eff_sorted.plot.line(".-", x="f")
plt.show()

fig, ax = plt.subplots(
    mode_spec.num_modes,
    len(freqs),
    tight_layout=True,
    figsize=(len(freqs), mode_spec.num_modes),
)
for j in range(mode_spec.num_modes):
    for i, freq in enumerate(freqs):
        ax[j, i].imshow(
            mode_data_sorted.Ex.isel(mode_index=j, f=i, drop=True).squeeze().abs.T,
            cmap="magma",
            origin="lower",
        )
        ax[j, i].axis("off")
plt.show()

../_images/notebooks_ModeSolver_62_0.png
../_images/notebooks_ModeSolver_62_1.png

Note that depending on particular situation the mode sorting, whether default or by hands, might not succesfully resolve tracking of all modes throughout the chosen frequency range. Among possible reasons for that: - A mode is physically disappearing at a certain frequency point. - A mode is not included in calculated number of modes, num_modes, for certain frequency points. In this case increasing num_modes might alleviate the issue. - Presence of degenerate modes at certain frequencies.

Group Index#

The mode solver can also be used to directly calculate the group index of the calculated mode. It is important to have mode tracking enabled, to avoid differentiating wrong mode pairs around mode crossings.

[34]:
# Use the previous mode solver with the appropriate settings
mode_spec = td.ModeSpec(num_modes=12, track_freq="central", group_index_step=True)
mode_solver = mode_solver.copy(
    update={"mode_spec": mode_spec}
)

mode_data = mode_solver.solve()
n_eff = mode_data.n_eff
n_group = mode_data.n_group

_, ax = plt.subplots(1, 2, figsize=(10, 5), tight_layout=True)

n_eff.plot.line(".-", x="f", ax=ax[0])
n_group.plot.line(".-", x="f", ax=ax[1])
ax[1].set_ylabel("Group Index")
plt.show()
../_images/notebooks_ModeSolver_65_0.png

Notes / Considerations#

  • This mode solver runs locally, which means it does not require credits to run.

  • It also means that the mode solver does not use subpixel-smoothing, even if this is specified in the simulation. On the other hand, when the modes are stored in a ModeSolverMonitor during a simulation run, the subpixel averaging is applied. Therefore, the latter results may not perfectly match those from a local run, and are more accurate.

  • Symmetries are applied if they are defined in the simulation and the mode plane center sits on the simulation center.

  • In cases when the number of modes/frequency points and/or the dimensions of mode solver plane are large and tracking of modes is not important then it might be beneficial to turn off the default mode sorting for increased computational efficiency.