The directional coupler as well as the source and mode monitor geometries are described by the GDSII file examples/coupler.gds. A snapshot of this file viewed using KLayout is shown below. The figure labels have been added in post processing. The design consists of two identical strip waveguides which are positioned close together via an adiabatic taper such that their modes couple evanescently. There is a source (labelled "Source") and four mode monitors (labelled "Port 1", etc.). The input pulse from Port 1 is split in two and exits through Ports 3 and 4. The design objective is to find the separation distance (labelled "d") which maximizes power in Port 4 at a wavelength of 1.55 μm. More generally, though not included in this example, it is possible to have two additional degrees of freedom: (1) the length of the straight waveguide section where the two waveguides are coupled and (2) the length of the tapered section (the taper profile is described by a hyperbolic tangent (tanh) function).
The GDSII file is adapted from the SiEPIC EBeam PDK with four major modifications:
the computational cell is centered at the origin of the xy plane and defined on layer 0
the source and four mode monitors are defined on layers 1-5
the lower and upper branches of the coupler are defined on layers 31 and 32
the straight waveguide sections are perfectly linear
Note that rather than being specified as part of the GDSII file, the volume regions of the source and flux monitors could have been specified in the simulation script.
import meep as mp
import numpy
import matplotlib.pyplot as plt
res = 50 # pixels/μm
three_d = False # 3d calculation?
d = 0.12 # branch separation
gdsII_file = 'coupler.gds'
CELL_LAYER = 0
PORT1_LAYER = 1
PORT2_LAYER = 2
PORT3_LAYER = 3
PORT4_LAYER = 4
SOURCE_LAYER = 5
UPPER_BRANCH_LAYER = 31
LOWER_BRANCH_LAYER = 32
default_d = 0.3
t_oxide = 1.0
t_Si = 0.22
t_air = 0.78
dpml = 1
cell_thickness = dpml+t_oxide+t_Si+t_air+dpml
si_zmin = 0
oxide = mp.Medium(epsilon=2.25)
silicon=mp.Medium(epsilon=12)
lcen = 1.55
fcen = 1/lcen
df = 0.2*fcen
cell_zmax = 0.5*cell_thickness if three_d else 0
cell_zmin = -0.5*cell_thickness if three_d else 0
si_zmax = t_Si if three_d else 0
# read cell size, volumes for source region and flux monitors,
# and coupler geometry from GDSII file
upper_branch = mp.get_GDSII_prisms(silicon, gdsII_file, UPPER_BRANCH_LAYER, si_zmin, si_zmax)
lower_branch = mp.get_GDSII_prisms(silicon, gdsII_file, LOWER_BRANCH_LAYER, si_zmin, si_zmax)
cell = mp.GDSII_vol(gdsII_file, CELL_LAYER, cell_zmin, cell_zmax)
p1 = mp.GDSII_vol(gdsII_file, PORT1_LAYER, si_zmin, si_zmax)
p2 = mp.GDSII_vol(gdsII_file, PORT2_LAYER, si_zmin, si_zmax)
p3 = mp.GDSII_vol(gdsII_file, PORT3_LAYER, si_zmin, si_zmax)
p4 = mp.GDSII_vol(gdsII_file, PORT4_LAYER, si_zmin, si_zmax)
src_vol = mp.GDSII_vol(gdsII_file, SOURCE_LAYER, si_zmin, si_zmax)
# displace upper and lower branches of coupler (as well as source and flux regions)
if d != default_d:
delta_y = 0.5*(d-default_d)
delta = mp.Vector3(y=delta_y)
p1.center += delta
p2.center -= delta
p3.center += delta
p4.center -= delta
src_vol.center += delta
cell.size += 2*delta
for np in range(len(lower_branch)):
lower_branch[np].center -= delta
for nv in range(len(lower_branch[np].vertices)):
lower_branch[np].vertices[nv] -= delta
for np in range(len(upper_branch)):
upper_branch[np].center += delta
for nv in range(len(upper_branch[np].vertices)):
upper_branch[np].vertices[nv] += delta
geometry = upper_branch+lower_branch
if three_d:
oxide_center = mp.Vector3(z=-0.5*t_oxide)
oxide_size = mp.Vector3(cell.size.x,cell.size.y,t_oxide)
oxide_layer = [mp.Block(material=oxide, center=oxide_center, size=oxide_size)]
geometry = geometry+oxide_layer
sources = [mp.EigenModeSource(src=mp.GaussianSource(fcen,fwidth=df),
size=src_vol.size,
center=src_vol.center,
eig_band=1,
eig_parity=mp.NO_PARITY if three_d else mp.EVEN_Y+mp.ODD_Z,
eig_match_freq=True)]
sim = mp.Simulation(resolution=res,
cell_size=cell.size,
boundary_layers=[mp.PML(dpml)],
sources=sources,
geometry=geometry)
mode1 = sim.add_mode_monitor(fcen, 0, 1, mp.ModeRegion(volume=p1))
mode2 = sim.add_mode_monitor(fcen, 0, 1, mp.ModeRegion(volume=p2))
mode3 = sim.add_mode_monitor(fcen, 0, 1, mp.ModeRegion(volume=p3))
mode4 = sim.add_mode_monitor(fcen, 0, 1, mp.ModeRegion(volume=p4))
sim.run(until_after_sources=100)
----------- Initializing structure... Working in 2D dimensions. Computational cell is 34.4 x 7.82 x 0 with resolution 50 prism, center = (-9.09425,1.32149,0) height 0, axis (0,0,1), 174 vertices: (-4,0.06,0) (-4.108,0.061,0) (-4.215,0.062,0) (-4.322,0.065,0) (-4.429,0.07,0) (-4.535,0.075,0) (-4.641,0.081,0) (-4.747,0.089,0) (-4.852,0.097,0) (-5.062,0.117,0) (-5.167,0.129,0) (-5.271,0.141,0) (-5.479,0.169,0) (-5.582,0.184,0) (-5.685,0.2,0) (-5.788,0.217,0) (-5.891,0.235,0) (-5.993,0.253,0) (-6.095,0.272,0) (-6.197,0.292,0) (-6.299,0.313,0) (-6.4,0.334,0) (-6.501,0.356,0) (-6.703,0.402,0) (-6.803,0.425,0) (-6.904,0.45,0) (-7.204,0.525,0) (-7.303,0.552,0) (-7.403,0.578,0) (-7.502,0.605,0) (-7.601,0.633,0) (-7.7,0.66,0) (-7.799,0.688,0) (-7.898,0.717,0) (-7.996,0.745,0) (-8.095,0.774,0) (-8.193,0.803,0) (-8.291,0.833,0) (-8.389,0.862,0) (-8.487,0.892,0) (-8.585,0.921,0) (-8.781,0.981,0) (-8.878,1.011,0) (-9.074,1.071,0) (-9.268,1.131,0) (-9.366,1.161,0) (-9.463,1.191,0) (-9.56,1.22,0) (-9.658,1.25,0) (-9.949,1.337,0) (-10.047,1.366,0) (-10.144,1.395,0) (-10.338,1.451,0) (-10.436,1.478,0) (-10.533,1.506,0) (-10.63,1.533,0) (-10.728,1.559,0) (-10.922,1.611,0) (-11.02,1.636,0) (-11.117,1.66,0) (-11.215,1.685,0) (-11.313,1.708,0) (-11.41,1.731,0) (-11.508,1.754,0) (-11.606,1.776,0) (-11.704,1.797,0) (-11.802,1.817,0) (-11.901,1.837,0) (-11.999,1.856,0) (-12.098,1.875,0) (-12.196,1.893,0) (-12.295,1.91,0) (-12.394,1.926,0) (-12.592,1.956,0) (-12.691,1.97,0) (-12.791,1.982,0) (-12.89,1.994,0) (-12.99,2.005,0) (-13.19,2.025,0) (-13.291,2.033,0) (-13.392,2.04,0) (-13.492,2.046,0) (-13.593,2.051,0) (-13.695,2.055,0) (-13.796,2.058,0) (-14,2.06,0) (-17.2,2.06,0) (-17.2,2.56,0) (-14,2.56,0) (-13.892,2.559,0) (-13.785,2.558,0) (-13.678,2.555,0) (-13.571,2.55,0) (-13.465,2.545,0) (-13.359,2.539,0) (-13.253,2.531,0) (-13.148,2.523,0) (-12.938,2.503,0) (-12.833,2.491,0) (-12.729,2.479,0) (-12.521,2.451,0) (-12.418,2.436,0) (-12.315,2.42,0) (-12.212,2.403,0) (-12.109,2.385,0) (-12.007,2.367,0) (-11.905,2.348,0) (-11.803,2.328,0) (-11.701,2.307,0) (-11.6,2.286,0) (-11.499,2.264,0) (-11.297,2.218,0) (-11.197,2.195,0) (-11.096,2.17,0) (-10.796,2.095,0) (-10.697,2.068,0) (-10.597,2.042,0) (-10.498,2.015,0) (-10.399,1.987,0) (-10.3,1.96,0) (-10.201,1.932,0) (-10.102,1.903,0) (-10.004,1.875,0) (-9.905,1.846,0) (-9.807,1.817,0) (-9.709,1.787,0) (-9.611,1.758,0) (-9.513,1.728,0) (-9.415,1.699,0) (-9.219,1.639,0) (-9.122,1.609,0) (-8.926,1.549,0) (-8.732,1.489,0) (-8.634,1.459,0) (-8.537,1.429,0) (-8.44,1.4,0) (-8.342,1.37,0) (-8.051,1.283,0) (-7.953,1.254,0) (-7.856,1.225,0) (-7.662,1.169,0) (-7.564,1.142,0) (-7.467,1.114,0) (-7.37,1.087,0) (-7.272,1.061,0) (-7.078,1.009,0) (-6.98,0.984,0) (-6.883,0.96,0) (-6.785,0.935,0) (-6.687,0.912,0) (-6.59,0.889,0) (-6.492,0.866,0) (-6.394,0.844,0) (-6.296,0.823,0) (-6.198,0.803,0) (-6.099,0.783,0) (-6.001,0.764,0) (-5.902,0.745,0) (-5.804,0.727,0) (-5.705,0.71,0) (-5.606,0.694,0) (-5.408,0.664,0) (-5.309,0.65,0) (-5.209,0.638,0) (-5.11,0.626,0) (-5.01,0.615,0) (-4.81,0.595,0) (-4.709,0.587,0) (-4.608,0.58,0) (-4.508,0.574,0) (-4.407,0.569,0) (-4.305,0.565,0) (-4.204,0.562,0) (-4,0.56,0) dielectric constant epsilon diagonal = (12,12,12) prism, center = (9.09425,1.32149,0) height 0, axis (0,0,1), 174 vertices: (4,0.06,0) (4,0.56,0) (4.204,0.562,0) (4.305,0.565,0) (4.407,0.569,0) (4.508,0.574,0) (4.608,0.58,0) (4.709,0.587,0) (4.81,0.595,0) (5.01,0.615,0) (5.11,0.626,0) (5.209,0.638,0) (5.309,0.65,0) (5.408,0.664,0) (5.606,0.694,0) (5.705,0.71,0) (5.804,0.727,0) (5.902,0.745,0) (6.001,0.764,0) (6.099,0.783,0) (6.198,0.803,0) (6.296,0.823,0) (6.394,0.844,0) (6.492,0.866,0) (6.59,0.889,0) (6.687,0.912,0) (6.785,0.935,0) (6.883,0.96,0) (6.98,0.984,0) (7.078,1.009,0) (7.272,1.061,0) (7.37,1.087,0) (7.467,1.114,0) (7.564,1.142,0) (7.662,1.169,0) (7.856,1.225,0) (7.953,1.254,0) (8.051,1.283,0) (8.342,1.37,0) (8.44,1.4,0) (8.537,1.429,0) (8.634,1.459,0) (8.732,1.489,0) (8.926,1.549,0) (9.122,1.609,0) (9.219,1.639,0) (9.415,1.699,0) (9.513,1.728,0) (9.611,1.758,0) (9.709,1.787,0) (9.807,1.817,0) (9.905,1.846,0) (10.004,1.875,0) (10.102,1.903,0) (10.201,1.932,0) (10.3,1.96,0) (10.399,1.987,0) (10.498,2.015,0) (10.597,2.042,0) (10.697,2.068,0) (10.796,2.095,0) (11.096,2.17,0) (11.197,2.195,0) (11.297,2.218,0) (11.499,2.264,0) (11.6,2.286,0) (11.701,2.307,0) (11.803,2.328,0) (11.905,2.348,0) (12.007,2.367,0) (12.109,2.385,0) (12.212,2.403,0) (12.315,2.42,0) (12.418,2.436,0) (12.521,2.451,0) (12.729,2.479,0) (12.833,2.491,0) (12.938,2.503,0) (13.148,2.523,0) (13.253,2.531,0) (13.359,2.539,0) (13.465,2.545,0) (13.571,2.55,0) (13.678,2.555,0) (13.785,2.558,0) (13.892,2.559,0) (14,2.56,0) (17.2,2.56,0) (17.2,2.06,0) (14,2.06,0) (13.796,2.058,0) (13.695,2.055,0) (13.593,2.051,0) (13.492,2.046,0) (13.392,2.04,0) (13.291,2.033,0) (13.19,2.025,0) (12.99,2.005,0) (12.89,1.994,0) (12.791,1.982,0) (12.691,1.97,0) (12.592,1.956,0) (12.394,1.926,0) (12.295,1.91,0) (12.196,1.893,0) (12.098,1.875,0) (11.999,1.856,0) (11.901,1.837,0) (11.802,1.817,0) (11.704,1.797,0) (11.606,1.776,0) (11.508,1.754,0) (11.41,1.731,0) (11.313,1.708,0) (11.215,1.685,0) (11.117,1.66,0) (11.02,1.636,0) (10.922,1.611,0) (10.728,1.559,0) (10.63,1.533,0) (10.533,1.506,0) (10.436,1.478,0) (10.338,1.451,0) (10.144,1.395,0) (10.047,1.366,0) (9.949,1.337,0) (9.658,1.25,0) (9.56,1.22,0) (9.463,1.191,0) (9.366,1.161,0) (9.268,1.131,0) (9.074,1.071,0) (8.878,1.011,0) (8.781,0.981,0) (8.585,0.921,0) (8.487,0.892,0) (8.389,0.862,0) (8.291,0.833,0) (8.193,0.803,0) (8.095,0.774,0) (7.996,0.745,0) (7.898,0.717,0) (7.799,0.688,0) (7.7,0.66,0) (7.601,0.633,0) (7.502,0.605,0) (7.403,0.578,0) (7.303,0.552,0) (7.204,0.525,0) (6.904,0.45,0) (6.803,0.425,0) (6.703,0.402,0) (6.501,0.356,0) (6.4,0.334,0) (6.299,0.313,0) (6.197,0.292,0) (6.095,0.272,0) (5.993,0.253,0) (5.891,0.235,0) (5.788,0.217,0) (5.685,0.2,0) (5.582,0.184,0) (5.479,0.169,0) (5.271,0.141,0) (5.167,0.129,0) (5.062,0.117,0) (4.852,0.097,0) (4.747,0.089,0) (4.641,0.081,0) (4.535,0.075,0) (4.429,0.07,0) (4.322,0.065,0) (4.215,0.062,0) (4.108,0.061,0) dielectric constant epsilon diagonal = (12,12,12) prism, center = (0,0.31,0) height 0, axis (0,0,1), 4 vertices: (-4,0.06,0) (-4,0.56,0) (4,0.56,0) (4,0.06,0) dielectric constant epsilon diagonal = (12,12,12) prism, center = (-9.09425,-1.32149,0) height 0, axis (0,0,1), 174 vertices: (-17.2,-2.56,0) (-17.2,-2.06,0) (-14,-2.06,0) (-13.796,-2.058,0) (-13.695,-2.055,0) (-13.593,-2.051,0) (-13.492,-2.046,0) (-13.392,-2.04,0) (-13.291,-2.033,0) (-13.19,-2.025,0) (-12.99,-2.005,0) (-12.89,-1.994,0) (-12.791,-1.982,0) (-12.691,-1.97,0) (-12.592,-1.956,0) (-12.394,-1.926,0) (-12.295,-1.91,0) (-12.196,-1.893,0) (-12.098,-1.875,0) (-11.999,-1.856,0) (-11.901,-1.837,0) (-11.802,-1.817,0) (-11.704,-1.797,0) (-11.606,-1.776,0) (-11.508,-1.754,0) (-11.41,-1.731,0) (-11.313,-1.708,0) (-11.215,-1.685,0) (-11.117,-1.66,0) (-11.02,-1.636,0) (-10.922,-1.611,0) (-10.728,-1.559,0) (-10.63,-1.533,0) (-10.533,-1.506,0) (-10.436,-1.478,0) (-10.338,-1.451,0) (-10.144,-1.395,0) (-10.047,-1.366,0) (-9.949,-1.337,0) (-9.658,-1.25,0) (-9.56,-1.22,0) (-9.463,-1.191,0) (-9.366,-1.161,0) (-9.268,-1.131,0) (-9.074,-1.071,0) (-8.878,-1.011,0) (-8.781,-0.981,0) (-8.585,-0.921,0) (-8.487,-0.892,0) (-8.389,-0.862,0) (-8.291,-0.833,0) (-8.193,-0.803,0) (-8.095,-0.774,0) (-7.996,-0.745,0) (-7.898,-0.717,0) (-7.799,-0.688,0) (-7.7,-0.66,0) (-7.601,-0.633,0) (-7.502,-0.605,0) (-7.403,-0.578,0) (-7.303,-0.552,0) (-7.204,-0.525,0) (-6.904,-0.45,0) (-6.803,-0.425,0) (-6.703,-0.402,0) (-6.501,-0.356,0) (-6.4,-0.334,0) (-6.299,-0.313,0) (-6.197,-0.292,0) (-6.095,-0.272,0) (-5.993,-0.253,0) (-5.891,-0.235,0) (-5.788,-0.217,0) (-5.685,-0.2,0) (-5.582,-0.184,0) (-5.479,-0.169,0) (-5.271,-0.141,0) (-5.167,-0.129,0) (-5.062,-0.117,0) (-4.852,-0.097,0) (-4.747,-0.089,0) (-4.641,-0.081,0) (-4.535,-0.075,0) (-4.429,-0.07,0) (-4.322,-0.065,0) (-4.215,-0.062,0) (-4.108,-0.061,0) (-4,-0.06,0) (-4,-0.56,0) (-4.204,-0.562,0) (-4.305,-0.565,0) (-4.407,-0.569,0) (-4.508,-0.574,0) (-4.608,-0.58,0) (-4.709,-0.587,0) (-4.81,-0.595,0) (-5.01,-0.615,0) (-5.11,-0.626,0) (-5.209,-0.638,0) (-5.309,-0.65,0) (-5.408,-0.664,0) (-5.606,-0.694,0) (-5.705,-0.71,0) (-5.804,-0.727,0) (-5.902,-0.745,0) (-6.001,-0.764,0) (-6.099,-0.783,0) (-6.198,-0.803,0) (-6.296,-0.823,0) (-6.394,-0.844,0) (-6.492,-0.866,0) (-6.59,-0.889,0) (-6.687,-0.912,0) (-6.785,-0.935,0) (-6.883,-0.96,0) (-6.98,-0.984,0) (-7.078,-1.009,0) (-7.272,-1.061,0) (-7.37,-1.087,0) (-7.467,-1.114,0) (-7.564,-1.142,0) (-7.662,-1.169,0) (-7.856,-1.225,0) (-7.953,-1.254,0) (-8.051,-1.283,0) (-8.342,-1.37,0) (-8.44,-1.4,0) (-8.537,-1.429,0) (-8.634,-1.459,0) (-8.732,-1.489,0) (-8.926,-1.549,0) (-9.122,-1.609,0) (-9.219,-1.639,0) (-9.415,-1.699,0) (-9.513,-1.728,0) (-9.611,-1.758,0) (-9.709,-1.787,0) (-9.807,-1.817,0) (-9.905,-1.846,0) (-10.004,-1.875,0) (-10.102,-1.903,0) (-10.201,-1.932,0) (-10.3,-1.96,0) (-10.399,-1.987,0) (-10.498,-2.015,0) (-10.597,-2.042,0) (-10.697,-2.068,0) (-10.796,-2.095,0) (-11.096,-2.17,0) (-11.197,-2.195,0) (-11.297,-2.218,0) (-11.499,-2.264,0) (-11.6,-2.286,0) (-11.701,-2.307,0) (-11.803,-2.328,0) (-11.905,-2.348,0) (-12.007,-2.367,0) (-12.109,-2.385,0) (-12.212,-2.403,0) (-12.315,-2.42,0) (-12.418,-2.436,0) (-12.521,-2.451,0) (-12.729,-2.479,0) (-12.833,-2.491,0) (-12.938,-2.503,0) (-13.148,-2.523,0) (-13.253,-2.531,0) (-13.359,-2.539,0) (-13.465,-2.545,0) (-13.571,-2.55,0) (-13.678,-2.555,0) (-13.785,-2.558,0) (-13.892,-2.559,0) (-14,-2.56,0) dielectric constant epsilon diagonal = (12,12,12) prism, center = (9.09425,-1.32149,0) height 0, axis (0,0,1), 174 vertices: (14,-2.56,0) (13.892,-2.559,0) (13.785,-2.558,0) (13.678,-2.555,0) (13.571,-2.55,0) (13.465,-2.545,0) (13.359,-2.539,0) (13.253,-2.531,0) (13.148,-2.523,0) (12.938,-2.503,0) (12.833,-2.491,0) (12.729,-2.479,0) (12.521,-2.451,0) (12.418,-2.436,0) (12.315,-2.42,0) (12.212,-2.403,0) (12.109,-2.385,0) (12.007,-2.367,0) (11.905,-2.348,0) (11.803,-2.328,0) (11.701,-2.307,0) (11.6,-2.286,0) (11.499,-2.264,0) (11.297,-2.218,0) (11.197,-2.195,0) (11.096,-2.17,0) (10.796,-2.095,0) (10.697,-2.068,0) (10.597,-2.042,0) (10.498,-2.015,0) (10.399,-1.987,0) (10.3,-1.96,0) (10.201,-1.932,0) (10.102,-1.903,0) (10.004,-1.875,0) (9.905,-1.846,0) (9.807,-1.817,0) (9.709,-1.787,0) (9.611,-1.758,0) (9.513,-1.728,0) (9.415,-1.699,0) (9.219,-1.639,0) (9.122,-1.609,0) (8.926,-1.549,0) (8.732,-1.489,0) (8.634,-1.459,0) (8.537,-1.429,0) (8.44,-1.4,0) (8.342,-1.37,0) (8.051,-1.283,0) (7.953,-1.254,0) (7.856,-1.225,0) (7.662,-1.169,0) (7.564,-1.142,0) (7.467,-1.114,0) (7.37,-1.087,0) (7.272,-1.061,0) (7.078,-1.009,0) (6.98,-0.984,0) (6.883,-0.96,0) (6.785,-0.935,0) (6.687,-0.912,0) (6.59,-0.889,0) (6.492,-0.866,0) (6.394,-0.844,0) (6.296,-0.823,0) (6.198,-0.803,0) (6.099,-0.783,0) (6.001,-0.764,0) (5.902,-0.745,0) (5.804,-0.727,0) (5.705,-0.71,0) (5.606,-0.694,0) (5.408,-0.664,0) (5.309,-0.65,0) (5.209,-0.638,0) (5.11,-0.626,0) (5.01,-0.615,0) (4.81,-0.595,0) (4.709,-0.587,0) (4.608,-0.58,0) (4.508,-0.574,0) (4.407,-0.569,0) (4.305,-0.565,0) (4.204,-0.562,0) (4,-0.56,0) (4,-0.06,0) (4.108,-0.061,0) (4.215,-0.062,0) (4.322,-0.065,0) (4.429,-0.07,0) (4.535,-0.075,0) (4.641,-0.081,0) (4.747,-0.089,0) (4.852,-0.097,0) (5.062,-0.117,0) (5.167,-0.129,0) (5.271,-0.141,0) (5.479,-0.169,0) (5.582,-0.184,0) (5.685,-0.2,0) (5.788,-0.217,0) (5.891,-0.235,0) (5.993,-0.253,0) (6.095,-0.272,0) (6.197,-0.292,0) (6.299,-0.313,0) (6.4,-0.334,0) (6.501,-0.356,0) (6.703,-0.402,0) (6.803,-0.425,0) (6.904,-0.45,0) (7.204,-0.525,0) (7.303,-0.552,0) (7.403,-0.578,0) (7.502,-0.605,0) (7.601,-0.633,0) (7.7,-0.66,0) (7.799,-0.688,0) (7.898,-0.717,0) (7.996,-0.745,0) (8.095,-0.774,0) (8.193,-0.803,0) (8.291,-0.833,0) (8.389,-0.862,0) (8.487,-0.892,0) (8.585,-0.921,0) (8.781,-0.981,0) (8.878,-1.011,0) (9.074,-1.071,0) (9.268,-1.131,0) (9.366,-1.161,0) (9.463,-1.191,0) (9.56,-1.22,0) (9.658,-1.25,0) (9.949,-1.337,0) (10.047,-1.366,0) (10.144,-1.395,0) (10.338,-1.451,0) (10.436,-1.478,0) (10.533,-1.506,0) (10.63,-1.533,0) (10.728,-1.559,0) (10.922,-1.611,0) (11.02,-1.636,0) (11.117,-1.66,0) (11.215,-1.685,0) (11.313,-1.708,0) (11.41,-1.731,0) (11.508,-1.754,0) (11.606,-1.776,0) (11.704,-1.797,0) (11.802,-1.817,0) (11.901,-1.837,0) (11.999,-1.856,0) (12.098,-1.875,0) (12.196,-1.893,0) (12.295,-1.91,0) (12.394,-1.926,0) (12.592,-1.956,0) (12.691,-1.97,0) (12.791,-1.982,0) (12.89,-1.994,0) (12.99,-2.005,0) (13.19,-2.025,0) (13.291,-2.033,0) (13.392,-2.04,0) (13.492,-2.046,0) (13.593,-2.051,0) (13.695,-2.055,0) (13.796,-2.058,0) (14,-2.06,0) (17.2,-2.06,0) (17.2,-2.56,0) dielectric constant epsilon diagonal = (12,12,12) prism, center = (0,-0.31,0) height 0, axis (0,0,1), 4 vertices: (-4,-0.56,0) (-4,-0.06,0) (4,-0.06,0) (4,-0.56,0) dielectric constant epsilon diagonal = (12,12,12) subpixel-averaging is 11.5011% done, 30.9631 s remaining subpixel-averaging is 20.4466% done, 15.8135 s remaining subpixel-averaging is 29.3921% done, 9.71826 s remaining subpixel-averaging is 58.5715% done, 2.83111 s remaining subpixel-averaging is 71.3508% done, 1.62492 s remaining subpixel-averaging is 80.2963% done, 0.998215 s remaining subpixel-averaging is 89.2418% done, 0.492406 s remaining subpixel-averaging is 11.5011% done, 30.9076 s remaining subpixel-averaging is 20.4466% done, 15.9068 s remaining subpixel-averaging is 29.3921% done, 9.80747 s remaining subpixel-averaging is 57.7195% done, 2.9313 s remaining subpixel-averaging is 71.1378% done, 1.66162 s remaining subpixel-averaging is 80.0833% done, 1.0166 s remaining subpixel-averaging is 88.8158% done, 0.504616 s remaining subpixel-averaging is 11.5011% done, 30.9943 s remaining subpixel-averaging is 20.4466% done, 15.8803 s remaining subpixel-averaging is 29.3921% done, 9.78632 s remaining subpixel-averaging is 56.8676% done, 3.03894 s remaining subpixel-averaging is 70.9248% done, 1.64887 s remaining subpixel-averaging is 79.8703% done, 1.02929 s remaining subpixel-averaging is 88.8158% done, 0.514541 s remaining time for set_epsilon = 102.547 s ----------- MPB solved for omega_1(2.2349,0,0) = 0.686107 after 12 iters MPB solved for omega_1(2.08827,0,0) = 0.645201 after 7 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 4 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 1 iters Meep progress: 5.5600000000000005/177.5 = 3.1% done in 4.0s, 123.7s to go on time step 556 (time=5.56), 0.00719569 s/step Meep progress: 10.59/177.5 = 6.0% done in 8.0s, 126.1s to go on time step 1059 (time=10.59), 0.00795647 s/step Meep progress: 15.25/177.5 = 8.6% done in 12.0s, 127.8s to go on time step 1525 (time=15.25), 0.0086014 s/step Meep progress: 19.76/177.5 = 11.1% done in 16.0s, 127.8s to go on time step 1976 (time=19.76), 0.00887541 s/step Meep progress: 24.02/177.5 = 13.5% done in 20.0s, 127.9s to go on time step 2402 (time=24.02), 0.00940361 s/step Meep progress: 29.82/177.5 = 16.8% done in 24.0s, 119.0s to go on time step 2982 (time=29.82), 0.00690337 s/step Meep progress: 35.64/177.5 = 20.1% done in 28.0s, 111.6s to go on time step 3564 (time=35.64), 0.00688025 s/step Meep progress: 41.43/177.5 = 23.3% done in 32.0s, 105.2s to go on time step 4143 (time=41.43), 0.00691613 s/step Meep progress: 47.25/177.5 = 26.6% done in 36.0s, 99.3s to go on time step 4725 (time=47.25), 0.0068737 s/step Meep progress: 52.99/177.5 = 29.9% done in 40.0s, 94.1s to go on time step 5299 (time=52.99), 0.0069787 s/step Meep progress: 58.68/177.5 = 33.1% done in 44.0s, 89.2s to go on time step 5868 (time=58.68), 0.00703139 s/step Meep progress: 64.46000000000001/177.5 = 36.3% done in 48.0s, 84.2s to go on time step 6446 (time=64.46), 0.0069234 s/step Meep progress: 70.27/177.5 = 39.6% done in 52.0s, 79.4s to go on time step 7027 (time=70.27), 0.00689556 s/step Meep progress: 76.01/177.5 = 42.8% done in 56.1s, 74.8s to go on time step 7601 (time=76.01), 0.00698002 s/step Meep progress: 81.79/177.5 = 46.1% done in 60.1s, 70.3s to go on time step 8179 (time=81.79), 0.00692945 s/step Meep progress: 87.48/177.5 = 49.3% done in 64.1s, 65.9s to go on time step 8748 (time=87.48), 0.00703933 s/step Meep progress: 93.25/177.5 = 52.5% done in 68.1s, 61.5s to go on time step 9325 (time=93.25), 0.00693846 s/step Meep progress: 99.02/177.5 = 55.8% done in 72.1s, 57.1s to go on time step 9902 (time=99.02), 0.00693694 s/step Meep progress: 104.81/177.5 = 59.0% done in 76.1s, 52.8s to go on time step 10481 (time=104.81), 0.00691733 s/step Meep progress: 110.57000000000001/177.5 = 62.3% done in 80.1s, 48.5s to go on time step 11057 (time=110.57), 0.00695527 s/step Meep progress: 116.32000000000001/177.5 = 65.5% done in 84.1s, 44.2s to go on time step 11632 (time=116.32), 0.00696429 s/step Meep progress: 122.05/177.5 = 68.8% done in 88.1s, 40.0s to go on time step 12205 (time=122.05), 0.00698096 s/step Meep progress: 127.81/177.5 = 72.0% done in 92.1s, 35.8s to go on time step 12781 (time=127.81), 0.00694492 s/step Meep progress: 133.59/177.5 = 75.3% done in 96.1s, 31.6s to go on time step 13359 (time=133.59), 0.00692692 s/step Meep progress: 139.34/177.5 = 78.5% done in 100.1s, 27.4s to go on time step 13934 (time=139.34), 0.00695738 s/step Meep progress: 145.11/177.5 = 81.8% done in 104.1s, 23.2s to go on time step 14511 (time=145.11), 0.00694059 s/step Meep progress: 150.94/177.5 = 85.0% done in 108.1s, 19.0s to go on time step 15094 (time=150.94), 0.00686312 s/step Meep progress: 156.67000000000002/177.5 = 88.3% done in 112.1s, 14.9s to go on time step 15667 (time=156.67), 0.00698277 s/step Meep progress: 162.42000000000002/177.5 = 91.5% done in 116.1s, 10.8s to go on time step 16242 (time=162.42), 0.00696221 s/step Meep progress: 168.18/177.5 = 94.7% done in 120.1s, 6.7s to go on time step 16818 (time=168.18), 0.00695355 s/step Meep progress: 173.85/177.5 = 97.9% done in 124.1s, 2.6s to go on time step 17385 (time=173.85), 0.00706013 s/step run 0 finished at t = 177.5 (17750 timesteps)
# S parameters
p1_coeff = sim.get_eigenmode_coefficients(mode1, [1], eig_parity=mp.NO_PARITY if three_d else mp.EVEN_Y+mp.ODD_Z).alpha[0,0,0]
p2_coeff = sim.get_eigenmode_coefficients(mode2, [1], eig_parity=mp.NO_PARITY if three_d else mp.EVEN_Y+mp.ODD_Z).alpha[0,0,1]
p3_coeff = sim.get_eigenmode_coefficients(mode3, [1], eig_parity=mp.NO_PARITY if three_d else mp.EVEN_Y+mp.ODD_Z).alpha[0,0,0]
p4_coeff = sim.get_eigenmode_coefficients(mode4, [1], eig_parity=mp.NO_PARITY if three_d else mp.EVEN_Y+mp.ODD_Z).alpha[0,0,0]
# transmittance
p2_trans = abs(p2_coeff)**2/abs(p1_coeff)**2
p3_trans = abs(p3_coeff)**2/abs(p1_coeff)**2
p4_trans = abs(p4_coeff)**2/abs(p1_coeff)**2
print("trans:, {:.2f}, {:.6f}, {:.6f}, {:.6f}".format(d,p2_trans,p3_trans,p4_trans))
MPB solved for omega_1(2.2349,0,0) = 0.686107 after 12 iters MPB solved for omega_1(2.08827,0,0) = 0.645201 after 7 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 4 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 1 iters Dominant planewave for band 1: (2.088124,-0.000000,0.000000) MPB solved for omega_1(2.2349,0,0) = 0.686107 after 12 iters MPB solved for omega_1(2.08827,0,0) = 0.645201 after 7 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 4 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 1 iters Dominant planewave for band 1: (2.088124,-0.000000,0.000000) MPB solved for omega_1(2.2349,0,0) = 0.686107 after 13 iters MPB solved for omega_1(2.08827,0,0) = 0.645201 after 7 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 4 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 1 iters Dominant planewave for band 1: (2.088124,-0.000000,0.000000) MPB solved for omega_1(2.2349,0,0) = 0.686107 after 14 iters MPB solved for omega_1(2.08827,0,0) = 0.645201 after 7 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 4 iters MPB solved for omega_1(2.08812,0,0) = 0.645161 after 1 iters Dominant planewave for band 1: (2.088124,-0.000000,0.000000) trans:, 0.12, 0.000002, 0.557065, 0.426816
For a given waveguide separation distance (d
), the simulation computes the transmittance of Ports 2, 3, and 4. The transmittance is the square of the S-parameter which is equivalent to the mode coefficient. There is an additional mode monitor at Port 1 to compute the input power from the adjacent eigenmode source; this is used for normalization when computing the transmittance. The eight layers of the GDSII file are each converted to a Simulation
object: the upper and lower branches of the coupler are defined as a collection of Prism
s, the rectilinear regions of the source and flux monitor as a Volume
and FluxRegion
. The size of the cell in the $y$ direction is dependent on d
. The default dimensionality is 2d. An optional input parameter (three_d
) converts the geometry to 3d by extruding the coupler geometry in the z direction and adding an oxide layer beneath similar to a silicon on insulator (SOI) substrate. A schematic of the coupler design in 3d generated using MayaVi is shown below.
The transmittance results are plotted in the figure below for a range of separation distances from 0.02 to 0.30 μm with increments of 0.02 μm. When the two waveguide branches are sufficiently separated (d
> 0.2 μm), practically all of the input power remains in the top branch and is transferred to Port 3. A small amount of the input power is lost due to scattering into radiative modes within the light cone in the tapered sections where the translational symmetry of the waveguide is broken. This is why the power in Port 3 never reaches exactly 100%. For separation distances of less than approximately 0.2 μm, evanescent coupling of the modes from the top to the lower branch begins to transfer some of the input power to Port 4. For d
of 0.13 μm, the input signal is split evenly into Ports 3 and 4. For d
of 0.06 μm, the input power is transferred completely to Port 4. Finally, for d
of less than 0.06 μm, the evanescent coupling becomes rapidly ineffective and the signal again remains mostly in Port 3. Note that there is never any power in Port 2 given its location relative to the input from Port 1.
These quantitative results can also be verified qualitatively using the field profiles shown below for d
of 0.06, 0.13, and 0.30 μm. To generate these images, the pulse source is replaced with a continuous wave (CW) and the fields are time stepped for a sufficiently long run time until they have reached steady state. The array slicing routines get_epsilon
and get_efield_z
are then used to obtain the dielectric and field data over the entire cell.
sim.reset_meep()
sources = [mp.EigenModeSource(src=mp.ContinuousSource(fcen,fwidth=df),
size=src_vol.size,
center=src_vol.center,
eig_band=1,
eig_parity=mp.EVEN_Y+mp.ODD_Z,
eig_match_freq=True)]
sim = mp.Simulation(resolution=res,
cell_size=cell.size,
boundary_layers=[mp.PML(dpml)],
sources=sources,
geometry=geometry)
sim.run(until=400) # arbitrary long run time to ensure that fields have reached steady state
eps_data = sim.get_epsilon()
ez_data = numpy.real(sim.get_efield_z())
plt.figure(dpi=200)
plt.imshow(numpy.transpose(eps_data), interpolation='spline36', cmap='binary')
plt.imshow(numpy.flipud(numpy.transpose(ez_data)), interpolation='spline36', cmap='RdBu', alpha=0.9)
plt.axis('off')
plt.show()
The field profiles confirm that for d
of 0.06 μm (Figure 1), the input signal in Port 1 of the top branch is almost completely transferred to Port 4 of the bottom branch. For d
of 0.13 μm (Figure 2), the input signal is split evenly between the two branches. Finally, for d
of 0.30 μm (Figure 3), there is no longer any evanescent coupling and the signal remains completely in the top branch. Note the absence of the fields in the PML regions of Ports 3 and 4.