Devices & Geometry¶
The classes defined in superscreen.device
, superscreen.geometry
,
and superscreen.parameters
are used to set up the inputs to a
SuperScreen
simulation, namely:
The geometry and penetration depth of all superconducting films
The spatial distribution of the applied magnetic field
Devices¶
Device¶
- class superscreen.Device(name, *, layers, films, holes=None, abstract_regions=None, length_units='um', solve_dtype='float64')[source]¶
An object representing a device composed of multiple layers of thin film superconductor.
- Parameters:
name (
str
) – Name of the device.layers (
Union
[List
[Layer
],Dict
[str
,Layer
]]) –Layers
making up the device.films (
Union
[List
[Polygon
],Dict
[str
,Polygon
]]) –Polygons
representing regions of superconductor.holes (
Union
[List
[Polygon
],Dict
[str
,Polygon
],None
]) –Holes
representing holes in superconducting films.abstract_regions (
Union
[List
[Polygon
],Dict
[str
,Polygon
],None
]) –Polygons
representing abstract regions in a device. Abstract regions will be meshed, and one can calculate the flux through them.length_units (
str
) – Distance units for the coordinate system.solve_dtype (
Union
[str
,dtype
]) – The float data type to use when solving the device.
- polygons_by_layer(polygon_type=None)[source]¶
Returns a dict of
{layer_name: list of polygons in layer}
.
- contains_points_by_layer(points, polygon_type=None, index=False, radius=0)[source]¶
For each layer, determines whether
points
lie within a polygon in that layer.- Parameters:
points (
ndarray
) – Shape(n, 2)
array of x, y coordinates.polygon_type (
Optional
[str
]) – One of ‘film’, ‘hole’, or ‘all’, specifying which types of polygons to include.index (
bool
) – If True, then return the indices of the points inpoints
that lie within the polygon. Otherwise, returns a shape(n, )
boolean array.radius (
float
) – An additional margin onpolygon.path
. Seematplotlib.path.Path.contains_points()
.
- Return type:
- Returns:
A dict of
{layer_name: contains_points}
. If index is True, thencontains_points
is an array of indices of the points inpoints
that lie within any polygon in the layer. Otherwise,contains_points
is a shape(n, )
boolean array indicating whether each point lies within a polygon in the layer.
- property poly_points: ndarray¶
Shape (n, 2) array of (x, y) coordinates of all Polygons in the Device.
- Return type:
- get_arrays(copy_arrays=False, dense=True)[source]¶
Returns a dict of the large arrays that belong to the device.
- set_arrays(arrays)[source]¶
Sets the Device’s large arrays from a dict like that returned by Device.get_arrays().
- 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 theorigin
.
- rotate(degrees, origin=(0, 0))[source]¶
Returns a new device with polygons rotated a given amount counterclockwise about specified origin.
- mirror_layers(about_z=0.0)[source]¶
Returns a new device with its layers mirrored about the plane
z = about_z
.- Parameters:
about_z (
float
) – The z-position of the plane (parallel to the x-y plane) about which to mirror the layers.- Return type:
- Returns:
The mirrored
superscreen.Device
.
- translate(dx=0, dy=0, dz=0, inplace=False)[source]¶
Translates the device polygons, layers, 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 layers 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:
- 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.
- make_mesh(bounding_polygon=None, compute_matrices=True, convex_hull=True, weight_method='half_cotangent', min_points=None, smooth=0, **meshpy_kwargs)[source]¶
Generates and optimizes the triangular mesh.
- Parameters:
compute_matrices (
bool
) – Whether to compute the field-independent matrices (weights, Q, Laplace operator) needed for Brandt simulations.convex_hull (
bool
) – If True, mesh the entire convex hull of the device’s polygons.weight_method (
str
) – Weight methods for computing the Laplace operator: one of “uniform”, “half_cotangent”, or “inv_euclidian”.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 iterations to perform.**meshpy_kwargs – Passed to meshpy.triangle.build().
bounding_polygon (str | None) –
- Return type:
- compute_matrices(weight_method='half_cotangent')[source]¶
Calculcates mesh weights, Laplace oeprator, and kernel functions.
- mutual_inductance_matrix(hole_polygon_mapping=None, units='pH', all_iterations=False, **solve_kwargs)[source]¶
Calculates the mutual inductance matrix \(\mathbf{M}\) for the Device.
\(\mathbf{M}\) is defined such that element \(M_{ij}\) is the mutual inductance \(\Phi^f_{S_i} / I_j\) between hole \(j\) and polygon \(S_{i}\), which encloses hole \(i\). \(\Phi^f_{S_i}\) is the fluxoid for polygon \(S_i\) and \(I_j\) is the current circulating around hole \(j\).
- Parameters:
hole_polygon_mapping (
Optional
[Dict
[str
,ndarray
]]) – A dict of{hole_name: polygon_coordinates}
specifying a mapping between holes in the device and polygons enclosing those holes, for which the fluxoid will be calculated. The length of this dict,n_holes
, will determine the dimension of the square mutual inductance matrix \(M\).units (
str
) – The units in which to report the mutual inductance.all_iterations (
bool
) – Whether to return mutual inductance matrices for alliterations + 1
solutions, or just the final solution.solve_kwargs – Keyword arguments passed to
superscreen.solve.solve()
, e.g.iterations
.
- Return type:
- Returns:
If all_iterations is False, returns a shape
(n_holes, n_holes)
mutual inductance matrix from the final iteration. Otherwise, returns a list of mutual inductance matrices, each with shape(n_holes, n_holes)
. The length of the list is1
if the device has a single layer, oriterations + 1
if the device has multiple layers.
- plot(ax=None, subplots=False, 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 (
bool
) – 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 toax.triplot()
ifmesh
is True.kwargs – Passed to
ax.plot()
for the polygon boundaries.
- Return type:
- Returns:
Matplotlib Figure and Axes
- patches()[source]¶
Returns a dict of
{layer_name: {film_name: PathPatch}}
for visualizing the device. `
- draw(ax=None, subplots=False, max_cols=3, legend=False, figsize=None, alpha=0.5, exclude=None, layer_order='increasing')[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.subplots (
bool
) – If True, plots each layer on a different subplot.max_cols (
int
) – The maximum number of columns to create ifsubplots
is True.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.layer_order (
str
) – If"increasing"
("decreasing"
) draw polygons in increasing (decreasing) order by layer heightlayer.z0
.
- Return type:
- Returns:
Matplotlib Figre and Axes.
- to_file(directory, save_mesh=True, compressed=True)[source]¶
Serializes the Device to disk.
- Parameters:
- Return type:
TransportDevice¶
- class superscreen.TransportDevice(name, *, layer, film, source_terminals=None, drain_terminal=None, holes=None, abstract_regions=None, length_units='um', solve_dtype='float64')[source]¶
Bases:
Device
An object representing a single-layer, single-film device to which a source-drain current bias can be applied.
There can be multiple source terminals, but only a single drain terminal. The drain terminal acts as a sink for currents from all source terminals.
- Parameters:
name (
str
) – Name of the device.layer (
Layer
) – TheLayer
making up the device.film (
Polygon
) – ThePolygon
representing superconducting film.source_terminals (
Optional
[List
[Polygon
]]) – A list of Polygons representing source terminals for a transport current. Any mesh vertices that are on the boundary of the film and lie inside a source terminal will have current source boundary conditions.drain_terminal (
Optional
[Polygon
]) – Polygon representing the sole current drain (or output) terminal of the device.holes (
Union
[List
[Polygon
],Dict
[str
,Polygon
],None
]) –Holes
representing holes in superconducting films.abstract_regions (
Union
[List
[Polygon
],Dict
[str
,Polygon
],None
]) –Polygons
representing abstract regions in a device. Abstract regions will be meshed, and one can calculate the flux through them.length_units (
str
) – Distance units for the coordinate system.solve_dtype (
Union
[str
,dtype
]) – The float data type to use when solving the device.
- copy(with_arrays=True, copy_arrays=False)[source]¶
Copy this Device to create a new one.
- Parameters:
- Return type:
- Returns:
A new Device instance, copied from self
- make_mesh(compute_matrices=True, weight_method='half_cotangent', min_points=None, smooth=0, **meshpy_kwargs)[source]¶
Generates and optimizes the triangular mesh.
- Parameters:
compute_matrices (
bool
) – Whether to compute the field-independent matrices (weights, Q, Laplace operator) needed for Brandt simulations.weight_method (
str
) – Weight methods for computing the Laplace operator: one of “uniform”, “half_cotangent”, or “inv_euclidian”.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:
Layer¶
- class superscreen.Layer(name, Lambda=None, london_lambda=None, thickness=None, z0=0)[source]¶
A single layer of a superconducting device.
You can provide either an effective penetration depth Lambda, or both a London penetration depth (lambda_london) and a layer thickness. Lambda and london_lambda can either be real numers or Parameters which compute the penetration depth as a function of position.
- Parameters:
name (
str
) – Name of the layer.Lambda (
Union
[float
,Parameter
,None
]) – The effective magnetic penetration depth of the superconducting film(s) in the layer.thickness (
Optional
[float
]) – Thickness of the superconducting film(s) located in the layer.london_lambda (
Union
[float
,Parameter
,None
]) – London penetration depth of the superconducting film(s) located in the layer.z0 (
float
) – Vertical location of the layer.
Polygon¶
- class superscreen.Polygon(name=None, *, layer=None, points, mesh=True)[source]¶
A polygonal region located in a Layer.
- Parameters:
layer (
Optional
[str
]) – Name of the layer in which the polygon is located.points (
Union
[ndarray
,LineString
,LinearRing
,Polygon
]) – The polygon vertices. This can be a shape(n, 2)
array of x, y coordinates or a shapelyLineString
,LinearRing
, orPolygon
.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.- Return type:
- property is_valid: bool¶
True if the
Polygon
has aname
andlayer
its geometry is valid.- Return type:
- property extents: Tuple[float, float]¶
Returns the total x, y extent of the polygon, (Delta_x, Delta_y).
- 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 inpoints
that lie within the polygon. Otherwise, returns a shape(n, )
boolean array.radius (
float
) – An additional margin onself.path
. Seematplotlib.path.Path.contains_points()
.
- Return type:
- 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:
- 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:
- 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:
- Returns:
The rotated polygon.
- scale(xfact=1.0, yfact=1.0, origin=(0, 0), inplace=False)[source]¶
Scales the polygon horizontally by
xfact
and vertically byyfact
.Negative
xfact
(yfact
) can be used to reflect the polygon horizontally (vertically) about theorigin
.- 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:
- Returns:
The scaled polygon.
- union(*others, name=None)[source]¶
Returns the union of the polygon and zero or more other polygons.
- Parameters:
others (
Union
[Polygon
,ndarray
,LineString
,LinearRing
,Polygon
]) – One or more objects with which to join the polygon.name (
Optional
[str
]) – A name for the resulting joined Polygon (defaults toself.name
.)
- Return type:
- Returns:
A new
Polygon
instance representing the union ofself
andothers
.
- intersection(*others, name=None)[source]¶
Returns the intersection of the polygon and zero or more other polygons.
- Parameters:
others (
Union
[Polygon
,ndarray
,LineString
,LinearRing
,Polygon
]) – One or more objects with which to join the polygon.name (
Optional
[str
]) – A name for the resulting joined Polygon (defaults toself.name
.)
- Return type:
- Returns:
A new
Polygon
instance representing the intersection ofself
andothers
.
- difference(*others, symmetric=False, name=None)[source]¶
Returns the difference of the polygon and zero more other polygons.
- Parameters:
others (
Union
[Polygon
,ndarray
,LineString
,LinearRing
,Polygon
]) – One or more objects with which to join the polygon.symmetric (
bool
) – Whether to join via a symmetric difference operation (aka XOR). See the shapely documentation.name (
Optional
[str
]) – A name for the resulting joined Polygon (defaults toself.name
.)
- Return type:
- Returns:
A new
Polygon
instance representing the difference ofself
andothers
.
- 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 givendistance
. Ifdistance > 0
this “inflates” the polygon, and ifdistance < 0
this shrinks the polygon.- Parameters:
distance (
float
) – The amount by which to inflate or deflate the polygon.join_style (
Union
[str
,int
]) – One of “round” (1), “mitre” (2), or “bevel” (3). See the shapely documentation.mitre_limit (
float
) – See the shapely documentation.single_sided (
bool
) – See the shapely documentation.as_polygon (
bool
) – If True, returns a newPolygon
instance, otherwise returns a shape(n, 2)
array of polygon vertices.
- Return type:
- Returns:
A new
Polygon
or an array of vertices offset bydistance
.
- 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. Ifnum_points
is None, the polygon is resampled tolen(self.points)
points. Ifnum_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 iterpolate. Defaults to 1 (linear spline).smooth (
float
) – Smoothing condition.
- Return type:
- classmethod from_union(items, *, name=None, layer=None, mesh=True)[source]¶
Creates a new
Polygon
from the union of a sequence of polygons.- Parameters:
items (
Iterable
[Union
[Polygon
,ndarray
,LineString
,LinearRing
,Polygon
]]) – A sequence of polygon-like objects to join.layer (
Optional
[str
]) – Name of the layer in which the polygon is located.mesh (
bool
) – Whether to include this polygon when computing a mesh.
- Return type:
- Returns:
A new
Polygon
.
- classmethod from_intersection(items, *, name=None, layer=None, mesh=True)[source]¶
Creates a new
Polygon
from the intersection of a sequence of polygons.- Parameters:
items (
Iterable
[Union
[Polygon
,ndarray
,LineString
,LinearRing
,Polygon
]]) – A sequence of polygon-like objects to join.layer (
Optional
[str
]) – Name of the layer in which the polygon is located.mesh (
bool
) – Whether to include this polygon when computing a mesh.
- Return type:
- Returns:
A new
Polygon
.
- classmethod from_difference(items, *, name=None, layer=None, mesh=True, symmetric=False)[source]¶
Creates a new
Polygon
from the difference of a sequence of polygons.- Parameters:
items (
Iterable
[Union
[Polygon
,ndarray
,LineString
,LinearRing
,Polygon
]]) – A sequence of polygon-like objects to join.layer (
Optional
[str
]) – Name of the layer in which the polygon is located.mesh (
bool
) – Whether to include this polygon when computing a mesh.symmetric (
bool
) – If True, creates a newPolygon
from the “symmetric difference” (aka XOR) of the inputs.
- Return type:
- Returns:
A new
Polygon
.
Meshing¶
- superscreen.device.generate_mesh(coords, min_points=None, convex_hull=False, boundary=None, **kwargs)[source]¶
Generates a Delaunay mesh for a given set of polygon vertex coordinates.
Additional keyword arguments are passed to
triangle.build()
.- Parameters:
coords (
ndarray
) – Shape(n, 2)
array of polygon(x, y)
coordinates.min_points (
Optional
[int
]) – The minimimum number of vertices in the resulting mesh.convex_hull (
bool
) – If True, then the entire convex hull of the polygon will be meshed. Otherwise, only the polygon interior is meshed.boundary (ndarray | None) –
- Return type:
- Returns:
Mesh vertex coordinates and triangle indices.
Parameters¶
Parameter¶
- class superscreen.Parameter(func, **kwargs)[source]¶
A callable object that computes a scalar or vector quantity as a function of position coordinates x, y (and optionally z).
Addition, subtraction, multiplication, and division between multiple Parameters and/or real numbers (ints and floats) is supported. The result of any of these operations is a
CompositeParameter
object.- Parameters:
func (
Callable
) – A callable/function that actually calculates the parameter’s value. The function must take x, y (and optionally z) as the first and only positional arguments, and all other arguments must be keyword arguments. Therefore func should have a signature likefunc(x, y, z, a=1, b=2, c=True)
,func(x, y, *, a, b, c)
,func(x, y, z, *, a, b, c)
, orfunc(x, y, z, *, a, b=None, c=3)
.kwargs – Keyword arguments for func.
CompositeParameter¶
- class superscreen.parameter.CompositeParameter(left, right, operator_)[source]¶
Bases:
Parameter
A callable object that behaves like a Parameter (i.e. it computes a scalar or vector quantity as a function of position coordinates x, y, z). A CompositeParameter object is created as a result of mathematical operations between Parameters, CompositeParameters, and/or real numbers.
Addition, subtraction, multiplication, division, and exponentiation between Parameters, CompositeParameters and real numbers (ints and floats) are supported. The result of any of these operations is a new CompositeParameter object.
- Parameters:
left (
Union
[int
,float
,Parameter
,CompositeParameter
]) – The object on the left-hand side of the operator.right (
Union
[int
,float
,Parameter
,CompositeParameter
]) – The object on the right-hand side of the operator.operator – The operator acting on left and right (or its string representation).
Constant¶
Geometry¶
- superscreen.geometry.rotate(coords, angle_degrees)[source]¶
Rotates an array of (x, y) coordinates counterclockwise by the specified angle.
- superscreen.geometry.translate(coords, dx, dy)[source]¶
Translates the given (x, y) coordinates in the palne by
(dx, dy)
.
- superscreen.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 lengthb (
float
) – Semi-major axis lengthpoints (
int
) – Number of points in the circlecenter (
Tuple
[float
,float
]) – Coordinates of the center of the circleangle (
float
) – Angle (in degrees) by which to rotate counterclockwise about (0, 0) before translating to the specified center.
- Returns:
A shape
(points, 2)
array of (x, y) coordinates
- superscreen.geometry.circle(radius, points=100, center=(0, 0))[source]¶
Returns the coordinates for a circle with a given radius, centered at the specified center.
- superscreen.geometry.box(width, height=None, points_per_side=25, 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_per_side (
int
) – Number of points on each side of the box.center (
Tuple
[float
,float
]) – Coordinates of the center of the rectangle.angle (
float
) – Angle (in degrees) by which to rotate counterclockwise about (0, 0) before translating to the specified center.
- Return type:
- Returns:
A shape
(4 * points_per_side, 2)
array of (x, y) coordinates