Mode solving anisotropic

Contents

Mode solving anisotropic#

[1]

Caution

This example is under construction As Julia-Dicts are not ordered, the mesh might become incorrect when adjusted (for now, better do the meshing in python)

using PyCall
np = pyimport("numpy")
shapely = pyimport("shapely")
shapely.affinity = pyimport("shapely.affinity")
clip_by_rect = pyimport("shapely.ops").clip_by_rect
OrderedDict = pyimport("collections").OrderedDict
mesh_from_OrderedDict = pyimport("femwell.mesh").mesh_from_OrderedDict

wg_width = 0.8
wg_thickness = 0.6076
core = shapely.geometry.box(-wg_width / 2, 0, +wg_width / 2, wg_thickness)
env = shapely.affinity.scale(core.buffer(5, resolution = 8), xfact = 0.5)

polygons = OrderedDict(
    core = core,
    box = clip_by_rect(env, -np.inf, -np.inf, np.inf, 0),
    clad = clip_by_rect(env, -np.inf, 0, np.inf, np.inf),
)

resolutions = Dict("core" => Dict("resolution" => 0.03, "distance" => 0.5))

mesh = mesh_from_OrderedDict(
    polygons,
    resolutions,
    default_resolution_max = 0.1,
    filename = "mesh.msh",
)
Hide code cell output
PyObject <meshio mesh object>
  Number of points: 10536
  Number of cells:
    line: 27
    line: 21
    line: 27
    line: 21
    line: 27
    line: 27
    triangle: 1442
    triangle: 7987
    triangle: 11363
  Cell sets: core___box, core___clad, box___clad, core, box, clad, gmsh:bounding_entities
  Point data: gmsh:dim_tags
  Cell data: gmsh:physical, gmsh:geometrical
  Field data: core___box, core___clad, box___clad, core, box, clad
using LinearAlgebra
using SparseArrays

using Gridap
using Gridap.Geometry
using Gridap.Visualization
using Gridap.ReferenceFEs
using GridapGmsh
using GridapMakie, CairoMakie

using Femwell.Maxwell.Waveguide

CairoMakie.inline!(true)
model = GmshDiscreteModel("mesh.msh")
Ω = Triangulation(model)

labels = get_face_labeling(model)

eye = diagonal_tensor(VectorValue([1.0, 1.0, 1.0]))
n2 = 2.302
Δ = 0.005
yig = TensorValue([n2^2 1im*Δ 0; -1im*Δ n2^2 0; 0 0 n2^2])
epsilons = ["core" => yig, "box" => 1.95^2 * eye, "clad" => 1.0^2 * eye]
ε(tag) = Dict(get_tag_from_name(labels, u) => v for (u, v) in epsilons)[tag]

τ = CellField(get_face_tag(labels, num_cell_dims(model)), Ω)

modes = calculate_modes(model, ε  τ, λ = 1.3, num = 2, order = 1)

plot_mode(modes[1])
plot_mode(modes[2])
modes
Info    : Reading 'mesh.msh'...
Info    : 85 entities
Info    : 10536 nodes
Info    : 20942 elements
Info    : Done reading 'mesh.msh'
../_images/66232bda16a5e1797e1de45cae699fe5d6aadd19392e89580f68350a69e3b0ac.png ../_images/4ec2afcd9b353158c05d341c2c619218897edac45360d904d737d185567de901.png
2-element Vector{Mode}:
 Mode(k=9.901995888561764 - 2.5013470641342013e-13im, n_eff=2.048737069782298 - 5.175322745389658e-14im)
 Mode(k=9.89343369466801 + 2.572585582948619e-13im, n_eff=2.046965539655825 + 5.3227162566918326e-14im)

Bibliography#

[1]

Arman B. Fallahkhair, Kai S. Li, and Thomas E. Murphy. Vector finite difference modesolver for anisotropic dielectric waveguides. Journal of Lightwave Technology, 26(11):1423–1431, June 2008. URL: https://doi.org/10.1109/jlt.2008.923643, doi:10.1109/jlt.2008.923643.