Coplanar waveguide - vary width

Contents

Coplanar waveguide - vary width#

In this example we calculate effective epsilon of the coplanar waveguides from [1]

Hide code cell source
from collections import OrderedDict

import matplotlib.pyplot as plt
import numpy as np
import scipy.constants
import shapely
import shapely.ops
from shapely.geometry import LineString, box
from skfem import Basis, ElementTriP0
from skfem.io.meshio import from_meshio
from tqdm import tqdm

from femwell.maxwell.waveguide import compute_modes
from femwell.mesh import mesh_from_OrderedDict
def mesh_waveguide(filename, wsim, hclad, hsi, wcore, hcore):
    core = box(-wcore / 2, -hcore / 2, wcore / 2, hcore / 2)
    clad = box(-wsim / 2, -hcore / 2, wsim / 2, -hcore / 2 + hclad)
    silicon = box(-wsim / 2, -hcore / 2, wsim / 2, -hcore / 2 - hsi)

    combined = shapely.ops.unary_union([clad, silicon])

    polygons = OrderedDict(
        surface=LineString(combined.exterior),
        core_interface=core.buffer(hcore / 10).exterior,
        core=core,
        clad=clad,
        silicon=silicon,
    )

    resolutions = dict(
        core={"resolution": min(wcore / 20, hcore / 4), "distance": 1, "SizeMax": 0.5},
        core_interface={
            "resolution": min(wcore / 20, hcore / 4),
            "distance": 1,
            "SizeMax": 0.5,
        },
        # silicon={"resolution": 0.2, "distance": 5},
    )

    return mesh_from_OrderedDict(polygons, resolutions, filename=filename, default_resolution_max=5)
frequencies = np.linspace(1e9, 16e9, 16)
# widths = (0.04, 0.6, 0.8, 0.1, 0.14, 0.18, 0.25, 0.3, 0.4, 0.5, 0.6, 0.8, 1, 1.2, 1.4, 1.6, 1.8, 2, 2.2, 2.6, 3,)
widths = (
    0.04,
    0.4,
    3,
)
epsilon_effs = np.zeros((len(widths), len(frequencies), 1), dtype=complex)
characteristic_impedances = np.zeros((len(widths), len(frequencies), 1), dtype=complex)

for i, width in enumerate(tqdm(widths)):
    mesh = from_meshio(
        mesh_waveguide(
            filename="mesh.msh",
            wsim=10,
            hclad=10,
            hsi=0.64,
            wcore=width,
            hcore=0.005,
        )
    )
    mesh = mesh.scaled((1e-3,) * 2)

    for j, frequency in enumerate(tqdm(frequencies, leave=False)):
        basis0 = Basis(mesh, ElementTriP0(), intorder=4)
        epsilon = basis0.zeros().astype(complex)
        epsilon[basis0.get_dofs(elements="silicon")] = 9.9 + 0.0005
        epsilon[basis0.get_dofs(elements="clad")] = 1.0
        epsilon[basis0.get_dofs(elements="core")] = 1 - 1j * 1 / (
            18e-6 * 1e-3
        ) / scipy.constants.epsilon_0 / (2 * np.pi * frequency)
        # basis0.plot(np.real(epsilon), colorbar=True).show()

        modes = compute_modes(
            basis0,
            epsilon,
            wavelength=scipy.constants.speed_of_light / frequency,
            num_modes=1,
            metallic_boundaries=True,
        )
        print(f"Width: {width}, Frequency: {frequency/1e9} GHz")
        print(f"Effective epsilons {modes.n_effs**2}")
        # modes[0].show("E", part="real", plot_vectors=True, colorbar=True)

        epsilon_effs[i, j] = modes.n_effs**2

        # In work

        conductors = ("core_interface",)

        from skfem import Functional
        from skfem.helpers import inner

        @Functional(dtype=np.complex64)
        def current_form(w):
            return inner(np.array([w.n[1], -w.n[0]]), w.H)

        currents = np.zeros((len(conductors), len(modes)))

        for mode_i, mode in enumerate(modes):
            # modes[0].show("H", part="real", plot_vectors=True, colorbar=True)

            (ht, ht_basis), (hz, hz_basis) = mode.basis.split(mode.H)
            for conductors_i, conductor in enumerate(conductors):
                facet_basis = ht_basis.boundary(facets=mesh.boundaries[conductor])
                current = abs(current_form.assemble(facet_basis, H=facet_basis.interpolate(ht)))
                currents[conductors_i, mode_i] = current

        characteristic_impedances[i, j] = np.linalg.inv(currents).T @ np.linalg.inv(currents)
        print("characteristic impedances", characteristic_impedances[i, j])
Hide code cell output
  0%|          | 0/3 [00:00<?, ?it/s]

  0%|          | 0/16 [00:00<?, ?it/s]


  6%|▋         | 1/16 [00:00<00:03,  4.03it/s]

Width: 0.04, Frequency: 1.0 GHz
Effective epsilons [5.87625225-0.15324565j]
characteristic impedances [136.33112309+0.j]

 12%|█▎        | 2/16 [00:00<00:03,  4.13it/s]

Width: 0.04, Frequency: 2.0 GHz
Effective epsilons [5.85860486-0.10264146j]
characteristic impedances [137.59970055+0.j]

 19%|█▉        | 3/16 [00:00<00:03,  4.18it/s]

Width: 0.04, Frequency: 3.0 GHz
Effective epsilons [5.8560037-0.08533542j]
characteristic impedances [138.22313856+0.j]

 25%|██▌       | 4/16 [00:00<00:02,  4.20it/s]

Width: 0.04, Frequency: 4.0 GHz
Effective epsilons [5.86064764-0.07579502j]
characteristic impedances [138.65591411+0.j]

 31%|███▏      | 5/16 [00:01<00:02,  4.17it/s]

Width: 0.04, Frequency: 5.0 GHz
Effective epsilons [5.87001526-0.06916168j]
characteristic impedances [139.03849787+0.j]

 38%|███▊      | 6/16 [00:01<00:02,  4.20it/s]

Width: 0.04, Frequency: 6.0 GHz
Effective epsilons [5.88262028-0.06403188j]
characteristic impedances [139.42204232+0.j]

 44%|████▍     | 7/16 [00:01<00:02,  4.22it/s]

Width: 0.04, Frequency: 7.0 GHz
Effective epsilons [5.89747783-0.05986636j]
characteristic impedances [139.82352233+0.j]

 50%|█████     | 8/16 [00:01<00:01,  4.24it/s]

Width: 0.04, Frequency: 8.0 GHz
Effective epsilons [5.91394572-0.05639351j]
characteristic impedances [140.24752605+0.j]

 56%|█████▋    | 9/16 [00:02<00:01,  4.25it/s]

Width: 0.04, Frequency: 9.0 GHz
Effective epsilons [5.93161637-0.05344512j]
characteristic impedances [140.69493616+0.j]

 62%|██████▎   | 10/16 [00:02<00:01,  4.25it/s]

Width: 0.04, Frequency: 10.0 GHz
Effective epsilons [5.95023457-0.05090389j]
characteristic impedances [141.1664123+0.j]

 69%|██████▉   | 11/16 [00:02<00:01,  4.24it/s]

Width: 0.04, Frequency: 11.0 GHz
Effective epsilons [5.96963997-0.04868326j]
characteristic impedances [141.66350401+0.j]

 75%|███████▌  | 12/16 [00:02<00:00,  4.21it/s]

Width: 0.04, Frequency: 12.0 GHz
Effective epsilons [5.98972972-0.04671755j]
characteristic impedances [142.18872407+0.j]

 81%|████████▏ | 13/16 [00:03<00:00,  4.22it/s]

Width: 0.04, Frequency: 13.0 GHz
Effective epsilons [6.01043545-0.0449563j]
characteristic impedances [142.74525716+0.j]

 88%|████████▊ | 14/16 [00:03<00:00,  4.22it/s]

Width: 0.04, Frequency: 14.0 GHz
Effective epsilons [6.0317094-0.04336037j]
characteristic impedances [143.33680844+0.j]

 94%|█████████▍| 15/16 [00:03<00:00,  4.20it/s]

Width: 0.04, Frequency: 15.0 GHz
Effective epsilons [6.05351647-0.04189926j]
characteristic impedances [143.96728503+0.j]

100%|██████████| 16/16 [00:03<00:00,  4.11it/s]


                                               

 33%|███▎      | 1/3 [00:04<00:08,  4.24s/it]
Width: 0.04, Frequency: 16.0 GHz
Effective epsilons [6.07582954-0.04054904j]
characteristic impedances [144.64074887+0.j]

  0%|          | 0/16 [00:00<?, ?it/s]


  6%|▋         | 1/16 [00:01<00:20,  1.39s/it]

Width: 0.4, Frequency: 1.0 GHz
Effective epsilons [6.38175462-0.04265092j]
characteristic impedances [68.07419793+0.j]

 12%|█▎        | 2/16 [00:02<00:19,  1.37s/it]

Width: 0.4, Frequency: 2.0 GHz
Effective epsilons [6.39140243-0.02904808j]
characteristic impedances [68.72273047+0.j]

 19%|█▉        | 3/16 [00:04<00:17,  1.36s/it]

Width: 0.4, Frequency: 3.0 GHz
Effective epsilons [6.41263982-0.0245779j]
characteristic impedances [69.23651174+0.j]

 25%|██▌       | 4/16 [00:05<00:16,  1.37s/it]

Width: 0.4, Frequency: 4.0 GHz
Effective epsilons [6.44126916-0.02206328j]
characteristic impedances [69.68543374+0.j]

 31%|███▏      | 5/16 [00:06<00:14,  1.36s/it]

Width: 0.4, Frequency: 5.0 GHz
Effective epsilons [6.47473968-0.02022274j]
characteristic impedances [70.0959054+0.j]

 38%|███▊      | 6/16 [00:08<00:13,  1.35s/it]

Width: 0.4, Frequency: 6.0 GHz
Effective epsilons [6.51125386-0.01873081j]
characteristic impedances [70.4853003+0.j]

 44%|████▍     | 7/16 [00:09<00:12,  1.35s/it]

Width: 0.4, Frequency: 7.0 GHz
Effective epsilons [6.54964206-0.01747792j]
characteristic impedances [70.86481591+0.j]

 50%|█████     | 8/16 [00:10<00:10,  1.35s/it]

Width: 0.4, Frequency: 8.0 GHz
Effective epsilons [6.58920107-0.01641021j]
characteristic impedances [71.24119986+0.j]

 56%|█████▋    | 9/16 [00:12<00:09,  1.35s/it]

Width: 0.4, Frequency: 9.0 GHz
Effective epsilons [6.62952953-0.0154915j]
characteristic impedances [71.6188911+0.j]

 62%|██████▎   | 10/16 [00:13<00:08,  1.35s/it]

Width: 0.4, Frequency: 10.0 GHz
Effective epsilons [6.67040198-0.01469381j]
characteristic impedances [72.0014128+0.j]

 69%|██████▉   | 11/16 [00:14<00:06,  1.35s/it]

Width: 0.4, Frequency: 11.0 GHz
Effective epsilons [6.71168706-0.01399464j]
characteristic impedances [72.39195995+0.j]

 75%|███████▌  | 12/16 [00:16<00:05,  1.36s/it]

Width: 0.4, Frequency: 12.0 GHz
Effective epsilons [6.75329951-0.01337589j]
characteristic impedances [72.79354843+0.j]

 81%|████████▏ | 13/16 [00:17<00:04,  1.36s/it]

Width: 0.4, Frequency: 13.0 GHz
Effective epsilons [6.79517441-0.01282303j]
characteristic impedances [73.20898252+0.j]

 88%|████████▊ | 14/16 [00:18<00:02,  1.36s/it]

Width: 0.4, Frequency: 14.0 GHz
Effective epsilons [6.83725459-0.01232447j]
characteristic impedances [73.64078+0.j]

 94%|█████████▍| 15/16 [00:20<00:01,  1.35s/it]

Width: 0.4, Frequency: 15.0 GHz
Effective epsilons [6.87948522-0.01187094j]
characteristic impedances [74.0911261+0.j]

100%|██████████| 16/16 [00:21<00:00,  1.35s/it]


                                               

 67%|██████▋   | 2/3 [00:26<00:15, 15.04s/it]
Width: 0.4, Frequency: 16.0 GHz
Effective epsilons [6.92181205-0.01145507j]
characteristic impedances [74.56189761+0.j]

  0%|          | 0/16 [00:00<?, ?it/s]


  6%|▋         | 1/16 [00:12<03:09, 12.61s/it]

Width: 3, Frequency: 1.0 GHz
Effective epsilons [7.72197706-0.02084239j]
characteristic impedances [19.06954516+0.j]

 12%|█▎        | 2/16 [00:24<02:53, 12.36s/it]

Width: 3, Frequency: 2.0 GHz
Effective epsilons [7.77089439-0.01451959j]
characteristic impedances [19.22120091+0.j]

 19%|█▉        | 3/16 [00:36<02:39, 12.24s/it]

Width: 3, Frequency: 3.0 GHz
Effective epsilons [7.84563215-0.01242895j]
characteristic impedances [19.41054896+0.j]

 25%|██▌       | 4/16 [00:49<02:28, 12.37s/it]

Width: 3, Frequency: 4.0 GHz
Effective epsilons [7.93428127-0.01127637j]
characteristic impedances [19.61995951+0.j]

 31%|███▏      | 5/16 [01:01<02:15, 12.32s/it]

Width: 3, Frequency: 5.0 GHz
Effective epsilons [8.02775497-0.01043376j]
characteristic impedances [19.83956641+0.j]

 38%|███▊      | 6/16 [01:13<02:02, 12.24s/it]

Width: 3, Frequency: 6.0 GHz
Effective epsilons [8.12026481-0.00974244j]
characteristic impedances [20.06260255+0.j]

 44%|████▍     | 7/16 [01:25<01:49, 12.17s/it]

Width: 3, Frequency: 7.0 GHz
Effective epsilons [8.20877938-0.00915282j]
characteristic impedances [20.28374168+0.j]

 50%|█████     | 8/16 [01:37<01:37, 12.17s/it]

Width: 3, Frequency: 8.0 GHz
Effective epsilons [8.29206805-0.00864227j]
characteristic impedances [20.49942836+0.j]

 56%|█████▋    | 9/16 [01:50<01:25, 12.20s/it]

Width: 3, Frequency: 9.0 GHz
Effective epsilons [8.36988911-0.00819566j]
characteristic impedances [20.70790646+0.j]

 62%|██████▎   | 10/16 [02:02<01:12, 12.14s/it]

Width: 3, Frequency: 10.0 GHz
Effective epsilons [8.44246521-0.00780105j]
characteristic impedances [20.908791+0.j]

 69%|██████▉   | 11/16 [02:15<01:01, 12.35s/it]

Width: 3, Frequency: 11.0 GHz
Effective epsilons [8.51019396-0.00744873j]
characteristic impedances [21.10250463+0.j]

 75%|███████▌  | 12/16 [02:28<00:50, 12.61s/it]

Width: 3, Frequency: 12.0 GHz
Effective epsilons [8.57350669-0.00713085j]
characteristic impedances [21.28983419+0.j]

 81%|████████▏ | 13/16 [02:41<00:38, 12.82s/it]

Width: 3, Frequency: 13.0 GHz
Effective epsilons [8.63280876-0.00684123j]
characteristic impedances [21.47164679+0.j]

 88%|████████▊ | 14/16 [02:54<00:25, 12.99s/it]

Width: 3, Frequency: 14.0 GHz
Effective epsilons [8.68846011-0.00657503j]
characteristic impedances [21.64872831+0.j]

 94%|█████████▍| 15/16 [03:08<00:13, 13.07s/it]

Width: 3, Frequency: 15.0 GHz
Effective epsilons [8.74077384-0.00632849j]
characteristic impedances [21.82173308+0.j]

100%|██████████| 16/16 [03:21<00:00, 13.09s/it]


                                               

100%|██████████| 3/3 [03:52<00:00, 102.09s/it]
100%|██████████| 3/3 [03:52<00:00, 77.51s/it] 
Width: 3, Frequency: 16.0 GHz
Effective epsilons [8.79002158-0.00609871j]
characteristic impedances [21.99117351+0.j]

Hide code cell source
plt.figure(figsize=(10, 14))
plt.xlabel("Frequency / GHz")
plt.ylabel("Effective dielectric constant")

for i, width in enumerate(widths):
    plt.plot(frequencies / 1e9, epsilon_effs[i, :, 0].real)
    plt.annotate(
        xy=(frequencies[-1] / 1e9, epsilon_effs[i, :, 0].real[-1]),
        text=str(width),
        va="center",
    )

plt.show()

plt.figure(figsize=(10, 14))
plt.xlabel("Frequency / GHz")
plt.ylabel("Characteristic Impedance / Ohm")

for i, width in enumerate(widths):
    plt.plot(frequencies / 1e9, characteristic_impedances[i, :, 0].real)
    plt.annotate(
        xy=(frequencies[-1] / 1e9, characteristic_impedances[i, :, 0].real[-1]),
        text=str(width),
        va="center",
    )

plt.show()
../../_images/d69ed87f19cf01fc88e428bae3121e91f3c991ee970a00b2dbccc1a2a60d50dc.png ../../_images/01ca9400f8a0da17c9a98bd026b60b55a196184d0a831028315100b14df62767.png

Bibliography#

[1]

R.H. Jansen. High-speed computation of single and coupled microstrip parameters including dispersion, high-order modes, loss and finite strip thickness. IEEE Transactions on Microwave Theory and Techniques, 26(2):75–82, February 1978. URL: http://dx.doi.org/10.1109/TMTT.1978.1129316, doi:10.1109/tmtt.1978.1129316.