Device Interface

The tdgl.device subpackage provides the following functionalities:

Overview

Here is a quick overview of the most useful functions and methods for creating and working with polygons and devices. For a demonstration of how all of these pieces work together in practice, see Working with polygons.

Layer

class tdgl.Layer(*, london_lambda, coherence_length, thickness, conductivity=None, u=5.79, gamma=10.0, z0=0)[source]

A superconducting thin film.

Parameters
  • london_lambda (float) – The London penetration depth of the film.

  • coherence_length (float) – The superconducting coherence length of the film.

  • thickness (float) – The thickness of the film.

  • conductivity (Optional[float]) – The normal state conductivity of the superconductor in Siemens / length_unit.

  • u (float) – The ratio of the relaxation times for the order parameter amplitude and phase. This value is 5.79 for dirty superconductors.

  • gamma (float) – This parameter quantifies the effect of inelastic phonon-electron scattering. \(\gamma\) is proportional to the inelastic scattering time and the size of the superconducting gap.

  • z0 (float) – Vertical location of the film.

property Lambda: float

Effective magnetic penetration depth, \(\Lambda=\lambda^2/d\).

copy()[source]

Create a deep copy of the tdgl.Layer.

Return type

Layer

to_hdf5(h5_group)[source]

Save the tdgl.Layer to an h5py.Group.

Parameters

h5_group (Group) – An open h5py.Group to which to save the layer.

Return type

None

static from_hdf5(h5_group)[source]

Load a tdgl.Layer from an h5py.Group.

Parameters

h5_group (Group) – An open h5py.Group from which to load the layer.

Return type

Layer

Returns

A new tdgl.Layer instance.

Polygon

class tdgl.Polygon(name=None, *, points, mesh=True)[source]

A polygonal region located in a tdgl.Layer.

Parameters
  • name (Optional[str]) – Name of the polygon.

  • points (Union[Polygon, ndarray, LineString, LinearRing, Polygon]) – The polygon vertices. This can be a shape (n, 2) array of x, y coordinates or a shapely LineString, LinearRing, or Polygon.

  • mesh (bool) – Whether to include this polygon when computing a mesh.

property points: ndarray

A shape (n, 2) array of counter-clockwise-oriented polygon vertices.

property is_valid: bool

True if the Polygon has a name and its geometry is valid.

property area: float

The area of the polygon.

property bbox: Tuple[Tuple[float, float], Tuple[float, float]]

Returns the coordinates of the lower left and upper right corners of the polygon’s bounding box.

property extents: Tuple[float, float]

Returns the total x, y extent of the polygon, (Delta_x, Delta_y).

property polygon: Polygon

A shapely Polygon representing the Polygon.

property path: Path

A matplotlib.path.Path representing the polygon boundary.

contains_points(points, index=False, radius=0)[source]

Determines whether points lie within the polygon.

Parameters
  • points (ndarray) – Shape (n, 2) array of x, y coordinates.

  • index (bool) – If True, then return the indices of the points in points that lie within the polygon. Otherwise, returns a shape (n, ) boolean array.

  • radius (float) – An additional margin on self.path. See matplotlib.path.Path.contains_points().

Return type

Union[bool, ndarray]

Returns

If index is True, returns the indices of the points in points that lie within the polygon. Otherwise, returns a shape (n, ) boolean array indicating whether each point lies within the polygon.

on_boundary(points, radius=0.001, index=False)[source]

Determines whether points lie within a given radius of the Polygon boundary.

Parameters
  • points (ndarray) – Shape (n, 2) array of x, y coordinates.

  • radius (float) – Points within radius of the boundary are considered to lie on the boundary.

  • index (bool) – If True, then return the indices of the points in points that lie on the boundary. Otherwise, returns a shape (n, ) boolean array.

Returns

If index is True, returns the indices of the points in points that lie within the polygon. Otherwise, returns a shape (n, ) boolean array indicating whether each point lies within the polygon.

make_mesh(min_points=None, smooth=0, **meshpy_kwargs)[source]

Returns the vertices and triangles of a Delaunay mesh covering the Polygon.

Parameters
  • min_points (Optional[int]) – Minimum number of vertices in the mesh. If None, then the number of vertices will be determined by meshpy_kwargs and the number of vertices in the underlying polygons.

  • smooth (int) – Number of Laplacian smoothing steps to perform.

  • **meshpy_kwargs – Passed to meshpy.triangle.build().

Return type

Mesh

Returns

Mesh vertex coordinates and triangle indices

rotate(degrees, origin=(0.0, 0.0), inplace=False)[source]

Rotates the polygon counterclockwise by a given angle.

Parameters
  • degrees (float) – The amount by which to rotate the polygon.

  • origin (Union[str, Tuple[float, float]]) – (x, y) coorindates about which to rotate, or the strings “center” (for the bounding box center) or “centroid” (for the polygon center of mass).

  • inplace (bool) – If True, modify the polygon in place. Otherwise, return a modified copy.

Return type

Polygon

Returns

The rotated polygon.

translate(dx=0.0, dy=0.0, inplace=False)[source]

Translates the polygon by a given distance.

Parameters
  • dx (float) – Distance by which to translate along the x-axis.

  • dy (float) – Distance by which to translate along the y-axis.

  • inplace (bool) – If True, modify the polygon in place. Otherwise, return a modified copy.

Return type

Polygon

Returns

The translated polygon.

scale(xfact=1.0, yfact=1.0, origin=(0, 0), inplace=False)[source]

Scales the polygon horizontally by xfact and vertically by yfact.

Negative xfact (yfact) can be used to reflect the polygon horizontally (vertically) about the origin.

Parameters
  • xfact (float) – Distance by which to translate along the x-axis.

  • yfact (float) – Distance by which to translate along the y-axis.

  • origin (Union[str, Tuple[float, float]]) – (x, y) coorindates for the scaling origin, or the strings “center” (for the bounding box center) or “centroid” (for the polygon center of mass).

  • inplace (bool) – If True, modify the polygon in place. Otherwise, return a modified copy.

Return type

Polygon

Returns

The scaled polygon.

union(*others, name=None)[source]

Returns the union of the polygon and zero or more other polygons.

Parameters
Return type

Polygon

Returns

A new Polygon instance representing the union of self and others.

intersection(*others, name=None)[source]

Returns the intersection of the polygon and zero or more other polygons.

Parameters
Return type

Polygon

Returns

A new Polygon instance representing the intersection of self and others.

difference(*others, name=None)[source]

Returns the difference of the polygon and zero more other polygons.

Parameters
Return type

Polygon

Returns

A new Polygon instance representing the difference of self and others.

buffer(distance, join_style='mitre', mitre_limit=5.0, single_sided=True, as_polygon=True)[source]

Returns polygon points or a new Polygon object with vertices offset from self.points by a given distance. If distance > 0 this “inflates” the polygon, and if distance < 0 this shrinks the polygon.

Parameters
Return type

Union[ndarray, Polygon]

Returns

A new Polygon or an array of vertices offset by distance.

resample(num_points=None, degree=1, smooth=0)[source]

Resample vertices so that they are approximately uniformly distributed along the polygon boundary.

Parameters
  • num_points (Optional[int]) – Number of points to interpolate to. If num_points is None, the polygon is resampled to len(self.points) points. If num_points is not None and has a boolean value of False, then an unaltered copy of the polygon is returned.

  • degree (int) – The degree of the spline with which to interpolate. Defaults to 1 (linear spline).

  • smooth (float) – Smoothing condition.

Return type

Polygon

set_name(name)[source]

Sets the Polygon’s name and returns self.

Return type

Polygon

Parameters

name (Optional[str]) –

plot(ax=None, **kwargs)[source]

Plots the Polygon’s vertices.

Parameters
  • ax (Optional[Axes]) – The matplotlib Axes on which to plot. If None is given, a new one is created.

  • kwargs – Passed to ax.plot().

Return type

Axes

Returns

The matplotlib Axes.

classmethod from_union(items, *, name=None, mesh=True)[source]

Creates a new Polygon from the union of a sequence of polygons.

Parameters
Return type

Polygon

Returns

A new Polygon.

classmethod from_intersection(items, *, name=None, mesh=True)[source]

Creates a new Polygon from the intersection of a sequence of polygons.

Parameters
Return type

Polygon

Returns

A new Polygon.

classmethod from_difference(items, *, name=None, mesh=True)[source]

Creates a new Polygon from the difference of a sequence of polygons.

Parameters
Return type

Polygon

Returns

A new Polygon.

to_hdf5(h5_group)[source]

Save the Polygon to an h5py.Group.

Return type

None

Parameters

h5_group (Group) –

classmethod from_hdf5(h5_group)[source]

Load a Polygon from an h5py.Group.

Return type

Polygon

Parameters

h5_group (Group) –

copy()[source]

Returns a deep copy of the Polygon

Return type

Polygon

Device

class tdgl.Device(name, *, layer, film, holes=None, terminals=None, probe_points=None, length_units='um')[source]

An object representing a thin film superconducting device.

Parameters
  • name (str) – Name of the device.

  • layer (Layer) – The superconducting Layer.

  • film (Polygon) – The Polygon representing the superconducting film.

  • holes (Optional[List[Polygon]]) – Polygons representing holes in the superconducting film.

  • terminals (Optional[List[Polygon]]) – A sequence of Polygons representing the current terminals. Any points that are on the boundary of the mesh and lie inside a terminal will have current source/sink boundary conditions.

  • probe_points (Optional[Sequence[Tuple[float, float]]]) – A shape (n, 2) sequence of floats, with each row representing the (x, y) position of a voltage probe.

  • length_units (str) – Distance units for the coordinate system.

property length_units: str

Length units used for the device geometry.

property coherence_length: Quantity

Ginzburg-Landau coherence length, \(\xi\)

property london_lambda: Quantity

London penetration depth, \(\lambda\)

property thickness: Quantity

Film thickness, \(d\)

property Lambda: Quantity

Effective magnetic penetration depth, \(\Lambda=\lambda^2/d\).

property conductivity: Optional[Quantity]

Film normal state conductivity, \(\sigma\)

property kappa: float

The Ginzburg-Landau parameter, \(\kappa=\lambda/\xi\).

property Bc2: Quantity

Upper critical field, \(B_{c2}=\Phi_0/(2\pi\xi^2)\).

property A0: Quantity

Scale for the magnetic vector potential, \(A_0=\xi B_{c2}\).

property K0: Quantity

Sheet current density scale (dimensions of current / length), \(K_0=4\xi B_{c2}/(\mu_0\Lambda)\).

tau0(conductivity=None)[source]

Time scale, \(\tau_0=\mu_0\sigma\lambda^2\).

Parameters

conductivity (Optional[Quantity]) – The normal state conductivity of the film, which defaults to device.layer.conductivity.

Returns

math:`tau_0=mu_0sigmalambda^2

Return type

The time scale,

V0(conductivity=None)[source]

Electric potential scale, \(\V_0=\xi J_0/\sigma\).

Parameters

conductivity (Optional[Quantity]) – The normal state conductivity of the film, which defaults to device.layer.conductivity.

Return type

Quantity

Returns

The electric potential scale, \(\V_0=\xi J_0/\sigma\)

property triangulation: Optional[Triangulation]

Matplotlib triangulation of the mesh.

terminal_info()[source]

Returns a tuple of TerminalInfo objects, one for each current terminal in the device.

Return type

Tuple[TerminalInfo, ...]

property polygons: Tuple[Polygon, ...]

Tuple of all Polygons in the Device.

property points: Optional[ndarray]

The mesh vertex coordinates in length_units (shape (n, 2), type float).

property triangles: Optional[ndarray]

Mesh triangle indices (shape (m, 3), type int).

property edges: Optional[ndarray]

Mesh edge indices (shape (p, 2), type int).

property edge_lengths: Optional[ndarray]

An array of the mesh vertex-to-vertex distances.

property areas: Optional[ndarray]

An array of the mesh Voronoi cell areas.

property probe_point_indices: Optional[List[int]]

A list of the mesh site indices for the probe points.

boundary_sites()[source]

Returns a dict of {polygon_name: boundary_indices}, where boundary_indices is an integer array of site indices for mesh sites on the boundary of each polygon.

The length of the returned dictionary will be the number of holes in the device plus one.

Return type

Optional[Dict[str, ndarray]]

Returns

{polygon_name: boundary_indices}

contains_points(points, index=False, radius=0)[source]

Determines whether points lie within the Device.

Parameters
  • points (ndarray) – Shape (n, 2) array of x, y coordinates.

  • index (bool) – If True, then return the indices of the points in points that lie within the polygon. Otherwise, returns a shape (n, ) boolean array.

  • radius (float) – An additional margin on self.path. See matplotlib.path.Path.contains_points().

Return type

ndarray

Returns

If index is True, returns the indices of the points in points that lie within the polygon. Otherwise, returns a shape (n, ) boolean array indicating whether each point lies within the polygon.

copy(with_mesh=True)[source]

Copy this Device to create a new one.

Parameters

with_mesh (bool) – Whether to copy the mesh.

Return type

Device

Returns

A new Device instance, copied from self.

scale(xfact=1, yfact=1, origin=(0, 0))[source]

Returns a new device with polygons scaled horizontally and/or vertically.

Negative xfact (yfact) can be used to reflect the device horizontally (vertically) about the origin.

Parameters
  • xfact (float) – Factor by which to scale the device horizontally.

  • yfact (float) – Factor by which to scale the device vertically.

  • origin (Tuple[float, float]) – (x, y) coorindates of the origin.

Return type

Device

Returns

The scaled Device.

rotate(degrees, origin=(0, 0))[source]

Returns a new device with polygons rotated a given amount counterclockwise about specified origin.

Parameters
  • degrees (float) – The amount by which to rotate the polygons.

  • origin (Tuple[float, float]) – (x, y) coorindates of the origin.

Return type

Device

Returns

The rotated Device.

translate(dx=0, dy=0, dz=0, inplace=False)[source]

Translates the device polygons, layer, and mesh in space by a given amount.

Parameters
  • dx (float) – Distance by which to translate along the x-axis.

  • dy (float) – Distance by which to translate along the y-axis.

  • dz (float) – Distance by which to translate layer along the z-axis.

  • inplace (bool) – If True, modifies the device (self) in-place and returns None, otherwise, creates a new device, translates it, and returns it.

Return type

Device

Returns

The translated device.

translation(dx, dy, dz=0)[source]

A context manager that temporarily translates a device in-place, then returns it to its original position.

Parameters
  • dx (float) – Distance by which to translate polygons along the x-axis.

  • dy (float) – Distance by which to translate polygons along the y-axis.

  • dz (float) – Distance by which to translate layers along the z-axis.

Return type

None

make_mesh(max_edge_length=None, min_points=None, smooth=0, **meshpy_kwargs)[source]

Generates and optimizes the triangular mesh.

Parameters
  • min_points (Optional[float]) – Minimum number of vertices in the mesh. If None, then the number of vertices will be determined by meshpy_kwargs, the number of vertices in the underlying polygons, and max_edge_length.

  • max_edge_length (Optional[float]) – The maximum distance between vertices in the mesh. Passing a value <= 0 means that the number of mesh points will be determined solely by the density of points in the Device’s film and holes. Defaults to 1.0 * self.coherence_length.

  • smooth (int) – Number of Laplacian smoothing iterations to perform.

  • **meshpy_kwargs – Passed to meshpy.triangle.build().

Return type

None

mesh_stats_dict()[source]

Returns a dictionary of information about the mesh.

Return type

Dict[str, Union[Real, str]]

mesh_stats(precision=3)[source]

When called with in Jupyter notebook, displays a table of information about the mesh.

Parameters

precision (int) – Number of digits after the decimal for float values.

Return type

HTML

Returns

An HTML table of mesh statistics.

plot(ax=None, legend=True, figsize=None, mesh=False, mesh_kwargs={'color': 'k', 'lw': 0.5}, **kwargs)[source]

Plot all of the device’s polygons.

Parameters
  • ax (Optional[Axes]) – matplotlib axis on which to plot. If None, a new figure is created.

  • subplots – If True, plots each layer on a different subplot.

  • legend (bool) – Whether to add a legend.

  • figsize (Optional[Tuple[float, float]]) – matplotlib figsize, only used if ax is None.

  • mesh (bool) – If True, plot the mesh.

  • mesh_kwargs (Dict[str, Any]) – Keyword arguments passed to ax.triplot() if mesh is True.

  • kwargs – Passed to ax.plot() for the polygon boundaries.

Return type

Tuple[Figure, Axes]

Returns

Matplotlib Figure and Axes

patches()[source]

Returns a dict of {polygon_name: PathPatch} for visualizing the device.

Return type

Dict[str, PathPatch]

draw(ax=None, legend=True, figsize=None, alpha=0.5, exclude=None)[source]

Draws all polygons in the device as matplotlib patches.

Parameters
  • ax (Optional[Axes]) – matplotlib axis on which to plot. If None, a new figure is created.

  • legend (bool) – Whether to add a legend.

  • figsize (Optional[Tuple[float, float]]) – matplotlib figsize, only used if ax is None.

  • alpha (float) – The alpha (opacity) value for the patches (0 <= alpha <= 1).

  • exclude (Union[str, List[str], None]) – A polygon name or list of polygon names to exclude from the figure.

Return type

Tuple[Figure, Union[Axes, ndarray]]

Returns

Matplotlib Figre and Axes.

to_hdf5(path_or_group, save_mesh=True)[source]

Serializes the Device to disk.

Parameters
  • path_or_group (Union[str, File, Group]) – A path to an HDF5 file, or an open HDF5 file or group into which to save the Device.

  • save_mesh (bool) – Whether to serialize the full mesh.

Return type

None

classmethod from_hdf5(path_or_group)[source]

Creates a new Device from one serialized to disk.

Parameters

path_or_group (Union[str, File, Group]) – A path to an HDF5 file, or an open HDF5 file or group containing the serialized Device.

Return type

Device

Returns

The loaded Device instance.

Geometry

The tdgl.geometry module contains functions for creating polygons approximating simple shapes, i.e., rectangles and ellipses, which can be combined using constructive solid geometry to create complex shapes.

tdgl.geometry.box(width, height=None, points=101, center=(0, 0), angle=0)[source]

Returns the coordinates for a rectangle with a given width and height, centered at the specified center.

Parameters
  • width (float) – Width of the rectangle (in the x direction).

  • height (Optional[float]) – Height of the rectangle (in the y direction). If None is given, then height is set to width and the function returns a square.

  • points (int) – The target number of points making up the box. The actual number of points may be slightly different than this value.

  • center (Tuple[float, float]) – Coordinates of the center of the rectangle.

  • angle (float) – Angle (in degrees) by which to rotate counterclockwise about (0, 0) after translating to the specified center.

Return type

ndarray

Returns

A shape (m, 2) or array of (x, y) coordinates.

tdgl.geometry.circle(radius, points=100, center=(0, 0))[source]

Returns the coordinates for a circle with a given radius, centered at the specified center.

Parameters
  • radius (float) – Radius of the circle

  • points (int) – Number of points in the circle

  • center (Tuple[float, float]) – Coordinates of the center of the circle

Return type

ndarray

Returns

A shape (points, 2) array of (x, y) coordinates

tdgl.geometry.ellipse(a, b, points=100, center=(0, 0), angle=0)[source]

Returns the coordinates for an ellipse with major axis a and semimajor axis b, rotated by the specified angle about (0, 0), then translated to the specified center.

Parameters
  • a (float) – Major axis length

  • b (float) – Semi-major axis length

  • points (int) – Number of points in the circle

  • center (Tuple[float, float]) – Coordinates of the center of the circle

  • angle (float) – Angle (in degrees) by which to rotate counterclockwise about (0, 0) after translating to the specified center.

Returns

A shape (points, 2) array of (x, y) coordinates

tdgl.geometry.rotate(coords, angle_degrees)[source]

Rotates an array of (x, y) coordinates counterclockwise by the specified angle.

Parameters
  • coords (ndarray) – Shape (n, 2) array of (x, y) coordinates.

  • angle_degrees (float) – The angle by which to rotate the coordinates.

Return type

ndarray

Returns

Shape (n, 2) array of rotated coordinates (x’, y’)

Mesh Generation

tdgl.generate_mesh(poly_coords, hole_coords=None, min_points=None, max_edge_length=None, convex_hull=False, boundary=None, min_angle=32.5, **kwargs)[source]

Generates a Delaunay mesh for a given set of polygon vertex coordinates.

Additional keyword arguments are passed to triangle.build().

Parameters
  • poly_coords (ndarray) – Shape (n, 2) array of polygon (x, y) coordinates.

  • hole_coords (Optional[List[ndarray]]) – A list of arrays of hole boundary coordinates.

  • min_points (Optional[int]) – The minimimum number of vertices in the resulting mesh.

  • max_edge_length (Optional[float]) – The maximum distance between vertices in the resulting mesh.

  • convex_hull (bool) – If True, then the entire convex hull of the polygon (minus holes) will be meshed. Otherwise, only the polygon interior is meshed.

  • boundary (Optional[ndarray]) – Shape (m, 2) (where m <= n) array of (x, y) coordinates for points on the boundary of the polygon.

  • min_angle (float) – The minimum angle in the mesh’s triangles. Setting a larger value will make the triangles closer to equilateral, but the mesh generation may fail if the value is too large.

Return type

Tuple[ndarray, ndarray]

Returns

Mesh vertex coordinates and triangle indices.