Solver

The superscreen.solve module contains the actual implementation of Brandt’s method, as described here.

Solve

superscreen.solve(device=None, *, model=None, applied_field=None, terminal_currents=None, circulating_currents=None, vortices=None, field_units='mT', current_units='uA', check_inversion=False, iterations=0, return_solutions=True, save_path=None, log_level=None, progress_bar=True, _solver='superscreen.solve')[source]

Computes the stream functions and magnetic fields for all layers in a Device.

The simulation strategy is:

1. Compute the stream functions and fields for each film given only the applied field.

2. If iterations > 1 and there are multiple films, then for each film, calculate the screening field from all other films and recompute the stream function and fields based on the sum of the applied field and the screening fields from all other films.

  1. Repeat step 2 (iterations - 1) times.

Parameters:
  • device (Optional[Device]) – The Device to simulate. Required if model is not provided.

  • model (Optional[FactorizedModel]) – A superscreen.FactorizedModel instance.

  • applied_field (Optional[Callable]) – A callable that computes the applied magnetic field as a function of x, y, z coordinates.

  • terminal_currents (Optional[Dict[str, Dict[str, Union[float, str, Quantity]]]]) – A dict like {film_name: {source_name: source_current}} for each film and terminal.

  • circulating_currents (Optional[Dict[str, Union[float, str, Quantity]]]) – A dict of {hole_name: circulating_current}. If circulating_current is a float, then it is assumed to be in units of current_units. If circulating_current is a string, then it is converted to a pint.Quantity.

  • vortices (Optional[Sequence[Vortex]]) – A list of Vortex objects located in the Device.

  • field_units (str) – Units of the applied field. Can either be magnetic field H or magnetic flux density B = mu0 * H.

  • current_units (str) – Units to use for current quantities. The applied field will be converted to units of [current_units / device.length_units].

  • check_inversion (bool) – Whether to verify the accuracy of the matrix inversion.

  • iterations (int) – Number of times to compute the interactions between layers.

  • return_solutions (bool) – Whether to return a list of Solution objects.

  • save_path (Optional[PathLike]) – Path to an HDF5 file in which to save the results.

  • log_level (Optional[int]) – Logging level to use, if any.

  • progress_bar (bool) – Show a progress bar for self-consistent iterations.

  • _solver (str) – Name of the solver method used.

Return type:

List[Solution]

Returns:

If return_solutions is True, returns a list of Solutions of length iterations + 1.

class superscreen.Vortex(x, y, film, nPhi0=1)[source]

A vortex located at position (x, y) in film containing a total flux nPhi0 in units of the flux quantum \(\Phi_0\).

Parameters:
  • x (float) – Vortex x-position.

  • y (float) – Vortex y-position.

  • film (str) – The name of the film in which the vortex is pinned.

  • nPhi0 (float) – The number of flux quanta contained in the vortex.

Pre-factorizing models

superscreen.factorize_model(*, device, current_units, terminal_currents=None, circulating_currents=None, vortices=None)[source]

Prepares the applied field-indepdendent portions of the model, including factorizing the linear system describing each film and hole.

Parameters:
  • device (Device) – The superscreen.Device to simulate.

  • current_units (str) – Units to use for current quantities. The applied field will be converted to units of [current_units / device.length_units].

  • terminal_currents (Optional[Dict[str, Dict[str, Union[float, str, Quantity]]]]) – A dict like {film_name: {source_name: source_current}} for each film and terminal.

  • circulating_currents (Optional[Dict[str, Union[float, str, Quantity]]]) – A dict of {hole_name: circulating_current}. If circulating_current is a float, then it is assumed to be in units of current_units. If circulating_current is a string, then it is converted to a pint.Quantity.

  • vortices (Optional[Sequence[Vortex]]) – A list of superscreen.Vortex objects located in the superscreen.Device.

Return type:

FactorizedModel

Returns:

A superscreen.FactorizedModel instance that can be used to solve the model.

class superscreen.FactorizedModel(device, film_info, film_systems, hole_systems, terminal_systems, terminal_currents, circulating_currents, vortices, current_units)[source]

A pre-factorized SuperScreen model.

Parameters:
to_hdf5(h5group)[source]

Save a superscreen.FactorizedModel to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group in which to save the model

Return type:

None

static from_hdf5(h5group)[source]

Load a superscreen.FactorizedModel from an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group from which to load the model

Return type:

FactorizedModel

Returns:

The loaded superscreen.FactorizedModel

set_circulating_currents(circulating_currents)[source]

Set the circulating currents for the model.

Parameters:

circulating_currents (Dict[str, float]) – A dict of {hole_name: current}, where current is a float in units of self.current_units.

Return type:

None

Solution

class superscreen.Solution(*, device, film_solutions, applied_field_func, field_units, current_units, circulating_currents=None, terminal_currents=None, vortices=None, solver='superscreen.solve')[source]

A container for the calculated stream functions and fields, with some convenient data processing methods.

Parameters:
  • device (Device) – The Device that was solved

  • film_solutions (Dict[str, FilmSolution]) – A dict of {film_name: film_solution} containing the raw simulation results in field_units, current_units, and device.length_units.

  • applied_field_func (Callable) – The function defining the applied field

  • field_units (str) – Units of the applied field

  • current_units (str) – Units used for current quantities.

  • circulating_currents (Optional[Dict[str, float]]) – A dict of {hole_name: circulating_current}.

  • terminal_currents (Optional[Dict[str, float]]) – A dict of {terminal_name: terminal_current}.

  • vortices (Optional[List[Vortex]]) – A list of Vortex objects located in the Device.

  • solver (str) – The solver method that generated the solution.

property field_units: str

The units in which magnetic fields are specified.

property current_units: str

The units in which currents are specified.

property solver: str

The solver method that generated the solution.

property time_created: datetime

The time at which the solution was originally created.

property version_info: Dict[str, str]

A dictionary of dependency versions.

interp_current_density(positions, *, film, method='linear', units=None, with_units=False)[source]

Interpolates the current density J = [dg/dy, -dg/dx] within a film.

Parameters:
  • positions (ndarray) – Shape (m, 2) array of x, y coordinates at which to evaluate the current density.

  • film (str) – The name of the film in which to interpolate current density.

  • method (Literal['linear', 'cubic']) – Interpolation method to use (“linear” or “cubic”).

  • units (Optional[str]) – The desired units for the current density. Defaults to self.current_units / self.device.length_units.

  • with_units (bool) – Whether to return arrays of pint.Quantities with units attached.

Return type:

ndarray

Returns:

The interpolated current density

current_through_path(path_coords, *, film, interp_method='linear', units=None, with_units=True)[source]

Calculates the total current crossing a given path.

Parameters:
  • path_coords (ndarray) – An (n, 2) array of (x, y) coordinates defining the path.

  • film (str) – The name of the film in which to interpolate current density.

  • interp_method (str) – Interpolation method to use (“linear” or “cubic”).

  • units (Optional[str]) – The current units to return.

  • with_units (bool) – Whether to return a pint.Quantity with units attached.

Return type:

Union[float, Quantity]

Returns:

The total current crossing the path as either a float or a pint.Quantity.

interp_field(positions, *, film, dataset='field', method='linear', units=None, with_units=False)[source]

Interpolates the z-component of the field within a film.

Parameters:
  • positions (ndarray) – Shape (m, 2) array of x, y coordinates at which to evaluate the fields.

  • film (str) – The name of the film in which to interpolate the field.

  • dataset (Literal['field', 'self_field', 'applied_field', 'field_from_other_films']) – The dataset to interpolate. One of ‘field’, ‘self_field’, ‘applied_field’, or ‘field_from_other_films’.

  • method (Literal['linear', 'cubic']) – Interpolation method to use: ‘linear’ or ‘cubic’.

  • units (Optional[str]) – The desired units for the current density. Defaults to self.field_units.

  • with_units (bool) – Whether to return arrays of pint.Quantities with units attached.

Returns:

The interpolated field

polygon_flux(name, units=None, with_units=True)[source]

Computes the flux through a given polygon.

Parameters:
  • name (str) – The name of the polygon for which to compute the flux.

  • units (Union[str, Unit, None]) – The flux units to use.

  • with_units (bool) – Whether to a return a pint.Quantity with units attached.

Return type:

Dict[str, Union[float, Quantity]]

Returns:

The polygon flux.

polygon_fluxoid(polygon_coords, *, film, interp_method='linear', units='Phi_0', with_units=True)[source]

Computes the Fluxoid (flux + supercurrent) for a given polygonal region.

The fluxoid for a closed region \(S\) with boundary \(\partial S\) is defined as:

\[\Phi^f_S = \underbrace{ \int_S \mu_0 H_z(\vec{r})\,\mathrm{d}^2r }_{\text{flux part}} + \underbrace{ \oint_{\partial S} \mu_0\Lambda(\vec{r})\vec{J}(\vec{r})\cdot\mathrm{d}\vec{r} }_{\text{supercurrent part}}\]
Parameters:
  • polygon_coords (Union[ndarray, Polygon]) – A shape (n, 2) array of (x, y) coordinates of polygon vertices defining the closed region \(S\), or a Polygon instance.

  • film (str) – The name of the film in which to evaluate the field and current.

  • interp_method (Literal['linear', 'cubic']) – Interpolation method to use.

  • units (Optional[str]) – The desired units for the current density. Defaults to \(\Phi_0\).

  • with_units (bool) – Whether to return values as pint.Quantities with units attached.

Return type:

Dict[str, Fluxoid]

Returns:

The Fluxoid for the given polygon.

hole_fluxoid(hole_name, points=None, interp_method='linear', units='Phi_0', with_units=True)[source]

Calculcates the fluxoid for a polygon enclosing the specified hole.

Parameters:
  • hole_name (str) – The name of the hole for which to calculate the fluxoid.

  • points (Optional[ndarray]) – The vertices of the polygon enclosing the hole. If None is given, a polygon is generated using supercreen.fluxoid.make_fluxoid_polygons().

  • interp_method (Literal['linear', 'cubic']) – Interpolation method to use.

  • units (Optional[str]) – The desired units for the current density. Defaults to \(\Phi_0\).

  • with_units (bool) – Whether to return values as pint.Quantities with units attached.

Return type:

Fluxoid

Returns:

The hole’s Fluxoid.

screening_field_at_position(positions, *, zs=None, vector=False, interp_method='linear', units=None, with_units=True, return_sum=True)[source]

Calculates the field due to currents in the device at any point(s) in space (excluding the applied field).

Parameters:
  • positions (ndarray) – Shape (m, 2) array of (x, y) coordinates, or (m, 3) array of (x, y, z) coordinates at which to calculate the magnetic field. A single sequence like [x, y] or [x, y, z] is also allowed.

  • zs (Union[float, ndarray, None]) – z coordinates at which to calculate the field. If positions has shape (m, 3), then this argument is not allowed. If zs is a scalar, then the fields are calculated in a plane parallel to the x-y plane. If zs is any array, then it must be same length as positions.

  • vector (bool) – Whether to return the full vector magnetic field or just the z component.

  • interp_method (Literal['linear', 'cubic']) – Interpolation method to use.

  • units (Optional[str]) – Units to which to convert the fields (can be either magnetic field H or magnetic flux density B = mu0 * H). If not given, then the fields are returned in units of self.field_units.

  • with_units (bool) – Whether to return the fields as pint.Quantity with units attached.

  • return_sum (bool) – Whether to return the sum of the fields from all layers in the device, or a dict of {layer_name: field_from_layer}.

Return type:

Union[ndarray, Dict[str, ndarray]]

Returns:

An np.ndarray if return_sum is True, otherwise a dict of {film_name: field_from_film}. If with_units is True, then the array(s) will contain pint.Quantities. field_from_film will have shape (m, ) if vector is False, or shape (m, 3) if vector is True.

field_at_position(positions, *, zs=None, interp_method='linear', units=None, with_units=True, return_sum=True)[source]

Calculates the field due to currents in the device at any point(s) in space.

Parameters:
  • positions (ndarray) – Shape (m, 2) array of (x, y) coordinates, or (m, 3) array of (x, y, z) coordinates at which to calculate the magnetic field. A single sequence like [x, y] or [x, y, z] is also allowed.

  • zs (Union[float, ndarray, None]) – z coordinates at which to calculate the field. If positions has shape (m, 3), then this argument is not allowed. If zs is a scalar, then the fields are calculated in a plane parallel to the x-y plane. If zs is any array, then it must be same length as positions.

  • interp_method (Literal['linear', 'cubic']) – Interpolation method to use.

  • units (Optional[str]) – Units to which to convert the fields (can be either magnetic field H or magnetic flux density B = mu0 * H). If not given, then the fields are returned in units of self.field_units.

  • with_units (bool) – Whether to return the fields as pint.Quantity with units attached.

  • return_sum (bool) – Whether to return the sum of the fields from all layers in the device, or a dict of {layer_name: field_from_layer}.

Return type:

Union[ndarray, Dict[str, ndarray]]

Returns:

An np.ndarray if return_sum is True, otherwise a dict of {film_name: field_from_film}. If with_units is True, then the array(s) will contain pint.Quantities.

vector_potential_at_position(positions, *, zs=None, units=None, with_units=True, return_sum=True)[source]

Calculates the vector potential due to currents in the device at any point(s) in space. Note that this only considers the vector potential due to currents in the device, so only represents the total vector potential in cases where the applied field is zero (e.g. models with only vortices and/or circulating currents).

The vector potential \(\vec{A}\) at position \(\vec{r}\) due to sheet current density \(\vec{J}(\vec{r}')\) flowing in a film with lateral geometry \(S\) is:

\[\vec{A}(\vec{r}) = \frac{\mu_0}{4\pi} \int_S\frac{\vec{J}(\vec{r}')}{|\vec{r}-\vec{r}'|}\mathrm{d}^2r'.\]
Parameters:
  • positions (ndarray) – Shape (m, 2) array of (x, y) coordinates, or (m, 3) array of (x, y, z) coordinates at which to calculate the vector potential. A single list like [x, y] or [x, y, z] is also allowed.

  • zs (Union[float, ndarray, None]) – z coordinates at which to calculate the potential. If positions has shape (m, 3), then this argument is not allowed. If zs is a scalar, then the fields are calculated in a plane parallel to the x-y plane. If zs is any array, then it must be same length as positions.

  • units (Optional[str]) – Units to which to convert the vector potential.

  • with_units (bool) – Whether to return the vector potential as a pint.Quantity with units attached.

  • return_sum (bool) – Whether to return the sum of the potential from all layers in the device, or a dict of {layer_name: potential_from_layer}.

Return type:

Union[ndarray, Dict[str, ndarray]]

Returns:

An np.ndarray if return_sum is True, otherwise a dict of {film_name: potential_from_film}. If with_units is True, then the array(s) will contain pint.Quantities. potential_from_film will have shape (m, 3).

to_hdf5(path_or_group, device_path=None, compress=True)[source]

Save the Solution to an HDF5 file.

Return type:

None

Parameters:
  • path_or_group (PathLike | Group) –

  • device_path (str | None) –

  • compress (bool) –

path_or_group: An HDF5 file path or an open h5py.Group in which to save

the Solution.

device_path: Path within the HDF5 file in which the Solution’s Device

is saved. If None, the Device will be saved at "/device".

compress: Save the mesh in a compressed format.

static from_hdf5(path_or_group)[source]

Load a Solution from and HDF5 file.

Parameters:

path_or_group (Union[PathLike, Group]) – An HDF5 file path or an open h5py.Group from which to load the Solution.

Return type:

Solution

Returns:

The loaded Solution

static save_solutions(solutions, path_or_group, compress=True)[source]

Save a series of Solutions to an HDF5 file.

Parameters:
  • solutions (Sequence[Solution]) – A series of Solutions to save.

  • path_or_group (Union[PathLike, Group]) – An HDF5 file path or an open h5py.Group in which to save the Solutions.

  • compress (bool) – Save the meshes in a compressed format.

Return type:

None

static load_solutions(path_or_group)[source]

Load a series of Solutions from an HDF5 file.

Parameters:

path_or_group (Union[PathLike, Group]) – An HDF5 file path or an open h5py.Group from which to load the Solutions.

Return type:

List[Solution]

Returns:

A list of loaded Solutions.

equals(other, require_same_timestamp=False)[source]

Checks whether two solutions are equal.

Parameters:
  • other (Any) – The Solution to compare for equality.

  • require_same_timestamp (bool) – If True, two solutions are only considered equal if they have the exact same time_created.

Return type:

bool

Returns:

A boolean indicating whether the two solutions are equal

plot_streams(**kwargs)[source]

Alias for superscreen.visualization.plot_streams().

Plots the stream function for multiple layers in a Device.

Additional keyword arguments are passed to plt.subplots().

Parameters:
  • films – Name(s) of films(s) for which to plot the stream function. By default, the stream function is plotted for all films in the Device.

  • units – Units in which to plot the stream function. Defaults to solution.current_units.

  • max_cols – Maximum number of columns in the grid of subplots.

  • cmap – Name of the matplotlib colormap to use.

  • colorbar – Whether to add a colorbar to each subplot.

  • shading – May be "flat" or "gouraud". The latter does some interpolation.

  • auto_range_cutoff – Cutoff percentile for auto_range_iqr.

  • share_color_scale – Whether to force all layers to use the same color scale.

  • symmetric_color_scale – Whether to use a symmetric color scale (vmin = -vmax).

  • vmin – Color scale minimum to use for all layers

  • vmax – Color scale maximum to use for all layers

Return type:

Tuple[Figure, ndarray]

Returns:

matplotlib figure and axes

plot_currents(**kwargs)[source]

Alias for superscreen.visualization.plot_currents().

Plots the sheet current density for one or more films in a Device.

Additional keyword arguments are passed to plt.subplots().

Parameters:
  • films – Name(s) of film(s) for which to plot the sheet current. By default, the sheet current is plotted for all films in the Device.

  • units – Units in which to plot the current density. Defaults to solution.current_units / solution.device.length_units.

  • max_cols – Maximum number of columns in the grid of subplots.

  • cmap – Name of the matplotlib colormap to use.

  • colorbar – Whether to add a colorbar to each subplot.

  • shading – May be "flat" or "gouraud". The latter does some interpolation.

  • auto_range_cutoff – Cutoff percentile for auto_range_iqr.

  • share_color_scale – Whether to force all layers to use the same color scale.

  • symmetric_color_scale – Whether to use a symmetric color scale (vmin = -vmax).

  • vmin – Color scale minimum to use for all layers (ignored if share_color_scale is True).

  • vmax – Color scale maximum to use for all layers (ignored if share_color_scale is True).

  • streamplot – Whether to overlay current streamlines on the plot.

  • grid_shape – Shape of the rectangular grid used to contruct streamlines. If a single integer N is given, the grid will be square, shape = (N, N).

  • min_stream_amp – Streamlines will not be drawn anywhere the current density is less than min_stream_amp * max(current_density). This avoids streamlines being drawn where there is no current flowing.

  • cross_section_coords – Shape (m, 2) array of (x, y) coordinates for a cross-section (or a list of such arrays).

Return type:

Tuple[Figure, ndarray]

Returns:

matplotlib figure and axes

plot_fields(**kwargs)[source]

Alias for superscreen.visualization.plot_fields().

Plots the fields for one or more films in a Device.

Additional keyword arguments are passed to plt.subplots().

Parameters:
  • films – Name(s) of films(s) for which to plot the fields. By default, the field is plotted for all films in the Device.

  • dataset – Which set of fields to plot, either “field”, “self_field”, “applied_field”, or “field_from_other_films”.

  • normalize – Whether to normalize the fields by the applied field.

  • units – Units in which to plot the fields. Defaults to solution.field_units. This argument is ignored if normalize is True.

  • shading – May be "flat" or "gouraud". The latter does some interpolation.

  • max_cols – Maximum number of columns in the grid of subplots.

  • cmap – Name of the matplotlib colormap to use.

  • colorbar – Whether to add a colorbar to each subplot.

  • auto_range_cutoff – Cutoff percentile for auto_range_iqr.

  • share_color_scale – Whether to force all layers to use the same color scale.

  • symmetric_color_scale – Whether to use a symmetric color scale (vmin = -vmax).

  • vmin – Color scale minimum to use for all layers

  • vmax – Color scale maximum to use for all layers

  • cross_section_coords – Shape (m, 2) array of (x, y) coordinates for a cross-section (or a list of such arrays).

Return type:

Tuple[Figure, ndarray]

Returns:

matplotlib figure and axes

plot_field_at_positions(points, **kwargs)[source]

Alias for superscreen.visualization.plot_field_at_positions().

Plots the total field (either all three components or just the z component) anywhere in space.

Additional keyword arguments are passed to plt.subplots().

Parameters:
  • positions – Shape (m, 2) array of (x, y) coordinates, or (m, 3) array of (x, y, z) coordinates at which to calculate the magnetic field. A single list like [x, y] or [x, y, z] is also allowed.

  • zs – z coordinates at which to calculate the field. If positions has shape (m, 3), then this argument is not allowed. If zs is a scalar, then the fields are calculated in a plane parallel to the x-y plane. If zs is any array, then it must be same length as positions.

  • units – Units in which to plot the fields. Defaults to solution.field_units. This argument is ignored if normalize is True.

  • shading – May be "flat" or "gouraud". The latter does some interpolation.

  • max_cols – Maximum number of columns in the grid of subplots.

  • cmap – Name of the matplotlib colormap to use.

  • colorbar – Whether to add a colorbar to each subplot.

  • auto_range_cutoff – Cutoff percentile for auto_range_iqr.

  • share_color_scale – Whether to force all layers to use the same color scale.

  • symmetric_color_scale – Whether to use a symmetric color scale (vmin = -vmax).

  • vmin – Color scale minimum to use for all layers

  • vmax – Color scale maximum to use for all layers

  • cross_section_coords – Shape (m, 2) array of (x, y) coordinates for a cross-section (or a list of such arrays).

  • points (ndarray) –

Return type:

Tuple[Figure, Axes]

Returns:

matplotlib figure and axes

class superscreen.FilmSolution(stream, current_density, applied_field, self_field, field_from_other_films=None)[source]

Raw solution data for a single film.

Parameters:
  • stream (ndarray) – The stream function

  • current_density (ndarray) – The sheet current density

  • applied_field (ndarray) – The applied field

  • self_field (ndarray) – The field due to screening currents in the film

  • field_from_other_films (Optional[ndarray]) – The field due to screening currents in other films

property total_field: ndarray

The total magnetic field in the film.

to_hdf5(h5group)[source]

Save the superscreen.FilmSolution to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group in which to save the superscreen.FilmSolution

Return type:

None

static from_hdf5(h5group)[source]

Load a superscreen.FilmSolution from an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group from which to load the superscreen.FilmSolution

Return type:

FilmSolution

Returns:

The loaded superscreen.FilmSolution

is_close(other, rtol=0.0001, atol=1e-07)[source]

Check whether two FilmSolutions are equal to within a tolerance.

Parameters:
Return type:

bool

Returns:

True if the two FilmSolutions are equal within the given tolerances.

Fluxoid

class superscreen.Fluxoid(flux_part, supercurrent_part)[source]

Bases: NamedTuple

The fluxoid for a closed region \(S\) with boundary \(\partial S\) is defined as:

\[\Phi^f_S = \underbrace{ \int_S \mu_0 H_z(\vec{r})\,\mathrm{d}^2r }_{\text{flux part}} + \underbrace{ \oint_{\partial S} \mu_0\Lambda(\vec{r})\vec{J}(\vec{r})\cdot\mathrm{d}\vec{r} }_{\text{supercurrent part}}\]
Parameters:
  • flux_part (float | Quantity) – \(\int_S \mu_0 H_z(\vec{r})\,\mathrm{d}^2r\).

  • supercurrent_part (float | Quantity) – \(\oint_{\partial S}\mu_0\Lambda(\vec{r})\vec{J}(\vec{r})\cdot\mathrm{d}\vec{r}\).

superscreen.fluxoid.make_fluxoid_polygons(device, holes=None, join_style='mitre', interp_points=None)[source]

Generates polygons enclosing the given holes to calculate the fluxoid.

Parameters:
  • device (Device) – The Device for which to generate polygons.

  • holes (Union[str, List[str], None]) – Name(s) of the hole(s) in the device for which to generate polygons. Defaults to all holes in the device.

  • join_style (str) – See superscreen.components.Polygon.buffer().

  • interp_points (Optional[int]) – If provided, the resulting polygons will be interpolated to interp_points vertices.

Return type:

Dict[str, ndarray]

Returns:

A dict of {hole_name: fluxoid_polygon}.

superscreen.fluxoid.find_fluxoid_solution(device, *, fluxoids=None, **solve_kwargs)[source]

Calculates the current(s) circulating around hole(s) in the device required to realize the specified fluxoid state.

Parameters:
  • device (Device) – The Device for which to find the given fluxoid solution.

  • fluxoids (Optional[Dict[str, float]]) – A dict of {hole_name: fluxoid_value}, where fluxoid_value is in units of \(\Phi_0\). The fluxoid for any holes not in this dict will default to 0.

  • solve_kwargs – Additional keyword arguments are passed to superscreen.solve.solve().

Return type:

Solution

Returns:

The optimized superscreen.Solution.

Supporting functions

superscreen.solver.solve_film(*, device, applied_field, film_info, film_system, hole_systems, field_conversion, vortex_flux, terminal_systems=None, field_from_other_films=None, check_inversion=False)[source]

Computes the stream function and magnetic field within a single film in a superscreen.Device.

Parameters:
Return type:

FilmSolution

Returns:

A superscreen.FilmSolution containing the results

superscreen.solver.solve_for_terminal_current_stream(device, film_info, terminal_systems, terminal_currents)[source]

Solves for the stream function associated with transport currents in a single film.

The algorithm is:

1. Solve for the stream function in the film assuming no applied field and ignoring the presence of any holes. 2. Set the stream function in each hole to the weighted average of the stream function found in step 1. 3. Re-solve for the stream function in the film with the new hole boundary conditions.

Parameters:
Return type:

ndarray

Returns:

The stream function associated with the transport current

superscreen.solver.factorize_linear_systems(device, film_info_dict)[source]

Build and factorize the linear systems for all films, holes, and terminals.

Parameters:
Return type:

Tuple[Dict[str, LinearSystem], Dict[str, Dict[str, LinearSystem]], Dict[str, TerminalSystems]]

Returns:

A dict of {film_name: film_system}, a dict of {film_name: {hole_name: hole_system}}, and a dict of {film_name: TerminalSystems}

superscreen.solver.convert_field(value, new_units, old_units=None, ureg=None, with_units=True)[source]

Converts a value between different field units, either magnetic field H [current] / [length] or flux density B = mu0 * H [mass] / ([curret] [time]^2)).

Parameters:
  • value (Union[ndarray, float, str, Quantity]) – The value to convert. It can either be a numpy array (no units), a float (no units), a string like “1 uA/um”, or a scalar or array pint.Quantity. If value is not a string with units or a pint.Quantity, then old_units must specify the units of the float or array.

  • new_units (Union[str, Unit]) – The units to convert to.

  • old_units (Union[str, Unit, None]) – The old units of value. This argument is required if value is not a string with units or a pint.Quantity.

  • ureg (Optional[UnitRegistry]) – The pint.UnitRegistry to use for conversion. If None is given, a new instance is created.

  • with_units (bool) – Whether to return a pint.Quantity with units attached.

Return type:

Union[Quantity, ndarray, float]

Returns:

The converted value, either a pint.Quantity (scalar or array with units), or an array or float without units, depending on the with_units argument.

superscreen.solver.field_conversion_factor(field_units, current_units, length_units='m', ureg=None)[source]

Returns a conversion factor from field_units to current_units / length_units.

Parameters:
  • field_units (str) – Magnetic field/flux unit to convert, having dimensionality either of magnetic field H (e.g. A / m or Oe) or of magnetic flux density B = mu0 * H (e.g. Tesla or Gauss).

  • current_units (str) – Current unit to use for the conversion.

  • length_units (str) – Lenght/distance unit to use for the conversion.

  • ureg (Optional[UnitRegistry]) – pint UnitRegistry to use for the conversion. If None is provided, a new UnitRegistry is created.

Return type:

Quantity

Returns:

Conversion factor as a pint.Quantity. conversion_factor.magnitude gives you the numerical value of the conversion factor.

Supporting classes

class superscreen.solver.FilmInfo(name, layer, lambda_info, vortices, interior_indices, boundary_indices, hole_indices, in_hole, circulating_currents, weights, kernel, laplacian, gradient=None, terminal_currents=None)[source]

A container for information about a single film required for the solver.

Parameters:
  • name (str) – The film name

  • layer (str) – The layer in which the film exists

  • lambda_info (LambdaInfo) – A superscreen.solver.LambdaInfo instance defining the effective penetration depth in the film

  • vortices (Tuple[Vortex]) – Any superscreen.Vortex instances located in the film

  • interior_indices (ndarray) – The indices of the film in its mesh

  • boundary_indices (ndarray) – Indices of the boundary vertices for the mesh, ordered counterclockwise

  • hole_indices (Dict[str, ndarray]) – A dict containing the indices of each hole in the film’s mesh

  • in_hole (ndarray) – A boolean array indicated which mesh sites lie inside a hole

  • circulating_currents (Dict[str, float]) – A dict of {hole_name, circulating_current}

  • weights (ndarray) – The mesh weights

  • kernel (ndarray) – The mesh self-field kernel \(\mathbf{Q}\)

  • laplacian (ndarray) – The mesh Laplacian \(\nabla^2\)

  • gradient (Optional[ndarray]) – The mesh gradient operator \(\vec{\nabla}\)

  • terminal_currents (Optional[Dict[str, float]]) – A dict of {terminal_name: terminal_current}

to_hdf5(h5group)[source]

Save the superscreen.solver.FilmInfo instance to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group in which to save the FilmInfo

Return type:

None

static from_hdf5(h5group)[source]

Load a superscreen.solver.FilmInfo instance from an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group from which to load the FilmInfo

Return type:

FilmInfo

Returns:

The loaded superscreen.solver.FilmInfo instance

class superscreen.solver.LambdaInfo(*, film, Lambda, london_lambda=None, thickness=None)[source]

An object containing information about the effective penetration depth for a given film.

Parameters:
  • film (str) – The name of the film

  • Lambda (ndarray) – The effective penetration depth at each mesh site

  • london_lambda (Optional[ndarray]) – The London penetration depth at each mesh site

  • thickness (Optional[float]) – The thickness of the layer in which the film exists

to_hdf5(h5group)[source]

Save the LambdaInfo instance to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group in which to save the LambdaInfo instance.

Return type:

None

static from_hdf5(h5group)[source]

Load a LambdaInfo instance from an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group from which to load the LambdaInfo instance.

Return type:

LambdaInfo

Returns:

The loaded LambdaInfo instance

class superscreen.solver.LinearSystem(A, indices, lu_piv=None, grad_Lambda_term=0.0)[source]

The linear system representing a given film or hole.

Parameters:
  • A (ndarray) – The matrix quantity to be inverted, \(\mathbf{Q}.\mathbf{w}^T-\mathbf{\Lambda}^T.\mathbf{\nabla}^2-\vec{\nabla}\mathbf{\Lambda}\cdot\vec{\nabla}\)

  • indices (ndarray) – The indices into the corresponding mesh

  • lu_piv (Optional[Tuple[ndarray, ndarray]]) – The LU factorization lu_piv = scipy.linalg.lu_factor(-A), see scipy.linalg.lu_factor()

  • grad_Lambda_term (Union[float, ndarray]) – The term corresponding to the gradient of the effective penetration depth.

to_hdf5(h5group)[source]

Save a superscreen.solver.LinearSystem to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group to which to save the LinearSystem

Return type:

None

static from_hdf5(h5group)[source]

Load a superscreen.solver.LinearSystem to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group from which to load the LinearSystem

Return type:

LinearSystem

Returns:

The loaded superscreen.solver.LinearSystem

class superscreen.solver.TerminalSystems(film, boundary, holes, film_without_boundary, film_without_boundary_or_holes=None)[source]

A container for the superscreen.solver.LinearSystem objects needed to calculate the stream function associated with transport currents.

Parameters:
to_hdf5(h5group)[source]

Save a superscreen.solver.TerminalSystems to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group to which to save the TerminalSystems

Return type:

None

static from_hdf5(h5group)[source]

Load a superscreen.solver.TerminalSystems to an h5py.Group.

Parameters:

h5group (Group) – The h5py.Group from which to load the TerminalSystems

Return type:

TerminalSystems

Returns:

The loaded superscreen.solver.TerminalSystems