4. Maxwell eigenmode extraction of a coplanar waveguide resonator with adaptive meshing

4.1. Requirements

4.1.1. Software components

  • QTCAD

  • Gmsh

4.1.2. Geometry file

  • qtcad/examples/tutorials/meshes/meandered_resonator.py

4.1.3. Python script

  • qtcad/examples/tutorials/maxwell_meandered_resonator.py

4.1.4. References

4.2. Briefing

In this example we extract the fundamental Maxwell eigenmode of a quarter-wavelength meandered coplanar waveguide resonator. The resonator is formed by the ground plane and a central strip, placed on a silicon substrate. The structure is enclosed inside a cavity.

4.3. Geometry of the problem

The geometry of the device is shown in the figure below:

Geometry of the meandered resonator device

The following table shows the physical groups setup in the geometry file:

Gmsh physical group

Dimension

Purpose

"gnd"

2D

The ground plane (modeled as a PEC [1] sheet)

"strip"

2D

The central strip of the coplanar waveguide (modeled as a PEC sheet)

"substrate"

3D

The silicon substrate on which the ground plane and strip lie

"envelope"

2D

The inner surface of the cavity enclosing the structure [2]

"air"

3D

The air that fills the space in the cavity above the silicon substrate

Layout of the meandered resonator

4.4. Setting up the device and finding the Maxwell eigenmodes

The procedure for finding the Maxwell eigenmodes is similar to that for the capacitance solver (see, for example, Capacitance extraction of a coaxial cable with adaptive meshing). You will also find it useful to look at Maxwell Eigenmodes and the API reference for the maxwell_eigenmode module.

4.4.1. Header, input parameters, and file paths

"""
Find the fundamental mode of a meandered quarter-wavelength coplanar waveguide
resonator inside a cavity.
"""

from pathlib import Path
import os
from time import time
import numpy as np

# import relevant modules of QTCAD
from qtcad.device.maxwell_eigenmode import Solver
from qtcad.device.maxwell_eigenmode import SolverParams
from qtcad.device.device import Device
from qtcad.device.mesh3d import Mesh
from qtcad.device import materials as mt
from qtcad.device import constants as ct

Similarly to capacitance, the maxwell_eigenmode module contains Solver and SolverParams classes. The Solver class contains algorithms needed to find the Maxwell eigenmodes and the SolverParams class can be used to control parameters of the solver.

# scale in the Gmsh files
scale = 1e-6

# the total length of the resonator in meters
length = 5950 * scale

# directories and file paths
script_dir = Path(__file__).parent.resolve()
input_dir = script_dir / "meshes"  # for mesh and raw geometry files
result_dir = script_dir / "output" / Path(__file__).stem  # for results
fpath_mesh = input_dir / "meandered_resonator.msh4" # mesh file
fpath_xao = input_dir / "meandered_resonator.xao" # raw geometry file

# code to check that the mesh and raw geometry files exist
if not os.path.isfile(fpath_mesh) or not os.path.isfile(fpath_xao):
    raise Exception(
        "Please run %s/meandered_resonator.py to generate the mesh and raw geometry files."
        % (input_dir)
    )

4.4.2. Loading the initial mesh and defining the device

As in the case of the capacitance solver, we first parse the mesh file, initialize the device,

######################################################################################
# setup the device
######################################################################################

# parse the mesh and initialize the device
mesh = Mesh(scale, fpath_mesh)
dvc = Device(mesh)

and assign the media to 3D physical groups and boundary conditions to 2D physical groups.

material_sub = mt.Si
material_air = mt.vacuum
# assign media to regions
dvc.new_region("substrate", material_sub)
dvc.new_region("air", material_air)
# assign perfect electric conductor boundary condition to all conductors
dvc.new_pec_bnd("gnd")
dvc.new_pec_bnd("strip")
dvc.new_pec_bnd("envelope")

In the third and forth line, we assign silicon to "substrate" and vacuum to "air". Next, we assign the perfect electric conductor (PEC) boundary condition to the conducting surfaces "gnd" and "strip" and to the inner surface of the cavity "envelope". See Boundary conditions for more details on the boundary conditions in the Maxwell eigensolver.

4.4.3. Creating and running the Maxwell eigenmode solver

Next, we define the parameters of the solver by creating a SolverParams object and modifying its attributes.

######################################################################################
# setup the solver
######################################################################################
params = SolverParams()
params.num_modes = 1 # number of modes to find
params.output_dir = result_dir # directory where output files will be stored
params.tol_rel = 0.05 # adaptive meshing tolerance on the frequency (relative)

# uncomment the line below in order to save all intermediate results at every iteration
# params.save_intermediate_results = True

# Number of consecutive iterations that must agree within the tolerance
# thresholds. Uncomment the line below to change the default behavior. Then, for each mode,
# the last result will be compared to 7 previous data points (i.e. 8 data points in total
# are compared for each mode).
# params.min_converged_iters = 7

# Do not terminate until min_iters refinements have been completed
params.min_iters = 35

# Uncomment the line below in order to impose an additional convergence criterion
# on the number of nodes
# params.min_nodes = 10000

In the code above, num_modes is the number of modes to be found (in this case only the fundamental mode). Parameter output_dir specifies the directory where the solver results will be saved, similarly to the capacitance solver. Parameter tol_rel is the relative tolerance on the frequencies.

Other useful parameters include save_intermediate_results and min_converged_iters. The former allows saving all intermediate results to the output_dir directory and is useful if the execution of the solver needs to be interrupted. Parameter min_converged_iters indicates how many of the last iterations are used to compute the current error in the eigenmode frequencies. It can be increased for improved robustness, but this would lead to longer simulations.

Parameter min_iters can be used to prevent the solver from terminating before the given number of refinements has been completed. It is useful in the cases where the error appears to have converged before the mesh is sufficiently refined. To avoid this issue, it is a good practice to study the final refined mesh and ensure that it is resolves all of important features in the device geometry. If the simulation suffers from the aforementioned false convevergence issues, min_iters can be used to force the solver not to terminate too early. See Studying the refined mesh for the comparison of the refined meshes with and without using this parameter.

Another parameter that can be used to control the convergence is min_nodes, which (if specified) prevents the algorithm from terminating before the refined mesh contains at least this amount of nodes.

Next, we initialize and execute the solver with the parameters params created above.

######################################################################################
# Run the solver (results are stored in the folder specified by params.output_dir)
######################################################################################
slv = Solver(dvc, params, geo_file=fpath_xao)
t0 = time()
slv.solve()
dt = time() - t0
print("Solution completed in %.2f s" % dt)

Similarly to the capacitance solver, the input geo_file is needed in order to use adaptive meshing solver instead of the static meshing. The mesh created by qtcad/examples/tutorials/meshes/meandered_resonator.py will then be used as the initial mesh.

The solver will display the following lines:

Convergence data will be saved in:
     [params.output_dir]\meandered_resonator-values.txt
     [params.output_dir]\meandered_resonator-delta.txt
     [params.output_dir]\meandered_resonator-delta_rel.txt

indicating the paths to files where one will find

  • mode frequencies at each iteration,

  • the corresponding mesh sizes,

  • timing information for each iteration and cumulative time since the start of the solver execution.

4.4.4. Studying the refined mesh

The refined mesh can be viewed via the following command:

gmsh [params.output_dir]\meandered_resonator-refined-mesh.msh4

The mesh is shown below for two simulations: with and without the line params.min_iters = 35 [3].

Refined mesh with and without the line ``params.min_iters = 35``

Setting the parameter min_iters to 35, we obtain a mesh with better refinement of the gap between the strip and the ground of the coplanar waveguide. We remark that the desired refinement strongly depends on the accuracy required for the particular application.

4.4.5. Extracting further solver results

In addition to the information in the convergence files, one can access the solver results via the Device object dvc. The mode frequencies are contained in the attribute maxwell_freqs, used in the following code to display the final computed frequency:

######################################################################################
# Extract additional information from the final results
######################################################################################

# Find the theoretical value

# NOTE:
#   The effective permittivity is computed for a coplanar waveguide placed on an
#   infinitely thick and wide substrate.
#
#       Simons, Rainee N. Coplanar waveguide circuits, components, and systems.
#       John Wiley & Sons, 2001.
#
#   Effects of the finite substrate size, or the cavity are neglected due to the
#   large distance between the resonator and the cavity (and consequently also the
#   edges of the substrate). Please refer to Simons (2001) for details on the
#   treatment of such effects. It should also be noted that the theoretical value
#   does not take into account the meandered layout of the resonator.

eps_eff = (material_air.eps + material_sub.eps) / 2
wavelength = length * 4
vphase = 1 / np.sqrt(eps_eff * ct.mu0)
freq_theor = vphase / wavelength
print("Frequency of the fundamental mode (analytical): %.3f GHz" % (freq_theor / 1e9))

# display the frequency of the fundamental mode and compare to theoretical value
err = abs(dvc.maxwell_freqs[0] - freq_theor) / freq_theor
print(
    "Frequency of the fundamental mode (computed): %.3f GHz (rel. err: %.3f)"
    % (dvc.maxwell_freqs[0] / 1e9, err)
)

with the results shown below:

Frequency of the fundamental mode (analytical): 4.979 GHz
Frequency of the fundamental mode (computed): 4.758 GHz (rel. err: 0.044)

Methods get_e_field_maxwell and get_b_field_maxwell can be used to find the electric field and magnetic flux density for each of the computed modes.

# find the electric and magnetic fields
e_field = dvc.get_e_field_maxwell()
b_field = dvc.get_b_field_maxwell()

The code below shows an example of how the resulting fields can be used to extract information about the energy stored in different parts of the device.

# fields in the fundamental mode
e_field_0 = e_field[..., 0]
b_field_0 = b_field[..., 0]

# Find the ratio between the energy stored in the substrate and the energy in the air
energy_sub = dvc.energy_e(e_field_0, "substrate") + dvc.energy_b(b_field_0, "substrate")
energy_air = dvc.energy_e(e_field_0, "air") + dvc.energy_b(b_field_0, "air")

print(
    "In the fundamental mode, the ratio of energy in the substrate and in the air is %.1f:1"
    % (energy_sub / energy_air)
)

The code produces the following result:

In the fundamental mode, the ratio of energy in the substrate and in the air is 2.5:1

In the params.output_dir directory, the solver also saves a .pos file, which can be opened in Gmsh, containing the magnitude of the electric field for each mode. To visualize the field in the fundamental mode, one can run the following command:

gmsh [params.output_dir]\meandered_resonator-e-magn-mode-0.pos

which opens the Gmsh GUI. One can then use the Clipping tool in Gmsh to view the fields at the height of the ground plane.

Clipping tool in Gmsh GUI

The resulting fields are shown in the figure below:

Magnitude of the electric field in the fundamental mode

4.5. Full code

__copyright__ = "Copyright 2022-2025, Nanoacademic Technologies Inc."

"""
Find the fundamental mode of a meandered quarter-wavelength coplanar waveguide 
resonator inside a cavity.
"""
from pathlib import Path
import os
from time import time
import numpy as np
# import relevant modules of QTCAD
from qtcad.device.maxwell_eigenmode import Solver
from qtcad.device.maxwell_eigenmode import SolverParams
from qtcad.device.device import Device
from qtcad.device.mesh3d import Mesh
from qtcad.device import materials as mt
from qtcad.device import constants as ct

# scale in the Gmsh files
scale = 1e-6

# the total length of the resonator in meters
length = 5950 * scale

# directories and file paths
script_dir = Path(__file__).parent.resolve()
input_dir = script_dir / "meshes"  # for mesh and raw geometry files
result_dir = script_dir / "output" / Path(__file__).stem  # for results
fpath_mesh = input_dir / "meandered_resonator.msh4" # mesh file
fpath_xao = input_dir / "meandered_resonator.xao" # raw geometry file

# code to check that the mesh and raw geometry files exist
if not os.path.isfile(fpath_mesh) or not os.path.isfile(fpath_xao):
    raise Exception(
        "Please run %s/meandered_resonator.py to generate the mesh and raw geometry files."
        % (input_dir)
    )


######################################################################################
# setup the device
######################################################################################

# parse the mesh and initialize the device
mesh = Mesh(scale, fpath_mesh)
dvc = Device(mesh)

material_sub = mt.Si
material_air = mt.vacuum
# assign media to regions
dvc.new_region("substrate", material_sub)
dvc.new_region("air", material_air)
# assign perfect electric conductor boundary condition to all conductors
dvc.new_pec_bnd("gnd")
dvc.new_pec_bnd("strip")
dvc.new_pec_bnd("envelope")

######################################################################################
# setup the solver
######################################################################################
params = SolverParams()
params.num_modes = 1 # number of modes to find
params.output_dir = result_dir # directory where output files will be stored
params.tol_rel = 0.05 # adaptive meshing tolerance on the frequency (relative)

# uncomment the line below in order to save all intermediate results at every iteration
# params.save_intermediate_results = True 

# Number of consecutive iterations that must agree within the tolerance
# thresholds. Uncomment the line below to change the default behavior. Then, for each mode, 
# the last result will be compared to 7 previous data points (i.e. 8 data points in total 
# are compared for each mode).
# params.min_converged_iters = 7 

# Do not terminate until min_iters refinements have been completed
params.min_iters = 35

# Uncomment the line below in order to impose an additional convergence criterion 
# on the number of nodes
# params.min_nodes = 10000

######################################################################################
# Run the solver (results are stored in the folder specified by params.output_dir)
######################################################################################
slv = Solver(dvc, params, geo_file=fpath_xao)
t0 = time()
slv.solve()
dt = time() - t0
print("Solution completed in %.2f s" % dt)

######################################################################################
# Extract additional information from the final results
######################################################################################

# Find the theoretical value

# NOTE:
#   The effective permittivity is computed for a coplanar waveguide placed on an 
#   infinitely thick and wide substrate. 
#
#       Simons, Rainee N. Coplanar waveguide circuits, components, and systems.
#       John Wiley & Sons, 2001.
#
#   Effects of the finite substrate size, or the cavity are neglected due to the 
#   large distance between the resonator and the cavity (and consequently also the 
#   edges of the substrate). Please refer to Simons (2001) for details on the 
#   treatment of such effects. It should also be noted that the theoretical value 
#   does not take into account the meandered layout of the resonator.

eps_eff = (material_air.eps + material_sub.eps) / 2
wavelength = length * 4
vphase = 1 / np.sqrt(eps_eff * ct.mu0)
freq_theor = vphase / wavelength
print("Frequency of the fundamental mode (analytical): %.3f GHz" % (freq_theor / 1e9))

# display the frequency of the fundamental mode and compare to theoretical value
err = abs(dvc.maxwell_freqs[0] - freq_theor) / freq_theor
print(
    "Frequency of the fundamental mode (computed): %.3f GHz (rel. err: %.3f)"
    % (dvc.maxwell_freqs[0] / 1e9, err)
)

# find the electric and magnetic fields
e_field = dvc.get_e_field_maxwell()
b_field = dvc.get_b_field_maxwell()

# fields in the fundamental mode
e_field_0 = e_field[..., 0]
b_field_0 = b_field[..., 0]

# Find the ratio between the energy stored in the substrate and the energy in the air
energy_sub = dvc.energy_e(e_field_0, "substrate") + dvc.energy_b(b_field_0, "substrate")
energy_air = dvc.energy_e(e_field_0, "air") + dvc.energy_b(b_field_0, "air")

print(
    "In the fundamental mode, the ratio of energy in the substrate and in the air is %.1f:1"
    % (energy_sub / energy_air)
)